GPUE  v1.0
GPU Gross-Pitaevskii Equation numerical solver for Bose-Einstein condensates
operators.h
Go to the documentation of this file.
1 //##############################################################################
13  //#############################################################################
14 
15 #ifndef OPERATORS_H
16 #define OPERATORS_H
17 
18 #include "../include/ds.h"
19 #include "../include/constants.h"
20 #include <sys/stat.h>
21 #include <unordered_map>
22 //#include <boost/math/special_functions.hpp>
23 
24 // function to find Bz, the curl of the gauge field
33 double *curl2d(Grid &par, double *Ax, double *Ay);
34 
44 double *curl3d_x(Grid &par, double *Ax, double *Ay, double *Az);
45 
55 double *curl3d_y(Grid &par, double *Ax, double *Ay, double *Az);
56 
66 double *curl3d_z(Grid &par, double *Ax, double *Ay, double *Az);
67 
77 double *curl3d_r(Grid &par, double *Bx, double *By);
78 
88 double *curl3d_phi(Grid &par, double *Bx, double *By);
89 
90 
96 std::string filecheck(std::string filename);
97 
105 void file_A(std::string filename, double *A, double omega);
106 
107 /*----------------------------------------------------------------------------//
108 * GPU KERNELS
109 *-----------------------------------------------------------------------------*/
110 
114 void generate_p_space(Grid &par);
115 
119 void generate_K(Grid &par);
120 
124 __global__ void simple_K(double *xp, double *yp, double *zp, double mass,
125  double *K);
126 
130 void generate_gauge(Grid &par);
131 
135 __global__ void kconstant_A(double *x, double *y, double *z,
136  double xMax, double yMax, double zMax,
137  double omegaX, double omegaY, double omegaZ,
138  double omega, double fudge, double *A);
139 
143 __global__ void krotation_Ax(double *x, double *y, double *z,
144  double xMax, double yMax, double zMax,
145  double omegaX, double omegaY, double omegaZ,
146  double omega, double fudge, double *A);
147 
151 __global__ void krotation_Ay(double *x, double *y, double *z,
152  double xMax, double yMax, double zMax,
153  double omegaX, double omegaY, double omegaZ,
154  double omega, double fudge, double *A);
155 
159 __global__ void kring_rotation_Ax(double *x, double *y, double *z,
160  double xMax, double yMax, double zMax,
161  double omegaX, double omegaY, double omegaZ,
162  double omega, double fudge, double *A);
163 
167 __global__ void kring_rotation_Ay(double *x, double *y, double *z,
168  double xMax, double yMax, double zMax,
169  double omegaX, double omegaY, double omegaZ,
170  double omega, double fudge, double *A);
171 
175 __global__ void kring_rotation_Az(double *x, double *y, double *z,
176  double xMax, double yMax, double zMax,
177  double omegaX, double omegaY, double omegaZ,
178  double omega, double fudge, double *A);
179 
180 
184 __global__ void ktest_Ay(double *x, double *y, double *z,
185  double xMax, double yMax, double zMax,
186  double omegaX, double omegaY, double omegaZ,
187  double omega, double fudge, double *A);
188 
192 __global__ void ktest_Ax(double *x, double *y, double *z,
193  double xMax, double yMax, double zMax,
194  double omegaX, double omegaY, double omegaZ,
195  double omega, double fudge, double *A);
196 
200 __global__ void kring_Az(double *x, double *y, double *z,
201  double xMax, double yMax, double zMax,
202  double omegaX, double omegaY, double omegaZ,
203  double omega, double fudge, double *A);
204 
205 
209 void generate_fields(Grid &par);
210 
214 __global__ void kharmonic_V(double *x, double *y, double *z, double *items,
215  double *Ax, double *Ay, double *Az, double *V);
216 
220 __global__ void ktorus_V(double *x, double *y, double *z, double *items,
221  double *Ax, double *Ay, double *Az, double *V);
222 
226 __global__ void kstd_wfc(double *x, double *y, double *z, double *items,
227  double winding, double *phi, double2 *wfc);
228 
232 __global__ void ktorus_wfc(double *x, double *y, double *z, double *items,
233  double winding, double *phi, double2 *wfc);
234 
239 __global__ void aux_fields(double *V, double *K, double gdt, double dt,
240  double* Ax, double *Ay, double* Az,
241  double *px, double *py, double *pz,
242  double* pAx, double* pAy, double* pAz,
243  double2* GV, double2* EV, double2* GK, double2* EK,
244  double2* GpAx, double2* GpAy, double2* GpAz,
245  double2* EpAx, double2* EpAy, double2* EpAz);
246 // Function to generate grid and treads
247 void generate_grid(Grid &par);
248 #endif
double * curl2d(Grid &par, double *Ax, double *Ay)
Finds Bz, the curl of the gauge field.
% % % starting wavefunction wfc
Definition: GPE_2d.m:52
double * curl3d_r(Grid &par, double *Bx, double *By)
Finds Br, the curl of the gauge field.
void generate_K(Grid &par)
This function calls the appropriate K kernel.
double * curl3d_phi(Grid &par, double *Bx, double *By)
Finds Bphi, the curl of the gauge field.
__global__ void kring_rotation_Ay(double *x, double *y, double *z, double xMax, double yMax, double zMax, double omegaX, double omegaY, double omegaZ, double omega, double fudge, double *A)
Kernel for simple triangular lattice of rings, Ay.
double * curl3d_y(Grid &par, double *Ax, double *Ay, double *Az)
Finds By, the curl of the gauge field.
__global__ void simple_K(double *xp, double *yp, double *zp, double mass, double *K)
Simple kernel for generating K.
__global__ void kring_Az(double *x, double *y, double *z, double xMax, double yMax, double zMax, double omegaX, double omegaY, double omegaZ, double omega, double fudge, double *A)
Kernel for testing Az.
V
Definition: en.py:68
void generate_grid(Grid &par)
tuple yMax
Definition: en.py:63
tuple omega
Definition: observables.py:73
double * curl3d_z(Grid &par, double *Ax, double *Ay, double *Az)
Finds Bz, the curl of the gauge field.
__global__ void ktorus_V(double *x, double *y, double *z, double *items, double *Ax, double *Ay, double *Az, double *V)
Kernel to generate toroidal V (3d)
A
Definition: GPE_2d.m:67
void generate_p_space(Grid &par)
Function to generate momentum grids.
__global__ void kring_rotation_Ax(double *x, double *y, double *z, double xMax, double yMax, double zMax, double omegaX, double omegaY, double omegaZ, double omega, double fudge, double *A)
Kernel for simple triangular lattice of rings, Ax.
__global__ void kring_rotation_Az(double *x, double *y, double *z, double xMax, double yMax, double zMax, double omegaX, double omegaY, double omegaZ, double omega, double fudge, double *A)
Kernel for simple triangular lattice of rings, Az.
__global__ void ktest_Ax(double *x, double *y, double *z, double xMax, double yMax, double zMax, double omegaX, double omegaY, double omegaZ, double omega, double fudge, double *A)
Kernel for testing Ax.
__global__ void ktorus_wfc(double *x, double *y, double *z, double *items, double winding, double *phi, double2 *wfc)
Kernel to generate toroidal wavefunction.
void generate_fields(Grid &par)
Function to generate V.
__global__ void kharmonic_V(double *x, double *y, double *z, double *items, double *Ax, double *Ay, double *Az, double *V)
Kernel to generate harmonic V.
__global__ void ktest_Ay(double *x, double *y, double *z, double xMax, double yMax, double zMax, double omegaX, double omegaY, double omegaZ, double omega, double fudge, double *A)
Kernel for testing Ay.
tuple omegaZ
Definition: observables.py:71
tuple xMax
Definition: en.py:62
void file_A(std::string filename, double *A, double omega)
Reads A from file.
__global__ void kstd_wfc(double *x, double *y, double *z, double *items, double winding, double *phi, double2 *wfc)
Kernel to generate simple gaussian wavefunction.
__global__ void aux_fields(double *V, double *K, double gdt, double dt, double *Ax, double *Ay, double *Az, double *px, double *py, double *pz, double *pAx, double *pAy, double *pAz, double2 *GV, double2 *EV, double2 *GK, double2 *EK, double2 *GpAx, double2 *GpAy, double2 *GpAz, double2 *EpAx, double2 *EpAy, double2 *EpAz)
Kernel to generate all auxiliary fields.
void generate_gauge(Grid &par)
Function to generate game fields.
__global__ void kconstant_A(double *x, double *y, double *z, double xMax, double yMax, double zMax, double omegaX, double omegaY, double omegaZ, double omega, double fudge, double *A)
constant Kernel A
tuple omegaX
Definition: observables.py:74
double * curl3d_x(Grid &par, double *Ax, double *Ay, double *Az)
Finds Bx, the curl of the gauge field.
tuple dt
Definition: en.py:61
__global__ void krotation_Ay(double *x, double *y, double *z, double xMax, double yMax, double zMax, double omegaX, double omegaY, double omegaZ, double omega, double fudge, double *A)
Kernel for simple rotational case, Ay.
K
Definition: GPE_2d.m:59
Class to hold the variable map and grid information.
Definition: ds.h:86
def par
Definition: plot.py:237
__global__ void krotation_Ax(double *x, double *y, double *z, double xMax, double yMax, double zMax, double omegaX, double omegaY, double omegaZ, double omega, double fudge, double *A)
Kernel for simple rotational case, Ax.
std::string filecheck(std::string filename)
Determines if file exists, requests new file if it does not.
tuple mass
Definition: observables.py:72