Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / src / turbulenceModels / incompressible / LES / dynLagrangian / dynLagrangian.C
blobf75df11c8c9811db9735140cd64dcdfd869c7aae
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2010-2010 OpenCFD Ltd.
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 the
13     Free Software Foundation; either version 3 of the License, or (at your
14     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, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "dynLagrangian.H"
28 #include "addToRunTimeSelectionTable.H"
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 namespace Foam
34 namespace incompressible
36 namespace LESModels
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 defineTypeNameAndDebug(dynLagrangian, 0);
42 addToRunTimeSelectionTable(LESModel, dynLagrangian, dictionary);
44 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
46 void dynLagrangian::updateSubGridScaleFields
48     const tmp<volTensorField>& gradU
51     nuSgs_ = (flm_/fmm_)*delta()*sqrt(k(gradU));
52     nuSgs_.correctBoundaryConditions();
56 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
58 dynLagrangian::dynLagrangian
60     const volVectorField& U,
61     const surfaceScalarField& phi,
62     transportModel& transport,
63     const word& turbulenceModelName,
64     const word& modelName
67     LESModel(modelName, U, phi, transport, turbulenceModelName),
68     GenEddyVisc(U, phi, transport),
70     flm_
71     (
72         IOobject
73         (
74             "flm",
75             runTime_.timeName(),
76             mesh_,
77             IOobject::MUST_READ,
78             IOobject::AUTO_WRITE
79         ),
80         mesh_
81     ),
82     fmm_
83     (
84         IOobject
85         (
86             "fmm",
87             runTime_.timeName(),
88             mesh_,
89             IOobject::MUST_READ,
90             IOobject::AUTO_WRITE
91         ),
92         mesh_
93     ),
94     theta_
95     (
96         dimensioned<scalar>::lookupOrAddToDict
97         (
98             "theta",
99             coeffDict_,
100             1.5
101         )
102     ),
103     simpleFilter_(U.mesh()),
104     filterPtr_(LESfilter::New(U.mesh(), coeffDict())),
105     filter_(filterPtr_()),
106     flm0_("flm0", flm_.dimensions(), 0.0),
107     fmm0_("fmm0", fmm_.dimensions(), VSMALL)
109     updateSubGridScaleFields(fvc::grad(U));
111     printCoeffs();
115 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
117 void dynLagrangian::correct(const tmp<volTensorField>& gradU)
119     LESModel::correct(gradU);
121     volSymmTensorField S(dev(symm(gradU())));
123     volScalarField magS(mag(S));
125     volVectorField Uf(filter_(U()));
127     volSymmTensorField Sf(dev(symm(fvc::grad(Uf))));
129     volScalarField magSf(mag(Sf));
131     volSymmTensorField L(dev(filter_(sqr(U())) - (sqr(filter_(U())))));
133     volSymmTensorField M(2.0*sqr(delta())*(filter_(magS*S) - 4.0*magSf*Sf));
135     volScalarField invT
136     (
137         (1.0/(theta_.value()*delta()))*pow(flm_*fmm_, 1.0/8.0)
138     );
140     volScalarField LM(L && M);
142     fvScalarMatrix flmEqn
143     (
144         fvm::ddt(flm_)
145       + fvm::div(phi(), flm_)
146      ==
147         invT*LM
148       - fvm::Sp(invT, flm_)
149     );
151     flmEqn.relax();
152     flmEqn.solve();
154     bound(flm_, flm0_);
156     volScalarField MM(M && M);
158     fvScalarMatrix fmmEqn
159     (
160         fvm::ddt(fmm_)
161       + fvm::div(phi(), fmm_)
162      ==
163         invT*MM
164       - fvm::Sp(invT, fmm_)
165     );
167     fmmEqn.relax();
168     fmmEqn.solve();
170     bound(fmm_, fmm0_);
172     updateSubGridScaleFields(gradU);
176 bool dynLagrangian::read()
178     if (GenEddyVisc::read())
179     {
180         filter_.read(coeffDict());
181         theta_.readIfPresent(coeffDict());
183         return true;
184     }
185     else
186     {
187         return false;
188     }
192 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
194 } // End namespace LESModels
195 } // End namespace incompressible
196 } // End namespace Foam
198 // ************************************************************************* //