From: Anthony Geay Date: Fri, 19 Dec 2014 15:22:45 +0000 (+0100) Subject: EDF9622: Correction of ELNOPoints part X-Git-Tag: V7_5_1b1~1 X-Git-Url: http://git.salome-platform.org/gitweb/?a=commitdiff_plain;h=c6ff7f6dfc471ef0fe7ab45e171e05a6335512a9;p=modules%2Fparavis.git EDF9622: Correction of ELNOPoints part --- diff --git a/src/Plugins/MEDReader/IO/vtkELNOFilter.cxx b/src/Plugins/MEDReader/IO/vtkELNOFilter.cxx index 8af02927..46d7a379 100644 --- a/src/Plugins/MEDReader/IO/vtkELNOFilter.cxx +++ b/src/Plugins/MEDReader/IO/vtkELNOFilter.cxx @@ -33,6 +33,7 @@ #include "vtkUnstructuredGrid.h" #include "MEDUtilities.hxx" +#include "InterpKernelAutoPtr.hxx" //vtkCxxRevisionMacro(vtkELNOFilter, "$Revision: 1.2.2.2 $"); vtkStandardNewMacro(vtkELNOFilter); @@ -158,6 +159,7 @@ int vtkELNOFilter::RequestData(vtkInformation *request, vtkInformationVector **i } } } + AttachCellFieldsOn(usgIn,pdOut->GetCellData(),pdOut->GetNumberOfCells()); return 1; } @@ -166,3 +168,38 @@ void vtkELNOFilter::PrintSelf(ostream& os, vtkIndent indent) this->Superclass::PrintSelf(os, indent); os << indent << "ShrinkFactor : " << this->ShrinkFactor << endl; } + +/*! + * This method attach fields on cell of \a inGrid and add it as a point data in \a outData. + */ +void vtkELNOFilter::AttachCellFieldsOn(vtkUnstructuredGrid *inGrid, vtkCellData *outData, int nbCellsOut) +{ + vtkCellData *cd(inGrid->GetCellData()); + int nbOfArrays(cd->GetNumberOfArrays()); + vtkIdType nbCells(inGrid->GetNumberOfCells()); + if(nbOfArrays==0) + return ; + INTERP_KERNEL::AutoPtr tmpPtr(new vtkIdType[nbCells]); + for(vtkIdType cellId=0;cellIdGetCell(cellId)); + tmpPtr[cellId]=cell->GetNumberOfPoints(); + } + for(int index=0;indexGetArray(index)); + vtkDataArray *newArray(data->NewInstance()); + newArray->SetName(data->GetName()); + outData->AddArray(newArray); + newArray->SetNumberOfComponents(data->GetNumberOfComponents()); + newArray->SetNumberOfTuples(nbCellsOut); + newArray->CopyComponentNames(data); + newArray->Delete(); + vtkIdType offset(0); + for(vtkIdType cellId=0;cellIdSetTuple(offset,cellId,data); + } + } +} diff --git a/src/Plugins/MEDReader/IO/vtkELNOFilter.h b/src/Plugins/MEDReader/IO/vtkELNOFilter.h index 87715f93..94319814 100644 --- a/src/Plugins/MEDReader/IO/vtkELNOFilter.h +++ b/src/Plugins/MEDReader/IO/vtkELNOFilter.h @@ -42,9 +42,8 @@ protected: vtkELNOFilter(); ~vtkELNOFilter(); - int RequestData(vtkInformation *, vtkInformationVector **, - vtkInformationVector *); - + int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *); + void AttachCellFieldsOn(vtkUnstructuredGrid *, vtkCellData *, int); double ShrinkFactor; private: diff --git a/src/Plugins/MEDReader/Test/CMakeLists.txt b/src/Plugins/MEDReader/Test/CMakeLists.txt index 42b017da..7f781895 100644 --- a/src/Plugins/MEDReader/Test/CMakeLists.txt +++ b/src/Plugins/MEDReader/Test/CMakeLists.txt @@ -52,3 +52,5 @@ ADD_TEST(testMEDReader15 ${PYTHON_EXECUTABLE} ${CMAKE_CURRENT_SOURCE_DIR}/testME SET_TESTS_PROPERTIES(testMEDReader15 PROPERTIES ENVIRONMENT "${tests_env}") ADD_TEST(testMEDReader16 ${PYTHON_EXECUTABLE} ${CMAKE_CURRENT_SOURCE_DIR}/testMEDReader16.py) SET_TESTS_PROPERTIES(testMEDReader16 PROPERTIES ENVIRONMENT "${tests_env}") +ADD_TEST(testMEDReader16 ${PYTHON_EXECUTABLE} ${CMAKE_CURRENT_SOURCE_DIR}/testMEDReader17.py) +SET_TESTS_PROPERTIES(testMEDReader17 PROPERTIES ENVIRONMENT "${tests_env}") diff --git a/src/Plugins/MEDReader/Test/testMEDReader1.py b/src/Plugins/MEDReader/Test/testMEDReader1.py index 49f033ba..69bde962 100644 --- a/src/Plugins/MEDReader/Test/testMEDReader1.py +++ b/src/Plugins/MEDReader/Test/testMEDReader1.py @@ -26,7 +26,6 @@ This test focused on ELNO. Here a 2 QUAD4 cells and a single ELNO field is defined. """ fname="testMEDReader1.med" -outImgName="testMEDReader1.png" coords=DataArrayDouble([(0,0,0),(2,1,0),(1,0,0),(1,1,0),(2,0,0),(0,1,0)]) m=MEDCouplingUMesh("mesh",2) ; m.setCoords(coords) diff --git a/src/Plugins/MEDReader/Test/testMEDReader16.py b/src/Plugins/MEDReader/Test/testMEDReader16.py index 012c5110..13ded9e4 100644 --- a/src/Plugins/MEDReader/Test/testMEDReader16.py +++ b/src/Plugins/MEDReader/Test/testMEDReader16.py @@ -24,7 +24,6 @@ from MEDLoader import * """ This test is a non regression test of EDF8662 : This bug revealed that ELNOMesh and ELNOPoints do not behave correctly after the call of ExtractGroup""" fname="testMEDReader16.med" -outImgName="testMEDReader16.png" arr=DataArrayDouble([0,1,2]) m=MEDCouplingCMesh() ; m.setCoords(arr,arr) ; m=m.buildUnstructured() ; m.setName("Mesh") diff --git a/src/Plugins/MEDReader/Test/testMEDReader17.py b/src/Plugins/MEDReader/Test/testMEDReader17.py new file mode 100644 index 00000000..e9db394a --- /dev/null +++ b/src/Plugins/MEDReader/Test/testMEDReader17.py @@ -0,0 +1,68 @@ +# -*- coding: iso-8859-1 -*- +# Copyright (C) 2007-2014 CEA/DEN, EDF R&D +# +# This library is free software; you can redistribute it and/or +# modify it under the terms of the GNU Lesser General Public +# License as published by the Free Software Foundation; either +# version 2.1 of the License. +# +# This library is distributed in the hope that it will be useful, +# but WITHOUT ANY WARRANTY; without even the implied warranty of +# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU +# Lesser General Public License for more details. +# +# You should have received a copy of the GNU Lesser General Public +# License along with this library; if not, write to the Free Software +# Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA +# +# See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com +# +# Author : Anthony Geay (EDF R&D) + +from MEDLoader import * + +""" This test is a non regression test of EDF9622.""" + +fname="testMEDReader17.med" + +arr=DataArrayDouble([0,1,2]) +m=MEDCouplingCMesh() ; m.setCoords(arr,arr) ; m=m.buildUnstructured() ; m.setName("Mesh") +mm=MEDFileUMesh() ; mm.setMeshAtLevel(0,m) +grp0=DataArrayInt([0,1]) ; grp0.setName("grp0") +grp1=DataArrayInt([2,3]) ; grp1.setName("grp1") +grp2=DataArrayInt([1,2]) ; grp2.setName("grp2") +grp3=DataArrayInt([0,1,2,3]) ; grp3.setName("grp3") +mm.setGroupsAtLevel(0,[grp0,grp1,grp2,grp3]) +f=MEDCouplingFieldDouble(ON_GAUSS_NE) ; f.setMesh(m) ; f.setName("MyField") ; f.setTime(0.,0,0) +arr2=DataArrayDouble(4*4*2) ; arr2.iota() ; arr2.rearrange(2) ; arr2.setInfoOnComponents(["aa","bbb"]) +f.setArray(arr2) ; arr2+=0.1 ; f.checkCoherency() +mm.write(fname,2) +MEDLoader.WriteFieldUsingAlreadyWrittenMesh(fname,f) +# +from paraview.simple import * +from paraview import servermanager + +paraview.simple._DisableFirstRenderCameraReset() +reader=MEDReader(FileName=fname) +ExpectedEntries=['TS0/Mesh/ComSup0/MyField@@][@@GSSNE','TS1/Mesh/ComSup0/Mesh@@][@@P0'] +assert(reader.GetProperty("FieldsTreeInfo")[::2]==ExpectedEntries) +reader.AllArrays=['TS0/Mesh/ComSup0/MyField@@][@@GSSNE'] +ExtractGroup1 = ExtractGroup(Input=reader) +ExtractGroup1.UpdatePipelineInformation() +ExtractGroup1.AllGroups=["GRP_grp1"] +ELNOPoints1=ELNOPoints(Input=ExtractGroup1) +ELNOPoints1.SelectSourceArray=['ELNO@MyField'] +# +ELNOPoints1=ELNOPoints(Input=reader) +ELNOPoints1.SelectSourceArray=['ELNO@MyField'] +ExtractGroup1 = ExtractGroup(Input=ELNOPoints1) +ExtractGroup1.UpdatePipelineInformation() +ExtractGroup1.AllGroups=["GRP_grp1"] +ExtractGroup1.UpdatePipeline() +res=servermanager.Fetch(ExtractGroup1,0) +assert(res.GetBlock(0).GetNumberOfCells()==8) +vtkArrToTest=res.GetBlock(0).GetPointData().GetArray("MyField") +assert(vtkArrToTest.GetNumberOfComponents()==2) +assert(vtkArrToTest.GetNumberOfTuples()==8) +vals=[vtkArrToTest.GetValue(i) for i in xrange(16)] +assert(DataArrayDouble([(16.1,17.1),(18.1,19.1),(20.1,21.1),(22.1,23.1),(24.1,25.1),(26.1,27.1),(28.1,29.1),(30.1,31.1)]).isEqual(DataArrayDouble(vals,8,2),1e-12))