diff --git a/src/postProcessing/functionObjects/utilities/ddt2/ddt2.C b/src/postProcessing/functionObjects/utilities/ddt2/ddt2.C
index 76975fcd070ed59041923cf1e63201c3442feb6f..44fc01630d49eb294466f13777ed55e7c71effc6 100644
--- a/src/postProcessing/functionObjects/utilities/ddt2/ddt2.C
+++ b/src/postProcessing/functionObjects/utilities/ddt2/ddt2.C
@@ -29,6 +29,7 @@ License
 #include "dictionary.H"
 #include "FieldFunctions.H"
 #include "fvcDdt.H"
+#include "steadyStateDdtScheme.H"
 
 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
 
@@ -39,76 +40,198 @@ namespace Foam
 
 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
 
-
-template<class FieldType>
-bool Foam::ddt2::calculate
-(
-    const fvMesh& mesh,
-    bool& done
-)
+bool Foam::ddt2::checkFormatName(const word& str)
 {
-    if (done)
+    if (str.find("@@") == string::npos)
     {
-        return true; // already done - skip
-    }
+        WarningInFunction
+            << "Bad result naming "
+            << "(no '@@' token found), deactivating."
+            << nl << endl;
 
-    done = mesh.foundObject<FieldType>(fieldName_);
-    if (!done)
+        return false;
+    }
+    else if (str == "@@")
     {
+        WarningInFunction
+            << "Bad result naming "
+            << "(only a '@@' token found), deactivating."
+            << nl
+            << endl;
+
         return false;
     }
+    else
+    {
+        return true;
+    }
+}
 
-    const FieldType& input =
-        mesh.lookupObject<FieldType>(fieldName_);
 
-    if (!mesh.foundObject<volScalarField>(resultName_))
+void Foam::ddt2::uniqWords(wordReList& lst)
+{
+    boolList retain(lst.size());
+    wordHashSet uniq;
+    forAll(lst, i)
     {
-        const dimensionSet dims
+        const wordRe& select = lst[i];
+
+        retain[i] =
         (
-            mag_
-          ? mag(input.dimensions()/dimTime)
-          : magSqr(input.dimensions()/dimTime)
+            select.isPattern()
+         || uniq.insert(static_cast<const word&>(select))
         );
+    }
+
+    inplaceSubset(retain, lst);
+}
+
+
+bool Foam::ddt2::accept(const word& fieldName) const
+{
+    // check input vs possible result names
+    // to avoid circular calculations
+    return !blacklist_.match(fieldName);
+}
+
+
+template<class FieldType>
+int Foam::ddt2::apply
+(
+    const fvMesh& mesh,
+    const word& inputName,
+    int& state
+)
+{
+    // state: return 0 (not-processed), -1 (skip), +1 ok
+
+    // already done, or not available
+    if (state || !mesh.foundObject<FieldType>(inputName))
+    {
+        return state;
+    }
+
+    const FieldType& input = mesh.lookupObject<FieldType>(inputName);
+
+    word outputName(resultName_);
+    outputName.replace("@@", inputName);
+
+    results_.set(outputName);
 
+    if (!mesh.foundObject<volScalarField>(outputName))
+    {
         mesh.objectRegistry::store
         (
             new volScalarField
             (
                 IOobject
                 (
-                    resultName_,
+                    outputName,
                     mesh.time().timeName(),
                     mesh,
                     IOobject::NO_READ,
-                    IOobject::NO_WRITE // OR IOobject::AUTO_WRITE
+                    IOobject::NO_WRITE
                 ),
                 mesh,
                 dimensionedScalar
                 (
-                    "zero",
-                    dims,
+                    "0",
+                    (
+                        mag_
+                      ? mag(input.dimensions()/dimTime)
+                      : magSqr(input.dimensions()/dimTime)
+                    ),
                     Zero
-                ),
-                emptyPolyPatch::typeName
+                )
             )
         );
     }
 
-    volScalarField& field = const_cast<volScalarField&>
+    volScalarField& output = const_cast<volScalarField&>
     (
-        mesh.lookupObject<volScalarField>(resultName_)
+        mesh.lookupObject<volScalarField>(outputName)
     );
 
     if (mag_)
     {
-        field = mag(fvc::ddt(input));
+        output = mag(fvc::ddt(input));
     }
     else
     {
-        field = magSqr(fvc::ddt(input));
+        output = magSqr(fvc::ddt(input));
+    }
+
+    if (log_)
+    {
+        // could add additional statistics here
+        Info<< type() << " " << name_
+            << " field " << output
+            << " average: " << gAverage(output) << endl;
+    }
+
+    state = +1;
+    return state;
+}
+
+
+int Foam::ddt2::process(const fvMesh& mesh, const word& fieldName)
+{
+    if (!accept(fieldName))
+    {
+        return -1;
+    }
+
+    int state = 0;
+
+    apply<volScalarField>(mesh, fieldName, state);
+    apply<volVectorField>(mesh, fieldName, state);
+
+    return state;
+}
+
+
+void Foam::ddt2::process(const fvMesh& mesh)
+{
+    results_.clear();
+
+    wordHashSet candidates = subsetStrings(selectFields_, mesh.names());
+    DynamicList<word> missing(selectFields_.size());
+    DynamicList<word> ignored(selectFields_.size());
+
+    // check exact matches first
+    forAll(selectFields_, i)
+    {
+        const wordRe& select = selectFields_[i];
+        if (!select.isPattern())
+        {
+            const word& fieldName = static_cast<const word&>(select);
+
+            if (!candidates.erase(fieldName))
+            {
+                missing.append(fieldName);
+            }
+            else if (process(mesh, fieldName) < 1)
+            {
+                ignored.append(fieldName);
+            }
+        }
     }
 
-    return done;
+    forAllConstIter(wordHashSet, candidates, iter)
+    {
+        process(mesh, iter.key());
+    }
+
+    if (missing.size())
+    {
+        WarningInFunction
+            << "Missing field " << missing << endl;
+    }
+    if (ignored.size())
+    {
+        WarningInFunction
+            << "Unprocessed field " << ignored << endl;
+    }
 }
 
 
@@ -125,13 +248,28 @@ Foam::ddt2::ddt2
     name_(name),
     obr_(obr),
     active_(true),
-    fieldName_("undefined-fieldName"),
+    selectFields_(),
     resultName_(word::null),
+    blacklist_(),
+    results_(),
     log_(true),
     mag_(dict.lookupOrDefault<Switch>("mag", false))
 {
-    // Check if the available mesh is an fvMesh, otherwise deactivate
-    if (!isA<fvMesh>(obr_))
+    // Check if the available mesh is an fvMesh and transient,
+    // otherwise deactivate
+    if (isA<fvMesh>(obr_))
+    {
+        const fvMesh& mesh = refCast<const fvMesh>(obr_);
+
+        if (word(mesh.ddtScheme("default")) == "steadyState")
+        {
+            active_ = false;
+            WarningInFunction
+                << "Steady-state, deactivating." << nl
+                << endl;
+        }
+    }
+    else
     {
         active_ = false;
         WarningInFunction
@@ -157,76 +295,62 @@ void Foam::ddt2::read(const dictionary& dict)
     {
         log_.readIfPresent("log", dict);
 
-        dict.lookup("fieldName") >> fieldName_;
-        dict.readIfPresent("resultName", resultName_);
+        dict.lookup("fields") >> selectFields_;
+        uniqWords(selectFields_);
+
+        resultName_ = dict.lookupOrDefault<word>
+        (
+            "result",
+            ( mag_ ? "mag(ddt(@@))" : "magSqr(ddt(@@))" )
+        );
 
-        if (resultName_.empty())
+        active_ = checkFormatName(resultName_);
+        if (active_)
         {
-            resultName_ =
+            blacklist_.set
             (
-                word(mag_ ? "mag" : "magSqr")
-              + "(ddt(" + fieldName_ + "))"
+                string::quotemeta<regExp>
+                (
+                    resultName_
+                ).replace("@@", "(.+)")
             );
         }
+        else
+        {
+            blacklist_.clear();
+        }
     }
 }
 
 
 void Foam::ddt2::execute()
 {
+    results_.clear();
     if (active_)
     {
-        const fvMesh& mesh = refCast<const fvMesh>(obr_);
-        bool done = false;
-
-        calculate<volScalarField>(mesh, done);
-        calculate<volVectorField>(mesh, done);
-
-        if (!done)
-        {
-            WarningInFunction
-                << "Unprocessed field " << fieldName_ << endl;
-        }
+        process(refCast<const fvMesh>(obr_));
     }
 }
 
 
-void Foam::ddt2::end()
-{
-    // Do nothing
-}
-
-
-void Foam::ddt2::timeSet()
-{
-    // Do nothing
-}
-
-
 void Foam::ddt2::write()
 {
     if (active_)
     {
-        if (obr_.foundObject<regIOobject>(resultName_))
+        // consistent output order
+        const wordList outputList = results_.sortedToc();
+
+        forAll(outputList, i)
         {
-            const regIOobject& io =
-                obr_.lookupObject<regIOobject>(resultName_);
+            const word& fieldName = outputList[i];
 
-            if (log_)
+            if (obr_.foundObject<regIOobject>(fieldName))
             {
-                const volScalarField& field = dynamicCast<const volScalarField&>
-                (
-                    io
-                );
+                const regIOobject& io =
+                    obr_.lookupObject<regIOobject>(fieldName);
 
-                // could add additional statistics here
-                Info<< type() << " " << name_
-                    << " output: writing field " << field.name()
-                    << " average: " << gAverage(field) << endl;
+                io.write();
             }
-
-            // could also add option to suppress writing?
-            io.write();
         }
     }
 }
diff --git a/src/postProcessing/functionObjects/utilities/ddt2/ddt2.H b/src/postProcessing/functionObjects/utilities/ddt2/ddt2.H
index 93ee8ea077f3e3137431eecfab3bf26473ea2dc5..acbcc5d47cc902394fa1aa4d61b96cffa2dbc7af 100644
--- a/src/postProcessing/functionObjects/utilities/ddt2/ddt2.H
+++ b/src/postProcessing/functionObjects/utilities/ddt2/ddt2.H
@@ -40,8 +40,8 @@ Description
     {
         type        ddt2;
         functionObjectLibs ("libutilityFunctionObjects.so");
-        fieldName   p;
-        resultName  dpdt2;
+        fields      (p);
+        result      d@@dt2;
         ...
     }
     \endverbatim
@@ -50,8 +50,8 @@ Description
     \table
         Property | Description                | Required  | Default value
         type     | type name: ddt2            | yes       |
-        fieldName | Name of field to process  | yes       |
-        resultName | Name of magnitude field  | no        | magSqr(ddt(fieldName))
+        fields   | Name of fields to process  | yes       |
+        result   | Name of results            | no        | magSqr(ddt(@@))
         log      | Log to standard output     | no        | yes
         mag      | Use 'mag' instead of 'magSqr' | no     | false
     \endtable
@@ -59,6 +59,13 @@ Description
     Note that the optional 'mag' entry cannot be changed during the simulation
     since it alters the dimensions of the output field.
 
+    A list of fields can contain exact names or regular expressions.
+    The token '\@\@' in the result name is replaced by the name of the source
+    field.
+
+    The function object will skip over fields that appear to have
+    already been processed (ie, their names are similar to the output names).
+
 SourceFiles
     ddt2.C
     IOddt2.H
@@ -69,9 +76,10 @@ SourceFiles
 #define ddt2_H
 
 #include "volFieldsFwd.H"
-#include "surfaceFieldsFwd.H"
-#include "pointFieldFwd.H"
 #include "OFstream.H"
+#include "wordReList.H"
+#include "regExp.H"
+#include "HashSet.H"
 #include "Switch.H"
 
 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
@@ -85,7 +93,6 @@ class dictionary;
 class fvMesh;
 class polyMesh;
 class mapPolyMesh;
-class dimensionSet;
 
 /*---------------------------------------------------------------------------*\
                             Class ddt2 Declaration
@@ -104,12 +111,18 @@ class ddt2
         //- On/off switch
         bool active_;
 
-        //- Name of field to process
-        word fieldName_;
+        //- Name of fields to process
+        wordReList selectFields_;
 
-        //- Name of result field
+        //- Formatting for the result fields.
         word resultName_;
 
+        //- Avoid processing the same field twice.
+        mutable regExp blacklist_;
+
+        //- Hashed names of result fields.
+        wordHashSet results_;
+
         //- Switch to send output to Info as well as to file
         Switch log_;
 
@@ -121,21 +134,25 @@ class ddt2
 
     // Private Member Functions
 
-        //- Create result field upon demand
-        volScalarField& getOrCreateResultField
-        (
-            const fvMesh&,
-            const dimensionSet& inputDims
-        );
+        //- Check that the word contains the appropriate substitution token(s).
+        static bool checkFormatName(const word&);
 
+        //- Eliminate duplicate 'word' entries
+        static void uniqWords(wordReList&);
 
-        //- Create result field upon demand
-        template<class FieldType>
-        bool calculate(const FieldType& input);
 
-        //- Create result field upon demand
+        //- Accept unless field name appears to have already been processed
+        bool accept(const word& fieldName) const;
+
+        //- Apply for the volume field type
         template<class FieldType>
-        bool calculate(const fvMesh&, bool& done);
+        int apply(const fvMesh& mesh, const word& inputName, int& state);
+
+        //- Process by trying to apply for various volume field types.
+        int process(const fvMesh& mesh, const word& inputName);
+
+        //- Calculate the ddt2 fields
+        void process(const fvMesh& mesh);
 
 
         //- Disallow default bitwise copy construct
@@ -179,16 +196,18 @@ public:
         //- Read the ddt2 specification
         virtual void read(const dictionary&);
 
-        //- Execute, currently does nothing
+        //- Calculate the ddt2 fields
         virtual void execute();
 
         //- Execute at the final time-loop, currently does nothing
-        virtual void end();
+        virtual void end()
+        {}
 
         //- Called when time was set at the end of the Time::operator++
-        virtual void timeSet();
+        virtual void timeSet()
+        {}
 
-        //- Calculate the ddt2 and write
+        //- Write the ddt fields
         virtual void write();
 
         //- Update for changes of mesh