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 template<class Type>
195 template<class GeoField>
196 scalar CrankNicolsonDdtScheme<Type>::coef_
197 (
198  const DDt0Field<GeoField>& ddt0
199 ) const
200 {
201  if (mesh().time().timeIndex() > ddt0.startTimeIndex())
202  {
203  return 1 + ocCoeff();
204  }
205  else
206  {
207  return 1;
208  }
209 }
210 
211 
212 template<class Type>
213 template<class GeoField>
214 scalar CrankNicolsonDdtScheme<Type>::coef0_
215 (
216  const DDt0Field<GeoField>& ddt0
217 ) const
218 {
219  if (mesh().time().timeIndex() > ddt0.startTimeIndex() + 1)
220  {
221  return 1 + ocCoeff();
222  }
223  else
224  {
225  return 1;
226  }
227 }
228 
229 
230 template<class Type>
231 template<class GeoField>
232 dimensionedScalar CrankNicolsonDdtScheme<Type>::rDtCoef_
233 (
234  const DDt0Field<GeoField>& ddt0
235 ) const
236 {
237  return coef_(ddt0)/mesh().time().deltaT();
238 }
239 
240 
241 template<class Type>
242 template<class GeoField>
243 dimensionedScalar CrankNicolsonDdtScheme<Type>::rDtCoef0_
244 (
245  const DDt0Field<GeoField>& ddt0
246 ) const
247 {
248  return coef0_(ddt0)/mesh().time().deltaT0();
249 }
250 
251 
252 template<class Type>
253 template<class GeoField>
254 tmp<GeoField> CrankNicolsonDdtScheme<Type>::offCentre_
255 (
256  const GeoField& ddt0
257 ) const
258 {
259  if (ocCoeff() < 1)
260  {
261  return ocCoeff()*ddt0;
262  }
263  else
264  {
265  return ddt0;
266  }
267 }
268 
269 
270 template<class Type>
272 (
274 )
275 {
276  return bf;
277 }
278 
279 
280 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
281 
282 template<class Type>
284 :
285  ddtScheme<Type>(mesh),
286  ocCoeff_(new Function1s::Constant<scalar>("ocCoeff", 1))
287 {
288  // Ensure the old-old-time cell volumes are available
289  // for moving meshes
290  if (mesh.moving())
291  {
292  mesh.V00();
293  }
294 }
295 
296 
297 template<class Type>
299 (
300  const fvMesh& mesh,
301  Istream& is
302 )
303 :
304  ddtScheme<Type>(mesh, is)
305 {
306  token firstToken(is);
307 
308  if (firstToken.isNumber())
309  {
310  const scalar ocCoeff = firstToken.number();
311 
312  if (ocCoeff < 0 || ocCoeff > 1)
313  {
315  << "Off-centreing coefficient = " << ocCoeff
316  << " should be >= 0 and <= 1"
317  << exit(FatalIOError);
318  }
319 
320  ocCoeff_ =
321  new Function1s::Constant<scalar>("ocCoeff", ocCoeff);
322  }
323  else
324  {
325  is.putBack(firstToken);
326  dictionary dict(is);
327  ocCoeff_ =
329  (
330  "ocCoeff",
331  mesh.time().userUnits(),
332  unitFraction,
333  dict
334  );
335  }
336 
337  // Ensure the old-old-time cell volumes are available
338  // for moving meshes
339  if (mesh.moving())
340  {
341  mesh.V00();
342  }
343 }
344 
345 
346 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
347 
348 template<class Type>
351 (
352  const dimensioned<Type>& dt
353 )
354 {
355  DDt0Field<VolField<Type>>& ddt0 =
356  ddt0_<VolField<Type>>
357  (
358  "ddt0(" + dt.name() + ')',
359  dt.dimensions()
360  );
361 
362  const word ddtName("ddt(" + dt.name() + ')');
363 
364  tmp<VolField<Type>> tdtdt
365  (
367  (
368  ddtName,
369  mesh(),
371  (
372  "0",
373  dt.dimensions()/dimTime,
374  Zero
375  )
376  )
377  );
378 
379  dimensionedScalar rDtCoef = rDtCoef_(ddt0);
380 
381  if (mesh().moving())
382  {
383  if (evaluate(ddt0))
384  {
385  dimensionedScalar rDtCoef0 = rDtCoef0_(ddt0);
386 
387  ddt0.internalFieldRef() =
388  (
389  (rDtCoef0*dt)*(mesh().V0() - mesh().V00())
390  - mesh().V00()*offCentre_(ddt0.internalField())
391  )/mesh().V0();
392  }
393 
394  tdtdt.ref().internalFieldRef() =
395  (
396  (rDtCoef*dt)*(mesh().V() - mesh().V0())
397  - mesh().V0()*offCentre_(ddt0.internalField())
398  )/mesh().V();
399  }
400 
401  return tdtdt;
402 }
403 
404 
405 template<class Type>
408 (
409  const VolField<Type>& vf
410 )
411 {
412  DDt0Field<VolField<Type>>& ddt0 =
413  ddt0_<VolField<Type>>
414  (
415  "ddt0(" + vf.name() + ')',
416  vf.dimensions()
417  );
418 
419  const word ddtName("ddt(" + vf.name() + ')');
420 
421  dimensionedScalar rDtCoef = rDtCoef_(ddt0);
422 
423  if (mesh().moving())
424  {
425  if (evaluate(ddt0))
426  {
427  const scalar rDtCoef0 = rDtCoef0_(ddt0).value();
428 
429  ddt0.primitiveFieldRef() =
430  (
431  rDtCoef0*
432  (
433  mesh().V0()*vf.oldTime().primitiveField()
434  - mesh().V00()*vf.oldTime().oldTime().primitiveField()
435  ) - mesh().V00()*offCentre_(ddt0.primitiveField())
436  )/mesh().V0();
437 
438  ddt0.boundaryFieldRef() =
439  (
440  rDtCoef0*
441  (
442  vf.oldTime().boundaryField()
443  - vf.oldTime().oldTime().boundaryField()
444  ) - offCentre_(ff(ddt0.boundaryField()))
445  );
446  }
447 
448  return VolField<Type>::New
449  (
450  ddtName,
451  (
452  rDtCoef*
453  (
454  mesh().V()*vf
455  - mesh().V0()*vf.oldTime()
456  ) - mesh().V0()*offCentre_(ddt0()())
457  )/mesh().V(),
458  rDtCoef.value()*
459  (
460  vf.boundaryField() - vf.oldTime().boundaryField()
461  ) - offCentre_(ff(ddt0.boundaryField()))
462  );
463  }
464  else
465  {
466  if (evaluate(ddt0))
467  {
468  ddt0 = rDtCoef0_(ddt0)*(vf.oldTime() - vf.oldTime().oldTime())
469  - offCentre_(ddt0());
470  }
471 
472  return VolField<Type>::New
473  (
474  ddtName,
475  rDtCoef*(vf - vf.oldTime()) - offCentre_(ddt0())
476  );
477  }
478 }
479 
480 
481 template<class Type>
484 (
485  const dimensionedScalar& rho,
486  const VolField<Type>& vf
487 )
488 {
489  DDt0Field<VolField<Type>>& ddt0 =
490  ddt0_<VolField<Type>>
491  (
492  "ddt0(" + rho.name() + ',' + vf.name() + ')',
493  rho.dimensions()*vf.dimensions()
494  );
495 
496  const word ddtName("ddt(" + rho.name() + ',' + vf.name() + ')');
497 
498  dimensionedScalar rDtCoef = rDtCoef_(ddt0);
499 
500  if (mesh().moving())
501  {
502  if (evaluate(ddt0))
503  {
504  const scalar rDtCoef0 = rDtCoef0_(ddt0).value();
505 
506  ddt0.primitiveFieldRef() =
507  (
508  rDtCoef0*rho.value()*
509  (
510  mesh().V0()*vf.oldTime().primitiveField()
511  - mesh().V00()*vf.oldTime().oldTime().primitiveField()
512  ) - mesh().V00()*offCentre_(ddt0.primitiveField())
513  )/mesh().V0();
514 
515  ddt0.boundaryFieldRef() =
516  (
517  rDtCoef0*rho.value()*
518  (
519  vf.oldTime().boundaryField()
520  - vf.oldTime().oldTime().boundaryField()
521  ) - offCentre_(ff(ddt0.boundaryField()))
522  );
523  }
524 
525  return VolField<Type>::New
526  (
527  ddtName,
528  (
529  rDtCoef*rho*
530  (
531  mesh().V()*vf()
532  - mesh().V0()*vf.oldTime()()
533  ) - mesh().V0()*offCentre_(ddt0())()
534  )/mesh().V(),
535  rDtCoef.value()*rho.value()*
536  (
537  vf.boundaryField() - vf.oldTime().boundaryField()
538  ) - offCentre_(ff(ddt0.boundaryField()))
539  );
540  }
541  else
542  {
543  if (evaluate(ddt0))
544  {
545  ddt0 = rDtCoef0_(ddt0)*rho*(vf.oldTime() - vf.oldTime().oldTime())
546  - offCentre_(ddt0());
547  }
548 
549  return VolField<Type>::New
550  (
551  ddtName,
552  rDtCoef*rho*(vf - vf.oldTime()) - offCentre_(ddt0())
553  );
554  }
555 }
556 
557 
558 template<class Type>
561 (
562  const volScalarField& rho,
563  const VolField<Type>& vf
564 )
565 {
566  DDt0Field<VolField<Type>>& ddt0 =
567  ddt0_<VolField<Type>>
568  (
569  "ddt0(" + rho.name() + ',' + vf.name() + ')',
570  rho.dimensions()*vf.dimensions()
571  );
572 
573  const word ddtName("ddt(" + rho.name() + ',' + vf.name() + ')');
574 
575  dimensionedScalar rDtCoef = rDtCoef_(ddt0);
576 
577  if (mesh().moving())
578  {
579  if (evaluate(ddt0))
580  {
581  const scalar rDtCoef0 = rDtCoef0_(ddt0).value();
582 
583  ddt0.primitiveFieldRef() =
584  (
585  rDtCoef0*
586  (
587  mesh().V0()*rho.oldTime().primitiveField()
588  *vf.oldTime().primitiveField()
589  - mesh().V00()*rho.oldTime().oldTime().primitiveField()
590  *vf.oldTime().oldTime().primitiveField()
591  ) - mesh().V00()*offCentre_(ddt0.primitiveField())
592  )/mesh().V0();
593 
594  ddt0.boundaryFieldRef() =
595  (
596  rDtCoef0*
597  (
598  rho.oldTime().boundaryField()
599  *vf.oldTime().boundaryField()
600  - rho.oldTime().oldTime().boundaryField()
601  *vf.oldTime().oldTime().boundaryField()
602  ) - offCentre_(ff(ddt0.boundaryField()))
603  );
604  }
605 
606  return VolField<Type>::New
607  (
608  ddtName,
609  (
610  rDtCoef*
611  (
612  mesh().V()*rho()*vf()
613  - mesh().V0()*rho.oldTime()()
614  *vf.oldTime()()
615  ) - mesh().V00()*offCentre_(ddt0())()
616  )/mesh().V(),
617  rDtCoef.value()*
618  (
619  rho.boundaryField()*vf.boundaryField()
620  - rho.oldTime().boundaryField()*vf.oldTime().boundaryField()
621  ) - offCentre_(ff(ddt0.boundaryField()))
622  );
623  }
624  else
625  {
626  if (evaluate(ddt0))
627  {
628  ddt0 = rDtCoef0_(ddt0)*
629  (
630  rho.oldTime()*vf.oldTime()
631  - rho.oldTime().oldTime()*vf.oldTime().oldTime()
632  ) - offCentre_(ddt0());
633  }
634 
635  return VolField<Type>::New
636  (
637  ddtName,
638  rDtCoef*(rho*vf - rho.oldTime()*vf.oldTime()) - offCentre_(ddt0())
639  );
640  }
641 }
642 
643 
644 template<class Type>
647 (
648  const volScalarField& alpha,
649  const volScalarField& rho,
650  const VolField<Type>& vf
651 )
652 {
653  DDt0Field<VolField<Type>>& ddt0 =
654  ddt0_<VolField<Type>>
655  (
656  "ddt0(" + alpha.name() + ',' + rho.name() + ',' + vf.name() + ')',
657  alpha.dimensions()*rho.dimensions()*vf.dimensions()
658  );
659 
660  const word ddtName
661  (
662  "ddt(" + alpha.name() + ',' + rho.name() + ',' + vf.name() + ')'
663  );
664 
665  dimensionedScalar rDtCoef = rDtCoef_(ddt0);
666 
667  if (mesh().moving())
668  {
669  if (evaluate(ddt0))
670  {
671  const scalar rDtCoef0 = rDtCoef0_(ddt0).value();
672 
673  ddt0.primitiveFieldRef() =
674  (
675  rDtCoef0*
676  (
677  mesh().V0()
678  *alpha.oldTime().primitiveField()
679  *rho.oldTime().primitiveField()
680  *vf.oldTime().primitiveField()
681 
682  - mesh().V00()
683  *alpha.oldTime().oldTime().primitiveField()
684  *rho.oldTime().oldTime().primitiveField()
685  *vf.oldTime().oldTime().primitiveField()
686  ) - mesh().V00()*offCentre_(ddt0.primitiveField())
687  )/mesh().V0();
688 
689  ddt0.boundaryFieldRef() =
690  (
691  rDtCoef0*
692  (
693  alpha.oldTime().boundaryField()
694  *rho.oldTime().boundaryField()
695  *vf.oldTime().boundaryField()
696 
697  - alpha.oldTime().oldTime().boundaryField()
698  *rho.oldTime().oldTime().boundaryField()
699  *vf.oldTime().oldTime().boundaryField()
700  ) - offCentre_(ff(ddt0.boundaryField()))
701  );
702  }
703 
704  return VolField<Type>::New
705  (
706  ddtName,
707  (
708  rDtCoef.value()*
709  (
710  mesh().V()*alpha()*rho()*vf()
711  - mesh().V0()*alpha.oldTime()()*rho.oldTime()()*vf.oldTime()()
712  ) - mesh().V00()*offCentre_(ddt0())()
713  )/mesh().V(),
714  rDtCoef.value()*
715  (
716  alpha.boundaryField()
717  *rho.boundaryField()
718  *vf.boundaryField()
719 
720  - alpha.oldTime().boundaryField()
721  *rho.oldTime().boundaryField()
722  *vf.oldTime().boundaryField()
723  ) - offCentre_(ff(ddt0.boundaryField()))
724  );
725  }
726  else
727  {
728  if (evaluate(ddt0))
729  {
730  ddt0 = rDtCoef0_(ddt0)*
731  (
732  alpha.oldTime()
733  *rho.oldTime()
734  *vf.oldTime()
735 
736  - alpha.oldTime().oldTime()
737  *rho.oldTime().oldTime()
738  *vf.oldTime().oldTime()
739  ) - offCentre_(ddt0());
740  }
741 
742  return VolField<Type>::New
743  (
744  ddtName,
745  rDtCoef
746  *(
747  alpha*rho*vf
748  - alpha.oldTime()*rho.oldTime()*vf.oldTime()
749  )
750  - offCentre_(ddt0())
751  );
752  }
753 }
754 
755 
756 template<class Type>
759 (
760  const VolField<Type>& vf
761 )
762 {
763  DDt0Field<VolField<Type>>& ddt0 =
764  ddt0_<VolField<Type>>
765  (
766  "ddt0(" + vf.name() + ')',
767  vf.dimensions()
768  );
769 
770  tmp<fvMatrix<Type>> tfvm
771  (
772  new fvMatrix<Type>
773  (
774  vf,
776  )
777  );
778 
779  fvMatrix<Type>& fvm = tfvm.ref();
780 
781  const scalar rDtCoef = rDtCoef_(ddt0).value();
782  fvm.diag() = rDtCoef*mesh().V();
783 
784  vf.oldTime().oldTime();
785 
786  if (mesh().moving())
787  {
788  if (evaluate(ddt0))
789  {
790  const scalar rDtCoef0 = rDtCoef0_(ddt0).value();
791 
792  ddt0.primitiveFieldRef() =
793  (
794  rDtCoef0*
795  (
796  mesh().V0()*vf.oldTime().primitiveField()
797  - mesh().V00()*vf.oldTime().oldTime().primitiveField()
798  )
799  - mesh().V00()*offCentre_(ddt0.primitiveField())
800  )/mesh().V0();
801 
802  ddt0.boundaryFieldRef() =
803  (
804  rDtCoef0*
805  (
806  vf.oldTime().boundaryField()
807  - vf.oldTime().oldTime().boundaryField()
808  )
809  - offCentre_(ff(ddt0.boundaryField()))
810  );
811  }
812 
813  fvm.source() =
814  (
815  rDtCoef*vf.oldTime().primitiveField()
816  + offCentre_(ddt0.primitiveField())
817  )*mesh().V0();
818  }
819  else
820  {
821  if (evaluate(ddt0))
822  {
823  ddt0 = rDtCoef0_(ddt0)*(vf.oldTime() - vf.oldTime().oldTime())
824  - offCentre_(ddt0());
825  }
826 
827  fvm.source() =
828  (
829  rDtCoef*vf.oldTime().primitiveField()
830  + offCentre_(ddt0.primitiveField())
831  )*mesh().V();
832  }
833 
834  return tfvm;
835 }
836 
837 
838 template<class Type>
841 (
842  const dimensionedScalar& rho,
843  const VolField<Type>& vf
844 )
845 {
846  DDt0Field<VolField<Type>>& ddt0 =
847  ddt0_<VolField<Type>>
848  (
849  "ddt0(" + rho.name() + ',' + vf.name() + ')',
850  rho.dimensions()*vf.dimensions()
851  );
852 
853  tmp<fvMatrix<Type>> tfvm
854  (
855  new fvMatrix<Type>
856  (
857  vf,
858  rho.dimensions()*vf.dimensions()*dimVolume/dimTime
859  )
860  );
861  fvMatrix<Type>& fvm = tfvm.ref();
862 
863  const scalar rDtCoef = rDtCoef_(ddt0).value();
864  fvm.diag() = rDtCoef*rho.value()*mesh().V();
865 
866  vf.oldTime().oldTime();
867 
868  if (mesh().moving())
869  {
870  if (evaluate(ddt0))
871  {
872  const scalar rDtCoef0 = rDtCoef0_(ddt0).value();
873 
874  ddt0.primitiveFieldRef() =
875  (
876  rDtCoef0*rho.value()*
877  (
878  mesh().V0()*vf.oldTime().primitiveField()
879  - mesh().V00()*vf.oldTime().oldTime().primitiveField()
880  )
881  - mesh().V00()*offCentre_(ddt0.primitiveField())
882  )/mesh().V0();
883 
884  ddt0.boundaryFieldRef() =
885  (
886  rDtCoef0*rho.value()*
887  (
888  vf.oldTime().boundaryField()
889  - vf.oldTime().oldTime().boundaryField()
890  )
891  - offCentre_(ff(ddt0.boundaryField()))
892  );
893  }
894 
895  fvm.source() =
896  (
897  rDtCoef*rho.value()*vf.oldTime().primitiveField()
898  + offCentre_(ddt0.primitiveField())
899  )*mesh().V0();
900  }
901  else
902  {
903  if (evaluate(ddt0))
904  {
905  ddt0 = rDtCoef0_(ddt0)*rho*(vf.oldTime() - vf.oldTime().oldTime())
906  - offCentre_(ddt0());
907  }
908 
909  fvm.source() =
910  (
911  rDtCoef*rho.value()*vf.oldTime().primitiveField()
912  + offCentre_(ddt0.primitiveField())
913  )*mesh().V();
914  }
915 
916  return tfvm;
917 }
918 
919 
920 template<class Type>
923 (
924  const volScalarField& rho,
925  const VolField<Type>& vf
926 )
927 {
928  DDt0Field<VolField<Type>>& ddt0 =
929  ddt0_<VolField<Type>>
930  (
931  "ddt0(" + rho.name() + ',' + vf.name() + ')',
932  rho.dimensions()*vf.dimensions()
933  );
934 
935  tmp<fvMatrix<Type>> tfvm
936  (
937  new fvMatrix<Type>
938  (
939  vf,
940  rho.dimensions()*vf.dimensions()*dimVolume/dimTime
941  )
942  );
943  fvMatrix<Type>& fvm = tfvm.ref();
944 
945  const scalar rDtCoef = rDtCoef_(ddt0).value();
946  fvm.diag() = rDtCoef*rho.primitiveField()*mesh().V();
947 
948  vf.oldTime().oldTime();
949  rho.oldTime().oldTime();
950 
951  if (mesh().moving())
952  {
953  if (evaluate(ddt0))
954  {
955  const scalar rDtCoef0 = rDtCoef0_(ddt0).value();
956 
957  ddt0.primitiveFieldRef() =
958  (
959  rDtCoef0*
960  (
961  mesh().V0()*rho.oldTime().primitiveField()
962  *vf.oldTime().primitiveField()
963  - mesh().V00()*rho.oldTime().oldTime().primitiveField()
964  *vf.oldTime().oldTime().primitiveField()
965  )
966  - mesh().V00()*offCentre_(ddt0.primitiveField())
967  )/mesh().V0();
968 
969  ddt0.boundaryFieldRef() =
970  (
971  rDtCoef0*
972  (
973  rho.oldTime().boundaryField()
974  *vf.oldTime().boundaryField()
975  - rho.oldTime().oldTime().boundaryField()
976  *vf.oldTime().oldTime().boundaryField()
977  )
978  - offCentre_(ff(ddt0.boundaryField()))
979  );
980  }
981 
982  fvm.source() =
983  (
984  rDtCoef*rho.oldTime().primitiveField()*vf.oldTime().primitiveField()
985  + offCentre_(ddt0.primitiveField())
986  )*mesh().V0();
987  }
988  else
989  {
990  if (evaluate(ddt0))
991  {
992  ddt0 = rDtCoef0_(ddt0)*
993  (
994  rho.oldTime()*vf.oldTime()
995  - rho.oldTime().oldTime()*vf.oldTime().oldTime()
996  ) - offCentre_(ddt0());
997  }
998 
999  fvm.source() =
1000  (
1001  rDtCoef*rho.oldTime().primitiveField()*vf.oldTime().primitiveField()
1002  + offCentre_(ddt0.primitiveField())
1003  )*mesh().V();
1004  }
1005 
1006  return tfvm;
1007 }
1008 
1009 
1010 template<class Type>
1013 (
1014  const volScalarField& alpha,
1015  const volScalarField& rho,
1016  const VolField<Type>& vf
1017 )
1018 {
1019  DDt0Field<VolField<Type>>& ddt0 =
1020  ddt0_<VolField<Type>>
1021  (
1022  "ddt0(" + alpha.name() + ',' + rho.name() + ',' + vf.name() + ')',
1023  alpha.dimensions()*rho.dimensions()*vf.dimensions()
1024  );
1025 
1026  tmp<fvMatrix<Type>> tfvm
1027  (
1028  new fvMatrix<Type>
1029  (
1030  vf,
1031  alpha.dimensions()
1032  *rho.dimensions()
1033  *vf.dimensions()
1034  *dimVolume
1035  /dimTime
1036  )
1037  );
1038  fvMatrix<Type>& fvm = tfvm.ref();
1039 
1040  const scalar rDtCoef = rDtCoef_(ddt0).value();
1041  fvm.diag() = rDtCoef*alpha.primitiveField()*rho.primitiveField()*mesh().V();
1042 
1043  vf.oldTime().oldTime();
1044  alpha.oldTime().oldTime();
1045  rho.oldTime().oldTime();
1046 
1047  if (mesh().moving())
1048  {
1049  if (evaluate(ddt0))
1050  {
1051  const scalar rDtCoef0 = rDtCoef0_(ddt0).value();
1052 
1053  ddt0.primitiveFieldRef() =
1054  (
1055  rDtCoef0*
1056  (
1057  mesh().V0()
1058  *alpha.oldTime().primitiveField()
1059  *rho.oldTime().primitiveField()
1060  *vf.oldTime().primitiveField()
1061 
1062  - mesh().V00()
1063  *alpha.oldTime().oldTime().primitiveField()
1064  *rho.oldTime().oldTime().primitiveField()
1065  *vf.oldTime().oldTime().primitiveField()
1066  )
1067  - mesh().V00()*offCentre_(ddt0.primitiveField())
1068  )/mesh().V0();
1069 
1070  ddt0.boundaryFieldRef() =
1071  (
1072  rDtCoef0*
1073  (
1074  alpha.oldTime().boundaryField()
1075  *rho.oldTime().boundaryField()
1076  *vf.oldTime().boundaryField()
1077 
1078  - alpha.oldTime().oldTime().boundaryField()
1079  *rho.oldTime().oldTime().boundaryField()
1080  *vf.oldTime().oldTime().boundaryField()
1081  )
1082  - offCentre_(ff(ddt0.boundaryField()))
1083  );
1084  }
1085 
1086  fvm.source() =
1087  (
1088  rDtCoef
1089  *alpha.oldTime().primitiveField()
1090  *rho.oldTime().primitiveField()
1091  *vf.oldTime().primitiveField()
1092  + offCentre_(ddt0.primitiveField())
1093  )*mesh().V0();
1094  }
1095  else
1096  {
1097  if (evaluate(ddt0))
1098  {
1099  ddt0 = rDtCoef0_(ddt0)*
1100  (
1101  alpha.oldTime()
1102  *rho.oldTime()
1103  *vf.oldTime()
1104 
1105  - alpha.oldTime().oldTime()
1106  *rho.oldTime().oldTime()
1107  *vf.oldTime().oldTime()
1108  ) - offCentre_(ddt0());
1109  }
1110 
1111  fvm.source() =
1112  (
1113  rDtCoef
1114  *alpha.oldTime().primitiveField()
1115  *rho.oldTime().primitiveField()
1116  *vf.oldTime().primitiveField()
1117  + offCentre_(ddt0.primitiveField())
1118  )*mesh().V();
1119  }
1120 
1121  return tfvm;
1122 }
1123 
1124 
1125 template<class Type>
1128 (
1129  const VolField<Type>& U,
1130  const SurfaceField<Type>& Uf
1131 )
1132 {
1133  DDt0Field<VolField<Type>>& ddt0 =
1134  ddt0_<VolField<Type>>
1135  (
1136  "ddtCorrDdt0(" + U.name() + ')',
1137  U.dimensions()
1138  );
1139 
1140  DDt0Field<SurfaceField<Type>>& dUfdt0 =
1141  ddt0_<SurfaceField<Type>>
1142  (
1143  "ddtCorrDdt0(" + Uf.name() + ')',
1144  Uf.dimensions()
1145  );
1146 
1147  dimensionedScalar rDtCoef = rDtCoef_(ddt0);
1148 
1149  if (evaluate(ddt0))
1150  {
1151  ddt0 =
1152  rDtCoef0_(ddt0)*(U.oldTime() - U.oldTime().oldTime())
1153  - offCentre_(ddt0());
1154  }
1155 
1156  if (evaluate(dUfdt0))
1157  {
1158  dUfdt0 =
1159  rDtCoef0_(dUfdt0)*(Uf.oldTime() - Uf.oldTime().oldTime())
1160  - offCentre_(dUfdt0());
1161  }
1162 
1163  return fluxFieldType::New
1164  (
1165  "ddtCorr(" + U.name() + ',' + Uf.name() + ')',
1166  this->fvcDdtPhiCoeff(U.oldTime(), mesh().Sf() & Uf.oldTime())
1167  *(
1168  mesh().Sf()
1169  & (
1170  (rDtCoef*Uf.oldTime() + offCentre_(dUfdt0()))
1171  - fvc::interpolate(rDtCoef*U.oldTime() + offCentre_(ddt0()))
1172  )
1173  )
1174  );
1175 }
1176 
1177 
1178 template<class Type>
1181 (
1182  const VolField<Type>& U,
1183  const fluxFieldType& phi
1184 )
1185 {
1186  DDt0Field<VolField<Type>>& ddt0 =
1187  ddt0_<VolField<Type>>
1188  (
1189  "ddtCorrDdt0(" + U.name() + ')',
1190  U.dimensions()
1191  );
1192 
1193  DDt0Field<fluxFieldType>& dphidt0 =
1194  ddt0_<fluxFieldType>
1195  (
1196  "ddtCorrDdt0(" + phi.name() + ')',
1197  phi.dimensions()
1198  );
1199 
1200  dimensionedScalar rDtCoef = rDtCoef_(ddt0);
1201 
1202  if (evaluate(ddt0))
1203  {
1204  ddt0 =
1205  rDtCoef0_(ddt0)*(U.oldTime() - U.oldTime().oldTime())
1206  - offCentre_(ddt0());
1207  }
1208 
1209  if (evaluate(dphidt0))
1210  {
1211  dphidt0 =
1212  rDtCoef0_(dphidt0)*(phi.oldTime() - phi.oldTime().oldTime())
1213  - offCentre_(dphidt0());
1214  }
1215 
1216  return fluxFieldType::New
1217  (
1218  "ddtCorr(" + U.name() + ',' + phi.name() + ')',
1219  this->fvcDdtPhiCoeff(U.oldTime(), phi.oldTime())
1220  *(
1221  (rDtCoef*phi.oldTime() + offCentre_(dphidt0()))
1223  (
1224  mesh().Sf(),
1225  rDtCoef*U.oldTime() + offCentre_(ddt0())
1226  )
1227  )
1228  );
1229 }
1230 
1231 
1232 template<class Type>
1235 (
1236  const volScalarField& rho,
1237  const VolField<Type>& U,
1238  const SurfaceField<Type>& rhoUf
1239 )
1240 {
1241  if
1242  (
1243  U.dimensions() == dimVelocity
1244  && rhoUf.dimensions() == rho.dimensions()*dimVelocity
1245  )
1246  {
1247  DDt0Field<VolField<Type>>& ddt0 =
1248  ddt0_<VolField<Type>>
1249  (
1250  "ddtCorrDdt0(" + rho.name() + ',' + U.name() + ')',
1251  rho.dimensions()*U.dimensions()
1252  );
1253 
1254  DDt0Field<SurfaceField<Type>>& drhoUfdt0 =
1255  ddt0_<SurfaceField<Type>>
1256  (
1257  "ddtCorrDdt0(" + rhoUf.name() + ')',
1258  rhoUf.dimensions()
1259  );
1260 
1261  dimensionedScalar rDtCoef = rDtCoef_(ddt0);
1262 
1263  VolField<Type> rhoU0
1264  (
1265  rho.oldTime()*U.oldTime()
1266  );
1267 
1268  if (evaluate(ddt0))
1269  {
1270  ddt0 =
1271  rDtCoef0_(ddt0)
1272  *(rhoU0 - rho.oldTime().oldTime()*U.oldTime().oldTime())
1273  - offCentre_(ddt0());
1274  }
1275 
1276  if (evaluate(drhoUfdt0))
1277  {
1278  drhoUfdt0 =
1279  rDtCoef0_(drhoUfdt0)
1280  *(rhoUf.oldTime() - rhoUf.oldTime().oldTime())
1281  - offCentre_(drhoUfdt0());
1282  }
1283 
1284  return fluxFieldType::New
1285  (
1286  "ddtCorr(" + rho.name() + ',' + U.name() + ',' + rhoUf.name() + ')',
1287  this->fvcDdtPhiCoeff
1288  (
1289  rhoU0,
1290  mesh().Sf() & rhoUf.oldTime(),
1291  rho.oldTime()
1292  )
1293  *(
1294  mesh().Sf()
1295  & (
1296  (rDtCoef*rhoUf.oldTime() + offCentre_(drhoUfdt0()))
1297  - fvc::interpolate(rDtCoef*rhoU0 + offCentre_(ddt0()))
1298  )
1299  )
1300  );
1301  }
1302  else if
1303  (
1304  U.dimensions() == rho.dimensions()*dimVelocity
1305  && rhoUf.dimensions() == rho.dimensions()*dimVelocity
1306  )
1307  {
1308  DDt0Field<VolField<Type>>& ddt0 =
1309  ddt0_<VolField<Type>>
1310  (
1311  "ddtCorrDdt0(" + U.name() + ')',
1312  U.dimensions()
1313  );
1314 
1315  DDt0Field<SurfaceField<Type>>& drhoUfdt0 =
1316  ddt0_<SurfaceField<Type>>
1317  (
1318  "ddtCorrDdt0(" + rhoUf.name() + ')',
1319  rhoUf.dimensions()
1320  );
1321 
1322  dimensionedScalar rDtCoef = rDtCoef_(ddt0);
1323 
1324  if (evaluate(ddt0))
1325  {
1326  ddt0 =
1327  rDtCoef0_(ddt0)*(U.oldTime() - U.oldTime().oldTime())
1328  - offCentre_(ddt0());
1329  }
1330 
1331  if (evaluate(drhoUfdt0))
1332  {
1333  drhoUfdt0 =
1334  rDtCoef0_(drhoUfdt0)
1335  *(rhoUf.oldTime() - rhoUf.oldTime().oldTime())
1336  - offCentre_(drhoUfdt0());
1337  }
1338 
1339  return fluxFieldType::New
1340  (
1341  "ddtCorr(" + U.name() + ',' + rhoUf.name() + ')',
1342  this->fvcDdtPhiCoeff
1343  (
1344  U.oldTime(),
1345  mesh().Sf() & rhoUf.oldTime(),
1346  rho.oldTime()
1347  )
1348  *(
1349  mesh().Sf()
1350  & (
1351  (rDtCoef*rhoUf.oldTime() + offCentre_(drhoUfdt0()))
1353  (
1354  rDtCoef*U.oldTime() + offCentre_(ddt0())
1355  )
1356  )
1357  )
1358  );
1359  }
1360  else
1361  {
1363  << "dimensions of rhoUf are not correct"
1364  << abort(FatalError);
1365 
1366  return fluxFieldType::null();
1367  }
1368 }
1369 
1370 
1371 template<class Type>
1374 (
1375  const volScalarField& rho,
1376  const VolField<Type>& U,
1377  const fluxFieldType& phi
1378 )
1379 {
1380  if
1381  (
1382  U.dimensions() == dimVelocity
1383  && phi.dimensions() == rho.dimensions()*dimVolumetricFlux
1384  )
1385  {
1386  DDt0Field<VolField<Type>>& ddt0 =
1387  ddt0_<VolField<Type>>
1388  (
1389  "ddtCorrDdt0(" + rho.name() + ',' + U.name() + ')',
1390  rho.dimensions()*U.dimensions()
1391  );
1392 
1393  DDt0Field<fluxFieldType>& dphidt0 =
1394  ddt0_<fluxFieldType>
1395  (
1396  "ddtCorrDdt0(" + phi.name() + ')',
1397  phi.dimensions()
1398  );
1399 
1400  dimensionedScalar rDtCoef = rDtCoef_(ddt0);
1401 
1402  VolField<Type> rhoU0
1403  (
1404  rho.oldTime()*U.oldTime()
1405  );
1406 
1407  if (evaluate(ddt0))
1408  {
1409  ddt0 =
1410  rDtCoef0_(ddt0)
1411  *(rhoU0 - rho.oldTime().oldTime()*U.oldTime().oldTime())
1412  - offCentre_(ddt0());
1413  }
1414 
1415  if (evaluate(dphidt0))
1416  {
1417  dphidt0 =
1418  rDtCoef0_(dphidt0)
1419  *(phi.oldTime() - phi.oldTime().oldTime())
1420  - offCentre_(dphidt0());
1421  }
1422 
1423  return fluxFieldType::New
1424  (
1425  "ddtCorr(" + rho.name() + ',' + U.name() + ',' + phi.name() + ')',
1426  this->fvcDdtPhiCoeff(rhoU0, phi.oldTime(), rho.oldTime())
1427  *(
1428  (rDtCoef*phi.oldTime() + offCentre_(dphidt0()))
1430  (
1431  mesh().Sf(),
1432  rDtCoef*rhoU0 + offCentre_(ddt0())
1433  )
1434  )
1435  );
1436  }
1437  else if
1438  (
1439  U.dimensions() == rho.dimensions()*dimVelocity
1440  && phi.dimensions() == rho.dimensions()*dimVolumetricFlux
1441  )
1442  {
1443  DDt0Field<VolField<Type>>& ddt0 =
1444  ddt0_<VolField<Type>>
1445  (
1446  "ddtCorrDdt0(" + U.name() + ')',
1447  U.dimensions()
1448  );
1449 
1450  DDt0Field<fluxFieldType>& dphidt0 =
1451  ddt0_<fluxFieldType>
1452  (
1453  "ddtCorrDdt0(" + phi.name() + ')',
1454  phi.dimensions()
1455  );
1456 
1457  dimensionedScalar rDtCoef = rDtCoef_(ddt0);
1458 
1459  if (evaluate(ddt0))
1460  {
1461  ddt0 =
1462  rDtCoef0_(ddt0)*(U.oldTime() - U.oldTime().oldTime())
1463  - offCentre_(ddt0());
1464  }
1465 
1466  if (evaluate(dphidt0))
1467  {
1468  dphidt0 =
1469  rDtCoef0_(dphidt0)*(phi.oldTime() - phi.oldTime().oldTime())
1470  - offCentre_(dphidt0());
1471  }
1472 
1473  return fluxFieldType::New
1474  (
1475  "ddtCorr(" + U.name() + ',' + phi.name() + ')',
1476  this->fvcDdtPhiCoeff(U.oldTime(), phi.oldTime(), rho.oldTime())
1477  *(
1478  (rDtCoef*phi.oldTime() + offCentre_(dphidt0()))
1480  (
1481  mesh().Sf(),
1482  rDtCoef*U.oldTime() + offCentre_(ddt0())
1483  )
1484  )
1485  );
1486  }
1487  else
1488  {
1490  << "dimensions of phi are not correct"
1491  << abort(FatalError);
1492 
1493  return fluxFieldType::null();
1494  }
1495 }
1496 
1497 
1498 template<class Type>
1500 (
1501  const VolField<Type>& vf
1502 )
1503 {
1504  DDt0Field<surfaceScalarField>& meshPhi0 =
1505  ddt0_<surfaceScalarField>("meshPhi0", dimVolume);
1506 
1507  if (evaluate(meshPhi0))
1508  {
1509  meshPhi0 =
1510  coef0_(meshPhi0)*mesh().phi().oldTime() - offCentre_(meshPhi0());
1511  }
1512 
1514  (
1515  mesh().phi().name(),
1516  coef_(meshPhi0)*mesh().phi() - offCentre_(meshPhi0())
1517  );
1518 }
1519 
1520 
1521 template<class Type>
1523 (
1524  const VolField<Type>& vf,
1525  const label patchi
1526 )
1527 {
1528  DDt0Field<surfaceScalarField>& meshPhi0 =
1529  ddt0_<surfaceScalarField>("meshPhi0", dimVolume);
1530 
1531  if (evaluate(meshPhi0))
1532  {
1533  meshPhi0 =
1534  coef0_(meshPhi0)*mesh().phi().oldTime() - offCentre_(meshPhi0());
1535  }
1536 
1537  return
1538  (
1539  coef_(meshPhi0)*mesh().phi().boundaryField()[patchi]
1540  - offCentre_
1541  (
1542  static_cast<const scalarField&>(meshPhi0().boundaryField()[patchi])
1543  )
1544  );
1545 }
1546 
1547 
1548 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1549 
1550 } // End namespace fv
1551 
1552 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1553 
1554 } // End namespace Foam
1555 
1556 // ************************************************************************* //
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, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, 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:310
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 FieldType & oldTime() const
Return the old time field.
Definition: OldTimeField.C:315
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:28
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:869
virtual dimensionedScalar startTime() const
Return start time.
Definition: Time.C:802
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
Definition: Time.C:641
virtual label startTimeIndex() const
Return start time index.
Definition: Time.C:785
A list of keyword definitions, which are a keyword followed by any number of values (e....
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:99
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:418
const DimensionedField< scalar, volMesh > & V00() const
Return old-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
const Time & time() const
Return time.
bool moving() const
Is mesh moving.
Definition: polyMesh.H:476
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:181
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:531
A class for handling words, derived from string.
Definition: word.H:62
#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))
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
const FieldField< fvPatchField, Type > & ff(const FieldField< fvPatchField, 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
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimVolumetricFlux
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet dimTime
const dimensionSet dimVolume
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:176
IOerror FatalIOError
const dimensionSet dimVelocity
void evaluate(GeometricField< Type, PatchField, GeoMesh > &result, const Function1< Type > &func, const GeometricField< Type, PatchField, GeoMesh > &x)
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