27 use,
intrinsic :: iso_fortran_env
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(:, :)
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(:, :)
100 type(type_t) :: func_type
101 integer :: d1, d2, d3
104 logical :: mix_spin_density_matrix
105 real(real64) :: last_residue
110 type(mixfield_t),
pointer :: p
113 integer,
parameter :: MAX_AUXMIXFIELD = 5
124 type(derivatives_t),
pointer :: der
127 real(real64),
allocatable :: coeff(:)
128 real(real64) :: adaptive_damping
129 real(real64) :: residual_coeff
132 integer,
public :: ns
133 integer,
public :: ns_restart
140 type(mixfield_t) :: mixfield
141 integer :: nauxmixfield
142 type(mixfield_ptr_t) :: auxmixfield(MAX_AUXMIXFIELD)
148 real(real64) :: kerker_factor
149 logical :: precondition
150 type(nl_operator_t) :: preconditioner
153 logical :: mix_spin_density_matrix
155 procedure :: compute_residuals_aux_field
174 subroutine mix_init(smix, namespace, space, der, d1, d2, def_, func_type_, prefix_)
175 type(mix_t),
intent(out) :: smix
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_
186 character(len=32) :: prefix
188 real(real64) :: coeff
194 def = option__mixingscheme__broyden
195 if (
present(def_)) def = def_
196 if (
present(func_type_))
then
197 func_type = func_type_
202 if (
present(prefix_)) prefix = prefix_
251 call parse_variable(namespace, trim(prefix)+
'MixingPreconditioner', .false., smix%precondition)
252 if (der%dim /= 3) smix%precondition = .false.
264 call parse_variable(namespace, trim(prefix)+
'Mixing', 0.3_real64, coeff)
268 safe_allocate(smix%coeff(1:d2))
270 smix%adaptive_damping =
m_one
283 call parse_variable(namespace, trim(prefix)+
'MixingSpinDensityMatrix', .false., smix%mix_spin_density_matrix)
293 call parse_variable(namespace, trim(prefix)+
'MixingMagnetization', 1.5_real64, coeff)
295 if (.not. smix%mix_spin_density_matrix) smix%coeff(2:d2) = coeff
297 smix%mix_spin_density_matrix = .
true.
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.')
322 if (smix%scheme /= option__mixingscheme__linear)
then
323 call parse_variable(namespace, trim(prefix)//
'MixNumberSteps', 4, smix%ns)
338 if (smix%scheme /= option__mixingscheme__linear)
then
339 call parse_variable(namespace, trim(prefix)//
'MixingRestart', 20, smix%ns_restart)
345 write(
message(1),
'(A,I4,A,I4,A)')
"Info: Mixing uses ", smix%ns,
" steps and restarts after ", &
346 smix%ns_restart,
" steps."
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')
375 call parse_variable(namespace, trim(prefix)+
'MixingKerker', .false., smix%kerker)
388 call parse_variable(namespace, trim(prefix)+
'MixingKerkerFactor', 1._real64, smix%kerker_factor)
390 smix%nauxmixfield = 0
391 do ii = 1,max_auxmixfield
392 nullify(smix%auxmixfield(ii)%p)
395 call mixfield_init(smix, smix%mixfield, d1, d2, smix%ns, func_type)
405 integer :: ns, maxp, ip, is
406 real(real64),
parameter :: weight = 50.0_real64
413 assert(.not. der%mesh%use_curvilinear)
418 call nl_operator_build(space, der%mesh, smix%preconditioner, der%mesh%np, const_w = .not. der%mesh%use_curvilinear)
420 ns = smix%preconditioner%stencil%size
422 if (smix%preconditioner%const_w)
then
431 select case (sum(abs(smix%preconditioner%stencil%points(1:der%dim, is))))
433 smix%preconditioner%w(is, ip) =
m_one + weight/8.0_real64
435 smix%preconditioner%w(is, ip) = weight/16.0_real64
437 smix%preconditioner%w(is, ip) = weight/32.0_real64
439 smix%preconditioner%w(is, ip) = weight/64.0_real64
458 type(
mix_t),
intent(inout) :: smix
466 smix%adaptive_damping =
m_one
474 type(
mix_t),
intent(inout) :: smix
484 smix%nauxmixfield = 0
485 do ii = 1,max_auxmixfield
486 nullify(smix%auxmixfield(ii)%p)
489 safe_deallocate_a(smix%coeff)
496 subroutine mix_dump(namespace, restart, smix, mesh, ierr)
499 type(
mix_t),
intent(in) :: smix
500 class(
mesh_t),
intent(in) :: mesh
501 integer,
intent(out) :: ierr
503 integer :: iunit, id2, id3, err, err2(4)
504 character(len=40) :: lines(7)
505 character(len=80) :: filename
511 if (restart%skip())
then
516 message(1) =
"Debug: Writing mixing restart."
520 assert(mesh%np == smix%mixfield%d1)
523 iunit = restart%open(
'mixing')
524 if(restart%do_i_write())
then
525 write(lines(1),
'(a11,i1)')
'scheme= ', smix%scheme
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)
536 call restart%close(iunit)
541 if (smix%scheme /= option__mixingscheme__linear)
then
542 do id2 = 1, smix%mixfield%d2
543 do id3 = 1, smix%mixfield%d3
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)
549 call restart%write_mesh_function(filename, mesh, smix%mixfield%zdf(1:mesh%np, id2, id3), err)
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)
557 call restart%write_mesh_function(filename, mesh, smix%mixfield%zdv(1:mesh%np, id2, id3), err)
559 err2(2) = err2(2) + err
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)
567 call restart%write_mesh_function(filename, mesh, smix%mixfield%zf_old(1:mesh%np, id2), err)
569 err2(3) = err2(3) + err
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)
575 call restart%write_mesh_function(filename, mesh, smix%mixfield%zvin_old(1:mesh%np, id2), err)
577 err2(4) = err2(4) + err
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
587 message(1) =
"Debug: Writing mixing restart done."
595 subroutine mix_load(namespace, restart, smix, mesh, ierr)
598 type(
mix_t),
intent(inout) :: smix
599 class(
mesh_t),
intent(in) :: mesh
600 integer,
intent(out) :: ierr
602 integer :: iunit, err, err2(4)
603 integer :: scheme, d1, d2, d3, ns
605 character(len=11) :: str
606 character(len=80) :: filename
607 character(len=256) :: lines(7)
613 if (restart%skip())
then
619 message(1) =
"Debug: Reading mixing restart."
623 iunit = restart%open(
'mixing')
624 call restart%read(iunit, lines, 7, err)
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
636 call restart%close(iunit)
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."
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."
660 if (smix%scheme /= option__mixingscheme__linear)
then
662 do id2 = 1, smix%mixfield%d2
663 do id3 = 1, smix%mixfield%d3
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)
669 call restart%read_mesh_function(filename, mesh, smix%mixfield%zdf(1:mesh%np, id2, id3), err)
671 if (err /= 0) err2(1) = err2(1) + 1
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)
677 call restart%read_mesh_function(filename, mesh, smix%mixfield%zdv(1:mesh%np, id2, id3), err)
679 if (err /= 0) err2(2) = err2(2) + 1
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)
687 call restart%read_mesh_function(filename, mesh, smix%mixfield%zf_old(1:mesh%np, id2), err)
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)
695 call restart%read_mesh_function(filename, mesh, smix%mixfield%zvin_old(1:mesh%np, id2), err)
697 if (err /= 0) err2(4) = err2(4) + 1
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
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
719 message(1) =
"Debug: Reading mixing restart done."
726 type(
mix_t),
intent(in) :: this
728 coefficient = this%coeff(1)
731 integer pure function
mix_scheme(this) result(scheme)
732 type(
mix_t),
intent(in) :: this
737 integer pure function
mix_d3(this)
738 type(
mix_t),
intent(in) :: this
744 type(
mix_t),
target,
intent(in) :: this
745 type(
mixfield_t),
pointer,
intent(out) :: mixfield
747 mixfield => this%mixfield
751 subroutine mixing(namespace, smix)
752 type(namespace_t),
intent(in) :: namespace
753 type(
mix_t),
intent(inout) :: smix
757 if (smix%mixfield%func_type == type_float)
then
758 call dmixing(namespace, smix, smix%mixfield%dvin, smix%mixfield%dvout, smix%mixfield%dvnew)
760 call zmixing(namespace, smix, smix%mixfield%zvin, smix%mixfield%zvout, smix%mixfield%zvnew)
767 type(namespace_t),
intent(in) :: namespace
768 type(
mix_t),
intent(inout) :: smix
769 type(
mixfield_t),
target,
intent(in) :: mixfield
773 smix%nauxmixfield = smix%nauxmixfield + 1
774 smix%auxmixfield(smix%nauxmixfield)%p => mixfield
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)
785 subroutine mixfield_init(smix, mixfield, d1, d2, d3, func_type)
786 type(
mix_t),
intent(in) :: smix
788 integer,
intent(in) :: d1, d2, d3
789 type(type_t),
intent(in) :: func_type
797 mixfield%func_type = func_type
799 mixfield%mix_spin_density_matrix = smix%mix_spin_density_matrix
800 mixfield%last_residue = m_zero
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))
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))
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))
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))
833 type(
mix_t),
intent(inout) :: smix
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)
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)
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)
859 safe_deallocate_a(mixfield%zvin)
860 safe_deallocate_a(mixfield%zvout)
861 safe_deallocate_a(mixfield%zvnew)
862 safe_deallocate_a(mixfield%zresidual)
870 integer,
intent(in) :: scheme
872 integer :: d1, d2, d3
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
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
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
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
920 class(
mix_t),
intent(inout) :: this
928 do i = 1, this%nauxmixfield
929 aux_field => this%auxmixfield(i)%p
930 dims = [aux_field%d1, aux_field%d2]
932 if (aux_field%func_type == type_float)
then
933 call dcompute_residual(this, dims, aux_field%dvin, aux_field%dvout, aux_field%dresidual)
935 call zcompute_residual(this, dims, aux_field%zvin, aux_field%zvout, aux_field%zresidual)
957 type(
mix_t),
intent(inout) :: smix
964 call smix%compute_residuals_aux_field()
966 do i = 1, smix%nauxmixfield
967 field => smix%auxmixfield(i)%p
969 if (field%func_type == type_float)
then
970 call dmixing_linear(field%d1, field%d2, smix%coeff, field%dvin, field%dresidual, field%dvnew)
972 call zmixing_linear(field%d1, field%d2, smix%coeff, field%zvin, field%zresidual, field%zvnew)
985#include "mix_inc.F90"
988#include "complex.F90"
990#include "mix_inc.F90"
subroutine init_preconditioner()
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
real(real64), parameter, public m_zero
real(real64), parameter, public m_one
This module is intended to contain "only mathematical" functions and procedures.
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
subroutine, public messages_warning(no_lines, all_nodes, namespace)
subroutine, public messages_obsolete_variable(namespace, name, rep)
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_input_error(namespace, var, details, row, column)
subroutine, public messages_experimental(name, namespace)
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
subroutine, public mixfield_end(smix, mixfield)
Deallocate all arrays of a mixfield instance.
integer pure function, public mix_scheme(this)
subroutine dmixfield_set_vin(mixfield, vin)
real(real64) pure function, public mix_coefficient(this)
subroutine, public mixing(namespace, smix)
Main entry-point to SCF mixer.
subroutine, public mixfield_init(smix, mixfield, d1, d2, d3, func_type)
Initialise all attributes of a mixfield instance.
subroutine zmixfield_get_vnew(mixfield, vnew)
subroutine dcompute_residual(smix, dims, vin, vout, residual)
Compute the residual of the potential (or density) for SCF mixing.
subroutine, public mix_get_field(this, mixfield)
subroutine compute_residuals_aux_field(this)
Compute the residuals of the auxilliary fields.
subroutine zmixfield_set_vout(mixfield, vout)
subroutine, public mix_add_auxmixfield(namespace, smix, mixfield)
subroutine dmixfield_get_vnew(mixfield, vnew)
subroutine dmixing_linear(d1, d2, coeff, vin, residual, vnew)
Linear mixing of the input and output potentials.
subroutine, public mix_dump(namespace, restart, smix, mesh, ierr)
subroutine, public mixfield_clear(scheme, mixfield)
Zero all potential and field attributes of a mixfield instance.
subroutine, public mix_init(smix, namespace, space, der, d1, d2, def_, func_type_, prefix_)
Initialise mix_t instance.
subroutine zmixfield_set_vin(mixfield, vin)
subroutine, public zmixing(namespace, smix, vin, vout, vnew)
Mix the input and output potentials (or densities) according to some scheme.
integer pure function, public mix_d3(this)
subroutine, public mix_load(namespace, restart, smix, mesh, ierr)
subroutine dmixfield_set_vout(mixfield, vout)
subroutine zcompute_residual(smix, dims, vin, vout, residual)
Compute the residual of the potential (or density) for SCF mixing.
subroutine, public dmixing(namespace, smix, vin, vout, vnew)
Mix the input and output potentials (or densities) according to some scheme.
subroutine linear_mixing_aux_field(smix)
Linear mixing of the auxilliary fields.
subroutine zmixing_linear(d1, d2, coeff, vin, residual, vnew)
Linear mixing of the input and output potentials.
subroutine, public mix_end(smix)
subroutine, public mix_clear(smix)
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.
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
class representing derivatives
Describes mesh distribution to nodes.
Quantities used in mixing: Input, output and new potentials, and the residuals.