2 % Copyright (C) 2009-2015 Peter Rakyta, Ph.D.
4 % This program is free software: you can redistribute it and/or modify
5 % it under the terms of the GNU General Public License as published by
6 % the Free Software Foundation, either version 3 of the License, or
7 % (at your option) any later version.
9 % This program is distributed in the hope that it will be useful,
10 % but WITHOUT ANY WARRANTY; without even the implied warranty of
11 % MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12 % GNU General Public License
for more details.
14 % You should have received a copy of the GNU General Public License
15 % along with
this program. If not, see http:
17 %> @addtogroup utilities Utilities
20 %> @brief A
class to calculate the DC Josephson current.
22 %> @brief A
class to calculate the DC Josephson current.
27 %> The magnetic field in Tesla (OBSOLETE)
33 methods (Access=
public)
35 %% Contructor of the
class 36 %> @brief Constructor of the
class.
38 %> @
param varargin Cell array of optional parameters. For details see #InputParsing.
39 %> @
return An instance of the
class 46 obj.scatterPotential = [];
47 obj.gfininvfromHamiltonian =
false;
49 obj.useSelfEnergy =
true;
51 if strcmpi(
class(obj),
'SNSJosephson')
52 obj.InputParsing( varargin{:} );
60 %% CurrentCalc_discrete
61 %> @brief Calculates the Josephson current
using the method in PRB 93, 224510 (2016)
62 %> @
param DeltaPhi_vec Array of phase diffrencies between the superconducting contacts.
63 %> @
param varargin Cell array of optional parameters (https:
64 %> @
param 'range' String. Set
'ABS' (
default)
for the Andreev bound states,
'NBS' for the normal bound states,
'CONT' for the scattering states, or
'ALL' for calculating all ranges at once.
65 %> @
param 'Edb' The number of the energy points over the contour integral. (
default value is 511)
66 %> @
param 'DeltaPhi' A parameter to control the incident angle of the integration contour near the real axis. (Default value is \f$\Delta\Phi=0.5\pi\f$).
67 %> @
param 'Evec' An array containing the complex energy points in
case of custom integration path.
68 %> @
param 'Emin' A minimum of the energy array on the real axis.
69 %> @
param 'Emax' A maximum of the energy array on the real axis.
70 %> @
param 'ProximityOn' Logical value. Set
true to recalculate the Green
operator of the scattering center
for every phase difference. Useful
for self-consistent modelling of the proximity effect.
71 %> @
return An array of the calculated current corresponding to the phase differencies DeltaPhi_vec
72 function currentvec = CurrentCalc_discrete( obj, DeltaPhi_vec, varargin )
75 p.addParameter(
'range',
'ABS' ); % NBS: normal bound states, ABS: andreev bound states, CONT: scattering range, ALL: all contributions at once
76 p.addParameter(
'Edb', 511 ); % number of energy points on the contour
77 p.addParameter(
'DeltaPhi', 0.5*pi );
78 p.addParameter(
'Evec', [] );
79 p.addParameter(
'Emin', [] );
80 p.addParameter(
'Emax', [] );
81 p.addParameter(
'ProximityOn',
false ); %New
if proxymity effect is taken into account, the current
operator must be recalculated
for each phase difference
84 range = p.Results.range;
86 Evec = p.Results.Evec;
87 Emin = p.Results.Emin;
88 Emax = p.Results.Emax;
89 DeltaPhi = p.Results.DeltaPhi;
90 ProximityOn = p.Results.ProximityOn;
92 DeltaPhi_vec = [0, DeltaPhi_vec];
93 DeltaPhi_vec = diff(DeltaPhi_vec);
95 % reset the phase difference to zero
100 obj.gfininvfromHamiltonian =
true;
104 phivec = (0.5 + 0.5*sin( -0.5*pi:pi/(Edb+1):0.5*pi )) * DeltaPhi + (pi-DeltaPhi)/2;
106 if strcmp( range,
'NBS' ) && isempty(Evec) && (isempty(Emax) && isempty(Emin))
107 obj.display([
'EQuUs:Utils:',
class(obj),
':CurrentCalc_discrete: Calculating current for the normal bound states.'], 1)
109 Emax = -obj.BandWidth.Lead.Emax;
110 Emin = -obj.BandWidth.Scatter.Emax;
111 % phivec = (0.5 + 0.5*sin( -0.5*pi:pi/(Edb+1):0.5*pi )) * pi;
112 %fact = -fact; %the contour integral is in the opposite direction as
for the ABS
113 elseif strcmp( range,
'ABS' ) && isempty(Evec) && (isempty(Emax) && isempty(Emin))
114 obj.display([
'EQuUs:Utils:',
class(obj),
':CurrentCalc_discrete: Calculating current for the Andreev bound states.'], 1)
116 param = obj.junction.FL_handles.Read(
'param' );
118 Emin = -pair_potential;
120 %phivec = (0.5 + 0.5*sin( -0.5*pi:pi/(Edb+1):0.5*pi )) * pi;
121 elseif strcmp( range,
'CONT' ) && isempty(Evec) && (isempty(Emax) && isempty(Emin))
122 obj.display([
'EQuUs:Utils:',
class(obj),
':CurrentCalc_discrete: Calculating current for the scattering states.'], 1)
124 param = obj.junction.FL_handles.Read(
'param' );
126 Emin = -obj.BandWidth.Lead.Emax;
127 Emax = -pair_potential;
128 %phivec = (0.5 + 0.5*sin( -0.5*pi:pi/(Edb+1):0.5*pi )) * pi;
129 elseif strcmp( range,
'ALL' ) && isempty(Evec) && (isempty(Emax) && isempty(Emin))
130 obj.display([
'EQuUs:Utils:',
class(obj),
':CurrentCalc_discrete: Calculating current for all the states.'], 1)
132 Emin = min([-obj.BandWidth.Lead.Emax,-obj.BandWidth.Scatter.Emax]);
134 %phivec = (0.5 + 0.5*sin( -0.5*pi:pi/(Edb+1):0.5*pi )) * pi;
135 elseif ~isempty(Evec)
136 obj.display([
'EQuUs:Utils:',
class(obj),
':CurrentCalc_discrete: Calculating current for custom energy contour.'], 1)
137 Emin = min(real(Evec));
138 Emax = max(real(Evec));
140 elseif (~isempty(Emax) && ~isempty(Emin))
142 err = MException([
'EQuUs:Utils:',
class(obj),
':CurrentCalc_discrete',
'Unrecognized energy range']);
143 save(
'Error_SNSJosephson_CurrentCalc_discrete.mat');
148 currentvec = zeros(size(DeltaPhi_vec));
152 % E -> -E symmetry
for both zero and nonzero temperature.
155 % creating the energy array
156 radius = abs( Emax-Emin) /2;
157 R = radius/sin(DeltaPhi/2);
159 Evec = R*(cos(phivec) + 1i*sin(phivec)) - radius + Emax - 1i*R*sin((pi-DeltaPhi)/2);
162 obj.display([
'EQuUs:Utils:',
class(obj),
':CurrentCalc_discrete: Calculating the current between Emin = ', num2str(Emin),
' and Emax = ',num2str(Emax),],1);
163 deltaEvec = [0, diff(Evec)];
164 currentsurf = complex(zeros( length(Evec), length(obj.T), length(DeltaPhi_vec) ));
166 % open the parallel pool
167 poolmanager =
Parallel( obj.junction.FL_handles.Read(
'Opt') );
170 % setting
function handles
for the parallel calculations
171 hdisplay = @obj.display;
172 hSetEnergy = @obj.SetEnergy;
173 hcreateCurrentOperator = @obj.createCurrentOperator;
174 hincreasePhaseDifference = @obj.increasePhaseDifference;
175 hcomplexCurrent = @complexCurrent;
177 junction_loc = obj.junction;
180 % starting parallel calculations
181 parfor (iidx = 1:length(Evec), poolmanager.NumWorkers)
182 %
for iidx = 1:length(Evec)
184 if isempty(junction_loc.getTransverseMomentum())
186 feval( hdisplay, [
'EQuUs:Utils:',
class(obj),
':CurrentCalc_discrete: Currently evaluating ', num2str(iidx),
' from ', num2str(length(Evec))], 1);
189 % in the first iteration the
self energies of all terminals should be recalculated
190 recalculateSurface = [1 2];
192 %
set the energy value in the current iteration
193 feval(hSetEnergy, Evec(iidx),
'junction', junction_loc );
196 % determine the current
operator and calculate the Green
operator of the scattering center
197 current_op = feval(hcreateCurrentOperator,
'junction', junction_loc);
199 currentsurf( iidx, :, : ) = NaN(length(obj.T), length(DeltaPhi_vec));
200 feval(hdisplay, [
'EQuUs:Utils:',
class(obj),
':CurrentCalc_discrete: Error in SNSJosephson:CurrentCalc_discrete at energy point E = ', num2str(Evec(iidx)),
' caused by:'],1);
201 feval(hdisplay, errCause.identifier, 1 );
202 feval(hdisplay, errCause.stack(1), 1 );
203 feval(hdisplay, [
'EQuUs:Utils:',
class(obj),
':CurrentCalc_discrete: Giving NaN for the current at the given energy point'], 1);
207 % preallocating array
for the calculated, phase resolved current values
208 currentvec2int = zeros( length( T ), length( DeltaPhi_vec ) );
209 for jjdx = 1:length( DeltaPhi_vec )
211 % increasing the phase different between the superconducting leads
212 feval(hincreasePhaseDifference, DeltaPhi_vec(jjdx), junction_loc );
213 if ProximityOn && jjdx > 1
214 %
if proximity effect is taken into account, the current
operator must be always recalculated
215 current_op = feval(hcreateCurrentOperator,
'junction', junction_loc);
218 % evaluate the current at a complex energy
219 currentvec2int( :, jjdx ) = feval(hcomplexCurrent, junction_loc, current_op, recalculateSurface);
221 % the phase different is increased by unitary transformation, so recalculating the leeds is not necessary
222 recalculateSurface = []; % see increasePhaseDifference
224 feval(hdisplay, [
'EQuUs:Utils:',
class(obj),
':CurrentCalc_discrete: Error in SNSJosephson:CurrentCalc_discrete at energy point E = ', num2str(Evec(iidx)),
'and at phase difference ', num2str(DeltaPhi_vec(jjdx)),
' caused by:'], 1);
225 feval(hdisplay, errCause.identifier, 1 );
226 for jjjdx = 1:length(errCause.stack)
227 feval(hdisplay, errCause.stack(jjjdx), 1 );
229 feval(hdisplay, [
'EQuUs:Utils:',
class(obj),
':CurrentCalc_discrete: Giving NaN for the calculated current at the given energy point and phase difference.'], 1);
230 currentvec2int( jjdx ) = NaN;
234 currentsurf( iidx, :, : ) = currentvec2int;
236 if isempty(junction_loc.getTransverseMomentum())
241 % closing the parallel pool
242 poolmanager.closePool();
245 % calculating the contour integrals
246 currentvec = cell( size(obj.T) );
247 for kkdx = 1:length( obj.T )
248 currentvec{kkdx} = zeros(size(DeltaPhi_vec));
249 for jjdx = 1:length( DeltaPhi_vec )
250 currentvec2integrate = currentsurf(:,kkdx,jjdx);
251 % omitting NaNs by substituting interpolates
252 indexes = isnan(currentvec2integrate);
253 indexesall = 1:length(currentvec2integrate);
255 currentvec2integrate(indexes) = interp1(indexesall(~indexes), currentvec2integrate(~indexes), indexesall(indexes),
'spline' );
256 current = fact*imag( obj.SimphsonInt( currentvec2integrate.*transpose(deltaEvec) ))/(pi);
257 currentvec{kkdx}( jjdx ) = currentvec{kkdx}( jjdx ) + current;
260 if length(obj.T) == 1
261 currentvec = currentvec{1};
264 obj.display([
'EQuUs:Utils:',
class(obj),
':CurrentCalc_discrete: Process finished.'], 1)
270 %------------------------------------------------
271 % calculates current
for complex energy
272 function current = complexCurrent( junction_loc, current_op, recalculateSurface )
275 if isfield(junction_loc,
'Ribbons')
276 junction_loc = junction_loc.Ribbons{1}; %
for the
CombineRibbons interface
280 G = junction_loc.CustomDysonFunc(
'recalculateSurface', recalculateSurface, ...
281 'decimate',
true,
'constant_channels',
false, ...
282 'SelfEnergy', obj.useSelfEnergy,
'keep_sites',
'scatter');
285 gfininv_size = size(gfininv,1);
286 if G_size ~= gfininv_size
287 warning( [
'EQuUs:Utils:',
class(obj),
':CurrentCalc_discrete: wrong size of Greens function'] )
292 fE = obj.Fermi(junction_loc.getEnergy());
294 scatter_sites = junction_loc.CreateH.Read(
'kulso_szabfokok');
295 if ~isempty(scatter_sites)
296 %
if the Green
operator is calculated from the whole Hamiltonian, the coordinates contain all
sites.
297 coordinates = junction_loc.CreateH.Read(
'coordinates');
298 BdG_indexes = coordinates.BdG_u(scatter_sites);
300 % obtaining only the surface points of the junction structure
301 coordinates = junction_loc.getCoordinates();
302 BdG_indexes = coordinates.BdG_u();
305 current = zeros(length( fE ), 1 );
307 for idx = 1:length( fE )
309 G_tmp( BdG_indexes, BdG_indexes) = G_tmp( BdG_indexes, BdG_indexes)*fE(idx);
310 G_tmp( ~BdG_indexes, ~BdG_indexes) = G_tmp( ~BdG_indexes, ~BdG_indexes)*(1-fE(idx));
312 current(idx) = trace( current_op*G_tmp );
318 % end nested functions
323 %% CurrentCalc_continuum
324 %> @brief Calculates the Josephson current of the continous scattering states.
325 %> @
param DeltaPhi_vec Array of phase diffrencies between the superconducting contacts.
326 %> @
param varargin Cell array of optional parameters (https:
327 %> @
param 'Edb' The number of the energy points over the contour integral. (
default value is 511)
328 %> @
return [1] An array of the calculated current corresponding to the phase differencies DeltaPhi_vec
329 %> @
return [2] Energy resolved Josephson current.
330 function [currentvec, current_surf] = CurrentCalc_continuum( obj, DeltaPhi_vec, varargin )
332 obj.display([
'EQuUs:Utils:',
class(obj),
':CurrentCalc_continuum: Calculating current in the continuum range from the scattering matrix.'], 1)
334 p.addParameter(
'Edb', 511 ); % number of energy points on the contour
335 p.parse(varargin{:});
338 DeltaPhi_vec = [0, DeltaPhi_vec];
339 DeltaPhi_vec = diff(DeltaPhi_vec);
340 currentvec = zeros(size(DeltaPhi_vec));
342 % reset the phase difference to zero
345 % determine the bandwidth
347 Emax = -obj.BandWidth.Lead.Emin;
348 Emin = -obj.BandWidth.Lead.Emax;
349 Evec = (0.5+0.5*sin((0:Edb-1)/(Edb-1)*pi - pi/2))*(Emax-Emin) + Emin;
350 deltaEvec = [0, diff(Evec) ];
352 current_surf = zeros(length(DeltaPhi_vec), length(Evec));
353 % E -> -E symmetry
for both zero and nonzero temperature
356 hCurrentCalc_continuum_forGivenEnergy = @CurrentCalc_continuum_forGivenEnergy;
357 ribbon_loc = obj.junction;%;.CreateClone();
358 for kdx = 1:length( Evec )
359 obj.display( [
'EQuUs:Utils:',
class(obj),
':CurrentCalc_continuum: Evaluating ', num2str(kdx),
' from ', num2str(length(Evec))], 1);
361 current_surf(:,kdx) = feval(hCurrentCalc_continuum_forGivenEnergy, DeltaPhi_vec, Evec(kdx), ribbon_loc );
365 for jdx = 1:length( DeltaPhi_vec )
366 currentvec2integrate = current_surf(jdx,:);
367 % omitting NaNs by substituting interpolates
368 indexes = isnan(currentvec2integrate);
369 indexesall = 1:length(currentvec2integrate);
371 currentvec2integrate(indexes) = interp1(indexesall(~indexes), currentvec2integrate(~indexes), indexesall(indexes),
'spline' );
372 current_surf(jdx,:) = currentvec2integrate;
373 currentvec(jdx) = currentvec(jdx) + fact*obj.SimphsonInt( currentvec2integrate.*(deltaEvec) )/(2*pi);
376 %reset the system
for a
new calculation
379 obj.display([
'EQuUs:Utils:',
class(obj),
':CurrentCalc_continuum: Process finished.'], 1)
381 %-----------------------------------------------------------
383 function current = CurrentCalc_continuum_forGivenEnergy( DeltaPhi_vec, E, junction_loc )
385 obj.SetEnergy( E,
'junction', junction_loc );
387 param = junction_loc.FL_handles.Read(
'param' );
390 current = NaN( size( DeltaPhi_vec ) );
391 if abs(obj.junction.getEnergy()) <= abs(pair_potential)
395 current_op = obj.createCurrentOperator(
'junction', junction_loc);
400 recalculateSurface = [1 2];
401 for idx = 1: length( DeltaPhi_vec )
404 % setting phi2 to phi2+delta_phi
405 obj.increasePhaseDifference( DeltaPhi_vec(idx), junction_loc );
406 %calculating the current
407 current_tmp = CurrentCalc_continuum_for_GivenDeltaPhi( junction_loc, current_op, recalculateSurface);
409 obj.display( [
'EQuUs:Utils:',
class(obj),
':CurrentCalc_continuum: Error in SNSJosephson:CurrentCalc_continuum_forGivenEnergy caused by:'])
410 obj.display( errCause.identifier, 1 );
411 obj.display( errCause.stack(1), 1 );
412 obj.display( [
'EQuUs:Utils:',
class(obj),
':CurrentCalc_continuum: Giving NaN for the calculated current'])
415 recalculateSurface = [];
416 current(idx) = current_tmp;
421 %-------------------------------------------------
422 function current = CurrentCalc_continuum_for_GivenDeltaPhi( ribbon_loc, current_op, recalculateSurface )
424 %Dysonfunc = @()CustomDysonFunc(
'recalculateSurface', recalculateSurface);
425 Dysonfunc = @()ribbon_loc.CustomDysonFunc(
'recalculateSurface', recalculateSurface, ...
426 'decimate',
false,
'constant_channels',
false,
'SelfEnergy', obj.useSelfEnergy,
'keep_sites',
'scatter' );
428 G = ribbon_loc.FL_handles.DysonEq(
'CustomDyson', Dysonfunc );
430 %[S,ny] = obj.junction.FL_handles.SmatrixCalc();
431 M = ribbon_loc.FL_handles.Read(
'M');
433 % The calculated leads
for calculate the scattering states
434 Leads = ribbon_loc.FL_handles.Read(
'Leads');
435 vcsop_p = cell(length(Leads),1);
436 expk_p = cell(length(Leads),1);
437 expk_m = cell(length(Leads),1);
438 modusmtx_p = cell(length(Leads),1);
439 Normamtx = cell(length(Leads),1);
440 modusmtx_m = cell(length(Leads),1);
441 d_modusmtx_m = cell(length(Leads),1);
442 for iidx = 1:length(Leads)
443 vcsop_p{iidx} = Leads{iidx}.Read(
'vcsop_p');
444 expk_p{iidx} = Leads{iidx}.Read(
'expk_p');
445 expk_m{iidx} = Leads{iidx}.Read(
'expk_m');
446 modusmtx_p{iidx} = Leads{iidx}.Read(
'modusmtx_p');
447 Normamtx{iidx} = Leads{iidx}.Read(
'Normamtx');
448 modusmtx_m{iidx} = Leads{iidx}.Read(
'modusmtx_m');
449 d_modusmtx_m{iidx} = Leads{iidx}.Read(
'd_modusmtx_m');
452 fE = obj.Fermi(ribbon_loc.getEnergy());
456 prop_indexes{1} = find(abs(abs(expk_p{1}) - 1) < tolerance);
457 prop_indexes{2} = find(abs(abs(expk_p{2}) - 1) < tolerance);
460 for cell_index = 1:length(Leads)
461 for iidx = 1:length(prop_indexes{cell_index})
462 current_tmp = CurrentFromGivenMode( cell_index, prop_indexes{cell_index}(iidx) );
463 current = current + current_tmp;
467 if imag(current) > 1e-6
468 warning([
'EQuUs:Utils:',
class(obj),
':CurrentCalc_continuum: Current is complex: ', num2str(current)])
470 current = real(current);
473 %--------------------------------------------------
474 % calculates current from given incomming mode
475 function current_trans = CurrentFromGivenMode( cell_index, idx )
477 modus_in = zeros(size(G,1), 1);
478 indexes = sum(M(1:cell_index-1)) + (1:M(cell_index));
479 modus_in(indexes) = Normamtx{cell_index}*modusmtx_p{cell_index}(:,idx)/sqrt(vcsop_p{cell_index}(idx));
481 modus_out = G*modus_in;
482 sc_wavefnc = modus_out(sum(M)+1:end);
484 coordinates = ribbon_loc.getCoordinates();
485 scatter_sites = ribbon_loc.CreateH.Read(
'kulso_szabfokok');
486 if ~isempty(scatter_sites) %
if the Green
operator is calculated from the full Hamiltonian, the coordinates contain all
sites.
487 BdG_indexes = coordinates.BdG_u(scatter_sites);
489 BdG_indexes = coordinates.BdG_u;
492 current_op_particle = zeros(size(current_op));
493 current_op_particle( BdG_indexes, BdG_indexes) = current_op( BdG_indexes, BdG_indexes)*fE;
494 current_op_particle( ~BdG_indexes, ~BdG_indexes) = current_op( ~BdG_indexes, ~BdG_indexes)*(1-fE);
497 %G_tmp = G(sum(M)+1:end, sum(M)+1:end);
498 % G_tmp( BdG_indexes, BdG_indexes) = G_tmp( BdG_indexes, BdG_indexes)*fE;
499 % G_tmp( ~BdG_indexes, ~BdG_indexes) = G_tmp( ~BdG_indexes, ~BdG_indexes)*(1-fE);
501 % current_trans_old = imag(trace( current_op*G_tmp ));
503 %G_tmp = G(sum(M)+1:end, sum(M)+1:end);
506 %current_trans = imag(trace( current_op_particle*G_tmp ));
510 current_trans = sc_wavefnc
'*current_op_particle*sc_wavefnc; 513 %transmitted quasiparticle current 514 quasi_current_trans = sc_wavefnc'*current_op*sc_wavefnc;
516 %H1sc = ribbon_loc.Scatter_UC.Read(
'H1');
517 %c1 = sc_wavefnc(1:M(cell_index));
518 %c2 = sc_wavefnc(M(cell_index)+1:2*M(cell_index));
519 %quasi_current_trans = 1/(1i)*( c2
'*H1sc'*c1 - c1
'*H1sc*c2) 521 % ******** unitarity test ********** 522 % coupling Hamiltonian to calculate the current in the incomming lead 523 H1 = Leads{cell_index}.Read('H1
'); 525 modus_in = modusmtx_p{cell_index}(:,idx)/sqrt(vcsop_p{cell_index}(idx)); 526 quasi_current_in = 1i*modus_in'*(H1*expk_p{cell_index}(idx)-H1
'/expk_p{cell_index}(idx))*modus_in; 528 modus_ref = modus_out(indexes) - modus_in; 529 r = d_modusmtx_m{cell_index}*modus_ref; 530 quasi_current_ref = 0; 531 for iiidx = 1:M(cell_index) 532 if abs(abs(expk_m{cell_index}(iiidx))-1) <= 1e-6 533 modus_ref_tmp = r(iiidx)*modusmtx_m{cell_index}(:,iiidx); 534 quasi_current_ref = quasi_current_ref + 1i*modus_ref_tmp'*(H1*expk_m{cell_index}(iiidx)-H1
'/expk_m{cell_index}(iiidx))*modus_ref_tmp; 537 quasi_current_ref = -quasi_current_ref; 539 delta_quasi_current = real(quasi_current_in - quasi_current_trans - quasi_current_ref); 541 if abs(delta_quasi_current) > tolerance 542 display(['EQuUs:Utils:
', class(obj), ':CurrentCalc_continuum: unitarity
test failed with \Delta I =
',num2str(delta_quasi_current)]); 544 % ******** unitarity test ********** 551 % end nested functions 554 %% createCurrentOperator 555 %> @brief Calculates the charge current operator from the inverse of the Greens function of the scattering region. 556 %> @param varargin Cell array of optional parameters (https://www.mathworks.com/help/matlab/ref/varargin.html): 557 %> @param 'junction
' An instance of class #NTerminal (or its subclass) describing the junction. 558 %> @return [1] The matrix operator of the current operator. 559 function current_op = createCurrentOperator( obj, varargin ) 562 p.addParameter('junction
', obj.junction); 563 p.parse(varargin{:}); 565 junction_loc = p.Results.junction; 567 obj.create_scatter_GreensFunction('junction
', junction_loc, 'gauge_trans
', true); 568 [~, gfininv ] = junction_loc.GetFiniteGreensFunction(); 574 if strcmpi( class(junction_loc), 'Ribbon' ) 577 current_op = zeros( size(gfininv) ); 578 % gfininv = E - Heff(E) !!!!! 579 current_op( 1:M, 2*M+1:3*M ) = -1/1i*(-1)*gfininv( 1:M, 2*M+1:3*M ); 580 current_op( M+1:2*M, 3*M+1:end ) = -1/1i*(-1)*gfininv( M+1:2*M, 3*M+1:end ); 581 current_op( 2*M+1:3*M, 1:M ) = 1/1i*(-1)*gfininv( 2*M+1:3*M, 1:M ); 582 current_op( 3*M+1:end, M+1:2*M ) = 1/1i*(-1)*gfininv(3*M+1:end, M+1:2*M ); 584 % gfininv = E - Heff(E) !!!!! 585 current_op = 1/1i*(triu(gfininv) - tril(gfininv)); 592 %> @brief Determines the band width of the leads and of the scattering region. 593 %> @return An instance of structure #BandWidth containing the bandwidths. 594 function ret = getBandWidth( obj ) 595 obj.display('Getting Bandwidths
', 1); 597 if ~isempty(obj.BandWidth) 601 poolmanager = Parallel( obj.Opt ); 602 poolmanager.openPool(); 604 Surface_1 = obj.junction.FL_handles.SurfaceGreenFunctionCalculator(1, 'Just_Create_Hamiltonians
', true, 'q
', obj.junction.getTransverseMomentum(), 'leadmodel
', obj.junction.leadmodel, ... 605 'CustomHamiltonian
', obj.junction.cCustom_Hamiltonians); 606 W = BandWidthCalculator( Surface_1 ); 608 param = obj.junction.FL_handles.Read( 'param' ); 609 pair_potential = min([abs(param.Leads{1}.pair_potential), abs(param.Leads{2}.pair_potential)]); 612 obj.BandWidth.Lead.Emin = abs(pair_potential); 613 obj.BandWidth.Lead.Emax = W; 616 if strcmpi( class(obj.junction), 'Ribbon' ) 617 obj.junction.CreateRibbon('justHamiltonians
', true); 618 W = BandWidthCalculator( obj.junction.Scatter_UC ); 619 elseif strcmpi( class(obj.junction), 'NTerminal' ) 620 obj.junction.CreateScatter(); 621 W = 1.3*BandWidthCalculatorMolecule( obj.junction.CreateH ); % multiplied by 1.3 just for sure to include normal bound states 626 obj.BandWidth.Scatter.Emin = 0; 627 obj.BandWidth.Scatter.Emax = W; 629 poolmanager.closePool(); 634 %-------------------------------------------- 635 function W = BandWidthCalculator( Scatter_UC ) 637 H0 = Scatter_UC.Read('H0
'); 638 if ~isempty( obj.junction.getTransverseMomentum() ) 639 H1_transverse = Scatter_UC.Read('H1_transverse
'); 640 H0 = H0 + H1_transverse*exp(1i*obj.junction.q) + H1_transverse'*exp(-1i*obj.junction.q);
642 H1 = Scatter_UC.Read(
'H1');
643 u = rand(size(H0,1));
646 hgetMaxEigenvalue = @getMaxEigenvalue;
649 arnoldi = zeros( length(k),1);
650 parfor (idx = 1:length(k), poolmanager.NumWorkers)
651 Heff = feval(hsecular_H,H0,H1,k(idx));
652 arnoldi(idx) = abs(feval(hgetMaxEigenvalue, Heff ));
657 % --------------------------------------
658 % Hamiltonian
for the secular equation
660 H = H0 + H1*exp(1i*k) + H1
'*exp(-1i*k); 663 %-------------------------------------------- 664 % getting the maximal eigenvalue with Arnoldi iteration 665 % Numerical Linear Algebra" by Trefethen and Bau. (p250) or 666 % http://en.wikipedia.org/wiki/Arnoldi_iteration 667 function ret = getMaxEigenvalue( Heff ) 668 %u = rand(size(H0,1)); 676 norm_new = norm(u_new); 677 lambda_new = norm_new/norm(u); 679 if abs((lambda - lambda_new)/lambda) < tolerance 690 %-------------------------------------------- 691 function W = BandWidthCalculatorMolecule( CreateH ) 692 Hscatter = CreateH.Read('Hscatter
'); 693 W = abs(eigs(Hscatter,1)); 696 % end nested functions 704 methods (Access=protected) 709 %% create_scatter_GreensFunction 710 %> @brief Calculates the surface Green operator of the scattering region. 711 %> The calculated Green operator is stored within the objec #junction (or in the object junction given as an optional parameter) 712 %> @param varargin Cell array of optional parameters (https://www.mathworks.com/help/matlab/ref/varargin.html): 713 %> @param 'gauge_trans
' Logical value. Set true (default) to perform gauge transformation on the Green's
function and on the
Hamiltonians, or
false otherwise.
714 %> @
param 'junction' An instance of
class #
NTerminal (or its subclass) describing the junction.
715 function create_scatter_GreensFunction( obj, varargin )
718 p.addParameter(
'gauge_trans',
true); % logical:
true if want to perform gauge transformation on the Green
's function and Hamiltonians 719 p.addParameter('junction
', obj.junction); %for parallel computations 720 p.parse(varargin{:}); 721 gauge_trans = p.Results.gauge_trans; 722 ribbon_loc = p.Results.junction; 724 if obj.gfininvfromHamiltonian 725 %calculating the Greens function from the Hamiltonian 726 ribbon_loc.CalcFiniteGreensFunctionFromHamiltonian('onlyGinv
', true, 'gauge_trans
', gauge_trans, 'PotInScatter
', obj.scatterPotential ); 728 % Calculating the Greens function by the fast way optimized for translational invariant systems 729 ribbon_loc.CalcFiniteGreensFunction('onlyGinv
', true, 'gauge_trans
', gauge_trans); 737 %% increasePhaseDifference 738 %> @brief Increase the phase difference between the Leads by delta_phi (apply gauge transformation on the second lead) 739 %> @param delta_phi The phase difference increment between the leads. 740 %> @param ribbon_loc An instance of class NTerminal (or its subclass) describing the junction. 741 function increasePhaseDifference( obj, delta_phi, junction_loc ) 744 if isfield(junction_loc, 'Ribbons
') 745 junction_loc = junction_loc.Ribbons{1}; %for the CombineRibbons interface 748 Leads = junction_loc.FL_handles.Read( 'Leads
' ); 753 gaugetrans( Leads{2} ) 756 % gauge transformation in the interface region 757 if isempty( junction_loc.Interface_Regions ) 760 gaugetrans( junction_loc.Interface_Regions{2} ) 765 %-------------------------------------- 766 % function that performs tha gauge transformation on the pair 768 function gaugetrans( Lead ) 769 if abs(delta_phi ) < 1e-16 773 coordinates = Lead.Read('coordinates
'); 774 regular_sites = Lead.Read('kulso_szabfokok
'); 776 if isempty(coordinates) 777 % Hamiltonians are not created yet 781 if isempty(coordinates.BdG_u) 782 % not a BdG calculation 783 warning(['EQuUs:Utils:
', class(obj), ':gaugetrans:Not a BdG calculation. Gauge transformation of the pair potential is skipped
']); 787 if isempty(regular_sites) 788 regular_sites = 1:Lead.Get_Neff(); 791 Umtx_diag = zeros(size(regular_sites)); 793 Umtx_diag( coordinates.BdG_u(regular_sites) ) = exp(1i*delta_phi/2); 794 Umtx_diag( ~coordinates.BdG_u(regular_sites) ) = exp(-1i*delta_phi/2); 795 Umtx = diag(Umtx_diag); 797 Lead.Unitary_Transform( Umtx ); 801 % end nested functions 807 %> @brief Resets the phase difference to the initial value. 808 %> @param varargin Cell array of optional parameters (https://www.mathworks.com/help/matlab/ref/varargin.html): 809 %> @param 'junction
' An instance of class #NTerminal (or its subclass) describing the junction. 810 function ResetPhase( obj, varargin ) 812 p.addParameter('junction
', obj.junction); 813 p.parse(varargin{:}); 815 junction_loc = p.Results.junction; 817 % retriving the param structure 818 param = junction_loc.getParam(); 820 % setting the phase of the second lead 821 junction_loc.setParam( param ); 830 %> @brief Parses the optional parameters for the class constructor. 831 %> @param varargin Cell array of optional parameters (https://www.mathworks.com/help/matlab/ref/varargin.html): 832 %> @param 'T
' The temperature in Kelvin (scalar or an array) 833 %> @param 'mu
' The Chemical potential in the same unit as other energy scales in the Hamiltonians. 834 %> @param 'scatterPotential
' A function handle pot=f( #coordinates ) or pot=f( #CreateHamiltonians, Energy) for the potential to be applied in the Hamiltonian (used when FiniteGreensFunctionFromHamiltonian=true). 835 %> @param 'gfininvfromHamiltonian
' Logical value. Set true calculate the surface Greens function of the scattering region from the Hamiltonaian of the scattering region, or false (default) to calculate it by the fast way (see Phys. Rev. B 90, 125428 (2014) for details). 836 %> @param 'useSelfEnergy
' Logical value. Set true to use the self energies of the leads in the Dyson equation, or false (default) to use the surface Green function instead. 837 %> @param 'junction
' An instance of class #NTerminal (or its derived class) representing the junction. 838 function InputParsing( obj, varargin ) 841 p.addParameter('T
', obj.T); 842 p.addParameter('mu
', obj.mu, @isscalar); 843 p.addParameter('scatterPotential
', obj.scatterPotential); 844 p.addParameter('gfininvfromHamiltonian
', obj.gfininvfromHamiltonian); 845 p.addParameter('ribbon
', obj.junction); % NEW OBSOLETE use parameter junction istead 846 p.addParameter('junction
', obj.junction); % NEW 847 p.addParameter('useSelfEnergy
', obj.useSelfEnergy); %NEW 849 p.parse(obj.varargin{:}); 851 junction = p.Results.junction; 852 if ~isempty(p.Results.ribbon) && isempty(junction) 853 junction = p.Results.ribbon; 857 InputParsing@UtilsBase( obj, 'T
', p.Results.T, ... 858 'mu
', p.Results.mu, ... 859 'scatterPotential
', p.Results.scatterPotential, ... 860 'gfininvfromHamiltonian
', p.Results.gfininvfromHamiltonian, ... 861 'useSelfEnergy
', p.Results.useSelfEnergy, ... 862 'junction
', junction); 867 obj.junction.FL_handles.LeadCalc('createCore
', true, 'q', obj.junction.getTransverseMomentum() ); 873 end % protected methods A class describing an N-terminal geometry for equilibrium calculations mostly in the zero temperature...
function test(arg1, arg2)
Brief description of the function.
lead_param Leads
A list of structures lead_param containing the physical parameters for the scattering region...
Property T
The temperature in Kelvin.
A class for controlling the parallel pool for paralell computations.
Structure Opt contains the basic computational parameters used in EQuUs.
A class for calculations on a ribbon of finite width for equilibrium calculations mostly in the zero ...
pair_potential
A physical parameter, see the individual lattice documentations for details.
function Transport(Energy, B)
creating the Ribbon class representing the twoterminal setup
function Hamiltonians(varargin)
Function to create the custom Hamiltonians for the 1D chain.
A class to calculate the DC Josephson current.
function secular_H(H0, H1, k)
Hamiltonian for the secular equation.
Property q
The transverse momentum quantum number.
Structure BandWidthLimits contains the bandwidth limits.
Structure param contains data structures describing the physical parameters of the scattering center ...
Structure sites contains data to identify the individual sites in a matrix.
function GetFiniteGreensFunction()
Query fo the calculated Green operator of the scattering center.
This class is a base class containing common properties and methods utilized in several other classes...
An object for combining multiple ribbon parts of equal width and from the same material in a two term...
function openPool()
Opens a parallel pool for parallel computations.