OGS
ParallelVectorMatrixAssembler.cpp
Go to the documentation of this file.
1
12
13#include <cstdlib>
14#include <fstream>
15#include <range/v3/view/iota.hpp>
16#include <vector>
17
19#include "BaseLib/StringTools.h"
26
27namespace
28{
30 const std::size_t mesh_item_id,
32 const NumLib::LocalToGlobalIndexMap& dof_table, const double t,
33 const double dt, const GlobalVector& x, const GlobalVector& x_prev,
34 std::vector<double>& local_b_data, std::vector<double>& local_Jac_data,
35 ProcessLib::AbstractJacobianAssembler& jacobian_assembler,
37{
38 std::vector<GlobalIndexType> const& indices =
39 NumLib::getIndices(mesh_item_id, dof_table);
40
41 local_b_data.clear();
42 local_Jac_data.clear();
43
44 auto const local_x = x.get(indices);
45 auto const local_x_prev = x_prev.get(indices);
46 jacobian_assembler.assembleWithJacobian(local_assembler, t, dt, local_x,
47 local_x_prev, local_b_data,
48 local_Jac_data);
49
50 if (local_Jac_data.empty())
51 {
53 "No Jacobian has been assembled! This might be due to "
54 "programming errors in the local assembler of the "
55 "current process.");
56 }
57
58 cache.add(local_b_data, local_Jac_data, indices);
59}
60
62 const std::size_t mesh_item_id,
64 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
65 const double t, const double dt, std::vector<GlobalVector*> const& x,
66 std::vector<GlobalVector*> const& x_prev, int const process_id,
67 std::vector<double>& local_b_data, std::vector<double>& local_Jac_data,
68 ProcessLib::AbstractJacobianAssembler& jacobian_assembler,
70{
71 std::vector<std::vector<GlobalIndexType>> indices_of_processes;
72 indices_of_processes.reserve(dof_tables.size());
73 transform(cbegin(dof_tables), cend(dof_tables),
74 back_inserter(indices_of_processes),
75 [&](auto const* dof_table)
76 { return NumLib::getIndices(mesh_item_id, *dof_table); });
77
78 auto local_coupled_xs =
79 ProcessLib::getCoupledLocalSolutions(x, indices_of_processes);
80 auto const local_x = MathLib::toVector(local_coupled_xs);
81
82 auto local_coupled_x_prevs =
83 ProcessLib::getCoupledLocalSolutions(x_prev, indices_of_processes);
84 auto const local_x_prev = MathLib::toVector(local_coupled_x_prevs);
85
86 std::vector<GlobalIndexType> const& indices =
87 indices_of_processes[process_id];
88
89 local_b_data.clear();
90 local_Jac_data.clear();
91
93 local_assembler, t, dt, local_x, local_x_prev, process_id, local_b_data,
94 local_Jac_data);
95
96 if (local_Jac_data.empty())
97 {
99 "No Jacobian has been assembled! This might be due to "
100 "programming errors in the local assembler of the "
101 "current process.");
102 }
103
104 cache.add(local_b_data, local_Jac_data, indices);
105}
106
108std::vector<
109 std::tuple<std::ptrdiff_t,
110 std::reference_wrapper<ProcessLib::LocalAssemblerInterface>>>
113 ProcessLib::LocalAssemblerInterface> const& local_assemblers,
114 std::vector<std::size_t> const& active_elements)
115{
116 auto create_ids_asm_pairs = [&](auto const& element_ids)
117 {
118 std::vector<std::tuple<
119 std::ptrdiff_t,
120 std::reference_wrapper<ProcessLib::LocalAssemblerInterface>>>
121 result;
122 result.reserve(static_cast<std::size_t>(element_ids.size()));
123 for (auto const id : element_ids)
124 {
125 result.push_back({id, local_assemblers[id]});
126 }
127 return result;
128 };
129
130 if (active_elements.empty())
131 {
132 return create_ids_asm_pairs(ranges::views::iota(
133 static_cast<std::size_t>(0), local_assemblers.size()));
134 }
135 return create_ids_asm_pairs(active_elements);
136}
137
139 std::vector<std::tuple<
140 std::ptrdiff_t,
141 std::reference_wrapper<ProcessLib::LocalAssemblerInterface>>> const&
142 ids_local_assemblers,
143 ThreadException& exception,
144 auto local_matrix_output,
145 auto assemble)
146{
147 std::ptrdiff_t n_elements =
148 static_cast<std::ptrdiff_t>(ids_local_assemblers.size());
149
150#pragma omp for nowait
151 for (std::ptrdiff_t i = 0; i < n_elements; ++i)
152 {
153 if (exception)
154 {
155 continue;
156 }
157 auto [element_id, loc_asm] = ids_local_assemblers[i];
158
159 try
160 {
161 assemble(element_id, loc_asm);
162 }
163 catch (...)
164 {
165 exception.capture();
166 continue;
167 }
168
169 local_matrix_output(element_id);
170 }
171}
172
173} // namespace
174
175namespace ProcessLib::Assembly
176{
178 AbstractJacobianAssembler& jacobian_assembler)
179 : jacobian_assembler_{jacobian_assembler},
180 num_threads_(BaseLib::getNumberOfAssemblyThreads())
181{
182 INFO("Threads used for ParallelVectorMatrixAssembler: {}.", num_threads_);
183}
184
187 LocalAssemblerInterface> const& local_assemblers,
188 std::vector<std::size_t> const& active_elements,
189 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
190 const double t, double const dt, std::vector<GlobalVector*> const& xs,
191 std::vector<GlobalVector*> const& x_prevs, int const process_id,
192 GlobalVector& b, GlobalMatrix& Jac)
193{
194 // checks //////////////////////////////////////////////////////////////////
195 if (dof_tables.size() != xs.size())
196 {
197 OGS_FATAL("Different number of DOF tables and solution vectors.");
198 }
199
200 std::size_t const number_of_processes = xs.size();
201 // algorithm ///////////////////////////////////////////////////////////////
202
204 ConcurrentMatrixView b_view(b);
205 ConcurrentMatrixView Jac_view(Jac);
206
207 ThreadException exception;
208#pragma omp parallel num_threads(num_threads_)
209 {
210#ifdef _OPENMP
211#pragma omp single nowait
212 {
213 INFO("Number of threads: {}", omp_get_num_threads());
214 }
215#endif
216
217 // temporary data only stored here in order to avoid frequent memory
218 // reallocations.
219 std::vector<double> local_b_data;
220 std::vector<double> local_Jac_data;
221
222 // copy to avoid concurrent access
223 auto const jac_asm = jacobian_assembler_.copy();
224 auto stats_this_thread = stats->clone();
225
226 MultiMatrixElementCache cache{b_view, Jac_view,
227 stats_this_thread->data};
228
229 auto local_matrix_output = [&](std::ptrdiff_t element_id)
230 {
231 local_matrix_output_(t, process_id, element_id, local_b_data,
232 local_Jac_data);
233 };
234
235 // TODO corner case: what if all elements on a submesh are deactivated?
236
237 // Monolithic scheme
238 if (number_of_processes == 1)
239 {
240 assert(process_id == 0);
241 auto const& dof_table = *dof_tables[0];
242 auto const& x = *xs[0];
243 auto const& x_prev = *x_prevs[0];
244
245 runAssembly(
246 collectActiveLocalAssemblers(local_assemblers, active_elements),
247 exception, local_matrix_output,
248 [&](auto element_id, auto& loc_asm)
249 {
250 assembleWithJacobianOneElement(
251 element_id, loc_asm, dof_table, t, dt, x, x_prev,
252 local_b_data, local_Jac_data, *jac_asm, cache);
253 });
254 }
255 else // Staggered scheme
256 {
257 runAssembly(
258 collectActiveLocalAssemblers(local_assemblers, active_elements),
259 exception, local_matrix_output,
260 [&](auto element_id, auto& loc_asm)
261 {
262 assembleWithJacobianForStaggeredSchemeOneElement(
263 element_id, loc_asm, dof_tables, t, dt, xs, x_prevs,
264 process_id, local_b_data, local_Jac_data, *jac_asm,
265 cache);
266 });
267 }
268 }
269
270 stats->print();
271
272 global_matrix_output_(t, process_id, b, Jac);
273 exception.rethrow();
274}
275} // namespace ProcessLib::Assembly
#define OGS_FATAL(...)
Definition Error.h:26
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:36
Definition of string helper functions.
Global vector based on Eigen vector.
Definition EigenVector.h:26
double get(IndexType rowId) const
get entry
Definition EigenVector.h:59
Base class for Jacobian assemblers.
virtual 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)=0
virtual std::unique_ptr< AbstractJacobianAssembler > copy() const =0
virtual void assembleWithJacobianForStaggeredScheme(LocalAssemblerInterface &, double const, double const, Eigen::VectorXd const &, Eigen::VectorXd const &, int const, std::vector< double > &, std::vector< double > &)
static std::shared_ptr< CumulativeStats< Data > > create()
void add(std::vector< double > const &local_b_data, std::vector< double > const &local_Jac_data, std::vector< GlobalIndexType > const &indices)
ParallelVectorMatrixAssembler(AbstractJacobianAssembler &jacobian_assembler)
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)
Eigen::Map< const Vector > toVector(std::vector< double > const &data, Eigen::VectorXd::Index size)
Creates an Eigen mapped vector from the given data vector.
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
std::vector< double > getCoupledLocalSolutions(std::vector< GlobalVector * > const &global_solutions, std::vector< std::vector< GlobalIndexType > > const &indices)
std::vector< std::tuple< std::ptrdiff_t, std::reference_wrapper< ProcessLib::LocalAssemblerInterface > > > collectActiveLocalAssemblers(BaseLib::PolymorphicRandomAccessContainerView< ProcessLib::LocalAssemblerInterface > const &local_assemblers, std::vector< std::size_t > const &active_elements)
Returns a vector of active element ids with corresponding local assemblers.
void assembleWithJacobianOneElement(const std::size_t mesh_item_id, ProcessLib::LocalAssemblerInterface &local_assembler, const NumLib::LocalToGlobalIndexMap &dof_table, const double t, const double dt, const GlobalVector &x, const GlobalVector &x_prev, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data, ProcessLib::AbstractJacobianAssembler &jacobian_assembler, ProcessLib::Assembly::MultiMatrixElementCache &cache)
void assembleWithJacobianForStaggeredSchemeOneElement(const std::size_t mesh_item_id, ProcessLib::LocalAssemblerInterface &local_assembler, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, const double t, const double dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data, ProcessLib::AbstractJacobianAssembler &jacobian_assembler, ProcessLib::Assembly::MultiMatrixElementCache &cache)
void runAssembly(std::vector< std::tuple< std::ptrdiff_t, std::reference_wrapper< ProcessLib::LocalAssemblerInterface > > > const &ids_local_assemblers, ThreadException &exception, auto local_matrix_output, auto assemble)