From 8165c462cd7fb184ff0d98d559ac1b8fc23d7d71 Mon Sep 17 00:00:00 2001
From: Henry Weller <http://openfoam.org>
Date: Fri, 23 Feb 2018 12:23:06 +0000
Subject: [PATCH] ENH: freestreamPressure, freestreamVelocity: New blended
 boundary conditions for the freestream

These BCs blend between typical inflow and outflow conditions based on the
velocity orientation.

airFoil2D tutorial updated to demonstrate these new BCs.
---
 src/finiteVolume/Make/files                   |   1 +
 .../freestream/freestreamFvPatchField.H       |   4 +-
 .../freestreamPressureFvPatchScalarField.C    |  95 ++++-----
 .../freestreamPressureFvPatchScalarField.H    |  40 ++--
 .../freestreamVelocityFvPatchVectorField.C    | 136 +++++++++++++
 .../freestreamVelocityFvPatchVectorField.H    | 184 ++++++++++++++++++
 .../incompressible/simpleFoam/airFoil2D/0/U   |   8 +-
 .../incompressible/simpleFoam/airFoil2D/0/p   |   2 +
 .../simpleFoam/airFoil2D/system/controlDict   |   2 +-
 .../simpleFoam/airFoil2D/system/fvSolution    |   2 -
 10 files changed, 388 insertions(+), 86 deletions(-)
 create mode 100644 src/finiteVolume/fields/fvPatchFields/derived/freestreamVelocity/freestreamVelocityFvPatchVectorField.C
 create mode 100644 src/finiteVolume/fields/fvPatchFields/derived/freestreamVelocity/freestreamVelocityFvPatchVectorField.H

diff --git a/src/finiteVolume/Make/files b/src/finiteVolume/Make/files
index f64ba15fe81..74dd85ebe42 100644
--- a/src/finiteVolume/Make/files
+++ b/src/finiteVolume/Make/files
@@ -154,6 +154,7 @@ $(derivedFvPatchFields)/flowRateInletVelocity/flowRateInletVelocityFvPatchVector
 $(derivedFvPatchFields)/flowRateOutletVelocity/flowRateOutletVelocityFvPatchVectorField.C
 $(derivedFvPatchFields)/fluxCorrectedVelocity/fluxCorrectedVelocityFvPatchVectorField.C
 $(derivedFvPatchFields)/freestream/freestreamFvPatchFields.C
+$(derivedFvPatchFields)/freestreamVelocity/freestreamVelocityFvPatchVectorField.C
 $(derivedFvPatchFields)/freestreamPressure/freestreamPressureFvPatchScalarField.C
 $(derivedFvPatchFields)/inletOutlet/inletOutletFvPatchFields.C
 $(derivedFvPatchFields)/inletOutletTotalTemperature/inletOutletTotalTemperatureFvPatchScalarField.C
diff --git a/src/finiteVolume/fields/fvPatchFields/derived/freestream/freestreamFvPatchField.H b/src/finiteVolume/fields/fvPatchFields/derived/freestream/freestreamFvPatchField.H
index 258c38ceefd..15757c90cf1 100644
--- a/src/finiteVolume/fields/fvPatchFields/derived/freestream/freestreamFvPatchField.H
+++ b/src/finiteVolume/fields/fvPatchFields/derived/freestream/freestreamFvPatchField.H
@@ -2,7 +2,7 @@
   =========                 |
   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
    \\    /   O peration     |
-    \\  /    A nd           | Copyright (C) 2011-2016 OpenFOAM Foundation
+    \\  /    A nd           | Copyright (C) 2011-2018 OpenFOAM Foundation
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -45,7 +45,7 @@ Usage
     <patchName>
     {
         type            freestream;
-        phi             phi;
+        freestreamValue uniform (300 0 0);
     }
     \endverbatim
 
diff --git a/src/finiteVolume/fields/fvPatchFields/derived/freestreamPressure/freestreamPressureFvPatchScalarField.C b/src/finiteVolume/fields/fvPatchFields/derived/freestreamPressure/freestreamPressureFvPatchScalarField.C
index 1986bc40d68..551b4d55e70 100644
--- a/src/finiteVolume/fields/fvPatchFields/derived/freestreamPressure/freestreamPressureFvPatchScalarField.C
+++ b/src/finiteVolume/fields/fvPatchFields/derived/freestreamPressure/freestreamPressureFvPatchScalarField.C
@@ -2,7 +2,7 @@
   =========                 |
   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
    \\    /   O peration     |
-    \\  /    A nd           | Copyright (C) 2011-2016 OpenFOAM Foundation
+    \\  /    A nd           | Copyright (C) 2011-2018 OpenFOAM Foundation
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -24,10 +24,7 @@ License
 \*---------------------------------------------------------------------------*/
 
 #include "freestreamPressureFvPatchScalarField.H"
-#include "freestreamFvPatchFields.H"
-#include "fvPatchFieldMapper.H"
 #include "volFields.H"
-#include "surfaceFields.H"
 #include "addToRunTimeSelectionTable.H"
 
 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
@@ -39,10 +36,8 @@ freestreamPressureFvPatchScalarField
     const DimensionedField<scalar, volMesh>& iF
 )
 :
-    zeroGradientFvPatchScalarField(p, iF),
-    UName_("U"),
-    phiName_("phi"),
-    rhoName_("rho")
+    mixedFvPatchScalarField(p, iF),
+    UName_("U")
 {}
 
 
@@ -54,11 +49,26 @@ freestreamPressureFvPatchScalarField
     const dictionary& dict
 )
 :
-    zeroGradientFvPatchScalarField(p, iF, dict),
-    UName_(dict.lookupOrDefault<word>("U", "U")),
-    phiName_(dict.lookupOrDefault<word>("phi", "phi")),
-    rhoName_(dict.lookupOrDefault<word>("rho", "rho"))
-{}
+    mixedFvPatchScalarField(p, iF),
+    UName_(dict.lookupOrDefault<word>("U", "U"))
+{
+    freestreamValue() = scalarField("freestreamValue", dict, p.size());
+
+    if (dict.found("value"))
+    {
+        fvPatchScalarField::operator=
+        (
+            scalarField("value", dict, p.size())
+        );
+    }
+    else
+    {
+        fvPatchScalarField::operator=(freestreamValue());
+    }
+
+    refGrad() = Zero;
+    valueFraction() = 0;
+}
 
 
 Foam::freestreamPressureFvPatchScalarField::
@@ -70,10 +80,8 @@ freestreamPressureFvPatchScalarField
     const fvPatchFieldMapper& mapper
 )
 :
-    zeroGradientFvPatchScalarField(ptf, p, iF, mapper),
-    UName_(ptf.UName_),
-    phiName_(ptf.phiName_),
-    rhoName_(ptf.rhoName_)
+    mixedFvPatchScalarField(ptf, p, iF, mapper),
+    UName_(ptf.UName_)
 {}
 
 
@@ -83,10 +91,8 @@ freestreamPressureFvPatchScalarField
     const freestreamPressureFvPatchScalarField& wbppsf
 )
 :
-    zeroGradientFvPatchScalarField(wbppsf),
-    UName_(wbppsf.UName_),
-    phiName_(wbppsf.phiName_),
-    rhoName_(wbppsf.rhoName_)
+    mixedFvPatchScalarField(wbppsf),
+    UName_(wbppsf.UName_)
 {}
 
 
@@ -97,10 +103,8 @@ freestreamPressureFvPatchScalarField
     const DimensionedField<scalar, volMesh>& iF
 )
 :
-    zeroGradientFvPatchScalarField(wbppsf, iF),
-    UName_(wbppsf.UName_),
-    phiName_(wbppsf.phiName_),
-    rhoName_(wbppsf.rhoName_)
+    mixedFvPatchScalarField(wbppsf, iF),
+    UName_(wbppsf.UName_)
 {}
 
 
@@ -113,43 +117,15 @@ void Foam::freestreamPressureFvPatchScalarField::updateCoeffs()
         return;
     }
 
-    const freestreamFvPatchVectorField& Up =
-        refCast<const freestreamFvPatchVectorField>
+    const Field<vector>& Up =
+        patch().template lookupPatchField<volVectorField, vector>
         (
-            patch().lookupPatchField<volVectorField, vector>(UName_)
+            UName_
         );
 
-    const surfaceScalarField& phi =
-        db().lookupObject<surfaceScalarField>(phiName_);
-
-    fvsPatchField<scalar>& phip =
-        const_cast<fvsPatchField<scalar>&>
-        (
-            patch().patchField<surfaceScalarField, scalar>(phi)
-        );
-
-    if (phi.dimensions() == dimVelocity*dimArea)
-    {
-        phip = patch().Sf() & Up.freestreamValue();
-    }
-    else if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
-    {
-        const fvPatchField<scalar>& rhop =
-            patch().lookupPatchField<volScalarField, scalar>(rhoName_);
-
-        phip = rhop*(patch().Sf() & Up.freestreamValue());
-    }
-    else
-    {
-        FatalErrorInFunction
-            << "dimensions of phi are not correct"
-            << "\n    on patch " << this->patch().name()
-            << " of field " << this->internalField().name()
-            << " in file " << this->internalField().objectPath()
-            << exit(FatalError);
-    }
+    valueFraction() = 0.5 + 0.5*(Up & patch().nf())/mag(Up);
 
-    zeroGradientFvPatchScalarField::updateCoeffs();
+    mixedFvPatchField<scalar>::updateCoeffs();
 }
 
 
@@ -157,8 +133,7 @@ void Foam::freestreamPressureFvPatchScalarField::write(Ostream& os) const
 {
     fvPatchScalarField::write(os);
     os.writeEntryIfDifferent<word>("U", "U", UName_);
-    os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
-    os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
+    freestreamValue().writeEntry("freestreamValue", os);
     writeEntry("value", os);
 }
 
diff --git a/src/finiteVolume/fields/fvPatchFields/derived/freestreamPressure/freestreamPressureFvPatchScalarField.H b/src/finiteVolume/fields/fvPatchFields/derived/freestreamPressure/freestreamPressureFvPatchScalarField.H
index be263b1af2a..8486e42805b 100644
--- a/src/finiteVolume/fields/fvPatchFields/derived/freestreamPressure/freestreamPressureFvPatchScalarField.H
+++ b/src/finiteVolume/fields/fvPatchFields/derived/freestreamPressure/freestreamPressureFvPatchScalarField.H
@@ -29,15 +29,16 @@ Group
 
 Description
     This boundary condition provides a free-stream condition for pressure.
-    It is a zero-gradient condition that constrains the flux across the patch
-    based on the free-stream velocity.
+
+    It is an outlet-inlet condition that uses the velocity orientation to
+    continuously blend between zero gradient for normal inlet and fixed value
+    for normal outlet flow.
 
 Usage
     \table
-        Property     | Description             | Required    | Default value
-        U            | velocity field name     | no          | U
-        phi          | flux field name         | no          | phi
-        rho          | density field name      | no          | none
+        Property        | Description             | Required    | Default value
+        U               | velocity field name     | no          | U
+        freestreamValue | freestream pressure     | yes         |
     \endtable
 
     Example of the boundary condition specification:
@@ -45,14 +46,15 @@ Usage
     <patchName>
     {
         type            freestreamPressure;
+        freestreamValue uniform 1e5;
     }
     \endverbatim
 
 Note
-    This condition is designed to operate with a freestream velocity condition
+    This condition is designed to operate with a freestreamVelocity condition
 
 See also
-    Foam::zeroGradientFvPatchField
+    Foam::mixedFvPatchField
     Foam::freestreamFvPatchField
 
 SourceFiles
@@ -64,7 +66,7 @@ SourceFiles
 #define freestreamPressureFvPatchScalarField_H
 
 #include "fvPatchFields.H"
-#include "zeroGradientFvPatchFields.H"
+#include "mixedFvPatchFields.H"
 
 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
 
@@ -77,20 +79,13 @@ namespace Foam
 
 class freestreamPressureFvPatchScalarField
 :
-    public zeroGradientFvPatchScalarField
+    public mixedFvPatchScalarField
 {
     // Private data
 
         //- Name of the velocity field
         word UName_;
 
-        //- Name of the flux transporting the field
-        word phiName_;
-
-        //- Name of the density field used to normalise the mass flux
-        //- if necessary
-        word rhoName_;
-
 
 public:
 
@@ -162,6 +157,17 @@ public:
 
     // Member functions
 
+            const scalarField& freestreamValue() const
+            {
+                return refValue();
+            }
+
+            scalarField& freestreamValue()
+            {
+                return refValue();
+            }
+
+
         // Evaluation functions
 
             //- Update the coefficients associated with the patch field
diff --git a/src/finiteVolume/fields/fvPatchFields/derived/freestreamVelocity/freestreamVelocityFvPatchVectorField.C b/src/finiteVolume/fields/fvPatchFields/derived/freestreamVelocity/freestreamVelocityFvPatchVectorField.C
new file mode 100644
index 00000000000..9990cb0c8e7
--- /dev/null
+++ b/src/finiteVolume/fields/fvPatchFields/derived/freestreamVelocity/freestreamVelocityFvPatchVectorField.C
@@ -0,0 +1,136 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | Copyright (C) 2018 OpenFOAM Foundation
+     \\/     M anipulation  |
+-------------------------------------------------------------------------------
+License
+    This file is part of OpenFOAM.
+
+    OpenFOAM is free software: you can redistribute it and/or modify it
+    under the terms of the GNU General Public License as published by
+    the Free Software Foundation, either version 3 of the License, or
+    (at your option) any later version.
+
+    OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
+    ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
+    FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
+    for more details.
+
+    You should have received a copy of the GNU General Public License
+    along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
+
+\*---------------------------------------------------------------------------*/
+
+#include "freestreamVelocityFvPatchVectorField.H"
+#include "addToRunTimeSelectionTable.H"
+
+// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
+
+Foam::freestreamVelocityFvPatchVectorField::freestreamVelocityFvPatchVectorField
+(
+    const fvPatch& p,
+    const DimensionedField<vector, volMesh>& iF
+)
+:
+    mixedFvPatchVectorField(p, iF)
+{}
+
+
+Foam::freestreamVelocityFvPatchVectorField::freestreamVelocityFvPatchVectorField
+(
+    const fvPatch& p,
+    const DimensionedField<vector, volMesh>& iF,
+    const dictionary& dict
+)
+:
+    mixedFvPatchVectorField(p, iF)
+{
+    freestreamValue() = vectorField("freestreamValue", dict, p.size());
+
+    if (dict.found("value"))
+    {
+        fvPatchVectorField::operator=
+        (
+            vectorField("value", dict, p.size())
+        );
+    }
+    else
+    {
+        fvPatchVectorField::operator=(freestreamValue());
+    }
+
+    refGrad() = Zero;
+    valueFraction() = 1;
+}
+
+
+Foam::freestreamVelocityFvPatchVectorField::freestreamVelocityFvPatchVectorField
+(
+    const freestreamVelocityFvPatchVectorField& ptf,
+    const fvPatch& p,
+    const DimensionedField<vector, volMesh>& iF,
+    const fvPatchFieldMapper& mapper
+)
+:
+    mixedFvPatchVectorField(ptf, p, iF, mapper)
+{}
+
+
+Foam::freestreamVelocityFvPatchVectorField::freestreamVelocityFvPatchVectorField
+(
+    const freestreamVelocityFvPatchVectorField& wbppsf
+)
+:
+    mixedFvPatchVectorField(wbppsf)
+{}
+
+
+Foam::freestreamVelocityFvPatchVectorField::freestreamVelocityFvPatchVectorField
+(
+    const freestreamVelocityFvPatchVectorField& wbppsf,
+    const DimensionedField<vector, volMesh>& iF
+)
+:
+    mixedFvPatchVectorField(wbppsf, iF)
+{}
+
+
+// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
+
+void Foam::freestreamVelocityFvPatchVectorField::updateCoeffs()
+{
+    if (updated())
+    {
+        return;
+    }
+
+    const Field<vector>& Up = *this;
+
+    valueFraction() = 0.5 - 0.5*(Up & patch().nf())/mag(Up);
+
+    mixedFvPatchField<vector>::updateCoeffs();
+}
+
+
+void Foam::freestreamVelocityFvPatchVectorField::write(Ostream& os) const
+{
+    fvPatchVectorField::write(os);
+    freestreamValue().writeEntry("freestreamValue", os);
+    writeEntry("value", os);
+}
+
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+    makePatchTypeField
+    (
+        fvPatchVectorField,
+        freestreamVelocityFvPatchVectorField
+    );
+}
+
+// ************************************************************************* //
diff --git a/src/finiteVolume/fields/fvPatchFields/derived/freestreamVelocity/freestreamVelocityFvPatchVectorField.H b/src/finiteVolume/fields/fvPatchFields/derived/freestreamVelocity/freestreamVelocityFvPatchVectorField.H
new file mode 100644
index 00000000000..9415256e495
--- /dev/null
+++ b/src/finiteVolume/fields/fvPatchFields/derived/freestreamVelocity/freestreamVelocityFvPatchVectorField.H
@@ -0,0 +1,184 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | Copyright (C) 2018 OpenFOAM Foundation
+     \\/     M anipulation  |
+-------------------------------------------------------------------------------
+License
+    This file is part of OpenFOAM.
+
+    OpenFOAM is free software: you can redistribute it and/or modify it
+    under the terms of the GNU General Public License as published by
+    the Free Software Foundation, either version 3 of the License, or
+    (at your option) any later version.
+
+    OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
+    ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
+    FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
+    for more details.
+
+    You should have received a copy of the GNU General Public License
+    along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
+
+Class
+    Foam::freestreamVelocityFvPatchVectorField
+
+Group
+    grpInletBoundaryConditions grpOutletBoundaryConditions
+
+Description
+    This boundary condition provides a free-stream condition for velocity.
+
+    It is an inlet-outlet condition that uses the velocity orientation to
+    continuously blend between fixed value for normal inlet and zero gradient
+    for normal outlet flow.
+
+Usage
+    \table
+        Property        | Description             | Required    | Default value
+        freestreamValue | freestream velocity     | yes         |
+    \endtable
+
+    Example of the boundary condition specification:
+    \verbatim
+    <patchName>
+    {
+        type            freestreamVelocity;
+        freestreamValue uniform (300 0 0);
+    }
+    \endverbatim
+
+Note
+    This condition is designed to operate with the freestreamPressure condition
+
+See also
+    Foam::mixedFvPatchField
+    Foam::freestreamFvPatchField
+
+SourceFiles
+    freestreamVelocityFvPatchVectorField.C
+
+\*---------------------------------------------------------------------------*/
+
+#ifndef freestreamVelocityFvPatchVectorField_H
+#define freestreamVelocityFvPatchVectorField_H
+
+#include "fvPatchFields.H"
+#include "mixedFvPatchFields.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+
+/*---------------------------------------------------------------------------*\
+             Class freestreamVelocityFvPatchVectorField Declaration
+\*---------------------------------------------------------------------------*/
+
+class freestreamVelocityFvPatchVectorField
+:
+    public mixedFvPatchVectorField
+{
+
+public:
+
+    //- Runtime type information
+    TypeName("freestreamVelocity");
+
+
+    // Constructors
+
+        //- Construct from patch and internal field
+        freestreamVelocityFvPatchVectorField
+        (
+            const fvPatch&,
+            const DimensionedField<vector, volMesh>&
+        );
+
+        //- Construct from patch, internal field and dictionary
+        freestreamVelocityFvPatchVectorField
+        (
+            const fvPatch&,
+            const DimensionedField<vector, volMesh>&,
+            const dictionary&
+        );
+
+        //- Construct by mapping given freestreamVelocityFvPatchVectorField onto
+        //  a new patch
+        freestreamVelocityFvPatchVectorField
+        (
+            const freestreamVelocityFvPatchVectorField&,
+            const fvPatch&,
+            const DimensionedField<vector, volMesh>&,
+            const fvPatchFieldMapper&
+        );
+
+        //- Construct as copy
+        freestreamVelocityFvPatchVectorField
+        (
+            const freestreamVelocityFvPatchVectorField&
+        );
+
+        //- Construct and return a clone
+        virtual tmp<fvPatchVectorField> clone() const
+        {
+            return tmp<fvPatchVectorField>
+            (
+                new freestreamVelocityFvPatchVectorField(*this)
+            );
+        }
+
+        //- Construct as copy setting internal field reference
+        freestreamVelocityFvPatchVectorField
+        (
+            const freestreamVelocityFvPatchVectorField&,
+            const DimensionedField<vector, volMesh>&
+        );
+
+        //- Construct and return a clone setting internal field reference
+        virtual tmp<fvPatchVectorField> clone
+        (
+            const DimensionedField<vector, volMesh>& iF
+        ) const
+        {
+            return tmp<fvPatchVectorField>
+            (
+                new freestreamVelocityFvPatchVectorField(*this, iF)
+            );
+        }
+
+
+    // Member functions
+
+            const vectorField& freestreamValue() const
+            {
+                return refValue();
+            }
+
+            vectorField& freestreamValue()
+            {
+                return refValue();
+            }
+
+
+        // Evaluation functions
+
+            //- Update the coefficients associated with the patch field
+            virtual void updateCoeffs();
+
+
+        //- Write
+        virtual void write(Ostream&) const;
+};
+
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+} // End namespace Foam
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#endif
+
+// ************************************************************************* //
diff --git a/tutorials/incompressible/simpleFoam/airFoil2D/0/U b/tutorials/incompressible/simpleFoam/airFoil2D/0/U
index ba119821d39..5121c5c3517 100644
--- a/tutorials/incompressible/simpleFoam/airFoil2D/0/U
+++ b/tutorials/incompressible/simpleFoam/airFoil2D/0/U
@@ -22,14 +22,14 @@ boundaryField
 {
     inlet
     {
-        type            freestream;
-        freestreamValue uniform (25.75 3.62 0);
+        type            freestreamVelocity;
+        freestreamValue $internalField;
     }
 
     outlet
     {
-        type            freestream;
-        freestreamValue uniform (25.75 3.62 0);
+        type            freestreamVelocity;
+        freestreamValue $internalField;
     }
 
     walls
diff --git a/tutorials/incompressible/simpleFoam/airFoil2D/0/p b/tutorials/incompressible/simpleFoam/airFoil2D/0/p
index 4494ee33c1f..75c231def0b 100644
--- a/tutorials/incompressible/simpleFoam/airFoil2D/0/p
+++ b/tutorials/incompressible/simpleFoam/airFoil2D/0/p
@@ -23,11 +23,13 @@ boundaryField
     inlet
     {
         type            freestreamPressure;
+        freestreamValue $internalField;
     }
 
     outlet
     {
         type            freestreamPressure;
+        freestreamValue $internalField;
     }
 
     walls
diff --git a/tutorials/incompressible/simpleFoam/airFoil2D/system/controlDict b/tutorials/incompressible/simpleFoam/airFoil2D/system/controlDict
index 3e4f6c5e380..6d0c7a03523 100644
--- a/tutorials/incompressible/simpleFoam/airFoil2D/system/controlDict
+++ b/tutorials/incompressible/simpleFoam/airFoil2D/system/controlDict
@@ -17,7 +17,7 @@ FoamFile
 
 application     simpleFoam;
 
-startFrom       latestTime;
+startFrom       startTime;
 
 startTime       0;
 
diff --git a/tutorials/incompressible/simpleFoam/airFoil2D/system/fvSolution b/tutorials/incompressible/simpleFoam/airFoil2D/system/fvSolution
index 1d22f18c896..0dc5cf88e4b 100644
--- a/tutorials/incompressible/simpleFoam/airFoil2D/system/fvSolution
+++ b/tutorials/incompressible/simpleFoam/airFoil2D/system/fvSolution
@@ -47,8 +47,6 @@ solvers
 SIMPLE
 {
     nNonOrthogonalCorrectors 0;
-    pRefCell        0;
-    pRefValue       0;
 
     residualControl
     {
-- 
GitLab