PeriDyno 1.0.0
Loading...
Searching...
No Matches
HeightField.cu
Go to the documentation of this file.
1#include "HeightField.h"
2#include <fstream>
3#include <iostream>
4#include <sstream>
5
6#include "Math/Lerp.h"
7
8namespace dyno
9{
10 IMPLEMENT_TCLASS(HeightField, TDataType)
11
12 template<typename TDataType>
13 HeightField<TDataType>::HeightField()
14 : TopologyModule()
15 {
16 mDisplacement.resize(128, 128);
17 mDisplacement.reset();
18
19 mGridSpacing = Real(0.1);
20
21 mOrigin = Coord(-mGridSpacing * 64, Real(0), -mGridSpacing * 64);
22 }
23
24 template<typename TDataType>
25 HeightField<TDataType>::~HeightField()
26 {
27 mDisplacement.clear();
28 }
29
30 template<typename TDataType>
31 void HeightField<TDataType>::setExtents(uint nx, uint ny)
32 {
33 mDisplacement.resize(nx, ny);
34 mHeights.resize(nx, ny);
35
36 mDisplacement.reset();
37 mHeights.reset();
38 }
39
40 template<typename TDataType>
41 uint HeightField<TDataType>::width()
42 {
43 return mDisplacement.nx();
44 }
45
46 template<typename TDataType>
47 uint HeightField<TDataType>::height()
48 {
49 return mDisplacement.ny();
50 }
51
52 template<typename TDataType>
53 void HeightField<TDataType>::copyFrom(HeightField<TDataType>& hf)
54 {
55 mOrigin = hf.mOrigin;
56 mGridSpacing = hf.mGridSpacing;
57 mDisplacement.assign(hf.mDisplacement);
58 mHeights.assign(hf.mHeights);
59 }
60
61 template <typename Real, typename Coord>
62 __global__ void PS_Scale(
63 DArray<Coord> vertex,
64 Real s)
65 {
66 int pId = threadIdx.x + (blockIdx.x * blockDim.x);
67 if (pId >= vertex.size()) return;
68 //return;
69 vertex[pId] = vertex[pId] * s;
70 }
71
72 template<typename TDataType>
73 void HeightField<TDataType>::scale(Real s)
74 {
75 //cuExecute(m_coords.size(), PS_Scale, m_coords, s);
76 }
77
78 template <typename Coord>
79 __global__ void PS_Scale(
80 DArray<Coord> vertex,
81 Coord s)
82 {
83 int pId = threadIdx.x + (blockIdx.x * blockDim.x);
84 if (pId >= vertex.size()) return;
85
86 Coord pos_i = vertex[pId];
87 vertex[pId] = Coord(pos_i[0] * s[0], pos_i[1] * s[1], pos_i[2] * s[2]);
88 }
89
90 template<typename TDataType>
91 void HeightField<TDataType>::scale(Coord s)
92 {
93 //cuExecute(m_coords.size(), PS_Scale, m_coords, s);
94 }
95
96 template <typename Coord>
97 __global__ void PS_Translate(
98 DArray<Coord> vertex,
99 Coord t)
100 {
101 int pId = threadIdx.x + (blockIdx.x * blockDim.x);
102 if (pId >= vertex.size()) return;
103
104 vertex[pId] = vertex[pId] + t;
105 }
106
107
108 template<typename TDataType>
109 void HeightField<TDataType>::translate(Coord t)
110 {
111 //cuExecute(m_coords.size(), PS_Translate, m_coords, t);
112
113// uint pDims = cudaGridSize(m_coords.size(), BLOCK_SIZE);
114//
115// PS_Translate << <pDims, BLOCK_SIZE >> > (
116// m_coords,
117// t);
118// cuSynchronize();
119 }
120
121 //Back tracing using the semi-Lagrangian scheme
122 template <typename Real, typename Coord>
123 __global__ void HF_RasterizeDisplacements(
124 DArray2D<Real> vertical,
125 DArray2D<Coord> displacements,
126 Coord origin,
127 Real h)
128 {
129 unsigned int i = blockIdx.x * blockDim.x + threadIdx.x;
130 unsigned int j = blockIdx.y * blockDim.y + threadIdx.y;
131
132 uint nx = displacements.nx();
133 uint ny = displacements.ny();
134
135 if (i < nx && j < ny)
136 {
137 Coord disp_ij = displacements(i, j);
138
139 Coord interp_ij = bilinear(displacements, i - disp_ij.x / h, j - disp_ij.z / h, LerpMode::REPEAT);
140
141 vertical(i, j) = interp_ij.y;
142 }
143 }
144
145 template<typename TDataType>
146 DArray2D<typename TDataType::Real>& HeightField<TDataType>::calculateHeightField()
147 {
148 uint nx = mDisplacement.nx();
149 uint ny = mDisplacement.ny();
150
151 if (nx != mHeights.nx() || ny != mHeights.ny()) {
152 mHeights.resize(nx, ny);
153 }
154
155 cuExecute2D(make_uint2(nx, ny),
156 HF_RasterizeDisplacements,
157 mHeights,
158 mDisplacement,
159 mOrigin,
160 mGridSpacing);
161
162 return mHeights;
163 }
164
165 DEFINE_CLASS(HeightField);
166}