|  | 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_isoBrittle.f90" 
    4 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/source_damage_isoBrittle.f90" 
   24   integer,                       
dimension(:),           
allocatable :: &
 
   32     character(len=pStringLen), 
allocatable, 
dimension(:) :: &
 
   54   integer :: ninstance,sourceoffset,nipcmyphase,p
 
   55   character(len=pStringLen) :: extmsg = 
'' 
   61     write(6,
'(a16,1x,i5,/)') 
'# instances:',ninstance
 
   65   allocate(
param(ninstance))
 
   82     prm%N                = 
config%getFloat(
'isobrittle_n')
 
   83     prm%critStrainEnergy = 
config%getFloat(
'isobrittle_criticalstrainenergy')
 
   86     if (prm%N                <= 0.0_preal) extmsg = trim(extmsg)//
' isobrittle_n' 
   87     if (prm%critStrainEnergy <= 0.0_preal) extmsg = trim(extmsg)//
' isobrittle_criticalstrainenergy' 
   91     sourcestate(p)%p(sourceoffset)%atol = 
config%getFloat(
'isobrittle_atol',defaultval=1.0e-3_preal)
 
   92     if(any(
sourcestate(p)%p(sourceoffset)%atol < 0.0_preal)) extmsg = trim(extmsg)//
' isobrittle_atol' 
  110   integer, 
intent(in) :: &
 
  111     ipc, &                                                                                          !< component-ID of integration point
 
  112     ip, &                                                                                           !< integration point
 
  114   real(preal),  
intent(in), 
dimension(3,3) :: &
 
  116   real(preal),  
intent(in), 
dimension(6,6) :: &
 
  123   real(preal), 
dimension(6) :: &
 
  128   phase = material_phaseat(ipc,el)                                                                  
 
  129   constituent = material_phasememberat(ipc,ip,el)                                                   
 
  132   strain = 0.5_preal*math_sym33to6(matmul(transpose(fe),fe)-math_i3)
 
  135   strainenergy = 2.0_preal*sum(strain*matmul(c,strain))/prm%critStrainEnergy
 
  138   if (strainenergy > sourcestate(phase)%p(sourceoffset)%subState0(1,constituent)) 
then 
  139     sourcestate(phase)%p(sourceoffset)%deltaState(1,constituent) = &
 
  140       strainenergy - sourcestate(phase)%p(sourceoffset)%state(1,constituent)
 
  142     sourcestate(phase)%p(sourceoffset)%deltaState(1,constituent) = &
 
  143       sourcestate(phase)%p(sourceoffset)%subState0(1,constituent) - &
 
  144       sourcestate(phase)%p(sourceoffset)%state(1,constituent)
 
  156   integer, 
intent(in) :: &
 
  159   real(preal),  
intent(in) :: &
 
  161   real(preal),  
intent(out) :: &
 
  171   localphidot = (1.0_preal - phi)**(prm%n - 1.0_preal) &
 
  172               - phi*sourcestate(phase)%p(sourceoffset)%state(1,constituent)
 
  173   dlocalphidot_dphi = - (prm%n - 1.0_preal)* (1.0_preal - phi)**max(0.0_preal,prm%n - 2.0_preal) &
 
  174                       - sourcestate(phase)%p(sourceoffset)%state(1,constituent)
 
  185   integer,          
intent(in) :: phase
 
  186   character(len=*), 
intent(in) :: group
 
  192   outputsloop: 
do o = 1,
size(prm%output)
 
  193     select case(trim(prm%output(o)))
 
  194       case (
'isobrittle_drivingforce')
 
  195         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
material subroutine incoprorating isotropic brittle damage source mechanism
integer, dimension(:), allocatable source_damage_isobrittle_offset
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
@, public source_damage_isobrittle_id
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.
Reads in the material configuration from file.
integer, dimension(:), allocatable source_damage_isobrittle_instance
setting precision for real and int type
character(len= *), parameter, public source_damage_isobrittle_label
integer, public, protected discretization_nip
subroutine, public source_damage_isobrittle_deltastate(C, Fe, ipc, ip, el)
calculates derived quantities from state
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
subroutine, public source_damage_isobrittle_results(phase, group)
writes results to HDF5 output file
subroutine, public source_damage_isobrittle_getrateanditstangent(localphiDot, dLocalphiDot_dPhi, phi, phase, constituent)
returns local part of nonlocal damage driving force
subroutine, public source_damage_isobrittle_init
module initialization
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
integer, dimension(:), allocatable, public, protected phase_nsources
number of source mechanisms active in each phase
container type for internal constitutive parameters
type(tparameters), dimension(:), allocatable param
containers of constitutive parameters (len Ninstance)