CrankNicolsonDdtScheme.H
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-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Class
27  Foam::fv::CrankNicolsonDdtScheme
28 
29 Group
30  grpFvDdtSchemes
31 
32 Description
33  Second-oder Crank-Nicolson implicit ddt using the current and
34  previous time-step fields as well as the previous time-step ddt.
35 
36  The Crank-Nicolson scheme is often unstable for complex flows in complex
37  geometries and it is necessary to "off-centre" the scheme to stabilize it
38  while retaining greater temporal accuracy than the first-order
39  Euler-implicit scheme. Off-centering is specified via the mandatory
40  coefficient \c ocCoeff in the range [0,1] following the scheme name e.g.
41  \verbatim
42  ddtSchemes
43  {
44  default CrankNicolson 0.9;
45  }
46  \endverbatim
47  or with an optional "ramp" function to transition from the Euler scheme to
48  Crank-Nicolson over a initial period to avoid start-up problems, e.g.
49  \verbatim
50  ddtSchemes
51  {
52  default CrankNicolson
53  ocCoeff
54  {
55  type scale;
56  scale linearRamp;
57  duration 0.01;
58  value 0.9;
59  };
60  }
61  \endverbatim
62  With a coefficient of 1 the scheme is fully centred and second-order,
63  with a coefficient of 0 the scheme is equivalent to Euler-implicit.
64  A coefficient of 0.9 has been found to be suitable for a range of cases for
65  which higher-order accuracy in time is needed and provides similar accuracy
66  and stability to the backward scheme. However, the backward scheme has
67  been found to be more robust and provides formal second-order accuracy in
68  time.
69 
70  The advantage of the Crank-Nicolson scheme over backward is that only the
71  new and old-time values are needed, the additional terms relating to the
72  fluxes and sources are evaluated at the mid-point of the time-step which
73  provides the opportunity to limit the fluxes in such a way as to ensure
74  boundedness while maintaining greater accuracy in time compared to the
75  Euler-implicit scheme. This approach is now used with MULES in the
76  interFoam family of solvers. Boundedness cannot be guaranteed with the
77  backward scheme.
78 
79 Note
80  The Crank-Nicolson coefficient for the implicit part of the RHS is related
81  to the off-centering coefficient by
82  \verbatim
83  cnCoeff = 1.0/(1.0 + ocCoeff);
84  \endverbatim
85 
86 See also
87  Foam::Euler
88  Foam::backward
89 
90 SourceFiles
91  CrankNicolsonDdtScheme.C
92 
93 \*---------------------------------------------------------------------------*/
94 
95 #ifndef CrankNicolsonDdtScheme_H
96 #define CrankNicolsonDdtScheme_H
97 
98 #include "ddtScheme.H"
99 #include "Function1.H"
100 
101 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
102 
103 namespace Foam
104 {
105 
106 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
107 
108 namespace fv
109 {
110 
111 /*---------------------------------------------------------------------------*\
112  Class CrankNicolsonDdtScheme Declaration
113 \*---------------------------------------------------------------------------*/
114 
115 template<class Type>
117 :
118  public fv::ddtScheme<Type>
119 {
120  // Private Data
121 
122  //- Class to store the ddt0 fields on the objectRegistry for use in the
123  // next time-step. The start-time index of the CN scheme is also
124  // stored to help handle the transition from Euler to CN
125  template<class GeoField>
126  class DDt0Field
127  :
128  public GeoField
129  {
130  label startTimeIndex_;
131 
132  public:
133 
134  //- Constructor from file for restart.
135  DDt0Field
136  (
137  const IOobject& io,
138  const fvMesh& mesh
139  );
140 
141  //- Constructor from components, initisalised to zero with given
142  // dimensions.
143  DDt0Field
144  (
145  const IOobject& io,
146  const fvMesh& mesh,
148  );
149 
150  //- Return the start-time index
151  label startTimeIndex() const;
152 
153  //- Cast to the underlying GeoField
154  GeoField& operator()();
155 
156  //- Assignment to a GeoField
157  void operator=(const GeoField& gf);
158  };
159 
160 
161  //- Off-centering coefficient function
162  // 1 -> CN, less than one blends with EI
163  autoPtr<Function1<scalar>> ocCoeff_;
164 
165 
166  // Private Member Functions
167 
168  //- No copy construct
170 
171  //- No copy assignment
172  void operator=(const CrankNicolsonDdtScheme&) = delete;
173 
174  template<class GeoField>
175  DDt0Field<GeoField>& ddt0_
176  (
177  const word& name,
178  const dimensionSet& dims
179  );
180 
181  //- Check if the ddt0 needs to be evaluated for this time-step
182  template<class GeoField>
183  bool evaluate(const DDt0Field<GeoField>& ddt0) const;
184 
185  //- Return the coefficient for Euler scheme for the first time-step
186  // for and CN thereafter
187  template<class GeoField>
188  scalar coef_(const DDt0Field<GeoField>&) const;
189 
190  //- Return the old time-step coefficient for Euler scheme for the
191  // second time-step and for CN thereafter
192  template<class GeoField>
193  scalar coef0_(const DDt0Field<GeoField>&) const;
194 
195  //- Return the reciprocal time-step coefficient for Euler for the
196  // first time-step and CN thereafter
197  template<class GeoField>
198  dimensionedScalar rDtCoef_(const DDt0Field<GeoField>&) const;
199 
200  //- Return the reciprocal old time-step coefficient for Euler for the
201  // second time-step and CN thereafter
202  template<class GeoField>
203  dimensionedScalar rDtCoef0_(const DDt0Field<GeoField>&) const;
204 
205  //- Return ddt0 multiplied by the off-centreing coefficient
206  template<class GeoField>
207  tmp<GeoField> offCentre_(const GeoField& ddt0) const;
208 
209 
210 public:
211 
212  //- Runtime type information
213  TypeName("CrankNicolson");
214 
215 
216  // Constructors
217 
218  //- Construct from mesh
220 
221  //- Construct from mesh and Istream
223 
224 
225  // Member Functions
226 
227  //- Return mesh reference
228  const fvMesh& mesh() const
229  {
230  return fv::ddtScheme<Type>::mesh();
231  }
232 
233  //- Return the current off-centreing coefficient
234  scalar ocCoeff() const
235  {
236  return ocCoeff_->value(mesh().time().value());
237  }
238 
240  (
241  const dimensioned<Type>&
242  );
243 
245  (
247  );
248 
250  (
251  const dimensionedScalar&,
253  );
254 
256  (
257  const volScalarField&,
259  );
260 
262  (
263  const volScalarField& alpha,
264  const volScalarField& rho,
266  );
267 
269  (
271  );
272 
274  (
275  const dimensionedScalar&,
277  );
278 
280  (
281  const volScalarField&,
283  );
284 
286  (
287  const volScalarField& alpha,
288  const volScalarField& rho,
290  );
291 
293 
295  (
298  );
299 
301  (
303  const fluxFieldType& phi
304  );
305 
307  (
308  const volScalarField& rho,
311  );
312 
314  (
315  const volScalarField& rho,
317  const fluxFieldType& phi
318  );
319 
320 
322  (
324  );
325 };
326 
327 
328 template<>
330 (
333 );
334 
335 template<>
337 (
338  const volScalarField& U,
339  const surfaceScalarField& phi
340 );
341 
342 template<>
344 (
345  const volScalarField& rho,
346  const volScalarField& U,
347  const surfaceScalarField& Uf
348 );
349 
350 template<>
352 (
353  const volScalarField& rho,
354  const volScalarField& U,
355  const surfaceScalarField& phi
356 );
357 
358 
359 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
360 
361 } // End namespace fv
362 
363 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
364 
365 } // End namespace Foam
366 
367 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
368 
369 #ifdef NoRepository
370  #include "CrankNicolsonDdtScheme.C"
371 #endif
372 
373 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
374 
375 #endif
376 
377 // ************************************************************************* //
Foam::fv::CrankNicolsonDdtScheme::mesh
const fvMesh & mesh() const
Return mesh reference.
Definition: CrankNicolsonDdtScheme.H:227
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Function1.H
Foam::fv::CrankNicolsonDdtScheme::fvcDdtUfCorr
tmp< fluxFieldType > fvcDdtUfCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
Definition: CrankNicolsonDdtScheme.C:1174
Uf
autoPtr< surfaceVectorField > Uf
Definition: createUfIfPresent.H:33
Foam::fv::CrankNicolsonDdtScheme::TypeName
TypeName("CrankNicolson")
Runtime type information.
CrankNicolsonDdtScheme.C
Foam::dimensionSet
Dimension set for the base types.
Definition: dimensionSet.H:65
Foam::fv::CrankNicolsonDdtScheme::fluxFieldType
ddtScheme< Type >::fluxFieldType fluxFieldType
Definition: CrankNicolsonDdtScheme.H:291
rho
rho
Definition: readInitialConditions.H:88
Foam::fv::CrankNicolsonDdtScheme::ocCoeff
scalar ocCoeff() const
Return the current off-centreing coefficient.
Definition: CrankNicolsonDdtScheme.H:233
Foam::fv::CrankNicolsonDdtScheme
Second-oder Crank-Nicolson implicit ddt using the current and previous time-step fields as well as th...
Definition: CrankNicolsonDdtScheme.H:115
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::fv::CrankNicolsonDdtScheme::fvmDdt
tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)
Definition: CrankNicolsonDdtScheme.C:809
Foam::fv::CrankNicolsonDdtScheme::fvcDdt
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const dimensioned< Type > &)
Definition: CrankNicolsonDdtScheme.C:345
Foam::fv::CrankNicolsonDdtScheme::fvcDdtPhiCorr
tmp< fluxFieldType > fvcDdtPhiCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
Definition: CrankNicolsonDdtScheme.C:1235
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:43
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::fv::CrankNicolsonDdtScheme::meshPhi
tmp< surfaceScalarField > meshPhi(const GeometricField< Type, fvPatchField, volMesh > &)
Definition: CrankNicolsonDdtScheme.C:1600
Foam::fv::ddtScheme::mesh
const fvMesh & mesh() const
Return mesh reference.
Definition: ddtScheme.H:162
ddtScheme.H
fv
labelList fv(nPoints)
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
U
U
Definition: pEqn.H:72
Foam::fv::ddtScheme
Abstract base class for ddt schemes.
Definition: ddtScheme.H:83
Foam::GeometricField< Type, fvPatchField, volMesh >
psi
const volScalarField & psi
Definition: createFieldRefs.H:1