correctUphiBCs.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) 2015-2016 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
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 "CorrectPhi.H"
29
30// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31
33(
36)
37{
38 const fvMesh& mesh = U.mesh();
39
40 if (mesh.changing())
41 {
42 volVectorField::Boundary& Ubf = U.boundaryFieldRef();
44 phi.boundaryFieldRef();
45
46 forAll(Ubf, patchi)
47 {
48 if (Ubf[patchi].fixesValue())
49 {
50 Ubf[patchi].initEvaluate();
51 }
52 }
53
54 forAll(Ubf, patchi)
55 {
56 if (Ubf[patchi].fixesValue())
57 {
58 Ubf[patchi].evaluate();
59
60 phibf[patchi] = Ubf[patchi] & mesh.Sf().boundaryField()[patchi];
61 }
62 }
63 }
64}
65
66
68(
69 const volScalarField& rho,
72)
73{
74 const fvMesh& mesh = U.mesh();
75
76 if (mesh.changing())
77 {
78 volVectorField::Boundary& Ubf = U.boundaryFieldRef();
80 phi.boundaryFieldRef();
81
82 forAll(Ubf, patchi)
83 {
84 if (Ubf[patchi].fixesValue())
85 {
86 Ubf[patchi].initEvaluate();
87 }
88 }
89
90 forAll(Ubf, patchi)
91 {
92 if (Ubf[patchi].fixesValue())
93 {
94 Ubf[patchi].evaluate();
95
96 phibf[patchi] =
97 rho.boundaryField()[patchi]
98 *(
99 Ubf[patchi]
100 & mesh.Sf().boundaryField()[patchi]
101 );
102 }
103 }
104 }
105}
106
107
108// ************************************************************************* //
surfaceScalarField & phi
void evaluate()
Evaluate boundary conditions.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
U
Definition: pEqn.H:72
dynamicFvMesh & mesh
void correctUphiBCs(volVectorField &U, surfaceScalarField &phi)
If the mesh is moving correct the velocity BCs on the moving walls to.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333