-
Notifications
You must be signed in to change notification settings - Fork 4
/
Copy pathpyruqt.py
437 lines (386 loc) · 21.7 KB
/
pyruqt.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
from ase import transport,Atoms,units
import matplotlib.pyplot as plt
import numpy as np
import ruqt
from numpy import linalg
import sys
class sie_negf:
def __init__(self, **kwargs):
self.input_parameters = {'output' : "pyruqt_results",
'exmol_prog' : "molcas",
'exmol_dir' : None,
'elec_prog' : "molcas",
'elec_dir' : None,
'elec2_dir' : None,
'run_molcas' : False,
'min_trans_energy' : -2.0,
'max_trans_energy' : 2.0,
'delta_energy' : 0.001,
'min_bias' : -2,
'max_bias' : 2,
'delta_bias' : 0.1,
'temp' : 1E-5,
'full_align' : True,
'dft_functional' : "pbe",
'basis_set' : None,
'ecp' : None,
'n_elec_units' : 2,
'exmol_geo' : None,
'elec_geo' : None,
'elec2_geo' : None,
'state_num' : 1,
'state_num_e' : 1,
'coupling_calc' : "none",
'coupled' : "molecule",
'spin_pol' : False,
'align_elec' : 0,
'dos_calc' : False,
'fd_change' : 0.001,
'ase_current' : False,
'elec_size' : [0,0],
'conv_tol' : 1E-7,
'max_iter' : 100,
'pyscf_pbc' : False,
'lattice_v' : None,
'meshnum' : None,
'verbosity' : 2,
'cell_dim' : 1,
'pbc_spin' : None,
'aux_basis' : None,
'pdos_states' : [],
'eigenchannels' : 0}
self.param_update(**kwargs)
def param_update(self,**kwargs):
inp=self.input_parameters
for key in kwargs:
if key in inp:
inp[key] = kwargs[key]
elif key not in inp:
raise KeyError('%r not a vaild keyword. Please check your input parameters.' % key)
def calc_setup(self):
inp=self.input_parameters
outputfile=open(inp['output']+".out",'w')
sys.stderr=open(inp['output']+".err",'w')
print("Performing non-self-consistent NEGF transport calculations using semi-infinite tight-binding electrodes",file=outputfile)
print("Using Atomic Simulation Environment to calculate electrode interactions and transport",file=outputfile)
print("Running Calculation using the following paramaters:",file=outputfile)
print(self.input_parameters,file=outputfile)
energies=np.arange(inp['min_trans_energy'],inp['max_trans_energy']+inp['delta_energy'],inp['delta_energy'])
bias=np.arange(inp['min_bias'],inp['max_bias']+inp['delta_bias'],inp['delta_bias'])
bias=np.where(abs(bias)<0.01, 0.0001, bias)
if inp['exmol_prog']=="molcas":
print("Using Molcas calculation at "+inp['exmol_dir']+" for extended molecular region",file=outputfile)
print("Using the effective Hamiltonian for electronic state "+str(inp['state_num'])+" of extended mol. region",file=outputfile)
elif inp['exmol_prog']=="pyscf":
print("Calculating extended molecular region using Pyscf with "+inp['dft_functional']+" in "+inp['basis_set']+" basis set",file=outputfile)
if inp['elec_prog']=="molcas":
print("Using Molcas calculation at "+inp['elec_dir']+" for left electrode",file=outputfile)
if inp['elec2_dir']!=None:
print("Using non-identical electrodes. Right electrode geometry taken from: "+inp['elec2_dir'],file=outputfile)
else:
print("Assuming symmetric electrodes",file=outputfile)
print("Using the effective Hamiltonian for electronic state "+str(inp['state_num'])+" of extended mol. region",file=outputfile)
if inp['exmol_prog']=="molcas":
h,s,norb,numelec,actorb,actelec,states=ruqt.esc_molcas2(inp['exmol_dir'],"MolEl.dat",inp['state_num'],outputfile)
#h,s=ruqt.esc_molcas(exmol_file,exmol_dir,exmol_molcasd,state_num,outputfile)
elif inp['exmol_prog']=="pyscf":
if inp['pyscf_pbc']==True:
h,s,norb,numelec=ruqt.esc_pyscf_pbc(inp['exmol_dir']+inp['exmol_geo'],inp['dft_functional'],inp['basis_set'],inp['ecp'],inp['conv_tol'],inp['max_iter'],inp['lattice_v'],inp['meshnum'],inp['verbosity'],inp['cell_dim'],inp['pbc_spin'],inp['aux_basis'])
else:
h,s,norb,numelec=ruqt.esc_pyscf(inp['exmol_dir']+inp['exmol_geo'],inp['dft_functional'],inp['basis_set'],inp['ecp'],inp['conv_tol'],inp['max_iter'])
if inp['elec_prog']=="molcas":
h1,s1,norb_le,numelec_le,actorb_le,actelec_le,states_le=ruqt.esc_molcas2(inp['elec_dir'],"MolEl.dat",inp['state_num'],outputfile)
if inp['elec2_dir']!=None:
h2,s2,norb_re,numelec_re,actorb_re,actelec_re,states_re=ruqt.esc_molcas2(inp['elec_dir'],"MolEl.dat",inp['state_num'],outputfile)
else:
h2=None
s2=None
elif inp['elec_prog']=="pyscf":
if inp['pyscf_pbc']==True:
h1,s1,norb_le,numelec_le=ruqt.esc_pyscf_pbc(inp['elec_dir']+inp['elec_geo'],inp['dft_functional'],inp['basis_set'],inp['ecp'],inp['conv_tol'],inp['max_iter'],inp['lattice_v'],inp['meshnum'],inp['verbosity'],inp['cell_dim'],inp['pbc_spin'],inp['aux_basis'])
else:
h1,s1,norb_le,numelec_le=ruqt.esc_pyscf(inp['elec_dir']+inp['elec_geo'],inp['dft_functional'],inp['basis_set'],inp['ecp'],inp['conv_tol'],inp['max_iter'])
if inp['elec2_geo']!=None:
if inp['pyscf_pbc']==True:
h2,s2,norb_re,numelec_re=ruqt.esc_pyscf_pbc(inp['elec_dir']+inp['elec_geo'],inp['dft_functional'],inp['basis_set'],inp['ecp'],inp['conv_tol'],inp['max_iter'],inp['lattice_v'],inp['meshnum'],inp['verbosity'],inp['cell_dim'],inp['pbc_spin'],inp['aux_basis'])
else:
h2,s2,norb_re,numelec_re=ruqt.esc_pyscf(inp['elec_dir']+inp['elec_geo'],inp['dft_functional'],inp['basis_set'],inp['ecp'],inp['conv_tol'],inp['max_iter'])
else:
h2=None
s2=None
elif inp['elec_prog']=="supercell":
l_elec=inp['elec_size'][0]
r_elec=inp['elec_size'][1]
h1=h[:l_elec,:l_elec]
h2=h[-r_elec:,-r_elec:]
s1=s[:l_elec,:l_elec]
s2=s[-r_elec:,-r_elec:]
if inp['coupling_calc']=="Fock_EX":
hc1,sc1,hc2,sc2=ruqt.calc_coupling(h,s,h1,h2,s1,s2,inp['coupled'],inp['n_elec_units'])
else:
hc1=None
sc1=None
hc2=None
sc2=None
return(energies,bias,outputfile,h,h1,h2,s,s1,s2,hc1,sc1,hc2,sc2)
def transmission(self):
energies,bias,outputfile,h,h1,h2,s,s1,s2,hc1,sc1,hc2,sc2=sie_negf.calc_setup(self)
inp=self.input_parameters
print("Calculating "+str(len(energies))+" transmission energies from: "+str(max(energies))+" eV to "+str(min(energies))+" eV",file=outputfile)
print("Final transmission values will be printed to "+inp['output']+".trans"+" in relative transmission vs eV",file=outputfile)
print("Performing NEGF Transport Calculations using the Atomic Simulation Environment",file=outputfile)
if inp['full_align']==True:
h=ruqt.full_align(h,h1,s,inp['n_elec_units'])
if inp['align_elec']>=1:
inp['align_elec']-=1
print("Aligning the "+str(inp['align_elec'])+" element of both electrode and extended molecule",file=outputfile)
calc = transport.TransportCalculator(h=h, h1=h1,h2=h2, s=s, s1=s1,s2=s2,hc1=hc1,sc1=sc1,hc2=hc2,sc2=sc2, energies=energies,dos=inp['dos_calc'],logfile=inp['output']+".trans",align_bf=inp['align_elec'],pdos=inp['pdos_states'],eigenchannels=inp['eigenchannels'])
elif inp['align_elec']<1:
calc=transport.TransportCalculator(h=h, h1=h1,h2=h2, s=s, s1=s1,s2=s2,hc1=hc1,sc1=sc1,hc2=hc2,sc2=sc2, energies=energies,dos=inp['dos_calc'],logfile=inp['output']+".trans",pdos=inp['pdos_states'],eigenchannels=inp['eigenchannels'])
T = calc.get_transmission()
t_plot=plt.plot(energies, T)
plt.xlabel('E-E(Fermi) (eV)')
plt.ylabel('Transmission (rel)')
plt.savefig(inp['output']+"_trans.png")
plt.clf()
#Automatic printing and plotting routines for Total and Projected Density of States
if inp['dos_calc']==True or inp['pdos_states']!=[]:
dosfile=open(inp['output']+".dos",'w')
if inp['dos_calc']==True:
print("Calculating Total DOS",file=outputfile)
d_plot=plt.plot(energies,calc.dos_e)
plt.xlabel('E-E(Fermi) (eV)')
plt.ylabel('Density of States')
plt.savefig(inp['output']+"_dos.png")
plt.clf()
print("Total DOS",len(calc.dos_e),file=dosfile)
temp=0
while temp <= len(energies)-1:
print(energies[temp],calc.dos_e[temp],file=dosfile)
temp+=1
if temp > len(energies):
break
if inp['pdos_states']!=[]:
pdos_len=len(inp['pdos_states'])
print("Plotting PDOS States:"+str(pdos_len),file=outputfile)
vec=0
pdos_ne_comb=np.zeros(len(energies))
while vec <= pdos_len-1:
#print(str(calc.pdos_ne.shape),file=outputfile)
d_plot=plt.plot(energies,calc.pdos_ne[vec,:])
plt.xlabel('E-E(Fermi) (eV)')
plt.ylabel('Projected Density of States')
plt.savefig(inp['output']+"_pdos"+str(vec+1)+".png")
plt.clf()
print("PDOS State",inp['pdos_states'][vec],len(energies),file=dosfile)
temp=0
while temp <= len(energies)-1:
print(energies[temp],calc.pdos_ne[vec,temp],file=dosfile)
pdos_ne_comb[temp]+=calc.pdos_ne[vec,temp]
temp+=1
if temp > len(energies):
break
vec+=1
if pdos_len < 0 or vec > pdos_len:
break
print("PDOS Combined",file=dosfile)
for i in range(0,len(energies)-1):
print(energies[i],pdos_ne_comb[i],file=dosfile)
d_plot=plt.plot(energies,pdos_ne_comb)
plt.xlabel('E-E(Fermi) (eV)')
plt.ylabel('Projected Density of States')
plt.savefig(inp['output']+"_pdos_combined"+".png")
plt.clf()
#if inp['eigenchannels'] > 0:
# d_plot=plt.plot(energies,calc.dos_e)
# plt.xlabel('E-E(Fermi) (eV)')
# plt.ylabel('Eigenchannels')
# plt.savefig(inp['output']+"_eigen.png")
# plt.clf()
return T,calc
def current(self):
energies,bias,outputfile,h,h1,h2,s,s1,s2,hc1,sc1,hc2,sc2=sie_negf.calc_setup(self)
inp=self.input_parameters
print("Performing a Landauer current calculation for the following bias voltage range(V): "+str(bias),file=outputfile)
print("Calculating "+str(len(energies))+" transmission energies from: "+str(max(energies))+" eV to "+str(min(energies))+" eV",file=outputfile)
print("Final transmission values will be printed to "+inp['output']+".trans"+" in relative transmission vs eV",file=outputfile)
print("Final current values will be printed to "+inp['output']+".iv"+" in volts vs ampheres",file=outputfile)
print("Final conductance values will be printed to "+inp['output']+".con"+" in volts vs G_0",file=outputfile)
T,calc=sie_negf.transmission(self)
if inp['ase_current']==True:
I = calc.get_current(bias,T=inp['temp'],E=energies,T_e=T,spinpol=inp['spin_pol'])
else:
I = ruqt.get_current(bias,T=inp['temp'],E=energies,T_e=T,spinpol=inp['spin_pol'],delta_e=inp['delta_energy'])
i_plot=plt.plot(bias,2.*units._e**2/units._hplanck*I)
plt.xlabel('Voltage (V)')
plt.ylabel('Current (A)')
plt.savefig(inp['output']+"_current.png")
b_range=len(bias)
cond=np.zeros(b_range)
for x in range(len(bias)):
cond[x]=I[x]/bias[x]
I=2.*units._e**2/units._hplanck*I
np.savetxt(inp['output']+".iv",np.c_[bias,I],fmt="%s")
plt.clf()
c_plot=plt.plot(bias,cond)
plt.xlabel('Voltage (V)')
plt.ylabel('Conductance (G_0)')
plt.savefig(inp['output']+"_conductance.png")
np.savetxt(inp['output']+".con",np.c_[bias,cond],fmt="%s")
def diff_conductance(self):
energies,bias,outputfile,h,h1,h2,s,s1,s2,hc1,sc1,hc2,sc2=sie_negf.calc_setup(self)
inp=self.input_parameters
print("Calculating differential conductance using numerical derivatives",file=outputfile)
print("Calculting each value using the +/-"+str(inp['fd_change'])+" voltage points around it.",file=outputfile)
print("Performing the diff. conductance calculation for the following bias voltage range(V): "+str(bias),file=outputfile)
print("Calculating "+str(len(energies))+" transmission energies from: "+str(max(energies))+" eV to "+str(min(energies))+" eV around electrode Fermi level.",file=outputfile)
print("Final transmission values will be printed to "+inp['output']+".trans"+" in relative transmission vs eV",file=outputfile)
print("Final diff. conductance values will be printed to "+inp['output']+".dcon"+" in volts vs G_0",file=outputfile)
T,calc=sie_negf.transmission(self)
DE=ruqt.get_diffcond(calc,bias,inp['temp'],energies,T,inp['fd_change'],inp['ase_current'],inp['delta_energy'])
c_plot=plt.plot(bias,DE)
plt.xlabel('Voltage (V)')
plt.ylabel('Diff. Conductance (G_0)')
plt.savefig(inp['output']+"_diffcon.png")
np.savetxt(inp['output']+".dcon",np.c_[bias,DE],fmt="%s")
class wbl_negf:
def __init__(self, **kwargs):
self.input_parameters = {'output' : "pyruqt.results",
'exmol_dir' : None,
'num_elec_atoms' : None,
'exmol_prog' : "molcas",
'run_molcas' : False,
'min_trans_energy' : -2.0,
'max_trans_energy' : 2.0,
'delta_energy' : 0.01,
'min_bias' : -2,
'max_bias' : 2,
'delta_bias' : 0.1,
'temp' : 1E-5,
'dft_functional' : None,
'basis_set' : None,
'ecp' : None,
'exmol_geo' : None,
'state_num' : 1,
'FermiE' :-5.30,
'FermiD' : 0.07,
'qc_method' : "dft",
'rdm_type' : 1,
'fort_data' : None,
'fort_trans' : False,
'fd_change' : 0.001,
'ase_current' : False,
'conv_tol' : 1E-7,
'max_iter' : 100}
self.param_update(**kwargs)
def param_update(self,**kwargs):
inp=self.input_parameters
for key in kwargs:
if key in inp:
inp[key] = kwargs[key]
elif key not in inp:
raise KeyError('%r not a vaild keyword. Please check your input parameters.' % key)
def calc_setup(self):
inp=self.input_parameters
outputfile=open(inp['output']+".out",'w')
sys.stderr=open(inp['output']+".err",'w')
print("Performing non-self-consistent NEGF calculations using metal wide band limit approximation for electrodes",file=outputfile)
print("Using RUQT-Fortran to calculate transport",file=outputfile)
print("Running Calculation using the following paramaters:",file=outputfile)
print(self.input_parameters,file=outputfile)
energies=np.arange(inp['min_trans_energy'],inp['max_trans_energy']+inp['delta_energy'],inp['delta_energy'])
bias=np.arange(inp['min_bias'],inp['max_bias']+inp['delta_bias'],inp['delta_bias'])
bias=np.where(abs(bias)<0.01, 0.01, bias)
if inp['exmol_prog']=="molcas":
print("Using Molcas calculation at "+inp['exmol_dir']+" for extended molecular region",file=outputfile)
print("Using the effective Hamiltonian for electronic state "+str(inp['state_num'])+" of extended mol. region",file=outputfile)
elif inp['exmol_prog']=="pyscf":
print("Calculating extended molecular region using Pyscf with "+inp['dft_functional']+" in "+inp['basis_set']+" basis set",file=outputfile)
if inp['exmol_prog']=="molcas":
h,s,norb,numelec,actorb,actelec,states=ruqt.esc_molcas2(inp['exmol_dir'],"MolEl.dat",inp['state_num'],outputfile)
elec_orb=0
#h,s=ruqt.esc_molcas(exmol_file,exmol_dir,exmol_molcasd,state_num,outputfile)
elif inp['exmol_prog']=="pyscf":
h,s,norb,numelec,elec_orb=ruqt.esc_pyscf_wbl(inp['exmol_dir']+inp['exmol_geo'],inp['dft_functional'],inp['basis_set'],inp['ecp'],inp['conv_tol'],inp['max_iter'],inp['num_elec_atoms'])
return(energies,bias,outputfile,h,s,norb,numelec,elec_orb)
def transmission(self):
inp=self.input_parameters
energies,bias,outputfile,h,s,norb,numelec,elec_orb=wbl_negf.calc_setup(self)
print("Calculating "+str(len(energies))+" transmission energies from: "+str(max(energies))+" eV to "+str(min(energies))+" eV",file=outputfile)
print("Final transmission values will be printed to "+inp['output']+".trans"+" in relative transmission vs eV",file=outputfile)
ruqt.fort_inputwrite("T",inp['FermiE'],inp['FermiD'],inp['temp'],inp['max_bias'],inp['min_bias'],inp['delta_bias'],inp['min_trans_energy'],inp['max_trans_energy'],inp['delta_energy'],inp['qc_method'],inp['rdm_type'],inp['exmol_dir'],inp['fort_data'],inp['exmol_prog'],inp['num_elec_atoms'],outputfile,inp['state_num'],norb,numelec,elec_orb)
T,I=ruqt.fort_calc("RUQT.x","fort_ruqt",energies,bias,"T",outputfile)
t_plot=plt.plot(energies, T)
plt.xlabel('E-E(Fermi) (eV)')
plt.ylabel('Transmission (rel)')
plt.savefig(inp['output']+"_trans.png")
plt.clf()
def current(self):
energies,bias,outputfile,h,s,norb,numelec,elec_orb=wbl_negf.calc_setup(self)
inp=self.input_parameters
print("Performing a Landauer current calculation for the following bias voltage range(V): "+str(bias),file=outputfile)
print("Calculating "+str(len(energies))+" transmission energies from: "+str(max(energies))+" eV to "+str(min(energies))+" eV",file=outputfile)
print("Final transmission values will be printed to "+inp['output']+".trans"+" in relative transmission vs eV",file=outputfile)
print("Final current values will be printed to "+inp['output']+".iv"+" in volts vs ampheres",file=outputfile)
print("Final conductance values will be printed to "+inp['output']+".con"+" in volts vs G_0",file=outputfile)
if inp['fort_trans']==False:
ruqt.fort_inputwrite("C",inp['FermiE'],inp['FermiD'],inp['temp'],inp['max_bias'],inp['min_bias'],inp['delta_bias'],inp['min_trans_energy'],inp['max_trans_energy'],inp['delta_energy'],inp['qc_method'],inp['rdm_type'],inp['exmol_dir'],inp['fort_data'],inp['exmol_prog'],inp['num_elec_atoms'],outputfile,inp['state_num'],norb,numelec,elec_orb)
T,I=ruqt.fort_calc("RUQT.x","fort_ruqt",energies,bias,"C",outputfile)
elif inp['fort_trans']==True:
print("Calculating current with ASE transport",file=outputfile)
ruqt.fort_inputwrite("T",inp['FermiE'],inp['FermiD'],inp['temp'],inp['max_bias'],inp['min_bias'],inp['delta_bias'],inp['min_trans_energy'],inp['max_trans_energy'],inp['delta_energy'],inp['qc_method'],inp['rdm_type'],inp['exmol_dir'],inp['fort_data'],inp['exmol_prog'],inp['num_elec_atoms'],outputfile,inp['state_num'],norb,numelec,elec_orb)
T,I=ruqt.fort_calc("RUQT.x","fort_ruqt",energies,bias,"T",outputfile)
h=np.zeros((2,2))
h1=np.zeros((2,2))
calc = transport.TransportCalculator(h=h, h1=h1, energies=energies,logfile="temp")
if inp['ase_current']==True:
I = calc.get_current(bias,T=inp['temp'],E=energies,T_e=T,spinpol=inp['spin_pol'])
else:
I = ruqt.get_current(bias,T=inp['temp'],E=energies,T_e=T,spinpol=inp['spin_pol'],delta_e=inp['delta_energy'])
I=2.*units._e**2/units._hplanck*I
t_plot=plt.plot(energies, T)
plt.xlabel('E-E(Fermi) (eV)')
plt.ylabel('Transmission (rel)')
plt.savefig(inp['output']+"_trans.png")
plt.clf()
i_plot=plt.plot(bias,I)
plt.xlabel('Voltage (V)')
plt.ylabel('Current (A)')
plt.savefig(inp['output']+"_current.png")
b_range=len(bias)
cond=np.zeros(b_range)
for x in range(len(bias)):
cond[x]=I[x]/(bias[x]*2.*units._e**2/units._hplanck)
np.savetxt(inp['output']+".iv",np.c_[bias,I],fmt="%s")
plt.clf()
c_plot=plt.plot(bias,cond)
plt.xlabel('Voltage (V)')
plt.ylabel('Conductance (G_0)')
plt.savefig(inp['output']+"_conductance.png")
np.savetxt(inp['output']+".con",np.c_[bias,cond],fmt="%s")
def diff_conductance(self):
energies,bias,outputfile,h,s,norb,numelec,elec_orb=wbl_negf.calc_setup(self)
inp=self.input_parameters
print("Calculating differential conductance using numerical derivatives",file=outputfile)
print("Calculting each value using the +/-"+str(inp['fd_change'])+" voltage points around it.",file=outputfile)
print("Performing the diff. conductance calculation for the following bias voltage range(V): "+str(bias),file=outputfile)
print("Calculating "+str(len(energies))+" transmission energies from: "+str(max(energies))+" eV to "+str(min(energies))+" eV around electrode Fermi level.",file=outputfile)
print("Final transmission values will be printed to "+inp['output']+".trans"+" in relative transmission vs eV",file=outputfile)
print("Final diff. conductance values will be printed to "+inp['output']+".dcon"+" in volts vs G_0",file=outputfile)
print("Not available in RUQT-Fortran. Using RUQT-Fortan transmission with pyRUQT DiffCond calculator.",file=outputfile)
ruqt.fort_inputwrite("T",inp['FermiE'],inp['FermiD'],inp['temp'],inp['max_bias'],inp['min_bias'],inp['delta_bias'],inp['min_trans_energy'],inp['max_trans_energy'],inp['delta_energy'],inp['qc_method'],inp['rdm_type'],inp['exmol_dir'],inp['fort_data'],inp['exmol_prog'],inp['num_elec_atoms'],outputfile,inp['state_num'],norb,numelec,elec_orb)
T,I=ruqt.fort_calc("RUQT.x","fort_ruqt",energies,bias,"T",outputfile)
t_plot=plt.plot(energies, T)
plt.xlabel('E-E(Fermi) (eV)')
plt.ylabel('Transmission (rel)')
plt.savefig(inp['output']+"_trans.png")
plt.clf()
h=np.zeros((2,2))
h1=np.zeros((2,2))
calc = transport.TransportCalculator(h=h, h1=h1, energies=energies,logfile="temp")
DE=ruqt.get_diffcond(calc,bias,inp['temp'],energies,T,inp['fd_change'],inp['ase_current'],inp['delta_energy'])
c_plot=plt.plot(bias,DE)
plt.xlabel('Voltage (V)')
plt.ylabel('Diff. Conductance (G_0)')
plt.savefig(inp['output']+"_diffcon.png")
np.savetxt(inp['output']+".dcon",np.c_[bias,DE],fmt="%s")