OGS
ProcessLib::Assembly::MatrixElementCache< Dim > Class Template Referencefinal

Detailed Description

template<std::size_t Dim>
class ProcessLib::Assembly::MatrixElementCache< Dim >

Definition at line 120 of file MatrixElementCache.h.

#include <MatrixElementCache.h>

Collaboration diagram for ProcessLib::Assembly::MatrixElementCache< Dim >:
[legend]

Public Member Functions

 MatrixElementCache (GlobalMatView &mat_or_vec, Stats &stats)
 
void add (std::vector< double > const &local_data, std::vector< GlobalIndexType > const &indices)
 
 ~MatrixElementCache ()
 

Private Types

using GlobalMatView = ConcurrentMatrixView<Dim>
 

Private Member Functions

void addToCache (std::vector< double > const &values, std::vector< GlobalIndexType > const &indices)
 
void addToCacheImpl (std::vector< double > const &values, std::vector< GlobalIndexType > const &indices, std::integral_constant< std::size_t, 1 >)
 
void addToCacheImpl (std::vector< double > const &values, std::vector< GlobalIndexType > const &indices, std::integral_constant< std::size_t, 2 >)
 
void ensureEnoughSpace (std::size_t const space_needed)
 
void addToGlobal ()
 

Private Attributes

std::vector< MatrixElementCacheEntry< Dim > > cache_
 
GlobalMatViewmat_or_vec_
 
Statsstats_
 

Static Private Attributes

static constexpr std::size_t cache_capacity = 1'000'000
 

Member Typedef Documentation

◆ GlobalMatView

template<std::size_t Dim>
using ProcessLib::Assembly::MatrixElementCache< Dim >::GlobalMatView = ConcurrentMatrixView<Dim>
private

Definition at line 123 of file MatrixElementCache.h.

Constructor & Destructor Documentation

◆ MatrixElementCache()

template<std::size_t Dim>
ProcessLib::Assembly::MatrixElementCache< Dim >::MatrixElementCache ( GlobalMatView & mat_or_vec,
Stats & stats )
inline

◆ ~MatrixElementCache()

Member Function Documentation

◆ add()

template<std::size_t Dim>
void ProcessLib::Assembly::MatrixElementCache< Dim >::add ( std::vector< double > const & local_data,
std::vector< GlobalIndexType > const & indices )
inline

Definition at line 134 of file MatrixElementCache.h.

136 {
137 addToCache(local_data, indices);
138 }
void addToCache(std::vector< double > const &values, std::vector< GlobalIndexType > const &indices)

References ProcessLib::Assembly::MatrixElementCache< Dim >::addToCache().

Referenced by ProcessLib::Assembly::MultiMatrixElementCache::add().

◆ addToCache()

template<std::size_t Dim>
void ProcessLib::Assembly::MatrixElementCache< Dim >::addToCache ( std::vector< double > const & values,
std::vector< GlobalIndexType > const & indices )
inlineprivate

Definition at line 143 of file MatrixElementCache.h.

145 {
146 if (values.empty())
147 {
148 return;
149 }
150
151 ensureEnoughSpace(values.size());
152
153 addToCacheImpl(values, indices,
154 std::integral_constant<std::size_t, Dim>{});
155 }
void addToCacheImpl(std::vector< double > const &values, std::vector< GlobalIndexType > const &indices, std::integral_constant< std::size_t, 1 >)
void ensureEnoughSpace(std::size_t const space_needed)

References ProcessLib::Assembly::MatrixElementCache< Dim >::addToCacheImpl(), and ProcessLib::Assembly::MatrixElementCache< Dim >::ensureEnoughSpace().

Referenced by ProcessLib::Assembly::MatrixElementCache< Dim >::add().

◆ addToCacheImpl() [1/2]

template<std::size_t Dim>
void ProcessLib::Assembly::MatrixElementCache< Dim >::addToCacheImpl ( std::vector< double > const & values,
std::vector< GlobalIndexType > const & indices,
std::integral_constant< std::size_t, 1 >  )
inlineprivate

Definition at line 158 of file MatrixElementCache.h.

161 {
162 auto const num_r_c = indices.size();
163
164 for (std::size_t r_local = 0; r_local < num_r_c; ++r_local)
165 {
166 ++stats_.count;
167 auto const value = values[r_local];
168
169 if (value == 0)
170 {
171 continue;
172 }
173 else
174 {
176 }
177
178 auto const r_global = indices[r_local];
179 cache_.emplace_back(std::array{r_global}, value);
180 }
181 }

References ProcessLib::Assembly::MatrixElementCache< Dim >::cache_, ProcessLib::Assembly::Stats::count, ProcessLib::Assembly::Stats::count_nonzero, and ProcessLib::Assembly::MatrixElementCache< Dim >::stats_.

Referenced by ProcessLib::Assembly::MatrixElementCache< Dim >::addToCache().

◆ addToCacheImpl() [2/2]

template<std::size_t Dim>
void ProcessLib::Assembly::MatrixElementCache< Dim >::addToCacheImpl ( std::vector< double > const & values,
std::vector< GlobalIndexType > const & indices,
std::integral_constant< std::size_t, 2 >  )
inlineprivate

Definition at line 184 of file MatrixElementCache.h.

187 {
188 auto const num_r_c = indices.size();
189
190 // Note: There is an implicit storage order assumption, here!
191 auto const local_mat = MathLib::toMatrix(values, num_r_c, num_r_c);
192
193 for (std::size_t r_local = 0; r_local < num_r_c; ++r_local)
194 {
195 auto const r_global = indices[r_local];
196
197 for (std::size_t c_local = 0; c_local < num_r_c; ++c_local)
198 {
199 ++stats_.count;
200 auto const value = local_mat(r_local, c_local);
201
202 // TODO skipping zero values sometimes does not work together
203 // with the Eigen SparseLU linear solver. See also
204 // https://gitlab.opengeosys.org/ogs/ogs/-/merge_requests/4556#note_125561
205#if 0
206 if (value == 0)
207 {
208 continue;
209 }
210#endif
212
213 auto const c_global = indices[c_local];
214 cache_.emplace_back(std::array{r_global, c_global}, value);
215 }
216 }
217 }
Eigen::Map< const Matrix > toMatrix(std::vector< double > const &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)

References ProcessLib::Assembly::MatrixElementCache< Dim >::cache_, ProcessLib::Assembly::Stats::count, ProcessLib::Assembly::Stats::count_nonzero, ProcessLib::Assembly::MatrixElementCache< Dim >::stats_, and MathLib::toMatrix().

◆ addToGlobal()

◆ ensureEnoughSpace()

template<std::size_t Dim>
void ProcessLib::Assembly::MatrixElementCache< Dim >::ensureEnoughSpace ( std::size_t const space_needed)
inlineprivate

Definition at line 219 of file MatrixElementCache.h.

220 {
221 auto const size_initial = cache_.size();
222 auto const cap_initial = cache_.capacity();
223
224 if (size_initial + space_needed <= cap_initial)
225 {
226 return;
227 }
228
229 addToGlobal();
230
231 // ensure again that there is enough capacity (corner case, if initial
232 // capacity is too small because of whatever arcane reason)
233 auto const size_after = cache_.size();
234 auto const cap_after = cache_.capacity();
235
236 if (size_after + space_needed > cap_after)
237 {
238 cache_.reserve(size_after + 2 * space_needed);
239 }
240 }

References ProcessLib::Assembly::MatrixElementCache< Dim >::addToGlobal(), and ProcessLib::Assembly::MatrixElementCache< Dim >::cache_.

Referenced by ProcessLib::Assembly::MatrixElementCache< Dim >::addToCache().

Member Data Documentation

◆ cache_

◆ cache_capacity

template<std::size_t Dim>
constexpr std::size_t ProcessLib::Assembly::MatrixElementCache< Dim >::cache_capacity = 1'000'000
staticconstexprprivate

◆ mat_or_vec_

template<std::size_t Dim>
GlobalMatView& ProcessLib::Assembly::MatrixElementCache< Dim >::mat_or_vec_
private

◆ stats_


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