moveMesh.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) 2022-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 "VoFSolver.H"
27 #include "fvCorrectPhi.H"
28 #include "fvcMeshPhi.H"
29 #include "fvcDiv.H"
30 
31 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
32 
34 {
36  {
37  if
38  (
40  && divergent()
41  && !divU.valid()
42  )
43  {
44  // Construct and register divU for correctPhi
45  divU = new volScalarField
46  (
47  "divU0",
49  );
50  }
51 
52  // Move the mesh
53  mesh_.move();
54 
55  if (mesh.changing())
56  {
58 
59  MRF.update();
60 
61  if (correctPhi || mesh.topoChanged())
62  {
63  // Calculate absolute flux
64  // from the mapped surface velocity
65  phi_ = mesh.Sf() & Uf();
66 
67  correctUphiBCs(U_, phi_, true);
68 
69  if (incompressible())
70  {
72  (
73  phi_,
74  U,
75  p_rgh,
76  rAU,
77  divU,
79  pimple
80  );
81  }
82  else
83  {
85  (
86  phi_,
87  p_rgh,
88  psiByRho(),
89  rAU,
90  divU(),
91  pimple
92  );
93  }
94 
95  // Make the fluxes relative to the mesh motion
97  }
98 
99  meshCourantNo();
100 
102  }
103 
104  divU.clear();
105  }
106 }
107 
108 
109 // ************************************************************************* //
void update()
Update MRFZone faces if the mesh topology changes.
Definition: MRFZoneList.C:366
const surfaceVectorField & Sf() const
Return cell face area vectors.
bool move()
Move the mesh.
Definition: fvMesh.C:698
bool moveMeshOuterCorrectors() const
Switch to move the mesh at the start of every PIMPLE.
bool firstIter() const
Flag to indicate the first iteration.
bool changing() const
Is mesh changing.
Definition: polyMesh.H:488
bool topoChanged() const
Has the mesh topology changed this time-step.
Definition: polyMesh.H:481
pimpleNoLoopControl pimple
PIMPLE inner-loop controls.
Definition: solver.H:100
fvMesh & mesh_
Region mesh.
Definition: solver.H:61
const fvMesh & mesh
Region mesh.
Definition: solver.H:94
volScalarField & p_rgh
Reference to the buoyant pressure for buoyant cases.
Definition: VoFSolver.H:107
const surfaceScalarField & phi
Reference to the mass-flux field.
Definition: VoFSolver.H:212
volVectorField U_
Velocity field.
Definition: VoFSolver.H:94
autoPtr< surfaceVectorField > Uf
Pointer to the surface momentum field.
Definition: VoFSolver.H:133
virtual void moveMesh()
Called at the start of the PIMPLE loop to move the mesh.
Definition: moveMesh.C:33
autoPtr< volScalarField > rAU
Inverse momentum equation diagonal.
Definition: VoFSolver.H:129
autoPtr< volScalarField > divU
Pointer to the momentum divergence field.
Definition: VoFSolver.H:138
const volVectorField & U
Reference to the velocity field.
Definition: VoFSolver.H:209
virtual const Foam::pressureReference & pressureReference() const =0
Return the pressure reference.
virtual bool divergent() const =0
Is the flow divergent?
IOMRFZoneList MRF
MRF zone list.
Definition: VoFSolver.H:122
virtual bool incompressible() const =0
Is the flow incompressible?
surfaceScalarField phi_
Volumetric flux field.
Definition: VoFSolver.H:97
virtual tmp< volScalarField > psiByRho() const =0
Return the mixture compressibility/density.
virtual void correctInterface()=0
Correct the interface properties following mesh-change.
Buoyancy related data for the Foam::solvers::isothermalFluid solver module when solving buoyant cases...
Definition: buoyancy.H:70
void moveMesh()
Update gh and ghf following mesh-motion.
Definition: buoyancy.C:124
bool correctPhi
Switch to correct the flux after mesh change.
Definition: fluidSolver.H:95
void meshCourantNo() const
Check mesh Courant numbers for moving mesh cases.
Definition: fluidSolver.C:67
Flux correction functions to ensure continuity.
Calculate the divergence of the given field.
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
void correctPhi(surfaceScalarField &phi, const volVectorField &U, const volScalarField &p, const autoPtr< volScalarField > &rAU, const autoPtr< volScalarField > &divU, const pressureReference &pressureReference, nonOrthogonalSolutionControl &pcorrControl)
Definition: fvCorrectPhi.C:32
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:202
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:89
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
void correctUphiBCs(volVectorField &U, surfaceScalarField &phi, const bool evaluateUBCs)
If the mesh is moving correct the velocity BCs on the moving walls to.