fvOptionList.H
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) 2011-2019 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 Class
25  Foam::fv::optionList
26 
27 Description
28  List of finite volume options
29 
30 SourceFile
31  optionList.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef fvOptionList_H
36 #define fvOptionList_H
37 
38 #include "fvOption.H"
39 #include "PtrList.H"
40 #include "GeometricField.H"
41 #include "geometricOneField.H"
42 #include "fvPatchField.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 
49 // Forward declaration of friend functions and operators
50 
51 namespace fv
52 {
53  class optionList;
54 }
55 
56 Ostream& operator<<(Ostream& os, const fv::optionList& options);
57 
58 namespace fv
59 {
60 
61 /*---------------------------------------------------------------------------*\
62  Class optionList Declaration
63 \*---------------------------------------------------------------------------*/
64 
65 class optionList
66 :
67  public PtrList<option>
68 {
69 protected:
70 
71  // Protected data
72 
73  //- Reference to the mesh database
74  const fvMesh& mesh_;
75 
76  //- Time index to check that all defined sources have been applied
77  label checkTimeIndex_;
78 
79 
80  // Protected Member Functions
81 
82  //- Return the "options" sub-dictionary if present otherwise return dict
83  const dictionary& optionsDict(const dictionary& dict) const;
84 
85  //- Read options dictionary
86  bool readOptions(const dictionary& dict);
87 
88  //- Check that all sources have been applied
89  void checkApplied() const;
90 
91  //- Return source for equation with specified name and dimensions
92  template<class Type>
93  tmp<fvMatrix<Type>> source
94  (
96  const word& fieldName,
97  const dimensionSet& ds
98  );
99 
100 
101 public:
102 
103  //- Runtime type information
104  TypeName("optionList");
105 
106 
107  // Constructors
108 
109  //- Construct null
110  optionList(const fvMesh& mesh);
111 
112  //- Construct from mesh and dictionary
113  optionList(const fvMesh& mesh, const dictionary& dict);
114 
115  //- Disallow default bitwise copy construction
116  optionList(const optionList&) = delete;
117 
118 
119  //- Destructor
120  virtual ~optionList()
121  {}
122 
123 
124  // Member Functions
125 
126  //- Reset the source list
127  void reset(const dictionary& dict);
128 
129  //- Return whether there is something to apply to the field
130  bool appliesToField(const word& fieldName) const;
131 
132 
133  // Sources
134 
135  //- Return source for equation
136  template<class Type>
137  tmp<fvMatrix<Type>> operator()
138  (
140  );
141 
142  //- Return source for equation with specified name
143  template<class Type>
144  tmp<fvMatrix<Type>> operator()
145  (
147  const word& fieldName
148  );
149 
150  //- Return source for equation
151  template<class Type>
152  tmp<fvMatrix<Type>> operator()
153  (
154  const volScalarField& rho,
156  );
157 
158  //- Return source for equation with specified name
159  template<class Type>
160  tmp<fvMatrix<Type>> operator()
161  (
162  const volScalarField& rho,
164  const word& fieldName
165  );
166 
167  //- Return source for equation
168  template<class Type>
169  tmp<fvMatrix<Type>> operator()
170  (
171  const volScalarField& alpha,
172  const volScalarField& rho,
174  );
175 
176  //- Return source for equation with specified name
177  template<class Type>
178  tmp<fvMatrix<Type>> operator()
179  (
180  const volScalarField& alpha,
181  const volScalarField& rho,
183  const word& fieldName
184  );
185 
186  //- Return source for equation
187  template<class Type>
188  tmp<fvMatrix<Type>> operator()
189  (
190  const volScalarField& alpha,
191  const geometricOneField& rho,
193  );
194 
195  //- Return source for equation
196  template<class Type>
197  tmp<fvMatrix<Type>> operator()
198  (
199  const geometricOneField& alpha,
200  const volScalarField& rho,
202  );
203 
204  //- Return source for equation
205  template<class Type>
206  tmp<fvMatrix<Type>> operator()
207  (
208  const geometricOneField& alpha,
209  const geometricOneField& rho,
211  );
212 
213  //- Return source for equation with second time derivative
214  template<class Type>
216  (
218  );
219 
220  //- Return source for equation with second time derivative
221  template<class Type>
223  (
225  const word& fieldName
226  );
227 
228 
229  // Constraints
230 
231  //- Apply constraints to equation
232  template<class Type>
233  void constrain(fvMatrix<Type>& eqn);
234 
235 
236  // Correction
237 
238  //- Apply correction to field
239  template<class Type>
241 
242 
243  // IO
244 
245  //- Read dictionary
246  virtual bool read(const dictionary& dict);
247 
248  //- Write data to Ostream
249  virtual bool writeData(Ostream& os) const;
250 
251  //- Ostream operator
252  friend Ostream& operator<<
253  (
254  Ostream& os,
255  const optionList& options
256  );
257 
258 
259  // Member Operators
260 
261  //- Disallow default bitwise assignment
262  void operator=(const optionList&) = delete;
263 };
264 
265 
266 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
267 
268 } // End namespace fv
269 } // End namespace Foam
270 
271 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
272 
273 #ifdef NoRepository
274  #include "fvOptionListTemplates.C"
275 #endif
276 
277 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
278 
279 #endif
280 
281 // ************************************************************************* //
dictionary dict
tmp< GeometricField< Type, fvPatchField, volMesh > > d2dt2(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcD2dt2.C:45
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
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
Generic GeometricField class.
Finite-volume options.
Definition: fvOptions.H:52
Dimension set for the base types.
Definition: dimensionSet.H:120
bool read(const char *, int32_t &)
Definition: int32IO.C:85
dynamicFvMesh & mesh
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:72
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:70
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-7/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:68
const bool writeData(readBool(pdfDictionary.lookup("writeData")))
Ostream & operator<<(Ostream &, const SemiImplicitSource< Type > &)
fvOptions constrain(rhoEqn)
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
List of finite volume options.
Definition: fvOptionList.H:64
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.