From f54b190f319a8a07fef208d25c64a9ad02f7f529 Mon Sep 17 00:00:00 2001 From: Andrew Heather <> Date: Mon, 17 May 2021 16:03:20 +0100 Subject: [PATCH] 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; } --- .../porousBafflePressureFvPatchField.C | 13 +++++--- .../derived/fan/fanFvPatchFields.C | 4 +-- .../derived/fixedJump/fixedJumpFvPatchField.C | 30 +++++++++++++++++++ .../derived/fixedJump/fixedJumpFvPatchField.H | 15 ++++++++-- .../swirlFanVelocityFvPatchField.C | 2 +- .../uniformJump/uniformJumpFvPatchField.C | 2 +- .../energyJump/energyJumpFvPatchScalarField.C | 2 +- 7 files changed, 56 insertions(+), 12 deletions(-) diff --git a/src/TurbulenceModels/turbulenceModels/derivedFvPatchFields/porousBafflePressure/porousBafflePressureFvPatchField.C b/src/TurbulenceModels/turbulenceModels/derivedFvPatchFields/porousBafflePressure/porousBafflePressureFvPatchField.C index 7fe61490348..4a84834687d 100644 --- a/src/TurbulenceModels/turbulenceModels/derivedFvPatchFields/porousBafflePressure/porousBafflePressureFvPatchField.C +++ b/src/TurbulenceModels/turbulenceModels/derivedFvPatchFields/porousBafflePressure/porousBafflePressureFvPatchField.C @@ -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() << ':' diff --git a/src/finiteVolume/fields/fvPatchFields/derived/fan/fanFvPatchFields.C b/src/finiteVolume/fields/fvPatchFields/derived/fan/fanFvPatchFields.C index c1d7c100544..53ec3fac7f6 100644 --- a/src/finiteVolume/fields/fvPatchFields/derived/fan/fanFvPatchFields.C +++ b/src/finiteVolume/fields/fvPatchFields/derived/fan/fanFvPatchFields.C @@ -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(); diff --git a/src/finiteVolume/fields/fvPatchFields/derived/fixedJump/fixedJumpFvPatchField.C b/src/finiteVolume/fields/fvPatchFields/derived/fixedJump/fixedJumpFvPatchField.C index 49fdedd3be4..d96fe621213 100644 --- a/src/finiteVolume/fields/fvPatchFields/derived/fixedJump/fixedJumpFvPatchField.C +++ b/src/finiteVolume/fields/fvPatchFields/derived/fixedJump/fixedJumpFvPatchField.C @@ -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); } diff --git a/src/finiteVolume/fields/fvPatchFields/derived/fixedJump/fixedJumpFvPatchField.H b/src/finiteVolume/fields/fvPatchFields/derived/fixedJump/fixedJumpFvPatchField.H index baef3236e66..8c5470782f5 100644 --- a/src/finiteVolume/fields/fvPatchFields/derived/fixedJump/fixedJumpFvPatchField.H +++ b/src/finiteVolume/fields/fvPatchFields/derived/fixedJump/fixedJumpFvPatchField.H @@ -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; diff --git a/src/finiteVolume/fields/fvPatchFields/derived/swirlFanVelocity/swirlFanVelocityFvPatchField.C b/src/finiteVolume/fields/fvPatchFields/derived/swirlFanVelocity/swirlFanVelocityFvPatchField.C index a7bcfd3a20f..b75e4b32225 100644 --- a/src/finiteVolume/fields/fvPatchFields/derived/swirlFanVelocity/swirlFanVelocityFvPatchField.C +++ b/src/finiteVolume/fields/fvPatchFields/derived/swirlFanVelocity/swirlFanVelocityFvPatchField.C @@ -104,7 +104,7 @@ void Foam::swirlFanVelocityFvPatchField::calcFanJump() // Calculate the tangential velocity const vectorField tangentialVelocity(magTangU*tanDir); - this->jump_ = tangentialVelocity; + this->setJump(tangentialVelocity); } } diff --git a/src/finiteVolume/fields/fvPatchFields/derived/uniformJump/uniformJumpFvPatchField.C b/src/finiteVolume/fields/fvPatchFields/derived/uniformJump/uniformJumpFvPatchField.C index 74ae08f31b8..75758640834 100644 --- a/src/finiteVolume/fields/fvPatchFields/derived/uniformJump/uniformJumpFvPatchField.C +++ b/src/finiteVolume/fields/fvPatchFields/derived/uniformJump/uniformJumpFvPatchField.C @@ -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(); diff --git a/src/thermophysicalModels/basic/derivedFvPatchFields/energyJump/energyJump/energyJumpFvPatchScalarField.C b/src/thermophysicalModels/basic/derivedFvPatchFields/energyJump/energyJump/energyJumpFvPatchScalarField.C index 7ccaecc93f5..5fd72a0f7fd 100644 --- a/src/thermophysicalModels/basic/derivedFvPatchFields/energyJump/energyJump/energyJumpFvPatchScalarField.C +++ b/src/thermophysicalModels/basic/derivedFvPatchFields/energyJump/energyJump/energyJumpFvPatchScalarField.C @@ -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(); -- GitLab