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), contiguous, 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), contiguous, 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 assert(istin_ > 0)
177 if (abs(weight(istin_)) <= m_epsilon) then
178 call profiling_out("CALC_STATE_DENSITY")
179 pop_sub(density_calc_state)
180 return
181 end if
182
183 !fix ist for the remaining code
184 ist = istin_
185
187 select case (psib%status())
188 case (batch_not_packed)
189 select case (this%st%d%ispin)
191 if (states_are_real(this%st)) then
192 !$omp parallel do simd schedule(static)
193 do ip = 1, this%gr%np
194 this%density(ip, ispin) = this%density(ip, ispin) + weight(ist)*psib%dff(ip, 1, ist)**2
195 end do
196 else
197 !$omp parallel do schedule(static)
198 do ip = 1, this%gr%np
199 this%density(ip, ispin) = this%density(ip, ispin) + weight(ist)* &
200 real(conjg(psib%zff(ip, 1, ist))*psib%zff(ip, 1, ist), real64)
201 end do
202 end if
203 case (spinors)
204 !$omp parallel do schedule(static) private(psi1, psi2, term)
205 do ip = 1, this%gr%np
206 psi1 = psib%zff(ip, 1, ist)
207 psi2 = psib%zff(ip, 2, ist)
208 this%density(ip, 1) = this%density(ip, 1) + weight(ist)*real(conjg(psi1)*psi1, real64)
209 this%density(ip, 2) = this%density(ip, 2) + weight(ist)*real(conjg(psi2)*psi2, real64)
210 term = weight(ist)*psi1*conjg(psi2)
211 this%density(ip, 3) = this%density(ip, 3) + real(term, real64)
212 this%density(ip, 4) = this%density(ip, 4) + aimag(term)
213 end do
214 end select
215
216 case (batch_packed)
217
218 select case (this%st%d%ispin)
220 if (states_are_real(this%st)) then
221 !$omp parallel do schedule(static)
222 do ip = 1, this%gr%np
223 this%density(ip, ispin) = this%density(ip, ispin) + weight(ist)*psib%dff_pack(ist, ip)**2
224 end do
225 else
226
227 !$omp parallel do schedule(static)
228 do ip = 1, this%gr%np
229 this%density(ip, ispin) = this%density(ip, ispin) + weight(ist)* &
230 real(conjg(psib%zff_pack(ist, ip))*psib%zff_pack(ist, ip), real64)
231 end do
232 end if
233 case (spinors)
234 assert(mod(psib%nst_linear, 2) == 0)
235 !$omp parallel do schedule(static) private(psi1, psi2, term)
236 do ip = 1, this%gr%np
237 psi1 = psib%zff_pack(2*ist - 1, ip)
238 psi2 = psib%zff_pack(2*ist, ip)
239 term = weight(ist)*psi1*conjg(psi2)
240
241 this%density(ip, 1) = this%density(ip, 1) + weight(ist)*real(conjg(psi1)*psi1, real64)
242 this%density(ip, 2) = this%density(ip, 2) + weight(ist)*real(conjg(psi2)*psi2, real64)
243 this%density(ip, 3) = this%density(ip, 3) + real(term, real64)
244 this%density(ip, 4) = this%density(ip, 4) + aimag(term)
245 end do
246 end select
247
249 if (.not. this%packed) call density_calc_pack(this, async=.true.)
250
251 call accel_create_buffer(buff_weight, accel_mem_read_only, type_float, psib%nst)
252 call accel_write_buffer(buff_weight, psib%nst, weight, async=.true.)
253
254 select case (this%st%d%ispin)
256 if (states_are_real(this%st)) then
257 kernel => kernel_density_real
258 else
259 kernel => kernel_density_complex
260 end if
261
262 call accel_set_kernel_arg(kernel, 0, psib%nst)
263 call accel_set_kernel_arg(kernel, 1, this%gr%np)
264 call accel_set_kernel_arg(kernel, 2, int(this%pnp*(ispin - 1), int32))
265 call accel_set_kernel_arg(kernel, 3, buff_weight)
266 call accel_set_kernel_arg(kernel, 4, psib%ff_device)
267 call accel_set_kernel_arg(kernel, 5, int(log2(psib%pack_size(1)), int32))
268 call accel_set_kernel_arg(kernel, 6, this%buff_density)
269
270 case (spinors)
271 kernel => kernel_density_spinors
272
273 call accel_set_kernel_arg(kernel, 0, psib%nst)
274 call accel_set_kernel_arg(kernel, 1, this%gr%np)
275 call accel_set_kernel_arg(kernel, 2, int(this%pnp, int32))
276 call accel_set_kernel_arg(kernel, 3, buff_weight)
277 call accel_set_kernel_arg(kernel, 4, psib%ff_device)
278 call accel_set_kernel_arg(kernel, 5, int(log2(psib%pack_size(1)), int32))
279 call accel_set_kernel_arg(kernel, 6, this%buff_density)
280 end select
281
282 wgsize = accel_kernel_workgroup_size(kernel)
283
284 call accel_kernel_run(kernel, (/pad(this%gr%np, wgsize)/), (/wgsize/))
285
286 call accel_finish()
287 call accel_release_buffer(buff_weight)
288 end select
289
290 safe_deallocate_a(weight)
291
292 call profiling_out("CALC_STATE_DENSITY")
293
294 pop_sub(density_calc_state)
295 end subroutine density_calc_state
296
297 ! ---------------------------------------------------
299 !
300 subroutine density_calc_accumulate(this, psib, async)
301 type(density_calc_t), intent(inout) :: this
302 type(wfs_elec_t), intent(in) :: psib
303 logical, optional, intent(in) :: async
304
305 integer :: ist, ip, ispin, wgsize
306 complex(real64) :: term, psi1, psi2
307 real(real64), allocatable :: weight(:)
308 type(accel_mem_t) :: buff_weight
309 type(accel_kernel_t), pointer :: kernel
310
312 call profiling_in("CALC_DENSITY")
313
314 ispin = this%st%d%get_spin_index(psib%ik)
315
316 safe_allocate(weight(1:psib%nst))
317 do ist = 1, psib%nst
318 weight(ist) = this%st%kweights(psib%ik)*this%st%occ(psib%ist(ist), psib%ik)
319 end do
320
321 select case (psib%status())
322 case (batch_not_packed)
323 select case (this%st%d%ispin)
325 if (states_are_real(this%st)) then
326 do ist = 1, psib%nst
327 if (abs(weight(ist)) <= m_epsilon) cycle
328 !$omp parallel do simd schedule(static)
329 do ip = 1, this%gr%np
330 this%density(ip, ispin) = this%density(ip, ispin) + weight(ist)*psib%dff(ip, 1, ist)**2
331 end do
332 end do
333 else
334 do ist = 1, psib%nst
335 if (abs(weight(ist)) <= m_epsilon) cycle
336 !$omp parallel do schedule(static)
337 do ip = 1, this%gr%np
338 this%density(ip, ispin) = this%density(ip, ispin) + weight(ist)* &
339 real(conjg(psib%zff(ip, 1, ist))*psib%zff(ip, 1, ist), real64)
340 end do
341 end do
342 end if
343 case (spinors)
344 do ist = 1, psib%nst
345 if (abs(weight(ist)) <= m_epsilon) cycle
346 !$omp parallel do schedule(static) private(psi1, psi2, term)
347 do ip = 1, this%gr%np
348 psi1 = psib%zff(ip, 1, ist)
349 psi2 = psib%zff(ip, 2, ist)
350 this%density(ip, 1) = this%density(ip, 1) + weight(ist)*real(conjg(psi1)*psi1, real64)
351 this%density(ip, 2) = this%density(ip, 2) + weight(ist)*real(conjg(psi2)*psi2, real64)
352 term = weight(ist)*psi1*conjg(psi2)
353 this%density(ip, 3) = this%density(ip, 3) + real(term, real64)
354 this%density(ip, 4) = this%density(ip, 4) + aimag(term)
355 end do
356 end do
357 end select
358
359 case (batch_packed)
360
361 select case (this%st%d%ispin)
363 if (states_are_real(this%st)) then
364 !$omp parallel do schedule(static) private(ist)
365 do ip = 1, this%gr%np
366 do ist = 1, psib%nst
367 this%density(ip, ispin) = this%density(ip, ispin) + weight(ist)*psib%dff_pack(ist, ip)**2
368 end do
369 end do
370 else
371 !$omp parallel do schedule(static) private(ist)
372 do ip = 1, this%gr%np
373 do ist = 1, psib%nst
374 this%density(ip, ispin) = this%density(ip, ispin) + weight(ist)* &
375 real(conjg(psib%zff_pack(ist, ip))*psib%zff_pack(ist, ip), real64)
376 end do
377 end do
378 end if
379 case (spinors)
380 assert(mod(psib%nst_linear, 2) == 0)
381 !$omp parallel do schedule(static) private(ist, psi1, psi2, term)
382 do ip = 1, this%gr%np
383 do ist = 1, psib%nst
384 psi1 = psib%zff_pack(2*ist - 1, ip)
385 psi2 = psib%zff_pack(2*ist, ip)
386 term = weight(ist)*psi1*conjg(psi2)
387
388 this%density(ip, 1) = this%density(ip, 1) + weight(ist)*real(conjg(psi1)*psi1, real64)
389 this%density(ip, 2) = this%density(ip, 2) + weight(ist)*real(conjg(psi2)*psi2, real64)
390 this%density(ip, 3) = this%density(ip, 3) + real(term, real64)
391 this%density(ip, 4) = this%density(ip, 4) + aimag(term)
392 end do
393 end do
394 end select
395
397 if (.not. this%packed) call density_calc_pack(this)
398
399 call accel_create_buffer(buff_weight, accel_mem_read_only, type_float, psib%nst)
400 call accel_write_buffer(buff_weight, psib%nst, weight, async=.true.)
401
402 select case (this%st%d%ispin)
404 if (states_are_real(this%st)) then
405 kernel => kernel_density_real
406 else
407 kernel => kernel_density_complex
408 end if
409
410 call accel_set_kernel_arg(kernel, 0, psib%nst)
411 call accel_set_kernel_arg(kernel, 1, this%gr%np)
412 call accel_set_kernel_arg(kernel, 2, int(this%pnp, int32)*(ispin - 1))
413 call accel_set_kernel_arg(kernel, 3, buff_weight)
414 call accel_set_kernel_arg(kernel, 4, psib%ff_device)
415 call accel_set_kernel_arg(kernel, 5, int(log2(psib%pack_size(1)), int32))
416 call accel_set_kernel_arg(kernel, 6, this%buff_density)
417
418 case (spinors)
419 kernel => kernel_density_spinors
420
421 call accel_set_kernel_arg(kernel, 0, psib%nst)
422 call accel_set_kernel_arg(kernel, 1, this%gr%np)
423 call accel_set_kernel_arg(kernel, 2, int(this%pnp, int32))
424 call accel_set_kernel_arg(kernel, 3, buff_weight)
425 call accel_set_kernel_arg(kernel, 4, psib%ff_device)
426 call accel_set_kernel_arg(kernel, 5, int(log2(psib%pack_size(1)), int32))
427 call accel_set_kernel_arg(kernel, 6, this%buff_density)
428 end select
429
430 wgsize = accel_kernel_workgroup_size(kernel)
431
432 call accel_kernel_run(kernel, (/pad(this%gr%np, wgsize)/), (/wgsize/))
433
434 call accel_release_buffer(buff_weight)
435
436 if(.not. optional_default(async, .false.)) call accel_finish()
437
438 end select
439
440 safe_deallocate_a(weight)
441
442
443 if (this%st%d%ispin /= spinors) then
444 if (states_are_real(this%st)) then
445 ! 1 add (1), 2 mult (1)
446 call profiling_count_operations(psib%nst*this%gr%np*(1 + 2))
447 else
448 ! 1 add (2), 1 complex mult (6), 1 real-complex mult (2)
449 call profiling_count_operations(psib%nst*this%gr%np*(2 + 6 + 2))
450 end if
451 else
452 ! 4 add (2), 3 complex mult (6), 3 real-complex mult (2)
453 call profiling_count_operations(psib%nst*this%gr%np*(4*2 + 3*6 + 3*2))
454 end if
455
456 call profiling_out("CALC_DENSITY")
457
459 end subroutine density_calc_accumulate
460
461 ! ---------------------------------------------------
466 !
467 subroutine density_calc_end(this, allreduce, symmetrize)
468 type(density_calc_t), intent(inout) :: this
469 logical, optional, intent(in) :: allreduce
470 logical, optional, intent(in) :: symmetrize
471
472 integer :: ispin, ip
473 real(real64), allocatable :: tmpdensity(:)
474
475 push_sub(density_calc_end)
476
477 if (this%packed) then
478 safe_allocate(tmpdensity(1:this%gr%np))
479
480 ! the density is in device memory
481 do ispin = 1, this%st%d%nspin
482 call accel_read_buffer(this%buff_density, int(this%gr%np, int64), tmpdensity, offset = (ispin - 1)*this%pnp)
483
484 !$omp parallel do simd
485 do ip = 1, this%gr%np
486 this%density(ip, ispin) = this%density(ip, ispin) + tmpdensity(ip)
487 end do
488 end do
489
490 this%packed = .false.
491 call accel_release_buffer(this%buff_density)
492 safe_deallocate_a(tmpdensity)
493 end if
494
495 ! reduce over states and k-points
496 if ((this%st%parallel_in_states .or. this%st%d%kpt%parallel) .and. optional_default(allreduce, .true.)) then
497 call profiling_in("DENSITY_REDUCE")
498 call comm_allreduce(this%st%st_kpt_mpi_grp, this%density, dim = (/this%gr%np, this%st%d%nspin/))
499 call profiling_out("DENSITY_REDUCE")
500 end if
501
502 if (this%st%symmetrize_density .and. optional_default(symmetrize, .true.)) then
503 do ispin = 1, this%st%d%nspin
504 call dgrid_symmetrize_scalar_field(this%gr, this%density(:, ispin))
505 end do
506 end if
507
508 pop_sub(density_calc_end)
509 end subroutine density_calc_end
510
511
512 ! ---------------------------------------------------------
517 !
518 subroutine density_calc(st, gr, density, istin)
519 type(states_elec_t), intent(in) :: st
520 type(grid_t), intent(in) :: gr
521 real(real64), contiguous,intent(out) :: density(:, :)
522 integer, optional, intent(in) :: istin
523
524 integer :: ik, ib
525 type(density_calc_t) :: dens_calc
526
527 push_sub(density_calc)
528
529 assert(ubound(density, dim = 1) == gr%np .or. ubound(density, dim = 1) == gr%np_part)
530
531 call density_calc_init(dens_calc, st, gr, density)
532
533 if (.not. present(istin)) then
534 ! ordinary density accumulate
535 do ik = st%d%kpt%start, st%d%kpt%end
536 do ib = st%group%block_start, st%group%block_end
537 call density_calc_accumulate(dens_calc, st%group%psib(ib, ik), async=.true.)
538 end do
539 end do
540 call accel_finish()
541 else
542 ! state-resolved density
543 do ik = st%d%kpt%start, st%d%kpt%end
544 do ib = st%group%block_start, st%group%block_end
545 if (st%group%psib(ib, ik)%ist(1) <= istin .and. st%group%psib(ib, ik)%ist(st%group%psib(ib, ik)%nst) >= istin) then
546 call density_calc_state(dens_calc, st%group%psib(ib, ik), istin)
547 end if
548 end do
549 end do
550
551 end if
552
553
554 call density_calc_end(dens_calc, allreduce = .not. present(istin))
555
556 pop_sub(density_calc)
557 end subroutine density_calc
558
559 ! ---------------------------------------------------------
561 subroutine states_elec_freeze_orbitals(st, namespace, space, gr, mc, kpoints, n, family_is_mgga)
562 type(states_elec_t), intent(inout) :: st
563 type(namespace_t), intent(in) :: namespace
564 class(space_t), intent(in) :: space
565 type(grid_t), intent(in) :: gr
566 type(multicomm_t), intent(in) :: mc
567 type(kpoints_t), intent(in) :: kpoints
568 integer, intent(in) :: n
569 logical, intent(in) :: family_is_mgga
570
571 integer :: ist, ik, ib, nblock
572 integer :: nodeto, nodefr
573 type(states_elec_t) :: staux
574 complex(real64), allocatable :: psi(:, :, :), rec_buffer(:,:)
575 type(wfs_elec_t) :: psib
576 type(density_calc_t) :: dens_calc
577
579
580 if (n >= st%nst) then
581 write(message(1),'(a)') 'Attempting to freeze a number of orbitals which is larger or equal to'
582 write(message(2),'(a)') 'the total number. The program has to stop.'
583 call messages_fatal(2, namespace=namespace)
584 end if
585
586 assert(states_are_complex(st))
587
588 if (.not. allocated(st%frozen_rho)) then
589 safe_allocate(st%frozen_rho(1:gr%np, 1:st%d%nspin))
590 end if
591
592 call density_calc_init(dens_calc, st, gr, st%frozen_rho)
593
594 do ik = st%d%kpt%start, st%d%kpt%end
595 if (n < st%st_start) cycle
596
597 do ib = st%group%block_start, st%group%block_end
598 !We can use the full batch
599 if (states_elec_block_max(st, ib) <= n) then
600
601 call density_calc_accumulate(dens_calc, st%group%psib(ib, ik))
602 if (states_elec_block_max(st, ib) == n) exit
603
604 else !Here we only use a part of this batch
605
606 nblock = n - states_elec_block_min(st, ib) + 1
607
608 safe_allocate(psi(1:gr%np, 1:st%d%dim, 1:nblock))
609
610 do ist = 1, nblock
611 call states_elec_get_state(st, gr, states_elec_block_min(st, ib) + ist - 1, ik, psi(:, :, ist))
612 end do
613
614 call wfs_elec_init(psib, st%d%dim, states_elec_block_min(st, ib), n, psi, ik)
615 call density_calc_accumulate(dens_calc, psib)
616 call psib%end()
617 safe_deallocate_a(psi)
618
619 exit
620
621 end if
622
623 end do
624 end do
625
626 call density_calc_end(dens_calc)
627
628 if (family_is_mgga) then
629 if (.not. allocated(st%frozen_tau)) then
630 safe_allocate(st%frozen_tau(1:gr%np, 1:st%d%nspin))
631 end if
632 if (.not. allocated(st%frozen_gdens)) then
633 safe_allocate(st%frozen_gdens(1:gr%np, 1:space%dim, 1:st%d%nspin))
634 end if
635 if (.not. allocated(st%frozen_ldens)) then
636 safe_allocate(st%frozen_ldens(1:gr%np, 1:st%d%nspin))
637 end if
638
639 call states_elec_calc_quantities(gr, st, kpoints, .true., kinetic_energy_density = st%frozen_tau, &
640 density_gradient = st%frozen_gdens, density_laplacian = st%frozen_ldens, st_end = n)
641 end if
642
643
644 call states_elec_copy(staux, st)
645
646 call states_elec_freeze_redistribute_states(st, namespace, gr, mc, n)
647
648 safe_allocate(psi(1:gr%np, 1:st%d%dim, 1))
649 safe_allocate(rec_buffer(1:gr%np, 1:st%d%dim))
650
651 do ik = st%d%kpt%start, st%d%kpt%end
652 ! loop over all frozen states
653 do ist = 1, st%nst
654 ! new distribution: frozen states
655 nodeto = st%node(ist)
656 ! old distribution: full states
657 nodefr = staux%node(ist+n)
658
659 if (st%mpi_grp%rank == nodeto .and. st%mpi_grp%rank == nodefr) then
660 ! do a simple copy on the same rank, also for serial version
661 call states_elec_get_state(staux, gr, ist+n, ik, psi(:, :, 1))
662 call states_elec_set_state(st, gr, ist, ik, psi(:, :, 1))
663 else
664 ! we can use blocking send/recv calls here because we always have different ranks
665 if (st%mpi_grp%rank == nodefr) then
666 call states_elec_get_state(staux, gr, ist+n, ik, psi(:, :, 1))
667 call st%mpi_grp%send(psi(1, 1, 1), gr%np*st%d%dim, mpi_double_complex, nodeto, ist)
668 end if
669 if (st%mpi_grp%rank == nodeto) then
670 call st%mpi_grp%recv(rec_buffer(1, 1), gr%np*st%d%dim, mpi_double_complex, nodefr, ist)
671 call states_elec_set_state(st, gr, ist, ik, rec_buffer(:, :))
672 end if
673 end if
674 end do
675 end do
676
677 safe_deallocate_a(psi)
678 safe_deallocate_a(rec_buffer)
679
680 do ik = 1, st%nik
681 do ist = 1, st%nst
682 st%occ(ist, ik) = staux%occ(n+ist, ik)
683 st%eigenval(ist, ik) = staux%eigenval(n+ist, ik)
684 end do
685 end do
686
687 call states_elec_end(staux)
689 end subroutine states_elec_freeze_orbitals
690
691 ! ---------------------------------------------------------
692 subroutine states_elec_freeze_redistribute_states(st, namespace, mesh, mc, nn)
693 type(states_elec_t), intent(inout) :: st
694 type(namespace_t), intent(in) :: namespace
695 class(mesh_t), intent(in) :: mesh
696 type(multicomm_t), intent(in) :: mc
697 integer, intent(in) :: nn
698
700
701 st%nst = st%nst - nn
702
704 call states_elec_distribute_nodes(st, namespace, mc)
705 call states_elec_allocate_wfns(st, mesh, type_cmplx, packed=.true.)
706
707 safe_deallocate_a(st%eigenval)
708 safe_allocate(st%eigenval(1:st%nst, 1:st%nik))
709 st%eigenval = huge(st%eigenval)
710
711 safe_deallocate_a(st%occ)
712 safe_allocate(st%occ (1:st%nst, 1:st%nik))
713 st%occ = m_zero
714
715
718
719 ! ---------------------------------------------------------
720 subroutine states_elec_freeze_adjust_qtot(st)
721 type(states_elec_t), intent(inout) :: st
722
723 integer :: ik, ist
724
726
727 ! Change the smearing method by fixing the occupations to
728 ! that of the ground-state such that the unfrozen states inherit
729 ! those values.
730 st%smear%method = smear_fixed_occ
731
732 ! Set total charge
733 st%qtot = m_zero
734 do ik = st%d%kpt%start, st%d%kpt%end
735 do ist = st%st_start, st%st_end
736 st%qtot = st%qtot + st%occ(ist, ik) * st%kweights(ik)
737 end do
738 end do
739
740 if (st%parallel_in_states .or. st%d%kpt%parallel) then
741 call comm_allreduce(st%st_kpt_mpi_grp, st%qtot)
742 end if
743
744
746 end subroutine states_elec_freeze_adjust_qtot
747
748
749 ! ---------------------------------------------------------
758 !
759 subroutine states_elec_total_density(st, mesh, total_rho)
760 type(states_elec_t), intent(in) :: st
761 class(mesh_t), intent(in) :: mesh
762 real(real64), contiguous, intent(out) :: total_rho(:,:)
763
764 integer :: is
765
767
768 call lalg_copy(mesh%np, st%d%nspin, st%rho, total_rho)
769
770 if (allocated(st%rho_core)) then
771 do is = 1, st%d%spin_channels
772 call lalg_axpy(mesh%np, m_one/st%d%spin_channels, st%rho_core, total_rho(:, is))
773 end do
774 end if
775
776 ! Add, if it exists, the frozen density from the inner orbitals.
777 if (allocated(st%frozen_rho)) then
778 call lalg_axpy(mesh%np, st%d%nspin, m_one, st%frozen_rho, total_rho)
779 end if
780
782 end subroutine states_elec_total_density
783
784#include "undef.F90"
785#include "real.F90"
786#include "density_inc.F90"
787
788#include "undef.F90"
789#include "complex.F90"
790#include "density_inc.F90"
791
792end module density_oct_m
793
794
795!! Local Variables:
796!! mode: f90
797!! coding: utf-8
798!! 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:1298
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:1248
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:1480
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:561
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:394
subroutine, public states_elec_freeze_adjust_qtot(st)
Definition: density.F90:814
subroutine, public states_elec_freeze_redistribute_states(st, namespace, mesh, mc, nn)
Definition: density.F90:786
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:948
subroutine, public zdensity_accumulate_grad(space, mesh, st, psib, grad_psib, grad_rho)
Accumulate within a batch.
Definition: density.F90:1134
subroutine, public states_elec_total_density(st, mesh, total_rho)
This routine calculates the total electronic density.
Definition: density.F90:853
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:655
subroutine, public density_calc(st, gr, density, istin)
Computes the density from the orbitals in st.
Definition: density.F90:612
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:161
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:415
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), parameter, public type_cmplx
Definition: types.F90:134
type(type_t), parameter, public type_float
Definition: types.F90:133
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)