2016-09-21 01:08:18 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2018-08-20 19:29:13 +05:30
|
|
|
|
!> @author Jaeyong Jung, Max-Planck-Institut für Eisenforschung GmbH
|
|
|
|
|
!> @author Pratheek Shanthraj, Max-Planck-Institut für Eisenforschung GmbH
|
2016-09-21 01:08:18 +05:30
|
|
|
|
!> @author Martin Diehl, Max-Planck-Institut für Eisenforschung GmbH
|
2013-01-02 22:32:12 +05:30
|
|
|
|
!> @author Philip Eisenlohr, Max-Planck-Institut für Eisenforschung GmbH
|
2018-08-20 19:29:13 +05:30
|
|
|
|
!> @brief Interfacing between the PETSc-based solvers and the material subroutines provided
|
2013-01-02 22:32:12 +05:30
|
|
|
|
!! by DAMASK
|
2018-08-20 19:29:13 +05:30
|
|
|
|
!> @details Interfacing between the PETSc-based solvers and the material subroutines provided
|
2019-03-06 20:17:48 +05:30
|
|
|
|
!> by DAMASK. Interpretating the command line arguments to get load case, geometry file,
|
2018-08-05 14:11:01 +05:30
|
|
|
|
!> and working directory.
|
2014-10-10 18:38:34 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-04-03 11:55:42 +05:30
|
|
|
|
#define GCC_MIN 6
|
2020-01-03 17:10:25 +05:30
|
|
|
|
#define INTEL_MIN 1700
|
2019-04-03 11:55:42 +05:30
|
|
|
|
#define PETSC_MAJOR 3
|
|
|
|
|
#define PETSC_MINOR_MIN 10
|
2019-10-07 23:53:17 +05:30
|
|
|
|
#define PETSC_MINOR_MAX 12
|
2019-05-28 15:36:21 +05:30
|
|
|
|
|
2012-03-06 20:22:48 +05:30
|
|
|
|
module DAMASK_interface
|
2019-05-28 15:36:21 +05:30
|
|
|
|
use, intrinsic :: iso_fortran_env
|
|
|
|
|
use PETScSys
|
|
|
|
|
|
2019-05-11 15:40:23 +05:30
|
|
|
|
use prec
|
|
|
|
|
use system_routines
|
|
|
|
|
|
2019-03-06 20:17:48 +05:30
|
|
|
|
implicit none
|
|
|
|
|
private
|
|
|
|
|
logical, public, protected :: &
|
2019-03-24 16:29:00 +05:30
|
|
|
|
SIGTERM, & !< termination signal
|
|
|
|
|
SIGUSR1, & !< 1. user-defined signal
|
|
|
|
|
SIGUSR2 !< 2. user-defined signal
|
2019-03-06 20:17:48 +05:30
|
|
|
|
integer, public, protected :: &
|
|
|
|
|
interface_restartInc = 0 !< Increment at which calculation starts
|
|
|
|
|
character(len=1024), public, protected :: &
|
|
|
|
|
geometryFile = '', & !< parameter given for geometry file
|
|
|
|
|
loadCaseFile = '' !< parameter given for load case file
|
|
|
|
|
|
|
|
|
|
public :: &
|
|
|
|
|
getSolverJobName, &
|
2019-03-24 16:29:00 +05:30
|
|
|
|
DAMASK_interface_init, &
|
|
|
|
|
setSIGTERM, &
|
|
|
|
|
setSIGUSR1, &
|
|
|
|
|
setSIGUSR2
|
2019-05-09 11:55:56 +05:30
|
|
|
|
|
2012-03-06 20:22:48 +05:30
|
|
|
|
contains
|
|
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2012-03-09 01:55:28 +05:30
|
|
|
|
!> @brief initializes the solver by interpreting the command line arguments. Also writes
|
2012-06-15 21:40:21 +05:30
|
|
|
|
!! information on computation to screen
|
2012-03-06 20:22:48 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-05-03 09:51:43 +05:30
|
|
|
|
subroutine DAMASK_interface_init
|
2018-05-17 15:34:21 +05:30
|
|
|
|
#include <petsc/finclude/petscsys.h>
|
2019-04-03 18:19:16 +05:30
|
|
|
|
#if defined(__GFORTRAN__) && __GNUC__<GCC_MIN
|
2018-10-15 08:33:53 +05:30
|
|
|
|
===================================================================================================
|
2019-04-03 11:55:42 +05:30
|
|
|
|
----- WRONG COMPILER VERSION ----- WRONG COMPILER VERSION ----- WRONG COMPILER VERSION -----
|
2018-10-15 08:33:53 +05:30
|
|
|
|
===================================================================================================
|
2019-04-03 11:55:42 +05:30
|
|
|
|
=============== THIS VERSION OF DAMASK REQUIRES A NEWER gfortran VERSION ======================
|
|
|
|
|
================== THIS VERSION OF DAMASK REQUIRES A NEWER gfortran VERSION ===================
|
|
|
|
|
===================== THIS VERSION OF DAMASK REQUIRES A NEWER gfortran VERSION ================
|
2018-10-15 08:33:53 +05:30
|
|
|
|
===================================================================================================
|
2019-04-03 11:55:42 +05:30
|
|
|
|
----- WRONG COMPILER VERSION ----- WRONG COMPILER VERSION ----- WRONG COMPILER VERSION -----
|
2018-10-15 08:33:53 +05:30
|
|
|
|
===================================================================================================
|
|
|
|
|
#endif
|
|
|
|
|
|
2019-04-03 18:19:16 +05:30
|
|
|
|
#if defined(__INTEL_COMPILER) && __INTEL_COMPILER<INTEL_MIN
|
2018-10-15 08:33:53 +05:30
|
|
|
|
===================================================================================================
|
2019-04-03 11:55:42 +05:30
|
|
|
|
----- WRONG COMPILER VERSION ----- WRONG COMPILER VERSION ----- WRONG COMPILER VERSION -----
|
2018-10-15 08:33:53 +05:30
|
|
|
|
===================================================================================================
|
2019-04-03 11:55:42 +05:30
|
|
|
|
================= THIS VERSION OF DAMASK REQUIRES A NEWER ifort VERSION =======================
|
|
|
|
|
==================== THIS VERSION OF DAMASK REQUIRES A NEWER ifort VERSION ====================
|
|
|
|
|
======================= THIS VERSION OF DAMASK REQUIRES A NEWER ifort VERSION =================
|
2018-10-15 08:33:53 +05:30
|
|
|
|
===================================================================================================
|
2019-04-03 11:55:42 +05:30
|
|
|
|
----- WRONG COMPILER VERSION ----- WRONG COMPILER VERSION ----- WRONG COMPILER VERSION -----
|
2018-10-15 08:33:53 +05:30
|
|
|
|
===================================================================================================
|
|
|
|
|
#endif
|
|
|
|
|
|
2019-04-03 18:19:16 +05:30
|
|
|
|
#if PETSC_VERSION_MAJOR!=3 || PETSC_VERSION_MINOR<PETSC_MINOR_MIN || PETSC_VERSION_MINOR>PETSC_MINOR_MAX
|
2018-08-20 21:27:15 +05:30
|
|
|
|
===================================================================================================
|
2019-04-03 15:23:59 +05:30
|
|
|
|
-- WRONG PETSc VERSION --- WRONG PETSc VERSION --- WRONG PETSc VERSION --- WRONG PETSc VERSION --
|
2018-08-20 21:27:15 +05:30
|
|
|
|
===================================================================================================
|
2019-04-03 15:23:59 +05:30
|
|
|
|
============ THIS VERSION OF DAMASK REQUIRES A DIFFERENT PETSc VERSION ========================
|
|
|
|
|
=============== THIS VERSION OF DAMASK REQUIRES A DIFFERENT PETSc VERSION =====================
|
|
|
|
|
================== THIS VERSION OF DAMASK REQUIRES A DIFFERENT PETSc VERSION ==================
|
2018-08-20 21:27:15 +05:30
|
|
|
|
===================================================================================================
|
2019-04-03 15:23:59 +05:30
|
|
|
|
-- WRONG PETSc VERSION --- WRONG PETSc VERSION --- WRONG PETSc VERSION --- WRONG PETSc VERSION --
|
2018-05-17 15:34:21 +05:30
|
|
|
|
===================================================================================================
|
|
|
|
|
#endif
|
2018-10-15 08:33:53 +05:30
|
|
|
|
|
2019-03-06 20:17:48 +05:30
|
|
|
|
character(len=1024) :: &
|
|
|
|
|
commandLine, & !< command line call as string
|
2019-03-24 15:18:46 +05:30
|
|
|
|
arg, & !< individual argument
|
|
|
|
|
loadCaseArg = '', & !< -l argument given to the executable
|
2019-03-06 20:17:48 +05:30
|
|
|
|
geometryArg = '', & !< -g argument given to the executable
|
|
|
|
|
workingDirArg = '', & !< -w argument given to the executable
|
|
|
|
|
userName !< name of user calling the executable
|
|
|
|
|
integer :: &
|
2019-03-24 15:18:46 +05:30
|
|
|
|
stat, &
|
2019-03-06 20:17:48 +05:30
|
|
|
|
i, &
|
2017-04-25 16:04:14 +05:30
|
|
|
|
#ifdef _OPENMP
|
2019-03-06 20:17:48 +05:30
|
|
|
|
threadLevel, &
|
2017-04-25 16:04:14 +05:30
|
|
|
|
#endif
|
2019-03-06 20:17:48 +05:30
|
|
|
|
worldrank = 0, &
|
2019-03-09 15:19:56 +05:30
|
|
|
|
worldsize = 0, &
|
|
|
|
|
typeSize
|
2019-03-06 20:17:48 +05:30
|
|
|
|
integer, dimension(8) :: &
|
|
|
|
|
dateAndTime
|
2019-03-09 15:19:56 +05:30
|
|
|
|
integer :: mpi_err
|
|
|
|
|
PetscErrorCode :: petsc_err
|
2019-03-06 20:17:48 +05:30
|
|
|
|
external :: &
|
|
|
|
|
quit
|
|
|
|
|
|
|
|
|
|
open(6, encoding='UTF-8') ! for special characters in output
|
2016-06-29 17:10:54 +05:30
|
|
|
|
|
2012-11-06 21:30:51 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
! PETSc Init
|
2018-02-12 19:33:35 +05:30
|
|
|
|
#ifdef _OPENMP
|
2019-03-06 20:17:48 +05:30
|
|
|
|
! If openMP is enabled, check if the MPI libary supports it and initialize accordingly.
|
|
|
|
|
! Otherwise, the first call to PETSc will do the initialization.
|
2019-03-09 15:19:56 +05:30
|
|
|
|
call MPI_Init_Thread(MPI_THREAD_FUNNELED,threadLevel,mpi_err)
|
|
|
|
|
if (mpi_err /= 0) call quit(1)
|
2019-03-06 20:17:48 +05:30
|
|
|
|
if (threadLevel<MPI_THREAD_FUNNELED) then
|
2019-05-04 21:17:52 +05:30
|
|
|
|
write(6,'(/,a)') ' ERROR: MPI library does not support OpenMP'
|
2019-03-06 20:17:48 +05:30
|
|
|
|
call quit(1)
|
|
|
|
|
endif
|
2016-06-27 15:53:42 +05:30
|
|
|
|
#endif
|
2019-05-03 09:51:43 +05:30
|
|
|
|
call PETScInitializeNoArguments(petsc_err) ! according to PETSc manual, that should be the first line in the code
|
2019-03-09 15:19:56 +05:30
|
|
|
|
CHKERRQ(petsc_err) ! this is a macro definition, it is case sensitive
|
|
|
|
|
|
|
|
|
|
call MPI_Comm_rank(PETSC_COMM_WORLD,worldrank,mpi_err)
|
|
|
|
|
if (mpi_err /= 0) call quit(1)
|
|
|
|
|
call MPI_Comm_size(PETSC_COMM_WORLD,worldsize,mpi_err)
|
|
|
|
|
if (mpi_err /= 0) call quit(1)
|
|
|
|
|
|
2019-03-06 20:17:48 +05:30
|
|
|
|
mainProcess: if (worldrank == 0) then
|
|
|
|
|
if (output_unit /= 6) then
|
2019-05-04 21:17:52 +05:30
|
|
|
|
write(output_unit,'(/,a)') ' ERROR: STDOUT != 6'
|
2019-03-06 20:17:48 +05:30
|
|
|
|
call quit(1)
|
|
|
|
|
endif
|
|
|
|
|
if (error_unit /= 0) then
|
2019-05-04 21:17:52 +05:30
|
|
|
|
write(output_unit,'(/,a)') ' ERROR: STDERR != 0'
|
2019-03-06 20:17:48 +05:30
|
|
|
|
call quit(1)
|
|
|
|
|
endif
|
|
|
|
|
else mainProcess
|
|
|
|
|
close(6) ! disable output for non-master processes (open 6 to rank specific file for debug)
|
|
|
|
|
open(6,file='/dev/null',status='replace') ! close(6) alone will leave some temp files in cwd
|
|
|
|
|
endif mainProcess
|
|
|
|
|
|
|
|
|
|
write(6,'(/,a)') ' <<<+- DAMASK_interface init -+>>>'
|
|
|
|
|
|
2019-03-14 02:53:09 +05:30
|
|
|
|
! http://patorjk.com/software/taag/#p=display&f=Lean&t=DAMASK
|
2019-03-13 10:46:31 +05:30
|
|
|
|
write(6,*) achar(27)//'[94m'
|
|
|
|
|
write(6,*) ' _/_/_/ _/_/ _/ _/ _/_/ _/_/_/ _/ _/'
|
|
|
|
|
write(6,*) ' _/ _/ _/ _/ _/_/ _/_/ _/ _/ _/ _/ _/'
|
|
|
|
|
write(6,*) ' _/ _/ _/_/_/_/ _/ _/ _/ _/_/_/_/ _/_/ _/_/'
|
|
|
|
|
write(6,*) ' _/ _/ _/ _/ _/ _/ _/ _/ _/ _/ _/'
|
|
|
|
|
write(6,*) ' _/_/_/ _/ _/ _/ _/ _/ _/ _/_/_/ _/ _/'
|
|
|
|
|
write(6,*) achar(27)//'[0m'
|
|
|
|
|
|
2019-03-13 03:26:09 +05:30
|
|
|
|
write(6,'(/,a)') ' Roters et al., Computational Materials Science 158:420–478, 2019'
|
2019-03-09 15:19:56 +05:30
|
|
|
|
write(6,'(a)') ' https://doi.org/10.1016/j.commatsci.2018.04.030'
|
|
|
|
|
|
|
|
|
|
write(6,'(/,a)') ' Version: '//DAMASKVERSION
|
2019-03-06 20:17:48 +05:30
|
|
|
|
|
|
|
|
|
! https://github.com/jeffhammond/HPCInfo/blob/master/docs/Preprocessor-Macros.md
|
2019-02-16 03:24:38 +05:30
|
|
|
|
#if defined(__GFORTRAN__) || __INTEL_COMPILER >= 1800
|
2019-03-13 10:46:31 +05:30
|
|
|
|
write(6,'(/,a)') ' Compiled with: '//compiler_version()
|
|
|
|
|
write(6,'(a)') ' Compiler options: '//compiler_options()
|
2019-02-16 03:24:38 +05:30
|
|
|
|
#elif defined(__INTEL_COMPILER)
|
2019-03-09 15:19:56 +05:30
|
|
|
|
write(6,'(/,a,i4.4,a,i8.8)') ' Compiled with Intel fortran version :', __INTEL_COMPILER,&
|
|
|
|
|
', build date :', __INTEL_COMPILER_BUILD_DATE
|
2019-02-16 03:24:38 +05:30
|
|
|
|
#elif defined(__PGI)
|
2019-03-09 15:19:56 +05:30
|
|
|
|
write(6,'(a,i4.4,a,i8.8)') ' Compiled with PGI fortran version :', __PGIC__,&
|
|
|
|
|
'.', __PGIC_MINOR__
|
2019-02-16 03:24:38 +05:30
|
|
|
|
#endif
|
|
|
|
|
|
2019-03-09 15:19:56 +05:30
|
|
|
|
write(6,'(/,a)') ' Compiled on: '//__DATE__//' at '//__TIME__
|
|
|
|
|
|
|
|
|
|
call date_and_time(values = dateAndTime)
|
|
|
|
|
write(6,'(/,a,2(i2.2,a),i4.4)') ' Date: ',dateAndTime(3),'/',dateAndTime(2),'/', dateAndTime(1)
|
|
|
|
|
write(6,'(a,2(i2.2,a),i2.2)') ' Time: ',dateAndTime(5),':', dateAndTime(6),':', dateAndTime(7)
|
2019-03-06 20:17:48 +05:30
|
|
|
|
|
2019-03-09 15:19:56 +05:30
|
|
|
|
call MPI_Type_size(MPI_INTEGER,typeSize,mpi_err)
|
|
|
|
|
if (mpi_err /= 0) call quit(1)
|
|
|
|
|
if (typeSize*8 /= bit_size(0)) then
|
|
|
|
|
write(6,'(a)') ' Mismatch between MPI and DAMASK integer'
|
|
|
|
|
call quit(1)
|
|
|
|
|
endif
|
|
|
|
|
|
|
|
|
|
call MPI_Type_size(MPI_DOUBLE,typeSize,mpi_err)
|
|
|
|
|
if (mpi_err /= 0) call quit(1)
|
|
|
|
|
if (typeSize*8 /= storage_size(0.0_pReal)) then
|
|
|
|
|
write(6,'(a)') ' Mismatch between MPI and DAMASK real'
|
|
|
|
|
call quit(1)
|
|
|
|
|
endif
|
2019-03-06 20:17:48 +05:30
|
|
|
|
|
2019-03-24 15:18:46 +05:30
|
|
|
|
do i = 1, command_argument_count()
|
|
|
|
|
call get_command_argument(i,arg)
|
|
|
|
|
select case(trim(arg)) ! extract key
|
2019-03-06 20:17:48 +05:30
|
|
|
|
case ('-h','--help')
|
|
|
|
|
write(6,'(a)') ' #######################################################################'
|
|
|
|
|
write(6,'(a)') ' DAMASK Command Line Interface:'
|
|
|
|
|
write(6,'(a)') ' For PETSc-based solvers for the Düsseldorf Advanced Material Simulation Kit'
|
|
|
|
|
write(6,'(a,/)')' #######################################################################'
|
|
|
|
|
write(6,'(a,/)')' Valid command line switches:'
|
|
|
|
|
write(6,'(a)') ' --geom (-g, --geometry)'
|
|
|
|
|
write(6,'(a)') ' --load (-l, --loadcase)'
|
|
|
|
|
write(6,'(a)') ' --workingdir (-w, --wd, --workingdirectory, -d, --directory)'
|
|
|
|
|
write(6,'(a)') ' --restart (-r, --rs)'
|
|
|
|
|
write(6,'(a)') ' --help (-h)'
|
|
|
|
|
write(6,'(/,a)')' -----------------------------------------------------------------------'
|
|
|
|
|
write(6,'(a)') ' Mandatory arguments:'
|
|
|
|
|
write(6,'(/,a)')' --geom PathToGeomFile/NameOfGeom'
|
|
|
|
|
write(6,'(a)') ' Specifies the location of the geometry definition file.'
|
|
|
|
|
write(6,'(/,a)')' --load PathToLoadFile/NameOfLoadFile'
|
|
|
|
|
write(6,'(a)') ' Specifies the location of the load case definition file.'
|
|
|
|
|
write(6,'(/,a)')' -----------------------------------------------------------------------'
|
|
|
|
|
write(6,'(a)') ' Optional arguments:'
|
|
|
|
|
write(6,'(/,a)')' --workingdirectory PathToWorkingDirectory'
|
|
|
|
|
write(6,'(a)') ' Specifies the working directory and overwrites the default ./'
|
|
|
|
|
write(6,'(a)') ' Make sure the file "material.config" exists in the working'
|
|
|
|
|
write(6,'(a)') ' directory.'
|
|
|
|
|
write(6,'(a)') ' For further configuration place "numerics.config"'
|
|
|
|
|
write(6,'(a)')' and "debug.config" in that directory.'
|
|
|
|
|
write(6,'(/,a)')' --restart XX'
|
|
|
|
|
write(6,'(a)') ' Reads in increment XX and continues with calculating'
|
|
|
|
|
write(6,'(a)') ' increment XX+1 based on this.'
|
|
|
|
|
write(6,'(a)') ' Appends to existing results file'
|
|
|
|
|
write(6,'(a)') ' "NameOfGeom_NameOfLoadFile".'
|
|
|
|
|
write(6,'(a)') ' Works only if the restart information for increment XX'
|
|
|
|
|
write(6,'(a)') ' is available in the working directory.'
|
|
|
|
|
write(6,'(/,a)')' -----------------------------------------------------------------------'
|
|
|
|
|
write(6,'(a)') ' Help:'
|
|
|
|
|
write(6,'(/,a)')' --help'
|
|
|
|
|
write(6,'(a,/)')' Prints this message and exits'
|
|
|
|
|
call quit(0) ! normal Termination
|
|
|
|
|
case ('-l', '--load', '--loadcase')
|
2019-03-24 15:18:46 +05:30
|
|
|
|
call get_command_argument(i+1,loadCaseArg)
|
2019-03-06 20:17:48 +05:30
|
|
|
|
case ('-g', '--geom', '--geometry')
|
2019-03-24 15:18:46 +05:30
|
|
|
|
call get_command_argument(i+1,geometryArg)
|
2019-03-06 20:17:48 +05:30
|
|
|
|
case ('-w', '-d', '--wd', '--directory', '--workingdir', '--workingdirectory')
|
2019-03-24 15:18:46 +05:30
|
|
|
|
call get_command_argument(i+1,workingDirArg)
|
2019-03-06 20:17:48 +05:30
|
|
|
|
case ('-r', '--rs', '--restart')
|
2019-03-24 15:18:46 +05:30
|
|
|
|
call get_command_argument(i+1,arg)
|
|
|
|
|
read(arg,*,iostat=stat) interface_restartInc
|
|
|
|
|
if (interface_restartInc < 0 .or. stat /=0) then
|
2019-05-04 21:17:52 +05:30
|
|
|
|
write(6,'(/,a)') ' ERROR: Could not parse restart increment: '//trim(arg)
|
2019-03-24 15:18:46 +05:30
|
|
|
|
call quit(1)
|
2019-03-06 20:17:48 +05:30
|
|
|
|
endif
|
|
|
|
|
end select
|
|
|
|
|
enddo
|
|
|
|
|
|
|
|
|
|
if (len_trim(loadcaseArg) == 0 .or. len_trim(geometryArg) == 0) then
|
2019-05-04 21:17:52 +05:30
|
|
|
|
write(6,'(/,a)') ' ERROR: Please specify geometry AND load case (-h for help)'
|
2019-03-06 20:17:48 +05:30
|
|
|
|
call quit(1)
|
|
|
|
|
endif
|
|
|
|
|
|
|
|
|
|
if (len_trim(workingDirArg) > 0) call setWorkingDirectory(trim(workingDirArg))
|
|
|
|
|
geometryFile = getGeometryFile(geometryArg)
|
|
|
|
|
loadCaseFile = getLoadCaseFile(loadCaseArg)
|
|
|
|
|
|
2019-03-24 15:18:46 +05:30
|
|
|
|
call get_command(commandLine)
|
2019-03-06 20:17:48 +05:30
|
|
|
|
call get_environment_variable('USER',userName)
|
|
|
|
|
! ToDo: https://stackoverflow.com/questions/8953424/how-to-get-the-username-in-c-c-in-linux
|
|
|
|
|
write(6,'(/,a,i4.1)') ' MPI processes: ',worldsize
|
|
|
|
|
write(6,'(a,a)') ' Host name: ', trim(getHostName())
|
|
|
|
|
write(6,'(a,a)') ' User name: ', trim(userName)
|
|
|
|
|
|
|
|
|
|
write(6,'(/a,a)') ' Command line call: ', trim(commandLine)
|
|
|
|
|
if (len(trim(workingDirArg)) > 0) &
|
|
|
|
|
write(6,'(a,a)') ' Working dir argument: ', trim(workingDirArg)
|
|
|
|
|
write(6,'(a,a)') ' Geometry argument: ', trim(geometryArg)
|
|
|
|
|
write(6,'(a,a)') ' Load case argument: ', trim(loadcaseArg)
|
|
|
|
|
write(6,'(a,a)') ' Working directory: ', trim(getCWD())
|
|
|
|
|
write(6,'(a,a)') ' Geometry file: ', trim(geometryFile)
|
|
|
|
|
write(6,'(a,a)') ' Loadcase file: ', trim(loadCaseFile)
|
|
|
|
|
write(6,'(a,a)') ' Solver job name: ', trim(getSolverJobName())
|
|
|
|
|
if (interface_restartInc > 0) &
|
|
|
|
|
write(6,'(a,i6.6)') ' Restart from increment: ', interface_restartInc
|
|
|
|
|
|
2019-03-24 16:29:00 +05:30
|
|
|
|
!call signalterm_c(c_funloc(catchSIGTERM))
|
|
|
|
|
call signalusr1_c(c_funloc(catchSIGUSR1))
|
|
|
|
|
call signalusr2_c(c_funloc(catchSIGUSR2))
|
|
|
|
|
call setSIGTERM(.false.)
|
|
|
|
|
call setSIGUSR1(.false.)
|
|
|
|
|
call setSIGUSR2(.false.)
|
2019-02-11 23:16:14 +05:30
|
|
|
|
|
|
|
|
|
|
2012-03-06 20:22:48 +05:30
|
|
|
|
end subroutine DAMASK_interface_init
|
2011-11-04 01:02:11 +05:30
|
|
|
|
|
2013-01-02 22:32:12 +05:30
|
|
|
|
|
2018-06-18 19:49:03 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
!> @brief extract working directory from given argument or from location of geometry file,
|
|
|
|
|
!! possibly converting relative arguments to absolut path
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2018-08-21 02:44:34 +05:30
|
|
|
|
subroutine setWorkingDirectory(workingDirectoryArg)
|
2019-03-06 20:17:48 +05:30
|
|
|
|
|
|
|
|
|
character(len=*), intent(in) :: workingDirectoryArg !< working directory argument
|
|
|
|
|
character(len=1024) :: workingDirectory !< working directory argument
|
|
|
|
|
logical :: error
|
|
|
|
|
external :: quit
|
|
|
|
|
|
|
|
|
|
absolutePath: if (workingDirectoryArg(1:1) == '/') then
|
|
|
|
|
workingDirectory = workingDirectoryArg
|
|
|
|
|
else absolutePath
|
|
|
|
|
workingDirectory = getCWD()
|
|
|
|
|
workingDirectory = trim(workingDirectory)//'/'//workingDirectoryArg
|
|
|
|
|
endif absolutePath
|
|
|
|
|
|
|
|
|
|
workingDirectory = trim(rectifyPath(workingDirectory))
|
|
|
|
|
error = setCWD(trim(workingDirectory))
|
|
|
|
|
if(error) then
|
2019-05-04 21:17:52 +05:30
|
|
|
|
write(6,'(/,a)') ' ERROR: Working directory "'//trim(workingDirectory)//'" does not exist'
|
2019-03-06 20:17:48 +05:30
|
|
|
|
call quit(1)
|
|
|
|
|
endif
|
2013-01-02 22:32:12 +05:30
|
|
|
|
|
2018-08-21 02:44:34 +05:30
|
|
|
|
end subroutine setWorkingDirectory
|
added fftw3 as fft(library will not versioned, should be in a linkable folder) , did some corrections on the code, splitted main file up (allows use of makefile), added makefile
changes on mpie_spectral.f90:
new structure, changed variable names, now using defgrad instead of disgrad, cleaned up, removed augmented Lagrange.
ToDo: Implement Augmented Lagrange again (but then a working version), implement Large strain, think about complex-to real-transform backwards, try to implement MP-support
2010-08-27 22:09:38 +05:30
|
|
|
|
|
2011-11-04 01:02:11 +05:30
|
|
|
|
|
2012-03-06 20:22:48 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2012-06-15 21:40:21 +05:30
|
|
|
|
!> @brief solver job name (no extension) as combination of geometry and load case name
|
2012-03-06 20:22:48 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
character(len=1024) function getSolverJobName()
|
added fftw3 as fft(library will not versioned, should be in a linkable folder) , did some corrections on the code, splitted main file up (allows use of makefile), added makefile
changes on mpie_spectral.f90:
new structure, changed variable names, now using defgrad instead of disgrad, cleaned up, removed augmented Lagrange.
ToDo: Implement Augmented Lagrange again (but then a working version), implement Large strain, think about complex-to real-transform backwards, try to implement MP-support
2010-08-27 22:09:38 +05:30
|
|
|
|
|
2019-03-06 20:17:48 +05:30
|
|
|
|
integer :: posExt,posSep
|
2012-02-21 21:34:16 +05:30
|
|
|
|
|
2019-05-09 12:03:12 +05:30
|
|
|
|
posExt = scan(geometryFile,'.',back=.true.)
|
|
|
|
|
posSep = scan(geometryFile,'/',back=.true.)
|
added fftw3 as fft(library will not versioned, should be in a linkable folder) , did some corrections on the code, splitted main file up (allows use of makefile), added makefile
changes on mpie_spectral.f90:
new structure, changed variable names, now using defgrad instead of disgrad, cleaned up, removed augmented Lagrange.
ToDo: Implement Augmented Lagrange again (but then a working version), implement Large strain, think about complex-to real-transform backwards, try to implement MP-support
2010-08-27 22:09:38 +05:30
|
|
|
|
|
2019-05-09 12:03:12 +05:30
|
|
|
|
getSolverJobName = geometryFile(posSep+1:posExt-1)
|
added fftw3 as fft(library will not versioned, should be in a linkable folder) , did some corrections on the code, splitted main file up (allows use of makefile), added makefile
changes on mpie_spectral.f90:
new structure, changed variable names, now using defgrad instead of disgrad, cleaned up, removed augmented Lagrange.
ToDo: Implement Augmented Lagrange again (but then a working version), implement Large strain, think about complex-to real-transform backwards, try to implement MP-support
2010-08-27 22:09:38 +05:30
|
|
|
|
|
2019-05-09 12:03:12 +05:30
|
|
|
|
posExt = scan(loadCaseFile,'.',back=.true.)
|
|
|
|
|
posSep = scan(loadCaseFile,'/',back=.true.)
|
2012-06-15 21:40:21 +05:30
|
|
|
|
|
2019-05-09 12:03:12 +05:30
|
|
|
|
getSolverJobName = trim(getSolverJobName)//'_'//loadCaseFile(posSep+1:posExt-1)
|
added fftw3 as fft(library will not versioned, should be in a linkable folder) , did some corrections on the code, splitted main file up (allows use of makefile), added makefile
changes on mpie_spectral.f90:
new structure, changed variable names, now using defgrad instead of disgrad, cleaned up, removed augmented Lagrange.
ToDo: Implement Augmented Lagrange again (but then a working version), implement Large strain, think about complex-to real-transform backwards, try to implement MP-support
2010-08-27 22:09:38 +05:30
|
|
|
|
|
2012-06-15 21:40:21 +05:30
|
|
|
|
end function getSolverJobName
|
added fftw3 as fft(library will not versioned, should be in a linkable folder) , did some corrections on the code, splitted main file up (allows use of makefile), added makefile
changes on mpie_spectral.f90:
new structure, changed variable names, now using defgrad instead of disgrad, cleaned up, removed augmented Lagrange.
ToDo: Implement Augmented Lagrange again (but then a working version), implement Large strain, think about complex-to real-transform backwards, try to implement MP-support
2010-08-27 22:09:38 +05:30
|
|
|
|
|
2011-02-07 20:05:42 +05:30
|
|
|
|
|
2012-03-06 20:22:48 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2012-06-15 21:40:21 +05:30
|
|
|
|
!> @brief basename of geometry file with extension from command line arguments
|
2012-03-06 20:22:48 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2012-06-15 21:40:21 +05:30
|
|
|
|
character(len=1024) function getGeometryFile(geometryParameter)
|
2011-02-07 20:05:42 +05:30
|
|
|
|
|
2019-03-06 20:17:48 +05:30
|
|
|
|
character(len=1024), intent(in) :: geometryParameter
|
|
|
|
|
logical :: file_exists
|
|
|
|
|
external :: quit
|
2011-02-07 20:05:42 +05:30
|
|
|
|
|
2019-03-06 20:17:48 +05:30
|
|
|
|
getGeometryFile = trim(geometryParameter)
|
|
|
|
|
if (scan(getGeometryFile,'/') /= 1) getGeometryFile = trim(getCWD())//'/'//trim(getGeometryFile)
|
|
|
|
|
getGeometryFile = makeRelativePath(trim(getCWD()), getGeometryFile)
|
2011-02-07 20:05:42 +05:30
|
|
|
|
|
2019-03-06 20:17:48 +05:30
|
|
|
|
inquire(file=trim(getGeometryFile), exist=file_exists)
|
|
|
|
|
if (.not. file_exists) then
|
2019-05-04 21:17:52 +05:30
|
|
|
|
write(6,'(/,a)') ' ERROR: Geometry file does not exists ('//trim(getGeometryFile)//')'
|
2019-03-06 20:17:48 +05:30
|
|
|
|
call quit(1)
|
|
|
|
|
endif
|
2011-02-07 20:05:42 +05:30
|
|
|
|
|
2012-06-15 21:40:21 +05:30
|
|
|
|
end function getGeometryFile
|
2011-02-07 20:05:42 +05:30
|
|
|
|
|
added fftw3 as fft(library will not versioned, should be in a linkable folder) , did some corrections on the code, splitted main file up (allows use of makefile), added makefile
changes on mpie_spectral.f90:
new structure, changed variable names, now using defgrad instead of disgrad, cleaned up, removed augmented Lagrange.
ToDo: Implement Augmented Lagrange again (but then a working version), implement Large strain, think about complex-to real-transform backwards, try to implement MP-support
2010-08-27 22:09:38 +05:30
|
|
|
|
|
2012-03-06 20:22:48 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
!> @brief relative path of loadcase from command line arguments
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2012-06-15 21:40:21 +05:30
|
|
|
|
character(len=1024) function getLoadCaseFile(loadCaseParameter)
|
added fftw3 as fft(library will not versioned, should be in a linkable folder) , did some corrections on the code, splitted main file up (allows use of makefile), added makefile
changes on mpie_spectral.f90:
new structure, changed variable names, now using defgrad instead of disgrad, cleaned up, removed augmented Lagrange.
ToDo: Implement Augmented Lagrange again (but then a working version), implement Large strain, think about complex-to real-transform backwards, try to implement MP-support
2010-08-27 22:09:38 +05:30
|
|
|
|
|
2019-03-06 20:17:48 +05:30
|
|
|
|
character(len=1024), intent(in) :: loadCaseParameter
|
|
|
|
|
logical :: file_exists
|
|
|
|
|
external :: quit
|
2013-12-28 01:33:28 +05:30
|
|
|
|
|
2019-03-06 20:17:48 +05:30
|
|
|
|
getLoadCaseFile = trim(loadCaseParameter)
|
|
|
|
|
if (scan(getLoadCaseFile,'/') /= 1) getLoadCaseFile = trim(getCWD())//'/'//trim(getLoadCaseFile)
|
|
|
|
|
getLoadCaseFile = makeRelativePath(trim(getCWD()), getLoadCaseFile)
|
2012-06-15 21:40:21 +05:30
|
|
|
|
|
2019-03-06 20:17:48 +05:30
|
|
|
|
inquire(file=trim(getLoadCaseFile), exist=file_exists)
|
|
|
|
|
if (.not. file_exists) then
|
2019-05-04 21:17:52 +05:30
|
|
|
|
write(6,'(/,a)') ' ERROR: Load case file does not exists ('//trim(getLoadCaseFile)//')'
|
2019-03-06 20:17:48 +05:30
|
|
|
|
call quit(1)
|
|
|
|
|
endif
|
2018-09-15 18:24:56 +05:30
|
|
|
|
|
2012-06-15 21:40:21 +05:30
|
|
|
|
end function getLoadCaseFile
|
added fftw3 as fft(library will not versioned, should be in a linkable folder) , did some corrections on the code, splitted main file up (allows use of makefile), added makefile
changes on mpie_spectral.f90:
new structure, changed variable names, now using defgrad instead of disgrad, cleaned up, removed augmented Lagrange.
ToDo: Implement Augmented Lagrange again (but then a working version), implement Large strain, think about complex-to real-transform backwards, try to implement MP-support
2010-08-27 22:09:38 +05:30
|
|
|
|
|
|
|
|
|
|
2012-03-06 20:22:48 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2018-07-10 13:53:21 +05:30
|
|
|
|
!> @brief remove ../, /./, and // from path.
|
2018-06-29 19:06:12 +05:30
|
|
|
|
!> @details works only if absolute path is given
|
2012-03-06 20:22:48 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
added fftw3 as fft(library will not versioned, should be in a linkable folder) , did some corrections on the code, splitted main file up (allows use of makefile), added makefile
changes on mpie_spectral.f90:
new structure, changed variable names, now using defgrad instead of disgrad, cleaned up, removed augmented Lagrange.
ToDo: Implement Augmented Lagrange again (but then a working version), implement Large strain, think about complex-to real-transform backwards, try to implement MP-support
2010-08-27 22:09:38 +05:30
|
|
|
|
function rectifyPath(path)
|
|
|
|
|
|
2019-03-09 15:19:56 +05:30
|
|
|
|
character(len=*) :: path
|
2019-03-06 20:17:48 +05:30
|
|
|
|
character(len=1024) :: rectifyPath
|
|
|
|
|
integer :: i,j,k,l
|
added fftw3 as fft(library will not versioned, should be in a linkable folder) , did some corrections on the code, splitted main file up (allows use of makefile), added makefile
changes on mpie_spectral.f90:
new structure, changed variable names, now using defgrad instead of disgrad, cleaned up, removed augmented Lagrange.
ToDo: Implement Augmented Lagrange again (but then a working version), implement Large strain, think about complex-to real-transform backwards, try to implement MP-support
2010-08-27 22:09:38 +05:30
|
|
|
|
|
2013-05-13 19:40:48 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
! remove /./ from path
|
2019-03-06 20:17:48 +05:30
|
|
|
|
rectifyPath = trim(path)
|
|
|
|
|
l = len_trim(rectifyPath)
|
|
|
|
|
do i = l,3,-1
|
|
|
|
|
if (rectifyPath(i-2:i) == '/./') rectifyPath(i-1:l) = rectifyPath(i+1:l)//' '
|
|
|
|
|
enddo
|
2018-07-10 13:53:21 +05:30
|
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
! remove // from path
|
2019-03-06 20:17:48 +05:30
|
|
|
|
l = len_trim(rectifyPath)
|
|
|
|
|
do i = l,2,-1
|
|
|
|
|
if (rectifyPath(i-1:i) == '//') rectifyPath(i-1:l) = rectifyPath(i:l)//' '
|
|
|
|
|
enddo
|
added fftw3 as fft(library will not versioned, should be in a linkable folder) , did some corrections on the code, splitted main file up (allows use of makefile), added makefile
changes on mpie_spectral.f90:
new structure, changed variable names, now using defgrad instead of disgrad, cleaned up, removed augmented Lagrange.
ToDo: Implement Augmented Lagrange again (but then a working version), implement Large strain, think about complex-to real-transform backwards, try to implement MP-support
2010-08-27 22:09:38 +05:30
|
|
|
|
|
2013-05-13 19:40:48 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
! remove ../ and corresponding directory from rectifyPath
|
2019-03-06 20:17:48 +05:30
|
|
|
|
l = len_trim(rectifyPath)
|
|
|
|
|
i = index(rectifyPath(i:l),'../')
|
|
|
|
|
j = 0
|
|
|
|
|
do while (i > j)
|
|
|
|
|
j = scan(rectifyPath(1:i-2),'/',back=.true.)
|
|
|
|
|
rectifyPath(j+1:l) = rectifyPath(i+3:l)//repeat(' ',2+i-j)
|
|
|
|
|
if (rectifyPath(j+1:j+1) == '/') then !search for '//' that appear in case of XXX/../../XXX
|
|
|
|
|
k = len_trim(rectifyPath)
|
|
|
|
|
rectifyPath(j+1:k-1) = rectifyPath(j+2:k)
|
|
|
|
|
rectifyPath(k:k) = ' '
|
|
|
|
|
endif
|
|
|
|
|
i = j+index(rectifyPath(j+1:l),'../')
|
|
|
|
|
enddo
|
|
|
|
|
if(len_trim(rectifyPath) == 0) rectifyPath = '/'
|
2011-08-01 23:40:55 +05:30
|
|
|
|
|
2012-03-06 20:22:48 +05:30
|
|
|
|
end function rectifyPath
|
added fftw3 as fft(library will not versioned, should be in a linkable folder) , did some corrections on the code, splitted main file up (allows use of makefile), added makefile
changes on mpie_spectral.f90:
new structure, changed variable names, now using defgrad instead of disgrad, cleaned up, removed augmented Lagrange.
ToDo: Implement Augmented Lagrange again (but then a working version), implement Large strain, think about complex-to real-transform backwards, try to implement MP-support
2010-08-27 22:09:38 +05:30
|
|
|
|
|
2019-03-06 20:17:48 +05:30
|
|
|
|
|
2012-03-06 20:22:48 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
!> @brief relative path from absolute a to absolute b
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
character(len=1024) function makeRelativePath(a,b)
|
added fftw3 as fft(library will not versioned, should be in a linkable folder) , did some corrections on the code, splitted main file up (allows use of makefile), added makefile
changes on mpie_spectral.f90:
new structure, changed variable names, now using defgrad instead of disgrad, cleaned up, removed augmented Lagrange.
ToDo: Implement Augmented Lagrange again (but then a working version), implement Large strain, think about complex-to real-transform backwards, try to implement MP-support
2010-08-27 22:09:38 +05:30
|
|
|
|
|
2019-03-06 20:17:48 +05:30
|
|
|
|
character (len=*), intent(in) :: a,b
|
|
|
|
|
character (len=1024) :: a_cleaned,b_cleaned
|
|
|
|
|
integer :: i,posLastCommonSlash,remainingSlashes
|
added fftw3 as fft(library will not versioned, should be in a linkable folder) , did some corrections on the code, splitted main file up (allows use of makefile), added makefile
changes on mpie_spectral.f90:
new structure, changed variable names, now using defgrad instead of disgrad, cleaned up, removed augmented Lagrange.
ToDo: Implement Augmented Lagrange again (but then a working version), implement Large strain, think about complex-to real-transform backwards, try to implement MP-support
2010-08-27 22:09:38 +05:30
|
|
|
|
|
2019-03-06 20:17:48 +05:30
|
|
|
|
posLastCommonSlash = 0
|
|
|
|
|
remainingSlashes = 0
|
|
|
|
|
a_cleaned = rectifyPath(trim(a)//'/')
|
|
|
|
|
b_cleaned = rectifyPath(b)
|
2012-02-21 21:34:16 +05:30
|
|
|
|
|
2019-03-06 20:17:48 +05:30
|
|
|
|
do i = 1, min(1024,len_trim(a_cleaned),len_trim(rectifyPath(b_cleaned)))
|
|
|
|
|
if (a_cleaned(i:i) /= b_cleaned(i:i)) exit
|
|
|
|
|
if (a_cleaned(i:i) == '/') posLastCommonSlash = i
|
|
|
|
|
enddo
|
|
|
|
|
do i = posLastCommonSlash+1,len_trim(a_cleaned)
|
|
|
|
|
if (a_cleaned(i:i) == '/') remainingSlashes = remainingSlashes + 1
|
|
|
|
|
enddo
|
2018-07-10 13:53:21 +05:30
|
|
|
|
|
2019-03-06 20:17:48 +05:30
|
|
|
|
makeRelativePath = repeat('..'//'/',remainingSlashes)//b_cleaned(posLastCommonSlash+1:len_trim(b_cleaned))
|
2011-08-01 23:40:55 +05:30
|
|
|
|
|
2012-03-06 20:22:48 +05:30
|
|
|
|
end function makeRelativePath
|
added fftw3 as fft(library will not versioned, should be in a linkable folder) , did some corrections on the code, splitted main file up (allows use of makefile), added makefile
changes on mpie_spectral.f90:
new structure, changed variable names, now using defgrad instead of disgrad, cleaned up, removed augmented Lagrange.
ToDo: Implement Augmented Lagrange again (but then a working version), implement Large strain, think about complex-to real-transform backwards, try to implement MP-support
2010-08-27 22:09:38 +05:30
|
|
|
|
|
2019-03-06 20:17:48 +05:30
|
|
|
|
|
2019-02-11 23:16:14 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-03-24 16:29:00 +05:30
|
|
|
|
!> @brief sets global variable SIGTERM to .true.
|
2019-02-11 23:16:14 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-03-24 16:29:00 +05:30
|
|
|
|
subroutine catchSIGTERM(signal) bind(C)
|
|
|
|
|
|
|
|
|
|
integer(C_INT), value :: signal
|
|
|
|
|
SIGTERM = .true.
|
|
|
|
|
|
|
|
|
|
write(6,'(a,i2.2,a)') ' received signal ',signal, ', set SIGTERM'
|
|
|
|
|
|
|
|
|
|
end subroutine catchSIGTERM
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
!> @brief sets global variable SIGTERM
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
subroutine setSIGTERM(state)
|
|
|
|
|
|
|
|
|
|
logical, intent(in) :: state
|
|
|
|
|
SIGTERM = state
|
|
|
|
|
|
|
|
|
|
end subroutine setSIGTERM
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
!> @brief sets global variable SIGUSR1 to .true.
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
subroutine catchSIGUSR1(signal) bind(C)
|
2019-03-06 20:17:48 +05:30
|
|
|
|
|
|
|
|
|
integer(C_INT), value :: signal
|
|
|
|
|
SIGUSR1 = .true.
|
|
|
|
|
|
2019-03-09 05:37:26 +05:30
|
|
|
|
write(6,'(a,i2.2,a)') ' received signal ',signal, ', set SIGUSR1'
|
2019-03-06 20:17:48 +05:30
|
|
|
|
|
2019-03-24 16:29:00 +05:30
|
|
|
|
end subroutine catchSIGUSR1
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
!> @brief sets global variable SIGUSR1
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
subroutine setSIGUSR1(state)
|
|
|
|
|
|
|
|
|
|
logical, intent(in) :: state
|
|
|
|
|
SIGUSR1 = state
|
|
|
|
|
|
2019-02-11 23:16:14 +05:30
|
|
|
|
end subroutine setSIGUSR1
|
2019-03-06 20:17:48 +05:30
|
|
|
|
|
|
|
|
|
|
2019-02-11 23:16:14 +05:30
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
!> @brief sets global variable SIGUSR2 to .true. if program receives SIGUSR2
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2019-03-24 16:29:00 +05:30
|
|
|
|
subroutine catchSIGUSR2(signal) bind(C)
|
2019-02-11 23:16:14 +05:30
|
|
|
|
|
2019-03-06 20:17:48 +05:30
|
|
|
|
integer(C_INT), value :: signal
|
|
|
|
|
SIGUSR2 = .true.
|
|
|
|
|
|
2019-03-09 05:37:26 +05:30
|
|
|
|
write(6,'(a,i2.2,a)') ' received signal ',signal, ', set SIGUSR2'
|
2019-02-11 23:16:14 +05:30
|
|
|
|
|
2019-03-24 16:29:00 +05:30
|
|
|
|
end subroutine catchSIGUSR2
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
!> @brief sets global variable SIGUSR2
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
|
subroutine setSIGUSR2(state)
|
|
|
|
|
|
|
|
|
|
logical, intent(in) :: state
|
|
|
|
|
SIGUSR2 = state
|
|
|
|
|
|
2019-02-11 23:16:14 +05:30
|
|
|
|
end subroutine setSIGUSR2
|
|
|
|
|
|
2019-03-24 16:29:00 +05:30
|
|
|
|
|
2019-02-16 03:24:38 +05:30
|
|
|
|
end module
|