1!--------------------------------------------------------------------------------------------------!
2!  DFTB+: general package for performing fast atomistic simulations                                !
3!  Copyright (C) 2006 - 2019  DFTB+ developers group                                               !
4!                                                                                                  !
5!  See the LICENSE file for terms of usage and distribution.                                       !
6!--------------------------------------------------------------------------------------------------!
7
8module libnegf_vars
9  use dftbp_accuracy, only : dp, mc, lc
10  use dftbp_commontypes
11  use dftbp_wrappedintr
12  use dftbp_xmlf90
13  implicit none
14  private
15
16  public :: TNEGFTunDos
17  public :: TNEGFGreenDensInfo
18  public :: TTransPar
19  public :: ContactInfo
20  public :: TElph
21
22
23
24  !> Options for electron-phonon model
25  type TElPh
26
27    !> True if filled up with info from an input block
28    logical :: defined = .false.
29
30    !> Specify which model in input, 1-3 elastic models: diagonal, block diagonal and overlap masked
31    integer :: model = 0
32
33    !> Coupling strength (?)
34    real(dp), allocatable :: coupling(:)
35
36    !> Iterations for self-consistent Born approximation
37    integer :: scba_niter = 0
38
39    !> List of orbital per atom for models = (2,3)
40    integer, allocatable :: orbsperatm(:)
41
42  end type TElPh
43
44
45  !Structure for contact information in a transport calculation
46  type ContactInfo
47
48    ! Beginning (1) and end (2) of contact in atoms (?)
49    integer :: idxrange(2)
50
51    !> Note: a contact id is specifically defined because, with multiple definition of contacts in
52    !> the input file, relying on contact ordering to assign an integer can be inconsistent
53    character(mc) :: name
54
55    !> Accuracy of rigid layer shift
56    real(dp) :: shiftAccuracy = 0.0_dp
57
58    integer :: dir = 0
59
60    real(dp) :: length = 0.0_dp
61
62    !> contact vector
63    real(dp) :: lattice(3)
64
65    real(dp) :: potential = 0.0_dp
66
67    !> for colinear spin we may need two Fermi levels (up and down)
68    real(dp) :: eFermi(2) = [0.0_dp, 0.0_dp]
69
70    !> contact temperature
71    real(dp) :: kbT = 0.0_dp
72
73    ! Is it a contact in the wide band approximation?
74    logical :: wideBand = .false.
75
76    real(dp) :: wideBandDos = 0.0_dp
77
78    !> Filename for contact infos (shiftcont_) TO BE MOVED?
79    character(lc) :: output
80
81    logical :: tWriteSelfEnergy = .false.
82    logical :: tReadSelfEnergy = .false.
83    logical :: tWriteSurfaceGF = .false.
84    logical :: tReadSurfaceGF = .false.
85
86  end type ContactInfo
87
88
89  !> Options for Landauer (Tunneling and DOS) calculation
90  type TNEGFTunDos
91
92    !> true only if filling block is defined
93    logical :: defined = .false.
94
95    !> verbosity level of the library
96    integer :: verbose
97
98    !> spin degeneracy (used in transmission and current integration)
99    integer :: gSpin
100
101    !> Min integration energy (possible to define them different for colinear spin calculation)
102    real(dp) :: emin
103
104    !> Max integration energy
105    real(dp) :: emax
106
107    !> Energy step
108    real(dp) :: estep
109
110    !> Delta for Green's function
111    real(dp) :: delta
112
113    !> An additional broadening delta for DOS and tunneling
114    real(dp) :: broadeningDelta
115
116    !> emitter contact(s)
117    integer, allocatable :: ni(:)
118
119    !> collector contact(s)
120    integer, allocatable :: nf(:)
121
122    !> Orbitals in regions
123    type(WrappedInt1), allocatable :: dosOrbitals(:)
124
125    !> Labels of regions for LDOS calculations
126    character(lc), allocatable :: dosLabels(:)
127
128    !> write DOS on separate files
129    logical :: writeLDOS = .false.
130
131    !> write tunneling on separate files
132    logical :: writeTunn = .false.
133
134    !> contact temperatures
135    real(dp), allocatable :: kbT(:)
136
137    !> Electron-phonon coupling
138    type(Telph) :: elph
139
140    !> Buttiker Probe for dephasing
141    type(Telph) :: bp
142
143  end type TNEGFTunDos
144
145
146  !> Information for Green's function charge density calculation
147  type TNEGFGreenDensInfo
148
149    !> true only if filling block is defined
150    logical :: defined = .false.
151
152    !> verbosity level of the library
153    integer :: verbose
154
155    !> Fermi level for closed system calculation. If a coliner spin calculation is defined, two
156    !> values are needed (up and down) unique Fermi closed systems
157    real(dp) :: oneFermi(2) = [0.0_dp, 0.0_dp]
158
159    !> delta function in G.F.
160    real(dp) :: delta
161
162    !> Number of points in contour
163    integer :: nP(3)
164
165    !> Lowest energy for contour int
166    real(dp) :: enLow
167
168    !> Number of kT for Fermi dist
169    integer :: nkT
170
171    !> Number of poles included in contour
172    integer :: nPoles
173
174    !> use or not Green solver
175    logical :: doGreenDens = .false.
176
177    !> save SGF in files
178    logical :: saveSGF
179
180    !> read SGF from files
181    logical :: readSGF
182
183    !> Calculate or not the local J. There is an independent definition of principal layers (pls),
184    !> since in a closed system Green's calculation a separate definition may be used
185    logical :: doLocalCurr = .false.
186
187    !> Number of principal layers
188    integer :: nPLs = 0
189
190    !> PL indices (starting atom)
191    integer, allocatable :: PL(:)
192
193    !> spin degeneracy (used in charge integration)
194    integer :: gSpin
195
196    !> contact temperatures
197    real(dp), allocatable :: kbT(:)
198
199    !> Electron-phonon coupling
200    type(Telph) :: elph
201
202    !> Buttiker Probe for dephasing
203    type(Telph) :: bp
204
205  end type TNEGFGreenDensInfo
206
207
208  !> Options from Transport section (geometry and task)
209  type TTransPar
210
211    !> True if the corresponding input block exists
212    logical :: defined = .false.
213
214    !> Contacts in the system
215    type(ContactInfo), allocatable :: contacts(:)
216
217    !> Number of contacts
218    integer :: ncont = 0
219
220    !> Start and end index of device region
221    integer :: idxdevice(2)
222
223    !> Number of principal layers
224    integer :: nPLs =1
225
226    !> PL indices (starting atom)
227    integer, allocatable, dimension(:) :: PL
228
229    !> False: run the full OBC calculation / True: upload contact phase
230    logical :: taskUpload = .false.
231
232    !> Index of contact for contact hamiltonian task, if any
233    integer :: taskContInd = 0
234
235    !> Not from input file
236    logical :: tPeriodic1D = .false.
237
238
239    !DAR begin - type TTransPar new items
240    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
241    logical :: tNoGeometry = .false.
242    logical :: tWriteDFTB = .false.
243    logical :: tReadDFTB = .false.
244    logical :: tOrthonormal = .false.
245    logical :: tOrthonormalDevice = .false.
246    logical :: tModel = .false.
247    logical :: tRead_negf_in = .false.
248    integer :: NumStates = 0
249    integer, allocatable :: cblk(:)
250    character(lc) :: FileName
251    logical :: tManyBody =.false.
252    logical :: tElastic =.true.
253    logical :: tDephasingVE = .false.
254    logical :: tDephasingBP = .false.
255    logical :: tZeroCurrent = .false.
256    integer :: MaxIter = 1000
257    logical :: tWriteDOS = .false.
258    logical :: tWrite_ldos = .false.
259    logical :: tWrite_negf_params = .false.
260    logical :: tDOSwithS =.true.
261    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
262    !DAR end
263
264  end type TTransPar
265
266contains
267
268  !> Copies contents of a Green's function density calculation structure
269  subroutine copyGreenDens(gIN, gOUT)
270
271    !> Original structure
272    type(TNEGFGreenDensInfo), intent(in) :: gIN
273
274    !> Duplicate
275    type(TNEGFGreenDensInfo), intent(out) :: gOUT
276
277    integer :: isz
278
279    ! copy greendens
280    gOUT%defined = gIN%defined
281    gOUT%verbose = gIN%verbose
282    gOUT%oneFermi = gIN%oneFermi
283    gOUT%delta = gIN%delta
284    gOUT%np = gIN%np
285    gOUT%enLow = gIN%enLow
286    gOUT%nkT = gIN%nkT
287    gOUT%nPoles = gIN%nPoles
288    gOUT%doGreenDens= gIN%doGreenDens
289    gOUT%saveSGF = gIN%saveSGF
290    gOUT%readSGF = gIN%readSGF
291    gOUT%doLocalCurr= gIN%doLocalCurr
292    gOUT%nPLs = gIN%nPLs
293    gOUT%gspin = gIN%gspin
294
295    if (allocated(gIN%kbT)) then
296      isz = size(gIN%kbT)
297      allocate(gOUT%kbT(isz))
298      gOUT%kbT = gIN%kbT
299    end if
300    if (allocated(gIN%PL)) then
301      isz = size(gIN%PL)
302      allocate(gOUT%PL(isz))
303      gOUT%PL = gIN%PL
304    end if
305
306    ! copy greendens%elph
307    if (allocated(gIN%elph%coupling)) then
308      isz = size(gIN%elph%coupling)
309      allocate(gOUT%elph%coupling(isz))
310      gOUT%elph%coupling = gIN%elph%coupling
311    end if
312    if (allocated(gIN%elph%orbsperatm)) then
313      isz = size(gIN%elph%orbsperatm)
314      allocate(gOUT%elph%orbsperatm(isz))
315      gOUT%elph%orbsperatm = gIN%elph%orbsperatm
316    end if
317    gOUT%elph%defined = gIN%elph%defined
318    gOUT%elph%model = gIN%elph%model
319    gOUT%elph%scba_niter = gIN%elph%scba_niter
320
321    ! copy greendens%bp
322    if (allocated(gIN%bp%coupling)) then
323      isz = size(gIN%bp%coupling)
324      allocate(gOUT%bp%coupling(isz))
325      gOUT%bp%coupling = gIN%bp%coupling
326    end if
327    if (allocated(gIN%bp%orbsperatm)) then
328      isz = size(gIN%bp%orbsperatm)
329      allocate(gOUT%bp%orbsperatm(isz))
330      gOUT%bp%orbsperatm = gIN%bp%orbsperatm
331    end if
332    gOUT%bp%defined = gIN%bp%defined
333    gOUT%bp%model = gIN%bp%model
334    gOUT%bp%scba_niter = gIN%bp%scba_niter
335
336  end subroutine copyGreenDens
337
338
339  !> Copies contents of a Landauer calculation structure
340  subroutine copyTunDos(gIN, gOUT)
341
342    !> Original structure
343    type(TNEGFTunDos), intent(in) :: gIN
344
345    !> Duplicate
346    type(TNEGFTunDos), intent(out) :: gOUT
347
348    integer :: isz
349
350    ! copy tundos
351    if (allocated(gIN%ni)) then
352      isz = size(gIN%ni)
353      allocate(gOUT%ni(isz))
354      gOUT%ni = gIN%ni
355    end if
356    if (allocated(gIN%nf)) then
357      isz = size(gIN%nf)
358      allocate(gOUT%nf(isz))
359      gOUT%nf = gIN%nf
360    end if
361    if (allocated(gIN%kbT)) then
362      isz = size(gIN%kbT)
363      allocate(gOUT%kbT(isz))
364      gOUT%kbT = gIN%kbT
365    end if
366    if (allocated(gIN%dosOrbitals)) then
367      isz = size(gIN%dosOrbitals)
368      allocate(gOUT%dosOrbitals(isz))
369      gOUT%dosOrbitals = gIN%dosOrbitals
370    end if
371    if (allocated(gIN%dosLabels)) then
372      isz = size(gIN%dosLabels)
373      allocate(gOUT%dosLabels(isz))
374      gOUT%dosLabels = gIN%dosLabels
375    end if
376
377    gOUT%defined = gIN%defined
378    gOUT%verbose = gIN%verbose
379    gOUT%gspin = gIN%gspin
380    gOUT%emin = gIN%emin
381    gOUT%emax = gIN%emax
382    gOUT%estep = gIN%estep
383    gOUT%delta = gIN%delta
384    gOUT%writeLDOS = gIN%writeLDOS
385    gOUT%writeTUNN = gIN%writeTUNN
386    gOUT%broadeningDelta = gIN%broadeningDelta
387
388    ! copy tundos%elph
389    if (allocated(gIN%elph%coupling)) then
390      isz = size(gIN%elph%coupling)
391      allocate(gOUT%elph%coupling(isz))
392      gOUT%elph%coupling = gIN%elph%coupling
393    end if
394    if (allocated(gIN%elph%orbsperatm)) then
395      isz = size(gIN%elph%orbsperatm)
396      allocate(gOUT%elph%orbsperatm(isz))
397      gOUT%elph%orbsperatm = gIN%elph%orbsperatm
398    end if
399    gOUT%elph%defined = gIN%elph%defined
400    gOUT%elph%model = gIN%elph%model
401    gOUT%elph%scba_niter = gIN%elph%scba_niter
402
403    ! copy tundos%bp
404    if (allocated(gIN%bp%coupling)) then
405      isz = size(gIN%bp%coupling)
406      allocate(gOUT%bp%coupling(isz))
407      gOUT%bp%coupling = gIN%bp%coupling
408    end if
409    if (allocated(gIN%bp%orbsperatm)) then
410      isz = size(gIN%bp%orbsperatm)
411      allocate(gOUT%bp%orbsperatm(isz))
412      gOUT%bp%orbsperatm = gIN%bp%orbsperatm
413    end if
414    gOUT%bp%defined = gIN%bp%defined
415    gOUT%bp%model = gIN%bp%model
416    gOUT%bp%scba_niter = gIN%bp%scba_niter
417
418  end subroutine copyTunDos
419
420
421  !> Copies contents of a transport calculation structure
422  subroutine copyTranspar(tIN, tOUT)
423
424    !> Original structure
425    type(TTransPar), intent(in) :: tIN
426
427    !> Duplicate
428    type(TTransPar), intent(out) :: tOUT
429
430    integer :: isz
431
432    tOUT%defined = tIN%defined
433    tOUT%ncont = tIN%ncont
434    tOUT%idxdevice = tIN%idxdevice
435    tOUT%nPls = tIN%nPls
436    tOUT%taskUpload = tIN%taskUpload
437    tOUT%taskContInd = tIN%taskContInd
438    tOUT%tPeriodic1D = tIN%tPeriodic1D
439    if (allocated(tIN%PL)) then
440      isz = size(tIN%PL)
441      allocate(tOUT%PL(isz))
442      tOUT%PL = tIN%PL
443    end if
444    if (allocated(tIN%contacts)) then
445      isz = size(tIN%contacts)
446      allocate(tOUT%contacts(isz))
447      tOUT%contacts = tIN%contacts
448    end if
449
450  end subroutine copyTranspar
451
452
453end module libnegf_vars
454