LagrangianModel.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 "LagrangianModel.H"
27 
28 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
34 }
35 
36 
37 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
38 
39 template<class Type>
41 (
42  const LagrangianSubScalarField& deltaT,
43  const LagrangianSubSubField<Type>& field,
45 ) const
46 {}
47 
48 
49 template<class Type>
51 (
52  const LagrangianSubScalarField& deltaT,
54  const LagrangianSubSubField<Type>& field,
56 ) const
57 {}
58 
59 
60 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
61 
63 (
64  const word& name,
65  const LagrangianMesh& mesh
66 )
67 :
68  name_(name),
69  mesh_(mesh)
70 {}
71 
72 
73 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
74 
76 (
77  const word& name,
78  const LagrangianMesh& mesh,
79  const dictionary& modelDict
80 )
81 {
82  const word type(modelDict.lookup("type"));
83 
84  Info<< "Selecting " << typeName
85  << " with name " << name
86  << " of type " << type << endl;
87 
88  if
89  (
90  !dictionaryConstructorTablePtr_
91  || dictionaryConstructorTablePtr_->find(type)
92  == dictionaryConstructorTablePtr_->end()
93  )
94  {
95  if (!libs.open(modelDict, "libs", dictionaryConstructorTablePtr_))
96  {
97  libs.open("lib" + type.remove(':') + ".so", false);
98  }
99 
100  if (!dictionaryConstructorTablePtr_)
101  {
103  << "Unknown " << typeName << " type "
104  << type << nl << nl
105  << "Table of " << typeName << "s is empty"
106  << exit(FatalError);
107  }
108  }
109 
110  dictionaryConstructorTable::iterator cstrIter =
111  dictionaryConstructorTablePtr_->find(type);
112 
113  if (cstrIter == dictionaryConstructorTablePtr_->end())
114  {
115  FatalIOErrorInFunction(modelDict)
116  << "Unknown " << typeName << " " << type << nl << nl
117  << "Valid " << typeName << "s are:" << nl
118  << dictionaryConstructorTablePtr_->sortedToc()
119  << exit(FatalIOError);
120  }
121 
122  return
123  cstrIter()
124  (
125  name,
126  mesh,
127  modelDict.optionalSubDict(type + "Coeffs"),
128  stateDict(name, mesh)
129  );
130 }
131 
132 
133 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
134 
136 {}
137 
138 
139 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
140 
142 {
143  return wordList::null();
144 }
145 
146 
147 bool Foam::LagrangianModel::addsSupToField(const word& fieldName) const
148 {
149  return findIndex(addSupFields(), fieldName) != -1;
150 }
151 
152 
154 {}
155 
156 
158 {}
159 
160 
162 (
163  const LagrangianMesh& mesh,
164  DynamicList<elementModification>& elementModifications
165 ) const
166 {}
167 
168 
170 (
172  const LagrangianSubMesh& modifiedMesh
173 ) const
174 {
175  return mesh_.subNone();
176 }
177 
178 
180 (
181  const LagrangianSubScalarField& deltaT,
182  const bool final
183 )
184 {}
185 
186 
188 (
189  const LagrangianSubScalarField& deltaT,
191 ) const
192 {}
193 
194 
196 
197 
199 
200 
201 void Foam::LagrangianModel::topoChange(const polyTopoChangeMap&)
202 {}
203 
204 
206 {}
207 
208 
210 {}
211 
212 
214 {
215  return true;
216 }
217 
218 
219 bool Foam::LagrangianModel::write(const bool write) const
220 {
221  return write;
222 }
223 
224 
225 // ************************************************************************* //
#define IMPLEMENT_LAGRANGIAN_MODEL_ADD_FIELD_SUP(Type, modelType)
#define IMPLEMENT_LAGRANGIAN_MODEL_ADD_M_FIELD_SUP(Type, modelType)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:78
This class stores the coefficients of a Lagrangian equation, and facilitates solving that equation an...
Definition: LagrangianEqn.H:56
Class containing Lagrangian geometry and topology.
Base class for Lagrangian models.
virtual ~LagrangianModel()
Destructor.
virtual wordList addSupFields() const
Return the list of fields for which the LagrangianModel adds.
virtual void addSup(const LagrangianSubScalarField &deltaT, LagrangianSubScalarField &S) const
Add a fractional source term.
virtual void postConstruct()
Do post construction steps which require access to other models.
void addSupType(const LagrangianSubScalarField &deltaT, const LagrangianSubSubField< Type > &field, LagrangianEqn< Type > &eqn) const
Add a source term to an equation.
virtual void preModify(const LagrangianMesh &mesh, DynamicList< elementModification > &elementModifications) const
Identify elements in the Lagrangian mesh which are to be.
virtual void correct()
Correct the LagrangianModel.
virtual bool addsSupToField(const word &) const
Return true if the LagrangianModel adds a source term to the.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
virtual bool read(const dictionary &modelDict)
Read dictionary.
virtual bool write(const bool write) const
Write data.
LagrangianModel(const word &name, const LagrangianMesh &mesh)
Construct from components.
virtual LagrangianSubMesh modify(LagrangianMesh &mesh, const LagrangianSubMesh &modifiedMesh) const
Instantaneously modify and/or create and remove elements in the.
virtual void calculate(const LagrangianSubScalarField &deltaT, const bool final)
Solve equations and/or update continually changing properties.
static autoPtr< LagrangianModel > New(const word &name, const LagrangianMesh &mesh, const dictionary &modelDict)
Selector.
Mesh that relates to a sub-section of a Lagrangian mesh. This is used to construct fields that relate...
static const List< word > & null()
Return a null List.
Definition: ListI.H:118
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:740
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:927
bool open(const fileName &libName, const bool verbose=true)
Open the named library, optionally with warnings if problems occur.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
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
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const VolField< Type > &)
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dlLibraryTable libs
Table of loaded dynamic libraries.
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
messageStream Info
defineTypeNameAndDebug(combustionModel, 0)
IOerror FatalIOError
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
error FatalError
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
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
dictionary dict