Salome HOME
d79200a2d68ce97ede5811bcaa7959377198521b
[modules/paravis.git] / src / Plugins / MEDReader / IO / vtkExtractCellType.cxx
1 // Copyright (C) 2010-2016  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 "vtkExtractCellType.h"
22 #include "MEDFileFieldRepresentationTree.hxx"
23 #include "MEDFileFieldOverView.hxx"
24
25 #include "vtkAdjacentVertexIterator.h"
26 #include "vtkDataArrayTemplate.h"
27 #include "vtkIntArray.h"
28 #include "vtkCellData.h"
29 #include "vtkPointData.h"
30
31 #include "vtkStreamingDemandDrivenPipeline.h"
32 #include "vtkUnstructuredGrid.h"
33 #include  "vtkMultiBlockDataSet.h"
34
35 #include "vtkInformationStringKey.h"
36 #include "vtkAlgorithmOutput.h"
37 #include "vtkObjectFactory.h"
38 #include "vtkMutableDirectedGraph.h"
39 #include "vtkMultiBlockDataSet.h"
40 #include "vtkDataSet.h"
41 #include "vtkInformationVector.h"
42 #include "vtkInformation.h"
43 #include "vtkDataArraySelection.h"
44 #include "vtkTimeStamp.h"
45 #include "vtkInEdgeIterator.h"
46 #include "vtkInformationDataObjectKey.h"
47 #include "vtkExecutive.h"
48 #include "vtkVariantArray.h"
49 #include "vtkStringArray.h"
50 #include "vtkDoubleArray.h"
51 #include "vtkCharArray.h"
52 #include "vtkUnsignedCharArray.h"
53 #include "vtkDataSetAttributes.h"
54 #include "vtkDemandDrivenPipeline.h"
55 #include "vtkDataObjectTreeIterator.h"
56 #include "vtkThreshold.h"
57
58 #include <map>
59 #include <deque>
60
61 vtkStandardNewMacro(vtkExtractCellType);
62
63 vtkCxxSetObjectMacro(vtkExtractCellType, SIL, vtkMutableDirectedGraph);
64
65 ///////////////////
66
67 class ExtractCellTypeStatus
68 {
69 public:
70   ExtractCellTypeStatus():_status(false),_vtkt(-1),_mct(INTERP_KERNEL::NORM_ERROR) { }
71   ExtractCellTypeStatus(int vtkt, INTERP_KERNEL::NormalizedCellType mct);
72   bool isSame(int vtkt, INTERP_KERNEL::NormalizedCellType mct) const { return _vtkt==vtkt && _mct==mct; }
73   bool getStatus() const { return _status; }
74   void setStatus(bool status) const { _status=status; }
75   void cpyStatusFrom(const ExtractCellTypeStatus& other) { _status=other._status; }
76   std::string getKey() const { return _type_str; }
77   const char *getKeyOfEntry() const { return _type_str.c_str(); }
78   int getVTKCellType() const { return _vtkt; }
79   void printMySelf(std::ostream& os) const;
80   bool isSameAs(const ExtractCellTypeStatus& other) const;
81   void feedSIL(vtkMutableDirectedGraph *sil, vtkIdType root, vtkVariantArray *childEdge, std::vector<std::string>& names) const;
82 protected:
83   mutable bool _status;
84   int _vtkt;
85   INTERP_KERNEL::NormalizedCellType _mct;
86   std::string _type_str;
87 };
88
89 class vtkExtractCellType::vtkExtractCellTypeInternal
90 {
91 public:
92   vtkExtractCellTypeInternal():_ref_mtime(0) { }
93   int getNumberOfEntries() const;
94   const char *getKeyOfEntry(int i) const;
95   bool getStatusOfEntryStr(const char *entry) const;
96   void setStatusOfEntryStr(const char *entry, bool status) const;
97   void feedSIL(vtkMutableDirectedGraph *sil) const;
98   std::vector<int> getIdsToKeep() const;
99   void printMySelf(std::ostream& os) const;
100   bool setRefTime(vtkObject *input) const;
101   // non const methods
102   void loadFrom(const std::map<int,INTERP_KERNEL::NormalizedCellType>& m);
103 private:
104   const ExtractCellTypeStatus& getEntry(const char *entry) const;
105   bool checkSame(const std::map<int,INTERP_KERNEL::NormalizedCellType>& m) const;
106 private:
107   std::vector<ExtractCellTypeStatus> _types;
108   mutable unsigned long _ref_mtime;
109 };
110
111 bool vtkExtractCellType::vtkExtractCellTypeInternal::setRefTime(vtkObject *input) const
112 {
113   unsigned long mtime(input->GetMTime());
114   if(mtime>_ref_mtime)
115     {
116       _ref_mtime=mtime;
117       return true;
118     }
119   else
120     return false;
121 }
122
123 std::vector<int> vtkExtractCellType::vtkExtractCellTypeInternal::getIdsToKeep() const
124 {
125   std::vector<int> ret;
126   for(std::vector<ExtractCellTypeStatus>::const_iterator it=_types.begin();it!=_types.end();it++)
127     {
128       if((*it).getStatus())
129         ret.push_back((*it).getVTKCellType());
130     }
131   return ret;
132 }
133
134 void vtkExtractCellType::vtkExtractCellTypeInternal::feedSIL(vtkMutableDirectedGraph *sil) const
135 {
136   vtkSmartPointer<vtkVariantArray> childEdge(vtkSmartPointer<vtkVariantArray>::New());
137   childEdge->InsertNextValue(0);
138   vtkSmartPointer<vtkVariantArray> crossEdge(vtkSmartPointer<vtkVariantArray>::New());
139   crossEdge->InsertNextValue(1);
140   // CrossEdge is an edge linking hierarchies.
141   vtkUnsignedCharArray* crossEdgesArray=vtkUnsignedCharArray::New();
142   crossEdgesArray->SetName("CrossEdges");
143   sil->GetEdgeData()->AddArray(crossEdgesArray);
144   crossEdgesArray->Delete();
145   std::vector<std::string> names;
146   // Add global fields root
147   vtkIdType root(sil->AddVertex());
148   names.push_back("CellTypesTree");
149   //
150   for(std::vector<ExtractCellTypeStatus>::const_iterator it=_types.begin();it!=_types.end();it++)
151     {
152       (*it).feedSIL(sil,root,childEdge,names);
153     }
154   // This array is used to assign names to nodes.
155   vtkStringArray *namesArray(vtkStringArray::New());
156   namesArray->SetName("Names");
157   namesArray->SetNumberOfTuples(sil->GetNumberOfVertices());
158   sil->GetVertexData()->AddArray(namesArray);
159   namesArray->Delete();
160   std::vector<std::string>::const_iterator iter;
161   vtkIdType cc;
162   for(cc=0, iter=names.begin(); iter!=names.end(); ++iter, ++cc)
163     namesArray->SetValue(cc,(*iter).c_str());
164 }
165
166 void vtkExtractCellType::vtkExtractCellTypeInternal::loadFrom(const std::map<int,INTERP_KERNEL::NormalizedCellType>& m)
167 {
168   if(checkSame(m))
169     return;
170   //
171   std::size_t sz(m.size()),ii(0);
172   _types.resize(sz);
173   for(std::map<int,INTERP_KERNEL::NormalizedCellType>::const_iterator it=m.begin();it!=m.end();it++,ii++)
174     {
175       ExtractCellTypeStatus elt((*it).first,(*it).second);
176       _types[ii]=elt;
177     }
178 }
179
180 int vtkExtractCellType::vtkExtractCellTypeInternal::getNumberOfEntries() const
181 {
182   return (int) _types.size();
183 }
184
185 const char *vtkExtractCellType::vtkExtractCellTypeInternal::getKeyOfEntry(int i) const
186 {
187   return _types[i].getKeyOfEntry();
188 }
189
190 bool vtkExtractCellType::vtkExtractCellTypeInternal::checkSame(const std::map<int,INTERP_KERNEL::NormalizedCellType>& m) const
191 {
192   std::size_t sz(m.size());
193   if(sz!=_types.size())
194     return false;
195   bool ret(true);
196   std::map<int,INTERP_KERNEL::NormalizedCellType>::const_iterator it(m.begin());
197   for(std::size_t i=0;i<sz && ret;i++,it++)
198     ret=_types[i].isSame((*it).first,(*it).second);
199   return ret;
200 }
201
202 const ExtractCellTypeStatus& vtkExtractCellType::vtkExtractCellTypeInternal::getEntry(const char *entry) const
203 {
204   std::string entryCpp(entry);
205   for(std::vector<ExtractCellTypeStatus>::const_iterator it0=_types.begin();it0!=_types.end();it0++)
206     if(entryCpp==(*it0).getKey())
207       return *it0;
208   std::ostringstream oss; oss << "vtkExtractCellTypeInternal::getEntry : no such entry \"" << entry << "\"!";
209   throw INTERP_KERNEL::Exception(oss.str().c_str());
210 }
211
212 bool vtkExtractCellType::vtkExtractCellTypeInternal::getStatusOfEntryStr(const char *entry) const
213 {
214   try 
215     {
216       const ExtractCellTypeStatus& elt(getEntry(entry));
217       return elt.getStatus();
218     }  
219   catch (INTERP_KERNEL::Exception e)
220     {      
221       //std::cerr << vtkDebugMacro"Exception has been thrown in vtkExtractCellType::vtkExtractCellTypeInternal::getStatusOfEntryStr : " << e.what() << std::endl;
222       return false;
223     }
224 }
225
226 void vtkExtractCellType::vtkExtractCellTypeInternal::setStatusOfEntryStr(const char *entry, bool status) const
227 {
228   try 
229     {
230       const ExtractCellTypeStatus& elt(getEntry(entry));
231       elt.setStatus(status);
232     }
233   catch (INTERP_KERNEL::Exception e)
234     {      
235       //std::cerr << "Exception has been thrown in vtkExtractCellType::vtkExtractCellTypeInternal::setStatusOfEntryStr : " << e.what() << std::endl;
236     }
237 }
238
239 void vtkExtractCellType::vtkExtractCellTypeInternal::printMySelf(std::ostream& os) const
240 {
241   for(std::vector<ExtractCellTypeStatus>::const_iterator it0=_types.begin();it0!=_types.end();it0++)
242     (*it0).printMySelf(os);
243 }
244
245 ExtractCellTypeStatus::ExtractCellTypeStatus(int vtkt, INTERP_KERNEL::NormalizedCellType mct):_status(false),_vtkt(vtkt),_mct(mct)
246 {
247   std::string name(INTERP_KERNEL::CellModel::GetCellModel(mct).getRepr());
248   _type_str=name.substr(5);//skip "NORM_"
249 }
250
251 void ExtractCellTypeStatus::printMySelf(std::ostream& os) const
252 {
253   os << "      -" << _type_str << "(";
254   if(_status)
255     os << "X";
256   else
257     os << " ";
258   os << ")" << std::endl;
259 }
260
261 bool ExtractCellTypeStatus::isSameAs(const ExtractCellTypeStatus& other) const
262 {
263   return _vtkt==other._vtkt && _mct==other._mct;
264 }
265
266 void ExtractCellTypeStatus::feedSIL(vtkMutableDirectedGraph *sil, vtkIdType root, vtkVariantArray *childEdge, std::vector<std::string>& names) const
267 {
268   vtkIdType InfoGeoType(sil->AddChild(root,childEdge));
269   names.push_back(_type_str);
270   vtkIdType InfoVTKID(sil->AddChild(InfoGeoType,childEdge));
271   std::ostringstream oss; oss << _vtkt;
272   names.push_back(oss.str());
273 }
274
275 ////////////////////
276
277 vtkExtractCellType::vtkExtractCellType():SIL(NULL),Internal(new vtkExtractCellTypeInternal),InsideOut(0)
278 {
279 }
280
281 vtkExtractCellType::~vtkExtractCellType()
282 {
283   if(this->SIL)
284     this->SIL->Delete();
285   delete this->Internal;
286 }
287
288 void vtkExtractCellType::SetInsideOut(int val)
289 {
290   if(this->InsideOut!=val)
291     {
292       this->InsideOut=val;
293       this->Modified();
294     }
295 }
296
297 int vtkExtractCellType::RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
298 {
299   try
300     {
301       //std::cerr << "########################################## vtkExtractCellType::RequestInformation ##########################################" << std::endl;
302       vtkInformation *outInfo(outputVector->GetInformationObject(0));
303       vtkInformation *inputInfo(inputVector[0]->GetInformationObject(0));
304       vtkDataSet *input(0);
305       {
306         vtkDataObject *inp(inputInfo->Get(vtkDataObject::DATA_OBJECT()));
307         if(vtkDataSet::SafeDownCast(inp))
308           input=vtkDataSet::SafeDownCast(inp);
309         else
310           {
311             vtkMultiBlockDataSet *inputTmp(vtkMultiBlockDataSet::SafeDownCast(inp));
312             if(inputTmp)
313               {
314                 if(inputTmp->GetNumberOfBlocks()!=1)
315                   {
316                     vtkDebugMacro("vtkExtractCellType::RequestInformation : input vtkMultiBlockDataSet must contain exactly 1 block !");
317                     return 0;
318                   }
319                 vtkDataSet *blk0(vtkDataSet::SafeDownCast(inputTmp->GetBlock(0)));
320                 if(!blk0)
321                   {
322                     vtkDebugMacro("vtkExtractCellType::RequestInformation : the single block in input vtkMultiBlockDataSet must be a vtkDataSet instance !");
323                     return 0;
324                   }
325                 input=blk0;
326               }
327             else
328               {
329                 vtkDebugMacro("vtkExtractCellType::RequestInformation : supported input are vtkDataSet or vtkMultiBlockDataSet !");
330                 return 0;
331               }
332           }
333       }
334       if(this->Internal->setRefTime(input))
335         {
336           vtkIdType nbOfCells(input->GetNumberOfCells());
337           std::map<int,INTERP_KERNEL::NormalizedCellType> m;
338           for(vtkIdType cellId=0;cellId<nbOfCells;cellId++)
339             {
340               int vtkCt(input->GetCellType(cellId));
341               const std::map<int,INTERP_KERNEL::NormalizedCellType>::const_iterator it(m.find(vtkCt));
342               if(it==m.end())
343                 {
344                   const unsigned char *pos(std::find(MEDCoupling::MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE,MEDCoupling::MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE+MEDCoupling::MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE_LGTH,vtkCt));
345                   if(pos==MEDCoupling::MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE+MEDCoupling::MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE_LGTH)
346                     {
347                       vtkDebugMacro("vtkExtractCellType::RequestInformation : cell #" << cellId << " has unrecognized type !");
348                       return 0;
349                     }
350                   m[vtkCt]=(INTERP_KERNEL::NormalizedCellType)std::distance(MEDCoupling::MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE,pos);
351             }
352             }
353           this->Internal->loadFrom(m);
354           if(this->SIL)
355             this->SIL->Delete();
356           this->SIL=vtkMutableDirectedGraph::New();
357           this->Internal->feedSIL(this->SIL);
358           //
359           outInfo->Set(vtkDataObject::SIL(),this->SIL);
360         }
361     }
362   catch(INTERP_KERNEL::Exception& e)
363     {
364       std::cerr << "Exception has been thrown in vtkExtractCellType::RequestInformation : " << e.what() << std::endl;
365       return 0;
366     }
367   return 1;
368 }
369
370 vtkDataSet *FilterFamilies(vtkDataSet *input, const std::vector<int>& idsToKeep, bool insideOut)
371 {
372   const int VTK_DATA_ARRAY_DELETE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_DELETE;
373   const char ZE_SELECTION_ARR_NAME[]="@@ZeSelection@@";
374   vtkDataSet *output(input->NewInstance());
375   output->ShallowCopy(input);
376   vtkSmartPointer<vtkThreshold> thres(vtkSmartPointer<vtkThreshold>::New());
377   thres->SetInputData(output);
378   vtkDataSetAttributes *dscIn(input->GetCellData()),*dscIn2(input->GetPointData());
379   vtkDataSetAttributes *dscOut(output->GetCellData()),*dscOut2(output->GetPointData());
380   //
381   double vMin(insideOut==0?1.:0.),vMax(insideOut==0?2.:1.);
382   thres->ThresholdBetween(vMin,vMax);
383   // OK for the output 
384   vtkIdType nbOfCells(input->GetNumberOfCells());
385   vtkCharArray *zeSelection(vtkCharArray::New());
386   zeSelection->SetName(ZE_SELECTION_ARR_NAME);
387   zeSelection->SetNumberOfComponents(1);
388   char *pt(new char[nbOfCells]);
389   zeSelection->SetArray(pt,nbOfCells,0,VTK_DATA_ARRAY_DELETE);
390   std::fill(pt,pt+nbOfCells,0);
391   std::vector<bool> pt2(nbOfCells,false);
392   for(std::vector<int>::const_iterator it=idsToKeep.begin();it!=idsToKeep.end();it++)
393     {
394       for(vtkIdType ii=0;ii<nbOfCells;ii++)
395         {
396           if(input->GetCellType(ii)==*it)
397             pt2[ii]=true;
398         }
399     }
400   for(int ii=0;ii<nbOfCells;ii++)
401     if(pt2[ii])
402       pt[ii]=2;
403   int idx(output->GetCellData()->AddArray(zeSelection));
404   output->GetCellData()->SetActiveAttribute(idx,vtkDataSetAttributes::SCALARS);
405   output->GetCellData()->CopyScalarsOff();
406   zeSelection->Delete();
407   //
408   thres->SetInputArrayToProcess(idx,0,0,"vtkDataObject::FIELD_ASSOCIATION_CELLS",ZE_SELECTION_ARR_NAME);
409   thres->Update();
410   vtkUnstructuredGrid *zeComputedOutput(thres->GetOutput());
411   zeComputedOutput->GetCellData()->RemoveArray(idx);
412   output->Delete();
413   zeComputedOutput->Register(0);
414   return zeComputedOutput;
415 }
416
417 int vtkExtractCellType::RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
418 {
419   try
420     {
421       //std::cerr << "########################################## vtkExtractCellType::RequestData        ##########################################" << std::endl;
422       vtkInformation* inputInfo=inputVector[0]->GetInformationObject(0);
423       vtkDataSet *input(vtkDataSet::SafeDownCast(inputInfo->Get(vtkDataObject::DATA_OBJECT())));
424       vtkInformation *info(input->GetInformation());
425       vtkInformation *outInfo(outputVector->GetInformationObject(0));
426       vtkDataSet *output(vtkDataSet::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT())));
427       std::vector<int> idsToKeep(this->Internal->getIdsToKeep());
428       vtkDataSet *tryOnCell(FilterFamilies(input,idsToKeep,this->InsideOut));
429       // first shrink the input
430       output->ShallowCopy(tryOnCell);
431       tryOnCell->Delete();
432     }
433   catch(INTERP_KERNEL::Exception& e)
434     {
435       std::cerr << "Exception has been thrown in vtkExtractCellType::RequestData : " << e.what() << std::endl;
436       return 0;
437     }
438   return 1;
439 }
440
441 int vtkExtractCellType::GetSILUpdateStamp()
442 {
443   return this->SILTime;
444 }
445
446 void vtkExtractCellType::PrintSelf(ostream& os, vtkIndent indent)
447 {
448   this->Superclass::PrintSelf(os, indent);
449 }
450
451 int vtkExtractCellType::GetNumberOfGeoTypesArrays()
452 {
453   int ret(this->Internal->getNumberOfEntries());
454   //std::cerr << "vtkExtractCellType::GetNumberOfGeoTypesArrays() -> " << ret << std::endl;
455   return ret;
456 }
457
458 const char *vtkExtractCellType::GetGeoTypesArrayName(int index)
459 {
460   const char *ret(this->Internal->getKeyOfEntry(index));
461   //std::cerr << "vtkExtractCellType::GetGeoTypesArrayName(" << index << ") -> " << ret << std::endl;
462   return ret;
463 }
464
465 int vtkExtractCellType::GetGeoTypesArrayStatus(const char *name)
466 {
467   int ret((int)this->Internal->getStatusOfEntryStr(name));
468   //std::cerr << "vtkExtractCellType::GetGeoTypesArrayStatus(" << name << ") -> " << ret << std::endl;
469   return ret;
470 }
471
472 void vtkExtractCellType::SetGeoTypesStatus(const char *name, int status)
473 {
474   //std::cerr << "vtkExtractCellType::SetGeoTypesStatus(" << name << "," << status << ")" << std::endl;
475   if (GetNumberOfGeoTypesArrays()<1)
476     return;
477   this->Internal->setStatusOfEntryStr(name,(bool)status);
478   if(std::string(name)==GetGeoTypesArrayName(GetNumberOfGeoTypesArrays()-1))
479     {
480       this->Modified();
481       //this->Internal->printMySelf(std::cerr);
482     }
483 }