Line data Source code
1 : // Copyright Amazon.com, Inc. or its affiliates. All Rights Reserved.
2 : // SPDX-License-Identifier: Apache-2.0
3 :
4 : #include "surfaceimpedanceoperator.hpp"
5 :
6 : #include <set>
7 : #include "models/materialoperator.hpp"
8 : #include "utils/communication.hpp"
9 : #include "utils/geodata.hpp"
10 : #include "utils/iodata.hpp"
11 : #include "utils/prettyprint.hpp"
12 :
13 : namespace palace
14 : {
15 :
16 17 : SurfaceImpedanceOperator::SurfaceImpedanceOperator(const IoData &iodata,
17 : const MaterialOperator &mat_op,
18 17 : const mfem::ParMesh &mesh)
19 17 : : mat_op(mat_op)
20 : {
21 : // Print out BC info for all impedance boundary attributes.
22 17 : SetUpBoundaryProperties(iodata, mesh);
23 17 : PrintBoundaryInfo(iodata, mesh);
24 17 : }
25 :
26 17 : void SurfaceImpedanceOperator::SetUpBoundaryProperties(const IoData &iodata,
27 : const mfem::ParMesh &mesh)
28 : {
29 : // Check that impedance boundary attributes have been specified correctly.
30 17 : int bdr_attr_max = mesh.bdr_attributes.Size() ? mesh.bdr_attributes.Max() : 0;
31 : mfem::Array<int> bdr_attr_marker;
32 17 : if (!iodata.boundaries.impedance.empty())
33 : {
34 : mfem::Array<int> impedance_marker(bdr_attr_max);
35 : bdr_attr_marker.SetSize(bdr_attr_max);
36 : bdr_attr_marker = 0;
37 : impedance_marker = 0;
38 0 : for (auto attr : mesh.bdr_attributes)
39 : {
40 0 : bdr_attr_marker[attr - 1] = 1;
41 : }
42 : std::set<int> bdr_warn_list;
43 0 : for (const auto &data : iodata.boundaries.impedance)
44 : {
45 0 : for (auto attr : data.attributes)
46 : {
47 0 : MFEM_VERIFY(
48 : !impedance_marker[attr - 1],
49 : "Multiple definitions of impedance boundary properties for boundary attribute "
50 : << attr << "!");
51 0 : impedance_marker[attr - 1] = 1;
52 : // MFEM_VERIFY(attr > 0 && attr <= bdr_attr_max,
53 : // "Impedance boundary attribute tags must be non-negative and
54 : // correspond " "to attributes in the mesh!");
55 : // MFEM_VERIFY(bdr_attr_marker[attr - 1],
56 : // "Unknown impedance boundary attribute " << attr << "!");
57 0 : if (attr <= 0 || attr > bdr_attr_max || !bdr_attr_marker[attr - 1])
58 : {
59 : bdr_warn_list.insert(attr);
60 : }
61 : }
62 : }
63 0 : if (!bdr_warn_list.empty())
64 : {
65 0 : Mpi::Print("\n");
66 0 : Mpi::Warning("Unknown impedance boundary attributes!\nSolver will just ignore them!");
67 0 : utils::PrettyPrint(bdr_warn_list, "Boundary attribute list:");
68 0 : Mpi::Print("\n");
69 : }
70 : }
71 :
72 : // Impedance boundaries are defined using the user provided impedance per square.
73 17 : boundaries.reserve(iodata.boundaries.impedance.size());
74 17 : for (const auto &data : iodata.boundaries.impedance)
75 : {
76 0 : MFEM_VERIFY(std::abs(data.Rs) + std::abs(data.Ls) + std::abs(data.Cs) > 0.0,
77 : "Impedance boundary has no Rs, Ls, or Cs defined!");
78 0 : auto &bdr = boundaries.emplace_back();
79 0 : bdr.Rs = data.Rs;
80 0 : bdr.Ls = data.Ls;
81 0 : bdr.Cs = data.Cs;
82 0 : bdr.attr_list.Reserve(static_cast<int>(data.attributes.size()));
83 0 : for (auto attr : data.attributes)
84 : {
85 0 : if (attr <= 0 || attr > bdr_attr_max || !bdr_attr_marker[attr - 1])
86 : {
87 0 : continue; // Can just ignore if wrong
88 : }
89 0 : bdr.attr_list.Append(attr);
90 : }
91 : }
92 17 : }
93 :
94 17 : void SurfaceImpedanceOperator::PrintBoundaryInfo(const IoData &iodata,
95 : const mfem::ParMesh &mesh)
96 : {
97 17 : if (boundaries.empty())
98 : {
99 17 : return;
100 : }
101 :
102 : fmt::memory_buffer buf{}; // Output buffer & buffer append lambda for cleaner code
103 0 : auto to = [&buf](auto fmt, auto &&...args)
104 0 : { fmt::format_to(std::back_inserter(buf), fmt, std::forward<decltype(args)>(args)...); };
105 :
106 : using VT = Units::ValueType;
107 :
108 0 : to("\nConfiguring Robin impedance BC at attributes:\n");
109 0 : for (const auto &bdr : boundaries)
110 : {
111 0 : for (auto attr : bdr.attr_list)
112 : {
113 0 : to(" {:d}:", attr);
114 0 : if (std::abs(bdr.Rs) > 0.0)
115 : {
116 0 : to(" Rs = {:.3e} Ω/sq,", iodata.units.Dimensionalize<VT::IMPEDANCE>(bdr.Rs));
117 : }
118 0 : if (std::abs(bdr.Ls) > 0.0)
119 : {
120 0 : to(" Ls = {:.3e} H/sq,", iodata.units.Dimensionalize<VT::INDUCTANCE>(bdr.Ls));
121 : }
122 0 : if (std::abs(bdr.Cs) > 0.0)
123 : {
124 0 : to(" Cs = {:.3e} F/sq,", iodata.units.Dimensionalize<VT::CAPACITANCE>(bdr.Cs));
125 : }
126 0 : to(" n = ({:+.1f})\n", fmt::join(mesh::GetSurfaceNormal(mesh, attr), ","));
127 : }
128 : }
129 0 : Mpi::Print("{}", fmt::to_string(buf));
130 : }
131 :
132 17 : mfem::Array<int> SurfaceImpedanceOperator::GetAttrList() const
133 : {
134 : mfem::Array<int> attr_list;
135 17 : for (const auto &bdr : boundaries)
136 : {
137 : attr_list.Append(bdr.attr_list);
138 : }
139 17 : return attr_list;
140 : }
141 :
142 17 : mfem::Array<int> SurfaceImpedanceOperator::GetRsAttrList() const
143 : {
144 : mfem::Array<int> attr_list;
145 17 : for (const auto &bdr : boundaries)
146 : {
147 0 : if (std::abs(bdr.Rs) > 0.0)
148 : {
149 : attr_list.Append(bdr.attr_list);
150 : }
151 : }
152 17 : return attr_list;
153 : }
154 :
155 17 : mfem::Array<int> SurfaceImpedanceOperator::GetLsAttrList() const
156 : {
157 : mfem::Array<int> attr_list;
158 17 : for (const auto &bdr : boundaries)
159 : {
160 0 : if (std::abs(bdr.Ls) > 0.0)
161 : {
162 : attr_list.Append(bdr.attr_list);
163 : }
164 : }
165 17 : return attr_list;
166 : }
167 :
168 0 : mfem::Array<int> SurfaceImpedanceOperator::GetCsAttrList() const
169 : {
170 : mfem::Array<int> attr_list;
171 0 : for (const auto &bdr : boundaries)
172 : {
173 0 : if (std::abs(bdr.Cs) > 0.0)
174 : {
175 : attr_list.Append(bdr.attr_list);
176 : }
177 : }
178 0 : return attr_list;
179 : }
180 :
181 0 : void SurfaceImpedanceOperator::AddStiffnessBdrCoefficients(double coeff,
182 : MaterialPropertyCoefficient &fb)
183 : {
184 : // Lumped inductor boundaries.
185 0 : for (const auto &bdr : boundaries)
186 : {
187 0 : if (std::abs(bdr.Ls) > 0.0)
188 : {
189 0 : fb.AddMaterialProperty(mat_op.GetCeedBdrAttributes(bdr.attr_list), coeff / bdr.Ls);
190 : }
191 : }
192 0 : }
193 :
194 0 : void SurfaceImpedanceOperator::AddDampingBdrCoefficients(double coeff,
195 : MaterialPropertyCoefficient &fb)
196 : {
197 : // Lumped resistor boundaries.
198 0 : for (const auto &bdr : boundaries)
199 : {
200 0 : if (std::abs(bdr.Rs) > 0.0)
201 : {
202 0 : fb.AddMaterialProperty(mat_op.GetCeedBdrAttributes(bdr.attr_list), coeff / bdr.Rs);
203 : }
204 : }
205 0 : }
206 :
207 0 : void SurfaceImpedanceOperator::AddMassBdrCoefficients(double coeff,
208 : MaterialPropertyCoefficient &fb)
209 : {
210 : // Lumped capacitor boundaries.
211 0 : for (const auto &bdr : boundaries)
212 : {
213 0 : if (std::abs(bdr.Cs) > 0.0)
214 : {
215 0 : fb.AddMaterialProperty(mat_op.GetCeedBdrAttributes(bdr.attr_list), coeff * bdr.Cs);
216 : }
217 : }
218 0 : }
219 :
220 : } // namespace palace
|