Blame view

stim/optics_old/halfspace.cuh 8.83 KB
ecfd14df   David Mayerich   implemented D fie...
1
2
3
4
  #ifndef	RTS_HALFSPACE_H
  #define	RTS_HALFSPACE_H
  
  #include "../math/plane.h"
559d0fcb   David Mayerich   wave interactions...
5
  
559d0fcb   David Mayerich   wave interactions...
6
  
8a86bd56   David Mayerich   changed rts names...
7
  namespace stim{
559d0fcb   David Mayerich   wave interactions...
8
  
ecfd14df   David Mayerich   implemented D fie...
9
  //GPU kernel to compute the electric field
559d0fcb   David Mayerich   wave interactions...
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
  template<typename T>
  __global__ void gpu_halfspace_pw2ef(complex<T>* X, complex<T>* Y, complex<T>* Z, unsigned int r0, unsigned int r1, 
  									 plane<T> P, planewave<T> w, rect<T> q, bool at_surface = false)
  {
      int iu = blockIdx.x * blockDim.x + threadIdx.x;
      int iv = blockIdx.y * blockDim.y + threadIdx.y;
  
      //make sure that the thread indices are in-bounds
      if(iu >= r0 || iv >= r1) return;
  
      //compute the index into the field
      int i = iv*r0 + iu;
  
  	//get the current position
  	vec<T> p = q( (T)iu/(T)r0, (T)iv/(T)r1 );
  
  	if(at_surface){
  		if(P.side(p) > 0)
  			return;
  	}
  	else{
  		if(P.side(p) >= 0)
  			return;
  	}
  
  	//if the current position is on the wrong side of the plane
  
  	//vec<T> r(p[0], p[1], p[2]);
  
  	complex<T> x( 0.0f, w.kvec().dot(p) );
ecfd14df   David Mayerich   implemented D fie...
40
  	//complex<T> phase( 0.0f, w.phase());
559d0fcb   David Mayerich   wave interactions...
41
42
43
44
45
46
  
      if(Y == NULL)                       //if this is a scalar simulation
          X[i] += w.E().len() * exp(x);    //use the vector magnitude as the plane wave amplitude
      else
      {
      	//X[i] = Y[i] = Z[i] = 1;
ecfd14df   David Mayerich   implemented D fie...
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
          X[i] += w.E()[0] * exp(x);// * exp(phase);
          Y[i] += w.E()[1] * exp(x);// * exp(phase);
          Z[i] += w.E()[2] * exp(x);// * exp(phase);
      }
  }
  
  //GPU kernel to compute the electric displacement field
  template<typename T>
  __global__ void gpu_halfspace_pw2df(complex<T>* X, complex<T>* Y, complex<T>* Z, unsigned int r0, unsigned int r1, 
  									 plane<T> P, planewave<T> w, rect<T> q, T n, bool at_surface = false)
  {
      int iu = blockIdx.x * blockDim.x + threadIdx.x;
      int iv = blockIdx.y * blockDim.y + threadIdx.y;
  
      //make sure that the thread indices are in-bounds
      if(iu >= r0 || iv >= r1) return;
  
      //compute the index into the field
      int i = iv*r0 + iu;
  
  	//get the current position
  	vec<T> p = q( (T)iu/(T)r0, (T)iv/(T)r1 );
  
  	if(at_surface){
  		if(P.side(p) > 0)
  			return;
  	}
  	else{
  		if(P.side(p) >= 0)
  			return;
  	}
  
  	//if the current position is on the wrong side of the plane
  
  	//vec<T> r(p[0], p[1], p[2]);
  
  	complex<T> x( 0.0f, w.kvec().dot(p) );
  	//complex<T> phase( 0.0f, w.phase());
  
  	//vec< complex<T> > testE(1, 0, 0);
  
  	
  
      if(Y == NULL)                       //if this is a scalar simulation
          X[i] += w.E().len() * exp(x);    //use the vector magnitude as the plane wave amplitude
      else
      {
      	plane< complex<T> > cplane = plane< complex<T>, 3>(P);
      	vec< complex<T> > E_para;// = cplane.parallel(w.E());
  		vec< complex<T> > E_perp;// = cplane.perpendicular(w.E()) * (n*n);
  		cplane.decompose(w.E(), E_para, E_perp);
  		T epsilon = n*n;
  
          X[i] += (E_para[0] + E_perp[0] * epsilon) * exp(x);
          Y[i] += (E_para[1] + E_perp[1] * epsilon) * exp(x);
          Z[i] += (E_para[2] + E_perp[2] * epsilon) * exp(x);
559d0fcb   David Mayerich   wave interactions...
103
104
      }
  }
ecfd14df   David Mayerich   implemented D fie...
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
  
  //computes a scalar field containing the refractive index of the half-space at each point
  template<typename T>
  __global__ void gpu_halfspace_n(T* n, unsigned int r0, unsigned int r1, rect<T> q, plane<T> P, T n0, T n1){
  
  	int iu = blockIdx.x * blockDim.x + threadIdx.x;
      int iv = blockIdx.y * blockDim.y + threadIdx.y;
  
      //make sure that the thread indices are in-bounds
      if(iu >= r0 || iv >= r1) return;
  
      //compute the index into the field
      int i = iv*r0 + iu;
  
  	//get the current position
  	vec<T> p = q( (T)iu/(T)r0, (T)iv/(T)r1 );
  
  	//set the appropriate refractive index
  	if(P.side(p) < 0) n[i] = n0;
  	else n[i] = n1;
559d0fcb   David Mayerich   wave interactions...
125
126
127
  }
  
  template<class T>
ecfd14df   David Mayerich   implemented D fie...
128
129
130
131
132
133
134
135
136
137
  class halfspace
  {
  private:
  	rts::plane<T> S;		//surface plane splitting the half space
  	rts::complex<T> n0;		//refractive index at the front of the plane
  	rts::complex<T> n1;		//refractive index at the back of the plane
  
  	//lists of waves in front (pw0) and behind (pw1) the plane
  	std::vector< rts::planewave<T> > w0;
  	std::vector< rts::planewave<T> > w1;
559d0fcb   David Mayerich   wave interactions...
138
  
ecfd14df   David Mayerich   implemented D fie...
139
140
141
  	//rts::planewave<T> pi;	//incident plane wave
  	//rts::planewave<T> pr;	//plane wave reflected from the surface
  	//rts::planewave<T> pt;	//plane wave transmitted through the surface
559d0fcb   David Mayerich   wave interactions...
142
  
ecfd14df   David Mayerich   implemented D fie...
143
144
145
146
  	void init(){
  		n0 = 1.0;
  		n1 = 1.5;
  	}
559d0fcb   David Mayerich   wave interactions...
147
  
ecfd14df   David Mayerich   implemented D fie...
148
149
150
151
152
153
154
  	//compute which side of the interface is hit by the incoming plane wave (0 = front, 1 = back)
  	bool facing(planewave<T> p){
  		if(p.kvec().dot(S.norm()) > 0)
  			return 1;
  		else
  			return 0;
  	}
559d0fcb   David Mayerich   wave interactions...
155
  
ecfd14df   David Mayerich   implemented D fie...
156
  	T calc_theta_i(vec<T> v){
559d0fcb   David Mayerich   wave interactions...
157
  
ecfd14df   David Mayerich   implemented D fie...
158
  		vec<T> v_hat = v.norm();
559d0fcb   David Mayerich   wave interactions...
159
  
ecfd14df   David Mayerich   implemented D fie...
160
161
162
163
  		//compute the cosine of the angle between k_hat and the plane normal
  		T cos_theta_i = v_hat.dot(S.norm());
  
  		return acos(abs(cos_theta_i));
559d0fcb   David Mayerich   wave interactions...
164
165
  	}
  
ecfd14df   David Mayerich   implemented D fie...
166
167
168
169
  	T calc_theta_t(T ni_nt, T theta_i){
  
  		T sin_theta_t = ni_nt * sin(theta_i);
  		return asin(sin_theta_t);
559d0fcb   David Mayerich   wave interactions...
170
  	}
559d0fcb   David Mayerich   wave interactions...
171
172
  
  
ecfd14df   David Mayerich   implemented D fie...
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
  public:
  
  	//constructors
  	halfspace(){
  		init();
  	}
  
  	halfspace(T na, T nb){
  		init();
  		n0 = na;
  		n1 = nb;
  	}
  
  	halfspace(T na, T nb, plane<T> p){
  		n0 = na;
  		n1 = nb;
  		S = p;
  	}
  
  	//compute the transmitted and reflective waves given the incident (vacuum) plane wave p
  	void incident(rts::planewave<T> p){
  
  		planewave<T> r, t;
  		p.scatter(S, n1.real()/n0.real(), r, t);
  
  		//std::cout<<"i+r: "<<p.r()[0] + r.r()[0]<<std::endl;
  		//std::cout<<"t:   "<<t.r()[0]<<std::endl;
  
  		if(facing(p)){
  			w1.push_back(p);
  
  			if(r.E().len() != 0)
  				w1.push_back(r);
  			if(t.E().len() != 0)
  				w0.push_back(t);
  		}
  		else{
  			w0.push_back(p);
  
  			if(r.E().len() != 0)
  				w0.push_back(r);
  			if(t.E().len() != 0)
  				w1.push_back(t);
  		}
  	}
  
  	void incident(rts::beam<T> b, unsigned int N = 10000){
  
  		//generate a plane wave decomposition for the beam
  		std::vector< planewave<T> > pw_list = b.mc(N);
  
  		//calculate the reflected and refracted waves for each incident wave
  		for(unsigned int w = 0; w < pw_list.size(); w++){
  			incident(pw_list[w]);
  		}
  	}
  
  	//return the electric field at the specified resolution and position
  	rts::efield<T> E(unsigned int r0, unsigned int r1, rect<T> R){
  		
  		int maxThreads = rts::maxThreadsPerBlock(); //compute the optimal block size
  		int SQRT_BLOCK = (int)std::sqrt((float)maxThreads);
  
  		//create one thread for each detector pixel
  		dim3 dimBlock(SQRT_BLOCK, SQRT_BLOCK);
  		dim3 dimGrid((r0 + SQRT_BLOCK -1)/SQRT_BLOCK, (r1 + SQRT_BLOCK - 1)/SQRT_BLOCK);
  
  		//create an electric field
  		rts::efield<T> ef(r0, r1);
  		ef.position(R);
559d0fcb   David Mayerich   wave interactions...
243
  
ecfd14df   David Mayerich   implemented D fie...
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
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
  		//render each plane wave
  
  		//plane waves at the surface front
  		for(unsigned int w = 0; w < w0.size(); w++){
  			//std::cout<<"w0 "<<w<<": "<<hs.w0[w].str()<<std::endl;
  			rts::gpu_halfspace_pw2ef<T> <<<dimGrid, dimBlock>>> (ef.x(), ef.y(), ef.z(), r0, r1, S, w0[w], ef.p());
  		}
  
  		//plane waves at the surface back
  		for(unsigned int w = 0; w < w1.size(); w++){
  			//std::cout<<"w1 "<<w<<": "<<hs.w1[w].str()<<std::endl;
  			rts::gpu_halfspace_pw2ef<T> <<<dimGrid, dimBlock>>> (ef.x(), ef.y(), ef.z(), r0, r1, -S, w1[w], ef.p(), true);
  		}
  
  		return ef;
  	}
  
  	//return the electric displacement at the specified resolution and position
  	rts::efield<T> D(unsigned int r0, unsigned int r1, rect<T> R){
  
  		int maxThreads = rts::maxThreadsPerBlock(); //compute the optimal block size
  		int SQRT_BLOCK = (int)std::sqrt((float)maxThreads);
  
  		//create one thread for each detector pixel
  		dim3 dimBlock(SQRT_BLOCK, SQRT_BLOCK);
  		dim3 dimGrid((r0 + SQRT_BLOCK -1)/SQRT_BLOCK, (r1 + SQRT_BLOCK - 1)/SQRT_BLOCK);
  		
  		//create a complex vector field
  		rts::efield<T> df(r0, r1);
  		df.position(R);
  		
  		//render each plane wave
  
  		//plane waves at the surface front
  		for(unsigned int w = 0; w < w0.size(); w++){
  			//std::cout<<"w0 "<<w<<": "<<hs.w0[w].str()<<std::endl;
  			rts::gpu_halfspace_pw2df<T> <<<dimGrid, dimBlock>>> (df.x(), df.y(), df.z(), r0, r1, S, w0[w], df.p(), n0.real());
  		}
  		
  		//plane waves at the surface back
  		for(unsigned int w = 0; w < w1.size(); w++){
  			//std::cout<<"w1 "<<w<<": "<<hs.w1[w].str()<<std::endl;
  			rts::gpu_halfspace_pw2df<T> <<<dimGrid, dimBlock>>> (df.x(), df.y(), df.z(), r0, r1, -S, w1[w], df.p(), n1.real(), true);
  		}
  		
  		return df;
  	}
  
  	realfield<T, 1, true> nfield(unsigned int Ru, unsigned int Rv, rect<T> p){
  
  		int maxThreads = rts::maxThreadsPerBlock(); //compute the optimal block size
  		int SQRT_BLOCK = (int)std::sqrt((float)maxThreads);
  
  		//create one thread for each detector pixel
  		dim3 dimBlock(SQRT_BLOCK, SQRT_BLOCK);
  		dim3 dimGrid((Ru + SQRT_BLOCK -1)/SQRT_BLOCK, (Rv + SQRT_BLOCK - 1)/SQRT_BLOCK);
  
  		realfield<T, 1, true> n(Ru, Rv);	//create a scalar field to store the result of n
  
  		rts::gpu_halfspace_n<T> <<<dimGrid, dimBlock>>> (n.ptr(), Ru, Rv, p, S, n0.real(), n1.real());
  
  		return n;
  	}
  
  	std::string str(){
  		std::stringstream ss;
  		ss<<"Half Space-------------"<<std::endl;
  		ss<<"n0: "<<n0<<std::endl;
  		ss<<"n1: "<<n1<<std::endl;
  		ss<<S.str();
  		
  
  		if(w0.size() != 0 || w1.size() != 0){
  			ss<<std::endl<<"Plane Waves:"<<std::endl;
  			for(unsigned int w = 0; w < w0.size(); w++){
  				ss<<"w0 = "<<w<<": "<<w0[w]<<std::endl;
  			}
  			for(unsigned int w = 0; w < w1.size(); w++){
  				ss<<"w1 = "<<w<<": "<<w1[w]<<std::endl;
  			}
  		}
  
  		return ss.str();
  	}
  	////////FRIENDSHIP
      //friend void operator<< <> (rts::efield<T> &ef, rts::halfspace<T> hs);
  
  };
  
  
  
  
  }
559d0fcb   David Mayerich   wave interactions...
337
338
  
  
8a86bd56   David Mayerich   changed rts names...
339
  #endif