PolyFEM
Loading...
Searching...
No Matches
MatrixCache.cpp
Go to the documentation of this file.
1#include "MatrixCache.hpp"
2
5
6namespace polyfem::utils
7{
9 {
10 init(size);
11 }
12
13 SparseMatrixCache::SparseMatrixCache(const size_t rows, const size_t cols)
14 {
15 init(rows, cols);
16 }
17
19 {
20 init(other);
21 }
22
24 const SparseMatrixCache &other, const bool copy_main_cache_ptr)
25 {
26 init(other, copy_main_cache_ptr);
27 }
28
29 void SparseMatrixCache::init(const size_t size)
30 {
31 assert(mapping().empty() || size_ == size);
32
33 size_ = size;
34 tmp_.resize(size_, size_);
35 mat_.resize(size_, size_);
36 mat_.setZero();
37 }
38
39 void SparseMatrixCache::init(const size_t rows, const size_t cols)
40 {
41 assert(mapping().empty());
42
43 size_ = rows == cols ? rows : 0;
44 tmp_.resize(rows, cols);
45 mat_.resize(rows, cols);
46 mat_.setZero();
47 }
48
50 {
51 assert(this != &other);
52 assert(&other == &dynamic_cast<const SparseMatrixCache &>(other));
53 init(dynamic_cast<const SparseMatrixCache &>(other));
54 }
55
57 const SparseMatrixCache &other, const bool copy_main_cache_ptr)
58 {
59 assert(this != &other);
60 if (copy_main_cache_ptr)
61 {
63 }
64 else if (main_cache_ == nullptr)
65 {
66 main_cache_ = other.main_cache();
67 // Only one level of cache
68 assert(main_cache_ != this && main_cache_ != nullptr && main_cache_->main_cache_ == nullptr);
69 }
70 size_ = other.size_;
71
72 values_.resize(other.values_.size());
73
74 tmp_.resize(other.mat_.rows(), other.mat_.cols());
75 mat_.resize(other.mat_.rows(), other.mat_.cols());
76 mat_.setZero();
77 std::fill(values_.begin(), values_.end(), 0);
78 }
79
81 {
82 tmp_.setZero();
83 mat_.setZero();
84
85 std::fill(values_.begin(), values_.end(), 0);
86 }
87
88 void SparseMatrixCache::add_value(const int e, const int i, const int j, const double value)
89 {
90 // caches have yet to be constructed (likely because the matrix has yet to be fully assembled)
91 if (mapping().empty())
92 {
93 // save entry so it can be added to the matrix later
94 entries_.emplace_back(i, j, value);
95
96 // save the index information so the cache can be built later
97 if (second_cache_entries_.size() <= e)
98 second_cache_entries_.resize(e + 1);
99 second_cache_entries_[e].emplace_back(i, j);
100 }
101 else
102 {
103 if (e != current_e_)
104 {
105 current_e_ = e;
107 }
108
109 // save entry directly to value buffer at the proper index
110 values_[second_cache()[e][current_e_index_]] += value;
112 }
113 }
114
116 {
117 // caches have yet to be constructed (likely because the matrix has yet to be fully assembled)
118 if (mapping().empty())
119 {
120 tmp_.setFromTriplets(entries_.begin(), entries_.end());
121 tmp_.makeCompressed();
122 mat_ += tmp_;
123
124 tmp_.setZero();
125 tmp_.data().squeeze();
126 mat_.makeCompressed();
127
128 entries_.clear();
129
130 mat_.makeCompressed();
131 }
132 }
133
135 {
136 prune();
137
138 // caches have yet to be constructed (likely because the matrix has yet to be fully assembled)
139 if (mapping().empty())
140 {
141 if (compute_mapping && size_ > 0)
142 {
143 assert(main_cache_ == nullptr);
144
145 values_.resize(mat_.nonZeros());
146 inner_index_.resize(mat_.nonZeros());
147 outer_index_.resize(mat_.rows() + 1);
148 mapping_.resize(mat_.rows());
149
150 // note: mat_ is column major
151 const auto inn_ptr = mat_.innerIndexPtr();
152 const auto out_ptr = mat_.outerIndexPtr();
153 inner_index_.assign(inn_ptr, inn_ptr + inner_index_.size());
154 outer_index_.assign(out_ptr, out_ptr + outer_index_.size());
155
156 size_t index = 0;
157 // loop over columns of the matrix
158 for (size_t i = 0; i < mat_.cols(); ++i)
159 {
160 const auto start = outer_index_[i];
161 const auto end = outer_index_[i + 1];
162
163 // loop over the nonzero elements of the given column
164 for (size_t ii = start; ii < end; ++ii)
165 {
166 // pick out current row
167 const auto j = inner_index_[ii];
168 auto &map = mapping_[j];
169 map.emplace_back(i, index);
170 ++index;
171 }
172 }
173
174 logger().trace("Cache computed");
175
176 second_cache_.clear();
178 // loop over each element
179 for (int e = 0; e < second_cache_entries_.size(); ++e)
180 {
181 // loop over each global index affected by the given element
182 for (const auto &p : second_cache_entries_[e])
183 {
184 const int i = p.first;
185 const int j = p.second;
186
187 // pick out column/sparse matrix index pairs for the given column
188 const auto &map = mapping()[i];
189 int index = -1;
190
191 // loop over column/sparse matrix index pairs
192 for (const auto &p : map)
193 {
194 // match columns
195 if (p.first == j)
196 {
197 assert(p.second < values_.size());
198 index = p.second;
199 break;
200 }
201 }
202 assert(index >= 0);
203
204 // save the sparse matrix index used by this element
205 second_cache_[e].emplace_back(index);
206 }
207 }
208
209 second_cache_entries_.resize(0);
210
211 logger().trace("Second cache computed");
212 }
213 }
214 else
215 {
216 assert(size_ > 0);
217 const auto &outer_index = main_cache()->outer_index_;
218 const auto &inner_index = main_cache()->inner_index_;
219 // directly write the values to the matrix
220 mat_ = Eigen::Map<const StiffnessMatrix>(
221 size_, size_, values_.size(), &outer_index[0], &inner_index[0], &values_[0]);
222
223 current_e_ = -1;
224 current_e_index_ = -1;
225
226 }
227 std::fill(values_.begin(), values_.end(), 0);
228 return mat_;
229 }
230
231 std::shared_ptr<MatrixCache> SparseMatrixCache::operator+(const MatrixCache &a) const
232 {
233 assert(&a == &dynamic_cast<const SparseMatrixCache &>(a));
234 return *this + dynamic_cast<const SparseMatrixCache &>(a);
235 }
236
237 std::shared_ptr<MatrixCache> SparseMatrixCache::operator+(const SparseMatrixCache &a) const
238 {
239 std::shared_ptr<SparseMatrixCache> out = std::make_shared<SparseMatrixCache>(a);
240
241 if (a.mapping().empty() || mapping().empty())
242 {
243 out->mat_ = a.mat_ + mat_;
244 const size_t this_e_size = second_cache_entries_.size();
245 const size_t a_e_size = a.second_cache_entries_.size();
246
247 out->second_cache_entries_.resize(std::max(this_e_size, a_e_size));
248 for (int e = 0; e < std::min(this_e_size, a_e_size); ++e)
249 {
250 assert(second_cache_entries_[e].size() == 0 || a.second_cache_entries_[e].size() == 0);
251 out->second_cache_entries_[e].insert(out->second_cache_entries_[e].end(), second_cache_entries_[e].begin(), second_cache_entries_[e].end());
252 out->second_cache_entries_[e].insert(out->second_cache_entries_[e].end(), a.second_cache_entries_[e].begin(), a.second_cache_entries_[e].end());
253 }
254
255 for (int e = std::min(this_e_size, a_e_size); e < std::max(this_e_size, a_e_size); ++e)
256 {
257 if (second_cache_entries_.size() < e)
258 out->second_cache_entries_[e].insert(out->second_cache_entries_[e].end(), second_cache_entries_[e].begin(), second_cache_entries_[e].end());
259 else
260 out->second_cache_entries_[e].insert(out->second_cache_entries_[e].end(), a.second_cache_entries_[e].begin(), a.second_cache_entries_[e].end());
261 }
262 }
263 else
264 {
265 const auto &outer_index = main_cache()->outer_index_;
266 const auto &inner_index = main_cache()->inner_index_;
267 const auto &aouter_index = a.main_cache()->outer_index_;
268 const auto &ainner_index = a.main_cache()->inner_index_;
269 assert(ainner_index.size() == inner_index.size());
270 assert(aouter_index.size() == outer_index.size());
271 assert(a.values_.size() == values_.size());
272
273 maybe_parallel_for(a.values_.size(), [&](int start, int end, int thread_id) {
274 for (int i = start; i < end; ++i)
275 {
276 out->values_[i] = a.values_[i] + values_[i];
277 }
278 });
279 }
280
281 return out;
282 }
283
285 {
286 assert(&o == &dynamic_cast<const SparseMatrixCache &>(o));
287 *this += dynamic_cast<const SparseMatrixCache &>(o);
288 }
289
291 {
292 if (mapping().empty() || o.mapping().empty())
293 {
294 mat_ += o.mat_;
295
296 const size_t this_e_size = second_cache_entries_.size();
297 const size_t o_e_size = o.second_cache_entries_.size();
298
299 second_cache_entries_.resize(std::max(this_e_size, o_e_size));
300 for (int e = 0; e < o_e_size; ++e)
301 {
302 assert(second_cache_entries_[e].size() == 0 || o.second_cache_entries_[e].size() == 0);
304 }
305 }
306 else
307 {
308 const auto &outer_index = main_cache()->outer_index_;
309 const auto &inner_index = main_cache()->inner_index_;
310 const auto &oouter_index = o.main_cache()->outer_index_;
311 const auto &oinner_index = o.main_cache()->inner_index_;
312 assert(inner_index.size() == oinner_index.size());
313 assert(outer_index.size() == oouter_index.size());
314 assert(values_.size() == o.values_.size());
315
316 maybe_parallel_for(o.values_.size(), [&](int start, int end, int thread_id) {
317 for (int i = start; i < end; ++i)
318 {
319 values_[i] += o.values_[i];
320 }
321 });
322 }
323 }
324
325 // ========================================================================
326
328 {
329 mat_.setZero(size, size);
330 }
331
332 DenseMatrixCache::DenseMatrixCache(const size_t rows, const size_t cols)
333 {
334 mat_.setZero(rows, cols);
335 }
336
338 {
339 init(other);
340 }
341
343 {
344 init(other);
345 }
346
347 void DenseMatrixCache::init(const size_t size)
348 {
349 mat_.setZero(size, size);
350 }
351
352 void DenseMatrixCache::init(const size_t rows, const size_t cols)
353 {
354 mat_.setZero(rows, cols);
355 }
356
358 {
359 init(dynamic_cast<const DenseMatrixCache &>(other));
360 }
361
363 {
364 mat_.setZero(other.mat_.rows(), other.mat_.cols());
365 }
366
368 {
369 mat_.setZero();
370 }
371
372 void DenseMatrixCache::add_value(const int e, const int i, const int j, const double value)
373 {
374 mat_(i, j) += value;
375 }
376
378
380 {
381 return mat_.sparseView();
382 }
383
384 std::shared_ptr<MatrixCache> DenseMatrixCache::operator+(const MatrixCache &a) const
385 {
386 return *this + dynamic_cast<const DenseMatrixCache &>(a);
387 }
388
389 std::shared_ptr<MatrixCache> DenseMatrixCache::operator+(const DenseMatrixCache &a) const
390 {
391 std::shared_ptr<DenseMatrixCache> out = std::make_shared<DenseMatrixCache>(a);
392 out->mat_ += mat_;
393 return out;
394 }
395
397 {
398 *this += dynamic_cast<const DenseMatrixCache &>(o);
399 }
400
402 {
403 mat_ += o.mat_;
404 }
405
406} // namespace polyfem::utils
std::shared_ptr< MatrixCache > operator+(const MatrixCache &a) const override
void init(const size_t size) override
StiffnessMatrix get_matrix(const bool compute_mapping=true) override
void add_value(const int e, const int i, const int j, const double value) override
void operator+=(const MatrixCache &o) override
abstract class used for caching
void operator+=(const MatrixCache &o) override
void add_value(const int e, const int i, const int j, const double value) override
e = element_index, i = global row_index, j = global column_index, value = value to add to matrix if t...
void set_zero() override
set matrix values to zero modifies tmp_, mat_, and values (setting all to zero)
std::vector< int > outer_index_
saves inner/outer indices for sparse matrix
void prune() override
if caches have yet to be constructed, add the saved triplets to mat_ modifies tmp_ and mat_,...
const SparseMatrixCache * main_cache_
std::shared_ptr< MatrixCache > operator+(const MatrixCache &a) const override
const std::vector< std::vector< std::pair< int, size_t > > > & mapping() const
std::vector< std::vector< std::pair< int, int > > > second_cache_entries_
maps element indices to global matrix indices
const SparseMatrixCache * main_cache() const
std::vector< double > values_
buffer for values (corresponds to inner/outer_index_ structure for sparse matrix)
const std::vector< std::vector< int > > & second_cache() const
std::vector< Eigen::Triplet< double > > entries_
contains global matrix indices and corresponding value
StiffnessMatrix get_matrix(const bool compute_mapping=true) override
if the cache is yet to be constructed, save the cached (ordered) indices in inner_index_ and outer_in...
void init(const size_t size) override
set matrix to be size x size
std::vector< std::vector< int > > second_cache_
maps element index to local index
std::vector< std::vector< std::pair< int, size_t > > > mapping_
maps row indices to column index/local index pairs
void maybe_parallel_for(int size, const std::function< void(int, int, int)> &partial_for)
spdlog::logger & logger()
Retrieves the current logger.
Definition Logger.cpp:42
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:22