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 "multiphaseEuler.H"
27 #include "fvcDiv.H"
28 #include "fvcMeshPhi.H"
29 
30 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
31 
33 {
34  if
35  (
36  pimple.flow()
38  )
39  {
40  if
41  (
43  // && divergent()
44  && !divU.valid()
45  )
46  {
47  // Construct and register divU for mapping
48  divU = new volScalarField
49  (
50  "divU0",
51  fvc::div
52  (
53  fvc::absolute(phi, fluid.movingPhases()[0].U())
54  )
55  );
56  }
57 
58  // Move the mesh
59  mesh_.move();
60 
61  if (mesh.changing())
62  {
64 
65  if (correctPhi || mesh.topoChanged())
66  {
68 
70  (
71  p_rgh,
72  divU,
74  pimple
75  );
76  }
77 
78  meshCourantNo();
79  }
80 
81  divU.clear();
82  }
83 }
84 
85 
86 // ************************************************************************* //
bool flow() const
Flag to indicate to solve for the flow.
bool move()
Move the mesh.
Definition: fvMesh.C:698
void correctPhi(const volScalarField &p_rgh, const autoPtr< volScalarField > &divU, const pressureReference &pressureReference, nonOrthogonalSolutionControl &pimple)
Definition: phaseSystem.C:730
virtual void meshUpdate()
Update the fluid properties for mesh changes.
Definition: phaseSystem.C:684
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
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
Solver module for steady or transient turbulent flow of compressible fluids with heat-transfer for HV...
Definition: fluid.H:70
const volVectorField & U
Velocity field.
volScalarField & p_rgh
Reference to the buoyant pressure for buoyant cases.
const surfaceScalarField & phi
Reference to the mass-flux field.
virtual void moveMesh()
Called at the start of the PIMPLE loop to move the mesh.
Definition: moveMesh.C:32
autoPtr< volScalarField > divU
Stored divU from the previous mesh so that it can be.
Calculate the divergence of the given field.
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
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
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61