From c07997d816b8dba1ea6926f59760c0b948c5ea16 Mon Sep 17 00:00:00 2001
From: Mattijs Janssens <ext-mjanssens@esi-group.com>
Date: Thu, 30 May 2024 10:25:24 +0000
Subject: [PATCH] Function1: added 'lookup' Function1

---
 src/OpenFOAM/matrices/solution/solution.C     |  33 ++-
 .../strainRateFunction/strainRateFunction.C   |   6 +
 src/engine/engineTime/freePiston/freePiston.C |   2 +-
 .../uniformFixedGradientFaPatchField.C        |   3 +-
 .../uniformMixed/uniformMixedFaPatchField.C   |   9 +-
 src/finiteVolume/Make/files                   |   3 +
 .../PatchFunction1/lookup/Lookup.C            |  89 ++++++++
 .../PatchFunction1/lookup/Lookup.H            | 142 +++++++++++++
 .../PatchFunction1/lookup/LookupField.C       | 198 ++++++++++++++++++
 .../PatchFunction1/lookup/LookupField.H       | 171 +++++++++++++++
 .../PatchFunction1/lookup/lookupBase.C        |  52 +++++
 .../PatchFunction1/lookup/lookupBase.H        |  97 +++++++++
 .../PatchFunction1/lookup/makeLookupFields.C  |  40 ++++
 .../PatchFunction1/lookup/makeLookups.C       |  46 ++++
 src/mesh/extrudeModel/radial/radial.C         |   2 -
 .../PatchFunction1/MappedFile/MappedFile.C    |  10 +-
 .../coordinateScaling/coordinateScaling.C     |   2 +-
 ...allBoilingWallFunctionFvPatchScalarField.C |   2 +-
 .../restraints/externalForce/externalForce.C  |   2 +-
 .../enthalpySorptionFvPatchScalarField.C      |   2 +-
 tutorials/mesh/blockMesh/pipe/0.orig/U        |  49 +++++
 tutorials/mesh/blockMesh/pipe/0.orig/epsilon  |  50 +++++
 tutorials/mesh/blockMesh/pipe/0.orig/k        |  50 +++++
 tutorials/mesh/blockMesh/pipe/0.orig/nuTilda  |  49 +++++
 tutorials/mesh/blockMesh/pipe/0.orig/nut      |  51 +++++
 tutorials/mesh/blockMesh/pipe/0.orig/omega    |  50 +++++
 tutorials/mesh/blockMesh/pipe/0.orig/p        |  48 +++++
 tutorials/mesh/blockMesh/pipe/Allclean        |   2 +-
 tutorials/mesh/blockMesh/pipe/Allrun          |   4 +
 .../pipe/constant/transportProperties         |  22 ++
 .../pipe/constant/turbulenceProperties        |  31 +++
 .../mesh/blockMesh/pipe/system/controlDict    |  73 ++++++-
 .../mesh/blockMesh/pipe/system/fvSchemes      |  39 +++-
 .../mesh/blockMesh/pipe/system/fvSolution     |  55 +++++
 34 files changed, 1456 insertions(+), 28 deletions(-)
 create mode 100644 src/finiteVolume/expressions/PatchFunction1/lookup/Lookup.C
 create mode 100644 src/finiteVolume/expressions/PatchFunction1/lookup/Lookup.H
 create mode 100644 src/finiteVolume/expressions/PatchFunction1/lookup/LookupField.C
 create mode 100644 src/finiteVolume/expressions/PatchFunction1/lookup/LookupField.H
 create mode 100644 src/finiteVolume/expressions/PatchFunction1/lookup/lookupBase.C
 create mode 100644 src/finiteVolume/expressions/PatchFunction1/lookup/lookupBase.H
 create mode 100644 src/finiteVolume/expressions/PatchFunction1/lookup/makeLookupFields.C
 create mode 100644 src/finiteVolume/expressions/PatchFunction1/lookup/makeLookups.C
 create mode 100644 tutorials/mesh/blockMesh/pipe/0.orig/U
 create mode 100644 tutorials/mesh/blockMesh/pipe/0.orig/epsilon
 create mode 100644 tutorials/mesh/blockMesh/pipe/0.orig/k
 create mode 100644 tutorials/mesh/blockMesh/pipe/0.orig/nuTilda
 create mode 100644 tutorials/mesh/blockMesh/pipe/0.orig/nut
 create mode 100644 tutorials/mesh/blockMesh/pipe/0.orig/omega
 create mode 100644 tutorials/mesh/blockMesh/pipe/0.orig/p
 create mode 100644 tutorials/mesh/blockMesh/pipe/constant/transportProperties
 create mode 100644 tutorials/mesh/blockMesh/pipe/constant/turbulenceProperties

diff --git a/src/OpenFOAM/matrices/solution/solution.C b/src/OpenFOAM/matrices/solution/solution.C
index e62f9ee6057..401b35393e7 100644
--- a/src/OpenFOAM/matrices/solution/solution.C
+++ b/src/OpenFOAM/matrices/solution/solution.C
@@ -36,6 +36,7 @@ License
 namespace Foam
 {
     defineDebugSwitchWithName(solution, "solution", 0);
+    registerDebugSwitchWithName(solution, solution, "solution");
 }
 
 // List of sub-dictionaries to rewrite
@@ -113,7 +114,7 @@ void Foam::solution::read(const dictionary& dict)
         {
             fieldRelaxDefault_.reset
             (
-                new Function1Types::Constant<scalar>("default", 0)
+                new Function1Types::Constant<scalar>("default", 0, &db())
             );
         }
 
@@ -127,7 +128,7 @@ void Foam::solution::read(const dictionary& dict)
         {
             eqnRelaxDefault_.reset
             (
-                new Function1Types::Constant<scalar>("default", 0)
+                new Function1Types::Constant<scalar>("default", 0, &db())
             );
         }
 
@@ -348,15 +349,29 @@ bool Foam::solution::relaxField(const word& name, scalar& factor) const
             &db()
         )().value(time().timeOutputValue());
 
+        DebugInfo
+            << "Field relaxation factor for " << name
+            << " is " << factor
+            << " (from Function1)" << endl;
+
         return true;
     }
     else if (fieldRelaxDict_.found("default") && fieldRelaxDefault_)
     {
         factor = fieldRelaxDefault_->value(time().timeOutputValue());
+
+        DebugInfo
+            << "Field relaxation factor for " << name
+            << " is " << factor
+            << " (from default " << eqnRelaxDefault_->type() << ')' << endl;
+
         return true;
     }
 
     // Fallthrough - nothing found
+
+    DebugInfo<< "No field relaxation factor for " << name << endl;
+
     return false;
 }
 
@@ -376,15 +391,29 @@ bool Foam::solution::relaxEquation(const word& name, scalar& factor) const
             &db()
         )().value(time().timeOutputValue());
 
+        DebugInfo
+            << "Equation relaxation factor for " << name
+            << " is " << factor
+            << " (from Function1)" << endl;
+
         return true;
     }
     else if (eqnRelaxDict_.found("default") && eqnRelaxDefault_)
     {
         factor = eqnRelaxDefault_->value(time().timeOutputValue());
+
+        DebugInfo
+            << "Equation relaxation factor for " << name
+            << " is " << factor
+            << " (from default " << eqnRelaxDefault_->type() << ')' << endl;
+
         return true;
     }
 
     // Fallthrough - nothing found
+
+    DebugInfo<< "No equation relaxation factor for " << name << endl;
+
     return false;
 }
 
diff --git a/src/TurbulenceModels/turbulenceModels/laminar/generalizedNewtonian/generalizedNewtonianViscosityModels/strainRateFunction/strainRateFunction.C b/src/TurbulenceModels/turbulenceModels/laminar/generalizedNewtonian/generalizedNewtonianViscosityModels/strainRateFunction/strainRateFunction.C
index 04d3792a9d8..a32f5e630b7 100644
--- a/src/TurbulenceModels/turbulenceModels/laminar/generalizedNewtonian/generalizedNewtonianViscosityModels/strainRateFunction/strainRateFunction.C
+++ b/src/TurbulenceModels/turbulenceModels/laminar/generalizedNewtonian/generalizedNewtonianViscosityModels/strainRateFunction/strainRateFunction.C
@@ -110,6 +110,12 @@ nu
         dimensionedScalar(dimViscosity, Zero)
     );
 
+    // Set the database - could not be done in constructor
+    const_cast<Function1<scalar>&>(strainRateFunction_()).resetDb
+    (
+        nu0.mesh().thisDb()
+    );
+
     tnu.ref().primitiveFieldRef() = strainRateFunction_->value(strainRate);
 
     volScalarField::Boundary& nuBf = tnu.ref().boundaryFieldRef();
diff --git a/src/engine/engineTime/freePiston/freePiston.C b/src/engine/engineTime/freePiston/freePiston.C
index ff490dcf110..e29ae83c5b8 100644
--- a/src/engine/engineTime/freePiston/freePiston.C
+++ b/src/engine/engineTime/freePiston/freePiston.C
@@ -59,7 +59,7 @@ Foam::freePiston::freePiston
     ),
     pistonPositionTime_
     (
-        Function1<scalar>::New("pistonPositionTime", dict_)
+        Function1<scalar>::New("pistonPositionTime", dict_, &db())
     )
 {}
 
diff --git a/src/finiteArea/fields/faPatchFields/derived/uniformFixedGradient/uniformFixedGradientFaPatchField.C b/src/finiteArea/fields/faPatchFields/derived/uniformFixedGradient/uniformFixedGradientFaPatchField.C
index 650cdbef1af..2ddafdbe8ce 100644
--- a/src/finiteArea/fields/faPatchFields/derived/uniformFixedGradient/uniformFixedGradientFaPatchField.C
+++ b/src/finiteArea/fields/faPatchFields/derived/uniformFixedGradient/uniformFixedGradientFaPatchField.C
@@ -69,7 +69,8 @@ Foam::uniformFixedGradientFaPatchField<Type>::uniformFixedGradientFaPatchField
         (
             /* p.patch(), */
             "uniformGradient",
-            dict
+            dict,
+            &iF.db()
         )
     )
 {
diff --git a/src/finiteArea/fields/faPatchFields/derived/uniformMixed/uniformMixedFaPatchField.C b/src/finiteArea/fields/faPatchFields/derived/uniformMixed/uniformMixedFaPatchField.C
index 8f29c0f425d..75b332323d2 100644
--- a/src/finiteArea/fields/faPatchFields/derived/uniformMixed/uniformMixedFaPatchField.C
+++ b/src/finiteArea/fields/faPatchFields/derived/uniformMixed/uniformMixedFaPatchField.C
@@ -74,7 +74,8 @@ Foam::uniformMixedFaPatchField<Type>::uniformMixedFaPatchField
         (
             /* p.patch(), */
             "uniformValue",
-            dict
+            dict,
+            &iF.db()
         )
     ),
     refGradFunc_
@@ -83,7 +84,8 @@ Foam::uniformMixedFaPatchField<Type>::uniformMixedFaPatchField
         (
             // p.patch(),
             "uniformGradient",
-            dict
+            dict,
+            &iF.db()
         )
     ),
     valueFractionFunc_(nullptr)
@@ -101,7 +103,8 @@ Foam::uniformMixedFaPatchField<Type>::uniformMixedFaPatchField
                 (
                     /* p.patch(), */
                     "uniformValueFraction",
-                    dict
+                    dict,
+                    &iF.db()
                 )
             );
         }
diff --git a/src/finiteVolume/Make/files b/src/finiteVolume/Make/files
index dca1ff2d9f8..36b69f77f08 100644
--- a/src/finiteVolume/Make/files
+++ b/src/finiteVolume/Make/files
@@ -295,6 +295,9 @@ $(patchExpr)/patchExprLemonParser.lyy-m4
 $(patchExpr)/patchExprScanner.cc
 
 $(expr)/PatchFunction1/makePatchFunction1Expression.C
+$(expr)/PatchFunction1/lookup/lookupBase.C
+$(expr)/PatchFunction1/lookup/makeLookups.C
+$(expr)/PatchFunction1/lookup/makeLookupFields.C
 
 volumeExpr = $(expr)/volume
 $(volumeExpr)/volumeExpr.C
diff --git a/src/finiteVolume/expressions/PatchFunction1/lookup/Lookup.C b/src/finiteVolume/expressions/PatchFunction1/lookup/Lookup.C
new file mode 100644
index 00000000000..9187277e3ee
--- /dev/null
+++ b/src/finiteVolume/expressions/PatchFunction1/lookup/Lookup.C
@@ -0,0 +1,89 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | www.openfoam.com
+     \\/     M anipulation  |
+-------------------------------------------------------------------------------
+    Copyright (C) 2024 OpenCFD Ltd.
+-------------------------------------------------------------------------------
+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 "uniformDimensionedFields.H"
+
+// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
+
+template<class Type>
+Foam::Function1Types::Lookup<Type>::Lookup
+(
+    const word& entryName,
+    const dictionary& dict,
+    const objectRegistry* obrPtr
+)
+:
+    lookupBase(dict),
+    Function1<Type>(entryName, dict, obrPtr)
+{}
+
+
+template<class Type>
+Foam::Function1Types::Lookup<Type>::Lookup(const Lookup<Type>& rhs)
+:
+    lookupBase(rhs),
+    Function1<Type>(rhs)
+{}
+
+
+// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
+
+template<class Type>
+inline Type Foam::Function1Types::Lookup<Type>::value(const scalar t) const
+{
+    const objectRegistry& db = function1Base::obr();
+    const auto& obj =
+        db.lookupObject<UniformDimensionedField<Type>>(name_, true);
+
+    return obj.value();
+}
+
+
+template<class Type>
+inline Type Foam::Function1Types::Lookup<Type>::integrate
+(
+    const scalar t1,
+    const scalar t2
+) const
+{
+    return (t2-t1)*value(0.5*(t1+t2));
+}
+
+
+template<class Type>
+void Foam::Function1Types::Lookup<Type>::writeData(Ostream& os) const
+{
+    Function1<Type>::writeData(os);
+    os.endEntry();
+
+    os.beginBlock(word(this->name() + "Coeffs"));
+    writeEntries(os);
+    os.endBlock();
+}
+
+
+// ************************************************************************* //
diff --git a/src/finiteVolume/expressions/PatchFunction1/lookup/Lookup.H b/src/finiteVolume/expressions/PatchFunction1/lookup/Lookup.H
new file mode 100644
index 00000000000..4a8c64c85a7
--- /dev/null
+++ b/src/finiteVolume/expressions/PatchFunction1/lookup/Lookup.H
@@ -0,0 +1,142 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | www.openfoam.com
+     \\/     M anipulation  |
+-------------------------------------------------------------------------------
+    Copyright (C) 2024 OpenCFD Ltd.
+-------------------------------------------------------------------------------
+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::Function1Types::Lookup
+
+Description
+    Function1 to lookup UniformDimensionedField from an objectregistry.
+
+    The dictionary specification would typically resemble this:
+    \verbatim
+    entry
+    {
+        type        lookup;
+        name        myScalar;
+    }
+    \endverbatim
+
+    where the entries mean:
+    \table
+        Property  | Description             | Type | Reqd | Default
+        type      | Function type: lookup   | word | yes |
+        name      | name of variable        | word | yes |
+    \endtable
+
+SourceFiles
+    Lookup.C
+
+\*---------------------------------------------------------------------------*/
+
+#ifndef Function1Types_Lookup_H
+#define Function1Types_Lookup_H
+
+#include "Function1.H"
+#include "lookupBase.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+namespace Function1Types
+{
+
+/*---------------------------------------------------------------------------*\
+                           Class Lookup Declaration
+\*---------------------------------------------------------------------------*/
+
+template<class Type>
+class Lookup
+:
+    public lookupBase,
+    public Function1<Type>
+{
+public:
+
+    // Runtime type information
+    TypeName("lookup");
+
+
+    // Generated Methods
+
+        //- No copy assignment
+        void operator=(const Lookup<Type>&) = delete;
+
+
+    // Constructors
+
+        //- Construct from entry name, dictionary and optional registry
+        Lookup
+        (
+            const word& entryName,
+            const dictionary& dict,
+            const objectRegistry* obrPtr = nullptr
+        );
+
+        //- Copy construct
+        explicit Lookup(const Lookup<Type>& rhs);
+
+        //- Construct and return a clone
+        virtual tmp<Function1<Type>> clone() const
+        {
+            return tmp<Function1<Type>>(new Lookup<Type>(*this));
+        }
+
+
+    //- Destructor
+    virtual ~Lookup() = default;
+
+
+    // Member Functions
+
+        //- Return value for time t
+        virtual Type value(const scalar t) const;
+
+        //- Integrate between two (scalar) values
+        virtual Type integrate(const scalar x1, const scalar x2) const;
+
+        //- Write in dictionary format
+        virtual void writeData(Ostream& os) const;
+
+        using lookupBase::writeEntries;
+};
+
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+} // End namespace Function1Types
+} // End namespace Foam
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#ifdef NoRepository
+    #include "Lookup.C"
+#endif
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#endif
+
+// ************************************************************************* //
diff --git a/src/finiteVolume/expressions/PatchFunction1/lookup/LookupField.C b/src/finiteVolume/expressions/PatchFunction1/lookup/LookupField.C
new file mode 100644
index 00000000000..2a202b6bfdc
--- /dev/null
+++ b/src/finiteVolume/expressions/PatchFunction1/lookup/LookupField.C
@@ -0,0 +1,198 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | www.openfoam.com
+     \\/     M anipulation  |
+-------------------------------------------------------------------------------
+    Copyright (C) 2024 OpenCFD Ltd.
+-------------------------------------------------------------------------------
+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 "volFields.H"
+#include "volMesh.H"
+#include "surfaceFields.H"
+#include "surfaceMesh.H"
+#include "pointFields.H"
+#include "pointMesh.H"
+#include "UniformDimensionedField.H"
+
+// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
+
+template<class Type>
+Foam::PatchFunction1Types::LookupField<Type>::LookupField
+(
+    const polyPatch& pp,
+    const word& redirectType,
+    const word& entryName,
+    const dictionary& dict,
+    const bool faceValues
+)
+:
+    PatchFunction1<Type>(pp, entryName, dict, faceValues),
+    lookupBase(dict)
+{}
+
+
+template<class Type>
+Foam::PatchFunction1Types::LookupField<Type>::LookupField
+(
+    const LookupField<Type>& rhs,
+    const polyPatch& pp
+)
+:
+    PatchFunction1<Type>(rhs, pp),
+    lookupBase(rhs)
+{}
+
+
+template<class Type>
+Foam::PatchFunction1Types::LookupField<Type>::LookupField
+(
+    const LookupField<Type>& rhs
+)
+:
+    LookupField<Type>(rhs, rhs.patch())
+{}
+
+
+// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
+
+template<class Type>
+Foam::tmp<Foam::Field<Type>>
+Foam::PatchFunction1Types::LookupField<Type>::value(const scalar x) const
+{
+    typedef UniformDimensionedField<Type> UType;
+
+    const objectRegistry& db = patchFunction1Base::obr();
+    const label patchi = patchFunction1Base::patch().index();
+
+    if (this->faceValues())
+    {
+        typedef GeometricField<Type, fvPatchField, volMesh> VType;
+        typedef GeometricField<Type, fvsPatchField, surfaceMesh> SType;
+
+        // Try:
+        // - as volField in local scope
+        // - as surfaceField in local scope
+        // - as UniformDimensionedField recursively
+
+        const regIOobject* ptr = db.cfindIOobject(name_, false);
+
+        if (ptr)
+        {
+            const auto* vPtr = dynamic_cast<const VType*>(ptr);
+            if (vPtr)
+            {
+                return vPtr->boundaryField()[patchi];
+            }
+
+            const auto* sPtr = dynamic_cast<const SType*>(ptr);
+            if (sPtr)
+            {
+                return sPtr->boundaryField()[patchi];
+            }
+
+            const auto* uPtr = dynamic_cast<const UType*>(ptr);
+            if (uPtr)
+            {
+                return Field<Type>(this->size(), uPtr->value());
+            }
+        }
+
+        // Done db level. Try recursion
+        ptr = db.parent().cfindIOobject(name_, true);
+
+        if (ptr)
+        {
+            const auto* uPtr = dynamic_cast<const UType*>(ptr);
+            if (uPtr)
+            {
+                return Field<Type>(this->size(), uPtr->value());
+            }
+        }
+
+        FatalErrorInFunction
+            << nl
+            << "    failed lookup of " << name_
+            << " (objectRegistry "
+            << db.name()
+            << ")\n    available objects of type " << VType::typeName
+            << ':' << nl
+            << db.names<VType>() << nl
+            << "    available objects of type " << SType::typeName
+            << ':' << nl
+            << db.names<SType>() << nl
+            << "    available objects of type " << UType::typeName
+            << ':' << nl
+            << db.names<UType>() << nl
+            << exit(FatalError);
+        return Field<Type>::null();
+    }
+    else
+    {
+        // Assume pointField
+        typedef GeometricField<Type, pointPatchField, pointMesh> PType;
+
+        const regIOobject* ptr = db.cfindIOobject(name_, false);
+
+        if (ptr)
+        {
+            const auto* pPtr = dynamic_cast<const PType*>(ptr);
+            if (pPtr)
+            {
+                return pPtr->boundaryField()[patchi].patchInternalField();
+            }
+        }
+
+        // Re-do as uniform field. Note: could repeat logic above
+        const auto& obj = db.lookupObject<UType>(name_, true);
+
+        return Field<Type>(this->size(), obj.value());
+    }
+}
+
+
+template<class Type>
+Foam::tmp<Foam::Field<Type>>
+Foam::PatchFunction1Types::LookupField<Type>::integrate
+(
+    const scalar x1,
+    const scalar x2
+) const
+{
+    return (x2-x1)*value(0.5*(x1+x2));
+}
+
+
+template<class Type>
+void Foam::PatchFunction1Types::LookupField<Type>::writeData
+(
+    Ostream& os
+) const
+{
+    PatchFunction1<Type>::writeData(os);
+    os.writeEntry(this->name(), type());
+    os.beginBlock(word(this->name() + "Coeffs"));
+    writeEntries(os);
+    os.endBlock();
+}
+
+
+// ************************************************************************* //
diff --git a/src/finiteVolume/expressions/PatchFunction1/lookup/LookupField.H b/src/finiteVolume/expressions/PatchFunction1/lookup/LookupField.H
new file mode 100644
index 00000000000..7fb5d25a6f0
--- /dev/null
+++ b/src/finiteVolume/expressions/PatchFunction1/lookup/LookupField.H
@@ -0,0 +1,171 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | www.openfoam.com
+     \\/     M anipulation  |
+-------------------------------------------------------------------------------
+    Copyright (C) 2024 OpenCFD Ltd.
+-------------------------------------------------------------------------------
+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::PatchFunction1Types::LookupField
+
+Description
+    PatchFunction1 to lookup volField/surfaceField or pointField from an
+    objectregistry and return its value on the patch.
+
+    The dictionary specification would typically resemble this:
+    \verbatim
+    entry
+    {
+        type        lookup;
+        name        myField;
+    }
+    \endverbatim
+
+    where the entries mean:
+    \table
+        Property  | Description                 | Type | Reqd | Default
+        type      | Function type: lookup       | word | yes |
+        name      | name of volField/pointField | word | yes |
+    \endtable
+
+SourceFiles
+    LookupField.C
+
+\*---------------------------------------------------------------------------*/
+
+#ifndef PatchFunction1Types_LookupField_H
+#define PatchFunction1Types_LookupField_H
+
+#include "PatchFunction1.H"
+#include "lookupBase.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+namespace PatchFunction1Types
+{
+
+/*---------------------------------------------------------------------------*\
+                            Class LookupField Declaration
+\*---------------------------------------------------------------------------*/
+
+template<class Type>
+class LookupField
+:
+    public PatchFunction1<Type>,
+    protected lookupBase
+{
+
+
+public:
+
+    // Runtime type information
+    TypeName("lookup");
+
+
+    // Generated Methods
+
+        //- No copy assignment
+        void operator=(const LookupField<Type>&) = delete;
+
+
+    // Constructors
+
+        //- Construct from entry name, dictionary and optional registry
+        LookupField
+        (
+            const polyPatch& pp,
+            const word& redirectType,
+            const word& entryName,
+            const dictionary& dict,
+            const bool faceValues = true
+        );
+
+        //- Copy construct, setting patch
+        explicit LookupField
+        (
+            const LookupField<Type>& rhs,
+            const polyPatch& pp
+        );
+
+        //- Copy construct
+        explicit LookupField(const LookupField<Type>& rhs);
+
+        //- Construct and return a clone
+        virtual tmp<PatchFunction1<Type>> clone() const
+        {
+            return tmp<PatchFunction1<Type>>(new LookupField<Type>(*this));
+        }
+
+        //- Return a clone, setting the patch
+        virtual tmp<PatchFunction1<Type>> clone(const polyPatch& pp) const
+        {
+            return PatchFunction1<Type>::Clone(*this, pp);
+        }
+
+
+    //- Destructor
+    virtual ~LookupField() = default;
+
+
+    // Member Functions
+
+        //- Is value uniform (i.e. independent of coordinate)
+        virtual inline bool uniform() const { return false; }
+
+
+        // Evaluation
+
+            //- Return Lookup value
+            virtual tmp<Field<Type>> value(const scalar x) const;
+
+            //- Integrate between two values
+            virtual tmp<Field<Type>> integrate
+            (
+                const scalar x1,
+                const scalar x2
+            ) const;
+
+
+        // I-O
+
+            //- Write in dictionary format
+            virtual void writeData(Ostream& os) const;
+};
+
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+} // End namespace PatchFunction1Types
+} // End namespace Foam
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#ifdef NoRepository
+    #include "LookupField.C"
+#endif
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#endif
+
+// ************************************************************************* //
diff --git a/src/finiteVolume/expressions/PatchFunction1/lookup/lookupBase.C b/src/finiteVolume/expressions/PatchFunction1/lookup/lookupBase.C
new file mode 100644
index 00000000000..4c3a10cee18
--- /dev/null
+++ b/src/finiteVolume/expressions/PatchFunction1/lookup/lookupBase.C
@@ -0,0 +1,52 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | www.openfoam.com
+     \\/     M anipulation  |
+-------------------------------------------------------------------------------
+    Copyright (C) 2024 OpenCFD Ltd.
+-------------------------------------------------------------------------------
+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 "lookupBase.H"
+
+// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
+
+Foam::lookupBase::lookupBase(const dictionary& dict)
+:
+    name_(dict.get<word>("name"))
+{}
+
+
+Foam::lookupBase::lookupBase(const lookupBase& rhs)
+:
+    name_(rhs.name_)
+{}
+
+
+// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
+
+void Foam::lookupBase::writeEntries(Ostream& os) const
+{
+    os.writeEntry<word>("name", name_);
+}
+
+
+// ************************************************************************* //
diff --git a/src/finiteVolume/expressions/PatchFunction1/lookup/lookupBase.H b/src/finiteVolume/expressions/PatchFunction1/lookup/lookupBase.H
new file mode 100644
index 00000000000..271466ff4be
--- /dev/null
+++ b/src/finiteVolume/expressions/PatchFunction1/lookup/lookupBase.H
@@ -0,0 +1,97 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | www.openfoam.com
+     \\/     M anipulation  |
+-------------------------------------------------------------------------------
+    Copyright (C) 2024 OpenCFD Ltd.
+-------------------------------------------------------------------------------
+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::lookupBase
+
+Description
+    Base level of Lookup and LookupField classes.
+
+SourceFiles
+    lookupBase.C
+
+\*---------------------------------------------------------------------------*/
+
+#ifndef lookupBase_H
+#define lookupBase_H
+
+#include "dictionary.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+
+/*---------------------------------------------------------------------------*\
+                         Class lookupBase Declaration
+\*---------------------------------------------------------------------------*/
+
+class lookupBase
+{
+protected:
+
+    // Protected Data
+
+        //- Key of object
+        word name_;
+
+
+public:
+
+    // Generated Methods
+
+        //- No copy assignment
+        void operator=(const lookupBase&) = delete;
+
+
+    // Constructors
+
+        //- Construct from dictionary
+        explicit lookupBase(const dictionary& dict);
+
+        //- Copy construct
+        explicit lookupBase(const lookupBase& rhs);
+
+
+    //- Destructor
+    virtual ~lookupBase() = default;
+
+
+    // Member Functions
+
+        //- Write coefficient entries in dictionary format
+        virtual void writeEntries(Ostream& os) const;
+};
+
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+} // End namespace Foam
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#endif
+
+// ************************************************************************* //
diff --git a/src/finiteVolume/expressions/PatchFunction1/lookup/makeLookupFields.C b/src/finiteVolume/expressions/PatchFunction1/lookup/makeLookupFields.C
new file mode 100644
index 00000000000..5989e59a8c8
--- /dev/null
+++ b/src/finiteVolume/expressions/PatchFunction1/lookup/makeLookupFields.C
@@ -0,0 +1,40 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | www.openfoam.com
+     \\/     M anipulation  |
+-------------------------------------------------------------------------------
+    Copyright (C) 2024 OpenCFD Ltd.
+-------------------------------------------------------------------------------
+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 "LookupField.H"
+#include "addToRunTimeSelectionTable.H"
+
+namespace Foam
+{
+    makePatchFunction1Type(LookupField, scalar);
+    makePatchFunction1Type(LookupField, vector);
+    makePatchFunction1Type(LookupField, sphericalTensor);
+    makePatchFunction1Type(LookupField, symmTensor);
+    makePatchFunction1Type(LookupField, tensor);
+}
+
+// ************************************************************************* //
diff --git a/src/finiteVolume/expressions/PatchFunction1/lookup/makeLookups.C b/src/finiteVolume/expressions/PatchFunction1/lookup/makeLookups.C
new file mode 100644
index 00000000000..e98d88563d8
--- /dev/null
+++ b/src/finiteVolume/expressions/PatchFunction1/lookup/makeLookups.C
@@ -0,0 +1,46 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | www.openfoam.com
+     \\/     M anipulation  |
+-------------------------------------------------------------------------------
+    Copyright (C) 2024 OpenCFD Ltd.
+-------------------------------------------------------------------------------
+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 "Lookup.H"
+#include "fieldTypes.H"
+#include "PatchFunction1.H"
+#include "addToRunTimeSelectionTable.H"
+#include "UniformValueField.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+    makeFunction1Type(Lookup, scalar);
+    makeFunction1Type(Lookup, vector);
+    makeFunction1Type(Lookup, sphericalTensor);
+    makeFunction1Type(Lookup, symmTensor);
+    makeFunction1Type(Lookup, tensor);
+}
+
+
+// ************************************************************************* //
diff --git a/src/mesh/extrudeModel/radial/radial.C b/src/mesh/extrudeModel/radial/radial.C
index 812cf1f613b..6e73f32390d 100644
--- a/src/mesh/extrudeModel/radial/radial.C
+++ b/src/mesh/extrudeModel/radial/radial.C
@@ -66,8 +66,6 @@ point radial::operator()
 
     scalar r = R_->value(layer);
 
-Pout<< "** for layer " << layer << " r:" << r << endl;
-
     return r*rsHat;
 }
 
diff --git a/src/meshTools/PatchFunction1/MappedFile/MappedFile.C b/src/meshTools/PatchFunction1/MappedFile/MappedFile.C
index 98f54c33342..16d7fbc850b 100644
--- a/src/meshTools/PatchFunction1/MappedFile/MappedFile.C
+++ b/src/meshTools/PatchFunction1/MappedFile/MappedFile.C
@@ -60,7 +60,15 @@ Foam::PatchFunction1Types::MappedFile<Type>::MappedFile
     sampleIndex_(-1, -1),
     sampleAverage_(Zero, Zero),
     sampleValues_(),
-    offset_(Function1<Type>::NewIfPresent("offset", dict))
+    offset_
+    (
+        Function1<Type>::NewIfPresent
+        (
+            "offset",
+            dict,
+            patchFunction1Base::whichDb()
+        )
+    )
 {
     if (fieldTableName_.empty())
     {
diff --git a/src/meshTools/PatchFunction1/coordinateScaling/coordinateScaling.C b/src/meshTools/PatchFunction1/coordinateScaling/coordinateScaling.C
index 6badd8d9e4e..e24f7ac804f 100644
--- a/src/meshTools/PatchFunction1/coordinateScaling/coordinateScaling.C
+++ b/src/meshTools/PatchFunction1/coordinateScaling/coordinateScaling.C
@@ -51,7 +51,7 @@ Foam::coordinateScaling<Type>::coordinateScaling
     {
         const word key("scale" + Foam::name(dir+1));
 
-        auto scaling = Function1<Type>::NewIfPresent(key, dict);
+        auto scaling = Function1<Type>::NewIfPresent(key, dict, &obr);
 
         if (scaling)
         {
diff --git a/src/phaseSystemModels/reactingEuler/multiphaseSystem/derivedFvPatchFields/alphatWallBoilingWallFunction/alphatWallBoilingWallFunctionFvPatchScalarField.C b/src/phaseSystemModels/reactingEuler/multiphaseSystem/derivedFvPatchFields/alphatWallBoilingWallFunction/alphatWallBoilingWallFunctionFvPatchScalarField.C
index c405a54f1a9..7adb2c0d55a 100644
--- a/src/phaseSystemModels/reactingEuler/multiphaseSystem/derivedFvPatchFields/alphatWallBoilingWallFunction/alphatWallBoilingWallFunctionFvPatchScalarField.C
+++ b/src/phaseSystemModels/reactingEuler/multiphaseSystem/derivedFvPatchFields/alphatWallBoilingWallFunction/alphatWallBoilingWallFunctionFvPatchScalarField.C
@@ -116,7 +116,7 @@ alphatWallBoilingWallFunctionFvPatchScalarField
     alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField(p, iF, dict),
     otherPhaseName_(dict.get<word>("otherPhase")),
     phaseType_(phaseTypeNames_.get("phaseType", dict)),
-    relax_(Function1<scalar>::New("relax", dict)),
+    relax_(Function1<scalar>::New("relax", dict, &db())),
     AbyV_(p.size(), 0),
     alphatConv_(p.size(), 0),
     dDep_(p.size(), 1e-5),
diff --git a/src/rigidBodyDynamics/restraints/externalForce/externalForce.C b/src/rigidBodyDynamics/restraints/externalForce/externalForce.C
index 26820cfc54b..ac7c57079e8 100644
--- a/src/rigidBodyDynamics/restraints/externalForce/externalForce.C
+++ b/src/rigidBodyDynamics/restraints/externalForce/externalForce.C
@@ -109,7 +109,7 @@ bool Foam::RBD::restraints::externalForce::read
 
     coeffs_.readEntry("location", location_);
 
-    externalForce_ = Function1<vector>::New("force", coeffs_);
+    externalForce_ = Function1<vector>::New("force", coeffs_, &model_.time());
 
     return true;
 }
diff --git a/src/thermophysicalModels/reactionThermo/derivedFvPatchFields/enthalpySorption/enthalpySorptionFvPatchScalarField.C b/src/thermophysicalModels/reactionThermo/derivedFvPatchFields/enthalpySorption/enthalpySorptionFvPatchScalarField.C
index a3caadc95aa..ca9bee458fa 100644
--- a/src/thermophysicalModels/reactionThermo/derivedFvPatchFields/enthalpySorption/enthalpySorptionFvPatchScalarField.C
+++ b/src/thermophysicalModels/reactionThermo/derivedFvPatchFields/enthalpySorption/enthalpySorptionFvPatchScalarField.C
@@ -89,7 +89,7 @@ Foam::enthalpySorptionFvPatchScalarField::enthalpySorptionFvPatchScalarField
         case enthalpyModelType::calculated:
         {
             enthalpyMassLoadPtr_ =
-                Function1<scalar>::New("enthalpyTable", dict);
+                Function1<scalar>::New("enthalpyTable", dict, &iF.db());
             break;
         }
         case enthalpyModelType::estimated:
diff --git a/tutorials/mesh/blockMesh/pipe/0.orig/U b/tutorials/mesh/blockMesh/pipe/0.orig/U
new file mode 100644
index 00000000000..62beac9a3d2
--- /dev/null
+++ b/tutorials/mesh/blockMesh/pipe/0.orig/U
@@ -0,0 +1,49 @@
+/*--------------------------------*- C++ -*----------------------------------*\
+| =========                 |                                                 |
+| \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox           |
+|  \\    /   O peration     | Version:  v2312                                 |
+|   \\  /    A nd           | Website:  www.openfoam.com                      |
+|    \\/     M anipulation  |                                                 |
+\*---------------------------------------------------------------------------*/
+FoamFile
+{
+    version     2.0;
+    format      ascii;
+    class       volVectorField;
+    object      U;
+}
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+dimensions      [0 1 -1 0 0 0 0];
+
+internalField   uniform (0 0 0);
+
+boundaryField
+{
+    #includeEtc "caseDicts/setConstraintTypes"
+
+    side
+    {
+        type            surfaceNormalFixedValue;
+        refValue        uniform -10;
+    }
+
+    inlet
+    {
+        type            surfaceNormalFixedValue;
+        refValue        uniform -10;
+    }
+
+    outlet
+    {
+        type            zeroGradient;
+    }
+
+    walls
+    {
+        type            noSlip;
+    }
+}
+
+
+// ************************************************************************* //
diff --git a/tutorials/mesh/blockMesh/pipe/0.orig/epsilon b/tutorials/mesh/blockMesh/pipe/0.orig/epsilon
new file mode 100644
index 00000000000..48a0a6dcf3d
--- /dev/null
+++ b/tutorials/mesh/blockMesh/pipe/0.orig/epsilon
@@ -0,0 +1,50 @@
+/*--------------------------------*- C++ -*----------------------------------*\
+| =========                 |                                                 |
+| \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox           |
+|  \\    /   O peration     | Version:  v2312                                 |
+|   \\  /    A nd           | Website:  www.openfoam.com                      |
+|    \\/     M anipulation  |                                                 |
+\*---------------------------------------------------------------------------*/
+FoamFile
+{
+    version     2.0;
+    format      ascii;
+    class       volScalarField;
+    object      epsilon;
+}
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+dimensions      [0 2 -3 0 0 0 0];
+
+internalField   uniform 14.855;
+
+boundaryField
+{
+    #includeEtc "caseDicts/setConstraintTypes"
+
+    side
+    {
+        type            fixedValue;
+        value           $internalField;
+    }
+
+    inlet
+    {
+        type            fixedValue;
+        value           $internalField;
+    }
+
+    outlet
+    {
+        type            zeroGradient;
+    }
+
+    walls
+    {
+        type            epsilonWallFunction;
+        value           $internalField;
+    }
+}
+
+
+// ************************************************************************* //
diff --git a/tutorials/mesh/blockMesh/pipe/0.orig/k b/tutorials/mesh/blockMesh/pipe/0.orig/k
new file mode 100644
index 00000000000..9c3ce53c15a
--- /dev/null
+++ b/tutorials/mesh/blockMesh/pipe/0.orig/k
@@ -0,0 +1,50 @@
+/*--------------------------------*- C++ -*----------------------------------*\
+| =========                 |                                                 |
+| \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox           |
+|  \\    /   O peration     | Version:  v2312                                 |
+|   \\  /    A nd           | Website:  www.openfoam.com                      |
+|    \\/     M anipulation  |                                                 |
+\*---------------------------------------------------------------------------*/
+FoamFile
+{
+    version     2.0;
+    format      ascii;
+    class       volScalarField;
+    object      k;
+}
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+dimensions      [0 2 -2 0 0 0 0];
+
+internalField   uniform 0.375;
+
+boundaryField
+{
+    #includeEtc "caseDicts/setConstraintTypes"
+
+    side
+    {
+        type            fixedValue;
+        value           $internalField;
+    }
+
+    inlet
+    {
+        type            fixedValue;
+        value           $internalField;
+    }
+
+    outlet
+    {
+        type            zeroGradient;
+    }
+
+    walls
+    {
+        type            kqRWallFunction;
+        value           $internalField;
+    }
+}
+
+
+// ************************************************************************* //
diff --git a/tutorials/mesh/blockMesh/pipe/0.orig/nuTilda b/tutorials/mesh/blockMesh/pipe/0.orig/nuTilda
new file mode 100644
index 00000000000..a7cbaa3c858
--- /dev/null
+++ b/tutorials/mesh/blockMesh/pipe/0.orig/nuTilda
@@ -0,0 +1,49 @@
+/*--------------------------------*- C++ -*----------------------------------*\
+| =========                 |                                                 |
+| \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox           |
+|  \\    /   O peration     | Version:  v2312                                 |
+|   \\  /    A nd           | Website:  www.openfoam.com                      |
+|    \\/     M anipulation  |                                                 |
+\*---------------------------------------------------------------------------*/
+FoamFile
+{
+    version     2.0;
+    format      ascii;
+    class       volScalarField;
+    object      nuTilda;
+}
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+dimensions      [0 2 -1 0 0 0 0];
+
+internalField   uniform 0;
+
+boundaryField
+{
+    #includeEtc "caseDicts/setConstraintTypes"
+
+    side
+    {
+        type            fixedValue;
+        value           uniform 0;
+    }
+
+    inlet
+    {
+        type            fixedValue;
+        value           uniform 0;
+    }
+
+    outlet
+    {
+        type            zeroGradient;
+    }
+
+    walls
+    {
+        type            zeroGradient;
+    }
+}
+
+
+// ************************************************************************* //
diff --git a/tutorials/mesh/blockMesh/pipe/0.orig/nut b/tutorials/mesh/blockMesh/pipe/0.orig/nut
new file mode 100644
index 00000000000..6d81e6a5d2e
--- /dev/null
+++ b/tutorials/mesh/blockMesh/pipe/0.orig/nut
@@ -0,0 +1,51 @@
+/*--------------------------------*- C++ -*----------------------------------*\
+| =========                 |                                                 |
+| \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox           |
+|  \\    /   O peration     | Version:  v2312                                 |
+|   \\  /    A nd           | Website:  www.openfoam.com                      |
+|    \\/     M anipulation  |                                                 |
+\*---------------------------------------------------------------------------*/
+FoamFile
+{
+    version     2.0;
+    format      ascii;
+    class       volScalarField;
+    object      nut;
+}
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+dimensions      [0 2 -1 0 0 0 0];
+
+internalField   uniform 0;
+
+boundaryField
+{
+    #includeEtc "caseDicts/setConstraintTypes"
+
+    side
+    {
+        type            calculated;
+        value           uniform 0;
+    }
+
+    inlet
+    {
+        type            calculated;
+        value           uniform 0;
+    }
+
+    outlet
+    {
+        type            calculated;
+        value           uniform 0;
+    }
+
+    walls
+    {
+        type            nutkWallFunction;
+        value           uniform 0;
+    }
+}
+
+
+// ************************************************************************* //
diff --git a/tutorials/mesh/blockMesh/pipe/0.orig/omega b/tutorials/mesh/blockMesh/pipe/0.orig/omega
new file mode 100644
index 00000000000..dd8bde00665
--- /dev/null
+++ b/tutorials/mesh/blockMesh/pipe/0.orig/omega
@@ -0,0 +1,50 @@
+/*--------------------------------*- C++ -*----------------------------------*\
+| =========                 |                                                 |
+| \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox           |
+|  \\    /   O peration     | Version:  v2312                                 |
+|   \\  /    A nd           | Website:  www.openfoam.com                      |
+|    \\/     M anipulation  |                                                 |
+\*---------------------------------------------------------------------------*/
+FoamFile
+{
+    version     2.0;
+    format      ascii;
+    class       volScalarField;
+    object      omega;
+}
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+dimensions      [0 0 -1 0 0 0 0];
+
+internalField   uniform 440.15;
+
+boundaryField
+{
+    #includeEtc "caseDicts/setConstraintTypes"
+
+    side
+    {
+        type            fixedValue;
+        value           $internalField;
+    }
+
+    inlet
+    {
+        type            fixedValue;
+        value           $internalField;
+    }
+
+    outlet
+    {
+        type            zeroGradient;
+    }
+
+    walls
+    {
+        type            omegaWallFunction;
+        value           $internalField;
+    }
+}
+
+
+// ************************************************************************* //
diff --git a/tutorials/mesh/blockMesh/pipe/0.orig/p b/tutorials/mesh/blockMesh/pipe/0.orig/p
new file mode 100644
index 00000000000..5deb51a1df2
--- /dev/null
+++ b/tutorials/mesh/blockMesh/pipe/0.orig/p
@@ -0,0 +1,48 @@
+/*--------------------------------*- C++ -*----------------------------------*\
+| =========                 |                                                 |
+| \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox           |
+|  \\    /   O peration     | Version:  v2312                                 |
+|   \\  /    A nd           | Website:  www.openfoam.com                      |
+|    \\/     M anipulation  |                                                 |
+\*---------------------------------------------------------------------------*/
+FoamFile
+{
+    version     2.0;
+    format      ascii;
+    class       volScalarField;
+    object      p;
+}
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+dimensions      [0 2 -2 0 0 0 0];
+
+internalField   uniform 0;
+
+boundaryField
+{
+    #includeEtc "caseDicts/setConstraintTypes"
+
+    side
+    {
+        type            zeroGradient;
+    }
+
+    inlet
+    {
+        type            zeroGradient;
+    }
+
+    outlet
+    {
+        type            fixedValue;
+        value           uniform 0;
+    }
+
+    walls
+    {
+        type            zeroGradient;
+    }
+}
+
+
+// ************************************************************************* //
diff --git a/tutorials/mesh/blockMesh/pipe/Allclean b/tutorials/mesh/blockMesh/pipe/Allclean
index d4f5975bc81..64435d02b14 100755
--- a/tutorials/mesh/blockMesh/pipe/Allclean
+++ b/tutorials/mesh/blockMesh/pipe/Allclean
@@ -5,6 +5,6 @@ cd "${0%/*}" || exit                                # Run from this directory
 
 cleanCase0
 
-rm -rf constant
+rm -rf constant/geometry
 
 #------------------------------------------------------------------------------
diff --git a/tutorials/mesh/blockMesh/pipe/Allrun b/tutorials/mesh/blockMesh/pipe/Allrun
index 9a4df7071b1..9c48a0769ce 100755
--- a/tutorials/mesh/blockMesh/pipe/Allrun
+++ b/tutorials/mesh/blockMesh/pipe/Allrun
@@ -11,4 +11,8 @@ cp -f \
 
 runApplication blockMesh
 
+restore0Dir
+
+runApplication $(getApplication)
+
 #------------------------------------------------------------------------------
diff --git a/tutorials/mesh/blockMesh/pipe/constant/transportProperties b/tutorials/mesh/blockMesh/pipe/constant/transportProperties
new file mode 100644
index 00000000000..4908cd4b363
--- /dev/null
+++ b/tutorials/mesh/blockMesh/pipe/constant/transportProperties
@@ -0,0 +1,22 @@
+/*--------------------------------*- C++ -*----------------------------------*\
+| =========                 |                                                 |
+| \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox           |
+|  \\    /   O peration     | Version:  v2312                                 |
+|   \\  /    A nd           | Website:  www.openfoam.com                      |
+|    \\/     M anipulation  |                                                 |
+\*---------------------------------------------------------------------------*/
+FoamFile
+{
+    version     2.0;
+    format      ascii;
+    class       dictionary;
+    object      transportProperties;
+}
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+transportModel  Newtonian;
+
+nu              1e-05;
+
+
+// ************************************************************************* //
diff --git a/tutorials/mesh/blockMesh/pipe/constant/turbulenceProperties b/tutorials/mesh/blockMesh/pipe/constant/turbulenceProperties
new file mode 100644
index 00000000000..00ce4ad5f38
--- /dev/null
+++ b/tutorials/mesh/blockMesh/pipe/constant/turbulenceProperties
@@ -0,0 +1,31 @@
+/*--------------------------------*- C++ -*----------------------------------*\
+| =========                 |                                                 |
+| \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox           |
+|  \\    /   O peration     | Version:  v2312                                 |
+|   \\  /    A nd           | Website:  www.openfoam.com                      |
+|    \\/     M anipulation  |                                                 |
+\*---------------------------------------------------------------------------*/
+FoamFile
+{
+    version     2.0;
+    format      ascii;
+    class       dictionary;
+    object      turbulenceProperties;
+}
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+simulationType      RAS;
+
+RAS
+{
+    // Tested with kEpsilon, realizableKE, kOmega, kOmegaSST,
+    // ShihQuadraticKE, LienCubicKE.
+    RASModel        kEpsilon;
+
+    turbulence      on;
+
+    printCoeffs     on;
+}
+
+
+// ************************************************************************* //
diff --git a/tutorials/mesh/blockMesh/pipe/system/controlDict b/tutorials/mesh/blockMesh/pipe/system/controlDict
index 170175045a4..328116a5f1e 100644
--- a/tutorials/mesh/blockMesh/pipe/system/controlDict
+++ b/tutorials/mesh/blockMesh/pipe/system/controlDict
@@ -18,12 +18,13 @@ libs            (blockMesh);
 
 DebugSwitches
 {
-//    project 1;
-//    searchableExtrudedCircle 1;
-//    projectCurve 1;
+    //project 1;
+    //searchableExtrudedCircle 1;
+    //projectCurve 1;
+    solution    1;
 }
 
-application     blockMesh;
+application     simpleFoam;
 
 startFrom       startTime;
 
@@ -31,13 +32,13 @@ startTime       0;
 
 stopAt          endTime;
 
-endTime         0;
+endTime         100;
 
-deltaT          0;
+deltaT          1;
 
 writeControl    timeStep;
 
-writeInterval   1;
+writeInterval   10;
 
 purgeWrite      0;
 
@@ -53,5 +54,63 @@ timePrecision   6;
 
 runTimeModifiable true;
 
+functions
+{
+    relaxationFactor
+    {
+        // Convoluted example - control relaxation factor
+        libs            (utilityFunctionObjects);
+        type            coded;
+
+        name            relaxationFactor;
+
+        codeRead
+        #{
+            const IOobject io
+            (
+                "banana",
+                mesh().time().constant(),
+                mesh()
+            );
+
+            auto* ptr = const_cast<objectRegistry&>(io.db()).
+                    findObject<uniformDimensionedScalarField>(io.name());
+            if (!ptr)
+            {
+                Info<< "Registering relaxation factor " << io.name() << endl;
+
+                ptr = new uniformDimensionedScalarField
+                (
+                    io,
+                    dimless,
+                    0.80
+                );
+                ptr->store();
+            }
+        #};
+
+        codeExecute
+        #{
+            const IOobject io
+            (
+                "banana",
+                mesh().time().constant(),
+                mesh()
+            );
+
+            auto& val =
+                io.db().lookupObjectRef<uniformDimensionedScalarField>
+                (
+                    io.name()
+                );
+
+            // Ramp a bit
+            val.value() = min(0.99, val.value()+0.01);
+
+            Info<< "Set relaxation factor to " << val << endl;
+        #};
+    }
+}
+
 
 // ************************************************************************* //
diff --git a/tutorials/mesh/blockMesh/pipe/system/fvSchemes b/tutorials/mesh/blockMesh/pipe/system/fvSchemes
index e81493017ea..34b7d87bd9c 100644
--- a/tutorials/mesh/blockMesh/pipe/system/fvSchemes
+++ b/tutorials/mesh/blockMesh/pipe/system/fvSchemes
@@ -15,22 +15,49 @@ FoamFile
 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
 
 ddtSchemes
-{}
+{
+    default         steadyState;
+}
 
 gradSchemes
-{}
+{
+    default         Gauss linear;
+}
 
 divSchemes
-{}
+{
+    default         none;
+
+    div(phi,U)      bounded Gauss linearUpwind grad(U);
+
+    turbulence      bounded Gauss limitedLinear 1;
+    div(phi,k)      $turbulence;
+    div(phi,epsilon) $turbulence;
+    div(phi,omega)  $turbulence;
+
+    div(nonlinearStress) Gauss linear;
+    div((nuEff*dev2(T(grad(U))))) Gauss linear;
+}
 
 laplacianSchemes
-{}
+{
+    default         Gauss linear corrected;
+}
 
 interpolationSchemes
-{}
+{
+    default         linear;
+}
 
 snGradSchemes
-{}
+{
+    default         corrected;
+}
+
+wallDist
+{
+    method          meshWave;
+}
 
 
 // ************************************************************************* //
diff --git a/tutorials/mesh/blockMesh/pipe/system/fvSolution b/tutorials/mesh/blockMesh/pipe/system/fvSolution
index 2ede5a9ce0b..bcc0bcc3ccd 100644
--- a/tutorials/mesh/blockMesh/pipe/system/fvSolution
+++ b/tutorials/mesh/blockMesh/pipe/system/fvSolution
@@ -14,5 +14,60 @@ FoamFile
 }
 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
 
+solvers
+{
+    p
+    {
+        solver          GAMG;
+        tolerance       1e-06;
+        relTol          0.1;
+        smoother        GaussSeidel;
+    }
+
+    "(U|k|epsilon|omega|f|v2)"
+    {
+        solver          smoothSolver;
+        smoother        symGaussSeidel;
+        tolerance       1e-05;
+        relTol          0.1;
+    }
+}
+
+SIMPLE
+{
+    nNonOrthogonalCorrectors 0;
+    consistent      yes;
+
+    residualControl
+    {
+        p               1e-2;
+        U               1e-3;
+        "(k|epsilon|omega|f|v2)" 1e-3;
+    }
+}
+
+relaxationFactors
+{
+    equations
+    {
+        //- Constant relaxation factor
+        //U               0.9; // 0.9 is more stable but 0.95 more convergent
+
+        //- Control relaxation factor through time-based table
+        //U               table
+        //(
+        //    (0  0.7)
+        //    (100 0.95)
+        //);
+
+        //- Control relaxation factor through registered variable
+        //- (type uniformDimensionedScalarField) controlled by functionObject
+        U               lookup;
+        name            banana;
+
+        ".*"            0.95; // 0.9 is more stable but 0.95 more convergent
+    }
+}
+
 
 // ************************************************************************* //
-- 
GitLab