initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / applications / utilities / postProcessing / graphics / PVFoamReader / vtkFoam / vtkFoam.C
blob6a1830b5d14dd60a851496247830f9a6942c983c
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 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 "vtkFoam.H"
29 #include "argList.H"
30 #include "Time.H"
31 #include "polyBoundaryMeshEntries.H"
32 #include "IOobjectList.H"
33 #include "wordList.H"
34 #include "fvMesh.H"
35 #include "volFields.H"
36 #include "pointMesh.H"
37 #include "volPointInterpolation.H"
39 #include "vtkFoamReader.h"
40 #include "vtkDataArraySelection.h"
41 #include "vtkUnstructuredGrid.h"
42 #include "vtkPointData.h"
43 #include "vtkCellData.h"
44 #include "vtkFloatArray.h"
45 #include "vtkCharArray.h"
47 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
49 defineTypeNameAndDebug(Foam::vtkFoam, 0);
51 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
53 #include "vtkFoamConvertFields.H"
55 void Foam::vtkFoam::SetName
57     vtkUnstructuredGrid* vtkMesh,
58     const char* name
61     vtkCharArray* nmArray =  vtkCharArray::New();
62     nmArray->SetName("Name");
63     size_t len = strlen(name);
64     nmArray->SetNumberOfTuples(static_cast<vtkIdType>(len)+1);
65     char* copy = nmArray->GetPointer(0);
66     memcpy(copy, name, len);
67     copy[len] = '\0';
68     vtkMesh->GetFieldData()->AddArray(nmArray);
69     nmArray->Delete();
73 Foam::string Foam::vtkFoam::padTimeString(const string& ts)
75     return ts + string("            ", max(label(12 - ts.size()), 0));
79 // Pad the patch name string in order to account for dynamic changes
80 // in patch names during topological changes
81 Foam::string Foam::vtkFoam::padPatchString(const string& ps)
83     label n = max(label(50 - ps.size()), 0);
84     return ps + string("                                                  ", n);
88 void Foam::vtkFoam::setSelectedTime
90     Time& runTime,
91     vtkFoamReader* reader
94     // Get times list
95     instantList Times = runTime.times();
96     int timeIndex = min(max(reader->GetTimeStep() + 1, 0), Times.size()-1);
98     // If this is the first call timeIndex will be 0 ("constant")
99     // so reset to the first time step if one exists and deselect every
100     // element of the selection array
101     if (timeIndex == 0)
102     {
103         timeIndex = min(1, Times.size()-1);
104         reader->GetTimeSelection()->DisableAllArrays();
105     }
107     label selectedTimeIndex = -1;
108     label nSelectedTimes = reader->GetTimeSelection()->GetNumberOfArrays();
110     for (label i=nSelectedTimes-1; i>=0; i--)
111     {
112         if(reader->GetTimeSelection()->GetArraySetting(i))
113         {
114             word timeName = string::validate<word>
115             (
116                 reader->GetTimeSelection()->GetArrayName(i)
117             );
119             forAll(Times, j)
120             {
121                 if (Times[j].name() == timeName)
122                 {
123                     selectedTimeIndex = j;
124                     break;
125                 }
126             }
127             break;
128         }
129     }
131     if (selectedTimeIndex != -1)
132     {
133         timeIndex = min(selectedTimeIndex, Times.size()-1);
134     }
136     if (debug)
137     {
138         Info<< "Selecting time " << Times[timeIndex].name() << endl;
139     }
141     runTime.setTime(Times[timeIndex], timeIndex);
143     Times = runTime.times();
145     reader->SetTimeStepRange(0, max(Times.size()-2, 0));
147     // reset the time steps ...
148     reader->GetTimeSelection()->RemoveAllArrays();
150     int* TimeStepLimits = reader->GetTimeStepLimits();
151     label maxStartTimes = min(Times.size(), TimeStepLimits[0]);
152     label maxNTimes = min(Times.size() - maxStartTimes, TimeStepLimits[1]);
154     for (label i=0; i<maxStartTimes; i++)
155     {
156         reader->GetTimeSelection()
157             ->AddArray(padTimeString(Times[i].name()).c_str());
158     }
160     if (Times.size() > TimeStepLimits[0] + TimeStepLimits[1])
161     {
162         reader->GetTimeSelection()->AddArray(padTimeString("...").c_str());
163     }
165     for (label i=Times.size() - maxNTimes; i<Times.size(); i++)
166     {
167         reader->GetTimeSelection()
168             ->AddArray(padTimeString(Times[i].name()).c_str());
169     }
171     // Disable all the time selections (which are all selected by default) ...
172     reader->GetTimeSelection()->DisableAllArrays();
174     // But maintain the selections made previously
175     if (selectedTimeIndex != -1 && selectedTimeIndex < Times.size())
176     {
177         reader->GetTimeSelection()->EnableArray
178             (padTimeString(Times[selectedTimeIndex].name()).c_str());
179     }
183 void Foam::vtkFoam::updateSelectedRegions()
185     if (debug)
186     {
187         Info<< "Foam::vtkFoam::updateSelectedRegions()" << endl;
188     }
190     label nRegions = reader_->GetRegionSelection()->GetNumberOfArrays();
192     selectedRegions_.setSize(nRegions);
194     // Read the selected patches and add to the region list
195     for (int i=0; i<nRegions; i++)
196     {
197         selectedRegions_[i] =
198             reader_->GetRegionSelection()->GetArraySetting(i);
199     }
203 void Foam::vtkFoam::convertMesh()
205     if (debug)
206     {
207         Info<< "Foam::vtkFoam::convertMesh()" << endl;
208     }
210     const fvMesh& mesh = *meshPtr_;
212     // Read the internal mesh as region 0 if selected
213     if (reader_->GetRegionSelection()->GetArraySetting(0))
214     {
215         selectedRegions_[0] = true;
216         addInternalMesh
217         (
218             mesh,
219             vtkUnstructuredGrid::SafeDownCast(reader_->GetOutput(0))
220         );
221     }
222     else
223     {
224         selectedRegions_[0] = false;
226         vtkUnstructuredGrid *vtkMesh =
227             vtkUnstructuredGrid::SafeDownCast(reader_->GetOutput(0));
229         vtkMesh->Initialize();
230         SetName(vtkMesh, "(Internal Mesh)");
231     }
234     // Read the selected patches and add to the region list
236     polyBoundaryMeshEntries patchEntries
237     (
238         IOobject
239         (
240             "boundary",
241             dbPtr_().findInstance(polyMesh::meshSubDir, "boundary"),
242             polyMesh::meshSubDir,
243             dbPtr_(),
244             IOobject::MUST_READ,
245             IOobject::NO_WRITE,
246             false
247         )
248     );
250     label regioni = 0;
251     label regioniLast = 0;
253     // Read in the number Outputs (patch regions) currently being used
254     label currNOutputs = reader_->GetNumberOfOutputs();
256     // Cycle through all the patches in the boundary file for the relevant
257     // time step
258     forAll(patchEntries, entryi)
259     {
260         // Number of faces in the current patch (Used to detect dummy patches
261         // of size zero)
262         label nFaces(readLabel(patchEntries[entryi].dict().lookup("nFaces")));
264         // Check to see if the patch is currently a part of the displayed list
265         if
266         (
267             reader_->GetRegionSelection()->ArrayExists
268             (
269                 padPatchString(patchEntries[entryi].keyword()).c_str()
270             )
271         )
272         {
273             if  (!nFaces)
274             {
275                 // Remove patch if it is only a dummy patch in the current
276                 // time step with zero faces
277                 reader_->GetRegionSelection()->RemoveArrayByName
278                 (
279                     padPatchString(patchEntries[entryi].keyword()).c_str()
280                 );
281             }
282             else
283             {
284                 // A patch already existent in the list and which
285                 // continues to exist found
286                 regioni++;
287             }
288         }
289         else
290         {
291             // A new patch so far not yet included into the list has been found
292             if  (nFaces)
293             {
294                 regioni++;
296                 // Add a new entry to the list of regions
297                 reader_->GetRegionSelection()->AddArray
298                 (
299                     padPatchString(patchEntries[entryi].keyword()).c_str()
300                 );
302                 // AddArray automatically enables a new array... disable
303                 // it manually
304                 reader_->GetRegionSelection()->DisableArray
305                 (
306                     padPatchString(patchEntries[entryi].keyword()).c_str()
307                 );
308             }
309         }
311         // Avoid Initialization of the same Output twice
312         if (regioni != regioniLast)
313         {
314             // Only setup an Output if it has not been setup before
315             if(regioni >= currNOutputs)
316             {
317                 vtkUnstructuredGrid* ugrid = vtkUnstructuredGrid::New();
318                 reader_->SetNthOutput(regioni,ugrid);
319                 ugrid->Delete();
320             }
321             // Initialize -> Delete memory used, and reset to zero state
322             reader_->GetOutput(regioni)->Initialize();
323             regioniLast = regioni;
324         }
325     }
327     // Initialize (reset to zero and free) any outputs which are not used
328     // anymore
329     if (regioni < currNOutputs)
330     {
331         for(label i = (regioni+1); i < currNOutputs;i++)
332         {
333             reader_->GetOutput(i)->Initialize();
334         }
335     }
337     selectedRegions_.setSize(regioni + 1);
339     regioni = 0;
341     const polyBoundaryMesh& patches = mesh.boundaryMesh();
343     forAll (patches, patchi)
344     {
345         if (patches[patchi].size())
346         {
347             regioni++;
349             if (reader_->GetRegionSelection()->GetArraySetting(regioni))
350             {
351                 selectedRegions_[regioni] = true;
352                 addPatch
353                 (
354                     patches[patchi],
355                     vtkUnstructuredGrid::SafeDownCast
356                     (
357                         reader_->GetOutput(regioni)
358                     )
359                 );
360             }
361             else
362             {
363                 selectedRegions_[regioni] = false;
365                 vtkUnstructuredGrid *vtkMesh =
366                     vtkUnstructuredGrid::SafeDownCast
367                     (
368                         reader_->GetOutput(regioni)
369                     );
371                 vtkMesh->Initialize();
372                 SetName
373                 (
374                     vtkMesh,
375                     ('(' + padPatchString(patches[patchi].name()) + ')').c_str()
376                 );
377             }
378         }
379     }
383 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
385 Foam::vtkFoam::vtkFoam(const char* const FileName, vtkFoamReader* reader)
387     reader_(reader),
388     argsPtr_(NULL),
389     dbPtr_(NULL),
390     meshPtr_(NULL)
392     fileName fullCasePath(fileName(FileName).path());
394     if (!isDir(fullCasePath))
395     {
396         return;
397     }
399     char* argvStrings[3];
400     argvStrings[0] = new char[9];
401     strcpy(argvStrings[0], "/vtkFoam");
402     argvStrings[1] = new char[6];
403     strcpy(argvStrings[1], "-case");
404     argvStrings[2] = new char[fullCasePath.size()+1];
405     strcpy(argvStrings[2], fullCasePath.c_str());
407     int argc = 3;
408     char** argv = &argvStrings[0];
409     argsPtr_.reset(new argList(argc, argv));
411     for(int i = 0; i < argc; i++)
412     {
413         delete[] argvStrings[i];
414     }
416     dbPtr_.reset
417     (
418         new Time
419         (
420             Time::controlDictName,
421             argsPtr_().rootPath(),
422             argsPtr_().caseName()
423         )
424     );
425     dbPtr_().functionObjects().off();
426     setSelectedTime(dbPtr_(), reader_);
428     if (debug)
429     {
430         Info<< "vtkFoam::ExecuteInformation: Initialising outputs" << endl;
431     }
433     reader_->GetRegionSelection()->AddArray("Internal Mesh");
435     vtkUnstructuredGrid* ugrid = vtkUnstructuredGrid::New();
436     reader_->SetNthOutput(0, ugrid);
437     ugrid->Delete();
438     reader_->GetOutput(0)->Initialize();
440     polyBoundaryMeshEntries patchEntries
441     (
442         IOobject
443         (
444             "boundary",
445             dbPtr_().findInstance(polyMesh::meshSubDir, "boundary"),
446             polyMesh::meshSubDir,
447             dbPtr_(),
448             IOobject::MUST_READ,
449             IOobject::NO_WRITE,
450             false
451         )
452     );
454     label regioni = 0;
455     forAll(patchEntries, entryi)
456     {
457         label nFaces(readLabel(patchEntries[entryi].dict().lookup("nFaces")));
459         if (nFaces)
460         {
461             regioni++;
463             reader_->GetRegionSelection()->AddArray
464             (
465                 padPatchString(patchEntries[entryi].keyword()).c_str()
466             );
468             vtkUnstructuredGrid* ugrid = vtkUnstructuredGrid::New();
469             reader_->SetNthOutput(regioni, ugrid);
470             ugrid->Delete();
471             reader_->GetOutput(regioni)->Initialize();
472         }
473     }
475     selectedRegions_.setSize(regioni + 1);
476     selectedRegions_ = true;
478     UpdateInformation();
482 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
484 Foam::vtkFoam::~vtkFoam()
486     // Do NOT delete meshPtr_ since still referenced somehow.
490 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
492 #include "vtkFoamAddFields.H"
494 void Foam::vtkFoam::UpdateInformation()
496     if (debug)
497     {
498         Info<< "TimeStep = " << reader_->GetTimeStep() << endl;
499     }
501     setSelectedTime(dbPtr_(), reader_);
503     // Search for list of objects for this time
504     IOobjectList objects(dbPtr_(), dbPtr_().timeName());
506     addFields<volScalarField>(reader_->GetVolFieldSelection(), objects);
507     addFields<volVectorField>(reader_->GetVolFieldSelection(), objects);
508     addFields<volSphericalTensorField>(reader_->GetVolFieldSelection(), objects);
509     addFields<volSymmTensorField>(reader_->GetVolFieldSelection(), objects);
510     addFields<volTensorField>(reader_->GetVolFieldSelection(), objects);
512     addFields<pointScalarField>(reader_->GetPointFieldSelection(), objects);
513     addFields<pointVectorField>(reader_->GetPointFieldSelection(), objects);
514     addFields<pointSphericalTensorField>(reader_->GetPointFieldSelection(), objects);
515     addFields<pointSymmTensorField>(reader_->GetPointFieldSelection(), objects);
516     addFields<pointTensorField>(reader_->GetPointFieldSelection(), objects);
520 void Foam::vtkFoam::Update()
522     if
523     (
524         !reader_->GetCacheMesh()
525      || reader_->GetTimeSelection()->GetArraySetting(0)
526     )
527     {
528         meshPtr_= NULL;
529     }
531     // Clear the current set of selected fields
533     for (label i=0; i<reader_->GetNumberOfOutputs(); i++)
534     {
535         vtkUnstructuredGrid *vtkMesh =
536             vtkUnstructuredGrid::SafeDownCast(reader_->GetOutput(i));
538         vtkCellData* cellData = vtkMesh->GetCellData();
539         int numberOfCellArrays = cellData->GetNumberOfArrays();
541         wordList cellFieldNames(numberOfCellArrays);
542         for (int j=0; j<numberOfCellArrays; j++)
543         {
544             cellFieldNames[j] = cellData->GetArrayName(j);
545         }
547         for (int j=0; j<numberOfCellArrays; j++)
548         {
549             cellData->RemoveArray(cellFieldNames[j].c_str());
550         }
552         vtkPointData* pointData = vtkMesh->GetPointData();
553         int numberOfPointArrays = pointData->GetNumberOfArrays();
555         wordList pointFieldNames(numberOfPointArrays);
556         for (int j=0; j<numberOfPointArrays; j++)
557         {
558             pointFieldNames[j] = pointData->GetArrayName(j);
559         }
561         for (int j=0; j<numberOfPointArrays; j++)
562         {
563             pointData->RemoveArray(pointFieldNames[j].c_str());
564         }
565     }
567     // Check to see if the mesh has been created
569     if (!meshPtr_)
570     {
571         if (debug)
572         {
573             Info<< "Reading Mesh" << endl;
574         }
575         meshPtr_ =
576             new fvMesh
577             (
578                 IOobject
579                 (
580                     fvMesh::defaultRegion,
581                     dbPtr_().timeName(),
582                     dbPtr_()
583                 )
584             );
585         convertMesh();
586     }
587     else
588     {
589         boolList oldSelectedRegions = selectedRegions_;
590         updateSelectedRegions();
591         if
592         (
593             meshPtr_->readUpdate() != fvMesh::UNCHANGED
594          || oldSelectedRegions != selectedRegions_
595         )
596         {
597             convertMesh();
598         }
599     }
601     if (debug)
602     {
603         Info<< "converting fields" << endl;
604     }
606     const fvMesh& mesh = *meshPtr_;
608     // Construct interpolation on the raw mesh
609     Foam::pointMesh pMesh(mesh);
611     Foam::volPointInterpolation pInterp(mesh, pMesh);
613     // Search for list of objects for this time
614     Foam::IOobjectList objects(mesh, dbPtr_().timeName());
616     convertVolFields<Foam::scalar>
617     (
618         mesh, pInterp, objects, reader_->GetVolFieldSelection()
619     );
620     convertVolFields<Foam::vector>
621     (
622         mesh, pInterp, objects, reader_->GetVolFieldSelection()
623     );
624     convertVolFields<Foam::sphericalTensor>
625     (
626         mesh, pInterp, objects, reader_->GetVolFieldSelection()
627     );
628     convertVolFields<Foam::symmTensor>
629     (
630         mesh, pInterp, objects, reader_->GetVolFieldSelection()
631     );
632     convertVolFields<Foam::tensor>
633     (
634         mesh, pInterp, objects, reader_->GetVolFieldSelection()
635     );
637     convertPointFields<Foam::scalar>
638     (
639         mesh, objects, reader_->GetPointFieldSelection()
640     );
641     convertPointFields<Foam::vector>
642     (
643         mesh, objects, reader_->GetPointFieldSelection()
644     );
645     convertPointFields<Foam::sphericalTensor>
646     (
647         mesh, objects, reader_->GetPointFieldSelection()
648     );
649     convertPointFields<Foam::symmTensor>
650     (
651         mesh, objects, reader_->GetPointFieldSelection()
652     );
653     convertPointFields<Foam::tensor>
654     (
655         mesh, objects, reader_->GetPointFieldSelection()
656     );
658     if (debug)
659     {
660         Info<< "done" << endl;
661     }
665 // ************************************************************************* //