faOption.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 "faOption.H"
29#include "areaFields.H"
30
31// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32
33namespace Foam
34{
35 namespace fa
36 {
39 }
40}
41
42
43// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
44
46{
48 applied_ = false;
49}
50
51
52// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53
55(
56 const word& name,
57 const word& modelType,
58 const dictionary& dict,
59 const fvPatch& patch
60)
61:
62 name_(name),
63 modelType_(modelType),
64 mesh_(patch.boundaryMesh().mesh()),
65 patch_(patch),
66 dict_(dict),
67 coeffs_(dict.optionalSubDict(modelType + "Coeffs")),
68 fieldNames_(),
69 applied_(),
70 regionName_(dict.get<word>("region")),
71 regionMeshPtr_(nullptr),
72 vsmPtr_(nullptr),
73 active_(dict.getOrDefault("active", true)),
74 log(true)
75{
76 Log << incrIndent << indent << "Source: " << name_ << endl << decrIndent;
77}
78
79
80// * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
81
83(
84 const word& name,
85 const dictionary& coeffs,
86 const fvPatch& patch
87)
88{
89 const word modelType(coeffs.get<word>("type"));
90
91 Info<< indent
92 << "Selecting finite area options type " << modelType << endl;
93
94 patch.boundaryMesh().mesh().time().libs().open
95 (
96 coeffs,
97 "libs",
98 dictionaryConstructorTablePtr_
99 );
100
101 auto* ctorPtr = dictionaryConstructorTable(modelType);
102
103 if (!ctorPtr)
104 {
106 (
107 coeffs,
108 "faOption",
109 modelType,
110 *dictionaryConstructorTablePtr_
111 ) << exit(FatalIOError);
112 }
113
114 return autoPtr<fa::option>(ctorPtr(name, modelType, coeffs, patch));
115}
116
117
118// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
119
121{
122 return active_;
123}
124
125
126Foam::label Foam::fa::option::applyToField(const word& fieldName) const
127{
128 return fieldNames_.find(fieldName);
129}
130
131
133{
134 forAll(applied_, i)
135 {
136 if (!applied_[i])
137 {
139 << "Source " << name_ << " defined for field "
140 << fieldNames_[i] << " but never used" << endl;
141 }
142 }
143}
144
145
147(
148 const areaScalarField& h,
149 faMatrix<scalar>& eqn,
150 const label fieldi
151)
152{}
153
154
156(
157 const areaScalarField& h,
158 faMatrix<vector>& eqn,
159 const label fieldi
160)
161{}
162
163
165(
166 const areaScalarField& h,
168 const label fieldi
169)
170{}
171
172
174(
175 const areaScalarField& h,
177 const label fieldi
178)
179{}
180
181
183(
184 const areaScalarField& h,
185 faMatrix<tensor>& eqn,
186 const label fieldi
187)
188{}
189
190
192(
193 const areaScalarField& h,
194 const areaScalarField& rho,
195 faMatrix<scalar>& eqn,
196 const label fieldi
197)
198{}
199
200
202(
203 const areaScalarField& h,
204 const areaScalarField& rho,
205 faMatrix<vector>& eqn,
206 const label fieldi
207)
208{}
209
210
212(
213 const areaScalarField& h,
214 const areaScalarField& rho,
216 const label fieldi
217)
218{}
219
220
222(
223 const areaScalarField& h,
224 const areaScalarField& rho,
226 const label fieldi
227)
228{}
229
230
232(
233 const areaScalarField& h,
234 const areaScalarField& rho,
235 faMatrix<tensor>& eqn,
236 const label fieldi
237)
238{}
239
240
241void Foam::fa::option::constrain(faMatrix<scalar>& eqn, const label fieldi)
242{}
243
244
245void Foam::fa::option::constrain(faMatrix<vector>& eqn, const label fieldi)
246{}
247
248
250(
252 const label fieldi
253)
254{}
255
256
258(
260 const label fieldi
261)
262{}
263
264
265void Foam::fa::option::constrain(faMatrix<tensor>& eqn, const label fieldi)
266{}
267
268
270{}
271
272
274{}
275
276
278{}
279
280
282{}
283
284
286{}
287
288
289// ************************************************************************* //
#define Log
Definition: PDRblock.C:35
void resize_nocopy(const label len)
Adjust allocated size of list without necessarily.
Definition: ListI.H:146
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:63
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
A special matrix type and solver, designed for finite area solutions of scalar equations....
Definition: faMatrix.H:76
Base abstract class for handling finite area options (i.e. faOption).
Definition: faOption.H:134
List< bool > applied_
Applied flag list - corresponds to each fieldNames_ entry.
Definition: faOption.H:167
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition: faOption.H:164
virtual void addSup(const areaScalarField &h, faMatrix< scalar > &eqn, const label fieldi)
Definition: faOption.C:147
virtual void checkApplied() const
Check that the source has been applied.
Definition: faOption.C:132
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
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition: faOption.C:45
virtual void constrain(faMatrix< scalar > &eqn, const label fieldi)
Definition: faOption.C:241
const word name_
Source name.
Definition: faOption.H:146
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
rDeltaTY field()
dynamicFvMesh & mesh
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:478
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
messageStream Info
Information stream (stdout output on master, null elsewhere)
dimensionedScalar log(const dimensionedScalar &ds)
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:349
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:342
IOerror FatalIOError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:356
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
volScalarField & h
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333