ENH: support rhoRef for derived surfMesh sampled fields (issue #898)

This commit is contained in:
Mark Olesen
2018-06-22 16:57:30 +02:00
parent 4fe6ed3c6a
commit e0db37f043
21 changed files with 1088 additions and 70 deletions

View File

@ -49,7 +49,6 @@ Description
// Scheme to obtain node values
// (only used if interpolate=true for the surfaces below)
interpolationScheme cell;
// Output surface format
@ -74,6 +73,22 @@ Description
}
\endverbatim
Entries:
\table
Property | Description | Required | Default
type | surfaces | yes |
surfaces | the list of sample surfaces | recommended |
fields | word/regex list of fields to sampled | yes |
sampleScheme | scheme to obtain face centre value | no | cell
interpolationScheme | scheme to obtain node values | yes |
surfaceFormat | output surface format | yes |
formatOptions | dictionary of format options | no |
\endtable
Note
The interpolationScheme is only used if interpolate=true is used by any
of the surfaces.
See also
Foam::surfMeshSamplers

View File

@ -95,6 +95,93 @@ void Foam::surfMeshSamplers::checkOutNames
// return tmpRegistry().subRegistry(subName, true, false);
// }
bool Foam::surfMeshSamplers::add_rhoU(const word& derivedName)
{
const objectRegistry& db = mesh_.thisDb();
const bool existed = db.foundObject<volVectorField>(derivedName);
if (existed)
{
return false; // Volume field already existed - not added.
}
// rhoU = rho * U
const auto* rhoPtr = mesh_.lookupObjectPtr<volScalarField>("rho");
const volVectorField& U = mesh_.lookupObject<volVectorField>("U");
tmp<volVectorField> tresult;
if (rhoPtr)
{
const auto& rho = *rhoPtr;
tresult = tmp<volVectorField>::New
(
derivedName, (rho * U)
);
}
else
{
const dimensionedScalar rho("rho", dimDensity, rhoRef_);
tresult = tmp<volVectorField>::New
(
derivedName, (rho * U)
);
}
db.store(tresult.ptr());
return !existed;
}
bool Foam::surfMeshSamplers::add_pTotal(const word& derivedName)
{
const objectRegistry& db = mesh_.thisDb();
const bool existed = db.foundObject<volVectorField>(derivedName);
if (existed)
{
return false; // Volume field already existed - not added.
}
// pTotal = p + rho * U^2 / 2
const auto* rhoPtr = mesh_.lookupObjectPtr<volScalarField>("rho");
const volScalarField& p = mesh_.lookupObject<volScalarField>("p");
const volVectorField& U = mesh_.lookupObject<volVectorField>("U");
tmp<volScalarField> tresult;
if (rhoPtr)
{
const auto& rho = *rhoPtr;
tresult = tmp<volScalarField>::New
(
derivedName, (p + 0.5 * rho * magSqr(U))
);
}
else
{
const dimensionedScalar rho("rho", dimDensity, rhoRef_);
tresult = tmp<volScalarField>::New
(
derivedName, (rho * (p + 0.5 * magSqr(U)))
);
}
db.store(tresult.ptr());
return !existed;
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
@ -110,7 +197,8 @@ Foam::surfMeshSamplers::surfMeshSamplers
mesh_(refCast<const fvMesh>(obr_)),
fieldSelection_(),
derivedNames_(),
sampleScheme_(word::null)
sampleScheme_(word::null),
rhoRef_(1.0)
{
read(dict);
}
@ -128,7 +216,8 @@ Foam::surfMeshSamplers::surfMeshSamplers
mesh_(refCast<const fvMesh>(obr)),
fieldSelection_(),
derivedNames_(),
sampleScheme_(word::null)
sampleScheme_(word::null),
rhoRef_(1.0)
{
read(dict);
}
@ -157,79 +246,23 @@ bool Foam::surfMeshSamplers::execute()
for (const word& derivedName : derivedNames_)
{
// This is a fairly ugly dispatch mechanism
if (derivedName == "rhoU")
{
added.append(derivedName);
if (!db.foundObject<volVectorField>(derivedName))
if (add_rhoU(derivedName))
{
cleanup.append(derivedName);
db.store
(
new volVectorField
(
derivedName,
// rhoU = rho * U
(
mesh_.lookupObject<volScalarField>("rho")
* mesh_.lookupObject<volVectorField>("U")
)
)
);
}
added.append(derivedName);
}
else if (derivedName == "pTotal")
{
added.append(derivedName);
if (!db.foundObject<volScalarField>(derivedName))
if (add_pTotal(derivedName))
{
cleanup.append(derivedName);
const volScalarField& p =
mesh_.lookupObject<volScalarField>("p");
if (p.dimensions() == dimPressure)
{
db.store
(
new volScalarField
(
derivedName,
// pTotal = p + rho U^2 / 2
(
p
+ 0.5
* mesh_.lookupObject<volScalarField>("rho")
* magSqr
(
mesh_.lookupObject<volVectorField>("U")
)
)
)
);
}
else
{
db.store
(
new volScalarField
(
derivedName,
// pTotal = p + U^2 / 2
(
p
+ 0.5
* magSqr
(
mesh_.lookupObject<volVectorField>("U")
)
)
)
);
}
}
added.append(derivedName);
}
else
{
@ -307,6 +340,8 @@ bool Foam::surfMeshSamplers::read(const dictionary& dict)
fieldSelection_.clear();
derivedNames_.clear();
rhoRef_ = dict.lookupOrDefault<scalar>("rhoRef", 1);
const bool createOnRead = dict.lookupOrDefault("createOnRead", false);
if (dict.found("surfaces"))

View File

@ -25,7 +25,8 @@ Class
Foam::surfMeshSamplers
Description
Set of surfaces to sample into a surfMesh/surfField.
Set of surfaces to sample from a volume field onto surfField that resides
on a surfMesh object.
The execute() method is used to sample, and the write() method to write.
It is fairly common to use for sampling only and have the write disabled.
@ -53,6 +54,9 @@ Description
// Optional: pre-defined derived fields to be sampled
derived (rhoU pTotal);
// Reference density for incompressible
rhoRef 1.25;
// Optional: create surface immediately on read
// The default is to create a placeholder without any faces.
createOnRead false;
@ -69,6 +73,23 @@ Description
}
\endverbatim
Entries:
\table
Property | Description | Required | Default
type | surfMeshes | yes |
surfaces | the list of sample surfaces | recommended |
fields | word/regex list of fields to sampled | yes |
derived | additional derived fields pTotal/rhoU | no |
rhoRef | reference density for derived fields (incompressible) | no | 1
sampleScheme | scheme to obtain face centre value | no | cell
createOnRead | Create surface immediately on read | no | false;
\endtable
Note
The default is to create a placeholder surMesh without any faces on
construction. This behaviour can be changed by the createOnRead option.
For incompressible cases, \c rhoRef can be specified for use in the
derived quantities. The default is 1.
See also
Foam::sampledSurfaces
@ -93,13 +114,13 @@ SourceFiles
namespace Foam
{
// Forward declaration of classes
// Forward declarations
class Time;
class fvMesh;
class dictionary;
/*---------------------------------------------------------------------------*\
Class surfMeshSamplers Declaration
Class surfMeshSamplers Declaration
\*---------------------------------------------------------------------------*/
class surfMeshSamplers
@ -119,7 +140,7 @@ class surfMeshSamplers
const fvMesh& mesh_;
// Read from dictonary
// Read from dictionary
//- Names of fields to sample
wordRes fieldSelection_;
@ -130,6 +151,9 @@ class surfMeshSamplers
//- Sample scheme to obtain face values
word sampleScheme_;
//- Reference density (to convert from kinematic to static pressure)
scalar rhoRef_;
// Private Member Functions
@ -140,6 +164,16 @@ class surfMeshSamplers
const UList<word>& names
);
//- Hard-coded derived field (rho * U)
// \return true if field did not previously exist
bool add_rhoU(const word& derivedName);
//- Hard-coded derived field (p + 1/2 * rho * U)
// \return true if field did not previously exist
bool add_pTotal(const word& derivedName);
//- Access the sampling surfaces
inline const PtrList<surfMeshSample>& surfaces() const
{