Salome HOME
Merge branch 'V7_main' of ssh://git.salome-platform.org/modules/paravis into V7_main
[modules/paravis.git] / src / Plugins / MEDReader / IO / vtkExtractGroup.cxx
1 // Copyright (C) 2010-2014  CEA/DEN, EDF R&D
2 //
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License.
7 //
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
11 // Lesser General Public License for more details.
12 //
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19 // 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
56 #include <map>
57 #include <deque>
58
59 vtkStandardNewMacro(vtkExtractGroup);
60
61 vtkCxxSetObjectMacro(vtkExtractGroup, SIL, vtkMutableDirectedGraph);
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 *getKeyOfEntry(int i) const;
115   bool getStatusOfEntryStr(const char *entry) const;
116   void setStatusOfEntryStr(const char *entry, bool status);
117   void printMySelf(std::ostream& os) const;
118   std::set<int> getIdsToKeep() const;
119   int getIdOfFamily(const std::string& famName) const;
120 private:
121   std::map<std::string,int> computeFamStrIdMap() const;
122   const ExtractGroupStatus& getEntry(const char *entry) const;
123   ExtractGroupStatus& getEntry(const char *entry);
124 private:
125   std::vector<ExtractGroupGrp> _groups;
126   std::vector<ExtractGroupFam> _fams;
127 };
128
129 const char ExtractGroupGrp::START[]="GRP_";
130
131 const char ExtractGroupFam::START[]="FAM_";
132
133 ExtractGroupStatus::ExtractGroupStatus(const char *name):_status(false),_name(name)
134 {
135 }
136
137 void ExtractGroupStatus::printMySelf(std::ostream& os) const
138 {
139   os << "      -" << _ze_key_name << "(";
140   if(_status)
141     os << "X";
142   else
143     os << " ";
144   os << ")" << std::endl;
145 }
146
147 bool ExtractGroupStatus::isSameAs(const ExtractGroupStatus& other) const
148 {
149   return _name==other._name && _ze_key_name==other._ze_key_name;
150 }
151
152 bool ExtractGroupGrp::isSameAs(const ExtractGroupGrp& other) const
153 {
154   bool ret(ExtractGroupStatus::isSameAs(other));
155   if(ret)
156     return _fams==other._fams;
157   else
158     return false;
159 }
160
161 void vtkExtractGroup::vtkExtractGroupInternal::loadFrom(vtkMutableDirectedGraph *sil)
162 {
163   std::vector<ExtractGroupGrp> oldGrps(_groups); _groups.clear();
164   std::vector<ExtractGroupFam> oldFams(_fams); _fams.clear();
165   int idNames(0);
166   vtkAbstractArray *verticesNames(sil->GetVertexData()->GetAbstractArray("Names",idNames));
167   vtkStringArray *verticesNames2(vtkStringArray::SafeDownCast(verticesNames));
168   vtkIdType id0;
169   bool found(false);
170   for(int i=0;i<verticesNames2->GetNumberOfValues();i++)
171     {
172       vtkStdString &st(verticesNames2->GetValue(i));
173       if(st=="MeshesFamsGrps")
174         {
175           id0=i;
176           found=true;
177         }
178     }
179   if(!found)
180     throw INTERP_KERNEL::Exception("There is an internal error ! The tree on server side has not the expected look !");
181   vtkAdjacentVertexIterator *it0(vtkAdjacentVertexIterator::New());
182   sil->GetAdjacentVertices(id0,it0);
183   int kk(0),ll(0);
184   while(it0->HasNext())
185     {
186       vtkIdType id1(it0->Next());
187       std::string meshName((const char *)verticesNames2->GetValue(id1));
188       vtkAdjacentVertexIterator *it1(vtkAdjacentVertexIterator::New());
189       sil->GetAdjacentVertices(id1,it1);
190       vtkIdType idZeGrps(it1->Next());//zeGroups
191       vtkAdjacentVertexIterator *itGrps(vtkAdjacentVertexIterator::New());
192       sil->GetAdjacentVertices(idZeGrps,itGrps);
193       while(itGrps->HasNext())
194         {
195           vtkIdType idg(itGrps->Next());
196           ExtractGroupGrp grp((const char *)verticesNames2->GetValue(idg));
197           vtkAdjacentVertexIterator *itGrps2(vtkAdjacentVertexIterator::New());
198           sil->GetAdjacentVertices(idg,itGrps2);
199           std::vector<std::string> famsOnGroup;
200           while(itGrps2->HasNext())
201             {
202               vtkIdType idgf(itGrps2->Next());
203               famsOnGroup.push_back(std::string((const char *)verticesNames2->GetValue(idgf)));
204             }
205           grp.setFamilies(famsOnGroup);
206           itGrps2->Delete();
207           _groups.push_back(grp);
208         }
209       itGrps->Delete();
210       vtkIdType idZeFams(it1->Next());//zeFams
211       it1->Delete();
212       vtkAdjacentVertexIterator *itFams(vtkAdjacentVertexIterator::New());
213       sil->GetAdjacentVertices(idZeFams,itFams);
214       while(itFams->HasNext())
215         {
216           vtkIdType idf(itFams->Next());
217           ExtractGroupFam fam((const char *)verticesNames2->GetValue(idf));
218           _fams.push_back(fam);
219         }
220       itFams->Delete();
221     }
222   it0->Delete(); 
223   //
224   std::size_t szg(_groups.size()),szf(_fams.size());
225   if(szg==oldGrps.size() && szf==oldFams.size())
226     {
227       bool isSame(true);
228       for(std::size_t i=0;i<szg && isSame;i++)
229         isSame=_groups[i].isSameAs(oldGrps[i]);
230       for(std::size_t i=0;i<szf && isSame;i++)
231         isSame=_fams[i].isSameAs(oldFams[i]);
232       if(isSame)
233         {
234           for(std::size_t i=0;i<szg;i++)
235             _groups[i].cpyStatusFrom(oldGrps[i]);
236           for(std::size_t i=0;i<szf;i++)
237             _fams[i].cpyStatusFrom(oldFams[i]);
238         }
239     }
240 }
241
242 int vtkExtractGroup::vtkExtractGroupInternal::getNumberOfEntries() const
243 {
244   std::size_t sz0(_groups.size()),sz1(_fams.size());
245   return (int)(sz0+sz1);
246 }
247
248 const char *vtkExtractGroup::vtkExtractGroupInternal::getKeyOfEntry(int i) const
249 {
250   int sz0((int)_groups.size());
251   if(i>=0 && i<sz0)
252     return _groups[i].getKeyOfEntry();
253   else
254     return _fams[i-sz0].getKeyOfEntry();
255 }
256
257 bool vtkExtractGroup::vtkExtractGroupInternal::getStatusOfEntryStr(const char *entry) const
258 {
259   const ExtractGroupStatus& elt(getEntry(entry));
260   return elt.getStatus();
261 }
262
263 void vtkExtractGroup::vtkExtractGroupInternal::setStatusOfEntryStr(const char *entry, bool status)
264 {
265   ExtractGroupStatus& elt(getEntry(entry));
266   elt.setStatus(status);
267 }
268
269 const ExtractGroupStatus& vtkExtractGroup::vtkExtractGroupInternal::getEntry(const char *entry) const
270 {
271   std::string entryCpp(entry);
272   for(std::vector<ExtractGroupGrp>::const_iterator it0=_groups.begin();it0!=_groups.end();it0++)
273     if(entryCpp==(*it0).getKeyOfEntry())
274       return *it0;
275   for(std::vector<ExtractGroupFam>::const_iterator it0=_fams.begin();it0!=_fams.end();it0++)
276     if(entryCpp==(*it0).getKeyOfEntry())
277       return *it0;
278   std::ostringstream oss; oss << "vtkExtractGroupInternal::getEntry : no such entry \"" << entry << "\"!";
279   throw INTERP_KERNEL::Exception(oss.str().c_str());
280 }
281
282 ExtractGroupStatus& vtkExtractGroup::vtkExtractGroupInternal::getEntry(const char *entry)
283 {
284   std::string entryCpp(entry);
285   for(std::vector<ExtractGroupGrp>::iterator it0=_groups.begin();it0!=_groups.end();it0++)
286     if(entryCpp==(*it0).getKeyOfEntry())
287       return *it0;
288   for(std::vector<ExtractGroupFam>::iterator it0=_fams.begin();it0!=_fams.end();it0++)
289     if(entryCpp==(*it0).getKeyOfEntry())
290       return *it0;
291   std::ostringstream oss; oss << "vtkExtractGroupInternal::getEntry : no such entry \"" << entry << "\"!";
292   throw INTERP_KERNEL::Exception(oss.str().c_str());
293 }
294
295 void vtkExtractGroup::vtkExtractGroupInternal::printMySelf(std::ostream& os) const
296 {
297   os << "Groups :" << std::endl;
298   for(std::vector<ExtractGroupGrp>::const_iterator it0=_groups.begin();it0!=_groups.end();it0++)
299     (*it0).printMySelf(os);
300   os << "Families :" << std::endl;
301   for(std::vector<ExtractGroupFam>::const_iterator it0=_fams.begin();it0!=_fams.end();it0++)
302     (*it0).printMySelf(os);
303 }
304
305 int vtkExtractGroup::vtkExtractGroupInternal::getIdOfFamily(const std::string& famName) const
306 {
307   for(std::vector<ExtractGroupFam>::const_iterator it=_fams.begin();it!=_fams.end();it++)
308     {
309       if((*it).getName()==famName)
310         return (*it).getId();
311     }
312 }
313
314 ExtractGroupFam::ExtractGroupFam(const char *name):ExtractGroupStatus(name),_id(0)
315 {
316   std::size_t pos(_name.find(MEDFileFieldRepresentationLeavesArrays::ZE_SEP));
317   std::string name0(_name.substr(0,pos)),name1(_name.substr(pos+strlen(MEDFileFieldRepresentationLeavesArrays::ZE_SEP)));
318   std::istringstream iss(name1);
319   iss >> _id;
320   std::ostringstream oss; oss << START << name; _ze_key_name=oss.str(); _name=name0;
321 }
322
323 bool ExtractGroupFam::isSameAs(const ExtractGroupFam& other) const
324 {
325   bool ret(ExtractGroupStatus::isSameAs(other));
326   if(ret)
327     return _id==other._id;
328   else
329     return false;
330 }
331
332 void ExtractGroupFam::printMySelf(std::ostream& os) const
333 {
334   os << "      -" << _ze_key_name << " famName : \"" << _name << "\" id : " << _id << " (";
335   if(_status)
336     os << "X";
337   else
338     os << " ";
339   os << ")" << std::endl;
340 }
341
342 void ExtractGroupFam::fillIdsToKeep(std::set<int>& s) const
343 {
344   s.insert(_id);
345 }
346
347 std::set<int> vtkExtractGroup::vtkExtractGroupInternal::getIdsToKeep() const
348 {
349   std::map<std::string,int> m(this->computeFamStrIdMap());
350   std::set<int> s;
351   for(std::vector<ExtractGroupGrp>::const_iterator it0=_groups.begin();it0!=_groups.end();it0++)
352     {
353       if((*it0).getStatus())
354         {
355           const std::vector<std::string>& fams((*it0).getFamiliesLyingOn());
356           for(std::vector<std::string>::const_iterator it1=fams.begin();it1!=fams.end();it1++)
357             {
358               std::map<std::string,int>::iterator it2(m.find((*it1)));
359               if(it2!=m.end())
360                 s.insert((*it2).second);
361             }
362         }
363      }
364   for(std::vector<ExtractGroupFam>::const_iterator it0=_fams.begin();it0!=_fams.end();it0++)
365     if((*it0).getStatus())
366       (*it0).fillIdsToKeep(s);
367   return s;
368 }
369
370 std::map<std::string,int> vtkExtractGroup::vtkExtractGroupInternal::computeFamStrIdMap() const
371 {
372   std::map<std::string,int> ret;
373   for(std::vector<ExtractGroupFam>::const_iterator it0=_fams.begin();it0!=_fams.end();it0++)
374     ret[(*it0).getName()]=(*it0).getId();
375   return ret;
376 }
377
378 ////////////////////
379
380 vtkExtractGroup::vtkExtractGroup():SIL(NULL),Internal(new vtkExtractGroupInternal),InsideOut(0)
381 {
382 }
383
384 vtkExtractGroup::~vtkExtractGroup()
385 {
386   delete this->Internal;
387 }
388
389 void vtkExtractGroup::SetInsideOut(int val)
390 {
391   if(this->InsideOut!=val)
392     {
393       this->InsideOut=val;
394       this->Modified();
395     }
396 }
397
398 int vtkExtractGroup::RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
399 {
400   try
401     {
402       //std::cerr << "########################################## vtkExtractGroup::RequestInformation ##########################################" << std::endl;
403       vtkInformation *outInfo(outputVector->GetInformationObject(0));
404       vtkInformation *inputInfo(inputVector[0]->GetInformationObject(0));//unfortunately inputInfo->Has(vtkDataObject::SIL) returns false... use executive to find it !
405       //
406       vtkExecutive *exe(GetExecutive());
407       vtkAlgorithm *alg(this);
408       vtkInformation *infoOnSIL(alg->GetOutputInformation(0));
409       while(!infoOnSIL->Has(vtkDataObject::SIL()))// skipping vtkPVPostFilter
410         {
411           if(exe->GetNumberOfInputConnections(0)<1)
412             {
413               vtkErrorMacro("No SIL Data available ! The source of this filter must be MEDReader !");
414               return 0;
415             }
416           vtkExecutive *exe2(exe->GetInputExecutive(0,0));
417           //
418           alg=exe2->GetAlgorithm(); exe=exe2; infoOnSIL=alg->GetOutputInformation(0);
419         }
420       //
421       this->SetSIL(vtkMutableDirectedGraph::SafeDownCast(infoOnSIL->Get(vtkDataObject::SIL())));
422       this->Internal->loadFrom(this->SIL);
423       //this->Internal->printMySelf(std::cerr);
424       outInfo->Set(vtkDataObject::SIL(),this->SIL);
425     }
426   catch(INTERP_KERNEL::Exception& e)
427     {
428       std::cerr << "Exception has been thrown in vtkExtractGroup::RequestInformation : " << e.what() << std::endl;
429       return 0;
430     }
431   return 1;
432 }
433
434 template<class CellPointExtractor>
435 vtkDataSet *FilterFamilies(vtkDataSet *input, const std::set<int>& idsToKeep, bool insideOut, const char *arrNameOfFamilyField,
436                            const char *associationForThreshold, bool& catchAll, bool& catchSmth)
437 {
438   static const int VTK_DATA_ARRAY_DELETE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_DELETE;
439   static const char ZE_SELECTION_ARR_NAME[]="@@ZeSelection@@";
440   vtkDataSet *output(input->NewInstance());
441   output->ShallowCopy(input);
442   vtkSmartPointer<vtkThreshold> thres(vtkSmartPointer<vtkThreshold>::New());
443   thres->SetInputData(output);
444   vtkDataSetAttributes *dscIn(input->GetCellData()),*dscIn2(input->GetPointData());
445   vtkDataSetAttributes *dscOut(output->GetCellData()),*dscOut2(output->GetPointData());
446   //
447   double vMin(insideOut==0?1.:0.),vMax(insideOut==0?2.:1.);
448   thres->ThresholdBetween(vMin,vMax);
449   // OK for the output 
450   //
451   CellPointExtractor cpe2(input);
452   vtkDataArray *da(cpe2.Get()->GetScalars(arrNameOfFamilyField));
453   if(!da)
454     return 0;
455   std::string daName(da->GetName());
456   vtkIntArray *dai(vtkIntArray::SafeDownCast(da));
457   if(daName!=arrNameOfFamilyField || !dai)
458     return 0;
459   //
460   int nbOfTuples(dai->GetNumberOfTuples());
461   vtkCharArray *zeSelection(vtkCharArray::New());
462   zeSelection->SetName(ZE_SELECTION_ARR_NAME);
463   zeSelection->SetNumberOfComponents(1);
464   char *pt(new char[nbOfTuples]);
465   zeSelection->SetArray(pt,nbOfTuples,0,VTK_DATA_ARRAY_DELETE);
466   const int *inPtr(dai->GetPointer(0));
467   std::fill(pt,pt+nbOfTuples,0);
468   catchAll=true; catchSmth=false;
469   std::vector<bool> pt2(nbOfTuples,false);
470   for(std::set<int>::const_iterator it=idsToKeep.begin();it!=idsToKeep.end();it++)
471     {
472       bool catchFid(false);
473       for(int i=0;i<nbOfTuples;i++)
474         if(inPtr[i]==*it)
475           { pt2[i]=true; catchFid=true; }
476       if(!catchFid)
477         catchAll=false;
478       else
479         catchSmth=true;
480     }
481   for(int ii=0;ii<nbOfTuples;ii++)
482     if(pt2[ii])
483       pt[ii]=2;
484   CellPointExtractor cpe3(output);
485   int idx(cpe3.Get()->AddArray(zeSelection));
486   cpe3.Get()->SetActiveAttribute(idx,vtkDataSetAttributes::SCALARS);
487   cpe3.Get()->CopyScalarsOff();
488   zeSelection->Delete();
489   //
490   thres->SetInputArrayToProcess(idx,0,0,associationForThreshold,ZE_SELECTION_ARR_NAME);
491   thres->Update();
492   vtkUnstructuredGrid *zeComputedOutput(thres->GetOutput());
493   CellPointExtractor cpe(zeComputedOutput);
494   cpe.Get()->RemoveArray(idx);
495   output->Delete();
496   zeComputedOutput->Register(0);
497   return zeComputedOutput;
498 }
499
500 class CellExtractor
501 {
502 public:
503   CellExtractor(vtkDataSet *ds):_ds(ds) { }
504   vtkDataSetAttributes *Get() { return _ds->GetCellData(); }
505 private:
506   vtkDataSet *_ds;
507 };
508
509 class PointExtractor
510 {
511 public:
512   PointExtractor(vtkDataSet *ds):_ds(ds) { }
513   vtkDataSetAttributes *Get() { return _ds->GetPointData(); }
514 private:
515   vtkDataSet *_ds;
516 };
517
518 int vtkExtractGroup::RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
519 {
520   try
521     {
522       //std::cerr << "########################################## vtkExtractGroup::RequestData        ##########################################" << std::endl;
523       vtkInformation* inputInfo=inputVector[0]->GetInformationObject(0);
524       vtkDataSet *input(vtkDataSet::SafeDownCast(inputInfo->Get(vtkDataObject::DATA_OBJECT())));
525       vtkInformation *info(input->GetInformation());
526       vtkInformation *outInfo(outputVector->GetInformationObject(0));
527       vtkDataSet *output(vtkDataSet::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT())));
528       std::set<int> idsToKeep(this->Internal->getIdsToKeep());
529       // first shrink the input
530       bool catchAll,catchSmth;
531       vtkDataSet *tryOnCell(FilterFamilies<CellExtractor>(input,idsToKeep,this->InsideOut,
532                                                           MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_CELL_NAME,"vtkDataObject::FIELD_ASSOCIATION_CELLS",catchAll,catchSmth));
533       if(tryOnCell)
534         {
535           if(catchAll)
536             {
537               output->ShallowCopy(tryOnCell);
538               tryOnCell->Delete();//
539               return 1;
540             }
541           else
542             {
543               if(catchSmth)
544                 {
545                   vtkDataSet *tryOnNode(FilterFamilies<PointExtractor>(tryOnCell,idsToKeep,this->InsideOut,
546                                                                        MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME,"vtkDataObject::FIELD_ASSOCIATION_POINTS",catchAll,catchSmth));
547                   if(tryOnNode && catchSmth)
548                     {
549                       output->ShallowCopy(tryOnNode);
550                       tryOnCell->Delete();
551                       tryOnNode->Delete();//
552                       return 1;
553                     }
554                   else
555                     {
556                       if(tryOnNode)
557                         tryOnNode->Delete();
558                       output->ShallowCopy(tryOnCell);
559                       tryOnCell->Delete();
560                       return 1;
561                     }
562                 }
563               else
564                 {
565                   vtkDataSet *tryOnNode(FilterFamilies<PointExtractor>(input,idsToKeep,this->InsideOut,
566                                                                        MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME,"vtkDataObject::FIELD_ASSOCIATION_POINTS",catchAll,catchSmth));
567                   if(tryOnNode)
568                     {
569                       tryOnCell->Delete();
570                       output->ShallowCopy(tryOnNode);
571                       tryOnNode->Delete();
572                       return 1;
573                     }
574                   else
575                     {
576                       output->ShallowCopy(tryOnCell);
577                       tryOnCell->Delete();
578                       return 0;
579                     }
580                 }
581             }
582         }
583       else
584         {
585           vtkDataSet *tryOnNode(FilterFamilies<PointExtractor>(tryOnCell,idsToKeep,this->InsideOut,
586                                                                MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME,"vtkDataObject::FIELD_ASSOCIATION_POINTS",catchAll,catchSmth));
587           if(tryOnNode)
588             {
589               output->ShallowCopy(tryOnNode);
590               tryOnNode->Delete();//
591               return 1;
592             }
593           else
594             {
595               std::ostringstream oss; oss << "vtkExtractGroup::RequestData : The integer array with name \""<< MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_CELL_NAME;
596               oss << "\" or \"" << MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME << "\" does not exist ! The extraction of group and/or family is not possible !";
597               if(this->HasObserver("ErrorEvent") )
598                 this->InvokeEvent("ErrorEvent",const_cast<char *>(oss.str().c_str()));
599               else
600                 vtkOutputWindowDisplayErrorText(const_cast<char *>(oss.str().c_str()));
601               vtkObject::BreakOnError();
602               return 0;
603             }
604         }
605     }
606   catch(INTERP_KERNEL::Exception& e)
607     {
608       std::cerr << "Exception has been thrown in vtkExtractGroup::RequestData : " << e.what() << std::endl;
609       return 0;
610     }
611 }
612
613 int vtkExtractGroup::GetSILUpdateStamp()
614 {
615   return this->SILTime;
616 }
617
618 void vtkExtractGroup::PrintSelf(ostream& os, vtkIndent indent)
619 {
620   this->Superclass::PrintSelf(os, indent);
621 }
622
623 int vtkExtractGroup::GetNumberOfGroupsFlagsArrays()
624 {
625   int ret(this->Internal->getNumberOfEntries());
626   //std::cerr << "vtkExtractGroup::GetNumberOfFieldsTreeArrays() -> " << ret << std::endl;
627   return ret;
628 }
629
630 const char *vtkExtractGroup::GetGroupsFlagsArrayName(int index)
631 {
632   const char *ret(this->Internal->getKeyOfEntry(index));
633   //std::cerr << "vtkExtractGroup::GetFieldsTreeArrayName(" << index << ") -> " << ret << std::endl;
634   return ret;
635 }
636
637 int vtkExtractGroup::GetGroupsFlagsArrayStatus(const char *name)
638 {
639   int ret((int)this->Internal->getStatusOfEntryStr(name));
640   //std::cerr << "vtkExtractGroup::GetGroupsFlagsArrayStatus(" << name << ") -> " << ret << std::endl;
641   return ret;
642 }
643
644 void vtkExtractGroup::SetGroupsFlagsStatus(const char *name, int status)
645 {
646   //std::cerr << "vtkExtractGroup::SetFieldsStatus(" << name << "," << status << ")" << std::endl;
647   this->Internal->setStatusOfEntryStr(name,(bool)status);
648   if(std::string(name)==GetGroupsFlagsArrayName(GetNumberOfGroupsFlagsArrays()-1))
649      this->Modified();
650   //this->Internal->printMySelf(std::cerr);
651 }