prghUniformDensityHydrostaticPressureFvPatchScalarField: New boundary condition for p_rgh

Sets the boundary values of p_rgh corresponding to a constant density hydrostatic
pressure distribution.

Description
    This boundary condition provides a hydrostatic pressure condition for p_rgh,
    calculated as:

        \f[
            p_{rgh} = p_{ref} - (\rho - \rho_0) g (h - h_{ref})
        \f]

    where
    \vartable
        p_{rgh}    | Pseudo hydrostatic pressure [Pa]
        p_{ref}    | Static pressure at hRef [Pa]
        h          | Height in the opposite direction to gravity
        h_{ref}    | Reference height in the opposite direction to gravity
        \rho       | Density field
        \rho_{ref} | Uniform reference density at boundary
        g          | Acceleration due to gravity [m/s^2]
    \endtable

Usage
    \table
        Property     | Description               | Required    | Default value
        pRef         | Reference static pressure | yes         |
        rhoRef       | Reference density         | yes         |
        rho          | Density field name        | no          | rho
    \endtable

    Example of the boundary condition specification:
    \verbatim
    <patchName>
    {
        type            prghUniformDensityHydrostaticPressure;
        rhoRef          1000;
        p               0;
        value           uniform 0; // optional initial value
    }
    \endverbatim
This commit is contained in:
Henry Weller
2018-03-08 21:59:29 +00:00
parent ba84383e26
commit 52a2ba84c2
4 changed files with 392 additions and 6 deletions

View File

@ -192,6 +192,7 @@ $(derivedFvPatchFields)/translatingWallVelocity/translatingWallVelocityFvPatchVe
$(derivedFvPatchFields)/turbulentInlet/turbulentInletFvPatchFields.C
$(derivedFvPatchFields)/turbulentIntensityKineticEnergyInlet/turbulentIntensityKineticEnergyInletFvPatchScalarField.C
$(derivedFvPatchFields)/uniformDensityHydrostaticPressure/uniformDensityHydrostaticPressureFvPatchScalarField.C
$(derivedFvPatchFields)/prghUniformDensityHydrostaticPressure/prghUniformDensityHydrostaticPressureFvPatchScalarField.C
$(derivedFvPatchFields)/uniformFixedGradient/uniformFixedGradientFvPatchFields.C
$(derivedFvPatchFields)/uniformFixedValue/uniformFixedValueFvPatchFields.C
$(derivedFvPatchFields)/uniformInletOutlet/uniformInletOutletFvPatchFields.C

View File

@ -0,0 +1,179 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / 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 "prghUniformDensityHydrostaticPressureFvPatchScalarField.H"
#include "addToRunTimeSelectionTable.H"
#include "fvPatchFieldMapper.H"
#include "volFields.H"
#include "uniformDimensionedFields.H"
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::prghUniformDensityHydrostaticPressureFvPatchScalarField::
prghUniformDensityHydrostaticPressureFvPatchScalarField
(
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF
)
:
fixedValueFvPatchScalarField(p, iF),
pRef_(0),
rhoRef_(0),
rhoName_("rho")
{}
Foam::prghUniformDensityHydrostaticPressureFvPatchScalarField::
prghUniformDensityHydrostaticPressureFvPatchScalarField
(
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF,
const dictionary& dict
)
:
fixedValueFvPatchScalarField(p, iF, dict, false),
pRef_(readScalar(dict.lookup("pRef"))),
rhoRef_(readScalar(dict.lookup("rhoRef"))),
rhoName_(dict.lookupOrDefault<word>("rho", "rho"))
{
if (dict.found("value"))
{
fvPatchScalarField::operator=
(
scalarField("value", dict, p.size())
);
}
else
{
fvPatchField<scalar>::operator=(pRef_);
}
}
Foam::prghUniformDensityHydrostaticPressureFvPatchScalarField::
prghUniformDensityHydrostaticPressureFvPatchScalarField
(
const prghUniformDensityHydrostaticPressureFvPatchScalarField& ptf,
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF,
const fvPatchFieldMapper& mapper
)
:
fixedValueFvPatchScalarField(ptf, p, iF, mapper),
pRef_(ptf.pRef_),
rhoRef_(ptf.rhoRef_),
rhoName_(ptf.rhoName_)
{}
Foam::prghUniformDensityHydrostaticPressureFvPatchScalarField::
prghUniformDensityHydrostaticPressureFvPatchScalarField
(
const prghUniformDensityHydrostaticPressureFvPatchScalarField& ptf
)
:
fixedValueFvPatchScalarField(ptf),
pRef_(ptf.pRef_),
rhoRef_(ptf.rhoRef_),
rhoName_(ptf.rhoName_)
{}
Foam::prghUniformDensityHydrostaticPressureFvPatchScalarField::
prghUniformDensityHydrostaticPressureFvPatchScalarField
(
const prghUniformDensityHydrostaticPressureFvPatchScalarField& ptf,
const DimensionedField<scalar, volMesh>& iF
)
:
fixedValueFvPatchScalarField(ptf, iF),
pRef_(ptf.pRef_),
rhoRef_(ptf.rhoRef_),
rhoName_(ptf.rhoName_)
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
void Foam::prghUniformDensityHydrostaticPressureFvPatchScalarField::
updateCoeffs()
{
if (updated())
{
return;
}
const scalarField& rhop = patch().lookupPatchField<volScalarField, scalar>
(
rhoName_
);
const uniformDimensionedVectorField& g =
db().lookupObject<uniformDimensionedVectorField>("g");
const uniformDimensionedScalarField& hRef =
db().lookupObject<uniformDimensionedScalarField>("hRef");
dimensionedScalar ghRef
(
mag(g.value()) > small
? g & (cmptMag(g.value())/mag(g.value()))*hRef
: dimensionedScalar("ghRef", g.dimensions()*dimLength, 0)
);
operator==
(
pRef_ - (rhop - rhoRef_)*((g.value() & patch().Cf()) - ghRef.value())
);
fixedValueFvPatchScalarField::updateCoeffs();
}
void Foam::prghUniformDensityHydrostaticPressureFvPatchScalarField::write
(
Ostream& os
) const
{
fvPatchScalarField::write(os);
os.writeKeyword("pRef") << pRef_ << token::END_STATEMENT << nl;
os.writeKeyword("rhoRef") << rhoRef_ << token::END_STATEMENT << nl;
writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
writeEntry("value", os);
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
makePatchTypeField
(
fvPatchScalarField,
prghUniformDensityHydrostaticPressureFvPatchScalarField
);
}
// ************************************************************************* //

View File

@ -0,0 +1,206 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / 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::prghUniformDensityHydrostaticPressureFvPatchScalarField
Group
grpGenericBoundaryConditions
Description
This boundary condition provides a hydrostatic pressure condition for p_rgh,
calculated as:
\f[
p_{rgh} = p_{ref} - (\rho - \rho_0) g (h - h_{ref})
\f]
where
\vartable
p_{rgh} | Pseudo hydrostatic pressure [Pa]
p_{ref} | Static pressure at hRef [Pa]
h | Height in the opposite direction to gravity
h_{ref} | Reference height in the opposite direction to gravity
\rho | Density field
\rho_{ref} | Uniform reference density at boundary
g | Acceleration due to gravity [m/s^2]
\endtable
Usage
\table
Property | Description | Required | Default value
pRef | Reference static pressure | yes |
rhoRef | Reference density | yes |
rho | Density field name | no | rho
\endtable
Example of the boundary condition specification:
\verbatim
<patchName>
{
type prghUniformDensityHydrostaticPressure;
rhoRef 1000;
p 0;
value uniform 0; // optional initial value
}
\endverbatim
See also
Foam::fixedValueFvPatchScalarField
SourceFiles
prghUniformDensityHydrostaticPressureFvPatchScalarField.C
\*---------------------------------------------------------------------------*/
#ifndef prghUniformDensityHydrostaticPressureFvPatchScalarField_H
#define prghUniformDensityHydrostaticPressureFvPatchScalarField_H
#include "fixedValueFvPatchFields.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
/*---------------------------------------------------------------------------*\
Class prghUniformDensityHydrostaticPressureFvPatchScalarField Declaration
\*---------------------------------------------------------------------------*/
class prghUniformDensityHydrostaticPressureFvPatchScalarField
:
public fixedValueFvPatchScalarField
{
protected:
// Protected data
//- Reference Static pressure
scalar pRef_;
//- Reference density
scalar rhoRef_;
//- Name of phase-fraction field
word rhoName_;
public:
//- Runtime type information
TypeName("prghUniformDensityHydrostaticPressure");
// Constructors
//- Construct from patch and internal field
prghUniformDensityHydrostaticPressureFvPatchScalarField
(
const fvPatch&,
const DimensionedField<scalar, volMesh>&
);
//- Construct from patch, internal field and dictionary
prghUniformDensityHydrostaticPressureFvPatchScalarField
(
const fvPatch&,
const DimensionedField<scalar, volMesh>&,
const dictionary&
);
//- Construct by mapping given
// prghUniformDensityHydrostaticPressureFvPatchScalarField
// onto a new patch
prghUniformDensityHydrostaticPressureFvPatchScalarField
(
const prghUniformDensityHydrostaticPressureFvPatchScalarField&,
const fvPatch&,
const DimensionedField<scalar, volMesh>&,
const fvPatchFieldMapper&
);
//- Construct as copy
prghUniformDensityHydrostaticPressureFvPatchScalarField
(
const prghUniformDensityHydrostaticPressureFvPatchScalarField&
);
//- Construct and return a clone
virtual tmp<fvPatchScalarField> clone() const
{
return tmp<fvPatchScalarField >
(
new prghUniformDensityHydrostaticPressureFvPatchScalarField
(
*this
)
);
}
//- Construct as copy setting internal field reference
prghUniformDensityHydrostaticPressureFvPatchScalarField
(
const prghUniformDensityHydrostaticPressureFvPatchScalarField&,
const DimensionedField<scalar, volMesh>&
);
//- Construct and return a clone setting internal field reference
virtual tmp<fvPatchScalarField> clone
(
const DimensionedField<scalar, volMesh>& iF
) const
{
return tmp<fvPatchScalarField>
(
new prghUniformDensityHydrostaticPressureFvPatchScalarField
(
*this,
iF
)
);
}
// Member functions
// Evaluation functions
//- Update the coefficients associated with the patch field
virtual void updateCoeffs();
//- Write
virtual void write(Ostream&) const;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //

View File

@ -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
@ -37,11 +37,11 @@ Description
where
\vartable
p_{hyd} | hyrostatic pressure [Pa]
p_{ref} | reference pressure [Pa]
x_{ref} | reference point in Cartesian co-ordinates
\rho | density (assumed uniform)
g | acceleration due to gravity [m/s2]
p_{hyd} | Hyrostatic pressure [Pa]
p_{ref} | Reference pressure [Pa]
x_{ref} | Reference point in Cartesian co-ordinates
\rho | Density (assumed uniform)
g | Acceleration due to gravity [m/s2]
\endtable
Usage