|  | GPUE
    v1.0
    GPU Gross-Pitaevskii Equation numerical solver for Bose-Einstein condensates | 
GPU kernel definitions. More...
#include <stdio.h>

Go to the source code of this file.
| Functions | |
| __device__ double2 | subtract (double2 a, double2 b) | 
| subtraction operation for 2 double2 values  More... | |
| __device__ double2 | add (double2 a, double2 b) | 
| addition operation for 2 double2 values  More... | |
| __device__ double2 | pow (double2 a, int b) | 
| power operation for a double2  More... | |
| __device__ double2 | mult (double2 a, double b) | 
| multiplication operation for double2 and double values  More... | |
| __device__ double2 | mult (double2 a, double2 b) | 
| multiplication operation for 2 double2 values  More... | |
| __global__ void | make_cufftDoubleComplex (double *in, double2 *out) | 
| transforms an array of doubles into double2's  More... | |
| __device__ unsigned int | getGid3d3d () | 
| Indexing of threads on grid.  More... | |
| __device__ unsigned int | getBid3d3d () | 
| Indexing of blocks on device.  More... | |
| __device__ unsigned int | getTid3d3d () | 
| Indexing of threads in a block on device.  More... | |
| __global__ void | is_eq (bool *a, bool *b, bool *ans) | 
| checks to arrays to see if they are equal  More... | |
| __device__ double | complexMagnitude (double2 in) | 
| Calculates magnitude of complex number. $|a + ib|$.  More... | |
| __global__ void | complexMultiply (double2 *in1, double2 *in2, double2 *out) | 
| Complex multiplication of two input arrays.  More... | |
| __host__ __device__ double2 | complexMultiply (double2 in1, double2 in2) | 
| Complex multiplication of two input values.  More... | |
| __device__ double2 | make_complex (double in, int evolution_type) | 
| Transforms field value into operator.  More... | |
| __global__ void | complexAbsSum (double2 *in1, double2 *in2, double *out) | 
| Sums the absolute value of two complex arrays.  More... | |
| __global__ void | complexMagnitude (double2 *in, double *out) | 
| Complex magnitude of a double2 array.  More... | |
| __device__ double | complexMagnitudeSquared (double2 in) | 
| Return the squared magnitude of a complex number. $|(a+{i}b)*(a-{i}b)|$.  More... | |
| __global__ void | complexMagnitudeSquared (double2 *in, double *out) | 
| Complex magnitude of a double2 array.  More... | |
| __global__ void | complexMagnitudeSquared (double2 *in, double2 *out) | 
| Complex magnitude of a double2 array.  More... | |
| __device__ double2 | conjugate (double2 in) | 
| Returns conjugate of the a complex number.  More... | |
| __device__ double2 | realCompMult (double scalar, double2 comp) | 
| Multiply real scalar by a complex number.  More... | |
| __global__ void | cMult (double2 *in1, double2 *in2, double2 *out) | 
| Kernel for complex multiplication.  More... | |
| __global__ void | cMultPhi (double2 *in1, double *in2, double2 *out) | 
| Kernel for multiplcation with real array and complex array.  More... | |
| __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.  More... | |
| __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.  More... | |
| __global__ void | pinVortex (double2 *in1, double2 *in2, double2 *out) | 
| Hold vortex at specified position. Not implemented. cMultPhi should implement required functionality.  More... | |
| __global__ void | vecMult (double2 *in, double *factor, double2 *out) | 
| Complex field scaling and renormalisation. Used mainly post-FFT.  More... | |
| __global__ void | l2_norm (double *in1, double *in2, double *in3, double *out) | 
| performs the l2 normalization of the provided terms  More... | |
| __global__ void | l2_norm (double2 *in1, double2 *in2, double2 *in3, double *out) | 
| performs the l2 normalization of the provided terms  More... | |
| __global__ void | l2_norm (double *in1, double *in2, double *out) | 
| performs the l2 normalization of the provided terms  More... | |
| __global__ void | l2_norm (double2 *in1, double2 *in2, double *out) | 
| performs the l2 normalization of the provided terms  More... | |
| __global__ void | scalarDiv (double2 *in, double factor, double2 *out) | 
| Complex field scaling and renormalisation. Used mainly post-FFT.  More... | |
| __global__ void | scalarDiv (double *in, double factor, double *out) | 
| Real field scaling and renormalisation. Used mainly post-FFT.  More... | |
| __global__ void | scalarMult (double2 *in, double factor, double2 *out) | 
| Complex field scaling and renormalisation. Used mainly post-FFT.  More... | |
| __global__ void | scalarPow (double2 *in, double param, double2 *out) | 
| Complex field raised to a power.  More... | |
| __global__ void | vecConjugate (double2 *in, double2 *out) | 
| Conjugate of double2*.  More... | |
| __global__ void | scalarDiv1D (double2 *, double2 *) | 
| Complex field scaling and renormalisation. Not implemented. Use scalarDiv.  More... | |
| __global__ void | scalarDiv2D (double2 *, double2 *) | 
| Complex field scaling and renormalisation. Not implemented. Use scalarDiv.  More... | |
| __global__ void | scalarDiv_wfcNorm (double2 *in, double dr, double *pSum, double2 *out) | 
| Used as part of multipass to renormalise the wavefucntion.  More... | |
| __global__ void | reduce (double2 *in, double *out) | 
| Not implemented.  More... | |
| __global__ void | thread_test (double *input, double *output) | 
| Performs wavefunction renormalisation using parallel summation and applying scalarDiv_wfcNorm.  More... | |
| __global__ void | multipass (double2 *input, double2 *output, int pass) | 
| Performs wavefunction renormalisation using parallel summation and applying scalarDiv_wfcNorm.  More... | |
| __global__ void | multipass (double *input, double *output) | 
| Performs parallel summation of double arrays.  More... | |
| __global__ void | angularOp (double omega, double dt, double2 *wfc, double *xpyypx, double2 *out) | 
| Calculates angular momentum. Not fully implemented. Handled in post-processing instead.  More... | |
| __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.  More... | |
| __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.  More... | |
| __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.  More... | |
| __device__ double2 | real_ast (double val, double dt) | 
| Function to find AST operator in real-time.  More... | |
| __device__ double2 | im_ast (double val, double dt) | 
| Function to find AST operator in imaginary-time.  More... | |
| __global__ void | zeros (bool *in, bool *out) | 
| Sets boolean array to 0.  More... | |
| __global__ void | set_eq (double *in1, double *in2) | 
| Sets in2 to be equal to in1.  More... | |
| __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.  More... | |
| __device__ double2 | braKetMult (double2 in1, double2 in2) | 
| Performs bra-ket state multiplication. Not fully implemented.  More... | |
| __global__ void | pSum (double *in1, double *output, int pass) | 
| Performs parallel sum. Not verified. I use multipass instead.  More... | |
GPU kernel definitions.
Kernel definitions for all CUDA-enabled routines for solving GPE.
Definition in file kernels.h.
| __device__ double2 add | ( | double2 | a, | 
| double2 | b | ||
| ) | 
addition operation for 2 double2 values
| __global__ void angularOp | ( | double | omega, | 
| double | dt, | ||
| double2 * | wfc, | ||
| double * | xpyypx, | ||
| double2 * | out | ||
| ) | 
Calculates angular momentum. Not fully implemented. Handled in post-processing instead.
| omega | Harmonic trap rotation frequency | 
| dt | Time-step for evolution | 
| wfc | Wavefunction | 
| xpyypx | L_z operator | 
| out | Output of calculation | 
| __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_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 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.
| 
 | inline | 
Performs bra-ket state multiplication. Not fully implemented.
| in1 | Bra | 
| in2 | Ket | 
| __global__ void cMult | ( | double2 * | in1, | 
| double2 * | in2, | ||
| double2 * | out | ||
| ) | 
Kernel for complex multiplication.
Multiplication for linear, non-linear and phase-imprinting of the condensate.
| in1 | Wavefunction input | 
| in2 | Evolution operator input | 
| out | Pass by reference output for multiplcation result | 
| __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.
| in1 | Wavefunction input | 
| in2 | Evolution operator input | 
| out | Pass by reference output for multiplication result | 
| dt | Timestep for evolution | 
| mass | Atomic species mass | 
| gState | If performing real (1) or imaginary (0) time evolution | 
| gDenConst | a constant for evolution | 
| __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.
| GPU | AST | 
| in | Wavefunction input | 
| out | Wavefunction output | 
| dx | |
| dy | |
| dz | |
| time | |
| element | number in AST | 
| dt | Timestep for evolution | 
| mass | Atomic species mass | 
| gState | If performing real (1) or imaginary (0) time evolution | 
| gDenConst | a constant for evolution | 
| __global__ void cMultPhi | ( | double2 * | in1, | 
| double * | in2, | ||
| double2 * | out | ||
| ) | 
Kernel for multiplcation with real array and complex array.
| in1 | Wavefunction input | 
| in2 | Evolution operator input | 
| out | Pass by reference output for multiplcation result | 
| __global__ void complexAbsSum | ( | double2 * | in1, | 
| double2 * | in2, | ||
| double * | out | ||
| ) | 
Sums the absolute value of two complex arrays.
| Array | 1 | 
| Array | 2 | 
| Output | 
| __device__ double complexMagnitude | ( | double2 | in | ) | 
Calculates magnitude of complex number. $|a + ib|$.
Helper functions for complex numbers
| in | Complex number | 
| __global__ void complexMagnitude | ( | double2 * | in, | 
| double * | out | ||
| ) | 
Complex magnitude of a double2 array.
| __device__ double complexMagnitudeSquared | ( | double2 | in | ) | 
Return the squared magnitude of a complex number. $|(a+{i}b)*(a-{i}b)|$.
| in | Complex number | 
| __global__ void complexMagnitudeSquared | ( | double2 * | in, | 
| double * | out | ||
| ) | 
Complex magnitude of a double2 array.
| __global__ void complexMagnitudeSquared | ( | double2 * | in, | 
| double2 * | out | ||
| ) | 
Complex magnitude of a double2 array.
| __global__ void complexMultiply | ( | double2 * | in1, | 
| double2 * | in2, | ||
| double2 * | out | ||
| ) | 
Complex multiplication of two input arrays.
| Input | 1 | 
| Input | 2 | 
| Output | 
| __host__ __device__ double2 complexMultiply | ( | double2 | in1, | 
| double2 | in2 | ||
| ) | 
Complex multiplication of two input values.
| Input | 1 | 
| Input | 2 | 
| __device__ double2 conjugate | ( | double2 | in | ) | 
Returns conjugate of the a complex number.
| in | Number to be conjugated | 
| __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.
| wfc | Wavefunction | 
| op | Operator to calculate energy for. | 
| dt | Time-step for evolution | 
| energy | Energy result output | 
| gnd_state | Wavefunction | 
| op_space | Check if position space with non-linear term or not. | 
| sqrt_omegaz_mass | sqrt(omegaZ/mass), part of the nonlin interaction term. | 
| __device__ unsigned int getBid3d3d | ( | ) | 
Indexing of blocks on device.
| __device__ unsigned int getGid3d3d | ( | ) | 
Indexing of threads on grid.
| __device__ unsigned int getTid3d3d | ( | ) | 
Indexing of threads in a block on device.
| __device__ double2 im_ast | ( | double | val, | 
| double | dt | ||
| ) | 
Function to find AST operator in imaginary-time.
| __global__ void is_eq | ( | bool * | a, | 
| bool * | b, | ||
| bool * | ans | ||
| ) | 
checks to arrays to see if they are equal
| __global__ void l2_norm | ( | double * | in1, | 
| double * | in2, | ||
| double * | in3, | ||
| double * | out | ||
| ) | 
performs the l2 normalization of the provided terms
| __global__ void l2_norm | ( | double2 * | in1, | 
| double2 * | in2, | ||
| double2 * | in3, | ||
| double * | out | ||
| ) | 
performs the l2 normalization of the provided terms
| __global__ void l2_norm | ( | double * | in1, | 
| double * | in2, | ||
| double * | out | ||
| ) | 
performs the l2 normalization of the provided terms
| __global__ void l2_norm | ( | double2 * | in1, | 
| double2 * | in2, | ||
| double * | out | ||
| ) | 
performs the l2 normalization of the provided terms
| __device__ double2 make_complex | ( | double | in, | 
| int | evolution_type | ||
| ) | 
Transforms field value into operator.
| Input | value | 
| Evolution | type (0 for imaginary, 1 for real) | 
| __global__ void make_cufftDoubleComplex | ( | double * | in, | 
| double2 * | out | ||
| ) | 
transforms an array of doubles into double2's
| __device__ double2 mult | ( | double2 | a, | 
| double | b | ||
| ) | 
multiplication operation for double2 and double values
| __device__ double2 mult | ( | double2 | a, | 
| double2 | b | ||
| ) | 
multiplication operation for 2 double2 values
| __global__ void multipass | ( | double2 * | input, | 
| double2 * | output, | ||
| int | pass | ||
| ) | 
Performs wavefunction renormalisation using parallel summation and applying scalarDiv_wfcNorm.
| input | Wavefunction to be renormalised | 
| output | Pass by reference return of renormalised wavefunction | 
| pass | Number of passes performed by routine | 
| __global__ void multipass | ( | double * | input, | 
| double * | output | ||
| ) | 
Performs parallel summation of double arrays.
| input | double array | 
| input | double array after summation | 
| __global__ void pinVortex | ( | double2 * | in1, | 
| double2 * | in2, | ||
| double2 * | out | ||
| ) | 
Hold vortex at specified position. Not implemented. cMultPhi should implement required functionality.
| in1 | Wavefunction input | 
| in2 | Evolution operator input | 
| out | Pass by reference output for multiplcation result | 
| __device__ double2 pow | ( | double2 | a, | 
| int | b | ||
| ) | 
power operation for a double2
| base | number | 
| power | 
Referenced by sepAvg().

| __global__ void pSum | ( | double * | in1, | 
| double * | output, | ||
| int | pass | ||
| ) | 
Performs parallel sum. Not verified. I use multipass instead.
| in1 | That which must be summed | 
| output | That which has been summed | 
| pass | Number of passes | 
| __device__ double2 real_ast | ( | double | val, | 
| double | dt | ||
| ) | 
Function to find AST operator in real-time.
| __device__ double2 realCompMult | ( | double | scalar, | 
| double2 | comp | ||
| ) | 
Multiply real scalar by a complex number.
| scalar | Scalar multiplier | 
| comp | Complex multiplicand | 
| __global__ void reduce | ( | double2 * | in, | 
| double * | out | ||
| ) | 
Not implemented.
| in | Input field | 
| out | Output values | 
| __global__ void scalarDiv | ( | double2 * | in, | 
| double | factor, | ||
| double2 * | out | ||
| ) | 
Complex field scaling and renormalisation. Used mainly post-FFT.
| in | Complex field to be scaled (divided, not multiplied) | 
| factor | Scaling factor to be used | 
| out | Pass by reference output for result | 
| __global__ void scalarDiv | ( | double * | in, | 
| double | factor, | ||
| double * | out | ||
| ) | 
Real field scaling and renormalisation. Used mainly post-FFT.
| in | Real field to be scaled (divided, not multiplied) | 
| factor | Scaling factor to be used | 
| out | Pass by reference output for result | 
| __global__ void scalarDiv1D | ( | double2 * | , | 
| double2 * | |||
| ) | 
Complex field scaling and renormalisation. Not implemented. Use scalarDiv.
| __global__ void scalarDiv2D | ( | double2 * | , | 
| double2 * | |||
| ) | 
Complex field scaling and renormalisation. Not implemented. Use scalarDiv.
| __global__ void scalarDiv_wfcNorm | ( | double2 * | in, | 
| double | dr, | ||
| double * | pSum, | ||
| double2 * | out | ||
| ) | 
Used as part of multipass to renormalise the wavefucntion.
| in | Complex field to be renormalised | 
| dr | Smallest area element of grid (dx*dy) | 
| pSum | GPU array used to store intermediate results during parallel summation | 
| __global__ void scalarMult | ( | double2 * | in, | 
| double | factor, | ||
| double2 * | out | ||
| ) | 
Complex field scaling and renormalisation. Used mainly post-FFT.
| in | Complex field to be scaled (multiplied, not divided) | 
| factor | Scaling factor to be used | 
| out | Pass by reference output for result | 
| __global__ void scalarPow | ( | double2 * | in, | 
| double | param, | ||
| double2 * | out | ||
| ) | 
Complex field raised to a power.
| in | Complex field to be scaled (multiplied, not divided) | 
| power | parameter | 
| out | Pass by reference output for result | 
| __global__ void set_eq | ( | double * | in1, | 
| double * | in2 | ||
| ) | 
Sets in2 to be equal to in1.
| __device__ double2 subtract | ( | double2 | a, | 
| double2 | b | ||
| ) | 
subtraction operation for 2 double2 values
| __global__ void thread_test | ( | double * | input, | 
| double * | output | ||
| ) | 
Performs wavefunction renormalisation using parallel summation and applying scalarDiv_wfcNorm.
| input | Wavefunction to be renormalised | 
| output | Pass by reference return of renormalised wavefunction | 
| pass | Number of passes performed by routine | 
| __global__ void vecConjugate | ( | double2 * | in, | 
| double2 * | out | ||
| ) | 
Conjugate of double2*.
| in | Complex field to be conjugated | 
| out | Pass by reference output for result | 
| __global__ void vecMult | ( | double2 * | in, | 
| double * | factor, | ||
| double2 * | out | ||
| ) | 
Complex field scaling and renormalisation. Used mainly post-FFT.
| in | Complex field to be scaled (divided, not multiplied) | 
| factor | Scaling vector to be used | 
| out | Pass by reference output for result | 
| __global__ void zeros | ( | bool * | in, | 
| bool * | out | ||
| ) | 
Sets boolean array to 0.