ENH: adding non-dimensional option to fanFvPatchField.C

This commit is contained in:
sergio 2017-11-15 12:54:09 -08:00 committed by Mark Olesen
parent ad00fed4aa
commit 219c520f0b
3 changed files with 93 additions and 11 deletions

View File

@ -3,7 +3,7 @@
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
\\/ M anipulation |
\\/ M anipulation | Copyright (C) 2017 OpenCFD Ltd.
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
@ -49,7 +49,10 @@ Foam::fanFvPatchField<Type>::fanFvPatchField
uniformJumpFvPatchField<Type>(p, iF),
phiName_("phi"),
rhoName_("rho"),
uniformJump_(false)
uniformJump_(false),
nonDimensional_(false),
rpm_(0.0),
dm_(0.0)
{}
@ -64,8 +67,17 @@ Foam::fanFvPatchField<Type>::fanFvPatchField
uniformJumpFvPatchField<Type>(p, iF, dict),
phiName_(dict.lookupOrDefault<word>("phi", "phi")),
rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
uniformJump_(dict.lookupOrDefault<bool>("uniformJump", false))
{}
uniformJump_(dict.lookupOrDefault<bool>("uniformJump", false)),
nonDimensional_(dict.lookupOrDefault<Switch>("nonDimensional", false)),
rpm_(dict.lookupOrDefault<scalar>("rpm", 0.0)),
dm_(dict.lookupOrDefault<scalar>("dm", 0.0))
{
if (nonDimensional_)
{
dict.lookup("rpm") >> rpm_;
dict.lookup("dm") >> dm_;
}
}
template<class Type>
@ -80,7 +92,10 @@ Foam::fanFvPatchField<Type>::fanFvPatchField
uniformJumpFvPatchField<Type>(ptf, p, iF, mapper),
phiName_(ptf.phiName_),
rhoName_(ptf.rhoName_),
uniformJump_(ptf.uniformJump_)
uniformJump_(ptf.uniformJump_),
nonDimensional_(ptf.nonDimensional_),
rpm_(ptf.rpm_),
dm_(ptf.dm_)
{}
@ -93,7 +108,10 @@ Foam::fanFvPatchField<Type>::fanFvPatchField
uniformJumpFvPatchField<Type>(ptf),
phiName_(ptf.phiName_),
rhoName_(ptf.rhoName_),
uniformJump_(ptf.uniformJump_)
uniformJump_(ptf.uniformJump_),
nonDimensional_(ptf.nonDimensional_),
rpm_(ptf.rpm_),
dm_(ptf.dm_)
{}
@ -107,7 +125,10 @@ Foam::fanFvPatchField<Type>::fanFvPatchField
uniformJumpFvPatchField<Type>(ptf, iF),
phiName_(ptf.phiName_),
rhoName_(ptf.rhoName_),
uniformJump_(ptf.uniformJump_)
uniformJump_(ptf.uniformJump_),
nonDimensional_(ptf.nonDimensional_),
rpm_(ptf.rpm_),
dm_(ptf.dm_)
{}
@ -138,6 +159,12 @@ void Foam::fanFvPatchField<Type>::write(Ostream& os) const
(
os, "uniformJump", false, uniformJump_
);
this->template writeEntryIfDifferent<bool>
(
os, "nonDimensional", false, nonDimensional_
);
this->template writeEntryIfDifferent<scalar>(os, "rpm", 0, rpm_);
this->template writeEntryIfDifferent<scalar>(os, "dm", 0, dm_);
}

View File

@ -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) 2017 OpenCFD Ltd.
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
@ -34,6 +34,24 @@ Description
The jump is specified as a \c Function1 type, to enable the use of, e.g.
contant, polynomial, table values.
The switch nonDimensional can be used for a non-dimensional table. It needs
inputs rpm and dm of the fan. It should be used with uniformJump true.
The nonDimensional flux for the table is calculate as :
phi = 4.0*mDot/(rho*sqr(PI)*dm^3*omega)
where:
dm is the mean diameter.
omega is rad/sec.
The nonDimensinal pressure :
Psi = 2 deltaP/(rho*(sqr(PI*omega*dm)))
where:
deltaP is the pressure drop
The non-dimensional table should be given as Psi = F(phi).
Usage
\table
Property | Description | Required | Default value
@ -43,6 +61,9 @@ Usage
rho | density field name | no | none
uniformJump | applies a uniform pressure based on the averaged
velocity | no | false
nonDimensional | uses non-dimensional table | no | false
rpm | fan rpm for non-dimensional table | no | 0.0
dm | mean diameter for non-dimensional table | no | 0.0
\endtable
Example of the boundary condition specification:
@ -114,6 +135,17 @@ class fanFvPatchField
//- Uniform pressure drop
bool uniformJump_;
//- Swtich for using non-dimensional curve
Switch nonDimensional_;
// Parameters for non-dimensional table
//- Fan rpm
scalar rpm_;
//- Fan mean diameter
scalar dm_;
// Private Member Functions

View File

@ -3,7 +3,7 @@
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2017 OpenFOAM Foundation
\\/ M anipulation |
\\/ M anipulation | Copyright (C) 2017 OpenCFD Ltd.
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
@ -48,6 +48,15 @@ void Foam::fanFvPatchField<Foam::scalar>::calcFanJump()
{
scalar area = gSum(patch().magSf());
Un = gSum(Un*patch().magSf())/area;
if (nonDimensional_)
{
// Create an adimensional volumetric flow rate
Un =
120.0*Un/pow3(constant::mathematical::pi)
* patch().magSf()
/ pow3(dm_);
}
}
if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
@ -55,9 +64,20 @@ void Foam::fanFvPatchField<Foam::scalar>::calcFanJump()
Un /= patch().lookupPatchField<volScalarField, scalar>(rhoName_);
}
if (nonDimensional_)
{
scalarField deltap = this->jumpTable_->value(Un);
// Convert adimensional deltap from curve into deltaP
scalarField pdFan =
deltap*pow4(constant::mathematical::pi)*rpm_*sqr(dm_)/1800;
this->jump_ = max(pdFan, scalar(0));
}
else
{
this->jump_ = max(this->jumpTable_->value(Un), scalar(0));
}
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
@ -73,7 +93,10 @@ Foam::fanFvPatchField<Foam::scalar>::fanFvPatchField
uniformJumpFvPatchField<scalar>(p, iF),
phiName_(dict.lookupOrDefault<word>("phi", "phi")),
rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
uniformJump_(dict.lookupOrDefault<bool>("uniformJump", false))
uniformJump_(dict.lookupOrDefault<bool>("uniformJump", false)),
nonDimensional_(dict.lookupOrDefault<Switch>("nonDimensional", false)),
rpm_(dict.lookupOrDefault<scalar>("rpm", 0.0)),
dm_(dict.lookupOrDefault<scalar>("dm", 0.0))
{
if (this->cyclicPatch().owner())
{