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 utilities Utilities
20 %> @brief A
class describing an N-terminal geometry
for equilibrium calculations mostly in the zero temperature limit.
21 %> @image html two_terminal_structure.jpg
22 %> @image latex two_terminal_structure.jpg
24 %> @brief A
class describing an N-terminal geometry for equilibrium calculations mostly in the zero temperature limit.
26 %> EQuUs v4.9 or later
27 %> <tr
class=
"heading"><td colspan=
"2"><h2
class=
"groupheader"><a name=
"avail"></a>Structure described by the
class</h2></td></tr>
28 %> @image html two_terminal_structure.jpg
29 %> @image latex two_terminal_structure.jpg
30 %> The drawing represents a two-terminal structure made of two leads, of a scattering region and of two
interface regions between the leads and scattering center.
31 %> However, additional lead can be added to the structure.
32 %> Each rectangle describes a unit cell including singular and non-singular
sites.
33 %> The scattering center is, on the other hand, represented by a larger rectangle corresponding to the Hamiltonian of the scattering region.
35 %> The orientation of the lead is +1
if the lead is terminated by the
interface region in the positive direction, and -1 if the lead is terminated by the interface region in the negative direction.
41 properties (Access =
protected )
42 %> An instance of strucutre #
param 44 %> Green
operator of the scattering region
46 %> The inverse of the Green
operator G
48 %> The energy used in the calculations
50 %> The Fermi energy. Attribute #E is measured from
this value. (Use
for equilibrium calculations in the zero temperature limit.)
52 %> The transverse momentum quantum number
55 CustomHamiltoniansHandle
59 properties (Access =
public)
60 %> An instance of
class #
CreateHamiltonians or its subclass to manipulate the Hamiltonian of the scattering region
64 %> A cell array of classes #
InterfaceRegion to describe the
interface region between the leads and the scattering region
66 %> An instance of
class #
Peierls object to describe the peirls substitution in the scattering region
67 PeierlsTransform_Scatter
68 %> An instance of
class #
Peierls object to describe the peirls substitution in the leads
69 PeierlsTransform_Leads
70 %> Function
handle S = f( x,y) of the
gauge transformation in the scattering center. (S is a N x 1 vector, where N is the number of the points given by the x and y coordinates.)
72 %>
function handle for individual physical model of the leads
74 %>
function handle for individual physical model of the
interface regions
76 %> Input filename containing the computational parameters. (Obsolete)
78 %> Output filename to export the computational parameters
80 %> A
string of the working directory
84 %>
if true, no output messages are print
89 methods ( Access =
public )
92 %% Contructor of the
class 93 %> @brief Constructor of the
class.
94 %> @
param varargin Cell array of optional parameters. For details see #InputParsing.
95 %> @
return An instance of the
class 98 obj.G = []; % Greens
function of the scattering region
99 obj.Ginv = []; % inverse of G
105 obj.Interface_Regions = [];
106 obj.PeierlsTransform_Scatter = [];
107 obj.PeierlsTransform_Leads = [];
110 obj.filenameIn = [pwd, filesep, 'Basic_Input_zigzag_leads_orig.xml'];
111 obj.filenameOut = [pwd, filesep, 'Basic_Input_running_parameters.xml'];
112 obj.WorkingDir = pwd;
120 % processing the optional parameters
121 obj.InputParsing( varargin{:} );
123 obj.cCustom_Hamiltonians = [];
125 obj.display([
'EQuUs:Utils:',
class(obj),
':Creating NTerminal class.'])
128 % exporting calculation parameters into an XML format
131 % create
class intances and initializing class attributes
134 % create Hamiltonian of the scattering region (or load from external source)
135 obj.CreateNTerminalHamiltonians();
146 %> @brief Calculates the transport through the two terminal setup. Use
for development pupose only.
147 %> @
param Energy The energy value.
148 %> @
param varargin Cell array of optional parameters (https:
149 %> @
param 'constant_channels' Logical value. Set true (
default) to keep constant the number of the open channels in the leads
for each energy value, or
false otherwise.
150 %> @
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).
151 %> @
param 'decimateDyson' Logical value. Set true (
default) to decimate the
sites of the scattering region in the Dyson equation.
152 %> @
param 'PotInScatter' Obsolete parameter. Use
'scatterPotential' instead.
153 %> @
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).
154 %> @
param 'selfEnergy' 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.
155 %> @
return [1] Conductance_matrix The conductance tensor
156 %> @
return [2] ny Array of the open channel in the leads.
157 %> @
return [3] S The scattering matrix.
158 function [Conductance_matrix,ny,S] =
Transport(obj, Energy, varargin)
161 p.addParameter(
'constant_channels',
false);
162 p.addParameter(
'gfininvfromHamiltonian',
false);
163 p.addParameter(
'decimateDyson',
true);
164 p.addParameter(
'PotInScatter', []) %OBSOLETE use scatterPotential instead
165 p.addParameter(
'scatterPotential', []) %NEW overrides optional argument
'PotInScatter' 166 p.addParameter(
'selfEnergy',
false)
167 p.parse(varargin{:});
168 constant_channels = p.Results.constant_channels;
169 gfininvfromHamiltonian = p.Results.gfininvfromHamiltonian;
171 scatterPotential = p.Results.PotInScatter;
172 if ~isempty( p.Results.scatterPotential )
173 scatterPotential = p.Results.scatterPotential;
176 decimateDyson = p.Results.decimateDyson;
177 selfEnergy = p.Results.selfEnergy;
179 obj.CalcSpectralFunction(Energy, 'constant_channels', constant_channels, 'gfininvfromHamiltonian', gfininvfromHamiltonian, ...
180 'decimateDyson', decimateDyson, 'scatterPotential', scatterPotential, 'SelfEnergy', selfEnergy );
183 [S,ny] = obj.FL_handles.SmatrixCalc();
185 norma = norm(S*S'-eye(sum(ny)));
188 obj.display( ['error of the unitarity of S-matrix: ',num2str(norma)] )
192 Conductance_matrix = obj.FL_handles.Conduktance();
197 %% CalcSpectralFunction
198 %> @brief Calculates the spectral density and the Green operator.
199 %> @
param Energy The energy value.
200 %> @
param varargin Cell array of optional parameters (https:
201 %> @
param 'constant_channels' Logical value. Set true (default) to keep constant the number of the open channels in the leads for each energy value, or false otherwise.
202 %> @
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).
203 %> @
param 'decimateDyson' Logical value. Set true (default) to decimate the
sites of the scattering region in the Dyson equation.
204 %> @
param 'PotInScatter' Obsolete parameter. Use 'scatterPotential' instead.
205 %> @
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). 206 %> @
param 'selfEnergy' 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.
207 %> @
return [1] The spectral density.
208 %> @
return [2] The Green
operator.
209 function [A,G] = CalcSpectralFunction(obj, Energy, varargin)
210 obj.setEnergy( Energy )
213 p.addParameter(
'constant_channels',
false);
214 p.addParameter(
'gfininvfromHamiltonian',
false);
215 p.addParameter(
'decimateDyson',
true);
216 p.addParameter(
'PotInScatter', []) %OBSOLETE use scatterPotential instead
217 p.addParameter(
'scatterPotential', []) %NEW overrides optional argument
'PotInScatter' 218 p.addParameter(
'selfEnergy',
false)
219 p.parse(varargin{:});
220 constant_channels = p.Results.constant_channels;
221 gfininvfromHamiltonian = p.Results.gfininvfromHamiltonian;
223 scatterPotential = p.Results.PotInScatter;
224 if ~isempty( p.Results.scatterPotential )
225 scatterPotential = p.Results.scatterPotential;
228 decimateDyson = p.Results.decimateDyson;
229 selfEnergy = p.Results.selfEnergy;
230 if ~gfininvfromHamiltonian
231 obj.CalcFiniteGreensFunction('gauge_trans', true, 'onlyGinv', true);
233 obj.CalcFiniteGreensFunctionFromHamiltonian('gauge_trans', true, 'scatterPotential', scatterPotential, 'onlyGinv', true);
236 Dysonfunc = @()(obj.CustomDysonFunc( 'constant_channels', constant_channels, 'decimate', decimateDyson, 'SelfEnergy', selfEnergy ));
238 G = obj.FL_handles.DysonEq( 'CustomDyson', Dysonfunc );
240 err = MException(['EQuUs:Utils:', class(obj),':CalcSpectralFunction'], 'Error occured in calculating the Greens
function');
241 err = addCause(err, errCause);
242 save(['Error_', class(obj), '_CalcSpectralFunction.mat']);
243 obj.display(errCause.identifier, 1 );
244 obj.display( errCause.stack(1), 1 );
248 A = -imag( trace(G))/pi;
254 %> @brief Gets the coordinates of the central region
256 %> @return [2]
Coordinates of the interface region.
257 function [coordinates, coordinates_interface] = getCoordinates( obj )
259 coordinates = obj.CreateH.Read('coordinates');
262 err = MException(['EQuUs:Utils:', class(obj), ':getCoordinates'], 'Error occured when retriving the coordinates');
263 err = addCause(err, errCause);
264 save('Error_Ribbon_getCoordinatesFromRibbon.mat');
270 if isempty( obj.Interface_Regions )
271 coordinates_interface = [];
275 coordinates_interface = cell( length(obj.Interface_Regions), 1);
276 for idx = 1:length( obj.Interface_Regions )
277 coordinates_interface{idx} = obj.Interface_Regions{idx}.Read(
'coordinates');
281 err = MException([
'EQuUs:Utils:',
class(obj),
':getCoordinates',
'Error occured when retriving the coordinates']);
282 err = addCause(err, errCause);
283 save(
'Error_Ribbon_getCoordinatesFromRibbon.mat');
291 %> @brief Creates an instance of
class #
CreateHamiltonians for storing and manipulate the Hamiltonian of the the scattering region. The created object is stored in attribute #CreateH.
292 function CreateH = CreateScatter( obj )
294 %> Create an instance of
class CreateHamiltonians to create the Hamiltonian of the scattering region
298 CreateH = obj.CreateH;
301 %y Create the Hamiltonian of the scattering region
if not created
302 if ~(obj.CreateH.Read(
'HamiltoniansCreated') && (~obj.CreateH.Read(
'HamiltoniansDecimated')))
303 obj.CreateH.
CreateScatterH(
'CustomHamiltonian', obj.cCustom_Hamiltonians );
306 %> apply magnetic field in scatter
307 if ~isempty( obj.PeierlsTransform_Scatter ) && isempty(obj.q) %
for finite
q vector potential must be parallel to
q, and perpendicular to the unit cell vector
308 obj.display(
'Applying magnetic field in the scattering region')
309 obj.PeierlsTransform_Scatter.PeierlsTransfor( obj.CreateH );
315 %> @brief Sets the energy
for the calculations
316 %> @
param Energy The value of the energy in the same units as the Hamiltonian.
317 function setEnergy( obj, Energy )
320 % resetting the
class describing the N-terminal geometry
321 if ~isempty( obj.FL_handles ) && strcmpi(
class(obj.FL_handles),
'Transport_Interface')
322 obj.FL_handles.setEnergy( obj.E );
325 % reset the class describing the scattering region
330 % create interface regions
331 for idx = 1:length(obj.Interface_Regions)
332 obj.Interface_Regions{idx}.Reset();
335 % recreate the Hamiltonian of the scattering region
336 if ~isempty( obj.CreateH ) && strcmpi(
class(obj),
'NTerminal' )
337 obj.CreateNTerminalHamiltonians();
344 %> @brief Retrives the energy value (attribute
#E) used in the calculations 345 %> @
return The energy value
346 function ret = getEnergy(obj)
352 %> @brief Retrives the Fermi level (attribute #EF) used in the calculations
353 %> @
return The Femi level
354 function ret = getFermiEnergy(obj)
360 %> @brief Custom Dyson
function for a general two terminal arrangement.
361 %> @
param varargin Cell array of optional parameters (https:
362 %> @
param 'gfininv' The inverse of the Greens
function of the scattering region. For
default the inverse of the attribute #G is used.
363 %> @
param 'constant_channels' Logical value. Set true (
default) to keep constant the number of the open channels in the leads
for each energy value, or
false otherwise.
364 %> @
param 'onlyGinverz' Logical value. Set
true to calculate only the inverse of the total Green
operator, or false (
default) to calculate #G as well.
365 %> @
param 'recalculateSurface' A vector of the identification numbers of the lead surfaces to be recalculated.
366 %> @
param 'decimate' Logical value. Set true (
default) to eliminate all inner
sites in the Greens
function and keep only the selected
sites. Set
false to omit the decimation procedure.
367 %> @
param 'kulso_szabfokok' Array of
sites to be kept after the decimation procedure. (Use parameter
'keep_sites' instead)
368 %> @
param 'selfEnergy' 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.
369 %> @
param 'keep_sites' Name of
sites to be kept in the resulted Green
function (POssible values are:
'scatter',
'interface',
'lead').
370 %> @
return [1] The calculated Greens
function.
371 %> @
return [2] The inverse of the Green
operator.
372 %> @
return [3] An instance of structure #
junction_sites describing the
sites in the calculated Green
operator.
373 function [Gret, Ginverz,
junction_sites] = CustomDysonFunc( obj, varargin )
376 p.addParameter(
'gfininv', []);
377 p.addParameter(
'constant_channels',
true);
378 p.addParameter(
'onlyGinverz',
false );
379 p.addParameter(
'recalculateSurface', 1:length(obj.param.Leads) );
380 p.addParameter(
'decimate',
true );
381 p.addParameter(
'kulso_szabfokok', []); %The list of
sites to be left after the decimation procedure
382 p.addParameter(
'SelfEnergy',
false); %set
true to calculate the Dyson equation with the
self energy
383 p.addParameter(
'keep_sites',
'lead'); %identification of
sites to be kept (scatter, interface, lead)
384 p.parse(varargin{:});
385 gfininv = p.Results.gfininv;
386 constant_channels = p.Results.constant_channels;
387 onlyGinverz = p.Results.onlyGinverz;
388 recalculateSurface = p.Results.recalculateSurface;
389 decimate = p.Results.decimate;
390 kulso_szabfokok = p.Results.kulso_szabfokok;
391 useSelfEnergy = p.Results.SelfEnergy;
392 keep_sites = p.Results.keep_sites;
395 if ~isempty( obj.Ginv )
398 rcond_G = rcond(obj.G);
399 if isnan(rcond_G) || abs( rcond_G ) < 1e-15
400 gfininv = obj.inv_SVD( obj.G );
402 gfininv = inv(obj.G);
407 if ~isempty(recalculateSurface)
409 % creating interfaces for the Leads
411 shiftLeads = ones(length(obj.
param.Leads),1)*obj.E;
413 shiftLeads = ones(length(obj.
param.Leads),1)*0;
416 % creating objects describing the Leads
417 Leads = obj.FL_handles.LeadCalc('shiftLeads', shiftLeads, ...
418 'SelfEnergy', useSelfEnergy, 'SurfaceGreensFunction', ~useSelfEnergy, '
gauge_field', obj.
gauge_field, 'leads', recalculateSurface, 'q', obj.q, ...
419 'leadmodel', obj.leadmodel, 'CustomHamiltonian', obj.cCustom_Hamiltonians);
421 Leads = obj.FL_handles.Read( 'Leads' );
425 Neffs = obj.FL_handles.Get_Neff();
429 Ginverz = CalcGinverzwithSelfEnergy();
431 Ginverz = CalcGinverz();
437 GetSitesForDecimation();
438 Ginverz = obj.DecimationFunction( kulso_szabfokok, Ginverz);
448 if issparse( Ginverz )
449 err = MException(['EQuUs:Utils:', class(obj), ':CustomDysonFunc', 'Ginverz is sparse. Calculation of the inverse of a sparse matrix is not supported at this point']);
450 save('Error_Ribbon_CustomDysonFunc.mat');
454 rcond_Ginverz = rcond(Ginverz);
455 if isnan( rcond_Ginverz )
456 err = MException(['EQuUs:Utils:', class(obj), ':CustomDysonFunc'], 'NaN is the condition number for Ginverz');
457 save('Error_Ribbon_CustomDysonFunc.mat');
461 if abs( rcond_Ginverz ) < 1e-15
462 obj.display( ['EQuUs:Utils:', class(obj), ':CustomDysonFunc: Regularizing Ginverz by SVD'],1);
463 Gret = obj.inv_SVD( Ginverz );
465 Gret = inv( Ginverz );
470 %-------------------------------
471 %> calculate the inverse Greens
function 474 Gsurfinverz = cell(length(Neffs),1);
475 for idx = 1:length(Neffs)
476 Gsurfinverz{idx} = Leads{idx}.Read(
'gsurfinv');
479 K_interface = cell(length(Neffs),1);
480 K1_sc = cell(length(Neffs),1);
481 K1adj_sc = cell(length(Neffs),1);
482 K1 = cell(length(Neffs),1);
483 K1adj = cell(length(Neffs),1);
484 for idx = 1:length(recalculateSurface)
485 obj.CreateInterface( recalculateSurface(idx) );
487 for idx = 1:length( obj.Interface_Regions )
488 [K_interface{idx}, K1{idx}, K1adj{idx}, K1_sc{idx}, K1adj_sc{idx}] = obj.Interface_Regions{idx}.Get_Effective_Hamiltonians();
489 if ~isempty( obj.q ) && ~obj.Interface_Regions{idx}.Read(
'HamiltoniansDecimated')
490 K1_transverse = obj.Interface_Regions{idx}.Read(
'K1_transverse');
491 K_interface{idx} = K_interface{idx} + K1_transverse*exp(1i*obj.q) + K1_transverse
'*exp(-1i*obj.q); 495 %> getting dimensions 496 rows_interface= zeros(size(Neffs)); 497 cols_interface= zeros(size(Neffs)); 498 for idx = 1:length(Neffs) 499 rows_interface(idx) = size(K_interface{idx},1); 500 cols_interface(idx) = size(K_interface{idx},2); 502 cols_scatter = size(gfininv,2); 504 % check whether gfininv is sparse 506 Ginverz = sparse( [], [], [], sum(Neffs)+sum(rows_interface)+size(gfininv,1), sum(Neffs)+sum(cols_interface)+size(gfininv,2) ); 508 Ginverz = zeros( sum(Neffs)+sum(rows_interface)+size(gfininv,1), sum(Neffs)+sum(cols_interface)+size(gfininv,2) ); 511 for idx = 1:length(Neffs) 512 % surface Green function of leads 513 row_indexes_lead = sum(Neffs(1:idx-1))+1:sum(Neffs(1:idx)); 514 col_indexes_lead = sum(Neffs(1:idx-1))+1:sum(Neffs(1:idx)); 515 Ginverz( row_indexes_lead, col_indexes_lead ) = Gsurfinverz{idx}; 517 %interface region and coupling to the leads 518 row_indexes_interface = sum(Neffs) + (sum(rows_interface(1:idx-1))+1:sum(rows_interface(1:idx))); 519 col_indexes_interface = sum(Neffs) + (sum(cols_interface(1:idx-1))+1:sum(cols_interface(1:idx))); 520 Ginverz( row_indexes_lead, col_indexes_interface ) = -K1{idx}; 521 Ginverz( row_indexes_interface, col_indexes_lead ) = -K1adj{idx}; 522 Ginverz( row_indexes_interface, col_indexes_interface ) = -K_interface{idx}; 524 %scattering region and coupling to the interface regions 525 col_indexes_scatter = sum(Neffs) + sum(cols_interface) + (1:cols_scatter); 526 Ginverz( row_indexes_interface, col_indexes_scatter ) = -K1_sc{idx}; 527 Ginverz( col_indexes_scatter, col_indexes_interface ) = -K1adj_sc{idx}; 530 Ginverz( end-size(gfininv,1)+1:end, end-size(gfininv,2)+1:end) = gfininv; 532 % [K0_lead, K1_lead, K1adj_lead] = Leads{1}.Get_Effective_Hamiltonians(); 533 % Ginverz = [Gsurfinverz{1}, -K1_lead; -K1adj_lead, Gsurfinverz{2}]; 537 %------------------------------- 538 % calculate the inverse Greens function with the self energy 539 function Ginverz = CalcGinverzwithSelfEnergy() 540 SelfEnergy = cell(length(Neffs),1); 541 for iidx = 1:length(Neffs) 542 SelfEnergy{iidx} = Leads{iidx}.Read('Sigma
'); 545 K_interface = cell(length(Neffs),1); 546 K1_sc = cell(length(Neffs),1); 547 K1adj_sc = cell(length(Neffs),1); 548 K1 = cell(length(Neffs),1); 549 K1adj = cell(length(Neffs),1); 551 for idx = 1:length(recalculateSurface) 552 obj.CreateInterface( recalculateSurface(idx) ); 554 for idx = 1:length( obj.Interface_Regions ) 555 [K_interface{idx}, K1{idx}, K1adj{idx}, K1_sc{idx}, K1adj_sc{idx}] = obj.Interface_Regions{idx}.Get_Effective_Hamiltonians(); 556 if length(obj.q) < 2 && ~isempty( obj.q ) && ~obj.Interface_Regions{idx}.Read('HamiltoniansDecimated
') 557 K1_transverse = obj.Interface_Regions{idx}.Read('K1_transverse
'); 558 K_interface{idx} = K_interface{idx} + K1_transverse*exp(1i*obj.q) + K1_transverse'*exp(-1i*obj.q);
563 rows_interface= zeros(size(Neffs));
564 cols_interface= zeros(size(Neffs));
565 for idx = 1:length(Neffs)
566 rows_interface(idx) = size(K_interface{idx},1);
567 cols_interface(idx) = size(K_interface{idx},2);
569 cols_scatter = size(gfininv,2);
571 % check whether gfininv is sparse
573 Ginverz = sparse( [], [], [], sum(rows_interface)+size(gfininv,1), sum(cols_interface)+size(gfininv,2) );
575 Ginverz = zeros( sum(rows_interface)+size(gfininv,1), sum(cols_interface)+size(gfininv,2) );
578 for idx = 1:length(Neffs)
579 %
interface region with the self energies of the leads
580 K_interface{idx}(1:Neffs(idx), 1:Neffs(idx)) = K_interface{idx}(1:Neffs(idx), 1:Neffs(idx)) + SelfEnergy{idx};
581 row_indexes_interface = (sum(rows_interface(1:idx-1))+1:sum(rows_interface(1:idx)));
582 col_indexes_interface = (sum(cols_interface(1:idx-1))+1:sum(cols_interface(1:idx)));
583 Ginverz( row_indexes_interface, col_indexes_interface ) = -K_interface{idx};
585 %scattering region and coupling to the
interface regions
586 col_indexes_scatter = sum(cols_interface) + (1:cols_scatter);
587 Ginverz( row_indexes_interface, col_indexes_scatter ) = -K1_sc{idx};
588 Ginverz( col_indexes_scatter, col_indexes_interface ) = -K1adj_sc{idx};
591 Ginverz( end-size(gfininv,1)+1:end, end-size(gfininv,2)+1:end) = gfininv;
593 % [K0, K1, K1adj] = Leads{1}.Get_Effective_Hamiltonians();
594 % %[K_interface{1}, K1_interface{1}, K1adj_interface{1}, K1_sc{1}, K1adj_sc{1}] = obj.Interface_Regions{1}.Get_Effective_Hamiltonians();
595 % Hscatter = obj.CreateH.Read(
'Hscatter');
596 % Sscatter = obj.CreateH.Read(
'Sscatter');
597 % Kscatter = Hscatter - Sscatter*obj.E;
599 % Ginverz2 = [-K_interface{1}, -K1, zeros(size(K_interface{2}));
600 % -K1
', gfininv, -K1; 601 % zeros(size(K_interface{2})), -K1', -K_interface{2}];
609 % creating site indexes corresponding to the elements of the calculated Green
operator 614 % determine the matrix sizes
615 size_K0_leads = zeros(length(Leads), 1);
616 size_K0_interface = zeros(length(Leads), 1);
617 for kdx = 1:length(Leads)
618 K0_lead = Leads{kdx}.Get_Effective_Hamiltonians();
619 size_K0_leads(kdx) = size(K0_lead,1);
621 K0_interface = obj.Interface_Regions{kdx}.Get_Effective_Hamiltonians();
622 size_K0_interface(kdx) = size(K0_interface,1);
625 % determine junction
sites corresponding to the scattering region
632 junction_sites.Scatter.site_indexes = sum(size_K0_leads) + sum(size_K0_interface) + (1:size(gfininv,1)); %including 2*
interface regions and 2 leads
635 % determine junction
sites corresponding to the interface
639 [~, coordinates_interface] = obj.getCoordinates();
640 for kdx = 1:length(Leads)
665 % Obtain the site indexes to be decimated out and removes the unnecesarry
sites from the structure
junction_sites according to the
666 % performed decimation
667 function GetSitesForDecimation()
669 if isempty(kulso_szabfokok)
670 if strcmpi(keep_sites, 'scatter')
671 kulso_szabfokok = get_scatter_sites();
677 elseif strcmpi(keep_sites, 'interface')
678 kulso_szabfokok = get_interface_sites();
688 elseif strcmpi(keep_sites,
'lead')
689 kulso_szabfokok = get_lead_sites();
699 error([
'EQuUs:Utils:',
class(obj),
':CustomDysonFunc'], [
'Undifined Sites: ', keep_sites]);
703 %
for custom given kulso_szabfokok
705 % determine
sites in the scattering center
709 % determine
sites in the
interface regions
715 % determine
sites in the
interface regions
723 %-------------------------------------------------------------
724 % an auxilary
function to select the
sites to be kept
725 function sites_ret = SelectSites(
sites )
727 start_index = find( min(site_indexes_tmp) <= kulso_szabfokok, 1,
'first');
728 end_index = find( max(site_indexes_tmp) >= kulso_szabfokok, 1,
'last');
730 if isempty(start_index)
738 % selecting site indexes
739 if isempty(end_index)
740 sites_ret.site_indexes = 1:length(kulso_szabfokok(start_index:end));
742 sites_ret.site_indexes = 1:length(kulso_szabfokok(start_index:end_index));
745 % selecting coordinate sof the
sites 747 for iidx = 1:length(names)
748 fieldname = names(iidx);
749 if strcmpi( fieldname,
'a') || strcmpi( fieldname,
'b')
753 field_tmp = sites_ret.coordinates.(fieldname);
755 if isempty(field_tmp)
759 if isempty(end_index)
760 field_tmp = field_tmp( kulso_szabfokok(start_index:end) - min(site_indexes_tmp) + 1 );
762 field_tmp = field_tmp( kulso_szabfokok(start_index:end_index) - min(site_indexes_tmp) + 1 );
764 sites_ret.coordinates.(fieldname) = field_tmp;
773 %-------------------------------
774 % determine the site indexes of the scattering region within Ginverz
775 function kulso_szabfokok = get_scatter_sites()
779 %-------------------------------
780 % determine the site indexes of the interface region within Ginverz
781 function kulso_szabfokok = get_interface_sites()
784 size_interface = size_interface + length(
junction_sites.Interface{idx}.site_indexes);
787 kulso_szabfokok = zeros( size_interface, 1);
790 indexes = size_interface + ( 1:length(
junction_sites.Interface{idx}.site_indexes) );
797 %-------------------------------
798 % determine the site indexes of the leads within Ginverz
799 function kulso_szabfokok = get_lead_sites()
803 size_leads = size_leads + length(
junction_sites.Leads{idx}.site_indexes);
806 kulso_szabfokok = zeros( size_leads, 1);
809 indexes = size_leads + ( 1:length(
junction_sites.Leads{idx}.site_indexes) );
817 % end nested functions
822 %% DecimationFunction
823 %> @brief Performs the decimation procedure on the inverse Green
operator.
824 %> @
param kulso_szabfokok Array of
sites to be kept after the decimation.
825 %> @
param ginv The inverse of the Green
operator.
826 %> @
param varargin Cell array of optional parameters (https:
827 %> @
param 'coordinates' An instance of the structure #coordinates containing the coordinates of the
sites.
828 %> @
return [1] The matrix of the decimated inverse Greens
operator.
829 %> @
return [2] An instance of structure #coordinates containing the
sites remained after the decimation.
830 function [ret, coordinates_ret] = DecimationFunction( obj, kulso_szabfokok, ginv, varargin)
832 p.addParameter(
'coordinates', [] );
833 p.parse(varargin{:});
838 Hamiltoni = eye(size(ginv))*obj.E - ginv;
841 CreateH.
Write(
'Hscatter', Hamiltoni);
844 CreateH.
Write(
'kulso_szabfokok', kulso_szabfokok);
845 CreateH.
Write(
'coordinates', coordinates);
846 CreateH.
Write(
'HamiltoniansCreated', 1);
847 CreateH.
Write(
'HamiltoniansDecimated', 0);
850 Decimation_Procedure.DecimationFunc(obj.E, CreateH,
'Hscatter',
'kulso_szabfokok');
852 Hamiltoni = CreateH.
Read(
'Hscatter');
853 CreateH.
Clear(
'Hscatter');
855 ret = eye(size(Hamiltoni))*obj.E - Hamiltoni;
859 coordinates_ret = CreateH.
Read(
'coordinates');
865 %% GetFiniteGreensFunction
866 %> @brief Query fo the calculated Green
operator of the scattering center.
867 %> @
return [1] The Green
operator.
868 %> @
return [2] The inverse Green
operator.
869 function [Gret, Ginvret] = GetFiniteGreensFunction( obj )
874 %% CalcFiniteGreensFunction
875 %> @brief Calculates the Green
operator of the scattering region.
876 %> @
param varargin Cell array of optional parameters (https:
877 %> @
param 'gauge_trans' Logical value. Set
true to perform gauge transformation on the Green
operator and on the
Hamiltonians.
878 %> @
param 'onlyGinv' Logical value. Set
true to calculate only the inverse of the surface Greens
function #Ginv, or false (
default) to calculate #G as well. In the latter
case the attribute #Ginv is set to empty at the end.
879 function CalcFiniteGreensFunction( obj, varargin )
882 p.addParameter(
'gauge_trans',
false); % logical:
true if want to perform gauge transformation on the Green
's function and Hamiltonians 883 p.addParameter('onlyGinv
', false); 884 p.parse(varargin{:}); 885 gauge_trans = p.Results.gauge_trans; 886 onlyGinv = p.Results.onlyGinv; 888 obj.CalcFiniteGreensFunctionFromHamiltonian('gauge_trans
', gauge_trans, 'onlyGinv
', onlyGinv); 893 %% CalcFiniteGreensFunctionFromHamiltonian 894 %> @brief Calculates the Green operator of the scattering region from the whole Hamiltonian. 895 %> @param varargin Cell array of optional parameters (https://www.mathworks.com/help/matlab/ref/varargin.html): 896 %> @param 'gauge_trans
' Logical value. Set true to perform gauge transformation on the Green operator and on the Hamiltonians. 897 %> @param 'onlyGinv
' Logical value. Set true to calculate only the inverse of the surface Greens function #Ginv, or false (default) to calculate #G as well. In the latter case the attribute #Ginv is set to empty at the end. 898 %> @param 'PotInScatter
' Obsolete parameter. Use 'scatterPotential
' instead. 899 %> @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). 900 function CalcFiniteGreensFunctionFromHamiltonian( obj, varargin ) 905 p.addParameter('gauge_trans
', false); % logical: true if want to perform gauge transformation on the Green's
function and
Hamiltonians 906 p.addParameter(
'onlyGinv',
false);
907 p.addParameter(
'PotInScatter', []) %OBSOLETE use scatterPotential instead
908 p.addParameter(
'scatterPotential', []) %NEW overrides optional argument
'PotInScatter' 909 p.parse(varargin{:});
910 gauge_trans = p.Results.gauge_trans;
911 onlyGinv = p.Results.onlyGinv;
914 scatterPotential = p.Results.PotInScatter;
915 if ~isempty( p.Results.scatterPotential )
916 scatterPotential = p.Results.scatterPotential;
921 %****************************************
922 CreateH = obj.CreateH;
923 Hscatter = CreateH.Read('Hscatter');
924 Sscatter = CreateH.Read('Sscatter');
925 coordinates = CreateH.Read('coordinates');
926 kulso_szabfokok = CreateH.Read( 'kulso_szabfokok' );
929 if ~isempty(scatterPotential)
930 obj.ApplyPotentialInScatter( CreateH, scatterPotential)
931 Hscatter = CreateH.Read('Hscatter');
935 q = obj.CreateH.Read('q');
936 if length(q)<2 && ~isempty( q ) && ~obj.CreateH.Read('HamiltoniansDecimated')
937 Hscatter_transverse = obj.CreateH.Read('Hscatter_transverse');
938 Hscatter = Hscatter + Hscatter_transverse*diag(exp(-1i*q)) + Hscatter_transverse'*diag(exp(1i*q));
939 obj.CreateH.Clear('Hscatter_transverse');
942 % Applying overlap matrices
943 if isempty( Sscatter )
944 Hscatter = (sparse(1:size(Hscatter,1), 1:size(Hscatter,1), obj.E, size(Hscatter,1),size(Hscatter,1)) - Hscatter);
946 Hscatter = obj.E*Sscatter - Hscatter;
949 indexes = false( size(Hscatter,1), 1);
950 indexes(kulso_szabfokok) = true;
951 Hscatter = [ Hscatter( ~indexes, ~indexes) , Hscatter(~indexes, indexes);
952 Hscatter(indexes, ~indexes), Hscatter(indexes, indexes)];
954 obj.G = obj.partialInv( Hscatter, length(kulso_szabfokok) );
957 % Hscatter = obj.CreateH.Read('Hscatter');
958 % Hscatter = (sparse(1:size(Hscatter,1), 1:size(Hscatter,1), obj.E, size(Hscatter,1),size(Hscatter,1)) - Hscatter);
959 % obj.G = inv(Hscatter);
960 %**************************************************
962 % gauge transformation of the vector potential in the effective
Hamiltonians 965 surface_coordinates = obj.getCoordinates();
966 % gauge transformation on Green's
function 967 if ~isempty(obj.G) && ~isempty(obj.PeierlsTransform_Scatter) && isempty(obj.q)
968 % gauge transformation on the inverse Green's
function 969 obj.G = obj.PeierlsTransform_Scatter.gaugeTransformation( obj.G, surface_coordinates, obj.
gauge_field );
972 err = MException(['EQuUs:Utils:', class(obj), ':CalcFiniteGreensFunctionFromHamiltonian'], 'Unable to perform gauge transformation');
973 err = addCause(err, errCause);
974 save('Error_Ribbon_CalcFiniteGreensFunctionFromHamiltonian.mat')
980 rcond_G = rcond(obj.G);
981 if isnan(rcond_G ) || abs( rcond_G ) < 1e-15
982 obj.display( ['EQuUs:Utils:', class(obj), ':CalcFiniteGreensFunctionFromHamiltonian: Regularizing Ginv by SVD'],1);
983 obj.Ginv = obj.inv_SVD( obj.G );
985 obj.Ginv = inv(obj.G);
997 %% setInterfaceRegions
998 %> @brief Replaces the attribute Interface_Region with the given value.
1000 function setInterfaceRegions( obj, Interface_Regions )
1001 if length( Interface_Regions ) ~= length(obj.
param.Leads)
1002 err = MException(['EQuUs:Utils:', class(obj), ':setInterfaceRegions'], ['The length of Interface_Regions must be ', num2str(length(obj.
param.Leads))]);
1003 save('Error_Ribbon_setInterfaceRegions.mat');
1007 if ~iscell( Interface_Regions )
1008 err = MException(['EQuUs:Utils:', class(obj), ':setInterfaceRegions'], 'The Interface_Regions must be a cell containing of instances of class
InterfaceRegion.');
1009 save('Error_Ribbon_setInterfaceRegions.mat');
1013 obj.Interface_Regions = Interface_Regions;
1018 %> @brief Creates the
Hamiltonians for the interface regions between the leads and scattering center.
1019 %> @
param idx Identification number of the interface region.
1020 function CreateInterface( obj, idx )
1022 Interface_Region = obj.Interface_Regions{idx} ;
1024 if ~obj.Interface_Regions{idx}.Read(
'HamiltoniansCreated' )
1025 obj.display([
'EQuUs:Utils:',
class(obj),
':CreateInterface: Creating interface region ', num2str(idx)])
1026 H0 = obj.cCustom_Hamiltonians.Read(
'H0' );
1027 S0 = obj.cCustom_Hamiltonians.Read(
'S0' );
1028 H1 = obj.cCustom_Hamiltonians.Read(
'H1' );
1029 S1 = obj.cCustom_Hamiltonians.Read(
'S1' );
1030 H_coupling = obj.cCustom_Hamiltonians.Read(
'Hcoupling' );
1031 S_coupling = obj.cCustom_Hamiltonians.Read(
'Scoupling' );
1032 coordinates = obj.cCustom_Hamiltonians.Read(
'coordinates' );
1033 Interface_Region.Write(
'H0', H0{idx} );
1035 Interface_Region.Write(
'S0', S0{idx} );
1038 Interface_Region.Write(
'H1', H1{idx} );
1039 Interface_Region.Write(
'H1adj', H1{idx}
' ); 1041 Interface_Region.Write( 'S1
', S1{idx} ); 1046 Interface_Region.Write('Hcoupling
', H_coupling{idx}); 1047 Interface_Region.Write('Hcouplingadj
', H_coupling{idx}');
1048 if ~isempty( S_coupling )
1049 Interface_Region.Write(
'Scoupling', S_coupling{idx});
1052 Interface_Region.Write(
'M', size(H0{idx},1) );
1053 Interface_Region.Write(
'coordinates', coordinates{idx} );
1054 Interface_Region.Write(
'coordinates2', obj.CreateH.Read(
'coordinates') );
1056 % Transforming the
Hamiltonians according to the BdG model
1057 Interface_Region.Transform2BdG();
1059 % truncating the coupling
Hamiltonians to fit the scattering region
1060 surface_points = obj.CreateH.Read(
'kulso_szabfokok' );
1061 H_coupling = Interface_Region.Read(
'Hcoupling');
1062 indexes =
false( size(H_coupling, 2), 1);
1063 indexes( surface_points ) =
true;
1064 H_coupling = H_coupling(:,indexes);
1065 Interface_Region.Write(
'Hcoupling', H_coupling);
1066 Interface_Region.Write(
'Hcouplingadj', H_coupling
'); 1068 S_coupling = Interface_Region.Read('Scoupling
'); 1069 if ~isempty( S_coupling ) 1070 S_coupling = S_coupling(:,surface_points); 1071 Interface_Region.Write('Scoupling
', S_coupling); 1074 % keeping data on the surface sites of the scattering region 1075 coordinates2 = obj.Interface_Regions{idx}.Read('coordinates2
'); 1076 coordinates2 = coordinates2.KeepSites( indexes ); 1077 obj.Interface_Regions{idx}.Write('coordinates2
', coordinates2); 1082 % apply magnetic field in the interface region 1083 if ~isempty( obj.PeierlsTransform_Leads ) 1084 % In superconducting lead one must not include nonzero magnetic 1086 % Hamiltonians in transverse computations must remain 1087 % traslational invariant. 1088 if ~Interface_Region.isSuperconducting() %&& isempty( obj.Interface_Regions{idx}.Read('q
') ) 1089 obj.display(['EQuUs:Utils:
', class(obj), 'CreateInterface: Applying magnetic field in the
Hamiltonians of interface region
', num2str(idx)]) 1090 obj.PeierlsTransform_Leads.PeierlsTransformLeads( Interface_Region ); 1091 else %&& isempty( obj.Interface_Regions{idx}.Read('q
') ) 1092 obj.display(['EQuUs:Utils:
', class(obj), 'CreateInterface: Applying gauge transformation in the
Hamiltonians of interface region
', num2str(idx)]) 1093 obj.PeierlsTransform_Leads.gaugeTransformationOnLead( Interface_Region, obj.gauge_field ); 1097 Interface_Region.ApplyOverlapMatrices( obj.E ); 1099 % method to adjust the interface region and coupling to the scattering region by an external function. 1100 if ~isempty( obj.interfacemodel ) 1101 obj.interfacemodel( Interface_Region ); 1104 % The regularization of the interface is performed according to the Leads 1105 Leads = obj.FL_handles.Read( 'Leads
' ); 1107 Interface_Region.Calc_Effective_Hamiltonians( obj.E, 'Lead', Lead ); 1115 %% setHandlesForMagneticField 1116 %> @brief Sets the function handles of the vector potential and gauge transformation. 1117 %> @param varargin Cell array of optional parameters (https://www.mathworks.com/help/matlab/ref/varargin.html): 1118 %> @param 'scatter
' Function handle A = f( x,y) of the vector potential to be used in the scattering region (A is a N x 2 vector, where N is the number of the points given by the x and y coordinates.) 1119 %> @param 'lead
' Function handle A = f( x,y) of the vector potential to be used in the leads. (A is a N x 2 vector, where N is the number of the points given by the x and y coordinates.) 1120 %> @param 'gauge_field' Function handle S = f( x,y) of the gauge transformation. (S is a N x 1 vector, where N is the number of the points given by the x and y coordinates.) 1121 function setHandlesForMagneticField( obj, varargin ) 1124 p.addParameter('scatter
', [] ); 1125 p.addParameter('lead
', [] ); 1128 p.parse(varargin{:}); 1131 hscatter = p.Results.scatter; 1132 hlead = p.Results.lead; 1133 obj.gauge_field = p.Results.gauge_field; 1136 obj.display(['EQuUs:Utils:
', class(obj), ':setHandlesForMagneticField: Adding handles of the magnetic field to the scattering region
']) 1137 if obj.Opt.magnetic_field && isempty( obj.PeierlsTransform_Scatter ) 1138 % setting the vector potential for the scattering region 1139 obj.PeierlsTransform_Scatter = Peierls(obj.Opt, 'Vectorpotential
', hscatter); 1140 elseif ~isempty( obj.PeierlsTransform_Scatter ) && ~isempty(hscatter) 1141 % setting the vector potential for the scattering region if the Peierls class was already constructed 1142 obj.PeierlsTransform_Scatter.setVectorPotential( hscatter ); 1145 if obj.Opt.magnetic_field && isempty( obj.PeierlsTransform_Leads ) 1146 % setting the vector potential for the leads 1147 obj.PeierlsTransform_Leads = Peierls(obj.Opt, 'Vectorpotential
', hlead); 1148 elseif ~isempty( obj.PeierlsTransform_Leads ) && ~isempty(hlead) 1149 % setting the vector potential for the leads if the Peierls class was already constructed 1150 obj.PeierlsTransform_Leads.setVectorPotential( hlead ); 1158 %> @brief Creates a clone of the present object. 1159 %> @return Returns with the cloned object. 1160 function ret = CreateClone( obj ) 1162 % cerating new instance of class NTerminal 1163 ret = NTerminal( ... 1164 'filenameIn
', obj.filenameIn, ... 1165 'filenameOut
', obj.filenameOut, ... 1169 'silent
', obj.silent, ... 1171 'param', obj.param, ... 1173 'leadmodel
', obj.leadmodel, ... 1174 'interfacemodel
', obj.interfacemodel ); 1176 %> cloning the individual attributes 1177 ret.CreateH = obj.CreateH.CreateClone(); 1178 ret.FL_handles = obj.FL_handles.CreateClone(); 1179 ret.Interface_Regions = cell(size(obj.Interface_Regions)); 1180 for idx = 1:length(obj.Interface_Regions) 1181 ret.Interface_Regions{idx} = obj.Interface_Regions{idx}.CreateClone(); 1183 if ~isempty( obj.PeierlsTransform_Scatter ) 1184 ret.PeierlsTransform_Scatter = obj.PeierlsTransform_Scatter.CreateClone(); 1187 if ~isempty( obj.PeierlsTransform_Leads ) 1188 ret.PeierlsTransform_Leads = obj.PeierlsTransform_Leads.CreateClone(); 1191 %> setting the gauge field 1192 ret.gauge_field = obj.gauge_field; % function handle for the scalar field to transform the vector potential from Landauy to Landaux 1199 %> @brief Sets the structure #param in the attributes. 1200 %> @param param An instance of structure #param. 1201 function setParam(obj, param) 1204 if ~isempty( obj.CreateH ) 1205 obj.CreateH.Write('param', param); 1208 if ~isempty( obj.FL_handles ) 1209 obj.FL_handles.Write('param', param); 1212 if ~isempty( obj.cCustom_Hamiltonians ) 1213 obj.cCustom_Hamiltonians.Write('param', param); 1216 for idx = 1:length( obj.Interface_Regions ) 1217 obj.Interface_Regions{idx}.Write('param', param); 1223 %> @brief Retrives a copy of the structure #param used in the current calculations. 1224 %> @return param An instance of structure #param. 1225 function ret = getParam(obj) 1229 %% getTransverseMomentum 1230 %> @brief Retrives the value of the transverse momentum quantum number. 1231 %> @return q The transverse momentum quantum number. 1232 function ret = getTransverseMomentum(obj) 1237 %% ApplyPotentialInScatter 1238 %> @brief Applies the potential in the scattering region 1239 %> @param CreateH An instance of class #CreateHamiltonians containing the Hamiltonian of the scattering center. 1240 %> @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). 1241 function ApplyPotentialInScatter( obj, CreateH, scatterPotential ) 1244 if nargin( scatterPotential ) == 1 1245 %> calculating the potential if the scattering potential has one input arguments 1246 %> obtaining coordinates of the scattering region 1247 coordinates = CreateH.Read('coordinates
'); 1248 %> calculating the potential from the coordinates 1249 potential = scatterPotential( coordinates ); 1250 elseif nargin( scatterPotential ) == 2 1251 %> calculating the potential if the scattering potential has two input arguments 1252 potential = scatterPotential( CreateH, obj.E ); 1253 %> obtaining coordinates of the scattering region that might have been changed inside the scatterPotential function 1254 coordinates = CreateH.Read('coordinates
'); % coordinates might be changed in function scatterPotential 1256 error('EQuUs:
Ribbon:ApplyPotentialInScatter
', 'To many input arguments in
function handle scatterPotential
'); 1259 % obtaining the scattering Hamiltonian 1260 Hscatter = CreateH.Read('Hscatter
'); 1262 if isvector(potential) && length(potential) == size(Hscatter,1) 1263 % applying potential in case of a diagonal on-site potential 1265 % transforming the potential for the case of the BdG model. 1266 if isprop(coordinates, 'BdG_u
') 1267 potential(~coordinates.BdG_u) = -potential(~coordinates.BdG_u); 1270 %creating diagonal sparse matrix from the potential 1271 potential = sparse(1:size(Hscatter,1), 1:size(Hscatter,1), potential, size(Hscatter,1), size(Hscatter,1)); 1273 elseif length(size(potential)) == 2 && norm( size(potential)-size(Hscatter) ) == 0 1274 % applying potential in case, when the potential is given in a matrix form 1276 % transforming the potential for the case of the BdG model. 1277 if isprop(coordinates, 'BdG_u
') 1278 potential(~coordinates.BdG_u, ~coordinates.BdG_u) = -potential(~coordinates.BdG_u, ~coordinates.BdG_u); 1281 error(['EQuUs:utils:
', class(obj), ':ApplyPotentialInScatter
'], 'Wrong output format of the
function handle scatterPotential
'); 1285 % Ensure not to overload the memory 1286 if issparse(Hscatter) && ~issparse(potential) 1287 error(['EQuUs:utils:
', class(obj), ':ApplyPotentialInScatter
'], 'If the Hamiltonian is sparse, potential should also be sparse
'); 1290 % adding potential to the scattering region 1291 Hscatter = Hscatter + potential; 1293 % save the Hamiltonain 1294 CreateH.Write('Hscatter
', Hscatter); 1301 methods ( Access = protected ) 1303 %% CreateNTerminalHamiltonians 1304 %> @brief Extracts the Hamiltonians from the external source. 1305 function CreateNTerminalHamiltonians( obj ) 1307 % Creating attribute Custom_Hamiltonian providing Hamiltonians from external source 1308 if ~strcmpi( class(obj.cCustom_Hamiltonians), 'Custom_Hamiltonian
' ) 1309 obj.cCustom_Hamiltonians = Custom_Hamiltonians( obj.Opt, obj.param, 'EF
', obj.EF, 'CustomHamiltoniansHandle
', obj.CustomHamiltoniansHandle ); 1312 % load Hamiltonians from the external source 1313 if ~obj.cCustom_Hamiltonians.Read( 'Hamiltonians_loaded
' ) 1314 obj.cCustom_Hamiltonians.LoadHamiltonians( 'WorkingDir
', obj.WorkingDir, 'q
', obj.q ); 1317 % cerating the Hamiltonain for the scattering region 1318 obj.CreateScatter(); 1323 %> @brief Initializes the attributes of the class. 1324 function CreateHandles( obj ) 1326 % creating classes for managing the magnetic field 1327 obj.setHandlesForMagneticField(); 1329 % create class storing the Hamiltonian of the scattering region 1330 obj.CreateH = CreateHamiltonians(obj.Opt, obj.param, 'q
', obj.q); 1332 % create class for transport calculations 1333 obj.FL_handles = Transport_Interface(obj.E, obj.Opt, obj.param, 'CreateH
', obj.CreateH, 'PeierlsTransform
', obj.PeierlsTransform_Leads ); 1335 % read out structure Opt containing the computational parameters 1336 obj.Opt = obj.FL_handles.Read( 'Opt' ); 1338 % creating classes of the interface regions 1339 obj.Interface_Regions = cell(length(obj.param.Leads),1); 1340 for idx = 1:length(obj.Interface_Regions) 1341 Lead_Orientation = obj.param.Leads{idx}.Lead_Orientation; 1342 obj.Interface_Regions{idx} = InterfaceRegion(obj.Opt, obj.param, 'Hanyadik_Lead
', idx, 'q
', obj.q, 'Lead_Orientation
', Lead_Orientation); 1348 end % protected methods 1350 methods (Access=protected) 1354 %> @brief Parses the optional parameters for the class constructor. 1355 %> @param varargin Cell array of optional parameters (https://www.mathworks.com/help/matlab/ref/varargin.html): 1356 %> @param 'filenameIn
' The input filename containing the computational parameters. (Use parameters 'Op
' and 'param' instead) 1357 %> @param 'filenameOut
' The output filename to export the computational parameters. 1358 %> @param 'WorkingDir
' The absolute path to the working directoy. 1359 %> @param 'CustomHamiltoniansHandle
' function handle for the custom Hamiltonians. Has the same inputs as #Custom_Hamiltonians.LoadHamiltonians and output values defined by the example #Hamiltonians. 1360 %> @param 'E
' The energy value used in the calculations (in the same units as the Hamiltonian). 1361 %> @param 'EF
' The Fermi energy in the same units as the Hamiltonian. Attribute #E is measured from this value. (Use for equilibrium calculations in the zero temperature limit. Overrides the one comming from the external source) 1362 %> @param 'silent
' Set true to suppress output messages. 1363 %> @param 'leadmodel
' A function handle #Lead=f( idx, E, varargin ) of the alternative lead model with equivalent inputs and return values as #Transport_Interface.SurfaceGreenFunctionCalculator and with E standing for the energy. 1364 %> @param 'interfacemodel
' A function handle f( #InterfaceRegion ) to manually adjus the interface regions. (Usefull when 'leadmodel
' is also given. For example see @InterfaceModel) 1365 %> @param 'Opt' An instance of the structure #Opt. 1366 %> @param 'param' An instance of the structure #param. 1367 %> @param 'q
' The transverse momentum quantum number. 1368 function InputParsing( obj, varargin) 1371 p.addParameter('filenameIn
', obj.filenameIn, @ischar); 1372 p.addParameter('filenameOut
', obj.filenameOut, @ischar); 1373 p.addParameter('WorkingDir
', obj.WorkingDir, @ischar); 1374 p.addParameter('CustomHamiltoniansHandle
', obj.CustomHamiltoniansHandle); %function handle for the custom Hamiltonians 1375 p.addParameter('E
', obj.E, @isscalar); 1376 p.addParameter('EF
', obj.EF); 1377 p.addParameter('silent
', obj.silent); 1378 p.addParameter('leadmodel
', obj.leadmodel); %individual physical model for the contacts 1379 p.addParameter('interfacemodel
', obj.interfacemodel); %individual physical model for the interface regions 1380 p.addParameter('Opt', obj.Opt); 1381 p.addParameter('param', obj.param); 1382 p.addParameter('q
', obj.q); 1384 p.parse(varargin{:}); 1387 obj.filenameIn = p.Results.filenameIn; 1388 obj.filenameOut = p.Results.filenameOut; 1389 obj.WorkingDir = p.Results.WorkingDir; 1390 obj.CustomHamiltoniansHandle = p.Results.CustomHamiltoniansHandle; 1391 obj.E = p.Results.E; 1392 obj.EF = p.Results.EF; 1393 obj.silent = p.Results.silent; 1394 obj.leadmodel = p.Results.leadmodel; 1395 obj.interfacemodel = p.Results.interfacemodel; 1396 obj.q = p.Results.q; 1397 obj.Opt = p.Results.Opt; 1398 obj.param = p.Results.param; 1401 if isempty(obj.Opt) || isempty(obj.param) 1402 [ obj.Opt, obj.param ] = parseInput( obj.filenameIn ); 1407 end % methods private A class describing an N-terminal geometry for equilibrium calculations mostly in the zero temperature...
function LoadHamiltonians(varargin)
Obtain the Hamiltonians from the external source.
A class containing methodes for displaying several standard messages.
function Write(MemberName, input)
Sets the value of an attribute in the class.
Structure Opt contains the basic computational parameters used in EQuUs.
Property gauge
String containing the name of the built-in gauge ('LandauX', 'LandauY').
Class to create and store Hamiltonian of the translational invariant leads.
A class for calculations on a ribbon of finite width for equilibrium calculations mostly in the zero ...
Property H1
The coupling Hamiltonian between the unit cells.
sites Interface
Structure sites describing the geometry of the interface regions.
function Transport(Energy, B)
Calculates the conductance at a given energy value.
A class to evaluate the Dyson equation and to calculate the scattering matrix for equilibrium process...
function Hamiltonians(varargin)
Function to create the custom Hamiltonians for the 1D chain.
function CreateScatterH(varargin)
Creates a Hamiltonian of a rectangle shaped area.
A class providing function handle to reduce the number of sites in the Hamiltonian via decimation pro...
function createOutput(filename, Opt, param)
This function creates an output file containing the running parameters.
function CreateHamiltonians(Opt, param, varargin)
Constructor of the class.
Property q
A transverse momentum.
A class to import custom Hamiltonians provided by other codes or created manually
A class containing common basic functionalities.
function Clear(MemberName)
Clears the value of an attribute in the class.
Property Hcoupling
Coupling Hamiltonian from the interface region to the scattering region.
A class describing the interface region between the scattering region and a lead.
A class to calculate the Green functions and self energies of a translational invariant lead The nota...
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 ...
site_indexes
An array containing the site indexes of the given system part.
Structure sites contains data to identify the individual sites in a matrix.
sites Leads
Structure sites describing the geometry of the leads.
A class responsible for the Peierls and gauge transformations.
function gauge_field(x, y, eta_B, Aconst, height, lattice_constant)
Scalar gauge field connecting Landaux and Landauy gauges.
coordinates
An instance of structure coordinates describing the geometry.
sites Scatter
Structure sites describing the geometry of the scattering region.
Structure containing the coordinates and other quantum number identifiers of the sites in the Hamilto...
function Read(MemberName)
Query for the value of an attribute in the class.
function structures(name)
Structure junction_sites contains data to identify the individual sites of the leads,...
A class to create and store Hamiltonian of the scattering region.