diff --git a/src/finiteArea/Make/files b/src/finiteArea/Make/files
index 039386be84fb4c13c1b97c5d86dca23e03a24925..803c8df34f99713d3e2dfbf378aa9e261228da83 100644
--- a/src/finiteArea/Make/files
+++ b/src/finiteArea/Make/files
@@ -68,7 +68,9 @@ $(derivedFaPatchFields)/outletInlet/outletInletFaPatchFields.C
 $(derivedFaPatchFields)/slip/slipFaPatchFields.C
 $(derivedFaPatchFields)/edgeNormalFixedValue/edgeNormalFixedValueFaPatchVectorField.C
 $(derivedFaPatchFields)/timeVaryingUniformFixedValue/timeVaryingUniformFixedValueFaPatchFields.C
+$(derivedFaPatchFields)/uniformFixedGradient/uniformFixedGradientFaPatchFields.C
 $(derivedFaPatchFields)/uniformFixedValue/uniformFixedValueFaPatchFields.C
+$(derivedFaPatchFields)/uniformMixed/uniformMixedFaPatchFields.C
 $(derivedFaPatchFields)/ignore/ignoreFaPatchFields.C
 $(derivedFaPatchFields)/clampedPlate/clampedPlateFaPatchFields.C
 
diff --git a/src/finiteArea/fields/faPatchFields/derived/uniformFixedGradient/uniformFixedGradientFaPatchField.C b/src/finiteArea/fields/faPatchFields/derived/uniformFixedGradient/uniformFixedGradientFaPatchField.C
new file mode 100644
index 0000000000000000000000000000000000000000..f74667880c15522bcf58086fa5436a1ac26bd460
--- /dev/null
+++ b/src/finiteArea/fields/faPatchFields/derived/uniformFixedGradient/uniformFixedGradientFaPatchField.C
@@ -0,0 +1,161 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | www.openfoam.com
+     \\/     M anipulation  |
+-------------------------------------------------------------------------------
+    Copyright (C) 2023 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 "uniformFixedGradientFaPatchField.H"
+
+// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
+
+template<class Type>
+Foam::uniformFixedGradientFaPatchField<Type>::uniformFixedGradientFaPatchField
+(
+    const faPatch& p,
+    const DimensionedField<Type, areaMesh>& iF
+)
+:
+    fixedGradientFaPatchField<Type>(p, iF),
+    refGradFunc_(nullptr)
+{}
+
+
+template<class Type>
+Foam::uniformFixedGradientFaPatchField<Type>::uniformFixedGradientFaPatchField
+(
+    const faPatch& p,
+    const DimensionedField<Type, areaMesh>& iF,
+    const Field<Type>& fld
+)
+:
+    fixedGradientFaPatchField<Type>(p, iF, fld),
+    refGradFunc_(nullptr)
+{}
+
+
+template<class Type>
+Foam::uniformFixedGradientFaPatchField<Type>::uniformFixedGradientFaPatchField
+(
+    const faPatch& p,
+    const DimensionedField<Type, areaMesh>& iF,
+    const dictionary& dict
+)
+:
+    fixedGradientFaPatchField<Type>(p, iF),  // Bypass dictionary constructor
+    refGradFunc_
+    (
+        Function1<Type>::New
+        (
+            /* p.patch(), */
+            "uniformGradient",
+            dict
+        )
+    )
+{
+    this->evaluate();
+}
+
+
+template<class Type>
+Foam::uniformFixedGradientFaPatchField<Type>::uniformFixedGradientFaPatchField
+(
+    const uniformFixedGradientFaPatchField<Type>& ptf,
+    const faPatch& p,
+    const DimensionedField<Type, areaMesh>& iF,
+    const faPatchFieldMapper& mapper
+)
+:
+    fixedGradientFaPatchField<Type>(ptf, p, iF, mapper),
+    refGradFunc_(ptf.refGradFunc_.clone(/*p.patch()*/))
+{}
+
+
+template<class Type>
+Foam::uniformFixedGradientFaPatchField<Type>::uniformFixedGradientFaPatchField
+(
+    const uniformFixedGradientFaPatchField<Type>& ptf
+)
+:
+    fixedGradientFaPatchField<Type>(ptf),
+    refGradFunc_(ptf.refGradFunc_.clone(/*p.patch()*/))
+{}
+
+
+template<class Type>
+Foam::uniformFixedGradientFaPatchField<Type>::uniformFixedGradientFaPatchField
+(
+    const uniformFixedGradientFaPatchField<Type>& ptf,
+    const DimensionedField<Type, areaMesh>& iF
+)
+:
+    fixedGradientFaPatchField<Type>(ptf, iF),
+    refGradFunc_(ptf.refGradFunc_.clone(/*this->patch().patch()*/))
+{
+    // Evaluate the profile if defined
+    if (refGradFunc_)
+    {
+        this->evaluate();
+    }
+}
+
+
+// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
+
+template<class Type>
+void Foam::uniformFixedGradientFaPatchField<Type>::updateCoeffs()
+{
+    if (this->updated())
+    {
+        return;
+    }
+
+    const scalar t = this->db().time().timeOutputValue();
+
+    // Extra safety on the evaluation
+    if (refGradFunc_)
+    {
+        this->gradient() = refGradFunc_->value(t);
+    }
+    else
+    {
+        this->gradient() = Zero;
+    }
+
+    fixedGradientFaPatchField<Type>::updateCoeffs();
+}
+
+
+template<class Type>
+void Foam::uniformFixedGradientFaPatchField<Type>::write(Ostream& os) const
+{
+    fixedGradientFaPatchField<Type>::write(os);
+    if (refGradFunc_)
+    {
+        refGradFunc_->writeData(os);
+    }
+    faPatchField<Type>::writeValueEntry(os);
+}
+
+
+// ************************************************************************* //
diff --git a/src/finiteArea/fields/faPatchFields/derived/uniformFixedGradient/uniformFixedGradientFaPatchField.H b/src/finiteArea/fields/faPatchFields/derived/uniformFixedGradient/uniformFixedGradientFaPatchField.H
new file mode 100644
index 0000000000000000000000000000000000000000..52e0f026b63304cb8af0115e9e05da48ca799d19
--- /dev/null
+++ b/src/finiteArea/fields/faPatchFields/derived/uniformFixedGradient/uniformFixedGradientFaPatchField.H
@@ -0,0 +1,188 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | www.openfoam.com
+     \\/     M anipulation  |
+-------------------------------------------------------------------------------
+    Copyright (C) 2023 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::uniformFixedGradientFaPatchField
+
+Group
+    grpGenericBoundaryConditions
+
+Description
+    This boundary condition provides a uniform fixed gradient condition.
+
+Usage
+    \table
+        Property     | Description                  | Required | Default
+        uniformGradient | uniform gradient          | yes |
+    \endtable
+
+    Example of the boundary condition specification:
+    \verbatim
+    <patchName>
+    {
+        type            uniformFixedGradient;
+        uniformGradient constant 0.2;
+    }
+    \endverbatim
+
+Note
+    The uniformGradient entry is a Function1 type.
+    The example above gives the usage for supplying a constant value.
+
+See also
+    Foam::Function1Types
+    Foam::fixedGradientFaPatchField
+
+SourceFiles
+    uniformFixedGradientFaPatchField.C
+
+\*---------------------------------------------------------------------------*/
+
+#ifndef Foam_uniformFixedGradientFaPatchField_H
+#define Foam_uniformFixedGradientFaPatchField_H
+
+#include "fixedGradientFaPatchField.H"
+#include "PatchFunction1.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+
+/*---------------------------------------------------------------------------*\
+              Class uniformFixedGradientFaPatchField Declaration
+\*---------------------------------------------------------------------------*/
+
+template<class Type>
+class uniformFixedGradientFaPatchField
+:
+    public fixedGradientFaPatchField<Type>
+{
+    // Private Data
+
+        //- Function providing the gradient
+        autoPtr<Function1<Type>> refGradFunc_;
+
+
+public:
+
+    //- Runtime type information
+    TypeName("uniformFixedGradient");
+
+
+    // Constructors
+
+        //- Construct from patch and internal field
+        uniformFixedGradientFaPatchField
+        (
+            const faPatch&,
+            const DimensionedField<Type, areaMesh>&
+        );
+
+        //- Construct from patch and internal field and patch field
+        uniformFixedGradientFaPatchField
+        (
+            const faPatch&,
+            const DimensionedField<Type, areaMesh>&,
+            const Field<Type>& fld
+        );
+
+        //- Construct from patch, internal field and dictionary
+        uniformFixedGradientFaPatchField
+        (
+            const faPatch&,
+            const DimensionedField<Type, areaMesh>&,
+            const dictionary&
+        );
+
+        //- Construct by mapping onto a new patch
+        uniformFixedGradientFaPatchField
+        (
+            const uniformFixedGradientFaPatchField<Type>&,
+            const faPatch&,
+            const DimensionedField<Type, areaMesh>&,
+            const faPatchFieldMapper&
+        );
+
+        //- Construct as copy
+        uniformFixedGradientFaPatchField
+        (
+            const uniformFixedGradientFaPatchField<Type>&
+        );
+
+        //- Construct and return a clone
+        virtual tmp<faPatchField<Type>> clone() const
+        {
+            return tmp<faPatchField<Type>>
+            (
+                new uniformFixedGradientFaPatchField<Type>(*this)
+            );
+        }
+
+        //- Construct as copy setting internal field reference
+        uniformFixedGradientFaPatchField
+        (
+            const uniformFixedGradientFaPatchField<Type>&,
+            const DimensionedField<Type, areaMesh>&
+        );
+
+        //- Construct and return a clone setting internal field reference
+        virtual tmp<faPatchField<Type>> clone
+        (
+            const DimensionedField<Type, areaMesh>& iF
+        ) const
+        {
+            return tmp<faPatchField<Type>>
+            (
+                new uniformFixedGradientFaPatchField<Type>(*this, iF)
+            );
+        }
+
+
+    // Member functions
+
+        //- Update the coefficients associated with the patch field
+        virtual void updateCoeffs();
+
+        //- Write
+        virtual void write(Ostream& os) const;
+};
+
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+} // End namespace Foam
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#ifdef NoRepository
+    #include "uniformFixedGradientFaPatchField.C"
+#endif
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#endif
+
+// ************************************************************************* //
diff --git a/src/finiteArea/fields/faPatchFields/derived/uniformFixedGradient/uniformFixedGradientFaPatchFields.C b/src/finiteArea/fields/faPatchFields/derived/uniformFixedGradient/uniformFixedGradientFaPatchFields.C
new file mode 100644
index 0000000000000000000000000000000000000000..2e6a854e5ce3b1eefc27680c61ce4117bc72e41e
--- /dev/null
+++ b/src/finiteArea/fields/faPatchFields/derived/uniformFixedGradient/uniformFixedGradientFaPatchFields.C
@@ -0,0 +1,40 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | www.openfoam.com
+     \\/     M anipulation  |
+-------------------------------------------------------------------------------
+    Copyright (C) 2023 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 "uniformFixedGradientFaPatchFields.H"
+#include "addToRunTimeSelectionTable.H"
+#include "areaFields.H"
+#include "edgeFields.H"
+
+// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
+
+namespace Foam
+{
+    makeFaPatchFields(uniformFixedGradient);
+}
+
+// ************************************************************************* //
diff --git a/src/finiteArea/fields/faPatchFields/derived/uniformFixedGradient/uniformFixedGradientFaPatchFields.H b/src/finiteArea/fields/faPatchFields/derived/uniformFixedGradient/uniformFixedGradientFaPatchFields.H
new file mode 100644
index 0000000000000000000000000000000000000000..7e7a2973d648f29d6a9583985af86f251df580e1
--- /dev/null
+++ b/src/finiteArea/fields/faPatchFields/derived/uniformFixedGradient/uniformFixedGradientFaPatchFields.H
@@ -0,0 +1,50 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | www.openfoam.com
+     \\/     M anipulation  |
+-------------------------------------------------------------------------------
+    Copyright (C) 2023 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/>.
+
+\*---------------------------------------------------------------------------*/
+
+#ifndef Foam_uniformFixedGradientFaPatchFields_H
+#define Foam_uniformFixedGradientFaPatchFields_H
+
+#include "uniformFixedGradientFaPatchField.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+makeFaPatchTypeFieldTypedefs(uniformFixedGradient);
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+} // End namespace Foam
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#endif
+
+// ************************************************************************* //
diff --git a/src/finiteArea/fields/faPatchFields/derived/uniformFixedValue/uniformFixedValueFaPatchField.C b/src/finiteArea/fields/faPatchFields/derived/uniformFixedValue/uniformFixedValueFaPatchField.C
index a02dd23823b610b981b5629d3a185b2899416235..804b049a5f52d4f61f87336e8da4f54da11bb3bf 100644
--- a/src/finiteArea/fields/faPatchFields/derived/uniformFixedValue/uniformFixedValueFaPatchField.C
+++ b/src/finiteArea/fields/faPatchFields/derived/uniformFixedValue/uniformFixedValueFaPatchField.C
@@ -73,11 +73,7 @@ Foam::uniformFixedValueFaPatchField<Type>::uniformFixedValueFaPatchField
         )
     )
 {
-    if (this->readValueEntry(dict))
-    {
-        // Full restart
-    }
-    else
+    if (!this->readValueEntry(dict))
     {
         this->evaluate();
     }
diff --git a/src/finiteArea/fields/faPatchFields/derived/uniformMixed/uniformMixedFaPatchField.C b/src/finiteArea/fields/faPatchFields/derived/uniformMixed/uniformMixedFaPatchField.C
new file mode 100644
index 0000000000000000000000000000000000000000..82ee838495f353225b430395911112e0d0d7ca64
--- /dev/null
+++ b/src/finiteArea/fields/faPatchFields/derived/uniformMixed/uniformMixedFaPatchField.C
@@ -0,0 +1,245 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | www.openfoam.com
+     \\/     M anipulation  |
+-------------------------------------------------------------------------------
+    Copyright (C) 2023 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 "uniformMixedFaPatchField.H"
+
+// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
+
+template<class Type>
+Foam::uniformMixedFaPatchField<Type>::uniformMixedFaPatchField
+(
+    const faPatch& p,
+    const DimensionedField<Type, areaMesh>& iF
+)
+:
+    mixedFaPatchField<Type>(p, iF),
+    refValueFunc_(nullptr),
+    refGradFunc_(nullptr),
+    valueFractionFunc_(nullptr)
+{}
+
+
+template<class Type>
+Foam::uniformMixedFaPatchField<Type>::uniformMixedFaPatchField
+(
+    const faPatch& p,
+    const DimensionedField<Type, areaMesh>& iF,
+    const Field<Type>& fld
+)
+:
+    mixedFaPatchField<Type>(p, iF, fld),
+    refValueFunc_(nullptr),
+    refGradFunc_(nullptr),
+    valueFractionFunc_(nullptr)
+{}
+
+
+template<class Type>
+Foam::uniformMixedFaPatchField<Type>::uniformMixedFaPatchField
+(
+    const faPatch& p,
+    const DimensionedField<Type, areaMesh>& iF,
+    const dictionary& dict
+)
+:
+    mixedFaPatchField<Type>(p, iF, dict, IOobjectOption::NO_READ),
+    refValueFunc_
+    (
+        Function1<Type>::NewIfPresent
+        (
+            /* p.patch(), */
+            "uniformValue",
+            dict
+        )
+    ),
+    refGradFunc_
+    (
+        Function1<Type>::NewIfPresent
+        (
+            // p.patch(),
+            "uniformGradient",
+            dict
+        )
+    ),
+    valueFractionFunc_(nullptr)
+{
+    if (refValueFunc_)
+    {
+        if (refGradFunc_)
+        {
+            // Both value + gradient: needs valueFraction
+            valueFractionFunc_.reset
+            (
+                Function1<scalar>::New
+                (
+                    /* p.patch(), */
+                    "uniformValueFraction",
+                    dict
+                )
+            );
+        }
+    }
+    else if (!refGradFunc_)
+    {
+        // Missing both value and gradient: FatalIOError
+        FatalIOErrorInFunction(dict)
+            << "For " << this->internalField().name() << " on "
+            << this->patch().name() << nl
+            << "Require either or both: uniformValue and uniformGradient"
+            << " (possibly uniformValueFraction as well)" << nl
+            << exit(FatalIOError);
+    }
+
+    // Use restart value if provided...
+    if (!this->readValueEntry(dict))
+    {
+        this->evaluate();
+    }
+}
+
+
+template<class Type>
+Foam::uniformMixedFaPatchField<Type>::uniformMixedFaPatchField
+(
+    const uniformMixedFaPatchField<Type>& ptf,
+    const faPatch& p,
+    const DimensionedField<Type, areaMesh>& iF,
+    const faPatchFieldMapper& mapper
+)
+:
+    mixedFaPatchField<Type>(ptf, p, iF, mapper),
+    refValueFunc_(ptf.refValueFunc_.clone(/*p.patch()*/)),
+    refGradFunc_(ptf.refGradFunc_.clone(/*p.patch()*/)),
+    valueFractionFunc_(ptf.valueFractionFunc_.clone(/*p.patch()*/))
+{}
+
+
+template<class Type>
+Foam::uniformMixedFaPatchField<Type>::uniformMixedFaPatchField
+(
+    const uniformMixedFaPatchField<Type>& ptf
+)
+:
+    mixedFaPatchField<Type>(ptf),
+    refValueFunc_(ptf.refValueFunc_.clone(/*this->patch().patch()*/)),
+    refGradFunc_(ptf.refGradFunc_.clone(/*this->patch().patch()*/)),
+    valueFractionFunc_(ptf.valueFractionFunc_.clone(/*this->patch().patch()*/))
+{}
+
+
+template<class Type>
+Foam::uniformMixedFaPatchField<Type>::uniformMixedFaPatchField
+(
+    const uniformMixedFaPatchField<Type>& ptf,
+    const DimensionedField<Type, areaMesh>& iF
+)
+:
+    mixedFaPatchField<Type>(ptf, iF),
+    refValueFunc_(ptf.refValueFunc_.clone(/*this->patch().patch()*/)),
+    refGradFunc_(ptf.refGradFunc_.clone(/*this->patch().patch()*/)),
+    valueFractionFunc_(ptf.valueFractionFunc_.clone(/*this->patch().patch()*/))
+{
+    // Evaluate the profile if defined
+    if (ptf.refValueFunc_ || ptf.refGradFunc_)
+    {
+        this->evaluate();
+    }
+}
+
+
+// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
+
+template<class Type>
+void Foam::uniformMixedFaPatchField<Type>::updateCoeffs()
+{
+    if (this->updated())
+    {
+        return;
+    }
+
+    const scalar t = this->db().time().timeOutputValue();
+
+    if (refValueFunc_)
+    {
+        this->refValue() = refValueFunc_->value(t);
+
+        if (refGradFunc_)
+        {
+            // Both value + gradient: has valueFraction too
+            this->valueFraction() = valueFractionFunc_->value(t);
+        }
+        else
+        {
+            // Has value only
+            this->valueFraction() = 1;
+        }
+    }
+    else
+    {
+        this->refValue() = Zero;
+        this->valueFraction() = 0;
+    }
+    if (refGradFunc_)
+    {
+        this->refGrad() = refGradFunc_->value(t);
+    }
+    else
+    {
+        this->refGrad() = Zero;
+    }
+
+    // Missing both value and gradient is caught as an error in
+    // dictionary constructor, but treated as zero-gradient here.
+
+    mixedFaPatchField<Type>::updateCoeffs();
+}
+
+
+template<class Type>
+void Foam::uniformMixedFaPatchField<Type>::write(Ostream& os) const
+{
+    faPatchField<Type>::write(os);
+
+    if (refValueFunc_)
+    {
+        refValueFunc_->writeData(os);
+    }
+    if (refGradFunc_)
+    {
+        refGradFunc_->writeData(os);
+    }
+    if (valueFractionFunc_)
+    {
+        valueFractionFunc_->writeData(os);
+    }
+
+    // Eg, for visualisation
+    faPatchField<Type>::writeValueEntry(os);
+}
+
+
+// ************************************************************************* //
diff --git a/src/finiteArea/fields/faPatchFields/derived/uniformMixed/uniformMixedFaPatchField.H b/src/finiteArea/fields/faPatchFields/derived/uniformMixed/uniformMixedFaPatchField.H
new file mode 100644
index 0000000000000000000000000000000000000000..bd2819f6f422d43a28f47f5ec812731d4cc59cdb
--- /dev/null
+++ b/src/finiteArea/fields/faPatchFields/derived/uniformMixed/uniformMixedFaPatchField.H
@@ -0,0 +1,210 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | www.openfoam.com
+     \\/     M anipulation  |
+-------------------------------------------------------------------------------
+    Copyright (C) 2023 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::uniformMixedFaPatchField
+
+Group
+    grpGenericBoundaryConditions
+
+Description
+    This boundary condition provides 'mixed' type boundary condition
+    that mix a \em uniform fixed value and a \em uniform patch-normal
+    gradient condition. The term "uniform" is a legacy name since the
+    prescribed values were previously spatially uniform across that patch.
+
+Usage
+    \table
+        Property     | Description                      | Required | Default
+        uniformValue | uniform value                    | partly | 0
+        uniformGradient | uniform gradient              | partly | 0
+        uniformValueFraction | uniform valueFraction    | partly | depends
+    \endtable
+
+    Example of the boundary condition specification:
+    \verbatim
+    <patchName>
+    {
+        type            uniformMixed;
+        uniformValue    constant 0.2;
+        uniformGradient constant 0.2;
+        uniformValueFraction
+        {
+            type  sine;
+            ...
+        }
+    }
+    \endverbatim
+
+Note
+    This boundary condition allows \em lazier definitions so that either
+    or both: \c uniformValue and \c uniformGradient must be defined.
+    If only of these entries is defined, the value fraction is automatically
+    treated appropriately (ie, 0 with \c uniformGradient and 1 with
+    uniformValue).
+    If both \c uniformValue and \c uniformGradient are defined,
+    the \c uniformValueFraction must also be defined.
+
+See also
+    Foam::Function1Types
+    Foam::mixedFaPatchField
+
+SourceFiles
+    uniformMixedFaPatchField.C
+
+\*---------------------------------------------------------------------------*/
+
+#ifndef Foam_uniformMixedFaPatchField_H
+#define Foam_uniformMixedFaPatchField_H
+
+#include "mixedFaPatchField.H"
+#include "Function1.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+
+/*---------------------------------------------------------------------------*\
+                  Class uniformMixedFaPatchField Declaration
+\*---------------------------------------------------------------------------*/
+
+template<class Type>
+class uniformMixedFaPatchField
+:
+    public mixedFaPatchField<Type>
+{
+    // Private Data
+
+        //- Function providing the value
+        autoPtr<Function1<Type>> refValueFunc_;
+
+        //- Function providing the gradient
+        autoPtr<Function1<Type>> refGradFunc_;
+
+        //- Function providing the value-fraction
+        autoPtr<Function1<scalar>> valueFractionFunc_;
+
+
+public:
+
+    //- Runtime type information
+    TypeName("uniformMixed");
+
+
+    // Constructors
+
+        //- Construct from patch and internal field
+        uniformMixedFaPatchField
+        (
+            const faPatch&,
+            const DimensionedField<Type, areaMesh>&
+        );
+
+        //- Construct from patch and internal field and patch field
+        uniformMixedFaPatchField
+        (
+            const faPatch&,
+            const DimensionedField<Type, areaMesh>&,
+            const Field<Type>& fld
+        );
+
+        //- Construct from patch, internal field and dictionary
+        uniformMixedFaPatchField
+        (
+            const faPatch&,
+            const DimensionedField<Type, areaMesh>&,
+            const dictionary&
+        );
+
+        //- Construct by mapping onto a new patch
+        uniformMixedFaPatchField
+        (
+            const uniformMixedFaPatchField<Type>&,
+            const faPatch&,
+            const DimensionedField<Type, areaMesh>&,
+            const faPatchFieldMapper&
+        );
+
+        //- Construct as copy
+        uniformMixedFaPatchField
+        (
+            const uniformMixedFaPatchField<Type>&
+        );
+
+        //- Construct and return a clone
+        virtual tmp<faPatchField<Type>> clone() const
+        {
+            return tmp<faPatchField<Type>>
+            (
+                new uniformMixedFaPatchField<Type>(*this)
+            );
+        }
+
+        //- Construct as copy setting internal field reference
+        uniformMixedFaPatchField
+        (
+            const uniformMixedFaPatchField<Type>&,
+            const DimensionedField<Type, areaMesh>&
+        );
+
+        //- Construct and return a clone setting internal field reference
+        virtual tmp<faPatchField<Type>> clone
+        (
+            const DimensionedField<Type, areaMesh>& iF
+        ) const
+        {
+            return tmp<faPatchField<Type>>
+            (
+                new uniformMixedFaPatchField<Type>(*this, iF)
+            );
+        }
+
+
+    // Member Functions
+
+        //- Update the coefficients associated with the patch field
+        virtual void updateCoeffs();
+
+        //- Write
+        virtual void write(Ostream& os) const;
+};
+
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+} // End namespace Foam
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#ifdef NoRepository
+    #include "uniformMixedFaPatchField.C"
+#endif
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#endif
+
+// ************************************************************************* //
diff --git a/src/finiteArea/fields/faPatchFields/derived/uniformMixed/uniformMixedFaPatchFields.C b/src/finiteArea/fields/faPatchFields/derived/uniformMixed/uniformMixedFaPatchFields.C
new file mode 100644
index 0000000000000000000000000000000000000000..d7f518c1a11f306f91658e48576e702fbab6cf47
--- /dev/null
+++ b/src/finiteArea/fields/faPatchFields/derived/uniformMixed/uniformMixedFaPatchFields.C
@@ -0,0 +1,40 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | www.openfoam.com
+     \\/     M anipulation  |
+-------------------------------------------------------------------------------
+    Copyright (C) 2023 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 "uniformMixedFaPatchFields.H"
+#include "addToRunTimeSelectionTable.H"
+#include "areaFields.H"
+#include "edgeFields.H"
+
+// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
+
+namespace Foam
+{
+    makeFaPatchFields(uniformMixed);
+}
+
+// ************************************************************************* //
diff --git a/src/finiteArea/fields/faPatchFields/derived/uniformMixed/uniformMixedFaPatchFields.H b/src/finiteArea/fields/faPatchFields/derived/uniformMixed/uniformMixedFaPatchFields.H
new file mode 100644
index 0000000000000000000000000000000000000000..4f91619e41d27e597f347783de5ea19fd4e4cdb8
--- /dev/null
+++ b/src/finiteArea/fields/faPatchFields/derived/uniformMixed/uniformMixedFaPatchFields.H
@@ -0,0 +1,50 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | www.openfoam.com
+     \\/     M anipulation  |
+-------------------------------------------------------------------------------
+    Copyright (C) 2023 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/>.
+
+\*---------------------------------------------------------------------------*/
+
+#ifndef Foam_uniformMixedFaPatchFields_H
+#define Foam_uniformMixedFaPatchFields_H
+
+#include "uniformMixedFaPatchField.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+makeFaPatchTypeFieldTypedefs(uniformMixed);
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+} // End namespace Foam
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#endif
+
+// ************************************************************************* //
diff --git a/src/finiteVolume/Make/files b/src/finiteVolume/Make/files
index aa58f5991cba79fb07efb5659f2f4b11f5620f8b..dca1ff2d9f89b4bc3e0b7e41cd8f2b83ffd6159c 100644
--- a/src/finiteVolume/Make/files
+++ b/src/finiteVolume/Make/files
@@ -227,6 +227,7 @@ $(derivedFvPatchFields)/turbulentIntensityKineticEnergyInlet/turbulentIntensityK
 $(derivedFvPatchFields)/uniformDensityHydrostaticPressure/uniformDensityHydrostaticPressureFvPatchScalarField.C
 $(derivedFvPatchFields)/uniformFixedGradient/uniformFixedGradientFvPatchFields.C
 $(derivedFvPatchFields)/uniformFixedValue/uniformFixedValueFvPatchFields.C
+$(derivedFvPatchFields)/uniformMixed/uniformMixedFvPatchFields.C
 $(derivedFvPatchFields)/uniformInletOutlet/uniformInletOutletFvPatchFields.C
 $(derivedFvPatchFields)/uniformJump/uniformJumpFvPatchFields.C
 $(derivedFvPatchFields)/uniformJumpAMI/uniformJumpAMIFvPatchFields.C
diff --git a/src/finiteVolume/fields/fvPatchFields/derived/uniformFixedGradient/uniformFixedGradientFvPatchField.C b/src/finiteVolume/fields/fvPatchFields/derived/uniformFixedGradient/uniformFixedGradientFvPatchField.C
index 9be8715829b4fb56e5f4c9f94fd588daf839c137..d265cadd197dc2d346ae59850e58ebc11827423d 100644
--- a/src/finiteVolume/fields/fvPatchFields/derived/uniformFixedGradient/uniformFixedGradientFvPatchField.C
+++ b/src/finiteVolume/fields/fvPatchFields/derived/uniformFixedGradient/uniformFixedGradientFvPatchField.C
@@ -6,7 +6,7 @@
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
     Copyright (C) 2013-2016 OpenFOAM Foundation
-    Copyright (C) 2017-2020 OpenCFD Ltd.
+    Copyright (C) 2017-2023 OpenCFD Ltd.
 -------------------------------------------------------------------------------
 License
     This file is part of OpenFOAM.
@@ -38,7 +38,7 @@ Foam::uniformFixedGradientFvPatchField<Type>::uniformFixedGradientFvPatchField
 )
 :
     fixedGradientFvPatchField<Type>(p, iF),
-    uniformGradient_(nullptr)
+    refGradFunc_(nullptr)
 {}
 
 
@@ -51,7 +51,7 @@ Foam::uniformFixedGradientFvPatchField<Type>::uniformFixedGradientFvPatchField
 )
 :
     fixedGradientFvPatchField<Type>(p, iF, fld),
-    uniformGradient_(nullptr)
+    refGradFunc_(nullptr)
 {}
 
 
@@ -64,7 +64,7 @@ Foam::uniformFixedGradientFvPatchField<Type>::uniformFixedGradientFvPatchField
 )
 :
     fixedGradientFvPatchField<Type>(p, iF),  // Bypass dictionary constructor
-    uniformGradient_
+    refGradFunc_
     (
         PatchFunction1<Type>::New(p.patch(), "uniformGradient", dict)
     )
@@ -84,7 +84,7 @@ Foam::uniformFixedGradientFvPatchField<Type>::uniformFixedGradientFvPatchField
 )
 :
     fixedGradientFvPatchField<Type>(ptf, p, iF, mapper),
-    uniformGradient_(ptf.uniformGradient_.clone())
+    refGradFunc_(ptf.refGradFunc_.clone())
 {}
 
 
@@ -95,7 +95,7 @@ Foam::uniformFixedGradientFvPatchField<Type>::uniformFixedGradientFvPatchField
 )
 :
     fixedGradientFvPatchField<Type>(ptf),
-    uniformGradient_(ptf.uniformGradient_.clone())
+    refGradFunc_(ptf.refGradFunc_.clone())
 {}
 
 
@@ -107,10 +107,10 @@ Foam::uniformFixedGradientFvPatchField<Type>::uniformFixedGradientFvPatchField
 )
 :
     fixedGradientFvPatchField<Type>(ptf, iF),
-    uniformGradient_(ptf.uniformGradient_.clone())
+    refGradFunc_(ptf.refGradFunc_.clone())
 {
     // Evaluate the profile if defined
-    if (ptf.uniformGradient_)
+    if (refGradFunc_)
     {
         this->evaluate();
     }
@@ -128,7 +128,15 @@ void Foam::uniformFixedGradientFvPatchField<Type>::updateCoeffs()
     }
 
     const scalar t = this->db().time().timeOutputValue();
-    this->gradient() = uniformGradient_->value(t);
+    // Extra safety on the evaluation
+    if (refGradFunc_)
+    {
+        this->gradient() = refGradFunc_->value(t);
+    }
+    else
+    {
+        this->gradient() = Zero;
+    }
 
     fixedGradientFvPatchField<Type>::updateCoeffs();
 }
@@ -138,7 +146,10 @@ template<class Type>
 void Foam::uniformFixedGradientFvPatchField<Type>::write(Ostream& os) const
 {
     fixedGradientFvPatchField<Type>::write(os);
-    uniformGradient_->writeData(os);
+    if (refGradFunc_)
+    {
+        refGradFunc_->writeData(os);
+    }
     fvPatchField<Type>::writeValueEntry(os);
 }
 
diff --git a/src/finiteVolume/fields/fvPatchFields/derived/uniformFixedGradient/uniformFixedGradientFvPatchField.H b/src/finiteVolume/fields/fvPatchFields/derived/uniformFixedGradient/uniformFixedGradientFvPatchField.H
index 0cce531299c893a1769099c3ff277bcae32dad90..34bdbb50a0cbe59b86f60ddb91283c264db0be16 100644
--- a/src/finiteVolume/fields/fvPatchFields/derived/uniformFixedGradient/uniformFixedGradientFvPatchField.H
+++ b/src/finiteVolume/fields/fvPatchFields/derived/uniformFixedGradient/uniformFixedGradientFvPatchField.H
@@ -6,7 +6,7 @@
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
     Copyright (C) 2013-2016 OpenFOAM Foundation
-    Copyright (C) 2020 OpenCFD Ltd.
+    Copyright (C) 2020-2023 OpenCFD Ltd.
 -------------------------------------------------------------------------------
 License
     This file is part of OpenFOAM.
@@ -62,8 +62,8 @@ SourceFiles
 
 \*---------------------------------------------------------------------------*/
 
-#ifndef uniformFixedGradientFvPatchField_H
-#define uniformFixedGradientFvPatchField_H
+#ifndef Foam_uniformFixedGradientFvPatchField_H
+#define Foam_uniformFixedGradientFvPatchField_H
 
 #include "fixedGradientFvPatchFields.H"
 #include "PatchFunction1.H"
@@ -84,8 +84,8 @@ class uniformFixedGradientFvPatchField
 {
     // Private Data
 
-        //- Gradient
-        autoPtr<PatchFunction1<Type>> uniformGradient_;
+        //- Function providing the gradient
+        autoPtr<PatchFunction1<Type>> refGradFunc_;
 
 
 public:
diff --git a/src/finiteVolume/fields/fvPatchFields/derived/uniformFixedValue/uniformFixedValueFvPatchField.C b/src/finiteVolume/fields/fvPatchFields/derived/uniformFixedValue/uniformFixedValueFvPatchField.C
index c1ba993d44631ecb2ac237d2fef50838916d2da9..4e876b1a07582fcdb5a975a1a865dce4ec3d2b9f 100644
--- a/src/finiteVolume/fields/fvPatchFields/derived/uniformFixedValue/uniformFixedValueFvPatchField.C
+++ b/src/finiteVolume/fields/fvPatchFields/derived/uniformFixedValue/uniformFixedValueFvPatchField.C
@@ -134,12 +134,12 @@ void Foam::uniformFixedValueFvPatchField<Type>::autoMap
     if (refValueFunc_)
     {
         refValueFunc_().autoMap(mapper);
-    }
 
-    if (refValueFunc_().constant())
-    {
-        // If mapper is not dependent on time we're ok to evaluate
-        this->evaluate();
+        if (refValueFunc_().constant())
+        {
+            // If mapper is not dependent on time we're ok to evaluate
+            this->evaluate();
+        }
     }
 }
 
@@ -156,7 +156,7 @@ void Foam::uniformFixedValueFvPatchField<Type>::rmap
     const uniformFixedValueFvPatchField& tiptf =
         refCast<const uniformFixedValueFvPatchField>(ptf);
 
-    if (refValueFunc_)
+    if (refValueFunc_ && tiptf.refValueFunc_)
     {
         refValueFunc_().rmap(tiptf.refValueFunc_(), addr);
     }
@@ -172,6 +172,7 @@ void Foam::uniformFixedValueFvPatchField<Type>::updateCoeffs()
     }
 
     const scalar t = this->db().time().timeOutputValue();
+
     fvPatchField<Type>::operator==(refValueFunc_->value(t));
     fixedValueFvPatchField<Type>::updateCoeffs();
 }
@@ -181,12 +182,10 @@ template<class Type>
 void Foam::uniformFixedValueFvPatchField<Type>::write(Ostream& os) const
 {
     fvPatchField<Type>::write(os);
-
     if (refValueFunc_)
     {
         refValueFunc_->writeData(os);
     }
-
     fvPatchField<Type>::writeValueEntry(os);
 }
 
diff --git a/src/finiteVolume/fields/fvPatchFields/derived/uniformFixedValue/uniformFixedValueFvPatchField.H b/src/finiteVolume/fields/fvPatchFields/derived/uniformFixedValue/uniformFixedValueFvPatchField.H
index f36be7f9765f74ea4cefb956102a79983366b910..fa5f8a344205b1eb998d0570beb12c91bdec7da5 100644
--- a/src/finiteVolume/fields/fvPatchFields/derived/uniformFixedValue/uniformFixedValueFvPatchField.H
+++ b/src/finiteVolume/fields/fvPatchFields/derived/uniformFixedValue/uniformFixedValueFvPatchField.H
@@ -6,7 +6,7 @@
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
     Copyright (C) 2011-2016 OpenFOAM Foundation
-    Copyright (C) 2018-2020 OpenCFD Ltd.
+    Copyright (C) 2018-2023 OpenCFD Ltd.
 -------------------------------------------------------------------------------
 License
     This file is part of OpenFOAM.
diff --git a/src/finiteVolume/fields/fvPatchFields/derived/uniformMixed/uniformMixedFvPatchField.C b/src/finiteVolume/fields/fvPatchFields/derived/uniformMixed/uniformMixedFvPatchField.C
new file mode 100644
index 0000000000000000000000000000000000000000..c73934f4ec1d718b60634c57a5fd7bd1db121491
--- /dev/null
+++ b/src/finiteVolume/fields/fvPatchFields/derived/uniformMixed/uniformMixedFvPatchField.C
@@ -0,0 +1,235 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | www.openfoam.com
+     \\/     M anipulation  |
+-------------------------------------------------------------------------------
+    Copyright (C) 2023 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 "uniformMixedFvPatchField.H"
+
+// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
+
+template<class Type>
+Foam::uniformMixedFvPatchField<Type>::uniformMixedFvPatchField
+(
+    const fvPatch& p,
+    const DimensionedField<Type, volMesh>& iF
+)
+:
+    mixedFvPatchField<Type>(p, iF),
+    refValueFunc_(nullptr),
+    refGradFunc_(nullptr),
+    valueFractionFunc_(nullptr)
+{}
+
+
+template<class Type>
+Foam::uniformMixedFvPatchField<Type>::uniformMixedFvPatchField
+(
+    const fvPatch& p,
+    const DimensionedField<Type, volMesh>& iF,
+    const Field<Type>& fld
+)
+:
+    mixedFvPatchField<Type>(p, iF, fld),
+    refValueFunc_(nullptr),
+    refGradFunc_(nullptr),
+    valueFractionFunc_(nullptr)
+{}
+
+
+template<class Type>
+Foam::uniformMixedFvPatchField<Type>::uniformMixedFvPatchField
+(
+    const fvPatch& p,
+    const DimensionedField<Type, volMesh>& iF,
+    const dictionary& dict
+)
+:
+    mixedFvPatchField<Type>(p, iF, dict, IOobjectOption::NO_READ),
+    refValueFunc_
+    (
+        PatchFunction1<Type>::NewIfPresent(p.patch(), "uniformValue", dict)
+    ),
+    refGradFunc_
+    (
+        PatchFunction1<Type>::NewIfPresent(p.patch(), "uniformGradient", dict)
+    ),
+    valueFractionFunc_(nullptr)
+{
+    if (refValueFunc_)
+    {
+        if (refGradFunc_)
+        {
+            // Both value + gradient: needs valueFraction
+            valueFractionFunc_.reset
+            (
+                PatchFunction1<scalar>::New
+                (
+                    p.patch(),
+                    "uniformValueFraction",
+                    dict
+                )
+            );
+        }
+    }
+    else if (!refGradFunc_)
+    {
+        // Missing both value and gradient: FatalIOError
+        FatalIOErrorInFunction(dict)
+            << "For " << this->internalField().name() << " on "
+            << this->patch().name() << nl
+            << "Require either or both: uniformValue and uniformGradient"
+            << " (possibly uniformValueFraction as well)" << nl
+            << exit(FatalIOError);
+    }
+
+    // Use restart value if provided...
+    if (!this->readValueEntry(dict))
+    {
+        this->evaluate();
+    }
+}
+
+
+template<class Type>
+Foam::uniformMixedFvPatchField<Type>::uniformMixedFvPatchField
+(
+    const uniformMixedFvPatchField<Type>& ptf,
+    const fvPatch& p,
+    const DimensionedField<Type, volMesh>& iF,
+    const fvPatchFieldMapper& mapper
+)
+:
+    mixedFvPatchField<Type>(ptf, p, iF, mapper),
+    refValueFunc_(ptf.refValueFunc_.clone(p.patch())),
+    refGradFunc_(ptf.refGradFunc_.clone(p.patch())),
+    valueFractionFunc_(ptf.valueFractionFunc_.clone(p.patch()))
+{}
+
+
+template<class Type>
+Foam::uniformMixedFvPatchField<Type>::uniformMixedFvPatchField
+(
+    const uniformMixedFvPatchField<Type>& ptf
+)
+:
+    mixedFvPatchField<Type>(ptf),
+    refValueFunc_(ptf.refValueFunc_.clone(this->patch().patch())),
+    refGradFunc_(ptf.refGradFunc_.clone(this->patch().patch())),
+    valueFractionFunc_(ptf.valueFractionFunc_.clone(this->patch().patch()))
+{}
+
+
+template<class Type>
+Foam::uniformMixedFvPatchField<Type>::uniformMixedFvPatchField
+(
+    const uniformMixedFvPatchField<Type>& ptf,
+    const DimensionedField<Type, volMesh>& iF
+)
+:
+    mixedFvPatchField<Type>(ptf, iF),
+    refValueFunc_(ptf.refValueFunc_.clone(this->patch().patch())),
+    refGradFunc_(ptf.refGradFunc_.clone(this->patch().patch())),
+    valueFractionFunc_(ptf.valueFractionFunc_.clone(this->patch().patch()))
+{
+    // Evaluate the profile if defined
+    if (ptf.refValueFunc_ || ptf.refGradFunc_)
+    {
+        this->evaluate();
+    }
+}
+
+
+// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
+
+template<class Type>
+void Foam::uniformMixedFvPatchField<Type>::updateCoeffs()
+{
+    if (this->updated())
+    {
+        return;
+    }
+
+    const scalar t = this->db().time().timeOutputValue();
+
+    if (refValueFunc_)
+    {
+        this->refValue() = refValueFunc_->value(t);
+
+        if (refGradFunc_)
+        {
+            // Both value + gradient: has valueFraction too
+            this->valueFraction() = valueFractionFunc_->value(t);
+        }
+        else
+        {
+            // Has value only
+            this->valueFraction() = 1;
+        }
+    }
+    else
+    {
+        this->refValue() = Zero;
+        this->valueFraction() = 0;
+    }
+    if (refGradFunc_)
+    {
+        this->refGrad() = refGradFunc_->value(t);
+    }
+    else
+    {
+        this->refGrad() = Zero;
+    }
+
+    // Missing both value and gradient is caught as an error in
+    // dictionary constructor, but treated as zero-gradient here.
+
+    mixedFvPatchField<Type>::updateCoeffs();
+}
+
+
+template<class Type>
+void Foam::uniformMixedFvPatchField<Type>::write(Ostream& os) const
+{
+    fvPatchField<Type>::write(os);
+
+    if (refValueFunc_)
+    {
+        refValueFunc_->writeData(os);
+    }
+    if (refGradFunc_)
+    {
+        refGradFunc_->writeData(os);
+    }
+    if (valueFractionFunc_)
+    {
+        valueFractionFunc_->writeData(os);
+    }
+
+    // Eg, for visualisation
+    fvPatchField<Type>::writeValueEntry(os);
+}
+
+
+// ************************************************************************* //
diff --git a/src/finiteVolume/fields/fvPatchFields/derived/uniformMixed/uniformMixedFvPatchField.H b/src/finiteVolume/fields/fvPatchFields/derived/uniformMixed/uniformMixedFvPatchField.H
new file mode 100644
index 0000000000000000000000000000000000000000..ccaa055b31af3332ba6973c8bf31d1331edad5eb
--- /dev/null
+++ b/src/finiteVolume/fields/fvPatchFields/derived/uniformMixed/uniformMixedFvPatchField.H
@@ -0,0 +1,213 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | www.openfoam.com
+     \\/     M anipulation  |
+-------------------------------------------------------------------------------
+    Copyright (C) 2023 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::uniformMixedFvPatchField
+
+Group
+    grpGenericBoundaryConditions
+
+Description
+    This boundary condition provides 'mixed' type boundary condition
+    that mix a \em uniform fixed value and a \em uniform patch-normal
+    gradient condition. The term "uniform" is a legacy name since the
+    prescribed values were previously spatially uniform across that patch.
+
+    In the meantime, a PatchFunction1 is used, which can have both
+    spatial and temporal dependencies.
+
+Usage
+    \table
+        Property     | Description                      | Required | Default
+        uniformValue | uniform value                    | partly | 0
+        uniformGradient | uniform gradient              | partly | 0
+        uniformValueFraction | uniform valueFraction    | partly | depends
+    \endtable
+
+    Example of the boundary condition specification:
+    \verbatim
+    <patchName>
+    {
+        type            uniformMixed;
+        uniformValue    constant 0.2;
+        uniformGradient constant 0.2;
+        uniformValueFraction
+        {
+            type  sine;
+            ...
+        }
+    }
+    \endverbatim
+
+Note
+    This boundary condition allows \em lazier definitions so that either
+    or both: \c uniformValue and \c uniformGradient must be defined.
+    If only of these entries is defined, the value fraction is automatically
+    treated appropriately (ie, 0 with \c uniformGradient and 1 with
+    uniformValue).
+    If both \c uniformValue and \c uniformGradient are defined,
+    the \c uniformValueFraction must also be defined.
+
+See also
+    Foam::Function1Types
+    Foam::mixedFvPatchField
+
+SourceFiles
+    uniformMixedFvPatchField.C
+
+\*---------------------------------------------------------------------------*/
+
+#ifndef Foam_uniformMixedFvPatchField_H
+#define Foam_uniformMixedFvPatchField_H
+
+#include "mixedFvPatchField.H"
+#include "PatchFunction1.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+
+/*---------------------------------------------------------------------------*\
+                  Class uniformMixedFvPatchField Declaration
+\*---------------------------------------------------------------------------*/
+
+template<class Type>
+class uniformMixedFvPatchField
+:
+    public mixedFvPatchField<Type>
+{
+    // Private Data
+
+        //- Function providing the value
+        autoPtr<PatchFunction1<Type>> refValueFunc_;
+
+        //- Function providing the gradient
+        autoPtr<PatchFunction1<Type>> refGradFunc_;
+
+        //- Function providing the value-fraction
+        autoPtr<PatchFunction1<scalar>> valueFractionFunc_;
+
+
+public:
+
+    //- Runtime type information
+    TypeName("uniformMixed");
+
+
+    // Constructors
+
+        //- Construct from patch and internal field
+        uniformMixedFvPatchField
+        (
+            const fvPatch&,
+            const DimensionedField<Type, volMesh>&
+        );
+
+        //- Construct from patch and internal field and patch field
+        uniformMixedFvPatchField
+        (
+            const fvPatch&,
+            const DimensionedField<Type, volMesh>&,
+            const Field<Type>& fld
+        );
+
+        //- Construct from patch, internal field and dictionary
+        uniformMixedFvPatchField
+        (
+            const fvPatch&,
+            const DimensionedField<Type, volMesh>&,
+            const dictionary&
+        );
+
+        //- Construct by mapping onto a new patch
+        uniformMixedFvPatchField
+        (
+            const uniformMixedFvPatchField<Type>&,
+            const fvPatch&,
+            const DimensionedField<Type, volMesh>&,
+            const fvPatchFieldMapper&
+        );
+
+        //- Construct as copy
+        uniformMixedFvPatchField
+        (
+            const uniformMixedFvPatchField<Type>&
+        );
+
+        //- Construct and return a clone
+        virtual tmp<fvPatchField<Type>> clone() const
+        {
+            return tmp<fvPatchField<Type>>
+            (
+                new uniformMixedFvPatchField<Type>(*this)
+            );
+        }
+
+        //- Construct as copy setting internal field reference
+        uniformMixedFvPatchField
+        (
+            const uniformMixedFvPatchField<Type>&,
+            const DimensionedField<Type, volMesh>&
+        );
+
+        //- Construct and return a clone setting internal field reference
+        virtual tmp<fvPatchField<Type>> clone
+        (
+            const DimensionedField<Type, volMesh>& iF
+        ) const
+        {
+            return tmp<fvPatchField<Type>>
+            (
+                new uniformMixedFvPatchField<Type>(*this, iF)
+            );
+        }
+
+
+    // Member Functions
+
+        //- Update the coefficients associated with the patch field
+        virtual void updateCoeffs();
+
+        //- Write
+        virtual void write(Ostream& os) const;
+};
+
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+} // End namespace Foam
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#ifdef NoRepository
+    #include "uniformMixedFvPatchField.C"
+#endif
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#endif
+
+// ************************************************************************* //
diff --git a/src/finiteVolume/fields/fvPatchFields/derived/uniformMixed/uniformMixedFvPatchFields.C b/src/finiteVolume/fields/fvPatchFields/derived/uniformMixed/uniformMixedFvPatchFields.C
new file mode 100644
index 0000000000000000000000000000000000000000..e307da0c6577790d688b9e365efba6b3eeb7cb13
--- /dev/null
+++ b/src/finiteVolume/fields/fvPatchFields/derived/uniformMixed/uniformMixedFvPatchFields.C
@@ -0,0 +1,39 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | www.openfoam.com
+     \\/     M anipulation  |
+-------------------------------------------------------------------------------
+    Copyright (C) 2023 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 "uniformMixedFvPatchFields.H"
+#include "addToRunTimeSelectionTable.H"
+#include "volFields.H"
+
+// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
+
+namespace Foam
+{
+    makePatchFields(uniformMixed);
+}
+
+// ************************************************************************* //
diff --git a/src/finiteVolume/fields/fvPatchFields/derived/uniformMixed/uniformMixedFvPatchFields.H b/src/finiteVolume/fields/fvPatchFields/derived/uniformMixed/uniformMixedFvPatchFields.H
new file mode 100644
index 0000000000000000000000000000000000000000..25eac16d840c66754631d85aa8925e2539c2e592
--- /dev/null
+++ b/src/finiteVolume/fields/fvPatchFields/derived/uniformMixed/uniformMixedFvPatchFields.H
@@ -0,0 +1,50 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | www.openfoam.com
+     \\/     M anipulation  |
+-------------------------------------------------------------------------------
+    Copyright (C) 2023 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/>.
+
+\*---------------------------------------------------------------------------*/
+
+#ifndef Foam_uniformMixedFvPatchFields_H
+#define Foam_uniformMixedFvPatchFields_H
+
+#include "uniformMixedFvPatchField.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+makePatchTypeFieldTypedefs(uniformMixed);
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+} // End namespace Foam
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#endif
+
+// ************************************************************************* //
diff --git a/src/meshTools/fields/pointPatchFields/uniformFixedValue/uniformFixedValuePointPatchField.C b/src/meshTools/fields/pointPatchFields/uniformFixedValue/uniformFixedValuePointPatchField.C
index 5708eab90088c4e33787327466531dc7e6f6bad4..2f291298c06021fe492dece7c9dcfb5de50e4487 100644
--- a/src/meshTools/fields/pointPatchFields/uniformFixedValue/uniformFixedValuePointPatchField.C
+++ b/src/meshTools/fields/pointPatchFields/uniformFixedValue/uniformFixedValuePointPatchField.C
@@ -152,12 +152,16 @@ void Foam::uniformFixedValuePointPatchField<Type>::autoMap
 )
 {
     fixedValuePointPatchField<Type>::autoMap(mapper);
-    refValueFunc_().autoMap(mapper);
 
-    if (refValueFunc_().constant())
+    if (refValueFunc_)
     {
-        // If mapper is not dependent on time we're ok to evaluate
-        this->evaluate();
+        refValueFunc_().autoMap(mapper);
+
+        if (refValueFunc_().constant())
+        {
+            // If mapper is not dependent on time we're ok to evaluate
+            this->evaluate();
+        }
     }
 }
 
@@ -174,7 +178,10 @@ void Foam::uniformFixedValuePointPatchField<Type>::rmap
     const uniformFixedValuePointPatchField& tiptf =
         refCast<const uniformFixedValuePointPatchField>(ptf);
 
-    refValueFunc_().rmap(tiptf.refValueFunc_(), addr);
+    if (refValueFunc_ && tiptf.refValueFunc_)
+    {
+        refValueFunc_().rmap(tiptf.refValueFunc_(), addr);
+    }
 }
 
 
@@ -186,6 +193,7 @@ void Foam::uniformFixedValuePointPatchField<Type>::updateCoeffs()
         return;
     }
     const scalar t = this->db().time().timeOutputValue();
+
     fixedValuePointPatchField<Type>::operator==(refValueFunc_->value(t));
     fixedValuePointPatchField<Type>::updateCoeffs();
 }
@@ -197,7 +205,10 @@ write(Ostream& os) const
 {
     // Note: write value
     fixedValuePointPatchField<Type>::write(os);
-    refValueFunc_->writeData(os);
+    if (refValueFunc_)
+    {
+        refValueFunc_->writeData(os);
+    }
 }
 
 
diff --git a/tutorials/compressible/rhoPimpleFoam/RAS/TJunctionAverage/0/U b/tutorials/compressible/rhoPimpleFoam/RAS/TJunctionAverage/0/U
index 6f72b041f00239179504b6ade26a5d8ab93ab932..10a031b68883b36972d14225d3f10de6445d7294 100644
--- a/tutorials/compressible/rhoPimpleFoam/RAS/TJunctionAverage/0/U
+++ b/tutorials/compressible/rhoPimpleFoam/RAS/TJunctionAverage/0/U
@@ -1,7 +1,7 @@
 /*--------------------------------*- C++ -*----------------------------------*\
 | =========                 |                                                 |
 | \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox           |
-|  \\    /   O peration     | Version:  v2212                                 |
+|  \\    /   O peration     | Version:  v2306                                 |
 |   \\  /    A nd           | Website:  www.openfoam.com                      |
 |    \\/     M anipulation  |                                                 |
 \*---------------------------------------------------------------------------*/
@@ -35,9 +35,19 @@ boundaryField
 
     outlet2
     {
-        type            inletOutlet;
-        inletValue      uniform (0 0 0);
+        // A do-it-yourself inletOutlet (for demonstration purposes)
+        // Need 'value' entry to avoid evaluation (with phi) on startup
+        type            uniformMixed;
         value           uniform (0 0 0);
+
+        uniformValue    zero;
+        uniformGradient zero;
+
+        uniformValueFraction
+        {
+            type        expression;
+            expression  "neg(phi)";
+        }
     }
 
     defaultFaces