36#include "ITHACAparameters.H"
38#pragma GCC diagnostic push
39#pragma GCC diagnostic ignored "-Wold-style-cast"
41#pragma GCC diagnostic pop
57 const word& hilbertSpacePOD,
58 const double& weightBC=0.,
const word& patchBC=
"inlet",
59 const double& weightH1 = 1.0);
71 const word& hilbertSpacePOD,
72 const double& weightBC=0.,
const word& patchBC=
"inlet",
73 const double& weightH1 = 1.0);
86 const word& hilbertSpacePOD,
87 const double& weightBC=0.,
const word& patchBC=
"inlet",
88 const double& weightH1 = 1.0);
103 const word& hilbertSpacePOD,
104 const double& weightBC,
106 const double& weightH1 = 1.0);
114double dot_product_L2(
const volVectorField& v,
const volVectorField& w);
123double dot_product_L2(
const volScalarField& v,
const volScalarField& w);
131double dot_product_L2(
const volTensorField& v,
const volTensorField& w);
141double dot_product_H1(
const T& v,
const T& w,
const double& weightH1 = 1.0);
154double dot_product_patch(
const Eigen::VectorXd& f1BC_i,
const Eigen::VectorXd& f2BC_i,
const scalarField& AreaFace,
const int& d);
179double dot_product_L2wBC(
const T& v,
const T& w,
const double& weightBC,
const word& patchBC=
"inlet");
192double norm_L2wBC(
const T& v,
const double& weightBC,
const word& patchBC=
"inlet");
207 const word& hilbertSpacePOD,
208 const double& weightBC=0.,
const word& patchBC=
"inlet",
209 const double& weightH1 = 1.0);
231double LinfNorm(GeometricField<T, fvPatchField, volMesh>& field);
243double H1Seminorm(GeometricField<T, fvPatchField, volMesh>& field);
254template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
255double frobNorm(GeometricField<Type, PatchField, GeoMesh>& field);
266double L2normOnPatch(fvMesh& mesh, volScalarField& field, word patch);
279 List<scalar>& field2, word patch);
290double LinfNormOnPatch(fvMesh& mesh, volScalarField& field, word patch);
301double integralOnPatch(fvMesh& mesh, volScalarField& field, word patch);
Header file of the ITHACAstream class, it contains the implementation of several methods for input ou...
Namespace to implement some useful assign operation of OF fields.
double norm_L2wBC(const T &v, const double &weightBC, const word &patchBC)
Compute the norm using the L2 dot product of v at the boundaries.
double LinfNormOnPatch(fvMesh &mesh, volScalarField &field, word patch)
Evaluate the Linf norm of a field on a boundary patch.
double dot_product_boundary(const T &v, const T &w, const word &patchBC)
Perform the L2 dot product of v and w at a given boundary.
Eigen::MatrixXd dot_product_POD(PtrList< T > &v, PtrList< T > &w, const word &hilbertSpacePOD, const double &weightBC, const word &patchBC, const double &weightH1)
Perform the dot product of two PtrList of type T.
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 frobNorm(GeometricField< Type, PatchField, GeoMesh > &field)
Evaluate the Frobenius norm of a field.
double dot_product_patch(const Eigen::VectorXd &f1BC_i, const Eigen::VectorXd &f2BC_i, const scalarField &AreaFace, const int &d)
Perform the L2 dot product of v and w at a given boundary patch.
double L2Norm(const T v)
Compute the L2 norm of v.
double dot_product_L2(const volVectorField &v, const volVectorField &w)
Perform the L2 dot product of v and w.
double norm_POD(const T &v, const word &hilbertSpacePOD, const double &weightBC, const word &patchBC, const double &weightH1)
Compute the norm of v.
double dot_product_H1(const T &v, const T &w, const double &weightH1)
Perform the H1 dot product of v and w.
double L2productOnPatch(fvMesh &mesh, List< scalar > &field1, List< scalar > &field2, word patch)
Evaluate the L2 inner product between two scalarLists.
double dot_product_L2wBC(const T &v, const T &w, const double &weightBC, const word &patchBC)
Perform the L2 dot product of v and w at the boundaries.