Salome HOME
Minimal CORBA mode: adapted PARAVIS' Python API to work in this mode.
[modules/paravis.git] / src / Plugins / IntegrationPoints / vtkELNOMeshFilter.cxx
1 // Copyright (C) 2010-2013  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.
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 #include "vtkELNOMeshFilter.h"
21 #include "vtkInformation.h"
22 #include "vtkInformationVector.h"
23 #include "vtkObjectFactory.h"
24 #include "vtkPolyDataAlgorithm.h"
25 #include "vtkPolyData.h"
26 #include "vtkIdTypeArray.h"
27 #include "vtkInformationQuadratureSchemeDefinitionVectorKey.h"
28 #include "vtkQuadratureSchemeDefinition.h"
29 #include "vtkUnstructuredGrid.h"
30 #include "vtkShrinkFilter.h"
31 #include "vtkSmartPointer.h"
32 #include "vtkPointData.h"
33 #include "vtkCellData.h"
34 #include "vtkIdList.h"
35
36 //vtkCxxRevisionMacro(vtkELNOMeshFilter, "$Revision$")
37 //;
38 vtkStandardNewMacro(vtkELNOMeshFilter)
39 ;
40
41 vtkELNOMeshFilter::vtkELNOMeshFilter()
42 {
43 }
44
45 vtkELNOMeshFilter::~vtkELNOMeshFilter()
46 {
47 }
48
49 int vtkELNOMeshFilter::RequestData(vtkInformation *request,
50     vtkInformationVector **input, vtkInformationVector *output)
51 {
52   vtkUnstructuredGrid *usgIn = vtkUnstructuredGrid::SafeDownCast(
53       input[0]->GetInformationObject(0)->Get(vtkDataObject::DATA_OBJECT()));
54
55   vtkUnstructuredGrid *usgOut = vtkUnstructuredGrid::SafeDownCast(
56       output->GetInformationObject(0)->Get(vtkDataObject::DATA_OBJECT()));
57
58   vtkIdTypeArray* usg_offsets = vtkIdTypeArray::SafeDownCast(
59       this->GetInputArrayToProcess(0, input));
60
61   if(usgIn == NULL || usg_offsets == NULL || usgOut == NULL)
62     {
63     vtkDebugMacro("vtkELNOMeshFilter no correctly configured : offsets = " << usg_offsets);
64     return 1;
65     }
66
67   vtkIdTypeArray* a = vtkIdTypeArray::SafeDownCast(
68       usgIn->GetCellData()->GetArray(usg_offsets->GetName()));
69
70   vtkInformationVector *inArrayVec =
71       this->Information->Get(INPUT_ARRAYS_TO_PROCESS());
72
73   // first shrink the input
74   vtkUnstructuredGrid* usgInClone = usgIn->NewInstance();
75
76   usgInClone->ShallowCopy(usgIn);
77
78   vtkSmartPointer<vtkShrinkFilter> shrink =
79       vtkSmartPointer<vtkShrinkFilter>::New();
80   shrink->SetInputData(usgInClone);
81   shrink->SetShrinkFactor(0.9999);
82   shrink->Update();
83   vtkUnstructuredGrid* shrinked = shrink->GetOutput();
84
85   usgInClone->Delete();
86
87   usgOut->ShallowCopy(shrinked);
88
89   vtkIdTypeArray* offsets = vtkIdTypeArray::SafeDownCast(
90       shrinked->GetCellData()->GetArray(usg_offsets->GetName()));
91
92   // now copy ELNO data. Start by verifying if it is possible to
93   // shallow copy the array.
94   vtkInformation *info = offsets->GetInformation();
95   vtkInformationQuadratureSchemeDefinitionVectorKey *key =
96       vtkQuadratureSchemeDefinition::DICTIONARY();
97   if(!key->Has(info))
98     {
99     vtkDebugMacro("Dictionary is not present in array " << offsets->GetName() << " " << offsets << " Aborting." );
100     return 0;
101     }
102   int dictSize = key->Size(info);
103   vtkQuadratureSchemeDefinition **dict =
104       new vtkQuadratureSchemeDefinition *[dictSize];
105   key->GetRange(info, dict, 0, 0, dictSize);
106
107   vtkIdType nVerts = shrinked->GetNumberOfPoints();
108   vtkIdType ncell = usgIn->GetNumberOfCells();
109   // first loop through all cells to check if a shallow copy is possible
110   bool shallowok = true;
111   vtkIdType previous = 0;
112
113   for(vtkIdType cellId = 0; cellId < ncell; cellId++)
114     {
115     vtkIdType offset = offsets->GetValue(cellId);
116     if(offset != previous)
117       {
118       shallowok = false;
119       break;
120       }
121     int cellType = usgIn->GetCellType(cellId);
122
123     if(dict[cellType] == NULL)
124       {
125       previous = offset;
126       }
127     else
128       {
129       previous = offset + dict[cellType]->GetNumberOfQuadraturePoints();
130       }
131     }
132   if(previous != nVerts)
133     {
134     shallowok = false;
135     }
136
137   vtkFieldData* fielddata = usgIn->GetFieldData();
138   for(int index = 0; index < fielddata->GetNumberOfArrays(); index++)
139     {
140     vtkDataArray* data = fielddata->GetArray(index);
141     if(data == NULL)
142       continue;
143
144     vtkInformation* info = data->GetInformation();
145     const char* arrayOffsetName = info->Get(
146         vtkQuadratureSchemeDefinition::QUADRATURE_OFFSET_ARRAY_NAME());
147
148     if(arrayOffsetName == NULL ||
149        strcmp(arrayOffsetName, offsets->GetName()) != 0)
150       {
151       shrinked->GetFieldData()->AddArray(data);
152       continue;
153       }
154
155     if(shallowok)
156       {
157       usgOut->GetPointData()->AddArray(data);
158       }
159     else
160       {
161       vtkDataArray* newArray = data->NewInstance();
162       newArray->SetName(data->GetName());
163       usgOut->GetPointData()->AddArray(newArray);
164       newArray->SetNumberOfComponents(data->GetNumberOfComponents());
165       newArray->SetNumberOfTuples(usgOut->GetNumberOfPoints());
166       newArray->CopyComponentNames(data);
167       newArray->Delete();
168       vtkIdList *ids = vtkIdList::New();
169
170       for(vtkIdType cellId = 0; cellId < ncell; cellId++)
171         {
172         vtkIdType offset = offsets->GetValue(cellId);
173         int cellType = shrinked->GetCellType(cellId);
174         shrinked->GetCellPoints(cellId, ids);
175         for(int id = 0; id < dict[cellType]->GetNumberOfQuadraturePoints(); id++)
176           {
177           const double * w = dict[cellType]->GetShapeFunctionWeights(id);
178           int j;
179           for(j = 0; j < dict[cellType]->GetNumberOfNodes(); j++)
180             {
181             if(w[j] == 1.0)
182               break;
183             }
184           if(j == dict[cellType]->GetNumberOfNodes())
185             {
186             j = id;
187             }
188           newArray->SetTuple(ids->GetId(id), offset + j, data);
189           }
190         }
191       ids->FastDelete();
192       }
193     }
194
195   delete[] dict;
196
197   return 1;
198 }
199
200 void vtkELNOMeshFilter::PrintSelf(ostream& os, vtkIndent indent)
201 {
202   this->Superclass::PrintSelf(os, indent);
203 }