Merge branch 'feature-curle-fo' into 'develop'

ENH: Updated Curle function object

See merge request Development/openfoam!373
This commit is contained in:
Andrew Heather
2020-06-18 21:56:16 +01:00
18 changed files with 2171 additions and 235 deletions

View File

@ -43,57 +43,12 @@ namespace functionObjects
}
}
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
bool Foam::functionObjects::Curle::calc()
{
if (foundObject<volScalarField>(fieldName_))
{
// Evaluate pressure force time derivative
const volScalarField& p = lookupObject<volScalarField>(fieldName_);
const volScalarField dpdt(scopedName("dpdt"), fvc::ddt(p));
const volScalarField::Boundary& dpdtBf = dpdt.boundaryField();
const surfaceVectorField::Boundary& SfBf = mesh_.Sf().boundaryField();
dimensionedVector dfdt("dfdt", p.dimensions()*dimArea/dimTime, Zero);
for (const label patchi : patchSet_)
{
dfdt.value() += sum(dpdtBf[patchi]*SfBf[patchi]);
}
reduce(dfdt.value(), sumOp<vector>());
// Construct and store Curle acoustic pressure
const volVectorField& C = mesh_.C();
auto tpDash = tmp<volScalarField>::New
(
IOobject
(
resultName_,
mesh_.time().timeName(),
mesh_,
IOobject::NO_READ,
IOobject::NO_WRITE
),
mesh_,
dimensionedScalar(p.dimensions(), Zero)
);
auto& pDash = tpDash.ref();
const volVectorField d(scopedName("d"), C - x0_);
pDash = (d/magSqr(d) & dfdt)/(4.0*mathematical::pi*c0_);
return store(resultName_, tpDash);
}
return false;
}
const Foam::Enum<Foam::functionObjects::Curle::modeType>
Foam::functionObjects::Curle::modeTypeNames_
({
{ modeType::point, "point" },
{ modeType::surface, "surface" },
});
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
@ -105,14 +60,17 @@ Foam::functionObjects::Curle::Curle
const dictionary& dict
)
:
fieldExpression(name, runTime, dict, "p"),
fvMeshFunctionObject(name, runTime, dict),
writeFile(mesh_, name),
pName_("p"),
patchSet_(),
x0_("x0", dimLength, Zero),
c0_("c0", dimVelocity, Zero)
observerPositions_(),
c0_(0),
rawFilePtrs_(),
inputSurface_(),
surfaceWriterPtr_(nullptr)
{
read(dict);
setResultName(typeName, fieldName_);
}
@ -120,57 +78,192 @@ Foam::functionObjects::Curle::Curle
bool Foam::functionObjects::Curle::read(const dictionary& dict)
{
if (fieldExpression::read(dict))
if (!(fvMeshFunctionObject::read(dict) && writeFile::read(dict)))
{
patchSet_ =
mesh_.boundaryMesh().patchSet
(
dict.get<wordRes>("patches")
);
if (patchSet_.empty())
{
WarningInFunction
<< "No patches defined"
<< endl;
return false;
}
// Read the reference speed of sound
dict.readEntry("c0", c0_);
if (c0_.value() < VSMALL)
{
FatalErrorInFunction
<< "Reference speed of sound = " << c0_
<< " cannot be negative or zero."
<< abort(FatalError);
}
// Set the location of the effective point source to the area-average
// of the patch face centres
const volVectorField::Boundary& Cbf = mesh_.C().boundaryField();
const surfaceScalarField::Boundary& magSfBf =
mesh_.magSf().boundaryField();
x0_.value() = vector::zero;
scalar sumMagSf = 0;
for (auto patchi : patchSet_)
{
x0_.value() += sum(Cbf[patchi]*magSfBf[patchi]);
sumMagSf += sum(magSfBf[patchi]);
}
reduce(x0_.value(), sumOp<vector>());
reduce(sumMagSf, sumOp<scalar>());
x0_.value() /= sumMagSf + ROOTVSMALL;
return true;
return false;
}
return false;
dict.readIfPresent("p", pName_);
patchSet_ = mesh_.boundaryMesh().patchSet(dict.get<wordRes>("patches"));
if (patchSet_.empty())
{
WarningInFunction
<< "No patches defined"
<< endl;
return false;
}
const modeType inputMode = modeTypeNames_.get("input", dict);
switch (inputMode)
{
case modeType::point:
{
observerPositions_ = dict.get<List<point>>("observerPositions");
break;
}
case modeType::surface:
{
const fileName fName = (dict.get<fileName>("surface")).expand();
inputSurface_ = MeshedSurface<face>(fName);
observerPositions_ = inputSurface_.Cf();
}
}
if (observerPositions_.empty())
{
WarningInFunction
<< "No observer positions defined"
<< endl;
return false;
}
const auto outputMode = modeTypeNames_.get("output", dict);
switch (outputMode)
{
case modeType::point:
{
rawFilePtrs_.setSize(observerPositions_.size());
forAll(rawFilePtrs_, filei)
{
rawFilePtrs_.set
(
filei,
createFile("observer" + Foam::name(filei))
);
if (rawFilePtrs_.set(filei))
{
OFstream& os = rawFilePtrs_[filei];
writeHeaderValue(os, "Observer", filei);
writeHeaderValue(os, "Position", observerPositions_[filei]);
writeCommented(os, "Time");
writeTabbed(os, "p(Curle)");
os << endl;
}
}
break;
}
case modeType::surface:
{
if (inputMode != modeType::surface)
{
FatalIOErrorInFunction(dict)
<< "Surface output is only available when input mode is "
<< "type '" << modeTypeNames_[modeType::surface] << "'"
<< abort(FatalIOError);
}
const word surfaceType(dict.get<word>("surfaceType"));
surfaceWriterPtr_ = surfaceWriter::New
(
surfaceType,
dict.subOrEmptyDict("formatOptions").subOrEmptyDict(surfaceType)
);
// Use outputDir/TIME/surface-name
surfaceWriterPtr_->useTimeDir() = true;
}
}
// Read the reference speed of sound
dict.readEntry("c0", c0_);
if (c0_ < VSMALL)
{
FatalErrorInFunction
<< "Reference speed of sound = " << c0_
<< " cannot be negative or zero."
<< abort(FatalError);
}
return true;
}
bool Foam::functionObjects::Curle::execute()
{
if (!foundObject<volScalarField>(pName_))
{
return false;
}
const volScalarField& p = lookupObject<volScalarField>(pName_);
const volScalarField::Boundary& pBf = p.boundaryField();
const volScalarField dpdt(scopedName("dpdt"), fvc::ddt(p));
const volScalarField::Boundary& dpdtBf = dpdt.boundaryField();
const surfaceVectorField::Boundary& SfBf = mesh_.Sf().boundaryField();
const surfaceVectorField::Boundary& CfBf = mesh_.Cf().boundaryField();
scalarField pDash(observerPositions_.size(), 0);
for (auto patchi : patchSet_)
{
const scalarField& pp = pBf[patchi];
const scalarField& dpdtp = dpdtBf[patchi];
const vectorField& Sfp = SfBf[patchi];
const vectorField& Cfp = CfBf[patchi];
forAll(observerPositions_, pointi)
{
const vectorField r(Cfp - observerPositions_[pointi]);
const scalarField invMagR(1/(mag(r) + ROOTVSMALL));
pDash[pointi] +=
sum((pp*sqr(invMagR) + invMagR/c0_*dpdtp)*(Sfp & r));
}
}
Pstream::listCombineGather(pDash, plusEqOp<scalar>());
Pstream::listCombineScatter(pDash);
if (surfaceWriterPtr_)
{
if (Pstream::master())
{
// Time-aware, with time spliced into the output path
surfaceWriterPtr_->beginTime(time_);
surfaceWriterPtr_->open
(
inputSurface_.points(),
inputSurface_.surfFaces(),
(baseFileDir()/name()/"surface"),
false // serial - already merged
);
surfaceWriterPtr_->write("p(Curle)", pDash);
surfaceWriterPtr_->endTime();
surfaceWriterPtr_->clear();
}
}
else
{
forAll(observerPositions_, pointi)
{
if (rawFilePtrs_.set(pointi))
{
OFstream& os = rawFilePtrs_[pointi];
writeCurrentTime(os);
os << pDash[pointi] << endl;
}
}
}
return true;
}
bool Foam::functionObjects::Curle::write()
{
return true;
}

View File

@ -35,15 +35,21 @@ Description
Curle's analogy is implemented as:
\f[
p' = \frac{1}{4 \pi c_0}\frac{\vec d}{|\vec d|^2}\cdot\frac{d\vec F}{dt}
p' = \frac{1}{4 \pi}
\frac{\vec{r}}{| \vec{r} | ^2} \cdot
\left(
\frac{\vec{F}}{| \vec{r} |}
+ \frac{1}{c_0}\frac{d(\vec{F})}{d(t)}
\right)
\f]
where
\vartable
p' | Curle's acoustic pressure [Pa] or [Pa \f$(m^3/\rho)\f$]
c_0 | Reference speed of sound [m/s]
\vec d | Distance vector to observer locations [m]
\vec F | Force [N] or [N (\f$m^3/\rho\f$)]
p' | Curle's acoustic pressure [Pa] or [Pa (m3/kg)]
c_0 | Reference speed of sound [m/s]
\vec{r} | Distance vector to observer locations [m]
\vec{F} | Force [N] or [N (m3/kg)]
t | time [s]
\endvartable
Operands:
@ -62,16 +68,27 @@ Usage
\verbatim
Curle1
{
// Mandatory entries (unmodifiable)
type Curle;
libs (fieldFunctionObjects);
// Mandatory entries (runtime modifiable)
patches (<patch1> <patch2> ... <patchN>)
c0 343;
// Optional (inherited) entries
libs ("libfieldFunctionObjects.so");
...
patches (surface1 surface2);
c0 330;
// Input - either points or surface
input points;
observerPositions ((0 0 0)(1 0 0));
//input surface;
//surface "inputSurface.obj"
// Output - either points or surface
output points;
//output surface;
//surfaceType ensight;
}
\endverbatim
@ -80,8 +97,14 @@ Usage
Property | Description | Type | Req'd | Dflt
type | Type name: Curle | word | yes | -
libs | Library name: fieldFunctionObjects | word | yes | -
patches | Names of the operand patches | wordList | yes | -
c0 | Reference speed of sound [m/s] | scalar | yes | -
p | Pressure field name | word | no | p
patches | Sound generation patch names | wordList | yes | -
c0 | Reference speed of sound | scalar | yes | -
input | Input type | word | yes | -
observerPositions | List of observer positions (x y z) | pointList | no |-
surface | Input surface file name | word | no | -
output | Output type | word | yes | -
surfaceType | Output surface type | word | no | -
\endtable
The inherited entries are elaborated in:
@ -93,7 +116,7 @@ Usage
See also
- Foam::functionObject
- Foam::functionObjects::fvMeshFunctionObject
- Foam::functionObjects::fieldExpression
- Foam::functionObjects::writeFile
- ExtendedCodeGuide::functionObjects::field::Curle
SourceFiles
@ -104,9 +127,11 @@ SourceFiles
#ifndef functionObjects_Curle_H
#define functionObjects_Curle_H
#include "fieldExpression.H"
#include "dimensionedScalar.H"
#include "dimensionedVector.H"
#include "fvMeshFunctionObject.H"
#include "writeFile.H"
#include "Enum.H"
#include "MeshedSurface.H"
#include "surfaceWriter.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
@ -121,26 +146,40 @@ namespace functionObjects
class Curle
:
public fieldExpression
public fvMeshFunctionObject,
public writeFile
{
enum class modeType
{
point,
surface
};
static const Enum<modeType> modeTypeNames_;
// Private Data
// Read from dictionary
//- Name of pressure field; default = p
word pName_;
//- Patches to integrate forces over
labelHashSet patchSet_;
//- Patches to integrate forces over
labelHashSet patchSet_;
//- Area-averaged centre of patch faces
dimensionedVector x0_;
//- Observer positions
List<point> observerPositions_;
//- Reference speed of sound
dimensionedScalar c0_;
//- Reference speed of sound
scalar c0_;
//- Output files when sampling points
PtrList<OFstream> rawFilePtrs_;
protected:
//- Input surface when sampling a surface
MeshedSurface<face> inputSurface_;
//- Calculate the acoustic pressure field and return true if successful
virtual bool calc();
//- Ouput surface when sampling a surface
autoPtr<surfaceWriter> surfaceWriterPtr_;
public:
@ -174,6 +213,10 @@ public:
//- Read the Curle data
virtual bool read(const dictionary&);
virtual bool execute();
virtual bool write();
};

View File

@ -104,7 +104,8 @@ void surfaceNoise::initialise(const fileName& fName)
deltaT_ = checkUniformTimeStep(times_);
// Read the surface geometry
const meshedSurface& surf = readerPtr_->geometry();
// Note: hard-coded to read mesh from first time index
const meshedSurface& surf = readerPtr_->geometry(0);
nFace_ = surf.size();
}
@ -208,7 +209,6 @@ void surfaceNoise::readSurfaceData
Info<< " time: " << times_[i] << endl;
const scalarField p(readerPtr_->field(timeI, pIndex_, scalar(0)));
forAll(p, faceI)
{
pData[faceI][i] = p[faceI]*rhoRef_;
@ -256,7 +256,8 @@ Foam::scalar surfaceNoise::writeSurfaceData
scalar areaAverage = 0;
if (Pstream::master())
{
const meshedSurface& surf = readerPtr_->geometry();
// Note: hard-coded to read mesh from first time index
const meshedSurface& surf = readerPtr_->geometry(0);
scalarField allData(surf.size());
@ -298,7 +299,7 @@ Foam::scalar surfaceNoise::writeSurfaceData
// TO BE VERIFIED: area-averaged values
// areaAverage = sum(allData*surf.magSf())/sum(surf.magSf());
areaAverage = sum(allData)/allData.size();
areaAverage = sum(allData)/(allData.size() + ROOTVSMALL);
}
Pstream::scatter(areaAverage);
@ -306,7 +307,8 @@ Foam::scalar surfaceNoise::writeSurfaceData
}
else
{
const meshedSurface& surf = readerPtr_->geometry();
// Note: hard-coded to read mesh from first time index
const meshedSurface& surf = readerPtr_->geometry(0);
if (writeSurface)
{
@ -329,7 +331,7 @@ Foam::scalar surfaceNoise::writeSurfaceData
// TO BE VERIFIED: area-averaged values
// return sum(data*surf.magSf())/sum(surf.magSf());
return sum(data)/data.size();
return sum(data)/(data.size() + ROOTVSMALL);
}
}
@ -357,7 +359,8 @@ Foam::scalar surfaceNoise::surfaceAverage
scalar areaAverage = 0;
if (Pstream::master())
{
const meshedSurface& surf = readerPtr_->geometry();
// Note: hard-coded to read mesh from first time index
const meshedSurface& surf = readerPtr_->geometry(0);
scalarField allData(surf.size());

View File

@ -123,6 +123,25 @@ void Foam::ensightSurfaceReader::debugSection
}
Foam::fileName Foam::ensightSurfaceReader::replaceMask
(
const fileName& fName,
const label timeIndex
) const
{
fileName result(fName);
std::ostringstream oss;
label nMask = stringOps::count(fName, '*');
const std::string maskStr(nMask, '*');
oss << std::setfill('0') << std::setw(nMask) << timeIndex;
const word indexStr = oss.str();
result.replace(maskStr, indexStr);
return result;
}
Foam::Pair<Foam::ensightSurfaceReader::idTypes>
Foam::ensightSurfaceReader::readGeometryHeader(ensightReadFile& is) const
{
@ -230,6 +249,8 @@ void Foam::ensightSurfaceReader::readCase(IFstream& is)
// - use the last entry
meshFileName_ = stringOps::splitSpace(buffer).last().str();
DebugInfo << "mesh file:" << meshFileName_ << endl;
debugSection("VARIABLE", is);
// Read the field description
@ -330,13 +351,17 @@ Foam::ensightSurfaceReader::ensightSurfaceReader(const fileName& fName)
// * * * * * * * * * * * * * Public Member Functions * * * * * * * * * * * //
const Foam::meshedSurface& Foam::ensightSurfaceReader::geometry()
const Foam::meshedSurface& Foam::ensightSurfaceReader::geometry
(
const label timeIndex
)
{
DebugInFunction << endl;
if (!surfPtr_.valid())
{
IFstream isBinary(baseDir_/meshFileName_, IOstream::BINARY);
fileName meshInstance(replaceMask(meshFileName_, timeIndex));
IFstream isBinary(baseDir_/meshInstance, IOstream::BINARY);
if (!isBinary.good())
{
@ -391,7 +416,7 @@ const Foam::meshedSurface& Foam::ensightSurfaceReader::geometry()
}
ensightReadFile is(baseDir_/meshFileName_, streamFormat_);
ensightReadFile is(baseDir_/meshInstance, streamFormat_);
DebugInfo
<< "File: " << is.name() << nl;

View File

@ -5,7 +5,7 @@
\\ / A nd | www.openfoam.com
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2015-2019 OpenCFD Ltd.
Copyright (C) 2015-2020 OpenCFD Ltd.
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
@ -112,6 +112,13 @@ protected:
//- Read and check a section header
void debugSection(const word& expected, IFstream& is) const;
//- Replace the mask with an 0 padded string
fileName replaceMask
(
const fileName& fName,
const label timeIndex
) const;
//- Read (and discard) geometry file header.
// \return information about node/element id handling
Pair<idTypes> readGeometryHeader(ensightReadFile& is) const;
@ -164,7 +171,7 @@ public:
// Member Functions
//- Return a reference to the surface geometry
virtual const meshedSurface& geometry();
virtual const meshedSurface& geometry(const label timeIndex);
//- Return a list of the available times
virtual instantList times() const;

View File

@ -70,24 +70,7 @@ Foam::tmp<Foam::Field<Type>> Foam::ensightSurfaceReader::readField
const word& fieldName(fieldNames_[fieldIndex]);
const label fileIndex = timeStartIndex_ + timeIndex*timeIncrement_;
fileName fieldFileName(fieldFileNames_[fieldIndex]);
std::ostringstream oss;
label nMask = 0;
for (size_t chari = 0; chari < fieldFileName.size(); ++chari)
{
if (fieldFileName[chari] == '*')
{
nMask++;
}
}
const std::string maskStr(nMask, '*');
oss << std::setfill('0') << std::setw(nMask) << fileIndex;
const word indexStr = oss.str();
fieldFileName.replace(maskStr, indexStr);
fileName fieldFileName(replaceMask(fieldFileNames_[fieldIndex], fileIndex));
ensightReadFile is(baseDir_/fieldFileName, streamFormat_);
if (!is.good())

View File

@ -101,7 +101,7 @@ public:
// Member Functions
//- Return a reference to the surface geometry
virtual const meshedSurface& geometry() = 0;
virtual const meshedSurface& geometry(const label timeIndex) = 0;
//- Return a list of the available times
virtual instantList times() const = 0;

View File

@ -16,19 +16,20 @@ FoamFile
dimensions [0 1 -1 0 0 0 0];
internalField uniform (40 0 0);
internalField uniform (0.025 0 0);
boundaryField
{
cylinder
{
type noSlip;
type fixedValue;
value uniform (0 0 0);
}
inlet
{
type fixedValue;
value $internalField;
value uniform (0.025 0 0);
}
outlet
@ -38,10 +39,14 @@ boundaryField
inletValue uniform (0 0 0);
}
"top|bottom"
"(top|bottom)"
{
type pressureInletOutletVelocity;
value uniform (0 0 0);
type slip;
}
frontAndBack
{
type empty;
}
}

View File

@ -17,7 +17,7 @@ FoamFile
dimensions [0 2 -1 0 0 0 0];
internalField uniform 1.51e-05;
internalField uniform 1e-05;
boundaryField
{
@ -33,12 +33,22 @@ boundaryField
value $internalField;
}
"outlet|top|bottom"
outlet
{
type inletOutlet;
inletValue $internalField;
value $internalField;
}
"(top|bottom)"
{
type zeroGradient;
}
frontAndBack
{
type empty;
}
}

View File

@ -33,11 +33,16 @@ boundaryField
value $internalField;
}
"outlet|top|bottom"
"(outlet|top|bottom)"
{
type calculated;
value $internalField;
}
frontAndBack
{
type empty;
}
}

View File

@ -30,15 +30,25 @@ boundaryField
type zeroGradient;
}
"outlet|top|bottom"
outlet
{
type totalPressure;
type fixedValue; // totalPressure;
rho none;
psi none;
p0 $internalField;
gamma 1;
value $internalField;
}
"(top|bottom)"
{
type zeroGradient;
}
frontAndBack
{
type empty;
}
}

View File

@ -1,3 +1,15 @@
- Acoustic pressure is calculated based on Curle's analogy.
- Such field is sampled for point (probe) and surface (cutting plane)
- noise utility is operated onto these samples to obtain frequency information
Summary
- Vortex street generated from flow over a sphere
- Acoustic pressure is calculated based on Curle's analogy
- noise utility is processes the acoustic pressure to obtain frequency information
Case detailis
- cylinder diameter: 0.04
- flow speed: 0.025
- kinematic viscosity: 1e-5
- Reynolds number based on cylinder diameter, Red = 100
- For this Red, St = 0.198(1 - 0.197/Red) = 0.159
- Shedding fequency, f = St*U/d = 0.1 Hz
Plotting the FFT of the pressure field shoud show a peak at the shedding
frequency

View File

@ -17,6 +17,6 @@ FoamFile
transportModel Newtonian;
nu 2e-05;
nu 1e-05;
// ************************************************************************* //

File diff suppressed because it is too large Load Diff

View File

@ -23,15 +23,15 @@ startTime 0;
stopAt endTime;
endTime 0.2;
endTime 200;
deltaT 1e-5;
deltaT 1e-2;
writeControl timeStep;
writeInterval 100; // every 0.001s
writeInterval 100; // every 1s
purgeWrite 200;
purgeWrite 50;
writeFormat binary;
@ -55,46 +55,81 @@ functions
fields (U p);
}
curle
curleSurface
{
// Mandatory entries
type Curle;
libs (fieldFunctionObjects);
patches (cylinder);
c0 343;
input surface;
surface surface.obj;
output surface;
surfaceType ensight;
formatOptions
{
ensight
{
format ascii; // needed for surfaceNoise ...
collateTimes true;
}
}
// Optional (inherited) entries
field p;
result Curle;
p p;
region region0;
enabled true;
log true;
timeStart 0.1;
timeStart 100;
timeEnd 1000;
executeControl timeStep;
executeInterval 1;
}
curlePoint
{
// Mandatory entries
type Curle;
libs (fieldFunctionObjects);
patches (cylinder);
c0 343;
input point;
output point;
observerPositions
(
(0.20 0.17 -0.01) // N
(0.22 0.15 -0.01) // E
(0.20 0.13 -0.01) // S
(0.18 0.15 -0.01) // W
);
// Optional (inherited) entries
p p;
region region0;
enabled true;
log true;
timeStart 100;
timeEnd 1000;
executeControl timeStep;
executeInterval 1;
writeControl writeTime;
writeInterval -1;
}
cuttingPlane
{
type surfaces;
libs (sampling);
writeControl timeStep;
writeInterval 2;
timeStart 0.1;
writeControl writeTime;
timeStart 100;
surfaceFormat ensight;
formatOptions
{
ensight
{
format binary;
collateTimes true;
format binary;
collateTimes true;
}
}
fields (p Curle);
fields (p U);
interpolationScheme cellPoint;
@ -106,8 +141,8 @@ functions
planeType pointAndNormal;
pointAndNormalDict
{
point (0 0 -0.01);
normal (0 0 1);
point (0 0 -0.01);
normal (0 0 1);
}
interpolate false;
}
@ -119,32 +154,14 @@ functions
type forces;
libs (forces);
writeControl writeTime;
timeStart 0.1;
timeStart 100;
patches (cylinder);
CofR (0.20 0.15 -0.01);
writeFields yes;
rho rhoInf;
rhoInf 1.205;
}
probes
{
type patchProbes;
libs (sampling);
writeControl timeStep;
timeStart 0.1;
patch cylinder;
probeLocations
(
(0.20 0.17 -0.01) // N
(0.22 0.15 -0.01) // E
(0.20 0.13 -0.01) // S
(0.18 0.15 -0.01) // W
);
fields (p forces:force Curle);
rhoInf 1;
}
}

View File

@ -17,13 +17,10 @@ FoamFile
solvers
{
p
{
solver GAMG;
smoother DICGaussSeidel;
// solver PCG;
// preconditioner DIC;
tolerance 1e-20;
relTol 0.05;
}
@ -37,10 +34,10 @@ solvers
"(U|nuTilda)"
{
// solver PBiCGStab;
// preconditioner DILU;
solver smoothSolver;
smoother symGaussSeidel;
//solver smoothSolver;
//smoother symGaussSeidel;
solver PBiCGStab;
preconditioner DILU;
tolerance 0;
relTol 0.1;
}
@ -64,9 +61,9 @@ PIMPLE
relaxationFactors
{
nuTilda 0.95;
U 0.95;
p 0.95;
nuTilda 0.8;
U 0.8;
p 0.8;
".*Final" 1.0;
}

View File

@ -34,8 +34,10 @@ pointNoiseCoeffs
files
(
"postProcessing/probes/0/p"
"postProcessing/probes/0.1/Curle"
"postProcessing/curlePoint/0/observer0.dat"
"postProcessing/curlePoint/0/observer1.dat"
"postProcessing/curlePoint/0/observer2.dat"
"postProcessing/curlePoint/0/observer3.dat"
);
nHeaderLine 6;
@ -48,7 +50,7 @@ pointNoiseCoeffs
// default = 2^16 (=65536)
N 4096;
rhoRef 1.205;
rhoRef 1;
}

View File

@ -34,7 +34,7 @@ surfaceNoiseCoeffs
// Input file
file "postProcessing/cuttingPlane/zNormal/zNormal.case";
file "postProcessing/curleSurface/surface/surface.case";
// Surface reader
reader ensight;
@ -53,13 +53,13 @@ surfaceNoiseCoeffs
}
// Pressure field
p Curle;
p pCurle;
// Number of samples in sampling window
// default = 2^16 (=65536)
N 1024;
N 4096;
rhoRef 1.205;
rhoRef 1;
}