DAMASK with MSC.Marc FEM solver  Revision: v2.0.3-2204-gdb1f2151
The Düsseldorf Advanced Material Simulation Kit with MSC.Marc
lattice Module Reference

contains lattice structure definitions including Schmid matrices for slip, twin, trans, More...

Data Types

interface  lattice_forestprojection_edge
 
interface  lattice_forestprojection_screw
 

Enumerations

enum  { lattice_undefined_id }
 

Functions/Subroutines

subroutine, public lattice_init
 Module initialization. More...
 
real(preal) function, dimension(sum(ntwin)), public lattice_characteristicshear_twin (Ntwin, structure, CoverA)
 Characteristic shear for twinning. More...
 
real(preal) function, dimension(6, 6, sum(ntwin)), public lattice_c66_twin (Ntwin, C66, structure, CoverA)
 Rotated elasticity matrices for twinning in 66-vector notation. More...
 
real(preal) function, dimension(6, 6, sum(ntrans)), public lattice_c66_trans (Ntrans, C_parent66, structure_target, cOverA_trans, a_bcc, a_fcc)
 Rotated elasticity matrices for transformation in 66-vector notation. More...
 
real(preal) function, dimension(1:3, 1:3, sum(nslip)), public lattice_nonschmidmatrix (Nslip, nonSchmidCoefficients, sense)
 Non-schmid projections for bcc with up to 6 coefficients. More...
 
real(preal) function, dimension(sum(nslip), sum(nslip)), public lattice_interaction_slipbyslip (Nslip, interactionValues, structure)
 Slip-slip interaction matrix details only active slip systems are considered. More...
 
real(preal) function, dimension(sum(ntwin), sum(ntwin)), public lattice_interaction_twinbytwin (Ntwin, interactionValues, structure)
 Twin-twin interaction matrix details only active twin systems are considered. More...
 
real(preal) function, dimension(sum(ntrans), sum(ntrans)), public lattice_interaction_transbytrans (Ntrans, interactionValues, structure)
 Trans-trans interaction matrix details only active trans systems are considered. More...
 
real(preal) function, dimension(sum(nslip), sum(ntwin)), public lattice_interaction_slipbytwin (Nslip, Ntwin, interactionValues, structure)
 Slip-twin interaction matrix details only active slip and twin systems are considered. More...
 
real(preal) function, dimension(sum(nslip), sum(ntrans)), public lattice_interaction_slipbytrans (Nslip, Ntrans, interactionValues, structure)
 Slip-trans interaction matrix details only active slip and trans systems are considered. More...
 
real(preal) function, dimension(sum(ntwin), sum(nslip)), public lattice_interaction_twinbyslip (Ntwin, Nslip, interactionValues, structure)
 Twin-slip interaction matrix details only active twin and slip systems are considered. More...
 
real(preal) function, dimension(3, 3, sum(nslip)), public lattice_schmidmatrix_slip (Nslip, structure, cOverA)
 Schmid matrix for slip details only active slip systems are considered. More...
 
real(preal) function, dimension(3, 3, sum(ntwin)), public lattice_schmidmatrix_twin (Ntwin, structure, cOverA)
 Schmid matrix for twinning details only active twin systems are considered. More...
 
real(preal) function, dimension(3, 3, sum(ntrans)), public lattice_schmidmatrix_trans (Ntrans, structure_target, cOverA, a_bcc, a_fcc)
 Schmid matrix for twinning details only active twin systems are considered. More...
 
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. More...
 
real(preal) function, dimension(3, sum(nslip)), public lattice_slip_direction (Nslip, structure, cOverA)
 Slip direction of slip systems (|| b) More...
 
real(preal) function, dimension(3, sum(nslip)), public lattice_slip_normal (Nslip, structure, cOverA)
 Normal direction of slip systems (|| n) More...
 
real(preal) function, dimension(3, sum(nslip)), public lattice_slip_transverse (Nslip, structure, cOverA)
 Transverse direction of slip systems ( || t = b x n) More...
 
character(len=:) function, dimension(:), allocatable, public lattice_labels_slip (Nslip, structure)
 Labels for slip systems details only active slip systems are considered. More...
 
real(preal) function, dimension(3, 3), public lattice_applylatticesymmetry33 (T, structure)
 Return 3x3 tensor with symmetry according to given crystal structure. More...
 
real(preal) function, dimension(6, 6) applylatticesymmetryc66 (C66, structure)
 Return stiffness matrix in 6x6 notation with symmetry according to given crystal structure. More...
 
character(len=:) function, dimension(:), allocatable, public lattice_labels_twin (Ntwin, structure)
 Labels for twin systems details only active twin systems are considered. More...
 
real(preal) function, dimension(sum(nslip), sum(nslip)) slipprojection_transverse (Nslip, structure, cOverA)
 Projection of the transverse direction onto the slip plane. More...
 
real(preal) function, dimension(sum(nslip), sum(nslip)) slipprojection_direction (Nslip, structure, cOverA)
 Projection of the slip direction onto the slip plane. More...
 
real(preal) function, dimension(3, 3, sum(nslip)) coordinatesystem_slip (Nslip, structure, cOverA)
 build a local coordinate system on slip systems More...
 
real(preal) function, dimension(sum(reacting_used), sum(acting_used)) buildinteraction (reacting_used, acting_used, reacting_max, acting_max, values, matrix)
 Populate reduced interaction matrix. More...
 
real(preal) function, dimension(3, 3, sum(active)) buildcoordinatesystem (active, potential, system, structure, cOverA)
 Build a local coordinate system on slip, twin, trans, cleavage systems. More...
 
subroutine buildtransformationsystem (Q, S, Ntrans, cOverA, a_fcc, a_bcc)
 Helper function to define transformation systems. More...
 
character(len=:) function, dimension(:), allocatable getlabels (active, potential, system)
 select active systems as strings More...
 
real(preal) function equivalent_nu (C, assumption)
 Equivalent Poisson's ratio (ν) More...
 
real(preal) function equivalent_mu (C, assumption)
 Equivalent shear modulus (μ) More...
 
subroutine unittest
 check correctness of some lattice functions More...
 

Variables

integer, dimension(2), parameter fcc_nslipsystem = [12, 6]
 
integer, dimension(1), parameter fcc_ntwinsystem = [12]
 
integer, dimension(1), parameter fcc_ntranssystem = [12]
 
integer, dimension(1), parameter fcc_ncleavagesystem = [3]
 
integer, parameter fcc_nslip = sum(FCC_NSLIPSYSTEM)
 total # of slip systems for fcc More...
 
integer, parameter fcc_ntwin = sum(FCC_NTWINSYSTEM)
 total # of twin systems for fcc More...
 
integer, parameter fcc_ntrans = sum(FCC_NTRANSSYSTEM)
 total # of transformation systems for fcc More...
 
integer, parameter fcc_ncleavage = sum(FCC_NCLEAVAGESYSTEM)
 total # of cleavage systems for fcc More...
 
real(preal), dimension(3+3, fcc_nslip), parameter fcc_systemslip = reshape(real([ 0, 1,-1, 1, 1, 1, -1, 0, 1, 1, 1, 1, 1,-1, 0, 1, 1, 1, 0,-1,-1, -1,-1, 1, 1, 0, 1, -1,-1, 1, -1, 1, 0, -1,-1, 1, 0,-1, 1, 1,-1,-1, -1, 0,-1, 1,-1,-1, 1, 1, 0, 1,-1,-1, 0, 1, 1, -1, 1,-1, 1, 0,-1, -1, 1,-1, -1,-1, 0, -1, 1,-1, 1, 1, 0, 1,-1, 0, 1,-1, 0, 1, 1, 0, 1, 0, 1, 1, 0,-1, 1, 0,-1, 1, 0, 1, 0, 1, 1, 0, 1,-1, 0, 1,-1, 0, 1, 1 ], pReal), shape(FCC_SYSTEMSLIP))
 Slip system <110>{111} directions. Sorted according to Eisenlohr & Hantcherli. More...
 
real(preal), dimension(3+3, fcc_ntwin), parameter fcc_systemtwin = reshape(real( [ -2, 1, 1, 1, 1, 1, 1,-2, 1, 1, 1, 1, 1, 1,-2, 1, 1, 1, 2,-1, 1, -1,-1, 1, -1, 2, 1, -1,-1, 1, -1,-1,-2, -1,-1, 1, -2,-1,-1, 1,-1,-1, 1, 2,-1, 1,-1,-1, 1,-1, 2, 1,-1,-1, 2, 1,-1, -1, 1,-1, -1,-2,-1, -1, 1,-1, -1, 1, 2, -1, 1,-1 ], pReal), shape(FCC_SYSTEMTWIN))
 Twin system <112>{111} directions. Sorted according to Eisenlohr & Hantcherli. More...
 
integer, dimension(2, fcc_ntwin), parameter, public lattice_fcc_twinnucleationslippair = reshape( [ 2,3, 1,3, 1,2, 5,6, 4,6, 4,5, 8,9, 7,9, 7,8, 11,12, 10,12, 10,11 ], shape(lattice_FCC_TWINNUCLEATIONSLIPPAIR))
 
real(preal), dimension(3+3, fcc_ncleavage), parameter fcc_systemcleavage = reshape(real([ 0, 1, 0, 1, 0, 0, 0, 0, 1, 0, 1, 0, 1, 0, 0, 0, 0, 1 ], pReal), shape(FCC_SYSTEMCLEAVAGE))
 
integer, dimension(2), parameter bcc_nslipsystem = [12, 12]
 
integer, dimension(1), parameter bcc_ntwinsystem = [12]
 
integer, dimension(1), parameter bcc_ncleavagesystem = [3]
 
integer, parameter bcc_nslip = sum(BCC_NSLIPSYSTEM)
 total # of slip systems for bcc More...
 
integer, parameter bcc_ntwin = sum(BCC_NTWINSYSTEM)
 total # of twin systems for bcc More...
 
integer, parameter bcc_ncleavage = sum(BCC_NCLEAVAGESYSTEM)
 total # of cleavage systems for bcc More...
 
real(preal), dimension(3+3, bcc_nslip), parameter bcc_systemslip = reshape(real([ 1,-1, 1, 0, 1, 1, -1,-1, 1, 0, 1, 1, 1, 1, 1, 0,-1, 1, -1, 1, 1, 0,-1, 1, -1, 1, 1, 1, 0, 1, -1,-1, 1, 1, 0, 1, 1, 1, 1, -1, 0, 1, 1,-1, 1, -1, 0, 1, -1, 1, 1, 1, 1, 0, -1, 1,-1, 1, 1, 0, 1, 1, 1, -1, 1, 0, 1, 1,-1, -1, 1, 0, -1, 1, 1, 2, 1, 1, 1, 1, 1, -2, 1, 1, 1, 1,-1, 2,-1, 1, 1,-1, 1, 2, 1,-1, 1,-1, 1, 1, 2, 1, 1, 1,-1, -1, 2, 1, 1, 1, 1, 1,-2, 1, -1, 1, 1, 1, 2,-1, 1, 1,-1, 1, 1, 2, 1,-1, 1, -1, 1, 2, -1, 1, 1, 1,-1, 2, 1, 1, 1, 1, 1,-2 ], pReal), shape(BCC_SYSTEMSLIP))
 
real(preal), dimension(3+3, bcc_ntwin), parameter bcc_systemtwin = reshape(real([ -1, 1, 1, 2, 1, 1, 1, 1, 1, -2, 1, 1, 1, 1,-1, 2,-1, 1, 1,-1, 1, 2, 1,-1, 1,-1, 1, 1, 2, 1, 1, 1,-1, -1, 2, 1, 1, 1, 1, 1,-2, 1, -1, 1, 1, 1, 2,-1, 1, 1,-1, 1, 1, 2, 1,-1, 1, -1, 1, 2, -1, 1, 1, 1,-1, 2, 1, 1, 1, 1, 1,-2 ], pReal), shape(BCC_SYSTEMTWIN))
 
real(preal), dimension(3+3, bcc_ncleavage), parameter bcc_systemcleavage = reshape(real([ 0, 1, 0, 1, 0, 0, 0, 0, 1, 0, 1, 0, 1, 0, 0, 0, 0, 1 ], pReal), shape(BCC_SYSTEMCLEAVAGE))
 
integer, dimension(6), parameter hex_nslipsystem = [3, 3, 3, 6, 12, 6]
 
integer, dimension(4), parameter hex_ntwinsystem = [6, 6, 6, 6]
 
integer, parameter hex_nslip = sum(HEX_NSLIPSYSTEM)
 total # of slip systems for hex More...
 
integer, parameter hex_ntwin = sum(HEX_NTWINSYSTEM)
 total # of twin systems for hex More...
 
real(preal), dimension(4+4, hex_nslip), parameter hex_systemslip = reshape(real([ 2, -1, -1, 0, 0, 0, 0, 1, -1, 2, -1, 0, 0, 0, 0, 1, -1, -1, 2, 0, 0, 0, 0, 1, 2, -1, -1, 0, 0, 1, -1, 0, -1, 2, -1, 0, -1, 0, 1, 0, -1, -1, 2, 0, 1, -1, 0, 0, -1, 1, 0, 0, 1, 1, -2, 0, 0, -1, 1, 0, -2, 1, 1, 0, 1, 0, -1, 0, 1, -2, 1, 0, -1, 2, -1, 0, 1, 0, -1, 1, -2, 1, 1, 0, 0, 1, -1, 1, -1, -1, 2, 0, -1, 1, 0, 1, 1, -2, 1, 0, -1, 0, 1, 1, 2, -1, -1, 0, 0, -1, 1, 1, 1, 1, -2, 0, 1, -1, 0, 1, -2, 1, 1, 3, 1, 0, -1, 1, -1, -1, 2, 3, 1, 0, -1, 1, -1, -1, 2, 3, 0, 1, -1, 1, 1, -2, 1, 3, 0, 1, -1, 1, 1, -2, 1, 3, -1, 1, 0, 1, 2, -1, -1, 3, -1, 1, 0, 1, 2, -1, -1, 3, -1, 0, 1, 1, 1, 1, -2, 3, -1, 0, 1, 1, 1, 1, -2, 3, 0, -1, 1, 1, -1, 2, -1, 3, 0, -1, 1, 1, -1, 2, -1, 3, 1, -1, 0, 1, -2, 1, 1, 3, 1, -1, 0, 1, -1, -1, 2, 3, 1, 1, -2, 2, 1, -2, 1, 3, -1, 2, -1, 2, 2, -1, -1, 3, -2, 1, 1, 2, 1, 1, -2, 3, -1, -1, 2, 2, -1, 2, -1, 3, 1, -2, 1, 2, -2, 1, 1, 3, 2, -1, -1, 2 ], pReal), shape(HEX_SYSTEMSLIP))
 slip systems for hex, sorted by P. Eisenlohr CCW around starting next to a_1 axis More...
 
real(preal), dimension(4+4, hex_ntwin), parameter hex_systemtwin = reshape(real([ -1, 0, 1, 1, 1, 0, -1, 2, 0, -1, 1, 1, 0, 1, -1, 2, 1, -1, 0, 1, -1, 1, 0, 2, 1, 0, -1, 1, -1, 0, 1, 2, 0, 1, -1, 1, 0, -1, 1, 2, -1, 1, 0, 1, 1, -1, 0, 2, -1, -1, 2, 6, 1, 1, -2, 1, 1, -2, 1, 6, -1, 2, -1, 1, 2, -1, -1, 6, -2, 1, 1, 1, 1, 1, -2, 6, -1, -1, 2, 1, -1, 2, -1, 6, 1, -2, 1, 1, -2, 1, 1, 6, 2, -1, -1, 1, 1, 0, -1, -2, 1, 0, -1, 1, 0, 1, -1, -2, 0, 1, -1, 1, -1, 1, 0, -2, -1, 1, 0, 1, -1, 0, 1, -2, -1, 0, 1, 1, 0, -1, 1, -2, 0, -1, 1, 1, 1, -1, 0, -2, 1, -1, 0, 1, 1, 1, -2, -3, 1, 1, -2, 2, -1, 2, -1, -3, -1, 2, -1, 2, -2, 1, 1, -3, -2, 1, 1, 2, -1, -1, 2, -3, -1, -1, 2, 2, 1, -2, 1, -3, 1, -2, 1, 2, 2, -1, -1, -3, 2, -1, -1, 2 ], pReal), shape(HEX_SYSTEMTWIN))
 twin systems for hex, sorted by P. Eisenlohr CCW around starting next to a_1 axis More...
 
integer, dimension(13), parameter bct_nslipsystem = [2, 2, 2, 4, 2, 4, 2, 2, 4, 8, 4, 8, 8 ]
 
integer, parameter bct_nslip = sum(BCT_NSLIPSYSTEM)
 total # of slip systems for bct More...
 
real(preal), dimension(3+3, bct_nslip), parameter bct_systemslip = reshape(real([ 0, 0, 1, 1, 0, 0, 0, 0, 1, 0, 1, 0, 0, 0, 1, 1, 1, 0, 0, 0, 1, -1, 1, 0, 0, 1, 0, 1, 0, 0, 1, 0, 0, 0, 1, 0, 1,-1, 1, 1, 1, 0, 1,-1,-1, 1, 1, 0, -1,-1,-1, -1, 1, 0, -1,-1, 1, -1, 1, 0, 1, -1, 0, 1, 1, 0, 1, 1, 0, 1,-1, 0, 0, 1, 1, 1, 0, 0, 0,-1, 1, 1, 0, 0, -1, 0, 1, 0, 1, 0, 1, 0, 1, 0, 1, 0, 0, 1, 0, 0, 0, 1, 1, 0, 0, 0, 0, 1, 1, 1, 0, 0, 0, 1, -1, 1, 0, 0, 0, 1, 0, 1,-1, 0, 1, 1, 0,-1,-1, 0,-1, 1, -1, 0,-1, -1, 0, 1, 1, 0,-1, 1, 0, 1, 1,-1, 1, 0, 1, 1, 1, 1,-1, 0, 1, 1, 1, 1, 1, 0, 1,-1, -1, 1, 1, 0, 1,-1, 1,-1,-1, 1, 0, 1, -1,-1, 1, 1, 0, 1, 1, 1, 1, 1, 0,-1, 1,-1, 1, 1, 0,-1, 1, 0, 0, 0, 1, 1, 1, 0, 0, 0, 1,-1, 0, 1, 0, 1, 0, 1, 0, 1, 0, 1, 0,-1, 0, 1,-1, 2, 1, 1, 0,-1,-1, 2,-1, 1, 1, 0,-1, 1, 2, 1, -1, 0,-1, -1, 2, 1, 0, 1,-1, -2, 1, 1, 0,-1,-1, -2,-1, 1, -1, 0,-1, -1,-2, 1, 1, 0,-1, 1,-2, 1, -1, 1, 1, 2, 1, 1, -1,-1, 1, 2,-1, 1, 1,-1, 1, 1, 2, 1, -1,-1, 1, -1, 2, 1, 1, 1, 1, -2, 1, 1, 1,-1, 1, -2,-1, 1, -1, 1, 1, -1,-2, 1, 1, 1, 1, 1,-2, 1 ], pReal), shape(BCT_SYSTEMSLIP))
 slip systems for bct sorted by Bieler More...
 
integer, dimension(3), parameter ort_ncleavagesystem = [1, 1, 1]
 
integer, parameter ort_ncleavage = sum(ORT_NCLEAVAGESYSTEM)
 total # of cleavage systems for ortho More...
 
real(preal), dimension(3+3, ort_ncleavage), parameter ort_systemcleavage = reshape(real([ 0, 1, 0, 1, 0, 0, 0, 0, 1, 0, 1, 0, 1, 0, 0, 0, 0, 1 ], pReal), shape(ORT_SYSTEMCLEAVAGE))
 
@, public lattice_iso_id
 
@, public lattice_fcc_id
 
@, public lattice_bcc_id
 
@, public lattice_bct_id
 
@, public lattice_hex_id
 
@, public lattice_ort_id
 
real(preal), dimension(:), allocatable, public, protected lattice_mu
 
real(preal), dimension(:), allocatable, public, protected lattice_nu
 
real(preal), dimension(:), allocatable, public, protected lattice_damagemobility
 
real(preal), dimension(:), allocatable, public, protected lattice_massdensity
 
real(preal), dimension(:), allocatable, public, protected lattice_specificheat
 
real(preal), dimension(:,:,:), allocatable, public, protected lattice_c66
 
real(preal), dimension(:,:,:), allocatable, public, protected lattice_thermalconductivity
 
real(preal), dimension(:,:,:), allocatable, public, protected lattice_damagediffusion
 
integer(kind(lattice_undefined_id)), dimension(:), allocatable, public, protected lattice_structure
 

Detailed Description

contains lattice structure definitions including Schmid matrices for slip, twin, trans,

Author
Franz Roters, Max-Planck-Institut für Eisenforschung GmbH
Philip Eisenlohr, Max-Planck-Institut für Eisenforschung GmbH
Pratheek Shanthraj, Max-Planck-Institut für Eisenforschung GmbH
Martin Diehl, Max-Planck-Institut für Eisenforschung GmbH

Enumeration Type Documentation

◆ anonymous enum

anonymous enum
private
Enumerator
lattice_undefined_id 

Definition at line 12334 of file DAMASK_marc.f90.

Function/Subroutine Documentation

◆ applylatticesymmetryc66()

real(preal) function, dimension(6,6) lattice::applylatticesymmetryc66 ( real(preal), dimension(6,6), intent(in)  C66,
character(len=*), intent(in)  structure 
)
private

Return stiffness matrix in 6x6 notation with symmetry according to given crystal structure.

J. A. Rayne and B. S. Chandrasekhar Phys. Rev. 120, 1658 Erratum Phys. Rev. 122, 1962

Definition at line 13676 of file DAMASK_marc.f90.

References io::io_error().

Referenced by lattice_c66_trans(), lattice_init(), and unittest().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ buildcoordinatesystem()

real(preal) function, dimension(3,3,sum(active)) lattice::buildcoordinatesystem ( integer, dimension(:), intent(in)  active,
integer, dimension(:), intent(in)  potential,
real(preal), dimension(:,:), intent(in)  system,
character(len=*), intent(in)  structure,
real(preal), intent(in)  cOverA 
)
private

Build a local coordinate system on slip, twin, trans, cleavage systems.

Order: Direction, plane (normal), and common perpendicular

Parameters
[in]potential# of potential systems per family
[in]structurelattice structure
Parameters
active# of active systems per family

Definition at line 13933 of file DAMASK_marc.f90.

References io::io_error(), and math::math_cross().

Referenced by coordinatesystem_slip(), lattice_c66_twin(), lattice_nonschmidmatrix(), lattice_schmidmatrix_cleavage(), lattice_schmidmatrix_slip(), lattice_schmidmatrix_twin(), and unittest().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ buildinteraction()

real(preal) function, dimension(sum(reacting_used),sum(acting_used)) lattice::buildinteraction ( integer, dimension(:), intent(in)  reacting_used,
integer, dimension(:), intent(in)  acting_used,
integer, dimension(:), intent(in)  reacting_max,
integer, dimension(:), intent(in)  acting_max,
real(preal), dimension(:), intent(in)  values,
integer, dimension(:,:), intent(in)  matrix 
)
private

Populate reduced interaction matrix.

Parameters
[in]acting_maxmax # of acting systems per family for given lattice
[in]valuesinteraction values
[in]matrixinteraction types
Parameters
reacting_used# of reacting systems per family as specified in material.config
acting_used# of acting systems per family as specified in material.config
reacting_maxmax # of reacting systems per family for given lattice

Definition at line 13889 of file DAMASK_marc.f90.

References io::io_error().

Referenced by lattice_interaction_slipbyslip(), lattice_interaction_slipbytrans(), lattice_interaction_slipbytwin(), lattice_interaction_transbytrans(), lattice_interaction_twinbyslip(), and lattice_interaction_twinbytwin().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ buildtransformationsystem()

subroutine lattice::buildtransformationsystem ( real(preal), dimension(3,3,sum(ntrans)), intent(out)  Q,
real(preal), dimension(3,3,sum(ntrans)), intent(out)  S,
integer, dimension(:), intent(in)  Ntrans,
real(preal), intent(in)  cOverA,
real(preal), intent(in)  a_fcc,
real(preal), intent(in)  a_bcc 
)
private

Helper function to define transformation systems.

Parameters
[out]sEigendeformation tensor for phase transformation
[in]a_fcclattice parameter a for parent fcc structure
Parameters
QTotal rotation: Q = R*B
cOverAc/a for target hex structure
a_bcclattice parameter a for target bcc structure

Definition at line 14004 of file DAMASK_marc.f90.

References prec::deq0(), io::io_error(), math::math_cross(), math::math_i3, and math::math_outer().

Referenced by lattice_c66_trans(), and lattice_schmidmatrix_trans().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ coordinatesystem_slip()

real(preal) function, dimension(3,3,sum(nslip)) lattice::coordinatesystem_slip ( integer, dimension(:), intent(in)  Nslip,
character(len=*), intent(in)  structure,
real(preal), intent(in)  cOverA 
)
private

build a local coordinate system on slip systems

Order: Direction, plane (normal), and common perpendicular

Parameters
[in]nslipnumber of active slip systems per family
[in]structurelattice structure
[in]coverac/a ratio

Definition at line 13846 of file DAMASK_marc.f90.

References bcc_nslipsystem, bcc_systemslip, bct_nslipsystem, bct_systemslip, buildcoordinatesystem(), fcc_nslipsystem, fcc_systemslip, hex_nslipsystem, hex_systemslip, and io::io_error().

Referenced by lattice_slip_direction(), lattice_slip_normal(), and lattice_slip_transverse().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ equivalent_mu()

real(preal) function lattice::equivalent_mu ( real(preal), dimension(6,6), intent(in)  C,
character(len=*), intent(in)  assumption 
)
private

Equivalent shear modulus (μ)

https://doi.org/10.1143/JPSJ.20.635

Parameters
[in]cStiffness tensor (Voigt notation)
[in]assumptionAssumption ('Voigt' = isostrain, 'Reuss' = isostress)

Definition at line 14217 of file DAMASK_marc.f90.

References io::io_error(), io::io_lc(), and math::math_invert().

Referenced by equivalent_nu(), lattice_init(), and unittest().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ equivalent_nu()

real(preal) function lattice::equivalent_nu ( real(preal), dimension(6,6), intent(in)  C,
character(len=*), intent(in)  assumption 
)
private

Equivalent Poisson's ratio (ν)

https://doi.org/10.1143/JPSJ.20.635

Parameters
[in]cStiffness tensor (Voigt notation)
[in]assumptionAssumption ('Voigt' = isostrain, 'Reuss' = isostress)

Definition at line 14185 of file DAMASK_marc.f90.

References equivalent_mu(), io::io_error(), io::io_lc(), and math::math_invert().

Referenced by lattice_init(), and unittest().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ getlabels()

character(len=:) function, dimension(:), allocatable lattice::getlabels ( integer, dimension(:), intent(in)  active,
integer, dimension(:), intent(in)  potential,
real(preal), dimension(:,:), intent(in)  system 
)
private

select active systems as strings

Parameters
[in]potential# of potential systems per family
Parameters
active# of active systems per family

Definition at line 14128 of file DAMASK_marc.f90.

Referenced by lattice_labels_slip(), and lattice_labels_twin().

+ Here is the caller graph for this function:

◆ lattice_applylatticesymmetry33()

real(preal) function, dimension(3,3), public lattice::lattice_applylatticesymmetry33 ( real(preal), dimension(3,3), intent(in)  T,
character(len=*), intent(in)  structure 
)

Return 3x3 tensor with symmetry according to given crystal structure.

Definition at line 13638 of file DAMASK_marc.f90.

References io::io_error().

Referenced by kinematics_thermal_expansion::kinematics_thermal_expansion_init(), and lattice_init().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ lattice_c66_trans()

real(preal) function, dimension(6,6,sum(ntrans)), public lattice::lattice_c66_trans ( integer, dimension(:), intent(in)  Ntrans,
real(preal), dimension(6,6), intent(in)  C_parent66,
character(len=*), intent(in)  structure_target,
real(preal)  cOverA_trans,
real(preal)  a_bcc,
real(preal)  a_fcc 
)

Rotated elasticity matrices for transformation in 66-vector notation.

Parameters
[in]ntransnumber of active twin systems per family
[in]structure_targetlattice structure

Definition at line 12611 of file DAMASK_marc.f90.

References applylatticesymmetryc66(), buildtransformationsystem(), io::io_error(), and prec::tol_math_check.

+ Here is the call graph for this function:

◆ lattice_c66_twin()

real(preal) function, dimension(6,6,sum(ntwin)), public lattice::lattice_c66_twin ( integer, dimension(:), intent(in)  Ntwin,
real(preal), dimension(6,6), intent(in)  C66,
character(len=*), intent(in)  structure,
real(preal), intent(in)  CoverA 
)

Rotated elasticity matrices for twinning in 66-vector notation.

Parameters
[in]ntwinnumber of active twin systems per family
[in]structurelattice structure
[in]c66unrotated parent stiffness matrix
[in]coverac/a ratio

Definition at line 12570 of file DAMASK_marc.f90.

References bcc_nslipsystem, bcc_systemtwin, buildcoordinatesystem(), fcc_nslipsystem, fcc_systemtwin, hex_nslipsystem, hex_systemtwin, io::io_error(), and math::pi.

+ Here is the call graph for this function:

◆ lattice_characteristicshear_twin()

real(preal) function, dimension(sum(ntwin)), public lattice::lattice_characteristicshear_twin ( integer, dimension(:), intent(in)  Ntwin,
character(len=*), intent(in)  structure,
real(preal), intent(in)  CoverA 
)

Characteristic shear for twinning.

Parameters
[in]ntwinnumber of active twin systems per family
[in]structurelattice structure
[in]coverac/a ratio

Definition at line 12493 of file DAMASK_marc.f90.

References hex_ntwin, hex_ntwinsystem, and io::io_error().

+ Here is the call graph for this function:

◆ lattice_init()

subroutine, public lattice::lattice_init

Module initialization.

Definition at line 12403 of file DAMASK_marc.f90.

References applylatticesymmetryc66(), config::config_phase, equivalent_mu(), equivalent_nu(), io::io_error(), lattice_applylatticesymmetry33(), lattice_bcc_id, lattice_bct_id, lattice_c66, lattice_damagediffusion, lattice_damagemobility, lattice_fcc_id, lattice_hex_id, lattice_iso_id, lattice_massdensity, lattice_mu, lattice_nu, lattice_ort_id, lattice_specificheat, lattice_structure, lattice_thermalconductivity, lattice_undefined_id, math::math_sym3333to66(), math::math_voigt66to3333(), prec::tol_math_check, and prec::unittest().

Referenced by cpfem_initall().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ lattice_interaction_slipbyslip()

real(preal) function, dimension(sum(nslip),sum(nslip)), public lattice::lattice_interaction_slipbyslip ( integer, dimension(:), intent(in)  Nslip,
real(preal), dimension(:), intent(in)  interactionValues,
character(len=*), intent(in)  structure 
)

Slip-slip interaction matrix details only active slip systems are considered.

Parameters
[in]nslipnumber of active slip systems per family
[in]interactionvaluesvalues for slip-slip interaction
[in]structurelattice structure

Definition at line 12721 of file DAMASK_marc.f90.

References bcc_nslipsystem, bct_nslipsystem, buildinteraction(), fcc_nslipsystem, hex_nslipsystem, and io::io_error().

+ Here is the call graph for this function:

◆ lattice_interaction_slipbytrans()

real(preal) function, dimension(sum(nslip),sum(ntrans)), public lattice::lattice_interaction_slipbytrans ( integer, dimension(:), intent(in)  Nslip,
integer, dimension(:), intent(in)  Ntrans,
real(preal), dimension(:), intent(in)  interactionValues,
character(len=*), intent(in)  structure 
)

Slip-trans interaction matrix details only active slip and trans systems are considered.

Parameters
[in]ntransnumber of active trans systems per family
[in]interactionvaluesvalues for slip-trans interaction
[in]structurelattice structure (parent crystal)
Parameters
Nslipnumber of active slip systems per family

Definition at line 13227 of file DAMASK_marc.f90.

References buildinteraction(), fcc_nslipsystem, fcc_ntranssystem, and io::io_error().

+ Here is the call graph for this function:

◆ lattice_interaction_slipbytwin()

real(preal) function, dimension(sum(nslip),sum(ntwin)), public lattice::lattice_interaction_slipbytwin ( integer, dimension(:), intent(in)  Nslip,
integer, dimension(:), intent(in)  Ntwin,
real(preal), dimension(:), intent(in)  interactionValues,
character(len=*), intent(in)  structure 
)

Slip-twin interaction matrix details only active slip and twin systems are considered.

Parameters
[in]ntwinnumber of active twin systems per family
[in]interactionvaluesvalues for slip-twin interaction
[in]structurelattice structure
Parameters
Nslipnumber of active slip systems per family

Definition at line 13087 of file DAMASK_marc.f90.

References bcc_nslipsystem, bcc_ntwinsystem, buildinteraction(), fcc_nslipsystem, fcc_ntwinsystem, hex_nslipsystem, hex_ntwinsystem, and io::io_error().

+ Here is the call graph for this function:

◆ lattice_interaction_transbytrans()

real(preal) function, dimension(sum(ntrans),sum(ntrans)), public lattice::lattice_interaction_transbytrans ( integer, dimension(:), intent(in)  Ntrans,
real(preal), dimension(:), intent(in)  interactionValues,
character(len=*), intent(in)  structure 
)

Trans-trans interaction matrix details only active trans systems are considered.

Parameters
[in]ntransnumber of active trans systems per family
[in]interactionvaluesvalues for trans-trans interaction
[in]structurelattice structure (parent crystal)

Definition at line 13042 of file DAMASK_marc.f90.

References buildinteraction(), fcc_ntranssystem, and io::io_error().

+ Here is the call graph for this function:

◆ lattice_interaction_twinbyslip()

real(preal) function, dimension(sum(ntwin),sum(nslip)), public lattice::lattice_interaction_twinbyslip ( integer, dimension(:), intent(in)  Ntwin,
integer, dimension(:), intent(in)  Nslip,
real(preal), dimension(:), intent(in)  interactionValues,
character(len=*), intent(in)  structure 
)

Twin-slip interaction matrix details only active twin and slip systems are considered.

Parameters
[in]nslipnumber of active slip systems per family
[in]interactionvaluesvalues for twin-twin interaction
[in]structurelattice structure
Parameters
Ntwinnumber of active twin systems per family

Definition at line 13283 of file DAMASK_marc.f90.

References bcc_nslipsystem, bcc_ntwinsystem, buildinteraction(), fcc_nslipsystem, fcc_ntwinsystem, hex_nslipsystem, hex_ntwinsystem, and io::io_error().

+ Here is the call graph for this function:

◆ lattice_interaction_twinbytwin()

real(preal) function, dimension(sum(ntwin),sum(ntwin)), public lattice::lattice_interaction_twinbytwin ( integer, dimension(:), intent(in)  Ntwin,
real(preal), dimension(:), intent(in)  interactionValues,
character(len=*), intent(in)  structure 
)

Twin-twin interaction matrix details only active twin systems are considered.

Parameters
[in]ntwinnumber of active twin systems per family
[in]interactionvaluesvalues for twin-twin interaction
[in]structurelattice structure

Definition at line 12941 of file DAMASK_marc.f90.

References bcc_ntwinsystem, buildinteraction(), fcc_ntwinsystem, hex_ntwinsystem, and io::io_error().

+ Here is the call graph for this function:

◆ lattice_labels_slip()

character(len=:) function, dimension(:), allocatable, public lattice::lattice_labels_slip ( integer, dimension(:), intent(in)  Nslip,
character(len=*), intent(in)  structure 
)

Labels for slip systems details only active slip systems are considered.

Parameters
[in]nslipnumber of active slip systems per family
[in]structurelattice structure

Definition at line 13595 of file DAMASK_marc.f90.

References bcc_nslipsystem, bcc_systemslip, bct_nslipsystem, bct_systemslip, fcc_nslipsystem, fcc_systemslip, getlabels(), hex_nslipsystem, hex_systemslip, and io::io_error().

+ Here is the call graph for this function:

◆ lattice_labels_twin()

character(len=:) function, dimension(:), allocatable, public lattice::lattice_labels_twin ( integer, dimension(:), intent(in)  Ntwin,
character(len=*), intent(in)  structure 
)

Labels for twin systems details only active twin systems are considered.

Parameters
[in]ntwinnumber of active slip systems per family
[in]structurelattice structure

Definition at line 13757 of file DAMASK_marc.f90.

References bcc_ntwinsystem, bcc_systemtwin, fcc_ntwinsystem, fcc_systemtwin, getlabels(), hex_ntwinsystem, hex_systemtwin, and io::io_error().

+ Here is the call graph for this function:

◆ lattice_nonschmidmatrix()

real(preal) function, dimension(1:3,1:3,sum(nslip)), public lattice::lattice_nonschmidmatrix ( integer, dimension(:), intent(in)  Nslip,
real(preal), dimension(:), intent(in)  nonSchmidCoefficients,
integer, intent(in)  sense 
)

Non-schmid projections for bcc with up to 6 coefficients.

Parameters
[in]nslipnumber of active slip systems per family
[in]nonschmidcoefficientsnon-Schmid coefficients for projections
[in]sensesense (-1,+1)

Definition at line 12674 of file DAMASK_marc.f90.

References bcc_nslipsystem, bcc_systemslip, buildcoordinatesystem(), io::io_error(), lattice_schmidmatrix_slip(), math::math_cross(), math::math_outer(), and prec::preal.

+ Here is the call graph for this function:

◆ lattice_schmidmatrix_cleavage()

real(preal) function, dimension(3,3,3,sum(ncleavage)), public lattice::lattice_schmidmatrix_cleavage ( integer, dimension(:), intent(in)  Ncleavage,
character(len=*), intent(in)  structure,
real(preal), intent(in)  cOverA 
)

Schmid matrix for cleavage details only active cleavage systems are considered.

Parameters
[in]ncleavagenumber of active cleavage systems per family
[in]structurelattice structure
[in]coverac/a ratio

Definition at line 13492 of file DAMASK_marc.f90.

References bcc_ncleavagesystem, bcc_systemcleavage, buildcoordinatesystem(), fcc_ncleavagesystem, fcc_systemcleavage, io::io_error(), math::math_outer(), ort_ncleavagesystem, and ort_systemcleavage.

Referenced by kinematics_cleavage_opening::kinematics_cleavage_opening_init(), and source_damage_anisobrittle::source_damage_anisobrittle_init().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ lattice_schmidmatrix_slip()

real(preal) function, dimension(3,3,sum(nslip)), public lattice::lattice_schmidmatrix_slip ( integer, dimension(:), intent(in)  Nslip,
character(len=*), intent(in)  structure,
real(preal), intent(in)  cOverA 
)

Schmid matrix for slip details only active slip systems are considered.

Parameters
[in]nslipnumber of active slip systems per family
[in]structurelattice structure

Definition at line 13361 of file DAMASK_marc.f90.

References bcc_nslipsystem, bcc_systemslip, bct_nslipsystem, bct_systemslip, buildcoordinatesystem(), fcc_nslipsystem, fcc_systemslip, hex_nslipsystem, hex_systemslip, io::io_error(), math::math_outer(), math::math_trace33(), and prec::tol_math_check.

Referenced by lattice_nonschmidmatrix().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ lattice_schmidmatrix_trans()

real(preal) function, dimension(3,3,sum(ntrans)), public lattice::lattice_schmidmatrix_trans ( integer, dimension(:), intent(in)  Ntrans,
character(len=*), intent(in)  structure_target,
real(preal), intent(in)  cOverA,
real(preal)  a_bcc,
real(preal)  a_fcc 
)

Schmid matrix for twinning details only active twin systems are considered.

Parameters
[in]ntransnumber of active twin systems per family
[in]structure_targetlattice structure
[in]coverac/a ratio

Definition at line 13462 of file DAMASK_marc.f90.

References buildtransformationsystem(), and io::io_error().

+ Here is the call graph for this function:

◆ lattice_schmidmatrix_twin()

real(preal) function, dimension(3,3,sum(ntwin)), public lattice::lattice_schmidmatrix_twin ( integer, dimension(:), intent(in)  Ntwin,
character(len=*), intent(in)  structure,
real(preal), intent(in)  cOverA 
)

Schmid matrix for twinning details only active twin systems are considered.

Parameters
[in]ntwinnumber of active twin systems per family
[in]structurelattice structure
[in]coverac/a ratio

Definition at line 13413 of file DAMASK_marc.f90.

References bcc_ntwinsystem, bcc_systemtwin, buildcoordinatesystem(), fcc_ntwinsystem, fcc_systemtwin, hex_ntwinsystem, hex_systemtwin, io::io_error(), math::math_outer(), math::math_trace33(), and prec::tol_math_check.

+ Here is the call graph for this function:

◆ lattice_slip_direction()

real(preal) function, dimension(3,sum(nslip)), public lattice::lattice_slip_direction ( integer, dimension(:), intent(in)  Nslip,
character(len=*), intent(in)  structure,
real(preal), intent(in)  cOverA 
)

Slip direction of slip systems (|| b)

Parameters
[in]nslipnumber of active slip systems per family
[in]structurelattice structure
[in]coverac/a ratio

Definition at line 13540 of file DAMASK_marc.f90.

References coordinatesystem_slip().

Referenced by kinematics_slipplane_opening::kinematics_slipplane_opening_init(), and slipprojection_direction().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ lattice_slip_normal()

real(preal) function, dimension(3,sum(nslip)), public lattice::lattice_slip_normal ( integer, dimension(:), intent(in)  Nslip,
character(len=*), intent(in)  structure,
real(preal), intent(in)  cOverA 
)

Normal direction of slip systems (|| n)

Parameters
[in]nslipnumber of active slip systems per family
[in]structurelattice structure
[in]coverac/a ratio

Definition at line 13558 of file DAMASK_marc.f90.

References coordinatesystem_slip().

Referenced by kinematics_slipplane_opening::kinematics_slipplane_opening_init(), slipprojection_direction(), and slipprojection_transverse().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ lattice_slip_transverse()

real(preal) function, dimension(3,sum(nslip)), public lattice::lattice_slip_transverse ( integer, dimension(:), intent(in)  Nslip,
character(len=*), intent(in)  structure,
real(preal), intent(in)  cOverA 
)

Transverse direction of slip systems ( || t = b x n)

Parameters
[in]nslipnumber of active slip systems per family
[in]structurelattice structure
[in]coverac/a ratio

Definition at line 13576 of file DAMASK_marc.f90.

References coordinatesystem_slip().

Referenced by kinematics_slipplane_opening::kinematics_slipplane_opening_init(), and slipprojection_transverse().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ slipprojection_direction()

real(preal) function, dimension(sum(nslip),sum(nslip)) lattice::slipprojection_direction ( integer, dimension(:), intent(in)  Nslip,
character(len=*), intent(in)  structure,
real(preal), intent(in)  cOverA 
)
private

Projection of the slip direction onto the slip plane.

: This projection is used to calculate forest hardening for screw dislocations

Parameters
[in]nslipnumber of active slip systems per family
[in]structurelattice structure
[in]coverac/a ratio

Definition at line 13822 of file DAMASK_marc.f90.

References lattice_slip_direction(), lattice_slip_normal(), and math::math_inner().

+ Here is the call graph for this function:

◆ slipprojection_transverse()

real(preal) function, dimension(sum(nslip),sum(nslip)) lattice::slipprojection_transverse ( integer, dimension(:), intent(in)  Nslip,
character(len=*), intent(in)  structure,
real(preal), intent(in)  cOverA 
)
private

Projection of the transverse direction onto the slip plane.

: This projection is used to calculate forest hardening for edge dislocations

Parameters
[in]nslipnumber of active slip systems per family
[in]structurelattice structure
[in]coverac/a ratio

Definition at line 13798 of file DAMASK_marc.f90.

References lattice_slip_normal(), lattice_slip_transverse(), and math::math_inner().

+ Here is the call graph for this function:

◆ unittest()

subroutine lattice::unittest
private

check correctness of some lattice functions

Definition at line 14245 of file DAMASK_marc.f90.

References applylatticesymmetryc66(), buildcoordinatesystem(), prec::dneq(), equivalent_mu(), equivalent_nu(), io::io_error(), and math::math_i3.

+ Here is the call graph for this function:

Variable Documentation

◆ bcc_ncleavage

integer, parameter lattice::bcc_ncleavage = sum(BCC_NCLEAVAGESYSTEM)
private

total # of cleavage systems for bcc

Definition at line 12071 of file DAMASK_marc.f90.

◆ bcc_ncleavagesystem

integer, dimension(1), parameter lattice::bcc_ncleavagesystem = [3]
private

of cleavage systems per family for bcc

Definition at line 12068 of file DAMASK_marc.f90.

Referenced by lattice_schmidmatrix_cleavage().

◆ bcc_nslip

integer, parameter lattice::bcc_nslip = sum(BCC_NSLIPSYSTEM)
private

total # of slip systems for bcc

Definition at line 12071 of file DAMASK_marc.f90.

◆ bcc_nslipsystem

integer, dimension(2), parameter lattice::bcc_nslipsystem = [12, 12]
private

◆ bcc_ntwin

integer, parameter lattice::bcc_ntwin = sum(BCC_NTWINSYSTEM)
private

total # of twin systems for bcc

Definition at line 12071 of file DAMASK_marc.f90.

◆ bcc_ntwinsystem

integer, dimension(1), parameter lattice::bcc_ntwinsystem = [12]
private

◆ bcc_systemcleavage

real(preal), dimension(3+3,bcc_ncleavage), parameter lattice::bcc_systemcleavage = reshape(real([ 0, 1, 0, 1, 0, 0, 0, 0, 1, 0, 1, 0, 1, 0, 0, 0, 0, 1 ], pReal), shape(BCC_SYSTEMCLEAVAGE))
private

Definition at line 12130 of file DAMASK_marc.f90.

Referenced by lattice_schmidmatrix_cleavage().

◆ bcc_systemslip

real(preal), dimension(3+3,bcc_nslip), parameter lattice::bcc_systemslip = reshape(real([ 1,-1, 1, 0, 1, 1, -1,-1, 1, 0, 1, 1, 1, 1, 1, 0,-1, 1, -1, 1, 1, 0,-1, 1, -1, 1, 1, 1, 0, 1, -1,-1, 1, 1, 0, 1, 1, 1, 1, -1, 0, 1, 1,-1, 1, -1, 0, 1, -1, 1, 1, 1, 1, 0, -1, 1,-1, 1, 1, 0, 1, 1, 1, -1, 1, 0, 1, 1,-1, -1, 1, 0, -1, 1, 1, 2, 1, 1, 1, 1, 1, -2, 1, 1, 1, 1,-1, 2,-1, 1, 1,-1, 1, 2, 1,-1, 1,-1, 1, 1, 2, 1, 1, 1,-1, -1, 2, 1, 1, 1, 1, 1,-2, 1, -1, 1, 1, 1, 2,-1, 1, 1,-1, 1, 1, 2, 1,-1, 1, -1, 1, 2, -1, 1, 1, 1,-1, 2, 1, 1, 1, 1, 1,-2 ], pReal), shape(BCC_SYSTEMSLIP))
private

◆ bcc_systemtwin

real(preal), dimension(3+3,bcc_ntwin), parameter lattice::bcc_systemtwin = reshape(real([ -1, 1, 1, 2, 1, 1, 1, 1, 1, -2, 1, 1, 1, 1,-1, 2,-1, 1, 1,-1, 1, 2, 1,-1, 1,-1, 1, 1, 2, 1, 1, 1,-1, -1, 2, 1, 1, 1, 1, 1,-2, 1, -1, 1, 1, 1, 2,-1, 1, 1,-1, 1, 1, 2, 1,-1, 1, -1, 1, 2, -1, 1, 1, 1,-1, 2, 1, 1, 1, 1, 1,-2 ], pReal), shape(BCC_SYSTEMTWIN))
private

◆ bct_nslip

integer, parameter lattice::bct_nslip = sum(BCT_NSLIPSYSTEM)
private

total # of slip systems for bct

Definition at line 12236 of file DAMASK_marc.f90.

◆ bct_nslipsystem

integer, dimension(13), parameter lattice::bct_nslipsystem = [2, 2, 2, 4, 2, 4, 2, 2, 4, 8, 4, 8, 8 ]
private

of slip systems per family for bct (Sn) Bieler J. Electr Mater 2009

Definition at line 12233 of file DAMASK_marc.f90.

Referenced by coordinatesystem_slip(), lattice_interaction_slipbyslip(), lattice_labels_slip(), and lattice_schmidmatrix_slip().

◆ bct_systemslip

real(preal), dimension(3+3,bct_nslip), parameter lattice::bct_systemslip = reshape(real([ 0, 0, 1, 1, 0, 0, 0, 0, 1, 0, 1, 0, 0, 0, 1, 1, 1, 0, 0, 0, 1, -1, 1, 0, 0, 1, 0, 1, 0, 0, 1, 0, 0, 0, 1, 0, 1,-1, 1, 1, 1, 0, 1,-1,-1, 1, 1, 0, -1,-1,-1, -1, 1, 0, -1,-1, 1, -1, 1, 0, 1, -1, 0, 1, 1, 0, 1, 1, 0, 1,-1, 0, 0, 1, 1, 1, 0, 0, 0,-1, 1, 1, 0, 0, -1, 0, 1, 0, 1, 0, 1, 0, 1, 0, 1, 0, 0, 1, 0, 0, 0, 1, 1, 0, 0, 0, 0, 1, 1, 1, 0, 0, 0, 1, -1, 1, 0, 0, 0, 1, 0, 1,-1, 0, 1, 1, 0,-1,-1, 0,-1, 1, -1, 0,-1, -1, 0, 1, 1, 0,-1, 1, 0, 1, 1,-1, 1, 0, 1, 1, 1, 1,-1, 0, 1, 1, 1, 1, 1, 0, 1,-1, -1, 1, 1, 0, 1,-1, 1,-1,-1, 1, 0, 1, -1,-1, 1, 1, 0, 1, 1, 1, 1, 1, 0,-1, 1,-1, 1, 1, 0,-1, 1, 0, 0, 0, 1, 1, 1, 0, 0, 0, 1,-1, 0, 1, 0, 1, 0, 1, 0, 1, 0, 1, 0,-1, 0, 1,-1, 2, 1, 1, 0,-1,-1, 2,-1, 1, 1, 0,-1, 1, 2, 1, -1, 0,-1, -1, 2, 1, 0, 1,-1, -2, 1, 1, 0,-1,-1, -2,-1, 1, -1, 0,-1, -1,-2, 1, 1, 0,-1, 1,-2, 1, -1, 1, 1, 2, 1, 1, -1,-1, 1, 2,-1, 1, 1,-1, 1, 1, 2, 1, -1,-1, 1, -1, 2, 1, 1, 1, 1, -2, 1, 1, 1,-1, 1, -2,-1, 1, -1, 1, 1, -1,-2, 1, 1, 1, 1, 1,-2, 1 ], pReal), shape(BCT_SYSTEMSLIP))
private

slip systems for bct sorted by Bieler

Definition at line 12243 of file DAMASK_marc.f90.

Referenced by coordinatesystem_slip(), lattice_labels_slip(), and lattice_schmidmatrix_slip().

◆ fcc_ncleavage

integer, parameter lattice::fcc_ncleavage = sum(FCC_NCLEAVAGESYSTEM)
private

total # of cleavage systems for fcc

Definition at line 11983 of file DAMASK_marc.f90.

◆ fcc_ncleavagesystem

integer, dimension(1), parameter lattice::fcc_ncleavagesystem = [3]
private

of cleavage systems per family for fcc

Definition at line 11980 of file DAMASK_marc.f90.

Referenced by lattice_schmidmatrix_cleavage().

◆ fcc_nslip

integer, parameter lattice::fcc_nslip = sum(FCC_NSLIPSYSTEM)
private

total # of slip systems for fcc

Definition at line 11983 of file DAMASK_marc.f90.

◆ fcc_nslipsystem

integer, dimension(2), parameter lattice::fcc_nslipsystem = [12, 6]
private

◆ fcc_ntrans

integer, parameter lattice::fcc_ntrans = sum(FCC_NTRANSSYSTEM)
private

total # of transformation systems for fcc

Definition at line 11983 of file DAMASK_marc.f90.

◆ fcc_ntranssystem

integer, dimension(1), parameter lattice::fcc_ntranssystem = [12]
private

of transformation systems per family for fcc

Definition at line 11977 of file DAMASK_marc.f90.

Referenced by lattice_interaction_slipbytrans(), and lattice_interaction_transbytrans().

◆ fcc_ntwin

integer, parameter lattice::fcc_ntwin = sum(FCC_NTWINSYSTEM)
private

total # of twin systems for fcc

Definition at line 11983 of file DAMASK_marc.f90.

◆ fcc_ntwinsystem

integer, dimension(1), parameter lattice::fcc_ntwinsystem = [12]
private

◆ fcc_systemcleavage

real(preal), dimension(3+3,fcc_ncleavage), parameter lattice::fcc_systemcleavage = reshape(real([ 0, 1, 0, 1, 0, 0, 0, 0, 1, 0, 1, 0, 1, 0, 0, 0, 0, 1 ], pReal), shape(FCC_SYSTEMCLEAVAGE))
private

Definition at line 12052 of file DAMASK_marc.f90.

Referenced by lattice_schmidmatrix_cleavage().

◆ fcc_systemslip

real(preal), dimension(3+3,fcc_nslip), parameter lattice::fcc_systemslip = reshape(real([ 0, 1,-1, 1, 1, 1, -1, 0, 1, 1, 1, 1, 1,-1, 0, 1, 1, 1, 0,-1,-1, -1,-1, 1, 1, 0, 1, -1,-1, 1, -1, 1, 0, -1,-1, 1, 0,-1, 1, 1,-1,-1, -1, 0,-1, 1,-1,-1, 1, 1, 0, 1,-1,-1, 0, 1, 1, -1, 1,-1, 1, 0,-1, -1, 1,-1, -1,-1, 0, -1, 1,-1, 1, 1, 0, 1,-1, 0, 1,-1, 0, 1, 1, 0, 1, 0, 1, 1, 0,-1, 1, 0,-1, 1, 0, 1, 0, 1, 1, 0, 1,-1, 0, 1,-1, 0, 1, 1 ], pReal), shape(FCC_SYSTEMSLIP))
private

Slip system <110>{111} directions. Sorted according to Eisenlohr & Hantcherli.

Definition at line 11996 of file DAMASK_marc.f90.

Referenced by coordinatesystem_slip(), lattice_labels_slip(), and lattice_schmidmatrix_slip().

◆ fcc_systemtwin

real(preal), dimension(3+3,fcc_ntwin), parameter lattice::fcc_systemtwin = reshape(real( [ -2, 1, 1, 1, 1, 1, 1,-2, 1, 1, 1, 1, 1, 1,-2, 1, 1, 1, 2,-1, 1, -1,-1, 1, -1, 2, 1, -1,-1, 1, -1,-1,-2, -1,-1, 1, -2,-1,-1, 1,-1,-1, 1, 2,-1, 1,-1,-1, 1,-1, 2, 1,-1,-1, 2, 1,-1, -1, 1,-1, -1,-2,-1, -1, 1,-1, -1, 1, 2, -1, 1,-1 ], pReal), shape(FCC_SYSTEMTWIN))
private

Twin system <112>{111} directions. Sorted according to Eisenlohr & Hantcherli.

Definition at line 12020 of file DAMASK_marc.f90.

Referenced by lattice_c66_twin(), lattice_labels_twin(), and lattice_schmidmatrix_twin().

◆ hex_nslip

integer, parameter lattice::hex_nslip = sum(HEX_NSLIPSYSTEM)
private

total # of slip systems for hex

Definition at line 12146 of file DAMASK_marc.f90.

◆ hex_nslipsystem

integer, dimension(6), parameter lattice::hex_nslipsystem = [3, 3, 3, 6, 12, 6]
private

◆ hex_ntwin

integer, parameter lattice::hex_ntwin = sum(HEX_NTWINSYSTEM)
private

total # of twin systems for hex

Definition at line 12146 of file DAMASK_marc.f90.

Referenced by lattice_characteristicshear_twin().

◆ hex_ntwinsystem

integer, dimension(4), parameter lattice::hex_ntwinsystem = [6, 6, 6, 6]
private

◆ hex_systemslip

real(preal), dimension(4+4,hex_nslip), parameter lattice::hex_systemslip = reshape(real([ 2, -1, -1, 0, 0, 0, 0, 1, -1, 2, -1, 0, 0, 0, 0, 1, -1, -1, 2, 0, 0, 0, 0, 1, 2, -1, -1, 0, 0, 1, -1, 0, -1, 2, -1, 0, -1, 0, 1, 0, -1, -1, 2, 0, 1, -1, 0, 0, -1, 1, 0, 0, 1, 1, -2, 0, 0, -1, 1, 0, -2, 1, 1, 0, 1, 0, -1, 0, 1, -2, 1, 0, -1, 2, -1, 0, 1, 0, -1, 1, -2, 1, 1, 0, 0, 1, -1, 1, -1, -1, 2, 0, -1, 1, 0, 1, 1, -2, 1, 0, -1, 0, 1, 1, 2, -1, -1, 0, 0, -1, 1, 1, 1, 1, -2, 0, 1, -1, 0, 1, -2, 1, 1, 3, 1, 0, -1, 1, -1, -1, 2, 3, 1, 0, -1, 1, -1, -1, 2, 3, 0, 1, -1, 1, 1, -2, 1, 3, 0, 1, -1, 1, 1, -2, 1, 3, -1, 1, 0, 1, 2, -1, -1, 3, -1, 1, 0, 1, 2, -1, -1, 3, -1, 0, 1, 1, 1, 1, -2, 3, -1, 0, 1, 1, 1, 1, -2, 3, 0, -1, 1, 1, -1, 2, -1, 3, 0, -1, 1, 1, -1, 2, -1, 3, 1, -1, 0, 1, -2, 1, 1, 3, 1, -1, 0, 1, -1, -1, 2, 3, 1, 1, -2, 2, 1, -2, 1, 3, -1, 2, -1, 2, 2, -1, -1, 3, -2, 1, 1, 2, 1, 1, -2, 3, -1, -1, 2, 2, -1, 2, -1, 3, 1, -2, 1, 2, -2, 1, 1, 3, 2, -1, -1, 2 ], pReal), shape(HEX_SYSTEMSLIP))
private

slip systems for hex, sorted by P. Eisenlohr CCW around starting next to a_1 axis

Definition at line 12155 of file DAMASK_marc.f90.

Referenced by coordinatesystem_slip(), lattice_labels_slip(), and lattice_schmidmatrix_slip().

◆ hex_systemtwin

real(preal), dimension(4+4,hex_ntwin), parameter lattice::hex_systemtwin = reshape(real([ -1, 0, 1, 1, 1, 0, -1, 2, 0, -1, 1, 1, 0, 1, -1, 2, 1, -1, 0, 1, -1, 1, 0, 2, 1, 0, -1, 1, -1, 0, 1, 2, 0, 1, -1, 1, 0, -1, 1, 2, -1, 1, 0, 1, 1, -1, 0, 2, -1, -1, 2, 6, 1, 1, -2, 1, 1, -2, 1, 6, -1, 2, -1, 1, 2, -1, -1, 6, -2, 1, 1, 1, 1, 1, -2, 6, -1, -1, 2, 1, -1, 2, -1, 6, 1, -2, 1, 1, -2, 1, 1, 6, 2, -1, -1, 1, 1, 0, -1, -2, 1, 0, -1, 1, 0, 1, -1, -2, 0, 1, -1, 1, -1, 1, 0, -2, -1, 1, 0, 1, -1, 0, 1, -2, -1, 0, 1, 1, 0, -1, 1, -2, 0, -1, 1, 1, 1, -1, 0, -2, 1, -1, 0, 1, 1, 1, -2, -3, 1, 1, -2, 2, -1, 2, -1, -3, -1, 2, -1, 2, -2, 1, 1, -3, -2, 1, 1, 2, -1, -1, 2, -3, -1, -1, 2, 2, 1, -2, 1, -3, 1, -2, 1, 2, 2, -1, -1, -3, 2, -1, -1, 2 ], pReal), shape(HEX_SYSTEMTWIN))
private

twin systems for hex, sorted by P. Eisenlohr CCW around starting next to a_1 axis

Definition at line 12199 of file DAMASK_marc.f90.

Referenced by lattice_c66_twin(), lattice_labels_twin(), and lattice_schmidmatrix_twin().

◆ lattice_bcc_id

@, public lattice::lattice_bcc_id

Definition at line 12341 of file DAMASK_marc.f90.

Referenced by crystallite_results(), and lattice_init().

◆ lattice_bct_id

@, public lattice::lattice_bct_id

Definition at line 12341 of file DAMASK_marc.f90.

Referenced by crystallite_results(), and lattice_init().

◆ lattice_c66

real(preal), dimension(:,:,:), allocatable, public, protected lattice::lattice_c66

Definition at line 12350 of file DAMASK_marc.f90.

Referenced by constitutive::constitutive_homogenizedc(), and lattice_init().

◆ lattice_damagediffusion

real(preal), dimension(:,:,:), allocatable, public, protected lattice::lattice_damagediffusion

Definition at line 12350 of file DAMASK_marc.f90.

Referenced by lattice_init().

◆ lattice_damagemobility

real(preal), dimension(:), allocatable, public, protected lattice::lattice_damagemobility

Definition at line 12345 of file DAMASK_marc.f90.

Referenced by lattice_init().

◆ lattice_fcc_id

@, public lattice::lattice_fcc_id

Definition at line 12341 of file DAMASK_marc.f90.

Referenced by crystallite_results(), and lattice_init().

◆ lattice_fcc_twinnucleationslippair

integer, dimension(2,fcc_ntwin), parameter, public lattice::lattice_fcc_twinnucleationslippair = reshape( [ 2,3, 1,3, 1,2, 5,6, 4,6, 4,5, 8,9, 7,9, 7,8, 11,12, 10,12, 10,11 ], shape(lattice_FCC_TWINNUCLEATIONSLIPPAIR))

Definition at line 12036 of file DAMASK_marc.f90.

◆ lattice_hex_id

@, public lattice::lattice_hex_id

Definition at line 12341 of file DAMASK_marc.f90.

Referenced by crystallite_results(), and lattice_init().

◆ lattice_iso_id

@, public lattice::lattice_iso_id

Definition at line 12341 of file DAMASK_marc.f90.

Referenced by crystallite_results(), and lattice_init().

◆ lattice_massdensity

real(preal), dimension(:), allocatable, public, protected lattice::lattice_massdensity

Definition at line 12345 of file DAMASK_marc.f90.

Referenced by lattice_init().

◆ lattice_mu

real(preal), dimension(:), allocatable, public, protected lattice::lattice_mu

Definition at line 12345 of file DAMASK_marc.f90.

Referenced by lattice_init().

◆ lattice_nu

real(preal), dimension(:), allocatable, public, protected lattice::lattice_nu

Definition at line 12345 of file DAMASK_marc.f90.

Referenced by lattice_init().

◆ lattice_ort_id

@, public lattice::lattice_ort_id

Definition at line 12341 of file DAMASK_marc.f90.

Referenced by crystallite_results(), and lattice_init().

◆ lattice_specificheat

real(preal), dimension(:), allocatable, public, protected lattice::lattice_specificheat

Definition at line 12345 of file DAMASK_marc.f90.

Referenced by lattice_init().

◆ lattice_structure

integer(kind(lattice_undefined_id)), dimension(:), allocatable, public, protected lattice::lattice_structure

Definition at line 12354 of file DAMASK_marc.f90.

Referenced by crystallite_results(), and lattice_init().

◆ lattice_thermalconductivity

real(preal), dimension(:,:,:), allocatable, public, protected lattice::lattice_thermalconductivity

Definition at line 12350 of file DAMASK_marc.f90.

Referenced by lattice_init().

◆ ort_ncleavage

integer, parameter lattice::ort_ncleavage = sum(ORT_NCLEAVAGESYSTEM)
private

total # of cleavage systems for ortho

Definition at line 12318 of file DAMASK_marc.f90.

◆ ort_ncleavagesystem

integer, dimension(3), parameter lattice::ort_ncleavagesystem = [1, 1, 1]
private

of cleavage systems per family for ortho

Definition at line 12315 of file DAMASK_marc.f90.

Referenced by lattice_schmidmatrix_cleavage().

◆ ort_systemcleavage

real(preal), dimension(3+3,ort_ncleavage), parameter lattice::ort_systemcleavage = reshape(real([ 0, 1, 0, 1, 0, 0, 0, 0, 1, 0, 1, 0, 1, 0, 0, 0, 0, 1 ], pReal), shape(ORT_SYSTEMCLEAVAGE))
private

Definition at line 12325 of file DAMASK_marc.f90.

Referenced by lattice_schmidmatrix_cleavage().