simplified multi processor (MPI) reporting

This commit is contained in:
Martin Diehl 2014-10-10 13:08:34 +00:00
parent 1657e0f7ba
commit d095c2484d
4 changed files with 216 additions and 224 deletions

View File

@ -9,16 +9,14 @@
!> by DAMASK. Interpretating the command line arguments or, in case of called from f2py,
!> the arguments parsed to the init routine to get load case, geometry file, working
!> directory, etc.
!-----------------------------------------------------------`---------------------------------------
!--------------------------------------------------------------------------------------------------
module DAMASK_interface
use prec, only: &
pInt
implicit none
private
#ifdef PETSc
#include <finclude/petscsys.h>
#endif
logical, public, protected :: appendToOutFile = .false. !< Append to existing spectralOut file (in case of restart, not in case of regridding)
integer(pInt), public, protected :: spectralRestartInc = 1_pInt !< Increment at which calculation starts
@ -65,31 +63,33 @@ subroutine DAMASK_interface_init(loadCaseParameterIn,geometryParameterIn)
userName, & !< name of user calling DAMASK_spectral.exe
tag
integer :: &
i
i, &
worldrank = 0
integer, parameter :: &
maxNchunks = 128 !< DAMASK_spectral + (l,g,w,r)*2 + h
integer, dimension(1+ 2* maxNchunks) :: &
positions
integer, dimension(8) :: &
dateAndTime ! type default integer
PetscErrorCode :: ierr
external :: &
quit
#ifdef PETSc
external :: &
quit,&
MPI_Comm_rank,&
PETScInitialize, &
MPI_abort
PetscErrorCode :: ierr
!--------------------------------------------------------------------------------------------------
! PETSc Init
call PetscInitialize(PETSC_NULL_CHARACTER,ierr) ! according to PETSc manual, that should be the first line in the code
CHKERRQ(ierr) ! this is a macro definition, it is case sensitive
#endif
open(6, encoding='UTF-8') ! modern fortran compilers (gfortran >4.4, ifort >11 support it)
call MPI_Comm_rank(PETSC_COMM_WORLD,worldrank,ierr);CHKERRQ(ierr)
mainProcess: if (worldrank == 0) then
write(6,'(/,a)') ' <<<+- DAMASK_spectral_interface init -+>>>'
write(6,'(a)') ' $Id$'
#include "compilation_info.f90"
endif mainProcess
if ( present(loadcaseParameterIn) .and. present(geometryParameterIn)) then ! both mandatory parameters given in function call
geometryArg = geometryParameterIn
@ -102,6 +102,7 @@ subroutine DAMASK_interface_init(loadCaseParameterIn,geometryParameterIn)
tag = IIO_lc(IIO_stringValue(commandLine,positions,i)) ! extract key
select case(tag)
case ('-h','--help')
mainProcess2: if (worldrank == 0) then
write(6,'(a)') ' #######################################################################'
write(6,'(a)') ' DAMASK_spectral:'
write(6,'(a)') ' The spectral method boundary value problem solver for'
@ -156,6 +157,7 @@ subroutine DAMASK_interface_init(loadCaseParameterIn,geometryParameterIn)
write(6,'(/,a)')' --help'
write(6,'(a,/)')' Prints this message and exits'
call quit(0_pInt) ! normal Termination
endif mainProcess2
case ('-l', '--load', '--loadcase')
loadcaseArg = IIO_stringValue(commandLine,positions,i+1_pInt)
case ('-g', '--geom', '--geometry')
@ -184,7 +186,7 @@ subroutine DAMASK_interface_init(loadCaseParameterIn,geometryParameterIn)
call get_environment_variable('HOSTNAME',hostName)
call get_environment_variable('USER',userName)
call date_and_time(values = dateAndTime)
mainProcess3: if (worldrank == 0) then
write(6,'(a,2(i2.2,a),i4.4)') ' Date: ',dateAndTime(3),'/',&
dateAndTime(2),'/',&
dateAndTime(1)
@ -203,9 +205,10 @@ subroutine DAMASK_interface_init(loadCaseParameterIn,geometryParameterIn)
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 (SpectralRestartInc > 1_pInt) write(6,'(a,i6.6)') &
' Restart at increment: ', spectralRestartInc
if (SpectralRestartInc > 1_pInt) &
write(6,'(a,i6.6)') ' Restart at increment: ', spectralRestartInc
write(6,'(a,l1,/)') ' Append to result file: ', appendToOutFile
endif mainProcess3
end subroutine DAMASK_interface_init

View File

@ -102,23 +102,22 @@ contains
!--------------------------------------------------------------------------------------------------
subroutine IO_init
use, intrinsic :: iso_fortran_env ! to get compiler_version and compiler_options (at least for gfortran 4.6 at the moment)
#ifdef FEM
implicit none
integer(pInt) :: worldrank
#ifdef PETSc
#include <finclude/petscsys.h>
PetscInt :: worldrank
PetscErrorCode :: ierr
#endif
#ifdef FEM
#ifdef PETSc
call MPI_Comm_rank(PETSC_COMM_WORLD,worldrank,ierr);CHKERRQ(ierr)
if (worldrank == 0) then
#endif
mainProcess: if (worldrank == 0) then
write(6,'(/,a)') ' <<<+- IO init -+>>>'
write(6,'(a)') ' $Id$'
write(6,'(a15,a)') ' Current time: ',IO_timeStamp()
#include "compilation_info.f90"
#ifdef FEM
endif
#endif
#ifdef HDF
call HDF5_createJobFile

View File

@ -12,7 +12,7 @@ module numerics
implicit none
private
#ifdef FEM
#ifdef PETSc
#include <finclude/petsc.h90>
#endif
character(len=64), parameter, private :: &
@ -27,7 +27,9 @@ module numerics
nState = 10_pInt, & !< state loop limit
nStress = 40_pInt, & !< stress loop limit
pert_method = 1_pInt, & !< method used in perturbation technique for tangent
fixedSeed = 0_pInt !< fixed seeding for pseudo-random number generator, Default 0: use random seed
fixedSeed = 0_pInt, & !< fixed seeding for pseudo-random number generator, Default 0: use random seed
worldrank = 0_pInt, & !< MPI worldrank (/=0 for MPI simulations only)
worldsize = 0_pInt !< MPI worldsize (/=0 for MPI simulations only)
integer, protected, public :: &
DAMASK_NumThreadsInt = 0 !< value stored in environment variable DAMASK_NUM_THREADS, set to zero if no OpenMP directive
integer(pInt), public :: &
@ -154,9 +156,7 @@ module numerics
integrationOrder = 1_pInt, &
structOrder = 2_pInt, &
thermalOrder = 2_pInt, &
damageOrder = 2_pInt, &
worldrank = 0_pInt, &
worldsize = 0_pInt
damageOrder = 2_pInt
#endif
public :: numerics_init
@ -201,18 +201,16 @@ subroutine numerics_init
line
!$ character(len=6) DAMASK_NumThreadsString ! environment variable DAMASK_NUM_THREADS
#ifdef FEM
#ifdef PETSc
call MPI_Comm_rank(PETSC_COMM_WORLD,worldrank,ierr);CHKERRQ(ierr)
call MPI_Comm_size(PETSC_COMM_WORLD,worldsize,ierr);CHKERRQ(ierr)
if (worldrank == 0) then
#endif
mainProcess: if (worldrank == 0) then
write(6,'(/,a)') ' <<<+- numerics init -+>>>'
write(6,'(a)') ' $Id$'
write(6,'(a15,a)') ' Current time: ',IO_timeStamp()
#include "compilation_info.f90"
#ifdef FEM
endif
#endif
endif mainProcess
!$ call GET_ENVIRONMENT_VARIABLE(NAME='DAMASK_NUM_THREADS',VALUE=DAMASK_NumThreadsString,STATUS=gotDAMASK_NUM_THREADS) ! get environment variable DAMASK_NUM_THREADS...
!$ if(gotDAMASK_NUM_THREADS /= 0) then ! could not get number of threads, set it to 1
@ -227,14 +225,10 @@ subroutine numerics_init
!--------------------------------------------------------------------------------------------------
! try to open the config file
fileExists: if(IO_open_file_stat(FILEUNIT,numerics_configFile)) then
#ifdef FEM
if (worldrank == 0) then
#endif
mainProcess2: if (worldrank == 0) then
write(6,'(a,/)') ' using values from config file'
flush(6)
#ifdef FEM
endif
#endif
endif mainProcess2
!--------------------------------------------------------------------------------------------------
! read variables from config file and overwrite default parameters if keyword is present
@ -422,9 +416,9 @@ subroutine numerics_init
case ('petsc_optionsfem')
petsc_optionsFEM = trim(line(positions(4):))
#else
case ('err_struct_tolabs','err_struct_tolrel','err_thermal_tol','err_damage_tol','residualstiffness',& ! found spectral parameter for FEM build
'itmaxfem', 'itminfem','maxcutbackfem','integrationorder','structorder','thermalorder', &
'damageorder','petsc_optionsfem')
case ('err_struct_tolabs','err_struct_tolrel','err_thermal_tol','err_damage_tol',& ! found FEM parameter for spectral/Abaqus/Marc build
'residualstiffness', 'itmaxfem', 'itminfem','maxcutbackfem','integrationorder',&
'structorder','thermalorder', 'damageorder','petsc_optionsfem')
call IO_warning(40_pInt,ext_msg=tag)
#endif
case default ! found unknown keyword
@ -462,11 +456,9 @@ subroutine numerics_init
numerics_timeSyncing = numerics_timeSyncing .and. all(numerics_integrator==2_pInt) ! timeSyncing only allowed for explicit Euler integrator
#ifdef FEM
if (worldrank == 0) then
#endif
!--------------------------------------------------------------------------------------------------
! writing parameters to output
mainProcess3: if (worldrank == 0) then
write(6,'(a24,1x,es8.1)') ' relevantStrain: ',relevantStrain
write(6,'(a24,1x,es8.1)') ' defgradTolerance: ',defgradTolerance
write(6,'(a24,1x,i8)') ' iJacoStiffness: ',iJacoStiffness
@ -565,9 +557,8 @@ subroutine numerics_init
write(6,'(a24,1x,es8.1)') ' residualStiffness: ',residualStiffness
write(6,'(a24,1x,a)') ' PETSc_optionsFEM: ',trim(petsc_optionsFEM)
#endif
#ifdef FEM
endif
#endif
endif mainProcess3
!--------------------------------------------------------------------------------------------------
! sanity checks

View File

@ -96,18 +96,19 @@ subroutine prec_init
use, intrinsic :: iso_fortran_env ! to get compiler_version and compiler_options (at least for gfortran 4.6 at the moment)
implicit none
#ifdef FEM
integer(pInt) :: worldrank = 0_pInt
#ifdef PETSc
#include <finclude/petscsys.h>
PetscInt :: worldrank
PetscErrorCode :: ierr
#endif
external :: &
quit
#ifdef FEM
#ifdef PETSc
call MPI_Comm_rank(PETSC_COMM_WORLD,worldrank,ierr);CHKERRQ(ierr)
if (worldrank == 0) then
#endif
mainProcess: if (worldrank == 0) then
write(6,'(/,a)') ' <<<+- prec init -+>>>'
write(6,'(a)') ' $Id$'
#include "compilation_info.f90"
@ -116,9 +117,7 @@ subroutine prec_init
write(6,'(a,i3)') ' Bytes for pLongInt: ',pLongInt
write(6,'(a,e10.3)') ' NaN: ', DAMASK_NaN
write(6,'(a,l3,/)') ' NaN /= NaN: ',DAMASK_NaN/=DAMASK_NaN
#ifdef FEM
endif
#endif
endif mainProcess
if (DAMASK_NaN == DAMASK_NaN) call quit(9000)