diff --git a/src/functionObjects/solvers/scalarTransport/scalarTransport.C b/src/functionObjects/solvers/scalarTransport/scalarTransport.C
index 954a5d6..4a1d996 100644
--- a/src/functionObjects/solvers/scalarTransport/scalarTransport.C
+++ b/src/functionObjects/solvers/scalarTransport/scalarTransport.C
@@ -82,42 +82,44 @@ Foam::tmp<Foam::volScalarField> Foam::functionObjects::scalarTransport::D
             )
         );
     }
-    else if (mesh_.foundObject<icoModel>(turbulenceModel::propertiesName))
+    else
     {
-        const icoModel& model = mesh_.lookupObject<icoModel>
-        (
-            turbulenceModel::propertiesName
-        );
+        if (phi.dimensions() == dimMass/dimTime)
+        {
+            const volScalarField& mu =
+                mesh_.lookupObject<volScalarField>(muName_);
 
-        return model.nuEff();
-    }
-    else if (mesh_.foundObject<cmpModel>(turbulenceModel::propertiesName))
-    {
-        const cmpModel& model = mesh_.lookupObject<cmpModel>
-        (
-            turbulenceModel::propertiesName
-        );
+            tmp<volScalarField> D(new volScalarField(Dname, mu*alphaD_));
 
-        return model.muEff();
-    }
-    else
-    {
-        return tmp<volScalarField>
-        (
-            new volScalarField
-            (
-                IOobject
+            if (mesh_.foundObject<cmpModel>(turbulenceModel::propertiesName))
+            {
+                const cmpModel& model = mesh_.lookupObject<cmpModel>
                 (
-                    Dname,
-                    mesh_.time().timeName(),
-                    mesh_.time(),
-                    IOobject::NO_READ,
-                    IOobject::NO_WRITE
-                ),
-                mesh_,
-                dimensionedScalar(Dname, phi.dimensions()/dimLength, 0.0)
-            )
-        );
+                    turbulenceModel::propertiesName
+                );
+                D = D + model.mut()*alphaDt_;
+            }
+
+            return D;
+        }
+        else
+        {
+            const volScalarField& nu =
+                mesh_.lookupObject<volScalarField>(nuName_);
+
+            tmp<volScalarField> D(new volScalarField(Dname, nu*alphaD_));
+
+            if (mesh_.foundObject<icoModel>(turbulenceModel::propertiesName))
+            {
+                const icoModel& model = mesh_.lookupObject<icoModel>
+                (
+                    turbulenceModel::propertiesName
+                );
+                D = D + model.nut()*alphaDt_;
+            }
+
+            return D;
+        }
     }
 }
 
@@ -167,11 +169,45 @@ bool Foam::functionObjects::scalarTransport::read(const dictionary& dict)
 
     phiName_ = dict.lookupOrDefault<word>("phi", "phi");
     rhoName_ = dict.lookupOrDefault<word>("rho", "rho");
+    muName_  = dict.lookupOrDefault<word>("mu", "thermo:mu");
+    nuName_  = dict.lookupOrDefault<word>("nu", "nu");
     schemesField_ = dict.lookupOrDefault<word>("schemesField", fieldName_);
 
+    bool alphaDPresent = false;
+    if (dict.readIfPresent("alphaD", alphaD_))
+    {
+        alphaDPresent = true;
+    }
+
+    if (dict.readIfPresent("alphaDt", alphaDt_))
+    {
+        if (!alphaDPresent)
+        {
+            FatalErrorInFunction
+                << "Laminar alphaD must be specified in dictionary "
+                << dict.name() << exit(FatalError);
+        }
+    }
+    // Set turbulent alphaDt to 0 if it is not specified
+    else
+    {
+        alphaDt_ = 0.0;
+    }
+
     constantD_ = false;
     if (dict.readIfPresent("D", D_))
     {
+        if (alphaDPresent)
+        {
+            FatalErrorInFunction
+                << "Both D and alphaD are specified in dictionary "
+                << dict.name() << exit(FatalError);
+        }
+        constantD_ = true;
+    }
+    // If D and alphaD are not specified then D=0 is applied
+    else if (!alphaDPresent)
+    {
         constantD_ = true;
     }
 
diff --git a/src/functionObjects/solvers/scalarTransport/scalarTransport.H b/src/functionObjects/solvers/scalarTransport/scalarTransport.H
index 86c92bd..7992419 100644
--- a/src/functionObjects/solvers/scalarTransport/scalarTransport.H
+++ b/src/functionObjects/solvers/scalarTransport/scalarTransport.H
@@ -2,7 +2,7 @@
   =========                 |
   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
    \\    /   O peration     |
-    \\  /    A nd           | Copyright (C) 2012-2016 OpenFOAM Foundation
+    \\  /    A nd           | Copyright (C) 2012-2017 OpenFOAM Foundation
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -33,8 +33,11 @@ Description
     - To specify the field name set the 'field' entry
     - To employ the same numerical schemes as another field set
       the 'schemesField' entry,
-    - The diffusivity can be set manually using the 'D' entry, or retrieved
-      from the turbulence model (if applicable).
+    - The constant diffusivity can be set manually using the 'D' entry
+    - The 'alphaD' and 'alphaDt' entries can be specified instead of 'D';
+      these are inverse Schmidt numbers (laminar and turbulent) defined as:
+      alphaD = 1/Sc, alphaDt = 1/Sct.
+      So diffusivity is evaluated as: D = nu*alphaD + nut*alphaDt
 
 See also
     Foam::functionObjects::fvMeshFunctionObject
@@ -77,12 +80,24 @@ class scalarTransport
         //- Name of density field (optional)
         word rhoName_;
 
+        //- Name of dynamic viscosity field (optional)
+        word muName_;
+
+        //- Name of kinematic viscosity field (optional)
+        word nuName_;
+
         //- Diffusion coefficient (optional)
         scalar D_;
 
         //- Flag to indicate whether a constant, uniform D_ is specified
         bool constantD_;
 
+        //- Inverse laminar Schmidt number alphaD = 1/Sc (optional)
+        scalar alphaD_;
+
+        //- Inverse turbulent Schmidt number alphaDt = 1/Sct (optional)
+        scalar alphaDt_;
+
         //- Number of corrector iterations (optional)
         label nCorr_;
 
