xref: /phasta/phastaIO/vtkPPhastaReader.cxx (revision 595995161822a203c8467e0e4a253d7bd7d6df32)
1 /*=========================================================================
2 
3   Program:   Visualization Toolkit
4   Module:    $RCSfile: vtkPPhastaReader.cxx,v $
5 
6   Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7   All rights reserved.
8   See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9 
10      This software is distributed WITHOUT ANY WARRANTY; without even
11      the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12      PURPOSE.  See the above copyright notice for more information.
13 
14 =========================================================================*/
15 #include "vtkPPhastaReader.h"
16 
17 #include "vtkCellData.h"
18 #include "vtkFieldData.h"
19 #include "vtkInformation.h"
20 #include "vtkInformationVector.h"
21 #include "vtkMultiBlockDataSet.h"
22 #include "vtkMultiPieceDataSet.h"
23 #include "vtkObjectFactory.h"
24 #include "vtkPointData.h"
25 #include "vtkPVXMLElement.h"
26 #include "vtkPVXMLParser.h"
27 #include "vtkPhastaReader.h"
28 #include "vtkSmartPointer.h"
29 #include "vtkStreamingDemandDrivenPipeline.h"
30 #include "vtkUnstructuredGrid.h"
31 
32 #include <vtksys/SystemTools.hxx>
33 
34 #include <vtkstd/map>
35 #include <vtksys/ios/sstream>
36 
37 int NUM_PIECES;
38 int NUM_FILES;
39 int TIME_STEP;
40 char * FILE_PATH;
41 int PART_ID;
42 int FILE_ID;
43 
44 double opentime_total = 0.0;
45 
46 /*
47  * Modified part is dealing with new phasta data format
48  * SyncIO and rbIO library, contact liun2@cs.rpi.edu
49  *
50  *               ------ Ning Liu
51  *               ------ Sept. 2010
52  */
53 
54 struct vtkPPhastaReaderInternal
55 {
56   struct TimeStepInfo
57   {
58     int GeomIndex;
59     int FieldIndex;
60     double TimeValue;
61 
TimeStepInfovtkPPhastaReaderInternal::TimeStepInfo62     TimeStepInfo() : GeomIndex(-1), FieldIndex(-1), TimeValue(0.0)
63       {
64       }
65   };
66 
67   typedef vtkstd::map<int, TimeStepInfo> TimeStepInfoMapType;
68   TimeStepInfoMapType TimeStepInfoMap;
69   typedef vtkstd::map<int, vtkSmartPointer<vtkUnstructuredGrid> >
70   CachedGridsMapType;
71   CachedGridsMapType CachedGrids;
72 };
73 
74 //----------------------------------------------------------------------------
75 vtkCxxRevisionMacro(vtkPPhastaReader, "$Revision: 1.6 $");
76 vtkStandardNewMacro(vtkPPhastaReader);
77 
78 //----------------------------------------------------------------------------
vtkPPhastaReader()79 vtkPPhastaReader::vtkPPhastaReader()
80 {
81   //this->DebugOn(); // comment out this line when in production
82   this->FileName = 0;
83 
84   this->TimeStepIndex = 0;
85   this->ActualTimeStep = 0;
86 
87   this->Reader = vtkPhastaReader::New();
88 
89   this->SetNumberOfInputPorts(0);
90 
91   this->Parser = 0;
92 
93   this->Internal = new vtkPPhastaReaderInternal;
94 
95   this->TimeStepRange[0] = 0;
96   this->TimeStepRange[1] = 0;
97 }
98 
99 //----------------------------------------------------------------------------
~vtkPPhastaReader()100 vtkPPhastaReader::~vtkPPhastaReader()
101 {
102   this->Reader->Delete();
103   this->SetFileName(0);
104 
105   if (this->Parser)
106     {
107     this->Parser->Delete();
108     }
109 
110   delete this->Internal;
111 }
112 
113 //----------------------------------------------------------------------------
RequestData(vtkInformation *,vtkInformationVector **,vtkInformationVector * outputVector)114 int vtkPPhastaReader::RequestData(vtkInformation*,
115                                   vtkInformationVector**,
116                                   vtkInformationVector* outputVector)
117 {
118   vtkDebugMacro("Entering PP RequestData()\n");
119   // get the data object
120   vtkInformation *outInfo =
121     outputVector->GetInformationObject(0);
122 
123   this->ActualTimeStep = this->TimeStepIndex;
124 
125   int tsLength =
126     outInfo->Length(vtkStreamingDemandDrivenPipeline::TIME_STEPS());
127   double* steps =
128     outInfo->Get(vtkStreamingDemandDrivenPipeline::TIME_STEPS());
129 
130   // Check if a particular time was requested.
131   if(outInfo->Has(vtkStreamingDemandDrivenPipeline::UPDATE_TIME_STEPS()))
132     {
133     // Get the requested time step. We only supprt requests of a single time
134     // step in this reader right now
135     double *requestedTimeSteps =
136       outInfo->Get(vtkStreamingDemandDrivenPipeline::UPDATE_TIME_STEPS());
137     double timeValue = requestedTimeSteps[0];
138 
139     // find the first time value larger than requested time value
140     // this logic could be improved
141     int cnt = 0;
142     while (cnt < tsLength-1 && steps[cnt] < timeValue)
143       {
144       cnt++;
145       }
146     this->ActualTimeStep = cnt;
147     }
148 
149   if (this->ActualTimeStep > this->TimeStepRange[1])
150     {
151     vtkErrorMacro("Timestep index too large.");
152     return 0;
153     }
154 
155   // get the current piece being requested
156   int piece =
157     outInfo->Get(vtkStreamingDemandDrivenPipeline::UPDATE_PIECE_NUMBER());
158 
159   int numProcPieces =
160     outInfo->Get(vtkStreamingDemandDrivenPipeline::UPDATE_NUMBER_OF_PIECES());
161 
162   // this is actually # of proc
163   vtkDebugMacro(<<"numProcPieces (i.e. # of proc) = "<< numProcPieces);
164 
165   if (!this->Parser)
166     {
167     vtkErrorMacro("No parser was created. Cannot read file");
168     return 0;
169     }
170 
171   vtkPVXMLElement* rootElement = this->Parser->GetRootElement();
172 
173   int numPieces;
174   int numFiles, timeStep;
175   if (!rootElement->GetScalarAttribute("number_of_pieces", &numPieces))
176     {
177     numPieces = 1;
178     }
179 
180   if (!rootElement->GetScalarAttribute("number_of_files", &numFiles))
181     {
182     numFiles = 1;
183     }
184 
185   NUM_PIECES=numPieces;
186   NUM_FILES=numFiles;
187 
188   vtkDebugMacro(<<"NEW PHT Parameter: number_of_files = "<< numFiles );
189 
190   vtkMultiBlockDataSet *output = vtkMultiBlockDataSet::SafeDownCast(
191     outInfo->Get(vtkDataObject::DATA_OBJECT()));
192   output->SetNumberOfBlocks(1);
193   vtkMultiPieceDataSet* MultiPieceDataSet = vtkMultiPieceDataSet::New();
194   // the following line was deleted in Ning's version
195   MultiPieceDataSet->SetNumberOfPieces(numPieces);
196   output->SetBlock(0, MultiPieceDataSet);
197   MultiPieceDataSet->Delete();
198 
199   const char* geometryPattern = 0;
200   int geomHasPiece = 0;
201   int geomHasTime = 0;
202   const char* fieldPattern = 0;
203   int fieldHasPiece = 0;
204   int fieldHasTime = 0;
205 
206   unsigned int numElements = rootElement->GetNumberOfNestedElements();
207   for (unsigned int i=0; i<numElements; i++)
208     {
209     vtkPVXMLElement* nested = rootElement->GetNestedElement(i);
210 
211     if (strcmp("GeometryFileNamePattern", nested->GetName()) == 0)
212       {
213       geometryPattern = nested->GetAttribute("pattern");
214       if (!nested->GetScalarAttribute("has_piece_entry", &geomHasPiece))
215         {
216         geomHasPiece = 0;
217         }
218       if (!nested->GetScalarAttribute("has_time_entry", &geomHasTime))
219         {
220         geomHasTime = 0;
221         }
222       }
223 
224     if (strcmp("FieldFileNamePattern", nested->GetName()) == 0)
225       {
226       fieldPattern = nested->GetAttribute("pattern");
227       if (!nested->GetScalarAttribute("has_piece_entry", &fieldHasPiece))
228         {
229         fieldHasPiece = 0;
230         }
231       if (!nested->GetScalarAttribute("has_time_entry", &fieldHasTime))
232         {
233         fieldHasTime = 0;
234         }
235       }
236     }
237 
238   if (!geometryPattern)
239     {
240     vtkErrorMacro("No geometry pattern was specified. Cannot load file");
241     return 0;
242     }
243 
244   if (!fieldPattern)
245     {
246     vtkErrorMacro("No field pattern was specified. Cannot load file");
247     return 0;
248     }
249 
250   char* geom_name = new char [ strlen(geometryPattern) + 60 ];
251   char* field_name = new char [ strlen(fieldPattern) + 60 ];
252 
253   //////////////////////
254   int numPiecesPerFile = numPieces/numFiles;
255   //////////////////////
256 
257   // now loop over all of the files that I should load
258   for(int loadingPiece=piece;loadingPiece<numPieces;loadingPiece+=numProcPieces)
259     {
260       TIME_STEP=this->Internal->TimeStepInfoMap[this->ActualTimeStep].FieldIndex;
261       FILE_ID = int(loadingPiece/numPiecesPerFile)+1; // this will be passed to PhastaReader by extern...
262       PART_ID = loadingPiece+1;
263 
264       vtkDebugMacro(<<"PP In loop, piece="<< piece <<", loadingPiece+1="<< loadingPiece +1 << ", numPieces="<<numPieces<<", FILE_ID=" << FILE_ID<<", numProcPieces=" << numProcPieces);
265 
266     if (geomHasTime && geomHasPiece)
267       {
268       sprintf(geom_name,
269               geometryPattern,
270               this->Internal->TimeStepInfoMap[this->ActualTimeStep].GeomIndex,
271               FILE_ID);
272       }
273     else if (geomHasPiece)
274       {
275         sprintf(geom_name, geometryPattern, FILE_ID);
276       }
277     else if (geomHasTime)
278       {
279       sprintf(geom_name,
280               geometryPattern,
281               this->Internal->TimeStepInfoMap[this->ActualTimeStep].GeomIndex);
282       }
283     else
284       {
285       strcpy(geom_name, geometryPattern);
286       }
287 
288     if (fieldHasTime && fieldHasPiece)
289       {
290       sprintf(field_name,
291               fieldPattern,
292               this->Internal->TimeStepInfoMap[this->ActualTimeStep].FieldIndex,
293   //          FILE_ID); // don't use file id, otherwise dup geom and field file id will make PhastaReader not update -- jingfu
294               loadingPiece+1);
295       }
296     else if (fieldHasPiece)
297       {
298       sprintf(field_name, fieldPattern, loadingPiece+1);
299   //FILE_ID);
300       }
301     else if (fieldHasTime)
302       {
303       sprintf(field_name,
304               fieldPattern,
305               this->Internal->TimeStepInfoMap[this->ActualTimeStep].FieldIndex);
306       }
307     else
308       {
309       strcpy(geom_name, fieldPattern);
310       }
311 
312     vtksys_ios::ostringstream geomFName;
313     vtkstd::string gpath = vtksys::SystemTools::GetFilenamePath(geom_name);
314     if (gpath.empty() || !vtksys::SystemTools::FileIsFullPath(gpath.c_str()))
315       {
316       vtkstd::string path = vtksys::SystemTools::GetFilenamePath(this->FileName);
317       if (!path.empty())
318         {
319         geomFName << path.c_str() << "/";
320         }
321       }
322     geomFName << geom_name << ends;
323     this->Reader->SetGeometryFileName(geomFName.str().c_str());
324 
325     vtksys_ios::ostringstream fieldFName;
326     // try to strip out the path of file, if it's a full path file name
327     vtkstd::string fpath = vtksys::SystemTools::GetFilenamePath(field_name);
328 
329     ///////////////////////////////////////////
330     FILE_PATH = new char[fpath.size()+1];
331     strcpy(FILE_PATH,fpath.c_str());
332     ///////////////////////////////////////////
333 
334     if (fpath.empty() || !vtksys::SystemTools::FileIsFullPath(fpath.c_str()))
335       {
336       vtkstd::string path = vtksys::SystemTools::GetFilenamePath(this->FileName); // FileName is the .pht file
337       if (!path.empty())
338         {
339           /////////////////////////////////////////
340           delete [] FILE_PATH;
341           FILE_PATH = new char[path.size()+1];
342           strcpy(FILE_PATH,path.c_str());
343           /////////////////////////////////////////
344           fieldFName << path.c_str() << "/";
345           //std::cout << "something might be wrong here, string path=" << path << ", fieldFName=" << fieldFName<< std::endl;
346         }
347       }
348     fieldFName << field_name << ends;
349     this->Reader->SetFieldFileName(fieldFName.str().c_str());
350 
351     vtkPPhastaReaderInternal::CachedGridsMapType::iterator CachedCopy =
352       this->Internal->CachedGrids.find(loadingPiece);
353 
354     // the following "if" was commented out in previous new version (tweaked by Ning)
355     // if there is a cached copy, use that
356     /*
357     if(CachedCopy != this->Internal->CachedGrids.end())
358       {
359       this->Reader->SetCachedGrid(CachedCopy->second);
360       printf("should use cached copy but can't compile\n");
361       }
362       */
363 
364     // In order to register etc, Reader need a new executative in every
365     // update call, otherwise it doesn't do anything
366     this->Reader->Update();
367     // the following "if" was commented out in previous new version (tweaked by Ning)
368 
369     /*
370     if(CachedCopy == this->Internal->CachedGrids.end())
371       {
372       vtkSmartPointer<vtkUnstructuredGrid> cached =
373         vtkSmartPointer<vtkUnstructuredGrid>::New();
374       cached->ShallowCopy(this->Reader->GetOutput());
375       cached->GetPointData()->Initialize();
376       cached->GetCellData()->Initialize();
377       cached->GetFieldData()->Initialize();
378       this->Internal->CachedGrids[loadingPiece] = cached;
379       }
380       */
381 
382     vtkSmartPointer<vtkUnstructuredGrid> copy =
383       vtkSmartPointer<vtkUnstructuredGrid>::New();
384     copy->ShallowCopy(this->Reader->GetOutput());
385     MultiPieceDataSet->SetPiece(loadingPiece, copy);
386     //MultiPieceDataSet->SetPiece(MultiPieceDataSet->GetNumberOfPieces(),copy); // Ning's version
387     }
388 
389   delete [] FILE_PATH;
390   delete [] geom_name;
391   delete [] field_name;
392 
393   if (steps)
394     {
395     output->GetInformation()->Set(vtkDataObject::DATA_TIME_STEPS(),
396                                   steps+this->ActualTimeStep, 1);
397     }
398 
399   vtkDebugMacro("End of PP RequestData()\n, total open time is " << opentime_total);
400   // if it's not too many printf, print it out
401   if( numProcPieces < 16) printf("total open time for sync-io is %lf (nf=%d, np=%d)\n", opentime_total, numFiles, numProcPieces);
402 
403   return 1;
404 }
405 
406 //----------------------------------------------------------------------------
RequestInformation(vtkInformation *,vtkInformationVector **,vtkInformationVector * outputVector)407 int vtkPPhastaReader::RequestInformation(vtkInformation*,
408                                        vtkInformationVector**,
409                                        vtkInformationVector* outputVector)
410 {
411   vtkDebugMacro(<<"In PP requestInformation() -- nothing modified in this func\n");
412   this->Internal->TimeStepInfoMap.clear();
413   this->Reader->ClearFieldInfo();
414 
415   vtkInformation *outInfo =
416     outputVector->GetInformationObject(0);
417 
418   if (!this->FileName)
419     {
420     vtkErrorMacro("FileName has to be specified.");
421     return 0;
422     }
423 
424   if (this->Parser)
425     {
426     this->Parser->Delete();
427     this->Parser = 0;
428     }
429 
430   vtkSmartPointer<vtkPVXMLParser> parser =
431     vtkSmartPointer<vtkPVXMLParser>::New();
432 
433   parser->SetFileName(this->FileName);
434   if (!parser->Parse())
435     {
436     return 0;
437     }
438 
439   vtkPVXMLElement* rootElement = parser->GetRootElement();
440   if (!rootElement)
441     {
442     vtkErrorMacro("Cannot parse file.");
443     return 0;
444     }
445 
446   if (strcmp(rootElement->GetName(), "PhastaMetaFile") != 0)
447     {
448     vtkErrorMacro("This is not a phasta metafile.");
449     return 0;
450     }
451 
452   this->Parser = parser;
453   parser->Register(this);
454 
455   int numTimeSteps=1;
456   int hasTimeValues = 0;
457 
458   unsigned int numElements = rootElement->GetNumberOfNestedElements();
459   for (unsigned int i=0; i<numElements; i++)
460     {
461     vtkPVXMLElement* nested = rootElement->GetNestedElement(i);
462     if (strcmp("TimeSteps", nested->GetName()) == 0)
463       {
464       if (!nested->GetScalarAttribute("number_of_steps", &numTimeSteps))
465         {
466         numTimeSteps = 1;
467         }
468       int autoGen;
469       int indexIncr;
470       int startIndex;
471       if (!nested->GetScalarAttribute("auto_generate_indices", &autoGen))
472         {
473         autoGen = 0;
474         }
475       if (!nested->GetScalarAttribute("increment_index_by", &indexIncr))
476         {
477         indexIncr = 1;
478         }
479       if (!nested->GetScalarAttribute("start_index", &startIndex))
480         {
481         startIndex = 0;
482         }
483       double startValue = 0.;
484       double valueIncr = 1.*indexIncr;
485       if (nested->GetScalarAttribute("start_value", &startValue))
486         {
487         hasTimeValues = 1;
488         }
489       if (nested->GetScalarAttribute("increment_value_by", &valueIncr))
490         {
491         hasTimeValues = 1;
492         }
493       if (autoGen)
494         {
495         for (int j=0; j<numTimeSteps; j++)
496           {
497           vtkPPhastaReaderInternal::TimeStepInfo& info =
498             this->Internal->TimeStepInfoMap[j];
499           info.GeomIndex = startIndex;
500           info.FieldIndex = startIndex;
501           info.TimeValue = startValue;
502           startIndex += indexIncr;
503           startValue += valueIncr;
504           }
505         }
506 
507       unsigned int numElements2 = nested->GetNumberOfNestedElements();
508       for (unsigned int j=0; j<numElements2; j++)
509         {
510         vtkPVXMLElement* nested2 = nested->GetNestedElement(j);
511         if (strcmp("TimeStep", nested2->GetName()) == 0)
512           {
513           int index;
514           if (nested2->GetScalarAttribute("index", &index))
515             {
516             if ( (index+1) > numTimeSteps )
517               {
518               numTimeSteps = index+1;
519               }
520             vtkPPhastaReaderInternal::TimeStepInfo& info =
521               this->Internal->TimeStepInfoMap[index];
522             int gIdx;
523             if (nested2->GetScalarAttribute("geometry_index",
524                                              &gIdx))
525               {
526               info.GeomIndex = gIdx;
527               }
528             int fIdx;
529             if (nested2->GetScalarAttribute("field_index",
530                                              &fIdx))
531               {
532               info.FieldIndex = fIdx;
533               }
534             double val;
535             if (nested2->GetScalarAttribute("value",
536                                             &val))
537               {
538               info.TimeValue = val;
539               hasTimeValues = 1;
540               }
541             }
542           }
543         }
544       break;
545       }
546     }
547 
548   int numberOfFields=0, numberOfFields2=0;
549   for (unsigned int i=0; i<numElements; i++)
550     {
551     vtkPVXMLElement* nested = rootElement->GetNestedElement(i);
552     if (strcmp("Fields", nested->GetName()) == 0)
553       {
554       if (!nested->GetScalarAttribute("number_of_fields", &numberOfFields))
555         {
556         numberOfFields = 1;
557         }
558 
559       numberOfFields2 = 0;
560       unsigned int numElements2 = nested->GetNumberOfNestedElements();
561       for (unsigned int j=0; j<numElements2; j++)
562         {
563         vtkPVXMLElement* nested2 = nested->GetNestedElement(j);
564         if (strcmp("Field", nested2->GetName()) == 0)
565           {
566           numberOfFields2++;
567           vtkstd::string paraviewFieldTagStr, dataTypeStr;
568           const char* paraviewFieldTag = 0;
569           paraviewFieldTag = nested2->GetAttribute("paraview_field_tag");
570           if (!paraviewFieldTag)
571             {
572             vtksys_ios::ostringstream paraviewFieldTagStrStream;
573             paraviewFieldTagStrStream << "Field " << numberOfFields2 << ends;
574             paraviewFieldTagStr = paraviewFieldTagStrStream.str();
575             paraviewFieldTag = paraviewFieldTagStr.c_str();
576             }
577           const char* phastaFieldTag = 0;
578           phastaFieldTag = nested2->GetAttribute("phasta_field_tag");
579           if (!phastaFieldTag)
580             {
581             vtkErrorMacro("No phasta field tag was specified");
582             return 0;
583             }
584           int index; // 0 as default (for safety)
585           if (!nested2->GetScalarAttribute("start_index_in_phasta_array", &index))
586             {
587             index = 0;
588             }
589           int numOfComps; // 1 as default (for safety)
590           if (!nested2->GetScalarAttribute("number_of_components", &numOfComps))
591             {
592             numOfComps = 1;
593             }
594           int dataDependency; // nodal as default
595           if (!nested2->GetScalarAttribute("data_dependency", &dataDependency))
596             {
597             dataDependency = 0;
598             }
599           const char* dataType = 0;
600           dataType = nested2->GetAttribute("data_type");
601           if (!dataType) // "double" as default
602             {
603             dataTypeStr = "double";
604             dataType = dataTypeStr.c_str();
605             }
606 
607           this->Reader->SetFieldInfo(paraviewFieldTag,phastaFieldTag,index,numOfComps,dataDependency,dataType);
608           }
609         }
610 
611       if (numberOfFields<numberOfFields2)
612         {
613         numberOfFields = numberOfFields2;
614         }
615 
616       break;
617       }
618     }
619 
620   if (!numberOfFields2) // by default take "solution" with only flow variables
621     {
622     numberOfFields = 3;
623     this->Reader->SetFieldInfo("pressure","solution",0,1,0,"double");
624     this->Reader->SetFieldInfo("velocity","solution",1,3,0,"double");
625     this->Reader->SetFieldInfo("temperature","solution",4,1,0,"double");
626     }
627 
628   int tidx;
629   // Make sure all indices are there
630   for (tidx=1; tidx<numTimeSteps; tidx++)
631     {
632     vtkPPhastaReaderInternal::TimeStepInfoMapType::iterator iter =
633       this->Internal->TimeStepInfoMap.find(tidx);
634     if (iter == this->Internal->TimeStepInfoMap.end())
635       {
636       vtkErrorMacro("Missing timestep, index=" << tidx);
637       return 0;
638       }
639     }
640 
641   if (hasTimeValues)
642     {
643     double* timeSteps = new double[numTimeSteps];
644     for (tidx=0; tidx<numTimeSteps; tidx++)
645       {
646       timeSteps[tidx] = this->Internal->TimeStepInfoMap[tidx].TimeValue;
647       }
648     outInfo->Set(vtkStreamingDemandDrivenPipeline::TIME_STEPS(),
649                  timeSteps,
650                  numTimeSteps);
651     double timeRange[2];
652     timeRange[0] = timeSteps[0];
653     timeRange[1] = timeSteps[numTimeSteps-1];
654     outInfo->Set(vtkStreamingDemandDrivenPipeline::TIME_RANGE(),
655                      timeRange, 2);
656     delete[] timeSteps;
657     }
658 
659   this->TimeStepRange[0] = 0;
660   this->TimeStepRange[1] = numTimeSteps-1;
661 
662   vtkInformation* info = outputVector->GetInformationObject(0);
663   info->Set(
664     vtkStreamingDemandDrivenPipeline::MAXIMUM_NUMBER_OF_PIECES(), -1);
665 
666   return 1;
667 }
668 
669 //-----------------------------------------------------------------------------
CanReadFile(const char * filename)670 int vtkPPhastaReader::CanReadFile(const char *filename)
671 {
672   vtkSmartPointer<vtkPVXMLParser> parser
673     = vtkSmartPointer<vtkPVXMLParser>::New();
674   parser->SuppressErrorMessagesOn();
675   parser->SetFileName(filename);
676 
677   // Make sure we can parse the XML metafile.
678   if (!parser->Parse()) return 0;
679 
680   // Make sure the XML file has a root element and it is of the right tag.
681   vtkPVXMLElement *rootElement = parser->GetRootElement();
682   if (!rootElement) return 0;
683   if (strcmp(rootElement->GetName(), "PhastaMetaFile") != 0) return 0;
684 
685   // The file clearly is supposed to be a Phasta file.
686   return 1;
687 }
688 
689 //----------------------------------------------------------------------------
PrintSelf(ostream & os,vtkIndent indent)690 void vtkPPhastaReader::PrintSelf(ostream& os, vtkIndent indent)
691 {
692   this->Superclass::PrintSelf(os,indent);
693 
694   os << indent << "FileName: "
695      << (this->FileName?this->FileName:"(none)")
696      << endl;
697   os << indent << "TimeStepIndex: " << this->TimeStepIndex << endl;
698   os << indent << "TimeStepRange: "
699      << this->TimeStepRange[0] << " " << this->TimeStepRange[1]
700      << endl;
701 }
702 
703