15 integer,
parameter ::
dp = kind(1.0d0)
19 integer :: norbs, seed
20 character(100) :: jobname
21 character(100) :: bml_type
27 real(
dp) :: dec, rcoeff
40 integer,
parameter :: nkey_char = 2, nkey_int = 2, nkey_re = 7, nkey_log = 1
41 character(len=*) :: filename
44 character(len=50),
parameter :: keyvector_char(nkey_char) = [character(len=50) :: &
45 'JobName=',
'BMLType=']
46 character(len=100) :: valvector_char(nkey_char) = [character(len=100) :: &
47 'GetModelHam',
'Dense' ]
49 character(len=50),
parameter :: keyvector_int(nkey_int) = [character(len=50) :: &
51 integer :: valvector_int(nkey_int) = (/ &
54 character(len=50),
parameter :: keyvector_re(nkey_re) = [character(len=50) :: &
55 'EpsilonA=',
'EpsilonB=',
'DeltaAB=',
'DeltaAiAj=',
'DeltaBiBj=',
'Decay=',
'RCoeff=']
56 real(
dp) :: valvector_re(nkey_re) = (/&
57 0.0, 0.0, -1.0, 0.0, -1.0, -100.0, 0.0 /)
59 character(len=50),
parameter :: keyvector_log(nkey_log) = [character(len=50) :: &
61 logical :: valvector_log(nkey_log) = (/&
65 character(len=50),
parameter :: startstop(2) = [character(len=50) :: &
69 ,keyvector_int,valvector_int,keyvector_re,valvector_re,&
70 keyvector_log,valvector_log,trim(filename),startstop)
73 mham%JobName = valvector_char(1)
75 if(valvector_char(2) ==
"Dense")
then
76 mham%bml_type = bml_matrix_dense
77 elseif(valvector_char(2) ==
"Ellpack")
then
78 mham%bml_type = bml_matrix_ellpack
79 elseif(valvector_char(2) ==
"CSR")
then
80 mham%bml_type = bml_matrix_csr
81 elseif(valvector_char(2) ==
"Ellblock")
then
82 mham%bml_type = bml_matrix_ellblock
86 mham%norbs = valvector_int(1)
87 mham%seed = valvector_int(2)
90 mham%ea = valvector_re(1)
91 mham%eb = valvector_re(2)
92 mham%dab = valvector_re(3)
93 mham%daiaj = valvector_re(4)
94 mham%dbibj = valvector_re(5)
95 mham%dec = valvector_re(6)
96 mham%rcoeff = valvector_re(7)
99 mham%reshuffle = valvector_log(1)
117 & seed, h_bml, verbose)
118 real(
dp),
intent(in) :: ea, eb, dab, daiaj, dbibj, rcoeff
119 integer,
intent(in) :: verbose, seed
120 integer,
allocatable :: seedin(:)
121 logical,
intent(in) :: reshuffle
122 type(bml_matrix_t),
intent(inout) :: h_bml
123 real(
dp),
allocatable :: diagonal(:), row(:), rowi(:), rowj(:)
124 type(bml_matrix_t) :: ht_bml
125 integer :: norbs, i, j, ssize
126 real(
dp) :: dec, dist, ran
128 norbs = bml_get_n(h_bml)
129 allocate(diagonal(norbs))
133 call random_seed(size=ssize)
134 allocate(seedin(ssize))
136 call random_seed(put=seedin)
139 if(mod(i,2) == 0)
then
140 call random_number(ran)
141 diagonal(i) = ea + rcoeff*(2.0_dp*ran - 1.0_dp)
143 call random_number(ran)
144 diagonal(i) = eb + rcoeff*(2.0_dp*ran - 1.0_dp)
150 if(abs(real(i-j,
dp)) <= norbs/2.0d0)
then
151 dist = max(abs(real(i-j,
dp))- 2.0_dp,0.0_dp)
153 dist = max((-abs(real(i-j,
dp))+norbs) - 2.0_dp,0.0_dp)
156 if((mod(i,2) .ne. 0) .and. (mod(j,2) .ne. 0))
then
157 call random_number(ran)
158 row(j) = (daiaj + rcoeff*(2.0_dp*ran - 1.0_dp))*exp(dec*dist)
160 elseif((mod(i,2) == 0) .and. (mod(j,2) == 0))
then
161 call random_number(ran)
162 row(j) = (dbibj + rcoeff*(2.0_dp*ran - 1.0_dp))*exp(dec*dist)
165 call random_number(ran)
166 row(j) = (dab + rcoeff*(2.0_dp*ran - 1.0_dp))*exp(dec*dist)
171 call bml_set_row(h_bml,i,row)
174 call bml_set_diagonal(h_bml,diagonal)
177 call bml_copy_new(h_bml,ht_bml)
178 call bml_transpose_new(h_bml,ht_bml)
180 call bml_print_matrix(
"h_bml",h_bml,0,10,0,10)
181 call bml_print_matrix(
"ht_bml",ht_bml,0,10,0,10)
183 call bml_add(h_bml,ht_bml,0.5d0,0.5d0,0.0d0)
185 call bml_deallocate(ht_bml)
188 allocate(rowj(norbs))
189 allocate(rowi(norbs))
191 call random_number(ran)
192 j = int(floor(ran*norbs+1))
193 call bml_get_row(h_bml,i,rowi)
194 call bml_get_row(h_bml,j,rowj)
195 call bml_set_row(h_bml,i,rowj)
196 call bml_set_row(h_bml,j,rowi)
218 & seed, h_bml, verbose)
219 real(
dp),
intent(in) :: ea, eb, dab, daiaj, dbibj, rcoeff
220 integer,
intent(in) :: verbose, seed
221 integer,
allocatable :: seedin(:)
222 logical,
intent(in) :: reshuffle
223 type(bml_matrix_t),
intent(inout) :: h_bml
224 real(
dp),
allocatable :: diagonal(:), row1(:), row2(:), rowi(:), rowj(:)
225 type(bml_matrix_t) :: ht_bml
226 integer :: norbs, i, j, ssize, i1, j1, k1, i2, j2, k2, n1d
227 real(
dp) :: dec, dist, di, dj, dk, ran, tmp
229 norbs = bml_get_n(h_bml)
231 allocate(diagonal(norbs))
232 allocate(row1(norbs))
233 allocate(row2(norbs))
242 call random_seed(size=ssize)
243 allocate(seedin(ssize))
245 call random_seed(put=seedin)
248 if(mod(i,2) == 0)
then
249 call random_number(ran)
250 diagonal(i) = ea + rcoeff*(2.0_dp*ran - 1.0_dp)
252 call random_number(ran)
253 diagonal(i) = eb + rcoeff*(2.0_dp*ran - 1.0_dp)
260 i = (i1-1) + n1d*(j1-1) + n1d*n1d*(k1-1)
264 j = (i2-1) + n1d*(j2-1) + n1d*n1d*(k2-1)
265 if(abs(real(i1-i2,
dp)) <= n1d/2.0d0)
then
266 di = max(abs(real(i1-i2,
dp))- 2.0_dp,0.0_dp)
268 di = max((-abs(real(i1-i2,
dp))+n1d)- 2.0_dp,0.0_dp)
270 if(abs(real(j1-j2,
dp)) <= n1d/2.0d0)
then
271 dj = max(abs(real(j1-j2,
dp))- 2.0_dp,0.0_dp)
273 dj = max((-abs(real(j1-j2,
dp))+n1d)- 2.0_dp,0.0_dp)
275 if(abs(real(k1-k2,
dp)) <= n1d/2.0d0)
then
276 dk = max(abs(real(k1-k2,
dp))- 2.0_dp,0.0_dp)
278 dk = max((-abs(real(k1-k2,
dp))+n1d)- 2.0_dp,0.0_dp)
280 dist = dsqrt(di+di+dj*dj+dk*dk)
284 call random_number(ran)
285 row1(2*j+1) = (daiaj + rcoeff*(2.0_dp*ran - 1.0_dp))*exp(dec*dist)
287 call random_number(ran)
288 row1(2*j+2) = (dab + rcoeff*(2.0_dp*ran - 1.0_dp))*exp(dec*dist)
289 row2(2*j+1) = (dab + rcoeff*(2.0_dp*ran - 1.0_dp))*exp(dec*dist)
291 call random_number(ran)
292 row2(2*j+2) = (dbibj + rcoeff*(2.0_dp*ran - 1.0_dp))*exp(dec*dist)
297 call bml_set_row(h_bml,2*i+1,row1)
298 call bml_set_row(h_bml,2*i+2,row2)
303 call bml_set_diagonal(h_bml,diagonal)
306 call bml_copy_new(h_bml,ht_bml)
307 call bml_transpose_new(h_bml,ht_bml)
309 call bml_print_matrix(
"h_bml",h_bml,0,10,0,10)
310 call bml_print_matrix(
"ht_bml",ht_bml,0,10,0,10)
312 call bml_add(h_bml,ht_bml,0.5d0,0.5d0,0.0d0)
314 call bml_deallocate(ht_bml)
317 allocate(rowj(norbs))
318 allocate(rowi(norbs))
320 call random_number(ran)
321 j = int(floor(ran*norbs+1))
322 call bml_get_row(h_bml,i,rowi)
323 call bml_get_row(h_bml,j,rowj)
324 call bml_set_row(h_bml,i,rowj)
325 call bml_set_row(h_bml,j,rowi)
Some general parsing functions.
subroutine, public prg_parsing_kernel(keyvector_char, valvector_char, keyvector_int, valvector_int, keyvector_re, valvector_re, keyvector_log, valvector_log, filename, startstop)
The general parsing function. It is used to vectorize a set of "keywords" "value" pairs as included i...
The prg_hamiltonian module.
subroutine, public prg_parse_mham(mham, filename)
Model Ham parse.
subroutine, public prg_twolevel_model(ea, eb, dab, daiaj, dbibj, dec, rcoeff, reshuffle, seed, h_bml, verbose)
Construct a two-level model Hamiltonian.
subroutine, public prg_twolevel_model3d(ea, eb, dab, daiaj, dbibj, dec, rcoeff, reshuffle, seed, h_bml, verbose)
Construct a two-level model Hamiltonian.
A module to read and handle chemical systems.