OGS
CompareJacobiansJacobianAssembler.cpp
Go to the documentation of this file.
1
12
13#include <sstream>
14
15#include "BaseLib/ConfigTree.h"
18
19namespace
20{
22void dump_py(std::ostream& fh, std::string const& var, double const val)
23{
24 fh << var << " = " << val << '\n';
25}
26
28void dump_py(std::ostream& fh, std::string const& var, std::size_t const val)
29{
30 fh << var << " = " << val << '\n';
31}
32
34template <typename Vec>
35void dump_py_vec(std::ostream& fh, std::string const& var, Vec const& val)
36{
37 fh << var << " = np.array([";
38 for (decltype(val.size()) i = 0; i < val.size(); ++i)
39 {
40 if (i != 0)
41 {
42 if (i % 8 == 0)
43 {
44 // Print at most eight entries on one line,
45 // indent with four spaces.
46 fh << ",\n ";
47 }
48 else
49 {
50 fh << ", ";
51 }
52 }
53 fh << val[i];
54 }
55 fh << "])\n";
56}
57
59void dump_py(std::ostream& fh, std::string const& var,
60 std::vector<double> const& val)
61{
62 dump_py_vec(fh, var, val);
63}
64
66template <typename Derived>
67void dump_py(std::ostream& fh, std::string const& var,
68 Eigen::ArrayBase<Derived> const& val,
69 std::integral_constant<int, 1> /*unused*/)
70{
71 dump_py_vec(fh, var, val);
72}
73
75template <typename Derived, int ColsAtCompileTime>
76void dump_py(std::ostream& fh, std::string const& var,
77 Eigen::ArrayBase<Derived> const& val,
78 std::integral_constant<int, ColsAtCompileTime> /*unused*/)
79{
80 fh << var << " = np.array([\n";
81 for (std::ptrdiff_t r = 0; r < val.rows(); ++r)
82 {
83 if (r != 0)
84 {
85 fh << ",\n";
86 }
87 fh << " [";
88 for (std::ptrdiff_t c = 0; c < val.cols(); ++c)
89 {
90 if (c != 0)
91 {
92 fh << ", ";
93 }
94 fh << val(r, c);
95 }
96 fh << "]";
97 }
98 fh << "])\n";
99}
100
102template <typename Derived>
103void dump_py(std::ostream& fh, std::string const& var,
104 Eigen::ArrayBase<Derived> const& val)
105{
106 dump_py(fh, var, val,
107 std::integral_constant<int, Derived::ColsAtCompileTime>{});
108}
109
111template <typename Derived>
112void dump_py(std::ostream& fh, std::string const& var,
113 Eigen::MatrixBase<Derived> const& val)
114{
115 dump_py(fh, var, val.array());
116}
117
119const std::string msg_fatal =
120 "The local matrices M or K or the local vectors b assembled with the two "
121 "different Jacobian assemblers differ.";
122
123} // anonymous namespace
124
125namespace ProcessLib
126{
128 LocalAssemblerInterface& local_assembler, double const t, double const dt,
129 std::vector<double> const& local_x, std::vector<double> const& local_x_prev,
130 std::vector<double>& local_b_data, std::vector<double>& local_Jac_data)
131{
132 ++_counter;
133
134 auto const num_dof = local_x.size();
135
136 // First assembly -- the one whose results will be added to the global
137 // equation system finally.
138 _asm1->assembleWithJacobian(local_assembler, t, dt, local_x, local_x_prev,
139 local_b_data, local_Jac_data);
140
141 auto const local_b1 = MathLib::toVector(local_b_data);
142
143 std::vector<double> local_b_data2;
144 std::vector<double> local_Jac_data2;
145
146 // Second assembly -- used for checking only.
147 _asm2->assembleWithJacobian(local_assembler, t, dt, local_x, local_x_prev,
148 local_b_data2, local_Jac_data2);
149
150 auto const local_b2 = MathLib::toVector(local_b_data2);
151
152 auto const local_Jac1 = MathLib::toMatrix(local_Jac_data, num_dof, num_dof);
153 auto const local_Jac2 =
154 MathLib::toMatrix(local_Jac_data2, num_dof, num_dof);
155
156 auto const abs_diff = (local_Jac2 - local_Jac1).array().eval();
157 auto const rel_diff =
158 (abs_diff == 0.0)
159 .select(abs_diff,
160 2. * abs_diff /
161 (local_Jac1.cwiseAbs() + local_Jac2.cwiseAbs()).array())
162 .eval();
163
164 auto const abs_diff_mask =
165 (abs_diff.abs() <= _abs_tol)
166 .select(decltype(abs_diff)::Zero(abs_diff.rows(), abs_diff.cols()),
167 decltype(abs_diff)::Ones(abs_diff.rows(), abs_diff.cols()))
168 .eval();
169 auto const rel_diff_mask =
170 (rel_diff.abs() <= _rel_tol)
171 .select(decltype(rel_diff)::Zero(rel_diff.rows(), rel_diff.cols()),
172 decltype(rel_diff)::Ones(rel_diff.rows(), rel_diff.cols()))
173 .eval();
174
175 auto const abs_diff_OK = !abs_diff_mask.any();
176 auto const rel_diff_OK = !rel_diff_mask.any();
177
178 std::ostringstream msg_tolerance;
179 bool tol_exceeded = true;
180 bool fatal_error = false;
181
182 if (abs_diff_OK)
183 {
184 tol_exceeded = false;
185 }
186 else
187 {
188 msg_tolerance << "absolute tolerance of " << _abs_tol << " exceeded";
189 }
190
191 if (rel_diff_OK)
192 {
193 tol_exceeded = false;
194 }
195 else
196 {
197 if (!msg_tolerance.str().empty())
198 {
199 msg_tolerance << " and ";
200 }
201
202 msg_tolerance << "relative tolerance of " << _rel_tol << " exceeded";
203 }
204
205 // basic consistency check if something went terribly wrong
206 auto check_equality = [&fatal_error](auto mat_or_vec1, auto mat_or_vec2)
207 {
208 if (mat_or_vec1.size() == 0 || mat_or_vec2.size() == 0)
209 {
210 return;
211 }
212 if (mat_or_vec1.rows() != mat_or_vec2.rows() ||
213 mat_or_vec1.cols() != mat_or_vec2.cols())
214 {
215 fatal_error = true;
216 }
217 else if (((mat_or_vec1 - mat_or_vec2).array().cwiseAbs() >
218 std::numeric_limits<double>::epsilon())
219 .any())
220 {
221 fatal_error = true;
222 }
223 };
224
225 check_equality(local_b1, local_b2);
226
227 Eigen::VectorXd res1 = Eigen::VectorXd::Zero(num_dof);
228 auto const x = MathLib::toVector(local_x);
229 auto const x_dot = ((x - MathLib::toVector(local_x_prev)) / dt).eval();
230 if (local_b1.size() != 0)
231 {
232 res1.noalias() -= local_b1;
233 }
234
235 Eigen::VectorXd res2 = Eigen::VectorXd::Zero(num_dof);
236 if (local_b2.size() != 0)
237 {
238 res2.noalias() -= local_b2;
239 }
240
241 check_equality(res1, res2);
242
243 if (tol_exceeded)
244 {
245 WARN("Compare Jacobians: {:s}", msg_tolerance.str());
246 }
247
248 bool const output = tol_exceeded || fatal_error;
249
250 if (output)
251 {
252 _log_file << "\n### counter: " << std::to_string(_counter)
253 << " (begin)\n";
254 }
255
256 if (fatal_error)
257 {
258 _log_file << '\n'
259 << "#######################################################\n"
260 << "# FATAL ERROR: " << msg_fatal << '\n'
261 << "# You cannot expect any meaningful insights "
262 "from the Jacobian data printed below!\n"
263 << "# The reason for the mentioned differences "
264 "might be\n"
265 << "# (a) that the assembly routine has side "
266 "effects or\n"
267 << "# (b) that the assembly routines for b "
268 "themselves differ.\n"
269 << "#######################################################\n"
270 << '\n';
271 }
272
273 if (tol_exceeded)
274 {
275 _log_file << "# " << msg_tolerance.str() << "\n\n";
276 }
277
278 if (output)
279 {
280 dump_py(_log_file, "counter", _counter);
281 dump_py(_log_file, "num_dof", num_dof);
282 dump_py(_log_file, "abs_tol", _abs_tol);
283 dump_py(_log_file, "rel_tol", _rel_tol);
284
285 _log_file << '\n';
286
287 dump_py(_log_file, "local_x", local_x);
288 dump_py(_log_file, "local_x_prev", local_x_prev);
289 dump_py(_log_file, "dt", dt);
290
291 _log_file << '\n';
292
293 dump_py(_log_file, "Jacobian_1", local_Jac1);
294 dump_py(_log_file, "Jacobian_2", local_Jac2);
295
296 _log_file << '\n';
297
298 _log_file << "# Jacobian_2 - Jacobian_1\n";
299 dump_py(_log_file, "abs_diff", abs_diff);
300 _log_file << "# Componentwise: 2 * abs_diff / (|Jacobian_1| + "
301 "|Jacobian_2|)\n";
302 dump_py(_log_file, "rel_diff", rel_diff);
303
304 _log_file << '\n';
305
306 _log_file << "# Masks: 0 ... tolerance met, 1 ... tolerance exceeded\n";
307 dump_py(_log_file, "abs_diff_mask", abs_diff_mask);
308 dump_py(_log_file, "rel_diff_mask", rel_diff_mask);
309
310 _log_file << '\n';
311
312 dump_py(_log_file, "b_1", local_b_data);
313 dump_py(_log_file, "b_2", local_b_data2);
314 if (fatal_error && local_b1.size() == local_b2.size())
315 {
316 dump_py(_log_file, "delta_b", local_b2 - local_b1);
317 _log_file << '\n';
318 }
319
320 dump_py(_log_file, "res_1", res1);
321 dump_py(_log_file, "res_2", res2);
322 if (fatal_error)
323 {
324 dump_py(_log_file, "delta_res", res2 - res1);
325 }
326
327 _log_file << '\n';
328
329 _log_file << "### counter: " << std::to_string(_counter) << " (end)\n";
330 }
331
332 if (fatal_error)
333 {
334 _log_file << std::flush;
335 OGS_FATAL("{:s}", msg_fatal);
336 }
337
338 if (tol_exceeded && _fail_on_error)
339 {
340 _log_file << std::flush;
341 OGS_FATAL(
342 "OGS failed, because the two Jacobian implementations returned "
343 "different results.");
344 }
345}
346
347std::unique_ptr<AbstractJacobianAssembler>
349{
350 OGS_FATAL(
351 "CompareJacobiansJacobianAssembler should not be copied. This class "
352 "logs to a file, which would most certainly break after copying "
353 "(concurrent file access) with the current implementation.");
354}
355
356std::unique_ptr<CompareJacobiansJacobianAssembler>
358{
359 // TODO doc script corner case: Parameter could occur at different
360 // locations.
362 config.checkConfigParameter("type", "CompareJacobians");
363
364 auto asm1 =
366 createJacobianAssembler(config.getConfigSubtree("jacobian_assembler"));
367
368 auto asm2 = createJacobianAssembler(
370 config.getConfigSubtree("reference_jacobian_assembler"));
371
373 auto const abs_tol = config.getConfigParameter<double>("abs_tol");
375 auto const rel_tol = config.getConfigParameter<double>("rel_tol");
376
378 auto const fail_on_error = config.getConfigParameter<bool>("fail_on_error");
379
381 auto const log_file = config.getConfigParameter<std::string>("log_file");
382
383 return std::make_unique<CompareJacobiansJacobianAssembler>(
384 std::move(asm1), std::move(asm2), abs_tol, rel_tol, fail_on_error,
385 log_file);
386}
387} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:26
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:40
T getConfigParameter(std::string const &param) const
ConfigTree getConfigSubtree(std::string const &root) const
void checkConfigParameter(std::string const &param, std::string_view const value) const
bool const _fail_on_error
Whether to abort if the tolerances are exceeded.
void assembleWithJacobian(LocalAssemblerInterface &local_assembler, double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) override
std::unique_ptr< AbstractJacobianAssembler > _asm2
std::unique_ptr< AbstractJacobianAssembler > copy() const override
std::unique_ptr< AbstractJacobianAssembler > _asm1
Eigen::Map< const Vector > toVector(std::vector< double > const &data, Eigen::VectorXd::Index size)
Creates an Eigen mapped vector from the given data vector.
Eigen::Map< const Matrix > toMatrix(std::vector< double > const &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
std::unique_ptr< CompareJacobiansJacobianAssembler > createCompareJacobiansJacobianAssembler(BaseLib::ConfigTree const &config)
std::unique_ptr< AbstractJacobianAssembler > createJacobianAssembler(std::optional< BaseLib::ConfigTree > const &config)
void dump_py_vec(std::ostream &fh, std::string const &var, Vec const &val)
Dumps an arbitrary vector as a Python script snippet.
void dump_py(std::ostream &fh, std::string const &var, double const val)
Dumps a double value as a Python script snippet.
const std::string msg_fatal
Will be printed if some consistency error is detected.