2Info <<
"*********************************************************" << endl;
3Info <<
"Performing test for the parameterized BC inverse solver" << endl;
5word outputFolder =
"./ITHACAoutput/parameterizedBCtest_RBFparameter/";
6volScalarField gTrueField = example_paramBC.list2Field(example_paramBC.gTrue);
16Eigen::VectorXd residualNorms;
17residualNorms.resize(rbfWidthTest_size);
18Eigen::VectorXd heatFluxL2norm(rbfWidthTest_size);
19Eigen::VectorXd heatFluxLinfNorm = heatFluxL2norm;
21Eigen::MatrixXd singVal;
23for (
int i = 0; i < rbfWidthTest_size; i++)
25 Info <<
"*********************************************************" << endl;
26 Info <<
"RBF parameter " << rbfWidth(i) << endl;
27 Info <<
"Test " << i + 1 << endl;
29 example_paramBC.set_gParametrized(
"rbf", rbfWidth(i));
30 example_paramBC.parameterizedBCoffline(1);
31 example_paramBC.parameterizedBC(
"fullPivLU");
32 volScalarField gParametrizedField = example_paramBC.list2Field(
35 std::to_string(i + 1),
38 volScalarField& T(example_paramBC._T());
40 std::to_string(i + 1),
43 residualNorms(i) = Foam::sqrt(
44 example_paramBC.residual.squaredNorm());
45 Eigen::MatrixXd A = example_paramBC.Theta.transpose() * example_paramBC.Theta;
46 Eigen::JacobiSVD<Eigen::MatrixXd> svd(A,
47 Eigen::ComputeThinU | Eigen::ComputeThinV);
48 Eigen::MatrixXd singularValues = svd.singularValues();
49 singVal.conservativeResize(singularValues.rows(), singVal.cols() + 1);
50 singVal.col(i) = singularValues;
51 double conditionNumber;
53 if (singularValues.minCoeff() > 0)
55 conditionNumber = singularValues.maxCoeff() / singularValues.minCoeff();
59 conditionNumber = DBL_MAX;
62 Info <<
"Condition number = " << conditionNumber << endl;
64 volScalarField gDiffField = (gParametrizedField - gTrueField).ref();
66 volScalarField relativeErrorField(gTrueField);
68 for (label i = 0; i < relativeErrorField.internalField().size(); i++)
70 if (std::abs(gTrueField.ref()[i]) < EPS)
72 relativeErrorField.ref()[i] = (std::abs(gDiffField.ref()[i])) / EPS;
76 relativeErrorField.ref()[i] = (std::abs(gDiffField.ref()[i])) /
82 std::to_string(i + 1), outputFolder,
83 "relativeErrorField");
100example_paramBC.postProcess(outputFolder,
"gParametrized");
101Info <<
"*********************************************************" << endl;
102Info <<
"*********************************************************" << endl;
Eigen::VectorXd ExpSpaced(double first, double last, int n)
Returns exponentially spaced vector.
T condNumber(Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > &A)
Conditioning number of a dense matrix.
void exportSolution(GeometricField< Type, PatchField, GeoMesh > &s, fileName subfolder, fileName folder, word fieldName)
Export a field to file in a certain folder and subfolder.
void exportMatrix(Eigen::Matrix< T, -1, dim > &matrix, word Name, word type, word folder)
Export the reduced matrices in numpy (type=python), matlab (type=matlab) and txt (type=eigen) format ...
double LinfNormOnPatch(fvMesh &mesh, volScalarField &field, word patch)
Evaluate the Linf norm of a field on a boundary patch.
double L2normOnPatch(fvMesh &mesh, volScalarField &field, word patch)
Evaluate the L2 norm of a field on a boundary patch.