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_M_data, std::vector<double>& local_K_data,
131 std::vector<double>& local_b_data, std::vector<double>& local_Jac_data)
132{
133 ++_counter;
134
135 auto const num_dof = local_x.size();
136 auto to_mat = [num_dof](std::vector<double> const& data)
137 {
138 if (data.empty())
139 {
140 return Eigen::Map<Eigen::Matrix<
141 double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor> const>(
142 nullptr, 0, 0);
143 }
144
145 return MathLib::toMatrix(data, num_dof, num_dof);
146 };
147
148 // First assembly -- the one whose results will be added to the global
149 // equation system finally.
150 _asm1->assembleWithJacobian(local_assembler, t, dt, local_x, local_x_prev,
151 local_M_data, local_K_data, local_b_data,
152 local_Jac_data);
153
154 auto const local_M1 = to_mat(local_M_data);
155 auto const local_K1 = to_mat(local_K_data);
156 auto const local_b1 = MathLib::toVector(local_b_data);
157
158 std::vector<double> local_M_data2;
159 std::vector<double> local_K_data2;
160 std::vector<double> local_b_data2;
161 std::vector<double> local_Jac_data2;
162
163 // Second assembly -- used for checking only.
164 _asm2->assembleWithJacobian(local_assembler, t, dt, local_x, local_x_prev,
165 local_M_data2, local_K_data2, local_b_data2,
166 local_Jac_data2);
167
168 auto const local_M2 = to_mat(local_M_data2);
169 auto const local_K2 = to_mat(local_K_data2);
170 auto const local_b2 = MathLib::toVector(local_b_data2);
171
172 auto const local_Jac1 = MathLib::toMatrix(local_Jac_data, num_dof, num_dof);
173 auto const local_Jac2 =
174 MathLib::toMatrix(local_Jac_data2, num_dof, num_dof);
175
176 auto const abs_diff = (local_Jac2 - local_Jac1).array().eval();
177 auto const rel_diff =
178 (abs_diff == 0.0)
179 .select(abs_diff,
180 2. * abs_diff /
181 (local_Jac1.cwiseAbs() + local_Jac2.cwiseAbs()).array())
182 .eval();
183
184 auto const abs_diff_mask =
185 (abs_diff.abs() <= _abs_tol)
186 .select(decltype(abs_diff)::Zero(abs_diff.rows(), abs_diff.cols()),
187 decltype(abs_diff)::Ones(abs_diff.rows(), abs_diff.cols()))
188 .eval();
189 auto const rel_diff_mask =
190 (rel_diff.abs() <= _rel_tol)
191 .select(decltype(rel_diff)::Zero(rel_diff.rows(), rel_diff.cols()),
192 decltype(rel_diff)::Ones(rel_diff.rows(), rel_diff.cols()))
193 .eval();
194
195 auto const abs_diff_OK = !abs_diff_mask.any();
196 auto const rel_diff_OK = !rel_diff_mask.any();
197
198 std::ostringstream msg_tolerance;
199 bool tol_exceeded = true;
200 bool fatal_error = false;
201
202 if (abs_diff_OK)
203 {
204 tol_exceeded = false;
205 }
206 else
207 {
208 msg_tolerance << "absolute tolerance of " << _abs_tol << " exceeded";
209 }
210
211 if (rel_diff_OK)
212 {
213 tol_exceeded = false;
214 }
215 else
216 {
217 if (!msg_tolerance.str().empty())
218 {
219 msg_tolerance << " and ";
220 }
221
222 msg_tolerance << "relative tolerance of " << _rel_tol << " exceeded";
223 }
224
225 // basic consistency check if something went terribly wrong
226 auto check_equality = [&fatal_error](auto mat_or_vec1, auto mat_or_vec2)
227 {
228 if (mat_or_vec1.size() == 0 || mat_or_vec2.size() == 0)
229 {
230 return;
231 }
232 if (mat_or_vec1.rows() != mat_or_vec2.rows() ||
233 mat_or_vec1.cols() != mat_or_vec2.cols())
234 {
235 fatal_error = true;
236 }
237 else if (((mat_or_vec1 - mat_or_vec2).array().cwiseAbs() >
238 std::numeric_limits<double>::epsilon())
239 .any())
240 {
241 fatal_error = true;
242 }
243 };
244
245 check_equality(local_M1, local_M2);
246 check_equality(local_K1, local_K2);
247 check_equality(local_b1, local_b2);
248
249 Eigen::VectorXd res1 = Eigen::VectorXd::Zero(num_dof);
250 auto const x = MathLib::toVector(local_x);
251 auto const x_dot = ((x - MathLib::toVector(local_x_prev)) / dt).eval();
252 if (local_M1.size() != 0)
253 {
254 res1.noalias() += local_M1 * x_dot;
255 }
256 if (local_K1.size() != 0)
257 {
258 res1.noalias() += local_K1 * x;
259 }
260 if (local_b1.size() != 0)
261 {
262 res1.noalias() -= local_b1;
263 }
264
265 Eigen::VectorXd res2 = Eigen::VectorXd::Zero(num_dof);
266 if (local_M2.size() != 0)
267 {
268 res2.noalias() += local_M2 * x_dot;
269 }
270 if (local_K2.size() != 0)
271 {
272 res2.noalias() += local_K2 * x;
273 }
274 if (local_b2.size() != 0)
275 {
276 res2.noalias() -= local_b2;
277 }
278
279 check_equality(res1, res2);
280
281 if (tol_exceeded)
282 {
283 WARN("Compare Jacobians: {:s}", msg_tolerance.str());
284 }
285
286 bool const output = tol_exceeded || fatal_error;
287
288 if (output)
289 {
290 _log_file << "\n### counter: " << std::to_string(_counter)
291 << " (begin)\n";
292 }
293
294 if (fatal_error)
295 {
296 _log_file << '\n'
297 << "#######################################################\n"
298 << "# FATAL ERROR: " << msg_fatal << '\n'
299 << "# You cannot expect any meaningful insights "
300 "from the Jacobian data printed below!\n"
301 << "# The reason for the mentioned differences "
302 "might be\n"
303 << "# (a) that the assembly routine has side "
304 "effects or\n"
305 << "# (b) that the assembly routines for M, K "
306 "and b themselves differ.\n"
307 << "#######################################################\n"
308 << '\n';
309 }
310
311 if (tol_exceeded)
312 {
313 _log_file << "# " << msg_tolerance.str() << "\n\n";
314 }
315
316 if (output)
317 {
318 dump_py(_log_file, "counter", _counter);
319 dump_py(_log_file, "num_dof", num_dof);
320 dump_py(_log_file, "abs_tol", _abs_tol);
321 dump_py(_log_file, "rel_tol", _rel_tol);
322
323 _log_file << '\n';
324
325 dump_py(_log_file, "local_x", local_x);
326 dump_py(_log_file, "local_x_prev", local_x_prev);
327 dump_py(_log_file, "dt", dt);
328
329 _log_file << '\n';
330
331 dump_py(_log_file, "Jacobian_1", local_Jac1);
332 dump_py(_log_file, "Jacobian_2", local_Jac2);
333
334 _log_file << '\n';
335
336 _log_file << "# Jacobian_2 - Jacobian_1\n";
337 dump_py(_log_file, "abs_diff", abs_diff);
338 _log_file << "# Componentwise: 2 * abs_diff / (|Jacobian_1| + "
339 "|Jacobian_2|)\n";
340 dump_py(_log_file, "rel_diff", rel_diff);
341
342 _log_file << '\n';
343
344 _log_file << "# Masks: 0 ... tolerance met, 1 ... tolerance exceeded\n";
345 dump_py(_log_file, "abs_diff_mask", abs_diff_mask);
346 dump_py(_log_file, "rel_diff_mask", rel_diff_mask);
347
348 _log_file << '\n';
349
350 dump_py(_log_file, "M_1", local_M1);
351 dump_py(_log_file, "M_2", local_M2);
352 if (fatal_error && local_M1.size() == local_M2.size())
353 {
354 dump_py(_log_file, "delta_M", local_M2 - local_M1);
355 _log_file << '\n';
356 }
357
358 dump_py(_log_file, "K_1", local_K1);
359 dump_py(_log_file, "K_2", local_K2);
360 if (fatal_error && local_K1.size() == local_K2.size())
361 {
362 dump_py(_log_file, "delta_K", local_K2 - local_K1);
363 _log_file << '\n';
364 }
365
366 dump_py(_log_file, "b_1", local_b_data);
367 dump_py(_log_file, "b_2", local_b_data2);
368 if (fatal_error && local_b1.size() == local_b2.size())
369 {
370 dump_py(_log_file, "delta_b", local_b2 - local_b1);
371 _log_file << '\n';
372 }
373
374 dump_py(_log_file, "res_1", res1);
375 dump_py(_log_file, "res_2", res2);
376 if (fatal_error)
377 {
378 dump_py(_log_file, "delta_res", res2 - res1);
379 }
380
381 _log_file << '\n';
382
383 _log_file << "### counter: " << std::to_string(_counter) << " (end)\n";
384 }
385
386 if (fatal_error)
387 {
388 _log_file << std::flush;
389 OGS_FATAL("{:s}", msg_fatal);
390 }
391
392 if (tol_exceeded && _fail_on_error)
393 {
394 _log_file << std::flush;
395 OGS_FATAL(
396 "OGS failed, because the two Jacobian implementations returned "
397 "different results.");
398 }
399}
400
401std::unique_ptr<AbstractJacobianAssembler>
403{
404 OGS_FATAL(
405 "CompareJacobiansJacobianAssembler should not be copied. This class "
406 "logs to a file, which would most certainly break after copying "
407 "(concurrent file access) with the current implementation.");
408}
409
410std::unique_ptr<CompareJacobiansJacobianAssembler>
412{
413 // TODO doc script corner case: Parameter could occur at different
414 // locations.
416 config.checkConfigParameter("type", "CompareJacobians");
417
418 auto asm1 =
420 createJacobianAssembler(config.getConfigSubtree("jacobian_assembler"));
421
422 auto asm2 = createJacobianAssembler(
424 config.getConfigSubtree("reference_jacobian_assembler"));
425
427 auto const abs_tol = config.getConfigParameter<double>("abs_tol");
429 auto const rel_tol = config.getConfigParameter<double>("rel_tol");
430
432 auto const fail_on_error = config.getConfigParameter<bool>("fail_on_error");
433
435 auto const log_file = config.getConfigParameter<std::string>("log_file");
436
437 return std::make_unique<CompareJacobiansJacobianAssembler>(
438 std::move(asm1), std::move(asm2), abs_tol, rel_tol, fail_on_error,
439 log_file);
440}
441} // 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.
std::unique_ptr< AbstractJacobianAssembler > _asm2
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_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) override
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.