2011-04-07 12:50:28 +05:30
|
|
|
! Copyright 2011 Max-Planck-Institut für Eisenforschung GmbH
|
2011-04-04 19:39:54 +05:30
|
|
|
!
|
|
|
|
! This file is part of DAMASK,
|
2011-04-07 12:50:28 +05:30
|
|
|
! the Düsseldorf Advanced MAterial Simulation Kit.
|
2011-04-04 19:39:54 +05:30
|
|
|
!
|
|
|
|
! DAMASK is free software: you can redistribute it and/or modify
|
|
|
|
! it under the terms of the GNU General Public License as published by
|
|
|
|
! the Free Software Foundation, either version 3 of the License, or
|
|
|
|
! (at your option) any later version.
|
|
|
|
!
|
|
|
|
! DAMASK is distributed in the hope that it will be useful,
|
|
|
|
! but WITHOUT ANY WARRANTY; without even the implied warranty of
|
|
|
|
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
|
|
|
! GNU General Public License for more details.
|
|
|
|
!
|
|
|
|
! You should have received a copy of the GNU General Public License
|
|
|
|
! along with DAMASK. If not, see <http://www.gnu.org/licenses/>.
|
|
|
|
!
|
|
|
|
!##############################################################
|
2009-08-31 20:39:15 +05:30
|
|
|
!* $Id$
|
2007-03-20 19:25:22 +05:30
|
|
|
!##############################################################
|
|
|
|
MODULE IO
|
|
|
|
!##############################################################
|
|
|
|
|
|
|
|
CONTAINS
|
|
|
|
!---------------------------
|
2010-07-13 15:56:07 +05:30
|
|
|
! function IO_abaqus_assembleInputFile
|
2007-03-20 19:25:22 +05:30
|
|
|
! function IO_open_file(unit,relPath)
|
2011-05-28 15:14:43 +05:30
|
|
|
! function IO_open_inputFile(unit, model)
|
|
|
|
! function IO_open_logFile(unit)
|
2007-03-21 22:34:10 +05:30
|
|
|
! function IO_hybridIA(Nast,ODFfileName)
|
|
|
|
! private function hybridIA_reps(dV_V,steps,C)
|
2007-04-10 16:51:34 +05:30
|
|
|
! function IO_stringPos(line,maxN)
|
2007-03-20 19:25:22 +05:30
|
|
|
! function IO_stringValue(line,positions,pos)
|
|
|
|
! function IO_floatValue(line,positions,pos)
|
|
|
|
! function IO_intValue(line,positions,pos)
|
2007-04-10 16:51:34 +05:30
|
|
|
! function IO_fixedStringValue(line,ends,pos)
|
|
|
|
! function IO_fixedFloatValue(line,ends,pos)
|
|
|
|
! function IO_fixedFloatNoEValue(line,ends,pos)
|
|
|
|
! function IO_fixedIntValue(line,ends,pos)
|
2009-04-02 18:32:25 +05:30
|
|
|
! function IO_continousIntValues(unit,maxN)
|
2007-04-10 16:51:34 +05:30
|
|
|
! function IO_lc(line)
|
|
|
|
! subroutine IO_lcInplace(line)
|
2007-03-20 19:25:22 +05:30
|
|
|
! subroutine IO_error(ID)
|
2009-03-31 14:51:57 +05:30
|
|
|
! subroutine IO_warning(ID)
|
2007-03-20 19:25:22 +05:30
|
|
|
!---------------------------
|
|
|
|
|
2009-08-31 20:39:15 +05:30
|
|
|
!********************************************************************
|
|
|
|
! output version number
|
|
|
|
!********************************************************************
|
2009-10-12 21:31:49 +05:30
|
|
|
subroutine IO_init ()
|
|
|
|
|
|
|
|
!$OMP CRITICAL (write2out)
|
|
|
|
write(6,*)
|
|
|
|
write(6,*) '<<<+- IO init -+>>>'
|
|
|
|
write(6,*) '$Id$'
|
|
|
|
write(6,*)
|
|
|
|
call flush(6)
|
|
|
|
!$OMP END CRITICAL (write2out)
|
|
|
|
|
2009-09-18 21:07:14 +05:30
|
|
|
endsubroutine
|
2007-03-20 19:25:22 +05:30
|
|
|
|
2010-07-13 15:56:07 +05:30
|
|
|
|
|
|
|
|
|
|
|
!********************************************************************
|
|
|
|
! AP: 12.07.10
|
|
|
|
! create a new input file for abaqus simulations
|
|
|
|
! by removing all comment lines and including "include"s
|
|
|
|
!********************************************************************
|
|
|
|
recursive function IO_abaqus_assembleInputFile(unit1,unit2) result(createSuccess)
|
|
|
|
use prec
|
2011-05-11 22:31:03 +05:30
|
|
|
use DAMASK_interface
|
2010-07-13 15:56:07 +05:30
|
|
|
implicit none
|
|
|
|
|
|
|
|
character(len=300) line,fname
|
|
|
|
integer(pInt), intent(in) :: unit1, unit2
|
|
|
|
logical createSuccess,fexist
|
2011-07-18 14:45:20 +05:30
|
|
|
integer(pInt), parameter :: maxNchunks = 6
|
|
|
|
integer(pInt), dimension(1+2*maxNchunks) :: positions
|
|
|
|
|
2010-07-13 15:56:07 +05:30
|
|
|
|
|
|
|
do
|
|
|
|
read(unit2,'(A300)',END=220) line
|
2011-07-18 14:45:20 +05:30
|
|
|
! line = IO_lc(trim(line))
|
|
|
|
! do not change the whole line to lower case, file names in Linux are case sensitive!
|
|
|
|
positions = IO_stringPos(line,maxNchunks)
|
|
|
|
|
2010-07-13 15:56:07 +05:30
|
|
|
! call IO_lcInPlace(line)
|
2011-07-18 14:45:20 +05:30
|
|
|
if (IO_lc(IO_StringValue(line,positions,1))=='*include') then
|
2010-07-13 15:56:07 +05:30
|
|
|
fname = trim(getSolverWorkingDirectoryName())//trim(line(9+scan(line(9:),'='):))
|
|
|
|
inquire(file=fname, exist=fexist)
|
|
|
|
if (.not.(fexist)) then
|
openmp parallelization working again (at least for j2 and nonlocal constitutive model).
In order to keep it like that, please follow these simple rules:
DON'T use implicit array subscripts:
example: real, dimension(3,3) :: A,B
A(:,2) = B(:,1) <--- DON'T USE
A(1:3,2) = B(1:3,1) <--- BETTER USE
In many cases the use of explicit array subscripts is inevitable for parallelization. Additionally, it is an easy means to prevent memory leaks.
Enclose all write statements with the following:
!$OMP CRITICAL (write2out)
<your write statement>
!$OMP END CRITICAL (write2out)
Whenever you change something in the code and are not sure if it affects parallelization and leads to nonconforming behavior, please ask me and/or Franz to check this.
2011-03-17 16:16:17 +05:30
|
|
|
!$OMP CRITICAL (write2out)
|
|
|
|
write(6,*)'ERROR: file does not exist error in IO_abaqus_assembleInputFile'
|
|
|
|
write(6,*)'filename: ', trim(fname)
|
|
|
|
!$OMP END CRITICAL (write2out)
|
2010-07-13 15:56:07 +05:30
|
|
|
createSuccess = .false.
|
|
|
|
return
|
|
|
|
endif
|
|
|
|
open(unit2+1,err=200,status='old',file=fname)
|
|
|
|
if (IO_abaqus_assembleInputFile(unit1,unit2+1)) then
|
|
|
|
createSuccess=.true.
|
|
|
|
close(unit2+1)
|
|
|
|
else
|
|
|
|
createSuccess=.false.
|
|
|
|
return
|
|
|
|
endif
|
|
|
|
else if (line(1:2) /= '**') then
|
|
|
|
write(unit1,'(A)') trim(line)
|
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
|
|
|
|
220 createSuccess = .true.
|
|
|
|
return
|
|
|
|
|
|
|
|
200 createSuccess =.false.
|
2011-09-13 21:24:06 +05:30
|
|
|
|
2010-07-13 15:56:07 +05:30
|
|
|
end function
|
|
|
|
|
|
|
|
!***********************************************************
|
|
|
|
! check if the input file for Abaqus contains part info
|
|
|
|
!***********************************************************
|
|
|
|
function IO_abaqus_hasNoPart(unit)
|
|
|
|
|
|
|
|
use prec, only: pInt
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
integer(pInt) unit
|
|
|
|
integer(pInt), parameter :: maxNchunks = 1
|
|
|
|
integer(pInt), dimension(1+2*maxNchunks) :: pos
|
|
|
|
logical IO_abaqus_hasNoPart
|
|
|
|
character(len=300) line
|
|
|
|
|
|
|
|
IO_abaqus_hasNoPart = .true.
|
|
|
|
|
|
|
|
610 FORMAT(A300)
|
|
|
|
rewind(unit)
|
|
|
|
do
|
|
|
|
read(unit,610,END=620) line
|
|
|
|
pos = IO_stringPos(line,maxNchunks)
|
|
|
|
if (IO_lc(IO_stringValue(line,pos,1)) == '*part' ) then
|
|
|
|
IO_abaqus_hasNoPart = .false.
|
|
|
|
exit
|
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
|
2011-10-18 14:52:33 +05:30
|
|
|
620 endfunction
|
2010-07-13 15:56:07 +05:30
|
|
|
|
|
|
|
|
|
|
|
|
2007-03-20 19:25:22 +05:30
|
|
|
!********************************************************************
|
|
|
|
! open existing file to given unit
|
|
|
|
! path to file is relative to working directory
|
|
|
|
!********************************************************************
|
2009-06-15 18:41:21 +05:30
|
|
|
logical function IO_open_file(unit,relPath)
|
2007-03-20 19:25:22 +05:30
|
|
|
|
|
|
|
use prec, only: pInt
|
2011-05-11 22:31:03 +05:30
|
|
|
use DAMASK_interface
|
2007-03-20 19:25:22 +05:30
|
|
|
implicit none
|
|
|
|
|
2010-08-04 05:17:00 +05:30
|
|
|
character(len=*), parameter :: pathSep = achar(47)//achar(92) ! forward and backward slash
|
2007-03-20 19:25:22 +05:30
|
|
|
character(len=*) relPath
|
|
|
|
integer(pInt) unit
|
|
|
|
|
2010-07-13 15:56:07 +05:30
|
|
|
IO_open_file = .false.
|
|
|
|
|
2010-05-11 12:27:15 +05:30
|
|
|
open(unit,status='old',err=100,file=trim(getSolverWorkingDirectoryName())//relPath)
|
2007-03-20 19:25:22 +05:30
|
|
|
IO_open_file = .true.
|
2010-07-13 15:56:07 +05:30
|
|
|
|
2011-10-18 14:52:33 +05:30
|
|
|
100 endfunction
|
2007-03-20 19:25:22 +05:30
|
|
|
|
|
|
|
|
|
|
|
!********************************************************************
|
|
|
|
! open FEM inputfile to given unit
|
2010-07-13 15:56:07 +05:30
|
|
|
! AP: 12.07.10
|
|
|
|
! : changed the function to open *.inp_assembly, which is basically
|
|
|
|
! the input file without comment lines and possibly assembled includes
|
2007-03-20 19:25:22 +05:30
|
|
|
!********************************************************************
|
2011-05-28 15:14:43 +05:30
|
|
|
logical function IO_open_inputFile(unit,model)
|
2007-03-20 19:25:22 +05:30
|
|
|
|
|
|
|
use prec, only: pReal, pInt
|
2011-05-11 22:31:03 +05:30
|
|
|
use DAMASK_interface
|
2007-03-20 19:25:22 +05:30
|
|
|
implicit none
|
|
|
|
|
2009-07-22 21:37:19 +05:30
|
|
|
integer(pInt), intent(in) :: unit
|
2011-05-28 15:14:43 +05:30
|
|
|
character(len=*) model
|
2010-05-10 20:32:59 +05:30
|
|
|
|
2010-07-13 15:56:07 +05:30
|
|
|
IO_open_inputFile = .false.
|
|
|
|
|
|
|
|
if (FEsolver == 'Abaqus') then
|
|
|
|
open(unit+1,status='old',err=100,&
|
|
|
|
file=trim(getSolverWorkingDirectoryName())//&
|
2011-05-28 15:14:43 +05:30
|
|
|
trim(model)//InputFileExtension)
|
2010-07-13 15:56:07 +05:30
|
|
|
open(unit,err=100,file=trim(getSolverWorkingDirectoryName())//&
|
2011-05-28 15:14:43 +05:30
|
|
|
trim(model)//InputFileExtension//'_assembly')
|
2010-07-13 15:56:07 +05:30
|
|
|
IO_open_inputFile = IO_abaqus_assembleInputFile(unit,unit+1) ! strip comments and concatenate any "include"s
|
|
|
|
close(unit+1)
|
|
|
|
else
|
|
|
|
open(unit,status='old',err=100,file=trim(getSolverWorkingDirectoryName())//&
|
2011-05-28 15:14:43 +05:30
|
|
|
trim(model)//InputFileExtension)
|
2010-07-13 15:56:07 +05:30
|
|
|
IO_open_inputFile = .true.
|
|
|
|
endif
|
|
|
|
|
2011-10-18 14:52:33 +05:30
|
|
|
100 endfunction
|
2007-03-20 19:25:22 +05:30
|
|
|
|
|
|
|
|
2010-11-03 20:09:18 +05:30
|
|
|
!********************************************************************
|
|
|
|
! open FEM logfile to given unit
|
|
|
|
!********************************************************************
|
|
|
|
logical function IO_open_logFile(unit)
|
|
|
|
|
|
|
|
use prec, only: pReal, pInt
|
2011-05-11 22:31:03 +05:30
|
|
|
use DAMASK_interface
|
2010-11-03 20:09:18 +05:30
|
|
|
implicit none
|
|
|
|
|
|
|
|
integer(pInt), intent(in) :: unit
|
|
|
|
|
|
|
|
IO_open_logFile = .false.
|
|
|
|
|
|
|
|
open(unit,status='old',err=100,file=trim(getSolverWorkingDirectoryName())//&
|
|
|
|
trim(getSolverJobName())//LogFileExtension)
|
|
|
|
IO_open_logFile = .true.
|
|
|
|
|
2011-10-18 14:52:33 +05:30
|
|
|
100 endfunction
|
2010-11-03 20:09:18 +05:30
|
|
|
|
|
|
|
|
2009-07-22 21:37:19 +05:30
|
|
|
!********************************************************************
|
|
|
|
! open (write) file related to current job
|
|
|
|
! but with different extension to given unit
|
|
|
|
!********************************************************************
|
|
|
|
logical function IO_open_jobFile(unit,newExt)
|
|
|
|
|
|
|
|
use prec, only: pReal, pInt
|
2011-05-11 22:31:03 +05:30
|
|
|
use DAMASK_interface
|
2009-07-22 21:37:19 +05:30
|
|
|
implicit none
|
|
|
|
|
|
|
|
integer(pInt), intent(in) :: unit
|
|
|
|
character(*), intent(in) :: newExt
|
|
|
|
|
2010-07-13 15:56:07 +05:30
|
|
|
IO_open_jobFile = .false.
|
|
|
|
|
2011-08-02 15:44:16 +05:30
|
|
|
open(unit,status='old',err=100,file=trim(getSolverWorkingDirectoryName())//&
|
2010-05-11 12:27:15 +05:30
|
|
|
trim(getSolverJobName())//'.'//newExt)
|
2009-07-22 21:37:19 +05:30
|
|
|
IO_open_jobFile = .true.
|
2010-07-13 15:56:07 +05:30
|
|
|
|
2011-10-18 14:52:33 +05:30
|
|
|
100 endfunction
|
2009-07-22 21:37:19 +05:30
|
|
|
|
|
|
|
|
2011-08-02 15:44:16 +05:30
|
|
|
!********************************************************************
|
|
|
|
! open (write) file related to current job
|
|
|
|
! but with different extension to given unit
|
|
|
|
!********************************************************************
|
|
|
|
logical function IO_write_jobFile(unit,newExt)
|
|
|
|
|
|
|
|
use prec, only: pReal, pInt
|
|
|
|
use DAMASK_interface
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
integer(pInt), intent(in) :: unit
|
|
|
|
character(*), intent(in) :: newExt
|
|
|
|
|
|
|
|
IO_write_jobFile = .false.
|
|
|
|
|
|
|
|
open(unit,status='replace',err=100,file=trim(getSolverWorkingDirectoryName())//&
|
|
|
|
trim(getSolverJobName())//'.'//newExt)
|
|
|
|
IO_write_jobFile = .true.
|
|
|
|
|
2011-10-18 14:52:33 +05:30
|
|
|
100 endfunction
|
2011-08-02 15:44:16 +05:30
|
|
|
|
|
|
|
|
2010-11-03 20:09:18 +05:30
|
|
|
!********************************************************************
|
|
|
|
! open (write) binary file related to current job
|
|
|
|
! but with different extension to given unit
|
|
|
|
!********************************************************************
|
|
|
|
logical function IO_write_jobBinaryFile(unit,newExt,recMultiplier)
|
|
|
|
|
|
|
|
use prec, only: pReal, pInt
|
2011-05-11 22:31:03 +05:30
|
|
|
use DAMASK_interface
|
2010-11-03 20:09:18 +05:30
|
|
|
implicit none
|
|
|
|
|
|
|
|
integer(pInt), intent(in) :: unit
|
|
|
|
integer(pInt), intent(in), optional :: recMultiplier
|
|
|
|
character(*), intent(in) :: newExt
|
|
|
|
|
|
|
|
IO_write_jobBinaryFile = .false.
|
|
|
|
if (present(recMultiplier)) then
|
|
|
|
open(unit,status='replace',form='unformatted',access='direct',recl=pReal*recMultiplier, &
|
|
|
|
err=100,file=trim(getSolverWorkingDirectoryName())//&
|
|
|
|
trim(getSolverJobName())//'.'//newExt)
|
|
|
|
else
|
|
|
|
open(unit,status='replace',form='unformatted',access='direct',recl=pReal, &
|
|
|
|
err=100,file=trim(getSolverWorkingDirectoryName())//&
|
|
|
|
trim(getSolverJobName())//'.'//newExt)
|
|
|
|
endif
|
|
|
|
IO_write_jobBinaryFile = .true.
|
|
|
|
|
2011-10-18 14:52:33 +05:30
|
|
|
100 endfunction
|
2010-11-03 20:09:18 +05:30
|
|
|
|
|
|
|
|
|
|
|
!********************************************************************
|
|
|
|
! open (read) binary file related to restored job
|
|
|
|
! and with different extension to given unit
|
|
|
|
!********************************************************************
|
|
|
|
logical function IO_read_jobBinaryFile(unit,newExt,jobName,recMultiplier)
|
|
|
|
|
|
|
|
use prec, only: pReal, pInt
|
2011-05-11 22:31:03 +05:30
|
|
|
use DAMASK_interface
|
2010-11-03 20:09:18 +05:30
|
|
|
implicit none
|
|
|
|
|
|
|
|
integer(pInt), intent(in) :: unit
|
|
|
|
integer(pInt), intent(in), optional :: recMultiplier
|
|
|
|
character(*), intent(in) :: newExt, jobName
|
|
|
|
|
|
|
|
IO_read_jobBinaryFile = .false.
|
|
|
|
if (present(recMultiplier)) then
|
|
|
|
open(unit,status='old',form='unformatted',access='direct',recl=pReal*recMultiplier, &
|
|
|
|
err=100,file=trim(getSolverWorkingDirectoryName())//&
|
|
|
|
trim(jobName)//'.'//newExt)
|
|
|
|
else
|
|
|
|
open(unit,status='old',form='unformatted',access='direct',recl=pReal, &
|
|
|
|
err=100,file=trim(getSolverWorkingDirectoryName())//&
|
|
|
|
trim(jobName)//'.'//newExt)
|
|
|
|
endif
|
|
|
|
IO_read_jobBinaryFile = .true.
|
|
|
|
|
2011-10-18 14:52:33 +05:30
|
|
|
100 endfunction
|
2010-11-03 20:09:18 +05:30
|
|
|
|
|
|
|
|
2007-03-21 18:02:15 +05:30
|
|
|
!********************************************************************
|
|
|
|
! hybrid IA repetition counter
|
|
|
|
!********************************************************************
|
2009-06-15 18:41:21 +05:30
|
|
|
function hybridIA_reps(dV_V,steps,C)
|
2007-03-21 18:02:15 +05:30
|
|
|
|
|
|
|
use prec, only: pReal, pInt
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
integer(pInt), intent(in), dimension(3) :: steps
|
|
|
|
integer(pInt) hybridIA_reps, phi1,Phi,phi2
|
|
|
|
real(pReal), intent(in), dimension(steps(3),steps(2),steps(1)) :: dV_V
|
|
|
|
real(pReal), intent(in) :: C
|
|
|
|
|
|
|
|
hybridIA_reps = 0_pInt
|
|
|
|
do phi1=1,steps(1)
|
|
|
|
do Phi =1,steps(2)
|
|
|
|
do phi2=1,steps(3)
|
|
|
|
hybridIA_reps = hybridIA_reps+nint(C*dV_V(phi2,Phi,phi1), pInt)
|
2009-06-15 18:41:21 +05:30
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
enddo
|
2007-03-21 18:02:15 +05:30
|
|
|
|
2009-06-15 18:41:21 +05:30
|
|
|
endfunction
|
2007-03-21 18:02:15 +05:30
|
|
|
|
|
|
|
|
|
|
|
!********************************************************************
|
|
|
|
! hybrid IA sampling of ODFfile
|
|
|
|
!********************************************************************
|
2009-06-15 18:41:21 +05:30
|
|
|
function IO_hybridIA(Nast,ODFfileName)
|
2007-03-21 18:02:15 +05:30
|
|
|
|
|
|
|
use prec, only: pReal, pInt
|
|
|
|
implicit none
|
|
|
|
|
2010-07-13 15:56:07 +05:30
|
|
|
real(pReal), parameter :: pi = 3.14159265358979323846264338327950288419716939937510_pReal
|
|
|
|
real(pReal), parameter :: inRad = pi/180.0_pReal
|
|
|
|
|
2007-03-21 18:02:15 +05:30
|
|
|
character(len=*) ODFfileName
|
|
|
|
character(len=80) line
|
|
|
|
character(len=*), parameter :: fileFormat = '(A80)'
|
2007-03-26 20:33:21 +05:30
|
|
|
integer(pInt) i,j,bin,Nast,NnonZero,Nset,Nreps,reps,phi1,Phi,phi2
|
|
|
|
integer(pInt), dimension(7) :: pos
|
2007-03-21 18:02:15 +05:30
|
|
|
integer(pInt), dimension(3) :: steps
|
|
|
|
integer(pInt), dimension(:), allocatable :: binSet
|
|
|
|
real(pReal) center,sum_dV_V,prob,dg_0,C,lowerC,upperC,rnd
|
|
|
|
real(pReal), dimension(3) :: limits,deltas
|
|
|
|
real(pReal), dimension(:,:,:), allocatable :: dV_V
|
2007-03-26 20:33:21 +05:30
|
|
|
real(pReal), dimension(3,Nast) :: IO_hybridIA
|
2007-03-21 18:02:15 +05:30
|
|
|
|
|
|
|
if (.not. IO_open_file(999,ODFfileName)) goto 100
|
|
|
|
|
|
|
|
!--- parse header of ODF file ---
|
|
|
|
!--- limits in phi1, Phi, phi2 ---
|
|
|
|
read(999,fmt=fileFormat,end=100) line
|
|
|
|
pos = IO_stringPos(line,3)
|
|
|
|
if (pos(1).ne.3) goto 100
|
|
|
|
do i=1,3
|
|
|
|
limits(i) = IO_intValue(line,pos,i)*inRad
|
2009-06-15 18:41:21 +05:30
|
|
|
enddo
|
2007-03-21 18:02:15 +05:30
|
|
|
|
|
|
|
!--- deltas in phi1, Phi, phi2 ---
|
|
|
|
read(999,fmt=fileFormat,end=100) line
|
|
|
|
pos = IO_stringPos(line,3)
|
|
|
|
if (pos(1).ne.3) goto 100
|
|
|
|
do i=1,3
|
|
|
|
deltas(i) = IO_intValue(line,pos,i)*inRad
|
2009-06-15 18:41:21 +05:30
|
|
|
enddo
|
2007-03-21 18:02:15 +05:30
|
|
|
steps = nint(limits/deltas,pInt)
|
|
|
|
allocate(dV_V(steps(3),steps(2),steps(1)))
|
|
|
|
|
|
|
|
!--- box boundary/center at origin? ---
|
|
|
|
read(999,fmt=fileFormat,end=100) line
|
|
|
|
if (index(IO_lc(line),'bound')>0) then
|
|
|
|
center = 0.5_pReal
|
|
|
|
else
|
|
|
|
center = 0.0_pReal
|
2009-06-15 18:41:21 +05:30
|
|
|
endif
|
2007-03-21 18:02:15 +05:30
|
|
|
|
|
|
|
!--- skip blank line ---
|
|
|
|
read(999,fmt=fileFormat,end=100) line
|
|
|
|
|
|
|
|
sum_dV_V = 0.0_pReal
|
|
|
|
dV_V = 0.0_pReal
|
|
|
|
dg_0 = deltas(1)*deltas(3)*2.0_pReal*sin(deltas(2)/2.0_pReal)
|
|
|
|
NnonZero = 0_pInt
|
|
|
|
|
|
|
|
do phi1=1,steps(1)
|
|
|
|
do Phi=1,steps(2)
|
|
|
|
do phi2=1,steps(3)
|
2008-03-12 19:23:00 +05:30
|
|
|
read(999,fmt=*,end=100) prob
|
2007-03-21 18:02:15 +05:30
|
|
|
if (prob > 0.0_pReal) then
|
|
|
|
NnonZero = NnonZero+1
|
|
|
|
sum_dV_V = sum_dV_V+prob
|
|
|
|
else
|
|
|
|
prob = 0.0_pReal
|
2009-06-15 18:41:21 +05:30
|
|
|
endif
|
2007-03-21 18:02:15 +05:30
|
|
|
dV_V(phi2,Phi,phi1) = prob*dg_0*sin((Phi-1.0_pReal+center)*deltas(2))
|
2009-06-15 18:41:21 +05:30
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
enddo
|
2009-08-13 18:51:22 +05:30
|
|
|
|
2007-03-21 18:02:15 +05:30
|
|
|
dV_V = dV_V/sum_dV_V ! normalize to 1
|
|
|
|
|
|
|
|
!--- now fix bounds ---
|
2009-10-14 18:51:03 +05:30
|
|
|
Nset = max(Nast,NnonZero) ! if less than non-zero voxel count requested, sample at least that much
|
2007-03-21 18:02:15 +05:30
|
|
|
lowerC = 0.0_pReal
|
|
|
|
upperC = real(Nset, pReal)
|
|
|
|
|
|
|
|
do while (hybridIA_reps(dV_V,steps,upperC) < Nset)
|
|
|
|
lowerC = upperC
|
|
|
|
upperC = upperC*2.0_pReal
|
2009-06-15 18:41:21 +05:30
|
|
|
enddo
|
2007-03-21 18:02:15 +05:30
|
|
|
!--- binary search for best C ---
|
|
|
|
do
|
|
|
|
C = (upperC+lowerC)/2.0_pReal
|
|
|
|
Nreps = hybridIA_reps(dV_V,steps,C)
|
|
|
|
if (abs(upperC-lowerC) < upperC*1.0e-14_pReal) then
|
|
|
|
C = upperC
|
|
|
|
Nreps = hybridIA_reps(dV_V,steps,C)
|
|
|
|
exit
|
|
|
|
elseif (Nreps < Nset) then
|
|
|
|
lowerC = C
|
|
|
|
elseif (Nreps > Nset) then
|
|
|
|
upperC = C
|
|
|
|
else
|
|
|
|
exit
|
2009-06-15 18:41:21 +05:30
|
|
|
endif
|
|
|
|
enddo
|
2009-08-13 18:51:22 +05:30
|
|
|
|
2007-03-21 18:02:15 +05:30
|
|
|
allocate(binSet(Nreps))
|
|
|
|
bin = 0 ! bin counter
|
|
|
|
i = 1 ! set counter
|
|
|
|
do phi1=1,steps(1)
|
|
|
|
do Phi=1,steps(2)
|
|
|
|
do phi2=1,steps(3)
|
|
|
|
reps = nint(C*dV_V(phi2,Phi,phi1), pInt)
|
|
|
|
binSet(i:i+reps-1) = bin
|
|
|
|
bin = bin+1 ! advance bin
|
|
|
|
i = i+reps ! advance set
|
2009-06-15 18:41:21 +05:30
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
enddo
|
2007-03-21 18:02:15 +05:30
|
|
|
|
|
|
|
do i=1,Nast
|
|
|
|
if (i < Nast) then
|
|
|
|
call random_number(rnd)
|
2009-08-13 18:51:22 +05:30
|
|
|
j = nint(rnd*(Nreps-i)+i+0.5_pReal,pInt)
|
2007-03-21 18:02:15 +05:30
|
|
|
else
|
|
|
|
j = i
|
2009-06-15 18:41:21 +05:30
|
|
|
endif
|
2007-03-21 18:02:15 +05:30
|
|
|
bin = binSet(j)
|
|
|
|
IO_hybridIA(1,i) = deltas(1)*(mod(bin/(steps(3)*steps(2)),steps(1))+center) ! phi1
|
|
|
|
IO_hybridIA(2,i) = deltas(2)*(mod(bin/ steps(3) ,steps(2))+center) ! Phi
|
|
|
|
IO_hybridIA(3,i) = deltas(3)*(mod(bin ,steps(3))+center) ! phi2
|
|
|
|
binSet(j) = binSet(i)
|
2009-06-15 18:41:21 +05:30
|
|
|
enddo
|
2007-03-21 18:02:15 +05:30
|
|
|
close(999)
|
|
|
|
return
|
|
|
|
|
|
|
|
! on error
|
2007-03-26 20:33:21 +05:30
|
|
|
100 IO_hybridIA = -1
|
2007-03-21 18:02:15 +05:30
|
|
|
close(999)
|
|
|
|
|
2009-06-15 18:41:21 +05:30
|
|
|
endfunction
|
2007-03-21 18:02:15 +05:30
|
|
|
|
|
|
|
|
2009-03-04 17:18:54 +05:30
|
|
|
!********************************************************************
|
|
|
|
! identifies lines without content
|
|
|
|
!********************************************************************
|
2009-06-15 18:41:21 +05:30
|
|
|
pure function IO_isBlank (line)
|
2009-03-04 17:18:54 +05:30
|
|
|
|
|
|
|
use prec, only: pInt
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
character(len=*), intent(in) :: line
|
|
|
|
character(len=*), parameter :: blank = achar(32)//achar(9)//achar(10)//achar(13) ! whitespaces
|
|
|
|
character(len=*), parameter :: comment = achar(35) ! comment id '#'
|
|
|
|
integer(pInt) posNonBlank, posComment
|
|
|
|
logical IO_isBlank
|
|
|
|
|
|
|
|
posNonBlank = verify(line,blank)
|
|
|
|
posComment = scan(line,comment)
|
|
|
|
IO_isBlank = posNonBlank == 0 .or. posNonBlank == posComment
|
|
|
|
|
2009-06-15 18:41:21 +05:30
|
|
|
endfunction
|
2009-03-04 17:18:54 +05:30
|
|
|
|
|
|
|
!********************************************************************
|
|
|
|
! get tagged content of line
|
|
|
|
!********************************************************************
|
2011-09-13 21:24:06 +05:30
|
|
|
pure function IO_getTag (line,openChar,closeChar)
|
2009-03-04 17:18:54 +05:30
|
|
|
|
|
|
|
use prec, only: pInt
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
character(len=*), intent(in) :: line,openChar,closeChar
|
|
|
|
character(len=*), parameter :: sep=achar(32)//achar(9)//achar(10)//achar(13) ! whitespaces
|
|
|
|
character(len=len_trim(line)) IO_getTag
|
|
|
|
integer(pInt) left,right
|
|
|
|
|
|
|
|
IO_getTag = ''
|
|
|
|
left = scan(line,openChar)
|
|
|
|
right = scan(line,closeChar)
|
|
|
|
|
|
|
|
if (left == verify(line,sep) .and. right > left) & ! openChar is first and closeChar occurs
|
|
|
|
IO_getTag = line(left+1:right-1)
|
|
|
|
|
2009-06-15 18:41:21 +05:30
|
|
|
endfunction
|
2009-03-04 17:18:54 +05:30
|
|
|
|
|
|
|
|
|
|
|
!*********************************************************************
|
2009-06-15 18:41:21 +05:30
|
|
|
function IO_countSections(file,part)
|
2009-03-04 17:18:54 +05:30
|
|
|
!*********************************************************************
|
|
|
|
use prec, only: pInt
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
!* Definition of variables
|
|
|
|
integer(pInt), intent(in) :: file
|
|
|
|
character(len=*), intent(in) :: part
|
|
|
|
integer(pInt) IO_countSections
|
|
|
|
character(len=1024) line
|
|
|
|
|
|
|
|
IO_countSections = 0
|
|
|
|
line = ''
|
|
|
|
rewind(file)
|
|
|
|
|
|
|
|
do while (IO_getTag(line,'<','>') /= part) ! search for part
|
|
|
|
read(file,'(a1024)',END=100) line
|
|
|
|
enddo
|
|
|
|
|
|
|
|
do
|
|
|
|
read(file,'(a1024)',END=100) line
|
|
|
|
if (IO_isBlank(line)) cycle ! skip empty lines
|
|
|
|
if (IO_getTag(line,'<','>') /= '') exit ! stop at next part
|
|
|
|
if (IO_getTag(line,'[',']') /= '') & ! found [section] identifier
|
|
|
|
IO_countSections = IO_countSections + 1
|
|
|
|
enddo
|
|
|
|
|
2011-10-18 14:52:33 +05:30
|
|
|
100 endfunction
|
2009-03-04 17:18:54 +05:30
|
|
|
|
|
|
|
|
|
|
|
!*********************************************************************
|
|
|
|
! return array of myTag counts within <part> for at most N[sections]
|
|
|
|
!*********************************************************************
|
2009-06-15 18:41:21 +05:30
|
|
|
function IO_countTagInPart(file,part,myTag,Nsections)
|
2009-03-04 17:18:54 +05:30
|
|
|
|
|
|
|
use prec, only: pInt
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
!* Definition of variables
|
|
|
|
integer(pInt), intent(in) :: file, Nsections
|
|
|
|
character(len=*), intent(in) :: part, myTag
|
|
|
|
integer(pInt), dimension(Nsections) :: IO_countTagInPart, counter
|
|
|
|
integer(pInt), parameter :: maxNchunks = 1
|
|
|
|
integer(pInt), dimension(1+2*maxNchunks) :: positions
|
|
|
|
integer(pInt) section
|
|
|
|
character(len=1024) line,tag
|
|
|
|
|
|
|
|
counter = 0_pInt
|
|
|
|
section = 0_pInt
|
|
|
|
line = ''
|
|
|
|
rewind(file)
|
|
|
|
|
|
|
|
do while (IO_getTag(line,'<','>') /= part) ! search for part
|
|
|
|
read(file,'(a1024)',END=100) line
|
|
|
|
enddo
|
|
|
|
|
|
|
|
do
|
|
|
|
read(file,'(a1024)',END=100) line
|
|
|
|
if (IO_isBlank(line)) cycle ! skip empty lines
|
|
|
|
if (IO_getTag(line,'<','>') /= '') exit ! stop at next part
|
|
|
|
if (IO_getTag(line,'[',']') /= '') & ! found [section] identifier
|
|
|
|
section = section + 1
|
|
|
|
if (section > 0) then
|
|
|
|
positions = IO_stringPos(line,maxNchunks)
|
|
|
|
tag = IO_lc(IO_stringValue(line,positions,1)) ! extract key
|
|
|
|
if (tag == myTag) & ! match
|
|
|
|
counter(section) = counter(section) + 1
|
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
|
|
|
|
100 IO_countTagInPart = counter
|
|
|
|
|
2009-06-15 18:41:21 +05:30
|
|
|
endfunction
|
2009-03-04 17:18:54 +05:30
|
|
|
|
|
|
|
|
2009-04-03 16:00:18 +05:30
|
|
|
!*********************************************************************
|
|
|
|
! return array of myTag presence within <part> for at most N[sections]
|
|
|
|
!*********************************************************************
|
2009-06-15 18:41:21 +05:30
|
|
|
function IO_spotTagInPart(file,part,myTag,Nsections)
|
2009-04-03 16:00:18 +05:30
|
|
|
|
|
|
|
use prec, only: pInt
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
!* Definition of variables
|
|
|
|
integer(pInt), intent(in) :: file, Nsections
|
|
|
|
character(len=*), intent(in) :: part, myTag
|
|
|
|
logical, dimension(Nsections) :: IO_spotTagInPart
|
|
|
|
integer(pInt), parameter :: maxNchunks = 1
|
|
|
|
integer(pInt), dimension(1+2*maxNchunks) :: positions
|
|
|
|
integer(pInt) section
|
|
|
|
character(len=1024) line,tag
|
|
|
|
|
|
|
|
IO_spotTagInPart = .false. ! assume to nowhere spot tag
|
|
|
|
section = 0_pInt
|
|
|
|
line = ''
|
|
|
|
rewind(file)
|
|
|
|
|
|
|
|
do while (IO_getTag(line,'<','>') /= part) ! search for part
|
|
|
|
read(file,'(a1024)',END=100) line
|
|
|
|
enddo
|
|
|
|
|
|
|
|
do
|
|
|
|
read(file,'(a1024)',END=100) line
|
|
|
|
if (IO_isBlank(line)) cycle ! skip empty lines
|
|
|
|
if (IO_getTag(line,'<','>') /= '') exit ! stop at next part
|
|
|
|
if (IO_getTag(line,'[',']') /= '') & ! found [section] identifier
|
|
|
|
section = section + 1
|
|
|
|
if (section > 0) then
|
|
|
|
positions = IO_stringPos(line,maxNchunks)
|
|
|
|
tag = IO_lc(IO_stringValue(line,positions,1)) ! extract key
|
|
|
|
if (tag == myTag) & ! match
|
|
|
|
IO_spotTagInPart(section) = .true.
|
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
|
2011-10-18 14:52:33 +05:30
|
|
|
100 endfunction
|
2009-04-03 16:00:18 +05:30
|
|
|
|
2007-03-20 19:25:22 +05:30
|
|
|
!********************************************************************
|
|
|
|
! locate at most N space-separated parts in line
|
2009-12-15 21:33:53 +05:30
|
|
|
! return array containing number of parts in line and
|
|
|
|
! the left/right positions of at most N to be used by IO_xxxVal
|
2007-03-20 19:25:22 +05:30
|
|
|
!********************************************************************
|
2009-12-15 21:33:53 +05:30
|
|
|
! pure function IO_stringPos (line,N)
|
|
|
|
function IO_stringPos (line,N)
|
2007-03-20 19:25:22 +05:30
|
|
|
|
|
|
|
use prec, only: pReal,pInt
|
|
|
|
implicit none
|
|
|
|
|
2009-01-16 20:55:37 +05:30
|
|
|
character(len=*), intent(in) :: line
|
2009-10-12 21:31:49 +05:30
|
|
|
character(len=*), parameter :: sep=achar(44)//achar(32)//achar(9)//achar(10)//achar(13) ! comma and whitespaces
|
2009-01-20 00:40:58 +05:30
|
|
|
integer(pInt), intent(in) :: N
|
2009-12-15 21:33:53 +05:30
|
|
|
integer(pInt) left,right
|
2007-03-20 19:25:22 +05:30
|
|
|
integer(pInt) IO_stringPos(1+N*2)
|
2007-04-25 20:08:22 +05:30
|
|
|
|
2007-03-20 19:25:22 +05:30
|
|
|
IO_stringPos = -1
|
|
|
|
IO_stringPos(1) = 0
|
2009-12-15 21:33:53 +05:30
|
|
|
right = 0
|
|
|
|
|
|
|
|
do while (verify(line(right+1:),sep)>0)
|
|
|
|
left = right + verify(line(right+1:),sep)
|
|
|
|
right = left + scan(line(left:),sep) - 2
|
2011-05-30 14:39:19 +05:30
|
|
|
if ( line(left:left) == '#' ) then
|
|
|
|
exit
|
|
|
|
endif
|
2009-12-15 21:33:53 +05:30
|
|
|
if ( IO_stringPos(1)<N ) then
|
|
|
|
IO_stringPos(1+IO_stringPos(1)*2+1) = left
|
|
|
|
IO_stringPos(1+IO_stringPos(1)*2+2) = right
|
|
|
|
endif
|
|
|
|
IO_stringPos(1) = IO_stringPos(1)+1
|
2009-06-15 18:41:21 +05:30
|
|
|
enddo
|
2009-12-15 21:33:53 +05:30
|
|
|
|
2011-10-18 14:52:33 +05:30
|
|
|
endfunction
|
2007-03-20 19:25:22 +05:30
|
|
|
|
|
|
|
|
|
|
|
!********************************************************************
|
|
|
|
! read string value at pos from line
|
|
|
|
!********************************************************************
|
2009-06-15 18:41:21 +05:30
|
|
|
pure function IO_stringValue (line,positions,pos)
|
2007-03-20 19:25:22 +05:30
|
|
|
|
|
|
|
use prec, only: pReal,pInt
|
|
|
|
implicit none
|
|
|
|
|
2009-01-16 20:55:37 +05:30
|
|
|
character(len=*), intent(in) :: line
|
|
|
|
integer(pInt), intent(in) :: positions(*),pos
|
2007-03-20 19:25:22 +05:30
|
|
|
character(len=1+positions(pos*2+1)-positions(pos*2)) IO_stringValue
|
2007-04-25 20:08:22 +05:30
|
|
|
|
|
|
|
if (positions(1) < pos) then
|
|
|
|
IO_stringValue = ''
|
2007-03-28 15:30:49 +05:30
|
|
|
else
|
2007-04-25 20:08:22 +05:30
|
|
|
IO_stringValue = line(positions(pos*2):positions(pos*2+1))
|
2007-03-28 15:30:49 +05:30
|
|
|
endif
|
2007-03-20 19:25:22 +05:30
|
|
|
|
2009-06-15 18:41:21 +05:30
|
|
|
endfunction
|
2007-03-20 19:25:22 +05:30
|
|
|
|
|
|
|
|
2007-03-21 20:15:03 +05:30
|
|
|
!********************************************************************
|
2007-03-21 20:19:21 +05:30
|
|
|
! read string value at pos from fixed format line
|
2007-03-21 20:15:03 +05:30
|
|
|
!********************************************************************
|
2009-06-15 18:41:21 +05:30
|
|
|
pure function IO_fixedStringValue (line,ends,pos)
|
2007-03-21 20:15:03 +05:30
|
|
|
|
|
|
|
use prec, only: pReal,pInt
|
|
|
|
implicit none
|
|
|
|
|
2009-01-16 20:55:37 +05:30
|
|
|
character(len=*), intent(in) :: line
|
|
|
|
integer(pInt), intent(in) :: ends(*),pos
|
2007-03-28 15:30:49 +05:30
|
|
|
character(len=ends(pos+1)-ends(pos)) IO_fixedStringValue
|
2007-03-21 20:15:03 +05:30
|
|
|
|
2007-03-28 15:30:49 +05:30
|
|
|
IO_fixedStringValue = line(ends(pos)+1:ends(pos+1))
|
2007-03-21 20:15:03 +05:30
|
|
|
|
2009-06-15 18:41:21 +05:30
|
|
|
endfunction
|
2007-03-21 20:15:03 +05:30
|
|
|
|
|
|
|
|
2007-03-20 19:25:22 +05:30
|
|
|
!********************************************************************
|
|
|
|
! read float value at pos from line
|
|
|
|
!********************************************************************
|
2011-10-18 14:52:33 +05:30
|
|
|
pure function IO_floatValue (line,positions,myPos)
|
2007-03-20 19:25:22 +05:30
|
|
|
|
|
|
|
use prec, only: pReal,pInt
|
|
|
|
implicit none
|
|
|
|
|
2009-01-16 20:55:37 +05:30
|
|
|
character(len=*), intent(in) :: line
|
2011-10-18 14:52:33 +05:30
|
|
|
integer(pInt), intent(in) :: positions(*),myPos
|
2009-01-20 00:40:58 +05:30
|
|
|
real(pReal) IO_floatValue
|
2007-03-20 19:25:22 +05:30
|
|
|
|
2011-10-18 14:52:33 +05:30
|
|
|
if (positions(1) < myPos) then
|
2009-06-19 12:39:39 +05:30
|
|
|
IO_floatValue = 0.0_pReal
|
|
|
|
else
|
2011-10-18 14:52:33 +05:30
|
|
|
read(UNIT=line(positions(myPos*2):positions(myPos*2+1)),ERR=100,FMT=*) IO_floatValue
|
2007-03-28 15:30:49 +05:30
|
|
|
endif
|
2009-06-19 12:39:39 +05:30
|
|
|
return
|
2007-04-26 18:10:06 +05:30
|
|
|
100 IO_floatValue = huge(1.0_pReal)
|
2007-03-20 19:25:22 +05:30
|
|
|
|
2009-06-15 18:41:21 +05:30
|
|
|
endfunction
|
2007-03-20 19:25:22 +05:30
|
|
|
|
|
|
|
|
2007-03-21 20:15:03 +05:30
|
|
|
!********************************************************************
|
2007-03-21 20:19:21 +05:30
|
|
|
! read float value at pos from fixed format line
|
2007-03-21 20:15:03 +05:30
|
|
|
!********************************************************************
|
2009-06-15 18:41:21 +05:30
|
|
|
pure function IO_fixedFloatValue (line,ends,pos)
|
2007-03-21 20:15:03 +05:30
|
|
|
|
|
|
|
use prec, only: pReal,pInt
|
|
|
|
implicit none
|
|
|
|
|
2009-01-16 20:55:37 +05:30
|
|
|
character(len=*), intent(in) :: line
|
|
|
|
integer(pInt), intent(in) :: ends(*),pos
|
2009-01-20 00:40:58 +05:30
|
|
|
real(pReal) IO_fixedFloatValue
|
2007-03-21 20:15:03 +05:30
|
|
|
|
2007-04-26 18:10:06 +05:30
|
|
|
read(UNIT=line(ends(pos-1)+1:ends(pos)),ERR=100,FMT=*) IO_fixedFloatValue
|
2007-03-21 20:15:03 +05:30
|
|
|
return
|
2007-04-26 18:10:06 +05:30
|
|
|
100 IO_fixedFloatValue = huge(1.0_pReal)
|
2007-03-21 20:15:03 +05:30
|
|
|
|
2009-06-15 18:41:21 +05:30
|
|
|
endfunction
|
2007-03-21 20:15:03 +05:30
|
|
|
|
|
|
|
|
|
|
|
!********************************************************************
|
2007-03-21 20:19:21 +05:30
|
|
|
! read float x.y+z value at pos from format line line
|
2007-03-21 20:15:03 +05:30
|
|
|
!********************************************************************
|
2009-06-15 18:41:21 +05:30
|
|
|
pure function IO_fixedNoEFloatValue (line,ends,pos)
|
2007-03-21 20:15:03 +05:30
|
|
|
|
|
|
|
use prec, only: pReal,pInt
|
|
|
|
implicit none
|
|
|
|
|
2009-01-16 20:55:37 +05:30
|
|
|
character(len=*), intent(in) :: line
|
2009-01-20 00:40:58 +05:30
|
|
|
integer(pInt), intent(in) :: ends(*),pos
|
|
|
|
integer(pInt) pos_exp,expon
|
|
|
|
real(pReal) IO_fixedNoEFloatValue,base
|
2007-03-21 20:15:03 +05:30
|
|
|
|
2007-03-28 15:30:49 +05:30
|
|
|
pos_exp = scan(line(ends(pos)+1:ends(pos+1)),'+-',back=.true.)
|
2007-03-21 20:15:03 +05:30
|
|
|
if (pos_exp > 1) then
|
2008-03-12 19:23:00 +05:30
|
|
|
read(UNIT=line(ends(pos)+1:ends(pos)+pos_exp-1),ERR=100,FMT=*) base
|
|
|
|
read(UNIT=line(ends(pos)+pos_exp:ends(pos+1)),ERR=100,FMT=*) expon
|
2007-03-21 20:15:03 +05:30
|
|
|
else
|
2007-04-26 18:10:06 +05:30
|
|
|
read(UNIT=line(ends(pos)+1:ends(pos+1)),ERR=100,FMT=*) base
|
2007-03-21 22:34:10 +05:30
|
|
|
expon = 0_pInt
|
2007-03-21 20:15:03 +05:30
|
|
|
endif
|
|
|
|
IO_fixedNoEFloatValue = base*10.0_pReal**expon
|
|
|
|
return
|
2007-04-26 18:10:06 +05:30
|
|
|
100 IO_fixedNoEFloatValue = huge(1.0_pReal)
|
2007-03-21 20:15:03 +05:30
|
|
|
|
2009-06-15 18:41:21 +05:30
|
|
|
endfunction
|
2007-03-21 20:15:03 +05:30
|
|
|
|
|
|
|
|
2007-03-20 19:25:22 +05:30
|
|
|
!********************************************************************
|
|
|
|
! read int value at pos from line
|
|
|
|
!********************************************************************
|
2009-06-15 18:41:21 +05:30
|
|
|
pure function IO_intValue (line,positions,pos)
|
2007-03-20 19:25:22 +05:30
|
|
|
|
|
|
|
use prec, only: pReal,pInt
|
|
|
|
implicit none
|
|
|
|
|
2009-01-16 20:55:37 +05:30
|
|
|
character(len=*), intent(in) :: line
|
|
|
|
integer(pInt), intent(in) :: positions(*),pos
|
2009-01-20 00:40:58 +05:30
|
|
|
integer(pInt) IO_intValue
|
2007-03-20 19:25:22 +05:30
|
|
|
|
2009-06-19 12:39:39 +05:30
|
|
|
if (positions(1) < pos) then
|
|
|
|
IO_intValue = 0_pInt
|
|
|
|
else
|
2008-03-12 19:23:00 +05:30
|
|
|
read(UNIT=line(positions(pos*2):positions(pos*2+1)),ERR=100,FMT=*) IO_intValue
|
2007-04-25 20:08:22 +05:30
|
|
|
endif
|
2009-06-19 12:39:39 +05:30
|
|
|
return
|
2007-04-26 18:10:06 +05:30
|
|
|
100 IO_intValue = huge(1_pInt)
|
2007-03-20 19:25:22 +05:30
|
|
|
|
2009-06-15 18:41:21 +05:30
|
|
|
endfunction
|
2007-03-20 19:25:22 +05:30
|
|
|
|
|
|
|
|
2007-03-21 20:15:03 +05:30
|
|
|
!********************************************************************
|
|
|
|
! read int value at pos from fixed format line
|
|
|
|
!********************************************************************
|
2009-06-15 18:41:21 +05:30
|
|
|
pure function IO_fixedIntValue (line,ends,pos)
|
2007-03-21 20:15:03 +05:30
|
|
|
|
|
|
|
use prec, only: pReal,pInt
|
|
|
|
implicit none
|
|
|
|
|
2009-01-16 20:55:37 +05:30
|
|
|
character(len=*), intent(in) :: line
|
|
|
|
integer(pInt), intent(in) :: ends(*),pos
|
2009-01-20 00:40:58 +05:30
|
|
|
integer(pInt) IO_fixedIntValue
|
2007-03-21 20:15:03 +05:30
|
|
|
|
2008-03-12 19:23:00 +05:30
|
|
|
read(UNIT=line(ends(pos)+1:ends(pos+1)),ERR=100,FMT=*) IO_fixedIntValue
|
2007-03-21 20:15:03 +05:30
|
|
|
return
|
2007-04-26 18:10:06 +05:30
|
|
|
100 IO_fixedIntValue = huge(1_pInt)
|
2007-03-21 20:15:03 +05:30
|
|
|
|
2009-06-15 18:41:21 +05:30
|
|
|
endfunction
|
2007-03-21 20:15:03 +05:30
|
|
|
|
|
|
|
|
2007-04-25 20:08:22 +05:30
|
|
|
!********************************************************************
|
|
|
|
! change character in line to lower case
|
|
|
|
!********************************************************************
|
2009-06-15 18:41:21 +05:30
|
|
|
pure function IO_lc (line)
|
2007-03-20 19:25:22 +05:30
|
|
|
|
|
|
|
use prec, only: pInt
|
|
|
|
implicit none
|
|
|
|
|
2009-01-16 20:55:37 +05:30
|
|
|
character (len=*), intent(in) :: line
|
2007-03-21 18:02:15 +05:30
|
|
|
character (len=len(line)) IO_lc
|
2007-03-20 19:25:22 +05:30
|
|
|
integer(pInt) i
|
|
|
|
|
2007-03-21 18:02:15 +05:30
|
|
|
IO_lc = line
|
2007-04-25 20:08:22 +05:30
|
|
|
do i=1,len(line)
|
|
|
|
if(64<iachar(line(i:i)) .and. iachar(line(i:i))<91) IO_lc(i:i)=achar(iachar(line(i:i))+32)
|
|
|
|
enddo
|
2007-03-20 19:25:22 +05:30
|
|
|
|
2009-06-15 18:41:21 +05:30
|
|
|
endfunction
|
2007-03-20 19:25:22 +05:30
|
|
|
|
|
|
|
|
|
|
|
!********************************************************************
|
2007-03-21 18:02:15 +05:30
|
|
|
! in place change of character in line to lower case
|
2007-03-20 19:25:22 +05:30
|
|
|
!********************************************************************
|
2009-06-15 18:41:21 +05:30
|
|
|
subroutine IO_lcInplace (line)
|
2007-03-20 19:25:22 +05:30
|
|
|
|
|
|
|
use prec, only: pInt
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
character (len=*) line
|
2007-03-28 14:01:12 +05:30
|
|
|
character (len=len(line)) IO_lc
|
2007-03-20 19:25:22 +05:30
|
|
|
integer(pInt) i
|
|
|
|
|
2007-03-28 14:01:12 +05:30
|
|
|
IO_lc = line
|
2007-04-25 20:08:22 +05:30
|
|
|
do i=1,len(line)
|
|
|
|
if(64<iachar(line(i:i)) .and. iachar(line(i:i))<91) IO_lc(i:i)=achar(iachar(line(i:i))+32)
|
|
|
|
enddo
|
2007-03-28 14:01:12 +05:30
|
|
|
line = IO_lc
|
2007-03-20 19:25:22 +05:30
|
|
|
|
2009-06-15 18:41:21 +05:30
|
|
|
endsubroutine
|
2007-03-20 19:25:22 +05:30
|
|
|
|
|
|
|
|
2009-04-03 12:34:31 +05:30
|
|
|
!********************************************************************
|
|
|
|
! read on in file to skip (at least) N chunks (may be over multiple lines)
|
|
|
|
!********************************************************************
|
2009-06-15 18:41:21 +05:30
|
|
|
subroutine IO_skipChunks (unit,N)
|
2009-04-03 12:34:31 +05:30
|
|
|
|
|
|
|
use prec, only: pReal,pInt
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
integer(pInt) remainingChunks,unit,N
|
|
|
|
integer(pInt), parameter :: maxNchunks = 64
|
|
|
|
integer(pInt), dimension(1+2*maxNchunks) :: pos
|
|
|
|
character(len=300) line
|
|
|
|
|
|
|
|
remainingChunks = N
|
|
|
|
do while (remainingChunks > 0)
|
|
|
|
read(unit,'(A300)',end=100) line
|
|
|
|
pos = IO_stringPos(line,maxNchunks)
|
|
|
|
remainingChunks = remainingChunks - pos(1)
|
2009-06-15 18:41:21 +05:30
|
|
|
enddo
|
2011-10-18 14:52:33 +05:30
|
|
|
100 endsubroutine
|
2009-04-03 12:34:31 +05:30
|
|
|
|
2009-10-12 21:31:49 +05:30
|
|
|
|
|
|
|
!********************************************************************
|
|
|
|
! extract value from key=value pair and check whether key matches
|
|
|
|
!********************************************************************
|
|
|
|
pure function IO_extractValue (line,key)
|
2009-04-03 12:34:31 +05:30
|
|
|
|
2009-10-12 21:31:49 +05:30
|
|
|
use prec, only: pReal,pInt
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
character(len=*), intent(in) :: line,key
|
|
|
|
character(len=*), parameter :: sep = achar(61) ! '='
|
|
|
|
integer(pInt) pos
|
|
|
|
character(len=300) IO_extractValue
|
|
|
|
|
|
|
|
IO_extractValue = ''
|
|
|
|
|
|
|
|
pos = scan(line,sep)
|
|
|
|
if (pos > 0 .and. line(:pos-1) == key(:pos-1)) & ! key matches expected key
|
|
|
|
IO_extractValue = line(pos+1:) ! extract value
|
|
|
|
|
|
|
|
endfunction
|
|
|
|
|
|
|
|
|
2007-04-25 20:08:22 +05:30
|
|
|
!********************************************************************
|
2009-10-12 21:31:49 +05:30
|
|
|
! count lines containig data up to next *keyword
|
2010-07-13 15:56:07 +05:30
|
|
|
! AP: changed the function to neglect comment lines between keyword definitions.
|
|
|
|
! : is not changed back to the original version since *.inp_assembly does not
|
|
|
|
! : contain any comment lines (12.07.2010)
|
2007-10-15 19:25:52 +05:30
|
|
|
!********************************************************************
|
2009-10-12 21:31:49 +05:30
|
|
|
function IO_countDataLines (unit)
|
2007-10-15 19:25:52 +05:30
|
|
|
|
|
|
|
use prec, only: pReal,pInt
|
|
|
|
implicit none
|
|
|
|
|
2009-10-12 21:31:49 +05:30
|
|
|
integer(pInt) IO_countDataLines,unit
|
2010-07-13 15:56:07 +05:30
|
|
|
integer(pInt), parameter :: maxNchunks = 1
|
2009-10-12 21:31:49 +05:30
|
|
|
integer(pInt), dimension(1+2*maxNchunks) :: pos
|
|
|
|
character(len=300) line,tmp
|
2007-10-15 19:25:52 +05:30
|
|
|
|
2009-10-12 21:31:49 +05:30
|
|
|
IO_countDataLines = 0
|
2007-10-15 19:25:52 +05:30
|
|
|
do
|
|
|
|
read(unit,'(A300)',end=100) line
|
2009-10-12 21:31:49 +05:30
|
|
|
pos = IO_stringPos(line,maxNchunks)
|
|
|
|
tmp = IO_lc(IO_stringValue(line,pos,1))
|
2010-07-13 15:56:07 +05:30
|
|
|
if (tmp(1:1) == '*' .and. tmp(2:2) /= '*') then ! found keyword
|
2009-10-12 21:31:49 +05:30
|
|
|
exit
|
2007-10-15 19:25:52 +05:30
|
|
|
else
|
2010-07-13 15:56:07 +05:30
|
|
|
if (tmp(2:2) /= '*') IO_countDataLines = IO_countDataLines + 1_pInt
|
2007-10-15 19:25:52 +05:30
|
|
|
endif
|
|
|
|
enddo
|
2009-10-12 21:31:49 +05:30
|
|
|
100 backspace(unit)
|
|
|
|
|
|
|
|
endfunction
|
|
|
|
|
|
|
|
|
|
|
|
!********************************************************************
|
|
|
|
! count items in consecutive lines
|
|
|
|
! Marc: ints concatenated by "c" as last char or range of values a "to" b
|
|
|
|
! Abaqus: triplet of start,stop,inc
|
|
|
|
!********************************************************************
|
|
|
|
function IO_countContinousIntValues (unit)
|
|
|
|
|
2011-05-11 22:31:03 +05:30
|
|
|
use DAMASK_interface
|
2009-10-12 21:31:49 +05:30
|
|
|
use prec, only: pReal,pInt
|
|
|
|
implicit none
|
|
|
|
|
2010-02-18 13:59:57 +05:30
|
|
|
integer(pInt) unit,l,count
|
2009-10-12 21:31:49 +05:30
|
|
|
integer(pInt) IO_countContinousIntValues
|
2012-01-12 22:31:24 +05:30
|
|
|
integer(pInt), parameter :: maxNchunks = 8192
|
2009-10-12 21:31:49 +05:30
|
|
|
integer(pInt), dimension(1+2*maxNchunks) :: pos
|
2012-01-12 22:31:24 +05:30
|
|
|
character(len=65536) line
|
2009-10-12 21:31:49 +05:30
|
|
|
|
|
|
|
IO_countContinousIntValues = 0_pInt
|
|
|
|
|
|
|
|
select case (FEsolver)
|
2012-01-12 22:31:24 +05:30
|
|
|
case ('Marc','Spectral')
|
2009-10-12 21:31:49 +05:30
|
|
|
|
|
|
|
do
|
|
|
|
read(unit,'(A300)',end=100) line
|
|
|
|
pos = IO_stringPos(line,maxNchunks)
|
|
|
|
if (IO_lc(IO_stringValue(line,pos,2)) == 'to' ) then ! found range indicator
|
|
|
|
IO_countContinousIntValues = 1 + IO_intValue(line,pos,3) - IO_intValue(line,pos,1)
|
|
|
|
exit ! only one single range indicator allowed
|
|
|
|
else
|
|
|
|
IO_countContinousIntValues = IO_countContinousIntValues+pos(1)-1 ! add line's count when assuming 'c'
|
|
|
|
if ( IO_lc(IO_stringValue(line,pos,pos(1))) /= 'c' ) then ! line finished, read last value
|
|
|
|
IO_countContinousIntValues = IO_countContinousIntValues+1
|
|
|
|
exit ! data ended
|
|
|
|
endif
|
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
|
|
|
|
case('Abaqus')
|
|
|
|
|
|
|
|
count = IO_countDataLines(unit)
|
|
|
|
do l = 1,count
|
|
|
|
backspace(unit)
|
|
|
|
enddo
|
|
|
|
|
|
|
|
do l = 1,count
|
|
|
|
read(unit,'(A300)',end=100) line
|
|
|
|
pos = IO_stringPos(line,maxNchunks)
|
|
|
|
IO_countContinousIntValues = IO_countContinousIntValues + 1 + & ! assuming range generation
|
|
|
|
(IO_intValue(line,pos,2)-IO_intValue(line,pos,1))/max(1,IO_intValue(line,pos,3))
|
|
|
|
enddo
|
|
|
|
|
|
|
|
endselect
|
|
|
|
|
2011-10-18 14:52:33 +05:30
|
|
|
100 endfunction
|
2007-10-15 19:25:52 +05:30
|
|
|
|
2009-10-12 21:31:49 +05:30
|
|
|
|
|
|
|
!********************************************************************
|
|
|
|
! return integer list corrsponding to items in consecutive lines
|
|
|
|
! Marc: ints concatenated by "c" as last char, range of a "to" b, or named set
|
|
|
|
! Abaqus: triplet of start,stop,inc or named set
|
|
|
|
!********************************************************************
|
2009-06-15 18:41:21 +05:30
|
|
|
function IO_continousIntValues (unit,maxN,lookupName,lookupMap,lookupMaxN)
|
2007-04-25 20:08:22 +05:30
|
|
|
|
2011-05-11 22:31:03 +05:30
|
|
|
use DAMASK_interface
|
2007-04-25 20:08:22 +05:30
|
|
|
use prec, only: pReal,pInt
|
|
|
|
implicit none
|
|
|
|
|
2009-10-12 21:31:49 +05:30
|
|
|
integer(pInt) unit,maxN,i,j,l,count,first,last
|
2007-04-25 20:08:22 +05:30
|
|
|
integer(pInt), dimension(1+maxN) :: IO_continousIntValues
|
2012-01-12 22:31:24 +05:30
|
|
|
integer(pInt), parameter :: maxNchunks = 8192_pInt
|
2009-10-12 21:31:49 +05:30
|
|
|
integer(pInt), dimension(1+2*maxNchunks) :: pos
|
2007-10-23 18:38:27 +05:30
|
|
|
character(len=64), dimension(:) :: lookupName
|
|
|
|
integer(pInt) :: lookupMaxN
|
|
|
|
integer(pInt), dimension(:,:) :: lookupMap
|
2012-01-12 22:31:24 +05:30
|
|
|
character(len=65536) line
|
2010-07-13 15:56:07 +05:30
|
|
|
logical rangeGeneration
|
2007-04-25 20:08:22 +05:30
|
|
|
|
2009-10-12 21:31:49 +05:30
|
|
|
IO_continousIntValues = 0
|
2010-07-13 15:56:07 +05:30
|
|
|
rangeGeneration = .false.
|
2009-10-12 21:31:49 +05:30
|
|
|
|
|
|
|
select case (FEsolver)
|
2012-01-12 22:31:24 +05:30
|
|
|
case ('Marc','Spectral')
|
2009-10-12 21:31:49 +05:30
|
|
|
|
|
|
|
do
|
2012-01-12 22:31:24 +05:30
|
|
|
read(unit,'(A65536)',end=100) line
|
2009-10-12 21:31:49 +05:30
|
|
|
pos = IO_stringPos(line,maxNchunks)
|
|
|
|
if (verify(IO_stringValue(line,pos,1),"0123456789") > 0) then ! a non-int, i.e. set name
|
|
|
|
do i = 1,lookupMaxN ! loop over known set names
|
|
|
|
if (IO_stringValue(line,pos,1) == lookupName(i)) then ! found matching name
|
|
|
|
IO_continousIntValues = lookupMap(:,i) ! return resp. entity list
|
|
|
|
exit
|
|
|
|
endif
|
|
|
|
enddo
|
2007-10-23 18:38:27 +05:30
|
|
|
exit
|
2012-01-12 22:31:24 +05:30
|
|
|
else if (pos(1) > 2_pInt .and. IO_lc(IO_stringValue(line,pos,2)) == 'to' ) then ! found range indicator
|
2009-10-12 21:31:49 +05:30
|
|
|
do i = IO_intValue(line,pos,1),IO_intValue(line,pos,3)
|
|
|
|
IO_continousIntValues(1) = IO_continousIntValues(1) + 1
|
|
|
|
IO_continousIntValues(1+IO_continousIntValues(1)) = i
|
|
|
|
enddo
|
|
|
|
exit
|
|
|
|
else
|
|
|
|
do i = 1,pos(1)-1 ! interpret up to second to last value
|
|
|
|
IO_continousIntValues(1) = IO_continousIntValues(1) + 1
|
|
|
|
IO_continousIntValues(1+IO_continousIntValues(1)) = IO_intValue(line,pos,i)
|
|
|
|
enddo
|
|
|
|
if ( IO_lc(IO_stringValue(line,pos,pos(1))) /= 'c' ) then ! line finished, read last value
|
|
|
|
IO_continousIntValues(1) = IO_continousIntValues(1)+1
|
|
|
|
IO_continousIntValues(1+IO_continousIntValues(1)) = IO_intValue(line,pos,pos(1))
|
|
|
|
exit
|
|
|
|
endif
|
2007-10-23 18:38:27 +05:30
|
|
|
endif
|
|
|
|
enddo
|
2009-10-12 21:31:49 +05:30
|
|
|
|
|
|
|
case('Abaqus')
|
|
|
|
|
|
|
|
count = IO_countDataLines(unit)
|
|
|
|
do l = 1,count
|
|
|
|
backspace(unit)
|
2007-04-25 20:08:22 +05:30
|
|
|
enddo
|
2009-10-12 21:31:49 +05:30
|
|
|
|
2010-07-13 15:56:07 +05:30
|
|
|
! check if the element values in the elset are auto generated
|
|
|
|
backspace(unit)
|
2012-01-12 22:31:24 +05:30
|
|
|
read(unit,'(A65536)',end=100) line
|
2010-07-13 15:56:07 +05:30
|
|
|
pos = IO_stringPos(line,maxNchunks)
|
|
|
|
do i = 1,pos(1)
|
|
|
|
if (IO_lc(IO_stringValue(line,pos,i)) == 'generate') rangeGeneration = .true.
|
|
|
|
enddo
|
|
|
|
|
2009-10-12 21:31:49 +05:30
|
|
|
do l = 1,count
|
2012-01-12 22:31:24 +05:30
|
|
|
read(unit,'(A65536)',end=100) line
|
2009-10-12 21:31:49 +05:30
|
|
|
pos = IO_stringPos(line,maxNchunks)
|
|
|
|
if (verify(IO_stringValue(line,pos,1),"0123456789") > 0) then ! a non-int, i.e. set names follow on this line
|
|
|
|
do i = 1,pos(1) ! loop over set names in line
|
|
|
|
do j = 1,lookupMaxN ! look thru known set names
|
|
|
|
if (IO_stringValue(line,pos,i) == lookupName(j)) then ! found matching name
|
|
|
|
first = 2 + IO_continousIntValues(1) ! where to start appending data
|
|
|
|
last = first + lookupMap(1,j) - 1 ! up to where to append data
|
|
|
|
IO_continousIntValues(first:last) = lookupMap(2:1+lookupMap(1,j),j) ! add resp. entity list
|
|
|
|
IO_continousIntValues(1) = IO_continousIntValues(1) + lookupMap(1,j) ! count them
|
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
enddo
|
2010-07-13 15:56:07 +05:30
|
|
|
else if (rangeGeneration) then ! range generation
|
2009-10-12 21:31:49 +05:30
|
|
|
do i = IO_intValue(line,pos,1),IO_intValue(line,pos,2),max(1,IO_intValue(line,pos,3))
|
|
|
|
IO_continousIntValues(1) = IO_continousIntValues(1) + 1
|
|
|
|
IO_continousIntValues(1+IO_continousIntValues(1)) = i
|
|
|
|
enddo
|
2010-07-13 15:56:07 +05:30
|
|
|
else ! read individual elem nums
|
|
|
|
do i = 1,pos(1)
|
|
|
|
! write(*,*)'IO_CIV-int',IO_intValue(line,pos,i)
|
|
|
|
IO_continousIntValues(1) = IO_continousIntValues(1) + 1
|
|
|
|
IO_continousIntValues(1+IO_continousIntValues(1)) = IO_intValue(line,pos,i)
|
|
|
|
enddo
|
2009-10-12 21:31:49 +05:30
|
|
|
endif
|
2007-04-25 20:08:22 +05:30
|
|
|
enddo
|
2009-10-12 21:31:49 +05:30
|
|
|
|
|
|
|
endselect
|
|
|
|
|
2011-10-18 14:52:33 +05:30
|
|
|
100 endfunction
|
2007-04-25 20:08:22 +05:30
|
|
|
|
|
|
|
|
2009-10-12 21:31:49 +05:30
|
|
|
|
2007-03-20 19:25:22 +05:30
|
|
|
!********************************************************************
|
|
|
|
! write error statements to standard out
|
|
|
|
! and terminate the Marc run with exit #9xxx
|
|
|
|
! in ABAQUS either time step is reduced or execution terminated
|
|
|
|
!********************************************************************
|
2011-11-02 20:08:42 +05:30
|
|
|
subroutine IO_error(error_ID,e,i,g,ext_msg)
|
2007-03-20 19:25:22 +05:30
|
|
|
|
2008-02-15 15:34:49 +05:30
|
|
|
use prec, only: pInt
|
2007-03-20 19:25:22 +05:30
|
|
|
implicit none
|
|
|
|
|
2011-11-02 20:08:42 +05:30
|
|
|
integer(pInt), intent(in) :: error_ID
|
2009-03-04 17:18:54 +05:30
|
|
|
integer(pInt), optional, intent(in) :: e,i,g
|
|
|
|
character(len=*), optional, intent(in) :: ext_msg
|
2011-11-02 20:08:42 +05:30
|
|
|
character(len=1024) msg
|
2007-03-20 19:25:22 +05:30
|
|
|
|
2011-11-02 20:08:42 +05:30
|
|
|
select case (error_ID)
|
2011-07-14 15:07:31 +05:30
|
|
|
case (30)
|
2011-10-18 14:52:33 +05:30
|
|
|
msg = 'could not open spectral loadcase'
|
2011-07-14 15:07:31 +05:30
|
|
|
case (31)
|
|
|
|
msg = 'mask consistency violated in spectral loadcase'
|
|
|
|
case (32)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'ill-defined L (each line should be either fully or not at all defined) in spectral loadcase'
|
2011-07-14 15:07:31 +05:30
|
|
|
case (34)
|
|
|
|
msg = 'negative time increment in spectral loadcase'
|
|
|
|
case (35)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'non-positive increments in spectral loadcase'
|
2011-07-14 15:07:31 +05:30
|
|
|
case (36)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'non-positive result frequency in spectral loadcase'
|
2011-08-26 19:36:37 +05:30
|
|
|
case (37)
|
|
|
|
msg = 'incomplete loadcase'
|
2011-09-02 19:20:05 +05:30
|
|
|
case (38)
|
|
|
|
msg = 'mixed boundary conditions allow rotation'
|
2011-11-07 16:34:57 +05:30
|
|
|
case (39)
|
|
|
|
msg = 'non-positive restart frequency in spectral loadcase'
|
2010-06-10 20:21:10 +05:30
|
|
|
case (40)
|
|
|
|
msg = 'path rectification error'
|
|
|
|
case (41)
|
|
|
|
msg = 'path too long'
|
|
|
|
case (42)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'missing header info in spectral mesh'
|
2010-06-10 20:21:10 +05:30
|
|
|
case (43)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'resolution in spectral mesh'
|
2010-06-10 20:21:10 +05:30
|
|
|
case (44)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'dimension in spectral mesh'
|
2011-10-18 14:52:33 +05:30
|
|
|
case (45)
|
|
|
|
msg = 'incomplete information in spectral mesh header'
|
2011-10-24 23:56:34 +05:30
|
|
|
case (46)
|
|
|
|
msg = 'not a rotation defined for loadcase rotation'
|
2012-01-25 14:30:40 +05:30
|
|
|
case (47)
|
|
|
|
msg = 'updating of gamma operator not possible if it is pre calculated'
|
2009-07-22 21:37:19 +05:30
|
|
|
case (50)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'writing constitutive output description'
|
2009-03-04 17:18:54 +05:30
|
|
|
case (100)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'opening material configuration'
|
2010-02-18 13:59:57 +05:30
|
|
|
case (101)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'opening input file'
|
2012-01-13 21:48:16 +05:30
|
|
|
case (102)
|
|
|
|
msg = 'precistion not suitable for FFTW'
|
2011-02-21 20:07:38 +05:30
|
|
|
case (103)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'odd resolution given'
|
2011-02-21 20:07:38 +05:30
|
|
|
case (104)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'initializing FFTW'
|
2009-03-04 17:18:54 +05:30
|
|
|
case (105)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'reading from ODF file'
|
2011-05-28 15:14:43 +05:30
|
|
|
case (106)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'reading info on old job'
|
2011-10-18 22:12:06 +05:30
|
|
|
case (107)
|
|
|
|
msg = 'writing spectralOut file'
|
2008-02-15 15:34:49 +05:30
|
|
|
case (110)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'no homogenization specified via State Variable 2'
|
2008-02-15 15:34:49 +05:30
|
|
|
case (120)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'no microstructure specified via State Variable 3'
|
2010-02-25 23:09:11 +05:30
|
|
|
case (125)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'no entries in config part'
|
2009-03-04 17:18:54 +05:30
|
|
|
case (130)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'homogenization index out of bounds'
|
2009-03-04 17:18:54 +05:30
|
|
|
case (140)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'microstructure index out of bounds'
|
2009-03-04 17:18:54 +05:30
|
|
|
case (150)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'crystallite index out of bounds'
|
2010-02-25 23:09:11 +05:30
|
|
|
case (155)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'phase index out of bounds'
|
2009-03-04 17:18:54 +05:30
|
|
|
case (160)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'texture index out of bounds'
|
2009-03-04 17:18:54 +05:30
|
|
|
case (170)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'sum of phase fractions differs from 1'
|
2012-01-12 22:31:24 +05:30
|
|
|
case (180)
|
|
|
|
msg = 'mismatch of microstructure count and a*b*c in geom file'
|
2007-03-20 19:25:22 +05:30
|
|
|
case (200)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'unknown constitution specified'
|
2009-03-04 17:18:54 +05:30
|
|
|
case (201)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'unknown homogenization specified'
|
2009-03-04 17:18:54 +05:30
|
|
|
case (205)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'unknown lattice structure encountered'
|
2009-07-22 21:37:19 +05:30
|
|
|
case (210)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'negative initial resistance'
|
2009-07-22 21:37:19 +05:30
|
|
|
case (211)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'non-positive reference shear rate'
|
2009-07-22 21:37:19 +05:30
|
|
|
case (212)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'non-positive stress exponent'
|
2009-07-22 21:37:19 +05:30
|
|
|
case (213)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'non-positive saturation stress'
|
2009-08-31 19:43:10 +05:30
|
|
|
case (214)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'zero hardening exponent'
|
2009-03-05 21:36:01 +05:30
|
|
|
case (220)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'negative initial dislocation density'
|
2009-03-05 21:36:01 +05:30
|
|
|
case (221)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'negative Burgers vector'
|
2009-03-05 21:36:01 +05:30
|
|
|
case (222)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'negative activation energy for edge dislocation glide'
|
2009-03-05 21:36:01 +05:30
|
|
|
case (223)
|
2011-09-26 15:25:08 +05:30
|
|
|
msg = 'zero stackin fault energy'
|
|
|
|
! case (224)
|
|
|
|
! msg = 'non-positive diffusion prefactor'
|
2009-07-22 21:37:19 +05:30
|
|
|
case (225)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'no slip systems specified'
|
2009-09-18 21:07:14 +05:30
|
|
|
case (226)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'non-positive prefactor for dislocation velocity'
|
2009-09-18 21:07:14 +05:30
|
|
|
case (227)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'non-positive prefactor for mean free path'
|
2009-09-18 21:07:14 +05:30
|
|
|
case (228)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'non-positive minimum stable dipole distance'
|
2009-09-18 21:07:14 +05:30
|
|
|
case (229)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'non-positive hardening interaction coefficients'
|
2009-09-18 21:07:14 +05:30
|
|
|
case (230)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'non-positive atomic volume'
|
2009-09-18 21:07:14 +05:30
|
|
|
case (231)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'non-positive prefactor for self-diffusion coefficient' ! what is the difference to 224 ??
|
2009-09-18 21:07:14 +05:30
|
|
|
case (232)
|
2011-09-26 15:25:08 +05:30
|
|
|
msg = 'non-positive activation energy for self-diffusion'
|
2009-09-18 21:07:14 +05:30
|
|
|
case (233)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'non-positive relevant dislocation density'
|
2011-09-26 15:25:08 +05:30
|
|
|
case (234)
|
|
|
|
msg = 'error in shear banding input'
|
2011-01-12 18:06:48 +05:30
|
|
|
case (235)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'material parameter in nonlocal constitutive phase out of bounds:'
|
2011-01-26 15:47:42 +05:30
|
|
|
case (236)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'unknown material parameter in nonlocal constitutive phase:'
|
2011-03-22 19:28:42 +05:30
|
|
|
case (237)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'singularity in internal stress calculation'
|
2009-03-05 22:37:32 +05:30
|
|
|
case (240)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'non-positive Taylor factor'
|
2009-07-22 21:37:19 +05:30
|
|
|
case (241)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'non-positive hardening exponent'
|
2009-09-18 21:07:14 +05:30
|
|
|
case (242)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'non-positive relevant slip resistance'
|
2009-06-15 18:41:21 +05:30
|
|
|
case (260)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'non-positive relevant strain'
|
2009-06-15 18:41:21 +05:30
|
|
|
case (261)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'frequency of stiffness update < 0'
|
2009-06-15 18:41:21 +05:30
|
|
|
case (262)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'frequency of Jacobian update of Lp residuum < 0'
|
2009-06-15 18:41:21 +05:30
|
|
|
case (263)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'non-positive perturbation value'
|
2009-06-15 18:41:21 +05:30
|
|
|
case (264)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'limit for homogenization loop too small'
|
2009-06-15 18:41:21 +05:30
|
|
|
case (265)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'limit for crystallite loop too small'
|
2009-06-15 18:41:21 +05:30
|
|
|
case (266)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'limit for state loop too small'
|
2009-06-15 18:41:21 +05:30
|
|
|
case (267)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'limit for stress loop too small'
|
2009-06-15 18:41:21 +05:30
|
|
|
case (268)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'non-positive minimum substep size'
|
2009-06-15 18:41:21 +05:30
|
|
|
case (269)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'non-positive relative state tolerance'
|
2009-06-15 18:41:21 +05:30
|
|
|
case (270)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'Non-positive relative stress tolerance'
|
2009-06-15 18:41:21 +05:30
|
|
|
case (271)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'Non-positive absolute stress tolerance'
|
2009-07-31 17:32:20 +05:30
|
|
|
|
|
|
|
!* Error messages related to RGC numerical parameters <<<updated 31.07.2009>>>
|
2009-06-15 18:41:21 +05:30
|
|
|
case (272)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'non-positive relative tolerance of residual in RGC'
|
2009-06-15 18:41:21 +05:30
|
|
|
case (273)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'non-positive absolute tolerance of residual in RGC'
|
2009-06-15 18:41:21 +05:30
|
|
|
case (274)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'non-positive relative maximum of residual in RGC'
|
2009-06-15 18:41:21 +05:30
|
|
|
case (275)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'non-positive absolute maximum of residual in RGC'
|
2009-07-31 17:32:20 +05:30
|
|
|
case (276)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'non-positive penalty perturbation in RGC'
|
2009-07-31 17:32:20 +05:30
|
|
|
case (277)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'non-positive relevant mismatch in RGC'
|
2009-11-17 19:12:38 +05:30
|
|
|
case (278)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'non-positive definite viscosity model in RGC'
|
2009-11-17 19:12:38 +05:30
|
|
|
case (288)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'non-positive maximum threshold of relaxation change in RGC'
|
2009-12-16 21:50:53 +05:30
|
|
|
case (289)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'non-positive definite volume discrepancy penalty in RGC'
|
2009-07-31 17:32:20 +05:30
|
|
|
|
2010-05-20 20:25:11 +05:30
|
|
|
case (294)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'non-positive tolerance for deformation gradient'
|
2010-05-20 20:25:11 +05:30
|
|
|
|
2010-10-01 17:48:49 +05:30
|
|
|
case (298)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'chosen integration method does not exist'
|
2009-11-10 19:06:27 +05:30
|
|
|
case (299)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'chosen perturbation method does not exist'
|
2009-11-10 19:06:27 +05:30
|
|
|
|
2007-03-20 19:25:22 +05:30
|
|
|
case (300)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'this material can only be used with elements with three direct stress components'
|
2007-03-20 19:25:22 +05:30
|
|
|
case (500)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'unknown lattice type specified'
|
2010-05-04 18:24:13 +05:30
|
|
|
case (550)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'unknown symmetry type specified'
|
2007-03-20 19:25:22 +05:30
|
|
|
case (600)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'convergence not reached'
|
2009-03-04 17:18:54 +05:30
|
|
|
case (610)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'stress loop not converged'
|
2009-10-16 01:32:52 +05:30
|
|
|
|
|
|
|
case (666)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'memory leak detected'
|
2011-05-04 21:32:18 +05:30
|
|
|
case (667)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'invalid materialpoint result requested'
|
2009-10-16 01:32:52 +05:30
|
|
|
|
2010-05-06 19:37:21 +05:30
|
|
|
case (670)
|
|
|
|
msg = 'math_check: quat -> axisAngle -> quat failed'
|
|
|
|
case (671)
|
|
|
|
msg = 'math_check: quat -> R -> quat failed'
|
|
|
|
case (672)
|
|
|
|
msg = 'math_check: quat -> euler -> quat failed'
|
|
|
|
case (673)
|
|
|
|
msg = 'math_check: R -> euler -> R failed'
|
|
|
|
|
2007-03-20 19:25:22 +05:30
|
|
|
case (700)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'singular matrix in stress iteration'
|
2010-07-13 15:56:07 +05:30
|
|
|
|
2011-02-21 20:07:38 +05:30
|
|
|
case (800)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'matrix inversion error'
|
2011-02-21 20:07:38 +05:30
|
|
|
|
2010-07-13 15:56:07 +05:30
|
|
|
! Error messages related to parsing of Abaqus input file
|
|
|
|
case (900)
|
|
|
|
msg = 'PARSE ERROR: Improper definition of nodes in input file (Nnodes < 2)'
|
|
|
|
case (901)
|
|
|
|
msg = 'PARSE ERROR: No Elements defined in input file (Nelems = 0)'
|
|
|
|
case (902)
|
|
|
|
msg = 'PARSE ERROR: No Element sets defined in input file (Atleast one *Elset must exist)'
|
|
|
|
case (903)
|
|
|
|
msg = 'PARSE ERROR: No Materials defined in input file (Look into section assigments)'
|
|
|
|
case (904)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'PARSE ERROR: No elements could be assigned for Elset: '
|
2010-07-13 15:56:07 +05:30
|
|
|
case (905)
|
|
|
|
msg = 'PARSE ERROR: Error in mesh_abaqus_map_materials'
|
|
|
|
case (906)
|
|
|
|
msg = 'PARSE ERROR: Error in mesh_abaqus_count_cpElements'
|
|
|
|
case (907)
|
|
|
|
msg = 'PARSE ERROR: Incorrect size of mesh_mapFEtoCPelem in mesh_abaqus_map_elements; Size cannot be zero'
|
|
|
|
case (908)
|
|
|
|
msg = 'PARSE ERROR: Incorrect size of mesh_mapFEtoCPnode in mesh_abaqus_map_nodes; Size cannot be zero'
|
|
|
|
case (909)
|
|
|
|
msg = 'PARSE ERROR: Incorrect size of mesh_node in mesh_abaqus_build_nodes; must be equal to mesh_Nnodes'
|
|
|
|
case(910)
|
2011-08-02 15:44:16 +05:30
|
|
|
msg = 'PARSE ERROR: Incorrect element type mapping in '
|
2010-07-13 15:56:07 +05:30
|
|
|
|
|
|
|
|
2007-03-20 19:25:22 +05:30
|
|
|
case default
|
2009-03-31 14:51:57 +05:30
|
|
|
msg = 'Unknown error number...'
|
2007-03-20 19:25:22 +05:30
|
|
|
end select
|
2009-01-20 00:40:58 +05:30
|
|
|
|
2008-05-26 18:41:25 +05:30
|
|
|
!$OMP CRITICAL (write2out)
|
2009-03-31 14:51:57 +05:30
|
|
|
write(6,*)
|
2011-11-02 20:08:42 +05:30
|
|
|
write(6,'(a38)') '+------------------------------------+'
|
|
|
|
write(6,'(a38)') '+ error +'
|
|
|
|
write(6,'(a17,i3,a18)') '+ ',error_ID,' +'
|
|
|
|
write(6,'(a38)') '+ +'
|
|
|
|
write(6,'(a2,a)') '+ ', trim(msg)
|
|
|
|
if (present(ext_msg)) write(6,'(a2,a)') '+ ', trim(ext_msg)
|
2009-03-04 17:18:54 +05:30
|
|
|
if (present(e)) then
|
|
|
|
if (present(i) .and. present(g)) then
|
2011-05-13 22:25:13 +05:30
|
|
|
write(6,'(a13,i6,a4,i2,a7,i4,a2)') '+ at element ',e,' IP ',i,' grain ',g,' +'
|
2009-03-04 17:18:54 +05:30
|
|
|
else
|
2011-05-13 22:25:13 +05:30
|
|
|
write(6,'(a18,i6,a14)') '+ at ',e,' +'
|
2009-03-04 17:18:54 +05:30
|
|
|
endif
|
|
|
|
endif
|
2009-05-28 22:08:40 +05:30
|
|
|
write(6,'(a38)') '+------------------------------------+'
|
2007-03-20 19:25:22 +05:30
|
|
|
call flush(6)
|
2011-11-02 20:08:42 +05:30
|
|
|
call quit(9000+error_ID)
|
2010-02-18 15:42:45 +05:30
|
|
|
!$OMP END CRITICAL (write2out)
|
2009-01-20 00:40:58 +05:30
|
|
|
|
2007-03-20 19:25:22 +05:30
|
|
|
! ABAQUS returns in some cases
|
|
|
|
|
2009-06-15 18:41:21 +05:30
|
|
|
endsubroutine
|
2007-03-20 19:25:22 +05:30
|
|
|
|
|
|
|
|
2009-03-31 14:51:57 +05:30
|
|
|
!********************************************************************
|
|
|
|
! write warning statements to standard out
|
|
|
|
!********************************************************************
|
2011-11-02 20:08:42 +05:30
|
|
|
subroutine IO_warning(warning_ID,e,i,g,ext_msg)
|
2009-03-31 14:51:57 +05:30
|
|
|
|
|
|
|
use prec, only: pInt
|
|
|
|
implicit none
|
|
|
|
|
2011-11-02 20:08:42 +05:30
|
|
|
integer(pInt), intent(in) :: warning_ID
|
2009-03-31 14:51:57 +05:30
|
|
|
integer(pInt), optional, intent(in) :: e,i,g
|
|
|
|
character(len=*), optional, intent(in) :: ext_msg
|
2011-11-02 20:08:42 +05:30
|
|
|
character(len=1024) msg
|
2009-03-31 14:51:57 +05:30
|
|
|
|
2011-11-02 20:08:42 +05:30
|
|
|
select case (warning_ID)
|
2011-12-06 22:28:17 +05:30
|
|
|
case (34)
|
|
|
|
msg = 'invalid restart increment given'
|
2011-11-15 23:24:18 +05:30
|
|
|
case (47_pInt)
|
|
|
|
msg = 'No valid parameter for FFTW given, using FFTW_PATIENT'
|
2011-11-02 20:08:42 +05:30
|
|
|
case (101_pInt)
|
2010-11-03 20:28:11 +05:30
|
|
|
msg = '+ crystallite debugging off... +'
|
2011-11-02 20:08:42 +05:30
|
|
|
case (600_pInt)
|
2010-11-03 20:28:11 +05:30
|
|
|
msg = '+ crystallite responds elastically +'
|
2011-11-02 20:08:42 +05:30
|
|
|
case (601_pInt)
|
2010-11-04 23:48:01 +05:30
|
|
|
msg = '+ stiffness close to zero +'
|
2011-11-02 20:08:42 +05:30
|
|
|
case (650_pInt)
|
2010-11-03 20:28:11 +05:30
|
|
|
msg = '+ polar decomposition failed +'
|
2011-11-02 20:08:42 +05:30
|
|
|
case (700_pInt)
|
2009-12-15 13:50:31 +05:30
|
|
|
msg = '+ unknown crystal symmetry +'
|
2009-03-31 14:51:57 +05:30
|
|
|
case default
|
2010-11-03 20:28:11 +05:30
|
|
|
msg = '+ unknown warning number... +'
|
2009-03-31 14:51:57 +05:30
|
|
|
end select
|
|
|
|
|
|
|
|
!$OMP CRITICAL (write2out)
|
|
|
|
write(6,*)
|
2011-11-02 20:08:42 +05:30
|
|
|
write(6,'(a38)') '+------------------------------------+'
|
|
|
|
write(6,'(a38)') '+ warning +'
|
|
|
|
write(6,'(a38)') '+ +'
|
|
|
|
write(6,'(a17,i3,a18)') '+ ',warning_ID,' +'
|
|
|
|
write(6,'(a2,a)') '+ ', trim(msg)
|
|
|
|
if (present(ext_msg)) write(6,'(a2,a)') '+ ', trim(ext_msg)
|
2009-03-31 14:51:57 +05:30
|
|
|
if (present(e)) then
|
2010-11-04 23:48:01 +05:30
|
|
|
if (present(i)) then
|
|
|
|
if (present(g)) then
|
|
|
|
write(6,'(a12,x,i6,x,a2,x,i2,x,a5,x,i4,a2)') '+ at element',e,'IP',i,'grain',g,' +'
|
|
|
|
else
|
|
|
|
write(6,'(a12,x,i6,x,a2,x,i2,a13)') '+ at element',e,'IP',i,' +'
|
|
|
|
endif
|
2009-03-31 14:51:57 +05:30
|
|
|
else
|
2010-11-04 23:48:01 +05:30
|
|
|
write(6,'(a12,x,i6,a19)') '+ at element',e,' +'
|
2009-03-31 14:51:57 +05:30
|
|
|
endif
|
|
|
|
endif
|
2009-05-28 22:08:40 +05:30
|
|
|
write(6,'(a38)') '+------------------------------------+'
|
constitutive_nonlocal:
- read in activation energy for dislocation glide from material.config
- changed naming of dDipMin/Max to dLower/dUpper
- added new outputs: rho_dot, rho_dot_dip, rho_dot_gen, rho_dot_sgl2dip, rho_dot_dip2sgl, rho_dot_ann_ath, rho_dot_ann_the, rho_dot_flux, d_upper_edge, d_upper_screw, d_upper_dot_edge, d_upper_dot_screw
- poisson's ratio is now calculated from elastic constants
- microstrucutre has state as first argument, since this is our output variable
- periodic boundary conditions are taken into account for fluxes and internal stresses. for the moment, flag has to be set in constitutive_nonlocal.
- corrected calculation for dipole formation by glide
- added terms for dipole formation/annihilation by stress decrease/increase
constitutive:
- passing of arguments is adapted for constitutive_nonlocal model
crystallite:
- in stiffness calculation: call to collect_dotState used wrong arguments
- crystallite_postResults uses own Tstar_v and temperature, no need for passing them from materialpoint_postResults
homogenization:
- crystallite_postResults uses own Tstar_v and temperature, no need for passing them from materialpoint_postResults
IO:
- changed error message 229
material.config:
- changed example for nonlocal constitution according to constitutive_nonlocal
all:
- added some flush statements
2009-10-20 20:06:03 +05:30
|
|
|
call flush(6)
|
2010-04-06 12:17:15 +05:30
|
|
|
!$OMP END CRITICAL (write2out)
|
2009-03-31 14:51:57 +05:30
|
|
|
|
2009-06-15 18:41:21 +05:30
|
|
|
endsubroutine
|
2009-03-31 14:51:57 +05:30
|
|
|
|
2007-03-20 19:25:22 +05:30
|
|
|
END MODULE IO
|