39#include "freestreamFvPatchField.H"
42#pragma GCC diagnostic push
43#pragma GCC diagnostic ignored "-Wold-style-cast"
45#pragma GCC diagnostic pop
48#include "polyMeshTools.H"
50#include "mixedFvPatchFields.H"
51#include "fvMeshSubset.H"
52using namespace std::placeholders;
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);
87double errorL2Rel(GeometricField<T, fvPatchField, volMesh>& field1,
88 GeometricField<T, fvPatchField, volMesh>& field2, List<label>* labels = NULL);
102double errorLinfRel(GeometricField<T, fvPatchField, volMesh>& field1,
103 GeometricField<T, fvPatchField, volMesh>& field2, List<label>* labels = NULL);
118double errorL2Abs(GeometricField<T, fvPatchField, volMesh>& field1,
119 GeometricField<T, fvPatchField, volMesh>& field2, volScalarField& Volumes);
133double errorL2Abs(GeometricField<T, fvPatchField, volMesh>& field1,
134 GeometricField<T, fvPatchField, volMesh>& field2, List<label>* labels = NULL);
148Eigen::MatrixXd
errorL2Rel(PtrList<GeometricField<T, fvPatchField, volMesh>>&
150 PtrList<GeometricField<T, fvPatchField, volMesh>>& fields2,
151 List<label>* labels = NULL);
164template<
class T,
template<
class>
class PatchField,
class GeoMesh>
165Eigen::MatrixXd
errorFrobRel(PtrList<GeometricField<T, PatchField, GeoMesh>>&
167 PtrList<GeometricField<T, PatchField, GeoMesh>>& fields2,
168 List<label>* labels = NULL);
184 PtrList<GeometricField<T, fvPatchField, volMesh>>& fields1,
185 PtrList<GeometricField<T, fvPatchField, volMesh>>& fields2,
186 PtrList<volScalarField>& Volumes);
200Eigen::MatrixXd
errorL2Abs(PtrList<GeometricField<T, fvPatchField, volMesh>>&
202 PtrList<GeometricField<T, fvPatchField, volMesh>>& fields2,
203 List<label>* labels = NULL);
215double L2Norm(GeometricField<T, fvPatchField, volMesh>& field);
227double LinfNorm(GeometricField<T, fvPatchField, volMesh>& field);
239double H1Seminorm(GeometricField<T, fvPatchField, volMesh>& field);
250template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
251double frobNorm(GeometricField<Type, PatchField, GeoMesh>& field);
275 List<scalar>& field2, word patch);
Header file of the Foam2Eigen class.
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)
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)
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.
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)