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 "isothermalFluid.H"
27 #include "fvCorrectPhi.H"
28 #include "fvcMeshPhi.H"
29 
30 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
31 
33 {
35  {
36  // Move the mesh
37  mesh_.move();
38 
39  // The rhoU and rhoU0 fields can be cleared following mesh-motion
40  // now the mesh has been re-stitched as necessary
41  rhoU.clear();
42  rhoU0.clear();
43 
44  if (mesh.changing())
45  {
46  if (buoyancy.valid())
47  {
48  buoyancy->moveMesh();
49  }
50 
51  MRF.update();
52 
53  if (correctPhi || mesh.topoChanged())
54  {
55  // Calculate absolute flux
56  // from the mapped surface velocity
57  phi_ = mesh.Sf() & rhoUf();
58 
59  correctUphiBCs(rho, U_, phi_, true);
60 
62  (
63  phi_,
64  buoyancy.valid() ? p_rgh : p,
65  thermo.psi(),
67  divrhoU(),
68  pimple
69  );
70 
71  // Make the fluxes relative to the mesh-motion
73  }
74 
75  meshCourantNo();
76  }
77  }
78 }
79 
80 
81 // ************************************************************************* //
void update()
Update MRFZone faces if the mesh topology changes.
Definition: MRFZoneList.C:366
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
virtual const volScalarField & psi() const =0
Compressibility [s^2/m^2].
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
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
autoPtr< volScalarField > divrhoU
Pointer to the vol momentum divergence field.
volScalarField & p_rgh
Reference to the buoyant pressure for buoyant cases.
volVectorField U_
Velocity field.
autoPtr< surfaceVectorField > rhoUf
Pointer to the surface momentum field.
autoPtr< volVectorField > rhoU
Pointer to the vol momentum field.
virtual void moveMesh()
Called at the start of the PIMPLE loop to move the mesh.
Definition: moveMesh.C:32
const volVectorField & U
Velocity field.
IOMRFZoneList MRF
MRF zone list.
autoPtr< volVectorField > rhoU0
Pointer to the old-time vol momentum field.
surfaceScalarField phi_
Mass-flux field.
const volScalarField & rho
Reference to the continuity density field.
const fluidThermo & thermo
Reference to the fluid thermophysical properties.
const volScalarField & p
Reference to the pressure field.
Flux correction functions to ensure continuity.
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
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:89
void correctUphiBCs(volVectorField &U, surfaceScalarField &phi, const bool evaluateUBCs)
If the mesh is moving correct the velocity BCs on the moving walls to.