https://www.cdslab.org/paramonte/fortran/2
Current view: top level - main - pm_arrayInsert@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 generic interfaces in [pm_arrayInsert](@ref pm_arrayInsert).
      19             : !>
      20             : !>  \author
      21             : !>  \AmirShahmoradi, September 1, 2017, 12:00 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin
      22             : 
      23             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      24             : 
      25             :         ! Determine assumed-length scalar character vs. array input arguments.
      26             : #if     D0_D0_ENABLED && SK_ENABLED
      27             : #define GET_SIZE len
      28             : #define ALL
      29             : #elif   D1_D0_ENABLED || D1_D1_ENABLED
      30             : #define GET_SIZE size
      31             : #define ALL all
      32             : #else
      33             : #error  "Unrecognized interface."
      34             : #endif
      35             :         ! Define logical vs. normal equivalence operators
      36             : #if     LK_ENABLED
      37             : #define IS_EQUAL .eqv.
      38             : #elif   SK_ENABLED || IK_ENABLED || CK_ENABLED || RK_ENABLED
      39             : #define IS_EQUAL ==
      40             : #else
      41             : #error  "Unrecognized interface."
      42             : #endif
      43             :         integer(IK)                 :: lenArray, lenArrayNew, lenIndex
      44        5072 :         integer(IK)                 :: IndexNew(size(index))
      45             :         integer(IK)                 :: start, stop, i
      46             :         logical(LK)                 :: negative
      47             :         logical(LK)                 :: unsorted
      48             :         ! Scalar vs. vector `insertion` argument.
      49             : #if     D1_D0_ENABLED
      50             :         integer(IK) , parameter     :: lenInsertion = 1_IK
      51             : #elif   D0_D0_ENABLED   || D1_D1_ENABLED
      52             :         integer(IK)                 :: lenInsertion
      53        1453 :         lenInsertion = GET_SIZE(insertion, kind = IK) ! \warning GET_SIZE is a preprocessor macro.
      54        2858 :         if (lenInsertion == 0_IK) then
      55         360 :             arrayNew = array
      56          18 :             return
      57             :         end if
      58             : #else
      59             : #error  "Unrecognized interface."
      60             : #endif
      61             :         lenIndex = size(index, kind = IK)
      62        4712 :         if (lenIndex > 0_IK) then
      63        1660 :             lenArray = GET_SIZE(array, kind = IK) ! \warning GET_SIZE is a preprocessor macro.
      64        3272 :             lenArrayNew = GET_SIZE(arrayNew, kind = IK) ! \warning GET_SIZE is a preprocessor macro.
      65       39328 :             CHECK_ASSERTION(__LINE__, All(abs(index) <= lenArray + 1_IK), \
      66             :             SK_": The condition `All(abs(index) <= lenArray + 1_IK)` must hold. lenArray, index = "\
      67             :             //getStr([lenArray,index])) ! fpp
      68             : #if         SK_ENABLED && (D1_D0_ENABLED || D1_D1_ENABLED)
      69         492 :             CHECK_ASSERTION(__LINE__, len(array, IK) == len(insertion, IK), \
      70             :             SK_": The condition `len(array) == len(insertion)` must hold. len(array), len(insertion) = "\
      71             :             //getStr([len(array, IK), len(insertion, IK)])) ! fpp
      72             : #endif
      73             : #if         SK_ENABLED && (D1_D0_ENABLED || D1_D1_ENABLED) && setInserted_ENABLED
      74         246 :             CHECK_ASSERTION(__LINE__, len(array, IK) == len(arrayNew, IK), \
      75             :             SK_": The condition `len(array) == len(arrayNew)` must hold. len(array) == len(arrayNew) = "\
      76             :             //getStr([len(array), len(arrayNew)])) ! fpp
      77             : #endif
      78             : #if         (D0_D0_ENABLED || D1_D1_ENABLED) && setInserted_ENABLED
      79        4475 :             CHECK_ASSERTION(__LINE__, lenArrayNew == lenArray + lenIndex * lenInsertion, \
      80             :             SK_": The condition `lenArrayNew /= lenArray + lenIndex * lenInsertion` must hold. lenArrayNew, lenArray, lenIndex, lenInsertion = "\
      81             :             //getStr([lenArrayNew, lenArray, lenIndex, lenInsertion])) ! fpp
      82             : #endif
      83             : #if         D1_D0_ENABLED || D1_D1_ENABLED
      84       15300 :             CHECK_ASSERTION(__LINE__, lenArrayNew == lenArray + lenIndex * lenInsertion, \
      85             :             SK_": The condition lenArrayNew == lenArray + lenIndex * lenInsertion must hold. lenArrayNew, lenArray, lenIndex, lenInsertion = "\
      86             :             //getStr([lenArrayNew, lenArray, lenIndex, lenInsertion])) ! fpp
      87             : #endif
      88             :             ! Set the positivity of index.
      89        3272 :             if (present(positive)) then
      90        1920 :                 negative = .not. positive
      91             :             else
      92             :                 negative = .true._LK
      93             :             end if
      94        1920 :             if (negative) then
      95       11720 :                 do i = 1_IK, lenIndex
      96       11720 :                     if (index(i) > 0_IK) then
      97        5384 :                         IndexNew(i) = index(i)
      98             :                     else
      99        3704 :                         IndexNew(i) = lenArray + index(i) + 1_IK
     100             :                     end if
     101             :                 end do
     102             :             else
     103        2480 :                 IndexNew(1_IK:lenIndex) = index
     104             :             end if
     105             :             ! Sort the index if needed.
     106        3272 :             if (present(sorted)) then
     107        2000 :                 unsorted = .not. sorted
     108             :             else
     109             :                 unsorted = .true._LK
     110             :             end if
     111        3272 :             if (unsorted) call setSorted(IndexNew)
     112             :             ! Insert the `insertion` into the proper locations.
     113        3272 :             if (IndexNew(1) > 1_IK) arrayNew(1_IK:IndexNew(1)-1_IK) = array(1_IK:IndexNew(1)-1_IK)
     114         376 :             start = IndexNew(1)
     115       10928 :             do i = 2_IK, lenIndex
     116        4052 :                 stop = start + lenInsertion - 1_IK
     117       18468 :                 arrayNew(start:stop) = insertion
     118        3604 :                 start = stop + 1_IK
     119        7656 :                 stop = start + IndexNew(i) - IndexNew(i-1_IK) - 1_IK
     120       13232 :                 arrayNew(start:stop) = array(IndexNew(i-1_IK):IndexNew(i)-1_IK)
     121        3272 :                 start = stop + 1_IK
     122             :             end do
     123        1742 :             stop = start + lenInsertion - 1_IK
     124        7862 :             arrayNew(start:stop) = insertion
     125        1530 :             start = stop + 1_IK
     126        3956 :             arrayNew(start:lenArrayNew) = array(IndexNew(lenIndex):lenArray)
     127             :         else ! index is empty.
     128        2808 :             arrayNew = array
     129          72 :             return
     130             :         end if
     131             : 
     132             : #undef  INDEXNEW
     133             : #undef  IS_EQUAL
     134             : #undef  GET_SIZE
     135             : #undef  ALL

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