functions only used within updatestate

This commit is contained in:
Martin Diehl 2018-11-03 21:11:43 +01:00
parent c16fdec749
commit 0aa21e507a
1 changed files with 322 additions and 323 deletions

View File

@ -74,11 +74,6 @@ module homogenization_RGC
homogenization_RGC_updateState, &
homogenization_RGC_postResults
private :: &
stressPenalty, &
volumePenalty, &
grainDeformation, &
surfaceCorrection, &
equivalentModuli, &
relaxationVector, &
interfaceNormal, &
getInterface, &
@ -798,95 +793,7 @@ function homogenization_RGC_updateState(P,F,F0,avgF,dt,dPdF,ip,el)
!$OMP END CRITICAL (write2out)
endif
end function homogenization_RGC_updateState
!--------------------------------------------------------------------------------------------------
!> @brief derive average stress and stiffness from constituent quantities
!--------------------------------------------------------------------------------------------------
subroutine homogenization_RGC_averageStressAndItsTangent(avgP,dAvgPdAvgF,P,dPdF,instance)
use material, only: &
homogenization_maxNgrains
implicit none
real(pReal), dimension (3,3), intent(out) :: avgP !< average stress at material point
real(pReal), dimension (3,3,3,3), intent(out) :: dAvgPdAvgF !< average stiffness at material point
real(pReal), dimension (3,3,homogenization_maxNgrains), intent(in) :: P !< array of current grain stresses
real(pReal), dimension (3,3,3,3,homogenization_maxNgrains), intent(in) :: dPdF !< array of current grain stiffnesses
integer(pInt), intent(in) :: instance
avgP = sum(P,3) /real(product(param(instance)%Nconstituents),pReal)
dAvgPdAvgF = sum(dPdF,5)/real(product(param(instance)%Nconstituents),pReal)
end subroutine homogenization_RGC_averageStressAndItsTangent
!--------------------------------------------------------------------------------------------------
!> @brief return array of homogenization results for post file inclusion
!--------------------------------------------------------------------------------------------------
pure function homogenization_RGC_postResults(ip,el) result(postResults)
use mesh, only: &
mesh_homogenizationAt
use material, only: &
homogenization_typeInstance,&
homogState, &
mappingHomogenization, &
homogenization_Noutput
implicit none
integer(pInt), intent(in) :: &
ip, & !< integration point number
el !< element number
integer(pInt) instance,o,c,nIntFaceTot
type(tParameters) :: prm
real(pReal), dimension(homogenization_RGC_sizePostResults(homogenization_typeInstance(mesh_homogenizationAt(el)))) :: &
postResults
instance = homogenization_typeInstance(mesh_homogenizationAt(el))
associate(prm => param(instance))
nIntFaceTot=(prm%Nconstituents(1)-1_pInt)*prm%Nconstituents(2)* prm%Nconstituents(3)&
+ prm%Nconstituents(1)* (prm%Nconstituents(2)-1_pInt)*prm%Nconstituents(3)&
+ prm%Nconstituents(1)* prm%Nconstituents(2)* (prm%Nconstituents(3)-1_pInt)
c = 0_pInt
postResults = 0.0_pReal
outputsLoop: do o = 1_pInt,size(prm%outputID)
select case(prm%outputID(o))
case (constitutivework_ID)
postResults(c+1) = homogState(mappingHomogenization(2,ip,el))% &
state(3*nIntFaceTot+1,mappingHomogenization(1,ip,el))
c = c + 1_pInt
case (magnitudemismatch_ID)
postResults(c+1) = homogState(mappingHomogenization(2,ip,el))% &
state(3*nIntFaceTot+2,mappingHomogenization(1,ip,el))
postResults(c+2) = homogState(mappingHomogenization(2,ip,el))% &
state(3*nIntFaceTot+3,mappingHomogenization(1,ip,el))
postResults(c+3) = homogState(mappingHomogenization(2,ip,el))% &
state(3*nIntFaceTot+4,mappingHomogenization(1,ip,el))
c = c + 3_pInt
case (penaltyenergy_ID)
postResults(c+1) = homogState(mappingHomogenization(2,ip,el))% &
state(3*nIntFaceTot+5,mappingHomogenization(1,ip,el))
c = c + 1_pInt
case (volumediscrepancy_ID)
postResults(c+1) = homogState(mappingHomogenization(2,ip,el))% &
state(3*nIntFaceTot+6,mappingHomogenization(1,ip,el))
c = c + 1_pInt
case (averagerelaxrate_ID)
postResults(c+1) = homogState(mappingHomogenization(2,ip,el))% &
state(3*nIntFaceTot+7,mappingHomogenization(1,ip,el))
c = c + 1_pInt
case (maximumrelaxrate_ID)
postResults(c+1) = homogState(mappingHomogenization(2,ip,el))% &
state(3*nIntFaceTot+8,mappingHomogenization(1,ip,el))
c = c + 1_pInt
end select
enddo outputsLoop
end associate
end function homogenization_RGC_postResults
contains
!--------------------------------------------------------------------------------------------------
!> @brief calculate stress-like penalty due to deformation mismatch
!--------------------------------------------------------------------------------------------------
@ -1168,6 +1075,138 @@ function equivalentModuli(grainID,ip,el)
end function equivalentModuli
!--------------------------------------------------------------------------------------------------
!> @brief calculating the grain deformation gradient (the same with
! homogenization_RGC_partionDeformation, but used only for perturbation scheme)
!--------------------------------------------------------------------------------------------------
subroutine grainDeformation(F, avgF, ip, el)
use mesh, only: &
mesh_homogenizationAt
use material, only: &
homogenization_maxNgrains,&
homogenization_typeInstance
implicit none
real(pReal), dimension (3,3,homogenization_maxNgrains), intent(out) :: F !< partioned F per grain
real(pReal), dimension (3,3), intent(in) :: avgF !<
integer(pInt), intent(in) :: &
el, & !< element number
ip !< integration point number
real(pReal), dimension (3) :: aVect,nVect
integer(pInt), dimension (4) :: intFace
integer(pInt), dimension (3) :: iGrain3
integer(pInt) :: instance, iGrain,iFace,i,j
!--------------------------------------------------------------------------------------------------
! compute the deformation gradient of individual grains due to relaxations
instance = homogenization_typeInstance(mesh_homogenizationAt(el))
F = 0.0_pReal
do iGrain = 1_pInt,product(param(instance)%Nconstituents)
iGrain3 = grain1to3(iGrain,instance)
do iFace = 1_pInt,6_pInt
intFace = getInterface(iFace,iGrain3)
aVect = relaxationVector(intFace,instance, ip, el)
nVect = interfaceNormal(intFace,ip,el)
forall (i=1_pInt:3_pInt,j=1_pInt:3_pInt) &
F(i,j,iGrain) = F(i,j,iGrain) + aVect(i)*nVect(j) ! effective relaxations
enddo
F(1:3,1:3,iGrain) = F(1:3,1:3,iGrain) + avgF ! relaxed deformation gradient
enddo
end subroutine grainDeformation
end function homogenization_RGC_updateState
!--------------------------------------------------------------------------------------------------
!> @brief derive average stress and stiffness from constituent quantities
!--------------------------------------------------------------------------------------------------
subroutine homogenization_RGC_averageStressAndItsTangent(avgP,dAvgPdAvgF,P,dPdF,instance)
use material, only: &
homogenization_maxNgrains
implicit none
real(pReal), dimension (3,3), intent(out) :: avgP !< average stress at material point
real(pReal), dimension (3,3,3,3), intent(out) :: dAvgPdAvgF !< average stiffness at material point
real(pReal), dimension (3,3,homogenization_maxNgrains), intent(in) :: P !< array of current grain stresses
real(pReal), dimension (3,3,3,3,homogenization_maxNgrains), intent(in) :: dPdF !< array of current grain stiffnesses
integer(pInt), intent(in) :: instance
avgP = sum(P,3) /real(product(param(instance)%Nconstituents),pReal)
dAvgPdAvgF = sum(dPdF,5)/real(product(param(instance)%Nconstituents),pReal)
end subroutine homogenization_RGC_averageStressAndItsTangent
!--------------------------------------------------------------------------------------------------
!> @brief return array of homogenization results for post file inclusion
!--------------------------------------------------------------------------------------------------
pure function homogenization_RGC_postResults(ip,el) result(postResults)
use mesh, only: &
mesh_homogenizationAt
use material, only: &
homogenization_typeInstance,&
homogState, &
mappingHomogenization, &
homogenization_Noutput
implicit none
integer(pInt), intent(in) :: &
ip, & !< integration point number
el !< element number
integer(pInt) instance,o,c,nIntFaceTot
type(tParameters) :: prm
real(pReal), dimension(homogenization_RGC_sizePostResults(homogenization_typeInstance(mesh_homogenizationAt(el)))) :: &
postResults
instance = homogenization_typeInstance(mesh_homogenizationAt(el))
associate(prm => param(instance))
nIntFaceTot=(prm%Nconstituents(1)-1_pInt)*prm%Nconstituents(2)* prm%Nconstituents(3)&
+ prm%Nconstituents(1)* (prm%Nconstituents(2)-1_pInt)*prm%Nconstituents(3)&
+ prm%Nconstituents(1)* prm%Nconstituents(2)* (prm%Nconstituents(3)-1_pInt)
c = 0_pInt
postResults = 0.0_pReal
outputsLoop: do o = 1_pInt,size(prm%outputID)
select case(prm%outputID(o))
case (constitutivework_ID)
postResults(c+1) = homogState(mappingHomogenization(2,ip,el))% &
state(3*nIntFaceTot+1,mappingHomogenization(1,ip,el))
c = c + 1_pInt
case (magnitudemismatch_ID)
postResults(c+1) = homogState(mappingHomogenization(2,ip,el))% &
state(3*nIntFaceTot+2,mappingHomogenization(1,ip,el))
postResults(c+2) = homogState(mappingHomogenization(2,ip,el))% &
state(3*nIntFaceTot+3,mappingHomogenization(1,ip,el))
postResults(c+3) = homogState(mappingHomogenization(2,ip,el))% &
state(3*nIntFaceTot+4,mappingHomogenization(1,ip,el))
c = c + 3_pInt
case (penaltyenergy_ID)
postResults(c+1) = homogState(mappingHomogenization(2,ip,el))% &
state(3*nIntFaceTot+5,mappingHomogenization(1,ip,el))
c = c + 1_pInt
case (volumediscrepancy_ID)
postResults(c+1) = homogState(mappingHomogenization(2,ip,el))% &
state(3*nIntFaceTot+6,mappingHomogenization(1,ip,el))
c = c + 1_pInt
case (averagerelaxrate_ID)
postResults(c+1) = homogState(mappingHomogenization(2,ip,el))% &
state(3*nIntFaceTot+7,mappingHomogenization(1,ip,el))
c = c + 1_pInt
case (maximumrelaxrate_ID)
postResults(c+1) = homogState(mappingHomogenization(2,ip,el))% &
state(3*nIntFaceTot+8,mappingHomogenization(1,ip,el))
c = c + 1_pInt
end select
enddo outputsLoop
end associate
end function homogenization_RGC_postResults
!--------------------------------------------------------------------------------------------------
!> @brief collect relaxation vectors of an interface
!--------------------------------------------------------------------------------------------------
@ -1394,44 +1433,4 @@ pure function interface1to4(iFace1D, instance)
end function interface1to4
!--------------------------------------------------------------------------------------------------
!> @brief calculating the grain deformation gradient (the same with
! homogenization_RGC_partionDeformation, but used only for perturbation scheme)
!--------------------------------------------------------------------------------------------------
subroutine grainDeformation(F, avgF, ip, el)
use mesh, only: &
mesh_homogenizationAt
use material, only: &
homogenization_maxNgrains,&
homogenization_typeInstance
implicit none
real(pReal), dimension (3,3,homogenization_maxNgrains), intent(out) :: F !< partioned F per grain
real(pReal), dimension (3,3), intent(in) :: avgF !<
integer(pInt), intent(in) :: &
el, & !< element number
ip !< integration point number
real(pReal), dimension (3) :: aVect,nVect
integer(pInt), dimension (4) :: intFace
integer(pInt), dimension (3) :: iGrain3
integer(pInt) :: instance, iGrain,iFace,i,j
!--------------------------------------------------------------------------------------------------
! compute the deformation gradient of individual grains due to relaxations
instance = homogenization_typeInstance(mesh_homogenizationAt(el))
F = 0.0_pReal
do iGrain = 1_pInt,product(param(instance)%Nconstituents)
iGrain3 = grain1to3(iGrain,instance)
do iFace = 1_pInt,6_pInt
intFace = getInterface(iFace,iGrain3)
aVect = relaxationVector(intFace,instance, ip, el)
nVect = interfaceNormal(intFace,ip,el)
forall (i=1_pInt:3_pInt,j=1_pInt:3_pInt) &
F(i,j,iGrain) = F(i,j,iGrain) + aVect(i)*nVect(j) ! effective relaxations
enddo
F(1:3,1:3,iGrain) = F(1:3,1:3,iGrain) + avgF ! relaxed deformation gradient
enddo
end subroutine grainDeformation
end module homogenization_RGC