faOptionListTemplates.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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2019-2021 OpenCFD Ltd.
9------------------------------------------------------------------------------
10License
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\*---------------------------------------------------------------------------*/
27
28#include "profiling.H"
29#include "areaFields.H"
30
31// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
32
33template<class Type>
35(
37 const areaScalarField& h,
38 const word& fieldName,
39 const dimensionSet& ds
40)
41{
43
45 faMatrix<Type>& mtx = tmtx.ref();
46
47 for (fa::option& source : *this)
48 {
49 const label fieldi = source.applyToField(fieldName);
50
51 if (fieldi != -1)
52 {
53 addProfiling(faopt, "faOption()." + source.name());
54
55 source.setApplied(fieldi);
56
57 const bool ok = source.isActive();
58
59 if (debug)
60 {
61 if (ok)
62 {
63 Info<< "Apply";
64 }
65 else
66 {
67 Info<< "(Inactive)";
68 }
69 Info<< " source " << source.name()
70 << " for field " << fieldName << endl;
71 }
72
73 if (ok)
74 {
75 source.addSup(h, mtx, fieldi);
76 }
77 }
78 }
79
80 return tmtx;
81}
82
83
84template<class Type>
86(
87 const areaScalarField& h,
89)
90{
91 return this->operator()(h, field, field.name());
92}
93
94
95template<class Type>
97(
98 const areaScalarField& h,
100 const word& fieldName
101)
102{
103 return source(field, h, fieldName, field.dimensions()/dimTime*dimArea);
104}
105
106
107template<class Type>
109(
110 const areaScalarField& h,
111 const areaScalarField& rho,
113)
114{
115 return this->operator()(h, rho, field, field.name());
116}
117
118
119template<class Type>
121(
122 const areaScalarField& h,
123 const areaScalarField& rho,
125 const word& fieldName
126)
127{
128 checkApplied();
129
130 const dimensionSet ds
131 (
132 rho.dimensions()*field.dimensions()/dimTime*dimArea
133 );
134
136 faMatrix<Type>& mtx = tmtx.ref();
137
138 for (fa::option& source : *this)
139 {
140 const label fieldi = source.applyToField(fieldName);
141
142 if (fieldi != -1)
143 {
144 addProfiling(faopt, "faOption()." + source.name());
145
146 source.setApplied(fieldi);
147
148 const bool ok = source.isActive();
149
150 if (debug)
151 {
152 if (ok)
153 {
154 Info<< "Apply";
155 }
156 else
157 {
158 Info<< "(Inactive)";
159 }
160 Info<< " source " << source.name()
161 << " for field " << fieldName << endl;
162 }
163
164 if (ok)
165 {
166 source.addSup(h, rho, mtx, fieldi);
167 }
168 }
169 }
170
171 return tmtx;
172}
173
174
175template<class Type>
177(
178 const areaScalarField& rho,
180 const dimensionSet& ds
181)
182{
183 checkApplied();
184
185 const dimensionSet dsMat(ds*dimArea);
186
187 tmp<faMatrix<Type>> tmtx(new faMatrix<Type>(field, dsMat));
188 faMatrix<Type>& mtx = tmtx.ref();
189
190 for (fa::option& source : *this)
191 {
192 const label fieldi = source.applyToField(field.name());
193
194 if (fieldi != -1)
195 {
196 addProfiling(faopt, "faOption()." + source.name());
197
198 source.setApplied(fieldi);
199
200 const bool ok = source.isActive();
201
202 if (debug)
203 {
204 if (ok)
205 {
206 Info<< "Apply";
207 }
208 else
209 {
210 Info<< "(Inactive)";
211 }
212 Info<< " source " << source.name()
213 << " for field " << field.name() << endl;
214 }
215
216 if (ok)
217 {
218 source.addSup(rho, mtx, fieldi);
219 }
220 }
221 }
222
223 return tmtx;
224}
225
226
227template<class Type>
229(
231)
232{
233 return this->d2dt2(field, field.name());
234}
235
236
237template<class Type>
239(
241 const word& fieldName
242)
243{
244 return source(field, fieldName, field.dimensions()/sqr(dimTime)*dimArea);
245}
246
247
248template<class Type>
250{
251 checkApplied();
252
253 for (fa::option& source : *this)
254 {
255 const label fieldi = source.applyToField(eqn.psi().name());
256
257 if (fieldi != -1)
258 {
259 addProfiling(faopt, "faOption::constrain." + eqn.psi().name());
260
261 source.setApplied(fieldi);
262
263 const bool ok = source.isActive();
264
265 if (debug)
266 {
267 if (ok)
268 {
269 Info<< "Constrain";
270 }
271 else
272 {
273 Info<< "(Inactive constrain)";
274 }
275 Info<< " source " << source.name()
276 << " for field " << eqn.psi().name() << endl;
277 }
278
279 if (ok)
280 {
281 source.constrain(eqn, fieldi);
282 }
283 }
284 }
285}
286
287
288template<class Type>
290(
292)
293{
294 const word& fieldName = field.name();
295
296 for (fa::option& source : *this)
297 {
298 const label fieldi = source.applyToField(fieldName);
299
300 if (fieldi != -1)
301 {
302 addProfiling(faopt, "faOption::correct." + source.name());
303
304 source.setApplied(fieldi);
305
306 const bool ok = source.isActive();
307
308 if (debug)
309 {
310 if (ok)
311 {
312 Info<< "Correct";
313 }
314 else
315 {
316 Info<< "(Inactive correct)";
317 }
318 Info<< " source " << source.name()
319 << " for field " << fieldName << endl;
320 }
321
322 if (ok)
323 {
324 source.correct(field);
325 }
326 }
327 }
328}
329
330
331// ************************************************************************* //
Generic GeometricField class.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:109
A special matrix type and solver, designed for finite area solutions of scalar equations....
Definition: faMatrix.H:76
const GeometricField< Type, faPatchField, areaMesh > & psi() const
Definition: faMatrix.H:270
friend Ostream & operator(Ostream &, const faMatrix< Type > &)
Field< Type > & source() noexcept
Definition: faMatrix.H:280
void checkApplied() const
Check that all sources have been applied.
Definition: faOptionList.C:68
void constrain(faMatrix< Type > &eqn)
Apply constraints to equation.
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.
tmp< faMatrix< Type > > d2dt2(GeometricField< Type, faPatchField, areaMesh > &field)
Return source for equation with second time derivative.
Base abstract class for handling finite area options (i.e. faOption).
Definition: faOption.H:134
const word & name() const noexcept
Return const access to the source name.
Definition: faOptionI.H:30
void setApplied(const label fieldi)
Set the applied flag to true for field index fieldi.
Definition: faOptionI.H:68
virtual void correct(areaScalarField &field)
Definition: faOption.C:269
virtual void addSup(const areaScalarField &h, faMatrix< scalar > &eqn, const label fieldi)
Definition: faOption.C:147
virtual label applyToField(const word &fieldName) const
Return index of field name if found in fieldNames list.
Definition: faOption.C:126
virtual bool isActive()
Is the source active?
Definition: faOption.C:120
virtual void constrain(faMatrix< scalar > &eqn, const label fieldi)
Definition: faOption.C:241
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
A class for handling words, derived from Foam::string.
Definition: word.H:68
rDeltaTY field()
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
#define addProfiling(name, descr)
Define profiling trigger with specified name and description string.
volScalarField & h