more control over initialisation of field values. specify initial field value in the homogenisation part of the material config file using the appropriate tags

This commit is contained in:
Pratheek Shanthraj 2015-07-24 14:53:50 +00:00
parent 87d42bf447
commit 7554647c8e
15 changed files with 87 additions and 47 deletions

View File

@ -68,6 +68,7 @@ subroutine damage_local_init(fileUnit)
damageState, &
damageMapping, &
damage, &
damage_initialPhi, &
material_partHomogenization
use numerics,only: &
worldrank
@ -162,9 +163,9 @@ subroutine damage_local_init(fileUnit)
sizeState = 1_pInt
damageState(homog)%sizeState = sizeState
damageState(homog)%sizePostResults = damage_local_sizePostResults(instance)
allocate(damageState(homog)%state0 (sizeState,NofMyHomog))
allocate(damageState(homog)%subState0(sizeState,NofMyHomog))
allocate(damageState(homog)%state (sizeState,NofMyHomog))
allocate(damageState(homog)%state0 (sizeState,NofMyHomog), source=damage_initialPhi(homog))
allocate(damageState(homog)%subState0(sizeState,NofMyHomog), source=damage_initialPhi(homog))
allocate(damageState(homog)%state (sizeState,NofMyHomog), source=damage_initialPhi(homog))
nullify(damageMapping(homog)%p)
damageMapping(homog)%p => mappingHomogenization(1,:,:)

View File

@ -46,12 +46,12 @@ subroutine damage_none_init()
NofMyHomog = count(material_homog == homog)
damageState(homog)%sizeState = 0_pInt
damageState(homog)%sizePostResults = 0_pInt
allocate(damageState(homog)%state0 (0_pInt,NofMyHomog), source=0.0_pReal)
allocate(damageState(homog)%subState0(0_pInt,NofMyHomog), source=0.0_pReal)
allocate(damageState(homog)%state (0_pInt,NofMyHomog), source=0.0_pReal)
allocate(damageState(homog)%state0 (0_pInt,NofMyHomog))
allocate(damageState(homog)%subState0(0_pInt,NofMyHomog))
allocate(damageState(homog)%state (0_pInt,NofMyHomog))
deallocate(damage(homog)%p)
allocate (damage(homog)%p(1), source=1.0_pReal)
allocate (damage(homog)%p(1), source=damage_initialPhi(homog))
endif myhomog
enddo initializeInstances

View File

@ -73,6 +73,7 @@ subroutine damage_nonlocal_init(fileUnit)
damageState, &
damageMapping, &
damage, &
damage_initialPhi, &
material_partHomogenization
use numerics,only: &
worldrank
@ -173,7 +174,7 @@ subroutine damage_nonlocal_init(fileUnit)
nullify(damageMapping(section)%p)
damageMapping(section)%p => mappingHomogenization(1,:,:)
deallocate(damage(section)%p)
allocate(damage(section)%p(NofMyHomog), source=1.0_pReal)
allocate(damage(section)%p(NofMyHomog), source=damage_initialPhi(section))
endif

View File

@ -158,11 +158,11 @@ subroutine homogenization_init(temperature_init)
if (.not. IO_open_jobFile_stat(FILEUNIT,material_localFileExt)) & ! no local material configuration present...
call IO_open_file(FILEUNIT,material_configFile) ! ... open material.config file
if (any(thermal_type == THERMAL_isothermal_ID)) &
call thermal_isothermal_init(temperature_init)
call thermal_isothermal_init()
if (any(thermal_type == THERMAL_adiabatic_ID)) &
call thermal_adiabatic_init(FILEUNIT,temperature_init)
call thermal_adiabatic_init(FILEUNIT)
if (any(thermal_type == THERMAL_conduction_ID)) &
call thermal_conduction_init(FILEUNIT,temperature_init)
call thermal_conduction_init(FILEUNIT)
close(FILEUNIT)
!--------------------------------------------------------------------------------------------------

View File

@ -79,6 +79,7 @@ subroutine hydrogenflux_cahnhilliard_init(fileUnit)
hydrogenfluxMapping, &
hydrogenConc, &
hydrogenConcRate, &
hydrogenflux_initialCh, &
material_partHomogenization, &
material_partPhase
use numerics,only: &
@ -187,7 +188,7 @@ subroutine hydrogenflux_cahnhilliard_init(fileUnit)
hydrogenfluxMapping(section)%p => mappingHomogenization(1,:,:)
deallocate(hydrogenConc (section)%p)
deallocate(hydrogenConcRate(section)%p)
allocate (hydrogenConc (section)%p(NofMyHomog), source=0.0_pReal)
allocate (hydrogenConc (section)%p(NofMyHomog), source=hydrogenflux_initialCh(section))
allocate (hydrogenConcRate(section)%p(NofMyHomog), source=0.0_pReal)
endif

View File

@ -52,7 +52,7 @@ subroutine hydrogenflux_isoconc_init()
deallocate(hydrogenConc (homog)%p)
deallocate(hydrogenConcRate(homog)%p)
allocate (hydrogenConc (homog)%p(1), source=0.0_pReal)
allocate (hydrogenConc (homog)%p(1), source=hydrogenflux_initialCh(homog))
allocate (hydrogenConcRate(homog)%p(1), source=0.0_pReal)
endif myhomog

View File

@ -190,15 +190,17 @@ module material
integer(pInt), dimension(:), allocatable, public, protected :: &
phase_Nsources, & !< number of source mechanisms active in each phase
phase_Nkinematics, & !< number of kinematic mechanisms active in each phase
phase_NstiffnessDegradations !< number of stiffness degradation mechanisms active in each phase
phase_NstiffnessDegradations, & !< number of stiffness degradation mechanisms active in each phase
phase_Noutput, & !< number of '(output)' items per phase
phase_elasticityInstance, & !< instance of particular elasticity of each phase
phase_plasticityInstance !< instance of particular plasticity of each phase
integer(pInt), dimension(:), allocatable, public, protected :: &
crystallite_Noutput !< number of '(output)' items per crystallite setting
integer(pInt), dimension(:), allocatable, public, protected :: &
homogenization_Ngrains, & !< number of grains in each homogenization
homogenization_Noutput, & !< number of '(output)' items per homogenization
phase_Noutput, & !< number of '(output)' items per phase
phase_elasticityInstance, & !< instance of particular elasticity of each phase
phase_plasticityInstance, & !< instance of particular plasticity of each phase
crystallite_Noutput, & !< number of '(output)' items per crystallite setting
homogenization_typeInstance, & !< instance of particular type of each homogenization
thermal_typeInstance, & !< instance of particular type of each thermal transport
damage_typeInstance, & !< instance of particular type of each nonlocal damage
@ -207,6 +209,13 @@ module material
hydrogenflux_typeInstance, & !< instance of particular type of each hydrogen flux
microstructure_crystallite !< crystallite setting ID of each microstructure
real(pReal), dimension(:), allocatable, public, protected :: &
thermal_initialT, & !< initial temperature per each homogenization
damage_initialPhi, & !< initial damage per each homogenization
vacancyflux_initialCv, & !< initial vacancy concentration per each homogenization
porosity_initialPhi, & !< initial posority per each homogenization
hydrogenflux_initialCh !< initial hydrogen concentration per each homogenization
integer(pInt), dimension(:,:,:), allocatable, public :: &
material_phase !< phase (index) of each grain,IP,element
integer(pInt), dimension(:,:), allocatable, public :: &
@ -511,11 +520,11 @@ subroutine material_init()
vacancyfluxMapping (myHomog)%p => mappingHomogenizationConst
porosityMapping (myHomog)%p => mappingHomogenizationConst
hydrogenfluxMapping(myHomog)%p => mappingHomogenizationConst
allocate(temperature (myHomog)%p(1), source=300.0_pReal)
allocate(damage (myHomog)%p(1), source=1.0_pReal)
allocate(vacancyConc (myHomog)%p(1), source=0.0_pReal)
allocate(porosity (myHomog)%p(1), source=1.0_pReal)
allocate(hydrogenConc (myHomog)%p(1), source=0.0_pReal)
allocate(temperature (myHomog)%p(1), source=thermal_initialT(myHomog))
allocate(damage (myHomog)%p(1), source=damage_initialPhi(myHomog))
allocate(vacancyConc (myHomog)%p(1), source=vacancyflux_initialCv(myHomog))
allocate(porosity (myHomog)%p(1), source=porosity_initialPhi(myHomog))
allocate(hydrogenConc (myHomog)%p(1), source=hydrogenflux_initialCh(myHomog))
allocate(temperatureRate (myHomog)%p(1), source=0.0_pReal)
allocate(vacancyConcRate (myHomog)%p(1), source=0.0_pReal)
allocate(hydrogenConcRate(myHomog)%p(1), source=0.0_pReal)
@ -539,6 +548,7 @@ subroutine material_parseHomogenization(fileUnit,myPart)
IO_isBlank, &
IO_stringValue, &
IO_intValue, &
IO_floatValue, &
IO_stringPos, &
IO_EOF
use mesh, only: &
@ -577,6 +587,11 @@ subroutine material_parseHomogenization(fileUnit,myPart)
allocate(homogenization_Ngrains(Nsections), source=0_pInt)
allocate(homogenization_Noutput(Nsections), source=0_pInt)
allocate(homogenization_active(Nsections), source=.false.) !!!!!!!!!!!!!!!
allocate(thermal_initialT(Nsections), source=300.0_pReal)
allocate(damage_initialPhi(Nsections), source=1.0_pReal)
allocate(vacancyflux_initialCv(Nsections), source=0.0_pReal)
allocate(porosity_initialPhi(Nsections), source=1.0_pReal)
allocate(hydrogenflux_initialCh(Nsections), source=0.0_pReal)
forall (s = 1_pInt:Nsections) homogenization_active(s) = any(mesh_element(3,:) == s) ! current homogenization used in model? Homogenization view, maximum operations depend on maximum number of homog schemes
homogenization_Noutput = IO_countTagInPart(fileUnit,myPart,'(output)',Nsections)
@ -677,6 +692,22 @@ subroutine material_parseHomogenization(fileUnit,myPart)
case ('nconstituents','ngrains')
homogenization_Ngrains(section) = IO_intValue(line,positions,2_pInt)
case ('initialtemperature','initialt')
thermal_initialT(section) = IO_floatValue(line,positions,2_pInt)
case ('initialdamage')
damage_initialPhi(section) = IO_floatValue(line,positions,2_pInt)
case ('initialvacancyconc','initialcv')
vacancyflux_initialCv(section) = IO_floatValue(line,positions,2_pInt)
case ('initialporosity')
porosity_initialPhi(section) = IO_floatValue(line,positions,2_pInt)
case ('initialhydrogenconc','initialch')
hydrogenflux_initialCh(section) = IO_floatValue(line,positions,2_pInt)
end select
endif
enddo

View File

@ -51,7 +51,7 @@ subroutine porosity_none_init()
allocate(porosityState(homog)%state (0_pInt,NofMyHomog), source=0.0_pReal)
deallocate(porosity(homog)%p)
allocate (porosity(homog)%p(1), source=1.0_pReal)
allocate (porosity(homog)%p(1), source=porosity_initialPhi(homog))
endif myhomog
enddo initializeInstances

View File

@ -75,6 +75,7 @@ subroutine porosity_phasefield_init(fileUnit)
porosityState, &
porosityMapping, &
porosity, &
porosity_initialPhi, &
material_partHomogenization, &
material_partPhase
use numerics,only: &
@ -176,7 +177,7 @@ subroutine porosity_phasefield_init(fileUnit)
nullify(porosityMapping(section)%p)
porosityMapping(section)%p => mappingHomogenization(1,:,:)
deallocate(porosity(section)%p)
allocate(porosity(section)%p(NofMyHomog), source=1.0_pReal)
allocate(porosity(section)%p(NofMyHomog), source=porosity_initialPhi(section))
endif

View File

@ -47,7 +47,7 @@ contains
!> @brief module initialization
!> @details reads in material parameters, allocates arrays, and does sanity checks
!--------------------------------------------------------------------------------------------------
subroutine thermal_adiabatic_init(fileUnit,temperature_init)
subroutine thermal_adiabatic_init(fileUnit)
use, intrinsic :: iso_fortran_env ! to get compiler_version and compiler_options (at least for gfortran 4.6 at the moment)
use IO, only: &
IO_read, &
@ -72,6 +72,7 @@ subroutine thermal_adiabatic_init(fileUnit,temperature_init)
mappingHomogenization, &
thermalState, &
thermalMapping, &
thermal_initialT, &
temperature, &
temperatureRate, &
material_partHomogenization
@ -79,7 +80,6 @@ subroutine thermal_adiabatic_init(fileUnit,temperature_init)
worldrank
implicit none
real(pReal), intent(in) :: temperature_init !< initial temperature
integer(pInt), intent(in) :: fileUnit
integer(pInt), parameter :: MAXNCHUNKS = 7_pInt
@ -168,9 +168,9 @@ subroutine thermal_adiabatic_init(fileUnit,temperature_init)
sizeState = 1_pInt
thermalState(section)%sizeState = sizeState
thermalState(section)%sizePostResults = thermal_adiabatic_sizePostResults(instance)
allocate(thermalState(section)%state0 (sizeState,NofMyHomog), source=temperature_init)
allocate(thermalState(section)%subState0(sizeState,NofMyHomog), source=temperature_init)
allocate(thermalState(section)%state (sizeState,NofMyHomog), source=temperature_init)
allocate(thermalState(section)%state0 (sizeState,NofMyHomog), source=thermal_initialT(section))
allocate(thermalState(section)%subState0(sizeState,NofMyHomog), source=thermal_initialT(section))
allocate(thermalState(section)%state (sizeState,NofMyHomog), source=thermal_initialT(section))
nullify(thermalMapping(section)%p)
thermalMapping(section)%p => mappingHomogenization(1,:,:)

View File

@ -48,7 +48,7 @@ contains
!> @brief module initialization
!> @details reads in material parameters, allocates arrays, and does sanity checks
!--------------------------------------------------------------------------------------------------
subroutine thermal_conduction_init(fileUnit,temperature_init)
subroutine thermal_conduction_init(fileUnit)
use, intrinsic :: iso_fortran_env ! to get compiler_version and compiler_options (at least for gfortran 4.6 at the moment)
use IO, only: &
IO_read, &
@ -73,6 +73,7 @@ subroutine thermal_conduction_init(fileUnit,temperature_init)
mappingHomogenization, &
thermalState, &
thermalMapping, &
thermal_initialT, &
temperature, &
temperatureRate, &
material_partHomogenization
@ -80,7 +81,6 @@ subroutine thermal_conduction_init(fileUnit,temperature_init)
worldrank
implicit none
real(pReal), intent(in) :: temperature_init !< initial temperature
integer(pInt), intent(in) :: fileUnit
integer(pInt), parameter :: MAXNCHUNKS = 7_pInt
@ -176,7 +176,7 @@ subroutine thermal_conduction_init(fileUnit,temperature_init)
nullify(thermalMapping(section)%p)
thermalMapping(section)%p => mappingHomogenization(1,:,:)
deallocate(temperature (section)%p)
allocate (temperature (section)%p(NofMyHomog), source=temperature_init)
allocate (temperature (section)%p(NofMyHomog), source=thermal_initialT(section))
deallocate(temperatureRate(section)%p)
allocate (temperatureRate(section)%p(NofMyHomog), source=0.0_pReal)

View File

@ -17,7 +17,7 @@ contains
!--------------------------------------------------------------------------------------------------
!> @brief allocates all neccessary fields, reads information from material configuration file
!--------------------------------------------------------------------------------------------------
subroutine thermal_isothermal_init(temperature_init)
subroutine thermal_isothermal_init()
use, intrinsic :: iso_fortran_env ! to get compiler_version and compiler_options (at least for gfortran 4.6 at the moment)
use prec, only: &
pReal, &
@ -29,7 +29,6 @@ subroutine thermal_isothermal_init(temperature_init)
worldrank
implicit none
real(pReal), intent(in) :: temperature_init !< initial temperature
integer(pInt) :: &
homog, &
NofMyHomog, &
@ -54,7 +53,7 @@ subroutine thermal_isothermal_init(temperature_init)
allocate(thermalState(homog)%state (sizeState,NofMyHomog), source=0.0_pReal)
deallocate(temperature (homog)%p)
allocate (temperature (homog)%p(1), source=temperature_init)
allocate (temperature (homog)%p(1), source=thermal_initialT(homog))
deallocate(temperatureRate(homog)%p)
allocate (temperatureRate(homog)%p(1), source=0.0_pReal)

View File

@ -89,6 +89,7 @@ subroutine vacancyflux_cahnhilliard_init(fileUnit)
vacancyfluxMapping, &
vacancyConc, &
vacancyConcRate, &
vacancyflux_initialCv, &
material_partHomogenization, &
material_partPhase
use numerics,only: &
@ -205,7 +206,7 @@ subroutine vacancyflux_cahnhilliard_init(fileUnit)
nullify(vacancyfluxMapping(section)%p)
vacancyfluxMapping(section)%p => mappingHomogenization(1,:,:)
deallocate(vacancyConc (section)%p)
allocate (vacancyConc (section)%p(NofMyHomog), source=0.0_pReal)
allocate (vacancyConc (section)%p(NofMyHomog), source=vacancyflux_initialCv(section))
deallocate(vacancyConcRate(section)%p)
allocate (vacancyConcRate(section)%p(NofMyHomog), source=0.0_pReal)
@ -510,8 +511,7 @@ subroutine vacancyflux_cahnhilliard_getChemPotAndItsTangent(ChemPot,dChemPot_dCv
VoidPhaseFrac = porosity(homog)%p(porosityMapping(homog)%p(ip,el))
kBT = vacancyflux_cahnhilliard_getEntropicCoeff(ip,el)
ChemPot = vacancyflux_cahnhilliard_getFormationEnergy(ip,el)* &
vacancyflux_cahnhilliard_thermalFluc(vacancyflux_typeInstance(homog))%p(mappingHomogenization(1,ip,el))
ChemPot = vacancyflux_cahnhilliard_getFormationEnergy(ip,el)
dChemPot_dCv = 0.0_pReal
do o = 1_pInt, vacancyPolyOrder
ChemPot = ChemPot + kBT*((2.0_pReal*Cv - 1.0_pReal)**real(2_pInt*o-1_pInt,pReal))/ &
@ -535,7 +535,12 @@ subroutine vacancyflux_cahnhilliard_getChemPotAndItsTangent(ChemPot,dChemPot_dCv
elseif (Cv > 1.0_pReal) then
ChemPot = ChemPot + 3.0_pReal*vacancyBoundPenalty*(1.0_pReal - Cv)*(1.0_pReal - Cv)
dChemPot_dCv = dChemPot_dCv - 6.0_pReal*vacancyBoundPenalty*(1.0_pReal - Cv)
endif
endif
ChemPot = ChemPot* &
vacancyflux_cahnhilliard_thermalFluc(vacancyflux_typeInstance(homog))%p(mappingHomogenization(1,ip,el))
dChemPot_dCv = dChemPot_dCv* &
vacancyflux_cahnhilliard_thermalFluc(vacancyflux_typeInstance(homog))%p(mappingHomogenization(1,ip,el))
end subroutine vacancyflux_cahnhilliard_getChemPotAndItsTangent

View File

@ -72,6 +72,7 @@ subroutine vacancyflux_isochempot_init(fileUnit)
vacancyfluxMapping, &
vacancyConc, &
vacancyConcRate, &
vacancyflux_initialCv, &
material_partHomogenization
use numerics,only: &
worldrank
@ -165,9 +166,9 @@ subroutine vacancyflux_isochempot_init(fileUnit)
sizeState = 1_pInt
vacancyfluxState(section)%sizeState = sizeState
vacancyfluxState(section)%sizePostResults = vacancyflux_isochempot_sizePostResults(instance)
allocate(vacancyfluxState(section)%state0 (sizeState,NofMyHomog), source=0.0_pReal)
allocate(vacancyfluxState(section)%subState0(sizeState,NofMyHomog), source=0.0_pReal)
allocate(vacancyfluxState(section)%state (sizeState,NofMyHomog), source=0.0_pReal)
allocate(vacancyfluxState(section)%state0 (sizeState,NofMyHomog), source=vacancyflux_initialCv(section))
allocate(vacancyfluxState(section)%subState0(sizeState,NofMyHomog), source=vacancyflux_initialCv(section))
allocate(vacancyfluxState(section)%state (sizeState,NofMyHomog), source=vacancyflux_initialCv(section))
nullify(vacancyfluxMapping(section)%p)
vacancyfluxMapping(section)%p => mappingHomogenization(1,:,:)

View File

@ -46,12 +46,12 @@ subroutine vacancyflux_isoconc_init()
NofMyHomog = count(material_homog == homog)
vacancyfluxState(homog)%sizeState = 0_pInt
vacancyfluxState(homog)%sizePostResults = 0_pInt
allocate(vacancyfluxState(homog)%state0 (0_pInt,NofMyHomog), source=0.0_pReal)
allocate(vacancyfluxState(homog)%subState0(0_pInt,NofMyHomog), source=0.0_pReal)
allocate(vacancyfluxState(homog)%state (0_pInt,NofMyHomog), source=0.0_pReal)
allocate(vacancyfluxState(homog)%state0 (0_pInt,NofMyHomog))
allocate(vacancyfluxState(homog)%subState0(0_pInt,NofMyHomog))
allocate(vacancyfluxState(homog)%state (0_pInt,NofMyHomog))
deallocate(vacancyConc (homog)%p)
allocate (vacancyConc (homog)%p(1), source=0.0_pReal)
allocate (vacancyConc (homog)%p(1), source=vacancyflux_initialCv(homog))
deallocate(vacancyConcRate(homog)%p)
allocate (vacancyConcRate(homog)%p(1), source=0.0_pReal)