diff --git a/src/postProcessing/functionObjects/utilities/Make/files b/src/postProcessing/functionObjects/utilities/Make/files
index b47cb83cc1685cd613dbd8b96461681923ad2e49..d5a78589b2fa13a72a42aa49b84babae6c3bb13e 100644
--- a/src/postProcessing/functionObjects/utilities/Make/files
+++ b/src/postProcessing/functionObjects/utilities/Make/files
@@ -9,6 +9,9 @@ dsmcFields/dsmcFieldsFunctionObject.C
 pressureTools/pressureTools.C
 pressureTools/pressureToolsFunctionObject.C
 
+scalarTransport/scalarTransport.C
+scalarTransport/scalarTransportFunctionObject.C
+
 timeActivatedFileUpdate/timeActivatedFileUpdate.C
 timeActivatedFileUpdate/timeActivatedFileUpdateFunctionObject.C
 
diff --git a/src/postProcessing/functionObjects/utilities/Make/options b/src/postProcessing/functionObjects/utilities/Make/options
index 61713e09134be2aac6a3f6db5eeb15df2d08953b..38f9d8a3a07c45599068eb6e7d96d01e32e5bbfc 100644
--- a/src/postProcessing/functionObjects/utilities/Make/options
+++ b/src/postProcessing/functionObjects/utilities/Make/options
@@ -1,5 +1,6 @@
 EXE_INC = \
     -I$(LIB_SRC)/finiteVolume/lnInclude \
+    -I$(LIB_SRC)/fieldSources/lnInclude \
     -I$(LIB_SRC)/lagrangian/basic/lnInclude \
     -I$(LIB_SRC)/lagrangian/dsmc/lnInclude \
     -I$(LIB_SRC)/meshTools/lnInclude \
@@ -13,6 +14,7 @@ EXE_INC = \
 
 LIB_LIBS = \
     -lfiniteVolume \
+    -lfieldSources \
     -lmeshTools \
     -lsampling \
     -llagrangian \
diff --git a/src/postProcessing/functionObjects/utilities/scalarTransport/IOscalarTransport.H b/src/postProcessing/functionObjects/utilities/scalarTransport/IOscalarTransport.H
new file mode 100644
index 0000000000000000000000000000000000000000..37d3ef0f3d61b31b74e254ba84838a06c007bdd6
--- /dev/null
+++ b/src/postProcessing/functionObjects/utilities/scalarTransport/IOscalarTransport.H
@@ -0,0 +1,49 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | Copyright (C) 2012 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/>.
+
+Typedef
+    Foam::IOscalarTransport
+
+Description
+    Instance of the generic IOOutputFilter for scalarTransport.
+
+\*---------------------------------------------------------------------------*/
+
+#ifndef IOscalarTransport_H
+#define IOscalarTransport_H
+
+#include "scalarTransport.H"
+#include "IOOutputFilter.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+    typedef IOOutputFilter<scalarTransport> IOscalarTransport;
+}
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#endif
+
+// ************************************************************************* //
diff --git a/src/postProcessing/functionObjects/utilities/scalarTransport/scalarTransport.C b/src/postProcessing/functionObjects/utilities/scalarTransport/scalarTransport.C
new file mode 100644
index 0000000000000000000000000000000000000000..2d1f866a2a17516358f551ad4a292c3c83d1ee8d
--- /dev/null
+++ b/src/postProcessing/functionObjects/utilities/scalarTransport/scalarTransport.C
@@ -0,0 +1,272 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | Copyright (C) 2012 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 "scalarTransport.H"
+#include "surfaceFields.H"
+#include "dictionary.H"
+#include "fixedValueFvPatchFields.H"
+#include "zeroGradientFvPatchFields.H"
+#include "fvScalarMatrix.H"
+#include "fvmDdt.H"
+#include "fvmDiv.H"
+#include "fvcDiv.H"
+#include "fvmLaplacian.H"
+#include "fvmSup.H"
+#include "incompressible/turbulenceModel/turbulenceModel.H"
+#include "compressible/turbulenceModel/turbulenceModel.H"
+
+// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
+
+defineTypeNameAndDebug(Foam::scalarTransport, 0);
+
+// * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
+
+Foam::wordList Foam::scalarTransport::boundaryTypes() const
+{
+    const volVectorField& U = mesh_.lookupObject<volVectorField>(UName_);
+
+    wordList bTypes(U.boundaryField().size());
+
+    forAll(bTypes, patchI)
+    {
+        const fvPatchField<vector>& pf = U.boundaryField()[patchI];
+        if (isA<fixedValueFvPatchVectorField>(pf))
+        {
+            bTypes[patchI] = fixedValueFvPatchScalarField::typeName;
+        }
+        else
+        {
+            bTypes[patchI] = zeroGradientFvPatchScalarField::typeName;
+        }
+    }
+
+    return bTypes;    
+}
+
+
+Foam::tmp<Foam::volScalarField> Foam::scalarTransport::DT
+(
+    const surfaceScalarField& phi
+) const
+{
+    typedef incompressible::turbulenceModel icoModel;
+    typedef compressible::turbulenceModel cmpModel;
+
+    if (userDT_)
+    {
+        return tmp<volScalarField>
+        (
+            new volScalarField
+            (
+                IOobject
+                (
+                    "DT",
+                    mesh_.time().timeName(),
+                    mesh_.time(),
+                    IOobject::NO_READ,
+                    IOobject::NO_WRITE
+                ),
+                mesh_,
+                dimensionedScalar("DT", phi.dimensions()/dimLength, DT_)
+            )
+        );
+    }
+    else if (mesh_.foundObject<icoModel>("turbulenceModel"))
+    {
+        const icoModel& model = mesh_.lookupObject<icoModel>("turbulenceModel");
+
+        return model.nuEff();
+    }
+    else if (mesh_.foundObject<cmpModel>("turbulenceModel"))
+    {
+        const cmpModel& model = mesh_.lookupObject<cmpModel>("turbulenceModel");
+
+        return model.muEff();
+    }
+    else
+    {
+        return tmp<volScalarField>
+        (
+            new volScalarField
+            (
+                IOobject
+                (
+                    "DT",
+                    mesh_.time().timeName(),
+                    mesh_.time(),
+                    IOobject::NO_READ,
+                    IOobject::NO_WRITE
+                ),
+                mesh_,
+                dimensionedScalar("DT", phi.dimensions()/dimLength, 0.0)
+            )
+        );
+    }
+}
+
+
+// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
+
+Foam::scalarTransport::scalarTransport
+(
+    const word& name,
+    const objectRegistry& obr,
+    const dictionary& dict,
+    const bool loadFromFiles
+)
+:
+    name_(name),
+    mesh_(refCast<const fvMesh>(obr)),
+    active_(true),
+    phiName_("phi"),
+    UName_("U"),
+    DT_(0.0),
+    userDT_(false),
+    resetOnStartUp_(false),
+    nCorr_(0),
+    autoSchemes_(false),
+    sources_(mesh_),
+    T_
+    (
+        IOobject
+        (
+            name,
+            mesh_.time().timeName(),
+            mesh_,
+            IOobject::READ_IF_PRESENT,
+            IOobject::AUTO_WRITE
+        ),
+        mesh_,
+        dimensionedScalar("zero", dimless, 0.0),
+        boundaryTypes()
+    )
+{
+    read(dict);
+
+    if (resetOnStartUp_)
+    {
+        T_ == dimensionedScalar("zero", dimless, 0.0);
+    }
+}
+
+
+// * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
+
+Foam::scalarTransport::~scalarTransport()
+{}
+
+
+// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
+
+void Foam::scalarTransport::read(const dictionary& dict)
+{
+    if (active_)
+    {
+        Info<< type() << ":" << nl;
+
+        phiName_ = dict.lookupOrDefault<word>("phiName", "phi");
+        UName_ = dict.lookupOrDefault<word>("UName", "U");
+
+        userDT_ = false;
+        if (dict.readIfPresent("DT", DT_))
+        {
+            userDT_ = true;
+        }
+
+        dict.lookup("resetOnStartUp") >> resetOnStartUp_;
+
+        dict.readIfPresent("nCorr", nCorr_);
+
+        dict.lookup("autoSchemes") >> autoSchemes_;
+
+        sources_.reset(dict.subDict("sources"));
+    }
+}
+
+
+void Foam::scalarTransport::execute()
+{
+    Info<< type() << " output:" << endl;
+
+    const surfaceScalarField& phi =
+        mesh_.lookupObject<surfaceScalarField>(phiName_);
+
+    // calculate the diffusivity
+    volScalarField DT(this->DT(phi));
+
+    // set schemes
+    word schemeVar = T_.name();
+    if (autoSchemes_)
+    {
+        schemeVar = UName_;
+    }
+    
+    word divScheme("div(phi," + schemeVar + ")");
+    word laplacianScheme("laplacian(" + DT.name() + "," + schemeVar + ")");
+
+    // set under-relaxation coeff
+    scalar relaxCoeff = 0.0;
+    if (mesh_.relaxEquation(schemeVar))
+    {
+        relaxCoeff = mesh_.equationRelaxationFactor(schemeVar);
+    }
+
+    // solve
+    for (label i = 0; i <= nCorr_; i++)
+    {
+        fvScalarMatrix TEqn
+        (
+            fvm::ddt(T_)
+          + fvm::div(phi, T_, divScheme)
+          - fvm::Sp(fvc::div(phi), T_)
+          - fvm::laplacian(DT, T_, laplacianScheme)
+         ==
+            sources_(T_)
+        );
+
+        TEqn.relax(relaxCoeff);
+
+        sources_.constrain(TEqn);
+
+        TEqn.solve(mesh_.solverDict(UName_));
+    }
+
+    Info<< endl;
+}
+
+
+void Foam::scalarTransport::end()
+{
+    // Do nothing
+}
+
+
+void Foam::scalarTransport::write()
+{
+    // Do nothing
+}
+
+
+// ************************************************************************* //
diff --git a/src/postProcessing/functionObjects/utilities/scalarTransport/scalarTransport.H b/src/postProcessing/functionObjects/utilities/scalarTransport/scalarTransport.H
new file mode 100644
index 0000000000000000000000000000000000000000..e08e40943a6de36f5a1686e87c1905020965267f
--- /dev/null
+++ b/src/postProcessing/functionObjects/utilities/scalarTransport/scalarTransport.H
@@ -0,0 +1,188 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | Copyright (C) 2012 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::scalarTransport
+
+Description
+    This function object evolves a passive scalar transport equation.  The
+    field in ininitially zero, to which sources are added.  The field name
+    is assigned the name of the function object.  Boundary conditions are
+    automatically applied, based on the velocity boundary conditions.
+
+    - the field can be zeroed on start-up using the resetOnStartUp flag
+    - to employ the same numerical schemes as the flow velocity, use the
+      autoSchemes flag
+    - the diffusivity can be set manually using the DT entry, or retrieved
+      from the turbulence model (if applicable)
+
+SourceFiles
+    scalarTransport.C
+    IOscalarTransport.H
+
+SeeAlso
+    Foam::basicSourceList
+
+\*---------------------------------------------------------------------------*/
+
+#ifndef scalarTransport_H
+#define scalarTransport_H
+
+#include "volFields.H"
+#include "surfaceFieldsFwd.H"
+#include "pointFieldFwd.H"
+#include "fvMatricesFwd.H"
+#include "basicSourceList.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+
+// Forward declaration of classes
+class objectRegistry;
+class dictionary;
+class mapPolyMesh;
+
+/*---------------------------------------------------------------------------*\
+                       Class scalarTransport Declaration
+\*---------------------------------------------------------------------------*/
+
+class scalarTransport
+{
+    // Private data
+
+        //- Name of this set of scalarTransport objects
+        word name_;
+
+        //- Reference to the mesh database
+        const fvMesh& mesh_;
+
+        //- On/off switch
+        bool active_;
+
+        //- Name of flux field (optional)
+        word phiName_;
+
+        //- Name of velocity field (optional)
+        word UName_;
+
+        //- Diffusion coefficient (optional)
+        scalar DT_;
+
+        //- Flag to indicate whether user DT_ is used
+        bool userDT_;
+
+        //- Flag to reset scalar field on start-up
+        bool resetOnStartUp_;
+
+        //- Number of corrector iterations (optional)
+        label nCorr_;
+
+        //- Flag to employ schemes for velocity for scalar transport
+        bool autoSchemes_;
+
+        //- Run-time selectable sources
+        basicSourceList sources_;
+
+        //- The scalar field
+        volScalarField T_;
+
+
+    // Private Member Functions
+
+        //- Return the boundary types for the scalar field
+        wordList boundaryTypes() const;
+
+        //- Return the diffusivity field
+        tmp<volScalarField> DT(const surfaceScalarField& phi) const;
+
+        //- Disallow default bitwise copy construct
+        scalarTransport(const scalarTransport&);
+
+        //- Disallow default bitwise assignment
+        void operator=(const scalarTransport&);
+
+
+public:
+
+    //- Runtime type information
+    TypeName("scalarTransport");
+
+
+    // Constructors
+
+        //- Construct for given objectRegistry and dictionary.
+        //  Allow the possibility to load fields from files
+        scalarTransport
+        (
+            const word& name,
+            const objectRegistry&,
+            const dictionary&,
+            const bool loadFromFiles = false
+        );
+
+
+    //- Destructor
+    virtual ~scalarTransport();
+
+
+    // Member Functions
+
+        //- Return name of the set of scalarTransport
+        virtual const word& name() const
+        {
+            return name_;
+        }
+
+        //- Read the scalarTransport data
+        virtual void read(const dictionary&);
+
+        //- Execute, currently does nothing
+        virtual void execute();
+
+        //- Execute at the final time-loop, currently does nothing
+        virtual void end();
+
+        //- Calculate the scalarTransport and write
+        virtual void write();
+
+        //- Update for changes of mesh
+        virtual void updateMesh(const mapPolyMesh&)
+        {}
+
+        //- Update for changes of mesh
+        virtual void movePoints(const pointField&)
+        {}
+};
+
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+} // End namespace Foam
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#endif
+
+// ************************************************************************* //
diff --git a/src/postProcessing/functionObjects/utilities/scalarTransport/scalarTransportFunctionObject.C b/src/postProcessing/functionObjects/utilities/scalarTransport/scalarTransportFunctionObject.C
new file mode 100644
index 0000000000000000000000000000000000000000..3f9ba1b3ad410d031a6e811ec28a16e6c6328a35
--- /dev/null
+++ b/src/postProcessing/functionObjects/utilities/scalarTransport/scalarTransportFunctionObject.C
@@ -0,0 +1,42 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | Copyright (C) 2012 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 "scalarTransportFunctionObject.H"
+
+// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
+
+namespace Foam
+{
+    defineNamedTemplateTypeNameAndDebug(scalarTransportFunctionObject, 0);
+
+    addToRunTimeSelectionTable
+    (
+        functionObject,
+        scalarTransportFunctionObject,
+        dictionary
+    );
+}
+
+// ************************************************************************* //
diff --git a/src/postProcessing/functionObjects/utilities/scalarTransport/scalarTransportFunctionObject.H b/src/postProcessing/functionObjects/utilities/scalarTransport/scalarTransportFunctionObject.H
new file mode 100644
index 0000000000000000000000000000000000000000..f5a417f4fdec04a7d24bf11989c1ba513ab478e7
--- /dev/null
+++ b/src/postProcessing/functionObjects/utilities/scalarTransport/scalarTransportFunctionObject.H
@@ -0,0 +1,54 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | Copyright (C) 2012 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/>.
+
+Typedef
+    Foam::scalarTransportFunctionObject
+
+Description
+    FunctionObject wrapper around scalarTransport to allow it to be
+    created via the functions entry within controlDict.
+
+SourceFiles
+    scalarTransportFunctionObject.C
+
+\*---------------------------------------------------------------------------*/
+
+#ifndef scalarTransportFunctionObject_H
+#define scalarTransportFunctionObject_H
+
+#include "scalarTransport.H"
+#include "OutputFilterFunctionObject.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+    typedef OutputFilterFunctionObject<scalarTransport>
+        scalarTransportFunctionObject;
+}
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#endif
+
+// ************************************************************************* //