2 % Copyright (C) 2016 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 basic Basic Functionalities
20 %> @brief Class to solve the eigenproblem of a translational invariant leads and calculate the group velocities.
22 %> @brief Class to solve the eigenproblem of a translational invariant leads and calculate the group velocities.
24 %> EQuUs v4.8 or later
29 properties ( Access =
protected )
30 %>
true for calculating the retarded Green
function, or
false for the advanced Green
function.
32 %> The unsorted right-sided eigenstates.
34 %> The unsorted left-sided eigenstates.
36 %> The unsorted wave numbers in form exp(1i*k).
38 %> The unsorted group velocities.
40 %> The wave numbers of the eigenstates, that propagates or decays in the positive direction.
42 %> The wave numbers of the eigenstates in form exp(1i*k), that propagates or decays in the negative direction.
44 %> The group velocities of the eigenstates in form exp(1i*k), that propagates or decays in the positive direction.
46 %> The group velocities of the eigenstates, that propagates or decays in the negative direction.
48 %> The right sided wave functions of the eigenstates, that propagates or decays in the positive direction.
50 %> The right-sided wave functions of the eigenstates, that propagates or decays in the negative direction.
52 %> The left-sided wave numbers of the eigenstates, that propagates or decays in the positive direction.
54 %> The left-sided wave functions of the eigenstates, that propagates or decays in the negative direction.
56 %> The dual basis of the right-sided wave functions of the eigenstates, that propagates or decays in the positive direction.
58 %> The dual basis of the right-sided wave functions of the eigenstates, that propagates or decays in the negative direction.
60 %> A real number corresponding to the tolerance used to sort the left and right moving (decaying) modes.
62 %> A structure #
open_channels containing info on the open channels
64 %> Logical matrix containing the degenerate k subspaces
65 degenerate_k_subspaces
66 %> The energy value
for which the #TrukkosSajatertekek eigenvalue problem was solved.
73 methods ( Access =
public )
74 %% constructorof the
class 75 %> @brief Constructor of the
class.
79 %> @
return An instance of the
class 84 %% Determine_Open_Channels
85 %> @brief Determines the open channels in the lead. The data are storen within the attribute #
open_channels.
86 function Determine_Open_Channels( obj )
87 if isempty(obj.expk_p)
88 obj.display('Right and left going states were not determined yet. Use method
EigenProblemLead::szetvalaszto first.', 1);
91 propagating_indexes_p = logical( abs( abs(obj.expk_p) - 1 ) < obj.sort_tolerance );
92 propagating_indexes_m = logical( abs( abs(obj.expk_m) - 1 ) < obj.sort_tolerance );
95 obj.open_channels.num = length(find(propagating_indexes_p));
96 obj.open_channels.open_channels_p = propagating_indexes_p;
97 obj.open_channels.open_channels_m = propagating_indexes_m;
100 if obj.isSuperconducting()
101 obj.open_channels.isSuperconducting = true;
103 obj.open_channels.BdG_u_incoming = logical(diag(obj.modusmtx_p(obj.coordinates.BdG_u,:)'*obj.modusmtx_p(obj.coordinates.BdG_u,:)));
104 obj.open_channels.BdG_u_outgoing = logical(diag(obj.modusmtx_m(obj.coordinates.BdG_u,:)'*obj.modusmtx_m(obj.coordinates.BdG_u,:)));
105 obj.open_channels.isSuperconducting = false;
108 obj.open_channels.isSuperconducting = false;
114 %> @brief Sorts the left and right propagating (decaying) modes.
115 %> @
param tolerance The tolerance to be used during the sort.
116 function szetvalaszto( obj, tolerance )
118 if isempty(obj.csoportseb)
123 if ~exist('tolerance', 'var')
125 elseif tolerance < 1e-16
126 err = MException(['EQuUs:', class(obj), ':szetvalaszto'], 'Failed to separate the left and right-going states');
127 save(['EQuUs_', class(obj), '_szetvalaszto.mat']);
132 indexes = (real(obj.csoportseb*obj.Lead_Orientation) > 0 & abs(abs(obj.expk)-1) < tolerance) | (1-abs(obj.expk))*obj.Lead_Orientation > tolerance;
134 indexes = (real(obj.csoportseb*obj.Lead_Orientation) < 0 & abs(abs(obj.expk)-1) < tolerance) | (1-abs(obj.expk))*obj.Lead_Orientation > tolerance;
136 obj.expk_p = obj.expk(indexes);
137 obj.vcsop_p = obj.csoportseb(indexes);
138 obj.modusmtx_p = obj.modusmtx(:,indexes);
139 if ~obj.
Opt.usingDualModes
140 obj.modusmtx_p_left = obj.modusmtx_left(indexes,:);
143 obj.expk_m = obj.expk(~indexes);
144 obj.vcsop_m = obj.csoportseb(~indexes);
145 obj.modusmtx_m = obj.modusmtx(:,~indexes);
146 if ~obj.
Opt.usingDualModes
147 obj.modusmtx_m_left = obj.modusmtx_left(~indexes,:);
150 obj.sort_tolerance = tolerance;
151 if length(obj.expk_p) ~= length(obj.expk_m)
152 obj.szetvalaszto( obj.sort_tolerance/10 );
164 %> @brief Calculates the group velocities corresponding to the propagating states. The calculated group velocities are stored within the class.
165 function Group_Velocity(obj)
167 usingDualModes = obj.
Opt.usingDualModes;
169 csoportseb = obj.csoportseb;
171 needs_to_be_diagonalized = false;
175 if needs_to_be_diagonalized
176 Diagonalize_Current_Operator();
180 obj.csoportseb = diag(csoportseb);
183 %------------------------------------------------------------------
185 function Diagonalize_Current_Operator()
187 [U, eigvals] = eig( full(csoportseb));
189 % rotate the wavefunctions
190 obj.modusmtx = obj.modusmtx*U;
192 obj.modusmtx_left = U'*obj.modusmtx_left;
196 % reordering the wavenumbers
197 obj.expk = diag(U'*diag(obj.expk)*U);
199 obj.degenerate_k_subspaces = logical(speye(length(obj.expk)));
202 %------------------------------------------------------------------
204 %------------------------------------------------------------------
205 function GroupVelocityCalc()
207 [~, K1_eff, K1adj_eff] = obj.Get_Effective_Hamiltonians();
208 Neff = obj.Get_Neff();
209 csoportseb = sparse([],[],[], 2*Neff, 2*Neff);
211 adj_modusmtx = obj.modusmtx';
213 jdx_indexes = obj.degenerate_k_subspaces(idx,:);
214 if ~needs_to_be_diagonalized && length(csoportseb(idx, jdx_indexes)) >1
215 needs_to_be_diagonalized = true;
217 %Eq (18) in GNew J. Phys 093029 (2014), (modes are already normalized with the overlap integrals)
218 csoportseb(idx, jdx_indexes) = 1i*adj_modusmtx(idx,:)*(K1_eff*obj.expk(idx)-K1adj_eff*(1/obj.expk(idx)) )*obj.modusmtx(:,jdx_indexes);
222 %equation (A4) in PRB 78, 035407 (difference in the evanescent group velocities)
224 jdx_indexes = obj.degenerate_k_subspaces(idx,:);
225 if ~needs_to_be_diagonalized && length(csoportseb(idx, jdx_indexes)) >1
226 needs_to_be_diagonalized = true;
228 csoportseb(idx, jdx_indexes) = 1i*obj.modusmtx_left(idx,:)*(K1_eff*obj.expk(idx)-K1adj_eff*(1/obj.expk(idx)) )*obj.modusmtx(:,jdx_indexes);
232 % end nested functions
233 %------------------------------------------------------------------
240 %% TrukkosSajatertekek
241 %> @brief Calculates the wave numbers corresponding to the propagating states at given energy. The calculated wave numbers are stored by attribute
#expk. 242 %> @
param E The energy value used in the calculations.
243 function TrukkosSajatertekek(obj, E)
246 obj.retarted =
false;
253 % create effective
Hamiltonians with SVD regularization
if necessary
254 obj.Calc_Effective_Hamiltonians( E );
255 [K0_eff, K1_eff, K1adj_eff] = obj.Get_Effective_Hamiltonians();
256 Neff = obj.Get_Neff();
258 usingDualModes = obj.Opt.usingDualModes;
260 if obj.isSuperconducting() || ~obj.Opt.BdG
261 %
if normal system, or
if BdG Hamiltonian with nonzero pair potential
262 [obj.modusmtx,obj.expk, obj.modusmtx_left] = calcTrickyEigenvalues( K0_eff, K1_eff, K1adj_eff );
264 %
if BdG Hamiltonian with nonzero pair potential
266 K0_e = K0_eff(1:M0, 1:M0);
267 K0_h = K0_eff(M0+1:Neff, M0+1:Neff);
268 K1_e = K1_eff(1:M0, 1:M0);
269 K1_h = K1_eff(M0+1:Neff, M0+1:Neff);
270 K1adj_e = K1adj_eff(1:M0, 1:M0);
271 K1adj_h = K1adj_eff(M0+1:Neff, M0+1:Neff);
273 % eigenvalues
for electron-like part
274 [modusmtx_e, expk_e, modusmtx_left_e] = calcTrickyEigenvalues( K0_e, K1_e, K1adj_e );
276 %eigenvalues
for hole-like part
277 [modusmtx_h, expk_h, modusmtx_left_h] = calcTrickyEigenvalues( K0_h, K1_h, K1adj_h );
279 obj.modusmtx = [ modusmtx_e, zeros( size(modusmtx_e,1), size(modusmtx_e,2)); ...
280 zeros(size(modusmtx_h,1), size(modusmtx_e,2)), modusmtx_h ];
282 obj.expk = [expk_e;expk_h];
285 obj.modusmtx_left = [ modusmtx_left_e, zeros( size(modusmtx_left_e,1), size(modusmtx_left_e,2)); ...
286 zeros(size(modusmtx_left_h,1), size(modusmtx_left_e,2)), modusmtx_left_h ];
291 obj.NormalizeModes();
292 obj.Determine_Degenerate_k_subspaces();
296 obj.d_modusmtx_p = [];
297 obj.d_modusmtx_m = [];
303 %------------------------------------------------------------------
305 % calculates the tricky eigenvalue problem
306 function [modusmtx,expk, modusmtx_left] = calcTrickyEigenvalues( H0, H1, H1adj )
308 [H_eff,B_eff] = Effective_H(H0,H1,H1adj);
311 [modusmtx,expk] = eig(H_eff,B_eff);
315 [modusmtx, modusmtx_left, expk] = obj.eig(H_eff,B_eff);
319 % Eq (A3) in PRB 78, 035407
320 %normfactors = diag(modusmtx_left*B_eff*modusmtx);
321 %modusmtx_left = diag(1./normfactors)*modusmtx_left;
322 modusmtx_left = modusmtx_left(:,1:M);
323 modusmtx_left = diag(sqrt(expk))*modusmtx_left;
327 % Eq (7) in PRB 78, 035407
328 modusmtx = modusmtx(1:M,:);
329 modusmtx = modusmtx*diag(1./sqrt(expk));
331 %------------------------------------------------------------------
334 %------------------------------------------------------------------
335 % Eq (6) PRB 78 035407
336 function [H_eff,B_eff] = Effective_H(H0,H1,H1adj)
338 %H1inv = H1\eye(size(H1));%inv(H1);
339 H_eff = full([-H0,full(-H1adj);eye(M),zeros(M)]);
340 B_eff = full([ H1, zeros(M); zeros(M), eye(M)]);
343 %------------------------------------------------------------------
344 % end nested functions
349 %> @brief Adds on-site potential to the Hamiltonian H0.
350 %> @
param V The potential calculated on the
sites.
351 function AddPotential( obj, V )
357 %> @brief Transforms the effective
Hamiltonians by a unitary transformation
358 %> @
param Umtx The matrix of the unitary transformation.
359 function Unitary_Transform(obj, Umtx)
363 if ~isempty(obj.modusmtx_p)
364 obj.modusmtx_p = Umtx*obj.modusmtx_p;
367 if ~isempty(obj.modusmtx_m)
368 obj.modusmtx_m = Umtx*obj.modusmtx_m;
371 if ~isempty(obj.modusmtx_p_left)
372 obj.modusmtx_p_left = Umtx*obj.modusmtx_p_left;
375 if ~isempty(obj.modusmtx_m_left)
376 obj.modusmtx_m_left = Umtx*obj.modusmtx_m_left;
379 if ~isempty(obj.d_modusmtx_p)
380 obj.d_modusmtx_p = obj.modusmtx_p*Umtx';
383 if ~isempty(obj.d_modusmtx_m)
384 obj.d_modusmtx_m = obj.modusmtx_m*Umtx';
391 %> @brief Creates a clone of the present
object.
392 %> @return Returns with the cloned class.
393 %> @
param varargin Cell array of optional parameters:
394 %> @
param 'empty' Set true to create an empty class, or false (default) to copy all the attributes.
395 %> @return Returns with an instance of the cloned class.
396 function ret = CreateClone( obj )
399 p.addParameter('empty', false );
401 p.parse(varargin{:});
402 empty = p.Results.empty;
405 'Hanyadik_Lead', obj.Hanyadik_Lead,...
406 'Lead_Orientation', obj.Lead_Orientation, ...
412 meta_data = metaclass(obj);
414 for idx = 1:length(meta_data.PropertyList)
415 prop_name = meta_data.PropertyList(idx).Name;
416 if strcmpi( prop_name, 'next_SVD_cycle' )
417 if ~isempty( obj.next_SVD_cycle )
418 ret.next_SVD_cycle = obj.next_SVD_cycle.CreateClone();
421 ret.Write( prop_name, obj.(prop_name));
428 %> @brief Resets all elements in the
object.
429 function Reset( obj )
432 meta_data = metaclass(obj);
434 for idx = 1:length(meta_data.PropertyList)
435 prop_name = meta_data.PropertyList(idx).Name;
436 if strcmp(prop_name, '
Opt') || strcmp( prop_name, '
param') || strcmp(prop_name, 'varargin')
439 obj.Clear( prop_name );
448 %> @brief Sets the value of an attribute in the class.
449 %> @
param MemberName The name of the attribute to be set.
450 %> @
param input The value to be set.
451 function Write(obj, MemberName, input)
453 if strcmp(MemberName, '
param') || strcmp(MemberName, 'params')
458 obj.(MemberName) = input;
460 error(['EQuUs:', class(obj), ':Read'], ['No property to write with name: ', MemberName]);
466 %> @brief Query for the value of an attribute in the class.
467 %> @
param MemberName The name of the attribute to be set.
468 %> @return Returns with the value of the attribute.
469 function ret = Read(obj, MemberName)
471 if strcmp(MemberName, 'Hamiltonian2Dec')
476 ret = obj.(MemberName);
478 error(['EQuUs:', class(obj), ':Read'], ['No property to read with name: ', MemberName]);
484 %> @brief Clears the value of an attribute in the class.
485 %> @
param MemberName The name of the attribute to be cleared.
486 function Clear(obj, MemberName)
489 obj.(MemberName) = [];
491 error(['EQuUs:', class(obj), ':Clear'], ['No property to clear with name: ', MemberName]);
500 methods (Access=private)
502 %% Determine_Degenerate_k_subspaces
503 %> @brief Determine the degenerate k subspaces for which the current operator needs to be diagonalized.
504 function Determine_Degenerate_k_subspaces(obj)
506 degeneracy_tolerance = 1e-12;
508 obj.degenerate_k_subspaces = obj.expk*transpose(1./obj.expk);
509 obj.degenerate_k_subspaces = abs(obj.degenerate_k_subspaces - 1) < degeneracy_tolerance;
515 %> @brief Renormalizes the wave functions with the overlap integrals (Eq (9) in PRB 78 035407
516 function NormalizeModes(obj)
517 usingDualModes = obj.
Opt.usingDualModes;
519 modusmtx_adj = obj.modusmtx';
520 [S0, S1, S1adj] = obj.Get_Effective_Overlaps();
522 for idx = 1:length(obj.expk)
523 expk = obj.expk(idx);
524 if ~isempty(obj.q) && ~obj.HamiltoniansDecimated
525 S0 = S0 + obj.S1_transverse*diag(exp(-1i*obj.q)) + obj.S1_transverse'*diag(exp(1i*obj.q));
529 S = S0 + S1*expk + S1adj*(1/expk);
531 l_n = modusmtx_adj(idx,:)*S*obj.modusmtx(:,idx);
532 obj.modusmtx(:,idx) = obj.modusmtx(:,idx)/sqrt(l_n);
534 % normalize the left sided eigenvectors
536 l_n = obj.modusmtx_left(idx,:)*S*obj.modusmtx(:,idx);
537 obj.modusmtx_left(idx,:) = obj.modusmtx_left(idx,:)/(l_n);
541 %no overlap integrals are present
542 normfactors = diag( modusmtx_adj*obj.modusmtx );
543 obj.modusmtx = obj.modusmtx*diag(1./sqrt(normfactors));
545 normfactors = diag( obj.modusmtx_left*obj.modusmtx );
546 obj.modusmtx_left = diag(1./(normfactors))*obj.modusmtx_left;
552 %% Extend_Wavefncs --- DEVELOPMENT
553 %> @brief Extend the reduced wave functions to the original basis before the SVD regularization (Eq (45) in PRB 78 035407) This function is in a development stage.
554 function Extend_Wavefncs(obj) %NEW
555 wavefnc_extended = zeros(obj.M, obj.Get_Neff());
556 for idx = 1:length(obj.expk)
557 wavefnc_extended(:,idx) = obj.Extend_Wavefnc( obj.modusmtx(:,idx), obj.expk(idx) );
560 obj.modusmtx = wavefnc_extended;
566 end % methods private
569 methods (Access=protected)
572 %> @brief Initializes
object properties.
573 function Initialize(obj)
575 obj.sort_tolerance = 1e-6;
Structure hole contains the data about the antidot used in class antidot.
Structure Opt contains the basic computational parameters used in EQuUs.
Structure open_channels describes the open channel in the individual leads.
Class to create and store Hamiltonian of the translational invariant leads.
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.
Class to regulerize the H1 problem of the leads by SVD decompozition.
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 SVDregularizationLead(Opt, param, varargin)
Constructor of the class.
Class to solve the eigenproblem of a translational invariant leads and calculate the group velocities...
function structures(name)