ATCModel.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) 2007-2021 PCOpt/NTUA
9 Copyright (C) 2013-2021 FOSS GP
10 Copyright (C) 2019-2021 OpenCFD Ltd.
11-------------------------------------------------------------------------------
12License
13 This file is part of OpenFOAM.
14
15 OpenFOAM is free software: you can redistribute it and/or modify it
16 under the terms of the GNU General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23 for more details.
24
25 You should have received a copy of the GNU General Public License
26 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27
28\*---------------------------------------------------------------------------*/
29
30#include "ATCModel.H"
31#include "localMin.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35namespace Foam
36{
37
38// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39
42
43// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
44
46{
48}
49
50
52{
54 DebugInfo<<
55 "max ATC mag " << gMax(ATC_) << endl;
56}
57
58
59// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
60
62(
63 const fvMesh& mesh,
64 const incompressibleVars& primalVars,
65 const incompressibleAdjointVars& adjointVars,
66 const dictionary& dict
67)
68:
70 (
72 (
73 "ATCModel" + adjointVars.solverName(),
74 mesh.time().timeName(),
75 mesh,
76 IOobject::NO_READ,
77 IOobject::NO_WRITE
78 )
79 ),
80 mesh_(mesh),
81 primalVars_(primalVars),
82 adjointVars_(adjointVars),
83 dict_(dict),
84 extraConvection_(dict_.getOrDefault<scalar>("extraConvection", Zero)),
85 extraDiffusion_(dict_.getOrDefault<scalar>("extraDiffusion", Zero)),
86 nSmooth_(dict_.getOrDefault<label>("nSmooth", 0)),
87 reconstructGradients_
88 (
89 dict_.getOrDefault("reconstructGradients", false)
90 ),
91 adjointSolverName_(adjointVars.solverName()),
92 zeroATCcells_(zeroATCcells::New(mesh, dict_)),
93 ATClimiter_
94 (
96 (
97 "ATClimiter" + adjointSolverName_,
98 mesh_.time().timeName(),
99 mesh_,
100 IOobject::NO_READ,
101 IOobject::NO_WRITE
102 ),
103 mesh_,
104 dimensionedScalar("limiter", dimless, 1.0),
105 zeroGradientFvPatchField<scalar>::typeName
106 ),
107 ATC_
108 (
110 (
111 "ATCField" + adjointSolverName_,
112 mesh_.time().timeName(),
113 mesh_,
114 IOobject::NO_READ,
115 IOobject::NO_WRITE
116 ),
117 mesh_,
118 dimensionedVector(dimensionSet(0, 1, -2, 0, 0), Zero)
119 )
120{
121 // Compute ATC limiter
123}
124
125
126// * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
127
129(
130 const fvMesh& mesh,
131 const incompressibleVars& primalVars,
132 const incompressibleAdjointVars& adjointVars,
133 const dictionary& dict
134)
135{
136 const word modelType(dict.get<word>("ATCModel"));
137
138 auto* ctorPtr = dictionaryConstructorTable(modelType);
139
140 Info<< "ATCModel type " << modelType << endl;
141
142 if (!ctorPtr)
143 {
145 (
146 dict,
147 "ATCModel",
148 modelType,
149 *dictionaryConstructorTablePtr_
150 ) << exit(FatalIOError);
151 }
152
153 return autoPtr<ATCModel>
154 (
155 ctorPtr(mesh, primalVars, adjointVars, dict)
156 );
157}
158
159
160// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
161
163{
164 return zeroATCcells_->getZeroATCcells();
165}
166
167
169{
170 return extraConvection_;
171}
172
173
175{
176 return extraDiffusion_;
177}
178
179
181{
182 return ATClimiter_;
183}
184
185
187(
189 const labelList& cells,
190 const label nSmooth
191)
192{
193 // Restore values
194 limiter.primitiveFieldRef() = 1;
195 limiter.correctBoundaryConditions();
196
197 // Set to zero in predefined cells
198 for (const label celli : cells)
199 {
200 limiter[celli] = Zero;
201 }
202
203 // Correct bcs to get the correct value for boundary faces
204 limiter.correctBoundaryConditions();
205
206 // Apply "laplacian" smoother
207 const fvMesh& mesh = limiter.mesh();
208 const localMin<scalar> scheme(mesh);
209 for (label iLimit = 0; iLimit < nSmooth; ++iLimit)
210 {
211 surfaceScalarField faceLimiter
212 (
213 scheme.interpolate(limiter)
214 );
215 limiter = fvc::average(faceLimiter);
216 }
217}
218
219
221(
222 const fvMesh& mesh,
223 const dictionary& dict
224)
225{
227 const labelList& zeroCells = zeroType->getZeroATCcells();
228 const label nSmooth = dict.getOrDefault<label>("nSmooth", 0);
229
230 tmp<volScalarField> tlimiter
231 (
233 (
235 (
236 "limiter",
237 mesh.time().timeName(),
238 mesh,
241 ),
242 mesh,
243 dimensionedScalar("limiter", dimless, 1.0),
245 )
246 );
247 volScalarField& limiter = tlimiter.ref();
248
250
251 return tlimiter;
252}
253
254
255bool ATCModel::writeData(Ostream&) const
256{
257 // Does nothing
258 return true;
259}
260
261
263{
264 // Does nothing in base
265}
266
267
268// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
269
270} // End namespace Foam
271
272// ************************************************************************* //
Base class for selecting the adjoint transpose convection model. Inherits from regIOobject to add loo...
Definition: ATCModel.H:63
static tmp< volScalarField > createLimiter(const fvMesh &mesh, const dictionary &dict)
Definition: ATCModel.C:221
const labelList & getZeroATCcells()
Get the list of cells on which to zero ATC.
Definition: ATCModel.C:162
const scalar extraDiffusion_
Definition: ATCModel.H:86
scalar getExtraDiffusionMultiplier()
Get the extra diffusion multiplier.
Definition: ATCModel.C:174
scalar getExtraConvectionMultiplier()
Get the extra convection multiplier.
Definition: ATCModel.C:168
const scalar extraConvection_
Definition: ATCModel.H:85
volVectorField ATC_
Definition: ATCModel.H:92
virtual void updatePrimalBasedQuantities()
Update quantities related with the primal fields.
Definition: ATCModel.C:262
const volScalarField & getLimiter() const
Get the list of cells on which to zero ATC.
Definition: ATCModel.C:180
void smoothATC()
Limit ATC field using ATClimiter_.
Definition: ATCModel.C:51
autoPtr< zeroATCcells > zeroATCcells_
Definition: ATCModel.H:90
volScalarField ATClimiter_
Definition: ATCModel.H:91
void computeLimiter()
Compute limiter based on the cells given by zeroATCcells.
Definition: ATCModel.C:45
const label nSmooth_
Definition: ATCModel.H:87
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
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
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:109
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
Class including all adjoint fields for incompressible flows.
Base class for solution control classes.
LocalMin-mean differencing scheme class.
Definition: localMin.H:62
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the face-interpolate of the given cell field.
Definition: localMin.H:123
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:76
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
Base class for selecting cells on which to zero the ATC term.
Definition: zeroATCcells.H:57
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
dynamicFvMesh & mesh
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:478
const cellShapeList & cells
word timeName
Definition: getTimeIndex.H:3
#define DebugInfo
Report an information message using Foam::Info.
tmp< GeometricField< Type, fvPatchField, volMesh > > average(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Area-weighted average a surfaceField creating a volField.
Definition: fvcAverage.C:48
Namespace for OpenFOAM.
const dimensionSet dimless
Dimensionless.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
IOerror FatalIOError
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Type gMax(const FieldField< Field, Type > &f)
tmp< areaScalarField > limiter(const areaScalarField &phi)
Definition: faNVDscheme.C:38
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
zeroCells(alpha, inletCells)