mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-12-28 03:37:59 +00:00
- has the selected values directly and use these lookup names to store directly into a hash. This replaces several parallel lists of decomp information etc and makes it easier.
884 lines
23 KiB
C
884 lines
23 KiB
C
/*---------------------------------------------------------------------------*\
|
|
========= |
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
\\ / O peration |
|
|
\\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
|
|
\\/ M anipulation | Copyright (C) 2017 OpenCFD Ltd.
|
|
-------------------------------------------------------------------------------
|
|
License
|
|
This file is part of OpenFOAM.
|
|
|
|
OpenFOAM is free software: you can redistribute it and/or modify it
|
|
under the terms of the GNU General Public License as published by
|
|
the Free Software Foundation, either version 3 of the License, or
|
|
(at your option) any later version.
|
|
|
|
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
|
|
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
|
|
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
|
|
for more details.
|
|
|
|
You should have received a copy of the GNU General Public License
|
|
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
#include "vtkPVFoam.H"
|
|
#include "vtkPVFoamReader.h"
|
|
|
|
// OpenFOAM includes
|
|
#include "fvMesh.H"
|
|
#include "Time.H"
|
|
#include "patchZones.H"
|
|
|
|
// VTK includes
|
|
#include "vtkDataArraySelection.h"
|
|
#include "vtkMultiBlockDataSet.h"
|
|
#include "vtkRenderer.h"
|
|
#include "vtkTextActor.h"
|
|
#include "vtkTextProperty.h"
|
|
#include "vtkSmartPointer.h"
|
|
|
|
// Templates (only needed here)
|
|
#include "vtkPVFoamUpdateTemplates.C"
|
|
|
|
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
|
|
|
|
namespace Foam
|
|
{
|
|
defineTypeNameAndDebug(vtkPVFoam, 0);
|
|
|
|
// file-scope
|
|
static word updateStateName(polyMesh::readUpdateState state)
|
|
{
|
|
switch (state)
|
|
{
|
|
case polyMesh::UNCHANGED: return "UNCHANGED";
|
|
case polyMesh::POINTS_MOVED: return "POINTS_MOVED";
|
|
case polyMesh::TOPO_CHANGE: return "TOPO_CHANGE";
|
|
case polyMesh::TOPO_PATCH_CHANGE: return "TOPO_PATCH_CHANGE";
|
|
};
|
|
|
|
return "UNKNOWN";
|
|
}
|
|
}
|
|
|
|
|
|
// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
|
|
|
|
namespace Foam
|
|
{
|
|
// file-scope
|
|
|
|
//- Create a text actor
|
|
vtkSmartPointer<vtkTextActor> createTextActor
|
|
(
|
|
const std::string& s,
|
|
const Foam::point& pt
|
|
)
|
|
{
|
|
vtkSmartPointer<vtkTextActor> txt =
|
|
vtkSmartPointer<vtkTextActor>::New();
|
|
|
|
txt->SetInput(s.c_str());
|
|
|
|
// Set text properties
|
|
vtkTextProperty* tprop = txt->GetTextProperty();
|
|
tprop->SetFontFamilyToArial();
|
|
tprop->BoldOn();
|
|
tprop->ShadowOff();
|
|
tprop->SetLineSpacing(1.0);
|
|
tprop->SetFontSize(14);
|
|
tprop->SetColor(1.0, 0.0, 1.0);
|
|
tprop->SetJustificationToCentered();
|
|
|
|
txt->GetPositionCoordinate()->SetCoordinateSystemToWorld();
|
|
txt->GetPositionCoordinate()->SetValue(pt.x(), pt.y(), pt.z());
|
|
|
|
return txt;
|
|
}
|
|
}
|
|
|
|
|
|
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
|
|
|
|
void Foam::vtkPVFoam::resetCounters()
|
|
{
|
|
// Reset array range information (ids and sizes)
|
|
rangeVolume_.reset();
|
|
rangePatches_.reset();
|
|
rangeLagrangian_.reset();
|
|
rangeCellZones_.reset();
|
|
rangeFaceZones_.reset();
|
|
rangePointZones_.reset();
|
|
rangeCellSets_.reset();
|
|
rangeFaceSets_.reset();
|
|
rangePointSets_.reset();
|
|
}
|
|
|
|
|
|
int Foam::vtkPVFoam::setTime(int nRequest, const double requestTimes[])
|
|
{
|
|
Time& runTime = dbPtr_();
|
|
|
|
// Get times list
|
|
instantList Times = runTime.times();
|
|
|
|
int nearestIndex = timeIndex_;
|
|
for (int requestI = 0; requestI < nRequest; ++requestI)
|
|
{
|
|
int index = Time::findClosestTimeIndex(Times, requestTimes[requestI]);
|
|
if (index >= 0 && index != timeIndex_)
|
|
{
|
|
nearestIndex = index;
|
|
break;
|
|
}
|
|
}
|
|
|
|
if (nearestIndex < 0)
|
|
{
|
|
nearestIndex = 0;
|
|
}
|
|
|
|
if (debug)
|
|
{
|
|
Info<< "<beg> setTime(";
|
|
for (int requestI = 0; requestI < nRequest; ++requestI)
|
|
{
|
|
if (requestI) Info<< ", ";
|
|
Info<< requestTimes[requestI];
|
|
}
|
|
Info<< ") - previousIndex = " << timeIndex_
|
|
<< ", nearestIndex = " << nearestIndex << endl;
|
|
}
|
|
|
|
// See what has changed
|
|
if (timeIndex_ != nearestIndex)
|
|
{
|
|
timeIndex_ = nearestIndex;
|
|
runTime.setTime(Times[nearestIndex], nearestIndex);
|
|
|
|
// When the changes, so do the fields
|
|
fieldsChanged_ = true;
|
|
meshState_ = meshPtr_ ? meshPtr_->readUpdate() : polyMesh::TOPO_CHANGE;
|
|
|
|
reader_->UpdateProgress(0.05);
|
|
|
|
// this seems to be needed for catching Lagrangian fields
|
|
updateInfo();
|
|
}
|
|
|
|
if (debug)
|
|
{
|
|
Info<< "<end> setTime() - selectedTime="
|
|
<< Times[nearestIndex].name() << " index=" << timeIndex_
|
|
<< "/" << Times.size()
|
|
<< " meshUpdateState=" << updateStateName(meshState_)
|
|
<< " fieldsChanged=" << Switch(fieldsChanged_) << endl;
|
|
}
|
|
|
|
return nearestIndex;
|
|
}
|
|
|
|
|
|
Foam::word Foam::vtkPVFoam::getPartName(const int partId) const
|
|
{
|
|
if (selectedPartIds_.found(partId))
|
|
{
|
|
return getFoamName(selectedPartIds_[partId]);
|
|
}
|
|
else
|
|
{
|
|
return word::null;
|
|
}
|
|
}
|
|
|
|
|
|
Foam::word Foam::vtkPVFoam::getReaderPartName(const int partId) const
|
|
{
|
|
return getFoamName(reader_->GetPartArrayName(partId));
|
|
}
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
|
|
|
Foam::vtkPVFoam::vtkPVFoam
|
|
(
|
|
const char* const FileName,
|
|
vtkPVFoamReader* reader
|
|
)
|
|
:
|
|
reader_(reader),
|
|
dbPtr_(nullptr),
|
|
meshPtr_(nullptr),
|
|
meshRegion_(polyMesh::defaultRegion),
|
|
meshDir_(polyMesh::meshSubDir),
|
|
timeIndex_(-1),
|
|
meshState_(polyMesh::TOPO_CHANGE),
|
|
fieldsChanged_(true),
|
|
rangeVolume_("unzoned"),
|
|
rangePatches_("patch"),
|
|
rangeLagrangian_("lagrangian"),
|
|
rangeCellZones_("cellZone"),
|
|
rangeFaceZones_("faceZone"),
|
|
rangePointZones_("pointZone"),
|
|
rangeCellSets_("cellSet"),
|
|
rangeFaceSets_("faceSet"),
|
|
rangePointSets_("pointSet")
|
|
{
|
|
if (debug)
|
|
{
|
|
Info<< "vtkPVFoam - " << FileName << endl;
|
|
printMemory();
|
|
}
|
|
|
|
// avoid argList and get rootPath/caseName directly from the file
|
|
fileName fullCasePath(fileName(FileName).path());
|
|
|
|
if (!isDir(fullCasePath))
|
|
{
|
|
return;
|
|
}
|
|
if (fullCasePath == ".")
|
|
{
|
|
fullCasePath = cwd();
|
|
}
|
|
|
|
// The name of the executable, unless already present in the environment
|
|
setEnv("FOAM_EXECUTABLE", "paraview", false);
|
|
|
|
// Set the case as an environment variable - some BCs might use this
|
|
if (fullCasePath.name().find("processor", 0) == 0)
|
|
{
|
|
const fileName globalCase = fullCasePath.path();
|
|
|
|
setEnv("FOAM_CASE", globalCase, true);
|
|
setEnv("FOAM_CASENAME", globalCase.name(), true);
|
|
}
|
|
else
|
|
{
|
|
setEnv("FOAM_CASE", fullCasePath, true);
|
|
setEnv("FOAM_CASENAME", fullCasePath.name(), true);
|
|
}
|
|
|
|
// look for 'case{region}.OpenFOAM'
|
|
// could be stringent and insist the prefix match the directory name...
|
|
// Note: cannot use fileName::name() due to the embedded '{}'
|
|
string caseName(fileName(FileName).lessExt());
|
|
string::size_type beg = caseName.find_last_of("/{");
|
|
string::size_type end = caseName.find('}', beg);
|
|
|
|
if
|
|
(
|
|
beg != string::npos && caseName[beg] == '{'
|
|
&& end != string::npos && end == caseName.size()-1
|
|
)
|
|
{
|
|
meshRegion_ = caseName.substr(beg+1, end-beg-1);
|
|
|
|
// some safety
|
|
if (meshRegion_.empty())
|
|
{
|
|
meshRegion_ = polyMesh::defaultRegion;
|
|
}
|
|
|
|
if (meshRegion_ != polyMesh::defaultRegion)
|
|
{
|
|
meshDir_ = meshRegion_/polyMesh::meshSubDir;
|
|
}
|
|
}
|
|
|
|
if (debug)
|
|
{
|
|
Info<< "fullCasePath=" << fullCasePath << nl
|
|
<< "FOAM_CASE=" << getEnv("FOAM_CASE") << nl
|
|
<< "FOAM_CASENAME=" << getEnv("FOAM_CASENAME") << nl
|
|
<< "region=" << meshRegion_ << endl;
|
|
}
|
|
|
|
// Create time object
|
|
dbPtr_.reset
|
|
(
|
|
new Time
|
|
(
|
|
Time::controlDictName,
|
|
fileName(fullCasePath.path()),
|
|
fileName(fullCasePath.name())
|
|
)
|
|
);
|
|
|
|
dbPtr_().functionObjects().off();
|
|
|
|
updateInfo();
|
|
}
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
|
|
|
|
Foam::vtkPVFoam::~vtkPVFoam()
|
|
{
|
|
if (debug)
|
|
{
|
|
Info<< "~vtkPVFoam" << endl;
|
|
}
|
|
|
|
delete meshPtr_;
|
|
}
|
|
|
|
|
|
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
|
|
|
void Foam::vtkPVFoam::updateInfo()
|
|
{
|
|
if (debug)
|
|
{
|
|
Info<< "<beg> updateInfo"
|
|
<< " [meshPtr=" << (meshPtr_ ? "set" : "nullptr") << "] timeIndex="
|
|
<< timeIndex_ << endl;
|
|
}
|
|
|
|
resetCounters();
|
|
|
|
vtkDataArraySelection* partSelection = reader_->GetPartSelection();
|
|
|
|
// there are two ways to ensure we have the correct list of parts:
|
|
// 1. remove everything and then set particular entries 'on'
|
|
// 2. build a 'char **' list and call SetArraysWithDefault()
|
|
//
|
|
// Nr. 2 has the potential advantage of not touching the modification
|
|
// time of the vtkDataArraySelection, but the qt/paraview proxy
|
|
// layer doesn't care about that anyhow.
|
|
|
|
HashSet<string> enabledEntries;
|
|
if (!partSelection->GetNumberOfArrays() && !meshPtr_)
|
|
{
|
|
// enable 'internalMesh' on the first call
|
|
enabledEntries = { "internalMesh" };
|
|
}
|
|
else
|
|
{
|
|
// preserve the enabled selections
|
|
enabledEntries = getSelectedArrayEntries(partSelection);
|
|
}
|
|
|
|
// Clear current mesh parts list
|
|
partSelection->RemoveAllArrays();
|
|
|
|
// Update mesh parts list - add Lagrangian at the bottom
|
|
updateInfoInternalMesh(partSelection);
|
|
updateInfoPatches(partSelection, enabledEntries);
|
|
updateInfoSets(partSelection);
|
|
updateInfoZones(partSelection);
|
|
updateInfoLagrangian(partSelection);
|
|
|
|
// Adjust/restore the enabled selections
|
|
setSelectedArrayEntries(partSelection, enabledEntries);
|
|
|
|
if (meshState_ != polyMesh::UNCHANGED)
|
|
{
|
|
fieldsChanged_ = true;
|
|
}
|
|
|
|
// Update volume, point and lagrangian fields
|
|
updateInfoFields<fvPatchField, volMesh>
|
|
(
|
|
reader_->GetVolFieldSelection()
|
|
);
|
|
updateInfoFields<pointPatchField, pointMesh>
|
|
(
|
|
reader_->GetPointFieldSelection()
|
|
);
|
|
updateInfoLagrangianFields
|
|
(
|
|
reader_->GetLagrangianFieldSelection()
|
|
);
|
|
|
|
if (debug)
|
|
{
|
|
Info<< "<end> updateInfo" << endl;
|
|
}
|
|
}
|
|
|
|
|
|
void Foam::vtkPVFoam::Update
|
|
(
|
|
vtkMultiBlockDataSet* output,
|
|
vtkMultiBlockDataSet* outputLagrangian
|
|
)
|
|
{
|
|
if (debug)
|
|
{
|
|
cout<< "<beg> Foam::vtkPVFoam::Update\n";
|
|
output->Print(cout);
|
|
if (outputLagrangian) outputLagrangian->Print(cout);
|
|
printMemory();
|
|
}
|
|
reader_->UpdateProgress(0.1);
|
|
|
|
// Set up mesh parts selection(s)
|
|
{
|
|
if (debug)
|
|
{
|
|
Info<< "<beg> updateMeshPartsStatus" << endl;
|
|
}
|
|
|
|
vtkDataArraySelection* selection = reader_->GetPartSelection();
|
|
const int n = selection->GetNumberOfArrays();
|
|
|
|
// All previously selected (enabled) names
|
|
HashSet<string> original;
|
|
forAllConstIters(selectedPartIds_, iter)
|
|
{
|
|
original.insert(iter.object());
|
|
}
|
|
|
|
selectedPartIds_.clear();
|
|
|
|
for (int id=0; id < n; ++id)
|
|
{
|
|
string str(selection->GetArrayName(id));
|
|
bool status = selection->GetArraySetting(id);
|
|
|
|
if (status)
|
|
{
|
|
selectedPartIds_.set(id, str); // id -> name
|
|
|
|
if (!original.erase(str))
|
|
{
|
|
// New part, or newly enabled
|
|
//? meshChanged_ = true;
|
|
}
|
|
}
|
|
else
|
|
{
|
|
if (original.erase(str))
|
|
{
|
|
// Part disappeared, or newly disabled
|
|
//? meshChanged_ = true;
|
|
}
|
|
}
|
|
|
|
if (debug > 1)
|
|
{
|
|
Info<< " part[" << id << "] = " << status
|
|
<< " : " << str << nl;
|
|
}
|
|
}
|
|
|
|
if (debug)
|
|
{
|
|
Info<< "<end> updateMeshPartsStatus" << endl;
|
|
}
|
|
}
|
|
|
|
reader_->UpdateProgress(0.15);
|
|
|
|
// Update the OpenFOAM mesh
|
|
{
|
|
if (debug)
|
|
{
|
|
Info<< "<beg> updateFoamMesh" << endl;
|
|
printMemory();
|
|
}
|
|
|
|
if (!reader_->GetCacheMesh())
|
|
{
|
|
delete meshPtr_;
|
|
meshPtr_ = nullptr;
|
|
}
|
|
|
|
// Check to see if the OpenFOAM mesh has been created
|
|
if (!meshPtr_)
|
|
{
|
|
if (debug)
|
|
{
|
|
Info<< "Creating OpenFOAM mesh for region " << meshRegion_
|
|
<< " at time=" << dbPtr_().timeName() << endl;
|
|
|
|
}
|
|
|
|
meshPtr_ = new fvMesh
|
|
(
|
|
IOobject
|
|
(
|
|
meshRegion_,
|
|
dbPtr_().timeName(),
|
|
dbPtr_(),
|
|
IOobject::MUST_READ
|
|
)
|
|
);
|
|
|
|
meshState_ = polyMesh::TOPO_CHANGE; // New mesh
|
|
}
|
|
else
|
|
{
|
|
if (debug)
|
|
{
|
|
Info<< "Using existing OpenFOAM mesh" << endl;
|
|
}
|
|
}
|
|
|
|
if (debug)
|
|
{
|
|
Info<< "<end> updateFoamMesh" << endl;
|
|
printMemory();
|
|
}
|
|
}
|
|
|
|
reader_->UpdateProgress(0.4);
|
|
|
|
// Update cached, saved, unneed values:
|
|
HashSet<string> nowActive;
|
|
forAllConstIters(selectedPartIds_, iter)
|
|
{
|
|
nowActive.insert(iter.object());
|
|
}
|
|
|
|
// Dispose of unneeded components
|
|
cachedVtp_.retain(nowActive);
|
|
cachedVtu_.retain(nowActive);
|
|
|
|
// Reset (expire) dataset ids
|
|
partDataset_.clear();
|
|
|
|
// Convert meshes - start port0 at block=0
|
|
int blockNo = 0;
|
|
|
|
convertMeshVolume(output, blockNo);
|
|
convertMeshPatches(output, blockNo);
|
|
reader_->UpdateProgress(0.6);
|
|
|
|
if (reader_->GetIncludeZones())
|
|
{
|
|
convertMeshCellZones(output, blockNo);
|
|
convertMeshFaceZones(output, blockNo);
|
|
convertMeshPointZones(output, blockNo);
|
|
reader_->UpdateProgress(0.65);
|
|
}
|
|
|
|
if (reader_->GetIncludeSets())
|
|
{
|
|
convertMeshCellSets(output, blockNo);
|
|
convertMeshFaceSets(output, blockNo);
|
|
convertMeshPointSets(output, blockNo);
|
|
reader_->UpdateProgress(0.7);
|
|
}
|
|
|
|
if (outputLagrangian)
|
|
{
|
|
// Lagrangian dual port - restart port1 at block=0
|
|
blockNo = 0;
|
|
convertMeshLagrangian(outputLagrangian, blockNo);
|
|
}
|
|
else
|
|
{
|
|
convertMeshLagrangian(output, blockNo);
|
|
}
|
|
|
|
reader_->UpdateProgress(0.8);
|
|
|
|
// Update fields
|
|
convertVolFields(output);
|
|
convertPointFields(output);
|
|
convertLagrangianFields(outputLagrangian ? outputLagrangian : output);
|
|
|
|
if (debug)
|
|
{
|
|
Info<< "done reader part" << nl << endl;
|
|
}
|
|
reader_->UpdateProgress(0.95);
|
|
|
|
fieldsChanged_ = false;
|
|
meshState_ = polyMesh::UNCHANGED;
|
|
|
|
// Standard memory cleanup
|
|
cachedVtp_.clear();
|
|
cachedVtu_.clear();
|
|
}
|
|
|
|
|
|
void Foam::vtkPVFoam::UpdateFinalize()
|
|
{
|
|
if (!reader_->GetCacheMesh())
|
|
{
|
|
delete meshPtr_;
|
|
meshPtr_ = nullptr;
|
|
}
|
|
|
|
reader_->UpdateProgress(1.0);
|
|
}
|
|
|
|
|
|
std::vector<double> Foam::vtkPVFoam::findTimes(const bool skipZero) const
|
|
{
|
|
std::vector<double> times;
|
|
|
|
if (dbPtr_.valid())
|
|
{
|
|
const Time& runTime = dbPtr_();
|
|
instantList timeLst = runTime.times();
|
|
|
|
// find the first time for which this mesh appears to exist
|
|
label begIndex = timeLst.size();
|
|
forAll(timeLst, timeI)
|
|
{
|
|
if
|
|
(
|
|
IOobject
|
|
(
|
|
"points",
|
|
timeLst[timeI].name(),
|
|
meshDir_,
|
|
runTime
|
|
).typeHeaderOk<pointIOField>(false, false)
|
|
)
|
|
{
|
|
begIndex = timeI;
|
|
break;
|
|
}
|
|
}
|
|
|
|
label nTimes = timeLst.size() - begIndex;
|
|
|
|
// skip "constant" time whenever possible
|
|
if (begIndex == 0 && nTimes > 1)
|
|
{
|
|
if (timeLst[begIndex].name() == runTime.constant())
|
|
{
|
|
++begIndex;
|
|
--nTimes;
|
|
}
|
|
}
|
|
|
|
// skip "0/" time if requested and possible
|
|
if (skipZero && nTimes > 1 && timeLst[begIndex].name() == "0")
|
|
{
|
|
++begIndex;
|
|
--nTimes;
|
|
}
|
|
|
|
times.reserve(nTimes);
|
|
while (nTimes-- > 0)
|
|
{
|
|
times.push_back(timeLst[begIndex++].value());
|
|
}
|
|
}
|
|
else
|
|
{
|
|
if (debug)
|
|
{
|
|
cout<< "no valid dbPtr:\n";
|
|
}
|
|
}
|
|
|
|
return times;
|
|
}
|
|
|
|
|
|
void Foam::vtkPVFoam::renderPatchNames
|
|
(
|
|
vtkRenderer* renderer,
|
|
const bool show
|
|
)
|
|
{
|
|
// always remove old actors first
|
|
|
|
forAll(patchTextActors_, patchi)
|
|
{
|
|
renderer->RemoveViewProp(patchTextActors_[patchi]);
|
|
}
|
|
patchTextActors_.clear();
|
|
|
|
if (show && meshPtr_)
|
|
{
|
|
// get the display patches, strip off any prefix/suffix
|
|
hashedWordList selectedPatches = getSelected
|
|
(
|
|
reader_->GetPartSelection(),
|
|
rangePatches_
|
|
);
|
|
|
|
if (selectedPatches.empty())
|
|
{
|
|
return;
|
|
}
|
|
|
|
const polyBoundaryMesh& pbMesh = meshPtr_->boundaryMesh();
|
|
|
|
// Find the total number of zones
|
|
// Each zone will take the patch name
|
|
// Number of zones per patch ... zero zones should be skipped
|
|
labelList nZones(pbMesh.size(), 0);
|
|
|
|
// Per global zone number the average face centre position
|
|
List<DynamicList<point>> zoneCentre(pbMesh.size());
|
|
|
|
|
|
// Loop through all patches to determine zones, and centre of each zone
|
|
forAll(pbMesh, patchi)
|
|
{
|
|
const polyPatch& pp = pbMesh[patchi];
|
|
|
|
// Only include the patch if it is selected
|
|
if (!selectedPatches.found(pp.name()))
|
|
{
|
|
continue;
|
|
}
|
|
|
|
const labelListList& edgeFaces = pp.edgeFaces();
|
|
const vectorField& n = pp.faceNormals();
|
|
|
|
boolList featEdge(pp.nEdges(), false);
|
|
|
|
forAll(edgeFaces, edgeI)
|
|
{
|
|
const labelList& eFaces = edgeFaces[edgeI];
|
|
|
|
if (eFaces.size() == 1)
|
|
{
|
|
// Note: could also do ones with > 2 faces but this gives
|
|
// too many zones for baffles
|
|
featEdge[edgeI] = true;
|
|
}
|
|
else if (mag(n[eFaces[0]] & n[eFaces[1]]) < 0.5)
|
|
{
|
|
featEdge[edgeI] = true;
|
|
}
|
|
}
|
|
|
|
// Do topological analysis of patch, find disconnected regions
|
|
patchZones pZones(pp, featEdge);
|
|
|
|
nZones[patchi] = pZones.nZones();
|
|
|
|
labelList zoneNFaces(pZones.nZones(), 0);
|
|
|
|
// Create storage for additional zone centres
|
|
forAll(zoneNFaces, zoneI)
|
|
{
|
|
zoneCentre[patchi].append(Zero);
|
|
}
|
|
|
|
// Do averaging per individual zone
|
|
forAll(pp, facei)
|
|
{
|
|
label zoneI = pZones[facei];
|
|
zoneCentre[patchi][zoneI] += pp[facei].centre(pp.points());
|
|
zoneNFaces[zoneI]++;
|
|
}
|
|
|
|
forAll(zoneCentre[patchi], zoneI)
|
|
{
|
|
zoneCentre[patchi][zoneI] /= zoneNFaces[zoneI];
|
|
}
|
|
}
|
|
|
|
// Count number of zones we're actually going to display.
|
|
// This is truncated to a max per patch
|
|
|
|
const label MAXPATCHZONES = 20;
|
|
|
|
label displayZoneI = 0;
|
|
|
|
forAll(pbMesh, patchi)
|
|
{
|
|
displayZoneI += min(MAXPATCHZONES, nZones[patchi]);
|
|
}
|
|
|
|
if (debug)
|
|
{
|
|
Info<< "displayed zone centres = " << displayZoneI << nl
|
|
<< "zones per patch = " << nZones << endl;
|
|
}
|
|
|
|
// Set the size of the patch labels to max number of zones
|
|
patchTextActors_.setSize(displayZoneI);
|
|
|
|
if (debug)
|
|
{
|
|
Info<< "constructing patch labels" << endl;
|
|
}
|
|
|
|
// Actor index
|
|
displayZoneI = 0;
|
|
|
|
forAll(pbMesh, patchi)
|
|
{
|
|
const polyPatch& pp = pbMesh[patchi];
|
|
|
|
// Only selected patches will have a non-zero number of zones
|
|
const label nDisplayZones = min(MAXPATCHZONES, nZones[patchi]);
|
|
label increment = 1;
|
|
if (nZones[patchi] >= MAXPATCHZONES)
|
|
{
|
|
increment = nZones[patchi]/MAXPATCHZONES;
|
|
}
|
|
|
|
label globalZoneI = 0;
|
|
for (label i = 0; i < nDisplayZones; ++i, globalZoneI += increment)
|
|
{
|
|
if (debug)
|
|
{
|
|
Info<< "patch name = " << pp.name() << nl
|
|
<< "anchor = " << zoneCentre[patchi][globalZoneI] << nl
|
|
<< "globalZoneI = " << globalZoneI << endl;
|
|
}
|
|
|
|
// Into a list for later removal
|
|
patchTextActors_[displayZoneI++] = createTextActor
|
|
(
|
|
pp.name(),
|
|
zoneCentre[patchi][globalZoneI]
|
|
);
|
|
}
|
|
}
|
|
|
|
// Resize the patch names list to the actual number of patch names added
|
|
patchTextActors_.setSize(displayZoneI);
|
|
}
|
|
|
|
// Add text to each renderer
|
|
forAll(patchTextActors_, actori)
|
|
{
|
|
renderer->AddViewProp(patchTextActors_[actori]);
|
|
}
|
|
}
|
|
|
|
|
|
void Foam::vtkPVFoam::PrintSelf(ostream& os, vtkIndent indent) const
|
|
{
|
|
os << indent << "Number of nodes: "
|
|
<< (meshPtr_ ? meshPtr_->nPoints() : 0) << "\n";
|
|
|
|
os << indent << "Number of cells: "
|
|
<< (meshPtr_ ? meshPtr_->nCells() : 0) << "\n";
|
|
|
|
os << indent << "Number of available time steps: "
|
|
<< (dbPtr_.valid() ? dbPtr_().times().size() : 0) << "\n";
|
|
|
|
os << indent << "mesh region: " << meshRegion_ << "\n";
|
|
}
|
|
|
|
|
|
void Foam::vtkPVFoam::printInfo() const
|
|
{
|
|
std::cout
|
|
<< "Region: " << meshRegion_ << "\n"
|
|
<< "nPoints: " << (meshPtr_ ? meshPtr_->nPoints() : 0) << "\n"
|
|
<< "nCells: " << (meshPtr_ ? meshPtr_->nCells() : 0) << "\n"
|
|
<< "nTimes: "
|
|
<< (dbPtr_.valid() ? dbPtr_().times().size() : 0) << "\n";
|
|
|
|
std::vector<double> times = this->findTimes(reader_->GetSkipZeroTime());
|
|
|
|
std::cout<<" " << times.size() << "(";
|
|
for (const double& val : times)
|
|
{
|
|
std::cout<< ' ' << val;
|
|
}
|
|
std::cout << " )" << std::endl;
|
|
}
|
|
|
|
|
|
// ************************************************************************* //
|