Eötvös Quantum Utilities  v5.0.144
Providing the Horsepowers in the Quantum Realm
Ribbon.m
Go to the documentation of this file.
1 %% Eotvos Quantum Transport Utilities - Ribbon
2 % Copyright (C) 2009-2015 Peter Rakyta, Ph.D.
3 %
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.
8 %
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.
13 %
14 % You should have received a copy of the GNU General Public License
15 % along with this program. If not, see http://www.gnu.org/licenses/.
16 %
17 %> @addtogroup utilities Utilities
18 %> @{
19 %> @file Ribbon.m
20 %> @brief A class for calculations on a ribbon of finite width for equilibrium calculations mostly in the zero temperature limit.
21 %> @image html Ribbon_structure.jpg
22 %> @image latex Ribbon_structure.jpg
23 %> @}
24 %> @brief A class for calculations on a ribbon of finite width for equilibrium calculations mostly in the zero temperature limit.
25 %> <tr class="heading"><td colspan="2"><h2 class="groupheader"><a name="avail"></a>Structure described by the class</h2></td></tr>
26 %> @image html Ribbon_structure.jpg
27 %> @image latex Ribbon_structure.jpg
28 %> 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.
29 %> Each rectangle describes a unit cell including singular and non-singular sites.
30 %> The scattering center is also described by a set of identical unit cells, but arbitrary potential can be used.
31 %> Arrows indicate the hopping direction stored by the attributes in the corresponding classes (see attributes #CreateLeadHamiltonians.H1 and #InterfaceRegion.Hcoupling for details).
32 %> 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.
33 %> (see attribute #CreateLeadHamiltonians.Lead_Orientation for details)
34 classdef Ribbon < NTerminal
35 
36 
37  properties ( Access = public )
38  %> An instance of class #Lead (or its subclass) describing the unit cell of the scattering region
39  Scatter_UC
40  %> A function handle pot = f( #Coordinates ) or pot=f( #CreateLeadHamiltonians, Energy) of the transverse potential applied in the lead. (Instead of #CreateLeadHamiltonians can be used its derived class)
41  transversepotential
42  %> width of the scattering region (number of the nonsingular atomic sites in the cross section)
43  width
44  %> height (length) of the scattering region (number of unit cells)
45  height
46  %> the shift of the coordinates of the sites (two component vector)
47  shift
48  end
49 
50 
51 methods ( Access = public )
52 %% Contructor of the class
53 %> @brief Constructor of the class.
54 %> @param varargin Cell array of optional parameters. For details see #InputParsing.
55 %> @return An instance of the class
56  function obj = Ribbon( varargin )
57  obj = obj@NTerminal();
58 
59 
60  obj.param = [];
61  obj.Scatter_UC = [];
62  obj.transversepotential = [];
63  obj.width = [];
64  obj.height = [];
65  obj.shift = 0;
66 
67 
68  if strcmpi( class(obj), 'Ribbon')
69  % processing the optional parameters
70  obj.InputParsing(varargin{:})
71 
72 
73  obj.display(['EQuUs:Utils:', class(obj), ':Ribbon: Creating a Ribbon object'])
74 
75  % create the shape of the scattering region
76  obj.createShape();
77 
78  % set the Fermi level
79  obj.setFermiEnergy();
80 
81  % exporting calculation parameters into an XML format
82  createOutput( obj.filenameOut, obj.Opt, obj.param );
83 
84  % create class intances and initializing class attributes
85  obj.CreateHandles();
86 
87  % create Hamiltonians and coordinates of the unit cells
88  obj.CreateRibbon();
89  end
90 
91 
92 
93  end
94 
95 %% Transport
96 %> @brief Calculates the transport through the two terminal setup on two dimensional lattices. Use for development pupose only.
97 %> @param Energy The energy value.
98 %> @param varargin Cell array of optional parameters (https://www.mathworks.com/help/matlab/ref/varargin.html):
99 %> @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.
100 %> @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).
101 %> @param 'decimateDyson' Logical value. Set true (default) to decimate the sites of the scattering region in the Dyson equation.
102 %> @param 'PotInScatter' Obsolete parameter. Use 'scatterPotential' instead.
103 %> @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).
104 %> @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.
105 %> @param 'Smatrix' Set true (default) to calculate the conductance by using the scattering matrix via #Transport_Interface.Conduktance, or false to use function #Transport_Interface.Conductance2. (The latter one also works with complex energies.)
106 %> @return [1] Conductivity The calculated conductivity.
107 %> @return [2] aspect_ratio The aspect ratio W/L of the junction.
108 %> @return [3] Conductance The conductance tensor
109 %> @return [4] ny Array of the open channel in the leads.
110 %> @return [5] DeltaC Error of the unitarity.
111 %> @return [6] S The scattering matrix.
112  function [Conductivity,aspect_ratio,Conductance,ny,DeltaC,S] = Transport(obj, Energy, varargin)
113 
114  p = inputParser;
115  p.addParameter('constant_channels', false);
116  p.addParameter('gfininvfromHamiltonian', false);
117  p.addParameter('decimateDyson', true);
118  p.addParameter('PotInScatter', []) %OBSOLETE use scatterPotential instead
119  p.addParameter('scatterPotential', []) %NEW overrides optional argument 'PotInScatter'
120  p.addParameter('selfEnergy', false)
121  p.addParameter('Smatrix', true) %logcal value to use the S-matrix method to calculate the conductance or not.
122  p.parse(varargin{:});
123  constant_channels = p.Results.constant_channels;
124  gfininvfromHamiltonian = p.Results.gfininvfromHamiltonian;
125 
126  scatterPotential = p.Results.PotInScatter;
127  if ~isempty( p.Results.scatterPotential )
128  scatterPotential = p.Results.scatterPotential;
129  end
130 
131  decimateDyson = p.Results.decimateDyson;
132  selfEnergy = p.Results.selfEnergy;
133  Smatrix = p.Results.Smatrix;
134 
135  obj.CalcSpectralFunction(Energy, 'constant_channels', constant_channels, 'gfininvfromHamiltonian', gfininvfromHamiltonian, ...
136  'decimateDyson', decimateDyson, 'scatterPotential', scatterPotential, 'SelfEnergy', selfEnergy);
137 
138  if strcmpi(obj.Opt.Lattice_Type, 'Graphene' ) && strcmp(obj.param.scatter.End_Type, 'A')
139  aspect_ratio = (obj.width*3/2)/((obj.height+3)*sqrt(3));
140  elseif strcmpi(obj.Opt.Lattice_Type, 'Graphene' ) && strcmp(obj.param.scatter.End_Type, 'Z')
141  aspect_ratio = ((obj.width-0.5)*sqrt(3))/((obj.height+3)*3);
142  elseif strcmpi(obj.Opt.Lattice_Type, 'Square' )
143  aspect_ratio = obj.width/obj.height;
144  else
145  aspect_ratio = 1;
146  end
147 
148  if Smatrix
149  [S,ny] = obj.FL_handles.SmatrixCalc();
150 
151  norma = norm(S*S'-eye(sum(ny)));
152 
153  if norma >= 1e-3
154  obj.display( ['error of the unitarity of S-matrix: ',num2str(norma)] )
155  end
156 
157  if ny(1) ~= ny(2)
158  obj.display( ['openchannels do not match: ',num2str(ny(1)), ' ', num2str(ny(2))] )
159  end
160 
161 
162 
163  conductance = obj.FL_handles.Conduktance();
164  conductance = abs([conductance(1,:), conductance(2,:)]);
165  DeltaC = std(conductance);
166 
167  C = mean(conductance);
168  Conductivity = C/aspect_ratio;
169  Conductance = C;
170 
171  obj.display( ['aspect ratio = ', num2str(aspect_ratio), ' conductance = ', num2str(C), ' open_channels= ', num2str(ny(1)), ' conductivity = ', num2str(Conductivity)])
172 
173  else
174  Conductance = obj.FL_handles.Conductance2();
175  DeltaC = NaN;
176  Conductivity = Conductance/aspect_ratio;
177  ny = NaN;
178  S = [];
179  obj.display( ['aspect ratio = ', num2str(aspect_ratio), ' conductance = ', num2str(Conductance), ' conductivity = ', num2str(Conductivity)])
180  end
181 
182  end
183 
184 %% getCoordinates
185 %> @brief Gets the coordinates of the central region
186 %> @return [1] Coordinates of the central region.
187 %> @return [2] Coordinates of the interface region.
188 function [coordinates, coordinates_interface] = getCoordinates( obj )
189  try
190  if isempty( obj.CreateH ) || ~obj.CreateH.Read('HamiltoniansCreated')
191  obj.CreateScatter()
192  end
193 
194  coordinates = obj.CreateH.Read('coordinates');
195  non_singular_sites = obj.CreateH.Read('kulso_szabfokok');
196  indexes = false( size(coordinates.x) );
197  indexes(non_singular_sites) = true;
198 
199  coordinates = coordinates.KeepSites( indexes );
200 
201  catch errCause
202  err = MException('EQuUs:Utils:Ribbon:getCoordinatesFromRibbon', 'Error occured when retriving the coordinates of the ribbon.');
203  err = addCause(err, errCause);
204  save('Error_Ribbon_getCoordinatesFromRibbon.mat');
205  throw( err );
206  end
207 
208 
209  try
210  if isempty( obj.Interface_Regions )
211  coordinates_interface = [];
212  return
213  end
214 
215  coordinates_interface = cell( length(obj.Interface_Regions), 1);
216  for idx = 1:length( obj.Interface_Regions )
217  coordinates_interface{idx} = obj.Interface_Regions{idx}.Read('coordinates');
218  end
219 
220  catch errCause
221  err = MException('EQuUs:Utils:Ribbon:getCoordinatesFromRibbon', 'Error occured when retriving the coordinates of the interface region from a Ribbon interface');
222  err = addCause(err, errCause);
223  save('Error_Ribbon_getCoordinatesFromRibbon.mat');
224  throw( err );
225  end
226 
227 end
228 
229 %% ShiftCoordinates
230 %> @brief Shifts the coordinates of the sites in the ribbon by an integer multiple of the lattice vector. The coordinates of the Leads are automatically adjusted later
231 %> @param shift An integer value.
232  function ShiftCoordinates( obj, shift )
233  obj.Scatter_UC.ShiftCoordinates( shift );
234  if ~isempty( obj.Interface_Regions )
235  for idx = 1:length(obj.Interface_Regions)
236  obj.Interface_Regions{idx}.ShiftCoordinates( shift );
237  end
238  end
239 
240  obj.shift = obj.shift + shift;
241  end
242 
243 %% CreateScatter
244 %> @brief Initializes class #CreateHamiltonians for storing and manipulate the Hamiltonian of the the scattering region. The created object is stored in attribute #CreateH.
245  function Scatter_UC = CreateScatter( obj )
246  Scatter_UC = obj.Scatter_UC.CreateClone('empty', true);
247 
248  % create Hamiltonian of one unit cell of the scattering region
249  Scatter_UC.CreateHamiltonians( 'toSave', 0);
250  Scatter_UC.ShiftCoordinates( obj.shift );
251 
252  %applying transverse potential
253  obj.ApplyTransversePotential( Scatter_UC )
254 
255 
256  % apply magnetic field in the unit cell of the scattering region
257  % can be applied if the vector potential is identical in each unit cells
258  if ~isempty( obj.PeierlsTransform_Scatter ) && obj.Opt.magnetic_field_trans_invariant %for finite q the vector potential must be parallel to q, and perpendicular to the unit cell vector
259  obj.display(['EQuUs:Utils:',class(obj),':CreateScatter: Applying magnetic field in the unit cell of the scattering region']);
260  obj.PeierlsTransform_Scatter.PeierlsTransformLeads( Scatter_UC );
261  end
262 
263  % Create the Hamiltonian of the scattering region
264  createH = CreateHamiltonians(obj.Opt, obj.param, 'q', obj.q);
265  createH.CreateScatterH( 'Scatter_UC', Scatter_UC );
266 
267  % apply magnetic field in the whole Hamiltonian of the scattering region
268  % can be applied for non-translational invariant vector potentials
269  if ~isempty( obj.PeierlsTransform_Scatter ) && ~obj.Opt.magnetic_field_trans_invariant
270  obj.display(['EQuUs:Utils:',class(obj),':CreateScatter: Applying magnetic field in the whole Hamiltonian of the scattering region']);
271  obj.PeierlsTransform_Scatter.PeierlsTransform( createH );
272  end
273 
274 
275  obj.CreateH = createH;
276 
277 
278  % obtaining the sites of the scattering region that are coupled to the leads
279  H0 = Scatter_UC.Read('H0');
280  H1 = Scatter_UC.Read('H1');
281 
282  [rows, cols] = find( H1 );
283  rows = unique(rows); % cols identical to non_singular_sites
284  cols = unique(cols); % cols identical to non_singular_sites
285 
286  % identify sites that would be directly connected to the leads
287  non_singular_sites = [reshape(cols, 1, length(cols)), (obj.height-1)*size(H0,1)+reshape(rows, 1, length(rows))];
288  createH.Write('kulso_szabfokok', non_singular_sites );
289 
290 
291  end
292 
293 %% setEnergy
294 %> @brief Sets the energy for the calculations
295 %> @param Energy The value of the energy in the same units as the Hamiltonian.
296  function setEnergy( obj, Energy )
297 
298  setEnergy@NTerminal( obj, Energy );
299 
300  if ~isempty( obj.Scatter_UC ) && strcmpi(class(obj.Scatter_UC), 'Lead')
301  obj.Scatter_UC.Reset();
302  end
303 
304  % recreate the Hamiltonian of the scattering region
305  if ~isempty(obj.Scatter_UC)
306  obj.CreateRibbon();
307  end
308 
309  end
310 
311 
312 %% CustomDysonFunc
313 %> @brief Custom Dyson function for a two terminal arrangement on a two dimensional lattice.
314 %> @param varargin Cell array of optional parameters (https://www.mathworks.com/help/matlab/ref/varargin.html):
315 %> @param 'gfininv' The inverse of the Greens function of the scattering region. For default the inverse of the attribute #G is used.
316 %> @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.
317 %> @param 'onlyGinverz' Logical value. Set true to calculate only the inverse of the total Green operator, or false (default) to calculate #G as well.
318 %> @param 'recalculateSurface' A vector of the identification numbers of the lead surfaces to be recalculated.
319 %> @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.
320 %> @param 'kulso_szabfokok' Array of sites to be kept after the decimation procedure. (Use parameter 'keep_sites' instead)
321 %> @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.
322 %> @param 'keep_sites' Name of sites to be kept in the resulted Green function (Possible values are: 'scatter', 'interface', 'lead').
323 %> @param 'UseHamiltonian' Set true if the interface region is matched to the whole Hamiltonian of the scattering center, or false (default) if the surface Green operator of the scattering center is used in the calculations.
324 %> @return [1] The calculated Greens function.
325 %> @return [2] The inverse of the Green operator.
326 %> @return [3] An instance of structure #junction_sites describing the sites in the calculated Green operator.
327  function [Gret, Ginverz, junction_sites] = CustomDysonFunc( obj, varargin ) %NEW output
328 
329  p = inputParser;
330  p.addParameter('gfininv', []);
331  p.addParameter('constant_channels', true);
332  p.addParameter('onlyGinverz', false );
333  p.addParameter('recalculateSurface', [1 2] );
334  p.addParameter('decimate', true );
335  p.addParameter('kulso_szabfokok', []); %The list of sites to be left after the decimation procedure
336  p.addParameter('SelfEnergy', false); %set true to calculate the Dyson equation with the self energy
337  p.addParameter('keep_sites', 'lead'); %Name of sites to be kept (scatter, interface, lead)
338  p.addParameter('UseHamiltonian', false); %true if the interface region is matched to the whole Hamiltonian of the scattering center, false if the surface Green operator of the scattering center is used in the calculations.
339  p.parse(varargin{:});
340  gfininv = p.Results.gfininv;
341  constant_channels = p.Results.constant_channels;
342  onlyGinverz = p.Results.onlyGinverz;
343  recalculateSurface = p.Results.recalculateSurface;
344  decimate = p.Results.decimate;
345  kulso_szabfokok = p.Results.kulso_szabfokok;
346  useSelfEnergy = p.Results.SelfEnergy;
347  keep_sites = p.Results.keep_sites;
348  UseHamiltonian = p.Results.UseHamiltonian;
349 
350 
351  if ~isempty(recalculateSurface)
352 
353  % creating interfaces for the Leads
354  if constant_channels
355  shiftLeads = ones(length(obj.param.Leads),1)*obj.E;
356  else
357  shiftLeads = ones(length(obj.param.Leads),1)*0;
358  end
359 
360  % creating Lead instaces and calculating the retarded surface Green operator/self-energy
361  coordinates_shift = [-2, obj.height+1] + obj.shift;
362  obj.FL_handles.LeadCalc('coordinates_shift', coordinates_shift, 'shiftLeads', shiftLeads, 'transversepotential', obj.transversepotential, ...
363  'SelfEnergy', useSelfEnergy, 'SurfaceGreensFunction', ~useSelfEnergy, 'gauge_field', obj.gauge_field, 'leads', recalculateSurface, 'q', obj.q, ...
364  'leadmodel', obj.leadmodel);
365 
366  for idx = 1:length(recalculateSurface)
367  obj.CreateInterface( recalculateSurface(idx), 'UseHamiltonian', UseHamiltonian );
368  end
369 
370  end
371 
372  [Gret, Ginverz, junction_sites] = CustomDysonFunc@NTerminal(obj, 'gfininv', gfininv, ...
373  'onlyGinverz', onlyGinverz, ...
374  'recalculateSurface', [], ...
375  'decimate', decimate, ...
376  'kulso_szabfokok', kulso_szabfokok, ...
377  'SelfEnergy', useSelfEnergy, ...
378  'keep_sites', keep_sites);
379 
380 
381  end
382 
383 %% CalcFiniteGreensFunction
384 %> @brief Calculates the Green operator of the scattering region by the fast way (see PRB 90, 125428 (2014)).
385 %> @param varargin Cell array of optional parameters identical to #NTerminal.CalcFiniteGreensFunction.
386  function CalcFiniteGreensFunction( obj, varargin )
387 
388  p = inputParser;
389  p.addParameter('gauge_trans', false); % logical: true if want to perform gauge transformation on the Green's function and Hamiltonians
390  p.addParameter('onlyGinv', false);
391  p.parse(varargin{:});
392  gauge_trans = p.Results.gauge_trans;
393  onlyGinv = p.Results.onlyGinv;
394 
395 
396  if obj.Opt.magnetic_field && ~obj.Opt.magnetic_field_trans_invariant
397  warning( ['EQuUs:Utils:', class(obj), ':CalcFiniteGreensFunction'], 'The vector potential should be translational invariant in this calculation!' )
398  end
399 
400  % recreate The Hamiltonian of the unit cell
401  obj.CreateRibbon();
402 
403  if ~obj.Scatter_UC.Read('OverlapApplied')
404  obj.Scatter_UC.ApplyOverlapMatrices(obj.E);
405  end
406 
407  % getting the Hamiltonian for the edge slabs
408  [K0, K1, K1adj] = obj.Scatter_UC.qDependentHamiltonians();
409 
410 
411  obj.display( ['EQuUs:Utils:', class(obj), ':CalcFiniteGreensFunction: Solving the eigenproblem in the scattering region'] )
412 
413  % Trukkos sajatertek
414  obj.display(['EQuUs:Utils:', class(obj), ':CalcFiniteGreensFunction: Eigenvalues of the scattering region'])
415  obj.Scatter_UC.TrukkosSajatertekek(obj.E);
416  % group velocity
417  obj.display(['EQuUs:Utils:', class(obj), ':CalcFiniteGreensFunction: Group velocity for the scattering region'])
418  obj.Scatter_UC.Group_Velocity();
419 
420 
421  % transforming the Hamiltonians by SVD if necessary
422  if obj.Scatter_UC.Read('is_SVD_transformed')
423  V = obj.Scatter_UC.Get_V();
424  K0 = V'*K0*V;
425  K1 = V'*K1*V;
426  K1adj = V'*K1adj*V;
427  end
428 
429 
430  %> the first two and last two slabs are added manualy at the end
431  if ~obj.Scatter_UC.Read('is_SVD_transformed')
432  z1 = 1;
433  z2 = obj.height-2;
434  else
435  z1 = 2;
436  z2 = obj.height-3;
437  end
438 
439  obj.G = [];
440  obj.Ginv = [];
441 
442  obj.Scatter_UC.FiniteGreenFunction(z1,z2, 'onlygfininv', true);
443  obj.Ginv = obj.Scatter_UC.Read( 'gfininv' );
444 
445  % adding to the scattering region the first set of transition layers
446  Neff = obj.Scatter_UC.Get_Neff();
447  non_singular_sites = obj.Scatter_UC.Read( 'kulso_szabfokok' );
448  if isempty( non_singular_sites )
449  non_singular_sites = 1:Neff;
450  end
451 
452  non_singular_sites_edges = [non_singular_sites, size(K0,1)+1:2*size(K0,1)];
453  ginv_edges = -[K0, K1; K1adj, K0];
454  ginv_edges = obj.DecimationFunction( non_singular_sites_edges, ginv_edges );
455 
456  %Ginv = [first_slab, H1, 0;
457  % H1', invG, H1;
458  % 0 , H1', last_slab];
459 
460  obj.Ginv = [ginv_edges(1:Neff,1:Neff), [ginv_edges(1:Neff, Neff+non_singular_sites), zeros(Neff, Neff+size(K0,2))]; ...
461  [ginv_edges(Neff+non_singular_sites, 1:Neff); zeros(Neff, Neff)], obj.Ginv, [zeros(Neff, size(ginv_edges,2)-Neff); ginv_edges(1:Neff, Neff+1:end)]; ...
462  [zeros(size(K0,2),2*Neff), ginv_edges(Neff+1:end, 1:Neff)], ginv_edges(Neff+1:end, Neff+1:end)];
463 
464  [rows, ~] = find( K1 );
465  rows = unique(rows); % cols identical to non_singular_sites
466 
467  %non_singular_sites_Ginv = [1:length(non_singular_sites), size(obj.Ginv,1)-size(K0,1)+1:size(obj.Ginv,1)];
468  non_singular_sites_Ginv = [1:length(non_singular_sites), size(obj.Ginv,1)-size(K0,1)+reshape(rows, 1, length(rows))];
469  obj.Ginv = obj.DecimationFunction( non_singular_sites_Ginv, obj.Ginv );
470 
471 
472  % terminate the scattering region by the first and last slabs and transform back to the normal space from the SVD representation
473  if obj.Scatter_UC.Read('is_SVD_transformed')
474  %> adding the first and last slab to
475  %Ginv = [first_slab, H1, 0;
476  % H1', invG, H1;
477  % 0 , H1', last_slab];
478 
479  obj.Ginv = [K0, [K1(:, 1:Neff), zeros(size(K0,1), 2*size(K0,2))]; ...
480  [K1adj(1:Neff,:); zeros(size(K0))], obj.Ginv, [zeros(Neff, size(K0,2)); K1]; ...
481  zeros(size(K0)), [zeros(size(K0,1), Neff), K1adj], K0];
482 
483 
484  non_singular_sites_Ginv = [1:size(K0,1), size(obj.Ginv,1)-size(K0,1)+1:size(obj.Ginv,1)];
485  obj.Ginv = obj.DecimationFunction( non_singular_sites_Ginv, obj.Ginv );
486 
487  % transform back to the normal space
488  V_tot = [V, zeros(size(K0)); zeros(size(K0)), V];
489  obj.Ginv = V_tot*obj.Ginv*V_tot';
490  end
491 
492 
493 %disp( obj.Ginv )
494  % gauge transformation of the vector potential in the effective Hamiltonians
495  if gauge_trans
496  try
497  % gauge transformation on Green's function
498  if ~isempty(obj.Ginv) && ~isempty(obj.PeierlsTransform_Scatter) && isempty(obj.q) && ~isempty(obj.gauge_field)
499  coordinates_scatter = obj.getCoordinates();
500  %> gauge transformation on the inverse Green's function
501  obj.Ginv = obj.PeierlsTransform_Scatter.gaugeTransformation( obj.Ginv, coordinates_scatter, obj.gauge_field );
502  end
503  catch errCause
504  err = MException(['EQuUs:Utils:', class(obj), ':CalcFiniteGreensFunction'], 'Unable to perform gauge transformation');
505  err = addCause(err, errCause);
506  save('Error_Ribbon_CalcFiniteGreensFunction.mat')
507  throw(err);
508  end
509  end
510 
511  if ~onlyGinv
512  rcond_Ginv = rcond(obj.Ginv);
513  if isnan(rcond_Ginv ) || abs( rcond_Ginv ) < 1e-15
514  obj.display( ['EQuUs:Utils:', class(obj), ':CalcFiniteGreensFunction: Regularizing Ginv by SVD'], 1);
515  obj.G = obj.inv_SVD( obj.Ginv );
516  else
517  obj.G = inv(obj.Ginv);
518  end
519  obj.Ginv = [];
520  else
521  obj.G = [];
522  end
523 
524  end
525 
526 %% CalcFiniteGreensFunctionFromHamiltonian
527 %> @brief Calculates the Green operator of the scattering region from the whole Hamiltonian.
528 %> @param varargin Cell array of optional parameters (https://www.mathworks.com/help/matlab/ref/varargin.html):
529 %> @param 'gauge_trans' Logical value. Set true to perform gauge transformation on the Green operator and on the Hamiltonians.
530 %> @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.
531 %> @param 'PotInScatter' Obsolete parameter. Use 'scatterPotential' instead.
532 %> @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).
533  function CalcFiniteGreensFunctionFromHamiltonian( obj, varargin )
534 
535  p = inputParser;
536  p.addParameter('gauge_trans', false); % logical: true if want to perform gauge transformation on the Green's function and Hamiltonians
537  p.addParameter('onlyGinv', false);
538  p.addParameter('PotInScatter', []) %OBSOLETE use scatterPotential instead
539  p.addParameter('scatterPotential', []) %NEW overrides optional argument 'PotInScatter', might be a cell array of function handles
540  p.parse(varargin{:});
541  gauge_trans = p.Results.gauge_trans;
542  onlyGinv = p.Results.onlyGinv;
543 
544  scatterPotential = p.Results.PotInScatter;
545  if ~isempty( p.Results.scatterPotential )
546  scatterPotential = p.Results.scatterPotential;
547  end
548 
549  % creating the Hamiltonian of the scattering region
550  obj.CreateScatter();
551  CreateH = obj.CreateH.CreateClone();
552 
553  % obtaining the Hamiltonian of the scattering region
554  Hscatter = CreateH.Read('Hscatter');
555  Hscatter_transverse = CreateH.Read('Hscatter_transverse');
556 
557  % apply custom potential in the scattering center
558  if ~isempty(scatterPotential)
559  if iscell( scatterPotential )
560  for idx = 1:length( scatterPotential )
561  obj.ApplyPotentialInScatter( CreateH, scatterPotential{idx} );
562  end
563  else
564  obj.ApplyPotentialInScatter( CreateH, scatterPotential);
565  end
566  Hscatter = CreateH.Read('Hscatter');
567  end
568 
569  non_singular_sites_Ginv = CreateH.Read('kulso_szabfokok');
570 
571  % apply the periodic boundary condition in the transverse direction
572  q = CreateH.Read('q');
573  if ~isempty( q ) && ~CreateH.Read('HamiltoniansDecimated')
574  Hscatter = Hscatter + Hscatter_transverse*diag(exp(1i*q)) + Hscatter_transverse'*diag(exp(-1i*q));
575  end
576 
577 
578  % reordering the sites of the scattering region to group the external point into the top corner of the matrix
579  obj.display(['EQuUs:Utils:', class(obj), ':CalcFiniteGreensFunctionFromHamiltonian: Calculating the surface Green function of the scattering region.'])
580  Hscatter = (sparse(1:size(Hscatter,1), 1:size(Hscatter,1), obj.E, size(Hscatter,1),size(Hscatter,1)) - Hscatter);
581  indexes = false( size(Hscatter,1), 1);
582  indexes(non_singular_sites_Ginv) = true;
583  Hscatter = [ Hscatter( ~indexes, ~indexes) , Hscatter(~indexes, indexes);
584  Hscatter(indexes, ~indexes), Hscatter(indexes, indexes)];
585 
586  obj.G = obj.partialInv( Hscatter, length(non_singular_sites_Ginv) );
587 
588  % recreate the Hamiltonian of the unit cell
589  obj.CreateRibbon();
590 
591 %disp( inv(obj.G) )
592  % gauge transformation of the vector potential in the effective Hamiltonians
593  if gauge_trans
594  try
595  % gauge transformation on Green's function
596  if ~isempty(obj.G) && ~isempty(obj.PeierlsTransform_Scatter) && isempty(obj.q) && ~isempty(obj.gauge_field)
597  surface_coordinates = obj.getCoordinates();
598  %> gauge transformation on the inverse Green's function
599  obj.G = obj.PeierlsTransform_Scatter.gaugeTransformation( obj.G, surface_coordinates, obj.gauge_field );
600  end
601  catch errCause
602  err = MException(['EQuUs:Utils:', class(obj), ':CalcFiniteGreensFunctionFromHamiltonian', 'Unable to perform gauge transformation']);
603  err = addCause(err, errCause);
604  save('Error_Ribbon_CalcFiniteGreensFunctionFromHamiltonian.mat')
605  throw(err);
606  end
607  end
608 
609  if onlyGinv
610  rcond_G = rcond(obj.G);
611  if isnan(rcond_G ) || abs( rcond_G ) < 1e-15
612  obj.display( ['EQuUs:Utils:', class(obj), ':CalcFiniteGreensFunctionFromHamiltonian: Regularizing Ginv by SVD'],1);
613  obj.Ginv = obj.inv_SVD( obj.G );
614  else
615  obj.Ginv = inv(obj.G);
616  end
617  obj.G = [];
618  else
619  obj.Ginv = [];
620  end
621  end
622 
623 %% CreateRibbon
624 %> @brief Creates the Hamiltonians of the unit cell of the ribbon shaped scattering region.
625  function CreateRibbon( obj )
626 
627  % create the Hamiltonians for a unit cell of the scattering region
628  if ~obj.Scatter_UC.Read( 'HamiltoniansCreated' )
629  obj.display( ['EQuUs:Utils:', class(obj), ':CreateRibbon: Creating ribbon Hamiltonian.'] )
630  obj.Scatter_UC.CreateHamiltonians( 'toSave', 0);
631  obj.Scatter_UC.ShiftCoordinates( obj.shift );
632 
633  %applying transverse potential
634  obj.ApplyTransversePotential( obj.Scatter_UC )
635  end
636 
637 
638 
639  % apply magnetic field in the unit cell of the scattering region
640  if ~isempty( obj.PeierlsTransform_Scatter ) && ~obj.Scatter_UC.Read('MagneticFieldApplied') && obj.Opt.magnetic_field_trans_invariant
641  obj.display(['EQuUs:Utils:', class(obj), ':CreateRibbon: Applying magnetic field in ribbon Hamiltonians'])
642  obj.PeierlsTransform_Scatter.PeierlsTransformLeads( obj.Scatter_UC );
643  end
644 
645 
646 
647  end
648 
649 
650 %% CreateInterface
651 %> @brief Creates the Hamiltonians for the interface regions between the leads and scattering center.
652 %> @param idx Identification number of the interface region.
653 %> @param varargin Cell array of optional parameters (https://www.mathworks.com/help/matlab/ref/varargin.html):
654 %> @param 'UseHamiltonian' Logical value. Set true if the interface region should be created to match to the whole Hamiltonian of the scattering center, false (default) if only the surface Green operator of the scattering center is used in the calculations.
655  function CreateInterface( obj, idx, varargin )
656 
657  p = inputParser;
658  p.addParameter('UseHamiltonian', false); %true if the interface region is matched to the whole Hamiltonian of the scattering center, false if the surface Green operator of the scattering center is used in the calculations.
659  p.parse(varargin{:});
660  UseHamiltonian = p.Results.UseHamiltonian;
661 
662  %> Hamiltoninans of the interface region
663  Interface_Region = obj.Interface_Regions{idx};
664 
665  % The regularization of the interface is performed according to the Leads
666  Leads = obj.FL_handles.Read( 'Leads' );
667  Lead = Leads{idx};
668 
669  Interface_Region.Write( 'coordinates2', obj.getCoordinates() );
670  Interface_Region.Write( 'K0', Lead.Read('K0'));
671  Interface_Region.Write( 'K1', Lead.Read('K1'));
672  Interface_Region.Write( 'K1adj', Lead.Read('K1adj'));
673  Interface_Region.Write( 'K1_transverse', Lead.Read('K1_transverse'));
674  Interface_Region.Write( 'K1_skew_left', Lead.Read('K1_skew_left'));
675  Interface_Region.Write( 'K1_skew_right', Lead.Read('K1_skew_right'));
676  Interface_Region.Write( 'coordinates', Lead.Read('coordinates'));
677  Interface_Region.Write( 'kulso_szabfokok', Lead.Read('kulso_szabfokok'));
678  Interface_Region.Write( 'OverlapApplied', true);
679 
680  coordinates_shift = [1, -1 ]; %relative to the leads
681  Interface_Region.ShiftCoordinates( coordinates_shift(idx) );
682 
683  %> coupling between the interface and the scattering region
684  Surface_sc = obj.createSurface_sc( idx );
685  Surface_sc.ApplyOverlapMatrices(obj.E);
686 
687 
688  Lead_Orientation = Surface_sc.Read('Lead_Orientation');
689  [~, K1, K1adj] = Surface_sc.qDependentHamiltonians();
690 
691  %K1 = Surface_sc.Read('K1');
692  %K1adj = Surface_sc.Read('K1adj');
693  [rows, cols] = find( K1 );
694  rows = unique(rows);
695  cols = unique(cols); % identical as non_singular_sites
696  if Lead_Orientation == 1
697 
698  if UseHamiltonian
699  Hscatter = obj.CreateH.Read('Hscatter');
700  if isempty( Hscatter)
701  error(['EQuUs:utils:', class(obj), ':CreateInterface: Hamiltonian of the scattering region needs to be constructed first.'])
702  end
703  Neff = size(Hscatter,1)-size(K1,2);
704  if ~Lead.Read('is_SVD_transformed')
705  % simple decimation
706  Kcoupling = [K1, sparse([], [], [], size(K1,1), Neff )];
707  Kcouplingadj = [K1adj; sparse([], [], [], Neff, size(K1,1))];
708  else
709  % regularization with SVD
710  Kcoupling = [K1, sparse([], [], [], size(K1,1), Neff )];
711  Kcouplingadj = [K1adj; sparse([], [], [], Neff, size(K1,1))];
712  end
713  else
714  Neff = length(rows); %number of coupling sites at Lead_Oriantation=-1
715  if ~Lead.Read('is_SVD_transformed')
716  % simple decimation
717  Kcoupling = [K1(:,cols), zeros( size(K1,1), Neff )];
718  Kcouplingadj = [K1adj(cols,:); zeros(Neff, size(K1,1))];
719  else
720  % regularization with SVD
721  Kcoupling = [K1(:,cols), zeros(size(K1,1), Neff )];
722  Kcouplingadj = [K1adj(cols,:); zeros(Neff, size(K1,1))];
723  end
724  end
725 
726 
727 
728  elseif Lead_Orientation == -1
729 
730  if UseHamiltonian
731  Hscatter = obj.CreateH.Read('Hscatter');
732  if isempty( Hscatter)
733  error('EQuUs:utils:Ribbon:CreateInterface: Hamiltonian of the scattering region needs to be constructed first.')
734  end
735  Neff = size(Hscatter,1)-size(K1adj,2);
736  if ~Lead.Read('is_SVD_transformed')
737  % simple decimation
738  Kcoupling = [sparse([], [], [], length(cols), Neff ), K1adj(cols,:)];
739  Kcouplingadj = [sparse([], [], [], Neff, length(cols)); K1(:,cols)];
740  else
741  % regularization with SVD
742  Kcoupling = [sparse([], [], [], size(K1adj,1), Neff ), K1adj];
743  Kcouplingadj = [sparse([], [], [], Neff, size(K1,1) ); K1];
744  end
745  else
746  Neff = length(cols); %number of coupling sites at Lead_Oriantation=1
747  if ~Lead.Read('is_SVD_transformed')
748  % simple decimation
749  Kcoupling = [zeros(length(cols), Neff), K1adj(cols,rows)];
750  Kcouplingadj = [zeros(Neff, length(cols)); K1(rows,cols)];
751  else
752  % regularization with SVD
753  Kcoupling = [zeros(size(K1,1), Neff), K1adj(:,rows)];
754  Kcouplingadj = [zeros(Neff, size(K1,2)); K1(rows,:)];
755  end
756  end
757 
758 
759  else
760  error('EQuUs:Utils:Ribbon:CreateInterface', 'Unknown lead orientation');
761  end
762 
763 
764 
765  Interface_Region.Write('Kcoupling', Kcoupling);
766  Interface_Region.Write('Kcouplingadj', Kcouplingadj);
767 
768  % method to adjust the interface region and coupling to the scattering region by an external function.
769  if ~isempty( obj.interfacemodel )
770  obj.interfacemodel( Interface_Region );
771  end
772 
773  Interface_Region.Calc_Effective_Hamiltonians( 0, 'Lead', Lead );
774 
775  end
776 
777 %% CreateClone
778 %> @brief Creates a clone of the present object.
779 %> @return Returns with the cloned object.
780  function ret = CreateClone( obj )
781 
782  ret = Ribbon( 'width', obj.width, ...
783  'height', obj.height, ...
784  'filenameIn', obj.filenameIn, ...
785  'filenameOut', obj.filenameOut, ...
786  'E', obj.E, ...
787  'EF', 0, ...
788  'phi', obj.phi, ...
789  'silent', obj.silent, ...
790  'transversepotential', obj.transversepotential, ...
791  'Opt', obj.Opt, ...
792  'param', obj.param, ...
793  'q', obj.q, ...
794  'leadmodel', obj.leadmodel, ...
795  'interfacemodel', obj.interfacemodel);
796 
797  ret.EF = obj.EF;
798  ret.CreateH = obj.CreateH.CreateClone();
799  ret.FL_handles = obj.FL_handles.CreateClone();
800  ret.Scatter_UC = obj.Scatter_UC.CreateClone();
801  ret.Interface_Regions = cell(size(obj.Interface_Regions));
802  for idx = 1:length(obj.Interface_Regions)
803  ret.Interface_Regions{idx} = obj.Interface_Regions{idx}.CreateClone();
804  end
805  if ~isempty( obj.PeierlsTransform_Scatter )
806  ret.PeierlsTransform_Scatter = obj.PeierlsTransform_Scatter.CreateClone();
807  end
808 
809  if ~isempty( obj.PeierlsTransform_Leads )
810  ret.PeierlsTransform_Leads = obj.PeierlsTransform_Leads.CreateClone();
811  end
812  ret.gauge_field = obj.gauge_field; % function handle for the scalar field to transform the vector potential from Landauy to Landaux
813 
814  end
815 
816 
817 end % methods public
818 
819 methods ( Access = protected )
820 
821 %% setFermiEnergy
822 %> @brief Sets the Fermi energy on the atomic sites for the calculations (use the same units as the elements of the Hamiltonian).
823  function setFermiEnergy( obj )
824  if ~isempty(obj.EF)
825  obj.param.scatter.epsilon = obj.param.scatter.epsilon - obj.EF;
826  for idx = 1:length(obj.param.Leads)
827  obj.param.Leads{idx}.epsilon = obj.param.Leads{idx}.epsilon - obj.EF;
828  end
829  end
830  end
831 
832 %% createSurface_sc
833 %> @brief Creates the copuling Hamiltonians between the scattering and interface region
834 %> @param idx The identification number of the interface region. (Integer value.)
835 %> @return An instance of class #Lead describing the copuling between the scattering and interface region
836  function Surface_sc = createSurface_sc( obj, idx )
837  Surface_sc = obj.Scatter_UC.CreateClone('empty', true);
838 
839  if ~isempty( obj.param.Leads{idx}.vargamma_sc )
840  params = Surface_sc.Read( 'params' );
841  params.vargamma = obj.param.Leads{idx}.vargamma_sc;
842  Surface_sc.Write( 'params', params );
843  end
844 
845  Surface_sc.CreateHamiltonians( 'toSave', 0);
846  if idx == 1
847  coordinates_shift = 0 + obj.shift ;
848  elseif idx == 2
849  coordinates_shift = obj.height-1 + obj.shift;
850  end
851  Surface_sc.ShiftCoordinates( coordinates_shift )
852  Surface_sc.Write('Hanyadik_Lead', idx);
853  Surface_sc.Write('Lead_Orientation', obj.Interface_Regions{idx}.Read('Lead_Orientation'));
854 
855  %> applying transverse potential
856  obj.ApplyTransversePotential( Surface_sc )
857 
858  %> applying magnetic field
859  if ~isempty( obj.PeierlsTransform_Leads )
860  %> In superconducting lead one must not include nonzero magnetic
861  %> field.
862  %> Hamiltonians in transverse computations must remain
863  %> traslational invariant.
864  if ~Surface_sc.isSuperconducting()
865  obj.display('EQuUs:Utils:Ribbon:createSurface_sc: Applying magnetic field in the Hamiltonians')
866  obj.PeierlsTransform_Leads.PeierlsTransformLeads( Surface_sc );
867  else
868  obj.display('EQuUs:Utils:Ribbon:createSurface_sc: Applying gauge transformation in the Hamiltonians')
869  obj.PeierlsTransform_Leads.gaugeTransformationOnLead( Surface_sc, obj.gauge_field );
870  end
871  end
872 
873 
874  end
875 
876 %% ApplyTransversePotential
877 %> @brief Apply the tranvesre potential in the Hamiltonians
878 %> @param Scatter_UC An instance of class #Lead containing the Hamiltonians.
879  function ApplyTransversePotential( obj, Scatter_UC )
880  if ~isempty(obj.transversepotential) && isempty(obj.q) %In transverse computations no transverse potential can be applied
881  coordinates = Scatter_UC.Read('coordinates');
882  if nargin( obj.transversepotential ) == 1
883  potential2apply = obj.transversepotential( coordinates );
884  elseif nargin( obj.transversepotential ) == 2
885  potential2apply = obj.transversepotential( Scatter_UC, obj.E );
886  else
887  error('EQuUs:Utils:Ribbon:ApplyTransversePotential', 'To many input arguments in function handle scatterpotential');
888  end
889 
890  if isprop(coordinates, 'BdG_u')
891  if size( potential2apply, 1) == 1 || size( potential2apply, 2) == 1
892  potential2apply(~coordinates.BdG_u) = -potential2apply(~coordinates.BdG_u);
893  else
894  potential2apply(~coordinates.BdG_u, ~coordinates.BdG_u) = -conj(potential2apply(~coordinates.BdG_u, ~coordinates.BdG_u));
895  end
896  end
897  Scatter_UC.AddPotential( potential2apply );
898  end
899  end
900 
901 
902 
903 
904 %% CreateHandles
905 %> @brief Initializes the attributes of the class.
906  function CreateHandles( obj )
907 
908  CreateHandles@NTerminal( obj )
909 
910  obj.Scatter_UC = obj.FL_handles.SurfaceGreenFunctionCalculator([], 'createCore', 1, 'q', obj.q);
911 
912  end
913 
914 %% calculate_lead_attach_points
915 %> @brief Determines the site indexes at which the leads are connected to the scattering center.
916  function calculate_lead_attach_points( obj )
917  for idx = 1:length(obj.param.Leads)
918  obj.param.Leads{idx}.M = obj.param.scatter.shape.width;
919  end
920  end
921 
922 %% createShape
923 %> @brief Creates the geometry data of the ribbon shaped scattering region.
924  function createShape( obj )
925 
926  if ~isempty( obj.width ) && ~isempty( obj.height )
927  obj.param.scatter.shape.width = obj.width;
928  obj.param.scatter.shape.height = obj.height;
929  end
930 
931  if ~isempty(obj.param.scatter.shape.width)
932  obj.width = obj.param.scatter.shape.width;
933  else
934  err = MException(['EQuUs:utils:', class(obj), ':createShape'], 'Shape is not given correctly, width is missing');
935  throw(err)
936  end
937 
938  if ~isempty(obj.param.scatter.shape.height)
939  obj.height = obj.param.scatter.shape.height;
940  else
941  err = MException(['EQuUs:utils:', class(obj), ':createShape'], 'Shape is not given correctly, height is missing');
942  throw(err)
943  end
944 
945  obj.calculate_lead_attach_points();
946 
947  end
948 
949 end % protected methods
950 
951 
952 methods (Access=protected)
953 
954 
955 %% InputParsing
956 %> @brief Parses the optional parameters for the class constructor.
957 %> @param varargin Cell array of optional parameters (https://www.mathworks.com/help/matlab/ref/varargin.html):
958 %> @param 'width' Integer. The number of the nonsingular atomic sites in the cross section of the ribbon.
959 %> @param 'height' Integer. The height of the ribbon in units of the lattice vector.
960 %> @param 'filenameIn' The input filename containing the computational parameters. (Use parameters 'Op' and 'param' instead)
961 %> @param 'filenameOut' The output filename to export the computational parameters.
962 %> @param 'WorkingDir' The absolute path to the working directoy.
963 %> @param 'E' The energy value used in the calculations (in the same units as the Hamiltonian).
964 %> @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)
965 %> @param 'silent' Set true to suppress output messages.
966 %> @param 'transversepotential' A function handle pot = f( #Coordinates ) or pot=f( #CreateLeadHamiltonians, Energy) of the transverse potential applied in the lead. (Instead of #CreateLeadHamiltonians can be used its derived class)
967 %> @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.
968 %> @param 'interfacemodel' A function handle f( #InterfaceRegion ) to manually adjus the interface regions. (Usefull when 'leadmodel' is also given. For example see @InterfaceModel)
969 %> @param 'Opt' An instance of the structure #Opt.
970 %> @param 'param' An instance of the structure #param.
971 %> @param 'q' The transverse momentum quantum number.
972  function InputParsing( obj, varargin )
973 
974  p = inputParser;
975  p.addParameter('width', obj.width);
976  p.addParameter('height', obj.height);
977  p.addParameter('filenameIn', obj.filenameIn, @ischar);
978  p.addParameter('filenameOut', obj.filenameOut, @ischar);
979  p.addParameter('WorkingDir', obj.WorkingDir, @ischar);
980  p.addParameter('E', obj.E, @isscalar);
981  p.addParameter('EF', obj.EF);
982  p.addParameter('silent', obj.silent);
983  p.addParameter('transversepotential', obj.transversepotential);
984  p.addParameter('leadmodel', obj.leadmodel); %individual physical model for the contacts
985  p.addParameter('interfacemodel', obj.interfacemodel); %individual physical model for the interface regions
986  p.addParameter('Opt', obj.Opt);
987  p.addParameter('param', obj.param);
988  p.addParameter('q', obj.q);
989 
990  p.parse(varargin{:});
991 
992  InputParsing@NTerminal( obj, 'filenameIn', p.Results.filenameIn, ...
993  'filenameOut', p.Results.filenameOut, ...
994  'WorkingDir', p.Results.WorkingDir, ...
995  'E', p.Results.E, ...
996  'EF', p.Results.EF, ...
997  'silent', p.Results.silent, ...
998  'leadmodel', p.Results.leadmodel, ...
999  'interfacemodel', p.Results.interfacemodel, ...
1000  'Opt', p.Results.Opt, ...
1001  'param', p.Results.param, ...
1002  'q', p.Results.q);
1003 
1004 
1005  obj.width = p.Results.width;
1006  obj.height = p.Results.height;
1007  obj.transversepotential = p.Results.transversepotential;
1008 
1009  end
1010 
1011 end % methdos protected
1012 
1013 
1014 
1015 end
A class describing an N-terminal geometry for equilibrium calculations mostly in the zero temperature...
Definition: NTerminal.m:38
function isSuperconducting()
Test, whether the lead is in the superconducting phase or not.
function InputParsing(varargin)
Parses the optional parameters for the class constructor.
lead_param Leads
A list of structures lead_param containing the physical parameters for the scattering region.
Definition: structures.m:49
Property q
The tranverse momentum for transverse computations.
function CreateClone(varargin)
Creates a clone of the present class.
Property params
An instance of the structure lead_param.
Structure Opt contains the basic computational parameters used in EQuUs.
Definition: structures.m:60
Structure open_channels describes the open channel in the individual leads.
Definition: structures.m:159
function Conduktance()
Calculates the conductance matrix from the scattering matrix.
Structure shape contains data about the geometry of the scattering region.
Definition: structures.m:106
function Write(MemberName, input)
Sets the value of an attribute in the class.
Class to create and store Hamiltonian of the translational invariant leads.
Property Interface_Regions
A cell array of classes InterfaceRegion to describe the interface region between the leads and the sc...
Definition: NTerminal.m:77
A class for calculations on a ribbon of finite width for equilibrium calculations mostly in the zero ...
Definition: Ribbon.m:34
Property H1
The coupling Hamiltonian between the unit cells.
function CalcFiniteGreensFunction(varargin)
Calculates the Green operator of the scattering region.
function Transport(Energy, B)
Calculates the conductance at a given energy value.
function createShape(theNode, param_structure)
This function creates structure shape for the scattering region.
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.
Property varargin
list of optional parameters (see http://www.mathworks.com/help/matlab/ref/varargin....
function createOutput(filename, Opt, param)
This function creates an output file containing the running parameters.
function Landauy(x, y, eta_B)
Vector potential in the Landau gauge parallel to the y direction.
function CreateHamiltonians(Opt, param, varargin)
Constructor of the class.
function CreateHamiltonians(varargin)
Creates the Hamiltonians H_0 and H_1 of the lead.
Property q
A transverse momentum.
kulso_szabfokok
An array containing the site indexes connected to other parts.
Definition: structures.m:193
function Reset()
Resets all elements 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...
Definition: Lead.m:29
Property varargin
list of optional parameters (see http://www.mathworks.com/help/matlab/ref/varargin....
function InterfaceModel(Interface_Region)
Method to adjust the Hamiltonians of the interface regions.
Structure param contains data structures describing the physical parameters of the scattering center ...
Definition: structures.m:45
Property Lead_Orientation
The orientation of the lead. Set +1 is the "incoming" direction of the propagating states is defined ...
scatter_param scatter
An instance of the structure scatter_param containing the physical parameters for the scattering regi...
Definition: structures.m:47
Structure sites contains data to identify the individual sites in a matrix.
Definition: structures.m:187
function ShiftCoordinates(shift)
Shifts the coordinates of the sites by an integer multiple of the lattice vector Coordinates....
function Read(MemberName)
Query for the value of an attribute in the class.
Property Opt
An instance of structure Opt.
Definition: Messages.m:31
function Landaux(x, y, eta_B, Aconst, height, lattice_constant)
Vector potential in the Landau gauge parallel to the x direction.
function CreateClone(varargin)
Creates a clone of the present Lead object.
function SurfaceGreenFunctionCalculator(idx, varargin)
Calculates the surface Green's function or the self energy of a Lead.
function gauge_field(x, y, eta_B, Aconst, height, lattice_constant)
Scalar gauge field connecting Landaux and Landauy gauges.
function Conductance2()
Conductance calculated by Eq (19) in PRB 73 085414.
Property E
The energy value for which the TrukkosSajatertekek eigenvalue problem was solved.
Structure containing the coordinates and other quantum number identifiers of the sites in the Hamilto...
Definition: Coordinates.m:24
function AddPotential(V)
Adds on-site potential to the Hamiltonian H0.
function Read(MemberName)
Query for the value of an attribute in the class.
Structure junction_sites contains data to identify the individual sites of the leads,...
Definition: structures.m:176
Property coordinates
An instance of the structure coordinates.
A class to create and store Hamiltonian of the scattering region.