localEulerDdtScheme.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "localEulerDdtScheme.H"
27 #include "surfaceInterpolate.H"
28 #include "fvMatrices.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace fv
38 {
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 template<class Type>
43 const volScalarField& localEulerDdtScheme<Type>::localRDeltaT() const
44 {
46 }
47 
48 
49 template<class Type>
50 const surfaceScalarField& localEulerDdtScheme<Type>::localRDeltaTf() const
51 {
53 }
54 
55 
56 template<class Type>
57 tmp<GeometricField<Type, fvPatchField, volMesh>>
59 (
60  const dimensioned<Type>& dt
61 )
62 {
63  const word ddtName("ddt(" + dt.name() + ')');
64 
66  (
68  (
69  ddtName,
70  mesh(),
72  (
73  "0",
74  dt.dimensions()/dimTime,
75  Zero
76  ),
78  )
79  );
80 }
81 
82 
83 template<class Type>
86 (
88 )
89 {
90  const volScalarField& rDeltaT = localRDeltaT();
91 
92  const word ddtName("ddt(" + vf.name() + ')');
93 
95  (
97  (
98  ddtName,
99  rDeltaT*(vf - vf.oldTime())
100  )
101  );
102 }
103 
104 
105 template<class Type>
108 (
109  const dimensionedScalar& rho,
111 )
112 {
113  const volScalarField& rDeltaT = localRDeltaT();
114 
115  const word ddtName("ddt(" + rho.name() + ',' + vf.name() + ')');
116 
118  (
120  (
121  ddtName,
122  rDeltaT*rho*(vf - vf.oldTime())
123  )
124  );
125 }
126 
127 
128 template<class Type>
131 (
132  const volScalarField& rho,
134 )
135 {
136  const volScalarField& rDeltaT = localRDeltaT();
137 
138  const word ddtName("ddt(" + rho.name() + ',' + vf.name() + ')');
139 
141  (
143  (
144  ddtName,
145  rDeltaT*(rho*vf - rho.oldTime()*vf.oldTime())
146  )
147  );
148 }
149 
150 
151 template<class Type>
154 (
155  const volScalarField& alpha,
156  const volScalarField& rho,
158 )
159 {
160  const volScalarField& rDeltaT = localRDeltaT();
161 
162  const word ddtName("ddt("+alpha.name()+','+rho.name()+','+vf.name()+')');
163 
165  (
167  (
168  ddtName,
169  rDeltaT
170  *(
171  alpha*rho*vf
172  - alpha.oldTime()*rho.oldTime()*vf.oldTime()
173  )
174  )
175  );
176 }
177 
178 
179 template<class Type>
182 (
184 )
185 {
186  const surfaceScalarField& rDeltaT = localRDeltaTf();
187 
188  const word ddtName("ddt("+sf.name()+')');
189 
191  (
192  ddtName,
193  rDeltaT*(sf - sf.oldTime())
194  );
195 }
196 
197 
198 template<class Type>
201 (
203 )
204 {
205  tmp<fvMatrix<Type>> tfvm
206  (
207  new fvMatrix<Type>
208  (
209  vf,
211  )
212  );
213 
214  fvMatrix<Type>& fvm = tfvm.ref();
215 
216  const scalarField& rDeltaT = localRDeltaT();
217 
218  fvm.diag() = rDeltaT*mesh().Vsc();
219  fvm.source() = rDeltaT*vf.oldTime().primitiveField()*mesh().Vsc();
220 
221  return tfvm;
222 }
223 
224 
225 template<class Type>
228 (
229  const dimensionedScalar& rho,
231 )
232 {
233  tmp<fvMatrix<Type>> tfvm
234  (
235  new fvMatrix<Type>
236  (
237  vf,
239  )
240  );
241  fvMatrix<Type>& fvm = tfvm.ref();
242 
243  const scalarField& rDeltaT = localRDeltaT();
244 
245  fvm.diag() = rDeltaT*rho.value()*mesh().Vsc();
246 
247  fvm.source() =
248  rDeltaT*rho.value()*vf.oldTime().primitiveField()*mesh().Vsc();
249 
250  return tfvm;
251 }
252 
253 
254 template<class Type>
257 (
258  const volScalarField& rho,
260 )
261 {
262  tmp<fvMatrix<Type>> tfvm
263  (
264  new fvMatrix<Type>
265  (
266  vf,
268  )
269  );
270  fvMatrix<Type>& fvm = tfvm.ref();
271 
272  const scalarField& rDeltaT = localRDeltaT();
273 
274  fvm.diag() = rDeltaT*rho.primitiveField()*mesh().Vsc();
275 
276  fvm.source() = rDeltaT
277  *rho.oldTime().primitiveField()
278  *vf.oldTime().primitiveField()*mesh().Vsc();
279 
280  return tfvm;
281 }
282 
283 
284 template<class Type>
287 (
288  const volScalarField& alpha,
289  const volScalarField& rho,
291 )
292 {
293  tmp<fvMatrix<Type>> tfvm
294  (
295  new fvMatrix<Type>
296  (
297  vf,
298  alpha.dimensions()*rho.dimensions()*vf.dimensions()*dimVol/dimTime
299  )
300  );
301  fvMatrix<Type>& fvm = tfvm.ref();
302 
303  const scalarField& rDeltaT = localRDeltaT();
304 
305  fvm.diag() =
306  rDeltaT*alpha.primitiveField()*rho.primitiveField()*mesh().Vsc();
307 
308  fvm.source() = rDeltaT
309  *alpha.oldTime().primitiveField()
310  *rho.oldTime().primitiveField()
311  *vf.oldTime().primitiveField()*mesh().Vsc();
312 
313  return tfvm;
314 }
315 
316 
317 /*
318 // Courant number limited formulation
319 template<class Type>
320 tmp<surfaceScalarField> localEulerDdtScheme<Type>::fvcDdtPhiCoeff
321 (
322  const GeometricField<Type, fvPatchField, volMesh>& U,
323  const fluxFieldType& phi,
324  const fluxFieldType& phiCorr
325 )
326 {
327  // Courant number limited formulation
328  tmp<surfaceScalarField> tddtCouplingCoeff = scalar(1)
329  - min
330  (
331  mag(phiCorr)*mesh().deltaCoeffs()
332  /(fvc::interpolate(localRDeltaT())*mesh().magSf()),
333  scalar(1)
334  );
335 
336  surfaceScalarField& ddtCouplingCoeff = tddtCouplingCoeff.ref();
337 
338  surfaceScalarField::Boundary& ccbf = ddtCouplingCoeff.boundaryFieldRef();
339 
340  forAll(U.boundaryField(), patchi)
341  {
342  if
343  (
344  U.boundaryField()[patchi].fixesValue()
345  || isA<cyclicAMIFvPatch>(mesh().boundary()[patchi])
346  )
347  {
348  ccbf[patchi] = 0.0;
349  }
350  }
351 
352  if (debug > 1)
353  {
354  InfoInFunction
355  << "ddtCouplingCoeff mean max min = "
356  << gAverage(ddtCouplingCoeff.primitiveField())
357  << " " << gMax(ddtCouplingCoeff.primitiveField())
358  << " " << gMin(ddtCouplingCoeff.primitiveField())
359  << endl;
360  }
361 
362  return tddtCouplingCoeff;
363 }
364 */
365 
366 
367 template<class Type>
370 (
373 )
374 {
375  const surfaceScalarField rDeltaT(fvc::interpolate(localRDeltaT()));
376 
377  fluxFieldType phiUf0(mesh().Sf() & Uf.oldTime());
378  fluxFieldType phiCorr
379  (
380  phiUf0 - fvc::dotInterpolate(mesh().Sf(), U.oldTime())
381  );
382 
383  return tmp<fluxFieldType>
384  (
385  new fluxFieldType
386  (
387  IOobject
388  (
389  "ddtCorr(" + U.name() + ',' + Uf.name() + ')',
390  mesh().time().timeName(),
391  mesh()
392  ),
393  this->fvcDdtPhiCoeff(U.oldTime(), phiUf0, phiCorr)
394  *rDeltaT*phiCorr
395  )
396  );
397 }
398 
399 
400 template<class Type>
403 (
405  const fluxFieldType& phi
406 )
407 {
408  const surfaceScalarField rDeltaT(fvc::interpolate(localRDeltaT()));
409 
410  fluxFieldType phiCorr
411  (
412  phi.oldTime() - fvc::dotInterpolate(mesh().Sf(), U.oldTime())
413  );
414 
415  return tmp<fluxFieldType>
416  (
417  new fluxFieldType
418  (
419  IOobject
420  (
421  "ddtCorr(" + U.name() + ',' + phi.name() + ')',
422  mesh().time().timeName(),
423  mesh()
424  ),
425  this->fvcDdtPhiCoeff(U.oldTime(), phi.oldTime(), phiCorr)
426  *rDeltaT*phiCorr
427  )
428  );
429 }
430 
431 
432 template<class Type>
435 (
436  const volScalarField& rho,
439 )
440 {
441  const surfaceScalarField rDeltaT(fvc::interpolate(localRDeltaT()));
442 
443  if
444  (
445  U.dimensions() == dimVelocity
447  )
448  {
450  (
451  rho.oldTime()*U.oldTime()
452  );
453 
454  fluxFieldType phiUf0(mesh().Sf() & Uf.oldTime());
455  fluxFieldType phiCorr(phiUf0 - fvc::dotInterpolate(mesh().Sf(), rhoU0));
456 
457  return tmp<fluxFieldType>
458  (
459  new fluxFieldType
460  (
461  IOobject
462  (
463  "ddtCorr("
464  + rho.name() + ',' + U.name() + ',' + Uf.name() + ')',
465  mesh().time().timeName(),
466  mesh()
467  ),
468  this->fvcDdtPhiCoeff(rhoU0, phiUf0, phiCorr, rho.oldTime())
469  *rDeltaT*phiCorr
470  )
471  );
472  }
473  else if
474  (
477  )
478  {
479  fluxFieldType phiUf0(mesh().Sf() & Uf.oldTime());
480  fluxFieldType phiCorr
481  (
482  phiUf0 - fvc::dotInterpolate(mesh().Sf(), U.oldTime())
483  );
484 
485  return tmp<fluxFieldType>
486  (
487  new fluxFieldType
488  (
489  IOobject
490  (
491  "ddtCorr("
492  + rho.name() + ',' + U.name() + ',' + Uf.name() + ')',
493  mesh().time().timeName(),
494  mesh()
495  ),
496  this->fvcDdtPhiCoeff
497  (
498  U.oldTime(),
499  phiUf0,
500  phiCorr,
501  rho.oldTime()
502  )*rDeltaT*phiCorr
503  )
504  );
505  }
506  else
507  {
509  << "dimensions of Uf are not correct"
510  << abort(FatalError);
511 
512  return fluxFieldType::null();
513  }
514 }
515 
516 
517 template<class Type>
520 (
521  const volScalarField& rho,
523  const fluxFieldType& phi
524 )
525 {
526  const surfaceScalarField rDeltaT(fvc::interpolate(localRDeltaT()));
527 
528  if
529  (
530  U.dimensions() == dimVelocity
531  && phi.dimensions() == rho.dimensions()*dimVelocity*dimArea
532  )
533  {
535  (
536  rho.oldTime()*U.oldTime()
537  );
538 
539  fluxFieldType phiCorr
540  (
541  phi.oldTime() - fvc::dotInterpolate(mesh().Sf(), rhoU0)
542  );
543 
544  return tmp<fluxFieldType>
545  (
546  new fluxFieldType
547  (
548  IOobject
549  (
550  "ddtCorr("
551  + rho.name() + ',' + U.name() + ',' + phi.name() + ')',
552  mesh().time().timeName(),
553  mesh()
554  ),
555  this->fvcDdtPhiCoeff
556  (
557  rhoU0,
558  phi.oldTime(),
559  phiCorr,
560  rho.oldTime()
561  )*rDeltaT*phiCorr
562  )
563  );
564  }
565  else if
566  (
567  U.dimensions() == rho.dimensions()*dimVelocity
568  && phi.dimensions() == rho.dimensions()*dimVelocity*dimArea
569  )
570  {
571  fluxFieldType phiCorr
572  (
573  phi.oldTime() - fvc::dotInterpolate(mesh().Sf(), U.oldTime())
574  );
575 
576  return tmp<fluxFieldType>
577  (
578  new fluxFieldType
579  (
580  IOobject
581  (
582  "ddtCorr("
583  + rho.name() + ',' + U.name() + ',' + phi.name() + ')',
584  mesh().time().timeName(),
585  mesh()
586  ),
587  this->fvcDdtPhiCoeff
588  (
589  U.oldTime(),
590  phi.oldTime(),
591  phiCorr,
592  rho.oldTime()
593  )*rDeltaT*phiCorr
594  )
595  );
596  }
597  else
598  {
600  << "dimensions of phi are not correct"
601  << abort(FatalError);
602 
603  return fluxFieldType::null();
604  }
605 }
606 
607 
608 template<class Type>
610 (
612 )
613 {
615  (
616  "meshPhi",
617  mesh(),
619  );
620 }
621 
622 
623 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
624 
625 } // End namespace fv
626 
627 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
628 
629 } // End namespace Foam
630 
631 // ************************************************************************* //
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
const word & name() const
Return name.
Definition: IOobject.H:303
static const surfaceScalarField & localRDeltaTf(const fvMesh &mesh)
Return the reciprocal of the local face time-step.
Definition: localEulerDdt.C:58
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
static tmp< GeometricField< typename innerProduct< vector, Type >::type, fvsPatchField, surfaceMesh > > dotInterpolate(const surfaceVectorField &Sf, const GeometricField< Type, fvPatchField, volMesh > &tvf)
Interpolate field onto faces.
virtual tmp< fluxFieldType > fvcDdtUfCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
virtual tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)
const dimensionSet dimVol(dimVolume)
Definition: dimensionSets.H:59
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
Generic GeometricField class.
Generic dimensioned Type class.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
const dimensionSet & dimensions() const
Return dimensions.
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
virtual tmp< fluxFieldType > fvcDdtPhiCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
This boundary condition is not designed to be evaluated; it is assumed that the value is assigned via...
const Type & value() const
Return const reference to value.
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:72
static const zero Zero
Definition: zero.H:97
static const volScalarField & localRDeltaT(const fvMesh &mesh)
Return the reciprocal of the local time-step.
Definition: localEulerDdt.C:46
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Field< Type > & source()
Definition: fvMatrix.H:292
virtual tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const dimensioned< Type > &)
const word & name() const
Return const reference to name.
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.
const dimensionSet dimDensity
virtual tmp< surfaceScalarField > meshPhi(const GeometricField< Type, fvPatchField, volMesh > &)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
A special matrix type and solver, designed for finite volume solutions of scalar equations.
const dimensionSet & dimensions() const
Return const reference to dimensions.
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
scalarField & diag()
Definition: lduMatrix.C:186
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.
const dimensionSet dimVelocity