source: XMLIO_V2/dev/common/src/xmlio/vtk/lscereader.cpp @ 219

Last change on this file since 219 was 219, checked in by hozdoba, 13 years ago

Préparation nouvelle arborescence

File size: 27.0 KB
Line 
1#include "lscereader.hpp"
2
3// BibliothÚque VTK
4#include "vtkObjectFactory.h"
5#include "vtkInformation.h"
6#include "vtkInformationVector.h"
7#include "vtkStreamingDemandDrivenPipeline.h"
8
9#include <vtkProperty.h>
10#include <vtkDataSetMapper.h>
11#include <vtkActor.h>
12#include <vtkRenderWindow.h>
13#include <vtkRenderer.h>
14#include <vtkRenderWindowInteractor.h>
15#include <vtkLegendScaleActor.h>
16
17#include "vtkMath.h"
18
19namespace xmlioserver
20{
21   namespace vtk
22   {
23      /// ////////////////////// Définitions ////////////////////// ///
24      vtkStandardNewMacro(vtkLSCEReader);
25
26      vtkLSCEReader::vtkLSCEReader(void)
27         : FileName(), CurGridType(RECTILINEAR), VarNames()
28         , A3D(true), ATemporal(true), ACell(true)
29         , Reader()
30      {
31         this->SetNumberOfInputPorts (0);
32         this->SetNumberOfOutputPorts(1);
33      }
34
35      vtkLSCEReader::~vtkLSCEReader(void)
36      {
37
38      }
39
40      ///---------------------------------------------------------------
41
42      void vtkLSCEReader::SetFileName(const StdString & fileName)
43      {
44         this->FileName = fileName;
45         this->Reader   = boost::shared_ptr<io::CINetCDF4>(new io::CINetCDF4(fileName));
46         this->Modified();
47      }
48
49      const StdString & vtkLSCEReader::GetFileName(void) const
50      {
51         return (this->FileName);
52      }
53
54      void vtkLSCEReader::PrintSelf(ostream& os, vtkIndent indent)
55      {
56        this->Superclass::PrintSelf(os,indent);
57      }
58
59      //----------------------------------------------------------------
60      void vtkLSCEReader::AddVariableToSelection(const StdString & varName)
61      {
62         if (this->Reader->isRectilinear(varName))
63            this->SetGridType (RECTILINEAR);
64         else if (this->Reader->isCurvilinear(varName))
65            this->SetGridType (CURVILINEAR);
66         else
67            this->SetGridType (UNSTRUCTURED);
68
69         this->VarNames.insert(varName);
70      }
71
72      void vtkLSCEReader::RemoveSelectedVariable(const StdString & varName)
73      {
74         std::set<StdString>::iterator it = this->VarNames.find(varName);
75         if (it != this->VarNames.end()) this->VarNames.erase(it);
76      }
77
78      void vtkLSCEReader::RemoveAllSelectedVariables(void)
79      {
80         this->VarNames.clear();
81      }
82
83      const std::set<StdString> &
84            vtkLSCEReader::GetSelectedVariables(void) const
85      {
86         return (this->VarNames);
87      }
88
89      bool vtkLSCEReader::HasSelectedVariable(void) const
90      {
91         return (this->VarNames.size() != 0);
92      }
93
94      //----------------------------------------------------------------
95
96      void vtkLSCEReader::SetGridType (GridType type)
97      {
98         if ((this->IsUnstructured() && (type == CURVILINEAR)) ||
99             (this->IsUnstructured() && (type == RECTILINEAR)) ||
100             (this->IsCurvilinear()  && (type == RECTILINEAR)))
101            this->RemoveAllSelectedVariables();
102         this->CurGridType = type;
103      }
104
105      bool vtkLSCEReader::IsUnstructured(void) const
106      {
107         return (this->CurGridType == UNSTRUCTURED);
108      }
109
110      bool vtkLSCEReader::IsCurvilinear(void) const
111      {
112         return (this->CurGridType == CURVILINEAR);
113      }
114
115      bool vtkLSCEReader::IsRectilinear(void) const
116      {
117         return (this->CurGridType == RECTILINEAR);
118      }
119
120      //----------------------------------------------------------------
121
122      void vtkLSCEReader::AcceptTemporalOnly(bool value)
123      {
124         this->ATemporal = value;
125      }
126
127      void vtkLSCEReader::Accept3DOnly(bool value)
128      {
129         this->A3D = value;
130      }
131
132      void vtkLSCEReader::AcceptCellOnly(bool value)
133      {
134         this->ACell = value;
135      }
136
137      //----------------------------------------------------------------
138
139      int vtkLSCEReader::ProcessRequest
140                           (vtkInformation       *  request,
141                            vtkInformationVector ** vtkNotUsed(inputVector),
142                            vtkInformationVector *  outputVector)
143      {
144         if(request->Has(vtkDemandDrivenPipeline::REQUEST_DATA()))
145            return (this->RequestData(NULL, NULL, outputVector));
146
147         //if(request->Has(vtkDemandDrivenPipeline::REQUEST_DATA_OBJECT()))
148         //   return (this->RequestDataObject(NULL, NULL, outputVector));
149
150         //if(request->Has(vtkDemandDrivenPipeline::REQUEST_INFORMATION()))
151         //   return (this->RequestInformation(NULL, NULL, outputVector));
152
153         if(request->Has(vtkStreamingDemandDrivenPipeline::REQUEST_UPDATE_EXTENT()))
154            return (this->RequestUpdateExtent(NULL, NULL, outputVector));
155
156         return (1);
157      }
158
159      //----------------------------------------------------------------
160
161      void vtkLSCEReader::ShowVariable(const vtkStdString & filename,
162                                       const vtkStdString & varname)
163      {
164         vtkSmartPointer<vtkLSCEReader> lscereader =
165               vtkSmartPointer<vtkLSCEReader>::New();
166
167         lscereader->SetFileName(filename);
168         lscereader->AddVariableToSelection(varname);
169         lscereader->Update();
170
171         // Création d'un mappeur et d'un acteur.
172         vtkSmartPointer<vtkDataSetMapper> mapper = vtkSmartPointer<vtkDataSetMapper>::New();
173         mapper->SetInputConnection(lscereader->GetOutput()->GetProducerPort());
174         mapper->SetScalarRange  (1000, 10000);
175
176         vtkSmartPointer<vtkActor> actor = vtkSmartPointer<vtkActor>::New();
177         actor->SetMapper(mapper);
178         actor->GetProperty()->SetRepresentationToWireframe();
179
180         // Création d'une fenêtre de rendu.
181         vtkSmartPointer<vtkRenderer> renderer = vtkSmartPointer<vtkRenderer>::New();
182         vtkSmartPointer<vtkRenderWindow> renderWindow = vtkSmartPointer<vtkRenderWindow>::New();
183         renderWindow->AddRenderer(renderer);
184         vtkSmartPointer<vtkRenderWindowInteractor> renderWindowInteractor =
185            vtkSmartPointer<vtkRenderWindowInteractor>::New();
186         renderWindowInteractor->SetRenderWindow(renderWindow);
187
188
189         vtkSmartPointer<vtkLegendScaleActor> legendScaleActor =
190            vtkSmartPointer<vtkLegendScaleActor>::New();
191
192         // Ajout de l'acteur à la scÚne.
193         renderer->AddActor(actor);
194         renderer->AddActor(legendScaleActor);
195
196         renderer->GradientBackgroundOn();
197         renderer->SetBackground(1,1,1);
198         renderer->SetBackground2(0,0,0);
199
200         // Rendu et intéraction
201         renderWindow->Render();
202         renderWindowInteractor->Start();
203      }
204
205      //----------------------------------------------------------------
206
207      void vtkLSCEReader::CreateRectilinearGrid(vtkRectilinearGrid * grid,
208                                                vtkInformation     * outInfo,
209                                                vtkFloatArray      * xspacing,
210                                                vtkFloatArray      * yspacing,
211                                                vtkFloatArray      * zspacing,
212                                                vtkIntArray        * dimensions)
213      {
214            int extent[6] = { 0, dimensions->GetValue(0) - 1
215                            , 0, dimensions->GetValue(1) - 1
216                            , 0, dimensions->GetValue(2) - 1};
217            outInfo->Set(vtkStreamingDemandDrivenPipeline::WHOLE_EXTENT(),extent,6);
218            grid->SetExtent (extent);
219            grid->SetDimensions (dimensions->GetValue(0),
220                                 dimensions->GetValue(1),
221                                 dimensions->GetValue(2));
222            grid->SetXCoordinates (xspacing);
223            grid->SetYCoordinates (yspacing);
224            grid->SetZCoordinates (zspacing);
225      }
226
227      void vtkLSCEReader::CreateStructuredGrid(vtkStructuredGrid * grid,
228                                               vtkInformation    * outInfo,
229                                               vtkPoints         * points,
230                                               vtkIntArray       * dimensions)
231      {
232            int extent[6] = { 0, dimensions->GetValue(0) - 1
233                            , 0, dimensions->GetValue(1) - 1
234                            , 0, dimensions->GetValue(2) - 1};
235            outInfo->Set(vtkStreamingDemandDrivenPipeline::WHOLE_EXTENT(),extent,6);
236            grid->SetExtent (extent);
237            grid->SetDimensions (dimensions->GetValue(0),
238                                 dimensions->GetValue(1),
239                                 dimensions->GetValue(2));
240            grid->SetPoints(points);
241      }
242
243      void vtkLSCEReader::CreateUnstructuredGrid(vtkUnstructuredGrid * grid,
244                                                 vtkInformation      * outInfo,
245                                                 vtkPoints           * points,
246                                                 vtkCellArray        * cells,
247                                                 int                   celltype)
248      {
249         grid->SetPoints(points);
250         grid->SetCells(celltype, cells);
251      }
252
253      //----------------------------------------------------------------
254
255      void vtkLSCEReader::GetRectilinearConnectivity(int dimZ, int dimY, int dimX, vtkCellArray * cells)
256      {
257         if (dimZ == 1 && dimY >  1 && dimX >  1) // 2D
258         {
259            vtkIdType pts[4] ;
260
261            for (int i = 0, k = 0; i < (dimX - 1); i++)
262               for (int j = 0; j < (dimY - 1); j++, k++)
263            {
264               pts[0] =  i * dimY + j;
265               pts[1] =  i * dimY + j + 1;
266               pts[2] =  (i+1) * dimY + j + 1;
267               pts[3] =  (i+1) * dimY + j;
268               cells->InsertNextCell(4,  pts);
269            }
270         }
271         else if (dimZ > 1 && dimY >  1 && dimX >  1) // 3D
272         {
273            vtkIdType pts[8];
274            for (int i = 0, l = 0; i < (dimX - 1); i++)
275               for (int j = 0; j < (dimY - 1); j++)
276                  for (int k = 0; k < (dimZ - 1); k++, l++)
277            {
278               pts[0] =  i * dimY * dimZ + j * dimZ + k;
279               pts[1] =  i * dimY * dimZ + j * dimZ + k + 1;
280               pts[2] =  (i+1) * dimY * dimZ + j * dimZ + k + 1;
281               pts[3] =  (i+1) * dimY * dimZ + j * dimZ + k;
282
283               pts[4] =  i * dimY * dimZ  + (j + 1) * dimZ + k ;
284               pts[5] =  i * dimY * dimZ  + (j + 1) * dimZ + k +1;
285               pts[6] =  (i+1) * dimY * dimZ + (j + 1) * dimZ + k + 1;
286               pts[7] =  (i+1) * dimY * dimZ + (j + 1) * dimZ + k ;
287               cells->InsertNextCell(8,  pts);
288            }
289         }
290      }
291
292      //----------------------------------------------------------------
293
294      void vtkLSCEReader::CreateSimpleGrid (int xi, int xf,
295                                            int yi, int yf,
296                                            int zi, int zf,
297                                            vtkFloatArray * xspacing,
298                                            vtkFloatArray * yspacing,
299                                            vtkFloatArray * zspacing,
300                                            vtkIntArray   * dimensions)
301      {
302         dimensions->InsertNextValue (xf - xi + 1); // X
303         dimensions->InsertNextValue (yf - yi + 1); // Y
304         dimensions->InsertNextValue (zf - zi + 1); // Z
305         for (int i = xi; i <= xf; i++)
306            xspacing->InsertNextValue (i);
307         for (int j = yi; j <= yf; j++)
308            yspacing->InsertNextValue (j);
309         for (int k = zi; k <= zf; k++)
310            zspacing->InsertNextValue (k);
311      }
312
313      void vtkLSCEReader::CreateSimpleGrid (int xi, int xf,
314                                            int yi, int yf,
315                                            int zi, int zf,
316                                            vtkPoints   * points,
317                                            vtkIntArray * dimensions)
318      {
319         dimensions->InsertNextValue (xf - xi + 1); // X
320         dimensions->InsertNextValue (yf - yi + 1); // Y
321         dimensions->InsertNextValue (zf - zi + 1); // Z
322         for (int i = xi; i <= xf; i++)
323            for (int j = yi; j <= yf; j++)
324               for (int k = zi; k <= zf; k++)
325                  points->InsertNextPoint(i, j, k);
326      }
327
328      void vtkLSCEReader::CreateSimpleGrid (int xi, int xf,
329                                            int yi, int yf,
330                                            int zi, int zf,
331                                            vtkPoints    * points,
332                                            vtkCellArray * cells,
333                                            vtkIntArray  * dimensions)
334      {
335         this->CreateSimpleGrid (xi, xf, yi, yf, zi, zf, points, dimensions);
336         int dimX = dimensions->GetValue(0),
337             dimY = dimensions->GetValue(1),
338             dimZ = dimensions->GetValue(2);
339         this->GetRectilinearConnectivity(dimZ, dimY, dimX, cells);
340      }
341
342      //----------------------------------------------------------------
343      void vtkLSCEReader::GetSpacings(const vtkStdString & coordinate,
344                                      bool bounds, vtkFloatArray * spacing)
345      {
346         StdSize size = (Reader->getDimensions(&coordinate).begin())->second;
347         ARRAY_CREATE(data_arr, float, 1, [1]);
348         if (Reader->hasVariable(coordinate))
349         {
350            if (bounds)
351            {
352               if (Reader->hasBounds(coordinate))
353               { // CellData, variable de coordonnée avec attributes bounds.
354                  StdString boundid = Reader->getBoundsId(coordinate);
355                  float f;
356                  Reader->getData(data_arr, boundid);
357                  spacing->InsertNextValue (f = (*data_arr)[0]);
358                  for (StdSize j = 1; j < (2 * size); j++) // << Super laid
359                     if (f != (*data_arr)[j])
360                        spacing->InsertNextValue (f = (*data_arr)[j]);
361               }
362               else
363               { // CellData, variable de coordonnée sans attributes bounds.
364                  Reader->getData(data_arr, coordinate);
365                  spacing->InsertNextValue ((*data_arr)[0] - ((*data_arr)[1] - (*data_arr)[0])/2);
366                  for (StdSize j = 0; j < (size-1); j++)
367                     spacing->InsertNextValue ((*data_arr)[j] + ((*data_arr)[j+1]-(*data_arr)[j])/2);
368                  spacing->InsertNextValue ((*data_arr)[size-2] + ((*data_arr)[size-1] - (*data_arr)[size-2])/2);
369               }
370            }
371            else
372            { // PointData, variable de coordonnée.
373               Reader->getData(data_arr, coordinate);
374               for (StdSize j = 0; j < (size-1); j++)
375                  spacing->InsertNextValue ((*data_arr)[j]);
376            }
377         }
378         else
379         {
380            if (bounds) // CellData, dimension de coordonnée.
381               for (float i = -0.5; i < ((float)size + 1.); i = i + 1.)
382                  spacing->InsertNextValue (i);
383            else // PointData, dimension de coordonnée.
384               for (float j = 0.; j < (float)size; j = j + 1.)
385                  spacing->InsertNextValue (j);
386         }
387      }
388      //----------------------------------------------------------------
389
390      void vtkLSCEReader::AddPoint(vtkPoints * points, float *  value, bool proj)
391      {
392         if (proj)
393         {
394            double h = 1.0;
395            double cartesianCoord[3];
396            value[0] =  0.017453292519943295 * value[0];
397            value[1] =  0.017453292519943295 * value[1];
398
399            cartesianCoord[0] = h * cos(value[0]) * cos(value[1]);
400            cartesianCoord[1] = h * sin(value[0]) * cos(value[1]);
401            cartesianCoord[2] = h * sin(value[1]);
402            points->InsertNextPoint (cartesianCoord);
403         }
404         else
405         {
406            points->InsertNextPoint (value);
407         }
408      }
409
410      //----------------------------------------------------------------
411      void vtkLSCEReader::GetPoints(const vtkStdString & xcoordinate,
412                                    const vtkStdString & ycoordinate,
413                                    const vtkStdString & zcoordinate,
414                                    bool bounds, bool proj,
415                                    vtkPoints * points, vtkIntArray * dimensions)
416      {
417         float value[3];
418         bool is2D = true;
419         ARRAY_CREATE(xvalue, float, 1, [1]);
420         ARRAY_CREATE(yvalue, float, 1, [1]);
421         vtkSmartPointer<vtkFloatArray> zspacing = vtkSmartPointer<vtkFloatArray>::New();
422         if (zcoordinate.size() > 0)
423         {
424            GetSpacings(zcoordinate, bounds, zspacing);
425            is2D = false;
426         }
427         else
428            zspacing->InsertNextValue(0);
429
430         dimensions->InsertNextValue
431            (Reader->getDimensions(&xcoordinate)[Reader->getLonCoordName(xcoordinate)]);
432         dimensions->InsertNextValue
433            (Reader->getDimensions(&xcoordinate)[Reader->getLatCoordName(xcoordinate)]);
434         dimensions->InsertNextValue(zspacing->GetNumberOfTuples());
435
436         if (!bounds)
437         {
438            Reader->getData(xvalue, xcoordinate);
439            Reader->getData(yvalue, ycoordinate);
440
441            for (int i = 0; i < dimensions->GetValue(0); i++)
442               for (int j = 0; j < dimensions->GetValue(1); j++)
443                  for (int k = 0; k < dimensions->GetValue(2); k++)
444                  {
445                     value[0] = (*xvalue)[i*dimensions->GetValue(1) + j];
446                     value[1] = (*yvalue)[i*dimensions->GetValue(1) + j];
447                     value[2] = zspacing->GetValue(k);
448                     AddPoint(points, value, proj);
449                  }
450            return;
451         }
452
453         dimensions->SetValue(0, dimensions->GetValue(0)+1);
454         dimensions->SetValue(1, dimensions->GetValue(1)+1);
455
456         if (Reader->hasBounds(xcoordinate))
457         {
458            StdString boundid = Reader->getBoundsId(xcoordinate);
459            Reader->getData(xvalue, boundid);
460         }
461         if (Reader->hasBounds(ycoordinate))
462         {
463            StdString boundid = Reader->getBoundsId(ycoordinate);
464            Reader->getData(yvalue, boundid);
465         }
466         int ydim = dimensions->GetValue(0) - 1 ,
467             xdim = dimensions->GetValue(1) - 1;
468
469         for (int i = 0; i <= xdim; i++)
470            for (int j = 0; j <= ydim; j++)
471               for (int k = 0; k < dimensions->GetValue(2); k++)
472               {
473                  if ((j != ydim) && (i != xdim))
474                  {
475                     value[0] = (*xvalue)[4 * (i * ydim + j) + 0];
476                     value[1] = (*yvalue)[4 * (i * ydim + j) + 0];
477                     value[2] = zspacing->GetValue(k);
478                     AddPoint(points, value, proj);
479
480                     continue;
481                  }
482                  else if ((j == ydim) && (i != xdim))
483                  {
484                     value[0] = (*xvalue)[4 * (i * ydim + j - 1) + 1];
485                     value[1] = (*yvalue)[4 * (i * ydim + j - 1) + 1];
486                     value[2] = zspacing->GetValue(k);
487                     AddPoint(points, value, proj);
488
489                     continue;
490                  }
491                  else if ((j != ydim) && (i == xdim))
492                  {
493                     value[0] = (*xvalue)[4 * ((i-1) * ydim + j) + 3];
494                     value[1] = (*yvalue)[4 * ((i-1) * ydim + j) + 3];
495                     value[2] = zspacing->GetValue(k);
496                     AddPoint(points, value, proj);
497
498                     continue;
499                  }
500                  else if ((j == ydim) && (i == xdim))
501                  {
502                     value[0] = (*xvalue)[4 * ((i-1) * ydim + j - 1) + 2];
503                     value[1] = (*yvalue)[4 * ((i-1) * ydim + j - 1) + 2];
504                     value[2] = zspacing->GetValue(k);
505                     AddPoint(points, value, proj);
506
507                     continue;
508                  }
509               }
510      }
511
512      //----------------------------------------------------------------
513
514      void vtkLSCEReader::GetCellsAndPoints(const vtkStdString & xcoordinate,
515                                            const vtkStdString & ycoordinate,
516                                            const vtkStdString & zcoordinate,
517                                            bool bounds, bool proj, StdSize nbvertex,
518                                            vtkCellArray * cells, vtkPoints * points,
519                                            vtkIntArray * dimensions)
520      {
521         float value[3];
522         double h = 1.0;
523         double cartesianCoord[3];
524         vtkIdType pts[nbvertex] ;
525         bool is2D = true;
526         ARRAY_CREATE(xvalue, float, 1, [1]);
527         ARRAY_CREATE(yvalue, float, 1, [1]);
528         vtkSmartPointer<vtkFloatArray> zspacing = vtkSmartPointer<vtkFloatArray>::New();
529         if (zcoordinate.size() > 0)
530         {
531            GetSpacings(zcoordinate, bounds, zspacing);
532            is2D = false;
533         }
534         else
535            zspacing->InsertNextValue(0);
536
537         StdSize size = (*Reader->getDimensions(&xcoordinate).begin()).second;
538         if (Reader->hasBounds(xcoordinate))
539         {
540            StdString boundid = Reader->getBoundsId(xcoordinate);
541            Reader->getData(xvalue, boundid);
542         }
543         if (Reader->hasBounds(ycoordinate))
544         {
545            StdString boundid = Reader->getBoundsId(ycoordinate);
546            Reader->getData(yvalue, boundid);
547         }
548
549         for (StdSize u = 0; u < size; u++)
550         {
551            for (int v = 0; v < zspacing->GetNumberOfTuples(); v++)
552            {
553               for (StdSize w = 0; w < nbvertex; w++)
554               {
555                  value[0] = (*xvalue)[nbvertex * u + w];
556                  value[1] = (*yvalue)[nbvertex * u + w];
557                  value[2] = zspacing->GetValue(v);
558                  pts[w] = nbvertex * u + w;
559
560                  if (proj)
561                  {
562                     value[0] = 0.017453292519943295 * value[0];
563                     value[1] = 0.017453292519943295 * value[1];
564
565                     cartesianCoord[0] = h * cos(value[0]) * cos(value[1]);
566                     cartesianCoord[1] = h * sin(value[0]) * cos(value[1]);
567                     cartesianCoord[2] = h * sin(value[1]);
568                     points->InsertPoint (pts[w], cartesianCoord);
569                  }
570                  else
571                  {
572                     points->InsertPoint (pts[w], value);
573                  }
574               }
575               cells->InsertNextCell(nbvertex,  pts);
576            }
577         }
578      }
579
580      //----------------------------------------------------------------
581
582      void vtkLSCEReader::AddScalarData(vtkDataSet * output,
583                                        const StdString & varname,
584                                        StdSize record, bool bounds)
585      {
586         vtkSmartPointer<vtkDoubleArray> data = vtkDoubleArray::New();
587         data->SetName(varname.c_str());
588         ARRAY_CREATE(data_arr, double, 1, [1]);
589         Reader->getData(data_arr, varname);
590
591         data->SetNumberOfValues (data_arr->size());
592         for (StdSize i = 0; i < data_arr->size(); i++)
593            data->SetValue (i, (*data_arr)[i]);
594
595         if (bounds)
596         {
597            output->GetCellData()->AddArray(data);
598            output->GetCellData()->SetActiveScalars (varname.c_str());
599         }
600         else
601         {
602            output->GetPointData()->AddArray(data);
603            output->GetPointData()->SetActiveScalars (varname.c_str());
604         }
605      }
606
607      //----------------------------------------------------------------
608
609      int vtkLSCEReader::RequestData
610                           (vtkInformation       *  vtkNotUsed(request),
611                            vtkInformationVector ** vtkNotUsed(inputVector),
612                            vtkInformationVector *  outputVector)
613      {
614         vtkInformation *  outInfo = outputVector->GetInformationObject(0);
615         vtkDataSet     *  output  = vtkDataSet::GetData(outInfo);
616
617         if (!output)
618         {
619            switch (this->CurGridType)
620            {
621               case(RECTILINEAR)  :
622                  output = vtkRectilinearGrid::New();
623                  break;
624               case(CURVILINEAR)  :
625                  output = vtkStructuredGrid::New();
626                  break;
627               default :
628                  output = vtkUnstructuredGrid::New();
629                  break;
630            }
631            output->SetPipelineInformation(outInfo);
632            output->Delete();
633         }
634
635         vtkUnstructuredGrid * uoutput  = vtkUnstructuredGrid::SafeDownCast(output);
636         vtkRectilinearGrid  * routput  = vtkRectilinearGrid::SafeDownCast(output);
637         vtkStructuredGrid   * soutput  = vtkStructuredGrid::SafeDownCast(output);
638
639         if (routput)
640         { // Grille rectiliniéaire
641            vtkSmartPointer<vtkFloatArray>
642                           xspacing = vtkSmartPointer<vtkFloatArray>::New(),
643                           yspacing = vtkSmartPointer<vtkFloatArray>::New(),
644                           zspacing = vtkSmartPointer<vtkFloatArray>::New();
645            vtkSmartPointer<vtkIntArray> dims = vtkSmartPointer<vtkIntArray>::New();
646
647            //CreateSimpleGrid (-1, 10, -1, 10, -1, 10, xspacing, yspacing, zspacing, dims);
648            GetSpacings(Reader->getLonCoordName(*VarNames.begin()), true, xspacing);
649            GetSpacings(Reader->getLatCoordName(*VarNames.begin()), true, yspacing);
650            GetSpacings(Reader->getVertCoordName(*VarNames.begin()), true, zspacing);
651            /*std::cout << zspacing->GetNumberOfTuples() << std::endl;*/
652
653            dims->InsertNextValue(xspacing->GetNumberOfTuples());
654            dims->InsertNextValue(yspacing->GetNumberOfTuples());
655            dims->InsertNextValue(zspacing->GetNumberOfTuples());
656
657            CreateRectilinearGrid(routput, outInfo, xspacing, yspacing, zspacing, dims);
658         }
659
660         if (soutput)
661         { // Grille structurée
662            vtkSmartPointer<vtkPoints>   pts  = vtkSmartPointer<vtkPoints>::New();
663            vtkSmartPointer<vtkIntArray> dims = vtkSmartPointer<vtkIntArray>::New();
664
665            //CreateSimpleGrid (-1, 10, -1, 10, -1, 10, pts, dims);
666            GetPoints(Reader->getLonCoordName(*VarNames.begin()),
667                      Reader->getLatCoordName(*VarNames.begin()),
668                      Reader->getVertCoordName(*VarNames.begin()), false, true, pts, dims);
669            CreateStructuredGrid(soutput, outInfo, pts, dims);
670         }
671
672         if (uoutput)
673         { // Grille non-structurée
674            vtkSmartPointer<vtkPoints>    pts   = vtkSmartPointer<vtkPoints>::New();
675            vtkSmartPointer<vtkIntArray>  dims  = vtkSmartPointer<vtkIntArray>::New();
676            vtkSmartPointer<vtkCellArray> cells = vtkSmartPointer<vtkCellArray>::New();
677            // VTK_QUAD, VTK_HEXAHEDRON, VTK_POLYGON
678
679            //CreateSimpleGrid (0, 12, 0, 12, 0, 12, pts, cells, dims);
680            GetCellsAndPoints(Reader->getLonCoordName(*VarNames.begin()),
681                              Reader->getLatCoordName(*VarNames.begin()),
682                              Reader->getVertCoordName(*VarNames.begin()), true, true,
683                              Reader->getNbVertex(*VarNames.begin()), cells, pts, dims);
684
685            AddScalarData(uoutput, *VarNames.begin(), 0, true);
686            CreateUnstructuredGrid(uoutput, outInfo, pts, cells, VTK_POLYGON);
687
688         }
689         return (1);
690      }
691
692      ///---------------------------------------------------------------
693
694
695   } // namespace vtk
696} // namespace xmlioserver
697
Note: See TracBrowser for help on using the repository browser.