All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
fvOptionList.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) 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 \*---------------------------------------------------------------------------*/
25 
26 #include "fvOptionList.H"
27 #include "surfaceFields.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace fv
34 {
35  defineTypeNameAndDebug(optionList, 0);
36 }
37 }
38 
39 
40 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
41 
43 (
44  const dictionary& dict
45 ) const
46 {
47  if (dict.found("options"))
48  {
49  return dict.subDict("options");
50  }
51  else
52  {
53  return dict;
54  }
55 }
56 
57 
59 {
60  checkTimeIndex_ = mesh_.time().timeIndex() + 2;
61 
62  bool allOk = true;
63  forAll(*this, i)
64  {
65  option& bs = this->operator[](i);
66  bool ok = bs.read(dict.subDict(bs.name()));
67  allOk = (allOk && ok);
68  }
69  return allOk;
70 }
71 
72 
74 {
75  if (mesh_.time().timeIndex() == checkTimeIndex_)
76  {
77  forAll(*this, i)
78  {
79  const option& bs = this->operator[](i);
80  bs.checkApplied();
81  }
82  }
83 }
84 
85 
86 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
87 
89 :
90  PtrList<option>(),
91  mesh_(mesh),
92  checkTimeIndex_(mesh_.time().startTimeIndex() + 2)
93 {
94  reset(optionsDict(dict));
95 }
96 
97 
99 :
100  PtrList<option>(),
101  mesh_(mesh),
102  checkTimeIndex_(mesh_.time().startTimeIndex() + 2)
103 {}
104 
105 
106 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
107 
109 {
110  // Count number of active fvOptions
111  label count = 0;
112  forAllConstIter(dictionary, dict, iter)
113  {
114  if (iter().isDict())
115  {
116  count++;
117  }
118  }
119 
120  this->setSize(count);
121  label i = 0;
122  forAllConstIter(dictionary, dict, iter)
123  {
124  if (iter().isDict())
125  {
126  const word& name = iter().keyword();
127  const dictionary& sourceDict = iter().dict();
128 
129  this->set
130  (
131  i++,
132  option::New(name, sourceDict, mesh_)
133  );
134  }
135  }
136 }
137 
138 
139 bool Foam::fv::optionList::appliesToField(const word& fieldName) const
140 {
141  forAll(*this, i)
142  {
143  const option& source = this->operator[](i);
144 
145  label fieldi = source.applyToField(fieldName);
146 
147  if (fieldi != -1)
148  {
149  return true;
150  }
151  }
152 
153  return false;
154 }
155 
156 
158 {
159  return readOptions(optionsDict(dict));
160 }
161 
162 
164 {
165  // Write list contents
166  forAll(*this, i)
167  {
168  os << nl;
169  this->operator[](i).writeHeader(os);
170  this->operator[](i).writeData(os);
171  this->operator[](i).writeFooter(os);
172  }
173 
174  // Check state of IOstream
175  return os.good();
176 }
177 
178 
179 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
180 
182 {
183  options.writeData(os);
184  return os;
185 }
186 
187 
188 // ************************************************************************* //
defineTypeNameAndDebug(option, 0)
Foam::surfaceFields.
bool readOptions(const dictionary &dict)
Read options dictionary.
Definition: fvOptionList.C:58
dictionary dict
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:667
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
const dictionary & optionsDict(const dictionary &dict) const
Return the "options" sub-dictionary if present otherwise return dict.
Definition: fvOptionList.C:43
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOptionList.H:73
const T & operator[](const label) const
Return element const reference.
Definition: UPtrListI.H:96
tmp< fvMatrix< Type > > source(GeometricField< Type, fvPatchField, volMesh > &field, const word &fieldName, const dimensionSet &ds)
Return source for equation with specified name and dimensions.
label checkTimeIndex_
Time index to check that all defined sources have been applied.
Definition: fvOptionList.H:76
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:333
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:934
static autoPtr< option > New(const word &name, const dictionary &dict, const fvMesh &mesh)
Return a reference to the selected fvOption model.
Definition: fvOption.C:67
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
optionList(const fvMesh &mesh)
Construct null.
Definition: fvOptionList.C:98
virtual void checkApplied() const
Check that the source has been applied.
Definition: fvOption.C:119
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvOptionIO.C:53
const word & name() const
Return const access to the source name.
Definition: fvOptionI.H:28
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
bool appliesToField(const word &fieldName) const
Return whether there is something to apply to the field.
Definition: fvOptionList.C:139
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
void reset(const dictionary &dict)
Reset the source list.
Definition: fvOptionList.C:108
static const char nl
Definition: Ostream.H:260
void checkApplied() const
Check that all sources have been applied.
Definition: fvOptionList.C:73
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
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
Ostream & operator<<(Ostream &, const ensightPart &)
List of finite volume options.
Definition: fvOptionList.H:64
virtual label applyToField(const word &fieldName) const
Return index of field name if found in fieldNames list.
Definition: fvOption.C:113
virtual bool writeData(Ostream &os) const
Write data to Ostream.
Definition: fvOptionList.C:163
Namespace for OpenFOAM.
Finite volume options abstract base class. Provides a base set of controls, e.g.: ...
Definition: fvOption.H:66
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: fvOptionList.C:157