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

ENH: Added under-relaxation to jump conditions

Currently only applied to the fanFvPatchField, e.g.

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

        // Optional under-relaxation
        relax           0.2;

        ...
    }
parent 30da494e
No related branches found
No related tags found
1 merge request!467ENH: Added under-relaxation to jump conditions
......@@ -6,7 +6,7 @@
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2011-2017 OpenFOAM Foundation
Copyright (C) 2017-2020 OpenCFD Ltd.
Copyright (C) 2017-2021 OpenCFD Ltd.
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
......@@ -80,6 +80,8 @@ void Foam::fanFvPatchField<Foam::scalar>::calcFanJump()
{
this->jump_ = this->jumpTable_->value(Un);
}
this->relax();
}
}
......@@ -94,7 +96,7 @@ Foam::fanFvPatchField<Foam::scalar>::fanFvPatchField
const dictionary& dict
)
:
uniformJumpFvPatchField<scalar>(p, iF),
uniformJumpFvPatchField<scalar>(p, iF, dict),
phiName_(dict.getOrDefault<word>("phi", "phi")),
rhoName_(dict.getOrDefault<word>("rho", "rho")),
uniformJump_(dict.getOrDefault("uniformJump", false)),
......
......@@ -6,6 +6,7 @@
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2011-2017 OpenFOAM Foundation
Copyright (C) 2021 OpenCFD Ltd.
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
......@@ -37,7 +38,10 @@ Foam::fixedJumpFvPatchField<Type>::fixedJumpFvPatchField
)
:
jumpCyclicFvPatchField<Type>(p, iF),
jump_(this->size(), Zero)
jump_(this->size(), Zero),
jump0_(this->size(), Zero),
relaxFactor_(-1),
timeIndex_(-1)
{}
......@@ -51,7 +55,10 @@ Foam::fixedJumpFvPatchField<Type>::fixedJumpFvPatchField
)
:
jumpCyclicFvPatchField<Type>(ptf, p, iF, mapper),
jump_(ptf.jump_, mapper)
jump_(ptf.jump_, mapper),
jump0_(ptf.jump0_, mapper),
relaxFactor_(ptf.relaxFactor_),
timeIndex_(ptf.timeIndex_)
{}
......@@ -63,12 +70,20 @@ Foam::fixedJumpFvPatchField<Type>::fixedJumpFvPatchField
const dictionary& dict
)
:
jumpCyclicFvPatchField<Type>(p, iF),
jump_(p.size(), Zero)
jumpCyclicFvPatchField<Type>(p, iF, dict),
jump_(p.size(), Zero),
jump0_(p.size(), Zero),
relaxFactor_(dict.getOrDefault<scalar>("relax", -1)),
timeIndex_(this->db().time().timeIndex())
{
if (this->cyclicPatch().owner())
{
jump_ = Field<Type>("jump", dict, p.size());
if (dict.found("jump0"))
{
jump0_ = Field<Type>("jump0", dict, p.size());
}
}
if (dict.found("value"))
......@@ -92,7 +107,10 @@ Foam::fixedJumpFvPatchField<Type>::fixedJumpFvPatchField
)
:
jumpCyclicFvPatchField<Type>(ptf),
jump_(ptf.jump_)
jump_(ptf.jump_),
jump0_(ptf.jump0_),
relaxFactor_(ptf.relaxFactor_),
timeIndex_(ptf.timeIndex_)
{}
......@@ -104,7 +122,10 @@ Foam::fixedJumpFvPatchField<Type>::fixedJumpFvPatchField
)
:
jumpCyclicFvPatchField<Type>(ptf, iF),
jump_(ptf.jump_)
jump_(ptf.jump_),
jump0_(ptf.jump0_),
relaxFactor_(ptf.relaxFactor_),
timeIndex_(ptf.timeIndex_)
{}
......@@ -127,6 +148,49 @@ Foam::tmp<Foam::Field<Type>> Foam::fixedJumpFvPatchField<Type>::jump() const
}
template<class Type>
Foam::tmp<Foam::Field<Type>> Foam::fixedJumpFvPatchField<Type>::jump0() const
{
if (this->cyclicPatch().owner())
{
return jump0_;
}
else
{
return refCast<const fixedJumpFvPatchField<Type>>
(
this->neighbourPatchField()
).jump0();
}
}
template<class Type>
Foam::scalar Foam::fixedJumpFvPatchField<Type>::relaxFactor() const
{
return relaxFactor_;
}
template<class Type>
void Foam::fixedJumpFvPatchField<Type>::relax()
{
if (!this->cyclicPatch().owner() || relaxFactor_ < 0)
{
return;
}
jump_ = relaxFactor_*jump_ + (1 - relaxFactor_)*jump0_;
if (timeIndex_ != this->db().time().timeIndex())
{
jump0_ = jump_;
timeIndex_ = this->db().time().timeIndex();
}
}
template<class Type>
void Foam::fixedJumpFvPatchField<Type>::autoMap
(
......@@ -135,6 +199,7 @@ void Foam::fixedJumpFvPatchField<Type>::autoMap
{
jumpCyclicFvPatchField<Type>::autoMap(m);
jump_.autoMap(m);
jump0_.autoMap(m);
}
......@@ -147,9 +212,9 @@ void Foam::fixedJumpFvPatchField<Type>::rmap
{
jumpCyclicFvPatchField<Type>::rmap(ptf, addr);
const fixedJumpFvPatchField<Type>& tiptf =
refCast<const fixedJumpFvPatchField<Type>>(ptf);
jump_.rmap(tiptf.jump_, addr);
const auto& fjptf = refCast<const fixedJumpFvPatchField<Type>>(ptf);
jump_.rmap(fjptf.jump_, addr);
jump0_.rmap(fjptf.jump0_, addr);
}
......@@ -157,11 +222,22 @@ template<class Type>
void Foam::fixedJumpFvPatchField<Type>::write(Ostream& os) const
{
fvPatchField<Type>::write(os);
os.writeEntry("patchType", this->interfaceFieldType());
// Write patchType if not done already by fvPatchField
if (!this->patchType().size())
{
os.writeEntry("patchType", this->interfaceFieldType());
}
if (this->cyclicPatch().owner())
{
jump_.writeEntry("jump", os);
if (relaxFactor_ > 0)
{
os.writeEntry("relax", relaxFactor_);
jump0_.writeEntry("jump0", os);
}
}
this->writeEntry("value", os);
......
......@@ -6,6 +6,7 @@
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2011-2016 OpenFOAM Foundation
Copyright (C) 2021 OpenCFD Ltd.
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
......@@ -41,6 +42,7 @@ Usage
Property | Description | Required | Default value
patchType | underlying patch type should be \c cyclic| yes |
jump | current jump value | yes |
relax | under-relaxation factor | no |
\endtable
Example of the boundary condition specification:
......@@ -93,6 +95,15 @@ protected:
//- "jump" field
Field<Type> jump_;
//- "jump" field at old time level
Field<Type> jump0_;
//- Under-relaxation factor
scalar relaxFactor_;
//- Time index
label timeIndex_;
public:
......@@ -168,6 +179,15 @@ public:
//- Return the "jump" across the patch
virtual tmp<Field<Type>> jump() const;
//- Return the old time "jump" across the patch
virtual tmp<Field<Type>> jump0() const;
//- Return the under-relaxation factor
virtual scalar relaxFactor() const;
//- Return the relaxed "jump" across the patch
virtual void relax();
// Mapping functions
......
......@@ -63,7 +63,7 @@ Foam::uniformJumpFvPatchField<Type>::uniformJumpFvPatchField
const dictionary& dict
)
:
fixedJumpFvPatchField<Type>(p, iF),
fixedJumpFvPatchField<Type>(p, iF, dict),
jumpTable_()
{
if (this->cyclicPatch().owner())
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment