OGS
CompareJacobiansJacobianAssembler.cpp
Go to the documentation of this file.
1 
12 
13 #include <sstream>
14 
17 
18 namespace
19 {
21 void dump_py(std::ostream& fh, std::string const& var, double const val)
22 {
23  fh << var << " = " << val << '\n';
24 }
25 
27 void dump_py(std::ostream& fh, std::string const& var, std::size_t const val)
28 {
29  fh << var << " = " << val << '\n';
30 }
31 
33 template <typename Vec>
34 void dump_py_vec(std::ostream& fh, std::string const& var, Vec const& val)
35 {
36  fh << var << " = np.array([";
37  for (decltype(val.size()) i = 0; i < val.size(); ++i)
38  {
39  if (i != 0)
40  {
41  if (i % 8 == 0)
42  {
43  // Print at most eight entries on one line,
44  // indent with four spaces.
45  fh << ",\n ";
46  }
47  else
48  {
49  fh << ", ";
50  }
51  }
52  fh << val[i];
53  }
54  fh << "])\n";
55 }
56 
58 void dump_py(std::ostream& fh, std::string const& var,
59  std::vector<double> const& val)
60 {
61  dump_py_vec(fh, var, val);
62 }
63 
65 template <typename Derived>
66 void dump_py(std::ostream& fh, std::string const& var,
67  Eigen::ArrayBase<Derived> const& val,
68  std::integral_constant<int, 1> /*unused*/)
69 {
70  dump_py_vec(fh, var, val);
71 }
72 
74 template <typename Derived, int ColsAtCompileTime>
75 void dump_py(std::ostream& fh, std::string const& var,
76  Eigen::ArrayBase<Derived> const& val,
77  std::integral_constant<int, ColsAtCompileTime> /*unused*/)
78 {
79  fh << var << " = np.array([\n";
80  for (std::ptrdiff_t r = 0; r < val.rows(); ++r)
81  {
82  if (r != 0)
83  {
84  fh << ",\n";
85  }
86  fh << " [";
87  for (std::ptrdiff_t c = 0; c < val.cols(); ++c)
88  {
89  if (c != 0)
90  {
91  fh << ", ";
92  }
93  fh << val(r, c);
94  }
95  fh << "]";
96  }
97  fh << "])\n";
98 }
99 
101 template <typename Derived>
102 void dump_py(std::ostream& fh, std::string const& var,
103  Eigen::ArrayBase<Derived> const& val)
104 {
105  dump_py(fh, var, val,
106  std::integral_constant<int, Derived::ColsAtCompileTime>{});
107 }
108 
110 template <typename Derived>
111 void dump_py(std::ostream& fh, std::string const& var,
112  Eigen::MatrixBase<Derived> const& val)
113 {
114  dump_py(fh, var, val.array());
115 }
116 
118 const std::string msg_fatal =
119  "The local matrices M or K or the local vectors b assembled with the two "
120  "different Jacobian assemblers differ.";
121 
122 } // anonymous namespace
123 
124 namespace ProcessLib
125 {
127  LocalAssemblerInterface& local_assembler, double const t, double const dt,
128  std::vector<double> const& local_x, std::vector<double> const& local_xdot,
129  const double dxdot_dx, const double dx_dx,
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_xdot,
151  dxdot_dx, dx_dx, local_M_data, local_K_data,
152  local_b_data, 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_xdot,
165  dxdot_dx, dx_dx, local_M_data2, local_K_data2,
166  local_b_data2, 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  if (local_M1.size() != 0)
251  {
252  res1.noalias() += local_M1 * MathLib::toVector(local_xdot);
253  }
254  if (local_K1.size() != 0)
255  {
256  res1.noalias() += local_K1 * MathLib::toVector(local_x);
257  }
258  if (local_b1.size() != 0)
259  {
260  res1.noalias() -= local_b1;
261  }
262 
263  Eigen::VectorXd res2 = Eigen::VectorXd::Zero(num_dof);
264  if (local_M2.size() != 0)
265  {
266  res2.noalias() += local_M2 * MathLib::toVector(local_xdot);
267  }
268  if (local_K2.size() != 0)
269  {
270  res2.noalias() += local_K2 * MathLib::toVector(local_x);
271  }
272  if (local_b2.size() != 0)
273  {
274  res2.noalias() -= local_b2;
275  }
276 
277  check_equality(res1, res2);
278 
279  if (tol_exceeded)
280  {
281  WARN("Compare Jacobians: {:s}", msg_tolerance.str());
282  }
283 
284  bool const output = tol_exceeded || fatal_error;
285 
286  if (output)
287  {
288  _log_file << "\n### counter: " << std::to_string(_counter)
289  << " (begin)\n";
290  }
291 
292  if (fatal_error)
293  {
294  _log_file << '\n'
295  << "#######################################################\n"
296  << "# FATAL ERROR: " << msg_fatal << '\n'
297  << "# You cannot expect any meaningful insights "
298  "from the Jacobian data printed below!\n"
299  << "# The reason for the mentioned differences "
300  "might be\n"
301  << "# (a) that the assembly routine has side "
302  "effects or\n"
303  << "# (b) that the assembly routines for M, K "
304  "and b themselves differ.\n"
305  << "#######################################################\n"
306  << '\n';
307  }
308 
309  if (tol_exceeded)
310  {
311  _log_file << "# " << msg_tolerance.str() << "\n\n";
312  }
313 
314  if (output)
315  {
316  dump_py(_log_file, "counter", _counter);
317  dump_py(_log_file, "num_dof", num_dof);
318  dump_py(_log_file, "abs_tol", _abs_tol);
319  dump_py(_log_file, "rel_tol", _rel_tol);
320 
321  _log_file << '\n';
322 
323  dump_py(_log_file, "local_x", local_x);
324  dump_py(_log_file, "local_x_dot", local_xdot);
325  dump_py(_log_file, "dxdot_dx", dxdot_dx);
326  dump_py(_log_file, "dx_dx", dx_dx);
327 
328  _log_file << '\n';
329 
330  dump_py(_log_file, "Jacobian_1", local_Jac1);
331  dump_py(_log_file, "Jacobian_2", local_Jac2);
332 
333  _log_file << '\n';
334 
335  _log_file << "# Jacobian_2 - Jacobian_1\n";
336  dump_py(_log_file, "abs_diff", abs_diff);
337  _log_file << "# Componentwise: 2 * abs_diff / (|Jacobian_1| + "
338  "|Jacobian_2|)\n";
339  dump_py(_log_file, "rel_diff", rel_diff);
340 
341  _log_file << '\n';
342 
343  _log_file << "# Masks: 0 ... tolerance met, 1 ... tolerance exceeded\n";
344  dump_py(_log_file, "abs_diff_mask", abs_diff_mask);
345  dump_py(_log_file, "rel_diff_mask", rel_diff_mask);
346 
347  _log_file << '\n';
348 
349  dump_py(_log_file, "M_1", local_M1);
350  dump_py(_log_file, "M_2", local_M2);
351  if (fatal_error && local_M1.size() == local_M2.size())
352  {
353  dump_py(_log_file, "delta_M", local_M2 - local_M1);
354  _log_file << '\n';
355  }
356 
357  dump_py(_log_file, "K_1", local_K1);
358  dump_py(_log_file, "K_2", local_K2);
359  if (fatal_error && local_K1.size() == local_K2.size())
360  {
361  dump_py(_log_file, "delta_K", local_K2 - local_K1);
362  _log_file << '\n';
363  }
364 
365  dump_py(_log_file, "b_1", local_b_data);
366  dump_py(_log_file, "b_2", local_b_data2);
367  if (fatal_error && local_b1.size() == local_b2.size())
368  {
369  dump_py(_log_file, "delta_b", local_b2 - local_b1);
370  _log_file << '\n';
371  }
372 
373  dump_py(_log_file, "res_1", res1);
374  dump_py(_log_file, "res_2", res2);
375  if (fatal_error)
376  {
377  dump_py(_log_file, "delta_res", res2 - res1);
378  }
379 
380  _log_file << '\n';
381 
382  _log_file << "### counter: " << std::to_string(_counter) << " (end)\n";
383  }
384 
385  if (fatal_error)
386  {
387  _log_file << std::flush;
388  OGS_FATAL("{:s}", msg_fatal);
389  }
390 
391  if (tol_exceeded && _fail_on_error)
392  {
393  _log_file << std::flush;
394  OGS_FATAL(
395  "OGS failed, because the two Jacobian implementations returned "
396  "different results.");
397  }
398 }
399 
400 std::unique_ptr<CompareJacobiansJacobianAssembler>
402 {
403  // TODO doc script corner case: Parameter could occur at different
404  // locations.
406  config.checkConfigParameter("type", "CompareJacobians");
407 
408  auto asm1 =
410  createJacobianAssembler(config.getConfigSubtree("jacobian_assembler"));
411 
412  auto asm2 = createJacobianAssembler(
414  config.getConfigSubtree("reference_jacobian_assembler"));
415 
417  auto const abs_tol = config.getConfigParameter<double>("abs_tol");
419  auto const rel_tol = config.getConfigParameter<double>("rel_tol");
420 
422  auto const fail_on_error = config.getConfigParameter<bool>("fail_on_error");
423 
425  auto const log_file = config.getConfigParameter<std::string>("log_file");
426 
427  return std::make_unique<CompareJacobiansJacobianAssembler>(
428  std::move(asm1), std::move(asm2), abs_tol, rel_tol, fail_on_error,
429  log_file);
430 }
431 
432 } // namespace ProcessLib
#define OGS_FATAL(...)
Definition: Error.h:26
void WARN(char const *fmt, Args const &... args)
Definition: Logging.h:37
void checkConfigParameter(std::string const &param, T const &value) const
T getConfigParameter(std::string const &param) const
ConfigTree getConfigSubtree(std::string const &root) const
Definition: ConfigTree.cpp:146
bool const _fail_on_error
Whether to abort if the tolerances are exceeded.
std::unique_ptr< AbstractJacobianAssembler > _asm2
std::unique_ptr< AbstractJacobianAssembler > _asm1
void assembleWithJacobian(LocalAssemblerInterface &local_assembler, double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_xdot, const double dxdot_dx, const double dx_dx, 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
Eigen::Map< const Vector > toVector(std::vector< double > const &data, Eigen::VectorXd::Index size)
Creates an Eigen mapped vector from the given data vector.
static const double r
Eigen::Map< const Matrix > toMatrix(std::vector< double > const &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
Definition: EigenMapTools.h:72
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, Eigen::MatrixBase< Derived > const &val)
Dumps an Eigen matrix as a Python script snippet.
const std::string msg_fatal
Will be printed if some consistency error is detected.