EulerFaDdtScheme.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2016-2017 Wikki Ltd
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "EulerFaDdtScheme.H"
29 #include "faMatrices.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace fa
39 {
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 template<class Type>
46 (
47  const dimensioned<Type> dt
48 )
49 {
50  dimensionedScalar rDeltaT = 1.0/mesh().time().deltaT();
51 
52  const IOobject ddtIOobject
53  (
54  mesh().thisDb().newIOobject
55  (
56  "ddt("+dt.name()+')',
58  )
59  );
60 
61  if (mesh().moving())
62  {
64  (
66  (
67  ddtIOobject,
68  mesh(),
70  )
71  );
72 
73  tdtdt.ref().primitiveFieldRef() =
74  rDeltaT.value()*dt.value()*(1.0 - mesh().S0()/mesh().S());
75 
76  return tdtdt;
77  }
78 
79 
81  (
82  ddtIOobject,
83  mesh(),
86  );
87 }
88 
89 
90 template<class Type>
93 (
94  const dimensioned<Type> dt
95 )
96 {
97  dimensionedScalar rDeltaT = 1.0/mesh().time().deltaT();
98 
99  const IOobject ddtIOobject
100  (
101  mesh().thisDb().newIOobject
102  (
103  "ddt("+dt.name()+')',
105  )
106  );
107 
109  (
111  (
112  ddtIOobject,
113  mesh(),
114  -rDeltaT*dt
115  )
116  );
117 
118  if (mesh().moving())
119  {
120  tdtdt0.ref().primitiveFieldRef() =
121  (-rDeltaT.value()*dt.value())*mesh().S0()/mesh().S();
122  }
123 
124  return tdtdt0;
125 }
126 
127 
128 template<class Type>
131 (
133 )
134 {
135  dimensionedScalar rDeltaT = 1.0/mesh().time().deltaT();
136 
137  const IOobject ddtIOobject
138  (
139  mesh().thisDb().newIOobject
140  (
141  "ddt("+vf.name()+')',
143  )
144  );
145 
146  if (mesh().moving())
147  {
149  (
151  (
152  ddtIOobject,
153  mesh(),
154  rDeltaT.dimensions()*vf.dimensions(),
155  rDeltaT.value()*
156  (
157  vf()
158  - vf.oldTime()()*mesh().S0()/mesh().S()
159  ),
160  rDeltaT.value()*
161  (
162  vf.boundaryField() - vf.oldTime().boundaryField()
163  )
164  )
165  );
166  }
167  else
168  {
169  return tmp<GeometricField<Type, faPatchField, areaMesh>>
170  (
171  new GeometricField<Type, faPatchField, areaMesh>
172  (
173  ddtIOobject,
174  rDeltaT*(vf - vf.oldTime())
175  )
176  );
177  }
178 }
179 
180 
181 template<class Type>
184 (
186 )
187 {
188  dimensionedScalar rDeltaT = 1.0/mesh().time().deltaT();
189 
190  const IOobject ddtIOobject
191  (
192  mesh().thisDb().newIOobject
193  (
194  "ddt0("+vf.name()+')',
196  )
197  );
198 
199  if (mesh().moving())
200  {
201  return tmp<GeometricField<Type, faPatchField, areaMesh>>
202  (
203  new GeometricField<Type, faPatchField, areaMesh>
204  (
205  ddtIOobject,
206  mesh(),
207  rDeltaT.dimensions()*vf.dimensions(),
208  (-rDeltaT.value())*vf.oldTime().internalField(),
209  (-rDeltaT.value())*vf.oldTime().boundaryField()
210  )
211  );
212  }
213  else
214  {
215  return tmp<GeometricField<Type, faPatchField, areaMesh>>
216  (
217  new GeometricField<Type, faPatchField, areaMesh>
218  (
219  ddtIOobject,
220  (-rDeltaT)*vf.oldTime()
221  )
222  );
223  }
224 }
225 
226 
227 template<class Type>
230 (
231  const dimensionedScalar& rho,
233 )
234 {
235  dimensionedScalar rDeltaT = 1.0/mesh().time().deltaT();
236 
237  const IOobject ddtIOobject
238  (
239  mesh().thisDb().newIOobject
240  (
241  "ddt("+rho.name()+','+vf.name()+')',
243  )
244  );
245 
246  if (mesh().moving())
247  {
248  return tmp<GeometricField<Type, faPatchField, areaMesh>>
249  (
250  new GeometricField<Type, faPatchField, areaMesh>
251  (
252  ddtIOobject,
253  mesh(),
254  rDeltaT.dimensions()*rho.dimensions()*vf.dimensions(),
255  rDeltaT.value()*rho.value()*
256  (
257  vf()
258  - vf.oldTime()()*mesh().S0()/mesh().S()
259  ),
260  rDeltaT.value()*rho.value()*
261  (
262  vf.boundaryField() - vf.oldTime().boundaryField()
263  )
264  )
265  );
266  }
267  else
268  {
269  return tmp<GeometricField<Type, faPatchField, areaMesh>>
270  (
271  new GeometricField<Type, faPatchField, areaMesh>
272  (
273  ddtIOobject,
274  rDeltaT*rho*(vf - vf.oldTime())
275  )
276  );
277  }
278 }
279 
280 template<class Type>
283 (
284  const dimensionedScalar& rho,
286 )
287 {
288  dimensionedScalar rDeltaT = 1.0/mesh().time().deltaT();
289 
290  const IOobject ddtIOobject
291  (
292  mesh().thisDb().newIOobject
293  (
294  "ddt0("+rho.name()+','+vf.name()+')',
296  )
297  );
298 
299  if (mesh().moving())
300  {
301  return tmp<GeometricField<Type, faPatchField, areaMesh>>
302  (
303  new GeometricField<Type, faPatchField, areaMesh>
304  (
305  ddtIOobject,
306  mesh(),
307  rDeltaT.dimensions()*rho.dimensions()*vf.dimensions(),
308  (-rDeltaT.value())*rho.value()*
309  vf.oldTime()()*mesh().S0()/mesh().S(),
310  (-rDeltaT.value())*rho.value()*
311  vf.oldTime().boundaryField()
312  )
313  );
314  }
315  else
316  {
317  return tmp<GeometricField<Type, faPatchField, areaMesh>>
318  (
319  new GeometricField<Type, faPatchField, areaMesh>
320  (
321  ddtIOobject,
322  (-rDeltaT)*rho*vf.oldTime()
323  )
324  );
325  }
326 }
327 
328 
329 template<class Type>
332 (
333  const areaScalarField& rho,
335 )
336 {
337  dimensionedScalar rDeltaT = 1.0/mesh().time().deltaT();
338 
339  const IOobject ddtIOobject
340  (
341  mesh().thisDb().newIOobject
342  (
343  "ddt("+rho.name()+','+vf.name()+')',
345  )
346  );
347 
348  if (mesh().moving())
349  {
350  return tmp<GeometricField<Type, faPatchField, areaMesh>>
351  (
352  new GeometricField<Type, faPatchField, areaMesh>
353  (
354  ddtIOobject,
355  mesh(),
356  rDeltaT.dimensions()*rho.dimensions()*vf.dimensions(),
357  rDeltaT.value()*
358  (
359  rho()*vf()
360  - rho.oldTime()()
361  *vf.oldTime()()*mesh().S0()/mesh().S()
362  ),
363  rDeltaT.value()*
364  (
365  rho.boundaryField()*vf.boundaryField()
366  - rho.oldTime().boundaryField()
367  *vf.oldTime().boundaryField()
368  )
369  )
370  );
371  }
372  else
373  {
374  return tmp<GeometricField<Type, faPatchField, areaMesh>>
375  (
376  new GeometricField<Type, faPatchField, areaMesh>
377  (
378  ddtIOobject,
379  rDeltaT*(rho*vf - rho.oldTime()*vf.oldTime())
380  )
381  );
382  }
383 }
384 
385 
386 template<class Type>
389 (
390  const areaScalarField& rho,
392 )
393 {
394  dimensionedScalar rDeltaT = 1.0/mesh().time().deltaT();
395 
396  const IOobject ddtIOobject
397  (
398  mesh().thisDb().newIOobject
399  (
400  "ddt0("+rho.name()+','+vf.name()+')',
402  )
403  );
404 
405  if (mesh().moving())
406  {
407  return tmp<GeometricField<Type, faPatchField, areaMesh>>
408  (
409  new GeometricField<Type, faPatchField, areaMesh>
410  (
411  ddtIOobject,
412  mesh(),
413  rDeltaT.dimensions()*rho.dimensions()*vf.dimensions(),
414  rDeltaT.value()*
415  (
416  - rho.oldTime()()
417  *vf.oldTime()()*mesh().S0()/mesh().S()
418  ),
419  rDeltaT.value()*
420  (
421  - rho.oldTime().boundaryField()
422  *vf.oldTime().boundaryField()
423  )
424  )
425  );
426  }
427  else
428  {
429  return tmp<GeometricField<Type, faPatchField, areaMesh>>
430  (
431  new GeometricField<Type, faPatchField, areaMesh>
432  (
433  ddtIOobject,
434  (-rDeltaT)*rho.oldTime()*vf.oldTime()
435  )
436  );
437  }
438 }
439 
440 template<class Type>
443 (
445 )
446 {
447  tmp<faMatrix<Type>> tfam
448  (
449  new faMatrix<Type>
450  (
451  vf,
453  )
454  );
455 
456  faMatrix<Type>& fam = tfam.ref();
457 
458  scalar rDeltaT = 1.0/mesh().time().deltaTValue();
459 
460  fam.diag() = rDeltaT*mesh().S();
461 
462  if (mesh().moving())
463  {
464  fam.source() = rDeltaT*vf.oldTime().primitiveField()*mesh().S0();
465  }
466  else
467  {
468  fam.source() = rDeltaT*vf.oldTime().primitiveField()*mesh().S();
469  }
470 
471  return tfam;
472 }
473 
474 
475 template<class Type>
478 (
479  const dimensionedScalar& rho,
481 )
482 {
483  tmp<faMatrix<Type>> tfam
484  (
485  new faMatrix<Type>
486  (
487  vf,
488  rho.dimensions()*vf.dimensions()*dimArea/dimTime
489  )
490  );
491  faMatrix<Type>& fam = tfam.ref();
492 
493  scalar rDeltaT = 1.0/mesh().time().deltaTValue();
494 
495  fam.diag() = rDeltaT*rho.value()*mesh().S();
496 
497  if (mesh().moving())
498  {
499  fam.source() = rDeltaT
500  *rho.value()*vf.oldTime().primitiveField()*mesh().S0();
501  }
502  else
503  {
504  fam.source() = rDeltaT
505  *rho.value()*vf.oldTime().primitiveField()*mesh().S();
506  }
507 
508  return tfam;
509 }
510 
511 
512 template<class Type>
515 (
516  const areaScalarField& rho,
518 )
519 {
520  tmp<faMatrix<Type>> tfam
521  (
522  new faMatrix<Type>
523  (
524  vf,
525  rho.dimensions()*vf.dimensions()*dimArea/dimTime
526  )
527  );
528  faMatrix<Type>& fam = tfam.ref();
529 
530  scalar rDeltaT = 1.0/mesh().time().deltaTValue();
531 
532  fam.diag() = rDeltaT*rho.primitiveField()*mesh().S();
533 
534  if (mesh().moving())
535  {
536  fam.source() = rDeltaT
537  *rho.oldTime().primitiveField()
538  *vf.oldTime().primitiveField()*mesh().S0();
539  }
540  else
541  {
542  fam.source() = rDeltaT
543  *rho.oldTime().primitiveField()
544  *vf.oldTime().primitiveField()*mesh().S();
545  }
546 
547  return tfam;
548 }
549 
550 
551 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
552 
553 } // End namespace fa
554 
555 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
556 
557 } // End namespace Foam
558 
559 // ************************************************************************* //
tmp< faMatrix< Type > > famDdt(const GeometricField< Type, faPatchField, areaMesh > &)
tmp< GeometricField< Type, faPatchField, areaMesh > > facDdt0(const dimensioned< Type >)
const Type & value() const noexcept
Return const reference to value.
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: tmpI.H:235
tmp< GeometricField< Type, faPatchField, areaMesh > > facDdt(const dimensioned< Type >)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Generic GeometricField class.
Generic dimensioned Type class.
dynamicFvMesh & mesh
const dimensionSet & dimensions() const noexcept
Return const reference to dimensions.
static const word & calculatedType() noexcept
The type name for calculated patch fields.
Definition: faPatchField.H:174
const word & name() const noexcept
Return const reference to name.
A special matrix type and solver, designed for finite area solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: faMatricesFwd.H:37
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
Calculate the matrix for the second temporal derivative.
Request registration (bool: true)
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.
const dimensionSet & dimensions() const noexcept
Return dimensions.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127