Added script to calculate VDK without cutoff for determining optical unit cell size

This commit is contained in:
Karthik 2025-03-14 01:38:42 +01:00
parent 8994dd76f8
commit 201e4a39c6
5 changed files with 38 additions and 10 deletions

View File

@ -5,7 +5,7 @@ OptionsStruct = struct;
OptionsStruct.NumberOfAtoms = 5E5; OptionsStruct.NumberOfAtoms = 5E5;
OptionsStruct.DipolarPolarAngle = 0; OptionsStruct.DipolarPolarAngle = 0;
OptionsStruct.DipolarAzimuthAngle = 0; OptionsStruct.DipolarAzimuthAngle = 0;
OptionsStruct.ScatteringLength = 81.0; OptionsStruct.ScatteringLength = 78.0;
AspectRatio = 2.8; AspectRatio = 2.8;
HorizontalTrapFrequency = 125; HorizontalTrapFrequency = 125;
@ -15,7 +15,7 @@ OptionsStruct.TrapPotentialType = 'Harmonic';
OptionsStruct.NumberOfGridPoints = [256, 256]; OptionsStruct.NumberOfGridPoints = [256, 256];
OptionsStruct.Dimensions = [18, 18]; OptionsStruct.Dimensions = [18, 18];
OptionsStruct.TimeStepSize = 0.001; % in s OptionsStruct.TimeStepSize = 1E-3; % in s
OptionsStruct.MinimumTimeStepSize = 1E-5; % in s OptionsStruct.MinimumTimeStepSize = 1E-5; % in s
OptionsStruct.TimeCutOff = 2E6; % in s OptionsStruct.TimeCutOff = 2E6; % in s
OptionsStruct.EnergyTolerance = 5E-10; OptionsStruct.EnergyTolerance = 5E-10;
@ -23,7 +23,7 @@ OptionsStruct.ResidualTolerance = 1E-04;
OptionsStruct.NoiseScaleFactor = 0.05; OptionsStruct.NoiseScaleFactor = 0.05;
OptionsStruct.MaxIterations = 10; OptionsStruct.MaxIterations = 10;
OptionsStruct.VariationalWidth = 1.8; OptionsStruct.VariationalWidth = 2.0;
OptionsStruct.WidthLowerBound = 0.01; OptionsStruct.WidthLowerBound = 0.01;
OptionsStruct.WidthUpperBound = 12; OptionsStruct.WidthUpperBound = 12;
OptionsStruct.WidthCutoff = 1e-2; OptionsStruct.WidthCutoff = 1e-2;
@ -43,14 +43,14 @@ solver.Potential = pot.trap();
%-% Run Solver %-% %-% Run Solver %-%
[Params, Transf, psi, V, VDk] = solver.run(); [Params, Transf, psi, V, VDk] = solver.run();
%% AR = 4 %% AR = 4.0
OptionsStruct = struct; OptionsStruct = struct;
OptionsStruct.NumberOfAtoms = 5E5; OptionsStruct.NumberOfAtoms = 5E5;
OptionsStruct.DipolarPolarAngle = 0; OptionsStruct.DipolarPolarAngle = 0;
OptionsStruct.DipolarAzimuthAngle = 0; OptionsStruct.DipolarAzimuthAngle = 0;
OptionsStruct.ScatteringLength = 81.0; OptionsStruct.ScatteringLength = 78.0;
AspectRatio = 4.0; AspectRatio = 4.0;
HorizontalTrapFrequency = 125; HorizontalTrapFrequency = 125;
@ -60,7 +60,7 @@ OptionsStruct.TrapPotentialType = 'Harmonic';
OptionsStruct.NumberOfGridPoints = [256, 256]; OptionsStruct.NumberOfGridPoints = [256, 256];
OptionsStruct.Dimensions = [18, 18]; OptionsStruct.Dimensions = [18, 18];
OptionsStruct.TimeStepSize = 0.001; % in s OptionsStruct.TimeStepSize = 1E-3; % in s
OptionsStruct.MinimumTimeStepSize = 1E-5; % in s OptionsStruct.MinimumTimeStepSize = 1E-5; % in s
OptionsStruct.TimeCutOff = 2E6; % in s OptionsStruct.TimeCutOff = 2E6; % in s
OptionsStruct.EnergyTolerance = 5E-10; OptionsStruct.EnergyTolerance = 5E-10;
@ -68,7 +68,7 @@ OptionsStruct.ResidualTolerance = 1E-04;
OptionsStruct.NoiseScaleFactor = 0.05; OptionsStruct.NoiseScaleFactor = 0.05;
OptionsStruct.MaxIterations = 10; OptionsStruct.MaxIterations = 10;
OptionsStruct.VariationalWidth = 1.5; OptionsStruct.VariationalWidth = 2.0;
OptionsStruct.WidthLowerBound = 0.01; OptionsStruct.WidthLowerBound = 0.01;
OptionsStruct.WidthUpperBound = 12; OptionsStruct.WidthUpperBound = 12;
OptionsStruct.WidthCutoff = 1e-2; OptionsStruct.WidthCutoff = 1e-2;

View File

@ -31,6 +31,7 @@ for Lx = 2:0.4:6
OptionsStruct.EnergyTolerance = 5E-10; OptionsStruct.EnergyTolerance = 5E-10;
OptionsStruct.ResidualTolerance = 1E-05; OptionsStruct.ResidualTolerance = 1E-05;
OptionsStruct.NoiseScaleFactor = 0.05; OptionsStruct.NoiseScaleFactor = 0.05;
OptionsStruct.IncludeDDICutOff = false;
OptionsStruct.MaxIterations = 10; OptionsStruct.MaxIterations = 10;
OptionsStruct.VariationalWidth = 1.00; OptionsStruct.VariationalWidth = 1.00;

View File

@ -0,0 +1,19 @@
function VDk = calculateVDkWithoutCutoff(~, Transf, Params, ell)
% == Calculating the DDI potential in Fourier space with appropriate cutoff == %
% Interaction in K space
QX = Transf.KX*ell/sqrt(2);
QY = Transf.KY*ell/sqrt(2);
Qsq = QX.^2 + QY.^2;
absQ = sqrt(Qsq);
QDsq = QX.^2*cos(Params.eta)^2 + QY.^2*sin(Params.eta)^2;
% Bare interaction
Fpar = -1 + 3*sqrt(pi)*QDsq.*erfcx(absQ)./absQ; % Scaled complementary error function erfcx(x) = e^(x^2) * erfc(x)
Fperp = 2 - 3*sqrt(pi).*absQ.*erfcx(absQ);
Fpar(absQ == 0) = -1;
% Full DDI
VDk = Fpar*sin(Params.theta)^2 + Fperp*cos(Params.theta)^2;
end

View File

@ -28,6 +28,7 @@ classdef DipolarGas < handle & matlab.mixin.Copyable
Calculator; Calculator;
SimulationParameters; SimulationParameters;
IncludeDDICutOff;
PlotLive; PlotLive;
JobNumber; JobNumber;
RunOnGPU; RunOnGPU;
@ -84,6 +85,8 @@ classdef DipolarGas < handle & matlab.mixin.Copyable
@(x) assert(isnumeric(x) && isscalar(x) && (x > 0))); @(x) assert(isnumeric(x) && isscalar(x) && (x > 0)));
addParameter(p, 'JobNumber', 0,... addParameter(p, 'JobNumber', 0,...
@(x) assert(isnumeric(x) && isscalar(x) && (x >= 0))); @(x) assert(isnumeric(x) && isscalar(x) && (x >= 0)));
addParameter(p, 'IncludeDDICutOff', true,...
@islogical);
addParameter(p, 'PlotLive', false,... addParameter(p, 'PlotLive', false,...
@islogical); @islogical);
addParameter(p, 'RunOnGPU', false,... addParameter(p, 'RunOnGPU', false,...
@ -119,6 +122,7 @@ classdef DipolarGas < handle & matlab.mixin.Copyable
this.VariationalWidthTolerance = p.Results.VariationalWidthTolerance; this.VariationalWidthTolerance = p.Results.VariationalWidthTolerance;
this.VariationalEnergyTolerance = p.Results.VariationalEnergyTolerance; this.VariationalEnergyTolerance = p.Results.VariationalEnergyTolerance;
this.IncludeDDICutOff = p.Results.IncludeDDICutOff;
this.PlotLive = p.Results.PlotLive; this.PlotLive = p.Results.PlotLive;
this.JobNumber = p.Results.JobNumber; this.JobNumber = p.Results.JobNumber;
this.RunOnGPU = p.Results.RunOnGPU; this.RunOnGPU = p.Results.RunOnGPU;

View File

@ -4,7 +4,11 @@ function [psi,V,VDk] = initialize(this,Params,VParams,Transf)
assert(~anynan(V), 'Potential not defined! Specify as <SimulatorObject>.Potential = <PotentialsObject>.trap() + <AdditionalTerms>.'); assert(~anynan(V), 'Potential not defined! Specify as <SimulatorObject>.Potential = <PotentialsObject>.trap() + <AdditionalTerms>.');
% == Calculating the DDIs == % % == Calculating the DDIs == %
if this.IncludeDDICutOff
VDk = this.Calculator.calculateVDkWithCutoff(Transf, Params, VParams.ell); VDk = this.Calculator.calculateVDkWithCutoff(Transf, Params, VParams.ell);
else
VDk = this.Calculator.calculateVDkWithoutCutoff(Transf, Params, VParams.ell);
end
% == Setting up the initial wavefunction == % % == Setting up the initial wavefunction == %
psi = this.setupWavefunction(Params,Transf); psi = this.setupWavefunction(Params,Transf);