1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | foam-extend: Open Source CFD
4 \\ / O peration | Version: 3.2
5 \\ / A nd | Web: http://www.foam-extend.org
6 \\/ M anipulation | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
9 This file is part of foam-extend.
11 foam-extend is free software: you can redistribute it and/or modify it
12 under the terms of the GNU General Public License as published by the
13 Free Software Foundation, either version 3 of the License, or (at your
14 option) any later version.
16 foam-extend is distributed in the hope that it will be useful, but
17 WITHOUT ANY WARRANTY; without even the implied warranty of
18 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 General Public License for more details.
21 You should have received a copy of the GNU General Public License
22 along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "vtkPV3Foam.H"
27 #include "vtkPV3FoamReader.h"
32 #include "patchZones.H"
35 #include "vtkDataArraySelection.h"
36 #include "vtkMultiBlockDataSet.h"
37 #include "vtkRenderer.h"
38 #include "vtkTextActor.h"
39 #include "vtkTextProperty.h"
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 defineTypeNameAndDebug(Foam::vtkPV3Foam, 0);
45 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 #include "vtkPV3FoamAddToSelection.H"
48 #include "vtkPV3FoamUpdateInfoFields.H"
50 void Foam::vtkPV3Foam::resetCounters()
52 // Reset mesh part ids and sizes
53 partInfoVolume_.reset();
54 partInfoPatches_.reset();
55 partInfoLagrangian_.reset();
56 partInfoCellZones_.reset();
57 partInfoFaceZones_.reset();
58 partInfoPointZones_.reset();
59 partInfoCellSets_.reset();
60 partInfoFaceSets_.reset();
61 partInfoPointSets_.reset();
65 void Foam::vtkPV3Foam::reduceMemory()
67 forAll(regionPolyDecomp_, i)
69 regionPolyDecomp_[i].clear();
72 forAll(zonePolyDecomp_, i)
74 zonePolyDecomp_[i].clear();
77 forAll(csetPolyDecomp_, i)
79 csetPolyDecomp_[i].clear();
82 if (!reader_->GetCacheMesh())
92 int Foam::vtkPV3Foam::setTime(int nRequest, const double requestTimes[])
96 Info<< "<beg> Foam::vtkPV3Foam::setTime(";
97 for (int requestI = 0; requestI < nRequest; ++requestI)
104 Info<< requestTimes[requestI];
106 Info << ") - previousIndex = " << timeIndex_ << endl;
109 Time& runTime = dbPtr_();
112 instantList Times = runTime.times();
114 int nearestIndex = timeIndex_;
115 for (int requestI = 0; requestI < nRequest; ++requestI)
117 int index = Time::findClosestTimeIndex(Times, requestTimes[requestI]);
118 if (index >= 0 && index != timeIndex_)
120 nearestIndex = index;
125 if (nearestIndex < 0)
131 // see what has changed
132 if (timeIndex_ != nearestIndex)
134 timeIndex_ = nearestIndex;
135 runTime.setTime(Times[nearestIndex], nearestIndex);
137 // the fields change each time
138 fieldsChanged_ = true;
142 if (meshPtr_->readUpdate() != polyMesh::UNCHANGED)
152 reader_->UpdateProgress(0.05);
154 // this seems to be needed for catching Lagrangian fields
160 Info<< "<end> Foam::vtkPV3Foam::setTime() - selectedTime="
161 << Times[nearestIndex].name() << " index=" << timeIndex_
162 << "/" << Times.size()
163 << " meshChanged=" << Switch(meshChanged_)
164 << " fieldsChanged=" << Switch(fieldsChanged_) << endl;
171 void Foam::vtkPV3Foam::updateMeshPartsStatus()
175 Info<< "<beg> Foam::vtkPV3Foam::updateMeshPartsStatus" << endl;
178 vtkDataArraySelection* selection = reader_->GetPartSelection();
179 label nElem = selection->GetNumberOfArrays();
181 if (partStatus_.size() != nElem)
183 partStatus_.setSize(nElem);
188 // this needs fixing if we wish to re-use the datasets
189 partDataset_.setSize(nElem);
192 // Read the selected mesh parts (zones, patches ...) and add to list
193 forAll(partStatus_, partId)
195 const int setting = selection->GetArraySetting(partId);
197 if (partStatus_[partId] != setting)
199 partStatus_[partId] = setting;
205 Info<< " part[" << partId << "] = "
206 << partStatus_[partId]
207 << " : " << selection->GetArrayName(partId) << endl;
212 Info<< "<end> Foam::vtkPV3Foam::updateMeshPartsStatus" << endl;
216 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
218 Foam::vtkPV3Foam::vtkPV3Foam
220 const char* const FileName,
221 vtkPV3FoamReader* reader
227 meshRegion_(polyMesh::defaultRegion),
228 meshDir_(polyMesh::meshSubDir),
231 fieldsChanged_(true),
232 partInfoVolume_("unzoned"),
233 partInfoPatches_("patches"),
234 partInfoLagrangian_("lagrangian"),
235 partInfoCellZones_("cellZone"),
236 partInfoFaceZones_("faceZone"),
237 partInfoPointZones_("pointZone"),
238 partInfoCellSets_("cellSet"),
239 partInfoFaceSets_("faceSet"),
240 partInfoPointSets_("pointSet")
244 Info<< "Foam::vtkPV3Foam::vtkPV3Foam - " << FileName << endl;
248 // avoid argList and get rootPath/caseName directly from the file
249 fileName fullCasePath(fileName(FileName).path());
251 if (!isDir(fullCasePath))
255 if (fullCasePath == ".")
257 fullCasePath = cwd();
260 // Set the case as an environment variable - some BCs might use this
261 if (fullCasePath.name().find("processor", 0) == 0)
263 const fileName globalCase = fullCasePath.path();
265 setEnv("FOAM_CASE", globalCase, true);
266 setEnv("FOAM_CASENAME", globalCase.name(), true);
270 setEnv("FOAM_CASE", fullCasePath, true);
271 setEnv("FOAM_CASENAME", fullCasePath.name(), true);
274 // look for 'case{region}.FOAM'
275 // could be stringent and insist the prefix match the directory name...
276 // Note: cannot use fileName::name() due to the embedded '{}'
277 string caseName(fileName(FileName).lessExt());
278 string::size_type beg = caseName.find_last_of("/{");
279 string::size_type end = caseName.find('}', beg);
283 beg != string::npos && caseName[beg] == '{'
284 && end != string::npos && end == caseName.size()-1
287 meshRegion_ = caseName.substr(beg+1, end-beg-1);
290 if (!meshRegion_.size())
292 meshRegion_ = polyMesh::defaultRegion;
295 if (meshRegion_ != polyMesh::defaultRegion)
297 meshDir_ = meshRegion_/polyMesh::meshSubDir;
303 Info<< "fullCasePath=" << fullCasePath << nl
304 << "FOAM_CASE=" << getEnv("FOAM_CASE") << nl
305 << "FOAM_CASENAME=" << getEnv("FOAM_CASENAME") << nl
306 << "region=" << meshRegion_ << endl;
309 // Create time object
314 Time::controlDictName,
315 fileName(fullCasePath.path()),
316 fileName(fullCasePath.name())
320 dbPtr_().functionObjects().off();
326 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
328 Foam::vtkPV3Foam::~vtkPV3Foam()
332 Info<< "<end> Foam::vtkPV3Foam::~vtkPV3Foam" << endl;
339 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
341 void Foam::vtkPV3Foam::updateInfo()
345 Info<< "<beg> Foam::vtkPV3Foam::updateInfo"
346 << " [meshPtr=" << (meshPtr_ ? "set" : "NULL") << "] timeIndex="
347 << timeIndex_ << endl;
352 vtkDataArraySelection* partSelection = reader_->GetPartSelection();
354 // enable 'internalMesh' on the first call
355 // or preserve the enabled selections
356 stringList enabledEntries;
357 if (!partSelection->GetNumberOfArrays() && !meshPtr_)
359 enabledEntries.setSize(1);
360 enabledEntries[0] = "internalMesh";
364 enabledEntries = getSelectedArrayEntries(partSelection);
367 // Clear current mesh parts list
368 partSelection->RemoveAllArrays();
370 // Update mesh parts list - add Lagrangian at the bottom
371 updateInfoInternalMesh();
375 updateInfoLagrangian();
377 // restore the enabled selections
378 setSelectedArrayEntries(partSelection, enabledEntries);
382 fieldsChanged_ = true;
385 // Update volume, point and lagrangian fields
386 updateInfoFields<fvPatchField, volMesh>
388 reader_->GetVolFieldSelection()
390 updateInfoFields<pointPatchField, pointMesh>
392 reader_->GetPointFieldSelection()
394 updateInfoLagrangianFields();
398 Info<< "<end> Foam::vtkPV3Foam::updateInfo" << endl;
404 void Foam::vtkPV3Foam::updateFoamMesh()
408 Info<< "<beg> Foam::vtkPV3Foam::updateFoamMesh" << endl;
412 if (!reader_->GetCacheMesh())
418 // Check to see if the FOAM mesh has been created
423 Info<< "Creating Foam mesh for region " << meshRegion_
424 << " at time=" << dbPtr_().timeName()
429 meshPtr_ = new fvMesh
446 Info<< "Using existing Foam mesh" << endl;
452 Info<< "<end> Foam::vtkPV3Foam::updateFoamMesh" << endl;
458 void Foam::vtkPV3Foam::Update
460 vtkMultiBlockDataSet* output,
461 vtkMultiBlockDataSet* lagrangianOutput
466 cout<< "<beg> Foam::vtkPV3Foam::Update - output with "
467 << output->GetNumberOfBlocks() << " and "
468 << lagrangianOutput->GetNumberOfBlocks() << " blocks\n";
470 lagrangianOutput->Print(cout);
473 reader_->UpdateProgress(0.1);
475 // Set up mesh parts selection(s)
476 updateMeshPartsStatus();
478 reader_->UpdateProgress(0.15);
480 // Update the Foam mesh
482 reader_->UpdateProgress(0.4);
484 // Convert meshes - start port0 at block=0
487 convertMeshVolume(output, blockNo);
488 convertMeshPatches(output, blockNo);
489 reader_->UpdateProgress(0.6);
491 if (reader_->GetIncludeZones())
493 convertMeshCellZones(output, blockNo);
494 convertMeshFaceZones(output, blockNo);
495 convertMeshPointZones(output, blockNo);
496 reader_->UpdateProgress(0.65);
499 if (reader_->GetIncludeSets())
501 convertMeshCellSets(output, blockNo);
502 convertMeshFaceSets(output, blockNo);
503 convertMeshPointSets(output, blockNo);
504 reader_->UpdateProgress(0.7);
507 #ifdef VTKPV3FOAM_DUALPORT
508 // restart port1 at block=0
511 convertMeshLagrangian(lagrangianOutput, blockNo);
513 reader_->UpdateProgress(0.8);
516 convertVolFields(output);
517 convertPointFields(output);
518 convertLagrangianFields(lagrangianOutput);
519 reader_->UpdateProgress(0.95);
521 meshChanged_ = fieldsChanged_ = false;
525 void Foam::vtkPV3Foam::CleanUp()
527 // reclaim some memory
529 reader_->UpdateProgress(1.0);
533 double* Foam::vtkPV3Foam::findTimes(int& nTimeSteps)
536 double* tsteps = NULL;
540 Time& runTime = dbPtr_();
541 instantList timeLst = runTime.times();
543 // find the first time for which this mesh appears to exist
545 for (; timeI < timeLst.size(); ++timeI)
547 const word& timeName = timeLst[timeI].name();
551 isFile(runTime.path()/timeName/meshDir_/"points")
552 && IOobject("points", timeName, meshDir_, runTime).headerOk()
559 nTimes = timeLst.size() - timeI;
561 // always skip "constant" time if possible
562 if (timeI == 0 && nTimes > 1)
570 tsteps = new double[nTimes];
571 for (label stepI = 0; stepI < nTimes; ++stepI, ++timeI)
573 tsteps[stepI] = timeLst[timeI].value();
581 cout<< "no valid dbPtr:\n";
585 // vector length returned via the parameter
592 void Foam::vtkPV3Foam::addPatchNames(vtkRenderer* renderer)
594 // Remove any patch names previously added to the renderer
595 removePatchNames(renderer);
597 // get the display patches, strip off any suffix
598 wordHashSet selectedPatches = getSelected
600 reader_->GetPartSelection(),
604 if (!selectedPatches.size())
611 Info<< "<beg> Foam::vtkPV3Foam::addPatchNames" << nl
612 <<"... add patches: " << selectedPatches << endl;
615 const polyBoundaryMesh& pbMesh = meshPtr_->boundaryMesh();
617 // Find the total number of zones
618 // Each zone will take the patch name
619 // Number of zones per patch ... zero zones should be skipped
620 labelList nZones(pbMesh.size(), 0);
622 // Per global zone number the average face centre position
623 DynamicList<point> zoneCentre(pbMesh.size());
627 Info<< "... determining patch zones" << endl;
630 // Loop through all patches to determine zones, and centre of each zone
631 forAll(pbMesh, patchI)
633 const polyPatch& pp = pbMesh[patchI];
635 // Only include the patch if it is selected
636 if (!selectedPatches.found(pp.name()))
641 const labelListList& edgeFaces = pp.edgeFaces();
642 const vectorField& n = pp.faceNormals();
644 boolList featEdge(pp.nEdges(), false);
646 forAll(edgeFaces, edgeI)
648 const labelList& eFaces = edgeFaces[edgeI];
650 if (eFaces.size() == 1)
652 // Note: could also do ones with > 2 faces but this gives
653 // too many zones for baffles
654 featEdge[edgeI] = true;
656 else if (mag(n[eFaces[0]] & n[eFaces[1]]) < 0.5)
658 featEdge[edgeI] = true;
662 // Do topological analysis of patch, find disconnected regions
663 patchZones pZones(pp, featEdge);
665 nZones[patchI] = pZones.nZones();
667 labelList zoneNFaces(pZones.nZones(), 0);
669 // Save start of information for current patch
670 label patchStart = zoneCentre.size();
672 // Create storage for additional zone centres
673 forAll(zoneNFaces, zoneI)
675 zoneCentre.append(vector::zero);
678 // Do averaging per individual zone
681 label zoneI = pZones[faceI];
682 zoneCentre[patchStart+zoneI] += pp[faceI].centre(pp.points());
686 for (label i=0; i<nZones[patchI]; i++)
688 zoneCentre[patchStart + i] /= zoneNFaces[i];
692 // Count number of zones we're actually going to display. This is truncated
693 // to a max per patch
695 const label MAXPATCHZONES = 20;
697 label displayZoneI = 0;
699 forAll(pbMesh, patchI)
701 displayZoneI += min(MAXPATCHZONES, nZones[patchI]);
709 Info<< "patch zone centres = " << zoneCentre << nl
710 << "displayed zone centres = " << displayZoneI << nl
711 << "zones per patch = " << nZones << endl;
714 // Set the size of the patch labels to max number of zones
715 patchTextActorsPtrs_.setSize(displayZoneI);
719 Info<< "constructing patch labels" << endl;
725 // Index in zone centres
726 label globalZoneI = 0;
728 forAll(pbMesh, patchI)
730 const polyPatch& pp = pbMesh[patchI];
732 // Only selected patches will have a non-zero number of zones
733 label nDisplayZones = min(MAXPATCHZONES, nZones[patchI]);
735 if (nZones[patchI] >= MAXPATCHZONES)
737 increment = nZones[patchI]/MAXPATCHZONES;
740 for (label i = 0; i < nDisplayZones; i++)
744 Info<< "patch name = " << pp.name() << nl
745 << "anchor = " << zoneCentre[globalZoneI] << nl
746 << "globalZoneI = " << globalZoneI << endl;
749 vtkTextActor* txt = vtkTextActor::New();
751 txt->SetInput(pp.name().c_str());
753 // Set text properties
754 vtkTextProperty* tprop = txt->GetTextProperty();
755 tprop->SetFontFamilyToArial();
758 tprop->SetLineSpacing(1.0);
759 tprop->SetFontSize(12);
760 tprop->SetColor(1.0, 0.0, 0.0);
761 tprop->SetJustificationToCentered();
763 // Set text to use 3-D world co-ordinates
764 txt->GetPositionCoordinate()->SetCoordinateSystemToWorld();
766 txt->GetPositionCoordinate()->SetValue
768 zoneCentre[globalZoneI].x(),
769 zoneCentre[globalZoneI].y(),
770 zoneCentre[globalZoneI].z()
773 // Add text to each renderer
774 renderer->AddViewProp(txt);
776 // Maintain a list of text labels added so that they can be
778 patchTextActorsPtrs_[displayZoneI] = txt;
780 globalZoneI += increment;
785 // Resize the patch names list to the actual number of patch names added
786 patchTextActorsPtrs_.setSize(displayZoneI);
790 Info<< "<end> Foam::vtkPV3Foam::addPatchNames" << endl;
795 void Foam::vtkPV3Foam::removePatchNames(vtkRenderer* renderer)
797 forAll(patchTextActorsPtrs_, patchI)
799 renderer->RemoveViewProp(patchTextActorsPtrs_[patchI]);
800 patchTextActorsPtrs_[patchI]->Delete();
802 patchTextActorsPtrs_.clear();
806 void Foam::vtkPV3Foam::PrintSelf(ostream& os, vtkIndent indent) const
808 os << indent << "Number of nodes: "
809 << (meshPtr_ ? meshPtr_->nPoints() : 0) << "\n";
811 os << indent << "Number of cells: "
812 << (meshPtr_ ? meshPtr_->nCells() : 0) << "\n";
814 os << indent << "Number of available time steps: "
815 << (dbPtr_.valid() ? dbPtr_().times().size() : 0) << endl;
817 os << indent << "mesh region: " << meshRegion_ << "\n";
820 // ************************************************************************* //