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
|
2018-08-05 14:11:01 +05:30
|
|
|
!> by DAMASK. Interpretating the command line arguments to get load case, geometry file,
|
|
|
|
!> and working directory.
|
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
|
2018-05-17 15:34:21 +05:30
|
|
|
|
2011-11-04 01:02:11 +05:30
|
|
|
implicit none
|
2012-03-06 20:22:48 +05:30
|
|
|
private
|
2018-08-21 02:44:34 +05:30
|
|
|
integer(pInt), public, protected :: &
|
|
|
|
interface_restartInc = 0_pInt !< Increment at which calculation starts
|
2013-01-02 22:32:12 +05:30
|
|
|
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-02-11 15:14:17 +05:30
|
|
|
public :: &
|
|
|
|
getSolverJobName, &
|
|
|
|
DAMASK_interface_init
|
|
|
|
private :: &
|
2018-07-10 11:54:45 +05:30
|
|
|
setWorkingDirectory, &
|
2013-02-11 15:14:17 +05:30
|
|
|
getGeometryFile, &
|
|
|
|
getLoadCaseFile, &
|
|
|
|
rectifyPath, &
|
|
|
|
makeRelativePath, &
|
|
|
|
IIO_stringValue, &
|
|
|
|
IIO_intValue, &
|
|
|
|
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()
|
2017-10-05 20:05:34 +05:30
|
|
|
use, intrinsic :: &
|
|
|
|
iso_fortran_env
|
2018-05-17 15:34:21 +05:30
|
|
|
#include <petsc/finclude/petscsys.h>
|
2018-10-15 08:33:53 +05:30
|
|
|
#if defined(__GFORTRAN__) && __GNUC__ < 5
|
|
|
|
===================================================================================================
|
|
|
|
5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0
|
|
|
|
===================================================================================================
|
|
|
|
================== THIS VERSION OF DAMASK REQUIRES gfortran > 5.0 ==============================
|
|
|
|
====================== THIS VERSION OF DAMASK REQUIRES gfortran > 5.0 ==========================
|
|
|
|
========================= THIS VERSION OF DAMASK REQUIRES gfortran > 5.0 =======================
|
|
|
|
===================================================================================================
|
|
|
|
5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0 5.0
|
|
|
|
===================================================================================================
|
|
|
|
#endif
|
|
|
|
|
|
|
|
#if defined(__INTEL_COMPILER) && __INTEL_COMPILER < 1600
|
|
|
|
===================================================================================================
|
|
|
|
16.0 16.0 16.0 16.0 16.0 16.0 16.0 16.0 16.0 16.0 16.0 16.0 16.0 16.0 16.0 16.0 16.0 16.0 16.0
|
|
|
|
===================================================================================================
|
|
|
|
================== THIS VERSION OF DAMASK REQUIRES ifort > 16.0 ================================
|
|
|
|
====================== THIS VERSION OF DAMASK REQUIRES ifort > 16.0 ===========================
|
|
|
|
========================= THIS VERSION OF DAMASK REQUIRES ifort > 16.0 ========================
|
|
|
|
===================================================================================================
|
|
|
|
16.0 16.0 16.0 16.0 16.0 16.0 16.0 16.0 16.0 16.0 16.0 16.0 16.0 16.0 16.0 16.0 16.0 16.0 16.0
|
|
|
|
===================================================================================================
|
|
|
|
#endif
|
|
|
|
|
2018-09-15 23:38:47 +05:30
|
|
|
#if PETSC_VERSION_MAJOR!=3 || PETSC_VERSION_MINOR!=10
|
2018-08-20 21:27:15 +05:30
|
|
|
===================================================================================================
|
2018-09-15 23:38:47 +05:30
|
|
|
3.10.x 3.10.x 3.10.x 3.10.x 3.10.x 3.10.x 3.10.x 3.10.x 3.10.x 3.10.x 3.10.x 3.10.x 3.10.x 3.10.x
|
2018-08-20 21:27:15 +05:30
|
|
|
===================================================================================================
|
2018-09-15 23:38:47 +05:30
|
|
|
=================== THIS VERSION OF DAMASK REQUIRES PETSc 3.10.x ==============================
|
|
|
|
====================== THIS VERSION OF DAMASK REQUIRES PETSc 3.10.x ===========================
|
|
|
|
========================= THIS VERSION OF DAMASK REQUIRES PETSc 3.10.x ========================
|
2018-08-20 21:27:15 +05:30
|
|
|
===================================================================================================
|
2018-09-15 23:38:47 +05:30
|
|
|
3.10.x 3.10.x 3.10.x 3.10.x 3.10.x 3.10.x 3.10.x 3.10.x 3.10.x 3.10.x 3.10.x 3.10.x 3.10.x 3.10.x
|
2018-05-17 15:34:21 +05:30
|
|
|
===================================================================================================
|
|
|
|
#endif
|
2018-10-15 08:33:53 +05:30
|
|
|
|
2018-05-17 15:34:21 +05:30
|
|
|
use PETScSys
|
2016-09-20 10:38:31 +05:30
|
|
|
use system_routines, only: &
|
2018-08-21 02:44:34 +05:30
|
|
|
getHostName, &
|
|
|
|
getCWD
|
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
|
2018-08-20 19:29:13 +05:30
|
|
|
loadcaseArg = '', & !< -l argument given to the executable
|
|
|
|
geometryArg = '', & !< -g argument given to the executable
|
|
|
|
workingDirArg = '', & !< -w argument given to the executable
|
2018-08-21 02:44:34 +05:30
|
|
|
userName !< name of user calling the executable
|
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
|
2013-02-13 23:24:56 +05:30
|
|
|
external :: &
|
2018-09-21 11:55:35 +05:30
|
|
|
quit
|
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
|
2018-02-12 19:33:35 +05:30
|
|
|
#ifdef _OPENMP
|
2018-02-12 16:16:01 +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.
|
|
|
|
call MPI_Init_Thread(MPI_THREAD_FUNNELED,threadLevel,ierr);CHKERRQ(ierr)
|
2016-06-27 15:53:42 +05:30
|
|
|
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
|
2018-05-17 15:34:21 +05:30
|
|
|
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
|
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)
|
2018-08-20 19:29:13 +05:30
|
|
|
write(6,'(/,a)') ' <<<+- DAMASK_interface init -+>>>'
|
2018-07-11 01:10:01 +05:30
|
|
|
write(6,'(a,/)') ' Roters et al., Computational Materials Science, 2018'
|
2016-06-29 17:10:54 +05:30
|
|
|
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
|
|
|
#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)
|
2018-07-11 01:10:01 +05:30
|
|
|
do i = 2_pInt, chunkPos(1)
|
|
|
|
select case(IIO_stringValue(commandLine,chunkPos,i)) ! extract key
|
2016-06-27 19:13:04 +05:30
|
|
|
case ('-h','--help')
|
2016-09-20 10:38:31 +05:30
|
|
|
write(6,'(a)') ' #######################################################################'
|
2018-08-20 19:29:13 +05:30
|
|
|
write(6,'(a)') ' DAMASK Command Line Interface:'
|
|
|
|
write(6,'(a)') ' For PETSc-based solvers for the Düsseldorf Advanced Material Simulation Kit'
|
2016-09-20 10:38:31 +05:30
|
|
|
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:'
|
2018-08-20 19:29:13 +05:30
|
|
|
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.'
|
2016-09-20 10:38:31 +05:30
|
|
|
write(6,'(/,a)')' -----------------------------------------------------------------------'
|
|
|
|
write(6,'(a)') ' Optional arguments:'
|
|
|
|
write(6,'(/,a)')' --workingdirectory PathToWorkingDirectory'
|
2018-08-20 19:29:13 +05:30
|
|
|
write(6,'(a)') ' Specifies the working directory and overwrites the default ./'
|
2016-09-20 10:38:31 +05:30
|
|
|
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'
|
2018-02-17 03:11:07 +05:30
|
|
|
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'
|
2018-08-20 19:29:13 +05:30
|
|
|
write(6,'(a)') ' "NameOfGeom_NameOfLoadFile".'
|
2018-02-17 03:11:07 +05:30
|
|
|
write(6,'(a)') ' Works only if the restart information for increment XX'
|
|
|
|
write(6,'(a)') ' is available in the working directory.'
|
2016-09-20 10:38:31 +05:30
|
|
|
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')
|
2018-07-11 01:10:01 +05:30
|
|
|
if ( i < chunkPos(1)) loadcaseArg = trim(IIO_stringValue(commandLine,chunkPos,i+1_pInt))
|
2016-06-27 19:13:04 +05:30
|
|
|
case ('-g', '--geom', '--geometry')
|
2018-07-11 01:10:01 +05:30
|
|
|
if (i < chunkPos(1)) geometryArg = trim(IIO_stringValue(commandLine,chunkPos,i+1_pInt))
|
2016-06-27 19:13:04 +05:30
|
|
|
case ('-w', '-d', '--wd', '--directory', '--workingdir', '--workingdirectory')
|
2018-07-11 01:10:01 +05:30
|
|
|
if (i < chunkPos(1)) workingDirArg = trim(IIO_stringValue(commandLine,chunkPos,i+1_pInt))
|
2016-06-27 19:13:04 +05:30
|
|
|
case ('-r', '--rs', '--restart')
|
2018-07-11 01:10:01 +05:30
|
|
|
if (i < chunkPos(1)) then
|
2018-08-20 19:29:13 +05:30
|
|
|
interface_restartInc = IIO_IntValue(commandLine,chunkPos,i+1_pInt)
|
2018-07-11 01:10:01 +05:30
|
|
|
endif
|
2016-06-27 19:13:04 +05:30
|
|
|
end select
|
|
|
|
enddo
|
2018-07-11 01:10:01 +05:30
|
|
|
|
|
|
|
if (len_trim(loadcaseArg) == 0 .or. len_trim(geometryArg) == 0) then
|
2013-01-02 22:32:12 +05:30
|
|
|
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
|
|
|
|
2018-08-21 02:44:34 +05:30
|
|
|
if (len_trim(workingDirArg) > 0) call setWorkingDirectory(trim(workingDirArg))
|
2018-06-26 19:39:37 +05:30
|
|
|
geometryFile = getGeometryFile(geometryArg)
|
|
|
|
loadCaseFile = getLoadCaseFile(loadCaseArg)
|
2012-06-15 21:40:21 +05:30
|
|
|
|
|
|
|
call get_environment_variable('USER',userName)
|
2018-08-21 02:44:34 +05:30
|
|
|
! ToDo: https://stackoverflow.com/questions/8953424/how-to-get-the-username-in-c-c-in-linux
|
|
|
|
write(6,'(a,a)') ' Host name: ', trim(getHostName())
|
2018-02-17 03:11:07 +05:30
|
|
|
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)
|
2018-08-21 02:44:34 +05:30
|
|
|
write(6,'(a,a)') ' Working directory: ', trim(getCWD())
|
2018-02-17 03:11:07 +05:30
|
|
|
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())
|
2018-08-20 19:29:13 +05:30
|
|
|
if (interface_restartInc > 0_pInt) &
|
|
|
|
write(6,'(a,i6.6)') ' Restart from increment: ', interface_restartInc
|
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
|
|
|
|
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)
|
2018-06-18 19:49:03 +05:30
|
|
|
use system_routines, only: &
|
2018-07-10 11:54:45 +05:30
|
|
|
getCWD, &
|
|
|
|
setCWD
|
2018-06-18 19:49:03 +05:30
|
|
|
|
|
|
|
implicit none
|
|
|
|
character(len=*), intent(in) :: workingDirectoryArg !< working directory argument
|
2018-08-21 02:44:34 +05:30
|
|
|
character(len=1024) :: workingDirectory !< working directory argument
|
|
|
|
logical :: error
|
2018-09-27 11:53:30 +05:30
|
|
|
external :: quit
|
2018-08-21 02:44:34 +05:30
|
|
|
|
|
|
|
absolutePath: if (workingDirectoryArg(1:1) == '/') then
|
|
|
|
workingDirectory = workingDirectoryArg
|
|
|
|
else absolutePath
|
|
|
|
workingDirectory = getCWD()
|
|
|
|
workingDirectory = trim(workingDirectory)//'/'//workingDirectoryArg
|
|
|
|
endif absolutePath
|
2016-03-12 01:29:14 +05:30
|
|
|
|
2018-08-21 02:44:34 +05:30
|
|
|
workingDirectory = trim(rectifyPath(workingDirectory))
|
|
|
|
error = setCWD(trim(workingDirectory))
|
2018-07-10 13:53:21 +05:30
|
|
|
if(error) then
|
2018-08-21 02:44:34 +05:30
|
|
|
write(6,'(a20,a,a16)') ' working directory "',trim(workingDirectory),'" does not exist'
|
2018-07-10 13:53:21 +05:30
|
|
|
call quit(1_pInt)
|
|
|
|
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
|
|
|
|
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)
|
2018-08-21 02:44:34 +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
|
2018-08-21 02:44:34 +05:30
|
|
|
character(len=1024), intent(in) :: geometryParameter
|
2018-09-15 18:24:56 +05:30
|
|
|
logical :: file_exists
|
2018-09-27 11:53:30 +05:30
|
|
|
external :: quit
|
2011-02-07 20:05:42 +05:30
|
|
|
|
2018-06-29 19:06:12 +05:30
|
|
|
getGeometryFile = trim(geometryParameter)
|
2018-08-21 02:44:34 +05:30
|
|
|
if (scan(getGeometryFile,'/') /= 1) getGeometryFile = trim(getCWD())//'/'//trim(getGeometryFile)
|
|
|
|
getGeometryFile = makeRelativePath(trim(getCWD()), getGeometryFile)
|
2011-02-07 20:05:42 +05:30
|
|
|
|
2018-09-15 18:24:56 +05:30
|
|
|
inquire(file=trim(getGeometryFile), exist=file_exists)
|
|
|
|
if (.not. file_exists) then
|
|
|
|
write(6,'(a)') ' Geometry file does not exists ('//trim(getGeometryFile)//')'
|
|
|
|
call quit(1_pInt)
|
|
|
|
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)
|
2018-08-21 02:44:34 +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
|
2018-08-21 02:44:34 +05:30
|
|
|
character(len=1024), intent(in) :: loadCaseParameter
|
2018-09-15 18:24:56 +05:30
|
|
|
logical :: file_exists
|
2018-09-27 11:53:30 +05:30
|
|
|
external :: quit
|
2013-12-28 01:33:28 +05:30
|
|
|
|
2018-06-29 19:06:12 +05:30
|
|
|
getLoadCaseFile = trim(loadCaseParameter)
|
2018-08-21 02:44:34 +05:30
|
|
|
if (scan(getLoadCaseFile,'/') /= 1) getLoadCaseFile = trim(getCWD())//'/'//trim(getLoadCaseFile)
|
|
|
|
getLoadCaseFile = makeRelativePath(trim(getCWD()), getLoadCaseFile)
|
2012-06-15 21:40:21 +05:30
|
|
|
|
2018-09-15 18:24:56 +05:30
|
|
|
inquire(file=trim(getLoadCaseFile), exist=file_exists)
|
|
|
|
if (.not. file_exists) then
|
|
|
|
write(6,'(a)') ' Geometry file does not exists ('//trim(getLoadCaseFile)//')'
|
|
|
|
call quit(1_pInt)
|
|
|
|
endif
|
|
|
|
|
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)
|
|
|
|
|
|
|
|
implicit none
|
2012-02-16 00:28:38 +05:30
|
|
|
character(len=*) :: path
|
2018-08-22 21:22:00 +05:30
|
|
|
character(len=1024) :: 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
|
2018-08-22 21:22:00 +05:30
|
|
|
rectifyPath = trim(path)
|
|
|
|
l = len_trim(rectifyPath)
|
2012-02-16 00:28:38 +05:30
|
|
|
do i = l,3,-1
|
2018-07-10 13:53:21 +05:30
|
|
|
if (rectifyPath(i-2:i) == '/./') rectifyPath(i-1:l) = rectifyPath(i+1:l)//' '
|
|
|
|
enddo
|
|
|
|
|
|
|
|
!--------------------------------------------------------------------------------------------------
|
|
|
|
! remove // from path
|
2018-08-22 21:22:00 +05:30
|
|
|
l = len_trim(rectifyPath)
|
2018-07-10 13:53:21 +05:30
|
|
|
do i = l,2,-1
|
|
|
|
if (rectifyPath(i-1:i) == '//') rectifyPath(i-1:l) = rectifyPath(i: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)
|
2018-06-29 19:06:12 +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
|
2018-06-29 19:06:12 +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
|
2018-07-10 13:53:21 +05:30
|
|
|
character (len=*), intent(in) :: a,b
|
|
|
|
character (len=1024) :: a_cleaned,b_cleaned
|
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
|
2018-07-10 13:53:21 +05:30
|
|
|
a_cleaned = rectifyPath(trim(a)//'/')
|
|
|
|
b_cleaned = rectifyPath(b)
|
2012-02-21 21:34:16 +05:30
|
|
|
|
2018-07-10 13:53:21 +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
|
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
|
2018-07-10 13:53:21 +05:30
|
|
|
do i = posLastCommonSlash+1,len_trim(a_cleaned)
|
|
|
|
if (a_cleaned(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
|
2018-07-10 13:53:21 +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
|
|
|
|
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
|
2018-07-11 01:10:01 +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=chunkPos(myChunk*2+1)-chunkPos(myChunk*2)+1) :: 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
|
|
|
|
2018-07-11 01:10:01 +05:30
|
|
|
IIO_stringValue = string(chunkPos(myChunk*2):chunkPos(myChunk*2+1))
|
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_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-03-06 20:22:48 +05:30
|
|
|
end module
|