LagrangianAccumulationScheme.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 
27 #include "LagrangianSubMesh.H"
28 
29 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
30 
31 template<class Type>
34 (
35  const LagrangianMesh& mesh,
36  Istream& is
37 )
38 {
39  if (is.eof())
40  {
42  << "Accumulation scheme not specified" << endl << endl
43  << "Valid accumulation schemes are :" << endl
44  << IstreamConstructorTablePtr_->sortedToc()
45  << exit(FatalIOError);
46  }
47 
48  const word schemeName(is);
49 
50  typename IstreamConstructorTable::iterator cstrIter =
51  IstreamConstructorTablePtr_->find(schemeName);
52 
53  if (cstrIter == IstreamConstructorTablePtr_->end())
54  {
56  << "Unknown accumulation scheme " << schemeName << nl << nl
57  << "Valid accumulation schemes are :" << endl
58  << IstreamConstructorTablePtr_->sortedToc()
59  << exit(FatalIOError);
60  }
61 
62  return cstrIter()(mesh, is);
63 }
64 
65 
66 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
67 
68 template<class Type>
70 {}
71 
72 
73 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
74 
75 template<class Type>
76 template<class CellMesh, template<class> class PrimitiveField>
79 (
81 )
82 {
83  typedef typename DimensionedField<Type, CellMesh>::Mesh resultMeshType;
84 
87  (
88  lPsi.name(),
89  refCast<const resultMeshType>(lPsi.mesh().mesh()),
91  );
92 
94  (
95  lPsi.mesh().subAll().sub(lPsi),
96  tResult.ref().primitiveFieldRef()
97  );
98 
99  return tResult;
100 }
101 
102 
103 template<class Type>
104 template<class CellMesh, template<class> class PrimitiveField>
106 (
109 )
110 {
111  vPsi.dimensions() += lPsi.dimensions();
112 
113  accumulate(toSubField(lPsi)(), vPsi.primitiveFieldRef());
114 }
115 
116 
117 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
GeoMesh::Mesh Mesh
Type of mesh on which this DimensionedField is instantiated.
const dimensionSet & dimensions() const
Return dimensions.
PrimitiveField< Type > & primitiveFieldRef()
Return a reference to the internal field.
const Mesh & mesh() const
Return mesh.
const word & name() const
Return name.
Definition: IOobject.H:307
bool eof() const
Return true if end of input seen.
Definition: IOstream.H:339
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
Class containing Lagrangian geometry and topology.
Basic cell Lagrangian accumulation scheme.
static tmp< accumulationScheme< Type > > New(const LagrangianMesh &mesh, Istream &is)
Return a pointer to a new accumulationScheme.
Generic dimensioned Type class.
Traits class for primitives.
Definition: pTraits.H:53
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 FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
tmp< DimensionedField< Type, CellMesh > > accumulate(const DimensionedField< Type, LagrangianMesh, PrimitiveField > &lPsi, const word &name)
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
IOerror FatalIOError
static const char nl
Definition: Ostream.H:267
tmp< DimensionedField< Type, GeoMesh, SubField > > toSubField(const DimensionedField< Type, GeoMesh, Field > &)
Return a temporary sub-field from a reference to a field.