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