DAMASK_EICMD/src/config.f90

822 lines
30 KiB
Fortran
Raw Normal View History

!-------------------------------------------------------------------------------------------------
2018-06-10 14:37:17 +05:30
!> @author Martin Diehl, Max-Planck-Institut für Eisenforschung GmbH
!> @brief Reads in the material configuration from file
!> @details Reads the material configuration file, where solverJobName.materialConfig takes
!! precedence over material.config. Stores the raw strings and the positions of delimiters for the
!! parts 'homogenization', 'crystallite', 'phase', 'texture', and 'microstucture'
!--------------------------------------------------------------------------------------------------
module config
use prec, only: &
2019-03-13 11:30:20 +05:30
pReal
implicit none
private
type, private :: tPartitionedString
character(len=:), allocatable :: val
2019-03-13 02:59:03 +05:30
integer, dimension(:), allocatable :: pos
end type tPartitionedString
2019-03-08 03:17:45 +05:30
type, private :: tPartitionedStringList
type(tPartitionedString) :: string
type(tPartitionedStringList), pointer :: next => null()
contains
procedure :: add => add
procedure :: show => show
procedure :: free => free
! currently, a finalize is needed for all shapes of tPartitionedStringList.
! with Fortran 2015, we can define one recursive elemental function
! https://software.intel.com/en-us/forums/intel-visual-fortran-compiler-for-windows/topic/543326
final :: finalize, &
finalizeArray
procedure :: keyExists => keyExists
procedure :: countKeys => countKeys
procedure :: getFloat => getFloat
procedure :: getInt => getInt
2018-06-22 03:18:34 +05:30
procedure :: getString => getString
procedure :: getFloats => getFloats
procedure :: getInts => getInts
2018-06-22 03:18:34 +05:30
procedure :: getStrings => getStrings
end type tPartitionedStringList
2018-08-22 18:00:51 +05:30
type(tPartitionedStringList), public, protected, allocatable, dimension(:) :: &
2018-06-27 00:24:54 +05:30
config_phase, &
config_microstructure, &
config_homogenization, &
config_texture, &
config_crystallite
type(tPartitionedStringList), public, protected :: &
config_numerics, &
config_debug
character(len=64), dimension(:), allocatable, public, protected :: &
phase_name, & !< name of each phase
homogenization_name, & !< name of each homogenization
crystallite_name, & !< name of each crystallite setting
microstructure_name, & !< name of each microstructure
texture_name !< name of each texture
2018-06-27 00:24:54 +05:30
! ToDo: Remove, use size(config_phase) etc
2019-03-13 02:59:03 +05:30
integer, public, protected :: &
material_Nphase, & !< number of phases
2019-03-09 05:37:26 +05:30
material_Nhomogenization !< number of homogenizations
public :: &
2018-08-22 18:00:51 +05:30
config_init, &
config_deallocate
contains
!--------------------------------------------------------------------------------------------------
!> @brief reads material.config and stores its content per part
!--------------------------------------------------------------------------------------------------
subroutine config_init
use prec, only: &
pStringLen
2018-07-16 15:12:58 +05:30
use DAMASK_interface, only: &
getSolverJobName
use IO, only: &
IO_read_ASCII, &
IO_error, &
IO_lc, &
2019-03-08 13:18:06 +05:30
IO_getTag
use debug, only: &
debug_level, &
debug_material, &
debug_levelBasic
implicit none
2019-03-13 02:59:03 +05:30
integer :: myDebug,i
character(len=pStringLen) :: &
line, &
part
character(len=pStringLen), dimension(:), allocatable :: fileContent
logical :: fileExists
write(6,'(/,a)') ' <<<+- config init -+>>>'
myDebug = debug_level(debug_material)
2019-03-08 13:18:06 +05:30
inquire(file=trim(getSolverJobName())//'.materialConfig',exist=fileExists)
if(fileExists) then
write(6,'(/,a)') ' reading '//trim(getSolverJobName())//'.materialConfig'; flush(6)
2019-03-13 02:18:33 +05:30
fileContent = read_materialConfig(trim(getSolverJobName())//'.materialConfig')
2018-07-16 15:12:58 +05:30
else
inquire(file='material.config',exist=fileExists)
2019-03-13 02:59:03 +05:30
if(.not. fileExists) call IO_error(100,ext_msg='material.config')
write(6,'(/,a)') ' reading material.config'; flush(6)
2019-03-13 02:18:33 +05:30
fileContent = read_materialConfig('material.config')
2018-07-16 15:12:58 +05:30
endif
2019-03-13 02:59:03 +05:30
do i = 1, size(fileContent)
2018-07-16 15:12:58 +05:30
line = trim(fileContent(i))
part = IO_lc(IO_getTag(line,'<','>'))
select case (trim(part))
2019-03-09 05:37:26 +05:30
case (trim('phase'))
2019-03-13 02:18:33 +05:30
call parse_materialConfig(phase_name,config_phase,line,fileContent(i+1:))
2019-03-09 05:37:26 +05:30
if (iand(myDebug,debug_levelBasic) /= 0) write(6,'(a)') ' Phase parsed'; flush(6)
2019-03-09 05:37:26 +05:30
case (trim('microstructure'))
2019-03-13 02:18:33 +05:30
call parse_materialConfig(microstructure_name,config_microstructure,line,fileContent(i+1:))
2019-03-09 05:37:26 +05:30
if (iand(myDebug,debug_levelBasic) /= 0) write(6,'(a)') ' Microstructure parsed'; flush(6)
2019-03-09 05:37:26 +05:30
case (trim('crystallite'))
2019-03-13 02:18:33 +05:30
call parse_materialConfig(crystallite_name,config_crystallite,line,fileContent(i+1:))
2019-03-09 05:37:26 +05:30
if (iand(myDebug,debug_levelBasic) /= 0) write(6,'(a)') ' Crystallite parsed'; flush(6)
2019-03-09 05:37:26 +05:30
case (trim('homogenization'))
2019-03-13 02:18:33 +05:30
call parse_materialConfig(homogenization_name,config_homogenization,line,fileContent(i+1:))
2019-03-09 05:37:26 +05:30
if (iand(myDebug,debug_levelBasic) /= 0) write(6,'(a)') ' Homogenization parsed'; flush(6)
2019-03-09 05:37:26 +05:30
case (trim('texture'))
2019-03-13 02:18:33 +05:30
call parse_materialConfig(texture_name,config_texture,line,fileContent(i+1:))
2019-03-09 05:37:26 +05:30
if (iand(myDebug,debug_levelBasic) /= 0) write(6,'(a)') ' Texture parsed'; flush(6)
end select
enddo
2019-03-08 13:18:06 +05:30
2018-06-27 00:24:54 +05:30
material_Nhomogenization = size(config_homogenization)
2019-03-08 13:18:06 +05:30
material_Nphase = size(config_phase)
2019-03-13 02:59:03 +05:30
if (material_Nhomogenization < 1) call IO_error(160,ext_msg='<homogenization>')
if (size(config_microstructure) < 1) call IO_error(160,ext_msg='<microstructure>')
if (size(config_crystallite) < 1) call IO_error(160,ext_msg='<crystallite>')
if (material_Nphase < 1) call IO_error(160,ext_msg='<phase>')
if (size(config_texture) < 1) call IO_error(160,ext_msg='<texture>')
inquire(file='numerics.config', exist=fileExists)
if (fileExists) then
write(6,'(/,a)') ' reading numerics.config'; flush(6)
fileContent = IO_read_ASCII('numerics.config')
call parse_debugAndNumericsConfig(config_numerics,fileContent)
endif
inquire(file='debug.config', exist=fileExists)
if (fileExists) then
write(6,'(/,a)') ' reading debug.config'; flush(6)
fileContent = IO_read_ASCII('debug.config')
call parse_debugAndNumericsConfig(config_debug,fileContent)
endif
2018-06-09 17:18:37 +05:30
2019-03-13 02:18:33 +05:30
contains
!--------------------------------------------------------------------------------------------------
!> @brief reads material.config
!! Recursion is triggered by "{path/to/inputfile}" in a line
!--------------------------------------------------------------------------------------------------
recursive function read_materialConfig(fileName,cnt) result(fileContent)
use IO, only: &
IO_warning
implicit none
character(len=*), intent(in) :: fileName
integer, intent(in), optional :: cnt !< recursion counter
character(len=pStringLen), dimension(:), allocatable :: fileContent !< file content, separated per lines
character(len=pStringLen), dimension(:), allocatable :: includedContent
character(len=pStringLen) :: line
character(len=pStringLen), parameter :: dummy = 'https://damask.mpie.de' !< to fill up remaining array
2019-03-13 02:18:33 +05:30
character(len=:), allocatable :: rawData
integer :: &
2019-03-13 02:18:33 +05:30
fileLength, &
fileUnit, &
startPos, endPos, &
myTotalLines, & !< # lines read from file without include statements
l,i, &
myStat
logical :: warned
if (present(cnt)) then
2019-03-13 02:59:03 +05:30
if (cnt>10) call IO_error(106,ext_msg=trim(fileName))
2019-03-13 02:18:33 +05:30
endif
!--------------------------------------------------------------------------------------------------
! read data as stream
inquire(file = fileName, size=fileLength)
if (fileLength == 0) then
allocate(fileContent(0))
return
endif
open(newunit=fileUnit, file=fileName, access='stream',&
status='old', position='rewind', action='read',iostat=myStat)
2019-03-13 02:59:03 +05:30
if(myStat /= 0) call IO_error(100,ext_msg=trim(fileName))
2019-03-13 02:18:33 +05:30
allocate(character(len=fileLength)::rawData)
read(fileUnit) rawData
close(fileUnit)
!--------------------------------------------------------------------------------------------------
! count lines to allocate string array
myTotalLines = 1
do l=1, len(rawData)
2019-03-13 02:18:33 +05:30
if (rawData(l:l) == new_line('')) myTotalLines = myTotalLines+1
enddo
allocate(fileContent(myTotalLines))
!--------------------------------------------------------------------------------------------------
! split raw data at end of line and handle includes
warned = .false.
startPos = 1
l = 1
2019-03-13 02:18:33 +05:30
do while (l <= myTotalLines)
endPos = merge(startPos + scan(rawData(startPos:),new_line('')) - 2,len(rawData),l /= myTotalLines)
if (endPos - startPos > pStringLen -1) then
line = rawData(startPos:startPos+pStringLen-1)
2019-03-13 02:18:33 +05:30
if (.not. warned) then
call IO_warning(207,ext_msg=trim(fileName),el=l)
2019-03-13 02:18:33 +05:30
warned = .true.
endif
else
line = rawData(startPos:endpos)
endif
startPos = endPos + 2 ! jump to next line start
2019-03-13 02:18:33 +05:30
recursion: if (scan(trim(adjustl(line)),'{') == 1 .and. scan(trim(line),'}') > 2) then
includedContent = read_materialConfig(trim(line(scan(line,'{')+1:scan(line,'}')-1)), &
merge(cnt,1,present(cnt))) ! to track recursion depth
fileContent = [ fileContent(1:l-1), includedContent, [(dummy,i=1,myTotalLines-l)] ] ! add content and grow array
myTotalLines = myTotalLines - 1 + size(includedContent)
l = l - 1 + size(includedContent)
2019-03-13 02:18:33 +05:30
else recursion
fileContent(l) = line
2019-03-13 02:59:03 +05:30
l = l + 1
2019-03-13 02:18:33 +05:30
endif recursion
enddo
end function read_materialConfig
!--------------------------------------------------------------------------------------------------
!> @brief parses the material.config file
!--------------------------------------------------------------------------------------------------
2019-03-13 02:18:33 +05:30
subroutine parse_materialConfig(sectionNames,part,line, &
fileContent)
implicit none
character(len=64), allocatable, dimension(:), intent(out) :: sectionNames
type(tPartitionedStringList), allocatable, dimension(:), intent(inout) :: part
character(len=pStringLen), intent(inout) :: line
character(len=pStringLen), dimension(:), intent(in) :: fileContent
2019-03-13 02:59:03 +05:30
integer, allocatable, dimension(:) :: partPosition ! position of [] tags + last line in section
integer :: i, j
logical :: echo
echo = .false.
2019-03-13 02:59:03 +05:30
if (allocated(part)) call IO_error(161,ext_msg=trim(line))
2018-08-22 18:00:51 +05:30
allocate(partPosition(0))
2019-03-13 02:59:03 +05:30
do i = 1, size(fileContent)
2018-07-16 15:12:58 +05:30
line = trim(fileContent(i))
if (IO_getTag(line,'<','>') /= '') exit
nextSection: if (IO_getTag(line,'[',']') /= '') then
2018-08-22 18:00:51 +05:30
partPosition = [partPosition, i]
cycle
endif nextSection
2019-03-13 02:59:03 +05:30
if (size(partPosition) < 1) &
echo = (trim(IO_getTag(line,'/','/')) == 'echo') .or. echo
enddo
2018-08-22 18:00:51 +05:30
allocate(sectionNames(size(partPosition)))
allocate(part(size(partPosition)))
partPosition = [partPosition, i] ! needed when actually storing content
2019-03-13 02:59:03 +05:30
do i = 1, size(partPosition) -1
sectionNames(i) = trim(adjustl(IO_getTag(fileContent(partPosition(i)),'[',']')))
2019-03-13 02:59:03 +05:30
do j = partPosition(i) + 1, partPosition(i+1) -1
call part(i)%add(trim(adjustl(fileContent(j))))
enddo
if (echo) then
write(6,*) 'section',i, '"'//trim(sectionNames(i))//'"'
2018-08-22 18:00:51 +05:30
call part(i)%show()
endif
enddo
2019-03-13 02:18:33 +05:30
end subroutine parse_materialConfig
!--------------------------------------------------------------------------------------------------
!> @brief parses the material.config file
!--------------------------------------------------------------------------------------------------
subroutine parse_debugAndNumericsConfig(config_list, &
fileContent)
implicit none
type(tPartitionedStringList), intent(out) :: config_list
character(len=pStringLen), dimension(:), intent(in) :: fileContent
integer :: i
2019-03-13 02:59:03 +05:30
do i = 1, size(fileContent)
call config_list%add(trim(adjustl(fileContent(i))))
enddo
end subroutine parse_debugAndNumericsConfig
2019-03-13 02:18:33 +05:30
end subroutine config_init
!--------------------------------------------------------------------------------------------------
!> @brief deallocates the linked lists that store the content of the configuration files
!--------------------------------------------------------------------------------------------------
subroutine config_deallocate(what)
use IO, only: &
IO_error
implicit none
character(len=*), intent(in) :: what
select case(trim(what))
case('material.config/phase')
deallocate(config_phase)
case('material.config/microstructure')
deallocate(config_microstructure)
case('material.config/crystallite')
deallocate(config_crystallite)
case('material.config/homogenization')
deallocate(config_homogenization)
case('material.config/texture')
deallocate(config_texture)
case('debug.config')
call config_debug%free
case('numerics.config')
call config_numerics%free
case default
2019-03-13 02:59:03 +05:30
call IO_error(0,ext_msg='config_deallocate')
end select
end subroutine config_deallocate
!##################################################################################################
! The folowing functions are part of the tPartitionedStringList object
!##################################################################################################
!--------------------------------------------------------------------------------------------------
!> @brief add element
!> @details Adds a string together with the start/end position of chunks in this string. The new
!! element is added at the end of the list. Empty strings are not added. All strings are converted
!! to lower case. The data is not stored in the new element but in the current.
!--------------------------------------------------------------------------------------------------
subroutine add(this,string)
use IO, only: &
IO_isBlank, &
IO_lc, &
IO_stringPos
implicit none
class(tPartitionedStringList), target, intent(in) :: this
character(len=*), intent(in) :: string
type(tPartitionedStringList), pointer :: new, temp
if (IO_isBlank(string)) return
allocate(new)
temp => this
do while (associated(temp%next))
temp => temp%next
enddo
temp%string%val = IO_lc (trim(string))
temp%string%pos = IO_stringPos(trim(string))
temp%next => new
end subroutine add
!--------------------------------------------------------------------------------------------------
!> @brief prints all elements
!> @details Strings are printed in order of insertion (FIFO)
!--------------------------------------------------------------------------------------------------
subroutine show(this)
implicit none
class(tPartitionedStringList), target, intent(in) :: this
type(tPartitionedStringList), pointer :: item
item => this
do while (associated(item%next))
write(6,'(a)') ' '//trim(item%string%val)
item => item%next
2019-01-06 04:25:10 +05:30
enddo
end subroutine show
!--------------------------------------------------------------------------------------------------
2018-08-22 18:00:51 +05:30
!> @brief empties list and frees associated memory
!> @details explicit interface to reset list. Triggers final statement (and following chain reaction)
!--------------------------------------------------------------------------------------------------
subroutine free(this)
implicit none
class(tPartitionedStringList), intent(inout) :: this
if(associated(this%next)) deallocate(this%next)
end subroutine free
!--------------------------------------------------------------------------------------------------
2018-08-22 18:00:51 +05:30
!> @brief empties list and frees associated memory
!> @details called when variable goes out of scope. Triggers chain reaction for list
!--------------------------------------------------------------------------------------------------
recursive subroutine finalize(this)
implicit none
type(tPartitionedStringList), intent(inout) :: this
if(associated(this%next)) deallocate(this%next)
2018-06-27 12:34:14 +05:30
end subroutine finalize
!--------------------------------------------------------------------------------------------------
2018-08-22 18:00:51 +05:30
!> @brief cleans entire array of linke lists
!> @details called when variable goes out of scope and deallocates the list at each array entry
!--------------------------------------------------------------------------------------------------
subroutine finalizeArray(this)
implicit none
integer :: i
type(tPartitionedStringList), intent(inout), dimension(:) :: this
2018-08-22 18:00:51 +05:30
type(tPartitionedStringList), pointer :: temp ! bug in Gfortran?
do i=1, size(this)
if (associated(this(i)%next)) then
temp => this(i)%next
!deallocate(this(i)) !internal compiler error: in gfc_build_final_call, at fortran/trans.c:975
deallocate(temp)
endif
enddo
end subroutine finalizeArray
!--------------------------------------------------------------------------------------------------
!> @brief reports wether a given key (string value at first position) exists in the list
!--------------------------------------------------------------------------------------------------
logical function keyExists(this,key)
use IO, only: &
IO_stringValue
implicit none
class(tPartitionedStringList), target, intent(in) :: this
character(len=*), intent(in) :: key
type(tPartitionedStringList), pointer :: item
keyExists = .false.
item => this
do while (associated(item%next) .and. .not. keyExists)
keyExists = trim(IO_stringValue(item%string%val,item%string%pos,1)) == trim(key)
item => item%next
2019-01-06 04:25:10 +05:30
enddo
end function keyExists
!--------------------------------------------------------------------------------------------------
!> @brief count number of key appearances
!> @details traverses list and counts each occurrence of specified key
!--------------------------------------------------------------------------------------------------
2019-03-13 02:59:03 +05:30
integer function countKeys(this,key)
use IO, only: &
IO_stringValue
implicit none
class(tPartitionedStringList), target, intent(in) :: this
character(len=*), intent(in) :: key
type(tPartitionedStringList), pointer :: item
2019-03-13 02:59:03 +05:30
countKeys = 0
item => this
do while (associated(item%next))
if (trim(IO_stringValue(item%string%val,item%string%pos,1)) == trim(key)) &
2019-03-13 02:59:03 +05:30
countKeys = countKeys + 1
item => item%next
2019-01-06 04:25:10 +05:30
enddo
end function countKeys
!--------------------------------------------------------------------------------------------------
2018-06-22 03:18:34 +05:30
!> @brief gets float value of for a given key from a linked list
!> @details gets the last value if the key occurs more than once. If key is not found exits with
!! error unless default is given
!--------------------------------------------------------------------------------------------------
real(pReal) function getFloat(this,key,defaultVal)
use IO, only : &
IO_error, &
IO_stringValue, &
IO_FloatValue
implicit none
class(tPartitionedStringList), target, intent(in) :: this
character(len=*), intent(in) :: key
real(pReal), intent(in), optional :: defaultVal
type(tPartitionedStringList), pointer :: item
logical :: found
found = present(defaultVal)
if (found) getFloat = defaultVal
item => this
do while (associated(item%next))
if (trim(IO_stringValue(item%string%val,item%string%pos,1)) == trim(key)) then
found = .true.
2019-03-13 02:59:03 +05:30
if (item%string%pos(1) < 2) call IO_error(143,ext_msg=key)
getFloat = IO_FloatValue(item%string%val,item%string%pos,2)
endif
item => item%next
2019-01-06 04:25:10 +05:30
enddo
2019-03-13 02:59:03 +05:30
if (.not. found) call IO_error(140,ext_msg=key)
end function getFloat
!--------------------------------------------------------------------------------------------------
2018-06-22 03:18:34 +05:30
!> @brief gets integer value of for a given key from a linked list
!> @details gets the last value if the key occurs more than once. If key is not found exits with
!! error unless default is given
!--------------------------------------------------------------------------------------------------
2019-03-13 02:59:03 +05:30
integer function getInt(this,key,defaultVal)
use IO, only: &
IO_error, &
IO_stringValue, &
IO_IntValue
implicit none
class(tPartitionedStringList), target, intent(in) :: this
character(len=*), intent(in) :: key
2019-03-13 02:59:03 +05:30
integer, intent(in), optional :: defaultVal
type(tPartitionedStringList), pointer :: item
logical :: found
found = present(defaultVal)
if (found) getInt = defaultVal
item => this
do while (associated(item%next))
if (trim(IO_stringValue(item%string%val,item%string%pos,1)) == trim(key)) then
found = .true.
2019-03-13 02:59:03 +05:30
if (item%string%pos(1) < 2) call IO_error(143,ext_msg=key)
getInt = IO_IntValue(item%string%val,item%string%pos,2)
endif
item => item%next
2019-01-06 04:25:10 +05:30
enddo
2019-03-13 02:59:03 +05:30
if (.not. found) call IO_error(140,ext_msg=key)
end function getInt
!--------------------------------------------------------------------------------------------------
2018-06-22 03:18:34 +05:30
!> @brief gets string value of for a given key from a linked list
!> @details gets the last value if the key occurs more than once. If key is not found exits with
!! error unless default is given. If raw is true, the the complete string is returned, otherwise
!! the individual chunks are returned
!--------------------------------------------------------------------------------------------------
character(len=65536) function getString(this,key,defaultVal,raw)
use IO, only: &
IO_error, &
IO_stringValue
implicit none
class(tPartitionedStringList), target, intent(in) :: this
character(len=*), intent(in) :: key
2019-03-23 11:22:26 +05:30
character(len=*), intent(in), optional :: defaultVal
logical, intent(in), optional :: raw
type(tPartitionedStringList), pointer :: item
logical :: found, &
whole
if (present(raw)) then
whole = raw
else
whole = .false.
endif
found = present(defaultVal)
if (found) then
2018-06-26 22:48:23 +05:30
getString = trim(defaultVal)
2019-03-23 11:22:26 +05:30
!if (len_trim(getString) /= len_trim(defaultVal)) call IO_error(0,ext_msg='getString')
endif
item => this
do while (associated(item%next))
if (trim(IO_stringValue(item%string%val,item%string%pos,1)) == trim(key)) then
found = .true.
2019-03-13 02:59:03 +05:30
if (item%string%pos(1) < 2) call IO_error(143,ext_msg=key)
if (whole) then
getString = trim(item%string%val(item%string%pos(4):)) ! raw string starting a second chunk
else
getString = IO_StringValue(item%string%val,item%string%pos,2)
endif
endif
item => item%next
2019-01-06 04:25:10 +05:30
enddo
2019-03-13 02:59:03 +05:30
if (.not. found) call IO_error(140,ext_msg=key)
end function getString
!--------------------------------------------------------------------------------------------------
2018-06-22 03:18:34 +05:30
!> @brief gets array of float values of for a given key from a linked list
!> @details for cumulative keys, "()", values from all occurrences are return. Otherwise only all
!! values from the last occurrence. If key is not found exits with error unless default is given.
!--------------------------------------------------------------------------------------------------
function getFloats(this,key,defaultVal,requiredSize)
2018-06-22 03:18:34 +05:30
use IO, only: &
IO_error, &
IO_stringValue, &
IO_FloatValue
implicit none
2018-06-22 03:18:34 +05:30
real(pReal), dimension(:), allocatable :: getFloats
class(tPartitionedStringList), target, intent(in) :: this
2018-06-22 03:18:34 +05:30
character(len=*), intent(in) :: key
real(pReal), dimension(:), intent(in), optional :: defaultVal
2019-03-13 02:59:03 +05:30
integer, intent(in), optional :: requiredSize
2018-06-22 03:18:34 +05:30
type(tPartitionedStringList), pointer :: item
2019-03-13 02:59:03 +05:30
integer :: i
logical :: found, &
2018-06-22 03:18:34 +05:30
cumulative
cumulative = (key(1:1) == '(' .and. key(len_trim(key):len_trim(key)) == ')')
found = .false.
2018-06-22 03:18:34 +05:30
allocate(getFloats(0))
item => this
do while (associated(item%next))
if (trim(IO_stringValue(item%string%val,item%string%pos,1)) == trim(key)) then
found = .true.
2018-06-26 22:48:23 +05:30
if (.not. cumulative) getFloats = [real(pReal)::]
2019-03-13 02:59:03 +05:30
if (item%string%pos(1) < 2) call IO_error(143,ext_msg=key)
do i = 2, item%string%pos(1)
2018-06-22 03:18:34 +05:30
getFloats = [getFloats,IO_FloatValue(item%string%val,item%string%pos,i)]
enddo
endif
item => item%next
2019-01-06 04:25:10 +05:30
enddo
if (.not. found) then
2019-03-13 02:59:03 +05:30
if (present(defaultVal)) then; getFloats = defaultVal; else; call IO_error(140,ext_msg=key); endif
endif
if (present(requiredSize)) then
if(requiredSize /= size(getFloats)) call IO_error(146,ext_msg=key)
endif
2018-06-22 03:18:34 +05:30
end function getFloats
!--------------------------------------------------------------------------------------------------
2018-06-22 03:18:34 +05:30
!> @brief gets array of integer values of for a given key from a linked list
!> @details for cumulative keys, "()", values from all occurrences are return. Otherwise only all
!! values from the last occurrence. If key is not found exits with error unless default is given.
!--------------------------------------------------------------------------------------------------
function getInts(this,key,defaultVal,requiredSize)
use IO, only: &
IO_error, &
IO_stringValue, &
IO_IntValue
implicit none
2019-03-13 02:59:03 +05:30
integer, dimension(:), allocatable :: getInts
class(tPartitionedStringList), target, intent(in) :: this
character(len=*), intent(in) :: key
2019-03-13 02:59:03 +05:30
integer, dimension(:), intent(in), optional :: defaultVal
integer, intent(in), optional :: requiredSize
type(tPartitionedStringList), pointer :: item
2019-03-13 02:59:03 +05:30
integer :: i
logical :: found, &
cumulative
cumulative = (key(1:1) == '(' .and. key(len_trim(key):len_trim(key)) == ')')
found = .false.
allocate(getInts(0))
item => this
do while (associated(item%next))
if (trim(IO_stringValue(item%string%val,item%string%pos,1)) == trim(key)) then
found = .true.
2019-03-13 02:59:03 +05:30
if (.not. cumulative) getInts = [integer::]
if (item%string%pos(1) < 2) call IO_error(143,ext_msg=key)
do i = 2, item%string%pos(1)
getInts = [getInts,IO_IntValue(item%string%val,item%string%pos,i)]
enddo
endif
item => item%next
2019-01-06 04:25:10 +05:30
enddo
if (.not. found) then
2019-03-13 02:59:03 +05:30
if (present(defaultVal)) then; getInts = defaultVal; else; call IO_error(140,ext_msg=key); endif
endif
if (present(requiredSize)) then
if(requiredSize /= size(getInts)) call IO_error(146,ext_msg=key)
endif
end function getInts
!--------------------------------------------------------------------------------------------------
2018-06-22 03:18:34 +05:30
!> @brief gets array of string values of for a given key from a linked list
!> @details for cumulative keys, "()", values from all occurrences are return. Otherwise only all
!! values from the last occurrence. If key is not found exits with error unless default is given.
!! If raw is true, the the complete string is returned, otherwise the individual chunks are returned
!--------------------------------------------------------------------------------------------------
function getStrings(this,key,defaultVal,raw)
use IO, only: &
IO_error, &
2018-06-22 03:18:34 +05:30
IO_StringValue
implicit none
2018-06-22 03:18:34 +05:30
character(len=65536),dimension(:), allocatable :: getStrings
class(tPartitionedStringList), target, intent(in) :: this
2018-06-22 03:18:34 +05:30
character(len=*), intent(in) :: key
character(len=65536),dimension(:), intent(in), optional :: defaultVal
logical, intent(in), optional :: raw
type(tPartitionedStringList), pointer :: item
character(len=65536) :: str
2019-03-13 02:59:03 +05:30
integer :: i
2018-06-22 03:18:34 +05:30
logical :: found, &
whole, &
2018-06-22 03:18:34 +05:30
cumulative
cumulative = (key(1:1) == '(' .and. key(len_trim(key):len_trim(key)) == ')')
if (present(raw)) then
whole = raw
else
whole = .false.
endif
found = .false.
item => this
do while (associated(item%next))
if (trim(IO_stringValue(item%string%val,item%string%pos,1)) == trim(key)) then
found = .true.
2018-06-22 03:18:34 +05:30
if (allocated(getStrings) .and. .not. cumulative) deallocate(getStrings)
2019-03-13 02:59:03 +05:30
if (item%string%pos(1) < 2) call IO_error(143,ext_msg=key)
2018-06-22 03:18:34 +05:30
notAllocated: if (.not. allocated(getStrings)) then
if (whole) then
str = item%string%val(item%string%pos(4):)
getStrings = [str]
else
2019-03-13 02:59:03 +05:30
str = IO_StringValue(item%string%val,item%string%pos,2)
2018-06-22 03:18:34 +05:30
allocate(getStrings(1),source=str)
2019-03-13 02:59:03 +05:30
do i=3,item%string%pos(1)
2018-06-22 03:18:34 +05:30
str = IO_StringValue(item%string%val,item%string%pos,i)
getStrings = [getStrings,str]
enddo
endif
else notAllocated
if (whole) then
2018-08-01 02:45:44 +05:30
str = item%string%val(item%string%pos(4):)
getStrings = [getStrings,str]
else
2019-03-13 02:59:03 +05:30
do i=2,item%string%pos(1)
2018-06-22 03:18:34 +05:30
str = IO_StringValue(item%string%val,item%string%pos,i)
getStrings = [getStrings,str]
enddo
endif
endif notAllocated
endif
item => item%next
2019-01-06 04:25:10 +05:30
enddo
if (.not. found) then
2019-03-13 02:59:03 +05:30
if (present(defaultVal)) then; getStrings = defaultVal; else; call IO_error(140,ext_msg=key); endif
endif
2018-06-22 03:18:34 +05:30
end function getStrings
end module config