Salome HOME
Merge remote branch 'origin/master'
[modules/paravis.git] / src / Plugins / MEDReader / IO / vtkExtractGroup.cxx
index 5eab0e9011aa5e435ae1036245fd9b960ca307f9..1db4ee4975a34415338d9e917fc951cafa492f14 100644 (file)
@@ -1,9 +1,9 @@
-// Copyright (C) 2010-2014  CEA/DEN, EDF R&D
+// Copyright (C) 2010-2015  CEA/DEN, EDF R&D
 //
 // This library is free software; you can redistribute it and/or
 // modify it under the terms of the GNU Lesser General Public
 // License as published by the Free Software Foundation; either
-// version 2.1 of the License.
+// version 2.1 of the License, or (at your option) any later version.
 //
 // This library is distributed in the hope that it will be useful,
 // but WITHOUT ANY WARRANTY; without even the implied warranty of
 #include "vtkDemandDrivenPipeline.h"
 #include "vtkDataObjectTreeIterator.h"
 #include "vtkThreshold.h"
+#include "vtkMultiBlockDataGroupFilter.h"
+#include "vtkCompositeDataToUnstructuredGridFilter.h"
 
 #include <map>
 #include <deque>
 
 vtkStandardNewMacro(vtkExtractGroup);
 
-vtkCxxSetObjectMacro(vtkExtractGroup, SIL, vtkMutableDirectedGraph);
-
 ///////////////////
 
 class ExtractGroupStatus
@@ -427,6 +427,7 @@ void vtkExtractGroup::SetInsideOut(int val)
 
 int vtkExtractGroup::RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
 {
+  vtkUnstructuredGridAlgorithm::RequestInformation(request,inputVector,outputVector);
   try
     {
       //std::cerr << "########################################## vtkExtractGroup::RequestInformation ##########################################" << std::endl;
@@ -461,15 +462,25 @@ int vtkExtractGroup::RequestInformation(vtkInformation *request, vtkInformationV
   return 1;
 }
 
+/*!
+ * Do not use vtkCxxSetObjectMacro macro because input mdg comes from an already managed in the pipeline just a ref on it.
+ */
+void vtkExtractGroup::SetSIL(vtkMutableDirectedGraph *mdg)
+{
+  if(this->SIL==mdg)
+    return ;
+  this->SIL=mdg;
+}
+
 template<class CellPointExtractor>
-vtkDataSet *FilterFamilies(vtkDataSet *input, const std::set<int>& idsToKeep, bool insideOut, const char *arrNameOfFamilyField,
+vtkDataSet *FilterFamilies(vtkSmartPointer<vtkThreshold>& thres,
+                           vtkDataSet *input, const std::set<int>& idsToKeep, bool insideOut, const char *arrNameOfFamilyField,
                            const char *associationForThreshold, bool& catchAll, bool& catchSmth)
 {
   static const int VTK_DATA_ARRAY_DELETE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_DELETE;
   static const char ZE_SELECTION_ARR_NAME[]="@@ZeSelection@@";
   vtkDataSet *output(input->NewInstance());
   output->ShallowCopy(input);
-  vtkSmartPointer<vtkThreshold> thres(vtkSmartPointer<vtkThreshold>::New());
   thres->SetInputData(output);
   vtkDataSetAttributes *dscIn(input->GetCellData()),*dscIn2(input->GetPointData());
   vtkDataSetAttributes *dscOut(output->GetCellData()),*dscOut2(output->GetPointData());
@@ -558,7 +569,8 @@ int vtkExtractGroup::RequestData(vtkInformation *request, vtkInformationVector *
       std::set<int> idsToKeep(this->Internal->getIdsToKeep());
       // first shrink the input
       bool catchAll,catchSmth;
-      vtkDataSet *tryOnCell(FilterFamilies<CellExtractor>(input,idsToKeep,this->InsideOut,
+      vtkSmartPointer<vtkThreshold> thres1(vtkSmartPointer<vtkThreshold>::New()),thres2(vtkSmartPointer<vtkThreshold>::New());
+      vtkDataSet *tryOnCell(FilterFamilies<CellExtractor>(thres1,input,idsToKeep,this->InsideOut,
                                                          MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_CELL_NAME,"vtkDataObject::FIELD_ASSOCIATION_CELLS",catchAll,catchSmth));
       if(tryOnCell)
        {
@@ -572,13 +584,20 @@ int vtkExtractGroup::RequestData(vtkInformation *request, vtkInformationVector *
            {
              if(catchSmth)
                {
-                 vtkDataSet *tryOnNode(FilterFamilies<PointExtractor>(tryOnCell,idsToKeep,this->InsideOut,
+                 vtkDataSet *tryOnNode(FilterFamilies<PointExtractor>(thres2,input,idsToKeep,this->InsideOut,
                                                                       MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME,"vtkDataObject::FIELD_ASSOCIATION_POINTS",catchAll,catchSmth));
                  if(tryOnNode && catchSmth)
                    {
-                     output->ShallowCopy(tryOnNode);
+                      vtkSmartPointer<vtkMultiBlockDataGroupFilter> mb(vtkSmartPointer<vtkMultiBlockDataGroupFilter>::New());
+                      vtkSmartPointer<vtkCompositeDataToUnstructuredGridFilter> cd(vtkSmartPointer<vtkCompositeDataToUnstructuredGridFilter>::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();//
+                     tryOnNode->Delete();
                      return 1;
                    }
                  else
@@ -592,7 +611,7 @@ int vtkExtractGroup::RequestData(vtkInformation *request, vtkInformationVector *
                }
              else
                {
-                 vtkDataSet *tryOnNode(FilterFamilies<PointExtractor>(input,idsToKeep,this->InsideOut,
+                 vtkDataSet *tryOnNode(FilterFamilies<PointExtractor>(thres1,input,idsToKeep,this->InsideOut,
                                                                       MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME,"vtkDataObject::FIELD_ASSOCIATION_POINTS",catchAll,catchSmth));
                  if(tryOnNode)
                    {
@@ -612,7 +631,7 @@ int vtkExtractGroup::RequestData(vtkInformation *request, vtkInformationVector *
        }
       else
        {
-         vtkDataSet *tryOnNode(FilterFamilies<PointExtractor>(tryOnCell,idsToKeep,this->InsideOut,
+         vtkDataSet *tryOnNode(FilterFamilies<PointExtractor>(thres1,input,idsToKeep,this->InsideOut,
                                                               MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME,"vtkDataObject::FIELD_ASSOCIATION_POINTS",catchAll,catchSmth));
          if(tryOnNode)
            {
@@ -674,6 +693,8 @@ int vtkExtractGroup::GetGroupsFlagsArrayStatus(const char *name)
 void vtkExtractGroup::SetGroupsFlagsStatus(const char *name, int status)
 {
   //std::cerr << "vtkExtractGroup::SetFieldsStatus(" << name << "," << status << ")" << std::endl;
+  if (GetNumberOfGroupsFlagsArrays()<1)
+    return;
   this->Internal->setStatusOfEntryStr(name,(bool)status);
   if(std::string(name)==GetGroupsFlagsArrayName(GetNumberOfGroupsFlagsArrays()-1))
      this->Modified();