Salome HOME
158d8369ade5a948d149d12482fe93bd8999eda8
[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->SetInput(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   vtkSmartPointer<vtkShrinkFilter> shrink=
88       vtkSmartPointer<vtkShrinkFilter>::New();
89   shrink->SetInput(surface);
90   shrink->SetShrinkFactor(0.9999);
91   shrink->Update();
92
93   vtkUnstructuredGrid* shrinked=shrink->GetOutput();
94
95   usgInClone->Delete();
96
97   usgOut->ShallowCopy(shrinked);
98
99   vtkIdTypeArray* offsets=vtkIdTypeArray::SafeDownCast(
100       shrinked->GetCellData()->GetArray(usg_offsets->GetName()));
101
102   // now copy ELNO data. Start by verifying if it is possible to
103   // shallow copy the array.
104   vtkInformation *info=offsets->GetInformation();
105   vtkInformationQuadratureSchemeDefinitionVectorKey *key=
106       vtkQuadratureSchemeDefinition::DICTIONARY();
107   if(!key->Has(info))
108     {
109     vtkDebugMacro("Dictionary is not present in array " << offsets->GetName()
110                   << " " << offsets << " Aborting." );
111     return 0;
112     }
113   int dictSize=key->Size(info);
114   vtkQuadratureSchemeDefinition **dict=
115       new vtkQuadratureSchemeDefinition *[dictSize];
116   key->GetRange(info, dict, 0, 0, dictSize);
117
118   vtkIdType ncell=shrinked->GetNumberOfCells();
119
120   vtkFieldData* fielddata=usgIn->GetFieldData();
121   vtkIdList *ids=vtkIdList::New();
122   vtkIdList *surfaceIds=vtkIdList::New();
123   vtkIdList *originalIds=vtkIdList::New();
124   for(int index=0; index<fielddata->GetNumberOfArrays(); index++)
125     {
126     vtkDataArray* data=fielddata->GetArray(index);
127     if(data==NULL)
128       continue;
129
130     vtkInformation* info=data->GetInformation();
131     const char* arrayOffsetName=info->Get(
132         vtkQuadratureSchemeDefinition::QUADRATURE_OFFSET_ARRAY_NAME());
133
134     if(arrayOffsetName == NULL ||
135        strcmp(arrayOffsetName, offsets->GetName())!=0)
136       {
137       usgOut->GetFieldData()->AddArray(data);
138
139       continue;
140       }
141
142     vtkDataArray* newArray=data->NewInstance();
143     newArray->SetName(data->GetName());
144     usgOut->GetPointData()->AddArray(newArray);
145     newArray->SetNumberOfComponents(data->GetNumberOfComponents());
146     newArray->SetNumberOfTuples(usgOut->GetNumberOfPoints());
147     newArray->CopyComponentNames(data);
148     newArray->Delete();
149
150     for(vtkIdType cellId=0; cellId<ncell; cellId++)
151       {
152       vtkIdType offset=offsets->GetValue(cellId);
153
154       vtkIdType originalCellId=originalCellIds->GetValue(cellId);
155       int originalCellType=usgIn->GetCellType(originalCellId);
156
157       shrinked->GetCellPoints(cellId, ids);
158       surface->GetCellPoints(cellId, surfaceIds);
159
160       for(int id=0; id<ids->GetNumberOfIds(); id++)
161         {
162         vtkIdType surfaceId=surfaceIds->GetId(id);
163         vtkIdType shrinkedId=ids->GetId(id);
164         vtkIdType originalPointId = originalPointIds->GetValue(surfaceId);
165
166         usgIn->GetCellPoints(originalCellId, originalIds);
167         int originalLocalId=-1;
168         for(int li=0; li<originalIds->GetNumberOfIds(); li++)
169           {
170           if(originalPointId==originalIds->GetId(li))
171             {
172             originalLocalId=li;
173             break;
174             }
175           }
176         if(originalLocalId==-1)
177           {
178           originalLocalId=0;
179           vtkErrorMacro("cannot find original id");
180           }
181
182         const double * w=dict[originalCellType]->GetShapeFunctionWeights(
183             originalLocalId);
184         int j;
185         for(j=0; j<dict[originalCellType]->GetNumberOfNodes(); j++)
186           {
187           if(w[j]==1.0)
188             break;
189           }
190         if(j==dict[originalCellType]->GetNumberOfNodes())
191           {
192           vtkErrorMacro("cannot find elno weigth.");
193           j=id;
194           }
195         newArray->SetTuple(shrinkedId, offset+j, data);
196         }
197       }
198     }
199
200   ids->FastDelete();
201   surfaceIds->FastDelete();
202   originalIds->FastDelete();
203   delete[] dict;
204
205   return 1;
206 }
207
208 void vtkELNOSurfaceFilter::PrintSelf(ostream& os, vtkIndent indent)
209 {
210   this->Superclass::PrintSelf(os, indent);
211 }