]> SALOME platform Git repositories - modules/paravis.git/commitdiff
Salome HOME
bos #23995: Porting to ParaView 5.9
authorvsr <vsr@opencascade.com>
Sat, 13 Mar 2021 09:53:43 +0000 (12:53 +0300)
committervsr <vsr@opencascade.com>
Sat, 13 Mar 2021 09:53:43 +0000 (12:53 +0300)
src/PVGUI/PVGUI_Module.cxx
src/Plugins/JSONReader/plugin/JSONReaderModule/vtkJSONReader.cxx
src/Plugins/MEDReader/plugin/MEDLoaderForPV/MEDFileFieldRepresentationTree.cxx
src/Plugins/MEDReader/plugin/MEDLoaderForPV/vtkGenerateVectors.cxx
src/Plugins/MEDReader/plugin/MEDReaderIO/vtkExtractCellType.cxx
src/Plugins/MEDReader/plugin/MEDReaderIO/vtkExtractGroup.cxx
src/Plugins/MEDReader/plugin/MEDReaderIO/vtkInformationGaussDoubleVectorKey.h

index 08fc02d20472f52902e9989b5ee7d80af07ef001..a472ce0998c6114901eb7ac11f3f7832472b3a45 100644 (file)
@@ -798,7 +798,7 @@ QString PVGUI_Module::getTraceString()
 
   vtkSMTrace* tracer = vtkSMTrace::GetActiveTracer();
   if ( tracer ) {
-    traceString = tracer->GetCurrentTrace();
+    traceString = tracer->GetCurrentTrace().c_str();
     // 'import pvsimple' is necessary to fix the first call to DisableFirstRenderCamera in the paraview trace
     // 'ShowParaviewView()' ensure there is an opened viewing window (otherwise SEGFAULT!)
     traceString = "import pvsimple" + end_line +
index ff72f761f8bf6a0468622ffebae779a8fa085974..082b322f6fd63f77b18dd8a677611d742e276f57 100644 (file)
@@ -34,7 +34,7 @@
 
 #include <stdexcept>
 #include <sstream>
-
+#include <fstream>
 // Key words
 #define MD  "_metadata"
 #define CMT "_comment"
@@ -183,7 +183,7 @@ int vtkJSONReader::CanParseFile(const char *fname, Json::Value &root)
     return 0;
   }
 
-  ifstream file;
+  std::ifstream file;
   std::ostringstream oss;
   bool parsedSuccess = true;
   Json::Reader reader;
index 61edea6cb5270e07fafb308034a2e2a8e995ce07..6b1857b66d1daaf2b39ab0c9bffe5bbe3b5caf16 100644 (file)
@@ -43,7 +43,7 @@
 #include "vtkInformationQuadratureSchemeDefinitionVectorKey.h"
 #include "vtkInformationIntegerKey.h"
 #include "vtkInformation.h"
-#include "vtkDataArrayTemplate.h"
+#include "vtkAOSDataArrayTemplate.h"
 #include "vtkIdTypeArray.h"
 #include "vtkDoubleArray.h"
 #include "vtkIntArray.h"
@@ -118,7 +118,7 @@ vtkIdTypeArray *ELGACmp::isExisting(const std::vector<std::string>& locsReallyUs
 template<class T>
 vtkIdTypeArray *ELGACmp::createNew(const MEDCoupling::MEDFileFieldGlobsReal *globs, const std::vector<std::string>& locsReallyUsed, vtkDataArray *vtkd, vtkDataSet *ds, ExportedTinyInfo *internalInfo) const
 {
-  const int VTK_DATA_ARRAY_DELETE=vtkDataArrayTemplate<T>::VTK_DATA_ARRAY_DELETE;
+  const int VTK_DATA_ARRAY_DELETE=vtkAOSDataArrayTemplate<T>::VTK_DATA_ARRAY_DELETE;
   std::vector< std::vector<std::string> > locNames(_loc_names);
   std::vector<vtkIdTypeArray *> elgas(_elgas);
   std::vector< std::pair< vtkQuadratureSchemeDefinition *, unsigned char > > defs;
@@ -263,9 +263,9 @@ template<class T>
 void AssignDataPointerToVTK(typename MEDFileVTKTraits<T>::VtkType *vtkTab, typename MEDFileVTKTraits<T>::MCType *mcTab, bool noCpyNumNodes)
 {
   if(noCpyNumNodes)
-    vtkTab->SetArray(mcTab->getPointer(),mcTab->getNbOfElems(),1,vtkDataArrayTemplate<T>::VTK_DATA_ARRAY_FREE);
+    vtkTab->SetArray(mcTab->getPointer(),mcTab->getNbOfElems(),1,vtkAOSDataArrayTemplate<T>::VTK_DATA_ARRAY_FREE);
  else
-   { vtkTab->SetArray(mcTab->getPointer(),mcTab->getNbOfElems(),0,vtkDataArrayTemplate<T>::VTK_DATA_ARRAY_FREE); mcTab->accessToMemArray().setSpecificDeallocator(0); }
+   { vtkTab->SetArray(mcTab->getPointer(),mcTab->getNbOfElems(),0,vtkAOSDataArrayTemplate<T>::VTK_DATA_ARRAY_FREE); mcTab->accessToMemArray().setSpecificDeallocator(0); }
 }
 
 // here copy is always assumed.
@@ -291,7 +291,7 @@ void AssignToFieldData(DataArray *vPtr, const MEDTimeReq *tr, vtkFieldData *att,
                        const std::vector<TypeOfField>& discs, const ELGACmp& elgaCmp, const MEDCoupling::MEDFileFieldGlobsReal *globs,
                        MEDFileAnyTypeField1TS *f1ts, vtkDataSet *ds, ExportedTinyInfo *internalInfo)
 {
-  const int VTK_DATA_ARRAY_DELETE=vtkDataArrayTemplate<T>::VTK_DATA_ARRAY_DELETE;
+  const int VTK_DATA_ARRAY_DELETE=vtkAOSDataArrayTemplate<T>::VTK_DATA_ARRAY_DELETE;
   typename MEDFileVTKTraits<T>::MCType *vi(static_cast<typename MEDFileVTKTraits<T>::MCType *>(vPtr));
   typename MEDFileVTKTraits<T>::VtkType *vtkd(MEDFileVTKTraits<T>::VtkType::New());
   vtkd->SetNumberOfComponents((int)vi->getNumberOfComponents());
index 56cc316466b9102a815c2680329ec0a396b2baa5..53f2f4c626137cdf7d6bce3a7dbd40cd4cc0f4c2 100644 (file)
@@ -19,7 +19,7 @@
 // Author : Anthony Geay
 
 #include "vtkGenerateVectors.h"
-#include "vtkDataArrayTemplate.h"
+#include "vtkAOSDataArrayTemplate.h"
 #include "vtkDoubleArray.h"
 #include "vtkInformation.h"
 #include "vtkUnstructuredGrid.h"
@@ -83,7 +83,7 @@ void vtkGenerateVectors::Operate(vtkFieldData *fd)
 
 vtkDoubleArray *vtkGenerateVectors::Operate2Compo(vtkDoubleArray *oldArr)
 {
-  const int VTK_DATA_ARRAY_FREE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_FREE;
+  const int VTK_DATA_ARRAY_FREE=vtkAOSDataArrayTemplate<double>::VTK_DATA_ARRAY_FREE;
   vtkDoubleArray *ret(vtkDoubleArray::New());
   vtkIdType nbOfTuples(oldArr->GetNumberOfTuples());
   const double *inPt(oldArr->GetPointer(0));
@@ -106,7 +106,7 @@ vtkDoubleArray *vtkGenerateVectors::Operate2Compo(vtkDoubleArray *oldArr)
 
 vtkDoubleArray *vtkGenerateVectors::OperateMoreThan3Compo(vtkDoubleArray *oldArr)
 {
-  const int VTK_DATA_ARRAY_FREE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_FREE;
+  const int VTK_DATA_ARRAY_FREE=vtkAOSDataArrayTemplate<double>::VTK_DATA_ARRAY_FREE;
   vtkDoubleArray *ret(vtkDoubleArray::New());
   int nbOfCompo(oldArr->GetNumberOfComponents());
   vtkIdType nbOfTuples(oldArr->GetNumberOfTuples());
index bb3fc7c992394cb6ccfbee2ab1ec0c00ab409854..03a0b8b6eb2616a1dbdba23d86ebbae3d126dd13 100644 (file)
@@ -23,7 +23,7 @@
 #include "MEDFileFieldOverView.hxx"
 
 #include "vtkAdjacentVertexIterator.h"
-#include "vtkDataArrayTemplate.h"
+#include "vtkAOSDataArrayTemplate.h"
 #include "vtkIntArray.h"
 #include "vtkCellData.h"
 #include "vtkPointData.h"
@@ -371,7 +371,7 @@ int vtkExtractCellType::RequestInformation(vtkInformation * /*request*/, vtkInfo
 
 vtkDataSet *FilterFamilies(vtkDataSet *input, const std::vector<int>& idsToKeep, bool insideOut)
 {
-  const int VTK_DATA_ARRAY_DELETE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_DELETE;
+  const int VTK_DATA_ARRAY_DELETE=vtkAOSDataArrayTemplate<double>::VTK_DATA_ARRAY_DELETE;
   const char ZE_SELECTION_ARR_NAME[]="@@ZeSelection@@";
   vtkDataSet *output(input->NewInstance());
   output->ShallowCopy(input);
index 5996b950c8adc5289fc1720c207d23fd13349685..40a76ea6ae9722e0d8121dbf38a60d2f31ac380c 100644 (file)
@@ -25,7 +25,7 @@
 #include "VTKMEDTraits.hxx"
 
 #include "vtkAdjacentVertexIterator.h"
-#include "vtkDataArrayTemplate.h"
+#include "vtkAOSDataArrayTemplate.h"
 #include "vtkIntArray.h"
 #include "vtkLongArray.h"
 #ifdef WIN32
@@ -135,7 +135,7 @@ vtkDataSet *FilterFamilies(vtkSmartPointer<vtkThreshold>& thres,
                            vtkDataSet *input, const std::set<int>& idsToKeep, bool insideOut, const char *arrNameOfFamilyField,
                            const char *associationForThreshold, bool& catchAll, bool& catchSmth)
 {
-  const int VTK_DATA_ARRAY_DELETE=vtkDataArrayTemplate<double>::VTK_DATA_ARRAY_DELETE;
+  const int VTK_DATA_ARRAY_DELETE=vtkAOSDataArrayTemplate<double>::VTK_DATA_ARRAY_DELETE;
   const char ZE_SELECTION_ARR_NAME[]="@@ZeSelection@@";
   vtkDataSet *output(input->NewInstance());
   output->ShallowCopy(input);
index 47429cd5874b30f98f230a27e6b08be106e9cc90..3c72138f6c9f144c5324fcf427da922c815f5bf4 100644 (file)
@@ -26,7 +26,7 @@ class VTK_EXPORT vtkInformationGaussDoubleVectorKey : public vtkInformationDoubl
 {
 public:
   vtkTypeMacro(vtkInformationGaussDoubleVectorKey, vtkInformationDoubleVectorKey)
-  void PrintSelf(ostream& /*os*/, vtkIndent /*indent*/) VTK_OVERRIDE{}
+  void PrintSelf(ostream& /*os*/, vtkIndent /*indent*/) override {};
 
   vtkInformationGaussDoubleVectorKey(const char* name, const char* location,
     int length = -1) : vtkInformationDoubleVectorKey(name, location, length) { }
@@ -48,7 +48,7 @@ public:
   */
   void CopyDefaultInformation(vtkInformation* /*request*/,
     vtkInformation* fromInfo,
-    vtkInformation* toInfo) VTK_OVERRIDE
+    vtkInformation* toInfo) override
   {
     this->ShallowCopy(fromInfo, toInfo);
   }