diff --git a/src/TurbulenceModels/compressible/turbulentFluidThermoModels/turbulentFluidThermoModels.C b/src/TurbulenceModels/compressible/turbulentFluidThermoModels/turbulentFluidThermoModels.C
index 519e15be4a..accb8ab287 100644
--- a/src/TurbulenceModels/compressible/turbulentFluidThermoModels/turbulentFluidThermoModels.C
+++ b/src/TurbulenceModels/compressible/turbulentFluidThermoModels/turbulentFluidThermoModels.C
@@ -72,6 +72,9 @@ makeRASModel(kOmegaSST);
#include "Smagorinsky.H"
makeLESModel(Smagorinsky);
+#include "WALE.H"
+makeLESModel(WALE);
+
#include "kEqn.H"
makeLESModel(kEqn);
diff --git a/src/TurbulenceModels/incompressible/turbulentTransportModels/turbulentTransportModels.C b/src/TurbulenceModels/incompressible/turbulentTransportModels/turbulentTransportModels.C
index c881910b36..508288db76 100644
--- a/src/TurbulenceModels/incompressible/turbulentTransportModels/turbulentTransportModels.C
+++ b/src/TurbulenceModels/incompressible/turbulentTransportModels/turbulentTransportModels.C
@@ -67,6 +67,9 @@ makeRASModel(kOmegaSST);
#include "Smagorinsky.H"
makeLESModel(Smagorinsky);
+#include "WALE.H"
+makeLESModel(WALE);
+
#include "kEqn.H"
makeLESModel(kEqn);
diff --git a/src/TurbulenceModels/turbulenceModels/LES/WALE/WALE.C b/src/TurbulenceModels/turbulenceModels/LES/WALE/WALE.C
new file mode 100644
index 0000000000..40cf764d97
--- /dev/null
+++ b/src/TurbulenceModels/turbulenceModels/LES/WALE/WALE.C
@@ -0,0 +1,205 @@
+/*---------------------------------------------------------------------------*\
+ ========= |
+ \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
+ \\ / O peration |
+ \\ / A nd | Copyright (C) 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 .
+
+\*---------------------------------------------------------------------------*/
+
+#include "WALE.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+namespace LESModels
+{
+
+// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
+
+template
+tmp WALE::Sd
+(
+ const volTensorField& gradU
+) const
+{
+ return dev(symm(gradU & gradU));
+}
+
+
+template
+tmp WALE::k
+(
+ const volTensorField& gradU
+) const
+{
+ volScalarField magSqrSd(magSqr(Sd(gradU)));
+
+ return tmp
+ (
+ new volScalarField
+ (
+ IOobject
+ (
+ IOobject::groupName("k", this->U_.group()),
+ this->runTime_.timeName(),
+ this->mesh_
+ ),
+ sqr(sqr(Cw_)*this->delta()/Ck_)*
+ (
+ pow3(magSqrSd)
+ /(
+ sqr
+ (
+ pow(magSqr(symm(gradU)), 5.0/2.0)
+ + pow(magSqrSd, 5.0/4.0)
+ )
+ + dimensionedScalar
+ (
+ "SMALL",
+ dimensionSet(0, 0, -10, 0, 0),
+ SMALL
+ )
+ )
+ )
+ )
+ );
+}
+
+
+template
+void WALE::correctNut()
+{
+ this->nut_ = Ck_*this->delta()*sqrt(this->k(fvc::grad(this->U_)));
+ this->nut_.correctBoundaryConditions();
+}
+
+
+// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
+
+template
+WALE::WALE
+(
+ 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
+ (
+ type,
+ alpha,
+ rho,
+ U,
+ alphaRhoPhi,
+ phi,
+ transport,
+ propertiesName
+ ),
+
+ Ck_
+ (
+ dimensioned::lookupOrAddToDict
+ (
+ "Ck",
+ this->coeffDict_,
+ 0.094
+ )
+ ),
+
+ Cw_
+ (
+ dimensioned::lookupOrAddToDict
+ (
+ "Cw",
+ this->coeffDict_,
+ 0.325
+ )
+ )
+{
+ if (type == typeName)
+ {
+ correctNut();
+ this->printCoeffs(type);
+ }
+}
+
+
+// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
+
+template
+bool WALE::read()
+{
+ if (LESeddyViscosity::read())
+ {
+ Ck_.readIfPresent(this->coeffDict());
+ Cw_.readIfPresent(this->coeffDict());
+
+ return true;
+ }
+ else
+ {
+ return false;
+ }
+}
+
+
+template
+tmp WALE::epsilon() const
+{
+ volScalarField k(this->k(fvc::grad(this->U_)));
+
+ return tmp
+ (
+ new volScalarField
+ (
+ IOobject
+ (
+ IOobject::groupName("epsilon", this->U_.group()),
+ this->runTime_.timeName(),
+ this->mesh_,
+ IOobject::NO_READ,
+ IOobject::NO_WRITE
+ ),
+ this->Ce_*k*sqrt(k)/this->delta()
+ )
+ );
+}
+
+
+template
+void WALE::correct()
+{
+ LESeddyViscosity::correct();
+ correctNut();
+}
+
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+} // End namespace LESModels
+} // End namespace Foam
+
+// ************************************************************************* //
diff --git a/src/TurbulenceModels/turbulenceModels/LES/WALE/WALE.H b/src/TurbulenceModels/turbulenceModels/LES/WALE/WALE.H
new file mode 100644
index 0000000000..7cae23f7da
--- /dev/null
+++ b/src/TurbulenceModels/turbulenceModels/LES/WALE/WALE.H
@@ -0,0 +1,173 @@
+/*---------------------------------------------------------------------------*\
+ ========= |
+ \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
+ \\ / O peration |
+ \\ / A nd | Copyright (C) 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 .
+
+Class
+ Foam::LESModels::WALE
+
+Group
+ grpLESTurbulence
+
+Description
+ The Wall-adapting local eddy-viscosity (WALE) SGS model.
+
+ Reference:
+ \verbatim
+ Nicoud, F., & Ducros, F. (1999).
+ Subgrid-scale stress modelling based on the square of the velocity
+ gradient tensor.
+ Flow, Turbulence and Combustion, 62(3), 183-200.
+ \endverbatim
+
+ The default model coefficients correspond to the following:
+ \verbatim
+ WALECoeffs
+ {
+ Ck 0.094;
+ Ce 1.048;
+ Cw 0.325;
+ }
+ \endverbatim
+
+SourceFiles
+ WALE.C
+
+\*---------------------------------------------------------------------------*/
+
+#ifndef WALE_H
+#define WALE_H
+
+#include "LESModel.H"
+#include "LESeddyViscosity.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+namespace LESModels
+{
+
+/*---------------------------------------------------------------------------*\
+ Class WALE Declaration
+\*---------------------------------------------------------------------------*/
+
+template
+class WALE
+:
+ public LESeddyViscosity
+{
+ // Private Member Functions
+
+ // Disallow default bitwise copy construct and assignment
+ WALE(const WALE&);
+ WALE& operator=(const WALE&);
+
+
+protected:
+
+ // Protected data
+
+ dimensionedScalar Ck_;
+ dimensionedScalar Cw_;
+
+
+ // Protected Member Functions
+
+ //- Return the deviatoric symmetric part of the square of the given
+ // velocity gradient field
+ tmp Sd(const volTensorField& gradU) const;
+
+ //- Return SGS kinetic energy
+ // calculated from the given velocity gradient
+ tmp k(const volTensorField& gradU) const;
+
+ //- Update the SGS eddy viscosity
+ virtual void correctNut();
+
+
+public:
+
+ typedef typename BasicTurbulenceModel::alphaField alphaField;
+ typedef typename BasicTurbulenceModel::rhoField rhoField;
+ typedef typename BasicTurbulenceModel::transportModel transportModel;
+
+
+ //- Runtime type information
+ TypeName("WALE");
+
+
+ // Constructors
+
+ //- Construct from components
+ WALE
+ (
+ 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 ~WALE()
+ {}
+
+
+ // Member Functions
+
+ //- Read model coefficients if they have changed
+ virtual bool read();
+
+ //- Return SGS kinetic energy
+ virtual tmp k() const
+ {
+ return k(fvc::grad(this->U_));
+ }
+
+ //- Return sub-grid disipation rate
+ virtual tmp epsilon() const;
+
+ //- Correct Eddy-Viscosity and related properties
+ virtual void correct();
+};
+
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+} // End namespace LESModels
+} // End namespace Foam
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#ifdef NoRepository
+# include "WALE.C"
+#endif
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#endif
+
+// ************************************************************************* //
diff --git a/tutorials/incompressible/pimpleFoam/channel395/constant/turbulenceProperties b/tutorials/incompressible/pimpleFoam/channel395/constant/turbulenceProperties
index b04390b6e2..1f1772f0e4 100644
--- a/tutorials/incompressible/pimpleFoam/channel395/constant/turbulenceProperties
+++ b/tutorials/incompressible/pimpleFoam/channel395/constant/turbulenceProperties
@@ -19,13 +19,13 @@ simulationType LES;
LES
{
- LESModel SpalartAllmarasDDES; //kEqn;
+ LESModel WALE;
turbulence on;
printCoeffs on;
- delta vanDriest;
+ delta cubeRootVol;
cubeRootVolCoeffs
{