store data where it is needed
avoid globals, even if they are read only
This commit is contained in:
parent
0bb7fea782
commit
e4792e56fb
|
@ -43,6 +43,23 @@ submodule(homogenization) homogenization_mech_RGC
|
||||||
orientation
|
orientation
|
||||||
end type tRGCdependentState
|
end type tRGCdependentState
|
||||||
|
|
||||||
|
type :: tNumerics
|
||||||
|
real(pReal) :: &
|
||||||
|
atol, & !< absolute tolerance of RGC residuum
|
||||||
|
rtol, & !< relative tolerance of RGC residuum
|
||||||
|
absMax, & !< absolute maximum of RGC residuum
|
||||||
|
relMax, & !< relative maximum of RGC residuum
|
||||||
|
pPert, & !< perturbation for computing RGC penalty tangent
|
||||||
|
xSmoo, & !< RGC penalty smoothing parameter (hyperbolic tangent)
|
||||||
|
viscPower, & !< power (sensitivity rate) of numerical viscosity in RGC scheme, Default 1.0e0: Newton viscosity (linear model)
|
||||||
|
viscModus, & !< stress modulus of RGC numerical viscosity, Default 0.0e0: No viscosity is applied
|
||||||
|
refRelaxRate, & !< reference relaxation rate in RGC viscosity
|
||||||
|
maxdRelax, & !< threshold of maximum relaxation vector increment (if exceed this then cutback)
|
||||||
|
maxVolDiscr, & !< threshold of maximum volume discrepancy allowed
|
||||||
|
volDiscrMod, & !< stiffness of RGC volume discrepancy (zero = without volume discrepancy constraint)
|
||||||
|
volDiscrPow !< powerlaw penalty for volume discrepancy
|
||||||
|
end type tNumerics
|
||||||
|
|
||||||
type(tparameters), dimension(:), allocatable :: &
|
type(tparameters), dimension(:), allocatable :: &
|
||||||
param
|
param
|
||||||
type(tRGCstate), dimension(:), allocatable :: &
|
type(tRGCstate), dimension(:), allocatable :: &
|
||||||
|
@ -50,6 +67,8 @@ submodule(homogenization) homogenization_mech_RGC
|
||||||
state0
|
state0
|
||||||
type(tRGCdependentState), dimension(:), allocatable :: &
|
type(tRGCdependentState), dimension(:), allocatable :: &
|
||||||
dependentState
|
dependentState
|
||||||
|
type(tNumerics) :: &
|
||||||
|
num ! numerics parameters. Better name?
|
||||||
|
|
||||||
contains
|
contains
|
||||||
|
|
||||||
|
@ -64,7 +83,6 @@ module subroutine mech_RGC_init
|
||||||
NofMyHomog, &
|
NofMyHomog, &
|
||||||
sizeState, nIntFaceTot
|
sizeState, nIntFaceTot
|
||||||
|
|
||||||
|
|
||||||
write(6,'(/,a)') ' <<<+- homogenization_'//HOMOGENIZATION_RGC_label//' init -+>>>'; flush(6)
|
write(6,'(/,a)') ' <<<+- homogenization_'//HOMOGENIZATION_RGC_label//' init -+>>>'; flush(6)
|
||||||
|
|
||||||
write(6,'(/,a)') ' Tjahjanto et al., International Journal of Material Forming 2(1):939–942, 2009'
|
write(6,'(/,a)') ' Tjahjanto et al., International Journal of Material Forming 2(1):939–942, 2009'
|
||||||
|
@ -82,6 +100,33 @@ module subroutine mech_RGC_init
|
||||||
allocate(state0(Ninstance))
|
allocate(state0(Ninstance))
|
||||||
allocate(dependentState(Ninstance))
|
allocate(dependentState(Ninstance))
|
||||||
|
|
||||||
|
num%atol = config_numerics%getFloat('atol_rgc', defaultVal=1.0e+4_pReal)
|
||||||
|
num%rtol = config_numerics%getFloat('rtol_rgc', defaultVal=1.0e-3_pReal)
|
||||||
|
num%absMax = config_numerics%getFloat('amax_rgc', defaultVal=1.0e+10_pReal)
|
||||||
|
num%relMax = config_numerics%getFloat('rmax_rgc', defaultVal=1.0e+2_pReal)
|
||||||
|
num%pPert = config_numerics%getFloat('perturbpenalty_rgc', defaultVal=1.0e-7_pReal)
|
||||||
|
num%xSmoo = config_numerics%getFloat('relvantmismatch_rgc', defaultVal=1.0e-5_pReal)
|
||||||
|
num%viscPower = config_numerics%getFloat('viscositypower_rgc', defaultVal=1.0e+0_pReal)
|
||||||
|
num%viscModus = config_numerics%getFloat('viscositymodulus_rgc', defaultVal=0.0e+0_pReal)
|
||||||
|
num%refRelaxRate = config_numerics%getFloat('refrelaxationrate_rgc',defaultVal=1.0e-3_pReal)
|
||||||
|
num%maxdRelax = config_numerics%getFloat('maxrelaxationrate_rgc',defaultVal=1.0e+0_pReal)
|
||||||
|
num%maxVolDiscr = config_numerics%getFloat('maxvoldiscrepancy_rgc',defaultVal=1.0e-5_pReal)
|
||||||
|
num%volDiscrMod = config_numerics%getFloat('voldiscrepancymod_rgc',defaultVal=1.0e+12_pReal)
|
||||||
|
num%volDiscrPow = config_numerics%getFloat('dicrepancypower_rgc', defaultVal=5.0_pReal)
|
||||||
|
|
||||||
|
if (num%atol <= 0.0_pReal) call IO_error(301,ext_msg='absTol_RGC')
|
||||||
|
if (num%rtol <= 0.0_pReal) call IO_error(301,ext_msg='relTol_RGC')
|
||||||
|
if (num%absMax <= 0.0_pReal) call IO_error(301,ext_msg='absMax_RGC')
|
||||||
|
if (num%relMax <= 0.0_pReal) call IO_error(301,ext_msg='relMax_RGC')
|
||||||
|
if (num%pPert <= 0.0_pReal) call IO_error(301,ext_msg='pPert_RGC')
|
||||||
|
if (num%xSmoo <= 0.0_pReal) call IO_error(301,ext_msg='xSmoo_RGC')
|
||||||
|
if (num%viscPower < 0.0_pReal) call IO_error(301,ext_msg='viscPower_RGC')
|
||||||
|
if (num%viscModus < 0.0_pReal) call IO_error(301,ext_msg='viscModus_RGC')
|
||||||
|
if (num%refRelaxRate <= 0.0_pReal) call IO_error(301,ext_msg='refRelaxRate_RGC')
|
||||||
|
if (num%maxdRelax <= 0.0_pReal) call IO_error(301,ext_msg='maxdRelax_RGC')
|
||||||
|
if (num%maxVolDiscr <= 0.0_pReal) call IO_error(301,ext_msg='maxVolDiscr_RGC')
|
||||||
|
if (num%volDiscrMod < 0.0_pReal) call IO_error(301,ext_msg='volDiscrMod_RGC')
|
||||||
|
if (num%volDiscrPow <= 0.0_pReal) call IO_error(301,ext_msg='volDiscrPw_RGC')
|
||||||
|
|
||||||
do h = 1, size(homogenization_type)
|
do h = 1, size(homogenization_type)
|
||||||
if (homogenization_type(h) /= HOMOGENIZATION_RGC_ID) cycle
|
if (homogenization_type(h) /= HOMOGENIZATION_RGC_ID) cycle
|
||||||
|
@ -97,6 +142,8 @@ module subroutine mech_RGC_init
|
||||||
endif
|
endif
|
||||||
#endif
|
#endif
|
||||||
|
|
||||||
|
prm%output = config%getStrings('(output)',defaultVal=emptyStringArray)
|
||||||
|
|
||||||
prm%Nconstituents = config%getInts('clustersize',requiredSize=3)
|
prm%Nconstituents = config%getInts('clustersize',requiredSize=3)
|
||||||
if (homogenization_Ngrains(h) /= product(prm%Nconstituents)) &
|
if (homogenization_Ngrains(h) /= product(prm%Nconstituents)) &
|
||||||
call IO_error(211,ext_msg='clustersize ('//HOMOGENIZATION_RGC_label//')')
|
call IO_error(211,ext_msg='clustersize ('//HOMOGENIZATION_RGC_label//')')
|
||||||
|
@ -107,8 +154,6 @@ module subroutine mech_RGC_init
|
||||||
prm%dAlpha = config%getFloats('grainsize', requiredSize=3)
|
prm%dAlpha = config%getFloats('grainsize', requiredSize=3)
|
||||||
prm%angles = config%getFloats('clusterorientation',requiredSize=3)
|
prm%angles = config%getFloats('clusterorientation',requiredSize=3)
|
||||||
|
|
||||||
prm%output = config%getStrings('(output)',defaultVal=emptyStringArray)
|
|
||||||
|
|
||||||
NofMyHomog = count(material_homogenizationAt == h)
|
NofMyHomog = count(material_homogenizationAt == h)
|
||||||
nIntFaceTot = 3*( (prm%Nconstituents(1)-1)*prm%Nconstituents(2)*prm%Nconstituents(3) &
|
nIntFaceTot = 3*( (prm%Nconstituents(1)-1)*prm%Nconstituents(2)*prm%Nconstituents(3) &
|
||||||
+ prm%Nconstituents(1)*(prm%Nconstituents(2)-1)*prm%Nconstituents(3) &
|
+ prm%Nconstituents(1)*(prm%Nconstituents(2)-1)*prm%Nconstituents(3) &
|
||||||
|
@ -296,7 +341,7 @@ module procedure mech_RGC_updateState
|
||||||
!--------------------------------------------------------------------------------------------------
|
!--------------------------------------------------------------------------------------------------
|
||||||
! compute the residual of traction at the interface (in local system, 4-dimensional index)
|
! compute the residual of traction at the interface (in local system, 4-dimensional index)
|
||||||
do i = 1,3
|
do i = 1,3
|
||||||
tract(iNum,i) = sign(viscModus_RGC*(abs(drelax(i+3*(iNum-1)))/(refRelaxRate_RGC*dt))**viscPower_RGC, &
|
tract(iNum,i) = sign(num%viscModus*(abs(drelax(i+3*(iNum-1)))/(num%refRelaxRate*dt))**num%viscPower, &
|
||||||
drelax(i+3*(iNum-1))) ! contribution from the relaxation viscosity
|
drelax(i+3*(iNum-1))) ! contribution from the relaxation viscosity
|
||||||
do j = 1,3
|
do j = 1,3
|
||||||
tract(iNum,i) = tract(iNum,i) + (P(i,j,iGrP) + R(i,j,iGrP) + D(i,j,iGrP))*normP(j) & ! contribution from material stress P, mismatch penalty R, and volume penalty D projected into the interface
|
tract(iNum,i) = tract(iNum,i) + (P(i,j,iGrP) + R(i,j,iGrP) + D(i,j,iGrP))*normP(j) & ! contribution from material stress P, mismatch penalty R, and volume penalty D projected into the interface
|
||||||
|
@ -337,7 +382,7 @@ module procedure mech_RGC_updateState
|
||||||
|
|
||||||
!--------------------------------------------------------------------------------------------------
|
!--------------------------------------------------------------------------------------------------
|
||||||
! If convergence reached => done and happy
|
! If convergence reached => done and happy
|
||||||
if (residMax < relTol_RGC*stresMax .or. residMax < absTol_RGC) then
|
if (residMax < num%rtol*stresMax .or. residMax < num%atol) then
|
||||||
mech_RGC_updateState = .true.
|
mech_RGC_updateState = .true.
|
||||||
#ifdef DEBUG
|
#ifdef DEBUG
|
||||||
if (iand(debug_level(debug_homogenization),debug_levelExtensive) /= 0 .and. prm%of_debug == of) &
|
if (iand(debug_level(debug_homogenization),debug_levelExtensive) /= 0 .and. prm%of_debug == of) &
|
||||||
|
@ -377,7 +422,7 @@ module procedure mech_RGC_updateState
|
||||||
|
|
||||||
!--------------------------------------------------------------------------------------------------
|
!--------------------------------------------------------------------------------------------------
|
||||||
! if residual blows-up => done but unhappy
|
! if residual blows-up => done but unhappy
|
||||||
elseif (residMax > relMax_RGC*stresMax .or. residMax > absMax_RGC) then ! try to restart when residual blows up exceeding maximum bound
|
elseif (residMax > num%relMax*stresMax .or. residMax > num%absMax) then ! try to restart when residual blows up exceeding maximum bound
|
||||||
mech_RGC_updateState = [.true.,.false.] ! with direct cut-back
|
mech_RGC_updateState = [.true.,.false.] ! with direct cut-back
|
||||||
|
|
||||||
#ifdef DEBUG
|
#ifdef DEBUG
|
||||||
|
@ -465,7 +510,7 @@ module procedure mech_RGC_updateState
|
||||||
|
|
||||||
do ipert = 1,3*nIntFaceTot
|
do ipert = 1,3*nIntFaceTot
|
||||||
p_relax = relax
|
p_relax = relax
|
||||||
p_relax(ipert) = relax(ipert) + pPert_RGC ! perturb the relaxation vector
|
p_relax(ipert) = relax(ipert) + num%pPert ! perturb the relaxation vector
|
||||||
stt%relaxationVector(:,of) = p_relax
|
stt%relaxationVector(:,of) = p_relax
|
||||||
call grainDeformation(pF,avgF,instance,of) ! rain deformation from perturbed state
|
call grainDeformation(pF,avgF,instance,of) ! rain deformation from perturbed state
|
||||||
call stressPenalty(pR,DevNull, avgF,pF,ip,el,instance,of) ! stress penalty due to interface mismatch from perturbed state
|
call stressPenalty(pR,DevNull, avgF,pF,ip,el,instance,of) ! stress penalty due to interface mismatch from perturbed state
|
||||||
|
@ -502,7 +547,7 @@ module procedure mech_RGC_updateState
|
||||||
+ (pD(i,j,iGrN) - D(i,j,iGrN))*normN(j)
|
+ (pD(i,j,iGrN) - D(i,j,iGrN))*normN(j)
|
||||||
enddo; enddo
|
enddo; enddo
|
||||||
enddo
|
enddo
|
||||||
pmatrix(:,ipert) = p_resid/pPert_RGC
|
pmatrix(:,ipert) = p_resid/num%pPert
|
||||||
enddo
|
enddo
|
||||||
|
|
||||||
#ifdef DEBUG
|
#ifdef DEBUG
|
||||||
|
@ -520,8 +565,8 @@ module procedure mech_RGC_updateState
|
||||||
! ... of the numerical viscosity traction "rmatrix"
|
! ... of the numerical viscosity traction "rmatrix"
|
||||||
allocate(rmatrix(3*nIntFaceTot,3*nIntFaceTot),source=0.0_pReal)
|
allocate(rmatrix(3*nIntFaceTot,3*nIntFaceTot),source=0.0_pReal)
|
||||||
do i=1,3*nIntFaceTot
|
do i=1,3*nIntFaceTot
|
||||||
rmatrix(i,i) = viscModus_RGC*viscPower_RGC/(refRelaxRate_RGC*dt)* & ! tangent due to numerical viscosity traction appears
|
rmatrix(i,i) = num%viscModus*num%viscPower/(num%refRelaxRate*dt)* & ! tangent due to numerical viscosity traction appears
|
||||||
(abs(drelax(i))/(refRelaxRate_RGC*dt))**(viscPower_RGC - 1.0_pReal) ! only in the main diagonal term
|
(abs(drelax(i))/(num%refRelaxRate*dt))**(num%viscPower - 1.0_pReal) ! only in the main diagonal term
|
||||||
enddo
|
enddo
|
||||||
|
|
||||||
#ifdef DEBUG
|
#ifdef DEBUG
|
||||||
|
@ -573,7 +618,7 @@ module procedure mech_RGC_updateState
|
||||||
drelax(i) = drelax(i) - jnverse(i,j)*resid(j) ! Calculate the correction for the state variable
|
drelax(i) = drelax(i) - jnverse(i,j)*resid(j) ! Calculate the correction for the state variable
|
||||||
enddo; enddo
|
enddo; enddo
|
||||||
stt%relaxationVector(:,of) = relax + drelax ! Updateing the state variable for the next iteration
|
stt%relaxationVector(:,of) = relax + drelax ! Updateing the state variable for the next iteration
|
||||||
if (any(abs(drelax) > maxdRelax_RGC)) then ! Forcing cutback when the incremental change of relaxation vector becomes too large
|
if (any(abs(drelax) > num%maxdRelax)) then ! Forcing cutback when the incremental change of relaxation vector becomes too large
|
||||||
mech_RGC_updateState = [.true.,.false.]
|
mech_RGC_updateState = [.true.,.false.]
|
||||||
!$OMP CRITICAL (write2out)
|
!$OMP CRITICAL (write2out)
|
||||||
write(6,'(1x,a,1x,i3,1x,a,1x,i3,1x,a)')'RGC_updateState: ip',ip,'| el',el,'enforces cutback'
|
write(6,'(1x,a,1x,i3,1x,a,1x,i3,1x,a)')'RGC_updateState: ip',ip,'| el',el,'enforces cutback'
|
||||||
|
@ -690,7 +735,7 @@ module procedure mech_RGC_updateState
|
||||||
*surfCorr(abs(intFace(1)))/prm%dAlpha(abs(intFace(1))) &
|
*surfCorr(abs(intFace(1)))/prm%dAlpha(abs(intFace(1))) &
|
||||||
*cosh(prm%ciAlpha*nDefNorm) &
|
*cosh(prm%ciAlpha*nDefNorm) &
|
||||||
*0.5_pReal*nVect(l)*nDef(i,k)/nDefNorm*math_LeviCivita(k,l,j) &
|
*0.5_pReal*nVect(l)*nDef(i,k)/nDefNorm*math_LeviCivita(k,l,j) &
|
||||||
*tanh(nDefNorm/xSmoo_RGC)
|
*tanh(nDefNorm/num%xSmoo)
|
||||||
enddo; enddo;enddo; enddo
|
enddo; enddo;enddo; enddo
|
||||||
enddo interfaceLoop
|
enddo interfaceLoop
|
||||||
#ifdef DEBUG
|
#ifdef DEBUG
|
||||||
|
@ -738,8 +783,8 @@ module procedure mech_RGC_updateState
|
||||||
! calculate the stress and penalty due to volume discrepancy
|
! calculate the stress and penalty due to volume discrepancy
|
||||||
vPen = 0.0_pReal
|
vPen = 0.0_pReal
|
||||||
do i = 1,nGrain
|
do i = 1,nGrain
|
||||||
vPen(:,:,i) = -1.0_pReal/real(nGrain,pReal)*volDiscrMod_RGC*volDiscrPow_RGC/maxVolDiscr_RGC* &
|
vPen(:,:,i) = -1.0_pReal/real(nGrain,pReal)*num%volDiscrMod*num%volDiscrPow/num%maxVolDiscr* &
|
||||||
sign((abs(vDiscrep)/maxVolDiscr_RGC)**(volDiscrPow_RGC - 1.0),vDiscrep)* &
|
sign((abs(vDiscrep)/num%maxVolDiscr)**(num%volDiscrPow - 1.0),vDiscrep)* &
|
||||||
gVol(i)*transpose(math_inv33(fDef(:,:,i)))
|
gVol(i)*transpose(math_inv33(fDef(:,:,i)))
|
||||||
|
|
||||||
#ifdef DEBUG
|
#ifdef DEBUG
|
||||||
|
|
|
@ -17,12 +17,12 @@ module numerics
|
||||||
private
|
private
|
||||||
|
|
||||||
integer, protected, public :: &
|
integer, protected, public :: &
|
||||||
iJacoStiffness = 1, & !< frequency of stiffness update
|
iJacoStiffness = 1, & !< frequency of stiffness update
|
||||||
nMPstate = 10, & !< materialpoint state loop limit
|
nMPstate = 10, & !< materialpoint state loop limit
|
||||||
randomSeed = 0, & !< fixed seeding for pseudo-random number generator, Default 0: use random seed
|
randomSeed = 0, & !< fixed seeding for pseudo-random number generator, Default 0: use random seed
|
||||||
worldrank = 0, & !< MPI worldrank (/=0 for MPI simulations only)
|
worldrank = 0, & !< MPI worldrank (/=0 for MPI simulations only)
|
||||||
worldsize = 1, & !< MPI worldsize (/=1 for MPI simulations only)
|
worldsize = 1, & !< MPI worldsize (/=1 for MPI simulations only)
|
||||||
numerics_integrator = 1 !< method used for state integration Default 1: fix-point iteration
|
numerics_integrator = 1 !< method used for state integration Default 1: fix-point iteration
|
||||||
integer(4), protected, public :: &
|
integer(4), protected, public :: &
|
||||||
DAMASK_NumThreadsInt = 0 !< value stored in environment variable DAMASK_NUM_THREADS, set to zero if no OpenMP directive
|
DAMASK_NumThreadsInt = 0 !< value stored in environment variable DAMASK_NUM_THREADS, set to zero if no OpenMP directive
|
||||||
real(pReal), protected, public :: &
|
real(pReal), protected, public :: &
|
||||||
|
@ -31,19 +31,6 @@ module numerics
|
||||||
subStepSizeHomog = 0.25_pReal, & !< size of first substep when cutback in homogenization
|
subStepSizeHomog = 0.25_pReal, & !< size of first substep when cutback in homogenization
|
||||||
stepIncreaseHomog = 1.5_pReal, & !< increase of next substep size when previous substep converged in homogenization
|
stepIncreaseHomog = 1.5_pReal, & !< increase of next substep size when previous substep converged in homogenization
|
||||||
numerics_unitlength = 1.0_pReal, & !< determines the physical length of one computational length unit
|
numerics_unitlength = 1.0_pReal, & !< determines the physical length of one computational length unit
|
||||||
absTol_RGC = 1.0e+4_pReal, & !< absolute tolerance of RGC residuum
|
|
||||||
relTol_RGC = 1.0e-3_pReal, & !< relative tolerance of RGC residuum
|
|
||||||
absMax_RGC = 1.0e+10_pReal, & !< absolute maximum of RGC residuum
|
|
||||||
relMax_RGC = 1.0e+2_pReal, & !< relative maximum of RGC residuum
|
|
||||||
pPert_RGC = 1.0e-7_pReal, & !< perturbation for computing RGC penalty tangent
|
|
||||||
xSmoo_RGC = 1.0e-5_pReal, & !< RGC penalty smoothing parameter (hyperbolic tangent)
|
|
||||||
viscPower_RGC = 1.0e+0_pReal, & !< power (sensitivity rate) of numerical viscosity in RGC scheme, Default 1.0e0: Newton viscosity (linear model)
|
|
||||||
viscModus_RGC = 0.0e+0_pReal, & !< stress modulus of RGC numerical viscosity, Default 0.0e0: No viscosity is applied
|
|
||||||
refRelaxRate_RGC = 1.0e-3_pReal, & !< reference relaxation rate in RGC viscosity
|
|
||||||
maxdRelax_RGC = 1.0e+0_pReal, & !< threshold of maximum relaxation vector increment (if exceed this then cutback)
|
|
||||||
maxVolDiscr_RGC = 1.0e-5_pReal, & !< threshold of maximum volume discrepancy allowed
|
|
||||||
volDiscrMod_RGC = 1.0e+12_pReal, & !< stiffness of RGC volume discrepancy (zero = without volume discrepancy constraint)
|
|
||||||
volDiscrPow_RGC = 5.0_pReal, & !< powerlaw penalty for volume discrepancy
|
|
||||||
charLength = 1.0_pReal, & !< characteristic length scale for gradient problems
|
charLength = 1.0_pReal, & !< characteristic length scale for gradient problems
|
||||||
residualStiffness = 1.0e-6_pReal !< non-zero residual damage
|
residualStiffness = 1.0e-6_pReal !< non-zero residual damage
|
||||||
logical, protected, public :: &
|
logical, protected, public :: &
|
||||||
|
@ -179,35 +166,6 @@ subroutine numerics_init
|
||||||
case ('unitlength')
|
case ('unitlength')
|
||||||
numerics_unitlength = IO_floatValue(line,chunkPos,2)
|
numerics_unitlength = IO_floatValue(line,chunkPos,2)
|
||||||
|
|
||||||
!--------------------------------------------------------------------------------------------------
|
|
||||||
! RGC parameters
|
|
||||||
case ('atol_rgc')
|
|
||||||
absTol_RGC = IO_floatValue(line,chunkPos,2)
|
|
||||||
case ('rtol_rgc')
|
|
||||||
relTol_RGC = IO_floatValue(line,chunkPos,2)
|
|
||||||
case ('amax_rgc')
|
|
||||||
absMax_RGC = IO_floatValue(line,chunkPos,2)
|
|
||||||
case ('rmax_rgc')
|
|
||||||
relMax_RGC = IO_floatValue(line,chunkPos,2)
|
|
||||||
case ('perturbpenalty_rgc')
|
|
||||||
pPert_RGC = IO_floatValue(line,chunkPos,2)
|
|
||||||
case ('relevantmismatch_rgc')
|
|
||||||
xSmoo_RGC = IO_floatValue(line,chunkPos,2)
|
|
||||||
case ('viscositypower_rgc')
|
|
||||||
viscPower_RGC = IO_floatValue(line,chunkPos,2)
|
|
||||||
case ('viscositymodulus_rgc')
|
|
||||||
viscModus_RGC = IO_floatValue(line,chunkPos,2)
|
|
||||||
case ('refrelaxationrate_rgc')
|
|
||||||
refRelaxRate_RGC = IO_floatValue(line,chunkPos,2)
|
|
||||||
case ('maxrelaxation_rgc')
|
|
||||||
maxdRelax_RGC = IO_floatValue(line,chunkPos,2)
|
|
||||||
case ('maxvoldiscrepancy_rgc')
|
|
||||||
maxVolDiscr_RGC = IO_floatValue(line,chunkPos,2)
|
|
||||||
case ('voldiscrepancymod_rgc')
|
|
||||||
volDiscrMod_RGC = IO_floatValue(line,chunkPos,2)
|
|
||||||
case ('discrepancypower_rgc')
|
|
||||||
volDiscrPow_RGC = IO_floatValue(line,chunkPos,2)
|
|
||||||
|
|
||||||
!--------------------------------------------------------------------------------------------------
|
!--------------------------------------------------------------------------------------------------
|
||||||
! random seeding parameter
|
! random seeding parameter
|
||||||
case ('random_seed','fixed_seed')
|
case ('random_seed','fixed_seed')
|
||||||
|
@ -280,8 +238,6 @@ subroutine numerics_init
|
||||||
#endif
|
#endif
|
||||||
end select
|
end select
|
||||||
enddo
|
enddo
|
||||||
|
|
||||||
|
|
||||||
else fileExists
|
else fileExists
|
||||||
write(6,'(a,/)') ' using standard values'
|
write(6,'(a,/)') ' using standard values'
|
||||||
flush(6)
|
flush(6)
|
||||||
|
@ -301,21 +257,6 @@ subroutine numerics_init
|
||||||
write(6,'(a24,1x,es8.1)') ' stepIncreaseHomog: ',stepIncreaseHomog
|
write(6,'(a24,1x,es8.1)') ' stepIncreaseHomog: ',stepIncreaseHomog
|
||||||
write(6,'(a24,1x,i8,/)') ' nMPstate: ',nMPstate
|
write(6,'(a24,1x,i8,/)') ' nMPstate: ',nMPstate
|
||||||
|
|
||||||
!--------------------------------------------------------------------------------------------------
|
|
||||||
! RGC parameters
|
|
||||||
write(6,'(a24,1x,es8.1)') ' aTol_RGC: ',absTol_RGC
|
|
||||||
write(6,'(a24,1x,es8.1)') ' rTol_RGC: ',relTol_RGC
|
|
||||||
write(6,'(a24,1x,es8.1)') ' aMax_RGC: ',absMax_RGC
|
|
||||||
write(6,'(a24,1x,es8.1)') ' rMax_RGC: ',relMax_RGC
|
|
||||||
write(6,'(a24,1x,es8.1)') ' perturbPenalty_RGC: ',pPert_RGC
|
|
||||||
write(6,'(a24,1x,es8.1)') ' relevantMismatch_RGC: ',xSmoo_RGC
|
|
||||||
write(6,'(a24,1x,es8.1)') ' viscosityrate_RGC: ',viscPower_RGC
|
|
||||||
write(6,'(a24,1x,es8.1)') ' viscositymodulus_RGC: ',viscModus_RGC
|
|
||||||
write(6,'(a24,1x,es8.1)') ' maxrelaxation_RGC: ',maxdRelax_RGC
|
|
||||||
write(6,'(a24,1x,es8.1)') ' maxVolDiscrepancy_RGC: ',maxVolDiscr_RGC
|
|
||||||
write(6,'(a24,1x,es8.1)') ' volDiscrepancyMod_RGC: ',volDiscrMod_RGC
|
|
||||||
write(6,'(a24,1x,es8.1,/)') ' discrepancyPower_RGC: ',volDiscrPow_RGC
|
|
||||||
|
|
||||||
!--------------------------------------------------------------------------------------------------
|
!--------------------------------------------------------------------------------------------------
|
||||||
! Random seeding parameter
|
! Random seeding parameter
|
||||||
write(6,'(a16,1x,i16,/)') ' random_seed: ',randomSeed
|
write(6,'(a16,1x,i16,/)') ' random_seed: ',randomSeed
|
||||||
|
@ -367,7 +308,6 @@ subroutine numerics_init
|
||||||
write(6,'(a24,1x,L8)') ' B-Bar stabilisation: ',BBarStabilisation
|
write(6,'(a24,1x,L8)') ' B-Bar stabilisation: ',BBarStabilisation
|
||||||
#endif
|
#endif
|
||||||
|
|
||||||
|
|
||||||
!--------------------------------------------------------------------------------------------------
|
!--------------------------------------------------------------------------------------------------
|
||||||
! sanity checks
|
! sanity checks
|
||||||
if (defgradTolerance <= 0.0_pReal) call IO_error(301,ext_msg='defgradTolerance')
|
if (defgradTolerance <= 0.0_pReal) call IO_error(301,ext_msg='defgradTolerance')
|
||||||
|
@ -379,19 +319,6 @@ subroutine numerics_init
|
||||||
if (numerics_integrator <= 0 .or. numerics_integrator >= 6) &
|
if (numerics_integrator <= 0 .or. numerics_integrator >= 6) &
|
||||||
call IO_error(301,ext_msg='integrator')
|
call IO_error(301,ext_msg='integrator')
|
||||||
if (numerics_unitlength <= 0.0_pReal) call IO_error(301,ext_msg='unitlength')
|
if (numerics_unitlength <= 0.0_pReal) call IO_error(301,ext_msg='unitlength')
|
||||||
if (absTol_RGC <= 0.0_pReal) call IO_error(301,ext_msg='absTol_RGC')
|
|
||||||
if (relTol_RGC <= 0.0_pReal) call IO_error(301,ext_msg='relTol_RGC')
|
|
||||||
if (absMax_RGC <= 0.0_pReal) call IO_error(301,ext_msg='absMax_RGC')
|
|
||||||
if (relMax_RGC <= 0.0_pReal) call IO_error(301,ext_msg='relMax_RGC')
|
|
||||||
if (pPert_RGC <= 0.0_pReal) call IO_error(301,ext_msg='pPert_RGC')
|
|
||||||
if (xSmoo_RGC <= 0.0_pReal) call IO_error(301,ext_msg='xSmoo_RGC')
|
|
||||||
if (viscPower_RGC < 0.0_pReal) call IO_error(301,ext_msg='viscPower_RGC')
|
|
||||||
if (viscModus_RGC < 0.0_pReal) call IO_error(301,ext_msg='viscModus_RGC')
|
|
||||||
if (refRelaxRate_RGC <= 0.0_pReal) call IO_error(301,ext_msg='refRelaxRate_RGC')
|
|
||||||
if (maxdRelax_RGC <= 0.0_pReal) call IO_error(301,ext_msg='maxdRelax_RGC')
|
|
||||||
if (maxVolDiscr_RGC <= 0.0_pReal) call IO_error(301,ext_msg='maxVolDiscr_RGC')
|
|
||||||
if (volDiscrMod_RGC < 0.0_pReal) call IO_error(301,ext_msg='volDiscrMod_RGC')
|
|
||||||
if (volDiscrPow_RGC <= 0.0_pReal) call IO_error(301,ext_msg='volDiscrPw_RGC')
|
|
||||||
if (residualStiffness < 0.0_pReal) call IO_error(301,ext_msg='residualStiffness')
|
if (residualStiffness < 0.0_pReal) call IO_error(301,ext_msg='residualStiffness')
|
||||||
if (itmax <= 1) call IO_error(301,ext_msg='itmax')
|
if (itmax <= 1) call IO_error(301,ext_msg='itmax')
|
||||||
if (itmin > itmax .or. itmin < 1) call IO_error(301,ext_msg='itmin')
|
if (itmin > itmax .or. itmin < 1) call IO_error(301,ext_msg='itmin')
|
||||||
|
|
Loading…
Reference in New Issue