https://mooseframework.inl.gov
Public Types | Public Member Functions | Static Public Member Functions | Public Attributes | Static Public Attributes | Protected Types | Protected Member Functions | Static Protected Member Functions | Protected Attributes | List of all members
DomainIntegralAction Class Reference

Action to set up all objects used in computation of fracture domain integrals. More...

#include <DomainIntegralAction.h>

Inheritance diagram for DomainIntegralAction:
[legend]

Public Types

typedef DataFileName DataFileParameterType
 

Public Member Functions

 DomainIntegralAction (const InputParameters &params)
 
 ~DomainIntegralAction ()
 
virtual void act () override
 
virtual void addRelationshipManagers (Moose::RelationshipManagerType input_rm_type) override
 
virtual void addRelationshipManagers (Moose::RelationshipManagerType when_type)
 
bool addRelationshipManagers (Moose::RelationshipManagerType when_type, const InputParameters &moose_object_pars)
 
void timedAct ()
 
MooseObjectName uniqueActionName () const
 
const std::string & specificTaskName () const
 
const std::set< std::string > & getAllTasks () const
 
void appendTask (const std::string &task)
 
MooseAppgetMooseApp () const
 
const std::string & type () const
 
const std::string & name () const
 
std::string typeAndName () const
 
MooseObjectParameterName uniqueParameterName (const std::string &parameter_name) const
 
MooseObjectName uniqueName () const
 
const InputParametersparameters () const
 
const hit::Node * getHitNode () const
 
bool hasBase () const
 
const std::string & getBase () const
 
const T & getParam (const std::string &name) const
 
std::vector< std::pair< T1, T2 > > getParam (const std::string &param1, const std::string &param2) const
 
const T * queryParam (const std::string &name) const
 
const T & getRenamedParam (const std::string &old_name, const std::string &new_name) const
 
getCheckedPointerParam (const std::string &name, const std::string &error_string="") const
 
bool isParamValid (const std::string &name) const
 
bool isParamSetByUser (const std::string &name) const
 
void connectControllableParams (const std::string &parameter, const std::string &object_type, const std::string &object_name, const std::string &object_parameter) const
 
void paramError (const std::string &param, Args... args) const
 
void paramWarning (const std::string &param, Args... args) const
 
void paramInfo (const std::string &param, Args... args) const
 
std::string messagePrefix (const bool hit_prefix=true) const
 
std::string errorPrefix (const std::string &) const
 
void mooseError (Args &&... args) const
 
void mooseDocumentedError (const std::string &repo_name, const unsigned int issue_num, Args &&... args) const
 
void mooseErrorNonPrefixed (Args &&... args) const
 
void mooseWarning (Args &&... args) const
 
void mooseWarningNonPrefixed (Args &&... args) const
 
void mooseDeprecated (Args &&... args) const
 
void mooseInfo (Args &&... args) const
 
void callMooseError (std::string msg, const bool with_prefix, const hit::Node *node=nullptr) const
 
std::string getDataFileName (const std::string &param) const
 
std::string getDataFileNameByName (const std::string &relative_path) const
 
std::string getDataFilePath (const std::string &relative_path) const
 
PerfGraphperfGraph ()
 
const Parallel::Communicator & comm () const
 
processor_id_type n_processors () const
 
processor_id_type processor_id () const
 

Static Public Member Functions

static InputParameters validParams ()
 
static void callMooseError (MooseApp *const app, const InputParameters &params, std::string msg, const bool with_prefix, const hit::Node *node)
 

Public Attributes

const ConsoleStream _console
 

Static Public Attributes

static const std::string unique_action_name_param
 
static const std::string type_param
 
static const std::string name_param
 
static const std::string unique_name_param
 
static const std::string app_param
 
static const std::string moose_base_param
 
static const std::string kokkos_object_param
 
static constexpr auto SYSTEM
 
static constexpr auto NAME
 

Protected Types

enum  INTEGRAL {
  J_INTEGRAL, C_INTEGRAL, K_FROM_J_INTEGRAL, INTERACTION_INTEGRAL_KI,
  INTERACTION_INTEGRAL_KII, INTERACTION_INTEGRAL_KIII, INTERACTION_INTEGRAL_T
}
 Enum used to select the type of integral to be performed. More...
 
enum  Q_FUNCTION_TYPE { GEOMETRY, TOPOLOGY }
 Enum used to select the method used to compute the Q function used in the fracture integrals. More...
 

Protected Member Functions

unsigned int calcNumCrackFrontPoints ()
 Compute the number of points on the crack front. More...
 
bool addRelationshipManagers (Moose::RelationshipManagerType when_type, const InputParameters &moose_object_pars)
 
void associateWithParameter (const std::string &param_name, InputParameters &params) const
 
void associateWithParameter (const InputParameters &from_params, const std::string &param_name, InputParameters &params) const
 
const T & getMeshProperty (const std::string &data_name, const std::string &prefix)
 
const T & getMeshProperty (const std::string &data_name)
 
bool hasMeshProperty (const std::string &data_name, const std::string &prefix) const
 
bool hasMeshProperty (const std::string &data_name, const std::string &prefix) const
 
bool hasMeshProperty (const std::string &data_name) const
 
bool hasMeshProperty (const std::string &data_name) const
 
std::string meshPropertyName (const std::string &data_name) const
 
PerfID registerTimedSection (const std::string &section_name, const unsigned int level) const
 
PerfID registerTimedSection (const std::string &section_name, const unsigned int level, const std::string &live_message, const bool print_dots=true) const
 
std::string timedSectionName (const std::string &section_name) const
 

Static Protected Member Functions

static std::string meshPropertyName (const std::string &data_name, const std::string &prefix)
 

Protected Attributes

std::set< INTEGRAL_integrals
 Container for enumerations describing the individual integrals computed. More...
 
const std::vector< BoundaryName > & _boundary_names
 Boundaries containing the crack front points. More...
 
std::vector< Point > _crack_front_points
 User-defined vector of crack front points. More...
 
bool _closed_loop
 Indicates whether the crack forms a closed loop. More...
 
UserObjectName _crack_front_points_provider
 Name of crack front points provider user object used to optionally define the crack points. More...
 
bool _use_crack_front_points_provider
 Whether to use a crack front points provider. More...
 
MooseEnum _direction_method_moose_enum
 Enum used to define the method to compute crack front direction. More...
 
MooseEnum _end_direction_method_moose_enum
 Enum used to define the method to compute crack front direction at ends of crack front. More...
 
const bool _have_crack_direction_vector
 Whether the crack direction vector has been provided. More...
 
const RealVectorValue _crack_direction_vector
 Vector optionally used to prescribe direction of crack extension. More...
 
const bool _have_crack_direction_vector_end_1
 Whether the crack direction vector at the 1st end of the crack has been provided. More...
 
const RealVectorValue _crack_direction_vector_end_1
 Vector optionally used to prescribe direction of crack extension at the 1st end of the crack. More...
 
const bool _have_crack_direction_vector_end_2
 Whether the crack direction vector at the 2nd end of the crack has been provided. More...
 
const RealVectorValue _crack_direction_vector_end_2
 Vector optionally used to prescribe direction of crack extension at the 2nd end of the crack. More...
 
std::vector< BoundaryName > _crack_mouth_boundary_names
 Names of boundaries optionally used to define the crack mouth location. More...
 
std::vector< BoundaryName > _intersecting_boundary_names
 Names of boundaries optionally used for improved computation of crack extension direction at ends of the crack where it intersects the prescribed boundaries. More...
 
bool _treat_as_2d
 Whether fracture computations for a 3D model should be treated as though it were a 2D model. More...
 
unsigned int _axis_2d
 Out-of-plane axis for 3D models treated as 2D. More...
 
unsigned int _ring_first
 Number of elements away from the crack tip to inside of inner ring with the topological q function. More...
 
unsigned int _ring_last
 Number of elements away from the crack tip to outside of outer ring with the topological q function. More...
 
std::vector< VariableName > _output_variables
 List of variables for which values are to be sampled and output at the crack front points. More...
 
Real _poissons_ratio
 Poisson's ratio of material. More...
 
Real _youngs_modulus
 Young's modulus of material. More...
 
std::vector< SubdomainName > _blocks
 Blocks for which the domain integrals are to be computed. More...
 
std::vector< VariableName > _displacements
 Vector of displacement variables. More...
 
VariableName _temp
 Temperature variable. More...
 
bool _has_symmetry_plane
 Whether the model has a symmetry plane passing through the plane of the crack. More...
 
unsigned int _symmetry_plane
 Identifier for which plane is the symmetry plane. More...
 
MooseEnum _position_type
 How the distance along the crack front is measured (angle or distance) More...
 
MooseEnum _q_function_type
 How the q function is evaluated (geometric distance from crack front or ring of elements) More...
 
bool _get_equivalent_k
 Whether to compute the equivalent K from the individual fracture integrals for mixed-mode fracture. More...
 
bool _use_displaced_mesh
 Whether to compute the fracture integrals on the displaced mesh. More...
 
bool _output_q
 Whether to ouput the q function as a set of AuxVariables. More...
 
std::vector< unsigned int_ring_vec
 Vector of ids for the individual rings on which the fracture integral is computed. More...
 
bool _incremental
 Whether the constitutive models for the mechanics calculations use an incremental form. More...
 
bool _convert_J_to_K
 Whether to convert the J-integral to a stress intensity factor (K) –deprecated. More...
 
bool _fgm_crack
 Whether the crack lives in a functionally-graded material. More...
 
MaterialPropertyName _functionally_graded_youngs_modulus_crack_dir_gradient
 Material property name for the Youngs modulus derivative for functionally graded materials. More...
 
MaterialPropertyName _functionally_graded_youngs_modulus
 Material property name for spatially-dependent Youngs modulus for functionally graded materials. More...
 
const bool _use_ad
 Whether to create automatic differentiation objects from the action. More...
 
std::string _registered_identifier
 
std::string _specific_task_name
 
std::set< std::string > _all_tasks
 
ActionWarehouse_awh
 
const std::string & _current_task
 
std::shared_ptr< MooseMesh > & _mesh
 
std::shared_ptr< MooseMesh > & _displaced_mesh
 
std::shared_ptr< FEProblemBase > & _problem
 
PerfID _act_timer
 
MooseApp_app
 
Factory_factory
 
ActionFactory_action_factory
 
const std::string & _type
 
const std::string & _name
 
const InputParameters_pars
 
MooseApp_pg_moose_app
 
const std::string _prefix
 
const Parallel::Communicator & _communicator
 
const std::string _order
 Order and family of the AuxVariables optionally created to output the values of q. More...
 
const std::string _family
 
std::vector< Real_radius_inner
 Sets of inner and outer radii of the rings used for the domain form of the fracture integrals. More...
 
std::vector< Real_radius_outer
 

Detailed Description

Action to set up all objects used in computation of fracture domain integrals.

Definition at line 22 of file DomainIntegralAction.h.

Member Enumeration Documentation

◆ INTEGRAL

◆ Q_FUNCTION_TYPE

Enum used to select the method used to compute the Q function used in the fracture integrals.

Enumerator
GEOMETRY 
TOPOLOGY 

Definition at line 51 of file DomainIntegralAction.h.

Constructor & Destructor Documentation

◆ DomainIntegralAction()

DomainIntegralAction::DomainIntegralAction ( const InputParameters params)

Definition at line 143 of file DomainIntegralAction.C.

144  : Action(params),
145  _boundary_names(getParam<std::vector<BoundaryName>>("boundary")),
146  _closed_loop(getParam<bool>("closed_loop")),
148  _order(getParam<std::string>("order")),
149  _family(getParam<std::string>("family")),
150  _direction_method_moose_enum(getParam<MooseEnum>("crack_direction_method")),
151  _end_direction_method_moose_enum(getParam<MooseEnum>("crack_end_direction_method")),
152  _have_crack_direction_vector(isParamValid("crack_direction_vector")),
154  _have_crack_direction_vector ? getParam<RealVectorValue>("crack_direction_vector") : 0.0),
155  _have_crack_direction_vector_end_1(isParamValid("crack_direction_vector_end_1")),
157  ? getParam<RealVectorValue>("crack_direction_vector_end_1")
158  : 0.0),
159  _have_crack_direction_vector_end_2(isParamValid("crack_direction_vector_end_2")),
161  ? getParam<RealVectorValue>("crack_direction_vector_end_2")
162  : 0.0),
163  _treat_as_2d(getParam<bool>("2d")),
164  _axis_2d(getParam<unsigned int>("axis_2d")),
165  _has_symmetry_plane(isParamValid("symmetry_plane")),
166  _symmetry_plane(_has_symmetry_plane ? getParam<unsigned int>("symmetry_plane")
167  : std::numeric_limits<unsigned int>::max()),
168  _position_type(getParam<MooseEnum>("position_type")),
169  _q_function_type(getParam<MooseEnum>("q_function_type")),
170  _get_equivalent_k(getParam<bool>("equivalent_k")),
171  _use_displaced_mesh(false),
172  _output_q(getParam<bool>("output_q")),
173  _incremental(getParam<bool>("incremental")),
174  _convert_J_to_K(isParamValid("convert_J_to_K") ? getParam<bool>("convert_J_to_K") : false),
175  _fgm_crack(false),
176  _use_ad(getParam<bool>("use_automatic_differentiation"))
177 {
178 
179  if (isParamValid("functionally_graded_youngs_modulus_crack_dir_gradient") !=
180  isParamValid("functionally_graded_youngs_modulus"))
181  paramError("functionally_graded_youngs_modulus_crack_dir_gradient",
182  "You have selected to compute the interaction integral for a crack in FGM. That "
183  "selection requires the user to provide a spatially varying elasticity modulus that "
184  "defines the transition of material properties (i.e. "
185  "'functionally_graded_youngs_modulus') and its "
186  "spatial derivative in the crack direction (i.e. "
187  "'functionally_graded_youngs_modulus_crack_dir_gradient').");
188 
189  if (isParamValid("functionally_graded_youngs_modulus_crack_dir_gradient") &&
190  isParamValid("functionally_graded_youngs_modulus"))
191  {
192  _fgm_crack = true;
194  getParam<MaterialPropertyName>("functionally_graded_youngs_modulus_crack_dir_gradient");
196  getParam<MaterialPropertyName>("functionally_graded_youngs_modulus");
197  }
198 
199  if (_q_function_type == GEOMETRY)
200  {
201  if (isParamValid("radius_inner") && isParamValid("radius_outer"))
202  {
203  _radius_inner = getParam<std::vector<Real>>("radius_inner");
204  _radius_outer = getParam<std::vector<Real>>("radius_outer");
205  }
206  else
207  mooseError("DomainIntegral error: must set radius_inner and radius_outer.");
208  for (unsigned int i = 0; i < _radius_inner.size(); ++i)
209  _ring_vec.push_back(i + 1);
210  }
211  else if (_q_function_type == TOPOLOGY)
212  {
213  if (isParamValid("ring_first") && isParamValid("ring_last"))
214  {
215  _ring_first = getParam<unsigned int>("ring_first");
216  _ring_last = getParam<unsigned int>("ring_last");
217  }
218  else
219  mooseError(
220  "DomainIntegral error: must set ring_first and ring_last if q_function_type = Topology.");
221  for (unsigned int i = _ring_first; i <= _ring_last; ++i)
222  _ring_vec.push_back(i);
223  }
224  else
225  paramError("q_function_type", "DomainIntegral error: invalid q_function_type.");
226 
227  if (isParamValid("crack_front_points"))
228  _crack_front_points = getParam<std::vector<Point>>("crack_front_points");
229 
230  if (isParamValid("crack_front_points_provider"))
231  {
232  if (!isParamValid("number_points_from_provider"))
233  paramError("number_points_from_provider",
234  "DomainIntegral error: when crack_front_points_provider is used, "
235  "number_points_from_provider must be provided.");
237  _crack_front_points_provider = getParam<UserObjectName>("crack_front_points_provider");
238  }
239  else if (isParamValid("number_points_from_provider"))
240  paramError("crack_front_points_provider",
241  "DomainIntegral error: number_points_from_provider is provided but "
242  "crack_front_points_provider cannot be found.");
243  if (isParamValid("crack_mouth_boundary"))
244  _crack_mouth_boundary_names = getParam<std::vector<BoundaryName>>("crack_mouth_boundary");
245  if (isParamValid("intersecting_boundary"))
246  _intersecting_boundary_names = getParam<std::vector<BoundaryName>>("intersecting_boundary");
247  if (_radius_inner.size() != _radius_outer.size())
248  mooseError("Number of entries in 'radius_inner' and 'radius_outer' must match.");
249 
250  bool youngs_modulus_set(false);
251  bool poissons_ratio_set(false);
252 
253  // All domain integral types block restrict the objects created by this action.
254  _blocks = getParam<std::vector<SubdomainName>>("block");
255 
256  MultiMooseEnum integral_moose_enums = getParam<MultiMooseEnum>("integrals");
257  for (unsigned int i = 0; i < integral_moose_enums.size(); ++i)
258  {
259  _displacements = getParam<std::vector<VariableName>>("displacements");
260 
261  if (_displacements.size() < 2)
262  paramError(
263  "displacements",
264  "DomainIntegral error: The size of the displacements vector should at least be 2.");
265 
266  if (integral_moose_enums[i] != "JIntegral" && integral_moose_enums[i] != "CIntegral" &&
267  integral_moose_enums[i] != "KFromJIntegral")
268  {
269  // Check that parameters required for interaction integrals are defined
270  if (!(isParamValid("poissons_ratio")) || !(isParamValid("youngs_modulus")))
271  mooseError(
272  "DomainIntegral error: must set Poisson's ratio and Young's modulus for integral: ",
273  integral_moose_enums[i]);
274 
275  _poissons_ratio = getParam<Real>("poissons_ratio");
276  poissons_ratio_set = true;
277  _youngs_modulus = getParam<Real>("youngs_modulus");
278  youngs_modulus_set = true;
279  }
280 
281  _integrals.insert(INTEGRAL(int(integral_moose_enums.get(i))));
282  }
283 
284  if ((_integrals.count(J_INTEGRAL) != 0 && _integrals.count(C_INTEGRAL) != 0) ||
285  (_integrals.count(J_INTEGRAL) != 0 && _integrals.count(K_FROM_J_INTEGRAL) != 0) ||
286  (_integrals.count(C_INTEGRAL) != 0 && _integrals.count(K_FROM_J_INTEGRAL) != 0))
287  paramError("integrals",
288  "JIntegral, CIntegral, and KFromJIntegral options are mutually exclusive");
289 
290  // Acommodate deprecated parameter convert_J_to_K
291  if (_convert_J_to_K && _integrals.count(K_FROM_J_INTEGRAL) != 0)
292  {
294  _integrals.erase(J_INTEGRAL);
295  }
296 
297  if (isParamValid("temperature"))
298  _temp = getParam<VariableName>("temperature");
299 
300  if (_temp != "" && !isParamValid("eigenstrain_names"))
301  paramError(
302  "eigenstrain_names",
303  "DomainIntegral error: must provide `eigenstrain_names` when temperature is coupled.");
304 
306  _integrals.count(INTERACTION_INTEGRAL_KII) == 0 ||
308  paramError("integrals",
309  "DomainIntegral error: must calculate KI, KII and KIII to get equivalent K.");
310 
311  if (isParamValid("output_variable"))
312  {
313  _output_variables = getParam<std::vector<VariableName>>("output_variable");
314  if (_crack_front_points.size() > 0)
315  paramError("output_variables",
316  "'output_variables' not yet supported with 'crack_front_points'");
317  }
318 
319  if (_integrals.count(K_FROM_J_INTEGRAL) != 0)
320  {
321  if (!isParamValid("youngs_modulus") || !isParamValid("poissons_ratio"))
322  mooseError("DomainIntegral error: must set Young's modulus and Poisson's ratio "
323  "if K_FROM_J_INTEGRAL is selected.");
324  if (!youngs_modulus_set)
325  _youngs_modulus = getParam<Real>("youngs_modulus");
326  if (!poissons_ratio_set)
327  _poissons_ratio = getParam<Real>("poissons_ratio");
328  }
329 
330  if (_integrals.count(J_INTEGRAL) != 0 || _integrals.count(C_INTEGRAL) != 0 ||
331  _integrals.count(K_FROM_J_INTEGRAL) != 0)
332  {
333  if (isParamValid("eigenstrain_gradient"))
334  paramError("eigenstrain_gradient",
335  "'eigenstrain_gradient' cannot be specified when the computed integrals include "
336  "JIntegral, CIntegral, or KFromJIntegral");
337  if (isParamValid("body_force"))
338  paramError("body_force",
339  "'body_force' cannot be specified when the computed integrals include JIntegral, "
340  "CIntegral, or KFromJIntegral");
341  }
342  if (isParamValid("eigenstrain_gradient") && (_temp != "" || isParamValid("eigenstrain_names")))
343  paramError("eigenstrain_gradient",
344  "'eigenstrain_gradient' cannot be specified together with 'temperature' or "
345  "'eigenstrain_names'. These are for separate, mutually exclusive systems for "
346  "including the effect of eigenstrains");
347 }
const RealVectorValue _crack_direction_vector
Vector optionally used to prescribe direction of crack extension.
unsigned int _symmetry_plane
Identifier for which plane is the symmetry plane.
VariableName _temp
Temperature variable.
std::vector< VariableName > _output_variables
List of variables for which values are to be sampled and output at the crack front points...
const std::vector< BoundaryName > & _boundary_names
Boundaries containing the crack front points.
UserObjectName _crack_front_points_provider
Name of crack front points provider user object used to optionally define the crack points...
const std::string _order
Order and family of the AuxVariables optionally created to output the values of q.
const bool _have_crack_direction_vector_end_2
Whether the crack direction vector at the 2nd end of the crack has been provided. ...
MaterialPropertyName _functionally_graded_youngs_modulus
Material property name for spatially-dependent Youngs modulus for functionally graded materials...
void paramError(const std::string &param, Args... args) const
const T & getParam(const std::string &name) const
Action(const InputParameters &parameters)
const bool _have_crack_direction_vector
Whether the crack direction vector has been provided.
std::vector< BoundaryName > _intersecting_boundary_names
Names of boundaries optionally used for improved computation of crack extension direction at ends of ...
const bool _use_ad
Whether to create automatic differentiation objects from the action.
std::vector< Real > _radius_inner
Sets of inner and outer radii of the rings used for the domain form of the fracture integrals...
bool _has_symmetry_plane
Whether the model has a symmetry plane passing through the plane of the crack.
bool _use_crack_front_points_provider
Whether to use a crack front points provider.
bool _treat_as_2d
Whether fracture computations for a 3D model should be treated as though it were a 2D model...
Real _poissons_ratio
Poisson&#39;s ratio of material.
const std::string _family
MooseEnum _direction_method_moose_enum
Enum used to define the method to compute crack front direction.
std::vector< VariableName > _displacements
Vector of displacement variables.
std::vector< BoundaryName > _crack_mouth_boundary_names
Names of boundaries optionally used to define the crack mouth location.
std::vector< SubdomainName > _blocks
Blocks for which the domain integrals are to be computed.
MooseEnum _end_direction_method_moose_enum
Enum used to define the method to compute crack front direction at ends of crack front.
std::vector< unsigned int > _ring_vec
Vector of ids for the individual rings on which the fracture integral is computed.
bool _incremental
Whether the constitutive models for the mechanics calculations use an incremental form...
std::vector< Real > _radius_outer
MooseEnum _q_function_type
How the q function is evaluated (geometric distance from crack front or ring of elements) ...
const bool _have_crack_direction_vector_end_1
Whether the crack direction vector at the 1st end of the crack has been provided. ...
bool _convert_J_to_K
Whether to convert the J-integral to a stress intensity factor (K) –deprecated.
MaterialPropertyName _functionally_graded_youngs_modulus_crack_dir_gradient
Material property name for the Youngs modulus derivative for functionally graded materials.
bool _use_displaced_mesh
Whether to compute the fracture integrals on the displaced mesh.
const RealVectorValue _crack_direction_vector_end_1
Vector optionally used to prescribe direction of crack extension at the 1st end of the crack...
INTEGRAL
Enum used to select the type of integral to be performed.
void mooseError(Args &&... args) const
bool _output_q
Whether to ouput the q function as a set of AuxVariables.
bool _get_equivalent_k
Whether to compute the equivalent K from the individual fracture integrals for mixed-mode fracture...
MooseEnum _position_type
How the distance along the crack front is measured (angle or distance)
bool isParamValid(const std::string &name) const
std::vector< Point > _crack_front_points
User-defined vector of crack front points.
unsigned int _ring_last
Number of elements away from the crack tip to outside of outer ring with the topological q function...
unsigned int _ring_first
Number of elements away from the crack tip to inside of inner ring with the topological q function...
bool _fgm_crack
Whether the crack lives in a functionally-graded material.
const RealVectorValue _crack_direction_vector_end_2
Vector optionally used to prescribe direction of crack extension at the 2nd end of the crack...
unsigned int _axis_2d
Out-of-plane axis for 3D models treated as 2D.
bool _closed_loop
Indicates whether the crack forms a closed loop.
Real _youngs_modulus
Young&#39;s modulus of material.
std::set< INTEGRAL > _integrals
Container for enumerations describing the individual integrals computed.

◆ ~DomainIntegralAction()

DomainIntegralAction::~DomainIntegralAction ( )

Definition at line 349 of file DomainIntegralAction.C.

349 {}

Member Function Documentation

◆ act()

void DomainIntegralAction::act ( )
overridevirtual

Implements Action.

Definition at line 352 of file DomainIntegralAction.C.

353 {
354  const std::string uo_name("crackFrontDefinition");
355  const std::string ak_base_name("q");
356  const std::string av_base_name("q");
357  const unsigned int num_crack_front_points = calcNumCrackFrontPoints();
358  const std::string aux_stress_base_name("aux_stress");
359  const std::string aux_grad_disp_base_name("aux_grad_disp");
360 
361  // checking if built with xfem and setting flags for vpps used by xfem
362  std::vector<std::string> xfem_exec_flags = {EXEC_XFEM_MARK, EXEC_TIMESTEP_END};
363 
364  std::string ad_prepend = "";
365  if (_use_ad)
366  ad_prepend = "AD";
367 
368  if (_current_task == "add_user_object")
369  {
370  const std::string uo_type_name("CrackFrontDefinition");
371 
372  InputParameters params = _factory.getValidParams(uo_type_name);
374  {
375  // The CrackFrontDefinition updates the vpps and MUST execute before them
376  params.set<int>("execution_order_group") = -1;
377  params.set<ExecFlagEnum>("execute_on") = xfem_exec_flags;
378  }
379  else
380  params.set<ExecFlagEnum>("execute_on") = {EXEC_INITIAL, EXEC_TIMESTEP_END};
381 
382  params.set<MooseEnum>("crack_direction_method") = _direction_method_moose_enum;
383  params.set<MooseEnum>("crack_end_direction_method") = _end_direction_method_moose_enum;
385  params.set<RealVectorValue>("crack_direction_vector") = _crack_direction_vector;
387  params.set<RealVectorValue>("crack_direction_vector_end_1") = _crack_direction_vector_end_1;
389  params.set<RealVectorValue>("crack_direction_vector_end_2") = _crack_direction_vector_end_2;
390  if (_crack_mouth_boundary_names.size() != 0)
391  params.set<std::vector<BoundaryName>>("crack_mouth_boundary") = _crack_mouth_boundary_names;
392  if (_intersecting_boundary_names.size() != 0)
393  params.set<std::vector<BoundaryName>>("intersecting_boundary") = _intersecting_boundary_names;
394  params.set<bool>("2d") = _treat_as_2d;
395  params.set<unsigned int>("axis_2d") = _axis_2d;
397  params.set<unsigned int>("symmetry_plane") = _symmetry_plane;
398  if (_boundary_names.size() != 0)
399  params.set<std::vector<BoundaryName>>("boundary") = _boundary_names;
400  if (_crack_front_points.size() != 0)
401  params.set<std::vector<Point>>("crack_front_points") = _crack_front_points;
403  params.applyParameters(parameters(),
404  {"crack_front_points_provider, number_points_from_provider"});
405  if (_closed_loop)
406  params.set<bool>("closed_loop") = _closed_loop;
407  params.set<bool>("use_displaced_mesh") = _use_displaced_mesh;
408  if (_integrals.count(INTERACTION_INTEGRAL_T) != 0)
409  {
410  params.set<VariableName>("disp_x") = _displacements[0];
411  params.set<VariableName>("disp_y") = _displacements[1];
412  if (_displacements.size() == 3)
413  params.set<VariableName>("disp_z") = _displacements[2];
414  params.set<bool>("t_stress") = true;
415  }
416 
417  unsigned int nrings = 0;
418  if (_q_function_type == TOPOLOGY)
419  {
420  params.set<bool>("q_function_rings") = true;
421  params.set<unsigned int>("last_ring") = _ring_last;
422  params.set<unsigned int>("first_ring") = _ring_first;
423  nrings = _ring_last - _ring_first + 1;
424  }
425  else if (_q_function_type == GEOMETRY)
426  {
427  params.set<std::vector<Real>>("j_integral_radius_inner") = _radius_inner;
428  params.set<std::vector<Real>>("j_integral_radius_outer") = _radius_outer;
429  nrings = _ring_vec.size();
430  }
431 
432  params.set<unsigned int>("nrings") = nrings;
433  params.set<MooseEnum>("q_function_type") = _q_function_type;
434 
435  _problem->addUserObject(uo_type_name, uo_name, params);
436  }
437  else if (_current_task == "add_aux_variable" && _output_q)
438  {
439  for (unsigned int ring_index = 0; ring_index < _ring_vec.size(); ++ring_index)
440  {
441  std::string aux_var_type;
442  if (_family == "LAGRANGE")
443  aux_var_type = "MooseVariable";
444  else if (_family == "MONOMIAL")
445  aux_var_type = "MooseVariableConstMonomial";
446  else if (_family == "SCALAR")
447  aux_var_type = "MooseVariableScalar";
448  else
449  mooseError("Unsupported finite element family in, " + name() +
450  ". Please use LAGRANGE, MONOMIAL, or SCALAR");
451 
452  auto params = _factory.getValidParams(aux_var_type);
453  params.set<MooseEnum>("order") = _order;
454  params.set<MooseEnum>("family") = _family;
455 
457  {
458  std::ostringstream av_name_stream;
459  av_name_stream << av_base_name << "_" << _ring_vec[ring_index];
460  _problem->addAuxVariable(aux_var_type, av_name_stream.str(), params);
461  }
462  else
463  {
464  for (unsigned int cfp_index = 0; cfp_index < num_crack_front_points; ++cfp_index)
465  {
466  std::ostringstream av_name_stream;
467  av_name_stream << av_base_name << "_" << cfp_index + 1 << "_" << _ring_vec[ring_index];
468  _problem->addAuxVariable(aux_var_type, av_name_stream.str(), params);
469  }
470  }
471  }
472  }
473 
474  else if (_current_task == "add_aux_kernel" && _output_q)
475  {
476  std::string ak_type_name;
477  unsigned int nrings = 0;
478  if (_q_function_type == GEOMETRY)
479  {
480  ak_type_name = "DomainIntegralQFunction";
481  nrings = _ring_vec.size();
482  }
483  else if (_q_function_type == TOPOLOGY)
484  {
485  ak_type_name = "DomainIntegralTopologicalQFunction";
486  nrings = _ring_last - _ring_first + 1;
487  }
488 
489  InputParameters params = _factory.getValidParams(ak_type_name);
490  params.set<ExecFlagEnum>("execute_on") = {EXEC_INITIAL, EXEC_TIMESTEP_END};
491  params.set<UserObjectName>("crack_front_definition") = uo_name;
492  params.set<bool>("use_displaced_mesh") = _use_displaced_mesh;
493 
494  for (unsigned int ring_index = 0; ring_index < nrings; ++ring_index)
495  {
496  if (_q_function_type == GEOMETRY)
497  {
498  params.set<Real>("j_integral_radius_inner") = _radius_inner[ring_index];
499  params.set<Real>("j_integral_radius_outer") = _radius_outer[ring_index];
500  }
501  else if (_q_function_type == TOPOLOGY)
502  {
503  params.set<unsigned int>("ring_index") = _ring_first + ring_index;
504  }
505 
507  {
508  std::ostringstream ak_name_stream;
509  ak_name_stream << ak_base_name << "_" << _ring_vec[ring_index];
510  std::ostringstream av_name_stream;
511  av_name_stream << av_base_name << "_" << _ring_vec[ring_index];
512  params.set<AuxVariableName>("variable") = av_name_stream.str();
513  _problem->addAuxKernel(ak_type_name, ak_name_stream.str(), params);
514  }
515  else
516  {
517  for (unsigned int cfp_index = 0; cfp_index < num_crack_front_points; ++cfp_index)
518  {
519  std::ostringstream ak_name_stream;
520  ak_name_stream << ak_base_name << "_" << cfp_index + 1 << "_" << _ring_vec[ring_index];
521  std::ostringstream av_name_stream;
522  av_name_stream << av_base_name << "_" << cfp_index + 1 << "_" << _ring_vec[ring_index];
523  params.set<AuxVariableName>("variable") = av_name_stream.str();
524  params.set<unsigned int>("crack_front_point_index") = cfp_index;
525  _problem->addAuxKernel(ak_type_name, ak_name_stream.str(), params);
526  }
527  }
528  }
529  }
530 
531  else if (_current_task == "add_postprocessor")
532  {
533  for (std::set<INTEGRAL>::iterator sit = _integrals.begin(); sit != _integrals.end(); ++sit)
534  {
535  std::string pp_base_name;
536  switch (*sit)
537  {
538  case J_INTEGRAL:
539  pp_base_name = "J";
540  break;
541 
542  case C_INTEGRAL:
543  pp_base_name = "C";
544  break;
545 
546  case K_FROM_J_INTEGRAL:
547  pp_base_name = "K";
548  break;
549 
551  pp_base_name = "II_KI";
552  break;
553 
555  pp_base_name = "II_KII";
556  break;
557 
559  pp_base_name = "II_KIII";
560  break;
561 
563  pp_base_name = "II_T";
564  break;
565  }
566  const std::string pp_type_name("VectorPostprocessorComponent");
567  InputParameters params = _factory.getValidParams(pp_type_name);
568  for (unsigned int ring_index = 0; ring_index < _ring_vec.size(); ++ring_index)
569  {
571  {
572  params.set<VectorPostprocessorName>("vectorpostprocessor") =
573  pp_base_name + "_2DVPP_" + Moose::stringify(_ring_vec[ring_index]);
574  std::string pp_name = pp_base_name + +"_" + Moose::stringify(_ring_vec[ring_index]);
575  params.set<unsigned int>("index") = 0;
576  params.set<std::string>("vector_name") =
577  pp_base_name + "_" + Moose::stringify(_ring_vec[ring_index]);
578  _problem->addPostprocessor(pp_type_name, pp_name, params);
579  }
580  else
581  {
582  for (unsigned int cfp_index = 0; cfp_index < num_crack_front_points; ++cfp_index)
583  {
584  params.set<VectorPostprocessorName>("vectorpostprocessor") =
585  pp_base_name + "_" + Moose::stringify(_ring_vec[ring_index]);
586  std::string pp_name = pp_base_name + "_" + Moose::stringify(cfp_index + 1) + "_" +
587  Moose::stringify(_ring_vec[ring_index]);
588  params.set<unsigned int>("index") = cfp_index;
589  params.set<std::string>("vector_name") =
590  pp_base_name + "_" + Moose::stringify(_ring_vec[ring_index]);
591  _problem->addPostprocessor(pp_type_name, pp_name, params);
592  }
593  }
594  }
595  }
596 
597  if (_get_equivalent_k)
598  {
599  std::string pp_base_name("Keq");
600  const std::string pp_type_name("VectorPostprocessorComponent");
601  InputParameters params = _factory.getValidParams(pp_type_name);
602  for (unsigned int ring_index = 0; ring_index < _ring_vec.size(); ++ring_index)
603  {
605  {
606  params.set<VectorPostprocessorName>("vectorpostprocessor") =
607  pp_base_name + "_2DVPP_" + Moose::stringify(_ring_vec[ring_index]);
608  std::string pp_name = pp_base_name + +"_" + Moose::stringify(_ring_vec[ring_index]);
609  params.set<unsigned int>("index") = 0;
610  params.set<std::string>("vector_name") =
611  pp_base_name + "_" + Moose::stringify(_ring_vec[ring_index]);
612  _problem->addPostprocessor(pp_type_name, pp_name, params);
613  }
614  else
615  {
616  for (unsigned int cfp_index = 0; cfp_index < num_crack_front_points; ++cfp_index)
617  {
618  params.set<VectorPostprocessorName>("vectorpostprocessor") =
619  pp_base_name + "_" + Moose::stringify(_ring_vec[ring_index]);
620  std::string pp_name = pp_base_name + "_" + Moose::stringify(cfp_index + 1) + "_" +
621  Moose::stringify(_ring_vec[ring_index]);
622  params.set<unsigned int>("index") = cfp_index;
623  params.set<std::string>("vector_name") =
624  pp_base_name + "_" + Moose::stringify(_ring_vec[ring_index]);
625  _problem->addPostprocessor(pp_type_name, pp_name, params);
626  }
627  }
628  }
629  }
630 
631  for (unsigned int i = 0; i < _output_variables.size(); ++i)
632  {
633  const std::string ov_base_name(_output_variables[i]);
634  const std::string pp_type_name("CrackFrontData");
635  InputParameters params = _factory.getValidParams(pp_type_name);
636  params.set<ExecFlagEnum>("execute_on") = EXEC_TIMESTEP_END;
637  params.set<UserObjectName>("crack_front_definition") = uo_name;
639  {
640  std::ostringstream pp_name_stream;
641  pp_name_stream << ov_base_name << "_crack";
642  params.set<VariableName>("variable") = _output_variables[i];
643  _problem->addPostprocessor(pp_type_name, pp_name_stream.str(), params);
644  }
645  else
646  {
647  for (unsigned int cfp_index = 0; cfp_index < num_crack_front_points; ++cfp_index)
648  {
649  std::ostringstream pp_name_stream;
650  pp_name_stream << ov_base_name << "_crack_" << cfp_index + 1;
651  params.set<VariableName>("variable") = _output_variables[i];
652  params.set<unsigned int>("crack_front_point_index") = cfp_index;
653  _problem->addPostprocessor(pp_type_name, pp_name_stream.str(), params);
654  }
655  }
656  }
657  }
658 
659  else if (_current_task == "add_vector_postprocessor")
660  {
661  if (_integrals.count(J_INTEGRAL) != 0 || _integrals.count(C_INTEGRAL) != 0 ||
662  _integrals.count(K_FROM_J_INTEGRAL) != 0)
663  {
664  std::string vpp_base_name;
665  std::string jintegral_selection = "JIntegral";
666 
667  if (_integrals.count(J_INTEGRAL) != 0)
668  {
669  vpp_base_name = "J";
670  jintegral_selection = "JIntegral";
671  }
672  else if (_integrals.count(K_FROM_J_INTEGRAL) != 0)
673  {
674  vpp_base_name = "K";
675  jintegral_selection = "KFromJIntegral";
676  }
677  else if (_integrals.count(C_INTEGRAL) != 0)
678  {
679  vpp_base_name = "C";
680  jintegral_selection = "CIntegral";
681  }
682 
684  vpp_base_name += "_2DVPP";
685 
686  const std::string vpp_type_name("JIntegral");
687  InputParameters params = _factory.getValidParams(vpp_type_name);
688  if (!getParam<bool>("output_vpp"))
689  params.set<std::vector<OutputName>>("outputs") = {"none"};
690 
691  params.set<ExecFlagEnum>("execute_on") = EXEC_TIMESTEP_END;
692  params.set<UserObjectName>("crack_front_definition") = uo_name;
693  params.set<std::vector<SubdomainName>>("block") = {_blocks};
694  params.set<MooseEnum>("position_type") = _position_type;
695 
696  if (_integrals.count(K_FROM_J_INTEGRAL) != 0)
697  {
698  params.set<Real>("youngs_modulus") = _youngs_modulus;
699  params.set<Real>("poissons_ratio") = _poissons_ratio;
700  }
701 
703  params.set<unsigned int>("symmetry_plane") = _symmetry_plane;
704 
705  // Select the integral type to be computed in JIntegral vector postprocessor
706  params.set<MooseEnum>("integral") = jintegral_selection;
707 
708  params.set<std::vector<VariableName>>("displacements") = _displacements;
709  params.set<bool>("use_displaced_mesh") = _use_displaced_mesh;
710  for (unsigned int ring_index = 0; ring_index < _ring_vec.size(); ++ring_index)
711  {
712  params.set<unsigned int>("ring_index") = _ring_vec[ring_index];
713  params.set<MooseEnum>("q_function_type") = _q_function_type;
714 
715  std::string vpp_name = vpp_base_name + "_" + Moose::stringify(_ring_vec[ring_index]);
716  _problem->addVectorPostprocessor(vpp_type_name, vpp_name, params);
717  }
718  }
719 
720  if (_integrals.count(INTERACTION_INTEGRAL_KI) != 0 ||
721  _integrals.count(INTERACTION_INTEGRAL_KII) != 0 ||
724  {
727  paramError("symmetry_plane",
728  "In DomainIntegral, symmetry_plane option cannot be used with mode-II or "
729  "mode-III interaction integral");
730 
731  std::string vpp_base_name;
732  std::string vpp_type_name(ad_prepend + "InteractionIntegral");
733 
734  InputParameters params = _factory.getValidParams(vpp_type_name);
735  if (!getParam<bool>("output_vpp"))
736  params.set<std::vector<OutputName>>("outputs") = {"none"};
737 
739  params.set<ExecFlagEnum>("execute_on") = xfem_exec_flags;
740  else
741  params.set<ExecFlagEnum>("execute_on") = {EXEC_TIMESTEP_END};
742 
743  if (isParamValid("additional_eigenstrain_00") && isParamValid("additional_eigenstrain_01") &&
744  isParamValid("additional_eigenstrain_11") && isParamValid("additional_eigenstrain_22"))
745  {
746  params.set<CoupledName>("additional_eigenstrain_00") =
747  getParam<CoupledName>("additional_eigenstrain_00");
748  params.set<CoupledName>("additional_eigenstrain_01") =
749  getParam<CoupledName>("additional_eigenstrain_01");
750  params.set<CoupledName>("additional_eigenstrain_11") =
751  getParam<CoupledName>("additional_eigenstrain_11");
752  params.set<CoupledName>("additional_eigenstrain_22") =
753  getParam<CoupledName>("additional_eigenstrain_22");
754  }
755 
756  params.set<UserObjectName>("crack_front_definition") = uo_name;
757  params.set<bool>("use_displaced_mesh") = _use_displaced_mesh;
758  params.set<std::vector<SubdomainName>>("block") = {_blocks};
759 
761  params.set<unsigned int>("symmetry_plane") = _symmetry_plane;
762 
763  if (_fgm_crack)
764  {
765  params.set<MaterialPropertyName>(
766  "functionally_graded_youngs_modulus_crack_dir_gradient") = {
768  params.set<MaterialPropertyName>("functionally_graded_youngs_modulus") = {
770  }
771 
772  params.set<Real>("poissons_ratio") = _poissons_ratio;
773  params.set<Real>("youngs_modulus") = _youngs_modulus;
774  params.set<std::vector<VariableName>>("displacements") = _displacements;
775  if (_temp != "")
776  params.set<std::vector<VariableName>>("temperature") = {_temp};
777 
778  if (parameters().isParamValid("eigenstrain_gradient"))
779  params.set<MaterialPropertyName>("eigenstrain_gradient") =
780  parameters().get<MaterialPropertyName>("eigenstrain_gradient");
781  if (parameters().isParamValid("body_force"))
782  params.set<MaterialPropertyName>("body_force") =
783  parameters().get<MaterialPropertyName>("body_force");
784 
785  for (std::set<INTEGRAL>::iterator sit = _integrals.begin(); sit != _integrals.end(); ++sit)
786  {
787  switch (*sit)
788  {
789  case J_INTEGRAL:
790  continue;
791 
792  case C_INTEGRAL:
793  continue;
794 
795  case K_FROM_J_INTEGRAL:
796  continue;
797 
799  vpp_base_name = "II_KI";
800  params.set<Real>("K_factor") =
801  0.5 * _youngs_modulus / (1.0 - std::pow(_poissons_ratio, 2.0));
802  params.set<MooseEnum>("sif_mode") = "KI";
803  break;
804 
806  vpp_base_name = "II_KII";
807  params.set<Real>("K_factor") =
808  0.5 * _youngs_modulus / (1.0 - std::pow(_poissons_ratio, 2.0));
809  params.set<MooseEnum>("sif_mode") = "KII";
810  break;
811 
813  vpp_base_name = "II_KIII";
814  params.set<Real>("K_factor") = 0.5 * _youngs_modulus / (1.0 + _poissons_ratio);
815  params.set<MooseEnum>("sif_mode") = "KIII";
816  break;
817 
819  vpp_base_name = "II_T";
820  params.set<Real>("K_factor") = _youngs_modulus / (1 - std::pow(_poissons_ratio, 2));
821  params.set<MooseEnum>("sif_mode") = "T";
822  break;
823  }
825  vpp_base_name += "_2DVPP";
826  for (unsigned int ring_index = 0; ring_index < _ring_vec.size(); ++ring_index)
827  {
828  params.set<unsigned int>("ring_index") = _ring_vec[ring_index];
829  params.set<MooseEnum>("q_function_type") = _q_function_type;
830 
831  std::string vpp_name = vpp_base_name + "_" + Moose::stringify(_ring_vec[ring_index]);
832  _problem->addVectorPostprocessor(vpp_type_name, vpp_name, params);
833  }
834  }
835  }
836 
837  if (_get_equivalent_k)
838  {
839  std::string vpp_base_name("Keq");
841  vpp_base_name += "_2DVPP";
842  const std::string vpp_type_name("MixedModeEquivalentK");
843  InputParameters params = _factory.getValidParams(vpp_type_name);
844  params.set<ExecFlagEnum>("execute_on") = EXEC_TIMESTEP_END;
845  params.set<Real>("poissons_ratio") = _poissons_ratio;
846 
847  for (unsigned int ring_index = 0; ring_index < _ring_vec.size(); ++ring_index)
848  {
849  std::string ki_name = "II_KI_";
850  std::string kii_name = "II_KII_";
851  std::string kiii_name = "II_KIII_";
852  params.set<unsigned int>("ring_index") = _ring_vec[ring_index];
854  {
855  params.set<VectorPostprocessorName>("KI_vectorpostprocessor") =
856  ki_name + "2DVPP_" + Moose::stringify(_ring_vec[ring_index]);
857  params.set<VectorPostprocessorName>("KII_vectorpostprocessor") =
858  kii_name + "2DVPP_" + Moose::stringify(_ring_vec[ring_index]);
859  params.set<VectorPostprocessorName>("KIII_vectorpostprocessor") =
860  kiii_name + "2DVPP_" + Moose::stringify(_ring_vec[ring_index]);
861  }
862  else
863  {
864  params.set<VectorPostprocessorName>("KI_vectorpostprocessor") =
865  ki_name + Moose::stringify(_ring_vec[ring_index]);
866  params.set<VectorPostprocessorName>("KII_vectorpostprocessor") =
867  kii_name + Moose::stringify(_ring_vec[ring_index]);
868  params.set<VectorPostprocessorName>("KIII_vectorpostprocessor") =
869  kiii_name + Moose::stringify(_ring_vec[ring_index]);
870  }
871  params.set<std::string>("KI_vector_name") =
872  ki_name + Moose::stringify(_ring_vec[ring_index]);
873  params.set<std::string>("KII_vector_name") =
874  kii_name + Moose::stringify(_ring_vec[ring_index]);
875  params.set<std::string>("KIII_vector_name") =
876  kiii_name + Moose::stringify(_ring_vec[ring_index]);
877  std::string vpp_name = vpp_base_name + "_" + Moose::stringify(_ring_vec[ring_index]);
878  _problem->addVectorPostprocessor(vpp_type_name, vpp_name, params);
879  }
880  }
881 
883  {
884  for (unsigned int i = 0; i < _output_variables.size(); ++i)
885  {
886  const std::string vpp_type_name("VectorOfPostprocessors");
887  InputParameters params = _factory.getValidParams(vpp_type_name);
888  params.set<ExecFlagEnum>("execute_on") = EXEC_TIMESTEP_END;
889  std::ostringstream vpp_name_stream;
890  vpp_name_stream << _output_variables[i] << "_crack";
891  std::vector<PostprocessorName> postprocessor_names;
892  for (unsigned int cfp_index = 0; cfp_index < num_crack_front_points; ++cfp_index)
893  {
894  std::ostringstream pp_name_stream;
895  pp_name_stream << vpp_name_stream.str() << "_" << cfp_index + 1;
896  postprocessor_names.push_back(pp_name_stream.str());
897  }
898  params.set<std::vector<PostprocessorName>>("postprocessors") = postprocessor_names;
899  _problem->addVectorPostprocessor(vpp_type_name, vpp_name_stream.str(), params);
900  }
901  }
902  }
903 
904  else if (_current_task == "add_material")
905  {
906  if (_temp != "")
907  {
908  std::string mater_name;
909  const std::string mater_type_name("ThermalFractureIntegral");
910  mater_name = "ThermalFractureIntegral";
911 
912  InputParameters params = _factory.getValidParams(mater_type_name);
913  params.set<std::vector<MaterialPropertyName>>("eigenstrain_names") =
914  getParam<std::vector<MaterialPropertyName>>("eigenstrain_names");
915  params.set<std::vector<VariableName>>("temperature") = {_temp};
916  params.set<std::vector<SubdomainName>>("block") = {_blocks};
917  _problem->addMaterial(mater_type_name, mater_name, params);
918  }
919  MultiMooseEnum integral_moose_enums = getParam<MultiMooseEnum>("integrals");
920  bool have_j_integral = false;
921  bool have_c_integral = false;
922 
923  for (auto ime : integral_moose_enums)
924  {
925  if (ime == "JIntegral" || ime == "CIntegral" || ime == "KFromJIntegral" ||
926  ime == "InteractionIntegralKI" || ime == "InteractionIntegralKII" ||
927  ime == "InteractionIntegralKIII" || ime == "InteractionIntegralT")
928  have_j_integral = true;
929 
930  if (ime == "CIntegral")
931  have_c_integral = true;
932  }
933  if (have_j_integral)
934  {
935  std::string mater_name;
936  const std::string mater_type_name(ad_prepend + "StrainEnergyDensity");
937  mater_name = ad_prepend + "StrainEnergyDensity";
938 
939  InputParameters params = _factory.getValidParams(mater_type_name);
940  _incremental = getParam<bool>("incremental");
941  params.set<bool>("incremental") = _incremental;
942  params.set<std::vector<SubdomainName>>("block") = {_blocks};
943  _problem->addMaterial(mater_type_name, mater_name, params);
944 
945  {
946  std::string mater_name;
947  const std::string mater_type_name(ad_prepend + "EshelbyTensor");
948  mater_name = ad_prepend + "EshelbyTensor";
949 
950  InputParameters params = _factory.getValidParams(mater_type_name);
951  _displacements = getParam<std::vector<VariableName>>("displacements");
952  params.set<std::vector<VariableName>>("displacements") = _displacements;
953  params.set<std::vector<SubdomainName>>("block") = {_blocks};
954 
955  if (have_c_integral)
956  params.set<bool>("compute_dissipation") = true;
957 
958  if (_temp != "")
959  params.set<std::vector<VariableName>>("temperature") = {_temp};
960 
961  _problem->addMaterial(mater_type_name, mater_name, params);
962  }
963  // Strain energy rate density needed for C(t)/C* integral
964  if (have_c_integral)
965  {
966  std::string mater_name;
967  const std::string mater_type_name(ad_prepend + "StrainEnergyRateDensity");
968  mater_name = ad_prepend + "StrainEnergyRateDensity";
969 
970  InputParameters params = _factory.getValidParams(mater_type_name);
971  params.set<std::vector<SubdomainName>>("block") = {_blocks};
972  params.set<std::vector<MaterialName>>("inelastic_models") =
973  getParam<std::vector<MaterialName>>("inelastic_models");
974 
975  _problem->addMaterial(mater_type_name, mater_name, params);
976  }
977  }
978  }
979 }
const RealVectorValue _crack_direction_vector
Vector optionally used to prescribe direction of crack extension.
unsigned int _symmetry_plane
Identifier for which plane is the symmetry plane.
VariableName _temp
Temperature variable.
std::vector< VariableName > _output_variables
List of variables for which values are to be sampled and output at the crack front points...
const std::vector< BoundaryName > & _boundary_names
Boundaries containing the crack front points.
const std::string _order
Order and family of the AuxVariables optionally created to output the values of q.
const bool _have_crack_direction_vector_end_2
Whether the crack direction vector at the 2nd end of the crack has been provided. ...
MaterialPropertyName _functionally_graded_youngs_modulus
Material property name for spatially-dependent Youngs modulus for functionally graded materials...
void paramError(const std::string &param, Args... args) const
const T & getParam(const std::string &name) const
std::vector< std::pair< R1, R2 > > get(const std::string &param1, const std::string &param2) const
const ExecFlagType EXEC_XFEM_MARK
Exec flag used to execute MooseObjects while elements are being marked for cutting by XFEM...
Factory & _factory
const InputParameters & parameters() const
T & set(const std::string &name, bool quiet_mode=false)
const bool _have_crack_direction_vector
Whether the crack direction vector has been provided.
std::vector< BoundaryName > _intersecting_boundary_names
Names of boundaries optionally used for improved computation of crack extension direction at ends of ...
InputParameters getValidParams(const std::string &name) const
unsigned int calcNumCrackFrontPoints()
Compute the number of points on the crack front.
void applyParameters(const InputParameters &common, const std::vector< std::string > &exclude={}, const bool allow_private=false)
const ExecFlagType EXEC_TIMESTEP_END
const bool _use_ad
Whether to create automatic differentiation objects from the action.
std::vector< Real > _radius_inner
Sets of inner and outer radii of the rings used for the domain form of the fracture integrals...
bool _has_symmetry_plane
Whether the model has a symmetry plane passing through the plane of the crack.
bool _use_crack_front_points_provider
Whether to use a crack front points provider.
bool _treat_as_2d
Whether fracture computations for a 3D model should be treated as though it were a 2D model...
Real _poissons_ratio
Poisson&#39;s ratio of material.
const std::string & name() const
const std::string _family
MooseEnum _direction_method_moose_enum
Enum used to define the method to compute crack front direction.
std::vector< VariableName > _displacements
Vector of displacement variables.
std::vector< BoundaryName > _crack_mouth_boundary_names
Names of boundaries optionally used to define the crack mouth location.
const std::string & _current_task
std::vector< SubdomainName > _blocks
Blocks for which the domain integrals are to be computed.
std::string stringify(const T &t)
MooseEnum _end_direction_method_moose_enum
Enum used to define the method to compute crack front direction at ends of crack front.
std::vector< unsigned int > _ring_vec
Vector of ids for the individual rings on which the fracture integral is computed.
bool _incremental
Whether the constitutive models for the mechanics calculations use an incremental form...
std::vector< Real > _radius_outer
MooseEnum _q_function_type
How the q function is evaluated (geometric distance from crack front or ring of elements) ...
const bool _have_crack_direction_vector_end_1
Whether the crack direction vector at the 1st end of the crack has been provided. ...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
MaterialPropertyName _functionally_graded_youngs_modulus_crack_dir_gradient
Material property name for the Youngs modulus derivative for functionally graded materials.
bool _use_displaced_mesh
Whether to compute the fracture integrals on the displaced mesh.
std::vector< VariableName > CoupledName
const RealVectorValue _crack_direction_vector_end_1
Vector optionally used to prescribe direction of crack extension at the 1st end of the crack...
void mooseError(Args &&... args) const
bool _output_q
Whether to ouput the q function as a set of AuxVariables.
std::shared_ptr< FEProblemBase > & _problem
bool _get_equivalent_k
Whether to compute the equivalent K from the individual fracture integrals for mixed-mode fracture...
MooseEnum _position_type
How the distance along the crack front is measured (angle or distance)
bool isParamValid(const std::string &name) const
std::vector< Point > _crack_front_points
User-defined vector of crack front points.
unsigned int _ring_last
Number of elements away from the crack tip to outside of outer ring with the topological q function...
unsigned int _ring_first
Number of elements away from the crack tip to inside of inner ring with the topological q function...
bool _fgm_crack
Whether the crack lives in a functionally-graded material.
const RealVectorValue _crack_direction_vector_end_2
Vector optionally used to prescribe direction of crack extension at the 2nd end of the crack...
MooseUnits pow(const MooseUnits &, int)
unsigned int _axis_2d
Out-of-plane axis for 3D models treated as 2D.
bool _closed_loop
Indicates whether the crack forms a closed loop.
Real _youngs_modulus
Young&#39;s modulus of material.
std::set< INTEGRAL > _integrals
Container for enumerations describing the individual integrals computed.
bool isParamValid(const std::string &name) const
const ExecFlagType EXEC_INITIAL

◆ addRelationshipManagers() [1/3]

virtual void Action::addRelationshipManagers

◆ addRelationshipManagers() [2/3]

bool Action::addRelationshipManagers

◆ addRelationshipManagers() [3/3]

void DomainIntegralAction::addRelationshipManagers ( Moose::RelationshipManagerType  input_rm_type)
overridevirtual

Reimplemented from Action.

Definition at line 1017 of file DomainIntegralAction.C.

1018 {
1019  if (_integrals.count(INTERACTION_INTEGRAL_T) != 0)
1020  {
1021  InputParameters params = _factory.getValidParams("CrackFrontDefinition");
1022  addRelationshipManagers(input_rm_type, params);
1023  }
1024 }
Factory & _factory
InputParameters getValidParams(const std::string &name) const
virtual void addRelationshipManagers(Moose::RelationshipManagerType input_rm_type) override
std::set< INTEGRAL > _integrals
Container for enumerations describing the individual integrals computed.

◆ calcNumCrackFrontPoints()

unsigned int DomainIntegralAction::calcNumCrackFrontPoints ( )
protected

Compute the number of points on the crack front.

This is either the number of points in the crack front nodeset, or the number of points from the crack front points provider.

Returns
Number of points on the crack front

Definition at line 982 of file DomainIntegralAction.C.

Referenced by act().

983 {
984  unsigned int num_points = 0;
985  if (_boundary_names.size() != 0)
986  {
987  std::vector<BoundaryID> bids = _mesh->getBoundaryIDs(_boundary_names, true);
988  std::set<unsigned int> nodes;
989 
990  ConstBndNodeRange & bnd_nodes = *_mesh->getBoundaryNodeRange();
991  for (ConstBndNodeRange::const_iterator nd = bnd_nodes.begin(); nd != bnd_nodes.end(); ++nd)
992  {
993  const BndNode * bnode = *nd;
994  BoundaryID boundary_id = bnode->_bnd_id;
995 
996  for (unsigned int ibid = 0; ibid < bids.size(); ++ibid)
997  {
998  if (boundary_id == bids[ibid])
999  {
1000  nodes.insert(bnode->_node->id());
1001  break;
1002  }
1003  }
1004  }
1005  num_points = nodes.size();
1006  }
1007  else if (_crack_front_points.size() != 0)
1008  num_points = _crack_front_points.size();
1010  num_points = getParam<unsigned int>("number_points_from_provider");
1011  else
1012  mooseError("Must define either 'boundary' or 'crack_front_points'");
1013  return num_points;
1014 }
const std::vector< BoundaryName > & _boundary_names
Boundaries containing the crack front points.
BoundaryID _bnd_id
bool _use_crack_front_points_provider
Whether to use a crack front points provider.
dof_id_type id() const
boundary_id_type BoundaryID
libMesh::Node * _node
const_iterator end() const
std::shared_ptr< MooseMesh > & _mesh
vec_type::const_iterator const_iterator
const_iterator begin() const
void mooseError(Args &&... args) const
std::vector< Point > _crack_front_points
User-defined vector of crack front points.

◆ validParams()

InputParameters DomainIntegralAction::validParams ( )
static

Definition at line 34 of file DomainIntegralAction.C.

35 {
38  MultiMooseEnum integral_vec(
39  "JIntegral CIntegral KFromJIntegral InteractionIntegralKI InteractionIntegralKII "
40  "InteractionIntegralKIII InteractionIntegralT");
41  params.addClassDescription(
42  "Creates the MOOSE objects needed to compute fraction domain integrals");
43  params.addRequiredParam<MultiMooseEnum>("integrals",
44  integral_vec,
45  "Domain integrals to calculate. Choices are: " +
46  integral_vec.getRawNames());
47  params.addParam<std::vector<BoundaryName>>(
48  "boundary", {}, "Boundary containing the crack front points");
49  params.addParam<std::vector<Point>>("crack_front_points", "Set of points to define crack front");
50  params.addParam<std::string>(
51  "order", "FIRST", "Specifies the order of the FE shape function to use for q AuxVariables");
52  params.addParam<std::string>(
53  "family", "LAGRANGE", "Specifies the family of FE shape functions to use for q AuxVariables");
54  params.addParam<std::vector<Real>>("radius_inner", "Inner radius for volume integral domain");
55  params.addParam<std::vector<Real>>("radius_outer", "Outer radius for volume integral domain");
56  params.addParam<unsigned int>("ring_first",
57  "The first ring of elements for volume integral domain");
58  params.addParam<unsigned int>("ring_last",
59  "The last ring of elements for volume integral domain");
60  params.addParam<std::vector<VariableName>>(
61  "output_variable", "Variable values to be reported along the crack front");
62  params.addParam<Real>("poissons_ratio", "Poisson's ratio");
63  params.addParam<Real>("youngs_modulus", "Young's modulus");
64  params.addParam<std::vector<SubdomainName>>(
65  "block", {}, "The block ids where integrals are defined");
66  params.addParam<std::vector<VariableName>>(
67  "displacements",
68  {},
69  "The displacements appropriate for the simulation geometry and coordinate system");
70  params.addParam<VariableName>("temperature", "", "The temperature");
71  params.addParam<MaterialPropertyName>(
72  "functionally_graded_youngs_modulus_crack_dir_gradient",
73  "Gradient of the spatially varying Young's modulus provided in "
74  "'functionally_graded_youngs_modulus' in the direction of crack extension.");
75  params.addParam<MaterialPropertyName>(
76  "functionally_graded_youngs_modulus",
77  "Spatially varying elasticity modulus variable. This input is required when "
78  "using the functionally graded material capability.");
79  params.addCoupledVar("additional_eigenstrain_00",
80  "Optional additional eigenstrain variable that will be accounted for in the "
81  "interaction integral (component 00 or XX).");
82  params.addCoupledVar("additional_eigenstrain_01",
83  "Optional additional eigenstrain variable that will be accounted for in the "
84  "interaction integral (component 01 or XY).");
85  params.addCoupledVar("additional_eigenstrain_11",
86  "Optional additional eigenstrain variable that will be accounted for in the "
87  "interaction integral (component 11 or YY).");
88  params.addCoupledVar("additional_eigenstrain_22",
89  "Optional additional eigenstrain variable that will be accounted for in the "
90  "interaction integral (component 22 or ZZ).");
91  params.addCoupledVar("additional_eigenstrain_02",
92  "Optional additional eigenstrain variable that will be accounted for in the "
93  "interaction integral (component 02 or XZ).");
94  params.addCoupledVar("additional_eigenstrain_12",
95  "Optional additional eigenstrain variable that will be accounted for in the "
96  "interaction integral (component 12 or XZ).");
97  params.addParamNamesToGroup(
98  "additional_eigenstrain_00 additional_eigenstrain_01 additional_eigenstrain_11 "
99  "additional_eigenstrain_22 additional_eigenstrain_02 additional_eigenstrain_12",
100  "Generic eigenstrains for the computation of the interaction integral.");
101  MooseEnum position_type("Angle Distance", "Distance");
102  params.addParam<MooseEnum>(
103  "position_type",
104  position_type,
105  "The method used to calculate position along crack front. Options are: " +
106  position_type.getRawNames());
107  MooseEnum q_function_type("Geometry Topology", "Geometry");
108  params.addParam<MooseEnum>("q_function_type",
109  q_function_type,
110  "The method used to define the integration domain. Options are: " +
111  q_function_type.getRawNames());
112  params.addParam<bool>(
113  "equivalent_k",
114  false,
115  "Calculate an equivalent K from KI, KII and KIII, assuming self-similar crack growth.");
116  params.addParam<bool>("output_q", true, "Output q");
117  params.addRequiredParam<bool>(
118  "incremental", "Flag to indicate whether an incremental or total model is being used.");
119  params.addParam<std::vector<MaterialPropertyName>>(
120  "eigenstrain_names", "List of eigenstrains applied in the strain calculation");
121  params.addDeprecatedParam<bool>("convert_J_to_K",
122  false,
123  "Convert J-integral to stress intensity factor K.",
124  "This input parameter is deprecated and will be removed soon. "
125  "Use 'integrals = KFromJIntegral' to request output of the "
126  "conversion from the J-integral to stress intensity factors");
127  params.addParam<std::vector<MaterialName>>(
128  "inelastic_models",
129  "The material objects to use to calculate the strain energy rate density.");
130  params.addParam<MaterialPropertyName>("eigenstrain_gradient",
131  "Material defining gradient of eigenstrain tensor");
132  params.addParam<MaterialPropertyName>("body_force", "Material defining body force");
133  params.addParam<bool>("use_automatic_differentiation",
134  false,
135  "Flag to use automatic differentiation (AD) objects when possible");
136  params.addParam<bool>("output_vpp",
137  true,
138  "Flag to control the vector postprocessor outputs. Select false to "
139  "suppress the redundant csv files for each time step and ring");
140  return params;
141 }
void addDeprecatedParam(const std::string &name, const T &value, const std::string &doc_string, const std::string &deprecation_message)
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
static void includeCrackFrontDefinitionParams(InputParameters &params)
used by Actions to add CrackFrontDefinitionParams
std::string getRawNames() const
void addRequiredParam(const std::string &name, const std::string &doc_string)
static InputParameters validParams()
void addCoupledVar(const std::string &name, const std::string &doc_string)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void addClassDescription(const std::string &doc_string)
void addParamNamesToGroup(const std::string &space_delim_names, const std::string group_name)

Member Data Documentation

◆ _axis_2d

unsigned int DomainIntegralAction::_axis_2d
protected

Out-of-plane axis for 3D models treated as 2D.

Definition at line 106 of file DomainIntegralAction.h.

Referenced by act().

◆ _blocks

std::vector<SubdomainName> DomainIntegralAction::_blocks
protected

Blocks for which the domain integrals are to be computed.

Definition at line 124 of file DomainIntegralAction.h.

Referenced by act(), and DomainIntegralAction().

◆ _boundary_names

const std::vector<BoundaryName>& DomainIntegralAction::_boundary_names
protected

Boundaries containing the crack front points.

Definition at line 68 of file DomainIntegralAction.h.

Referenced by act(), and calcNumCrackFrontPoints().

◆ _closed_loop

bool DomainIntegralAction::_closed_loop
protected

Indicates whether the crack forms a closed loop.

Definition at line 72 of file DomainIntegralAction.h.

Referenced by act().

◆ _convert_J_to_K

bool DomainIntegralAction::_convert_J_to_K
protected

Whether to convert the J-integral to a stress intensity factor (K) –deprecated.

Definition at line 148 of file DomainIntegralAction.h.

Referenced by DomainIntegralAction().

◆ _crack_direction_vector

const RealVectorValue DomainIntegralAction::_crack_direction_vector
protected

Vector optionally used to prescribe direction of crack extension.

Definition at line 89 of file DomainIntegralAction.h.

Referenced by act().

◆ _crack_direction_vector_end_1

const RealVectorValue DomainIntegralAction::_crack_direction_vector_end_1
protected

Vector optionally used to prescribe direction of crack extension at the 1st end of the crack.

Definition at line 93 of file DomainIntegralAction.h.

Referenced by act().

◆ _crack_direction_vector_end_2

const RealVectorValue DomainIntegralAction::_crack_direction_vector_end_2
protected

Vector optionally used to prescribe direction of crack extension at the 2nd end of the crack.

Definition at line 97 of file DomainIntegralAction.h.

Referenced by act().

◆ _crack_front_points

std::vector<Point> DomainIntegralAction::_crack_front_points
protected

User-defined vector of crack front points.

Definition at line 70 of file DomainIntegralAction.h.

Referenced by act(), calcNumCrackFrontPoints(), and DomainIntegralAction().

◆ _crack_front_points_provider

UserObjectName DomainIntegralAction::_crack_front_points_provider
protected

Name of crack front points provider user object used to optionally define the crack points.

Definition at line 74 of file DomainIntegralAction.h.

Referenced by DomainIntegralAction().

◆ _crack_mouth_boundary_names

std::vector<BoundaryName> DomainIntegralAction::_crack_mouth_boundary_names
protected

Names of boundaries optionally used to define the crack mouth location.

Definition at line 99 of file DomainIntegralAction.h.

Referenced by act(), and DomainIntegralAction().

◆ _direction_method_moose_enum

MooseEnum DomainIntegralAction::_direction_method_moose_enum
protected

Enum used to define the method to compute crack front direction.

Definition at line 83 of file DomainIntegralAction.h.

Referenced by act().

◆ _displacements

std::vector<VariableName> DomainIntegralAction::_displacements
protected

Vector of displacement variables.

Definition at line 126 of file DomainIntegralAction.h.

Referenced by act(), and DomainIntegralAction().

◆ _end_direction_method_moose_enum

MooseEnum DomainIntegralAction::_end_direction_method_moose_enum
protected

Enum used to define the method to compute crack front direction at ends of crack front.

Definition at line 85 of file DomainIntegralAction.h.

Referenced by act().

◆ _family

const std::string DomainIntegralAction::_family
protected

Definition at line 80 of file DomainIntegralAction.h.

Referenced by act().

◆ _fgm_crack

bool DomainIntegralAction::_fgm_crack
protected

Whether the crack lives in a functionally-graded material.

Definition at line 150 of file DomainIntegralAction.h.

Referenced by act(), and DomainIntegralAction().

◆ _functionally_graded_youngs_modulus

MaterialPropertyName DomainIntegralAction::_functionally_graded_youngs_modulus
protected

Material property name for spatially-dependent Youngs modulus for functionally graded materials.

Definition at line 154 of file DomainIntegralAction.h.

Referenced by act(), and DomainIntegralAction().

◆ _functionally_graded_youngs_modulus_crack_dir_gradient

MaterialPropertyName DomainIntegralAction::_functionally_graded_youngs_modulus_crack_dir_gradient
protected

Material property name for the Youngs modulus derivative for functionally graded materials.

Definition at line 152 of file DomainIntegralAction.h.

Referenced by act(), and DomainIntegralAction().

◆ _get_equivalent_k

bool DomainIntegralAction::_get_equivalent_k
protected

Whether to compute the equivalent K from the individual fracture integrals for mixed-mode fracture.

Definition at line 138 of file DomainIntegralAction.h.

Referenced by act(), and DomainIntegralAction().

◆ _has_symmetry_plane

bool DomainIntegralAction::_has_symmetry_plane
protected

Whether the model has a symmetry plane passing through the plane of the crack.

Definition at line 130 of file DomainIntegralAction.h.

Referenced by act().

◆ _have_crack_direction_vector

const bool DomainIntegralAction::_have_crack_direction_vector
protected

Whether the crack direction vector has been provided.

Definition at line 87 of file DomainIntegralAction.h.

Referenced by act().

◆ _have_crack_direction_vector_end_1

const bool DomainIntegralAction::_have_crack_direction_vector_end_1
protected

Whether the crack direction vector at the 1st end of the crack has been provided.

Definition at line 91 of file DomainIntegralAction.h.

Referenced by act().

◆ _have_crack_direction_vector_end_2

const bool DomainIntegralAction::_have_crack_direction_vector_end_2
protected

Whether the crack direction vector at the 2nd end of the crack has been provided.

Definition at line 95 of file DomainIntegralAction.h.

Referenced by act().

◆ _incremental

bool DomainIntegralAction::_incremental
protected

Whether the constitutive models for the mechanics calculations use an incremental form.

Definition at line 146 of file DomainIntegralAction.h.

Referenced by act().

◆ _integrals

std::set<INTEGRAL> DomainIntegralAction::_integrals
protected

Container for enumerations describing the individual integrals computed.

Definition at line 66 of file DomainIntegralAction.h.

Referenced by act(), addRelationshipManagers(), and DomainIntegralAction().

◆ _intersecting_boundary_names

std::vector<BoundaryName> DomainIntegralAction::_intersecting_boundary_names
protected

Names of boundaries optionally used for improved computation of crack extension direction at ends of the crack where it intersects the prescribed boundaries.

Definition at line 102 of file DomainIntegralAction.h.

Referenced by act(), and DomainIntegralAction().

◆ _order

const std::string DomainIntegralAction::_order
protected

Order and family of the AuxVariables optionally created to output the values of q.

Definition at line 79 of file DomainIntegralAction.h.

Referenced by act().

◆ _output_q

bool DomainIntegralAction::_output_q
protected

Whether to ouput the q function as a set of AuxVariables.

Definition at line 142 of file DomainIntegralAction.h.

Referenced by act().

◆ _output_variables

std::vector<VariableName> DomainIntegralAction::_output_variables
protected

List of variables for which values are to be sampled and output at the crack front points.

Definition at line 118 of file DomainIntegralAction.h.

Referenced by act(), and DomainIntegralAction().

◆ _poissons_ratio

Real DomainIntegralAction::_poissons_ratio
protected

Poisson's ratio of material.

Definition at line 120 of file DomainIntegralAction.h.

Referenced by act(), and DomainIntegralAction().

◆ _position_type

MooseEnum DomainIntegralAction::_position_type
protected

How the distance along the crack front is measured (angle or distance)

Definition at line 134 of file DomainIntegralAction.h.

Referenced by act().

◆ _q_function_type

MooseEnum DomainIntegralAction::_q_function_type
protected

How the q function is evaluated (geometric distance from crack front or ring of elements)

Definition at line 136 of file DomainIntegralAction.h.

Referenced by act(), and DomainIntegralAction().

◆ _radius_inner

std::vector<Real> DomainIntegralAction::_radius_inner
protected

Sets of inner and outer radii of the rings used for the domain form of the fracture integrals.

These are defined in corresponding pairs for each ring.

Definition at line 110 of file DomainIntegralAction.h.

Referenced by act(), and DomainIntegralAction().

◆ _radius_outer

std::vector<Real> DomainIntegralAction::_radius_outer
protected

Definition at line 111 of file DomainIntegralAction.h.

Referenced by act(), and DomainIntegralAction().

◆ _ring_first

unsigned int DomainIntegralAction::_ring_first
protected

Number of elements away from the crack tip to inside of inner ring with the topological q function.

Definition at line 114 of file DomainIntegralAction.h.

Referenced by act(), and DomainIntegralAction().

◆ _ring_last

unsigned int DomainIntegralAction::_ring_last
protected

Number of elements away from the crack tip to outside of outer ring with the topological q function.

Definition at line 116 of file DomainIntegralAction.h.

Referenced by act(), and DomainIntegralAction().

◆ _ring_vec

std::vector<unsigned int> DomainIntegralAction::_ring_vec
protected

Vector of ids for the individual rings on which the fracture integral is computed.

Definition at line 144 of file DomainIntegralAction.h.

Referenced by act(), and DomainIntegralAction().

◆ _symmetry_plane

unsigned int DomainIntegralAction::_symmetry_plane
protected

Identifier for which plane is the symmetry plane.

Definition at line 132 of file DomainIntegralAction.h.

Referenced by act().

◆ _temp

VariableName DomainIntegralAction::_temp
protected

Temperature variable.

Definition at line 128 of file DomainIntegralAction.h.

Referenced by act(), and DomainIntegralAction().

◆ _treat_as_2d

bool DomainIntegralAction::_treat_as_2d
protected

Whether fracture computations for a 3D model should be treated as though it were a 2D model.

Definition at line 104 of file DomainIntegralAction.h.

Referenced by act().

◆ _use_ad

const bool DomainIntegralAction::_use_ad
protected

Whether to create automatic differentiation objects from the action.

Definition at line 156 of file DomainIntegralAction.h.

Referenced by act().

◆ _use_crack_front_points_provider

bool DomainIntegralAction::_use_crack_front_points_provider
protected

Whether to use a crack front points provider.

Definition at line 76 of file DomainIntegralAction.h.

Referenced by act(), calcNumCrackFrontPoints(), and DomainIntegralAction().

◆ _use_displaced_mesh

bool DomainIntegralAction::_use_displaced_mesh
protected

Whether to compute the fracture integrals on the displaced mesh.

Definition at line 140 of file DomainIntegralAction.h.

Referenced by act().

◆ _youngs_modulus

Real DomainIntegralAction::_youngs_modulus
protected

Young's modulus of material.

Definition at line 122 of file DomainIntegralAction.h.

Referenced by act(), and DomainIntegralAction().


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