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