Eötvös Quantum Utilities  v4.9.146
Providing the Horsepowers in the Quantum Realm
Transport_Interface_Keldysh.m
Go to the documentation of this file.
1 %% Eotvos Quantum Transport Utilities - Transport_Interface_Keldysh
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 basic Basic Functionalities
18 %> @{
20 %> @brief A class to evaluate the Dyson equation derived from class #Transport_Interface with specific modifications for steady state non-equilibrium processes.
21 %> @}
22 %> @brief A class to evaluate the Dyson equation derived from class #Transport_Interface with specific modifications for steady state non-equilibrium processes.
23 %%
25 
26 
27  properties ( Access = protected )
28 
29  end
30 
31 
32 
33 
34 methods
35 
36 %% Contructor of the class
37 %> @brief Constructor of the class.
38 %> @param Opt An instance of the structure #Opt.
39 %> @param param An instance of the structure #param.
40 %> @param varargin Cell array of optional parameters. For details see #InputParsing.
41 %> @return An instance of the class
42 function obj = Transport_Interface_Keldysh(E, Opt, param, varargin)
43  obj = obj@Transport_Interface( E, Opt, param, varargin{:} );
44 
45  obj.Initialize();
46 
47 
48 end
49 
50 
51 %% LeadCalc
52 %> @brief Invokes the function #SurfaceGreenFunctionCalculator for each lead.
53 %> @param varargin Cell array of optional parameters (https://www.mathworks.com/help/matlab/ref/varargin.html):
54 %> @param 'shiftLeads' A real number. If given, the on-site potentials of the sites in the lead are shifted by this value.
55 %> @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)
56 %> @param 'coordinates_shift' An integer. If given, the coordinates of the sites in the lead are shifted by coordinates_shift*lattice vector.
57 %> @param 'gauge_field' Function handle S = f( x,y) of the gauge transformation on the Hamiltonians of the leads. (S is a N x 1 vector, where N is the number of the points given by the x and y coordinates.)
58 %> @param 'createCore' Set 1 for creating the class instances #Lead but omitting any further calculations, or false (default) to continue with the calculations.
59 %> @param 'SelfEnergy' Logical value. Set true to calculate the self energies of the leads. (default is false)
60 %> @param 'SurfaceGreensFunction' Logical value. Set true to calculate the surface Green functions of the leads. (default is true)
61 %> @param 'leads' Array of lead identification numbers for which the calculations should be performed. (Default is a list of all leads)
62 %> @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.
63 %> @param 'CustomHamiltonian' An instance of class #Custom_Hamiltonians or its derived class to create custom Hamiltonians in the leads.
64 %> @param 'Just_Create_Hamiltonians' Logical value. Set true to create the Hamiltonians of the leads, but stop further calculations. Default value is false.
65 %> @param 'q' The tranverse momentum for transverse computations.
66 %> @param 'T' The absolute temperature in Kelvins.
67 %> @param 'mu_leads' An array containing the chemical potentials of the leads in the same unit as other energy scales in the Hamiltonians.
68 %> @return Returns with an array of the created #Lead instances.
69 function lead_ret = LeadCalc( obj, varargin )
70 
71  p = inputParser;
72  p.addParameter('shiftLeads', zeros(length(obj.param.Leads),1));
73  p.addParameter('transversepotential', []);
74  p.addParameter('coordinates_shift', zeros(length(obj.param.Leads),1) );
75  p.addParameter('gauge_field', [] );
76  p.addParameter('createCore', 0);
77  p.addParameter('SelfEnergy',false); % set true to calculate the self energy of the semi-infinite lead
78  p.addParameter('SurfaceGreensFunction', true );% set true to calculate the surface Greens function of the semi-infinite lead
79  p.addParameter('leads', 1:length(obj.param.Leads) );% list of leads to be calculated
80  p.addParameter('leadmodel', []); %function handle for an individual physical model for the contacts
81  p.addParameter('CustomHamiltonian', []);
82  p.addParameter('Just_Create_Hamiltonians', 0);
83  p.addParameter('q', [] ) %transverse momentum
84  p.addParameter('T', 0 ) %transverse momentum
85  p.addParameter('mu_leads', zeros( size(obj.param.Leads) ) ) %chemical potentials in the leads
86  p.addParameter('bias_leads', zeros( size(obj.param.Leads) ) ) %array containign the bias of the leads
87  p.parse(varargin{:});
88 
89  shiftLeads = p.Results.shiftLeads;
90  transversepotential = p.Results.transversepotential;
91  coordinates_shift = p.Results.coordinates_shift; % shifts coordinates by "coordinates_shift" times of a lattice vector.
92  gauge_field = p.Results.gauge_field; % if a valid function handle is given, it is used to perform gauge transformation on Lead Hamiltonians
93  createCore = p.Results.createCore;
94  SelfEnergy = p.Results.SelfEnergy;
95  SurfaceGreensFunction = p.Results.SurfaceGreensFunction;
96  leads = p.Results.leads;
97  Just_Create_Hamiltonians = p.Results.Just_Create_Hamiltonians;
98  q = p.Results.q;
99  CustomHamiltonian = p.Results.CustomHamiltonian;
100  leadmodel = p.Results.leadmodel;
101  T = p.Results.T;
102  mu_leads = p.Results.mu_leads;
103  bias_leads = p.Results.bias_leads;
104 
105  if isempty(obj.Leads)
106  obj.Leads = cell(length(obj.param.Leads),1);
107  end
108  obj.display(' ')
109  obj.display('Creating surface_Greens_function interfaces:')
110  obj.display(' ')
111  idx = 1;
112  while idx <= length(leads)
113  obj.display('------------------------------')
114  obj.display(['Lead: ',num2str(leads(idx))])
115  obj.Leads{leads(idx)} = obj.SurfaceGreenFunctionCalculator(leads(idx), 'shiftLead', shiftLeads(leads(idx)), 'transversepotential', transversepotential, 'coordinates_shift', coordinates_shift(leads(idx)),...
116  'gauge_field', gauge_field, 'createCore', createCore, ...
117  'SelfEnergy', SelfEnergy, 'SurfaceGreensFunction', SurfaceGreensFunction, 'q', q, 'leadmodel', leadmodel, 'CustomHamiltonian', CustomHamiltonian, 'Just_Create_Hamiltonians', Just_Create_Hamiltonians, ...
118  'T', T, 'mu', mu_leads(idx), 'bias', bias_leads(idx) );
119  obj.display(' ')
120  idx = idx + 1;
121  end
122 
123  lead_ret = obj.Leads;
124 
125 end
126 
127 %% setTemperature
128 %> @brief Sets the temperature in the calculations
129 %> @param T The value of the temperature in Kelvins.
130  function setTemperature( obj, T )
131 
132  % set the temperature in the leads
133  if ~isempty( obj.Leads )
134  for idx = 1:length( obj.Leads )
135  obj.Leads{idx}.setTemperature(T);
136  end
137  end
138 
139 
140  end
141 
142 %% SurfaceGreenFunctionCalculator
143 %> @brief Calculates the surface Green's function or the self energy of a Lead.
144 %> @param idx the identification number of the lead in the attribute Leads (see attribute #CreateLeadHamiltonians.Hanyadik_Lead).
145 %> @param varargin Cell array of optional parameters (https://www.mathworks.com/help/matlab/ref/varargin.html):
146 %> @param 'shiftLeads' A real number. If given, the on-site potentials of the sites in the lead are shifted by this value.
147 %> @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)
148 %> @param 'Lead' An instance of class #Lead_Keldysh (or its subclass) that is intended to be used in the calculations.
149 %> @param 'coordinates_shift' An integer. If given, the coordinates of the sites in the lead are shifted by coordinates_shift*lattice vector.
150 %> @param 'gauge_field' Function handle S = f( x,y) of the gauge transformation on the Hamiltonians of the leads. (S is a N x 1 vector, where N is the number of the points given by the x and y coordinates.)
151 %> @param 'createCore' Set 1 for creating the class instances #Lead but omitting any further calculations, or false (default) to continue with the calculations.
152 %> @param 'SelfEnergy' Logical value. Set true to calculate the self energies of the leads. (default is true)
153 %> @param 'SurfaceGreensFunction' Logical value. Set true to calculate the surface Green functions of the leads. (default is false)
154 %> @param 'leads' Array of lead identification numbers for which the calculations should be performed. (Default is a list of all leads)
155 %> @param 'leadmodel' A function handle #Lead_keldysh=f( idx, E, varargin ) of the alternative lead model with equivalent inputs and return values as #Transport_Interface_Keldysh.SurfaceGreenFunctionCalculator and with E standing for the energy.
156 %> @param 'CustomHamiltonian' An instance of class #Custom_Hamiltonians or its derived class to create custom Hamiltonians in the leads.
157 %> @param 'Just_Create_Hamiltonians' Logical value. Set true to create the Hamiltonians of the leads, but stop further calculations. Default value is false.
158 %> @param 'q' The tranverse momentum for transverse computations.
159 %> @param 'T' The absolute temperature in Kelvins.
160 %> @param 'mu' The Chemical potential in the same unit as other energy scales in the Hamiltonians.
161 %> @return Returns with the created #Lead_Keldysh instance.
162 function Lead_ret = SurfaceGreenFunctionCalculator( obj, idx, varargin)
163 
164  p = inputParser;
165  p.addParameter('createCore', 0);
166  p.addParameter('Just_Create_Hamiltonians', 0);
167  p.addParameter('shiftLead', 0);
168  p.addParameter('coordinates_shift', 0);
169  p.addParameter('transversepotential',[]);
170  p.addParameter('Lead',[]);
171  p.addParameter('gauge_field', [] );% gauge field for performing gauge transformation
172  p.addParameter('SelfEnergy',true); % set true to calculate the self energy of the semi-infinite lead
173  p.addParameter('SurfaceGreensFunction', false );% set true to calculate the surface Greens function of the semi-infinite lead
174  p.addParameter('leadmodel', []); %function handle for an individual physical model for the contacts
175  p.addParameter('CustomHamiltonian', []);
176  p.addParameter('q', [] ) %transverse momentum
177  p.addParameter('T', 0 ) %transverse momentum
178  p.addParameter('mu', 0 ) %chemical potential in the lead
179  p.addParameter('bias', 0 ) %electrical bias in the lead
180  p.parse(varargin{:});
181 
182  createCore = p.Results.createCore;
183  Lead_ret = p.Results.Lead;
184  q = p.Results.q;
185  SurfaceGreensFunction = p.Results.SurfaceGreensFunction;
186  SelfEnergy = p.Results.SelfEnergy;
187  mu = p.Results.mu;
188  bias = p.Results.bias;
189  T = p.Results.T;
190 
191 
192  if ~isempty( Lead_ret )
193  % do not create a new class instance if Lead_ret was given
194  Lead_ret.Reset();
195 
196  elseif ~isempty(idx) && ~isempty(obj.param.Leads{idx}.Lead_Orientation)
197  Lead_ret = Lead_Keldysh(obj.Opt, obj.param,...
198  'Hanyadik_Lead', idx,....
199  'Lead_Orientation', obj.param.Leads{idx}.Lead_Orientation, ...
200  'q', q, ...
201  'mu', mu, ...
202  'bias', bias, ...
203  'T', T);
204  else
205  Lead_ret = Lead_Keldysh(obj.Opt, obj.param,...
206  'Hanyadik_Lead', idx, ...
207  'mu', mu, ...
208  'bias', bias, ...
209  'T', T, ...
210  'q', q);
211  end
212 
213  % calculate the retarded surface Green operator/self-energy
214  SurfaceGreenFunctionCalculator@Transport_Interface( obj, idx, 'createCore', createCore, ...
215  'Just_Create_Hamiltonians', p.Results.Just_Create_Hamiltonians, ...
216  'shiftLead', p.Results.shiftLead, ...
217  'coordinates_shift', p.Results.coordinates_shift, ...
218  'transversepotential', p.Results.transversepotential, ...
219  'Lead', Lead_ret, ...
220  'gauge_field', p.Results.gauge_field, ...
221  'SelfEnergy', SelfEnergy, ...
222  'SurfaceGreensFunction', SurfaceGreensFunction, ...
223  'q', p.Results.q, ...
224  'CustomHamiltonian', p.Results.CustomHamiltonian, ...
225  'leadmodel', p.Results.leadmodel);
226 
227  if createCore
228  return
229  end
230 
231 
232  % lesser surface Green operator
233  if SurfaceGreensFunction
234  Lead_ret.LesserSurfaceGreenFunction();
235  end
236 
237  % lesser SelfEnergy
238  if SelfEnergy
239  Lead_ret.LesserSelfEnergy();
240  end
241 
242 
243 end
244 
245 
246 
247 %% replaceLead
248 %> @brief Adds/replaces a lead to/in the system.
249 %> @param Lead_tmp An instance of class #Lead_Keldysh or its subclass
250 %> @param idx Identification number of the lead. (see attribute #CreateLeadHamiltonians.Hanyadik_Lead).
251  function replaceLead( obj, Lead_tmp, idx )
252 
253  if ~strcmpi( class(obj.junction), 'Lead_Keldysh' )
254  supClasses = superclasses(Lead_tmp);
255  if sum( strcmp( supClasses, 'Lead_Keldysh' ) ) == 0
256  error(['EQuUs:', class(obj), ':replaceLead'], 'Input Lead_tmp is not valid.');
257  end
258  end
259 
260  obj.Leads{ idx } = Lead_tmp;
261  obj.Opt.NofLeads = length(obj.Leads);
262  end
263 
264 
265 %% CreateClone
266 %> @brief Creates a clone of the present class.
267 %> @return Returns with the cloned object.
268  function ret = CreateClone( obj )
269 
270  if isempty( obj.CreateH )
271  CreateH_loc = [];
272  else
273  CreateH_loc = obj.CreateH.CreateClone();
274  end
275 
276  if isempty( obj.PeierlsTransform )
277  PeierlsTransform_loc = [];
278  else
279  PeierlsTransform_loc = obj.PeierlsTransform.CreateClone();
280  end
281 
282  ret = Transport_Interface_Keldysh(obj.E, obj.Opt, obj.param, ...
283  'CreateH', CreateH_loc, ...
284  'PeierlsTransform', PeierlsTransform_loc );
285 
286  for idx = 1:length(obj.Leads)
287  ret.replaceLead( obj.Leads{idx}.CreateClone(), idx );
288  end
289 
290  end
291 
292 %% Write
293 %> @brief Sets the value of an attribute in the class.
294 %> @param MemberName The name of the attribute to be set.
295 %> @param input The value to be set.
296  function Write(obj, MemberName, input)
297 
298  if strcmp(MemberName, 'param')
299  obj.param = input;
300  obj.Reset();
301  return
302  elseif strcmp(MemberName, 'E')
303  obj.setEnergy( input );
304  elseif strcmp(MemberName, 'T')
305  obj.setTemperature( input );
306  else
307  try
308  obj.(MemberName) = input;
309  catch
310  error(['EQuUs:', class(obj), ':Read'], ['No property to write with name: ', MemberName]);
311  end
312  end
313 
314  end
315 
316 %% Reset
317 %> @brief Resets all elements in the class.
318  function Reset(obj)
319 
320  if strcmpi( class(obj), 'Transport_Interface_Keldysh' )
321  meta_data = metaclass(obj);
322 
323  for idx = 1:length(meta_data.PropertyList)
324  prop_name = meta_data.PropertyList(idx).Name;
325  if strcmp(prop_name, 'Opt') || strcmp( prop_name, 'param') || strcmp(prop_name, 'varargin') || strcmp(prop_name, 'E')
326  continue
327  end
328  obj.Clear( prop_name );
329  end
330  end
331 
332  obj.Initialize();
333 
334 
335  end
336 
337 
338 end % method public
339 
340 
341 methods ( Access = protected )
342 
343 %% Initialize
344 %> @brief Initializes object attributes.
345  function Initialize(obj)
346 
347  Initialize@Transport_Interface(obj);
348 
349  if strcmpi( class(obj), 'Transport_Interface_Keldysh')
350  obj.InputParsing( obj.varargin{:});
351  end
352 
353  end
354 
355 %% InputParsing
356 %> @brief Parses the optional parameters for the class constructor.
357 %> @param varargin Cell array of optional parameters (https://www.mathworks.com/help/matlab/ref/varargin.html):
358 %> @param 'CreateH' An instance of class #CreateHamiltonians storing the Hamiltonian of the scattering center.
359 %> @param 'PeierlsTransform' An instance of class #Peierls for the Peierls transformations.
360  function InputParsing( obj, varargin )
361 
362  p = inputParser;
363  p.addParameter('CreateH', []);
364  p.addParameter('PeierlsTransform', []);
365 
366  p.parse(varargin{:});
367 
368  InputParsing@Transport_Interface( obj, 'PeierlsTransform', p.Results.PeierlsTransform, ...
369  'CreateH', p.Results.CreateH);
370 
371  end
372 
373 end % methods protected
374 
375 end
Structure Opt contains the basic computational parameters used in EQuUs.
Definition: structures.m:60
Class to create and store Hamiltonian of the translational invariant leads.
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.
A class to evaluate the Dyson equation derived from class Transport_Interface with specific modificat...
A class to import custom Hamiltonians provided by other codes or created manually
A class to calculate the Green functions and self energies of a translational invariant lead The nota...
Definition: Lead.m:29
function Initialize()
Initializes object attributes.
A class to calculate the Keldysh lesser and greater Green functions and self energies of a translatio...
Definition: Lead_Keldysh.m:29
Structure param contains data structures describing the physical parameters of the scattering center ...
Definition: structures.m:45
Structure sites contains data to identify the individual sites in a matrix.
Definition: structures.m:187
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.