initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / applications / utilities / postProcessing / graphics / PV3FoamReader / vtkPV3Foam / vtkPV3Foam.C
blob8af700450c0e31c10a2a2544fae5603bd61b9992
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
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
19     for more details.
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"
29 // Foam includes
30 #include "Time.H"
31 #include "fvMesh.H"
32 #include "IOobjectList.H"
33 #include "patchZones.H"
34 #include "vtkPV3FoamReader.h"
35 #include "IFstream.H"
37 // VTK includes
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,
65     vtkDataSet* dataset,
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)
74     {
75         FatalErrorIn("Foam::vtkPV3Foam::AddToBlock")
76             << "Block already has a vtkDataSet assigned to it" << nl << endl;
77         return;
78     }
80     if (!block)
81     {
82         block = vtkMultiBlockDataSet::New();
83         output->SetBlock(blockNo, block);
84         block->Delete();
85     }
87     if (block)
88     {
89         if (debug)
90         {
91             Info<< "block[" << blockNo << "] has "
92                 << block->GetNumberOfBlocks()
93                 <<  " datasets prior to adding set " << datasetNo
94                 <<  " with name: " << blockName << endl;
95         }
97         // when assigning dataset 0, also name the parent block
98         if (!datasetNo && selector.name())
99         {
100             output->GetMetaData(blockNo)->Set
101             (
102                 vtkCompositeDataSet::NAME(),
103                 selector.name()
104             );
105         }
106     }
109     block->SetBlock(datasetNo, dataset);
111     if (blockName.size())
112     {
113         block->GetMetaData(datasetNo)->Set
114         (
115             vtkCompositeDataSet::NAME(), blockName.c_str()
116         );
117     }
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);
133     if (block)
134     {
135         return vtkDataSet::SafeDownCast(block->GetBlock(datasetNo));
136     }
138     return 0;
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);
152     if (block)
153     {
154         return block->GetNumberOfBlocks();
155     }
157     return 0;
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)
178     if (debug)
179     {
180         Info<< "<beg> Foam::vtkPV3Foam::setTime(" << requestedTime << ")"
181             << endl;
182     }
184     Time& runTime = dbPtr_();
186     // Get times list
187     instantList Times = runTime.times();
189     bool found = false;
190     int nearestIndex = Time::findClosestTimeIndex(Times, requestedTime);
192     if (nearestIndex == -1)
193     {
194         nearestIndex = 0;
195         found = false;
196     }
197     else
198     {
199         found = true;
200     }
202     // see what has changed
203     if (timeIndex_ != nearestIndex)
204     {
205         timeIndex_ = nearestIndex;
206         runTime.setTime(Times[nearestIndex], nearestIndex);
208         // the fields change each time
209         fieldsChanged_ = true;
211         if (meshPtr_)
212         {
213             if (meshPtr_->readUpdate() != polyMesh::UNCHANGED)
214             {
215                 meshChanged_ = true;
216                 // patches, zones etc might have changed
217                 UpdateInformation();
218             }
219         }
220         else
221         {
222             meshChanged_ = true;
223         }
224     }
226     if (debug)
227     {
228         Info<< "<end> Foam::vtkPV3Foam::setTime() - selected time "
229             << Times[nearestIndex].name() << " index=" << nearestIndex
230             << " meshChanged=" << meshChanged_
231             << " fieldsChanged=" << fieldsChanged_ << endl;
232     }
234     return found;
238 void Foam::vtkPV3Foam::updateSelectedRegions()
240     if (debug)
241     {
242         Info<< "<beg> Foam::vtkPV3Foam::updateSelectedRegions" << endl;
243     }
245     vtkDataArraySelection* arraySelection = reader_->GetRegionSelection();
247     const label nSelect = arraySelection->GetNumberOfArrays();
249     if (selectedRegions_.size() != nSelect)
250     {
251         selectedRegions_.setSize(nSelect);
252         selectedRegions_ = 0;
253         meshChanged_ = true;
254     }
256     selectedRegionDatasetIds_.setSize(nSelect);
258     // Read the selected cell regions, zones, patches and add to region list
259     forAll (selectedRegions_, regionId)
260     {
261         int setting = arraySelection->GetArraySetting(regionId);
263         if (selectedRegions_[regionId] != setting)
264         {
265             selectedRegions_[regionId] = setting;
266             meshChanged_ = true;
267         }
269         selectedRegionDatasetIds_[regionId] = -1;
271         if (debug)
272         {
273             Info<< "  region[" << regionId << "] = "
274                 << selectedRegions_[regionId]
275                 << " : " << arraySelection->GetArrayName(regionId) << endl;
276         }
277     }
278     if (debug)
279     {
280         Info<< "<end> Foam::vtkPV3Foam::updateSelectedRegions" << endl;
281     }
285 Foam::stringList Foam::vtkPV3Foam::getSelectedArrayEntries
287     vtkDataArraySelection* arraySelection,
288     const bool firstWord
291     stringList selections(arraySelection->GetNumberOfArrays());
292     label nElem = 0;
294     if (debug)
295     {
296         Info<< "available(";
297         forAll (selections, elemI)
298         {
299             Info<< " \"" << arraySelection->GetArrayName(elemI) << "\"";
300         }
301         Info<< " )\n"
302             << "selected(";
303     }
305     forAll (selections, elemI)
306     {
307         if (arraySelection->GetArraySetting(elemI))
308         {
309             if (firstWord)
310             {
311                 selections[nElem] = getFirstWord
312                 (
313                     arraySelection->GetArrayName(elemI)
314                 );
315             }
316             else
317             {
318                 selections[nElem] = arraySelection->GetArrayName(elemI);
319             }
321             if (debug)
322             {
323                 Info<< " " << selections[nElem];
324             }
326             ++nElem;
327         }
328     }
330     if (debug)
331     {
332         Info<< " )" << endl;
333     }
335     selections.setSize(nElem);
336     return selections;
340 Foam::stringList Foam::vtkPV3Foam::getSelectedArrayEntries
342     vtkDataArraySelection* arraySelection,
343     const selectionInfo& selector,
344     const bool firstWord
347     stringList selections(selector.size());
348     label nElem = 0;
350     if (debug)
351     {
352         Info<< "available(";
353         for
354         (
355             int elemI = selector.start();
356             elemI < selector.end();
357             ++elemI
358         )
359         {
360             Info<< " \"" << arraySelection->GetArrayName(elemI) << "\"";
361         }
363         Info<< " )\n"
364             << "selected(";
365     }
367     for
368     (
369         int elemI = selector.start();
370         elemI < selector.end();
371         ++elemI
372     )
373     {
374         if (arraySelection->GetArraySetting(elemI))
375         {
376             if (firstWord)
377             {
378                 selections[nElem] = getFirstWord
379                 (
380                     arraySelection->GetArrayName(elemI)
381                 );
382             }
383             else
384             {
385                 selections[nElem] = arraySelection->GetArrayName(elemI);
386             }
388             if (debug)
389             {
390                 Info<< " " << selections[nElem];
391             }
393             ++nElem;
394         }
395     }
397     if (debug)
398     {
399         Info<< " )" << endl;
400     }
402     selections.setSize(nElem);
403     return selections;
407 void Foam::vtkPV3Foam::setSelectedArrayEntries
409     vtkDataArraySelection* arraySelection,
410     const stringList& selections
413     if (debug > 1)
414     {
415         Info<< "<beg> Foam::vtkPV3Foam::setSelectedArrayEntries" << endl;
416     }
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)
424     {
425         if (debug > 1)
426         {
427             Info<< "selections[" << elemI << "] = " << selections[elemI]
428                 << endl;
429         }
431         for (label i=0; i<nEntries; i++)
432         {
433             string arrayName = arraySelection->GetArrayName(i);
435             if (arrayName == selections[elemI])
436             {
437                 if (debug > 1)
438                 {
439                     Info<< "enabling array: " << arrayName << " Index = "
440                         << i
441                         << endl;
442                 }
444                 arraySelection->EnableArray(arrayName.c_str());
445                 break;
446             }
447         }
448     }
449     if (debug > 1)
450     {
451         Info<< "<end> Foam::vtkPV3Foam::setSelectedArrayEntries" << endl;
452     }
456 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
458 Foam::vtkPV3Foam::vtkPV3Foam
460     const char* const FileName,
461     vtkPV3FoamReader* reader
464     reader_(reader),
465     dbPtr_(NULL),
466     meshPtr_(NULL),
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),
477     nMesh_(0),
478     timeIndex_(-1),
479     meshChanged_(true),
480     fieldsChanged_(true)
482     if (debug)
483     {
484         Info<< "Foam::vtkPV3Foam::vtkPV3Foam - " << FileName << endl;
485         printMemory();
486     }
488     // avoid argList and get rootPath/caseName directly from the file
489     fileName fullCasePath(fileName(FileName).path());
491     if (!dir(fullCasePath))
492     {
493         return;
494     }
495     if (fullCasePath == ".")
496     {
497         fullCasePath = cwd();
498     }
500     // Set the case as an environment variable - some BCs might use this
501     if (fullCasePath.name().find("processor", 0) == 0)
502     {
503         setEnv("FOAM_CASE", fullCasePath.path(), true);
504     }
505     else
506     {
507         setEnv("FOAM_CASE", fullCasePath, true);
508     }
510     if (debug)
511     {
512         Info<< "fullCasePath=" << fullCasePath << nl
513             << "FOAM_CASE=" << getEnv("FOAM_CASE") << endl;
514     }
516     // Create time object
517     dbPtr_.reset
518     (
519         new Time
520         (
521             Time::controlDictName,
522             fileName(fullCasePath.path()),
523             fileName(fullCasePath.name())
524         )
525     );
527     dbPtr_().functionObjects().off();
529     // Set initial cloud name
530     // TODO - TEMPORARY MEASURE UNTIL CAN PROCESS MULTIPLE CLOUDS
531     cloudName_ = "";
533     UpdateInformation();
537 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
539 Foam::vtkPV3Foam::~vtkPV3Foam()
541     if (debug)
542     {
543         Info<< "<end> Foam::vtkPV3Foam::~vtkPV3Foam" << endl;
544     }
546     if (meshPtr_)
547     {
548         delete meshPtr_;
549         meshPtr_ = NULL;
550     }
554 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
556 void Foam::vtkPV3Foam::UpdateInformation()
558     if (debug)
559     {
560         Info<< "<beg> Foam::vtkPV3Foam::UpdateInformation"
561             << " [meshPtr=" << (meshPtr_ ? "set" : "NULL") << "] TimeStep="
562             << reader_->GetTimeStep() << endl;
563     }
565     resetCounters();
567     vtkDataArraySelection* arraySelection = reader_->GetRegionSelection();
569     stringList selectedEntries;
570     // enable 'internalMesh' on the first call
571     if (arraySelection->GetNumberOfArrays() == 0 && !meshPtr_)
572     {
573         selectedEntries.setSize(1);
574         selectedEntries[0] = "internalMesh";
575     }
576     else
577     {
578         // preserve the currently selected values
579         selectedEntries = getSelectedArrayEntries
580         (
581             arraySelection
582         );
583     }
585     // Clear current region list/array
586     arraySelection->RemoveAllArrays();
588     // Update region array
589     updateInformationInternalMesh();
590     updateInformationLagrangian();
591     updateInformationPatches();
593     if (reader_->GetIncludeSets())
594     {
595         updateInformationSets();
596     }
598     if (reader_->GetIncludeZones())
599     {
600         updateInformationZones();
601     }
603     // restore the currently enabled values
604     setSelectedArrayEntries
605     (
606         arraySelection,
607         selectedEntries
608     );
610     if (meshChanged_)
611     {
612         fieldsChanged_ = true;
613     }
615     // Update volField array
616     updateInformationFields<fvPatchField, volMesh>
617     (
618         reader_->GetVolFieldSelection()
619     );
621     // Update pointField array
622     updateInformationFields<pointPatchField, pointMesh>
623     (
624         reader_->GetPointFieldSelection()
625     );
627     // Update lagrangian field array
628     updateInformationLagrangianFields();
630     if (debug)
631     {
632         Info<< "<end> Foam::vtkPV3Foam::UpdateInformation" << endl;
633     }
638 void Foam::vtkPV3Foam::Update
640     vtkMultiBlockDataSet* output
643     if (debug)
644     {
645         cout<< "<beg> Foam::vtkPV3Foam::Update" << nl
646             <<"Update\n";
647         output->Print(cout);
649         cout<<"Internally:\n";
650         output->Print(cout);
652         cout<< " has " << output->GetNumberOfBlocks() << " blocks\n";
653         printMemory();
654     }
656     // Set up region selection(s)
657     updateSelectedRegions();
659     // Update the Foam mesh
660     updateFoamMesh();
661     reader_->UpdateProgress(0.2);
663     // Convert meshes
664     convertMeshVolume(output);
665     convertMeshLagrangian(output);
666     convertMeshPatches(output);
667     reader_->UpdateProgress(0.4);
669     if (reader_->GetIncludeZones())
670     {
671         convertMeshCellZones(output);
672         convertMeshFaceZones(output);
673         convertMeshPointZones(output);
674     }
676     if (reader_->GetIncludeSets())
677     {
678         convertMeshCellSets(output);
679         convertMeshFaceSets(output);
680         convertMeshPointSets(output);
681     }
682     reader_->UpdateProgress(0.8);
684     // Update fields
685     updateVolFields(output);
686     updatePointFields(output);
687     updateLagrangianFields(output);
688     reader_->UpdateProgress(1.0);
690     if (debug)
691     {
692         Info<< "Number of data sets after update" << nl
693             << "  VOLUME = "
694             << GetNumberOfDataSets(output, selectInfoVolume_) << nl
695             << "  PATCHES = "
696             << GetNumberOfDataSets(output, selectInfoPatches_) << nl
697             << "  LAGRANGIAN = "
698             << GetNumberOfDataSets(output, selectInfoLagrangian_) << nl
699             << "  CELLZONE = "
700             << GetNumberOfDataSets(output, selectInfoCellZones_) << nl
701             << "  FACEZONE = "
702             << GetNumberOfDataSets(output, selectInfoFaceZones_) << nl
703             << "  POINTZONE = "
704             << GetNumberOfDataSets(output, selectInfoPointZones_) << nl
705             << "  CELLSET = "
706             << GetNumberOfDataSets(output, selectInfoCellSets_) << nl
707             << "  FACESET = "
708             << GetNumberOfDataSets(output, selectInfoFaceSets_) << nl
709             << "  POINTSET = "
710             << GetNumberOfDataSets(output, selectInfoPointSets_) << nl;
712         // traverse blocks:
713         cout<< "nBlocks = " << output->GetNumberOfBlocks() << "\n";
714         cout<< "done Update\n";
715         output->Print(cout);
716         cout<< " has " << output->GetNumberOfBlocks() << " blocks\n";
717         output->GetInformation()->Print(cout);
719         cout<<"ShouldIReleaseData :" << output->ShouldIReleaseData() << "\n";
720         printMemory();
721     }
723     meshChanged_ = fieldsChanged_ = false;
727 double* Foam::vtkPV3Foam::findTimes(int& nTimeSteps)
729     int nTimes = 0;
730     double* tsteps = NULL;
732     if (dbPtr_.valid())
733     {
734         Time& runTime = dbPtr_();
735         instantList timeLst = runTime.times();
737         // always skip "constant" time, unless there are no other times
738         nTimes = timeLst.size();
739         label timeI = 0;
741         if (nTimes > 1)
742         {
743             timeI = 1;
744             --nTimes;
745         }
747         if (nTimes)
748         {
749             tsteps = new double[nTimes];
750             for (label stepI = 0; stepI < nTimes; ++stepI, ++timeI)
751             {
752                 tsteps[stepI] = timeLst[timeI].value();
753             }
754         }
755     }
756     else
757     {
758         if (debug)
759         {
760             cout<< "no valid dbPtr:\n";
761         }
762     }
764     // vector length returned via the parameter
765     nTimeSteps = nTimes;
767     return tsteps;
771 void Foam::vtkPV3Foam::addPatchNames(vtkRenderer* renderer)
773     // Remove any patch names previously added to the renderer
774     removePatchNames(renderer);
776     if (debug)
777     {
778         Info<< "<beg> Foam::vtkPV3Foam::addPatchNames" << endl;
779     }
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
787     (
788         reader_->GetRegionSelection(),
789         selector,
790         true
791     );
793     if (debug)
794     {
795         Info<<"... add patches: " << selectedPatches << endl;
796     }
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());
807     if (debug)
808     {
809         Info<< "... determining patch zones" << endl;
810     }
812     // Loop through all patches to determine zones, and centre of each zone
813     forAll(pbMesh, patchI)
814     {
815         const polyPatch& pp = pbMesh[patchI];
817         // Only include the patch if it is selected
818         bool isSelected = false;
819         forAll(selectedPatches, elemI)
820         {
821             if (pp.name() == selectedPatches[elemI])
822             {
823                 isSelected = true;
824                 break;
825             }
826         }
828         if (isSelected)
829         {
830             const labelListList& edgeFaces = pp.edgeFaces();
831             const vectorField& n = pp.faceNormals();
833             boolList featEdge(pp.nEdges(), false);
835             forAll(edgeFaces, edgeI)
836             {
837                 const labelList& eFaces = edgeFaces[edgeI];
839                 if (eFaces.size() != 2)
840                 {
841                     featEdge[edgeI] = true;
842                 }
843                 else if (mag(n[eFaces[0]] & n[eFaces[1]]) < 0.5)
844                 {
845                     featEdge[edgeI] = true;
846                 }
847             }
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)
861             {
862                 zoneCentre.append(vector::zero);
863             }
865             // Do averaging per individual zone
867             forAll(pp, faceI)
868             {
869                 label zoneI = pZones[faceI];
870                 zoneCentre[patchStart+zoneI] += pp[faceI].centre(pp.points());
871                 zoneNFaces[zoneI]++;
872             }
874             for (label i=0; i<nZones[patchI]; i++)
875             {
876                 zoneCentre[patchStart + i] /= zoneNFaces[i];
877             }
878         }
879     }
880     zoneCentre.shrink();
882     if (debug)
883     {
884         Info<< "patch zone centres = " << zoneCentre << nl
885             << "zones per patch = " << nZones << endl;
886     }
888     // Set the size of the patch labels to max number of zones
889     patchTextActorsPtrs_.setSize(zoneCentre.size());
891     if (debug)
892     {
893         Info<< "constructing patch labels" << endl;
894     }
896     label globalZoneI = 0;
897     forAll(pbMesh, patchI)
898     {
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++)
903         {
904             if (debug)
905             {
906                 Info<< "patch name = " << pp.name() << nl
907                     << "anchor = " << zoneCentre[globalZoneI] << nl
908                     << "globalZoneI = " << globalZoneI << endl;
909             }
911             vtkTextActor* txt = vtkTextActor::New();
913             txt->SetInput(pp.name().c_str());
915             // Set text properties
916             vtkTextProperty* tprop = txt->GetTextProperty();
917             tprop->SetFontFamilyToArial();
918             tprop->BoldOff();
919             tprop->ShadowOff();
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
929             (
930                 zoneCentre[globalZoneI].x(),
931                 zoneCentre[globalZoneI].y(),
932                 zoneCentre[globalZoneI].z()
933             );
935             // Add text to each renderer
936             renderer->AddViewProp(txt);
938             // Maintain a list of text labels added so that they can be
939             // removed later
940             patchTextActorsPtrs_[globalZoneI] = txt;
942             globalZoneI++;
943         }
944     }
946     // Resize the patch names list to the actual number of patch names added
947     patchTextActorsPtrs_.setSize(globalZoneI);
949     if (debug)
950     {
951         Info<< "<end> Foam::vtkPV3Foam::addPatchNames)" << endl;
952     }
956 void Foam::vtkPV3Foam::removePatchNames(vtkRenderer* renderer)
958     if (debug)
959     {
960         Info<< "removePatchNames()" << endl;
961     }
963     forAll(patchTextActorsPtrs_, patchI)
964     {
965         renderer->RemoveViewProp(patchTextActorsPtrs_[patchI]);
966         patchTextActorsPtrs_[patchI]->Delete();
967     }
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";
997     if (exists(meminfo))
998     {
999         IFstream is(meminfo);
1000         label memTotal = 0;
1001         label memFree = 0;
1003         string line;
1005         while (is.getLine(line).good())
1006         {
1007             char tag[32];
1008             int value;
1010             if (sscanf(line.c_str(), "%30s %d", tag, &value) == 2)
1011             {
1012                 if (!strcmp(tag, "MemTotal:"))
1013                 {
1014                     memTotal = value;
1015                 }
1016                 else if (!strcmp(tag, "MemFree:"))
1017                 {
1018                     memFree = value;
1019                 }
1020             }
1021         }
1023         Info << "memUsed: " << (memTotal - memFree) << " kB\n";
1024     }
1027 // ************************************************************************* //