1 // Copyright (C) 2010-2015 CEA/DEN, EDF R&D
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.
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.
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
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
19 // Author : Anthony Geay
21 #include "vtkExtractGroup.h"
22 #include "MEDFileFieldRepresentationTree.hxx"
23 #include "vtkMEDReader.h"
25 #include "vtkAdjacentVertexIterator.h"
26 #include "vtkIntArray.h"
27 #include "vtkCellData.h"
28 #include "vtkPointData.h"
30 #include "vtkStreamingDemandDrivenPipeline.h"
31 #include "vtkUnstructuredGrid.h"
32 #include "vtkMultiBlockDataSet.h"
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"
63 vtkStandardNewMacro(vtkExtractGroup);
67 class ExtractGroupStatus
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;
82 std::string _ze_key_name;
85 class ExtractGroupGrp : public ExtractGroupStatus
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;
93 static const char START[];
94 std::vector<std::string> _fams;
97 class ExtractGroupFam : public ExtractGroupStatus
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;
106 static const char START[];
111 class vtkExtractGroup::vtkExtractGroupInternal
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);
125 std::map<std::string,int> computeFamStrIdMap() const;
126 const ExtractGroupStatus& getEntry(const char *entry) const;
127 ExtractGroupStatus& getEntry(const char *entry);
129 std::vector<ExtractGroupGrp> _groups;
130 std::vector<ExtractGroupFam> _fams;
131 std::string _mesh_name;
134 const char ExtractGroupGrp::START[]="GRP_";
136 const char ExtractGroupFam::START[]="FAM_";
138 ExtractGroupStatus::ExtractGroupStatus(const char *name):_status(false),_name(name)
142 void ExtractGroupStatus::printMySelf(std::ostream& os) const
144 os << " -" << _ze_key_name << "(";
149 os << ")" << std::endl;
152 bool ExtractGroupStatus::isSameAs(const ExtractGroupStatus& other) const
154 return _name==other._name && _ze_key_name==other._ze_key_name;
157 bool ExtractGroupGrp::isSameAs(const ExtractGroupGrp& other) const
159 bool ret(ExtractGroupStatus::isSameAs(other));
161 return _fams==other._fams;
166 bool vtkExtractGroup::vtkExtractGroupInternal::IsInformationOK(vtkInformation *info)
168 // Check the information contain meta data key
169 if(!info->Has(vtkMEDReader::META_DATA()))
173 vtkMutableDirectedGraph *sil(vtkMutableDirectedGraph::SafeDownCast(info->Get(vtkMEDReader::META_DATA())));
177 vtkAbstractArray *verticesNames(sil->GetVertexData()->GetAbstractArray("Names",idNames));
178 vtkStringArray *verticesNames2(vtkStringArray::SafeDownCast(verticesNames));
181 for(int i=0;i<verticesNames2->GetNumberOfValues();i++)
183 vtkStdString &st(verticesNames2->GetValue(i));
184 if(st=="MeshesFamsGrps")
190 const char* vtkExtractGroup::GetGrpStart()
192 return ExtractGroupGrp::START;
195 const char* vtkExtractGroup::GetFamStart()
197 return ExtractGroupFam::START;
200 const char *vtkExtractGroup::vtkExtractGroupInternal::getMeshName() const
202 return this->_mesh_name.c_str();
205 void vtkExtractGroup::vtkExtractGroupInternal::loadFrom(vtkMutableDirectedGraph *sil)
207 std::vector<ExtractGroupGrp> oldGrps(_groups); _groups.clear();
208 std::vector<ExtractGroupFam> oldFams(_fams); _fams.clear();
210 vtkAbstractArray *verticesNames(sil->GetVertexData()->GetAbstractArray("Names",idNames));
211 vtkStringArray *verticesNames2(vtkStringArray::SafeDownCast(verticesNames));
214 for(int i=0;i<verticesNames2->GetNumberOfValues();i++)
216 vtkStdString &st(verticesNames2->GetValue(i));
217 if(st=="MeshesFamsGrps")
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);
228 while(it0->HasNext())
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())
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())
247 vtkIdType idgf(itGrps2->Next());
248 famsOnGroup.push_back(std::string((const char *)verticesNames2->GetValue(idgf)));
250 grp.setFamilies(famsOnGroup);
252 _groups.push_back(grp);
255 vtkIdType idZeFams(it1->Next());//zeFams
257 vtkAdjacentVertexIterator *itFams(vtkAdjacentVertexIterator::New());
258 sil->GetAdjacentVertices(idZeFams,itFams);
259 while(itFams->HasNext())
261 vtkIdType idf(itFams->Next());
262 ExtractGroupFam fam((const char *)verticesNames2->GetValue(idf));
263 _fams.push_back(fam);
269 std::size_t szg(_groups.size()),szf(_fams.size());
270 if(szg==oldGrps.size() && szf==oldFams.size())
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]);
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]);
287 int vtkExtractGroup::vtkExtractGroupInternal::getNumberOfEntries() const
289 std::size_t sz0(_groups.size()),sz1(_fams.size());
290 return (int)(sz0+sz1);
293 const char *vtkExtractGroup::vtkExtractGroupInternal::getKeyOfEntry(int i) const
295 int sz0((int)_groups.size());
297 return _groups[i].getKeyOfEntry();
299 return _fams[i-sz0].getKeyOfEntry();
302 bool vtkExtractGroup::vtkExtractGroupInternal::getStatusOfEntryStr(const char *entry) const
304 const ExtractGroupStatus& elt(getEntry(entry));
305 return elt.getStatus();
308 void vtkExtractGroup::vtkExtractGroupInternal::setStatusOfEntryStr(const char *entry, bool status)
310 ExtractGroupStatus& elt(getEntry(entry));
311 elt.setStatus(status);
314 const ExtractGroupStatus& vtkExtractGroup::vtkExtractGroupInternal::getEntry(const char *entry) const
316 std::string entryCpp(entry);
317 for(std::vector<ExtractGroupGrp>::const_iterator it0=_groups.begin();it0!=_groups.end();it0++)
318 if(entryCpp==(*it0).getKeyOfEntry())
320 for(std::vector<ExtractGroupFam>::const_iterator it0=_fams.begin();it0!=_fams.end();it0++)
321 if(entryCpp==(*it0).getKeyOfEntry())
323 std::ostringstream oss; oss << "vtkExtractGroupInternal::getEntry : no such entry \"" << entry << "\"!";
324 throw INTERP_KERNEL::Exception(oss.str().c_str());
327 ExtractGroupStatus& vtkExtractGroup::vtkExtractGroupInternal::getEntry(const char *entry)
329 std::string entryCpp(entry);
330 for(std::vector<ExtractGroupGrp>::iterator it0=_groups.begin();it0!=_groups.end();it0++)
331 if(entryCpp==(*it0).getKeyOfEntry())
333 for(std::vector<ExtractGroupFam>::iterator it0=_fams.begin();it0!=_fams.end();it0++)
334 if(entryCpp==(*it0).getKeyOfEntry())
336 std::ostringstream oss; oss << "vtkExtractGroupInternal::getEntry : no such entry \"" << entry << "\"!";
337 throw INTERP_KERNEL::Exception(oss.str().c_str());
340 void vtkExtractGroup::vtkExtractGroupInternal::printMySelf(std::ostream& os) const
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);
350 int vtkExtractGroup::vtkExtractGroupInternal::getIdOfFamily(const std::string& famName) const
352 for(std::vector<ExtractGroupFam>::const_iterator it=_fams.begin();it!=_fams.end();it++)
354 if((*it).getName()==famName)
355 return (*it).getId();
359 ExtractGroupFam::ExtractGroupFam(const char *name):ExtractGroupStatus(name),_id(0)
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);
365 std::ostringstream oss; oss << START << name; _ze_key_name=oss.str(); _name=name0;
368 bool ExtractGroupFam::isSameAs(const ExtractGroupFam& other) const
370 bool ret(ExtractGroupStatus::isSameAs(other));
372 return _id==other._id;
377 void ExtractGroupFam::printMySelf(std::ostream& os) const
379 os << " -" << _ze_key_name << " famName : \"" << _name << "\" id : " << _id << " (";
384 os << ")" << std::endl;
387 void ExtractGroupFam::fillIdsToKeep(std::set<int>& s) const
392 std::set<int> vtkExtractGroup::vtkExtractGroupInternal::getIdsToKeep() const
394 std::map<std::string,int> m(this->computeFamStrIdMap());
396 for(std::vector<ExtractGroupGrp>::const_iterator it0=_groups.begin();it0!=_groups.end();it0++)
398 if((*it0).getStatus())
400 const std::vector<std::string>& fams((*it0).getFamiliesLyingOn());
401 for(std::vector<std::string>::const_iterator it1=fams.begin();it1!=fams.end();it1++)
403 std::map<std::string,int>::iterator it2(m.find((*it1)));
405 s.insert((*it2).second);
409 for(std::vector<ExtractGroupFam>::const_iterator it0=_fams.begin();it0!=_fams.end();it0++)
410 if((*it0).getStatus())
411 (*it0).fillIdsToKeep(s);
415 std::map<std::string,int> vtkExtractGroup::vtkExtractGroupInternal::computeFamStrIdMap() const
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();
425 vtkExtractGroup::vtkExtractGroup():SIL(NULL),Internal(new vtkExtractGroupInternal),InsideOut(0)
429 vtkExtractGroup::~vtkExtractGroup()
431 delete this->Internal;
434 void vtkExtractGroup::SetInsideOut(int val)
436 if(this->InsideOut!=val)
443 int vtkExtractGroup::RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
445 // vtkUnstructuredGridAlgorithm::RequestInformation(request,inputVector,outputVector);
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))
454 vtkErrorMacro("No SIL Data available ! The source of this filter must be MEDReader !");
458 this->SetSIL(vtkMutableDirectedGraph::SafeDownCast(inputInfo->Get(vtkMEDReader::META_DATA())));
459 this->Internal->loadFrom(this->SIL);
460 //this->Internal->printMySelf(std::cerr);
462 catch(INTERP_KERNEL::Exception& e)
464 std::cerr << "Exception has been thrown in vtkExtractGroup::RequestInformation : " << e.what() << std::endl;
471 * Do not use vtkCxxSetObjectMacro macro because input mdg comes from an already managed in the pipeline just a ref on it.
473 void vtkExtractGroup::SetSIL(vtkMutableDirectedGraph *mdg)
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)
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());
493 double vMin(insideOut==0?1.:0.),vMax(insideOut==0?2.:1.);
494 thres->ThresholdBetween(vMin,vMax);
497 CellPointExtractor cpe2(input);
498 vtkDataArray *da(cpe2.Get()->GetScalars(arrNameOfFamilyField));
501 std::string daName(da->GetName());
502 vtkIntArray *dai(vtkIntArray::SafeDownCast(da));
503 if(daName!=arrNameOfFamilyField || !dai)
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++)
518 bool catchFid(false);
519 for(int i=0;i<nbOfTuples;i++)
521 { pt2[i]=true; catchFid=true; }
527 for(int ii=0;ii<nbOfTuples;ii++)
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();
536 thres->SetInputArrayToProcess(idx,0,0,associationForThreshold,ZE_SELECTION_ARR_NAME);
538 vtkUnstructuredGrid *zeComputedOutput(thres->GetOutput());
539 CellPointExtractor cpe(zeComputedOutput);
540 cpe.Get()->RemoveArray(idx);
542 zeComputedOutput->Register(0);
543 return zeComputedOutput;
549 CellExtractor(vtkDataSet *ds):_ds(ds) { }
550 vtkDataSetAttributes *Get() { return _ds->GetCellData(); }
558 PointExtractor(vtkDataSet *ds):_ds(ds) { }
559 vtkDataSetAttributes *Get() { return _ds->GetPointData(); }
564 int vtkExtractGroup::RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
568 // std::cerr << "########################################## vtkExtractGroup::RequestData ##########################################" << std::endl;
569 // request->Print(cout);
570 vtkInformation* inputInfo=inputVector[0]->GetInformationObject(0);
571 vtkDataSet *input(vtkDataSet::SafeDownCast(inputInfo->Get(vtkDataObject::DATA_OBJECT())));
572 vtkInformation *info(input->GetInformation());
573 vtkInformation *outInfo(outputVector->GetInformationObject(0));
574 vtkDataSet *output(vtkDataSet::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT())));
575 std::set<int> idsToKeep(this->Internal->getIdsToKeep());
576 // first shrink the input
577 bool catchAll,catchSmth;
578 vtkSmartPointer<vtkThreshold> thres1(vtkSmartPointer<vtkThreshold>::New()),thres2(vtkSmartPointer<vtkThreshold>::New());
579 vtkDataSet *tryOnCell(FilterFamilies<CellExtractor>(thres1,input,idsToKeep,this->InsideOut,
580 MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_CELL_NAME,"vtkDataObject::FIELD_ASSOCIATION_CELLS",catchAll,catchSmth));
585 output->ShallowCopy(tryOnCell);
586 tryOnCell->Delete();//
593 vtkDataSet *tryOnNode(FilterFamilies<PointExtractor>(thres2,input,idsToKeep,this->InsideOut,
594 MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME,"vtkDataObject::FIELD_ASSOCIATION_POINTS",catchAll,catchSmth));
595 if(tryOnNode && catchSmth)
597 vtkSmartPointer<vtkMultiBlockDataGroupFilter> mb(vtkSmartPointer<vtkMultiBlockDataGroupFilter>::New());
598 vtkSmartPointer<vtkCompositeDataToUnstructuredGridFilter> cd(vtkSmartPointer<vtkCompositeDataToUnstructuredGridFilter>::New());
599 mb->AddInputConnection(thres1->GetOutputPort());
600 mb->AddInputConnection(thres2->GetOutputPort());
601 cd->SetInputConnection(mb->GetOutputPort());
602 cd->SetMergePoints(0);
604 output->ShallowCopy(cd->GetOutput());
613 output->ShallowCopy(tryOnCell);
620 vtkDataSet *tryOnNode(FilterFamilies<PointExtractor>(thres1,input,idsToKeep,this->InsideOut,
621 MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME,"vtkDataObject::FIELD_ASSOCIATION_POINTS",catchAll,catchSmth));
625 output->ShallowCopy(tryOnNode);
631 output->ShallowCopy(tryOnCell);
640 vtkDataSet *tryOnNode(FilterFamilies<PointExtractor>(thres1,input,idsToKeep,this->InsideOut,
641 MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME,"vtkDataObject::FIELD_ASSOCIATION_POINTS",catchAll,catchSmth));
644 output->ShallowCopy(tryOnNode);
645 tryOnNode->Delete();//
650 std::ostringstream oss; oss << "vtkExtractGroup::RequestData : The integer array with name \""<< MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_CELL_NAME;
651 oss << "\" or \"" << MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME << "\" does not exist ! The extraction of group and/or family is not possible !";
652 if(this->HasObserver("ErrorEvent") )
653 this->InvokeEvent("ErrorEvent",const_cast<char *>(oss.str().c_str()));
655 vtkOutputWindowDisplayErrorText(const_cast<char *>(oss.str().c_str()));
656 vtkObject::BreakOnError();
661 catch(INTERP_KERNEL::Exception& e)
663 std::cerr << "Exception has been thrown in vtkExtractGroup::RequestData : " << e.what() << std::endl;
668 int vtkExtractGroup::GetSILUpdateStamp()
670 return this->SILTime;
673 void vtkExtractGroup::PrintSelf(ostream& os, vtkIndent indent)
675 this->Superclass::PrintSelf(os, indent);
678 int vtkExtractGroup::GetNumberOfGroupsFlagsArrays()
680 int ret(this->Internal->getNumberOfEntries());
681 //std::cerr << "vtkExtractGroup::GetNumberOfFieldsTreeArrays() -> " << ret << std::endl;
685 const char *vtkExtractGroup::GetGroupsFlagsArrayName(int index)
687 const char *ret(this->Internal->getKeyOfEntry(index));
688 // std::cerr << "vtkExtractGroup::GetFieldsTreeArrayName(" << index << ") -> " << ret << std::endl;
692 int vtkExtractGroup::GetGroupsFlagsArrayStatus(const char *name)
694 int ret((int)this->Internal->getStatusOfEntryStr(name));
695 // std::cerr << "vtkExtractGroup::GetGroupsFlagsArrayStatus(" << name << ") -> " << ret << std::endl;
699 void vtkExtractGroup::SetGroupsFlagsStatus(const char *name, int status)
701 //std::cerr << "vtkExtractGroup::SetFieldsStatus(" << name << "," << status << ")" << std::endl;
702 if (GetNumberOfGroupsFlagsArrays()<1)
704 this->Internal->setStatusOfEntryStr(name,(bool)status);
705 if(std::string(name)==GetGroupsFlagsArrayName(GetNumberOfGroupsFlagsArrays()-1))
707 //this->Internal->printMySelf(std::cerr);
710 const char *vtkExtractGroup::GetMeshName()
712 return this->Internal->getMeshName();