Skip to content

Commit ce0f171

Browse files
authored
Add error computation for openfoam solver (#449)
1 parent 5321b03 commit ce0f171

File tree

2 files changed

+113
-3
lines changed

2 files changed

+113
-3
lines changed

partitioned-heat-conduction/solver-openfoam/heatTransfer.C

Lines changed: 67 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -37,6 +37,46 @@
3737
#include "simpleControl.H"
3838

3939
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40+
// Some helper functions in order to investigate this case
41+
namespace Functions {
42+
double get_rhs(const double x, const double y, const double alpha, const double beta)
43+
{
44+
return beta - 2 - 2 * alpha;
45+
}
46+
double get_solution(const double x, const double y,
47+
const double alpha, const double beta, const double time)
48+
{
49+
return 1 + std::pow(x, 2) + (alpha * std::pow(y, 2)) + (beta * time);
50+
}
51+
52+
void compute_and_print_errors(const Foam::fvMesh &mesh, const double alpha,
53+
const double beta, const double time)
54+
{
55+
double error = 0;
56+
const Foam::volScalarField *T_(&mesh.lookupObject<volScalarField>("T"));
57+
58+
double max_error = std::numeric_limits<double>::min();
59+
60+
// Get the locations of the volume centered mesh vertices
61+
const vectorField &CellCenters = mesh.C();
62+
unsigned int numDataLocations = CellCenters.size();
63+
for (int i = 0; i < CellCenters.size(); i++) {
64+
const double coord_x = CellCenters[i].x();
65+
const double coord_y = CellCenters[i].y();
66+
const double exact_solution = Functions::get_solution(coord_x, coord_y, alpha, beta, time);
67+
const auto result = (exact_solution - T_->internalField()[i]) *
68+
(exact_solution - T_->internalField()[i]);
69+
error += result;
70+
max_error = std::max(result, max_error);
71+
}
72+
73+
Info << "\nError metrics at t = " << time << "s:\n";
74+
Info << "l2 error (sqrt(sum(err^2)/n)):\t\t" << std::sqrt(error / numDataLocations) << endl;
75+
Info << "Maximum absolute error (max(err)):\t" << std::sqrt(max_error) << endl;
76+
Info << "Global absolute error (sum(err^2)):\t" << error << "\n";
77+
Info << endl;
78+
}
79+
} // namespace Functions
4080

4181
int main(int argc, char *argv[])
4282
{
@@ -61,8 +101,8 @@ int main(int argc, char *argv[])
61101

62102
const double alpha = 3;
63103
const double beta = 1.2;
64-
const double rhs = beta - 2 - 2 * alpha;
65104

105+
// Initialize the RHS with zero
66106
volScalarField f(
67107
IOobject(
68108
"RHS",
@@ -74,7 +114,32 @@ int main(int argc, char *argv[])
74114
dimensionedScalar(
75115
"Tdim",
76116
dimensionSet(0, 0, -1, 1, 0, 0, 0),
77-
Foam::scalar(rhs)));
117+
Foam::scalar(0)));
118+
119+
// Now assign the respective values to the RHS
120+
//(strictly speaking not required here, only for space-dependent RHS functions)
121+
// Get the locations of the volume centered mesh vertices
122+
const vectorField &CellCenters = mesh.C();
123+
for (int i = 0; i < CellCenters.size(); i++) {
124+
const double coord_x = CellCenters[i].x();
125+
const double coord_y = CellCenters[i].y();
126+
f.ref()[i] = Functions::get_rhs(coord_x, coord_y, alpha, beta);
127+
}
128+
129+
for (int j = 0; j < mesh.boundaryMesh().size() - 1; ++j) {
130+
// Get the face centers of the current patch
131+
const vectorField faceCenters =
132+
mesh.boundaryMesh()[j].faceCentres();
133+
134+
// Assign the (x,y,z) locations to the vertices
135+
for (int i = 0; i < faceCenters.size(); i++) {
136+
const double coord_x = faceCenters[i].x();
137+
const double coord_y = faceCenters[i].y();
138+
f.boundaryFieldRef()[j][i] = Functions::get_rhs(coord_x, coord_y, alpha, beta);
139+
}
140+
}
141+
142+
Functions::compute_and_print_errors(mesh, alpha, beta, runTime.value());
78143

79144
while (simple.loop()) {
80145
Info << "Time = " << runTime.timeName() << nl << endl;

partitioned-heat-conduction/solver-openfoam/write.H

Lines changed: 46 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -1,4 +1,7 @@
1-
if (runTime.writeTime()) {
1+
// We need to disable the runTime.writeTime if condition for implicit coupling
2+
// as the condition only returns true in the very first coupling iteration
3+
// if (runTime.writeTime())
4+
{
25
volVectorField gradT(fvc::grad(T));
36

47
volScalarField gradTx(
@@ -36,6 +39,48 @@ if (runTime.writeTime()) {
3639
IOobject::NO_READ,
3740
IOobject::AUTO_WRITE),
3841
DT * gradT);
42+
volScalarField error_total(
43+
IOobject(
44+
"error",
45+
runTime.timeName(),
46+
mesh,
47+
IOobject::NO_READ,
48+
IOobject::AUTO_WRITE),
49+
DT * 0);
50+
51+
// Print error metrics
52+
Functions::compute_and_print_errors(mesh, alpha, beta, runTime.value());
53+
54+
// compute and store the error also in a paraView field
55+
{
56+
const Foam::volScalarField *T_(&mesh.lookupObject<volScalarField>("T"));
57+
58+
// Get the locations of the volume centered mesh vertices
59+
const vectorField &CellCenters = mesh.C();
60+
61+
for (int i = 0; i < CellCenters.size(); i++) {
62+
const double coord_x = CellCenters[i].x();
63+
const double coord_y = CellCenters[i].y();
64+
const double exact_solution = Functions::get_solution(coord_x, coord_y, alpha, beta, runTime.value());
65+
66+
error_total.ref()[i] = std::abs((exact_solution - T_->internalField()[i]));
67+
}
68+
69+
T.correctBoundaryConditions();
70+
for (int j = 0; j < mesh.boundaryMesh().size() - 1; ++j) {
71+
// Get the face centers of the current patch
72+
const vectorField faceCenters =
73+
mesh.boundaryMesh()[j].faceCentres();
74+
75+
// Assign the (x,y,z) locations to the vertices
76+
for (int i = 0; i < faceCenters.size(); i++) {
77+
const double coord_x = faceCenters[i].x();
78+
const double coord_y = faceCenters[i].y();
79+
const double exact_solution = Functions::get_solution(coord_x, coord_y, alpha, beta, runTime.value());
80+
error_total.boundaryFieldRef()[j][i] = std::abs(exact_solution - T_->boundaryField()[j][i]);
81+
}
82+
}
83+
}
3984

4085
runTime.write();
4186
}

0 commit comments

Comments
 (0)