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  Copyright (C) 2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "constrainPressure.H"
30 #include "volFields.H"
31 #include "surfaceFields.H"
32 #include "geometricOneField.H"
33 #include "updateableSnGrad.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 template<class RhoType, class RAUType, class MRFType>
39 (
41  const RhoType& rho,
42  const volVectorField& U,
44  const RAUType& rhorAU,
45  const MRFType& MRF
46 )
47 {
48  const fvMesh& mesh = p.mesh();
49 
50  volScalarField::Boundary& pBf = p.boundaryFieldRef();
51 
52  const volVectorField::Boundary& UBf = U.boundaryField();
53  const surfaceScalarField::Boundary& phiHbyABf =
54  phiHbyA.boundaryField();
55  const typename RAUType::Boundary& rhorAUBf = rhorAU.boundaryField();
56  const surfaceVectorField::Boundary& SfBf = mesh.Sf().boundaryField();
57  const surfaceScalarField::Boundary& magSfBf =
58  mesh.magSf().boundaryField();
59 
60  forAll(pBf, patchi)
61  {
62  typedef updateablePatchTypes::updateableSnGrad snGradType;
63  const auto* snGradPtr = isA<snGradType>(pBf[patchi]);
64 
65  if (snGradPtr)
66  {
67  const_cast<snGradType&>(*snGradPtr).updateSnGrad
68  (
69  (
70  phiHbyABf[patchi]
71  - rho.boundaryField()[patchi]
72  *MRF.relative(SfBf[patchi] & UBf[patchi], patchi)
73  )
74  /(magSfBf[patchi]*rhorAUBf[patchi])
75  );
76  }
77  }
78 }
79 
80 
81 template<class RAUType>
83 (
85  const volScalarField& rho,
86  const volVectorField& U,
88  const RAUType& rAU
89 )
90 {
92 }
93 
94 
95 template<class RAUType, class MRFType>
97 (
99  const volVectorField& U,
101  const RAUType& rAU,
102  const MRFType& MRF
103 )
104 {
106 }
107 
108 
109 template<class RAUType>
111 (
112  volScalarField& p,
113  const volVectorField& U,
115  const RAUType& rAU
116 )
117 {
119 }
120 
121 
122 // ************************************************************************* //
Foam::constrainPressure
void constrainPressure(volScalarField &p, const RhoType &rho, const volVectorField &U, const surfaceScalarField &phiHbyA, const RAUType &rhorAU, const MRFType &MRF)
Definition: constrainPressure.C:39
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
Foam::updateablePatchTypes::updateableSnGrad
Definition: updateableSnGrad.H:55
MRF
IOMRFZoneList & MRF
Definition: setRegionFluidFields.H:22
phiHbyA
phiHbyA
Definition: pcEqn.H:73
updateableSnGrad.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
U
U
Definition: pEqn.H:72
geometricOneField.H
Foam::GeometricField< scalar, fvPatchField, volMesh >