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