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