Octopus
mix.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2014 M. Marques, A. Castro, A. Rubio, G. Bertsch, M. Oliveira
2!! Copyright (C) 2024 A. Buccheri
3!!
4!! This program is free software; you can redistribute it and/or modify
5!! it under the terms of the GNU General Public License as published by
6!! the Free Software Foundation; either version 2, or (at your option)
7!! any later version.
8!!
9!! This program is distributed in the hope that it will be useful,
10!! but WITHOUT ANY WARRANTY; without even the implied warranty of
11!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12!! GNU General Public License for more details.
13!!
14!! You should have received a copy of the GNU General Public License
15!! along with this program; if not, write to the Free Software
16!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
17!! 02110-1301, USA.
18!!
19
20#include "global.h"
21
22module mix_oct_m
23 use comm_oct_m
24 use debug_oct_m
26 use global_oct_m
27 use, intrinsic :: iso_fortran_env
30 use math_oct_m
31 use mesh_oct_m
36 use parser_oct_m
41 use space_oct_m
43 use types_oct_m
45
46 implicit none
47
48 private
49 public :: &
50 mix_t, &
51 mix_init, &
52 mix_clear, &
53 mix_end, &
54 mix_dump, &
55 mix_load, &
56 mixing, &
57 dmixing, &
58 zmixing, &
60 mix_scheme, &
61 mix_d3, &
63 mixfield_t, &
71
72 ! A better design would be to have mixfield_t and mixfield_ptr_t separate from mix_t:
73 ! mix_t would contain constants, and should be immutable with only getters.
74 ! mixfield_type encapsulates the data used in mixing, and has bound procedures for
75 ! operating on this data.
76
78 type mixfield_t
79 private
80 ! Real versions of the data
81 real(real64), allocatable :: ddf(:, :, :)
82 real(real64), allocatable :: ddv(:, :, :)
83 real(real64), allocatable :: df_old(:, :)
84 real(real64), allocatable :: dvin_old(:, :)
85 real(real64), allocatable :: dvin(:, :)
86 real(real64), allocatable :: dvout(:, :)
87 real(real64), allocatable :: dvnew(:, :)
88 real(real64), allocatable :: dresidual(:, :)
89
90 ! Complex versions of the data
91 complex(real64), allocatable :: zdf(:, :, :)
92 complex(real64), allocatable :: zdv(:, :, :)
93 complex(real64), allocatable :: zf_old(:, :)
94 complex(real64), allocatable :: zvin_old(:, :)
95 complex(real64), allocatable :: zvin(:, :)
96 complex(real64), allocatable :: zvout(:, :)
97 complex(real64), allocatable :: zvnew(:, :)
98 complex(real64), allocatable :: zresidual(:, :)
99
100 type(type_t) :: func_type
101 integer :: d1, d2, d3
102 ! d1 is mesh%np, d2 vary depending on where this called
103 ! d3 is for SCF history of a quantity.
104 logical :: mix_spin_density_matrix
105 real(real64) :: last_residue
106 end type mixfield_t
107
108 type mixfield_ptr_t
109 private
110 type(mixfield_t), pointer :: p
111 end type mixfield_ptr_t
112
113 integer, parameter :: MAX_AUXMIXFIELD = 5
114
122 type mix_t
123 private
124 type(derivatives_t), pointer :: der
125
126 integer :: scheme
127 real(real64), allocatable :: coeff(:)
128 real(real64) :: adaptive_damping
129 real(real64) :: residual_coeff
130 integer :: iter
131 ! !< include the iterations done in previous calculations.
132 integer, public :: ns
133 integer, public :: ns_restart
134 integer :: interval
135
136 integer :: last_ipos
137 ! This should have been an attribute of type(mixfield_t)
138
139 ! Fields
140 type(mixfield_t) :: mixfield
141 integer :: nauxmixfield
142 type(mixfield_ptr_t) :: auxmixfield(MAX_AUXMIXFIELD)
143 integer :: ipos
144 ! This should have been an attribute of type(mixfield_ptr_t)
145
146 ! Preconditioning
147 logical :: kerker
148 real(real64) :: kerker_factor
149 logical :: precondition
150 type(nl_operator_t) :: preconditioner
151 ! [GPAW](https://wiki.fysik.dtu.dk/gpaw/documentation/densitymix/densitymix.html)
152
153 logical :: mix_spin_density_matrix
154 contains
155 procedure :: compute_residuals_aux_field
156
157 end type mix_t
158
159 interface mixfield_set_vin
160 module procedure dmixfield_set_vin, zmixfield_set_vin
161 end interface mixfield_set_vin
162
163 interface mixfield_set_vout
164 module procedure dmixfield_set_vout, zmixfield_set_vout
165 end interface mixfield_set_vout
166
167 interface mixfield_get_vnew
168 module procedure dmixfield_get_vnew, zmixfield_get_vnew
169 end interface mixfield_get_vnew
170
171contains
172
174 subroutine mix_init(smix, namespace, space, der, d1, d2, def_, func_type_, prefix_)
175 type(mix_t), intent(out) :: smix
176 type(namespace_t), intent(in) :: namespace
177 class(space_t), intent(in) :: space
178 type(derivatives_t), target, intent(in) :: der
179 integer, intent(in) :: d1
180 integer, intent(in) :: d2
181 integer, optional, intent(in) :: def_
182 type(type_t), optional, intent(in) :: func_type_
183 character(len=*), optional, intent(in) :: prefix_
184
185 integer :: def, ii
186 character(len=32) :: prefix
187 type(type_t) :: func_type
188 real(real64) :: coeff
190 push_sub(mix_init)
192 smix%der => der
194 def = option__mixingscheme__broyden
195 if (present(def_)) def = def_
196 if (present(func_type_)) then
197 func_type = func_type_
198 else
199 func_type = type_float
200 end if
201 prefix = ''
202 if (present(prefix_)) prefix = prefix_
204 call messages_obsolete_variable(namespace, 'TypeOfMixing', 'MixingScheme')
206 !%Variable MixingScheme
207 !%Type integer
208 !%Default broyden
209 !%Section SCF::Mixing
210 !%Description
211 !% The scheme used to produce, at each iteration in the self-consistent cycle
212 !% that attempts to solve the Kohn-Sham equations, the input density from the value
213 !% of the input and output densities of previous iterations.
214 !%Option linear 0
215 !% Simple linear mixing.
216 !%Option broyden 2
217 !% Broyden scheme [C. G Broyden, <i>Math. Comp.</i> <b>19</b>, 577 (1965);
218 !% D. D. Johnson, <i>Phys. Rev. B</i> <b>38</b>, 12807 (1988)].
219 !% The scheme is slightly adapted, see the comments in the code.
220 !% For complex functions (e.g. Sternheimer with <tt>EMEta</tt> > 0), we use the generalization
221 !% with a complex dot product.
222 !%Option diis 9
223 !% Direct inversion in the iterative subspace (diis)
224 !% scheme [P. Pulay, <i>Chem. Phys. Lett.</i>, <b>73</b>, 393
225 !% (1980)] as described in [G. Kresse, and J. Furthmueller,
226 !% <i>Phys. Rev. B</i> <b>54</b>, 11169 (1996)].
227 !%Option broyden_adaptive 11
228 !% Broyden mixing with adaptive damping. The Broyden update follows
229 !% [C. G. Broyden, <i>Math. Comp.</i> <b>19</b>, 577 (1965);
230 !% D. D. Johnson, <i>Phys. Rev. B</i> <b>38</b>, 12807 (1988)].
231 !% The adaptive damping logic is inspired by
232 !% [M. F. Herbst and A. Levitt, <i>J. Comput. Phys.</i> <b>459</b>, 111127 (2022)]
233 !% and by the adaptive-coefficient analysis of
234 !% [M. Novak et al., <i>Comput. Phys. Commun.</i> <b>292</b>, 108865 (2023)].
235 !%End
236 call parse_variable(namespace, trim(prefix)//'MixingScheme', def, smix%scheme)
237 if (.not. varinfo_valid_option('MixingScheme', smix%scheme)) then
238 call messages_input_error(namespace, 'MixingScheme', 'invalid option')
239 end if
240 call messages_print_var_option("MixingScheme", smix%scheme, namespace=namespace)
241
242 !%Variable MixingPreconditioner
243 !%Type logical
244 !%Default false
245 !%Section SCF::Mixing
246 !%Description
247 !% (Experimental) If set to yes, Octopus will use a preconditioner
248 !% for the mixing operator.
249 !% This preconditioner is disabled for systems with dimension other than 3.
250 !%End
251 call parse_variable(namespace, trim(prefix)+'MixingPreconditioner', .false., smix%precondition)
252 if (der%dim /= 3) smix%precondition = .false.
253 if (smix%precondition) call messages_experimental('MixingPreconditioner', namespace=namespace)
255 !%Variable Mixing
256 !%Type float
257 !%Default 0.3
258 !%Section SCF::Mixing
259 !%Description
260 !% The linear, Broyden, adaptive Broyden, and DIIS scheme depend on a "mixing parameter", set by this variable.
261 !% Must be 0 < <tt>Mixing</tt> <= 1.
262 !% For adaptive Broyden this is the maximum damping, which may be reduced automatically.
263 !%End
264 call parse_variable(namespace, trim(prefix)+'Mixing', 0.3_real64, coeff)
265 if (coeff <= m_zero .or. coeff > m_one) then
266 call messages_input_error(namespace, 'Mixing', 'Value should be positive and smaller than one')
267 end if
268 safe_allocate(smix%coeff(1:d2))
269 smix%coeff = coeff
270 smix%adaptive_damping = m_one
271
272 if (d2 > 1) then
273 !%Variable MixingSpinDensityMatrix
274 !%Type logical
275 !%Default false
276 !%Section SCF::Mixing
277 !%Description
278 !% If set to yes, Octopus mixes the components up, down, and real and imaginary parts of
279 !% the up down (for spinor calculations) of the density or potential.
280 !% By default, Octopus mixes the magnetization and the charge density with a different weight,
281 !% where the magnetization mixing coefficient is given by <tt>MixingMagnetization</tt>
282 !%End
283 call parse_variable(namespace, trim(prefix)+'MixingSpinDensityMatrix', .false., smix%mix_spin_density_matrix)
284
285 !%Variable MixingMagnetization
286 !%Type float
287 !%Default 1.5
288 !%Section SCF::Mixing
289 !%Description
290 !% Same as <tt>Mixing</tt>, but for the "magnetization".
291 !% The default value is taken from D. D. Johnson, <i>Phys. Rev. B</i> <b>38</b>, 12807 (1988).
292 !%End
293 call parse_variable(namespace, trim(prefix)+'MixingMagnetization', 1.5_real64, coeff)
294
295 if (.not. smix%mix_spin_density_matrix) smix%coeff(2:d2) = coeff
296 else
297 smix%mix_spin_density_matrix = .true.
298 end if
299
300 !%Variable MixingResidual
301 !%Type float
302 !%Default 0.05
303 !%Section SCF::Mixing
304 !%Description
305 !% In the DIIS mixing it is benefitial to include a bit of
306 !% residual into the mixing. This parameter controls this amount.
307 !%End
308 call parse_variable(namespace, trim(prefix)+'MixingResidual', 0.05_real64, smix%residual_coeff)
309 if (smix%residual_coeff <= m_zero .or. smix%residual_coeff > m_one) then
310 call messages_input_error(namespace, 'MixingResidual', 'Value should be positive and smaller than one.')
311 end if
312
313 !%Variable MixNumberSteps
314 !%Type integer
315 !%Default 4
316 !%Section SCF::Mixing
317 !%Description
318 !% In the Broyden scheme, the new input density or potential is constructed
319 !% from the values of the densities/potentials of a given number of previous iterations.
320 !% This number is set by this variable. Must be greater than 1.
321 !%End
322 if (smix%scheme /= option__mixingscheme__linear) then
323 call parse_variable(namespace, trim(prefix)//'MixNumberSteps', 4, smix%ns)
324 if (smix%ns <= 1) call messages_input_error(namespace, 'MixNumberSteps')
325 else
326 smix%ns = 0
327 end if
328
329 !%Variable MixingRestart
330 !%Type integer
331 !%Default 20
332 !%Section SCF::Mixing
333 !%Description
334 !% In the Broyden scheme, the mixing is restarted after
335 !% the number of iterations given by this variable.
336 !% Set this to zero to disable restarting the mixing.
337 !%End
338 if (smix%scheme /= option__mixingscheme__linear) then
339 call parse_variable(namespace, trim(prefix)//'MixingRestart', 20, smix%ns_restart)
340 if (smix%ns_restart < 0) call messages_input_error(namespace, 'MixingRestart')
341 else
342 smix%ns_restart = 0
343 end if
344
345 write(message(1), '(A,I4,A,I4,A)') "Info: Mixing uses ", smix%ns, " steps and restarts after ", &
346 smix%ns_restart, " steps."
347 call messages_info(1, namespace=namespace)
348
349 !%Variable MixInterval
350 !%Type integer
351 !%Default 1
352 !%Section SCF::Mixing
353 !%Description
354 !% When this variable is set to a value different than 1 (the
355 !% default) a combined mixing scheme will be used, with MixInterval
356 !% - 1 steps of linear mixing followed by 1 step of the selected
357 !% mixing. For the moment this variable only works with DIIS mixing.
358 !%End
359 call parse_variable(namespace, trim(prefix)//'MixInterval', 1, smix%interval)
360 if (smix%interval < 1) call messages_input_error(namespace, 'MixInterval', 'MixInterval must be larger or equal than 1')
361
362 smix%iter = 0
363
364 !%Variable MixingKerker
365 !%Type logical
366 !%Default false
367 !%Section SCF::Mixing
368 !%Description
369 !% (Experimental) If set to yes, Octopus will use the Kerker preconditioner
370 !% for the mixing operator applied in the linear and Broyden mixing schemes.
371 !% When using preconditioning with the Broyden scheme, the user should increase the <tt>Mixing</tt> from its default
372 !% to between 0.5 and 0.8, for optimal results.
373 !% The implementation follows Yoshinori Shiihara et al., Modelling Simul. Mater. Sci. Eng. 16 (2008), 035004
374 !%End
375 call parse_variable(namespace, trim(prefix)+'MixingKerker', .false., smix%kerker)
376 if(smix%kerker) call messages_experimental('MixingKerker')
377
378 !%Variable MixingKerkerFactor
379 !%Type float
380 !%Default 1.0
381 !%Section SCF::Mixing
382 !%Description
383 !% The screening factor, $q_0$ in the Kerker preconditioner in units of inverse length.
384 !% For small wave vectors, q, screening behaves like $\frac{Aq}{q_0}$, where A is the mixing coefficient.
385 !% An optimal value is determined empirically for a given system and mixing coeffient, however the default
386 !% is a reasonable choice in many cases.
387 !%End
388 call parse_variable(namespace, trim(prefix)+'MixingKerkerFactor', 1._real64, smix%kerker_factor)
389
390 smix%nauxmixfield = 0
391 do ii = 1,max_auxmixfield
392 nullify(smix%auxmixfield(ii)%p)
393 end do
394
395 call mixfield_init(smix, smix%mixfield, d1, d2, smix%ns, func_type)
396 call mix_clear(smix)
397 if (smix%precondition) call init_preconditioner()
398
399 pop_sub(mix_init)
400
401 contains
402
403 subroutine init_preconditioner()
404
405 integer :: ns, maxp, ip, is
406 real(real64), parameter :: weight = 50.0_real64
407
408 ! This the mixing preconditioner from GPAW:
409 !
410 ! https://wiki.fysik.dtu.dk/gpaw/documentation/densitymix/densitymix.html
411 !
412
413 assert(.not. der%mesh%use_curvilinear)
414 assert(der%dim == 3)
415
416 call nl_operator_init(smix%preconditioner, "Mixing preconditioner")
417 call stencil_cube_get_lapl(smix%preconditioner%stencil, der%dim, 1)
418 call nl_operator_build(space, der%mesh, smix%preconditioner, der%mesh%np, const_w = .not. der%mesh%use_curvilinear)
419
420 ns = smix%preconditioner%stencil%size
421
422 if (smix%preconditioner%const_w) then
423 maxp = 1
424 else
425 maxp = der%mesh%np
426 end if
427
428 do ip = 1, maxp
429
430 do is = 1, ns
431 select case (sum(abs(smix%preconditioner%stencil%points(1:der%dim, is))))
432 case (0)
433 smix%preconditioner%w(is, ip) = m_one + weight/8.0_real64
434 case (1)
435 smix%preconditioner%w(is, ip) = weight/16.0_real64
436 case (2)
437 smix%preconditioner%w(is, ip) = weight/32.0_real64
438 case (3)
439 smix%preconditioner%w(is, ip) = weight/64.0_real64
440 case default
441 assert(.false.)
442 end select
443
444 end do
445 end do
446
447 call nl_operator_allocate_gpu_buffers(smix%preconditioner)
448 call nl_operator_update_gpu_buffers(smix%preconditioner)
449 call nl_operator_output_weights(smix%preconditioner)
450
451 end subroutine init_preconditioner
452
453 end subroutine mix_init
454
455
456 ! ---------------------------------------------------------
457 subroutine mix_clear(smix)
458 type(mix_t), intent(inout) :: smix
459
460 push_sub(mix_clear)
461
462 call mixfield_clear(smix%scheme, smix%mixfield)
463
464 smix%iter = 0
465 smix%last_ipos = 0
466 smix%adaptive_damping = m_one
467
468 pop_sub(mix_clear)
469 end subroutine mix_clear
470
471
472 ! ---------------------------------------------------------
473 subroutine mix_end(smix)
474 type(mix_t), intent(inout) :: smix
475
476 integer :: ii
477
478 push_sub(mix_end)
479
480 if (smix%precondition) call nl_operator_end(smix%preconditioner)
481
482 call mixfield_end(smix, smix%mixfield)
483
484 smix%nauxmixfield = 0
485 do ii = 1,max_auxmixfield
486 nullify(smix%auxmixfield(ii)%p)
487 end do
488
489 safe_deallocate_a(smix%coeff)
490
491 pop_sub(mix_end)
492 end subroutine mix_end
493
494
495 ! ---------------------------------------------------------
496 subroutine mix_dump(namespace, restart, smix, mesh, ierr)
497 type(namespace_t), intent(in) :: namespace
498 class(restart_t), intent(in) :: restart
499 type(mix_t), intent(in) :: smix
500 class(mesh_t), intent(in) :: mesh
501 integer, intent(out) :: ierr
502
503 integer :: iunit, id2, id3, err, err2(4)
504 character(len=40) :: lines(7)
505 character(len=80) :: filename
506
507 push_sub(mix_dump)
508
509 ierr = 0
510
511 if (restart%skip()) then
512 pop_sub(mix_dump)
513 return
514 end if
515
516 message(1) = "Debug: Writing mixing restart."
517 call messages_info(1, namespace=namespace, debug_only=.true.)
518
519 ! functions to be written need to be compatible with the mesh
520 assert(mesh%np == smix%mixfield%d1)
521
522 ! First we write some information about the mixing
523 iunit = restart%open('mixing')
524 if(restart%do_i_write()) then
525 write(lines(1), '(a11,i1)') 'scheme= ', smix%scheme
526 ! Number of global mesh points have to be written, not only smix%d1
527 write(lines(2), '(a11,i10)') 'd1= ', mesh%np_global
528 write(lines(3), '(a11,i10)') 'd2= ', smix%mixfield%d2
529 write(lines(4), '(a11,i10)') 'd3= ', smix%mixfield%d3
530 write(lines(5), '(a11,i10)') 'iter= ', smix%iter
531 write(lines(6), '(a11,i10)') 'ns= ', smix%ns
532 write(lines(7), '(a11,i10)') 'last_ipos= ', smix%last_ipos
533 call restart%write(iunit, lines, 7, err)
534 ierr = ierr + err
535 end if
536 call restart%close(iunit)
537
538 ! Now we write the different functions.
539 ! These are not needed when using linear mixing, so we will make sure we skip this step in that case.
540 err2 = 0
541 if (smix%scheme /= option__mixingscheme__linear) then
542 do id2 = 1, smix%mixfield%d2
543 do id3 = 1, smix%mixfield%d3
544
545 write(filename,'(a3,i2.2,i2.2)') 'df_', id2, id3
546 if (smix%mixfield%func_type == type_float) then
547 call restart%write_mesh_function(filename, mesh, smix%mixfield%ddf(1:mesh%np, id2, id3), err)
548 else
549 call restart%write_mesh_function(filename, mesh, smix%mixfield%zdf(1:mesh%np, id2, id3), err)
550 end if
551 err2(1) = err2(1) + err
553 write(filename,'(a3,i2.2,i2.2)') 'dv_', id2, id3
554 if (smix%mixfield%func_type == type_float) then
555 call restart%write_mesh_function(filename, mesh, smix%mixfield%ddv(1:mesh%np, id2, id3), err)
556 else
557 call restart%write_mesh_function(filename, mesh, smix%mixfield%zdv(1:mesh%np, id2, id3), err)
558 end if
559 err2(2) = err2(2) + err
560
561 end do
562
563 write(filename,'(a6,i2.2)') 'f_old_', id2
564 if (smix%mixfield%func_type == type_float) then
565 call restart%write_mesh_function(filename, mesh, smix%mixfield%df_old(1:mesh%np, id2), err)
566 else
567 call restart%write_mesh_function(filename, mesh, smix%mixfield%zf_old(1:mesh%np, id2), err)
568 end if
569 err2(3) = err2(3) + err
570
571 write(filename,'(a8,i2.2)') 'vin_old_', id2
572 if (smix%mixfield%func_type == type_float) then
573 call restart%write_mesh_function(filename, mesh, smix%mixfield%dvin_old(1:mesh%np, id2), err)
574 else
575 call restart%write_mesh_function(filename, mesh, smix%mixfield%zvin_old(1:mesh%np, id2), err)
576 end if
577 err2(4) = err2(4) + err
578
579 end do
580
581 if (err2(1) /= 0) ierr = ierr + 2
582 if (err2(2) /= 0) ierr = ierr + 4
583 if (err2(3) /= 0) ierr = ierr + 8
584 if (err2(4) /= 0) ierr = ierr + 16
585 end if
586
587 message(1) = "Debug: Writing mixing restart done."
588 call messages_info(1, namespace=namespace, debug_only=.true.)
589
590 pop_sub(mix_dump)
591 end subroutine mix_dump
592
593
594 !---------------------------------------------------------
595 subroutine mix_load(namespace, restart, smix, mesh, ierr)
596 type(namespace_t), intent(in) :: namespace
597 type(restart_t), intent(in) :: restart
598 type(mix_t), intent(inout) :: smix
599 class(mesh_t), intent(in) :: mesh
600 integer, intent(out) :: ierr
601
602 integer :: iunit, err, err2(4)
603 integer :: scheme, d1, d2, d3, ns
604 integer :: id2, id3
605 character(len=11) :: str
606 character(len=80) :: filename
607 character(len=256) :: lines(7)
608
609 push_sub(mix_load)
610
611 ierr = 0
612
613 if (restart%skip()) then
614 ierr = -1
615 pop_sub(mix_load)
616 return
617 end if
618
619 message(1) = "Debug: Reading mixing restart."
620 call messages_info(1, namespace=namespace, debug_only=.true.)
621
622 ! First we read some information about the mixing
623 iunit = restart%open('mixing')
624 call restart%read(iunit, lines, 7, err)
625 if (err /= 0) then
626 ierr = ierr + 1
627 else
628 read(lines(1), *) str, scheme
629 read(lines(2), *) str, d1
630 read(lines(3), *) str, d2
631 read(lines(4), *) str, d3
632 read(lines(5), *) str, smix%iter
633 read(lines(6), *) str, ns
634 read(lines(7), *) str, smix%last_ipos
635 end if
636 call restart%close(iunit)
637
638
639 if (ierr == 0) then
640 ! We can only use the restart information if the mixing scheme and the number of steps used remained the same
641 if (scheme /= smix%scheme .or. ns /= smix%ns) then
642 message(1) = "The mixing scheme from the restart data is not the same as the one used in the current calculation."
643 call messages_warning(1, namespace=namespace)
644 ierr = ierr + 2
645 end if
646
647 ! Check the dimensions of the arrays to be read
648 if (mesh%np_global /= d1 .or. mesh%np /= smix%mixfield%d1 .or. d2 /= smix%mixfield%d2) then
649 message(1) = "The dimensions of the arrays from the mixing restart data"
650 message(2) = "are not the same as the ones used in this calculation."
651 call messages_warning(2, namespace=namespace)
652 ierr = ierr + 4
653 end if
654 end if
655
656
657 ! Now we read the different functions.
658 ! Note that we may have more or less functions than the ones needed (d3 /= smix%d3)
659 if (ierr == 0) then
660 if (smix%scheme /= option__mixingscheme__linear) then
661 err2 = 0
662 do id2 = 1, smix%mixfield%d2
663 do id3 = 1, smix%mixfield%d3
664
665 write(filename,'(a3,i2.2,i2.2)') 'df_', id2, id3
666 if (smix%mixfield%func_type == type_float) then
667 call restart%read_mesh_function(filename, mesh, smix%mixfield%ddf(1:mesh%np, id2, id3), err)
668 else
669 call restart%read_mesh_function(filename, mesh, smix%mixfield%zdf(1:mesh%np, id2, id3), err)
670 end if
671 if (err /= 0) err2(1) = err2(1) + 1
672
673 write(filename,'(a3,i2.2,i2.2)') 'dv_', id2, id3
674 if (smix%mixfield%func_type == type_float) then
675 call restart%read_mesh_function(filename, mesh, smix%mixfield%ddv(1:mesh%np, id2, id3), err)
676 else
677 call restart%read_mesh_function(filename, mesh, smix%mixfield%zdv(1:mesh%np, id2, id3), err)
678 end if
679 if (err /= 0) err2(2) = err2(2) + 1
680
681 end do
682
683 write(filename,'(a6,i2.2)') 'f_old_', id2
684 if (smix%mixfield%func_type == type_float) then
685 call restart%read_mesh_function(filename, mesh, smix%mixfield%df_old(1:mesh%np, id2), err)
686 else
687 call restart%read_mesh_function(filename, mesh, smix%mixfield%zf_old(1:mesh%np, id2), err)
688 end if
689 if (err /= 0) err2(3) = err2(3) + 1
691 write(filename,'(a8,i2.2)') 'vin_old_', id2
692 if (smix%mixfield%func_type == type_float) then
693 call restart%read_mesh_function(filename, mesh, smix%mixfield%dvin_old(1:mesh%np, id2), err)
694 else
695 call restart%read_mesh_function(filename, mesh, smix%mixfield%zvin_old(1:mesh%np, id2), err)
696 end if
697 if (err /= 0) err2(4) = err2(4) + 1
698
699 end do
700
701 if (err2(1) /= 0) ierr = ierr + 8
702 if (err2(2) /= 0) ierr = ierr + 16
703 if (err2(3) /= 0) ierr = ierr + 32
704 if (err2(4) /= 0) ierr = ierr + 64
705 end if
706 end if
707
708 if (ierr /= 0) then
709 ! Something went wront, so make sure we start from scratch
710 call mix_clear(smix)
711 end if
712
713 if (ierr == 0 .and. smix%scheme == option__mixingscheme__broyden_adaptive) then
714 smix%adaptive_damping = m_one
715 smix%mixfield%last_residue = m_zero
716 end if
717
718
719 message(1) = "Debug: Reading mixing restart done."
720 call messages_info(1, namespace=namespace, debug_only=.true.)
721
722 pop_sub(mix_load)
723 end subroutine mix_load
724
725 real(real64) pure function mix_coefficient(this) result(coefficient)
726 type(mix_t), intent(in) :: this
727
728 coefficient = this%coeff(1)
729 end function mix_coefficient
730
731 integer pure function mix_scheme(this) result(scheme)
732 type(mix_t), intent(in) :: this
733
734 scheme = this%scheme
735 end function mix_scheme
736
737 integer pure function mix_d3(this)
738 type(mix_t), intent(in) :: this
739
740 mix_d3 = this%mixfield%d3
741 end function mix_d3
742
743 subroutine mix_get_field(this, mixfield)
744 type(mix_t), target, intent(in) :: this
745 type(mixfield_t), pointer, intent(out) :: mixfield
746
747 mixfield => this%mixfield
748 end subroutine mix_get_field
749
751 subroutine mixing(namespace, smix)
752 type(namespace_t), intent(in) :: namespace
753 type(mix_t), intent(inout) :: smix
754
755 push_sub(mixing)
756
757 if (smix%mixfield%func_type == type_float) then
758 call dmixing(namespace, smix, smix%mixfield%dvin, smix%mixfield%dvout, smix%mixfield%dvnew)
759 else
760 call zmixing(namespace, smix, smix%mixfield%zvin, smix%mixfield%zvout, smix%mixfield%zvnew)
761 end if
762
763 pop_sub(mixing)
764 end subroutine mixing
765
766 subroutine mix_add_auxmixfield(namespace, smix, mixfield)
767 type(namespace_t), intent(in) :: namespace
768 type(mix_t), intent(inout) :: smix
769 type(mixfield_t), target, intent(in) :: mixfield
770
771 push_sub(mix_add_auxmixfield)
772
773 smix%nauxmixfield = smix%nauxmixfield + 1
774 smix%auxmixfield(smix%nauxmixfield)%p => mixfield
775
776 if (smix%scheme == option__mixingscheme__diis) then
777 message(1) = 'Mixing scheme DIIS is not implemented for auxiliary mixing fields'
778 call messages_fatal(1, namespace=namespace)
779 end if
780
781 pop_sub(mix_add_auxmixfield)
782 end subroutine mix_add_auxmixfield
783
785 subroutine mixfield_init(smix, mixfield, d1, d2, d3, func_type)
786 type(mix_t), intent(in) :: smix
787 type(mixfield_t), intent(inout) :: mixfield
788 integer, intent(in) :: d1, d2, d3
789 type(type_t), intent(in) :: func_type
790
791 push_sub(mixfield_init)
792
793 mixfield%d1 = d1
794 mixfield%d2 = d2
795 mixfield%d3 = d3
796
797 mixfield%func_type = func_type
798
799 mixfield%mix_spin_density_matrix = smix%mix_spin_density_matrix
800 mixfield%last_residue = m_zero
801
802 if (smix%scheme /= option__mixingscheme__linear) then
803 if (mixfield%func_type == type_float) then
804 safe_allocate( mixfield%ddf(1:d1, 1:d2, 1:d3))
805 safe_allocate( mixfield%ddv(1:d1, 1:d2, 1:d3))
806 safe_allocate(mixfield%dvin_old(1:d1, 1:d2))
807 safe_allocate( mixfield%df_old(1:d1, 1:d2))
808 else
809 safe_allocate( mixfield%zdf(1:d1, 1:d2, 1:d3))
810 safe_allocate( mixfield%zdv(1:d1, 1:d2, 1:d3))
811 safe_allocate(mixfield%zvin_old(1:d1, 1:d2))
812 safe_allocate( mixfield%zf_old(1:d1, 1:d2))
813 end if
814 end if
815
816 if (mixfield%func_type == type_float) then
817 safe_allocate(mixfield%dvin(1:d1, 1:d2))
818 safe_allocate(mixfield%dvout(1:d1, 1:d2))
819 safe_allocate(mixfield%dvnew(1:d1, 1:d2))
820 safe_allocate(mixfield%dresidual(1:d1, 1:d2))
821 else
822 safe_allocate(mixfield%zvin(1:d1, 1:d2))
823 safe_allocate(mixfield%zvout(1:d1, 1:d2))
824 safe_allocate(mixfield%zvnew(1:d1, 1:d2))
825 safe_allocate(mixfield%zresidual(1:d1, 1:d2))
826 end if
827
828 pop_sub(mixfield_init)
829 end subroutine mixfield_init
830
832 subroutine mixfield_end(smix, mixfield)
833 type(mix_t), intent(inout) :: smix
834 type(mixfield_t), intent(inout) :: mixfield
835
836 push_sub(mixfield_end)
837
838 ! Arrays got allocated for all mixing schemes, except linear mixing
839 if (smix%scheme /= option__mixingscheme__linear) then
840 if (mixfield%func_type == type_float) then
841 safe_deallocate_a(mixfield%ddf)
842 safe_deallocate_a(mixfield%ddv)
843 safe_deallocate_a(mixfield%dvin_old)
844 safe_deallocate_a(mixfield%df_old)
845 else
846 safe_deallocate_a(mixfield%zdf)
847 safe_deallocate_a(mixfield%zdv)
848 safe_deallocate_a(mixfield%zvin_old)
849 safe_deallocate_a(mixfield%zf_old)
850 end if
851 end if
852
853 if (mixfield%func_type == type_float) then
854 safe_deallocate_a(mixfield%dvin)
855 safe_deallocate_a(mixfield%dvout)
856 safe_deallocate_a(mixfield%dvnew)
857 safe_deallocate_a(mixfield%dresidual)
858 else
859 safe_deallocate_a(mixfield%zvin)
860 safe_deallocate_a(mixfield%zvout)
861 safe_deallocate_a(mixfield%zvnew)
862 safe_deallocate_a(mixfield%zresidual)
863 end if
864
865 pop_sub(mixfield_end)
866 end subroutine mixfield_end
867
869 subroutine mixfield_clear(scheme, mixfield)
870 integer, intent(in) :: scheme
871 type(mixfield_t), intent(inout) :: mixfield
872 integer :: d1, d2, d3
873
874 push_sub(mixfield_clear)
875
876 d1 = mixfield%d1
877 d2 = mixfield%d2
878 d3 = mixfield%d3
879
880 if (scheme /= option__mixingscheme__linear) then
881 if (mixfield%func_type == type_float) then
882 assert(allocated(mixfield%ddf))
883 mixfield%ddf(1:d1, 1:d2, 1:d3) = m_zero
884 mixfield%ddv(1:d1, 1:d2, 1:d3) = m_zero
885 mixfield%dvin_old(1:d1, 1:d2) = m_zero
886 mixfield%df_old(1:d1, 1:d2) = m_zero
887 else
888 assert(allocated(mixfield%zdf))
889 mixfield%zdf(1:d1, 1:d2, 1:d3) = m_z0
890 mixfield%zdv(1:d1, 1:d2, 1:d3) = m_z0
891 mixfield%zvin_old(1:d1, 1:d2) = m_z0
892 mixfield%zf_old(1:d1, 1:d2) = m_z0
893 end if
894 end if
895
896 if (mixfield%func_type == type_float) then
897 mixfield%dvin(1:d1, 1:d2) = m_zero
898 mixfield%dvout(1:d1, 1:d2) = m_zero
899 mixfield%dvnew(1:d1, 1:d2) = m_zero
900 mixfield%dresidual(1:d1, 1:d2) = m_zero
901 else
902 mixfield%zvin(1:d1, 1:d2) = m_z0
903 mixfield%zvout(1:d1, 1:d2) = m_z0
904 mixfield%zvnew(1:d1, 1:d2) = m_z0
905 mixfield%zresidual(1:d1, 1:d2) = m_z0
906 end if
907
908 pop_sub(mixfield_clear)
909 end subroutine mixfield_clear
910
919 subroutine compute_residuals_aux_field(this)
920 class(mix_t), intent(inout) :: this
921
922 type(mixfield_t), pointer :: aux_field
923 integer :: dims(2)
924 integer :: i
925
928 do i = 1, this%nauxmixfield
929 aux_field => this%auxmixfield(i)%p
930 dims = [aux_field%d1, aux_field%d2]
931
932 if (aux_field%func_type == type_float) then
933 call dcompute_residual(this, dims, aux_field%dvin, aux_field%dvout, aux_field%dresidual)
934 else
935 call zcompute_residual(this, dims, aux_field%zvin, aux_field%zvout, aux_field%zresidual)
936 end if
937
938 end do
939 nullify(aux_field)
940
942
943 end subroutine compute_residuals_aux_field
944
956 subroutine linear_mixing_aux_field(smix)
957 type(mix_t), intent(inout) :: smix
958
959 type(mixfield_t), pointer :: field
960 integer :: i
961
963
964 call smix%compute_residuals_aux_field()
965
966 do i = 1, smix%nauxmixfield
967 field => smix%auxmixfield(i)%p
968
969 if (field%func_type == type_float) then
970 call dmixing_linear(field%d1, field%d2, smix%coeff, field%dvin, field%dresidual, field%dvnew)
971 else
972 call zmixing_linear(field%d1, field%d2, smix%coeff, field%zvin, field%zresidual, field%zvnew)
973 end if
974
975 end do
976
977 nullify(field)
978
980 end subroutine linear_mixing_aux_field
981
982#include "undef.F90"
983#include "real.F90"
984
985#include "mix_inc.F90"
986
987#include "undef.F90"
988#include "complex.F90"
989
990#include "mix_inc.F90"
991
992end module mix_oct_m
993
994!! Local Variables:
995!! mode: f90
996!! coding: utf-8
997!! End:
subroutine init_preconditioner()
Definition: mix.F90:499
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
real(real64), parameter, public m_zero
Definition: global.F90:191
real(real64), parameter, public m_one
Definition: global.F90:192
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:117
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:524
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:999
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:161
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:690
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1039
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:593
subroutine, public mixfield_end(smix, mixfield)
Deallocate all arrays of a mixfield instance.
Definition: mix.F90:928
integer pure function, public mix_scheme(this)
Definition: mix.F90:827
subroutine dmixfield_set_vin(mixfield, vin)
Definition: mix.F90:1941
real(real64) pure function, public mix_coefficient(this)
Definition: mix.F90:821
subroutine, public mixing(namespace, smix)
Main entry-point to SCF mixer.
Definition: mix.F90:847
subroutine, public mixfield_init(smix, mixfield, d1, d2, d3, func_type)
Initialise all attributes of a mixfield instance.
Definition: mix.F90:881
subroutine zmixfield_get_vnew(mixfield, vnew)
Definition: mix.F90:2947
subroutine dcompute_residual(smix, dims, vin, vout, residual)
Compute the residual of the potential (or density) for SCF mixing.
Definition: mix.F90:1240
subroutine, public mix_get_field(this, mixfield)
Definition: mix.F90:839
subroutine compute_residuals_aux_field(this)
Compute the residuals of the auxilliary fields.
Definition: mix.F90:1015
subroutine zmixfield_set_vout(mixfield, vout)
Definition: mix.F90:2928
subroutine, public mix_add_auxmixfield(namespace, smix, mixfield)
Definition: mix.F90:862
subroutine dmixfield_get_vnew(mixfield, vnew)
Definition: mix.F90:2007
subroutine dmixing_linear(d1, d2, coeff, vin, residual, vnew)
Linear mixing of the input and output potentials.
Definition: mix.F90:1208
subroutine, public mix_dump(namespace, restart, smix, mesh, ierr)
Definition: mix.F90:592
subroutine, public mixfield_clear(scheme, mixfield)
Zero all potential and field attributes of a mixfield instance.
Definition: mix.F90:965
subroutine, public mix_init(smix, namespace, space, der, d1, d2, def_, func_type_, prefix_)
Initialise mix_t instance.
Definition: mix.F90:270
subroutine zmixfield_set_vin(mixfield, vin)
Definition: mix.F90:2909
subroutine, public zmixing(namespace, smix, vin, vout, vnew)
Mix the input and output potentials (or densities) according to some scheme.
Definition: mix.F90:2119
integer pure function, public mix_d3(this)
Definition: mix.F90:833
subroutine, public mix_load(namespace, restart, smix, mesh, ierr)
Definition: mix.F90:691
subroutine dmixfield_set_vout(mixfield, vout)
Definition: mix.F90:1974
subroutine zcompute_residual(smix, dims, vin, vout, residual)
Compute the residual of the potential (or density) for SCF mixing.
Definition: mix.F90:2208
subroutine, public dmixing(namespace, smix, vin, vout, vnew)
Mix the input and output potentials (or densities) according to some scheme.
Definition: mix.F90:1151
subroutine linear_mixing_aux_field(smix)
Linear mixing of the auxilliary fields.
Definition: mix.F90:1052
subroutine zmixing_linear(d1, d2, coeff, vin, residual, vnew)
Linear mixing of the input and output potentials.
Definition: mix.F90:2176
subroutine, public mix_end(smix)
Definition: mix.F90:569
subroutine, public mix_clear(smix)
Definition: mix.F90:553
This module defines non-local operators.
subroutine, public nl_operator_init(op, label, symm)
initialize an instance of a non-local operator by setting the label
subroutine, public nl_operator_update_gpu_buffers(op)
subroutine, public nl_operator_output_weights(this)
subroutine, public nl_operator_end(op)
subroutine, public nl_operator_build(space, mesh, op, np, const_w, regenerate)
Creates the nonlocal operators for the stencils used for finite differences.
subroutine, public nl_operator_allocate_gpu_buffers(op)
This module is intended to contain "only mathematical" functions and procedures.
Definition: solvers.F90:117
This module defines routines, generating operators for a cubic stencil.
subroutine, public stencil_cube_get_lapl(this, dim, order)
type(type_t), public type_float
Definition: types.F90:135
class representing derivatives
Describes mesh distribution to nodes.
Definition: mesh.F90:187
God class for mixing.
Definition: mix.F90:217
Quantities used in mixing: Input, output and new potentials, and the residuals.
Definition: mix.F90:173
int true(void)