unified style

This commit is contained in:
Martin Diehl 2019-04-15 15:53:46 +02:00
parent 18f9deef1a
commit 4793f964f8
2 changed files with 430 additions and 433 deletions

View File

@ -27,8 +27,7 @@ module grid_mech_spectral_basic
type(tSolutionParams), private :: params type(tSolutionParams), private :: params
type, private :: tNumerics type, private :: tNumerics
logical :: & logical :: update_gamma !< update gamma operator with current stiffness
update_gamma !< update gamma operator with current stiffness
end type tNumerics end type tNumerics
type(tNumerics) :: num ! numerics parameters. Better name? type(tNumerics) :: num ! numerics parameters. Better name?
@ -42,8 +41,8 @@ module grid_mech_spectral_basic
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! common pointwise data ! common pointwise data
real(pReal), private, dimension(:,:,:,:,:), allocatable :: & real(pReal), private, dimension(:,:,:,:,:), allocatable :: &
F_lastInc, & F_lastInc, & !< field of previous compatible deformation gradients
Fdot Fdot !< field of assumed rate of compatible deformation gradient
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! stress, stiffness and compliance average etc. ! stress, stiffness and compliance average etc.
@ -101,8 +100,7 @@ subroutine grid_mech_spectral_basic_init
use spectral_utilities, only: & use spectral_utilities, only: &
utilities_constitutiveResponse, & utilities_constitutiveResponse, &
utilities_updateGamma, & utilities_updateGamma, &
utilities_updateIPcoords, & utilities_updateIPcoords
wgt
use mesh, only: & use mesh, only: &
grid, & grid, &
grid3 grid3
@ -179,6 +177,7 @@ subroutine grid_mech_spectral_basic_init
write(rankStr,'(a1,i0)')'_',worldrank write(rankStr,'(a1,i0)')'_',worldrank
fileHandle = HDF5_openFile(trim(getSolverJobName())//trim(rankStr)//'.hdf5') fileHandle = HDF5_openFile(trim(getSolverJobName())//trim(rankStr)//'.hdf5')
call HDF5_read(fileHandle,F_aim, 'F_aim') call HDF5_read(fileHandle,F_aim, 'F_aim')
call HDF5_read(fileHandle,F_aim_lastInc,'F_aim_lastInc') call HDF5_read(fileHandle,F_aim_lastInc,'F_aim_lastInc')
call HDF5_read(fileHandle,F_aimDot, 'F_aimDot') call HDF5_read(fileHandle,F_aimDot, 'F_aimDot')
@ -226,7 +225,6 @@ function grid_mech_spectral_basic_solution(incInfoIn,timeinc,timeinc_old,stress_
terminallyIll terminallyIll
implicit none implicit none
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! input data for solution ! input data for solution
character(len=*), intent(in) :: & character(len=*), intent(in) :: &
@ -248,8 +246,8 @@ function grid_mech_spectral_basic_solution(incInfoIn,timeinc,timeinc_old,stress_
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! update stiffness (and gamma operator) ! update stiffness (and gamma operator)
S = Utilities_maskedCompliance(rotation_BC,stress_BC%maskLogical,C_volAvg) S = utilities_maskedCompliance(rotation_BC,stress_BC%maskLogical,C_volAvg)
if (num%update_gamma) call Utilities_updateGamma(C_minMaxAvg,restartWrite) if (num%update_gamma) call utilities_updateGamma(C_minMaxAvg,restartWrite)
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! set module wide available data ! set module wide available data
@ -329,7 +327,7 @@ subroutine grid_mech_spectral_basic_forward(guess,timeinc,timeinc_old,loadCaseTi
else else
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! restart information for spectral solver ! restart information for spectral solver
if (restartWrite) then ! QUESTION: where is this logical properly set? if (restartWrite) then
write(6,'(/,a)') ' writing converged results for restart';flush(6) write(6,'(/,a)') ' writing converged results for restart';flush(6)
write(rankStr,'(a1,i0)')'_',worldrank write(rankStr,'(a1,i0)')'_',worldrank
@ -344,6 +342,7 @@ subroutine grid_mech_spectral_basic_forward(guess,timeinc,timeinc_old,loadCaseTi
call HDF5_write(fileHandle,C_volAvg, 'C_volAvg') call HDF5_write(fileHandle,C_volAvg, 'C_volAvg')
call HDF5_write(fileHandle,C_volAvgLastInc,'C_volAvgLastInc') call HDF5_write(fileHandle,C_volAvgLastInc,'C_volAvgLastInc')
call HDF5_write(fileHandle,C_minMaxAvg, 'C_minMaxAvg') call HDF5_write(fileHandle,C_minMaxAvg, 'C_minMaxAvg')
call HDF5_closeFile(fileHandle) call HDF5_closeFile(fileHandle)
endif endif
@ -390,7 +389,7 @@ end subroutine grid_mech_spectral_basic_forward
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief convergence check !> @brief convergence check
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
subroutine converged(snes_local,PETScIter,xnorm,snorm,fnorm,reason,dummy,ierr) subroutine converged(snes_local,PETScIter,devNull1,devNull2,devNull3,reason,dummy,ierr)
use numerics, only: & use numerics, only: &
itmax, & itmax, &
itmin, & itmin, &
@ -403,11 +402,11 @@ subroutine converged(snes_local,PETScIter,xnorm,snorm,fnorm,reason,dummy,ierr)
implicit none implicit none
SNES :: snes_local SNES :: snes_local
PetscInt :: PETScIter PetscInt, intent(in) :: PETScIter
PetscReal :: & PetscReal, intent(in) :: &
xnorm, & ! not used devNull1, &
snorm, & ! not used devNull2, &
fnorm ! not used devNull3
SNESConvergedReason :: reason SNESConvergedReason :: reason
PetscObject :: dummy PetscObject :: dummy
PetscErrorCode :: ierr PetscErrorCode :: ierr
@ -506,7 +505,7 @@ subroutine formResidual(in, F, &
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! evaluate constitutive response ! evaluate constitutive response
call Utilities_constitutiveResponse(residuum, & ! "residuum" gets field of first PK stress (to save memory) call utilities_constitutiveResponse(residuum, & ! "residuum" gets field of first PK stress (to save memory)
P_av,C_volAvg,C_minMaxAvg, & P_av,C_volAvg,C_minMaxAvg, &
F,params%timeinc,params%rotation_BC) F,params%timeinc,params%rotation_BC)
call MPI_Allreduce(MPI_IN_PLACE,terminallyIll,1,MPI_LOGICAL,MPI_LOR,PETSC_COMM_WORLD,ierr) call MPI_Allreduce(MPI_IN_PLACE,terminallyIll,1,MPI_LOGICAL,MPI_LOR,PETSC_COMM_WORLD,ierr)

View File

@ -27,8 +27,7 @@ module grid_mech_spectral_polarisation
type(tSolutionParams), private :: params type(tSolutionParams), private :: params
type, private :: tNumerics type, private :: tNumerics
logical :: & logical :: update_gamma !< update gamma operator with current stiffness
update_gamma !< update gamma operator with current stiffness
end type tNumerics end type tNumerics
type(tNumerics) :: num ! numerics parameters. Better name? type(tNumerics) :: num ! numerics parameters. Better name?
@ -145,10 +144,10 @@ subroutine grid_mech_spectral_polarisation_init
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! allocate global fields ! allocate global fields
allocate (F_lastInc (3,3,grid(1),grid(2),grid3),source = 0.0_pReal) allocate(F_lastInc (3,3,grid(1),grid(2),grid3),source = 0.0_pReal)
allocate (Fdot (3,3,grid(1),grid(2),grid3),source = 0.0_pReal) allocate(Fdot (3,3,grid(1),grid(2),grid3),source = 0.0_pReal)
allocate (F_tau_lastInc(3,3,grid(1),grid(2),grid3),source = 0.0_pReal) allocate(F_tau_lastInc(3,3,grid(1),grid(2),grid3),source = 0.0_pReal)
allocate (F_tauDot (3,3,grid(1),grid(2),grid3),source = 0.0_pReal) allocate(F_tauDot (3,3,grid(1),grid(2),grid3),source = 0.0_pReal)
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! initialize solver specific parts of PETSc ! initialize solver specific parts of PETSc
@ -194,7 +193,7 @@ subroutine grid_mech_spectral_polarisation_init
call HDF5_read(fileHandle,F, 'F') call HDF5_read(fileHandle,F, 'F')
call HDF5_read(fileHandle,F_lastInc, 'F_lastInc') call HDF5_read(fileHandle,F_lastInc, 'F_lastInc')
call HDF5_read(fileHandle,F_tau, 'F_tau') call HDF5_read(fileHandle,F_tau, 'F_tau')
call HDF5_read(fileHandle,F_tau_lastInc, 'F_tau_lastInc') call HDF5_read(fileHandle,F_tau_lastInc,'F_tau_lastInc')
elseif (restartInc == 0) then restart elseif (restartInc == 0) then restart
F_lastInc = spread(spread(spread(math_I3,3,grid(1)),4,grid(2)),5,grid3) ! initialize to identity F_lastInc = spread(spread(spread(math_I3,3,grid(1)),4,grid(2)),5,grid3) ! initialize to identity
@ -216,6 +215,7 @@ subroutine grid_mech_spectral_polarisation_init
call HDF5_read(fileHandle,C_volAvg, 'C_volAvg') call HDF5_read(fileHandle,C_volAvg, 'C_volAvg')
call HDF5_read(fileHandle,C_volAvgLastInc,'C_volAvgLastInc') call HDF5_read(fileHandle,C_volAvgLastInc,'C_volAvgLastInc')
call HDF5_closeFile(fileHandle) call HDF5_closeFile(fileHandle)
fileUnit = IO_open_jobFile_binary('C_ref') fileUnit = IO_open_jobFile_binary('C_ref')
read(fileUnit) C_minMaxAvg; close(fileUnit) read(fileUnit) C_minMaxAvg; close(fileUnit)
endif restartRead endif restartRead
@ -242,7 +242,6 @@ function grid_mech_spectral_polarisation_solution(incInfoIn,timeinc,timeinc_old,
terminallyIll terminallyIll
implicit none implicit none
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! input data for solution ! input data for solution
character(len=*), intent(in) :: & character(len=*), intent(in) :: &
@ -264,7 +263,7 @@ function grid_mech_spectral_polarisation_solution(incInfoIn,timeinc,timeinc_old,
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! update stiffness (and gamma operator) ! update stiffness (and gamma operator)
S = Utilities_maskedCompliance(rotation_BC,stress_BC%maskLogical,C_volAvg) S = utilities_maskedCompliance(rotation_BC,stress_BC%maskLogical,C_volAvg)
if (num%update_gamma) then if (num%update_gamma) then
call utilities_updateGamma(C_minMaxAvg,restartWrite) call utilities_updateGamma(C_minMaxAvg,restartWrite)
C_scale = C_minMaxAvg C_scale = C_minMaxAvg
@ -324,7 +323,6 @@ subroutine grid_mech_spectral_polarisation_forward(guess,timeinc,timeinc_old,loa
use FEsolving, only: & use FEsolving, only: &
restartWrite restartWrite
implicit none
logical, intent(in) :: & logical, intent(in) :: &
guess guess
real(pReal), intent(in) :: & real(pReal), intent(in) :: &
@ -354,25 +352,24 @@ subroutine grid_mech_spectral_polarisation_forward(guess,timeinc,timeinc_old,loa
else else
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! restart information for spectral solver ! restart information for spectral solver
if (restartWrite) then ! QUESTION: where is this logical properly set? if (restartWrite) then
write(6,'(/,a)') ' writing converged results for restart';flush(6) write(6,'(/,a)') ' writing converged results for restart';flush(6)
write(rankStr,'(a1,i0)')'_',worldrank write(rankStr,'(a1,i0)')'_',worldrank
fileHandle = HDF5_openFile(trim(getSolverJobName())//trim(rankStr)//'.hdf5','w') fileHandle = HDF5_openFile(trim(getSolverJobName())//trim(rankStr)//'.hdf5','w')
call HDF5_write(fileHandle,C_volAvg, 'C_volAvg')
call HDF5_write(fileHandle,C_volAvgLastInc,'C_volAvgLastInc')
call HDF5_write(fileHandle,F_aim, 'F_aim') call HDF5_write(fileHandle,F_aim, 'F_aim')
call HDF5_write(fileHandle,F_aimDot, 'F_aimDot')
call HDF5_write(fileHandle,F_aim_lastInc, 'F_aim_lastInc') call HDF5_write(fileHandle,F_aim_lastInc, 'F_aim_lastInc')
call HDF5_write(fileHandle,F_aimDot, 'F_aimDot')
call HDF5_write(fileHandle,F, 'F') call HDF5_write(fileHandle,F, 'F')
call HDF5_write(fileHandle,F_lastInc, 'F_lastInc') call HDF5_write(fileHandle,F_lastInc, 'F_lastInc')
call HDF5_write(fileHandle,F_tau, 'F_tau') call HDF5_write(fileHandle,F_tau, 'F_tau')
call HDF5_write(fileHandle,F_tau_lastInc, 'F_tau_lastInc') call HDF5_write(fileHandle,F_tau_lastInc, 'F_tau_lastInc')
call HDF5_closeFile(fileHandle)
call HDF5_write(fileHandle,C_volAvg, 'C_volAvg')
call HDF5_write(fileHandle,C_volAvgLastInc,'C_volAvgLastInc')
call HDF5_closeFile(fileHandle)
endif endif
call CPFEM_age ! age state and kinematics call CPFEM_age ! age state and kinematics
@ -438,7 +435,7 @@ end subroutine grid_mech_spectral_polarisation_forward
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief convergence check !> @brief convergence check
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
subroutine converged(snes_local,PETScIter,xnorm,snorm,fnorm,reason,dummy,ierr) subroutine converged(snes_local,PETScIter,devNull1,devNull2,devNull3,reason,dummy,ierr)
use numerics, only: & use numerics, only: &
itmax, & itmax, &
itmin, & itmin, &
@ -453,11 +450,11 @@ subroutine converged(snes_local,PETScIter,xnorm,snorm,fnorm,reason,dummy,ierr)
implicit none implicit none
SNES :: snes_local SNES :: snes_local
PetscInt :: PETScIter PetscInt, intent(in) :: PETScIter
PetscReal :: & PetscReal, intent(in) :: &
xnorm, & ! not used devNull1, &
snorm, & ! not used devNull2, &
fnorm ! not used devNull3
SNESConvergedReason :: reason SNESConvergedReason :: reason
PetscObject :: dummy PetscObject :: dummy
PetscErrorCode :: ierr PetscErrorCode :: ierr
@ -622,6 +619,7 @@ subroutine formResidual(in, FandF_tau, &
tensorField_real(1:3,1:3,1:grid(1),1:grid(2),1:grid3) = residual_F !< stress field in disguise tensorField_real(1:3,1:3,1:grid(1),1:grid(2),1:grid3) = residual_F !< stress field in disguise
call utilities_FFTtensorForward call utilities_FFTtensorForward
err_div = Utilities_divergenceRMS() !< root mean squared error in divergence of stress err_div = Utilities_divergenceRMS() !< root mean squared error in divergence of stress
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
! constructing residual ! constructing residual
e = 0 e = 0