Go to the documentation of this file.
11 #include "RankFourTensor.h"
12 #include "libmesh/utility.h"
22 params.addRangeCheckedParam<
unsigned>(
"max_iterations",
25 "Maximum iterations for custom MeanCapTC return map");
26 params.addParam<
bool>(
27 "use_custom_returnMap",
true,
"Whether to use the custom MeanCapTC returnMap algorithm.");
28 params.addParam<
bool>(
"use_custom_cto",
30 "Whether to use the custom consistent tangent operator computations.");
31 params.addRequiredParam<UserObjectName>(
"tensile_strength",
32 "A TensorMechanicsHardening UserObject that defines "
33 "hardening of the mean-cap tensile strength (it will "
34 "typically be positive). Yield function = trace(stress) "
35 "- tensile_strength for trace(stress)>tensile_strength.");
36 params.addRequiredParam<UserObjectName>(
37 "compressive_strength",
38 "A TensorMechanicsHardening UserObject that defines hardening of the "
39 "mean-cap compressive strength. This should always be less than "
40 "tensile_strength (it will typically be negative). Yield function = "
41 "- (trace(stress) - compressive_strength) for "
42 "trace(stress)<compressive_strength.");
43 params.addClassDescription(
44 "Associative mean-cap tensile and compressive plasticity with hardening/softening");
51 _max_iters(getParam<unsigned>(
"max_iterations")),
52 _use_custom_returnMap(getParam<bool>(
"use_custom_returnMap")),
53 _use_custom_cto(getParam<bool>(
"use_custom_cto")),
60 mooseError(
"MeanCapTC: tensile strength (which is usually positive) must not be less than "
61 "compressive strength (which is usually negative)");
67 const Real tr = stress.trace();
79 return (c_str - t_str) / libMesh::pi * std::sin(libMesh::pi * (tr - c_str) / (t_str - c_str));
93 const Real tr = stress.trace();
104 return (dc - dt) / libMesh::pi * std::sin(libMesh::pi * (tr - c_str) / (t_str - c_str)) +
105 1.0 / (t_str - c_str) * std::cos(libMesh::pi * (tr - c_str) / (t_str - c_str)) *
106 ((tr - c_str) * dt - (tr - t_str) * dc);
112 const Real tr = stress.trace();
115 return stress.dtrace();
119 return -stress.dtrace();
121 return -std::cos(libMesh::pi * (tr - c_str) / (t_str - c_str)) * stress.dtrace();
134 const Real tr = stress.trace();
143 return libMesh::pi / (t_str - c_str) * std::sin(libMesh::pi * (tr - c_str) / (t_str - c_str)) *
144 stress.dtrace().outerProduct(stress.dtrace());
151 const Real tr = stress.trace();
162 return std::sin(libMesh::pi * (tr - c_str) / (t_str - c_str)) * stress.dtrace() * libMesh::pi /
163 Utility::pow<2>(t_str - c_str) * ((tr - t_str) * dc - (tr - c_str) * dt);
170 const Real tr = stress.trace();
182 return std::cos(libMesh::pi * (tr - c_str) /
190 const Real tr = stress.trace();
199 return -std::sin(libMesh::pi * (tr - c_str) / (t_str - c_str)) * libMesh::pi / (t_str - c_str) *
207 const Real tr = stress.trace();
218 return -std::sin(libMesh::pi * (tr - c_str) / (t_str - c_str)) * libMesh::pi /
219 Utility::pow<2>(t_str - c_str) * ((tr - t_str) * dc - (tr - c_str) * dt);
251 std::vector<bool> & act,
254 act.assign(1,
false);
258 returned_stress = stress;
262 const Real tr = stress.trace();
278 for (
unsigned i = 0; i < 3; ++i)
279 for (
unsigned j = 0; j < 3; ++j)
280 for (
unsigned k = 0; k < 3; ++k)
281 n(i, j) += dirn * Eijkl(i, j, k, k);
287 Real gamma = (stress.trace() - str) / n.trace();
289 for (
unsigned i = 0; i < 3; ++i)
290 for (
unsigned j = 0; j < 3; ++j)
291 returned_stress(i, j) = stress(i, j) - gamma * n(i, j);
300 Real ep_plastic_tolerance,
302 Real & returned_intnl,
303 std::vector<Real> & dpm,
305 std::vector<Real> & yf,
306 bool & trial_stress_inadmissible)
const
312 ep_plastic_tolerance,
318 trial_stress_inadmissible);
329 trial_stress_inadmissible =
false;
333 trial_stress_inadmissible =
true;
350 const bool tensile_failure = (trial_stress.trace() >=
tensile_strength(intnl_old));
351 const Real dirn = (tensile_failure ? 1.0 : -1.0);
354 for (
unsigned i = 0; i < 3; ++i)
355 for (
unsigned j = 0; j < 3; ++j)
356 for (
unsigned k = 0; k < 3; ++k)
357 n(i, j) += dirn * E_ijkl(i, j, k, k);
358 const Real n_trace = n.trace();
366 Real trial_trace = trial_stress.trace();
370 unsigned int iter = 0;
375 residual = trial_trace -
tensile_strength(intnl_old + dpm[0]) - dpm[0] * n_trace;
383 dpm[0] += -residual / jac;
391 returned_intnl = intnl_old + dirn * dpm[0];
392 returned_stress = trial_stress - dpm[0] * n;
393 delta_dp = dpm[0] * dirn * returned_stress.dtrace();
405 const std::vector<Real> & cumulative_pm)
const
409 trial_stress, intnl_old, stress, intnl, E_ijkl, cumulative_pm);
425 for (
unsigned int i = 0; i < 3; ++i)
426 for (
unsigned int j = 0; j < 3; ++j)
427 for (
unsigned int k = 0; k < 3; ++k)
428 elas(i, j) += E_ijkl(i, j, k, k);
430 const Real hw = -df_dq + alpha * elas.trace();
432 return E_ijkl - alpha / hw * elas.outerProduct(elas);
virtual bool returnMap(const RankTwoTensor &trial_stress, Real intnl_old, const RankFourTensor &E_ijkl, Real ep_plastic_tolerance, RankTwoTensor &returned_stress, Real &returned_intnl, std::vector< Real > &dpm, RankTwoTensor &delta_dp, std::vector< Real > &yf, bool &trial_stress_inadmissible) const override
Performs a custom return-map.
TensorMechanicsPlasticMeanCapTC(const InputParameters ¶meters)
static InputParameters validParams()
virtual RankFourTensor consistentTangentOperator(const RankTwoTensor &trial_stress, Real intnl_old, const RankTwoTensor &stress, Real intnl, const RankFourTensor &E_ijkl, const std::vector< Real > &cumulative_pm) const override
Calculates a custom consistent tangent operator.
const bool _use_custom_returnMap
Whether to use the custom return-map algorithm.
Rate-independent associative mean-cap tensile AND compressive failure with hardening/softening of the...
const TensorMechanicsHardeningModel & _strength
the tensile strength
const TensorMechanicsHardeningModel & _c_strength
the compressive strength
RankTwoTensor dyieldFunction_dstress(const RankTwoTensor &stress, Real intnl) const override
The derivative of yield function with respect to stress.
Real hardPotential(const RankTwoTensor &stress, Real intnl) const override
The hardening potential.
registerMooseObject("TensorMechanicsApp", TensorMechanicsPlasticMeanCapTC)
const unsigned _max_iters
max iters for custom return map loop
const Real _f_tol
Tolerance on yield function.
RankTwoTensor df_dsig(const RankTwoTensor &stress, Real intnl) const
Derivative of the yield function with respect to stress.
virtual bool useCustomCTO() const override
Returns false. You will want to override this in your derived class if you write a custom consistent ...
virtual RankTwoTensor dhardPotential_dstress(const RankTwoTensor &stress, Real intnl) const override
The derivative of the hardening potential with respect to stress.
virtual Real compressive_strength(const Real internal_param) const
compressive strength as a function of residual value, rate, and internal_param
virtual Real derivative(Real intnl) const
const bool _use_custom_cto
Whether to use the custom consistent tangent operator algorithm.
virtual bool returnMap(const RankTwoTensor &trial_stress, Real intnl_old, const RankFourTensor &E_ijkl, Real ep_plastic_tolerance, RankTwoTensor &returned_stress, Real &returned_intnl, std::vector< Real > &dpm, RankTwoTensor &delta_dp, std::vector< Real > &yf, bool &trial_stress_inadmissible) const
Performs a custom return-map.
virtual Real value(Real intnl) const
RankTwoTensorTempl< Real > RankTwoTensor
virtual std::string modelName() const override
virtual Real dhardPotential_dintnl(const RankTwoTensor &stress, Real intnl) const override
The derivative of the hardening potential with respect to the internal parameter.
Real yieldFunction(const RankTwoTensor &stress, Real intnl) const override
The following functions are what you should override when building single-plasticity models.
virtual RankFourTensor consistentTangentOperator(const RankTwoTensor &trial_stress, Real intnl_old, const RankTwoTensor &stress, Real intnl, const RankFourTensor &E_ijkl, const std::vector< Real > &cumulative_pm) const
Calculates a custom consistent tangent operator.
RankTwoTensor dflowPotential_dintnl(const RankTwoTensor &stress, Real intnl) const override
The derivative of the flow potential with respect to the internal parameter.
defineLegacyParams(TensorMechanicsPlasticMeanCapTC)
static InputParameters validParams()
virtual Real dtensile_strength(const Real internal_param) const
d(tensile strength)/d(internal_param) as a function of residual value, rate, and internal_param
Plastic Model base class The virtual functions written below must be over-ridden in derived classes t...
RankTwoTensor flowPotential(const RankTwoTensor &stress, Real intnl) const override
The flow potential.
virtual void activeConstraints(const std::vector< Real > &f, const RankTwoTensor &stress, Real intnl, const RankFourTensor &Eijkl, std::vector< bool > &act, RankTwoTensor &returned_stress) const override
The active yield surfaces, given a vector of yield functions.
RankFourTensor dflowPotential_dstress(const RankTwoTensor &stress, Real intnl) const override
The derivative of the flow potential with respect to stress.
Real dyieldFunction_dintnl(const RankTwoTensor &stress, Real intnl) const override
The derivative of yield function with respect to the internal parameter.
virtual bool useCustomReturnMap() const override
Returns false. You will want to override this in your derived class if you write a custom returnMap f...
virtual Real tensile_strength(const Real internal_param) const
tensile strength as a function of residual value, rate, and internal_param
RankFourTensorTempl< Real > RankFourTensor
virtual Real dcompressive_strength(const Real internal_param) const
d(compressive strength)/d(internal_param) as a function of residual value, rate, and internal_param
Hardening Model base class.