NEML2 1.4.0
Loading...
Searching...
No Matches
AssociativeIsotropicPlasticHardening.cxx
1// Copyright 2023, UChicago Argonne, LLC
2// All Rights Reserved
3// Software Name: NEML2 -- the New Engineering material Model Library, version 2
4// By: Argonne National Laboratory
5// OPEN SOURCE LICENSE (MIT)
6//
7// Permission is hereby granted, free of charge, to any person obtaining a copy
8// of this software and associated documentation files (the "Software"), to deal
9// in the Software without restriction, including without limitation the rights
10// to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
11// copies of the Software, and to permit persons to whom the Software is
12// furnished to do so, subject to the following conditions:
13//
14// The above copyright notice and this permission notice shall be included in
15// all copies or substantial portions of the Software.
16//
17// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
18// IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
19// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
20// AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
21// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
22// OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN
23// THE SOFTWARE.
24
25#include "neml2/models/solid_mechanics/AssociativeIsotropicPlasticHardening.h"
26
27namespace neml2
28{
30
33{
35 options.doc() += " This object calculates the rate of equivalent plastic strain following "
36 "associative flow rule, i.e. \\f$ \\dot{\\varepsilon}_p = - \\dot{\\gamma} "
37 "\\frac{\\partial f}{\\partial k} \\f$, where \\f$ \\dot{\\varepsilon}_p \\f$ "
38 "is the equivalent plastic strain, \\f$ \\dot{\\gamma} \\f$ is the flow rate, "
39 "\\f$ f \\f$ is the yield function, and \\f$ k \\f$ is the isotropic hardening.";
40
41 options.set<VariableName>("isotropic_hardening_direction") =
42 VariableName("state", "internal", "Nk");
43 options.set("isotropic_hardening_direction").doc() =
44 "Direction of associative isotropic hardening which can be calculated using Normality.";
45
46 options.set<VariableName>("equivalent_plastic_strain_rate") =
47 VariableName("state", "internal", "ep_rate");
48 options.set("equivalent_plastic_strain_rate").doc() = "Rate of equivalent plastic strain";
49
50 return options;
51}
52
54 const OptionSet & options)
55 : FlowRule(options),
56 _Nk(declare_input_variable<Scalar>("isotropic_hardening_direction")),
57 _ep_dot(declare_output_variable<Scalar>("equivalent_plastic_strain_rate"))
58{
59}
60
61void
63{
64 // For associative flow,
65 // ep_dot = - gamma_dot * Nk
66 // Nk = df/dk
67
68 if (out)
70
71 if (dout_din)
72 {
75 }
76
77 if (d2out_din2)
78 {
79 // I don't know when this will be useful, but since it's easy...
83 }
84}
85} // namespace neml2
Definition AssociativeIsotropicPlasticHardening.h:32
Variable< Scalar > & _ep_dot
Rate of equivalent plastic strain.
Definition AssociativeIsotropicPlasticHardening.h:51
const Variable< Scalar > & _Nk
Isotropic hardening direction.
Definition AssociativeIsotropicPlasticHardening.h:48
AssociativeIsotropicPlasticHardening(const OptionSet &options)
Definition AssociativeIsotropicPlasticHardening.cxx:53
static OptionSet expected_options()
Definition AssociativeIsotropicPlasticHardening.cxx:32
void set_value(bool out, bool dout_din, bool d2out_din2) override
The map between input -> output, and optionally its derivatives.
Definition AssociativeIsotropicPlasticHardening.cxx:62
The wrapper (decorator) for cross-referencing unresolved values at parse time.
Definition CrossRef.h:52
Definition FlowRule.h:32
const Variable< Scalar > & _gamma_dot
Definition FlowRule.h:39
static OptionSet expected_options()
Definition FlowRule.cxx:30
The accessor containing all the information needed to access an item in a LabeledAxis.
Definition LabeledAxisAccessor.h:44
const torch::TensorOptions & options() const
This model's tensor options.
Definition Model.h:116
A custom map-like data structure. The keys are strings, and the values can be nonhomogeneously typed.
Definition OptionSet.h:59
The (logical) scalar.
Definition Scalar.h:38
static Scalar identity_map(const torch::TensorOptions &options=default_tensor_options())
The derivative of a Scalar with respect to itself.
Definition Scalar.cxx:35
Derivative d(const VariableBase &x)
Create a wrapper representing the derivative dy/dx.
Definition Variable.cxx:102
Definition CrossRef.cxx:32
LabeledAxisAccessor VariableName
Definition Variable.h:35