Salome HOME
2b2af12ddb0ff9e046037d3303ae78451e783060
[modules/paravis.git] / src / Plugins / MEDReader / IO / vtkExtractGroup.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 // Author : Anthony Geay
20
21 #include "vtkExtractGroup.h"
22 #include "MEDFileFieldRepresentationTree.hxx"
23 #include "vtkMEDReader.h"
24 #include "VTKMEDTraits.hxx"
25
26 #include "vtkAdjacentVertexIterator.h"
27 #include "vtkDataArrayTemplate.h"
28 #include "vtkIntArray.h"
29 #include "vtkLongArray.h"
30 #include "vtkCellData.h"
31 #include "vtkPointData.h"
32
33 #include "vtkStreamingDemandDrivenPipeline.h"
34 #include "vtkUnstructuredGrid.h"
35 #include  "vtkMultiBlockDataSet.h"
36
37 #include "vtkInformationStringKey.h"
38 #include "vtkAlgorithmOutput.h"
39 #include "vtkObjectFactory.h"
40 #include "vtkMutableDirectedGraph.h"
41 #include "vtkMultiBlockDataSet.h"
42 #include "vtkDataSet.h"
43 #include "vtkInformationVector.h"
44 #include "vtkInformation.h"
45 #include "vtkDataArraySelection.h"
46 #include "vtkTimeStamp.h"
47 #include "vtkInEdgeIterator.h"
48 #include "vtkInformationDataObjectKey.h"
49 #include "vtkExecutive.h"
50 #include "vtkVariantArray.h"
51 #include "vtkStringArray.h"
52 #include "vtkDoubleArray.h"
53 #include "vtkCharArray.h"
54 #include "vtkUnsignedCharArray.h"
55 #include "vtkDataSetAttributes.h"
56 #include "vtkDemandDrivenPipeline.h"
57 #include "vtkDataObjectTreeIterator.h"
58 #include "vtkThreshold.h"
59 #include "vtkMultiBlockDataGroupFilter.h"
60 #include "vtkCompositeDataToUnstructuredGridFilter.h"
61 #include "vtkInformationDataObjectMetaDataKey.h"
62
63 #include <map>
64 #include <deque>
65
66 vtkStandardNewMacro(vtkExtractGroup);
67
68 ///////////////////
69
70 class ExtractGroupStatus
71 {
72 public:
73   ExtractGroupStatus():_status(false) { }
74   ExtractGroupStatus(const char *name);
75   bool getStatus() const { return _status; }
76   void setStatus(bool status) const { _status=status; }
77   void cpyStatusFrom(const ExtractGroupStatus& other) { _status=other._status; }
78   std::string getName() const { return _name; }
79   void resetStatus() const { _status=false; }
80   const char *getKeyOfEntry() const { return _ze_key_name.c_str(); }
81   virtual void printMySelf(std::ostream& os) const;
82   virtual bool isSameAs(const ExtractGroupStatus& other) const;
83 protected:
84 mutable bool _status;
85 std::string _name;
86 std::string _ze_key_name;
87 };
88
89 class ExtractGroupGrp : public ExtractGroupStatus
90 {
91 public:
92   ExtractGroupGrp(const char *name):ExtractGroupStatus(name) { std::ostringstream oss; oss << START << name; _ze_key_name=oss.str(); }
93   void setFamilies(const std::vector<std::string>& fams) { _fams=fams; }
94   const std::vector<std::string>& getFamiliesLyingOn() const { return _fams; }
95   bool isSameAs(const ExtractGroupGrp& other) const;
96 public:
97   static const char START[];
98   std::vector<std::string> _fams;
99 };
100
101 class ExtractGroupFam : public ExtractGroupStatus
102 {
103 public:
104   ExtractGroupFam(const char *name);
105   void printMySelf(std::ostream& os) const;
106   void fillIdsToKeep(std::set<int>& s) const;
107   int getId() const { return _id; }
108   bool isSameAs(const ExtractGroupFam& other) const;
109 public:
110   static const char START[];
111 private:
112   int _id;
113 };
114
115 class vtkExtractGroup::vtkExtractGroupInternal
116 {
117 public:
118   void loadFrom(vtkMutableDirectedGraph *sil);
119   int getNumberOfEntries() const;
120   const char *getMeshName() const;
121   const char *getKeyOfEntry(int i) const;
122   bool getStatusOfEntryStr(const char *entry) const;
123   void setStatusOfEntryStr(const char *entry, bool status);
124   void printMySelf(std::ostream& os) const;
125   std::set<int> getIdsToKeep() const;
126   void clearSelection() const;
127   int getIdOfFamily(const std::string& famName) const;
128   static bool IsInformationOK(vtkInformation *info);
129 private:
130   std::map<std::string,int> computeFamStrIdMap() const;
131   const ExtractGroupStatus& getEntry(const char *entry) const;
132   ExtractGroupStatus& getEntry(const char *entry);
133 private:
134   std::vector<ExtractGroupGrp> _groups;
135   std::vector<ExtractGroupFam> _fams;
136   mutable std::vector< std::pair<std::string,bool> > _selection;
137   std::string _mesh_name;
138 };
139
140 const char ExtractGroupGrp::START[]="GRP_";
141
142 const char ExtractGroupFam::START[]="FAM_";
143
144 ExtractGroupStatus::ExtractGroupStatus(const char *name):_status(false),_name(name)
145 {
146 }
147
148 void ExtractGroupStatus::printMySelf(std::ostream& os) const
149 {
150   os << "      -" << _ze_key_name << "(";
151   if(_status)
152     os << "X";
153   else
154     os << " ";
155   os << ")" << std::endl;
156 }
157
158 bool ExtractGroupStatus::isSameAs(const ExtractGroupStatus& other) const
159 {
160   return _name==other._name && _ze_key_name==other._ze_key_name;
161 }
162
163 bool ExtractGroupGrp::isSameAs(const ExtractGroupGrp& other) const
164 {
165   bool ret(ExtractGroupStatus::isSameAs(other));
166   if(ret)
167     return _fams==other._fams;
168   else
169     return false;
170 }
171
172 bool vtkExtractGroup::vtkExtractGroupInternal::IsInformationOK(vtkInformation *info)
173 {
174   // Check the information contain meta data key
175   if(!info->Has(vtkMEDReader::META_DATA()))
176     return false;
177
178   // Recover Meta Data
179   vtkMutableDirectedGraph *sil(vtkMutableDirectedGraph::SafeDownCast(info->Get(vtkMEDReader::META_DATA())));
180   if(!sil)
181     return false;
182   int idNames(0);
183   vtkAbstractArray *verticesNames(sil->GetVertexData()->GetAbstractArray("Names",idNames));
184   vtkStringArray *verticesNames2(vtkStringArray::SafeDownCast(verticesNames));
185   if(!verticesNames2)
186     return false;
187   for(int i=0;i<verticesNames2->GetNumberOfValues();i++)
188     {
189       vtkStdString &st(verticesNames2->GetValue(i));
190       if(st=="MeshesFamsGrps")
191         return true;
192     }
193   return false;
194 }
195
196 const char* vtkExtractGroup::GetGrpStart()
197 {
198   return ExtractGroupGrp::START;
199 }
200
201 const char* vtkExtractGroup::GetFamStart()
202 {
203   return ExtractGroupFam::START;
204 }
205
206 const char *vtkExtractGroup::vtkExtractGroupInternal::getMeshName() const
207 {
208   return this->_mesh_name.c_str();
209 }
210
211 void vtkExtractGroup::vtkExtractGroupInternal::loadFrom(vtkMutableDirectedGraph *sil)
212 {
213   std::vector<ExtractGroupGrp> oldGrps(_groups); _groups.clear();
214   std::vector<ExtractGroupFam> oldFams(_fams); _fams.clear();
215   int idNames(0);
216   vtkAbstractArray *verticesNames(sil->GetVertexData()->GetAbstractArray("Names",idNames));
217   vtkStringArray *verticesNames2(vtkStringArray::SafeDownCast(verticesNames));
218   vtkIdType id0;
219   bool found(false);
220   for(int i=0;i<verticesNames2->GetNumberOfValues();i++)
221     {
222       vtkStdString &st(verticesNames2->GetValue(i));
223       if(st=="MeshesFamsGrps")
224         {
225           id0=i;
226           found=true;
227         }
228     }
229   if(!found)
230     throw INTERP_KERNEL::Exception("There is an internal error ! The tree on server side has not the expected look !");
231   vtkAdjacentVertexIterator *it0(vtkAdjacentVertexIterator::New());
232   sil->GetAdjacentVertices(id0,it0);
233   int kk(0),ll(0);
234   while(it0->HasNext())
235     {
236       vtkIdType id1(it0->Next());
237       std::string meshName((const char *)verticesNames2->GetValue(id1));
238       this->_mesh_name=meshName;
239       vtkAdjacentVertexIterator *it1(vtkAdjacentVertexIterator::New());
240       sil->GetAdjacentVertices(id1,it1);
241       vtkIdType idZeGrps(it1->Next());//zeGroups
242       vtkAdjacentVertexIterator *itGrps(vtkAdjacentVertexIterator::New());
243       sil->GetAdjacentVertices(idZeGrps,itGrps);
244       while(itGrps->HasNext())
245         {
246           vtkIdType idg(itGrps->Next());
247           ExtractGroupGrp grp((const char *)verticesNames2->GetValue(idg));
248           vtkAdjacentVertexIterator *itGrps2(vtkAdjacentVertexIterator::New());
249           sil->GetAdjacentVertices(idg,itGrps2);
250           std::vector<std::string> famsOnGroup;
251           while(itGrps2->HasNext())
252             {
253               vtkIdType idgf(itGrps2->Next());
254               famsOnGroup.push_back(std::string((const char *)verticesNames2->GetValue(idgf)));
255             }
256           grp.setFamilies(famsOnGroup);
257           itGrps2->Delete();
258           _groups.push_back(grp);
259         }
260       itGrps->Delete();
261       vtkIdType idZeFams(it1->Next());//zeFams
262       it1->Delete();
263       vtkAdjacentVertexIterator *itFams(vtkAdjacentVertexIterator::New());
264       sil->GetAdjacentVertices(idZeFams,itFams);
265       while(itFams->HasNext())
266         {
267           vtkIdType idf(itFams->Next());
268           ExtractGroupFam fam((const char *)verticesNames2->GetValue(idf));
269           _fams.push_back(fam);
270         }
271       itFams->Delete();
272     }
273   it0->Delete();
274   //
275   std::size_t szg(_groups.size()),szf(_fams.size());
276   if(szg==oldGrps.size() && szf==oldFams.size())
277     {
278       bool isSame(true);
279       for(std::size_t i=0;i<szg && isSame;i++)
280         isSame=_groups[i].isSameAs(oldGrps[i]);
281       for(std::size_t i=0;i<szf && isSame;i++)
282         isSame=_fams[i].isSameAs(oldFams[i]);
283       if(isSame)
284         {
285           for(std::size_t i=0;i<szg;i++)
286             _groups[i].cpyStatusFrom(oldGrps[i]);
287           for(std::size_t i=0;i<szf;i++)
288             _fams[i].cpyStatusFrom(oldFams[i]);
289         }
290     }
291 }
292
293 int vtkExtractGroup::vtkExtractGroupInternal::getNumberOfEntries() const
294 {
295   std::size_t sz0(_groups.size()),sz1(_fams.size());
296   return (int)(sz0+sz1);
297 }
298
299 const char *vtkExtractGroup::vtkExtractGroupInternal::getKeyOfEntry(int i) const
300 {
301   int sz0((int)_groups.size());
302   if(i>=0 && i<sz0)
303     return _groups[i].getKeyOfEntry();
304   else
305     return _fams[i-sz0].getKeyOfEntry();
306 }
307
308 bool vtkExtractGroup::vtkExtractGroupInternal::getStatusOfEntryStr(const char *entry) const
309 {
310   const ExtractGroupStatus& elt(getEntry(entry));
311   return elt.getStatus();
312 }
313
314 void vtkExtractGroup::vtkExtractGroupInternal::setStatusOfEntryStr(const char *entry, bool status)
315 {
316   _selection.emplace_back(entry,status);
317 }
318
319 const ExtractGroupStatus& vtkExtractGroup::vtkExtractGroupInternal::getEntry(const char *entry) const
320 {
321   std::string entryCpp(entry);
322   for(std::vector<ExtractGroupGrp>::const_iterator it0=_groups.begin();it0!=_groups.end();it0++)
323     if(entryCpp==(*it0).getKeyOfEntry())
324       return *it0;
325   for(std::vector<ExtractGroupFam>::const_iterator it0=_fams.begin();it0!=_fams.end();it0++)
326     if(entryCpp==(*it0).getKeyOfEntry())
327       return *it0;
328   std::ostringstream oss; oss << "vtkExtractGroupInternal::getEntry : no such entry \"" << entry << "\"!";
329   throw INTERP_KERNEL::Exception(oss.str().c_str());
330 }
331
332 ExtractGroupStatus& vtkExtractGroup::vtkExtractGroupInternal::getEntry(const char *entry)
333 {
334   std::string entryCpp(entry);
335   for(std::vector<ExtractGroupGrp>::iterator it0=_groups.begin();it0!=_groups.end();it0++)
336     if(entryCpp==(*it0).getKeyOfEntry())
337       return *it0;
338   for(std::vector<ExtractGroupFam>::iterator it0=_fams.begin();it0!=_fams.end();it0++)
339     if(entryCpp==(*it0).getKeyOfEntry())
340       return *it0;
341   std::ostringstream oss; oss << "vtkExtractGroupInternal::getEntry : no such entry \"" << entry << "\"!";
342   throw INTERP_KERNEL::Exception(oss.str().c_str());
343 }
344
345 void vtkExtractGroup::vtkExtractGroupInternal::printMySelf(std::ostream& os) const
346 {
347   os << "Groups :" << std::endl;
348   for(std::vector<ExtractGroupGrp>::const_iterator it0=_groups.begin();it0!=_groups.end();it0++)
349     (*it0).printMySelf(os);
350   os << "Families :" << std::endl;
351   for(std::vector<ExtractGroupFam>::const_iterator it0=_fams.begin();it0!=_fams.end();it0++)
352     (*it0).printMySelf(os);
353 }
354
355 int vtkExtractGroup::vtkExtractGroupInternal::getIdOfFamily(const std::string& famName) const
356 {
357   for(std::vector<ExtractGroupFam>::const_iterator it=_fams.begin();it!=_fams.end();it++)
358     {
359       if((*it).getName()==famName)
360         return (*it).getId();
361     }
362   return std::numeric_limits<int>::max();
363 }
364
365 ExtractGroupFam::ExtractGroupFam(const char *name):ExtractGroupStatus(name),_id(0)
366 {
367   std::size_t pos(_name.find(MEDFileFieldRepresentationLeavesArrays::ZE_SEP));
368   std::string name0(_name.substr(0,pos)),name1(_name.substr(pos+strlen(MEDFileFieldRepresentationLeavesArrays::ZE_SEP)));
369   std::istringstream iss(name1);
370   iss >> _id;
371   std::ostringstream oss; oss << START << name; _ze_key_name=oss.str(); _name=name0;
372 }
373
374 bool ExtractGroupFam::isSameAs(const ExtractGroupFam& other) const
375 {
376   bool ret(ExtractGroupStatus::isSameAs(other));
377   if(ret)
378     return _id==other._id;
379   else
380     return false;
381 }
382
383 void ExtractGroupFam::printMySelf(std::ostream& os) const
384 {
385   os << "      -" << _ze_key_name << " famName : \"" << _name << "\" id : " << _id << " (";
386   if(_status)
387     os << "X";
388   else
389     os << " ";
390   os << ")" << std::endl;
391 }
392
393 void ExtractGroupFam::fillIdsToKeep(std::set<int>& s) const
394 {
395   s.insert(_id);
396 }
397
398 std::set<int> vtkExtractGroup::vtkExtractGroupInternal::getIdsToKeep() const
399 {
400   for(auto it: _selection)
401     {
402       const ExtractGroupStatus& elt(getEntry(it.first.c_str()));
403       elt.setStatus(it.second);
404     }
405   std::map<std::string,int> m(this->computeFamStrIdMap());
406   std::set<int> s;
407   for(std::vector<ExtractGroupGrp>::const_iterator it0=_groups.begin();it0!=_groups.end();it0++)
408     {
409       if((*it0).getStatus())
410         {
411           const std::vector<std::string>& fams((*it0).getFamiliesLyingOn());
412           for(std::vector<std::string>::const_iterator it1=fams.begin();it1!=fams.end();it1++)
413             {
414               std::map<std::string,int>::iterator it2(m.find((*it1)));
415               if(it2!=m.end())
416                 s.insert((*it2).second);
417             }
418         }
419      }
420   for(std::vector<ExtractGroupFam>::const_iterator it0=_fams.begin();it0!=_fams.end();it0++)
421     if((*it0).getStatus())
422       (*it0).fillIdsToKeep(s);
423   return s;
424 }
425
426 void vtkExtractGroup::vtkExtractGroupInternal::clearSelection() const
427 {
428   _selection.clear();
429   for(auto it : _groups)
430     it.resetStatus();
431   for(auto it : _fams)
432     it.resetStatus();
433 }
434
435 std::map<std::string,int> vtkExtractGroup::vtkExtractGroupInternal::computeFamStrIdMap() const
436 {
437   std::map<std::string,int> ret;
438   for(std::vector<ExtractGroupFam>::const_iterator it0=_fams.begin();it0!=_fams.end();it0++)
439     ret[(*it0).getName()]=(*it0).getId();
440   return ret;
441 }
442
443 ////////////////////
444
445 vtkExtractGroup::vtkExtractGroup():SIL(NULL),Internal(new vtkExtractGroupInternal),InsideOut(0)
446 {
447 }
448
449 vtkExtractGroup::~vtkExtractGroup()
450 {
451   delete this->Internal;
452 }
453
454 void vtkExtractGroup::SetInsideOut(int val)
455 {
456   if(this->InsideOut!=val)
457     {
458       this->InsideOut=val;
459       this->Modified();
460     }
461 }
462
463 int vtkExtractGroup::RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
464 {
465 //  vtkUnstructuredGridAlgorithm::RequestInformation(request,inputVector,outputVector);
466   try
467     {
468 //      std::cerr << "########################################## vtkExtractGroup::RequestInformation ##########################################" << std::endl;
469 //      request->Print(cout);
470       vtkInformation *outInfo(outputVector->GetInformationObject(0));
471       vtkInformation *inputInfo(inputVector[0]->GetInformationObject(0));
472       if(!vtkExtractGroup::vtkExtractGroupInternal::IsInformationOK(inputInfo))
473         {
474         vtkErrorMacro("No SIL Data available ! The source of this filter must be MEDReader !");
475         return 0;
476         }
477
478       this->SetSIL(vtkMutableDirectedGraph::SafeDownCast(inputInfo->Get(vtkMEDReader::META_DATA())));
479       this->Internal->loadFrom(this->SIL);
480       //this->Internal->printMySelf(std::cerr);
481     }
482   catch(INTERP_KERNEL::Exception& e)
483     {
484       std::cerr << "Exception has been thrown in vtkExtractGroup::RequestInformation : " << e.what() << std::endl;
485       return 0;
486     }
487   return 1;
488 }
489
490 /*!
491  * Do not use vtkCxxSetObjectMacro macro because input mdg comes from an already managed in the pipeline just a ref on it.
492  */
493 void vtkExtractGroup::SetSIL(vtkMutableDirectedGraph *mdg)
494 {
495   if(this->SIL==mdg)
496     return ;
497   this->SIL=mdg;
498 }
499
500 template<class CellPointExtractor>
501 vtkDataSet *FilterFamilies(vtkSmartPointer<vtkThreshold>& thres,
502                            vtkDataSet *input, const std::set<int>& idsToKeep, bool insideOut, const char *arrNameOfFamilyField,
503                            const char *associationForThreshold, bool& catchAll, bool& catchSmth)
504 {
505   const int VTK_DATA_ARRAY_DELETE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_DELETE;
506   const char ZE_SELECTION_ARR_NAME[]="@@ZeSelection@@";
507   vtkDataSet *output(input->NewInstance());
508   output->ShallowCopy(input);
509   thres->SetInputData(output);
510   vtkDataSetAttributes *dscIn(input->GetCellData()),*dscIn2(input->GetPointData());
511   vtkDataSetAttributes *dscOut(output->GetCellData()),*dscOut2(output->GetPointData());
512   //
513   double vMin(insideOut==0?1.:0.),vMax(insideOut==0?2.:1.);
514   thres->ThresholdBetween(vMin,vMax);
515   // OK for the output
516   //
517   CellPointExtractor cpe2(input);
518   vtkDataArray *da(cpe2.Get()->GetScalars(arrNameOfFamilyField));
519   if(!da)
520     return 0;
521   std::string daName(da->GetName());
522   typedef MEDFileVTKTraits<mcIdType>::VtkType vtkMCIdTypeArray;
523   vtkMCIdTypeArray *dai(vtkMCIdTypeArray::SafeDownCast(da));
524   if(daName!=arrNameOfFamilyField || !dai)
525     return 0;
526   //
527   int nbOfTuples(dai->GetNumberOfTuples());
528   vtkCharArray *zeSelection(vtkCharArray::New());
529   zeSelection->SetName(ZE_SELECTION_ARR_NAME);
530   zeSelection->SetNumberOfComponents(1);
531   char *pt(new char[nbOfTuples]);
532   zeSelection->SetArray(pt,nbOfTuples,0,VTK_DATA_ARRAY_DELETE);
533   const mcIdType *inPtr(dai->GetPointer(0));
534   std::fill(pt,pt+nbOfTuples,0);
535   catchAll=true; catchSmth=false;
536   std::vector<bool> pt2(nbOfTuples,false);
537   for(std::set<int>::const_iterator it=idsToKeep.begin();it!=idsToKeep.end();it++)
538     {
539       bool catchFid(false);
540       for(int i=0;i<nbOfTuples;i++)
541         if(inPtr[i]==*it)
542           { pt2[i]=true; catchFid=true; }
543       if(!catchFid)
544         catchAll=false;
545       else
546         catchSmth=true;
547     }
548   for(int ii=0;ii<nbOfTuples;ii++)
549     if(pt2[ii])
550       pt[ii]=2;
551   CellPointExtractor cpe3(output);
552   int idx(cpe3.Get()->AddArray(zeSelection));
553   cpe3.Get()->SetActiveAttribute(idx,vtkDataSetAttributes::SCALARS);
554   cpe3.Get()->CopyScalarsOff();
555   zeSelection->Delete();
556   //
557   thres->SetInputArrayToProcess(idx,0,0,associationForThreshold,ZE_SELECTION_ARR_NAME);
558   thres->Update();
559   vtkUnstructuredGrid *zeComputedOutput(thres->GetOutput());
560   CellPointExtractor cpe(zeComputedOutput);
561   cpe.Get()->RemoveArray(idx);
562   output->Delete();
563   zeComputedOutput->Register(0);
564   return zeComputedOutput;
565 }
566
567 class CellExtractor
568 {
569 public:
570   CellExtractor(vtkDataSet *ds):_ds(ds) { }
571   vtkDataSetAttributes *Get() { return _ds->GetCellData(); }
572 private:
573   vtkDataSet *_ds;
574 };
575
576 class PointExtractor
577 {
578 public:
579   PointExtractor(vtkDataSet *ds):_ds(ds) { }
580   vtkDataSetAttributes *Get() { return _ds->GetPointData(); }
581 private:
582   vtkDataSet *_ds;
583 };
584 int vtkExtractGroup::RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
585 {
586   try
587     {
588       // std::cerr << "########################################## vtkExtractGroup::RequestData        ##########################################" << std::endl;
589       // request->Print(cout);
590       vtkInformation* inputInfo=inputVector[0]->GetInformationObject(0);
591       vtkMultiBlockDataSet *inputMB(vtkMultiBlockDataSet::SafeDownCast(inputInfo->Get(vtkDataObject::DATA_OBJECT())));
592       if(inputMB->GetNumberOfBlocks()!=1)
593         {
594           std::ostringstream oss; oss << "vtkExtractGroup::RequestData : input has not the right number of parts ! Expected 1 !";
595           if(this->HasObserver("ErrorEvent") )
596             this->InvokeEvent("ErrorEvent",const_cast<char *>(oss.str().c_str()));
597           else
598             vtkOutputWindowDisplayErrorText(const_cast<char *>(oss.str().c_str()));
599           vtkObject::BreakOnError();
600           return 0;
601         }
602       vtkDataSet *input(vtkDataSet::SafeDownCast(inputMB->GetBlock(0)));
603       vtkInformation *info(input->GetInformation());
604       vtkInformation *outInfo(outputVector->GetInformationObject(0));
605       vtkMultiBlockDataSet *output(vtkMultiBlockDataSet::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT())));
606       std::set<int> idsToKeep(this->Internal->getIdsToKeep());
607       this->Internal->clearSelection();
608       // first shrink the input
609       bool catchAll,catchSmth;
610       vtkSmartPointer<vtkThreshold> thres1(vtkSmartPointer<vtkThreshold>::New()),thres2(vtkSmartPointer<vtkThreshold>::New());
611       vtkDataSet *tryOnCell(FilterFamilies<CellExtractor>(thres1,input,idsToKeep,this->InsideOut,
612                                                           MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_CELL_NAME,"vtkDataObject::FIELD_ASSOCIATION_CELLS",catchAll,catchSmth));
613       if(tryOnCell)
614         {
615           if(catchAll)
616             {
617               output->SetBlock(0,tryOnCell);
618               tryOnCell->Delete();//
619               return 1;
620             }
621           else
622             {
623               if(catchSmth)
624                 {
625                   vtkDataSet *tryOnNode(FilterFamilies<PointExtractor>(thres2,input,idsToKeep,this->InsideOut,
626                                                                        MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME,"vtkDataObject::FIELD_ASSOCIATION_POINTS",catchAll,catchSmth));
627                   if(tryOnNode && catchSmth)
628                     {
629                       output->SetBlock(0,tryOnCell);
630                       output->SetBlock(1,tryOnNode);
631                       tryOnCell->Delete();
632                       tryOnNode->Delete();
633                       return 1;
634                     }
635                   else
636                     {
637                       if(tryOnNode)
638                         tryOnNode->Delete();
639                       output->SetBlock(0,tryOnCell);
640                       tryOnCell->Delete();
641                       return 1;
642                     }
643                 }
644               else
645                 {
646                   vtkDataSet *tryOnNode(FilterFamilies<PointExtractor>(thres1,input,idsToKeep,this->InsideOut,
647                                                                        MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME,"vtkDataObject::FIELD_ASSOCIATION_POINTS",catchAll,catchSmth));
648                   if(tryOnNode)
649                     {
650                       tryOnCell->Delete();
651                       output->SetBlock(0,tryOnNode);
652                       tryOnNode->Delete();
653                       return 1;
654                     }
655                   else
656                     {
657                       output->SetBlock(0,tryOnNode);
658                       tryOnCell->Delete();
659                       return 0;
660                     }
661                 }
662             }
663         }
664       else
665         {
666           vtkDataSet *tryOnNode(FilterFamilies<PointExtractor>(thres1,input,idsToKeep,this->InsideOut,
667                                                                MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME,"vtkDataObject::FIELD_ASSOCIATION_POINTS",catchAll,catchSmth));
668           if(tryOnNode)
669             {
670               output->ShallowCopy(tryOnNode);
671               tryOnNode->Delete();//
672               return 1;
673             }
674           else
675             {
676               std::ostringstream oss; oss << "vtkExtractGroup::RequestData : The integer array with name \""<< MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_CELL_NAME;
677               oss << "\" or \"" << MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME << "\" does not exist ! The extraction of group and/or family is not possible !";
678               if(this->HasObserver("ErrorEvent") )
679                 this->InvokeEvent("ErrorEvent",const_cast<char *>(oss.str().c_str()));
680               else
681                 vtkOutputWindowDisplayErrorText(const_cast<char *>(oss.str().c_str()));
682               vtkObject::BreakOnError();
683               return 0;
684             }
685         }
686     }
687   catch(INTERP_KERNEL::Exception& e)
688     {
689       std::cerr << "Exception has been thrown in vtkExtractGroup::RequestData : " << e.what() << std::endl;
690       return 0;
691     }
692 }
693
694 int vtkExtractGroup::GetSILUpdateStamp()
695 {
696   return (int)this->SILTime;
697 }
698
699 void vtkExtractGroup::PrintSelf(ostream& os, vtkIndent indent)
700 {
701   this->Superclass::PrintSelf(os, indent);
702 }
703
704 int vtkExtractGroup::GetNumberOfGroupsFlagsArrays()
705 {
706   int ret(this->Internal->getNumberOfEntries());
707   //std::cerr << "vtkExtractGroup::GetNumberOfFieldsTreeArrays() -> " << ret << std::endl;
708   return ret;
709 }
710
711 const char *vtkExtractGroup::GetGroupsFlagsArrayName(int index)
712 {
713   const char *ret(this->Internal->getKeyOfEntry(index));
714 //  std::cerr << "vtkExtractGroup::GetFieldsTreeArrayName(" << index << ") -> " << ret << std::endl;
715   return ret;
716 }
717
718 int vtkExtractGroup::GetGroupsFlagsArrayStatus(const char *name)
719 {
720   int ret((int)this->Internal->getStatusOfEntryStr(name));
721 //  std::cerr << "vtkExtractGroup::GetGroupsFlagsArrayStatus(" << name << ") -> " << ret << std::endl;
722   return ret;
723 }
724
725 void vtkExtractGroup::SetGroupsFlagsStatus(const char *name, int status)
726 {
727   //std::cerr << "vtkExtractGroup::SetFieldsStatus(" << name << "," << status << ")" << std::endl;
728   this->Internal->setStatusOfEntryStr(name,(bool)status);
729   this->Modified();
730   //this->Internal->printMySelf(std::cerr);
731 }
732
733 const char *vtkExtractGroup::GetMeshName()
734 {
735   return this->Internal->getMeshName();
736 }