Salome HOME
307004c5feacdb5bdd4447ac021f34e5f0931533
[modules/paravis.git] / src / Plugins / MEDReader / plugin / MEDReaderIO / vtkMEDQuadraturePointsGenerator.cxx
1 // Copyright (C) 2010-2022  CEA/DEN, EDF R&D
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, or (at your option) any later version.
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 // Author : Roman NIKOLAEV
20
21 //Local includes
22 #include "vtkMEDQuadraturePointsGenerator.h"
23 #include "MEDFileFieldRepresentationTree.hxx"
24
25 //VTK includes
26 #include <vtkObjectFactory.h>
27 #include <vtkInformation.h>
28 #include <vtkInformationVector.h>
29 #include <vtkUnstructuredGrid.h>
30 #include <vtkCellData.h>
31 #include <vtkPointData.h>
32 #include <vtkInformationQuadratureSchemeDefinitionVectorKey.h>
33 #include <vtkQuadratureSchemeDefinition.h>
34
35 //-----------------------------------------------------------------------------
36 vtkStandardNewMacro(vtkMEDQuadraturePointsGenerator)
37
38 //-----------------------------------------------------------------------------
39 vtkMEDQuadraturePointsGenerator::vtkMEDQuadraturePointsGenerator()
40 {
41 }
42
43 //-----------------------------------------------------------------------------
44 vtkMEDQuadraturePointsGenerator::~vtkMEDQuadraturePointsGenerator()
45 {}
46
47
48 //-----------------------------------------------------------------------------
49 int vtkMEDQuadraturePointsGenerator::RequestData(
50       vtkInformation* request,
51       vtkInformationVector **input,
52       vtkInformationVector *output)
53 {
54   if (this->Superclass::RequestData(request, input, output) == 0 )
55     {
56       return 0;
57     }    
58
59   //Fill MED internal array
60   vtkDataObject *tmpDataObj;
61   
62   // Get the input.
63   tmpDataObj = input[0]->GetInformationObject(0)->Get(vtkDataObject::DATA_OBJECT());
64   vtkUnstructuredGrid *usgIn = vtkUnstructuredGrid::SafeDownCast(tmpDataObj);
65   // Get the output.
66   tmpDataObj = output->GetInformationObject(0)->Get(vtkDataObject::DATA_OBJECT());
67   vtkPolyData *pdOut = vtkPolyData::SafeDownCast(tmpDataObj);
68   vtkDataArray* offsets = this->GetInputArrayToProcess(0, input);
69
70   if (usgIn == NULL || pdOut == NULL || offsets == NULL )
71     {
72     vtkErrorMacro("Filter data has not been configured correctly. Aborting.");
73     return 1;
74     }
75   
76   vtkInformation *info = offsets->GetInformation();
77   vtkInformationQuadratureSchemeDefinitionVectorKey *key 
78     = vtkQuadratureSchemeDefinition::DICTIONARY();
79   if (!key->Has(info))
80     {
81     vtkErrorMacro(
82     << "Dictionary is not present in array "
83     << offsets->GetName() << " " << offsets
84     << " Aborting.");
85     return 0;
86     }
87
88   vtkIdType nCells = usgIn->GetNumberOfCells();
89   int dictSize = key->Size(info);
90   vtkQuadratureSchemeDefinition **dict
91     = new vtkQuadratureSchemeDefinition *[dictSize];
92   key->GetRange(info, dict, 0, 0, dictSize);  
93
94   // Loop over all fields to map the internal MED cell array to the points array
95   int nCArrays = usgIn->GetCellData()->GetNumberOfArrays();
96   for (int i = 0; i<nCArrays; ++i)
97     {
98     vtkDataArray* array = usgIn->GetCellData()->GetArray(i);
99     if ( !array )
100       {
101       continue;
102       }
103     std::string arrName = array->GetName();
104     if ( arrName == MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_CELL_NAME ) 
105       { 
106         vtkDataArray *out_id_cells = array->NewInstance();
107         out_id_cells->SetName(MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME);
108         out_id_cells->SetNumberOfComponents(array->GetNumberOfComponents());
109         out_id_cells->CopyComponentNames( array );
110         for (int cellId = 0; cellId < nCells; cellId++)
111           {
112           int cellType = usgIn->GetCellType(cellId);
113
114           // a simple check to see if a scheme really exists for this cell type.
115           // should not happen if the cell type has not been modified.
116           if (dict[cellType] == NULL)
117             {
118             continue;
119             }
120
121           int np = dict[cellType]->GetNumberOfQuadraturePoints();
122           for (int id = 0; id < np; id++)
123             {
124             out_id_cells->InsertNextTuple(cellId, array);
125             }
126           }
127         out_id_cells->Squeeze();
128         pdOut->GetPointData()->AddArray(out_id_cells);
129         out_id_cells->Delete();
130       }   
131     }
132   delete[] dict;  
133   return 1;
134 }