Alien  1.3.0
Developer documentation
Loading...
Searching...
No Matches
StreamVBlockMatrixBuilderT.h
1// -*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
2//-----------------------------------------------------------------------------
3// Copyright 2000-2026 CEA (www.cea.fr) IFPEN (www.ifpenergiesnouvelles.com)
4// See the top-level COPYRIGHT file for details.
5// SPDX-License-Identifier: Apache-2.0
6//-----------------------------------------------------------------------------
7
8
9#pragma once
10
11/*---------------------------------------------------------------------------*/
12/*---------------------------------------------------------------------------*/
13
17#include <alien/kernels/simple_csr/CSRStructInfo.h>
18#include <alien/kernels/simple_csr/SimpleCSRInternal.h>
19#include <alien/kernels/simple_csr/SimpleCSRMatrix.h>
20
21#define USE_VMAP
22
23#ifdef USE_VMAP
24#include <alien/utils/VMap.h>
25#undef KEY_OF
26#undef VALUE_OF
27#define KEY_OF(i) (i).key()
28#define VALUE_OF(i) (i).value()
29#else /* USE_VMAP */
30#include <algorithm>
31#include <cstdlib>
32#include <map>
33#include <utility>
34
35#undef KEY_OF
36#undef VALUE_OF
37#define KEY_OF(i) (i)->first
38#define VALUE_OF(i) (i)->second
39#endif /* USE_VMAP */
40
41/*---------------------------------------------------------------------------*/
42
43namespace Alien
44{
45
46/*---------------------------------------------------------------------------*/
48template <typename ValueT>
50VBlockMatrix& matrix, bool init_and_start)
51: m_matrix(matrix)
52, m_matrix_impl(NULL)
53, m_local_size(0)
54, m_global_size(0)
55, m_local_offset(0)
56, m_ghost_size(-1)
57, m_block_ghost_size(-1)
58, m_matrix_size(0)
59, m_block_matrix_size(0)
60, m_myrank(-1)
61, m_nproc(-1)
62, m_col_ordering(eUndef)
63, m_parallel_mng(NULL)
64, m_trace(NULL)
65, m_state(eNone)
66{
67 if (init_and_start) {
68 init();
69 start();
70 }
71}
72
73/*---------------------------------------------------------------------------*/
75template <typename ValueT>
80
81template <typename ValueT>
82void StreamVBlockMatrixBuilderT<ValueT>::finalize()
83{
84 bool has_error = false;
85 for (auto iter = m_inserters.begin(); iter != m_inserters.end(); ++iter) {
86 Filler& filler = **iter;
87 if (!filler.isEnd()) {
88 if (m_trace)
89 m_trace->warning() << "Inserter #" << filler.getId() << " not at final index";
90 has_error = true;
91 }
92 }
93 if (has_error)
94 throw FatalErrorException(A_FUNCINFO, "All inserters are not at final index");
95
96 m_matrix_impl->updateTimestamp();
97 m_matrix.impl()->unlock();
98}
99
100template <typename ValueT>
101void StreamVBlockMatrixBuilderT<ValueT>::end()
102{
103 _freeInserters();
104}
105/*---------------------------------------------------------------------------*/
106
107template <typename ValueT>
109StreamVBlockMatrixBuilderT<ValueT>::getNewInserter()
110{
111 Inserter* inserter = new Inserter(this, m_inserters.size());
112 m_inserters.add(inserter);
113 return *inserter;
114}
115
116/*---------------------------------------------------------------------------*/
117
118template <typename ValueT>
120StreamVBlockMatrixBuilderT<ValueT>::getInserter(Integer id)
121{
122 if ((id < 0) || (id >= (Integer)m_inserters.size()))
123 throw FatalErrorException(A_FUNCINFO, "Bad Inserter id");
124 return *m_inserters[id];
125}
126
127/*---------------------------------------------------------------------------*/
128
129template <typename ValueT>
130void StreamVBlockMatrixBuilderT<ValueT>::allocate()
131{
132 computeProfile();
133}
134
135/*---------------------------------------------------------------------------*/
136
137template <typename ValueT>
138void StreamVBlockMatrixBuilderT<ValueT>::_freeInserters()
139{
140 for (auto iter = m_inserters.begin(); iter != m_inserters.end(); ++iter)
141 delete *iter;
142 m_inserters.resize(0);
143}
144
145/*---------------------------------------------------------------------------*/
146
147template <typename ValueT>
148void StreamVBlockMatrixBuilderT<ValueT>::init()
149{
150 if (m_state == eInit)
151 return;
152
153 //const ISpace& space = m_matrix.rowSpace();
154 // if (space != m_matrix.colSpace())
155 // throw FatalErrorException(
156 // "stream matrix builder must be used with square matrix");
157
158 const MatrixDistribution& dist = m_matrix.impl()->distribution();
159 const VBlock* vblock = m_matrix.impl()->vblock();
160
161 if (!vblock)
162 throw FatalErrorException(
163 A_FUNCINFO, "Space is not variable block size - Can't use block builder");
164
165 // m_trace = space.traceMng();
166 m_parallel_mng = dist.parallelMng();
167
168 m_matrix.impl()->lock();
169 m_matrix_impl = &m_matrix.impl()->template get<BackEnd::tag::simplecsr>(true);
170 m_matrix_impl->free();
171
172 m_local_size = dist.localRowSize();
173 m_global_size = dist.globalRowSize();
174 m_local_offset = dist.rowOffset();
175 m_col_ordering = eUndef;
176
177 m_state = eInit;
178}
179
180/*---------------------------------------------------------------------------*/
181
182template <typename ValueT>
183void StreamVBlockMatrixBuilderT<ValueT>::start()
184{
185 m_matrix_impl->free();
186 m_state = ePrepared;
187}
188
189/*---------------------------------------------------------------------------*/
190
191template <typename ValueT>
192void StreamVBlockMatrixBuilderT<ValueT>::fillZero()
193{
194 m_matrix_impl->internal()->getValues().fill(0.);
195}
196
197/*---------------------------------------------------------------------------*/
198
199template <typename ValueT>
200const VBlock*
201StreamVBlockMatrixBuilderT<ValueT>::vblock() const
202{
203 return m_matrix.impl()->vblock();
204}
205
206/*---------------------------------------------------------------------------*/
207
208template <typename ValueT>
210{
217 // ALIEN_ASSERT((m_state == ePrepared),("Unexpected state: %d vs
218 // %d",m_state,ePrepared));
219
220 // InternalSpace::StructInfo const& info = m_matrix.space().structInfo();
221
222 // Attributes from classes relocated into this method (never used elsewhere)
223
224#ifdef USE_VMAP
225 typedef VMap<Integer, Integer> RowCols;
226#else /* USE_VMAP */
227 typedef std::map<Integer, Integer> RowCols;
228#endif /* USE_VMAP */
229 UniqueArray<RowCols> row_cols;
230 UniqueArray<RowCols> block_row_cols;
231
232 if (m_parallel_mng == NULL) {
233 m_myrank = 0;
234 m_nproc = 1;
235 }
236 else {
237 m_myrank = m_parallel_mng->commRank();
238 m_nproc = m_parallel_mng->commSize();
239 }
240
241 const VBlock* block_sizes = this->vblock();
242
243 m_ghost_size = 0;
244 m_block_ghost_size = 0;
245 m_offset.resize(m_nproc + 1);
246 Arccore::MessagePassing::mpAllGather(m_parallel_mng,
247 ConstArrayView<Integer>(1, &m_local_offset), m_offset.subView(0, m_nproc));
248 m_offset[m_nproc] = m_global_size;
249
250 // Utilitaire local (en attendant les lambda-functions)
251 class IsLocal
252 {
253 public:
254 IsLocal(const ConstArrayView<Integer> offset, const Integer myrank)
255 : m_offset(offset)
256 , m_myrank(myrank)
257 {}
258 bool operator()(Arccore::Integer col) const
259 {
260 return (col >= m_offset[m_myrank]) && (col < m_offset[m_myrank + 1]);
261 }
262
263 private:
264 const ConstArrayView<Integer> m_offset;
265 const Integer m_myrank;
266 } isLocal(m_offset, m_myrank);
267
268 SimpleCSRInternal::CSRStructInfo& profile = m_matrix_impl->internal()->getCSRProfile();
269 profile.init(m_local_size);
270
271 m_row_size.resize(m_local_size);
272 m_ghost_row_size.resize(m_local_size);
273 m_row_size.fill(1);
274 m_ghost_row_size.fill(0);
275 row_cols.resize(m_local_size);
276 block_row_cols.resize(m_local_size);
277 for (Integer row = 0; row < m_local_size; ++row) {
278 row_cols[row][m_local_offset + row] = 0;
279 block_row_cols[row][m_local_offset + row] = 0;
280 }
281
282 // Pour les blocs
283 m_block_row_size.resize(m_local_size);
284 m_block_ghost_row_size.resize(m_local_size);
285 m_block_row_size.fill(0);
286 m_block_ghost_row_size.fill(0);
287 for (Integer row = 0; row < m_local_size; ++row) {
288 const Integer id = m_local_offset + row;
289 const Integer block_size = block_sizes->size(id);
290 m_block_row_size[row] += block_size * block_size;
291 }
292
293 UniqueArray<UniqueArray<Integer>> m_block_data_index(m_inserters.size());
294
295 Integer nb_ins = 0;
296 // LOOP on inserter to compute filling position
297 for (auto iter = m_inserters.begin(); iter != m_inserters.end(); ++iter) {
298 BaseInserter* ins = *iter;
299 ins->m_data_index.resize(ins->count());
300 m_block_data_index[nb_ins].resize(ins->count());
301 ins->m_n.add(0); // position de end
302 ins->m_block_size_row.add(0);
303 ins->m_block_size_col.add(0);
304 if (m_trace)
305 m_trace->info() << "Inserter id=" << ins->getId() << " count=" << ins->count();
306 for (Integer i = 0; i < ins->count(); ++i) {
307 const Integer row = ins->m_row_index[i] - m_local_offset;
308 const Integer col = ins->m_col_index[i];
309 if (row == m_local_size) {
310 // case of ghost row
311 // save -1, latter will be turn to the null position of the CSR structure
312 ins->m_data_index[i] = -1;
313 m_block_data_index[nb_ins][i] = -1;
314 // if(m_trace) m_trace->info()<<"Ghost Row : "<<i<<" "<<ins->m_row_index[i];
315 }
316 else {
317#ifdef USE_VMAP
318 std::pair<RowCols::iterator, bool> finder = row_cols[row].insert(col);
319 std::pair<RowCols::iterator, bool> block_finder = block_row_cols[row].insert(col);
320#else
321 std::pair<RowCols::iterator, bool> finder =
322 row_cols[row].insert(std::pair<Integer, Integer>(col, 0));
323#endif
324 RowCols::iterator inner_iter = finder.first;
325 RowCols::iterator block_iter = block_finder.first;
326 if (finder.second) {
327 const Integer k = m_row_size[row];
328 const Integer gk = m_ghost_row_size[row];
329 const Integer block_k = m_block_row_size[row];
330 const Integer block_gk = m_block_ghost_row_size[row];
331 if (isLocal(col)) {
332 VALUE_OF(inner_iter) = k - gk;
333 VALUE_OF(block_iter) = block_k - block_gk;
334 ins->m_data_index[i] = k - gk;
335 m_block_data_index[nb_ins][i] = block_k - block_gk;
336 }
337 else {
338 // be careful, position will increment with ghost_offset after
339 VALUE_OF(inner_iter) = gk;
340 VALUE_OF(block_iter) = block_gk;
341 ins->m_data_index[i] = gk;
342 m_block_data_index[nb_ins][i] = block_gk;
343 ++m_ghost_row_size[row];
344 m_block_ghost_row_size[row] += block_sizes->size(ins->m_row_index[i]) * block_sizes->size(ins->m_col_index[i]);
345 }
346 ++m_row_size[row];
347 m_block_row_size[row] += block_sizes->size(ins->m_row_index[i]) * block_sizes->size(ins->m_col_index[i]);
348 }
349 else {
350 ins->m_data_index[i] = VALUE_OF(inner_iter);
351 m_block_data_index[nb_ins][i] = VALUE_OF(block_iter);
352 }
353 }
354 }
355 nb_ins++;
356 }
357
358 ArrayView<Integer> m_row_offsets =
359 m_matrix_impl->internal()->getCSRProfile().getRowOffset();
360 ArrayView<Integer> m_block_row_offsets =
361 m_matrix_impl->internal()->getCSRProfile().getBlockRowOffset();
362 m_matrix_size = 0;
363 m_block_matrix_size = 0;
364 for (Integer row = 0; row < m_local_size; ++row) {
365 m_row_offsets[row] = m_matrix_size;
366 m_matrix_size += m_row_size[row];
367 m_block_row_offsets[row] = m_block_matrix_size;
368 m_block_matrix_size += m_block_row_size[row];
369 }
370 m_row_offsets[m_local_size] = m_matrix_size;
371 m_block_row_offsets[m_local_size] = m_block_matrix_size;
372
373 profile.allocate();
374 ArrayView<Integer> m_cols = profile.getCols();
375
376 //Integer icount = 0;
377 Integer offset = 0;
378 for (Integer row = 0; row < m_local_size; ++row) {
379 const Integer ghost_offset = m_row_size[row] - m_ghost_row_size[row];
380 for (RowCols::iterator iter = row_cols[row].begin(); iter != row_cols[row].end();
381 ++iter) {
382 const Integer col_uid = KEY_OF(iter);
383 // increment position for ghost cols
384 if (!isLocal(col_uid))
385 iter.value() += ghost_offset;
386 m_cols[offset + VALUE_OF(iter)] = col_uid;
387 //++icount;
388 }
389 offset += m_row_size[row];
390 }
391
392 // ALIEN_ASSERT((icount==m_matrix_size),("matrix total size problem"));
393
394 ArrayView<Integer> m_block_cols = profile.getBlockCols();
395
396 Integer previous_block_size = 0, index = 0;
397 for (Integer irow = 0; irow < m_local_size; ++irow) { // Attention, c'est local !!!!
398 for (Integer j = m_row_offsets[irow]; j < m_row_offsets[irow + 1]; ++j) {
399 const Integer col = m_cols[j];
400 if (index == 0) {
401 m_block_cols[index] = 0;
402 }
403 else {
404 m_block_cols[index] = m_block_cols[index - 1] + previous_block_size;
405 }
406 previous_block_size =
407 block_sizes->size(col) * block_sizes->size(irow + m_local_offset);
408 index++;
409 }
410 }
411
412 m_matrix_impl->allocate();
413 m_matrix_impl->internal()->getValues().fill(0.);
414
415 nb_ins = 0;
416 for (auto iter = m_inserters.begin(); iter != m_inserters.end(); ++iter) {
417 BaseInserter* ins = *iter;
418 for (Integer i = 0; i < ins->count(); ++i) {
419 const Integer row = ins->m_row_index[i] - m_local_offset;
420 if (row == m_local_size) {
421 // Ghost row
422 ins->m_data_index[i] = m_block_matrix_size; // equivalent to the null position
423 }
424 else {
425 const Integer col = ins->m_col_index[i];
426 const Integer ghost_offset =
427 (isLocal(col) ? 0 : m_block_row_size[row] - m_block_ghost_row_size[row]);
428 ins->m_data_index[i] =
429 m_block_data_index[nb_ins][i] + m_block_row_offsets[row] + ghost_offset;
430 }
431 }
432 nb_ins++;
433 ins->setMatrixValues(m_matrix_impl->internal()->getDataPtr());
434 ins->m_col_index.dispose();
435 ins->m_row_index.dispose();
436 }
437
438 if (m_nproc > 1)
439 m_matrix_impl->parallelStart(m_offset, m_parallel_mng);
440 else
441 m_matrix_impl->sequentialStart();
442 m_col_ordering = eOwnAndGhost;
443 m_state = eStart;
444
445 m_row_size.resize(0);
446 m_ghost_row_size.resize(0);
447 m_block_row_size.resize(0);
448 m_block_ghost_row_size.resize(0);
449}
450
451/*---------------------------------------------------------------------------*/
452
453} // namespace Alien
454
455/*---------------------------------------------------------------------------*/
456/*---------------------------------------------------------------------------*/
IBlockBuilder.h.
MultiMatrixImpl.h.
VBlock.h.
Computes a matrix distribution.
Integer getId() const
Identifiant de l'inserter dans son StreamVBlockMatrixBuilder.
UniqueArray< Integer > m_data_index
positon of entry in the Matrix CSR structure
StreamVBlockMatrixBuilderT(VBlockMatrix &matrix, bool init_and_start=true)
Variable size block elements for block matrices.
Definition VBlock.h:46
Arccore::Integer size(Arccore::Integer index) const
Get the size of a block.
Definition VBlock.cc:162
-- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature --
Definition BackEnd.h:17