OGS
ProcessLib::CompareJacobiansJacobianAssembler Class Referencefinal

Detailed Description

Assembles the Jacobian matrix using two different Jacobian assemblers and compares the assembled local Jacobian matrices.

If the provided tolerances are exceeded, debugging information is logged in the form of a Python script.

Definition at line 31 of file CompareJacobiansJacobianAssembler.h.

#include <CompareJacobiansJacobianAssembler.h>

Inheritance diagram for ProcessLib::CompareJacobiansJacobianAssembler:
[legend]
Collaboration diagram for ProcessLib::CompareJacobiansJacobianAssembler:
[legend]

Public Member Functions

 CompareJacobiansJacobianAssembler (std::unique_ptr< AbstractJacobianAssembler > &&asm1, std::unique_ptr< AbstractJacobianAssembler > &&asm2, double abs_tol, double rel_tol, bool fail_on_error, std::string const &log_file_path)
 
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< AbstractJacobianAssemblercopy () const override
 
- Public Member Functions inherited from ProcessLib::AbstractJacobianAssembler
virtual void assembleWithJacobianForStaggeredScheme (LocalAssemblerInterface &, double const, double const, Eigen::VectorXd const &, Eigen::VectorXd const &, int const, std::vector< double > &, std::vector< double > &)
 
virtual ~AbstractJacobianAssembler ()=default
 

Private Attributes

std::unique_ptr< AbstractJacobianAssembler_asm1
 
std::unique_ptr< AbstractJacobianAssembler_asm2
 
double const _abs_tol
 
double const _rel_tol
 
bool const _fail_on_error
 Whether to abort if the tolerances are exceeded.
 
std::ofstream _log_file
 
std::size_t _counter = 0
 

Constructor & Destructor Documentation

◆ CompareJacobiansJacobianAssembler()

ProcessLib::CompareJacobiansJacobianAssembler::CompareJacobiansJacobianAssembler ( std::unique_ptr< AbstractJacobianAssembler > && asm1,
std::unique_ptr< AbstractJacobianAssembler > && asm2,
double abs_tol,
double rel_tol,
bool fail_on_error,
std::string const & log_file_path )
inline

Definition at line 34 of file CompareJacobiansJacobianAssembler.h.

41 : _asm1(std::move(asm1)),
42 _asm2(std::move(asm2)),
43 _abs_tol(abs_tol),
44 _rel_tol(rel_tol),
45 _fail_on_error(fail_on_error),
46 _log_file(log_file_path)
47 {
48 _log_file.precision(std::numeric_limits<double>::max_digits10);
49 _log_file << "#!/usr/bin/env python\n"
50 "import numpy as np\n"
51 "from numpy import nan\n"
52 << std::endl;
53 }
bool const _fail_on_error
Whether to abort if the tolerances are exceeded.
std::unique_ptr< AbstractJacobianAssembler > _asm2
std::unique_ptr< AbstractJacobianAssembler > _asm1

References _log_file.

Member Function Documentation

◆ assembleWithJacobian()

void ProcessLib::CompareJacobiansJacobianAssembler::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 )
overridevirtual

Assembles the Jacobian, the matrices \(M\) and \(K\), and the vector \(b\).

Implements ProcessLib::AbstractJacobianAssembler.

Definition at line 127 of file CompareJacobiansJacobianAssembler.cpp.

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}
#define OGS_FATAL(...)
Definition Error.h:26
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:40
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)
auto eval(Function &f, Tuples &... ts) -> typename detail::GetFunctionReturnType< decltype(&Function::eval)>::type
Definition Apply.h:274
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.

References _abs_tol, _asm1, _asm2, _counter, _fail_on_error, _log_file, _rel_tol, OGS_FATAL, MathLib::toMatrix(), MathLib::toVector(), and WARN().

◆ copy()

std::unique_ptr< AbstractJacobianAssembler > ProcessLib::CompareJacobiansJacobianAssembler::copy ( ) const
overridevirtual

Implements ProcessLib::AbstractJacobianAssembler.

Definition at line 348 of file CompareJacobiansJacobianAssembler.cpp.

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}

References OGS_FATAL.

Member Data Documentation

◆ _abs_tol

double const ProcessLib::CompareJacobiansJacobianAssembler::_abs_tol
private

Definition at line 69 of file CompareJacobiansJacobianAssembler.h.

Referenced by assembleWithJacobian().

◆ _asm1

std::unique_ptr<AbstractJacobianAssembler> ProcessLib::CompareJacobiansJacobianAssembler::_asm1
private

Definition at line 65 of file CompareJacobiansJacobianAssembler.h.

Referenced by assembleWithJacobian().

◆ _asm2

std::unique_ptr<AbstractJacobianAssembler> ProcessLib::CompareJacobiansJacobianAssembler::_asm2
private

Definition at line 66 of file CompareJacobiansJacobianAssembler.h.

Referenced by assembleWithJacobian().

◆ _counter

std::size_t ProcessLib::CompareJacobiansJacobianAssembler::_counter = 0
private

Counter used for identifying blocks in the _log_file. It is incremented upon each call of the assembly routine, i.e., for each element in each iteration etc.

Definition at line 82 of file CompareJacobiansJacobianAssembler.h.

Referenced by assembleWithJacobian().

◆ _fail_on_error

bool const ProcessLib::CompareJacobiansJacobianAssembler::_fail_on_error
private

Whether to abort if the tolerances are exceeded.

Definition at line 73 of file CompareJacobiansJacobianAssembler.h.

Referenced by assembleWithJacobian().

◆ _log_file

std::ofstream ProcessLib::CompareJacobiansJacobianAssembler::_log_file
private

Path where a Python script will be placed, which contains information about exceeded tolerances and assembled local matrices.

Definition at line 77 of file CompareJacobiansJacobianAssembler.h.

Referenced by CompareJacobiansJacobianAssembler(), and assembleWithJacobian().

◆ _rel_tol

double const ProcessLib::CompareJacobiansJacobianAssembler::_rel_tol
private

Definition at line 70 of file CompareJacobiansJacobianAssembler.h.

Referenced by assembleWithJacobian().


The documentation for this class was generated from the following files: