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