PeriDyno 1.0.0
Loading...
Searching...
No Matches
Array3D.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() || m_nz != src.size()) {
7 this->resize(src.nx(), src.ny(), src.nz());
8 }
9
10 cuSafeCall(cudaMemcpy2D(m_data.data(), sizeof(T) * m_nx, src.begin(), src.pitch(), sizeof(T) * src.nx(), src.ny() * src.nz(), cudaMemcpyDeviceToHost));
11 }
12
13 template<typename T>
14 class Array3D<T, DeviceType::GPU>
15 {
16 public:
17 Array3D()
18 {};
19
20 Array3D(uint nx, uint ny, uint nz)
21 {
22 this->resize(nx, ny, nz);
23 };
24
28 ~Array3D() { };
29
30 void resize(const uint nx, const uint ny, const uint nz);
31
32 void reset();
33
34 void clear();
35
36 inline T* begin() const { return m_data; }
37
38 DYN_FUNC inline uint nx() const { return m_nx; }
39 DYN_FUNC inline uint ny() const { return m_ny; }
40 DYN_FUNC inline uint nz() const { return m_nz; }
41 DYN_FUNC inline uint pitch() const { return m_pitch_x; }
42
43 DYN_FUNC inline T operator () (const int i, const int j, const int k) const
44 {
45 char* addr = (char*)m_data;
46 addr += (j * m_pitch_x + k * m_nxy);
47 return ((T*)addr)[i];
48 }
49
50 DYN_FUNC inline T& operator () (const int i, const int j, const int k)
51 {
52 char* addr = (char*)m_data;
53 addr += (j * m_pitch_x + k * m_nxy);
54 return ((T*)addr)[i];
55 }
56
57 DYN_FUNC inline T operator [] (const int id) const
58 {
59 return m_data[id];
60 }
61
62 DYN_FUNC inline T& operator [] (const int id)
63 {
64 return m_data[id];
65 }
66
67 DYN_FUNC inline size_t index(const uint i, const uint j, const uint k) const
68 {
69 return i + j * m_nx + k * m_nx * m_ny;
70 }
71
72 DYN_FUNC inline size_t size() const { return m_nx * m_ny * m_nz; }
73 DYN_FUNC inline bool isCPU() const { return false; }
74 DYN_FUNC inline bool isGPU() const { return true; }
75
76 void assign(const Array3D<T, DeviceType::GPU>& src);
77 void assign(const Array3D<T, DeviceType::CPU>& src);
78
79 private:
80 uint m_nx = 0;
81 uint m_pitch_x = 0;
82
83 uint m_ny = 0;
84 uint m_nz = 0;
85 uint m_nxy = 0;
86 T* m_data = nullptr;
87 };
88
89 template<typename T>
91
95
96
97 template<typename T>
98 void Array3D<T, DeviceType::GPU>::resize(const uint nx, const uint ny, const uint nz)
99 {
100 if (NULL != m_data) clear();
101
102 cuSafeCall(cudaMallocPitch((void**)&m_data, (size_t*)&m_pitch_x, (size_t)sizeof(T) * nx, (size_t)ny*nz));
103
104 //TODO: check whether it has problem when m_pitch_x is not divisible by sizeof(T)
105 m_nx = nx; m_ny = ny; m_nz = nz;
106 m_nxy = m_pitch_x * m_ny;
107 }
108
109 template<typename T>
110 void Array3D<T, DeviceType::GPU>::reset()
111 {
112 cuSafeCall(cudaMemset(m_data, 0, m_nxy * m_nz));
113 }
114
115 template<typename T>
117 {
118 if(m_data != nullptr) cuSafeCall(cudaFree(m_data));
119
120 m_data = nullptr;
121 m_nx = 0;
122 m_ny = 0;
123 m_nz = 0;
124 m_nxy = 0;
125 }
126
127 template<typename T>
129 {
130 if (m_nx != src.nx() || m_ny != src.ny() || m_nz != src.nz()) {
131 this->resize(src.nx(), src.ny(), src.nz());
132 }
133
134 cuSafeCall(cudaMemcpy2D(m_data, m_pitch_x, src.begin(), src.pitch(), sizeof(T) *src.nx(), src.ny()*src.nz(), cudaMemcpyDeviceToDevice));
135 }
136
137 template<typename T>
139 {
140 if (m_nx != src.nx() || m_ny != src.ny() || m_nz != src.nz()) {
141 this->resize(src.nx(), src.ny(), src.nz());
142 }
143
144 cuSafeCall(cudaMemcpy2D(m_data, m_pitch_x, src.begin(), sizeof(T) *src.nx(), sizeof(T) *src.nx(), src.ny()*src.nz(), cudaMemcpyHostToDevice));
145 }
146}
#define T(t)
This is an implementation of AdditiveCCD based on peridyno.
Definition Array.h:25
DArray3D< float3 > Grid3f
Definition Array3D.inl:93
Array3D< T, DeviceType::GPU > DArray3D
Definition Array3D.inl:90
DArray3D< bool > Grid1b
Definition Array3D.inl:94
DArray3D< float > Grid1f
Definition Array3D.inl:92
unsigned int uint
Definition VkProgram.h:14