rhoCentralDyMFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2013 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 Application
25  rhoCentralDyMFoam
26 
27 Description
28  Density-based compressible flow solver based on central-upwind schemes of
29  Kurganov and Tadmor
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #include "fvCFD.H"
34 #include "psiThermo.H"
35 #include "turbulenceModel.H"
38 #include "motionSolver.H"
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 int main(int argc, char *argv[])
43 {
44  #include "setRootCase.H"
45 
46  #include "createTime.H"
47  #include "createMesh.H"
48  #include "createFields.H"
49  #include "readTimeControls.H"
50 
51  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53  #include "readFluxScheme.H"
54 
55  dimensionedScalar v_zero("v_zero", dimVolume/dimTime, 0.0);
56 
57  Info<< "\nStarting time loop\n" << endl;
58 
59  autoPtr<Foam::motionSolver> motionPtr = motionSolver::New(mesh);
60 
61  while (runTime.run())
62  {
63  // --- upwind interpolation of primitive fields on faces
64 
65  surfaceScalarField rho_pos
66  (
67  fvc::interpolate(rho, pos, "reconstruct(rho)")
68  );
69  surfaceScalarField rho_neg
70  (
71  fvc::interpolate(rho, neg, "reconstruct(rho)")
72  );
73 
74  surfaceVectorField rhoU_pos
75  (
76  fvc::interpolate(rhoU, pos, "reconstruct(U)")
77  );
78  surfaceVectorField rhoU_neg
79  (
80  fvc::interpolate(rhoU, neg, "reconstruct(U)")
81  );
82 
83  volScalarField rPsi(1.0/psi);
84  surfaceScalarField rPsi_pos
85  (
86  fvc::interpolate(rPsi, pos, "reconstruct(T)")
87  );
88  surfaceScalarField rPsi_neg
89  (
90  fvc::interpolate(rPsi, neg, "reconstruct(T)")
91  );
92 
93  surfaceScalarField e_pos
94  (
95  fvc::interpolate(e, pos, "reconstruct(T)")
96  );
97  surfaceScalarField e_neg
98  (
99  fvc::interpolate(e, neg, "reconstruct(T)")
100  );
101 
102  surfaceVectorField U_pos(rhoU_pos/rho_pos);
103  surfaceVectorField U_neg(rhoU_neg/rho_neg);
104 
105  surfaceScalarField p_pos(rho_pos*rPsi_pos);
106  surfaceScalarField p_neg(rho_neg*rPsi_neg);
107 
108  surfaceScalarField phiv_pos(U_pos & mesh.Sf());
109  surfaceScalarField phiv_neg(U_neg & mesh.Sf());
110 
111  fvc::makeRelative(phiv_pos, U);
112  fvc::makeRelative(phiv_neg, U);
113 
114  volScalarField c(sqrt(thermo.Cp()/thermo.Cv()*rPsi));
115  surfaceScalarField cSf_pos
116  (
117  fvc::interpolate(c, pos, "reconstruct(T)")*mesh.magSf()
118  );
119  surfaceScalarField cSf_neg
120  (
121  fvc::interpolate(c, neg, "reconstruct(T)")*mesh.magSf()
122  );
123 
124  surfaceScalarField ap
125  (
126  max(max(phiv_pos + cSf_pos, phiv_neg + cSf_neg), v_zero)
127  );
128  surfaceScalarField am
129  (
130  min(min(phiv_pos - cSf_pos, phiv_neg - cSf_neg), v_zero)
131  );
132 
133  surfaceScalarField a_pos(ap/(ap - am));
134 
135  surfaceScalarField amaxSf("amaxSf", max(mag(am), mag(ap)));
136 
137  surfaceScalarField aSf(am*a_pos);
138 
139  if (fluxScheme == "Tadmor")
140  {
141  aSf = -0.5*amaxSf;
142  a_pos = 0.5;
143  }
144 
145  surfaceScalarField a_neg(1.0 - a_pos);
146 
147  phiv_pos *= a_pos;
148  phiv_neg *= a_neg;
149 
150  surfaceScalarField aphiv_pos(phiv_pos - aSf);
151  surfaceScalarField aphiv_neg(phiv_neg + aSf);
152 
153  // Reuse amaxSf for the maximum positive and negative fluxes
154  // estimated by the central scheme
155  amaxSf = max(mag(aphiv_pos), mag(aphiv_neg));
156 
157  #include "compressibleCourantNo.H"
158  #include "readTimeControls.H"
159  #include "setDeltaT.H"
160 
161  runTime++;
162 
163  Info<< "Time = " << runTime.timeName() << nl << endl;
164 
165  mesh.movePoints(motionPtr->newPoints());
166 
167  phi = aphiv_pos*rho_pos + aphiv_neg*rho_neg;
168 
169  surfaceVectorField phiUp
170  (
171  (aphiv_pos*rhoU_pos + aphiv_neg*rhoU_neg)
172  + (a_pos*p_pos + a_neg*p_neg)*mesh.Sf()
173  );
174 
175  surfaceScalarField phiEp
176  (
177  aphiv_pos*(rho_pos*(e_pos + 0.5*magSqr(U_pos)) + p_pos)
178  + aphiv_neg*(rho_neg*(e_neg + 0.5*magSqr(U_neg)) + p_neg)
179  + mesh.phi()*(a_pos*p_pos + a_neg*p_neg)
180  + aSf*p_pos - aSf*p_neg
181  );
182 
183  volScalarField muEff(turbulence->muEff());
184  volTensorField tauMC("tauMC", muEff*dev2(Foam::T(fvc::grad(U))));
185 
186  // --- Solve density
188 
189  // --- Solve momentum
190  solve(fvm::ddt(rhoU) + fvc::div(phiUp));
191 
192  U.dimensionedInternalField() =
193  rhoU.dimensionedInternalField()
194  /rho.dimensionedInternalField();
195  U.correctBoundaryConditions();
196  rhoU.boundaryField() = rho.boundaryField()*U.boundaryField();
197 
198  if (!inviscid)
199  {
200  solve
201  (
202  fvm::ddt(rho, U) - fvc::ddt(rho, U)
204  - fvc::div(tauMC)
205  );
206  rhoU = rho*U;
207  }
208 
209  // --- Solve energy
210  surfaceScalarField sigmaDotU
211  (
212  (
214  + (mesh.Sf() & fvc::interpolate(tauMC))
215  )
216  & (a_pos*U_pos + a_neg*U_neg)
217  );
218 
219  solve
220  (
221  fvm::ddt(rhoE)
222  + fvc::div(phiEp)
223  - fvc::div(sigmaDotU)
224  );
225 
226  e = rhoE/rho - 0.5*magSqr(U);
227  e.correctBoundaryConditions();
228  thermo.correct();
229  rhoE.boundaryField() =
230  rho.boundaryField()*
231  (
232  e.boundaryField() + 0.5*magSqr(U.boundaryField())
233  );
234 
235  if (!inviscid)
236  {
237  solve
238  (
239  fvm::ddt(rho, e) - fvc::ddt(rho, e)
240  - fvm::laplacian(turbulence->alphaEff(), e)
241  );
242  thermo.correct();
243  rhoE = rho*(e + 0.5*magSqr(U));
244  }
245 
246  p.dimensionedInternalField() =
247  rho.dimensionedInternalField()
248  /psi.dimensionedInternalField();
249  p.correctBoundaryConditions();
250  rho.boundaryField() = psi.boundaryField()*p.boundaryField();
251 
252  turbulence->correct();
253 
254  runTime.write();
255 
256  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
257  << " ClockTime = " << runTime.elapsedClockTime() << " s"
258  << nl << endl;
259  }
260 
261  Info<< "End\n" << endl;
262 
263  return 0;
264 }
265 
266 // ************************************************************************* //