fvOptionList.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2019 OpenCFD Ltd.
10 Copyright (C) 2020 PCOpt/NTUA
11 Copyright (C) 2020 FOSS GP
12-------------------------------------------------------------------------------
13License
14 This file is part of OpenFOAM.
15
16 OpenFOAM is free software: you can redistribute it and/or modify it
17 under the terms of the GNU General Public License as published by
18 the Free Software Foundation, either version 3 of the License, or
19 (at your option) any later version.
20
21 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
22 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
23 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
24 for more details.
25
26 You should have received a copy of the GNU General Public License
27 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
28
29Class
30 Foam::fv::optionList
31
32Description
33 List of finite volume options
34
35SourceFile
36 optionList.C
37
38\*---------------------------------------------------------------------------*/
39
40#ifndef Foam_fvOptionList_H
41#define Foam_fvOptionList_H
42
43#include "fvOption.H"
44#include "PtrList.H"
45#include "GeometricField.H"
46#include "geometricOneField.H"
47#include "fvPatchField.H"
48
49// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50
51namespace Foam
52{
53
54// Forward Declarations
55
56namespace fv
57{
58 class optionList;
59}
60
61Ostream& operator<<(Ostream& os, const fv::optionList& options);
62
63namespace fv
64{
65
66/*---------------------------------------------------------------------------*\
67 Class optionList Declaration
68\*---------------------------------------------------------------------------*/
70class optionList
71:
72 public PtrList<fv::option>
73{
74protected:
75
76 // Protected Data
77
78 //- Reference to the mesh database
79 const fvMesh& mesh_;
80
81 //- Time index to check that all defined sources have been applied
82 label checkTimeIndex_;
83
84
85 // Protected Member Functions
86
87 //- Return "options" sub-dictionary (if present) or return dict
88 static const dictionary& optionsDict(const dictionary& dict);
89
90 //- Read options dictionary
91 bool readOptions(const dictionary& dict);
92
93 //- Check that all sources have been applied
94 void checkApplied() const;
95
96 //- Return source for equation with specified name and dimensions
97 template<class Type>
99 (
101 const word& fieldName,
102 const dimensionSet& ds
103 );
104
105 //- No copy construct
106 optionList(const optionList&) = delete;
107
108 //- No copy assignment
109 void operator=(const optionList&) = delete;
110
111
112public:
113
114 //- Runtime type information
115 TypeName("optionList");
116
117
118 // Constructors
119
120 //- Default construct from mesh
121 explicit optionList(const fvMesh& mesh);
122
123 //- Construct from mesh and dictionary
124 optionList(const fvMesh& mesh, const dictionary& dict);
125
126
127 //- Destructor
128 virtual ~optionList() = default;
129
130
131 // Member Functions
132
133 //- Reset the source list
134 void reset(const dictionary& dict);
135
136 //- Return whether there is something to apply to the field
137 bool appliesToField(const word& fieldName) const;
138
139
140 // Sources
141
142 //- Return source for equation
143 template<class Type>
145 (
147 );
148
149 //- Return source for equation with specified name
150 template<class Type>
152 (
154 const word& fieldName
155 );
156
157 //- Return source for equation
158 template<class Type>
160 (
161 const volScalarField& rho,
163 );
164
165 //- Return source for equation with specified name
166 template<class Type>
168 (
169 const volScalarField& rho,
171 const word& fieldName
172 );
173
174 //- Return source for equation
175 template<class Type>
177 (
178 const volScalarField& alpha,
179 const volScalarField& rho,
181 );
182
183 //- Return source for equation with specified name
184 template<class Type>
186 (
187 const volScalarField& alpha,
188 const volScalarField& rho,
190 const word& fieldName
191 );
192
193 //- Return source for equation
194 template<class Type>
196 (
197 const volScalarField& alpha,
198 const geometricOneField& rho,
200 );
201
202 //- Return source for equation
203 template<class Type>
205 (
207 const volScalarField& rho,
209 );
210
211 //- Return source for equation
212 template<class Type>
214 (
216 const geometricOneField& rho,
218 );
219
220 //- Return source for equation with second time derivative
221 template<class Type>
223 (
225 );
226
227 //- Return source for equation with second time derivative
228 template<class Type>
230 (
232 const word& fieldName
233 );
234
235
236 // Constraints
237
238 //- Apply constraints to equation
239 template<class Type>
240 void constrain(fvMatrix<Type>& eqn);
241
242
243 // Correction
244
245 //- Apply correction to field
246 template<class Type>
248
249
250 //- Post process sensitivity field related to the fvOption
251 template<class Type>
252 void postProcessSens
253 (
254 Field<Type>& sensField,
255 const word& fieldName = word::null,
256 const word& designVariablesName = word::null
257 );
258
259
260 // IO
261
262 //- Read dictionary
263 virtual bool read(const dictionary& dict);
264
265 //- Write data to Ostream
266 virtual bool writeData(Ostream& os) const;
267
268 //- Ostream operator
269 friend Ostream& operator<<
270 (
271 Ostream& os,
272 const optionList& options
273 );
274};
275
276
277// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
278
279} // End namespace fv
280} // End namespace Foam
281
282// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
283
284#ifdef NoRepository
285 #include "fvOptionListTemplates.C"
286#endif
287
288// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
289
290#endif
291
292// ************************************************************************* //
Generic templated field type.
Definition: Field.H:82
Generic GeometricField class.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
friend Ostream & operator(Ostream &os, const UPtrList< T > &list)
Write UPtrList to Ostream.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:109
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
List of finite volume options.
Definition: fvOptionList.H:72
void reset(const dictionary &dict)
Reset the source list.
Definition: fvOptionList.C:104
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOptionList.H:78
virtual ~optionList()=default
Destructor.
virtual bool writeData(Ostream &os) const
Write data to Ostream.
Definition: fvOptionList.C:158
void checkApplied() const
Check that all sources have been applied.
Definition: fvOptionList.C:70
tmp< fvMatrix< Type > > d2dt2(GeometricField< Type, fvPatchField, volMesh > &field)
Return source for equation with second time derivative.
bool readOptions(const dictionary &dict)
Read options dictionary.
Definition: fvOptionList.C:56
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: fvOptionList.C:152
void operator=(const optionList &)=delete
No copy assignment.
TypeName("optionList")
Runtime type information.
tmp< fvMatrix< Type > > d2dt2(GeometricField< Type, fvPatchField, volMesh > &field, const word &fieldName)
Return source for equation with second time derivative.
label checkTimeIndex_
Time index to check that all defined sources have been applied.
Definition: fvOptionList.H:81
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
optionList(const optionList &)=delete
No copy construct.
void postProcessSens(Field< Type > &sensField, const word &fieldName=word::null, const word &designVariablesName=word::null)
Post process sensitivity field related to the fvOption.
tmp< fvMatrix< Type > > source(GeometricField< Type, fvPatchField, volMesh > &field, const word &fieldName, const dimensionSet &ds)
Return source for equation with specified name and dimensions.
bool appliesToField(const word &fieldName) const
Return whether there is something to apply to the field.
Definition: fvOptionList.C:136
static const dictionary & optionsDict(const dictionary &dict)
Return "options" sub-dictionary (if present) or return dict.
Definition: fvOptionList.C:46
Finite-volume options.
Definition: fvOptions.H:59
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
thermo correct()
rDeltaTY field()
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
labelList fv(nPoints)
volScalarField & alpha
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73