Loading...
Searching...
No Matches
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-2023 PCOpt/NTUA
9 Copyright (C) 2013-2023 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// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35namespace Foam
36{
37
38// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39
40defineTypeNameAndDebug(ATCModel, 0);
42
43// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
46{
48}
49
50
52{
55 "max ATC mag " << gMax(ATC_) << endl;
56}
57
58
59// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
60
61ATCModel::ATCModel
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 scalar(1),
105 dimless,
106 fvPatchFieldBase::zeroGradientType()
107 ),
108 ATC_
109 (
111 (
112 "ATCField" + adjointSolverName_,
113 mesh_.time().timeName(),
114 mesh_,
115 IOobject::NO_READ,
116 IOobject::NO_WRITE
117 ),
118 mesh_,
119 dimensionedVector(dimensionSet(0, 1, -2, 0, 0), Zero)
120 )
121{
122 // Compute ATC limiter
124}
125
126
127// * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
128
130(
131 const fvMesh& mesh,
132 const incompressibleVars& primalVars,
133 const incompressibleAdjointVars& adjointVars,
134 const dictionary& dict
135)
136{
137 const word modelType(dict.get<word>("ATCModel"));
138
139 auto* ctorPtr = dictionaryConstructorTable(modelType);
140
141 Info<< "ATCModel type " << modelType << endl;
142
143 if (!ctorPtr)
144 {
146 (
147 dict,
148 "ATCModel",
149 modelType,
150 *dictionaryConstructorTablePtr_
151 ) << exit(FatalIOError);
152 }
153
154 return autoPtr<ATCModel>
155 (
156 ctorPtr(mesh, primalVars, adjointVars, dict)
157 );
158}
159
160
161// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
164{
165 return zeroATCcells_->getZeroATCcells();
166}
167
172}
173
178}
179
182{
183 return ATClimiter_;
184}
185
186
188(
190 const labelList& cells,
191 const label nSmooth
192)
193{
194 // Restore values
195 limiter.primitiveFieldRef() = 1;
196 limiter.correctBoundaryConditions();
197
198 // Set to zero in predefined cells
199 for (const label celli : cells)
200 {
201 limiter[celli] = Zero;
202 }
203
204 // Correct bcs to get the correct value for boundary faces
205 limiter.correctBoundaryConditions();
206
207 // Apply "laplacian" smoother
208 const fvMesh& mesh = limiter.mesh();
209 const localMin<scalar> scheme(mesh);
210 for (label iLimit = 0; iLimit < nSmooth; ++iLimit)
211 {
212 surfaceScalarField faceLimiter
213 (
214 scheme.interpolate(limiter)
215 );
216 limiter = fvc::average(faceLimiter);
217 limiter.correctBoundaryConditions();
218 }
219}
220
221
222tmp<volScalarField> ATCModel::createLimiter
223(
224 const fvMesh& mesh,
225 const dictionary& dict
226)
227{
228 autoPtr<zeroATCcells> zeroType(zeroATCcells::New(mesh, dict));
229 const labelList& zeroCells = zeroType->getZeroATCcells();
230 const label nSmooth = dict.getOrDefault<label>("nSmooth", 0);
231
232 auto tlimiter = volScalarField::New
233 (
234 "limiter",
236 mesh,
237 scalar(1),
238 dimless,
240 );
241 auto& limiter = tlimiter.ref();
244
245 return tlimiter;
246}
247
248
250{
251 // Does nothing
252 return true;
253}
254
255
257{
258 // Does nothing in base
259}
260
261
262// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
263
264} // End namespace Foam
265
266// ************************************************************************* //
Base class for selecting the adjoint transpose convection model. Inherits from regIOobject to add loo...
Definition ATCModel.H:59
static tmp< volScalarField > createLimiter(const fvMesh &mesh, const dictionary &dict)
Return a limiter based on given cells. For use with classes other than ATC which employ the same smoo...
Definition ATCModel.C:216
const labelList & getZeroATCcells()
Get the list of cells on which to zero ATC.
Definition ATCModel.C:156
const incompressibleAdjointVars & adjointVars_
Definition ATCModel.H:81
const scalar extraDiffusion_
Definition ATCModel.H:86
const fvMesh & mesh_
Definition ATCModel.H:79
virtual bool writeData(Ostream &) const
Dummy writeData function required from regIOobject.
Definition ATCModel.C:242
scalar getExtraDiffusionMultiplier()
Get the extra diffusion multiplier.
Definition ATCModel.C:168
scalar getExtraConvectionMultiplier()
Get the extra convection multiplier.
Definition ATCModel.C:162
const scalar extraConvection_
Definition ATCModel.H:85
volVectorField ATC_
Definition ATCModel.H:92
const dictionary & dict_
Definition ATCModel.H:84
word adjointSolverName_
Definition ATCModel.H:89
virtual void updatePrimalBasedQuantities()
Update quantities related with the primal fields.
Definition ATCModel.C:249
const bool reconstructGradients_
Definition ATCModel.H:88
const volScalarField & getLimiter() const
Get the list of cells on which to zero ATC.
Definition ATCModel.C:174
void smoothATC()
Limit ATC field using ATClimiter_.
Definition ATCModel.C:44
const incompressibleVars & primalVars_
Definition ATCModel.H:80
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:38
const label nSmooth_
Definition ATCModel.H:87
static autoPtr< ATCModel > New(const fvMesh &mesh, const incompressibleVars &primalVars, const incompressibleAdjointVars &adjointVars, const dictionary &dict)
Return a reference to the selected turbulence model.
Definition ATCModel.C:123
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
@ 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 Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Template invariant parts for fvPatchField.
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Class including all adjoint fields for incompressible flows.
Base class for solution control classes.
LocalMin-mean differencing scheme class.
Definition localMin.H:58
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:133
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition regIOobject.H:71
regIOobject(const IOobject &io, const bool isTimeObject=false)
Construct from IOobject. The optional flag adds special handling if the object is the top-level regIO...
Definition regIOobject.C:43
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
Base class for selecting cells on which to zero the ATC term.
static autoPtr< zeroATCcells > New(const fvMesh &mesh, const dictionary &dict)
Return a reference to the selected turbulence model.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition error.H:637
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:41
Namespace for OpenFOAM.
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.
List< label > labelList
A List of labels.
Definition List.H:62
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
Type gMax(const FieldField< Field, Type > &f)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
tmp< areaScalarField > limiter(const areaScalarField &phi)
Definition faNVDscheme.C:31
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
zeroCells(alpha, inletCells)