https://www.cdslab.org/paramonte/fortran/2
Current view: top level - test - test_pm_arrayCopy@routines.inc.F90 (source / functions) Hit Total Coverage
Test: ParaMonte 2.0.0 :: Serial Fortran - Code Coverage Report Lines: 67 67 100.0 %
Date: 2024-04-08 03:18:57 Functions: 80 80 100.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 procedure implementations of [test_pm_arrayCopy](@ref test_pm_arrayCopy).
      19             : !>
      20             : !>  \fintest
      21             : !>
      22             : !>  \author
      23             : !>  \AmirShahmoradi, September 1, 2017, 11:35 PM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin
      24             : 
      25             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      26             : 
      27             : #if     LK_ENABLED
      28             : #define IS_EQUAL .eqv.
      29             : #else
      30             : #define IS_EQUAL ==
      31             : #endif
      32             : 
      33             : #if     SK_ENABLED && D0_ENABLED
      34             : #define GET_SIZE(x) len(x, kind = IK)
      35             : #define GET_INDEX(i) i:i
      36             : #else
      37             : #define GET_SIZE(x) size(x, kind = IK)
      38             : #define GET_INDEX(i) i
      39             : #endif
      40             : 
      41             : #if     SK_ENABLED && D0_ENABLED
      42             : #define ALL
      43             :         use pm_str, only: getCharSeq, getCharVec
      44           4 :         character(:,SKC), allocatable   :: From, To, To_ref
      45             :         character(1,SKC), parameter     :: low = SKC_"a"
      46             :         character(1,SKC), parameter     :: upp = SKC_"z"
      47             : #elif   SK_ENABLED && D1_ENABLED
      48             :         character(2,SKC), dimension(:), allocatable :: From, To, To_ref
      49             :         character(2,SKC), parameter                 :: low = SKC_"aa"
      50             :         character(2,SKC), parameter                 :: upp = SKC_"zz"
      51             : #elif   IK_ENABLED && D1_ENABLED
      52             :         integer(IKC)    , dimension(:), allocatable :: From, To, To_ref
      53             :         integer(IKC)    , parameter                 :: low = -huge(1_IKC)
      54             :         integer(IKC)    , parameter                 :: upp = +huge(1_IKC)
      55             : #elif   LK_ENABLED && D1_ENABLED
      56             :         logical(LKC)    , dimension(:), allocatable :: From, To, To_ref
      57             :         logical(LKC)    , parameter                 :: low = .false._LKC
      58             :         logical(LKC)    , parameter                 :: upp = .true._LKC
      59             : #elif   CK_ENABLED && D1_ENABLED
      60             :         complex(CKC)    , dimension(:), allocatable :: From, To, To_ref
      61             :         complex(CKC)    , parameter                 :: low = -cmplx(huge(0._CKC), huge(0._CKC), kind = CKC)
      62             :         complex(CKC)    , parameter                 :: upp = +cmplx(huge(0._CKC), huge(0._CKC), kind = CKC)
      63             : #elif   RK_ENABLED && D1_ENABLED
      64             :         real(RKC)       , dimension(:), allocatable :: From, To, To_ref
      65             :         real(RKC)       , parameter                 :: low = -huge(0._RKC)
      66             :         real(RKC)       , parameter                 :: upp = +huge(0._RKC)
      67             : #else
      68             : #error  "Unrecognized interface."
      69             : #endif
      70             :         integer(IK) :: itest, lenFrom, lenTo, lenCopy
      71             : 
      72             :         !%%%%%%%%%%%%%%%%%%%%%
      73             : #if     setCopyIndexed_ENABLED
      74             :         !%%%%%%%%%%%%%%%%%%%%%
      75             : 
      76             :         integer(IK), allocatable :: indexF(:), indexT(:)
      77             : 
      78          42 :         assertion = .true._LK
      79        4244 :         do itest = 1, 100
      80        4200 :             lenTo = getUnifRand(1_IK, 10_IK)
      81        4200 :             lenFrom = getUnifRand(1_IK, 10_IK)
      82        4200 :             lenCopy = getUnifRand(0_IK, max(lenFrom, lenTo))
      83       23779 :             indexF = getUnifRand(1_IK, lenFrom, lenCopy)
      84       23779 :             indexT = getUnifRand(1_IK, lenTo, lenCopy)
      85        4200 :             call setResized(From, lenFrom)
      86        4200 :             call setResized(To, lenTo)
      87       26412 :             call setUnifRand(From)
      88       26197 :             call setUnifRand(To)
      89       33830 :             To_ref = To
      90        4200 :             call setCopyIndexed_ref()
      91        4200 :             call setCopyIndexed(From, To, indexF, indexT)
      92        4242 :             call report()
      93             :         end do
      94             : 
      95             :     contains
      96             : 
      97        4200 :         subroutine setCopyIndexed_ref()
      98             :             integer(IK) :: i
      99       19579 :             do i = 1_IK, lenCopy
     100       19579 :                 To_ref(GET_INDEX(indexT(i))) = From(GET_INDEX(indexF(i)))
     101             :             end do
     102        4200 :         end subroutine
     103             : 
     104        4200 :         subroutine report()
     105             :             use pm_io, only: display_type
     106        4200 :             type(display_type) :: disp
     107       26197 :             assertion = assertion .and. logical(ALL(To IS_EQUAL To_ref), LK)
     108        4200 :             if (test%traceable .and. .not. assertion) then
     109             :                 ! LCOV_EXCL_START
     110             :                 call disp%skip()
     111             :                 call disp%show("From")
     112             :                 call disp%show( From )
     113             :                 call disp%show("indexF")
     114             :                 call disp%show( indexF )
     115             :                 call disp%show("indexT")
     116             :                 call disp%show( indexT )
     117             :                 call disp%show("To_ref")
     118             :                 call disp%show( To_ref )
     119             :                 call disp%show("To")
     120             :                 call disp%show( To )
     121             :                 ! LCOV_EXCL_STOP
     122             :             end if
     123        4200 :             call test%assert(assertion, SK_"The array must be copied correctly.", int(__LINE__, IK))
     124        4200 :         end subroutine
     125             : 
     126             :         !%%%%%%%%%%%%%%%%%%%%%
     127             : #elif   setCopyStrided_ENABLED
     128             :         !%%%%%%%%%%%%%%%%%%%%%
     129             : 
     130             :         integer(IK) :: incf, inct
     131             : 
     132          42 :         assertion = .true._LK
     133        4244 :         do itest = 1, 100
     134             :             do
     135        4232 :                 incf = getUnifRand(-5_IK, 5_IK)
     136        4232 :                 inct = getUnifRand(-5_IK, 5_IK)
     137        4232 :                 if (incf == 0_IK .and. inct == 0_IK) cycle
     138          32 :                 exit
     139             :             end do
     140        4200 :             lenCopy = getUnifRand(1_IK, 10_IK)
     141        4200 :             lenFrom = merge(lenCopy, (lenCopy - 1_IK) * abs(incf) + 1_IK, incf == 0_IK)
     142        4200 :             lenTo   = merge(lenCopy, (lenCopy - 1_IK) * abs(inct) + 1_IK, inct == 0_IK)
     143        4200 :             call setResized(From, lenFrom)
     144        4200 :             call setResized(To, lenTo)
     145       59509 :             call setUnifRand(From)
     146       59000 :             call setUnifRand(To)
     147       66819 :             To_ref = To
     148             : 
     149        4200 :             call setCopyStrided_ref()
     150        4200 :             call setCopyStrided(From, To, incf, inct)
     151        4242 :             call report()
     152             :         end do
     153             :         
     154             :     contains
     155             :         
     156        4200 :         subroutine setCopyStrided_ref()
     157             :             integer(IK) :: ifrom, ito
     158        4200 :             if (incf > 0_IK) then
     159             :                 ito = 1_IK
     160        1917 :                 if (inct < 0_IK) ito = GET_SIZE(To_ref)
     161        1917 :                 do ifrom = 1_IK, GET_SIZE(From), incf
     162       10452 :                     To_ref(GET_INDEX(ito)) = From(GET_INDEX(ifrom))
     163       10452 :                     ito = ito + inct
     164             :                 end do
     165        2283 :             elseif (incf < 0_IK) then
     166             :                 ito = 1_IK
     167        1930 :                 if (inct < 0_IK) ito = GET_SIZE(To_ref)
     168        1930 :                 do ifrom = GET_SIZE(From), 1_IK, incf
     169       10847 :                     To_ref(GET_INDEX(ito)) = From(GET_INDEX(ifrom))
     170       10847 :                     ito = ito + inct
     171             :                 end do
     172         353 :             elseif (inct > 0_IK) then
     173         173 :                 do concurrent(ito = 1_IK : GET_SIZE(To_ref) : inct)
     174        1115 :                     To_ref(GET_INDEX(ito)) = From(GET_INDEX(1_IK))
     175             :                 end do
     176         180 :             elseif (inct < 0_IK) then
     177         171 :                 do concurrent(ito = GET_SIZE(To_ref) : 1_IK : inct)
     178        1239 :                     To_ref(GET_INDEX(ito)) = From(GET_INDEX(1_IK))
     179             :                 end do
     180             :             end if
     181        4200 :         end subroutine
     182             :         
     183        4200 :         subroutine report()
     184             :             use pm_io, only: display_type
     185        4200 :             type(display_type) :: disp
     186       59000 :             assertion = assertion .and. logical(ALL(To IS_EQUAL To_ref), LK)
     187        4200 :             if (test%traceable .and. .not. assertion) then
     188             :                 ! LCOV_EXCL_START
     189             :                 call disp%skip()
     190             :                 call disp%show("From")
     191             :                 call disp%show( From )
     192             :                 call disp%show("incf")
     193             :                 call disp%show( incf )
     194             :                 call disp%show("inct")
     195             :                 call disp%show( inct )
     196             :                 call disp%show("To_ref")
     197             :                 call disp%show( To_ref )
     198             :                 call disp%show("To")
     199             :                 call disp%show( To )
     200             :                 ! LCOV_EXCL_STOP
     201             :             end if
     202        4200 :             call test%assert(assertion, SK_"The array must be copied correctly.", int(__LINE__, IK))
     203        4200 :         end subroutine
     204             : 
     205             : #else
     206             :         !%%%%%%%%%%%%%%%%%%%%%%%%
     207             : #error  "Unrecognized interface."
     208             :         !%%%%%%%%%%%%%%%%%%%%%%%%
     209             : #endif
     210             : 
     211             : #undef GET_INDEX
     212             : #undef GET_SIZE
     213             : #undef IS_EQUAL
     214             : #undef ALL

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