Salome HOME
Merge from BR_PORTING_VTK6 01/03/2013
[modules/paravis.git] / src / Plugins / IntegrationPoints / vtkELNOSurfaceFilter.cxx
1 // Copyright (C) 2010-2012  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 "vtkELNOSurfaceFilter.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 "vtkPVGeometryFilter.h"
31 #include "vtkShrinkFilter.h"
32 #include "vtkSmartPointer.h"
33 #include "vtkPointData.h"
34 #include "vtkCellData.h"
35 #include "vtkIdList.h"
36
37 //vtkCxxRevisionMacro(vtkELNOSurfaceFilter, "$Revision$")
38 //;
39 vtkStandardNewMacro(vtkELNOSurfaceFilter)
40 ;
41
42 vtkELNOSurfaceFilter::vtkELNOSurfaceFilter()
43 {
44 }
45
46 vtkELNOSurfaceFilter::~vtkELNOSurfaceFilter()
47 {
48 }
49
50 int vtkELNOSurfaceFilter::RequestData(vtkInformation *request,
51     vtkInformationVector **input, vtkInformationVector *output)
52 {
53   vtkUnstructuredGrid *usgIn=vtkUnstructuredGrid::SafeDownCast(
54       input[0]->GetInformationObject(0)->Get(vtkDataObject::DATA_OBJECT()));
55
56   vtkUnstructuredGrid *usgOut=vtkUnstructuredGrid::SafeDownCast(
57       output->GetInformationObject(0)->Get(vtkDataObject::DATA_OBJECT()));
58
59   vtkIdTypeArray* usg_offsets=vtkIdTypeArray::SafeDownCast(
60       this->GetInputArrayToProcess(0, input));
61
62   if(usgIn==NULL||usg_offsets==NULL||usgOut==NULL)
63     {
64     vtkDebugMacro("vtkELNOSurfaceFilter no correctly configured : offsets = " << usg_offsets);
65     return 1;
66     }
67
68   // first shrink the input
69   vtkUnstructuredGrid* usgInClone=usgIn->NewInstance();
70
71   usgInClone->ShallowCopy(usgIn);
72
73   vtkSmartPointer<vtkPVGeometryFilter> geomFilter=vtkSmartPointer<
74       vtkPVGeometryFilter>::New();
75   geomFilter->SetInputData(usgInClone);
76   geomFilter->SetPassThroughCellIds(1);
77   geomFilter->SetPassThroughPointIds(1);
78   geomFilter->SetUseOutline(0);
79   geomFilter->Update();
80
81   vtkPolyData* surface=vtkPolyData::SafeDownCast(geomFilter->GetOutput());
82   vtkIdTypeArray* originalCellIds=vtkIdTypeArray::SafeDownCast(
83       surface->GetCellData()->GetArray("vtkOriginalCellIds"));
84   vtkIdTypeArray* originalPointIds=vtkIdTypeArray::SafeDownCast(
85       surface->GetPointData()->GetArray("vtkOriginalPointIds"));
86
87   if( originalCellIds == NULL )
88   {
89     vtkErrorMacro("vtkPVGeometryFilter return NULL 'vtkOriginalCellIds' array");
90     return 0;
91   }
92
93   if(originalPointIds==NULL)
94   {
95     vtkErrorMacro("vtkPVGeometryFilter return NULL 'vtkOriginalPointIds' array");
96     return 0;
97   }
98
99   vtkSmartPointer<vtkShrinkFilter> shrink=
100       vtkSmartPointer<vtkShrinkFilter>::New();
101   shrink->SetInputConnection(geomFilter->GetOutputPort(0));
102   shrink->SetShrinkFactor(0.9999);
103   shrink->Update();
104
105   vtkUnstructuredGrid* shrinked=shrink->GetOutput();
106
107   usgInClone->Delete();
108
109   usgOut->ShallowCopy(shrinked);
110
111   vtkIdTypeArray* offsets=vtkIdTypeArray::SafeDownCast(
112       shrinked->GetCellData()->GetArray(usg_offsets->GetName()));
113
114   // now copy ELNO data. Start by verifying if it is possible to
115   // shallow copy the array.
116   vtkInformation *info=offsets->GetInformation();
117   vtkInformationQuadratureSchemeDefinitionVectorKey *key=
118       vtkQuadratureSchemeDefinition::DICTIONARY();
119   if(!key->Has(info))
120     {
121     vtkDebugMacro("Dictionary is not present in array " << offsets->GetName()
122                   << " " << offsets << " Aborting." );
123     return 0;
124     }
125   int dictSize=key->Size(info);
126   vtkQuadratureSchemeDefinition **dict=
127       new vtkQuadratureSchemeDefinition *[dictSize];
128   key->GetRange(info, dict, 0, 0, dictSize);
129
130   vtkIdType ncell=shrinked->GetNumberOfCells();
131
132   vtkFieldData* fielddata=usgIn->GetFieldData();
133   vtkIdList *ids=vtkIdList::New();
134   vtkIdList *surfaceIds=vtkIdList::New();
135   vtkIdList *originalIds=vtkIdList::New();
136   for(int index=0; index<fielddata->GetNumberOfArrays(); index++)
137     {
138     vtkDataArray* data=fielddata->GetArray(index);
139     if(data==NULL)
140       continue;
141
142     vtkInformation* info=data->GetInformation();
143     const char* arrayOffsetName=info->Get(
144         vtkQuadratureSchemeDefinition::QUADRATURE_OFFSET_ARRAY_NAME());
145
146     if(arrayOffsetName == NULL ||
147        strcmp(arrayOffsetName, offsets->GetName())!=0)
148       {
149       usgOut->GetFieldData()->AddArray(data);
150
151       continue;
152       }
153
154     vtkDataArray* newArray=data->NewInstance();
155     newArray->SetName(data->GetName());
156     usgOut->GetPointData()->AddArray(newArray);
157     newArray->SetNumberOfComponents(data->GetNumberOfComponents());
158     newArray->SetNumberOfTuples(usgOut->GetNumberOfPoints());
159     newArray->CopyComponentNames(data);
160     newArray->Delete();
161
162     for(vtkIdType cellId=0; cellId<ncell; cellId++)
163       {
164       vtkIdType offset=offsets->GetValue(cellId);
165
166       vtkIdType originalCellId=originalCellIds->GetValue(cellId);
167       int originalCellType=usgIn->GetCellType(originalCellId);
168
169       shrinked->GetCellPoints(cellId, ids);
170       surface->GetCellPoints(cellId, surfaceIds);
171
172       for(int id=0; id<ids->GetNumberOfIds(); id++)
173         {
174         vtkIdType surfaceId=surfaceIds->GetId(id);
175         vtkIdType shrinkedId=ids->GetId(id);
176         vtkIdType originalPointId = originalPointIds->GetValue(surfaceId);
177
178         usgIn->GetCellPoints(originalCellId, originalIds);
179         int originalLocalId=-1;
180         for(int li=0; li<originalIds->GetNumberOfIds(); li++)
181           {
182           if(originalPointId==originalIds->GetId(li))
183             {
184             originalLocalId=li;
185             break;
186             }
187           }
188         if(originalLocalId==-1)
189           {
190           originalLocalId=0;
191           vtkErrorMacro("cannot find original id");
192           }
193
194         const double * w=dict[originalCellType]->GetShapeFunctionWeights(
195             originalLocalId);
196         int j;
197         for(j=0; j<dict[originalCellType]->GetNumberOfNodes(); j++)
198           {
199           if(w[j]==1.0)
200             break;
201           }
202         if(j==dict[originalCellType]->GetNumberOfNodes())
203           {
204           vtkErrorMacro("cannot find elno weigth.");
205           j=id;
206           }
207         newArray->SetTuple(shrinkedId, offset+j, data);
208         }
209       }
210     }
211
212   ids->FastDelete();
213   surfaceIds->FastDelete();
214   originalIds->FastDelete();
215   delete[] dict;
216
217   return 1;
218 }
219
220 void vtkELNOSurfaceFilter::PrintSelf(ostream& os, vtkIndent indent)
221 {
222   this->Superclass::PrintSelf(os, indent);
223 }