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 ! !!dimensions (1:np, 1:d%nspin)
80 type(states_elec_t), pointer :: st
81 type(grid_t), pointer :: gr
82 type(accel_mem_t) :: buff_density
83 integer(int64) :: pnp
84 logical :: packed
85 end type density_calc_t
86
87contains
88
93 !
94 subroutine density_calc_init(this, st, gr, density)
95 type(density_calc_t), intent(out) :: this
96 type(states_elec_t), target, intent(in) :: st
97 type(grid_t), target, intent(in) :: gr
98 real(real64), target, intent(out) :: density(:, :)
99
100 integer :: ip, ispin
101 logical :: correct_size
102
103 push_sub(density_calc_init)
104
105 this%st => st
106 this%gr => gr
107
108 this%density => density
109 !$omp parallel private(ip)
110 do ispin = 1, st%d%nspin
111 !$omp do
112 do ip = 1, gr%np
113 this%density(ip, ispin) = m_zero
114 end do
115 !$omp end do nowait
116 end do
117 !$omp end parallel
118
119 this%packed = .false.
120
121 correct_size = ubound(this%density, dim = 1) == this%gr%np .or. &
122 ubound(this%density, dim = 1) == this%gr%np_part
123 assert(correct_size)
124
125 pop_sub(density_calc_init)
126 end subroutine density_calc_init
127
128 ! ---------------------------------------------------
130 !
131 subroutine density_calc_pack(this, async)
132 type(density_calc_t), intent(inout) :: this
133 logical, optional, intent(in) :: async
134
135 push_sub(density_calc_pack)
136
137 this%packed = .true.
138 this%pnp = accel_padded_size(this%gr%np)
139 call accel_create_buffer(this%buff_density, accel_mem_read_write, type_float, this%pnp*this%st%d%nspin)
140
141 ! set to zero
142 call accel_set_buffer_to_zero(this%buff_density, type_float, this%pnp*this%st%d%nspin, async=async)
143
144 pop_sub(density_calc_pack)
145 end subroutine density_calc_pack
146
147 ! ---------------------------------------------------
151 !
152 subroutine density_calc_state(this, psib, istin)
153 type(density_calc_t), intent(inout) :: this
154 type(wfs_elec_t), intent(in) :: psib
155 integer, intent(in) :: istin
156
157 integer :: ist, ip, ispin, istin_, bsize, gsize
158 complex(real64) :: term, psi1, psi2
159 real(real64), allocatable :: weight(:)
160 type(accel_mem_t) :: buff_weight
161 type(accel_kernel_t), pointer :: kernel
162
163 push_sub(density_calc_state)
164 call profiling_in("CALC_STATE_DENSITY")
165
166 ispin = this%st%d%get_spin_index(psib%ik)
167
168 istin_= 0
169
170 safe_allocate(weight(1:psib%nst))
171 do ist = 1, psib%nst
172 weight(ist) = this%st%kweights(psib%ik)*this%st%occ(psib%ist(ist), psib%ik)
173 if (psib%ist(ist) == istin) istin_ = ist
174 end do
176 if (abs(weight(istin_)) <= m_epsilon) then
177 return
179 end if
180
181 !fix ist for the remaining code
182 ist = istin_
183
184
185 select case (psib%status())
186 case (batch_not_packed)
187 select case (this%st%d%ispin)
189 if (states_are_real(this%st)) then
190 !$omp parallel do simd schedule(static)
191 do ip = 1, this%gr%np
192 this%density(ip, ispin) = this%density(ip, ispin) + weight(ist)*psib%dff(ip, 1, ist)**2
193 end do
194 else
195 !$omp parallel do schedule(static)
196 do ip = 1, this%gr%np
197 this%density(ip, ispin) = this%density(ip, ispin) + weight(ist)* &
198 real(conjg(psib%zff(ip, 1, ist))*psib%zff(ip, 1, ist), real64)
199 end do
200 end if
201 case (spinors)
202 !$omp parallel do schedule(static) private(psi1, psi2, term)
203 do ip = 1, this%gr%np
204 psi1 = psib%zff(ip, 1, ist)
205 psi2 = psib%zff(ip, 2, ist)
206 this%density(ip, 1) = this%density(ip, 1) + weight(ist)*real(conjg(psi1)*psi1, real64)
207 this%density(ip, 2) = this%density(ip, 2) + weight(ist)*real(conjg(psi2)*psi2, real64)
208 term = weight(ist)*psi1*conjg(psi2)
209 this%density(ip, 3) = this%density(ip, 3) + real(term, real64)
210 this%density(ip, 4) = this%density(ip, 4) + aimag(term)
211 end do
212 end select
213
214 case (batch_packed)
215
216 select case (this%st%d%ispin)
218 if (states_are_real(this%st)) then
219 !$omp parallel do schedule(static)
220 do ip = 1, this%gr%np
221 this%density(ip, ispin) = this%density(ip, ispin) + weight(ist)*psib%dff_pack(ist, ip)**2
222 end do
223 else
224
225 !$omp parallel do schedule(static)
226 do ip = 1, this%gr%np
227 this%density(ip, ispin) = this%density(ip, ispin) + weight(ist)* &
228 real(conjg(psib%zff_pack(ist, ip))*psib%zff_pack(ist, ip), real64)
229 end do
230 end if
231 case (spinors)
232 assert(mod(psib%nst_linear, 2) == 0)
233 !$omp parallel do schedule(static) private(ist, psi1, psi2, term)
234 do ip = 1, this%gr%np
235 psi1 = psib%zff_pack(2*ist - 1, ip)
236 psi2 = psib%zff_pack(2*ist, ip)
237 term = weight(ist)*psi1*conjg(psi2)
238
239 this%density(ip, 1) = this%density(ip, 1) + weight(ist)*real(conjg(psi1)*psi1, real64)
240 this%density(ip, 2) = this%density(ip, 2) + weight(ist)*real(conjg(psi2)*psi2, real64)
241 this%density(ip, 3) = this%density(ip, 3) + real(term, real64)
242 this%density(ip, 4) = this%density(ip, 4) + aimag(term)
243 end do
244 end select
245
247 if (.not. this%packed) call density_calc_pack(this, async=.true.)
248
249 call accel_create_buffer(buff_weight, accel_mem_read_only, type_float, psib%nst)
250 call accel_write_buffer(buff_weight, psib%nst, weight, async=.true.)
251
252 select case (this%st%d%ispin)
254 if (states_are_real(this%st)) then
255 kernel => kernel_density_real
256 else
257 kernel => kernel_density_complex
258 end if
259
260 call accel_set_kernel_arg(kernel, 0, psib%nst)
261 call accel_set_kernel_arg(kernel, 1, this%gr%np)
262 call accel_set_kernel_arg(kernel, 2, int(this%pnp*(ispin - 1), int32))
263 call accel_set_kernel_arg(kernel, 3, buff_weight)
264 call accel_set_kernel_arg(kernel, 4, psib%ff_device)
265 call accel_set_kernel_arg(kernel, 5, int(log2(psib%pack_size(1)), int32))
266 call accel_set_kernel_arg(kernel, 6, this%buff_density)
267
268 case (spinors)
269 kernel => kernel_density_spinors
270
271 call accel_set_kernel_arg(kernel, 0, psib%nst)
272 call accel_set_kernel_arg(kernel, 1, this%gr%np)
273 call accel_set_kernel_arg(kernel, 2, int(this%pnp, int32))
274 call accel_set_kernel_arg(kernel, 3, buff_weight)
275 call accel_set_kernel_arg(kernel, 4, psib%ff_device)
276 call accel_set_kernel_arg(kernel, 5, int(log2(psib%pack_size(1)), int32))
277 call accel_set_kernel_arg(kernel, 6, this%buff_density)
278 end select
279
280 ! Compute the grid
281 bsize = accel_kernel_block_size(kernel)
282 call accel_grid_size(this%gr%np, bsize, gsize)
283
284 call accel_kernel_run(kernel, (/gsize/), (/bsize/))
285
286 call accel_finish()
287 call accel_free_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, bsize, gsize
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
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 ! Compute the grid
431 bsize = accel_kernel_block_size(kernel)
432 call accel_grid_size(this%gr%np, bsize, gsize)
433
434 call accel_kernel_run(kernel, (/gsize/), (/bsize/))
435
436 call accel_free_buffer(buff_weight)
437
438 if(.not. optional_default(async, .false.)) call accel_finish()
439
440 end select
441
442 safe_deallocate_a(weight)
443
444
445 if (this%st%d%ispin /= spinors) then
446 if (states_are_real(this%st)) then
447 ! 1 add (1), 2 mult (1)
448 call profiling_count_operations(psib%nst*this%gr%np*(1 + 2))
449 else
450 ! 1 add (2), 1 complex mult (6), 1 real-complex mult (2)
451 call profiling_count_operations(psib%nst*this%gr%np*(2 + 6 + 2))
452 end if
453 else
454 ! 4 add (2), 3 complex mult (6), 3 real-complex mult (2)
455 call profiling_count_operations(psib%nst*this%gr%np*(4*2 + 3*6 + 3*2))
456 end if
457
458 call profiling_out("CALC_DENSITY")
459
461 end subroutine density_calc_accumulate
462
463 ! ---------------------------------------------------
468 !
469 subroutine density_calc_end(this, allreduce, symmetrize)
470 type(density_calc_t), intent(inout) :: this
471 logical, optional, intent(in) :: allreduce
472 logical, optional, intent(in) :: symmetrize
473
474 integer :: ispin, ip
475 real(real64), allocatable :: tmpdensity(:)
476
477 push_sub(density_calc_end)
478
479 if (this%packed) then
480 safe_allocate(tmpdensity(1:this%gr%np))
481
482 ! the density is in device memory
483 do ispin = 1, this%st%d%nspin
484 call accel_read_buffer(this%buff_density, int(this%gr%np, int64), tmpdensity, offset = (ispin - 1)*this%pnp)
485
486 !$omp parallel do simd
487 do ip = 1, this%gr%np
488 this%density(ip, ispin) = this%density(ip, ispin) + tmpdensity(ip)
489 end do
490 end do
491
492 this%packed = .false.
493 call accel_free_buffer(this%buff_density)
494 safe_deallocate_a(tmpdensity)
495 end if
496
497 ! reduce over states and k-points
498 if ((this%st%parallel_in_states .or. this%st%d%kpt%parallel) .and. optional_default(allreduce, .true.)) then
499 call profiling_in("DENSITY_REDUCE")
500 call comm_allreduce(this%st%st_kpt_mpi_grp, this%density, dim = (/this%gr%np, this%st%d%nspin/))
501 call profiling_out("DENSITY_REDUCE")
502 end if
503
504 if (this%st%symmetrize_density .and. optional_default(symmetrize, .true.)) then
505 do ispin = 1, this%st%d%nspin
506 call dgrid_symmetrize_scalar_field(this%gr, this%density(:, ispin))
507 end do
508 end if
509
510 pop_sub(density_calc_end)
511 end subroutine density_calc_end
512
513
514 ! ---------------------------------------------------------
519 !
520 subroutine density_calc(st, gr, density, istin)
521 type(states_elec_t), intent(in) :: st
522 type(grid_t), intent(in) :: gr
523 real(real64), intent(out) :: density(:, :)
524 integer, optional, intent(in) :: istin
525
526 integer :: ik, ib
527 type(density_calc_t) :: dens_calc
528
529 push_sub(density_calc)
530
531 assert(ubound(density, dim = 1) == gr%np .or. ubound(density, dim = 1) == gr%np_part)
532
533 call density_calc_init(dens_calc, st, gr, density)
534
535 if (.not. present(istin)) then
536 ! ordinary density accumulate
537 do ik = st%d%kpt%start, st%d%kpt%end
538 do ib = st%group%block_start, st%group%block_end
539 call density_calc_accumulate(dens_calc, st%group%psib(ib, ik), async=.true.)
540 end do
541 end do
542 call accel_finish()
543 else
544 ! state-resolved density
545 do ik = st%d%kpt%start, st%d%kpt%end
546 do ib = st%group%block_start, st%group%block_end
547 if (st%group%psib(ib, ik)%ist(1) <= istin .and. st%group%psib(ib, ik)%ist(st%group%psib(ib, ik)%nst) >= istin) then
548 call density_calc_state(dens_calc, st%group%psib(ib, ik), istin)
549 end if
550 end do
551 end do
552
553 end if
554
555
556 call density_calc_end(dens_calc, allreduce = .not. present(istin))
557
558 pop_sub(density_calc)
559 end subroutine density_calc
560
561 ! ---------------------------------------------------------
563 subroutine states_elec_freeze_orbitals(st, namespace, space, gr, mc, kpoints, n, family_is_mgga)
564 type(states_elec_t), intent(inout) :: st
565 type(namespace_t), intent(in) :: namespace
566 class(space_t), intent(in) :: space
567 type(grid_t), intent(in) :: gr
568 type(multicomm_t), intent(in) :: mc
569 type(kpoints_t), intent(in) :: kpoints
570 integer, intent(in) :: n
571 logical, intent(in) :: family_is_mgga
572
573 integer :: ist, ik, ib, nblock
574 integer :: nodeto, nodefr
575 type(states_elec_t) :: staux
576 complex(real64), allocatable :: psi(:, :, :), rec_buffer(:,:)
577 type(wfs_elec_t) :: psib
578 type(density_calc_t) :: dens_calc
579
581
582 if (n >= st%nst) then
583 write(message(1),'(a)') 'Attempting to freeze a number of orbitals which is larger or equal to'
584 write(message(2),'(a)') 'the total number. The program has to stop.'
585 call messages_fatal(2, namespace=namespace)
586 end if
587
588 assert(states_are_complex(st))
589
590 if (.not. allocated(st%frozen_rho)) then
591 safe_allocate(st%frozen_rho(1:gr%np, 1:st%d%nspin))
592 end if
593
594 call density_calc_init(dens_calc, st, gr, st%frozen_rho)
595
596 do ik = st%d%kpt%start, st%d%kpt%end
597 if (n < st%st_start) cycle
598
599 do ib = st%group%block_start, st%group%block_end
600 !We can use the full batch
601 if (states_elec_block_max(st, ib) <= n) then
602
603 call density_calc_accumulate(dens_calc, st%group%psib(ib, ik))
604 if (states_elec_block_max(st, ib) == n) exit
605
606 else !Here we only use a part of this batch
607
608 nblock = n - states_elec_block_min(st, ib) + 1
609
610 safe_allocate(psi(1:gr%np, 1:st%d%dim, 1:nblock))
611
612 do ist = 1, nblock
613 call states_elec_get_state(st, gr, states_elec_block_min(st, ib) + ist - 1, ik, psi(:, :, ist))
614 end do
616 call wfs_elec_init(psib, st%d%dim, states_elec_block_min(st, ib), n, psi, ik)
617 call density_calc_accumulate(dens_calc, psib)
618 call psib%end()
619 safe_deallocate_a(psi)
620
621 exit
622
623 end if
624
625 end do
626 end do
627
628 call density_calc_end(dens_calc)
629
630 if (family_is_mgga) then
631 if (.not. allocated(st%frozen_tau)) then
632 safe_allocate(st%frozen_tau(1:gr%np, 1:st%d%nspin))
633 end if
634 if (.not. allocated(st%frozen_gdens)) then
635 safe_allocate(st%frozen_gdens(1:gr%np, 1:space%dim, 1:st%d%nspin))
636 end if
637 if (.not. allocated(st%frozen_ldens)) then
638 safe_allocate(st%frozen_ldens(1:gr%np, 1:st%d%nspin))
639 end if
640
641 call states_elec_calc_quantities(gr, st, kpoints, .true., kinetic_energy_density = st%frozen_tau, &
642 density_gradient = st%frozen_gdens, density_laplacian = st%frozen_ldens, st_end = n)
643 end if
644
645
646 call states_elec_copy(staux, st)
647
648 call states_elec_freeze_redistribute_states(st, namespace, gr, mc, n)
649
650 safe_allocate(psi(1:gr%np, 1:st%d%dim, 1))
651 safe_allocate(rec_buffer(1:gr%np, 1:st%d%dim))
652
653 do ik = st%d%kpt%start, st%d%kpt%end
654 ! loop over all frozen states
655 do ist = 1, st%nst
656 ! new distribution: frozen states
657 nodeto = st%node(ist)
658 ! old distribution: full states
659 nodefr = staux%node(ist+n)
660
661 if (st%mpi_grp%rank == nodeto .and. st%mpi_grp%rank == nodefr) then
662 ! do a simple copy on the same rank, also for serial version
663 call states_elec_get_state(staux, gr, ist+n, ik, psi(:, :, 1))
664 call states_elec_set_state(st, gr, ist, ik, psi(:, :, 1))
665 else
666 ! we can use blocking send/recv calls here because we always have different ranks
667 if (st%mpi_grp%rank == nodefr) then
668 call states_elec_get_state(staux, gr, ist+n, ik, psi(:, :, 1))
669 call st%mpi_grp%send(psi(1, 1, 1), gr%np*st%d%dim, mpi_double_complex, nodeto, ist)
670 end if
671 if (st%mpi_grp%rank == nodeto) then
672 call st%mpi_grp%recv(rec_buffer(1, 1), gr%np*st%d%dim, mpi_double_complex, nodefr, ist)
673 call states_elec_set_state(st, gr, ist, ik, rec_buffer(:, :))
674 end if
675 end if
676 end do
677 end do
678
679 safe_deallocate_a(psi)
680 safe_deallocate_a(rec_buffer)
681
682 do ik = 1, st%nik
683 do ist = 1, st%nst
684 st%occ(ist, ik) = staux%occ(n+ist, ik)
685 st%eigenval(ist, ik) = staux%eigenval(n+ist, ik)
686 end do
687 end do
688
689 call states_elec_end(staux)
691 end subroutine states_elec_freeze_orbitals
692
693 ! ---------------------------------------------------------
694 subroutine states_elec_freeze_redistribute_states(st, namespace, mesh, mc, nn)
695 type(states_elec_t), intent(inout) :: st
696 type(namespace_t), intent(in) :: namespace
697 class(mesh_t), intent(in) :: mesh
698 type(multicomm_t), intent(in) :: mc
699 integer, intent(in) :: nn
700
702
703 st%nst = st%nst - nn
704
706 call states_elec_distribute_nodes(st, namespace, mc)
707 call states_elec_allocate_wfns(st, mesh, type_cmplx, packed=.true.)
708
709 safe_deallocate_a(st%eigenval)
710 safe_allocate(st%eigenval(1:st%nst, 1:st%nik))
711 st%eigenval = huge(st%eigenval)
712
713 safe_deallocate_a(st%occ)
714 safe_allocate(st%occ (1:st%nst, 1:st%nik))
715 st%occ = m_zero
716
717
720
721 ! ---------------------------------------------------------
722 subroutine states_elec_freeze_adjust_qtot(st)
723 type(states_elec_t), intent(inout) :: st
724
725 integer :: ik, ist
726
728
729 ! Change the smearing method by fixing the occupations to
730 ! that of the ground-state such that the unfrozen states inherit
731 ! those values.
732 st%smear%method = smear_fixed_occ
733
734 ! Set total charge
735 st%qtot = m_zero
736 do ik = st%d%kpt%start, st%d%kpt%end
737 do ist = st%st_start, st%st_end
738 st%qtot = st%qtot + st%occ(ist, ik) * st%kweights(ik)
739 end do
740 end do
741
742 if (st%parallel_in_states .or. st%d%kpt%parallel) then
743 call comm_allreduce(st%st_kpt_mpi_grp, st%qtot)
744 end if
745
746
748 end subroutine states_elec_freeze_adjust_qtot
749
750
751 ! ---------------------------------------------------------
760 !
761 subroutine states_elec_total_density(st, mesh, total_rho)
762 type(states_elec_t), intent(in) :: st
763 class(mesh_t), intent(in) :: mesh
764 real(real64), contiguous, intent(out) :: total_rho(:,:)
765
766 integer :: is
767
769
770 call lalg_copy(mesh%np, st%d%nspin, st%rho, total_rho)
771
772 if (allocated(st%rho_core)) then
773 do is = 1, st%d%spin_channels
774 call lalg_axpy(mesh%np, m_one/st%d%spin_channels, st%rho_core, total_rho(:, is))
775 end do
776 end if
777
778 ! Add, if it exists, the frozen density from the inner orbitals.
779 if (allocated(st%frozen_rho)) then
780 call lalg_axpy(mesh%np, st%d%nspin, m_one, st%frozen_rho, total_rho)
781 end if
782
784 end subroutine states_elec_total_density
785
786#include "undef.F90"
787#include "real.F90"
788#include "density_inc.F90"
790#include "undef.F90"
791#include "complex.F90"
792#include "density_inc.F90"
793
794end module density_oct_m
795
796
797!! Local Variables:
798!! mode: f90
799!! coding: utf-8
800!! 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
type(accel_kernel_t), target, save, public kernel_density_real
Definition: accel.F90:270
integer function, public accel_kernel_block_size(kernel)
Definition: accel.F90:1167
subroutine, public accel_free_buffer(this, async)
Definition: accel.F90:1004
subroutine, public accel_finish()
Definition: accel.F90:1077
type(accel_kernel_t), target, save, public kernel_density_spinors
Definition: accel.F90:272
integer, parameter, public accel_mem_read_write
Definition: accel.F90:184
type(accel_kernel_t), target, save, public kernel_density_complex
Definition: accel.F90:271
integer, parameter, public accel_mem_read_only
Definition: accel.F90:184
This module implements batches of mesh functions.
Definition: batch.F90:135
integer, parameter, public batch_not_packed
functions are stored in CPU memory, unpacked order
Definition: batch.F90:286
integer, parameter, public batch_device_packed
functions are stored in device memory in packed order
Definition: batch.F90:286
integer, parameter, public batch_packed
functions are stored in CPU memory, in transposed (packed) order
Definition: batch.F90:286
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:118
This module implements a calculator for the density and defines related functions.
Definition: density.F90:122
subroutine, public density_calc_end(this, allreduce, symmetrize)
Finalize the density calculation.
Definition: density.F90:565
subroutine density_calc_state(this, psib, istin)
Calculate contribution to the density from state istin.
Definition: density.F90:248
subroutine, public density_calc_accumulate(this, psib, async)
Accumulate weighted orbital densities for a batch psib.
Definition: density.F90:396
subroutine, public states_elec_freeze_adjust_qtot(st)
Definition: density.F90:818
subroutine, public states_elec_freeze_redistribute_states(st, namespace, mesh, mc, nn)
Definition: density.F90:790
subroutine, public density_calc_init(this, st, gr, density)
initialize the density calculator
Definition: density.F90:190
subroutine density_calc_pack(this, async)
prepare the density calculator for GPU use
Definition: density.F90:227
subroutine, public ddensity_accumulate_grad(space, mesh, st, psib, grad_psib, grad_rho)
Accumulate within a batch.
Definition: density.F90:952
subroutine, public zdensity_accumulate_grad(space, mesh, st, psib, grad_psib, grad_rho)
Accumulate within a batch.
Definition: density.F90:1141
subroutine, public states_elec_total_density(st, mesh, total_rho)
This routine calculates the total electronic density.
Definition: density.F90:857
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:659
subroutine, public density_calc(st, gr, density, istin)
Computes the density from the orbitals in st.
Definition: density.F90:616
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:191
real(real64), parameter, public m_epsilon
Definition: global.F90:207
real(real64), parameter, public m_one
Definition: global.F90:192
This module implements the underlying real-space grid.
Definition: grid.F90:119
subroutine, public dgrid_symmetrize_scalar_field(gr, field, suppress_warning)
Definition: grid.F90:672
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_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:410
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:147
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:631
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:554
integer, parameter, public smear_fixed_occ
Definition: smear.F90:173
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:135
type(type_t), public type_cmplx
Definition: types.F90:136
The calculator object.
Definition: density.F90:171
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:171
Describes mesh distribution to nodes.
Definition: mesh.F90:187
Stores all communicators and groups.
Definition: multicomm.F90:208
The states_elec_t class contains all electronic wave functions.
batches of electronic states
Definition: wfs_elec.F90:141
int true(void)