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
|
|
|
|
!> @brief Interfacing between the spectral solver and the material subroutines provided
|
|
|
|
!! by DAMASK
|
|
|
|
!> @details Interfacing between the spectral solver and the material subroutines provided
|
|
|
|
!> 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.
|
2014-10-10 18:38:34 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2012-03-06 20:22:48 +05:30
|
|
|
module DAMASK_interface
|
2012-06-15 21:40:21 +05:30
|
|
|
use prec, only: &
|
|
|
|
pInt
|
2011-11-04 01:02:11 +05:30
|
|
|
implicit none
|
2012-03-06 20:22:48 +05:30
|
|
|
private
|
2015-08-04 20:34:53 +05:30
|
|
|
#include <petsc/finclude/petscsys.h>
|
2013-01-02 22:32:12 +05:30
|
|
|
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
|
|
|
|
character(len=1024), public, protected :: &
|
2012-06-15 21:40:21 +05:30
|
|
|
geometryFile = '', & !< parameter given for geometry file
|
|
|
|
loadCaseFile = '' !< parameter given for load case file
|
2013-01-02 22:32:12 +05:30
|
|
|
character(len=1024), private :: workingDirectory !< accessed by getSolverWorkingDirectoryName for compatibility reasons
|
2012-06-15 21:40:21 +05:30
|
|
|
|
2013-02-11 15:14:17 +05:30
|
|
|
public :: &
|
|
|
|
getSolverWorkingDirectoryName, &
|
|
|
|
getSolverJobName, &
|
|
|
|
DAMASK_interface_init
|
|
|
|
private :: &
|
|
|
|
storeWorkingDirectory, &
|
|
|
|
getGeometryFile, &
|
|
|
|
getLoadCaseFile, &
|
|
|
|
rectifyPath, &
|
|
|
|
makeRelativePath, &
|
|
|
|
IIO_stringValue, &
|
|
|
|
IIO_intValue, &
|
|
|
|
IIO_lc, &
|
|
|
|
IIO_stringPos
|
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
|
|
|
!--------------------------------------------------------------------------------------------------
|
2016-06-27 19:13:04 +05:30
|
|
|
subroutine DAMASK_interface_init()
|
2012-03-06 20:22:48 +05:30
|
|
|
use, intrinsic :: iso_fortran_env ! to get compiler_version and compiler_options (at least for gfortran 4.6 at the moment)
|
2016-09-20 10:38:31 +05:30
|
|
|
use system_routines, only: &
|
|
|
|
getHostName
|
2011-11-04 01:02:11 +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
|
|
|
implicit none
|
2012-04-11 22:58:08 +05:30
|
|
|
character(len=1024) :: &
|
|
|
|
commandLine, & !< command line call as string
|
2013-01-02 22:32:12 +05:30
|
|
|
loadCaseArg ='', & !< -l argument given to DAMASK_spectral.exe
|
|
|
|
geometryArg ='', & !< -g argument given to DAMASK_spectral.exe
|
|
|
|
workingDirArg ='', & !< -w argument given to DAMASK_spectral.exe
|
|
|
|
hostName, & !< name of machine on which DAMASK_spectral.exe is execute (might require export HOSTNAME)
|
|
|
|
userName, & !< name of user calling DAMASK_spectral.exe
|
2012-06-15 21:40:21 +05:30
|
|
|
tag
|
2012-04-11 22:58:08 +05:30
|
|
|
integer :: &
|
2014-10-10 18:38:34 +05:30
|
|
|
i, &
|
2017-04-25 16:04:14 +05:30
|
|
|
#ifdef _OPENMP
|
2016-06-27 15:53:42 +05:30
|
|
|
threadLevel, &
|
2017-04-25 16:04:14 +05:30
|
|
|
#endif
|
2016-09-03 17:55:45 +05:30
|
|
|
worldrank = 0, &
|
|
|
|
worldsize = 0
|
2015-08-28 13:08:48 +05:30
|
|
|
integer, allocatable, dimension(:) :: &
|
|
|
|
chunkPos
|
2012-04-11 22:58:08 +05:30
|
|
|
integer, dimension(8) :: &
|
|
|
|
dateAndTime ! type default integer
|
2014-10-10 18:38:34 +05:30
|
|
|
PetscErrorCode :: ierr
|
2016-09-20 10:38:31 +05:30
|
|
|
logical :: error
|
2013-02-13 23:24:56 +05:30
|
|
|
external :: &
|
2014-10-10 18:38:34 +05:30
|
|
|
quit,&
|
|
|
|
MPI_Comm_rank,&
|
2016-09-03 17:55:45 +05:30
|
|
|
MPI_Comm_size,&
|
2013-02-13 23:24:56 +05:30
|
|
|
PETScInitialize, &
|
2016-06-27 15:53:42 +05:30
|
|
|
MPI_Init_Thread, &
|
2013-02-13 23:24:56 +05:30
|
|
|
MPI_abort
|
2013-01-02 22:32:12 +05:30
|
|
|
|
2016-06-29 17:10:54 +05:30
|
|
|
open(6, encoding='UTF-8') ! for special characters in output
|
|
|
|
|
2012-11-06 21:30:51 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! PETSc Init
|
2016-06-27 15:53:42 +05:30
|
|
|
#ifdef _OPENMP
|
|
|
|
call MPI_Init_Thread(MPI_THREAD_FUNNELED,threadLevel,ierr);CHKERRQ(ierr) ! in case of OpenMP, don't rely on PETScInitialize doing MPI init
|
|
|
|
if (threadLevel<MPI_THREAD_FUNNELED) then
|
2016-09-03 17:55:45 +05:30
|
|
|
write(6,'(a)') ' MPI library does not support OpenMP'
|
2016-06-27 15:53:42 +05:30
|
|
|
call quit(1_pInt)
|
|
|
|
endif
|
|
|
|
#endif
|
|
|
|
call PetscInitialize(PETSC_NULL_CHARACTER,ierr) ! according to PETSc manual, that should be the first line in the code
|
2012-11-06 21:30:51 +05:30
|
|
|
CHKERRQ(ierr) ! this is a macro definition, it is case sensitive
|
2014-10-10 18:38:34 +05:30
|
|
|
call MPI_Comm_rank(PETSC_COMM_WORLD,worldrank,ierr);CHKERRQ(ierr)
|
2016-10-24 14:03:01 +05:30
|
|
|
call MPI_Comm_size(PETSC_COMM_WORLD,worldsize,ierr);CHKERRQ(ierr)
|
2014-10-10 18:38:34 +05:30
|
|
|
mainProcess: if (worldrank == 0) then
|
2016-06-28 03:06:09 +05:30
|
|
|
if (output_unit /= 6) then
|
2016-09-03 17:55:45 +05:30
|
|
|
write(output_unit,'(a)') ' STDOUT != 6'
|
2016-06-28 03:06:09 +05:30
|
|
|
call quit(1_pInt)
|
|
|
|
endif
|
2017-02-13 03:29:14 +05:30
|
|
|
if (error_unit /= 0) then
|
2017-04-10 21:57:07 +05:30
|
|
|
write(output_unit,'(a)') ' STDERR != 0'
|
|
|
|
call quit(1_pInt)
|
|
|
|
endif
|
|
|
|
if (PETSC_VERSION_MAJOR /= 3 .or. PETSC_VERSION_MINOR /= 7) then
|
|
|
|
write(6,'(a,2(i1.1,a))') 'PETSc ',PETSC_VERSION_MAJOR,'.',PETSC_VERSION_MINOR,'.x not supported'
|
2017-02-13 03:29:14 +05:30
|
|
|
call quit(1_pInt)
|
|
|
|
endif
|
2016-06-28 02:45:41 +05:30
|
|
|
else mainProcess
|
2016-06-29 17:10:54 +05:30
|
|
|
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
|
2014-10-10 18:38:34 +05:30
|
|
|
endif mainProcess
|
2016-06-29 17:10:54 +05:30
|
|
|
|
|
|
|
call date_and_time(values = dateAndTime)
|
|
|
|
write(6,'(/,a)') ' <<<+- DAMASK_spectral -+>>>'
|
|
|
|
write(6,'(/,a)') ' Version: '//DAMASKVERSION
|
|
|
|
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)
|
2016-09-03 17:55:45 +05:30
|
|
|
write(6,'(/,a,i4.1)') ' MPI processes: ',worldsize
|
2016-06-29 17:10:54 +05:30
|
|
|
write(6,'(/,a)') ' <<<+- DAMASK_interface init -+>>>'
|
|
|
|
#include "compilation_info.f90"
|
2016-07-01 18:37:10 +05:30
|
|
|
|
2016-06-27 19:13:04 +05:30
|
|
|
call get_command(commandLine)
|
|
|
|
chunkPos = IIO_stringPos(commandLine)
|
|
|
|
do i = 1, chunkPos(1)
|
|
|
|
tag = IIO_lc(IIO_stringValue(commandLine,chunkPos,i)) ! extract key
|
|
|
|
select case(tag)
|
|
|
|
case ('-h','--help')
|
2016-09-20 10:38:31 +05:30
|
|
|
write(6,'(a)') ' #######################################################################'
|
|
|
|
write(6,'(a)') ' DAMASK_spectral:'
|
|
|
|
write(6,'(a)') ' The spectral method boundary value problem solver for'
|
|
|
|
write(6,'(a)') ' 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.geom'
|
|
|
|
write(6,'(a)') ' Specifies the location of the geometry definition file,'
|
|
|
|
write(6,'(a)') ' if no extension is given, .geom will be appended.'
|
|
|
|
write(6,'(a)') ' "PathToGeomFile" will be the working directory if not specified'
|
|
|
|
write(6,'(a)') ' via --workingdir.'
|
|
|
|
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 "numerics.config" in that directory.'
|
|
|
|
write(6,'(/,a)')' --load PathToLoadFile/NameOfLoadFile.load'
|
|
|
|
write(6,'(a)') ' Specifies the location of the load case definition file,'
|
|
|
|
write(6,'(a)') ' if no extension is given, .load will be appended.'
|
|
|
|
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)') ' "PathToGeomFile".'
|
|
|
|
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"'
|
2016-10-19 01:54:23 +05:30
|
|
|
write(6,'(a)')' and "debug.config" in that directory.'
|
2016-09-20 10:38:31 +05:30
|
|
|
write(6,'(/,a)')' --restart XX'
|
|
|
|
write(6,'(a)') ' Reads in total increment No. XX-1 and continues to'
|
|
|
|
write(6,'(a)') ' calculate total increment No. XX.'
|
|
|
|
write(6,'(a)') ' Appends to existing results file '
|
|
|
|
write(6,'(a)') ' "NameOfGeom_NameOfLoadFile.spectralOut".'
|
|
|
|
write(6,'(a)') ' Works only if the restart information for total increment'
|
|
|
|
write(6,'(a)') ' No. XX-1 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'
|
2017-04-10 21:34:26 +05:30
|
|
|
call quit(0_pInt) ! normal Termination
|
2016-06-27 19:13:04 +05:30
|
|
|
case ('-l', '--load', '--loadcase')
|
|
|
|
loadcaseArg = IIO_stringValue(commandLine,chunkPos,i+1_pInt)
|
|
|
|
case ('-g', '--geom', '--geometry')
|
|
|
|
geometryArg = IIO_stringValue(commandLine,chunkPos,i+1_pInt)
|
|
|
|
case ('-w', '-d', '--wd', '--directory', '--workingdir', '--workingdirectory')
|
|
|
|
workingDirArg = IIO_stringValue(commandLine,chunkPos,i+1_pInt)
|
|
|
|
case ('-r', '--rs', '--restart')
|
|
|
|
spectralRestartInc = IIO_IntValue(commandLine,chunkPos,i+1_pInt)
|
|
|
|
appendToOutFile = .true.
|
|
|
|
end select
|
|
|
|
enddo
|
2012-06-15 21:40:21 +05:30
|
|
|
|
2013-01-02 22:32:12 +05:30
|
|
|
if (len(trim(loadcaseArg)) == 0 .or. len(trim(geometryArg)) == 0) then
|
|
|
|
write(6,'(a)') ' Please specify geometry AND load case (-h for help)'
|
2012-06-15 21:40:21 +05:30
|
|
|
call quit(1_pInt)
|
2012-02-21 21:34:16 +05:30
|
|
|
endif
|
2012-01-30 19:22:41 +05:30
|
|
|
|
2016-05-05 18:41:28 +05:30
|
|
|
workingDirectory = trim(storeWorkingDirectory(trim(workingDirArg),trim(geometryArg)))
|
2013-01-02 22:32:12 +05:30
|
|
|
geometryFile = getGeometryFile(geometryArg)
|
|
|
|
loadCaseFile = getLoadCaseFile(loadCaseArg)
|
2012-06-15 21:40:21 +05:30
|
|
|
|
|
|
|
call get_environment_variable('USER',userName)
|
2016-09-20 10:38:31 +05:30
|
|
|
error = getHostName(hostName)
|
|
|
|
write(6,'(a,a)') ' Host name: ', trim(hostName)
|
|
|
|
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)') ' Loadcase argument: ', trim(loadcaseArg)
|
|
|
|
write(6,'(a,a)') ' Working directory: ', trim(getSolverWorkingDirectoryName())
|
|
|
|
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
|
|
|
|
write(6,'(a,l1,/)') ' Append to result file: ', appendToOutFile
|
2011-11-04 01:02:11 +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
|
|
|
|
2012-03-06 20:22:48 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2013-09-14 16:29:35 +05:30
|
|
|
!> @brief extract working directory from given argument or from location of geometry file,
|
|
|
|
!! possibly converting relative arguments to absolut path
|
|
|
|
!> @todo change working directory with call chdir(storeWorkingDirectory)?
|
2012-03-06 20:22:48 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2013-01-02 22:32:12 +05:30
|
|
|
character(len=1024) function storeWorkingDirectory(workingDirectoryArg,geometryArg)
|
2016-03-12 01:29:14 +05:30
|
|
|
use system_routines, only: &
|
|
|
|
isDirectory, &
|
2016-05-05 16:30:46 +05:30
|
|
|
getCWD
|
2011-11-04 01:02:11 +05:30
|
|
|
|
|
|
|
implicit none
|
2013-09-14 16:29:35 +05:30
|
|
|
character(len=*), intent(in) :: workingDirectoryArg !< working directory argument
|
|
|
|
character(len=*), intent(in) :: geometryArg !< geometry argument
|
2013-01-02 22:32:12 +05:30
|
|
|
character(len=1024) :: cwd
|
2016-03-12 01:29:14 +05:30
|
|
|
logical :: error
|
2013-05-13 19:40:48 +05:30
|
|
|
external :: quit
|
2016-03-12 01:29:14 +05:30
|
|
|
|
2016-03-22 01:39:45 +05:30
|
|
|
wdGiven: if (len(workingDirectoryArg)>0) then
|
2017-04-10 21:34:26 +05:30
|
|
|
absolutePath: if (workingDirectoryArg(1:1) == '/') then
|
2013-01-02 22:32:12 +05:30
|
|
|
storeWorkingDirectory = workingDirectoryArg
|
2016-05-05 18:41:28 +05:30
|
|
|
else absolutePath
|
|
|
|
error = getCWD(cwd)
|
2016-05-05 19:46:21 +05:30
|
|
|
if (error) call quit(1_pInt)
|
2017-04-10 21:34:26 +05:30
|
|
|
storeWorkingDirectory = trim(cwd)//'/'//workingDirectoryArg
|
2016-05-05 18:41:28 +05:30
|
|
|
endif absolutePath
|
2017-04-10 21:34:26 +05:30
|
|
|
if (storeWorkingDirectory(len(trim(storeWorkingDirectory)):len(trim(storeWorkingDirectory))) /= '/') &
|
|
|
|
storeWorkingDirectory = trim(storeWorkingDirectory)//'/' ! if path seperator is not given, append it
|
2016-03-22 01:39:45 +05:30
|
|
|
else wdGiven
|
2017-04-10 21:34:26 +05:30
|
|
|
if (geometryArg(1:1) == '/') then ! absolute path given as command line argument
|
|
|
|
storeWorkingDirectory = geometryArg(1:scan(geometryArg,'/',back=.true.))
|
2013-01-02 22:32:12 +05:30
|
|
|
else
|
2016-05-18 11:22:19 +05:30
|
|
|
error = getCWD(cwd) ! relative path given as command line argument
|
2016-05-05 19:46:21 +05:30
|
|
|
if (error) call quit(1_pInt)
|
2017-04-10 21:34:26 +05:30
|
|
|
storeWorkingDirectory = trim(cwd)//'/'//geometryArg(1:scan(geometryArg,'/',back=.true.))
|
2013-01-02 22:32:12 +05:30
|
|
|
endif
|
2016-03-22 01:39:45 +05:30
|
|
|
endif wdGiven
|
2016-05-18 11:22:19 +05:30
|
|
|
|
2016-05-05 18:41:28 +05:30
|
|
|
storeWorkingDirectory = trim(rectifyPath(storeWorkingDirectory))
|
2016-05-18 11:22:19 +05:30
|
|
|
if(.not. isDirectory(trim(storeWorkingDirectory))) then ! check if the directory exists
|
|
|
|
write(6,'(a20,a,a16)') ' working directory "',trim(storeWorkingDirectory),'" does not exist'
|
|
|
|
call quit(1_pInt)
|
|
|
|
endif
|
2016-03-12 01:29:14 +05:30
|
|
|
|
2013-01-02 22:32:12 +05:30
|
|
|
end function storeWorkingDirectory
|
|
|
|
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief simply returns the private string workingDir
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
character(len=1024) function getSolverWorkingDirectoryName()
|
|
|
|
|
|
|
|
implicit none
|
2016-09-20 10:38:31 +05:30
|
|
|
|
2013-01-02 22:32:12 +05:30
|
|
|
getSolverWorkingDirectoryName = workingDirectory
|
|
|
|
|
2012-03-06 20:22:48 +05:30
|
|
|
end function getSolverWorkingDirectoryName
|
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
|
|
|
|
2012-03-06 20:22:48 +05:30
|
|
|
implicit none
|
2012-02-16 00:28:38 +05:30
|
|
|
integer :: posExt,posSep
|
2012-06-15 21:40:21 +05:30
|
|
|
character(len=1024) :: tempString
|
2012-02-21 21:34:16 +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-06-15 21:40:21 +05:30
|
|
|
tempString = geometryFile
|
|
|
|
posExt = scan(tempString,'.',back=.true.)
|
2017-04-10 21:34:26 +05:30
|
|
|
posSep = scan(tempString,'/',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
|
|
|
|
2012-06-15 21:40:21 +05:30
|
|
|
getSolverJobName = tempString(posSep+1:posExt-1)
|
|
|
|
|
|
|
|
tempString = loadCaseFile
|
|
|
|
posExt = scan(tempString,'.',back=.true.)
|
2017-04-10 21:34:26 +05:30
|
|
|
posSep = scan(tempString,'/',back=.true.)
|
2012-06-15 21:40:21 +05:30
|
|
|
|
|
|
|
getSolverJobName = trim(getSolverJobName)//'_'//tempString(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)
|
2016-05-18 11:22:19 +05:30
|
|
|
use system_routines, only: &
|
|
|
|
getCWD
|
2011-02-07 20:05:42 +05:30
|
|
|
|
2012-03-06 20:22:48 +05:30
|
|
|
implicit none
|
2012-06-15 21:40:21 +05:30
|
|
|
character(len=1024), intent(in) :: &
|
|
|
|
geometryParameter
|
2016-05-18 11:22:19 +05:30
|
|
|
character(len=1024) :: &
|
|
|
|
cwd
|
2012-06-15 21:40:21 +05:30
|
|
|
integer :: posExt, posSep
|
2016-05-18 11:22:19 +05:30
|
|
|
logical :: error
|
2016-06-27 15:53:42 +05:30
|
|
|
external :: quit
|
2011-02-07 20:05:42 +05:30
|
|
|
|
2012-06-15 21:40:21 +05:30
|
|
|
getGeometryFile = geometryParameter
|
|
|
|
posExt = scan(getGeometryFile,'.',back=.true.)
|
2017-04-10 21:34:26 +05:30
|
|
|
posSep = scan(getGeometryFile,'/',back=.true.)
|
2012-06-15 21:40:21 +05:30
|
|
|
|
|
|
|
if (posExt <= posSep) getGeometryFile = trim(getGeometryFile)//('.geom') ! no extension present
|
2017-04-10 21:34:26 +05:30
|
|
|
if (scan(getGeometryFile,'/') /= 1) then ! relative path given as command line argument
|
2016-05-18 11:22:19 +05:30
|
|
|
error = getcwd(cwd)
|
2016-06-27 15:53:42 +05:30
|
|
|
if (error) call quit(1_pInt)
|
2017-04-10 21:34:26 +05:30
|
|
|
getGeometryFile = rectifyPath(trim(cwd)//'/'//getGeometryFile)
|
2012-06-15 21:40:21 +05:30
|
|
|
else
|
|
|
|
getGeometryFile = rectifyPath(getGeometryFile)
|
|
|
|
endif
|
2011-02-07 20:05:42 +05:30
|
|
|
|
2012-06-15 21:40:21 +05:30
|
|
|
getGeometryFile = makeRelativePath(getSolverWorkingDirectoryName(), getGeometryFile)
|
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)
|
2016-05-18 11:22:19 +05:30
|
|
|
use system_routines, only: &
|
|
|
|
getCWD
|
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
|
|
|
|
|
|
|
implicit none
|
2012-06-15 21:40:21 +05:30
|
|
|
character(len=1024), intent(in) :: &
|
|
|
|
loadCaseParameter
|
2016-05-18 11:22:19 +05:30
|
|
|
character(len=1024) :: &
|
|
|
|
cwd
|
2016-03-12 01:29:14 +05:30
|
|
|
integer :: posExt, posSep
|
2016-05-18 11:22:19 +05:30
|
|
|
logical :: error
|
2016-06-27 15:53:42 +05:30
|
|
|
external :: quit
|
2013-12-28 01:33:28 +05:30
|
|
|
|
2012-06-15 21:40:21 +05:30
|
|
|
getLoadCaseFile = loadcaseParameter
|
|
|
|
posExt = scan(getLoadCaseFile,'.',back=.true.)
|
2017-04-10 21:34:26 +05:30
|
|
|
posSep = scan(getLoadCaseFile,'/',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
|
|
|
|
2012-06-15 21:40:21 +05:30
|
|
|
if (posExt <= posSep) getLoadCaseFile = trim(getLoadCaseFile)//('.load') ! no extension present
|
2017-04-10 21:34:26 +05:30
|
|
|
if (scan(getLoadCaseFile,'/') /= 1) then ! relative path given as command line argument
|
2016-05-18 11:22:19 +05:30
|
|
|
error = getcwd(cwd)
|
2016-06-27 15:53:42 +05:30
|
|
|
if (error) call quit(1_pInt)
|
2017-04-10 21:34:26 +05:30
|
|
|
getLoadCaseFile = rectifyPath(trim(cwd)//'/'//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
|
|
|
else
|
2012-06-15 21:40:21 +05:30
|
|
|
getLoadCaseFile = rectifyPath(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
|
|
|
endif
|
|
|
|
|
2012-06-15 21:40:21 +05:30
|
|
|
getLoadCaseFile = makeRelativePath(getSolverWorkingDirectoryName(), getLoadCaseFile)
|
|
|
|
|
|
|
|
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
|
|
|
!--------------------------------------------------------------------------------------------------
|
2013-05-13 19:40:48 +05:30
|
|
|
!> @brief remove ../ and /./ from path
|
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)
|
|
|
|
|
|
|
|
implicit none
|
2012-02-16 00:28:38 +05:30
|
|
|
character(len=*) :: path
|
|
|
|
character(len=len_trim(path)) :: rectifyPath
|
2013-05-13 19:40:48 +05:30
|
|
|
integer :: i,j,k,l ! no pInt
|
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
|
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
|
|
|
l = len_trim(path)
|
|
|
|
rectifyPath = path
|
2012-02-16 00:28:38 +05:30
|
|
|
do i = l,3,-1
|
2017-04-10 21:34:26 +05:30
|
|
|
if (rectifyPath(i-2:i) == '/'//'.'//'/') &
|
2013-05-08 21:22:29 +05:30
|
|
|
rectifyPath(i-1:l) = rectifyPath(i+1: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
|
|
|
enddo
|
|
|
|
|
2013-05-13 19:40:48 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! remove ../ and corresponding directory from 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
|
|
|
l = len_trim(rectifyPath)
|
2017-04-10 21:34:26 +05:30
|
|
|
i = index(rectifyPath(i:l),'..'//'/')
|
2012-02-16 00:28:38 +05:30
|
|
|
j = 0
|
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
|
|
|
do while (i > j)
|
2017-04-10 21:34:26 +05:30
|
|
|
j = scan(rectifyPath(1:i-2),'/',back=.true.)
|
2012-02-16 00:28:38 +05:30
|
|
|
rectifyPath(j+1:l) = rectifyPath(i+3:l)//repeat(' ',2+i-j)
|
2017-04-10 21:34:26 +05:30
|
|
|
if (rectifyPath(j+1:j+1) == '/') then !search for '//' that appear in case of XXX/../../XXX
|
2012-01-13 20:52:42 +05:30
|
|
|
k = len_trim(rectifyPath)
|
2012-02-16 00:28:38 +05:30
|
|
|
rectifyPath(j+1:k-1) = rectifyPath(j+2:k)
|
2012-01-13 20:52:42 +05:30
|
|
|
rectifyPath(k:k) = ' '
|
|
|
|
endif
|
2017-04-10 21:34:26 +05:30
|
|
|
i = j+index(rectifyPath(j+1: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
|
|
|
enddo
|
2017-04-10 21:34:26 +05:30
|
|
|
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
|
|
|
|
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
|
|
|
|
|
|
|
implicit none
|
|
|
|
character (len=*) :: a,b
|
2012-02-16 00:28:38 +05:30
|
|
|
integer :: i,posLastCommonSlash,remainingSlashes !no pInt
|
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-02-16 00:28:38 +05:30
|
|
|
posLastCommonSlash = 0
|
|
|
|
remainingSlashes = 0
|
2012-02-21 21:34:16 +05:30
|
|
|
|
2012-02-16 00:28:38 +05:30
|
|
|
do i = 1, min(1024,len_trim(a),len_trim(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
|
|
|
if (a(i:i) /= b(i:i)) exit
|
2017-04-10 21:34:26 +05:30
|
|
|
if (a(i:i) == '/') posLastCommonSlash = i
|
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
|
|
|
enddo
|
2012-02-16 00:28:38 +05:30
|
|
|
do i = posLastCommonSlash+1,len_trim(a)
|
2017-04-10 21:34:26 +05:30
|
|
|
if (a(i:i) == '/') remainingSlashes = remainingSlashes + 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
|
|
|
enddo
|
2017-04-10 21:34:26 +05:30
|
|
|
makeRelativePath = repeat('..'//'/',remainingSlashes)//b(posLastCommonSlash+1:len_trim(b))
|
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
|
|
|
|
2012-02-21 21:34:16 +05:30
|
|
|
|
2013-01-02 22:32:12 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief taken from IO, check IO_stringValue for documentation
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2015-08-28 13:08:48 +05:30
|
|
|
pure function IIO_stringValue(string,chunkPos,myChunk)
|
2012-06-15 21:40:21 +05:30
|
|
|
|
|
|
|
implicit none
|
2015-08-28 13:08:48 +05:30
|
|
|
integer(pInt), dimension(:), intent(in) :: chunkPos !< positions of start and end of each tag/chunk in given string
|
|
|
|
integer(pInt), intent(in) :: myChunk !< position number of desired chunk
|
|
|
|
character(len=1+chunkPos(myChunk*2+1)-chunkPos(myChunk*2)) :: IIO_stringValue
|
|
|
|
character(len=*), intent(in) :: string !< raw input with known start and end of each chunk
|
2012-06-15 21:40:21 +05:30
|
|
|
|
2015-08-28 13:08:48 +05:30
|
|
|
|
|
|
|
valuePresent: if (myChunk > chunkPos(1) .or. myChunk < 1_pInt) then
|
2013-01-02 22:32:12 +05:30
|
|
|
IIO_stringValue = ''
|
2015-08-28 13:08:48 +05:30
|
|
|
else valuePresent
|
|
|
|
IIO_stringValue = string(chunkPos(myChunk*2):chunkPos(myChunk*2+1))
|
|
|
|
endif valuePresent
|
2012-06-15 21:40:21 +05:30
|
|
|
|
2013-01-02 22:32:12 +05:30
|
|
|
end function IIO_stringValue
|
|
|
|
|
2012-06-15 21:40:21 +05:30
|
|
|
|
2013-01-02 22:32:12 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2015-08-28 13:08:48 +05:30
|
|
|
!> @brief taken from IO, check IO_intValue for documentation
|
2013-01-02 22:32:12 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
2015-08-28 13:08:48 +05:30
|
|
|
integer(pInt) pure function IIO_intValue(string,chunkPos,myChunk)
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
character(len=*), intent(in) :: string !< raw input with known start and end of each chunk
|
|
|
|
integer(pInt), intent(in) :: myChunk !< position number of desired sub string
|
|
|
|
integer(pInt), dimension(:), intent(in) :: chunkPos !< positions of start and end of each tag/chunk in given string
|
2012-06-15 21:40:21 +05:30
|
|
|
|
|
|
|
|
2015-08-28 13:08:48 +05:30
|
|
|
valuePresent: if (myChunk > chunkPos(1) .or. myChunk < 1_pInt) then
|
2013-01-02 22:32:12 +05:30
|
|
|
IIO_intValue = 0_pInt
|
2015-08-28 13:08:48 +05:30
|
|
|
else valuePresent
|
|
|
|
read(UNIT=string(chunkPos(myChunk*2):chunkPos(myChunk*2+1)),ERR=100,FMT=*) IIO_intValue
|
|
|
|
endif valuePresent
|
2012-06-15 21:40:21 +05:30
|
|
|
return
|
2013-01-02 22:32:12 +05:30
|
|
|
100 IIO_intValue = huge(1_pInt)
|
|
|
|
|
|
|
|
end function IIO_intValue
|
2012-06-15 21:40:21 +05:30
|
|
|
|
|
|
|
|
2013-01-02 22:32:12 +05:30
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief taken from IO, check IO_lc for documentation
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2015-08-28 13:08:48 +05:30
|
|
|
pure function IIO_lc(string)
|
2012-06-15 21:40:21 +05:30
|
|
|
|
|
|
|
implicit none
|
2015-08-28 13:08:48 +05:30
|
|
|
character(len=*), intent(in) :: string !< string to convert
|
|
|
|
character(len=len(string)) :: IIO_lc
|
|
|
|
|
|
|
|
character(26), parameter :: LOWER = 'abcdefghijklmnopqrstuvwxyz'
|
|
|
|
character(26), parameter :: UPPER = 'ABCDEFGHIJKLMNOPQRSTUVWXYZ'
|
2012-06-15 21:40:21 +05:30
|
|
|
|
2013-01-02 22:32:12 +05:30
|
|
|
integer :: i,n ! no pInt (len returns default integer)
|
2012-06-15 21:40:21 +05:30
|
|
|
|
2015-08-28 13:08:48 +05:30
|
|
|
IIO_lc = string
|
|
|
|
do i=1,len(string)
|
|
|
|
n = index(UPPER,IIO_lc(i:i))
|
|
|
|
if (n/=0) IIO_lc(i:i) = LOWER(n:n)
|
2012-06-15 21:40:21 +05:30
|
|
|
enddo
|
|
|
|
|
2013-01-02 22:32:12 +05:30
|
|
|
end function IIO_lc
|
2012-06-15 21:40:21 +05:30
|
|
|
|
2013-01-02 22:32:12 +05:30
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
!> @brief taken from IO, check IO_stringPos for documentation
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
2015-08-28 13:08:48 +05:30
|
|
|
pure function IIO_stringPos(string)
|
2012-06-15 21:40:21 +05:30
|
|
|
|
|
|
|
implicit none
|
2015-08-28 13:08:48 +05:30
|
|
|
integer(pInt), dimension(:), allocatable :: IIO_stringPos
|
|
|
|
character(len=*), intent(in) :: string !< string in which chunks are searched for
|
2012-06-15 21:40:21 +05:30
|
|
|
|
2015-08-06 14:54:56 +05:30
|
|
|
character(len=*), parameter :: SEP=achar(44)//achar(32)//achar(9)//achar(10)//achar(13) ! comma and whitespaces
|
2013-05-13 19:40:48 +05:30
|
|
|
integer :: left, right ! no pInt (verify and scan return default integer)
|
2012-06-15 21:40:21 +05:30
|
|
|
|
2015-08-28 13:08:48 +05:30
|
|
|
allocate(IIO_stringPos(1), source=0_pInt)
|
2012-06-15 21:40:21 +05:30
|
|
|
right = 0
|
|
|
|
|
2015-08-06 14:54:56 +05:30
|
|
|
do while (verify(string(right+1:),SEP)>0)
|
|
|
|
left = right + verify(string(right+1:),SEP)
|
|
|
|
right = left + scan(string(left:),SEP) - 2
|
|
|
|
if ( string(left:left) == '#' ) exit
|
2015-08-28 13:08:48 +05:30
|
|
|
IIO_stringPos = [IIO_stringPos,int(left, pInt), int(right, pInt)]
|
2013-01-02 22:32:12 +05:30
|
|
|
IIO_stringPos(1) = IIO_stringPos(1)+1_pInt
|
2012-06-15 21:40:21 +05:30
|
|
|
enddo
|
|
|
|
|
2013-01-02 22:32:12 +05:30
|
|
|
end function IIO_stringPos
|
|
|
|
|
2012-02-21 21:34:16 +05:30
|
|
|
|
2012-03-06 20:22:48 +05:30
|
|
|
end module
|