LagrangiancAccumulate.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 "LagrangiancAccumulate.H"
28 #include "LagrangianMesh.H"
29 #include "LagrangianSubFields.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 template<class CellMesh, class Type, template<class> class PrimitiveField>
36 (
37  const DimensionedField<Type, LagrangianMesh, PrimitiveField>& lPsi
38 )
39 {
40  return accumulate<CellMesh>(lPsi, lPsi.name());
41 }
42 
43 
44 template<class CellMesh, class Type, template<class> class PrimitiveField>
47 (
48  const DimensionedField<Type, LagrangianMesh, PrimitiveField>& lPsi,
49  const word& name
50 )
51 {
52  const LagrangianMesh& mesh = lPsi.mesh();
53 
54  return
56  (
57  mesh,
58  mesh.schemes().accumulation(name)
59  ).ref().template accumulate<CellMesh>(lPsi);
60 }
61 
62 
63 template<class CellMesh, class Type, template<class> class PrimitiveField>
65 (
68 )
69 {
70  accumulate<CellMesh>(lPsi, vPsi, lPsi.name());
71 }
72 
73 
74 template<class CellMesh, class Type, template<class> class PrimitiveField>
76 (
79  const word& name
80 )
81 {
82  const LagrangianMesh& mesh = lPsi.mesh().mesh();
83 
85  (
86  mesh,
87  mesh.schemes().accumulation(name)
88  ).ref().template accumulate<CellMesh>(lPsi, vPsi);
89 }
90 
91 
92 // ************************************************************************* //
Functions for accumulating Lagrangian values into cell fields.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const Mesh & mesh() const
Return mesh.
const word & name() const
Return name.
Definition: IOobject.H:307
Class containing Lagrangian geometry and topology.
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1799
const polyMesh & mesh() const
Return reference to polyMesh.
Definition: fvMesh.H:443
A class for managing temporary objects.
Definition: tmp.H:55
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)
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)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.