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-2023 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 {
45  return localEulerDdt::localRDeltaT(mesh());
46 }
47 
48 
49 template<class Type>
50 const surfaceScalarField& localEulerDdtScheme<Type>::localRDeltaTf() const
51 {
52  return localEulerDdt::localRDeltaTf(mesh());
53 }
54 
55 
56 template<class Type>
57 tmp<VolField<Type>>
59 (
60  const dimensioned<Type>& dt
61 )
62 {
63  const word ddtName("ddt(" + dt.name() + ')');
64 
65  return tmp<VolField<Type>>
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 (
87  const VolField<Type>& vf
88 )
89 {
90  const volScalarField& rDeltaT = localRDeltaT();
91 
92  const word ddtName("ddt(" + vf.name() + ')');
93 
94  return tmp<VolField<Type>>
95  (
97  (
98  ddtName,
99  rDeltaT*(vf - vf.oldTime())
100  )
101  );
102 }
103 
104 
105 template<class Type>
108 (
109  const dimensionedScalar& rho,
110  const VolField<Type>& vf
111 )
112 {
113  const volScalarField& rDeltaT = localRDeltaT();
114 
115  const word ddtName("ddt(" + rho.name() + ',' + vf.name() + ')');
116 
117  return tmp<VolField<Type>>
118  (
120  (
121  ddtName,
122  rDeltaT*rho*(vf - vf.oldTime())
123  )
124  );
125 }
126 
127 
128 template<class Type>
131 (
132  const volScalarField& rho,
133  const VolField<Type>& vf
134 )
135 {
136  const volScalarField& rDeltaT = localRDeltaT();
137 
138  const word ddtName("ddt(" + rho.name() + ',' + vf.name() + ')');
139 
140  return tmp<VolField<Type>>
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,
157  const VolField<Type>& vf
158 )
159 {
160  const volScalarField& rDeltaT = localRDeltaT();
161 
162  const word ddtName("ddt("+alpha.name()+','+rho.name()+','+vf.name()+')');
163 
164  return tmp<VolField<Type>>
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 (
183  const SurfaceField<Type>& sf
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 (
202  const VolField<Type>& vf
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,
230  const VolField<Type>& vf
231 )
232 {
233  tmp<fvMatrix<Type>> tfvm
234  (
235  new fvMatrix<Type>
236  (
237  vf,
238  rho.dimensions()*vf.dimensions()*dimVol/dimTime
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,
259  const VolField<Type>& vf
260 )
261 {
262  tmp<fvMatrix<Type>> tfvm
263  (
264  new fvMatrix<Type>
265  (
266  vf,
267  rho.dimensions()*vf.dimensions()*dimVol/dimTime
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,
290  const VolField<Type>& vf
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 VolField<Type>& 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 (U.boundaryField()[patchi].fixesValue())
343  {
344  ccbf[patchi] = 0;
345  }
346  }
347 
348  if (debug > 1)
349  {
350  InfoInFunction
351  << "ddtCouplingCoeff mean max min = "
352  << gAverage(ddtCouplingCoeff.primitiveField())
353  << " " << gMax(ddtCouplingCoeff.primitiveField())
354  << " " << gMin(ddtCouplingCoeff.primitiveField())
355  << endl;
356  }
357 
358  return tddtCouplingCoeff;
359 }
360 */
361 
362 
363 template<class Type>
366 (
367  const VolField<Type>& U,
368  const SurfaceField<Type>& Uf
369 )
370 {
371  const surfaceScalarField rDeltaT(fvc::interpolate(localRDeltaT()));
372 
373  fluxFieldType phiUf0(mesh().Sf() & Uf.oldTime());
374  fluxFieldType phiCorr
375  (
376  phiUf0 - fvc::dotInterpolate(mesh().Sf(), U.oldTime())
377  );
378 
379  return tmp<fluxFieldType>
380  (
381  new fluxFieldType
382  (
383  IOobject
384  (
385  "ddtCorr(" + U.name() + ',' + Uf.name() + ')',
386  mesh().time().name(),
387  mesh()
388  ),
389  this->fvcDdtPhiCoeff(U.oldTime(), phiUf0, phiCorr)
390  *rDeltaT*phiCorr
391  )
392  );
393 }
394 
395 
396 template<class Type>
399 (
400  const VolField<Type>& U,
401  const fluxFieldType& phi
402 )
403 {
404  const surfaceScalarField rDeltaT(fvc::interpolate(localRDeltaT()));
405 
406  fluxFieldType phiCorr
407  (
408  phi.oldTime() - fvc::dotInterpolate(mesh().Sf(), U.oldTime())
409  );
410 
411  return tmp<fluxFieldType>
412  (
413  new fluxFieldType
414  (
415  IOobject
416  (
417  "ddtCorr(" + U.name() + ',' + phi.name() + ')',
418  mesh().time().name(),
419  mesh()
420  ),
421  this->fvcDdtPhiCoeff(U.oldTime(), phi.oldTime(), phiCorr)
422  *rDeltaT*phiCorr
423  )
424  );
425 }
426 
427 
428 template<class Type>
431 (
432  const volScalarField& rho,
433  const VolField<Type>& U,
434  const SurfaceField<Type>& rhoUf
435 )
436 {
437  const surfaceScalarField rDeltaT(fvc::interpolate(localRDeltaT()));
438 
439  if
440  (
441  U.dimensions() == dimVelocity
442  && rhoUf.dimensions() == dimDensity*dimVelocity
443  )
444  {
445  VolField<Type> rhoU0
446  (
447  rho.oldTime()*U.oldTime()
448  );
449 
450  fluxFieldType phiUf0(mesh().Sf() & rhoUf.oldTime());
451  fluxFieldType phiCorr(phiUf0 - fvc::dotInterpolate(mesh().Sf(), rhoU0));
452 
453  return tmp<fluxFieldType>
454  (
455  new fluxFieldType
456  (
457  IOobject
458  (
459  "ddtCorr("
460  + rho.name() + ',' + U.name() + ',' + rhoUf.name() + ')',
461  mesh().time().name(),
462  mesh()
463  ),
464  this->fvcDdtPhiCoeff(rhoU0, phiUf0, phiCorr, rho.oldTime())
465  *rDeltaT*phiCorr
466  )
467  );
468  }
469  else if
470  (
471  U.dimensions() == dimDensity*dimVelocity
472  && rhoUf.dimensions() == dimDensity*dimVelocity
473  )
474  {
475  fluxFieldType phiUf0(mesh().Sf() & rhoUf.oldTime());
476  fluxFieldType phiCorr
477  (
478  phiUf0 - fvc::dotInterpolate(mesh().Sf(), U.oldTime())
479  );
480 
481  return tmp<fluxFieldType>
482  (
483  new fluxFieldType
484  (
485  IOobject
486  (
487  "ddtCorr("
488  + rho.name() + ',' + U.name() + ',' + rhoUf.name() + ')',
489  mesh().time().name(),
490  mesh()
491  ),
492  this->fvcDdtPhiCoeff
493  (
494  U.oldTime(),
495  phiUf0,
496  phiCorr,
497  rho.oldTime()
498  )*rDeltaT*phiCorr
499  )
500  );
501  }
502  else
503  {
505  << "dimensions of rhoUf are not correct"
506  << abort(FatalError);
507 
508  return fluxFieldType::null();
509  }
510 }
511 
512 
513 template<class Type>
516 (
517  const volScalarField& rho,
518  const VolField<Type>& U,
519  const fluxFieldType& phi
520 )
521 {
522  const surfaceScalarField rDeltaT(fvc::interpolate(localRDeltaT()));
523 
524  if
525  (
526  U.dimensions() == dimVelocity
527  && phi.dimensions() == rho.dimensions()*dimFlux
528  )
529  {
530  VolField<Type> rhoU0
531  (
532  rho.oldTime()*U.oldTime()
533  );
534 
535  fluxFieldType phiCorr
536  (
537  phi.oldTime() - fvc::dotInterpolate(mesh().Sf(), rhoU0)
538  );
539 
540  return tmp<fluxFieldType>
541  (
542  new fluxFieldType
543  (
544  IOobject
545  (
546  "ddtCorr("
547  + rho.name() + ',' + U.name() + ',' + phi.name() + ')',
548  mesh().time().name(),
549  mesh()
550  ),
551  this->fvcDdtPhiCoeff
552  (
553  rhoU0,
554  phi.oldTime(),
555  phiCorr,
556  rho.oldTime()
557  )*rDeltaT*phiCorr
558  )
559  );
560  }
561  else if
562  (
563  U.dimensions() == rho.dimensions()*dimVelocity
564  && phi.dimensions() == rho.dimensions()*dimFlux
565  )
566  {
567  fluxFieldType phiCorr
568  (
569  phi.oldTime() - fvc::dotInterpolate(mesh().Sf(), U.oldTime())
570  );
571 
572  return tmp<fluxFieldType>
573  (
574  new fluxFieldType
575  (
576  IOobject
577  (
578  "ddtCorr("
579  + rho.name() + ',' + U.name() + ',' + phi.name() + ')',
580  mesh().time().name(),
581  mesh()
582  ),
583  this->fvcDdtPhiCoeff
584  (
585  U.oldTime(),
586  phi.oldTime(),
587  phiCorr,
588  rho.oldTime()
589  )*rDeltaT*phiCorr
590  )
591  );
592  }
593  else
594  {
596  << "dimensions of phi are not correct"
597  << abort(FatalError);
598 
599  return fluxFieldType::null();
600  }
601 }
602 
603 
604 template<class Type>
607 (
608  const volScalarField& alpha,
609  const volScalarField& rho,
610  const VolField<Type>& U,
611  const SurfaceField<Type>& Uf
612 )
613 {
614  const surfaceScalarField rDeltaT(fvc::interpolate(localRDeltaT()));
615 
616  if (U.dimensions() == dimVelocity && Uf.dimensions() == dimVelocity)
617  {
618  const volScalarField alphaRho0(alpha.oldTime()*rho.oldTime());
619 
620  return fluxFieldType::New
621  (
622  "ddtCorr("
623  + alpha.name() + rho.name() + ',' + U.name() + ',' + Uf.name()
624  + ')',
625  this->fvcDdtPhiCoeff(U.oldTime(), mesh().Sf() & Uf.oldTime())
626  *rDeltaT
627  *(
628  mesh().Sf()
629  & (
630  fvc::interpolate(alphaRho0)*Uf.oldTime()
631  - fvc::interpolate(alphaRho0*U.oldTime())
632  )
633  )
634  );
635  }
636  else
637  {
639  << "dimensions of Uf are not correct"
640  << abort(FatalError);
641 
642  return fluxFieldType::null();
643  }
644 }
645 
646 
647 template<class Type>
650 (
651  const volScalarField& alpha,
652  const volScalarField& rho,
653  const VolField<Type>& U,
654  const fluxFieldType& phi
655 )
656 {
657  const surfaceScalarField rDeltaT(fvc::interpolate(localRDeltaT()));
658 
659  if (U.dimensions() == dimVelocity && phi.dimensions() == dimFlux)
660  {
661  const volScalarField alphaRho0(alpha.oldTime()*rho.oldTime());
662 
663  return fluxFieldType::New
664  (
665  "ddtCorr("
666  + alpha.name() + rho.name() + ',' + U.name() + ',' + phi.name()
667  + ')',
668  this->fvcDdtPhiCoeff(U.oldTime(), phi.oldTime())
669  *rDeltaT
670  *(
671  fvc::interpolate(alphaRho0)*phi.oldTime()
672  -fvc::dotInterpolate(mesh().Sf(), alphaRho0*U.oldTime())
673  )
674  );
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 (
690  const VolField<Type>&
691 )
692 {
694  (
695  "meshPhi",
696  mesh(),
698  );
699 }
700 
701 
702 template<class Type>
704 (
705  const VolField<Type>&,
706  const label patchi
707 )
708 {
709  return tmp<scalarField>
710  (
711  new scalarField(mesh().boundary()[patchi].size(), 0)
712  );
713 }
714 
715 
716 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
717 
718 } // End namespace fv
719 
720 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
721 
722 } // End namespace Foam
723 
724 // ************************************************************************* //
const dimensionSet & dimensions() const
Return dimensions.
Generic GeometricField class.
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const word & name() const
Return name.
Definition: IOobject.H:310
This boundary condition is not designed to be evaluated; it is assumed that the value is assigned via...
Generic dimensioned Type class.
const dimensionSet & dimensions() const
Return const reference to dimensions.
const word & name() const
Return const reference to name.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Field< Type > & source()
Definition: fvMatrix.H:307
virtual tmp< fluxFieldType > fvcDdtUfCorr(const VolField< Type > &U, const SurfaceField< Type > &Uf)
virtual tmp< fvMatrix< Type > > fvmDdt(const VolField< Type > &)
virtual tmp< surfaceScalarField > meshPhi(const VolField< Type > &)
virtual tmp< fluxFieldType > fvcDdtPhiCorr(const VolField< Type > &U, const fluxFieldType &phi)
virtual tmp< VolField< Type > > fvcDdt(const dimensioned< Type > &)
static const surfaceScalarField & localRDeltaTf(const fvMesh &mesh)
Return the reciprocal of the local face time-step.
Definition: localEulerDdt.C:58
static const volScalarField & localRDeltaT(const fvMesh &mesh)
Return the reciprocal of the local time-step.
Definition: localEulerDdt.C:46
scalarField & diag()
Definition: lduMatrix.C:186
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
A special matrix type and solver, designed for finite volume solutions of scalar equations.
label patchi
volScalarField sf(fieldObject, mesh)
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
static tmp< SurfaceField< typename innerProduct< vector, Type >::type > > dotInterpolate(const surfaceVectorField &Sf, const VolField< Type > &tvf)
Interpolate field onto faces.
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
static const zero Zero
Definition: zero.H:97
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
errorManip< error > abort(error &err)
Definition: errorManip.H:131
SurfaceField< scalar > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const dimensionSet dimTime
const dimensionSet dimVol
const dimensionSet dimDensity
const dimensionSet dimVolume
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
const dimensionSet dimVelocity
error FatalError
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
const dimensionSet dimFlux
faceListList boundary(nPatches)
labelList fv(nPoints)