Calculations/DipolarGasSimulator/+Simulator/@Solver/Initialize.m

64 lines
2.1 KiB
Matlab

function [psi,V,VDk] = Initialize(Params,Transf)
format long
X = Transf.X; Y = Transf.Y; Z = Transf.Z;
Zcutoff = Params.Lz/2;
% == Potential == %
V = 0.5*(Params.gx.*X.^2+Params.gy.*Y.^2+Params.gz*Z.^2);
% == Calculating the DDIs == %
% For a cylindrical cutoff, we first construct a kr grid based on the 3D parameters using Bessel quadrature
loadDDI = 1;
if loadDDI == 1
VDk = load(sprintf('./Data/VDk_M.mat'));
VDk = VDk.VDk;
else
Params.Lr = 0.5*min(Params.Lx,Params.Ly);
Params.Nr = max(Params.Nx,Params.Ny);
[TransfRad] = Simulator.SetupSpaceRadial(Params); %morder really doesn't matter
VDk = Simulator.VDcutoff(TransfRad.kr,TransfRad.kz,TransfRad.Rmax,Zcutoff);
disp('Calculated radial grid and cutoff')
% VDk = interp2(DDI.kz,DDI.kr,DDI.VDk,Transf.kz,Transf.kr,'spline');
fullkr = [-flip(TransfRad.kr)',TransfRad.kr'];
[KR,KZ] = ndgrid(fullkr,TransfRad.kz);
[KX3D,KY3D,KZ3D] = ndgrid(ifftshift(Transf.kx),ifftshift(Transf.ky),ifftshift(Transf.kz));
KR3D = sqrt(KX3D.^2 + KY3D.^2);
fullVDK = [flip(VDk',2),VDk']';
VDk = interpn(KR,KZ,fullVDK,KR3D,KZ3D,'spline',-1/3); %Interpolating the radial VDk onto a new grid
VDk = fftshift(VDk);
save(sprintf('./Data/VDk_M.mat'),'VDk');
end
disp('Finished DDI')
% == Setting up the initial wavefunction == %
ellx = sqrt(Params.hbar/(Params.m*Params.wx))/Params.l0;
elly = sqrt(Params.hbar/(Params.m*Params.wy))/Params.l0;
ellz = sqrt(Params.hbar/(Params.m*Params.wz))/Params.l0;
Rx = 4*sqrt(2)*ellx;
Ry = 4*sqrt(2)*elly;
Rz = sqrt(2)*ellz;
X0 = 0.0*Transf.Xmax;
Y0 = 0.0*Transf.Ymax;
Z0 = 0*Transf.Zmax;
psiz = exp(-(Z-Z0).^2/Rz^2)/sqrt(ellz*sqrt(pi));
psi2d = load(sprintf('./Data/Seed/psi_2d_SS.mat'),'psiseed_2d'); psi2d = psi2d.psiseed_2d;
psi = psiz.*repmat(psi2d,[1 1 length(Transf.z)]);
% Add some noise
r = normrnd(0,1,size(X));
theta = rand(size(X));
noise = r.*exp(2*pi*1i*theta);
psi = psi + 0.00*noise;
Norm = trapz(abs(psi(:)).^2)*Transf.dx*Transf.dy*Transf.dz;
psi = sqrt(Params.N)*psi/sqrt(Norm);
end