OGS
ComponentTransportProcess.cpp
Go to the documentation of this file.
1
12
13#include <cassert>
14#include <range/v3/algorithm/copy.hpp>
15#include <range/v3/view/drop.hpp>
16
17#include "BaseLib/RunTime.h"
35
36namespace ProcessLib
37{
38namespace ComponentTransport
39{
41 std::string name,
42 MeshLib::Mesh& mesh,
43 std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
44 std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
45 unsigned const integration_order,
46 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
47 process_variables,
48 ComponentTransportProcessData&& process_data,
49 SecondaryVariableCollection&& secondary_variables,
50 bool const use_monolithic_scheme,
51 std::unique_ptr<ProcessLib::SurfaceFluxData>&& surfaceflux,
52 std::unique_ptr<ChemistryLib::ChemicalSolverInterface>&&
53 chemical_solver_interface,
54 bool const is_linear,
55 bool const ls_compute_only_upon_timestep_change)
56 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
57 integration_order, std::move(process_variables),
58 std::move(secondary_variables), use_monolithic_scheme),
59 _process_data(std::move(process_data)),
60 _surfaceflux(std::move(surfaceflux)),
61 _chemical_solver_interface(std::move(chemical_solver_interface)),
62 _asm_mat_cache{is_linear, use_monolithic_scheme},
63 _ls_compute_only_upon_timestep_change{
64 ls_compute_only_upon_timestep_change}
65{
66 if (ls_compute_only_upon_timestep_change)
67 {
68 // TODO move this feature to some common location for all processes.
69 if (!is_linear)
70 {
72 "Using the linear solver compute() method only upon timestep "
73 "change only makes sense for linear model equations.");
74 }
75
76 WARN(
77 "You specified that the ComponentTransport linear solver will do "
78 "the compute() step only upon timestep change. This is an expert "
79 "option. It is your responsibility to ensure that "
80 "the conditions for the correct use of this feature are met! "
81 "Otherwise OGS might compute garbage without being recognized. "
82 "There is no safety net!");
83
84 WARN(
85 "You specified that the ComponentTransport linear solver will do "
86 "the compute() step only upon timestep change. This option will "
87 "only be used by the Picard non-linear solver. The Newton-Raphson "
88 "non-linear solver will silently ignore this setting, i.e., it "
89 "won't do any harm, there, but you won't observe the speedup you "
90 "probably expect.");
91 }
92
93 auto residuum_name = [&](auto const& pv) -> std::string
94 {
95 std::string const& pv_name = pv.getName();
96 if (pv_name == "pressure")
97 {
98 return "LiquidMassFlowRate";
99 }
100 if (pv_name == "temperature")
101 {
102 return "HeatFlowRate";
103 }
104 return pv_name + "FlowRate";
105 };
106
108 {
109 int const process_id = 0;
110 for (auto const& pv : _process_variables[process_id])
111 {
112 _residua.push_back(MeshLib::getOrCreateMeshProperty<double>(
113 mesh, residuum_name(pv.get()), MeshLib::MeshItemType::Node, 1));
114 }
115 }
116 else
117 {
118 for (auto const& pv : _process_variables)
119 {
120 _residua.push_back(MeshLib::getOrCreateMeshProperty<double>(
121 mesh, residuum_name(pv[0].get()), MeshLib::MeshItemType::Node,
122 1));
123 }
124 }
125}
126
128 NumLib::LocalToGlobalIndexMap const& dof_table,
129 MeshLib::Mesh const& mesh,
130 unsigned const integration_order)
131{
132 int const mesh_space_dimension = _process_data.mesh_space_dimension;
133
134 _process_data.mesh_prop_velocity = MeshLib::getOrCreateMeshProperty<double>(
135 const_cast<MeshLib::Mesh&>(mesh), "velocity",
136 MeshLib::MeshItemType::Cell, mesh_space_dimension);
137
138 _process_data.mesh_prop_porosity = MeshLib::getOrCreateMeshProperty<double>(
139 const_cast<MeshLib::Mesh&>(mesh), "porosity_avg",
141
142 std::vector<std::reference_wrapper<ProcessLib::ProcessVariable>>
143 transport_process_variables;
145 {
146 const int process_id = 0;
147 for (auto pv_iter = std::next(_process_variables[process_id].begin());
148 pv_iter != _process_variables[process_id].end();
149 ++pv_iter)
150 {
151 transport_process_variables.push_back(*pv_iter);
152 }
153 }
154 else
155 {
156 // All process variables but the pressure and optionally the temperature
157 // are transport variables.
158 for (auto const& pv :
160 ranges::views::drop(_process_data.isothermal ? 1 : 2))
161 {
162 transport_process_variables.push_back(pv[0]);
163 }
164 }
165
166 ProcessLib::createLocalAssemblers<LocalAssemblerData>(
167 mesh_space_dimension, mesh.getElements(), dof_table, _local_assemblers,
168 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
169 _process_data, transport_process_variables);
170
172 {
176
177 _chemical_solver_interface->initialize();
178 }
179
181 "darcy_velocity",
183 mesh_space_dimension, getExtrapolator(), _local_assemblers,
185
186 for (std::size_t component_id = 0;
187 component_id < transport_process_variables.size();
188 ++component_id)
189 {
190 auto const& pv = transport_process_variables[component_id].get();
191
192 auto integration_point_values_accessor = [component_id](
194 const double t,
195 std::vector<GlobalVector*> const& x,
196 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
197 std::vector<double>& cache) -> auto const&
198 {
199 return loc_asm.getIntPtMolarFlux(t, x, dof_tables, cache,
200 component_id);
201 };
202
204 pv.getName() + "Flux",
205 makeExtrapolator(mesh_space_dimension, getExtrapolator(),
207 integration_point_values_accessor));
208 }
209}
210
212 std::vector<GlobalVector*>& x, double const t, int const process_id)
213{
215 {
216 return;
217 }
218
219 if (process_id != static_cast<int>(x.size() - 1))
220 {
221 return;
222 }
223
224 std::for_each(
225 x.begin(), x.end(),
226 [](auto const process_solution)
227 { MathLib::LinAlg::setLocalAccessibleVector(*process_solution); });
228
229 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
230 dof_tables.reserve(x.size());
231 std::generate_n(std::back_inserter(dof_tables), x.size(),
232 [&]() { return _local_to_global_index_map.get(); });
233
237 dof_tables, x, t);
238}
239
241 const double t, double const dt, std::vector<GlobalVector*> const& x,
242 std::vector<GlobalVector*> const& x_prev, int const process_id,
244{
245 DBUG("Assemble ComponentTransportProcess.");
246
247 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
249 {
250 dof_tables.push_back(_local_to_global_index_map.get());
251 }
252 else
253 {
254 std::generate_n(std::back_inserter(dof_tables),
255 _process_variables.size(),
256 [&]() { return _local_to_global_index_map.get(); });
257 }
258
259 _asm_mat_cache.assemble(t, dt, x, x_prev, process_id, &M, &K, &b,
262
263 // TODO (naumov) What about temperature? A test with FCT would reveal any
264 // problems.
265 if (process_id == _process_data.hydraulic_process_id)
266 {
267 return;
268 }
269 auto const matrix_specification = getMatrixSpecifications(process_id);
270
272 x_prev, process_id,
273 matrix_specification, M, K, b);
274}
275
277 const double t, double const dt, std::vector<GlobalVector*> const& x,
278 std::vector<GlobalVector*> const& x_prev, int const process_id,
279 GlobalVector& b, GlobalMatrix& Jac)
280{
281 DBUG("AssembleWithJacobian ComponentTransportProcess.");
282
284 {
285 OGS_FATAL(
286 "The AssembleWithJacobian() for ComponentTransportProcess is "
287 "implemented only for staggered scheme.");
288 }
289
290 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
291
292 std::generate_n(std::back_inserter(dof_tables), _process_variables.size(),
293 [&]() { return _local_to_global_index_map.get(); });
294
295 // Call global assembler for each local assembly item.
298 _local_assemblers, getActiveElementIDs(), dof_tables, t, dt, x, x_prev,
299 process_id, &b, &Jac);
300
301 // b is the negated residumm used in the Newton's method.
302 // Here negating b is to recover the primitive residuum.
303 transformVariableFromGlobalVector(b, 0, *_local_to_global_index_map,
304 *_residua[process_id],
305 std::negate<double>());
306}
307
309 std::size_t const element_id,
310 MathLib::Point3d const& p,
311 double const t,
312 std::vector<GlobalVector*> const& x) const
313{
314 std::vector<GlobalIndexType> indices_cache;
315 auto const r_c_indices = NumLib::getRowColumnIndices(
316 element_id, *_local_to_global_index_map, indices_cache);
317
318 std::vector<std::vector<GlobalIndexType>> indices_of_all_coupled_processes{
319 x.size(), r_c_indices.rows};
320 auto const local_xs =
321 getCoupledLocalSolutions(x, indices_of_all_coupled_processes);
322
323 return _local_assemblers[element_id]->getFlux(p, t, local_xs);
324}
325
327 std::vector<GlobalVector*>& x, std::vector<GlobalVector*> const& x_prev,
328 double const t, double const dt, NumLib::EquationSystem& ode_sys,
329 int const process_id)
330{
331 // todo (renchao): move chemical calculation to elsewhere.
332 if (_process_data.lookup_table && process_id == 0)
333 {
334 INFO("Update process solutions via the look-up table approach");
335 _process_data.lookup_table->lookup(x, x_prev, _mesh.getNumberOfNodes());
336
337 return;
338 }
339
341 {
342 return;
343 }
344
345 // Sequential non-iterative approach applied here to split the reactive
346 // transport process into the transport stage followed by the reaction
347 // stage.
348 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
349 dof_tables.reserve(x.size());
350 std::generate_n(std::back_inserter(dof_tables), x.size(),
351 [&]() { return _local_to_global_index_map.get(); });
352
353 if (process_id == 0)
354 {
358 dof_tables, x, t, dt);
359
360 BaseLib::RunTime time_phreeqc;
361 time_phreeqc.start();
362
363 _chemical_solver_interface->setAqueousSolutionsPrevFromDumpFile();
364
365 _chemical_solver_interface->executeSpeciationCalculation(dt);
366
367 INFO("[time] Phreeqc took {:g} s.", time_phreeqc.elapsed());
368
371 postSpeciationCalculation,
373 dt);
374
375 return;
376 }
377
378 auto const matrix_specification = getMatrixSpecifications(process_id);
379
380 std::size_t matrix_id = 0u;
382 matrix_specification, matrix_id);
384 matrix_specification, matrix_id);
385 auto& b =
387
388 M.setZero();
389 K.setZero();
390 b.setZero();
391
394 _local_assemblers, getActiveElementIDs(), dof_tables, x, t, dt, M, K, b,
395 process_id);
396
397 // todo (renchao): incorporate Neumann boundary condition
401
403 matrix_specification, matrix_id);
404 auto& rhs =
406
407 A.setZero();
408 rhs.setZero();
409
412
413 // compute A
415 MathLib::LinAlg::aypx(A, 1.0 / dt, K);
416
417 // compute rhs
418 MathLib::LinAlg::matMult(M, *x[process_id], rhs);
419 MathLib::LinAlg::aypx(rhs, 1.0 / dt, b);
420
421 using Tag = NumLib::NonlinearSolverTag;
423 auto& equation_system = static_cast<EqSys&>(ode_sys);
424 equation_system.applyKnownSolutionsPicard(A, rhs, *x[process_id]);
425
426 auto& linear_solver =
429 linear_solver.solve(A, rhs, *x[process_id]);
430
436}
437
439 double const t,
440 double const dt,
441 std::vector<GlobalVector*> const& x,
442 GlobalVector const& x_prev,
443 int const process_id)
444{
445 if (process_id != 0)
446 {
447 return;
448 }
449
450 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
451 dof_tables.reserve(x.size());
452 std::generate_n(std::back_inserter(dof_tables), x.size(),
453 [&]() { return _local_to_global_index_map.get(); });
454
457 _local_assemblers, getActiveElementIDs(), dof_tables, t, dt, x, x_prev,
458 process_id);
459
461 {
462 return;
463 }
464
467 computeReactionRelatedSecondaryVariable,
469}
470
472 std::vector<GlobalVector*> const& x,
473 std::vector<GlobalVector*> const& x_prev,
474 const double t,
475 const double dt,
476 int const process_id)
477{
478 if (process_id != 0)
479 {
480 return;
481 }
482
483 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
484 dof_tables.reserve(x.size());
485 std::generate_n(std::back_inserter(dof_tables), x.size(),
486 [&]() { return _local_to_global_index_map.get(); });
487
490 _local_assemblers, getActiveElementIDs(), dof_tables, x, x_prev, t, dt,
491 process_id);
492
493 if (!_surfaceflux) // computing the surfaceflux is optional
494 {
495 return;
496 }
497 _surfaceflux->integrate(x, t, *this, process_id, _integration_order, _mesh,
499}
500
502 const double t,
503 double const dt,
504 std::vector<GlobalVector*> const& x,
505 std::vector<GlobalVector*> const& x_prev,
506 int const process_id)
507{
508 auto const matrix_specification = getMatrixSpecifications(process_id);
509
511 matrix_specification);
513 matrix_specification);
515 matrix_specification);
516
517 M->setZero();
518 K->setZero();
519 b->setZero();
520
521 assembleConcreteProcess(t, dt, x, x_prev, process_id, *M, *K, *b);
522
523 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
525 {
526 dof_tables.push_back(_local_to_global_index_map.get());
527 }
528 else
529 {
530 std::generate_n(std::back_inserter(dof_tables),
531 _process_variables.size(),
532 [&]() { return _local_to_global_index_map.get(); });
533 }
534
535 BaseLib::RunTime time_residuum;
536 time_residuum.start();
537
539 {
540 auto const residuum =
541 computeResiduum(dt, *x[0], *x_prev[0], *M, *K, *b);
542 for (std::size_t variable_id = 0; variable_id < _residua.size();
543 ++variable_id)
544 {
545 transformVariableFromGlobalVector(
546 residuum, variable_id, *dof_tables[0], *_residua[variable_id],
547 std::negate<double>());
548 }
549 }
550 else
551 {
552 auto const residuum = computeResiduum(dt, *x[process_id],
553 *x_prev[process_id], *M, *K, *b);
554 transformVariableFromGlobalVector(residuum, 0, *dof_tables[process_id],
555 *_residua[process_id],
556 std::negate<double>());
557 }
558
559 INFO("[time] Computing residuum flow rates took {:g} s",
560 time_residuum.elapsed());
561}
562} // namespace ComponentTransport
563} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:26
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:35
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:40
Definition of the RunTime class.
Count the running time.
Definition RunTime.h:29
double elapsed() const
Get the elapsed time in seconds.
Definition RunTime.h:42
void start()
Start the timer.
Definition RunTime.h:32
Global vector based on Eigen vector.
Definition EigenVector.h:25
bool isAxiallySymmetric() const
Definition Mesh.h:137
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition Mesh.h:109
std::size_t getNumberOfNodes() const
Get the number of nodes.
Definition Mesh.h:100
virtual void releaseMatrix(GlobalMatrix const &A)=0
virtual GlobalMatrix & getMatrix(std::size_t &id)=0
Get an uninitialized matrix with the given id.
virtual GlobalVector & getVector(std::size_t &id)=0
Get an uninitialized vector with the given id.
virtual void releaseVector(GlobalVector const &x)=0
virtual std::vector< double > const & getIntPtDarcyVelocity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
void initializeChemicalSystem(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t)
void setChemicalSystem(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, double const dt)
void assembleReactionEquation(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, double const dt, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, int const process_id)
virtual std::vector< double > const & getIntPtMolarFlux(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache, int const component_id) const =0
void setInitialConditionsConcreteProcess(std::vector< GlobalVector * > &x, double const t, int const process_id) override
void preOutputConcreteProcess(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id) override
void assembleConcreteProcess(const double 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) override
void initializeConcreteProcess(NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
Process specific initialization called by initialize().
std::unique_ptr< ChemistryLib::ChemicalSolverInterface > _chemical_solver_interface
std::vector< MeshLib::PropertyVector< double > * > _residua
void computeSecondaryVariableConcrete(double const, double const, std::vector< GlobalVector * > const &x, GlobalVector const &, int const) override
void solveReactionEquation(std::vector< GlobalVector * > &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, NumLib::EquationSystem &ode_sys, int const process_id) override
void assembleWithJacobianConcreteProcess(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalVector &b, GlobalMatrix &Jac) override
void postTimestepConcreteProcess(std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, const double dt, int const process_id) override
std::vector< std::unique_ptr< ComponentTransportLocalAssemblerInterface > > _local_assemblers
std::unique_ptr< ProcessLib::SurfaceFluxData > _surfaceflux
ComponentTransportProcess(std::string name, MeshLib::Mesh &mesh, std::unique_ptr< ProcessLib::AbstractJacobianAssembler > &&jacobian_assembler, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, unsigned const integration_order, std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > &&process_variables, ComponentTransportProcessData &&process_data, SecondaryVariableCollection &&secondary_variables, bool const use_monolithic_scheme, std::unique_ptr< ProcessLib::SurfaceFluxData > &&surfaceflux, std::unique_ptr< ChemistryLib::ChemicalSolverInterface > &&chemical_solver_interface, bool const is_linear, bool const ls_compute_only_upon_timestep_change)
Eigen::Vector3d getFlux(std::size_t const element_id, MathLib::Point3d const &p, double const t, std::vector< GlobalVector * > const &x) const override
virtual void postTimestep(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
virtual void computeSecondaryVariable(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id)
MeshLib::Mesh & _mesh
Definition Process.h:365
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > _process_variables
Definition Process.h:400
std::vector< std::size_t > const & getActiveElementIDs() const
Definition Process.h:167
SecondaryVariableCollection _secondary_variables
Definition Process.h:370
VectorMatrixAssembler _global_assembler
Definition Process.h:377
unsigned const _integration_order
Definition Process.h:384
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:368
MathLib::MatrixSpecifications getMatrixSpecifications(const int process_id) const override
Definition Process.cpp:211
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:208
const bool _use_monolithic_scheme
Definition Process.h:379
Handles configuration of several secondary variables from the project file.
void addSecondaryVariable(std::string const &internal_name, SecondaryVariableFunctions &&fcts)
void assembleWithJacobian(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalVector *b, GlobalMatrix *Jac)
NonlinearSolverTag
Tag used to specify which nonlinear solver will be used.
Definition Types.h:20
void copy(PETScVector const &x, PETScVector &y)
Definition LinAlg.cpp:37
void setLocalAccessibleVector(PETScVector const &x)
Definition LinAlg.cpp:27
void matMult(PETScMatrix const &A, PETScVector const &x, PETScVector &y)
Definition LinAlg.cpp:149
void aypx(PETScVector &y, PetscScalar const a, PETScVector const &x)
Definition LinAlg.cpp:50
void finalizeVectorAssembly(VEC_T &)
General function to finalize the vector assembly.
bool finalizeMatrixAssembly(MAT_T &)
void computeFluxCorrectedTransport(NumericalStabilization const &stabilizer, const double t, const double dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, const MathLib::MatrixSpecifications &matrix_specification, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b)
NumLib::LocalToGlobalIndexMap::RowColumnIndices getRowColumnIndices(std::size_t const id, NumLib::LocalToGlobalIndexMap const &dof_table, std::vector< GlobalIndexType > &indices)
GlobalVector computeResiduum(double const dt, GlobalVector const &x, GlobalVector const &x_prev, GlobalMatrix const &M, GlobalMatrix const &K, GlobalVector const &b)
std::vector< double > getCoupledLocalSolutions(std::vector< GlobalVector * > const &global_solutions, std::vector< std::vector< GlobalIndexType > > const &indices)
SecondaryVariableFunctions makeExtrapolator(const unsigned num_components, NumLib::Extrapolator &extrapolator, LocalAssemblerCollection const &local_assemblers, typename NumLib::ExtrapolatableLocalAssemblerCollection< LocalAssemblerCollection >::IntegrationPointValuesMethod integration_point_values_method)
static NUMLIB_EXPORT MatrixProvider & provider
static NUMLIB_EXPORT VectorProvider & provider
static void executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)
static void executeSelectedMemberDereferenced(Object &object, Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)
void assemble(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix *const M, GlobalMatrix *const K, GlobalVector *const b, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, VectorMatrixAssembler &global_assembler, VectorOfLocalAssemblers const &local_assemblers, std::vector< std::size_t > const &active_element_ids)