https://www.cdslab.org/paramonte/fortran/2
Current view: top level - main - pm_complexMinMax@routines.inc.F90 (source / functions) Hit Total Coverage
Test: ParaMonte 2.0.0 :: Serial Fortran - Code Coverage Report Lines: 39 39 100.0 %
Date: 2024-04-08 03:18:57 Functions: 0 0 -
Legend: Lines: hit not hit

          Line data    Source code
       1             : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
       2             : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
       3             : !!!!                                                                                                                            !!!!
       4             : !!!!    ParaMonte: Parallel Monte Carlo and Machine Learning Library.                                                           !!!!
       5             : !!!!                                                                                                                            !!!!
       6             : !!!!    Copyright (C) 2012-present, The Computational Data Science Lab                                                          !!!!
       7             : !!!!                                                                                                                            !!!!
       8             : !!!!    This file is part of the ParaMonte library.                                                                             !!!!
       9             : !!!!                                                                                                                            !!!!
      10             : !!!!    LICENSE                                                                                                                 !!!!
      11             : !!!!                                                                                                                            !!!!
      12             : !!!!       https://github.com/cdslaborg/paramonte/blob/main/LICENSE.md                                                          !!!!
      13             : !!!!                                                                                                                            !!!!
      14             : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
      15             : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
      16             : 
      17             : !>  \brief
      18             : !>  This file contains the implementation details of the routines under the
      19             : !>  generic interfaces of module [pm_complexAbs](@ref pm_complexAbs).
      20             : !>
      21             : !>  \finmain
      22             : !>
      23             : !>  \author
      24             : !>  \FatemehBagheri, Wednesday 12:20 AM, October 13, 2021, Dallas, TX
      25             : 
      26             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      27             : 
      28             :         ! Define the initial values and comparison operations.
      29             : #if     minval_ENABLED || minloc_ENABLED
      30             : #define INIT+huge(0._CKC)
      31             : #define COMPARES_WITH>
      32             : #define GETLOC minloc
      33             : #define GETVAL minval
      34             : #elif   maxval_ENABLED || maxloc_ENABLED
      35             : #define INIT-huge(0._CKC)
      36             : #define COMPARES_WITH<
      37             : #define GETLOC maxloc
      38             : #define GETVAL maxval
      39             : #endif
      40             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      41             : #if     min_ENABLED && D0_ENABLED && CK_ENABLED
      42             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      43             : 
      44         125 :         val%re = min(a1%re, a2%re)
      45         125 :         val%im = min(a1%im, a2%im)
      46             : 
      47             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      48             : #elif   max_ENABLED && D0_ENABLED && CK_ENABLED
      49             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      50             : 
      51         121 :         val%re = max(a1%re, a2%re)
      52         121 :         val%im = max(a1%im, a2%im)
      53             : 
      54             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      55             : #elif   (minval_ENABLED || maxval_ENABLED) && D1_ENABLED && CK_ENABLED
      56             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      57             : 
      58             :         integer(IK) :: iell
      59       61020 :         val%re = INIT
      60       61020 :         val%im = INIT
      61      473246 :         do iell = 1, size(array, 1, IK)
      62      392106 :             if (val%re COMPARES_WITH array(iell)%re) val%re = array(iell)%re
      63      473246 :             if (val%im COMPARES_WITH array(iell)%im) val%im = array(iell)%im
      64             :         end do
      65             : 
      66             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      67             : #elif   (minval_ENABLED || maxval_ENABLED) && D2_ENABLED && ALL_ENABLED && CK_ENABLED
      68             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      69             : 
      70             :         complex(CKC), pointer :: parray(:)
      71        2410 :         parray(1 : size(array, kind = IK)) => array
      72        2410 :         val = GETVAL(parray)
      73        2410 :         nullify(parray)
      74             : 
      75             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      76             : #elif   (minval_ENABLED || maxval_ENABLED) && D2_ENABLED && DIM_ENABLED && CK_ENABLED
      77             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      78             : 
      79             :         integer(IK) :: idim, ndim, nsam
      80       64932 :         CHECK_ASSERTION(__LINE__, 1 <= dim .and. dim <= rank(array), SK_": The condition `1 <= dim .and. dim <= rank(array)` must hold. dim, rank(array) = "//getStr([dim, int(rank(array), IK)]))
      81       21644 :         ndim = size(array, 3 - dim, IK)
      82       21644 :         nsam = size(array, dim, IK)
      83       21644 :         if (dim == 1_IK) then
      84             :             do concurrent(idim = 1 : ndim)
      85       29208 :                 val(idim) = GETVAL(array(1 : nsam, idim))
      86             :             end do
      87       14110 :         elseif (dim == 2_IK) then
      88             :             do concurrent(idim = 1 : ndim)
      89       56720 :                 val(idim) = GETVAL(array(idim, 1 : nsam))
      90             :             end do
      91             :         end if
      92             : 
      93             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      94             : #elif   (minloc_ENABLED || maxloc_ENABLED) && D1_ENABLED && CK_ENABLED
      95             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      96             : 
      97             :         complex(CKC) :: val
      98             :         integer(IK) :: iell
      99         216 :         loc = 0_IK
     100          48 :         val%re = INIT
     101          48 :         val%im = INIT
     102         342 :         do iell = 1, size(array, 1, IK)
     103         270 :             if (val%re COMPARES_WITH array(iell)%re) then
     104         109 :                 val%re = array(iell)%re
     105         109 :                 loc(1) = iell
     106             :             end if
     107         342 :             if (val%im COMPARES_WITH array(iell)%im) then
     108         126 :                 val%im = array(iell)%im
     109         126 :                 loc(2) = iell
     110             :             end if
     111             :         end do
     112             : 
     113             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     114             : #elif   (minloc_ENABLED || maxloc_ENABLED) && D2_ENABLED && ALL_ENABLED && CK_ENABLED
     115             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     116             : 
     117             :         complex(CKC), pointer :: parray(:)
     118           6 :         parray(1 : size(array, kind = IK)) => array
     119           6 :         loc = GETLOC(parray)
     120           6 :         nullify(parray)
     121             : 
     122             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     123             : #elif   (minloc_ENABLED || maxloc_ENABLED) && D2_ENABLED && DIM_ENABLED && CK_ENABLED
     124             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     125             : 
     126             :         integer(IK) :: idim, ndim, nsam
     127          36 :         CHECK_ASSERTION(__LINE__, 1 <= dim .and. dim <= rank(array), SK_": The condition `1 <= dim .and. dim <= rank(array)` must hold. dim, rank(array) = "//getStr([dim, int(rank(array), IK)]))
     128          12 :         ndim = size(array, 3 - dim, IK)
     129          12 :         nsam = size(array, dim, IK)
     130          12 :         if (dim == 1_IK) then
     131             :             do concurrent(idim = 1 : ndim)
     132          36 :                 loc(1 : 2, idim) = GETLOC(array(1 : nsam, idim))
     133             :             end do
     134           6 :         elseif (dim == 2_IK) then
     135             :             do concurrent(idim = 1 : ndim)
     136          24 :                 loc(1 : 2, idim) = GETLOC(array(idim, 1 : nsam))
     137             :             end do
     138             :         end if
     139             : 
     140             : #else
     141             :         !%%%%%%%%%%%%%%%%%%%%%%%%
     142             : #error  "Unrecognized interface."
     143             :         !%%%%%%%%%%%%%%%%%%%%%%%%
     144             : #endif
     145             : #undef  COMPARES_WITH
     146             : #undef  GETLOC
     147             : #undef  GETVAL
     148             : #undef  INIT

ParaMonte: Parallel Monte Carlo and Machine Learning Library 
The Computational Data Science Lab
© Copyright 2012 - 2024