diff --git a/src/TurbulenceModels/compressible/turbulentFluidThermoModels/turbulentFluidThermoModels.C b/src/TurbulenceModels/compressible/turbulentFluidThermoModels/turbulentFluidThermoModels.C
index 10bacefd988db14ad57c2a1cb9eb328e20697e53..048ea5ee57ca3847b318b6a5ee78238b5e31aa0e 100644
--- a/src/TurbulenceModels/compressible/turbulentFluidThermoModels/turbulentFluidThermoModels.C
+++ b/src/TurbulenceModels/compressible/turbulentFluidThermoModels/turbulentFluidThermoModels.C
@@ -94,6 +94,9 @@ makeLESModel(WALE);
 #include "kEqn.H"
 makeLESModel(kEqn);
 
+#include "dynamicKEqn.H"
+makeLESModel(dynamicKEqn);
+
 #include "SpalartAllmarasDES.H"
 makeLESModel(SpalartAllmarasDES);
 
diff --git a/src/TurbulenceModels/incompressible/turbulentTransportModels/turbulentTransportModels.C b/src/TurbulenceModels/incompressible/turbulentTransportModels/turbulentTransportModels.C
index 786038a8df93fed6f1b4c9297a5306229a021c6f..ac9c460de2334df56164fa52218122b3cc2f25ae 100644
--- a/src/TurbulenceModels/incompressible/turbulentTransportModels/turbulentTransportModels.C
+++ b/src/TurbulenceModels/incompressible/turbulentTransportModels/turbulentTransportModels.C
@@ -89,6 +89,9 @@ makeLESModel(WALE);
 #include "kEqn.H"
 makeLESModel(kEqn);
 
+#include "dynamicKEqn.H"
+makeLESModel(dynamicKEqn);
+
 #include "SpalartAllmarasDES.H"
 makeLESModel(SpalartAllmarasDES);
 
diff --git a/src/TurbulenceModels/turbulenceModels/LES/dynamicKEqn/dynamicKEqn.C b/src/TurbulenceModels/turbulenceModels/LES/dynamicKEqn/dynamicKEqn.C
new file mode 100644
index 0000000000000000000000000000000000000000..9f98bea02479c0f44efeef19992ff5505395a088
--- /dev/null
+++ b/src/TurbulenceModels/turbulenceModels/LES/dynamicKEqn/dynamicKEqn.C
@@ -0,0 +1,285 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | Copyright (C) 2011-2015 OpenFOAM Foundation
+     \\/     M anipulation  |
+-------------------------------------------------------------------------------
+License
+    This file is part of OpenFOAM.
+
+    OpenFOAM is free software: you can redistribute it and/or modify it
+    under the terms of the GNU General Public License as published by
+    the Free Software Foundation, either version 3 of the License, or
+    (at your option) any later version.
+
+    OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
+    ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
+    FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
+    for more details.
+
+    You should have received a copy of the GNU General Public License
+    along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
+
+\*---------------------------------------------------------------------------*/
+
+#include "dynamicKEqn.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+namespace LESModels
+{
+
+// * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
+
+template<class BasicTurbulenceModel>
+volScalarField dynamicKEqn<BasicTurbulenceModel>::Ck
+(
+    const volSymmTensorField& D,
+    const volScalarField& KK
+) const
+{
+    const volSymmTensorField LL
+    (
+        simpleFilter_(dev(filter_(sqr(this->U_)) - (sqr(filter_(this->U_)))))
+    );
+
+    const volSymmTensorField MM
+    (
+        simpleFilter_(-2.0*this->delta()*sqrt(KK)*filter_(D))
+    );
+
+    const volScalarField Ck
+    (
+        simpleFilter_(0.5*(LL && MM))
+       /(
+            simpleFilter_(magSqr(MM))
+          + dimensionedScalar("small", sqr(MM.dimensions()), VSMALL)
+        )
+    );
+
+    tmp<volScalarField> tfld = 0.5*(mag(Ck) + Ck);
+    return tfld();
+}
+
+
+template<class BasicTurbulenceModel>
+volScalarField dynamicKEqn<BasicTurbulenceModel>::Ce
+(
+    const volSymmTensorField& D,
+    const volScalarField& KK
+) const
+{
+    const volScalarField Ce
+    (
+        simpleFilter_(this->nuEff()*(filter_(magSqr(D)) - magSqr(filter_(D))))
+       /simpleFilter_(pow(KK, 1.5)/(2.0*this->delta()))
+    );
+
+    tmp<volScalarField> tfld = 0.5*(mag(Ce) + Ce);
+    return tfld();
+}
+
+
+template<class BasicTurbulenceModel>
+volScalarField dynamicKEqn<BasicTurbulenceModel>::Ce() const
+{
+    const volSymmTensorField D(dev(symm(fvc::grad(this->U_))));
+
+    volScalarField KK
+    (
+        0.5*(filter_(magSqr(this->U_)) - magSqr(filter_(this->U_)))
+    );
+    KK.max(dimensionedScalar("small", KK.dimensions(), SMALL));
+
+    return Ce(D, KK);
+}
+
+
+template<class BasicTurbulenceModel>
+void dynamicKEqn<BasicTurbulenceModel>::correctNut
+(
+    const volSymmTensorField& D,
+    const volScalarField& KK
+)
+{
+    this->nut_ = Ck(D, KK)*sqrt(k_)*this->delta();
+    this->nut_.correctBoundaryConditions();
+}
+
+
+template<class BasicTurbulenceModel>
+void dynamicKEqn<BasicTurbulenceModel>::correctNut()
+{
+    const volScalarField KK
+    (
+        0.5*(filter_(magSqr(this->U_)) - magSqr(filter_(this->U_)))
+    );
+
+    correctNut(symm(fvc::grad(this->U_)), KK);
+}
+
+
+template<class BasicTurbulenceModel>
+tmp<fvScalarMatrix> dynamicKEqn<BasicTurbulenceModel>::kSource() const
+{
+    return tmp<fvScalarMatrix>
+    (
+        new fvScalarMatrix
+        (
+            k_,
+            dimVolume*this->rho_.dimensions()*k_.dimensions()
+            /dimTime
+        )
+    );
+}
+
+
+// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
+
+template<class BasicTurbulenceModel>
+dynamicKEqn<BasicTurbulenceModel>::dynamicKEqn
+(
+    const alphaField& alpha,
+    const rhoField& rho,
+    const volVectorField& U,
+    const surfaceScalarField& alphaRhoPhi,
+    const surfaceScalarField& phi,
+    const transportModel& transport,
+    const word& propertiesName,
+    const word& type
+)
+:
+    LESeddyViscosity<BasicTurbulenceModel>
+    (
+        type,
+        alpha,
+        rho,
+        U,
+        alphaRhoPhi,
+        phi,
+        transport,
+        propertiesName
+    ),
+
+    k_
+    (
+        IOobject
+        (
+            IOobject::groupName("k", this->U_.group()),
+            this->runTime_.timeName(),
+            this->mesh_,
+            IOobject::MUST_READ,
+            IOobject::AUTO_WRITE
+        ),
+        this->mesh_
+    ),
+
+    simpleFilter_(this->mesh_),
+    filterPtr_(LESfilter::New(this->mesh_, this->coeffDict())),
+    filter_(filterPtr_())
+{
+    bound(k_, this->kMin_);
+
+    if (type == typeName)
+    {
+        correctNut();
+        this->printCoeffs(type);
+    }
+}
+
+
+// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
+
+template<class BasicTurbulenceModel>
+bool dynamicKEqn<BasicTurbulenceModel>::read()
+{
+    if (LESeddyViscosity<BasicTurbulenceModel>::read())
+    {
+        filter_.read(this->coeffDict());
+
+        return true;
+    }
+    else
+    {
+        return false;
+    }
+}
+
+
+template<class BasicTurbulenceModel>
+tmp<volScalarField> dynamicKEqn<BasicTurbulenceModel>::epsilon() const
+{
+    return tmp<volScalarField>
+    (
+        new volScalarField
+        (
+            IOobject
+            (
+                IOobject::groupName("epsilon", this->U_.group()),
+                this->runTime_.timeName(),
+                this->mesh_,
+                IOobject::NO_READ,
+                IOobject::NO_WRITE
+            ),
+            Ce()*k()*sqrt(k())/this->delta()
+        )
+    );
+}
+
+
+template<class BasicTurbulenceModel>
+void dynamicKEqn<BasicTurbulenceModel>::correct()
+{
+    if (!this->turbulence_)
+    {
+        return;
+    }
+
+    // Local references
+    const alphaField& alpha = this->alpha_;
+    const rhoField& rho = this->rho_;
+    const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
+    const volVectorField& U = this->U_;
+    volScalarField& nut = this->nut_;
+
+    LESeddyViscosity<BasicTurbulenceModel>::correct();
+
+    volScalarField divU(fvc::div(fvc::absolute(this->phi(), U)));
+
+    tmp<volTensorField> tgradU(fvc::grad(U));
+    const volSymmTensorField D(dev(symm(tgradU())));
+    const volScalarField G(this->GName(), 2.0*nut*(tgradU() && D));
+    tgradU.clear();
+
+    volScalarField KK(0.5*(filter_(magSqr(U)) - magSqr(filter_(U))));
+    KK.max(dimensionedScalar("small", KK.dimensions(), SMALL));
+
+    tmp<fvScalarMatrix> kEqn
+    (
+        fvm::ddt(alpha, rho, k_)
+      + fvm::div(alphaRhoPhi, k_)
+      - fvm::laplacian(alpha*rho*DkEff(), k_)
+    ==
+        alpha*rho*G
+      - fvm::SuSp((2.0/3.0)*alpha*rho*divU, k_)
+      - fvm::Sp(Ce(D, KK)*alpha*rho*sqrt(k_)/this->delta(), k_)
+      + kSource()
+    );
+
+    kEqn().relax();
+    kEqn().solve();
+    bound(k_, this->kMin_);
+
+    correctNut(D, KK);
+}
+
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+} // End namespace LESModels
+} // End namespace Foam
+
+// ************************************************************************* //
diff --git a/src/TurbulenceModels/turbulenceModels/LES/dynamicKEqn/dynamicKEqn.H b/src/TurbulenceModels/turbulenceModels/LES/dynamicKEqn/dynamicKEqn.H
new file mode 100644
index 0000000000000000000000000000000000000000..4c2ec87683bd803320bc8dde374d435a9e55a68f
--- /dev/null
+++ b/src/TurbulenceModels/turbulenceModels/LES/dynamicKEqn/dynamicKEqn.H
@@ -0,0 +1,210 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | Copyright (C) 2011-2015 OpenFOAM Foundation
+     \\/     M anipulation  |
+-------------------------------------------------------------------------------
+License
+    This file is part of OpenFOAM.
+
+    OpenFOAM is free software: you can redistribute it and/or modify it
+    under the terms of the GNU General Public License as published by
+    the Free Software Foundation, either version 3 of the License, or
+    (at your option) any later version.
+
+    OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
+    ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
+    FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
+    for more details.
+
+    You should have received a copy of the GNU General Public License
+    along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
+
+Class
+     Foam::LESModels::dynamicKEqn
+
+Group
+    grpLESTurbulence
+
+Description
+    Dynamic one equation eddy-viscosity model
+
+    Eddy viscosity SGS model using a modeled balance equation to simulate
+    the behaviour of k in which a dynamic procedure is applied to evaluate the
+    coefficients.
+
+    Reference:
+    \verbatim
+        Kim, W and Menon, S. (1995).
+        A new dynamic one-equation subgrid-scale model for
+        large eddy simulation.
+        In 33rd Aerospace Sciences Meeting and Exhibit, Reno, NV, 1995.
+    \endverbatim
+
+    There are no default model coefficients but the filter used for KK must be
+    supplied, e.g.
+    \verbatim
+        dynamicKEqnCoeffs
+        {
+            filter simple;
+        }
+    \endverbatim
+
+SourceFiles
+    dynamicKEqn.C
+
+\*---------------------------------------------------------------------------*/
+
+#ifndef dynamicKEqn_H
+#define dynamicKEqn_H
+
+#include "LESeddyViscosity.H"
+#include "simpleFilter.H"
+#include "LESfilter.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+namespace LESModels
+{
+
+/*---------------------------------------------------------------------------*\
+                           Class dynamicKEqn Declaration
+\*---------------------------------------------------------------------------*/
+
+template<class BasicTurbulenceModel>
+class dynamicKEqn
+:
+    public LESeddyViscosity<BasicTurbulenceModel>
+{
+    // Private Member Functions
+
+        // Disallow default bitwise copy construct and assignment
+        dynamicKEqn(const dynamicKEqn&);
+        dynamicKEqn& operator=(const dynamicKEqn&);
+
+
+protected:
+
+    // Protected data
+
+        // Fields
+
+            volScalarField k_;
+
+
+        // Filters
+
+            simpleFilter simpleFilter_;
+            autoPtr<LESfilter> filterPtr_;
+            LESfilter& filter_;
+
+
+    // Protected Member Functions
+
+        //- Calculate Ck by filtering the velocity field U
+        volScalarField Ck
+        (
+            const volSymmTensorField& D,
+            const volScalarField& KK
+        ) const;
+
+        //- Calculate Ce by filtering the velocity field U
+        volScalarField Ce
+        (
+            const volSymmTensorField& D,
+            const volScalarField& KK
+        ) const;
+
+        volScalarField Ce() const;
+
+        //- Update sub-grid eddy-viscosity
+        void correctNut
+        (
+            const volSymmTensorField& D,
+            const volScalarField& KK
+        );
+
+        virtual void correctNut();
+
+        virtual tmp<fvScalarMatrix> kSource() const;
+
+
+public:
+
+    typedef typename BasicTurbulenceModel::alphaField alphaField;
+    typedef typename BasicTurbulenceModel::rhoField rhoField;
+    typedef typename BasicTurbulenceModel::transportModel transportModel;
+
+
+    //- Runtime type information
+    TypeName("dynamicKEqn");
+
+
+    // Constructors
+
+        //- Construct from components
+        dynamicKEqn
+        (
+            const alphaField& alpha,
+            const rhoField& rho,
+            const volVectorField& U,
+            const surfaceScalarField& alphaRhoPhi,
+            const surfaceScalarField& phi,
+            const transportModel& transport,
+            const word& propertiesName = turbulenceModel::propertiesName,
+            const word& type = typeName
+        );
+
+
+    //- Destructor
+    virtual ~dynamicKEqn()
+    {}
+
+
+    // Member Functions
+
+        //- Read model coefficients if they have changed
+        virtual bool read();
+
+        //- Return SGS kinetic energy
+        virtual tmp<volScalarField> k() const
+        {
+            return k_;
+        }
+
+        //- Return sub-grid disipation rate
+        virtual tmp<volScalarField> epsilon() const;
+
+        //- Return the effective diffusivity for k
+        tmp<volScalarField> DkEff() const
+        {
+            return tmp<volScalarField>
+            (
+                new volScalarField("DkEff", this->nut_ + this->nu())
+            );
+        }
+
+        //- Correct Eddy-Viscosity and related properties
+        virtual void correct();
+};
+
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+} // End namespace LESModels
+} // End namespace Foam
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#ifdef NoRepository
+#   include "dynamicKEqn.C"
+#endif
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#endif
+
+// ************************************************************************* //
diff --git a/src/TurbulenceModels/turbulenceModels/LES/kEqn/kEqn.C b/src/TurbulenceModels/turbulenceModels/LES/kEqn/kEqn.C
index 168cb477e46c05f6e8e6dd3b8a883181256afb7c..eda73f69294578b6dd9e44720a45daa4e58c1e79 100644
--- a/src/TurbulenceModels/turbulenceModels/LES/kEqn/kEqn.C
+++ b/src/TurbulenceModels/turbulenceModels/LES/kEqn/kEqn.C
@@ -109,6 +109,8 @@ kEqn<BasicTurbulenceModel>::kEqn
         )
     )
 {
+    bound(k_, this->kMin_);
+
     if (type == typeName)
     {
         correctNut();
diff --git a/src/TurbulenceModels/turbulenceModels/LES/kEqn/kEqn.H b/src/TurbulenceModels/turbulenceModels/LES/kEqn/kEqn.H
index bd27700afa6e8e9b48fa30dc1ff9d89407173139..6220fc790f6180dbfc5c219fe4ad61cc89e21bf0 100644
--- a/src/TurbulenceModels/turbulenceModels/LES/kEqn/kEqn.H
+++ b/src/TurbulenceModels/turbulenceModels/LES/kEqn/kEqn.H
@@ -28,7 +28,7 @@ Group
     grpLESTurbulence
 
 Description
-    One Equation Eddy Viscosity Model
+    One equation eddy-viscosity model
 
     Eddy viscosity SGS model using a modeled balance equation to simulate the
     behaviour of k.
diff --git a/tutorials/incompressible/pisoFoam/les/pitzDaily/constant/turbulenceProperties b/tutorials/incompressible/pisoFoam/les/pitzDaily/constant/turbulenceProperties
index 19b9695576877f971ba9f5084614db04322241e8..9c4722d37cfa8963fddc16c6b91364792b161fcb 100644
--- a/tutorials/incompressible/pisoFoam/les/pitzDaily/constant/turbulenceProperties
+++ b/tutorials/incompressible/pisoFoam/les/pitzDaily/constant/turbulenceProperties
@@ -19,7 +19,7 @@ simulationType  LES;
 
 LES
 {
-    LESModel        kEqn;
+    LESModel        dynamicKEqn;
 
     turbulence      on;
 
@@ -27,6 +27,11 @@ LES
 
     delta           cubeRootVol;
 
+    dynamicKEqnCoeffs
+    {
+        filter simple;
+    }
+
     cubeRootVolCoeffs
     {
         deltaCoeff      1;
@@ -89,4 +94,5 @@ LES
     }
 }
 
+
 // ************************************************************************* //
diff --git a/tutorials/incompressible/pisoFoam/les/pitzDaily/system/fvSchemes b/tutorials/incompressible/pisoFoam/les/pitzDaily/system/fvSchemes
index 752f20f660a2d10471ef576996e8adb2c6e7789a..59c8e98676b3a9d7e47ab45ffc8ad39a83591d19 100644
--- a/tutorials/incompressible/pisoFoam/les/pitzDaily/system/fvSchemes
+++ b/tutorials/incompressible/pisoFoam/les/pitzDaily/system/fvSchemes
@@ -28,11 +28,8 @@ gradSchemes
 divSchemes
 {
     default         none;
-    div(phi,U)      Gauss linear;
+    div(phi,U)      Gauss LUST grad(U);
     div(phi,k)      Gauss limitedLinear 1;
-    div(phi,B)      Gauss limitedLinear 1;
-    div(phi,nuTilda) Gauss limitedLinear 1;
-    div(B)          Gauss linear;
     div((nuEff*dev2(T(grad(U))))) Gauss linear;
 }