Loading...
Searching...
No Matches
postProcess.H
Go to the documentation of this file.
1Info << "Computing errors" << endl;
2gField.resize(0);
3ITHACAstream::read_fields(gField, heatFluxFieldName,
4 folder);
5volScalarField gTrueField = list2Field(gTrue, innerField);
6label Nsolutions = gField.size();
7Eigen::MatrixXd heatFluxL2norm;
8heatFluxL2norm.resize(Nsolutions, 1);
9Eigen::MatrixXd heatFluxLinfNorm = heatFluxL2norm;
10Eigen::VectorXd relativeErrorNorm_L2(Nsolutions);
11Eigen::VectorXd relativeErrorNorm_Linf(Nsolutions);
12
13forAll(gField, solutionI)
14{
15 volScalarField gDiffField = (gField[solutionI] - gTrueField).ref();
17 std::to_string(solutionI + 1), folder,
18 "gDiffField");
20 std::to_string(solutionI + 1), folder,
21 "gTrueField");
22 relativeErrorNorm_L2(solutionI) = ITHACAutilities::L2normOnPatch(mesh,
23 gDiffField,
24 "hotSide") / ITHACAutilities::L2normOnPatch(mesh, gTrueField, "hotSide");
25 relativeErrorNorm_Linf(solutionI) = ITHACAutilities::LinfNormOnPatch(mesh,
26 gDiffField,
27 "hotSide") / ITHACAutilities::LinfNormOnPatch(mesh, gTrueField, "hotSide");
28 scalar EPS = 1e-6;
29 volScalarField relativeErrorField(gDiffField);
30
31 for (label i = 0; i < relativeErrorField.internalField().size(); i++)
32 {
33 if (std::abs(gTrueField.ref()[i]) < EPS)
34 {
35 relativeErrorField.ref()[i] = (std::abs(gDiffField.ref()[i])) / EPS;
36 }
37 else
38 {
39 relativeErrorField.ref()[i] = (std::abs(gDiffField.ref()[i])) /
40 gTrueField.ref()[i];
41 }
42 }
43
44 ITHACAstream::exportSolution(relativeErrorField,
45 std::to_string(solutionI + 1), folder,
46 "relativeErrorField");
47 heatFluxL2norm(solutionI) = ITHACAutilities::L2normOnPatch(mesh,
48 relativeErrorField,
49 "hotSide");
50 heatFluxLinfNorm(solutionI) = ITHACAutilities::LinfNormOnPatch(mesh,
51 relativeErrorField,
52 "hotSide");
53}
54ITHACAstream::exportMatrix(heatFluxL2norm, "relError_L2norm", "eigen",
55 folder);
56ITHACAstream::exportMatrix(heatFluxLinfNorm, "relError_LinfNorm", "eigen",
57 folder);
58ITHACAstream::exportMatrix(relativeErrorNorm_L2, "relativeErrorNorm_L2",
59 "eigen",
60 folder);
61ITHACAstream::exportMatrix(relativeErrorNorm_Linf, "relativeErrorNorm_Linf",
62 "eigen",
63 folder);
forAll(example_CG.gList, solutionI)
Definition CGtest.H:21
Foam::fvMesh & mesh
Definition createMesh.H:47
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 ...
void read_fields(PtrList< GeometricField< Type, PatchField, GeoMesh > > &Lfield, word Name, fileName casename, int first_snap, int n_snap)
Function to read a list of fields from the name of the field and casename.
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.
label i
Definition pEqn.H:46