PeriDyno 1.0.0
Loading...
Searching...
No Matches
Array2D.inl
Go to the documentation of this file.
1namespace dyno {
2
3 template<typename T>
5 {
6 if (m_nx != src.size() || m_ny != src.size()) {
7 this->resize(src.nx(), src.ny());
8 }
9
10 cuSafeCall(cudaMemcpy2D(m_data.data(), sizeof(T) * m_nx, src.begin(), src.pitch(), sizeof(T) * src.nx(), src.ny(), cudaMemcpyDeviceToHost));
11 }
12
13 template<typename T>
14 class Array2D<T, DeviceType::GPU>
15 {
16 public:
17 Array2D() {};
18
19 Array2D(uint nx, uint ny)
20 {
21 this->resize(nx, ny);
22 };
23
27 ~Array2D() {};
28
29 void resize(uint nx, uint ny);
30
31 void reset();
32
33 void clear();
34
35 inline T* begin() const { return m_data; }
36
37 DYN_FUNC inline uint nx() const { return m_nx; }
38 DYN_FUNC inline uint ny() const { return m_ny; }
39 DYN_FUNC inline uint pitch() const { return m_pitch; }
40
41 GPU_FUNC inline T operator () (const uint i, const uint j) const
42 {
43 char* addr = (char*)m_data;
44 addr += j * m_pitch;
45
46 return ((T*)addr)[i];
47 //return m_data[i + j* m_pitch];
48 }
49
50 GPU_FUNC inline T& operator () (const uint i, const uint j)
51 {
52 char* addr = (char*)m_data;
53 addr += j * m_pitch;
54
55 return ((T*)addr)[i];
56
57 //return m_data[i + j* m_pitch];
58 }
59
60 DYN_FUNC inline int index(const uint i, const uint j) const
61 {
62 return i + j * m_nx;
63 }
64
65 GPU_FUNC inline T operator [] (const uint id) const
66 {
67 return m_data[id];
68 }
69
70 GPU_FUNC inline T& operator [] (const uint id)
71 {
72 return m_data[id];
73 }
74
75 DYN_FUNC inline uint size() const { return m_nx * m_ny; }
76 DYN_FUNC inline bool isCPU() const { return false; }
77 DYN_FUNC inline bool isGPU() const { return true; }
78
79 void assign(const Array2D<T, DeviceType::GPU>& src);
80 void assign(const Array2D<T, DeviceType::CPU>& src);
81
82 private:
83 uint m_nx = 0;
84 uint m_ny = 0;
85 uint m_pitch = 0;
86 T* m_data = nullptr;
87 };
88
89 template<typename T>
91
92 template<typename T>
94 {
95 if (nullptr != m_data) clear();
96
97 cuSafeCall(cudaMallocPitch((void**)&m_data, (size_t*)&m_pitch, (size_t)sizeof(T) * nx, (size_t)ny));
98
99 m_nx = nx;
100 m_ny = ny;
101 }
102
103 template<typename T>
104 void Array2D<T, DeviceType::GPU>::reset()
105 {
106 cuSafeCall(cudaMemset((void*)m_data, 0, m_pitch * m_ny));
107 }
108
109 template<typename T>
111 {
112 if (m_data != nullptr)
113 cuSafeCall(cudaFree((void*)m_data));
114
115 m_nx = 0;
116 m_ny = 0;
117 m_pitch = 0;
118 m_data = nullptr;
119 }
120
121 template<typename T>
123 {
124 if (m_nx != src.nx() || m_ny != src.ny()){
125 this->resize(src.nx(), src.ny());
126 }
127
128 cuSafeCall(cudaMemcpy2D(m_data, m_pitch, src.begin(), src.pitch(), sizeof(T) * src.nx(), src.ny(), cudaMemcpyDeviceToDevice));
129 }
130
131 template<typename T>
133 {
134 if (m_nx != src.nx() || m_ny != src.ny()) {
135 this->resize(src.nx(), src.ny());
136 }
137
138 cuSafeCall(cudaMemcpy2D(m_data, m_pitch, src.begin(), sizeof(T) *src.nx(), sizeof(T) * src.nx(), src.ny(), cudaMemcpyHostToDevice));
139 }
140}
#define T(t)
This is an implementation of AdditiveCCD based on peridyno.
Definition Array.h:25
Array2D< T, DeviceType::GPU > DArray2D
Definition Array2D.inl:90
unsigned int uint
Definition VkProgram.h:14