EulerDdtScheme.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 | Copyright (C) 2011-2015 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 "EulerDdtScheme.H"
27 #include "surfaceInterpolate.H"
28 #include "fvcDiv.H"
29 #include "fvMatrices.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace fv
39 {
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 template<class Type>
44 tmp<GeometricField<Type, fvPatchField, volMesh> >
46 (
47  const dimensioned<Type>& dt
48 )
49 {
50  dimensionedScalar rDeltaT = 1.0/mesh().time().deltaT();
51 
52  IOobject ddtIOobject
53  (
54  "ddt("+dt.name()+')',
55  mesh().time().timeName(),
56  mesh()
57  );
58 
59  if (mesh().moving())
60  {
62  (
64  (
65  ddtIOobject,
66  mesh(),
68  (
69  "0",
70  dt.dimensions()/dimTime,
72  )
73  )
74  );
75 
76  tdtdt().internalField() =
77  rDeltaT.value()*dt.value()*(1.0 - mesh().Vsc0()/mesh().Vsc());
78 
79  return tdtdt;
80  }
81  else
82  {
84  (
86  (
87  ddtIOobject,
88  mesh(),
90  (
91  "0",
92  dt.dimensions()/dimTime,
94  ),
96  )
97  );
98  }
99 }
100 
101 
102 template<class Type>
105 (
107 )
108 {
109  dimensionedScalar rDeltaT = 1.0/mesh().time().deltaT();
110 
111  IOobject ddtIOobject
112  (
113  "ddt("+vf.name()+')',
114  mesh().time().timeName(),
115  mesh()
116  );
117 
118  if (mesh().moving())
119  {
121  (
123  (
124  ddtIOobject,
125  mesh(),
126  rDeltaT.dimensions()*vf.dimensions(),
127  rDeltaT.value()*
128  (
129  vf.internalField()
130  - vf.oldTime().internalField()*mesh().Vsc0()/mesh().Vsc()
131  ),
132  rDeltaT.value()*
133  (
134  vf.boundaryField() - vf.oldTime().boundaryField()
135  )
136  )
137  );
138  }
139  else
140  {
142  (
144  (
145  ddtIOobject,
146  rDeltaT*(vf - vf.oldTime())
147  )
148  );
149  }
150 }
151 
152 
153 template<class Type>
156 (
157  const dimensionedScalar& rho,
159 )
160 {
161  dimensionedScalar rDeltaT = 1.0/mesh().time().deltaT();
162 
163  IOobject ddtIOobject
164  (
165  "ddt("+rho.name()+','+vf.name()+')',
166  mesh().time().timeName(),
167  mesh()
168  );
169 
170  if (mesh().moving())
171  {
173  (
175  (
176  ddtIOobject,
177  mesh(),
178  rDeltaT.dimensions()*rho.dimensions()*vf.dimensions(),
179  rDeltaT.value()*rho.value()*
180  (
181  vf.internalField()
182  - vf.oldTime().internalField()*mesh().Vsc0()/mesh().Vsc()
183  ),
184  rDeltaT.value()*rho.value()*
185  (
186  vf.boundaryField() - vf.oldTime().boundaryField()
187  )
188  )
189  );
190  }
191  else
192  {
194  (
196  (
197  ddtIOobject,
198  rDeltaT*rho*(vf - vf.oldTime())
199  )
200  );
201  }
202 }
203 
204 
205 template<class Type>
208 (
209  const volScalarField& rho,
211 )
212 {
213  dimensionedScalar rDeltaT = 1.0/mesh().time().deltaT();
214 
215  IOobject ddtIOobject
216  (
217  "ddt("+rho.name()+','+vf.name()+')',
218  mesh().time().timeName(),
219  mesh()
220  );
221 
222  if (mesh().moving())
223  {
225  (
227  (
228  ddtIOobject,
229  mesh(),
230  rDeltaT.dimensions()*rho.dimensions()*vf.dimensions(),
231  rDeltaT.value()*
232  (
233  rho.internalField()*vf.internalField()
234  - rho.oldTime().internalField()
235  *vf.oldTime().internalField()*mesh().Vsc0()/mesh().Vsc()
236  ),
237  rDeltaT.value()*
238  (
239  rho.boundaryField()*vf.boundaryField()
240  - rho.oldTime().boundaryField()
241  *vf.oldTime().boundaryField()
242  )
243  )
244  );
245  }
246  else
247  {
249  (
251  (
252  ddtIOobject,
253  rDeltaT*(rho*vf - rho.oldTime()*vf.oldTime())
254  )
255  );
256  }
257 }
258 
259 
260 template<class Type>
263 (
264  const volScalarField& alpha,
265  const volScalarField& rho,
267 )
268 {
269  dimensionedScalar rDeltaT = 1.0/mesh().time().deltaT();
270 
271  IOobject ddtIOobject
272  (
273  "ddt("+alpha.name()+','+rho.name()+','+vf.name()+')',
274  mesh().time().timeName(),
275  mesh()
276  );
277 
278  if (mesh().moving())
279  {
281  (
283  (
284  ddtIOobject,
285  mesh(),
286  rDeltaT.dimensions()
287  *alpha.dimensions()*rho.dimensions()*vf.dimensions(),
288  rDeltaT.value()*
289  (
290  alpha.internalField()
291  *rho.internalField()
292  *vf.internalField()
293 
294  - alpha.oldTime().internalField()
295  *rho.oldTime().internalField()
296  *vf.oldTime().internalField()*mesh().Vsc0()/mesh().Vsc()
297  ),
298  rDeltaT.value()*
299  (
300  alpha.boundaryField()
301  *rho.boundaryField()
302  *vf.boundaryField()
303 
304  - alpha.oldTime().boundaryField()
305  *rho.oldTime().boundaryField()
306  *vf.oldTime().boundaryField()
307  )
308  )
309  );
310  }
311  else
312  {
314  (
316  (
317  ddtIOobject,
318  rDeltaT
319  *(
320  alpha*rho*vf
321  - alpha.oldTime()*rho.oldTime()*vf.oldTime()
322  )
323  )
324  );
325  }
326 }
327 
328 
329 template<class Type>
332 (
334 )
335 {
336  tmp<fvMatrix<Type> > tfvm
337  (
338  new fvMatrix<Type>
339  (
340  vf,
341  vf.dimensions()*dimVol/dimTime
342  )
343  );
344 
345  fvMatrix<Type>& fvm = tfvm();
346 
347  scalar rDeltaT = 1.0/mesh().time().deltaTValue();
348 
349  fvm.diag() = rDeltaT*mesh().Vsc();
350 
351  if (mesh().moving())
352  {
353  fvm.source() = rDeltaT*vf.oldTime().internalField()*mesh().Vsc0();
354  }
355  else
356  {
357  fvm.source() = rDeltaT*vf.oldTime().internalField()*mesh().Vsc();
358  }
359 
360  return tfvm;
361 }
362 
363 
364 template<class Type>
367 (
368  const dimensionedScalar& rho,
370 )
371 {
372  tmp<fvMatrix<Type> > tfvm
373  (
374  new fvMatrix<Type>
375  (
376  vf,
377  rho.dimensions()*vf.dimensions()*dimVol/dimTime
378  )
379  );
380  fvMatrix<Type>& fvm = tfvm();
381 
382  scalar rDeltaT = 1.0/mesh().time().deltaTValue();
383 
384  fvm.diag() = rDeltaT*rho.value()*mesh().Vsc();
385 
386  if (mesh().moving())
387  {
388  fvm.source() = rDeltaT
389  *rho.value()*vf.oldTime().internalField()*mesh().Vsc0();
390  }
391  else
392  {
393  fvm.source() = rDeltaT
394  *rho.value()*vf.oldTime().internalField()*mesh().Vsc();
395  }
396 
397  return tfvm;
398 }
399 
400 
401 template<class Type>
404 (
405  const volScalarField& rho,
407 )
408 {
409  tmp<fvMatrix<Type> > tfvm
410  (
411  new fvMatrix<Type>
412  (
413  vf,
414  rho.dimensions()*vf.dimensions()*dimVol/dimTime
415  )
416  );
417  fvMatrix<Type>& fvm = tfvm();
418 
419  scalar rDeltaT = 1.0/mesh().time().deltaTValue();
420 
421  fvm.diag() = rDeltaT*rho.internalField()*mesh().Vsc();
422 
423  if (mesh().moving())
424  {
425  fvm.source() = rDeltaT
426  *rho.oldTime().internalField()
427  *vf.oldTime().internalField()*mesh().Vsc0();
428  }
429  else
430  {
431  fvm.source() = rDeltaT
432  *rho.oldTime().internalField()
433  *vf.oldTime().internalField()*mesh().Vsc();
434  }
435 
436  return tfvm;
437 }
438 
439 
440 template<class Type>
443 (
444  const volScalarField& alpha,
445  const volScalarField& rho,
447 )
448 {
449  tmp<fvMatrix<Type> > tfvm
450  (
451  new fvMatrix<Type>
452  (
453  vf,
454  alpha.dimensions()*rho.dimensions()*vf.dimensions()*dimVol/dimTime
455  )
456  );
457  fvMatrix<Type>& fvm = tfvm();
458 
459  scalar rDeltaT = 1.0/mesh().time().deltaTValue();
460 
461  fvm.diag() = rDeltaT*alpha.internalField()*rho.internalField()*mesh().Vsc();
462 
463  if (mesh().moving())
464  {
465  fvm.source() = rDeltaT
466  *alpha.oldTime().internalField()
467  *rho.oldTime().internalField()
468  *vf.oldTime().internalField()*mesh().Vsc0();
469  }
470  else
471  {
472  fvm.source() = rDeltaT
473  *alpha.oldTime().internalField()
474  *rho.oldTime().internalField()
475  *vf.oldTime().internalField()*mesh().Vsc();
476  }
477 
478  return tfvm;
479 }
480 
481 
482 template<class Type>
485 (
488 )
489 {
490  dimensionedScalar rDeltaT = 1.0/mesh().time().deltaT();
491 
492  fluxFieldType phiCorr
493  (
494  mesh().Sf() & (Uf.oldTime() - fvc::interpolate(U.oldTime()))
495  );
496 
497  return tmp<fluxFieldType>
498  (
499  new fluxFieldType
500  (
501  IOobject
502  (
503  "ddtCorr(" + U.name() + ',' + Uf.name() + ')',
504  mesh().time().timeName(),
505  mesh()
506  ),
507  this->fvcDdtPhiCoeff
508  (
509  U.oldTime(),
510  mesh().Sf() & Uf.oldTime(),
511  phiCorr
512  )
513  *rDeltaT*phiCorr
514  )
515  );
516 }
517 
518 
519 template<class Type>
522 (
524  const fluxFieldType& phi
525 )
526 {
527  dimensionedScalar rDeltaT = 1.0/mesh().time().deltaT();
528 
529  fluxFieldType phiCorr
530  (
531  phi.oldTime() - (mesh().Sf() & fvc::interpolate(U.oldTime()))
532  );
533 
534  return tmp<fluxFieldType>
535  (
536  new fluxFieldType
537  (
538  IOobject
539  (
540  "ddtCorr(" + U.name() + ',' + phi.name() + ')',
541  mesh().time().timeName(),
542  mesh()
543  ),
544  this->fvcDdtPhiCoeff(U.oldTime(), phi.oldTime(), phiCorr)
545  *rDeltaT*phiCorr
546  )
547  );
548 }
549 
550 
551 template<class Type>
554 (
555  const volScalarField& rho,
558 )
559 {
560  if
561  (
562  U.dimensions() == dimVelocity
563  && Uf.dimensions() == rho.dimensions()*dimVelocity
564  )
565  {
566  dimensionedScalar rDeltaT = 1.0/mesh().time().deltaT();
567 
569  (
570  rho.oldTime()*U.oldTime()
571  );
572 
573  fluxFieldType phiCorr
574  (
575  mesh().Sf() & (Uf.oldTime() - fvc::interpolate(rhoU0))
576  );
577 
578  return tmp<fluxFieldType>
579  (
580  new fluxFieldType
581  (
582  IOobject
583  (
584  "ddtCorr("
585  + rho.name() + ',' + U.name() + ',' + Uf.name() + ')',
586  mesh().time().timeName(),
587  mesh()
588  ),
589  this->fvcDdtPhiCoeff
590  (
591  rhoU0,
592  mesh().Sf() & Uf.oldTime(),
593  phiCorr
594  )
595  *rDeltaT*phiCorr
596  )
597  );
598  }
599  else if
600  (
601  U.dimensions() == rho.dimensions()*dimVelocity
602  && Uf.dimensions() == rho.dimensions()*dimVelocity
603  )
604  {
605  return fvcDdtUfCorr(U, Uf);
606  }
607  else
608  {
610  << "dimensions of Uf are not correct"
611  << abort(FatalError);
612 
613  return fluxFieldType::null();
614  }
615 }
616 
617 
618 template<class Type>
621 (
622  const volScalarField& rho,
624  const fluxFieldType& phi
625 )
626 {
627  if
628  (
629  U.dimensions() == dimVelocity
630  && phi.dimensions() == rho.dimensions()*dimVelocity*dimArea
631  )
632  {
633  dimensionedScalar rDeltaT = 1.0/mesh().time().deltaT();
634 
636  (
637  rho.oldTime()*U.oldTime()
638  );
639 
640  fluxFieldType phiCorr
641  (
642  phi.oldTime() - (mesh().Sf() & fvc::interpolate(rhoU0))
643  );
644 
645  return tmp<fluxFieldType>
646  (
647  new fluxFieldType
648  (
649  IOobject
650  (
651  "ddtCorr("
652  + rho.name() + ',' + U.name() + ',' + phi.name() + ')',
653  mesh().time().timeName(),
654  mesh()
655  ),
656  this->fvcDdtPhiCoeff(rhoU0, phi.oldTime(), phiCorr)
657  *rDeltaT*phiCorr
658  )
659  );
660  }
661  else if
662  (
663  U.dimensions() == rho.dimensions()*dimVelocity
664  && phi.dimensions() == rho.dimensions()*dimVelocity*dimArea
665  )
666  {
667  return fvcDdtPhiCorr(U, phi);
668  }
669  else
670  {
672  << "dimensions of phi are not correct"
673  << abort(FatalError);
674 
675  return fluxFieldType::null();
676  }
677 }
678 
679 
680 template<class Type>
682 (
684 )
685 {
686  return mesh().phi();
687 }
688 
689 
690 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
691 
692 } // End namespace fv
693 
694 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
695 
696 } // End namespace Foam
697 
698 // ************************************************************************* //
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::fv::EulerDdtScheme::fvcDdtPhiCorr
tmp< fluxFieldType > fvcDdtPhiCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
Definition: EulerDdtScheme.C:522
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::fv::EulerDdtScheme::meshPhi
tmp< surfaceScalarField > meshPhi(const GeometricField< Type, fvPatchField, volMesh > &)
Definition: EulerDdtScheme.C:682
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
Foam::fvc::interpolate
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Definition: surfaceInterpolate.C:110
Foam::GeometricField::oldTime
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
Definition: GeometricField.C:804
fvcDiv.H
Calculate the divergence of the given field.
Foam::dimensioned::name
const word & name() const
Return const reference to name.
Definition: dimensionedType.C:235
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::fv::EulerDdtScheme::fvcDdtUfCorr
tmp< fluxFieldType > fvcDdtUfCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
Definition: EulerDdtScheme.C:485
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:261
U
U
Definition: pEqn.H:46
fvMatrices.H
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
EulerDdtScheme.H
Foam::GeometricField::internalField
InternalField & internalField()
Return internal field.
Definition: GeometricField.C:724
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::dimensioned< Type >
surfaceInterpolate.H
Surface Interpolation.
Foam::calculatedFvPatchField
This boundary condition is not designed to be evaluated; it is assmued that the value is assigned via...
Definition: calculatedFvPatchField.H:65
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
rho
rho
Definition: pEqn.H:3
fv
labelList fv(nPoints)
Foam::fv::EulerDdtScheme::fvmDdt
tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)
Definition: EulerDdtScheme.C:332
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:50
Foam::fvMatrix::source
Field< Type > & source()
Definition: fvMatrix.H:291
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::fvMatrix< Type >
Foam::dimVol
const dimensionSet dimVol(dimVolume)
Definition: dimensionSets.H:59
Foam::dimensioned::dimensions
const dimensionSet & dimensions() const
Return const reference to dimensions.
Definition: dimensionedType.C:248
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::fv::EulerDdtScheme::fvcDdt
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const dimensioned< Type > &)
Definition: EulerDdtScheme.C:46
Uf
Uf
Definition: pEqn.H:78