OGS
AssemblyMixin.h
Go to the documentation of this file.
1
11#pragma once
12
13#include "BaseLib/MPI.h"
16#include "NumLib/Exceptions.h"
20
21namespace ProcessLib
22{
26{
27 explicit SubmeshAssemblyData(
28 MeshLib::Mesh const& mesh,
29 std::vector<std::vector<
30 std::reference_wrapper<MeshLib::PropertyVector<double>>>>&&
32
35 std::vector<std::size_t> active_element_ids;
36 std::vector<
37 std::vector<std::reference_wrapper<MeshLib::PropertyVector<double>>>>
39};
40
44{
51
52protected:
53 explicit AssemblyMixinBase(AbstractJacobianAssembler& jacobian_assembler)
54 : pvma_{jacobian_assembler} {};
55
57 MeshLib::Mesh& bulk_mesh,
58 std::vector<std::reference_wrapper<MeshLib::Mesh>> const& submeshes,
59 std::vector<std::vector<std::string>> const& residuum_names,
60 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>> const&
61 pvs);
62
63 void updateActiveElements(ProcessLib::Process const& process);
64
66 GlobalVector const& rhs,
67 NumLib::LocalToGlobalIndexMap const& local_to_global_index_map,
68 std::vector<std::reference_wrapper<MeshLib::PropertyVector<double>>>
69 residuum_vectors);
70
72 int const process_id,
73 GlobalVector const& rhs,
74 NumLib::LocalToGlobalIndexMap const& local_to_global_index_map,
75 SubmeshAssemblyData const& sad);
76
77private:
78 void updateActiveElementsImpl(Process const& process);
79
80protected:
81 std::vector<SubmeshAssemblyData> submesh_assembly_data_;
82 std::vector<
83 std::vector<std::reference_wrapper<MeshLib::PropertyVector<double>>>>
85
87 std::size_t b_submesh_id_ = 0;
88
90
91private:
93};
94
100template <typename Process>
102{
103 // Enforce correct use of CRTP, i.e., that Process is derived from
104 // AssemblyMixin<Process>.
106 friend Process;
107
108public:
131 std::vector<std::reference_wrapper<MeshLib::Mesh>> const& submeshes,
132 std::vector<std::vector<std::string>> const& residuum_names)
133 {
135 derived().getMesh(), submeshes, residuum_names,
136 derived().getProcessVariables());
137 }
138
143
144 // cppcheck-suppress functionStatic
145 void assemble(double const t, double const dt,
146 std::vector<GlobalVector*> const& /*x*/,
147 std::vector<GlobalVector*> const& /*x_prev*/,
148 int const process_id, GlobalMatrix& /*M*/,
149 GlobalMatrix& /*K*/, GlobalVector& /*b*/)
150 {
151 DBUG("AssemblyMixin assemble(t={}, dt={}, process_id={}).", t, dt,
152 process_id);
153
157 OGS_FATAL("AssemblyMixin for Picard scheme is not yet implemented.");
158 }
159
160 void assembleWithJacobian(double const t, double const dt,
161 std::vector<GlobalVector*> const& x,
162 std::vector<GlobalVector*> const& x_prev,
163 int const process_id, GlobalVector& b,
164 GlobalMatrix& Jac)
165 {
166 DBUG("AssemblyMixin assembleWithJacobian(t={}, dt={}, process_id={}).",
167 t, dt, process_id);
168
169 std::vector<NumLib::LocalToGlobalIndexMap const*> const dof_tables =
170 derived().getDOFTables(x.size());
171
172 auto const& loc_asms = derived().local_assemblers_;
173
174 std::exception_ptr exception = nullptr;
175 if (!submesh_assembly_data_.empty())
176 {
178 b, b_submesh_id_);
179
180 for (auto const& sad : submesh_assembly_data_)
181 {
182 b_submesh.setZero();
183
184 try
185 {
186 pvma_.assembleWithJacobian(loc_asms, sad.active_element_ids,
187 dof_tables, t, dt, x, x_prev,
188 process_id, b_submesh, Jac);
189 }
190 catch (NumLib::AssemblyException const&)
191 {
192 exception = std::current_exception();
193 }
194
195 MathLib::LinAlg::axpy(b, 1.0, b_submesh);
196
198 process_id, b_submesh, *(dof_tables[process_id]), sad);
199 }
200
202 }
203 else
204 {
205 try
206 {
208 loc_asms, derived().getActiveElementIDs(), dof_tables, t,
209 dt, x, x_prev, process_id, b, Jac);
210 }
211 catch (NumLib::AssemblyException const&)
212 {
213 exception = std::current_exception();
214 }
215 }
216
219
220 if (BaseLib::MPI::anyOf(exception != nullptr))
221 {
222 if (exception) // Only the rank with the exception rethrows...
223 {
224 std::rethrow_exception(exception);
225 }
226 // ... but all ranks quit.
227 return;
228 }
229
231 b, *(dof_tables[process_id]), residuum_vectors_bulk_[process_id]);
232 }
233
234private:
235 Process& derived() { return static_cast<Process&>(*this); }
236 Process const& derived() const
237 {
238 return static_cast<Process const&>(*this);
239 }
240};
241
242} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:26
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
Global vector based on Eigen vector.
Definition EigenVector.h:25
virtual GlobalVector & getVector(std::size_t &id)=0
Get an uninitialized vector with the given id.
virtual void releaseVector(GlobalVector const &x)=0
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)
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)
void updateActiveElementsImpl(Process const &process)
std::size_t b_submesh_id_
ID of the b vector on submeshes, cf. NumLib::VectorProvider.
ActiveElementIDsState ids_state_
void updateActiveElements(ProcessLib::Process const &process)
static void copyResiduumVectorsToSubmesh(int const process_id, GlobalVector const &rhs, NumLib::LocalToGlobalIndexMap const &local_to_global_index_map, SubmeshAssemblyData const &sad)
std::vector< SubmeshAssemblyData > submesh_assembly_data_
std::vector< std::vector< std::reference_wrapper< MeshLib::PropertyVector< double > > > > residuum_vectors_bulk_
void assemble(double const t, double const dt, std::vector< GlobalVector * > const &, std::vector< GlobalVector * > const &, int const process_id, GlobalMatrix &, GlobalMatrix &, GlobalVector &)
Process const & derived() 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)
void assembleWithJacobian(BaseLib::PolymorphicRandomAccessContainerView< LocalAssemblerInterface > const &local_assemblers, std::vector< std::size_t > const &active_elements, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, const double t, double const dt, std::vector< GlobalVector * > const &xs, std::vector< GlobalVector * > const &x_prevs, int const process_id, GlobalVector &b, GlobalMatrix &Jac)
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables(int const number_of_processes) const
Definition Process.cpp:384
static bool anyOf(bool const val, Mpi const &mpi=Mpi{MPI_COMM_WORLD})
Definition MPI.h:187
void finalizeAssembly(PETScMatrix &A)
Definition LinAlg.cpp:198
void axpy(PETScVector &y, PetscScalar const a, PETScVector const &x)
Definition LinAlg.cpp:57
OGSMesh getMesh(std::string const &name)
static NUMLIB_EXPORT VectorProvider & provider
std::vector< std::vector< std::reference_wrapper< MeshLib::PropertyVector< double > > > > residuum_vectors
MeshLib::PropertyVector< std::size_t > const & bulk_node_ids
MeshLib::PropertyVector< std::size_t > const & bulk_element_ids
std::vector< std::size_t > active_element_ids
SubmeshAssemblyData(MeshLib::Mesh const &mesh, std::vector< std::vector< std::reference_wrapper< MeshLib::PropertyVector< double > > > > &&residuum_vectors)