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