Function1LagrangianFieldSource.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 "LagrangianMesh.H"
28 #include "LagrangianModel.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class Type>
33 template<class Derived>
35 (
36  const Derived& field
37 )
38 :
39  field_(static_cast<const LagrangianFieldSource<Type>&>(field))
40 {}
41 
42 
43 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
44 
45 template<class Type>
46 template<class OtherType>
49 (
50  const LagrangianModel& model,
51  const LagrangianSubMesh& subMesh,
52  const dimensionSet& dims,
53  const Function1<OtherType>& function
54 )
55 {
56  const objectRegistry& db = subMesh.mesh();
57 
58  const scalar t1 = db.time().value();
59  const scalar t0 = t1 - db.time().deltaTValue();
60 
61  const word name = model.name() + ":" + function.name();
62 
63  if
64  (
65  db.template foundObject<LagrangianScalarInternalDynamicField>
66  (
68  )
69  )
70  {
71  const LagrangianSubScalarSubField fractionSf
72  (
73  subMesh.sub
74  (
75  db.template lookupObject<LagrangianScalarInternalDynamicField>
76  (
78  )
79  )
80  );
81 
82  return
84  (
85  name,
86  subMesh,
87  dims,
88  function.value(t0 + fractionSf*(t1 - t0))
89  );
90  }
91  else
92  {
93  return
95  (
96  name,
97  subMesh,
98  dimensioned<OtherType>(dims, function.value(t0))
99  );
100  }
101 }
102 
103 
104 template<class Type>
107 (
108  const LagrangianModel& model,
109  const LagrangianSubMesh& subMesh,
110  const Function1<Type>& function
111 ) const
112 {
113  return value(model, subMesh, field_.internalDimensions(), function);
114 }
115 
116 
117 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
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,.
static tmp< LagrangianSubField< OtherType > > value(const LagrangianModel &model, const LagrangianSubMesh &subMesh, const dimensionSet &dims, const Function1< OtherType > &function)
Return the source value.
Function1LagrangianFieldSource(const Derived &)
Construct with reference to the derived field source.
Run-time selectable general function of one variable.
Definition: Function1.H:125
Base class for Lagrangian source conditions.
static const word fractionName
Name of the tracked fraction field.
Base class for Lagrangian models.
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...
SubList< Type > sub(const List< Type > &list) const
Return a sub-list corresponding to this sub-mesh.
const LagrangianMesh & mesh() const
Return the mesh.
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:34
Dimension set for the base types.
Definition: dimensionSet.H:125
Generic dimensioned Type class.
const Type & value() const
Return const reference to value.
Registry of regIOobjects.
const Time & time() const
Return time.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.