fvOptionListTemplates.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2016 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 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
27 
28 template<class Type>
30 (
32 )
33 {
34  return this->operator()(field, field.name());
35 }
36 
37 
38 template<class Type>
40 (
42  const word& fieldName
43 )
44 {
45  checkApplied();
46 
47  const dimensionSet ds = field.dimensions()/dimTime*dimVolume;
48 
49  tmp<fvMatrix<Type>> tmtx(new fvMatrix<Type>(field, ds));
50  fvMatrix<Type>& mtx = tmtx.ref();
51 
52  forAll(*this, i)
53  {
54  option& source = this->operator[](i);
55 
56  label fieldi = source.applyToField(fieldName);
57 
58  if (fieldi != -1)
59  {
60  source.setApplied(fieldi);
61 
62  if (source.isActive())
63  {
64  if (debug)
65  {
66  Info<< "Applying source " << source.name() << " to field "
67  << fieldName << endl;
68  }
69 
70  source.addSup(mtx, fieldi);
71  }
72  }
73  }
74 
75  return tmtx;
76 }
77 
78 
79 template<class Type>
81 (
82  const volScalarField& rho,
84 )
85 {
86  return this->operator()(rho, field, field.name());
87 }
88 
89 
90 template<class Type>
92 (
93  const volScalarField& rho,
95  const word& fieldName
96 )
97 {
98  checkApplied();
99 
100  const dimensionSet ds
101  (
102  rho.dimensions()*field.dimensions()/dimTime*dimVolume
103  );
104 
105  tmp<fvMatrix<Type>> tmtx(new fvMatrix<Type>(field, ds));
106  fvMatrix<Type>& mtx = tmtx.ref();
107 
108  forAll(*this, i)
109  {
110  option& source = this->operator[](i);
111 
112  label fieldi = source.applyToField(fieldName);
113 
114  if (fieldi != -1)
115  {
116  source.setApplied(fieldi);
117 
118  if (source.isActive())
119  {
120  if (debug)
121  {
122  Info<< "Applying source " << source.name() << " to field "
123  << fieldName << endl;
124  }
125 
126  source.addSup(rho, mtx, fieldi);
127  }
128  }
129  }
130 
131  return tmtx;
132 }
133 
134 
135 template<class Type>
137 (
138  const volScalarField& alpha,
139  const volScalarField& rho,
141 )
142 {
143  return this->operator()(alpha, rho, field, field.name());
144 }
145 
146 
147 template<class Type>
149 (
150  const volScalarField& alpha,
151  const volScalarField& rho,
153  const word& fieldName
154 )
155 {
156  checkApplied();
157 
158  const dimensionSet ds
159  (
160  alpha.dimensions()*rho.dimensions()*field.dimensions()
162  );
163 
164  tmp<fvMatrix<Type>> tmtx(new fvMatrix<Type>(field, ds));
165  fvMatrix<Type>& mtx = tmtx.ref();
166 
167  forAll(*this, i)
168  {
169  option& source = this->operator[](i);
170 
171  label fieldi = source.applyToField(fieldName);
172 
173  if (fieldi != -1)
174  {
175  source.setApplied(fieldi);
176 
177  if (source.isActive())
178  {
179  if (debug)
180  {
181  Info<< "Applying source " << source.name() << " to field "
182  << fieldName << endl;
183  }
184 
185  source.addSup(alpha, rho, mtx, fieldi);
186  }
187  }
188  }
189 
190  return tmtx;
191 }
192 
193 
194 template<class Type>
196 (
197  const geometricOneField& alpha,
198  const geometricOneField& rho,
200 )
201 {
202  return this->operator()(field, field.name());
203 }
204 
205 
206 template<class Type>
208 (
209  const volScalarField& alpha,
210  const geometricOneField& rho,
212 )
213 {
215  (
216  IOobject
217  (
218  "one",
219  this->mesh_.time().timeName(),
220  this->mesh_,
221  IOobject::NO_READ,
222  IOobject::NO_WRITE,
223  false
224  ),
225  this->mesh_,
226  dimensionedScalar("one", dimless, 1.0)
227  );
228 
229  return this->operator()(alpha, one, field, field.name());
230 }
231 
232 
233 template<class Type>
235 (
236  const geometricOneField& alpha,
237  const volScalarField& rho,
239 )
240 {
241  return this->operator()(rho, field, field.name());
242 }
243 
244 
245 template<class Type>
247 {
248  checkApplied();
249 
250  forAll(*this, i)
251  {
252  option& source = this->operator[](i);
253 
254  label fieldi = source.applyToField(eqn.psi().name());
255 
256  if (fieldi != -1)
257  {
258  source.setApplied(fieldi);
259 
260  if (source.isActive())
261  {
262  if (debug)
263  {
264  Info<< "Applying constraint " << source.name()
265  << " to field " << eqn.psi().name() << endl;
266  }
267 
268  source.constrain(eqn, fieldi);
269  }
270  }
271  }
272 }
273 
274 
275 template<class Type>
277 (
279 )
280 {
281  const word& fieldName = field.name();
282 
283  forAll(*this, i)
284  {
285  option& source = this->operator[](i);
286 
287  label fieldi = source.applyToField(fieldName);
288 
289  if (fieldi != -1)
290  {
291  source.setApplied(fieldi);
292 
293  if (source.isActive())
294  {
295  if (debug)
296  {
297  Info<< "Correcting source " << source.name()
298  << " for field " << fieldName << endl;
299  }
300 
301  source.correct(field);
302  }
303  }
304  }
305 }
306 
307 
308 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
virtual void correct(volScalarField &field)
Definition: fvOption.C:306
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
const word & name() const
Return name.
Definition: IOobject.H:291
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:284
const T & operator[](const label) const
Return element const reference.
Definition: UPtrListI.H:103
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Generic GeometricField class.
void setApplied(const label fieldi)
Set the applied flag to true for field index fieldi.
Definition: fvOptionI.H:52
friend Ostream & operator(Ostream &, const UPtrList< T > &)
Write UPtrList to Ostream.
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
const dimensionSet & dimensions() const
Return dimensions.
Dimension set for the base types.
Definition: dimensionSet.H:120
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
virtual void addSup(fvMatrix< scalar > &eqn, const label fieldi)
Definition: fvOption.C:134
A class for handling words, derived from string.
Definition: word.H:59
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
const word & name() const
Return const access to the source name.
Definition: fvOptionI.H:28
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
void checkApplied() const
Check that all sources have been applied.
Definition: fvOptionList.C:73
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
virtual void constrain(fvMatrix< scalar > &eqn, const label fieldi)
Definition: fvOption.C:278
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
virtual bool isActive()
Is the source active?
Definition: fvOption.C:107
virtual label applyToField(const word &fieldName) const
Return index of field name if found in fieldNames list.
Definition: fvOption.C:113
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
A class representing the concept of 1 (scalar(1)) used to avoid unnecessary manipulations for objects...
Definition: one.H:50
Finite volume options abstract base class. Provides a base set of controls, e.g.: ...
Definition: fvOption.H:66