constrainPressure.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 OpenFOAM Foundation
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 "constrainPressure.H"
29 #include "volFields.H"
30 #include "surfaceFields.H"
31 #include "geometricOneField.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 template<class RhoType, class RAUType, class MRFType>
38 (
40  const RhoType& rho,
41  const volVectorField& U,
43  const RAUType& rhorAU,
44  const MRFType& MRF
45 )
46 {
47  const fvMesh& mesh = p.mesh();
48 
49  volScalarField::Boundary& pBf = p.boundaryFieldRef();
50 
51  const volVectorField::Boundary& UBf = U.boundaryField();
52  const surfaceScalarField::Boundary& phiHbyABf =
53  phiHbyA.boundaryField();
54  const typename RAUType::Boundary& rhorAUBf =
55  rhorAU.boundaryField();
56  const surfaceVectorField::Boundary& SfBf =
57  mesh.Sf().boundaryField();
58  const surfaceScalarField::Boundary& magSfBf =
59  mesh.magSf().boundaryField();
60 
61  forAll(pBf, patchi)
62  {
63  if (isA<fixedFluxPressureFvPatchScalarField>(pBf[patchi]))
64  {
65  refCast<fixedFluxPressureFvPatchScalarField>
66  (
67  pBf[patchi]
68  ).updateSnGrad
69  (
70  (
71  phiHbyABf[patchi]
72  - rho.boundaryField()[patchi]
73  *MRF.relative(SfBf[patchi] & UBf[patchi], patchi)
74  )
75  /(magSfBf[patchi]*rhorAUBf[patchi])
76  );
77  }
78  }
79 }
80 
81 
82 template<class RAUType>
84 (
86  const volScalarField& rho,
87  const volVectorField& U,
89  const RAUType& rAU
90 )
91 {
93 }
94 
95 
96 template<class RAUType, class MRFType>
98 (
100  const volVectorField& U,
102  const RAUType& rAU,
103  const MRFType& MRF
104 )
105 {
107 }
108 
109 
110 template<class RAUType>
112 (
113  volScalarField& p,
114  const volVectorField& U,
116  const RAUType& rAU
117 )
118 {
120 }
121 
122 
123 // ************************************************************************* //
Foam::constrainPressure
void constrainPressure(volScalarField &p, const RhoType &rho, const volVectorField &U, const surfaceScalarField &phiHbyA, const RAUType &rhorAU, const MRFType &MRF)
Definition: constrainPressure.C:38
volFields.H
p
volScalarField & p
Definition: createFieldRefs.H:8
rAU
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1
Foam::NullMRF
Definition: constrainPressure.H:53
constrainPressure.H
Foam::geometricOneField
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Definition: geometricOneField.H:55
surfaceFields.H
Foam::surfaceFields.
rho
rho
Definition: readInitialConditions.H:88
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
MRF
IOMRFZoneList & MRF
Definition: setRegionFluidFields.H:22
phiHbyA
phiHbyA
Definition: pcEqn.H:73
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:83
U
U
Definition: pEqn.H:72
geometricOneField.H
fixedFluxPressureFvPatchScalarField.H
Foam::GeometricField< scalar, fvPatchField, volMesh >