24 #ifndef ASLVTKCASTERS_H 25 #define ASLVTKCASTERS_H 30 #include <vtkSmartPointer.h> 31 #include <acl/aclHardware.h> 32 #include <aslGenerators.h> 33 #include <math/aslVectors.h> 62 unsigned int save = 0,
63 const std::string &name =
"");
69 const std::string &name =
"");
76 const std::string &name =
"");
81 const std::string &name);
90 const std::string &name =
"");
93 const std::string &name =
"");
95 vtkSmartPointer<vtkImageData>
castVTKData(
const Block & b);
97 vtkSmartPointer<vtkImageData>
castVTKData(
double *d,
99 unsigned int save = 0,
100 const std::string &name =
"");
102 vtkSmartPointer<vtkImageData>
castVTKData(
double *d1,
105 const std::string &name =
"");
107 vtkSmartPointer<vtkImageData>
castVTKData(
double *d1,
111 const std::string &name =
"");
113 void putToVTKData(
double *d, vtkSmartPointer<vtkImageData> target);
114 void putToVTKData(
double *d1,
double *d2, vtkSmartPointer<vtkImageData> target);
115 void putToVTKData(
double *d1,
double *d2,
double *d3, vtkSmartPointer<vtkImageData> target);
118 vtkSmartPointer<vtkImageData>
castVTKData(
const AbstractData & d,
119 const std::vector<std::string> &names = std::vector<std::string>(0));
122 unsigned int arrayNum = 0,
129 const std::string &name =
"")
136 const std::string &name =
"")
142 template<
typename T>
inline vtkSmartPointer<vtkDataArray>
castVTKDataArray(std::vector<T> & d,
144 const std::string &name =
"")
150 vtkSmartPointer<vtkImageData>
inline castVTKData(std::vector<double> & d,
152 const std::string &name =
"")
166 unsigned int nComponents = 3);
174 #endif // ASLVTKCASTERS_H vtkSmartPointer< vtkIdTypeArray > castVTKIdTypeArray(unsigned int *d0, unsigned int *d1, unsigned int *d2, unsigned int *d3, unsigned int np, const std::string &name="")
creates VTKDataArray with 3 component d1, d2 and d3 and length np and name
Advanced Simulation Library.
Advanced Computational Language.
std::shared_ptr< DataWithGhostNodesACLData > SPDataWithGhostNodesACLData
void putToVTKData(double *d1, double *d2, double *d3, vtkSmartPointer< vtkImageData > target)
void combineArrays(T *d1, T *d2, T *d3, unsigned int size, T *dTarget, unsigned int nComponents=3)
SPDataWithGhostNodesACLData makeData(vtkSmartPointer< vtkImageData > image, unsigned int arrayNum=0, acl::CommandQueue queue=acl::hardware.defaultQueue)
AVec< T > castVTKVector(AVec< T > a, T fill=0)
vtkSmartPointer< vtkDataArray > castVTKDataArray(std::vector< T > &d, unsigned int np, const std::string &name="")
std::shared_ptr< Block > makeBlock(vtkSmartPointer< vtkImageData > image)
std::shared_ptr< ElementBase > Element
vtkSmartPointer< vtkImageData > castVTKData(std::vector< double > &d, const Block &b, const std::string &name="")
vtkSmartPointer< vtkDataArray > castVTKDataArray2in3(T *d1, T *d2, unsigned int np, const std::string &name)
creates VTKDataArray with 3 component d2, d1 and 0 and length np and name