Skip to content

Commit

Permalink
fvPatchFields/derived/fixedProfile: New BC which applies the specifie…
Browse files Browse the repository at this point in the history
…d 1D profile

This is useful when applying an experimentally obtained profile as an
inlet condition:

    Example of the boundary condition specification:
    \verbatim
    myPatch
    {
        type            fixedProfile;
        profile    csvFile;

        profileCoeffs
        {
            nHeaderLine         0;          // Number of header lines
            refColumn           0;          // Reference column index
            componentColumns    (1 2 3);    // Component column indices
            separator           ",";        // Optional (defaults to ",")
            mergeSeparators     no;         // Merge multiple separators
            fileName            "Uprofile.csv";  // name of csv data file
            outOfBounds         clamp;      // Optional out-of-bounds handling
            interpolationScheme linear;     // Optional interpolation scheme
        }
        direction        (0 1 0);
        origin           0;
    }
    \endverbatim

or a simple polynomial profile:

    Example setting a parabolic inlet profile for the PitzDaily case:
    \verbatim
    inlet
    {
        type            fixedProfile;

        profile         polynomial
        (
            ((1 0 0)        (0 0 0))
            ((-6200 0 0)    (2 0 0))
        );
        direction       (0 1 0);
        origin          0.0127;
    }
    \endverbatim

Based on code provided by Hassan Kassem:
http://www.openfoam.org/mantisbt/view.php?id=1922
  • Loading branch information
Henry Weller committed Nov 24, 2015
1 parent 9168a4a commit 528cc56
Show file tree
Hide file tree
Showing 6 changed files with 540 additions and 0 deletions.
1 change: 1 addition & 0 deletions src/finiteVolume/Make/files
Expand Up @@ -200,6 +200,7 @@ $(derivedFvPatchFields)/waveSurfacePressure/waveSurfacePressureFvPatchScalarFiel
$(derivedFvPatchFields)/interstitialInletVelocity/interstitialInletVelocityFvPatchVectorField.C $(derivedFvPatchFields)/interstitialInletVelocity/interstitialInletVelocityFvPatchVectorField.C
$(derivedFvPatchFields)/prghPressure/prghPressureFvPatchScalarField.C $(derivedFvPatchFields)/prghPressure/prghPressureFvPatchScalarField.C
$(derivedFvPatchFields)/prghTotalPressure/prghTotalPressureFvPatchScalarField.C $(derivedFvPatchFields)/prghTotalPressure/prghTotalPressureFvPatchScalarField.C
$(derivedFvPatchFields)/fixedProfile/fixedProfileFvPatchFields.C


fvsPatchFields = fields/fvsPatchFields fvsPatchFields = fields/fvsPatchFields
$(fvsPatchFields)/fvsPatchField/fvsPatchFields.C $(fvsPatchFields)/fvsPatchField/fvsPatchFields.C
Expand Down
@@ -0,0 +1,167 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2015 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 "fixedProfileFvPatchField.H"

// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //

template<class Type>
Foam::fixedProfileFvPatchField<Type>::fixedProfileFvPatchField
(
const fvPatch& p,
const DimensionedField<Type, volMesh>& iF
)
:
fixedValueFvPatchField<Type>(p, iF),
profile_(),
dir_(pTraits<vector>::zero),
origin_(0)
{}


template<class Type>
Foam::fixedProfileFvPatchField<Type>::fixedProfileFvPatchField
(
const fvPatch& p,
const DimensionedField<Type, volMesh>& iF,
const Field<Type>& fld
)
:
fixedValueFvPatchField<Type>(p, iF, fld),
profile_(),
dir_(pTraits<vector>::zero),
origin_(0)
{}


template<class Type>
Foam::fixedProfileFvPatchField<Type>::fixedProfileFvPatchField
(
const fvPatch& p,
const DimensionedField<Type, volMesh>& iF,
const dictionary& dict
)
:
fixedValueFvPatchField<Type>(p, iF),
profile_(DataEntry<Type>::New("profile", dict)),
dir_(dict.lookup("direction")),
origin_(readScalar(dict.lookup("origin")))
{
if (mag(dir_) < SMALL)
{
FatalErrorInFunction
<< "magnitude Direction must be greater than zero"
<< abort(FatalError);
}

// Ensure direction vector is normalized
dir_ /= mag(dir_);

// Evaluate profile
this->evaluate();
}


template<class Type>
Foam::fixedProfileFvPatchField<Type>::fixedProfileFvPatchField
(
const fixedProfileFvPatchField<Type>& ptf,
const fvPatch& p,
const DimensionedField<Type, volMesh>& iF,
const fvPatchFieldMapper& mapper
)
:
fixedValueFvPatchField<Type>(p, iF), // Don't map
profile_(ptf.profile_, false),
dir_(ptf.dir_),
origin_(ptf.origin_)
{
// Evaluate profile since value not mapped
this->evaluate();
}


template<class Type>
Foam::fixedProfileFvPatchField<Type>::fixedProfileFvPatchField
(
const fixedProfileFvPatchField<Type>& ptf
)
:
fixedValueFvPatchField<Type>(ptf),
profile_(ptf.profile_, false),
dir_(ptf.dir_),
origin_(ptf.origin_)
{}


template<class Type>
Foam::fixedProfileFvPatchField<Type>::fixedProfileFvPatchField
(
const fixedProfileFvPatchField<Type>& ptf,
const DimensionedField<Type, volMesh>& iF
)
:
fixedValueFvPatchField<Type>(ptf, iF),
profile_(ptf.profile_, false),
dir_(ptf.dir_),
origin_(ptf.origin_)
{
// Evaluate the profile if defined
if (ptf.profile_.valid())
{
this->evaluate();
}
}


// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //

template<class Type>
void Foam::fixedProfileFvPatchField<Type>::updateCoeffs()
{
if (this->updated())
{
return;
}

const scalarField dirCmpt((dir_ & this->patch().Cf()) - origin_);
fvPatchField<Type>::operator==(profile_->value(dirCmpt));

fixedValueFvPatchField<Type>::updateCoeffs();
}


template<class Type>
void Foam::fixedProfileFvPatchField<Type>::write(Ostream& os) const
{
fvPatchField<Type>::write(os);
profile_->writeData(os);
os.writeKeyword("direction") << dir_ << token::END_STATEMENT << nl;
os.writeKeyword("origin") << origin_ << token::END_STATEMENT << nl;
this->writeEntry("value", os);
}


// ************************************************************************* //

0 comments on commit 528cc56

Please sign in to comment.