13 USE mctc_env,
ONLY: error_type, &
16 USE tblite_lapack,
ONLY: getrf, &
18 USE tblite_scf_mixer_type,
ONLY: mixer_type
26 PUBLIC :: new_cp2k_tblite_mixer
28 TYPE,
EXTENDS(mixer_type) :: cp2k_tblite_broyden_mixer_type
35 REAL(wp) :: damp = 0.0_wp
36 REAL(wp) :: max_weight = 0.0_wp
37 REAL(wp) :: min_weight = 0.0_wp
38 REAL(wp) :: omega0 = 0.0_wp
39 REAL(wp) :: weight_factor = 0.0_wp
40 REAL(wp),
ALLOCATABLE,
DIMENSION(:, :) :: a, df, u
41 REAL(wp),
ALLOCATABLE,
DIMENSION(:) :: dq, dqlast, omega, q_in, qlast_in
43 PROCEDURE :: diff_1d => cp2k_tblite_mixer_diff_1d
44 PROCEDURE :: get_1d => cp2k_tblite_mixer_get_1d
45 PROCEDURE :: get_error => cp2k_tblite_mixer_get_error
46 PROCEDURE :: next => cp2k_tblite_mixer_next
47 PROCEDURE :: set_1d => cp2k_tblite_mixer_set_1d
48 END TYPE cp2k_tblite_broyden_mixer_type
65 SUBROUTINE new_cp2k_tblite_mixer(mixer, memory, ndim, damping, omega0, min_weight, max_weight, &
67 CLASS(mixer_type),
ALLOCATABLE,
INTENT(OUT) :: mixer
68 INTEGER,
INTENT(IN) :: memory, ndim
69 REAL(wp),
INTENT(IN) :: damping, omega0, min_weight, max_weight, &
72 TYPE(cp2k_tblite_broyden_mixer_type),
ALLOCATABLE :: broyden
76 broyden%memory = memory
77 broyden%damp = damping
78 broyden%omega0 = omega0
79 broyden%min_weight = min_weight
80 broyden%max_weight = max_weight
81 broyden%weight_factor = weight_factor
82 ALLOCATE (broyden%a(memory, memory))
83 ALLOCATE (broyden%df(ndim, memory))
84 ALLOCATE (broyden%u(ndim, memory))
85 ALLOCATE (broyden%dq(ndim))
86 ALLOCATE (broyden%dqlast(ndim))
87 ALLOCATE (broyden%omega(memory))
88 ALLOCATE (broyden%q_in(ndim))
89 ALLOCATE (broyden%qlast_in(ndim))
90 CALL move_alloc(broyden, mixer)
92 END SUBROUTINE new_cp2k_tblite_mixer
99 SUBROUTINE cp2k_tblite_mixer_set_1d(self, qvec)
100 CLASS(cp2k_tblite_broyden_mixer_type),
INTENT(INOUT) :: self
101 REAL(wp),
DIMENSION(:),
INTENT(IN) :: qvec
103 self%q_in(self%iset + 1:self%iset +
SIZE(qvec)) = qvec
104 self%iset = self%iset +
SIZE(qvec)
106 END SUBROUTINE cp2k_tblite_mixer_set_1d
113 SUBROUTINE cp2k_tblite_mixer_diff_1d(self, qvec)
114 CLASS(cp2k_tblite_broyden_mixer_type),
INTENT(INOUT) :: self
115 REAL(wp),
DIMENSION(:),
INTENT(IN) :: qvec
117 self%dq(self%idif + 1:self%idif +
SIZE(qvec)) = &
118 qvec - self%q_in(self%idif + 1:self%idif +
SIZE(qvec))
119 self%idif = self%idif +
SIZE(qvec)
121 END SUBROUTINE cp2k_tblite_mixer_diff_1d
128 SUBROUTINE cp2k_tblite_mixer_next(self, error)
129 CLASS(cp2k_tblite_broyden_mixer_type),
INTENT(INOUT) :: self
130 TYPE(error_type),
ALLOCATABLE,
INTENT(OUT) :: error
137 self%iter = self%iter + 1
138 CALL cp2k_tblite_broyden(self%ndim, self%q_in, self%qlast_in, self%dq, self%dqlast, &
139 self%iter, self%memory, self%damp, self%omega0, &
140 self%min_weight, self%max_weight, self%weight_factor, &
141 self%omega, self%df, self%u, self%a, info)
142 IF (info /= 0)
CALL fatal_error(error,
"Broyden mixing failed to obtain next iteration")
144 END SUBROUTINE cp2k_tblite_mixer_next
151 SUBROUTINE cp2k_tblite_mixer_get_1d(self, qvec)
152 CLASS(cp2k_tblite_broyden_mixer_type),
INTENT(INOUT) :: self
153 REAL(wp),
DIMENSION(:),
INTENT(OUT) :: qvec
155 qvec(:) = self%q_in(self%iget + 1:self%iget +
SIZE(qvec))
156 self%iget = self%iget +
SIZE(qvec)
158 END SUBROUTINE cp2k_tblite_mixer_get_1d
165 PURE FUNCTION cp2k_tblite_mixer_get_error(self)
RESULT(error)
166 CLASS(cp2k_tblite_broyden_mixer_type),
INTENT(IN) :: self
169 error = sqrt(sum(self%dq**2)/
SIZE(self%dq))
171 END FUNCTION cp2k_tblite_mixer_get_error
193 SUBROUTINE cp2k_tblite_broyden(n, q, qlast, dq, dqlast, iter, memory, alpha, omega0, minw, &
194 maxw, wfac, omega, df, u, a, info)
195 INTEGER,
INTENT(IN) :: n
196 REAL(wp),
DIMENSION(n),
INTENT(INOUT) :: q, qlast
197 REAL(wp),
DIMENSION(n),
INTENT(IN) :: dq
198 REAL(wp),
DIMENSION(n),
INTENT(INOUT) :: dqlast
199 INTEGER,
INTENT(IN) :: iter, memory
200 REAL(wp),
INTENT(IN) :: alpha, omega0, minw, maxw, wfac
201 REAL(wp),
DIMENSION(memory),
INTENT(INOUT) :: omega
202 REAL(wp),
DIMENSION(n, memory),
INTENT(INOUT) :: df, u
203 REAL(wp),
DIMENSION(memory, memory),
INTENT(INOUT) :: a
204 INTEGER,
INTENT(OUT) :: info
206 INTEGER :: i, it1, itn, j
208 REAL(wp),
ALLOCATABLE,
DIMENSION(:, :) :: beta, c
212 it1 = mod(itn - 1, memory) + 1
221 ALLOCATE (beta(min(memory, itn), min(memory, itn)), c(min(memory, itn), 1))
223 omega(it1) = sqrt(dot_product(dq, dq))
224 IF (omega(it1) > (wfac/maxw))
THEN
225 omega(it1) = wfac/omega(it1)
229 IF (omega(it1) < minw) omega(it1) = minw
231 df(:, it1) = dq - dqlast
232 inv = max(sqrt(dot_product(df(:, it1), df(:, it1))), epsilon(1.0_wp))
234 df(:, it1) = inv*df(:, it1)
236 DO j = max(1, itn - memory + 1), itn
237 i = mod(j - 1, memory) + 1
238 a(i, it1) = dot_product(df(:, i), df(:, it1))
239 a(it1, i) = a(i, it1)
240 c(i, 1) = omega(i)*dot_product(df(:, i), dq)
243 DO j = max(1, itn - memory + 1), itn
244 i = mod(j - 1, memory) + 1
245 beta(:it1, i) = omega(:it1)*omega(i)*a(:it1, i)
246 beta(i, i) = beta(i, i) + omega0*omega0
249 CALL cp2k_tblite_lineq(beta, c, info)
250 IF (info /= 0)
RETURN
252 u(:, it1) = alpha*df(:, it1) + inv*(q - qlast)
257 DO j = max(1, itn - memory + 1), itn
258 i = mod(j - 1, memory) + 1
259 q(:) = q - omega(i)*c(i, 1)*u(:, i)
262 END SUBROUTINE cp2k_tblite_broyden
270 SUBROUTINE cp2k_tblite_lineq(a, c, info)
271 REAL(wp),
DIMENSION(:, :),
INTENT(INOUT) :: a, c
272 INTEGER,
INTENT(OUT) :: info
274 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: ipiv
276 ALLOCATE (ipiv(
SIZE(a, 1)))
277 CALL getrf(a, ipiv, info)
278 IF (info == 0)
CALL getrs(a, c, ipiv, info, trans=
"t")
280 END SUBROUTINE cp2k_tblite_lineq
CP2K-side tblite-compatible SCC Broyden mixer.