Loading...
Searching...
No Matches
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// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
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)
99 return nullptr;
100}
101
102
103template<class Type>
105(
106 const volScalarField& alpha,
107 const volScalarField& rho,
109)
112 return nullptr;
113}
114
115
116template<class Type>
118(
120)
123 return nullptr;
124}
125
126
127template<class Type>
129(
131 const fluxFieldType& phi,
132 const fluxFieldType& phiCorr
133)
134{
135 if (fv::debug)
136 {
137 InfoInFunction << "Using standard version" << endl;
138 }
139
140 tmp<surfaceScalarField> tddtCouplingCoeff
141 (
143 (
145 (
146 "ddtCouplingCoeff",
147 U.mesh().time().timeName(),
148 U.mesh().thisDb()
149 ),
150 U.mesh(),
151 dimensionedScalar("one", dimless, 1.0)
152 )
153 );
154
155 surfaceScalarField& ddtCouplingCoeff = tddtCouplingCoeff.ref();
156
157 if (ddtPhiCoeff_ < 0)
158 {
159 // v1712 and earlier
160 ddtCouplingCoeff -= min
161 (
162 mag(phiCorr)
163 /(mag(phi) + dimensionedScalar("small", phi.dimensions(), SMALL)),
164 scalar(1)
165 );
166 }
167 else
168 {
169 ddtCouplingCoeff =
170 dimensionedScalar("ddtPhiCoeff", dimless, ddtPhiCoeff_);
171 }
172
173 surfaceScalarField::Boundary& ccbf = ddtCouplingCoeff.boundaryFieldRef();
174
175 forAll(U.boundaryField(), patchi)
176 {
177 if
178 (
179 U.boundaryField()[patchi].fixesValue()
180 || isA<cyclicAMIFvPatch>(mesh().boundary()[patchi])
181 )
182 {
183 ccbf[patchi] = 0.0;
184 }
185 }
186
187 if (debug > 1)
188 {
189 auto limits = gMinMax(ddtCouplingCoeff.primitiveField());
190 auto avg = gAverage(ddtCouplingCoeff.primitiveField());
191
193 << "ddtCouplingCoeff mean max min = "
194 << avg << ' ' << limits.max() << ' ' << limits.min() << endl;
196
197 return tddtCouplingCoeff;
198}
199
200
201template<class Type>
203(
205 const fluxFieldType& phi,
206 const fluxFieldType& phiCorr
207)
208{
209 if (fv::debug)
210 {
211 InfoInFunction << "Using experimental version" << endl;
212 }
213
214 tmp<surfaceScalarField> tddtCouplingCoeff
215 (
217 (
218 IOobject
219 (
220 "ddtCouplingCoeff",
221 U.mesh().time().timeName(),
222 U.mesh().thisDb()
223 ),
224 U.mesh(),
225 dimensionedScalar("one", dimless, 1.0)
226 )
227 );
228
229 surfaceScalarField& ddtCouplingCoeff = tddtCouplingCoeff.ref();
230
231 if (ddtPhiCoeff_ < 0)
232 {
233 // See note below re: commented code
234 ddtCouplingCoeff -= min
235 (
236 // mag(phiCorr)
237 // *mesh().time().deltaT()*mag(mesh().deltaCoeffs())/mesh().magSf(),
238 // scalar(1)
239 mag(phiCorr)
240 *mesh().time().deltaT()*mesh().deltaCoeffs()/mesh().magSf(),
241 scalar(1)
242 );
243
244 // Note: setting oriented to false to avoid having to use mag(deltaCoeffs)
245 // - the deltaCoeffs field is always positive (scalars)
246 ddtCouplingCoeff.setOriented(false);
247 }
248 else
249 {
250 ddtCouplingCoeff =
251 dimensionedScalar("ddtPhiCoeff", dimless, ddtPhiCoeff_);
252 }
253
254 surfaceScalarField::Boundary& ccbf = ddtCouplingCoeff.boundaryFieldRef();
255
256 forAll(U.boundaryField(), patchi)
257 {
258 if
259 (
260 U.boundaryField()[patchi].fixesValue()
261 || isA<cyclicAMIFvPatch>(mesh().boundary()[patchi])
262 )
263 {
264 ccbf[patchi] = 0.0;
265 }
266 }
267
268 if (debug > 1)
269 {
270 auto limits = gMinMax(ddtCouplingCoeff.primitiveField());
271 auto avg = gAverage(ddtCouplingCoeff.primitiveField());
272
274 << "ddtCouplingCoeff mean max min = "
275 << avg << ' ' << limits.max() << ' ' << limits.min() << endl;
277
278 return tddtCouplingCoeff;
279}
280
281
282template<class Type>
284(
286 const fluxFieldType& phi,
287 const fluxFieldType& phiCorr,
288 const volScalarField& rho
289)
290{
291 if (experimentalDdtCorr)
292 {
293 return
294 fvcDdtPhiCoeffExperimental
295 (
296 U,
297 phi,
298 phiCorr/fvc::interpolate(rho)
299 );
300 }
301 else
303 return fvcDdtPhiCoeff(U, phi, phiCorr);
304 }
305}
306
307
308template<class Type>
310(
312 const fluxFieldType& phi
313)
314{
315 if (experimentalDdtCorr)
316 {
317 return
318 fvcDdtPhiCoeffExperimental
319 (
320 U,
321 phi,
322 phi - fvc::dotInterpolate(mesh().Sf(), U)
323 );
324 }
325 else
326 {
327 return
328 fvcDdtPhiCoeff
329 (
330 U,
331 phi,
333 );
334 }
335}
336
337
338template<class Type>
340(
342 const fluxFieldType& phi,
343 const volScalarField& rho
344)
345{
346 if (experimentalDdtCorr)
347 {
348 return fvcDdtPhiCoeffExperimental
349 (
350 rhoU,
351 phi,
352 (phi - fvc::dotInterpolate(mesh().Sf(), rhoU))
354 );
355 }
356 else
357 {
358 return fvcDdtPhiCoeff
359 (
360 rhoU,
361 phi,
362 (phi - fvc::dotInterpolate(mesh().Sf(), rhoU))
363 );
364 }
365}
366
367
368// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
369
370} // End namespace fv
371
372// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
373
374} // End namespace Foam
375
376// ************************************************************************* //
void setOriented(bool on=true) noexcept
Set the oriented flag: on/off.
Generic GeometricField class.
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef().
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
GeometricBoundaryField< scalar, fvsPatchField, surfaceMesh > Boundary
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
bool eof() const noexcept
True if end of input seen.
Definition IOstream.H:289
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
static bool experimentalDdtCorr
Flag to use experimental ddtCorr from org version Default is off for backwards compatibility.
Definition ddtScheme.H:73
scalar ddtPhiCoeff_
Input for fvcDdtPhiCoeff.
Definition ddtScheme.H:99
static tmp< ddtScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return a pointer to a new ddtScheme created on freestore.
Definition ddtScheme.C:43
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:196
const fvMesh & mesh() const
Return mesh reference.
Definition ddtScheme.H:179
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:122
GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMesh > fluxFieldType
Definition ddtScheme.H:247
A class for managing temporary objects.
Definition tmp.H:75
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
auto limits
Definition setRDeltaT.H:186
faceListList boundary
dynamicFvMesh & mesh
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition error.H:637
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition error.H:688
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
#define InfoInFunction
Report an information message using Foam::Info.
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for finite-volume.
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.
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
const dimensionSet dimless
Dimensionless.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
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:26
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
volScalarField & alpha
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299