faOptionList.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) 2019-2022 OpenCFD Ltd.
9 ------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Class
27  Foam::fa::optionList
28 
29 Description
30  List of finite volume options
31 
32 SourceFile
33  optionList.C
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef Foam_faOptionList_H
38 #define Foam_faOptionList_H
39 
40 #include "faOption.H"
41 #include "PtrList.H"
42 #include "GeometricField.H"
43 #include "geometricOneField.H"
44 #include "faPatchField.H"
45 #include "fvMesh.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 
52 // Forward Declarations
53 
54 namespace fa
55 {
56  class optionList;
57 }
58 
59 Ostream& operator<<(Ostream& os, const fa::optionList& options);
60 
61 namespace fa
62 {
63 
64 /*---------------------------------------------------------------------------*\
65  Class optionList Declaration
66 \*---------------------------------------------------------------------------*/
67 
68 class optionList
69 :
70  public PtrList<fa::option>
71 {
72 protected:
73 
74  // Protected Data
75 
76  //- Reference to the mesh database
77  const fvMesh& mesh_;
78 
79  //- Time index to check that all defined sources have been applied
80  label checkTimeIndex_;
81 
82 
83  // Protected Member Functions
84 
85  //- Return "options" sub-dictionary (if present) or return dict
86  static const dictionary& optionsDict(const dictionary& dict);
87 
88  //- Read options dictionary
89  bool readOptions(const dictionary& dict);
90 
91  //- Check that all sources have been applied
92  void checkApplied() const;
93 
94  //- Return source for equation with specified name and dimensions
95  template<class Type>
97  (
99  const areaScalarField& h,
100  const word& fieldName,
101  const dimensionSet& ds
102  );
103 
104  //- No copy construct
105  optionList(const optionList&) = delete;
106 
107  //- No copy assignment
108  void operator=(const optionList&) = delete;
109 
110 
111 public:
112 
113  //- Runtime type information
114  TypeName("optionList");
115 
116 
117  // Constructors
118 
119  //- Default construct from mesh
120  explicit optionList(const fvMesh& mesh);
121 
122  //- Construct from mesh and dictionary
123  optionList(const fvMesh& mesh, const dictionary& dict);
124 
125 
126  //- Destructor
127  virtual ~optionList() = default;
128 
129 
130  // Member Functions
131 
132  //- Reset the source list
133  void reset(const dictionary& dict);
134 
135  //- Return whether there is something to apply to the field
136  bool appliesToField(const word& fieldName) const;
137 
138 
139  // Sources
140 
141  //- Return source for equation
142  template<class Type>
144  (
145  const areaScalarField& h,
147  );
148 
149  //- Return source for equation with specified name
150  template<class Type>
152  (
153  const areaScalarField& h,
155  const word& fieldName
156  );
157 
158  //- Return source for equation
159  template<class Type>
161  (
162  const areaScalarField& h,
163  const areaScalarField& rho,
165  );
166 
167  //- Return source for equation with specified name
168  template<class Type>
170  (
171  const areaScalarField& h,
172  const areaScalarField& rho,
174  const word& fieldName
175  );
176 
177 
178  //- Return source for equation with specified name and dimensios
179  template<class Type>
181  (
182  const areaScalarField& rho,
184  const dimensionSet& ds
185  );
186 
187 
188  //- Return source for equation with second time derivative
189  template<class Type>
191  (
193  );
194 
195  //- Return source for equation with second time derivative
196  template<class Type>
198  (
200  const word& fieldName
201  );
202 
203 
204  // Constraints
205 
206  //- Apply constraints to equation
207  template<class Type>
208  void constrain(faMatrix<Type>& eqn);
209 
210 
211  // Correction
212 
213  //- Apply correction to field
214  template<class Type>
216 
217 
218  // IO
219 
220  //- Read dictionary
221  virtual bool read(const dictionary& dict);
222 
223  //- Write data to Ostream
224  virtual bool writeData(Ostream& os) const;
225 
226  //- Ostream operator
227  friend Ostream& operator<<
228  (
229  Ostream& os,
230  const optionList& options
231  );
232 };
233 
234 
235 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
236 
237 } // End namespace fa
238 } // End namespace Foam
239 
240 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
241 
242 #ifdef NoRepository
243  #include "faOptionListTemplates.C"
244 #endif
245 
246 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
247 
248 #endif
249 
250 // ************************************************************************* //
label checkTimeIndex_
Time index to check that all defined sources have been applied.
Definition: faOptionList.H:77
dictionary dict
void operator=(const optionList &)=delete
No copy assignment.
rDeltaTY field()
void constrain(faMatrix< Type > &eqn)
Apply constraints to equation.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
void reset(const dictionary &dict)
Reset the source list.
Definition: faOptionList.C:93
List of finite volume options.
Definition: faOptionList.H:61
bool appliesToField(const word &fieldName) const
Return whether there is something to apply to the field.
Definition: faOptionList.C:125
static const dictionary & optionsDict(const dictionary &dict)
Return "options" sub-dictionary (if present) or return dict.
Definition: faOptionList.C:37
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
void correct(GeometricField< Type, faPatchField, areaMesh > &field)
Apply correction to field.
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:105
dynamicFvMesh & mesh
Finite-area options.
Definition: faOptions.H:50
A class for handling words, derived from Foam::string.
Definition: word.H:63
tmp< faMatrix< Type > > source(GeometricField< Type, faPatchField, areaMesh > &field, const areaScalarField &h, const word &fieldName, const dimensionSet &ds)
Return source for equation with specified name and dimensions.
friend Ostream & operator(Ostream &os, const UPtrList< T > &list)
Write UPtrList to Ostream.
optionList(const optionList &)=delete
No copy construct.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
OBJstream os(runTime.globalPath()/outputName)
tmp< faMatrix< Type > > d2dt2(GeometricField< Type, faPatchField, areaMesh > &field)
Return source for equation with second time derivative.
const dimensionedScalar h
Planck constant.
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:77
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A special matrix type and solver, designed for finite area solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: faMatricesFwd.H:37
virtual bool writeData(Ostream &os) const
Write data to Ostream.
Definition: faOptionList.C:147
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: faOptionList.C:141
const fvMesh & mesh_
Reference to the mesh database.
Definition: faOptionList.H:72
A class for managing temporary objects.
Definition: HashPtrTable.H:50
bool readOptions(const dictionary &dict)
Read options dictionary.
Definition: faOptionList.C:47
virtual ~optionList()=default
Destructor.
TypeName("optionList")
Runtime type information.
Namespace for OpenFOAM.
void checkApplied() const
Check that all sources have been applied.
Definition: faOptionList.C:61