GPUE  v1.0
GPU Gross-Pitaevskii Equation numerical solver for Bose-Einstein condensates
vis Namespace Reference

Functions

def delaunay (dataName, dataType, value)
 
def voronoi (dataName, dataType, value)
 
def laplacian (density, name, imgdpi)
 
def struct_fact (density, name, imgdpi)
 
def opPot (dataName, imgdpi)
 
def hist_gen (name, value, num_bins)
 
def image_gen (dataName, initValue, finalValue, increment, imgdpi)
 
def image_gen_single (dataName, value, imgdpi, opmode)
 
def vort_traj (name, imgdpi)
 
def scaleAxis (data, dataName, label, value, imgdpi)
 
def overlap (dataName, initValue, finalValue, increment)
 

Variables

 CPUs = os.environ['SLURM_JOB_CPUS_PER_NODE']
 
 prec
 
 c = ConfigParser.ConfigParser()
 
 xDim = int(c.getfloat('Params','xDim'))
 
 yDim = int(c.getfloat('Params','yDim'))
 
 gndMaxVal = int(c.getfloat('Params','gsteps'))
 
 evMaxVal = int(c.getfloat('Params','esteps'))
 
 incr = int(c.getfloat('Params','print_out'))
 
tuple sep = (c.getfloat('Params','dx'))
 
tuple dx = (c.getfloat('Params','dx'))
 
tuple dt = (c.getfloat('Params','dt'))
 
tuple xMax = (c.getfloat('Params','xMax'))
 
tuple yMax = (c.getfloat('Params','yMax'))
 
int num_vort = 0
 
 data = numpy.ndarray(shape=(xDim,yDim))
 
list gndImgList = []
 
list evImgList = []
 
list gnd_proc = []
 
list ev_proc = []
 
list i = gndImgList.pop()
 
list proc = gnd_proc + ev_proc
 
list p = proc.pop()
 

Detailed Description

vis.py - GPUE: Split Operator based GPU solver for Nonlinear 
Schrodinger Equation, Copyright (C) 2011-2015, Lee J. O'Riordan 
<loriordan@gmail.com>, Tadhg Morgan, Neil Crowley. All rights reserved.

Redistribution and use in source and binary forms, with or without 
modification, are permitted provided that the following conditions are 
met:

1. Redistributions of source code must retain the above copyright 
notice, this list of conditions and the following disclaimer.

2. Redistributions in binary form must reproduce the above copyright 
notice, this list of conditions and the following disclaimer in the 
documentation and/or other materials provided with the distribution.

3. Neither the name of the copyright holder nor the names of its 
contributors may be used to endorse or promote products derived from 
this software without specific prior written permission.

THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS 
"AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT 
LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A 
PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT 
HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, 
SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED
TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR 
PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF 
LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING 
NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS 
SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.

Function Documentation

◆ delaunay()

def vis.delaunay (   dataName,
  dataType,
  value 
)

Definition at line 76 of file vis.py.

76 def delaunay(dataName,dataType,value):
77  v_arr=genfromtxt(dataName + str(value) + dataType,delimiter=',' )
78  data = np.array([[row[0],row[1]] for row in v_arr])
79  dln = sp.spatial.Delaunay(data)
80  plt.triplot(data[:,0],data[:,1],dln.simplices.copy(),linewidth=0.5,color='b',marker='.')
81  plt.xlim(300,700);plt.ylim(300,700);
82  plt.savefig('delaun_' + str(value) + '.png',dpi=200)
83  print 'Saved Delaunay @ t=' + str(value)
84 
def delaunay(dataName, dataType, value)
Definition: vis.py:76

◆ hist_gen()

def vis.hist_gen (   name,
  value,
  num_bins 
)

Definition at line 131 of file vis.py.

131 def hist_gen(name,value,num_bins):
132  v_arr=genfromtxt('vort_arr_' + str(value),delimiter=',' )
133  H=[]
134  count=0
135 
136  for i1 in range(0,v_arr.size/2):
137  for i2 in range(i1,v_arr.size/2):
138  H.append(m.sqrt( abs(v_arr[i1][0]*sep - v_arr[i2][0]*sep)**2 + abs(v_arr[i1][1]*sep - v_arr[i2][1]*sep)**2 ))
139  count = count + 1
140  plt.title('Vortex lattice @ t=' + str(value*dt))
141  plt.ticklabel_format(style='scientific')
142  plt.ticklabel_format(style='scientific',axis='x', scilimits=(0,0))
143  h = plt.hist(H, bins=num_bins)
144  plt.savefig(name + "_" + str(value) + ".pdf")
145  plt.close()
146 
def hist_gen(name, value, num_bins)
Definition: vis.py:131

◆ image_gen()

def vis.image_gen (   dataName,
  initValue,
  finalValue,
  increment,
  imgdpi 
)

Definition at line 147 of file vis.py.

147 def image_gen(dataName, initValue, finalValue, increment,imgdpi):
148  for i in range(initValue,finalValue,increment):
149  if not os.path.exists(dataName+"r_"+str(i)+"_abspsi2.png"):
150  real=open(dataName + '_' + str(i)).read().splitlines()
151  img=open(dataName + 'i_' + str(i)).read().splitlines()
152  a_r = numpy.asanyarray(real,dtype='f8') #64-bit double
153  a_i = numpy.asanyarray(img,dtype='f8') #64-bit double
154  a = a_r[:] + 1j*a_i[:]
155  b = np.reshape(a,(xDim,yDim))
156  f = plt.imshow(abs(b)**2)
157  plt.jet()
158  plt.gca().invert_yaxis()
159  plt.savefig(dataName+"r_"+str(i)+"_abspsi2.png",dpi=imgdpi)
160  plt.close()
161  g = plt.imshow(np.angle(b))
162  plt.gca().invert_yaxis()
163  plt.savefig(dataName+"r_"+str(i)+"_phi.png",dpi=imgdpi)
164  plt.close()
165  f = plt.imshow(abs(np.fft.fftshift(np.fft.fft2(b)))**2)
166  plt.gca().invert_yaxis()
167  plt.jet()
168  plt.savefig(dataName+"p_"+str(i)+"_abspsi2.png",dpi=imgdpi)
169  plt.close()
170  g = plt.imshow(np.angle(np.fft.fftshift(np.fft.fft2(b))))
171  plt.gca().invert_yaxis()
172  plt.savefig(dataName+"p_"+str(i)+"_phi.png",dpi=imgdpi)
173  plt.close()
174  print "Saved figure: " + str(i) + ".png"
175  plt.close()
176  else:
177  print "File(s) " + str(i) +".png already exist."
178 
def image_gen(dataName, initValue, finalValue, increment, imgdpi)
Definition: vis.py:147

◆ image_gen_single()

def vis.image_gen_single (   dataName,
  value,
  imgdpi,
  opmode 
)

Definition at line 179 of file vis.py.

References laplacian(), and struct_fact().

179 def image_gen_single(dataName, value, imgdpi,opmode):
180  real=open(dataName + '_' + str(0)).read().splitlines()
181  img=open(dataName + 'i_' + str(0)).read().splitlines()
182  a1_r = numpy.asanyarray(real,dtype='f8') #128-bit complex
183  a1_i = numpy.asanyarray(img,dtype='f8') #128-bit complex
184  a1 = a1_r[:] + 1j*a1_i[:]
185  b1 = np.reshape(a1,(xDim,yDim))
186 
187  if not os.path.exists(dataName+"r_"+str(value)+"_abspsi2.png"):
188  real=open(dataName + '_' + str(value)).read().splitlines()
189  img=open(dataName + 'i_' + str(value)).read().splitlines()
190  a_r = numpy.asanyarray(real,dtype='f8') #128-bit complex
191  a_i = numpy.asanyarray(img,dtype='f8') #128-bit complex
192  a = a_r[:] + 1j*a_i[:]
193  b = np.reshape(a,(xDim,yDim))
194  m_val=np.max(np.abs(b)**2)
195  #scaleAxis(b,dataName,"_abspsi2",value,imgdpi)
196  if opmode & 0b100000 > 0:
197 # fig, ax = plt.subplots()
198 # #plt.rc('text',usetex=True)
199 # #plt.rc('font',family='serif')
200 # f = plt.imshow((abs(b)**2 - abs(b1)**2),cmap='gnuplot2',vmin=-6,vmax=6)
201 # plt.title(r'$\left(\rho( r,t ) - \rho( r,t_0 )\right),t=$' + str(value*dt))
202 # cbar = fig.colorbar(f)
203 # plt.gca().set_xlabel('x '+ str((dx)))
204 # plt.gca().set_ylabel('x '+ str(dx))
205 # plt.gca().invert_yaxis()
206 # plt.savefig(dataName+"r_"+str(value)+"_diffabspsi2.png",dpi=imgdpi)
207 # plt.close()
208 # #plt.rc('text',usetex=True)
209 # #plt.rc('font',family='serif')
210 
211  fig, ax = plt.subplots()
212  f = plt.imshow((abs(b)**2),cmap='gnuplot2',vmin=0,vmax=1e7)
213  plt.title('rho(r) @ t=' + str(value*dt))
214  # plt.title(r'$\\rho \left( r,t \right),\,t=$' + str(value*dt))
215 
216  #plugins.connect(fig, plugins.MousePosition(fontsize=14))
217 
218  cbar = fig.colorbar(f)
219  plt.gca().set_xlabel('x '+ str((dx)))
220  plt.gca().set_ylabel('x '+ str(dx))
221  plt.gca().invert_yaxis()
222  plt.savefig(dataName+"r_"+str(value)+"_abspsi2.png",dpi=imgdpi)
223  plt.axis('off')
224  plt.savefig(dataName+"r_"+str(value)+"_abspsi2_axis0.pdf",bbox_inches='tight',dpi=imgdpi)
225  plt.close()
226 
227  if opmode & 0b010000 > 0:
228  fig, ax = plt.subplots()
229  g = plt.imshow(np.angle(b))
230  cbar = fig.colorbar(g)
231  plt.gca().invert_yaxis()
232  plt.title('theta(r) @ t=' + str(value*dt))
233  plt.savefig(dataName+"r_"+str(value)+"_phi.png",dpi=imgdpi)
234  plt.close()
235 
236  if opmode & 0b001000 > 0:
237  fig, ax = plt.subplots()
238  f = plt.imshow(abs(np.fft.fftshift(np.fft.fft2(b)))**2)
239  cbar = fig.colorbar(f)
240  plt.gca().invert_yaxis()
241  plt.jet()
242  plt.title('rho(p) @ t=' + str(value*dt))
243  plt.savefig(dataName+"p_"+str(value)+"_abspsi2.png",dpi=imgdpi)
244  plt.close()
245 
246  if opmode & 0b000100 > 0:
247  fig, ax = plt.subplots()
248  g = plt.imshow(np.angle(np.fft.fftshift(np.fft.fft2(b))))
249  cbar = fig.colorbar(g)
250  plt.gca().invert_yaxis()
251  plt.title('theta(p) @ t=' + str(value*dt))
252  plt.savefig(dataName+"p_"+str(value)+"_phi.png",dpi=imgdpi)
253  plt.close()
254 
255  if opmode & 0b000010 > 0:
256  struct_fact(abs(b)**2,dataName+"_" + str(value),imgdpi)
257 
258  if opmode & 0b000001 > 0:
259  laplacian(abs(b)**2,dataName+"_" + str(value),imgdpi)
260 
261  print "Saved figure: " + str(value) + ".png"
262  plt.close()
263  else:
264  print "File(s) " + str(value) +".png already exist."
265 
def image_gen_single(dataName, value, imgdpi, opmode)
Definition: vis.py:179
def struct_fact(density, name, imgdpi)
Definition: vis.py:109
def laplacian(density, name, imgdpi)
Definition: vis.py:94
Here is the call graph for this function:

◆ laplacian()

def vis.laplacian (   density,
  name,
  imgdpi 
)

Definition at line 94 of file vis.py.

Referenced by image_gen_single().

94 def laplacian(density,name,imgdpi):
95  gx,gy = np.gradient(density)
96  g2x,gxgy = np.gradient(gx)
97  gygx,g2y = np.gradient(gy)
98  fig, ax = plt.subplots()
99  #f = plt.quiver(gx,gy)
100  f = plt.imshow((g2x**2 + g2y**2),cmap=plt.get_cmap('spectral'))
101  cbar = fig.colorbar(f)
102  plt.savefig(name + "_laplacian.png",dpi=imgdpi)
103  plt.close()
104  f = plt.imshow((gxgy - gygx),cmap=plt.get_cmap('spectral'))
105  cbar = fig.colorbar(f)
106  plt.savefig(name + "_dxdy.png",dpi=imgdpi)
107  plt.close()
108 
def laplacian(density, name, imgdpi)
Definition: vis.py:94
Here is the caller graph for this function:

◆ opPot()

def vis.opPot (   dataName,
  imgdpi 
)

Definition at line 119 of file vis.py.

Referenced by overlap().

119 def opPot(dataName,imgdpi):
120  data = open(dataName).read().splitlines()
121  a = numpy.asanyarray(data,dtype='f8')
122  b = np.reshape(a,(xDim,yDim))
123  fig, ax = plt.subplots()
124  f = plt.imshow((b))
125  plt.gca().invert_yaxis()
126  cbar = fig.colorbar(f)
127  plt.jet()
128  plt.savefig(dataName + ".png",dpi=imgdpi)
129  plt.close()
130 
def opPot(dataName, imgdpi)
Definition: vis.py:119
Here is the caller graph for this function:

◆ overlap()

def vis.overlap (   dataName,
  initValue,
  finalValue,
  increment 
)

Definition at line 307 of file vis.py.

References opPot().

307 def overlap(dataName, initValue, finalValue, increment):
308  real=open(dataName + '_' + str(0)).read().splitlines()
309  img=open(dataName + 'i_' + str(0)).read().splitlines()
310  a_r = numpy.asanyarray(real,dtype='f8') #128-bit complex
311  a_i = numpy.asanyarray(img,dtype='f8') #128-bit complex
312  wfc0 = a_r[:] + 1j*a_i[:]
313  for i in range(initValue,finalValue,increment):
314  real=open(dataName + '_' + str(value)).read().splitlines()
315  img=open(dataName + 'i_' + str(value)).read().splitlines()
316  a_r = numpy.asanyarray(real,dtype='f8') #128-bit complex
317  a_i = numpy.asanyarray(img,dtype='f8') #128-bit complex
318  a = a_r[:] + 1j*a_i[:]
319  b = np.dot(wfc0,a)
320  print i, np.sum(b)
321 
def overlap(dataName, initValue, finalValue, increment)
Definition: vis.py:307
Here is the call graph for this function:

◆ scaleAxis()

def vis.scaleAxis (   data,
  dataName,
  label,
  value,
  imgdpi 
)

Definition at line 296 of file vis.py.

296 def scaleAxis(data,dataName,label,value,imgdpi):
297  fig, ax = plt.subplots()
298  ax.xaxis.set_major_locator(ScaledLocator(dx=dx))
299  ax.xaxis.set_major_formatter(ScaledLocator(dx=dx))
300  f = plt.imshow(abs(data)**2)
301  cbar = fig.colorbar(f)
302  plt.gca().invert_yaxis()
303  plt.jet()
304  plt.savefig(dataName+"r_"+str(value)+"_"+label +".png",dpi=imgdpi)
305  plt.close()
306 
def scaleAxis(data, dataName, label, value, imgdpi)
Definition: vis.py:296

◆ struct_fact()

def vis.struct_fact (   density,
  name,
  imgdpi 
)

Definition at line 109 of file vis.py.

Referenced by image_gen_single().

109 def struct_fact(density,name,imgdpi):
110  fig, ax = plt.subplots()
111  #f = plt.quiver(gx,gy)
112  f = plt.imshow((np.abs(np.fft.fftshift(np.fft.fft2(density)))),cmap=plt.get_cmap('prism'))
113  cbar = fig.colorbar(f)
114  cbar.set_clim(1e6,1e11)
115  plt.jet()
116  plt.savefig(name + "_struct_log10.png",dpi=imgdpi)
117  plt.close()
118 
def struct_fact(density, name, imgdpi)
Definition: vis.py:109
Here is the caller graph for this function:

◆ voronoi()

def vis.voronoi (   dataName,
  dataType,
  value 
)

Definition at line 85 of file vis.py.

85 def voronoi(dataName,dataType,value):
86  v_arr=genfromtxt(dataName + str(value) + dataType,delimiter=',' )
87  data = [[row[0],row[1]] for row in v_arr]
88  vor = Voronoi(data)
89  voronoi_plot_2d(vor)
90  plt.xlim(300,700);plt.ylim(300,700);
91  plt.savefig('voronoi_' + str(value) + '.png',dpi=200)
92  print 'Saved Voronoi @ t=' + str(value)
93 
def voronoi(dataName, dataType, value)
Definition: vis.py:85

◆ vort_traj()

def vis.vort_traj (   name,
  imgdpi 
)

Definition at line 266 of file vis.py.

266 def vort_traj(name,imgdpi):
267  evMaxVal_l = evMaxVal
268  H=genfromtxt('vort_arr_0',delimiter=',' )
269  count=0
270  for i1 in range(incr,evMaxVal_l,incr):
271  try:
272  v_arr=genfromtxt('vort_lsq_' + str(i1) + '.csv',delimiter=',' )
273  H=np.column_stack((H,v_arr))
274  except:
275  evMaxVal_l = i1
276  break
277  X=np.zeros((evMaxVal_l/incr),dtype=np.float64)
278  Y=np.zeros((evMaxVal_l/incr),dtype=np.float64)
279  H=np.reshape(H,([num_vort,2,evMaxVal_l/incr]),order='F')
280  for i1 in range(0, num_vort):
281  for i2 in range(0,evMaxVal_l/incr):
282  X[i2]=(H[i1,0,i2]*dx) - xMax
283  Y[i2]=(H[i1,1,i2]*dx) - yMax
284  h = plt.plot(X,Y,color=(r.random(),r.random(),r.random(),0.85),linewidth=0.1)
285  plt.axis('equal')
286  plt.title('Vort(x,y) from t=0 to t='+str(evMaxVal_l*dt)+" s")
287 
288  plt.axis((-xMax/2.0, xMax/2.0, -yMax/2.0, yMax/2.0))
289  plt.ticklabel_format(style='scientific')
290  plt.ticklabel_format(style='scientific',axis='x', scilimits=(0,0))
291  plt.ticklabel_format(style='scientific',axis='y', scilimits=(0,0))
292  plt.savefig(name +".pdf")
293  plt.close()
294  print "Trajectories plotted."
295 
def vort_traj(name, imgdpi)
Definition: vis.py:266

Variable Documentation

◆ c

vis.c = ConfigParser.ConfigParser()

Definition at line 57 of file vis.py.

◆ CPUs

vis.CPUs = os.environ['SLURM_JOB_CPUS_PER_NODE']

Definition at line 34 of file vis.py.

◆ data

vis.data = numpy.ndarray(shape=(xDim,yDim))

Definition at line 74 of file vis.py.

◆ dt

tuple vis.dt = (c.getfloat('Params','dt'))

Definition at line 69 of file vis.py.

◆ dx

tuple vis.dx = (c.getfloat('Params','dx'))

Definition at line 68 of file vis.py.

◆ ev_proc

list vis.ev_proc = []

Definition at line 333 of file vis.py.

◆ evImgList

list vis.evImgList = []

Definition at line 327 of file vis.py.

◆ evMaxVal

vis.evMaxVal = int(c.getfloat('Params','esteps'))

Definition at line 65 of file vis.py.

◆ gnd_proc

list vis.gnd_proc = []

Definition at line 332 of file vis.py.

◆ gndImgList

list vis.gndImgList = []

Definition at line 326 of file vis.py.

◆ gndMaxVal

vis.gndMaxVal = int(c.getfloat('Params','gsteps'))

Definition at line 64 of file vis.py.

◆ i

list vis.i = gndImgList.pop()

Definition at line 335 of file vis.py.

◆ incr

vis.incr = int(c.getfloat('Params','print_out'))

Definition at line 66 of file vis.py.

◆ num_vort

int vis.num_vort = 0

Definition at line 72 of file vis.py.

◆ p

list vis.p = proc.pop()

Definition at line 347 of file vis.py.

◆ prec

vis.prec

Definition at line 56 of file vis.py.

◆ proc

list vis.proc = gnd_proc + ev_proc

Definition at line 341 of file vis.py.

◆ sep

tuple vis.sep = (c.getfloat('Params','dx'))

Definition at line 67 of file vis.py.

◆ xDim

vis.xDim = int(c.getfloat('Params','xDim'))

Definition at line 62 of file vis.py.

◆ xMax

tuple vis.xMax = (c.getfloat('Params','xMax'))

Definition at line 70 of file vis.py.

◆ yDim

vis.yDim = int(c.getfloat('Params','yDim'))

Definition at line 63 of file vis.py.

◆ yMax

tuple vis.yMax = (c.getfloat('Params','yMax'))

Definition at line 71 of file vis.py.