27#include "bitpit_LA.hpp"
28#include "bitpit_operators.hpp"
30#include "bitpit_private_cblas.hpp"
31#include "bitpit_private_lapacke.hpp"
33#include "reconstruction.hpp"
58const bool ReconstructionPolynomial::ENABLE_FAST_PATH_OPTIMIZATIONS =
true;
63const int ReconstructionPolynomial::MAX_STACK_WORKSPACE_SIZE = 10;
68const std::vector<std::vector<uint16_t>> ReconstructionPolynomial::m_countCoefficientCache = generateCountCoefficientCache();
74const std::vector<std::vector<uint16_t>> ReconstructionPolynomial::m_countDegreeCoefficientCache = generateCountDegreeCoefficientCache();
79std::vector<std::vector<uint16_t>> ReconstructionPolynomial::generateCountCoefficientCache()
82 for (uint8_t dimensions = 0; dimensions <=
MAX_DIMENSIONS; ++dimensions) {
83 for (uint8_t degree = 0; degree <=
MAX_DEGREE; ++degree) {
95std::vector<std::vector<uint16_t>> ReconstructionPolynomial::generateCountDegreeCoefficientCache()
98 for (uint8_t dimensions = 0; dimensions <=
MAX_DIMENSIONS; ++dimensions) {
99 for (uint8_t degree = 0; degree <=
MAX_DEGREE; ++degree) {
121 return m_countCoefficientCache[dimensions][degree];
136 return m_countCoefficientCache[dimensions];
150 uint16_t nCoeffs = 0;
151 for (
int i = 0; i <= degree; ++i) {
172 return m_countDegreeCoefficientCache[dimensions][degree];
187 return m_countDegreeCoefficientCache[dimensions];
201 if (dimensions == 0) {
219 const std::array<double, 3> &point,
double *csi)
227 const std::array<double, 3> distance = point - origin;
230 for (
int i = 0; i < dimensions; ++i) {
231 csi[offset++] = distance[i];
236 for (
int i = 0; i < dimensions; ++i) {
237 csi[offset++] = 0.5 * distance[i] * distance[i];
240 if (dimensions >= 2) {
241 csi[offset++] = distance[0] * distance[1];
243 if (dimensions >= 3) {
244 csi[offset++] = distance[0] * distance[2];
245 csi[offset++] = distance[1] * distance[2];
268 const std::array<double, 3> &point,
const std::array<double, 3> &direction,
279 for (
int i = 0; i < dimensions; ++i) {
280 dcsi[offset++] = direction[i];
285 const std::array<double, 3> distance = point - origin;
287 for (
int i = 0; i < dimensions; ++i) {
288 dcsi[offset++] = distance[i] * direction[i];
291 if (dimensions >= 2) {
292 dcsi[offset++] = distance[0] * direction[1] + distance[1] * direction[0];
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];
339 const Cell &cell,
const std::array<double, 3> *vertexCoords,
356 bool cellTypeSupported;
357 if (cellType == ElementType::PIXEL) {
358 cellTypeSupported = (degree <= 2);
359 }
else if (cellType == ElementType::VOXEL) {
360 cellTypeSupported = (degree <= 2);
362 cellTypeSupported = (degree <= 1);
365 if (!cellTypeSupported) {
366 throw std::runtime_error(
"Cell type not supported.");
375 const std::array<double, 3> distance = cell.
evalCentroid(vertexCoords) - origin;
378 for (
int i = 0; i < dimensions; ++i) {
379 csi[offset++] = distance[i];
384 double cellSize = cell.
evalSize(vertexCoords);
386 for (
int i = 0; i < dimensions; ++i) {
387 csi[offset++] = 0.5 * (distance[i] * distance[i] + cellSize * cellSize / 12.);
390 if (dimensions >= 2) {
391 csi[offset++] = distance[0] * distance[1];
393 if (dimensions >= 3) {
394 csi[offset++] = distance[0] * distance[2];
395 csi[offset++] = distance[1] * distance[2];
409 : m_degree(0), m_dimensions(0), m_nFields(0), m_nCoeffs(0)
422 const std::array<double, 3> &origin,
int nFields)
423 : m_degree(0), m_dimensions(0), m_nFields(0), m_nCoeffs(0)
425 initialize(degree, dimensions, origin, nFields);
438 int nWeights = m_nCoeffs * m_nFields;
439 std::copy(other.m_coeffs.get(), other.m_coeffs.get() + nWeights, m_coeffs.get());
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);
486 const std::array<double, 3> &origin,
487 int nFields,
bool release)
493 m_dimensions = dimensions;
496 int currentStorageSize = m_nCoeffs * m_nFields;
501 int storageSize = m_nCoeffs * m_nFields;
505 reallocate = (currentStorageSize != storageSize);
507 reallocate = (currentStorageSize < storageSize);
511 m_coeffs = std::unique_ptr<double[]>(
new double[storageSize]);
593 return m_coeffs.get();
603 return m_coeffs.get();
614 const double *coefficients = m_coeffs.get() + computeFieldCoefficientsOffset(0, field);
627 double *coefficients = m_coeffs.get() + computeFieldCoefficientsOffset(0, field);
643 const double *coefficients = m_coeffs.get() + computeFieldCoefficientsOffset(degree, field);
659 double *coefficients = m_coeffs.get() + computeFieldCoefficientsOffset(degree, field);
673std::size_t ReconstructionPolynomial::computeFieldCoefficientsOffset(uint8_t degree,
int field)
const
675 std::size_t offset = field * getFieldCoefficientsStride() + degree;
685std::size_t ReconstructionPolynomial::getFieldCoefficientsStride()
const
698 int field,
double *value)
const
710 double *values)
const
723 int nFields,
double *values)
const
737 int nFields,
int offset,
double *values)
const
751 int field,
double *values)
const
764 double *values)
const
778 int nFields,
double *values)
const
793 int nFields,
int offset,
double *values)
const
803 if (m_dimensions == 0) {
808 const std::size_t fieldValuesStride = 1;
809 double *fieldValue = values;
810 const double *fieldValueEnd = fieldValue + fieldValuesStride * nFields;
813 std::size_t fieldCoeffsStride = getFieldCoefficientsStride();
814 const double *fieldCoeffs = m_coeffs.get() + fieldCoeffsStride * offset;
817 if (ENABLE_FAST_PATH_OPTIMIZATIONS && degree == 0) {
819 *fieldValue = fieldCoeffs[0];
821 fieldCoeffs += fieldCoeffsStride;
822 fieldValue += fieldValuesStride;
823 }
while (fieldValue != fieldValueEnd);
833 *fieldValue = fieldCoeffs[0] * csi[0];
834 for (
int i = 1; i < m_nCoeffs; ++i) {
835 *fieldValue += fieldCoeffs[i] * csi[i];
838 fieldCoeffs += fieldCoeffsStride;
839 fieldValue += fieldValuesStride;
840 }
while (fieldValue != fieldValueEnd);
857 const double *limiters,
int field,
858 double *values)
const
876 const double *limiters,
877 double *values)
const
896 const double *limiters,
int nFields,
897 double *values)
const
917 const double *limiters,
int nFields,
int offset,
918 double *values)
const
938 const double *limiters,
int field,
939 double *values)
const
958 const double *limiters,
959 double *values)
const
979 const double *limiters,
int nFields,
980 double *values)
const
1001 const double *limiters,
int nFields,
int offset,
1002 double *values)
const
1012 if (m_dimensions == 0) {
1017 std::size_t fieldValuesStride = 1;
1018 double *fieldValue = values;
1019 double *fieldValueEnd = fieldValue + fieldValuesStride * nFields;
1022 std::size_t fieldCoeffsStride = getFieldCoefficientsStride();
1023 const double *fieldCoeffs = m_coeffs.get() + fieldCoeffsStride * offset;
1026 if (ENABLE_FAST_PATH_OPTIMIZATIONS && degree == 0) {
1028 *fieldValue = fieldCoeffs[0];
1030 fieldCoeffs += fieldCoeffsStride;
1031 fieldValue += fieldValuesStride;
1032 }
while (fieldValue != fieldValueEnd);
1038 std::size_t fieldLimitersStride = degree;
1039 const double *fieldLimiters = limiters;
1049 *fieldValue = fieldCoeffs[0] * csi[0];
1052 int coeffEnd = nDegreeCoeffs[0];
1053 for (
int n = 1; n <= degree; ++n) {
1054 int coeffBegin = coeffEnd;
1055 coeffEnd = coeffBegin + nDegreeCoeffs[n];
1057 double fieldDegreeValue = 0;
1058 for (
int i = coeffBegin; i < coeffEnd; ++i) {
1059 fieldDegreeValue += fieldCoeffs[i] * csi[i];
1062 *fieldValue += fieldDegreeValue * fieldLimiters[n - 1];
1065 fieldCoeffs += fieldCoeffsStride;
1066 fieldValue += fieldValuesStride;
1067 fieldLimiters += fieldLimitersStride;
1068 }
while (fieldValue != fieldValueEnd);
1080 const std::array<double, 3> &direction,
1081 int field,
double *derivative)
const
1094 const std::array<double, 3> &direction,
1095 double *derivatives)
const
1109 const std::array<double, 3> &direction,
1110 int nFields,
double *derivatives)
const
1125 const std::array<double, 3> &direction,
1126 int nFields,
int offset,
double *derivatives)
const
1141 const std::array<double, 3> &direction,
1142 int field,
double *derivative)
const
1156 const std::array<double, 3> &direction,
1157 double *derivatives)
const
1172 const std::array<double, 3> &direction,
1173 int nFields,
double *derivatives)
const
1189 const std::array<double, 3> &direction,
1190 int nFields,
int offset,
double *derivatives)
const
1200 if (m_dimensions == 0) {
1205 std::size_t fieldDerivativeStride = 1;
1206 double *fieldDerivative = derivatives;
1207 double *fieldDerivativeEnd = fieldDerivative + fieldDerivativeStride * nFields;
1210 if (ENABLE_FAST_PATH_OPTIMIZATIONS && degree == 0) {
1212 *fieldDerivative = 0.;
1214 fieldDerivative += fieldDerivativeStride;
1215 }
while (fieldDerivative != fieldDerivativeEnd);
1221 std::size_t fieldCoeffsStride = getFieldCoefficientsStride();
1222 const double *fieldCoeffs = m_coeffs.get() + fieldCoeffsStride * offset;
1225 if (ENABLE_FAST_PATH_OPTIMIZATIONS && degree == 1) {
1227 *fieldDerivative = fieldCoeffs[1] * direction[0];
1228 for (
int d = 1; d < m_dimensions; ++d) {
1229 *fieldDerivative += fieldCoeffs[d + 1] * direction[d];
1232 fieldCoeffs += fieldCoeffsStride;
1233 fieldDerivative += fieldDerivativeStride;
1234 }
while (fieldDerivative != fieldDerivativeEnd);
1244 *fieldDerivative = fieldCoeffs[0] * csi[0];
1245 for (
int i = 1; i < m_nCoeffs; ++i) {
1246 *fieldDerivative += fieldCoeffs[i] * csi[i];
1249 fieldCoeffs += fieldCoeffsStride;
1250 fieldDerivative += fieldDerivativeStride;
1251 }
while (fieldDerivative != fieldDerivativeEnd);
1269 const std::array<double, 3> &direction,
1270 const double *limiters,
int field,
1271 double *derivative)
const
1290 const std::array<double, 3> &direction,
1291 const double *limiters,
1292 double *derivatives)
const
1312 const std::array<double, 3> &direction,
1313 const double *limiters,
int nFields,
1314 double *derivatives)
const
1335 const std::array<double, 3> &direction,
1336 const double *limiters,
int nFields,
int offset,
1337 double *derivatives)
const
1358 const std::array<double, 3> &direction,
1359 const double *limiters,
int field,
1360 double *derivative)
const
1380 const std::array<double, 3> &direction,
1381 const double *limiters,
1382 double *derivatives)
const
1403 const std::array<double, 3> &direction,
1404 const double *limiters,
int nFields,
1405 double *derivatives)
const
1427 const std::array<double, 3> &direction,
1428 const double *limiters,
int nFields,
int offset,
1429 double *derivatives)
const
1439 if (m_dimensions == 0) {
1444 std::size_t fieldDerivativeStride = 1;
1445 double *fieldDerivative = derivatives;
1446 double *fieldDerivativeEnd = fieldDerivative + fieldDerivativeStride * nFields;
1449 if (ENABLE_FAST_PATH_OPTIMIZATIONS && degree == 0) {
1451 *fieldDerivative = 0.;
1453 fieldDerivative += fieldDerivativeStride;
1454 }
while (fieldDerivative != fieldDerivativeEnd);
1460 std::size_t fieldCoeffsStride = getFieldCoefficientsStride();
1461 const double *fieldCoeffs = m_coeffs.get() + fieldCoeffsStride * offset;
1464 std::size_t fieldLimitersStride = degree;
1465 const double *fieldLimiters = limiters;
1468 if (ENABLE_FAST_PATH_OPTIMIZATIONS && degree == 1) {
1470 *fieldDerivative = fieldCoeffs[1] * direction[0];
1471 for (
int d = 1; d < m_dimensions; ++d) {
1472 *fieldDerivative += fieldCoeffs[d + 1] * direction[d];
1475 *fieldDerivative *= fieldLimiters[0];
1477 fieldCoeffs += fieldCoeffsStride;
1478 fieldDerivative += fieldDerivativeStride;
1479 fieldLimiters += fieldLimitersStride;
1480 }
while (fieldDerivative != fieldDerivativeEnd);
1492 *fieldDerivative = fieldCoeffs[0] * dcsi[0];
1494 int coeffEnd = nDegreeCoeffs[0];
1495 for (
int n = 1; n <= degree; ++n) {
1496 int coeffBegin = coeffEnd;
1497 coeffEnd = coeffBegin + nDegreeCoeffs[n];
1499 double fieldDegreeDerivative = 0;
1500 for (
int i = coeffBegin; i < coeffEnd; ++i) {
1501 fieldDegreeDerivative += fieldCoeffs[i] * dcsi[i];
1504 *fieldDerivative += fieldLimiters[n - 1] * fieldDegreeDerivative;
1507 fieldCoeffs += fieldCoeffsStride;
1508 fieldDerivative += fieldDerivativeStride;
1509 fieldLimiters += fieldLimitersStride;
1510 }
while (fieldDerivative != fieldDerivativeEnd);
1521 int field, std::array<double, 3> *gradient)
const
1533 std::array<double, 3> *gradients)
const
1546 int nFields, std::array<double, 3> *gradients)
const
1560 int nFields,
int offset,
1561 std::array<double, 3> *gradients)
const
1575 int field, std::array<double, 3> *gradient)
const
1588 std::array<double, 3> *gradients)
const
1602 int nFields, std::array<double, 3> *gradients)
const
1617 int nFields,
int offset,
1618 std::array<double, 3> *gradients)
const
1628 if (m_dimensions == 0) {
1633 const std::size_t fieldGradientStride = 1;
1634 std::array<double, 3> *fieldGradient = gradients;
1635 const std::array<double, 3> *fieldGradientEnd = fieldGradient + fieldGradientStride * nFields;
1638 if (ENABLE_FAST_PATH_OPTIMIZATIONS && degree == 0) {
1640 *fieldGradient = {{0., 0., 0.}};
1642 fieldGradient += fieldGradientStride;
1643 }
while (fieldGradient != fieldGradientEnd);
1649 const std::size_t fieldCoeffsStride = getFieldCoefficientsStride();
1650 const double *fieldCoeffs = m_coeffs.get() + fieldCoeffsStride * offset;
1653 if (ENABLE_FAST_PATH_OPTIMIZATIONS && degree == 1) {
1656 for (
int d = 0; d < m_dimensions; ++d) {
1657 (*fieldGradient)[d] = fieldCoeffs[1 + d];
1662 (*fieldGradient)[d] = 0.;
1666 fieldCoeffs += fieldCoeffsStride;
1667 fieldGradient += fieldGradientStride;
1668 }
while (fieldGradient != fieldGradientEnd);
1674 BITPIT_CREATE_WORKSPACE(dcsi,
double,
static_cast<std::size_t
>(m_dimensions * m_nCoeffs), 3 * MAX_STACK_WORKSPACE_SIZE);
1677 if (m_dimensions == 3) {
1683 for (
int d = 0; d < m_dimensions; ++d) {
1684 const double *dcsi_dimension = dcsi + d * m_nCoeffs;
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];
1694 (*fieldGradient)[d] = 0.;
1698 fieldCoeffs += fieldCoeffsStride;
1699 fieldGradient += fieldGradientStride;
1700 }
while (fieldGradient != fieldGradientEnd);
1717 const double *limiters,
int field,
1718 std::array<double, 3> *gradient)
const
1736 const double *limiters,
1737 std::array<double, 3> *gradients)
const
1756 const double *limiters,
int nFields,
1757 std::array<double, 3> *gradients)
const
1777 const double *limiters,
int nFields,
int offset,
1778 std::array<double, 3> *gradients)
const
1798 const double *limiters,
int field,
1799 std::array<double, 3> *gradient)
const
1818 const double *limiters,
1819 std::array<double, 3> *gradients)
const
1839 const double *limiters,
int nFields,
1840 std::array<double, 3> *gradients)
const
1861 const double *limiters,
int nFields,
int offset,
1862 std::array<double, 3> *gradients)
const
1872 if (m_dimensions == 0) {
1877 const std::size_t fieldGradientStride = 1;
1878 std::array<double, 3> *fieldGradient = gradients;
1879 const std::array<double, 3> *fieldGradientEnd = fieldGradient + fieldGradientStride * nFields;
1882 if (ENABLE_FAST_PATH_OPTIMIZATIONS && degree == 0) {
1884 *fieldGradient = {{0., 0., 0.}};
1886 fieldGradient += fieldGradientStride;
1887 }
while (fieldGradient != fieldGradientEnd);
1893 const std::size_t fieldCoeffsStride = getFieldCoefficientsStride();
1894 const double *fieldCoeffs = m_coeffs.get() + fieldCoeffsStride * offset;
1897 const std::size_t fieldLimitersStride = degree;
1898 const double *fieldLimiters = limiters;
1901 if (ENABLE_FAST_PATH_OPTIMIZATIONS && degree == 1) {
1904 for (
int d = 0; d < m_dimensions; ++d) {
1905 (*fieldGradient)[d] = fieldLimiters[0] * fieldCoeffs[1 + d];
1910 (*fieldGradient)[d] = 0.;
1914 fieldCoeffs += fieldCoeffsStride;
1915 fieldGradient += fieldGradientStride;
1916 fieldLimiters += fieldLimitersStride;
1917 }
while (fieldGradient != fieldGradientEnd);
1925 BITPIT_CREATE_WORKSPACE(dcsi,
double,
static_cast<std::size_t
>(m_dimensions * nDegreeCoeffs[degree]), 3 * MAX_STACK_WORKSPACE_SIZE);
1928 if (m_dimensions == 3) {
1934 for (
int d = 0; d < m_dimensions; ++d) {
1935 const double *dcsi_dimension = dcsi + d * nDegreeCoeffs[degree];
1937 (*fieldGradient)[d] = fieldCoeffs[0] * dcsi_dimension[0];
1940 int coeffEnd = nDegreeCoeffs[0];
1941 for (
int n = 1; n <= degree; ++n) {
1942 double fieldsLimiter = fieldLimiters[n - 1];
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];
1950 (*fieldGradient)[d] += fieldsLimiter * fieldCoeffs[i] * dcsi_dimension[i];
1957 (*fieldGradient)[d] = 0.;
1961 fieldCoeffs += fieldCoeffsStride;
1962 fieldGradient += fieldGradientStride;
1963 fieldLimiters += fieldLimitersStride;
1964 }
while (fieldGradient != fieldGradientEnd);
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 <<
" : " ;
1983 for (
int n = 0; n < nDegreeCoeffs; ++n) {
1984 out << degreeCoeffs[n];
1985 if (n != nDegreeCoeffs - 1) {
2006const bool ReconstructionKernel::ENABLE_FAST_PATH_OPTIMIZATIONS =
true;
2011const int ReconstructionKernel::MAX_STACK_WORKSPACE_SIZE = 10;
2017 : m_nEquations(0), m_nCoeffs(0), m_degree(0), m_dimensions(0)
2030 : m_nEquations(0), m_nCoeffs(0), m_degree(0), m_dimensions(0)
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());
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);
2096 m_dimensions = dimensions;
2098 int currentStorageSize = m_nCoeffs * m_nEquations;
2101 m_nEquations = nEquations;
2103 int storageSize = m_nCoeffs * m_nEquations;
2107 reallocate = (currentStorageSize != storageSize);
2109 reallocate = (currentStorageSize < storageSize);
2113 m_weights = std::unique_ptr<double[]>(
new double[storageSize]);
2134 return m_dimensions;
2154 return m_nEquations;
2166 return m_weights.get();
2178 return m_weights.get();
2194 const double *values,
2201 polynomial->
initialize(degree, dimensions, origin, 1,
true);
2226 int nFields,
const double **values,
2233 polynomial->
initialize(degree, dimensions, origin, nFields,
true);
2256 const double *values,
2262 polynomial->
initialize(degree, dimensions, origin, 1,
true);
2289 int nFields,
const double **values,
2295 polynomial->
initialize(degree, dimensions, origin, nFields,
true);
2364 if (ENABLE_FAST_PATH_OPTIMIZATIONS && degree == 0) {
2370 polynomialCoeffs[0] = 0.;
2371 for (
int j = 0; j < nEquations; ++j) {
2372 polynomialCoeffs[0] += values[j] * constantPolynomialWeights[j];
2384 cblas_dgemv(CBLAS_ORDER::CblasColMajor, CBLAS_TRANSPOSE::CblasTrans,
2385 nEquations, nCoeffs, 1., polynomialWeights, nEquations,
2386 values, 1, 0, polynomialCoeffs, 1);
2417 if (ENABLE_FAST_PATH_OPTIMIZATIONS && degree == 0) {
2422 for (
int k = 0; k < nFields; ++k) {
2425 fieldCoeffs[0] = 0.;
2426 for (
int j = 0; j < nEquations; ++j) {
2427 fieldCoeffs[0] += values[j][k] * constantPolynomialWeights[j];
2442 for (
int k = 0; k < nFields; ++k) {
2443 for (
int j = 0; j < nEquations; ++j) {
2444 fieldValues[j] = values[j][k];
2449 cblas_dgemv(CBLAS_ORDER::CblasColMajor, CBLAS_TRANSPOSE::CblasTrans,
2450 nEquations, nCoeffs, 1., polynomialWeights, nEquations,
2451 fieldValues, 1, 0, fieldCoeffs, 1);
2470 const std::array<double, 3> &point,
double *valueWeights)
const
2494 const std::array<double, 3> &point,
double *valueWeights)
const
2520 const std::array<double, 3> &point,
double *valueWeights)
const
2547 const double *limiters,
double *valueWeights)
const
2575 const std::array<double, 3> &point,
const double *limiters,
2576 double *valueWeights)
const
2608 const std::array<double, 3> &point,
const double *limiters,
2609 double *valueWeights)
const
2615 if (ENABLE_FAST_PATH_OPTIMIZATIONS && degree == 0) {
2620 std::copy_n(constantPolynomialWeights, nEquations, valueWeights);
2639 cblas_dgemv(CBLAS_ORDER::CblasColMajor, CBLAS_TRANSPOSE::CblasNoTrans,
2640 nEquations, nCoeffs, 1., polynomialWeights, polynomialWeightsRowCount,
2641 csi, 1, 0, valueWeights, 1);
2661 const std::array<double, 3> &point,
const std::array<double, 3> &direction,
2662 double *derivativeWeights)
const
2688 const std::array<double, 3> &point,
const std::array<double, 3> &direction,
2689 double *derivativeWeights)
const
2716 const std::array<double, 3> &point,
const std::array<double, 3> &direction,
2717 double *derivativeWeights)
const
2746 const std::array<double, 3> &point,
const std::array<double, 3> &direction,
2747 const double *limiters,
2748 double *derivativeWeights)
const
2781 const std::array<double, 3> &point,
const std::array<double, 3> &direction,
2782 const double *limiters,
2783 double *derivativeWeights)
const
2817 const std::array<double, 3> &point,
const std::array<double, 3> &direction,
2818 const double *limiters,
2819 double *derivativeWeights)
const
2825 if (ENABLE_FAST_PATH_OPTIMIZATIONS && degree == 0) {
2826 std::fill_n(derivativeWeights, nEquations, 0.);
2845 cblas_dgemv(CBLAS_ORDER::CblasColMajor, CBLAS_TRANSPOSE::CblasNoTrans,
2846 nEquations, nCoeffs, 1., polynomialWeights, polynomialWeightsRowCount,
2847 dcsi, 1, 0, derivativeWeights, 1);
2865 const std::array<double, 3> &point,
2866 std::array<double, 3> *gradientWeights)
const
2890 const std::array<double, 3> &point,
2891 std::array<double, 3> *gradientWeights)
const
2916 const std::array<double, 3> &point,
2917 std::array<double, 3> *gradientWeights)
const
2943 const std::array<double, 3> &point,
const double *limiters,
2944 std::array<double, 3> *gradientWeights)
const
2974 const std::array<double, 3> &point,
const double *limiters,
2975 std::array<double, 3> *gradientWeights)
const
3006 const std::array<double, 3> &point,
const double *limiters,
3007 std::array<double, 3> *gradientWeights)
const
3013 if (ENABLE_FAST_PATH_OPTIMIZATIONS && degree == 0) {
3014 std::fill_n(gradientWeights, nEquations, std::array<double, 3>{0., 0., 0.});
3027 for (
int d = 0; d < dimensions; ++d) {
3029 std::array<double, 3> direction = {{0., 0., 0.}};
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);
3047 for (
int j = 0; j < nEquations; ++j) {
3049 gradientWeights[j][d] = 0.;
3076 int coeffEnd = nDegreeCoeffs[0];
3077 for (
int n = 1; n <= degree; ++n) {
3078 double limiterScaleFactor = limiters[n - 1];
3080 int coeffBegin = coeffEnd;
3081 coeffEnd = coeffBegin + nDegreeCoeffs[n];
3082 for (
int k = coeffBegin; k < coeffEnd; ++k) {
3083 coeffs[k] *= limiterScaleFactor;
3098 for (
int i = 0; i < nCoeffs; ++i) {
3100 for (
int j = 0; j < m_nEquations; ++j) {
3102 if (std::abs(weigth) < tolerance) {
3106 out <<
"(" << j <<
"," << weigth <<
") ";
3124const double ReconstructionAssembler::SVD_ZERO_THRESHOLD = 1e-14;
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);
3186 m_dimensions = dimensions;
3201 m_constraintsOrder.clear();
3202 m_leastSquaresOrder.clear();
3203 m_leastSquaresScaleFactors.clear();
3212 m_SVDWorkspace.resize(1);
3216 m_constraintsOrder.shrink_to_fit();
3217 m_leastSquaresOrder.shrink_to_fit();
3218 m_leastSquaresScaleFactors.shrink_to_fit();
3220 m_A.shrink_to_fit();
3221 m_C.shrink_to_fit();
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();
3249 return m_dimensions;
3269 return m_constraintsOrder.size();
3279 return m_leastSquaresOrder.size();
3291 int nEquations = nConstraints + nLeastSquares;
3305 const std::array<double, 3> &origin,
3306 const std::array<double, 3> &point,
3309 double *equationCoeffs = _addEquation(type, scaleFactor);
3323 const std::array<double, 3> &origin,
3324 const std::array<double, 3> &point,
3325 const std::array<double, 3> &direction,
3328 double *equationCoeffs = _addEquation(type, scaleFactor);
3343 const std::array<double, 3> &origin,
3344 const std::array<double, 3> *vertexCoords,
3347 double *equationCoeffs = _addEquation(type, scaleFactor);
3357double * ReconstructionAssembler::_addEquation(ReconstructionType type,
double scaleFactor)
3363 case TYPE_CONSTRAINT:
3364 m_constraintsOrder.emplace_back(nEquations);
3367 case TYPE_LEAST_SQUARE:
3368 m_leastSquaresOrder.emplace_back(nEquations);
3369 m_leastSquaresScaleFactors.emplace_back(scaleFactor);
3377 double *equationCoeffsStorage =
nullptr;
3380 case TYPE_CONSTRAINT:
3381 m_C.resize(m_C.size() + nCoeffs);
3382 equationCoeffsStorage = m_C.data() + m_C.size() - nCoeffs;
3385 case TYPE_LEAST_SQUARE:
3386 m_A.resize(m_A.size() + nCoeffs);
3387 equationCoeffsStorage = m_A.data() + m_A.size() - nCoeffs;
3392 return equationCoeffsStorage;
3411 kernel->
initialize(degree, dimensions, nEquations,
true);
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);
3441 m_w.resize(nLeastSquares);
3442 for (
int k = 0; k < nLeastSquares; ++k) {
3443 m_w[k] = m_leastSquaresScaleFactors[k] / maxLeastSquareScaleFactor;
3466 int nUnknowns = nCoeffs + nConstraints;
3468 m_S.resize(nUnknowns * nUnknowns);
3469 for (
int i = 0; i < nCoeffs; ++i) {
3470 for (
int j = i; j < nCoeffs; ++j) {
3473 for (
int k = 0; k < nLeastSquares; ++k) {
3476 ATA_ij += m_A[A_ki_idx] * m_A[A_kj_idx] * m_w[k];
3486 for (
int j = nCoeffs; j < nUnknowns; ++j) {
3489 m_S[l] = m_C[C_idx];
3500 computePseudoInverse(nUnknowns, nUnknowns, SVD_ZERO_THRESHOLD, m_S.data());
3516 for (
int j = 0; j < nEquations; ++j) {
3518 if (j < nLeastSquares) {
3519 equation = m_leastSquaresOrder[j];
3521 equation = m_constraintsOrder[j - nLeastSquares];
3524 for (
int i = 0; i < nCoeffs; ++i) {
3526 for (
int k = 0; k < nUnknowns; ++k) {
3528 if (k < nCoeffs && j < nLeastSquares) {
3530 value += m_S[l] * m_A[A_jk_idx] * m_w[j];
3531 }
else if ((k - nCoeffs) == (j - nLeastSquares)) {
3537 weights[weightLineraIndex] = value;
3554void ReconstructionAssembler::computePseudoInverse(
int m,
int n,
double zeroThreshold,
double *A)
const
3559 int k = std::min(m,n);
3568 int workspaceSize = -1;
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);
3575 workspaceSize =
static_cast<int>(m_SVDWorkspace[0]);
3576 if (workspaceSize > (
int) m_SVDWorkspace.size()) {
3577 m_SVDWorkspace.resize(workspaceSize);
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);
3584 log::cout() <<
"SVD failed in ReconstructionAssembler::computePseudoInverse()" <<std::endl;
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);
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);
The Cell class defines the cells.
ElementType getType() const
double evalSize(const std::array< double, 3 > *coordinates) const
std::array< double, 3 > evalCentroid(const std::array< double, 3 > *coordinates) const
The ReconstructionAssembler class allows to define reconstruction polynomial.
uint8_t getDegree() const
uint16_t getCoefficientCount() const
void swap(ReconstructionAssembler &other) noexcept
void clear(bool release=true)
uint8_t getDimensions() const
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)
int countConstraints() const
int countEquations() const
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.)
ReconstructionAssembler()
void updateKernel(ReconstructionKernel *kernel) const
int countLeastSquares() 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
int getEquationCount() const
uint8_t getDegree() 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)
uint8_t getDimensions() const
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
uint16_t getCoefficientCount() 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
static const uint8_t MAX_DIMENSIONS
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
uint16_t getCoefficientCount() const
void display(std::ostream &out) const
uint8_t getDegree() const
void computeValuesLimited(const std::array< double, 3 > &point, const double *limiters, double *values) const
void clear(bool release=true)
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)
int getFieldCount() const
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
static const uint8_t MAX_DEGREE
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
uint8_t getDimensions() 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)
ReconstructionPolynomial()
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)