https://www.cdslab.org/paramonte/fortran/2
Current view: top level - main - pm_paramonte.F90 (source / functions) Hit Total Coverage
Test: ParaMonte 2.0.0 :: Serial Fortran - Code Coverage Report Lines: 36 36 100.0 %
Date: 2024-04-08 03:18:57 Functions: 1 1 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             : ! Define a quoting macro for better illustration in the output file.
      18             : !> \cond excluded
      19             : #ifdef  __GFORTRAN__
      20             : # define QSTART(X) "&
      21             : # define QEND(X) &X"
      22             : #else   /* default quoting macro */
      23             : # define GET_QUOTED(X) #X
      24             : # define QSTART(X) &
      25             : # define QEND(X) GET_QUOTED(X)
      26             : #endif
      27             : !> \endcond excluded
      28             : 
      29             : !>  \brief
      30             : !>  This module contains procedures and data that provide general information about the ParaMonte library, its interfaces, and its build.
      31             : !>
      32             : !>  \finmain
      33             : !>
      34             : !>  \author
      35             : !>  \AmirShahmoradi, Monday 00:01 AM, January 1, 2018, Institute for Computational Engineering and Sciences, University of Texas Austin
      36             : 
      37             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      38             : 
      39             : module pm_paramonte
      40             : 
      41             :     use pm_kind, only: SK, IK, LK
      42             :     use pm_str, only: NLC, UNDEFINED
      43             :     use iso_fortran_env, only: compiler_options, compiler_version
      44             : 
      45             :     implicit none
      46             : 
      47             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      48             : 
      49             :     !>  \brief
      50             :     !>  The scalar constant of type `character` of default kind \SK,
      51             :     !>  containing the web portal address for the ParaMonte library repository.<br>
      52             :     !>
      53             :     !>  \see
      54             :     !>  [paramonte](@ref pm_paramonte::paramonte)<br>
      55             :     !>  [PARAMONTE_WEB_DOC](@ref pm_paramonte::PARAMONTE_WEB_DOC)<br>
      56             :     !>  [PARAMONTE_WEB_ISSUES](@ref pm_paramonte::PARAMONTE_WEB_ISSUES)<br>
      57             :     !>  [PARAMONTE_WEB_REPOSITORY](@ref pm_paramonte::PARAMONTE_WEB_REPOSITORY)<br>
      58             :     !>  [PARAMONTE_COMPILER_OPTIONS](@ref pm_paramonte::PARAMONTE_COMPILER_OPTIONS)<br>
      59             :     !>  [PARAMONTE_COMPILER_VERSION](@ref pm_paramonte::PARAMONTE_COMPILER_VERSION)<br>
      60             :     !>  [paramonte_type](@ref pm_paramonte::paramonte_type)<br>
      61             :     !>  [getBaseName](@ref pm_sysPath::getBaseName)<br>
      62             :     !>  [getDirName](@ref pm_sysPath::getDirName)<br>
      63             :     !>
      64             :     !>  \finmain{PARAMONTE_WEB_REPOSITORY}
      65             :     !>
      66             :     !>  \author
      67             :     !>  \AmirShahmoradi, September 1, 2017, 11:35 PM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin
      68             :     character(*, SK), parameter :: PARAMONTE_WEB_REPOSITORY = SK_"https://github.com/cdslaborg/paramonte"
      69             : 
      70             :     !>  \brief
      71             :     !>  The scalar constant of type `character` of default kind \SK,
      72             :     !>  containing the web portal address for reporting the ParaMonte library issues.<br>
      73             :     !>
      74             :     !>  \see
      75             :     !>  [paramonte](@ref pm_paramonte::paramonte)<br>
      76             :     !>  [PARAMONTE_WEB_DOC](@ref pm_paramonte::PARAMONTE_WEB_DOC)<br>
      77             :     !>  [PARAMONTE_WEB_ISSUES](@ref pm_paramonte::PARAMONTE_WEB_ISSUES)<br>
      78             :     !>  [PARAMONTE_WEB_REPOSITORY](@ref pm_paramonte::PARAMONTE_WEB_REPOSITORY)<br>
      79             :     !>  [PARAMONTE_COMPILER_OPTIONS](@ref pm_paramonte::PARAMONTE_COMPILER_OPTIONS)<br>
      80             :     !>  [PARAMONTE_COMPILER_VERSION](@ref pm_paramonte::PARAMONTE_COMPILER_VERSION)<br>
      81             :     !>  [paramonte_type](@ref pm_paramonte::paramonte_type)<br>
      82             :     !>  [getBaseName](@ref pm_sysPath::getBaseName)<br>
      83             :     !>  [getDirName](@ref pm_sysPath::getDirName)<br>
      84             :     !>
      85             :     !>  \finmain{PARAMONTE_WEB_ISSUES}
      86             :     !>
      87             :     !>  \author
      88             :     !>  \AmirShahmoradi, September 1, 2017, 11:35 PM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin
      89             :     character(*, SK), parameter :: PARAMONTE_WEB_ISSUES = SK_"https://github.com/cdslaborg/paramonte/issues"
      90             : 
      91             :     !>  \brief
      92             :     !>  The scalar constant of type `character` of default kind \SK,
      93             :     !>  containing the web portal address for reporting the ParaMonte library issues.<br>
      94             :     !>
      95             :     !>  \see
      96             :     !>  [paramonte](@ref pm_paramonte::paramonte)<br>
      97             :     !>  [PARAMONTE_WEB_DOC](@ref pm_paramonte::PARAMONTE_WEB_DOC)<br>
      98             :     !>  [PARAMONTE_WEB_ISSUES](@ref pm_paramonte::PARAMONTE_WEB_ISSUES)<br>
      99             :     !>  [PARAMONTE_WEB_REPOSITORY](@ref pm_paramonte::PARAMONTE_WEB_REPOSITORY)<br>
     100             :     !>  [PARAMONTE_COMPILER_OPTIONS](@ref pm_paramonte::PARAMONTE_COMPILER_OPTIONS)<br>
     101             :     !>  [PARAMONTE_COMPILER_VERSION](@ref pm_paramonte::PARAMONTE_COMPILER_VERSION)<br>
     102             :     !>  [paramonte_type](@ref pm_paramonte::paramonte_type)<br>
     103             :     !>  [getBaseName](@ref pm_sysPath::getBaseName)<br>
     104             :     !>  [getDirName](@ref pm_sysPath::getDirName)<br>
     105             :     !>
     106             :     !>  \finmain{PARAMONTE_WEB_DOC}
     107             :     !>
     108             :     !>  \author
     109             :     !>  \AmirShahmoradi, September 1, 2017, 11:35 PM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin
     110             :     character(*, SK), parameter :: PARAMONTE_WEB_DOC = SK_"https://www.cdslab.org/paramonte/"
     111             : 
     112             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     113             : 
     114             :     !>  \brief
     115             :     !>  The scalar constant of type `character` of default kind \SK,
     116             :     !>  containing the version of the compiler with which the ParaMonte library has been built.<br>
     117             :     !>
     118             :     !>  \see
     119             :     !>  [paramonte](@ref pm_paramonte::paramonte)<br>
     120             :     !>  [PARAMONTE_WEB_DOC](@ref pm_paramonte::PARAMONTE_WEB_DOC)<br>
     121             :     !>  [PARAMONTE_WEB_ISSUES](@ref pm_paramonte::PARAMONTE_WEB_ISSUES)<br>
     122             :     !>  [PARAMONTE_WEB_REPOSITORY](@ref pm_paramonte::PARAMONTE_WEB_REPOSITORY)<br>
     123             :     !>  [PARAMONTE_COMPILER_OPTIONS](@ref pm_paramonte::PARAMONTE_COMPILER_OPTIONS)<br>
     124             :     !>  [PARAMONTE_COMPILER_VERSION](@ref pm_paramonte::PARAMONTE_COMPILER_VERSION)<br>
     125             :     !>  [paramonte_type](@ref pm_paramonte::paramonte_type)<br>
     126             :     !>  [getBaseName](@ref pm_sysPath::getBaseName)<br>
     127             :     !>  [getDirName](@ref pm_sysPath::getDirName)<br>
     128             :     !>
     129             :     !>  \finmain{PARAMONTE_COMPILER_VERSION}
     130             :     !>
     131             :     !>  \author
     132             :     !>  \AmirShahmoradi, September 1, 2017, 11:35 PM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin
     133             :     character(*, SK), parameter :: PARAMONTE_COMPILER_VERSION = compiler_version()
     134             : 
     135             :     !>  \brief
     136             :     !>  The scalar constant of type `character` of default kind \SK,
     137             :     !>  containing the compiler options used to build the ParaMonte library.<br>
     138             :     !>
     139             :     !>  \see
     140             :     !>  [paramonte](@ref pm_paramonte::paramonte)<br>
     141             :     !>  [PARAMONTE_WEB_DOC](@ref pm_paramonte::PARAMONTE_WEB_DOC)<br>
     142             :     !>  [PARAMONTE_WEB_ISSUES](@ref pm_paramonte::PARAMONTE_WEB_ISSUES)<br>
     143             :     !>  [PARAMONTE_WEB_REPOSITORY](@ref pm_paramonte::PARAMONTE_WEB_REPOSITORY)<br>
     144             :     !>  [PARAMONTE_COMPILER_OPTIONS](@ref pm_paramonte::PARAMONTE_COMPILER_OPTIONS)<br>
     145             :     !>  [PARAMONTE_COMPILER_VERSION](@ref pm_paramonte::PARAMONTE_COMPILER_VERSION)<br>
     146             :     !>  [paramonte_type](@ref pm_paramonte::paramonte_type)<br>
     147             :     !>  [getBaseName](@ref pm_sysPath::getBaseName)<br>
     148             :     !>  [getDirName](@ref pm_sysPath::getDirName)<br>
     149             :     !>
     150             :     !>  \finmain{PARAMONTE_COMPILER_OPTIONS}
     151             :     !>
     152             :     !>  \author
     153             :     !>  \AmirShahmoradi, September 1, 2017, 11:35 PM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin
     154             :     character(*, SK), parameter :: PARAMONTE_COMPILER_OPTIONS = compiler_options()
     155             : 
     156             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     157             : 
     158             :     !>  \brief
     159             :     !>  This a derived type devoid of any components or methods whose instantiated objects are used within the ParaMonte library
     160             :     !>  to signify the use of the ParaMonte style (vs. alternative approaches) in performing various actions within the library.<br>
     161             :     !>
     162             :     !>  \see
     163             :     !>  [paramonte](@ref pm_paramonte::paramonte)<br>
     164             :     !>  [PARAMONTE_WEB_DOC](@ref pm_paramonte::PARAMONTE_WEB_DOC)<br>
     165             :     !>  [PARAMONTE_WEB_ISSUES](@ref pm_paramonte::PARAMONTE_WEB_ISSUES)<br>
     166             :     !>  [PARAMONTE_WEB_REPOSITORY](@ref pm_paramonte::PARAMONTE_WEB_REPOSITORY)<br>
     167             :     !>  [PARAMONTE_COMPILER_OPTIONS](@ref pm_paramonte::PARAMONTE_COMPILER_OPTIONS)<br>
     168             :     !>  [PARAMONTE_COMPILER_VERSION](@ref pm_paramonte::PARAMONTE_COMPILER_VERSION)<br>
     169             :     !>  [paramonte_type](@ref pm_paramonte::paramonte_type)<br>
     170             :     !>  [getBaseName](@ref pm_sysPath::getBaseName)<br>
     171             :     !>  [getDirName](@ref pm_sysPath::getDirName)<br>
     172             :     !>
     173             :     !>  \finmain{paramonte_type}
     174             :     !>
     175             :     !>  \author
     176             :     !>  \AmirShahmoradi, September 1, 2017, 11:35 PM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin
     177             :     type :: paramonte_type
     178             :     end type
     179             : 
     180             :     !>  \brief
     181             :     !>  The scalar module constant of type [paramonte_type](@ref pm_paramonte::paramonte_type) that can be used within the ParaMonte library
     182             :     !>  to signify the use of the ParaMonte style (vs. alternative approaches) in performing various actions within the library.<br>
     183             :     !>
     184             :     !>  \see
     185             :     !>  [paramonte](@ref pm_paramonte::paramonte)<br>
     186             :     !>  [PARAMONTE_WEB_DOC](@ref pm_paramonte::PARAMONTE_WEB_DOC)<br>
     187             :     !>  [PARAMONTE_WEB_ISSUES](@ref pm_paramonte::PARAMONTE_WEB_ISSUES)<br>
     188             :     !>  [PARAMONTE_WEB_REPOSITORY](@ref pm_paramonte::PARAMONTE_WEB_REPOSITORY)<br>
     189             :     !>  [PARAMONTE_COMPILER_OPTIONS](@ref pm_paramonte::PARAMONTE_COMPILER_OPTIONS)<br>
     190             :     !>  [PARAMONTE_COMPILER_VERSION](@ref pm_paramonte::PARAMONTE_COMPILER_VERSION)<br>
     191             :     !>  [paramonte_type](@ref pm_paramonte::paramonte_type)<br>
     192             :     !>  [getBaseName](@ref pm_sysPath::getBaseName)<br>
     193             :     !>  [getDirName](@ref pm_sysPath::getDirName)<br>
     194             :     !>
     195             :     !>  \finmain{paramonte}
     196             :     !>
     197             :     !>  \author
     198             :     !>  \AmirShahmoradi, September 1, 2017, 11:35 PM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin
     199             :     type(paramonte_type), parameter :: paramonte = paramonte_type()
     200             : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
     201             :     !DIR$ ATTRIBUTES DLLEXPORT :: paramonte
     202             : #endif
     203             : 
     204             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     205             : 
     206             :     !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     207             :     ! The ParaMonte programming language interface.
     208             :     !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     209             : 
     210             : !>  \cond excluded
     211             : #if     C_ENABLED
     212             : 
     213             : #define PROGLANG c
     214             : #define PROGNAME"C"
     215             : 
     216             : #elif   COBOL_ENABLED
     217             : 
     218             : #define PROGLANG cobol
     219             : #define PROGNAME"COBOL"
     220             : 
     221             : #elif   CPP_ENABLED
     222             : 
     223             : #define PROGLANG cpp
     224             : #define PROGNAME"C++"
     225             : 
     226             : #elif   CSHARP_ENABLED
     227             : 
     228             : #define PROGLANG csharp
     229             : #define PROGNAME"C#"
     230             : 
     231             : #elif   GO_ENABLED
     232             : 
     233             : #define PROGLANG go
     234             : #define PROGNAME"Go"
     235             : 
     236             : #elif   FORTRAN_ENABLED
     237             : 
     238             : #define PROGLANG fortran
     239             : #define PROGNAME"Fortran"
     240             : 
     241             : #elif   JAVA_ENABLED
     242             : 
     243             : #define PROGLANG java
     244             : #define PROGNAME"Java"
     245             : 
     246             : #elif   JAVASCRIPT_ENABLED
     247             : 
     248             : #define PROGLANG javascript
     249             : #define PROGNAME"JavaScript"
     250             : 
     251             : #elif   JULIA_ENABLED
     252             : 
     253             : #define PROGLANG julia
     254             : #define PROGNAME"Julia"
     255             : 
     256             : #elif   MATLAB_ENABLED
     257             : 
     258             : #define PROGLANG matlab
     259             : #define PROGNAME"MATLAB"
     260             : 
     261             : #elif   MATTHEMATICA_ENABLED
     262             : 
     263             : #define PROGLANG mathematica
     264             : #define PROGNAME"Wolfram Mathematica"
     265             : 
     266             : #elif   PYTHON_ENABLED
     267             : 
     268             : #define PROGLANG python
     269             : #define PROGNAME"Python"
     270             : 
     271             : #elif   R_ENABLED
     272             : 
     273             : #define PROGLANG r
     274             : #define PROGNAME"R"
     275             : 
     276             : #elif   RUST_ENABLED
     277             : 
     278             : #define PROGLANG rust
     279             : #define PROGNAME"Rust"
     280             : 
     281             : #elif   SAS_ENABLED
     282             : 
     283             : #define PROGLANG sas
     284             : #define PROGNAME"SAS"
     285             : 
     286             : #elif   SWIFT_ENABLED
     287             : 
     288             : #define PROGLANG swift
     289             : #define PROGNAME"Swift"
     290             : 
     291             : #else
     292             : #error  "Unrecognized interface."
     293             : #endif
     294             : !>  \endcond excluded
     295             : 
     296             :     !>  \brief
     297             :     !>  The scalar constant of type `character` of default kind \SK containing the full generic name
     298             :     !>  of the programming language environment for which the ParaMonte library has been configured.
     299             :     !>
     300             :     !>  \see
     301             :     !>  [envis](@ref pm_paramonte::envis)<br>
     302             :     !>  [envis_type](@ref pm_paramonte::envis_type)<br>
     303             :     !>  [envname](@ref pm_paramonte::envname)<br>
     304             :     !>
     305             :     !>  \finmain{filext_type}
     306             :     !>
     307             :     !>  \author
     308             :     !>  \AmirShahmoradi, September 1, 2017, 11:35 PM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin
     309             :     character(*, SK), parameter :: envname = PROGNAME//" programming language"
     310             : 
     311             :     !>  \brief
     312             :     !>  This is the derived type for generating objects containing scalar components of type `logical` named after
     313             :     !>  different programming language environments that are currently or will be supported for access to the ParaMonte library.
     314             :     !>
     315             :     !>  \details
     316             :     !>  This derived type is of minimal usage outside the ParaMonte library internal routines.<br>
     317             :     !>  If needed, the constant object instances of this derived type [envis](@ref pm_paramonte::envis) can be used.<br>
     318             :     !>
     319             :     !>  \see
     320             :     !>  [envis](@ref pm_paramonte::envis)<br>
     321             :     !>  [envis_type](@ref pm_paramonte::envis_type)<br>
     322             :     !>  [envname](@ref pm_paramonte::envname)<br>
     323             :     !>
     324             :     !>  \finmain{envis_type}
     325             :     !>
     326             :     !>  \author
     327             :     !>  \AmirShahmoradi, September 1, 2017, 11:35 PM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin
     328             :     type :: envis_type
     329             :         logical(LK) :: c            = .false._LK
     330             :         logical(LK) :: cobol        = .false._LK
     331             :         logical(LK) :: cpp          = .false._LK
     332             :         logical(LK) :: csharp       = .false._LK
     333             :         logical(LK) :: go           = .false._LK
     334             :         logical(LK) :: fortran      = .false._LK
     335             :         logical(LK) :: java         = .false._LK
     336             :         logical(LK) :: javascript   = .false._LK
     337             :         logical(LK) :: julia        = .false._LK
     338             :         logical(LK) :: matlab       = .false._LK
     339             :         logical(LK) :: mathematica  = .false._LK
     340             :         logical(LK) :: python       = .false._LK
     341             :         logical(LK) :: r            = .false._LK
     342             :         logical(LK) :: rust         = .false._LK
     343             :         logical(LK) :: sas          = .false._LK
     344             :         logical(LK) :: swift        = .false._LK
     345             :     end type
     346             : 
     347             :     !>  \brief
     348             :     !>  The scalar constant object of type [envis_type](@ref pm_paramonte::envis_type) whose components
     349             :     !>  are all `.false.` except the environment for which the ParaMonte library is currently configured.
     350             :     !>
     351             :     !>  \see
     352             :     !>  [envis](@ref pm_paramonte::envis)<br>
     353             :     !>  [envis_type](@ref pm_paramonte::envis_type)<br>
     354             :     !>  [envname](@ref pm_paramonte::envname)<br>
     355             :     !>
     356             :     !>  \finmain{envis}
     357             :     !>
     358             :     !>  \author
     359             :     !>  \AmirShahmoradi, September 1, 2017, 11:35 PM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin
     360             :     type(envis_type), parameter :: envis = envis_type(PROGLANG = .true._LK)
     361             : 
     362             : !>  \cond excluded
     363             : #undef PROGLANG
     364             : #undef PROGNAME
     365             : !>  \endcond excluded
     366             : 
     367             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     368             : 
     369             :     !>  \brief
     370             :     !>  The `public` scalar `character` constant of default kind \SK containing the ParaMonte library build date.
     371             :     !>
     372             :     !>  \details
     373             :     !>  The generated date is,
     374             :     !>  <ol>
     375             :     !>      <li>    the value of the predefined compiler preprocessor macro `__TIMESTAMP__` if the compiler suite is Intel or GNU.<br>
     376             :     !>      <li>    the value of CMake build generator `TIMESTAMP` with format `"%a %b %d %H:%M:%S %Y"` if the build is configured with CMake.<br>
     377             :     !>  </ol>
     378             :     !>
     379             :     !>  \interface{PARAMONTE_BUILD_DATE}
     380             :     !>  \code{.F90}
     381             :     !>
     382             :     !>      use pm_paramonte, only: PARAMONTE_BUILD_DATE
     383             :     !>
     384             :     !>      print *, PARAMONTE_BUILD_DATE
     385             :     !>
     386             :     !>  \endcode
     387             :     !>
     388             :     !>  \see
     389             :     !>  [PARAMONTE_SPLASH](@ref pm_paramonte::PARAMONTE_SPLASH)<br>
     390             :     !>  [PARAMONTE_VERSION](@ref pm_paramonte::PARAMONTE_VERSION)<br>
     391             :     !>  [PARAMONTE_BUILD_DATE](@ref pm_paramonte::PARAMONTE_BUILD_DATE)<br>
     392             :     !>  [getParaMonteSplash](@ref pm_paramonte::getParaMonteSplash)<br>
     393             :     !>
     394             :     !>  \example{PARAMONTE_BUILD_DATE}
     395             :     !>  \include{lineno} example/pm_paramonte/PARAMONTE_BUILD_DATE/main.F90
     396             :     !>  \compilef{PARAMONTE_BUILD_DATE}
     397             :     !>  \output{PARAMONTE_BUILD_DATE}
     398             :     !>  \include{lineno} example/pm_paramonte/PARAMONTE_BUILD_DATE/main.out.F90
     399             :     !>
     400             :     !>  \finmain{PARAMONTE_BUILD_DATE}
     401             :     !>
     402             :     !>  \author
     403             :     !>  \AmirShahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
     404             : #if INTEL_ENABLED || GNU_ENABLED
     405             :     character(*, SK), parameter :: PARAMONTE_BUILD_DATE = __TIMESTAMP__
     406             :     !>  \cond excluded
     407             : #elif defined PARAMONTE_BUILD_TIMESTAMP
     408             :     character(*, SK), parameter :: PARAMONTE_BUILD_DATE = PARAMONTE_BUILD_TIMESTAMP
     409             : #else
     410             :     character(*, SK), parameter :: PARAMONTE_BUILD_DATE = UNDEFINED
     411             : #endif
     412             :     !>  \endcond excluded
     413             : 
     414             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     415             : 
     416             :     !>  \brief
     417             :     !>  The `public` scalar `character` constant of default kind \SK containing the ParaMonte library version.
     418             :     !>
     419             :     !>  \details
     420             :     !>  The generated version is the value reported in the
     421             :     !>  [VERSION.md](https://raw.githubusercontent.com/cdslaborg/paramonte/main/src/fortran/VERSION.md)
     422             :     !>  in the root directory of the ParaMonte library git repository.
     423             :     !>
     424             :     !>  \interface{PARAMONTE_VERSION}
     425             :     !>  \code{.F90}
     426             :     !>
     427             :     !>      use pm_paramonte, only: PARAMONTE_VERSION
     428             :     !>
     429             :     !>      print *, PARAMONTE_VERSION
     430             :     !>
     431             :     !>  \endcode
     432             :     !>
     433             :     !>  \see
     434             :     !>  [PARAMONTE_SPLASH](@ref pm_paramonte::PARAMONTE_SPLASH)<br>
     435             :     !>  [PARAMONTE_VERSION](@ref pm_paramonte::PARAMONTE_VERSION)<br>
     436             :     !>  [PARAMONTE_BUILD_DATE](@ref pm_paramonte::PARAMONTE_BUILD_DATE)<br>
     437             :     !>  [getParaMonteSplash](@ref pm_paramonte::getParaMonteSplash)<br>
     438             :     !>
     439             :     !>  \example{PARAMONTE_VERSION}
     440             :     !>  \include{lineno} example/pm_paramonte/PARAMONTE_VERSION/main.F90
     441             :     !>  \compilef{PARAMONTE_VERSION}
     442             :     !>  \output{PARAMONTE_VERSION}
     443             :     !>  \include{lineno} example/pm_paramonte/PARAMONTE_VERSION/main.out.F90
     444             :     !>
     445             :     !>  \finmain{PARAMONTE_VERSION}
     446             :     !>
     447             :     !>  \author
     448             :     !>  \AmirShahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
     449             : #if defined __PARAMONTE_VERSION__
     450             :     character(*, SK), parameter :: PARAMONTE_VERSION = __PARAMONTE_VERSION__
     451             :     !>  \cond excluded
     452             : #else
     453             :     ! WARNING: The following ParaMonte library version tag will be taken from the declaration
     454             :     ! WARNING: that is generated by the preprocessing build script of the ParaMonte library.
     455             :     ! WARNING: This is superior to the above method of using the compiler Preprocessor
     456             :     ! WARNING: since CMAKE triggers a complete lengthy rebuild of the library when the
     457             :     ! WARNING: only the library version has changed.
     458             : #include "pm_paramonte@version.inc.F90"
     459             : #endif
     460             :     !>  \endcond excluded
     461             : 
     462             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     463             : 
     464             :     !>  \brief
     465             :     !>  The `public` scalar `character` constant of default kind \SK containing the ParaMonte splash information.
     466             :     !>
     467             :     !>  \details
     468             :     !>  The splash information includes the following items,
     469             :     !>  <ol>
     470             :     !>      <li>    The library name.
     471             :     !>      <li>    The library version.
     472             :     !>      <li>    The library build date.
     473             :     !>      <li>    The developing organization name.
     474             :     !>      <li>    The developing organization affiliations.
     475             :     !>      <li>    The project PI contact information.
     476             :     !>      <li>    The project GitHub web address.
     477             :     !>  </ol>
     478             :     !>
     479             :     !>  \interface{PARAMONTE_SPLASH}
     480             :     !>  \code{.F90}
     481             :     !>
     482             :     !>      use pm_paramonte, only: PARAMONTE_SPLASH
     483             :     !>
     484             :     !>      print *, PARAMONTE_SPLASH
     485             :     !>
     486             :     !>  \endcode
     487             :     !>
     488             :     !>  \see
     489             :     !>  [PARAMONTE_SPLASH](@ref pm_paramonte::PARAMONTE_SPLASH)<br>
     490             :     !>  [PARAMONTE_VERSION](@ref pm_paramonte::PARAMONTE_VERSION)<br>
     491             :     !>  [PARAMONTE_BUILD_DATE](@ref pm_paramonte::PARAMONTE_BUILD_DATE)<br>
     492             :     !>  [getParaMonteSplash](@ref pm_paramonte::getParaMonteSplash)<br>
     493             :     !>
     494             :     !>  \example{PARAMONTE_SPLASH}
     495             :     !>  \include{lineno} example/pm_paramonte/PARAMONTE_SPLASH/main.F90
     496             :     !>  \compilef{PARAMONTE_SPLASH}
     497             :     !>  \output{PARAMONTE_SPLASH}
     498             :     !>  \include{lineno} example/pm_paramonte/PARAMONTE_SPLASH/main.out.F90
     499             :     !>
     500             :     !>  \finmain{PARAMONTE_SPLASH}
     501             :     !>
     502             :     !>  \author
     503             :     !>  \AmirShahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
     504             :     character(*, SK), parameter :: PARAMONTE_SPLASH =   NLC//NLC// &
     505             :                                                         SK_"ParaMonte"//NLC// &
     506             :                                                         SK_"Parallel Monte Carlo and"//NLC// &
     507             :                                                         SK_"Machine Learning Library"//NLC// &
     508             :                                                         NLC// &
     509             :                                                         SK_"Version "//PARAMONTE_VERSION//NLC// &
     510             :                                                         NLC// &
     511             :                                                         SK_"Build: "//PARAMONTE_BUILD_DATE//NLC// &
     512             :                                                         NLC// &
     513             :                                                         SK_"Developed by"//NLC// &
     514             :                                                         NLC// &
     515             :                                                         SK_"The computational Data Science Lab"//NLC// &
     516             :                                                         NLC// &
     517             :                                                         SK_"at"//NLC// &
     518             :                                                         NLC// &
     519             :                                                         SK_"Department of Physics"//NLC// &
     520             :                                                         SK_"Data Science Program, College of Science"//NLC// &
     521             :                                                         SK_"The University of Texas at Arlington"//NLC// &
     522             :                                                         NLC// &
     523             :                                                         SK_"Multiscale Modeling Group"//NLC// &
     524             :                                                         SK_"Center for Computational Oncology"//NLC// &
     525             :                                                         SK_"Oden Institute for Computational Engineering and Sciences"//NLC// &
     526             :                                                         SK_"Department of Aerospace Engineering and Engineering Mechanics"//NLC// &
     527             :                                                         SK_"Department of Neurology, Dell-Seton Medical School"//NLC// &
     528             :                                                         SK_"Department of Biomedical Engineering"//NLC// &
     529             :                                                         SK_"The University of Texas at Austin"//NLC// &
     530             :                                                         NLC// &
     531             :                                                         SK_"For questions and further information, please contact the PI:"//NLC// &
     532             :                                                         NLC// &
     533             :                                                         SK_"Amir Shahmoradi"//NLC// &
     534             :                                                         !NLC// &
     535             :                                                         SK_"shahmoradi@utexas.edu"//NLC// &
     536             :                                                         !SK_"amir.shahmoradi@uta.edu"//NLC// &
     537             :                                                         !SK_"ashahmoradi@gmail.com"//NLC// &
     538             :                                                         !SK_"amir@physics.utexas.edu"//NLC// &
     539             :                                                         !SK_"amir@austin.utexas.edu"//NLC// &
     540             :                                                         !SK_"amir@ph.utexas.edu"//NLC// &
     541             :                                                         !NLC// &
     542             :                                                         !NLC// &
     543             :                                                         !SK_"https://www.cdslab.org/pm"//NLC// &
     544             :                                                         SK_"cdslab.org/pm"//NLC// &
     545             :                                                         NLC//NLC
     546             : 
     547             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     548             : 
     549             : contains
     550             : 
     551             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     552             : 
     553             :     !>  \brief
     554             :     !>  Generate and return an `allocatable` scalar `character` of default kind \SK
     555             :     !>  containing the ParaMonte library splash information centered with the requested line width, margins, and fillings.
     556             :     !>
     557             :     !>  \details
     558             :     !>  See the documentation of [PARAMONTE_SPLASH](@ref pm_paramonte::PARAMONTE_SPLASH) for further splash information details.<br>
     559             :     !>  The output splash begins and ends with a newline character returned by the intrinsic `new_line("a")` for the requested `character` kind.<br>
     560             :     !>
     561             :     !>  \param[in]      width   :   The input scalar `integer` of default kind \IK, containing the width of the splash **including** the lengths of the left and right margins.<br>
     562             :     !>                              (**optional**, default = `132`.)
     563             :     !>  \param[in]      lwsize  :   The input scalar `integer` of default kind \IK representing the width of the left wrapper margin of each line of the output `splash`.<br>
     564             :     !>                              (**optional**, default = `4`.)
     565             :     !>  \param[in]      twsize  :   The input scalar `integer` of default kind \IK representing the width of the top wrapper margin of the output `splash`.<br>
     566             :     !>                              (**optional**, default = `2`.)
     567             :     !>  \param[in]      rwsize  :   The input scalar `integer` of default kind \IK representing the width of the right wrapper margin of each line of the output `splash`.<br>
     568             :     !>                              (**optional**, default = `4`.)
     569             :     !>  \param[in]      bwsize  :   The input scalar `integer` of default kind \IK representing the width of the bottom wrapper margin of the output `splash`.<br>
     570             :     !>                              (**optional**, default = `2`.)
     571             :     !>  \param[in]      fill    :   The input scalar of the same type and kind as the output `slpash`, of length type parameter `1`,
     572             :     !>                              containing the value to fill the new elements (if any, excluding the margins) surrounding the splash text in each line of the output splash.<br>
     573             :     !>                              (**optional**, default = `" "`, the whitespace character.)
     574             :     !>  \param[in]      lwfill  :   The input scalar of the same type and kind as the output `slpash`, of length type parameter `1`,
     575             :     !>                              containing the value to fill the left wrapper margin (if any) of each line of the output `splash`.<br>
     576             :     !>                              (**optional**, default = `"%"`.)
     577             :     !>  \param[in]      twfill  :   The input scalar of the same type and kind as the output `slpash`, of length type parameter `1`,
     578             :     !>                              containing the value to fill the top wrapper margin (if any) of the output `splash`.<br>
     579             :     !>                              (**optional**, default = `"%"`.)
     580             :     !>  \param[in]      rwfill  :   The input scalar of the same type and kind as the output `slpash`, of length type parameter `1`,
     581             :     !>                              containing the value to fill the right wrapper margin (if any) of each line of the output `splash`.<br>
     582             :     !>                              (**optional**, default = `"%"`.)
     583             :     !>  \param[in]      bwfill  :   The input scalar of the same type and kind as the output `slpash`, of length type parameter `1`,
     584             :     !>                              containing the value to fill the bottom wrapper margin (if any) of the output `splash`.<br>
     585             :     !>                              (**optional**, default = `"%"`.)
     586             :     !>
     587             :     !>  \return
     588             :     !>  `splash`                :   The output `allocatable` scalar of type `character` of default kind \SK, containing the ParaMonte splash.<br>
     589             :     !>                              Lines within the `splash` are separated via the newline character returned by the intrinsic `new_line(SK_"a")`.<br>
     590             :     !>
     591             :     !>  \interface{getParaMonteSplash}
     592             :     !>  \code{.F90}
     593             :     !>
     594             :     !>      use pm_paramonte, only: getParaMonteSplash
     595             :     !>      character(:, SK), allocatable :: splash
     596             :     !>
     597             :     !>      splash = getParaMonteSplash(width = width, lwsize = lwsize, twsize = twsize, rwsize = rwsize, bwsize = bwsize, fill = fill, lwfill = lwfill, twfill = twfill, rwfill = rwfill, bwfill = bwfill)
     598             :     !>
     599             :     !>  \endcode
     600             :     !>
     601             :     !>  \warnpure
     602             :     !>
     603             :     !>  \see
     604             :     !>  [PARAMONTE_SPLASH](@ref pm_paramonte::PARAMONTE_SPLASH)<br>
     605             :     !>  [PARAMONTE_VERSION](@ref pm_paramonte::PARAMONTE_VERSION)<br>
     606             :     !>  [PARAMONTE_BUILD_DATE](@ref pm_paramonte::PARAMONTE_BUILD_DATE)<br>
     607             :     !>  [getParaMonteSplash](@ref pm_paramonte::getParaMonteSplash)<br>
     608             :     !>
     609             :     !>  \example{getParaMonteSplash}
     610             :     !>  \include{lineno} example/pm_paramonte/getParaMonteSplash/main.F90
     611             :     !>  \compilef{getParaMonteSplash}
     612             :     !>  \output{getParaMonteSplash}
     613             :     !>  \include{lineno} example/pm_paramonte/getParaMonteSplash/main.out.F90
     614             :     !>
     615             :     !>  \finmain{getParaMonteSplash}
     616             :     !>
     617             :     !>  \author
     618             :     !>  \AmirShahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
     619          30 :     PURE function getParaMonteSplash(width, lwsize, twsize, rwsize, bwsize, fill, lwfill, twfill, rwfill, bwfill) result(splash)
     620             : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
     621             :         !DEC$ ATTRIBUTES DLLEXPORT :: getParaMonteSplash
     622             : #endif
     623             :         use pm_option, only: getOption
     624             :         use pm_arrayFind, only: setLoc
     625             :         use pm_arrayResize, only: setResized
     626             :         use pm_arrayCenter, only: setCentered
     627             :         integer(IK)     , intent(in)    , optional  :: width, lwsize, twsize, rwsize, bwsize
     628             :         character(1, SK), intent(in)    , optional  :: fill , lwfill, twfill, rwfill, bwfill
     629             :         character(:, SK), allocatable               :: splash
     630             :         integer(IK)     , allocatable               :: locNLC(:)
     631             :         character(1, SK)                            :: fill_def , lwfill_def, twfill_def, rwfill_def, bwfill_def
     632             :         integer(IK)                                 :: splashLen, iline, ibeg, iend, sbeg, lenNLC, linewidth, lenLocNLC
     633             :         integer(IK)                                 :: lwsize_def, twsize_def, rwsize_def, bwsize_def
     634          30 :         lwsize_def = getOption(4_IK, lwsize)
     635          30 :         twsize_def = getOption(2_IK, twsize)
     636          30 :         rwsize_def = getOption(4_IK, rwsize)
     637          30 :         bwsize_def = getOption(2_IK, bwsize)
     638          86 :         lwfill_def = getOption(SK_"%", lwfill)
     639          88 :         twfill_def = getOption(SK_"%", twfill)
     640          86 :         rwfill_def = getOption(SK_"%", rwfill)
     641          88 :         bwfill_def = getOption(SK_"%", bwfill)
     642          86 :         fill_def = getOption(SK_" ", fill)
     643          30 :         lenNLC = len(NLC, IK)
     644          30 :         call setResized(locNLC, 47_IK)
     645          30 :         call setLoc(locNLC, lenLocNLC, PARAMONTE_SPLASH, NLC, blindness = lenNLC)
     646             :         ! Ensure the width is minimally set to the maximum line length (~61) of the splash message.
     647             :         !print *, getOption(132_IK, width)
     648        1050 :         linewidth = max(maxval(locNLC(2 : lenLocNLC) - locNLC(1 : lenLocNLC - 1)) - lenNLC, getOption(132_IK, width))! + lwsize_def + rwsize_def
     649          30 :         splashLen = lenNLC + (twsize_def + lenLocNLC + bwsize_def) * (linewidth + lenNLC)
     650          30 :         allocate(character(splashLen, SK) :: splash)
     651          30 :         splash(1 : lenNLC) = NLC
     652          30 :         ibeg = lenNLC + 1_IK
     653          89 :         do iline = 1, twsize_def
     654          59 :             iend = ibeg + linewidth - 1_IK
     655        7475 :             splash(ibeg : iend) = repeat(twfill_def, linewidth)
     656          59 :             ibeg = iend + lenNLC + 1
     657          89 :             splash(iend + 1 : ibeg - 1) = NLC
     658             :         end do
     659          30 :         if (lenLocNLC > 0_IK) then
     660             :             sbeg = 1_IK ! splash line beginning.
     661        1080 :             do iline = 1, lenLocNLC
     662        1050 :                 iend = ibeg + linewidth - 1_IK
     663             :                 call setCentered( splash(ibeg : iend) & ! LCOV_EXCL_LINE
     664             :                                 , PARAMONTE_SPLASH(sbeg : locNLC(iline) - 1) & ! LCOV_EXCL_LINE
     665             :                                 , lmsize = lwsize_def & ! LCOV_EXCL_LINE
     666             :                                 , rmsize = rwsize_def & ! LCOV_EXCL_LINE
     667             :                                 , lmfill = lwfill_def & ! LCOV_EXCL_LINE
     668             :                                 , rmfill = rwfill_def & ! LCOV_EXCL_LINE
     669             :                                 , fill = fill_def & ! LCOV_EXCL_LINE
     670        1050 :                                 )
     671        1050 :                 sbeg = locNLC(iline) + lenNLC
     672        1050 :                 ibeg = iend + lenNLC + 1_IK
     673        1080 :                 splash(iend + 1 : ibeg - 1) = NLC
     674             :             end do
     675             :         else
     676             :             call setCentered( splash & ! LCOV_EXCL_LINE
     677             :                             , PARAMONTE_SPLASH & ! LCOV_EXCL_LINE
     678             :                             , lmsize = lwsize_def & ! LCOV_EXCL_LINE
     679             :                             , rmsize = rwsize_def & ! LCOV_EXCL_LINE
     680             :                             , lmfill = lwfill_def & ! LCOV_EXCL_LINE
     681             :                             , rmfill = rwfill_def & ! LCOV_EXCL_LINE
     682             :                             , fill = fill_def & ! LCOV_EXCL_LINE
     683             :                             ) ! LCOV_EXCL_LINE
     684             :         end if
     685          91 :         do iline = 1, bwsize_def
     686          61 :             iend = ibeg + linewidth - 1_IK
     687        7657 :             splash(ibeg : iend) = repeat(bwfill_def, linewidth)
     688          61 :             ibeg = iend + lenNLC + 1
     689          91 :             splash(iend + 1 : ibeg - 1) = NLC
     690             :         end do
     691          30 :     end function
     692             : 
     693             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     694             : 
     695             : end module pm_paramonte ! LCOV_EXCL_LINE

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