1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright held by original author
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
11 OpenFOAM 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 2 of the License, or (at your
14 option) any later version.
16 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM; if not, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "vtkPV3Foam.H"
28 #include "vtkPV3FoamReader.h"
33 #include "patchZones.H"
36 #include "vtkDataArraySelection.h"
37 #include "vtkMultiBlockDataSet.h"
38 #include "vtkRenderer.h"
39 #include "vtkTextActor.h"
40 #include "vtkTextProperty.h"
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44 defineTypeNameAndDebug(Foam::vtkPV3Foam, 0);
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 #include "vtkPV3FoamAddToSelection.H"
49 #include "vtkPV3FoamUpdateInfoFields.H"
51 void Foam::vtkPV3Foam::resetCounters()
53 // Reset mesh part ids and sizes
54 partInfoVolume_.reset();
55 partInfoPatches_.reset();
56 partInfoLagrangian_.reset();
57 partInfoCellZones_.reset();
58 partInfoFaceZones_.reset();
59 partInfoPointZones_.reset();
60 partInfoCellSets_.reset();
61 partInfoFaceSets_.reset();
62 partInfoPointSets_.reset();
66 void Foam::vtkPV3Foam::reduceMemory()
68 forAll(regionPolyDecomp_, i)
70 regionPolyDecomp_[i].clear();
73 forAll(zonePolyDecomp_, i)
75 zonePolyDecomp_[i].clear();
78 forAll(csetPolyDecomp_, i)
80 csetPolyDecomp_[i].clear();
83 if (!reader_->GetCacheMesh())
93 int Foam::vtkPV3Foam::setTime(int nRequest, const double requestTimes[])
97 Info<< "<beg> Foam::vtkPV3Foam::setTime(";
98 for (int requestI = 0; requestI < nRequest; ++requestI)
105 Info<< requestTimes[requestI];
107 Info << ") - previousIndex = " << timeIndex_ << endl;
110 Time& runTime = dbPtr_();
113 instantList Times = runTime.times();
115 int nearestIndex = timeIndex_;
116 for (int requestI = 0; requestI < nRequest; ++requestI)
118 int index = Time::findClosestTimeIndex(Times, requestTimes[requestI]);
119 if (index >= 0 && index != timeIndex_)
121 nearestIndex = index;
126 if (nearestIndex < 0)
132 // see what has changed
133 if (timeIndex_ != nearestIndex)
135 timeIndex_ = nearestIndex;
136 runTime.setTime(Times[nearestIndex], nearestIndex);
138 // the fields change each time
139 fieldsChanged_ = true;
143 if (meshPtr_->readUpdate() != polyMesh::UNCHANGED)
153 reader_->UpdateProgress(0.05);
155 // this seems to be needed for catching Lagrangian fields
161 Info<< "<end> Foam::vtkPV3Foam::setTime() - selectedTime="
162 << Times[nearestIndex].name() << " index=" << timeIndex_
163 << "/" << Times.size()
164 << " meshChanged=" << Switch(meshChanged_)
165 << " fieldsChanged=" << Switch(fieldsChanged_) << endl;
172 void Foam::vtkPV3Foam::updateMeshPartsStatus()
176 Info<< "<beg> Foam::vtkPV3Foam::updateMeshPartsStatus" << endl;
179 vtkDataArraySelection* selection = reader_->GetPartSelection();
180 label nElem = selection->GetNumberOfArrays();
182 if (partStatus_.size() != nElem)
184 partStatus_.setSize(nElem);
189 // this needs fixing if we wish to re-use the datasets
190 partDataset_.setSize(nElem);
193 // Read the selected mesh parts (zones, patches ...) and add to list
194 forAll(partStatus_, partId)
196 const int setting = selection->GetArraySetting(partId);
198 if (partStatus_[partId] != setting)
200 partStatus_[partId] = setting;
206 Info<< " part[" << partId << "] = "
207 << partStatus_[partId]
208 << " : " << selection->GetArrayName(partId) << endl;
213 Info<< "<end> Foam::vtkPV3Foam::updateMeshPartsStatus" << endl;
217 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
219 Foam::vtkPV3Foam::vtkPV3Foam
221 const char* const FileName,
222 vtkPV3FoamReader* reader
228 meshRegion_(polyMesh::defaultRegion),
229 meshDir_(polyMesh::meshSubDir),
232 fieldsChanged_(true),
233 partInfoVolume_("unzoned"),
234 partInfoPatches_("patches"),
235 partInfoLagrangian_("lagrangian"),
236 partInfoCellZones_("cellZone"),
237 partInfoFaceZones_("faceZone"),
238 partInfoPointZones_("pointZone"),
239 partInfoCellSets_("cellSet"),
240 partInfoFaceSets_("faceSet"),
241 partInfoPointSets_("pointSet")
245 Info<< "Foam::vtkPV3Foam::vtkPV3Foam - " << FileName << endl;
249 // avoid argList and get rootPath/caseName directly from the file
250 fileName fullCasePath(fileName(FileName).path());
252 if (!isDir(fullCasePath))
256 if (fullCasePath == ".")
258 fullCasePath = cwd();
261 // Set the case as an environment variable - some BCs might use this
262 if (fullCasePath.name().find("processor", 0) == 0)
264 const fileName globalCase = fullCasePath.path();
266 setEnv("FOAM_CASE", globalCase, true);
267 setEnv("FOAM_CASENAME", globalCase.name(), true);
271 setEnv("FOAM_CASE", fullCasePath, true);
272 setEnv("FOAM_CASENAME", fullCasePath.name(), true);
275 // look for 'case{region}.OpenFOAM'
276 // could be stringent and insist the prefix match the directory name...
277 // Note: cannot use fileName::name() due to the embedded '{}'
278 string caseName(fileName(FileName).lessExt());
279 string::size_type beg = caseName.find_last_of("/{");
280 string::size_type end = caseName.find('}', beg);
284 beg != string::npos && caseName[beg] == '{'
285 && end != string::npos && end == caseName.size()-1
288 meshRegion_ = caseName.substr(beg+1, end-beg-1);
291 if (!meshRegion_.size())
293 meshRegion_ = polyMesh::defaultRegion;
296 if (meshRegion_ != polyMesh::defaultRegion)
298 meshDir_ = meshRegion_/polyMesh::meshSubDir;
304 Info<< "fullCasePath=" << fullCasePath << nl
305 << "FOAM_CASE=" << getEnv("FOAM_CASE") << nl
306 << "FOAM_CASENAME=" << getEnv("FOAM_CASENAME") << nl
307 << "region=" << meshRegion_ << endl;
310 // Create time object
315 Time::controlDictName,
316 fileName(fullCasePath.path()),
317 fileName(fullCasePath.name())
321 dbPtr_().functionObjects().off();
327 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
329 Foam::vtkPV3Foam::~vtkPV3Foam()
333 Info<< "<end> Foam::vtkPV3Foam::~vtkPV3Foam" << endl;
340 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
342 void Foam::vtkPV3Foam::updateInfo()
346 Info<< "<beg> Foam::vtkPV3Foam::updateInfo"
347 << " [meshPtr=" << (meshPtr_ ? "set" : "NULL") << "] timeIndex="
348 << timeIndex_ << endl;
353 vtkDataArraySelection* partSelection = reader_->GetPartSelection();
355 // enable 'internalMesh' on the first call
356 // or preserve the enabled selections
357 stringList enabledEntries;
358 if (!partSelection->GetNumberOfArrays() && !meshPtr_)
360 enabledEntries.setSize(1);
361 enabledEntries[0] = "internalMesh";
365 enabledEntries = getSelectedArrayEntries(partSelection);
368 // Clear current mesh parts list
369 partSelection->RemoveAllArrays();
371 // Update mesh parts list - add Lagrangian at the bottom
372 updateInfoInternalMesh();
376 updateInfoLagrangian();
378 // restore the enabled selections
379 setSelectedArrayEntries(partSelection, enabledEntries);
383 fieldsChanged_ = true;
386 // Update volume, point and lagrangian fields
387 updateInfoFields<fvPatchField, volMesh>
389 reader_->GetVolFieldSelection()
391 updateInfoFields<pointPatchField, pointMesh>
393 reader_->GetPointFieldSelection()
395 updateInfoLagrangianFields();
399 Info<< "<end> Foam::vtkPV3Foam::updateInfo" << endl;
405 void Foam::vtkPV3Foam::updateFoamMesh()
409 Info<< "<beg> Foam::vtkPV3Foam::updateFoamMesh" << endl;
413 if (!reader_->GetCacheMesh())
419 // Check to see if the FOAM mesh has been created
424 Info<< "Creating Foam mesh for region " << meshRegion_
425 << " at time=" << dbPtr_().timeName()
430 meshPtr_ = new fvMesh
447 Info<< "Using existing Foam mesh" << endl;
453 Info<< "<end> Foam::vtkPV3Foam::updateFoamMesh" << endl;
459 void Foam::vtkPV3Foam::Update
461 vtkMultiBlockDataSet* output,
462 vtkMultiBlockDataSet* lagrangianOutput
467 cout<< "<beg> Foam::vtkPV3Foam::Update - output with "
468 << output->GetNumberOfBlocks() << " and "
469 << lagrangianOutput->GetNumberOfBlocks() << " blocks\n";
471 lagrangianOutput->Print(cout);
474 reader_->UpdateProgress(0.1);
476 // Set up mesh parts selection(s)
477 updateMeshPartsStatus();
479 reader_->UpdateProgress(0.15);
481 // Update the Foam mesh
483 reader_->UpdateProgress(0.4);
485 // Convert meshes - start port0 at block=0
488 convertMeshVolume(output, blockNo);
489 convertMeshPatches(output, blockNo);
490 reader_->UpdateProgress(0.6);
492 if (reader_->GetIncludeZones())
494 convertMeshCellZones(output, blockNo);
495 convertMeshFaceZones(output, blockNo);
496 convertMeshPointZones(output, blockNo);
497 reader_->UpdateProgress(0.65);
500 if (reader_->GetIncludeSets())
502 convertMeshCellSets(output, blockNo);
503 convertMeshFaceSets(output, blockNo);
504 convertMeshPointSets(output, blockNo);
505 reader_->UpdateProgress(0.7);
508 #ifdef VTKPV3FOAM_DUALPORT
509 // restart port1 at block=0
512 convertMeshLagrangian(lagrangianOutput, blockNo);
514 reader_->UpdateProgress(0.8);
517 convertVolFields(output);
518 convertPointFields(output);
519 convertLagrangianFields(lagrangianOutput);
520 reader_->UpdateProgress(0.95);
522 meshChanged_ = fieldsChanged_ = false;
526 void Foam::vtkPV3Foam::CleanUp()
528 // reclaim some memory
530 reader_->UpdateProgress(1.0);
534 double* Foam::vtkPV3Foam::findTimes(int& nTimeSteps)
537 double* tsteps = NULL;
541 Time& runTime = dbPtr_();
542 instantList timeLst = runTime.times();
544 // find the first time for which this mesh appears to exist
546 for (; timeI < timeLst.size(); ++timeI)
548 const word& timeName = timeLst[timeI].name();
552 isFile(runTime.path()/timeName/meshDir_/"points")
553 && IOobject("points", timeName, meshDir_, runTime).headerOk()
560 nTimes = timeLst.size() - timeI;
562 // always skip "constant" time if possible
563 if (timeI == 0 && nTimes > 1)
571 tsteps = new double[nTimes];
572 for (label stepI = 0; stepI < nTimes; ++stepI, ++timeI)
574 tsteps[stepI] = timeLst[timeI].value();
582 cout<< "no valid dbPtr:\n";
586 // vector length returned via the parameter
593 void Foam::vtkPV3Foam::addPatchNames(vtkRenderer* renderer)
595 // Remove any patch names previously added to the renderer
596 removePatchNames(renderer);
598 // get the display patches, strip off any suffix
599 wordHashSet selectedPatches = getSelected
601 reader_->GetPartSelection(),
605 if (!selectedPatches.size())
612 Info<< "<beg> Foam::vtkPV3Foam::addPatchNames" << nl
613 <<"... add patches: " << selectedPatches << endl;
616 const polyBoundaryMesh& pbMesh = meshPtr_->boundaryMesh();
618 // Find the total number of zones
619 // Each zone will take the patch name
620 // Number of zones per patch ... zero zones should be skipped
621 labelList nZones(pbMesh.size(), 0);
623 // Per global zone number the average face centre position
624 DynamicList<point> zoneCentre(pbMesh.size());
628 Info<< "... determining patch zones" << endl;
631 // Loop through all patches to determine zones, and centre of each zone
632 forAll(pbMesh, patchI)
634 const polyPatch& pp = pbMesh[patchI];
636 // Only include the patch if it is selected
637 if (!selectedPatches.found(pp.name()))
642 const labelListList& edgeFaces = pp.edgeFaces();
643 const vectorField& n = pp.faceNormals();
645 boolList featEdge(pp.nEdges(), false);
647 forAll(edgeFaces, edgeI)
649 const labelList& eFaces = edgeFaces[edgeI];
651 if (eFaces.size() == 1)
653 // Note: could also do ones with > 2 faces but this gives
654 // too many zones for baffles
655 featEdge[edgeI] = true;
657 else if (mag(n[eFaces[0]] & n[eFaces[1]]) < 0.5)
659 featEdge[edgeI] = true;
663 // Do topological analysis of patch, find disconnected regions
664 patchZones pZones(pp, featEdge);
666 nZones[patchI] = pZones.nZones();
668 labelList zoneNFaces(pZones.nZones(), 0);
670 // Save start of information for current patch
671 label patchStart = zoneCentre.size();
673 // Create storage for additional zone centres
674 forAll(zoneNFaces, zoneI)
676 zoneCentre.append(vector::zero);
679 // Do averaging per individual zone
682 label zoneI = pZones[faceI];
683 zoneCentre[patchStart+zoneI] += pp[faceI].centre(pp.points());
687 for (label i=0; i<nZones[patchI]; i++)
689 zoneCentre[patchStart + i] /= zoneNFaces[i];
693 // Count number of zones we're actually going to display. This is truncated
694 // to a max per patch
696 const label MAXPATCHZONES = 20;
698 label displayZoneI = 0;
700 forAll(pbMesh, patchI)
702 displayZoneI += min(MAXPATCHZONES, nZones[patchI]);
710 Info<< "patch zone centres = " << zoneCentre << nl
711 << "displayed zone centres = " << displayZoneI << nl
712 << "zones per patch = " << nZones << endl;
715 // Set the size of the patch labels to max number of zones
716 patchTextActorsPtrs_.setSize(displayZoneI);
720 Info<< "constructing patch labels" << endl;
726 // Index in zone centres
727 label globalZoneI = 0;
729 forAll(pbMesh, patchI)
731 const polyPatch& pp = pbMesh[patchI];
733 // Only selected patches will have a non-zero number of zones
734 label nDisplayZones = min(MAXPATCHZONES, nZones[patchI]);
736 if (nZones[patchI] >= MAXPATCHZONES)
738 increment = nZones[patchI]/MAXPATCHZONES;
741 for (label i = 0; i < nDisplayZones; i++)
745 Info<< "patch name = " << pp.name() << nl
746 << "anchor = " << zoneCentre[globalZoneI] << nl
747 << "globalZoneI = " << globalZoneI << endl;
750 vtkTextActor* txt = vtkTextActor::New();
752 txt->SetInput(pp.name().c_str());
754 // Set text properties
755 vtkTextProperty* tprop = txt->GetTextProperty();
756 tprop->SetFontFamilyToArial();
759 tprop->SetLineSpacing(1.0);
760 tprop->SetFontSize(12);
761 tprop->SetColor(1.0, 0.0, 0.0);
762 tprop->SetJustificationToCentered();
764 // Set text to use 3-D world co-ordinates
765 txt->GetPositionCoordinate()->SetCoordinateSystemToWorld();
767 txt->GetPositionCoordinate()->SetValue
769 zoneCentre[globalZoneI].x(),
770 zoneCentre[globalZoneI].y(),
771 zoneCentre[globalZoneI].z()
774 // Add text to each renderer
775 renderer->AddViewProp(txt);
777 // Maintain a list of text labels added so that they can be
779 patchTextActorsPtrs_[displayZoneI] = txt;
781 globalZoneI += increment;
786 // Resize the patch names list to the actual number of patch names added
787 patchTextActorsPtrs_.setSize(displayZoneI);
791 Info<< "<end> Foam::vtkPV3Foam::addPatchNames" << endl;
796 void Foam::vtkPV3Foam::removePatchNames(vtkRenderer* renderer)
798 forAll(patchTextActorsPtrs_, patchI)
800 renderer->RemoveViewProp(patchTextActorsPtrs_[patchI]);
801 patchTextActorsPtrs_[patchI]->Delete();
803 patchTextActorsPtrs_.clear();
807 void Foam::vtkPV3Foam::PrintSelf(ostream& os, vtkIndent indent) const
809 os << indent << "Number of nodes: "
810 << (meshPtr_ ? meshPtr_->nPoints() : 0) << "\n";
812 os << indent << "Number of cells: "
813 << (meshPtr_ ? meshPtr_->nCells() : 0) << "\n";
815 os << indent << "Number of available time steps: "
816 << (dbPtr_.valid() ? dbPtr_().times().size() : 0) << endl;
818 os << indent << "mesh region: " << meshRegion_ << "\n";
821 // ************************************************************************* //