OGS
ProcessLib::Assembly::ParallelVectorMatrixAssembler Class Reference

Detailed Description

Definition at line 20 of file ParallelVectorMatrixAssembler.h.

#include <ParallelVectorMatrixAssembler.h>

Collaboration diagram for ProcessLib::Assembly::ParallelVectorMatrixAssembler:
[legend]

Public Member Functions

 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)

Private Attributes

AbstractJacobianAssemblerjacobian_assembler_
LocalMatrixOutput local_matrix_output_
GlobalMatrixOutput global_matrix_output_
int const num_threads_

Constructor & Destructor Documentation

◆ ParallelVectorMatrixAssembler()

ProcessLib::Assembly::ParallelVectorMatrixAssembler::ParallelVectorMatrixAssembler ( AbstractJacobianAssembler & jacobian_assembler)
explicit

Definition at line 177 of file ParallelVectorMatrixAssembler.cpp.

179 : jacobian_assembler_{jacobian_assembler},
181{
182 INFO("Threads used for ParallelVectorMatrixAssembler: {}.", num_threads_);
183}
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:36
int getNumberOfAssemblyThreads()

References INFO(), jacobian_assembler_, and num_threads_.

Member Function Documentation

◆ assembleWithJacobian()

void ProcessLib::Assembly::ParallelVectorMatrixAssembler::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 )

Definition at line 185 of file ParallelVectorMatrixAssembler.cpp.

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
246 collectActiveLocalAssemblers(local_assemblers, active_elements),
247 exception, local_matrix_output,
248 [&](auto element_id, auto& loc_asm)
249 {
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 {
258 collectActiveLocalAssemblers(local_assemblers, active_elements),
259 exception, local_matrix_output,
260 [&](auto element_id, auto& loc_asm)
261 {
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}
#define OGS_FATAL(...)
Definition Error.h:26
static std::shared_ptr< CumulativeStats< Data > > create()
ConcurrentMatrixView(GlobalVector &) -> ConcurrentMatrixView< 1 >
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)

References ProcessLib::Assembly::CumulativeStats< Data >::create(), global_matrix_output_, INFO(), jacobian_assembler_, local_matrix_output_, OGS_FATAL, and ThreadException::rethrow().

Member Data Documentation

◆ global_matrix_output_

GlobalMatrixOutput ProcessLib::Assembly::ParallelVectorMatrixAssembler::global_matrix_output_
private

Definition at line 38 of file ParallelVectorMatrixAssembler.h.

Referenced by assembleWithJacobian().

◆ jacobian_assembler_

AbstractJacobianAssembler& ProcessLib::Assembly::ParallelVectorMatrixAssembler::jacobian_assembler_
private

◆ local_matrix_output_

LocalMatrixOutput ProcessLib::Assembly::ParallelVectorMatrixAssembler::local_matrix_output_
private

Definition at line 37 of file ParallelVectorMatrixAssembler.h.

Referenced by assembleWithJacobian().

◆ num_threads_

int const ProcessLib::Assembly::ParallelVectorMatrixAssembler::num_threads_
private

Definition at line 40 of file ParallelVectorMatrixAssembler.h.

Referenced by ParallelVectorMatrixAssembler().


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