]> SALOME platform Git repositories - tools/paravisaddons_common.git/blob - src/ExtractThreeD/plugin/ExtractThreeDModule/vtkExtractThreeD.cxx
Salome HOME
nabil tuleap23995
[tools/paravisaddons_common.git] / src / ExtractThreeD / plugin / ExtractThreeDModule / vtkExtractThreeD.cxx
1 // Copyright (C) 2021  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 "vtkExtractThreeD.h"
22
23 #include <NormalizedGeometricTypes>
24 #include <CellModel.hxx>
25 #include <InterpKernelException.hxx>
26 #include <MEDFileFieldOverView.hxx>
27
28 #include <vtkAdjacentVertexIterator.h>
29 #include <vtkAlgorithmOutput.h>
30 #include <vtkCellData.h>
31 #include <vtkCharArray.h>
32 #include <vtkDataArraySelection.h>
33 #include <vtkAOSDataArrayTemplate.h>
34 #include <vtkDataObjectTreeIterator.h>
35 #include <vtkDataSet.h>
36 #include <vtkDataSetAttributes.h>
37 #include <vtkDemandDrivenPipeline.h>
38 #include <vtkDoubleArray.h>
39 #include <vtkExecutive.h>
40 #include <vtkInEdgeIterator.h>
41 #include <vtkInformation.h>
42 #include <vtkInformationDataObjectKey.h>
43 #include <vtkInformationStringKey.h>
44 #include <vtkInformationVector.h>
45 #include <vtkIntArray.h>
46 #include <vtkMultiBlockDataSet.h>
47 #include <vtkMutableDirectedGraph.h>
48 #include <vtkObjectFactory.h>
49 #include <vtkPointData.h>
50 #include <vtkStreamingDemandDrivenPipeline.h>
51 #include <vtkStringArray.h>
52 #include <vtkThreshold.h>
53 #include <vtkTimeStamp.h>
54 #include <vtkUnsignedCharArray.h>
55 #include <vtkUnstructuredGrid.h>
56 #include <vtkVariantArray.h>
57
58 #include <deque>
59 #include <map>
60 #include <sstream>
61
62 vtkStandardNewMacro(vtkExtractThreeD);
63
64 ///////////////////
65
66 class vtkExtractThreeD::vtkExtractThreeDInternal
67 {
68 public:
69   vtkExtractThreeDInternal() {}
70   std::vector<int> getIdsToKeep() const { return _types; }
71   void setIdsToKeep(const std::vector<int> &types) { _types = types; }
72
73 private:
74   std::vector<int> _types;
75 };
76
77 vtkExtractThreeD::vtkExtractThreeD() : Internal(new vtkExtractThreeDInternal)
78 {
79 }
80
81 vtkExtractThreeD::~vtkExtractThreeD()
82 {
83   delete this->Internal;
84 }
85
86 int vtkExtractThreeD::RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
87 {
88   try
89   {
90     //std::cerr << "########################################## vtkExtractThreeD::RequestInformation ##########################################" << std::endl;
91     vtkInformation *outInfo(outputVector->GetInformationObject(0));
92     vtkInformation *inputInfo(inputVector[0]->GetInformationObject(0));
93     vtkDataSet *input(0);
94     {
95       vtkDataObject *inp(inputInfo->Get(vtkDataObject::DATA_OBJECT()));
96       if (vtkDataSet::SafeDownCast(inp))
97         input = vtkDataSet::SafeDownCast(inp);
98       else
99       {
100         vtkMultiBlockDataSet *inputTmp(vtkMultiBlockDataSet::SafeDownCast(inp));
101         if (inputTmp)
102         {
103           if (inputTmp->GetNumberOfBlocks() != 1)
104           {
105             vtkDebugMacro("vtkExtractThreeD::RequestInformation : input vtkMultiBlockDataSet must contain exactly 1 block !");
106             return 0;
107           }
108           vtkDataSet *blk0(vtkDataSet::SafeDownCast(inputTmp->GetBlock(0)));
109           if (!blk0)
110           {
111             vtkDebugMacro("vtkExtractThreeD::RequestInformation : the single block in input vtkMultiBlockDataSet must be a vtkDataSet instance !");
112             return 0;
113           }
114           input = blk0;
115         }
116         else
117         {
118           vtkDebugMacro("vtkExtractThreeD::RequestInformation : supported input are vtkDataSet or vtkMultiBlockDataSet !");
119           return 0;
120         }
121       }
122     }
123     {
124       vtkIdType nbOfCells(input->GetNumberOfCells());
125       std::map<int, INTERP_KERNEL::NormalizedCellType> m;
126       std::set<int> typesToKeep;
127       for (vtkIdType cellId = 0; cellId < nbOfCells; cellId++)
128       {
129         int vtkCt(input->GetCellType(cellId));
130         const std::map<int, INTERP_KERNEL::NormalizedCellType>::const_iterator it(m.find(vtkCt));
131         if (it == m.end())
132         {
133           const unsigned char *pos(std::find(MEDCoupling::MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE, MEDCoupling::MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE + MEDCoupling::MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE_LGTH, vtkCt));
134           if (pos == MEDCoupling::MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE + MEDCoupling::MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE_LGTH)
135           {
136             vtkDebugMacro("vtkExtractThreeD::RequestInformation : cell #" << cellId << " has unrecognized type !");
137             return 0;
138           }
139           INTERP_KERNEL::NormalizedCellType mcCtype((INTERP_KERNEL::NormalizedCellType)std::distance(MEDCoupling::MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE, pos));
140           const INTERP_KERNEL::CellModel &cm(INTERP_KERNEL::CellModel::GetCellModel(mcCtype));
141           if (cm.getDimension() == 3)
142             typesToKeep.insert(vtkCt);
143         }
144       }
145       std::vector<int> typesToKeep2(typesToKeep.begin(), typesToKeep.end());
146       this->Internal->setIdsToKeep(typesToKeep2);
147     }
148   }
149   catch (INTERP_KERNEL::Exception &e)
150   {
151     std::cerr << "Exception has been thrown in vtkExtractThreeD::RequestInformation : " << e.what() << std::endl;
152     return 0;
153   }
154   return 1;
155 }
156
157 vtkDataSet *FilterFamilies(vtkDataSet *input, const std::vector<int> &idsToKeep)
158 {
159   const int VTK_DATA_ARRAY_DELETE = vtkAOSDataArrayTemplate<double>::VTK_DATA_ARRAY_DELETE;
160   const char ZE_SELECTION_ARR_NAME[] = "@@ZeSelection@@";
161   vtkDataSet *output(input->NewInstance());
162   output->ShallowCopy(input);
163   vtkSmartPointer<vtkThreshold> thres(vtkSmartPointer<vtkThreshold>::New());
164   thres->SetInputData(output);
165   vtkDataSetAttributes *dscIn(input->GetCellData()), *dscIn2(input->GetPointData());
166   vtkDataSetAttributes *dscOut(output->GetCellData()), *dscOut2(output->GetPointData());
167   //
168   double vMin(1.), vMax(2.);
169   thres->ThresholdBetween(vMin, vMax);
170   // OK for the output
171   vtkIdType nbOfCells(input->GetNumberOfCells());
172   vtkCharArray *zeSelection(vtkCharArray::New());
173   zeSelection->SetName(ZE_SELECTION_ARR_NAME);
174   zeSelection->SetNumberOfComponents(1);
175   char *pt(new char[nbOfCells]);
176   zeSelection->SetArray(pt, nbOfCells, 0, VTK_DATA_ARRAY_DELETE);
177   std::fill(pt, pt + nbOfCells, 0);
178   std::vector<bool> pt2(nbOfCells, false);
179   for (std::vector<int>::const_iterator it = idsToKeep.begin(); it != idsToKeep.end(); it++)
180   {
181     for (vtkIdType ii = 0; ii < nbOfCells; ii++)
182     {
183       if (input->GetCellType(ii) == *it)
184         pt2[ii] = true;
185     }
186   }
187   for (int ii = 0; ii < nbOfCells; ii++)
188     if (pt2[ii])
189       pt[ii] = 2;
190   int idx(output->GetCellData()->AddArray(zeSelection));
191   output->GetCellData()->SetActiveAttribute(idx, vtkDataSetAttributes::SCALARS);
192   output->GetCellData()->CopyScalarsOff();
193   zeSelection->Delete();
194   //
195   thres->SetInputArrayToProcess(idx, 0, 0, "vtkDataObject::FIELD_ASSOCIATION_CELLS", ZE_SELECTION_ARR_NAME);
196   thres->Update();
197   vtkUnstructuredGrid *zeComputedOutput(thres->GetOutput());
198   zeComputedOutput->GetCellData()->RemoveArray(idx);
199   output->Delete();
200   zeComputedOutput->Register(0);
201   return zeComputedOutput;
202 }
203
204 int vtkExtractThreeD::RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
205 {
206   try
207   {
208     //std::cerr << "########################################## vtkExtractThreeD::RequestData        ##########################################" << std::endl;
209     vtkInformation *inputInfo = inputVector[0]->GetInformationObject(0);
210     vtkDataSet *input(vtkDataSet::SafeDownCast(inputInfo->Get(vtkDataObject::DATA_OBJECT())));
211     vtkInformation *info(input->GetInformation());
212     vtkInformation *outInfo(outputVector->GetInformationObject(0));
213     vtkDataSet *output(vtkDataSet::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT())));
214     std::vector<int> idsToKeep(this->Internal->getIdsToKeep());
215     vtkDataSet *tryOnCell(FilterFamilies(input, idsToKeep));
216     // first shrink the input
217     output->ShallowCopy(tryOnCell);
218     tryOnCell->Delete();
219   }
220   catch (INTERP_KERNEL::Exception &e)
221   {
222     std::cerr << "Exception has been thrown in vtkExtractThreeD::RequestData : " << e.what() << std::endl;
223     return 0;
224   }
225   return 1;
226 }
227
228 void vtkExtractThreeD::PrintSelf(ostream &os, vtkIndent indent)
229 {
230   this->Superclass::PrintSelf(os, indent);
231 }