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