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-2024 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 }
56 
57 
59 {
60  if (pimple.firstIter() || pimple.moveMeshOuterCorrectors())
61  {
62  if (mesh.changing())
63  {
65 
66  MRF.update();
67 
68  if (correctPhi || mesh.topoChanged())
69  {
70  // Calculate absolute flux
71  // from the mapped surface velocity
72  phi_ = mesh.Sf() & Uf();
73 
74  correctUphiBCs(U_, phi_, true);
75 
76  if (incompressible())
77  {
79  (
80  phi_,
81  U,
82  p_rgh,
83  rAU,
84  divU,
86  pimple
87  );
88  }
89  else
90  {
92  (
93  phi_,
94  p_rgh,
95  psiByRho(),
96  rAU,
97  divU(),
98  pimple
99  );
100  }
101 
102  // Make the fluxes relative to the mesh motion
103  MRF.makeRelative(phi_);
104  fvc::makeRelative(phi_, U);
105  }
106 
107  meshCourantNo();
108 
109  correctInterface();
110  }
111 
112  divU.clear();
113  }
114 }
115 
116 
117 // ************************************************************************* //
bool move()
Move the mesh.
Definition: fvMesh.C:653
bool moveMeshOuterCorrectors() const
Switch to move the mesh at the start of every PIMPLE.
bool firstIter() const
Flag to indicate the first iteration.
bool topoChanged() const
Has the mesh topology changed this time-step.
Definition: polyMesh.H:482
Provides controls for the pressure reference in closed-volume simulations.
pimpleNoLoopControl pimple
PIMPLE inner-loop controls.
Definition: solver.H:107
fvMesh & mesh_
Region mesh.
Definition: solver.H:62
const fvMesh & mesh
Region mesh.
Definition: solver.H:101
const surfaceScalarField & phi
Reference to the mass-flux field.
Definition: VoFSolver.H:212
virtual void moveMesh()
Called at the start of the PIMPLE loop to move the mesh.
Definition: moveMesh.C:33
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 bool divergent() const =0
Is the flow divergent?
virtual void motionCorrector()
Corrections that follow mesh motion.
Definition: moveMesh.C:58
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:101
IOMRFZoneList MRF(mesh)
pimpleControl pimple(mesh)
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.
volScalarField rAU(1.0/UEqn.A())
U
Definition: pEqn.H:72
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:64
void correctUphiBCs(volVectorField &U, surfaceScalarField &phi, const bool evaluateUBCs)
If the mesh is moving correct the velocity BCs on the moving walls to.
correctPhi