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