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