Loading...
Searching...
No Matches
ITHACAnorm.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 ITHACAnorm_H
31#define ITHACAnorm_H
32
35
36#include "ITHACAparameters.H"
37#include "fvCFD.H"
38#pragma GCC diagnostic push
39#pragma GCC diagnostic ignored "-Wold-style-cast"
40#include <Eigen/Eigen>
41#pragma GCC diagnostic pop
42#include "ITHACAstream.H"
43
44namespace ITHACAutilities
45{
46
47//------------------------------------------------------------------------------
49
56double dot_product_POD(const volVectorField& v, const volVectorField& w,
57 const word& hilbertSpacePOD,
58 const double& weightBC=0., const word& patchBC="inlet",
59 const double& weightH1 = 1.0);
60
61//------------------------------------------------------------------------------
63
70double dot_product_POD(const volScalarField& v, const volScalarField& w,
71 const word& hilbertSpacePOD,
72 const double& weightBC=0., const word& patchBC="inlet",
73 const double& weightH1 = 1.0);
74
75
76//------------------------------------------------------------------------------
78
85double dot_product_POD(const volTensorField& v, const volTensorField& w,
86 const word& hilbertSpacePOD,
87 const double& weightBC=0., const word& patchBC="inlet",
88 const double& weightH1 = 1.0);
89
90
91//------------------------------------------------------------------------------
93
101template<typename T>
102Eigen::MatrixXd dot_product_POD(PtrList<T>& v, PtrList<T>& w,
103 const word& hilbertSpacePOD,
104 const double& weightBC,
105 const word& patchBC,
106 const double& weightH1 = 1.0);
107
108//------------------------------------------------------------------------------
110
114double dot_product_L2(const volVectorField& v, const volVectorField& w);
115
116
117//------------------------------------------------------------------------------
119
123double dot_product_L2(const volScalarField& v, const volScalarField& w);
124
125//------------------------------------------------------------------------------
127
131double dot_product_L2(const volTensorField& v, const volTensorField& w);
132
133//------------------------------------------------------------------------------
135
140template<typename T>
141double dot_product_H1(const T& v, const T& w, const double& weightH1 = 1.0);
142
143//------------------------------------------------------------------------------
145
152
154double dot_product_patch(const Eigen::VectorXd& f1BC_i,const Eigen::VectorXd& f2BC_i, const scalarField& AreaFace, const int& d);
155template<typename T>
156
157//------------------------------------------------------------------------------
159
164
166double dot_product_boundary(const T& v,const T& w, const word& patchBC);
167
168//------------------------------------------------------------------------------
170
176
178template<typename T>
179double dot_product_L2wBC(const T& v, const T& w, const double& weightBC, const word& patchBC="inlet");
180
181//------------------------------------------------------------------------------
183
189
191template<typename T>
192double norm_L2wBC(const T& v, const double& weightBC, const word& patchBC="inlet");
193
194//------------------------------------------------------------------------------
196
203
205template<typename T>
206double norm_POD(const T& v,
207 const word& hilbertSpacePOD,
208 const double& weightBC=0., const word& patchBC="inlet",
209 const double& weightH1 = 1.0);
210
211//------------------------------------------------------------------------------
213
216
218template<typename T>
219double L2Norm(const T v);
220
221//------------------------------------------------------------------------------
230template<class T>
231double LinfNorm(GeometricField<T, fvPatchField, volMesh>& field);
232
233//------------------------------------------------------------------------------
242template<class T>
243double H1Seminorm(GeometricField<T, fvPatchField, volMesh>& field);
244
245//------------------------------------------------------------------------------
254template<class Type, template<class> class PatchField, class GeoMesh>
255double frobNorm(GeometricField<Type, PatchField, GeoMesh>& field);
256
257//--------------------------------------------------------------------------
266double L2normOnPatch(fvMesh& mesh, volScalarField& field, word patch);
267
268//--------------------------------------------------------------------------
278double L2productOnPatch(fvMesh& mesh, List<scalar>& field1,
279 List<scalar>& field2, word patch);
280
281//--------------------------------------------------------------------------
290double LinfNormOnPatch(fvMesh& mesh, volScalarField& field, word patch);
291
292//--------------------------------------------------------------------------
301double integralOnPatch(fvMesh& mesh, volScalarField& field, word patch);
302
303//--------------------------------------------------------------------------
312double integralOnPatch(fvMesh& mesh, List<scalar> field, word patch);
313
314}
315
316#endif
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.
Definition ITHACAnorm.C:281
double LinfNormOnPatch(fvMesh &mesh, volScalarField &field, word patch)
Evaluate the Linf norm of a field on a boundary patch.
Definition ITHACAnorm.C:406
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.
Definition ITHACAnorm.C:220
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.
Definition ITHACAnorm.C:41
double integralOnPatch(fvMesh &mesh, volScalarField &field, word patch)
Evaluate the integral on a patch.
Definition ITHACAnorm.C:432
double L2normOnPatch(fvMesh &mesh, volScalarField &field, word patch)
Evaluate the L2 norm of a field on a boundary patch.
Definition ITHACAnorm.C:366
double frobNorm(GeometricField< Type, PatchField, GeoMesh > &field)
Evaluate the Frobenius norm of a field.
Definition ITHACAnorm.C:355
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.
Definition ITHACAnorm.C:200
double L2Norm(const T v)
Compute the L2 norm of v.
Definition ITHACAnorm.C:300
double dot_product_L2(const volVectorField &v, const volVectorField &w)
Perform the L2 dot product of v and w.
Definition ITHACAnorm.C:182
double norm_POD(const T &v, const word &hilbertSpacePOD, const double &weightBC, const word &patchBC, const double &weightH1)
Compute the norm of v.
Definition ITHACAnorm.C:291
double dot_product_H1(const T &v, const T &w, const double &weightH1)
Perform the H1 dot product of v and w.
Definition ITHACAnorm.C:174
double L2productOnPatch(fvMesh &mesh, List< scalar > &field1, List< scalar > &field2, word patch)
Evaluate the L2 inner product between two scalarLists.
Definition ITHACAnorm.C:386
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.
Definition ITHACAnorm.C:264