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-2019 PCOpt/NTUA
9  Copyright (C) 2013-2019 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 
40 defineTypeNameAndDebug(ATCModel, 0);
41 defineRunTimeSelectionTable(ATCModel, dictionary);
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,
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_,
102  ),
103  mesh_,
104  dimensionedScalar("limiter", dimless, 1.0),
106  ),
107  ATC_
108  (
109  IOobject
110  (
111  "ATCField" + adjointSolverName_,
112  mesh_.time().timeName(),
113  mesh_,
116  ),
117  mesh_,
118  dimensionedVector(dimensionSet(0, 1, -2, 0, 0), Zero)
119  )
120 {
121  // Compute ATC limiter
122  computeLimiter();
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();
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  (
232  new volScalarField
233  (
234  IOobject
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 
249  computeLimiter(limiter, zeroCells, nSmooth);
250 
251  return tlimiter;
252 }
253 
254 
256 {
257  // Does nothing
258  return true;
259 }
260 
261 
262 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
263 
264 } // End namespace Foam
265 
266 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::variablesSet::solverName
const word & solverName() const
Return solver name.
Definition: variablesSet.C:84
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
ATCModel.H
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::localMin
LocalMin-mean differencing scheme class.
Definition: localMin.H:59
Foam::defineRunTimeSelectionTable
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Foam::ATCModel::writeData
virtual bool writeData(Ostream &) const
Dummy writeData function required from regIOobject.
Definition: ATCModel.C:255
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::ATCModel::getExtraConvectionMultiplier
scalar getExtraConvectionMultiplier()
Get the extra convection multiplier.
Definition: ATCModel.C:168
Foam::ATCModel::computeLimiter
void computeLimiter()
Compute limiter based on the cells given by zeroATCcells.
Definition: ATCModel.C:45
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:107
Foam::dimensionSet
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:108
Foam::fvc::average
tmp< GeometricField< Type, fvPatchField, volMesh > > average(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Area-weighted average a surfaceField creating a volField.
Definition: fvcAverage.C:48
Foam::incompressibleAdjointVars
Class including all adjoint fields for incompressible flows.
Definition: incompressibleAdjointVars.H:52
Foam::ATCModel::extraDiffusion_
const scalar extraDiffusion_
Definition: ATCModel.H:86
Foam::zeroGradientFvPatchField
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
Definition: zeroGradientFvPatchField.H:64
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:50
zeroCells
zeroCells(alpha, inletCells)
FatalIOErrorInLookup
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:478
Foam::ATCModel::getZeroATCcells
const labelList & getZeroATCcells()
Get the list of cells on which to zero ATC.
Definition: ATCModel.C:162
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::ATCModel::getLimiter
const volScalarField & getLimiter() const
Get the list of cells on which to zero ATC.
Definition: ATCModel.C:180
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::ATCModel::zeroATCcells_
autoPtr< zeroATCcells > zeroATCcells_
Definition: ATCModel.H:90
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::ATCModel::getExtraDiffusionMultiplier
scalar getExtraDiffusionMultiplier()
Get the extra diffusion multiplier.
Definition: ATCModel.C:174
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
localMin.H
Foam::regIOobject
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:73
Foam::zeroATCcells::New
static autoPtr< zeroATCcells > New(const fvMesh &mesh, const dictionary &dict)
Return a reference to the selected turbulence model.
Definition: zeroATCcells.C:86
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:382
Foam::ATCModel::nSmooth_
const label nSmooth_
Definition: ATCModel.H:87
Foam::zeroATCcells::getZeroATCcells
const labelList & getZeroATCcells()
Get the zeroATCcells.
Definition: zeroATCcells.H:125
Foam::List< label >
Foam::ATCModel::ATClimiter_
volScalarField ATClimiter_
Definition: ATCModel.H:91
Foam::ATCModel::createLimiter
static tmp< volScalarField > createLimiter(const fvMesh &mesh, const dictionary &dict)
Definition: ATCModel.C:221
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:280
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::ATCModel::smoothATC
void smoothATC()
Limit ATC field using ATClimiter_.
Definition: ATCModel.C:51
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:148
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::ATCModel::ATC_
volVectorField ATC_
Definition: ATCModel.H:92
Foam::limiter
tmp< areaScalarField > limiter(const areaScalarField &phi)
Definition: faNVDscheme.C:37
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::fac::scheme
static tmp< edgeInterpolationScheme< Type > > scheme(const edgeScalarField &faceFlux, Istream &schemeData)
Return weighting factors for scheme given from Istream.
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592
Foam::ATCModel::extraConvection_
const scalar extraConvection_
Definition: ATCModel.H:85
Foam::incompressibleVars
Base class for solution control classes.
Definition: incompressibleVars.H:54
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189
Foam::ATCModel::New
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:129