ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_paramonte.F90
Go to the documentation of this file.
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.
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
28
36
37!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
38
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
68 character(*, SK), parameter :: PARAMONTE_WEB_REPOSITORY = SK_"https://github.com/cdslaborg/paramonte"
69
89 character(*, SK), parameter :: PARAMONTE_WEB_ISSUES = SK_"https://github.com/cdslaborg/paramonte/issues"
90
110 character(*, SK), parameter :: PARAMONTE_WEB_DOC = SK_"https://www.cdslab.org/paramonte/"
111
112!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
113
133 character(*, SK), parameter :: PARAMONTE_COMPILER_VERSION = compiler_version()
134
154 character(*, SK), parameter :: PARAMONTE_COMPILER_OPTIONS = compiler_options()
155
156!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
157
178 end type
179
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
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
295
309 character(*, SK), parameter :: envname = PROGNAME//" programming language"
310
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
360 type(envis_type), parameter :: envis = envis_type(PROGLANG = .true._LK)
361
363#undef PROGLANG
364#undef PROGNAME
365
366
367!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
368
404#if INTEL_ENABLED || GNU_ENABLED
405 character(*, SK), parameter :: PARAMONTE_BUILD_DATE = __TIMESTAMP__
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
413
414!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
415
449#if defined __PARAMONTE_VERSION__
450 character(*, SK), parameter :: PARAMONTE_VERSION = __PARAMONTE_VERSION__
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
461
462!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
463
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
549contains
550
551!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
552
619 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 lwsize_def = getOption(4_IK, lwsize)
635 twsize_def = getOption(2_IK, twsize)
636 rwsize_def = getOption(4_IK, rwsize)
637 bwsize_def = getOption(2_IK, bwsize)
638 lwfill_def = getOption(SK_"%", lwfill)
639 twfill_def = getOption(SK_"%", twfill)
640 rwfill_def = getOption(SK_"%", rwfill)
641 bwfill_def = getOption(SK_"%", bwfill)
642 fill_def = getOption(SK_" ", fill)
643 lenNLC = len(NLC, IK)
644 call setResized(locNLC, 47_IK)
645 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 linewidth = max(maxval(locNLC(2 : lenLocNLC) - locNLC(1 : lenLocNLC - 1)) - lenNLC, getOption(132_IK, width))! + lwsize_def + rwsize_def
649 splashLen = lenNLC + (twsize_def + lenLocNLC + bwsize_def) * (linewidth + lenNLC)
650 allocate(character(splashLen, SK) :: splash)
651 splash(1 : lenNLC) = NLC
652 ibeg = lenNLC + 1_IK
653 do iline = 1, twsize_def
654 iend = ibeg + linewidth - 1_IK
655 splash(ibeg : iend) = repeat(twfill_def, linewidth)
656 ibeg = iend + lenNLC + 1
657 splash(iend + 1 : ibeg - 1) = NLC
658 end do
659 if (lenLocNLC > 0_IK) then
660 sbeg = 1_IK ! splash line beginning.
661 do iline = 1, lenLocNLC
662 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 )
671 sbeg = locNLC(iline) + lenNLC
672 ibeg = iend + lenNLC + 1_IK
673 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 do iline = 1, bwsize_def
686 iend = ibeg + linewidth - 1_IK
687 splash(ibeg : iend) = repeat(bwfill_def, linewidth)
688 ibeg = iend + lenNLC + 1
689 splash(iend + 1 : ibeg - 1) = NLC
690 end do
691 end function
692
693!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
694
695end module pm_paramonte ! LCOV_EXCL_LINE
Center the contents of the input array within the output arrayCentered while filling the newly added ...
Return an allocatable array containing the indices of the locations within the input array where the ...
Allocate or resize (shrink or expand) an input allocatable scalar string or array of rank 1....
Generate and return the value of the optional input argument if it is present, otherwise,...
Definition: pm_option.F90:135
This module contains procedures and generic interfaces for resizing an input array and centering the ...
This module contains procedures and generic interfaces for finding locations of a pattern in arrays o...
This module contains procedures and generic interfaces for resizing allocatable arrays of various typ...
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
Definition: pm_kind.F90:268
integer, parameter LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
Definition: pm_kind.F90:541
integer, parameter IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
Definition: pm_kind.F90:540
integer, parameter SK
The default character kind in the ParaMonte library: kind("a") in Fortran, c_char in C-Fortran Intero...
Definition: pm_kind.F90:539
This module contains procedures, generic interfaces, and types for generating default values for opti...
Definition: pm_option.F90:34
This module contains procedures and data that provide general information about the ParaMonte library...
character(*, SK), parameter PARAMONTE_SPLASH
The public scalar character constant of default kind SK containing the ParaMonte splash information.
type(envis_type), parameter envis
The scalar constant object of type envis_type whose components are all .false. except the environment...
PURE character(:, SK) function, allocatable getParaMonteSplash(width, lwsize, twsize, rwsize, bwsize, fill, lwfill, twfill, rwfill, bwfill)
Generate and return an allocatable scalar character of default kind SK containing the ParaMonte libra...
character(*, SK), parameter PARAMONTE_VERSION
The public scalar character constant of default kind SK containing the ParaMonte library version.
character(*, SK), parameter PARAMONTE_BUILD_DATE
The public scalar character constant of default kind SK containing the ParaMonte library build date.
type(paramonte_type), parameter paramonte
The scalar module constant of type paramonte_type that can be used within the ParaMonte library to si...
character(*, SK), parameter PARAMONTE_COMPILER_OPTIONS
The scalar constant of type character of default kind SK, containing the compiler options used to bui...
character(*, SK), parameter envname
The scalar constant of type character of default kind SK containing the full generic name of the prog...
character(*, SK), parameter PARAMONTE_WEB_REPOSITORY
The scalar constant of type character of default kind SK, containing the web portal address for the P...
character(*, SK), parameter PARAMONTE_WEB_DOC
The scalar constant of type character of default kind SK, containing the web portal address for repor...
character(*, SK), parameter PARAMONTE_WEB_ISSUES
The scalar constant of type character of default kind SK, containing the web portal address for repor...
character(*, SK), parameter PARAMONTE_COMPILER_VERSION
The scalar constant of type character of default kind SK, containing the version of the compiler with...
This module contains classes and procedures for various string manipulations and inquiries.
Definition: pm_str.F90:49
character(*, SK), parameter UNDEFINED
The scalar character parameter of default kind SK containing "UNDEFINED".
Definition: pm_str.F90:93
character(*, SK), parameter NLC
The newline character of default kind SK as returned by new_line(SK_"a").
Definition: pm_str.F90:88
This is the derived type for generating objects containing scalar components of type logical named af...
This a derived type devoid of any components or methods whose instantiated objects are used within th...