]> SALOME platform Git repositories - modules/paravis.git/blob - src/Plugins/MEDReader/IO/vtkELNOFilter.cxx
Salome HOME
8af029279e19edbd389e9d1bd8f5518d6126a85a
[modules/paravis.git] / src / Plugins / MEDReader / IO / vtkELNOFilter.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
20 #include "vtkELNOFilter.h"
21 #include "vtkInformation.h"
22 #include "vtkInformationVector.h"
23 #include "vtkObjectFactory.h"
24 #include "vtkPolyDataAlgorithm.h"
25 #include "vtkPolyData.h"
26 #include "vtkIdTypeArray.h"
27 #include "vtkFieldData.h"
28 #include "vtkCellData.h"
29 #include "vtkPointData.h"
30 #include "vtkCell.h"
31 #include "vtkInformationQuadratureSchemeDefinitionVectorKey.h"
32 #include "vtkQuadratureSchemeDefinition.h"
33 #include "vtkUnstructuredGrid.h"
34
35 #include "MEDUtilities.hxx"
36
37 //vtkCxxRevisionMacro(vtkELNOFilter, "$Revision: 1.2.2.2 $");
38 vtkStandardNewMacro(vtkELNOFilter);
39
40 vtkELNOFilter::vtkELNOFilter()
41 {
42   this->ShrinkFactor = 0.5;
43 }
44
45 vtkELNOFilter::~vtkELNOFilter()
46 {
47 }
48
49 int vtkELNOFilter::RequestData(vtkInformation *request, vtkInformationVector **input, vtkInformationVector *output)
50 {
51   vtkUnstructuredGrid *usgIn(vtkUnstructuredGrid::SafeDownCast( input[0]->GetInformationObject(0)->Get(vtkDataObject::DATA_OBJECT())));
52   vtkPolyData *pdOut(vtkPolyData::SafeDownCast(output->GetInformationObject(0)->Get(vtkDataObject::DATA_OBJECT())));
53
54   vtkDataArray *array(this->GetInputArrayToProcess(0, input));
55   vtkIdTypeArray* offsets(vtkIdTypeArray::SafeDownCast(this->GetInputArrayToProcess(0, input)));
56
57   if(usgIn == NULL || offsets == NULL || pdOut == NULL)
58     {
59       vtkDebugMacro("vtkELNOFilter no correctly configured : offsets = " << offsets);
60       return 1;
61     }
62
63   vtkInformation *info(offsets->GetInformation());
64   vtkInformationQuadratureSchemeDefinitionVectorKey *key(vtkQuadratureSchemeDefinition::DICTIONARY());
65   if(!key->Has(info))
66     {
67       vtkDebugMacro("Dictionary is not present in array " << offsets->GetName() << " " << offsets << " Aborting." );
68       return 1;
69     }
70
71   int res(this->Superclass::RequestData(request, input, output));
72   if(res == 0)
73     return 0;
74   
75   int dictSize(key->Size(info));
76   vtkQuadratureSchemeDefinition **dict(new vtkQuadratureSchemeDefinition *[dictSize]);
77   key->GetRange(info, dict, 0, 0, dictSize);
78
79   vtkIdType ncell(usgIn->GetNumberOfCells());
80   vtkPoints *points(pdOut->GetPoints());
81   vtkIdType start(0);
82   for(vtkIdType cellId = 0; cellId < ncell; cellId++)
83     {
84       vtkIdType offset(offsets->GetValue(cellId));
85       int cellType(usgIn->GetCellType(cellId));
86       // a simple check to see if a scheme really exists for this cell type.
87       // should not happen if the cell type has not been modified.
88       if(dict[cellType] == NULL)
89         continue;
90       int np = dict[cellType]->GetNumberOfQuadraturePoints();
91       double center[3] = {0, 0, 0};
92       for(int id = start; id < start + np; id++)
93         {
94           double *position = points->GetPoint(id);
95           center[0] += position[0];
96           center[1] += position[1];
97           center[2] += position[2];
98         }
99       center[0] /= np;
100       center[1] /= np;
101       center[2] /= np;
102       for(int id = start; id < start + np; id++)
103         {
104           double *position = points->GetPoint(id);
105           double newpos[3];
106           newpos[0] = position[0] * this->ShrinkFactor + center[0] * (1 - this->ShrinkFactor);
107           newpos[1] = position[1] * this->ShrinkFactor + center[1] * (1 - this->ShrinkFactor);
108           newpos[2] = position[2] * this->ShrinkFactor + center[2] * (1 - this->ShrinkFactor);
109           points->SetPoint(id, newpos);
110         }
111       start += np;
112     }
113   //// bug EDF 8407 PARAVIS - mantis 22610
114   vtkFieldData *fielddata(usgIn->GetFieldData());
115   for(int index=0;index<fielddata->GetNumberOfArrays();index++)
116     {
117       vtkDataArray *data(fielddata->GetArray(index));
118       vtkInformation *info(data->GetInformation());
119       const char *arrayOffsetName(info->Get(vtkQuadratureSchemeDefinition::QUADRATURE_OFFSET_ARRAY_NAME()));
120       vtkIdTypeArray *offData(0);
121       bool isELNO(false);
122       if(arrayOffsetName)
123         {
124           vtkCellData *cellData(usgIn->GetCellData());
125           vtkDataArray *offDataTmp(cellData->GetArray(arrayOffsetName));
126           isELNO=offDataTmp->GetInformation()->Get(MEDUtilities::ELNO())==1;
127           offData=dynamic_cast<vtkIdTypeArray *>(offDataTmp);
128         }
129       if(isELNO && offData)
130         {
131           vtkIdType nbCellsInput(usgIn->GetNumberOfCells());
132           if(nbCellsInput==0)
133             continue ;//no cells -> no fields
134           // First trying to detected if we are in the special case where data can be used directly. To detect that look at offData. If offData.front()==0 && offData->back()+NbOfNodesInLastCell==data->GetNumberOfTuples() OK.
135           vtkCell *cell(usgIn->GetCell(nbCellsInput-1));
136           bool statement0(offData->GetTuple1(0)==0);
137           bool statement1(offData->GetTuple1(nbCellsInput-1)+cell->GetNumberOfPoints()==data->GetNumberOfTuples());
138           if(statement0 && statement1)
139             pdOut->GetPointData()->AddArray(data);//We are lucky ! No extraction needed.
140           else
141             {//not lucky ! Extract values from data. A previous threshold has been done... Bug EDF8662
142               vtkDataArray *newArray(data->NewInstance());
143               newArray->SetName(data->GetName());
144               pdOut->GetPointData()->AddArray(newArray);
145               newArray->SetNumberOfComponents(data->GetNumberOfComponents());
146               newArray->SetNumberOfTuples(pdOut->GetNumberOfPoints());
147               newArray->CopyComponentNames(data);
148               newArray->Delete();
149               vtkIdType *offsetPtr(offData->GetPointer(0));
150               vtkIdType zeId(0);
151               for(vtkIdType cellId=0;cellId<nbCellsInput;cellId++)
152                 {
153                   vtkCell *cell(usgIn->GetCell(cellId));
154                   vtkIdType nbPoints(cell->GetNumberOfPoints()),offset(offsetPtr[cellId]);
155                   for(vtkIdType j=0;j<nbPoints;j++,zeId++)
156                     newArray->SetTuple(zeId,offsetPtr[cellId]+j,data);
157                 }
158             }
159         }
160     }
161   return 1;
162 }
163
164 void vtkELNOFilter::PrintSelf(ostream& os, vtkIndent indent)
165 {
166   this->Superclass::PrintSelf(os, indent);
167   os << indent << "ShrinkFactor : " << this->ShrinkFactor << endl;
168 }