CrankNicolsonDdtScheme.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-2024 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 "CrankNicolsonDdtScheme.H"
27 #include "surfaceInterpolate.H"
28 #include "fvcDiv.H"
29 #include "fvMatrices.H"
30 #include "Constant.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
39 namespace fv
40 {
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 template<class Type>
45 template<class GeoField>
46 CrankNicolsonDdtScheme<Type>::DDt0Field<GeoField>::DDt0Field
47 (
48  const IOobject& io,
49  const fvMesh& mesh
50 )
51 :
52  GeoField(io, mesh),
53  startTimeIndex_(-2) // This field is for a restart and thus correct so set
54  // the start time-index to correspond to a previous run
55 {
56  // Set the time-index to the beginning of the run to ensure the field
57  // is updated during the first time-step
58  this->timeIndex() = mesh.time().startTimeIndex();
59 }
60 
61 
62 template<class Type>
63 template<class GeoField>
64 CrankNicolsonDdtScheme<Type>::DDt0Field<GeoField>::DDt0Field
65 (
66  const IOobject& io,
67  const fvMesh& mesh,
68  const dimensioned<typename GeoField::value_type>& dimType
69 )
70 :
71  GeoField(io, mesh, dimType),
72  startTimeIndex_(mesh.time().timeIndex())
73 {}
74 
75 
76 template<class Type>
77 template<class GeoField>
78 label CrankNicolsonDdtScheme<Type>::DDt0Field<GeoField>::
79 startTimeIndex() const
80 {
81  return startTimeIndex_;
82 }
83 
84 
85 template<class Type>
86 template<class GeoField>
87 GeoField& CrankNicolsonDdtScheme<Type>::DDt0Field<GeoField>::
88 operator()()
89 {
90  return *this;
91 }
92 
93 
94 template<class Type>
95 template<class GeoField>
97 operator=(const GeoField& gf)
98 {
99  GeoField::operator=(gf);
100 }
101 
102 
103 template<class Type>
104 template<class GeoField>
105 typename CrankNicolsonDdtScheme<Type>::template DDt0Field<GeoField>&
107 (
108  const word& unTypedName,
109  const dimensionSet& dims
110 )
111 {
112  const word name = typedName(unTypedName);
113 
114  if (!mesh().objectRegistry::template foundObject<GeoField>(name))
115  {
116  const Time& runTime = mesh().time();
117  word startTimeName = runTime.timeName(runTime.startTime().value());
118 
119  if
120  (
121  (
122  runTime.timeIndex() == runTime.startTimeIndex()
123  || runTime.timeIndex() == runTime.startTimeIndex() + 1
124  )
125  && typeIOobject<DDt0Field<GeoField>>
126  (
127  name,
128  startTimeName,
129  mesh()
130  ).headerOk()
131  )
132  {
134  (
135  new DDt0Field<GeoField>
136  (
137  IOobject
138  (
139  name,
140  startTimeName,
141  mesh(),
144  ),
145  mesh()
146  )
147  );
148  }
149  else
150  {
152  (
153  new DDt0Field<GeoField>
154  (
155  IOobject
156  (
157  name,
158  mesh().time().name(),
159  mesh(),
162  ),
163  mesh(),
165  (
166  "0",
167  dims/dimTime,
168  Zero
169  )
170  )
171  );
172  }
173  }
174 
175  DDt0Field<GeoField>& ddt0 = static_cast<DDt0Field<GeoField>&>
176  (
177  mesh().objectRegistry::template lookupObjectRef<GeoField>(name)
178  );
179 
180  return ddt0;
181 }
182 
183 
184 template<class Type>
185 template<class GeoField>
187 (
188  const DDt0Field<GeoField>& ddt0
189 ) const
190 {
191  return ddt0.timeIndex() != mesh().time().timeIndex();
192 }
193 
194 
195 template<class Type>
196 template<class GeoField>
197 scalar CrankNicolsonDdtScheme<Type>::coef_
198 (
199  const DDt0Field<GeoField>& ddt0
200 ) const
201 {
202  if (mesh().time().timeIndex() > ddt0.startTimeIndex())
203  {
204  return 1 + ocCoeff();
205  }
206  else
207  {
208  return 1;
209  }
210 }
211 
212 
213 template<class Type>
214 template<class GeoField>
215 scalar CrankNicolsonDdtScheme<Type>::coef0_
216 (
217  const DDt0Field<GeoField>& ddt0
218 ) const
219 {
220  if (mesh().time().timeIndex() > ddt0.startTimeIndex() + 1)
221  {
222  return 1 + ocCoeff();
223  }
224  else
225  {
226  return 1;
227  }
228 }
229 
230 
231 template<class Type>
232 template<class GeoField>
233 dimensionedScalar CrankNicolsonDdtScheme<Type>::rDtCoef_
234 (
235  const DDt0Field<GeoField>& ddt0
236 ) const
237 {
238  return coef_(ddt0)/mesh().time().deltaT();
239 }
240 
241 
242 template<class Type>
243 template<class GeoField>
244 dimensionedScalar CrankNicolsonDdtScheme<Type>::rDtCoef0_
245 (
246  const DDt0Field<GeoField>& ddt0
247 ) const
248 {
249  return coef0_(ddt0)/mesh().time().deltaT0();
250 }
251 
252 
253 template<class Type>
254 template<class GeoField>
255 tmp<GeoField> CrankNicolsonDdtScheme<Type>::offCentre_
256 (
257  const GeoField& ddt0
258 ) const
259 {
260  if (ocCoeff() < 1)
261  {
262  return ocCoeff()*ddt0;
263  }
264  else
265  {
266  return ddt0;
267  }
268 }
269 
270 
271 template<class Type>
273 (
275 )
276 {
277  return bf;
278 }
279 
280 
281 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
282 
283 template<class Type>
285 :
286  ddtScheme<Type>(mesh),
287  ocCoeff_(new Function1s::Constant<scalar>("ocCoeff", 1))
288 {
289  // Ensure the old-old-time cell volumes are available
290  // for moving meshes
291  if (mesh.moving())
292  {
293  mesh.V00();
294  }
295 }
296 
297 
298 template<class Type>
300 (
301  const fvMesh& mesh,
302  Istream& is
303 )
304 :
305  ddtScheme<Type>(mesh, is)
306 {
307  token firstToken(is);
308 
309  if (firstToken.isNumber())
310  {
311  const scalar ocCoeff = firstToken.number();
312 
313  if (ocCoeff < 0 || ocCoeff > 1)
314  {
316  << "Off-centreing coefficient = " << ocCoeff
317  << " should be >= 0 and <= 1"
318  << exit(FatalIOError);
319  }
320 
321  ocCoeff_ =
322  new Function1s::Constant<scalar>("ocCoeff", ocCoeff);
323  }
324  else
325  {
326  is.putBack(firstToken);
327  dictionary dict(is);
328  ocCoeff_ =
330  (
331  "ocCoeff",
332  mesh.time().userUnits(),
333  unitFraction,
334  dict
335  );
336  }
337 
338  // Ensure the old-old-time cell volumes are available
339  // for moving meshes
340  if (mesh.moving())
341  {
342  mesh.V00();
343  }
344 }
345 
346 
347 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
348 
349 template<class Type>
352 (
353  const dimensioned<Type>& dt
354 )
355 {
356  DDt0Field<VolField<Type>>& ddt0 =
357  ddt0_<VolField<Type>>
358  (
359  "ddt0(" + dt.name() + ')',
360  dt.dimensions()
361  );
362 
363  const word ddtName("ddt(" + dt.name() + ')');
364 
365  tmp<VolField<Type>> tdtdt
366  (
368  (
369  ddtName,
370  mesh(),
372  (
373  "0",
374  dt.dimensions()/dimTime,
375  Zero
376  )
377  )
378  );
379 
380  dimensionedScalar rDtCoef = rDtCoef_(ddt0);
381 
382  if (mesh().moving())
383  {
384  if (evaluate(ddt0))
385  {
386  dimensionedScalar rDtCoef0 = rDtCoef0_(ddt0);
387 
388  ddt0.internalFieldRef() =
389  (
390  (rDtCoef0*dt)*(mesh().V0() - mesh().V00())
391  - mesh().V00()*offCentre_(ddt0.internalField())
392  )/mesh().V0();
393  }
394 
395  tdtdt.ref().internalFieldRef() =
396  (
397  (rDtCoef*dt)*(mesh().V() - mesh().V0())
398  - mesh().V0()*offCentre_(ddt0.internalField())
399  )/mesh().V();
400  }
401 
402  return tdtdt;
403 }
404 
405 
406 template<class Type>
409 (
410  const VolField<Type>& vf
411 )
412 {
413  DDt0Field<VolField<Type>>& ddt0 =
414  ddt0_<VolField<Type>>
415  (
416  "ddt0(" + vf.name() + ')',
417  vf.dimensions()
418  );
419 
420  const word ddtName("ddt(" + vf.name() + ')');
421 
422  dimensionedScalar rDtCoef = rDtCoef_(ddt0);
423 
424  if (mesh().moving())
425  {
426  if (evaluate(ddt0))
427  {
428  const scalar rDtCoef0 = rDtCoef0_(ddt0).value();
429 
430  ddt0.primitiveFieldRef() =
431  (
432  rDtCoef0*
433  (
434  mesh().V0()*vf.oldTime().primitiveField()
435  - mesh().V00()*vf.oldTime().oldTime().primitiveField()
436  ) - mesh().V00()*offCentre_(ddt0.primitiveField())
437  )/mesh().V0();
438 
439  ddt0.boundaryFieldRef() =
440  (
441  rDtCoef0*
442  (
443  vf.oldTime().boundaryField()
444  - vf.oldTime().oldTime().boundaryField()
445  ) - offCentre_(ff(ddt0.boundaryField()))
446  );
447  }
448 
449  return VolField<Type>::New
450  (
451  ddtName,
452  (
453  rDtCoef*
454  (
455  mesh().V()*vf
456  - mesh().V0()*vf.oldTime()
457  ) - mesh().V0()*offCentre_(ddt0()())
458  )/mesh().V(),
459  rDtCoef.value()*
460  (
461  vf.boundaryField() - vf.oldTime().boundaryField()
462  ) - offCentre_(ff(ddt0.boundaryField()))
463  );
464  }
465  else
466  {
467  if (evaluate(ddt0))
468  {
469  ddt0 = rDtCoef0_(ddt0)*(vf.oldTime() - vf.oldTime().oldTime())
470  - offCentre_(ddt0());
471  }
472 
473  return VolField<Type>::New
474  (
475  ddtName,
476  rDtCoef*(vf - vf.oldTime()) - offCentre_(ddt0())
477  );
478  }
479 }
480 
481 
482 template<class Type>
485 (
486  const dimensionedScalar& rho,
487  const VolField<Type>& vf
488 )
489 {
490  DDt0Field<VolField<Type>>& ddt0 =
491  ddt0_<VolField<Type>>
492  (
493  "ddt0(" + rho.name() + ',' + vf.name() + ')',
494  rho.dimensions()*vf.dimensions()
495  );
496 
497  const word ddtName("ddt(" + rho.name() + ',' + vf.name() + ')');
498 
499  dimensionedScalar rDtCoef = rDtCoef_(ddt0);
500 
501  if (mesh().moving())
502  {
503  if (evaluate(ddt0))
504  {
505  const scalar rDtCoef0 = rDtCoef0_(ddt0).value();
506 
507  ddt0.primitiveFieldRef() =
508  (
509  rDtCoef0*rho.value()*
510  (
511  mesh().V0()*vf.oldTime().primitiveField()
512  - mesh().V00()*vf.oldTime().oldTime().primitiveField()
513  ) - mesh().V00()*offCentre_(ddt0.primitiveField())
514  )/mesh().V0();
515 
516  ddt0.boundaryFieldRef() =
517  (
518  rDtCoef0*rho.value()*
519  (
520  vf.oldTime().boundaryField()
521  - vf.oldTime().oldTime().boundaryField()
522  ) - offCentre_(ff(ddt0.boundaryField()))
523  );
524  }
525 
526  return VolField<Type>::New
527  (
528  ddtName,
529  (
530  rDtCoef*rho*
531  (
532  mesh().V()*vf()
533  - mesh().V0()*vf.oldTime()()
534  ) - mesh().V0()*offCentre_(ddt0())()
535  )/mesh().V(),
536  rDtCoef.value()*rho.value()*
537  (
538  vf.boundaryField() - vf.oldTime().boundaryField()
539  ) - offCentre_(ff(ddt0.boundaryField()))
540  );
541  }
542  else
543  {
544  if (evaluate(ddt0))
545  {
546  ddt0 = rDtCoef0_(ddt0)*rho*(vf.oldTime() - vf.oldTime().oldTime())
547  - offCentre_(ddt0());
548  }
549 
550  return VolField<Type>::New
551  (
552  ddtName,
553  rDtCoef*rho*(vf - vf.oldTime()) - offCentre_(ddt0())
554  );
555  }
556 }
557 
558 
559 template<class Type>
562 (
563  const volScalarField& rho,
564  const VolField<Type>& vf
565 )
566 {
567  DDt0Field<VolField<Type>>& ddt0 =
568  ddt0_<VolField<Type>>
569  (
570  "ddt0(" + rho.name() + ',' + vf.name() + ')',
571  rho.dimensions()*vf.dimensions()
572  );
573 
574  const word ddtName("ddt(" + rho.name() + ',' + vf.name() + ')');
575 
576  dimensionedScalar rDtCoef = rDtCoef_(ddt0);
577 
578  if (mesh().moving())
579  {
580  if (evaluate(ddt0))
581  {
582  const scalar rDtCoef0 = rDtCoef0_(ddt0).value();
583 
584  ddt0.primitiveFieldRef() =
585  (
586  rDtCoef0*
587  (
588  mesh().V0()*rho.oldTime().primitiveField()
589  *vf.oldTime().primitiveField()
590  - mesh().V00()*rho.oldTime().oldTime().primitiveField()
591  *vf.oldTime().oldTime().primitiveField()
592  ) - mesh().V00()*offCentre_(ddt0.primitiveField())
593  )/mesh().V0();
594 
595  ddt0.boundaryFieldRef() =
596  (
597  rDtCoef0*
598  (
599  rho.oldTime().boundaryField()
600  *vf.oldTime().boundaryField()
601  - rho.oldTime().oldTime().boundaryField()
602  *vf.oldTime().oldTime().boundaryField()
603  ) - offCentre_(ff(ddt0.boundaryField()))
604  );
605  }
606 
607  return VolField<Type>::New
608  (
609  ddtName,
610  (
611  rDtCoef*
612  (
613  mesh().V()*rho()*vf()
614  - mesh().V0()*rho.oldTime()()
615  *vf.oldTime()()
616  ) - mesh().V00()*offCentre_(ddt0())()
617  )/mesh().V(),
618  rDtCoef.value()*
619  (
620  rho.boundaryField()*vf.boundaryField()
621  - rho.oldTime().boundaryField()*vf.oldTime().boundaryField()
622  ) - offCentre_(ff(ddt0.boundaryField()))
623  );
624  }
625  else
626  {
627  if (evaluate(ddt0))
628  {
629  ddt0 = rDtCoef0_(ddt0)*
630  (
631  rho.oldTime()*vf.oldTime()
632  - rho.oldTime().oldTime()*vf.oldTime().oldTime()
633  ) - offCentre_(ddt0());
634  }
635 
636  return VolField<Type>::New
637  (
638  ddtName,
639  rDtCoef*(rho*vf - rho.oldTime()*vf.oldTime()) - offCentre_(ddt0())
640  );
641  }
642 }
643 
644 
645 template<class Type>
648 (
649  const volScalarField& alpha,
650  const volScalarField& rho,
651  const VolField<Type>& vf
652 )
653 {
654  DDt0Field<VolField<Type>>& ddt0 =
655  ddt0_<VolField<Type>>
656  (
657  "ddt0(" + alpha.name() + ',' + rho.name() + ',' + vf.name() + ')',
658  alpha.dimensions()*rho.dimensions()*vf.dimensions()
659  );
660 
661  const word ddtName
662  (
663  "ddt(" + alpha.name() + ',' + rho.name() + ',' + vf.name() + ')'
664  );
665 
666  dimensionedScalar rDtCoef = rDtCoef_(ddt0);
667 
668  if (mesh().moving())
669  {
670  if (evaluate(ddt0))
671  {
672  const scalar rDtCoef0 = rDtCoef0_(ddt0).value();
673 
674  ddt0.primitiveFieldRef() =
675  (
676  rDtCoef0*
677  (
678  mesh().V0()
679  *alpha.oldTime().primitiveField()
680  *rho.oldTime().primitiveField()
681  *vf.oldTime().primitiveField()
682 
683  - mesh().V00()
684  *alpha.oldTime().oldTime().primitiveField()
685  *rho.oldTime().oldTime().primitiveField()
686  *vf.oldTime().oldTime().primitiveField()
687  ) - mesh().V00()*offCentre_(ddt0.primitiveField())
688  )/mesh().V0();
689 
690  ddt0.boundaryFieldRef() =
691  (
692  rDtCoef0*
693  (
694  alpha.oldTime().boundaryField()
695  *rho.oldTime().boundaryField()
696  *vf.oldTime().boundaryField()
697 
698  - alpha.oldTime().oldTime().boundaryField()
699  *rho.oldTime().oldTime().boundaryField()
700  *vf.oldTime().oldTime().boundaryField()
701  ) - offCentre_(ff(ddt0.boundaryField()))
702  );
703  }
704 
705  return VolField<Type>::New
706  (
707  ddtName,
708  (
709  rDtCoef.value()*
710  (
711  mesh().V()*alpha()*rho()*vf()
712  - mesh().V0()*alpha.oldTime()()*rho.oldTime()()*vf.oldTime()()
713  ) - mesh().V00()*offCentre_(ddt0())()
714  )/mesh().V(),
715  rDtCoef.value()*
716  (
717  alpha.boundaryField()
718  *rho.boundaryField()
719  *vf.boundaryField()
720 
721  - alpha.oldTime().boundaryField()
722  *rho.oldTime().boundaryField()
723  *vf.oldTime().boundaryField()
724  ) - offCentre_(ff(ddt0.boundaryField()))
725  );
726  }
727  else
728  {
729  if (evaluate(ddt0))
730  {
731  ddt0 = rDtCoef0_(ddt0)*
732  (
733  alpha.oldTime()
734  *rho.oldTime()
735  *vf.oldTime()
736 
737  - alpha.oldTime().oldTime()
738  *rho.oldTime().oldTime()
739  *vf.oldTime().oldTime()
740  ) - offCentre_(ddt0());
741  }
742 
743  return VolField<Type>::New
744  (
745  ddtName,
746  rDtCoef
747  *(
748  alpha*rho*vf
749  - alpha.oldTime()*rho.oldTime()*vf.oldTime()
750  )
751  - offCentre_(ddt0())
752  );
753  }
754 }
755 
756 
757 template<class Type>
760 (
761  const VolField<Type>& vf
762 )
763 {
764  DDt0Field<VolField<Type>>& ddt0 =
765  ddt0_<VolField<Type>>
766  (
767  "ddt0(" + vf.name() + ')',
768  vf.dimensions()
769  );
770 
771  tmp<fvMatrix<Type>> tfvm
772  (
773  new fvMatrix<Type>
774  (
775  vf,
777  )
778  );
779 
780  fvMatrix<Type>& fvm = tfvm.ref();
781 
782  const scalar rDtCoef = rDtCoef_(ddt0).value();
783  fvm.diag() = rDtCoef*mesh().V();
784 
785  vf.oldTime().oldTime();
786 
787  if (mesh().moving())
788  {
789  if (evaluate(ddt0))
790  {
791  const scalar rDtCoef0 = rDtCoef0_(ddt0).value();
792 
793  ddt0.primitiveFieldRef() =
794  (
795  rDtCoef0*
796  (
797  mesh().V0()*vf.oldTime().primitiveField()
798  - mesh().V00()*vf.oldTime().oldTime().primitiveField()
799  )
800  - mesh().V00()*offCentre_(ddt0.primitiveField())
801  )/mesh().V0();
802 
803  ddt0.boundaryFieldRef() =
804  (
805  rDtCoef0*
806  (
807  vf.oldTime().boundaryField()
808  - vf.oldTime().oldTime().boundaryField()
809  )
810  - offCentre_(ff(ddt0.boundaryField()))
811  );
812  }
813 
814  fvm.source() =
815  (
816  rDtCoef*vf.oldTime().primitiveField()
817  + offCentre_(ddt0.primitiveField())
818  )*mesh().V0();
819  }
820  else
821  {
822  if (evaluate(ddt0))
823  {
824  ddt0 = rDtCoef0_(ddt0)*(vf.oldTime() - vf.oldTime().oldTime())
825  - offCentre_(ddt0());
826  }
827 
828  fvm.source() =
829  (
830  rDtCoef*vf.oldTime().primitiveField()
831  + offCentre_(ddt0.primitiveField())
832  )*mesh().V();
833  }
834 
835  return tfvm;
836 }
837 
838 
839 template<class Type>
842 (
843  const dimensionedScalar& rho,
844  const VolField<Type>& vf
845 )
846 {
847  DDt0Field<VolField<Type>>& ddt0 =
848  ddt0_<VolField<Type>>
849  (
850  "ddt0(" + rho.name() + ',' + vf.name() + ')',
851  rho.dimensions()*vf.dimensions()
852  );
853 
854  tmp<fvMatrix<Type>> tfvm
855  (
856  new fvMatrix<Type>
857  (
858  vf,
859  rho.dimensions()*vf.dimensions()*dimVolume/dimTime
860  )
861  );
862  fvMatrix<Type>& fvm = tfvm.ref();
863 
864  const scalar rDtCoef = rDtCoef_(ddt0).value();
865  fvm.diag() = rDtCoef*rho.value()*mesh().V();
866 
867  vf.oldTime().oldTime();
868 
869  if (mesh().moving())
870  {
871  if (evaluate(ddt0))
872  {
873  const scalar rDtCoef0 = rDtCoef0_(ddt0).value();
874 
875  ddt0.primitiveFieldRef() =
876  (
877  rDtCoef0*rho.value()*
878  (
879  mesh().V0()*vf.oldTime().primitiveField()
880  - mesh().V00()*vf.oldTime().oldTime().primitiveField()
881  )
882  - mesh().V00()*offCentre_(ddt0.primitiveField())
883  )/mesh().V0();
884 
885  ddt0.boundaryFieldRef() =
886  (
887  rDtCoef0*rho.value()*
888  (
889  vf.oldTime().boundaryField()
890  - vf.oldTime().oldTime().boundaryField()
891  )
892  - offCentre_(ff(ddt0.boundaryField()))
893  );
894  }
895 
896  fvm.source() =
897  (
898  rDtCoef*rho.value()*vf.oldTime().primitiveField()
899  + offCentre_(ddt0.primitiveField())
900  )*mesh().V0();
901  }
902  else
903  {
904  if (evaluate(ddt0))
905  {
906  ddt0 = rDtCoef0_(ddt0)*rho*(vf.oldTime() - vf.oldTime().oldTime())
907  - offCentre_(ddt0());
908  }
909 
910  fvm.source() =
911  (
912  rDtCoef*rho.value()*vf.oldTime().primitiveField()
913  + offCentre_(ddt0.primitiveField())
914  )*mesh().V();
915  }
916 
917  return tfvm;
918 }
919 
920 
921 template<class Type>
924 (
925  const volScalarField& rho,
926  const VolField<Type>& vf
927 )
928 {
929  DDt0Field<VolField<Type>>& ddt0 =
930  ddt0_<VolField<Type>>
931  (
932  "ddt0(" + rho.name() + ',' + vf.name() + ')',
933  rho.dimensions()*vf.dimensions()
934  );
935 
936  tmp<fvMatrix<Type>> tfvm
937  (
938  new fvMatrix<Type>
939  (
940  vf,
941  rho.dimensions()*vf.dimensions()*dimVolume/dimTime
942  )
943  );
944  fvMatrix<Type>& fvm = tfvm.ref();
945 
946  const scalar rDtCoef = rDtCoef_(ddt0).value();
947  fvm.diag() = rDtCoef*rho.primitiveField()*mesh().V();
948 
949  vf.oldTime().oldTime();
950  rho.oldTime().oldTime();
951 
952  if (mesh().moving())
953  {
954  if (evaluate(ddt0))
955  {
956  const scalar rDtCoef0 = rDtCoef0_(ddt0).value();
957 
958  ddt0.primitiveFieldRef() =
959  (
960  rDtCoef0*
961  (
962  mesh().V0()*rho.oldTime().primitiveField()
963  *vf.oldTime().primitiveField()
964  - mesh().V00()*rho.oldTime().oldTime().primitiveField()
965  *vf.oldTime().oldTime().primitiveField()
966  )
967  - mesh().V00()*offCentre_(ddt0.primitiveField())
968  )/mesh().V0();
969 
970  ddt0.boundaryFieldRef() =
971  (
972  rDtCoef0*
973  (
974  rho.oldTime().boundaryField()
975  *vf.oldTime().boundaryField()
976  - rho.oldTime().oldTime().boundaryField()
977  *vf.oldTime().oldTime().boundaryField()
978  )
979  - offCentre_(ff(ddt0.boundaryField()))
980  );
981  }
982 
983  fvm.source() =
984  (
985  rDtCoef*rho.oldTime().primitiveField()*vf.oldTime().primitiveField()
986  + offCentre_(ddt0.primitiveField())
987  )*mesh().V0();
988  }
989  else
990  {
991  if (evaluate(ddt0))
992  {
993  ddt0 = rDtCoef0_(ddt0)*
994  (
995  rho.oldTime()*vf.oldTime()
996  - rho.oldTime().oldTime()*vf.oldTime().oldTime()
997  ) - offCentre_(ddt0());
998  }
999 
1000  fvm.source() =
1001  (
1002  rDtCoef*rho.oldTime().primitiveField()*vf.oldTime().primitiveField()
1003  + offCentre_(ddt0.primitiveField())
1004  )*mesh().V();
1005  }
1006 
1007  return tfvm;
1008 }
1009 
1010 
1011 template<class Type>
1014 (
1015  const volScalarField& alpha,
1016  const volScalarField& rho,
1017  const VolField<Type>& vf
1018 )
1019 {
1020  DDt0Field<VolField<Type>>& ddt0 =
1021  ddt0_<VolField<Type>>
1022  (
1023  "ddt0(" + alpha.name() + ',' + rho.name() + ',' + vf.name() + ')',
1024  alpha.dimensions()*rho.dimensions()*vf.dimensions()
1025  );
1026 
1027  tmp<fvMatrix<Type>> tfvm
1028  (
1029  new fvMatrix<Type>
1030  (
1031  vf,
1032  alpha.dimensions()
1033  *rho.dimensions()
1034  *vf.dimensions()
1035  *dimVolume
1036  /dimTime
1037  )
1038  );
1039  fvMatrix<Type>& fvm = tfvm.ref();
1040 
1041  const scalar rDtCoef = rDtCoef_(ddt0).value();
1042  fvm.diag() = rDtCoef*alpha.primitiveField()*rho.primitiveField()*mesh().V();
1043 
1044  vf.oldTime().oldTime();
1045  alpha.oldTime().oldTime();
1046  rho.oldTime().oldTime();
1047 
1048  if (mesh().moving())
1049  {
1050  if (evaluate(ddt0))
1051  {
1052  const scalar rDtCoef0 = rDtCoef0_(ddt0).value();
1053 
1054  ddt0.primitiveFieldRef() =
1055  (
1056  rDtCoef0*
1057  (
1058  mesh().V0()
1059  *alpha.oldTime().primitiveField()
1060  *rho.oldTime().primitiveField()
1061  *vf.oldTime().primitiveField()
1062 
1063  - mesh().V00()
1064  *alpha.oldTime().oldTime().primitiveField()
1065  *rho.oldTime().oldTime().primitiveField()
1066  *vf.oldTime().oldTime().primitiveField()
1067  )
1068  - mesh().V00()*offCentre_(ddt0.primitiveField())
1069  )/mesh().V0();
1070 
1071  ddt0.boundaryFieldRef() =
1072  (
1073  rDtCoef0*
1074  (
1075  alpha.oldTime().boundaryField()
1076  *rho.oldTime().boundaryField()
1077  *vf.oldTime().boundaryField()
1078 
1079  - alpha.oldTime().oldTime().boundaryField()
1080  *rho.oldTime().oldTime().boundaryField()
1081  *vf.oldTime().oldTime().boundaryField()
1082  )
1083  - offCentre_(ff(ddt0.boundaryField()))
1084  );
1085  }
1086 
1087  fvm.source() =
1088  (
1089  rDtCoef
1090  *alpha.oldTime().primitiveField()
1091  *rho.oldTime().primitiveField()
1092  *vf.oldTime().primitiveField()
1093  + offCentre_(ddt0.primitiveField())
1094  )*mesh().V0();
1095  }
1096  else
1097  {
1098  if (evaluate(ddt0))
1099  {
1100  ddt0 = rDtCoef0_(ddt0)*
1101  (
1102  alpha.oldTime()
1103  *rho.oldTime()
1104  *vf.oldTime()
1105 
1106  - alpha.oldTime().oldTime()
1107  *rho.oldTime().oldTime()
1108  *vf.oldTime().oldTime()
1109  ) - offCentre_(ddt0());
1110  }
1111 
1112  fvm.source() =
1113  (
1114  rDtCoef
1115  *alpha.oldTime().primitiveField()
1116  *rho.oldTime().primitiveField()
1117  *vf.oldTime().primitiveField()
1118  + offCentre_(ddt0.primitiveField())
1119  )*mesh().V();
1120  }
1121 
1122  return tfvm;
1123 }
1124 
1125 
1126 template<class Type>
1129 (
1130  const VolField<Type>& U,
1131  const SurfaceField<Type>& Uf
1132 )
1133 {
1134  DDt0Field<VolField<Type>>& ddt0 =
1135  ddt0_<VolField<Type>>
1136  (
1137  "ddtCorrDdt0(" + U.name() + ')',
1138  U.dimensions()
1139  );
1140 
1141  DDt0Field<SurfaceField<Type>>& dUfdt0 =
1142  ddt0_<SurfaceField<Type>>
1143  (
1144  "ddtCorrDdt0(" + Uf.name() + ')',
1145  Uf.dimensions()
1146  );
1147 
1148  dimensionedScalar rDtCoef = rDtCoef_(ddt0);
1149 
1150  if (evaluate(ddt0))
1151  {
1152  ddt0 =
1153  rDtCoef0_(ddt0)*(U.oldTime() - U.oldTime().oldTime())
1154  - offCentre_(ddt0());
1155  }
1156 
1157  if (evaluate(dUfdt0))
1158  {
1159  dUfdt0 =
1160  rDtCoef0_(dUfdt0)*(Uf.oldTime() - Uf.oldTime().oldTime())
1161  - offCentre_(dUfdt0());
1162  }
1163 
1164  return fluxFieldType::New
1165  (
1166  "ddtCorr(" + U.name() + ',' + Uf.name() + ')',
1167  this->fvcDdtPhiCoeff(U.oldTime(), mesh().Sf() & Uf.oldTime())
1168  *(
1169  mesh().Sf()
1170  & (
1171  (rDtCoef*Uf.oldTime() + offCentre_(dUfdt0()))
1172  - fvc::interpolate(rDtCoef*U.oldTime() + offCentre_(ddt0()))
1173  )
1174  )
1175  );
1176 }
1177 
1178 
1179 template<class Type>
1182 (
1183  const VolField<Type>& U,
1184  const fluxFieldType& phi
1185 )
1186 {
1187  DDt0Field<VolField<Type>>& ddt0 =
1188  ddt0_<VolField<Type>>
1189  (
1190  "ddtCorrDdt0(" + U.name() + ')',
1191  U.dimensions()
1192  );
1193 
1194  DDt0Field<fluxFieldType>& dphidt0 =
1195  ddt0_<fluxFieldType>
1196  (
1197  "ddtCorrDdt0(" + phi.name() + ')',
1198  phi.dimensions()
1199  );
1200 
1201  dimensionedScalar rDtCoef = rDtCoef_(ddt0);
1202 
1203  if (evaluate(ddt0))
1204  {
1205  ddt0 =
1206  rDtCoef0_(ddt0)*(U.oldTime() - U.oldTime().oldTime())
1207  - offCentre_(ddt0());
1208  }
1209 
1210  if (evaluate(dphidt0))
1211  {
1212  dphidt0 =
1213  rDtCoef0_(dphidt0)*(phi.oldTime() - phi.oldTime().oldTime())
1214  - offCentre_(dphidt0());
1215  }
1216 
1217  return fluxFieldType::New
1218  (
1219  "ddtCorr(" + U.name() + ',' + phi.name() + ')',
1220  this->fvcDdtPhiCoeff(U.oldTime(), phi.oldTime())
1221  *(
1222  (rDtCoef*phi.oldTime() + offCentre_(dphidt0()))
1224  (
1225  mesh().Sf(),
1226  rDtCoef*U.oldTime() + offCentre_(ddt0())
1227  )
1228  )
1229  );
1230 }
1231 
1232 
1233 template<class Type>
1236 (
1237  const volScalarField& rho,
1238  const VolField<Type>& U,
1239  const SurfaceField<Type>& rhoUf
1240 )
1241 {
1242  if
1243  (
1244  U.dimensions() == dimVelocity
1245  && rhoUf.dimensions() == rho.dimensions()*dimVelocity
1246  )
1247  {
1248  DDt0Field<VolField<Type>>& ddt0 =
1249  ddt0_<VolField<Type>>
1250  (
1251  "ddtCorrDdt0(" + rho.name() + ',' + U.name() + ')',
1252  rho.dimensions()*U.dimensions()
1253  );
1254 
1255  DDt0Field<SurfaceField<Type>>& drhoUfdt0 =
1256  ddt0_<SurfaceField<Type>>
1257  (
1258  "ddtCorrDdt0(" + rhoUf.name() + ')',
1259  rhoUf.dimensions()
1260  );
1261 
1262  dimensionedScalar rDtCoef = rDtCoef_(ddt0);
1263 
1264  VolField<Type> rhoU0
1265  (
1266  rho.oldTime()*U.oldTime()
1267  );
1268 
1269  if (evaluate(ddt0))
1270  {
1271  ddt0 =
1272  rDtCoef0_(ddt0)
1273  *(rhoU0 - rho.oldTime().oldTime()*U.oldTime().oldTime())
1274  - offCentre_(ddt0());
1275  }
1276 
1277  if (evaluate(drhoUfdt0))
1278  {
1279  drhoUfdt0 =
1280  rDtCoef0_(drhoUfdt0)
1281  *(rhoUf.oldTime() - rhoUf.oldTime().oldTime())
1282  - offCentre_(drhoUfdt0());
1283  }
1284 
1285  return fluxFieldType::New
1286  (
1287  "ddtCorr(" + rho.name() + ',' + U.name() + ',' + rhoUf.name() + ')',
1288  this->fvcDdtPhiCoeff
1289  (
1290  rhoU0,
1291  mesh().Sf() & rhoUf.oldTime(),
1292  rho.oldTime()
1293  )
1294  *(
1295  mesh().Sf()
1296  & (
1297  (rDtCoef*rhoUf.oldTime() + offCentre_(drhoUfdt0()))
1298  - fvc::interpolate(rDtCoef*rhoU0 + offCentre_(ddt0()))
1299  )
1300  )
1301  );
1302  }
1303  else if
1304  (
1305  U.dimensions() == rho.dimensions()*dimVelocity
1306  && rhoUf.dimensions() == rho.dimensions()*dimVelocity
1307  )
1308  {
1309  DDt0Field<VolField<Type>>& ddt0 =
1310  ddt0_<VolField<Type>>
1311  (
1312  "ddtCorrDdt0(" + U.name() + ')',
1313  U.dimensions()
1314  );
1315 
1316  DDt0Field<SurfaceField<Type>>& drhoUfdt0 =
1317  ddt0_<SurfaceField<Type>>
1318  (
1319  "ddtCorrDdt0(" + rhoUf.name() + ')',
1320  rhoUf.dimensions()
1321  );
1322 
1323  dimensionedScalar rDtCoef = rDtCoef_(ddt0);
1324 
1325  if (evaluate(ddt0))
1326  {
1327  ddt0 =
1328  rDtCoef0_(ddt0)*(U.oldTime() - U.oldTime().oldTime())
1329  - offCentre_(ddt0());
1330  }
1331 
1332  if (evaluate(drhoUfdt0))
1333  {
1334  drhoUfdt0 =
1335  rDtCoef0_(drhoUfdt0)
1336  *(rhoUf.oldTime() - rhoUf.oldTime().oldTime())
1337  - offCentre_(drhoUfdt0());
1338  }
1339 
1340  return fluxFieldType::New
1341  (
1342  "ddtCorr(" + U.name() + ',' + rhoUf.name() + ')',
1343  this->fvcDdtPhiCoeff
1344  (
1345  U.oldTime(),
1346  mesh().Sf() & rhoUf.oldTime(),
1347  rho.oldTime()
1348  )
1349  *(
1350  mesh().Sf()
1351  & (
1352  (rDtCoef*rhoUf.oldTime() + offCentre_(drhoUfdt0()))
1354  (
1355  rDtCoef*U.oldTime() + offCentre_(ddt0())
1356  )
1357  )
1358  )
1359  );
1360  }
1361  else
1362  {
1364  << "dimensions of rhoUf are not correct"
1365  << abort(FatalError);
1366 
1367  return fluxFieldType::null();
1368  }
1369 }
1370 
1371 
1372 template<class Type>
1375 (
1376  const volScalarField& rho,
1377  const VolField<Type>& U,
1378  const fluxFieldType& phi
1379 )
1380 {
1381  if
1382  (
1383  U.dimensions() == dimVelocity
1384  && phi.dimensions() == rho.dimensions()*dimVolumetricFlux
1385  )
1386  {
1387  DDt0Field<VolField<Type>>& ddt0 =
1388  ddt0_<VolField<Type>>
1389  (
1390  "ddtCorrDdt0(" + rho.name() + ',' + U.name() + ')',
1391  rho.dimensions()*U.dimensions()
1392  );
1393 
1394  DDt0Field<fluxFieldType>& dphidt0 =
1395  ddt0_<fluxFieldType>
1396  (
1397  "ddtCorrDdt0(" + phi.name() + ')',
1398  phi.dimensions()
1399  );
1400 
1401  dimensionedScalar rDtCoef = rDtCoef_(ddt0);
1402 
1403  VolField<Type> rhoU0
1404  (
1405  rho.oldTime()*U.oldTime()
1406  );
1407 
1408  if (evaluate(ddt0))
1409  {
1410  ddt0 =
1411  rDtCoef0_(ddt0)
1412  *(rhoU0 - rho.oldTime().oldTime()*U.oldTime().oldTime())
1413  - offCentre_(ddt0());
1414  }
1415 
1416  if (evaluate(dphidt0))
1417  {
1418  dphidt0 =
1419  rDtCoef0_(dphidt0)
1420  *(phi.oldTime() - phi.oldTime().oldTime())
1421  - offCentre_(dphidt0());
1422  }
1423 
1424  return fluxFieldType::New
1425  (
1426  "ddtCorr(" + rho.name() + ',' + U.name() + ',' + phi.name() + ')',
1427  this->fvcDdtPhiCoeff(rhoU0, phi.oldTime(), rho.oldTime())
1428  *(
1429  (rDtCoef*phi.oldTime() + offCentre_(dphidt0()))
1431  (
1432  mesh().Sf(),
1433  rDtCoef*rhoU0 + offCentre_(ddt0())
1434  )
1435  )
1436  );
1437  }
1438  else if
1439  (
1440  U.dimensions() == rho.dimensions()*dimVelocity
1441  && phi.dimensions() == rho.dimensions()*dimVolumetricFlux
1442  )
1443  {
1444  DDt0Field<VolField<Type>>& ddt0 =
1445  ddt0_<VolField<Type>>
1446  (
1447  "ddtCorrDdt0(" + U.name() + ')',
1448  U.dimensions()
1449  );
1450 
1451  DDt0Field<fluxFieldType>& dphidt0 =
1452  ddt0_<fluxFieldType>
1453  (
1454  "ddtCorrDdt0(" + phi.name() + ')',
1455  phi.dimensions()
1456  );
1457 
1458  dimensionedScalar rDtCoef = rDtCoef_(ddt0);
1459 
1460  if (evaluate(ddt0))
1461  {
1462  ddt0 =
1463  rDtCoef0_(ddt0)*(U.oldTime() - U.oldTime().oldTime())
1464  - offCentre_(ddt0());
1465  }
1466 
1467  if (evaluate(dphidt0))
1468  {
1469  dphidt0 =
1470  rDtCoef0_(dphidt0)*(phi.oldTime() - phi.oldTime().oldTime())
1471  - offCentre_(dphidt0());
1472  }
1473 
1474  return fluxFieldType::New
1475  (
1476  "ddtCorr(" + U.name() + ',' + phi.name() + ')',
1477  this->fvcDdtPhiCoeff(U.oldTime(), phi.oldTime(), rho.oldTime())
1478  *(
1479  (rDtCoef*phi.oldTime() + offCentre_(dphidt0()))
1481  (
1482  mesh().Sf(),
1483  rDtCoef*U.oldTime() + offCentre_(ddt0())
1484  )
1485  )
1486  );
1487  }
1488  else
1489  {
1491  << "dimensions of phi are not correct"
1492  << abort(FatalError);
1493 
1494  return fluxFieldType::null();
1495  }
1496 }
1497 
1498 
1499 template<class Type>
1501 (
1502  const VolField<Type>& vf
1503 )
1504 {
1505  DDt0Field<surfaceScalarField>& meshPhi0 =
1506  ddt0_<surfaceScalarField>("meshPhi0", dimVolume);
1507 
1508  if (evaluate(meshPhi0))
1509  {
1510  meshPhi0 =
1511  coef0_(meshPhi0)*mesh().phi().oldTime() - offCentre_(meshPhi0());
1512  }
1513 
1515  (
1516  mesh().phi().name(),
1517  coef_(meshPhi0)*mesh().phi() - offCentre_(meshPhi0())
1518  );
1519 }
1520 
1521 
1522 template<class Type>
1524 (
1525  const VolField<Type>& vf,
1526  const label patchi
1527 )
1528 {
1529  DDt0Field<surfaceScalarField>& meshPhi0 =
1530  ddt0_<surfaceScalarField>("meshPhi0", dimVolume);
1531 
1532  if (evaluate(meshPhi0))
1533  {
1534  meshPhi0 =
1535  coef0_(meshPhi0)*mesh().phi().oldTime() - offCentre_(meshPhi0());
1536  }
1537 
1538  return
1539  (
1540  coef_(meshPhi0)*mesh().phi().boundaryField()[patchi]
1541  - offCentre_
1542  (
1543  static_cast<const scalarField&>(meshPhi0().boundaryField()[patchi])
1544  )
1545  );
1546 }
1547 
1548 
1549 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1550 
1551 } // End namespace fv
1552 
1553 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1554 
1555 } // End namespace Foam
1556 
1557 // ************************************************************************* //
const dimensionSet & dimensions() const
Return dimensions.
Generic field type.
Definition: FieldField.H:77
static autoPtr< Function1< Type > > New(const word &name, const Function1s::unitConversions &units, const dictionary &dict)
Select from dictionary.
Definition: Function1New.C:32
Templated function that returns a constant value.
Definition: Constant.H:60
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
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:307
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
void putBack(const token &)
Put back token.
Definition: Istream.C:30
const Field0Type & oldTime() const
Return the old-time field.
Definition: OldTimeField.C:322
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:28
dimensionedScalar deltaT0() const
Return old time step.
Definition: TimeStateI.H:52
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:46
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
const unitConversion & userUnits() const
Return the user-time unit conversion.
Definition: Time.C:858
virtual dimensionedScalar startTime() const
Return start time.
Definition: Time.C:791
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
Definition: Time.C:630
virtual label startTimeIndex() const
Return start time index.
Definition: Time.C:774
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Dimension set for the base types.
Definition: dimensionSet.H:125
Generic dimensioned Type class.
const dimensionSet & dimensions() const
Return const reference to dimensions.
const Type & value() const
Return const reference to value.
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
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const DimensionedField< scalar, volMesh > & V00() const
Return old-old-time cell volumes.
const surfaceScalarField & phi() const
Return cell face motion fluxes.
const DimensionedField< scalar, volMesh > & V0() const
Return old-time cell volumes.
Second-oder Crank-Nicolson implicit ddt using the current and previous time-step fields as well as th...
virtual tmp< fvMatrix< Type > > fvmDdt(const VolField< Type > &)
const fvMesh & mesh() const
Return mesh reference.
scalar ocCoeff() const
Return the current off-centering coefficient.
virtual tmp< fluxFieldType > fvcDdtPhiCorr(const VolField< Type > &U, const fluxFieldType &phi)
virtual tmp< surfaceScalarField > meshPhi(const VolField< Type > &)
virtual tmp< fluxFieldType > fvcDdtUfCorr(const VolField< Type > &U, const SurfaceField< Type > &Uf)
void operator=(const CrankNicolsonDdtScheme &)=delete
Disallow default bitwise assignment.
virtual tmp< VolField< Type > > fvcDdt(const dimensioned< Type > &)
CrankNicolsonDdtScheme(const fvMesh &mesh)
Construct from mesh.
Abstract base class for ddt schemes.
Definition: ddtScheme.H:68
scalarField & diag()
Definition: lduMatrix.C:186
bool moving() const
Is mesh moving.
Definition: polyMesh.H:473
void store()
Transfer ownership of this object to its registry.
Definition: regIOobjectI.H:40
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:197
A token holds items read from Istream.
Definition: token.H:73
bool isNumber() const
Definition: tokenI.H:745
scalar number() const
Definition: tokenI.H:755
Templated form of IOobject providing type information for file reading and header type checking.
Definition: IOobject.H:530
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Calculate the divergence of the given field.
label patchi
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
tmp< LagrangianEqn< Type > > ddt0(const LagrangianSubScalarField &deltaT, const LagrangianSubSubField< Type > &psi)
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
const FieldField< volMesh::PatchField, Type > & ff(const FieldField< volMesh::PatchField, Type > &bf)
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.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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
const dimensionSet dimVolumetricFlux
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet dimTime
void evaluate(GeometricField< Type, GeoMesh > &result, const Function1< Type > &func, const GeometricField< Type, GeoMesh > &x)
const dimensionSet dimVolume
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:181
IOerror FatalIOError
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
const dimensionSet dimVelocity
error FatalError
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const unitConversion unitFraction
label timeIndex
Definition: getTimeIndex.H:4
labelList fv(nPoints)
dictionary dict