https://www.cdslab.org/paramonte/fortran/2
Current view: top level - main - pm_matrixInit@routines.inc.F90 (source / functions) Hit Total Coverage
Test: ParaMonte 2.0.0 :: Serial Fortran - Code Coverage Report Lines: 166 171 97.1 %
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 the generic interface [pm_matrixInit](@ref pm_matrixInit).
      19             : !>
      20             : !>  \finmain
      21             : !>
      22             : !>  \author
      23             : !>  \AmirShahmoradi, Thursday 01:00 AM, September 23, 2021, Dallas, TX
      24             : 
      25             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      26             : 
      27             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      28             : #if     getMatInit_ENABLED && XXD_ENABLED
      29             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      30             : 
      31             :         integer(IK) :: roff_def, coff_def
      32             : 
      33             :         ! set the row offset.
      34           9 :         if (present(roff)) then
      35           4 :             roff_def = roff
      36             :         else
      37           5 :             roff_def = 0_IK
      38             :         end if
      39             : 
      40             :         ! set the column offset.
      41           9 :         if (present(coff)) then
      42           4 :             coff_def = coff
      43             :         else
      44           5 :             coff_def = 0_IK
      45             :         end if
      46             : 
      47             : #if     D2XX0_ENABLED
      48             :         ! set the number of the diagonal elements initialized.
      49           6 :         call setMatInit(mat, subset, vdia, getOption(min(shape(1) - roff_def, shape(2) - coff_def), ndia), roff_def, coff_def)
      50             : #elif   D2XX1_ENABLED
      51           3 :         call setMatInit(mat, subset, vdia, roff_def, coff_def)
      52             : #else
      53             : #error  "Unrecognized interface."
      54             : #endif
      55             : 
      56             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      57             : #elif   getMatInit_ENABLED && XLD_ENABLED
      58             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      59             : 
      60             :         integer(IK) :: roff_def, coff_def
      61             :         integer(IK) :: nrow_def, ncol_def
      62             : 
      63             :         ! set the row offset.
      64          10 :         if (present(roff)) then
      65           6 :             roff_def = roff
      66             :         else
      67           4 :             roff_def = 0_IK
      68             :         end if
      69             : 
      70             :         ! set the column offset.
      71          10 :         if (present(coff)) then
      72           5 :             coff_def = coff
      73             :         else
      74           5 :             coff_def = 0_IK
      75             :         end if
      76             : 
      77             :         ! set the number of the rows initialized.
      78          10 :         if (present(nrow)) then
      79           7 :             nrow_def = nrow
      80             :         else
      81           3 :             nrow_def = shape(1) - roff_def
      82             :         end if
      83             : 
      84             :         ! set the number of the columns initialized.
      85             :         if (present(ncol)) then
      86             :             ncol_def = ncol
      87             :         else
      88             :             ncol_def = shape(2) - coff_def
      89             :         end if
      90             : 
      91             :         ! set the number of the columns initialized.
      92          10 :         if (present(ncol)) then
      93           7 :             ncol_def = ncol
      94             :         else
      95           3 :             ncol_def = shape(2) - coff_def
      96             :         end if
      97             : 
      98          10 :         call setMatInit(mat, subset, vlow, vdia, nrow_def, ncol_def, roff_def, coff_def, doff)
      99             : 
     100             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     101             : #elif   getMatInit_ENABLED && UXD_ENABLED
     102             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     103             : 
     104             :         integer(IK) :: roff_def, coff_def
     105             :         integer(IK) :: nrow_def, ncol_def
     106             : 
     107             :         ! set the row offset.
     108          10 :         if (present(roff)) then
     109           6 :             roff_def = roff
     110             :         else
     111           4 :             roff_def = 0_IK
     112             :         end if
     113             : 
     114             :         ! set the column offset.
     115          10 :         if (present(coff)) then
     116           5 :             coff_def = coff
     117             :         else
     118           5 :             coff_def = 0_IK
     119             :         end if
     120             : 
     121             :         ! set the number of the rows initialized.
     122          10 :         if (present(nrow)) then
     123           7 :             nrow_def = nrow
     124             :         else
     125           3 :             nrow_def = shape(1) - roff_def
     126             :         end if
     127             : 
     128             :         ! set the number of the columns initialized.
     129             :         if (present(ncol)) then
     130             :             ncol_def = ncol
     131             :         else
     132             :             ncol_def = shape(2) - coff_def
     133             :         end if
     134             : 
     135             :         ! set the number of the columns initialized.
     136          10 :         if (present(ncol)) then
     137           7 :             ncol_def = ncol
     138             :         else
     139           3 :             ncol_def = shape(2) - coff_def
     140             :         end if
     141             : 
     142          10 :         call setMatInit(mat, subset, vupp, vdia, nrow_def, ncol_def, roff_def, coff_def, doff)
     143             : 
     144             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     145             : #elif   getMatInit_ENABLED && ULX_ENABLED
     146             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     147             : 
     148             :         integer(IK) :: roff_def, coff_def
     149             :         integer(IK) :: nrow_def, ncol_def
     150             : 
     151             :         ! set the row offset.
     152          11 :         if (present(roff)) then
     153           7 :             roff_def = roff
     154             :         else
     155           4 :             roff_def = 0_IK
     156             :         end if
     157             : 
     158             :         ! set the column offset.
     159          11 :         if (present(coff)) then
     160           5 :             coff_def = coff
     161             :         else
     162           6 :             coff_def = 0_IK
     163             :         end if
     164             : 
     165             :         ! set the number of the rows initialized.
     166          11 :         if (present(nrow)) then
     167           8 :             nrow_def = nrow
     168             :         else
     169           3 :             nrow_def = shape(1) - roff_def
     170             :         end if
     171             : 
     172             :         ! set the number of the columns initialized.
     173             :         if (present(ncol)) then
     174             :             ncol_def = ncol
     175             :         else
     176             :             ncol_def = shape(2) - coff_def
     177             :         end if
     178             : 
     179             :         ! set the number of the columns initialized.
     180          11 :         if (present(ncol)) then
     181           8 :             ncol_def = ncol
     182             :         else
     183           3 :             ncol_def = shape(2) - coff_def
     184             :         end if
     185             : 
     186          11 :         call setMatInit(mat, subset, vupp, vlow, nrow_def, ncol_def, roff_def, coff_def, doff)
     187             : 
     188             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     189             : #elif   getMatInit_ENABLED && ULD_ENABLED
     190             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     191             : 
     192             :         integer(IK) :: roff_def, coff_def
     193             :         integer(IK) :: nrow_def, ncol_def
     194             : 
     195             :         ! set the row offset.
     196        3335 :         if (present(roff)) then
     197           7 :             roff_def = roff
     198             :         else
     199        3328 :             roff_def = 0_IK
     200             :         end if
     201             : 
     202             :         ! set the column offset.
     203        3335 :         if (present(coff)) then
     204           5 :             coff_def = coff
     205             :         else
     206        3330 :             coff_def = 0_IK
     207             :         end if
     208             : 
     209             :         ! set the number of the rows initialized.
     210        3335 :         if (present(nrow)) then
     211           8 :             nrow_def = nrow
     212             :         else
     213        3327 :             nrow_def = shape(1) - roff_def
     214             :         end if
     215             : 
     216             :         ! set the number of the columns initialized.
     217             :         if (present(ncol)) then
     218             :             ncol_def = ncol
     219             :         else
     220             :             ncol_def = shape(2) - coff_def
     221             :         end if
     222             : 
     223             :         ! set the number of the columns initialized.
     224        3335 :         if (present(ncol)) then
     225           8 :             ncol_def = ncol
     226             :         else
     227        3327 :             ncol_def = shape(2) - coff_def
     228             :         end if
     229             : 
     230        3335 :         call setMatInit(mat, subset, vupp, vlow, vdia, nrow_def, ncol_def, roff_def, coff_def, doff)
     231             : 
     232             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     233             : #elif   setMatInit_ENABLED && XXD_ENABLED
     234             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     235             : 
     236             :         integer(IK) :: i
     237             : #if     IMP_ENABLED
     238             :         integer(IK), parameter :: roff = 0_IK, coff = 0_IK
     239             : #elif   !EXP_ENABLED
     240             : #error  "Unrecognized interface."
     241             : #endif
     242             : #if     D2XXF_ENABLED
     243             :         integer(IK) :: ndia
     244        4400 :         ndia = min(size(mat, 1, IK) - roff, size(mat, 2, IK) - coff)
     245             : #define GET_DIA(i) vdia
     246             : #elif   D2XX0_ENABLED
     247             : #define NDIA ndia
     248             : #define GET_DIA(i) vdia
     249         160 :         CHECK_ASSERTION(__LINE__, ndia <= minval(shape(mat, IK) - [roff, coff]), SK_"@setMatInit(): The condition `ndia <= minval(shape(mat, IK) - [roff, coff])` must hold. ndia, shape(mat), roff, coff = "//getStr([ndia, shape(mat, kind = IK), roff, coff]))
     250             : #elif   D2XX1_ENABLED
     251             : #define NDIA size(vdia, kind = IK)
     252             : #define GET_DIA(i) vdia(i)
     253        8060 :         CHECK_ASSERTION(__LINE__, NDIA <= minval(shape(mat, IK) - [roff, coff]), SK_"@setMatInit(): The condition `size(vdia) <= minval(shape(mat, IK) - [roff, coff])` must hold. size(vdia), shape(mat), roff, coff = "//getStr([NDIA, shape(mat, kind = IK), roff, coff]))
     254             : #else
     255             : #error  "Unrecognized interface."
     256             : #endif
     257             : #if     SK_ENABLED
     258          21 :         CHECK_ASSERTION(__LINE__, len(vdia, IK) <= len(mat, IK), \
     259             :         SK_"@setMatInit(): The condition `len(vdia) <= len(mat)` must hold. len(vdia), len(mat) = "//getStr([len(vdia, IK), len(mat, IK)]))
     260             : #endif
     261        5222 :         CHECK_ASSERTION(__LINE__, 0_IK <= roff, SK_"@setMatInit(): The condition `0_IK <= roff` must hold. roff = "//getStr(roff))
     262        5222 :         CHECK_ASSERTION(__LINE__, 0_IK <= coff, SK_"@setMatInit(): The condition `0_IK <= coff` must hold. coff = "//getStr(coff))
     263             :         do concurrent(i = 1 : NDIA)
     264       21700 :             mat(i + roff, i + coff) = GET_DIA(i)
     265             :         end do
     266             : 
     267             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     268             : #elif   setMatInit_ENABLED && (XLX_ENABLED || XLD_ENABLED)
     269             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     270             : 
     271             :         integer(IK) :: i, lenDia, doff_def
     272             : #if     IMP_ENABLED
     273             :         integer(IK), parameter :: roff = 0_IK, coff = 0_IK
     274             :         integer(IK) :: nrow, ncol
     275        3300 :         nrow = size(mat, 1, IK)
     276        3300 :         ncol = size(mat, 2, IK)
     277             : #elif   EXP_ENABLED
     278          18 :         CHECK_ASSERTION(__LINE__, 0_IK <= roff, SK_"@setMatInit(): The condition `0_IK <= roff` must hold. roff = "//getStr(roff))
     279          18 :         CHECK_ASSERTION(__LINE__, 0_IK <= coff, SK_"@setMatInit(): The condition `0_IK <= coff` must hold. coff = "//getStr(coff))
     280          72 :         CHECK_ASSERTION(__LINE__, 0_IK <= nrow .and. nrow + roff <= size(mat, 1, IK), SK_"@setMatInit(): The condition `0_IK <= nrow .and. nrow + roff <= size(mat,1)` must hold. nrow, roff, size(mat,1) = "//getStr([nrow, roff, size(mat, 1, IK)]))
     281          72 :         CHECK_ASSERTION(__LINE__, 0_IK <= ncol .and. ncol + coff <= size(mat, 2, IK), SK_"@setMatInit(): The condition `0_IK <= ncol .and. ncol + coff <= size(mat,2)` must hold. ncol, coff, size(mat,2) = "//getStr([ncol, coff, size(mat, 2, IK)]))
     282             : #else
     283             : #error  "Unrecognized interface."
     284             : #endif
     285             : #if     SK_ENABLED
     286          27 :         CHECK_ASSERTION(__LINE__, len(vlow, IK) <= len(mat, IK), SK_"@setMatInit(): The condition `len(vlow) <= len(mat)` must hold. len(vlow), len(mat) = "//getStr([len(vlow, IK), len(mat, IK)]))
     287             : #if     XLD_ENABLED
     288          27 :         CHECK_ASSERTION(__LINE__, len(vdia, IK) <= len(mat, IK), SK_"@setMatInit(): The condition `len(vdia) <= len(mat)` must hold. len(vdia), len(mat) = "//getStr([len(vdia, IK), len(mat, IK)]))
     289             : #endif
     290             : #endif
     291        3318 :         if (nrow < 1_IK .or. ncol < 1_IK) return ! this must be here.
     292        3310 :         if (present(doff)) then
     293          12 :             CHECK_ASSERTION(__LINE__, 0 <= doff .and. doff < ncol, SK_"@setMatInit(): The condition `0 <= doff .and. doff < ncol` must hold. doff, ncol = "//getStr([doff, ncol]))
     294             :             doff_def = doff
     295             :         else
     296             :             doff_def = 0_IK
     297             :         endif
     298       13240 :         CHECK_ASSERTION(__LINE__, ncol - doff_def <= nrow, SK_"@setMatInit(): The condition `ncol - doff <= nrow` must hold. You cannot ask for more columns than possible. ncol, doff, nrow = "//getStr([ncol, doff_def, nrow]))
     299             : #if     D2X00_ENABLED || D2X0X_ENABLED
     300             : #define GET_DIA(i) vdia
     301             :         lenDia = ncol - doff_def ! min(ncol - doff_def, nrow)
     302             : #elif   D2X01_ENABLED
     303             : #define GET_DIA(i) vdia(i)
     304           2 :         lenDia = size(vdia, 1, IK)
     305          10 :         CHECK_ASSERTION(__LINE__, lenDia == ncol - doff_def, SK_"@setMatInit(): The condition `size(vdia) == ncol - doff` must hold. size(vdia), nrow, doff, ncol = "//getStr([lenDia, nrow, doff_def, ncol]))
     306             : #else
     307             : #error  "Unrecognized interface."
     308             : #endif
     309             :         !check_assertion(__LINE__, all(0_IK < [nrow, ncol]), SK_"@setMatInit(): The condition `all(0 < [nrow, ncol])` must hold. nrow, ncol = "//getStr([[nrow, ncol]]))
     310             :         do concurrent(i = 1 : doff_def)
     311        3366 :             mat(1 : nrow, i) = vlow
     312             :         end do
     313             : #if     XLD_ENABLED
     314           9 :         mat(1 , 1 + doff_def) = GET_DIA(1)
     315             : #endif
     316       10476 :         mat(2 : nrow , 1 + doff_def) = vlow
     317       10439 :         do i = 2_IK, lenDia
     318             : #if         XLD_ENABLED
     319          26 :             mat(i , i + doff_def) = GET_DIA(i)
     320             : #endif
     321       18681 :             mat(i + 1 : nrow , i + doff_def) = vlow
     322             :         end do
     323             : 
     324             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     325             : #elif   setMatInit_ENABLED && (UXX_ENABLED || UXD_ENABLED)
     326             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     327             : 
     328             :         integer(IK) :: i, lenDia, doff_def
     329             : #if     IMP_ENABLED
     330             :         integer(IK), parameter :: roff = 0_IK, coff = 0_IK
     331             :         integer(IK) :: nrow, ncol
     332        3690 :         nrow = size(mat, 1, IK)
     333        3690 :         ncol = size(mat, 2, IK)
     334             : #elif   EXP_ENABLED
     335          18 :         CHECK_ASSERTION(__LINE__, 0_IK <= roff, SK_"@setMatInit(): The condition `0_IK <= roff` must hold. roff = "//getStr(roff))
     336          18 :         CHECK_ASSERTION(__LINE__, 0_IK <= coff, SK_"@setMatInit(): The condition `0_IK <= coff` must hold. coff = "//getStr(coff))
     337          72 :         CHECK_ASSERTION(__LINE__, 0_IK <= nrow .and. nrow + roff <= size(mat, 1, IK), SK_"@setMatInit(): The condition `0_IK <= nrow .and. nrow + roff <= size(mat,1)` must hold. nrow, roff, size(mat,1) = "//getStr([nrow, roff, size(mat, 1, IK)]))
     338          72 :         CHECK_ASSERTION(__LINE__, 0_IK <= ncol .and. ncol + coff <= size(mat, 2, IK), SK_"@setMatInit(): The condition `0_IK <= ncol .and. ncol + coff <= size(mat,2)` must hold. ncol, coff, size(mat,2) = "//getStr([ncol, coff, size(mat, 2, IK)]))
     339             : #else
     340             : #error  "Unrecognized interface."
     341             : #endif
     342             : #if     SK_ENABLED
     343          27 :         CHECK_ASSERTION(__LINE__, len(vupp, IK) <= len(mat, IK), SK_"@setMatInit(): The condition `len(vupp) <= len(mat)` must hold. len(vupp), len(mat) = "//getStr([len(vupp, IK), len(mat, IK)]))
     344             : #if     UXD_ENABLED
     345          27 :         CHECK_ASSERTION(__LINE__, len(vdia, IK) <= len(mat, IK), SK_"@setMatInit(): The condition `len(vdia) <= len(mat)` must hold. len(vdia), len(mat) = "//getStr([len(vdia, IK), len(mat, IK)]))
     346             : #endif
     347             : #endif
     348             :         !check_assertion(__LINE__, all(0_IK < [nrow, ncol]), SK_"@setMatInit(): The condition `all(0 < [nrow, ncol])` must hold. nrow, ncol = "//getStr([[nrow, ncol]]))
     349        3708 :         if (nrow < 1_IK .or. ncol < 1_IK) return ! this must be here.
     350        3705 :         if (present(doff)) then
     351          12 :             CHECK_ASSERTION(__LINE__, -nrow < doff .and. doff <= 0_IK, SK_"@setMatInit(): The condition `-nrow < doff .and. doff <= 0_IK` must hold. nrow, doff = "//getStr([nrow, doff]))
     352           4 :             doff_def = -doff
     353             :         else
     354             :             doff_def = 0_IK
     355             :         endif
     356       14820 :         CHECK_ASSERTION(__LINE__, nrow - doff_def <= ncol, SK_"@setMatInit(): The condition `nrow + doff <= ncol` must hold. You cannot ask for more rows than possible. nrow, doff, ncol = "//getStr([nrow, -doff_def, ncol]))
     357             : #if     D20X0_ENABLED || D20XX_ENABLED
     358             : #define GET_DIA(i) vdia
     359             :         lenDia = nrow - doff_def ! min(nrow - doff_def, ncol)
     360             : #elif   D20X1_ENABLED
     361          12 :         lenDia = size(vdia, 1, IK)
     362          48 :         CHECK_ASSERTION(__LINE__, lenDia == nrow - doff_def, SK_"@setMatInit(): The condition `size(vdia) == nrow + doff` must hold. size(vdia), nrow, doff = "//getStr([lenDia, nrow, doff_def]))
     363             : #define GET_DIA(i) vdia(i)
     364             : #else
     365             : #error  "Unrecognized interface."
     366             : #endif
     367       15609 :         do i = 1_IK, lenDia
     368       31147 :             mat(1 : i + doff_def - 1, i) = vupp
     369             : #if         UXD_ENABLED
     370       10074 :             mat(i + doff_def , i) = GET_DIA(i)
     371             : #endif
     372             :         end do
     373        3759 :         do i = lenDia + 1, ncol
     374        3982 :             mat(1 : nrow, i) = vupp
     375             :         end do
     376             : 
     377             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     378             : #elif   setMatInit_ENABLED && ULX_ENABLED
     379             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     380             : 
     381             :         integer(IK) :: i, lenDia, doff_def
     382             : #if     IMP_ENABLED
     383             :         integer(IK), parameter :: roff = 0_IK, coff = 0_IK
     384             :         integer(IK) :: nrow, ncol
     385           0 :         nrow = size(mat, 1, IK)
     386           0 :         ncol = size(mat, 2, IK)
     387             : #elif   EXP_ENABLED
     388          21 :         CHECK_ASSERTION(__LINE__, 0_IK <= roff, SK_"@setMatInit(): The condition `0_IK <= roff` must hold. roff = "//getStr(roff))
     389          21 :         CHECK_ASSERTION(__LINE__, 0_IK <= coff, SK_"@setMatInit(): The condition `0_IK <= coff` must hold. coff = "//getStr(coff))
     390          84 :         CHECK_ASSERTION(__LINE__, 0_IK <= nrow .and. nrow + roff <= size(mat, 1, IK), SK_"@setMatInit(): The condition `0_IK <= nrow .and. nrow + roff <= size(mat,1)` must hold. nrow, roff, size(mat,1) = "//getStr([nrow, roff, size(mat, 1, IK)]))
     391          84 :         CHECK_ASSERTION(__LINE__, 0_IK <= ncol .and. ncol + coff <= size(mat, 2, IK), SK_"@setMatInit(): The condition `0_IK <= ncol .and. ncol + coff <= size(mat,2)` must hold. ncol, coff, size(mat,2) = "//getStr([ncol, coff, size(mat, 2, IK)]))
     392             : #else
     393             : #error  "Unrecognized interface."
     394             : #endif
     395             : #if     SK_ENABLED
     396          36 :         CHECK_ASSERTION(__LINE__, len(vupp, IK) <= len(mat, IK), SK_"@setMatInit(): The condition `len(vupp) <= len(mat)` must hold. len(vupp), len(mat) = "//getStr([len(vupp, IK), len(mat, IK)]))
     397          36 :         CHECK_ASSERTION(__LINE__, len(vlow, IK) <= len(mat, IK), SK_"@setMatInit(): The condition `len(vlow) <= len(mat)` must hold. len(vlow), len(mat) = "//getStr([len(vlow, IK), len(mat, IK)]))
     398             : #endif
     399             :         !check_assertion(__LINE__, all(0_IK < [nrow, ncol]), SK_"@setMatInit(): The condition `all(0 < [nrow, ncol])` must hold. nrow, ncol = "//getStr([[nrow, ncol]]))
     400          21 :         if (nrow < 1_IK .or. ncol < 1_IK) return ! this must be here.
     401          21 :         if (present(doff)) then
     402          24 :             CHECK_ASSERTION(__LINE__, -nrow < doff .and. doff < ncol, SK_"@setMatInit(): The condition `-nrow < doff .and. doff < ncol` must hold. nrow, doff, ncol = "//getStr([nrow, doff, ncol]))
     403             :             doff_def = doff
     404           6 :             if (doff_def < 0_IK) then
     405           1 :                 lenDia = min(nrow + doff_def, ncol)
     406             :             else
     407           5 :                 lenDia = min(ncol - doff_def, nrow)
     408             :             end if
     409             :         else
     410             :             doff_def = 0_IK
     411          15 :             lenDia = min(ncol, nrow)
     412             :         endif
     413           6 :         if (doff_def > 0_IK) then
     414          17 :             do i = 1, doff_def
     415          69 :                 mat(1 : nrow, i) = vlow
     416             :             end do
     417          24 :             mat(2 : nrow , doff_def + 1_IK) = vlow
     418          18 :             do i = 2_IK, lenDia
     419          38 :                 mat(1 : i - 1, i + doff_def) = vupp
     420          42 :                 mat(i + 1 : nrow , i + doff_def) = vlow
     421             :             end do
     422           5 :             do i = lenDia + doff_def + 1, ncol
     423           5 :                 mat(1 : nrow, i) = vupp
     424             :             end do
     425             :         else
     426          96 :             do i = 1_IK, lenDia
     427         299 :                 mat(1 : i - doff_def - 1, i) = vupp
     428         423 :                 mat(i - doff_def + 1 : nrow , i) = vlow
     429             :             end do
     430          22 :             do i = lenDia + 1, ncol
     431          46 :                 mat(1 : nrow, i) = vupp
     432             :             end do
     433             :         end if
     434             : 
     435             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     436             : #elif   setMatInit_ENABLED && ULD_ENABLED
     437             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     438             : 
     439             :         integer(IK) :: i, lenDia, doff_def
     440             : #if     IMP_ENABLED
     441             :         integer(IK), parameter :: roff = 0_IK, coff = 0_IK
     442             :         integer(IK) :: nrow, ncol
     443           0 :         nrow = size(mat, 1, IK)
     444           0 :         ncol = size(mat, 2, IK)
     445             : #elif   EXP_ENABLED
     446        3347 :         CHECK_ASSERTION(__LINE__, 0_IK <= roff, SK_"@setMatInit(): The condition `0_IK <= roff` must hold. roff = "//getStr(roff))
     447        3347 :         CHECK_ASSERTION(__LINE__, 0_IK <= coff, SK_"@setMatInit(): The condition `0_IK <= coff` must hold. coff = "//getStr(coff))
     448       13388 :         CHECK_ASSERTION(__LINE__, 0_IK <= nrow .and. nrow + roff <= size(mat, 1, IK), SK_"@setMatInit(): The condition `0_IK <= nrow .and. nrow + roff <= size(mat,1)` must hold. nrow, roff, size(mat,1) = "//getStr([nrow, roff, size(mat, 1, IK)]))
     449       13388 :         CHECK_ASSERTION(__LINE__, 0_IK <= ncol .and. ncol + coff <= size(mat, 2, IK), SK_"@setMatInit(): The condition `0_IK <= ncol .and. ncol + coff <= size(mat,2)` must hold. ncol, coff, size(mat,2) = "//getStr([ncol, coff, size(mat, 2, IK)]))
     450             : #else
     451             : #error  "Unrecognized interface."
     452             : #endif
     453             : #if     SK_ENABLED
     454          66 :         CHECK_ASSERTION(__LINE__, len(vupp, IK) <= len(mat, IK), SK_"@setMatInit(): The condition `len(vupp) <= len(mat)` must hold. len(vupp), len(mat) = "//getStr([len(vupp, IK), len(mat, IK)]))
     455          66 :         CHECK_ASSERTION(__LINE__, len(vlow, IK) <= len(mat, IK), SK_"@setMatInit(): The condition `len(vlow) <= len(mat)` must hold. len(vlow), len(mat) = "//getStr([len(vlow, IK), len(mat, IK)]))
     456          66 :         CHECK_ASSERTION(__LINE__, len(vdia, IK) <= len(mat, IK), SK_"@setMatInit(): The condition `len(vdia) <= len(mat)` must hold. len(vdia), len(mat) = "//getStr([len(vdia, IK), len(mat, IK)]))
     457             : #endif
     458             :         !check_assertion(__LINE__, all(0_IK < [nrow, ncol]), SK_"@setMatInit(): The condition `all(0 < [nrow, ncol])` must hold. nrow, ncol = "//getStr([[nrow, ncol]]))
     459        3347 :         if (nrow < 1_IK .or. ncol < 1_IK) return ! this must be here.
     460        3347 :         if (present(doff)) then
     461          24 :             CHECK_ASSERTION(__LINE__, -nrow < doff .and. doff < ncol, SK_"@setMatInit(): The condition `-nrow < doff .and. doff < ncol` must hold. nrow, doff, ncol = "//getStr([nrow, doff, ncol]))
     462             :             doff_def = doff
     463             :         else
     464             :             doff_def = 0_IK
     465             :         endif
     466             : #if     D2000_ENABLED
     467             : #define GET_DIA(i) vdia
     468           5 :         if (doff_def < 0_IK) then
     469           0 :             lenDia = min(nrow + doff_def, ncol)
     470             :         else
     471        2540 :             lenDia = min(ncol - doff_def, nrow)
     472             :         end if
     473             : #elif   D2001_ENABLED
     474             : #define GET_DIA(i) vdia(i)
     475         807 :         lenDia = size(vdia, 1, IK)
     476        4841 :         CHECK_ASSERTION(__LINE__, lenDia == merge(min(nrow + doff_def, ncol), min(nrow, ncol - doff_def), doff_def < 0_IK), \
     477             :         SK_"@setMatInit(): The condition `size(vdia) == merge(min(nrow + doff, ncol), min(nrow, ncol - doff), doff < 0)` must hold. size(vdia), nrow, ncol, doff = "//\
     478             :         getStr([size(vdia, 1, IK), nrow, ncol, doff_def]))
     479             : #else
     480             : #error  "Unrecognized interface."
     481             : #endif
     482             :         !check_assertion(__LINE__, all(0_IK < [nrow, ncol]), SK_"@setMatInit(): The condition `all(0 < [nrow, ncol])` must hold. nrow, ncol = "//getStr([[nrow, ncol]]))
     483             :         if (nrow < 1_IK .or. ncol < 1_IK) return ! this must be here.
     484        3347 :         if (doff_def > 0_IK) then
     485          17 :             do i = 1, doff_def
     486          69 :                 mat(1 : nrow, i) = vlow
     487             :             end do
     488           5 :             i = doff_def + 1_IK
     489           4 :             mat(1 , i) = GET_DIA(1)
     490          24 :             mat(2 : nrow , i) = vlow
     491          18 :             do i = 2_IK, lenDia
     492          38 :                 mat(1 : i - 1, i + doff_def) = vupp
     493          11 :                 mat(i , i + doff_def) = GET_DIA(i)
     494          42 :                 mat(i + 1 : nrow , i + doff_def) = vlow
     495             :             end do
     496           5 :             do i = lenDia + doff_def + 1, ncol
     497           5 :                 mat(1 : nrow, i) = vupp
     498             :             end do
     499             :         else
     500       11756 :             do i = 1_IK, lenDia
     501       17777 :                 mat(1 : i - doff_def - 1, i) = vupp
     502          58 :                 mat(i - doff_def , i) = GET_DIA(i)
     503       21266 :                 mat(i - doff_def + 1 : nrow , i) = vlow
     504             :             end do
     505        3362 :             do i = lenDia + 1, ncol
     506        3419 :                 mat(1 : nrow, i) = vupp
     507             :             end do
     508             :         end if
     509             : #else
     510             :         !%%%%%%%%%%%%%%%%%%%%%%%%
     511             : #error  "Unrecognized interface."
     512             :         !%%%%%%%%%%%%%%%%%%%%%%%%
     513             : #endif
     514             : #undef  GET_DIA
     515             : #undef  NDIA

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