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-2018 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  //- Disallow default bitwise copy construct
101  optionList(const optionList&);
102 
103  //- Disallow default bitwise assignment
104  void operator=(const optionList&);
105 
106 
107 public:
108 
109  //- Runtime type information
110  TypeName("optionList");
111 
112 
113  // Constructors
114 
115  //- Construct null
116  optionList(const fvMesh& mesh);
117 
118  //- Construct from mesh and dictionary
119  optionList(const fvMesh& mesh, const dictionary& dict);
120 
121 
122  //- Destructor
123  virtual ~optionList()
124  {}
125 
126 
127  // Member Functions
128 
129  //- Reset the source list
130  void reset(const dictionary& dict);
131 
132  //- Return whether there is something to apply to the field
133  bool appliesToField(const word& fieldName) const;
134 
135 
136  // Sources
137 
138  //- Return source for equation
139  template<class Type>
140  tmp<fvMatrix<Type>> operator()
141  (
143  );
144 
145  //- Return source for equation with specified name
146  template<class Type>
147  tmp<fvMatrix<Type>> operator()
148  (
150  const word& fieldName
151  );
152 
153  //- Return source for equation
154  template<class Type>
155  tmp<fvMatrix<Type>> operator()
156  (
157  const volScalarField& rho,
159  );
160 
161  //- Return source for equation with specified name
162  template<class Type>
163  tmp<fvMatrix<Type>> operator()
164  (
165  const volScalarField& rho,
167  const word& fieldName
168  );
169 
170  //- Return source for equation
171  template<class Type>
172  tmp<fvMatrix<Type>> operator()
173  (
174  const volScalarField& alpha,
175  const volScalarField& rho,
177  );
178 
179  //- Return source for equation with specified name
180  template<class Type>
181  tmp<fvMatrix<Type>> operator()
182  (
183  const volScalarField& alpha,
184  const volScalarField& rho,
186  const word& fieldName
187  );
188 
189  //- Return source for equation
190  template<class Type>
191  tmp<fvMatrix<Type>> operator()
192  (
193  const volScalarField& alpha,
194  const geometricOneField& rho,
196  );
197 
198  //- Return source for equation
199  template<class Type>
200  tmp<fvMatrix<Type>> operator()
201  (
202  const geometricOneField& alpha,
203  const volScalarField& rho,
205  );
206 
207  //- Return source for equation
208  template<class Type>
209  tmp<fvMatrix<Type>> operator()
210  (
211  const geometricOneField& alpha,
212  const geometricOneField& rho,
214  );
215 
216  //- Return source for equation with second time derivative
217  template<class Type>
219  (
221  );
222 
223  //- Return source for equation with second time derivative
224  template<class Type>
226  (
228  const word& fieldName
229  );
230 
231 
232  // Constraints
233 
234  //- Apply constraints to equation
235  template<class Type>
236  void constrain(fvMatrix<Type>& eqn);
237 
238 
239  // Correction
240 
241  //- Apply correction to field
242  template<class Type>
244 
245 
246  // IO
247 
248  //- Read dictionary
249  virtual bool read(const dictionary& dict);
250 
251  //- Write data to Ostream
252  virtual bool writeData(Ostream& os) const;
253 
254  //- Ostream operator
255  friend Ostream& operator<<
256  (
257  Ostream& os,
258  const optionList& options
259  );
260 };
261 
262 
263 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
264 
265 } // End namespace fv
266 } // End namespace Foam
267 
268 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
269 
270 #ifdef NoRepository
271  #include "fvOptionListTemplates.C"
272 #endif
273 
274 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
275 
276 #endif
277 
278 // ************************************************************************* //
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:137
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(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil) - psi) *pSat, rhoMin);# 1 "/home/ubuntu/OpenFOAM-6/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:72
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:63
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
List of finite volume options.
Definition: fvOptionList.H:64
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Namespace for OpenFOAM.