important bugfix for reading in results in case of restart

This commit is contained in:
Martin Diehl 2011-11-17 22:11:05 +00:00
parent 75e20dffb7
commit dc6c29a910
1 changed files with 427 additions and 461 deletions

View File

@ -53,7 +53,7 @@ program DAMASK_spectral
use mesh, only: mesh_ipCenterOfGravity use mesh, only: mesh_ipCenterOfGravity
use CPFEM, only: CPFEM_general, CPFEM_initAll use CPFEM, only: CPFEM_general, CPFEM_initAll
use FEsolving, only: restartWrite, restartReadSpectral, restartReadStep use FEsolving, only: restartWrite, restartReadSpectral, restartReadStep
use numerics, only: err_div_tol, err_stress_tol, err_stress_tolrel , rotation_tol,& use numerics, only: err_div_tol, err_stress_tolrel , rotation_tol,&
itmax, memory_efficient, DAMASK_NumThreadsInt, divergence_correction, & itmax, memory_efficient, DAMASK_NumThreadsInt, divergence_correction, &
fftw_planner_flag, fftw_timelimit fftw_planner_flag, fftw_timelimit
use homogenization, only: materialpoint_sizeResults, materialpoint_results use homogenization, only: materialpoint_sizeResults, materialpoint_results
@ -129,13 +129,13 @@ program DAMASK_spectral
! loop variables, convergence etc. ! loop variables, convergence etc.
real(pReal) :: time = 0.0_pReal, time0 = 0.0_pReal, timeinc ! elapsed time, begin of interval, time interval real(pReal) :: time = 0.0_pReal, time0 = 0.0_pReal, timeinc ! elapsed time, begin of interval, time interval
real(pReal) :: guessmode, err_div, err_stress, p_hat_avg real(pReal) :: guessmode, err_div, err_stress, err_stress_tol, p_hat_avg
complex(pReal), parameter :: img = cmplx(0.0_pReal,1.0_pReal) complex(pReal), parameter :: img = cmplx(0.0_pReal,1.0_pReal)
real(pReal), dimension(3,3), parameter :: ones = 1.0_pReal, zeroes = 0.0_pReal real(pReal), dimension(3,3), parameter :: ones = 1.0_pReal, zeroes = 0.0_pReal
complex(pReal), dimension(3,3) :: temp33_Complex complex(pReal), dimension(3,3) :: temp33_Complex
real(pReal), dimension(3,3) :: temp33_Real real(pReal), dimension(3,3) :: temp33_Real
integer(pInt) :: i, j, k, l, m, n, p integer(pInt) :: i, j, k, l, m, n, p
integer(pInt) :: N_Loadcases, loadcase, step, iter, ielem, CPFEM_mode, stepZero=1_pInt, & integer(pInt) :: N_Loadcases, loadcase, step, iter, ielem, CPFEM_mode, &
ierr, notConvergedCounter = 0_pInt, totalStepsCounter = 0_pInt ierr, notConvergedCounter = 0_pInt, totalStepsCounter = 0_pInt
logical :: errmatinv, regrid = .false. logical :: errmatinv, regrid = .false.
real(pReal) :: defgradDet, defgradDetMax, defgradDetMin real(pReal) :: defgradDet, defgradDetMax, defgradDetMin
@ -146,16 +146,13 @@ program DAMASK_spectral
real(pReal) :: p_real_avg, err_div_max, err_real_div_avg, err_real_div_max real(pReal) :: p_real_avg, err_div_max, err_real_div_avg, err_real_div_max
logical :: debugGeneral = .false., debugDivergence = .false., debugRestart = .false. logical :: debugGeneral = .false., debugDivergence = .false., debugRestart = .false.
!Initializing ! Initializing model size independed parameters
!$ call omp_set_num_threads(DAMASK_NumThreadsInt) ! set number of threads for parallel execution set by DAMASK_NUM_THREADS !$ call omp_set_num_threads(DAMASK_NumThreadsInt) ! set number of threads for parallel execution set by DAMASK_NUM_THREADS
if (.not.(command_argument_count()==4 .or. command_argument_count()==6)) call IO_error(error_ID=102_pInt) ! check for correct number of given arguments if (.not.(command_argument_count()==4 .or. command_argument_count()==6)) &! check for correct number of given arguments
call IO_error(error_ID=102_pInt)
call DAMASK_interface_init() call DAMASK_interface_init()
if (iand(spectral_debug_verbosity,1_pInt)==1_pInt) debugGeneral = .true.
if (iand(spectral_debug_verbosity,2_pInt)==2_pInt) debugDivergence = .true.
if (iand(spectral_debug_verbosity,4_pInt)==4_pInt) debugRestart = .true.
!$OMP CRITICAL (write2out) !$OMP CRITICAL (write2out)
print '(a)', '' print '(a)', ''
print '(a,a)', ' <<<+- DAMASK_spectral init -+>>>' print '(a,a)', ' <<<+- DAMASK_spectral init -+>>>'
@ -166,7 +163,7 @@ program DAMASK_spectral
print '(a)', '' print '(a)', ''
!$OMP END CRITICAL (write2out) !$OMP END CRITICAL (write2out)
! Reading the loadcase file and allocate variables ! Reading the loadcase file and allocate variables for loadcases
path = getLoadcaseName() path = getLoadcaseName()
if (.not. IO_open_file(myUnit,path)) call IO_error(error_ID=30_pInt,ext_msg = trim(path)) if (.not. IO_open_file(myUnit,path)) call IO_error(error_ID=30_pInt,ext_msg = trim(path))
rewind(myUnit) rewind(myUnit)
@ -194,6 +191,7 @@ program DAMASK_spectral
allocate (bc(N_Loadcases)) allocate (bc(N_Loadcases))
! Reading the loadcase and assign values to the allocated data structure
rewind(myUnit) rewind(myUnit)
loadcase = 0_pInt loadcase = 0_pInt
do do
@ -269,8 +267,8 @@ program DAMASK_spectral
101 close(myUnit) 101 close(myUnit)
!read header of geom file to get the information needed before the complete geom file is intepretated by mesh.f90 !read header of geom file to get the information needed before the complete geom file is intepretated by mesh.f90
path = getModelName() path = getModelName()
if (.not. IO_open_file(myUnit,trim(path)//InputFileExtension))& if (.not. IO_open_file(myUnit,trim(path)//InputFileExtension))&
call IO_error(error_ID=101_pInt,ext_msg = trim(path)//InputFileExtension) call IO_error(error_ID=101_pInt,ext_msg = trim(path)//InputFileExtension)
rewind(myUnit) rewind(myUnit)
@ -326,10 +324,15 @@ program DAMASK_spectral
mod(resolution(2),2_pInt)/=0_pInt .or.& mod(resolution(2),2_pInt)/=0_pInt .or.&
(mod(resolution(3),2_pInt)/=0_pInt .and. resolution(3)/= 1_pInt)) call IO_error(error_ID=103_pInt) (mod(resolution(3),2_pInt)/=0_pInt .and. resolution(3)/= 1_pInt)) call IO_error(error_ID=103_pInt)
! Initialization of CPFEM_general (= constitutive law) and of deformation gradient field ! Initialization of CPFEM_general (= constitutive law)
call CPFEM_initAll(bc(1)%temperature,1_pInt,1_pInt) call CPFEM_initAll(bc(1)%temperature,1_pInt,1_pInt)
!Output of geom file ! Get debugging parameters
if (iand(spectral_debug_verbosity,1_pInt)==1_pInt) debugGeneral = .true.
if (iand(spectral_debug_verbosity,2_pInt)==2_pInt) debugDivergence = .true.
if (iand(spectral_debug_verbosity,4_pInt)==4_pInt) debugRestart = .true.
!Output of geometry
!$OMP CRITICAL (write2out) !$OMP CRITICAL (write2out)
print '(a)', '' print '(a)', ''
print '(a)', '#############################################################' print '(a)', '#############################################################'
@ -350,6 +353,7 @@ program DAMASK_spectral
call IO_warning(warning_ID=33_pInt) ! cannot guess along trajectory for first step of first loadcase call IO_warning(warning_ID=33_pInt) ! cannot guess along trajectory for first step of first loadcase
bc(1)%followFormerTrajectory = .false. bc(1)%followFormerTrajectory = .false.
endif endif
! consistency checks and output of loadcase ! consistency checks and output of loadcase
do loadcase = 1_pInt, N_Loadcases do loadcase = 1_pInt, N_Loadcases
!$OMP CRITICAL (write2out) !$OMP CRITICAL (write2out)
@ -418,9 +422,7 @@ program DAMASK_spectral
call dfftw_plan_with_nthreads(DAMASK_NumThreadsInt) call dfftw_plan_with_nthreads(DAMASK_NumThreadsInt)
endif endif
#endif #endif
!call dfftw_timelimit(fftw_timelimit) ! is not working, have to fix it in FFTW source file !call dfftw_timelimit(fftw_timelimit) ! is not working, have to fix it in FFTW source file
select case(IO_lc(fftw_planner_flag)) ! setting parameters for the plan creation of FFTW. Basically a translation from fftw3.f select case(IO_lc(fftw_planner_flag)) ! setting parameters for the plan creation of FFTW. Basically a translation from fftw3.f
case('estimate','fftw_estimate') ! ordered from slow execution (but fast plan creation) to fast execution case('estimate','fftw_estimate') ! ordered from slow execution (but fast plan creation) to fast execution
fftw_flag = 64 fftw_flag = 64
@ -434,53 +436,11 @@ program DAMASK_spectral
call IO_warning(warning_ID=47_pInt,ext_msg=trim(IO_lc(fftw_planner_flag))) call IO_warning(warning_ID=47_pInt,ext_msg=trim(IO_lc(fftw_planner_flag)))
fftw_flag = 32 fftw_flag = 32
end select end select
if (.not. restartReadSpectral) then ! start at first step of first loadcase
loadcase = 1_pInt
step = 1_pInt
else ! going forwarnd and use old values
do i = 1_pInt, N_Loadcases ! looping over ALL loadcases
time0 = time ! loadcase start time
timeinc = bc(i)%timeIncrement/bc(i)%steps ! only valid for given linear time scale. will be overwritten later in case loglinear scale is used
do j = 1_pInt, bc(i)%steps ! looping over ALL steps in current loadcase
if (totalStepsCounter < restartReadStep) then ! forwarding to restart step
totalStepsCounter = totalStepsCounter + 1_pInt
if (bc(i)%logscale == 1_pInt) then ! loglinear scale
if (i == 1_pInt) then ! 1st loadcase of loglinear scale
if (j == 1_pInt) then ! 1st step of 1st loadcase of loglinear scale
timeinc = bc(1)%timeIncrement*(2.0_pReal**real( 1_pInt-bc(1)%steps ,pReal)) ! assume 1st step is equal to 2nd
else ! not-1st step of 1st loadcase of loglinear scale
timeinc = bc(1)%timeIncrement*(2.0_pReal**real(j-1_pInt-bc(1)%steps ,pReal))
endif
else ! not-1st loadcase of loglinear scale
timeinc = time0 *( (1.0_pReal + bc(i)%timeIncrement/time0 )**real( j/bc(i)%steps ,pReal) &
-(1.0_pReal + bc(i)%timeIncrement/time0 )**real( (j-1_pInt)/bc(i)%steps ,pReal) ) !ToDo: correct? how should the float casting be done
endif
endif
time = time + timeinc
endif
enddo
enddo
do i = 1_pInt, N_Loadcases ! looping over ALL loadcases
do j = 1_pInt, bc(i)%steps ! looping over ALL steps in current loadcase
if (totalStepsCounter -1_pInt < restartReadStep) then ! forwarding to restart step
step = j
loadcase = i
endif
enddo
enddo
endif
print*, totalStepsCounter
print*, loadcase
print*, step
pause
!************************************************************* !*************************************************************
! Loop over loadcases defined in the loadcase file ! Loop over loadcases defined in the loadcase file
do loadcase = loadcase, N_Loadcases do loadcase = 1_pInt, N_Loadcases
!************************************************************* !*************************************************************
time0 = time ! loadcase start time time0 = time ! loadcase start time
if (bc(loadcase)%followFormerTrajectory) then ! continue to guess along former trajectory where applicable if (bc(loadcase)%followFormerTrajectory) then ! continue to guess along former trajectory where applicable
guessmode = 1.0_pReal guessmode = 1.0_pReal
else else
@ -495,16 +455,33 @@ pause
timeinc = bc(loadcase)%timeIncrement/bc(loadcase)%steps ! only valid for given linear time scale. will be overwritten later in case loglinear scale is used timeinc = bc(loadcase)%timeIncrement/bc(loadcase)%steps ! only valid for given linear time scale. will be overwritten later in case loglinear scale is used
fDot = bc(loadcase)%deformation ! only valid for given fDot. will be overwritten later in case L is given fDot = bc(loadcase)%deformation ! only valid for given fDot. will be overwritten later in case L is given
!************************************************************* !*************************************************************
! loop oper steps defined in input file for current loadcase ! loop oper steps defined in input file for current loadcase
do step = step, bc(loadcase)%steps do step = 1_pInt, bc(loadcase)%steps
!************************************************************* !*************************************************************
! forwarding time
if (bc(loadcase)%logscale == 1_pInt) then ! loglinear scale
if (loadcase == 1_pInt) then ! 1st loadcase of loglinear scale
if (step == 1_pInt) then ! 1st step of 1st loadcase of loglinear scale
timeinc = bc(1)%timeIncrement*(2.0_pReal**real( 1_pInt-bc(1)%steps ,pReal)) ! assume 1st step is equal to 2nd
else ! not-1st step of 1st loadcase of loglinear scale
timeinc = bc(1)%timeIncrement*(2.0_pReal**real(step-1_pInt-bc(1)%steps ,pReal))
endif
else ! not-1st loadcase of loglinear scale
timeinc = time0 *( (1.0_pReal + bc(loadcase)%timeIncrement/time0 )**real( step/bc(loadcase)%steps ,pReal) &
-(1.0_pReal + bc(loadcase)%timeIncrement/time0 )**real( (step-1_pInt)/bc(loadcase)%steps ,pReal) )
endif
endif
time = time + timeinc
totalStepsCounter = totalStepsCounter + 1_pInt
!************************************************************* !*************************************************************
! Initialization Start ! Initialization Start
!************************************************************* !*************************************************************
if(stepZero==1_pInt) then ! we start if(totalStepsCounter >= restartReadStep) then ! Do calculations (otherwise just forwarding)
stepZero = 0_pInt if (regrid==.true. ) then ! 'DeInitialize' the values changing in case of regridding
if (regrid==.true. ) then ! 'real' start vs. regrid start regrid = .false.
call dfftw_destroy_plan(fftw_plan(1)); call dfftw_destroy_plan(fftw_plan(2)) call dfftw_destroy_plan(fftw_plan(1)); call dfftw_destroy_plan(fftw_plan(2))
if(debugDivergence) call dfftw_destroy_plan(fftw_plan(3)) if(debugDivergence) call dfftw_destroy_plan(fftw_plan(3))
deallocate (defgrad) deallocate (defgrad)
@ -513,8 +490,10 @@ pause
deallocate (temperature) deallocate (temperature)
deallocate (xi) deallocate (xi)
deallocate (workfft) deallocate (workfft)
! here we have to create the new geometry and assign the values from the previous step !ToDo: here we have to create the new geometry and assign the values from the previous step
endif endif
if(totalStepsCounter == restartReadStep) then ! Initialize values
guessmode = 0.0_pReal ! change of load case, homogeneous guess for the first step
allocate (defgrad ( resolution(1),resolution(2),resolution(3),3,3)); defgrad = 0.0_pReal allocate (defgrad ( resolution(1),resolution(2),resolution(3),3,3)); defgrad = 0.0_pReal
allocate (defgradold ( resolution(1),resolution(2),resolution(3),3,3)); defgradold = 0.0_pReal allocate (defgradold ( resolution(1),resolution(2),resolution(3),3,3)); defgradold = 0.0_pReal
allocate (coordinates(3,resolution(1),resolution(2),resolution(3))); coordinates = 0.0_pReal allocate (coordinates(3,resolution(1),resolution(2),resolution(3))); coordinates = 0.0_pReal
@ -539,14 +518,13 @@ pause
write (6,*) 'FFTW initialized' write (6,*) 'FFTW initialized'
!$OMP END CRITICAL (write2out) !$OMP END CRITICAL (write2out)
endif endif
if (restartReadStep==1_pInt) then ! no deformation at the beginning
if (.not. restartReadSpectral) then ! no deformation at the beginning
do k = 1_pInt, resolution(3); do j = 1_pInt, resolution(2); do i = 1_pInt, resolution(1) do k = 1_pInt, resolution(3); do j = 1_pInt, resolution(2); do i = 1_pInt, resolution(1)
defgrad(i,j,k,1:3,1:3) = math_I3 defgrad(i,j,k,1:3,1:3) = math_I3
defgradold(i,j,k,1:3,1:3) = math_I3 defgradold(i,j,k,1:3,1:3) = math_I3
enddo; enddo; enddo enddo; enddo; enddo
else ! using old values else ! using old values
if (IO_read_jobBinaryFile(777,'convergedSpectralDefgrad',trim(getModelName()),size(defgrad))) then if (IO_read_jobBinaryFile(777,'convergedSpectralDefgrad',trim(getSolverJobName()),size(defgrad))) then
read (777,rec=1) defgrad read (777,rec=1) defgrad
close (777) close (777)
endif endif
@ -557,8 +535,8 @@ pause
enddo; enddo; enddo enddo; enddo; enddo
defgradAim = defgradAim * wgt defgradAim = defgradAim * wgt
defgradAimOld = defgradAim defgradAimOld = defgradAim
guessmode=0.0_pInt
endif endif
ielem = 0_pInt ielem = 0_pInt
do k = 1_pInt, resolution(3); do j = 1_pInt, resolution(2); do i = 1_pInt, resolution(1) do k = 1_pInt, resolution(3); do j = 1_pInt, resolution(2); do i = 1_pInt, resolution(1)
ielem = ielem + 1_pInt ielem = ielem + 1_pInt
@ -627,13 +605,13 @@ pause
write(538), 'dimension', geomdimension write(538), 'dimension', geomdimension
write(538), 'materialpoint_sizeResults', materialpoint_sizeResults write(538), 'materialpoint_sizeResults', materialpoint_sizeResults
write(538), 'loadcases', N_Loadcases write(538), 'loadcases', N_Loadcases
write(538), 'logscale', bc(loadcase)%logscale ! one entry per loadcase (0: linear, 1: log) write(538), 'logscale', bc(1:N_Loadcases)%logscale ! one entry per loadcase (0: linear, 1: log)
write(538), 'frequencies', bc(loadcase)%outputfrequency ! one entry per loadcase write(538), 'frequencies', bc(1:N_Loadcases)%outputfrequency ! one entry per loadcase
write(538), 'times', bc(loadcase)%timeIncrement ! one entry per loadcase write(538), 'times', bc(1:N_Loadcases)%timeIncrement ! one entry per loadcase
bc(1)%steps= bc(1)%steps + 1_pInt bc(1)%steps= bc(1)%steps + 1_pInt
write(538), 'increments', bc(loadcase)%steps ! one entry per loadcase ToDo: rename keyword to steps write(538), 'increments', bc(1:N_Loadcases)%steps ! one entry per loadcase ToDo: rename keyword to steps
bc(1)%steps= bc(1)%steps - 1_pInt bc(1)%steps= bc(1)%steps - 1_pInt
write(538), 'startingIncrement', totalStepsCounter write(538), 'startingIncrement', restartReadStep -1_pInt ! start with writing out the previous step
write(538), 'eoh' ! end of header write(538), 'eoh' ! end of header
write(538), materialpoint_results(:,1,:) ! initial (non-deformed) results !ToDo: define array size write(538), materialpoint_results(:,1,:) ! initial (non-deformed) results !ToDo: define array size
!$OMP END CRITICAL (write2out) !$OMP END CRITICAL (write2out)
@ -641,31 +619,13 @@ pause
!************************************************************* !*************************************************************
! Initialization End ! Initialization End
!************************************************************* !*************************************************************
totalStepsCounter = totalStepsCounter + 1_pInt
if (mod(step - 1_pInt,bc(loadcase)%restartFrequency)==0_pInt) then ! at frequency of writing restart information if (mod(step - 1_pInt,bc(loadcase)%restartFrequency)==0_pInt) then ! at frequency of writing restart information
restartWrite = .true. ! setting restart parameter for FEsolving (first call to CPFEM_general will write ToDo: true?) restartWrite = .true. ! setting restart parameter for FEsolving (first call to CPFEM_general will write ToDo: true?)
if(IO_write_jobBinaryFile(777,'convergedSpectralDefgrad',size(defgrad))) then ! and writing deformation gradient field to file
write (777,rec=1) defgrad
close (777)
endif
else else
restartWrite = .false. restartWrite = .false.
endif endif
if (bc(loadcase)%logscale == 1_pInt) then ! loglinear scale
if (loadcase == 1_pInt) then ! 1st loadcase of loglinear scale
if (step == 1_pInt) then ! 1st step of 1st loadcase of loglinear scale
timeinc = bc(1)%timeIncrement*(2.0_pReal**real( 1_pInt-bc(1)%steps ,pReal)) ! assume 1st step is equal to 2nd
else ! not-1st step of 1st loadcase of loglinear scale
timeinc = bc(1)%timeIncrement*(2.0_pReal**real(step-1_pInt-bc(1)%steps ,pReal))
endif
else ! not-1st loadcase of loglinear scale
timeinc = time0 *( (1.0_pReal + bc(loadcase)%timeIncrement/time0 )**real( step/bc(loadcase)%steps ,pReal) &
-(1.0_pReal + bc(loadcase)%timeIncrement/time0 )**real( (step-1_pInt)/bc(loadcase)%steps ,pReal) )
endif
endif
time = time + timeinc
if (bc(loadcase)%velGradApplied) & ! calculate fDot from given L and current F if (bc(loadcase)%velGradApplied) & ! calculate fDot from given L and current F
fDot = math_mul33x33(bc(loadcase)%deformation, defgradAim) fDot = math_mul33x33(bc(loadcase)%deformation, defgradAim)
@ -737,9 +697,14 @@ pause
!$OMP CRITICAL (write2out) !$OMP CRITICAL (write2out)
print '(a)', '#############################################################' print '(a)', '#############################################################'
print '(A,I5.5,A,es12.6)', 'Increment ', totalStepsCounter, ' Time ',time print '(A,I5.5,A,es12.6)', 'Increment ', totalStepsCounter, ' Time ',time
if (restartWrite .eq. .true. ) print '(A)', 'Writing converged Results of previous Step for Restart' if (restartWrite ) then
print '(A)', 'Writing converged Results of previous Step for Restart'
if(IO_write_jobBinaryFile(777,'convergedSpectralDefgrad',size(defgrad))) then ! and writing deformation gradient field to file
write (777,rec=1) defgrad
close (777)
endif
endif
!$OMP END CRITICAL (write2out) !$OMP END CRITICAL (write2out)
!************************************************************* !*************************************************************
! convergence loop ! convergence loop
do while(iter < itmax .and. & do while(iter < itmax .and. &
@ -917,6 +882,7 @@ pause
write(538), materialpoint_results(:,1,:) ! write result to file write(538), materialpoint_results(:,1,:) ! write result to file
endif endif
!$OMP END CRITICAL (write2out) !$OMP END CRITICAL (write2out)
endif
enddo ! end looping over steps in current loadcase enddo ! end looping over steps in current loadcase
deallocate(c_reduced) deallocate(c_reduced)
deallocate(s_reduced) deallocate(s_reduced)