Superimposed viscoplastic flow rule¶
Overview¶
This model sums the contribution of multiple individual viscoplastic hardening rules according to the relations:

where
here is all of
,
, and
and
represents concatenation. The internal variables is then
the set of all variables maintained by each individual flow rule (and these
variables cannot overlap).
Parameters¶
Parameter |
Object type |
Description |
Default |
|---|---|---|---|
|
|
List of individual flow rules |
No |
Class description¶
-
class SuperimposedViscoPlasticFlowRule : public neml::ViscoPlasticFlowRule¶
Superimpose multiple flow rules into a composite model.
Public Functions
-
SuperimposedViscoPlasticFlowRule(ParameterSet ¶ms)¶
-
size_t nmodels() const¶
Number of individual models being summed.
-
virtual void y(const double *const s, const double *const alpha, double T, double &yv) const¶
Scalar flow rate.
-
virtual void dy_ds(const double *const s, const double *const alpha, double T, double *const dyv) const¶
Derivative of scalar flow wrt stress.
-
virtual void dy_da(const double *const s, const double *const alpha, double T, double *const dyv) const¶
Derivative of scalar flow wrt history.
-
virtual void g(const double *const s, const double *const alpha, double T, double *const gv) const¶
Contribution towards the flow proportional to the scalar inelastic strain rate
-
virtual void dg_ds(const double *const s, const double *const alpha, double T, double *const dgv) const¶
Derivative of g wrt stress.
-
virtual void dg_da(const double *const s, const double *const alpha, double T, double *const dgv) const¶
Derivative of g wrt history.
-
virtual void g_time(const double *const s, const double *const alpha, double T, double *const gv) const¶
Contribution towards the flow proportional directly to time.
-
virtual void dg_ds_time(const double *const s, const double *const alpha, double T, double *const dgv) const¶
Derivative of g_time wrt stress.
-
virtual void dg_da_time(const double *const s, const double *const alpha, double T, double *const dgv) const¶
Derivative of g_time wrt history.
-
virtual void g_temp(const double *const s, const double *const alpha, double T, double *const gv) const¶
Contribution towards the flow proportional to the temperature rate.
-
virtual void dg_ds_temp(const double *const s, const double *const alpha, double T, double *const dgv) const¶
Derivative of g_temp wrt stress.
-
virtual void dg_da_temp(const double *const s, const double *const alpha, double T, double *const dgv) const¶
Derivative of g_temp wrt history.
-
virtual void h(const double *const s, const double *const alpha, double T, double *const hv) const¶
Hardening rate proportional to the scalar inelastic strain rate.
-
virtual void dh_ds(const double *const s, const double *const alpha, double T, double *const dhv) const¶
Derivative of h wrt stress.
-
virtual void dh_da(const double *const s, const double *const alpha, double T, double *const dhv) const¶
Derivative of h wrt history.
-
virtual void h_time(const double *const s, const double *const alpha, double T, double *const hv) const¶
Hardening rate proportional directly to time.
-
virtual void dh_ds_time(const double *const s, const double *const alpha, double T, double *const dhv) const¶
Derivative of h_time wrt stress.
-
virtual void dh_da_time(const double *const s, const double *const alpha, double T, double *const dhv) const¶
Derivative of h_time wrt history.
-
virtual void h_temp(const double *const s, const double *const alpha, double T, double *const hv) const¶
Hardening rate proportional to the temperature rate.
-
virtual void dh_ds_temp(const double *const s, const double *const alpha, double T, double *const dhv) const¶
Derivative of h_temp wrt. stress.
-
virtual void dh_da_temp(const double *const s, const double *const alpha, double T, double *const dhv) const¶
Derivative of h_temp wrt history.
Public Static Functions
-
static std::string type()¶
String type for the object system.
-
static std::unique_ptr<NEMLObject> initialize(ParameterSet ¶ms)¶
Default parameters.
-
static ParameterSet parameters()¶
Initialize from parameters.
-
SuperimposedViscoPlasticFlowRule(ParameterSet ¶ms)¶