Line data Source code
1 : // Copyright 2024, 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/chemical_reactions/EffectiveVolume.h"
26 :
27 : #include "neml2/misc/assertions.h"
28 :
29 : namespace neml2
30 : {
31 : register_NEML2_object(EffectiveVolume);
32 : OptionSet
33 2 : EffectiveVolume::expected_options()
34 : {
35 2 : OptionSet options = Model::expected_options();
36 2 : options.doc() =
37 : "Calculate the total volume of a control mass. The volume has the form of \\f$ V = "
38 : "\\dfrac{M}{1-\\phi_{o}} \\sum_i \\frac{\\omega_i}{\\rho_i} \\f$, where \\f$ \\omega_i "
39 : "\\f$ and \\f$ \\rho_i \\f$ are respectively the mass fraction and the density of each "
40 : "component; \\f$ \\phi_{o} \\f$ is the volume fraction accounting for leakage from the "
41 2 : "control mass; \\f$ M \\f$ is the reference mass of the composite.";
42 :
43 4 : options.set_parameter<TensorName<Scalar>>("reference_mass");
44 4 : options.set("reference_mass").doc() = "Reference mass of the composite";
45 :
46 4 : options.set_input("open_volume_fraction");
47 4 : options.set("open_volume_fraction").doc() = "Open volume fraction accounting for leakage";
48 :
49 4 : options.set<std::vector<VariableName>>("mass_fractions");
50 4 : options.set("mass_fractions").doc() = "Mass fractions of the components in the composite";
51 :
52 4 : options.set_parameter<std::vector<TensorName<Scalar>>>("densities");
53 2 : options.set("densities").doc() = "Densities of the components in the composite";
54 :
55 6 : options.set_output("composite_volume") = VariableName(STATE, "V");
56 2 : options.set("composite_volume").doc() = "Volume of the composite";
57 :
58 2 : return options;
59 0 : }
60 :
61 1 : EffectiveVolume::EffectiveVolume(const OptionSet & options)
62 : : Model(options),
63 4 : _M(declare_parameter<Scalar>("M", "reference_mass")),
64 1 : _phio(options.get("open_volume_fraction").user_specified()
65 1 : ? &declare_input_variable<Scalar>("open_volume_fraction")
66 : : nullptr),
67 2 : _V(declare_output_variable<Scalar>("composite_volume"))
68 : {
69 6 : for (const auto & v : options.get<std::vector<VariableName>>("mass_fractions"))
70 5 : _ws.push_back(&declare_input_variable<Scalar>(v));
71 :
72 : // densities
73 1 : const auto rho_refs = options.get<std::vector<TensorName<Scalar>>>("densities");
74 :
75 : // sizes must match
76 1 : neml_assert(rho_refs.size() == _ws.size(),
77 : "Number of mass fractions (",
78 1 : _ws.size(),
79 : ") does not match number of densities (",
80 1 : rho_refs.size(),
81 : ").");
82 :
83 : // declare densities as parameters
84 1 : _rhos.resize(_ws.size());
85 5 : for (std::size_t i = 0; i < _ws.size(); i++)
86 4 : _rhos[i] = &declare_parameter<Scalar>(
87 8 : "rho_" + std::to_string(i), rho_refs[i], /*allow_nonlinear=*/false);
88 1 : }
89 :
90 : void
91 2 : EffectiveVolume::set_value(bool out, bool dout_din, bool /*d2out_din2*/)
92 : {
93 2 : auto sum = Scalar::zeros_like(*_ws[0]);
94 :
95 10 : for (std::size_t i = 0; i < _ws.size(); i++)
96 8 : sum = sum + *_ws[i] / *_rhos[i];
97 :
98 2 : auto coef = _phio ? _M / (1 - *_phio) : _M;
99 :
100 2 : if (out)
101 1 : _V = coef * sum;
102 :
103 2 : if (dout_din)
104 : {
105 5 : for (std::size_t i = 0; i < _ws.size(); i++)
106 4 : if (_ws[i]->is_dependent())
107 4 : _V.d(*_ws[i]) = coef / *_rhos[i];
108 :
109 1 : if (_phio && _phio->is_dependent())
110 1 : _V.d(*_phio) = _M * sum / ((1 - *_phio) * (1 - *_phio));
111 : }
112 2 : }
113 : }
|