OGS
anonymous_namespace{ParallelVectorMatrixAssembler.cpp} Namespace Reference

Functions

void assembleOneElement (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_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data, ProcessLib::Assembly::MultiMatrixElementCache< 2 > &cache)
void assembleForStaggeredSchemeOneElement (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_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data, ProcessLib::Assembly::MultiMatrixElementCache< 2 > &cache)
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< 1 > &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< 1 > &cache)
void runAssemblyImpl (BaseLib::PolymorphicRandomAccessContainerView< ProcessLib::LocalAssemblerInterface > const &local_assemblers, ThreadException &exception, auto local_matrix_output, auto assemble, std::ptrdiff_t const num_active_local_assemblers, auto get_element_id_callback)
 Runs the passed assemble functor for each active local assembler.
void runAssembly (BaseLib::PolymorphicRandomAccessContainerView< ProcessLib::LocalAssemblerInterface > const &local_assemblers, std::vector< std::size_t > const *const active_elements, ThreadException &exception, auto local_matrix_output, auto assemble)

Function Documentation

◆ assembleForStaggeredSchemeOneElement()

void anonymous_namespace{ParallelVectorMatrixAssembler.cpp}::assembleForStaggeredSchemeOneElement ( 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_M_data,
std::vector< double > & local_K_data,
std::vector< double > & local_b_data,
ProcessLib::Assembly::MultiMatrixElementCache< 2 > & cache )

Definition at line 48 of file ParallelVectorMatrixAssembler.cpp.

57{
58 std::vector<std::vector<GlobalIndexType>> indices_of_processes;
59 indices_of_processes.reserve(dof_tables.size());
60 transform(cbegin(dof_tables), cend(dof_tables),
61 back_inserter(indices_of_processes), [&](auto const* dof_table)
62 { return NumLib::getIndices(mesh_item_id, *dof_table); });
63
64 auto local_coupled_xs =
65 ProcessLib::getCoupledLocalSolutions(x, indices_of_processes);
66 auto const local_x = MathLib::toVector(local_coupled_xs);
67
68 auto local_coupled_x_prevs =
69 ProcessLib::getCoupledLocalSolutions(x_prev, indices_of_processes);
70 auto const local_x_prev = MathLib::toVector(local_coupled_x_prevs);
71
72 local_M_data.clear();
73 local_K_data.clear();
74 local_b_data.clear();
75
76 local_assembler.assembleForStaggeredScheme(t, dt, local_x, local_x_prev,
77 process_id, local_M_data,
78 local_K_data, local_b_data);
79
80 auto const& indices = indices_of_processes[process_id];
81 cache.add(local_M_data, local_K_data, local_b_data, indices);
82}
void add(std::vector< double > const &local_M_data, std::vector< double > const &local_K_data, std::vector< double > const &local_b_data, std::vector< GlobalIndexType > const &indices)
virtual void assembleForStaggeredScheme(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, int const process_id, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data)
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)

References ProcessLib::LocalAssemblerInterface::assembleForStaggeredScheme(), ProcessLib::getCoupledLocalSolutions(), NumLib::getIndices(), and MathLib::toVector().

◆ assembleOneElement()

void anonymous_namespace{ParallelVectorMatrixAssembler.cpp}::assembleOneElement ( 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_M_data,
std::vector< double > & local_K_data,
std::vector< double > & local_b_data,
ProcessLib::Assembly::MultiMatrixElementCache< 2 > & cache )

Definition at line 23 of file ParallelVectorMatrixAssembler.cpp.

32{
33 std::vector<GlobalIndexType> indices =
34 NumLib::getIndices(mesh_item_id, dof_table);
35
36 local_M_data.clear();
37 local_K_data.clear();
38 local_b_data.clear();
39
40 auto const local_x = x.get(indices);
41 auto const local_x_prev = x_prev.get(indices);
42 local_assembler.assemble(t, dt, local_x, local_x_prev, local_M_data,
43 local_K_data, local_b_data);
44
45 cache.add(local_M_data, local_K_data, local_b_data, std::move(indices));
46}
double get(IndexType rowId) const
get entry
Definition EigenVector.h:52
virtual void assemble(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data)

References ProcessLib::LocalAssemblerInterface::assemble(), MathLib::EigenVector::get(), and NumLib::getIndices().

◆ assembleWithJacobianForStaggeredSchemeOneElement()

void anonymous_namespace{ParallelVectorMatrixAssembler.cpp}::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< 1 > & cache )

Definition at line 116 of file ParallelVectorMatrixAssembler.cpp.

125{
126 std::vector<std::vector<GlobalIndexType>> indices_of_processes;
127 indices_of_processes.reserve(dof_tables.size());
128 transform(cbegin(dof_tables), cend(dof_tables),
129 back_inserter(indices_of_processes),
130 [&](auto const* dof_table)
131 { return NumLib::getIndices(mesh_item_id, *dof_table); });
132
133 auto local_coupled_xs =
134 ProcessLib::getCoupledLocalSolutions(x, indices_of_processes);
135 auto const local_x = MathLib::toVector(local_coupled_xs);
136
137 auto local_coupled_x_prevs =
138 ProcessLib::getCoupledLocalSolutions(x_prev, indices_of_processes);
139 auto const local_x_prev = MathLib::toVector(local_coupled_x_prevs);
140
141 std::vector<GlobalIndexType> const& indices =
142 indices_of_processes[process_id];
143
144 local_b_data.clear();
145 local_Jac_data.clear();
146
148 local_assembler, t, dt, local_x, local_x_prev, process_id, local_b_data,
149 local_Jac_data);
150
151 if (local_Jac_data.empty())
152 {
153 OGS_FATAL(
154 "No Jacobian has been assembled! This might be due to "
155 "programming errors in the local assembler of the "
156 "current process.");
157 }
158
159 cache.add(local_b_data, local_Jac_data, indices);
160}
#define OGS_FATAL(...)
Definition Error.h:19
virtual void assembleWithJacobianForStaggeredScheme(LocalAssemblerInterface &, double const, double const, Eigen::VectorXd const &, Eigen::VectorXd const &, int const, std::vector< double > &, std::vector< double > &)
void add(std::vector< double > const &local_b_data, std::vector< double > const &local_Jac_data, std::vector< GlobalIndexType > const &indices)

References ProcessLib::AbstractJacobianAssembler::assembleWithJacobianForStaggeredScheme(), ProcessLib::getCoupledLocalSolutions(), NumLib::getIndices(), OGS_FATAL, and MathLib::toVector().

◆ assembleWithJacobianOneElement()

void anonymous_namespace{ParallelVectorMatrixAssembler.cpp}::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< 1 > & cache )

Definition at line 84 of file ParallelVectorMatrixAssembler.cpp.

92{
93 std::vector<GlobalIndexType> const& indices =
94 NumLib::getIndices(mesh_item_id, dof_table);
95
96 local_b_data.clear();
97 local_Jac_data.clear();
98
99 auto const local_x = x.get(indices);
100 auto const local_x_prev = x_prev.get(indices);
101 jacobian_assembler.assembleWithJacobian(local_assembler, t, dt, local_x,
102 local_x_prev, local_b_data,
103 local_Jac_data);
104
105 if (local_Jac_data.empty())
106 {
107 OGS_FATAL(
108 "No Jacobian has been assembled! This might be due to "
109 "programming errors in the local assembler of the "
110 "current process.");
111 }
112
113 cache.add(local_b_data, local_Jac_data, indices);
114}
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

References ProcessLib::AbstractJacobianAssembler::assembleWithJacobian(), MathLib::EigenVector::get(), NumLib::getIndices(), and OGS_FATAL.

◆ runAssembly()

void anonymous_namespace{ParallelVectorMatrixAssembler.cpp}::runAssembly ( BaseLib::PolymorphicRandomAccessContainerView< ProcessLib::LocalAssemblerInterface > const & local_assemblers,
std::vector< std::size_t > const *const active_elements,
ThreadException & exception,
auto local_matrix_output,
auto assemble )

Runs the passed assemble functor for each active local assembler.

Forwards all calls to runAssemblyImpl().

Definition at line 200 of file ParallelVectorMatrixAssembler.cpp.

207{
208 [[likely]] if (!active_elements)
209 {
210 std::ptrdiff_t const n_elements =
211 static_cast<std::ptrdiff_t>(local_assemblers.size());
212
213 runAssemblyImpl(local_assemblers, exception, local_matrix_output,
214 assemble, n_elements, [](auto i) { return i; });
215 }
216 else
217 {
218 std::ptrdiff_t const n_elements =
219 static_cast<std::ptrdiff_t>(active_elements->size());
220
221 runAssemblyImpl(local_assemblers, exception, local_matrix_output,
222 assemble, n_elements, [active_elements](auto i)
223 { return (*active_elements)[i]; });
224 }
225}
void runAssemblyImpl(BaseLib::PolymorphicRandomAccessContainerView< ProcessLib::LocalAssemblerInterface > const &local_assemblers, ThreadException &exception, auto local_matrix_output, auto assemble, std::ptrdiff_t const num_active_local_assemblers, auto get_element_id_callback)
Runs the passed assemble functor for each active local assembler.

References runAssemblyImpl().

◆ runAssemblyImpl()

void anonymous_namespace{ParallelVectorMatrixAssembler.cpp}::runAssemblyImpl ( BaseLib::PolymorphicRandomAccessContainerView< ProcessLib::LocalAssemblerInterface > const & local_assemblers,
ThreadException & exception,
auto local_matrix_output,
auto assemble,
std::ptrdiff_t const num_active_local_assemblers,
auto get_element_id_callback )

Runs the passed assemble functor for each active local assembler.

Definition at line 163 of file ParallelVectorMatrixAssembler.cpp.

171{
172#pragma omp for nowait
173 for (std::ptrdiff_t i = 0; i < num_active_local_assemblers; ++i)
174 {
175 [[unlikely]] if (exception)
176 {
177 continue;
178 }
179
180 auto const element_id = get_element_id_callback(i);
181 auto& loc_asm = local_assemblers[element_id];
182
183 try
184 {
185 assemble(element_id, loc_asm);
186 }
187 catch (...)
188 {
189 exception.capture();
190 continue;
191 }
192
193 local_matrix_output(element_id);
194 }
195}

References ThreadException::capture().

Referenced by runAssembly().