From: Anthony Geay Date: Wed, 30 Dec 2015 15:45:30 +0000 (+0100) Subject: ExtractGroup is now a multiblock filter and not an unstructured like MEDReader. For... X-Git-Tag: V7_8_0a1~9 X-Git-Url: http://git.salome-platform.org/gitweb/?a=commitdiff_plain;h=b705d58ad8589dde0c4d7fc673f1fb723c3f2ea5;p=modules%2Fparavis.git ExtractGroup is now a multiblock filter and not an unstructured like MEDReader. For users mixing cells and nodes fetch the two parts are in different blocks. --- diff --git a/src/Plugins/MEDReader/IO/vtkExtractGroup.cxx b/src/Plugins/MEDReader/IO/vtkExtractGroup.cxx index 06ceb7a9..47c5e1ce 100644 --- a/src/Plugins/MEDReader/IO/vtkExtractGroup.cxx +++ b/src/Plugins/MEDReader/IO/vtkExtractGroup.cxx @@ -452,7 +452,7 @@ int vtkExtractGroup::RequestInformation(vtkInformation *request, vtkInformationV if(!vtkExtractGroup::vtkExtractGroupInternal::IsInformationOK(inputInfo)) { vtkErrorMacro("No SIL Data available ! The source of this filter must be MEDReader !"); - return 0; + return 0; } this->SetSIL(vtkMutableDirectedGraph::SafeDownCast(inputInfo->Get(vtkMEDReader::META_DATA()))); @@ -568,95 +568,100 @@ int vtkExtractGroup::RequestData(vtkInformation *request, vtkInformationVector * // std::cerr << "########################################## vtkExtractGroup::RequestData ##########################################" << std::endl; // request->Print(cout); vtkInformation* inputInfo=inputVector[0]->GetInformationObject(0); - vtkDataSet *input(vtkDataSet::SafeDownCast(inputInfo->Get(vtkDataObject::DATA_OBJECT()))); + vtkMultiBlockDataSet *inputMB(vtkMultiBlockDataSet::SafeDownCast(inputInfo->Get(vtkDataObject::DATA_OBJECT()))); + if(inputMB->GetNumberOfBlocks()!=1) + { + std::ostringstream oss; oss << "vtkExtractGroup::RequestData : input has not the right number of parts ! Expected 1 !"; + if(this->HasObserver("ErrorEvent") ) + this->InvokeEvent("ErrorEvent",const_cast(oss.str().c_str())); + else + vtkOutputWindowDisplayErrorText(const_cast(oss.str().c_str())); + vtkObject::BreakOnError(); + return 0; + } + vtkDataSet *input(vtkDataSet::SafeDownCast(inputMB->GetBlock(0))); vtkInformation *info(input->GetInformation()); vtkInformation *outInfo(outputVector->GetInformationObject(0)); - vtkDataSet *output(vtkDataSet::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT()))); + vtkMultiBlockDataSet *output(vtkMultiBlockDataSet::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT()))); std::set idsToKeep(this->Internal->getIdsToKeep()); // first shrink the input bool catchAll,catchSmth; vtkSmartPointer thres1(vtkSmartPointer::New()),thres2(vtkSmartPointer::New()); vtkDataSet *tryOnCell(FilterFamilies(thres1,input,idsToKeep,this->InsideOut, - MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_CELL_NAME,"vtkDataObject::FIELD_ASSOCIATION_CELLS",catchAll,catchSmth)); + MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_CELL_NAME,"vtkDataObject::FIELD_ASSOCIATION_CELLS",catchAll,catchSmth)); if(tryOnCell) - { - if(catchAll) - { - output->ShallowCopy(tryOnCell); - tryOnCell->Delete();// - return 1; - } - else - { - if(catchSmth) - { - vtkDataSet *tryOnNode(FilterFamilies(thres2,input,idsToKeep,this->InsideOut, - MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME,"vtkDataObject::FIELD_ASSOCIATION_POINTS",catchAll,catchSmth)); - if(tryOnNode && catchSmth) - { - vtkSmartPointer mb(vtkSmartPointer::New()); - vtkSmartPointer cd(vtkSmartPointer::New()); - mb->AddInputConnection(thres1->GetOutputPort()); - mb->AddInputConnection(thres2->GetOutputPort()); - cd->SetInputConnection(mb->GetOutputPort()); - cd->SetMergePoints(0); - cd->Update(); - output->ShallowCopy(cd->GetOutput()); - tryOnCell->Delete(); - tryOnNode->Delete(); - return 1; - } - else - { - if(tryOnNode) - tryOnNode->Delete(); - output->ShallowCopy(tryOnCell); - tryOnCell->Delete(); - return 1; - } - } - else - { - vtkDataSet *tryOnNode(FilterFamilies(thres1,input,idsToKeep,this->InsideOut, - MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME,"vtkDataObject::FIELD_ASSOCIATION_POINTS",catchAll,catchSmth)); - if(tryOnNode) - { - tryOnCell->Delete(); - output->ShallowCopy(tryOnNode); - tryOnNode->Delete(); - return 1; - } - else - { - output->ShallowCopy(tryOnCell); - tryOnCell->Delete(); - return 0; - } - } - } - } + { + if(catchAll) + { + output->SetBlock(0,tryOnCell); + tryOnCell->Delete();// + return 1; + } + else + { + if(catchSmth) + { + vtkDataSet *tryOnNode(FilterFamilies(thres2,input,idsToKeep,this->InsideOut, + MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME,"vtkDataObject::FIELD_ASSOCIATION_POINTS",catchAll,catchSmth)); + if(tryOnNode && catchSmth) + { + output->SetBlock(0,tryOnCell); + output->SetBlock(1,tryOnNode); + tryOnCell->Delete(); + tryOnNode->Delete(); + return 1; + } + else + { + if(tryOnNode) + tryOnNode->Delete(); + output->SetBlock(0,tryOnCell); + tryOnCell->Delete(); + return 1; + } + } + else + { + vtkDataSet *tryOnNode(FilterFamilies(thres1,input,idsToKeep,this->InsideOut, + MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME,"vtkDataObject::FIELD_ASSOCIATION_POINTS",catchAll,catchSmth)); + if(tryOnNode) + { + tryOnCell->Delete(); + output->SetBlock(0,tryOnNode); + tryOnNode->Delete(); + return 1; + } + else + { + output->SetBlock(0,tryOnNode); + tryOnCell->Delete(); + return 0; + } + } + } + } else - { - vtkDataSet *tryOnNode(FilterFamilies(thres1,input,idsToKeep,this->InsideOut, - MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME,"vtkDataObject::FIELD_ASSOCIATION_POINTS",catchAll,catchSmth)); - if(tryOnNode) - { - output->ShallowCopy(tryOnNode); - tryOnNode->Delete();// - return 1; - } - else - { - std::ostringstream oss; oss << "vtkExtractGroup::RequestData : The integer array with name \""<< MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_CELL_NAME; - oss << "\" or \"" << MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME << "\" does not exist ! The extraction of group and/or family is not possible !"; - if(this->HasObserver("ErrorEvent") ) - this->InvokeEvent("ErrorEvent",const_cast(oss.str().c_str())); - else - vtkOutputWindowDisplayErrorText(const_cast(oss.str().c_str())); - vtkObject::BreakOnError(); - return 0; - } - } + { + vtkDataSet *tryOnNode(FilterFamilies(thres1,input,idsToKeep,this->InsideOut, + MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME,"vtkDataObject::FIELD_ASSOCIATION_POINTS",catchAll,catchSmth)); + if(tryOnNode) + { + output->ShallowCopy(tryOnNode); + tryOnNode->Delete();// + return 1; + } + else + { + std::ostringstream oss; oss << "vtkExtractGroup::RequestData : The integer array with name \""<< MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_CELL_NAME; + oss << "\" or \"" << MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME << "\" does not exist ! The extraction of group and/or family is not possible !"; + if(this->HasObserver("ErrorEvent") ) + this->InvokeEvent("ErrorEvent",const_cast(oss.str().c_str())); + else + vtkOutputWindowDisplayErrorText(const_cast(oss.str().c_str())); + vtkObject::BreakOnError(); + return 0; + } + } } catch(INTERP_KERNEL::Exception& e) { diff --git a/src/Plugins/MEDReader/IO/vtkExtractGroup.h b/src/Plugins/MEDReader/IO/vtkExtractGroup.h index d3cca44b..b7b06938 100644 --- a/src/Plugins/MEDReader/IO/vtkExtractGroup.h +++ b/src/Plugins/MEDReader/IO/vtkExtractGroup.h @@ -21,15 +21,15 @@ #ifndef vtkExtractGroup_h__ #define vtkExtractGroup_h__ -#include "vtkUnstructuredGridAlgorithm.h" +#include "vtkMultiBlockDataSetAlgorithm.h" class vtkMutableDirectedGraph; -class VTK_EXPORT vtkExtractGroup: public vtkUnstructuredGridAlgorithm +class VTK_EXPORT vtkExtractGroup: public vtkMultiBlockDataSetAlgorithm { public: static vtkExtractGroup* New(); - vtkTypeMacro(vtkExtractGroup, vtkUnstructuredGridAlgorithm) + vtkTypeMacro(vtkExtractGroup, vtkMultiBlockDataSetAlgorithm) void PrintSelf(ostream& os, vtkIndent indent); virtual int GetNumberOfGroupsFlagsArrays(); const char *GetGroupsFlagsArrayName(int index); diff --git a/src/Plugins/MEDReader/Test/testMEDReader18.py b/src/Plugins/MEDReader/Test/testMEDReader18.py index 235a4042..f4f48a5d 100644 --- a/src/Plugins/MEDReader/Test/testMEDReader18.py +++ b/src/Plugins/MEDReader/Test/testMEDReader18.py @@ -49,4 +49,6 @@ ExtractGroup1 = ExtractGroup(Input=reader) ExtractGroup1.AllGroups=["GRP_grp0","GRP_grp1"] ExtractGroup1.UpdatePipelineInformation() res=servermanager.Fetch(ExtractGroup1,0) -assert(res.GetBlock(0).GetNumberOfCells()==3) +assert(res.GetNumberOfBlocks()==2) +assert(res.GetBlock(1).GetNumberOfCells()==1) +assert(res.GetBlock(0).GetNumberOfCells()==2)