Salome HOME
Update copyrights
[modules/paravis.git] / src / Plugins / MEDReader / IO / vtkELNOSurfaceFilter.cxx
1 // Copyright (C) 2010-2019  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
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, 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("vtkELNOSurfaceFilter no correctly configured : offsets = " << usg_offsets);
64     return 1;
65     }
66
67   // first shrink the input
68   vtkUnstructuredGrid* usgInClone=usgIn->NewInstance();
69
70   usgInClone->ShallowCopy(usgIn);
71
72   vtkSmartPointer<vtkPVGeometryFilter> geomFilter=vtkSmartPointer<
73       vtkPVGeometryFilter>::New();
74   geomFilter->SetInputData(usgInClone);
75   geomFilter->SetPassThroughCellIds(1);
76   geomFilter->SetPassThroughPointIds(1);
77   geomFilter->SetUseOutline(0);
78   geomFilter->Update();
79
80   vtkPolyData* surface=vtkPolyData::SafeDownCast(geomFilter->GetOutput());
81   vtkIdTypeArray* originalCellIds=vtkIdTypeArray::SafeDownCast(
82       surface->GetCellData()->GetArray("vtkOriginalCellIds"));
83   vtkIdTypeArray* originalPointIds=vtkIdTypeArray::SafeDownCast(
84       surface->GetPointData()->GetArray("vtkOriginalPointIds"));
85
86   if( originalCellIds == NULL )
87   {
88     vtkErrorMacro("vtkPVGeometryFilter return NULL 'vtkOriginalCellIds' array");
89     return 0;
90   }
91
92   if(originalPointIds==NULL)
93   {
94     vtkErrorMacro("It appears that your dataset is not reduced using vtkPVGeometryFilter (NULL 'vtkOriginalPointIds).\n==================================================================================================\nProbably your dataset is not 3D.\nIf it is not a 3D dataset you are expected to use ELNO Mesh filter instead of ELNO Surface filter.\n==================================================================================================\n");
95     return 0;
96   }
97
98   vtkSmartPointer<vtkShrinkFilter> shrink=
99       vtkSmartPointer<vtkShrinkFilter>::New();
100   shrink->SetInputConnection(geomFilter->GetOutputPort(0));
101   shrink->SetShrinkFactor(0.9999);
102   shrink->Update();
103
104   vtkUnstructuredGrid* shrinked=shrink->GetOutput();
105
106   usgInClone->Delete();
107
108   usgOut->ShallowCopy(shrinked);
109
110   vtkIdTypeArray* offsets=vtkIdTypeArray::SafeDownCast(
111       shrinked->GetCellData()->GetArray(usg_offsets->GetName()));
112
113   // now copy ELNO data. Start by verifying if it is possible to
114   // shallow copy the array.
115   vtkInformation *info=offsets->GetInformation();
116   vtkInformationQuadratureSchemeDefinitionVectorKey *key=
117       vtkQuadratureSchemeDefinition::DICTIONARY();
118   if(!key->Has(info))
119     {
120     vtkDebugMacro("Dictionary is not present in array " << offsets->GetName()
121                   << " " << offsets << " Aborting." );
122     return 0;
123     }
124   int dictSize=key->Size(info);
125   vtkQuadratureSchemeDefinition **dict=
126       new vtkQuadratureSchemeDefinition *[dictSize];
127   key->GetRange(info, dict, 0, 0, dictSize);
128
129   vtkIdType ncell=shrinked->GetNumberOfCells();
130
131   vtkFieldData* fielddata=usgIn->GetFieldData();
132   vtkIdList *ids=vtkIdList::New();
133   vtkIdList *surfaceIds=vtkIdList::New();
134   vtkIdList *originalIds=vtkIdList::New();
135   for(int index=0; index<fielddata->GetNumberOfArrays(); index++)
136     {
137     vtkDataArray* data=fielddata->GetArray(index);
138     if(data==NULL)
139       continue;
140
141     vtkInformation* info=data->GetInformation();
142     const char* arrayOffsetName=info->Get(
143         vtkQuadratureSchemeDefinition::QUADRATURE_OFFSET_ARRAY_NAME());
144
145     if(arrayOffsetName == NULL ||
146        strcmp(arrayOffsetName, offsets->GetName())!=0)
147       {
148       usgOut->GetFieldData()->AddArray(data);
149
150       continue;
151       }
152
153     vtkDataArray* newArray=data->NewInstance();
154     newArray->SetName(data->GetName());
155     usgOut->GetPointData()->AddArray(newArray);
156     newArray->SetNumberOfComponents(data->GetNumberOfComponents());
157     newArray->SetNumberOfTuples(usgOut->GetNumberOfPoints());
158     newArray->CopyComponentNames(data);
159     newArray->Delete();
160
161     for(vtkIdType cellId=0; cellId<ncell; cellId++)
162       {
163       vtkIdType offset=offsets->GetValue(cellId);
164
165       vtkIdType originalCellId=originalCellIds->GetValue(cellId);
166       int originalCellType=usgIn->GetCellType(originalCellId);
167
168       shrinked->GetCellPoints(cellId, ids);
169       surface->GetCellPoints(cellId, surfaceIds);
170
171       for(int id=0; id<ids->GetNumberOfIds(); id++)
172         {
173         vtkIdType surfaceId=surfaceIds->GetId(id);
174         vtkIdType shrinkedId=ids->GetId(id);
175         vtkIdType originalPointId = originalPointIds->GetValue(surfaceId);
176
177         usgIn->GetCellPoints(originalCellId, originalIds);
178         int originalLocalId=-1;
179         for(int li=0; li<originalIds->GetNumberOfIds(); li++)
180           {
181           if(originalPointId==originalIds->GetId(li))
182             {
183             originalLocalId=li;
184             break;
185             }
186           }
187         if(originalLocalId==-1)
188           {
189           originalLocalId=0;
190           vtkErrorMacro("cannot find original id");
191           }
192
193         const double * w=dict[originalCellType]->GetShapeFunctionWeights(
194             originalLocalId);
195         int j;
196         for(j=0; j<dict[originalCellType]->GetNumberOfNodes(); j++)
197           {
198           if(w[j]==1.0)
199             break;
200           }
201         if(j==dict[originalCellType]->GetNumberOfNodes())
202           {
203             //vtkErrorMacro("cannot find elno weigth.");
204           j=id;
205           }
206         newArray->SetTuple(shrinkedId, offset+j, data);
207         }
208       }
209     }
210
211   ids->FastDelete();
212   surfaceIds->FastDelete();
213   originalIds->FastDelete();
214   delete[] dict;
215
216   return 1;
217 }
218
219 void vtkELNOSurfaceFilter::PrintSelf(ostream& os, vtkIndent indent)
220 {
221   this->Superclass::PrintSelf(os, indent);
222 }