24 __device__ double2
subtract(double2
a, double2 b);
30 __device__ double2
add(double2
a, double2 b);
38 __device__ double2
pow(double2
a,
int b);
44 __device__ double2
mult(double2
a,
double b);
50 __device__ double2
mult(double2
a, double2 b);
80 __global__
void is_eq(
bool *
a,
bool *b,
bool *ans);
103 __global__
void complexMultiply(double2 *in1, double2 *in2, double2 *out);
121 __device__ double2
make_complex(
double in,
int evolution_type);
130 __global__
void complexAbsSum(double2 *in1, double2 *in2,
double *out);
164 __device__ double2
conjugate(double2 in);
173 __device__ double2
realCompMult(
double scalar, double2 comp);
188 __global__
void cMult(double2* in1, double2* in2, double2* out);
197 __global__
void cMultPhi(double2* in1,
double* in2, double2* out);
210 __global__
void cMultDensity(double2* in1, double2* in2, double2* out,
double dt,
double mass,
int gstate,
double gDenConst);
230 double dx,
double dy,
double dz,
double time,
231 int e_num,
double dt,
double mass,
int gstate,
242 __global__
void pinVortex(double2* in1, double2* in2, double2* out);
252 __global__
void vecMult(double2 *in,
double *factor, double2 *out);
258 __global__
void l2_norm(
double *in1,
double *in2,
double *in3,
double *out);
264 __global__
void l2_norm(double2 *in1, double2 *in2, double2 *in3,
double *out);
270 __global__
void l2_norm(
double *in1,
double *in2,
double *out);
276 __global__
void l2_norm(double2 *in1, double2 *in2,
double *out);
285 __global__
void scalarDiv(double2* in,
double factor, double2* out);
294 __global__
void scalarDiv(
double* in,
double factor,
double* out);
303 __global__
void scalarMult(double2* in,
double factor, double2* out);
312 __global__
void scalarPow(double2* in,
double param, double2* out);
320 __global__
void vecConjugate(double2 *in, double2 *out);
351 __global__
void reduce(double2* in,
double* out);
360 __global__
void thread_test(
double* input,
double* output);
369 __global__
void multipass(double2* input, double2* output,
int pass);
377 __global__
void multipass(
double* input,
double* output);
390 __global__
void angularOp(
double omega,
double dt, double2*
wfc,
double* xpyypx, double2* out);
397 double dx,
double dy,
double dz,
double time,
404 double dx,
double dy,
double dz,
double time,
410 __global__
void ast_op_mult(double2 *array, double2 *array_out,
412 double dx,
double dy,
double dz,
double time,
413 int element_num,
int evolution_type,
double dt);
431 __global__
void zeros(
bool *in,
bool *out);
437 __global__
void set_eq(
double *in1,
double *in2);
450 __global__
void energyCalc(double2 *
wfc, double2 *op,
double dt, double2 *energy,
int gnd_state,
int op_space,
double sqrt_omegaz_mass,
double gDenConst);
459 inline __device__ double2
braKetMult(double2 in1, double2 in2);
469 __global__
void pSum(
double* in1,
double* output,
int pass);
__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
__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.
__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.
__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.
__global__ void scalarDiv(double2 *in, double factor, double2 *out)
Complex field scaling and renormalisation. Used mainly post-FFT.
__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.
__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's
__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
__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.
__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.
__global__ void scalarDiv1D(double2 *, double2 *)
Complex field scaling and renormalisation. Not implemented. Use scalarDiv.