LagrangianAverage.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "LagrangianAverage.H"
27 #include "LagrangianFields.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
31 template<class Type>
33 (
34  const word& name,
35  const LagrangianMesh& mesh,
37 )
38 :
39  name_(name),
40  mesh_(mesh),
41  dimensions_(dimensions)
42 {}
43 
44 
45 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
46 
47 template<class Type>
48 template<template<class> class WeightPF, template<class> class PsiPF>
51 (
52  const word& type,
53  const word& name,
56 )
57 {
58  typename dictionaryConstructorTable::iterator cstrIter =
59  dictionaryConstructorTablePtr_->find(type);
60 
61  if (cstrIter == dictionaryConstructorTablePtr_->end())
62  {
64  << "Unknown Lagrangian average type " << type
65  << " for field " << psi.name() << nl << nl
66  << "Valid average types : " << endl
67  << dictionaryConstructorTablePtr_->sortedToc()
68  << exit(FatalError);
69  }
70 
72  cstrIter()
73  (
74  name,
75  psi.mesh(),
76  psi.dimensions(),
77  NullObjectRef<Field<scalar>>()
78  );
79 
80  result->add
81  (
82  weight.mesh().subAll().sub(weight),
83  psi.mesh().subAll().sub(psi),
84  false
85  );
86 
87  return result;
88 }
89 
90 
91 template<class Type>
92 template<class CellMesh, template<class> class PsiPF>
95 (
96  const word& type,
97  const word& name,
98  const DimensionedField<scalar, CellMesh>& cellWeightSum,
100 )
101 {
102  typename dictionaryConstructorTable::iterator cstrIter =
103  dictionaryConstructorTablePtr_->find(type);
104 
105  if (cstrIter == dictionaryConstructorTablePtr_->end())
106  {
108  << "Unknown Lagrangian average type " << type
109  << " for field " << weightPsi.name() << nl << nl
110  << "Valid average types : " << endl
111  << dictionaryConstructorTablePtr_->sortedToc()
112  << exit(FatalError);
113  }
114 
116  cstrIter()
117  (
118  name,
119  weightPsi.mesh(),
120  weightPsi.dimensions()/cellWeightSum.dimensions(),
121  cellWeightSum
122  );
123 
124  result->add
125  (
127  weightPsi.mesh().subAll().sub(weightPsi),
128  false
129  );
130 
131  return result;
132 }
133 
134 
135 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
136 
137 template<class Type>
139 {}
140 
141 
142 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
143 
144 template<class Type>
146 (
147  const LagrangianSubSubField<Type>& weightPsi
148 )
149 {
150  remove(NullObjectRef<LagrangianSubSubField<scalar>>(), weightPsi);
151 }
152 
153 
154 template<class Type>
156 (
157  const LagrangianSubSubField<Type>& weightPsi,
158  const bool cache
159 )
160 {
161  add(NullObjectRef<LagrangianSubSubField<scalar>>(), weightPsi, cache);
162 }
163 
164 
165 template<class Type>
167 (
168  const LagrangianSubSubField<Type>& weightPsi,
169  const bool cache
170 )
171 {
173 }
174 
175 
176 template<class Type>
179 (
180  const LagrangianSubMesh& subMesh
181 ) const
182 {
184  LagrangianSubField<Type>::New(name_, subMesh, dimensions_);
185 
186  interpolate(tresult.ref());
187 
188  return tresult;
189 }
190 
191 
192 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
const Mesh & mesh() const
Return mesh.
static tmp< DimensionedField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Mesh &mesh, const dimensionSet &, const PrimitiveField< Type > &)
Return a temporary field constructed from name, mesh,.
const word & name() const
Return name.
Definition: IOobject.H:307
Basic Lagrangian averaging process in which values are averaged in the cells and then interpolated as...
virtual void correct(const LagrangianSubSubField< scalar > &weight, const LagrangianSubSubField< Type > &psi, const bool cache)=0
Correct weighted values in the average.
LagrangianAverage(const word &name, const LagrangianMesh &mesh, const dimensionSet &dimensions)
Construct with a name, for a mesh and with given dimensions.
virtual void add(const LagrangianSubSubField< scalar > &weight, const LagrangianSubSubField< Type > &psi, const bool cache)=0
Add weighted values to the average.
virtual void remove(const LagrangianSubSubField< scalar > &weight, const LagrangianSubSubField< Type > &psi)=0
Remove weighted values from the average.
virtual ~LagrangianAverage()
Destructor.
static autoPtr< LagrangianAverage< Type > > New(const word &type, const word &name, const DimensionedField< scalar, LagrangianMesh, WeightPF > &weight, const DimensionedField< Type, LagrangianMesh, PsiPF > &psi)
Select to average a field.
Class containing Lagrangian geometry and topology.
Mesh that relates to a sub-section of a Lagrangian mesh. This is used to construct fields that relate...
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Dimension set for the base types.
Definition: dimensionSet.H:125
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
const volScalarField & psi
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp)
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
Definition: dimensionSets.C:96
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
const T & NullObjectRef()
Return const reference to the nullObject of type T.
Definition: nullObjectI.H:27
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
error FatalError
void add(LagrangianPatchField< typename typeOfSum< Type1, Type2 >::type > &f, const LagrangianPatchField< Type1 > &f1, const LagrangianPatchField< Type2 > &f2)
static const char nl
Definition: Ostream.H:267
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488