https://www.cdslab.org/paramonte/fortran/2
Current view: top level - main - pm_matrixCopy@routines.inc.F90 (source / functions) Hit Total Coverage
Test: ParaMonte 2.0.0 :: Serial Fortran - Code Coverage Report Lines: 141 178 79.2 %
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 include file contains procedure implementation of [pm_matrixCopy](@ref pm_matrixCopy).
      19             : !>
      20             : !>  \finmain
      21             : !>
      22             : !>  \author
      23             : !>  \FatemehBagheri, Sunday 11:23 PM, September 19, 2021, Dallas, TX
      24             : 
      25             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      26             : 
      27             :         ! Set the conjugation rule for triangular rfpack formats and Hermitian transpose.
      28             : #if     CK_ENABLED && (THO_ENABLED || (TSO_ENABLED && (RFP_RDP_ENABLED || RDP_RFP_ENABLED)))
      29             : #define GET_CONJG(X)conjg(X)
      30             : #elif   SK_ENABLED || IK_ENABLED || LK_ENABLED || CK_ENABLED || RK_ENABLED
      31             : #define GET_CONJG(X)X
      32             : #else
      33             : #error  "Unrecognized interface."
      34             : #endif
      35             : !       ! Set the conjugate equation rules for triangular rfpack format.
      36             : !#if     RFP_RDP_ENABLED
      37             : !#define EQUATE_CONJG(OBJ1, OBJ2, IND1, IND2) OBJ1 IND1 = GET_CONJG(OBJ2 IND2)
      38             : !#elif   RDP_RFP_ENABLED
      39             : !#define EQUATE_CONJG(OBJ1, OBJ2, IND1, IND2) OBJ1 IND2 = GET_CONJG(OBJ2 IND1)
      40             : !#endif
      41             :         ! Set the slicing rules for different operations.
      42             : #if     AIO_ENABLED
      43             : #define GET_SLICE(I, J)I, J
      44             : #elif   TSO_ENABLED || THO_ENABLED
      45             : #define GET_SLICE(I, J)J, I
      46             : #else
      47             : #error  "Unrecognized interface."
      48             : #endif
      49             :         ! Set the runtime check for string length compatibility.
      50             : #if     SK_ENABLED
      51             : #define CHECK_STRLEN CHECK_ASSERTION(__LINE__, len(source,IK) <= len(destin,IK), \
      52             : SK_"@setMatCopy(): The condition `len(source) <= len(destin)` must hold. len(source), len(destin) = "//getStr([len(source,IK), len(destin,IK)])); ! fpp
      53             : #else
      54             : #define CHECK_STRLEN
      55             : #endif
      56             :         ! Set the runtime check for `doff` given the subset.
      57             : #if     UXD_ENABLED || UXX_ENABLED
      58             : #define CHECK_DOFF CHECK_ASSERTION(__LINE__, -nrow <= doff .and. doff <= 0_IK, \
      59             : SK_"@setMatCopy(): The condition `-nrow <= doff .and. doff <= 0` must hold. nrow, doff = "//getStr([nrow, doff])); ! fpp
      60             : #elif   XLD_ENABLED || XLX_ENABLED
      61             : #define CHECK_DOFF CHECK_ASSERTION(__LINE__, 0_IK <= doff .and. doff <= ncol, \
      62             : SK_"@setMatCopy(): The condition `0 <= doff .and. doff <= ncol` must hold. doff, ncol = "//getStr([doff, ncol])); ! fpp
      63             : #elif   ULX_ENABLED || XXD_ENABLED
      64             : #define CHECK_DOFF CHECK_ASSERTION(__LINE__, -nrow <= doff .and. doff <= ncol, \
      65             : SK_"@setMatCopy(): The condition `-nrow <= doff .and. doff <= ncol` must hold. nrow, doff, ncol = "//getStr([nrow, doff, ncol])); ! fpp
      66             : #else
      67             : #error  "Unrecognized interface."
      68             : #endif
      69             :         ! Define value setting rules for LFP format.
      70             : #if     RDP_LFP_ENABLED
      71             : #define MATRDP destin
      72             : #define MATLFP source
      73             : #define EQUATE(OBJ1, OBJ2, IND1, IND2) OBJ1 IND1 = GET_CONJG(OBJ2 IND2)
      74             : #elif   LFP_RDP_ENABLED
      75             : #define MATLFP destin
      76             : #define MATRDP source
      77             : #define EQUATE(OBJ1, OBJ2, IND1, IND2) OBJ1 IND2 = GET_CONJG(OBJ2 IND1)
      78             : #elif   !(RDP_RDP_ENABLED || RDP_RFP_ENABLED || RFP_RDP_ENABLED)
      79             : #error  "Unrecognized interface."
      80             : #endif
      81             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      82             : #if     getMatCopy_ENABLED && RDP_RDP_ENABLED
      83             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      84             : 
      85      190298 :         if (present(init)) destin = init
      86             : #if     AIO_ENABLED
      87       14680 :         call setMatCopy(destin, rdpack, source, rdpack, subset, doff)
      88             : #elif   TSO_ENABLED || THO_ENABLED
      89          80 :         call setMatCopy(destin, rdpack, source, rdpack, subset, operation, doff)
      90             : #else
      91             : #error  "Unrecognized interface."
      92             : #endif
      93             : 
      94             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      95             : #elif   getMatCopy_ENABLED && RDP_LFP_ENABLED
      96             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      97             : 
      98             :         integer(IK) :: shape_def(2), doff_def
      99           9 :         if (present(doff) .and. present(shape)) then
     100           8 :             shape_def = shape
     101             :             doff_def = doff
     102             :         else
     103           3 :             CHECK_ASSERTION(__LINE__, .not. (present(shape) .or. present(doff)), \
     104             :             SK_"@setMatCopy(): The condition `.not. (present(shape) .or. present(doff))` must hold. present(shape), present(doff) = "\
     105             :             //getStr([present(shape), present(doff)])) ! fpp
     106             : #if         XXD_ENABLED
     107           1 :             shape_def(1) = size(source, 1, IK)
     108             : #elif       UXD_ENABLED || XLD_ENABLED
     109           0 :             shape_def(1) = (getSqrt(size(source, 1, IK) * 8 + 1) - 1) / 2
     110             : #endif
     111           1 :             shape_def(2) = shape_def(1)
     112             :             doff_def = 0_IK
     113             :         end if
     114          18 :         allocate(destin(shape_def(1), shape_def(2)))
     115           9 :         if (present(init)) destin = init
     116             : #if     AIO_ENABLED
     117           6 :         call setMatCopy(destin, dpack, source, spack, subset, doff)
     118             : #elif   TSO_ENABLED || THO_ENABLED
     119           3 :         call setMatCopy(destin, dpack, source, spack, subset, operation, doff)
     120             : #else
     121             : #error  "Unrecognized interface."
     122             : #endif
     123             : 
     124             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     125             : #elif   getMatCopy_ENABLED && LFP_RDP_ENABLED
     126             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     127             : 
     128             : #if     XXD_ENABLED
     129             :         integer(IK) :: doff_def
     130         980 :         doff_def = 0_IK; if (present(doff)) doff_def = doff
     131           3 :         if (doff_def < 0_IK) then
     132           1 :             allocate(destin(min(size(source, 1, IK) + doff_def, size(source, 2, IK))))
     133             :         else
     134         979 :             allocate(destin(min(size(source, 1, IK), size(source, 2, IK) - doff_def)))
     135             :         end if
     136             : #elif   UXD_ENABLED || XLD_ENABLED
     137             :         integer(IK) :: dsize
     138             :         dsize = 0_IK
     139       20465 :         if (all(0_IK < shape(source, IK))) dsize = getMatIndex(dpack, rdpack, shape(source, IK), subset, shape(source, IK), doff)
     140        7867 :         allocate(destin(dsize))
     141             : #else
     142             : #error  "Unrecognized interface."
     143             : #endif
     144        8295 :         if (present(init)) destin = init
     145             : #if     AIO_ENABLED
     146        8292 :         call setMatCopy(destin, dpack, source, spack, subset, doff)
     147             : #elif   TSO_ENABLED || THO_ENABLED
     148           3 :         call setMatCopy(destin, dpack, source, spack, subset, operation, doff)
     149             : #else
     150             : #error  "Unrecognized interface."
     151             : #endif
     152             : 
     153             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     154             : #elif   getMatCopy_ENABLED && RDP_RFP_ENABLED
     155             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     156             : 
     157             :         integer(IK) :: ndim
     158           2 :         ndim = size(source, 1, IK)
     159           2 :         if (mod(size(source, 2, IK), 2_IK) == 0_IK) ndim = ndim - 1
     160           4 :         allocate(destin(ndim, ndim))
     161           2 :         if (present(init)) destin = init
     162             : #if     AIO_ENABLED
     163           2 :         call setMatCopy(destin, dpack, source, spack, subset)
     164             : #elif   TSO_ENABLED || THO_ENABLED
     165             : #error  "Unrecognized interface."
     166             :         !call setMatCopy(destin, dpack, source, spack, subset, operation)
     167             : #else
     168             : #error  "Unrecognized interface."
     169             : #endif
     170             : 
     171             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     172             : #elif   getMatCopy_ENABLED && RFP_RDP_ENABLED
     173             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     174             : 
     175             :         integer(IK) :: ndim
     176           3 :         ndim = size(source, 1, IK)
     177           3 :         if (mod(ndim, 2_IK) == 0_IK) then
     178           0 :             allocate(destin(ndim + 1, ndim / 2))
     179             :         else
     180           6 :             allocate(destin(ndim, ndim / 2 + 1))
     181             :         end if
     182           3 :         if (present(init)) destin = init
     183             : #if     AIO_ENABLED
     184           3 :         call setMatCopy(destin, dpack, source, spack, subset)
     185             : #elif   TSO_ENABLED || THO_ENABLED
     186             : #error  "Unrecognized interface."
     187             :        !call setMatCopy(destin, dpack, source, spack, subset, operation)
     188             : #else
     189             : #error  "Unrecognized interface."
     190             : #endif
     191             : 
     192             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     193             : #elif   setMatCopy_ENABLED && RDP_RDP_ENABLED && (UXX_ENABLED || XLX_ENABLED || XXD_ENABLED || UXD_ENABLED || XLD_ENABLED || ULX_ENABLED)
     194             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     195             : 
     196             :         integer(IK) :: nrow, ncol, icol, doff_def
     197       75981 :         nrow = size(source, 1, IK)
     198       75981 :         ncol = size(source, 2, IK)
     199             :         ! Check the positivity of the lower bounds of `source`.
     200             : #if     AIO_ENABLED
     201      168750 :         CHECK_ASSERTION(__LINE__, all([nrow, ncol] == shape(destin, IK)), \
     202             :         SK_"@setMatCopy(): The condition `all([size(source,1), size(source,2)] == shape(destin, IK))` must hold. shape(source), shape(destin) = "\
     203             :         //getStr([nrow, ncol, shape(destin, IK)])) ! fpp
     204             : #elif   TSO_ENABLED || THO_ENABLED
     205      515079 :         CHECK_ASSERTION(__LINE__, all([ncol, nrow] == shape(destin, IK)), \
     206             :         SK_"@setMatCopy(): The condition `all([size(source,2), size(source,1)] == shape(destin, IK))` must hold. shape(source), shape(destin) = "\
     207             :         //getStr([nrow, ncol, shape(destin, IK)])) ! fpp
     208             : #else
     209             : #error  "Unrecognized interface."
     210             : #endif
     211         636 :         CHECK_STRLEN
     212       75981 :         if (present(doff)) then
     213         540 :             CHECK_DOFF
     214             :             doff_def = doff
     215             :         else
     216             :             doff_def = 0_IK
     217             :         end if
     218             : #if     UXX_ENABLED || XLX_ENABLED
     219             : #define INCREMENT(X,SIGN)X SIGN 1
     220             : #elif   XXD_ENABLED || UXD_ENABLED || XLD_ENABLED || ULX_ENABLED
     221             : #define INCREMENT(X,SIGN)X
     222             : #else
     223             : #error  "Unrecognized interface."
     224             : #endif
     225             :         ! Start the copy action.
     226             :         ! The following loops must remain loop due to a macro expansion limitation.
     227             : #if     UXD_ENABLED || UXX_ENABLED
     228             :         block
     229             :             integer(IK) :: limit
     230       38845 :             limit = min(nrow + doff_def, ncol)
     231             :             do concurrent(icol = 1 : limit)
     232      704627 :                 destin(GET_SLICE(1 : icol - INCREMENT(doff_def,-), icol)) = GET_CONJG(source(1 : icol - INCREMENT(doff_def,-), icol))
     233             :             end do
     234             :             do concurrent(icol = limit + 1 : ncol)
     235       39157 :                 destin(GET_SLICE(1 : nrow, icol)) = GET_CONJG(source(1 : nrow, icol))
     236             :             end do
     237             :         end block
     238             : #elif   XLD_ENABLED || XLX_ENABLED
     239             :         do concurrent(icol = 1 : doff_def)
     240       37242 :             destin(GET_SLICE(1 : nrow, icol)) = GET_CONJG(source(1 : nrow, icol))
     241             :         end do
     242       36972 :         do concurrent(icol = 1 : ncol - doff_def)
     243      688212 :             destin(GET_SLICE(INCREMENT(icol,+) : nrow, icol + doff_def)) = GET_CONJG(source(INCREMENT(icol,+) : nrow, icol + doff_def))
     244             :         end do
     245             : #elif   ULX_ENABLED
     246             :         block
     247             :             integer(IK) :: limit
     248          36 :             if (0_IK == doff_def) then
     249          18 :                 limit = min(nrow, ncol)
     250             :                 do concurrent(icol = 1 : limit)
     251         162 :                     destin(GET_SLICE(1 : icol - 1, icol)) = GET_CONJG(source(1 : icol - 1, icol))
     252         216 :                     destin(GET_SLICE(icol + 1 : nrow, icol)) = GET_CONJG(source(icol + 1 : nrow, icol))
     253             :                 end do
     254             :                 do concurrent(icol = limit + 1 : ncol)
     255          66 :                     destin(GET_SLICE(1 : nrow, icol)) = GET_CONJG(source(1 : nrow, icol))
     256             :                 end do
     257          30 :             elseif (0_IK < doff_def) then
     258             :                 do concurrent(icol = 1 : doff_def)
     259         204 :                     destin(GET_SLICE(1 : nrow, icol)) = GET_CONJG(source(1 : nrow, icol))
     260             :                 end do
     261          14 :                 do concurrent(icol = 1 : ncol - doff_def)
     262          80 :                     destin(GET_SLICE(1 : icol - 1, icol + doff_def)) = GET_CONJG(source(1 : icol - 1, icol + doff_def))
     263         130 :                     destin(GET_SLICE(icol + 1 : nrow, icol + doff_def)) = GET_CONJG(source(icol + 1 : nrow, icol + doff_def))
     264             :                 end do
     265             :             else
     266          16 :                 limit = min(nrow + doff_def, ncol)
     267             :                 do concurrent(icol = 1 : limit)
     268         138 :                     destin(GET_SLICE(1 : icol - doff_def - 1, icol)) = GET_CONJG(source(1 : icol - doff_def - 1, icol))
     269          64 :                     destin(GET_SLICE(icol - doff_def + 1 : nrow, icol)) = GET_CONJG(source(icol - doff_def + 1 : nrow, icol))
     270             :                 end do
     271             :                 do concurrent(icol = limit + 1 : ncol)
     272         296 :                     destin(GET_SLICE(1 : nrow, icol)) = GET_CONJG(source(1 : nrow, icol))
     273             :                 end do
     274             :             end if
     275             :         end block
     276             : #elif   XXD_ENABLED
     277          66 :         if (0_IK == doff_def) then
     278          54 :             do concurrent(icol = 1 : min(nrow, ncol))
     279         252 :                 destin(GET_SLICE(icol, icol)) = GET_CONJG(source(icol, icol))
     280             :             end do
     281          62 :         elseif (0_IK < doff_def) then
     282          30 :             do concurrent(icol = 1 : ncol - doff_def)
     283         112 :                 destin(GET_SLICE(icol, icol + doff_def)) = GET_CONJG(source(icol, icol + doff_def))
     284             :             end do
     285             :         else
     286          32 :             do concurrent(icol = 1 : min(nrow + doff_def, ncol))
     287         102 :                 destin(GET_SLICE(icol - doff_def, icol)) = GET_CONJG(source(icol - doff_def, icol))
     288             :             end do
     289             :         end if
     290             : #else
     291             : #error  "Unrecognized interface."
     292             : #endif
     293             : #undef  INCREMENT
     294             : 
     295             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     296             : #elif   setMatCopy_ENABLED && (RDP_LFP_ENABLED || LFP_RDP_ENABLED) && XXD_ENABLED
     297             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     298             : 
     299             :         integer(IK) :: iell
     300          48 :         CHECK_STRLEN
     301         992 :         if (present(doff)) then
     302          12 :             if (doff < 0_IK) then
     303          32 :                 CHECK_ASSERTION(__LINE__, size(MATLFP, 1, IK) == min(size(MATRDP, 1, IK) + doff, size(MATRDP, 2, IK)), \
     304             :                 SK_"@setMatCopy(): The condition `The number of specified diagonal elements in the two matrix subsets must match. doff, shape(destin), shape(source) = "//\
     305             :                 getStr([doff, shape(destin, IK), shape(source, IK)])) ! fpp
     306             :                 do concurrent(iell = 1 : size(MATLFP, 1, IK))
     307          20 :                     EQUATE(destin, source, (GET_SLICE(iell - doff, iell)), (iell))
     308             :                 end do
     309           4 :                 return
     310           8 :             elseif (0_IK < doff) then
     311          64 :                 CHECK_ASSERTION(__LINE__, \
     312             :                 size(MATLFP, 1, IK) == min(size(MATRDP, 1, IK), size(MATRDP, 2, IK) - doff), \
     313             :                 SK_"@setMatCopy(): The condition `The number of specified diagonal elements in the two matrix subsets must match. doff, shape(destin), shape(source) = "//\
     314             :                 getStr([doff, shape(destin, IK), shape(source, IK)])) ! fpp
     315             :                 do concurrent(iell = 1 : size(MATLFP, 1, IK))
     316          32 :                     EQUATE(destin, source, (GET_SLICE(iell, iell + doff)), (iell))
     317             :                 end do
     318           8 :                 return
     319             :             end if
     320             :         else
     321        8820 :             CHECK_ASSERTION(__LINE__, size(MATLFP, 1, IK) == minval(shape(MATRDP, IK)), \
     322             :             SK_"@setMatCopy(): The condition `The number of specified diagonal elements in the two matrix subsets must match. shape(destin), shape(source) = "//\
     323             :             getStr([shape(destin, IK), shape(source, IK)])) ! fpp
     324             :             do concurrent(iell = 1 : size(MATLFP, 1, IK))
     325        5044 :                 EQUATE(destin, source, (GET_SLICE(iell, iell)), (iell))
     326             :             end do
     327             :         end if
     328             : 
     329             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     330             : #elif   setMatCopy_ENABLED && (RDP_LFP_ENABLED || LFP_RDP_ENABLED)
     331             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     332             : 
     333             :         integer(IK) :: nrow, ncol, irow, icol, nell, doff_def
     334             : #if     AIO_ENABLED
     335        7317 :         nrow = size(MATRDP, 1, IK)
     336        7317 :         ncol = size(MATRDP, 2, IK)
     337             : #elif   TSO_ENABLED || THO_ENABLED
     338          12 :         nrow = size(MATRDP, 2, IK)
     339          12 :         ncol = size(MATRDP, 1, IK)
     340             : #else
     341             : #error  "Unrecognized interface."
     342             : #endif
     343          54 :         CHECK_STRLEN
     344        7329 :         if (present(doff)) then
     345          54 :             CHECK_DOFF
     346          18 :             doff_def = abs(doff)
     347             :         else
     348             :             doff_def = 0_IK
     349             :         end if
     350       58632 :         CHECK_ASSERTION(__LINE__, nrow * ncol - (nrow - doff_def) * (nrow - doff_def - 1) / 2 <= size(MATRDP, kind = IK), \
     351             :         SK_"@setMatCopy(): The condition `The number of specified elements in the two matrix subsets must match. doff, shape(destin), shape(source) = "//\
     352             :         getStr([doff_def, shape(destin, IK), shape(source, IK)])) ! fpp
     353             :         irow = 0_IK
     354             : #if     UXD_ENABLED
     355       21554 :         do icol = 1, min(nrow - doff_def, ncol)
     356       17900 :             nell = icol + doff_def
     357       89528 :             EQUATE(destin, source, (GET_SLICE(1 : nell, icol)), (irow + 1 : irow + nell))
     358        3654 :             irow = irow + nell
     359             :         end do
     360        3668 :         do icol = icol, ncol
     361          14 :             nell = irow + nrow
     362          92 :             EQUATE(destin, source, (GET_SLICE(1 : nrow, icol)), (irow + 1 : nell))
     363        3654 :             irow = nell
     364             :         end do
     365             : #elif   XLD_ENABLED && AIO_ENABLED
     366        3671 :         do icol = 1, doff_def
     367           6 :             nell = irow + nrow
     368          36 :             EQUATE(destin, source, (1 : nrow, icol), (irow + 1 : nell))
     369        3665 :             irow = nell
     370             :         end do
     371        3665 :         irow = irow + 1
     372       21615 :         do icol = 1, ncol - doff_def
     373       89692 :             EQUATE(destin, source, (icol : nrow, icol + doff_def), (irow : irow + nrow - icol))
     374       21615 :             irow = irow + nrow - icol + 1_IK
     375             :         end do
     376             : #elif   XLD_ENABLED && (TSO_ENABLED || THO_ENABLED)
     377          30 :         do icol = 1, nrow - doff_def
     378          20 :             nell = doff_def + icol
     379         110 :             EQUATE(destin, source, (icol, 1 : nell), (irow + 1 : irow + nell))
     380          10 :             irow = irow + nell
     381             :         end do
     382          60 :         do icol = icol, ncol
     383         300 :             EQUATE(destin, source, (icol, 1 : nrow), (irow + 1 : irow + nrow))
     384          10 :             irow = irow + nrow
     385             :         end do
     386             : #else
     387             : #error  "Unrecognized interface."
     388             : #endif
     389             : 
     390             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     391             : #elif   setMatCopy_ENABLED && RFP_RDP_ENABLED && AIO_ENABLED
     392             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     393             : 
     394             :         logical(LK) :: isNotTransMatB
     395             :         integer(IK) :: ndim, ndimHalf, ndimHalfP1 !, irow, icol
     396           6 :         ndim = size(source, 1, IK)
     397             : 
     398          18 :         CHECK_STRLEN
     399             : 
     400           6 :         CHECK_ASSERTION(__LINE__, size(source, 1, IK) == size(source, 2, IK), \
     401             :         SK_"@setMatCopy(): The condition `size(source, 1) == size(source, 2)` must hold. shape(source) = "//getStr(shape(source, IK))) ! fpp
     402             : 
     403             :         ! quick return if possible.
     404             : 
     405           6 :         if(ndim <= 1_IK) then
     406           0 :             CHECK_ASSERTION(__LINE__, all(shape(destin, IK) == shape(source, IK)), \
     407             :             SK_"@setMatCopy(): The condition `all(shape(destin) == shape(source))` must hold. shape(destin), shape(source) = "//getStr([shape(destin, IK), shape(source, IK)])) ! fpp
     408           0 :             if(ndim == 1_IK) destin(1, 1) = GET_CONJG(source(1, 1))
     409           0 :             return
     410             :         end if
     411             : 
     412           6 :         ndimHalf = ndim / 2_IK
     413           6 :         ndimHalfP1 = ndimHalf + 1_IK
     414           6 :         isNotTransMatB = logical(size(destin, 2, IK) < size(destin, 1, IK), LK)
     415           6 :         if (ndim < ndimHalf + ndimHalfP1) then ! ndim is even.
     416             : 
     417           0 :             if (isNotTransMatB) then ! `destin` is in normal form (no transposition).
     418           0 :                 CHECK_ASSERTION(__LINE__, all(shape(destin, IK) == [size(source, 1, IK) + 1_IK, size(source, 1, IK) / 2_IK]), \
     419             :                 SK_"@setMatCopy(): The condition `all(shape(destin) == [size(source, 1) + 1, size(source, 1) / 2]))` must hold. shape(destin), shape(source) = "\
     420             :                 //getStr([shape(destin, IK), shape(source, IK)])) ! fpp
     421             : #if             UXD_ENABLED
     422           0 :                 call setMatCopy(destin(1 : ndim, 1 : ndimHalf)                , rdpack, source(1 : ndim, ndimHalfP1 : ndim) , rdpack, subset = uppDia, doff = -ndimHalf)
     423           0 :                 call setMatCopy(destin(ndimHalf + 2 : ndim + 1, 1 : ndimHalf) , rdpack, source(1 : ndimHalf, 1 : ndimHalf)  , rdpack, subset = uppDia, operation = transHerm)
     424             : #elif           XLD_ENABLED
     425           0 :                 call setMatCopy(destin(2 : ndim + 1, 1 : ndimHalf), rdpack, source(1 : ndim, 1 : ndimHalf)              , rdpack, subset = lowDia)
     426           0 :                 call setMatCopy(destin(1 : ndimHalf, 1 : ndimHalf), rdpack, source(ndimHalfP1 : ndim, ndimHalfP1 : ndim), rdpack, subset = lowDia, operation = transHerm)
     427             : #else
     428             : #error          "Unrecognized interface."
     429             : #endif
     430             :             else ! `destin` is in LAPACK-style transposed form. xxx \todo Is the second transposition Hermitian? must be checked with lapack.
     431           0 :                 CHECK_ASSERTION(__LINE__, all(shape(destin, IK) == [size(source, 1, IK) / 2_IK, size(source, 1, IK) + 1_IK]), \
     432             :                 SK_"@setMatCopy(): The condition `all(shape(destin) == [size(source, 1) / 2, size(source, 1) + 1])` must hold. shape(destin), shape(source) = "\
     433             :                 //getStr([shape(destin, IK), shape(source, IK)])) ! fpp
     434             : #if             UXD_ENABLED
     435           0 :                 call setMatCopy(destin(1 : ndimHalf, 1 : ndim)                , rdpack, source(1 : ndim, ndimHalfP1 : ndim) , rdpack, subset = uppDia, operation = transHerm, doff = -ndimHalf)
     436           0 :                 call setMatCopy(destin(1 : ndimHalf, ndimHalf + 2 : ndim + 1) , rdpack, source(1 : ndimHalf, 1 : ndimHalf)  , rdpack, subset = uppDia)
     437             : #elif           XLD_ENABLED
     438           0 :                 call setMatCopy(destin(1 : ndimHalf, 2 : ndim + 1), rdpack, source(1 : ndim, 1 : ndimHalf)              , rdpack, subset = lowDia, operation = transHerm)
     439           0 :                 call setMatCopy(destin(1 : ndimHalf, 1 : ndimHalf), rdpack, source(ndimHalfP1 : ndim, ndimHalfP1 : ndim), rdpack, subset = lowDia)
     440             : #else
     441             : #error          "Unrecognized interface."
     442             : #endif
     443             :             end if
     444             : 
     445             :         else ! ndim is odd.
     446             : 
     447           6 :             if (isNotTransMatB) then ! `destin` is in normal form (no transposition).
     448          55 :                 CHECK_ASSERTION(__LINE__, all(shape(destin, IK) == [size(source, 1, IK), size(source, 1, IK) / 2 + 1]), \
     449             :                 SK_"@setMatCopy(): The condition `all(shape(destin) == [size(source, 1), size(source, 1) / 2 + 1]))` must hold. shape(destin), shape(source) = "\
     450             :                 //getStr([shape(destin, IK), shape(source, IK)])) ! fpp
     451             : #if             UXD_ENABLED
     452           3 :                 call setMatCopy(destin(1 : ndim, 1 : ndimHalfP1)          , rdpack, source(1 : ndim, ndimHalfP1 : ndim) , rdpack, subset = uppDia, doff = -ndimHalf)
     453           3 :                 call setMatCopy(destin(ndimHalf + 2 : ndim, 1 : ndimHalf) , rdpack, source(1 : ndimHalf, 1 : ndimHalf)  , rdpack, subset = uppDia, operation = transHerm)
     454             : #elif           XLD_ENABLED
     455           2 :                 call setMatCopy(destin(1 : ndim, 1 : ndimHalfP1)      , rdpack, source(1 : ndim, 1 : ndimHalfP1)                    , rdpack, subset = lowDia)
     456           2 :                 call setMatCopy(destin(1 : ndimHalf, 2 : ndimHalfP1)  , rdpack, source(ndimHalfP1 + 1 : ndim, ndimHalfP1 + 1 : ndim), rdpack, subset = lowDia, operation = transHerm)
     457             : #else
     458             : #error          "Unrecognized interface."
     459             : #endif
     460             :             else ! `destin` is in LAPACK-style transposed form. xxx \todo Is the second transposition Hermitian? must be checked with lapack.
     461          11 :                 CHECK_ASSERTION(__LINE__, all(shape(destin, IK) == [size(source, 1, IK) / 2 + 1, size(source, 1, IK)]), \
     462             :                 SK_"@setMatCopy(): The condition `all(shape(destin) == [size(source, 1) / 2 + 1, size(source, 1)]))` must hold. shape(destin), shape(source) = "\
     463             :                 //getStr([shape(destin, IK), shape(source, IK)])) ! fpp
     464             : #if             UXD_ENABLED
     465           1 :                 call setMatCopy(destin(1 : ndimHalfP1, 1 : ndim)          , rdpack, source(1 : ndim, ndimHalfP1 : ndim) , rdpack, subset = uppDia, operation = transHerm, doff = -ndimHalf)
     466           1 :                 call setMatCopy(destin(1 : ndimHalf, ndimHalf + 2 : ndim) , rdpack, source(1 : ndimHalf, 1 : ndimHalf)  , rdpack, subset = uppDia)
     467             : #elif           XLD_ENABLED
     468           0 :                 call setMatCopy(destin(1 : ndimHalfP1, 1 : ndim)      , rdpack, source(1 : ndim, 1 : ndimHalfP1)            , rdpack, subset = lowDia, operation = transHerm)
     469           0 :                 call setMatCopy(destin(2 : ndimHalfP1, 1 : ndimHalf)  , rdpack, source(ndimHalfP1 : ndim, ndimHalfP1 : ndim), rdpack, subset = lowDia)
     470             : #else
     471             : #error          "Unrecognized interface."
     472             : #endif
     473             :             end if
     474             : 
     475             :         end if
     476             : 
     477             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     478             : #elif   setMatCopy_ENABLED && RDP_RFP_ENABLED && AIO_ENABLED
     479             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     480             : 
     481             :         logical(LK) :: isNotTransMatA
     482             :         integer(IK) :: ndim, ndimHalf, ndimHalfP1
     483           4 :         ndim = size(destin, 1, IK)
     484             : 
     485          12 :         CHECK_STRLEN
     486           4 :         CHECK_ASSERTION(__LINE__, size(destin, 1, IK) == size(destin, 2, IK), \
     487             :         SK_"@setMatCopy(): The condition `size(destin, 1) == size(destin, 2)` must hold. shape(destin) = "//getStr(shape(destin, IK))) ! fpp
     488             : 
     489             :         ! quick return if possible.
     490             : 
     491           4 :         if(ndim <= 1_IK) then
     492           0 :             CHECK_ASSERTION(__LINE__, all(shape(destin, IK) == shape(source, IK)), \
     493             :             SK_"@setMatCopy(): The condition `all(shape(destin) == shape(source))` must hold. shape(destin), shape(source) = "//getStr([shape(destin, IK), shape(source, IK)])) ! fpp
     494           0 :             if(ndim == 1_IK) destin(1, 1) = GET_CONJG(source(1, 1))
     495           0 :             return
     496             :         end if
     497             : 
     498           4 :         ndimHalf = ndim / 2_IK
     499           4 :         ndimHalfP1 = ndimHalf + 1_IK
     500           4 :         isNotTransMatA = logical(size(source, 2, IK) < size(source, 1, IK), LK)
     501           4 :         if (ndim < ndimHalf + ndimHalfP1) then ! ndim is even.
     502             : 
     503           0 :             if (isNotTransMatA) then ! `source` is in normal form (no transposition).
     504           0 :                 CHECK_ASSERTION(__LINE__, all(shape(source, IK) == [size(destin, 1, IK) + 1_IK, size(destin, 1, IK) / 2_IK]), \
     505             :                 SK_"@setMatCopy(): The condition `all(shape(source) == [size(destin, 1) + 1, size(destin, 1) / 2]))` must hold. shape(source), shape(destin) = "\
     506             :                 //getStr([shape(source, IK), shape(destin, IK)])) ! fpp
     507             : #if             UXD_ENABLED
     508           0 :                 call setMatCopy(destin(1 : ndim, ndimHalfP1 : ndim)   , rdpack, source(1 : ndim, 1 : ndimHalf)                  , rdpack, subset = uppDia, doff = -ndimHalf)
     509           0 :                 call setMatCopy(destin(1 : ndimHalf, 1 : ndimHalf)    , rdpack, source(ndimHalf + 2 : ndim + 1, 1 : ndimHalf)   , rdpack, subset = uppDia, operation = transHerm)
     510             : #elif           XLD_ENABLED
     511           0 :                 call setMatCopy(destin(1 : ndim, 1 : ndimHalf)                , rdpack, source(2 : ndim + 1, 1 : ndimHalf), rdpack, subset = lowDia)
     512           0 :                 call setMatCopy(destin(ndimHalfP1 : ndim, ndimHalfP1 : ndim)  , rdpack, source(1 : ndimHalf, 1 : ndimHalf), rdpack, subset = lowDia, operation = transHerm)
     513             : #else
     514             : #error          "Unrecognized interface."
     515             : #endif
     516             :             else ! `source` is in LAPACK-style transposed form. xxx \todo Is the second transposition Hermitian? must be checked with lapack.
     517           0 :                 CHECK_ASSERTION(__LINE__, all(shape(source, IK) == [size(destin, 1, IK) / 2_IK, size(destin, 1, IK) + 1_IK]), \
     518             :                 SK_"@setMatCopy(): The condition `all(shape(source) == [size(destin, 1) / 2, size(destin, 1) + 1])` must hold. shape(source), shape(destin) = "\
     519             :                 //getStr([shape(source, IK), shape(destin, IK)])) ! fpp
     520             : #if             UXD_ENABLED
     521           0 :                 call setMatCopy(destin(1 : ndim, ndimHalfP1 : ndim)   , rdpack, source(1 : ndimHalf, 1 : ndim)                  , rdpack, subset = uppDia, operation = transHerm, doff = -ndimHalf)
     522           0 :                 call setMatCopy(destin(1 : ndimHalf, 1 : ndimHalf)    , rdpack, source(1 : ndimHalf, ndimHalf + 2 : ndim + 1)   , rdpack, subset = uppDia)
     523             : #elif           XLD_ENABLED
     524           0 :                 call setMatCopy(destin(1 : ndim, 1 : ndimHalf)                , rdpack, source(1 : ndimHalf, 2 : ndim + 1), rdpack, subset = lowDia, operation = transHerm)
     525           0 :                 call setMatCopy(destin(ndimHalfP1 : ndim, ndimHalfP1 : ndim)  , rdpack, source(1 : ndimHalf, 1 : ndimHalf), rdpack, subset = lowDia)
     526             : #else
     527             : #error          "Unrecognized interface."
     528             : #endif
     529             :             end if
     530             : 
     531             :         else ! ndim is odd.
     532             : 
     533           4 :             if (isNotTransMatA) then ! `source` is in normal form (no transposition).
     534          44 :                 CHECK_ASSERTION(__LINE__, all(shape(source, IK) == [size(destin, 1, IK), size(destin, 1, IK) / 2 + 1]), \
     535             :                 SK_"@setMatCopy(): The condition `all(shape(source) == [size(destin, 1), size(destin, 1) / 2 + 1]))` must hold. shape(source), shape(destin) = "\
     536             :                 //getStr([shape(source, IK), shape(destin, IK)])) ! fpp
     537             : #if             UXD_ENABLED
     538           2 :                 call setMatCopy(destin(1 : ndim, ndimHalfP1 : ndim)   , rdpack, source(1 : ndim, 1 : ndimHalfP1)            , rdpack, subset = uppDia, doff = -ndimHalf)
     539           2 :                 call setMatCopy(destin(1 : ndimHalf, 1 : ndimHalf)    , rdpack, source(ndimHalf + 2 : ndim, 1 : ndimHalf)   , rdpack, subset = lowDia, operation = transHerm)
     540             : #elif           XLD_ENABLED
     541           2 :                 call setMatCopy(destin(1 : ndim, 1 : ndimHalfP1)                      , rdpack, source(1 : ndim, 1 : ndimHalfP1)    , rdpack, subset = lowDia)
     542           2 :                 call setMatCopy(destin(ndimHalfP1 + 1 : ndim, ndimHalfP1 + 1 : ndim)  , rdpack, source(1 : ndimHalf, 2 : ndimHalfP1), rdpack, subset = uppDia, operation = transHerm)
     543             : #else
     544             : #error          "Unrecognized interface."
     545             : #endif
     546             :             else ! `source` is in LAPACK-style transposed form. xxx \todo Is the second transposition Hermitian? must be checked with lapack.
     547           0 :                 CHECK_ASSERTION(__LINE__, all(shape(source, IK) == [size(destin, 1, IK) / 2 + 1, size(destin, 1, IK)]), \
     548             :                 SK_"@setMatCopy(): The condition `all(shape(source) == [size(destin, 1) / 2 + 1, size(destin, 1)]))` must hold. shape(source), shape(destin) = "\
     549             :                 //getStr([shape(source, IK), shape(destin, IK)])) ! fpp
     550             : #if             UXD_ENABLED
     551           0 :                 call setMatCopy(destin(1 : ndim, ndimHalfP1 : ndim)   , rdpack, source(1 : ndimHalfP1, 1 : ndim)            , rdpack, subset = uppDia, operation = transHerm, doff = -ndimHalf)
     552           0 :                 call setMatCopy(destin(1 : ndimHalf, 1 : ndimHalf)    , rdpack, source(1 : ndimHalf, ndimHalf + 2 : ndim)   , rdpack, subset = uppDia)
     553             : #elif           XLD_ENABLED
     554           0 :                 call setMatCopy(destin(1 : ndim, 1 : ndimHalfP1)              , rdpack, source(1 : ndimHalfP1, 1 : ndim)    , rdpack, subset = lowDia, operation = transHerm)
     555           0 :                 call setMatCopy(destin(ndimHalfP1 : ndim, ndimHalfP1 : ndim)  , rdpack, source(2 : ndimHalfP1, 1 : ndimHalf), rdpack, subset = lowDia)
     556             : #else
     557             : #error          "Unrecognized interface."
     558             : #endif
     559             :             end if
     560             : 
     561             :         end if
     562             : 
     563             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     564             : #elif   setMatCopy_ENABLED && RFP_RDP_ENABLED && (TSO_ENABLED || THO_ENABLED)
     565             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     566             : 
     567             : #error  "Unrecognized interface."
     568             : 
     569             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     570             : #elif   setMatCopy_ENABLED && RDP_RFP_ENABLED && (TSO_ENABLED || THO_ENABLED)
     571             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     572             : 
     573             : #error  "Unrecognized interface."
     574             : 
     575             : #else
     576             :         !%%%%%%%%%%%%%%%%%%%%%%%%
     577             : #error  "Unrecognized interface."
     578             :         !%%%%%%%%%%%%%%%%%%%%%%%%
     579             : #endif
     580             : 
     581             : #undef CHECK_STRLEN
     582             : #undef CHECK_DOFF
     583             : #undef GET_CONJG
     584             : #undef GET_SLICE
     585             : #undef MATLFP
     586             : #undef MATRDP
     587             : #undef EQUATE

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