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