make results self-contained for reproducibility

ToDo: same functionality for load and geom
This commit is contained in:
Martin Diehl 2021-07-27 07:28:35 +02:00
parent f20302e0ee
commit 3b06498c2f
2 changed files with 36 additions and 11 deletions

View File

@ -8,7 +8,8 @@ module config
use IO
use YAML_parse
use YAML_types
use results
use parallelization
implicit none
private
@ -31,6 +32,7 @@ subroutine config_init
print'(/,a)', ' <<<+- config init -+>>>'; flush(IO_STDOUT)
call parse_material
call parse_numerics
call parse_debug
@ -41,14 +43,21 @@ end subroutine config_init
!--------------------------------------------------------------------------------------------------
!> @brief Read material.yaml or <jobname>.yaml.
!--------------------------------------------------------------------------------------------------
subroutine parse_material
subroutine parse_material()
logical :: fileExists
character(len=:), allocatable :: fileContent
inquire(file='material.yaml',exist=fileExists)
if(.not. fileExists) call IO_error(100,ext_msg='material.yaml')
print*, 'reading material.yaml'; flush(IO_STDOUT)
fileContent = IO_read('material.yaml')
if (worldrank == 0) then
call results_openJobFile(parallel=.false.)
call results_writeDataset_str(fileContent,'setup','material.yaml','DAMASK main configuration')
call results_closeJobFile
endif
config_material => YAML_parse_file('material.yaml')
end subroutine parse_material
@ -57,15 +66,22 @@ end subroutine parse_material
!--------------------------------------------------------------------------------------------------
!> @brief Read numerics.yaml.
!--------------------------------------------------------------------------------------------------
subroutine parse_numerics
subroutine parse_numerics()
logical :: fexist
logical :: fileExists
character(len=:), allocatable :: fileContent
config_numerics => emptyDict
inquire(file='numerics.yaml', exist=fexist)
if (fexist) then
inquire(file='numerics.yaml', exist=fileExists)
if (fileExists) then
print*, 'reading numerics.yaml'; flush(IO_STDOUT)
fileContent = IO_read('numerics.yaml')
if (worldrank == 0) then
call results_openJobFile(parallel=.false.)
call results_writeDataset_str(fileContent,'setup','numerics.yaml','numerics configuration (optional)')
call results_closeJobFile
endif
config_numerics => YAML_parse_file('numerics.yaml')
endif
@ -75,17 +91,24 @@ end subroutine parse_numerics
!--------------------------------------------------------------------------------------------------
!> @brief Read debug.yaml.
!--------------------------------------------------------------------------------------------------
subroutine parse_debug
subroutine parse_debug()
logical :: fexist
logical :: fileExists
character(len=:), allocatable :: fileContent
config_debug => emptyDict
inquire(file='debug.yaml', exist=fexist)
fileExists: if (fexist) then
inquire(file='debug.yaml', exist=fileExists)
if (fileExists) then
print*, 'reading debug.yaml'; flush(IO_STDOUT)
fileContent = IO_read('debug.yaml')
if (worldrank == 0) then
call results_openJobFile(parallel=.false.)
call results_writeDataset_str(fileContent,'setup','debug.yaml','debug configuration (optional)')
call results_closeJobFile
endif
config_debug => YAML_parse_file('debug.yaml')
endif fileExists
endif
end subroutine parse_debug

View File

@ -82,6 +82,8 @@ subroutine results_init(restart)
call results_addAttribute('call',trim(commandLine))
call results_closeGroup(results_addGroup('cell_to'))
call results_addAttribute('description','mappings to place data in space','cell_to')
call results_closeGroup(results_addGroup('setup'))
call results_addAttribute('description','input data used to run the simulation','setup')
call results_closeJobFile
endif