OGS
anonymous_namespace{ParallelVectorMatrixAssembler.cpp} Namespace Reference

Functions

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)
 
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 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)
 
int getNumberOfThreads ()
 

Function Documentation

◆ 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 & cache )

Definition at line 60 of file ParallelVectorMatrixAssembler.cpp.

69{
70 std::vector<std::vector<GlobalIndexType>> indices_of_processes;
71 indices_of_processes.reserve(dof_tables.size());
72 transform(cbegin(dof_tables), cend(dof_tables),
73 back_inserter(indices_of_processes),
74 [&](auto const* dof_table)
75 { return NumLib::getIndices(mesh_item_id, *dof_table); });
76
77 auto local_coupled_xs =
78 ProcessLib::getCoupledLocalSolutions(x, indices_of_processes);
79 auto const local_x = MathLib::toVector(local_coupled_xs);
80
81 auto local_coupled_x_prevs =
82 ProcessLib::getCoupledLocalSolutions(x_prev, indices_of_processes);
83 auto const local_x_prev = MathLib::toVector(local_coupled_x_prevs);
84
85 std::vector<GlobalIndexType> const& indices =
86 indices_of_processes[process_id];
87
88 local_b_data.clear();
89 local_Jac_data.clear();
90
92 local_assembler, t, dt, local_x, local_x_prev, process_id, local_b_data,
93 local_Jac_data);
94
95 if (local_Jac_data.empty())
96 {
98 "No Jacobian has been assembled! This might be due to "
99 "programming errors in the local assembler of the "
100 "current process.");
101 }
102
103 cache.add(local_b_data, local_Jac_data, indices);
104}
#define OGS_FATAL(...)
Definition Error.h:26
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)
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::Assembly::MultiMatrixElementCache::add(), 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 & cache )

Definition at line 28 of file ParallelVectorMatrixAssembler.cpp.

36{
37 std::vector<GlobalIndexType> const& indices =
38 NumLib::getIndices(mesh_item_id, dof_table);
39
40 local_b_data.clear();
41 local_Jac_data.clear();
42
43 auto const local_x = x.get(indices);
44 auto const local_x_prev = x_prev.get(indices);
45 jacobian_assembler.assembleWithJacobian(local_assembler, t, dt, local_x,
46 local_x_prev, local_b_data,
47 local_Jac_data);
48
49 if (local_Jac_data.empty())
50 {
52 "No Jacobian has been assembled! This might be due to "
53 "programming errors in the local assembler of the "
54 "current process.");
55 }
56
57 cache.add(local_b_data, local_Jac_data, indices);
58}
double get(IndexType rowId) const
get entry
Definition EigenVector.h:58
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::Assembly::MultiMatrixElementCache::add(), ProcessLib::AbstractJacobianAssembler::assembleWithJacobian(), MathLib::EigenVector::get(), NumLib::getIndices(), and OGS_FATAL.

◆ collectActiveLocalAssemblers()

std::vector< std::tuple< std::ptrdiff_t, std::reference_wrapper< ProcessLib::LocalAssemblerInterface > > > anonymous_namespace{ParallelVectorMatrixAssembler.cpp}::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.

Definition at line 110 of file ParallelVectorMatrixAssembler.cpp.

114{
115 auto create_ids_asm_pairs = [&](auto const& element_ids)
116 {
117 std::vector<std::tuple<
118 std::ptrdiff_t,
119 std::reference_wrapper<ProcessLib::LocalAssemblerInterface>>>
120 result;
121 result.reserve(static_cast<std::size_t>(element_ids.size()));
122 for (auto const id : element_ids)
123 {
124 result.push_back({id, local_assemblers[id]});
125 }
126 return result;
127 };
128
129 if (active_elements.empty())
130 {
131 return create_ids_asm_pairs(ranges::views::iota(
132 static_cast<std::size_t>(0), local_assemblers.size()));
133 }
134 return create_ids_asm_pairs(active_elements);
135}

◆ getNumberOfThreads()

int anonymous_namespace{ParallelVectorMatrixAssembler.cpp}::getNumberOfThreads ( )

Definition at line 172 of file ParallelVectorMatrixAssembler.cpp.

173{
174 char const* const num_threads_env = std::getenv("OGS_ASM_THREADS");
175
176 if (!num_threads_env)
177 {
178 return 1;
179 }
180
181 if (std::strlen(num_threads_env) == 0)
182 {
183 OGS_FATAL("The environment variable OGS_ASM_THREADS is set but empty.");
184 }
185
186 std::string num_threads_str{num_threads_env};
187 BaseLib::trim(num_threads_str);
188
189 std::istringstream num_threads_iss{num_threads_str};
190 int num_threads = -1;
191
192 num_threads_iss >> num_threads;
193
194 if (!num_threads_iss)
195 {
196 OGS_FATAL("Error parsing OGS_ASM_THREADS (= \"{}\").", num_threads_env);
197 }
198
199 if (!num_threads_iss.eof())
200 {
201 OGS_FATAL(
202 "Error parsing OGS_ASM_THREADS (= \"{}\"): not read entirely, the "
203 "remainder is \"{}\"",
204 num_threads_env,
205 num_threads_iss.str().substr(num_threads_iss.tellg()));
206 }
207
208 if (num_threads < 1)
209 {
210 OGS_FATAL(
211 "You asked (via OGS_ASM_THREADS) to assemble with {} < 1 thread.",
212 num_threads);
213 }
214
215 return num_threads;
216}
void trim(std::string &str, char ch)

References OGS_FATAL, and BaseLib::trim().

◆ runAssembly()

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

Definition at line 137 of file ParallelVectorMatrixAssembler.cpp.

145{
146 std::ptrdiff_t n_elements =
147 static_cast<std::ptrdiff_t>(ids_local_assemblers.size());
148
149#pragma omp for nowait
150 for (std::ptrdiff_t i = 0; i < n_elements; ++i)
151 {
152 if (exception)
153 {
154 continue;
155 }
156 auto [element_id, loc_asm] = ids_local_assemblers[i];
157
158 try
159 {
160 assemble(element_id, loc_asm);
161 }
162 catch (...)
163 {
164 exception.capture();
165 continue;
166 }
167
168 local_matrix_output(element_id);
169 }
170}

References ThreadException::capture().