www.mooseframework.org
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Private Types | Private Attributes | List of all members
InteractionIntegralSM Class Reference

This postprocessor computes the Interaction Integral. More...

#include <InteractionIntegralSM.h>

Inheritance diagram for InteractionIntegralSM:
[legend]

Public Member Functions

 InteractionIntegralSM (const InputParameters &parameters)
 
virtual Real getValue ()
 

Static Public Member Functions

static MooseEnum qFunctionType ()
 
static MooseEnum sifModeType ()
 

Protected Member Functions

virtual void initialSetup ()
 
virtual Real computeQpIntegral ()
 
virtual Real computeIntegral ()
 
void computeAuxFields (RankTwoTensor &aux_stress, RankTwoTensor &grad_disp)
 
void computeTFields (RankTwoTensor &aux_stress, RankTwoTensor &grad_disp)
 

Protected Attributes

unsigned int _ndisp
 
const CrackFrontDefinition *const _crack_front_definition
 
bool _has_crack_front_point_index
 
const unsigned int _crack_front_point_index
 
bool _treat_as_2d
 
const MaterialProperty< SymmTensor > & _stress
 
const MaterialProperty< SymmTensor > & _strain
 
std::vector< const VariableGradient * > _grad_disp
 
const bool _has_temp
 
const VariableGradient & _grad_temp
 
const MaterialProperty< Real > * _current_instantaneous_thermal_expansion_coef
 
Real _K_factor
 
bool _has_symmetry_plane
 
Real _poissons_ratio
 
Real _youngs_modulus
 
unsigned int _ring_index
 
std::vector< Real > _q_curr_elem
 
const std::vector< std::vector< Real > > * _phi_curr_elem
 
const std::vector< std::vector< RealGradient > > * _dphi_curr_elem
 
Real _kappa
 
Real _shear_modulus
 
Real _r
 
Real _theta
 

Private Types

enum  QMethod { QMethod::Geometry, QMethod::Topology }
 
enum  SifMethod { SifMethod::KI, SifMethod::KII, SifMethod::KIII, SifMethod::T }
 

Private Attributes

const QMethod _q_function_type
 
const SifMethod _sif_mode
 

Detailed Description

This postprocessor computes the Interaction Integral.

Definition at line 28 of file InteractionIntegralSM.h.

Member Enumeration Documentation

◆ QMethod

enum InteractionIntegralSM::QMethod
strongprivate
Enumerator
Geometry 
Topology 

Definition at line 70 of file InteractionIntegralSM.h.

71  {
72  Geometry,
73  Topology
74  };

◆ SifMethod

Enumerator
KI 
KII 
KIII 

Definition at line 78 of file InteractionIntegralSM.h.

79  {
80  KI,
81  KII,
82  KIII,
83  T
84  };

Constructor & Destructor Documentation

◆ InteractionIntegralSM()

InteractionIntegralSM::InteractionIntegralSM ( const InputParameters &  parameters)

Definition at line 73 of file InteractionIntegralSM.C.

74  : ElementIntegralPostprocessor(parameters),
75  _ndisp(coupledComponents("displacements")),
76  _crack_front_definition(&getUserObject<CrackFrontDefinition>("crack_front_definition")),
77  _has_crack_front_point_index(isParamValid("crack_front_point_index")),
79  _has_crack_front_point_index ? getParam<unsigned int>("crack_front_point_index") : 0),
80  _treat_as_2d(false),
81  _stress(getMaterialPropertyByName<SymmTensor>("stress")),
82  _strain(getMaterialPropertyByName<SymmTensor>("elastic_strain")),
83  _grad_disp(3),
84  _has_temp(isCoupled("temperature")),
85  _grad_temp(_has_temp ? coupledGradient("temperature") : _grad_zero),
87  hasMaterialProperty<Real>("current_instantaneous_thermal_expansion_coef")
88  ? &getMaterialProperty<Real>("current_instantaneous_thermal_expansion_coef")
89  : NULL),
90  _K_factor(getParam<Real>("K_factor")),
91  _has_symmetry_plane(isParamValid("symmetry_plane")),
92  _poissons_ratio(getParam<Real>("poissons_ratio")),
93  _youngs_modulus(getParam<Real>("youngs_modulus")),
94  _ring_index(getParam<unsigned int>("ring_index")),
95  _q_function_type(getParam<MooseEnum>("q_function_type").getEnum<QMethod>()),
96  _sif_mode(getParam<MooseEnum>("sif_mode").getEnum<SifMethod>())
97 {
99  mooseError("To include thermal strain term in interaction integral, must both couple "
100  "temperature in DomainIntegral block and compute thermal expansion property in "
101  "material model using compute_InteractionIntegral = true.");
102 
103  // plane strain
104  _kappa = 3.0 - 4.0 * _poissons_ratio;
105  _shear_modulus = _youngs_modulus / (2.0 * (1.0 + _poissons_ratio));
106 
107  // Checking for consistency between mesh size and length of the provided displacements vector
108  if (_ndisp != _mesh.dimension())
109  mooseError(
110  "The number of variables supplied in 'displacements' must match the mesh dimension.");
111  // fetch gradients of coupled variables
112  for (unsigned int i = 0; i < _ndisp; ++i)
113  _grad_disp[i] = &coupledGradient("displacements", i);
114 
115  // set unused dimensions to zero
116  for (unsigned i = _ndisp; i < 3; ++i)
117  _grad_disp[i] = &_grad_zero;
118 }
const unsigned int _crack_front_point_index
const CrackFrontDefinition *const _crack_front_definition
const MaterialProperty< SymmTensor > & _stress
const MaterialProperty< Real > * _current_instantaneous_thermal_expansion_coef
std::vector< const VariableGradient * > _grad_disp
const VariableGradient & _grad_temp
const MaterialProperty< SymmTensor > & _strain

Member Function Documentation

◆ computeAuxFields()

void InteractionIntegralSM::computeAuxFields ( RankTwoTensor &  aux_stress,
RankTwoTensor &  grad_disp 
)
protected

Definition at line 291 of file InteractionIntegralSM.C.

Referenced by computeQpIntegral().

292 {
293  RealVectorValue k(0.0);
294  if (_sif_mode == SifMethod::KI)
295  k(0) = 1.0;
296  else if (_sif_mode == SifMethod::KII)
297  k(1) = 1.0;
298  else if (_sif_mode == SifMethod::KIII)
299  k(2) = 1.0;
300 
301  Real t = _theta;
302  Real t2 = _theta / 2.0;
303  Real tt2 = 3.0 * _theta / 2.0;
304  Real st = std::sin(t);
305  Real ct = std::cos(t);
306  Real st2 = std::sin(t2);
307  Real ct2 = std::cos(t2);
308  Real stt2 = std::sin(tt2);
309  Real ctt2 = std::cos(tt2);
310  Real ct2sq = Utility::pow<2>(ct2);
311  Real ct2cu = Utility::pow<3>(ct2);
312  Real sqrt2PiR = std::sqrt(2.0 * libMesh::pi * _r);
313 
314  // Calculate auxiliary stress tensor
315  aux_stress.zero();
316 
317  aux_stress(0, 0) =
318  1.0 / sqrt2PiR * (k(0) * ct2 * (1.0 - st2 * stt2) - k(1) * st2 * (2.0 + ct2 * ctt2));
319  aux_stress(1, 1) = 1.0 / sqrt2PiR * (k(0) * ct2 * (1.0 + st2 * stt2) + k(1) * st2 * ct2 * ctt2);
320  aux_stress(0, 1) = 1.0 / sqrt2PiR * (k(0) * ct2 * st2 * ctt2 + k(1) * ct2 * (1.0 - st2 * stt2));
321  aux_stress(0, 2) = -1.0 / sqrt2PiR * k(2) * st2;
322  aux_stress(1, 2) = 1.0 / sqrt2PiR * k(2) * ct2;
323  // plane stress
324  // Real s33 = 0;
325  // plane strain
326  aux_stress(2, 2) = _poissons_ratio * (aux_stress(0, 0) + aux_stress(1, 1));
327 
328  aux_stress(1, 0) = aux_stress(0, 1);
329  aux_stress(2, 0) = aux_stress(0, 2);
330  aux_stress(2, 1) = aux_stress(1, 2);
331 
332  // Calculate x1 derivative of auxiliary displacements
333  grad_disp.zero();
334 
335  grad_disp(0, 0) = k(0) / (4.0 * _shear_modulus * sqrt2PiR) *
336  (ct * ct2 * _kappa + ct * ct2 - 2.0 * ct * ct2cu + st * st2 * _kappa +
337  st * st2 - 6.0 * st * st2 * ct2sq) +
338  k(1) / (4.0 * _shear_modulus * sqrt2PiR) *
339  (ct * st2 * _kappa + ct * st2 + 2.0 * ct * st2 * ct2sq - st * ct2 * _kappa +
340  3.0 * st * ct2 - 6.0 * st * ct2cu);
341 
342  grad_disp(0, 1) = k(0) / (4.0 * _shear_modulus * sqrt2PiR) *
343  (ct * st2 * _kappa + ct * st2 - 2.0 * ct * st2 * ct2sq - st * ct2 * _kappa -
344  5.0 * st * ct2 + 6.0 * st * ct2cu) +
345  k(1) / (4.0 * _shear_modulus * sqrt2PiR) *
346  (-ct * ct2 * _kappa + 3.0 * ct * ct2 - 2.0 * ct * ct2cu -
347  st * st2 * _kappa + 3.0 * st * st2 - 6.0 * st * st2 * ct2sq);
348 
349  grad_disp(0, 2) = k(2) / (_shear_modulus * sqrt2PiR) * (st2 * ct - ct2 * st);
350 }

◆ computeIntegral()

Real InteractionIntegralSM::computeIntegral ( )
protectedvirtual

Definition at line 252 of file InteractionIntegralSM.C.

253 {
254  Real sum = 0.0;
255 
256  // calculate phi and dphi for this element
257  FEType fe_type(Utility::string_to_enum<Order>("first"),
258  Utility::string_to_enum<FEFamily>("lagrange"));
259  const unsigned int dim = _current_elem->dim();
260  std::unique_ptr<FEBase> fe(FEBase::build(dim, fe_type));
261  fe->attach_quadrature_rule(_qrule);
262  _phi_curr_elem = &fe->get_phi();
263  _dphi_curr_elem = &fe->get_dphi();
264  fe->reinit(_current_elem);
265 
266  // calculate q for all nodes in this element
267  _q_curr_elem.clear();
268  unsigned int ring_base = (_q_function_type == QMethod::Topology) ? 0 : 1;
269 
270  for (unsigned int i = 0; i < _current_elem->n_nodes(); ++i)
271  {
272  Node * this_node = _current_elem->get_node(i);
273  Real q_this_node;
274 
277  _crack_front_point_index, _ring_index - ring_base, this_node);
280  _crack_front_point_index, _ring_index - ring_base, this_node);
281 
282  _q_curr_elem.push_back(q_this_node);
283  }
284 
285  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
286  sum += _JxW[_qp] * _coord[_qp] * computeQpIntegral();
287  return sum;
288 }
const unsigned int _crack_front_point_index
const CrackFrontDefinition *const _crack_front_definition
const std::vector< std::vector< Real > > * _phi_curr_elem
Real DomainIntegralTopologicalQFunction(unsigned int crack_front_point_index, unsigned int ring_index, const Node *const current_node) const
std::vector< Real > _q_curr_elem
const std::vector< std::vector< RealGradient > > * _dphi_curr_elem
Real DomainIntegralQFunction(unsigned int crack_front_point_index, unsigned int ring_index, const Node *const current_node) const

◆ computeQpIntegral()

Real InteractionIntegralSM::computeQpIntegral ( )
protectedvirtual

Definition at line 140 of file InteractionIntegralSM.C.

Referenced by computeIntegral().

141 {
142  Real scalar_q = 0.0;
143  RealVectorValue grad_q(0.0, 0.0, 0.0);
144 
145  const std::vector<std::vector<Real>> & phi_curr_elem = *_phi_curr_elem;
146  const std::vector<std::vector<RealGradient>> & dphi_curr_elem = *_dphi_curr_elem;
147 
148  for (unsigned int i = 0; i < _current_elem->n_nodes(); ++i)
149  {
150  scalar_q += phi_curr_elem[i][_qp] * _q_curr_elem[i];
151 
152  for (unsigned int j = 0; j < _current_elem->dim(); ++j)
153  grad_q(j) += dphi_curr_elem[i][_qp](j) * _q_curr_elem[i];
154  }
155 
156  // In the crack front coordinate system, the crack direction is (1,0,0)
157  RealVectorValue crack_direction(0.0);
158  crack_direction(0) = 1.0;
159 
160  // Calculate (r,theta) position of qp relative to crack front
161  Point p(_q_point[_qp]);
163 
164  RankTwoTensor aux_stress;
165  RankTwoTensor aux_du;
166 
168  computeAuxFields(aux_stress, aux_du);
169  else if (_sif_mode == SifMethod::T)
170  computeTFields(aux_stress, aux_du);
171 
172  RankTwoTensor stress;
173  stress(0, 0) = _stress[_qp].xx();
174  stress(0, 1) = _stress[_qp].xy();
175  stress(0, 2) = _stress[_qp].xz();
176  stress(1, 0) = _stress[_qp].xy();
177  stress(1, 1) = _stress[_qp].yy();
178  stress(1, 2) = _stress[_qp].yz();
179  stress(2, 0) = _stress[_qp].xz();
180  stress(2, 1) = _stress[_qp].yz();
181  stress(2, 2) = _stress[_qp].zz();
182 
183  RankTwoTensor strain;
184  strain(0, 0) = _strain[_qp].xx();
185  strain(0, 1) = _strain[_qp].xy();
186  strain(0, 2) = _strain[_qp].xz();
187  strain(1, 0) = _strain[_qp].xy();
188  strain(1, 1) = _strain[_qp].yy();
189  strain(1, 2) = _strain[_qp].yz();
190  strain(2, 0) = _strain[_qp].xz();
191  strain(2, 1) = _strain[_qp].yz();
192  strain(2, 2) = _strain[_qp].zz();
193 
194  RankTwoTensor grad_disp((*_grad_disp[0])[_qp], (*_grad_disp[1])[_qp], (*_grad_disp[2])[_qp]);
195 
196  // Rotate stress, strain, displacement and temperature to crack front coordinate system
197  RealVectorValue grad_q_cf =
199  RankTwoTensor grad_disp_cf =
201  RankTwoTensor stress_cf =
203  RankTwoTensor strain_cf =
205  RealVectorValue grad_temp_cf =
207 
208  RankTwoTensor dq;
209  dq(0, 0) = crack_direction(0) * grad_q_cf(0);
210  dq(0, 1) = crack_direction(0) * grad_q_cf(1);
211  dq(0, 2) = crack_direction(0) * grad_q_cf(2);
212 
213  // Calculate interaction integral terms
214 
215  // Term1 = stress * x1-derivative of aux disp * dq
216  RankTwoTensor tmp1 = dq * stress_cf;
217  Real term1 = aux_du.doubleContraction(tmp1);
218 
219  // Term2 = aux stress * x1-derivative of disp * dq
220  RankTwoTensor tmp2 = dq * aux_stress;
221  Real term2 = grad_disp_cf(0, 0) * tmp2(0, 0) + grad_disp_cf(1, 0) * tmp2(0, 1) +
222  grad_disp_cf(2, 0) * tmp2(0, 2);
223 
224  // Term3 = aux stress * strain * dq_x (= stress * aux strain * dq_x)
225  Real term3 = dq(0, 0) * aux_stress.doubleContraction(strain_cf);
226 
227  // Term4 (thermal strain term) = q * aux_stress * alpha * dtheta_x
228  // - the term including the derivative of alpha is not implemented
229  Real term4 = 0.0;
230  if (_has_temp)
231  term4 = scalar_q * aux_stress.trace() * (*_current_instantaneous_thermal_expansion_coef)[_qp] *
232  grad_temp_cf(0);
233 
234  Real q_avg_seg = 1.0;
236  {
237  q_avg_seg =
240  2.0;
241  }
242 
243  Real eq = term1 + term2 - term3 + term4;
244 
246  eq *= 2.0;
247 
248  return eq / q_avg_seg;
249 }
const unsigned int _crack_front_point_index
void calculateRThetaToCrackFront(const Point qp, const unsigned int point_index, Real &r, Real &theta) const
calculate r and theta in the crack front polar cooridnate
Real getCrackFrontBackwardSegmentLength(const unsigned int point_index) const
void computeAuxFields(RankTwoTensor &aux_stress, RankTwoTensor &grad_disp)
const CrackFrontDefinition *const _crack_front_definition
void computeTFields(RankTwoTensor &aux_stress, RankTwoTensor &grad_disp)
const MaterialProperty< SymmTensor > & _stress
const std::vector< std::vector< Real > > * _phi_curr_elem
RealVectorValue rotateToCrackFrontCoords(const RealVectorValue vector, const unsigned int point_index) const
Real getCrackFrontForwardSegmentLength(const unsigned int point_index) const
std::vector< const VariableGradient * > _grad_disp
std::vector< Real > _q_curr_elem
const std::vector< std::vector< RealGradient > > * _dphi_curr_elem
const VariableGradient & _grad_temp
const MaterialProperty< SymmTensor > & _strain

◆ computeTFields()

void InteractionIntegralSM::computeTFields ( RankTwoTensor &  aux_stress,
RankTwoTensor &  grad_disp 
)
protected

Definition at line 353 of file InteractionIntegralSM.C.

Referenced by computeQpIntegral().

354 {
355 
356  Real t = _theta;
357  Real st = std::sin(t);
358  Real ct = std::cos(t);
359  Real stsq = Utility::pow<2>(st);
360  Real ctsq = Utility::pow<2>(ct);
361  Real ctcu = Utility::pow<3>(ct);
362  Real oneOverPiR = 1.0 / (libMesh::pi * _r);
363 
364  aux_stress.zero();
365  aux_stress(0, 0) = -oneOverPiR * ctcu;
366  aux_stress(0, 1) = -oneOverPiR * st * ctsq;
367  aux_stress(1, 0) = -oneOverPiR * st * ctsq;
368  aux_stress(1, 1) = -oneOverPiR * ct * stsq;
369  aux_stress(2, 2) = -oneOverPiR * _poissons_ratio * (ctcu + ct * stsq);
370 
371  grad_disp.zero();
372  grad_disp(0, 0) = oneOverPiR / (4.0 * _youngs_modulus) *
373  (ct * (4.0 * Utility::pow<2>(_poissons_ratio) - 3.0 + _poissons_ratio) -
374  std::cos(3.0 * t) * (1.0 + _poissons_ratio));
375  grad_disp(0, 1) = -oneOverPiR / (4.0 * _youngs_modulus) *
376  (st * (4.0 * Utility::pow<2>(_poissons_ratio) - 3.0 + _poissons_ratio) +
377  std::sin(3.0 * t) * (1.0 + _poissons_ratio));
378 }

◆ getValue()

Real InteractionIntegralSM::getValue ( )
virtual

Definition at line 127 of file InteractionIntegralSM.C.

128 {
129  gatherSum(_integral_value);
130 
132  _integral_value +=
135 
136  return _K_factor * _integral_value;
137 }
const unsigned int _crack_front_point_index
const CrackFrontDefinition *const _crack_front_definition
Real getCrackFrontTangentialStrain(const unsigned int node_index) const

◆ initialSetup()

void InteractionIntegralSM::initialSetup ( )
protectedvirtual

Definition at line 121 of file InteractionIntegralSM.C.

122 {
124 }
const CrackFrontDefinition *const _crack_front_definition

◆ qFunctionType()

MooseEnum InteractionIntegralSM::qFunctionType ( )
static

Definition at line 20 of file InteractionIntegralSM.C.

Referenced by validParams< InteractionIntegralSM >().

21 {
22  return MooseEnum("Geometry Topology", "Geometry");
23 }

◆ sifModeType()

MooseEnum InteractionIntegralSM::sifModeType ( )
static

Definition at line 26 of file InteractionIntegralSM.C.

Referenced by validParams< InteractionIntegralSM >().

27 {
28  return MooseEnum("KI KII KIII T", "KI");
29 }

Member Data Documentation

◆ _crack_front_definition

const CrackFrontDefinition* const InteractionIntegralSM::_crack_front_definition
protected

Definition at line 46 of file InteractionIntegralSM.h.

Referenced by computeIntegral(), computeQpIntegral(), getValue(), and initialSetup().

◆ _crack_front_point_index

const unsigned int InteractionIntegralSM::_crack_front_point_index
protected

Definition at line 48 of file InteractionIntegralSM.h.

Referenced by computeIntegral(), computeQpIntegral(), and getValue().

◆ _current_instantaneous_thermal_expansion_coef

const MaterialProperty<Real>* InteractionIntegralSM::_current_instantaneous_thermal_expansion_coef
protected

Definition at line 55 of file InteractionIntegralSM.h.

Referenced by InteractionIntegralSM().

◆ _dphi_curr_elem

const std::vector<std::vector<RealGradient> >* InteractionIntegralSM::_dphi_curr_elem
protected

Definition at line 63 of file InteractionIntegralSM.h.

Referenced by computeIntegral(), and computeQpIntegral().

◆ _grad_disp

std::vector<const VariableGradient *> InteractionIntegralSM::_grad_disp
protected

Definition at line 52 of file InteractionIntegralSM.h.

Referenced by computeQpIntegral(), and InteractionIntegralSM().

◆ _grad_temp

const VariableGradient& InteractionIntegralSM::_grad_temp
protected

Definition at line 54 of file InteractionIntegralSM.h.

Referenced by computeQpIntegral().

◆ _has_crack_front_point_index

bool InteractionIntegralSM::_has_crack_front_point_index
protected

Definition at line 47 of file InteractionIntegralSM.h.

◆ _has_symmetry_plane

bool InteractionIntegralSM::_has_symmetry_plane
protected

Definition at line 57 of file InteractionIntegralSM.h.

Referenced by computeQpIntegral().

◆ _has_temp

const bool InteractionIntegralSM::_has_temp
protected

Definition at line 53 of file InteractionIntegralSM.h.

Referenced by computeQpIntegral(), and InteractionIntegralSM().

◆ _K_factor

Real InteractionIntegralSM::_K_factor
protected

Definition at line 56 of file InteractionIntegralSM.h.

Referenced by getValue().

◆ _kappa

Real InteractionIntegralSM::_kappa
protected

Definition at line 64 of file InteractionIntegralSM.h.

Referenced by computeAuxFields(), and InteractionIntegralSM().

◆ _ndisp

unsigned int InteractionIntegralSM::_ndisp
protected

Definition at line 45 of file InteractionIntegralSM.h.

Referenced by InteractionIntegralSM().

◆ _phi_curr_elem

const std::vector<std::vector<Real> >* InteractionIntegralSM::_phi_curr_elem
protected

Definition at line 62 of file InteractionIntegralSM.h.

Referenced by computeIntegral(), and computeQpIntegral().

◆ _poissons_ratio

Real InteractionIntegralSM::_poissons_ratio
protected

◆ _q_curr_elem

std::vector<Real> InteractionIntegralSM::_q_curr_elem
protected

Definition at line 61 of file InteractionIntegralSM.h.

Referenced by computeIntegral(), and computeQpIntegral().

◆ _q_function_type

const QMethod InteractionIntegralSM::_q_function_type
private

Definition at line 76 of file InteractionIntegralSM.h.

Referenced by computeIntegral().

◆ _r

Real InteractionIntegralSM::_r
protected

Definition at line 66 of file InteractionIntegralSM.h.

Referenced by computeAuxFields(), computeQpIntegral(), and computeTFields().

◆ _ring_index

unsigned int InteractionIntegralSM::_ring_index
protected

Definition at line 60 of file InteractionIntegralSM.h.

Referenced by computeIntegral().

◆ _shear_modulus

Real InteractionIntegralSM::_shear_modulus
protected

Definition at line 65 of file InteractionIntegralSM.h.

Referenced by computeAuxFields(), and InteractionIntegralSM().

◆ _sif_mode

const SifMethod InteractionIntegralSM::_sif_mode
private

Definition at line 86 of file InteractionIntegralSM.h.

Referenced by computeAuxFields(), computeQpIntegral(), and getValue().

◆ _strain

const MaterialProperty<SymmTensor>& InteractionIntegralSM::_strain
protected

Definition at line 51 of file InteractionIntegralSM.h.

Referenced by computeQpIntegral().

◆ _stress

const MaterialProperty<SymmTensor>& InteractionIntegralSM::_stress
protected

Definition at line 50 of file InteractionIntegralSM.h.

Referenced by computeQpIntegral().

◆ _theta

Real InteractionIntegralSM::_theta
protected

Definition at line 67 of file InteractionIntegralSM.h.

Referenced by computeAuxFields(), computeQpIntegral(), and computeTFields().

◆ _treat_as_2d

bool InteractionIntegralSM::_treat_as_2d
protected

Definition at line 49 of file InteractionIntegralSM.h.

Referenced by getValue(), and initialSetup().

◆ _youngs_modulus

Real InteractionIntegralSM::_youngs_modulus
protected

Definition at line 59 of file InteractionIntegralSM.h.

Referenced by computeTFields(), and InteractionIntegralSM().


The documentation for this class was generated from the following files: