1 // Copyright (C) 2020-2021 CEA/DEN, EDF R&D
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.
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.
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
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
19 // Author : Anthony Geay (EDF R&D)
21 #include "vtkUgSelectCellIds.h"
23 #include "vtkCellArray.h"
24 #include "vtkUnstructuredGrid.h"
25 #include "vtkInformationVector.h"
26 #include "vtkUnsignedCharArray.h"
27 #include "vtkInformation.h"
28 #include "vtkCellData.h"
29 #include "vtkPointData.h"
31 vtkStandardNewMacro(vtkUgSelectCellIds)
33 void vtkUgSelectCellIds::SetIds(vtkIdTypeArray *ids)
38 int vtkUgSelectCellIds::RequestData(vtkInformation *vtkNotUsed(request), vtkInformationVector **inputVector, vtkInformationVector *outputVector)
40 vtkInformation* inputInfo=inputVector[0]->GetInformationObject(0);
41 vtkUnstructuredGrid *ds(vtkUnstructuredGrid::SafeDownCast(inputInfo->Get(vtkDataObject::DATA_OBJECT())));
44 vtkErrorMacro("vtkUgSelectCellIds::RequestData : input is expected to be an unstructuredgrid !");
47 if( ! this->_ids.Get() )
49 vtkErrorMacro("vtkUgSelectCellIds::RequestData : not initialized array !");
52 vtkInformation *outInfo(outputVector->GetInformationObject(0));
53 vtkUnstructuredGrid *output(vtkUnstructuredGrid::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT())));
54 output->SetPoints(ds->GetPoints());
55 vtkNew<vtkCellArray> outCellArray;
56 vtkNew<vtkIdTypeArray> outConn;
57 vtkCellData *inputCellData(ds->GetCellData());
58 vtkPointData *inputPointData(ds->GetPointData());
59 // Finish the parlote : feed data arrays
60 vtkCellData *outCellData(output->GetCellData());
61 vtkPointData *outPointData(output->GetPointData());
62 vtkIdType inpNbCells( ds->GetNumberOfCells() );
63 vtkIdType outputNbCells( _ids->GetNumberOfTuples() );
64 vtkNew<vtkIdTypeArray> outNodalConn;
65 vtkNew<vtkUnsignedCharArray> outCellTypes;
66 outCellTypes->SetNumberOfComponents(1); outCellTypes->SetNumberOfTuples(outputNbCells);
67 vtkNew<vtkIdTypeArray> outCellLocations;
68 outCellLocations->SetNumberOfComponents(1); outCellLocations->SetNumberOfTuples(outputNbCells);
69 vtkIdType outConnLgth(0);
70 for( const vtkIdType *cellId = _ids->GetPointer(0) ; cellId != _ids->GetPointer(outputNbCells) ; ++cellId )
72 if( *cellId >=0 && *cellId < inpNbCells )
76 ds->GetCellPoints(*cellId, npts, pts);
77 outConnLgth += 1+npts;
81 vtkErrorMacro(<< "vtkUgSelectCellIds::RequestData : presence of " << *cellId << " in array must be in [0," << inpNbCells << "[ !");
85 outNodalConn->SetNumberOfComponents(1); outNodalConn->SetNumberOfTuples(outConnLgth);
86 vtkIdType *outConnPt(outNodalConn->GetPointer(0)),*outCellLocPt(outCellLocations->GetPointer(0));
87 vtkIdType outCurCellLoc = 0;
88 unsigned char *outCellTypePt(outCellTypes->GetPointer(0));
89 for( const vtkIdType *cellId = _ids->GetPointer(0) ; cellId != _ids->GetPointer(outputNbCells) ; ++cellId )
91 outCellLocPt[0] = outCurCellLoc;
94 ds->GetCellPoints(*cellId, npts, pts);
96 outConnPt = std::copy(pts,pts+npts,outConnPt);
97 *outCellTypePt++ = ds->GetCellType(*cellId);
98 outCurCellLoc += npts+1;//not sure
100 for( int cellFieldId = 0 ; cellFieldId < inputCellData->GetNumberOfArrays() ; ++cellFieldId )
102 vtkDataArray *array( inputCellData->GetArray(cellFieldId) );
103 vtkSmartPointer<vtkDataArray> outArray;
104 outArray.TakeReference( array->NewInstance() );
105 outArray->SetNumberOfComponents(array->GetNumberOfComponents()); outArray->SetNumberOfTuples(outputNbCells);
106 outArray->SetName(array->GetName());
108 for( const vtkIdType *cellId = _ids->GetPointer(0) ; cellId != _ids->GetPointer(outputNbCells) ; ++cellId )
110 outArray->SetTuple(pos++,*cellId,array);
112 outCellData->AddArray(outArray);
114 for( int pointFieldId = 0 ; pointFieldId < inputPointData->GetNumberOfArrays() ; ++pointFieldId )
116 vtkDataArray *array( inputPointData->GetArray(pointFieldId) );
117 vtkSmartPointer<vtkDataArray> outArray;
118 outArray.TakeReference( array->NewInstance() );
119 outArray->ShallowCopy(array);
120 outPointData->AddArray(outArray);
123 outCellArray->SetCells(outputNbCells,outNodalConn);
124 output->SetCells(outCellTypes,outCellLocations,outCellArray);