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 "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 }
40 
41 
43 {
44  if (pimple.firstIter() || pimple.moveMeshOuterCorrectors())
45  {
46  // The rhoU field can be cleared following mesh-motion
47  rhoU.clear();
48 
49  if (mesh.changing())
50  {
51  if (buoyancy.valid())
52  {
53  buoyancy->moveMesh();
54  }
55 
56  MRF.update();
57 
58  if (correctPhi || mesh.topoChanged())
59  {
60  // Calculate absolute flux
61  // from the mapped surface velocity
62  phi_ = mesh.Sf() & rhoUf();
63 
64  correctUphiBCs(rho, U_, phi_, true);
65 
67  (
68  phi_,
69  buoyancy.valid() ? p_rgh : p,
70  thermo.psi(),
72  divrhoU(),
73  pimple
74  );
75 
76  // Make the fluxes relative to the mesh-motion
77  MRF.makeRelative(phi_);
78  fvc::makeRelative(phi_, rho, U);
79  }
80 
81  meshCourantNo();
82  }
83  }
84 }
85 
86 
87 // ************************************************************************* //
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
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.
pimpleNoLoopControl pimple
PIMPLE inner-loop controls.
Definition: solver.H:107
fvMesh & mesh_
Region mesh.
Definition: solver.H:62
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
virtual void moveMesh()
Called at the start of the PIMPLE loop to move the mesh.
Definition: moveMesh.C:32
virtual void motionCorrector()
Corrections that follow mesh motion.
Definition: moveMesh.C:42
IOMRFZoneList MRF(mesh)
pimpleControl pimple(mesh)
Flux correction functions to ensure continuity.
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
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
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.
correctPhi
volScalarField & p
fluidMulticomponentThermo & thermo
Definition: createFields.H:31