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