X-Git-Url: http://git.salome-platform.org/gitweb/?a=blobdiff_plain;f=src%2FPlugins%2FMEDReader%2FIO%2FvtkGenerateVectors.cxx;h=a3df978df03350889c0c990032aca83fede6678f;hb=f24ad457e027761c8df0671ba283a27af1ae3511;hp=12327c159a0920d5a60ad06a49b60c63ea45c75b;hpb=866822387420d8da61654fcabf49a9e51308c507;p=modules%2Fparavis.git diff --git a/src/Plugins/MEDReader/IO/vtkGenerateVectors.cxx b/src/Plugins/MEDReader/IO/vtkGenerateVectors.cxx index 12327c15..a3df978d 100644 --- a/src/Plugins/MEDReader/IO/vtkGenerateVectors.cxx +++ b/src/Plugins/MEDReader/IO/vtkGenerateVectors.cxx @@ -1,4 +1,4 @@ -// Copyright (C) 2010-2014 CEA/DEN, EDF R&D +// Copyright (C) 2010-2016 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 @@ -19,18 +19,30 @@ // Author : Anthony Geay #include "vtkGenerateVectors.h" +#include "vtkDataArrayTemplate.h" #include "vtkDoubleArray.h" #include "vtkInformation.h" #include "vtkQuadratureSchemeDefinition.h" +#include "vtkInformationQuadratureSchemeDefinitionVectorKey.h" #include "MEDUtilities.hxx" -#define protected public #include "vtkFieldData.h" +#include + +const char vtkGenerateVectors::VECTOR_SUFFIX[]="_Vector"; + +std::string vtkGenerateVectors::SuffixFieldName(const std::string& name) +{ + std::ostringstream oss; oss << name << VECTOR_SUFFIX; + return oss.str(); +} + void vtkGenerateVectors::Operate(vtkFieldData *fd) { if(!fd) return ; const int nbOfArrs(fd->GetNumberOfArrays()); + std::vector daToAppend; for(int i=0;iGetArray(i)); @@ -43,15 +55,23 @@ void vtkGenerateVectors::Operate(vtkFieldData *fd) if(nbOfCompo<=1 || nbOfCompo==3) continue; if(nbOfCompo==2) - fd->SetArray(i,Operate2Compo(arrc)); + daToAppend.push_back(Operate2Compo(arrc)); else - fd->SetArray(i,OperateMoreThan3Compo(arrc)); + daToAppend.push_back(OperateMoreThan3Compo(arrc)); + } + for(std::vector::const_iterator it=daToAppend.begin();it!=daToAppend.end();it++) + { + vtkDoubleArray *elt(*it); + if(!elt) + continue; + fd->AddArray(elt); + elt->Delete(); } } vtkDoubleArray *vtkGenerateVectors::Operate2Compo(vtkDoubleArray *oldArr) { - static const int VTK_DATA_ARRAY_FREE=vtkDataArrayTemplate::VTK_DATA_ARRAY_FREE; + const int VTK_DATA_ARRAY_FREE=vtkDataArrayTemplate::VTK_DATA_ARRAY_FREE; vtkDoubleArray *ret(vtkDoubleArray::New()); std::size_t nbOfTuples(oldArr->GetNumberOfTuples()); const double *inPt(oldArr->GetPointer(0)); @@ -63,7 +83,8 @@ vtkDoubleArray *vtkGenerateVectors::Operate2Compo(vtkDoubleArray *oldArr) pt[3*i+2]=0.; } ret->SetNumberOfComponents(3); - ret->SetName(oldArr->GetName()); + std::string newName(SuffixFieldName(oldArr->GetName())); + ret->SetName(newName.c_str()); ret->SetComponentName(0,oldArr->GetComponentName(0)); ret->SetComponentName(1,oldArr->GetComponentName(1)); ret->SetArray(pt,3*nbOfTuples,0,VTK_DATA_ARRAY_FREE); @@ -73,7 +94,7 @@ vtkDoubleArray *vtkGenerateVectors::Operate2Compo(vtkDoubleArray *oldArr) vtkDoubleArray *vtkGenerateVectors::OperateMoreThan3Compo(vtkDoubleArray *oldArr) { - static const int VTK_DATA_ARRAY_FREE=vtkDataArrayTemplate::VTK_DATA_ARRAY_FREE; + const int VTK_DATA_ARRAY_FREE=vtkDataArrayTemplate::VTK_DATA_ARRAY_FREE; vtkDoubleArray *ret(vtkDoubleArray::New()); int nbOfCompo(oldArr->GetNumberOfComponents()); std::size_t nbOfTuples(oldArr->GetNumberOfTuples()); @@ -86,7 +107,8 @@ vtkDoubleArray *vtkGenerateVectors::OperateMoreThan3Compo(vtkDoubleArray *oldArr pt[3*i+2]=inPt[nbOfCompo*i+2]; } ret->SetNumberOfComponents(3); - ret->SetName(oldArr->GetName()); + std::string newName(SuffixFieldName(oldArr->GetName())); + ret->SetName(newName.c_str()); ret->SetComponentName(0,oldArr->GetComponentName(0)); ret->SetComponentName(1,oldArr->GetComponentName(1)); ret->SetComponentName(2,oldArr->GetComponentName(2)); @@ -103,6 +125,15 @@ void vtkGenerateVectors::UpdateInformationOfArray(vtkDoubleArray *oldArr, vtkDou } if(oldArr->GetInformation()->Get(MEDUtilities::ELGA())) arr->GetInformation()->Set(MEDUtilities::ELGA(),1); + vtkInformationQuadratureSchemeDefinitionVectorKey *key(vtkQuadratureSchemeDefinition::DICTIONARY()); + if(key->Has(oldArr->GetInformation())) + { + int dictSize(key->Size(oldArr->GetInformation())); + vtkQuadratureSchemeDefinition **dict = new vtkQuadratureSchemeDefinition *[dictSize]; + key->GetRange(oldArr->GetInformation(),dict,0,0,dictSize); + key->SetRange(arr->GetInformation(),dict,0,0,dictSize); + delete [] dict; + } if(oldArr->GetInformation()->Get(MEDUtilities::ELNO())) arr->GetInformation()->Set(MEDUtilities::ELNO(),1); }