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
 
virtual const std::string & name () const
 
std::string typeAndName () const
 
std::string errorPrefix (const std::string &error_type) const
 
void callMooseError (std::string msg, const bool with_prefix) const
 
MooseObjectParameterName uniqueParameterName (const std::string &parameter_name) const
 
const InputParametersparameters () const
 
MooseObjectName uniqueName () 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 &nm) 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
 
void connectControllableParams (const std::string &parameter, const std::string &object_type, const std::string &object_name, const std::string &object_parameter) const
 
void mooseError (Args &&... args) const
 
void mooseErrorNonPrefixed (Args &&... args) const
 
void mooseDocumentedError (const std::string &repo_name, const unsigned int issue_num, Args &&... args) const
 
void mooseWarning (Args &&... args) const
 
void mooseWarningNonPrefixed (Args &&... args) const
 
void mooseDeprecated (Args &&... args) const
 
void mooseInfo (Args &&... args) 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 ()
 

Public Attributes

const ConsoleStream _console
 

Static Public Attributes

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
 
const std::string _type
 
const std::string _name
 
const InputParameters_pars
 
Factory_factory
 
ActionFactory_action_factory
 
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...
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 isParamValid(const std::string &name) const
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.
const T & getParam(const std::string &name) const
void paramError(const std::string &param, Args... args) const
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)
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;
363  for (const auto & item :
365  {
366  if (item.name() == "XFEM_MARK")
367  {
368  // TODO the xfem_flags should be something like this:
369  // xfem_exec_flags = {item, EXEC_NONLINEAR, EXEC_TIMESTEP_END};
370  // but the item=XFEM_MARK flag causes xfem tests to diverge
371  xfem_exec_flags = {EXEC_INITIAL, EXEC_NONLINEAR, EXEC_TIMESTEP_END};
372  }
373  }
374 
375  std::string ad_prepend = "";
376  if (_use_ad)
377  ad_prepend = "AD";
378 
379  if (_current_task == "add_user_object")
380  {
381  const std::string uo_type_name("CrackFrontDefinition");
382 
383  InputParameters params = _factory.getValidParams(uo_type_name);
384  if (_use_crack_front_points_provider && !xfem_exec_flags.empty())
385  {
386  // The CrackFrontDefinition updates the vpps and MUST execute before them
387  params.set<int>("execution_order_group") = -1;
388  params.set<ExecFlagEnum>("execute_on") = xfem_exec_flags;
389  }
390  else
391  params.set<ExecFlagEnum>("execute_on") = {EXEC_INITIAL, EXEC_TIMESTEP_END};
392 
393  params.set<MooseEnum>("crack_direction_method") = _direction_method_moose_enum;
394  params.set<MooseEnum>("crack_end_direction_method") = _end_direction_method_moose_enum;
396  params.set<RealVectorValue>("crack_direction_vector") = _crack_direction_vector;
398  params.set<RealVectorValue>("crack_direction_vector_end_1") = _crack_direction_vector_end_1;
400  params.set<RealVectorValue>("crack_direction_vector_end_2") = _crack_direction_vector_end_2;
401  if (_crack_mouth_boundary_names.size() != 0)
402  params.set<std::vector<BoundaryName>>("crack_mouth_boundary") = _crack_mouth_boundary_names;
403  if (_intersecting_boundary_names.size() != 0)
404  params.set<std::vector<BoundaryName>>("intersecting_boundary") = _intersecting_boundary_names;
405  params.set<bool>("2d") = _treat_as_2d;
406  params.set<unsigned int>("axis_2d") = _axis_2d;
408  params.set<unsigned int>("symmetry_plane") = _symmetry_plane;
409  if (_boundary_names.size() != 0)
410  params.set<std::vector<BoundaryName>>("boundary") = _boundary_names;
411  if (_crack_front_points.size() != 0)
412  params.set<std::vector<Point>>("crack_front_points") = _crack_front_points;
414  params.applyParameters(parameters(),
415  {"crack_front_points_provider, number_points_from_provider"});
416  if (_closed_loop)
417  params.set<bool>("closed_loop") = _closed_loop;
418  params.set<bool>("use_displaced_mesh") = _use_displaced_mesh;
419  if (_integrals.count(INTERACTION_INTEGRAL_T) != 0)
420  {
421  params.set<VariableName>("disp_x") = _displacements[0];
422  params.set<VariableName>("disp_y") = _displacements[1];
423  if (_displacements.size() == 3)
424  params.set<VariableName>("disp_z") = _displacements[2];
425  params.set<bool>("t_stress") = true;
426  }
427 
428  unsigned int nrings = 0;
429  if (_q_function_type == TOPOLOGY)
430  {
431  params.set<bool>("q_function_rings") = true;
432  params.set<unsigned int>("last_ring") = _ring_last;
433  params.set<unsigned int>("first_ring") = _ring_first;
434  nrings = _ring_last - _ring_first + 1;
435  }
436  else if (_q_function_type == GEOMETRY)
437  {
438  params.set<std::vector<Real>>("j_integral_radius_inner") = _radius_inner;
439  params.set<std::vector<Real>>("j_integral_radius_outer") = _radius_outer;
440  nrings = _ring_vec.size();
441  }
442 
443  params.set<unsigned int>("nrings") = nrings;
444  params.set<MooseEnum>("q_function_type") = _q_function_type;
445 
446  _problem->addUserObject(uo_type_name, uo_name, params);
447  }
448  else if (_current_task == "add_aux_variable" && _output_q)
449  {
450  for (unsigned int ring_index = 0; ring_index < _ring_vec.size(); ++ring_index)
451  {
452  std::string aux_var_type;
453  if (_family == "LAGRANGE")
454  aux_var_type = "MooseVariable";
455  else if (_family == "MONOMIAL")
456  aux_var_type = "MooseVariableConstMonomial";
457  else if (_family == "SCALAR")
458  aux_var_type = "MooseVariableScalar";
459  else
460  mooseError("Unsupported finite element family in, " + name() +
461  ". Please use LAGRANGE, MONOMIAL, or SCALAR");
462 
463  auto params = _factory.getValidParams(aux_var_type);
464  params.set<MooseEnum>("order") = _order;
465  params.set<MooseEnum>("family") = _family;
466 
468  {
469  std::ostringstream av_name_stream;
470  av_name_stream << av_base_name << "_" << _ring_vec[ring_index];
471  _problem->addAuxVariable(aux_var_type, av_name_stream.str(), params);
472  }
473  else
474  {
475  for (unsigned int cfp_index = 0; cfp_index < num_crack_front_points; ++cfp_index)
476  {
477  std::ostringstream av_name_stream;
478  av_name_stream << av_base_name << "_" << cfp_index + 1 << "_" << _ring_vec[ring_index];
479  _problem->addAuxVariable(aux_var_type, av_name_stream.str(), params);
480  }
481  }
482  }
483  }
484 
485  else if (_current_task == "add_aux_kernel" && _output_q)
486  {
487  std::string ak_type_name;
488  unsigned int nrings = 0;
489  if (_q_function_type == GEOMETRY)
490  {
491  ak_type_name = "DomainIntegralQFunction";
492  nrings = _ring_vec.size();
493  }
494  else if (_q_function_type == TOPOLOGY)
495  {
496  ak_type_name = "DomainIntegralTopologicalQFunction";
497  nrings = _ring_last - _ring_first + 1;
498  }
499 
500  InputParameters params = _factory.getValidParams(ak_type_name);
501  params.set<ExecFlagEnum>("execute_on") = {EXEC_INITIAL, EXEC_TIMESTEP_END};
502  params.set<UserObjectName>("crack_front_definition") = uo_name;
503  params.set<bool>("use_displaced_mesh") = _use_displaced_mesh;
504 
505  for (unsigned int ring_index = 0; ring_index < nrings; ++ring_index)
506  {
507  if (_q_function_type == GEOMETRY)
508  {
509  params.set<Real>("j_integral_radius_inner") = _radius_inner[ring_index];
510  params.set<Real>("j_integral_radius_outer") = _radius_outer[ring_index];
511  }
512  else if (_q_function_type == TOPOLOGY)
513  {
514  params.set<unsigned int>("ring_index") = _ring_first + ring_index;
515  }
516 
518  {
519  std::ostringstream ak_name_stream;
520  ak_name_stream << ak_base_name << "_" << _ring_vec[ring_index];
521  std::ostringstream av_name_stream;
522  av_name_stream << av_base_name << "_" << _ring_vec[ring_index];
523  params.set<AuxVariableName>("variable") = av_name_stream.str();
524  _problem->addAuxKernel(ak_type_name, ak_name_stream.str(), params);
525  }
526  else
527  {
528  for (unsigned int cfp_index = 0; cfp_index < num_crack_front_points; ++cfp_index)
529  {
530  std::ostringstream ak_name_stream;
531  ak_name_stream << ak_base_name << "_" << cfp_index + 1 << "_" << _ring_vec[ring_index];
532  std::ostringstream av_name_stream;
533  av_name_stream << av_base_name << "_" << cfp_index + 1 << "_" << _ring_vec[ring_index];
534  params.set<AuxVariableName>("variable") = av_name_stream.str();
535  params.set<unsigned int>("crack_front_point_index") = cfp_index;
536  _problem->addAuxKernel(ak_type_name, ak_name_stream.str(), params);
537  }
538  }
539  }
540  }
541 
542  else if (_current_task == "add_postprocessor")
543  {
544  for (std::set<INTEGRAL>::iterator sit = _integrals.begin(); sit != _integrals.end(); ++sit)
545  {
546  std::string pp_base_name;
547  switch (*sit)
548  {
549  case J_INTEGRAL:
550  pp_base_name = "J";
551  break;
552 
553  case C_INTEGRAL:
554  pp_base_name = "C";
555  break;
556 
557  case K_FROM_J_INTEGRAL:
558  pp_base_name = "K";
559  break;
560 
562  pp_base_name = "II_KI";
563  break;
564 
566  pp_base_name = "II_KII";
567  break;
568 
570  pp_base_name = "II_KIII";
571  break;
572 
574  pp_base_name = "II_T";
575  break;
576  }
577  const std::string pp_type_name("VectorPostprocessorComponent");
578  InputParameters params = _factory.getValidParams(pp_type_name);
579  for (unsigned int ring_index = 0; ring_index < _ring_vec.size(); ++ring_index)
580  {
582  {
583  params.set<VectorPostprocessorName>("vectorpostprocessor") =
584  pp_base_name + "_2DVPP_" + Moose::stringify(_ring_vec[ring_index]);
585  std::string pp_name = pp_base_name + +"_" + Moose::stringify(_ring_vec[ring_index]);
586  params.set<unsigned int>("index") = 0;
587  params.set<std::string>("vector_name") =
588  pp_base_name + "_" + Moose::stringify(_ring_vec[ring_index]);
589  _problem->addPostprocessor(pp_type_name, pp_name, params);
590  }
591  else
592  {
593  for (unsigned int cfp_index = 0; cfp_index < num_crack_front_points; ++cfp_index)
594  {
595  params.set<VectorPostprocessorName>("vectorpostprocessor") =
596  pp_base_name + "_" + Moose::stringify(_ring_vec[ring_index]);
597  std::string pp_name = pp_base_name + "_" + Moose::stringify(cfp_index + 1) + "_" +
598  Moose::stringify(_ring_vec[ring_index]);
599  params.set<unsigned int>("index") = cfp_index;
600  params.set<std::string>("vector_name") =
601  pp_base_name + "_" + Moose::stringify(_ring_vec[ring_index]);
602  _problem->addPostprocessor(pp_type_name, pp_name, params);
603  }
604  }
605  }
606  }
607 
608  if (_get_equivalent_k)
609  {
610  std::string pp_base_name("Keq");
611  const std::string pp_type_name("VectorPostprocessorComponent");
612  InputParameters params = _factory.getValidParams(pp_type_name);
613  for (unsigned int ring_index = 0; ring_index < _ring_vec.size(); ++ring_index)
614  {
616  {
617  params.set<VectorPostprocessorName>("vectorpostprocessor") =
618  pp_base_name + "_2DVPP_" + Moose::stringify(_ring_vec[ring_index]);
619  std::string pp_name = pp_base_name + +"_" + Moose::stringify(_ring_vec[ring_index]);
620  params.set<unsigned int>("index") = 0;
621  params.set<std::string>("vector_name") =
622  pp_base_name + "_" + Moose::stringify(_ring_vec[ring_index]);
623  _problem->addPostprocessor(pp_type_name, pp_name, params);
624  }
625  else
626  {
627  for (unsigned int cfp_index = 0; cfp_index < num_crack_front_points; ++cfp_index)
628  {
629  params.set<VectorPostprocessorName>("vectorpostprocessor") =
630  pp_base_name + "_" + Moose::stringify(_ring_vec[ring_index]);
631  std::string pp_name = pp_base_name + "_" + Moose::stringify(cfp_index + 1) + "_" +
632  Moose::stringify(_ring_vec[ring_index]);
633  params.set<unsigned int>("index") = cfp_index;
634  params.set<std::string>("vector_name") =
635  pp_base_name + "_" + Moose::stringify(_ring_vec[ring_index]);
636  _problem->addPostprocessor(pp_type_name, pp_name, params);
637  }
638  }
639  }
640  }
641 
642  for (unsigned int i = 0; i < _output_variables.size(); ++i)
643  {
644  const std::string ov_base_name(_output_variables[i]);
645  const std::string pp_type_name("CrackFrontData");
646  InputParameters params = _factory.getValidParams(pp_type_name);
647  params.set<ExecFlagEnum>("execute_on") = EXEC_TIMESTEP_END;
648  params.set<UserObjectName>("crack_front_definition") = uo_name;
650  {
651  std::ostringstream pp_name_stream;
652  pp_name_stream << ov_base_name << "_crack";
653  params.set<VariableName>("variable") = _output_variables[i];
654  _problem->addPostprocessor(pp_type_name, pp_name_stream.str(), params);
655  }
656  else
657  {
658  for (unsigned int cfp_index = 0; cfp_index < num_crack_front_points; ++cfp_index)
659  {
660  std::ostringstream pp_name_stream;
661  pp_name_stream << ov_base_name << "_crack_" << cfp_index + 1;
662  params.set<VariableName>("variable") = _output_variables[i];
663  params.set<unsigned int>("crack_front_point_index") = cfp_index;
664  _problem->addPostprocessor(pp_type_name, pp_name_stream.str(), params);
665  }
666  }
667  }
668  }
669 
670  else if (_current_task == "add_vector_postprocessor")
671  {
672  if (_integrals.count(J_INTEGRAL) != 0 || _integrals.count(C_INTEGRAL) != 0 ||
673  _integrals.count(K_FROM_J_INTEGRAL) != 0)
674  {
675  std::string vpp_base_name;
676  std::string jintegral_selection = "JIntegral";
677 
678  if (_integrals.count(J_INTEGRAL) != 0)
679  {
680  vpp_base_name = "J";
681  jintegral_selection = "JIntegral";
682  }
683  else if (_integrals.count(K_FROM_J_INTEGRAL) != 0)
684  {
685  vpp_base_name = "K";
686  jintegral_selection = "KFromJIntegral";
687  }
688  else if (_integrals.count(C_INTEGRAL) != 0)
689  {
690  vpp_base_name = "C";
691  jintegral_selection = "CIntegral";
692  }
693 
695  vpp_base_name += "_2DVPP";
696 
697  const std::string vpp_type_name("JIntegral");
698  InputParameters params = _factory.getValidParams(vpp_type_name);
699  if (!getParam<bool>("output_vpp"))
700  params.set<std::vector<OutputName>>("outputs") = {"none"};
701 
702  params.set<ExecFlagEnum>("execute_on") = EXEC_TIMESTEP_END;
703  params.set<UserObjectName>("crack_front_definition") = uo_name;
704  params.set<std::vector<SubdomainName>>("block") = {_blocks};
705  params.set<MooseEnum>("position_type") = _position_type;
706 
707  if (_integrals.count(K_FROM_J_INTEGRAL) != 0)
708  {
709  params.set<Real>("youngs_modulus") = _youngs_modulus;
710  params.set<Real>("poissons_ratio") = _poissons_ratio;
711  }
712 
714  params.set<unsigned int>("symmetry_plane") = _symmetry_plane;
715 
716  // Select the integral type to be computed in JIntegral vector postprocessor
717  params.set<MooseEnum>("integral") = jintegral_selection;
718 
719  params.set<std::vector<VariableName>>("displacements") = _displacements;
720  params.set<bool>("use_displaced_mesh") = _use_displaced_mesh;
721  for (unsigned int ring_index = 0; ring_index < _ring_vec.size(); ++ring_index)
722  {
723  params.set<unsigned int>("ring_index") = _ring_vec[ring_index];
724  params.set<MooseEnum>("q_function_type") = _q_function_type;
725 
726  std::string vpp_name = vpp_base_name + "_" + Moose::stringify(_ring_vec[ring_index]);
727  _problem->addVectorPostprocessor(vpp_type_name, vpp_name, params);
728  }
729  }
730 
731  if (_integrals.count(INTERACTION_INTEGRAL_KI) != 0 ||
732  _integrals.count(INTERACTION_INTEGRAL_KII) != 0 ||
735  {
738  paramError("symmetry_plane",
739  "In DomainIntegral, symmetry_plane option cannot be used with mode-II or "
740  "mode-III interaction integral");
741 
742  std::string vpp_base_name;
743  std::string vpp_type_name(ad_prepend + "InteractionIntegral");
744 
745  InputParameters params = _factory.getValidParams(vpp_type_name);
746  if (!getParam<bool>("output_vpp"))
747  params.set<std::vector<OutputName>>("outputs") = {"none"};
748 
749  if (_use_crack_front_points_provider && !xfem_exec_flags.empty())
750  params.set<ExecFlagEnum>("execute_on") = xfem_exec_flags;
751  else
752  params.set<ExecFlagEnum>("execute_on") = {EXEC_TIMESTEP_END};
753 
754  if (isParamValid("additional_eigenstrain_00") && isParamValid("additional_eigenstrain_01") &&
755  isParamValid("additional_eigenstrain_11") && isParamValid("additional_eigenstrain_22"))
756  {
757  params.set<CoupledName>("additional_eigenstrain_00") =
758  getParam<CoupledName>("additional_eigenstrain_00");
759  params.set<CoupledName>("additional_eigenstrain_01") =
760  getParam<CoupledName>("additional_eigenstrain_01");
761  params.set<CoupledName>("additional_eigenstrain_11") =
762  getParam<CoupledName>("additional_eigenstrain_11");
763  params.set<CoupledName>("additional_eigenstrain_22") =
764  getParam<CoupledName>("additional_eigenstrain_22");
765  }
766 
767  params.set<UserObjectName>("crack_front_definition") = uo_name;
768  params.set<bool>("use_displaced_mesh") = _use_displaced_mesh;
769  params.set<std::vector<SubdomainName>>("block") = {_blocks};
770 
772  params.set<unsigned int>("symmetry_plane") = _symmetry_plane;
773 
774  if (_fgm_crack)
775  {
776  params.set<MaterialPropertyName>(
777  "functionally_graded_youngs_modulus_crack_dir_gradient") = {
779  params.set<MaterialPropertyName>("functionally_graded_youngs_modulus") = {
781  }
782 
783  params.set<Real>("poissons_ratio") = _poissons_ratio;
784  params.set<Real>("youngs_modulus") = _youngs_modulus;
785  params.set<std::vector<VariableName>>("displacements") = _displacements;
786  if (_temp != "")
787  params.set<std::vector<VariableName>>("temperature") = {_temp};
788 
789  if (parameters().isParamValid("eigenstrain_gradient"))
790  params.set<MaterialPropertyName>("eigenstrain_gradient") =
791  parameters().get<MaterialPropertyName>("eigenstrain_gradient");
792  if (parameters().isParamValid("body_force"))
793  params.set<MaterialPropertyName>("body_force") =
794  parameters().get<MaterialPropertyName>("body_force");
795 
796  for (std::set<INTEGRAL>::iterator sit = _integrals.begin(); sit != _integrals.end(); ++sit)
797  {
798  switch (*sit)
799  {
800  case J_INTEGRAL:
801  continue;
802 
803  case C_INTEGRAL:
804  continue;
805 
806  case K_FROM_J_INTEGRAL:
807  continue;
808 
810  vpp_base_name = "II_KI";
811  params.set<Real>("K_factor") =
812  0.5 * _youngs_modulus / (1.0 - std::pow(_poissons_ratio, 2.0));
813  params.set<MooseEnum>("sif_mode") = "KI";
814  break;
815 
817  vpp_base_name = "II_KII";
818  params.set<Real>("K_factor") =
819  0.5 * _youngs_modulus / (1.0 - std::pow(_poissons_ratio, 2.0));
820  params.set<MooseEnum>("sif_mode") = "KII";
821  break;
822 
824  vpp_base_name = "II_KIII";
825  params.set<Real>("K_factor") = 0.5 * _youngs_modulus / (1.0 + _poissons_ratio);
826  params.set<MooseEnum>("sif_mode") = "KIII";
827  break;
828 
830  vpp_base_name = "II_T";
831  params.set<Real>("K_factor") = _youngs_modulus / (1 - std::pow(_poissons_ratio, 2));
832  params.set<MooseEnum>("sif_mode") = "T";
833  break;
834  }
836  vpp_base_name += "_2DVPP";
837  for (unsigned int ring_index = 0; ring_index < _ring_vec.size(); ++ring_index)
838  {
839  params.set<unsigned int>("ring_index") = _ring_vec[ring_index];
840  params.set<MooseEnum>("q_function_type") = _q_function_type;
841 
842  std::string vpp_name = vpp_base_name + "_" + Moose::stringify(_ring_vec[ring_index]);
843  _problem->addVectorPostprocessor(vpp_type_name, vpp_name, params);
844  }
845  }
846  }
847 
848  if (_get_equivalent_k)
849  {
850  std::string vpp_base_name("Keq");
852  vpp_base_name += "_2DVPP";
853  const std::string vpp_type_name("MixedModeEquivalentK");
854  InputParameters params = _factory.getValidParams(vpp_type_name);
855  params.set<ExecFlagEnum>("execute_on") = EXEC_TIMESTEP_END;
856  params.set<Real>("poissons_ratio") = _poissons_ratio;
857 
858  for (unsigned int ring_index = 0; ring_index < _ring_vec.size(); ++ring_index)
859  {
860  std::string ki_name = "II_KI_";
861  std::string kii_name = "II_KII_";
862  std::string kiii_name = "II_KIII_";
863  params.set<unsigned int>("ring_index") = _ring_vec[ring_index];
865  {
866  params.set<VectorPostprocessorName>("KI_vectorpostprocessor") =
867  ki_name + "2DVPP_" + Moose::stringify(_ring_vec[ring_index]);
868  params.set<VectorPostprocessorName>("KII_vectorpostprocessor") =
869  kii_name + "2DVPP_" + Moose::stringify(_ring_vec[ring_index]);
870  params.set<VectorPostprocessorName>("KIII_vectorpostprocessor") =
871  kiii_name + "2DVPP_" + Moose::stringify(_ring_vec[ring_index]);
872  }
873  else
874  {
875  params.set<VectorPostprocessorName>("KI_vectorpostprocessor") =
876  ki_name + Moose::stringify(_ring_vec[ring_index]);
877  params.set<VectorPostprocessorName>("KII_vectorpostprocessor") =
878  kii_name + Moose::stringify(_ring_vec[ring_index]);
879  params.set<VectorPostprocessorName>("KIII_vectorpostprocessor") =
880  kiii_name + Moose::stringify(_ring_vec[ring_index]);
881  }
882  params.set<std::string>("KI_vector_name") =
883  ki_name + Moose::stringify(_ring_vec[ring_index]);
884  params.set<std::string>("KII_vector_name") =
885  kii_name + Moose::stringify(_ring_vec[ring_index]);
886  params.set<std::string>("KIII_vector_name") =
887  kiii_name + Moose::stringify(_ring_vec[ring_index]);
888  std::string vpp_name = vpp_base_name + "_" + Moose::stringify(_ring_vec[ring_index]);
889  _problem->addVectorPostprocessor(vpp_type_name, vpp_name, params);
890  }
891  }
892 
894  {
895  for (unsigned int i = 0; i < _output_variables.size(); ++i)
896  {
897  const std::string vpp_type_name("VectorOfPostprocessors");
898  InputParameters params = _factory.getValidParams(vpp_type_name);
899  params.set<ExecFlagEnum>("execute_on") = EXEC_TIMESTEP_END;
900  std::ostringstream vpp_name_stream;
901  vpp_name_stream << _output_variables[i] << "_crack";
902  std::vector<PostprocessorName> postprocessor_names;
903  for (unsigned int cfp_index = 0; cfp_index < num_crack_front_points; ++cfp_index)
904  {
905  std::ostringstream pp_name_stream;
906  pp_name_stream << vpp_name_stream.str() << "_" << cfp_index + 1;
907  postprocessor_names.push_back(pp_name_stream.str());
908  }
909  params.set<std::vector<PostprocessorName>>("postprocessors") = postprocessor_names;
910  _problem->addVectorPostprocessor(vpp_type_name, vpp_name_stream.str(), params);
911  }
912  }
913  }
914 
915  else if (_current_task == "add_material")
916  {
917  if (_temp != "")
918  {
919  std::string mater_name;
920  const std::string mater_type_name("ThermalFractureIntegral");
921  mater_name = "ThermalFractureIntegral";
922 
923  InputParameters params = _factory.getValidParams(mater_type_name);
924  params.set<std::vector<MaterialPropertyName>>("eigenstrain_names") =
925  getParam<std::vector<MaterialPropertyName>>("eigenstrain_names");
926  params.set<std::vector<VariableName>>("temperature") = {_temp};
927  params.set<std::vector<SubdomainName>>("block") = {_blocks};
928  _problem->addMaterial(mater_type_name, mater_name, params);
929  }
930  MultiMooseEnum integral_moose_enums = getParam<MultiMooseEnum>("integrals");
931  bool have_j_integral = false;
932  bool have_c_integral = false;
933 
934  for (auto ime : integral_moose_enums)
935  {
936  if (ime == "JIntegral" || ime == "CIntegral" || ime == "KFromJIntegral" ||
937  ime == "InteractionIntegralKI" || ime == "InteractionIntegralKII" ||
938  ime == "InteractionIntegralKIII" || ime == "InteractionIntegralT")
939  have_j_integral = true;
940 
941  if (ime == "CIntegral")
942  have_c_integral = true;
943  }
944  if (have_j_integral)
945  {
946  std::string mater_name;
947  const std::string mater_type_name(ad_prepend + "StrainEnergyDensity");
948  mater_name = ad_prepend + "StrainEnergyDensity";
949 
950  InputParameters params = _factory.getValidParams(mater_type_name);
951  _incremental = getParam<bool>("incremental");
952  params.set<bool>("incremental") = _incremental;
953  params.set<std::vector<SubdomainName>>("block") = {_blocks};
954  _problem->addMaterial(mater_type_name, mater_name, params);
955 
956  {
957  std::string mater_name;
958  const std::string mater_type_name(ad_prepend + "EshelbyTensor");
959  mater_name = ad_prepend + "EshelbyTensor";
960 
961  InputParameters params = _factory.getValidParams(mater_type_name);
962  _displacements = getParam<std::vector<VariableName>>("displacements");
963  params.set<std::vector<VariableName>>("displacements") = _displacements;
964  params.set<std::vector<SubdomainName>>("block") = {_blocks};
965 
966  if (have_c_integral)
967  params.set<bool>("compute_dissipation") = true;
968 
969  if (_temp != "")
970  params.set<std::vector<VariableName>>("temperature") = {_temp};
971 
972  _problem->addMaterial(mater_type_name, mater_name, params);
973  }
974  // Strain energy rate density needed for C(t)/C* integral
975  if (have_c_integral)
976  {
977  std::string mater_name;
978  const std::string mater_type_name(ad_prepend + "StrainEnergyRateDensity");
979  mater_name = ad_prepend + "StrainEnergyRateDensity";
980 
981  InputParameters params = _factory.getValidParams(mater_type_name);
982  params.set<std::vector<SubdomainName>>("block") = {_blocks};
983  params.set<std::vector<MaterialName>>("inelastic_models") =
984  getParam<std::vector<MaterialName>>("inelastic_models");
985 
986  _problem->addMaterial(mater_type_name, mater_name, params);
987  }
988  }
989  }
990 }
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...
std::vector< std::pair< R1, R2 > > get(const std::string &param1, const std::string &param2) 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...
static ExecFlagRegistry & getExecFlagRegistry()
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.
virtual const std::string & name() const
bool isParamValid(const std::string &name) const
Factory & _factory
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.
const T & getParam(const std::string &name) const
const std::string & _current_task
void paramError(const std::string &param, Args... args) const
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...
const ExecFlagType EXEC_NONLINEAR
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
const InputParameters & parameters() const
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)
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 1028 of file DomainIntegralAction.C.

1029 {
1030  if (_integrals.count(INTERACTION_INTEGRAL_T) != 0)
1031  {
1032  InputParameters params = _factory.getValidParams("CrackFrontDefinition");
1033  addRelationshipManagers(input_rm_type, params);
1034  }
1035 }
InputParameters getValidParams(const std::string &name) const
Factory & _factory
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 993 of file DomainIntegralAction.C.

Referenced by act().

994 {
995  unsigned int num_points = 0;
996  if (_boundary_names.size() != 0)
997  {
998  std::vector<BoundaryID> bids = _mesh->getBoundaryIDs(_boundary_names, true);
999  std::set<unsigned int> nodes;
1000 
1001  ConstBndNodeRange & bnd_nodes = *_mesh->getBoundaryNodeRange();
1002  for (ConstBndNodeRange::const_iterator nd = bnd_nodes.begin(); nd != bnd_nodes.end(); ++nd)
1003  {
1004  const BndNode * bnode = *nd;
1005  BoundaryID boundary_id = bnode->_bnd_id;
1006 
1007  for (unsigned int ibid = 0; ibid < bids.size(); ++ibid)
1008  {
1009  if (boundary_id == bids[ibid])
1010  {
1011  nodes.insert(bnode->_node->id());
1012  break;
1013  }
1014  }
1015  }
1016  num_points = nodes.size();
1017  }
1018  else if (_crack_front_points.size() != 0)
1019  num_points = _crack_front_points.size();
1021  num_points = getParam<unsigned int>("number_points_from_provider");
1022  else
1023  mooseError("Must define either 'boundary' or 'crack_front_points'");
1024  return num_points;
1025 }
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 addCrackFrontDefinitionParams(InputParameters &params)
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
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: