Loading...
Searching...
No Matches
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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2016-2021 OpenCFD Ltd.
10 Copyright (C) 2020,2023 PCOpt/NTUA
11 Copyright (C) 2020,2023 FOSS GP
12-------------------------------------------------------------------------------
13License
14 This file is part of OpenFOAM.
15
16 OpenFOAM is free software: you can redistribute it and/or modify it
17 under the terms of the GNU General Public License as published by
18 the Free Software Foundation, either version 3 of the License, or
19 (at your option) any later version.
20
21 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
22 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
23 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
24 for more details.
25
26 You should have received a copy of the GNU General Public License
27 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
28
29\*---------------------------------------------------------------------------*/
30
31#include "profiling.H"
32
33// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
34
35template<class Type>
37(
39 const word& fieldName,
40 const dimensionSet& ds
41)
42{
44
45 auto tmtx = tmp<fvMatrix<Type>>::New(field, ds);
46 auto& mtx = tmtx.ref();
47
48 for (fv::option& source : *this)
49 {
50 const label fieldi = source.applyToField(fieldName);
51
52 if (fieldi != -1)
53 {
54 addProfiling(fvopt, "fvOption().", source.name());
55
56 source.setApplied(fieldi);
57
58 const bool ok = source.isActive();
59
60 if (debug)
61 {
62 if (ok)
63 {
64 Info<< "Apply";
65 }
66 else
67 {
68 Info<< "(Inactive)";
69 }
70 Info<< " source " << source.name()
71 << " for field " << fieldName << endl;
72 }
73
74 if (ok)
75 {
76 source.addSup(mtx, fieldi);
77 }
78 }
79 }
80
81 return tmtx;
82}
83
84
85template<class Type>
86Foam::tmp<Foam::fvMatrix<Type>> Foam::fv::optionList::operator()
87(
90{
91 return this->operator()(field, field.name());
92}
93
94
95template<class Type>
96Foam::tmp<Foam::fvMatrix<Type>> Foam::fv::optionList::operator()
97(
99 const word& fieldName
101{
102 return source(field, fieldName, field.dimensions()/dimTime*dimVolume);
103}
104
105
106template<class Type>
107Foam::tmp<Foam::fvMatrix<Type>> Foam::fv::optionList::operator()
108(
109 const volScalarField& rho,
112{
113 return this->operator()(rho, field, field.name());
114}
115
116
117template<class Type>
118Foam::tmp<Foam::fvMatrix<Type>> Foam::fv::optionList::operator()
119(
120 const volScalarField& rho,
122 const word& fieldName
123)
124{
125 checkApplied();
126
127 const dimensionSet ds
128 (
129 rho.dimensions()*field.dimensions()/dimTime*dimVolume
130 );
131
132 auto tmtx = tmp<fvMatrix<Type>>::New(field, ds);
133 auto& mtx = tmtx.ref();
134
135 for (fv::option& source : *this)
136 {
137 const label fieldi = source.applyToField(fieldName);
138
139 if (fieldi != -1)
140 {
141 addProfiling(fvopt, "fvOption().", source.name());
142
143 source.setApplied(fieldi);
144
145 const bool ok = source.isActive();
146
147 if (debug)
148 {
149 if (ok)
150 {
151 Info<< "Apply";
152 }
153 else
154 {
155 Info<< "(Inactive)";
156 }
157 Info<< " source " << source.name()
158 << " for field " << fieldName << endl;
159 }
160
161 if (ok)
162 {
163 source.addSup(rho, mtx, fieldi);
164 }
165 }
167
168 return tmtx;
169}
170
171
172template<class Type>
173Foam::tmp<Foam::fvMatrix<Type>> Foam::fv::optionList::operator()
174(
175 const volScalarField& alpha,
176 const volScalarField& rho,
179{
180 return this->operator()(alpha, rho, field, field.name());
181}
182
183
184template<class Type>
185Foam::tmp<Foam::fvMatrix<Type>> Foam::fv::optionList::operator()
186(
187 const volScalarField& alpha,
188 const volScalarField& rho,
190 const word& fieldName
191)
192{
193 checkApplied();
194
195 const dimensionSet ds
196 (
197 alpha.dimensions()*rho.dimensions()*field.dimensions()
199 );
200
201 auto tmtx = tmp<fvMatrix<Type>>::New(field, ds);
202 auto& mtx = tmtx.ref();
203
204 for (fv::option& source : *this)
205 {
206 const label fieldi = source.applyToField(fieldName);
207
208 if (fieldi != -1)
209 {
210 addProfiling(fvopt, "fvOption().", source.name());
211
212 source.setApplied(fieldi);
213
214 const bool ok = source.isActive();
215
216 if (debug)
217 {
218 if (ok)
219 {
220 Info<< "Apply";
221 }
222 else
223 {
224 Info<< "(Inactive)";
225 }
226 Info<< " source " << source.name()
227 << " for field " << fieldName << endl;
228 }
229
230 if (ok)
231 {
232 source.addSup(alpha, rho, mtx, fieldi);
233 }
234 }
236
237 return tmtx;
238}
239
240
241template<class Type>
242Foam::tmp<Foam::fvMatrix<Type>> Foam::fv::optionList::operator()
243(
245 const geometricOneField& rho,
248{
249 return this->operator()(field, field.name());
250}
251
252
253template<class Type>
254Foam::tmp<Foam::fvMatrix<Type>> Foam::fv::optionList::operator()
255(
256 const volScalarField& alpha,
257 const geometricOneField& rho,
259)
260{
262 (
264 (
265 "one",
266 this->mesh_.time().timeName(),
267 this->mesh_,
271 ),
272 this->mesh_,
273 dimensionedScalar("one", dimless, scalar(1))
274 );
275
276 return this->operator()(alpha, one, field, field.name());
277}
278
279
280template<class Type>
281Foam::tmp<Foam::fvMatrix<Type>> Foam::fv::optionList::operator()
282(
284 const volScalarField& rho,
287{
288 return this->operator()(rho, field, field.name());
289}
290
291
292template<class Type>
294(
297{
298 return this->d2dt2(field, field.name());
299}
300
301
302template<class Type>
304(
306 const word& fieldName
308{
309 return source(field, fieldName, field.dimensions()/sqr(dimTime)*dimVolume);
310}
311
312
313template<class Type>
315{
316 checkApplied();
317
318 for (fv::option& source : *this)
319 {
320 const label fieldi = source.applyToField(eqn.psi().name());
321
322 if (fieldi != -1)
323 {
324 addProfiling(fvopt, "fvOption::constrain.", eqn.psi().name());
325
326 source.setApplied(fieldi);
327
328 const bool ok = source.isActive();
329
330 if (debug)
331 {
332 if (ok)
333 {
334 Info<< "Constrain";
335 }
336 else
337 {
338 Info<< "(Inactive constrain)";
339 }
340 Info<< " source " << source.name()
341 << " for field " << eqn.psi().name() << endl;
342 }
343
344 if (ok)
345 {
346 source.constrain(eqn, fieldi);
348 }
349 }
350}
351
352
353template<class Type, template<class> class PatchField, class GeoMesh>
355(
356 GeometricField<Type, PatchField, GeoMesh>& field
357)
358{
359 const word& fieldName = field.name();
360
361 for (fv::option& source : *this)
362 {
363 const label fieldi = source.applyToField(fieldName);
364
365 if (fieldi != -1)
366 {
367 addProfiling(fvopt, "fvOption::correct.", source.name());
368
369 source.setApplied(fieldi);
370
371 const bool ok = source.isActive();
372
373 if (debug)
374 {
375 if (ok)
376 {
377 Info<< "Correct";
378 }
379 else
380 {
381 Info<< "(Inactive correct)";
382 }
383 Info<< " source " << source.name()
384 << " for field " << fieldName << endl;
385 }
386
387 if (ok)
388 {
389 source.correct(field);
391 }
392 }
393}
394
395
396template<class Type>
398(
399 Field<Type>& sensField,
400 const word& fieldName,
401 const word& designVariablesName
402)
403{
404 for (fv::option& source : *this)
405 {
406 const label fieldi = source.applyToField(fieldName);
407
408 if (fieldi != -1)
409 {
410 addProfiling(fvopt, "fvOption::postProcessSens.", source.name());
411
412 const bool ok = source.isActive();
413
414 if (debug && ok)
415 {
416 Info<< "Post processing sensitivity source "
417 << source.name() << " for field " << fieldName << endl;
418 }
419
420 if (ok)
421 {
422 source.postProcessSens
423 (
424 sensField,
425 fieldName,
426 designVariablesName
427 );
429 }
430 }
431}
432
433
434template<class Type, template<class> class PatchField, class GeoMesh>
436(
437 const GeometricField<Type, PatchField, GeoMesh>& primal,
438 const GeometricField<Type, PatchField, GeoMesh>& adjoint,
439 scalarField& sensField,
440 const word& fieldName
441)
442{
443 for (fv::option& source : *this)
444 {
445 const label fieldi = source.applyToField(fieldName);
446
447 if (fieldi != -1)
448 {
450 (
451 fvopt,
452 "fvOption::postProcessAuxSens.", source.name()
453 );
454
455 const bool ok = source.isActive();
456
457 if (debug && ok)
458 {
459 Info<< "Post processing sensitivity source "
460 << source.name() << " for field " << fieldName << endl;
461 }
462
463 if (ok)
464 {
465 source.postProcessAuxSens
466 (
467 primal,
468 adjoint,
469 sensField,
470 fieldName
471 );
472 }
473 }
474 }
475}
476
477
478// ************************************************************************* //
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition GeoMesh.H:46
Generic GeometricField class.
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
const fvMesh & mesh_
Reference to the mesh database.
void checkApplied() const
Check that all sources have been applied.
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.
tmp< faMatrix< Type > > operator()(const areaScalarField &h, GeometricField< Type, faPatchField, areaMesh > &field)
Return source for equation.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition fvMatrix.H:118
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition fvMatrix.H:487
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
void checkApplied() const
Check that all sources have been applied.
tmp< fvMatrix< Type > > d2dt2(GeometricField< Type, fvPatchField, volMesh > &field)
Return source for equation with second time derivative.
void postProcessAuxSens(const GeometricField< Type, PatchField, GeoMesh > &primal, const GeometricField< Type, PatchField, GeoMesh > &adjoint, scalarField &sensField, const word &fieldName=word::null)
Post process auxiliary sensitivity field related to the fvOption.
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
void postProcessSens(Field< Type > &sensField, const word &fieldName=word::null, const word &designVariablesName=word::null)
Post process sensitivity field related to the fvOption.
tmp< fvMatrix< Type > > source(GeometricField< Type, fvPatchField, volMesh > &field, const word &fieldName, const dimensionSet &ds)
Return source for equation with specified name and dimensions.
Base abstract class for handling finite volume options (i.e. fvOption).
Definition fvOption.H:124
const word & name() const noexcept
Return const access to the source name.
Definition fvOptionI.H:24
void setApplied(const label fieldi)
Set the applied flag to true for field index fieldi.
Definition fvOptionI.H:56
virtual label applyToField(const word &fieldName) const
Return index of field name if found in fieldNames list.
Definition fvOption.C:121
virtual void addSup(fvMatrix< scalar > &eqn, const label fieldi)
Definition fvOption.C:142
virtual void constrain(fvMatrix< scalar > &eqn, const label fieldi)
Definition fvOption.C:286
virtual bool isActive()
Is the source active?
Definition fvOption.C:115
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition one.H:57
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
rDeltaTY field()
Namespace for handling debugging switches.
Definition debug.C:45
int debug
Static debugging option.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const dimensionSet dimVolume(pow3(dimLength))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
#define addProfiling(Name,...)
Define profiling trigger with specified name and description string. The description is generated by ...
volScalarField & alpha