|
DAMASK with grid solvers
Revision: v2.0.3-2204-gdb1f2151
The Düsseldorf Advanced Material Simulation Kit with Grid Solvers
|
|
Go to the documentation of this file. 1 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/source_damage_anisoDuctile.f90"
4 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/source_damage_anisoDuctile.f90"
24 integer,
dimension(:),
allocatable :: &
31 real(
preal),
dimension(:),
allocatable :: &
33 character(len=pStringLen),
allocatable,
dimension(:) :: &
55 integer :: ninstance,sourceoffset,nipcmyphase,p
56 integer,
dimension(:),
allocatable :: n_sl
57 character(len=pStringLen) :: extmsg =
''
63 write(6,
'(a16,1x,i5,/)')
'# instances:',ninstance
67 allocate(
param(ninstance))
85 prm%n =
config%getFloat(
'anisoductile_ratesensitivity')
86 prm%critPlasticStrain =
config%getFloats(
'anisoductile_criticalplasticstrain',requiredsize=
size(n_sl))
89 prm%critPlasticStrain =
math_expand(prm%critPlasticStrain,n_sl)
92 if (prm%n <= 0.0_preal) extmsg = trim(extmsg)//
' anisoductile_ratesensitivity'
93 if (any(prm%critPlasticStrain < 0.0_preal)) extmsg = trim(extmsg)//
' anisoductile_criticalplasticstrain'
97 sourcestate(p)%p(sourceoffset)%atol =
config%getFloat(
'anisoductile_atol',defaultval=1.0e-3_preal)
98 if(any(
sourcestate(p)%p(sourceoffset)%atol < 0.0_preal)) extmsg = trim(extmsg)//
' anisoductile_atol'
116 integer,
intent(in) :: &
117 ipc, & !< component-ID of integration point
118 ip, & !< integration point
128 phase = material_phaseat(ipc,el)
129 constituent = material_phasememberat(ipc,ip,el)
131 homog = material_homogenizationat(el)
132 damageoffset = damagemapping(homog)%p(ip,el)
135 sourcestate(phase)%p(sourceoffset)%dotState(1,constituent) &
136 = sum(plasticstate(phase)%slipRate(:,constituent)/(damage(homog)%p(damageoffset)**prm%n)/prm%critPlasticStrain)
147 integer,
intent(in) :: &
150 real(preal),
intent(in) :: &
152 real(preal),
intent(out) :: &
161 dlocalphidot_dphi = -sourcestate(phase)%p(sourceoffset)%state(1,constituent)
163 localphidot = 1.0_preal &
164 + dlocalphidot_dphi*phi
174 integer,
intent(in) :: phase
175 character(len=*),
intent(in) :: group
181 outputsloop:
do o = 1,
size(prm%output)
182 select case(trim(prm%output(o)))
183 case (
'anisoductile_drivingforce')
184 call results_writedataset(group,stt,
'tbd',
'driving force',
'tbd')
type(tsourcestate), dimension(:), allocatable, public sourcestate
integer, dimension(:,:), allocatable, public, protected material_phaseat
phase ID of each element
pure real(preal) function, dimension(sum(how)) math_expand(what, how)
vector expansion
integer, dimension(:), allocatable source_damage_anisoductile_instance
instance of damage source mechanism
subroutine, public io_error(error_ID, el, ip, g, instance, ext_msg)
write error statements to standard out and terminate the Marc/spectral run with exit #9xxx
type(tpartitionedstringlist), dimension(:), allocatable, public, protected config_phase
integer, dimension(debug_maxntype+2), public, protected debug_level
character(len=pstringlen), dimension(0), parameter emptystringarray
Parses material config file, either solverJobName.materialConfig or material.config.
subroutine, public source_damage_anisoductile_results(phase, group)
writes results to HDF5 output file
Reads in the material configuration from file.
setting precision for real and int type
subroutine, public source_damage_anisoductile_init
module initialization
integer, public, protected discretization_nip
type(tparameters), dimension(:), allocatable, private param
containers of constitutive parameters (len Ninstance)
input/output functions, partly depending on chosen solver
material subroutine incorporating anisotropic ductile damage source mechanism
integer, parameter, public debug_levelbasic
integer, parameter preal
number with 15 significant digits, up to 1e+-307 (typically 64 bit)
Reading in and interpretating the debugging settings for the various modules.
integer, dimension(:), allocatable source_damage_anisoductile_offset
which source is my current damage mechanism?
subroutine, public material_allocatesourcestate(phase, of, NipcMyPhase, sizeState, sizeDotState, sizeDeltaState)
allocates the source state of a phase
integer, dimension(0), parameter emptyintarray
subroutine, public source_damage_anisoductile_getrateanditstangent(localphiDot, dLocalphiDot_dPhi, phi, phase, constituent)
returns local part of nonlocal damage driving force
Mathematical library, including random number generation and tensor representations.
character(len= *), parameter, public source_damage_anisoductile_label
integer, parameter, public debug_constitutive
stores debug level for constitutive part of DAMASK bitwise coded
integer(kind(source_undefined_id)), dimension(:,:), allocatable, public, protected phase_source
active sources mechanisms of each phase
integer, dimension(:), allocatable, public, protected phase_nsources
number of source mechanisms active in each phase
container type for internal constitutive parameters
subroutine, public source_damage_anisoductile_dotstate(ipc, ip, el)
calculates derived quantities from state
@, public source_damage_anisoductile_id