GPUE  v1.0
GPU Gross-Pitaevskii Equation numerical solver for Bose-Einstein condensates
observables.py
Go to the documentation of this file.
1 '''
2 observables.py - GPUE: Split Operator based GPU solver for Nonlinear
3 Schrodinger Equation, Copyright (C) 2011-2015, Lee J. O'Riordan
4 <loriordan@gmail.com>, Tadhg Morgan, Neil Crowley. All rights reserved.
5 
6 Redistribution and use in source and binary forms, with or without
7 modification, are permitted provided that the following conditions are
8 met:
9 
10 1. Redistributions of source code must retain the above copyright
11 notice, this list of conditions and the following disclaimer.
12 
13 2. Redistributions in binary form must reproduce the above copyright
14 notice, this list of conditions and the following disclaimer in the
15 documentation and/or other materials provided with the distribution.
16 
17 3. Neither the name of the copyright holder nor the names of its
18 contributors may be used to endorse or promote products derived from
19 this software without specific prior written permission.
20 
21 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
22 "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
23 LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A
24 PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT
25 HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
26 SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED
27 TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
28 PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
29 LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
30 NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
31 SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
32 '''
33 import os
34 from numpy import genfromtxt
35 import math as m
36 import matplotlib as mpl
37 import numpy as np
38 import scipy as sp
39 import numpy.matlib
40 mpl.use('Agg')
41 import multiprocessing as mp
42 from multiprocessing import Pool
43 from multiprocessing import Process
44 from matplotlib.ticker import ScalarFormatter
45 import matplotlib.pyplot as plt
46 import ConfigParser
47 import random as r
48 from decimal import *
49 from scipy.spatial import Delaunay
50 
51 HBAR = 1.05457148e-34
52 PI = 3.141592653589793
53 
54 getcontext().prec = 4
55 c = ConfigParser.ConfigParser()
56 c.readfp(open(r'Params.dat'))
57 
58 xDim = int(c.getfloat('Params','xDim'))
59 yDim = int(c.getfloat('Params','yDim'))
60 gndMaxVal = int(c.getfloat('Params','gsteps'))
61 evMaxVal = int(c.getfloat('Params','esteps'))
62 incr = int(c.getfloat('Params','print_out'))
63 #sep = (c.getfloat('Params','dx'))
64 dx = (c.getfloat('Params','dx'))
65 dy = (c.getfloat('Params','dy'))
66 dkx = (c.getfloat('Params','dpx'))
67 dky = (c.getfloat('Params','dpy'))
68 dt = (c.getfloat('Params','dt'))
69 xMax = (c.getfloat('Params','xMax'))
70 yMax = (c.getfloat('Params','yMax'))
71 omegaZ = (c.getfloat('Params','omegaZ'))
72 mass = (c.getfloat('Params','Mass'))
73 omega = (c.getfloat('Params','omega'))
74 omegaX = (c.getfloat('Params','omegaX'))
75 
76 try:
77  num_vort = int(c.getfloat('Params','Num_vort'))
78 except:
79  print '!num_vort undefined!'
80 N = int(c.getfloat('Params','atoms'))
81 
82 data = numpy.ndarray(shape=(xDim,yDim))
83 
84 x=np.asarray(open('x_0').read().splitlines(),dtype='f8')
85 y=np.asarray(open('y_0').read().splitlines(),dtype='f8')
86 #kx=np.asarray(open('px_0').read().splitlines(),dtype='f8')
87 #ky=np.asarray(open('py_0').read().splitlines(),dtype='f8')
88 
89 kx = np.reshape( np.array( [np.linspace( 0, (xDim/2-1)*dkx, xDim/2), np.linspace( (-xDim/2-1)*dkx, -dkx, xDim/2)]), (xDim,1) )
90 ky = np.reshape( np.array( [np.linspace( 0, (yDim/2-1)*dky, yDim/2), np.linspace( (-yDim/2-1)*dky, -dky, yDim/2)]), (yDim,1) )
91 kxm, kym = np.meshgrid(kx,ky)
92 km_mag = np.sqrt( kxm**2 + kym**2 )
93 k_mag = np.sqrt( kx**2 + ky**2 )
94 kMax = max(max(k_mag))
95 
96 hbar = 1.05457e-34
97 m = 1.4431607e-25
98 
99 
104 def kinertrum(Psi, dx, i, quOn):
105 
106  kMax = np.max(np.max(kx))
107  Psi[np.where(Psi==0)] = 1e-100
108  n_r = np.abs(Psi)**2
109  n_r[np.where(n_r==0)] = 1e-100
110  cPsi = np.conj(Psi)
111  phi = np.angle(Psi)
112 
113  ph1 = np.unwrap(phi, axis=0)
114  ph2 = np.unwrap(phi, axis=1)
115 
116  vel_ph1_x, vel_ph1_y = np.gradient(ph1,dx,dy)
117  vel_ph2_x, vel_ph2_y = np.gradient(ph2,dx,dy)
118 
119  v_x = (hbar/m)*vel_ph1_x;
120  v_y = (hbar/m)*vel_ph2_y;
121  v_x[np.where(v_x==0)] = 1e-100
122  v_y[np.where(v_y==0)] = 1e-100
123 
124  u_x = np.multiply(np.abs(Psi),v_x)
125  u_y = np.multiply(np.abs(Psi),v_y)
126 
127  if quOn:
128  u_x = np.multiply(u_x,np.exp(1j*np.angle(Psi)))
129  u_y = np.multiply(u_y,np.exp(1j*np.angle(Psi)))
130 
131  u_kx = np.fft.fftn(u_x)
132  u_ky = np.fft.fftn(u_y)
133 
134  uc_kx = ( kxm**2*u_kx + kxm*kym*u_ky ) / ( km_mag**2 + 1e-100 )
135  uc_ky = ( kym*kxm*u_kx + kym**2*u_ky ) / ( km_mag**2 + 1e-100 )
136 
137  ui_kx = u_kx - uc_kx
138  ui_ky = u_ky - uc_ky
139 
140  uc_x = np.fft.ifftn(uc_kx)
141  uc_y = np.fft.ifftn(uc_ky)
142  ui_x = np.fft.ifftn(ui_kx)
143  ui_y = np.fft.ifftn(ui_ky)
144 
145  Ec = 0.5*np.abs(np.square(uc_x) + np.square(uc_y))
146  Ei = 0.5*np.abs(np.square(ui_x) + np.square(ui_y))
147 
148  fig, ax = plt.subplots()
149  f = plt.imshow((Ec),cmap=plt.get_cmap('gnuplot2'))
150  cbar = fig.colorbar(f)
151  plt.gca().invert_yaxis()
152  plt.savefig("Ec_" + str(i/incr) + ".png",dpi=200)
153  plt.close()
154  fig, ax = plt.subplots()
155  f = plt.imshow((Ei),cmap=plt.get_cmap('gnuplot2'))
156  cbar = fig.colorbar(f)
157  plt.gca().invert_yaxis()
158  plt.savefig("Ei_" + str(i/incr) + ".png",dpi=200)
159  plt.close()
160 
161  print Ec
162  #exit()
163  ekc = np.zeros((xDim/2-1,1))
164  eki = np.zeros((xDim/2-1,1))
165  for i1 in np.arange(0,np.size(k_mag)/2 -2):
166  iX = np.array(np.where(np.logical_and( k_mag[i1] >= km_mag, k_mag[i1+1] < km_mag )))
167 # Ei_kx = np.sum(np.sum(np.abs(ui_kx[iX]**2*k[iX]))
168 # Ei_ky = np.sum(np.sum(np.abs(ui_ky[iX]**2*k[iX]))
169  ekc[i1] = (0.5*m*k_mag[i1]) * (np.sum(np.abs(uc_kx[iX]**2 + uc_ky[iX]**2)))/len(iX)
170  eki[i1] = (0.5*m*k_mag[i1]) * (np.sum(np.abs(ui_kx[iX]**2 + ui_ky[iX]**2)))/len(iX)
171  print i1
172  np.savetxt('ekc_' + str(i) + '.csv',ekc,delimiter=',')
173  np.savetxt('eki_' + str(i) + '.csv',eki,delimiter=',')
174  fig, ax = plt.subplots()
175  print eki[0:(xDim/2-2)]
176  f = plt.loglog(np.ravel(k_mag[0:(xDim/2 -2)]),eki[0:(xDim/2-2)])
177  plt.savefig("eki_" + str(i) + ".png",dpi=200)
178  f = plt.loglog(np.ravel(k_mag[0:(xDim/2 -2)]),np.ravel(ekc[0:(xDim/2-2)]))
179  plt.savefig("ekc_" + str(i) + ".png",dpi=200)
180  plt.close()
181 
182 
183 def kinertrum_loop(dataName, initValue, finalValue, incr):
184  for i in range(initValue,incr*(finalValue/incr),incr):
185  if os.path.exists(dataName + '_' + str(i)):
186  real=open(dataName + '_' + str(i)).read().splitlines()
187  img=open(dataName + 'i_' + str(i)).read().splitlines()
188  a_r = numpy.asanyarray(real,dtype='f8') #64-bit double
189  a_i = numpy.asanyarray(img,dtype='f8') #64-bit double
190  a = a_r[:] + 1j*a_i[:]
191 
192  kinertrum(np.reshape(a,(xDim,yDim)),dx,i,1)
193 
194 def dens_struct_fact(dataName, initValue, finalValue,incr):
195  n_k=np.zeros(finalValue/incr)
196  n_k_t=np.zeros((finalValue/incr,xDim,yDim),dtype=np.complex128)
197  for i in range(initValue,incr*(finalValue/incr),incr):
198  if os.path.exists(dataName + '_' + str(i)):
199  real=open(dataName + '_' + str(i)).read().splitlines()
200  img=open(dataName + 'i_' + str(i)).read().splitlines()
201  a_r = numpy.asanyarray(real,dtype='f8') #64-bit double
202  a_i = numpy.asanyarray(img,dtype='f8') #64-bit double
203  a = a_r[:] + 1j*a_i[:]
204  n = np.abs(a)**2
205 
206  kinertrum(np.reshape(a,(xDim,yDim)),dx,i,0)
207  sf = np.fft.fftshift(np.fft.fft2(np.reshape(n,(xDim,yDim))))
208  n_k_t[i/incr][:][:] = sf[:][:];
209  n_k[i/incr]=(abs(np.sum(np.sum(sf))*dkx**2))
210 
211  fig, ax = plt.subplots()
212  f = plt.imshow(np.log10(abs(sf)),cmap=plt.get_cmap('gnuplot2'))
213  cbar = fig.colorbar(f)
214  plt.gca().invert_yaxis()
215  plt.savefig("struct_" + str(i/incr) + ".png",vmin=0,vmax=12,dpi=200)
216  plt.close()
217  print i/incr
218 
219  np.savetxt('Struct' + '.csv',n_k,delimiter=',')
220  plt.plot(range(initValue,finalValue,incr),n_k)
221  sp.io.savemat('Struct_t.mat',mdict={'n_k_t',n_k_t})
222  plt.savefig("Struct.pdf",dpi=200)
223  plt.close()
224 
225 V = np.array(open('V_0').read().splitlines(),dtype='f8')
226 V = np.reshape(V,(xDim,yDim))
227 K = np.array(open('K_0').read().splitlines(),dtype='f8')
228 K = np.reshape(K,(xDim,yDim))
229 xPy = np.array(open('xPy_0').read().splitlines(),dtype='f8')
230 xPy = np.reshape(xPy,(xDim,yDim))
231 yPx = np.array(open('yPx_0').read().splitlines(),dtype='f8')
232 yPx = np.reshape(yPx,(xDim,yDim))
233 g = (0.5*N)*4.0*HBAR*HBAR*PI*(4.67e-9/mass)*np.sqrt(mass*omegaZ/(2.0*PI*HBAR))
234 
235 def energy_total(dataName, initValue, finalValue, increment):
236  E=np.zeros((finalValue,1))
237  E_k=np.zeros((finalValue,1))
238  E_vi=np.zeros((finalValue,1))
239  E_l=np.zeros((finalValue,1))
240  for i in range(initValue,incr*(finalValue/incr),incr):
241  if os.path.exists(dataName + '_' + str(i)):
242  real=open(dataName + '_' + str(i)).read().splitlines()
243  img=open(dataName + 'i_' + str(i)).read().splitlines()
244  a_r = np.array(real,dtype='f8') #64-bit double
245  a_i = np.array(img,dtype='f8') #64-bit double
246  wfcr = np.reshape(a_r[:] + 1j*a_i[:],(xDim,yDim))
247  wfcp = np.array(np.fft.fft2(wfcr))
248  wfcr_c = np.conj(wfcr)
249 
250  E1 = np.fft.ifft2(K*wfcp)
251  E2 = (V + 0.5*g*np.abs(wfcr)**2)*wfcr
252  E3 = -(omega*omegaX)*(np.fft.ifft(xPy*np.fft.fft(wfcr,axis=0),axis=0) - np.fft.ifft(yPx*np.fft.fft(wfcr,axis=1),axis=1) )
253 
254  E_k[i/incr] = np.trapz(np.trapz(wfcr_c*E1))*dx*dy
255  E_vi[i/incr] = np.trapz(np.trapz(wfcr_c*E2))*dx*dy
256  E_l[i/incr] = np.trapz(np.trapz(wfcr_c*E3))*dx*dy
257  E[i/incr] = E_k[i/incr] + E_vi[i/incr] + E_l[i/incr]
258  print (i/float(evMaxVal))
259  np.savetxt('E_'+ str(i) + '.csv',E,delimiter=',')
260  np.savetxt('E_k_'+ str(i) + '.csv',E_k,delimiter=',')
261  np.savetxt('E_vi_'+ str(i) + '.csv',E_vi,delimiter=',')
262  np.savetxt('E_l_'+ str(i) + '.csv',E_l,delimiter=',')
263  t = np.array(range(initValue,finalValue,incr))/dt
264  plt.plot(t,E,'r-',t,E_k,'g-',t,E_vi,'b-',t,E_l,'y-')
265  plt.savefig("EnergyVst.pdf",dpi=200)
266  plt.close()
267 
268 def energy_kinetic(dataName, initValue, finalValue, increment):
269  px1 = np.fft.fftshift(px)
270  py1 = np.fft.fftshift(py)
271  dk=[]
272  dk2[:] = (px1[:]**2 + py1[:]**2)
273  Lz = np.zeros( (finalValue/incr))
274  for i in range(initValue,incr*(finalValue/incr),incr):
275  if os.path.exists(dataName + '_' + str(i)):
276  real=open(dataName + '_' + str(i)).read().splitlines()
277  img=open(dataName + 'i_' + str(i)).read().splitlines()
278  a_r = numpy.asanyarray(real,dtype='f8') #64-bit double
279  a_i = numpy.asanyarray(img,dtype='f8') #64-bit double
280  a = a_r[:] + 1j*a_i[:]
281  wfcp = np.fft.fft2(np.reshape(a,(xDim,yDim)))
282  conjwfcp = np.conj(wfcp)
283  E_k = np.zeros(len(px1))
284  for ii in range(0,len(px1)):
285  E_k[ii] = np.sum( np.sum( np.multiply(wfcp,conjwfcp) ) )*dk2[ii]
286 
287  np.savetxt('E_k_' + str(i) + '.csv',E_k,delimiter=',')
288  print i
289 
290 def energy_potential(dataName, initValue, finalValue, increment):
291  print 'energy'
292 
293 def ang_mom(dataName, initValue, finalValue, incr, ev_type, imgdpi):
294  xm, ym = np.meshgrid(x,y)
295  pxm, pym = np.meshgrid(px,py)
296  dx2=dx**2
297  Lz = np.zeros( (finalValue/incr))
298  for i in range(initValue,incr*(finalValue/incr),incr):
299  if os.path.exists(dataName + '_' + str(i)):
300  real=open(dataName + '_' + str(i)).read().splitlines()
301  img=open(dataName + 'i_' + str(i)).read().splitlines()
302  a_r = numpy.asanyarray(real,dtype='f8') #64-bit double
303  a_i = numpy.asanyarray(img,dtype='f8') #64-bit double
304  a = a_r[:] + 1j*a_i[:]
305  wfc = np.reshape(a,(xDim,yDim))
306  conjwfc = np.conj(wfc)
307 
308  wfc_ypx = np.multiply(ym,np.fft.ifft(np.multiply(pxm,np.fft.fft(wfc,axis=1)),axis=1))
309  wfc_xpy = np.multiply(xm,np.fft.ifft(np.multiply(pym,np.fft.fft(wfc,axis=0)),axis=0))
310  result = np.sum( np.sum( np.multiply(conjwfc,wfc_xpy - wfc_ypx) ) )*dx2
311  else:
312  print "Skipped " + dataName + "_"+ str(i)
313  result = np.nan
314 
315  print i, incr
316  Lz[(i/incr)] = np.real(result)
317  type=""
318  if ev_type == 0:
319  type = "gnd"
320  else:
321  type = "ev"
322  np.savetxt('Lz.csv',Lz,delimiter=',')
323 
324  plt.plot(Lz)
325  plt.savefig("Lz_"+type+".pdf",dpi=imgdpi)
326  plt.axis('off')
327  plt.savefig("Lz_"+type+"_axis0.pdf",bbox_inches='tight',dpi=imgdpi)
328  plt.close()
329 
330 def expec_val_monopole(dataName, initValue, finalValue, incr):
331  x=np.asarray(open('x_0').read().splitlines(),dtype='f8')
332  y=np.asarray(open('y_0').read().splitlines(),dtype='f8')
333 # px=open('px_0')
334 # py=open('py_0')
335  xm, ym = np.meshgrid(x, y)
336  result = []
337  for i in range(initValue,finalValue,incr):
338  if not os.path.exists(dataName):
339  real=open(dataName + '_' + str(i)).read().splitlines()
340  img=open(dataName + 'i_' + str(i)).read().splitlines()
341  a_r = numpy.asanyarray(real,dtype='f8') #64-bit double
342  a_i = numpy.asanyarray(img,dtype='f8') #64-bit double
343  a = a_r[:] + 1j*a_i[:]
344  wfc = np.reshape(a,(xDim,yDim))
345  conjwfc = np.conj(wfc)
346 
347  d1 = np.multiply( np.square(xm) + np.square(ym), wfc )
348  d2 = np.multiply( conjwfc, d1)
349  result.append( np.real( np.sum( np.sum( d2 ) ) )*dx*dx )
350  print str(100*float(i)/finalValue) + '%'
351  np.savetxt('monopole.csv',result,delimiter=',')
352  plt.plot(range(initValue,finalValue,incr),result)
353  plt.savefig("Monopole.png",dpi=200)
354  plt.close()
355 
356 def expec_val_quadrupole(dataName, initValue, finalValue, incr):
357  x=np.asarray(open('x_0').read().splitlines(),dtype='f8')
358  y=np.asarray(open('y_0').read().splitlines(),dtype='f8')
359 # px=open('px_0')
360 # py=open('py_0')
361  xm, ym = np.meshgrid(x, y)
362  result = []
363  for i in range(initValue,finalValue,incr):
364  if not os.path.exists(dataName):
365  real=open(dataName + '_' + str(i)).read().splitlines()
366  img=open(dataName + 'i_' + str(i)).read().splitlines()
367  a_r = numpy.asanyarray(real,dtype='f8') #64-bit double
368  a_i = numpy.asanyarray(img,dtype='f8') #64-bit double
369  a = a_r[:] + 1j*a_i[:]
370  wfc = np.reshape(a,(xDim,yDim))
371  conjwfc = np.conj(wfc)
372 
373  d1 = np.multiply( np.square(xm) - np.square(ym), wfc )
374  d2 = np.multiply( conjwfc, d1)
375  result.append( np.real( np.sum( np.sum( d2 ) ) )*dx*dx )
376  print str(100*float(i)/finalValue) + '%'
377  np.savetxt('quadrupole.csv',result,delimiter=',')
378  plt.plot(range(initValue,finalValue,incr),result)
379  plt.savefig("Quadrupole.png",dpi=200)
380  plt.close()
381 
382 def expec_val_(quant_name, quantity, dataName, initValue, finalValue, incr):
383  x=np.asarray(open('x_0').read().splitlines(),dtype='f8')
384  y=np.asarray(open('y_0').read().splitlines(),dtype='f8')
385  xm, ym = np.meshgrid(x, y)
386  result = []
387  for i in range(initValue,finalValue,incr):
388  if not os.path.exists(dataName):
389  real=open(dataName + '_' + str(i)).read().splitlines()
390  img=open(dataName + 'i_' + str(i)).read().splitlines()
391  a_r = numpy.asanyarray(real,dtype='f8') #64-bit double
392  a_i = numpy.asanyarray(img,dtype='f8') #64-bit double
393  a = a_r[:] + 1j*a_i[:]
394  wfc = np.reshape(a,(xDim,yDim))
395  conjwfc = np.conj(wfc)
396 
397  d1 = np.multiply( quantity, wfc )
398  d2 = np.multiply( conjwfc, d1)
399  result.append( np.real( np.sum( np.sum( d2 ) ) )*dx*dx )
400  print str(100*float(i)/finalValue) + '%'
401  np.savetxt(quant_name + '.csv',result,delimiter=',')
402  plt.plot(range(initValue,finalValue,incr),result)
403  plt.savefig(quant_name + ".pdf",dpi=200)
404  plt.close()
405 
406 if __name__ == '__main__':
407  kinertrum_loop('wfc_ev', 0, evMaxVal, incr)
408  energy_total('wfc_ev',0,evMaxVal,incr)
409  dens_struct_fact('wfc_ev', 0, evMaxVal, 500)
410 
411  energy_kinetic('wfc_ev', 0, evMaxVal, 200)
412  ang_mom('wfc_ev', 0, evMaxVal, incr, 1, 200)
413  expec_val_monopole('wfc_ev',0,evMaxVal,incr)
414  expec_val_quadrupole('wfc_ev',0,evMaxVal,incr)
def energy_total(dataName, initValue, finalValue, increment)
Definition: observables.py:235
def energy_potential(dataName, initValue, finalValue, increment)
Definition: observables.py:290
def dens_struct_fact(dataName, initValue, finalValue, incr)
Definition: observables.py:194
def expec_val_monopole(dataName, initValue, finalValue, incr)
Definition: observables.py:330
def energy_kinetic(dataName, initValue, finalValue, increment)
Definition: observables.py:268
def kinertrum_loop(dataName, initValue, finalValue, incr)
Definition: observables.py:183
def expec_val_quadrupole(dataName, initValue, finalValue, incr)
Definition: observables.py:356
def ang_mom(dataName, initValue, finalValue, incr, ev_type, imgdpi)
Definition: observables.py:293
def kinertrum(Psi, dx, i, quOn)
Kinetic energy spectrum = kinertrum.
Definition: observables.py:104
def expec_val_(quant_name, quantity, dataName, initValue, finalValue, incr)
Definition: observables.py:382