localEulerDdtScheme.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 "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>
44 {
46 }
47 
48 
49 template<class Type>
52 (
53  const dimensioned<Type>& dt
54 )
55 {
56  const volScalarField& rDeltaT = localRDeltaT();
57 
58  IOobject ddtIOobject
59  (
60  "ddt(" + dt.name() + ')',
61  mesh().time().timeName(),
62  mesh()
63  );
64 
65  if (mesh().moving())
66  {
68  (
70  (
71  ddtIOobject,
72  mesh(),
74  (
75  "0",
76  dt.dimensions()/dimTime,
78  )
79  )
80  );
81 
82  tdtdt().internalField() =
83  rDeltaT.internalField()*dt.value()
84  *(1.0 - mesh().Vsc0()/mesh().Vsc());
85 
86  return tdtdt;
87  }
88  else
89  {
91  (
93  (
94  ddtIOobject,
95  mesh(),
97  (
98  "0",
99  dt.dimensions()/dimTime,
101  ),
103  )
104  );
105  }
106 }
107 
108 
109 template<class Type>
112 (
114 )
115 {
116  const volScalarField& rDeltaT = localRDeltaT();
117 
118  IOobject ddtIOobject
119  (
120  "ddt(" + vf.name() + ')',
121  mesh().time().timeName(),
122  mesh()
123  );
124 
125  if (mesh().moving())
126  {
128  (
130  (
131  ddtIOobject,
132  mesh(),
133  rDeltaT.dimensions()*vf.dimensions(),
134  rDeltaT.internalField()*
135  (
136  vf.internalField()
137  - vf.oldTime().internalField()*mesh().Vsc0()/mesh().Vsc()
138  ),
139  rDeltaT.boundaryField()*
140  (
141  vf.boundaryField() - vf.oldTime().boundaryField()
142  )
143  )
144  );
145  }
146  else
147  {
149  (
151  (
152  ddtIOobject,
153  rDeltaT*(vf - vf.oldTime())
154  )
155  );
156  }
157 }
158 
159 
160 template<class Type>
163 (
164  const dimensionedScalar& rho,
166 )
167 {
168  const volScalarField& rDeltaT = localRDeltaT();
169 
170  IOobject ddtIOobject
171  (
172  "ddt(" + rho.name() + ',' + vf.name() + ')',
173  mesh().time().timeName(),
174  mesh()
175  );
176 
177  if (mesh().moving())
178  {
180  (
182  (
183  ddtIOobject,
184  mesh(),
185  rDeltaT.dimensions()*rho.dimensions()*vf.dimensions(),
186  rDeltaT.internalField()*rho.value()*
187  (
188  vf.internalField()
189  - vf.oldTime().internalField()*mesh().Vsc0()/mesh().Vsc()
190  ),
191  rDeltaT.boundaryField()*rho.value()*
192  (
193  vf.boundaryField() - vf.oldTime().boundaryField()
194  )
195  )
196  );
197  }
198  else
199  {
201  (
203  (
204  ddtIOobject,
205  rDeltaT*rho*(vf - vf.oldTime())
206  )
207  );
208  }
209 }
210 
211 
212 template<class Type>
215 (
216  const volScalarField& rho,
218 )
219 {
220  const volScalarField& rDeltaT = localRDeltaT();
221 
222  IOobject ddtIOobject
223  (
224  "ddt(" + rho.name() + ',' + vf.name() + ')',
225  mesh().time().timeName(),
226  mesh()
227  );
228 
229  if (mesh().moving())
230  {
232  (
234  (
235  ddtIOobject,
236  mesh(),
237  rDeltaT.dimensions()*rho.dimensions()*vf.dimensions(),
238  rDeltaT.internalField()*
239  (
240  rho.internalField()*vf.internalField()
241  - rho.oldTime().internalField()
242  *vf.oldTime().internalField()*mesh().Vsc0()/mesh().Vsc()
243  ),
244  rDeltaT.boundaryField()*
245  (
246  rho.boundaryField()*vf.boundaryField()
247  - rho.oldTime().boundaryField()
248  *vf.oldTime().boundaryField()
249  )
250  )
251  );
252  }
253  else
254  {
256  (
258  (
259  ddtIOobject,
260  rDeltaT*(rho*vf - rho.oldTime()*vf.oldTime())
261  )
262  );
263  }
264 }
265 
266 
267 template<class Type>
270 (
271  const volScalarField& alpha,
272  const volScalarField& rho,
274 )
275 {
276  const volScalarField& rDeltaT = localRDeltaT();
277 
278  IOobject ddtIOobject
279  (
280  "ddt("+alpha.name()+','+rho.name()+','+vf.name()+')',
281  mesh().time().timeName(),
282  mesh()
283  );
284 
285  if (mesh().moving())
286  {
288  (
290  (
291  ddtIOobject,
292  mesh(),
293  rDeltaT.dimensions()
294  *alpha.dimensions()*rho.dimensions()*vf.dimensions(),
295  rDeltaT.internalField()*
296  (
297  alpha.internalField()
298  *rho.internalField()
299  *vf.internalField()
300 
301  - alpha.oldTime().internalField()
302  *rho.oldTime().internalField()
303  *vf.oldTime().internalField()*mesh().Vsc0()/mesh().Vsc()
304  ),
305  rDeltaT.boundaryField()*
306  (
307  alpha.boundaryField()
308  *rho.boundaryField()
309  *vf.boundaryField()
310 
311  - alpha.oldTime().boundaryField()
312  *rho.oldTime().boundaryField()
313  *vf.oldTime().boundaryField()
314  )
315  )
316  );
317  }
318  else
319  {
321  (
323  (
324  ddtIOobject,
325  rDeltaT
326  *(
327  alpha*rho*vf
328  - alpha.oldTime()*rho.oldTime()*vf.oldTime()
329  )
330  )
331  );
332  }
333 }
334 
335 
336 template<class Type>
339 (
341 )
342 {
343  tmp<fvMatrix<Type> > tfvm
344  (
345  new fvMatrix<Type>
346  (
347  vf,
348  vf.dimensions()*dimVol/dimTime
349  )
350  );
351 
352  fvMatrix<Type>& fvm = tfvm();
353 
354  const scalarField& rDeltaT = localRDeltaT().internalField();
355 
356  fvm.diag() = rDeltaT*mesh().Vsc();
357 
358  if (mesh().moving())
359  {
360  fvm.source() = rDeltaT*vf.oldTime().internalField()*mesh().Vsc0();
361  }
362  else
363  {
364  fvm.source() = rDeltaT*vf.oldTime().internalField()*mesh().Vsc();
365  }
366 
367  return tfvm;
368 }
369 
370 
371 template<class Type>
374 (
375  const dimensionedScalar& rho,
377 )
378 {
379  tmp<fvMatrix<Type> > tfvm
380  (
381  new fvMatrix<Type>
382  (
383  vf,
384  rho.dimensions()*vf.dimensions()*dimVol/dimTime
385  )
386  );
387  fvMatrix<Type>& fvm = tfvm();
388 
389  const scalarField& rDeltaT = localRDeltaT().internalField();
390 
391  fvm.diag() = rDeltaT*rho.value()*mesh().Vsc();
392 
393  if (mesh().moving())
394  {
395  fvm.source() = rDeltaT
396  *rho.value()*vf.oldTime().internalField()*mesh().Vsc0();
397  }
398  else
399  {
400  fvm.source() = rDeltaT
401  *rho.value()*vf.oldTime().internalField()*mesh().Vsc();
402  }
403 
404  return tfvm;
405 }
406 
407 
408 template<class Type>
411 (
412  const volScalarField& rho,
414 )
415 {
416  tmp<fvMatrix<Type> > tfvm
417  (
418  new fvMatrix<Type>
419  (
420  vf,
421  rho.dimensions()*vf.dimensions()*dimVol/dimTime
422  )
423  );
424  fvMatrix<Type>& fvm = tfvm();
425 
426  const scalarField& rDeltaT = localRDeltaT().internalField();
427 
428  fvm.diag() = rDeltaT*rho.internalField()*mesh().Vsc();
429 
430  if (mesh().moving())
431  {
432  fvm.source() = rDeltaT
433  *rho.oldTime().internalField()
434  *vf.oldTime().internalField()*mesh().Vsc0();
435  }
436  else
437  {
438  fvm.source() = rDeltaT
439  *rho.oldTime().internalField()
440  *vf.oldTime().internalField()*mesh().Vsc();
441  }
442 
443  return tfvm;
444 }
445 
446 
447 template<class Type>
450 (
451  const volScalarField& alpha,
452  const volScalarField& rho,
454 )
455 {
456  tmp<fvMatrix<Type> > tfvm
457  (
458  new fvMatrix<Type>
459  (
460  vf,
461  alpha.dimensions()*rho.dimensions()*vf.dimensions()*dimVol/dimTime
462  )
463  );
464  fvMatrix<Type>& fvm = tfvm();
465 
466  const scalarField& rDeltaT = localRDeltaT().internalField();
467 
468  fvm.diag() = rDeltaT*alpha.internalField()*rho.internalField()*mesh().Vsc();
469 
470  if (mesh().moving())
471  {
472  fvm.source() = rDeltaT
473  *alpha.oldTime().internalField()
474  *rho.oldTime().internalField()
475  *vf.oldTime().internalField()*mesh().Vsc0();
476  }
477  else
478  {
479  fvm.source() = rDeltaT
480  *alpha.oldTime().internalField()
481  *rho.oldTime().internalField()
482  *vf.oldTime().internalField()*mesh().Vsc();
483  }
484 
485  return tfvm;
486 }
487 
488 
489 template<class Type>
492 (
495 )
496 {
497  const surfaceScalarField rDeltaT(fvc::interpolate(localRDeltaT()));
498 
499  fluxFieldType phiCorr
500  (
501  mesh().Sf() & (Uf.oldTime() - fvc::interpolate(U.oldTime()))
502  );
503 
504  return tmp<fluxFieldType>
505  (
506  new fluxFieldType
507  (
508  IOobject
509  (
510  "ddtCorr(" + U.name() + ',' + Uf.name() + ')',
511  mesh().time().timeName(),
512  mesh()
513  ),
514  this->fvcDdtPhiCoeff
515  (
516  U.oldTime(),
517  (mesh().Sf() & Uf.oldTime()),
518  phiCorr
519  )
520  *rDeltaT*phiCorr
521  )
522  );
523 }
524 
525 
526 template<class Type>
529 (
531  const fluxFieldType& phi
532 )
533 {
534  const surfaceScalarField rDeltaT(fvc::interpolate(localRDeltaT()));
535 
536  fluxFieldType phiCorr
537  (
538  phi.oldTime() - (mesh().Sf() & fvc::interpolate(U.oldTime()))
539  );
540 
541  return tmp<fluxFieldType>
542  (
543  new fluxFieldType
544  (
545  IOobject
546  (
547  "ddtCorr(" + U.name() + ',' + phi.name() + ')',
548  mesh().time().timeName(),
549  mesh()
550  ),
551  this->fvcDdtPhiCoeff(U.oldTime(), phi.oldTime(), phiCorr)
552  *rDeltaT*phiCorr
553  )
554  );
555 }
556 
557 
558 template<class Type>
561 (
562  const volScalarField& rho,
565 )
566 {
567  if
568  (
569  U.dimensions() == dimVelocity
570  && Uf.dimensions() == dimDensity*dimVelocity
571  )
572  {
573  const surfaceScalarField rDeltaT(fvc::interpolate(localRDeltaT()));
574 
576  (
577  rho.oldTime()*U.oldTime()
578  );
579 
580  fluxFieldType phiCorr
581  (
582  mesh().Sf() & (Uf.oldTime() - fvc::interpolate(rhoU0))
583  );
584 
585  return tmp<fluxFieldType>
586  (
587  new fluxFieldType
588  (
589  IOobject
590  (
591  "ddtCorr("
592  + rho.name() + ',' + U.name() + ',' + Uf.name() + ')',
593  mesh().time().timeName(),
594  mesh()
595  ),
596  this->fvcDdtPhiCoeff
597  (
598  rhoU0,
599  mesh().Sf() & Uf.oldTime(),
600  phiCorr
601  )
602  *rDeltaT*phiCorr
603  )
604  );
605  }
606  else if
607  (
608  U.dimensions() == dimDensity*dimVelocity
609  && Uf.dimensions() == dimDensity*dimVelocity
610  )
611  {
612  return fvcDdtUfCorr(U, Uf);
613  }
614  else
615  {
617  << "dimensions of Uf are not correct"
618  << abort(FatalError);
619 
620  return fluxFieldType::null();
621  }
622 }
623 
624 
625 template<class Type>
628 (
629  const volScalarField& rho,
631  const fluxFieldType& phi
632 )
633 {
634  if
635  (
636  U.dimensions() == dimVelocity
637  && phi.dimensions() == rho.dimensions()*dimVelocity*dimArea
638  )
639  {
640  dimensionedScalar rDeltaT = 1.0/mesh().time().deltaT();
641 
643  (
644  rho.oldTime()*U.oldTime()
645  );
646 
647  fluxFieldType phiCorr
648  (
649  phi.oldTime() - (mesh().Sf() & fvc::interpolate(rhoU0))
650  );
651 
652  return tmp<fluxFieldType>
653  (
654  new fluxFieldType
655  (
656  IOobject
657  (
658  "ddtCorr("
659  + rho.name() + ',' + U.name() + ',' + phi.name() + ')',
660  mesh().time().timeName(),
661  mesh()
662  ),
663  this->fvcDdtPhiCoeff(rhoU0, phi.oldTime(), phiCorr)
664  *rDeltaT*phiCorr
665  )
666  );
667  }
668  else if
669  (
670  U.dimensions() == rho.dimensions()*dimVelocity
671  && phi.dimensions() == rho.dimensions()*dimVelocity*dimArea
672  )
673  {
674  return fvcDdtPhiCorr(U, phi);
675  }
676  else
677  {
679  << "dimensions of phi are not correct"
680  << abort(FatalError);
681 
682  return fluxFieldType::null();
683  }
684 }
685 
686 
687 template<class Type>
689 (
691 )
692 {
694  (
696  (
697  IOobject
698  (
699  "meshPhi",
700  mesh().time().timeName(),
701  mesh(),
704  false
705  ),
706  mesh(),
708  )
709  );
710 }
711 
712 
713 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
714 
715 } // End namespace fv
716 
717 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
718 
719 } // End namespace Foam
720 
721 // ************************************************************************* //
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::fv::localEulerDdt::localRDeltaT
static const volScalarField & localRDeltaT(const fvMesh &mesh)
Return the reciprocal of the local time-step.
Definition: localEulerDdt.C:45
Foam::fv::localEulerDdtScheme::fvcDdtUfCorr
tmp< fluxFieldType > fvcDdtUfCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
Definition: localEulerDdtScheme.C:492
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::dimDensity
const dimensionSet dimDensity
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
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::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:261
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
U
U
Definition: pEqn.H:46
Foam::fv::localEulerDdtScheme::fvcDdt
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const dimensioned< Type > &)
Definition: localEulerDdtScheme.C:52
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
localEulerDdtScheme.H
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::GeometricField::internalField
InternalField & internalField()
Return internal field.
Definition: GeometricField.C:724
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
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::fv::localEulerDdtScheme::meshPhi
tmp< surfaceScalarField > meshPhi(const GeometricField< Type, fvPatchField, volMesh > &)
Definition: localEulerDdtScheme.C:689
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
rho
rho
Definition: pEqn.H:3
fv
labelList fv(nPoints)
Foam::fv::localEulerDdtScheme::fvmDdt
tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)
Definition: localEulerDdtScheme.C:339
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
Foam::fv::localEulerDdtScheme::fvcDdtPhiCorr
tmp< fluxFieldType > fvcDdtPhiCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
Definition: localEulerDdtScheme.C:529
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::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Uf
Uf
Definition: pEqn.H:78
Foam::fv::localEulerDdtScheme::localRDeltaT
const volScalarField & localRDeltaT() const
Return the reciprocal of the local time-step.
Definition: localEulerDdtScheme.C:43