correctUphiBCs.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) 2015-2023 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 "fvCorrectPhi.H"
27 #include "volFields.H"
28 #include "surfaceFields.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
33 (
35  surfaceScalarField& phi,
36  const bool evaluateUBCs
37 )
38 {
39  const fvMesh& mesh = U.mesh();
40 
41  if (mesh.changing())
42  {
43  volVectorField::Boundary& Ubf = U.boundaryFieldRef();
45 
46  if (evaluateUBCs)
47  {
48  forAll(Ubf, patchi)
49  {
50  if (Ubf[patchi].fixesValue())
51  {
52  Ubf[patchi].initEvaluate();
53  }
54  }
55  }
56 
57  forAll(Ubf, patchi)
58  {
59  if (Ubf[patchi].fixesValue())
60  {
61  if (evaluateUBCs)
62  {
63  Ubf[patchi].evaluate();
64  }
65 
66  phibf[patchi] = Ubf[patchi] & mesh.Sf().boundaryField()[patchi];
67  }
68  }
69  }
70 }
71 
72 
74 (
75  const volScalarField& rho,
77  surfaceScalarField& phi,
78  const bool evaluateUBCs
79 )
80 {
81  const fvMesh& mesh = U.mesh();
82 
83  if (mesh.changing())
84  {
85  volVectorField::Boundary& Ubf = U.boundaryFieldRef();
87 
88  if (evaluateUBCs)
89  {
90  forAll(Ubf, patchi)
91  {
92  if (Ubf[patchi].fixesValue())
93  {
94  Ubf[patchi].initEvaluate();
95  }
96  }
97  }
98 
99  forAll(Ubf, patchi)
100  {
101  if (Ubf[patchi].fixesValue())
102  {
103  if (evaluateUBCs)
104  {
105  Ubf[patchi].evaluate();
106  }
107 
108  phibf[patchi] =
109  rho.boundaryField()[patchi]
110  *(
111  Ubf[patchi]
112  & mesh.Sf().boundaryField()[patchi]
113  );
114  }
115  }
116  }
117 }
118 
119 
120 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Generic GeometricBoundaryField class.
void evaluate()
Evaluate boundary conditions.
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const surfaceVectorField & Sf() const
Return cell face area vectors.
bool changing() const
Is mesh changing.
Definition: polyMesh.H:490
Flux correction functions to ensure continuity.
label patchi
U
Definition: pEqn.H:72
void correctUphiBCs(volVectorField &U, surfaceScalarField &phi, const bool evaluateUBCs)
If the mesh is moving correct the velocity BCs on the moving walls to.
Foam::surfaceFields.