OGS
AssemblyMixin.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
4#pragma once
5
6#include "BaseLib/MPI.h"
7#include "BaseLib/RunTime.h"
9#include "NumLib/Exceptions.h"
17
18namespace ProcessLib
19{
23{
30
31protected:
32 explicit AssemblyMixinBase(AbstractJacobianAssembler& jacobian_assembler,
33 bool const is_linear,
34 bool const use_monolithic_scheme);
35
37 MeshLib::Mesh& bulk_mesh,
38 std::vector<std::reference_wrapper<MeshLib::Mesh>> const& submeshes,
39 std::vector<std::vector<std::string>> const& residuum_names,
40 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>> const&
41 pvs);
42
43 void updateActiveElements(ProcessLib::Process const& process);
44
45 bool isLinear() const { return is_linear_; }
46
48 GlobalVector const& rhs,
49 NumLib::LocalToGlobalIndexMap const& local_to_global_index_map,
50 std::vector<std::reference_wrapper<MeshLib::PropertyVector<double>>>
51 residuum_vectors);
52
54 int const process_id,
55 GlobalVector const& rhs,
56 NumLib::LocalToGlobalIndexMap const& local_to_global_index_map,
58
59private:
60 void updateActiveElementsImpl(Process const& process);
61
62protected:
64 std::vector<ProcessLib::Assembly::SubmeshAssemblyData>
66
68 std::vector<AssembledMatrixCache> submesh_matrix_cache_;
69
71 std::optional<ProcessLib::Assembly::BulkMeshAssemblyData>
73
75
77
81
82 std::optional<NumLib::NonlinearSolverTag> last_assembly_was_;
83
84private:
86};
87
93template <typename Process>
95{
96 // Enforce correct use of CRTP, i.e., that Process is derived from
97 // AssemblyMixin<Process>.
99 friend Process;
100
101public:
124 std::vector<std::reference_wrapper<MeshLib::Mesh>> const& submeshes,
125 std::vector<std::vector<std::string>> const& residuum_names)
126 {
128 derived().getMesh(), submeshes, residuum_names,
129 derived().getProcessVariables());
130 }
131
136
137 // sorted_element_subset can be used to restrict assembly to a custom
138 // element subset. That's used at the moment (Jan 2026) for some assembly
139 // optimizations in the HeatTransportBHE process.
141 double const t, double const dt, std::vector<GlobalVector*> const& x,
142 std::vector<GlobalVector*> const& x_prev, int const process_id,
144 std::vector<std::size_t> const* const sorted_element_subset = nullptr,
145 bool const copy_residua_to_mesh = false)
146 {
147 DBUG("AssemblyMixin assemble(t={:g}, dt={:g}, process_id={}).", t, dt,
148 process_id);
149
151
152 std::vector<NumLib::LocalToGlobalIndexMap const*> const dof_tables =
153 derived().getDOFTables(x.size());
154
155 std::exception_ptr const exception =
157 ? assembleOnBulkMesh(t, dt, x, x_prev, dof_tables, process_id,
158 M, K, b, sorted_element_subset)
159 : assembleOnSubmeshes(t, dt, x, x_prev, dof_tables, process_id,
160 M, K, b, sorted_element_subset,
161 copy_residua_to_mesh);
162
163 [[unlikely]] if (BaseLib::MPI::anyOf(exception != nullptr))
164 {
165 if (exception) // Only the rank with the exception rethrows...
166 {
167 std::rethrow_exception(exception);
168 }
169 // TODO should other ranks rather throw a new exception?
170 // ... but all ranks quit.
171 return;
172 }
173
174 if (copy_residua_to_mesh && bulk_mesh_assembly_data_)
175 {
176 // Note: computeResiduum() uses bwd/fwd Euler scheme.
177 GlobalVector const neg_res_bulkmesh = ProcessLib::computeResiduum(
178 dt, *x[process_id], *x_prev[process_id], M, K, b);
179
181 neg_res_bulkmesh, *(dof_tables[process_id]),
182 bulk_mesh_assembly_data_->residuum_vectors[process_id]);
183 }
184 }
185
186 // sorted_element_subset can be used to restrict assembly to a custom
187 // element subset. That's used at the moment (Jan 2026) for some assembly
188 // optimizations in the HeatTransportBHE process.
190 double const t, double const dt, std::vector<GlobalVector*> const& x,
191 std::vector<GlobalVector*> const& x_prev, int const process_id,
192 GlobalVector& b, GlobalMatrix& Jac,
193 std::vector<std::size_t> const* const sorted_element_subset = nullptr,
194 bool const copy_residua_to_mesh = false)
195 {
196 DBUG(
197 "AssemblyMixin assembleWithJacobian(t={:g}, dt={:g}, "
198 "process_id={}).",
199 t, dt, process_id);
200
202
203 // Note: We do not cache the Jacobian or b vector. It does not make
204 // sense, because even in the linear case the b vector changes if
205 // the solution changes (unless steady state is reached).
206 // Todo: If we separated Jacobian and residuum assembly, we could
207 // cache at least the Jacobian and only update the residuum.
208 if (is_linear_)
209 {
210 DBUG(
211 "You specified that this process is linear. The assembly for "
212 "the Newton-Raphson method will be run anyways.");
213 }
214
215 std::vector<NumLib::LocalToGlobalIndexMap const*> const dof_tables =
216 derived().getDOFTables(x.size());
217
218 std::exception_ptr const exception =
221 *bulk_mesh_assembly_data_, t, dt, x, x_prev, dof_tables,
222 process_id, b, Jac, sorted_element_subset)
224 t, dt, x, x_prev, dof_tables, process_id, b, Jac,
225 sorted_element_subset, copy_residua_to_mesh);
226
227 [[unlikely]] if (BaseLib::MPI::anyOf(exception != nullptr))
228 {
229 if (exception) // Only the rank with the exception rethrows...
230 {
231 std::rethrow_exception(exception);
232 }
233 // ... but all ranks quit.
234 return;
235 }
236
238
239 // TODO better for submesh case, too?
240 // TODO check copy_residua_to_mesh, too. But that needs a refactoring of
241 // postTimestep(), computeSecondaryVariable() etc. hooks.
242 if (/*copy_residua_to_mesh &&*/ bulk_mesh_assembly_data_)
243 {
245 b, *(dof_tables[process_id]),
246 bulk_mesh_assembly_data_->residuum_vectors[process_id]);
247 }
248 }
249
250 void preOutput(double const t,
251 double const dt,
252 std::vector<GlobalVector*> const& x,
253 std::vector<GlobalVector*> const& x_prev,
254 int const process_id)
255 {
257 {
258 WARN("AssemblyMixin. Skipping preOutput() step.");
259 return;
260 }
261
262 switch (*last_assembly_was_)
263 {
265 preOutputPicard(t, dt, x, x_prev, process_id);
266 return;
268 preOutputNewton(t, dt, x, x_prev, process_id);
269 return;
270 }
271 }
272
273private:
274 Process& derived() { return static_cast<Process&>(*this); }
275 Process const& derived() const
276 {
277 return static_cast<Process const&>(*this);
278 }
279
281 [[nodiscard]] std::exception_ptr assembleOnBulkMeshOrOnSubmeshCommon(
282 Assembly::CommonAssemblyData const& assembly_data,
283 AssembledMatrixCache& mat_cache, double const t, double const dt,
284 std::vector<GlobalVector*> const& x,
285 std::vector<GlobalVector*> const& x_prev,
286 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
287 int const process_id, GlobalMatrix& M, GlobalMatrix& K, GlobalVector& b,
288 std::vector<std::size_t> const* const sorted_element_subset)
289 {
290 std::exception_ptr exception = nullptr;
291
292 try
293 {
294 auto const& loc_asms = derived().local_assemblers_;
295
296 pvma_.assemble(
297 loc_asms,
298 assembly_data.activeElementIDsSorted(sorted_element_subset)
299 .get(),
300 dof_tables, t, dt, x, x_prev, process_id, M, K, b);
301 }
302 catch (NumLib::AssemblyException const&)
303 {
304 exception = std::current_exception();
305 }
306
310
311 if (is_linear_)
312 {
313 DBUG("Saving global M, K, b for later reuse.");
314
315 BaseLib::RunTime time_save;
316 time_save.start();
317
318 mat_cache.storeMKb(M, K, b);
319
320 INFO("[time] Saving global M, K, b took {:g} s",
321 time_save.elapsed());
322 }
323
324 return exception;
325 }
326
327 std::exception_ptr assembleOnBulkMesh(
328 double const t, double const dt, std::vector<GlobalVector*> const& x,
329 std::vector<GlobalVector*> const& x_prev,
330 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
331 int const process_id, GlobalMatrix& M, GlobalMatrix& K, GlobalVector& b,
332 std::vector<std::size_t> const* const sorted_element_subset)
333 {
334 if (!bulk_mesh_matrix_cache_.hasMKb())
335 {
338 x_prev, dof_tables, process_id, M, K, b, sorted_element_subset);
339 }
340
341 DBUG("Reusing saved global M, K, b.");
342
343 BaseLib::RunTime time_restore;
344 time_restore.start();
345
346 auto const [M_, K_, b_] = bulk_mesh_matrix_cache_.MKb();
350
351 INFO("[time] Restoring global M, K, b took {:g} s",
352 time_restore.elapsed());
353
354 return nullptr;
355 }
356
357 [[nodiscard]] std::exception_ptr assembleOnSubmeshes(
358 double const t, double const dt, std::vector<GlobalVector*> const& x,
359 std::vector<GlobalVector*> const& x_prev,
360 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
361 int const process_id, GlobalMatrix& M, GlobalMatrix& K, GlobalVector& b,
362 std::vector<std::size_t> const* const sorted_element_subset,
363 bool const copy_residua_to_mesh)
364 {
365 auto const mat_spec = derived().getMatrixSpecifications(process_id);
366
367 // Temporary matrices for assembly (reused across submeshes).
368 // They will remain nullptr if no assembly is necessary (linear case &
369 // cache populated).
370 std::unique_ptr<GlobalMatrix> M_tmp;
371 std::unique_ptr<GlobalMatrix> K_tmp;
372 std::unique_ptr<GlobalVector> b_tmp;
373
374 auto getMKb = [&M_tmp, &K_tmp, &b_tmp,
375 &mat_spec](AssembledMatrixCache const& cache)
376 {
377 if (cache.hasMKb())
378 {
379 return cache.MKb();
380 }
381
382 if (!(M_tmp && K_tmp && b_tmp))
383 {
385 mat_spec);
387 mat_spec);
389 mat_spec);
390 }
391
392 return std::tie(*M_tmp, *K_tmp, *b_tmp);
393 };
394
395 std::exception_ptr exception = nullptr;
396 for (auto const& [sad, mat_cache] :
397 ranges::zip_view(submesh_assembly_data_ | ranges::views::all,
398 submesh_matrix_cache_ | ranges::views::all))
399 {
400 auto [M_submesh, K_submesh, b_submesh] = getMKb(mat_cache);
401
402 if (!mat_cache.hasMKb())
403 {
404 M_submesh.setZero();
405 K_submesh.setZero();
406 b_submesh.setZero();
407
409 sad, mat_cache, t, dt, x, x_prev, dof_tables, process_id,
410 M_submesh, K_submesh, b_submesh, sorted_element_subset);
411 }
412
413 // TODO: This operates on the entire global vector b once
414 // for each submesh. If the number of submeshes is high,
415 // that might not perform well.
416 MathLib::LinAlg::axpy(M, 1.0, M_submesh);
417 MathLib::LinAlg::axpy(K, 1.0, K_submesh);
418 MathLib::LinAlg::axpy(b, 1.0, b_submesh);
419
420 if (copy_residua_to_mesh)
421 {
422 // Note: computeResiduum() uses bwd/fwd Euler scheme.
423 GlobalVector const neg_res_submesh =
424 computeResiduum(dt, *x[process_id], *x_prev[process_id],
425 M_submesh, K_submesh, b_submesh);
426
428 process_id, neg_res_submesh, *(dof_tables[process_id]),
429 sad);
430 }
431 }
432
433 return exception;
434 }
435
438 [[nodiscard]] std::exception_ptr
440 Assembly::CommonAssemblyData const& assembly_data, double const t,
441 double const dt, std::vector<GlobalVector*> const& x,
442 std::vector<GlobalVector*> const& x_prev,
443 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
444 int const process_id, GlobalVector& b, GlobalMatrix& Jac,
445 std::vector<std::size_t> const* const sorted_element_subset)
446 {
447 auto const& loc_asms = derived().local_assemblers_;
448 std::exception_ptr exception = nullptr;
449
450 try
451 {
452 pvma_.assembleWithJacobian(
453 loc_asms,
454 assembly_data.activeElementIDsSorted(sorted_element_subset)
455 .get(),
456 dof_tables, t, dt, x, x_prev, process_id, b, Jac);
457 }
458 catch (NumLib::AssemblyException const&)
459 {
460 exception = std::current_exception();
461 }
462
464
465 return exception;
466 }
467
468 [[nodiscard]] std::exception_ptr assembleWithJacobianOnSubmeshes(
469 double const t, double const dt, std::vector<GlobalVector*> const& x,
470 std::vector<GlobalVector*> const& x_prev,
471 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
472 int const process_id, GlobalVector& b, GlobalMatrix& Jac,
473 std::vector<std::size_t> const* const sorted_element_subset,
474 bool const /*copy_residua_to_mesh*/)
475 {
476 std::exception_ptr exception = nullptr;
477
478 auto const& mat_spec = derived().getMatrixSpecifications(process_id);
479 auto b_submesh =
481
482 for (auto const& sad : submesh_assembly_data_)
483 {
484 b_submesh->setZero();
485
487 sad, t, dt, x, x_prev, dof_tables, process_id, *b_submesh, Jac,
488 sorted_element_subset);
489
490 MathLib::LinAlg::axpy(b, 1.0, *b_submesh);
491
492 // TODO enabling this check needs a refactoring of postTimestep(),
493 // computeSecondaryVariable() etc. hooks.
494 if (/*copy_residua_to_mesh*/ true)
495 {
497 process_id, *b_submesh, *(dof_tables[process_id]), sad);
498 }
499 }
500
501 return exception;
502 }
503
504 void preOutputPicard(double const t,
505 double const dt,
506 std::vector<GlobalVector*> const& x,
507 std::vector<GlobalVector*> const& x_prev,
508 int const process_id)
509 {
510 auto const matrix_specification =
511 derived().getMatrixSpecifications(process_id);
512
514 matrix_specification);
516 matrix_specification);
518 matrix_specification);
519
520 M->setZero();
521 K->setZero();
522 b->setZero();
523
524 assemble(t, dt, x, x_prev, process_id, *M, *K, *b, nullptr, true);
525 }
526
527 void preOutputNewton(double const /*t*/,
528 double const /*dt*/,
529 std::vector<GlobalVector*> const& /*x*/,
530 std::vector<GlobalVector*> const& /*x_prev*/,
531 int const /*process_id*/)
532 {
533 // TODO enabling this needs a refactoring of postTimestep(),
534 // computeSecondaryVariable() etc. hooks.
535#if 0
536 auto const matrix_specification =
537 derived().getMatrixSpecifications(process_id);
538
540 matrix_specification);
542 matrix_specification);
543
544 Jac->setZero();
545 b->setZero();
546
547 assembleWithJacobian(t, dt, x, x_prev, process_id, *b, *Jac, nullptr,
548 true);
549#endif
550 }
551};
552
553} // namespace ProcessLib
MathLib::EigenMatrix GlobalMatrix
MathLib::EigenVector GlobalVector
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:28
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:22
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:34
Count the running time.
Definition RunTime.h:18
double elapsed() const
Get the elapsed time in seconds.
Definition RunTime.h:31
void start()
Start the timer.
Definition RunTime.h:21
Base class for Jacobian assemblers.
void initializeAssemblyOnSubmeshes(MeshLib::Mesh &bulk_mesh, std::vector< std::reference_wrapper< MeshLib::Mesh > > const &submeshes, std::vector< std::vector< std::string > > const &residuum_names, std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > const &pvs)
AssemblyMixinBase(AbstractJacobianAssembler &jacobian_assembler, bool const is_linear, bool const use_monolithic_scheme)
Assembly::ParallelVectorMatrixAssembler pvma_
static void copyResiduumVectorsToBulkMesh(GlobalVector const &rhs, NumLib::LocalToGlobalIndexMap const &local_to_global_index_map, std::vector< std::reference_wrapper< MeshLib::PropertyVector< double > > > residuum_vectors)
std::optional< NumLib::NonlinearSolverTag > last_assembly_was_
void updateActiveElementsImpl(Process const &process)
std::vector< AssembledMatrixCache > submesh_matrix_cache_
AssembledMatrixCache for each submesh.
static void copyResiduumVectorsToSubmesh(int const process_id, GlobalVector const &rhs, NumLib::LocalToGlobalIndexMap const &local_to_global_index_map, ProcessLib::Assembly::SubmeshAssemblyData const &sad)
ActiveElementIDsState ids_state_
void updateActiveElements(ProcessLib::Process const &process)
std::vector< ProcessLib::Assembly::SubmeshAssemblyData > submesh_assembly_data_
SubmeshAssemblyData for each submesh.
std::optional< ProcessLib::Assembly::BulkMeshAssemblyData > bulk_mesh_assembly_data_
Empty if submesh assembly is used.
AssembledMatrixCache bulk_mesh_matrix_cache_
AssemblyMixinBase(AbstractJacobianAssembler &jacobian_assembler, bool const is_linear, bool const use_monolithic_scheme)
void assemble(double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, std::vector< std::size_t > const *const sorted_element_subset=nullptr, bool const copy_residua_to_mesh=false)
std::exception_ptr assembleOnSubmeshes(double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, std::vector< std::size_t > const *const sorted_element_subset, bool const copy_residua_to_mesh)
std::exception_ptr assembleOnBulkMeshOrOnSubmeshCommon(Assembly::CommonAssemblyData const &assembly_data, AssembledMatrixCache &mat_cache, double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, std::vector< std::size_t > const *const sorted_element_subset)
Common code for Picard assembly on the bulk mesh or on submeshes.
Process const & derived() const
std::exception_ptr assembleWithJacobianOnBulkMeshOrOnSubmeshCommon(Assembly::CommonAssemblyData const &assembly_data, double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, int const process_id, GlobalVector &b, GlobalMatrix &Jac, std::vector< std::size_t > const *const sorted_element_subset)
void preOutputPicard(double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id)
std::exception_ptr assembleOnBulkMesh(double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, std::vector< std::size_t > const *const sorted_element_subset)
void preOutput(double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id)
void preOutputNewton(double const, double const, std::vector< GlobalVector * > const &, std::vector< GlobalVector * > const &, int const)
void initializeAssemblyOnSubmeshes(std::vector< std::reference_wrapper< MeshLib::Mesh > > const &submeshes, std::vector< std::vector< std::string > > const &residuum_names)
void assembleWithJacobian(double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalVector &b, GlobalMatrix &Jac, std::vector< std::size_t > const *const sorted_element_subset=nullptr, bool const copy_residua_to_mesh=false)
std::exception_ptr assembleWithJacobianOnSubmeshes(double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, int const process_id, GlobalVector &b, GlobalMatrix &Jac, std::vector< std::size_t > const *const sorted_element_subset, bool const)
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables(int const number_of_processes) const
Definition Process.cpp:376
MathLib::MatrixSpecifications getMatrixSpecifications(const int process_id) const override
Definition Process.cpp:203
static bool anyOf(bool const val, Mpi const &mpi=Mpi{OGS_COMM_WORLD})
Definition MPI.h:185
void finalizeAssembly(PETScMatrix &A)
Definition LinAlg.cpp:191
void copy(PETScVector const &x, PETScVector &y)
Definition LinAlg.cpp:30
void axpy(PETScVector &y, PetscScalar const a, PETScVector const &x)
Definition LinAlg.cpp:50
GlobalVector computeResiduum(double const dt, GlobalVector const &x, GlobalVector const &x_prev, GlobalMatrix const &M, GlobalMatrix const &K, GlobalVector const &b)
Stores assembled global M, K, b matrices/vectors for later reuse.
void storeMKb(GlobalMatrix const &M, GlobalMatrix const &K, GlobalVector const &b)
Store data.
std::shared_ptr< std::vector< std::size_t > const > activeElementIDsSorted(std::vector< std::size_t > const *const sorted_element_subset) const