1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2008 OpenCFD Ltd.
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"
32 #include "IOobjectList.H"
33 #include "patchZones.H"
34 #include "vtkPV3FoamReader.h"
38 #include "vtkCharArray.h"
39 #include "vtkDataArraySelection.h"
40 #include "vtkDataSet.h"
41 #include "vtkFieldData.h"
42 #include "vtkMultiBlockDataSet.h"
43 #include "vtkRenderer.h"
44 #include "vtkTextActor.h"
45 #include "vtkTextProperty.h"
46 #include "vtkPolyData.h"
47 #include "vtkUnstructuredGrid.h"
48 #include "vtkInformation.h"
50 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
52 defineTypeNameAndDebug(Foam::vtkPV3Foam, 0);
54 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
56 #include "vtkPV3FoamAddToSelection.H"
57 #include "vtkPV3FoamUpdateInformationFields.H"
60 void Foam::vtkPV3Foam::AddToBlock
62 vtkMultiBlockDataSet* output,
63 const selectionInfo& selector,
64 const label datasetNo,
66 const string& blockName
69 const int blockNo = selector.block();
71 vtkDataObject* blockDO = output->GetBlock(blockNo);
72 vtkMultiBlockDataSet* block = vtkMultiBlockDataSet::SafeDownCast(blockDO);
73 if (blockDO && !block)
75 FatalErrorIn("Foam::vtkPV3Foam::AddToBlock")
76 << "Block already has a vtkDataSet assigned to it" << nl << endl;
82 block = vtkMultiBlockDataSet::New();
83 output->SetBlock(blockNo, block);
91 Info<< "block[" << blockNo << "] has "
92 << block->GetNumberOfBlocks()
93 << " datasets prior to adding set " << datasetNo
94 << " with name: " << blockName << endl;
97 // when assigning dataset 0, also name the parent block
98 if (!datasetNo && selector.name())
100 output->GetMetaData(blockNo)->Set
102 vtkCompositeDataSet::NAME(),
109 block->SetBlock(datasetNo, dataset);
111 if (blockName.size())
113 block->GetMetaData(datasetNo)->Set
115 vtkCompositeDataSet::NAME(), blockName.c_str()
122 vtkDataSet* Foam::vtkPV3Foam::GetDataSetFromBlock
124 vtkMultiBlockDataSet* output,
125 const selectionInfo& selector,
126 const label datasetNo
129 const int blockNo = selector.block();
131 vtkDataObject* blockDO = output->GetBlock(blockNo);
132 vtkMultiBlockDataSet* block = vtkMultiBlockDataSet::SafeDownCast(blockDO);
135 return vtkDataSet::SafeDownCast(block->GetBlock(datasetNo));
142 Foam::label Foam::vtkPV3Foam::GetNumberOfDataSets
144 vtkMultiBlockDataSet* output,
145 const selectionInfo& selector
148 const int blockNo = selector.block();
150 vtkDataObject* blockDO = output->GetBlock(blockNo);
151 vtkMultiBlockDataSet* block = vtkMultiBlockDataSet::SafeDownCast(blockDO);
154 return block->GetNumberOfBlocks();
161 void Foam::vtkPV3Foam::resetCounters()
163 // Reset region ids and sizes
164 selectInfoVolume_.reset();
165 selectInfoPatches_.reset();
166 selectInfoLagrangian_.reset();
167 selectInfoCellZones_.reset();
168 selectInfoFaceZones_.reset();
169 selectInfoPointZones_.reset();
170 selectInfoCellSets_.reset();
171 selectInfoFaceSets_.reset();
172 selectInfoPointSets_.reset();
176 bool Foam::vtkPV3Foam::setTime(const double& requestedTime)
180 Info<< "<beg> Foam::vtkPV3Foam::setTime(" << requestedTime << ")"
184 Time& runTime = dbPtr_();
187 instantList Times = runTime.times();
190 int nearestIndex = Time::findClosestTimeIndex(Times, requestedTime);
192 if (nearestIndex == -1)
202 // see what has changed
203 if (timeIndex_ != nearestIndex)
205 timeIndex_ = nearestIndex;
206 runTime.setTime(Times[nearestIndex], nearestIndex);
208 // the fields change each time
209 fieldsChanged_ = true;
213 if (meshPtr_->readUpdate() != polyMesh::UNCHANGED)
216 // patches, zones etc might have changed
228 Info<< "<end> Foam::vtkPV3Foam::setTime() - selected time "
229 << Times[nearestIndex].name() << " index=" << nearestIndex
230 << " meshChanged=" << meshChanged_
231 << " fieldsChanged=" << fieldsChanged_ << endl;
238 void Foam::vtkPV3Foam::updateSelectedRegions()
242 Info<< "<beg> Foam::vtkPV3Foam::updateSelectedRegions" << endl;
245 vtkDataArraySelection* arraySelection = reader_->GetRegionSelection();
247 const label nSelect = arraySelection->GetNumberOfArrays();
249 if (selectedRegions_.size() != nSelect)
251 selectedRegions_.setSize(nSelect);
252 selectedRegions_ = 0;
256 selectedRegionDatasetIds_.setSize(nSelect);
258 // Read the selected cell regions, zones, patches and add to region list
259 forAll (selectedRegions_, regionId)
261 int setting = arraySelection->GetArraySetting(regionId);
263 if (selectedRegions_[regionId] != setting)
265 selectedRegions_[regionId] = setting;
269 selectedRegionDatasetIds_[regionId] = -1;
273 Info<< " region[" << regionId << "] = "
274 << selectedRegions_[regionId]
275 << " : " << arraySelection->GetArrayName(regionId) << endl;
280 Info<< "<end> Foam::vtkPV3Foam::updateSelectedRegions" << endl;
285 Foam::stringList Foam::vtkPV3Foam::getSelectedArrayEntries
287 vtkDataArraySelection* arraySelection,
291 stringList selections(arraySelection->GetNumberOfArrays());
297 forAll (selections, elemI)
299 Info<< " \"" << arraySelection->GetArrayName(elemI) << "\"";
305 forAll (selections, elemI)
307 if (arraySelection->GetArraySetting(elemI))
311 selections[nElem] = getFirstWord
313 arraySelection->GetArrayName(elemI)
318 selections[nElem] = arraySelection->GetArrayName(elemI);
323 Info<< " " << selections[nElem];
335 selections.setSize(nElem);
340 Foam::stringList Foam::vtkPV3Foam::getSelectedArrayEntries
342 vtkDataArraySelection* arraySelection,
343 const selectionInfo& selector,
347 stringList selections(selector.size());
355 int elemI = selector.start();
356 elemI < selector.end();
360 Info<< " \"" << arraySelection->GetArrayName(elemI) << "\"";
369 int elemI = selector.start();
370 elemI < selector.end();
374 if (arraySelection->GetArraySetting(elemI))
378 selections[nElem] = getFirstWord
380 arraySelection->GetArrayName(elemI)
385 selections[nElem] = arraySelection->GetArrayName(elemI);
390 Info<< " " << selections[nElem];
402 selections.setSize(nElem);
407 void Foam::vtkPV3Foam::setSelectedArrayEntries
409 vtkDataArraySelection* arraySelection,
410 const stringList& selections
415 Info<< "<beg> Foam::vtkPV3Foam::setSelectedArrayEntries" << endl;
417 const label nEntries = arraySelection->GetNumberOfArrays();
419 // Reset all current entries to 'not selected'
420 arraySelection->DisableAllArrays();
422 // Loop through entries, setting values from selectedEntries
423 forAll (selections, elemI)
427 Info<< "selections[" << elemI << "] = " << selections[elemI]
431 for (label i=0; i<nEntries; i++)
433 string arrayName = arraySelection->GetArrayName(i);
435 if (arrayName == selections[elemI])
439 Info<< "enabling array: " << arrayName << " Index = "
444 arraySelection->EnableArray(arrayName.c_str());
451 Info<< "<end> Foam::vtkPV3Foam::setSelectedArrayEntries" << endl;
456 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
458 Foam::vtkPV3Foam::vtkPV3Foam
460 const char* const FileName,
461 vtkPV3FoamReader* reader
467 selectInfoVolume_(VOLUME, "unzoned"),
468 selectInfoPatches_(PATCHES, "patches"),
469 selectInfoLagrangian_(LAGRANGIAN, "lagrangian"),
470 selectInfoCellZones_(CELLZONE, "cellZone"),
471 selectInfoFaceZones_(FACEZONE, "faceZone"),
472 selectInfoPointZones_(POINTZONE, "pointZone"),
473 selectInfoCellSets_(CELLSET, "cellSet"),
474 selectInfoFaceSets_(FACESET, "faceSet"),
475 selectInfoPointSets_(POINTSET, "pointSet"),
476 patchTextActorsPtrs_(0),
484 Info<< "Foam::vtkPV3Foam::vtkPV3Foam - " << FileName << endl;
488 // avoid argList and get rootPath/caseName directly from the file
489 fileName fullCasePath(fileName(FileName).path());
491 if (!dir(fullCasePath))
495 if (fullCasePath == ".")
497 fullCasePath = cwd();
500 // Set the case as an environment variable - some BCs might use this
501 if (fullCasePath.name().find("processor", 0) == 0)
503 setEnv("FOAM_CASE", fullCasePath.path(), true);
507 setEnv("FOAM_CASE", fullCasePath, true);
512 Info<< "fullCasePath=" << fullCasePath << nl
513 << "FOAM_CASE=" << getEnv("FOAM_CASE") << endl;
516 // Create time object
521 Time::controlDictName,
522 fileName(fullCasePath.path()),
523 fileName(fullCasePath.name())
527 dbPtr_().functionObjects().off();
529 // Set initial cloud name
530 // TODO - TEMPORARY MEASURE UNTIL CAN PROCESS MULTIPLE CLOUDS
537 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
539 Foam::vtkPV3Foam::~vtkPV3Foam()
543 Info<< "<end> Foam::vtkPV3Foam::~vtkPV3Foam" << endl;
554 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
556 void Foam::vtkPV3Foam::UpdateInformation()
560 Info<< "<beg> Foam::vtkPV3Foam::UpdateInformation"
561 << " [meshPtr=" << (meshPtr_ ? "set" : "NULL") << "] TimeStep="
562 << reader_->GetTimeStep() << endl;
567 vtkDataArraySelection* arraySelection = reader_->GetRegionSelection();
569 stringList selectedEntries;
570 // enable 'internalMesh' on the first call
571 if (arraySelection->GetNumberOfArrays() == 0 && !meshPtr_)
573 selectedEntries.setSize(1);
574 selectedEntries[0] = "internalMesh";
578 // preserve the currently selected values
579 selectedEntries = getSelectedArrayEntries
585 // Clear current region list/array
586 arraySelection->RemoveAllArrays();
588 // Update region array
589 updateInformationInternalMesh();
590 updateInformationLagrangian();
591 updateInformationPatches();
593 if (reader_->GetIncludeSets())
595 updateInformationSets();
598 if (reader_->GetIncludeZones())
600 updateInformationZones();
603 // restore the currently enabled values
604 setSelectedArrayEntries
612 fieldsChanged_ = true;
615 // Update volField array
616 updateInformationFields<fvPatchField, volMesh>
618 reader_->GetVolFieldSelection()
621 // Update pointField array
622 updateInformationFields<pointPatchField, pointMesh>
624 reader_->GetPointFieldSelection()
627 // Update lagrangian field array
628 updateInformationLagrangianFields();
632 Info<< "<end> Foam::vtkPV3Foam::UpdateInformation" << endl;
638 void Foam::vtkPV3Foam::Update
640 vtkMultiBlockDataSet* output
645 cout<< "<beg> Foam::vtkPV3Foam::Update" << nl
649 cout<<"Internally:\n";
652 cout<< " has " << output->GetNumberOfBlocks() << " blocks\n";
656 // Set up region selection(s)
657 updateSelectedRegions();
659 // Update the Foam mesh
661 reader_->UpdateProgress(0.2);
664 convertMeshVolume(output);
665 convertMeshLagrangian(output);
666 convertMeshPatches(output);
667 reader_->UpdateProgress(0.4);
669 if (reader_->GetIncludeZones())
671 convertMeshCellZones(output);
672 convertMeshFaceZones(output);
673 convertMeshPointZones(output);
676 if (reader_->GetIncludeSets())
678 convertMeshCellSets(output);
679 convertMeshFaceSets(output);
680 convertMeshPointSets(output);
682 reader_->UpdateProgress(0.8);
685 updateVolFields(output);
686 updatePointFields(output);
687 updateLagrangianFields(output);
688 reader_->UpdateProgress(1.0);
692 Info<< "Number of data sets after update" << nl
694 << GetNumberOfDataSets(output, selectInfoVolume_) << nl
696 << GetNumberOfDataSets(output, selectInfoPatches_) << nl
698 << GetNumberOfDataSets(output, selectInfoLagrangian_) << nl
700 << GetNumberOfDataSets(output, selectInfoCellZones_) << nl
702 << GetNumberOfDataSets(output, selectInfoFaceZones_) << nl
704 << GetNumberOfDataSets(output, selectInfoPointZones_) << nl
706 << GetNumberOfDataSets(output, selectInfoCellSets_) << nl
708 << GetNumberOfDataSets(output, selectInfoFaceSets_) << nl
710 << GetNumberOfDataSets(output, selectInfoPointSets_) << nl;
713 cout<< "nBlocks = " << output->GetNumberOfBlocks() << "\n";
714 cout<< "done Update\n";
716 cout<< " has " << output->GetNumberOfBlocks() << " blocks\n";
717 output->GetInformation()->Print(cout);
719 cout<<"ShouldIReleaseData :" << output->ShouldIReleaseData() << "\n";
723 meshChanged_ = fieldsChanged_ = false;
727 double* Foam::vtkPV3Foam::findTimes(int& nTimeSteps)
730 double* tsteps = NULL;
734 Time& runTime = dbPtr_();
735 instantList timeLst = runTime.times();
737 // always skip "constant" time, unless there are no other times
738 nTimes = timeLst.size();
749 tsteps = new double[nTimes];
750 for (label stepI = 0; stepI < nTimes; ++stepI, ++timeI)
752 tsteps[stepI] = timeLst[timeI].value();
760 cout<< "no valid dbPtr:\n";
764 // vector length returned via the parameter
771 void Foam::vtkPV3Foam::addPatchNames(vtkRenderer* renderer)
773 // Remove any patch names previously added to the renderer
774 removePatchNames(renderer);
778 Info<< "<beg> Foam::vtkPV3Foam::addPatchNames" << endl;
781 const polyBoundaryMesh& pbMesh = meshPtr_->boundaryMesh();
783 const selectionInfo& selector = selectInfoPatches_;
785 // the currently selected patches, strip off any suffix
786 const stringList selectedPatches = getSelectedArrayEntries
788 reader_->GetRegionSelection(),
795 Info<<"... add patches: " << selectedPatches << endl;
798 // Find the total number of zones
799 // Each zone will take the patch name
801 // Number of zones per patch ... zero zones should be skipped
802 labelList nZones(pbMesh.size(), 0);
804 // Per global zone number the average face centre position
805 DynamicList<point> zoneCentre(pbMesh.size());
809 Info<< "... determining patch zones" << endl;
812 // Loop through all patches to determine zones, and centre of each zone
813 forAll(pbMesh, patchI)
815 const polyPatch& pp = pbMesh[patchI];
817 // Only include the patch if it is selected
818 bool isSelected = false;
819 forAll(selectedPatches, elemI)
821 if (pp.name() == selectedPatches[elemI])
830 const labelListList& edgeFaces = pp.edgeFaces();
831 const vectorField& n = pp.faceNormals();
833 boolList featEdge(pp.nEdges(), false);
835 forAll(edgeFaces, edgeI)
837 const labelList& eFaces = edgeFaces[edgeI];
839 if (eFaces.size() != 2)
841 featEdge[edgeI] = true;
843 else if (mag(n[eFaces[0]] & n[eFaces[1]]) < 0.5)
845 featEdge[edgeI] = true;
849 // Do topological analysis of patch. Determine disconnected regions
850 patchZones pZones(pp, featEdge);
852 nZones[patchI] = pZones.nZones();
854 labelList zoneNFaces(pZones.nZones(), 0);
856 // Save start of information for current patch
857 label patchStart = zoneCentre.size();
859 // Create storage for additional zone centres
860 forAll(zoneNFaces, zoneI)
862 zoneCentre.append(vector::zero);
865 // Do averaging per individual zone
869 label zoneI = pZones[faceI];
870 zoneCentre[patchStart+zoneI] += pp[faceI].centre(pp.points());
874 for (label i=0; i<nZones[patchI]; i++)
876 zoneCentre[patchStart + i] /= zoneNFaces[i];
884 Info<< "patch zone centres = " << zoneCentre << nl
885 << "zones per patch = " << nZones << endl;
888 // Set the size of the patch labels to max number of zones
889 patchTextActorsPtrs_.setSize(zoneCentre.size());
893 Info<< "constructing patch labels" << endl;
896 label globalZoneI = 0;
897 forAll(pbMesh, patchI)
899 const polyPatch& pp = pbMesh[patchI];
901 // Only selected patches will have a non-zero number of zones
902 for (label i=0; i<nZones[patchI]; i++)
906 Info<< "patch name = " << pp.name() << nl
907 << "anchor = " << zoneCentre[globalZoneI] << nl
908 << "globalZoneI = " << globalZoneI << endl;
911 vtkTextActor* txt = vtkTextActor::New();
913 txt->SetInput(pp.name().c_str());
915 // Set text properties
916 vtkTextProperty* tprop = txt->GetTextProperty();
917 tprop->SetFontFamilyToArial();
920 tprop->SetLineSpacing(1.0);
921 tprop->SetFontSize(12);
922 tprop->SetColor(1.0, 0.0, 0.0);
923 tprop->SetJustificationToCentered();
925 // Set text to use 3-D world co-ordinates
926 txt->GetPositionCoordinate()->SetCoordinateSystemToWorld();
928 txt->GetPositionCoordinate()->SetValue
930 zoneCentre[globalZoneI].x(),
931 zoneCentre[globalZoneI].y(),
932 zoneCentre[globalZoneI].z()
935 // Add text to each renderer
936 renderer->AddViewProp(txt);
938 // Maintain a list of text labels added so that they can be
940 patchTextActorsPtrs_[globalZoneI] = txt;
946 // Resize the patch names list to the actual number of patch names added
947 patchTextActorsPtrs_.setSize(globalZoneI);
951 Info<< "<end> Foam::vtkPV3Foam::addPatchNames)" << endl;
956 void Foam::vtkPV3Foam::removePatchNames(vtkRenderer* renderer)
960 Info<< "removePatchNames()" << endl;
963 forAll(patchTextActorsPtrs_, patchI)
965 renderer->RemoveViewProp(patchTextActorsPtrs_[patchI]);
966 patchTextActorsPtrs_[patchI]->Delete();
968 patchTextActorsPtrs_.setSize(0);
972 void Foam::vtkPV3Foam::PrintSelf(ostream& os, vtkIndent indent) const
974 os << indent << "Number of meshes: " << nMesh_ << "\n";
975 os << indent << "Number of nodes: "
976 << (meshPtr_ ? meshPtr_->nPoints() : 0) << "\n";
978 os << indent << "Number of cells: "
979 << (meshPtr_ ? meshPtr_->nCells() : 0) << "\n";
981 os << indent << "Number of available time steps: "
982 << (dbPtr_.valid() ? dbPtr_().times().size() : 0) << endl;
986 // parse these bits of info from /proc/meminfo (Linux)
988 // MemTotal: 2062660 kB
989 // MemFree: 1124400 kB
991 // used = MemTotal - MemFree is what the free(1) uses.
993 void Foam::vtkPV3Foam::printMemory()
995 const char* meminfo = "/proc/meminfo";
999 IFstream is(meminfo);
1005 while (is.getLine(line).good())
1010 if (sscanf(line.c_str(), "%30s %d", tag, &value) == 2)
1012 if (!strcmp(tag, "MemTotal:"))
1016 else if (!strcmp(tag, "MemFree:"))
1023 Info << "memUsed: " << (memTotal - memFree) << " kB\n";
1027 // ************************************************************************* //