less public variables, clearer names

This commit is contained in:
Martin Diehl 2020-02-29 12:57:19 +01:00
parent ca76014e45
commit 2bc36121b2
6 changed files with 70 additions and 87 deletions

View File

@ -29,7 +29,7 @@ module damage_nonlocal
public :: &
damage_nonlocal_init, &
damage_nonlocal_getSourceAndItsTangent, &
damage_nonlocal_getDiffusion33, &
damage_nonlocal_getDiffusion, &
damage_nonlocal_getMobility, &
damage_nonlocal_putNonLocalDamage, &
damage_nonlocal_results
@ -130,28 +130,28 @@ end subroutine damage_nonlocal_getSourceAndItsTangent
!--------------------------------------------------------------------------------------------------
!> @brief returns homogenized non local damage diffusion tensor in reference configuration
!--------------------------------------------------------------------------------------------------
function damage_nonlocal_getDiffusion33(ip,el)
function damage_nonlocal_getDiffusion(ip,el)
integer, intent(in) :: &
ip, & !< integration point number
el !< element number
real(pReal), dimension(3,3) :: &
damage_nonlocal_getDiffusion33
damage_nonlocal_getDiffusion
integer :: &
homog, &
grain
homog = material_homogenizationAt(el)
damage_nonlocal_getDiffusion33 = 0.0_pReal
damage_nonlocal_getDiffusion = 0.0_pReal
do grain = 1, homogenization_Ngrains(homog)
damage_nonlocal_getDiffusion33 = damage_nonlocal_getDiffusion33 + &
crystallite_push33ToRef(grain,ip,el,lattice_DamageDiffusion33(1:3,1:3,material_phaseAt(grain,el)))
damage_nonlocal_getDiffusion = damage_nonlocal_getDiffusion + &
crystallite_push33ToRef(grain,ip,el,lattice_DamageDiffusion(1:3,1:3,material_phaseAt(grain,el)))
enddo
damage_nonlocal_getDiffusion33 = &
charLength**2*damage_nonlocal_getDiffusion33/real(homogenization_Ngrains(homog),pReal)
damage_nonlocal_getDiffusion = &
charLength**2*damage_nonlocal_getDiffusion/real(homogenization_Ngrains(homog),pReal)
end function damage_nonlocal_getDiffusion33
end function damage_nonlocal_getDiffusion
!--------------------------------------------------------------------------------------------------

View File

@ -249,8 +249,8 @@ subroutine formResidual(in,x_scal,f_scal,dummy,ierr)
cell = 0
do k = 1, grid3; do j = 1, grid(2); do i = 1,grid(1)
cell = cell + 1
vectorField_real(1:3,i,j,k) = matmul(damage_nonlocal_getDiffusion33(1,cell) - K_ref, &
vectorField_real(1:3,i,j,k))
vectorField_real(1:3,i,j,k) = matmul(damage_nonlocal_getDiffusion(1,cell) - K_ref, &
vectorField_real(1:3,i,j,k))
enddo; enddo; enddo
call utilities_FFTvectorForward
call utilities_fourierVectorDivergence !< calculate damage divergence in fourier field
@ -294,7 +294,7 @@ subroutine updateReference
mu_ref = 0.0_pReal
do k = 1, grid3; do j = 1, grid(2); do i = 1,grid(1)
cell = cell + 1
K_ref = K_ref + damage_nonlocal_getDiffusion33(1,cell)
K_ref = K_ref + damage_nonlocal_getDiffusion(1,cell)
mu_ref = mu_ref + damage_nonlocal_getMobility(1,cell)
enddo; enddo; enddo
K_ref = K_ref*wgt

View File

@ -252,8 +252,8 @@ subroutine formResidual(in,x_scal,f_scal,dummy,ierr)
cell = 0
do k = 1, grid3; do j = 1, grid(2); do i = 1,grid(1)
cell = cell + 1
vectorField_real(1:3,i,j,k) = matmul(thermal_conduction_getConductivity33(1,cell) - K_ref, &
vectorField_real(1:3,i,j,k))
vectorField_real(1:3,i,j,k) = matmul(thermal_conduction_getConductivity(1,cell) - K_ref, &
vectorField_real(1:3,i,j,k))
enddo; enddo; enddo
call utilities_FFTvectorForward
call utilities_fourierVectorDivergence !< calculate temperature divergence in fourier field
@ -294,9 +294,8 @@ subroutine updateReference
mu_ref = 0.0_pReal
do k = 1, grid3; do j = 1, grid(2); do i = 1,grid(1)
cell = cell + 1
K_ref = K_ref + thermal_conduction_getConductivity33(1,cell)
mu_ref = mu_ref + thermal_conduction_getMassDensity(1,cell)* &
thermal_conduction_getSpecificHeat(1,cell)
K_ref = K_ref + thermal_conduction_getConductivity(1,cell)
mu_ref = mu_ref + thermal_conduction_getMassDensity(1,cell)* thermal_conduction_getSpecificHeat(1,cell)
enddo; enddo; enddo
K_ref = K_ref*wgt
call MPI_Allreduce(MPI_IN_PLACE,K_ref,9,MPI_DOUBLE,MPI_SUM,PETSC_COMM_WORLD,ierr)

View File

@ -69,6 +69,9 @@ subroutine kinematics_thermal_expansion_init
prm%expansion(2,2,1:size(temp)) = temp
temp = config%getFloats('thermal_expansion33',defaultVal=[(0.0_pReal, i=1,size(temp))],requiredSize=size(temp))
prm%expansion(3,3,1:size(temp)) = temp
do i=1, size(prm%expansion,3)
prm%expansion(1:3,1:3,i) = lattice_symmetrize33(prm%expansion(1:3,1:3,i),config%getString('lattice_structure'))
enddo
end associate
enddo
@ -90,12 +93,9 @@ pure function kinematics_thermal_expansion_initialStrain(homog,phase,offset)
associate(prm => param(kinematics_thermal_expansion_instance(phase)))
kinematics_thermal_expansion_initialStrain = &
(temperature(homog)%p(offset) - prm%T_ref)**1 / 1. * &
lattice_thermalExpansion33(1:3,1:3,1,phase) + & ! constant coefficient
(temperature(homog)%p(offset) - prm%T_ref)**2 / 2. * &
lattice_thermalExpansion33(1:3,1:3,2,phase) + & ! linear coefficient
(temperature(homog)%p(offset) - prm%T_ref)**3 / 3. * &
lattice_thermalExpansion33(1:3,1:3,3,phase) ! quadratic coefficient
(temperature(homog)%p(offset) - prm%T_ref)**1 / 1. * prm%expansion(1:3,1:3,1) + & ! constant coefficient
(temperature(homog)%p(offset) - prm%T_ref)**2 / 2. * prm%expansion(1:3,1:3,2) + & ! linear coefficient
(temperature(homog)%p(offset) - prm%T_ref)**3 / 3. * prm%expansion(1:3,1:3,3) ! quadratic coefficient
end associate
end function kinematics_thermal_expansion_initialStrain
@ -128,14 +128,14 @@ subroutine kinematics_thermal_expansion_LiAndItsTangent(Li, dLi_dTstar, ipc, ip,
associate(prm => param(kinematics_thermal_expansion_instance(phase)))
Li = TDot * ( &
lattice_thermalExpansion33(1:3,1:3,1,phase)*(T - prm%T_ref)**0 & ! constant coefficient
+ lattice_thermalExpansion33(1:3,1:3,2,phase)*(T - prm%T_ref)**1 & ! linear coefficient
+ lattice_thermalExpansion33(1:3,1:3,3,phase)*(T - prm%T_ref)**2 & ! quadratic coefficient
prm%expansion(1:3,1:3,1)*(T - prm%T_ref)**0 & ! constant coefficient
+ prm%expansion(1:3,1:3,2)*(T - prm%T_ref)**1 & ! linear coefficient
+ prm%expansion(1:3,1:3,3)*(T - prm%T_ref)**2 & ! quadratic coefficient
) / &
(1.0_pReal &
+ lattice_thermalExpansion33(1:3,1:3,1,phase)*(T - prm%T_ref)**1 / 1. &
+ lattice_thermalExpansion33(1:3,1:3,2,phase)*(T - prm%T_ref)**2 / 2. &
+ lattice_thermalExpansion33(1:3,1:3,3,phase)*(T - prm%T_ref)**3 / 3. &
+ prm%expansion(1:3,1:3,1)*(T - prm%T_ref)**1 / 1. &
+ prm%expansion(1:3,1:3,2)*(T - prm%T_ref)**2 / 2. &
+ prm%expansion(1:3,1:3,3)*(T - prm%T_ref)**3 / 3. &
)
end associate
dLi_dTstar = 0.0_pReal

View File

@ -387,11 +387,9 @@ module lattice
! SHOULD NOT BE PART OF LATTICE BEGIN
real(pReal), dimension(:), allocatable, public, protected :: &
lattice_mu, lattice_nu
real(pReal), dimension(:,:,:,:), allocatable, public, protected :: & ! with higher-order parameters (e.g. temperature-dependent)
lattice_thermalExpansion33
real(pReal), dimension(:,:,:), allocatable, public, protected :: &
lattice_thermalConductivity33, &
lattice_damageDiffusion33
lattice_thermalConductivity, &
lattice_damageDiffusion
real(pReal), dimension(:), allocatable, public, protected :: &
lattice_damageMobility, &
lattice_massDensity, &
@ -399,13 +397,13 @@ module lattice
! SHOULD NOT BE PART OF LATTICE END
enum, bind(c)
enumerator :: LATTICE_undefined_ID, &
LATTICE_iso_ID, &
LATTICE_fcc_ID, &
LATTICE_bcc_ID, &
LATTICE_hex_ID, &
LATTICE_bct_ID, &
LATTICE_ort_ID
enumerator :: LATTICE_UNDEFINED_ID, &
LATTICE_ISO_ID, &
LATTICE_FCC_ID, &
LATTICE_BCC_ID, &
LATTICE_BCT_ID, &
LATTICE_HEX_ID, &
LATTICE_ORT_ID
end enum
integer(kind(LATTICE_undefined_ID)), dimension(:), allocatable, public, protected :: &
@ -421,12 +419,13 @@ module lattice
public :: &
lattice_init, &
LATTICE_iso_ID, &
LATTICE_fcc_ID, &
LATTICE_bcc_ID, &
LATTICE_bct_ID, &
LATTICE_hex_ID, &
LATTICE_ort_ID, &
LATTICE_ISO_ID, &
LATTICE_FCC_ID, &
LATTICE_BCC_ID, &
LATTICE_BCT_ID, &
LATTICE_HEX_ID, &
LATTICE_ORT_ID, &
lattice_symmetrize33, &
lattice_SchmidMatrix_slip, &
lattice_SchmidMatrix_twin, &
lattice_SchmidMatrix_trans, &
@ -458,8 +457,6 @@ subroutine lattice_init
integer :: Nphases, p,i
character(len=pStringLen) :: structure = ''
real(pReal), dimension(:), allocatable :: &
temp
write(6,'(/,a)') ' <<<+- lattice init -+>>>'; flush(6)
@ -468,12 +465,11 @@ subroutine lattice_init
allocate(lattice_structure(Nphases),source = LATTICE_undefined_ID)
allocate(lattice_C66(6,6,Nphases), source=0.0_pReal)
allocate(lattice_thermalExpansion33 (3,3,3,Nphases), source=0.0_pReal) ! constant, linear, quadratic coefficients
allocate(lattice_thermalConductivity33 (3,3,Nphases), source=0.0_pReal)
allocate(lattice_damageDiffusion33 (3,3,Nphases), source=0.0_pReal)
allocate(lattice_damageMobility ( Nphases), source=0.0_pReal)
allocate(lattice_massDensity ( Nphases), source=0.0_pReal)
allocate(lattice_specificHeat ( Nphases), source=0.0_pReal)
allocate(lattice_thermalConductivity (3,3,Nphases), source=0.0_pReal)
allocate(lattice_damageDiffusion (3,3,Nphases), source=0.0_pReal)
allocate(lattice_damageMobility ( Nphases), source=0.0_pReal)
allocate(lattice_massDensity ( Nphases), source=0.0_pReal)
allocate(lattice_specificHeat ( Nphases), source=0.0_pReal)
allocate(lattice_mu(Nphases), source=0.0_pReal)
allocate(lattice_nu(Nphases), source=0.0_pReal)
@ -523,33 +519,21 @@ subroutine lattice_init
! should not be part of lattice
lattice_thermalConductivity33(1,1,p) = config_phase(p)%getFloat('thermal_conductivity11',defaultVal=0.0_pReal)
lattice_thermalConductivity33(2,2,p) = config_phase(p)%getFloat('thermal_conductivity22',defaultVal=0.0_pReal)
lattice_thermalConductivity33(3,3,p) = config_phase(p)%getFloat('thermal_conductivity33',defaultVal=0.0_pReal)
lattice_thermalConductivity(1,1,p) = config_phase(p)%getFloat('thermal_conductivity11',defaultVal=0.0_pReal)
lattice_thermalConductivity(2,2,p) = config_phase(p)%getFloat('thermal_conductivity22',defaultVal=0.0_pReal)
lattice_thermalConductivity(3,3,p) = config_phase(p)%getFloat('thermal_conductivity33',defaultVal=0.0_pReal)
lattice_thermalConductivity(1:3,1:3,p) = lattice_symmetrize33(lattice_thermalConductivity(1:3,1:3,p),structure)
temp = config_phase(p)%getFloats('thermal_expansion11',defaultVal=[0.0_pReal]) ! read up to three parameters (constant, linear, quadratic with T)
lattice_thermalExpansion33(1,1,1:size(temp),p) = temp
temp = config_phase(p)%getFloats('thermal_expansion22',defaultVal=[0.0_pReal]) ! read up to three parameters (constant, linear, quadratic with T)
lattice_thermalExpansion33(2,2,1:size(temp),p) = temp
temp = config_phase(p)%getFloats('thermal_expansion33',defaultVal=[0.0_pReal]) ! read up to three parameters (constant, linear, quadratic with T)
lattice_thermalExpansion33(3,3,1:size(temp),p) = temp
lattice_specificHeat(p) = config_phase(p)%getFloat('specific_heat',defaultVal=0.0_pReal)
lattice_massDensity(p) = config_phase(p)%getFloat('mass_density', defaultVal=0.0_pReal)
lattice_specificHeat(p) = config_phase(p)%getFloat( 'specific_heat',defaultVal=0.0_pReal)
lattice_massDensity(p) = config_phase(p)%getFloat( 'mass_density',defaultVal=0.0_pReal)
lattice_DamageDiffusion33(1,1,p) = config_phase(p)%getFloat( 'damage_diffusion11',defaultVal=0.0_pReal)
lattice_DamageDiffusion33(2,2,p) = config_phase(p)%getFloat( 'damage_diffusion22',defaultVal=0.0_pReal)
lattice_DamageDiffusion33(3,3,p) = config_phase(p)%getFloat( 'damage_diffusion33',defaultVal=0.0_pReal)
lattice_DamageDiffusion(1,1,p) = config_phase(p)%getFloat('damage_diffusion11',defaultVal=0.0_pReal)
lattice_DamageDiffusion(2,2,p) = config_phase(p)%getFloat('damage_diffusion22',defaultVal=0.0_pReal)
lattice_DamageDiffusion(3,3,p) = config_phase(p)%getFloat('damage_diffusion33',defaultVal=0.0_pReal)
lattice_DamageDiffusion(1:3,1:3,p) = lattice_symmetrize33(lattice_DamageDiffusion(1:3,1:3,p),structure)
lattice_DamageMobility(p) = config_phase(p)%getFloat( 'damage_mobility',defaultVal=0.0_pReal)
do i = 1,3
lattice_thermalExpansion33 (1:3,1:3,i,p) = lattice_symmetrize33(lattice_thermalExpansion33 (1:3,1:3,i,p), &
structure)
enddo
lattice_thermalConductivity33 (1:3,1:3,p) = lattice_symmetrize33(lattice_thermalConductivity33(1:3,1:3,p), &
structure)
lattice_DamageDiffusion33 (1:3,1:3,p) = lattice_symmetrize33(lattice_DamageDiffusion33 (1:3,1:3,p), &
structure)
enddo
end subroutine lattice_init

View File

@ -26,7 +26,7 @@ module thermal_conduction
public :: &
thermal_conduction_init, &
thermal_conduction_getSourceAndItsTangent, &
thermal_conduction_getConductivity33, &
thermal_conduction_getConductivity, &
thermal_conduction_getSpecificHeat, &
thermal_conduction_getMassDensity, &
thermal_conduction_putTemperatureAndItsRate, &
@ -137,27 +137,27 @@ end subroutine thermal_conduction_getSourceAndItsTangent
!--------------------------------------------------------------------------------------------------
!> @brief returns homogenized thermal conductivity in reference configuration
!--------------------------------------------------------------------------------------------------
function thermal_conduction_getConductivity33(ip,el)
function thermal_conduction_getConductivity(ip,el)
integer, intent(in) :: &
ip, & !< integration point number
el !< element number
real(pReal), dimension(3,3) :: &
thermal_conduction_getConductivity33
thermal_conduction_getConductivity
integer :: &
grain
thermal_conduction_getConductivity33 = 0.0_pReal
thermal_conduction_getConductivity = 0.0_pReal
do grain = 1, homogenization_Ngrains(material_homogenizationAt(el))
thermal_conduction_getConductivity33 = thermal_conduction_getConductivity33 + &
crystallite_push33ToRef(grain,ip,el,lattice_thermalConductivity33(:,:,material_phaseAt(grain,el)))
thermal_conduction_getConductivity = thermal_conduction_getConductivity + &
crystallite_push33ToRef(grain,ip,el,lattice_thermalConductivity(:,:,material_phaseAt(grain,el)))
enddo
thermal_conduction_getConductivity33 = thermal_conduction_getConductivity33 &
/ real(homogenization_Ngrains(material_homogenizationAt(el)),pReal)
thermal_conduction_getConductivity = thermal_conduction_getConductivity &
/ real(homogenization_Ngrains(material_homogenizationAt(el)),pReal)
end function thermal_conduction_getConductivity33
end function thermal_conduction_getConductivity
!--------------------------------------------------------------------------------------------------