1 // Copyright (C) 2010-2019 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"
24 #include "VTKMEDTraits.hxx"
26 #include "vtkAdjacentVertexIterator.h"
27 #include "vtkDataArrayTemplate.h"
28 #include "vtkIntArray.h"
29 #include "vtkLongArray.h"
30 #include "vtkCellData.h"
31 #include "vtkPointData.h"
33 #include "vtkStreamingDemandDrivenPipeline.h"
34 #include "vtkUnstructuredGrid.h"
35 #include "vtkMultiBlockDataSet.h"
37 #include "vtkInformationStringKey.h"
38 #include "vtkAlgorithmOutput.h"
39 #include "vtkObjectFactory.h"
40 #include "vtkMutableDirectedGraph.h"
41 #include "vtkMultiBlockDataSet.h"
42 #include "vtkDataSet.h"
43 #include "vtkInformationVector.h"
44 #include "vtkInformation.h"
45 #include "vtkDataArraySelection.h"
46 #include "vtkTimeStamp.h"
47 #include "vtkInEdgeIterator.h"
48 #include "vtkInformationDataObjectKey.h"
49 #include "vtkExecutive.h"
50 #include "vtkVariantArray.h"
51 #include "vtkStringArray.h"
52 #include "vtkDoubleArray.h"
53 #include "vtkCharArray.h"
54 #include "vtkUnsignedCharArray.h"
55 #include "vtkDataSetAttributes.h"
56 #include "vtkDemandDrivenPipeline.h"
57 #include "vtkDataObjectTreeIterator.h"
58 #include "vtkThreshold.h"
59 #include "vtkMultiBlockDataGroupFilter.h"
60 #include "vtkCompositeDataToUnstructuredGridFilter.h"
61 #include "vtkInformationDataObjectMetaDataKey.h"
66 vtkStandardNewMacro(vtkExtractGroup);
70 class ExtractGroupStatus
73 ExtractGroupStatus():_status(false) { }
74 ExtractGroupStatus(const char *name);
75 bool getStatus() const { return _status; }
76 void setStatus(bool status) const { _status=status; }
77 void cpyStatusFrom(const ExtractGroupStatus& other) { _status=other._status; }
78 std::string getName() const { return _name; }
79 void resetStatus() const { _status=false; }
80 const char *getKeyOfEntry() const { return _ze_key_name.c_str(); }
81 virtual void printMySelf(std::ostream& os) const;
82 virtual bool isSameAs(const ExtractGroupStatus& other) const;
86 std::string _ze_key_name;
89 class ExtractGroupGrp : public ExtractGroupStatus
92 ExtractGroupGrp(const char *name):ExtractGroupStatus(name) { std::ostringstream oss; oss << START << name; _ze_key_name=oss.str(); }
93 void setFamilies(const std::vector<std::string>& fams) { _fams=fams; }
94 const std::vector<std::string>& getFamiliesLyingOn() const { return _fams; }
95 bool isSameAs(const ExtractGroupGrp& other) const;
97 static const char START[];
98 std::vector<std::string> _fams;
101 class ExtractGroupFam : public ExtractGroupStatus
104 ExtractGroupFam(const char *name);
105 void printMySelf(std::ostream& os) const;
106 void fillIdsToKeep(std::set<int>& s) const;
107 int getId() const { return _id; }
108 bool isSameAs(const ExtractGroupFam& other) const;
110 static const char START[];
115 class vtkExtractGroup::vtkExtractGroupInternal
118 void loadFrom(vtkMutableDirectedGraph *sil);
119 int getNumberOfEntries() const;
120 const char *getMeshName() const;
121 const char *getKeyOfEntry(int i) const;
122 bool getStatusOfEntryStr(const char *entry) const;
123 void setStatusOfEntryStr(const char *entry, bool status);
124 void printMySelf(std::ostream& os) const;
125 std::set<int> getIdsToKeep() const;
126 void clearSelection() const;
127 int getIdOfFamily(const std::string& famName) const;
128 static bool IsInformationOK(vtkInformation *info);
130 std::map<std::string,int> computeFamStrIdMap() const;
131 const ExtractGroupStatus& getEntry(const char *entry) const;
132 ExtractGroupStatus& getEntry(const char *entry);
134 std::vector<ExtractGroupGrp> _groups;
135 std::vector<ExtractGroupFam> _fams;
136 mutable std::vector< std::pair<std::string,bool> > _selection;
137 std::string _mesh_name;
140 const char ExtractGroupGrp::START[]="GRP_";
142 const char ExtractGroupFam::START[]="FAM_";
144 ExtractGroupStatus::ExtractGroupStatus(const char *name):_status(false),_name(name)
148 void ExtractGroupStatus::printMySelf(std::ostream& os) const
150 os << " -" << _ze_key_name << "(";
155 os << ")" << std::endl;
158 bool ExtractGroupStatus::isSameAs(const ExtractGroupStatus& other) const
160 return _name==other._name && _ze_key_name==other._ze_key_name;
163 bool ExtractGroupGrp::isSameAs(const ExtractGroupGrp& other) const
165 bool ret(ExtractGroupStatus::isSameAs(other));
167 return _fams==other._fams;
172 bool vtkExtractGroup::vtkExtractGroupInternal::IsInformationOK(vtkInformation *info)
174 // Check the information contain meta data key
175 if(!info->Has(vtkMEDReader::META_DATA()))
179 vtkMutableDirectedGraph *sil(vtkMutableDirectedGraph::SafeDownCast(info->Get(vtkMEDReader::META_DATA())));
183 vtkAbstractArray *verticesNames(sil->GetVertexData()->GetAbstractArray("Names",idNames));
184 vtkStringArray *verticesNames2(vtkStringArray::SafeDownCast(verticesNames));
187 for(int i=0;i<verticesNames2->GetNumberOfValues();i++)
189 vtkStdString &st(verticesNames2->GetValue(i));
190 if(st=="MeshesFamsGrps")
196 const char* vtkExtractGroup::GetGrpStart()
198 return ExtractGroupGrp::START;
201 const char* vtkExtractGroup::GetFamStart()
203 return ExtractGroupFam::START;
206 const char *vtkExtractGroup::vtkExtractGroupInternal::getMeshName() const
208 return this->_mesh_name.c_str();
211 void vtkExtractGroup::vtkExtractGroupInternal::loadFrom(vtkMutableDirectedGraph *sil)
213 std::vector<ExtractGroupGrp> oldGrps(_groups); _groups.clear();
214 std::vector<ExtractGroupFam> oldFams(_fams); _fams.clear();
216 vtkAbstractArray *verticesNames(sil->GetVertexData()->GetAbstractArray("Names",idNames));
217 vtkStringArray *verticesNames2(vtkStringArray::SafeDownCast(verticesNames));
220 for(int i=0;i<verticesNames2->GetNumberOfValues();i++)
222 vtkStdString &st(verticesNames2->GetValue(i));
223 if(st=="MeshesFamsGrps")
230 throw INTERP_KERNEL::Exception("There is an internal error ! The tree on server side has not the expected look !");
231 vtkAdjacentVertexIterator *it0(vtkAdjacentVertexIterator::New());
232 sil->GetAdjacentVertices(id0,it0);
234 while(it0->HasNext())
236 vtkIdType id1(it0->Next());
237 std::string meshName((const char *)verticesNames2->GetValue(id1));
238 this->_mesh_name=meshName;
239 vtkAdjacentVertexIterator *it1(vtkAdjacentVertexIterator::New());
240 sil->GetAdjacentVertices(id1,it1);
241 vtkIdType idZeGrps(it1->Next());//zeGroups
242 vtkAdjacentVertexIterator *itGrps(vtkAdjacentVertexIterator::New());
243 sil->GetAdjacentVertices(idZeGrps,itGrps);
244 while(itGrps->HasNext())
246 vtkIdType idg(itGrps->Next());
247 ExtractGroupGrp grp((const char *)verticesNames2->GetValue(idg));
248 vtkAdjacentVertexIterator *itGrps2(vtkAdjacentVertexIterator::New());
249 sil->GetAdjacentVertices(idg,itGrps2);
250 std::vector<std::string> famsOnGroup;
251 while(itGrps2->HasNext())
253 vtkIdType idgf(itGrps2->Next());
254 famsOnGroup.push_back(std::string((const char *)verticesNames2->GetValue(idgf)));
256 grp.setFamilies(famsOnGroup);
258 _groups.push_back(grp);
261 vtkIdType idZeFams(it1->Next());//zeFams
263 vtkAdjacentVertexIterator *itFams(vtkAdjacentVertexIterator::New());
264 sil->GetAdjacentVertices(idZeFams,itFams);
265 while(itFams->HasNext())
267 vtkIdType idf(itFams->Next());
268 ExtractGroupFam fam((const char *)verticesNames2->GetValue(idf));
269 _fams.push_back(fam);
275 std::size_t szg(_groups.size()),szf(_fams.size());
276 if(szg==oldGrps.size() && szf==oldFams.size())
279 for(std::size_t i=0;i<szg && isSame;i++)
280 isSame=_groups[i].isSameAs(oldGrps[i]);
281 for(std::size_t i=0;i<szf && isSame;i++)
282 isSame=_fams[i].isSameAs(oldFams[i]);
285 for(std::size_t i=0;i<szg;i++)
286 _groups[i].cpyStatusFrom(oldGrps[i]);
287 for(std::size_t i=0;i<szf;i++)
288 _fams[i].cpyStatusFrom(oldFams[i]);
293 int vtkExtractGroup::vtkExtractGroupInternal::getNumberOfEntries() const
295 std::size_t sz0(_groups.size()),sz1(_fams.size());
296 return (int)(sz0+sz1);
299 const char *vtkExtractGroup::vtkExtractGroupInternal::getKeyOfEntry(int i) const
301 int sz0((int)_groups.size());
303 return _groups[i].getKeyOfEntry();
305 return _fams[i-sz0].getKeyOfEntry();
308 bool vtkExtractGroup::vtkExtractGroupInternal::getStatusOfEntryStr(const char *entry) const
310 const ExtractGroupStatus& elt(getEntry(entry));
311 return elt.getStatus();
314 void vtkExtractGroup::vtkExtractGroupInternal::setStatusOfEntryStr(const char *entry, bool status)
316 _selection.emplace_back(entry,status);
319 const ExtractGroupStatus& vtkExtractGroup::vtkExtractGroupInternal::getEntry(const char *entry) const
321 std::string entryCpp(entry);
322 for(std::vector<ExtractGroupGrp>::const_iterator it0=_groups.begin();it0!=_groups.end();it0++)
323 if(entryCpp==(*it0).getKeyOfEntry())
325 for(std::vector<ExtractGroupFam>::const_iterator it0=_fams.begin();it0!=_fams.end();it0++)
326 if(entryCpp==(*it0).getKeyOfEntry())
328 std::ostringstream oss; oss << "vtkExtractGroupInternal::getEntry : no such entry \"" << entry << "\"!";
329 throw INTERP_KERNEL::Exception(oss.str().c_str());
332 ExtractGroupStatus& vtkExtractGroup::vtkExtractGroupInternal::getEntry(const char *entry)
334 std::string entryCpp(entry);
335 for(std::vector<ExtractGroupGrp>::iterator it0=_groups.begin();it0!=_groups.end();it0++)
336 if(entryCpp==(*it0).getKeyOfEntry())
338 for(std::vector<ExtractGroupFam>::iterator it0=_fams.begin();it0!=_fams.end();it0++)
339 if(entryCpp==(*it0).getKeyOfEntry())
341 std::ostringstream oss; oss << "vtkExtractGroupInternal::getEntry : no such entry \"" << entry << "\"!";
342 throw INTERP_KERNEL::Exception(oss.str().c_str());
345 void vtkExtractGroup::vtkExtractGroupInternal::printMySelf(std::ostream& os) const
347 os << "Groups :" << std::endl;
348 for(std::vector<ExtractGroupGrp>::const_iterator it0=_groups.begin();it0!=_groups.end();it0++)
349 (*it0).printMySelf(os);
350 os << "Families :" << std::endl;
351 for(std::vector<ExtractGroupFam>::const_iterator it0=_fams.begin();it0!=_fams.end();it0++)
352 (*it0).printMySelf(os);
355 int vtkExtractGroup::vtkExtractGroupInternal::getIdOfFamily(const std::string& famName) const
357 for(std::vector<ExtractGroupFam>::const_iterator it=_fams.begin();it!=_fams.end();it++)
359 if((*it).getName()==famName)
360 return (*it).getId();
362 return std::numeric_limits<int>::max();
365 ExtractGroupFam::ExtractGroupFam(const char *name):ExtractGroupStatus(name),_id(0)
367 std::size_t pos(_name.find(MEDFileFieldRepresentationLeavesArrays::ZE_SEP));
368 std::string name0(_name.substr(0,pos)),name1(_name.substr(pos+strlen(MEDFileFieldRepresentationLeavesArrays::ZE_SEP)));
369 std::istringstream iss(name1);
371 std::ostringstream oss; oss << START << name; _ze_key_name=oss.str(); _name=name0;
374 bool ExtractGroupFam::isSameAs(const ExtractGroupFam& other) const
376 bool ret(ExtractGroupStatus::isSameAs(other));
378 return _id==other._id;
383 void ExtractGroupFam::printMySelf(std::ostream& os) const
385 os << " -" << _ze_key_name << " famName : \"" << _name << "\" id : " << _id << " (";
390 os << ")" << std::endl;
393 void ExtractGroupFam::fillIdsToKeep(std::set<int>& s) const
398 std::set<int> vtkExtractGroup::vtkExtractGroupInternal::getIdsToKeep() const
400 for(auto it: _selection)
402 const ExtractGroupStatus& elt(getEntry(it.first.c_str()));
403 elt.setStatus(it.second);
405 std::map<std::string,int> m(this->computeFamStrIdMap());
407 for(std::vector<ExtractGroupGrp>::const_iterator it0=_groups.begin();it0!=_groups.end();it0++)
409 if((*it0).getStatus())
411 const std::vector<std::string>& fams((*it0).getFamiliesLyingOn());
412 for(std::vector<std::string>::const_iterator it1=fams.begin();it1!=fams.end();it1++)
414 std::map<std::string,int>::iterator it2(m.find((*it1)));
416 s.insert((*it2).second);
420 for(std::vector<ExtractGroupFam>::const_iterator it0=_fams.begin();it0!=_fams.end();it0++)
421 if((*it0).getStatus())
422 (*it0).fillIdsToKeep(s);
426 void vtkExtractGroup::vtkExtractGroupInternal::clearSelection() const
429 for(auto it : _groups)
435 std::map<std::string,int> vtkExtractGroup::vtkExtractGroupInternal::computeFamStrIdMap() const
437 std::map<std::string,int> ret;
438 for(std::vector<ExtractGroupFam>::const_iterator it0=_fams.begin();it0!=_fams.end();it0++)
439 ret[(*it0).getName()]=(*it0).getId();
445 vtkExtractGroup::vtkExtractGroup():SIL(NULL),Internal(new vtkExtractGroupInternal),InsideOut(0)
449 vtkExtractGroup::~vtkExtractGroup()
451 delete this->Internal;
454 void vtkExtractGroup::SetInsideOut(int val)
456 if(this->InsideOut!=val)
463 int vtkExtractGroup::RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
465 // vtkUnstructuredGridAlgorithm::RequestInformation(request,inputVector,outputVector);
468 // std::cerr << "########################################## vtkExtractGroup::RequestInformation ##########################################" << std::endl;
469 // request->Print(cout);
470 vtkInformation *outInfo(outputVector->GetInformationObject(0));
471 vtkInformation *inputInfo(inputVector[0]->GetInformationObject(0));
472 if(!vtkExtractGroup::vtkExtractGroupInternal::IsInformationOK(inputInfo))
474 vtkErrorMacro("No SIL Data available ! The source of this filter must be MEDReader !");
478 this->SetSIL(vtkMutableDirectedGraph::SafeDownCast(inputInfo->Get(vtkMEDReader::META_DATA())));
479 this->Internal->loadFrom(this->SIL);
480 //this->Internal->printMySelf(std::cerr);
482 catch(INTERP_KERNEL::Exception& e)
484 std::cerr << "Exception has been thrown in vtkExtractGroup::RequestInformation : " << e.what() << std::endl;
491 * Do not use vtkCxxSetObjectMacro macro because input mdg comes from an already managed in the pipeline just a ref on it.
493 void vtkExtractGroup::SetSIL(vtkMutableDirectedGraph *mdg)
500 template<class CellPointExtractor>
501 vtkDataSet *FilterFamilies(vtkSmartPointer<vtkThreshold>& thres,
502 vtkDataSet *input, const std::set<int>& idsToKeep, bool insideOut, const char *arrNameOfFamilyField,
503 const char *associationForThreshold, bool& catchAll, bool& catchSmth)
505 const int VTK_DATA_ARRAY_DELETE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_DELETE;
506 const char ZE_SELECTION_ARR_NAME[]="@@ZeSelection@@";
507 vtkDataSet *output(input->NewInstance());
508 output->ShallowCopy(input);
509 thres->SetInputData(output);
510 vtkDataSetAttributes *dscIn(input->GetCellData()),*dscIn2(input->GetPointData());
511 vtkDataSetAttributes *dscOut(output->GetCellData()),*dscOut2(output->GetPointData());
513 double vMin(insideOut==0?1.:0.),vMax(insideOut==0?2.:1.);
514 thres->ThresholdBetween(vMin,vMax);
517 CellPointExtractor cpe2(input);
518 vtkDataArray *da(cpe2.Get()->GetScalars(arrNameOfFamilyField));
521 std::string daName(da->GetName());
522 typedef MEDFileVTKTraits<mcIdType>::VtkType vtkMCIdTypeArray;
523 vtkMCIdTypeArray *dai(vtkMCIdTypeArray::SafeDownCast(da));
524 if(daName!=arrNameOfFamilyField || !dai)
527 int nbOfTuples(dai->GetNumberOfTuples());
528 vtkCharArray *zeSelection(vtkCharArray::New());
529 zeSelection->SetName(ZE_SELECTION_ARR_NAME);
530 zeSelection->SetNumberOfComponents(1);
531 char *pt(new char[nbOfTuples]);
532 zeSelection->SetArray(pt,nbOfTuples,0,VTK_DATA_ARRAY_DELETE);
533 const mcIdType *inPtr(dai->GetPointer(0));
534 std::fill(pt,pt+nbOfTuples,0);
535 catchAll=true; catchSmth=false;
536 std::vector<bool> pt2(nbOfTuples,false);
537 for(std::set<int>::const_iterator it=idsToKeep.begin();it!=idsToKeep.end();it++)
539 bool catchFid(false);
540 for(int i=0;i<nbOfTuples;i++)
542 { pt2[i]=true; catchFid=true; }
548 for(int ii=0;ii<nbOfTuples;ii++)
551 CellPointExtractor cpe3(output);
552 int idx(cpe3.Get()->AddArray(zeSelection));
553 cpe3.Get()->SetActiveAttribute(idx,vtkDataSetAttributes::SCALARS);
554 cpe3.Get()->CopyScalarsOff();
555 zeSelection->Delete();
557 thres->SetInputArrayToProcess(idx,0,0,associationForThreshold,ZE_SELECTION_ARR_NAME);
559 vtkUnstructuredGrid *zeComputedOutput(thres->GetOutput());
560 CellPointExtractor cpe(zeComputedOutput);
561 cpe.Get()->RemoveArray(idx);
563 zeComputedOutput->Register(0);
564 return zeComputedOutput;
570 CellExtractor(vtkDataSet *ds):_ds(ds) { }
571 vtkDataSetAttributes *Get() { return _ds->GetCellData(); }
579 PointExtractor(vtkDataSet *ds):_ds(ds) { }
580 vtkDataSetAttributes *Get() { return _ds->GetPointData(); }
584 int vtkExtractGroup::RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
588 // std::cerr << "########################################## vtkExtractGroup::RequestData ##########################################" << std::endl;
589 // request->Print(cout);
590 vtkInformation* inputInfo=inputVector[0]->GetInformationObject(0);
591 vtkMultiBlockDataSet *inputMB(vtkMultiBlockDataSet::SafeDownCast(inputInfo->Get(vtkDataObject::DATA_OBJECT())));
592 if(inputMB->GetNumberOfBlocks()!=1)
594 std::ostringstream oss; oss << "vtkExtractGroup::RequestData : input has not the right number of parts ! Expected 1 !";
595 if(this->HasObserver("ErrorEvent") )
596 this->InvokeEvent("ErrorEvent",const_cast<char *>(oss.str().c_str()));
598 vtkOutputWindowDisplayErrorText(const_cast<char *>(oss.str().c_str()));
599 vtkObject::BreakOnError();
602 vtkDataSet *input(vtkDataSet::SafeDownCast(inputMB->GetBlock(0)));
603 vtkInformation *info(input->GetInformation());
604 vtkInformation *outInfo(outputVector->GetInformationObject(0));
605 vtkMultiBlockDataSet *output(vtkMultiBlockDataSet::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT())));
606 std::set<int> idsToKeep(this->Internal->getIdsToKeep());
607 this->Internal->clearSelection();
608 // first shrink the input
609 bool catchAll,catchSmth;
610 vtkSmartPointer<vtkThreshold> thres1(vtkSmartPointer<vtkThreshold>::New()),thres2(vtkSmartPointer<vtkThreshold>::New());
611 vtkDataSet *tryOnCell(FilterFamilies<CellExtractor>(thres1,input,idsToKeep,this->InsideOut,
612 MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_CELL_NAME,"vtkDataObject::FIELD_ASSOCIATION_CELLS",catchAll,catchSmth));
617 output->SetBlock(0,tryOnCell);
618 tryOnCell->Delete();//
625 vtkDataSet *tryOnNode(FilterFamilies<PointExtractor>(thres2,input,idsToKeep,this->InsideOut,
626 MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME,"vtkDataObject::FIELD_ASSOCIATION_POINTS",catchAll,catchSmth));
627 if(tryOnNode && catchSmth)
629 output->SetBlock(0,tryOnCell);
630 output->SetBlock(1,tryOnNode);
639 output->SetBlock(0,tryOnCell);
646 vtkDataSet *tryOnNode(FilterFamilies<PointExtractor>(thres1,input,idsToKeep,this->InsideOut,
647 MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME,"vtkDataObject::FIELD_ASSOCIATION_POINTS",catchAll,catchSmth));
651 output->SetBlock(0,tryOnNode);
657 output->SetBlock(0,tryOnNode);
666 vtkDataSet *tryOnNode(FilterFamilies<PointExtractor>(thres1,input,idsToKeep,this->InsideOut,
667 MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME,"vtkDataObject::FIELD_ASSOCIATION_POINTS",catchAll,catchSmth));
670 output->ShallowCopy(tryOnNode);
671 tryOnNode->Delete();//
676 std::ostringstream oss; oss << "vtkExtractGroup::RequestData : The integer array with name \""<< MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_CELL_NAME;
677 oss << "\" or \"" << MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME << "\" does not exist ! The extraction of group and/or family is not possible !";
678 if(this->HasObserver("ErrorEvent") )
679 this->InvokeEvent("ErrorEvent",const_cast<char *>(oss.str().c_str()));
681 vtkOutputWindowDisplayErrorText(const_cast<char *>(oss.str().c_str()));
682 vtkObject::BreakOnError();
687 catch(INTERP_KERNEL::Exception& e)
689 std::cerr << "Exception has been thrown in vtkExtractGroup::RequestData : " << e.what() << std::endl;
694 int vtkExtractGroup::GetSILUpdateStamp()
696 return (int)this->SILTime;
699 void vtkExtractGroup::PrintSelf(ostream& os, vtkIndent indent)
701 this->Superclass::PrintSelf(os, indent);
704 int vtkExtractGroup::GetNumberOfGroupsFlagsArrays()
706 int ret(this->Internal->getNumberOfEntries());
707 //std::cerr << "vtkExtractGroup::GetNumberOfFieldsTreeArrays() -> " << ret << std::endl;
711 const char *vtkExtractGroup::GetGroupsFlagsArrayName(int index)
713 const char *ret(this->Internal->getKeyOfEntry(index));
714 // std::cerr << "vtkExtractGroup::GetFieldsTreeArrayName(" << index << ") -> " << ret << std::endl;
718 int vtkExtractGroup::GetGroupsFlagsArrayStatus(const char *name)
720 int ret((int)this->Internal->getStatusOfEntryStr(name));
721 // std::cerr << "vtkExtractGroup::GetGroupsFlagsArrayStatus(" << name << ") -> " << ret << std::endl;
725 void vtkExtractGroup::SetGroupsFlagsStatus(const char *name, int status)
727 //std::cerr << "vtkExtractGroup::SetFieldsStatus(" << name << "," << status << ")" << std::endl;
728 this->Internal->setStatusOfEntryStr(name,(bool)status);
730 //this->Internal->printMySelf(std::cerr);
733 const char *vtkExtractGroup::GetMeshName()
735 return this->Internal->getMeshName();