mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
MRFZone: Added function to convert between absolute and relative velocities.
This commit is contained in:
@ -375,6 +375,45 @@ void Foam::MRFZone::relativeVelocity(volVectorField& U) const
|
||||
}
|
||||
|
||||
|
||||
void Foam::MRFZone::absoluteVelocity(volVectorField& U) const
|
||||
{
|
||||
const volVectorField& C = mesh_.C();
|
||||
|
||||
const vector& origin = origin_.value();
|
||||
const vector& Omega = Omega_.value();
|
||||
|
||||
const labelList& cells = mesh_.cellZones()[cellZoneID_];
|
||||
|
||||
forAll(cells, i)
|
||||
{
|
||||
label celli = cells[i];
|
||||
U[celli] += (Omega ^ (C[celli] - origin));
|
||||
}
|
||||
|
||||
// Included patches
|
||||
forAll(includedFaces_, patchi)
|
||||
{
|
||||
forAll(includedFaces_[patchi], i)
|
||||
{
|
||||
label patchFacei = includedFaces_[patchi][i];
|
||||
U.boundaryField()[patchi][patchFacei] =
|
||||
(Omega ^ (C.boundaryField()[patchi][patchFacei] - origin));
|
||||
}
|
||||
}
|
||||
|
||||
// Excluded patches
|
||||
forAll(excludedFaces_, patchi)
|
||||
{
|
||||
forAll(excludedFaces_[patchi], i)
|
||||
{
|
||||
label patchFacei = excludedFaces_[patchi][i];
|
||||
U.boundaryField()[patchi][patchFacei] +=
|
||||
(Omega ^ (C.boundaryField()[patchi][patchFacei] - origin));
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
void Foam::MRFZone::relativeFlux(surfaceScalarField& phi) const
|
||||
{
|
||||
relativeRhoFlux(geometricOneField(), phi);
|
||||
|
||||
@ -175,6 +175,9 @@ public:
|
||||
//- Make the given absolute velocity relative within the MRF region
|
||||
void relativeVelocity(volVectorField& U) const;
|
||||
|
||||
//- Make the given relative velocity absolute within the MRF region
|
||||
void absoluteVelocity(volVectorField& U) const;
|
||||
|
||||
//- Make the given absolute flux relative within the MRF region
|
||||
void relativeFlux(surfaceScalarField& phi) const;
|
||||
|
||||
|
||||
@ -86,6 +86,15 @@ void Foam::MRFZones::relativeVelocity(volVectorField& U) const
|
||||
}
|
||||
|
||||
|
||||
void Foam::MRFZones::absoluteVelocity(volVectorField& U) const
|
||||
{
|
||||
forAll(*this, i)
|
||||
{
|
||||
operator[](i).absoluteVelocity(U);
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
void Foam::MRFZones::relativeFlux(surfaceScalarField& phi) const
|
||||
{
|
||||
forAll(*this, i)
|
||||
|
||||
@ -81,6 +81,9 @@ public:
|
||||
//- Make the given absolute velocity relative within the MRF region
|
||||
void relativeVelocity(volVectorField& U) const;
|
||||
|
||||
//- Make the given relative velocity absolute within the MRF region
|
||||
void absoluteVelocity(volVectorField& U) const;
|
||||
|
||||
//- Make the given absolute flux relative within the MRF region
|
||||
void relativeFlux(surfaceScalarField& phi) const;
|
||||
|
||||
|
||||
Reference in New Issue
Block a user