www.mooseframework.org
Functions
RankTwoScalarTools Namespace Reference

Functions

MooseEnum scalarOptions ()
 
template<typename T >
component (const RankTwoTensorTempl< T > &r2tensor, unsigned int i, unsigned int j)
 
template<typename T >
component (const RankTwoTensorTempl< T > &r2tensor, unsigned int i, unsigned int j, Point &direction)
 
template<typename T >
vonMisesStress (const RankTwoTensorTempl< T > &stress)
 
template<typename T >
effectiveStrain (const RankTwoTensorTempl< T > &strain)
 
template<typename T >
hydrostatic (const RankTwoTensorTempl< T > &r2tensor)
 
template<typename T >
L2norm (const RankTwoTensorTempl< T > &r2tensor)
 
template<typename T >
volumetricStrain (const RankTwoTensorTempl< T > &strain)
 
template<typename T >
firstInvariant (const RankTwoTensorTempl< T > &r2tensor)
 
template<typename T >
secondInvariant (const RankTwoTensorTempl< T > &r2tensor)
 
template<typename T >
thirdInvariant (const RankTwoTensorTempl< T > &r2tensor)
 
template<typename T >
calcEigenValuesEigenVectors (const RankTwoTensorTempl< T > &r2tensor, unsigned int index, Point &eigenvec)
 
template<typename T >
maxPrincipal (const RankTwoTensorTempl< T > &r2tensor, Point &direction)
 
template<typename T >
midPrincipal (const RankTwoTensorTempl< T > &r2tensor, Point &direction)
 
template<typename T >
minPrincipal (const RankTwoTensorTempl< T > &r2tensor, Point &direction)
 
template<typename T >
axialStress (const RankTwoTensorTempl< T > &stress, const Point &point1, const Point &point2, Point &direction)
 
void normalPositionVector (const Point &point1, const Point &point2, const Point &curr_point, Point &normalPosition)
 
template<typename T >
hoopStress (const RankTwoTensorTempl< T > &stress, const Point &point1, const Point &point2, const Point &curr_point, Point &direction)
 
template<typename T >
radialStress (const RankTwoTensorTempl< T > &stress, const Point &point1, const Point &point2, const Point &curr_point, Point &direction)
 
template<typename T >
directionValueTensor (const RankTwoTensorTempl< T > &r2tensor, Point &direction)
 
template<typename T >
triaxialityStress (const RankTwoTensorTempl< T > &stress)
 
template<typename T >
maxShear (const RankTwoTensorTempl< T > &stress)
 
template<typename T >
stressIntensity (const RankTwoTensorTempl< T > &stress)
 
template<typename T >
getQuantity (const RankTwoTensorTempl< T > &tensor, const MooseEnum &scalar_type, const Point &point1, const Point &point2, const Point &curr_point, Point &direction)
 

Function Documentation

◆ axialStress()

template<typename T >
T RankTwoScalarTools::axialStress ( const RankTwoTensorTempl< T > &  stress,
const Point &  point1,
const Point &  point2,
Point &  direction 
)

Definition at line 258 of file RankTwoScalarTools.h.

262 {
263  Point axis = point2 - point1;
264  axis /= axis.norm();
265 
266  // Calculate the stress in the direction of the axis specifed by the user
267  T axial_stress = 0.0;
268  for (unsigned int i = 0; i < 3; ++i)
269  for (unsigned int j = 0; j < 3; ++j)
270  axial_stress += axis(j) * stress(j, i) * axis(i);
271 
272  direction = axis;
273 
274  return axial_stress;
275 }

Referenced by RankTwoBasedFailureCriteriaNOSPD::computeFailureCriterionValue(), and getQuantity().

◆ calcEigenValuesEigenVectors()

template<typename T >
T RankTwoScalarTools::calcEigenValuesEigenVectors ( const RankTwoTensorTempl< T > &  r2tensor,
unsigned int  index,
Point &  eigenvec 
)

Definition at line 195 of file RankTwoScalarTools.h.

198 {
199  std::vector<T> eigenval(LIBMESH_DIM);
200  RankTwoTensorTempl<T> eigvecs;
201  r2tensor.symmetricEigenvaluesEigenvectors(eigenval, eigvecs);
202 
203  T val = eigenval[index];
204  eigenvec = eigvecs.column(index);
205 
206  return val;
207 }

Referenced by maxPrincipal(), midPrincipal(), and minPrincipal().

◆ component() [1/2]

template<typename T >
T RankTwoScalarTools::component ( const RankTwoTensorTempl< T > &  r2tensor,
unsigned int  i,
unsigned int  j 
)

◆ component() [2/2]

template<typename T >
T RankTwoScalarTools::component ( const RankTwoTensorTempl< T > &  r2tensor,
unsigned int  i,
unsigned int  j,
Point &  direction 
)

Definition at line 42 of file RankTwoScalarTools.h.

43 {
44  direction.zero();
45  if (i == j)
46  direction(i) = 1.0;
47  else
48  {
49  direction(i) = std::sqrt(0.5);
50  direction(j) = std::sqrt(0.5);
51  }
52 
53  return r2tensor(i, j);
54 }

◆ directionValueTensor()

template<typename T >
T RankTwoScalarTools::directionValueTensor ( const RankTwoTensorTempl< T > &  r2tensor,
Point &  direction 
)

Definition at line 366 of file RankTwoScalarTools.h.

367 {
368  T tensor_value_in_direction = 0.0;
369  for (unsigned int i = 0; i < 3; ++i)
370  for (unsigned int j = 0; j < 3; ++j)
371  tensor_value_in_direction += direction(j) * r2tensor(j, i) * direction(i);
372 
373  return tensor_value_in_direction;
374 }

Referenced by getQuantity().

◆ effectiveStrain()

template<typename T >
T RankTwoScalarTools::effectiveStrain ( const RankTwoTensorTempl< T > &  strain)

Definition at line 77 of file RankTwoScalarTools.h.

78 {
79  return std::sqrt(2.0 / 3.0 * strain.doubleContraction(strain));
80 }

Referenced by RankTwoBasedFailureCriteriaNOSPD::computeFailureCriterionValue(), and getQuantity().

◆ firstInvariant()

template<typename T >
T RankTwoScalarTools::firstInvariant ( const RankTwoTensorTempl< T > &  r2tensor)

Definition at line 139 of file RankTwoScalarTools.h.

140 {
141  return r2tensor.trace();
142 }

Referenced by getQuantity().

◆ getQuantity()

template<typename T >
T RankTwoScalarTools::getQuantity ( const RankTwoTensorTempl< T > &  tensor,
const MooseEnum &  scalar_type,
const Point &  point1,
const Point &  point2,
const Point &  curr_point,
Point &  direction 
)

Definition at line 420 of file RankTwoScalarTools.h.

426 {
427  switch (scalar_type)
428  {
429  case 0:
430  return vonMisesStress(tensor);
431  case 1:
432  return effectiveStrain(tensor);
433  case 2:
434  return hydrostatic(tensor);
435  case 3:
436  return L2norm(tensor);
437  case 4:
438  return maxPrincipal(tensor, direction);
439  case 5:
440  return midPrincipal(tensor, direction);
441  case 6:
442  return minPrincipal(tensor, direction);
443  case 7:
444  return volumetricStrain(tensor);
445  case 8:
446  return firstInvariant(tensor);
447  case 9:
448  return secondInvariant(tensor);
449  case 10:
450  return thirdInvariant(tensor);
451  case 11:
452  return axialStress(tensor, point1, point2, direction);
453  case 12:
454  return hoopStress(tensor, point1, point2, curr_point, direction);
455  case 13:
456  return radialStress(tensor, point1, point2, curr_point, direction);
457  case 14:
458  return triaxialityStress(tensor);
459  case 15:
460  return directionValueTensor(tensor, direction);
461  case 16:
462  return maxShear(tensor);
463  case 17:
464  return stressIntensity(tensor);
465  default:
466  mooseError("RankTwoScalarAux Error: Pass valid scalar type - " +
467  scalarOptions().getRawNames());
468  }
469 }

Referenced by NodalRankTwoPD::computeValue(), RankTwoScalarAux::computeValue(), XFEMRankTwoTensorMarkerUserObject::doesElementCrack(), NodalRankTwoScalarPD::gatherWeightedValue(), and LineMaterialRankTwoScalarSampler::getScalarFromProperty().

◆ hoopStress()

template<typename T >
T RankTwoScalarTools::hoopStress ( const RankTwoTensorTempl< T > &  stress,
const Point &  point1,
const Point &  point2,
const Point &  curr_point,
Point &  direction 
)

Definition at line 303 of file RankTwoScalarTools.h.

308 {
309  // Calculate the cross of the normal to the axis of rotation from the current
310  // location and the axis of rotation
311  Point xp;
312  normalPositionVector(point1, point2, curr_point, xp);
313  Point axis_rotation = point2 - point1;
314  Point yp = axis_rotation / axis_rotation.norm();
315  Point zp = xp.cross(yp);
316 
317  // Calculate the scalar value of the hoop stress
318  T hoop_stress = 0.0;
319  for (unsigned int i = 0; i < 3; ++i)
320  for (unsigned int j = 0; j < 3; ++j)
321  hoop_stress += zp(j) * stress(j, i) * zp(i);
322 
323  direction = zp;
324 
325  return hoop_stress;
326 }

Referenced by getQuantity().

◆ hydrostatic()

template<typename T >
T RankTwoScalarTools::hydrostatic ( const RankTwoTensorTempl< T > &  r2tensor)

Definition at line 88 of file RankTwoScalarTools.h.

89 {
90  return r2tensor.trace() / 3.0;
91 }

Referenced by getQuantity(), and triaxialityStress().

◆ L2norm()

template<typename T >
T RankTwoScalarTools::L2norm ( const RankTwoTensorTempl< T > &  r2tensor)

◆ maxPrincipal()

template<typename T >
T RankTwoScalarTools::maxPrincipal ( const RankTwoTensorTempl< T > &  r2tensor,
Point &  direction 
)

Definition at line 217 of file RankTwoScalarTools.h.

218 {
219  return calcEigenValuesEigenVectors(r2tensor, (LIBMESH_DIM - 1), direction);
220 }

Referenced by RankTwoBasedFailureCriteriaNOSPD::computeFailureCriterionValue(), getQuantity(), and maxShear().

◆ maxShear()

template<typename T >
T RankTwoScalarTools::maxShear ( const RankTwoTensorTempl< T > &  stress)

Definition at line 392 of file RankTwoScalarTools.h.

393 {
394  Point dummy;
395  return (maxPrincipal(stress, dummy) - minPrincipal(stress, dummy)) / 2.;
396 }

Referenced by RankTwoBasedFailureCriteriaNOSPD::computeFailureCriterionValue(), getQuantity(), and stressIntensity().

◆ midPrincipal()

template<typename T >
T RankTwoScalarTools::midPrincipal ( const RankTwoTensorTempl< T > &  r2tensor,
Point &  direction 
)

Definition at line 230 of file RankTwoScalarTools.h.

231 {
232  return calcEigenValuesEigenVectors(r2tensor, 1, direction);
233 }

Referenced by getQuantity().

◆ minPrincipal()

template<typename T >
T RankTwoScalarTools::minPrincipal ( const RankTwoTensorTempl< T > &  r2tensor,
Point &  direction 
)

Definition at line 243 of file RankTwoScalarTools.h.

244 {
245  return calcEigenValuesEigenVectors(r2tensor, 0, direction);
246 }

Referenced by getQuantity(), and maxShear().

◆ normalPositionVector()

void RankTwoScalarTools::normalPositionVector ( const Point &  point1,
const Point &  point2,
const Point &  curr_point,
Point &  normalPosition 
)

Definition at line 25 of file RankTwoScalarTools.C.

29 {
30  // Find the nearest point on the axis of rotation (defined by point2 - point1)
31  // to the current position, e.g. the normal to the axis of rotation at the
32  // current position
33  Point axis_rotation = point2 - point1;
34  Point positionWRTpoint1 = point1 - curr_point;
35  Real projection = (axis_rotation * positionWRTpoint1) / axis_rotation.norm_sq();
36  Point normal = point1 - projection * axis_rotation;
37 
38  // Calculate the direction normal to the plane formed by the axis of rotation
39  // and the normal to the axis of rotation from the current position.
40  normalPosition = curr_point - normal;
41  normalPosition /= normalPosition.norm();
42 }

Referenced by hoopStress(), and radialStress().

◆ radialStress()

template<typename T >
T RankTwoScalarTools::radialStress ( const RankTwoTensorTempl< T > &  stress,
const Point &  point1,
const Point &  point2,
const Point &  curr_point,
Point &  direction 
)

Definition at line 339 of file RankTwoScalarTools.h.

344 {
345  Point radial_norm;
346  normalPositionVector(point1, point2, curr_point, radial_norm);
347 
348  // Compute the scalar stress component in the direction of the normal vector from the
349  // user-defined axis of rotation.
350  T radial_stress = 0.0;
351  for (unsigned int i = 0; i < 3; ++i)
352  for (unsigned int j = 0; j < 3; ++j)
353  radial_stress += radial_norm(j) * stress(j, i) * radial_norm(i);
354 
355  direction = radial_norm;
356 
357  return radial_stress;
358 }

Referenced by getQuantity().

◆ scalarOptions()

MooseEnum RankTwoScalarTools::scalarOptions ( )

Definition at line 16 of file RankTwoScalarTools.C.

17 {
18  return MooseEnum("VonMisesStress EffectiveStrain Hydrostatic L2norm MaxPrincipal "
19  "MidPrincipal MinPrincipal VolumetricStrain FirstInvariant SecondInvariant "
20  "ThirdInvariant AxialStress HoopStress RadialStress TriaxialityStress "
21  "Direction MaxShear StressIntensity");
22 }

Referenced by getQuantity(), RankTwoScalarAux::validParams(), LineMaterialRankTwoScalarSampler::validParams(), validParams< NodalRankTwoPD >(), validParams< NodalRankTwoScalarPD >(), and validParams< XFEMRankTwoTensorMarkerUserObject >().

◆ secondInvariant()

template<typename T >
T RankTwoScalarTools::secondInvariant ( const RankTwoTensorTempl< T > &  r2tensor)

Definition at line 153 of file RankTwoScalarTools.h.

154 {
155  T val = 0.0;
156  for (unsigned int i = 0; i < 2; ++i)
157  for (unsigned int j = i + 1; j < 3; ++j)
158  {
159  val += r2tensor(i, i) * r2tensor(j, j);
160  val -= (r2tensor(i, j) * r2tensor(i, j) + r2tensor(j, i) * r2tensor(j, i)) * 0.5;
161  }
162 
163  return val;
164 }

Referenced by getQuantity().

◆ stressIntensity()

template<typename T >
T RankTwoScalarTools::stressIntensity ( const RankTwoTensorTempl< T > &  stress)

Definition at line 402 of file RankTwoScalarTools.h.

403 {
404  return 2. * maxShear(stress);
405 }

Referenced by getQuantity().

◆ thirdInvariant()

template<typename T >
T RankTwoScalarTools::thirdInvariant ( const RankTwoTensorTempl< T > &  r2tensor)

Definition at line 172 of file RankTwoScalarTools.h.

173 {
174  T val = 0.0;
175  val = r2tensor(0, 0) * r2tensor(1, 1) * r2tensor(2, 2) -
176  r2tensor(0, 0) * r2tensor(1, 2) * r2tensor(2, 1) +
177  r2tensor(0, 1) * r2tensor(1, 2) * r2tensor(2, 0) -
178  r2tensor(0, 1) * r2tensor(1, 0) * r2tensor(2, 2) +
179  r2tensor(0, 2) * r2tensor(1, 0) * r2tensor(2, 1) -
180  r2tensor(0, 2) * r2tensor(1, 1) * r2tensor(2, 0);
181 
182  return val;
183 }

Referenced by getQuantity().

◆ triaxialityStress()

template<typename T >
T RankTwoScalarTools::triaxialityStress ( const RankTwoTensorTempl< T > &  stress)

Definition at line 381 of file RankTwoScalarTools.h.

382 {
383  return hydrostatic(stress) / vonMisesStress(stress);
384 }

Referenced by getQuantity().

◆ volumetricStrain()

template<typename T >
T RankTwoScalarTools::volumetricStrain ( const RankTwoTensorTempl< T > &  strain)

Definition at line 128 of file RankTwoScalarTools.h.

129 {
130  return std::exp(strain(0, 0)) * std::exp(strain(1, 1)) * std::exp(strain(2, 2)) - 1.0;
131 }

Referenced by getQuantity().

◆ vonMisesStress()

template<typename T >
T RankTwoScalarTools::vonMisesStress ( const RankTwoTensorTempl< T > &  stress)

Definition at line 64 of file RankTwoScalarTools.h.

65 {
66  RankTwoTensorTempl<T> dev_stress = stress.deviatoric();
67 
68  return std::sqrt(1.5 * dev_stress.doubleContraction(dev_stress));
69 }

Referenced by RankTwoBasedFailureCriteriaNOSPD::computeFailureCriterionValue(), getQuantity(), and triaxialityStress().

RankTwoScalarTools::firstInvariant
T firstInvariant(const RankTwoTensorTempl< T > &r2tensor)
Definition: RankTwoScalarTools.h:139
RankTwoScalarTools::hydrostatic
T hydrostatic(const RankTwoTensorTempl< T > &r2tensor)
Definition: RankTwoScalarTools.h:88
RankTwoScalarTools::triaxialityStress
T triaxialityStress(const RankTwoTensorTempl< T > &stress)
Definition: RankTwoScalarTools.h:381
RankTwoScalarTools::effectiveStrain
T effectiveStrain(const RankTwoTensorTempl< T > &strain)
Definition: RankTwoScalarTools.h:77
RankTwoScalarTools::maxShear
T maxShear(const RankTwoTensorTempl< T > &stress)
Definition: RankTwoScalarTools.h:392
RankTwoScalarTools::radialStress
T radialStress(const RankTwoTensorTempl< T > &stress, const Point &point1, const Point &point2, const Point &curr_point, Point &direction)
Definition: RankTwoScalarTools.h:339
RankTwoScalarTools::axialStress
T axialStress(const RankTwoTensorTempl< T > &stress, const Point &point1, const Point &point2, Point &direction)
Definition: RankTwoScalarTools.h:258
RankTwoScalarTools::secondInvariant
T secondInvariant(const RankTwoTensorTempl< T > &r2tensor)
Definition: RankTwoScalarTools.h:153
RankTwoScalarTools::normalPositionVector
void normalPositionVector(const Point &point1, const Point &point2, const Point &curr_point, Point &normalPosition)
Definition: RankTwoScalarTools.C:25
RankTwoScalarTools::midPrincipal
T midPrincipal(const RankTwoTensorTempl< T > &r2tensor, Point &direction)
Definition: RankTwoScalarTools.h:230
RankTwoScalarTools::hoopStress
T hoopStress(const RankTwoTensorTempl< T > &stress, const Point &point1, const Point &point2, const Point &curr_point, Point &direction)
Definition: RankTwoScalarTools.h:303
RankTwoScalarTools::directionValueTensor
T directionValueTensor(const RankTwoTensorTempl< T > &r2tensor, Point &direction)
Definition: RankTwoScalarTools.h:366
RankTwoScalarTools::maxPrincipal
T maxPrincipal(const RankTwoTensorTempl< T > &r2tensor, Point &direction)
Definition: RankTwoScalarTools.h:217
RankTwoScalarTools::thirdInvariant
T thirdInvariant(const RankTwoTensorTempl< T > &r2tensor)
Definition: RankTwoScalarTools.h:172
RankTwoScalarTools::calcEigenValuesEigenVectors
T calcEigenValuesEigenVectors(const RankTwoTensorTempl< T > &r2tensor, unsigned int index, Point &eigenvec)
Definition: RankTwoScalarTools.h:195
RankTwoScalarTools::stressIntensity
T stressIntensity(const RankTwoTensorTempl< T > &stress)
Definition: RankTwoScalarTools.h:402
RankTwoScalarTools::volumetricStrain
T volumetricStrain(const RankTwoTensorTempl< T > &strain)
Definition: RankTwoScalarTools.h:128
RankTwoScalarTools::vonMisesStress
T vonMisesStress(const RankTwoTensorTempl< T > &stress)
Definition: RankTwoScalarTools.h:64
RankTwoTensorTempl
Definition: ACGrGrElasticDrivingForce.h:17
RankTwoScalarTools::scalarOptions
MooseEnum scalarOptions()
Definition: RankTwoScalarTools.C:16
RankTwoScalarTools::L2norm
T L2norm(const RankTwoTensorTempl< T > &r2tensor)
Definition: RankTwoScalarTools.h:98
RankTwoScalarTools::minPrincipal
T minPrincipal(const RankTwoTensorTempl< T > &r2tensor, Point &direction)
Definition: RankTwoScalarTools.h:243