Skip to content
Snippets Groups Projects
Commit f54b190f authored by Andrew Heather's avatar Andrew Heather Committed by Andrew Heather
Browse files

ENH: Jump conditions - added optional minimum jump value

Example:

    plane
    {
        type            fan;
        patchType       cyclic;
        jump            uniform 0;
        value           uniform 0;
        uniformJump     false;

        relax           0.2;

        jumpTable       ( (1000 0) (-200 1) (20 2) (-0.75 3) );

        // Optional minimum jump value
        minJump         0;
    }
parent 05e0749d
No related branches found
No related tags found
1 merge request!467ENH: Added under-relaxation to jump conditions
......@@ -162,21 +162,26 @@ void Foam::porousBafflePressureFvPatchField::updateCoeffs()
const scalar D = D_->value(t);
const scalar I = I_->value(t);
jump_ =
setJump
(
-sign(Un)
*(
D*turbModel.nu(patch().index())
+ I*0.5*magUn
)*magUn*length_;
)*magUn*length_
);
if (internalField().dimensions() == dimPressure)
{
jump_ *= patch().lookupPatchField<volScalarField, scalar>(rhoName_);
setJump
(
jump()*patch().lookupPatchField<volScalarField, scalar>(rhoName_)
);
}
if (debug)
{
scalar avePressureJump = gAverage(jump_);
scalar avePressureJump = gAverage(jump());
scalar aveVelocity = gAverage(Un);
Info<< patch().boundaryMesh().mesh().name() << ':'
......
......@@ -74,11 +74,11 @@ void Foam::fanFvPatchField<Foam::scalar>::calcFanJump()
deltap*pow4(constant::mathematical::pi)*sqr(dm_*rpm_)/1800.0
);
this->jump_ = pdFan;
this->setJump(pdFan);
}
else
{
this->jump_ = this->jumpTable_->value(Un);
this->setJump(jumpTable_->value(Un));
}
this->relax();
......
......@@ -40,6 +40,7 @@ Foam::fixedJumpFvPatchField<Type>::fixedJumpFvPatchField
jumpCyclicFvPatchField<Type>(p, iF),
jump_(this->size(), Zero),
jump0_(this->size(), Zero),
minJump_(pTraits<Type>::min),
relaxFactor_(-1),
timeIndex_(-1)
{}
......@@ -57,6 +58,7 @@ Foam::fixedJumpFvPatchField<Type>::fixedJumpFvPatchField
jumpCyclicFvPatchField<Type>(ptf, p, iF, mapper),
jump_(ptf.jump_, mapper),
jump0_(ptf.jump0_, mapper),
minJump_(ptf.minJump_),
relaxFactor_(ptf.relaxFactor_),
timeIndex_(ptf.timeIndex_)
{}
......@@ -73,6 +75,7 @@ Foam::fixedJumpFvPatchField<Type>::fixedJumpFvPatchField
jumpCyclicFvPatchField<Type>(p, iF, dict),
jump_(p.size(), Zero),
jump0_(p.size(), Zero),
minJump_(dict.getOrDefault<Type>("minJump", pTraits<Type>::min)),
relaxFactor_(dict.getOrDefault<scalar>("relax", -1)),
timeIndex_(this->db().time().timeIndex())
{
......@@ -109,6 +112,7 @@ Foam::fixedJumpFvPatchField<Type>::fixedJumpFvPatchField
jumpCyclicFvPatchField<Type>(ptf),
jump_(ptf.jump_),
jump0_(ptf.jump0_),
minJump_(ptf.minJump_),
relaxFactor_(ptf.relaxFactor_),
timeIndex_(ptf.timeIndex_)
{}
......@@ -124,6 +128,7 @@ Foam::fixedJumpFvPatchField<Type>::fixedJumpFvPatchField
jumpCyclicFvPatchField<Type>(ptf, iF),
jump_(ptf.jump_),
jump0_(ptf.jump0_),
minJump_(ptf.minJump_),
relaxFactor_(ptf.relaxFactor_),
timeIndex_(ptf.timeIndex_)
{}
......@@ -131,6 +136,26 @@ Foam::fixedJumpFvPatchField<Type>::fixedJumpFvPatchField
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
template<class Type>
void Foam::fixedJumpFvPatchField<Type>::setJump(const Field<Type>& jump)
{
if (this->cyclicPatch().owner())
{
jump_ = max(jump, minJump_);
}
}
template<class Type>
void Foam::fixedJumpFvPatchField<Type>::setJump(const Type& jump)
{
if (this->cyclicPatch().owner())
{
jump_ = max(jump, minJump_);
}
}
template<class Type>
Foam::tmp<Foam::Field<Type>> Foam::fixedJumpFvPatchField<Type>::jump() const
{
......@@ -240,6 +265,11 @@ void Foam::fixedJumpFvPatchField<Type>::write(Ostream& os) const
}
}
if (minJump_ != pTraits<Type>::min)
{
os.writeEntry("minJump", minJump_);
}
this->writeEntry("value", os);
}
......
......@@ -43,6 +43,7 @@ Usage
patchType | underlying patch type should be \c cyclic| yes |
jump | current jump value | yes |
relax | under-relaxation factor | no |
minJump | Minimum jump value | no |
\endtable
Example of the boundary condition specification:
......@@ -88,9 +89,7 @@ class fixedJumpFvPatchField
public jumpCyclicFvPatchField<Type>
{
protected:
// Protected data
// Private data
//- "jump" field
Field<Type> jump_;
......@@ -98,6 +97,10 @@ protected:
//- "jump" field at old time level
Field<Type> jump0_;
//- Minimum allowable jump value
Type minJump_;
//- Under-relaxation factor
scalar relaxFactor_;
......@@ -176,6 +179,12 @@ public:
// Access
//- Set the jump field
virtual void setJump(const Field<Type>& jump);
//- Set the jump field (uniform value)
virtual void setJump(const Type& jump);
//- Return the "jump" across the patch
virtual tmp<Field<Type>> jump() const;
......
......@@ -104,7 +104,7 @@ void Foam::swirlFanVelocityFvPatchField::calcFanJump()
// Calculate the tangential velocity
const vectorField tangentialVelocity(magTangU*tanDir);
this->jump_ = tangentialVelocity;
this->setJump(tangentialVelocity);
}
}
......
......@@ -120,7 +120,7 @@ void Foam::uniformJumpFvPatchField<Type>::updateCoeffs()
if (this->cyclicPatch().owner())
{
this->jump_ = jumpTable_->value(this->db().time().value());
this->setJump(jumpTable_->value(this->db().time().value()));
}
fixedJumpFvPatchField<Type>::updateCoeffs();
......
......@@ -125,7 +125,7 @@ void Foam::energyJumpFvPatchScalarField::updateCoeffs()
const labelUList& faceCells = this->patch().faceCells();
jump_ = thermo.he(pp, Tbp.jump(), faceCells);
setJump(thermo.he(pp, Tbp.jump(), faceCells));
}
fixedJumpFvPatchField<scalar>::updateCoeffs();
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment