|  | 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_anisoBrittle.f90" 
    4 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/source_damage_anisoBrittle.f90" 
   25   integer,                       
dimension(:),           
allocatable :: &
 
   33     real(
preal), 
dimension(:), 
allocatable :: &
 
   36     real(
preal), 
dimension(:,:,:,:), 
allocatable :: &
 
   40     character(len=pStringLen), 
allocatable, 
dimension(:) :: &
 
   62   integer :: ninstance,sourceoffset,nipcmyphase,p
 
   63   integer, 
dimension(:), 
allocatable :: n_cl
 
   64   character(len=pStringLen) :: extmsg = 
'' 
   70     write(6,
'(a16,1x,i5,/)') 
'# instances:',ninstance
 
   74   allocate(
param(ninstance))
 
   92     prm%sum_N_cl = sum(abs(n_cl))
 
   94     prm%n         = 
config%getFloat(
'anisobrittle_ratesensitivity')
 
   95     prm%sdot_0    = 
config%getFloat(
'anisobrittle_sdot0')
 
   97     prm%critDisp  = 
config%getFloats(
'anisobrittle_criticaldisplacement',requiredsize=
size(n_cl))
 
   98     prm%critLoad  = 
config%getFloats(
'anisobrittle_criticalload',        requiredsize=
size(n_cl))
 
  101                                                          config%getFloat(
'c/a',defaultval=0.0_preal))
 
  108     if (prm%n            <= 0.0_preal)  extmsg = trim(extmsg)//
' anisobrittle_n' 
  109     if (prm%sdot_0       <= 0.0_preal)  extmsg = trim(extmsg)//
' anisobrittle_sdot0' 
  110     if (any(prm%critLoad <  0.0_preal)) extmsg = trim(extmsg)//
' anisobrittle_critLoad' 
  111     if (any(prm%critDisp <  0.0_preal)) extmsg = trim(extmsg)//
' anisobrittle_critDisp' 
  115     sourcestate(p)%p(sourceoffset)%atol = 
config%getFloat(
'anisobrittle_atol',defaultval=1.0e-3_preal)
 
  116     if(any(
sourcestate(p)%p(sourceoffset)%atol < 0.0_preal)) extmsg = trim(extmsg)//
' anisobrittle_atol' 
  134   integer, 
intent(in) :: &
 
  135     ipc, &                                                                                          !< component-ID of integration point
 
  136     ip, &                                                                                           !< integration point
 
  138   real(preal),  
intent(in), 
dimension(3,3) :: &
 
  149     traction_d, traction_t, traction_n, traction_crit
 
  151   phase = material_phaseat(ipc,el)
 
  152   constituent = material_phasememberat(ipc,ip,el)
 
  154   homog = material_homogenizationat(el)
 
  155   damageoffset = damagemapping(homog)%p(ip,el)
 
  158   sourcestate(phase)%p(sourceoffset)%dotState(1,constituent) = 0.0_preal
 
  159   do i = 1, prm%sum_N_cl
 
  160     traction_d    = math_tensordot(s,prm%cleavage_systems(1:3,1:3,1,i))
 
  161     traction_t    = math_tensordot(s,prm%cleavage_systems(1:3,1:3,2,i))
 
  162     traction_n    = math_tensordot(s,prm%cleavage_systems(1:3,1:3,3,i))
 
  164     traction_crit = prm%critLoad(i)*damage(homog)%p(damageoffset)**2.0_preal
 
  166     sourcestate(phase)%p(sourceoffset)%dotState(1,constituent) &
 
  167     = sourcestate(phase)%p(sourceoffset)%dotState(1,constituent) &
 
  168     + prm%sdot_0 / prm%critDisp(i) &
 
  169       * ((max(0.0_preal, abs(traction_d) - traction_crit)/traction_crit)**prm%n + &
 
  170          (max(0.0_preal, abs(traction_t) - traction_crit)/traction_crit)**prm%n + &
 
  171          (max(0.0_preal, abs(traction_n) - traction_crit)/traction_crit)**prm%n)
 
  183   integer, 
intent(in) :: &
 
  186   real(preal),  
intent(in) :: &
 
  188   real(preal),  
intent(out) :: &
 
  197   dlocalphidot_dphi = -sourcestate(phase)%p(sourceoffset)%state(1,constituent)
 
  199   localphidot = 1.0_preal &
 
  200               + dlocalphidot_dphi*phi
 
  210   integer,          
intent(in) :: phase
 
  211   character(len=*), 
intent(in) :: group
 
  217   outputsloop: 
do o = 1,
size(prm%output)
 
  218     select case(trim(prm%output(o)))
 
  219       case (
'anisobrittle_drivingforce')
 
  220         call results_writedataset(group,stt,
'tbd',
'driving force',
'tbd')
 
 
 
type(tsourcestate), dimension(:), allocatable, public sourcestate
subroutine, public source_damage_anisobrittle_getrateanditstangent(localphiDot, dLocalphiDot_dPhi, phi, phase, constituent)
returns local part of nonlocal damage driving force
integer, dimension(:,:), allocatable, public, protected material_phaseat
phase ID of each element
integer, dimension(:), allocatable source_damage_anisobrittle_offset
which source is my current source mechanism?
pure real(preal) function, dimension(sum(how)) math_expand(what, how)
vector expansion
integer, dimension(:), allocatable source_damage_anisobrittle_instance
instance of 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
subroutine, public source_damage_anisobrittle_init
module initialization
Parses material config file, either solverJobName.materialConfig or material.config.
Reads in the material configuration from file.
setting precision for real and int type
integer, public, protected discretization_nip
input/output functions, partly depending on chosen solver
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.
subroutine, public material_allocatesourcestate(phase, of, NipcMyPhase, sizeState, sizeDotState, sizeDeltaState)
allocates the source state of a phase
@, public source_damage_anisobrittle_id
material subroutine incorporating anisotropic brittle damage source mechanism
character(len= *), parameter, public source_damage_anisobrittle_label
type(tparameters), dimension(:), allocatable param
containers of constitutive parameters (len Ninstance)
contains lattice structure definitions including Schmid matrices for slip, twin, trans,
integer, dimension(0), parameter emptyintarray
real(preal) function, dimension(3, 3, 3, sum(ncleavage)), public lattice_schmidmatrix_cleavage(Ncleavage, structure, cOverA)
Schmid matrix for cleavage details only active cleavage systems are considered.
Mathematical library, including random number generation and tensor representations.
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
subroutine, public source_damage_anisobrittle_dotstate(S, ipc, ip, el)
calculates derived quantities from state
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_anisobrittle_results(phase, group)
writes results to HDF5 output file