noneStateLagrangianLabelFieldSource.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 
28 
29 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
30 
33 {}
34 
35 
36 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
37 
40 (
41  const LagrangianInjection& injection,
42  const LagrangianSubMesh& subMesh
43 ) const
44 {
45  return
47  (
48  internalField().name() + ":" + injection.name(),
49  subMesh,
51  (
53  internalDimensions(),
54  static_cast<label>(LagrangianState::none)
55  )
56  );
57 }
58 
59 
60 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61 
62 namespace Foam
63 {
65  (
66  LagrangianLabelFieldSource,
68  );
69 }
70 
71 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
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,.
Base class for Lagrangian injections. Minimal wrapper over LagrangianSource. Implements some utility ...
const word & name() const
The source name.
Mesh that relates to a sub-section of a Lagrangian mesh. This is used to construct fields that relate...
Generic dimensioned Type class.
This source condition provides a label corresponding to the LagrangianState::none enumeration.
virtual tmp< LagrangianSubLabelField > value(const LagrangianInjection &, const LagrangianSubMesh &) const
Return the value for an instantaneous injection.
A class for managing temporary objects.
Definition: tmp.H:55
Namespace for OpenFOAM.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
makeNullConstructableLagrangianTypeFieldSource(LagrangianScalarFieldSource, massLagrangianScalarFieldSource)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.