Loading...
Searching...
No Matches
ITHACAfieldsOperations.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 ITHACAfieldsOperations_H
31#define ITHACAfieldsOperations_H
32
35
36#include "fvCFD.H"
37#pragma GCC diagnostic push
38#pragma GCC diagnostic ignored "-Wold-style-cast"
39#include <Eigen/Eigen>
40#pragma GCC diagnostic pop
41using namespace std::placeholders;
42#include "ITHACAstream.H"
43#include "ITHACAparameters.H"
44#include "ITHACAgeometry.H"
45#include "ITHACAsystem.H"
46#include "ITHACAassign.H"
47
48
49namespace ITHACAutilities
50{
51
52//------------------------------------------------------------------------------
59template<typename T>
60void multField(T& f1, double alpha);
61
62//------------------------------------------------------------------------------
69template<typename T>
70void multField(T &f1, const Eigen::VectorXd alphaVec);
71
72//------------------------------------------------------------------------------
80template<typename T>
81void multField(PtrList<T> &f1, const Eigen::VectorXd alphaVec);
82
83
84//------------------------------------------------------------------------------
94template<typename T>
95void addFields(T& f1, const T& f2c, double alpha = 1.0);
96
97
98//------------------------------------------------------------------------------
105template<typename T>
106void subtractFields(T& f1, const T& f2);
107
108
109//------------------------------------------------------------------------------
115volTensorField tensorFieldProduct(const volScalarField& coef, const volTensorField& S);
116
117//------------------------------------------------------------------------------
123volTensorField tensorFieldProduct(const volTensorField& coef, const volTensorField& S);
124
125
126//------------------------------------------------------------------------------
130int dimensionField(const volTensorField& v);
131
132//------------------------------------------------------------------------------
136int dimensionField(const volVectorField& v);
137
138//------------------------------------------------------------------------------
142int dimensionField(const volScalarField& v);
143
144
145//-----------------------------------------------------------------------------
158template<class TypeField>
159PtrList<TypeField> averageSubtract(PtrList<TypeField>
160 fields, Eigen::MatrixXd ind, PtrList<TypeField>& ave);
161
162
163//-----------------------------------------------------------------------------
172template<class TypeField>
173TypeField computeAverage(PtrList<TypeField>& fields);
174
175
176//--------------------------------------------------------------------------
183template<typename Type>
184void normalizeFields(
185 PtrList<GeometricField<Type, fvPatchField, volMesh >>& fields);
186
187
188//------------------------------------------------------------------------------
199template<typename Type>
200Eigen::MatrixXd getValues(GeometricField<Type, fvPatchField,
201 volMesh>& field, labelList& indices, labelList* xyz = NULL);
202
203
204//------------------------------------------------------------------------------
215template<typename Type>
216Eigen::MatrixXd getValues(PtrList<GeometricField<Type, fvPatchField,
217 volMesh >>& field, labelList& indices, labelList* xyz = NULL);
218
219};
220
221#endif
Header file of the ITHACAassign file.
Header file of the geometry namespace.
Header file of the ITHACAstream class, it contains the implementation of several methods for input ou...
Header file of the ITHACAsystem file.
Namespace to implement some useful assign operation of OF fields.
int dimensionField(const volTensorField &v)
Return the dimension of a volTensorField.
volTensorField tensorFieldProduct(const volScalarField &coef, const volTensorField &S)
Tensor field product between a volScalarField and a volTensorField.
void multField(T &f1, double alpha)
Multiplication between a field of type vol[Scalar|Vector|Tensor]Field and a double.
TypeField computeAverage(PtrList< TypeField > &fields)
Calculates the average of a list of fields.
void normalizeFields(PtrList< GeometricField< Type, fvPatchField, volMesh > > &fields)
Normalize list of Geometric fields.
void addFields(T &f1, const T &f2c, double alpha)
Perform the following operation f1 + f2 * alpha with f1 and f2 being of type vol[Scalar|Vector|Tensor...
PtrList< TypeField > averageSubtract(PtrList< TypeField > fields, Eigen::MatrixXd ind, PtrList< TypeField > &ave)
A function to compute time-averaged fields for a set of different parameter samples and also the fiel...
void subtractFields(T &f1, const T &f2)
Perform the substraction (f1 - f2) between two fields of type vol[Scalar|Vector|Tensor]Field and alph...