Loading...
Searching...
No Matches
stencil_solver.tpp
1/*---------------------------------------------------------------------------*\
2 *
3 * bitpit
4 *
5 * Copyright (C) 2015-2021 OPTIMAD engineering Srl
6 *
7 * -------------------------------------------------------------------------
8 * License
9 * This file is part of bitpit.
10 *
11 * bitpit is free software: you can redistribute it and/or modify it
12 * under the terms of the GNU Lesser General Public License v3 (LGPL)
13 * as published by the Free Software Foundation.
14 *
15 * bitpit is distributed in the hope that it will be useful, but WITHOUT
16 * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
17 * FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
18 * License for more details.
19 *
20 * You should have received a copy of the GNU Lesser General Public License
21 * along with bitpit. If not, see <http://www.gnu.org/licenses/>.
22 *
23\*---------------------------------------------------------------------------*/
24
25#include "stencil_solver.hpp"
26
27namespace bitpit {
28
51template<typename stencil_t, typename stencil_container_t>
54 m_stencils(stencils), m_stride(stride)
55{
56}
57
63template<typename stencil_t, typename stencil_container_t>
65 : DiscretizationStencilProxyBaseStorage<stencil_t, stencil_container_t>(stencils, 1)
67}
68
74template<typename stencil_t, typename stencil_container_t>
76{
77 return this->m_stencils->size();
78}
79
86template<typename stencil_t, typename stencil_container_t>
88{
89 return (*(this->m_stencils))[rowIndex];
90}
91
98template<typename stencil_t, typename stencil_container_t>
100{
101 return (*(this->m_stencils))[rowRawIndex];
102}
103
111template<typename stencil_t, typename stencil_container_t>
112const stencil_t & DiscretizationStencilProxyStorage<stencil_t, stencil_container_t>::at(long blockIndex, int componentIdx) const
113{
114 return (*(this->m_stencils))[blockIndex * this->m_stride + componentIdx];
115}
116
124template<typename stencil_t, typename stencil_container_t>
125const stencil_t & DiscretizationStencilProxyStorage<stencil_t, stencil_container_t>::rawAt(std::size_t blockRawIndex, int componentIdx) const
126{
127 return (*(this->m_stencils))[blockRawIndex * this->m_stride + componentIdx];
128}
129
135template<typename stencil_t>
137 : DiscretizationStencilProxyBaseStorage<stencil_t, PiercedStorage<stencil_t>>(stencils, stencils->getFieldCount())
138{
139}
140
146template<typename stencil_t>
148{
149 return this->m_stencils->getKernel()->size() * this->m_stride;
150}
151
158template<typename stencil_t>
160{
161 return this->m_stencils->at(rowIndex);
162}
163
170template<typename stencil_t>
171const stencil_t & DiscretizationStencilProxyStorage<stencil_t, PiercedStorage<stencil_t>>::rawAt(std::size_t rowRawIndex) const
172{
173 return this->m_stencils->rawAt(rowRawIndex);
174}
175
183template<typename stencil_t>
184const stencil_t & DiscretizationStencilProxyStorage<stencil_t, PiercedStorage<stencil_t>>::at(long blockIndex, int componentIdx) const
185{
186 return this->m_stencils->at(blockIndex, componentIdx);
187}
188
196template<typename stencil_t>
197const stencil_t & DiscretizationStencilProxyStorage<stencil_t, PiercedStorage<stencil_t>>::rawAt(std::size_t blockRawIndex, int componentIdx) const
198{
199 return this->m_stencils->rawAt(blockRawIndex, componentIdx);
200}
201
210#if BITPIT_ENABLE_MPI==1
218template<typename stencil_t, typename solver_kernel_t>
219template<typename stencil_container_t, typename... AssemblerKernelArgs>
221 AssemblerKernelArgs&&... assemblerKernelArgs)
222 : DiscretizationStencilSolverAssembler(MPI_COMM_SELF, false, stencils, std::forward<AssemblerKernelArgs>(assemblerKernelArgs)...)
223{
224}
225
235template<typename stencil_t, typename solver_kernel_t>
236template<typename stencil_container_t, typename... AssemblerKernelArgs>
238 const stencil_container_t *stencils,
239 AssemblerKernelArgs&&... assemblerKernelArgs)
240 : DiscretizationStencilSolverAssembler(communicator, partitioned,
241 std::unique_ptr<DiscretizationStencilStorageInterface<stencil_t>>(new DiscretizationStencilProxyStorage<stencil_t, stencil_container_t>(stencils)),
242 std::forward<AssemblerKernelArgs>(assemblerKernelArgs)...)
243{
244}
245#else
253template<typename stencil_t, typename solver_kernel_t>
254template<typename stencil_container_t, typename... AssemblerKernelArgs>
256 AssemblerKernelArgs&&... assemblerKernelArgs)
257 : DiscretizationStencilSolverAssembler(std::unique_ptr<DiscretizationStencilStorageInterface<stencil_t>>(new DiscretizationStencilProxyStorage<stencil_t, stencil_container_t>(stencils)),
258 std::forward<AssemblerKernelArgs>(assemblerKernelArgs)...)
259{
260}
261#endif
262
263#if BITPIT_ENABLE_MPI==1
271template<typename stencil_t, typename solver_kernel_t>
272template<typename... AssemblerKernelArgs>
274 AssemblerKernelArgs&&... assemblerKernelArgs)
275 : DiscretizationStencilSolverAssembler(MPI_COMM_SELF, false, std::move(stencils), std::forward<AssemblerKernelArgs>(assemblerKernelArgs)...)
276{
277}
278
288template<typename stencil_t, typename solver_kernel_t>
289template<typename... AssemblerKernelArgs>
291 std::unique_ptr<DiscretizationStencilStorageInterface<stencil_t>> &&stencils,
292 AssemblerKernelArgs&&... assemblerKernelArgs)
293 : DiscretizationStencilSolverAssembler(communicator, partitioned, std::forward<AssemblerKernelArgs>(assemblerKernelArgs)...)
294#else
302template<typename stencil_t, typename solver_kernel_t>
303template<typename... AssemblerKernelArgs>
304DiscretizationStencilSolverAssembler<stencil_t, solver_kernel_t>::DiscretizationStencilSolverAssembler(std::unique_ptr<DiscretizationStencilStorageInterface<stencil_t>> &&stencils,
305 AssemblerKernelArgs&&... assemblerKernelArgs)
306 : DiscretizationStencilSolverAssembler(std::forward<AssemblerKernelArgs>(assemblerKernelArgs)...)
307#endif
308{
309 setStencils(std::move(stencils));
311 setBlockSize();
313}
314
315#if BITPIT_ENABLE_MPI==1
322template<typename stencil_t, typename solver_kernel_t>
323template<typename... AssemblerKernelArgs>
325 : DiscretizationStencilSolverAssembler(MPI_COMM_SELF, false, std::forward<AssemblerKernelArgs>(assemblerKernelArgs)...)
326{
327}
328
337template<typename stencil_t, typename solver_kernel_t>
338template<typename... AssemblerKernelArgs>
340 AssemblerKernelArgs&&... assemblerKernelArgs)
341 : solver_kernel_type::Assembler(std::forward<AssemblerKernelArgs>(assemblerKernelArgs)...),
342 m_partitioned(partitioned), m_communicator(communicator)
343{
344}
345#else
352template<typename stencil_t, typename solver_kernel_t>
353template<typename... AssemblerKernelArgs>
355 : solver_kernel_type::Assembler(std::forward<AssemblerKernelArgs>(assemblerKernelArgs)...)
356{
357}
358#endif
359
360#if BITPIT_ENABLE_MPI==1
366template<typename stencil_t, typename solver_kernel_t>
371
377template<typename stencil_t, typename solver_kernel_t>
379{
380 return m_communicator;
381}
382#endif
383
389template<typename stencil_t, typename solver_kernel_t>
390typename DiscretizationStencilSolverAssembler<stencil_t, solver_kernel_t>::assembly_options_type DiscretizationStencilSolverAssembler<stencil_t, solver_kernel_t>::getOptions() const
391{
392 assembly_options_type options;
393 options.full = true;
394 options.sorted = false;
395
396 return options;
397}
398
402template<typename stencil_t, typename solver_kernel_t>
403template<typename W, typename V, typename std::enable_if<std::is_fundamental<W>::value>::type *>
408
412template<typename stencil_t, typename solver_kernel_t>
413template<typename W, typename V, std::size_t D, typename std::enable_if<std::is_same<std::array<V, D>, W>::value>::type *>
415{
416 setBlockSize(sizeof(typename StencilVector::weight_type) / sizeof(typename StencilVector::weight_type::value_type));
417}
418
431template<typename stencil_t, typename solver_kernel_t>
432template<typename W, typename V, typename std::enable_if<std::is_same<std::vector<V>, W>::value>::type *>
434{
435 // Get block size
436 //
437 // Block size is evaluated from the first weight.
438 //
439 // The block size is set equal to the square root of the weight size;
440 // if the square root of the weight type is not an integer number, an
441 // exception is throw.
442 m_blockSize = - 1;
443 for (long i = 0; i < getRowCount(); ++i) {
444 const StencilBlock &stencil = getRowStencil(i);
445 std::size_t stencilSize = stencil.size();
446 if (stencilSize == 0) {
447 break;
448 }
449
450 const StencilBlock::weight_type *weightData = stencil.weightData();
451 double blockSize = std::sqrt(weightData[0].size());
452 if (blockSize != std::sqrt(weightData[0].size())) {
453 throw std::runtime_error("Weights size should be a square.");
454 }
455 setBlockSize(static_cast<int>(blockSize));
456 break;
457 }
458
459 if (m_blockSize == -1) {
460 throw std::runtime_error("All weights should have a size greater than zero.");
461 }
462
463#ifdef DEBUG
464 // Validate block size
465 //
466 // All weight sizes should match
467 for (long i = 0; i < getRowCount(); ++i) {
468 const StencilBlock &stencil = getRowStencil(i);
469 const StencilBlock::weight_type *weightData = stencil.weightData();
470 std::size_t stencilSize = stencil.size();
471
472 for (std::size_t k = 0; k < stencilSize; ++k) {
473 if (weightData[k].size() != m_blockSize)) {
474 throw std::runtime_error("All stencils weights should have the same size.");
475 }
476 }
477
478 if (stencil.getConstant().size() != m_blockSize) {
479 throw std::runtime_error("The stencils constant should have the same size of the stencil weights.");
480 }
481 }
482#endif
483}
484
490template<typename stencil_t, typename solver_kernel_t>
492{
493 m_blockSize = blockSize;
494}
495
501template<typename stencil_t, typename solver_kernel_t>
506
510template<typename stencil_t, typename solver_kernel_t>
512{
513 setMatrixSizes(m_stencils->size(), m_stencils->size());
514}
515
522template<typename stencil_t, typename solver_kernel_t>
524{
525 // Set system sizes
526 m_nRows = nRows;
527 m_nCols = nCols;
528
529#if BITPIT_ENABLE_MPI==1
530 // Global system sizes
531 m_nGlobalRows = nRows;
532 m_nGlobalCols = nCols;
533 if (m_partitioned) {
534 MPI_Allreduce(MPI_IN_PLACE, &m_nGlobalRows, 1, MPI_LONG, MPI_SUM, m_communicator);
535 MPI_Allreduce(MPI_IN_PLACE, &m_nGlobalCols, 1, MPI_LONG, MPI_SUM, m_communicator);
536 }
537
538 // Global offsets
539 m_globalRowOffset = 0;
540 m_globalColOffset = 0;
541 if (m_partitioned) {
542 int nProcessors;
543 MPI_Comm_size(m_communicator, &nProcessors);
544
545 std::vector<long> nRankRows(nProcessors);
546 MPI_Allgather(&m_nRows, 1, MPI_LONG, nRankRows.data(), 1, MPI_LONG, m_communicator);
547
548 std::vector<long> nRankCols(nProcessors);
549 MPI_Allgather(&m_nCols, 1, MPI_LONG, nRankCols.data(), 1, MPI_LONG, m_communicator);
550
551 int rank;
552 MPI_Comm_rank(m_communicator, &rank);
553 for (int i = 0; i < rank; ++i) {
554 m_globalRowOffset += nRankRows[i];
555 m_globalColOffset += nRankCols[i];
556 }
557 }
558#endif
559}
560
564template<typename stencil_t, typename solver_kernel_t>
566{
567 long maxRowNZ = 0;
568 for (long n = 0; n < getRowCount(); ++n) {
569 maxRowNZ = std::max(getRowNZCount(n), maxRowNZ);
570 }
571
572 setMaximumRowNZ(maxRowNZ);
573}
574
580template<typename stencil_t, typename solver_kernel_t>
582{
583 m_maxRowNZ = maxRowNZ;
584}
585
591template<typename stencil_t, typename solver_kernel_t>
596
606template<typename stencil_t, typename solver_kernel_t>
611
621template<typename stencil_t, typename solver_kernel_t>
626
635template<typename stencil_t, typename solver_kernel_t>
637{
638 long nRowElements = getBlockSize() * getRowCount();
639
640 return nRowElements;
641}
642
651template<typename stencil_t, typename solver_kernel_t>
653{
654 long nColElements = getBlockSize() * getColCount();
655
656 return nColElements;
657}
658
659#if BITPIT_ENABLE_MPI==1
669template<typename stencil_t, typename solver_kernel_t>
674
684template<typename stencil_t, typename solver_kernel_t>
689
698template<typename stencil_t, typename solver_kernel_t>
700{
701 long nElements = getBlockSize() * getRowGlobalCount();
702
703 return nElements;
704}
705
714template<typename stencil_t, typename solver_kernel_t>
716{
717 long nElements = getBlockSize() * getColGlobalCount();
718
719 return nElements;
720}
721
727template<typename stencil_t, typename solver_kernel_t>
732
738template<typename stencil_t, typename solver_kernel_t>
743
752template<typename stencil_t, typename solver_kernel_t>
754{
755 long offset = getBlockSize() * getRowGlobalOffset();
756
757 return offset;
758}
759
768template<typename stencil_t, typename solver_kernel_t>
770{
771 long offset = getBlockSize() * getColGlobalOffset();
772
773 return offset;
774}
775#endif
776
783template<typename stencil_t, typename solver_kernel_t>
785{
786 const stencil_t &stencil = getRowStencil(rowIndex);
787 std::size_t stencilSize = stencil.size();
788
789 return stencilSize;
790}
791
797template<typename stencil_t, typename solver_kernel_t>
802
814template<typename stencil_t, typename solver_kernel_t>
816{
817 // Get stencil information
818 const stencil_t &stencil = getRowStencil(rowIndex);
819
820 // Get pattern
821 getPattern(stencil, pattern);
822}
823
830template<typename stencil_t, typename solver_kernel_t>
832{
833 std::size_t stencilSize = stencil.size();
834
835 const long *patternData = stencil.patternData();
836 pattern->set(patternData, stencilSize);
837}
838
852template<typename stencil_t, typename solver_kernel_t>
854{
855 // Get stencil information
856 const stencil_t &stencil = getRowStencil(rowIndex);
857
858 // Get values
859 getValues(stencil, values);
860}
861
870template<typename stencil_t, typename solver_kernel_t>
871template<typename U, typename std::enable_if<std::is_fundamental<U>::value>::type *>
873{
874 values->set(stencil.weightData(), m_blockSize * stencil.size());
875}
876
885template<typename stencil_t, typename solver_kernel_t>
886template<typename U, typename std::enable_if<!std::is_fundamental<U>::value>::type *>
888{
889 std::size_t stencilSize = stencil.size();
890 const stencil_weight_type *stencilWeightData = stencil.weightData();
891
892 int nBlockElements = m_blockSize * m_blockSize;
893 std::size_t nRowValues = m_blockSize * stencilSize;
894 std::size_t nValues = nBlockElements * stencilSize;
896 ConstProxyVector<double>::storage_pointer expandedValuesStorage = values->storedData();
897
898
899 for (std::size_t k = 0; k < stencilSize; ++k) {
900 const double *weightData = stencilWeightData[k].data();
901 for (int i = 0; i < m_blockSize; ++i) {
902 int weightOffset = linearalgebra::linearIndexRowMajor(i, 0, m_blockSize, m_blockSize);
903 int valuesOffset = linearalgebra::linearIndexRowMajor(i, m_blockSize * k, m_blockSize, nRowValues);
904
905 std::copy_n(weightData + weightOffset, m_blockSize, expandedValuesStorage + valuesOffset);
906 }
907 }
908}
909
919template<typename stencil_t, typename solver_kernel_t>
921{
922 // Get stencil information
923 const stencil_t &stencil = getRowStencil(rowIndex);
924
925 // Get pattern
926 getPattern(stencil, pattern);
927
928 // Get values
929 getValues(stencil, values);
930}
931
940template<typename stencil_t, typename solver_kernel_t>
942{
943 // Get stencil information
944 const stencil_t &stencil = getRowStencil(rowIndex);
945
946 // Get constant
947 getConstant(stencil, constant);
948}
949
958template<typename stencil_t, typename solver_kernel_t>
959template<typename U, typename std::enable_if<std::is_fundamental<U>::value>::type *>
961{
962 const stencil_weight_type &stencilConstant = stencil.getConstant();
963
964 constant->set(ConstProxyVector<double>::INTERNAL_STORAGE, m_blockSize);
965 ConstProxyVector<double>::storage_pointer constantStorage = constant->storedData();
966 std::copy_n(&stencilConstant, m_blockSize, constantStorage);
967}
968
977template<typename stencil_t, typename solver_kernel_t>
978template<typename U, typename std::enable_if<!std::is_fundamental<U>::value>::type *>
980{
981 const stencil_weight_type &stencilConstant = stencil.getConstant();
982
983 constant->set(ConstProxyVector<double>::INTERNAL_STORAGE, m_blockSize);
984 ConstProxyVector<double>::storage_pointer constantStorage = constant->storedData();
985 for (int i = 0; i < m_blockSize; ++i) {
986 int offset_i = linearalgebra::linearIndexRowMajor(i, 0, m_blockSize, m_blockSize);
987
988 constantStorage[i] = 0;
989 for (int j = 0; j < m_blockSize; ++j) {
990 constantStorage[i] += stencil.getWeightManager().at(stencilConstant, offset_i + j);
991 }
992 }
993}
994
1001template<typename stencil_t, typename solver_kernel_t>
1003{
1004 return m_stencils->at(rowIndex);
1005}
1006
1018template<typename stencil_t, typename solver_kernel_t>
1020{
1021 solver_kernel_t::clear();
1022
1023 std::vector<double>().swap(m_constants);
1024}
1025
1026#if BITPIT_ENABLE_MPI==1
1034template<typename stencil_t, typename solver_kernel_t>
1035template<typename stencil_container_t>
1037{
1038 assembly(MPI_COMM_SELF, false, stencils);
1039}
1040
1048template<typename stencil_t, typename solver_kernel_t>
1049template<typename stencil_container_t>
1050void DiscretizationStencilSolver<stencil_t, solver_kernel_t>::assembly(MPI_Comm communicator, bool partitioned, const stencil_container_t &stencils)
1051#else
1057template<typename stencil_t, typename solver_kernel_t>
1058template<typename stencil_container_t>
1059void DiscretizationStencilSolver<stencil_t, solver_kernel_t>::assembly(const stencil_container_t &stencils)
1060#endif
1061{
1062 // Create the assembler
1063#if BITPIT_ENABLE_MPI==1
1064 DiscretizationStencilSolverAssembler<stencil_t, solver_kernel_t> assembler(communicator, partitioned, &stencils);
1065#else
1067#endif
1068
1069 // Assembly the system
1070 solver_kernel_t::template assembly<DiscretizationStencilSolver<stencil_t, solver_kernel_t>>(assembler, NaturalSystemMatrixOrdering());
1071}
1072
1081template<typename stencil_t, typename solver_kernel_t>
1083{
1084 solver_kernel_t::template assembly<DiscretizationStencilSolver<stencil_t, solver_kernel_t>>(assembler, NaturalSystemMatrixOrdering());
1085}
1086
1092template<typename stencil_t, typename solver_kernel_t>
1094{
1095 // Assemble matrix
1096 solver_kernel_t::matrixAssembly(assembler);
1097
1098 // Assemble constants
1099 assembleConstants(assembler);
1100}
1101
1118template<typename stencil_t, typename solver_kernel_t>
1120{
1121 // Update the system
1122 solver_kernel_t::matrixUpdate(nRows, rows, assembler);
1123
1124 // Update the constants
1125 updateConstants(nRows, rows, assembler);
1126}
1127
1136template<typename stencil_t, typename solver_kernel_t>
1137template<typename stencil_container_t>
1139{
1140 update(this->getRowCount(), nullptr, stencils);
1141}
1142
1152template<typename stencil_t, typename solver_kernel_t>
1153template<typename stencil_container_t>
1154void DiscretizationStencilSolver<stencil_t, solver_kernel_t>::update(const std::vector<long> &rows, const stencil_container_t &stencils)
1155{
1156 update(rows.size(), rows.data(), stencils);
1157}
1158
1171template<typename stencil_t, typename solver_kernel_t>
1172template<typename stencil_container_t>
1173void DiscretizationStencilSolver<stencil_t, solver_kernel_t>::update(std::size_t nRows, const long *rows, const stencil_container_t &stencils)
1174{
1175#if BITPIT_ENABLE_MPI==1
1176 DiscretizationStencilSolverAssembler<stencil_t, solver_kernel_t> assembler(this->getCommunicator(), this->isPartitioned(), &stencils);
1177#else
1179#endif
1180
1181 // Update the system
1182 solver_kernel_t::template update<DiscretizationStencilSolver<stencil_t, solver_kernel_t>>(nRows, rows, assembler);
1183}
1184
1185
1196template<typename stencil_t, typename solver_kernel_t>
1197void DiscretizationStencilSolver<stencil_t, solver_kernel_t>::update(long nRows, const long *rows, const Assembler &assembler)
1198{
1199 solver_kernel_t::template update<DiscretizationStencilSolver<stencil_t, solver_kernel_t>>(nRows, rows, assembler);
1200}
1201
1207template<typename stencil_t, typename solver_kernel_t>
1209{
1210 long nRows = assembler.getRowCount();
1211 int blockSize = assembler.getBlockSize();
1212
1213 // Create the storage
1214 m_constants.resize(nRows * blockSize);
1215
1216 // Update the constants
1217 updateConstants(nRows, nullptr, assembler);
1218}
1219
1229template<typename stencil_t, typename solver_kernel_t>
1231 const Assembler &assembler)
1232{
1233 int blockSize = assembler.getBlockSize();
1234 ConstProxyVector<double> rowConstant(blockSize);
1235 for (std::size_t n = 0; n < nRows; ++n) {
1236 long row;
1237 if (rows) {
1238 row = rows[n];
1239 } else {
1240 row = n;
1241 }
1242
1243 assembler.getRowConstant(n, &rowConstant);
1244 std::copy_n(rowConstant.data(), blockSize, m_constants.data() + row * blockSize);
1245 }
1246}
1247
1251template<typename stencil_t, typename solver_kernel_t>
1253{
1254 // Check if the stencil solver is assembled
1255 if (!this->isAssembled()) {
1256 throw std::runtime_error("Unable to solve the system. The stencil solver is not yet assembled.");
1257 }
1258
1259 // Subtract constant terms to the RHS
1260 long nUnknowns = this->getBlockSize() * this->getRowCount();
1261 double *raw_rhs = this->getRHSRawPtr();
1262 for (long i = 0; i < nUnknowns; ++i) {
1263 raw_rhs[i] -= m_constants[i];
1264 }
1265 this->restoreRHSRawPtr(raw_rhs);
1266
1267 // Solve the system
1268 solver_kernel_t::solve();
1269}
1270
1271}
Metafunction for generating a discretization stencil.
Definition stencil.hpp:59
weight_t * weightData()
Definition stencil.tpp:380
std::size_t size() const
Definition stencil.tpp:259
The DiscretizationStencilProxyBaseStorage class defines a proxy for stencil storage.
DiscretizationStencilProxyBaseStorage(const stencil_container_t *stencils, int stride)
const stencil_t & rawAt(std::size_t rowRawIndex) const override
const stencil_t & at(long rowIndex) const override
DiscretizationStencilProxyStorage(const stencil_container_t *stencils)
The DiscretizationStencilSolverAssembler class defines an assembler for building the stencil solver.
long getRowNZCount(long rowIndex) const override
virtual void getRowConstant(long rowIndex, bitpit::ConstProxyVector< double > *constant) const
void getRowValues(long rowIndex, ConstProxyVector< double > *values) const override
void getRowPattern(long rowIndex, ConstProxyVector< long > *pattern) const override
void setStencils(std::unique_ptr< DiscretizationStencilStorageInterface< stencil_t > > &&stencils)
DiscretizationStencilSolverAssembler(const stencil_container_t *stencils, AssemblerKernelArgs &&... assemblerKernelArgs)
void getValues(const stencil_t &stencil, ConstProxyVector< double > *values) const
void getRowData(long rowIndex, ConstProxyVector< long > *pattern, ConstProxyVector< double > *values) const override
const MPI_Comm & getCommunicator() const override
void getPattern(const stencil_t &stencil, ConstProxyVector< long > *pattern) const
assembly_options_type getOptions() const override
virtual const stencil_t & getRowStencil(long rowIndex) const
void getConstant(const stencil_t &stencil, bitpit::ConstProxyVector< double > *constant) const
void matrixAssembly(const Assembler &assembler)
void assembleConstants(const Assembler &assembler)
void matrixUpdate(long nRows, const long *rows, const Assembler &assembler)
void assembly(const stencil_container_t &stencils)
void update(const stencil_container_t &stencils)
void updateConstants(std::size_t nRows, const long *rows, const Assembler &assembler)
The DiscretizationStencilStorageInterface class defines the interface for stencil storage.
The NaturalSystemMatrixOrdering class defines allows to use a matrix natural ordering.
Metafunction for generating a pierced storage.
Metafunction for generating a list of elements that can be either stored in an external vectror or,...
container_type::pointer storage_pointer
__PXV_POINTER__ data() noexcept
__PXV_STORAGE_POINTER__ storedData() noexcept
__PXV_REFERENCE__ at(std::size_t n)
void set(__PXV_POINTER__ data, std::size_t size)
int linearIndexRowMajor(int row, int col, int nRows, int nCols)
--- layout: doxygen_footer ---