(git:d5c4d39)
Loading...
Searching...
No Matches
tblite_scc_mixer.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2026 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
9!> \brief CP2K-side tblite-compatible SCC Broyden mixer.
10! **************************************************************************************************
12#if defined(__TBLITE)
13 USE mctc_env, ONLY: error_type, &
14 fatal_error, &
15 wp
16 USE tblite_lapack, ONLY: getrf, &
17 getrs
18 USE tblite_scf_mixer_type, ONLY: mixer_type
19#endif
20
21 IMPLICIT NONE
22
23 PRIVATE
24
25#if defined(__TBLITE)
26 PUBLIC :: new_cp2k_tblite_mixer
27
28 TYPE, EXTENDS(mixer_type) :: cp2k_tblite_broyden_mixer_type
29 INTEGER :: idif = 0
30 INTEGER :: iget = 0
31 INTEGER :: iset = 0
32 INTEGER :: iter = 0
33 INTEGER :: memory = 0
34 INTEGER :: ndim = 0
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
42 CONTAINS
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
49#endif
50
51CONTAINS
52
53#if defined(__TBLITE)
54! **************************************************************************************************
55!> \brief Create a CP2K-side tblite-compatible Broyden mixer.
56!> \param mixer mixer instance
57!> \param memory Broyden history length
58!> \param ndim number of SCC variables
59!> \param damping first-step damping
60!> \param omega0 Broyden regularization weight
61!> \param min_weight minimum dynamic Broyden weight
62!> \param max_weight maximum dynamic Broyden weight
63!> \param weight_factor residual-to-weight scaling factor
64! **************************************************************************************************
65 SUBROUTINE new_cp2k_tblite_mixer(mixer, memory, ndim, damping, omega0, min_weight, max_weight, &
66 weight_factor)
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, &
70 weight_factor
71
72 TYPE(cp2k_tblite_broyden_mixer_type), ALLOCATABLE :: broyden
73
74 ALLOCATE (broyden)
75 broyden%ndim = ndim
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)
91
92 END SUBROUTINE new_cp2k_tblite_mixer
93
94! **************************************************************************************************
95!> \brief Set input SCC variables.
96!> \param self mixer
97!> \param qvec SCC-variable vector
98! **************************************************************************************************
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
102
103 self%q_in(self%iset + 1:self%iset + SIZE(qvec)) = qvec
104 self%iset = self%iset + SIZE(qvec)
105
106 END SUBROUTINE cp2k_tblite_mixer_set_1d
107
108! **************************************************************************************************
109!> \brief Set SCC residual.
110!> \param self mixer
111!> \param qvec output SCC-variable vector before mixing
112! **************************************************************************************************
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
116
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)
120
121 END SUBROUTINE cp2k_tblite_mixer_diff_1d
122
123! **************************************************************************************************
124!> \brief Apply the next Broyden update.
125!> \param self mixer
126!> \param error error handling
127! **************************************************************************************************
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
131
132 INTEGER :: info
133
134 self%iset = 0
135 self%idif = 0
136 self%iget = 0
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")
143
144 END SUBROUTINE cp2k_tblite_mixer_next
145
146! **************************************************************************************************
147!> \brief Get mixed SCC variables.
148!> \param self mixer
149!> \param qvec mixed SCC-variable vector
150! **************************************************************************************************
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
154
155 qvec(:) = self%q_in(self%iget + 1:self%iget + SIZE(qvec))
156 self%iget = self%iget + SIZE(qvec)
157
158 END SUBROUTINE cp2k_tblite_mixer_get_1d
159
160! **************************************************************************************************
161!> \brief Get RMS SCC residual.
162!> \param self mixer
163!> \return residual
164! **************************************************************************************************
165 PURE FUNCTION cp2k_tblite_mixer_get_error(self) RESULT(error)
166 CLASS(cp2k_tblite_broyden_mixer_type), INTENT(IN) :: self
167 REAL(wp) :: error
168
169 error = sqrt(sum(self%dq**2)/SIZE(self%dq))
170
171 END FUNCTION cp2k_tblite_mixer_get_error
172
173! **************************************************************************************************
174!> \brief Modified Broyden update matching tblite's default algorithm, with explicit constants.
175!> \param n ...
176!> \param q ...
177!> \param qlast ...
178!> \param dq ...
179!> \param dqlast ...
180!> \param iter ...
181!> \param memory ...
182!> \param alpha ...
183!> \param omega0 ...
184!> \param minw ...
185!> \param maxw ...
186!> \param wfac ...
187!> \param omega ...
188!> \param df ...
189!> \param u ...
190!> \param a ...
191!> \param info ...
192! **************************************************************************************************
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
205
206 INTEGER :: i, it1, itn, j
207 REAL(wp) :: inv
208 REAL(wp), ALLOCATABLE, DIMENSION(:, :) :: beta, c
209
210 info = 0
211 itn = iter - 1
212 it1 = mod(itn - 1, memory) + 1
213
214 IF (iter == 1) THEN
215 dqlast(:) = dq
216 qlast(:) = q
217 q(:) = q + alpha*dq
218 RETURN
219 END IF
220
221 ALLOCATE (beta(min(memory, itn), min(memory, itn)), c(min(memory, itn), 1))
222
223 omega(it1) = sqrt(dot_product(dq, dq))
224 IF (omega(it1) > (wfac/maxw)) THEN
225 omega(it1) = wfac/omega(it1)
226 ELSE
227 omega(it1) = maxw
228 END IF
229 IF (omega(it1) < minw) omega(it1) = minw
230
231 df(:, it1) = dq - dqlast
232 inv = max(sqrt(dot_product(df(:, it1), df(:, it1))), epsilon(1.0_wp))
233 inv = 1.0_wp/inv
234 df(:, it1) = inv*df(:, it1)
235
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)
241 END DO
242
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
247 END DO
248
249 CALL cp2k_tblite_lineq(beta, c, info)
250 IF (info /= 0) RETURN
251
252 u(:, it1) = alpha*df(:, it1) + inv*(q - qlast)
253 dqlast(:) = dq
254 qlast(:) = q
255 q(:) = q + alpha*dq
256
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)
260 END DO
261
262 END SUBROUTINE cp2k_tblite_broyden
263
264! **************************************************************************************************
265!> \brief Solve the small Broyden linear system.
266!> \param a ...
267!> \param c ...
268!> \param info ...
269! **************************************************************************************************
270 SUBROUTINE cp2k_tblite_lineq(a, c, info)
271 REAL(wp), DIMENSION(:, :), INTENT(INOUT) :: a, c
272 INTEGER, INTENT(OUT) :: info
273
274 INTEGER, ALLOCATABLE, DIMENSION(:) :: ipiv
275
276 ALLOCATE (ipiv(SIZE(a, 1)))
277 CALL getrf(a, ipiv, info)
278 IF (info == 0) CALL getrs(a, c, ipiv, info, trans="t")
279
280 END SUBROUTINE cp2k_tblite_lineq
281#endif
282
283END MODULE tblite_scc_mixer
CP2K-side tblite-compatible SCC Broyden mixer.