Octopus
species_pot.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
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
24 use debug_oct_m
25 use global_oct_m
27 use index_oct_m
28 use, intrinsic :: iso_fortran_env
33 use logrid_oct_m
35 use mesh_oct_m
37 use mpi_oct_m
39 use parser_oct_m
41 use ps_oct_m
44 use space_oct_m
49 use unit_oct_m
51 use volume_oct_m
52
53 implicit none
54
55 private
56 public :: &
64
65 type(mesh_t), pointer :: mesh_p
66 real(real64), allocatable :: rho_p(:)
67 real(real64), allocatable :: grho_p(:, :)
68 real(real64) :: alpha2_p
69 real(real64), pointer :: pos_p(:)
70
71contains
72
73
74 ! ---------------------------------------------------------
75 subroutine species_atom_density(species, namespace, space, latt, pos, mesh, spin_channels, rho)
76 class(species_t), target, intent(in) :: species
77 type(namespace_t), intent(in) :: namespace
78 class(space_t), intent(in) :: space
79 type(lattice_vectors_t), intent(in) :: latt
80 real(real64), intent(in) :: pos(1:space%dim)
81 type(mesh_t), intent(in) :: mesh
82 integer, intent(in) :: spin_channels
83 real(real64), intent(inout) :: rho(:, :)
84
85 integer :: isp, ip, in_points, icell
86 real(real64) :: rr, x, pos_pc(space%dim), nrm, rmax
87 real(real64) :: xx(space%dim), yy(space%dim), rerho, imrho
88 real(real64), allocatable :: dorbital(:)
89 type(ps_t), pointer :: ps
90 type(volume_t) :: volume
91 integer :: in_points_red
92 type(lattice_iterator_t) :: latt_iter
93 integer :: iorb, ii, nn, ll, mm
94 real(real64) :: radius, density
95 type(submesh_t) :: sphere
96
97 push_sub(species_atom_density)
98
99 assert(spin_channels == 1 .or. spin_channels == 2)
100
101 rho = m_zero
102
103 ! build density ...
104 select type (species)
105 type is(species_from_file_t)
107
110
111 type is(soft_coulomb_t)
113
114 class is(allelectron_t)
115
116 do isp = 1, spin_channels
117 do iorb = 1, species%get_niwfs()
118 call species%get_iwf_ilm(iorb, isp, ii, ll, mm)
119 ! For all-electron species, we want to use the principal quantum number
120 call species%get_iwf_n(iorb, isp, nn)
121
122 radius = species%get_iwf_radius(nn, isp)
123 ! make sure that if the spacing is too large, the orbitals fit in a few points at least
124 radius = max(radius, m_two*maxval(mesh%spacing))
125
126 call submesh_init(sphere, space, mesh, latt, pos, radius)
127 safe_allocate(dorbital(1:sphere%np))
128
129 call datomic_orbital_get_submesh(species, sphere, nn, ll, mm, isp, dorbital)
130 ! The occupations are for one type of orbitals, e.g. 2p gets 6 electrons
131 ! So we normalize them by (2*l+1) such that they get distributed evenly
132 ! for each value of m
133 do ip = 1, sphere%np
134 dorbital(ip) = species%conf%occ(ii, isp)/real(2*ll+1, real64) *dorbital(ip)*dorbital(ip)
135 end do
136 call submesh_add_to_mesh(sphere, dorbital, rho(:, isp))
137 safe_deallocate_a(dorbital)
138
139 call submesh_end(sphere)
140 end do
141 end do
142
143 type is (jellium_charge_t)
144 ! We put, for the electron density, the same as the positive density that
145 ! creates the external potential.
146 ! This code is repeated in get_density, and should therefore be cleaned!!!!!
147
148 call volume_init(volume)
149 call volume_read_from_block(volume, namespace, trim(species%rho_string()))
150
151 rmax = latt%max_length()
152 latt_iter = lattice_iterator_t(latt, rmax)
153 rho = m_zero
154 do icell = 1, latt_iter%n_cells
155 yy = latt_iter%get(icell)
156 do ip = 1, mesh%np
157 call mesh_r(mesh, ip, rr, origin = pos, coords = xx)
158 xx = xx + yy
159 rr = norm2(xx)
161 rerho = m_zero
162 if (volume_in_volume(space, volume, xx)) rerho = m_one
163 rho(ip, 1) = rho(ip, 1) + rerho
164 end do
165 end do
166
167 call volume_end(volume)
168
169 if (spin_channels > 1) then
170 rho(:, 1) = m_half*rho(:, 1)
171 rho(:, 2) = rho(:, 1)
172 end if
173
174 ! rescale to match the valence charge
175 do isp = 1, spin_channels
176 x = species%get_zval() / dmf_integrate(mesh, rho(:, isp))
177 rho(1:mesh%np, isp) = x * rho(1:mesh%np, isp)
178 end do
179
181 ! We put, for the electron density, the same as the positive density that
182 ! creates the external potential.
183 ! This code is repeated in get_density, and should therefore be cleaned!!!!!
184
185 rmax = latt%max_length()
186 latt_iter = lattice_iterator_t(latt, rmax)
187 rho = m_zero
188 do icell = 1, latt_iter%n_cells
189 yy = latt_iter%get(icell)
190 do ip = 1, mesh%np
191 call mesh_r(mesh, ip, rr, origin = pos, coords = xx)
192 xx = xx + yy
193 rr = norm2(xx)
194
195 rerho = m_zero
196 call parse_expression(rerho, imrho, space%dim, xx, rr, m_zero, trim(species%rho_string()))
197 rho(ip, 1) = rho(ip, 1) + rerho
198 end do
199 end do
200
201 if (spin_channels > 1) then
202 rho(:, 1) = m_half*rho(:, 1)
203 rho(:, 2) = rho(:, 1)
204 end if
205
206 ! rescale to match the valence charge
207 do isp = 1, spin_channels
208 x = species%get_zval() / dmf_integrate(mesh, rho(:, isp))
209 rho(1:mesh%np, isp) = x * rho(1:mesh%np, isp)
210 end do
211
212
213 type is (jellium_sphere_t) ! ... from jellium
214 in_points = 0
215 do ip = 1, mesh%np
216 call mesh_r(mesh, ip, rr, origin = pos)
217 if (rr <= species%radius()) then
218 in_points = in_points + 1
219 end if
220 end do
221
222 if (mesh%parallel_in_domains) then
223 call mesh%mpi_grp%allreduce(in_points, in_points_red, 1, mpi_integer, mpi_sum)
224 in_points = in_points_red
225 end if
226
227 if (in_points > 0) then
228 ! This probably should be done inside the mesh_function_oct_m module.
229
230 if (mesh%use_curvilinear) then
231 do ip = 1, mesh%np
232 call mesh_r(mesh, ip, rr, origin = pos)
233 if (rr <= species%radius()) then
234 rho(ip, 1:spin_channels) = species%get_zval() / &
235 (mesh%vol_pp(ip) * real(in_points*spin_channels, real64) )
236 end if
237 end do
238 else
239 do ip = 1, mesh%np
240 call mesh_r(mesh, ip, rr, origin = pos)
241 if (rr <= species%radius()) then
242 rho(ip, 1:spin_channels) = species%get_zval() / &
243 (mesh%vol_pp(1) * real(in_points * spin_channels, real64) )
244 end if
245 end do
246 end if
247 end if
248
249 type is (jellium_slab_t) ! ... from jellium slab
250 density = species%get_density(mesh%box%bounding_box_l) / spin_channels
251
252 do ip = 1, mesh%np
253 rr = abs(mesh%x(3, ip) - pos(3))
254 if (rr <= species%thickness() / m_two) then
255 rho(ip, 1:spin_channels) = density
256 end if
257 end do
258
259 class is (pseudopotential_t)
260 ! ...from pseudopotentials
261
262 ps => species%ps
263
264 if (ps_has_density(ps)) then
265
266 assert(allocated(ps%density))
267
268 rmax = m_zero
269 do isp = 1, spin_channels
270 rmax = max(rmax, ps%density(isp)%x_threshold)
271 end do
272
273 latt_iter = lattice_iterator_t(latt, rmax)
274 do icell = 1, latt_iter%n_cells
275 pos_pc = pos + latt_iter%get(icell)
276 do ip = 1, mesh%np
277 call mesh_r(mesh, ip, rr, origin = pos_pc)
278 rr = max(rr, r_small)
279
280 do isp = 1, spin_channels
281 if (rr >= spline_range_max(ps%density(isp))) cycle
282 rho(ip, isp) = rho(ip, isp) + spline_eval(ps%density(isp), rr)
283 end do
284
285 end do
286 end do
287
288 else
289
290 !we use the square root of the short-range local potential, just to put something that looks like a density
291
292 latt_iter = lattice_iterator_t(latt, ps%vl%x_threshold)
293 do icell = 1, latt_iter%n_cells
294 pos_pc = pos + latt_iter%get(icell)
295 do ip = 1, mesh%np
296 call mesh_r(mesh, ip, rr, origin = pos_pc)
297 rr = max(rr, r_small)
298
299 if (rr >= spline_range_max(ps%vl)) cycle
300
301 do isp = 1, spin_channels
302 rho(ip, isp) = rho(ip, isp) + sqrt(abs(spline_eval(ps%vl, rr)))
303 end do
304
305 end do
306 end do
307
308 ! normalize
309 nrm = m_zero
310 do isp = 1, spin_channels
311 nrm = nrm + dmf_integrate(mesh, rho(:, isp))
312 end do
313
314 rho(1:mesh%np, 1:spin_channels) = rho(1:mesh%np, 1:spin_channels)*species%get_zval()/nrm
315
316 end if
317 class default
318 assert(.false.)
319 end select
320
321 pop_sub(species_atom_density)
322 contains
323 subroutine generate_uniform_density()
324 do isp = 1, spin_channels
325 rho(1:mesh%np, isp) = m_one
326 x = (species%get_zval()/real(spin_channels, real64) ) / dmf_integrate(mesh, rho(:, isp))
327 rho(1:mesh%np, isp) = x * rho(1:mesh%np, isp)
328 end do
329 end subroutine generate_uniform_density
330 end subroutine species_atom_density
331
332 ! ---------------------------------------------------------
333 ! A non periodized version of the routine species_atom_density
334 ! This is used for the Hirshfeld routines
335 ! TODO: implement it for other approaches than pseudo potentials.
336 subroutine species_atom_density_np(species, namespace, pos, mesh, spin_channels, rho)
337 class(species_t), target, intent(in) :: species
338 type(namespace_t), intent(in) :: namespace
339 real(real64), intent(in) :: pos(:)
340 type(mesh_t), intent(in) :: mesh
341 integer, intent(in) :: spin_channels
342 real(real64), intent(inout) :: rho(:, :)
343
344 integer :: isp, ip
345 real(real64) :: rr, nrm
346 type(ps_t), pointer :: ps
347
349
350 call profiling_in("SPECIES_ATOM_DEN_NP")
351
352 rho = m_zero
353 select type(species)
354 class is(pseudopotential_t)
355 ! ...from pseudopotentials
356
357 ps => species%ps
358 if (ps_has_density(ps)) then
359
360 assert(allocated(ps%density))
361
362 !$omp parallel private(ip, rr)
363 do isp = 1, spin_channels
364 !$omp do
365 do ip = 1, mesh%np
366 call mesh_r(mesh, ip, rr, origin = pos)
367 if (rr >= spline_range_max(ps%density(isp))) cycle
368 rr = max(rr, r_small)
369 rho(ip, isp) = rho(ip, isp) + spline_eval(ps%density(isp), rr)
370 end do
371 !$omp end do nowait
372 end do
373 !$omp end parallel
374
375 else
376
377 !we use the square root of the short-range local potential, just to put something that looks like a density
378
379 do ip = 1, mesh%np
380 call mesh_r(mesh, ip, rr, origin = pos)
381 rr = max(rr, r_small)
382
383 if (rr >= spline_range_max(ps%vl)) cycle
384
385 do isp = 1, spin_channels
386 rho(ip, isp) = rho(ip, isp) + sqrt(abs(spline_eval(ps%vl, rr)))
387 end do
388
389 end do
390
391 ! normalize
392 nrm = m_zero
393 do isp = 1, spin_channels
394 nrm = nrm + dmf_integrate(mesh, rho(:, isp))
395 end do
396
397 rho(1:mesh%np, 1:spin_channels) = rho(1:mesh%np, 1:spin_channels)*species%get_zval()/nrm
398
399 end if
400 class default
401 call messages_not_implemented('species_atom_density_np for non-pseudopotential species', namespace=namespace)
402
403 end select
404
405 call profiling_out("SPECIES_ATOM_DEN_NP")
406
408 end subroutine species_atom_density_np
409
410 ! ---------------------------------------------------------
411 !! Non-periodic version of the above routine
412 subroutine species_atom_density_derivative_np(species, namespace, pos, mesh, spin_channels, drho)
413 class(species_t), target, intent(in) :: species
414 type(namespace_t), intent(in) :: namespace
415 real(real64), intent(in) :: pos(:)
416 type(mesh_t), intent(in) :: mesh
417 integer, intent(in) :: spin_channels
418 real(real64), intent(inout) :: drho(:, :)
419
420 integer :: isp, ip
421 real(real64) :: rr
422 type(ps_t), pointer :: ps
423
425
426 call profiling_in("SPECIES_ATOM_DEN_DER_NP")
427
428 select type(species)
429 class is(pseudopotential_t)
430 ps => species%ps
432 if (ps_has_density(ps)) then
433 !$omp parallel private(ip, rr, isp)
434 do isp = 1, spin_channels
435 !$omp do
436 do ip = 1, mesh%np
437 call mesh_r(mesh, ip, rr, origin = pos)
438 if (rr >= spline_range_max(ps%density_der(isp))) cycle
439 rr = max(rr, r_small)
440 drho(ip, isp) = drho(ip, isp) + spline_eval(ps%density_der(isp), rr)
441 end do
442 !$omp end do nowait
443 end do
444 !$omp end parallel
445
446 else
447 call messages_write('The pseudopotential for')
448 call messages_write(species%get_label())
449 call messages_write(' does not contain the density.')
450 call messages_fatal(namespace=namespace)
451 end if
452 class default
453 assert(.false.)
454 end select
455
456 call profiling_out("SPECIES_ATOM_DEN_DER_NP")
457
460
461
462 ! ---------------------------------------------------------
463 !! Gradient of the atomic density, if available
464 subroutine species_atom_density_grad(species, namespace, space, latt, pos, mesh, spin_channels, drho)
465 class(species_t), target, intent(in) :: species
466 type(namespace_t), intent(in) :: namespace
467 class(space_t), intent(in) :: space
468 type(lattice_vectors_t), intent(in) :: latt
469 real(real64), intent(in) :: pos(1:space%dim)
470 type(mesh_t), intent(in) :: mesh
471 integer, intent(in) :: spin_channels
472 real(real64), intent(inout) :: drho(:, :, :)
473
474 integer :: isp, ip, icell, idir
475 real(real64) :: rr, pos_pc(space%dim), range, spline
476 type(ps_t), pointer :: ps
477 type(lattice_iterator_t) :: latt_iter
478
480
481 assert(spin_channels == 1 .or. spin_channels == 2)
482
483 drho = m_zero
484
485 ! build density ...
486 select type(species)
487 class is(pseudopotential_t)
488 ps => species%ps
489 ! ...from pseudopotentials
490
491 if (ps_has_density(ps)) then
492
493 range = ps%density_der(1)%x_threshold
494 if (spin_channels == 2) range = max(range, ps%density_der(2)%x_threshold)
495 latt_iter = lattice_iterator_t(latt, range)
496
497 do icell = 1, latt_iter%n_cells
498 pos_pc = pos + latt_iter%get(icell)
499
500 do ip = 1, mesh%np
501 call mesh_r(mesh, ip, rr, origin = pos_pc)
502 rr = max(rr, r_small)
503
504 do isp = 1, spin_channels
505 if (rr >= spline_range_max(ps%density_der(isp))) cycle
506 spline = spline_eval(ps%density_der(isp), rr)
508 if(abs(spline) < 1e-150_real64) cycle
509
510 do idir = 1, space%dim
511 drho(ip, isp, idir) = drho(ip, isp, idir) - spline*(mesh%x(idir, ip) - pos_pc(idir))/rr
512 end do
513 end do
514 end do
515 end do
516
517 else
518 call messages_write('The pseudopotential for')
519 call messages_write(species%get_label())
520 call messages_write(' does not contain the density.')
521 call messages_fatal(namespace=namespace)
522 end if
523
524 class default
525 call messages_not_implemented('species_atom_density_grad for non-pseudopotential species', namespace=namespace)
526
527 end select
528
530 end subroutine species_atom_density_grad
531
532 ! ---------------------------------------------------------
533
534 subroutine species_get_long_range_density(species, namespace, space, latt, pos, mesh, rho, sphere_inout, nlr_x)
535 class(species_t), target, intent(in) :: species
536 type(namespace_t), intent(in) :: namespace
537 class(space_t), intent(in) :: space
538 type(lattice_vectors_t), intent(in) :: latt
539 real(real64), target, intent(in) :: pos(1:space%dim)
540 class(mesh_t), target, intent(in) :: mesh
541 real(real64), intent(out) :: rho(:)
542 type(submesh_t), optional, target, intent(inout) :: sphere_inout
543 real(real64), optional, intent(inout) :: nlr_x(:,:)
544
545 type(root_solver_t) :: rs
546 logical :: conv
547 real(real64) :: startval(space%dim)
548 real(real64) :: delta, alpha, xx(space%dim), yy(space%dim), rr, imrho1, rerho
549 real(real64) :: dist2_min
550 integer :: icell, ipos, ip, idir, rankmin
551 type(lattice_iterator_t) :: latt_iter
552 type(ps_t), pointer :: ps
553 type(volume_t) :: volume
554 type(submesh_t), target :: sphere_local
555 type(submesh_t), pointer :: sphere
556 logical :: have_point
557 real(real64), allocatable :: rho_sphere(:)
558 real(real64), parameter :: threshold = 1e-6_real64
559 real(real64) :: norm_factor, range, radius, radius_nlr, radius_vl
560
562
563 call profiling_in("SPECIES_LR_DENSITY")
564
565 if(present(nlr_x)) then
566 assert(species%is_ps())
567 end if
568
569 select type (species)
570 type is(pseudopotential_t)
571 ps => species%ps
572 radius_nlr = spline_x_threshold(ps%nlr, threshold)
573 if (present(sphere_inout)) then
574 radius_vl = ps%vl%x_threshold*1.05_real64
575 radius = max(radius_nlr, radius_vl)
576 call submesh_init(sphere_inout, space, mesh, latt, pos, radius)
577 sphere => sphere_inout
578 else
579 radius = radius_nlr
580 call submesh_init(sphere_local, space, mesh, latt, pos, radius)
581 sphere => sphere_local
582 endif
583
584 safe_allocate(rho_sphere(1:sphere%np))
585 if (.not. present(sphere_inout) .and. sphere%np > 0) then
586 call lalg_copy(sphere%np, sphere%r, rho_sphere)
587 call spline_eval_vec(ps%nlr, sphere%np, rho_sphere)
588 else
589 do ip = 1, sphere%np
590 if(sphere%r(ip) <= radius_nlr) then
591 rho_sphere(ip) = spline_eval(ps%nlr, sphere%r(ip))
592 else
593 rho_sphere(ip) = m_zero
594 endif
595 end do
596 end if
597
598 rho(1:mesh%np) = m_zero
599
600 ! A small amount of charge is missing with the cutoff, we
601 ! renormalize so that the long range potential is exact
602 norm_factor = abs(species%get_zval()/dsm_integrate(mesh, sphere, rho_sphere))
603 do ip = 1, sphere%np
604 rho(sphere%map(ip)) = rho(sphere%map(ip)) + norm_factor*rho_sphere(ip)
605 end do
606
607 if (present(nlr_x)) then
608 do idir = 1, space%dim
609 do ip = 1, sphere%np
610 nlr_x(sphere%map(ip), idir) = nlr_x(sphere%map(ip), idir) + norm_factor*rho_sphere(ip)*sphere%rel_x(idir, ip)
611 end do
612 end do
613 end if
614
615 safe_deallocate_a(rho_sphere)
616 nullify(ps)
617 if ( .not. present(sphere_inout) ) then
618 call submesh_end(sphere)
619 end if
620 nullify(sphere)
621
622 type is (full_delta_t)
623
624 rho(1:mesh%np) = m_zero
625
626 ipos = mesh_nearest_point(mesh, pos, dist2_min, rankmin)
627 have_point = .true.
628 if (mesh%mpi_grp%rank /= rankmin) have_point = .false.
630 if (have_point) then
631 if (mesh%use_curvilinear) then
632 rho(ipos) = -species%get_z()/mesh%vol_pp(ipos)
633 else
634 rho(ipos) = -species%get_z()/mesh%vol_pp(1)
635 end if
636 end if
637
638 write(message(1), '(3a,f5.2,3a)') &
639 "Info: species_full_delta species ", trim(species%get_label()), &
640 " atom displaced ", units_from_atomic(units_out%length, sqrt(dist2_min)), &
641 " [ ", trim(units_abbrev(units_out%length)), " ]"
642 call messages_info(1, namespace=namespace)
643
644 type is (full_gaussian_t)
645
646 ! periodic copies are not considered in this routine
647 if (space%is_periodic()) then
648 call messages_not_implemented("species_full_gaussian for periodic systems", namespace=namespace)
649 end if
650
651 ! We need to work with \xi and \xi_0, not x(\xi) and x(\xi_0) as we do now
652 ! for the argument of the Gaussian
653 if (mesh%use_curvilinear) then
654 call messages_not_implemented("species_full_gaussian with curvilinear coordinates", namespace=namespace)
655 end if
656
657 ! --------------------------------------------------------------
658 ! Constructs density for an all-electron atom with the procedure
659 ! sketched in Modine et al. [Phys. Rev. B 55, 10289 (1997)],
660 ! section II.B
661 ! --------------------------------------------------------------
662
663 safe_allocate(rho_p(1:mesh%np))
664 safe_allocate(grho_p(1:mesh%np, 1:space%dim))
665
666 mesh_p => mesh
667 pos_p => pos
668
669 ! Initial guess.
670 delta = mesh%spacing(1)
671 alpha = sqrt(m_two)*species%get_sigma()*delta
672 alpha2_p = alpha**2 ! global copy of alpha
673
674 ! the dim variables are the position of the delta function
675 startval(1:space%dim) = pos
676
677 ! solve equation
678 ! Setting a tolerance such that the distance to the first moment is smaller than 1e-5 Bohr
679 call root_solver_init(rs, namespace, space%dim, solver_type=root_newton, maxiter=500, abs_tolerance=1.0e-10_real64)
680 call droot_solver_run(rs, func, xx, conv, startval=startval)
681
682 if (.not. conv) then
683 write(message(1),'(a)') 'Root finding in species_get_density did not converge.'
684 call messages_fatal(1, namespace=namespace)
685 end if
686
687 if(debug%info .and. space%dim == 3) then
688 write(message(1),'(a,3(f6.3,a))') 'Debug: Gaussian charge position (', xx(1), ', ', xx(2), ', ', xx(3), ')'
689 call messages_info(1, namespace=namespace)
690 end if
691
692 ! we want a charge of -Z
693 rho = -species%get_z()*rho_p
694
695 nullify(mesh_p)
696 nullify(pos_p)
697 safe_deallocate_a(grho_p)
698 safe_deallocate_a(rho_p)
699
700 type is (full_anc_t)
701
702 rho = m_zero
703
704 type is(jellium_charge_t)
705
706 call volume_init(volume)
707 call volume_read_from_block(volume, namespace, trim(species%rho_string()))
708
709 range = latt%max_length()
710 latt_iter = lattice_iterator_t(latt, range)
711
712 rho = m_zero
713 do icell = 1, latt_iter%n_cells
714 yy = latt_iter%get(icell)
715 do ip = 1, mesh%np
716 call mesh_r(mesh, ip, rr, origin = pos, coords = xx)
717 xx = xx + yy
718 rr = norm2(xx)
719
720 rerho = m_zero
721 if (volume_in_volume(space, volume, xx)) rerho = m_one
722 rho(ip) = rho(ip) - rerho
723 end do
724 end do
725
726 call volume_end(volume)
727
729
730 range = latt%max_length()
731 latt_iter = lattice_iterator_t(latt, range)
732
733 rho = m_zero
734 do icell = 1, latt_iter%n_cells
735 yy = latt_iter%get(icell)
736 do ip = 1, mesh%np
737 call mesh_r(mesh, ip, rr, origin = pos, coords = xx)
738 xx = xx + yy
739 rr = norm2(xx)
740
741 rerho = m_zero
742 call parse_expression(rerho, imrho1, space%dim, xx, rr, m_zero, trim(species%rho_string()))
743 rho(ip) = rho(ip) - rerho
744 end do
745 end do
746
747 rr = species%get_zval() / abs(dmf_integrate(mesh, rho(:)))
748 rho(1:mesh%np) = rr * rho(1:mesh%np)
749
750 class default
751 assert(.false.)
752 end select
753
754 call profiling_out("SPECIES_LR_DENSITY")
756 end subroutine species_get_long_range_density
757
758
759 ! ---------------------------------------------------------
760 subroutine func(xin, ff, jacobian)
761 real(real64), intent(in) :: xin(:)
762 real(real64), intent(out) :: ff(:), jacobian(:,:)
763
764 real(real64), allocatable :: xrho(:)
765 integer :: idir, jdir, dim, ip
766
767 push_sub(func)
768
769 dim = mesh_p%box%dim
770
771 call getrho(dim, xin)
772 safe_allocate(xrho(1:mesh_p%np))
773
774 ! First, we calculate the function ff.
775 do idir = 1, dim
776 !$omp parallel do simd
777 do ip = 1, mesh_p%np
778 xrho(ip) = rho_p(ip) * mesh_p%x_t(ip, idir)
779 end do
780 ff(idir) = dmf_integrate(mesh_p, xrho) - pos_p(idir)
781 end do
782
783 ! Now the jacobian.
784 do idir = 1, dim
785 do jdir = 1, dim
786 !$omp parallel do simd
787 do ip = 1, mesh_p%np
788 xrho(ip) = grho_p(ip, jdir) * mesh_p%x_t(ip, idir)
789 end do
790 jacobian(idir, jdir) = dmf_integrate(mesh_p, xrho)
791 end do
792 end do
793
794 safe_deallocate_a(xrho)
795 pop_sub(func)
796 end subroutine func
797
798 ! ---------------------------------------------------------
799 subroutine species_get_nlcc(species, space, latt, pos, mesh, rho_core, accumulate)
800 class(species_t), target, intent(in) :: species
801 class(space_t), intent(in) :: space
802 type(lattice_vectors_t), intent(in) :: latt
803 real(real64), intent(in) :: pos(1:space%dim)
804 class(mesh_t), intent(in) :: mesh
805 real(real64), intent(inout) :: rho_core(:)
806 logical, optional, intent(in) :: accumulate
807
808 real(real64) :: center(space%dim), rr
809 integer :: icell, ip
810 type(lattice_iterator_t) :: latt_iter
811 type(ps_t), pointer :: ps
812
813 push_sub(species_get_nlcc)
814
815 ! only for 3D pseudopotentials, please
816 select type(species)
817 class is(pseudopotential_t)
818 ps => species%ps
819 if (.not. optional_default(accumulate, .false.)) rho_core = m_zero
820
821 latt_iter = lattice_iterator_t(latt, ps%core%x_threshold)
822 do icell = 1, latt_iter%n_cells
823 center = pos + latt_iter%get(icell)
824 do ip = 1, mesh%np
825 rr = norm2(mesh%x(1:space%dim, ip) - center)
826 if (rr < spline_range_max(ps%core)) then
827 rho_core(ip) = rho_core(ip) + spline_eval(ps%core, rr)
828 end if
829 end do
830 end do
831 class default
832 if (.not. optional_default(accumulate, .false.)) rho_core = m_zero
833 end select
834
835 pop_sub(species_get_nlcc)
836 end subroutine species_get_nlcc
837
838 ! ---------------------------------------------------------
839 ! Return the density of a normalized Gaussian centered on xin
840 ! as well as its gradient with respect to the central position
841 subroutine getrho(dim, xin)
842 integer, intent(in) :: dim
843 real(real64), intent(in) :: xin(1:dim)
844
845 integer :: ip, idir
846 real(real64) :: r2, chi(dim), norm, threshold
847
848 push_sub(getrho)
849
850 ! We set here a threshold of 0.001 for the tail of the Gaussian, similar to what we do for the
851 ! pseudopotentials. The value of the threshold corresponds to the default for pseudopotentials
852 threshold = -log(0.001_real64)*alpha2_p
853
854 do ip = 1, mesh_p%np
855 ! This is not correct for curvilinear meshes
856 chi(1:dim) = mesh_p%x(1:dim, ip)
857 r2 = sum((chi - xin(1:dim))**2)
858
859 if (r2 < threshold) then
860 rho_p(ip) = exp(-r2/alpha2_p)
861 else
862 rho_p(ip) = m_zero
863 end if
864
865 do idir = 1, dim
866 grho_p(ip, idir) = (chi(idir) - xin(idir)) * rho_p(ip)
867 end do
868 end do
869
870 norm = dmf_integrate(mesh_p, rho_p)
871 call lalg_scal(mesh_p%np, m_one/norm, rho_p)
872 call lalg_scal(mesh_p%np, dim, m_two/alpha2_p/norm, grho_p)
873
874 pop_sub(getrho)
875 end subroutine getrho
876
877
878 ! ---------------------------------------------------------
880 subroutine species_get_local(species, namespace, space, latt, pos, mesh, vl)
881 class(species_t), target, intent(in) :: species
882 type(namespace_t), intent(in) :: namespace
883 class(space_t), intent(in) :: space
884 type(lattice_vectors_t), intent(in) :: latt
885 real(real64), intent(in) :: pos(1:space%dim)
886 type(mesh_t), intent(in) :: mesh
887 real(real64), intent(out) :: vl(:)
888
889 real(real64) :: a1, a2, Rb2, range, density ! for jellium
890 real(real64) :: xx(space%dim), pos_pc(space%dim), r, r2, threshold
891 integer :: ip, err, icell
892 complex(real64) :: zpot
893 type(lattice_iterator_t) :: latt_iter
894 real(real64) :: aa, bb
895
896 push_sub_with_profile(species_get_local)
897
898 select type(species)
899
900 type is (soft_coulomb_t)
901
902 call parse_variable(namespace, 'SpeciesProjectorSphereThreshold', 0.001_real64, threshold)
903
904 !Assuming that we want to take the contribution from all replica that contributes up to 0.001
905 ! to the center of the cell, we arrive to a range of 1000 a.u..
906 latt_iter = lattice_iterator_t(latt, species%get_zval() / threshold)
907 vl = m_zero
908 do icell = 1, latt_iter%n_cells
909 pos_pc = pos + latt_iter%get(icell)
910 do ip = 1, mesh%np
911 call mesh_r(mesh, ip, r, origin = pos_pc)
912 r2 = r*r
913 vl(ip) = vl(ip) -species%get_zval()/sqrt(r2+species%get_softening2())
914 end do
915 end do
916
917 type is (species_user_defined_t)
918 !TODO: we should control the value of 5 by a variable.
919 range = 5.0_real64 * latt%max_length()
920 latt_iter = lattice_iterator_t(latt, range)
921 vl = m_zero
922 do icell = 1, latt_iter%n_cells
923 pos_pc = pos + latt_iter%get(icell)
924 do ip = 1, mesh%np
925 call mesh_r(mesh, ip, r, origin = pos_pc, coords = xx)
926
927 zpot = species%user_pot(space%dim, xx, r)
928 vl(ip) = vl(ip) + real(zpot, real64)
929 end do
930 end do
931
932 type is(species_from_file_t)
933
934 call dio_function_input(trim(species%get_filename()), namespace, space, mesh, vl, err)
935 if (err /= 0) then
936 write(message(1), '(a)') 'Error loading file '//trim(species%get_filename())//'.'
937 write(message(2), '(a,i4)') 'Error code returned = ', err
938 call messages_fatal(2, namespace=namespace)
939 end if
940
941 type is(jellium_sphere_t)
942
943 assert(.not. space%is_periodic())
944
945 a1 = species%get_z()/(m_two*species%radius()**3)
946 a2 = species%get_z()/species%radius()
947 rb2= species%radius()**2
948
949 do ip = 1, mesh%np
950
951 xx = mesh%x(:, ip) - pos(1:space%dim)
952 r = norm2(xx)
953
954 if (r <= species%radius()) then
955 vl(ip) = (a1*(r*r - rb2) - a2)
956 else
957 vl(ip) = -species%get_z()/r
958 end if
959
960 end do
961
962 type is (jellium_slab_t)
963
964 ! Electrostatic potential from an infinite slab of thickness species%thickness
965 ! Potential and electric fields are continuous at +/- L/2
966 density = species%get_density(mesh%box%bounding_box_l)
967 a1 = m_four * m_pi * density * species%thickness() / m_two
968
969 do ip = 1, mesh%np
970
971 r = abs(mesh%x(3, ip) - pos(3))
972
973 if (r <= species%thickness()/m_two) then
974 vl(ip) = a1 * (r * r / species%thickness() + species%thickness() / m_four)
975 else
976 vl(ip) = a1 * r
977 end if
978
979 end do
980
981 class is (pseudopotential_t)
982
983 assert(.not. space%is_periodic())
984
985 !$omp parallel do private(r)
986 do ip = 1, mesh%np
987 r = norm2(mesh%x(:, ip) - pos)
988 vl(ip) = long_range_potential(r, species%ps%sigma_erf, species%ps%z_val)
989 end do
990
991 type is (full_anc_t)
992 ! periodic copies are not considered in this routine
993 if (space%is_periodic()) then
994 call messages_experimental("species_full_anc for periodic systems", namespace=namespace)
995 end if
996
997 aa = species%a()
998 bb = species%b()
999 assert(bb < m_zero) ! To be sure it was computed
1000
1001 ! Evaluation of the scaled potential, see Eq. 19
1002 do ip = 1, mesh%np
1003 r2 = sum((mesh%x(:, ip) - pos)**2)*(species%get_z()*aa)**2
1004 if(r2 > r_small**2) then
1005 r = sqrt(r2)
1006 vl(ip) = -m_half &
1007 - (erf(r) + m_two*(aa*bb + m_one/sqrt(m_pi))*r*exp(-r2))/r*aa &
1008 + (erf(r) + m_two*(aa*bb + m_one/sqrt(m_pi))*r*exp(-r2))**2*m_half &
1009 + (-m_two*aa**2*bb - m_four*aa/sqrt(m_pi) &
1010 + m_four*aa*(aa*bb + m_one/sqrt(m_pi))*r2)*exp(-r2)*m_half
1011 else ! Eq. 10
1012 vl(ip) = -m_half - m_three * aa**2*bb - 6.0_real64*aa/sqrt(m_pi)
1013 end if
1014 vl(ip) = vl(ip) * (species%get_z())**2
1015 end do
1016
1017 class default
1018 vl(1:mesh%np) = m_zero
1019 end select
1020
1021 pop_sub_with_profile(species_get_local)
1022 end subroutine species_get_local
1023
1024end module species_pot_oct_m
1025
1026!! Local Variables:
1027!! mode: f90
1028!! coding: utf-8
1029!! End:
Copies a vector x, to a vector y.
Definition: lalg_basic.F90:188
scales a vector by a constant
Definition: lalg_basic.F90:159
Both the filling of the function, and the retrieval of the values may be done using single- or double...
Definition: splines.F90:166
double log(double __x) __attribute__((__nothrow__
double exp(double __x) __attribute__((__nothrow__
subroutine, public datomic_orbital_get_submesh(species, submesh, ii, ll, mm, ispin, phi, derivative)
type(debug_t), save, public debug
Definition: debug.F90:158
real(real64), parameter, public m_two
Definition: global.F90:193
real(real64), parameter, public r_small
Definition: global.F90:183
real(real64), parameter, public m_zero
Definition: global.F90:191
real(real64), parameter, public m_four
Definition: global.F90:195
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:189
real(real64), parameter, public m_half
Definition: global.F90:197
real(real64), parameter, public m_one
Definition: global.F90:192
real(real64), parameter, public m_three
Definition: global.F90:194
This module implements the index, used for the mesh points.
Definition: index.F90:124
subroutine, public dio_function_input(filename, namespace, space, mesh, ff, ierr, map)
Reads a mesh function from file filename, and puts it into ff. If the map argument is passed,...
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
integer function, public mesh_nearest_point(mesh, pos, dmin, rankmin)
Returns the index of the point which is nearest to a given vector position pos.
Definition: mesh.F90:386
pure subroutine, public mesh_r(mesh, ip, rr, origin, coords)
return the distance to the origin for a given grid point
Definition: mesh.F90:342
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1067
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:409
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1039
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:593
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
Definition: ps.F90:116
pure logical function, public ps_has_density(ps)
Definition: ps.F90:1634
real(real64) pure function, public long_range_potential(r, sigma, z_val)
Evaluate the long-range potential at a given distance.
Definition: ps.F90:750
integer, parameter, public root_newton
subroutine, public root_solver_init(rs, namespace, dimensionality, solver_type, maxiter, rel_tolerance, abs_tolerance)
subroutine, public droot_solver_run(rs, func, root, success, startval)
subroutine, public species_get_local(species, namespace, space, latt, pos, mesh, vl)
used when the density is not available, or otherwise the Poisson eqn would be used instead
subroutine func(xin, ff, jacobian)
subroutine, public species_atom_density_np(species, namespace, pos, mesh, spin_channels, rho)
subroutine, public species_get_long_range_density(species, namespace, space, latt, pos, mesh, rho, sphere_inout, nlr_x)
subroutine, public species_atom_density_derivative_np(species, namespace, pos, mesh, spin_channels, drho)
subroutine, public species_atom_density_grad(species, namespace, space, latt, pos, mesh, spin_channels, drho)
subroutine getrho(dim, xin)
subroutine, public species_get_nlcc(species, space, latt, pos, mesh, rho_core, accumulate)
subroutine, public species_atom_density(species, namespace, space, latt, pos, mesh, spin_channels, rho)
real(real64) function, public spline_x_threshold(spl, threshold)
Determines the largest value of x for which the spline values are above the threshold.
Definition: splines.F90:1052
real(real64) function, public spline_eval(spl, x)
Definition: splines.F90:441
real(real64) pure function, public spline_range_max(this)
Definition: splines.F90:1091
real(real64) function, public dsm_integrate(mesh, sm, ff, reduce)
Definition: submesh.F90:1091
subroutine, public submesh_end(this)
Definition: submesh.F90:678
subroutine, public submesh_init(this, space, mesh, latt, center, rc)
Definition: submesh.F90:227
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:134
character(len=20) pure function, public units_abbrev(this)
Definition: unit.F90:225
This module defines the unit system, used for input and output.
type(unit_system_t), public units_out
subroutine, public volume_read_from_block(vol, namespace, block_name)
Definition: volume.F90:159
logical function, public volume_in_volume(space, vol, xx)
Definition: volume.F90:227
subroutine, public volume_end(vol)
Definition: volume.F90:151
subroutine, public volume_init(vol)
Definition: volume.F90:145
subroutine generate_uniform_density()
An abstract type for all electron species.
The following class implements a lattice iterator. It allows one to loop over all cells that are with...
Describes mesh distribution to nodes.
Definition: mesh.F90:187
A type storing the information and data about a pseudopotential.
Definition: ps.F90:188
An abstract class for species. Derived classes include jellium, all electron, and pseudopotential spe...
Definition: species.F90:147
A submesh is a type of mesh, used for the projectors in the pseudopotentials It contains points on a ...
Definition: submesh.F90:175
int true(void)