Loading...
Searching...
No Matches
ITHACAerror.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ██╗████████╗██╗ ██╗ █████╗ ██████╗ █████╗ ███████╗██╗ ██╗
3 ██║╚══██╔══╝██║ ██║██╔══██╗██╔════╝██╔══██╗ ██╔════╝██║ ██║
4 ██║ ██║ ███████║███████║██║ ███████║█████╗█████╗ ██║ ██║
5 ██║ ██║ ██╔══██║██╔══██║██║ ██╔══██║╚════╝██╔══╝ ╚██╗ ██╔╝
6 ██║ ██║ ██║ ██║██║ ██║╚██████╗██║ ██║ ██║ ╚████╔╝
7 ╚═╝ ╚═╝ ╚═╝ ╚═╝╚═╝ ╚═╝ ╚═════╝╚═╝ ╚═╝ ╚═╝ ╚═══╝
8
9 * In real Time Highly Advanced Computational Applications for Finite Volumes
10 * Copyright (C) 2017 by the ITHACA-FV authors
11-------------------------------------------------------------------------------
12
13License
14 This file is part of ITHACA-FV
15
16 ITHACA-FV is free software: you can redistribute it and/or modify
17 it under the terms of the GNU Lesser General Public License as published by
18 the Free Software Foundation, either version 3 of the License, or
19 (at your option) any later version.
20
21 ITHACA-FV is distributed in the hope that it will be useful,
22 but WITHOUT ANY WARRANTY; without even the implied warranty of
23 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
24 GNU Lesser General Public License for more details.
25
26 You should have received a copy of the GNU Lesser General Public License
27 along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
28\*---------------------------------------------------------------------------*/
29
30#ifndef ITHACAerror_H
31#define ITHACAerror_H
32
35
36
37#include "fvCFD.H"
38#include "IOmanip.H"
39#include "freestreamFvPatchField.H"
40#include <sys/stat.h>
41#include <unistd.h>
42#pragma GCC diagnostic push
43#pragma GCC diagnostic ignored "-Wold-style-cast"
44#include <Eigen/Eigen>
45#pragma GCC diagnostic pop
46#include <functional>
47#include "./colormod.H"
48#include "polyMeshTools.H"
49#include <chrono>
50#include "mixedFvPatchFields.H"
51#include "fvMeshSubset.H"
52using namespace std::placeholders;
53#include "Foam2Eigen.H"
54
55namespace ITHACAutilities
56{
57
58//------------------------------------------------------------------------------
71template<class Type, template<class> class PatchField, class GeoMesh>
72double errorFrobRel(GeometricField<Type, PatchField, GeoMesh>& field1,
73 GeometricField<Type, PatchField, GeoMesh>& field2, List<label>* labels = NULL);
74
75//--------------------------------------------------------------------------
86template<typename T>
87double errorL2Rel(GeometricField<T, fvPatchField, volMesh>& field1,
88 GeometricField<T, fvPatchField, volMesh>& field2, List<label>* labels = NULL);
89
90//--------------------------------------------------------------------------
101template<typename T>
102double errorLinfRel(GeometricField<T, fvPatchField, volMesh>& field1,
103 GeometricField<T, fvPatchField, volMesh>& field2, List<label>* labels = NULL);
104
105//--------------------------------------------------------------------------
117template<typename T>
118double errorL2Abs(GeometricField<T, fvPatchField, volMesh>& field1,
119 GeometricField<T, fvPatchField, volMesh>& field2, volScalarField& Volumes);
120
121//--------------------------------------------------------------------------
132template<typename T>
133double errorL2Abs(GeometricField<T, fvPatchField, volMesh>& field1,
134 GeometricField<T, fvPatchField, volMesh>& field2, List<label>* labels = NULL);
135
136//--------------------------------------------------------------------------
147template<typename T>
148Eigen::MatrixXd errorL2Rel(PtrList<GeometricField<T, fvPatchField, volMesh>>&
149 fields1,
150 PtrList<GeometricField<T, fvPatchField, volMesh>>& fields2,
151 List<label>* labels = NULL);
152
153//--------------------------------------------------------------------------
164template<class T, template<class> class PatchField, class GeoMesh>
165Eigen::MatrixXd errorFrobRel(PtrList<GeometricField<T, PatchField, GeoMesh>>&
166 fields1,
167 PtrList<GeometricField<T, PatchField, GeoMesh>>& fields2,
168 List<label>* labels = NULL);
169//--------------------------------------------------------------------------
170
182template<typename T>
183Eigen::MatrixXd errorL2Abs(
184 PtrList<GeometricField<T, fvPatchField, volMesh>>& fields1,
185 PtrList<GeometricField<T, fvPatchField, volMesh>>& fields2,
186 PtrList<volScalarField>& Volumes);
187
188//--------------------------------------------------------------------------
199template<typename T>
200Eigen::MatrixXd errorL2Abs(PtrList<GeometricField<T, fvPatchField, volMesh>>&
201 fields1,
202 PtrList<GeometricField<T, fvPatchField, volMesh>>& fields2,
203 List<label>* labels = NULL);
204
205//------------------------------------------------------------------------------
214template<class T>
215double L2Norm(GeometricField<T, fvPatchField, volMesh>& field);
216
217//------------------------------------------------------------------------------
226template<class T>
227double LinfNorm(GeometricField<T, fvPatchField, volMesh>& field);
228
229//------------------------------------------------------------------------------
238template<class T>
239double H1Seminorm(GeometricField<T, fvPatchField, volMesh>& field);
240
241//------------------------------------------------------------------------------
250template<class Type, template<class> class PatchField, class GeoMesh>
251double frobNorm(GeometricField<Type, PatchField, GeoMesh>& field);
252
253//--------------------------------------------------------------------------
262double L2normOnPatch(fvMesh& mesh, volScalarField& field, word patch);
263
264//--------------------------------------------------------------------------
274double L2productOnPatch(fvMesh& mesh, List<scalar>& field1,
275 List<scalar>& field2, word patch);
276
277//--------------------------------------------------------------------------
286double LinfNormOnPatch(fvMesh& mesh, volScalarField& field, word patch);
287
288//--------------------------------------------------------------------------
297double integralOnPatch(fvMesh& mesh, volScalarField& field, word patch);
298
299//--------------------------------------------------------------------------
308double integralOnPatch(fvMesh& mesh, List<scalar> field, word patch);
309
310}
311
312
313#endif
Header file of the Foam2Eigen class.
Foam::fvMesh & mesh
Definition createMesh.H:47
Simple header and source file of the Color::Modifier class to change color to the output stream.
Namespace to implement some useful assign operation of OF fields.
double LinfNorm(GeometricField< scalar, fvPatchField, volMesh > &field)
Definition ITHACAerror.C:57
double errorLinfRel(GeometricField< T, fvPatchField, volMesh > &field1, GeometricField< T, fvPatchField, volMesh > &field2, List< label > *labels)
Computes the relative error between two geometric Fields in Linf norm.
double LinfNormOnPatch(fvMesh &mesh, volScalarField &field, word patch)
Evaluate the Linf norm of a field on a boundary patch.
double H1Seminorm(GeometricField< scalar, fvPatchField, volMesh > &field)
double L2Norm(GeometricField< scalar, fvPatchField, volMesh > &field)
Definition ITHACAerror.C:41
double integralOnPatch(fvMesh &mesh, volScalarField &field, word patch)
Evaluate the integral on a patch.
double L2normOnPatch(fvMesh &mesh, volScalarField &field, word patch)
Evaluate the L2 norm of a field on a boundary patch.
double errorL2Rel(GeometricField< T, fvPatchField, volMesh > &field1, GeometricField< T, fvPatchField, volMesh > &field2, List< label > *labels)
Computes the relative error between two geometric Fields in L2 norm.
double frobNorm(GeometricField< Type, PatchField, GeoMesh > &field)
Evaluate the Frobenius norm of a field.
double errorFrobRel(GeometricField< Type, PatchField, GeoMesh > &field1, GeometricField< Type, PatchField, GeoMesh > &field2, List< label > *labels)
Computes the relative error between two Fields in the Frobenius norm.
Definition ITHACAerror.C:78
double L2productOnPatch(fvMesh &mesh, List< scalar > &field1, List< scalar > &field2, word patch)
Evaluate the L2 inner product between two scalarLists.
double errorL2Abs(GeometricField< vector, fvPatchField, volMesh > &field1, GeometricField< vector, fvPatchField, volMesh > &field2, volScalarField &Volumes)