fvcMeshPhiTemplates.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) 2022 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 template<class MRFType>
28 (
30  const volVectorField& U,
31  const surfaceScalarField& phi,
32  const MRFType& MRF
33 )
34 {
35  const fvMesh& mesh = U.mesh();
36 
37  if (mesh.dynamic())
38  {
39  Uf() = fvc::interpolate(U);
40  surfaceVectorField n(mesh.Sf()/mesh.magSf());
41  Uf() += n*(MRF.absolute(phi)/mesh.magSf() - (n & Uf()));
42  }
43 }
44 
45 
46 template<class MRFType>
48 (
50  const volScalarField& rho,
51  const volVectorField& U,
52  const surfaceScalarField& phi,
53  const MRFType& MRF
54 )
55 {
56  const fvMesh& mesh = U.mesh();
57 
58  if (mesh.dynamic())
59  {
60  rhoUf() = fvc::interpolate(rho*U);
61  surfaceVectorField n(mesh.Sf()/mesh.magSf());
62  rhoUf() +=
63  n*(
64  MRF.absolute(fvc::absolute(phi, rho, U))/mesh.magSf()
65  - (n & rhoUf())
66  );
67  }
68 }
69 
70 
71 // ************************************************************************* //
const surfaceVectorField & Sf() const
Return cell face area vectors.
fvMesh & mesh
autoPtr< surfaceVectorField > Uf
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
const Mesh & mesh() const
Return mesh.
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
void correctRhoUf(autoPtr< surfaceVectorField > &rhoUf, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const MRFType &MRF)
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:202
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
bool dynamic() const
Is this mesh dynamic?
Definition: fvMesh.C:578
void correctUf(autoPtr< surfaceVectorField > &Uf, const volVectorField &U, const surfaceScalarField &phi, const MRFType &MRF)
label n
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
autoPtr< surfaceVectorField > rhoUf