Salome HOME
Copyright update 2022
[modules/paravis.git] / src / Plugins / MEDReader / plugin / MEDLoaderForPV / vtkGenerateVectors.cxx
1 // Copyright (C) 2010-2022  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 "vtkGenerateVectors.h"
22 #include "vtkAOSDataArrayTemplate.h"
23 #include "vtkDoubleArray.h"
24 #include "vtkInformation.h"
25 #include "vtkUnstructuredGrid.h"
26 #include "vtkQuadratureSchemeDefinition.h"
27 #include "vtkInformationQuadratureSchemeDefinitionVectorKey.h"
28 #include "MEDUtilities.hxx"
29 #include "vtkFieldData.h"
30
31 #include <sstream>
32
33 const char vtkGenerateVectors::VECTOR_SUFFIX[]="_Vector";
34
35 std::string vtkGenerateVectors::SuffixFieldName(const std::string& name)
36 {
37   std::ostringstream oss; oss << name << VECTOR_SUFFIX;
38   return oss.str();
39 }
40
41 /*!
42  * This method forces MeshMTime modification. To do so, points are declared as modified.
43  */
44 void vtkGenerateVectors::ChangeMeshTimeToUpdateCache(vtkDataSet *dataset)
45 {
46   vtkUnstructuredGrid *ug(vtkUnstructuredGrid::SafeDownCast(dataset));
47   if(!ug)
48     return ;
49   ug->GetPoints()->Modified();
50 }
51
52 void vtkGenerateVectors::Operate(vtkFieldData *fd)
53 {
54   if(!fd)
55     return ;
56   const int nbOfArrs(fd->GetNumberOfArrays());
57   std::vector<vtkDoubleArray *> daToAppend;
58   for(int i=0;i<nbOfArrs;i++)
59     {
60       vtkDataArray *arr(fd->GetArray(i));
61       if(!arr)
62         continue;
63       vtkDoubleArray *arrc(vtkDoubleArray::SafeDownCast(arr));
64       if(!arrc)
65         continue;
66       int nbOfCompo(arrc->GetNumberOfComponents());
67       if(nbOfCompo<=1 || nbOfCompo==3)
68         continue;
69       if(nbOfCompo==2)
70         daToAppend.push_back(Operate2Compo(arrc));
71       else
72         daToAppend.push_back(OperateMoreThan3Compo(arrc));
73     }
74   for(std::vector<vtkDoubleArray *>::const_iterator it=daToAppend.begin();it!=daToAppend.end();it++)
75     {
76       vtkDoubleArray *elt(*it);
77       if(!elt)
78         continue;
79       fd->AddArray(elt);
80       elt->Delete();
81     }
82 }
83
84 vtkDoubleArray *vtkGenerateVectors::Operate2Compo(vtkDoubleArray *oldArr)
85 {
86   const int VTK_DATA_ARRAY_FREE=vtkAOSDataArrayTemplate<double>::VTK_DATA_ARRAY_FREE;
87   vtkDoubleArray *ret(vtkDoubleArray::New());
88   vtkIdType nbOfTuples(oldArr->GetNumberOfTuples());
89   const double *inPt(oldArr->GetPointer(0));
90   double *pt((double *)malloc(nbOfTuples*3*sizeof(double)));
91   for(vtkIdType i=0;i<nbOfTuples;i++)
92     {
93       pt[3*i+0]=inPt[2*i+0];
94       pt[3*i+1]=inPt[2*i+1];
95       pt[3*i+2]=0.;
96     }
97   ret->SetNumberOfComponents(3);
98   std::string newName(SuffixFieldName(oldArr->GetName()));
99   ret->SetName(newName.c_str());
100   ret->SetComponentName(0,oldArr->GetComponentName(0));
101   ret->SetComponentName(1,oldArr->GetComponentName(1));
102   ret->SetArray(pt,3*nbOfTuples,0,VTK_DATA_ARRAY_FREE);
103   UpdateInformationOfArray(oldArr,ret);
104   return ret;
105 }
106
107 vtkDoubleArray *vtkGenerateVectors::OperateMoreThan3Compo(vtkDoubleArray *oldArr)
108 {
109   const int VTK_DATA_ARRAY_FREE=vtkAOSDataArrayTemplate<double>::VTK_DATA_ARRAY_FREE;
110   vtkDoubleArray *ret(vtkDoubleArray::New());
111   int nbOfCompo(oldArr->GetNumberOfComponents());
112   vtkIdType nbOfTuples(oldArr->GetNumberOfTuples());
113   const double *inPt(oldArr->GetPointer(0));
114   double *pt((double *)malloc(nbOfTuples*3*sizeof(double)));
115   for(vtkIdType i=0;i<nbOfTuples;i++)
116     {
117       pt[3*i+0]=inPt[nbOfCompo*i+0];
118       pt[3*i+1]=inPt[nbOfCompo*i+1];
119       pt[3*i+2]=inPt[nbOfCompo*i+2];
120     }
121   ret->SetNumberOfComponents(3);
122   std::string newName(SuffixFieldName(oldArr->GetName()));
123   ret->SetName(newName.c_str());
124   ret->SetComponentName(0,oldArr->GetComponentName(0));
125   ret->SetComponentName(1,oldArr->GetComponentName(1));
126   ret->SetComponentName(2,oldArr->GetComponentName(2));
127   ret->SetArray(pt,3*nbOfTuples,0,VTK_DATA_ARRAY_FREE);
128   UpdateInformationOfArray(oldArr,ret);
129   return ret;
130 }
131
132 void vtkGenerateVectors::UpdateInformationOfArray(vtkDoubleArray *oldArr, vtkDoubleArray *arr)
133 {
134   if(oldArr->GetInformation()->Has(vtkQuadratureSchemeDefinition::QUADRATURE_OFFSET_ARRAY_NAME()))
135     {
136       arr->GetInformation()->Set(vtkQuadratureSchemeDefinition::QUADRATURE_OFFSET_ARRAY_NAME(),oldArr->GetInformation()->Get((vtkQuadratureSchemeDefinition::QUADRATURE_OFFSET_ARRAY_NAME())));
137     }
138   if(oldArr->GetInformation()->Get(MEDUtilities::ELGA()))
139     arr->GetInformation()->Set(MEDUtilities::ELGA(),1);
140   vtkInformationQuadratureSchemeDefinitionVectorKey *key(vtkQuadratureSchemeDefinition::DICTIONARY());
141   if(key->Has(oldArr->GetInformation()))
142     {
143       int dictSize(key->Size(oldArr->GetInformation()));
144       vtkQuadratureSchemeDefinition **dict = new vtkQuadratureSchemeDefinition *[dictSize];
145       key->GetRange(oldArr->GetInformation(),dict,0,0,dictSize);
146       key->SetRange(arr->GetInformation(),dict,0,0,dictSize);
147       delete [] dict;
148     }
149   if(oldArr->GetInformation()->Get(MEDUtilities::ELNO()))
150     arr->GetInformation()->Set(MEDUtilities::ELNO(),1);
151 }