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
30 use global_oct_m
33 use, intrinsic :: iso_fortran_env
35 use math_oct_m
36 use mesh_oct_m
39 use mpi_oct_m
40 use mix_oct_m
46 use space_oct_m
48 use types_oct_m
49 use unit_oct_m
52
53 implicit none
54
55 private
56 public :: &
63
64 integer, public, parameter :: &
65 INDEPENDENT_PARTICLES = 2, &
66 hartree = 1, &
67 hartree_fock = 3, &
68 kohn_sham_dft = 4, &
70 rdmft = 7
71
72
74 private
75 integer :: np, nspin
76 integer :: theory_level = independent_particles
77 logical :: needs_vtau = .false.
78
79 real(real64), public, allocatable :: vhartree(:)
80 real(real64), public, allocatable :: vxc(:,:)
81 real(real64), public, allocatable :: vhxc(:,:)
82 ! !! other contributions: photon xc potential, constrained DFT, ...
83 ! !! Check v_ks_calc for more details
84 real(real64), allocatable :: vtau(:,:)
85
86 type(accel_mem_t) :: vtau_accel
87
88 contains
89 procedure :: init => ks_potential_init
90 procedure :: end => ks_potential_end
91 procedure :: set_vtau => ks_potential_set_vtau
92 procedure :: add_vhxc => ks_potential_add_vhxc
93 procedure :: dmult_vhxc => dks_potential_mult_vhxc
94 procedure :: zmult_vhxc => zks_potential_mult_vhxc
95 procedure :: load => ks_potential_load_vhxc
96 procedure :: dump => ks_potential_dump_vhxc
97 procedure :: init_interpolation => ks_potential_init_interpolation
98 procedure :: run_zero_iter => ks_potential_run_zero_iter
99 procedure :: interpolation_new => ks_potential_interpolation_new
100 procedure :: get_interpolated_potentials => ks_potential_get_interpolated_potentials
101 procedure :: set_interpolated_potentials => ks_potential_set_interpolated_potentials
102 procedure :: interpolate_potentials => ks_potential_interpolate_potentials
103 procedure :: dapply_vtau_psi => dks_potential_apply_vtau_psi
104 procedure :: zapply_vtau_psi => zks_potential_apply_vtau_psi
105 procedure :: dcurrent_mass_renormalization => dks_potential_current_mass_renormalization
106 procedure :: zcurrent_mass_renormalization => zks_potential_current_mass_renormalization
107 procedure :: output_potentials => ks_potential_output_potentials
108 procedure :: store_potentials => ks_potential_store_copy
109 procedure :: restore_potentials => ks_potential_restore_copy
110 procedure :: check_convergence => ks_potential_check_convergence
111 procedure :: perform_interpolation => ks_potential_perform_interpolation
112 procedure :: mix_potentials => ks_potential_mix_potentials
113 end type ks_potential_t
114
115 ! We sometimes need to store a copy of the potentials, say for interpolation
116 ! This is used to store the corresponding data, still keeping them private
118 private
119 real(real64), public, allocatable :: vhxc(:,:)
120 real(real64), allocatable :: 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, np, np_part, nspin, theory_level, needs_vtau)
130 class(ks_potential_t), intent(inout) :: this
131 integer, intent(in) :: np, np_part
132 integer, intent(in) :: nspin
133 integer, intent(in) :: theory_level
134 logical, intent(in) :: needs_vtau
135
136 push_sub(ks_potential_init)
137
138 this%theory_level = theory_level
139 this%needs_vtau = needs_vtau
140 this%np = np
141 this%nspin = nspin
142
143 ! In the case of spinors, vxc_11 = hm%vxc(:, 1), vxc_22 = hm%vxc(:, 2), Re(vxc_12) = hm%vxc(:. 3);
144 ! Im(vxc_12) = hm%vxc(:, 4)
145 safe_allocate(this%vhxc(1:np, 1:nspin))
146 this%vhxc(1:np, 1:nspin) = m_zero
147
148 if (theory_level /= independent_particles) then
149
150 safe_allocate(this%vhartree(1:np_part))
151 this%vhartree=m_zero
152
153 safe_allocate(this%vxc(1:np, 1:nspin))
154 this%vxc=m_zero
155
156 if (needs_vtau) then
157 safe_allocate(this%vtau(1:np, 1:nspin))
158 this%vtau=m_zero
159 if (accel_is_enabled()) then
161 end if
162
163 end if
164 end if
165
167 end subroutine ks_potential_init
170 subroutine ks_potential_end(this)
171 class(ks_potential_t), intent(inout) :: this
175 safe_deallocate_a(this%vhxc)
176 safe_deallocate_a(this%vhartree)
177 safe_deallocate_a(this%vxc)
178 safe_deallocate_a(this%vtau)
180 if (this%needs_vtau .and. accel_is_enabled()) then
181 call accel_release_buffer(this%vtau_accel)
182 end if
185 end subroutine ks_potential_end
188 subroutine ks_potential_set_vtau(this, vtau)
189 class(ks_potential_t), intent(inout) :: this
190 real(real64), contiguous, intent(in) :: vtau(:,:)
194 assert(this%needs_vtau)
196 call lalg_copy(this%np, this%nspin, vtau, this%vtau)
205 class(ks_potential_t), intent(inout) :: this
206
207 integer :: offset, ispin
208
211 assert(this%needs_vtau)
214 offset = 0
215 do ispin = 1, this%nspin
216 call accel_write_buffer(this%vtau_accel, this%np, this%vtau(:, ispin), offset = offset)
217 offset = offset + accel_padded_size(this%np)
218 end do
219 end if
220
223
225 subroutine ks_potential_add_vhxc(this, pot, nspin)
226 class(ks_potential_t), intent(in) :: this
227 real(real64), contiguous, intent(inout) :: pot(:,:)
228 integer, optional, intent(in) :: nspin
229
230 integer :: nspin_
231
232 push_sub(ks_potential_add_vhxc)
233
234 assert(not_in_openmp())
235 nspin_ = optional_default(nspin, this%nspin)
236 assert(size(pot, dim=1) >= this%np)
237 assert(size(pot, dim=2) >= nspin_)
238
239 !TODO: We could skip this for IPA
240
241 call lalg_axpy(this%np, nspin_, m_one, this%vhxc, pot)
242
243 pop_sub(ks_potential_add_vhxc)
244 end subroutine ks_potential_add_vhxc
245
246 ! -----------------------------------------------------------------
248 subroutine ks_potential_dump_vhxc(this, restart, space, mesh, ierr)
249 class(ks_potential_t), intent(in) :: this
250 type(restart_t), intent(in) :: restart
251 class(space_t), intent(in) :: space
252 class(mesh_t), intent(in) :: mesh
253 integer, intent(out) :: ierr
254
255 integer :: iunit, err, err2(2), isp
256 character(len=MAX_PATH_LEN) :: filename
257 character(len=100) :: lines(2)
258
259 push_sub(ks_potential_dump_vhxc)
260
261 ierr = 0
262
263 if (restart_skip(restart) .or. this%theory_level == independent_particles) then
265 return
266 end if
267
268 if (debug%info) then
269 message(1) = "Debug: Writing Vhxc restart."
270 call messages_info(1)
271 end if
272
273 !write the different components of the Hartree+XC potential
274 iunit = restart_open(restart, 'vhxc')
275 lines(1) = '# #spin #nspin filename'
276 lines(2) = '%vhxc'
277 call restart_write(restart, iunit, lines, 2, err)
278 if (err /= 0) ierr = ierr + 1
279
280 err2 = 0
281 do isp = 1, this%nspin
282 filename = get_filename_with_spin('vhxc', this%nspin, isp)
283 write(lines(1), '(i8,a,i8,a)') isp, ' | ', this%nspin, ' | "'//trim(adjustl(filename))//'"'
284 call restart_write(restart, iunit, lines, 1, err)
285 if (err /= 0) err2(1) = err2(1) + 1
286
287 call restart_write_mesh_function(restart, space, filename, mesh, this%vhxc(:,isp), err)
288 if (err /= 0) err2(2) = err2(2) + 1
289
290 end do
291 if (err2(1) /= 0) ierr = ierr + 2
292 if (err2(2) /= 0) ierr = ierr + 4
293
294 lines(1) = '%'
295 call restart_write(restart, iunit, lines, 1, err)
296 if (err /= 0) ierr = ierr + 4
298 ! MGGAs and hybrid MGGAs have an extra term that also needs to be dumped
299 if (this%needs_vtau) then
300 lines(1) = '# #spin #nspin filename'
301 lines(2) = '%vtau'
302 call restart_write(restart, iunit, lines, 2, err)
303 if (err /= 0) ierr = ierr + 8
304
305 err2 = 0
306 do isp = 1, this%nspin
307 filename = get_filename_with_spin('vtau', this%nspin, isp)
308 write(lines(1), '(i8,a,i8,a)') isp, ' | ', this%nspin, ' | "'//trim(adjustl(filename))//'"'
309 call restart_write(restart, iunit, lines, 1, err)
310 if (err /= 0) err2(1) = err2(1) + 16
311
312 call restart_write_mesh_function(restart, space, filename, mesh, this%vtau(:,isp), err)
313 if (err /= 0) err2(1) = err2(1) + 1
314
315 end do
316 if (err2(1) /= 0) ierr = ierr + 32
317 if (err2(2) /= 0) ierr = ierr + 64
319 lines(1) = '%'
320 call restart_write(restart, iunit, lines, 1, err)
321 if (err /= 0) ierr = ierr + 128
322 end if
323
324 call restart_close(restart, iunit)
325
326 if (debug%info) then
327 message(1) = "Debug: Writing Vhxc restart done."
328 call messages_info(1)
329 end if
330
332 end subroutine ks_potential_dump_vhxc
333
334
335 ! ---------------------------------------------------------
337 subroutine ks_potential_load_vhxc(this, restart, space, mesh, ierr)
338 class(ks_potential_t), intent(inout) :: this
339 type(restart_t), intent(in) :: restart
340 class(space_t), intent(in) :: space
341 class(mesh_t), intent(in) :: mesh
342 integer, intent(out) :: ierr
343
344 integer :: err, err2, isp
345 character(len=MAX_PATH_LEN) :: filename
346
347 push_sub(ks_potential_load_vhxc)
348
349 ierr = 0
350
351 if (restart_skip(restart) .or. this%theory_level == independent_particles) then
352 ierr = -1
354 return
355 end if
356
357 if (debug%info) then
358 message(1) = "Debug: Reading Vhxc restart."
359 call messages_info(1)
360 end if
361
362 err2 = 0
363 do isp = 1, this%nspin
364 filename = get_filename_with_spin('vhxc', this%nspin, isp)
365
366 call restart_read_mesh_function(restart, space, filename, mesh, this%vhxc(:,isp), err)
367 if (err /= 0) err2 = err2 + 1
368
369 end do
370 if (err2 /= 0) ierr = ierr + 1
371
372 ! MGGAs and hybrid MGGAs have an extra term that also needs to be read
373 err2 = 0
374 if (this%needs_vtau) then
375 do isp = 1, this%nspin
376 filename = get_filename_with_spin('vtau', this%nspin, isp)
377
378 call restart_read_mesh_function(restart, space, filename, mesh, this%vtau(:,isp), err)
379 if (err /= 0) err2 = err2 + 1
380
381 end do
382
383 if (err2 /= 0) ierr = ierr + 2
384 end if
385
386 if (debug%info) then
387 message(1) = "Debug: Reading Vhxc restart done."
388 call messages_info(1)
389 end if
390
392 end subroutine ks_potential_load_vhxc
393
395 subroutine ks_potential_init_interpolation(this, vksold, order)
396 class(ks_potential_t), intent(in) :: this
397 type(potential_interpolation_t), intent(inout) :: vksold
398 integer, optional, intent(in) :: order
399
401
402 call potential_interpolation_init(vksold, this%np, this%nspin, this%needs_vtau, order=order)
403
406
408 subroutine ks_potential_run_zero_iter(this, vksold)
409 class(ks_potential_t), intent(in) :: this
410 type(potential_interpolation_t), intent(inout) :: vksold
411
413
414 call potential_interpolation_run_zero_iter(vksold, this%np, this%nspin, this%vhxc, vtau=this%vtau)
415
417 end subroutine ks_potential_run_zero_iter
418
420 subroutine ks_potential_interpolation_new(this, vksold, current_time, dt)
421 class(ks_potential_t), intent(inout) :: this
422 type(potential_interpolation_t), intent(inout) :: vksold
423 real(real64), intent(in) :: current_time
424 real(real64), intent(in) :: dt
425
427
428 call potential_interpolation_new(vksold, this%np, this%nspin, current_time, dt, this%vhxc, vtau=this%vtau)
429
431 end subroutine ks_potential_interpolation_new
432
436 subroutine ks_potential_get_interpolated_potentials(this, vksold, history, storage)
437 class(ks_potential_t), intent(inout) :: this
438 type(potential_interpolation_t), intent(inout) :: vksold
439 integer, intent(in) :: history
440 type(xc_copied_potentials_t), optional, intent(inout) :: storage
441
442 if (this%theory_level == independent_particles) return
443
445
446 if (.not. present(storage)) then
447 if (this%needs_vtau) then
448 call potential_interpolation_get(vksold, this%np, this%nspin, history, this%vhxc, vtau = this%vtau)
450 else
451 call potential_interpolation_get(vksold, this%np, this%nspin, history, this%vhxc)
452 end if
453 else
454 call ks_potential_storage_allocate(this, storage)
455 if (this%needs_vtau) then
456 call potential_interpolation_get(vksold, this%np, this%nspin, history, storage%vhxc, vtau = storage%vtau)
457 else
458 call potential_interpolation_get(vksold, this%np, this%nspin, history, storage%vhxc)
459 end if
460 end if
461
464
466 subroutine ks_potential_set_interpolated_potentials(this, vksold, history)
467 class(ks_potential_t), intent(inout) :: this
468 type(potential_interpolation_t), intent(inout) :: vksold
469 integer, intent(in) :: history
470
471 if (this%theory_level == independent_particles) return
472
474
475 if (this%needs_vtau) then
476 call potential_interpolation_set(vksold, this%np, this%nspin, history, this%vhxc, vtau = this%vtau)
477 else
478 call potential_interpolation_set(vksold, this%np, this%nspin, history, this%vhxc)
479 end if
480
483
484
486 subroutine ks_potential_interpolate_potentials(this, vksold, order, current_time, dt, interpolation_time)
487 class(ks_potential_t), intent(inout) :: this
488 type(potential_interpolation_t), intent(inout) :: vksold
489 integer, intent(in) :: order
490 real(real64), intent(in) :: current_time
491 real(real64), intent(in) :: dt
492 real(real64), intent(in) :: interpolation_time
493
494 if (this%theory_level == independent_particles) return
495
497
498 if (this%needs_vtau) then
499 call potential_interpolation_interpolate(vksold, order, current_time, dt, interpolation_time, &
500 this%vhxc, vtau = this%vtau)
502 else
503 call potential_interpolation_interpolate(vksold, order, current_time, dt, interpolation_time, &
504 this%vhxc)
505 end if
506
509
510 ! ---------------------------------------------------------
511 subroutine vtau_set_vout(field, this)
512 type(mixfield_t), intent(inout) :: field
513 type(ks_potential_t), intent(in) :: this
514
515 push_sub(vtau_set_vout)
516
517 call mixfield_set_vout(field, this%vtau)
518
519 pop_sub(vtau_set_vout)
520 end subroutine vtau_set_vout
521
522 ! ---------------------------------------------------------
523 subroutine vtau_set_vin(field, this)
524 type(mixfield_t), intent(inout) :: field
525 type(ks_potential_t), intent(in) :: this
526
527 push_sub(vtau_set_vin)
528
529 call mixfield_set_vin(field, this%vtau)
530
531 pop_sub(vtau_set_vin)
532 end subroutine vtau_set_vin
533
534 ! ---------------------------------------------------------
535 subroutine vtau_get_vnew(field, this)
536 type(mixfield_t), intent(in) :: field
537 type(ks_potential_t), intent(inout) :: this
538
539 push_sub(vtau_get_vnew)
540
541 call mixfield_get_vnew(field, this%vtau)
542
544
545 pop_sub(vtau_get_vnew)
546 end subroutine vtau_get_vnew
547
549 subroutine ks_potential_output_potentials(this, namespace, how, dir, space, mesh, pos, atoms, grp)
550 class(ks_potential_t), intent(in) :: this
551 type(namespace_t), intent(in) :: namespace
552 integer(int64), intent(in) :: how
553 character(len=*), intent(in) :: dir
554 class(space_t), intent(in) :: space
555 class(mesh_t), intent(in) :: mesh
556 real(real64), optional, intent(in) :: pos(:,:)
557 type(atom_t), optional, intent(in) :: atoms(:)
558 type(mpi_grp_t), optional, intent(in) :: grp
560 character(len=MAX_PATH_LEN) :: fname
561 integer :: is, err
562
564
565 call dio_function_output(how, dir, 'vh', namespace, &
566 space, mesh, this%vhartree, units_out%energy, err, pos=pos, atoms=atoms, grp = grp)
567
568 do is = 1, this%nspin
569 fname = get_filename_with_spin('vxc', this%nspin, is)
570 call dio_function_output(how, dir, fname, namespace, space, &
571 mesh, this%vxc(:, is), units_out%energy, err, pos=pos, atoms=atoms, grp = grp)
572
573 if (this%needs_vtau) then
574 fname = get_filename_with_spin('vtau', this%nspin, is)
575 call dio_function_output(how, dir, fname, namespace, space, &
576 mesh, this%vtau(:, is), units_out%energy, err, pos=pos, atoms=atoms, grp = grp)
577 end if
578 end do
581 end subroutine ks_potential_output_potentials
582
584 subroutine ks_potential_storage_allocate(this, copy)
585 class(ks_potential_t), intent(in) :: this
586 type(xc_copied_potentials_t), intent(inout) :: copy
587
589
590 if (.not. allocated(copy%vhxc)) then
591 safe_allocate(copy%vhxc(1:this%np, 1:this%nspin))
592 end if
593 if (this%needs_vtau) then
594 if (.not. allocated(copy%vtau)) then
595 safe_allocate(copy%vtau(1:this%np, 1:this%nspin))
596 end if
597 end if
598
600 end subroutine ks_potential_storage_allocate
601
603 subroutine ks_potential_store_copy(this, copy)
604 class(ks_potential_t), intent(in) :: this
605 type(xc_copied_potentials_t), intent(inout) :: copy
606
607 if (this%theory_level == independent_particles) return
608
610
611 call ks_potential_storage_allocate(this, copy)
612 call lalg_copy(this%np, this%nspin, this%vhxc, copy%vhxc)
613 if (this%needs_vtau) then
614 call lalg_copy(this%np, this%nspin, this%vtau, copy%vtau)
615 end if
618 end subroutine ks_potential_store_copy
619
621 subroutine ks_potential_restore_copy(this, copy)
622 class(ks_potential_t), intent(inout) :: this
623 type(xc_copied_potentials_t), intent(in) :: copy
624
625 if (this%theory_level == independent_particles) return
626
629 assert(allocated(copy%vhxc))
630 call lalg_copy(this%np, this%nspin, copy%vhxc, this%vhxc)
631 if (this%needs_vtau) then
632 assert(allocated(copy%vtau))
633 call this%set_vtau(copy%vtau)
634 end if
635
637 end subroutine ks_potential_restore_copy
638
640 subroutine xc_copied_potentials_copy_vhxc_to_buffer(this, np, nspin, pnp, buffer)
641 class(xc_copied_potentials_t), intent(in) :: this
642 integer(int64), intent(in) :: np
643 integer, intent(in) :: nspin
644 integer(int64), intent(in) :: pnp
645 type(accel_mem_t), intent(inout) :: buffer
646
647 integer :: ispin
648
650
651 do ispin = 1, nspin
652 call accel_write_buffer(buffer, np, this%vhxc(:, ispin), offset = (ispin - 1)*pnp)
653 end do
654
657
659 subroutine xc_copied_potentials_end(this)
660 type(xc_copied_potentials_t), intent(inout) :: this
661
662 safe_deallocate_a(this%vhxc)
663 safe_deallocate_a(this%vtau)
664 end subroutine xc_copied_potentials_end
665
672 real(real64) function ks_potential_check_convergence(this, copy, mesh, rho, qtot) result(diff)
673 class(ks_potential_t), intent(in) :: this
674 type(xc_copied_potentials_t), intent(in) :: copy
675 class(mesh_t), intent(in) :: mesh
676 real(real64), intent(in) :: rho(:,:)
677 real(real64), intent(in) :: qtot
678
679 real(real64), allocatable :: diff_pot(:)
680 integer :: ip, is
681
682 push_sub(ks_potential_check_convergence)
683
684 safe_allocate(diff_pot(1:this%np))
685
686 !$omp parallel
687 !$omp do
688 do ip = 1, this%np
689 diff_pot(ip) = abs(copy%vhxc(ip, 1) - this%vhxc(ip, 1))*rho(ip, 1)
690 end do
691 do is = 2, this%nspin
692 !$omp do
693 do ip = 1, this%np
694 diff_pot(ip) = diff_pot(ip) + abs(copy%vhxc(ip, is) - this%vhxc(ip, is))*rho(ip, is)
695 end do
696 end do
697 !$omp end parallel
698 diff = dmf_integrate(mesh, diff_pot) / qtot
699
700 safe_deallocate_a(diff_pot)
701
702 pop_sub(ks_potential_check_convergence)
704
706 subroutine ks_potential_perform_interpolation(this, vksold, times, current_time)
707 class(ks_potential_t), intent(inout) :: this
708 type(potential_interpolation_t), intent(inout) :: vksold
709 real(real64), intent(in) :: times(:)
710 real(real64), intent(in) :: current_time
711
713
714 assert(size(times) == 3)
715
716 call interpolate(times, vksold%v_old(:, :, 1:3), current_time, vksold%v_old(:, :, 0))
717 if (this%needs_vtau) then
718 call interpolate(times, vksold%vtau_old(:, :, 1:3), current_time, vksold%vtau_old(:, :, 0))
719 end if
720
723
725 subroutine ks_potential_mix_potentials(this, vold, dt)
726 class(ks_potential_t), intent(in) :: this
727 type(xc_copied_potentials_t), intent(inout) :: vold
728 real(real64), intent(in) :: dt
729
730 integer :: ispin, ip
731
734 assert(not_in_openmp())
735
736 !$omp parallel private(ip)
737 do ispin = 1, this%nspin
738 !$omp do simd
739 do ip = 1, this%np
740 vold%vhxc(ip, ispin) = m_half*dt*(this%vhxc(ip, ispin) - vold%vhxc(ip, ispin))
741 end do
742
743 if (this%needs_vtau) then
744 !$omp do simd
745 do ip = 1, this%np
746 vold%vtau(ip, ispin) = m_half*dt*(this%vtau(ip, ispin) - vold%vtau(ip, ispin))
747 end do
748 end if
749 end do
750 !$omp end parallel
751
753 end subroutine ks_potential_mix_potentials
754
755
756#include "undef.F90"
757#include "complex.F90"
758#include "ks_potential_inc.F90"
759
760#include "undef.F90"
761#include "real.F90"
762#include "ks_potential_inc.F90"
763
764end module ks_potential_oct_m
766!! Local Variables:
767!! mode: f90
768!! coding: utf-8
769!! 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:1248
pure logical function, public accel_is_enabled()
Definition: accel.F90:401
integer, parameter, public accel_mem_read_only
Definition: accel.F90:183
type(debug_t), save, public debug
Definition: debug.F90:156
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 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.
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_init(this, np, np_part, nspin, theory_level, needs_vtau)
Allocate the memory for the KS potentials.
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.
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:952
subroutine, public restart_write(restart, iunit, lines, nlines, ierr)
Definition: restart.F90:908
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:969
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:861
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