mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
BUG: Corrected sampling - handling of field regExps and field surface output
This commit is contained in:
@ -2,7 +2,7 @@
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 2011-2012 OpenFOAM Foundation
|
||||
\\ / A nd | Copyright (C) 2011-2013 OpenFOAM Foundation
|
||||
\\/ M anipulation |
|
||||
-------------------------------------------------------------------------------
|
||||
License
|
||||
@ -259,6 +259,7 @@ Foam::tmp<Foam::vectorField> Foam::sampledSurface::sample
|
||||
return tmp<vectorField>(NULL);
|
||||
}
|
||||
|
||||
|
||||
Foam::tmp<Foam::sphericalTensorField> Foam::sampledSurface::sample
|
||||
(
|
||||
const surfaceSphericalTensorField& sField
|
||||
|
||||
@ -171,17 +171,19 @@ void Foam::sampledSurfaces::write()
|
||||
writeGeometry();
|
||||
}
|
||||
|
||||
sampleAndWrite<volScalarField>();
|
||||
sampleAndWrite<volVectorField>();
|
||||
sampleAndWrite<volSphericalTensorField>();
|
||||
sampleAndWrite<volSymmTensorField>();
|
||||
sampleAndWrite<volTensorField>();
|
||||
const IOobjectList objects(mesh_, mesh_.time().timeName());
|
||||
|
||||
sampleAndWrite<surfaceScalarField>();
|
||||
sampleAndWrite<surfaceVectorField>();
|
||||
sampleAndWrite<surfaceSphericalTensorField>();
|
||||
sampleAndWrite<surfaceSymmTensorField>();
|
||||
sampleAndWrite<surfaceTensorField>();
|
||||
sampleAndWrite<volScalarField>(objects);
|
||||
sampleAndWrite<volVectorField>(objects);
|
||||
sampleAndWrite<volSphericalTensorField>(objects);
|
||||
sampleAndWrite<volSymmTensorField>(objects);
|
||||
sampleAndWrite<volTensorField>(objects);
|
||||
|
||||
sampleAndWrite<surfaceScalarField>(objects);
|
||||
sampleAndWrite<surfaceVectorField>(objects);
|
||||
sampleAndWrite<surfaceSphericalTensorField>(objects);
|
||||
sampleAndWrite<surfaceSymmTensorField>(objects);
|
||||
sampleAndWrite<surfaceTensorField>(objects);
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
@ -89,6 +89,7 @@ class sampledSurfaces
|
||||
//- Tolerance for merging points (fraction of mesh bounding box)
|
||||
static scalar mergeTol_;
|
||||
|
||||
|
||||
// Private data
|
||||
|
||||
//- Name of this set of surfaces,
|
||||
@ -113,6 +114,7 @@ class sampledSurfaces
|
||||
//- Interpolation scheme to use
|
||||
word interpolationScheme_;
|
||||
|
||||
|
||||
// surfaces
|
||||
|
||||
//- Information for merging surfaces
|
||||
@ -159,7 +161,7 @@ class sampledSurfaces
|
||||
);
|
||||
|
||||
//- Sample and write all sampled fields
|
||||
template<class Type> void sampleAndWrite();
|
||||
template<class Type> void sampleAndWrite(const IOobjectList& objects);
|
||||
|
||||
//- Disallow default bitwise copy construct and assignment
|
||||
sampledSurfaces(const sampledSurfaces&);
|
||||
@ -233,7 +235,6 @@ public:
|
||||
|
||||
//- Update for changes of mesh due to readUpdate - expires the surfaces
|
||||
virtual void readUpdate(const polyMesh::readUpdateState state);
|
||||
|
||||
};
|
||||
|
||||
|
||||
|
||||
@ -166,36 +166,26 @@ void Foam::sampledSurfaces::sampleAndWrite
|
||||
|
||||
|
||||
template<class GeoField>
|
||||
void Foam::sampledSurfaces::sampleAndWrite()
|
||||
void Foam::sampledSurfaces::sampleAndWrite(const IOobjectList& objects)
|
||||
{
|
||||
DynamicList<wordRe> fields;
|
||||
bool found = false;
|
||||
forAll (fieldSelection_, fieldI)
|
||||
if (loadFromFiles_)
|
||||
{
|
||||
const wordRe field = fieldSelection_[fieldI];
|
||||
if (mesh_.thisDb().foundObject<GeoField>(field))
|
||||
{
|
||||
found = true;
|
||||
fields.append(field);
|
||||
}
|
||||
}
|
||||
if (fields.size() && found)
|
||||
{
|
||||
forAll (fields, fieldI)
|
||||
{
|
||||
const wordRe field = fields[fieldI];
|
||||
if ((Pstream::master()) && verbose_)
|
||||
{
|
||||
Pout<< "sampleAndWrite: " << field << endl;
|
||||
}
|
||||
IOobjectList fieldObjects(objects.lookupClass(GeoField::typeName));
|
||||
|
||||
if (loadFromFiles_)
|
||||
forAll(fieldSelection_, fieldI)
|
||||
{
|
||||
const wordRe fieldNameRe = fieldSelection_[fieldI];
|
||||
IOobjectList fieldIO = fieldObjects.lookupRe(fieldNameRe);
|
||||
|
||||
forAllConstIter(IOobjectList, fieldIO, iter)
|
||||
{
|
||||
const GeoField geoField
|
||||
const word& fieldName = iter()->name();
|
||||
|
||||
const GeoField fld
|
||||
(
|
||||
IOobject
|
||||
(
|
||||
field,
|
||||
fieldName,
|
||||
mesh_.time().timeName(),
|
||||
mesh_,
|
||||
IOobject::MUST_READ
|
||||
@ -203,13 +193,38 @@ void Foam::sampledSurfaces::sampleAndWrite()
|
||||
mesh_
|
||||
);
|
||||
|
||||
sampleAndWrite(geoField);
|
||||
if ((Pstream::master()) && verbose_)
|
||||
{
|
||||
Pout<< "sampleAndWrite: " << fieldName << endl;
|
||||
}
|
||||
|
||||
sampleAndWrite(fld);
|
||||
}
|
||||
else
|
||||
}
|
||||
}
|
||||
else
|
||||
{
|
||||
forAll(fieldSelection_, fieldI)
|
||||
{
|
||||
const wordRe& fieldNameRe = fieldSelection_[fieldI];
|
||||
|
||||
const wordList dbFields
|
||||
(
|
||||
mesh_.thisDb().foundObjectRe<GeoField>(fieldNameRe)
|
||||
);
|
||||
|
||||
forAll(dbFields, i)
|
||||
{
|
||||
const word& fieldName = dbFields[i];
|
||||
|
||||
if ((Pstream::master()) && verbose_)
|
||||
{
|
||||
Pout<< "sampleAndWrite: " << fieldName << endl;
|
||||
}
|
||||
|
||||
sampleAndWrite
|
||||
(
|
||||
mesh_.thisDb().lookupObject<GeoField>(field)
|
||||
mesh_.thisDb().lookupObject<GeoField>(fieldName)
|
||||
);
|
||||
}
|
||||
}
|
||||
|
||||
Reference in New Issue
Block a user