]> SALOME platform Git repositories - tools/paravisaddons_common.git/blob - src/RosetteCIH/plugin/RosetteCIHFilters/vtkRosetteCIH.cxx
Salome HOME
[EDF22802] : Rosette filter improvement
[tools/paravisaddons_common.git] / src / RosetteCIH / plugin / RosetteCIHFilters / vtkRosetteCIH.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
20 #include "vtkRosetteCIH.h"
21
22 #include <vtkCellArray.h>
23 #include <vtkCellData.h>
24 #include <vtkCompositeDataToUnstructuredGridFilter.h>
25 #include <vtkDataSetSurfaceFilter.h>
26 #include <vtkDoubleArray.h>
27 #include <vtkInformation.h>
28 #include <vtkInformationVector.h>
29 #include <vtkLineSource.h>
30 #include <vtkMultiBlockDataGroupFilter.h>
31 #include <vtkMultiBlockDataSet.h>
32 #include <vtkNew.h>
33 #include <vtkObjectFactory.h>
34 #include <vtkRibbonFilter.h>
35 #include <vtkPVGlyphFilter.h>
36 #include <vtkPointData.h>
37 #include <vtkPolyData.h>
38 #include <vtkPolyDataNormals.h>
39 #include <vtkTable.h>
40 #include <vtkTessellatorFilter.h>
41 #include <vtkUnstructuredGrid.h>
42 #include <vtkVariant.h>
43 #include <vtkVariantArray.h>
44
45 //-----------------------------------------------------------------------------
46 void vtkRosetteCIH::ExtractInfo(
47   vtkInformationVector* inputVector, vtkSmartPointer<vtkUnstructuredGrid>& usgIn)
48 {
49   vtkInformation* inputInfo(inputVector->GetInformationObject(0));
50   vtkDataSet* input(0);
51   vtkDataSet* input0(vtkDataSet::SafeDownCast(inputInfo->Get(vtkDataObject::DATA_OBJECT())));
52   vtkMultiBlockDataSet* input1(
53     vtkMultiBlockDataSet::SafeDownCast(inputInfo->Get(vtkDataObject::DATA_OBJECT())));
54   if (input0)
55   {
56     input = input0;
57   }
58   else
59   {
60     if (!input1)
61     {
62       vtkErrorMacro("Input dataSet must be a DataSet or single elt multi block dataset expected !");
63       return;
64     }
65     if (input1->GetNumberOfBlocks() != 1)
66     {
67       vtkErrorMacro("Input dataSet is a multiblock dataset with not exactly one block ! Use "
68                     "MergeBlocks or ExtractBlocks filter before calling this filter !");
69       return;
70     }
71     vtkDataObject* input2(input1->GetBlock(0));
72     if (!input2)
73     {
74       vtkErrorMacro("Input dataSet is a multiblock dataset with exactly one block but this single "
75                     "element is NULL !");
76       return;
77     }
78     vtkDataSet* input2c(vtkDataSet::SafeDownCast(input2));
79     if (!input2c)
80     {
81       vtkErrorMacro(
82         "Input dataSet is a multiblock dataset with exactly one block but this single element is "
83         "not a dataset ! Use MergeBlocks or ExtractBlocks filter before calling this filter !");
84       return;
85     }
86     input = input2c;
87   }
88
89   if (!input)
90   {
91     vtkErrorMacro("Input data set is NULL !");
92     return;
93   }
94
95   usgIn = vtkUnstructuredGrid::SafeDownCast(input);
96   if (!usgIn)
97   {
98     if (!input1)
99     {
100       vtkNew<vtkMultiBlockDataGroupFilter> mb;
101       vtkNew<vtkCompositeDataToUnstructuredGridFilter> cd;
102       mb->AddInputData(input);
103       cd->SetInputConnection(mb->GetOutputPort());
104       cd->SetMergePoints(0);
105       cd->Update();
106       usgIn = cd->GetOutput();
107     }
108     else
109     {
110       vtkNew<vtkCompositeDataToUnstructuredGridFilter> filter;
111       filter->SetMergePoints(0);
112       filter->SetInputData(input1);
113       filter->Update();
114       usgIn = filter->GetOutput();
115     }
116   }
117 }
118
119 //-----------------------------------------------------------------------------
120 vtkStandardNewMacro(vtkRosetteCIH);
121
122 //-----------------------------------------------------------------------------
123 vtkSmartPointer<vtkDataSet> vtkRosetteCIH::GenerateGlyphLinesFor(vtkUnstructuredGrid* usgIn,
124   const char* keyPoint, const char COMPRESS_TRACTION[])
125 {
126   vtkFieldData* dsa(usgIn->GetCellData());
127   std::string arrayForGlyph(this->RetrieveFieldForGlyph(usgIn, keyPoint));
128   int compoId(-1);
129   vtkDoubleArray* arrayForPosNeg2(RetrieveFieldForPost(usgIn, keyPoint, compoId));
130   // vtkAbstractArray *arrayForPosNeg(dsa->GetAbstractArray(FIELD_NAME_2));
131   // vtkDoubleArray *arrayForPosNeg2(vtkDoubleArray::SafeDownCast(arrayForPosNeg));
132   int nbCompo(arrayForPosNeg2->GetNumberOfComponents());
133   vtkIdType nbTuples(arrayForPosNeg2->GetNumberOfTuples());
134   vtkNew<vtkDoubleArray> compressionOrTraction;
135   compressionOrTraction->SetNumberOfComponents(1);
136   compressionOrTraction->SetNumberOfTuples(nbTuples);
137   compressionOrTraction->SetName(COMPRESS_TRACTION);
138   const double* pt(arrayForPosNeg2->GetPointer(0));
139   double* ptOut(compressionOrTraction->GetPointer(0));
140   for (vtkIdType i = 0; i < nbTuples; i++)
141   {
142     if (pt[i * nbCompo + compoId] > 0.)
143       ptOut[i] = 1.;
144     else
145       ptOut[i] = -1.;
146   }
147   int arrId(dsa->AddArray(compressionOrTraction));
148   //
149   vtkNew<vtkPVGlyphFilter> glyph;
150   glyph->SetInputData(usgIn);
151   glyph->SetGlyphMode(0);       // vtkPVGlyphFilter::ALL_POINTS
152   glyph->SetVectorScaleMode(0); // vtkPVGlyphFilter::SCALE_BY_MAGNITUDE
153
154   //
155   vtkNew<vtkLineSource> arrow;
156   glyph->SetSourceConnection(arrow->GetOutputPort());
157   // idx,port,connection,fieldAssociation,name
158   glyph->SetInputArrayToProcess(
159     0, 0, 0, vtkDataObject::FIELD_ASSOCIATION_CELLS, arrayForGlyph.c_str()); // idx==0 -> scaleArray
160   glyph->SetInputArrayToProcess(1, 0, 0, vtkDataObject::FIELD_ASSOCIATION_CELLS,
161     arrayForGlyph.c_str()); // idx==1 -> orientationArray
162   glyph->SetScaleFactor(this->ScaleFactor);
163   glyph->Update();
164
165   return vtkSmartPointer<vtkDataSet>(glyph->GetOutput());
166 }
167
168 //-----------------------------------------------------------------------------
169 void vtkRosetteCIH::PostTraitementT1etT2(
170   vtkUnstructuredGrid* usgIn, vtkUnstructuredGrid* output)
171 {
172   constexpr char COMPRESS_TRACTION[] = "CompressionOrTraction";
173   // "RESUNL__SIRO_ELEM_T1_Vector" , "RESUNL__SIRO_ELEM_T1"
174   vtkSmartPointer<vtkDataSet> gl1 =
175     this->GenerateGlyphLinesFor(usgIn, "T1", COMPRESS_TRACTION);
176   vtkSmartPointer<vtkDataSet> gl2 =
177     this->GenerateGlyphLinesFor(usgIn, "T2", COMPRESS_TRACTION);
178   //
179   vtkNew<vtkDataSetSurfaceFilter> surface;
180   surface->SetNonlinearSubdivisionLevel(0);
181   surface->SetInputData(usgIn);
182   surface->Update();
183   vtkNew<vtkPolyData> surfaceCpy;
184   surfaceCpy->ShallowCopy(surface->GetOutput());
185   vtkNew<vtkDoubleArray> compressionOrTraction;
186   auto nbOfTuples(surface->GetOutput()->GetNumberOfPoints());
187   compressionOrTraction->SetNumberOfComponents(1);
188   compressionOrTraction->SetNumberOfTuples(nbOfTuples);
189   compressionOrTraction->SetName(COMPRESS_TRACTION);
190   compressionOrTraction->Fill(NAN);
191   surfaceCpy->GetPointData()->AddArray(compressionOrTraction);
192   //
193   vtkNew<vtkMultiBlockDataGroupFilter> mb;
194   vtkNew<vtkCompositeDataToUnstructuredGridFilter> cd;
195   mb->AddInputData(surfaceCpy);
196   mb->AddInputData(gl1);
197   mb->AddInputData(gl2);
198   cd->SetInputConnection(mb->GetOutputPort());
199   cd->SetMergePoints(0);
200   cd->Update();
201   //
202   output->ShallowCopy(cd->GetOutput());
203   //
204   vtkFieldData* dsa(output->GetPointData());
205   int nbOfArrays(dsa->GetNumberOfArrays());
206   for (int i = nbOfArrays - 1; i >= 0; i--)
207   {
208     const char* arrName(dsa->GetArrayName(i));
209     if (std::string(arrName) != COMPRESS_TRACTION)
210     {
211       dsa->RemoveArray(i);
212     }
213   }
214   output->GetPointData()->SetActiveAttribute(0, vtkDataSetAttributes::SCALARS);
215 }
216
217 //-----------------------------------------------------------------------------
218 int vtkRosetteCIH::ComponentIdOfArray(vtkAbstractArray* array, const std::string& compoName)
219 {
220   int nbCompo(array->GetNumberOfComponents());
221   int ret(-1);
222   for (int i = 0; i < nbCompo; i++)
223   {
224     if (compoName == array->GetComponentName(i))
225     {
226       if (ret != -1)
227       {
228         vtkErrorMacro("ComponentIdOfArray : already found !");
229         return ret;
230       }
231       ret = i;
232     }
233   }
234   if (ret == -1)
235   {
236     vtkErrorMacro(
237       "ComponentIdOfArray : component " << compoName << " in array " << array->GetName() << " !");
238   }
239   return ret;
240 }
241
242 //-----------------------------------------------------------------------------
243 std::string vtkRosetteCIH::GenerateAValidFieldForGlyph(
244   vtkUnstructuredGrid* usgInCpy, const std::string& arrayName, const char* keyPoint)
245 {
246   vtkFieldData* dsa(usgInCpy->GetCellData());
247   vtkAbstractArray* array(dsa->GetAbstractArray(arrayName.c_str()));
248   vtkDoubleArray* array2(vtkDoubleArray::SafeDownCast(array));
249   //
250   std::string ret(arrayName);
251   ret += std::string("_vveeccttoorr");
252   int compoIds[3] = { -1, -1, -1 };
253   for (int i = 0; i < 3; i++)
254   {
255     std::string compoName("SIG_");
256     compoName += keyPoint;
257     compoName += 'X' + i;
258     compoIds[i] = ComponentIdOfArray(array, compoName);
259   }
260   //
261   vtkIdType nbTuples(array2->GetNumberOfTuples());
262   int nbCompo(array2->GetNumberOfComponents());
263   vtkNew<vtkDoubleArray> vect;
264   vect->SetNumberOfComponents(3);
265   vect->SetNumberOfTuples(nbTuples);
266   vect->SetName(ret.c_str());
267
268   const double* pt(array2->GetPointer(0));
269   double* ptOut(vect->GetPointer(0));
270   for (vtkIdType i = 0; i < nbTuples; i++)
271   {
272     for (int j = 0; j < 3; j++)
273     {
274       ptOut[3 * i + j] = pt[nbCompo * i + compoIds[j]];
275     }
276   }
277   //
278   dsa->AddArray(vect);
279   return ret;
280 }
281
282 //-----------------------------------------------------------------------------
283 bool vtkRosetteCIH::EndWith(const std::string& arrayName, const std::string& end)
284 {
285   std::size_t lenOfLastChance(end.length());
286   if (arrayName.length() < lenOfLastChance)
287   {
288     return false;
289   }
290
291   std::string endOfArrayName(arrayName.substr(arrayName.length() - lenOfLastChance));
292   return endOfArrayName == end;
293 }
294
295 //-----------------------------------------------------------------------------
296 bool vtkRosetteCIH::IsFirstChance(const std::string& arrayName, const char* keyPoint)
297 {
298   std::string PATTERN("SIRO_ELEM");
299   PATTERN += std::string("_") + keyPoint;
300   return this->EndWith(arrayName, PATTERN);
301 }
302
303 //-----------------------------------------------------------------------------
304 bool vtkRosetteCIH::IsLastChanceArray(const std::string& arrayName)
305 {
306   return this->EndWith(arrayName, "SIRO_ELEM");
307 }
308
309 //-----------------------------------------------------------------------------
310 std::string vtkRosetteCIH::GetFieldName(vtkUnstructuredGrid* usgInCpy, const char* keyPoint)
311 {
312   vtkFieldData* dsa(usgInCpy->GetCellData());
313   std::string arrayNameOK;
314   int nbOfArrays(dsa->GetNumberOfArrays());
315   bool found(false);
316   for (int i = 0; i < nbOfArrays; ++i)
317   {
318     vtkAbstractArray* arrayAbstract(dsa->GetAbstractArray(i));
319     std::string arrayName(arrayAbstract->GetName());
320     if (this->IsFirstChance(arrayName, keyPoint) || this->IsLastChanceArray(arrayName))
321     {
322       if (found)
323       {
324         vtkErrorMacro("GetFieldName : already found !");
325       }
326       arrayNameOK = arrayName;
327       found = true;
328     }
329   }
330   if (!found)
331   {
332     vtkErrorMacro("GetFieldName : Impossible to find a valid array !");
333   }
334   return arrayNameOK;
335 }
336
337 //-----------------------------------------------------------------------------
338 std::string vtkRosetteCIH::RetrieveFieldForGlyph(
339   vtkUnstructuredGrid* usgInCpy, const char* keyPoint)
340 {
341   std::string arrayNameOK(this->GetFieldName(usgInCpy, keyPoint));
342   return this->GenerateAValidFieldForGlyph(usgInCpy, arrayNameOK, keyPoint);
343 }
344
345 //-----------------------------------------------------------------------------
346 vtkDoubleArray* vtkRosetteCIH::RetrieveFieldForPost(
347   vtkUnstructuredGrid* usgInCpy, const char* keyPoint, int& compId)
348 {
349   std::string FIELD_NAME_2(this->GetFieldName(usgInCpy, keyPoint));
350   vtkFieldData* dsa(usgInCpy->GetCellData());
351   vtkAbstractArray* arrayForPosNeg(dsa->GetAbstractArray(FIELD_NAME_2.c_str()));
352   vtkDoubleArray* arrayForPosNeg2(vtkDoubleArray::SafeDownCast(arrayForPosNeg));
353   std::string compoToFind("SIG_");
354   compoToFind += keyPoint;
355   compId = this->ComponentIdOfArray(arrayForPosNeg, compoToFind);
356   return arrayForPosNeg2;
357 }
358
359 //-----------------------------------------------------------------------------
360 void vtkRosetteCIH::PostTraitementOnlyOneCompo(vtkUnstructuredGrid* usgIn,
361   vtkUnstructuredGrid* output, const char* keyPoint,
362   const char* COMPRESS_TRACTION)
363 {
364   vtkNew<vtkUnstructuredGrid> usgInCpy;
365   usgInCpy->DeepCopy(usgIn);
366   //
367   vtkFieldData* dsa(usgInCpy->GetCellData());
368   int compId(-1);
369   // vtkAbstractArray *arrayForPosNeg(dsa->GetAbstractArray(FIELD_NAME_2));
370   std::string arrayForGlyph(this->RetrieveFieldForGlyph(usgInCpy, keyPoint));
371   vtkDoubleArray* arrayForPosNeg2(this->RetrieveFieldForPost(usgInCpy, keyPoint, compId));
372   // vtkDoubleArray::SafeDownCast(arrayForPosNeg);
373   int nbCompo(arrayForPosNeg2->GetNumberOfComponents());
374   vtkIdType nbTuples(arrayForPosNeg2->GetNumberOfTuples());
375
376   vtkNew<vtkDoubleArray> compressionOrTraction;
377   compressionOrTraction->SetNumberOfComponents(1);
378   compressionOrTraction->SetNumberOfTuples(nbTuples);
379   compressionOrTraction->SetName(COMPRESS_TRACTION);
380
381   const double* pt(arrayForPosNeg2->GetPointer(0));
382   double* ptOut(compressionOrTraction->GetPointer(0));
383   double valMin(std::numeric_limits<double>::max());
384   double valMax(-std::numeric_limits<double>::max());
385
386   for (vtkIdType i = 0; i < nbTuples; i++)
387   {
388     double val(pt[i * nbCompo + compId]);
389     valMin = std::min(valMin, val);
390     valMax = std::max(valMax, val);
391     ptOut[i] = val;
392   }
393   //
394   for (int i = dsa->GetNumberOfArrays() - 1; i >= 0; i--)
395   {
396     if (arrayForGlyph != dsa->GetAbstractArray(i)->GetName())
397     {
398       dsa->RemoveArray(i);
399     }
400   }
401   int arrId(dsa->AddArray(compressionOrTraction));
402
403   vtkNew<vtkLineSource> arrow;
404
405   vtkNew<vtkDataSetSurfaceFilter> surface;
406   surface->SetNonlinearSubdivisionLevel(0);
407   surface->SetInputData(usgInCpy);
408
409   vtkNew<vtkPolyDataNormals> normals;
410   normals->ComputeCellNormalsOn();
411   normals->ComputePointNormalsOff();
412   normals->SplittingOff();
413   normals->SetInputConnection(surface->GetOutputPort());
414   normals->Update();
415
416   // for some reasons, the glyph filter removes scalars and normals, we have to duplicate them
417   vtkDataArray* normalsArray = normals->GetOutput()->GetCellData()->GetNormals();
418
419   vtkSmartPointer<vtkDataArray> savedNormalsArray;
420   savedNormalsArray.TakeReference(normalsArray->NewInstance());
421   savedNormalsArray->DeepCopy(normalsArray);
422   savedNormalsArray->SetName("CellNormals");
423
424   normals->GetOutput()->GetCellData()->AddArray(savedNormalsArray);
425
426   vtkNew<vtkPVGlyphFilter> glyph;
427   glyph->SetInputConnection(normals->GetOutputPort());
428   glyph->SetGlyphMode(0);       // vtkPVGlyphFilter::ALL_POINTS
429   glyph->SetVectorScaleMode(0); // vtkPVGlyphFilter::SCALE_BY_MAGNITUDE
430   glyph->SetSourceConnection(arrow->GetOutputPort());
431   // idx,port,connection,fieldAssociation,name
432   glyph->SetInputArrayToProcess(
433     0, 0, 0, vtkDataObject::FIELD_ASSOCIATION_CELLS, arrayForGlyph.c_str()); // idx==0 -> scaleArray
434   glyph->SetInputArrayToProcess(1, 0, 0, vtkDataObject::FIELD_ASSOCIATION_CELLS,
435     arrayForGlyph.c_str()); // idx==1 -> orientationArray
436   glyph->SetScaleFactor(this->ScaleFactor);
437
438   vtkNew<vtkRibbonFilter> ribbon;
439   ribbon->SetWidth(this->WidthFactor);
440   ribbon->VaryWidthOff();
441   ribbon->UseDefaultNormalOff();
442   ribbon->SetInputArrayToProcess(1, 0, 0, vtkDataObject::FIELD_ASSOCIATION_POINTS, "CellNormals");
443   ribbon->SetInputConnection(glyph->GetOutputPort());
444   ribbon->Update();
445
446   vtkDataSet* ret = ribbon->GetOutput();
447
448   vtkFieldData* fieldData = ret->GetPointData();
449   for (int i = fieldData->GetNumberOfArrays() - 1; i >= 0; i--)
450   {
451     fieldData->RemoveArray(i);
452   }
453
454   fieldData = ret->GetCellData();
455   for (int i = fieldData->GetNumberOfArrays() - 1; i >= 0; i--)
456   {
457     fieldData->RemoveArray(i);
458   }
459
460   vtkNew<vtkDoubleArray> compressionOrTractionNaN;
461   compressionOrTractionNaN->SetNumberOfComponents(1);
462   compressionOrTractionNaN->SetNumberOfTuples(ret->GetNumberOfCells());
463   compressionOrTractionNaN->SetName(COMPRESS_TRACTION);
464   compressionOrTractionNaN->Fill(NAN);
465   fieldData->AddArray(compressionOrTractionNaN);
466
467   vtkNew<vtkMultiBlockDataGroupFilter> mb;
468   mb->AddInputData(usgInCpy);
469   mb->AddInputData(ret);
470
471   vtkNew<vtkCompositeDataToUnstructuredGridFilter> cd;
472   cd->SetInputConnection(mb->GetOutputPort());
473   cd->SetMergePoints(0);
474   cd->Update();
475
476   output->ShallowCopy(cd->GetOutput());
477
478   int arrayId;
479   output->GetCellData()->GetAbstractArray(COMPRESS_TRACTION, arrayId);
480   output->GetCellData()->SetActiveAttribute(arrayId, vtkDataSetAttributes::SCALARS);
481 }
482
483 //-----------------------------------------------------------------------------
484 void vtkRosetteCIH::PostTraitementT1(vtkUnstructuredGrid* usgIn, vtkUnstructuredGrid* output)
485 {
486   // constexpr char FIELD_NAME[]="RESUNL__SIRO_ELEM_T1_Vector";
487   // constexpr char FIELD_NAME_2[]="RESUNL__SIRO_ELEM_T1";
488   constexpr char COMPRESS_TRACTION[] = "Contrainte specifique 1";
489   this->PostTraitementOnlyOneCompo(usgIn, output, "T1", COMPRESS_TRACTION);
490 }
491
492 //-----------------------------------------------------------------------------
493 void vtkRosetteCIH::PostTraitementT2(vtkUnstructuredGrid* usgIn, vtkUnstructuredGrid* output)
494 {
495   // constexpr char FIELD_NAME[]="RESUNL__SIRO_ELEM_T2_Vector";
496   // constexpr char FIELD_NAME_2[]="RESUNL__SIRO_ELEM_T2";
497   constexpr char COMPRESS_TRACTION[] = "Contrainte specifique 3";
498   this->PostTraitementOnlyOneCompo(usgIn, output, "T2", COMPRESS_TRACTION);
499 }
500
501 //-----------------------------------------------------------------------------
502 int vtkRosetteCIH::RequestData(vtkInformation* vtkNotUsed(request),
503   vtkInformationVector** inputVector, vtkInformationVector* outputVector)
504 {
505   vtkInformation* outInfo(outputVector->GetInformationObject(0));
506   vtkUnstructuredGrid* output(
507     vtkUnstructuredGrid::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT())));
508   //
509   vtkSmartPointer<vtkUnstructuredGrid> usgIn;
510   this->ExtractInfo(inputVector[0], usgIn);
511   switch (this->TypeOfDisplay)
512   {
513     case 0:
514       this->PostTraitementT1etT2(usgIn, output);
515       break;
516     case 1:
517       this->PostTraitementT1(usgIn, output);
518       break;
519     case 2:
520       this->PostTraitementT2(usgIn, output);
521       break;
522     default:
523       vtkErrorMacro("GetFieldName : Impossible to find a valid array !");
524   }
525
526   return 1;
527 }
528
529 //-----------------------------------------------------------------------------
530 void vtkRosetteCIH::PrintSelf(ostream& os, vtkIndent indent)
531 {
532   this->Superclass::PrintSelf(os, indent);
533 }