Salome HOME
Merge from V6_main 01/04/2013
[modules/paravis.git] / src / Plugins / MedReader / IO / vtkExtractGroup.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 "vtkExtractGroup.h"
21
22 #include "vtkObjectFactory.h"
23 #include "vtkMutableDirectedGraph.h"
24 #include "vtkMultiBlockDataSet.h"
25 #include "vtkInformationVector.h"
26 #include "vtkInformation.h"
27 #include "vtkDataArraySelection.h"
28 #include "vtkMedUtilities.h"
29 #include "vtkTimeStamp.h"
30 #include "vtkInEdgeIterator.h"
31 #include "vtkMedReader.h"
32 #include "vtkInformationDataObjectKey.h"
33 #include "vtkExecutive.h"
34 #include "vtkVariantArray.h"
35 #include "vtkStringArray.h"
36 #include "vtkUnsignedCharArray.h"
37 #include "vtkDataSetAttributes.h"
38 #include "vtkDemandDrivenPipeline.h"
39 //#include "vtkCompositeDataIterator.h"
40 #include "vtkDataObjectTreeIterator.h"
41
42 #include <map>
43 #include <deque>
44 //#include <vtkstd/map>
45 //#include <vtkstd/deque>
46
47 // vtkCxxRevisionMacro(vtkExtractGroup, "$Revision$");
48 vtkStandardNewMacro(vtkExtractGroup);
49
50 vtkCxxSetObjectMacro(vtkExtractGroup, SIL, vtkMutableDirectedGraph);
51
52 vtkExtractGroup::vtkExtractGroup()
53 {
54   this->SIL=NULL;
55   this->Entities=vtkDataArraySelection::New();
56   this->Families=vtkDataArraySelection::New();
57   this->Groups=vtkDataArraySelection::New();
58   this->PruneOutput=0;
59 }
60
61 vtkExtractGroup::~vtkExtractGroup()
62 {
63   this->Entities->Delete();
64   this->Families->Delete();
65   this->Groups->Delete();
66 }
67
68 int vtkExtractGroup::ModifyRequest(vtkInformation* request, int when)
69 {
70   request->Set(vtkDemandDrivenPipeline::REQUEST_REGENERATE_INFORMATION(), 1);
71   return this->Superclass::ModifyRequest(request, when);
72 }
73
74 int vtkExtractGroup::RequestInformation(vtkInformation *request,
75     vtkInformationVector **inputVector, vtkInformationVector *outputVector)
76 {
77   vtkInformation* outInfo=outputVector->GetInformationObject(0);
78
79   vtkInformation* inputInfo=inputVector[0]->GetInformationObject(0);
80
81   vtkMutableDirectedGraph* old_SIL=this->GetSIL();
82
83   if(inputInfo->Has(vtkDataObject::SIL()))
84     {
85     this->SetSIL(vtkMutableDirectedGraph::SafeDownCast(inputInfo->Get(
86         vtkDataObject::SIL())));
87     }
88   else
89     {
90     vtkMutableDirectedGraph* sil=vtkMutableDirectedGraph::New();
91     this->BuildDefaultSIL(sil);
92     this->SetSIL(sil);
93     sil->Delete();
94     }
95
96   if(this->GetSIL()!=old_SIL||this->GetSIL()->GetMTime()>this->SILTime)
97     {
98     this->ClearSelections();
99     this->SILTime.Modified();
100     outInfo->Set(vtkDataObject::SIL(), this->GetSIL());
101     }
102
103   return 1;
104 }
105
106 vtkIdType vtkExtractGroup::FindVertex(const char* name)
107 {
108   vtkStringArray* names=vtkStringArray::SafeDownCast(
109       this->GetSIL()->GetVertexData()->GetAbstractArray("Names"));
110
111   return names->LookupValue(name);
112 }
113
114 void vtkExtractGroup::ClearSelections()
115 {
116   this->Families->RemoveAllArrays();
117   this->Entities->RemoveAllArrays();
118   this->Groups->RemoveAllArrays();
119 }
120
121 void vtkExtractGroup::BuildDefaultSIL(vtkMutableDirectedGraph* sil)
122 {
123   sil->Initialize();
124
125   vtkSmartPointer<vtkVariantArray> childEdge=
126       vtkSmartPointer<vtkVariantArray>::New();
127   childEdge->InsertNextValue(0);
128
129   vtkSmartPointer<vtkVariantArray> crossEdge=
130       vtkSmartPointer<vtkVariantArray>::New();
131   crossEdge->InsertNextValue(1);
132
133   // CrossEdge is an edge linking hierarchies.
134   vtkUnsignedCharArray* crossEdgesArray=vtkUnsignedCharArray::New();
135   crossEdgesArray->SetName("CrossEdges");
136   sil->GetEdgeData()->AddArray(crossEdgesArray);
137   crossEdgesArray->Delete();
138   std::deque<std::string> names;
139
140   // Now build the hierarchy.
141   vtkIdType rootId=sil->AddVertex();
142   names.push_back("SIL");
143
144   // Add a global entry to encode global names for the families
145   vtkIdType globalFamilyRoot=sil->AddChild(rootId, childEdge);
146   names.push_back("Families");
147
148   // Add a global entry to encode global names for the families
149   vtkIdType globalGroupRoot=sil->AddChild(rootId, childEdge);
150   names.push_back("Groups");
151
152   // Add the groups subtree
153   vtkIdType groupsRoot=sil->AddChild(rootId, childEdge);
154   names.push_back("GroupTree");
155
156   // Add the attributes subtree
157   vtkIdType attributesRoot=sil->AddChild(rootId, childEdge);
158   names.push_back("Attributes");
159
160   // Add a global entry to encode names for the cell types
161   vtkIdType globalEntityRoot=sil->AddChild(rootId, childEdge);
162   names.push_back("Entity");
163
164   // Add the cell types subtree
165   vtkIdType entityTypesRoot=sil->AddChild(rootId, childEdge);
166   names.push_back("EntityTree");
167
168   // This array is used to assign names to nodes.
169   vtkStringArray* namesArray=vtkStringArray::New();
170   namesArray->SetName("Names");
171   namesArray->SetNumberOfTuples(sil->GetNumberOfVertices());
172   sil->GetVertexData()->AddArray(namesArray);
173   namesArray->Delete();
174   std::deque<std::string>::iterator iter;
175   vtkIdType cc;
176   for(cc=0, iter=names.begin(); iter!=names.end(); ++iter, ++cc)
177     {
178     namesArray->SetValue(cc, (*iter).c_str());
179     }
180
181 }
182
183 int vtkExtractGroup::RequestData(vtkInformation *request,
184     vtkInformationVector **inputVector, vtkInformationVector *outputVector)
185
186 {
187   vtkInformation* inputInfo=inputVector[0]->GetInformationObject(0);
188   vtkMultiBlockDataSet* inmb=vtkMultiBlockDataSet::SafeDownCast(inputInfo->Get(
189       vtkDataObject::DATA_OBJECT()));
190
191   if(inmb==NULL)
192     return 0;
193
194   vtkMultiBlockDataSet* outmb=this->GetOutput();
195
196   outmb->CopyStructure(inmb);
197
198   //vtkCompositeDataIterator* iterator = inmb->NewIterator();
199   vtkDataObjectTreeIterator* iterator = vtkDataObjectTreeIterator::SafeDownCast( inmb->NewIterator() );
200   iterator->SetVisitOnlyLeaves(true);
201   iterator->InitTraversal();
202   while(!iterator->IsDoneWithTraversal())
203     {
204     vtkDataObject* indo = iterator->GetCurrentDataObject();
205     if(indo == NULL)
206       continue;
207
208     if(indo->GetFieldData()->HasArray("BLOCK_NAME"))
209       {
210
211       vtkStringArray* path = vtkStringArray::SafeDownCast(
212           indo->GetFieldData()->GetAbstractArray("BLOCK_NAME"));
213
214       if(this->IsBlockSelected(path))
215         {
216         vtkMultiBlockDataSet* parent = vtkMedUtilities::GetParent(outmb, path);
217         int nb = parent->GetNumberOfBlocks();
218         parent->SetNumberOfBlocks(nb+1);
219         vtkDataObject *outdo = indo->NewInstance();
220         outdo->ShallowCopy(indo);
221         parent->SetBlock(nb, outdo);
222         outdo->Delete();
223         }
224       }
225     iterator->GoToNextItem();
226     }
227
228   if(PruneOutput)
229     {
230     this->PruneEmptyBlocks(outmb);
231     }
232   return 1;
233 }
234
235 void vtkExtractGroup::SetGroupStatus(const char* key, int flag)
236 {
237   vtkIdType index=this->Groups->GetArrayIndex(key);
238   if(index==-1)
239     {
240     index = this->Groups->AddArray(key);
241     this->Modified();
242     }
243   int status=this->Groups->GetArraySetting(index);
244   if(status!=flag)
245     {
246     if(flag)
247       {
248       this->Groups->EnableArray(key);
249       }
250     else
251       {
252       this->Groups->DisableArray(key);
253       }
254     this->Modified();
255     }
256   this->GroupSelectionTime.Modified();
257 }
258
259 void vtkExtractGroup::PruneEmptyBlocks(vtkMultiBlockDataSet* mb)
260 {
261   if(mb==NULL)
262     return;
263   vtkIdType nn=0;
264   while(nn<mb->GetNumberOfBlocks())
265     {
266     bool remove=false;
267     vtkDataObject* dataObj=mb->GetBlock(nn);
268     if(dataObj==NULL)
269       {
270       remove=true;
271       }
272     else
273       {
274       vtkMultiBlockDataSet* child=vtkMultiBlockDataSet::SafeDownCast(dataObj);
275       if(child!=NULL)
276         {
277         this->PruneEmptyBlocks(child);
278         if(child->GetNumberOfBlocks()==0)
279           {
280           remove=true;
281           }
282         }
283       }
284     if(remove)
285       {
286       mb->RemoveBlock(nn);
287       }
288     else
289       {
290       nn++;
291       }
292     }
293 }
294
295 int vtkExtractGroup::IsBlockSelected(vtkStringArray* path)
296 {
297   const char* meshName = (path->GetNumberOfValues()>0?
298                           path->GetValue(0) : NULL);
299   const char* cellOrPoint = (path->GetNumberOfValues()>1?
300                              path->GetValue(1) : NULL);
301   const char* familyName = (path->GetNumberOfValues()>2?
302                             path->GetValue(2) : NULL);
303
304   if(!this->IsFamilySelected(meshName, cellOrPoint, familyName))
305     {
306     return 0;
307     }
308
309   bool isOnPoint = (strcmp(cellOrPoint, vtkMedUtilities::OnPointName)==0);
310
311   const char* entityName = (isOnPoint || path->GetNumberOfValues()<=3 ? NULL :
312                             path->GetValue(3));
313
314   if(isOnPoint)
315     return true;
316
317   return IsEntitySelected(entityName);
318
319 }
320
321 int vtkExtractGroup::IsEntitySelected(const char* entityKey)
322 {
323   return this->Entities->GetArraySetting(entityKey);
324 }
325
326 int vtkExtractGroup::IsFamilySelected(const char* meshName,
327     const char* pointOrCellKey, const char* familyName)
328 {
329   if(this->FamilySelectionTime <= this->GroupSelectionTime)
330     {
331     this->SelectFamiliesFromGroups();
332     }
333
334   int
335       pointOrCell= (strcmp(vtkMedUtilities::OnPointName, pointOrCellKey)==0?
336                     vtkMedUtilities::OnPoint
337                     : vtkMedUtilities::OnCell);
338
339   std::string name=
340       vtkMedUtilities::FamilyKey(meshName, pointOrCell, familyName);
341
342   return this->Families->GetArraySetting(name.c_str());
343 }
344
345 void vtkExtractGroup::SelectFamiliesFromGroups()
346 {
347   this->Families->DisableAllArrays();
348   vtkStringArray* names=vtkStringArray::SafeDownCast(
349       this->GetSIL()->GetVertexData()->GetAbstractArray("Names"));
350
351   for(int index = 0; index < this->Groups->GetNumberOfArrays(); index++)
352     {
353     if(this->Groups->GetArraySetting(index) == 0)
354       continue;
355
356     const char* name = this->Groups->GetArrayName(index);
357     vtkIdType silindex = this->FindVertex(name);
358
359     vtkInEdgeIterator* it = vtkInEdgeIterator::New();
360
361     this->GetSIL()->GetInEdges(silindex, it);
362     while(it->HasNext())
363       {
364       vtkIdType famId = it->Next().Source;
365       vtkStdString famName = names->GetValue(famId);
366       if(strncmp(famName, "FAMILY", 6)==0)
367         {
368         this->Families->EnableArray(famName.c_str());
369         }
370       }
371     it->Delete();
372     }
373
374   this->FamilySelectionTime.Modified();
375 }
376
377 void vtkExtractGroup::SetEntityStatus(const char* key, int flag)
378 {
379   vtkIdType index=this->Entities->GetArrayIndex(key);
380   if(index==-1)
381     {
382     index = this->Entities->AddArray(key);
383     this->Modified();
384     }
385   int status=this->Entities->GetArraySetting(index);
386   if(status!=flag)
387     {
388     if(flag)
389       {
390       this->Entities->EnableArray(key);
391       }
392     else
393       {
394       this->Entities->DisableArray(key);
395       }
396     this->Modified();
397     }
398 }
399
400 void vtkExtractGroup::SetFamilyStatus(const char* key, int flag)
401 {
402   vtkIdType index=this->Families->GetArrayIndex(key);
403   if(index==-1)
404     {
405     return;
406     }
407   int status=this->Families->GetArraySetting(index);
408   if(status!=flag)
409     {
410     if(flag)
411       {
412       this->Families->EnableArray(key);
413       }
414     else
415       {
416       this->Families->DisableArray(key);
417       }
418     }
419 }
420
421 int vtkExtractGroup::GetSILUpdateStamp()
422 {
423   return this->SILTime;
424 }
425
426 void vtkExtractGroup::PrintSelf(ostream& os, vtkIndent indent)
427 {
428   this->Superclass::PrintSelf(os, indent);
429 }