1 /*=========================================================================
2 
3   Program:   Visualization Toolkit
4   Module:    vtkTecplotReader.cxx
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 
16 /*****************************************************************************
17 *
18 * Copyright (c) 2000 - 2009, Lawrence Livermore National Security, LLC
19 * Produced at the Lawrence Livermore National Laboratory
20 * LLNL-CODE-400124
21 * All rights reserved.
22 *
23 * This file was adapted from the ASCII Tecplot reader of VisIt. For  details,
24 * see https://visit.llnl.gov/.  The full copyright notice is contained in the
25 * file COPYRIGHT located at the root of the VisIt distribution or at
26 * http://www.llnl.gov/visit/copyright.html.
27 *
28 *****************************************************************************/
29 
30 #include "vtkTecplotReader.h"
31 
32 #include "vtkInformation.h"
33 #include "vtkObjectFactory.h"
34 #include "vtkCallbackCommand.h"
35 #include "vtkInformationVector.h"
36 #include "vtkMultiBlockDataSet.h"
37 #include "vtkCompositeDataPipeline.h"
38 
39 #include "vtkPoints.h"
40 #include "vtkCellData.h"
41 #include "vtkPointData.h"
42 #include "vtkCellArray.h"
43 #include "vtkFloatArray.h"
44 #include "vtkIdTypeArray.h"
45 #include "vtkStructuredGrid.h"
46 #include "vtkUnstructuredGrid.h"
47 #include "vtkUnsignedCharArray.h"
48 #include "vtkDataArraySelection.h"
49 
50 #include "vtk_zlib.h"
51 #include <vtksys/SystemTools.hxx>
52 
53 #include <fstream>
54 
55 #include <ctype.h> // for isspace(), isalnum()
56 
57 vtkStandardNewMacro( vtkTecplotReader );
58 
59 
60 // ============================================================================
61 class FileStreamReader
62 {
63 public:
64   FileStreamReader();
65   ~FileStreamReader();
66 
67   bool open(const char* fileName);
is_open() const68   bool is_open()const {return Open;}
eof() const69   bool eof()const {return Eof;}
70 
71   void rewind();
72   void close();
73   int get();
74   bool operator !() const;
75 
76 protected:
77   bool Open;
78   bool Eof;
79   static const unsigned int BUFF_SIZE = 2048;
80   char buff[BUFF_SIZE];
81   unsigned int Pos;
82   unsigned int BuffEnd;
83   gzFile file;
84   std::string FileName;
85 
86 };
87 
88 // ----------------------------------------------------------------------------
FileStreamReader()89 FileStreamReader::FileStreamReader()
90   : Open(false),Eof(true),Pos(BUFF_SIZE),BuffEnd(BUFF_SIZE),FileName()
91 {
92 
93 }
94 
95 // ----------------------------------------------------------------------------
~FileStreamReader()96 FileStreamReader::~FileStreamReader()
97 {
98   this->close();
99 }
100 
101 // ----------------------------------------------------------------------------
open(const char * fileName)102 bool FileStreamReader::open( const char* fileName )
103   {
104   if ( !this->Open )
105     {
106     this->FileName = std::string(fileName);
107     //zlib handles both compressed and uncompressed file
108     //we just have peak into the file and see if it has the magic
109     //flags or not
110     unsigned char magic[2];
111     FILE *ff = fopen(fileName,"rb");
112     size_t count = fread(magic,1,2,ff);
113     fclose(ff);
114 
115     // only continue if fread succeeded
116     if (count == 2)
117       {
118       const char* mode = (magic[0] == 0x1f && magic[1] == 0x8b) ? "rb" : "r";
119       this->file = gzopen(fileName,mode);
120 
121       this->Eof = (this->file == 0);
122       this->Open = (this->file != 0);
123       this->Pos = BUFF_SIZE;
124       }
125     }
126   return this->Open;
127   }
128 // ----------------------------------------------------------------------------
get()129 int FileStreamReader::get()
130 {
131   if (!this->is_open() || this->eof() )
132     {
133     return this->eof();
134     }
135 
136   //when reading uncompressed data, zlib will return if it hits
137   //and eol character
138 
139   if (this->Pos >= this->BuffEnd)
140     {
141     this->Pos = 0;
142     //read the first buffer
143     this->BuffEnd = gzread(this->file,this->buff,this->BUFF_SIZE);
144     //assign EOF to what gzread returned
145     this->Eof = (this->BuffEnd <= 0);
146     if (this->Eof)
147       {
148       return this->Eof;
149       }
150     }
151   return this->buff[this->Pos++];
152 }
153 
154 // ----------------------------------------------------------------------------
rewind()155 void FileStreamReader::rewind( )
156 {
157   if ( this->Open )
158     {
159     //we don't want to use gzrewind as it rewinds to not the start of the
160     //file, but to start of the data in the file, meaning we are past any
161     //comments or headers.
162     std::string fn = this->FileName;
163     this->close();
164     this->open(fn.c_str());
165     }
166 }
167 
168 // ----------------------------------------------------------------------------
close()169 void FileStreamReader::close()
170 {
171   if ( this->Open )
172     {
173     this->Open = false;
174     this->Eof = false;
175     this->Pos = this->BUFF_SIZE;
176     this->BuffEnd = this->BUFF_SIZE;
177     this->FileName = std::string();
178 
179     gzclose(this->file);
180     }
181 }
182 
183 // ----------------------------------------------------------------------------
operator !() const184 bool FileStreamReader::operator!() const
185 {
186   return this->Eof;
187 }
188 
189 // ==========================================================================//
190 class vtkTecplotReaderInternal
191 {
192 public:
vtkTecplotReaderInternal()193   vtkTecplotReaderInternal()  { this->Init(); }
~vtkTecplotReaderInternal()194   ~vtkTecplotReaderInternal() { this->Init(); }
195 
196   int     XIdInList;
197   int     YIdInList;
198   int     ZIdInList;
199   int     Completed;
200   int     GeometryDim;
201   int     TopologyDim;
202   char    TheNextChar;
203   bool    NextCharEOF;
204   bool    NextCharEOL;
205   bool    NextCharValid;
206   bool    TokenIsString;
207   bool    IsCompressed;
208   FileStreamReader ASCIIStream;
209   std::string TokenBackup;
210 
211 public:
Init()212   void Init()
213   {
214     this->Completed =  0;
215     this->XIdInList = -1;
216     this->YIdInList = -1;
217     this->ZIdInList = -1;
218 
219     this->TopologyDim   = 0;
220     this->GeometryDim   = 1;
221     this->TheNextChar   = '\0';
222     this->TokenBackup   = "";
223     this->NextCharEOF   = false;
224     this->NextCharEOL   = false;
225     this->NextCharValid = false;
226     this->TokenIsString = false;
227     this->IsCompressed  = false;
228   }
229 
230   // This functions obtains the next token from the ASCII stream.
231   // Note that it is assumed that the ASCII stream is ready and no
232   // reading error occurs.
GetNextToken()233   std::string GetNextToken()
234   {
235     // this is where we take a one-token lookahead
236     if ( !this->TokenBackup.empty() )
237       {
238       std::string  retval = this->TokenBackup;
239       this->TokenBackup = "";
240       return  retval;
241       }
242 
243     // oops!  we hit EOF and someone still wants more.
244     if ( this->NextCharEOF )
245       {
246       return "";
247       }
248 
249     this->NextCharEOL   = false;
250     this->TokenIsString = false;
251 
252     std::string  retval = "";
253     if ( !this->NextCharValid )
254       {
255       this->TheNextChar   = this->ASCIIStream.get();
256       this->NextCharValid = true;
257 
258       if ( !this->ASCIIStream )
259         {
260         this->NextCharEOF = true;
261         }
262       }
263 
264     //if the token is a comment token, skip the entire line
265     if (!this->NextCharEOF && this->TheNextChar == '#')
266       {
267       while ( (!this->NextCharEOF) &&
268                 (this->TheNextChar != '\n') &&
269                 (this->TheNextChar != '\r') )
270         {
271         this->TheNextChar = this->ASCIIStream.get();
272         if ( this->TheNextChar == '\n'||
273              this->TheNextChar == '\r' )
274           {
275           this->NextCharEOL = true;
276           }
277         }
278       }
279 
280     // skip inter-token whitespace
281     while ( ! this->NextCharEOF &&
282             ( this->TheNextChar == ' '  ||
283               this->TheNextChar == '\n' ||
284               this->TheNextChar == '\r' ||
285               this->TheNextChar == '\t' ||
286               this->TheNextChar == '='  ||
287               this->TheNextChar == '('  ||
288               this->TheNextChar == ')'  ||
289               this->TheNextChar == ','
290             )
291           )
292       {
293       if ( this->TheNextChar == '\n' || this->TheNextChar == '\r' )
294         {
295         this->NextCharEOL = true;
296         }
297 
298       this->TheNextChar = this->ASCIIStream.get();
299       if ( !this->ASCIIStream )
300         {
301         this->NextCharEOF = true;
302         }
303 
304       // Ignore blank lines since they don't return a token
305       if ( this->NextCharEOL )
306         {
307         return this->GetNextToken();
308         }
309       }
310 
311     if ( this->TheNextChar == '\"' )
312       {
313       this->TokenIsString = true;
314       this->TheNextChar = this->ASCIIStream.get();
315       if ( !this->ASCIIStream )
316         {
317         this->NextCharEOF = true;
318         }
319 
320       while ( !this->NextCharEOF && this->TheNextChar != '\"' )
321         {
322         retval += this->TheNextChar;
323         this->TheNextChar = this->ASCIIStream.get();
324 
325         if ( !this->ASCIIStream )
326           {
327           this->NextCharEOF = true;
328           }
329         }
330 
331       this->TheNextChar = this->ASCIIStream.get();
332       if ( !this->ASCIIStream )
333         {
334         this->NextCharEOF = true;
335         }
336       }
337     else
338       {
339       // handle a normal token
340       while ( ! this->NextCharEOF &&
341               ( this->TheNextChar != ' '  &&
342                 this->TheNextChar != '\n' &&
343                 this->TheNextChar != '\r' &&
344                 this->TheNextChar != '\t' &&
345                 this->TheNextChar != '='  &&
346                 this->TheNextChar != '('  &&
347                 this->TheNextChar != ')'  &&
348                 this->TheNextChar != ','
349               )
350             )
351         {
352         if ( this->TheNextChar >= 'a' && this->TheNextChar <= 'z' )
353           {
354           this->TheNextChar += (  int( 'A' ) - int( 'a' )  );
355           }
356 
357         retval += this->TheNextChar;
358         this->TheNextChar = this->ASCIIStream.get();
359 
360         if ( !this->ASCIIStream )
361           {
362           this->NextCharEOF = true;
363           }
364         }
365       }
366 
367     // skip whitespace to EOL
368     while ( ! this->NextCharEOF &&
369             ( this->TheNextChar == ' '  ||
370               this->TheNextChar == '\n' ||
371               this->TheNextChar == '\r' ||
372               this->TheNextChar == '\t' ||
373               this->TheNextChar == '='  ||
374               this->TheNextChar == '('  ||
375               this->TheNextChar == ')'  ||
376               this->TheNextChar == ','
377             )
378           )
379       {
380       if ( this->TheNextChar == '\n' || this->TheNextChar == '\r' )
381         {
382         this->NextCharEOL = true;
383         }
384 
385       this->TheNextChar = this->ASCIIStream.get();
386       if ( !this->ASCIIStream )
387         {
388         this->NextCharEOF = true;
389         }
390 
391       if ( this->NextCharEOL )
392         {
393         break;
394         }
395       }
396     return retval;
397   }
398 
399 private:
400 
401   vtkTecplotReaderInternal( const vtkTecplotReaderInternal & );  // Not implemented.
402   void operator = ( const vtkTecplotReaderInternal & );          // Not implemented.
403 };
404 // ==========================================================================//
405 
406 // ----------------------------------------------------------------------------
407 //                         Supporting Functions (begin)
408 // ----------------------------------------------------------------------------
409 
410 #ifndef MAX
411 #define MAX( a, b ) (  (a) > (b) ? (a) : (b)  )
412 #endif
413 
414 // ----------------------------------------------------------------------------
GetCoord(const std::string & theToken)415 static int GetCoord( const std::string & theToken )
416 {
417   if ( theToken == "X" || theToken == "x" || theToken == "I" )
418     {
419     return 0;
420     }
421 
422   if ( theToken == "Y" || theToken == "y" || theToken == "J" )
423     {
424     return 1;
425     }
426 
427   if ( theToken == "Z" || theToken == "z" || theToken == "K" )
428     {
429     return 2;
430     }
431 
432   return -1;
433 }
434 
435 // ----------------------------------------------------------------------------
GuessCoord(const std::string & theToken)436 static int GuessCoord( const std::string & theToken )
437 {
438   int  guessVal = GetCoord( theToken );
439 
440   if ( theToken.length() >= 3 )
441     {
442     // do match: "x[m]" or "x (m)", etc. don't match: "x velocity"
443     if (  (  !isspace( theToken[1] ) && !isalnum( theToken[1] )  ) ||
444           (   isspace( theToken[1]   && !isalnum( theToken[2] )  )   )
445        )
446       {
447       guessVal = GetCoord(  theToken.substr( 0, 1 )  );
448       }
449     }
450 
451   return guessVal;
452 }
453 
454 // ----------------------------------------------------------------------------
SimplifyWhitespace(const std::string & s)455 static std::string SimplifyWhitespace( const std::string & s )
456 {
457   int headIndx = 0;
458   int tailIndx = int( s.length() ) - 1;
459 
460   while (  headIndx < tailIndx && ( s[headIndx] == ' ' || s[headIndx] == '\t' )  )
461     {
462     headIndx ++;
463     }
464 
465   while (  tailIndx > headIndx && ( s[tailIndx]  == ' ' || s[tailIndx]  == '\t' )  )
466     {
467     tailIndx --;
468     }
469 
470   return s.substr( headIndx, tailIndx - headIndx + 1 );
471 }
472 
473 // ----------------------------------------------------------------------------
474 //                         Supporting Functions ( end )
475 // ----------------------------------------------------------------------------
476 
477 // ----------------------------------------------------------------------------
vtkTecplotReader()478 vtkTecplotReader::vtkTecplotReader()
479 {
480   this->SelectionObserver  = vtkCallbackCommand::New();
481   this->SelectionObserver->SetClientData( this );
482   this->SelectionObserver->SetCallback
483         ( &vtkTecplotReader::SelectionModifiedCallback );
484 
485   this->DataArraySelection = vtkDataArraySelection::New();
486   this->DataArraySelection->AddObserver( vtkCommand::ModifiedEvent,
487                                          this->SelectionObserver );
488 
489   this->FileName = NULL;
490   this->Internal = new vtkTecplotReaderInternal;
491   this->SetNumberOfInputPorts( 0 );
492 
493   this->Init();
494 }
495 
496 // ----------------------------------------------------------------------------
~vtkTecplotReader()497 vtkTecplotReader::~vtkTecplotReader()
498 {
499   this->Init();
500 
501   delete [] this->FileName;
502 
503   delete this->Internal;
504   this->Internal = NULL;
505 
506   this->DataArraySelection->RemoveAllArrays();
507   this->DataArraySelection->RemoveObserver( this->SelectionObserver );
508   this->DataArraySelection->Delete();
509   this->DataArraySelection = NULL;
510 
511   this->SelectionObserver->SetClientData( NULL );
512   this->SelectionObserver->SetCallback( NULL );
513   this->SelectionObserver->Delete();
514   this->SelectionObserver  = NULL;
515 }
516 
517 // ----------------------------------------------------------------------------
Init()518 void vtkTecplotReader::Init()
519 {
520   // do NOT address this->FileName in this function !!!
521 
522   this->DataTitle         = "";
523   this->NumberOfVariables = 0;
524   this->CellBased.clear();
525   this->ZoneNames.clear();
526   this->Variables.clear();
527 
528   this->Internal->Init();
529 }
530 
531 // ----------------------------------------------------------------------------
SetFileName(const char * fileName)532 void vtkTecplotReader::SetFileName( const char * fileName )
533 {
534   if (    fileName
535        && strcmp( fileName, "" )
536        && ( ( this->FileName == NULL ) || strcmp( fileName, this->FileName ) )
537      )
538     {
539     delete [] this->FileName;
540     this->FileName = new char[  strlen( fileName ) + 1  ];
541     strcpy( this->FileName, fileName );
542     this->FileName[ strlen( fileName ) ] = '\0';
543 
544     this->Modified();
545     this->Internal->Completed = 0;
546     }
547 }
548 
549 //----------------------------------------------------------------------------
SelectionModifiedCallback(vtkObject *,unsigned long,void * tpReader,void *)550 void vtkTecplotReader::SelectionModifiedCallback( vtkObject *, unsigned long,
551                                                   void * tpReader, void * )
552 {
553   static_cast< vtkTecplotReader * >( tpReader )->Modified();
554 }
555 
556 // ----------------------------------------------------------------------------
FillOutputPortInformation(int vtkNotUsed (port),vtkInformation * info)557 int vtkTecplotReader::FillOutputPortInformation
558   (  int vtkNotUsed( port ),  vtkInformation * info  )
559 {
560   info->Set( vtkDataObject::DATA_TYPE_NAME(), "vtkMultiBlockDataSet" );
561   return 1;
562 }
563 
564 // ----------------------------------------------------------------------------
RequestInformation(vtkInformation * request,vtkInformationVector ** inputVector,vtkInformationVector * outputVector)565 int vtkTecplotReader::RequestInformation( vtkInformation * request,
566                                           vtkInformationVector ** inputVector,
567                                           vtkInformationVector  * outputVector )
568 {
569   if(  !this->Superclass::RequestInformation
570         ( request, inputVector, outputVector )
571     )
572     {
573     return 0;
574     }
575 
576   this->GetDataArraysList();
577 
578   return 1;
579 }
580 
581 // ----------------------------------------------------------------------------
RequestData(vtkInformation * vtkNotUsed (request),vtkInformationVector ** vtkNotUsed (inputVector),vtkInformationVector * outputVector)582 int vtkTecplotReader::RequestData( vtkInformation *        vtkNotUsed( request ),
583                                    vtkInformationVector ** vtkNotUsed( inputVector ),
584                                    vtkInformationVector *  outputVector )
585 {
586   vtkInformation *       outInf = outputVector->GetInformationObject( 0 );
587   vtkMultiBlockDataSet * output = vtkMultiBlockDataSet::SafeDownCast
588                                   (  outInf->Get( vtkDataObject::DATA_OBJECT() )  );
589 
590   this->Internal->Completed = 0;
591   this->ReadFile( output );
592   outInf = NULL;
593   output = NULL;
594 
595   return 1;
596 }
597 
598 // ----------------------------------------------------------------------------
GetDataTitle()599 const char * vtkTecplotReader::GetDataTitle()
600 {
601   return this->DataTitle.c_str();
602 }
603 
604 // ----------------------------------------------------------------------------
GetNumberOfBlocks()605 int   vtkTecplotReader::GetNumberOfBlocks()
606 {
607   return int( this->ZoneNames.size() );
608 }
609 
610 // ----------------------------------------------------------------------------
GetBlockName(int blockIdx)611 const char * vtkTecplotReader::GetBlockName( int blockIdx )
612 {
613   if (  blockIdx < 0  ||  blockIdx >= int( this->ZoneNames.size() )  )
614     {
615     return NULL;
616     }
617 
618   return this->ZoneNames[ blockIdx ].c_str();
619 }
620 
621 // ----------------------------------------------------------------------------
GetNumberOfDataAttributes()622 int  vtkTecplotReader::GetNumberOfDataAttributes()
623 {
624   return
625     this->NumberOfVariables - (   !(  !( this->Internal->XIdInList + 1 )  )   )
626                             - (   !(  !( this->Internal->YIdInList + 1 )  )   )
627                             - (   !(  !( this->Internal->ZIdInList + 1 )  )   );
628 }
629 
630 // ----------------------------------------------------------------------------
GetDataAttributeName(int attrIndx)631 const char * vtkTecplotReader::GetDataAttributeName( int attrIndx )
632 {
633   if ( attrIndx < 0 && attrIndx >= this->GetNumberOfDataAttributes() )
634     {
635     return NULL;
636     }
637 
638   return this->Variables[ attrIndx +
639                           this->Variables.size() -
640                           this->GetNumberOfDataAttributes()
641                         ].c_str();
642 }
643 
644 // ----------------------------------------------------------------------------
IsDataAttributeCellBased(int attrIndx)645 int   vtkTecplotReader::IsDataAttributeCellBased( int attrIndx )
646 {
647   int  cellBasd  =-1;
648   if ( attrIndx >= 0 && attrIndx < this->GetNumberOfDataAttributes() )
649     {
650     // the if-statement ensures that this->CellBased has been ready
651     cellBasd = this->CellBased[ attrIndx +
652                                 this->CellBased.size() -
653                                 this->GetNumberOfDataAttributes()
654                               ];
655     }
656 
657   return cellBasd;
658 }
659 
660 // ----------------------------------------------------------------------------
IsDataAttributeCellBased(const char * attrName)661 int   vtkTecplotReader::IsDataAttributeCellBased( const char * attrName )
662 {
663   int  cellBasd = -1;
664   int  varIndex = -1;
665 
666   if ( attrName )
667     {
668     for ( unsigned int i = 0; i < this->Variables.size(); i ++ )
669       {
670       if (  strcmp( this->Variables[i].c_str(), attrName ) == 0  )
671         {
672         varIndex = i;
673         break;
674         }
675       }
676 
677     cellBasd = ( varIndex == -1 ) ? -1: this->CellBased[ varIndex ];
678     }
679 
680   return cellBasd;
681 }
682 
683 //-----------------------------------------------------------------------------
GetNumberOfDataArrays()684 int vtkTecplotReader::GetNumberOfDataArrays()
685 {
686   return this->DataArraySelection->GetNumberOfArrays();
687 }
688 
689 //-----------------------------------------------------------------------------
GetDataArrayName(int arrayIdx)690 const char* vtkTecplotReader::GetDataArrayName( int arrayIdx )
691 {
692   return this->DataArraySelection->GetArrayName( arrayIdx );
693 }
694 
695 //-----------------------------------------------------------------------------
GetDataArrayStatus(const char * arayName)696 int vtkTecplotReader::GetDataArrayStatus( const char * arayName )
697 {
698   return this->DataArraySelection->ArrayIsEnabled( arayName );
699 }
700 
701 //-----------------------------------------------------------------------------
SetDataArrayStatus(const char * arayName,int bChecked)702 void vtkTecplotReader::SetDataArrayStatus( const char * arayName, int bChecked )
703 {
704   vtkDebugMacro( "Set cell array \"" << arayName
705                  << "\" status to: " << bChecked );
706 
707   if( bChecked )
708     {
709     this->DataArraySelection->EnableArray( arayName );
710     }
711   else
712     {
713     this->DataArraySelection->DisableArray( arayName );
714     }
715 }
716 
717 // ----------------------------------------------------------------------------
PrintSelf(ostream & os,vtkIndent indent)718 void vtkTecplotReader::PrintSelf( ostream & os, vtkIndent indent )
719 {
720   this->Superclass::PrintSelf( os, indent );
721 
722   os << indent << "DataTitle: "          << this->DataTitle         << endl;
723   os << indent << "Size of CellBased: "  << this->CellBased.size()  << endl;
724   os << indent << "Size of ZoneNames: "  << this->ZoneNames.size()  << endl;
725   os << indent << "Size of Variables: "  << this->Variables.size()  << endl;
726   os << indent << "NumberOfVariables: "  << this->NumberOfVariables << endl;
727 }
728 
729 // ----------------------------------------------------------------------------
GetArraysFromPointPackingZone(int numNodes,vtkPoints * theNodes,vtkPointData * nodeData)730 void vtkTecplotReader::GetArraysFromPointPackingZone
731   ( int numNodes, vtkPoints * theNodes, vtkPointData * nodeData )
732 {
733   // NOTE: The Tecplot ASCII file format mandates that cell data of any zone be
734   // stored in block-packing mode (VARLOCATION, pp. 158, Tecplot 360 Data Format
735   // Guide 2009). Thus we do not need to consider any cell data in this function.
736 
737   if ( !theNodes || !nodeData ||
738     !this->Internal->ASCIIStream.is_open()
739      )
740     {
741     vtkErrorMacro( << "File not open, errors with reading, or NULL vtkPoints /"
742                    << "vtkPointData." );
743     return;
744     }
745 
746   int     n,  v;
747   int     zArrayId;        // indexing zoneData
748   int     cordBase;        // offset of a 3D-coordinate triple in cordsPtr
749   int     isXcoord;
750   int     isYcoord;
751   int     isZcoord;
752   int   * anyCoord = NULL; // is any coordinate?
753   int   * coordIdx = NULL; // index of the coordinate array, just in case
754   int   * selected = NULL; // is a selected data array?
755   float * cordsPtr = NULL;
756   float * arrayPtr = NULL;
757   float   theValue;
758   vtkFloatArray * theArray = NULL;
759   std::vector< float * > pointers;
760   std::vector< vtkFloatArray * > zoneData;
761 
762   pointers.clear();
763   zoneData.clear();
764 
765   // geoemtry: 3D point coordinates (note that this array must be initialized
766   // since only 2D coordinates might be provided by a Tecplot file)
767   theNodes->SetNumberOfPoints( numNodes );
768   cordsPtr = static_cast< float * > (  theNodes->GetVoidPointer( 0 )  );
769   memset( cordsPtr, 0, sizeof( float ) * 3 * numNodes );
770 
771   // three arrays used to determine the role of each variable (including
772   // the coordinate arrays)
773   anyCoord = new int[ this->NumberOfVariables ];
774   coordIdx = new int[ this->NumberOfVariables ];
775   selected = new int[ this->NumberOfVariables ];
776 
777   // allocate arrays only if necessary to load the zone data
778   for ( v = 0; v < this->NumberOfVariables; v ++ )
779     {
780     isXcoord    = int(  !( v - this->Internal->XIdInList )  );
781     isYcoord    = int(  !( v - this->Internal->YIdInList )  );
782     isZcoord    = int(  !( v - this->Internal->ZIdInList )  );
783     anyCoord[v] = isXcoord + isYcoord + isZcoord;
784     coordIdx[v] = isYcoord + ( isZcoord << 1 );
785     selected[v] = this->DataArraySelection
786                       ->ArrayIsEnabled( this->Variables[v].c_str() );
787 
788     if ( anyCoord[v] + selected[v] )
789       {
790       theArray = vtkFloatArray::New();
791       theArray->SetNumberOfTuples( numNodes );
792       theArray->SetName( this->Variables[v].c_str() );
793       zoneData.push_back( theArray );
794       arrayPtr = static_cast< float * > (  theArray->GetVoidPointer( 0 )  );
795       pointers.push_back( arrayPtr );
796       arrayPtr = NULL;
797       theArray = NULL;
798       }
799     }
800 
801   // load the zone data (number of tuples <= number of points / nodes)
802   for ( n = 0; n < numNodes; n ++ )
803     {
804     cordBase = ( n << 1 ) + n;
805 
806     zArrayId = 0;
807     for ( v = 0; v < this->NumberOfVariables; v ++ )
808       {
809       // obtain a value that is either a coordinate or a selected attribute
810       if ( anyCoord[v] || selected[v] )
811         {
812         theValue = atof( this->Internal->GetNextToken().c_str() );
813         pointers[ zArrayId ++ ][n] = theValue;
814 
815         // collect the coordinate
816         if ( anyCoord[v] )
817           {
818           cordsPtr[ cordBase + coordIdx[v] ] = theValue;
819           }
820         }
821       else
822         {
823         // a value of an un-selected data array
824         this->Internal->GetNextToken();
825         }
826       }
827     }
828   cordsPtr = NULL;
829 
830   // attach the node-based data attributes to the grid
831   zArrayId = 0;
832   for ( v = 0; v < this->NumberOfVariables; v ++ )
833     {
834     if ( !anyCoord[v] && selected[v] )
835       {
836       nodeData->AddArray( zoneData[zArrayId] );
837       }
838 
839     zArrayId += int(    !(  !( anyCoord[v] + selected[v] )  )    );
840     }
841 
842   pointers.clear();
843 
844   //remove all the float arrays from vector so they won't leak
845   for ( unsigned int i=0; i < zoneData.size(); ++i)
846     {
847     vtkFloatArray *fa = zoneData.at(i);
848     if ( fa )
849       {
850       fa->FastDelete();
851       }
852     }
853   zoneData.clear();
854 
855   delete [] anyCoord;
856   delete [] coordIdx;
857   delete [] selected;
858   anyCoord = NULL;
859   coordIdx = NULL;
860   selected = NULL;
861 }
862 
863 // ----------------------------------------------------------------------------
GetArraysFromBlockPackingZone(int numNodes,int numCells,vtkPoints * theNodes,vtkPointData * nodeData,vtkCellData * cellData)864 void vtkTecplotReader::GetArraysFromBlockPackingZone( int numNodes, int numCells,
865      vtkPoints * theNodes, vtkPointData * nodeData, vtkCellData * cellData )
866 {
867   // NOTE: The Tecplot ASCII file format states that a block-packing zone may
868   // contain point data or cell data (VARLOCATION, pp. 158, Tecplot 360 Data
869   // Format Guide 2009). Thus we need to consider both cases in this function.
870 
871   if ( !theNodes || !nodeData || !cellData ||
872       !this->Internal->ASCIIStream.is_open()
873      )
874     {
875     vtkErrorMacro( << "File not open, errors with reading, or NULL vtkPoints /"
876                    << "vtkPointData / vtkCellData." );
877     return;
878     }
879 
880   int     v;
881   int     zArrayId;        // indexing zoneData
882   int     arraySiz;
883   int     isXcoord;
884   int     isYcoord;
885   int     isZcoord;
886   int   * anyCoord = NULL; // is any coordinate?
887   int   * selected = NULL; // is a selected data attribute?
888   float * cordsPtr = NULL;
889   float * arrayPtr = NULL;
890   vtkFloatArray * theArray = NULL;
891   std::vector< vtkFloatArray * > zoneData;
892   vtkDataSetAttributes * attribut[2] = { nodeData, cellData };
893 
894   zoneData.clear();
895 
896   // geoemtry: 3D point coordinates (note that this array must be initialized
897   // since only 2D coordinates might be provided by a Tecplot file)
898   theNodes->SetNumberOfPoints( numNodes );
899   cordsPtr = static_cast< float * > (  theNodes->GetVoidPointer( 0 )  );
900   memset( cordsPtr, 0, sizeof( float ) * 3 * numNodes );
901 
902   // two arrays used to determine the role of each variable (including
903   // the coordinate arrays)
904   anyCoord = new int[ this->NumberOfVariables ];
905   selected = new int[ this->NumberOfVariables ];
906 
907   for ( v = 0; v < this->NumberOfVariables; v ++ )
908     {
909     // check if this variable refers to a coordinate array
910     isXcoord = int(  !( v - this->Internal->XIdInList )  );
911     isYcoord = int(  !( v - this->Internal->YIdInList )  );
912     isZcoord = int(  !( v - this->Internal->ZIdInList )  );
913     anyCoord[v] = isXcoord + isYcoord + isZcoord;
914 
915     // in case of a data attribute, is it selected by the user?
916     selected[v] = this->DataArraySelection
917                       ->ArrayIsEnabled( this->Variables[v].c_str() );
918 
919     // obtain the size of the block
920     arraySiz = ( this->CellBased[v] ? numCells : numNodes );
921 
922     if ( anyCoord[v] || selected[v] )
923       {
924       // parse the block to extract either coordinates or data attribute
925       // values
926 
927       // extract the variable array throughout a block
928       theArray = vtkFloatArray::New();
929       theArray->SetNumberOfTuples( arraySiz );
930       theArray->SetName( this->Variables[v].c_str() );
931       zoneData.push_back( theArray );
932 
933       arrayPtr = static_cast< float * > (  theArray->GetVoidPointer( 0 )  );
934       for (int i = 0; i < arraySiz; i ++ )
935         {
936         arrayPtr[i] = atof( this->Internal->GetNextToken().c_str() );
937         }
938       theArray = NULL;
939 
940       // three special arrays are 'combined' to fill the 3D coord array
941       if ( anyCoord[v] )
942         {
943         float * coordPtr = cordsPtr + isYcoord + ( isZcoord << 1 );
944         for (int i = 0; i < arraySiz; i ++, coordPtr += 3 )
945           {
946           *coordPtr = arrayPtr[i];
947           }
948         coordPtr = NULL;
949         }
950 
951       arrayPtr = NULL;
952       }
953     else
954       {
955       // this block contains an un-selected data attribute and we
956       // need to read but ignore the values
957       for (int i = 0; i < arraySiz; i ++ )
958         {
959         this->Internal->GetNextToken();
960         }
961       }
962     }
963   cordsPtr = NULL;
964 
965   // attach the dataset attributes (node-based and cell-based) to the grid
966   // NOTE: zoneData[] and this->Variables (and this->CellBased) may differ
967   // in the number of the maintained arrays
968   zArrayId = 0;
969   for ( v = 0; v < this->NumberOfVariables; v ++ )
970     {
971     if ( !anyCoord[v] && selected[v] )
972       {
973       attribut[ this->CellBased[v] ]->AddArray( zoneData[zArrayId] );
974       }
975 
976     zArrayId += int(    !(  !( anyCoord[v] + selected[v] )  )    );
977     }
978 
979   //remove all the float arrays from vector so they won't leak
980   for ( unsigned int i=0; i < zoneData.size(); ++i)
981     {
982     vtkFloatArray *fa = zoneData.at(i);
983     if ( fa )
984       {
985       fa->FastDelete();
986       }
987     }
988   zoneData.clear();
989   attribut[0] = attribut[1] = NULL;
990 
991   delete [] anyCoord;
992   delete [] selected;
993   anyCoord = NULL;
994   selected = NULL;
995 }
996 
997 // ----------------------------------------------------------------------------
GetStructuredGridFromBlockPackingZone(int iDimSize,int jDimSize,int kDimSize,int zoneIndx,const char * zoneName,vtkMultiBlockDataSet * multZone)998 void vtkTecplotReader::GetStructuredGridFromBlockPackingZone
999   ( int iDimSize, int jDimSize, int kDimSize,
1000     int zoneIndx, const char * zoneName, vtkMultiBlockDataSet * multZone )
1001 {
1002   if ( !zoneName || !multZone )
1003     {
1004     vtkErrorMacro( "Zone name un-specified or NULL vtkMultiBlockDataSet." );
1005     return;
1006     }
1007 
1008   // determine the topological dimension
1009   if ( jDimSize == 1 && kDimSize == 1 )
1010     {
1011     this->Internal->TopologyDim = MAX( this->Internal->TopologyDim, 1 );
1012     }
1013   else
1014   if ( kDimSize == 1 )
1015     {
1016     this->Internal->TopologyDim = MAX( this->Internal->TopologyDim, 2 );
1017     }
1018   else
1019     {
1020     this->Internal->TopologyDim = MAX( this->Internal->TopologyDim, 3 );
1021     }
1022 
1023   // number of points, number of cells, and dimensionality
1024   int numNodes    = iDimSize * jDimSize * kDimSize;
1025   int numCells    = (  ( iDimSize <= 1 ) ? 1 : ( iDimSize - 1 )  ) *
1026                     (  ( jDimSize <= 1 ) ? 1 : ( jDimSize - 1 )  ) *
1027                     (  ( kDimSize <= 1 ) ? 1 : ( kDimSize - 1 )  );
1028   int gridDims[3] = { iDimSize, jDimSize, kDimSize };
1029 
1030   // Create vtkPoints and vtkStructuredGrid and associate them
1031   vtkPoints *         pntCords = vtkPoints::New();
1032   vtkStructuredGrid * strcGrid = vtkStructuredGrid::New();
1033   this->GetArraysFromBlockPackingZone( numNodes, numCells,
1034         pntCords, strcGrid->GetPointData(), strcGrid->GetCellData() );
1035   strcGrid->SetDimensions( gridDims );
1036   strcGrid->SetPoints( pntCords );
1037   pntCords->Delete();
1038   pntCords = NULL;
1039 
1040   if (    (    this->Internal->TopologyDim == 2
1041             || this->Internal->TopologyDim == 3
1042           )
1043        || (    this->Internal->TopologyDim == 0
1044             && this->Internal->GeometryDim >  1
1045           )
1046      )
1047     {
1048     multZone->SetBlock( zoneIndx, strcGrid );
1049     multZone->GetMetaData( zoneIndx )
1050             ->Set( vtkCompositeDataSet::NAME(), zoneName );
1051     }
1052   strcGrid->Delete();
1053   strcGrid = NULL;
1054 }
1055 
1056 // ----------------------------------------------------------------------------
GetStructuredGridFromPointPackingZone(int iDimSize,int jDimSize,int kDimSize,int zoneIndx,const char * zoneName,vtkMultiBlockDataSet * multZone)1057 void vtkTecplotReader::GetStructuredGridFromPointPackingZone
1058   ( int iDimSize, int jDimSize, int kDimSize,
1059     int zoneIndx, const char * zoneName, vtkMultiBlockDataSet * multZone )
1060 {
1061   if ( !zoneName || !multZone )
1062     {
1063     vtkErrorMacro( "Zone name un-specified or NULL vtkMultiBlockDataSet." );
1064     return;
1065     }
1066 
1067   if ( jDimSize == 1 && kDimSize == 1 )
1068     {
1069     this->Internal->TopologyDim = MAX( this->Internal->TopologyDim, 1 );
1070     }
1071   else if ( kDimSize == 1 )
1072     {
1073     this->Internal->TopologyDim = MAX( this->Internal->TopologyDim, 2 );
1074     }
1075   else
1076     {
1077     this->Internal->TopologyDim = MAX( this->Internal->TopologyDim, 3 );
1078     }
1079 
1080   // number of points, number of cells, and dimensionality
1081   int numNodes    = iDimSize * jDimSize * kDimSize;
1082   int gridDims[3] = { iDimSize, jDimSize, kDimSize };
1083 
1084   // Create vtkPoints and vtkStructuredGrid and associate them
1085   vtkPoints *         pntCords = vtkPoints::New();
1086   vtkStructuredGrid * strcGrid = vtkStructuredGrid::New();
1087   this->GetArraysFromPointPackingZone
1088         ( numNodes, pntCords, strcGrid->GetPointData() );
1089   strcGrid->SetDimensions( gridDims );
1090   strcGrid->SetPoints( pntCords );
1091   pntCords->Delete();
1092   pntCords = NULL;
1093 
1094   if (    (    this->Internal->TopologyDim == 2
1095             || this->Internal->TopologyDim == 3
1096           )
1097        || (    this->Internal->TopologyDim == 0
1098             && this->Internal->GeometryDim >  1
1099           )
1100      )
1101     {
1102     multZone->SetBlock( zoneIndx, strcGrid );
1103     multZone->GetMetaData( zoneIndx )
1104             ->Set( vtkCompositeDataSet::NAME(), zoneName );
1105     }
1106   strcGrid->Delete();
1107   strcGrid = NULL;
1108 }
1109 
1110 // ----------------------------------------------------------------------------
GetUnstructuredGridFromBlockPackingZone(int numNodes,int numCells,const char * cellType,int zoneIndx,const char * zoneName,vtkMultiBlockDataSet * multZone)1111 void vtkTecplotReader::GetUnstructuredGridFromBlockPackingZone
1112   ( int numNodes, int numCells, const char * cellType,
1113     int zoneIndx, const char * zoneName, vtkMultiBlockDataSet * multZone )
1114 {
1115   if ( !cellType || !zoneName || !multZone )
1116     {
1117     vtkErrorMacro( << "Zone name / cell type un-specified, or NULL "
1118                    << "vtkMultiBlockDataSet object." );
1119     return;
1120     }
1121 
1122   vtkPoints *           gridPnts = vtkPoints::New();
1123   vtkUnstructuredGrid * unstruct = vtkUnstructuredGrid::New();
1124   this->GetArraysFromBlockPackingZone( numNodes, numCells,
1125         gridPnts, unstruct->GetPointData(), unstruct->GetCellData() );
1126   this->GetUnstructuredGridCells( numCells, cellType, unstruct );
1127   unstruct->SetPoints( gridPnts );
1128   gridPnts->Delete();
1129   gridPnts = NULL;
1130 
1131   if (    (    this->Internal->TopologyDim == 2
1132             || this->Internal->TopologyDim == 3
1133           )
1134        || (    this->Internal->TopologyDim == 0
1135             && this->Internal->GeometryDim >  1
1136           )
1137      )
1138     {
1139     multZone->SetBlock( zoneIndx, unstruct );
1140     multZone->GetMetaData( zoneIndx )
1141             ->Set( vtkCompositeDataSet::NAME(), zoneName );
1142     }
1143   unstruct->Delete();
1144   unstruct = NULL;
1145 }
1146 
1147 // ----------------------------------------------------------------------------
GetUnstructuredGridFromPointPackingZone(int numNodes,int numCells,const char * cellType,int zoneIndx,const char * zoneName,vtkMultiBlockDataSet * multZone)1148 void vtkTecplotReader::GetUnstructuredGridFromPointPackingZone
1149   ( int numNodes, int numCells, const char * cellType,
1150     int zoneIndx, const char * zoneName, vtkMultiBlockDataSet * multZone )
1151 {
1152   if ( !cellType || !zoneName || !multZone )
1153     {
1154     vtkErrorMacro( << "Zone name / cell type un-specified, or NULL "
1155                    << "vtkMultiBlockDataSet object." );
1156     return;
1157     }
1158 
1159   vtkPoints *           gridPnts = vtkPoints::New();
1160   vtkUnstructuredGrid * unstruct = vtkUnstructuredGrid::New();
1161   this->GetArraysFromPointPackingZone
1162         ( numNodes, gridPnts, unstruct->GetPointData() );
1163   this->GetUnstructuredGridCells( numCells, cellType, unstruct );
1164   unstruct->SetPoints( gridPnts );
1165   gridPnts->Delete();
1166   gridPnts = NULL;
1167 
1168   if (    (    this->Internal->TopologyDim == 2
1169             || this->Internal->TopologyDim == 3
1170           )
1171        || (    this->Internal->TopologyDim == 0
1172             && this->Internal->GeometryDim >  1
1173           )
1174      )
1175     {
1176     multZone->SetBlock( zoneIndx, unstruct );
1177     multZone->GetMetaData( zoneIndx )
1178             ->Set( vtkCompositeDataSet::NAME(), zoneName );
1179     }
1180   unstruct->Delete();
1181   unstruct = NULL;
1182 }
1183 
1184 // ----------------------------------------------------------------------------
GetUnstructuredGridCells(int numberCells,const char * cellTypeStr,vtkUnstructuredGrid * unstrctGrid)1185 void vtkTecplotReader::GetUnstructuredGridCells( int numberCells,
1186   const char * cellTypeStr, vtkUnstructuredGrid * unstrctGrid )
1187 {
1188   if ( !cellTypeStr || !unstrctGrid )
1189     {
1190     vtkErrorMacro( << "Cell type (connectivity type) unspecified or NULL "
1191                    << "vtkUnstructuredGrid object." );
1192     return;
1193     }
1194 
1195   // determine the number of points per cell and the cell type
1196   int numCellPnts = -1;
1197   int theCellType = -1;
1198 
1199   if (  strcmp( cellTypeStr, "BRICK" ) == 0  )
1200     {
1201     numCellPnts = 8;
1202     theCellType = VTK_HEXAHEDRON;
1203     this->Internal->TopologyDim = MAX( this->Internal->TopologyDim, 3 );
1204     }
1205   else if (  strcmp( cellTypeStr, "TRIANGLE" ) == 0  )
1206     {
1207     numCellPnts = 3;
1208     theCellType = VTK_TRIANGLE;
1209     this->Internal->TopologyDim = MAX( this->Internal->TopologyDim, 2 );
1210     }
1211   else if (  strcmp( cellTypeStr, "QUADRILATERAL" ) == 0  )
1212     {
1213     numCellPnts = 4;
1214     theCellType = VTK_QUAD;
1215     this->Internal->TopologyDim = MAX( this->Internal->TopologyDim, 2 );
1216     }
1217   else if ( strcmp( cellTypeStr, "TETRAHEDRON" ) == 0  )
1218     {
1219     numCellPnts = 4;
1220     theCellType = VTK_TETRA;
1221     this->Internal->TopologyDim = MAX( this->Internal->TopologyDim, 3 );
1222     }
1223   else if (  strcmp( cellTypeStr, "POINT" ) == 0 || strcmp( cellTypeStr, "" ) == 0  )
1224     {
1225     numCellPnts = 1;
1226     theCellType = VTK_VERTEX;
1227     this->Internal->TopologyDim = MAX( this->Internal->TopologyDim, 0 );
1228     }
1229   else
1230     {
1231     vtkErrorMacro( << this->FileName << ": Unknown cell type for a zone." );
1232     }
1233 
1234   // the storage of each cell begins with the number of points per cell,
1235   // followed by a list of point ids representing the cell
1236   vtkIdTypeArray * cellInfoList = vtkIdTypeArray::New();
1237   cellInfoList->SetNumberOfValues(  ( numCellPnts + 1 ) * numberCells  );
1238   vtkIdType *      cellInforPtr = cellInfoList->GetPointer( 0 );
1239 
1240   // type of each cell
1241   vtkUnsignedCharArray * cellTypeList = vtkUnsignedCharArray::New();
1242   cellTypeList->SetNumberOfValues( numberCells );
1243   unsigned char *        cellTypesPtr = cellTypeList->GetPointer( 0 );
1244 
1245   // location of each cell in support of fast random (non-sequential) access
1246   vtkIdTypeArray * cellLocArray = vtkIdTypeArray::New();
1247   cellLocArray->SetNumberOfValues( numberCells );
1248   vtkIdType *      cellLocatPtr = cellLocArray->GetPointer( 0 );
1249 
1250   // fill the three arrays
1251   int locateOffset = 0;
1252   for ( int c = 0; c < numberCells; c ++ )
1253     {
1254     *cellTypesPtr ++ = theCellType;
1255     *cellInforPtr ++ = numCellPnts;
1256 
1257     // 1-origin connectivity array
1258     for ( int j = 0; j < numCellPnts; j ++ )
1259       {
1260       *cellInforPtr ++ = (   theCellType == VTK_VERTEX
1261                            ? c
1262                            : atoi( this->Internal->GetNextToken().c_str() ) - 1
1263                          );
1264       }
1265 
1266     *cellLocatPtr ++ = locateOffset;
1267     locateOffset    += numCellPnts + 1;
1268     }
1269   cellInforPtr = NULL;
1270   cellTypesPtr = NULL;
1271   cellLocatPtr = NULL;
1272 
1273   // create a cell array object to accept the cell info
1274   vtkCellArray * theCellArray = vtkCellArray::New();
1275   theCellArray->SetCells( numberCells, cellInfoList );
1276   cellInfoList->Delete();
1277   cellInfoList = NULL;
1278 
1279   // create a vtkUnstructuredGrid object and attach the 3 arrays (types, locations,
1280   // and cells) to it for export.
1281   unstrctGrid->SetCells( cellTypeList, cellLocArray, theCellArray );
1282   theCellArray->Delete();
1283   cellTypeList->Delete();
1284   cellLocArray->Delete();
1285   theCellArray = NULL;
1286   cellTypeList = NULL;
1287   cellLocArray = NULL;
1288 }
1289 
1290 // ----------------------------------------------------------------------------
GetDataArraysList()1291 void vtkTecplotReader::GetDataArraysList()
1292 {
1293   if (    ( this->Internal->Completed == 1 )
1294        || ( this->DataArraySelection->GetNumberOfArrays() > 0 )
1295        || ( this->FileName == NULL ) || (  strcmp( this->FileName, "" ) == 0  )
1296      )
1297     {
1298     return;
1299     }
1300 
1301   #define READ_UNTIL_TITLE_OR_VARIABLES !this->Internal->NextCharEOF &&\
1302                                         theTpToken != "TITLE"        &&\
1303                                         theTpToken != "VARIABLES"
1304   int             i;
1305   int             tpTokenLen = 0;
1306   int             guessedXid = -1;
1307   int             guessedYid = -1;
1308   int             guessedZid = -1;
1309   bool            tokenReady = false;
1310   std::string  theTpToken = "";
1311   std::string  noSpaceTok = "";
1312 
1313   this->Variables.clear();
1314   this->NumberOfVariables = 0;
1315 
1316   this->Internal->Init();
1317   this->Internal->ASCIIStream.open( this->FileName );
1318   theTpToken = this->Internal->GetNextToken();
1319 
1320   while ( !this->Internal->NextCharEOF )
1321     {
1322     tokenReady = false;
1323 
1324     if ( theTpToken == "" )
1325       {
1326       // whitespace: do nothing
1327       }
1328     else if ( theTpToken == "TITLE" )
1329       {
1330       this->Internal->GetNextToken();
1331       }
1332     else if ( theTpToken == "VARIABLES" )
1333       {
1334       theTpToken = this->Internal->GetNextToken();
1335 
1336       while ( this->Internal->TokenIsString )
1337         {
1338         tpTokenLen = int( theTpToken.length() );
1339         for ( i = 0; i < tpTokenLen; i ++ )
1340           {
1341           if ( theTpToken[i] == '(' )
1342             {
1343             theTpToken[i] = '[';
1344             }
1345           else if ( theTpToken[i] == ')' )
1346             {
1347             theTpToken[i] = ']';
1348             }
1349           else if ( theTpToken[i] == '/' )
1350             {
1351             theTpToken[i] = '_';
1352             }
1353           }
1354 
1355         noSpaceTok = SimplifyWhitespace( theTpToken );
1356 
1357         switch (  GetCoord( noSpaceTok )  )
1358           {
1359           case 0:  this->Internal->XIdInList = this->NumberOfVariables; break;
1360           case 1:  this->Internal->YIdInList = this->NumberOfVariables; break;
1361           case 2:  this->Internal->ZIdInList = this->NumberOfVariables; break;
1362           default: break;
1363           }
1364 
1365         switch (  GuessCoord( noSpaceTok )  )
1366           {
1367           case 0:  guessedXid = this->NumberOfVariables; break;
1368           case 1:  guessedYid = this->NumberOfVariables; break;
1369           case 2:  guessedZid = this->NumberOfVariables; break;
1370           default: break;
1371           }
1372 
1373         this->Variables.push_back( theTpToken );
1374         this->NumberOfVariables ++;
1375         theTpToken = this->Internal->GetNextToken();
1376         }
1377 
1378       if ( this->NumberOfVariables == 0 )
1379         {
1380         while ( true )
1381           {
1382           noSpaceTok = SimplifyWhitespace( theTpToken );
1383 
1384           switch (  GetCoord( noSpaceTok )  )
1385             {
1386             case 0:  this->Internal->XIdInList = this->NumberOfVariables; break;
1387             case 1:  this->Internal->YIdInList = this->NumberOfVariables; break;
1388             case 2:  this->Internal->ZIdInList = this->NumberOfVariables; break;
1389             default: break;
1390             }
1391 
1392           switch (  GuessCoord( noSpaceTok )  )
1393             {
1394             case 0:  guessedXid = this->NumberOfVariables; break;
1395             case 1:  guessedYid = this->NumberOfVariables; break;
1396             case 2:  guessedZid = this->NumberOfVariables; break;
1397             default: break;
1398             }
1399 
1400           this->Variables.push_back( theTpToken );
1401           this->NumberOfVariables ++;
1402 
1403           if ( this->Internal->NextCharEOL )
1404             {
1405             break;
1406             }
1407           theTpToken = this->Internal->GetNextToken();
1408           }
1409         }
1410 
1411       // in case there is not an exact match for coordinate axis vars
1412       this->Internal->XIdInList = ( this->Internal->XIdInList < 0 )
1413                                   ? guessedXid
1414                                   : this->Internal->XIdInList;
1415       this->Internal->YIdInList = ( this->Internal->YIdInList < 0 )
1416                                   ? guessedYid
1417                                   : this->Internal->YIdInList;
1418       this->Internal->ZIdInList = ( this->Internal->ZIdInList < 0 )
1419                                   ? guessedZid
1420                                   : this->Internal->ZIdInList;
1421 
1422       break;
1423       }
1424     else
1425       {
1426       do
1427         {
1428         theTpToken = this->Internal->GetNextToken();
1429         } while ( READ_UNTIL_TITLE_OR_VARIABLES );
1430 
1431       tokenReady = true;
1432       }
1433 
1434     if ( !tokenReady )
1435       {
1436       theTpToken = this->Internal->GetNextToken();
1437       }
1438 
1439     }
1440 
1441   this->Internal->ASCIIStream.rewind();
1442 
1443   // register the data arrays
1444   for ( i = 0; i < this->GetNumberOfDataAttributes(); i ++ )
1445     {
1446     // all data arrays are selected here by default
1447     this->DataArraySelection->EnableArray(  this->GetDataAttributeName( i )  );
1448     }
1449 }
1450 
1451 // ----------------------------------------------------------------------------
ReadFile(vtkMultiBlockDataSet * multZone)1452 void vtkTecplotReader::ReadFile( vtkMultiBlockDataSet * multZone )
1453 {
1454   if (    ( this->Internal->Completed == 1 )
1455        || ( this->FileName == NULL ) || (  strcmp( this->FileName, "" ) == 0  )
1456      )
1457     {
1458     return;
1459     }
1460 
1461   if ( multZone == NULL )
1462     {
1463     vtkErrorMacro( "vtkMultiBlockDataSet multZone NULL!" );
1464     return;
1465     }
1466 
1467   #define READ_UNTIL_LINE_END !this->Internal->NextCharEOF &&\
1468                               tok != "TITLE"               &&\
1469                               tok != "VARIABLES"           &&\
1470                               tok != "ZONE"                &&\
1471                               tok != "GEOMETRY"            &&\
1472                               tok != "TEXT"                &&\
1473                               tok != "DATASETAUXDATA"
1474   int  zoneIndex  = 0;
1475   bool firstToken = true;
1476   bool tokenReady = false;
1477 
1478   this->Init();
1479   this->Internal->ASCIIStream.open( this->FileName );
1480   std::string tok = this->Internal->GetNextToken();
1481 
1482   while ( !this->Internal->NextCharEOF )
1483     {
1484     tokenReady = false;
1485     if ( tok == "" )
1486       {
1487       // whitespace: do nothing
1488       }
1489     else if ( tok == "TITLE" )
1490       {
1491       this->DataTitle = this->Internal->GetNextToken();
1492       }
1493     else if ( tok == "GEOMETRY" )
1494       {
1495       // unsupported
1496       tok = this->Internal->GetNextToken();
1497       while ( READ_UNTIL_LINE_END )
1498         {
1499         // skipping token
1500         tok = this->Internal->GetNextToken();
1501         }
1502       tokenReady = true;
1503       }
1504     else if ( tok == "TEXT" )
1505       {
1506       // unsupported
1507       tok = this->Internal->GetNextToken();
1508       while ( READ_UNTIL_LINE_END )
1509         {
1510         // Skipping token
1511         tok = this->Internal->GetNextToken();
1512         }
1513       tokenReady = true;
1514       }
1515     else if ( tok == "VARIABLES" )
1516       {
1517       int guessedXindex = -1;
1518       int guessedYindex = -1;
1519       int guessedZindex = -1;
1520 
1521       // variable lists
1522       tok = this->Internal->GetNextToken();
1523       while ( this->Internal->TokenIsString )
1524         {
1525         int tokLen = int( tok.length() );
1526         for ( int i = 0; i < tokLen; i ++ )
1527           {
1528           if ( tok[i] == '(' )
1529             {
1530             tok[i] = '[';
1531             }
1532           else
1533           if ( tok[i] == ')' )
1534             {
1535             tok[i] = ']';
1536             }
1537           else
1538           if ( tok[i] == '/' )
1539             {
1540             tok[i] = '_';
1541             }
1542           }
1543 
1544         std::string   tok_nw = SimplifyWhitespace( tok );
1545 
1546         switch (  GetCoord( tok_nw )  )
1547           {
1548           case 0:  this->Internal->XIdInList = this->NumberOfVariables; break;
1549           case 1:  this->Internal->YIdInList = this->NumberOfVariables; break;
1550           case 2:  this->Internal->ZIdInList = this->NumberOfVariables; break;
1551           default: break;
1552           }
1553 
1554         switch (  GuessCoord( tok_nw )  )
1555           {
1556           case 0:  guessedXindex = this->NumberOfVariables; break;
1557           case 1:  guessedYindex = this->NumberOfVariables; break;
1558           case 2:  guessedZindex = this->NumberOfVariables; break;
1559           default: break;
1560           }
1561 
1562         this->Variables.push_back( tok );
1563         this->NumberOfVariables ++;
1564         tok = this->Internal->GetNextToken();
1565         }
1566 
1567       if ( this->NumberOfVariables == 0 )
1568         {
1569         while ( true )
1570           {
1571           std::string      tok_nw = SimplifyWhitespace( tok );
1572 
1573           switch (  GetCoord( tok_nw )  )
1574             {
1575             case 0:  this->Internal->XIdInList = this->NumberOfVariables; break;
1576             case 1:  this->Internal->YIdInList = this->NumberOfVariables; break;
1577             case 2:  this->Internal->ZIdInList = this->NumberOfVariables; break;
1578             default: break;
1579             }
1580 
1581           switch (  GuessCoord( tok_nw )  )
1582             {
1583             case 0:  guessedXindex = this->NumberOfVariables; break;
1584             case 1:  guessedYindex = this->NumberOfVariables; break;
1585             case 2:  guessedZindex = this->NumberOfVariables; break;
1586             default: break;
1587             }
1588 
1589           this->Variables.push_back( tok );
1590           this->NumberOfVariables ++;
1591 
1592           if ( this->Internal->NextCharEOL )
1593             {
1594             tok = this->Internal->GetNextToken();
1595             break;
1596             }
1597           else
1598             {
1599             tok = this->Internal->GetNextToken();
1600             }
1601           }
1602         }
1603 
1604       // Default the centering to nodal
1605       this->CellBased.clear();
1606       this->CellBased.resize( this->NumberOfVariables, 0 );
1607 
1608       // If we didn't find an exact match for coordinate axis vars, guess
1609       if ( this->Internal->XIdInList < 0 )
1610         {
1611         this->Internal->XIdInList = guessedXindex;
1612         }
1613       if ( this->Internal->YIdInList < 0 )
1614         {
1615         this->Internal->YIdInList = guessedYindex;
1616         }
1617       if ( this->Internal->ZIdInList < 0 )
1618         {
1619         this->Internal->ZIdInList = guessedZindex;
1620         }
1621 
1622       // Based on how many spatial coords we got, guess the spatial dim
1623       if ( this->Internal->XIdInList >= 0 )
1624         {
1625         this->Internal->GeometryDim = 1;
1626         if ( this->Internal->YIdInList >= 0 )
1627           {
1628           this->Internal->GeometryDim = 2;
1629           if ( this->Internal->ZIdInList >= 0 )
1630             {
1631             this->Internal->GeometryDim = 3;
1632             }
1633           }
1634         }
1635 
1636       tokenReady = true;
1637       }
1638     else if ( tok == "ZONE" )
1639       {
1640       int      numI = 1;
1641       int      numJ = 1;
1642       int      numK = 1;
1643       int      numNodes    = 0;
1644       int      numElements = 0;
1645       char     untitledZoneName[40];
1646       sprintf( untitledZoneName, "zone%05d", zoneIndex );
1647 
1648       std::string format    = "";
1649       std::string elemType  = "";
1650       std::string ZoneName = untitledZoneName;
1651 
1652       tok = this->Internal->GetNextToken();
1653       while (  !( tok != "T"  &&
1654                   tok != "I"  &&
1655                   tok != "J"  &&
1656                   tok != "K"  &&
1657                   tok != "N"  &&
1658                   tok != "E"  &&
1659                   tok != "ET" &&
1660                   tok != "F"  &&
1661                   tok != "D"  &&
1662                   tok != "DT" &&
1663                   tok != "STRANDID"     &&
1664                   tok != "SOLUTIONTIME" &&
1665                   tok != "DATAPACKING"  &&
1666                   tok != "VARLOCATION"
1667                 )
1668             )
1669         {
1670         if ( tok == "T" )
1671           {
1672           ZoneName = this->Internal->GetNextToken();
1673           if ( !this->Internal->TokenIsString )
1674             {
1675             vtkErrorMacro( << this->FileName << ": Zone titles MUST be "
1676                            << "quoted." );
1677             }
1678           }
1679         else if ( tok == "I" )
1680           {
1681           numI = atoi( this->Internal->GetNextToken().c_str() );
1682           }
1683         else if ( tok == "J" )
1684           {
1685           numJ = atoi( this->Internal->GetNextToken().c_str() );
1686           }
1687         else if ( tok == "K" )
1688           {
1689           numK = atoi( this->Internal->GetNextToken().c_str() );
1690           }
1691         else if ( tok == "N" )
1692           {
1693           numNodes = atoi( this->Internal->GetNextToken().c_str() );
1694           }
1695         else if ( tok == "E" )
1696           {
1697           numElements = atoi( this->Internal->GetNextToken().c_str() );
1698           }
1699         else if ( tok == "ET" )
1700           {
1701           elemType = this->Internal->GetNextToken();
1702           }
1703         else if ( tok == "F" || tok == "DATAPACKING" )
1704           {
1705           format = this->Internal->GetNextToken();
1706           }
1707         else if ( tok == "VARLOCATION" )
1708           {
1709           std::string  centering;
1710           this->CellBased.clear();
1711           this->CellBased.resize( this->NumberOfVariables, 0 );
1712 
1713           //read token to ascertain VARLOCATION syntax usage
1714           std::string var_format_type = this->Internal->GetNextToken();
1715 
1716           //if each variable will have data type specified explicitly (as is handled in old Tecplot reader),
1717           //else a range is specified for CELLCENTERED only, with NODAL values assumed implicitly
1718           if ( var_format_type == "NODAL" || var_format_type == "CELLCENTERED" )
1719             {
1720             if ( var_format_type == "CELLCENTERED" )
1721               {
1722               this->CellBased[0] = 1;
1723               }
1724             for ( int i = 1; i < this->NumberOfVariables; i ++ )
1725               {
1726               centering = this->Internal->GetNextToken();
1727               if ( centering == "CELLCENTERED" )
1728                 {
1729                 this->CellBased[i] = 1;
1730                 }
1731               }
1732             }
1733           else
1734             {
1735             do
1736               {
1737               //remove left square bracket, if it exists
1738               size_t brack_pos = var_format_type.find("[");
1739               if ( brack_pos != std::string::npos )
1740                 var_format_type.erase(brack_pos, brack_pos+1);
1741 
1742               //remove right square bracket, if it exists
1743               brack_pos = var_format_type.find("]");
1744               if ( brack_pos != std::string::npos )
1745                 var_format_type.erase(brack_pos,brack_pos+1);
1746 
1747               //if a range is defined, then split again, convert to int and set to cell data
1748               //else if a single value is defined, then just set the flag directly
1749               if ( var_format_type.find("-") != std::string::npos )
1750                 {
1751                 std::vector<std::string> var_range;
1752                 vtksys::SystemTools::Split(var_format_type.c_str(), var_range, '-');
1753 
1754                 int cell_start = atoi(var_range[0].c_str()) - 1;
1755                 int cell_end = atoi(var_range[1].c_str());
1756                 for ( int i = cell_start; i != cell_end; ++i )
1757                   {
1758                   this->CellBased[i] = 1;
1759                   }
1760                 }
1761               else
1762                 {
1763                 int index = atoi(var_format_type.c_str())-1;
1764                 this->CellBased[index] = 1;
1765                 }
1766 
1767               //get next value
1768               var_format_type = this->Internal->GetNextToken();
1769 
1770               //continue until the CELLCENTERED keyword is found
1771               } while ( var_format_type != "CELLCENTERED" );
1772             }
1773           }
1774         else if ( tok == "DT" )
1775           {
1776           for ( int i = 0; i < this->NumberOfVariables; i ++ )
1777             {
1778             this->Internal->GetNextToken();
1779             }
1780           }
1781         else if ( tok == "D" )
1782           {
1783           vtkErrorMacro( << this->FileName << "; Tecplot zone record parameter "
1784                          << "'D' is currently unsupported." );
1785           }
1786         else if ( tok == "STRANDID" )
1787           {
1788           vtkErrorMacro( << this->FileName << "; Tecplot zone record parameter "
1789                          << "'STRANDID' is currently unsupported." );
1790           this->Internal->GetNextToken();
1791           }
1792         else if ( tok == "SOLUTIONTIME" )
1793           {
1794           vtkErrorMacro( << this->FileName << "; Tecplot zone record parameter "
1795                          << "'SOLUTIONTIME' is currently unsupported." );
1796           this->Internal->GetNextToken();
1797           }
1798         tok = this->Internal->GetNextToken();
1799         }
1800 
1801       this->Internal->TokenBackup = tok;
1802 
1803       this->ZoneNames.push_back( ZoneName );
1804 
1805       if ( format == "FEBLOCK" )
1806         {
1807         this->GetUnstructuredGridFromBlockPackingZone( numNodes, numElements,
1808               elemType.c_str(), zoneIndex, ZoneName.c_str(),  multZone );
1809         }
1810       else if ( format == "FEPOINT" )
1811         {
1812         this->GetUnstructuredGridFromPointPackingZone( numNodes, numElements,
1813               elemType.c_str(), zoneIndex, ZoneName.c_str(),  multZone );
1814         }
1815       else if ( format == "BLOCK" )
1816         {
1817         this->GetStructuredGridFromBlockPackingZone
1818               ( numI, numJ, numK, zoneIndex, ZoneName.c_str(),  multZone );
1819         }
1820       else if ( format == "POINT" )
1821         {
1822         this->GetStructuredGridFromPointPackingZone
1823               ( numI, numJ, numK, zoneIndex, ZoneName.c_str(),  multZone );
1824         }
1825       else if ( format == "" )
1826         {
1827         // No format given; we will assume we got a POINT format
1828         this->GetStructuredGridFromPointPackingZone
1829               ( numI, numJ, numK, zoneIndex, ZoneName.c_str(),  multZone );
1830         }
1831       else
1832         {
1833         // UNKNOWN FORMAT
1834         vtkErrorMacro( << this->FileName << ": The format " << format.c_str()
1835                        << " found in the file is unknown." );
1836         }
1837 
1838       zoneIndex ++;
1839       }
1840     else if( tok == "DATASETAUXDATA" )
1841       {
1842       int   tokIndex       = 0;
1843       bool  haveVectorExpr = false;
1844       tok = this->Internal->GetNextToken();
1845 
1846       while ( READ_UNTIL_LINE_END )
1847         {
1848         if( tokIndex == 0 )
1849           {
1850           haveVectorExpr = ( tok == "VECTOR" );
1851           }
1852         else if( tokIndex == 1 )
1853           {
1854           if( haveVectorExpr )
1855             {
1856             // Remove spaces
1857             std::string::size_type pos = tok.find( " " );
1858             while( pos != std::string::npos )
1859               {
1860               tok.replace( pos, 1, "" );
1861               pos = tok.find( " " );
1862               }
1863 
1864             // Look for '('
1865             pos = tok.find( "(" );
1866             if( pos != std::string::npos )
1867               {
1868               std::string  exprName(  tok.substr( 0, pos )  );
1869               std::string  exprDef (  tok.substr( pos, tok.size() - pos )  );
1870 
1871               exprDef.replace( 0, 1, "{" );
1872 
1873               // Replace ')' with '}'
1874               pos = exprDef.find( ")" );
1875               if( pos != std::string::npos )
1876                 {
1877                 exprDef.replace( pos, 1, "}" );
1878                 vtkDebugMacro( "Expr name = " << exprName.c_str()
1879                                << ", Expr def = " << exprDef.c_str() );
1880                 }
1881               }
1882             }
1883           }
1884 
1885         // Skipping token
1886         tok = this->Internal->GetNextToken();
1887         tokIndex ++;
1888         }
1889 
1890       tokenReady = true;
1891       }
1892     else if ( firstToken && this->Internal->TokenIsString )
1893       {
1894       // Robust: assume it's a title
1895       this->DataTitle = tok;
1896       }
1897     else
1898       {
1899       // UNKNOWN RECORD TYPE
1900       vtkErrorMacro( << this->FileName << ": The record type " << tok.c_str()
1901                      << " found in the file is unknown." );
1902       }
1903 
1904     firstToken = false;
1905     if ( !tokenReady )
1906       {
1907       tok = this->Internal->GetNextToken();
1908       }
1909     }
1910   this->Internal->ASCIIStream.close();
1911 
1912   if ( this->Internal->TopologyDim > this->Internal->GeometryDim )
1913     {
1914     this->Internal->TopologyDim = this->Internal->GeometryDim;
1915     }
1916 
1917   this->Internal->Completed = 1;
1918 }
1919