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 "incompressibleFluid.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  if (mesh.changing())
40  {
41  MRF.update();
42 
43  if (correctPhi || mesh.topoChanged())
44  {
45  // Calculate absolute flux
46  // from the mapped surface velocity
47  phi_ = mesh.Sf() & Uf();
48 
49  correctUphiBCs(U_, phi_, true);
50 
52  (
53  phi_,
54  U,
55  p,
59  pimple
60  );
61 
62  // Make the flux relative to the mesh motion
64  }
65 
66  meshCourantNo();
67  }
68  }
69 }
70 
71 
72 // ************************************************************************* //
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
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
Provides controls for the pressure reference in closed-volume simulations.
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
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
volVectorField U_
Velocity field.
autoPtr< surfaceVectorField > Uf
Pointer to the surface momentum field.
virtual void moveMesh()
Called at the start of the PIMPLE loop to move the mesh.
Definition: moveMesh.C:32
const volVectorField & U
Reference to the velocity field.
IOMRFZoneList MRF
MRF zone list.
surfaceScalarField phi_
Volumetric-flux field.
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.