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