mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
ENH: rotor momentum source - trim model - enable target forces as well as force coeffs
This commit is contained in:
@ -70,21 +70,29 @@ Foam::vector Foam::targetCoeffTrim::calcCoeffs
|
||||
{
|
||||
label cellI = cells[i];
|
||||
|
||||
// create normalisation coefficient
|
||||
scalar radius = x[i].x();
|
||||
scalar coeff2 = coeff1*pow4(radius);
|
||||
if (compressible)
|
||||
{
|
||||
coeff2 *= trho()[cellI];
|
||||
}
|
||||
|
||||
vector fc = force[cellI];
|
||||
vector mc = fc^(C[cellI] - origin);
|
||||
|
||||
// add to coefficient vector
|
||||
cf[0] += (fc & yawAxis)/(coeff2 + ROOTVSMALL);
|
||||
cf[1] += (mc & pitchAxis)/(coeff2*radius + ROOTVSMALL);
|
||||
cf[2] += (mc & rollAxis)/(coeff2*radius + ROOTVSMALL);
|
||||
if (useCoeffs_)
|
||||
{
|
||||
scalar radius = x[i].x();
|
||||
scalar coeff2 = coeff1*pow4(radius);
|
||||
if (compressible)
|
||||
{
|
||||
coeff2 *= trho()[cellI];
|
||||
}
|
||||
|
||||
// add to coefficient vector
|
||||
cf[0] += (fc & yawAxis)/(coeff2 + ROOTVSMALL);
|
||||
cf[1] += (mc & pitchAxis)/(coeff2*radius + ROOTVSMALL);
|
||||
cf[2] += (mc & rollAxis)/(coeff2*radius + ROOTVSMALL);
|
||||
}
|
||||
else
|
||||
{
|
||||
cf[0] += fc & yawAxis;
|
||||
cf[1] += mc & pitchAxis;
|
||||
cf[2] += mc & rollAxis;
|
||||
}
|
||||
}
|
||||
|
||||
reduce(cf, sumOp<vector>());
|
||||
@ -103,6 +111,7 @@ Foam::targetCoeffTrim::targetCoeffTrim
|
||||
:
|
||||
trimModel(rotor, dict, typeName),
|
||||
calcFrequency_(-1),
|
||||
useCoeffs_(true),
|
||||
target_(vector::zero),
|
||||
theta_(vector::zero),
|
||||
nIter_(50),
|
||||
@ -128,9 +137,16 @@ void Foam::targetCoeffTrim::read(const dictionary& dict)
|
||||
trimModel::read(dict);
|
||||
|
||||
const dictionary& targetDict(coeffs_.subDict("target"));
|
||||
target_[0] = readScalar(targetDict.lookup("thrustCoeff"));
|
||||
target_[1] = readScalar(targetDict.lookup("pitchCoeff"));
|
||||
target_[2] = readScalar(targetDict.lookup("rollCoeff"));
|
||||
useCoeffs_ = targetDict.lookupOrDefault<bool>("useCoeffs", true);
|
||||
word ext = "";
|
||||
if (useCoeffs_)
|
||||
{
|
||||
ext = "Coeff";
|
||||
}
|
||||
|
||||
target_[0] = readScalar(targetDict.lookup("thrust" + ext));
|
||||
target_[1] = readScalar(targetDict.lookup("pitch" + ext));
|
||||
target_[2] = readScalar(targetDict.lookup("roll" + ext));
|
||||
|
||||
const dictionary& pitchAngleDict(coeffs_.subDict("pitchAngles"));
|
||||
theta_[0] = degToRad(readScalar(pitchAngleDict.lookup("theta0Ini")));
|
||||
@ -173,8 +189,14 @@ void Foam::targetCoeffTrim::correct(const vectorField& U, vectorField& force)
|
||||
{
|
||||
if (rotor_.mesh().time().timeIndex() % calcFrequency_ == 0)
|
||||
{
|
||||
word calcType = "forces";
|
||||
if (useCoeffs_)
|
||||
{
|
||||
calcType = "coefficients";
|
||||
}
|
||||
|
||||
Info<< type() << ":" << nl
|
||||
<< " solving for target trim coefficients" << nl;
|
||||
<< " solving for target trim " << calcType << nl;
|
||||
|
||||
const scalar rhoRef = rotor_.rhoRef();
|
||||
|
||||
@ -237,7 +259,7 @@ void Foam::targetCoeffTrim::correct(const vectorField& U, vectorField& force)
|
||||
<< "), iterations = " << iter << endl;
|
||||
}
|
||||
|
||||
Info<< " current and target coefficients:" << nl
|
||||
Info<< " current and target " << calcType << nl
|
||||
<< " thrust = " << old[0]*rhoRef << ", " << target_[0] << nl
|
||||
<< " pitch = " << old[1]*rhoRef << ", " << target_[1] << nl
|
||||
<< " roll = " << old[2]*rhoRef << ", " << target_[2] << nl
|
||||
|
||||
@ -25,7 +25,7 @@ Class
|
||||
Foam::targetCoeffTrim
|
||||
|
||||
Description
|
||||
Target trim coefficients
|
||||
Target trim forces/coefficients
|
||||
|
||||
Solves:
|
||||
|
||||
@ -43,8 +43,8 @@ Description
|
||||
Newton-Raphson iterative method. The solver tolerance can be user-input,
|
||||
using the 'tol' entry.
|
||||
|
||||
|
||||
The coefficients are determined using:
|
||||
If coefficients are requested (useCoeffs = true), the force and moments
|
||||
are normalised using:
|
||||
|
||||
force
|
||||
c = ---------------------------------
|
||||
@ -97,6 +97,9 @@ protected:
|
||||
//- Number of iterations between calls to 'correct'
|
||||
label calcFrequency_;
|
||||
|
||||
//- Flag to indicate whether to solve coeffs (true) or forces (false)
|
||||
bool useCoeffs_;
|
||||
|
||||
//- Target coefficient vector (thrust force, roll moment, pitch moment)
|
||||
vector target_;
|
||||
|
||||
|
||||
Reference in New Issue
Block a user