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 representing a two-terminal structure defined on a preprogrammed lattices
for steady state non-equilibrium calculations.
21 %> @image html Ribbon_structure.jpg
22 %> @image latex Ribbon_structure.jpg
24 %> @brief A
class representing a two-terminal structure defined on a preprogrammed lattices for steady state non-equilibrium calculations.
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 Ribbon_structure.jpg
29 %> @image latex Ribbon_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 %> Each rectangle describes a unit cell including singular and non-singular
sites.
32 %> The scattering center is also described by a set of identical unit cells, but arbitrary potential can be used.
34 %> 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.
39 properties (Access =
protected)
42 properties (Access =
public)
47 methods ( Access = public )
50 %% Contructor of the class
51 %> @brief Constructor of the class.
52 %> @
param varargin Cell array of optional parameters. For details see
#InputParsing. 53 %> @
return An instance of the
class 58 if strcmpi(
class(obj),
'Ribbon_Keldysh')
59 % processing the optional parameters
60 obj.InputParsing(varargin{:})
63 obj.display([
'EQuUs:Utils:',
class(obj),
':Ribbon_Keldysh: Creating a Ribbon_Keldysh object'])
65 % create the
shape of the scattering region
71 % exporting calculation parameters into an XML format
74 % create
class intances and initializing class attributes
78 obj.CreateRibbon(
'justHamiltonians',
true);
86 %> @brief Custom Dyson
function for a two terminal arrangement on a two dimensional lattice.
87 %> @
param varargin Cell array of optional parameters (https:
88 %> @
param 'gfininv' The inverse of the Greens
function of the scattering region. For
default the inverse of the attribute #G is used.
89 %> @
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.
90 %> @
param 'onlyGinverz' Logical value. Set
true to calculate only the inverse of the total Green
operator, or false (
default) to calculate #G as well.
91 %> @
param 'recalculateSurface' A vector of the identification numbers of the lead surfaces to be recalculated.
92 %> @
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.
93 %> @
param 'kulso_szabfokok' Array of
sites to be kept after the decimation procedure. (Use parameter
'keep_sites' instead)
94 %> @
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.
95 %> @
param 'keep_sites' Name of
sites to be kept in the resulted Green
function (Possible values are:
'scatter',
'interface',
'lead').
96 %> @
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.
97 %> @
return [1] The calculated Greens
function.
98 %> @
return [2] The inverse of the Green
operator.
99 %> @
return [3] An instance of structure #
junction_sites describing the
sites in the calculated Green
operator.
100 function [Gret, Ginverz,
junction_sites] = CustomDysonFunc( obj, varargin ) %NEW output
103 p.addParameter(
'gfininv', []);
104 p.addParameter(
'constant_channels',
true);
105 p.addParameter(
'onlyGinverz',
false );
106 p.addParameter(
'recalculateSurface', [1 2] );
107 p.addParameter(
'decimate',
true );
108 p.addParameter(
'kulso_szabfokok', []); %The list of
sites to be left after the decimation procedure
109 p.addParameter(
'SelfEnergy',
false); %set
true to calculate the Dyson equation with the
self energy
110 p.addParameter(
'keep_sites',
'lead'); %Name of
sites to be kept (scatter, interface, lead)
111 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.
112 p.parse(varargin{:});
113 gfininv = p.Results.gfininv;
114 constant_channels = p.Results.constant_channels;
115 onlyGinverz = p.Results.onlyGinverz;
116 recalculateSurface = p.Results.recalculateSurface;
117 decimate = p.Results.decimate;
118 kulso_szabfokok = p.Results.kulso_szabfokok;
119 useSelfEnergy = p.Results.SelfEnergy;
120 keep_sites = p.Results.keep_sites;
121 UseHamiltonian = p.Results.UseHamiltonian;
124 if ~isempty(recalculateSurface)
126 % creating interfaces
for the Leads
128 shiftLeads = ones(length(obj.param.Leads),1)*obj.E;
130 shiftLeads = ones(length(obj.param.Leads),1)*0;
133 % creating
Lead instaces and calculating the retarded surface Green
operator/
self-energy
134 coordinates_shift = [-2, obj.height+1] + obj.shift;
135 obj.FL_handles.LeadCalc(
'coordinates_shift', coordinates_shift,
'shiftLeads', shiftLeads,
'transversepotential', obj.transversepotential, ...
136 'SelfEnergy', useSelfEnergy,
'SurfaceGreensFunction', ~useSelfEnergy,
'gauge_field', obj.gauge_field,
'leads', recalculateSurface,
'q', obj.q, ...
137 'leadmodel', obj.leadmodel,
'bias_leads', obj.bias_leads,
'T', obj.T);
139 for idx = 1:length(recalculateSurface)
140 obj.CreateInterface( recalculateSurface(idx),
'UseHamiltonian', UseHamiltonian );
146 'onlyGinverz', onlyGinverz, ...
147 'recalculateSurface', [], ...
148 'decimate', decimate, ...
149 'kulso_szabfokok', kulso_szabfokok, ...
150 'SelfEnergy', useSelfEnergy, ...
151 'keep_sites', keep_sites);
157 %> @brief Sets the energy
for the calculations
158 %> @
param Energy The value of the energy in the same units as the Hamiltonian.
159 function setEnergy( obj, Energy )
163 if ~isempty( obj.Scatter_UC ) && strcmpi(
class(obj.Scatter_UC),
'Lead_Keldysh')
164 obj.Scatter_UC.Reset();
167 if ~isempty(obj.Scatter_UC)
168 obj.CreateRibbon('justHamiltonians', true);
174 end % end methods public
178 methods (Access=protected)
181 %> @brief Initializes the attributes of the class.
186 obj.Scatter_UC = obj.FL_handles.SurfaceGreenFunctionCalculator([], 'createCore', 1, 'q', obj.q, 'T', obj.T, 'mu', obj.mu);
192 %> @brief Parses the optional parameters for the class constructor.
193 %> @
param varargin Cell array of optional parameters (https:
194 %> @
param 'filenameIn' The input filename containing the computational parameters. (Use parameters 'Op' and '
param' instead)
195 %> @
param 'filenameOut' The output filename to export the computational parameters.
196 %> @
param 'WorkingDir' The absolute path to the working directoy.
197 %> @
param 'CustomHamiltoniansHandle'
function handle for the custom
Hamiltonians. Has the same inputs as
#Custom_Hamiltonians.LoadHamiltonians and output values defined by the example #Hamiltonians. 198 %> @
param 'E' The energy value used in the calculations (in the same units as the Hamiltonian).
199 %> @
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)
200 %> @
param 'silent' Set
true to suppress output messages.
203 %> @
param 'Opt' An instance of the structure #
Opt.
204 %> @
param 'param' An instance of the structure #
param.
205 %> @
param 'q' The transverse momentum quantum number.
206 %> @
param 'mu' An array containing the chemical potentials of the leads.
207 function InputParsing( obj, varargin)
210 p.addParameter(
'filenameIn', obj.filenameIn, @ischar);
211 p.addParameter(
'filenameOut', obj.filenameOut, @ischar);
212 p.addParameter(
'WorkingDir', obj.WorkingDir, @ischar);
213 p.addParameter(
'E', obj.E, @isscalar);
214 p.addParameter(
'silent', obj.silent);
215 p.addParameter(
'leadmodel', obj.leadmodel); %individual physical model
for the contacts
216 p.addParameter(
'interfacemodel', obj.interfacemodel); %individual physical model
for the
interface regions
217 p.addParameter('
Opt', obj.
Opt);
218 p.addParameter(
'param', obj.param);
219 p.addParameter(
'q', obj.q);
220 p.addParameter(
'T', obj.T);
221 p.addParameter(
'EF', obj.EF);
222 p.addParameter(
'bias_leads', obj.bias_leads);
224 p.addParameter(
'width', obj.width);
225 p.addParameter(
'height', obj.height);
226 p.addParameter(
'transversepotential', obj.transversepotential);
228 p.parse(varargin{:});
232 'filenameOut', p.Results.filenameOut, ...
233 'WorkingDir', p.Results.WorkingDir, ...
234 'E', p.Results.E, ...
235 'silent', p.Results.silent, ...
236 'leadmodel', p.Results.leadmodel, ...
237 'interfacemodel', p.Results.interfacemodel, ...
238 'q', p.Results.q, ...
239 'Opt', p.Results.Opt, ...
240 'param', p.Results.param, ...
241 'EF', p.Results.EF, ...
242 'bias_leads', p.Results.bias_leads, ...
245 InputParsing@
Ribbon( obj,
'width', p.Results.width, ...
246 'height', p.Results.height, ...
247 'transversepotential', p.Results.transversepotential);
252 end % methods
private A class representing a two-terminal structure defined on a preprogrammed lattices for steady state no...
Structure Opt contains the basic computational parameters used in EQuUs.
Structure shape contains data about the geometry of the scattering region.
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.
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 createOutput(filename, Opt, param)
This function creates an output file containing the running parameters.
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...
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 ...
Property Lead_Orientation
The orientation of the lead. Set +1 is the "incoming" direction of the propagating states is defined ...
A class describing an N-terminal geometry for steady state non-equilibrium calculations.
Structure sites contains data to identify the individual sites in a matrix.
function SurfaceGreenFunctionCalculator(idx, varargin)
Calculates the surface Green's function or the self energy of a Lead.
Structure junction_sites contains data to identify the individual sites of the leads,...