1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308
|
/*=========================================================================
Program: ParaView
Module: vtkCPPythonScriptPipeline.cxx
Copyright (c) Kitware, Inc.
All rights reserved.
See Copyright.txt or http://www.paraview.org/HTML/Copyright.html for details.
This software is distributed WITHOUT ANY WARRANTY; without even
the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
PURPOSE. See the above copyright notice for more information.
=========================================================================*/
#include "vtkCPPythonScriptPipeline.h"
#include "vtkCPDataDescription.h"
#include "vtkDataObject.h"
#include "vtkMultiProcessController.h"
#include "vtkNew.h"
#include "vtkObjectFactory.h"
#include "vtkProcessModule.h"
#include "vtkPVConfig.h"
#include "vtkPVPythonOptions.h"
#include "vtkPythonInterpreter.h"
#include "vtkSMObject.h"
#include "vtkSMProxyManager.h"
#include <string>
#include <vtksys/SystemTools.hxx>
#include <sstream>
extern "C" {
void vtkPVInitializePythonModules();
}
namespace
{
//----------------------------------------------------------------------------
void InitializePython()
{
static bool initialized = false;
if (initialized)
{
return;
}
initialized = true;
// register callback to initialize modules statically. The callback is
// empty when BUILD_SHARED_LIBS is ON.
vtkPVInitializePythonModules();
vtkPythonInterpreter::Initialize();
std::ostringstream loadPythonModules;
loadPythonModules
<< "import sys\n"
<< "import paraview\n"
<< "f1 = paraview.print_error\n"
<< "f2 = paraview.print_debug_info\n"
<< "def print_dummy(text):\n"
<< " pass\n"
<< "paraview.print_error = print_dummy\n"
<< "paraview.print_debug_info = print_dummy\n"
<< "paraview.print_error = f1\n"
<< "paraview.print_debug_info = f2\n"
<< "import vtkPVCatalystPython\n";
vtkPythonInterpreter::RunSimpleString(loadPythonModules.str().c_str());
}
//----------------------------------------------------------------------------
// for things like programmable filters that have a '\n' in their strings,
// we need to fix them to have \\n so that everything works smoothly
void fixEOL(std::string& str)
{
const std::string from = "\\n";
const std::string to = "\\\\n";
size_t start_pos = 0;
while((start_pos = str.find(from, start_pos)) != std::string::npos)
{
str.replace(start_pos, from.length(), to);
start_pos += to.length();
}
return;
}
}
vtkStandardNewMacro(vtkCPPythonScriptPipeline);
//----------------------------------------------------------------------------
vtkCPPythonScriptPipeline::vtkCPPythonScriptPipeline()
{
this->PythonScriptName = 0;
}
//----------------------------------------------------------------------------
vtkCPPythonScriptPipeline::~vtkCPPythonScriptPipeline()
{
this->SetPythonScriptName(0);
}
//----------------------------------------------------------------------------
int vtkCPPythonScriptPipeline::Initialize(const char* fileName)
{
// only process 0 checks if the file exists and broadcasts that information
// to the other processes
int fileExists = 0;
vtkMultiProcessController* controller =
vtkMultiProcessController::GetGlobalController();
if(controller->GetLocalProcessId()==0)
{
fileExists = vtksys::SystemTools::FileExists(fileName, true);
}
controller->Broadcast(&fileExists, 1, 0);
if(fileExists == 0)
{
vtkErrorMacro("Could not find file " << fileName);
return 0;
}
InitializePython();
// for now do not check on filename extension:
//vtksys::SystemTools::GetFilenameLastExtension(FileName) == ".py" == 0)
std::string fileNamePath = vtksys::SystemTools::GetFilenamePath(fileName);
std::string fileNameName = vtksys::SystemTools::GetFilenameWithoutExtension(
vtksys::SystemTools::GetFilenameName(fileName));
// need to save the script name as it is used as the name of the module
this->SetPythonScriptName(fileNameName.c_str());
// only process 0 reads the actual script and then broadcasts it out
char* scriptText = NULL;
// we need to add the script path to PYTHONPATH
char* scriptPath = NULL;
int rank = controller->GetLocalProcessId();
int scriptSizes[2] = {0, 0};
if(rank == 0)
{
std::string line;
std::ifstream myfile (fileName);
std::string desiredString;
if (myfile.is_open())
{
while ( getline (myfile,line) )
{
fixEOL(line);
desiredString.append(line).append("\n");
}
myfile.close();
}
if(fileNamePath.empty())
{
fileNamePath = ".";
}
scriptSizes[0] = static_cast<int>(fileNamePath.size()+1);
scriptPath = new char[scriptSizes[0]];
memcpy(scriptPath, fileNamePath.c_str(), sizeof(char)*scriptSizes[0]);
scriptSizes[1] = static_cast<int>(desiredString.size()+1);
scriptText = new char[scriptSizes[1]];
memcpy(scriptText, desiredString.c_str(), sizeof(char)*scriptSizes[1]);
}
controller->Broadcast(scriptSizes, 2, 0);
if (rank != 0)
{
scriptPath = new char[scriptSizes[0]];
scriptText = new char[scriptSizes[1]];
}
controller->Broadcast(scriptPath, scriptSizes[0], 0);
controller->Broadcast(scriptText, scriptSizes[1], 0);
vtkPythonInterpreter::PrependPythonPath(scriptPath);
// The code below creates a module from the scriptText string.
// This requires the manual creation of a module object like this:
//
// import types
// _foo = types.ModuleType('foo')
// _foo.__file__ = 'foo.pyc'
// import sys
// sys.module['foo'] = _foo
// _source= scriptText
// _code = compile(_source, 'foo.py', 'exec')
// exec _code in _foo.__dict__
// del _source
// del _code
// import foo
std::ostringstream loadPythonModules;
loadPythonModules << "import types" << std::endl;
loadPythonModules << "_" << fileNameName << " = types.ModuleType('" << fileNameName << "')" << std::endl;
loadPythonModules << "_" << fileNameName << ".__file__ = '" << fileNameName << ".pyc'" << std::endl;
loadPythonModules << "import sys" << std::endl;
loadPythonModules << "sys.modules['" << fileNameName << "'] = _" << fileNameName << std::endl;
loadPythonModules << "_source = \"\"\"" << std::endl;
loadPythonModules << scriptText;
loadPythonModules << "\"\"\"" << std::endl;
loadPythonModules << "_code = compile(_source, \"" << fileNameName << ".py\", \"exec\")" << std::endl;
loadPythonModules << "exec _code in _" << fileNameName << ".__dict__" << std::endl;
loadPythonModules << "del _source" << std::endl;
loadPythonModules << "del _code" << std::endl;
loadPythonModules << "import " << fileNameName << std::endl;
delete[] scriptPath;
delete[] scriptText;
vtkPythonInterpreter::RunSimpleString(loadPythonModules.str().c_str());
return 1;
}
//----------------------------------------------------------------------------
int vtkCPPythonScriptPipeline::RequestDataDescription(
vtkCPDataDescription* dataDescription)
{
if(!dataDescription)
{
vtkWarningMacro("dataDescription is NULL.");
return 0;
}
InitializePython();
// check the script to see if it should be run...
vtkStdString dataDescriptionString = this->GetPythonAddress(dataDescription);
std::ostringstream pythonInput;
pythonInput << "dataDescription = vtkPVCatalystPython.vtkCPDataDescription('"
<< dataDescriptionString << "')\n"
<< this->PythonScriptName << ".RequestDataDescription(dataDescription)\n";
vtkPythonInterpreter::RunSimpleString(pythonInput.str().c_str());
return dataDescription->GetIfAnyGridNecessary()? 1: 0;
}
//----------------------------------------------------------------------------
int vtkCPPythonScriptPipeline::CoProcess(
vtkCPDataDescription* dataDescription)
{
if(!dataDescription)
{
vtkWarningMacro("DataDescription is NULL.");
return 0;
}
InitializePython();
vtkStdString dataDescriptionString = this->GetPythonAddress(dataDescription);
std::ostringstream pythonInput;
pythonInput
<< "dataDescription = vtkPVCatalystPython.vtkCPDataDescription('"
<< dataDescriptionString << "')\n"
<< this->PythonScriptName << ".DoCoProcessing(dataDescription)\n";
vtkPythonInterpreter::RunSimpleString(pythonInput.str().c_str());
return 1;
}
//----------------------------------------------------------------------------
int vtkCPPythonScriptPipeline::Finalize()
{
InitializePython();
std::ostringstream pythonInput;
pythonInput
<< "if hasattr(" << this->PythonScriptName << ", 'Finalize'):\n"
<< " " << this->PythonScriptName << ".Finalize()\n";
vtkPythonInterpreter::RunSimpleString(pythonInput.str().c_str());
return 1;
}
//----------------------------------------------------------------------------
vtkStdString vtkCPPythonScriptPipeline::GetPythonAddress(void* pointer)
{
char addressOfPointer[1024];
#ifdef COPROCESSOR_WIN32_BUILD
sprintf_s(addressOfPointer, "%p", pointer);
#else
sprintf(addressOfPointer, "%p", pointer);
#endif
char *aplus = addressOfPointer;
if ((addressOfPointer[0] == '0') &&
((addressOfPointer[1] == 'x') || addressOfPointer[1] == 'X'))
{
aplus += 2; //skip over "0x"
}
vtkStdString value = aplus;
return value;
}
//----------------------------------------------------------------------------
void vtkCPPythonScriptPipeline::PrintSelf(ostream& os, vtkIndent indent)
{
this->Superclass::PrintSelf(os, indent);
os << indent << "PythonScriptName: " << this->PythonScriptName << "\n";
}
|