https://www.cdslab.org/paramonte/fortran/2
Current view: top level - main - pm_test.F90 (source / functions) Hit Total Coverage
Test: ParaMonte 2.0.0 :: Serial Fortran - Code Coverage Report Lines: 103 120 85.8 %
Date: 2024-04-08 03:18:57 Functions: 6 6 100.0 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
       2             : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
       3             : !!!!                                                                                                                            !!!!
       4             : !!!!    ParaMonte: Parallel Monte Carlo and Machine Learning Library.                                                           !!!!
       5             : !!!!                                                                                                                            !!!!
       6             : !!!!    Copyright (C) 2012-present, The Computational Data Science Lab                                                          !!!!
       7             : !!!!                                                                                                                            !!!!
       8             : !!!!    This file is part of the ParaMonte library.                                                                             !!!!
       9             : !!!!                                                                                                                            !!!!
      10             : !!!!    LICENSE                                                                                                                 !!!!
      11             : !!!!                                                                                                                            !!!!
      12             : !!!!       https://github.com/cdslaborg/paramonte/blob/main/LICENSE.md                                                          !!!!
      13             : !!!!                                                                                                                            !!!!
      14             : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
      15             : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
      16             : 
      17             : !>  \brief
      18             : !>  This module contains a simple unit-testing framework for the Fortran libraries, including the ParaMonte library.
      19             : !>
      20             : !>  \details
      21             : !>  Tests within the ParaMonte testing framework are based on modules.<br>
      22             : !>  Each module `pm_*` has a corresponding test module `test_pm_*` that implements tests of all objects and routines of the target module.<br>
      23             : !>  Each test module must contain a public `subroutine` named `setTest()` with the following interface,
      24             : !>  \code{.F90}
      25             : !>      subroutine setTest()
      26             : !>      end subroutine
      27             : !>  \endcode
      28             : !>  within which all test are implemented as desired.<br>
      29             : !>  These module subroutines are subsequently called in the ParaMonte `main.F90` file.<br>
      30             : !>  See the documentation of [test_type](@ref pm_test::test_type) for further details on the testing approach and assertion verification.<br>
      31             : !>
      32             : !>  \see
      33             : !>  [pm_err](@ref pm_err)<br>
      34             : !>
      35             : !>  \finmain
      36             : !>
      37             : !>  \author
      38             : !>  \AmirShahmoradi, September 1, 2017, 12:00 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin
      39             : 
      40             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      41             : 
      42             : module pm_test
      43             : 
      44             :     use pm_io, only: display_type
      45             :     use pm_kind, only: SK, IK, LK, RKC => RK
      46             :     use, intrinsic  :: iso_fortran_env, only: output_unit
      47             :     use pm_strANSI, only: bright, fcyan, fgreen, fyellow, fred, underlined, reset
      48             :     use pm_parallelism, only: image_type
      49             :     use pm_arrayResize, only: setResized
      50             :     use pm_container, only: css_type
      51             :     use pm_timer, only: timer_type
      52             :     use pm_val2str, only: getStr
      53             :     use pm_err, only: err_type
      54             : 
      55             :     implicit none
      56             : 
      57             :     private
      58             :     public :: SK, IK, LK, test_type, setSummary
      59             : 
      60             :     !>  \cond excluded
      61             :     character(*, SK)        , parameter     :: NLC = new_line(SK_"a")
      62             :     character(*, SK)        , parameter     :: MODULE_NAME = SK_"pm_test"
      63             :    !integer(IK)             , parameter     :: TEST_COUNTER_LEN = 8_IK
      64             :     integer(IK)             , parameter     :: TIME_FIELD_LEN = 9_IK
      65             :     integer(IK)                             :: mv_testCounterOld
      66             :     integer(IK)                             :: mv_testCounter
      67             :     integer(IK)                             :: mv_nfail
      68             :     integer(IK)                             :: mv_npass
      69             :     type(css_type)          , allocatable   :: mv_failedTestFuncName(:)
      70             :     character(:, SK)        , allocatable   :: mc_passedString
      71             :     character(:, SK)        , allocatable   :: mc_failedString
      72             :     logical(LK)                             :: mv_uninit = .true._LK
      73             :     class(timer_type)       , allocatable   :: mv_timer !<  \public This timer is exclusively used to time tests at the global level.
      74             :     type(image_type)                        :: mv_image
      75             :     !>  \endcond excluded
      76             : 
      77             :     !>  \brief
      78             :     !>  This is the module `private` derived type for constructing objects that contain the input and output directory paths for tests.
      79             :     !>
      80             :     !>  \details
      81             :     !>  This derived type is meant to be used internally within the parent module.<br>
      82             :     !>
      83             :     !>  \see
      84             :     !>  [dir_type](@ref pm_test::dir_type)<br>
      85             :     !>  [file_type](@ref pm_test::file_type)<br>
      86             :     !>  [test_type](@ref pm_test::test_type)<br>
      87             :     !>
      88             :     !>  \finmain{dir_type}
      89             :     !>
      90             :     !>  \author
      91             :     !>  \AmirShahmoradi, September 1, 2017, 12:00 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin
      92             :     type :: dir_type
      93             :         character(:, SK)    , allocatable   :: inp !<  \public The scalar `allocatable` of type `character` of default kind \SK containing the path to the input  directory for the current test.
      94             :         character(:, SK)    , allocatable   :: out !<  \public The scalar `allocatable` of type `character` of default kind \SK containing the path to the output directory for the current test.
      95             :     end type
      96             : 
      97             : 
      98             :     !>  \brief
      99             :     !>  This is the module `private` derived type for constructing objects that contain the path to a file and its unit.
     100             :     !>
     101             :     !>  \details
     102             :     !>  This derived type is meant to be used internally within the parent module.<br>
     103             :     !>
     104             :     !>  \see
     105             :     !>  [dir_type](@ref pm_test::dir_type)<br>
     106             :     !>  [file_type](@ref pm_test::file_type)<br>
     107             :     !>  [test_type](@ref pm_test::test_type)<br>
     108             :     !>
     109             :     !>  \finmain{file_type}
     110             :     !>
     111             :     !>  \author
     112             :     !>  \AmirShahmoradi, September 1, 2017, 12:00 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin
     113             :     type :: file_type
     114             :         integer(IK)                         :: unit
     115             :         character(:, SK)    , allocatable   :: path
     116             :     end type
     117             : 
     118             :     !type                                    :: strSplit_type
     119             :     !    integer(IK)         , allocatable   :: sindex(:,:)
     120             :     !    character(:, SK)    , allocatable   :: val
     121             :     !end type
     122             : 
     123             :     !type                                    :: subject_type
     124             :     !    type(strSplit_type)                 :: name             !<  \public An object containing information about the name and its parts (separated by underscore) of the current test subject (a procedure, type, ...).
     125             :     !end type
     126             : 
     127             :     !type                                    :: testFunc_type
     128             :     !    character(:, SK)    , allocatable   :: id               !<  \public The current number (ID) of the test function for a given object to be tested.
     129             :     !    character(:, SK)    , allocatable   :: name             !<  \public The name of the test function for the currently given object that is to be tested.
     130             :     !    class(timer_type)   , allocatable   :: timer            !<  \public This timer is exclusively used to time test cases within a given test internally.
     131             :     !    integer(IK)                         :: counter          !<  \public The counter for the test cases within a given test function for a given test subject.
     132             :     !end type
     133             : 
     134             :     !type                                    :: current_type
     135             :     !    integer(IK)                         :: counter          !<  \public The `protected` scalar object of type `integer` of default kind \IK, containing count of individual assertions made within a specified test function or test module.<br>
     136             :     !                                                            !!          The procedure name is automatically extracted from the input procedure name to the [run](@ref pm_test::setTestFunc) method of the object of type [test_type](@ref pm_test::test_type).<br>
     137             :     !end type
     138             : 
     139             :     type                                    :: scope_type
     140             :         class(timer_type)   , allocatable   :: timer            !<  \public This timer is exclusively used to time tests at the target level/scope (that is, to time all tests within a given module, or within a given test function).<br>
     141             :         character(:, SK)    , allocatable   :: name             !<  \public The `public` scalar `allocatable` of type `character` of default kind \SK, containing the name of the current scope of testing (e.g., host name of the procedure being tested, or the procedure name).<br>
     142             :     end type
     143             : 
     144             :     type, extends(scope_type)               :: func_type
     145             :         integer(IK)                         :: counter = 0_IK   !<  \private    The `private` scalar object of type `integer` of default kind \IK, containing count of individual assertions made within a
     146             :                                                                 !!              specified test function passed, via a single call, to the `run()` method of an object of type [test_type](@ref pm_test::testt_type).<br>
     147             :         character(:, SK)    , allocatable   :: id               !<  \public     The `public` scalar `allocatable` of type `character` of default kind \SK, containing the ID of the current test for a given procedure to be tested.<br>
     148             :                                                                 !!              The presence of this component is historical when individual tests for a single target procedure were implemented in separate test functions.<br>
     149             :                                                                 !!              In such a case, the test functions of a given procedure can be suffixed with a unique ID as in `test_*_ID()` where `ID` is replaced with the test function number.<br>
     150             :                                                                 !!              Test function IDs are not necessary anymore as of ParaMonvte V.2 because individual test functions can now implement as many tests as desired.<br>
     151             :                                                                 !!              If a test function does not have an ID, the `id` component of the parent object of type [test_type](@ref pm_test::test_type) is set to the default value of `1`.<br>
     152             :     end type
     153             : 
     154             :     !>  \brief
     155             :     !>  This is the derived type [test_type](@ref pm_test::test_type) for generating objects
     156             :     !>  that facilitate testing of a series of procedures or concepts within a unified testing framework.
     157             :     !>
     158             :     !>  \details
     159             :     !>  This procedure is a constructor of [test_type](@ref pm_test::test_type).<br>
     160             :     !>  See also the documentation of [constructTest](@ref pm_test::constructTest), the constructor of [test_type](@ref pm_test::test_type).<br>
     161             :     !>
     162             :     !>  \param[in]  host        :   The input scalar of type `character` of default kind \SK containing the name of the host of the procedure or code that is being tested.<br>
     163             :     !>                              Frequently, in modern Fortran, this is the name of the module containing the procedure being tested.<br>
     164             :     !>                              This name is used only for information display purposes.<br>
     165             :     !>                              (**optional**, default = `SK_""`, that is, the target procedures to be tested are assumed to be `external` and not in a module or in any specific host scope.)
     166             :     !>  \param[in]  inp         :   The input scalar of type `character` of arbitrary length type parameter, containing the input directory
     167             :     !>                              for the current testing suite (i.e., all tests performed with the output object of this constructor).<br>
     168             :     !>                              The specified value for this argument filles the corresponding `inp` component of the `dir` component of the output `test` object.<br>
     169             :     !>                              (**optional**, default = `SK_"./input"`)
     170             :     !>  \param[in]  out         :   The input scalar of type `character` of arbitrary length type parameter, containing the output directory
     171             :     !>                              for the current testing suite (i.e., all tests performed with the output object of this constructor).<br>
     172             :     !>                              The specified value for this argument filles the corresponding `out` component of the `out` component of the output `test` object.<br>
     173             :     !>                              (**optional**, default = `SK_"./output"`)
     174             :     !>  \param[in]  traceable   :   The input scalar of type `logical` of default kind \LK, that can be used to control the display of debugging information for failed tests.<br>
     175             :     !>                              Specifying this argument directly sets the corresponding `traceable` component of the output `test` object.<br>
     176             :     !>                              By convention, setting `traceable` to `.true.` is used within the test implementations to allow displaying information about assertion failures should any occur.<br>
     177             :     !>                              The ability to display descriptive information about the source of a test failure and its origins if often desired at the time of debugging.<br>
     178             :     !>                              Therefore, this component is set to `.true.` when runtime checks are enabled (i.e., FFP macro `CHECK_ENABLED=1` is set, frequently corresponding to the debug mode)
     179             :     !>                              and `.false.` otherwise (frequently corresponding to the release mode).<br>
     180             :     !>                              This optional input argument overrides the default behavior.<br>
     181             :     !>                              **Note that if `traceable` is set to `.true._LK`, the program will `error stop` by encountering the first test failure**.<br>
     182             :     !>                              (**optional**, default = `.true._LK` if the library is built with preprocessor macro `CHECK_ENABLED=1`, otherwise `.false._LK`)
     183             :     !>
     184             :     !>  \return
     185             :     !>  `test`                  :   The output scalar object of type [test_type](@ref pm_test::test_type) containing the specifics of the runtime system shell.<br>
     186             :     !>
     187             :     !>  \interface{constructTest}
     188             :     !>  \code{.F90}
     189             :     !>
     190             :     !>      use pm_test, only: test_type
     191             :     !>      type(test_type) :: test
     192             :     !>
     193             :     !>      test = test_type(host = host, inp = inp, out = out, traceable = traceable)
     194             :     !>
     195             :     !>  \endcode
     196             :     !>
     197             :     !>  \warning
     198             :     !>  This procedure must be called by all images/processes in parallel mode as it contains a call to `image_type()%%sync`.<br>
     199             :     !>
     200             :     !>  \impure
     201             :     !>
     202             :     !>  \see
     203             :     !>  [test_type](@ref pm_test::test_type)<br>
     204             :     !>
     205             :     !>  \finmain{constructTest}
     206             :     !>
     207             :     !>  \author
     208             :     !>  \AmirShahmoradi, Tuesday March 7, 2017, 3:50 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin
     209             :     type :: test_type
     210             :         integer(IK)                         , private   :: counter = 0_IK           !<  \private    The `private` scalar object of type `integer` of default kind \IK, containing count of individual assertions made by the `run()` method of an object of type [test_type](@ref pm_test::test_type).<br>
     211             :         logical(LK)                         , private   :: asserted = .true._LK     !<  \private    The `private` scalar of type `logical` of default kind \LK that is `.true.` <b>if and only if</b> all tests within the current test suite performed by an object of [test_type](@ref pm_test::test_type) have passed successfully. It is used for test summary.<br>
     212             :         logical(LK)                         , public    :: traceable                !<  \public     The `protected` scalar object of type `logical` of default kind \LK, that can be used to control the display of debugging information for failed tests. See also the corresponding `optional` argument of the [constructor](@ref pm_test::test_type).<br>
     213             :         type(dir_type)                      , public    :: dir                      !<  \public     The `public` scalar object of type [dir_type](@ref pm_test::dir_type) containing the input and output directories for the test data.<br>
     214             :         type(func_type)                     , public    :: func                     !<  \public     The scalar object of type [func_type](@ref pm_test::func_type) containing the name of the current procedure being tested and the timer starting at the beginning of the current test (corresponding to the most recent call to the `run()` method of the object of type [test_type](@ref pm_test::test_type).<br>
     215             :         type(scope_type)                    , public    :: host                     !<  \public     The scalar object of type [scope_type](@ref pm_test::scope_type) containing the name of the current host (scoping unit) being tested and the timer starting at the beginning of the creating the object of type [test_type](@ref pm_test::test_type).<br>
     216             :         type(display_type)                  , public    :: disp                     !<  \public     The `public` scalar of type [display_type](@ref pm_io::display_type) containing tools for displaying information in the desired output, for example, debugging information when a test fails.<br>
     217             :         type(image_type)                    , public    :: image                    !<  \public     The `public` scalar of type [image_type](@ref pm_parallelism::image_type) containing information about the parallel processes in parallel mode.<br>
     218             :         type(file_type)                     , public    :: file                     !<  \public     The `public` scalar of type [file_type](@ref pm_test::file_type) that serves as convenience to store the unit number and file path of any working file during the testing.<br>
     219             :        !type(err_type)                      , private   :: err                      !<  \private    The `private` scalar of type [err_type](@ref pm_err::err_type) containing information about any errors hat may occur during the testing.<br>
     220             :     contains
     221             :        !procedure                           , pass      :: openFile
     222             :         procedure                           , pass      :: summarize => setTestSummary
     223             :         procedure                           , pass      :: assert => setTestAsserted
     224             :         procedure                           , pass      :: run => setTestFunc
     225             :     end type
     226             : 
     227             :     interface test_type
     228             :         module procedure :: constructTest
     229             :     end interface
     230             : 
     231             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     232             : 
     233             : contains
     234             : 
     235             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     236             : 
     237             :     !>  \brief
     238             :     !>  Initialize the global module variables for testing.<br>
     239             :     !>
     240             :     !>  \details
     241             :     !>  This subroutine is called internally by the procedures of this module and methods of objects of type [test_type](@ref pm_test::test_type)
     242             :     !>  at the beginning of the first test, before any testing is done to ensure information for
     243             :     !>  the final global summary by [setSummary](@ref pm_test::setSummary) is properly collected.
     244             :     !>
     245             :     !>  \interface{setInitial}
     246             :     !>  \code{.F90}
     247             :     !>
     248             :     !>      call setInitial()
     249             :     !>
     250             :     !>  \endcode
     251             :     !>
     252             :     !>  \see
     253             :     !>  [test_type](@ref pm_test::test_type)<br>
     254             :     !>  [setSummary](@ref pm_test::setSummary)<br>
     255             :     !>
     256             :     !>  \finmain{setInitial}
     257             :     !>
     258             :     !>  \author
     259             :     !>  \AmirShahmoradi, September 1, 2017, 12:00 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin
     260          46 :     subroutine setInitial()
     261             : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
     262             :         !DEC$ ATTRIBUTES DLLEXPORT :: setInitial
     263             : #endif
     264             :         ! set the mode of error handling to testing mode
     265             :         !mv_isTestingMode = .true._LK
     266          46 :         mv_uninit = .false._LK
     267          46 :         mv_image = image_type()
     268          46 :         mv_timer = timer_type()
     269          46 :         mc_passedString = bright//fgreen//SK_"passed"//reset
     270          46 :         mc_failedString = bright//fred  //SK_"FAILED"//reset
     271          46 :         mv_testCounterOld = 0_IK
     272          46 :         mv_testCounter = 0_IK
     273          46 :         mv_npass = 0_IK
     274          46 :         mv_nfail = 0_IK
     275             :         ! preallocate the names of failed tests.
     276          46 :         call setResized(mv_failedTestFuncName, mv_nfail + 1_IK)
     277          46 :         write(output_unit, "(*(g0,:,' '))")
     278          46 :     end subroutine
     279             : 
     280             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     281             : 
     282             :     !>  \brief
     283             :     !>  Summarize the collection of all tests performed on all modules (or scoping units).<br>
     284             :     !>
     285             :     !>  \details
     286             :     !>  Calling this procedure at the end of a series of tests performed by a collection of objects of type [test_type](@ref pm_test::test_type)
     287             :     !>  will lead to display of relevant summary of all passed and failed tests performed and further details about the failed tests and test timing.<br>
     288             :     !>
     289             :     !>  \interface{setSummary}
     290             :     !>  \code{.F90}
     291             :     !>
     292             :     !>      call setSummary()
     293             :     !>
     294             :     !>  \endcode
     295             :     !>
     296             :     !>  \warning
     297             :     !>  This procedure must be called by all parallel images/processes/threads.<br>
     298             :     !>
     299             :     !>  \see
     300             :     !>  [test_type](@ref pm_test::test_type)<br>
     301             :     !>  [setInitial](@ref pm_test::setInitial)<br>
     302             :     !>
     303             :     !>  \finmain{setSummary}
     304             :     !>
     305             :     !>  \author
     306             :     !>  \AmirShahmoradi, September 1, 2017, 12:00 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin
     307          46 :     subroutine setSummary()
     308             : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
     309             :         !DEC$ ATTRIBUTES DLLEXPORT :: setSummary
     310             : #endif
     311             :         use pm_io, only: TABEQV
     312             :         use pm_err, only: setAborted
     313             :         use pm_mathNumSys, only: getCountDigit
     314             :         real(RKC) :: percentageTestPassed
     315             :         character(:, SK), allocatable :: msg, color
     316             :         character(*), parameter :: format = "(*(g0,:,' '))"
     317             :         integer(IK) :: ntotal, ifail, ndigit
     318          46 :         if (mv_uninit) call setInitial()
     319          46 :         ntotal = mv_npass + mv_nfail
     320          46 :         percentageTestPassed = 0._RKC
     321          46 :         if (ntotal /= 0_IK) percentageTestPassed = 100_IK * mv_npass / real(ntotal, RKC)
     322          46 :         if (mv_image%is%first) then
     323          46 :             msg = bright//fyellow//getStr(ntotal)//SK_" tests performed. "//reset
     324          46 :             if (nint(percentageTestPassed, IK) == 0_IK) then
     325          39 :                 color = fred
     326             :             else
     327           7 :                 color = fgreen
     328             :             end if
     329          46 :             msg = msg//bright//color//getStr(percentageTestPassed,SK_"(f0.2)")//SK_"% of "//getStr(ntotal)//SK_" tests passed. "//reset
     330          46 :             if (0_IK < mv_nfail) then
     331           0 :                 color = fred
     332             :             else
     333          46 :                 color = fgreen
     334             :             end if
     335          46 :             msg = msg//bright//color//getStr(mv_nfail, SK_"(g0)")//SK_" tests failed. "//reset
     336          46 :             msg = msg//bright//fyellow//SK_"The total elapsed wall-clock time: "//getStr(mv_timer%time(since = mv_timer%start), SK_"(f0.6)")//SK_" seconds."//reset
     337          46 :             write(output_unit, format) NLC//msg//NLC
     338          46 :             ndigit = getCountDigit(mv_nfail)
     339          46 :             if (mv_nfail > 0_IK) then
     340           0 :                 write(output_unit, format) NLC//bright//fred//SK_"The following tests FAILED:"//reset//NLC
     341           0 :                 do ifail = 1, mv_nfail
     342           0 :                     write(output_unit, format) bright//fred//TABEQV//adjustr(getStr(ifail, length = ndigit))//SK_") "//mv_failedTestFuncName(ifail)%val//reset
     343             :                 end do
     344           0 :                 write(output_unit, format) NLC//bright//fred//SK_"Errors occurred while running the ParaMonte tests."
     345           0 :                 write(output_unit, format) SK_"Please report this issue at: "//bright//fcyan//underlined//SK_"https://github.com/cdslaborg/paramonte/issues"//reset//NLC//NLC
     346             : #if             !CHECK_ENABLED
     347             :                 write(output_unit, format) NLC//bright//fred//SK_"To get information about the cause of failure, compile the library with FPP runtime checks enabled and rerun the tests."//reset//NLC
     348             : #endif
     349             :             end if
     350             :         end if
     351             :         !if (mv_nfail > 0_IK) !error stop !call setAborted()
     352          46 :         call mv_image%finalize()
     353          46 :         mv_uninit = .true._LK
     354          46 :     end subroutine
     355             : 
     356             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     357             : 
     358             :     !>  \brief
     359             :     !>  Generate and return an object of type [test_type](@ref pm_test::test_type).
     360             :     !>
     361             :     !>  \details
     362             :     !>  This procedure is a constructor of [test_type](@ref pm_test::test_type).<br>
     363             :     !>  See also the documentation of [test_type](@ref pm_test::test_type).<br>
     364             :     !>
     365             :     !>  \param[in]  host        :   The input scalar of type `character` of default kind \SK containing the name of the host of the procedure or code that is being tested.<br>
     366             :     !>                              Frequently, in modern Fortran, this is the name of the module containing the procedure being tested.<br>
     367             :     !>                              This name is used only for information display purposes.<br>
     368             :     !>                              (**optional**, default = `SK_"@unknown_scope"`, that is, the target procedures to be tested are assumed to be `external` and not in a module or in any specific host scope.)
     369             :     !>  \param[in]  inp         :   The input scalar of type `character` of arbitrary length type parameter, containing the input directory
     370             :     !>                              for the current testing suite (i.e., all tests performed with the output object of this constructor).<br>
     371             :     !>                              The specified value for this argument filles the corresponding `inp` component of the `dir` component of the output `test` object.<br>
     372             :     !>                              (**optional**, default = `SK_"./input"`)
     373             :     !>  \param[in]  out         :   The input scalar of type `character` of arbitrary length type parameter, containing the output directory
     374             :     !>                              for the current testing suite (i.e., all tests performed with the output object of this constructor).<br>
     375             :     !>                              The specified value for this argument filles the corresponding `out` component of the `out` component of the output `test` object.<br>
     376             :     !>                              (**optional**, default = `SK_"./output"`)
     377             :     !>  \param[in]  traceable   :   The input scalar of type `logical` of default kind \LK, that can be used to control the display of debugging information for failed tests.<br>
     378             :     !>                              Specifying this argument directly sets the corresponding `traceable` component of the output `test` object.<br>
     379             :     !>                              By convention, setting `traceable` to `.true.` is used within the test implementations to allow displaying information about assertion failures should any occur.<br>
     380             :     !>                              The ability to display descriptive information about the source of a test failure and its origins if often desired at the time of debugging.<br>
     381             :     !>                              Therefore, this component is set to `.true.` when runtime checks are enabled (i.e., FFP macro `CHECK_ENABLED=1` is set, frequently corresponding to the debug mode)
     382             :     !>                              and `.false.` otherwise (frequently corresponding to the release mode).<br>
     383             :     !>                              This optional input argument overrides the default behavior.<br>
     384             :     !>                              **Note that if `traceable` is set to `.true._LK`, the program will `error stop` by encountering the first test failure**.<br>
     385             :     !>                              (**optional**, default = `.true._LK` if the library is built with preprocessor macro `CHECK_ENABLED=1`, otherwise `.false._LK`)
     386             :     !>
     387             :     !>  \return
     388             :     !>  `test`                  :   The output scalar object of type [test_type](@ref pm_test::test_type) containing the specifics of the runtime system shell.<br>
     389             :     !>
     390             :     !>  \interface{constructTest}
     391             :     !>  \code{.F90}
     392             :     !>
     393             :     !>      use pm_test, only: test_type
     394             :     !>      type(test_type) :: test
     395             :     !>
     396             :     !>      test = test_type(host = host, inp = inp, out = out, traceable = traceable)
     397             :     !>
     398             :     !>  \endcode
     399             :     !>
     400             :     !>  \warning
     401             :     !>  This procedure must be called by all images/processes in parallel mode as it contains a call to [image_type()%sync()](@ref pm_parallelism::image_type)`.<br>
     402             :     !>
     403             :     !>  \impure
     404             :     !>
     405             :     !>  \see
     406             :     !>  [test_type](@ref pm_test::test_type)<br>
     407             :     !>
     408             :     !>  \finmain{constructTest}
     409             :     !>
     410             :     !>  \author
     411             :     !>  \AmirShahmoradi, Tuesday March 7, 2017, 3:50 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin
     412          74 :     function constructTest(host, inp, out, traceable) result(test)
     413             : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
     414             :         !DEC$ ATTRIBUTES DLLEXPORT :: constructTest
     415             : #endif
     416             :         use pm_err, only: err_type, setAborted
     417             :         use pm_sysPath, only: isFailedMakeDir, isDir
     418             :         character(*, SK), intent(in), optional :: host, inp, out
     419             :         logical(LK), intent(in), optional :: traceable
     420             :         type(test_type) :: test
     421           7 :         if (mv_uninit) call setInitial()
     422          74 :         test%image = mv_image
     423          74 :         test%disp = display_type()
     424          74 :         test%func%timer = timer_type()
     425          74 :         test%host%timer = timer_type()
     426          74 :         if (present(traceable)) then
     427           0 :             test%traceable = traceable
     428             :         else
     429             : #if         CHECK_ENABLED
     430          74 :             test%traceable = .true._LK
     431             : #else
     432             :             test%traceable = .false._LK
     433             : #endif
     434             :         end if
     435          74 :         if (present(host)) then
     436          74 :             test%host%name = trim(adjustl(host))
     437             :         else
     438           0 :             test%host%name = SK_"@unknown_scope"
     439             :         end if
     440          74 :         if (present(inp)) then
     441           0 :             test%dir%inp = inp
     442             :         else
     443          74 :             test%dir%inp = SK_"./input"
     444             :         end if
     445          74 :         if (present(out)) then
     446           0 :             test%dir%out = out
     447             :         else
     448          74 :             test%dir%out = SK_"./output"
     449             :         end if
     450             :         ! mkdir the output directory if it does not exists.
     451          74 :         if (test%image%is%first .and. .not. isDir(test%dir%out)) then
     452           1 :             if (isFailedMakeDir(test%dir%out)) error stop MODULE_NAME//SK_"@constructTest(): Failed to generate the test module output directory: "//test%dir%out
     453             :         end if
     454          74 :         call test%image%sync()
     455          74 :     end function
     456             : 
     457             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     458             : 
     459             :     !>  \brief
     460             :     !>  Run the input test function and verify the assertion returned by the test function.
     461             :     !>
     462             :     !>  \details
     463             :     !>  This procedure is a method of [test_type](@ref pm_test::test_type).<br>
     464             :     !>  See also the documentation of [test_type](@ref pm_test::test_type).<br>
     465             :     !>
     466             :     !>  \param[inout]   self            :   The input/output object of class [test_type](@ref pm_test::test_type) that is passed implicitly.
     467             :     !>  \param          getAssertion    :   The `external` user-specified test function that takes no input arguments and returns a scalar of type `logical` of default kind \LK.<br>
     468             :     !>                                      The returned value must be `.true.` <b>if and only if</b> the test assertion is verified within the specified test function.<br>
     469             :     !>                                      The following illustrates the generic interface of `getAssertion()`,
     470             :     !>                                      \code{.F90}
     471             :     !>                                          function getAssertion() result(assertion)
     472             :     !>                                              use pm_kind, only: LK
     473             :     !>                                              logical(LK) :: assertion
     474             :     !>                                          end function
     475             :     !>                                      \endcode
     476             :     !>  \param[in]      name            :   The input scalar of type `character` of default kind \SK containing the name of the test function as is,
     477             :     !>                                      that is, the actual name of target of the dummy argument `getAssertion()`.<br>
     478             :     !>                                      The test function name is assumed to follow the template `test_*()` or `test_*_ID()` where,<br>
     479             :     !>                                      <ol>
     480             :     !>                                          <li>    `*` is taken as the name of the target procedure that is being tested.<br>
     481             :     !>                                                  This name is used for displaying testing details on the specified output.<br>
     482             :     !>                                          <li>    `ID` is the optional integer in test function name representing the ID of the current test function for the target procedure being tested.<br>
     483             :     !>                                                  The default value for `ID` is `1`, if missing.<br>
     484             :     !>                                      </ol>
     485             :     !>                                      If the specified function name is not prefixed with `test_`, the whole name (excluding any `ID` suffix) is taken as the name of the procedure being tested.<br>
     486             :     !>                                      (**optional**, default = `SK_"@test_unknown"`.)
     487             :     !>
     488             :     !>  \interface{setTestFunc}
     489             :     !>  \code{.F90}
     490             :     !>
     491             :     !>      use pm_test, only: test_type
     492             :     !>      type(test_type) :: test
     493             :     !>
     494             :     !>      test = test_type()
     495             :     !>      call test%run(getAssertion, name = name)
     496             :     !>
     497             :     !>  \endcode
     498             :     !>
     499             :     !>  \impure
     500             :     !>
     501             :     !>  \see
     502             :     !>  [test_type](@ref pm_test::test_type)<br>
     503             :     !>
     504             :     !>  \finmain{setTestFunc}
     505             :     !>
     506             :     !>  \author
     507             :     !>  \AmirShahmoradi, Tuesday March 7, 2017, 3:50 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin
     508        3195 :     subroutine setTestFunc(self, getAssertion, name)
     509             : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
     510             :         !DEC$ ATTRIBUTES DLLEXPORT :: setTestFunc
     511             : #endif
     512             :         use pm_strASCII, only: isStrInteger
     513             :         procedure(logical(LK)) :: getAssertion
     514             :         class(test_type), intent(inout) :: self
     515             :         character(*, SK), intent(in), optional :: name
     516             :         logical(LK) :: assertion
     517        3195 :         if (present(name)) then
     518        3195 :             self%func%name = trim(adjustl(name))
     519        3195 :             if (4_IK < len(self%func%name, IK)) then
     520        3195 :                 if (self%func%name(1:5) == SK_"test_") self%func%name = self%func%name(6:)
     521             :             end if
     522             :             block
     523             :                 use pm_arrayFind, only: setLoc
     524             :                 integer(IK), allocatable :: loc(:)
     525             :                 integer(IK) :: nloc
     526        3195 :                 call setResized(loc, 2_IK)
     527        3195 :                 call setLoc(loc, nloc, self%func%name, SK_"_", blindness = 1_IK)
     528        3195 :                 if (nloc == 0_IK) then
     529          26 :                     self%func%id = SK_"#1"
     530        3169 :                 elseif (isStrInteger(self%func%name(loc(nloc) + 1 :))) then
     531        1439 :                     self%func%id = SK_"#"//self%func%name(loc(nloc) + 1 :)
     532        1439 :                     self%func%name = self%func%name(1 : loc(nloc) - 1)
     533             :                 else
     534        1730 :                     self%func%id = SK_"#1"
     535        1730 :                     self%func%name = name
     536             :                 end if
     537             :             end block
     538             :         else
     539           0 :             self%func%id = SK_"#1"
     540           0 :             self%func%name = SK_"unknown"
     541             :         end if
     542        3195 :         self%counter = 0_IK
     543             :        !self%func%name = SK_"Test_"//self%host%name(2:)//SK_"@"//trim(adjustl(name))
     544        3195 :         self%func%name = self%host%name//SK_"@"//trim(adjustl(self%func%name))
     545        3195 :         self%func%timer = timer_type()
     546        3195 :         assertion = getAssertion()
     547        3195 :         self%asserted = self%asserted .and. assertion
     548        3195 :         call self%assert(assertion)
     549        3195 :         self%func%counter = 0_IK
     550        3195 :     end subroutine
     551             : 
     552             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     553             : 
     554             :     !>  \brief
     555             :     !>  Test the validity of the input `assertion` and if it does not hold,<br>
     556             :     !>  <ol>
     557             :     !>      <li>    in `debug` and `testing` modes, call `error stop` after displaying the assertion description `desc`.<br>
     558             :     !>      <li>    in `release` mode, report the test failure and continue with the rest of the tests.<br>
     559             :     !>  </ol>
     560             :     !>
     561             :     !>  \details
     562             :     !>  This procedure is a method of the class [test_type](@ref pm_test::test_type).<br>
     563             :     !>
     564             :     !>  \param[inout]   self        :   The input/output object of class [test_type](@ref pm_test::test_type) that is passed implicitly.
     565             :     !>  \param[in]      assertion   :   The input scalar of type `logical` of default kind \LK containing the assertion to be verified.
     566             :     !>  \param[in]      desc        :   The input scalar of type `character` of default kind \SK of arbitrary length type parameter containing
     567             :     !>                                  a description of assertion for display if the assertion fails and `self%traceable` is `.true.`.<br>
     568             :     !>                                  (**optional**, default = `"unavailable"`)
     569             :     !>  \param[in]      line        :   The input scalar of type `integer` of default kind \IK,
     570             :     !>                                  containing the line number for the source of the assertion within the test.<br>
     571             :     !>                                  Supplying this information helps trace the source of the test failure should it occur.<br>
     572             :     !>                                  (**optional**)
     573             :     !>
     574             :     !>  \interface{setTestAsserted}
     575             :     !>  \code{.F90}
     576             :     !>
     577             :     !>      use pm_test, only: test_type
     578             :     !>      type(test_type) :: test
     579             :     !>
     580             :     !>      test = test_type()
     581             :     !>      call test%run(assertion = assertion, desc = desc, line = line)
     582             :     !>
     583             :     !>  \endcode
     584             :     !>
     585             :     !>  \see
     586             :     !>  [test_type](@ref pm_test::test_type)<br>
     587             :     !>
     588             :     !>  \finmain{setTestAsserted}
     589             :     !>
     590             :     !>  \author
     591             :     !>  \AmirShahmoradi, April 20, 2015, 1:11 AM, National Institute for Fusion Studies, The University of Texas at Austin
     592     2419120 :     subroutine setTestAsserted(self, assertion, desc, line)
     593             : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
     594             :         !DEC$ ATTRIBUTES DLLEXPORT :: setTestAsserted
     595             : #endif
     596             :         use pm_err, only: setAborted
     597             :         use pm_arrayInit, only: getCoreHalo
     598             :         logical(LK), intent(in) :: assertion
     599             :         class(test_type), intent(inout) :: self
     600             :         integer(IK), intent(in), optional :: line
     601             :         character(*, SK), intent(in), optional :: desc
     602             :         character(:, SK), allocatable :: statusMsg
     603             :         character(:, SK), allocatable :: desc_def
     604             :         character(:, SK), allocatable :: traceback
     605             :         character(:, SK), allocatable :: assertionID
     606             :         character(:, SK), allocatable :: errmsg
     607             :         integer(IK) :: iid
     608     4838240 :         self%func%timer%delta = self%func%timer%time(since = self%func%timer%start) - self%func%timer%clock
     609     2419120 :         self%func%timer%clock = self%func%timer%time(since = self%func%timer%start)
     610     2419120 :         self%func%counter = self%func%counter + 1_IK
     611     2419120 :         if (present(desc)) then
     612     2415888 :             desc_def = desc
     613             :         else
     614        3232 :             desc_def = SK_"unavailable"
     615             :         end if
     616     2419120 :         if (present(line)) then
     617     2283495 :             traceback = SK_"traceback test source line: "//getStr(line)
     618             :         else
     619      135625 :             traceback = SK_""
     620             :         end if
     621     2419120 :         assertionID = self%func%name//self%func%id//SK_"-"//getStr(self%func%counter)//SK_" "
     622     2419120 :         if (assertion) then
     623     2419120 :             mv_npass = mv_npass + 1_IK
     624     2419120 :             statusMsg = mc_passedString
     625             :         else
     626           0 :             mv_nfail = mv_nfail + 1_IK
     627           0 :             statusMsg = mc_failedString
     628           0 :             if (size(mv_failedTestFuncName, 1, IK) < mv_nfail) call setResized(mv_failedTestFuncName)
     629           0 :             mv_failedTestFuncName(mv_nfail)%val = assertionID
     630             :         end if
     631     2419120 :         mv_testCounter = mv_testCounter + 1_IK
     632     4838240 :         do iid = 1, self%image%count
     633     4838240 :             if (iid == self%image%id) then
     634             :                 write(self%disp%unit, "(*(A))") & ! LCOV_EXCL_LINE
     635             :                 SK_"["//adjustr(getStr(mv_testCounter))//SK_"] testing "// & ! LCOV_EXCL_LINE
     636             :                 getCoreHalo(79_IK, assertionID, SK_".", 0_IK)//SK_" "//statusMsg// & ! LCOV_EXCL_LINE
     637             :                 SK_" in "//adjustr(getStr(self%func%timer%delta, SK_"(f0.4)", TIME_FIELD_LEN))// & ! LCOV_EXCL_LINE
     638             :                 SK_" out of "//adjustr(getStr(self%func%timer%clock, SK_"(f0.4)", TIME_FIELD_LEN))// & ! LCOV_EXCL_LINE
     639     2419120 :                 SK_" seconds on image "//getStr(self%image%id)
     640             :             end if
     641             : #if         CAF_ENABLED || MPI_ENABLED
     642             :             block
     643             :                 integer :: stat
     644             :                 call execute_command_line(" ", cmdstat = stat)
     645             :                 flush(output_unit)
     646             :                 call self%image%sync()
     647             :             end block
     648             : #endif
     649             :         end do
     650     2419120 :         call self%func%timer%wait(0.0001_RKC)
     651     2419120 :         if (self%traceable .and. .not. assertion) then
     652             :             errmsg = SK_"The test assertion is FALSE."//NLC//SK_"The assertion description: "//trim(adjustl(desc_def))//NLC//traceback ! LCOV_EXCL_LINE
     653             :             call setAborted(msg = errmsg, prefix = SK_"ParaMonteTest"//self%func%name) ! LCOV_EXCL_LINE
     654             :             error stop ! LCOV_EXCL_LINE
     655             :         end if
     656     4838240 :     end subroutine
     657             : 
     658             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     659             : 
     660             :     !>  \brief
     661             :     !>  Summarize the suite of tests performed by the parent object of type [test_type](@ref pm_test::test_type) of this method.
     662             :     !>
     663             :     !>  \details
     664             :     !>  This procedure is a method of [test_type](@ref pm_test::test_type).<br>
     665             :     !>  See also the documentation of [test_type](@ref pm_test::test_type).<br>
     666             :     !>  This procedure can be effectively considered as the destructor of an object of type [test_type](@ref pm_test::test_type).<br>
     667             :     !>
     668             :     !>  \param[inout]   self    :   The input/output scalar object of class [test_type](@ref pm_test::test_type) that is to be destructed.<br>
     669             :     !>
     670             :     !>  \warning
     671             :     !>  This method must be called by all images/processes/threads as it contains a global synchronization upon exit.<br>
     672             :     !>
     673             :     !>  \impure
     674             :     !>
     675             :     !>  \see
     676             :     !>  [test_type](@ref pm_test::test_type)<br>
     677             :     !>
     678             :     !>  \finmain{setTestSummary}
     679             :     !>
     680             :     !>  \author
     681             :     !>  \AmirShahmoradi, Tuesday March 7, 2017, 3:50 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin
     682          74 :     subroutine setTestSummary(self)
     683             : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
     684             :         !DEC$ ATTRIBUTES DLLEXPORT :: setTestSummary
     685             : #endif
     686             :         use pm_arrayInit, only: getCoreHalo
     687             :         use pm_mathNumSys, only: getCountDigit
     688             :         class(test_type), intent(inout) :: self
     689             :         character(:, SK), allocatable :: statusMsg
     690          74 :         if (self%asserted) then
     691          74 :             statusMsg = mc_passedString
     692             :         else
     693           0 :             statusMsg = mc_failedString
     694             :         end if
     695          74 :         self%host%timer%delta = self%host%timer%time(since = self%host%timer%start)
     696             :         write(self%disp%unit, "(*(g0,:,' '))") bright//fyellow// & ! LCOV_EXCL_LINE
     697             :         SK_"["//adjustr(getStr(mv_testCounter - mv_testCounterOld, length = getCountDigit(mv_testCounter)))//SK_"] testing "// & ! LCOV_EXCL_LINE
     698             :         getCoreHalo(79_IK, self%host%name//SK_" ", SK_".", 0_IK)//SK_" "//statusMsg// & ! LCOV_EXCL_LINE
     699             :         SK_" in "//adjustr(getStr(self%host%timer%delta, SK_"(f0.4)", TIME_FIELD_LEN))// & ! LCOV_EXCL_LINE
     700             :         SK_" out of "//adjustr(getStr(mv_timer%time(since = mv_timer%start), SK_"(f0.4)", TIME_FIELD_LEN))// & ! LCOV_EXCL_LINE
     701          74 :         SK_" seconds on image "//getStr(self%image%id)//reset
     702          74 :         mv_testCounterOld = mv_testCounter
     703          74 :         call self%image%sync()
     704          74 :     end subroutine
     705             : 
     706             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     707             : 
     708             : #if 0
     709             :     !>  \cond excluded
     710             :     function openFile(self, path, label, status, position) result(file)
     711             :         class(test_type), intent(in) :: self
     712             :         character(*, SK), intent(in), optional :: path, label, status, position
     713             :         character(:, SK), allocatable :: prefix_def, status_def, position_def
     714             :         type(file_type) :: file
     715             :         if (present(position)) then
     716             :             position_def = position
     717             :         else
     718             :             position_def = SK_"asis"
     719             :         end if
     720             :         if (present(status)) then
     721             :             status_def = status
     722             :         else
     723             :             status_def = SK_"unknown"
     724             :         end if
     725             :         if (present(path)) then
     726             :             file%path = path
     727             :         else
     728             :             if (present(label)) then
     729             :                 prefix_def = SK_"@"//label
     730             :             else
     731             :                 prefix_def = SK_""
     732             :             end if
     733             :             file%path = self%dir%out//SK_"/"//self%func%name//prefix_def//SK_"@"//getStr(self%image%id)//SK_".txt"
     734             :         end if
     735             :         !if (getStrLower(position_def)=="append") status_def = SK_"old"
     736             : #if     INTEL_ENABLED && WINDOWS_ENABLED
     737             : #define SHARED, shared
     738             : #else
     739             : #define SHARED
     740             : #endif
     741             :         open(file = file%path, newunit = file%unit, status = status_def, position = position_def SHARED)
     742             : #undef  SHARED
     743             :     end function
     744             : 
     745             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     746             : 
     747             :     !>  \brief
     748             :     ! Return the negative natural logarithm of MVN distribution evaluated at the input vector `Point` of length `ndim`.
     749             :     function getLogFuncMVN(ndim,Point) result(logFunc) BIND_C
     750             :         implicit none
     751             :         integer(IK) , intent(in)        :: ndim
     752             :         real(RKC)   , intent(in)        :: Point(ndim)
     753             :         real(RKC)                       :: logFunc
     754             : #if     CFI_ENABLED
     755             :         value                           :: ndim
     756             : #endif
     757             :         real(RKC)   , parameter         :: LOG_INVERSE_SQRT_TWO_PI = log(1._RKC / sqrt(2._RKC * acos(-1._RKC)))
     758             : 
     759             :         !block
     760             :         !    use pm_sysShell, only: sleep
     761             :         !    use pm_err, only: err_type
     762             :         !    type(err_type) :: err
     763             :         !    call sleep(seconds=5000.e-6_RKC,err=err)
     764             :         !end block
     765             : 
     766             :         !block
     767             :         !    real(RKC), allocatable :: unifrnd(:,:)
     768             :         !    allocate(unifrnd(200,20))
     769             :         !    call random_number(unifrnd)
     770             :         !    logFunc = sum(unifrnd) - 0.5_RKC * sum(Point**2) - sum(unifrnd)
     771             :         !    deallocate(unifrnd)
     772             :         !end block
     773             : 
     774             :         logFunc = ndim * LOG_INVERSE_SQRT_TWO_PI - 0.5_RKC * sum(Point**2_IK)
     775             : 
     776             :         !block
     777             :         !    integer(IK), parameter :: nmode = 2_IK
     778             :         !    real(RKC):: LogAmplitude(nmode), mean(nmode), invCov(nmode), logSqrtDetInvCovMat(nmode)
     779             :         !    LogAmplitude           = [1._RKC, 1._RKC]
     780             :         !    mean                   = [0._RKC, 7._RKC]
     781             :         !    invCov              = [1._RKC,1._RKC]
     782             :         !    logSqrtDetInvCovMat    = [1._RKC,1._RKC]
     783             :         !    logFunc = getLogProbGausMix ( nmode = 2_IK &
     784             :         !                                , nd = 1_IK &
     785             :         !                                , np = 1_IK &
     786             :         !                                , LogAmplitude = LogAmplitude &
     787             :         !                                , mean = mean &
     788             :         !                                , invCov = invCov &
     789             :         !                                , logSqrtDetInvCovMat = logSqrtDetInvCovMat &
     790             :         !                                , Point = Point(1) &
     791             :         !                                )
     792             :         !end block
     793             : 
     794             :     end function getLogFuncMVN
     795             : 
     796             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     797             : 
     798             :     !>  \brief
     799             :     !>  Return the negative natural logarithm of MVN distribution evaluated at the input vector `Point` of length `ndim`.
     800             :     function getLogFuncBanana2D(ndim,Point) result(logFunc) BIND_C
     801             :         use pm_distMultiNorm, only: getMultiNormLogPDFNF
     802             :         use pm_distMultiNorm, only: getMultiNormLogPDF
     803             :         use pm_mathLogSumExp, only: getLogSumExp
     804             :         implicit none
     805             :         integer(IK) , intent(in)        :: ndim
     806             :         real(RKC)   , intent(in)        :: Point(ndim)
     807             :         real(RKC)                       :: logFunc
     808             : #if     CFI_ENABLED
     809             :         value                           :: ndim
     810             : #endif
     811             :         integer(IK) , parameter         :: NPAR = 2_IK ! sum(Banana,gaussian) normalization factor
     812             :         real(RKC)   , parameter         :: normfac = 0.3_RKC ! 0.6_RKC ! sum(Banana,gaussian) normalization factor
     813             :         real(RKC)   , parameter         :: lognormfac = log(normfac) ! sum(Banana,gaussian) normalization factor
     814             :         real(RKC)   , parameter         :: a = 0.7_RKC, b = 1.5_RKC ! parameters  of Banana function
     815             :         real(RKC)   , parameter         :: MeanB(NPAR) = [ -5.0_RKC , 0._RKC ] ! mean vector of Banana function
     816             :         real(RKC)   , parameter         :: MeanG(NPAR) = [  3.5_RKC , 0._RKC ] ! mean vector of Gaussian function
     817             :         real(RKC)   , parameter         :: CovMatB(NPAR,NPAR) = reshape([0.25_RKC,0._RKC,0._RKC,0.81_RKC],shape=shape(CovMatB)) ! Covariance matrix of Banana function
     818             :         real(RKC)   , parameter         :: CovMatG(NPAR,NPAR) = reshape([0.15_RKC,0._RKC,0._RKC,0.15_RKC],shape=shape(CovMatB)) ! Covariance matrix of Gaussian function
     819             :         real(RKC)   , parameter         :: InvCovMatB(NPAR,NPAR) = reshape([4._RKC,0._RKC,0._RKC,1.23456790123457_RKC],shape=shape(InvCovMatB)) ! Inverse Covariance matrix of Banana function
     820             :         real(RKC)   , parameter         :: InvCovMatG(NPAR,NPAR) = reshape([6.66666666666667_RKC,0._RKC,0._RKC,6.66666666666667_RKC],shape=shape(InvCovMatG)) ! Inverse Covariance matrix of Gaussian function
     821             :         real(RKC)   , parameter         :: logSqrtDetInvCovB = log(sqrt(4.93827160493827_RKC)) ! Determinant of the Inverse Covariance matrix of Banana function
     822             :         real(RKC)   , parameter         :: logSqrtDetInvCovG = log(sqrt(44.4444444444445_RKC)) ! Determinant of the Inverse Covariance matrix of Gaussian function
     823             :         real(RKC)                       :: PointSkewed(NPAR) ! transformed parameters that transform the Gaussian to the Banana function
     824             :         real(RKC)                       :: LogProb(2) ! logProbMVN, logProbBanana
     825             : 
     826             :         PointSkewed(1) = -Point(1)
     827             :         PointSkewed(2) = +Point(2)
     828             : 
     829             :         ! Gaussian function
     830             : 
     831             :         LogProb(1) = lognormfac + getMultiNormLogPDF(PointSkewed, MeanG, InvCovMatG, getMultiNormLogPDFNF(ndim, logSqrtDetInvCovG)) ! logProbMVN
     832             : 
     833             :         ! Do variable transformations for the Skewed-Gaussian (banana) function.
     834             : 
     835             :         PointSkewed(2)  = a * PointSkewed(2)
     836             :         PointSkewed(1)  = PointSkewed(1)/a - b*(PointSkewed(2)**2 + a**2)
     837             : 
     838             :         ! Banana function
     839             : 
     840             :         LogProb(2) = lognormfac + getMultiNormLogPDF(PointSkewed, MeanB, InvCovMatB, getMultiNormLogPDFNF(ndim, logSqrtDetInvCovB)) ! logProbBanana
     841             : 
     842             :         !MeanBnew(2) = a * MeanBnew(2)
     843             :         !MeanBnew(1) = MeanBnew(1)/a - b*(MeanBnew(2)**2 + a**2)
     844             : 
     845             :         !logFunc = LogProb(2)
     846             :         !logFunc = lognormfac + LogProb(1)
     847             :         logFunc = getLogSumExp(LogProb, maxval(LogProb))
     848             : 
     849             :     end function getLogFuncBanana2D
     850             : 
     851             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     852             : 
     853             :     !>  \brief
     854             :     !>  Return the 2D EggBox function in the cubic unit domain [0,1].
     855             :     !>  \remark
     856             :     !>  The log(integral) of this function in the unit domain is `235.88`.
     857             :     function getLogFuncEggBox2D(ndim,Point) result(logFunc) BIND_C
     858             :         integer(IK) , intent(in)        :: ndim
     859             :         real(RKC)   , intent(in)        :: Point(ndim)
     860             :         real(RKC)                       :: logFunc
     861             :         real(RKC)   , parameter         :: PI = acos(-1._RKC)
     862             : #if     CFI_ENABLED
     863             :         value :: ndim
     864             : #endif
     865             :         !logFunc = (2._RKC + cos(0.5_RKC*Point(1)) * cos(0.5_RKC*Point(2)) )**5_IK
     866             :         logFunc = (2._RKC + cos(5_IK * PI * Point(1) - 2.5_RKC * PI) * cos(5_IK * PI * Point(2) - 2.5_RKC * PI)) ** 5.0
     867             :     end function getLogFuncEggBox2D
     868             :     !>  \endcond excluded
     869             : 
     870             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     871             : #endif
     872             : 
     873             : end module pm_test ! LCOV_EXCL_LINE

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