Octopus
density.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2010 M. Marques, A. Castro, A. Rubio, G. Bertsch, X. Andrade
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
26!
27module density_oct_m
28 use accel_oct_m
29 use batch_oct_m
31 use iso_c_binding
32 use comm_oct_m
33 use debug_oct_m
36 use global_oct_m
37 use grid_oct_m
38 use, intrinsic :: iso_fortran_env
41 use math_oct_m
42 use mesh_oct_m
45 use mpi_oct_m
48 use parser_oct_m
50 use smear_oct_m
51 use space_oct_m
55 use types_oct_m
57
58 implicit none
59
60 private
61
62 public :: &
74
77 private
78 real(real64), pointer :: density(:, :)
79 type(states_elec_t), pointer :: st
80 type(grid_t), pointer :: gr
81 type(accel_mem_t) :: buff_density
82 integer(int64) :: pnp
83 logical :: packed
84 end type density_calc_t
85
86contains
87
92 !
93 subroutine density_calc_init(this, st, gr, density)
94 type(density_calc_t), intent(out) :: this
95 type(states_elec_t), target, intent(in) :: st
96 type(grid_t), target, intent(in) :: gr
97 real(real64), target, intent(out) :: density(:, :)
98
99 integer :: ip, ispin
100 logical :: correct_size
101
102 push_sub(density_calc_init)
103
104 this%st => st
105 this%gr => gr
106
107 this%density => density
108 !$omp parallel private(ip)
109 do ispin = 1, st%d%nspin
110 !$omp do
111 do ip = 1, gr%np
112 this%density(ip, ispin) = m_zero
113 end do
114 !$omp end do nowait
115 end do
116 !$omp end parallel
117
118 this%packed = .false.
119
120 correct_size = ubound(this%density, dim = 1) == this%gr%np .or. &
121 ubound(this%density, dim = 1) == this%gr%np_part
122 assert(correct_size)
123
124 pop_sub(density_calc_init)
125 end subroutine density_calc_init
126
127 ! ---------------------------------------------------
129 !
130 subroutine density_calc_pack(this, async)
131 type(density_calc_t), intent(inout) :: this
132 logical, optional, intent(in) :: async
133
134 push_sub(density_calc_pack)
135
136 this%packed = .true.
137 this%pnp = accel_padded_size(this%gr%np)
138 call accel_create_buffer(this%buff_density, accel_mem_read_write, type_float, this%pnp*this%st%d%nspin)
139
140 ! set to zero
141 call accel_set_buffer_to_zero(this%buff_density, type_float, this%pnp*this%st%d%nspin, async=async)
142
143 pop_sub(density_calc_pack)
144 end subroutine density_calc_pack
145
146 ! ---------------------------------------------------
150 !
151 subroutine density_calc_state(this, psib, istin)
152 type(density_calc_t), intent(inout) :: this
153 type(wfs_elec_t), intent(in) :: psib
154 integer, intent(in) :: istin
155
156 integer :: ist, ip, ispin, istin_, wgsize
157 complex(real64) :: term, psi1, psi2
158 real(real64), allocatable :: weight(:)
159 type(accel_mem_t) :: buff_weight
160 type(accel_kernel_t), pointer :: kernel
161
162 push_sub(density_calc_state)
163 call profiling_in("CALC_STATE_DENSITY")
164
165 ispin = this%st%d%get_spin_index(psib%ik)
166
167 istin_= 0
168
169 safe_allocate(weight(1:psib%nst))
170 do ist = 1, psib%nst
171 weight(ist) = this%st%kweights(psib%ik)*this%st%occ(psib%ist(ist), psib%ik)
172 if (psib%ist(ist) == istin) istin_ = ist
173 end do
175 if (abs(weight(istin_)) <= m_epsilon) then
176 return
177 pop_sub(density_calc_state)
178 end if
179
180 !fix ist for the remaining code
181 ist = istin_
182
183
184 select case (psib%status())
185 case (batch_not_packed)
186 select case (this%st%d%ispin)
188 if (states_are_real(this%st)) then
189 !$omp parallel do simd schedule(static)
190 do ip = 1, this%gr%np
191 this%density(ip, ispin) = this%density(ip, ispin) + weight(ist)*psib%dff(ip, 1, ist)**2
192 end do
193 else
194 !$omp parallel do schedule(static)
195 do ip = 1, this%gr%np
196 this%density(ip, ispin) = this%density(ip, ispin) + weight(ist)* &
197 real(conjg(psib%zff(ip, 1, ist))*psib%zff(ip, 1, ist), real64)
198 end do
199 end if
200 case (spinors)
201 !$omp parallel do schedule(static) private(psi1, psi2, term)
202 do ip = 1, this%gr%np
203 psi1 = psib%zff(ip, 1, ist)
204 psi2 = psib%zff(ip, 2, ist)
205 this%density(ip, 1) = this%density(ip, 1) + weight(ist)*real(conjg(psi1)*psi1, real64)
206 this%density(ip, 2) = this%density(ip, 2) + weight(ist)*real(conjg(psi2)*psi2, real64)
207 term = weight(ist)*psi1*conjg(psi2)
208 this%density(ip, 3) = this%density(ip, 3) + real(term, real64)
209 this%density(ip, 4) = this%density(ip, 4) + aimag(term)
210 end do
211 end select
212
213 case (batch_packed)
214
215 select case (this%st%d%ispin)
217 if (states_are_real(this%st)) then
218 !$omp parallel do schedule(static)
219 do ip = 1, this%gr%np
220 this%density(ip, ispin) = this%density(ip, ispin) + weight(ist)*psib%dff_pack(ist, ip)**2
221 end do
222 else
224 !$omp parallel do schedule(static)
225 do ip = 1, this%gr%np
226 this%density(ip, ispin) = this%density(ip, ispin) + weight(ist)* &
227 real(conjg(psib%zff_pack(ist, ip))*psib%zff_pack(ist, ip), real64)
228 end do
229 end if
230 case (spinors)
231 assert(mod(psib%nst_linear, 2) == 0)
232 !$omp parallel do schedule(static) private(ist, psi1, psi2, term)
233 do ip = 1, this%gr%np
234 psi1 = psib%zff_pack(2*ist - 1, ip)
235 psi2 = psib%zff_pack(2*ist, ip)
236 term = weight(ist)*psi1*conjg(psi2)
237
238 this%density(ip, 1) = this%density(ip, 1) + weight(ist)*real(conjg(psi1)*psi1, real64)
239 this%density(ip, 2) = this%density(ip, 2) + weight(ist)*real(conjg(psi2)*psi2, real64)
240 this%density(ip, 3) = this%density(ip, 3) + real(term, real64)
241 this%density(ip, 4) = this%density(ip, 4) + aimag(term)
242 end do
243 end select
246 if (.not. this%packed) call density_calc_pack(this, async=.true.)
247
248 call accel_create_buffer(buff_weight, accel_mem_read_only, type_float, psib%nst)
249 call accel_write_buffer(buff_weight, psib%nst, weight, async=.true.)
250
251 select case (this%st%d%ispin)
253 if (states_are_real(this%st)) then
254 kernel => kernel_density_real
255 else
256 kernel => kernel_density_complex
257 end if
258
259 call accel_set_kernel_arg(kernel, 0, psib%nst)
260 call accel_set_kernel_arg(kernel, 1, this%gr%np)
261 call accel_set_kernel_arg(kernel, 2, int(this%pnp*(ispin - 1), int32))
262 call accel_set_kernel_arg(kernel, 3, buff_weight)
263 call accel_set_kernel_arg(kernel, 4, psib%ff_device)
264 call accel_set_kernel_arg(kernel, 5, int(log2(psib%pack_size(1)), int32))
265 call accel_set_kernel_arg(kernel, 6, this%buff_density)
266
267 case (spinors)
268 kernel => kernel_density_spinors
269
270 call accel_set_kernel_arg(kernel, 0, psib%nst)
271 call accel_set_kernel_arg(kernel, 1, this%gr%np)
272 call accel_set_kernel_arg(kernel, 2, int(this%pnp, int32))
273 call accel_set_kernel_arg(kernel, 3, buff_weight)
274 call accel_set_kernel_arg(kernel, 4, psib%ff_device)
275 call accel_set_kernel_arg(kernel, 5, int(log2(psib%pack_size(1)), int32))
276 call accel_set_kernel_arg(kernel, 6, this%buff_density)
277 end select
278
279 wgsize = accel_kernel_workgroup_size(kernel)
280
281 call accel_kernel_run(kernel, (/pad(this%gr%np, wgsize)/), (/wgsize/))
282
283 call accel_finish()
284 call accel_release_buffer(buff_weight)
285 end select
286
287 safe_deallocate_a(weight)
288
289 call profiling_out("CALC_STATE_DENSITY")
290
291 pop_sub(density_calc_state)
292 end subroutine density_calc_state
293
294 ! ---------------------------------------------------
296 !
297 subroutine density_calc_accumulate(this, psib, async)
298 type(density_calc_t), intent(inout) :: this
299 type(wfs_elec_t), intent(in) :: psib
300 logical, optional, intent(in) :: async
301
302 integer :: ist, ip, ispin, wgsize
303 complex(real64) :: term, psi1, psi2
304 real(real64), allocatable :: weight(:)
305 type(accel_mem_t) :: buff_weight
306 type(accel_kernel_t), pointer :: kernel
307
309 call profiling_in("CALC_DENSITY")
310
311 ispin = this%st%d%get_spin_index(psib%ik)
312
313 safe_allocate(weight(1:psib%nst))
314 do ist = 1, psib%nst
315 weight(ist) = this%st%kweights(psib%ik)*this%st%occ(psib%ist(ist), psib%ik)
316 end do
317
318 select case (psib%status())
319 case (batch_not_packed)
320 select case (this%st%d%ispin)
322 if (states_are_real(this%st)) then
323 do ist = 1, psib%nst
324 if (abs(weight(ist)) <= m_epsilon) cycle
325 !$omp parallel do simd schedule(static)
326 do ip = 1, this%gr%np
327 this%density(ip, ispin) = this%density(ip, ispin) + weight(ist)*psib%dff(ip, 1, ist)**2
328 end do
329 end do
330 else
331 do ist = 1, psib%nst
332 if (abs(weight(ist)) <= m_epsilon) cycle
333 !$omp parallel do schedule(static)
334 do ip = 1, this%gr%np
335 this%density(ip, ispin) = this%density(ip, ispin) + weight(ist)* &
336 real(conjg(psib%zff(ip, 1, ist))*psib%zff(ip, 1, ist), real64)
337 end do
338 end do
339 end if
340 case (spinors)
341 do ist = 1, psib%nst
342 if (abs(weight(ist)) <= m_epsilon) cycle
343 !$omp parallel do schedule(static) private(psi1, psi2, term)
344 do ip = 1, this%gr%np
345 psi1 = psib%zff(ip, 1, ist)
346 psi2 = psib%zff(ip, 2, ist)
347 this%density(ip, 1) = this%density(ip, 1) + weight(ist)*real(conjg(psi1)*psi1, real64)
348 this%density(ip, 2) = this%density(ip, 2) + weight(ist)*real(conjg(psi2)*psi2, real64)
349 term = weight(ist)*psi1*conjg(psi2)
350 this%density(ip, 3) = this%density(ip, 3) + real(term, real64)
351 this%density(ip, 4) = this%density(ip, 4) + aimag(term)
352 end do
353 end do
354 end select
355
356 case (batch_packed)
357
358 select case (this%st%d%ispin)
360 if (states_are_real(this%st)) then
361 !$omp parallel do schedule(static) private(ist)
362 do ip = 1, this%gr%np
363 do ist = 1, psib%nst
364 this%density(ip, ispin) = this%density(ip, ispin) + weight(ist)*psib%dff_pack(ist, ip)**2
365 end do
366 end do
367 else
368 !$omp parallel do schedule(static) private(ist)
369 do ip = 1, this%gr%np
370 do ist = 1, psib%nst
371 this%density(ip, ispin) = this%density(ip, ispin) + weight(ist)* &
372 real(conjg(psib%zff_pack(ist, ip))*psib%zff_pack(ist, ip), real64)
373 end do
374 end do
375 end if
376 case (spinors)
377 assert(mod(psib%nst_linear, 2) == 0)
378 !$omp parallel do schedule(static) private(ist, psi1, psi2, term)
379 do ip = 1, this%gr%np
380 do ist = 1, psib%nst
381 psi1 = psib%zff_pack(2*ist - 1, ip)
382 psi2 = psib%zff_pack(2*ist, ip)
383 term = weight(ist)*psi1*conjg(psi2)
384
385 this%density(ip, 1) = this%density(ip, 1) + weight(ist)*real(conjg(psi1)*psi1, real64)
386 this%density(ip, 2) = this%density(ip, 2) + weight(ist)*real(conjg(psi2)*psi2, real64)
387 this%density(ip, 3) = this%density(ip, 3) + real(term, real64)
388 this%density(ip, 4) = this%density(ip, 4) + aimag(term)
389 end do
390 end do
391 end select
392
394 if (.not. this%packed) call density_calc_pack(this)
395
396 call accel_create_buffer(buff_weight, accel_mem_read_only, type_float, psib%nst)
397 call accel_write_buffer(buff_weight, psib%nst, weight, async=.true.)
398
399 select case (this%st%d%ispin)
401 if (states_are_real(this%st)) then
402 kernel => kernel_density_real
403 else
404 kernel => kernel_density_complex
405 end if
406
407 call accel_set_kernel_arg(kernel, 0, psib%nst)
408 call accel_set_kernel_arg(kernel, 1, this%gr%np)
409 call accel_set_kernel_arg(kernel, 2, int(this%pnp, int32)*(ispin - 1))
410 call accel_set_kernel_arg(kernel, 3, buff_weight)
411 call accel_set_kernel_arg(kernel, 4, psib%ff_device)
412 call accel_set_kernel_arg(kernel, 5, int(log2(psib%pack_size(1)), int32))
413 call accel_set_kernel_arg(kernel, 6, this%buff_density)
414
415 case (spinors)
416 kernel => kernel_density_spinors
417
418 call accel_set_kernel_arg(kernel, 0, psib%nst)
419 call accel_set_kernel_arg(kernel, 1, this%gr%np)
420 call accel_set_kernel_arg(kernel, 2, int(this%pnp, int32))
421 call accel_set_kernel_arg(kernel, 3, buff_weight)
422 call accel_set_kernel_arg(kernel, 4, psib%ff_device)
423 call accel_set_kernel_arg(kernel, 5, int(log2(psib%pack_size(1)), int32))
424 call accel_set_kernel_arg(kernel, 6, this%buff_density)
425 end select
426
427 wgsize = accel_kernel_workgroup_size(kernel)
428
429 call accel_kernel_run(kernel, (/pad(this%gr%np, wgsize)/), (/wgsize/))
430
431 call accel_release_buffer(buff_weight)
432
433 if(.not. optional_default(async, .false.)) call accel_finish()
434
435 end select
436
437 safe_deallocate_a(weight)
438
439
440 if (this%st%d%ispin /= spinors) then
441 if (states_are_real(this%st)) then
442 ! 1 add (1), 2 mult (1)
443 call profiling_count_operations(psib%nst*this%gr%np*(1 + 2))
444 else
445 ! 1 add (2), 1 complex mult (6), 1 real-complex mult (2)
446 call profiling_count_operations(psib%nst*this%gr%np*(2 + 6 + 2))
447 end if
448 else
449 ! 4 add (2), 3 complex mult (6), 3 real-complex mult (2)
450 call profiling_count_operations(psib%nst*this%gr%np*(4*2 + 3*6 + 3*2))
451 end if
452
453 call profiling_out("CALC_DENSITY")
454
456 end subroutine density_calc_accumulate
457
458 ! ---------------------------------------------------
463 !
464 subroutine density_calc_end(this, allreduce, symmetrize)
465 type(density_calc_t), intent(inout) :: this
466 logical, optional, intent(in) :: allreduce
467 logical, optional, intent(in) :: symmetrize
468
469 integer :: ispin, ip
470 real(real64), allocatable :: tmpdensity(:)
471
472 push_sub(density_calc_end)
473
474 if (this%packed) then
475 safe_allocate(tmpdensity(1:this%gr%np))
476
477 ! the density is in device memory
478 do ispin = 1, this%st%d%nspin
479 call accel_read_buffer(this%buff_density, int(this%gr%np, int64), tmpdensity, offset = (ispin - 1)*this%pnp)
480
481 !$omp parallel do simd
482 do ip = 1, this%gr%np
483 this%density(ip, ispin) = this%density(ip, ispin) + tmpdensity(ip)
484 end do
485 end do
486
487 this%packed = .false.
488 call accel_release_buffer(this%buff_density)
489 safe_deallocate_a(tmpdensity)
490 end if
491
492 ! reduce over states and k-points
493 if ((this%st%parallel_in_states .or. this%st%d%kpt%parallel) .and. optional_default(allreduce, .true.)) then
494 call profiling_in("DENSITY_REDUCE")
495 call comm_allreduce(this%st%st_kpt_mpi_grp, this%density, dim = (/this%gr%np, this%st%d%nspin/))
496 call profiling_out("DENSITY_REDUCE")
497 end if
498
499 if (this%st%symmetrize_density .and. optional_default(symmetrize, .true.)) then
500 do ispin = 1, this%st%d%nspin
501 call dgrid_symmetrize_scalar_field(this%gr, this%density(:, ispin))
502 end do
503 end if
504
505 pop_sub(density_calc_end)
506 end subroutine density_calc_end
507
508
509 ! ---------------------------------------------------------
514 !
515 subroutine density_calc(st, gr, density, istin)
516 type(states_elec_t), intent(in) :: st
517 type(grid_t), intent(in) :: gr
518 real(real64), intent(out) :: density(:, :)
519 integer, optional, intent(in) :: istin
520
521 integer :: ik, ib
522 type(density_calc_t) :: dens_calc
523
524 push_sub(density_calc)
525
526 assert(ubound(density, dim = 1) == gr%np .or. ubound(density, dim = 1) == gr%np_part)
527
528 call density_calc_init(dens_calc, st, gr, density)
529
530 if (.not. present(istin)) then
531 ! ordinary density accumulate
532 do ik = st%d%kpt%start, st%d%kpt%end
533 do ib = st%group%block_start, st%group%block_end
534 call density_calc_accumulate(dens_calc, st%group%psib(ib, ik), async=.true.)
535 end do
536 end do
537 call accel_finish()
538 else
539 ! state-resolved density
540 do ik = st%d%kpt%start, st%d%kpt%end
541 do ib = st%group%block_start, st%group%block_end
542 if (st%group%psib(ib, ik)%ist(1) <= istin .and. st%group%psib(ib, ik)%ist(st%group%psib(ib, ik)%nst) >= istin) then
543 call density_calc_state(dens_calc, st%group%psib(ib, ik), istin)
544 end if
545 end do
546 end do
547
548 end if
549
550
551 call density_calc_end(dens_calc, allreduce = .not. present(istin))
552
553 pop_sub(density_calc)
554 end subroutine density_calc
555
556 ! ---------------------------------------------------------
558 subroutine states_elec_freeze_orbitals(st, namespace, space, gr, mc, kpoints, n, family_is_mgga)
559 type(states_elec_t), intent(inout) :: st
560 type(namespace_t), intent(in) :: namespace
561 class(space_t), intent(in) :: space
562 type(grid_t), intent(in) :: gr
563 type(multicomm_t), intent(in) :: mc
564 type(kpoints_t), intent(in) :: kpoints
565 integer, intent(in) :: n
566 logical, intent(in) :: family_is_mgga
567
568 integer :: ist, ik, ib, nblock
569 integer :: nodeto, nodefr
570 type(states_elec_t) :: staux
571 complex(real64), allocatable :: psi(:, :, :), rec_buffer(:,:)
572 type(wfs_elec_t) :: psib
573 type(density_calc_t) :: dens_calc
574
576
577 if (n >= st%nst) then
578 write(message(1),'(a)') 'Attempting to freeze a number of orbitals which is larger or equal to'
579 write(message(2),'(a)') 'the total number. The program has to stop.'
580 call messages_fatal(2, namespace=namespace)
581 end if
582
583 assert(states_are_complex(st))
584
585 if (.not. allocated(st%frozen_rho)) then
586 safe_allocate(st%frozen_rho(1:gr%np, 1:st%d%nspin))
587 end if
588
589 call density_calc_init(dens_calc, st, gr, st%frozen_rho)
590
591 do ik = st%d%kpt%start, st%d%kpt%end
592 if (n < st%st_start) cycle
593
594 do ib = st%group%block_start, st%group%block_end
595 !We can use the full batch
596 if (states_elec_block_max(st, ib) <= n) then
597
598 call density_calc_accumulate(dens_calc, st%group%psib(ib, ik))
599 if (states_elec_block_max(st, ib) == n) exit
600
601 else !Here we only use a part of this batch
602
603 nblock = n - states_elec_block_min(st, ib) + 1
604
605 safe_allocate(psi(1:gr%np, 1:st%d%dim, 1:nblock))
606
607 do ist = 1, nblock
608 call states_elec_get_state(st, gr, states_elec_block_min(st, ib) + ist - 1, ik, psi(:, :, ist))
609 end do
610
611 call wfs_elec_init(psib, st%d%dim, states_elec_block_min(st, ib), n, psi, ik)
612 call density_calc_accumulate(dens_calc, psib)
613 call psib%end()
614 safe_deallocate_a(psi)
615
616 exit
617
618 end if
619
620 end do
621 end do
622
623 call density_calc_end(dens_calc)
624
625 if (family_is_mgga) then
626 if (.not. allocated(st%frozen_tau)) then
627 safe_allocate(st%frozen_tau(1:gr%np, 1:st%d%nspin))
628 end if
629 if (.not. allocated(st%frozen_gdens)) then
630 safe_allocate(st%frozen_gdens(1:gr%np, 1:space%dim, 1:st%d%nspin))
631 end if
632 if (.not. allocated(st%frozen_ldens)) then
633 safe_allocate(st%frozen_ldens(1:gr%np, 1:st%d%nspin))
634 end if
635
636 call states_elec_calc_quantities(gr, st, kpoints, .true., kinetic_energy_density = st%frozen_tau, &
637 density_gradient = st%frozen_gdens, density_laplacian = st%frozen_ldens, st_end = n)
638 end if
639
640
641 call states_elec_copy(staux, st)
642
643 call states_elec_freeze_redistribute_states(st, namespace, gr, mc, n)
644
645 safe_allocate(psi(1:gr%np, 1:st%d%dim, 1))
646 safe_allocate(rec_buffer(1:gr%np, 1:st%d%dim))
647
648 do ik = st%d%kpt%start, st%d%kpt%end
649 ! loop over all frozen states
650 do ist = 1, st%nst
651 ! new distribution: frozen states
652 nodeto = st%node(ist)
653 ! old distribution: full states
654 nodefr = staux%node(ist+n)
655
656 if (st%mpi_grp%rank == nodeto .and. st%mpi_grp%rank == nodefr) then
657 ! do a simple copy on the same rank, also for serial version
658 call states_elec_get_state(staux, gr, ist+n, ik, psi(:, :, 1))
659 call states_elec_set_state(st, gr, ist, ik, psi(:, :, 1))
660 else
661 ! we can use blocking send/recv calls here because we always have different ranks
662 if (st%mpi_grp%rank == nodefr) then
663 call states_elec_get_state(staux, gr, ist+n, ik, psi(:, :, 1))
664 call st%mpi_grp%send(psi(1, 1, 1), gr%np*st%d%dim, mpi_double_complex, nodeto, ist)
665 end if
666 if (st%mpi_grp%rank == nodeto) then
667 call st%mpi_grp%recv(rec_buffer(1, 1), gr%np*st%d%dim, mpi_double_complex, nodefr, ist)
668 call states_elec_set_state(st, gr, ist, ik, rec_buffer(:, :))
669 end if
670 end if
671 end do
672 end do
673
674 safe_deallocate_a(psi)
675 safe_deallocate_a(rec_buffer)
676
677 do ik = 1, st%nik
678 do ist = 1, st%nst
679 st%occ(ist, ik) = staux%occ(n+ist, ik)
680 st%eigenval(ist, ik) = staux%eigenval(n+ist, ik)
681 end do
682 end do
683
684 call states_elec_end(staux)
686 end subroutine states_elec_freeze_orbitals
687
688 ! ---------------------------------------------------------
689 subroutine states_elec_freeze_redistribute_states(st, namespace, mesh, mc, nn)
690 type(states_elec_t), intent(inout) :: st
691 type(namespace_t), intent(in) :: namespace
692 class(mesh_t), intent(in) :: mesh
693 type(multicomm_t), intent(in) :: mc
694 integer, intent(in) :: nn
695
697
698 st%nst = st%nst - nn
699
701 call states_elec_distribute_nodes(st, namespace, mc)
702 call states_elec_allocate_wfns(st, mesh, type_cmplx, packed=.true.)
703
704 safe_deallocate_a(st%eigenval)
705 safe_allocate(st%eigenval(1:st%nst, 1:st%nik))
706 st%eigenval = huge(st%eigenval)
707
708 safe_deallocate_a(st%occ)
709 safe_allocate(st%occ (1:st%nst, 1:st%nik))
710 st%occ = m_zero
711
712
715
716 ! ---------------------------------------------------------
717 subroutine states_elec_freeze_adjust_qtot(st)
718 type(states_elec_t), intent(inout) :: st
719
720 integer :: ik, ist
721
723
724 ! Change the smearing method by fixing the occupations to
725 ! that of the ground-state such that the unfrozen states inherit
726 ! those values.
727 st%smear%method = smear_fixed_occ
728
729 ! Set total charge
730 st%qtot = m_zero
731 do ik = st%d%kpt%start, st%d%kpt%end
732 do ist = st%st_start, st%st_end
733 st%qtot = st%qtot + st%occ(ist, ik) * st%kweights(ik)
734 end do
735 end do
736
737 if (st%parallel_in_states .or. st%d%kpt%parallel) then
738 call comm_allreduce(st%st_kpt_mpi_grp, st%qtot)
739 end if
740
741
743 end subroutine states_elec_freeze_adjust_qtot
744
745
746 ! ---------------------------------------------------------
755 !
756 subroutine states_elec_total_density(st, mesh, total_rho)
757 type(states_elec_t), intent(in) :: st
758 class(mesh_t), intent(in) :: mesh
759 real(real64), contiguous, intent(out) :: total_rho(:,:)
760
761 integer :: is
762
764
765 call lalg_copy(mesh%np, st%d%nspin, st%rho, total_rho)
766
767 if (allocated(st%rho_core)) then
768 do is = 1, st%d%spin_channels
769 call lalg_axpy(mesh%np, m_one/st%d%spin_channels, st%rho_core, total_rho(:, is))
770 end do
771 end if
772
773 ! Add, if it exists, the frozen density from the inner orbitals.
774 if (allocated(st%frozen_rho)) then
775 call lalg_axpy(mesh%np, st%d%nspin, m_one, st%frozen_rho, total_rho)
776 end if
777
779 end subroutine states_elec_total_density
780
781#include "undef.F90"
782#include "real.F90"
783#include "density_inc.F90"
784
785#include "undef.F90"
786#include "complex.F90"
787#include "density_inc.F90"
788
789end module density_oct_m
790
791
792!! Local Variables:
793!! mode: f90
794!! coding: utf-8
795!! 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
type(accel_kernel_t), target, save, public kernel_density_real
Definition: accel.F90:285
subroutine, public accel_finish()
Definition: accel.F90:1296
type(accel_kernel_t), target, save, public kernel_density_spinors
Definition: accel.F90:287
integer, parameter, public accel_mem_read_write
Definition: accel.F90:183
subroutine, public accel_release_buffer(this)
Definition: accel.F90:1246
type(accel_kernel_t), target, save, public kernel_density_complex
Definition: accel.F90:286
integer function, public accel_kernel_workgroup_size(kernel)
Definition: accel.F90:1478
integer, parameter, public accel_mem_read_only
Definition: accel.F90:183
This module implements batches of mesh functions.
Definition: batch.F90:133
integer, parameter, public batch_not_packed
functions are stored in CPU memory, unpacked order
Definition: batch.F90:276
integer, parameter, public batch_device_packed
functions are stored in device memory in packed order
Definition: batch.F90:276
integer, parameter, public batch_packed
functions are stored in CPU memory, in transposed (packed) order
Definition: batch.F90:276
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:116
This module implements a calculator for the density and defines related functions.
Definition: density.F90:120
subroutine, public density_calc_end(this, allreduce, symmetrize)
Finalize the density calculation.
Definition: density.F90:558
subroutine density_calc_state(this, psib, istin)
Calculate contribution to the density from state istin.
Definition: density.F90:245
subroutine, public density_calc_accumulate(this, psib, async)
Accumulate weighted orbital densities for a batch psib.
Definition: density.F90:391
subroutine, public states_elec_freeze_adjust_qtot(st)
Definition: density.F90:811
subroutine, public states_elec_freeze_redistribute_states(st, namespace, mesh, mc, nn)
Definition: density.F90:783
subroutine, public density_calc_init(this, st, gr, density)
initialize the density calculator
Definition: density.F90:187
subroutine density_calc_pack(this, async)
prepare the density calculator for GPU use
Definition: density.F90:224
subroutine, public ddensity_accumulate_grad(space, mesh, st, psib, grad_psib, grad_rho)
Accumulate within a batch.
Definition: density.F90:945
subroutine, public zdensity_accumulate_grad(space, mesh, st, psib, grad_psib, grad_rho)
Accumulate within a batch.
Definition: density.F90:1131
subroutine, public states_elec_total_density(st, mesh, total_rho)
This routine calculates the total electronic density.
Definition: density.F90:850
subroutine, public states_elec_freeze_orbitals(st, namespace, space, gr, mc, kpoints, n, family_is_mgga)
Calculate partial density for frozen orbitals.
Definition: density.F90:652
subroutine, public density_calc(st, gr, density, istin)
Computes the density from the orbitals in st.
Definition: density.F90:609
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
integer, parameter, public unpolarized
Parameters...
integer, parameter, public spinors
integer, parameter, public spin_polarized
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public m_epsilon
Definition: global.F90:204
real(real64), parameter, public m_one
Definition: global.F90:189
This module implements the underlying real-space grid.
Definition: grid.F90:117
subroutine, public dgrid_symmetrize_scalar_field(gr, field, suppress_warning)
Definition: grid.F90:661
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_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:414
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:145
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:623
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:552
integer, parameter, public smear_fixed_occ
Definition: smear.F90:171
pure logical function, public states_are_complex(st)
pure logical function, public states_are_real(st)
This module handles spin dimensions of the states and the k-point distribution.
integer pure function, public states_elec_block_max(st, ib)
return index of last state in block ib
subroutine, public states_elec_distribute_nodes(st, namespace, mc)
@Brief. Distribute states over the processes for states parallelization
subroutine, public states_elec_calc_quantities(gr, st, kpoints, nlcc, kinetic_energy_density, paramagnetic_current, density_gradient, density_laplacian, gi_kinetic_energy_density, st_end)
calculated selected quantities
subroutine, public states_elec_end(st)
finalize the states_elec_t object
subroutine, public states_elec_deallocate_wfns(st)
Deallocates the KS wavefunctions defined within a states_elec_t structure.
subroutine, public states_elec_allocate_wfns(st, mesh, wfs_type, skip, packed)
Allocates the KS wavefunctions defined within a states_elec_t structure.
subroutine, public states_elec_copy(stout, stin, exclude_wfns, exclude_eigenval, special)
make a (selective) copy of a states_elec_t object
integer pure function, public states_elec_block_min(st, ib)
return index of first state in block ib
type(type_t), public type_float
Definition: types.F90:133
type(type_t), public type_cmplx
Definition: types.F90:134
The calculator object.
Definition: density.F90:169
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:168
Describes mesh distribution to nodes.
Definition: mesh.F90:186
Stores all communicators and groups.
Definition: multicomm.F90:206
The states_elec_t class contains all electronic wave functions.
batches of electronic states
Definition: wfs_elec.F90:138
int true(void)