Octopus
ks_potential.F90
Go to the documentation of this file.
1!! Copyright (C) 2024 N. Tancogne-Dejean
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
27 use accel_oct_m
28 use atom_oct_m
29 use debug_oct_m
31 use global_oct_m
34 use, intrinsic :: iso_fortran_env
36 use math_oct_m
37 use mesh_oct_m
40 use mpi_oct_m
41 use mix_oct_m
47 use space_oct_m
49 use types_oct_m
50 use unit_oct_m
53
54 implicit none
55
56 private
57 public :: &
64
65 integer, public, parameter :: &
66 INDEPENDENT_PARTICLES = 2, &
67 hartree = 1, &
68 hartree_fock = 3, &
69 kohn_sham_dft = 4, &
71 rdmft = 7
72
73
75 private
76 integer :: np, nspin
77 integer :: theory_level = independent_particles
78 logical :: needs_vtau = .false.
79
80 real(real64), public, allocatable :: vhartree(:)
81 real(real64), public, allocatable :: vxc(:,:)
82 real(real64), public, allocatable :: vhxc(:,:)
83 ! !! other contributions: photon xc potential, constrained DFT, ...
84 ! !! Check v_ks_calc for more details
85 real(real64), allocatable :: vtau(:,:)
86 real(real64), allocatable :: lapl_vtau(:,:)
87
88 type(accel_mem_t) :: vtau_accel
89 type(accel_mem_t) :: lapl_vtau_accel
90
91 type(derivatives_t), pointer :: der => null()
92
93 contains
94 procedure :: init => ks_potential_init
95 procedure :: end => ks_potential_end
96 procedure :: set_vtau => ks_potential_set_vtau
97 procedure :: add_vhxc => ks_potential_add_vhxc
98 procedure :: dmult_vhxc => dks_potential_mult_vhxc
99 procedure :: zmult_vhxc => zks_potential_mult_vhxc
100 procedure :: load => ks_potential_load_vhxc
101 procedure :: dump => ks_potential_dump_vhxc
102 procedure :: init_interpolation => ks_potential_init_interpolation
103 procedure :: run_zero_iter => ks_potential_run_zero_iter
104 procedure :: interpolation_new => ks_potential_interpolation_new
105 procedure :: get_interpolated_potentials => ks_potential_get_interpolated_potentials
106 procedure :: set_interpolated_potentials => ks_potential_set_interpolated_potentials
107 procedure :: interpolate_potentials => ks_potential_interpolate_potentials
108 procedure :: dapply_vtau_psi => dks_potential_apply_vtau_psi
109 procedure :: zapply_vtau_psi => zks_potential_apply_vtau_psi
110 procedure :: dapply_lapl_vtau_psi => dks_potential_apply_lapl_vtau_psi
111 procedure :: zapply_lapl_vtau_psi => zks_potential_apply_lapl_vtau_psi
112 procedure :: dcurrent_mass_renormalization => dks_potential_current_mass_renormalization
113 procedure :: zcurrent_mass_renormalization => zks_potential_current_mass_renormalization
114 procedure :: output_potentials => ks_potential_output_potentials
115 procedure :: store_potentials => ks_potential_store_copy
116 procedure :: restore_potentials => ks_potential_restore_copy
117 procedure :: check_convergence => ks_potential_check_convergence
118 procedure :: perform_interpolation => ks_potential_perform_interpolation
119 procedure :: mix_potentials => ks_potential_mix_potentials
120 end type ks_potential_t
121
122 ! We sometimes need to store a copy of the potentials, say for interpolation
123 ! This is used to store the corresponding data, still keeping them private
125 private
126 real(real64), public, allocatable :: vhxc(:,:)
127 real(real64), allocatable :: vtau(:,:)
128 real(real64), allocatable :: lapl_vtau(:,:)
129 contains
130 procedure :: copy_vhxc_to_buffer => xc_copied_potentials_copy_vhxc_to_buffer
133
134contains
135
137 subroutine ks_potential_init(this, der, np, np_part, nspin, theory_level, needs_vtau)
138 class(ks_potential_t), intent(inout) :: this
139 type(derivatives_t), target, intent(in) :: der
140 integer, intent(in) :: np, np_part
141 integer, intent(in) :: nspin
142 integer, intent(in) :: theory_level
143 logical, intent(in) :: needs_vtau
144
145 push_sub(ks_potential_init)
146
147 this%theory_level = theory_level
148 this%needs_vtau = needs_vtau
149 this%np = np
150 this%nspin = nspin
151 this%der => der
152
153 ! In the case of spinors, vxc_11 = hm%vxc(:, 1), vxc_22 = hm%vxc(:, 2), Re(vxc_12) = hm%vxc(:. 3);
154 ! Im(vxc_12) = hm%vxc(:, 4)
155 safe_allocate(this%vhxc(1:np, 1:nspin))
156 this%vhxc(1:np, 1:nspin) = m_zero
157
158 if (theory_level /= independent_particles) then
159
160 safe_allocate(this%vhartree(1:np_part))
161 this%vhartree=m_zero
162
163 safe_allocate(this%vxc(1:np, 1:nspin))
164 this%vxc=m_zero
165
166 if (needs_vtau) then
167 safe_allocate(this%vtau(1:np_part, 1:nspin))
168 this%vtau = m_zero
169 safe_allocate(this%lapl_vtau(1:np, 1:nspin))
170 this%lapl_vtau = m_zero
172 call accel_create_buffer(this%vtau_accel, accel_mem_read_only, type_float, &
176 end if
177
178 end if
179 end if
180
182 end subroutine ks_potential_init
183
185 subroutine ks_potential_end(this)
186 class(ks_potential_t), intent(inout) :: this
190 safe_deallocate_a(this%vhxc)
191 safe_deallocate_a(this%vhartree)
192 safe_deallocate_a(this%vxc)
193 safe_deallocate_a(this%vtau)
194 safe_deallocate_a(this%lapl_vtau)
196 if (this%needs_vtau .and. accel_is_enabled()) then
197 call accel_release_buffer(this%vtau_accel)
198 call accel_release_buffer(this%lapl_vtau_accel)
199 end if
202 end subroutine ks_potential_end
205 subroutine ks_potential_set_vtau(this, vtau)
206 class(ks_potential_t), intent(inout) :: this
207 real(real64), contiguous, intent(in) :: vtau(:,:)
211 assert(this%needs_vtau)
213 call lalg_copy(this%np, this%nspin, vtau, this%vtau)
214
216
218 end subroutine ks_potential_set_vtau
222 class(ks_potential_t), intent(inout) :: this
224 integer :: offset, ispin
225
227
228 assert(this%needs_vtau)
229
230 do ispin = 1, this%nspin
231 call dderivatives_lapl(this%der, this%vtau(:, ispin), this%lapl_vtau(:, ispin))
232 end do
233
234 if (accel_is_enabled()) then
235 offset = 0
236 do ispin = 1, this%nspin
237 call accel_write_buffer(this%vtau_accel, this%np, this%vtau(:, ispin), offset = offset)
238 call accel_write_buffer(this%lapl_vtau_accel, this%np, this%lapl_vtau(:, ispin), offset = offset)
239 offset = offset + accel_padded_size(this%np)
240 end do
241 end if
242
245
247 subroutine ks_potential_add_vhxc(this, pot, nspin)
248 class(ks_potential_t), intent(in) :: this
249 real(real64), contiguous, intent(inout) :: pot(:,:)
250 integer, optional, intent(in) :: nspin
251
252 integer :: nspin_
253
254 push_sub(ks_potential_add_vhxc)
255
256 assert(not_in_openmp())
257 nspin_ = optional_default(nspin, this%nspin)
258 assert(size(pot, dim=1) >= this%np)
259 assert(size(pot, dim=2) >= nspin_)
260
261 !TODO: We could skip this for IPA
262
263 call lalg_axpy(this%np, nspin_, m_one, this%vhxc, pot)
264
265 pop_sub(ks_potential_add_vhxc)
266 end subroutine ks_potential_add_vhxc
267
268 ! -----------------------------------------------------------------
270 subroutine ks_potential_dump_vhxc(this, restart, space, mesh, ierr)
271 class(ks_potential_t), intent(in) :: this
272 type(restart_t), intent(in) :: restart
273 class(space_t), intent(in) :: space
274 class(mesh_t), intent(in) :: mesh
275 integer, intent(out) :: ierr
276
277 integer :: iunit, err, err2(2), isp
278 character(len=MAX_PATH_LEN) :: filename
279 character(len=100) :: lines(2)
280
281 push_sub(ks_potential_dump_vhxc)
282
283 ierr = 0
284
285 if (restart_skip(restart) .or. this%theory_level == independent_particles) then
287 return
288 end if
289
290 if (debug%info) then
291 message(1) = "Debug: Writing Vhxc restart."
292 call messages_info(1)
293 end if
294
295 !write the different components of the Hartree+XC potential
296 iunit = restart_open(restart, 'vhxc')
297 lines(1) = '# #spin #nspin filename'
298 lines(2) = '%vhxc'
299 call restart_write(restart, iunit, lines, 2, err)
300 if (err /= 0) ierr = ierr + 1
301
302 err2 = 0
303 do isp = 1, this%nspin
304 filename = get_filename_with_spin('vhxc', this%nspin, isp)
305 write(lines(1), '(i8,a,i8,a)') isp, ' | ', this%nspin, ' | "'//trim(adjustl(filename))//'"'
306 call restart_write(restart, iunit, lines, 1, err)
307 if (err /= 0) err2(1) = err2(1) + 1
308
309 call restart_write_mesh_function(restart, space, filename, mesh, this%vhxc(:,isp), err)
310 if (err /= 0) err2(2) = err2(2) + 1
311
312 end do
313 if (err2(1) /= 0) ierr = ierr + 2
314 if (err2(2) /= 0) ierr = ierr + 4
315
316 lines(1) = '%'
317 call restart_write(restart, iunit, lines, 1, err)
318 if (err /= 0) ierr = ierr + 4
319
320 ! MGGAs and hybrid MGGAs have an extra term that also needs to be dumped
321 if (this%needs_vtau) then
322 lines(1) = '# #spin #nspin filename'
323 lines(2) = '%vtau'
324 call restart_write(restart, iunit, lines, 2, err)
325 if (err /= 0) ierr = ierr + 8
326
327 err2 = 0
328 do isp = 1, this%nspin
329 filename = get_filename_with_spin('vtau', this%nspin, isp)
330 write(lines(1), '(i8,a,i8,a)') isp, ' | ', this%nspin, ' | "'//trim(adjustl(filename))//'"'
331 call restart_write(restart, iunit, lines, 1, err)
332 if (err /= 0) err2(1) = err2(1) + 16
333
334 call restart_write_mesh_function(restart, space, filename, mesh, this%vtau(:,isp), err)
335 if (err /= 0) err2(1) = err2(1) + 1
336
337 end do
338 if (err2(1) /= 0) ierr = ierr + 32
339 if (err2(2) /= 0) ierr = ierr + 64
341 lines(1) = '%'
342 call restart_write(restart, iunit, lines, 1, err)
343 if (err /= 0) ierr = ierr + 128
344 end if
345
346 call restart_close(restart, iunit)
347
348 if (debug%info) then
349 message(1) = "Debug: Writing Vhxc restart done."
350 call messages_info(1)
351 end if
352
354 end subroutine ks_potential_dump_vhxc
355
356
357 ! ---------------------------------------------------------
359 subroutine ks_potential_load_vhxc(this, restart, space, mesh, ierr)
360 class(ks_potential_t), intent(inout) :: this
361 type(restart_t), intent(in) :: restart
362 class(space_t), intent(in) :: space
363 class(mesh_t), intent(in) :: mesh
364 integer, intent(out) :: ierr
365
366 integer :: err, err2, isp
367 character(len=MAX_PATH_LEN) :: filename
368
369 push_sub(ks_potential_load_vhxc)
370
371 ierr = 0
372
373 if (restart_skip(restart) .or. this%theory_level == independent_particles) then
374 ierr = -1
376 return
377 end if
378
379 if (debug%info) then
380 message(1) = "Debug: Reading Vhxc restart."
381 call messages_info(1)
382 end if
383
384 err2 = 0
385 do isp = 1, this%nspin
386 filename = get_filename_with_spin('vhxc', this%nspin, isp)
387
388 call restart_read_mesh_function(restart, space, filename, mesh, this%vhxc(:,isp), err)
389 if (err /= 0) err2 = err2 + 1
390
391 end do
392 if (err2 /= 0) ierr = ierr + 1
393
394 ! MGGAs and hybrid MGGAs have an extra term that also needs to be read
395 err2 = 0
396 if (this%needs_vtau) then
397 do isp = 1, this%nspin
398 filename = get_filename_with_spin('vtau', this%nspin, isp)
399
400 call restart_read_mesh_function(restart, space, filename, mesh, this%vtau(:,isp), err)
401 if (err /= 0) err2 = err2 + 1
402
403 end do
404
405 if (err2 /= 0) ierr = ierr + 2
406 end if
407
408 if (debug%info) then
409 message(1) = "Debug: Reading Vhxc restart done."
410 call messages_info(1)
411 end if
412
414 end subroutine ks_potential_load_vhxc
415
417 subroutine ks_potential_init_interpolation(this, vksold, order)
418 class(ks_potential_t), intent(in) :: this
419 type(potential_interpolation_t), intent(inout) :: vksold
420 integer, optional, intent(in) :: order
421
423
424 call potential_interpolation_init(vksold, this%np, this%nspin, this%needs_vtau, order=order)
425
428
430 subroutine ks_potential_run_zero_iter(this, vksold)
431 class(ks_potential_t), intent(in) :: this
432 type(potential_interpolation_t), intent(inout) :: vksold
433
435
436 call potential_interpolation_run_zero_iter(vksold, this%np, this%nspin, this%vhxc, vtau=this%vtau)
437
439 end subroutine ks_potential_run_zero_iter
440
442 subroutine ks_potential_interpolation_new(this, vksold, current_time, dt)
443 class(ks_potential_t), intent(in) :: this
444 type(potential_interpolation_t), intent(inout) :: vksold
445 real(real64), intent(in) :: current_time
446 real(real64), intent(in) :: dt
447
449
450 call potential_interpolation_new(vksold, this%np, this%nspin, current_time, dt, this%vhxc, vtau=this%vtau)
451
453 end subroutine ks_potential_interpolation_new
454
458 subroutine ks_potential_get_interpolated_potentials(this, vksold, history, storage)
459 class(ks_potential_t), intent(inout) :: this
460 type(potential_interpolation_t), intent(inout) :: vksold
461 integer, intent(in) :: history
462 type(xc_copied_potentials_t), optional, intent(inout) :: storage
463
464 if (this%theory_level == independent_particles) return
465
467
468 if (.not. present(storage)) then
469 if (this%needs_vtau) then
470 call potential_interpolation_get(vksold, this%np, this%nspin, history, this%vhxc, vtau = this%vtau)
472 else
473 call potential_interpolation_get(vksold, this%np, this%nspin, history, this%vhxc)
474 end if
475 else
476 call ks_potential_storage_allocate(this, storage)
477 if (this%needs_vtau) then
478 call potential_interpolation_get(vksold, this%np, this%nspin, history, storage%vhxc, vtau = storage%vtau)
479 else
480 call potential_interpolation_get(vksold, this%np, this%nspin, history, storage%vhxc)
481 end if
482 end if
483
486
488 subroutine ks_potential_set_interpolated_potentials(this, vksold, history)
489 class(ks_potential_t), intent(in) :: this
490 type(potential_interpolation_t), intent(inout) :: vksold
491 integer, intent(in) :: history
492
493 if (this%theory_level == independent_particles) return
494
496
497 if (this%needs_vtau) then
498 call potential_interpolation_set(vksold, this%np, this%nspin, history, this%vhxc, vtau = this%vtau)
499 else
500 call potential_interpolation_set(vksold, this%np, this%nspin, history, this%vhxc)
501 end if
502
505
506
508 subroutine ks_potential_interpolate_potentials(this, vksold, order, current_time, dt, interpolation_time)
509 class(ks_potential_t), intent(inout) :: this
510 type(potential_interpolation_t), intent(inout) :: vksold
511 integer, intent(in) :: order
512 real(real64), intent(in) :: current_time
513 real(real64), intent(in) :: dt
514 real(real64), intent(in) :: interpolation_time
515
516 if (this%theory_level == independent_particles) return
517
519
520 if (this%needs_vtau) then
521 call potential_interpolation_interpolate(vksold, order, current_time, dt, interpolation_time, &
522 this%vhxc, vtau = this%vtau)
524 else
525 call potential_interpolation_interpolate(vksold, order, current_time, dt, interpolation_time, &
526 this%vhxc)
527 end if
528
531
532 ! ---------------------------------------------------------
533 subroutine vtau_set_vout(field, this)
534 type(mixfield_t), intent(inout) :: field
535 type(ks_potential_t), intent(in) :: this
536
537 push_sub(vtau_set_vout)
538
539 call mixfield_set_vout(field, this%vtau)
540
541 pop_sub(vtau_set_vout)
542 end subroutine vtau_set_vout
543
544 ! ---------------------------------------------------------
545 subroutine vtau_set_vin(field, this)
546 type(mixfield_t), intent(inout) :: field
547 type(ks_potential_t), intent(in) :: this
548
549 push_sub(vtau_set_vin)
550
551 call mixfield_set_vin(field, this%vtau)
552
553 pop_sub(vtau_set_vin)
554 end subroutine vtau_set_vin
555
556 ! ---------------------------------------------------------
557 subroutine vtau_get_vnew(field, this)
558 type(mixfield_t), intent(in) :: field
559 type(ks_potential_t), intent(inout) :: this
560
561 push_sub(vtau_get_vnew)
562
563 call mixfield_get_vnew(field, this%vtau)
564
566
567 pop_sub(vtau_get_vnew)
568 end subroutine vtau_get_vnew
569
571 subroutine ks_potential_output_potentials(this, namespace, how, dir, space, mesh, pos, atoms, grp)
572 class(ks_potential_t), intent(in) :: this
573 type(namespace_t), intent(in) :: namespace
574 integer(int64), intent(in) :: how
575 character(len=*), intent(in) :: dir
576 class(space_t), intent(in) :: space
577 class(mesh_t), intent(in) :: mesh
578 real(real64), optional, intent(in) :: pos(:,:)
579 type(atom_t), optional, intent(in) :: atoms(:)
580 type(mpi_grp_t), optional, intent(in) :: grp
582 character(len=MAX_PATH_LEN) :: fname
583 integer :: is, err
584
586
587 call dio_function_output(how, dir, 'vh', namespace, &
588 space, mesh, this%vhartree, units_out%energy, err, pos=pos, atoms=atoms, grp = grp)
589
590 do is = 1, this%nspin
591 fname = get_filename_with_spin('vxc', this%nspin, is)
592 call dio_function_output(how, dir, fname, namespace, space, &
593 mesh, this%vxc(:, is), units_out%energy, err, pos=pos, atoms=atoms, grp = grp)
594
595 if (this%needs_vtau) then
596 fname = get_filename_with_spin('vtau', this%nspin, is)
597 call dio_function_output(how, dir, fname, namespace, space, &
598 mesh, this%vtau(:, is), units_out%energy, err, pos=pos, atoms=atoms, grp = grp)
599 end if
600 end do
603 end subroutine ks_potential_output_potentials
604
606 subroutine ks_potential_storage_allocate(this, copy)
607 class(ks_potential_t), intent(in) :: this
608 type(xc_copied_potentials_t), intent(inout) :: copy
609
611
612 if (.not. allocated(copy%vhxc)) then
613 safe_allocate(copy%vhxc(1:this%np, 1:this%nspin))
614 end if
615 if (this%needs_vtau) then
616 if (.not. allocated(copy%vtau)) then
617 safe_allocate(copy%vtau(1:this%np, 1:this%nspin))
618 end if
619 end if
620
622 end subroutine ks_potential_storage_allocate
623
625 subroutine ks_potential_store_copy(this, copy)
626 class(ks_potential_t), intent(in) :: this
627 type(xc_copied_potentials_t), intent(inout) :: copy
628
629 if (this%theory_level == independent_particles) return
630
632
633 call ks_potential_storage_allocate(this, copy)
634 call lalg_copy(this%np, this%nspin, this%vhxc, copy%vhxc)
635 if (this%needs_vtau) then
636 call lalg_copy(this%np, this%nspin, this%vtau, copy%vtau)
637 end if
640 end subroutine ks_potential_store_copy
641
643 subroutine ks_potential_restore_copy(this, copy)
644 class(ks_potential_t), intent(inout) :: this
645 type(xc_copied_potentials_t), intent(in) :: copy
646
647 if (this%theory_level == independent_particles) return
648
651 assert(allocated(copy%vhxc))
652 call lalg_copy(this%np, this%nspin, copy%vhxc, this%vhxc)
653 if (this%needs_vtau) then
654 assert(allocated(copy%vtau))
655 call this%set_vtau(copy%vtau)
656 end if
657
659 end subroutine ks_potential_restore_copy
660
662 subroutine xc_copied_potentials_copy_vhxc_to_buffer(this, np, nspin, pnp, buffer)
663 class(xc_copied_potentials_t), intent(in) :: this
664 integer(int64), intent(in) :: np
665 integer, intent(in) :: nspin
666 integer(int64), intent(in) :: pnp
667 type(accel_mem_t), intent(inout) :: buffer
668
669 integer :: ispin
670
672
673 do ispin = 1, nspin
674 call accel_write_buffer(buffer, np, this%vhxc(:, ispin), offset = (ispin - 1)*pnp)
675 end do
676
679
681 subroutine xc_copied_potentials_end(this)
682 type(xc_copied_potentials_t), intent(inout) :: this
683
684 safe_deallocate_a(this%vhxc)
685 safe_deallocate_a(this%vtau)
686 end subroutine xc_copied_potentials_end
687
694 real(real64) function ks_potential_check_convergence(this, copy, mesh, rho, qtot) result(diff)
695 class(ks_potential_t), intent(in) :: this
696 type(xc_copied_potentials_t), intent(in) :: copy
697 class(mesh_t), intent(in) :: mesh
698 real(real64), intent(in) :: rho(:,:)
699 real(real64), intent(in) :: qtot
700
701 real(real64), allocatable :: diff_pot(:)
702 integer :: ip, is
703
704 push_sub(ks_potential_check_convergence)
705
706 safe_allocate(diff_pot(1:this%np))
707
708 !$omp parallel
709 !$omp do
710 do ip = 1, this%np
711 diff_pot(ip) = abs(copy%vhxc(ip, 1) - this%vhxc(ip, 1))*rho(ip, 1)
712 end do
713 do is = 2, this%nspin
714 !$omp do
715 do ip = 1, this%np
716 diff_pot(ip) = diff_pot(ip) + abs(copy%vhxc(ip, is) - this%vhxc(ip, is))*rho(ip, is)
717 end do
718 end do
719 !$omp end parallel
720 diff = dmf_integrate(mesh, diff_pot) / qtot
721
722 safe_deallocate_a(diff_pot)
723
724 pop_sub(ks_potential_check_convergence)
726
728 subroutine ks_potential_perform_interpolation(this, vksold, times, current_time)
729 class(ks_potential_t), intent(inout) :: this
730 type(potential_interpolation_t), intent(inout) :: vksold
731 real(real64), intent(in) :: times(:)
732 real(real64), intent(in) :: current_time
733
735
736 assert(size(times) == 3)
737
738 call interpolate(times, vksold%v_old(:, :, 1:3), current_time, vksold%v_old(:, :, 0))
739 if (this%needs_vtau) then
740 call interpolate(times, vksold%vtau_old(:, :, 1:3), current_time, vksold%vtau_old(:, :, 0))
741 end if
742
745
747 subroutine ks_potential_mix_potentials(this, vold, dt)
748 class(ks_potential_t), intent(in) :: this
749 type(xc_copied_potentials_t), intent(inout) :: vold
750 real(real64), intent(in) :: dt
751
752 integer :: ispin, ip
753
756 assert(not_in_openmp())
757
758 !$omp parallel private(ip)
759 do ispin = 1, this%nspin
760 !$omp do simd
761 do ip = 1, this%np
762 vold%vhxc(ip, ispin) = m_half*dt*(this%vhxc(ip, ispin) - vold%vhxc(ip, ispin))
763 end do
764
765 if (this%needs_vtau) then
766 !$omp do simd
767 do ip = 1, this%np
768 vold%vtau(ip, ispin) = m_half*dt*(this%vtau(ip, ispin) - vold%vtau(ip, ispin))
769 end do
770 end if
771 end do
772 !$omp end parallel
773
775 end subroutine ks_potential_mix_potentials
776
777
778#include "undef.F90"
779#include "complex.F90"
780#include "ks_potential_inc.F90"
781
782#include "undef.F90"
783#include "real.F90"
784#include "ks_potential_inc.F90"
785
786end module ks_potential_oct_m
788!! Local Variables:
789!! mode: f90
790!! coding: utf-8
791!! End:
constant times a vector plus a vector
Definition: lalg_basic.F90:171
Copies a vector x, to a vector y.
Definition: lalg_basic.F90:186
subroutine, public accel_release_buffer(this)
Definition: accel.F90:1250
pure logical function, public accel_is_enabled()
Definition: accel.F90:402
integer, parameter, public accel_mem_read_only
Definition: accel.F90:183
type(debug_t), save, public debug
Definition: debug.F90:156
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
subroutine, public dderivatives_lapl(der, ff, op_ff, ghost_update, set_bc, factor)
apply the Laplacian to a mesh function
real(real64), parameter, public m_zero
Definition: global.F90:188
logical pure function, public not_in_openmp()
Definition: global.F90:454
real(real64), parameter, public m_one
Definition: global.F90:189
subroutine, public dio_function_output(how, dir, fname, namespace, space, mesh, ff, unit, ierr, pos, atoms, grp, root)
Top-level IO routine for functions defined on the mesh.
A module to handle KS potential, without the external potential.
subroutine ks_potential_set_interpolated_potentials(this, vksold, history)
Set the interpolated potentials to history.
subroutine zks_potential_current_mass_renormalization(this, gpsi, space_dim, ndim, ispin)
Nonlocal contribution of vtau for the current.
subroutine dks_potential_apply_lapl_vtau_psi(this, mesh, d, ispin, psib, vpsib)
Wrapper to hamiltonian_elec_base_local_sub to hide the data of lapl_vtau.
subroutine ks_potential_add_vhxc(this, pot, nspin)
Adds vHxc to the potential.
subroutine dks_potential_current_mass_renormalization(this, gpsi, space_dim, ndim, ispin)
Nonlocal contribution of vtau for the current.
subroutine ks_potential_load_vhxc(this, restart, space, mesh, ierr)
Loads the vhxc potentials.
subroutine, public vtau_set_vin(field, this)
subroutine ks_potential_get_interpolated_potentials(this, vksold, history, storage)
Get the interpolated potentials from history.
subroutine zks_potential_apply_vtau_psi(this, mesh, d, ispin, psib, vpsib)
Wrapper to hamiltonian_elec_base_local_sub to hide the data of vtau.
integer, parameter, public rdmft
subroutine, public vtau_set_vout(field, this)
subroutine zks_potential_mult_vhxc(this, mf, ispin)
Multiply a mesh function by vHxc.
subroutine ks_potential_init(this, der, np, np_part, nspin, theory_level, needs_vtau)
Allocate the memory for the KS potentials.
integer, parameter, public hartree
subroutine ks_potential_interpolate_potentials(this, vksold, order, current_time, dt, interpolation_time)
Interpolate potentials to a new time.
integer, parameter, public hartree_fock
subroutine ks_potential_update_vtau_buffer(this)
Update vtau GPU buffer.
subroutine ks_potential_run_zero_iter(this, vksold)
Run zero iter for the interpolation.
subroutine ks_potential_store_copy(this, copy)
Copy the potentials to a storage object.
subroutine ks_potential_mix_potentials(this, vold, dt)
Replace vold potentials by 0.5*dt(vold + vhxc)
subroutine ks_potential_perform_interpolation(this, vksold, times, current_time)
Perform a time interpolation of the potentials.
subroutine, public vtau_get_vnew(field, this)
subroutine ks_potential_init_interpolation(this, vksold, order)
Initialize the potential interpolation.
real(real64) function ks_potential_check_convergence(this, copy, mesh, rho, qtot)
Check the convergence of a vhxc for predictor-corrector.
subroutine ks_potential_restore_copy(this, copy)
Copy the potentials from a storage object.
subroutine ks_potential_storage_allocate(this, copy)
Copy the potentials to a storage object.
subroutine ks_potential_interpolation_new(this, vksold, current_time, dt)
New interpolation point for the interpolation.
integer, parameter, public generalized_kohn_sham_dft
subroutine ks_potential_end(this)
Releases the memory for the KS potentials.
subroutine dks_potential_mult_vhxc(this, mf, ispin)
Multiply a mesh function by vHxc.
subroutine, public xc_copied_potentials_end(this)
Finalizer for the copied potentials.
integer, parameter, public kohn_sham_dft
subroutine ks_potential_dump_vhxc(this, restart, space, mesh, ierr)
Dumps the vhxc potentials.
subroutine dks_potential_apply_vtau_psi(this, mesh, d, ispin, psib, vpsib)
Wrapper to hamiltonian_elec_base_local_sub to hide the data of vtau.
subroutine ks_potential_set_vtau(this, vtau)
Set vtau and update the corresponding GPU buffer.
subroutine xc_copied_potentials_copy_vhxc_to_buffer(this, np, nspin, pnp, buffer)
Copy the vhxc potential to a gpu buffer.
subroutine ks_potential_output_potentials(this, namespace, how, dir, space, mesh, pos, atoms, grp)
Outputs vh, vxc, and vtau potentials.
subroutine zks_potential_apply_lapl_vtau_psi(this, mesh, d, ispin, psib, vpsib)
Wrapper to hamiltonian_elec_base_local_sub to hide the data of lapl_vtau.
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:616
this module contains the low-level part of the output system
Definition: output_low.F90:115
character(len=max_path_len) function, public get_filename_with_spin(output, nspin, spin_index)
Returns the filame as output, or output-spX is spin polarized.
Definition: output_low.F90:232
subroutine, public potential_interpolation_interpolate(potential_interpolation, order, time, dt, t, vhxc, vtau)
subroutine, public potential_interpolation_set(potential_interpolation, np, nspin, i, vhxc, vtau)
subroutine, public potential_interpolation_new(potential_interpolation, np, nspin, time, dt, vhxc, vtau)
subroutine, public potential_interpolation_run_zero_iter(potential_interpolation, np, nspin, vhxc, vtau)
subroutine, public potential_interpolation_get(potential_interpolation, np, nspin, i, vhxc, vtau)
subroutine, public potential_interpolation_init(potential_interpolation, np, nspin, mgga_with_exc, order)
subroutine, public restart_close(restart, iunit)
Close a file previously opened with restart_open.
Definition: restart.F90:953
subroutine, public restart_write(restart, iunit, lines, nlines, ierr)
Definition: restart.F90:909
logical pure function, public restart_skip(restart)
Returns true if the restart information should neither be read nor written. This might happen because...
Definition: restart.F90:970
integer function, public restart_open(restart, filename, status, position, silent)
Open file "filename" found inside the current restart directory. Depending on the type of restart,...
Definition: restart.F90:862
This module handles spin dimensions of the states and the k-point distribution.
type(type_t), public type_float
Definition: types.F90:133
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
This module defines the unit system, used for input and output.
type(unit_system_t), public units_out
Describes mesh distribution to nodes.
Definition: mesh.F90:186
Quantities used in mixing: Input, output and new potentials, and the residuals.
Definition: mix.F90:171
This is defined even when running serial.
Definition: mpi.F90:142