]> SALOME platform Git repositories - modules/visu.git/blob - src/CONVERTOR/VISU_TableReader.cxx
Salome HOME
Fix for the "0051899: curves are not shown in opened study" issue.
[modules/visu.git] / src / CONVERTOR / VISU_TableReader.cxx
1 // Copyright (C) 2007-2013  CEA/DEN, EDF R&D, OPEN CASCADE
2 //
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License.
7 //
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
11 // Lesser General Public License for more details.
12 //
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19
20 //  VISU OBJECT : interactive object for VISU entities implementation
21 //  File:
22 //  Author:
23 //  Module : VISU
24 //
25 #include "VISU_TableReader.hxx"
26
27 #include <QFileInfo>
28 #include <QString>
29 #include <QRegExp>
30 #include <QFile>
31 #include <QStringList>
32
33 #include <fstream>
34 #include <iostream>
35 #include <sstream>
36
37 #include <vtkPoints.h>
38 #include <vtkDoubleArray.h>
39 #include <vtkPointData.h>
40 #include <vtkCellData.h>
41 #include <vtkPolyData.h>
42 #include <vtkPassThroughFilter.h>
43
44 #include <vtkStructuredGrid.h>
45 #include <vtkStructuredGridGeometryFilter.h>
46
47 #include <Basics_Utils.hxx>
48
49 #ifdef _DEBUG_
50 static int MYDEBUG = 0;
51 #else
52 static int MYDEBUG = 0;
53 #endif
54
55
56 //---------------------------------------------------------------
57 int
58 VISU::TTable2D
59 ::Check()
60 {
61   if ( myRows.empty() ) 
62     return 0;
63
64   int iEnd = myRows[0].myValues.size();
65   if ( iEnd == 0 )
66     return 0;
67
68   if ( myColumnTitles.size() != iEnd ) 
69     myColumnTitles.resize( iEnd );
70
71   if ( myColumnUnits.size() != iEnd )
72     myColumnUnits.resize( iEnd );
73
74   int jEnd = myRows.size();
75   for ( int j = 0; j < jEnd; j++ )
76     if ( myRows[j].myValues.size() != iEnd )
77       return 0;
78
79   return 1;
80 }
81
82
83 //---------------------------------------------------------------
84 void
85 VISU::TTable2D
86 ::getColumns(VISU::TTable2D& theTable2D) const
87 {
88   TRows& aRows = theTable2D.myRows;
89   aRows.clear();
90   if ( myRows.empty() )
91     return;
92
93   int jEnd = myRows.size();
94
95   //Define Titles & Units
96   theTable2D.myColumnTitles.resize(jEnd);
97   theTable2D.myColumnUnits.resize(jEnd);
98   for ( int j = 0; j < jEnd; j++ ) {
99     theTable2D.myColumnTitles[j] = myRows[j].myTitle;
100     theTable2D.myColumnUnits[j] = myRows[j].myUnit;
101   }
102
103   //Define Rows
104   int iEnd = myRows[0].myValues.size();
105   for ( int i = 0; i < iEnd; i++ ) {
106     TRow aNewRow;
107     aNewRow.myTitle = myColumnTitles[i];
108     aNewRow.myUnit = myColumnUnits[i];
109     aNewRow.myValues.resize(jEnd);
110     for ( int j = 0; j < jEnd; j++ ) {
111       aNewRow.myValues[j] = myRows[j].myValues[i];
112     }
113     aRows.push_back(aNewRow);
114   }
115 }
116
117
118 //---------------------------------------------------------------
119 namespace
120 {
121   int getLine( std::ifstream& theStmIn, QString& theString )
122   {
123     char tmp;
124     std::ostringstream aStrOut;
125
126     while ( theStmIn.get( tmp ) ) {
127       aStrOut<<tmp;
128       if ( tmp == '\n' ) 
129         break;
130     }
131
132     aStrOut<<std::ends;
133     theString = aStrOut.str().c_str();
134
135     return !theStmIn.eof();
136   }
137
138   //=======================================================================
139   //function : findNextCell
140   //purpose  : auxilary for ImportCSVTable
141   //=======================================================================
142   bool findNextCell(std::ifstream& aStmIn, QString& aStr,
143                    QString& aCell, const char theSeparator)
144   {
145     aCell = "";
146     int index, tmpind = 0;
147     if( aStr.at(0) == theSeparator ) {
148       aStr.remove(0,1);
149       aStr = aStr.trimmed();
150       if(aStr.size()==0) return true;
151     }
152     QString aTmp = aStr;
153     if( aTmp.at(0) == '"' ) {
154       // find closed "
155       while( !aStmIn.eof() ) {
156         tmpind = aTmp.indexOf('"',1);
157         if( tmpind < 0 ) {
158           while( !aStmIn.eof() ) {
159             aCell.push_back(aTmp);
160             aCell.push_back('\n');
161             getLine( aStmIn, aTmp );
162             tmpind = aTmp.indexOf('"',1);
163             if( tmpind >= 0 ) {
164               break;
165             }
166           }
167         }
168         if(tmpind < 0)
169           return false;
170         // read next symbol
171         if( aTmp.at(tmpind+1) == '"' ) {
172           // "" is found => need to continue
173           aCell.push_back(aTmp.left(tmpind+1));
174           aTmp = aTmp.mid(tmpind+2);
175         }
176         else if( aTmp.at(tmpind+1) == theSeparator ) {
177           index = tmpind+1;
178           break;
179         }
180         else {
181           return false;
182         }
183       }
184     }
185     else {
186       // find index for separator
187       index = aTmp.indexOf( theSeparator );
188     }
189     if( index < 0 ) {
190       aCell += aTmp;
191       aStr = "";
192     }
193     else {
194       if(index>0) aCell += aTmp.left(index);
195       aStr = aTmp.mid(index).trimmed();
196     }
197     if( aCell.size()>0 && aCell.at(0) == '"' ) {
198       // remove first and last symbols
199       int last = aCell.size()-1;
200       aCell.remove(last,1);
201       aCell.remove(0,1);
202     }
203     // replace "" to " in aCell
204     QChar ctmp = '"';
205     QString tmp(ctmp);
206     tmp += ctmp;
207     index = aCell.indexOf(tmp);
208     while(index>=0) {
209       aCell.remove(index,1);
210       index = aCell.indexOf(tmp);
211     }
212     return true;
213   }
214
215 }
216
217
218 //---------------------------------------------------------------
219 void 
220 VISU::ImportTables( const char* theFileName, TTableContainer& theContainer,
221                     bool theFirstStrAsTitle)
222 {
223   std::ifstream aStmIn;
224   QFileInfo aFileInfo( theFileName );
225   if( !aFileInfo.isFile() || !aFileInfo.isReadable() || !aFileInfo.size() )
226     return;
227
228   QString tmp(theFileName);
229   tmp = tmp.trimmed();
230   tmp = tmp.right(3).trimmed();
231
232   if( tmp == QString("csv") ) {
233     const char separator = ',';
234     ImportCSVTable(theFileName, theContainer, theFirstStrAsTitle, separator);
235     return;
236   }
237
238   aStmIn.open( theFileName );
239   QString aTmp;
240   do {
241     // find beginning of table (tables are separated by empty lines)
242     while( getLine( aStmIn, aTmp ) && aTmp.trimmed() == "" );
243
244     PTableIDMapper aTableIDMapper( new TTableIDMapper() );
245     TTable2D& aTable2D = *aTableIDMapper;
246     if(MYDEBUG) std::cout << "New table is found" << std::endl;
247
248     bool IsFirst = true;
249     while( !aStmIn.eof() && aTmp.trimmed() != "" ){
250       QString data = aTmp.trimmed();
251       QString cmt = "";
252       QString keyword = "";
253       // split string to data and comment (comment starts from '#' symbol)
254       int index = aTmp.indexOf( "#" );
255       if ( index >= 0 ) {
256         data = aTmp.left( index ).trimmed();
257         cmt = aTmp.mid( index+1 ).trimmed();
258       }
259       // if comment is not empty, try to get keyword from it (separated by ':' symbol)
260       if ( !cmt.isEmpty() ) {
261         int index1 = cmt.indexOf( ":" );
262         if ( index1 >= 0 ) {
263           QString tmpstr = cmt.left( index1 ).trimmed();
264           if ( tmpstr == QString( "TITLE" ) ||
265                tmpstr == QString( "COLUMN_TITLES" ) ||
266                tmpstr == QString( "COLUMN_UNITS" ) ||
267                tmpstr == QString( "COMMENT" ) ) {
268             keyword = tmpstr;
269             cmt = cmt.mid( index1+1 ).trimmed();
270           }
271         }
272       }
273       // if data is empty, process only comment
274       if ( data.isEmpty() ) {
275         // if keyword is found, try to process it
276         // elsewise it is a simple comment, just ignore it
277         if ( !keyword.isEmpty() ) {
278           if ( keyword == QString( "TITLE" ) ) {
279             QString title = cmt;
280             if ( aTable2D.myTitle != "" )
281               title = QString( aTable2D.myTitle.c_str() ) + QString( " " ) + title;
282             if(MYDEBUG) std::cout << "...Table TITLE is: " << title.toLatin1().constData() << std::endl;
283             aTable2D.myTitle = title.toLatin1().constData();
284           }
285           else if ( keyword == QString( "COLUMN_TITLES" ) ) {
286             // comment may contain column headers
287             QStringList aStrList = cmt.split( "|", QString::SkipEmptyParts );
288             if(MYDEBUG) std::cout << "...Column TITLES are: ";
289             for ( int i = 0; i < aStrList.count(); i++ ) {
290               QString tmpstr = aStrList[ i ].trimmed();
291               if(MYDEBUG) std::cout << tmpstr.toLatin1().constData() << " ";
292               aTable2D.myColumnTitles.push_back( tmpstr.toLatin1().constData() );
293             }
294             if(MYDEBUG) std::cout << std::endl;
295           }
296           else if ( keyword == QString( "COLUMN_UNITS" ) ) {
297             // comment may contain column units
298             QStringList aStrList = cmt.split( " ", QString::SkipEmptyParts );
299             if(MYDEBUG) std::cout << "...Column UNITS are: ";
300             for ( int i = 0; i < aStrList.count(); i++ ) {
301               QString tmpstr = aStrList[ i ].trimmed();
302               if(MYDEBUG) std::cout << tmpstr.toLatin1().constData() << " ";
303               aTable2D.myColumnUnits.push_back( tmpstr.toLatin1().constData() );
304             }
305             if(MYDEBUG) std::cout << std::endl;
306           }
307           else if ( keyword == QString( "COMMENT" ) ) {
308             // keyword 'COMMENT' processing can be here
309             // currently it is ignored
310             if(MYDEBUG) std::cout << "...COMMENT: " << cmt.toLatin1().constData() << std::endl;
311           }
312         }
313         else {
314           if(MYDEBUG) std::cout << "...comment: " << cmt.toLatin1().constData() << std::endl;
315           // simple comment processing can be here
316           // currently it is ignored
317         }
318       }
319       // if data is not empty, try to process it
320       else {
321         TTable2D::TRow aRow;
322         if(MYDEBUG) std::cout << "...New row is found: " << std::endl;
323         QString datar1 = data.replace(QRegExp("\t"), " ");
324         QStringList aValList = datar1.split( " ", QString::SkipEmptyParts );
325         if( aTable2D.myColumnTitles.size()==0 && IsFirst && theFirstStrAsTitle ) {
326           for ( int i = 0; i < aValList.count(); i++ ) {
327             QString tmpstr = aValList[ i ].trimmed();
328             aTable2D.myColumnTitles.push_back( tmpstr.toLatin1().constData() );
329           }
330         }
331         else {
332           if ( !cmt.isEmpty() ) {
333             aRow.myTitle = cmt.toLatin1().constData();
334             if(MYDEBUG) std::cout << "......ROW TITLE is: " << cmt.toLatin1().constData() << std::endl;
335           }
336           //QString datar1 = data.replace(QRegExp("\t"), " ");
337           //QStringList aValList = datar1.split( " ", QString::SkipEmptyParts );
338           for ( int i = 0; i < aValList.count(); i++ ) {
339             if ( aValList[i].trimmed() != "" ) {
340               TTable2D::TValue aVal = aValList[i].trimmed().toLatin1().constData();
341               aRow.myValues.push_back( aVal );
342             }
343           }
344           if( aRow.myValues.size() > 0 )
345             aTable2D.myRows.push_back( aRow );
346         }
347         IsFirst = false;
348         // ************** OLD CODE ******************
349         /*
350         TValue aVal;
351         istringstream aStream( data );
352         aStream.precision( STRPRECISION );
353         while( aStream >> aVal ) {
354           aRow.myValues.push_back( aVal );
355         }
356         if( aRow.myValues.size() > 0 )
357           aTable2D.myRows.push_back( aRow );
358         */
359         // ************** OLD CODE ******************
360       }
361       getLine( aStmIn, aTmp );
362     }
363     if( aTable2D.Check() ) {
364       if(MYDEBUG) std::cout << "aTable2D is checked OK " << aTable2D.myTitle << std::endl;
365       theContainer.push_back( aTableIDMapper );
366     }
367   } while ( !aStmIn.eof() );
368   aStmIn.close();
369
370   if(MYDEBUG) std::cout << "After close" << std::endl;
371 }
372
373
374 //=======================================================================
375 //function : ImportCSVTable
376 //purpose  : 
377 //=======================================================================
378 void VISU::ImportCSVTable(const char* theFileName, TTableContainer& theContainer,
379                           bool theFirstStrAsTitle, const char theSeparator)
380 {
381   std::ifstream aStmIn;
382   QFileInfo aFileInfo( theFileName );
383   if( !aFileInfo.isFile() || !aFileInfo.isReadable() || !aFileInfo.size() )
384     return;
385   aStmIn.open( theFileName );
386   QString aTmp;
387   do {
388     // find beginning of table (tables are separated by empty lines)
389     while( getLine( aStmIn, aTmp ) && aTmp.trimmed() == "" );
390
391     PTableIDMapper aTableIDMapper( new TTableIDMapper() );
392     TTable2D& aTable2D = *aTableIDMapper;
393     if(MYDEBUG) std::cout << "New table is found" << std::endl;
394
395     bool IsFirst = true;
396     QStringList aValList;
397
398     while( !aStmIn.eof() ) {
399       QString aCell = "";
400       if( !( findNextCell(aStmIn, aTmp, aCell, theSeparator)) ) {
401         return;
402       }
403       if( aTmp.size()==0 ) {
404         // make table row
405         aValList.push_back(aCell);
406         if( IsFirst && theFirstStrAsTitle ) {
407           for ( int i = 0; i < aValList.count(); i++ ) {
408             aTable2D.myColumnTitles.push_back( aValList[i].trimmed().toLatin1().constData() );
409           }
410         }
411         else {
412           TTable2D::TRow aRow;
413           for ( int i = 0; i < aValList.count(); i++ ) {
414             if ( aValList[i].trimmed() != "" ) {
415               TTable2D::TValue aVal = aValList[i].trimmed().toLatin1().constData();
416               aRow.myValues.push_back( aVal );
417             }
418             else {
419               aRow.myValues.push_back( "Empty" );
420             }
421           }
422           if( aRow.myValues.size() > 0 ) {
423             aTable2D.myRows.push_back( aRow );
424           }
425         }
426         // clear list of values and read next string
427         aValList.clear();
428         getLine( aStmIn, aTmp );
429         IsFirst = false;
430       }
431       else {
432         // put value to table cell
433         aValList.push_back(aCell);
434       }
435     }
436
437     if( aTable2D.Check() ) {
438       if(MYDEBUG) std::cout << "aTable2D is checked OK " << aTable2D.myTitle << std::endl;
439       theContainer.push_back( aTableIDMapper );
440     }
441
442   } while ( !aStmIn.eof() );
443   aStmIn.close();
444
445   if(MYDEBUG) std::cout << "After close" << std::endl;
446 }
447
448
449 //---------------------------------------------------------------
450 VISU::TTableIDMapper
451 ::TTableIDMapper():
452   myOutput( vtkPolyData::New() ),
453   myFilter( vtkPassThroughFilter::New() ),
454   myXAxisPosition( -1 )
455 {
456 }
457
458 VISU::TTableIDMapper
459 ::~TTableIDMapper()
460 {
461   myOutput->Delete();
462   myFilter->Delete();
463 }
464
465 vtkPolyData*
466 VISU::TTableIDMapper
467 ::GetPolyDataOutput()
468 {
469   if ( myXAxisPosition == -1 )
470     SetXAxisPosition( 0 );
471
472   return myOutput;
473 }
474
475 vtkAlgorithmOutput*
476 VISU::TTableIDMapper
477 ::GetOutputPort()
478 {
479   return myFilter->GetOutputPort();
480 }
481
482 long unsigned int
483 VISU::TTableIDMapper
484 ::GetMemorySize()
485 {
486   return myOutput->GetActualMemorySize() * 1024;
487 }
488
489 void
490 VISU::TTableIDMapper
491 ::SetXAxisPosition( vtkIdType theAxisPosition )
492 {
493   // Set "C" numeric locale to import numbers correctly
494   Kernel_Utils::Localizer loc;
495
496   if ( myXAxisPosition == theAxisPosition || !Check() )
497     return;
498
499   myOutput->Initialize();
500
501   if ( !Check() )
502     return;
503
504   TTable2D aTable2D;
505   getColumns( aTable2D );
506   
507   vtkIdType aXSize = aTable2D.myRows[0].myValues.size();
508
509   // It is necessary to decrease the size at 1 take intoa account X axis
510   vtkIdType anYSize = aTable2D.myRows.size() - 1; 
511
512   vtkIdType aNbPoints = aXSize * anYSize;
513
514   std::vector<double> anXAxis(aXSize);
515   const TTable2D::TValues& aValues = aTable2D.myRows[theAxisPosition].myValues;
516   for ( vtkIdType aX = 0; aX < aXSize; aX++ )
517     anXAxis[aX] = atof( aValues[aX].c_str() );
518
519   double aXRange = anXAxis[aXSize - 1] - anXAxis[0];
520   double anYDelta = aXRange / anYSize;
521   std::vector<double> anYAxis(anYSize);
522   for ( vtkIdType anY = 0; anY < anYSize; anY++ )
523     anYAxis[anY] = anY * anYDelta;
524
525   vtkPoints* aPoints = vtkPoints::New();
526   aPoints->SetNumberOfPoints( aNbPoints );
527
528   vtkIntArray *aPointsIDMapper = vtkIntArray::New();
529   aPointsIDMapper->SetName("VISU_POINTS_MAPPER");
530   aPointsIDMapper->SetNumberOfComponents(2);
531   aPointsIDMapper->SetNumberOfTuples(aNbPoints);
532   int *aPointsIDMapperPtr = aPointsIDMapper->GetPointer(0);
533
534   //vtkIntArray *aCellIDMapper = vtkIntArray::New();
535   //aCellIDMapper->SetName("VISU_POINTS_MAPPER");
536   //aCellIDMapper->SetNumberOfComponents(2);
537   //aCellIDMapper->SetNumberOfTuples(aNbPoints);
538   //int *aCellIDMapperPtr = aCellIDMapper->GetPointer(0);
539
540   for ( vtkIdType aY = 0, aPntId = 0; aY < anYSize; aY++ ) {
541     for ( vtkIdType aX = 0; aX < aXSize; aX++, aPntId++ ) {
542       aPoints->SetPoint( aPntId, anXAxis[aX], anYAxis[aY], 0.0 );
543
544       *aPointsIDMapperPtr++ = aPntId;
545       *aPointsIDMapperPtr++ = 0;
546
547       //*aCellIDMapperPtr++ = aPntId;
548       //*aCellIDMapperPtr++ = 0;
549     }
550   }
551
552   std::vector<TValues> anYData;
553   for ( vtkIdType anY = 0; anY < anYSize + 1; anY++ ) {
554     if ( anY == theAxisPosition )
555       continue;
556     anYData.push_back( aTable2D.myRows[anY].myValues );
557   }
558
559   vtkDoubleArray* aScalars = vtkDoubleArray::New();
560   aScalars->SetNumberOfComponents( 1 );
561   aScalars->SetNumberOfTuples( aNbPoints );
562   double *aScalarsPtr = aScalars->GetPointer(0);
563   for ( vtkIdType anY = 0; anY < anYSize; anY++ ) {
564     const TTable2D::TValues& aValues = anYData[anY];
565     for ( vtkIdType aX = 0; aX < aXSize; aX++ ) {
566       double aValue = atof( aValues[aX].c_str() );
567       *aScalarsPtr++ = aValue;
568     }
569   }
570
571   vtkStructuredGrid* aStructuredGrid = vtkStructuredGrid::New();
572   aStructuredGrid->SetPoints( aPoints );
573   aPoints->Delete();
574
575   aStructuredGrid->SetDimensions( aXSize, anYSize, 1 );
576
577   aStructuredGrid->GetPointData()->AddArray( aPointsIDMapper );
578   aPointsIDMapper->Delete();
579
580   //aStructuredGrid->GetCellData()->AddArray( aCellIDMapper );
581   //aCellIDMapper->Delete();
582
583   aStructuredGrid->GetPointData()->SetScalars( aScalars );
584   aScalars->Delete();
585
586   vtkStructuredGridGeometryFilter* aFilter = vtkStructuredGridGeometryFilter::New();
587   aFilter->SetInputData( aStructuredGrid );
588   aFilter->Update();
589   myOutput->ShallowCopy( aFilter->GetOutput() );
590   aFilter->Delete();
591
592   myFilter->SetInputData( myOutput );
593 }
594