2 % Copyright (C) 2009-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 A
class describing the interface region between the scattering region and a lead.
22 %> @brief A
class describing the interface region between the scattering region and a lead.
24 %> EQuUs v4.8 or later
29 properties ( Access =
protected )
30 %> Effective Hamiltonian of the
interface region
32 %> Hinterface -E*Sinterface
34 %> Hamiltonian of the
interface region
36 %> Hcoupling-E*Scoupling
38 %> Hcouplingadj-E*Scoupling
' 40 %> Coupling Hamiltonian from the interface region to the scattering region 42 %> Coupling Hamiltonian from the scattering region to the interface region 44 %> Overlap integrals corresponding to Hcoupling 46 %> Overlap integrals corresponding to Hinterface 48 %> Coordinates of the surface sites of the scattering region to be connected to. 56 methods ( Access = public ) 57 %% constructorof the class 58 %> @brief Constructor of the class. 59 %> @param Opt An instance of the structure #Opt. 60 %> @param param An instance of structure #param. 61 %> @param varargin Cell array of optional parameters identical to #SVDregularizationLead.SVDregularizationLead. 62 %> @return An instance of the class 63 function obj = InterfaceRegion(Opt, param, varargin) 64 obj = obj@SVDregularizationLead( Opt, param, varargin{:} ); 72 %% ApplyOverlapMatrices 73 %> @brief Applies the overlap matrices to the Hamiltonians: K = H-ES 74 %> @param E The energy value. 75 function ApplyOverlapMatrices(obj, E) 78 obj.display('EQuUs:
InterfaceRegion:ApplyOverlapMatrices: Overlap matrices were already applied.
'); 82 ApplyOverlapMatrices@CreateLeadHamiltonians( obj, E ); 84 if ~isempty( obj.Scoupling ) 85 obj.Kcoupling = obj.Hcoupling - E*obj.Scoupling; 86 obj.Kcouplingadj = obj.Hcouplingadj - E*obj.Scoupling';
89 obj.Kcouplingadj = obj.Hcouplingadj;
92 if ~isempty( obj.Sinterface )
93 obj.Kinterface = obj.Hinterface - E*obj.Sinterface;
95 obj.Kinterface = obj.Hinterface;
102 %> @brief Transforms the
Hamiltonians and the overlap matrices into the BdG model.
103 function Transform2BdG( obj )
105 if isempty(obj.
Opt.BdG) || ~obj.
Opt.BdG
109 if ~isempty( obj.coordinates.BdG_u )
110 % already transformed into BdG model
116 if isempty( obj.coordinates2 )
117 fnames = fieldnames( obj.coordinates2 );
118 for idx = 1:length(fnames)
120 if strcmp(fname,
'a') || strcmp(fname,
'b')
123 obj.coordinates2.(fname) = [ obj.coordinates2.(fname); obj.coordinates2.(fname) ];
125 obj.coordinates2.BdG_u = [
true(size(obj.H0,1),1);
false(size(obj.H0,1),1)];
129 pair_potential = obj.
params.pair_potential;
132 if isempty(obj.Scoupling)
133 Scoupling = sparse([],[],[], size(obj.Hcoupling,1), size(obj.Hcoupling,2));
135 Scoupling = obj.Scoupling;
138 obj.Hcoupling = [obj.Hcoupling, Scoupling*pair_potential; Scoupling*conj(pair_potential), -conj(obj.Hcoupling)];
139 obj.Hcouplingadj = [obj.Hcouplingadj, Scoupling
'*pair_potential; Scoupling'*conj(pair_potential), -conj(obj.Hcouplingadj)];
141 if ~isempty( obj.Scoupling )
142 obj.Scoupling = [obj.Scoupling, sparse([], [], [], size(obj.Scoupling,1), size(obj.Scoupling,2)); sparse([], [], [], size(obj.Scoupling,1), size(obj.Scoupling,2)), obj.Scoupling];
148 %% Calc_Effective_Hamiltonians
149 %> @brief Calculates the effective
Hamiltonians according to Eq (48) of of PRB 78, 035407
150 %> @
param E The energy value
151 %> @
param varargin Cell array of optional parameters (https:
153 function Calc_Effective_Hamiltonians( obj, E, varargin )
156 p.addParameter(
'Lead', []);
157 p.parse(varargin{:});
160 obj.ApplyOverlapMatrices(E);
163 if ~isempty(obj.kulso_szabfokok) && length(obj.kulso_szabfokok) ~= size(obj.K0,1)
164 obj.Decimate_Hamiltonians();
166 elseif
Lead.
Read(
'is_SVD_transformed')
167 obj.SVD_transform(
Lead );
170 % in
case no SVD regularization or simple decimation is needed
171 obj.Neff = size(obj.K0,1);
182 %% Get_Effective_Hamiltonians
183 %> @brief Gets the effective
Hamiltonians of the
interface region
184 %> @
return [1] The effective Hamiltonian of the unit cell
185 %> @
return [2] The effective coupling between the unit cells
186 %> @
return [3] The adjungate of the effective coupling between the unit cells
187 %> @
return [4] The effective coupling between the lead and the scattering region
188 %> @
return [5] The adjungate of the effective coupling between the lead and the scattering region
189 function [K0_eff, K1_eff, K1adj_eff, Kcoupling, Kcouplingadj] = Get_Effective_Hamiltonians(obj)
190 if isempty(obj.next_SVD_cycle)
192 if obj.Lead_Orientation == 1
194 K1adj_eff = obj.K1adj;
195 elseif obj.Lead_Orientation == -1
199 Kcoupling = obj.Kcoupling;
200 Kcouplingadj = obj.Kcouplingadj;
201 elseif strcmpi(
class(obj.next_SVD_cycle),
'InterfaceRegion' )
202 [K0_eff, K1_eff, K1adj_eff, Kcoupling, Kcouplingadj] = obj.next_SVD_cycle.Get_Effective_Hamiltonians();
204 error(
'EQuUs:InterfaceRegion:Get_Effective_Hamiltonians',
'Unrecognized type of attribute next_SVD_cycle');
211 %> @brief Transforms the effective
Hamiltonians by a unitary transformation
212 %> @
param Umtx The matrix of the unitary transformation.
213 function Unitary_Transform(obj, Umtx)
215 if isempty(obj.next_SVD_cycle)
216 obj.K00 = Umtx*obj.K00*Umtx
'; 217 obj.K1 = Umtx*obj.K1*Umtx';
218 obj.K1adj = Umtx*obj.K1adj*Umtx
'; 219 %obj.Kcoupling = Umtx*obj.Kcoupling; 220 %obj.Kcouplingadj = obj.Kcouplingadj*Umtx';
222 elseif strcmpi(
class(obj.next_SVD_cycle),
'InterfaceRegion' )
223 obj.next_SVD_cycle.Unitary_Transform(Umtx);
225 error(
'EQuUs:InterfaceRegion:Unitary_Transform',
'Unrecognized type of attribute next_SVD_cycle');
231 %> @brief Transforms the
Hamiltonians of the
interface region according to the SVD regularization in the lead.
233 function SVD_transform( obj, Lead )
236 obj.Neff = Lead.
Read(
'Neff');
238 if ~Lead.
Read(
'is_SVD_transformed')
243 obj.
V = Lead.
Read(
'V');
246 obj.
U = Lead.
Read(
'U');
251 %> Performing SVD regularization
257 K0 = Lead.Read('K0'); 260 % Eqs (48) in PRB 78, 035407 261 % here we use diiferent transformation for the leads with diffferent orientation 266 K1(1:Neff,1:Neff) =
U'*obj.K1(1:Neff,1:Neff)*U; 267 K1adj(1:Neff,1:Neff) = U'*obj.
K1adj(1:Neff,1:Neff)*
U;
273 K00(1:Neff, 1:Neff) =
U'*obj.K00(1:Neff, 1:Neff)*U; 274 K00(1:Neff, end-Neff+1:end) = U'*obj.K00(1:Neff, end-Neff+1:end);
275 K00(end-Neff+1:end, 1:Neff) = obj.K00(end-Neff+1:end, 1:Neff)*
U;
278 S = abs(Lead.
Read(
'S'));
279 non_singular_sites_slab = transpose(find(
S > obj.
tolerance*max(
S)));
283 non_singular_sites = [non_singular_sites_slab, size(
K0,1)+1:size(
K0,1)+size(K00,1)];
286 non_singular_sites = [1:size(K00,1), size(K00,1) + non_singular_sites_slab];
293 CreateH.Write(
'Hscatter', K);
294 CreateH.Write(
'kulso_szabfokok', non_singular_sites);
295 CreateH.Write(
'HamiltoniansCreated',
true);
296 CreateH.Write(
'HamiltoniansDecimated',
false);
299 Decimation_Procedure.DecimationFunc(0, CreateH,
'Hscatter',
'kulso_szabfokok');
301 K = CreateH.Read(
'Hscatter');
302 CreateH.Clear(
'Hscatter');
304 Neff_new = length(non_singular_sites_slab);
307 K00 = K(Neff_new+1:end, Neff_new+1:end);
308 K1 = K(1:Neff_new, Neff_new+1:end);
309 K1adj = K(Neff_new+1:end, 1:Neff_new);
311 K00 = K(1:end-Neff_new, 1:end-Neff_new);
312 K1 = K(1:end-Neff_new, end-Neff_new+1:end);
313 K1adj = K(end-Neff_new+1:end, 1:end-Neff_new);
320 Kcoupling = obj.Kcoupling;
321 Kcouplingadj = obj.Kcouplingadj;
323 Kcoupling(1:Neff,:) =
U'*Kcoupling(1:Neff,:); 324 Kcouplingadj(:,1:Neff) = Kcouplingadj(:,1:Neff)*U; 326 Kcoupling(1:Neff,:) = U'*Kcoupling(1:Neff,:);
327 Kcouplingadj(:,1:Neff) = Kcouplingadj(:,1:Neff)*
U;
346 error(
'EQuUs:InterfaceRegion:Calc_Effective_Hamiltonians',
'Unrecognized type of attribute next_SVD_cycle');
351 %% Decimate_Interface
352 %> @brief Decimates the
Hamiltonians of the
interface region.
353 function Decimate_Hamiltonians( obj )
359 K = [K0loc,obj.K1;obj.K1adj,K0loc]; 361 non_singular_sites_slab = obj.kulso_szabfokok; 362 Neff_new = length(non_singular_sites_slab); 365 non_singular_sites = [non_singular_sites_slab, size(K,1)/2+1:size(K,1)]; 367 CreateH = CreateHamiltonians(obj.Opt, obj.param); 368 CreateH.Write('Hscatter
', K); 369 CreateH.Write('kulso_szabfokok
', non_singular_sites); 370 CreateH.Write('HamiltoniansCreated
', true); 371 CreateH.Write('HamiltoniansDecimated
', false); 373 Decimation_Procedure = Decimation(obj.Opt); 374 Decimation_Procedure.DecimationFunc(0, CreateH, 'Hscatter
', 'kulso_szabfokok
'); 376 K = CreateH.Read('Hscatter
'); 377 CreateH.Clear('Hscatter
'); 379 cCoordinates = obj.coordinates; 380 coord_fieldnames = fieldnames(cCoordinates); 381 if obj.Lead_Orientation == 1 382 K00 = K(Neff_new+1:end, Neff_new+1:end); 383 K1 = K(1:Neff_new, Neff_new+1:end); 384 K1adj = K(Neff_new+1:end, 1:Neff_new); 386 % reordering the sites to get the non-singular sites at the top corner 387 indexes = false( size(K00,1), 1); 388 indexes(non_singular_sites_slab) = true; 389 K00 = [K00(indexes, indexes), K00(indexes, ~indexes); K00(~indexes, indexes), K00(~indexes, ~indexes)]; 390 K1 = [K1(:, indexes), K1(:, ~indexes)]; 391 K1adj = [K1adj(indexes, :); K1adj(~indexes, :)]; 392 Kcoupling = obj.Kcoupling; 393 Kcoupling = [Kcoupling(indexes, :); Kcoupling(~indexes, :)]; 394 Kcouplingadj = obj.Kcouplingadj; 395 Kcouplingadj = [Kcouplingadj(:,indexes), Kcouplingadj(:, ~indexes)]; 397 % reordering the sites to get the non-singular sites at the top corner 398 for idx = 1:length(coord_fieldnames) 399 fieldname = coord_fieldnames{idx}; 400 if strcmpi(fieldname, 'a
') || strcmpi(fieldname, 'b
') || strcmpi(fieldname, 'LatticeConstant
') 402 elseif ~isempty(cCoordinates.(fieldname)) 403 tmp = cCoordinates.(fieldname); 404 cCoordinates.(fieldname) = [tmp(indexes); tmp(~indexes)]; 408 elseif obj.Lead_Orientation == -1 409 % keeping only regular sites 410 K00 = K(1:Neff_new, 1:Neff_new); 411 K1 = K(1:Neff_new, Neff_new+non_singular_sites_slab); 412 K1adj = K(Neff_new+non_singular_sites_slab, 1:Neff_new); 413 Kcoupling = obj.Kcoupling; 414 Kcouplingadj = obj.Kcouplingadj; 416 % keeping only regular sites 417 cCoordinates = cCoordinates.KeepSites(non_singular_sites_slab); 419 error('EQuUs:
InterfaceRegion:Decimate_Interface
', 'Unknown lead orientation
'); 422 obj.next_SVD_cycle = InterfaceRegion(obj.Opt, obj.param, obj.varargin{:}); 423 obj.next_SVD_cycle.Write('K00
', K00); 424 obj.next_SVD_cycle.Write('K1
', K1); 425 obj.next_SVD_cycle.Write('K1adj
', K1adj); 426 obj.next_SVD_cycle.Write('Kcoupling
', Kcoupling); 427 obj.next_SVD_cycle.Write('Kcouplingadj
', Kcouplingadj); 428 obj.next_SVD_cycle.Write('HamiltoniansDecimated
', true); 429 obj.next_SVD_cycle.Write('Neff
', Neff_new); 430 obj.next_SVD_cycle.Write('coordinates
', cCoordinates); 433 obj.HamiltoniansDecimated = true; 434 obj.Neff = size(obj.K0,1); 440 %> @brief Gets the effective number of sites after the elimination of the singular values. 441 %> @return Returns with the effective number of sites 442 function Neff = Get_Neff(obj) 443 if isempty(obj.next_SVD_cycle) 446 Neff = obj.next_SVD_cycle.Get_Neff(); 448 error('EQuUs:
InterfaceRegion:Get_Neff
', 'Unrecognized type of attribute next_SVD_cycle
'); 455 %> @brief Creates a clone of the present class. 456 %> @return Returns with the cloned object. 457 function ret = CreateClone( obj ) 460 ret = InterfaceRegion(obj.Opt, obj.param,... 461 'Hanyadik_Lead
', obj.Hanyadik_Lead,... 462 'Lead_Orientation
', obj.Lead_Orientation, ... 465 meta_data = metaclass(obj); 467 for idx = 1:length(meta_data.PropertyList) 468 prop_name = meta_data.PropertyList(idx).Name; 469 if strcmpi( prop_name, 'next_SVD_cycle
' ) 470 if ~isempty( obj.next_SVD_cycle ) 471 ret.next_SVD_cycle = obj.next_SVD_cycle.CreateClone(); 474 ret.Write( prop_name, obj.(prop_name)); 481 %> @brief Resets all elements in the object. 482 function Reset( obj ) 485 meta_data = metaclass(obj); 487 for idx = 1:length(meta_data.PropertyList) 488 prop_name = meta_data.PropertyList(idx).Name; 489 if strcmp(prop_name, 'Opt') || strcmp( prop_name, 'param') || strcmp(prop_name, 'varargin
') 492 obj.Clear( prop_name ); 496 Reset@SVDregularizationLead( obj ); 501 %> @brief Sets the value of an attribute in the interface. 502 %> @param MemberName The name of the attribute to be set. 503 %> @param input The value to be set. 504 function Write(obj, MemberName, input) 507 if strcmp(MemberName, 'param') || strcmp(MemberName, 'params
') 508 Write@CreateLeadHamiltonians( obj, MemberName, input ); 512 obj.(MemberName) = input; 514 error(['EQuUs:
', class(obj), ':Read
'], ['No
property to write with name:
', MemberName]); 520 %> @brief Query for the value of an attribute in the interface. 521 %> @param MemberName The name of the attribute to be set. 522 %> @return Returns with the value of the attribute. 523 function ret = Read(obj, MemberName) 526 ret = obj.(MemberName); 528 error(['EQuUs:
', class(obj), ':Read
'], ['No
property to read with name:
', MemberName]); 533 %> @brief Clears the value of an attribute in the interface. 534 %> @param MemberName The name of the attribute to be cleared. 535 function Clear(obj, MemberName) 537 obj.(MemberName) = []; 539 error(['EQuUs:
', class(obj), ':Clear
'], ['No
property to clear with name:
', MemberName]); 546 %------------------------------------------------------------------ Property next_SVD_cycle
Somethimes it is needed to perform another SVD cycle to regularize the H1 matrix. ...
Property q
The tranverse momentum for transverse computations.
Property params
An instance of the structure lead_param.
Structure Opt contains the basic computational parameters used in EQuUs.
Class to create and store Hamiltonian of the translational invariant leads.
Property V
V matrix from the SVD decompozition, see Eq (41) of PRB 78, 035407.
function Transport(Energy, B)
creating the Ribbon class representing the twoterminal setup
function Read(MemberName)
Query for the value of an attribute in the interface.
function Hamiltonians(varargin)
Function to create the custom Hamiltonians for the 1D chain.
A class providing function handle to reduce the number of sites in the Hamiltonian via decimation pro...
Property varargin
list of optional parameters (see http://www.mathworks.com/help/matlab/ref/varargin.html for details)
function CreateHamiltonians(varargin)
Creates the Hamiltonians H_0 and H_1 of the lead.
Property is_SVD_transformed
true if the Hamiltonians were SVD transformed, false otherwise
Class to regulerize the H1 problem of the leads by SVD decompozition.
A class containing common basic functionalities.
A class describing the interface region between the scattering region and a lead. ...
Property Kcoupling
Hcoupling-E*Scoupling.
A class to calculate the Green functions and self energies of a translational invariant lead The nota...
Property param
An instance of the structure param.
Structure param contains data structures describing the physical parameters of the scattering center ...
Property Lead_Orientation
The orientation of the lead. Set +1 is the "incoming" direction of the propagating states is defined ...
Property K1
K1=H1-E*S1, see Eq (4) of PRB 78, 035407.
function Read(MemberName)
Query for the value of an attribute in the class.
Property Opt
An instance of structure Opt.
Property S
S matrix from the SVD decompozition, see Eq (41) of PRB 78, 035407.
function Lead(Opt, param, varargin)
Constructor of the class.
Property tolerance
SVD tolerance to identify singular values.
Property U
U matrix from the SVD decompozition, see Eq (41) of PRB 78, 035407.
Property K1adj
K1adj=H1adj-E*S1', see Eq (4) of PRB 78, 035407.
Property HamiltoniansDecimated
A logical value. True if the Hamiltonians were decimated, false otherwise.
Property K0
K0=H0-E*S0, see Eq (4) of PRB 78, 035407.
Property K1_transverse
K1_transverse=H1_transverse-E*S1_transverse.