following paper

This commit is contained in:
Martin Diehl 2021-04-08 13:31:21 +02:00
parent d59051f576
commit c4765d3742
5 changed files with 45 additions and 64 deletions

View File

@ -259,7 +259,7 @@ subroutine formResidual(in,x_scal,f_scal,dummy,ierr)
PetscObject :: dummy
PetscErrorCode :: ierr
integer :: i, j, k, ce
real(pReal) :: phiDot, mobility
real(pReal) :: mobility
phi_current = x_scal
!--------------------------------------------------------------------------------------------------
@ -272,7 +272,7 @@ subroutine formResidual(in,x_scal,f_scal,dummy,ierr)
ce = 0
do k = 1, grid3; do j = 1, grid(2); do i = 1,grid(1)
ce = ce + 1
vectorField_real(1:3,i,j,k) = matmul(damage_nonlocal_getDiffusion(ce) - K_ref, &
vectorField_real(1:3,i,j,k) = matmul(homogenization_K_phi(ce) - K_ref, &
vectorField_real(1:3,i,j,k))
enddo; enddo; enddo
call utilities_FFTvectorForward
@ -281,9 +281,8 @@ subroutine formResidual(in,x_scal,f_scal,dummy,ierr)
ce = 0
do k = 1, grid3; do j = 1, grid(2); do i = 1,grid(1)
ce = ce + 1
call damage_nonlocal_getSourceAndItsTangent(phiDot, phi_current(i,j,k),ce)
mobility = damage_nonlocal_getMobility(ce)
scalarField_real(i,j,k) = params%timeinc*(scalarField_real(i,j,k) + phiDot) &
mobility = homogenization_mu_phi(ce)
scalarField_real(i,j,k) = params%timeinc*(scalarField_real(i,j,k) + homogenization_f_phi(phi_current(i,j,k),ce)) &
+ mobility*(phi_lastInc(i,j,k) - phi_current(i,j,k)) &
+ mu_ref*phi_current(i,j,k)
enddo; enddo; enddo
@ -317,8 +316,8 @@ subroutine updateReference
mu_ref = 0.0_pReal
do k = 1, grid3; do j = 1, grid(2); do i = 1,grid(1)
ce = ce + 1
K_ref = K_ref + damage_nonlocal_getDiffusion(ce)
mu_ref = mu_ref + damage_nonlocal_getMobility(ce)
K_ref = K_ref + homogenization_K_phi(ce)
mu_ref = mu_ref + homogenization_mu_phi(ce)
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

@ -161,18 +161,16 @@ module homogenization
real(pReal) :: f_T
end function homogenization_f_T
module function damage_nonlocal_getMobility(ce) result(M)
module function homogenization_mu_phi(ce) result(M)
integer, intent(in) :: ce
real(pReal) :: M
end function damage_nonlocal_getMobility
end function homogenization_mu_phi
module subroutine damage_nonlocal_getSourceAndItsTangent(phiDot, phi, ce)
module function homogenization_f_phi(phi,ce) result(f_phi)
integer, intent(in) :: ce
real(pReal), intent(in) :: &
phi
real(pReal), intent(out) :: &
phiDot
end subroutine damage_nonlocal_getSourceAndItsTangent
real(pReal), intent(in) :: phi
real(pReal) :: f_phi
end function homogenization_f_phi
module subroutine homogenization_set_phi(phi,ce)
integer, intent(in) :: ce
@ -188,8 +186,8 @@ module homogenization
homogenization_thermal_mu_T, &
homogenization_K, &
homogenization_f_T, &
damage_nonlocal_getMobility, &
damage_nonlocal_getSourceAndItsTangent, &
homogenization_mu_phi, &
homogenization_f_phi, &
homogenization_set_phi, &
homogenization_thermal_setfield, &
homogenization_T, &
@ -201,7 +199,7 @@ module homogenization
DAMAGE_NONLOCAL_ID
public :: &
damage_nonlocal_getDiffusion
homogenization_K_phi
contains
@ -448,27 +446,27 @@ end subroutine homogenization_restartRead
!--------------------------------------------------------------------------------------------------
!> @brief returns homogenized non local damage diffusion tensor in reference configuration
!--------------------------------------------------------------------------------------------------
function damage_nonlocal_getDiffusion(ce)
function homogenization_K_phi(ce)
integer, intent(in) :: ce
real(pReal), dimension(3,3) :: &
damage_nonlocal_getDiffusion
homogenization_K_phi
integer :: &
ho, &
co
ho = material_homogenizationID(ce)
damage_nonlocal_getDiffusion = 0.0_pReal
homogenization_K_phi = 0.0_pReal
do co = 1, homogenization_Nconstituents(ho)
damage_nonlocal_getDiffusion = damage_nonlocal_getDiffusion + &
homogenization_K_phi = homogenization_K_phi + &
crystallite_push33ToRef(co,ce,lattice_D(1:3,1:3,material_phaseID(co,ce)))
enddo
damage_nonlocal_getDiffusion = &
num_damage%charLength**2*damage_nonlocal_getDiffusion/real(homogenization_Nconstituents(ho),pReal)
homogenization_K_phi = &
num_damage%charLength**2*homogenization_K_phi/real(homogenization_Nconstituents(ho),pReal)
end function damage_nonlocal_getDiffusion
end function homogenization_K_phi
!--------------------------------------------------------------------------------------------------

View File

@ -100,7 +100,7 @@ module subroutine damage_partition(ce)
if(damageState_h(material_homogenizationID(ce))%sizeState < 1) return
phi = damagestate_h(material_homogenizationID(ce))%state(1,material_homogenizationEntry(ce))
call phase_damage_set_phi(phi,1,ce)
call phase_set_phi(phi,1,ce)
end subroutine damage_partition
@ -109,30 +109,29 @@ end subroutine damage_partition
!--------------------------------------------------------------------------------------------------
!> @brief Returns homogenized nonlocal damage mobility
!--------------------------------------------------------------------------------------------------
module function damage_nonlocal_getMobility(ce) result(M)
module function homogenization_mu_phi(ce) result(M)
integer, intent(in) :: ce
real(pReal) :: M
M = lattice_M(material_phaseID(1,ce))
end function damage_nonlocal_getMobility
end function homogenization_mu_phi
!--------------------------------------------------------------------------------------------------
!> @brief calculates homogenized damage driving forces
!--------------------------------------------------------------------------------------------------
module subroutine damage_nonlocal_getSourceAndItsTangent(phiDot, phi, ce)
module function homogenization_f_phi(phi,ce) result(f_phi)
integer, intent(in) :: ce
real(pReal), intent(in) :: &
phi
real(pReal), intent(out) :: &
phiDot
real(pReal) :: f_phi
phiDot = phase_damage_phi_dot(phi, 1, ce)
f_phi = phase_f_phi(phi, 1, ce)
end subroutine damage_nonlocal_getSourceAndItsTangent
end function homogenization_f_phi
!--------------------------------------------------------------------------------------------------

View File

@ -145,26 +145,22 @@ module phase
real(pReal), dimension(3,3) :: L_p
end function mechanical_L_p
module function phase_F(co,ce) result(F)
integer, intent(in) :: co, ce
real(pReal), dimension(3,3) :: F
end function phase_F
module function mechanical_F_e(ph,me) result(F_e)
integer, intent(in) :: ph,me
real(pReal), dimension(3,3) :: F_e
end function mechanical_F_e
module function phase_F(co,ce) result(F)
integer, intent(in) :: co, ce
real(pReal), dimension(3,3) :: F
end function phase_F
module function phase_P(co,ce) result(P)
integer, intent(in) :: co, ce
real(pReal), dimension(3,3) :: P
end function phase_P
module function phase_damage_get_phi(co,ip,el) result(phi)
integer, intent(in) :: co, ip, el
real(pReal) :: phi
end function phase_damage_get_phi
module function thermal_T(ph,me) result(T)
integer, intent(in) :: ph,me
real(pReal) :: T
@ -191,10 +187,10 @@ module phase
integer, intent(in) :: ce, co
end subroutine phase_thermal_setField
module subroutine phase_damage_set_phi(phi,co,ce)
module subroutine phase_set_phi(phi,co,ce)
real(pReal), intent(in) :: phi
integer, intent(in) :: co, ce
end subroutine phase_damage_set_phi
end subroutine phase_set_phi
! == cleaned:end ===================================================================================
@ -227,13 +223,13 @@ module phase
end function phase_homogenizedC
module function phase_damage_phi_dot(phi,co,ce) result(phi_dot)
module function phase_f_phi(phi,co,ce) result(phi_dot)
integer, intent(in) :: ce,co
real(pReal), intent(in) :: &
phi !< damage parameter
real(pReal) :: &
phi_dot
end function phase_damage_phi_dot
end function phase_f_phi
module function phase_f_T(ph,me) result(f_T)
integer, intent(in) :: ph, me
@ -299,7 +295,7 @@ module phase
public :: &
phase_init, &
phase_homogenizedC, &
phase_damage_phi_dot, &
phase_f_phi, &
phase_f_T, &
phase_results, &
phase_allocateState, &
@ -317,8 +313,7 @@ module phase
phase_restartRead, &
integrateDamageState, &
phase_thermal_setField, &
phase_damage_set_phi, &
phase_damage_get_phi, &
phase_set_phi, &
phase_P, &
phase_set_F, &
phase_F

View File

@ -142,7 +142,7 @@ end subroutine damage_init
!----------------------------------------------------------------------------------------------
!< @brief returns local part of nonlocal damage driving force
!----------------------------------------------------------------------------------------------
module function phase_damage_phi_dot(phi,co,ce) result(phi_dot)
module function phase_f_phi(phi,co,ce) result(phi_dot)
integer, intent(in) :: ce,co
real(pReal), intent(in) :: &
@ -165,7 +165,7 @@ module function phase_damage_phi_dot(phi,co,ce) result(phi_dot)
phi_dot = 0.0_pReal
end select
end function phase_damage_phi_dot
end function phase_f_phi
@ -411,7 +411,7 @@ end function source_active
!----------------------------------------------------------------------------------------------
!< @brief Set damage parameter
!----------------------------------------------------------------------------------------------
module subroutine phase_damage_set_phi(phi,co,ce)
module subroutine phase_set_phi(phi,co,ce)
real(pReal), intent(in) :: phi
integer, intent(in) :: ce, co
@ -419,17 +419,7 @@ module subroutine phase_damage_set_phi(phi,co,ce)
current(material_phaseID(co,ce))%phi(material_phaseEntry(co,ce)) = phi
end subroutine phase_damage_set_phi
module function phase_damage_get_phi(co,ip,el) result(phi)
integer, intent(in) :: co, ip, el
real(pReal) :: phi
phi = current(material_phaseAt(co,el))%phi(material_phaseMemberAt(co,ip,el))
end function phase_damage_get_phi
end subroutine phase_set_phi
module function damage_phi(ph,me) result(phi)