20 integer,
parameter ::
dp = kind(1.0d0)
47 mu, beta, occErrLimit, threshold, tol,SCF_IT, occiter, totns)
51 type(bml_matrix_t),
intent(in) :: h_bml
52 type(bml_matrix_t),
intent(inout) :: p_bml, inv_bml(nsteps)
53 integer,
intent(in) :: nsteps, scf_it
54 real(
dp),
intent(in) :: nocc, threshold
55 real(
dp),
intent(in) :: tol
56 real(
dp),
intent(in) :: occerrlimit, beta
57 real(
dp),
intent(inout) :: mu
58 integer,
intent(inout) :: occiter, totns
60 type(bml_matrix_t) :: w_bml, y_bml, d_bml, aux_bml, p2_bml, i_bml, ai_bml
61 real(
dp) :: trdpdmu, trp0, occerr, alpha, newerr
62 real(
dp) :: cnst, ofactor, mustep, preverr
63 real(
dp),
allocatable :: trace(:), gbnd(:)
64 character(20) :: bml_type
65 integer :: n, m, i, iter, muadj, prev, maxiter, nsiter, nsdm
67 bml_type = bml_get_type(h_bml)
72 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, p2_bml)
73 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, d_bml)
74 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, w_bml)
75 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, aux_bml)
76 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, y_bml)
77 call bml_identity_matrix(bml_type, bml_element_real,
dp, n, m, i_bml)
78 call bml_identity_matrix(bml_type, bml_element_real,
dp, n, m, ai_bml)
87 cnst = beta/(1.0_dp*2**(nsteps+2))
89 if (scf_it .eq. 1)
then
93 call bml_copy(h_bml, p_bml)
96 call bml_multiply_x2(p_bml, p2_bml, threshold, trace)
98 call bml_copy(p2_bml, y_bml)
99 call bml_add(y_bml, p_bml, 1.0_dp, -1.0_dp, threshold)
100 call bml_scale_add_identity(y_bml, 2.0_dp, 1.0_dp, threshold)
101 call prg_conjgrad(y_bml, ai_bml, i_bml, aux_bml, d_bml, w_bml, 0.0001_dp, threshold)
103 call bml_copy(i_bml, inv_bml(i))
107 do while ((occerr .gt. occerrlimit .or. muadj .eq. 1) .and. iter < maxiter)
110 write(*,*)
'mu =', mu
113 call bml_copy(h_bml, p_bml)
118 call bml_multiply_x2(p_bml, p2_bml, threshold, trace)
120 call bml_copy(p2_bml, y_bml)
121 call bml_add(y_bml, p_bml, 1.0_dp, -1.0_dp, threshold)
122 call bml_scale_add_identity(y_bml, 2.0_dp, 1.0_dp, threshold)
126 if (iter .eq. 1)
then
127 call prg_conjgrad(y_bml, inv_bml(i), i_bml, aux_bml, d_bml, w_bml, 0.001_dp, threshold)
129 call prg_newtonschulz(y_bml, inv_bml(i), d_bml, w_bml, aux_bml, i_bml, tol, threshold, nsiter)
131 call bml_multiply(inv_bml(i), p2_bml, p_bml, 1.0_dp, 0.0_dp, threshold)
135 write(*,*)
'Number of Newton-Schulz iterations:',nsdm
138 trp0 = bml_trace(p_bml)
139 trdpdmu = beta*(trp0 - bml_sum_squares(p_bml))
140 occerr = abs(trp0 - nocc)
141 write(*,*)
'occerr =', occerr
144 if (occerr > 1.0_dp)
then
145 if (nocc-trp0 < 0.0_dp .and. prev .eq. -1)
then
147 else if (nocc-trp0 > 0.0_dp .and. prev .eq. 1)
then
149 else if (nocc-trp0 > 0.0_dp .and. prev .eq. -1)
then
160 else if (occerr .gt. occerrlimit)
then
161 mustep = (nocc -trp0)/trdpdmu
168 if (iter .ge. maxiter)
then
169 write(*,*)
'Could not converge chemical potential in prg_impplicit_fermi_save_inverse'
181 occiter = occiter + iter
182 call bml_scale(2.0_dp,p_bml)
185 call bml_deallocate(p2_bml)
186 call bml_deallocate(w_bml)
187 call bml_deallocate(d_bml)
188 call bml_deallocate(y_bml)
189 call bml_deallocate(aux_bml)
190 call bml_deallocate(ai_bml)
191 call bml_deallocate(i_bml)
210 mu, beta, method, osteps, occErrLimit, threshold, tol)
214 type(bml_matrix_t),
intent(in) :: h_bml
215 type(bml_matrix_t),
intent(inout) :: p_bml
216 integer,
intent(in) :: osteps, nsteps, method, k
217 real(
dp),
intent(in) :: nocc, threshold
218 real(
dp),
intent(in) :: tol
219 real(
dp),
intent(in) :: occerrlimit, beta
220 real(
dp),
intent(inout) :: mu
222 type(bml_matrix_t) :: w_bml, y_bml, d_bml, p2_bml, aux1_bml, aux2_bml, i_bml, ai_bml
223 real(
dp) :: trdpdmu, trp0, occerr
224 real(
dp) :: cnst, ofactor
225 real(
dp),
allocatable :: trace(:), gbnd(:)
226 character(20) :: bml_type
227 integer :: n, m, i, iter, exp_order
229 bml_type = bml_get_type(h_bml)
234 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, p2_bml)
235 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, d_bml)
236 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, w_bml)
237 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, y_bml)
239 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, aux1_bml)
240 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, aux2_bml)
242 if (method .eq. 1)
then
243 call bml_identity_matrix(bml_type, bml_element_real,
dp, n, m, i_bml)
244 call bml_identity_matrix(bml_type, bml_element_real,
dp, n, m, ai_bml)
249 exp_order = k**nsteps
250 cnst = beta/(4*exp_order)
253 do while ((osteps .eq. 0 .and. occerr .gt. occerrlimit) .or. &
254 (osteps .gt. 0 .and. iter .lt. osteps))
259 call bml_copy(h_bml, p_bml)
262 if (method .eq. 0)
then
263 write(*,*)
"Doing CG"
267 call bml_multiply_x2(p_bml, p2_bml, threshold, trace)
270 call bml_copy(p2_bml, y_bml)
271 call bml_add(y_bml, p_bml, 1.0_dp, -1.0_dp, threshold)
272 call bml_scale_add_identity(y_bml, 2.0_dp, 1.0_dp, threshold)
274 call prg_setup_linsys(p_bml, y_bml, p2_bml, d_bml, w_bml, aux1_bml, aux2_bml, k, threshold)
276 call prg_conjgrad(y_bml, p_bml, p2_bml, d_bml, aux1_bml, w_bml, tol, threshold)
279 write(*,*)
"Doing NS"
283 call bml_multiply_x2(p_bml, p2_bml, threshold, trace)
286 call bml_copy(p2_bml, y_bml)
287 call bml_add(y_bml, p_bml, 1.0_dp, -1.0_dp, threshold)
288 call bml_scale_add_identity(y_bml, 2.0_dp, 1.0_dp, threshold)
290 call prg_setup_linsys(p_bml, y_bml, p2_bml, d_bml, w_bml, aux1_bml, aux2_bml, k, threshold)
293 call prg_conjgrad(y_bml, ai_bml, i_bml, aux1_bml, d_bml, w_bml, 0.9_dp, threshold)
296 call bml_multiply(ai_bml, p2_bml, p_bml, 1.0_dp, 0.0_dp, threshold)
300 trdpdmu = bml_trace(p_bml)
302 trdpdmu = trdpdmu - bml_sum_squares(p_bml)
303 trdpdmu = beta * trdpdmu
304 occerr = abs(trp0 - nocc)
305 if (occerr .gt. occerrlimit)
then
306 mu = mu + (nocc - trp0)/trdpdmu
308 write(*,*)
"mu =", mu
313 call bml_copy(p_bml, d_bml)
314 call bml_scale_add_identity(d_bml, -1.0_dp, 1.0_dp, threshold)
316 call bml_multiply(p_bml, d_bml, w_bml, 1.0_dp, 0.0_dp, threshold)
317 ofactor = ((nocc - trp0)/trdpdmu) * beta
323 call bml_deallocate(p2_bml)
324 call bml_deallocate(w_bml)
325 call bml_deallocate(d_bml)
326 call bml_deallocate(y_bml)
328 call bml_deallocate(aux1_bml)
329 call bml_deallocate(aux2_bml)
331 if (method .eq. 1)
then
332 call bml_deallocate(ai_bml)
333 call bml_deallocate(i_bml)
351 type(bml_matrix_t),
intent(in) :: h_bml
352 type(bml_matrix_t),
intent(inout) :: p_bml
353 integer,
intent(in) :: nsteps, method
354 real(
dp),
intent(in) :: mu, threshold
355 real(
dp),
intent(inout),
optional :: tol
357 type(bml_matrix_t) :: w_bml, y_bml, c_bml, d_bml, p2_bml, aux1_bml, aux2_bml, i_bml, ai_bml
359 real(
dp),
allocatable :: trace(:), gbnd(:)
360 character(20) :: bml_type
363 bml_type = bml_get_type(h_bml)
369 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, p2_bml)
370 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, d_bml)
371 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, w_bml)
372 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, y_bml)
373 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, c_bml)
374 if (method .eq. 1)
then
375 call bml_identity_matrix(bml_type, bml_element_real,
dp, n, m, i_bml)
376 call bml_identity_matrix(bml_type, bml_element_real,
dp, n, m, ai_bml)
379 call bml_copy(h_bml, p_bml)
380 call bml_gershgorin(p_bml, gbnd)
381 cnst = 0.5*min(1/(mu-gbnd(1)),1/(gbnd(2)-mu))
384 if (method .eq. 0)
then
385 write(*,*)
"Doing CG"
388 call bml_multiply_x2(p_bml, p2_bml, threshold, trace)
391 call bml_copy(p2_bml, y_bml)
392 call bml_add(y_bml, p_bml, 1.0_dp, -1.0_dp, threshold)
393 call bml_scale_add_identity(y_bml, 2.0_dp, 1.0_dp, threshold)
394 call prg_conjgrad(y_bml, p_bml, p2_bml, d_bml, w_bml, c_bml, tol, threshold)
397 write(*,*)
"Doing NS"
401 call bml_copy(p2_bml, y_bml)
402 call bml_add(y_bml, p_bml, 1.0_dp, -1.0_dp, threshold)
403 call bml_scale_add_identity(y_bml, 2.0_dp, 1.0_dp, threshold)
405 call prg_conjgrad(y_bml, ai_bml, i_bml, c_bml, d_bml, w_bml, 0.9_dp, threshold)
408 call bml_multiply(ai_bml, p2_bml, p_bml, 1.0_dp, 0.0_dp, threshold)
415 call bml_deallocate(p2_bml)
416 call bml_deallocate(w_bml)
417 call bml_deallocate(d_bml)
418 call bml_deallocate(y_bml)
419 call bml_deallocate(c_bml)
420 if (method .eq. 1)
then
421 call bml_deallocate(ai_bml)
422 call bml_deallocate(i_bml)
440 Inv_bml, nsteps, mu0, beta, nocc, threshold)
444 type(bml_matrix_t),
intent(in) :: h0_bml, h1_bml, inv_bml(nsteps)
445 type(bml_matrix_t),
intent(inout) :: p0_bml,p1_bml
446 real(
dp),
intent(in) :: mu0, threshold
448 real(
dp),
intent(in) :: beta, nocc
449 integer,
intent(in) :: nsteps
450 type(bml_matrix_t) :: b0_bml, b_bml, c_bml, c0_bml
451 character(20) :: bml_type
452 real(
dp) :: p1_trace, dpdmu_trace, p1b_trace, mu1b, cnst
453 integer :: n, m, i, j, k
455 bml_type = bml_get_type(h0_bml)
456 n = bml_get_n(h0_bml)
457 m = bml_get_m(h0_bml)
459 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, b0_bml)
460 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, b_bml)
461 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, c_bml)
462 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, c0_bml)
464 cnst = beta/(2**(2+nsteps))
467 call bml_copy(h0_bml, p0_bml)
471 call bml_copy(h1_bml, p1_bml)
472 call bml_scale(-1.0_dp*cnst, p1_bml)
477 call bml_multiply(p0_bml, p0_bml, c0_bml, 1.0_dp, 0.0_dp, threshold)
479 call bml_multiply(p0_bml, p1_bml, c_bml, 1.0_dp, 0.0_dp, threshold)
480 call bml_multiply(p1_bml, p0_bml, c_bml, 1.0_dp, 1.0_dp, threshold)
481 call bml_copy(p1_bml, b_bml)
482 call bml_add(b_bml, c_bml, 2.0_dp, -2.0_dp, threshold)
484 call bml_multiply(inv_bml(i), c0_bml, p0_bml, 1.0_dp, 0.0_dp, threshold)
487 call bml_multiply(b_bml, p0_bml, c_bml, 1.0_dp, 1.0_dp, threshold)
488 call bml_multiply(inv_bml(i), c_bml, p1_bml, 1.0_dp, 0.0_dp, threshold)
513 call bml_copy(p0_bml, b_bml)
514 call bml_scale_add_identity(b_bml, -1.0_dp, 1.0_dp, threshold)
515 call bml_multiply(p0_bml, b_bml, c_bml, beta, 0.0_dp, threshold)
516 dpdmu_trace = bml_trace(c_bml)
517 p1_trace = bml_trace(p1_bml)
518 mu1 = - p1_trace/dpdmu_trace
519 if (abs(dpdmu_trace) > 1e-8)
then
520 call bml_add(p1_bml,c_bml,1.0_dp,mu1,threshold)
523 call bml_deallocate(b_bml)
524 call bml_deallocate(b0_bml)
525 call bml_deallocate(c_bml)
526 call bml_deallocate(c0_bml)
548 nsteps, mu0, mu, beta, nocc, occ_tol, lin_tol, order, threshold)
552 type(bml_matrix_t),
intent(in) :: h0_bml, h1_bml, h2_bml, h3_bml
553 type(bml_matrix_t),
intent(inout) :: p0_bml, p1_bml, p2_bml, p3_bml
554 real(
dp),
intent(inout) :: mu0
555 real(
dp),
allocatable,
intent(inout) :: mu(:)
556 real(
dp),
intent(in) :: beta, occ_tol, lin_tol, nocc
557 integer,
intent(in) :: nsteps
558 type(bml_matrix_t) :: i_bml, tmp1_bml, tmp2_bml, tmp3_bml, c0_bml, t_bml, ti_bml
559 type(bml_matrix_t),
allocatable :: b_bml(:), p_bml(:), c_bml(:), h_bml(:)
560 real(
dp),
allocatable :: p_trace(:), trace(:)
561 character(20) :: bml_type
562 real(
dp) :: occ_err, p0_trace, pmu_trace, cnst, threshold, tol, lambda, h
563 integer :: n, m, order, i, j, k
567 allocate(p_trace(order))
568 allocate(b_bml(order))
569 allocate(c_bml(order))
570 allocate(p_bml(order))
571 allocate(h_bml(order))
573 bml_type = bml_get_type(h0_bml)
574 n = bml_get_n(h0_bml)
575 m = bml_get_m(h0_bml)
579 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, b_bml(i))
580 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, c_bml(i))
583 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, tmp1_bml)
584 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, tmp3_bml)
585 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, tmp2_bml)
586 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, c0_bml)
587 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, t_bml)
588 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, ti_bml)
589 call bml_identity_matrix(bml_type, bml_element_real,
dp, n, m, i_bml)
593 if (order .gt. 1)
then
597 if (order .gt. 2)
then
602 cnst = beta/(2**(2+nsteps))
604 do while (occ_err .gt. occ_tol)
607 call bml_copy(h0_bml, p0_bml)
612 call bml_copy(h_bml(j), p_bml(j))
614 call bml_scale_add_identity(p_bml(j), 1.0_dp, -0.5_dp, threshold)
622 call bml_multiply(p0_bml, p0_bml, c0_bml, 1.0_dp, 0.0_dp, threshold)
624 call bml_multiply(p0_bml, p_bml(1), c_bml(1), 1.0_dp, 0.0_dp, threshold)
625 call bml_multiply(p_bml(1), p0_bml, c_bml(1), 1.0_dp, 1.0_dp, threshold)
626 call bml_copy(p_bml(1), b_bml(1))
627 call bml_add(b_bml(1), c_bml(1), 2.0_dp, -2.0_dp, threshold)
630 call bml_multiply(p_bml(1), p_bml(1), c_bml(2), 1.0_dp, 0.0_dp, threshold)
631 call bml_multiply(p0_bml, p_bml(2), c_bml(2), 1.0_dp, 1.0_dp, threshold)
632 call bml_multiply(p_bml(2), p0_bml, c_bml(2), 1.0_dp, 1.0_dp, threshold)
633 call bml_copy(p_bml(2), b_bml(2))
634 call bml_add(b_bml(2), c_bml(2), 2.0_dp, -2.0_dp, threshold)
638 call bml_multiply(p_bml(1), p_bml(2), c_bml(3), 1.0_dp, 0.0_dp, threshold)
639 call bml_multiply(p_bml(2), p_bml(1), c_bml(3), 1.0_dp, 1.0_dp, threshold)
640 call bml_multiply(p0_bml, p_bml(3), c_bml(3), 1.0_dp, 1.0_dp, threshold)
641 call bml_multiply(p_bml(3), p0_bml, c_bml(3), 1.0_dp, 1.0_dp, threshold)
642 call bml_copy(p_bml(3), b_bml(3))
643 call bml_add(b_bml(3), c_bml(3), 2.0_dp, -2.0_dp, threshold)
646 call bml_copy(c0_bml, t_bml)
647 call bml_add(t_bml, p0_bml, 1.0_dp, -1.0_dp, threshold)
648 call bml_scale_add_identity(t_bml, 2.0_dp, 1.0_dp, threshold)
651 call prg_conjgrad(t_bml, ti_bml, i_bml, tmp1_bml, tmp2_bml, tmp3_bml,0.01_dp, threshold)
652 call bml_identity_matrix(bml_type, bml_element_real,
dp, n, m, i_bml)
656 call bml_multiply(ti_bml, c0_bml, p0_bml, 1.0_dp, 0.0_dp, threshold)
658 call bml_multiply(b_bml(1), p0_bml, c_bml(1), 1.0_dp, 1.0_dp, threshold)
659 call bml_multiply(ti_bml, c_bml(1), p_bml(1), 1.0_dp, 0.0_dp, threshold)
662 call bml_multiply(b_bml(2), p0_bml, c_bml(2), 1.0_dp, 1.0_dp, threshold)
663 call bml_multiply(b_bml(1), p_bml(1), c_bml(2), 1.0_dp, 1.0_dp, threshold)
664 call bml_multiply(ti_bml, c_bml(2), p_bml(2), 1.0_dp, 0.0_dp, threshold)
668 call bml_multiply(b_bml(3), p0_bml, c_bml(3), 1.0_dp, 1.0_dp, threshold)
669 call bml_multiply(b_bml(2), p_bml(1), c_bml(3), 1.0_dp, 1.0_dp, threshold)
670 call bml_multiply(b_bml(1), p_bml(2), c_bml(3), 1.0_dp, 1.0_dp, threshold)
671 call bml_multiply(ti_bml, c_bml(3), p_bml(3), 1.0_dp, 0.0_dp, threshold)
676 call bml_copy(p0_bml, tmp1_bml)
677 call bml_scale_add_identity(tmp1_bml, -1.0_dp, 1.0_dp, threshold)
678 call bml_multiply(p0_bml, tmp1_bml, tmp2_bml, beta, 0.0_dp, threshold)
680 pmu_trace = bml_trace(tmp2_bml)
681 p0_trace = bml_trace(p0_bml)
682 occ_err = abs(p0_trace-nocc)
683 mu0 = mu0 + (nocc - p0_trace)/pmu_trace
685 p_trace(i) = bml_trace(p_bml(i))
686 mu(i) = mu(i) - p_trace(i)/pmu_trace
687 occ_err = occ_err + abs(p_trace(i))
690 write(*,*)
"occ_err =", occ_err
692 write(*,*)
"Chemical potential is not converging"
699 call bml_deallocate(b_bml(i))
700 call bml_deallocate(c_bml(i))
703 call bml_deallocate(i_bml)
704 call bml_deallocate(tmp1_bml)
705 call bml_deallocate(tmp2_bml)
706 call bml_deallocate(t_bml)
707 call bml_deallocate(ti_bml)
726 subroutine prg_finite_diff(H0_bml, H_list, mu0, mu_list, beta, order, lambda, h, threshold)
730 type(bml_matrix_t),
intent(in) :: h0_bml
731 real(
dp),
intent(in) :: mu0
732 type(bml_matrix_t),
allocatable,
intent(in) :: h_list(:)
733 real(
dp),
allocatable,
intent(in) :: mu_list(:)
734 real(
dp),
intent(in) :: lambda, beta, threshold, h
735 integer,
intent(in) :: order
736 character(20) :: bml_type
737 real(
dp) :: mu_1minus, mu_1plus, mu_2minus, mu_2plus, mu_central
738 real(
dp) :: lambda_f, lambda_b, lambda_2f, lambda_2b
740 type(bml_matrix_t) :: d0_bml, d1minus_bml, d1plus_bml, d2plus_bml, d2minus_bml, tmp1_bml
742 bml_type = bml_get_type(h0_bml)
743 n = bml_get_n(h0_bml)
744 m = bml_get_m(h0_bml)
746 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, tmp1_bml)
747 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, d0_bml)
748 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, d1minus_bml)
749 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, d1plus_bml)
750 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, d2plus_bml)
751 call bml_zero_matrix(bml_type, bml_element_real,
dp, n, m, d2minus_bml)
755 lambda_2f = lambda+2*h
756 lambda_2b = lambda-2*h
761 call bml_copy(h0_bml, tmp1_bml)
764 call bml_add(tmp1_bml, h_list(i), 1.0_dp, lambda**i, threshold)
765 mu_central = mu_central + lambda**i*mu_list(i)
769 call bml_copy(h0_bml, tmp1_bml)
772 call bml_add(tmp1_bml, h_list(i), 1.0_dp, lambda_f**i, threshold)
773 mu_1plus = mu_1plus + lambda_f**i*mu_list(i)
776 call bml_copy(d1plus_bml, tmp1_bml)
777 call bml_add(tmp1_bml, d0_bml, 1.0_dp/h, -1.0_dp/h, threshold)
779 call bml_scale(1000.0_dp, tmp1_bml)
780 call bml_print_matrix(
"Finite diff - Order 1 * 1000", tmp1_bml, 0,10,0,10)
782 if (order .gt. 1)
then
784 call bml_copy(h0_bml, tmp1_bml)
787 call bml_add(tmp1_bml, h_list(i), 1.0_dp, lambda_b**i, threshold)
788 mu_1minus = mu_1minus + lambda_b**i*mu_list(i)
791 call bml_copy(d0_bml, tmp1_bml)
792 call bml_add(tmp1_bml, d1minus_bml, -1.0_dp/(h*h), 0.5_dp/(h*h), threshold)
793 call bml_add(tmp1_bml, d1plus_bml, 1.0_dp, 0.5_dp/(h*h))
795 call bml_scale(1000.0_dp, tmp1_bml)
796 call bml_print_matrix(
"Finite diff - Order 2 * 1000", tmp1_bml, 0,10,0,10)
799 if (order .gt. 2)
then
801 call bml_copy(h0_bml, tmp1_bml)
804 call bml_add(tmp1_bml, h_list(i), 1.0_dp, lambda_2b**i, threshold)
805 mu_2minus = mu_2minus + lambda_2b**i*mu_list(i)
808 call bml_copy(h0_bml, tmp1_bml)
811 call bml_add(tmp1_bml, h_list(i), 1.0_dp, lambda_2f**i, threshold)
812 mu_2plus = mu_2plus + lambda_2f**i*mu_list(i)
815 call bml_copy(d2plus_bml, tmp1_bml)
816 call bml_add(tmp1_bml, d2minus_bml, 1.0_dp/(12.0_dp*h**3), -1.0_dp/(12.0_dp*h**3), threshold)
817 call bml_add(tmp1_bml, d1plus_bml, 1.0_dp, -1.0/(6.0_dp*h**3), threshold)
818 call bml_add(tmp1_bml, d1minus_bml, 1.0_dp, 1.0/(6.0_dp*h**3), threshold)
820 call bml_scale(1000.0_dp, tmp1_bml)
821 call bml_print_matrix(
"Finite diff - Order 3 * 1000", tmp1_bml, 0,10,0,10)
824 call bml_deallocate(tmp1_bml)
825 call bml_deallocate(d0_bml)
826 call bml_deallocate(d1minus_bml)
827 call bml_deallocate(d1plus_bml)
828 call bml_deallocate(d2plus_bml)
829 call bml_deallocate(d2minus_bml)
844 aux1_bml, k, threshold)
848 type(bml_matrix_t),
intent(inout) :: A_bml, b_bml, p2_bml, y_bml, aux_bml, aux1_bml
849 type(bml_matrix_t),
intent(in) :: p_bml
850 real(dp),
intent(in) :: threshold
851 integer,
intent(in) :: k
852 character(20) :: bml_type
856 call bml_multiply(p_bml, p_bml, b_bml, 1.0_dp, 0.0_dp, threshold)
857 call bml_copy(b_bml, a_bml)
858 call bml_add(a_bml, p_bml, 2.0_dp, -2.0_dp, threshold)
859 call bml_scale_add_identity(a_bml, 1.0_dp, 1.0_dp, threshold)
862 call bml_multiply(p_bml, p_bml, p2_bml, 1.0_dp, 0.0_dp, threshold)
863 call bml_copy(p2_bml, y_bml)
864 call bml_add(y_bml, p_bml, 1.0_dp, -2.0_dp, threshold)
865 call bml_scale_add_identity(y_bml, 1.0_dp, 1.0_dp, threshold)
867 call bml_copy(p2_bml, b_bml)
868 call bml_copy(y_bml, a_bml)
870 call bml_multiply(b_bml, p2_bml, aux_bml, 1.0_dp, 0.0_dp, threshold)
871 call bml_multiply(a_bml, y_bml, aux1_bml, 1.0_dp, 0.0_dp, threshold)
872 call bml_copy(aux_bml, b_bml)
873 call bml_copy(aux1_bml, a_bml)
875 call bml_add(a_bml, b_bml, 1.0_dp, 1.0_dp, threshold)
888 subroutine prg_newtonschulz(a_bml, ai_bml, r_bml, tmp_bml, d_bml, I_bml, tol, threshold, num_iter)
892 type(bml_matrix_t),
intent(inout) :: ai_bml, r_bml, tmp_bml, d_bml
893 type(bml_matrix_t),
intent(in) :: a_bml, I_bml
894 real(dp),
intent(in) :: threshold, tol
895 integer,
intent(out) :: num_iter
896 real(dp) :: err,prev_err,scaled_tol
905 call bml_copy(ai_bml, tmp_bml)
906 call bml_multiply(a_bml, ai_bml, r_bml, 1.0_dp, 0.0_dp, threshold)
907 call bml_scale_add_identity(r_bml, -1.0_dp, 1.0_dp, threshold)
909 err = bml_fnorm(r_bml)
912 if (prev_err+0.01 < err)
then
913 write(*,*)
'NS did not converge, calling conjugate gradient'
914 call prg_conjgrad(a_bml, ai_bml, i_bml, r_bml, tmp_bml, d_bml, 0.00001_dp, threshold)
916 call bml_multiply(tmp_bml, r_bml, ai_bml, 1.0_dp, 1.0_dp, threshold)
934 subroutine prg_pcg(A_bml, p_bml, p2_bml, d_bml, wtmp_bml, cg_tol, threshold)
938 type(bml_matrix_t),
intent(in) :: A_bml
939 type(bml_matrix_t),
intent(inout) :: p_bml, p2_bml, d_bml, wtmp_bml
940 real(dp),
intent(in) :: cg_tol, threshold
942 type(bml_matrix_t) :: M_bml, z_bml
943 real(dp),
allocatable :: diagonal(:)
944 real(dp) :: alpha, beta
945 character(20) :: bml_type
947 real(dp) :: r_norm_old, r_norm_new
949 bml_type = bml_get_type(p_bml)
953 call bml_zero_matrix(bml_type, bml_element_real, dp, n, m, z_bml)
954 call bml_zero_matrix(bml_type, bml_element_real, dp, n, m, m_bml)
956 allocate(diagonal(n))
957 call bml_get_diagonal(a_bml, diagonal)
959 diagonal(k) = 1.0_dp/diagonal(k)
961 call bml_set_diagonal(m_bml, diagonal)
963 call bml_multiply(a_bml, p_bml, p2_bml, -1.0_dp, 1.0_dp, threshold)
964 call bml_multiply(m_bml, p2_bml, z_bml, 1.0_dp, 0.0_dp, threshold)
965 r_norm_new = bml_trace_mult(z_bml, p2_bml)
966 call bml_copy(z_bml, d_bml)
969 do while (bml_sum_squares(p2_bml) .gt. cg_tol)
971 write(*,*)
"r_norm", bml_sum_squares(p2_bml)
974 beta = r_norm_new/r_norm_old
975 call bml_add(d_bml, z_bml, beta, 1.0_dp, threshold)
978 call bml_multiply(a_bml, d_bml, wtmp_bml, 1.0_dp, 0.0_dp, threshold)
979 alpha = bml_trace_mult(p2_bml,z_bml)/bml_trace_mult(d_bml, wtmp_bml)
980 call bml_add(p_bml, d_bml, 1.0_dp, alpha, threshold)
981 call bml_add(p2_bml, wtmp_bml, 1.0_dp, -alpha, threshold)
982 call bml_multiply(m_bml, p2_bml, z_bml, 1.0_dp, 0.0_dp, threshold)
983 r_norm_old = r_norm_new
984 r_norm_new = bml_trace_mult(p2_bml,z_bml)
986 write(*,*)
"PCG is not converging"
990 write(*,*)
"Number of iterations:", k
992 call bml_deallocate(z_bml)
993 call bml_deallocate(m_bml)
1006 subroutine prg_conjgrad(A_bml, p_bml, p2_bml, tmp_bml, d_bml, w_bml, cg_tol, threshold)
1010 type(bml_matrix_t),
intent(in) :: A_bml, p2_bml
1011 type(bml_matrix_t),
intent(inout) :: p_bml, tmp_bml, d_bml, w_bml
1012 real(dp),
intent(in) :: cg_tol, threshold
1014 real(dp) :: alpha, beta
1016 real(dp) :: r_norm_old, r_norm_new
1018 call bml_copy(p2_bml,tmp_bml)
1019 call bml_multiply(a_bml, p_bml, tmp_bml, -1.0_dp, 1.0_dp, threshold)
1020 r_norm_new = bml_sum_squares(tmp_bml)
1023 do while (r_norm_new .gt. cg_tol)
1028 call bml_copy(tmp_bml, d_bml)
1030 beta = r_norm_new/r_norm_old
1031 call bml_add(d_bml, tmp_bml, beta, 1.0_dp, threshold)
1034 call bml_multiply(a_bml, d_bml, w_bml, 1.0_dp, 0.0_dp, threshold)
1035 alpha = r_norm_new/bml_trace_mult(d_bml, w_bml)
1037 call bml_add(p_bml, d_bml, 1.0_dp, alpha, threshold)
1038 call bml_add(tmp_bml, w_bml, 1.0_dp, -alpha, threshold)
1039 r_norm_old = r_norm_new
1040 r_norm_new = bml_sum_squares(tmp_bml)
1042 write(*,*)
"Conjugate gradient is not converging"
1060 type(bml_matrix_t),
intent(in) :: ham_bml
1061 type(bml_matrix_t),
intent(inout) :: p_bml
1062 real(dp),
intent(in) :: beta, threshold
1063 real(dp),
intent(in) :: mu
1064 character(20) :: bml_type
1066 real(dp),
allocatable :: eigenvalues(:)
1067 type(bml_matrix_t) :: eigenvectors_bml,occupation_bml,aux_bml,aux1_bml,i_bml
1069 bml_type = bml_get_type(p_bml)
1070 n = bml_get_n(p_bml)
1071 m = bml_get_m(p_bml)
1073 allocate(eigenvalues(n))
1075 call bml_zero_matrix(bml_type,bml_element_real,dp,n,m,eigenvectors_bml)
1076 call bml_zero_matrix(bml_type,bml_element_real,dp,n,m,occupation_bml)
1077 call bml_zero_matrix(bml_type,bml_element_real,dp,n,m,aux_bml)
1078 call bml_zero_matrix(bml_type,bml_element_real,dp,n,m,aux1_bml)
1079 call bml_identity_matrix(bml_type,bml_element_real,dp,n,m,i_bml)
1081 call bml_diagonalize(ham_bml,eigenvalues,eigenvectors_bml)
1084 eigenvalues(i) =
fermi(eigenvalues(i),mu,beta)
1087 call bml_set_diagonal(occupation_bml, eigenvalues)
1088 call bml_multiply(eigenvectors_bml, occupation_bml, aux_bml, 1.0_dp, 0.0_dp,threshold)
1089 call bml_transpose(eigenvectors_bml, aux1_bml)
1090 call bml_multiply(aux_bml, aux1_bml, p_bml, 1.0_dp, 0.0_dp, threshold)
1092 call bml_deallocate(eigenvectors_bml)
1093 call bml_deallocate(occupation_bml)
1094 call bml_deallocate(aux_bml)
1095 call bml_deallocate(aux1_bml)
1096 call bml_deallocate(i_bml)
1098 deallocate(eigenvalues)
1116 type(bml_matrix_t),
intent(in) :: ham_bml
1117 type(bml_matrix_t),
intent(inout) :: p_bml
1118 real(
dp),
intent(in) :: beta, nocc, occerrlimit, threshold
1119 real(
dp),
intent(inout) :: mu
1120 integer,
intent(in) :: osteps
1121 character(20) :: bml_type
1122 integer :: n, m, i, iter
1123 real(
dp) :: trdpdmu, trp0, ofactor, occerr
1124 real(
dp),
allocatable :: eigenvalues(:)
1125 type(bml_matrix_t) :: eigenvectors_bml,occupation_bml,aux_bml,aux1_bml,i_bml
1127 bml_type = bml_get_type(p_bml)
1128 n = bml_get_n(p_bml)
1129 m = bml_get_m(p_bml)
1131 allocate(eigenvalues(n))
1133 call bml_zero_matrix(bml_type,bml_element_real,
dp,n,m,eigenvectors_bml)
1134 call bml_zero_matrix(bml_type,bml_element_real,
dp,n,m,occupation_bml)
1135 call bml_zero_matrix(bml_type,bml_element_real,
dp,n,m,aux_bml)
1136 call bml_zero_matrix(bml_type,bml_element_real,
dp,n,m,aux1_bml)
1137 call bml_identity_matrix(bml_type,bml_element_real,
dp,n,m,i_bml)
1142 do while ((osteps .eq. 0 .and. occerr .gt. occerrlimit) .or. &
1143 (osteps .gt. 0 .and. iter .lt. osteps))
1146 call bml_diagonalize(ham_bml,eigenvalues,eigenvectors_bml)
1149 eigenvalues(i) =
fermi(eigenvalues(i),mu,beta)
1152 call bml_set_diagonal(occupation_bml, eigenvalues)
1153 call bml_multiply(eigenvectors_bml, occupation_bml, aux_bml, 1.0_dp, 0.0_dp,threshold)
1154 call bml_transpose(eigenvectors_bml, aux1_bml)
1155 call bml_multiply(aux_bml, aux1_bml, p_bml, 1.0_dp, 0.0_dp, threshold)
1156 call bml_print_matrix(
"test density",p_bml,0,10,0,10)
1157 trdpdmu = bml_trace(p_bml)
1159 trdpdmu = trdpdmu - bml_sum_squares(p_bml)
1160 trdpdmu = beta * trdpdmu
1161 occerr = abs(trp0 - nocc)
1162 if (occerr .gt. occerrlimit)
then
1163 mu = mu + (nocc - trp0)/trdpdmu
1170 call bml_copy(p_bml, aux_bml)
1171 call bml_scale_add_identity(aux_bml, -1.0_dp, 1.0_dp, threshold)
1173 call bml_multiply(p_bml, aux_bml, aux1_bml, 1.0_dp, 0.0_dp, threshold)
1174 ofactor = ((nocc - trp0)/trdpdmu) * beta
1178 call bml_deallocate(eigenvectors_bml)
1179 call bml_deallocate(occupation_bml)
1180 call bml_deallocate(aux_bml)
1181 call bml_deallocate(aux1_bml)
1182 call bml_deallocate(i_bml)
1184 deallocate(eigenvalues)
1194 real(
dp),
intent(in) :: e, mu, beta
1196 fermi = 1.0_dp/(1.0_dp+exp(beta*(e-mu)))
Module to obtain the density matrix by diagonalizing an orthogonalized Hamiltonian.
real(dp) function fermi(e, ef, kbt)
Gives the Fermi distribution value for energy e.
real(dp) function fermi(e, mu, beta)
Gives the Fermi distribution value for energy e.
subroutine, public prg_implicit_fermi(h_bml, p_bml, nsteps, k, nocc, mu, beta, method, osteps, occErrLimit, threshold, tol)
Recursive Implicit Fermi Dirac for finite temperature.
subroutine, public prg_implicit_fermi_save_inverse(Inv_bml, h_bml, p_bml, nsteps, nocc, mu, beta, occErrLimit, threshold, tol, SCF_IT, occiter, totns)
Recursive Implicit Fermi Dirac for finite temperature.
subroutine prg_pcg(A_bml, p_bml, p2_bml, d_bml, wtmp_bml, cg_tol, threshold)
Solve the system AX = B with conjugate gradient.
subroutine, public prg_implicit_fermi_first_order_response(H0_bml, H1_bml, P0_bml, P1_bml, Inv_bml, nsteps, mu0, beta, nocc, threshold)
Calculate first order density matrix response to perturbations using Implicit Fermi Dirac.
subroutine prg_get_density_matrix(ham_bml, p_bml, beta, mu, threshold)
Calculate the density matrix with diagonalization.
subroutine, public prg_implicit_fermi_response(H0_bml, H1_bml, H2_bml, H3_bml, P0_bml, P1_bml, P2_bml, P3_bml, nsteps, mu0, mu, beta, nocc, occ_tol, lin_tol, order, threshold)
Calculate density matrix response to perturbations using Implicit Fermi Dirac.
subroutine prg_setup_linsys(p_bml, A_bml, b_bml, p2_bml, y_bml, aux_bml, aux1_bml, k, threshold)
Set up linear system for Implicit Fermi Dirac.
subroutine, public prg_finite_diff(H0_bml, H_list, mu0, mu_list, beta, order, lambda, h, threshold)
Calculate density matrix response from perturbations in the Hamiltonian.
subroutine, public prg_test_density_matrix(ham_bml, p_bml, beta, mu, nocc, osteps, occErrLimit, threshold)
Calculate the density matrix with diagonalization and converge chemical.
subroutine prg_conjgrad(A_bml, p_bml, p2_bml, tmp_bml, d_bml, w_bml, cg_tol, threshold)
Solve the system AX = B with conjugate gradient.
subroutine prg_newtonschulz(a_bml, ai_bml, r_bml, tmp_bml, d_bml, I_bml, tol, threshold, num_iter)
Find the inverse of the matrix A with Newton-Schulz iteration.
subroutine, public prg_implicit_fermi_zero(h_bml, p_bml, nsteps, mu, method, threshold, tol)
Recursive Implicit Fermi Dirac for zero temperature.
The prg_normalize module.
subroutine, public prg_normalize_implicit_fermi(h_bml, cnst, mu)
Normalize a Hamiltonian matrix prior to running the implicit fermi dirac algorithm.