GPUE  v1.0
GPU Gross-Pitaevskii Equation numerical solver for Bose-Einstein condensates
kernels.h
Go to the documentation of this file.
1 //##############################################################################
14 //##############################################################################
15 
16 #ifndef KERNELS_H
17 #define KERNELS_H
18 #include<stdio.h>
19 
24 __device__ double2 subtract(double2 a, double2 b);
25 
30 __device__ double2 add(double2 a, double2 b);
31 
38 __device__ double2 pow(double2 a, int b);
39 
44 __device__ double2 mult(double2 a, double b);
45 
50 __device__ double2 mult(double2 a, double2 b);
51 
56 __global__ void make_cufftDoubleComplex(double *in, double2 *out);
57 
62 __device__ unsigned int getGid3d3d();
63 
68 __device__ unsigned int getBid3d3d();
69 
74 __device__ unsigned int getTid3d3d();
75 
80 __global__ void is_eq(bool *a, bool *b, bool *ans);
81 
82 //##############################################################################
86 //##############################################################################
87 
94 __device__ double complexMagnitude(double2 in);
95 
103 __global__ void complexMultiply(double2 *in1, double2 *in2, double2 *out);
104 
112 __host__ __device__ double2 complexMultiply(double2 in1, double2 in2);
113 
121 __device__ double2 make_complex(double in, int evolution_type);
122 
130 __global__ void complexAbsSum(double2 *in1, double2 *in2, double *out);
131 
136 __global__ void complexMagnitude(double2 *in, double *out);
137 
144 __device__ double complexMagnitudeSquared(double2 in);
145 
150 __global__ void complexMagnitudeSquared(double2 *in, double *out);
151 
156 __global__ void complexMagnitudeSquared(double2 *in, double2 *out);
157 
164 __device__ double2 conjugate(double2 in);
165 
173 __device__ double2 realCompMult(double scalar, double2 comp);
174 
175 //##############################################################################
179 //##############################################################################
180 
188 __global__ void cMult(double2* in1, double2* in2, double2* out);
189 
197 __global__ void cMultPhi(double2* in1, double* in2, double2* out);
198 
210 __global__ void cMultDensity(double2* in1, double2* in2, double2* out, double dt, double mass, int gstate, double gDenConst);
211 
229 __global__ void cMultDensity_ast(EqnNode_gpu *eqn, double2* in, double2* out,
230  double dx, double dy, double dz, double time,
231  int e_num, double dt, double mass, int gstate,
232  double gDenConst);
233 
234 
242 __global__ void pinVortex(double2* in1, double2* in2, double2* out);
243 
244 
252 __global__ void vecMult(double2 *in, double *factor, double2 *out);
253 
258 __global__ void l2_norm(double *in1, double *in2, double *in3, double *out);
259 
264 __global__ void l2_norm(double2 *in1, double2 *in2, double2 *in3, double *out);
265 
270 __global__ void l2_norm(double *in1, double *in2, double *out);
271 
276 __global__ void l2_norm(double2 *in1, double2 *in2, double *out);
277 
285 __global__ void scalarDiv(double2* in, double factor, double2* out);
286 
294 __global__ void scalarDiv(double* in, double factor, double* out);
295 
303 __global__ void scalarMult(double2* in, double factor, double2* out);
304 
312 __global__ void scalarPow(double2* in, double param, double2* out);
313 
320 __global__ void vecConjugate(double2 *in, double2 *out);
321 
326 __global__ void scalarDiv1D(double2*, double2*);
327 
332 __global__ void scalarDiv2D(double2*, double2*);
333 
341 __global__ void scalarDiv_wfcNorm(double2* in, double dr, double* pSum, double2* out);
342 
343 //##############################################################################
344 
351 __global__ void reduce(double2* in, double* out);
352 
360 __global__ void thread_test(double* input, double* output);
361 
369 __global__ void multipass(double2* input, double2* output, int pass);
370 
377 __global__ void multipass(double* input, double* output);
378 
379 //##############################################################################
380 
390 __global__ void angularOp(double omega, double dt, double2* wfc, double* xpyypx, double2* out);
391 
396 __global__ void ast_mult(double *array, double *array_out, EqnNode_gpu *eqn,
397  double dx, double dy, double dz, double time,
398  int element_num);
403 __global__ void ast_cmult(double2 *array, double2 *array_out, EqnNode_gpu *eqn,
404  double dx, double dy, double dz, double time,
405  int element_num);
410 __global__ void ast_op_mult(double2 *array, double2 *array_out,
411  EqnNode_gpu *eqn,
412  double dx, double dy, double dz, double time,
413  int element_num, int evolution_type, double dt);
414 
419 __device__ double2 real_ast(double val, double dt);
420 
425 __device__ double2 im_ast(double val, double dt);
426 
431 __global__ void zeros(bool *in, bool *out);
432 
437 __global__ void set_eq(double *in1, double *in2);
438 
450 __global__ void energyCalc(double2 *wfc, double2 *op, double dt, double2 *energy, int gnd_state, int op_space, double sqrt_omegaz_mass, double gDenConst);
451 
459 inline __device__ double2 braKetMult(double2 in1, double2 in2);
460 
461 
469 __global__ void pSum(double* in1, double* output, int pass);
470 
471 
472 #endif
__device__ double complexMagnitude(double2 in)
Calculates magnitude of complex number. $|a + ib|$.
__device__ double2 braKetMult(double2 in1, double2 in2)
Performs bra-ket state multiplication. Not fully implemented.
% % % starting wavefunction wfc
Definition: GPE_2d.m:52
__device__ double2 conjugate(double2 in)
Returns conjugate of the a complex number.
__device__ double2 pow(double2 a, int b)
power operation for a double2
__device__ double2 add(double2 a, double2 b)
addition operation for 2 double2 values
__device__ double2 realCompMult(double scalar, double2 comp)
Multiply real scalar by a complex number.
tuple dy
Definition: en.py:60
__global__ void set_eq(double *in1, double *in2)
Sets in2 to be equal to in1.
__device__ double2 make_complex(double in, int evolution_type)
Transforms field value into operator.
__global__ void cMultPhi(double2 *in1, double *in2, double2 *out)
Kernel for multiplcation with real array and complex array.
__global__ void complexAbsSum(double2 *in1, double2 *in2, double *out)
Sums the absolute value of two complex arrays.
tuple dx
Definition: en.py:59
__global__ void cMultDensity(double2 *in1, double2 *in2, double2 *out, double dt, double mass, int gstate, double gDenConst)
Kernel for complex multiplication with nonlinear density term.
__global__ void scalarMult(double2 *in, double factor, double2 *out)
Complex field scaling and renormalisation. Used mainly post-FFT.
__global__ void cMult(double2 *in1, double2 *in2, double2 *out)
Kernel for complex multiplication.
end % idx is for time
__global__ void scalarDiv(double2 *in, double factor, double2 *out)
Complex field scaling and renormalisation. Used mainly post-FFT.
tuple omega
Definition: observables.py:73
__global__ void reduce(double2 *in, double *out)
Not implemented.
__global__ void cMultDensity_ast(EqnNode_gpu *eqn, double2 *in, double2 *out, double dx, double dy, double dz, double time, int e_num, double dt, double mass, int gstate, double gDenConst)
Kernel for complex multiplication with nonlinear density term.
__global__ void scalarDiv_wfcNorm(double2 *in, double dr, double *pSum, double2 *out)
Used as part of multipass to renormalise the wavefucntion.
__device__ unsigned int getTid3d3d()
Indexing of threads in a block on device.
__global__ void vecConjugate(double2 *in, double2 *out)
Conjugate of double2*.
__device__ double complexMagnitudeSquared(double2 in)
Return the squared magnitude of a complex number. $|(a+{i}b)*(a-{i}b)|$.
__global__ void ast_cmult(double2 *array, double2 *array_out, EqnNode_gpu *eqn, double dx, double dy, double dz, double time, int element_num)
Complex multiplication of array with AST.
__global__ void ast_op_mult(double2 *array, double2 *array_out, EqnNode_gpu *eqn, double dx, double dy, double dz, double time, int element_num, int evolution_type, double dt)
Multiplication of array with AST Operator.
__global__ void thread_test(double *input, double *output)
Performs wavefunction renormalisation using parallel summation and applying scalarDiv_wfcNorm.
__device__ unsigned int getBid3d3d()
Indexing of blocks on device.
__global__ void energyCalc(double2 *wfc, double2 *op, double dt, double2 *energy, int gnd_state, int op_space, double sqrt_omegaz_mass, double gDenConst)
Calculates energy of the current state during evolution. Not implemented.
__device__ double2 subtract(double2 a, double2 b)
subtraction operation for 2 double2 values
Struct to hold the node information for the AST on the GPU.
Definition: ds.h:71
__global__ void scalarPow(double2 *in, double param, double2 *out)
Complex field raised to a power.
__global__ void make_cufftDoubleComplex(double *in, double2 *out)
transforms an array of doubles into double2&#39;s
tic for a
Definition: GPE_2d.m:90
__global__ void scalarDiv2D(double2 *, double2 *)
Complex field scaling and renormalisation. Not implemented. Use scalarDiv.
__global__ void l2_norm(double *in1, double *in2, double *in3, double *out)
performs the l2 normalization of the provided terms
int val
Definition: vort.py:104
__device__ double2 im_ast(double val, double dt)
Function to find AST operator in imaginary-time.
__global__ void angularOp(double omega, double dt, double2 *wfc, double *xpyypx, double2 *out)
Calculates angular momentum. Not fully implemented. Handled in post-processing instead.
__global__ void is_eq(bool *a, bool *b, bool *ans)
checks to arrays to see if they are equal
__device__ unsigned int getGid3d3d()
Indexing of threads on grid.
tuple dt
Definition: en.py:61
__global__ void multipass(double2 *input, double2 *output, int pass)
Performs wavefunction renormalisation using parallel summation and applying scalarDiv_wfcNorm.
__device__ double2 real_ast(double val, double dt)
Function to find AST operator in real-time.
__device__ double2 mult(double2 a, double b)
multiplication operation for double2 and double values
__global__ void pinVortex(double2 *in1, double2 *in2, double2 *out)
Hold vortex at specified position. Not implemented. cMultPhi should implement required functionality...
__global__ void pSum(double *in1, double *output, int pass)
Performs parallel sum. Not verified. I use multipass instead.
__global__ void zeros(bool *in, bool *out)
Sets boolean array to 0.
__global__ void ast_mult(double *array, double *array_out, EqnNode_gpu *eqn, double dx, double dy, double dz, double time, int element_num)
Multiplication of array with AST.
__global__ void complexMultiply(double2 *in1, double2 *in2, double2 *out)
Complex multiplication of two input arrays.
__global__ void vecMult(double2 *in, double *factor, double2 *out)
Complex field scaling and renormalisation. Used mainly post-FFT.
tuple mass
Definition: observables.py:72
__global__ void scalarDiv1D(double2 *, double2 *)
Complex field scaling and renormalisation. Not implemented. Use scalarDiv.