Salome HOME
Merge branch 'V7_5_BR'
[modules/paravis.git] / src / Plugins / MEDReader / IO / vtkExtractCellType.cxx
1 // Copyright (C) 2010-2014  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   const ExtractCellTypeStatus& elt(getEntry(entry));
214   return elt.getStatus();
215 }
216
217 void vtkExtractCellType::vtkExtractCellTypeInternal::setStatusOfEntryStr(const char *entry, bool status) const
218 {
219   const ExtractCellTypeStatus& elt(getEntry(entry));
220   elt.setStatus(status);
221 }
222
223 void vtkExtractCellType::vtkExtractCellTypeInternal::printMySelf(std::ostream& os) const
224 {
225   for(std::vector<ExtractCellTypeStatus>::const_iterator it0=_types.begin();it0!=_types.end();it0++)
226     (*it0).printMySelf(os);
227 }
228
229 ExtractCellTypeStatus::ExtractCellTypeStatus(int vtkt, INTERP_KERNEL::NormalizedCellType mct):_status(false),_vtkt(vtkt),_mct(mct)
230 {
231   std::string name(INTERP_KERNEL::CellModel::GetCellModel(mct).getRepr());
232   _type_str=name.substr(5);//skip "NORM_"
233 }
234
235 void ExtractCellTypeStatus::printMySelf(std::ostream& os) const
236 {
237   os << "      -" << _type_str << "(";
238   if(_status)
239     os << "X";
240   else
241     os << " ";
242   os << ")" << std::endl;
243 }
244
245 bool ExtractCellTypeStatus::isSameAs(const ExtractCellTypeStatus& other) const
246 {
247   return _vtkt==other._vtkt && _mct==other._mct;
248 }
249
250 void ExtractCellTypeStatus::feedSIL(vtkMutableDirectedGraph *sil, vtkIdType root, vtkVariantArray *childEdge, std::vector<std::string>& names) const
251 {
252   vtkIdType InfoGeoType(sil->AddChild(root,childEdge));
253   names.push_back(_type_str);
254   vtkIdType InfoVTKID(sil->AddChild(InfoGeoType,childEdge));
255   std::ostringstream oss; oss << _vtkt;
256   names.push_back(oss.str());
257 }
258
259 ////////////////////
260
261 vtkExtractCellType::vtkExtractCellType():SIL(NULL),Internal(new vtkExtractCellTypeInternal),InsideOut(0)
262 {
263 }
264
265 vtkExtractCellType::~vtkExtractCellType()
266 {
267   if(this->SIL)
268     this->SIL->Delete();
269   delete this->Internal;
270 }
271
272 void vtkExtractCellType::SetInsideOut(int val)
273 {
274   if(this->InsideOut!=val)
275     {
276       this->InsideOut=val;
277       this->Modified();
278     }
279 }
280
281 int vtkExtractCellType::RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
282 {
283   try
284     {
285       //std::cerr << "########################################## vtkExtractCellType::RequestInformation ##########################################" << std::endl;
286       vtkInformation *outInfo(outputVector->GetInformationObject(0));
287       vtkInformation *inputInfo(inputVector[0]->GetInformationObject(0));
288       vtkDataSet *input(0);
289       {
290         vtkDataObject *inp(inputInfo->Get(vtkDataObject::DATA_OBJECT()));
291         if(vtkDataSet::SafeDownCast(inp))
292           input=vtkDataSet::SafeDownCast(inp);
293         else
294           {
295             vtkMultiBlockDataSet *inputTmp(vtkMultiBlockDataSet::SafeDownCast(inp));
296             if(inputTmp)
297               {
298                 if(inputTmp->GetNumberOfBlocks()!=1)
299                   {
300                     vtkDebugMacro("vtkExtractCellType::RequestInformation : input vtkMultiBlockDataSet must contain exactly 1 block !");
301                     return 0;
302                   }
303                 vtkDataSet *blk0(vtkDataSet::SafeDownCast(inputTmp->GetBlock(0)));
304                 if(!blk0)
305                   {
306                     vtkDebugMacro("vtkExtractCellType::RequestInformation : the single block in input vtkMultiBlockDataSet must be a vtkDataSet instance !");
307                     return 0;
308                   }
309                 input=blk0;
310               }
311             else
312               {
313                 vtkDebugMacro("vtkExtractCellType::RequestInformation : supported input are vtkDataSet or vtkMultiBlockDataSet !");
314                 return 0;
315               }
316           }
317       }
318       if(this->Internal->setRefTime(input))
319         {
320           vtkIdType nbOfCells(input->GetNumberOfCells());
321           std::map<int,INTERP_KERNEL::NormalizedCellType> m;
322           for(vtkIdType cellId=0;cellId<nbOfCells;cellId++)
323             {
324               int vtkCt(input->GetCellType(cellId));
325               const std::map<int,INTERP_KERNEL::NormalizedCellType>::const_iterator it(m.find(vtkCt));
326               if(it==m.end())
327                 {
328                   const unsigned char *pos(std::find(ParaMEDMEM::MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE,ParaMEDMEM::MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE+ParaMEDMEM::MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE_LGTH,vtkCt));
329                   if(pos==ParaMEDMEM::MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE+ParaMEDMEM::MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE_LGTH)
330                     {
331                       vtkDebugMacro("vtkExtractCellType::RequestInformation : cell #" << cellId << " has unrecognized type !");
332                       return 0;
333                     }
334                   m[vtkCt]=(INTERP_KERNEL::NormalizedCellType)std::distance(ParaMEDMEM::MEDMeshMultiLev::PARAMEDMEM_2_VTKTYPE,pos);
335             }
336             }
337           this->Internal->loadFrom(m);
338           if(this->SIL)
339             this->SIL->Delete();
340           this->SIL=vtkMutableDirectedGraph::New();
341           this->Internal->feedSIL(this->SIL);
342           //
343           outInfo->Set(vtkDataObject::SIL(),this->SIL);
344         }
345     }
346   catch(INTERP_KERNEL::Exception& e)
347     {
348       std::cerr << "Exception has been thrown in vtkExtractCellType::RequestInformation : " << e.what() << std::endl;
349       return 0;
350     }
351   return 1;
352 }
353
354 vtkDataSet *FilterFamilies(vtkDataSet *input, const std::vector<int>& idsToKeep, bool insideOut)
355 {
356   static const int VTK_DATA_ARRAY_DELETE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_DELETE;
357   static const char ZE_SELECTION_ARR_NAME[]="@@ZeSelection@@";
358   vtkDataSet *output(input->NewInstance());
359   output->ShallowCopy(input);
360   vtkSmartPointer<vtkThreshold> thres(vtkSmartPointer<vtkThreshold>::New());
361   thres->SetInputData(output);
362   vtkDataSetAttributes *dscIn(input->GetCellData()),*dscIn2(input->GetPointData());
363   vtkDataSetAttributes *dscOut(output->GetCellData()),*dscOut2(output->GetPointData());
364   //
365   double vMin(insideOut==0?1.:0.),vMax(insideOut==0?2.:1.);
366   thres->ThresholdBetween(vMin,vMax);
367   // OK for the output 
368   vtkIdType nbOfCells(input->GetNumberOfCells());
369   vtkCharArray *zeSelection(vtkCharArray::New());
370   zeSelection->SetName(ZE_SELECTION_ARR_NAME);
371   zeSelection->SetNumberOfComponents(1);
372   char *pt(new char[nbOfCells]);
373   zeSelection->SetArray(pt,nbOfCells,0,VTK_DATA_ARRAY_DELETE);
374   std::fill(pt,pt+nbOfCells,0);
375   std::vector<bool> pt2(nbOfCells,false);
376   for(std::vector<int>::const_iterator it=idsToKeep.begin();it!=idsToKeep.end();it++)
377     {
378       for(vtkIdType ii=0;ii<nbOfCells;ii++)
379         {
380           if(input->GetCellType(ii)==*it)
381             pt2[ii]=true;
382         }
383     }
384   for(int ii=0;ii<nbOfCells;ii++)
385     if(pt2[ii])
386       pt[ii]=2;
387   int idx(output->GetCellData()->AddArray(zeSelection));
388   output->GetCellData()->SetActiveAttribute(idx,vtkDataSetAttributes::SCALARS);
389   output->GetCellData()->CopyScalarsOff();
390   zeSelection->Delete();
391   //
392   thres->SetInputArrayToProcess(idx,0,0,"vtkDataObject::FIELD_ASSOCIATION_CELLS",ZE_SELECTION_ARR_NAME);
393   thres->Update();
394   vtkUnstructuredGrid *zeComputedOutput(thres->GetOutput());
395   zeComputedOutput->GetCellData()->RemoveArray(idx);
396   output->Delete();
397   zeComputedOutput->Register(0);
398   return zeComputedOutput;
399 }
400
401 int vtkExtractCellType::RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
402 {
403   try
404     {
405       //std::cerr << "########################################## vtkExtractCellType::RequestData        ##########################################" << std::endl;
406       vtkInformation* inputInfo=inputVector[0]->GetInformationObject(0);
407       vtkDataSet *input(vtkDataSet::SafeDownCast(inputInfo->Get(vtkDataObject::DATA_OBJECT())));
408       vtkInformation *info(input->GetInformation());
409       vtkInformation *outInfo(outputVector->GetInformationObject(0));
410       vtkDataSet *output(vtkDataSet::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT())));
411       std::vector<int> idsToKeep(this->Internal->getIdsToKeep());
412       vtkDataSet *tryOnCell(FilterFamilies(input,idsToKeep,this->InsideOut));
413       // first shrink the input
414       output->ShallowCopy(tryOnCell);
415       tryOnCell->Delete();
416     }
417   catch(INTERP_KERNEL::Exception& e)
418     {
419       std::cerr << "Exception has been thrown in vtkExtractCellType::RequestData : " << e.what() << std::endl;
420       return 0;
421     }
422   return 1;
423 }
424
425 int vtkExtractCellType::GetSILUpdateStamp()
426 {
427   return this->SILTime;
428 }
429
430 void vtkExtractCellType::PrintSelf(ostream& os, vtkIndent indent)
431 {
432   this->Superclass::PrintSelf(os, indent);
433 }
434
435 int vtkExtractCellType::GetNumberOfGeoTypesArrays()
436 {
437   int ret(this->Internal->getNumberOfEntries());
438   //std::cerr << "vtkExtractCellType::GetNumberOfGeoTypesArrays() -> " << ret << std::endl;
439   return ret;
440 }
441
442 const char *vtkExtractCellType::GetGeoTypesArrayName(int index)
443 {
444   const char *ret(this->Internal->getKeyOfEntry(index));
445   //std::cerr << "vtkExtractCellType::GetGeoTypesArrayName(" << index << ") -> " << ret << std::endl;
446   return ret;
447 }
448
449 int vtkExtractCellType::GetGeoTypesArrayStatus(const char *name)
450 {
451   int ret((int)this->Internal->getStatusOfEntryStr(name));
452   //std::cerr << "vtkExtractCellType::GetGeoTypesArrayStatus(" << name << ") -> " << ret << std::endl;
453   return ret;
454 }
455
456 void vtkExtractCellType::SetGeoTypesStatus(const char *name, int status)
457 {
458   //std::cerr << "vtkExtractCellType::SetGeoTypesStatus(" << name << "," << status << ")" << std::endl;
459   this->Internal->setStatusOfEntryStr(name,(bool)status);
460   if(std::string(name)==GetGeoTypesArrayName(GetNumberOfGeoTypesArrays()-1))
461     {
462       this->Modified();
463       //this->Internal->printMySelf(std::cerr);
464     }
465 }