ddtScheme.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) 2011-2018 OpenFOAM Foundation
9 Copyright (C) 2017-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "fv.H"
30#include "HashTable.H"
31#include "surfaceInterpolate.H"
32#include "fvMatrix.H"
33#include "cyclicAMIFvPatch.H"
34#include "registerSwitch.H"
35
36// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37
38namespace Foam
39{
40
41// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42
43namespace fv
44{
45
46// * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
47
48template<class Type>
50(
51 const fvMesh& mesh,
52 Istream& schemeData
53)
54{
55 if (fv::debug)
56 {
57 InfoInFunction << "Constructing ddtScheme<Type>" << endl;
58 }
59
60 if (schemeData.eof())
61 {
62 FatalIOErrorInFunction(schemeData)
63 << "Ddt scheme not specified" << endl << endl
64 << "Valid ddt schemes are :" << endl
65 << IstreamConstructorTablePtr_->sortedToc()
67 }
68
69 const word schemeName(schemeData);
70
71 auto* ctorPtr = IstreamConstructorTable(schemeName);
72
73 if (!ctorPtr)
74 {
76 (
77 schemeData,
78 "ddt",
79 schemeName,
80 *IstreamConstructorTablePtr_
81 ) << exit(FatalIOError);
82 }
83
84 return ctorPtr(mesh, schemeData);
85}
86
87
88// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
89
90template<class Type>
92(
93 const volScalarField& alpha,
94 const volScalarField& rho,
96)
97{
99
101 (
103 );
104}
105
106
107template<class Type>
109(
110 const volScalarField& alpha,
111 const volScalarField& rho,
113)
114{
116
118 (
119 vf,
120 alpha.dimensions()*rho.dimensions()
122 );
123}
124
125
126template<class Type>
128(
130)
131{
133
135 (
137 );
138}
139
140
141
142template<class Type>
144(
146 const fluxFieldType& phi,
147 const fluxFieldType& phiCorr
148)
149{
150 if (fv::debug)
151 {
152 InfoInFunction << "Using standard version" << endl;
153 }
154
155 tmp<surfaceScalarField> tddtCouplingCoeff
156 (
158 (
160 (
161 "ddtCouplingCoeff",
162 U.mesh().time().timeName(),
163 U.mesh()
164 ),
165 U.mesh(),
166 dimensionedScalar("one", dimless, 1.0)
167 )
168 );
169
170 surfaceScalarField& ddtCouplingCoeff = tddtCouplingCoeff.ref();
171
172 if (ddtPhiCoeff_ < 0)
173 {
174 // v1712 and earlier
175 ddtCouplingCoeff -= min
176 (
177 mag(phiCorr)
178 /(mag(phi) + dimensionedScalar("small", phi.dimensions(), SMALL)),
179 scalar(1)
180 );
181 }
182 else
183 {
184 ddtCouplingCoeff =
185 dimensionedScalar("ddtPhiCoeff", dimless, ddtPhiCoeff_);
186 }
187
188 surfaceScalarField::Boundary& ccbf = ddtCouplingCoeff.boundaryFieldRef();
189
190 forAll(U.boundaryField(), patchi)
191 {
192 if
193 (
194 U.boundaryField()[patchi].fixesValue()
195 || isA<cyclicAMIFvPatch>(mesh().boundary()[patchi])
196 )
197 {
198 ccbf[patchi] = 0.0;
199 }
200 }
201
202 if (debug > 1)
203 {
205 << "ddtCouplingCoeff mean max min = "
206 << gAverage(ddtCouplingCoeff.primitiveField())
207 << " " << gMax(ddtCouplingCoeff.primitiveField())
208 << " " << gMin(ddtCouplingCoeff.primitiveField())
209 << endl;
210 }
211
212 return tddtCouplingCoeff;
213}
214
215
216template<class Type>
218(
220 const fluxFieldType& phi,
221 const fluxFieldType& phiCorr
222)
223{
224 if (fv::debug)
225 {
226 InfoInFunction << "Using experimental version" << endl;
227 }
228
229 tmp<surfaceScalarField> tddtCouplingCoeff
230 (
232 (
234 (
235 "ddtCouplingCoeff",
236 U.mesh().time().timeName(),
237 U.mesh()
238 ),
239 U.mesh(),
240 dimensionedScalar("one", dimless, 1.0)
241 )
242 );
243
244 surfaceScalarField& ddtCouplingCoeff = tddtCouplingCoeff.ref();
245
246 if (ddtPhiCoeff_ < 0)
247 {
248 // See note below re: commented code
249 ddtCouplingCoeff -= min
250 (
251 // mag(phiCorr)
252 // *mesh().time().deltaT()*mag(mesh().deltaCoeffs())/mesh().magSf(),
253 // scalar(1)
254 mag(phiCorr)
255 *mesh().time().deltaT()*mesh().deltaCoeffs()/mesh().magSf(),
256 scalar(1)
257 );
258
259 // Note: setting oriented to false to avoid having to use mag(deltaCoeffs)
260 // - the deltaCoeffs field is always positive (scalars)
261 ddtCouplingCoeff.setOriented(false);
262 }
263 else
264 {
265 ddtCouplingCoeff =
266 dimensionedScalar("ddtPhiCoeff", dimless, ddtPhiCoeff_);
267 }
268
269 surfaceScalarField::Boundary& ccbf = ddtCouplingCoeff.boundaryFieldRef();
270
271 forAll(U.boundaryField(), patchi)
272 {
273 if
274 (
275 U.boundaryField()[patchi].fixesValue()
276 || isA<cyclicAMIFvPatch>(mesh().boundary()[patchi])
277 )
278 {
279 ccbf[patchi] = 0.0;
280 }
281 }
282
283 if (debug > 1)
284 {
286 << "ddtCouplingCoeff mean max min = "
287 << gAverage(ddtCouplingCoeff.primitiveField())
288 << " " << gMax(ddtCouplingCoeff.primitiveField())
289 << " " << gMin(ddtCouplingCoeff.primitiveField())
290 << endl;
291 }
292
293 return tddtCouplingCoeff;
294}
295
296
297template<class Type>
299(
301 const fluxFieldType& phi,
302 const fluxFieldType& phiCorr,
303 const volScalarField& rho
304)
305{
306 if (experimentalDdtCorr)
307 {
308 return
309 fvcDdtPhiCoeffExperimental
310 (
311 U,
312 phi,
313 phiCorr/fvc::interpolate(rho)
314 );
315 }
316 else
317 {
318 return fvcDdtPhiCoeff(U, phi, phiCorr);
319 }
320}
321
322
323template<class Type>
325(
327 const fluxFieldType& phi
328)
329{
330 if (experimentalDdtCorr)
331 {
332 return
333 fvcDdtPhiCoeffExperimental
334 (
335 U,
336 phi,
337 phi - fvc::dotInterpolate(mesh().Sf(), U)
338 );
339 }
340 else
341 {
342 return
343 fvcDdtPhiCoeff
344 (
345 U,
346 phi,
347 phi - fvc::dotInterpolate(mesh().Sf(), U)
348 );
349 }
350}
351
352
353template<class Type>
355(
357 const fluxFieldType& phi,
358 const volScalarField& rho
359)
360{
361 if (experimentalDdtCorr)
362 {
363 return fvcDdtPhiCoeffExperimental
364 (
365 rhoU,
366 phi,
367 (phi - fvc::dotInterpolate(mesh().Sf(), rhoU))
369 );
370 }
371 else
372 {
373 return fvcDdtPhiCoeff
374 (
375 rhoU,
376 phi,
377 (phi - fvc::dotInterpolate(mesh().Sf(), rhoU))
378 );
379 }
380}
381
382
383// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
384
385} // End namespace fv
386
387// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
388
389} // End namespace Foam
390
391// ************************************************************************* //
surfaceScalarField & phi
const dimensionSet & dimensions() const
Return dimensions.
void setOriented(const bool oriented=true) noexcept
Set the oriented flag.
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
static const GeometricField< Type, PatchField, GeoMesh > & null()
Return a null geometric field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
bool eof() const noexcept
True if end of input seen.
Definition: IOstream.H:239
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
virtual tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)=0
tmp< surfaceScalarField > fvcDdtPhiCoeffExperimental(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi, const fluxFieldType &phiCorr)
Definition: ddtScheme.C:218
virtual tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const dimensioned< Type > &)=0
tmp< surfaceScalarField > fvcDdtPhiCoeff(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi, const fluxFieldType &phiCorr)
Definition: ddtScheme.C:144
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
faceListList boundary
dynamicFvMesh & mesh
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:478
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:517
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
#define InfoInFunction
Report an information message using Foam::Info.
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
static tmp< GeometricField< typename innerProduct< vector, Type >::type, fvsPatchField, surfaceMesh > > dotInterpolate(const surfaceVectorField &Sf, const GeometricField< Type, fvPatchField, volMesh > &tvf)
Interpolate field onto faces.
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimless
Dimensionless.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
IOerror FatalIOError
Type gAverage(const FieldField< Field, Type > &f)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Type gMin(const FieldField< Field, Type > &f)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Type gMax(const FieldField< Field, Type > &f)
const dimensionSet dimVol(dimVolume)
Older spelling for dimVolume.
Definition: dimensionSets.H:64
labelList fv(nPoints)
volScalarField & alpha
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333