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 :: &
65
67 private
68 integer :: np, nspin
69 integer :: theory_level = independent_particles
70 logical :: needs_vtau = .false.
71
72 real(real64), public, allocatable :: vhartree(:)
73 real(real64), public, allocatable :: vxc(:,:)
74 real(real64), public, allocatable :: vhxc(:,:)
75 ! !! other contributions: photon xc potential, constrained DFT, ...
76 ! !! Check v_ks_calc for more details
77 real(real64), allocatable :: vtau(:,:)
78 real(real64), allocatable :: lapl_vtau(:,:)
79
80 type(accel_mem_t) :: vtau_accel
81 type(accel_mem_t) :: lapl_vtau_accel
82
83 type(derivatives_t), pointer :: der => null()
84
85 contains
86 procedure :: init => ks_potential_init
87 procedure :: end => ks_potential_end
88 procedure :: set_vtau => ks_potential_set_vtau
89 procedure :: add_vhxc => ks_potential_add_vhxc
90 procedure :: dmult_vhxc => dks_potential_mult_vhxc
91 procedure :: zmult_vhxc => zks_potential_mult_vhxc
92 procedure :: load => ks_potential_load_vhxc
93 procedure :: dump => ks_potential_dump_vhxc
94 procedure :: init_interpolation => ks_potential_init_interpolation
95 procedure :: run_zero_iter => ks_potential_run_zero_iter
96 procedure :: interpolation_new => ks_potential_interpolation_new
97 procedure :: get_interpolated_potentials => ks_potential_get_interpolated_potentials
98 procedure :: set_interpolated_potentials => ks_potential_set_interpolated_potentials
99 procedure :: interpolate_potentials => ks_potential_interpolate_potentials
100 procedure :: dapply_vtau_psi => dks_potential_apply_vtau_psi
101 procedure :: zapply_vtau_psi => zks_potential_apply_vtau_psi
102 procedure :: dapply_lapl_vtau_psi => dks_potential_apply_lapl_vtau_psi
103 procedure :: zapply_lapl_vtau_psi => zks_potential_apply_lapl_vtau_psi
104 procedure :: dcurrent_mass_renormalization => dks_potential_current_mass_renormalization
105 procedure :: zcurrent_mass_renormalization => zks_potential_current_mass_renormalization
106 procedure :: output_potentials => ks_potential_output_potentials
107 procedure :: store_potentials => ks_potential_store_copy
108 procedure :: restore_potentials => ks_potential_restore_copy
109 procedure :: check_convergence => ks_potential_check_convergence
110 procedure :: perform_interpolation => ks_potential_perform_interpolation
111 procedure :: mix_potentials => ks_potential_mix_potentials
112 end type ks_potential_t
113
114 ! We sometimes need to store a copy of the potentials, say for interpolation
115 ! This is used to store the corresponding data, still keeping them private
117 private
118 real(real64), public, allocatable :: vhxc(:,:)
119 real(real64), allocatable :: vtau(:,:)
120 real(real64), allocatable :: lapl_vtau(:,:)
121 contains
122 procedure :: copy_vhxc_to_buffer => xc_copied_potentials_copy_vhxc_to_buffer
125
126contains
127
129 subroutine ks_potential_init(this, der, np, np_part, nspin, theory_level, needs_vtau)
130 class(ks_potential_t), intent(inout) :: this
131 type(derivatives_t), target, intent(in) :: der
132 integer, intent(in) :: np, np_part
133 integer, intent(in) :: nspin
134 integer, intent(in) :: theory_level
135 logical, intent(in) :: needs_vtau
136
137 push_sub(ks_potential_init)
138
139 this%theory_level = theory_level
140 this%needs_vtau = needs_vtau
141 this%np = np
142 this%nspin = nspin
143 this%der => der
144
145 ! In the case of spinors, vxc_11 = hm%vxc(:, 1), vxc_22 = hm%vxc(:, 2), Re(vxc_12) = hm%vxc(:. 3);
146 ! Im(vxc_12) = hm%vxc(:, 4)
147 safe_allocate(this%vhxc(1:np, 1:nspin))
148 this%vhxc(1:np, 1:nspin) = m_zero
149
150 if (theory_level /= independent_particles) then
151
152 safe_allocate(this%vhartree(1:np_part))
153 this%vhartree=m_zero
154
155 safe_allocate(this%vxc(1:np, 1:nspin))
156 this%vxc=m_zero
157
158 if (needs_vtau) then
159 safe_allocate(this%vtau(1:np_part, 1:nspin))
160 this%vtau = m_zero
161 safe_allocate(this%lapl_vtau(1:np, 1:nspin))
162 this%lapl_vtau = m_zero
166 call accel_create_buffer(this%lapl_vtau_accel, accel_mem_read_only, type_float, &
168 end if
170 end if
171 end if
174 end subroutine ks_potential_init
177 subroutine ks_potential_end(this)
178 class(ks_potential_t), intent(inout) :: this
179
180 push_sub(ks_potential_end)
182 safe_deallocate_a(this%vhxc)
183 safe_deallocate_a(this%vhartree)
184 safe_deallocate_a(this%vxc)
185 safe_deallocate_a(this%vtau)
186 safe_deallocate_a(this%lapl_vtau)
188 if (this%needs_vtau .and. accel_is_enabled()) then
189 call accel_free_buffer(this%vtau_accel)
190 call accel_free_buffer(this%lapl_vtau_accel)
191 end if
194 end subroutine ks_potential_end
197 subroutine ks_potential_set_vtau(this, vtau)
198 class(ks_potential_t), intent(inout) :: this
199 real(real64), contiguous, intent(in) :: vtau(:,:)
203 assert(this%needs_vtau)
205 call lalg_copy(this%np, this%nspin, vtau, this%vtau)
208
209 pop_sub(ks_potential_set_vtau)
210 end subroutine ks_potential_set_vtau
214 class(ks_potential_t), intent(inout) :: this
216 integer :: offset, ispin
219
220 assert(this%needs_vtau)
221
222 do ispin = 1, this%nspin
223 call dderivatives_lapl(this%der, this%vtau(:, ispin), this%lapl_vtau(:, ispin))
224 end do
225
226 if (accel_is_enabled()) then
227 offset = 0
228 do ispin = 1, this%nspin
229 call accel_write_buffer(this%vtau_accel, this%np, this%vtau(:, ispin), offset = offset)
230 call accel_write_buffer(this%lapl_vtau_accel, this%np, this%lapl_vtau(:, ispin), offset = offset)
231 offset = offset + accel_padded_size(this%np)
232 end do
233 end if
234
237
239 subroutine ks_potential_accel_rebuild(this)
240 class(ks_potential_t), intent(inout) :: this
241
242 integer :: offset, ispin
243
245
246 if (.not. accel_is_enabled()) then
248 return
249 end if
250
251 call accel_detach_buffer(this%vtau_accel)
252 call accel_detach_buffer(this%lapl_vtau_accel)
253
254 if (.not. this%needs_vtau) then
256 return
257 end if
258
259 assert(allocated(this%vtau))
260 assert(allocated(this%lapl_vtau))
261
262 call accel_create_buffer(this%vtau_accel, accel_mem_read_only, type_float, &
263 accel_padded_size(this%np)*this%nspin)
264 call accel_create_buffer(this%lapl_vtau_accel, accel_mem_read_only, type_float, &
265 accel_padded_size(this%np)*this%nspin)
266
267 offset = 0
268 do ispin = 1, this%nspin
269 call accel_write_buffer(this%vtau_accel, this%np, this%vtau(:, ispin), offset = offset)
270 call accel_write_buffer(this%lapl_vtau_accel, this%np, this%lapl_vtau(:, ispin), offset = offset)
271 offset = offset + accel_padded_size(this%np)
272 end do
273
275 end subroutine ks_potential_accel_rebuild
276
278 subroutine ks_potential_add_vhxc(this, pot, nspin)
279 class(ks_potential_t), intent(in) :: this
280 real(real64), contiguous, intent(inout) :: pot(:,:)
281 integer, optional, intent(in) :: nspin
282
283 integer :: nspin_
284
285 push_sub(ks_potential_add_vhxc)
286
287 assert(not_in_openmp())
288 nspin_ = optional_default(nspin, this%nspin)
289 assert(size(pot, dim=1) >= this%np)
290 assert(size(pot, dim=2) >= nspin_)
291
292 !TODO: We could skip this for IPA
293
294 call lalg_axpy(this%np, nspin_, m_one, this%vhxc, pot)
295
296 pop_sub(ks_potential_add_vhxc)
297 end subroutine ks_potential_add_vhxc
298
299 ! -----------------------------------------------------------------
301 subroutine ks_potential_dump_vhxc(this, restart, mesh, ierr)
302 class(ks_potential_t), intent(in) :: this
303 type(restart_t), intent(in) :: restart
304 class(mesh_t), intent(in) :: mesh
305 integer, intent(out) :: ierr
306
307 integer :: iunit, err, err2(2), isp
308 character(len=MAX_PATH_LEN) :: filename
309 character(len=100) :: lines(2)
310
311 push_sub(ks_potential_dump_vhxc)
312
313 ierr = 0
314
315 if (restart%skip() .or. this%theory_level == independent_particles) then
317 return
318 end if
319
320 if (debug%info) then
321 message(1) = "Debug: Writing Vhxc restart."
322 call messages_info(1)
323 end if
324
325 !write the different components of the Hartree+XC potential
326 iunit = restart%open('vhxc')
327 lines(1) = '# #spin #nspin filename'
328 lines(2) = '%vhxc'
329 call restart%write(iunit, lines, 2, err)
330 if (err /= 0) ierr = ierr + 1
331
332 err2 = 0
333 do isp = 1, this%nspin
334 filename = get_filename_with_spin('vhxc', this%nspin, isp)
335 write(lines(1), '(i8,a,i8,a)') isp, ' | ', this%nspin, ' | "'//trim(adjustl(filename))//'"'
336 call restart%write(iunit, lines, 1, err)
337 if (err /= 0) err2(1) = err2(1) + 1
338
339 call restart%write_mesh_function(filename, mesh, this%vhxc(:,isp), err)
340 if (err /= 0) err2(2) = err2(2) + 1
341
342 end do
343 if (err2(1) /= 0) ierr = ierr + 2
344 if (err2(2) /= 0) ierr = ierr + 4
345
346 lines(1) = '%'
347 call restart%write(iunit, lines, 1, err)
348 if (err /= 0) ierr = ierr + 4
349
350 ! MGGAs and hybrid MGGAs have an extra term that also needs to be dumped
351 if (this%needs_vtau) then
352 lines(1) = '# #spin #nspin filename'
353 lines(2) = '%vtau'
354 call restart%write(iunit, lines, 2, err)
355 if (err /= 0) ierr = ierr + 8
356
357 err2 = 0
358 do isp = 1, this%nspin
359 filename = get_filename_with_spin('vtau', this%nspin, isp)
360 write(lines(1), '(i8,a,i8,a)') isp, ' | ', this%nspin, ' | "'//trim(adjustl(filename))//'"'
361 call restart%write(iunit, lines, 1, err)
362 if (err /= 0) err2(1) = err2(1) + 16
363
364 call restart%write_mesh_function(filename, mesh, this%vtau(:,isp), err)
365 if (err /= 0) err2(1) = err2(1) + 1
366
367 end do
368 if (err2(1) /= 0) ierr = ierr + 32
369 if (err2(2) /= 0) ierr = ierr + 64
370
371 lines(1) = '%'
372 call restart%write(iunit, lines, 1, err)
373 if (err /= 0) ierr = ierr + 128
374 end if
375
376 call restart%close(iunit)
377
378 if (debug%info) then
379 message(1) = "Debug: Writing Vhxc restart done."
380 call messages_info(1)
381 end if
382
384 end subroutine ks_potential_dump_vhxc
385
386
387 ! ---------------------------------------------------------
389 subroutine ks_potential_load_vhxc(this, restart, mesh, ierr)
390 class(ks_potential_t), intent(inout) :: this
391 type(restart_t), intent(in) :: restart
392 class(mesh_t), intent(in) :: mesh
393 integer, intent(out) :: ierr
394
395 integer :: err, err2, isp
396 character(len=MAX_PATH_LEN) :: filename
397
398 push_sub(ks_potential_load_vhxc)
399
400 ierr = 0
401
402 if (restart%skip() .or. this%theory_level == independent_particles) then
403 ierr = -1
405 return
406 end if
407
408 if (debug%info) then
409 message(1) = "Debug: Reading Vhxc restart."
410 call messages_info(1)
411 end if
412
413 err2 = 0
414 do isp = 1, this%nspin
415 filename = get_filename_with_spin('vhxc', this%nspin, isp)
416
417 call restart%read_mesh_function(filename, mesh, this%vhxc(:,isp), err)
418 if (err /= 0) err2 = err2 + 1
419
420 end do
421 if (err2 /= 0) ierr = ierr + 1
422
423 ! MGGAs and hybrid MGGAs have an extra term that also needs to be read
424 err2 = 0
425 if (this%needs_vtau) then
426 do isp = 1, this%nspin
427 filename = get_filename_with_spin('vtau', this%nspin, isp)
428
429 call restart%read_mesh_function(filename, mesh, this%vtau(:,isp), err)
430 if (err /= 0) err2 = err2 + 1
431
432 end do
433
434 if (err2 /= 0) ierr = ierr + 2
435 end if
436
437 if (debug%info) then
438 message(1) = "Debug: Reading Vhxc restart done."
439 call messages_info(1)
440 end if
441
443 end subroutine ks_potential_load_vhxc
444
446 subroutine ks_potential_init_interpolation(this, vksold, order)
447 class(ks_potential_t), intent(in) :: this
448 type(potential_interpolation_t), intent(inout) :: vksold
449 integer, optional, intent(in) :: order
450
452
453 call potential_interpolation_init(vksold, this%np, this%nspin, this%needs_vtau, order=order)
454
457
459 subroutine ks_potential_run_zero_iter(this, vksold)
460 class(ks_potential_t), intent(in) :: this
461 type(potential_interpolation_t), intent(inout) :: vksold
462
464
465 call potential_interpolation_run_zero_iter(vksold, this%np, this%nspin, this%vhxc, vtau=this%vtau)
466
468 end subroutine ks_potential_run_zero_iter
469
471 subroutine ks_potential_interpolation_new(this, vksold, current_time, dt)
472 class(ks_potential_t), intent(in) :: this
473 type(potential_interpolation_t), intent(inout) :: vksold
474 real(real64), intent(in) :: current_time
475 real(real64), intent(in) :: dt
476
478
479 call potential_interpolation_new(vksold, this%np, this%nspin, current_time, dt, this%vhxc, vtau=this%vtau)
480
482 end subroutine ks_potential_interpolation_new
483
487 subroutine ks_potential_get_interpolated_potentials(this, vksold, history, storage)
488 class(ks_potential_t), intent(inout) :: this
489 type(potential_interpolation_t), intent(inout) :: vksold
490 integer, intent(in) :: history
491 type(xc_copied_potentials_t), optional, intent(inout) :: storage
492
493 if (this%theory_level == independent_particles) return
494
496
497 if (.not. present(storage)) then
498 if (this%needs_vtau) then
499 call potential_interpolation_get(vksold, this%np, this%nspin, history, this%vhxc, vtau = this%vtau)
501 else
502 call potential_interpolation_get(vksold, this%np, this%nspin, history, this%vhxc)
503 end if
504 else
505 call ks_potential_storage_allocate(this, storage)
506 if (this%needs_vtau) then
507 call potential_interpolation_get(vksold, this%np, this%nspin, history, storage%vhxc, vtau = storage%vtau)
508 else
509 call potential_interpolation_get(vksold, this%np, this%nspin, history, storage%vhxc)
510 end if
511 end if
512
515
517 subroutine ks_potential_set_interpolated_potentials(this, vksold, history)
518 class(ks_potential_t), intent(in) :: this
519 type(potential_interpolation_t), intent(inout) :: vksold
520 integer, intent(in) :: history
521
522 if (this%theory_level == independent_particles) return
523
525
526 if (this%needs_vtau) then
527 call potential_interpolation_set(vksold, this%np, this%nspin, history, this%vhxc, vtau = this%vtau)
528 else
529 call potential_interpolation_set(vksold, this%np, this%nspin, history, this%vhxc)
530 end if
531
534
535
537 subroutine ks_potential_interpolate_potentials(this, vksold, order, current_time, dt, interpolation_time)
538 class(ks_potential_t), intent(inout) :: this
539 type(potential_interpolation_t), intent(inout) :: vksold
540 integer, intent(in) :: order
541 real(real64), intent(in) :: current_time
542 real(real64), intent(in) :: dt
543 real(real64), intent(in) :: interpolation_time
544
545 if (this%theory_level == independent_particles) return
546
548
549 if (this%needs_vtau) then
550 call potential_interpolation_interpolate(vksold, order, current_time, dt, interpolation_time, &
551 this%vhxc, vtau = this%vtau)
553 else
554 call potential_interpolation_interpolate(vksold, order, current_time, dt, interpolation_time, &
555 this%vhxc)
556 end if
557
560
561 ! ---------------------------------------------------------
562 subroutine vtau_set_vout(field, this)
563 type(mixfield_t), intent(inout) :: field
564 type(ks_potential_t), intent(in) :: this
565
567
568 call mixfield_set_vout(field, this%vtau)
569
570 pop_sub(vtau_set_vout)
571 end subroutine vtau_set_vout
572
573 ! ---------------------------------------------------------
574 subroutine vtau_set_vin(field, this)
575 type(mixfield_t), intent(inout) :: field
576 type(ks_potential_t), intent(in) :: this
577
578 push_sub(vtau_set_vin)
579
580 call mixfield_set_vin(field, this%vtau)
581
583 end subroutine vtau_set_vin
584
585 ! ---------------------------------------------------------
586 subroutine vtau_get_vnew(field, this)
587 type(mixfield_t), intent(in) :: field
588 type(ks_potential_t), intent(inout) :: this
589
590 push_sub(vtau_get_vnew)
591
592 call mixfield_get_vnew(field, this%vtau)
593
595
596 pop_sub(vtau_get_vnew)
597 end subroutine vtau_get_vnew
598
600 subroutine ks_potential_output_potentials(this, namespace, how, dir, space, mesh, pos, atoms, grp)
601 class(ks_potential_t), intent(in) :: this
602 type(namespace_t), intent(in) :: namespace
603 integer(int64), intent(in) :: how
604 character(len=*), intent(in) :: dir
605 class(space_t), intent(in) :: space
606 class(mesh_t), intent(in) :: mesh
607 real(real64), optional, intent(in) :: pos(:,:)
608 type(atom_t), optional, intent(in) :: atoms(:)
609 type(mpi_grp_t), optional, intent(in) :: grp
610
611 character(len=MAX_PATH_LEN) :: fname
612 integer :: is, err
613
615
616 call dio_function_output(how, dir, 'vh', namespace, &
617 space, mesh, this%vhartree, units_out%energy, err, pos=pos, atoms=atoms, grp = grp)
618
619 do is = 1, this%nspin
620 fname = get_filename_with_spin('vxc', this%nspin, is)
621 call dio_function_output(how, dir, fname, namespace, space, &
622 mesh, this%vxc(:, is), units_out%energy, err, pos=pos, atoms=atoms, grp = grp)
623
624 fname = get_filename_with_spin('vhxc', this%nspin, is)
625 call dio_function_output(how, dir, fname, namespace, space, &
626 mesh, this%vhxc(:, is), units_out%energy, err, pos=pos, atoms=atoms, grp = grp)
627
628 if (this%needs_vtau) then
629 fname = get_filename_with_spin('vtau', this%nspin, is)
630 call dio_function_output(how, dir, fname, namespace, space, &
631 mesh, this%vtau(:, is), units_out%energy, err, pos=pos, atoms=atoms, grp = grp)
632 end if
633 end do
634
636 end subroutine ks_potential_output_potentials
637
639 subroutine ks_potential_storage_allocate(this, copy)
640 class(ks_potential_t), intent(in) :: this
641 type(xc_copied_potentials_t), intent(inout) :: copy
642
644
645 if (.not. allocated(copy%vhxc)) then
646 safe_allocate(copy%vhxc(1:this%np, 1:this%nspin))
647 end if
648 if (this%needs_vtau) then
649 if (.not. allocated(copy%vtau)) then
650 safe_allocate(copy%vtau(1:this%np, 1:this%nspin))
651 end if
652 end if
653
655 end subroutine ks_potential_storage_allocate
656
658 subroutine ks_potential_store_copy(this, copy)
659 class(ks_potential_t), intent(in) :: this
660 type(xc_copied_potentials_t), intent(inout) :: copy
661
662 if (this%theory_level == independent_particles) return
663
665
666 call ks_potential_storage_allocate(this, copy)
667 call lalg_copy(this%np, this%nspin, this%vhxc, copy%vhxc)
668 if (this%needs_vtau) then
669 call lalg_copy(this%np, this%nspin, this%vtau, copy%vtau)
670 end if
671
673 end subroutine ks_potential_store_copy
674
676 subroutine ks_potential_restore_copy(this, copy)
677 class(ks_potential_t), intent(inout) :: this
678 type(xc_copied_potentials_t), intent(in) :: copy
679
680 if (this%theory_level == independent_particles) return
683
684 assert(allocated(copy%vhxc))
685 call lalg_copy(this%np, this%nspin, copy%vhxc, this%vhxc)
686 if (this%needs_vtau) then
687 assert(allocated(copy%vtau))
688 call this%set_vtau(copy%vtau)
689 end if
690
692 end subroutine ks_potential_restore_copy
693
695 subroutine xc_copied_potentials_copy_vhxc_to_buffer(this, np, nspin, pnp, buffer)
696 class(xc_copied_potentials_t), intent(in) :: this
697 integer(int64), intent(in) :: np
698 integer, intent(in) :: nspin
699 integer(int64), intent(in) :: pnp
700 type(accel_mem_t), intent(inout) :: buffer
701
702 integer :: ispin
703
705
706 do ispin = 1, nspin
707 call accel_write_buffer(buffer, np, this%vhxc(:, ispin), offset = (ispin - 1)*pnp)
708 end do
709
712
714 subroutine xc_copied_potentials_end(this)
715 type(xc_copied_potentials_t), intent(inout) :: this
716
717 safe_deallocate_a(this%vhxc)
718 safe_deallocate_a(this%vtau)
719 end subroutine xc_copied_potentials_end
720
727 real(real64) function ks_potential_check_convergence(this, copy, mesh, rho, qtot) result(diff)
728 class(ks_potential_t), intent(in) :: this
729 type(xc_copied_potentials_t), intent(in) :: copy
730 class(mesh_t), intent(in) :: mesh
731 real(real64), intent(in) :: rho(:,:)
732 real(real64), intent(in) :: qtot
733
734 real(real64), allocatable :: diff_pot(:)
735 integer :: ip, is
736
737 push_sub(ks_potential_check_convergence)
738
739 safe_allocate(diff_pot(1:this%np))
740
741 !$omp parallel
742 !$omp do
743 do ip = 1, this%np
744 diff_pot(ip) = abs(copy%vhxc(ip, 1) - this%vhxc(ip, 1))*rho(ip, 1)
745 end do
746 do is = 2, this%nspin
747 !$omp do
748 do ip = 1, this%np
749 diff_pot(ip) = diff_pot(ip) + abs(copy%vhxc(ip, is) - this%vhxc(ip, is))*rho(ip, is)
750 end do
751 end do
752 !$omp end parallel
753 diff = dmf_integrate(mesh, diff_pot) / qtot
754
755 safe_deallocate_a(diff_pot)
756
757 pop_sub(ks_potential_check_convergence)
759
761 subroutine ks_potential_perform_interpolation(this, vksold, times, current_time)
762 class(ks_potential_t), intent(inout) :: this
763 type(potential_interpolation_t), intent(inout) :: vksold
764 real(real64), intent(in) :: times(:)
765 real(real64), intent(in) :: current_time
766
768
769 assert(size(times) == 3)
770
771 call interpolate(times, vksold%v_old(:, :, 1:3), current_time, vksold%v_old(:, :, 0))
772 if (this%needs_vtau) then
773 call interpolate(times, vksold%vtau_old(:, :, 1:3), current_time, vksold%vtau_old(:, :, 0))
774 end if
775
778
780 subroutine ks_potential_mix_potentials(this, vold, dt)
781 class(ks_potential_t), intent(in) :: this
782 type(xc_copied_potentials_t), intent(inout) :: vold
783 real(real64), intent(in) :: dt
784
785 integer :: ispin, ip
786
788
789 assert(not_in_openmp())
791 !$omp parallel private(ip)
792 do ispin = 1, this%nspin
793 !$omp do simd
794 do ip = 1, this%np
795 vold%vhxc(ip, ispin) = m_half*dt*(this%vhxc(ip, ispin) - vold%vhxc(ip, ispin))
796 end do
797
798 if (this%needs_vtau) then
799 !$omp do simd
800 do ip = 1, this%np
801 vold%vtau(ip, ispin) = m_half*dt*(this%vtau(ip, ispin) - vold%vtau(ip, ispin))
802 end do
803 end if
804 end do
805 !$omp end parallel
806
808 end subroutine ks_potential_mix_potentials
810
811#include "undef.F90"
812#include "complex.F90"
813#include "ks_potential_inc.F90"
814
815#include "undef.F90"
816#include "real.F90"
817#include "ks_potential_inc.F90"
818
819end module ks_potential_oct_m
820
821!! Local Variables:
822!! mode: f90
823!! coding: utf-8
824!! End:
constant times a vector plus a vector
Definition: lalg_basic.F90:173
Copies a vector x, to a vector y.
Definition: lalg_basic.F90:188
subroutine, public accel_free_buffer(this, async)
Definition: accel.F90:1005
subroutine, public accel_detach_buffer(this)
Clear a buffer handle without freeing device memory.
Definition: accel.F90:1049
pure logical function, public accel_is_enabled()
Definition: accel.F90:402
integer, parameter, public accel_mem_read_only
Definition: accel.F90:185
type(debug_t), save, public debug
Definition: debug.F90:158
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:191
integer, parameter, public independent_particles
Theory level.
Definition: global.F90:241
logical pure function, public not_in_openmp()
Definition: global.F90:530
real(real64), parameter, public m_one
Definition: global.F90:192
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, 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.
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.
subroutine ks_potential_load_vhxc(this, restart, mesh, ierr)
Loads the vhxc potentials.
subroutine ks_potential_interpolate_potentials(this, vksold, order, current_time, dt, interpolation_time)
Interpolate potentials to a new time.
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, public ks_potential_accel_rebuild(this)
Rebuild accelerator buffers after an intrinsic copy.
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_dump_vhxc(this, restart, mesh, ierr)
Dumps the vhxc potentials.
subroutine ks_potential_interpolation_new(this, vksold, current_time, dt)
New interpolation point for the interpolation.
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.
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:117
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:594
this module contains the low-level part of the output system
Definition: output_low.F90:117
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:233
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)
This module handles spin dimensions of the states and the k-point distribution.
type(type_t), public type_float
Definition: types.F90:135
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:134
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:187
Quantities used in mixing: Input, output and new potentials, and the residuals.
Definition: mix.F90:173
This is defined even when running serial.
Definition: mpi.F90:144