Skip to content
Snippets Groups Projects
Commit dada835d authored by sergio's avatar sergio
Browse files

ENH: adding swirlFanVelocity BC for jumpCyclic types

It creates a swirling velocity across a baffle using jumpCyclic for vectors
parent 4e47be34
Branches
Tags
No related merge requests found
......@@ -3,7 +3,7 @@
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
\\/ M anipulation |
\\/ M anipulation | Copyright (C) 2018 OpenCFD Ltd.
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
......@@ -159,6 +159,18 @@ void jumpCyclicFvPatchField<scalar>::updateInterfaceMatrix
) const;
template<>
void jumpCyclicFvPatchField<vector>::updateInterfaceMatrix
(
scalarField& result,
const bool add,
const scalarField& psiInternal,
const scalarField& coeffs,
const direction cmpt,
const Pstream::commsTypes commsType
) const;
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
......
......@@ -3,7 +3,7 @@
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
\\/ M anipulation |
\\/ M anipulation | Copyright (C) 2018 OpenCFD Ltd.
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
......@@ -83,4 +83,54 @@ void Foam::jumpCyclicFvPatchField<Foam::scalar>::updateInterfaceMatrix
}
template<>
void Foam::jumpCyclicFvPatchField<Foam::vector>::updateInterfaceMatrix
(
scalarField& result,
const bool add,
const scalarField& psiInternal,
const scalarField& coeffs,
const direction cmpt,
const Pstream::commsTypes
) const
{
scalarField pnf(this->size());
const labelUList& nbrFaceCells =
this->cyclicPatch().neighbFvPatch().faceCells();
const Field<vector>& iField = this->primitiveField();
// only apply jump to original field
if (&psiInternal == &(iField.component(cmpt).ref()))
{
Field<vector> jf(this->jump());
if (!this->cyclicPatch().owner())
{
jf *= -1.0;
}
forAll(*this, facei)
{
pnf[facei] =
psiInternal[nbrFaceCells[facei]]
- jf[facei].component(cmpt);
}
}
else
{
forAll(*this, facei)
{
pnf[facei] = psiInternal[nbrFaceCells[facei]];
}
}
// Transform according to the transformation tensors
this->transformCoupleField(pnf, cmpt);
// Multiply the field by coefficients and add into the result
this->addToInternalField(result, !add, coeffs, pnf);
}
// ************************************************************************* //
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2018 OpenCFD Ltd.
\\/ 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 "swirlFanVelocityFvPatchField.H"
#include "addToRunTimeSelectionTable.H"
#include "fvPatchFieldMapper.H"
#include "volFields.H"
#include "surfaceFields.H"
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
void Foam::swirlFanVelocityFvPatchField::calcFanJump()
{
if (this->cyclicPatch().owner())
{
const surfaceScalarField& phi =
db().lookupObject<surfaceScalarField>(phiName_);
const fvPatchField<scalar>& pOwner =
patch().lookupPatchField<volScalarField, scalar>(pName_);
const label nbrIndex = this->cyclicPatch().neighbPatchID();
const fvPatch& nbrPatch = patch().boundaryMesh()[nbrIndex];
const fvPatchField<scalar>& pSlave =
nbrPatch.lookupPatchField<volScalarField, scalar>(pName_);
scalarField deltaP(mag(pOwner - pSlave));
if (phi.dimensions() == dimMass/dimTime)
{
deltaP /=
patch().lookupPatchField<volScalarField, scalar>(rhoName_);
}
const vector axisHat =
gSum(patch().nf()*patch().magSf())/gSum(patch().magSf());
vectorField tanDir
(
axisHat ^ (patch().Cf() - origin_)
);
tanDir /= (mag(tanDir) + SMALL);
scalarField magTangU(patch().size(), 0.0);
if (useRealRadius_)
{
const scalarField rMag(mag(patch().Cf() - origin_));
forAll (rMag, i)
{
if (rMag[i] > rInner_ && rMag[i] < rOuter_)
{
magTangU =
deltaP
/rMag[i]
/fanEff_
/rpm_*constant::mathematical::pi/30.0;
}
}
}
else
{
magTangU =
deltaP/rEff_/fanEff_/rpm_*constant::mathematical::pi/30.0;
}
// U upstream from fan
const vectorField Uup(this->patchInternalField());
// U normal to fan
const vectorField Un((Uup & patch().nf())*patch().nf());
// U tangential into fan
const vectorField UtanIn(Uup - Un);
// Calculate the tangential velocity
const vectorField tangentialVelocity(magTangU*tanDir);
this->jump_ = tangentialVelocity + UtanIn;
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::swirlFanVelocityFvPatchField::swirlFanVelocityFvPatchField
(
const fvPatch& p,
const DimensionedField<vector, volMesh>& iF
)
:
fixedJumpFvPatchField<vector>(p, iF),
phiName_("phi"),
pName_("p"),
rhoName_("rho"),
origin_(),
rpm_(0.0),
rEff_(0.0),
fanEff_(1.0),
useRealRadius_(false),
rInner_(0.0),
rOuter_(0.0)
{}
Foam::swirlFanVelocityFvPatchField::swirlFanVelocityFvPatchField
(
const fvPatch& p,
const DimensionedField<vector, volMesh>& iF,
const dictionary& dict
)
:
fixedJumpFvPatchField<vector>(p, iF, dict),
phiName_(dict.lookupOrDefault<word>("phi", "phi")),
pName_(dict.lookupOrDefault<word>("p", "p")),
rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
origin_
(
dict.lookupOrDefault
(
"origin",
patch().size()
? gSum(patch().Cf()*patch().magSf())/gSum(patch().magSf())
: Zero
)
),
rpm_(dict.lookupOrDefault<scalar>("rpm", 0.0)),
rEff_(dict.lookupOrDefault<scalar>("rEff", 0.0)),
fanEff_(dict.lookupOrDefault<scalar>("fanEff", 1.0)),
useRealRadius_(dict.lookupOrDefault<Switch>("useRealRadius", false)),
rInner_(dict.lookupOrDefault<scalar>("rInner", 0.0)),
rOuter_(dict.lookupOrDefault<scalar>("rOuter", 0.0))
{
if (dict.found("value"))
{
fvPatchVectorField::operator=
(
vectorField("value", dict, p.size())
);
}
else
{
this->evaluate(Pstream::commsTypes::blocking);
}
}
Foam::swirlFanVelocityFvPatchField::swirlFanVelocityFvPatchField
(
const swirlFanVelocityFvPatchField& ptf,
const fvPatch& p,
const DimensionedField<vector, volMesh>& iF,
const fvPatchFieldMapper& mapper
)
:
fixedJumpFvPatchField<vector>(ptf, p, iF, mapper),
phiName_(ptf.phiName_),
pName_(ptf.pName_),
rhoName_(ptf.rhoName_),
origin_(ptf.origin_),
rpm_(ptf.rpm_),
rEff_(ptf.rEff_),
fanEff_(ptf.fanEff_),
useRealRadius_(ptf.useRealRadius_),
rInner_(ptf.rInner_),
rOuter_(ptf.rOuter_)
{}
Foam::swirlFanVelocityFvPatchField::swirlFanVelocityFvPatchField
(
const swirlFanVelocityFvPatchField& ptf
)
:
fixedJumpFvPatchField<vector>(ptf),
phiName_(ptf.phiName_),
rhoName_(ptf.rhoName_),
origin_(ptf.origin_),
rpm_(ptf.rpm_),
rEff_(ptf.rEff_),
useRealRadius_(ptf.useRealRadius_),
rInner_(ptf.rInner_),
rOuter_(ptf.rOuter_)
{}
Foam::swirlFanVelocityFvPatchField::swirlFanVelocityFvPatchField
(
const swirlFanVelocityFvPatchField& ptf,
const DimensionedField<vector, volMesh>& iF
)
:
fixedJumpFvPatchField<vector>(ptf, iF),
phiName_(ptf.phiName_),
rhoName_(ptf.rhoName_),
origin_(ptf.origin_),
rpm_(ptf.rpm_),
rEff_(ptf.rEff_),
useRealRadius_(ptf.useRealRadius_),
rInner_(ptf.rInner_),
rOuter_(ptf.rOuter_)
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
void Foam::swirlFanVelocityFvPatchField::updateCoeffs()
{
if (this->updated())
{
return;
}
calcFanJump();
}
void Foam::swirlFanVelocityFvPatchField::write(Ostream& os) const
{
fixedJumpFvPatchField<vector>::write(os);
if (this->cyclicPatch().owner())
{
os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
os.writeEntryIfDifferent<word>("p", "p", pName_);
os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
os.writeEntry("origin", origin_);
os.writeEntry("rpm", rpm_);
os.writeEntryIfDifferent<scalar>("rEff", 0.0, rEff_);
os.writeEntryIfDifferent<bool>("useRealRadius", false, useRealRadius_);
os.writeEntryIfDifferent<scalar>("rInner", 0.0, rInner_);
os.writeEntryIfDifferent<scalar>("rOuter", 0.0, rOuter_);
}
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
makePatchTypeField
(
fvPatchVectorField,
swirlFanVelocityFvPatchField
);
}
// ************************************************************************* //
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2018 OpenCFD Ltd.
\\/ 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::swirlFanVelocityFvPatchField
Group
grpCoupledBoundaryConditions
Description
This boundary condition provides a jump condition for U across a
cyclic pressure jump condition and applies a transformation to U.
The U-jump is specified with a swirl component as follows:
Utan = deltaP/rEff/fanEff/rpm/pi/30.0;
where:
deltaP : pressure drop across the cyclic.
rEff : effective radius
fanEff : fan efficiency coefficient
rmp : RPM of the fan
Alternatively an inner and outer radii can be used instead of rEff. The
Utan is as follow for r > rInner and r < rOuter
Utan = deltaP/r/fanEff/rpm/pi/30.0;
where:
r : p - origin, p is the face center
Outside rInner and rOuter, Utan = 0. The input for this mode is:
useRealRadius true;
rInner 0.005;
rOuter 0.01;
The radial velocity is zero in the present model.
Usage
\table
Property | Description | Required | Default value
patchType | underlying patch type should be \c cyclic| yes |
phi | flux field name | no | phi
rho | density field name | no | rho
p | pressure field name | no | p
origin | fan centre | no | calculated
rpm | RPM of the fan | yes
rEff | Effective radius | no | 0.0
fanEff | Fan efficiency | no | 1.0
useRealRadius| Use inner/outer radii | no | false
rInner | Inner radius | no | 0.0
rOuter | Outer radius | no | 0.0
\endtable
Example of the boundary condition specification:
\verbatim
<patchName>
{
cyclicFaces_master
{
type swirlFanVelocity;
patchType cyclic;
jump uniform (0 0 0);
value uniform (0 0 0);
rpm 1000;
rEff 0.01;
}
}
\endverbatim
SourceFiles
swirlFanVelocityFvPatchField.C
\*---------------------------------------------------------------------------*/
#ifndef swirlFanVelocityFvPatchField_H
#define swirlFanVelocityFvPatchField_H
#include "fixedJumpFvPatchField.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
/*---------------------------------------------------------------------------*\
Class swirlFanVelocityFvPatchField Declaration
\*---------------------------------------------------------------------------*/
class swirlFanVelocityFvPatchField
:
public fixedJumpFvPatchField<vector>
{
// Private data
//- Name of the flux field
const word phiName_;
//- Name of the pressure field
const word pName_;
//- Name of the rho field
const word rhoName_;
//- Origin of the rotation
const vector origin_;
//- Fan rpm
scalar rpm_;
//- Effective fan radius
scalar rEff_;
//- Fan efficiency
scalar fanEff_;
//- Switch to use effective radius or inner and outer radius
Switch useRealRadius_;
//- Inner radius
scalar rInner_;
//- Outer radius
scalar rOuter_;
// Private Member Functions
//- Calculate the fan pressure jump
void calcFanJump();
public:
//- Runtime type information
TypeName("swirlFanVelocity");
// Constructors
//- Construct from patch and internal field
swirlFanVelocityFvPatchField
(
const fvPatch&,
const DimensionedField<vector, volMesh>&
);
//- Construct from patch, internal field and dictionary
swirlFanVelocityFvPatchField
(
const fvPatch&,
const DimensionedField<vector, volMesh>&,
const dictionary&
);
//- Construct by mapping given swirlFanVelocityFvPatchField onto a new patch
swirlFanVelocityFvPatchField
(
const swirlFanVelocityFvPatchField&,
const fvPatch&,
const DimensionedField<vector, volMesh>&,
const fvPatchFieldMapper&
);
//- Construct as copy
swirlFanVelocityFvPatchField
(
const swirlFanVelocityFvPatchField&
);
//- Construct and return a clone
virtual tmp<fvPatchField<vector>> clone() const
{
return tmp<fvPatchField<vector>>
(
new swirlFanVelocityFvPatchField(*this)
);
}
//- Construct as copy setting internal field reference
swirlFanVelocityFvPatchField
(
const swirlFanVelocityFvPatchField&,
const DimensionedField<vector, volMesh>&
);
//- Construct and return a clone setting internal field reference
virtual tmp<fvPatchField<vector>> clone
(
const DimensionedField<vector, volMesh>& iF
) const
{
return tmp<fvPatchField<vector>>
(
new swirlFanVelocityFvPatchField(*this, iF)
);
}
// Member functions
//- Update the coefficients associated with the patch field
virtual void updateCoeffs();
//- Write
virtual void write(Ostream&) const;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
#endif
// ************************************************************************* //
0% or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment