added shell scripts to run tests from crontab, improved numbering in spectral compile tests and added some first files for testing restart capabilities of spectral solver.

added configuration file for generation of code documentation using doxygen
restart of spectral solver is fixed and seem to work now also for restart at significant deformation.
spectral solver now gives exit code 2 if some increments did not converge
This commit is contained in:
Martin Diehl 2012-06-18 15:27:01 +00:00
parent c36a73c142
commit 8537e87b7e
4 changed files with 106 additions and 85 deletions

View File

@ -35,9 +35,6 @@
! MPI fuer Eisenforschung, Duesseldorf ! MPI fuer Eisenforschung, Duesseldorf
#include "spectral_quit.f90" #include "spectral_quit.f90"
!#ifdef PETSC
!#include "finclude/petscdef.h"
!#endif
program DAMASK_spectral program DAMASK_spectral
use, intrinsic :: iso_fortran_env ! to get compiler_version and compiler_options (at least for gfortran >4.6 at the moment) use, intrinsic :: iso_fortran_env ! to get compiler_version and compiler_options (at least for gfortran >4.6 at the moment)
@ -47,7 +44,8 @@ program DAMASK_spectral
loadCaseFile, & loadCaseFile, &
geometryFile, & geometryFile, &
getSolverWorkingDirectoryName, & getSolverWorkingDirectoryName, &
getSolverJobName getSolverJobName, &
appendToOutFile
use prec, only: & use prec, only: &
pInt, & pInt, &
@ -177,7 +175,7 @@ program DAMASK_spectral
C_ref = 0.0_pReal, & C_ref = 0.0_pReal, &
C = 0.0_pReal, & C = 0.0_pReal, &
S_lastInc, & S_lastInc, &
C_lastInc ! stiffness and compliance C_lastInc ! stiffness and compliance
real(pReal), dimension(6) :: sigma ! cauchy stress real(pReal), dimension(6) :: sigma ! cauchy stress
real(pReal), dimension(6,6) :: dsde real(pReal), dimension(6,6) :: dsde
@ -213,7 +211,7 @@ program DAMASK_spectral
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, errorID integer(pInt) :: i, j, k, l, m, n, p, errorID
integer(pInt) :: N_Loadcases, loadcase = 0_pInt, inc, iter, ielem, CPFEM_mode, & integer(pInt) :: N_Loadcases, loadcase = 0_pInt, inc, iter, ielem, CPFEM_mode=1_pInt, &
ierr, totalIncsCounter = 0_pInt,& ierr, totalIncsCounter = 0_pInt,&
notConvergedCounter = 0_pInt, convergedCounter = 0_pInt notConvergedCounter = 0_pInt, convergedCounter = 0_pInt
logical :: errmatinv logical :: errmatinv
@ -527,18 +525,51 @@ program DAMASK_spectral
enddo; enddo; enddo enddo; enddo; enddo
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! get reference material stiffness and init fields to no deformation ! in case of no restart get reference material stiffness and init fields to no deformation
ielem = 0_pInt if (restartInc == 1_pInt) then
do k = 1_pInt, res(3); do j = 1_pInt, res(2); do i = 1_pInt, res(1) ielem = 0_pInt
ielem = ielem + 1_pInt do k = 1_pInt, res(3); do j = 1_pInt, res(2); do i = 1_pInt, res(1)
F(i,j,k,1:3,1:3) = math_I3 ielem = ielem + 1_pInt
F_lastInc(i,j,k,1:3,1:3) = math_I3 F(i,j,k,1:3,1:3) = math_I3
coordinates(i,j,k,1:3) = geomdim/real(res,pReal)*real([i,j,k],pReal) - geomdim/real(2_pInt*res,pReal) F_lastInc(i,j,k,1:3,1:3) = math_I3
call CPFEM_general(2_pInt,coordinates(i,j,k,1:3),math_I3,math_I3,temperature(i,j,k),& coordinates(i,j,k,1:3) = geomdim/real(res,pReal)*real([i,j,k],pReal) - geomdim/real(2_pInt*res,pReal)
0.0_pReal,ielem,1_pInt,sigma,dsde,P_real(i,j,k,1:3,1:3),dPdF) call CPFEM_general(2_pInt,coordinates(i,j,k,1:3),math_I3,math_I3,temperature(i,j,k),&
C = C + dPdF 0.0_pReal,ielem,1_pInt,sigma,dsde,P_real(i,j,k,1:3,1:3),dPdF)
enddo; enddo; enddo C = C + dPdF
C_ref = C * wgt enddo; enddo; enddo
C_ref = C * wgt
call IO_write_jobBinaryFile(777,'C_ref',size(C_ref))
write (777,rec=1) C_ref
close(777)
call IO_write_jobBinaryFile(777,'C',size(C))
write (777,rec=1) C
close(777)
!--------------------------------------------------------------------------------------------------
! restore deformation gradient and stiffness from saved state
elseif (restartInc > 1_pInt) then ! using old values from file
if (debugRestart) write(6,'(a,i6,a)') 'Reading values of increment ',&
restartInc - 1_pInt,' from file'
call IO_read_jobBinaryFile(777,'convergedSpectralDefgrad',&
trim(getSolverJobName()),size(F))
read (777,rec=1) F
close (777)
F_lastInc = F
F_aim = 0.0_pReal
do k = 1_pInt, res(3); do j = 1_pInt, res(2); do i = 1_pInt, res(1)
F_aim = F_aim + F(i,j,k,1:3,1:3) ! calculating old average deformation
enddo; enddo; enddo
F_aim = F_aim * wgt
F_aim_lastInc = F_aim
coordinates = 0.0 ! change it later!!!
call IO_read_jobBinaryFile(777,'C_ref',trim(getSolverJobName()),size(C_ref))
read (777,rec=1) C_ref
close (777)
call IO_read_jobBinaryFile(777,'C',trim(getSolverJobName()),size(C))
read (777,rec=1) C
close (777)
CPFEM_mode = 2_pInt
endif
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! calculate the gamma operator ! calculate the gamma operator
@ -560,43 +591,31 @@ C_ref = C * wgt
gamma_hat(1,1,1, 1:3,1:3,1:3,1:3) = 0.0_pReal ! singular point at xi=(0.0,0.0,0.0) i.e. i=j=k=1 gamma_hat(1,1,1, 1:3,1:3,1:3,1:3) = 0.0_pReal ! singular point at xi=(0.0,0.0,0.0) i.e. i=j=k=1
endif endif
!--------------------------------------------------------------------------------------------------
! possible restore deformation gradient from saved state
if (restartInc > 1_pInt) then ! using old values from file
if (debugRestart) write(6,'(a,i6,a)') 'Reading values of increment ',&
restartInc - 1_pInt,' from file'
call IO_read_jobBinaryFile(777,'convergedSpectralDefgrad',&
trim(getSolverJobName()),size(F))
read (777,rec=1) F
close (777)
F_lastInc = F
F_aim = 0.0_pReal
do k = 1_pInt, res(3); do j = 1_pInt, res(2); do i = 1_pInt, res(1)
F_aim = F_aim + F(i,j,k,1:3,1:3) ! calculating old average deformation
enddo; enddo; enddo
F_aim = F_aim * wgt
F_aim_lastInc = F_aim
endif
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! write header of output file ! write header of output file
open(538,file=trim(getSolverWorkingDirectoryName())//trim(getSolverJobName())& if (appendToOutFile) then
//'.spectralOut',form='UNFORMATTED',status='REPLACE') open(538,file=trim(getSolverWorkingDirectoryName())//trim(getSolverJobName())//'.spectralOut',&
write(538) 'load', trim(loadCaseFile) form='UNFORMATTED',status='REPLACE')
write(538) 'workingdir', trim(getSolverWorkingDirectoryName()) else
write(538) 'geometry', trim(geometryFile) open(538,file=trim(getSolverWorkingDirectoryName())//trim(getSolverJobName())//'.spectralOut',&
write(538) 'resolution', res form='UNFORMATTED',status='REPLACE',position='APPEND')
write(538) 'dimension', geomdim write(538) 'load', trim(loadCaseFile)
write(538) 'materialpoint_sizeResults', materialpoint_sizeResults write(538) 'workingdir', trim(getSolverWorkingDirectoryName())
write(538) 'loadcases', N_Loadcases write(538) 'geometry', trim(geometryFile)
write(538) 'frequencies', bc(1:N_Loadcases)%outputfrequency ! one entry per loadcase write(538) 'resolution', res
write(538) 'times', bc(1:N_Loadcases)%time ! one entry per loadcase write(538) 'dimension', geomdim
write(538) 'logscales', bc(1:N_Loadcases)%logscale write(538) 'materialpoint_sizeResults', materialpoint_sizeResults
write(538) 'increments', bc(1:N_Loadcases)%incs ! one entry per loadcase write(538) 'loadcases', N_Loadcases
write(538) 'startingIncrement', restartInc - 1_pInt ! start with writing out the previous inc write(538) 'frequencies', bc(1:N_Loadcases)%outputfrequency ! one entry per loadcase
write(538) 'eoh' ! end of header write(538) 'times', bc(1:N_Loadcases)%time ! one entry per loadcase
write(538) materialpoint_results(1_pInt:materialpoint_sizeResults,1,1_pInt:Npoints) ! initial (non-deformed or read-in) results write(538) 'logscales', bc(1:N_Loadcases)%logscale
if (debugGeneral) write(6,'(a)') 'Header of result file written out' write(538) 'increments', bc(1:N_Loadcases)%incs ! one entry per loadcase
write(538) 'startingIncrement', restartInc - 1_pInt ! start with writing out the previous inc
write(538) 'eoh' ! end of header
write(538) materialpoint_results(1_pInt:materialpoint_sizeResults,1,1_pInt:Npoints) ! initial (non-deformed or read-in) results
if (debugGeneral) write(6,'(a)') 'Header of result file written out'
endif
flush(538)
!################################################################################################## !##################################################################################################
! Loop over loadcases defined in the loadcase file ! Loop over loadcases defined in the loadcase file
@ -712,7 +731,6 @@ C_ref = C * wgt
write(6,'(A,I5.5,A,es12.5)') 'Increment ', totalIncsCounter, ' Time ',time write(6,'(A,I5.5,A,es12.5)') 'Increment ', totalIncsCounter, ' Time ',time
guessmode = 1.0_pReal ! keep guessing along former trajectory during same loadcase guessmode = 1.0_pReal ! keep guessing along former trajectory during same loadcase
CPFEM_mode = 1_pInt ! winding forward
iter = 0_pInt iter = 0_pInt
err_div = huge(err_div_tol) ! go into loop err_div = huge(err_div_tol) ! go into loop
@ -1015,6 +1033,7 @@ C_ref = C * wgt
enddo ! end looping when convergency is achieved enddo ! end looping when convergency is achieved
CPFEM_mode = 1_pInt ! winding forward
write(6,'(a)') '' write(6,'(a)') ''
write(6,'(a)') '==================================================================' write(6,'(a)') '=================================================================='
if(err_div > err_div_tol .or. err_stress > err_stress_tol) then if(err_div > err_div_tol .or. err_stress > err_stress_tol) then
@ -1029,6 +1048,7 @@ C_ref = C * wgt
write(6,'(a)') '' write(6,'(a)') ''
write(6,'(a)') '... writing results to file ......................................' write(6,'(a)') '... writing results to file ......................................'
write(538) materialpoint_results(1_pInt:materialpoint_sizeResults,1,1_pInt:Npoints) ! write result to file write(538) materialpoint_results(1_pInt:materialpoint_sizeResults,1,1_pInt:Npoints) ! write result to file
flush(538)
endif endif
if( bc(loadcase)%restartFrequency > 0_pInt .and. & if( bc(loadcase)%restartFrequency > 0_pInt .and. &
@ -1065,5 +1085,6 @@ C_ref = C * wgt
call fftw_destroy_plan(plan_scalarField_forth) call fftw_destroy_plan(plan_scalarField_forth)
call fftw_destroy_plan(plan_scalarField_back) call fftw_destroy_plan(plan_scalarField_back)
endif endif
if (notConvergedCounter > 0_pInt) call quit(2_pInt)
call quit(0_pInt) call quit(0_pInt)
end program DAMASK_spectral end program DAMASK_spectral

View File

@ -73,25 +73,25 @@ module IO
contains contains
!******************************************************************** !--------------------------------------------------------------------------------------------------
! output version number !> @brief only output of revision number
!******************************************************************** !--------------------------------------------------------------------------------------------------
subroutine IO_init subroutine IO_init
use, intrinsic :: iso_fortran_env ! to get compiler_version and compiler_options (at least for gfortran 4.6 at the moment) use, intrinsic :: iso_fortran_env ! to get compiler_version and compiler_options (at least for gfortran 4.6 at the moment)
!$OMP CRITICAL (write2out)
write(6,*) write(6,*)
write(6,*) '<<<+- IO init -+>>>' write(6,*) '<<<+- IO init -+>>>'
write(6,*) '$Id$' write(6,*) '$Id$'
#include "compilation_info.f90" #include "compilation_info.f90"
flush(6) flush(6)
!$OMP END CRITICAL (write2out)
end subroutine IO_init end subroutine IO_init
!********************************************************************
! checks if myUnit is opened for reading and rewinds !--------------------------------------------------------------------------------------------------
!******************************************************************** !> @brief Checks if unit is opened for reading, if true rewinds. Otherwise stops with
!! error message 102
!--------------------------------------------------------------------------------------------------
subroutine IO_checkAndRewind(myUnit) subroutine IO_checkAndRewind(myUnit)
implicit none implicit none
@ -104,13 +104,13 @@ implicit none
end subroutine IO_checkAndRewind end subroutine IO_checkAndRewind
!********************************************************************
! open existing file to given myUnit
! path to file is relative to working directory
!********************************************************************
logical function IO_open_file_stat(myUnit,relPath)
use DAMASK_interface, only: getSolverWorkingDirectoryName !--------------------------------------------------------------------------------------------------
!> @brief Open existing file to given unit path to file is relative to working directory
!--------------------------------------------------------------------------------------------------
logical function IO_open_file_stat(myUnit,relPath)
use DAMASK_interface, &
only: getSolverWorkingDirectoryName
implicit none implicit none
integer(pInt), intent(in) :: myUnit integer(pInt), intent(in) :: myUnit
@ -126,14 +126,13 @@ logical function IO_open_file_stat(myUnit,relPath)
end function IO_open_file_stat end function IO_open_file_stat
!******************************************************************** !--------------------------------------------------------------------------------------------------
! open (write) file related to current job !> @brief Open (write) file related to current job but with different extension to given unit
! but with different extension to given myUnit !--------------------------------------------------------------------------------------------------
!********************************************************************
logical function IO_open_jobFile_stat(myUnit,newExt) logical function IO_open_jobFile_stat(myUnit,newExt)
use DAMASK_interface, only: &
use DAMASK_interface, only: getSolverWorkingDirectoryName, & getSolverWorkingDirectoryName, &
getSolverJobName getSolverJobName
implicit none implicit none
integer(pInt), intent(in) :: myUnit integer(pInt), intent(in) :: myUnit
@ -149,10 +148,9 @@ logical function IO_open_jobFile_stat(myUnit,newExt)
end function IO_open_JobFile_stat end function IO_open_JobFile_stat
!******************************************************************** !--------------------------------------------------------------------------------------------------
! open existing file to given unit !> @brief Open existing file to given unit path to file is relative to working directory
! path to file is relative to working directory !--------------------------------------------------------------------------------------------------
!********************************************************************
subroutine IO_open_file(myUnit,relPath) subroutine IO_open_file(myUnit,relPath)
use DAMASK_interface, only: getSolverWorkingDirectoryName use DAMASK_interface, only: getSolverWorkingDirectoryName
@ -171,10 +169,9 @@ subroutine IO_open_file(myUnit,relPath)
end subroutine IO_open_file end subroutine IO_open_file
!******************************************************************** !--------------------------------------------------------------------------------------------------
! open (write) file related to current job !> @brief Open (write) file related to current job but with different extension to given unit
! but with different extension to given unit !--------------------------------------------------------------------------------------------------
!********************************************************************
subroutine IO_open_jobFile(myUnit,newExt) subroutine IO_open_jobFile(myUnit,newExt)
use DAMASK_interface, only: getSolverWorkingDirectoryName, & use DAMASK_interface, only: getSolverWorkingDirectoryName, &

View File

@ -43,11 +43,13 @@ LAPACKROOT :=/usr
F90 ?=ifort F90 ?=ifort
COMPILERNAME ?= $(F90) COMPILERNAME ?= $(F90)
INCLUDE_DIRS +=-I$(DAMASK_ROOT)/lib INCLUDE_DIRS +=-I$(DAMASK_ROOT)/lib
ifdef PETSC_DIR ifdef PETSC_DIR
INCLUDE_DIRS +=-I$(PETSC_DIR)/include INCLUDE_DIRS +=-I$(PETSC_DIR)/include
INCLUDE_DIRS +=-I$(PETSC_DIR)/$(PETSC_ARCH)/include INCLUDE_DIRS +=-I$(PETSC_DIR)/$(PETSC_ARCH)/include
INCLUDE_DIRS +=-DPETSC # just for the moment, as long as PETSC is non standard # just for the moment, as long as PETSC is non standard
INCLUDE_DIRS +=-DPETSC
endif endif
ifeq "$(FASTBUILD)" "YES" ifeq "$(FASTBUILD)" "YES"

View File

@ -54,6 +54,7 @@ subroutine quit(stop_id)
dateAndTime(6),':',& dateAndTime(6),':',&
dateAndTime(7) dateAndTime(7)
if (stop_id == 0_pInt) stop 0 ! normal termination if (stop_id == 0_pInt) stop 0 ! normal termination
if (stop_id == 2_pInt) stop 2 ! not all steps converged
if (stop_id < 0_pInt) then ! trigger regridding if (stop_id < 0_pInt) then ! trigger regridding
write(0,'(a,i6)') 'restart a', stop_id*(-1_pInt) write(0,'(a,i6)') 'restart a', stop_id*(-1_pInt)
stop 2 stop 2