Loading...
Searching...
No Matches
reconstruction.cpp
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 <cassert>
26#include <limits>
27#include "bitpit_LA.hpp"
28#include "bitpit_operators.hpp"
29
30#include "bitpit_private_cblas.hpp"
31#include "bitpit_private_lapacke.hpp"
32
33#include "reconstruction.hpp"
34
35namespace bitpit {
36
49
54
58const bool ReconstructionPolynomial::ENABLE_FAST_PATH_OPTIMIZATIONS = true;
59
63const int ReconstructionPolynomial::MAX_STACK_WORKSPACE_SIZE = 10;
64
68const std::vector<std::vector<uint16_t>> ReconstructionPolynomial::m_countCoefficientCache = generateCountCoefficientCache();
69
74const std::vector<std::vector<uint16_t>> ReconstructionPolynomial::m_countDegreeCoefficientCache = generateCountDegreeCoefficientCache();
75
79std::vector<std::vector<uint16_t>> ReconstructionPolynomial::generateCountCoefficientCache()
80{
81 std::vector<std::vector<uint16_t>> cache(MAX_DIMENSIONS + 1, std::vector<uint16_t>(MAX_DEGREE + 1, 0));
82 for (uint8_t dimensions = 0; dimensions <= MAX_DIMENSIONS; ++dimensions) {
83 for (uint8_t degree = 0; degree <= MAX_DEGREE; ++degree) {
84 cache[dimensions][degree] = countCoefficients(degree, dimensions);
85 }
86 }
87
88 return cache;
89}
90
95std::vector<std::vector<uint16_t>> ReconstructionPolynomial::generateCountDegreeCoefficientCache()
96{
97 std::vector<std::vector<uint16_t>> cache(MAX_DIMENSIONS + 1, std::vector<uint16_t>(MAX_DEGREE + 1, 0));
98 for (uint8_t dimensions = 0; dimensions <= MAX_DIMENSIONS; ++dimensions) {
99 for (uint8_t degree = 0; degree <= MAX_DEGREE; ++degree) {
100 cache[dimensions][degree] = countDegreeCoefficients(degree, dimensions);
101 }
102 }
103
104 return cache;
105}
106
116uint16_t ReconstructionPolynomial::getCoefficientCount(uint8_t degree, uint8_t dimensions)
117{
118 assert(degree <= ReconstructionPolynomial::MAX_DEGREE);
119 assert(dimensions <= ReconstructionPolynomial::MAX_DIMENSIONS);
120
121 return m_countCoefficientCache[dimensions][degree];
122}
123
132const std::vector<uint16_t> & ReconstructionPolynomial::getCoefficientsCount(uint8_t dimensions)
133{
134 assert(dimensions <= ReconstructionPolynomial::MAX_DIMENSIONS);
135
136 return m_countCoefficientCache[dimensions];
137}
138
148uint16_t ReconstructionPolynomial::countCoefficients(uint8_t degree, uint8_t dimensions)
149{
150 uint16_t nCoeffs = 0;
151 for (int i = 0; i <= degree; ++i) {
152 nCoeffs += countDegreeCoefficients(i, dimensions);
153 }
154
155 return nCoeffs;
156}
157
167uint16_t ReconstructionPolynomial::getDegreeCoefficientCount(uint8_t degree, uint8_t dimensions)
168{
169 assert(degree <= ReconstructionPolynomial::MAX_DEGREE);
170 assert(dimensions <= ReconstructionPolynomial::MAX_DIMENSIONS);
171
172 return m_countDegreeCoefficientCache[dimensions][degree];
173}
174
183const std::vector<uint16_t> & ReconstructionPolynomial::getDegreeCoefficientsCount(uint8_t dimensions)
184{
185 assert(dimensions <= ReconstructionPolynomial::MAX_DIMENSIONS);
186
187 return m_countDegreeCoefficientCache[dimensions];
188}
189
199uint16_t ReconstructionPolynomial::countDegreeCoefficients(uint8_t degree, uint8_t dimensions)
200{
201 if (dimensions == 0) {
202 return 0;
203 }
204
205 return static_cast<uint16_t>(utils::factorial(dimensions - 1 + degree) / utils::factorial(dimensions - 1) / utils::factorial(degree));
206}
207
218void ReconstructionPolynomial::evalPointBasisValues(uint8_t degree, uint8_t dimensions, const std::array<double, 3> &origin,
219 const std::array<double, 3> &point, double *csi)
220{
221 // Set 0-th degree coefficients
222 csi[0] = 1.;
223
224 // Set high degree coefficients
225 if (degree >= 1) {
226 int offset = 1;
227 const std::array<double, 3> distance = point - origin;
228
229 // Set 1-st degree coefficients
230 for (int i = 0; i < dimensions; ++i) {
231 csi[offset++] = distance[i];
232 }
233
234 // Set 2-nd degree coefficients
235 if (degree >= 2) {
236 for (int i = 0; i < dimensions; ++i) {
237 csi[offset++] = 0.5 * distance[i] * distance[i];
238 }
239
240 if (dimensions >= 2) {
241 csi[offset++] = distance[0] * distance[1];
242
243 if (dimensions >= 3) {
244 csi[offset++] = distance[0] * distance[2];
245 csi[offset++] = distance[1] * distance[2];
246 }
247 }
248 }
249
250 // Check if all coefficients have been set
251 assert(offset == ReconstructionPolynomial::getCoefficientCount(degree, dimensions));
252 }
253}
254
267void ReconstructionPolynomial::evalPointBasisDerivatives(uint8_t degree, uint8_t dimensions, const std::array<double, 3> &origin,
268 const std::array<double, 3> &point, const std::array<double, 3> &direction,
269 double *dcsi)
270{
271 // Set 0-th degree coefficients
272 dcsi[0] = 0.;
273
274 // Set high degree coefficients
275 if (degree >= 1) {
276 int offset = 1;
277
278 // Set 1-st degree coefficients
279 for (int i = 0; i < dimensions; ++i) {
280 dcsi[offset++] = direction[i];
281 }
282
283 // Set 2-nd degree coefficients
284 if (degree >= 2) {
285 const std::array<double, 3> distance = point - origin;
286
287 for (int i = 0; i < dimensions; ++i) {
288 dcsi[offset++] = distance[i] * direction[i];
289 }
290
291 if (dimensions >= 2) {
292 dcsi[offset++] = distance[0] * direction[1] + distance[1] * direction[0];
293
294 if (dimensions >= 3) {
295 dcsi[offset++] = distance[0] * direction[2] + distance[2] * direction[0];
296 dcsi[offset++] = distance[1] * direction[2] + distance[2] * direction[1];
297 }
298 }
299 }
300
301 // Check if all coefficients have been set
302 assert(offset == ReconstructionPolynomial::getCoefficientCount(degree, dimensions));
303 }
304}
305
338void ReconstructionPolynomial::evalCellBasisValues(uint8_t degree, uint8_t dimensions, const std::array<double, 3> &origin,
339 const Cell &cell, const std::array<double, 3> *vertexCoords,
340 double *csi)
341{
342 // Check if cell type is supported
343 //
344 // Polynomials of degree 0 and 1 are supported on all type of cells, for
345 // pixels and voxels also degree 2 is supported.
346 //
347 // Coefficients are evaluated as the volume average of the Taylor series
348 // expansion coefficients. In order to support higher order degrees on
349 // all type of cells, we need to implement numerical integration of the
350 // expansions terms over the volume of the cell. This can be done using
351 // Gauss quadrature rules and evaluating the terms on the integration
352 // points using the function that evaluate the basis function on a
353 // specified point.
354 ElementType cellType = cell.getType();
355
356 bool cellTypeSupported;
357 if (cellType == ElementType::PIXEL) {
358 cellTypeSupported = (degree <= 2);
359 } else if (cellType == ElementType::VOXEL) {
360 cellTypeSupported = (degree <= 2);
361 } else {
362 cellTypeSupported = (degree <= 1);
363 }
364
365 if (!cellTypeSupported) {
366 throw std::runtime_error("Cell type not supported.");
367 }
368
369 // Set 0-th degree coefficients
370 csi[0] = 1.;
371
372 // Set high degree coefficients
373 if (degree >= 1) {
374 int offset = 1;
375 const std::array<double, 3> distance = cell.evalCentroid(vertexCoords) - origin;
376
377 // Set 1-st degree coefficients
378 for (int i = 0; i < dimensions; ++i) {
379 csi[offset++] = distance[i];
380 }
381
382 // Set 2-nd degree coefficients
383 if (degree >= 2) {
384 double cellSize = cell.evalSize(vertexCoords);
385
386 for (int i = 0; i < dimensions; ++i) {
387 csi[offset++] = 0.5 * (distance[i] * distance[i] + cellSize * cellSize / 12.);
388 }
389
390 if (dimensions >= 2) {
391 csi[offset++] = distance[0] * distance[1];
392
393 if (dimensions >= 3) {
394 csi[offset++] = distance[0] * distance[2];
395 csi[offset++] = distance[1] * distance[2];
396 }
397 }
398 }
399
400 // Check if all coefficients have been set
401 assert(offset == ReconstructionPolynomial::getCoefficientCount(degree, dimensions));
402 }
403}
404
409 : m_degree(0), m_dimensions(0), m_nFields(0), m_nCoeffs(0)
410{
411 initialize(0, 0, {{0., 0., 0.}}, 0);
412}
413
421ReconstructionPolynomial::ReconstructionPolynomial(uint8_t degree, uint8_t dimensions,
422 const std::array<double, 3> &origin, int nFields)
423 : m_degree(0), m_dimensions(0), m_nFields(0), m_nCoeffs(0)
424{
425 initialize(degree, dimensions, origin, nFields);
426}
427
435 : ReconstructionPolynomial(other.getDegree(), other.getDimensions(), other.m_origin, other.m_nFields)
436{
437 if (m_nFields > 0) {
438 int nWeights = m_nCoeffs * m_nFields;
439 std::copy(other.m_coeffs.get(), other.m_coeffs.get() + nWeights, m_coeffs.get());
440 }
441}
442
450{
451 ReconstructionPolynomial tmp(other);
452 swap(tmp);
453
454 return *this;
455}
456
465{
466 std::swap(other.m_degree, m_degree);
467 std::swap(other.m_dimensions, m_dimensions);
468 std::swap(other.m_nCoeffs, m_nCoeffs);
469 std::swap(other.m_nFields, m_nFields);
470 std::swap(other.m_coeffs, m_coeffs);
471 std::swap(other.m_origin, m_origin);
472}
473
485void ReconstructionPolynomial::initialize(uint8_t degree, uint8_t dimensions,
486 const std::array<double, 3> &origin,
487 int nFields, bool release)
488{
489 assert(degree <= ReconstructionPolynomial::MAX_DEGREE);
490 m_degree = degree;
491
492 assert(dimensions <= ReconstructionPolynomial::MAX_DIMENSIONS);
493 m_dimensions = dimensions;
494
495 if (nFields > 0) {
496 int currentStorageSize = m_nCoeffs * m_nFields;
497
498 m_nFields = nFields;
499 m_nCoeffs = ReconstructionPolynomial::getCoefficientCount(m_degree, m_dimensions);
500
501 int storageSize = m_nCoeffs * m_nFields;
502
503 bool reallocate;
504 if (release) {
505 reallocate = (currentStorageSize != storageSize);
506 } else {
507 reallocate = (currentStorageSize < storageSize);
508 }
509
510 if (reallocate) {
511 m_coeffs = std::unique_ptr<double[]>(new double[storageSize]);
512 }
513 } else {
514 clear(release);
515 }
516
517 m_origin = origin;
518}
519
527{
528 m_nFields = 0;
529
530 if (release) {
531 m_nCoeffs = 0;
532 m_coeffs.reset();
533 }
534}
535
542{
543 return m_degree;
544}
545
552{
553 return m_dimensions;
554}
555
561const std::array<double, 3> & ReconstructionPolynomial::getOrigin() const
562{
563 return m_origin;
564}
565
572{
573 return m_nCoeffs;
574}
575
582{
583 return m_nFields;
584}
585
592{
593 return m_coeffs.get();
594}
595
602{
603 return m_coeffs.get();
604}
605
612const double * ReconstructionPolynomial::getCoefficients(int field) const
613{
614 const double *coefficients = m_coeffs.get() + computeFieldCoefficientsOffset(0, field);
615
616 return coefficients;
617}
618
626{
627 double *coefficients = m_coeffs.get() + computeFieldCoefficientsOffset(0, field);
628
629 return coefficients;
630}
631
641const double * ReconstructionPolynomial::getDegreeCoefficients(uint8_t degree, int field) const
642{
643 const double *coefficients = m_coeffs.get() + computeFieldCoefficientsOffset(degree, field);
644
645 return coefficients;
646}
647
657double * ReconstructionPolynomial::getDegreeCoefficients(uint8_t degree, int field)
658{
659 double *coefficients = m_coeffs.get() + computeFieldCoefficientsOffset(degree, field);
660
661 return coefficients;
662}
663
673std::size_t ReconstructionPolynomial::computeFieldCoefficientsOffset(uint8_t degree, int field) const
674{
675 std::size_t offset = field * getFieldCoefficientsStride() + degree;
676
677 return offset;
678}
679
685std::size_t ReconstructionPolynomial::getFieldCoefficientsStride() const
686{
687 return m_nCoeffs;
688}
689
697void ReconstructionPolynomial::computeValue(const std::array<double, 3> &point,
698 int field, double *value) const
699{
700 computeValues(m_degree, point, 1, field, value);
701}
702
709void ReconstructionPolynomial::computeValues(const std::array<double, 3> &point,
710 double *values) const
711{
712 computeValues(m_degree, point, m_nFields, 0, values);
713}
714
722void ReconstructionPolynomial::computeValues(const std::array<double, 3> &point,
723 int nFields, double *values) const
724{
725 computeValues(m_degree, point, nFields, 0, values);
726}
727
736void ReconstructionPolynomial::computeValues(const std::array<double, 3> &point,
737 int nFields, int offset, double *values) const
738{
739 computeValues(m_degree, point, nFields, offset, values);
740}
741
750void ReconstructionPolynomial::computeValue(int degree, const std::array<double, 3> &point,
751 int field, double *values) const
752{
753 computeValues(degree, point, 1, field, values);
754}
755
763void ReconstructionPolynomial::computeValues(int degree, const std::array<double, 3> &point,
764 double *values) const
765{
766 computeValues(degree, point, m_nFields, 0, values);
767}
768
777void ReconstructionPolynomial::computeValues(int degree, const std::array<double, 3> &point,
778 int nFields, double *values) const
779{
780 computeValues(degree, point, nFields, 0, values);
781}
782
792void ReconstructionPolynomial::computeValues(int degree, const std::array<double, 3> &point,
793 int nFields, int offset, double *values) const
794{
795 assert(degree <= getDegree());
796
797 // Early return if there are no field to process
798 if (nFields == 0) {
799 return;
800 }
801
802 // Early return if the polynomial is not initialized
803 if (m_dimensions == 0) {
804 return;
805 }
806
807 // Get values
808 const std::size_t fieldValuesStride = 1;
809 double *fieldValue = values;
810 const double *fieldValueEnd = fieldValue + fieldValuesStride * nFields;
811
812 // Get coefficients
813 std::size_t fieldCoeffsStride = getFieldCoefficientsStride();
814 const double *fieldCoeffs = m_coeffs.get() + fieldCoeffsStride * offset;
815
816 // Constant polynomial
817 if (ENABLE_FAST_PATH_OPTIMIZATIONS && degree == 0) {
818 do {
819 *fieldValue = fieldCoeffs[0];
820
821 fieldCoeffs += fieldCoeffsStride;
822 fieldValue += fieldValuesStride;
823 } while (fieldValue != fieldValueEnd);
824
825 return;
826 }
827
828 // Generic polynomial
829 BITPIT_CREATE_WORKSPACE(csi, double, m_nCoeffs, MAX_STACK_WORKSPACE_SIZE);
830 ReconstructionPolynomial::evalPointBasisValues(degree, m_dimensions, m_origin, point, csi);
831
832 do {
833 *fieldValue = fieldCoeffs[0] * csi[0];
834 for (int i = 1; i < m_nCoeffs; ++i) {
835 *fieldValue += fieldCoeffs[i] * csi[i];
836 }
837
838 fieldCoeffs += fieldCoeffsStride;
839 fieldValue += fieldValuesStride;
840 } while (fieldValue != fieldValueEnd);
841}
842
856void ReconstructionPolynomial::computeValueLimited(const std::array<double, 3> &point,
857 const double *limiters, int field,
858 double *values) const
859{
860 computeValuesLimited(m_degree, point, limiters, 1, field, values);
861}
862
875void ReconstructionPolynomial::computeValuesLimited(const std::array<double, 3> &point,
876 const double *limiters,
877 double *values) const
878{
879 computeValuesLimited(m_degree, point, limiters, m_nFields, 0, values);
880}
881
895void ReconstructionPolynomial::computeValuesLimited(const std::array<double, 3> &point,
896 const double *limiters, int nFields,
897 double *values) const
898{
899 computeValuesLimited(m_degree, point, limiters, nFields, 0, values);
900}
901
916void ReconstructionPolynomial::computeValuesLimited(const std::array<double, 3> &point,
917 const double *limiters, int nFields, int offset,
918 double *values) const
919{
920 computeValuesLimited(m_degree, point, limiters, nFields, offset, values);
921}
922
937void ReconstructionPolynomial::computeValueLimited(int degree, const std::array<double, 3> &point,
938 const double *limiters, int field,
939 double *values) const
940{
941 computeValuesLimited(degree, point, limiters, 1, field, values);
942}
943
957void ReconstructionPolynomial::computeValuesLimited(int degree, const std::array<double, 3> &point,
958 const double *limiters,
959 double *values) const
960{
961 computeValuesLimited(degree, point, limiters, m_nFields, 0, values);
962}
963
978void ReconstructionPolynomial::computeValuesLimited(int degree, const std::array<double, 3> &point,
979 const double *limiters, int nFields,
980 double *values) const
981{
982 computeValuesLimited(degree, point, limiters, nFields, 0, values);
983}
984
1000void ReconstructionPolynomial::computeValuesLimited(int degree, const std::array<double, 3> &point,
1001 const double *limiters, int nFields, int offset,
1002 double *values) const
1003{
1004 assert(degree <= getDegree());
1005
1006 // Early return if there are no field to process
1007 if (nFields == 0) {
1008 return;
1009 }
1010
1011 // Early return if the polynomial is not initialized
1012 if (m_dimensions == 0) {
1013 return;
1014 }
1015
1016 // Get values
1017 std::size_t fieldValuesStride = 1;
1018 double *fieldValue = values;
1019 double *fieldValueEnd = fieldValue + fieldValuesStride * nFields;
1020
1021 // Get coefficients
1022 std::size_t fieldCoeffsStride = getFieldCoefficientsStride();
1023 const double *fieldCoeffs = m_coeffs.get() + fieldCoeffsStride * offset;
1024
1025 // Constant polynomial
1026 if (ENABLE_FAST_PATH_OPTIMIZATIONS && degree == 0) {
1027 do {
1028 *fieldValue = fieldCoeffs[0];
1029
1030 fieldCoeffs += fieldCoeffsStride;
1031 fieldValue += fieldValuesStride;
1032 } while (fieldValue != fieldValueEnd);
1033
1034 return;
1035 }
1036
1037 // Get limiters
1038 std::size_t fieldLimitersStride = degree;
1039 const double *fieldLimiters = limiters;
1040
1041 // Generic polynomial
1042 const std::vector<uint16_t> &nDegreeCoeffs = ReconstructionPolynomial::getDegreeCoefficientsCount(m_dimensions);
1043
1044 BITPIT_CREATE_WORKSPACE(csi, double, nDegreeCoeffs[degree], MAX_STACK_WORKSPACE_SIZE);
1045 ReconstructionPolynomial::evalPointBasisValues(degree, m_dimensions, m_origin, point, csi);
1046
1047 do {
1048 // Degree 0
1049 *fieldValue = fieldCoeffs[0] * csi[0];
1050
1051 // Degrees greater than 0
1052 int coeffEnd = nDegreeCoeffs[0];
1053 for (int n = 1; n <= degree; ++n) {
1054 int coeffBegin = coeffEnd;
1055 coeffEnd = coeffBegin + nDegreeCoeffs[n];
1056
1057 double fieldDegreeValue = 0;
1058 for (int i = coeffBegin; i < coeffEnd; ++i) {
1059 fieldDegreeValue += fieldCoeffs[i] * csi[i];
1060 }
1061
1062 *fieldValue += fieldDegreeValue * fieldLimiters[n - 1];
1063 }
1064
1065 fieldCoeffs += fieldCoeffsStride;
1066 fieldValue += fieldValuesStride;
1067 fieldLimiters += fieldLimitersStride;
1068 } while (fieldValue != fieldValueEnd);
1069}
1070
1079void ReconstructionPolynomial::computeDerivative(const std::array<double, 3> &point,
1080 const std::array<double, 3> &direction,
1081 int field, double *derivative) const
1082{
1083 computeDerivatives(m_degree, point, direction, 1, field, derivative);
1084}
1085
1093void ReconstructionPolynomial::computeDerivatives(const std::array<double, 3> &point,
1094 const std::array<double, 3> &direction,
1095 double *derivatives) const
1096{
1097 computeDerivatives(m_degree, point, direction, m_nFields, 0, derivatives);
1098}
1099
1108void ReconstructionPolynomial::computeDerivatives(const std::array<double, 3> &point,
1109 const std::array<double, 3> &direction,
1110 int nFields, double *derivatives) const
1111{
1112 computeDerivatives(m_degree, point, direction, nFields, 0, derivatives);
1113}
1114
1124void ReconstructionPolynomial::computeDerivatives(const std::array<double, 3> &point,
1125 const std::array<double, 3> &direction,
1126 int nFields, int offset, double *derivatives) const
1127{
1128 computeDerivatives(m_degree, point, direction, nFields, offset, derivatives);
1129}
1130
1140void ReconstructionPolynomial::computeDerivative(int degree, const std::array<double, 3> &point,
1141 const std::array<double, 3> &direction,
1142 int field, double *derivative) const
1143{
1144 computeDerivatives(degree, point, direction, 1, field, derivative);
1145}
1146
1155void ReconstructionPolynomial::computeDerivatives(int degree, const std::array<double, 3> &point,
1156 const std::array<double, 3> &direction,
1157 double *derivatives) const
1158{
1159 computeDerivatives(degree, point, direction, m_nFields, 0, derivatives);
1160}
1161
1171void ReconstructionPolynomial::computeDerivatives(int degree, const std::array<double, 3> &point,
1172 const std::array<double, 3> &direction,
1173 int nFields, double *derivatives) const
1174{
1175 computeDerivatives(degree, point, direction, nFields, 0, derivatives);
1176}
1177
1188void ReconstructionPolynomial::computeDerivatives(int degree, const std::array<double, 3> &point,
1189 const std::array<double, 3> &direction,
1190 int nFields, int offset, double *derivatives) const
1191{
1192 assert(degree <= getDegree());
1193
1194 // Early return if there are no field to process
1195 if (nFields == 0) {
1196 return;
1197 }
1198
1199 // Early return if the polynomial is not initialized
1200 if (m_dimensions == 0) {
1201 return;
1202 }
1203
1204 // Get derivatives
1205 std::size_t fieldDerivativeStride = 1;
1206 double *fieldDerivative = derivatives;
1207 double *fieldDerivativeEnd = fieldDerivative + fieldDerivativeStride * nFields;
1208
1209 // Constant polynomial
1210 if (ENABLE_FAST_PATH_OPTIMIZATIONS && degree == 0) {
1211 do {
1212 *fieldDerivative = 0.;
1213
1214 fieldDerivative += fieldDerivativeStride;
1215 } while (fieldDerivative != fieldDerivativeEnd);
1216
1217 return;
1218 }
1219
1220 // Get coefficients
1221 std::size_t fieldCoeffsStride = getFieldCoefficientsStride();
1222 const double *fieldCoeffs = m_coeffs.get() + fieldCoeffsStride * offset;
1223
1224 // Linear polynomial
1225 if (ENABLE_FAST_PATH_OPTIMIZATIONS && degree == 1) {
1226 do {
1227 *fieldDerivative = fieldCoeffs[1] * direction[0];
1228 for (int d = 1; d < m_dimensions; ++d) {
1229 *fieldDerivative += fieldCoeffs[d + 1] * direction[d];
1230 }
1231
1232 fieldCoeffs += fieldCoeffsStride;
1233 fieldDerivative += fieldDerivativeStride;
1234 } while (fieldDerivative != fieldDerivativeEnd);
1235
1236 return;
1237 }
1238
1239 // Generic polynomial
1240 BITPIT_CREATE_WORKSPACE(csi, double, m_nCoeffs, MAX_STACK_WORKSPACE_SIZE);
1241 ReconstructionPolynomial::evalPointBasisDerivatives(degree, m_dimensions, m_origin, point, direction, csi);
1242
1243 do {
1244 *fieldDerivative = fieldCoeffs[0] * csi[0];
1245 for (int i = 1; i < m_nCoeffs; ++i) {
1246 *fieldDerivative += fieldCoeffs[i] * csi[i];
1247 }
1248
1249 fieldCoeffs += fieldCoeffsStride;
1250 fieldDerivative += fieldDerivativeStride;
1251 } while (fieldDerivative != fieldDerivativeEnd);
1252}
1253
1268void ReconstructionPolynomial::computeDerivativeLimited(const std::array<double, 3> &point,
1269 const std::array<double, 3> &direction,
1270 const double *limiters, int field,
1271 double *derivative) const
1272{
1273 computeDerivativesLimited(m_degree, point, direction, limiters, 1, field, derivative);
1274}
1275
1289void ReconstructionPolynomial::computeDerivativesLimited(const std::array<double, 3> &point,
1290 const std::array<double, 3> &direction,
1291 const double *limiters,
1292 double *derivatives) const
1293{
1294 computeDerivativesLimited(m_degree, point, direction, limiters, m_nFields, 0, derivatives);
1295}
1296
1311void ReconstructionPolynomial::computeDerivativesLimited(const std::array<double, 3> &point,
1312 const std::array<double, 3> &direction,
1313 const double *limiters, int nFields,
1314 double *derivatives) const
1315{
1316 computeDerivativesLimited(m_degree, point, direction, limiters, nFields, 0, derivatives);
1317}
1318
1334void ReconstructionPolynomial::computeDerivativesLimited(const std::array<double, 3> &point,
1335 const std::array<double, 3> &direction,
1336 const double *limiters, int nFields, int offset,
1337 double *derivatives) const
1338{
1339 computeDerivativesLimited(m_degree, point, direction, limiters, nFields, offset, derivatives);
1340}
1341
1357void ReconstructionPolynomial::computeDerivativeLimited(int degree, const std::array<double, 3> &point,
1358 const std::array<double, 3> &direction,
1359 const double *limiters, int field,
1360 double *derivative) const
1361{
1362 computeDerivativesLimited(degree, point, direction, limiters, 1, field, derivative);
1363}
1364
1379void ReconstructionPolynomial::computeDerivativesLimited(int degree, const std::array<double, 3> &point,
1380 const std::array<double, 3> &direction,
1381 const double *limiters,
1382 double *derivatives) const
1383{
1384 computeDerivativesLimited(degree, point, direction, limiters, m_nFields, 0, derivatives);
1385}
1386
1402void ReconstructionPolynomial::computeDerivativesLimited(int degree, const std::array<double, 3> &point,
1403 const std::array<double, 3> &direction,
1404 const double *limiters, int nFields,
1405 double *derivatives) const
1406{
1407 computeDerivativesLimited(degree, point, direction, limiters, nFields, 0, derivatives);
1408}
1409
1426void ReconstructionPolynomial::computeDerivativesLimited(int degree, const std::array<double, 3> &point,
1427 const std::array<double, 3> &direction,
1428 const double *limiters, int nFields, int offset,
1429 double *derivatives) const
1430{
1431 assert(degree <= getDegree());
1432
1433 // Early return if there are no field to process
1434 if (nFields == 0) {
1435 return;
1436 }
1437
1438 // Early return if the polynomial is not initialized
1439 if (m_dimensions == 0) {
1440 return;
1441 }
1442
1443 // Get derivatives
1444 std::size_t fieldDerivativeStride = 1;
1445 double *fieldDerivative = derivatives;
1446 double *fieldDerivativeEnd = fieldDerivative + fieldDerivativeStride * nFields;
1447
1448 // Constant polynomial
1449 if (ENABLE_FAST_PATH_OPTIMIZATIONS && degree == 0) {
1450 do {
1451 *fieldDerivative = 0.;
1452
1453 fieldDerivative += fieldDerivativeStride;
1454 } while (fieldDerivative != fieldDerivativeEnd);
1455
1456 return;
1457 }
1458
1459 // Get coefficients
1460 std::size_t fieldCoeffsStride = getFieldCoefficientsStride();
1461 const double *fieldCoeffs = m_coeffs.get() + fieldCoeffsStride * offset;
1462
1463 // Get limiters
1464 std::size_t fieldLimitersStride = degree;
1465 const double *fieldLimiters = limiters;
1466
1467 // Linear polynomial
1468 if (ENABLE_FAST_PATH_OPTIMIZATIONS && degree == 1) {
1469 do {
1470 *fieldDerivative = fieldCoeffs[1] * direction[0];
1471 for (int d = 1; d < m_dimensions; ++d) {
1472 *fieldDerivative += fieldCoeffs[d + 1] * direction[d];
1473 }
1474
1475 *fieldDerivative *= fieldLimiters[0];
1476
1477 fieldCoeffs += fieldCoeffsStride;
1478 fieldDerivative += fieldDerivativeStride;
1479 fieldLimiters += fieldLimitersStride;
1480 } while (fieldDerivative != fieldDerivativeEnd);
1481
1482 return;
1483 }
1484
1485 // Generic polynomial
1486 const std::vector<uint16_t> &nDegreeCoeffs = ReconstructionPolynomial::getDegreeCoefficientsCount(m_dimensions);
1487
1488 BITPIT_CREATE_WORKSPACE(dcsi, double, nDegreeCoeffs[degree], MAX_STACK_WORKSPACE_SIZE);
1489 ReconstructionPolynomial::evalPointBasisDerivatives(degree, m_dimensions, m_origin, point, direction, dcsi);
1490
1491 do {
1492 *fieldDerivative = fieldCoeffs[0] * dcsi[0];
1493
1494 int coeffEnd = nDegreeCoeffs[0];
1495 for (int n = 1; n <= degree; ++n) {
1496 int coeffBegin = coeffEnd;
1497 coeffEnd = coeffBegin + nDegreeCoeffs[n];
1498
1499 double fieldDegreeDerivative = 0;
1500 for (int i = coeffBegin; i < coeffEnd; ++i) {
1501 fieldDegreeDerivative += fieldCoeffs[i] * dcsi[i];
1502 }
1503
1504 *fieldDerivative += fieldLimiters[n - 1] * fieldDegreeDerivative;
1505 }
1506
1507 fieldCoeffs += fieldCoeffsStride;
1508 fieldDerivative += fieldDerivativeStride;
1509 fieldLimiters += fieldLimitersStride;
1510 } while (fieldDerivative != fieldDerivativeEnd);
1511}
1512
1520void ReconstructionPolynomial::computeGradient(const std::array<double, 3> &point,
1521 int field, std::array<double, 3> *gradient) const
1522{
1523 computeGradients(m_degree, point, 1, field, gradient);
1524}
1525
1532void ReconstructionPolynomial::computeGradients(const std::array<double, 3> &point,
1533 std::array<double, 3> *gradients) const
1534{
1535 computeGradients(m_degree, point, m_nFields, 0, gradients);
1536}
1537
1545void ReconstructionPolynomial::computeGradients(const std::array<double, 3> &point,
1546 int nFields, std::array<double, 3> *gradients) const
1547{
1548 computeGradients(m_degree, point, nFields, 0, gradients);
1549}
1550
1559void ReconstructionPolynomial::computeGradients(const std::array<double, 3> &point,
1560 int nFields, int offset,
1561 std::array<double, 3> *gradients) const
1562{
1563 computeGradients(m_degree, point, nFields, offset, gradients);
1564}
1565
1574void ReconstructionPolynomial::computeGradient(int degree, const std::array<double, 3> &point,
1575 int field, std::array<double, 3> *gradient) const
1576{
1577 computeGradients(degree, point, 1, field, gradient);
1578}
1579
1587void ReconstructionPolynomial::computeGradients(int degree, const std::array<double, 3> &point,
1588 std::array<double, 3> *gradients) const
1589{
1590 computeGradients(degree, point, m_nFields, 0, gradients);
1591}
1592
1601void ReconstructionPolynomial::computeGradients(int degree, const std::array<double, 3> &point,
1602 int nFields, std::array<double, 3> *gradients) const
1603{
1604 computeGradients(degree, point, nFields, 0, gradients);
1605}
1606
1616void ReconstructionPolynomial::computeGradients(int degree, const std::array<double, 3> &point,
1617 int nFields, int offset,
1618 std::array<double, 3> *gradients) const
1619{
1620 assert(degree <= getDegree());
1621
1622 // Early return if there are no field to process
1623 if (nFields == 0) {
1624 return;
1625 }
1626
1627 // Early return if the polynomial is not initialized
1628 if (m_dimensions == 0) {
1629 return;
1630 }
1631
1632 // Get gradients
1633 const std::size_t fieldGradientStride = 1;
1634 std::array<double, 3> *fieldGradient = gradients;
1635 const std::array<double, 3> *fieldGradientEnd = fieldGradient + fieldGradientStride * nFields;
1636
1637 // Constant polynomial
1638 if (ENABLE_FAST_PATH_OPTIMIZATIONS && degree == 0) {
1639 do {
1640 *fieldGradient = {{0., 0., 0.}};
1641
1642 fieldGradient += fieldGradientStride;
1643 } while (fieldGradient != fieldGradientEnd);
1644
1645 return;
1646 }
1647
1648 // Get coefficients
1649 const std::size_t fieldCoeffsStride = getFieldCoefficientsStride();
1650 const double *fieldCoeffs = m_coeffs.get() + fieldCoeffsStride * offset;
1651
1652 // Linear polynomial
1653 if (ENABLE_FAST_PATH_OPTIMIZATIONS && degree == 1) {
1654 do {
1655 // Evaluate gradients
1656 for (int d = 0; d < m_dimensions; ++d) {
1657 (*fieldGradient)[d] = fieldCoeffs[1 + d];
1658 }
1659
1660 // Explicitly zero unused components
1661 for (int d = m_dimensions; d < ReconstructionPolynomial::MAX_DIMENSIONS; ++d) {
1662 (*fieldGradient)[d] = 0.;
1663 }
1664
1665 // Advance to the next field
1666 fieldCoeffs += fieldCoeffsStride;
1667 fieldGradient += fieldGradientStride;
1668 } while (fieldGradient != fieldGradientEnd);
1669
1670 return;
1671 }
1672
1673 // Generic polynomial
1674 BITPIT_CREATE_WORKSPACE(dcsi, double, static_cast<std::size_t>(m_dimensions * m_nCoeffs), 3 * MAX_STACK_WORKSPACE_SIZE);
1675 ReconstructionPolynomial::evalPointBasisDerivatives(degree, m_dimensions, m_origin, point, {{1., 0., 0.}}, dcsi);
1676 ReconstructionPolynomial::evalPointBasisDerivatives(degree, m_dimensions, m_origin, point, {{0., 1., 0.}}, dcsi + m_nCoeffs);
1677 if (m_dimensions == 3) {
1678 ReconstructionPolynomial::evalPointBasisDerivatives(degree, m_dimensions, m_origin, point, {{0., 0., 1.}}, dcsi + 2 * m_nCoeffs);
1679 }
1680
1681 do {
1682 // Evaluate gradient
1683 for (int d = 0; d < m_dimensions; ++d) {
1684 const double *dcsi_dimension = dcsi + d * m_nCoeffs;
1685
1686 (*fieldGradient)[d] = fieldCoeffs[0] * dcsi_dimension[0];
1687 for (int i = 1; i < m_nCoeffs; ++i) {
1688 (*fieldGradient)[d] += fieldCoeffs[i] * dcsi_dimension[i];
1689 }
1690 }
1691
1692 // Explicitly zero unused components
1693 for (int d = m_dimensions; d < ReconstructionPolynomial::MAX_DIMENSIONS; ++d) {
1694 (*fieldGradient)[d] = 0.;
1695 }
1696
1697 // Advance to the next field
1698 fieldCoeffs += fieldCoeffsStride;
1699 fieldGradient += fieldGradientStride;
1700 } while (fieldGradient != fieldGradientEnd);
1701}
1702
1716void ReconstructionPolynomial::computeGradientLimited(const std::array<double, 3> &point,
1717 const double *limiters, int field,
1718 std::array<double, 3> *gradient) const
1719{
1720 computeGradientsLimited(m_degree, point, limiters, 1, field, gradient);
1721}
1722
1735void ReconstructionPolynomial::computeGradientsLimited(const std::array<double, 3> &point,
1736 const double *limiters,
1737 std::array<double, 3> *gradients) const
1738{
1739 computeGradientsLimited(m_degree, point, limiters, m_nFields, 0, gradients);
1740}
1741
1755void ReconstructionPolynomial::computeGradientsLimited(const std::array<double, 3> &point,
1756 const double *limiters, int nFields,
1757 std::array<double, 3> *gradients) const
1758{
1759 computeGradientsLimited(m_degree, point, limiters, nFields, 0, gradients);
1760}
1761
1776void ReconstructionPolynomial::computeGradientsLimited(const std::array<double, 3> &point,
1777 const double *limiters, int nFields, int offset,
1778 std::array<double, 3> *gradients) const
1779{
1780 computeGradientsLimited(m_degree, point, limiters, nFields, offset, gradients);
1781}
1782
1797void ReconstructionPolynomial::computeGradientLimited(int degree, const std::array<double, 3> &point,
1798 const double *limiters, int field,
1799 std::array<double, 3> *gradient) const
1800{
1801 computeGradientsLimited(degree, point, limiters, 1, field, gradient);
1802}
1803
1817void ReconstructionPolynomial::computeGradientsLimited(int degree, const std::array<double, 3> &point,
1818 const double *limiters,
1819 std::array<double, 3> *gradients) const
1820{
1821 computeGradientsLimited(degree, point, limiters, m_nFields, 0, gradients);
1822}
1823
1838void ReconstructionPolynomial::computeGradientsLimited(int degree, const std::array<double, 3> &point,
1839 const double *limiters, int nFields,
1840 std::array<double, 3> *gradients) const
1841{
1842 computeGradientsLimited(degree, point, limiters, nFields, 0, gradients);
1843}
1844
1860void ReconstructionPolynomial::computeGradientsLimited(int degree, const std::array<double, 3> &point,
1861 const double *limiters, int nFields, int offset,
1862 std::array<double, 3> *gradients) const
1863{
1864 assert(degree <= getDegree());
1865
1866 // Early return if there are no field to process
1867 if (nFields == 0) {
1868 return;
1869 }
1870
1871 // Early return if the polynomial is not initialized
1872 if (m_dimensions == 0) {
1873 return;
1874 }
1875
1876 // Get gradients
1877 const std::size_t fieldGradientStride = 1;
1878 std::array<double, 3> *fieldGradient = gradients;
1879 const std::array<double, 3> *fieldGradientEnd = fieldGradient + fieldGradientStride * nFields;
1880
1881 // Constant polynomial
1882 if (ENABLE_FAST_PATH_OPTIMIZATIONS && degree == 0) {
1883 do {
1884 *fieldGradient = {{0., 0., 0.}};
1885
1886 fieldGradient += fieldGradientStride;
1887 } while (fieldGradient != fieldGradientEnd);
1888
1889 return;
1890 }
1891
1892 // Get coefficients
1893 const std::size_t fieldCoeffsStride = getFieldCoefficientsStride();
1894 const double *fieldCoeffs = m_coeffs.get() + fieldCoeffsStride * offset;
1895
1896 // Get limiters
1897 const std::size_t fieldLimitersStride = degree;
1898 const double *fieldLimiters = limiters;
1899
1900 // Linear polynomial
1901 if (ENABLE_FAST_PATH_OPTIMIZATIONS && degree == 1) {
1902 do {
1903 // Evaluate gradients
1904 for (int d = 0; d < m_dimensions; ++d) {
1905 (*fieldGradient)[d] = fieldLimiters[0] * fieldCoeffs[1 + d];
1906 }
1907
1908 // Explicitly zero unused components
1909 for (int d = m_dimensions; d < ReconstructionPolynomial::MAX_DIMENSIONS; ++d) {
1910 (*fieldGradient)[d] = 0.;
1911 }
1912
1913 // Advance to the next field
1914 fieldCoeffs += fieldCoeffsStride;
1915 fieldGradient += fieldGradientStride;
1916 fieldLimiters += fieldLimitersStride;
1917 } while (fieldGradient != fieldGradientEnd);
1918
1919 return;
1920 }
1921
1922 // Generic polynomial
1923 const std::vector<uint16_t> &nDegreeCoeffs = ReconstructionPolynomial::getDegreeCoefficientsCount(m_dimensions);
1924
1925 BITPIT_CREATE_WORKSPACE(dcsi, double, static_cast<std::size_t>(m_dimensions * nDegreeCoeffs[degree]), 3 * MAX_STACK_WORKSPACE_SIZE);
1926 ReconstructionPolynomial::evalPointBasisDerivatives(degree, m_dimensions, m_origin, point, {{1., 0., 0.}}, dcsi);
1927 ReconstructionPolynomial::evalPointBasisDerivatives(degree, m_dimensions, m_origin, point, {{0., 1., 0.}}, dcsi + nDegreeCoeffs[degree]);
1928 if (m_dimensions == 3) {
1929 ReconstructionPolynomial::evalPointBasisDerivatives(degree, m_dimensions, m_origin, point, {{0., 0., 1.}}, dcsi + 2 * nDegreeCoeffs[degree]);
1930 }
1931
1932 do {
1933 // Evaluate gradients
1934 for (int d = 0; d < m_dimensions; ++d) {
1935 const double *dcsi_dimension = dcsi + d * nDegreeCoeffs[degree];
1936
1937 (*fieldGradient)[d] = fieldCoeffs[0] * dcsi_dimension[0];
1938 }
1939
1940 int coeffEnd = nDegreeCoeffs[0];
1941 for (int n = 1; n <= degree; ++n) {
1942 double fieldsLimiter = fieldLimiters[n - 1];
1943
1944 int coeffBegin = coeffEnd;
1945 coeffEnd = coeffBegin + nDegreeCoeffs[n];
1946 for (int i = coeffBegin; i < coeffEnd; ++i) {
1947 for (int d = 0; d < m_dimensions; ++d) {
1948 const double *dcsi_dimension = dcsi + d * nDegreeCoeffs[degree];
1949
1950 (*fieldGradient)[d] += fieldsLimiter * fieldCoeffs[i] * dcsi_dimension[i];
1951 }
1952 }
1953 }
1954
1955 // Explicitly zero unused components
1956 for (int d = m_dimensions; d < ReconstructionPolynomial::MAX_DIMENSIONS; ++d) {
1957 (*fieldGradient)[d] = 0.;
1958 }
1959
1960 // Advance to the next field
1961 fieldCoeffs += fieldCoeffsStride;
1962 fieldGradient += fieldGradientStride;
1963 fieldLimiters += fieldLimitersStride;
1964 } while (fieldGradient != fieldGradientEnd);
1965}
1966
1972void ReconstructionPolynomial::display(std::ostream &out) const
1973{
1974 uint8_t dimensions = getDimensions();
1975
1976 for (int k = 0; k < m_nFields; ++k) {
1977 out << " field " << k << "\n";
1978 for (int i = 0; i <= m_degree; ++i) {
1979 out << " degree = " << i << " : " ;
1980
1981 uint16_t nDegreeCoeffs = ReconstructionPolynomial::getDegreeCoefficientCount(i, dimensions);
1982 const double *degreeCoeffs = getDegreeCoefficients(i, k);
1983 for (int n = 0; n < nDegreeCoeffs; ++n) {
1984 out << degreeCoeffs[n];
1985 if (n != nDegreeCoeffs - 1) {
1986 out << " , ";
1987 }
1988 }
1989
1990 out << "\n";
1991 }
1992 }
1993}
1994
2006const bool ReconstructionKernel::ENABLE_FAST_PATH_OPTIMIZATIONS = true;
2007
2011const int ReconstructionKernel::MAX_STACK_WORKSPACE_SIZE = 10;
2012
2017 : m_nEquations(0), m_nCoeffs(0), m_degree(0), m_dimensions(0)
2018{
2019 initialize(0, 0, 0);
2020}
2021
2029ReconstructionKernel::ReconstructionKernel(uint8_t degree, uint8_t dimensions, int nEquations)
2030 : m_nEquations(0), m_nCoeffs(0), m_degree(0), m_dimensions(0)
2031{
2032 initialize(degree, dimensions, nEquations);
2033}
2034
2042 : ReconstructionKernel(other.getDegree(), other.getDimensions(), other.m_nEquations)
2043{
2044 if (m_nEquations > 0) {
2045 int nWeights = m_nCoeffs * m_nEquations;
2046 std::copy(other.m_weights.get(), other.m_weights.get() + nWeights, m_weights.get());
2047 }
2048}
2049
2057{
2058 ReconstructionKernel tmp(other);
2059 swap(tmp);
2060
2061 return *this;
2062}
2063
2072{
2073 std::swap(other.m_degree, m_degree);
2074 std::swap(other.m_dimensions, m_dimensions);
2075 std::swap(other.m_nCoeffs, m_nCoeffs);
2076 std::swap(other.m_nEquations, m_nEquations);
2077 std::swap(other.m_weights, m_weights);
2078}
2079
2090void ReconstructionKernel::initialize(uint8_t degree, uint8_t dimensions, int nEquations, bool release)
2091{
2092 assert(degree <= ReconstructionPolynomial::MAX_DEGREE);
2093 m_degree = degree;
2094
2095 assert(dimensions <= ReconstructionPolynomial::MAX_DIMENSIONS);
2096 m_dimensions = dimensions;
2097
2098 int currentStorageSize = m_nCoeffs * m_nEquations;
2099
2100 m_nCoeffs = ReconstructionPolynomial::getCoefficientCount(m_degree, m_dimensions);
2101 m_nEquations = nEquations;
2102
2103 int storageSize = m_nCoeffs * m_nEquations;
2104
2105 bool reallocate;
2106 if (release) {
2107 reallocate = (currentStorageSize != storageSize);
2108 } else {
2109 reallocate = (currentStorageSize < storageSize);
2110 }
2111
2112 if (reallocate) {
2113 m_weights = std::unique_ptr<double[]>(new double[storageSize]);
2114 }
2115}
2116
2123{
2124 return m_degree;
2125}
2126
2133{
2134 return m_dimensions;
2135}
2136
2143{
2144 return m_nCoeffs;
2145}
2146
2153{
2154 return m_nEquations;
2155}
2156
2165{
2166 return m_weights.get();
2167}
2168
2177{
2178 return m_weights.get();
2179}
2180
2193void ReconstructionKernel::assemblePolynomial(const std::array<double, 3> &origin,
2194 const double *values,
2195 ReconstructionPolynomial *polynomial) const
2196{
2197 uint8_t degree = getDegree();
2198 uint8_t dimensions = getDimensions();
2199
2200 // Initialize the polynomial
2201 polynomial->initialize(degree, dimensions, origin, 1, true);
2202
2203 // Update the polynomial
2204 updatePolynomial(degree, values, polynomial);
2205}
2206
2225void ReconstructionKernel::assemblePolynomial(const std::array<double, 3> &origin,
2226 int nFields, const double **values,
2227 ReconstructionPolynomial *polynomial) const
2228{
2229 uint8_t degree = getDegree();
2230 uint8_t dimensions = getDimensions();
2231
2232 // Initialize the polynomial
2233 polynomial->initialize(degree, dimensions, origin, nFields, true);
2234
2235 // Update the polynomial
2236 updatePolynomial(degree, nFields, values, polynomial);
2237}
2238
2255void ReconstructionKernel::assemblePolynomial(uint8_t degree, const std::array<double, 3> &origin,
2256 const double *values,
2257 ReconstructionPolynomial *polynomial) const
2258{
2259 uint8_t dimensions = getDimensions();
2260
2261 // Initialize the polynomial
2262 polynomial->initialize(degree, dimensions, origin, 1, true);
2263
2264 // Update the polynomial
2265 updatePolynomial(degree, values, polynomial);
2266}
2267
2288void ReconstructionKernel::assemblePolynomial(uint8_t degree, const std::array<double, 3> &origin,
2289 int nFields, const double **values,
2290 ReconstructionPolynomial *polynomial) const
2291{
2292 uint8_t dimensions = getDimensions();
2293
2294 // Initialize the polynomial
2295 polynomial->initialize(degree, dimensions, origin, nFields, true);
2296
2297 // Update the polynomial
2298 updatePolynomial(degree, nFields, values, polynomial);
2299}
2300
2313 ReconstructionPolynomial *polynomial) const
2314{
2315 updatePolynomial(getDegree(), values, polynomial);
2316}
2317
2333void ReconstructionKernel::updatePolynomial(int nFields, const double **values,
2334 ReconstructionPolynomial *polynomial) const
2335{
2336 updatePolynomial(getDegree(), nFields, values, polynomial);
2337}
2338
2354void ReconstructionKernel::updatePolynomial(uint8_t degree, const double *values,
2355 ReconstructionPolynomial *polynomial) const
2356{
2357 assert(degree <= getDegree());
2358
2359 int nEquations = getEquationCount();
2360
2361 double *polynomialCoeffs = polynomial->getCoefficients();
2362
2363 // Constant polynomial
2364 if (ENABLE_FAST_PATH_OPTIMIZATIONS && degree == 0) {
2365 // Since the polynomial weights are stored using a col-major order,
2366 // the weights for the constant polynomial are the first nEquation
2367 // values.
2368 const double *constantPolynomialWeights = getPolynomialWeights();
2369
2370 polynomialCoeffs[0] = 0.;
2371 for (int j = 0; j < nEquations; ++j) {
2372 polynomialCoeffs[0] += values[j] * constantPolynomialWeights[j];
2373 }
2374
2375 return;
2376 }
2377
2378 // Generic polynomial
2379 uint8_t dimensions = getDimensions();
2380 int nCoeffs = ReconstructionPolynomial::getCoefficientCount(degree, dimensions);
2381
2382 const double *polynomialWeights = getPolynomialWeights();
2383
2384 cblas_dgemv(CBLAS_ORDER::CblasColMajor, CBLAS_TRANSPOSE::CblasTrans,
2385 nEquations, nCoeffs, 1., polynomialWeights, nEquations,
2386 values, 1, 0, polynomialCoeffs, 1);
2387}
2388
2409void ReconstructionKernel::updatePolynomial(uint8_t degree, int nFields, const double **values,
2410 ReconstructionPolynomial *polynomial) const
2411{
2412 assert(degree <= getDegree());
2413
2414 int nEquations = getEquationCount();
2415
2416 // Constant polynomial
2417 if (ENABLE_FAST_PATH_OPTIMIZATIONS && degree == 0) {
2418 // Since the polynomial weights are stored using a col-major order,
2419 // the weights for the constant polynomial are the first nEquation
2420 // values.
2421 const double *constantPolynomialWeights = getPolynomialWeights();
2422 for (int k = 0; k < nFields; ++k) {
2423 double *fieldCoeffs = polynomial->getCoefficients(k);
2424
2425 fieldCoeffs[0] = 0.;
2426 for (int j = 0; j < nEquations; ++j) {
2427 fieldCoeffs[0] += values[j][k] * constantPolynomialWeights[j];
2428 }
2429 }
2430
2431 return;
2432 }
2433
2434 // Generic polynomial
2435 BITPIT_CREATE_WORKSPACE(fieldValues, double, nEquations, MAX_STACK_WORKSPACE_SIZE);
2436
2437 uint8_t dimensions = getDimensions();
2438 int nCoeffs = ReconstructionPolynomial::getCoefficientCount(degree, dimensions);
2439
2440 const double *polynomialWeights = getPolynomialWeights();
2441
2442 for (int k = 0; k < nFields; ++k) {
2443 for (int j = 0; j < nEquations; ++j) {
2444 fieldValues[j] = values[j][k];
2445 }
2446
2447 double *fieldCoeffs = polynomial->getCoefficients(k);
2448
2449 cblas_dgemv(CBLAS_ORDER::CblasColMajor, CBLAS_TRANSPOSE::CblasTrans,
2450 nEquations, nCoeffs, 1., polynomialWeights, nEquations,
2451 fieldValues, 1, 0, fieldCoeffs, 1);
2452 }
2453}
2454
2469void ReconstructionKernel::computeValueWeights(const std::array<double, 3> &origin,
2470 const std::array<double, 3> &point, double *valueWeights) const
2471{
2472 computeValueLimitedWeights(getEquationCount(), getDegree(), origin, point, nullptr, valueWeights);
2473}
2474
2493void ReconstructionKernel::computeValueWeights(uint8_t degree, const std::array<double, 3> &origin,
2494 const std::array<double, 3> &point, double *valueWeights) const
2495{
2496 computeValueLimitedWeights(getEquationCount(), degree, origin, point, nullptr, valueWeights);
2497}
2498
2499
2519void ReconstructionKernel::computeValueWeights(int nEquations, uint8_t degree, const std::array<double, 3> &origin,
2520 const std::array<double, 3> &point, double *valueWeights) const
2521{
2522 computeValueLimitedWeights(nEquations, degree, origin, point, nullptr, valueWeights);
2523}
2524
2546void ReconstructionKernel::computeValueLimitedWeights(const std::array<double, 3> &origin, const std::array<double, 3> &point,
2547 const double *limiters, double *valueWeights) const
2548{
2549 computeValueLimitedWeights(getEquationCount(), getDegree(), origin, point, limiters, valueWeights);
2550}
2551
2574void ReconstructionKernel::computeValueLimitedWeights(uint8_t degree, const std::array<double, 3> &origin,
2575 const std::array<double, 3> &point, const double *limiters,
2576 double *valueWeights) const
2577{
2578 computeValueLimitedWeights(getEquationCount(), degree, origin, point, limiters, valueWeights);
2579}
2580
2607void ReconstructionKernel::computeValueLimitedWeights(int nEquations, uint8_t degree, const std::array<double, 3> &origin,
2608 const std::array<double, 3> &point, const double *limiters,
2609 double *valueWeights) const
2610{
2611 assert(nEquations <= getEquationCount());
2612 assert(degree <= getDegree());
2613
2614 // Constant polynomial
2615 if (ENABLE_FAST_PATH_OPTIMIZATIONS && degree == 0) {
2616 // Since the polynomial weights are stored using a col-major order,
2617 // the weights for the constant polynomial are the first nEquation
2618 // values.
2619 const double *constantPolynomialWeights = getPolynomialWeights();
2620 std::copy_n(constantPolynomialWeights, nEquations, valueWeights);
2621
2622 return;
2623 }
2624
2625 // Evaluate basis
2626 uint8_t dimensions = getDimensions();
2627 int nCoeffs = ReconstructionPolynomial::getCoefficientCount(degree, dimensions);
2628
2629 BITPIT_CREATE_WORKSPACE(csi, double, nCoeffs, MAX_STACK_WORKSPACE_SIZE);
2630 ReconstructionPolynomial::evalPointBasisValues(degree, dimensions, origin, point, csi);
2631 if (limiters) {
2632 applyLimiter(degree, limiters, csi);
2633 }
2634
2635 // Generic polynomial
2636 int polynomialWeightsRowCount = getEquationCount();
2637 const double *polynomialWeights = getPolynomialWeights();
2638
2639 cblas_dgemv(CBLAS_ORDER::CblasColMajor, CBLAS_TRANSPOSE::CblasNoTrans,
2640 nEquations, nCoeffs, 1., polynomialWeights, polynomialWeightsRowCount,
2641 csi, 1, 0, valueWeights, 1);
2642}
2643
2660void ReconstructionKernel::computeDerivativeWeights(const std::array<double, 3> &origin,
2661 const std::array<double, 3> &point, const std::array<double, 3> &direction,
2662 double *derivativeWeights) const
2663{
2664 computeDerivativeLimitedWeights(getEquationCount(), getDegree(), origin, point, direction, nullptr, derivativeWeights);
2665}
2666
2687void ReconstructionKernel::computeDerivativeWeights(uint8_t degree, const std::array<double, 3> &origin,
2688 const std::array<double, 3> &point, const std::array<double, 3> &direction,
2689 double *derivativeWeights) const
2690{
2691 computeDerivativeLimitedWeights(getEquationCount(), degree, origin, point, direction, nullptr, derivativeWeights);
2692}
2693
2715void ReconstructionKernel::computeDerivativeWeights(int nEquations, uint8_t degree, const std::array<double, 3> &origin,
2716 const std::array<double, 3> &point, const std::array<double, 3> &direction,
2717 double *derivativeWeights) const
2718{
2719 computeDerivativeLimitedWeights(nEquations, degree, origin, point, direction, nullptr, derivativeWeights);
2720}
2721
2745void ReconstructionKernel::computeDerivativeLimitedWeights(const std::array<double, 3> &origin,
2746 const std::array<double, 3> &point, const std::array<double, 3> &direction,
2747 const double *limiters,
2748 double *derivativeWeights) const
2749{
2750 computeDerivativeLimitedWeights(getEquationCount(), getDegree(), origin, point, direction, limiters, derivativeWeights);
2751}
2752
2780void ReconstructionKernel::computeDerivativeLimitedWeights(uint8_t degree, const std::array<double, 3> &origin,
2781 const std::array<double, 3> &point, const std::array<double, 3> &direction,
2782 const double *limiters,
2783 double *derivativeWeights) const
2784{
2785 computeDerivativeLimitedWeights(getEquationCount(), degree, origin, point, direction, limiters, derivativeWeights);
2786}
2787
2816void ReconstructionKernel::computeDerivativeLimitedWeights(int nEquations, uint8_t degree, const std::array<double, 3> &origin,
2817 const std::array<double, 3> &point, const std::array<double, 3> &direction,
2818 const double *limiters,
2819 double *derivativeWeights) const
2820{
2821 assert(nEquations <= getEquationCount());
2822 assert(degree <= getDegree());
2823
2824 // Constant polynomial
2825 if (ENABLE_FAST_PATH_OPTIMIZATIONS && degree == 0) {
2826 std::fill_n(derivativeWeights, nEquations, 0.);
2827
2828 return;
2829 }
2830
2831 // Evaluate basis
2832 uint8_t dimensions = getDimensions();
2833 int nCoeffs = ReconstructionPolynomial::getCoefficientCount(degree, dimensions);
2834
2835 BITPIT_CREATE_WORKSPACE(dcsi, double, nCoeffs, MAX_STACK_WORKSPACE_SIZE);
2836 ReconstructionPolynomial::evalPointBasisDerivatives(degree, dimensions, origin, point, direction, dcsi);
2837 if (limiters) {
2838 applyLimiter(degree, limiters, dcsi);
2839 }
2840
2841 // Generic polynomial
2842 int polynomialWeightsRowCount = getEquationCount();
2843 const double *polynomialWeights = getPolynomialWeights();
2844
2845 cblas_dgemv(CBLAS_ORDER::CblasColMajor, CBLAS_TRANSPOSE::CblasNoTrans,
2846 nEquations, nCoeffs, 1., polynomialWeights, polynomialWeightsRowCount,
2847 dcsi, 1, 0, derivativeWeights, 1);
2848}
2849
2864void ReconstructionKernel::computeGradientWeights(const std::array<double, 3> &origin,
2865 const std::array<double, 3> &point,
2866 std::array<double, 3> *gradientWeights) const
2867{
2868 computeGradientLimitedWeights(getEquationCount(), getDegree(), origin, point, nullptr, gradientWeights);
2869}
2870
2889void ReconstructionKernel::computeGradientWeights(uint8_t degree, const std::array<double, 3> &origin,
2890 const std::array<double, 3> &point,
2891 std::array<double, 3> *gradientWeights) const
2892{
2893 computeGradientLimitedWeights(getEquationCount(), degree, origin, point, nullptr, gradientWeights);
2894}
2895
2915void ReconstructionKernel::computeGradientWeights(int nEquations, uint8_t degree, const std::array<double, 3> &origin,
2916 const std::array<double, 3> &point,
2917 std::array<double, 3> *gradientWeights) const
2918{
2919 computeGradientLimitedWeights(nEquations, degree, origin, point, nullptr, gradientWeights);
2920}
2921
2942void ReconstructionKernel::computeGradientLimitedWeights(const std::array<double, 3> &origin,
2943 const std::array<double, 3> &point, const double *limiters,
2944 std::array<double, 3> *gradientWeights) const
2945{
2946 computeGradientLimitedWeights(getEquationCount(), getDegree(), origin, point, limiters, gradientWeights);
2947}
2948
2973void ReconstructionKernel::computeGradientLimitedWeights(uint8_t degree, const std::array<double, 3> &origin,
2974 const std::array<double, 3> &point, const double *limiters,
2975 std::array<double, 3> *gradientWeights) const
2976{
2977 computeGradientLimitedWeights(getEquationCount(), degree, origin, point, limiters, gradientWeights);
2978}
2979
3005void ReconstructionKernel::computeGradientLimitedWeights(int nEquations, uint8_t degree, const std::array<double, 3> &origin,
3006 const std::array<double, 3> &point, const double *limiters,
3007 std::array<double, 3> *gradientWeights) const
3008{
3009 assert(nEquations <= getEquationCount());
3010 assert(degree <= getDegree());
3011
3012 // Constant polynomial
3013 if (ENABLE_FAST_PATH_OPTIMIZATIONS && degree == 0) {
3014 std::fill_n(gradientWeights, nEquations, std::array<double, 3>{0., 0., 0.});
3015
3016 return;
3017 }
3018
3019 // Generic polynomial
3020 uint8_t dimensions = getDimensions();
3021 int nCoeffs = ReconstructionPolynomial::getCoefficientCount(degree, dimensions);
3022
3023 int polynomialWeightsRowCount = getEquationCount();
3024 const double *polynomialWeights = getPolynomialWeights();
3025
3026 BITPIT_CREATE_WORKSPACE(dcsi, double, nCoeffs * dimensions, 3 * MAX_STACK_WORKSPACE_SIZE);
3027 for (int d = 0; d < dimensions; ++d) {
3028 // Select derivative direction
3029 std::array<double, 3> direction = {{0., 0., 0.}};
3030 direction[d] = 1.;
3031
3032 // Evaluate basis derivatives
3033 int offset = linearalgebra::linearIndexColMajor(0, d, nCoeffs, dimensions);
3034
3035 ReconstructionPolynomial::evalPointBasisDerivatives(degree, dimensions, origin, point, direction, dcsi + offset);
3036 if (limiters) {
3037 applyLimiter(degree, limiters, dcsi + offset);
3038 }
3039 }
3040
3041 cblas_dgemm(CBLAS_ORDER::CblasColMajor, CBLAS_TRANSPOSE::CblasTrans, CBLAS_TRANSPOSE::CblasTrans,
3042 dimensions, nEquations, nCoeffs, 1., dcsi, nCoeffs, polynomialWeights, polynomialWeightsRowCount,
3043 0., gradientWeights[0].data(), 3);
3044
3045 // Explicitly zero unused components
3046 if (dimensions != ReconstructionPolynomial::MAX_DIMENSIONS) {
3047 for (int j = 0; j < nEquations; ++j) {
3048 for (int d = dimensions; d < ReconstructionPolynomial::MAX_DIMENSIONS; ++d) {
3049 gradientWeights[j][d] = 0.;
3050 }
3051 }
3052 }
3053}
3054
3067void ReconstructionKernel::applyLimiter(uint8_t degree, const double *limiters, double *coeffs) const
3068{
3069 if (degree < 1) {
3070 return;
3071 }
3072
3073 uint8_t dimensions = getDimensions();
3074 const std::vector<uint16_t> &nDegreeCoeffs = ReconstructionPolynomial::getDegreeCoefficientsCount(dimensions);
3075
3076 int coeffEnd = nDegreeCoeffs[0];
3077 for (int n = 1; n <= degree; ++n) {
3078 double limiterScaleFactor = limiters[n - 1];
3079
3080 int coeffBegin = coeffEnd;
3081 coeffEnd = coeffBegin + nDegreeCoeffs[n];
3082 for (int k = coeffBegin; k < coeffEnd; ++k) {
3083 coeffs[k] *= limiterScaleFactor;
3084 }
3085 }
3086}
3087
3095void ReconstructionKernel::display(std::ostream &out, double tolerance) const
3096{
3097 int nCoeffs = getCoefficientCount();
3098 for (int i = 0; i < nCoeffs; ++i) {
3099 out << i << " ";
3100 for (int j = 0; j < m_nEquations; ++j) {
3101 double weigth = m_weights[linearalgebra::linearIndexColMajor(j, i, m_nEquations, m_nCoeffs)];
3102 if (std::abs(weigth) < tolerance) {
3103 continue;
3104 }
3105
3106 out << "(" << j << "," << weigth << ") ";
3107 }
3108
3109 out << std::endl;
3110 }
3111}
3112
3124const double ReconstructionAssembler::SVD_ZERO_THRESHOLD = 1e-14;
3125
3126
3134
3141ReconstructionAssembler::ReconstructionAssembler(uint8_t degree, uint8_t dimensions)
3142{
3143 initialize(degree, dimensions, false);
3144}
3145
3154{
3155 std::swap(other.m_degree, m_degree);
3156 std::swap(other.m_dimensions, m_dimensions);
3157 std::swap(other.m_nCoeffs, m_nCoeffs);
3158 std::swap(other.m_constraintsOrder, m_constraintsOrder);
3159 std::swap(other.m_leastSquaresOrder, m_leastSquaresOrder);
3160 std::swap(other.m_leastSquaresScaleFactors, m_leastSquaresScaleFactors);
3161 std::swap(other.m_A, m_A);
3162 std::swap(other.m_C, m_C);
3163 std::swap(other.m_sigma, m_sigma);
3164 std::swap(other.m_U, m_U);
3165 std::swap(other.m_S, m_S);
3166 std::swap(other.m_Vt, m_Vt);
3167 std::swap(other.m_SVDWorkspace, m_SVDWorkspace);
3168 std::swap(other.m_w, m_w);
3169}
3170
3180void ReconstructionAssembler::initialize(uint8_t degree, uint8_t dimensions, bool release)
3181{
3182 assert(degree <= ReconstructionPolynomial::MAX_DEGREE);
3183 m_degree = degree;
3184
3185 assert(dimensions <= ReconstructionPolynomial::MAX_DIMENSIONS);
3186 m_dimensions = dimensions;
3187
3188 m_nCoeffs = ReconstructionPolynomial::getCoefficientCount(m_degree, m_dimensions);
3189
3190 clear(release);
3191}
3192
3200{
3201 m_constraintsOrder.clear();
3202 m_leastSquaresOrder.clear();
3203 m_leastSquaresScaleFactors.clear();
3204
3205 m_A.clear();
3206 m_C.clear();
3207
3208 m_sigma.clear();
3209 m_U.clear();
3210 m_S.clear();
3211 m_Vt.clear();
3212 m_SVDWorkspace.resize(1);
3213 m_w.clear();
3214
3215 if (release) {
3216 m_constraintsOrder.shrink_to_fit();
3217 m_leastSquaresOrder.shrink_to_fit();
3218 m_leastSquaresScaleFactors.shrink_to_fit();
3219
3220 m_A.shrink_to_fit();
3221 m_C.shrink_to_fit();
3222
3223 m_sigma.shrink_to_fit();
3224 m_U.shrink_to_fit();
3225 m_S.shrink_to_fit();
3226 m_Vt.shrink_to_fit();
3227 m_SVDWorkspace.shrink_to_fit();
3228 m_w.shrink_to_fit();
3229 }
3230}
3231
3238{
3239 return m_degree;
3240}
3241
3248{
3249 return m_dimensions;
3250}
3251
3258{
3259 return m_nCoeffs;
3260}
3261
3268{
3269 return m_constraintsOrder.size();
3270}
3271
3278{
3279 return m_leastSquaresOrder.size();
3280}
3281
3288{
3289 int nConstraints = countConstraints();
3290 int nLeastSquares = countLeastSquares();
3291 int nEquations = nConstraints + nLeastSquares;
3292
3293 return nEquations;
3294}
3295
3305 const std::array<double, 3> &origin,
3306 const std::array<double, 3> &point,
3307 double scaleFactor)
3308{
3309 double *equationCoeffs = _addEquation(type, scaleFactor);
3310 ReconstructionPolynomial::evalPointBasisValues(getDegree(), getDimensions(), origin, point, equationCoeffs);
3311}
3312
3323 const std::array<double, 3> &origin,
3324 const std::array<double, 3> &point,
3325 const std::array<double, 3> &direction,
3326 double scaleFactor)
3327{
3328 double *equationCoeffs = _addEquation(type, scaleFactor);
3329 ReconstructionPolynomial::evalPointBasisDerivatives(getDegree(), getDimensions(), origin, point, direction, equationCoeffs);
3330}
3331
3342 const Cell &cell,
3343 const std::array<double, 3> &origin,
3344 const std::array<double, 3> *vertexCoords,
3345 double scaleFactor)
3346{
3347 double *equationCoeffs = _addEquation(type, scaleFactor);
3348 ReconstructionPolynomial::evalCellBasisValues(getDegree(), getDimensions(), origin, cell, vertexCoords, equationCoeffs);
3349}
3350
3357double * ReconstructionAssembler::_addEquation(ReconstructionType type, double scaleFactor)
3358{
3359 // Update equation information
3360 int nEquations = countEquations();
3361 switch (type) {
3362
3363 case TYPE_CONSTRAINT:
3364 m_constraintsOrder.emplace_back(nEquations);
3365 break;
3366
3367 case TYPE_LEAST_SQUARE:
3368 m_leastSquaresOrder.emplace_back(nEquations);
3369 m_leastSquaresScaleFactors.emplace_back(scaleFactor);
3370 break;
3371
3372 }
3373
3374 // Prepare storage for equation coefficients
3375 int nCoeffs = getCoefficientCount();
3376
3377 double *equationCoeffsStorage = nullptr;
3378 switch (type) {
3379
3380 case TYPE_CONSTRAINT:
3381 m_C.resize(m_C.size() + nCoeffs);
3382 equationCoeffsStorage = m_C.data() + m_C.size() - nCoeffs;
3383 break;
3384
3385 case TYPE_LEAST_SQUARE:
3386 m_A.resize(m_A.size() + nCoeffs);
3387 equationCoeffsStorage = m_A.data() + m_A.size() - nCoeffs;
3388 break;
3389
3390 }
3391
3392 return equationCoeffsStorage;
3393}
3394
3404{
3405 // Initialize reconstruction kernel
3406 uint8_t degree = getDegree();
3407 uint8_t dimensions = getDimensions();
3408
3409 int nEquations = countEquations();
3410
3411 kernel->initialize(degree, dimensions, nEquations, true);
3412
3413 // Update the kernel
3414 updateKernel(kernel);
3415}
3416
3425{
3426 // Get the number of equations
3427 int nConstraints = countConstraints();
3428 int nLeastSquares = countLeastSquares();
3429 int nEquations = countEquations();
3430
3431 // Get the number of polynomial coefficients
3432 int nCoeffs = getCoefficientCount();
3433
3434 // Evaluate normalized least square scale factors
3435 if (nLeastSquares > 0) {
3436 double maxLeastSquareScaleFactor = std::abs(m_leastSquaresScaleFactors[0]);
3437 for (int k = 1; k < nLeastSquares; ++k) {
3438 maxLeastSquareScaleFactor = std::max(std::abs(m_leastSquaresScaleFactors[k]), maxLeastSquareScaleFactor);
3439 }
3440
3441 m_w.resize(nLeastSquares);
3442 for (int k = 0; k < nLeastSquares; ++k) {
3443 m_w[k] = m_leastSquaresScaleFactors[k] / maxLeastSquareScaleFactor;
3444 }
3445 }
3446
3447 // The linear-constrained are introduced in the least-squares problem
3448 // through Lagrange multipliers. The resulting linear system is:
3449 //
3450 // | A^t A w C^t | |x | = |A^t w b|
3451 // | C 0 | |lambda| |d |
3452 //
3453 // with A and C the least-squares and constraints equations respectively,
3454 // and b and d their corresponding RHSs. x are the coefficients of
3455 // the polynomial and lambda the lagrange multipliers. w are the normalized
3456 // east square scale factors.
3457 //
3458 // This system is denoted by S:
3459 //
3460 // | x | |A^t w 0| |b|
3461 // |S| | | = | | | |
3462 // |lambda| |0 I| |d|
3463 //
3464 // The matrices S and S^-1 are symmetric and only the upper portions are
3465 // computed
3466 int nUnknowns = nCoeffs + nConstraints;
3467
3468 m_S.resize(nUnknowns * nUnknowns);
3469 for (int i = 0; i < nCoeffs; ++i) {
3470 for (int j = i; j < nCoeffs; ++j) {
3471 // Compute A^t A on the fly
3472 double ATA_ij = 0.;
3473 for (int k = 0; k < nLeastSquares; ++k) {
3474 int A_ki_idx = linearalgebra::linearIndexRowMajor(k, i, nLeastSquares, nCoeffs);
3475 int A_kj_idx = linearalgebra::linearIndexRowMajor(k, j, nLeastSquares, nCoeffs);
3476 ATA_ij += m_A[A_ki_idx] * m_A[A_kj_idx] * m_w[k];
3477 }
3478
3479 int l = linearalgebra::linearIndexColMajor(i, j, nUnknowns, nUnknowns);
3480 m_S[l] = ATA_ij;
3481
3482 int m = linearalgebra::linearIndexColMajor(j, i, nUnknowns, nUnknowns);
3483 m_S[m] = ATA_ij;
3484 }
3485
3486 for (int j = nCoeffs; j < nUnknowns; ++j) {
3487 int l = linearalgebra::linearIndexColMajor(i, j, nUnknowns, nUnknowns);
3488 int C_idx = linearalgebra::linearIndexRowMajor(j - nCoeffs, i, nConstraints, nCoeffs);
3489 m_S[l] = m_C[C_idx];
3490
3491 int m = linearalgebra::linearIndexColMajor(j, i, nUnknowns, nUnknowns);
3492 m_S[m] = m_S[l];
3493 }
3494 }
3495
3496 // Compute inverse S matrix
3497 // Since S may me be rank-deficit (eg if not enough neighbours are available)
3498 // the pseudo-inverse is used. This corresponds of computing the least-norm
3499 // solution of the problem.
3500 computePseudoInverse(nUnknowns, nUnknowns, SVD_ZERO_THRESHOLD, m_S.data());
3501
3502 // Weights needed to evaluate the polynomial coefficients come from the
3503 // following equation:
3504 //
3505 // | x | |A^t w 0| |b| |b|
3506 // | | = S^-1 | | | | = S^-1 Q | |
3507 // |lambda| |0 I| |d| |d|
3508 //
3509 // Since we are interested only in x (the polynomial coefficients) only
3510 // the first nCoeffs rows of the matrix S^-1 Q are computed. Those values
3511 // are the polynomial weights.
3512 //
3513 // Weights are stored according the order in which the equations have been
3514 // added.
3515 double *weights = kernel->getPolynomialWeights();
3516 for (int j = 0; j < nEquations; ++j) {
3517 int equation;
3518 if (j < nLeastSquares) {
3519 equation = m_leastSquaresOrder[j];
3520 } else {
3521 equation = m_constraintsOrder[j - nLeastSquares];
3522 }
3523
3524 for (int i = 0; i < nCoeffs; ++i) {
3525 double value = 0;
3526 for (int k = 0; k < nUnknowns; ++k) {
3527 int l = linearalgebra::linearIndexColMajorSymmetric(i, k, nUnknowns, nUnknowns, 'U');
3528 if (k < nCoeffs && j < nLeastSquares) {
3529 int A_jk_idx = linearalgebra::linearIndexRowMajor(j, k, nLeastSquares, nCoeffs);
3530 value += m_S[l] * m_A[A_jk_idx] * m_w[j];
3531 } else if ((k - nCoeffs) == (j - nLeastSquares)) {
3532 value += m_S[l];
3533 }
3534 }
3535
3536 int weightLineraIndex = linearalgebra::linearIndexColMajor(equation, i, nEquations, nCoeffs);
3537 weights[weightLineraIndex] = value;
3538 }
3539 }
3540}
3554void ReconstructionAssembler::computePseudoInverse(int m, int n, double zeroThreshold, double *A) const
3555{
3556 // Compute SVD
3557 //
3558 // A = U * Sigma * Vt (Equation 21)
3559 int k = std::min(m,n);
3560
3561 m_sigma.resize(k);
3562 m_U.resize(m * k);
3563 m_Vt.resize(k * n);
3564
3565 char jobU = 'S';
3566 char jobVT = 'S';
3567
3568 int workspaceSize = -1;
3569
3570 int info;
3571
3572 info = LAPACKE_dgesvd_work(LAPACK_COL_MAJOR, jobU, jobVT, m, n, A, m, m_sigma.data(), m_U.data(), m, m_Vt.data(), k,
3573 m_SVDWorkspace.data(), workspaceSize);
3574
3575 workspaceSize = static_cast<int>(m_SVDWorkspace[0]);
3576 if (workspaceSize > (int) m_SVDWorkspace.size()) {
3577 m_SVDWorkspace.resize(workspaceSize);
3578 }
3579
3580 info = LAPACKE_dgesvd_work(LAPACK_COL_MAJOR, jobU, jobVT, m, n, A, m, m_sigma.data(), m_U.data(), m, m_Vt.data(), k,
3581 m_SVDWorkspace.data(), workspaceSize);
3582
3583 if (info > 0) {
3584 log::cout() << "SVD failed in ReconstructionAssembler::computePseudoInverse()" <<std::endl;
3585 exit(1);
3586 }
3587
3588 // Inv(A) = V * Sigma^ + *U^T (Equation 22)
3589 //
3590 // u = sigma^ + *U
3591 // and is stored in U
3592 //
3593 // Sigma are singular values, hence they are non-negative by definition.
3594 // To check if a signular value is non zero, we only need to check if it's
3595 // greater than the defined threshold.
3596 for (int i = 0; i < k; ++i) {
3597 double sigma_plus = (m_sigma[i] > zeroThreshold) ? (1. / m_sigma[i]) : 0.;
3598 cblas_dscal(m, sigma_plus, &m_U[i*m], 1);
3599 }
3600
3601 // Inv(A) = (Vt)^T * u^T
3602 cblas_dgemm(CBLAS_ORDER::CblasColMajor, CBLAS_TRANSPOSE::CblasTrans, CBLAS_TRANSPOSE::CblasTrans,
3603 n, m, k, 1., m_Vt.data(), k, m_U.data(), m, 0., A, n);
3604}
3605
3620Reconstruction::Reconstruction(uint8_t degree, uint8_t dimensions)
3621 : ReconstructionAssembler(degree, dimensions)
3622{
3623}
3624
3633{
3636}
3637
3647void Reconstruction::initialize(uint8_t degree, uint8_t dimensions, bool release)
3648{
3649 ReconstructionAssembler::initialize(degree, dimensions, release);
3650}
3651
3659void Reconstruction::clear(bool release)
3660{
3662
3663 if (release) {
3664 ReconstructionKernel().swap(*this);
3665 }
3666}
3667
3681
3682}
The Cell class defines the cells.
Definition cell.hpp:42
ElementType getType() const
Definition element.cpp:551
double evalSize(const std::array< double, 3 > *coordinates) const
Definition element.cpp:1549
std::array< double, 3 > evalCentroid(const std::array< double, 3 > *coordinates) const
Definition element.cpp:1521
The ReconstructionAssembler class allows to define reconstruction polynomial.
void swap(ReconstructionAssembler &other) noexcept
void addPointValueEquation(ReconstructionType type, const std::array< double, 3 > &origin, const std::array< double, 3 > &point, double scaleFactor=1.)
void assembleKernel(ReconstructionKernel *kernel) const
void initialize(uint8_t degree, uint8_t dimensions, bool release=true)
void addPointDerivativeEquation(ReconstructionType type, const std::array< double, 3 > &origin, const std::array< double, 3 > &point, const std::array< double, 3 > &direction, double scaleFactor=1.)
void addCellAverageEquation(ReconstructionType type, const Cell &cell, const std::array< double, 3 > &origin, const std::array< double, 3 > *vertexCoords, double scaleFactor=1.)
void updateKernel(ReconstructionKernel *kernel) const
The ReconstructionKernel class allows to evaluate the weight of a reconstruction polynomial previousl...
void computeGradientWeights(const std::array< double, 3 > &origin, const std::array< double, 3 > &point, std::array< double, 3 > *gradientWeights) const
void computeValueWeights(const std::array< double, 3 > &origin, const std::array< double, 3 > &point, double *valueWeights) const
void computeGradientLimitedWeights(const std::array< double, 3 > &origin, const std::array< double, 3 > &point, const double *limiters, std::array< double, 3 > *gradientWeights) const
void updatePolynomial(const double *values, ReconstructionPolynomial *polynomial) const
void computeDerivativeLimitedWeights(const std::array< double, 3 > &origin, const std::array< double, 3 > &point, const std::array< double, 3 > &direction, const double *limiters, double *derivativeWeights) const
void swap(ReconstructionKernel &other) noexcept
void initialize(uint8_t degree, uint8_t dimensions, int nEquations, bool release=true)
void display(std::ostream &out, double tolerance=1.e-10) const
void computeValueLimitedWeights(const std::array< double, 3 > &origin, const std::array< double, 3 > &point, const double *limiters, double *valueWeights) const
void computeDerivativeWeights(const std::array< double, 3 > &origin, const std::array< double, 3 > &point, const std::array< double, 3 > &direction, double *derivativeWeights) const
ReconstructionKernel & operator=(const ReconstructionKernel &other)
void assemblePolynomial(const std::array< double, 3 > &origin, const double *values, ReconstructionPolynomial *polynomial) const
const double * getPolynomialWeights() const
void applyLimiter(uint8_t degree, const double *limiters, double *coeffs) const
The ReconstructionPolynomial class allows to apply a reconstruction polynomial previously assembled.
static uint16_t countCoefficients(uint8_t degree, uint8_t dimensions)
static void evalPointBasisDerivatives(uint8_t degree, uint8_t dimensions, const std::array< double, 3 > &origin, const std::array< double, 3 > &point, const std::array< double, 3 > &direction, double *dcsi)
void computeDerivativeLimited(const std::array< double, 3 > &point, const std::array< double, 3 > &direction, const double *limiters, int field, double *derivative) const
void computeValue(const std::array< double, 3 > &point, int field, double *values) const
const double * getDegreeCoefficients(uint8_t degree, int field=0) const
void computeDerivativesLimited(const std::array< double, 3 > &point, const std::array< double, 3 > &direction, const double *limiters, double *derivatives) const
void display(std::ostream &out) const
void computeValuesLimited(const std::array< double, 3 > &point, const double *limiters, double *values) const
void initialize(uint8_t degree, uint8_t dimensions, const std::array< double, 3 > &origin, int nFields=1, bool release=true)
void swap(ReconstructionPolynomial &other) noexcept
const double * getCoefficients() const
static void evalPointBasisValues(uint8_t degree, uint8_t dimensions, const std::array< double, 3 > &origin, const std::array< double, 3 > &point, double *csi)
void computeGradientsLimited(const std::array< double, 3 > &point, const double *limiters, std::array< double, 3 > *gradients) const
static const std::vector< uint16_t > & getCoefficientsCount(uint8_t dimensions)
void computeGradients(const std::array< double, 3 > &point, std::array< double, 3 > *gradients) const
void computeDerivative(const std::array< double, 3 > &point, const std::array< double, 3 > &direction, int field, double *derivative) const
void computeGradientLimited(const std::array< double, 3 > &point, const double *limiters, int field, std::array< double, 3 > *gradient) const
static void evalCellBasisValues(uint8_t degree, uint8_t dimensions, const std::array< double, 3 > &origin, const Cell &cell, const std::array< double, 3 > *vertexCoords, double *csi)
static uint16_t getDegreeCoefficientCount(uint8_t degree, uint8_t dimensions)
void computeValueLimited(const std::array< double, 3 > &point, const double *limiters, int field, double *value) const
void computeDerivatives(const std::array< double, 3 > &point, const std::array< double, 3 > &direction, double *derivatives) const
void computeValues(const std::array< double, 3 > &point, double *values) const
ReconstructionPolynomial & operator=(const ReconstructionPolynomial &other)
const std::array< double, 3 > & getOrigin() const
static uint16_t countDegreeCoefficients(uint8_t degree, uint8_t dimensions)
static const std::vector< uint16_t > & getDegreeCoefficientsCount(uint8_t dimensions)
void computeGradient(const std::array< double, 3 > &point, int field, std::array< double, 3 > *gradient) const
The Reconstruction class allows to build and apply a polynomial reconstructions.
void clear(bool release=true)
void initialize(uint8_t degree, uint8_t dimensions, bool release=true)
Reconstruction(uint8_t degree, uint8_t dimensions)
void swap(Reconstruction &other) noexcept
#define BITPIT_CREATE_WORKSPACE(workspace, item_type, size, stack_size)
unsigned long factorial(unsigned long n)
int linearIndexRowMajor(int row, int col, int nRows, int nCols)
int linearIndexColMajorSymmetric(int row, int col, int nRows, int nCols, char uplo)
int linearIndexColMajor(int row, int col, int nRows, int nCols)
Logger & cout(log::Level defaultSeverity, log::Visibility defaultVisibility)
Definition logger.cpp:1705
--- layout: doxygen_footer ---