ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_timer Module Reference

This module contains the timer procedures and derived types to facilitate timing applications at runtime.
More...

Data Types

interface  getTime_proc
 This is the abstract interface of the setTime static type-bound procedure component of timer_type abstract type that performs the timing since a user-specified or a processor-dependent origin. More...
 
interface  setIdle_proc
 This is the abstract interface of the wait static type-bound procedure component of timer_type abstract type that keeps the system waiting for the specified amount of seconds.
More...
 
interface  timer_type
 This is the abstract base derived type that serves as a simple container template for other timer classes in the ParaMonte library. More...
 
interface  timerCPU_type
 This is the timerCPU_type class, containing attributes and static methods for setting up a timer based on the CPU-clock, using the Fortran intrinsic cpu_time(). More...
 
interface  timerDAT_type
 This is the timerDAT_type class, containing attributes and static methods for setting up a timer based on the system-clock, using the Fortran intrinsic date_and_time(). More...
 
interface  timerMPI_type
 This is the timerMPI_type class, containing attributes and static methods for setting up a timer based on the MPI intrinsic MPI_Wtime(). More...
 
interface  timerOMP_type
 This is the timerMPI_type class, containing attributes and static methods for setting up a timer based on the OpenMP intrinsic omp_get_wtime(). More...
 
interface  timerSYS_type
 This is the timerSYS_type class, containing attributes and static methods for setting up a timer based on the system-clock, using the Fortran intrinsic system_clock(). More...
 

Functions/Subroutines

type(timerMPI_type) type(timerOMP_type) type(timerSYS_type) function timer_typer ()
 This is the constructor of the timer_type class.
More...
 
type(timerCPU_type) function timerCPU_typer ()
 This is the constructor of the timerCPU_type class.
More...
 
type(timerDAT_type) function timerDAT_typer ()
 This is the constructor of the timerDAT_type class.
More...
 
type(timerMPI_type) function timerMPI_typer ()
 This is the constructor of the timerMPI_type class.
More...
 
type(timerOMP_type) function timerOMP_typer ()
 This is the constructor of the timerOMP_type class.
More...
 
type(timerSYS_type) function timerSYS_typer ()
 This is the constructor of the timerSYS_type class.
More...
 
real(RKD) function getTime (since)
 Generate and return the time in units of seconds since the specified time since or since an arbitrary processor-dependent time. More...
 
real(RKD) function getTimeCPU (since)
 Generate and return the CPU time in units of seconds since the specified time since or since an arbitrary processor-dependent time. More...
 
real(RKD) function getTimeDAT (since)
 Generate and return the calendrical clock time in units of seconds since the specified time since or since the Gregorian calendrical time origin. More...
 
real(RKD) function getTimeMPI (since)
 Generate and return the MPI clock time in units of seconds since the specified time since or since an arbitrary time origin set by the MPI library. More...
 
real(RKD) function getTimeOMP (since)
 Generate and return the OpenMP clock time in units of seconds since the specified time since or since an arbitrary time origin set by the OpenMP library. More...
 
real(RKD) function getTimeSYS (since)
 Generate and return the system clock time in units of seconds since the specified time since or since an arbitrary processor-dependent time origin. More...
 
real(RKD) function getResTimer ()
 Generate and return the time resolution by the build-time default timer of the ParaMonte library as used in timer_type, in units of seconds. More...
 
real(RKD) function getResTimerCPU ()
 Generate and return the time resolution of the Fortran intrinsic cpu_time(), used in timerCPU_type, in units of seconds. More...
 
pure real(RKD) function getResTimerDAT ()
 Generate and return the time resolution of the Fortran intrinsic date_and_time(), used in timerDAT_type, in units of seconds. More...
 
real(RKD) function getResTimerMPI ()
 Generate and return the time resolution of the MPI intrinsic timer, used in timerMPI_type, in units of seconds. More...
 
real(RKD) function getResTimerOMP ()
 Generate and return the time resolution of the OpenMP intrinsic timer, used in timerOMP_type, in units of seconds. More...
 
real(RKD) function getResTimerSYS ()
 Generate and return the time resolution of the Fortran intrinsic system_clock(), used in timerSYS_type, in units of seconds. More...
 
subroutine setIdle (seconds)
 Set the processor in sleep mode for the input requested number of seconds via the default timer of the ParaMonte library. More...
 
subroutine setIdleCPU (seconds)
 Set the processor in sleep mode for the input requested number of seconds via the Fortran intrinsic cpu_time(). More...
 
subroutine setIdleDAT (seconds)
 Set the processor in sleep mode for the input requested number of seconds via the Fortran intrinsic date_and_time(). More...
 
subroutine setIdleMPI (seconds)
 Set the processor in sleep mode for the input requested number of seconds via the MPI intrinsic timer MPI_Wtime(). More...
 
subroutine setIdleOMP (seconds)
 Set the processor in sleep mode for the input requested number of seconds via the OpenMP intrinsic timer omp_get_wtime(). More...
 
subroutine setIdleSYS (seconds)
 Set the processor in sleep mode for the input requested number of seconds via the Fortran intrinsic system_clock(). More...
 

Variables

character(*, SK), parameter MODULE_NAME = "@pm_timer"
 

Detailed Description

This module contains the timer procedures and derived types to facilitate timing applications at runtime.

Currently, five types of timers are available in this module in addition to a generic timer class.

  1. timer_type
    1. This is a generic timer class which defaults to the most appropriate timer type based on the ParaMonte library build.
      1. If the ParaMonte library is built for serial applications, then timerSYS_type is used.
      2. If the ParaMonte library is built for MPI-parallel applications, then timerMPI_type is used.
      3. If the ParaMonte library is built for OpenMP-parallel applications, then timerOMP_type is used.
    2. Use this timer class if you are not sure which timer is suitable for the task.
  2. timerCPU_type
    1. This timer relies on the Fortran intrinsic cpu_time().
    2. This timer is not ideal for timing parallel or multi-threaded applications because cpu_time() returns the sum of the total time spent by all processes.
    3. The timing precision of the Fortran intrinsic cpu_time() is processor-dependent.
    4. There is currently no direct way of measuring the time resolution of cpu_time().
    5. This module uses a simple cycling method to compute the resolution of the time returned by cpu_time().
      It computes the resolution as the least noticeable non-zero time difference between any two time points returned by cpu_time().
    6. On processors with no clock, the time returned by cpu_time() is negative to indicate the lack of a processor clock.
      The clock existence is verified only if the library is built with the preprocessor macro CHECK_ENABLED=1.
  3. timerDAT_type
    1. This timer relies on the Fortran intrinsic date_and_time().
    2. This timer is not ideal for high resolution time measurements because the resolution is limited to milliseconds.
  4. timerMPI_type
    1. This timer relies on the MPI library intrinsic mpi_wtime().
    2. The clock time resolution is computed by calling mpi_wtick().
    3. This timer is ideal for high resolution (e.g., nanseconds) time measurements in distributed parallel applications.
    4. This timer is available only when the ParaMonte library is built with MPI support enabled.
    5. This timer requires the MPI library to have been already initialized via mpi_init() and not have been finalized via mpi_finalize().
      These conditions are verified only if the library is built with the preprocessor macro CHECK_ENABLED=1.
  5. timerOMP_type
    1. This timer relies on the OpenMP library intrinsic omp_get_wtime().
    2. The clock time resolution is computed by calling omp_get_wtick().
    3. This timer is ideal for high resolution time measurements in multithreaded parallel applications.
    4. This timer is available only when the ParaMonte library is built with OpenMP support enabled.
    5. This timer measures elapsed wall clock time in seconds.
    6. The time is measured per thread.
    7. No guarantee can be made that two distinct threads measure the same time.
    8. Time is measured from some time in the past, which is an arbitrary time guaranteed not to change during the execution of the program.
  6. timerSYS_type
    1. This timer relies on the Fortran intrinsic system_clock().
    2. This timer can be used in most parallel or serial timing scenarios.
    3. The resolution of this timer is processor-dependent.
    4. Time is measured from some time in the past, which is an arbitrary time guaranteed not to change during the execution of the program.
    5. The time is measured per processor.
    6. On processors with no clock, sysctem_clock() returns a negative clock count to indicate the lack of a processor clock.
      The clock existence is verified only if the library is built with the preprocessor macro CHECK_ENABLED=1.
Test:
test_pm_timer


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, March 22, 2012, 2:21 PM, National Institute for Fusion Studies, The University of Texas Austin

Function/Subroutine Documentation

◆ getResTimer()

real(RKD) function pm_timer::getResTimer

Generate and return the time resolution by the build-time default timer of the ParaMonte library as used in timer_type, in units of seconds.

This function relies on either of the following lower-level module functions depending on the ParaMonte library build:

  1. getResTimerMPI for MPI-enabled parallel builds of the ParaMonte library.
  2. getResTimerOMP for OpenMP-enabled parallel builds of the ParaMonte library.
  3. getResTimerSYS for all other parallel or serial builds of the ParaMonte library.
Returns
resolution : The output scalar of type real of kind double precision RKD, containing the time resolution of the clock in units of seconds.


Possible calling interfaces

real(RKD) :: resolution
resolution = getResTimer()
This module contains the timer procedures and derived types to facilitate timing applications at runt...
Definition: pm_timer.F90:99
real(RKD) function getResTimer()
Generate and return the time resolution by the build-time default timer of the ParaMonte library as u...
Definition: pm_timer.F90:1261
Remarks
The procedures under discussion are impure.
See also
timer_type


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK, RKD
4 use pm_io, only: display_type
5 use pm_timer, only: getResTimer
6
7 implicit none
8
9 type(display_type) :: disp
10 disp = display_type(file = "main.out.F90")
11
12 call disp%skip()
13 call disp%show("getResTimer()")
14 call disp%show( getResTimer() )
15 call disp%skip()
16
17end program example
This is a generic method of the derived type display_type with pass attribute.
Definition: pm_io.F90:11726
This is a generic method of the derived type display_type with pass attribute.
Definition: pm_io.F90:11508
This module contains classes and procedures for input/output (IO) or generic display operations on st...
Definition: pm_io.F90:252
type(display_type) disp
This is a scalar module variable an object of type display_type for general display.
Definition: pm_io.F90:11393
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 RKD
The double precision real kind in Fortran mode. On most platforms, this is an 64-bit real kind.
Definition: pm_kind.F90:568
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
Generate and return an object of type display_type.
Definition: pm_io.F90:10282

Example Unix compile command via Intel ifort compiler
1#!/usr/bin/env sh
2rm main.exe
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example Windows Batch compile command via Intel ifort compiler
1del main.exe
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
4main.exe

Example Unix / MinGW compile command via GNU gfortran compiler
1#!/usr/bin/env sh
2rm main.exe
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example output
1
3+0.10000000000000001E-8
4
5
Test:
test_pm_timer


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas Austin

Definition at line 1260 of file pm_timer.F90.

References getResTimerMPI(), getResTimerOMP(), and getResTimerSYS().

Here is the call graph for this function:

◆ getResTimerCPU()

real(RKD) function pm_timer::getResTimerCPU

Generate and return the time resolution of the Fortran intrinsic cpu_time(), used in timerCPU_type, in units of seconds.

See also the documentation details of pm_timer.

Returns
resolution : The output scalar of type real of kind double precision RKD, containing the time resolution of the CPU clock in units of seconds.


Possible calling interfaces

real(RKD) :: resolution
resolution = getResTimerCPU()
real(RKD) function getResTimerCPU()
Generate and return the time resolution of the Fortran intrinsic cpu_time(), used in timerCPU_type,...
Definition: pm_timer.F90:1316
Remarks
The procedures under discussion are impure.
See also
getResTimerSYS


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK, RKD
4 use pm_io, only: display_type
5 use pm_timer, only: getResTimerCPU
6
7 implicit none
8
9 type(display_type) :: disp
10 disp = display_type(file = "main.out.F90")
11
12 call disp%skip()
13 call disp%show("getResTimerCPU()")
14 call disp%show( getResTimerCPU() )
15 call disp%skip()
16
17end program example

Example Unix compile command via Intel ifort compiler
1#!/usr/bin/env sh
2rm main.exe
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example Windows Batch compile command via Intel ifort compiler
1del main.exe
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
4main.exe

Example Unix / MinGW compile command via GNU gfortran compiler
1#!/usr/bin/env sh
2rm main.exe
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example output
1
3+0.99999999999991589E-6
4
5
Test:
test_pm_timer


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas Austin

Definition at line 1315 of file pm_timer.F90.

Referenced by pm_timer::timerCPU_type::timerCPU_typer().

Here is the caller graph for this function:

◆ getResTimerDAT()

pure real(RKD) function pm_timer::getResTimerDAT

Generate and return the time resolution of the Fortran intrinsic date_and_time(), used in timerDAT_type, in units of seconds.

See also the documentation details of pm_timer.

Returns
resolution : The output scalar of type real of kind double precision RKD, containing the time resolution of date_and_time() in units of seconds.


Possible calling interfaces

real(RKD) :: resolution
resolution = getResTimerDAT()
pure real(RKD) function getResTimerDAT()
Generate and return the time resolution of the Fortran intrinsic date_and_time(), used in timerDAT_ty...
Definition: pm_timer.F90:1373
Remarks
The procedures under discussion are pure.
See also
getResTimerSYS


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK, RKD
4 use pm_io, only: display_type
5 use pm_timer, only: getResTimerDAT
6
7 implicit none
8
9 type(display_type) :: disp
10 disp = display_type(file = "main.out.F90")
11
12 call disp%skip()
13 call disp%show("getResTimerDAT()")
14 call disp%show( getResTimerDAT() )
15 call disp%skip()
16
17end program example

Example Unix compile command via Intel ifort compiler
1#!/usr/bin/env sh
2rm main.exe
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example Windows Batch compile command via Intel ifort compiler
1del main.exe
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
4main.exe

Example Unix / MinGW compile command via GNU gfortran compiler
1#!/usr/bin/env sh
2rm main.exe
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example output
1
3+0.10000000000000000E-2
4
5
Test:
test_pm_timer


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas Austin

Definition at line 1372 of file pm_timer.F90.

◆ getResTimerMPI()

real(RKD) function pm_timer::getResTimerMPI

Generate and return the time resolution of the MPI intrinsic timer, used in timerMPI_type, in units of seconds.

See also the documentation details of pm_timer.

Returns
resolution : The output scalar of type real of kind double precision RKD, containing the time resolution of mpi_wtime() in units of seconds.


Possible calling interfaces

real(RKD) :: resolution
resolution = getResTimerMPI()
real(RKD) function getResTimerMPI()
Generate and return the time resolution of the MPI intrinsic timer, used in timerMPI_type,...
Definition: pm_timer.F90:1424
Warning
This procedure exists only if the library is built with the preprocessor macro MPI_ENABLED=1.
See also
getResTimerSYS


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK, RKD
4 use pm_io, only: display_type
5
6 implicit none
7
8 type(display_type) :: disp
9 disp = display_type(file = "main.out.F90")
10
11#if defined MPI_VERSION
12 block
13 use pm_timer, only: getResTimerMPI
14 call disp%skip()
15 call disp%show("getResTimerMPI()")
16 call disp%show( getResTimerMPI() )
17 call disp%skip()
18 end block
19#else
20 call disp%show("The MPI timer is unavailable, because this example has not been linked with an MPI library.")
21#endif
22
23end program example

Example Unix compile command via Intel ifort compiler
1#!/usr/bin/env sh
2rm main.exe
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example Windows Batch compile command via Intel ifort compiler
1del main.exe
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
4main.exe

Example Unix / MinGW compile command via GNU gfortran compiler
1#!/usr/bin/env sh
2rm main.exe
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example output
1The MPI timer is unavailable, because this example has not been linked with an MPI library.
2
Test:
test_pm_timer


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas Austin

Definition at line 1423 of file pm_timer.F90.

Referenced by getResTimer().

Here is the caller graph for this function:

◆ getResTimerOMP()

real(RKD) function pm_timer::getResTimerOMP

Generate and return the time resolution of the OpenMP intrinsic timer, used in timerOMP_type, in units of seconds.

See also the documentation details of pm_timer.

Returns
resolution : The output scalar of type real of kind double precision RKD, containing the time resolution of mpi_wtime() in units of seconds.


Possible calling interfaces

real(RKD) :: resolution
resolution = getResTimerOMP()
real(RKD) function getResTimerOMP()
Generate and return the time resolution of the OpenMP intrinsic timer, used in timerOMP_type,...
Definition: pm_timer.F90:1477
Warning
This procedure exists only if the library is built with the preprocessor macro OMP_ENABLED=1.
See also
getResTimerSYS


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK, RKD
4 use pm_io, only: display_type
5
6 implicit none
7
8 type(display_type) :: disp
9 disp = display_type(file = "main.out.F90")
10
11#if defined _OPENMP
12 block
13 use pm_timer, only: getResTimerOMP
14 call disp%skip()
15 call disp%show("getResTimerOMP()")
16 call disp%show( getResTimerOMP() )
17 call disp%skip()
18 end block
19#else
20 call disp%show("The OMP timer is unavailable, because this example has not been linked with an OpenMP library.")
21#endif
22
23end program example

Example Unix compile command via Intel ifort compiler
1#!/usr/bin/env sh
2rm main.exe
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example Windows Batch compile command via Intel ifort compiler
1del main.exe
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
4main.exe

Example Unix / MinGW compile command via GNU gfortran compiler
1#!/usr/bin/env sh
2rm main.exe
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example output
1The OMP timer is unavailable, because this example has not been linked with an OpenMP library.
2
Test:
test_pm_timer


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas Austin

Definition at line 1476 of file pm_timer.F90.

Referenced by getResTimer().

Here is the caller graph for this function:

◆ getResTimerSYS()

real(RKD) function pm_timer::getResTimerSYS

Generate and return the time resolution of the Fortran intrinsic system_clock(), used in timerSYS_type, in units of seconds.

See also the documentation details of pm_timer.

Returns
resolution : The output scalar of type real of kind double precision RKD, containing the time resolution of the system clock in units of seconds.


Possible calling interfaces

real(RKD) :: resolution
resolution = getResTimerSYS()
real(RKD) function getResTimerSYS()
Generate and return the time resolution of the Fortran intrinsic system_clock(), used in timerSYS_typ...
Definition: pm_timer.F90:1528
Remarks
The procedures under discussion are impure.
See also
timerSYS_type


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK, RKD
4 use pm_io, only: display_type
5 use pm_timer, only: getResTimerSYS
6
7 implicit none
8
9 type(display_type) :: disp
10 disp = display_type(file = "main.out.F90")
11
12 call disp%skip()
13 call disp%show("getResTimerSYS()")
14 call disp%show( getResTimerSYS() )
15 call disp%skip()
16
17end program example

Example Unix compile command via Intel ifort compiler
1#!/usr/bin/env sh
2rm main.exe
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example Windows Batch compile command via Intel ifort compiler
1del main.exe
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
4main.exe

Example Unix / MinGW compile command via GNU gfortran compiler
1#!/usr/bin/env sh
2rm main.exe
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example output
1
3+0.10000000000000001E-8
4
5
Test:
test_pm_timer


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas Austin

Definition at line 1527 of file pm_timer.F90.

Referenced by getResTimer(), and pm_timer::timerSYS_type::timerSYS_typer().

Here is the caller graph for this function:

◆ getTime()

real(RKD) function pm_timer::getTime ( real(RKD), intent(in), optional  since)

Generate and return the time in units of seconds since the specified time since or since an arbitrary processor-dependent time.

This function relies on either of the following lower-level module functions depending on the ParaMonte library build:

  1. getTimeMPI for MPI-enabled parallel builds of the ParaMonte library.
  2. getTimeOMP for OpenMP-enabled parallel builds of the ParaMonte library.
  3. getTimeSYS for all other parallel or serial builds of the ParaMonte library.
Parameters
[in]since: The input scalar of type real of kind double precision RKD representing the time origin (epoch) in seconds.
(optional. If not present, its value is processor-dependent but constant at runtime.)
Returns
timeInSec : The output scalar of type real of kind double precision RKD, containing the time in seconds since a specified time since or since an arbitrary processor-dependent time.


Possible calling interfaces

use pm_timer, only: getTime
real(RKD) :: timeInSec
logical(LK) :: failed
timeInSec = getTime(since = since)
real(RKD) function getTime(since)
Generate and return the time in units of seconds since the specified time since or since an arbitrary...
Definition: pm_timer.F90:858
Remarks
The procedures under discussion are impure.
See also
timer_type


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK, RKD
4 use pm_io, only: display_type
5 use pm_timer, only: getTime
6 use pm_timer, only: setIdle
7
8 implicit none
9
10 real(RKD) :: since
11
12 type(display_type) :: disp
13 disp = display_type(file = "main.out.F90")
14
15 call disp%skip()
16 call disp%show("since = getTime()")
17 since = getTime()
18 call disp%show("call setIdle(seconds = 0.1_RKD)")
19 call setIdle(seconds = 0.1_RKD)
20 call disp%show("getTime(since)")
21 call disp%show( getTime(since) )
22 call disp%skip()
23
24end program example
subroutine setIdle(seconds)
Set the processor in sleep mode for the input requested number of seconds via the default timer of th...
Definition: pm_timer.F90:1593

Example Unix compile command via Intel ifort compiler
1#!/usr/bin/env sh
2rm main.exe
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example Windows Batch compile command via Intel ifort compiler
1del main.exe
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
4main.exe

Example Unix / MinGW compile command via GNU gfortran compiler
1#!/usr/bin/env sh
2rm main.exe
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example output
1
2since = getTime()
3call setIdle(seconds = 0.1_RKD)
4getTime(since)
5+0.10002790705766529
6
7
Test:
test_pm_timer


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas Austin

Definition at line 857 of file pm_timer.F90.

References getTimeMPI(), getTimeOMP(), and getTimeSYS().

Here is the call graph for this function:

◆ getTimeCPU()

real(RKD) function pm_timer::getTimeCPU ( real(RKD), intent(in), optional  since)

Generate and return the CPU time in units of seconds since the specified time since or since an arbitrary processor-dependent time.

See also the documentation details of pm_timer.

Parameters
[in]since: The input scalar of type real of kind double precision RKD representing the time origin (epoch) in seconds.
(optional. If not present, its value is processor-dependent but constant at runtime.)
Returns
timeInSec : The output scalar of type real of kind double precision RKD, containing the time in seconds since a specified time since or since an arbitrary processor-dependent time.


Possible calling interfaces

use pm_timer, only: getTimeCPU
real(RKD) :: timeInSec
logical(LK) :: failed
timeInSec = getTimeCPU(since = since)
real(RKD) function getTimeCPU(since)
Generate and return the CPU time in units of seconds since the specified time since or since an arbit...
Definition: pm_timer.F90:918
Remarks
The procedures under discussion are impure.
See also
timerCPU_type


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK, RKD
4 use pm_io, only: display_type
5 use pm_timer, only: getTimeCPU
6 use pm_timer, only: setIdleCPU
7
8 implicit none
9
10 real(RKD) :: since
11
12 type(display_type) :: disp
13 disp = display_type(file = "main.out.F90")
14
15 call disp%skip()
16 call disp%show("since = getTimeCPU()")
17 since = getTimeCPU()
18 call disp%show("call setIdleCPU(seconds = 0.1_RKD)")
19 call setIdleCPU(seconds = 0.1_RKD)
20 call disp%show("getTimeCPU(since)")
21 call disp%show( getTimeCPU(since) )
22 call disp%skip()
23
24end program example
subroutine setIdleCPU(seconds)
Set the processor in sleep mode for the input requested number of seconds via the Fortran intrinsic c...
Definition: pm_timer.F90:1658

Example Unix compile command via Intel ifort compiler
1#!/usr/bin/env sh
2rm main.exe
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example Windows Batch compile command via Intel ifort compiler
1del main.exe
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
4main.exe

Example Unix / MinGW compile command via GNU gfortran compiler
1#!/usr/bin/env sh
2rm main.exe
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example output
1
2since = getTimeCPU()
3call setIdleCPU(seconds = 0.1_RKD)
4getTimeCPU(since)
5+0.10003200000000000
6
7
Test:
test_pm_timer


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas Austin

Definition at line 917 of file pm_timer.F90.

Referenced by setIdleCPU().

Here is the caller graph for this function:

◆ getTimeDAT()

real(RKD) function pm_timer::getTimeDAT ( real(RKD), intent(in), optional  since)

Generate and return the calendrical clock time in units of seconds since the specified time since or since the Gregorian calendrical time origin.

See also the documentation details of pm_timer.

Parameters
[in]since: The input scalar of type real of kind double precision RKD representing the time origin (epoch) in seconds.
(optional. If not present, its value is the Gregorian calendrical time origin.)
Returns
timeInSec : The output scalar of type real of kind double precision RKD, containing the time in seconds since the specified time origin.


Possible calling interfaces

use pm_timer, only: getTimeDAT
real(RKD) :: timeInSec
timeInSec = getTimeDAT(since = since)
real(RKD) function getTimeDAT(since)
Generate and return the calendrical clock time in units of seconds since the specified time since or ...
Definition: pm_timer.F90:973
Remarks
The procedures under discussion are impure.
See also
timerDAT_type


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK, RKD
4 use pm_io, only: display_type
5 use pm_timer, only: getTimeDAT
6 use pm_timer, only: setIdleDAT
7
8 implicit none
9
10 real(RKD) :: since
11
12 type(display_type) :: disp
13 disp = display_type(file = "main.out.F90")
14
15 call disp%skip()
16 call disp%show("since = getTimeDAT()")
17 since = getTimeDAT()
18 call disp%show("call setIdleDAT(seconds = 0.1_RKD)")
19 call setIdleDAT(seconds = 0.1_RKD)
20 call disp%show("getTimeDAT(since)")
21 call disp%show( getTimeDAT(since) )
22 call disp%skip()
23
24end program example
subroutine setIdleDAT(seconds)
Set the processor in sleep mode for the input requested number of seconds via the Fortran intrinsic d...
Definition: pm_timer.F90:1721

Example Unix compile command via Intel ifort compiler
1#!/usr/bin/env sh
2rm main.exe
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example Windows Batch compile command via Intel ifort compiler
1del main.exe
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
4main.exe

Example Unix / MinGW compile command via GNU gfortran compiler
1#!/usr/bin/env sh
2rm main.exe
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example output
1
2since = getTimeDAT()
3call setIdleDAT(seconds = 0.1_RKD)
4getTimeDAT(since)
5+0.10098266601562500
6
7
Test:
test_pm_timer


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas Austin

Definition at line 972 of file pm_timer.F90.

References pm_kind::IK, pm_kind::RKD, and pm_dateTime::SECONDS_PER_DAY.

Referenced by setIdleDAT().

Here is the caller graph for this function:

◆ getTimeMPI()

real(RKD) function pm_timer::getTimeMPI ( real(RKD), intent(in), optional  since)

Generate and return the MPI clock time in units of seconds since the specified time since or since an arbitrary time origin set by the MPI library.

See also the documentation details of pm_timer.

Parameters
[in]since: The input scalar of type real of kind double precision RKD representing the time origin (epoch) in seconds.
(optional. If not present, its value is an arbitrary time origin set by the MPI library.)
Returns
timeInSec : The output scalar of type real of kind double precision RKD, containing the time in seconds since the specified time origin.


Possible calling interfaces

use pm_timer, only: getTimeMPI
real(RKD) :: timeInSec
timeInSec = getTimeMPI(since = since)
real(RKD) function getTimeMPI(since)
Generate and return the MPI clock time in units of seconds since the specified time since or since an...
Definition: pm_timer.F90:1036
Warning
This procedure exists only if the library is built with the preprocessor macro MPI_ENABLED=1.
Remarks
The procedures under discussion are impure.
See also
timerMPI_type


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK, RKD
4 use pm_io, only: display_type
5
6 implicit none
7
8 real(RKD) :: since
9
10 type(display_type) :: disp
11 disp = display_type(file = "main.out.F90")
12
13#if defined MPI_VERSION
14 block
15 use pm_timer, only: getTimeMPI
16 use pm_timer, only: setIdleMPI
17 call disp%skip()
18 call disp%show("since = getTimeMPI()")
19 since = getTimeMPI()
20 call disp%show("call setIdleMPI(seconds = 0.1_RKD)")
21 call setIdleMPI(seconds = 0.1_RKD)
22 call disp%show("getTimeMPI(since)")
23 call disp%show( getTimeMPI(since) )
24 call disp%skip()
25 end block
26#else
27 call disp%show("The MPI timer is unavailable, because this example has not been linked with an MPI library.")
28#endif
29
30end program example
subroutine setIdleMPI(seconds)
Set the processor in sleep mode for the input requested number of seconds via the MPI intrinsic timer...
Definition: pm_timer.F90:1787

Example Unix compile command via Intel ifort compiler
1#!/usr/bin/env sh
2rm main.exe
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example Windows Batch compile command via Intel ifort compiler
1del main.exe
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
4main.exe

Example Unix / MinGW compile command via GNU gfortran compiler
1#!/usr/bin/env sh
2rm main.exe
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example output
1The MPI timer is unavailable, because this example has not been linked with an MPI library.
2
Test:
test_pm_timer


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas Austin

Definition at line 1035 of file pm_timer.F90.

References pm_kind::IK, and pm_kind::RKD.

Referenced by getTime(), and setIdleMPI().

Here is the caller graph for this function:

◆ getTimeOMP()

real(RKD) function pm_timer::getTimeOMP ( real(RKD), intent(in), optional  since)

Generate and return the OpenMP clock time in units of seconds since the specified time since or since an arbitrary time origin set by the OpenMP library.

See also the documentation details of pm_timer.

Parameters
[in]since: The input scalar of type real of kind double precision RKD representing the time origin (epoch) in seconds.
(optional. If not present, its value is an arbitrary time origin set by the OpenMP library.)
Returns
timeInSec : The output scalar of type real of kind double precision RKD, containing the time in seconds since the specified time origin.


Possible calling interfaces

use pm_timer, only: getTimeOMP
real(RKD) :: timeInSec
timeInSec = getTimeOMP(since = since)
real(RKD) function getTimeOMP(since)
Generate and return the OpenMP clock time in units of seconds since the specified time since or since...
Definition: pm_timer.F90:1097
Warning
This procedure exists only if the library is built with the preprocessor macro OMP_ENABLED=1.
Remarks
The procedures under discussion are impure.
See also
timerOMP_type


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK, RKD
4 use pm_io, only: display_type
5
6 implicit none
7
8 real(RKD) :: since
9
10 type(display_type) :: disp
11 disp = display_type(file = "main.out.F90")
12
13#if defined _OPENMP
14 block
15 use pm_timer, only: getTimeOMP
16 use pm_timer, only: setIdleOMP
17 call disp%skip()
18 call disp%show("since = getTimeOMP()")
19 since = getTimeOMP()
20 call disp%show("call setIdleOMP(seconds = 0.1_RKD)")
21 call setIdleOMP(seconds = 0.1_RKD)
22 call disp%show("getTimeOMP(since)")
23 call disp%show( getTimeOMP(since) )
24 call disp%skip()
25 end block
26#else
27 call disp%show("The OMP timer is unavailable, because this example has not been linked with an OpenMP library.")
28#endif
29
30end program example
subroutine setIdleOMP(seconds)
Set the processor in sleep mode for the input requested number of seconds via the OpenMP intrinsic ti...
Definition: pm_timer.F90:1850

Example Unix compile command via Intel ifort compiler
1#!/usr/bin/env sh
2rm main.exe
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example Windows Batch compile command via Intel ifort compiler
1del main.exe
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
4main.exe

Example Unix / MinGW compile command via GNU gfortran compiler
1#!/usr/bin/env sh
2rm main.exe
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example output
1The OMP timer is unavailable, because this example has not been linked with an OpenMP library.
2
Test:
test_pm_timer


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas Austin

Definition at line 1096 of file pm_timer.F90.

References pm_kind::IK, and pm_kind::RKD.

Referenced by getTime(), and setIdleOMP().

Here is the caller graph for this function:

◆ getTimeSYS()

real(RKD) function pm_timer::getTimeSYS ( real(RKD), intent(in), optional  since)

Generate and return the system clock time in units of seconds since the specified time since or since an arbitrary processor-dependent time origin.

See also the documentation details of pm_timer.

Parameters
[in]since: The input scalar of type real of kind double precision RKD representing the time origin (epoch) in seconds.
(optional. If not present, its value is processor-dependent but constant at runtime.)
Returns
timeInSec : The output scalar of type real of kind double precision RKD, containing the time in seconds since a processor-defined origin of time.


Possible calling interfaces

use pm_timer, only: getTimeSYS
real(RKD) :: timeInSec
timeInSec = getTimeSYS(since = since)
real(RKD) function getTimeSYS(since)
Generate and return the system clock time in units of seconds since the specified time since or since...
Definition: pm_timer.F90:1198
Remarks
The procedures under discussion are impure.
See also
timerSYS_type


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK, RKD
4 use pm_io, only: display_type
5 use pm_timer, only: getTimeSYS
6 use pm_timer, only: setIdleSYS
7
8 implicit none
9
10 real(RKD) :: since
11
12 type(display_type) :: disp
13 disp = display_type(file = "main.out.F90")
14
15 call disp%skip()
16 call disp%show("since = getTimeSYS()")
17 since = getTimeSYS()
18 call disp%show("call setIdleSYS(seconds = 0.1_RKD)")
19 call setIdleSYS(seconds = 0.1_RKD)
20 call disp%show("getTimeSYS(since)")
21 call disp%show( getTimeSYS(since) )
22 call disp%skip()
23
24end program example
subroutine setIdleSYS(seconds)
Set the processor in sleep mode for the input requested number of seconds via the Fortran intrinsic s...
Definition: pm_timer.F90:1914

Example Unix compile command via Intel ifort compiler
1#!/usr/bin/env sh
2rm main.exe
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example Windows Batch compile command via Intel ifort compiler
1del main.exe
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
4main.exe

Example Unix / MinGW compile command via GNU gfortran compiler
1#!/usr/bin/env sh
2rm main.exe
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example output
1
2since = getTimeSYS()
3call setIdleSYS(seconds = 0.1_RKD)
4getTimeSYS(since)
5+0.10003740899264812
6
7
Test:
test_pm_timer
Todo:
Normal Priority: The performance of this procedure could be improved by avoiding the unnecessary retrieval of the count rate and computing its inverse to get the period (which is fixed in the entire runtime).
This requires redefining this procedure as a type-bound procedure of a parent timer type.
However, preliminary benchmarks indicate that the division (instead of multiplication) slows down the computation by at most 25 percent. On modern architecture, the extra cost is likely on the of 25 CPU cycles, approximately 25 nanoseconds. This is likely negligible in most practical scenarios. Here is a benchmark script for the cost of division vs. multiplication,
implicit none
double precision :: t0, t1, x(2), summ
integer :: v0(8), v1(8)
integer :: i
call cpu_time(t0)
!call date_and_time(values = v0)
do i = 1, 10**7
call random_number(x)
summ = summ + x(2) * x(1)
end do
call cpu_time(t1)
print *, summ, "prod", t1 - t0
!call date_and_time(values = v1)
!print *, summ, "prod", v1(7) - v0(7) + (v1(8) - v0(8)) / 1.d3
summ = 0.d0
call cpu_time(t0)
!call date_and_time(values = v0)
do i = 1, 10**7
call random_number(x)
summ = summ + x(2) / x(1)
end do
!call date_and_time(values = v1)
call cpu_time(t1)
print *, summ, t1 - t0
!print *, summ, "divi", v1(7) - v0(7) + (v1(8) - v0(8)) / 1.d3
end

The above code yields highly similar time costs for multiplication vs. division operation in most benchmarks.


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas Austin

Definition at line 1197 of file pm_timer.F90.

Referenced by getTime(), and setIdleSYS().

Here is the caller graph for this function:

◆ setIdle()

subroutine pm_timer::setIdle ( real(RKD), intent(in)  seconds)

Set the processor in sleep mode for the input requested number of seconds via the default timer of the ParaMonte library.

This function relies on either of the following lower-level module functions depending on the ParaMonte library build:

  1. setIdleMPI for MPI-enabled parallel builds of the ParaMonte library.
  2. setIdleOMP for OpenMP-enabled parallel builds of the ParaMonte library.
  3. setIdleSYS for all other parallel or serial builds of the ParaMonte library.
Parameters
[in]seconds: The input scalar of type real of kind double precision RKD representing the number of seconds to put the system into sleep.


Possible calling interfaces

use pm_timer, only: setIdle
real(RKD) :: seconds
call setIdle(seconds)
Remarks
The procedures under discussion are impure.
See also
setIdleCPU
setIdleDAT
setIdleMPI
setIdleOMP
setIdleSYS
getTimeCPU
getTimeDAT
getTimeMPI
getTimeOMP
getTimeSYS


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK, RKD
4 use pm_io, only: display_type
5 use pm_timer, only: getTime
6 use pm_timer, only: setIdle
7
8 implicit none
9
10 real(RKD) :: since
11
12 type(display_type) :: disp
13 disp = display_type(file = "main.out.F90")
14
15 call disp%skip()
16 call disp%show("since = getTime()")
17 since = getTime()
18 call disp%show("call setIdle(seconds = 0.2_RKD)")
19 call setIdle(seconds = 0.2_RKD)
20 call disp%show("getTime(since)")
21 call disp%show( getTime(since) )
22 call disp%skip()
23
24end program example

Example Unix compile command via Intel ifort compiler
1#!/usr/bin/env sh
2rm main.exe
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example Windows Batch compile command via Intel ifort compiler
1del main.exe
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
4main.exe

Example Unix / MinGW compile command via GNU gfortran compiler
1#!/usr/bin/env sh
2rm main.exe
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example output
1
2since = getTime()
3call setIdle(seconds = 0.2_RKD)
4getTime(since)
5+0.20003350905608386
6
7
Test:
test_pm_timer


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas Austin

Definition at line 1592 of file pm_timer.F90.

References pm_kind::RKD, setIdleMPI(), setIdleOMP(), and setIdleSYS().

Here is the call graph for this function:

◆ setIdleCPU()

subroutine pm_timer::setIdleCPU ( real(RKD), intent(in)  seconds)

Set the processor in sleep mode for the input requested number of seconds via the Fortran intrinsic cpu_time().

Parameters
[in]seconds: The input scalar of type real of kind double precision RKD representing the number of seconds to put the system into sleep.


Possible calling interfaces

use pm_timer, only: setIdleCPU
real(RKD) :: seconds
call setIdleCPU(seconds)
Warning
If the system does not possess a clock, the results are unpredictable and the procedure may enter an indefinite loop.
However, this condition rarely occurs and can be readily check only once on each new platform before a production run.
The lack of a system clock will lead to runtime error only if the library is built with the preprocessor macro CHECK_ENABLED=1.
Remarks
The procedures under discussion are impure.
See also
setIdleCPU
setIdleDAT
setIdleMPI
setIdleOMP
setIdleSYS
getTimeCPU
getTimeDAT
getTimeMPI
getTimeOMP
getTimeSYS


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK, RKD
4 use pm_io, only: display_type
5 use pm_timer, only: getTimeCPU
6 use pm_timer, only: setIdleCPU
7
8 implicit none
9
10 real(RKD) :: since
11
12 type(display_type) :: disp
13 disp = display_type(file = "main.out.F90")
14
15 call disp%skip()
16 call disp%show("since = getTimeCPU()")
17 since = getTimeCPU()
18 call disp%show("call setIdleCPU(seconds = 0.1_RKD)")
19 call setIdleCPU(seconds = 0.1_RKD)
20 call disp%show("getTimeCPU(since)")
21 call disp%show( getTimeCPU(since) )
22 call disp%skip()
23
24end program example

Example Unix compile command via Intel ifort compiler
1#!/usr/bin/env sh
2rm main.exe
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example Windows Batch compile command via Intel ifort compiler
1del main.exe
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
4main.exe

Example Unix / MinGW compile command via GNU gfortran compiler
1#!/usr/bin/env sh
2rm main.exe
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example output
1
2since = getTimeCPU()
3call setIdleCPU(seconds = 0.1_RKD)
4getTimeCPU(since)
5+0.10003300000000000
6
7
Test:
test_pm_timer


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas Austin

Definition at line 1657 of file pm_timer.F90.

References getTimeCPU(), and pm_kind::RKD.

Here is the call graph for this function:

◆ setIdleDAT()

subroutine pm_timer::setIdleDAT ( real(RKD), intent(in)  seconds)

Set the processor in sleep mode for the input requested number of seconds via the Fortran intrinsic date_and_time().

Parameters
[in]seconds: The input scalar of type real of kind double precision RKD representing the number of seconds to put the system into sleep.


Possible calling interfaces

use pm_timer, only: setIdleDAT
real(RKD) :: seconds
call setIdleDAT(seconds)
Warning
If the system does not possess a clock, the results are unpredictable and the procedure may enter an indefinite loop.
However, this condition rarely occurs and can be readily check only once on each new platform before a production run.
The lack of a system clock will lead to runtime error only if the library is built with the preprocessor macro CHECK_ENABLED=1.
Remarks
The procedures under discussion are impure.
See also
setIdleCPU
setIdleDAT
setIdleMPI
setIdleOMP
setIdleSYS
getTimeCPU
getTimeDAT
getTimeMPI
getTimeOMP
getTimeSYS


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK, RKD
4 use pm_io, only: display_type
5 use pm_timer, only: getTimeDAT
6 use pm_timer, only: setIdleDAT
7
8 implicit none
9
10 real(RKD) :: since
11
12 type(display_type) :: disp
13 disp = display_type(file = "main.out.F90")
14
15 call disp%skip()
16 call disp%show("since = getTimeDAT()")
17 since = getTimeDAT()
18 call disp%show("call setIdleDAT(seconds = 0.2_RKD)")
19 call setIdleDAT(seconds = 0.2_RKD)
20 call disp%show("getTimeDAT(since)")
21 call disp%show( getTimeDAT(since) )
22 call disp%skip()
23
24end program example

Example Unix compile command via Intel ifort compiler
1#!/usr/bin/env sh
2rm main.exe
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example Windows Batch compile command via Intel ifort compiler
1del main.exe
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
4main.exe

Example Unix / MinGW compile command via GNU gfortran compiler
1#!/usr/bin/env sh
2rm main.exe
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example output
1
2since = getTimeDAT()
3call setIdleDAT(seconds = 0.2_RKD)
4getTimeDAT(since)
5+0.20001220703125000
6
7
Test:
test_pm_timer


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas Austin

Definition at line 1720 of file pm_timer.F90.

References getTimeDAT(), and pm_kind::RKD.

Here is the call graph for this function:

◆ setIdleMPI()

subroutine pm_timer::setIdleMPI ( real(RKD), intent(in)  seconds)

Set the processor in sleep mode for the input requested number of seconds via the MPI intrinsic timer MPI_Wtime().

Parameters
[in]seconds: The input scalar of type real of kind double precision RKD representing the number of seconds to put the system into sleep.


Possible calling interfaces

use pm_timer, only: setIdleMPI
real(RKD) :: seconds
call setIdleMPI(seconds)
Warning
This procedure exists only if the library is built with the preprocessor macro MPI_ENABLED=1.
Calling this procedure requires the MPI library to have been initialized and not finalized.
This condition is verified only if the library is built with the preprocessor macro CHECK_ENABLED=1.
Remarks
The procedures under discussion are impure.
See also
setIdleCPU
setIdleDAT
setIdleMPI
setIdleOMP
setIdleSYS
getTimeCPU
getTimeDAT
getTimeMPI
getTimeOMP
getTimeSYS


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK, RKD
4 use pm_io, only: display_type
5
6 implicit none
7
8
9 type(display_type) :: disp
10 disp = display_type(file = "main.out.F90")
11
12#if defined MPI_VERSION
13 block
14 use pm_timer, only: getTimeMPI
15 use pm_timer, only: setIdleMPI
16 real(RKD) :: since
17 call disp%skip()
18 call disp%show("since = getTimeMPI()")
19 since = getTimeMPI()
20 call disp%show("call setIdleMPI(seconds = 0.2_RKD)")
21 call setIdleMPI(seconds = 0.2_RKD)
22 call disp%show("getTimeMPI(since)")
23 call disp%show( getTimeMPI(since) )
24 call disp%skip()
25 end block
26#else
27 call disp%show("The MPI timer is unavailable, because this example has not been linked with an MPI library.")
28#endif
29
30end program example

Example Unix compile command via Intel ifort compiler
1#!/usr/bin/env sh
2rm main.exe
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example Windows Batch compile command via Intel ifort compiler
1del main.exe
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
4main.exe

Example Unix / MinGW compile command via GNU gfortran compiler
1#!/usr/bin/env sh
2rm main.exe
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example output
1The MPI timer is unavailable, because this example has not been linked with an MPI library.
2
Test:
test_pm_timer


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas Austin

Definition at line 1786 of file pm_timer.F90.

References getTimeMPI(), and pm_kind::RKD.

Referenced by setIdle().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ setIdleOMP()

subroutine pm_timer::setIdleOMP ( real(RKD), intent(in)  seconds)

Set the processor in sleep mode for the input requested number of seconds via the OpenMP intrinsic timer omp_get_wtime().

Parameters
[in]seconds: The input scalar of type real of kind double precision RKD representing the number of seconds to put the system into sleep.


Possible calling interfaces

use pm_timer, only: setIdleOMP
real(RKD) :: seconds
call setIdleOMP(seconds)
Warning
This procedure exists only if the library is built with the preprocessor macro OMP_ENABLED=1.
Remarks
The procedures under discussion are impure.
See also
setIdleCPU
setIdleDAT
setIdleMPI
setIdleOMP
setIdleSYS
getTimeCPU
getTimeDAT
getTimeMPI
getTimeOMP
getTimeSYS


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK, RKD
4 use pm_io, only: display_type
5
6 implicit none
7
8 type(display_type) :: disp
9 disp = display_type(file = "main.out.F90")
10
11#if defined _OPENMP
12 block
13 use pm_timer, only: getTimeOMP
14 use pm_timer, only: setIdleOMP
15 real(RKD) :: since
16 call disp%skip()
17 call disp%show("since = getTimeOMP()")
18 since = getTimeOMP()
19 call disp%show("call setIdleOMP(seconds = 0.2_RKD)")
20 call setIdleOMP(seconds = 0.2_RKD)
21 call disp%show("getTimeOMP(since)")
22 call disp%show( getTimeOMP(since) )
23 call disp%skip()
24 end block
25#else
26 call disp%show("The OMP timer is unavailable, because this example has not been linked with an OpenMP library.")
27#endif
28
29end program example

Example Unix compile command via Intel ifort compiler
1#!/usr/bin/env sh
2rm main.exe
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example Windows Batch compile command via Intel ifort compiler
1del main.exe
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
4main.exe

Example Unix / MinGW compile command via GNU gfortran compiler
1#!/usr/bin/env sh
2rm main.exe
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example output
1The OMP timer is unavailable, because this example has not been linked with an OpenMP library.
2
Test:
test_pm_timer


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas Austin

Definition at line 1849 of file pm_timer.F90.

References getTimeOMP(), and pm_kind::RKD.

Referenced by setIdle().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ setIdleSYS()

subroutine pm_timer::setIdleSYS ( real(RKD), intent(in)  seconds)

Set the processor in sleep mode for the input requested number of seconds via the Fortran intrinsic system_clock().

Parameters
[in]seconds: The input scalar of type real of kind double precision RKD representing the number of seconds to put the system into sleep.


Possible calling interfaces

use pm_timer, only: setIdleSYS
real(RKD) :: seconds
call setIdleSYS(seconds)
Warning
If the system does not possess a clock, the results are unpredictable and the procedure may enter an indefinite loop.
However, this condition rarely occurs and can be readily check only once on each new platform before a production run.
The lack of a system clock will lead to runtime error only if the library is built with the preprocessor macro CHECK_ENABLED=1.
Remarks
The procedures under discussion are impure.
See also
setIdleCPU
setIdleDAT
setIdleMPI
setIdleOMP
setIdleSYS
getTimeCPU
getTimeDAT
getTimeMPI
getTimeOMP
getTimeSYS


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK, RKD
4 use pm_io, only: display_type
5 use pm_timer, only: getTimeSYS
6 use pm_timer, only: setIdleSYS
7
8 implicit none
9
10 real(RKD) :: since
11
12 type(display_type) :: disp
13 disp = display_type(file = "main.out.F90")
14
15 call disp%skip()
16 call disp%show("since = getTimeSYS()")
17 since = getTimeSYS()
18 call disp%show("call setIdleSYS(seconds = 0.2_RKD)")
19 call setIdleSYS(seconds = 0.2_RKD)
20 call disp%show("getTimeSYS(since)")
21 call disp%show( getTimeSYS(since) )
22 call disp%skip()
23
24end program example

Example Unix compile command via Intel ifort compiler
1#!/usr/bin/env sh
2rm main.exe
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example Windows Batch compile command via Intel ifort compiler
1del main.exe
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
4main.exe

Example Unix / MinGW compile command via GNU gfortran compiler
1#!/usr/bin/env sh
2rm main.exe
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example output
1
2since = getTimeSYS()
3call setIdleSYS(seconds = 0.2_RKD)
4getTimeSYS(since)
5+0.20003396400716156
6
7
Test:
test_pm_timer


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas Austin

Definition at line 1913 of file pm_timer.F90.

References getTimeSYS(), and pm_kind::RKD.

Referenced by setIdle().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ timer_typer()

type(timerMPI_type) type(timerOMP_type) type(timerSYS_type) function pm_timer::timer_typer

This is the constructor of the timer_type class.

Upon return, the constructor initializes all components of the timer object.
See also the documentation details of pm_timer.

Returns
timer : The output scalar object of class timer_type.
Specifically, the output of type,
  1. timerMPI_type if the ParaMonte library is built with preprocessor macro MPI_ENABLED defined.
  2. timerOMP_type if the ParaMonte library is built with preprocessor macro OMP_ENABLED defined.
  3. timerSYS_type if the ParaMonte library is built is serial mode.


Possible calling interfaces

use pm_timer, only: timer_type
class(timer_type), allocatable :: timer
timer = timer_type()
This is the abstract base derived type that serves as a simple container template for other timer cla...
Definition: pm_timer.F90:212
Remarks
See the documentation of timer_type for example usage.


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas Austin

Definition at line 556 of file pm_timer.F90.

◆ timerCPU_typer()

type(timerCPU_type) function pm_timer::timerCPU_typer

This is the constructor of the timerCPU_type class.

Upon return, the constructor initializes all components of the timer object.
See also the documentation details of pm_timer.

Returns
timer : The output scalar object of class timerCPU_type.


Possible calling interfaces

type(timerCPU_type) :: timer
timer = timerCPU_type()
This is the timerCPU_type class, containing attributes and static methods for setting up a timer base...
Definition: pm_timer.F90:271
Remarks
See the documentation of timerCPU_type for example usage.


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas Austin

Definition at line 601 of file pm_timer.F90.

◆ timerDAT_typer()

type(timerDAT_type) function pm_timer::timerDAT_typer

This is the constructor of the timerDAT_type class.

Upon return, the constructor initializes all components of the timer object.
See also the documentation details of pm_timer.


Possible calling interfaces

type(timerDAT_type) :: timer
timer = timerDAT_type()
This is the timerDAT_type class, containing attributes and static methods for setting up a timer base...
Definition: pm_timer.F90:322
Returns
timer : The output scalar object of class timerDAT_type.
Remarks
See the documentation of timerDAT_type for example usage.


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas Austin

Definition at line 641 of file pm_timer.F90.

◆ timerMPI_typer()

type(timerMPI_type) function pm_timer::timerMPI_typer

This is the constructor of the timerMPI_type class.

Upon return, the constructor initializes all components of the timer object.
See also the documentation details of pm_timer.


Possible calling interfaces

type(timerMPI_type) :: timer
timer = timerMPI_type()
This is the timerMPI_type class, containing attributes and static methods for setting up a timer base...
Definition: pm_timer.F90:377
Returns
timer : The output scalar object of class timerMPI_type.
Warning
This constructor requires the MPI library to have been initialized and not finalized.
This condition is verified only if the library is built with the preprocessor macro CHECK_ENABLED=1.
This constructor is available only if the library is built with the preprocessor macro MPI_ENABLED=1.
Remarks
See the documentation of timerMPI_type for example usage.


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, March 22, 2012, 02:51 AM, National Institute for Fusion Studies, The University of Texas Austin

Definition at line 689 of file pm_timer.F90.

◆ timerOMP_typer()

type(timerOMP_type) function pm_timer::timerOMP_typer

This is the constructor of the timerOMP_type class.

Upon return, the constructor initializes all components of the timer object.
See also the documentation details of pm_timer.


Possible calling interfaces

type(timerOMP_type) :: timer
timer = timerOMP_type()
This is the timerMPI_type class, containing attributes and static methods for setting up a timer base...
Definition: pm_timer.F90:433
Returns
timer : The output scalar object of class timerOMP_type.
Warning
This constructor is available only if the library is built with the preprocessor macro OMP_ENABLED=1.
Remarks
See the documentation of timerOMP_type for example usage.


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, March 22, 2012, 00:00 AM, National Institute for Fusion Studies, The University of Texas Austin

Definition at line 752 of file pm_timer.F90.

◆ timerSYS_typer()

type(timerSYS_type) function pm_timer::timerSYS_typer

This is the constructor of the timerSYS_type class.

Upon return, the constructor initializes all components of the timer object.
See also the documentation details of pm_timer.


Possible calling interfaces

type(timerSYS_type) :: timer
timer = timerSYS_type()
This is the timerSYS_type class, containing attributes and static methods for setting up a timer base...
Definition: pm_timer.F90:502
Returns
timer : The output scalar object of class timerSYS_type.
Remarks
See the documentation of timerSYS_type for example usage.


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, March 22, 2012, 03:00 AM, National Institute for Fusion Studies, The University of Texas Austin

Definition at line 794 of file pm_timer.F90.

Variable Documentation

◆ MODULE_NAME

character(*, SK), parameter pm_timer::MODULE_NAME = "@pm_timer"

Definition at line 117 of file pm_timer.F90.