Line data Source code
1 : MODULE m_broyden 2 : USE m_juDFT 3 : !################################################################ 4 : ! IMIX = 3 : BROYDEN'S FIRST METHOD 5 : ! IMIX = 5 : BROYDEN'S SECOND METHOD 6 : ! IMIX = 7 : GENERALIZED ANDERSEN METHOD 7 : ! sm : input charge density of iteration m 8 : ! afterwards update rho(m+1) 9 : ! fm : output minus input charge density of iteration m 10 : ! sm1 : input charge density of iteration m-1 11 : ! fm1 : output minus inputcharge density of iteration m-1 12 : !################################################################ 13 : CONTAINS 14 654 : SUBROUTINE broyden(alpha,fm,sm,l_dfpt) 15 : USE m_types 16 : USE m_types_mixvector 17 : IMPLICIT NONE 18 : 19 : real,INTENT(IN) :: alpha 20 : TYPE(t_mixvector),INTENT(IN) :: fm(:) 21 : TYPE(t_mixvector),INTENT(INOUT) :: sm(:) 22 : LOGICAL, INTENT(IN) :: l_dfpt 23 : 24 : ! Locals 25 : INTEGER :: n,it,hlen 26 : REAL :: fmvm,vmnorm 27 654 : REAL,ALLOCATABLE :: am(:),dfivi(:) 28 654 : TYPE(t_mixvector) :: fm1,sm1,ui,um,vi,vm 29 654 : TYPE(t_mixvector),allocatable :: u_store(:),v_store(:) 30 : 31 654 : hlen=size(fm) 32 654 : IF (hlen<2) THEN !Do a simple mixing step 33 144 : sm(hlen)=sm(hlen)+alpha*fm(hlen) 34 144 : RETURN 35 : ENDIF 36 : 37 8200 : ALLOCATE(u_store(hlen-2),v_store(hlen-2)) 38 3590 : do it=1,hlen-2 39 3080 : call u_store(it)%alloc() 40 3590 : call v_store(it)%alloc() 41 : enddo 42 : 43 510 : CALL fm1%alloc() 44 510 : CALL sm1%alloc() 45 510 : CALL ui%alloc() 46 510 : CALL um%alloc() 47 510 : CALL vi%alloc() 48 510 : CALL vm%alloc() 49 : 50 2040 : ALLOCATE (am(hlen-1),dfivi(hlen-1)) 51 4100 : dfivi = 0.0 52 4100 : am = 0.0 53 510 : call timestart("Broyden-loop") 54 4100 : DO n=2,hlen 55 3590 : sm1 = sm(n) - sm(n-1) 56 3590 : fm1 = fm(n) - fm(n-1) 57 : ! |vi> = w|vi> 58 : ! loop to generate um : um = sm1 + alpha*fm1 - \sum <fm1|w|vi> ui 59 3590 : um = alpha * fm1 + sm1 60 3590 : call timestart("Broyden-1.loop") 61 20824 : DO it = n-2,1,-1 62 17234 : ui=u_store(it) 63 17234 : vi=v_store(it) 64 : 65 17234 : am(it) = vi.dot.fm1 66 : ! calculate um(:) = -am(it)*ui(:) + um(:) 67 20824 : um=um-am(it)*ui 68 : !WRITE(oUnit,FMT='(5x,"<vi|w|Fm> for it",i2,5x,f10.6)')it,am(it) 69 : END DO 70 3590 : call timestop("Broyden-1.loop") 71 : 72 : ! calculate vm = alpha*wfm1 -\sum <fm1|w|vi> <fi1|w|vi><vi| 73 : ! convolute fm1 with the metrik and store in vm 74 3590 : vm=fm1%apply_metric(l_dfpt) 75 3590 : call timestart("Broyden-2.loop") 76 20824 : DO it = n-2,1,-1 77 17234 : vi=v_store(it) 78 : ! calculate vm(:) = -am(it)*dfivi*vi(:) + vm 79 20824 : vm=vm-am(it)*dfivi(it)*vi 80 : END DO 81 3590 : call timestop("Broyden-2.loop") 82 : 83 3590 : vmnorm=fm1.dot.vm 84 : ! if (vmnorm.lt.tol_10) stop 85 : 86 : ! calculate vm(:) = (1.0/vmnorm)*vm(:) 87 3590 : vm=(1.0/vmnorm)*vm 88 : 89 : ! save dfivi(mit) for next iteration 90 3590 : dfivi(n-1) = vmnorm 91 3590 : IF (n<hlen) u_store(n-1)=um 92 4100 : IF (n<hlen) v_store(n-1)=vm 93 : enddo 94 510 : call timestop("Broyden-loop") 95 : ! update rho(m+1) 96 : ! calculate <fm|w|vm> 97 510 : fmvm = vm.dot.fm(hlen) 98 : ! calculate sm(:) = (1.0-fmvm)*um(:) + sm 99 510 : sm(hlen)=sm(hlen)+(1.0-fmvm)*um 100 : 101 6814 : END SUBROUTINE broyden 102 : END MODULE m_broyden