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 -------------------------------------------------------------------------------
12 License
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 
35 namespace Foam
36 {
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
42 
43 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
44 
46 {
47  computeLimiter(ATClimiter_, zeroATCcells_->getZeroATCcells(), nSmooth_);
48 }
49 
50 
52 {
53  ATC_ *= ATClimiter_;
54  DebugInfo<<
55  "max ATC mag " << gMax(ATC_) << endl;
56 }
57 
58 
59 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
60 
61 ATCModel::ATCModel
62 (
63  const fvMesh& mesh,
64  const incompressibleVars& primalVars,
65  const incompressibleAdjointVars& adjointVars,
66  const dictionary& dict
67 )
68 :
70  (
71  IOobject
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  (
95  IOobject
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,
107  ),
108  ATC_
109  (
110  IOobject
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 
170 {
171  return extraConvection_;
172 }
173 
176 {
177  return extraDiffusion_;
178 }
179 
180 
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 
222 tmp<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  tmp<volScalarField> tlimiter
233  (
234  new volScalarField
235  (
236  IOobject
237  (
238  "limiter",
239  mesh.time().timeName(),
240  mesh,
243  ),
244  mesh,
245  scalar(1),
246  dimless,
248  )
249  );
250  volScalarField& limiter = tlimiter.ref();
252  computeLimiter(limiter, zeroCells, nSmooth);
253 
254  return tlimiter;
255 }
256 
257 
259 {
260  // Does nothing
261  return true;
262 }
263 
264 
266 {
267  // Does nothing in base
268 }
269 
270 
271 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
272 
273 } // End namespace Foam
274 
275 // ************************************************************************* //
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: fvPatchField.H:221
void computeLimiter()
Compute limiter based on the cells given by zeroATCcells.
Definition: ATCModel.C:38
virtual void updatePrimalBasedQuantities()
Update quantities related with the primal fields.
Definition: ATCModel.C:258
dictionary dict
volScalarField ATClimiter_
Definition: ATCModel.H:91
const word zeroGradientType
A zeroGradient patch field type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
void smoothATC()
Limit ATC field using ATClimiter_.
Definition: ATCModel.C:44
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
autoPtr< zeroATCcells > zeroATCcells_
Definition: ATCModel.H:90
static tmp< edgeInterpolationScheme< Type > > scheme(const edgeScalarField &faceFlux, Istream &schemeData)
Return weighting factors for scheme given from Istream.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
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.
Template invariant parts for fvPatchField.
Definition: fvPatchField.H:77
zeroCells(alpha, inletCells)
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
tmp< areaScalarField > limiter(const areaScalarField &phi)
Definition: faNVDscheme.C:31
scalar getExtraDiffusionMultiplier()
Get the extra diffusion multiplier.
Definition: ATCModel.C:168
Class including all adjoint fields for incompressible flows.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
static autoPtr< zeroATCcells > New(const fvMesh &mesh, const dictionary &dict)
Return a reference to the selected turbulence model.
Definition: zeroATCcells.C:79
word timeName
Definition: getTimeIndex.H:3
Base class for solution control classes.
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:105
dynamicFvMesh & mesh
const cellShapeList & cells
const label nSmooth_
Definition: ATCModel.H:87
virtual bool writeData(Ostream &) const
Dummy writeData function required from regIOobject.
Definition: ATCModel.C:251
A class for handling words, derived from Foam::string.
Definition: word.H:63
Base class for selecting the adjoint transpose convection model. Inherits from regIOobject to add loo...
Definition: ATCModel.H:56
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
volVectorField ATC_
Definition: ATCModel.H:92
#define DebugInfo
Report an information message using Foam::Info.
const scalar extraConvection_
Definition: ATCModel.H:85
Base class for selecting cells on which to zero the ATC term.
Definition: zeroATCcells.H:52
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
scalar getExtraConvectionMultiplier()
Get the extra convection multiplier.
Definition: ATCModel.C:162
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
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
Type gMax(const FieldField< Field, Type > &f)
defineTypeNameAndDebug(combustionModel, 0)
const labelList & getZeroATCcells()
Get the list of cells on which to zero ATC.
Definition: ATCModel.C:156
const scalar extraDiffusion_
Definition: ATCModel.H:86
const volScalarField & getLimiter() const
Get the list of cells on which to zero ATC.
Definition: ATCModel.C:174
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Nothing to be read.
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:66
messageStream Info
Information stream (stdout output on master, null elsewhere)
tmp< GeometricField< Type, fvPatchField, volMesh > > average(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Area-weighted average a surfaceField creating a volField.
Definition: fvcAverage.C:41
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
List< label > labelList
A List of labels.
Definition: List.H:62
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:635
Namespace for OpenFOAM.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127