Octopus
atomic.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
21module atomic_oct_m
22 use debug_oct_m
23 use global_oct_m
24 use, intrinsic :: iso_fortran_env
25 use logrid_oct_m
29 use sort_oct_m
30 use xc_f03_lib_m
31#ifdef HAVE_LIBXC_FUNCS
32 use xc_f03_funcs_m
33#endif
34 use iso_c_binding
35 implicit none
36
37 private
38 public :: &
39 atomxc, &
40 atomhxc, &
42
44 public :: &
45 valconf_t, &
50
51 integer, parameter :: VALCONF_STRING_LENGTH = 80
52
53 type valconf_t
54 ! Components are public by default
55 integer :: z = 0
56 character(len=3) :: symbol = ""
57 integer :: type = 0
58 integer :: p = 0
59 integer :: n(15) = 0
60 integer :: l(15) = 0
61 real(real64) :: occ(15,2) = m_zero
62 real(real64) :: j(15) = m_zero
63 end type valconf_t
64
65
66contains
67
68 ! ---------------------------------------------------------
69 subroutine valconf_copy(cout, cin)
70 type(valconf_t), intent(out) :: cout
71 type(valconf_t), intent(in) :: cin
72
73 push_sub(valconf_copy)
74
75 cout%z = cin%z
76 cout%symbol = cin%symbol
77 cout%type = cin%type
78 cout%p = cin%p
79 cout%n = cin%n
80 cout%l = cin%l
81 cout%occ = cin%occ
82 cout%j = cin%j
83
84 pop_sub(valconf_copy)
85 end subroutine valconf_copy
86
87 subroutine valconf_sort_by_l(c, lmax)
88 type(valconf_t), intent(inout) :: c
89 integer, intent(in) :: lmax
90
91 integer :: ll
92 integer, allocatable :: order(:)
93 type(valconf_t) :: tmp_c
94
95 push_sub(valconf_sort_by_l)
96
97 call valconf_copy(tmp_c, c)
98
99 ! In order to work in the following, we need to sort the occupations by angular momentum
100 safe_allocate(order(1:lmax+1))
101 call sort(c%l(1:lmax+1), order)
102 do ll = 0, lmax
103 c%occ(ll+1, 1:2) = tmp_c%occ(order(ll+1), 1:2)
104 c%n(ll+1) = tmp_c%n(order(ll+1))
105 c%j(ll+1) = tmp_c%j(order(ll+1))
106 end do
107 safe_deallocate_a(order)
108
109 pop_sub(valconf_sort_by_l)
110 end subroutine valconf_sort_by_l
111
112 !In case of spin-polarized calculations, we properly distribute the electrons
113 subroutine valconf_unpolarized_to_polarized(conf)
114 type(valconf_t), intent(inout) :: conf
115
116 integer :: ii, ll
117 real(real64) :: xx
118
120
121 do ii = 1, conf%p
122 ll = conf%l(ii)
123 xx = conf%occ(ii, 1)
124 conf%occ(ii, 1) = min(xx, real(2*ll+1, real64) )
125 conf%occ(ii, 2) = xx - conf%occ(ii,1)
126 end do
127
130
131
132 ! ---------------------------------------------------------
133 subroutine atomhxc(namespace, functl, g, nspin, dens, v, extra)
134 type(namespace_t), intent(in) :: namespace
135 character(len=*), intent(in) :: functl
136 type(logrid_t), intent(in) :: g
137 integer, intent(in) :: nspin
138 real(real64), intent(in) :: dens(g%nrval, nspin)
139 real(real64), intent(out) :: v(g%nrval, nspin)
140 real(real64), intent(in), optional :: extra(g%nrval)
141
142 character(len=5) :: xcfunc, xcauth
143 integer :: is, ir
144 real(real64), allocatable :: xc(:, :), ve(:, :), rho(:, :)
145 real(real64) :: r2, ex, ec, dx, dc
147 push_sub(atomhxc)
149 safe_allocate( ve(1:g%nrval, 1:nspin))
150 safe_allocate( xc(1:g%nrval, 1:nspin))
151 safe_allocate(rho(1:g%nrval, 1:nspin))
154 rho = m_zero
156 ! To calculate the Hartree term, we put all the density in one variable.
157 do is = 1, nspin
158 rho(:, 1) = rho(:, 1) + dens(:, is)
159 end do
160 ve = m_zero
161 do is = 1, nspin
162 call vhrtre(rho(:, 1), ve(:, is), g%rofi, g%drdi, g%s, g%nrval, g%a)
163 end do
164 ve(1, 1:nspin) = ve(2, 1:nspin)
165
166 select case (functl)
167 case ('LDA')
168 xcfunc = 'LDA'
169 xcauth = 'PZ'
170 case ('GGA')
171 xcfunc = 'GGA'
172 xcauth = 'PBE'
173 case default
174 message(1) = 'Internal Error in atomhxc: unknown functl'
175 call messages_fatal(1, namespace=namespace)
176 end select
177
178 rho(:, :) = dens(:, :)
179 do is = 1, nspin
180 if (present(extra)) rho(:, is) = rho(:, is) + extra(:)/nspin
181 rho(2:g%nrval, is) = rho(2:g%nrval, is)/(m_four*m_pi*g%r2ofi(2:g%nrval))
182 end do
183 r2 = g%rofi(2)/(g%rofi(3)-g%rofi(2))
184 rho(1, 1:nspin) = rho(2, 1:nspin) - (rho(3, 1:nspin)-rho(2, 1:nspin))*r2
185
186 do is = 1, nspin
187 do ir = 1, g%nrval
188 if (rho(ir, is) < m_epsilon) rho(ir, is) = m_zero
189 end do
190 end do
191
192 call atomxc(namespace, xcfunc, xcauth, g%nrval, g%nrval, g%rofi, nspin, rho, ex, ec, dx, dc, xc)
193 v = ve + xc
194
195 safe_deallocate_a(ve)
196 safe_deallocate_a(xc)
197 safe_deallocate_a(rho)
198
199 pop_sub(atomhxc)
200 end subroutine atomhxc
201
202
203!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
242 subroutine atomxc(namespace, FUNCTL, AUTHOR, NR, MAXR, RMESH, NSPIN, DENS, EX, EC, DX, DC, V_XC)
243 type(namespace_t), intent(in) :: namespace
244 character(len=*), intent(in) :: functl, author
245 integer, intent(in) :: nr, maxr
246 real(real64), intent(in) :: rmesh(maxr)
247 integer, intent(in) :: nspin
248 real(real64), intent(in) :: dens(maxr,nspin)
249 real(real64), intent(out) :: dc, dx, ec, ex, v_xc(maxr,nspin)
250
251!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
252! Internal parameters !
253! NN : order of the numerical derivatives: the number of radial !
254! points used is 2*NN+1 !
255!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
256
257 integer, parameter :: nn = 5
258
259!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
260! Fix energy unit: EUNIT=1.0 => Hartrees, !
261! EUNIT=0.5 => Rydbergs, !
262! EUNIT=0.03674903 => eV !
263!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
264
265 real(real64), parameter :: eunit = m_half
266
267!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
268! DVMIN is added to differential of volume to avoid division by zero !
269!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
270
271 real(real64), parameter :: dvmin = 1.0e-12_real64
272
273!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
274! Local variables and arrays !
275!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
276
277 logical :: gga
278 integer :: in, in1, in2, ir, is, jn
279 real(real64) :: aux(maxr), d(nspin), &
280 dgdm(-nn:nn), dgidfj(-nn:nn), drdm, dvol, &
281 f1, f2, gd(3,nspin), sigma(3), vxsigma(3), vcsigma(3), &
282 epsx(1), epsc(1), &
283 dexdd(nspin), dexdgd(3,nspin), decdd(nspin), decdgd(3,nspin)
284 type(xc_f03_func_t) :: x_conf, c_conf
285
286 push_sub(atomxc)
287
288 ! sanity check
289 assert(nspin == 1 .or. nspin == 2)
290
291 ! Set GGA switch
292 if (functl .eq. 'LDA' .or. functl .eq. 'lda' .or. &
293 functl .eq. 'LSD' .or. functl .eq. 'lsd') then
294 gga = .false.
295 else if (functl .eq. 'GGA' .or. functl .eq. 'gga') then
296 gga = .true.
297 else
298 gga = .false.
299 write(message(1),'(a,a)') 'atomxc: Unknown functional ', functl
300 call messages_fatal(1, namespace=namespace)
301 end if
302
303 ! initialize xc functional
304 if (gga) then
305 call xc_f03_func_init(x_conf, xc_gga_x_pbe, nspin)
306 call xc_f03_func_init(c_conf, xc_gga_c_pbe, nspin)
307 else
308 call xc_f03_func_init(x_conf, xc_lda_x, nspin)
309 if (author .eq. 'CA' .or. author .eq. 'ca' .or. &
310 author .eq. 'PZ' .or. author .eq. 'pz') then
311 call xc_f03_func_init(c_conf, xc_lda_c_pz, nspin)
312 else if (author .eq. 'PW92' .or. author .eq. 'pw92') then
313 call xc_f03_func_init(c_conf, xc_lda_c_pw, nspin)
314 else
315 write(message(1),'(a,a)') 'LDAXC: Unknown author ', author
316 call messages_fatal(1, namespace=namespace)
317 end if
318 end if
319
320!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
321! Initialize output !
322!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
323
324 ex = m_zero
325 ec = m_zero
326 dx = m_zero
327 dc = m_zero
328 do is = 1,nspin
329 do ir = 1,nr
330 v_xc(ir,is) = m_zero
331 end do
332 end do
333
334!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
335! Loop on mesh points !
336!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
338 do ir = 1,nr
339
340 ! Find interval of neighbour points to calculate derivatives
341 in1 = max( 1, ir-nn ) - ir
342 in2 = min( nr, ir+nn ) - ir
343
344 ! Find weights of numerical derivative from Lagrange
345 ! interpolation formula
346 do in = in1,in2
347 if (in == 0) then
348 dgdm(in) = 0
349 do jn = in1,in2
350 if (jn /= 0) dgdm(in) = dgdm(in) + m_one / (0 - jn)
351 end do
352 else
353 f1 = 1
354 f2 = 1
355 do jn = in1,in2
356 if (jn /= in .and. jn /= 0) f1 = f1 * (0 - jn)
357 if (jn /= in) f2 = f2 * (in - jn)
358 end do
359 dgdm(in) = f1 / f2
360 end if
361 end do
362
363 ! Find dr/dmesh
364 drdm = 0
365 do in = in1,in2
366 drdm = drdm + rmesh(ir+in) * dgdm(in)
367 end do
368
369 ! Find differential of volume. Use trapezoidal integration rule
370 dvol = 4 * m_pi * rmesh(ir)**2 * drdm
371
372 ! DVMIN is a small number added to avoid a division by zero
373 dvol = safe_tol(dvol, dvmin)
374 if (ir == 1 .or. ir == nr) dvol = dvol / 2
375 if (gga) aux(ir) = dvol
376
377 ! Find the weights for the derivative d(gradF(i))/d(F(j)), of
378 ! the gradient at point i with respect to the value at point j
379 if (gga) then
380 do in = in1,in2
381 dgidfj(in) = dgdm(in) / drdm
382 end do
383 end if
384
385 ! Find density and gradient of density at this point
386 do is = 1,nspin
387 d(is) = dens(ir,is)
388 end do
389 if (gga) then
390 do is = 1,nspin
391 gd(:,is) = m_zero
392 do in = in1,in2
393 gd(3,is) = gd(3,is) + dgidfj(in) * dens(ir+in,is)
394 end do
395 end do
396 else
397 gd = m_zero
398 end if
399
400 ! The xc_f03_gga routines need as input these combinations of
401 ! the gradient of the densities.
402 sigma(1) = gd(3, 1)*gd(3, 1)
403 if (nspin > 1) then
404 sigma(2) = gd(3, 1) * gd(3, 2)
405 sigma(3) = gd(3, 2) * gd(3, 2)
406 end if
407
408 ! Find exchange and correlation energy densities and their
409 ! derivatives with respect to density and density gradient
410 if (gga) then
411 call xc_f03_gga_exc_vxc(x_conf, int(1, c_size_t), d, sigma, epsx, dexdd, vxsigma)
412 call xc_f03_gga_exc_vxc(c_conf, int(1, c_size_t), d, sigma, epsc, decdd, vcsigma)
413 else
414 call xc_f03_lda_exc_vxc(x_conf, int(1, c_size_t), d, epsx, dexdd)
415 call xc_f03_lda_exc_vxc(c_conf, int(1, c_size_t), d, epsc, decdd)
416 end if
417
418 ! The derivatives of the exchange and correlation energies per particle with
419 ! respect to the density gradient have to be recovered from the derivatives
420 ! with respect to the sigma values defined above.
421 if (gga) then
422 dexdgd(:, 1) = m_two * vxsigma(1) * gd(:, 1)
423 decdgd(:, 1) = m_two * vcsigma(1) * gd(:, 1)
424 if (nspin == 2) then
425 dexdgd(:, 1) = dexdgd(:, 1) + vxsigma(2)*gd(:, 2)
426 dexdgd(:, 2) = m_two * vxsigma(3) * gd(:, 2) + vxsigma(2) * gd(:, 1)
427 decdgd(:, 1) = decdgd(:, 1) + vcsigma(2)*gd(:, 2)
428 decdgd(:, 2) = m_two * vcsigma(3) * gd(:, 2) + vcsigma(2) * gd(:, 1)
429 end if
430 end if
431
432 ! Add contributions to exchange-correlation energy and its
433 ! derivatives with respect to density at all points
434 do is = 1, nspin
435 ex = ex + dvol * d(is) * epsx(1)
436 ec = ec + dvol * d(is) * epsc(1)
437 dx = dx + dvol * d(is) * (epsx(1) - dexdd(is))
438 dc = dc + dvol * d(is) * (epsc(1) - decdd(is))
439 if (gga) then
440 v_xc(ir,is) = v_xc(ir,is) + dvol * (dexdd(is) + decdd(is))
441 do in = in1,in2
442 dx = dx - dvol * dens(ir+in,is) * dexdgd(3,is) * dgidfj(in)
443 dc = dc - dvol * dens(ir+in,is) * decdgd(3,is) * dgidfj(in)
444 v_xc(ir+in,is) = v_xc(ir+in,is) + dvol * &
445 (dexdgd(3,is) + decdgd(3,is)) * dgidfj(in)
446 end do
447 else
448 v_xc(ir,is) = dexdd(is) + decdd(is)
449 end if
450 end do
451
452 end do
453
454!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
455! Divide by volume element to obtain the potential (per electron) !
456!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
457
458 if (gga) then
459 do is = 1,nspin
460 do ir = 1,nr
461 dvol = aux(ir)
462 v_xc(ir,is) = v_xc(ir,is) / dvol
463 end do
464 end do
465 end if
466
467!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
468! Divide by energy unit !
469!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
470
471 ex = ex / eunit
472 ec = ec / eunit
473 dx = dx / eunit
474 dc = dc / eunit
475 do is = 1,nspin
476 do ir = 1,nr
477 v_xc(ir,is) = v_xc(ir,is) / eunit
478 end do
479 end do
480
481 call xc_f03_func_end(x_conf)
482 call xc_f03_func_end(c_conf)
483
484 pop_sub(atomxc)
485 end subroutine atomxc
486
487
488!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
507 subroutine vhrtre(rho, v, r, drdi, srdrdi, nr, a)
508 real(real64), intent(in) :: rho(:)
509 real(real64), intent(out) :: v(:)
510 real(real64), intent(in) :: r(:), drdi(:), srdrdi(:)
511 integer, intent(in) :: nr
512 real(real64), intent(in) :: a
513
514 real(real64) :: x,y,q,a2by4,ybyq,qbyy,qpartc,v0,qt,dz,t,beta,dv
515 integer :: nrm1,nrm2,ir
516
517 push_sub(vhrtre)
518
519 nrm1 = nr - 1
520 nrm2 = nr - 2
521 a2by4 = a*a/m_four
522 ybyq = m_one - a*a/48.0_real64
523 qbyy = m_one/ybyq
524
525!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
526! SIMPSON`S RULE IS USED TO PERFORM TWO INTEGRALS OVER THE ELECTRON !
527! DENSITY. THE TOTAL CHARGE QT IS USED TO FIX THE POTENTIAL AT R=R(NR) !
528! AND V0 (THE INTEGRAL OF THE ELECTRON DENSITY DIVIDED BY R) FIXES !
529! THE ELECTROSTATIC POTENTIAL AT THE ORIGIN !
530! Modified to be consistent with pseudopotential generation (use the !
531! trapezoidal rule for integration). A. Rubio. Jan. 2000 !
532!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
533
534 v0 = m_zero
535 qt = m_zero
536 do ir = 2, nrm1, 2
537 dz = drdi(ir)*rho(ir)
538 qt = qt + dz
539 v0 = v0 + dz/r(ir)
540 end do
541
542 do ir=3, nrm2, 2
543 dz = drdi(ir)*rho(ir)
544 qt = qt + dz
545 v0 = v0 + dz/r(ir)
546 end do
547 dz = drdi(nr)*rho(nr)
548
549 qt = (qt + qt + dz)*m_half
550 v0 = (v0 + v0 + dz/r(nr))
551 v(1) = m_two*v0
552
553!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
554! THE ELECTROSTATIC POTENTIAL AT R=0 IS SET EQUAL TO !
555! THE AVERAGE VALUE OF RHO(R)/R !
556! BEGIN CONSTRUCTION OF THE POTENTIAL AT FINITE !
557!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
558
559 ir = 2
560 t = srdrdi(ir)/r(ir)
561 beta = drdi(ir)*t*rho(ir)
562
563!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
564! THE NEXT 4 STATEMENTS INDICATE THAT WE FIRST FIND THE PARTICULAR !
565! SOLUTION TO THE INHOMOGENEOUS EQN. FOR WHICH Q(2)=0, WE THEN !
566! ADD TO THIS PARTICULAR SOLUTION A SOLUTION OF THE HOMOGENEOUS EQN. !
567! (A CONSTANT IN V OR A Q PROPORTIONAL TO R) !
568! WHICH WHEN DIVIDED BY R IN GOING FROM Q TO V GIVES !
569! THE POTENTIAL THE DESIRED COULOMB TAIL OUTSIDE THE ELECTRON DENSITY. !
570! THE SIGNIFICANCE OF THE SOLUTION VANISHING AT THE SECOND RADIAL !
571! MESH POINT IS THAT, SINCE ALL REGULAR SOLUTIONS OF THE EQUATION !
572! FOR Q=R*V VANISH AT THE ORIGIN, THE KNOWLEDGE OF THE SOLUTION !
573! VALUE AT THE SECOND MESH POINT PROVIDES THE TWO SOLUTION VALUES !
574! REQUIRED TO START THE NUMEROV PROCEDURE. !
575!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
576
577 x = m_zero
578 y = m_zero
579 q = (y - beta / 12.0_real64) * qbyy
580 v(ir) = m_two * t * q
581
582!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
583! BEGINNING OF THE NUMEROV ALGORITHM !
584!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
585
586 do ir = 2, nr
587 x = x + a2by4 * q - beta
588 y = y + x
589 t = srdrdi(ir) / r(ir)
590 beta = t * drdi(ir) * rho(ir)
591 q = (y - beta / 12.0_real64) * qbyy
592 v(ir) = m_two * t * q
593 end do
594
595!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
596! END OF THE NUMEROV ALGORITHM !
597! !
598! WE HAVE NOW FOUND A PARTICULAR SOLUTION TO THE INHOMOGENEOUS EQN. !
599! FOR WHICH Q(R) AT THE SECOND RADIAL MESH POINT EQUALS ZERO. !
600! NOTE THAT ALL REGULAR SOLUTIONS TO THE EQUATION FOR Q=R*V !
601! VANISH AT THE ORIGIN. !
602!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
603
604 qpartc = r(nr)*v(nr)/m_two
605 dz = qt - qpartc
606 dv = m_two*dz/r(nr)
607
608!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
609! THE LOOP FOLLOWING ADDS THE CONSTANT SOLUTION OF THE HOMOGENEOUS !
610! EQN TO THE PARTICULAR SOLUTION OF THE INHOMOGENEOUS EQN. !
611! NOTE THAT V(1) IS CONSTRUCTED INDEPENDENTLY !
612!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
613
614 do ir = 2, nr
615 v(ir) = v(ir) + dv
616 end do
617
618 pop_sub(vhrtre)
619 end subroutine vhrtre
620
621!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
622! FROM NOW, ON, ALL SUBROUTINES ARE IN DOUBLE PRECISION, NO MATTER IF THE CODE
623! IS COMPILED WITH SINGLE PRECISION OR NOT. APPARENTLY DOUBLE PRECISION IS
624! NEEDED FOR THESE PROCEDURES.
625!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
626
627!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
645 subroutine egofv(namespace, h, s, n, e, g, y, l, z, a, b, rmax, nprin, nnode, dr, ierr)
646 type(namespace_t), intent(in) :: namespace
647 real(real64), intent(in) :: h(:), s(:)
648 integer, intent(in) :: n
649 real(real64), intent(inout) :: e
650 real(real64), intent(out) :: g(:), y(:)
651 integer, intent(in) :: l
652 real(real64), intent(in) :: z, a, b, rmax
653 integer, intent(in) :: nprin, nnode
654 real(real64), intent(in) :: dr
655 integer, intent(out) :: ierr
656
657 integer :: i,ncor,n1,n2,niter,nt
658 real(real64) :: e1, e2, del, de, et, t
659 real(real64), parameter :: tol = 1.0e-5_real64
660
661 push_sub(egofv)
662
663 ncor = nprin-l-1
664 n1 = nnode
665 n2 = nnode-1
666 e1 = e
667 e2 = e
668 ierr = 1
669
670!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
671! !
672! the labels 1 and 2 refer to the bisection process, defining the !
673! range in which the desired solution is located. the initial !
674! settings of n1, n2, e1 and e2 are not consistent with the bisection !
675! algorithm; they are set to consistent values when the desired !
676! energy interval has been located. !
677! !
678!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
679
680 nt = 0
681 del = m_half
682 de = m_zero
683
684 do niter = 0, 40
685
686 et = e + de
687 ! the following line is the fundamental "bisection"
688 e = m_half*(e1+e2)
689
690!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
691! !
692! the following concatenation of logical ors ensures that node !
693! counting is used unless the previous integration of the radial !
694! eq produced both the correct number of nodes and a sensible !
695! prediction for the energy. !
696! !
697! sensible means that et must be greater than e1 and less than e2 !
698! correct number of nodes means that nt = nnode or nnode-1. !
699! !
700! leaving e set to its bisection value, and transfering to !
701! the call to yofe means that we are performing bisection, !
702! whereas setting e to et is use of the variational estimate. !
703! !
704!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
705
706
707 if (.not. (et <= e1 .or. et >= e2 .or. nt < nnode-1 .or. nt > nnode)) then
708 e = et
709 if (abs(de) < tol) then
710 ierr = 0
711 exit
712 end if
713 end if
714 call yofe(e,de,dr,rmax,h,s,y,n,l,ncor,nt,z,a,b)
715
716
717!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
718! !
719! yofe integrates the schro eq.; now the bisection logic !
720! !
721!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
722
723 if (nt < nnode) then
724 ! too few nodes; set e1 and n1
725 e1 = e
726 n1 = nt
727
728!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
729! !
730! at this point, we have just set the bottom of the bisection range; !
731! if the top is also set, we proceed. if the top of the range has not !
732! been set, it means that we have yet to find an e greater than the !
733! desired energy. the upper end of the range is extended. !
734! !
735!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
736
737 if (n2 >= nnode) cycle
738 del=del*m_two
739 e2=e1+del
740 cycle
741 end if
742
743 ! too many nodes; set e2 and n2
744 e2=e
745 n2=nt
746
747
748!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
749! !
750! at this point, we have just set the top of the bisection range; !
751! if the top is also set, we proceed. if the top of the range has !
752! not been set, it means that we have yet to find an e less than the !
753! desired energy. the lower end of the range is extended. !
754! !
755!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
756
757 if (n1 < nnode) cycle
758 del=del*m_two
759 e1=e2-del
760
761 end do
762
763!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
764! !
765! the Numerov method uses a transformation of the radial wavefunction. !
766! that we call "y". having located the eigenenergy, we transform !
767! y to "g", from which the density is easily constructed. !
768! finally, the call to "nrmlzg" normalizes g to one electron. !
769! !
770!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
771
772 if (ierr /= 1) then
773 g(1) = m_zero
774 do i = 2, n
775 t = h(i) - e * s(i)
776 g(i) = y(i) / (m_one - t / 12._real64)
777 end do
778 call nrmlzg(namespace, g, s, n)
779 end if
780
781 pop_sub(egofv)
782 end subroutine egofv
783
784
785 subroutine yofe(e,de,dr,rmax,h,s,y,nmax,l,ncor,nnode,z,a,b)
786
787
788!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
805
806 real(real64), intent(inout) :: e
807 real(real64), intent(out) :: de
808 real(real64), intent(in) :: dr, rmax
809 real(real64), intent(in) :: h(:), s(:)
810 real(real64), intent(out) :: y(:)
811 integer, intent(in) :: nmax, l
812 integer, intent(in) :: ncor
813 integer, intent(out) :: nnode
814 real(real64), intent(in) :: z, a, b
815
816 real(real64) :: y2, g, gsg, x, gin, gsgin, xin
817 integer :: n,knk,nndin,i
818 real(real64) :: zdr, yn, ratio, t
819
820 ! No PUSH SUB, called too often.
821
822 zdr = z*a*b
823
824 do n = nmax, 1, -1
825 if (h(n)-e*s(n) < m_one) then
826 knk = n
827 exit
828 end if
829 y(n) = m_zero
830 end do
831
832 call bcorgn(e,h,s,l,zdr,y2)
833
834
835!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
836! !
837! bcorgn computes y2, which embodies the boundary condition !
838! satisfied by the radial wavefunction at the origin !
839!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
840
841 call numout(e,h,s,y,ncor,knk,nnode,y2,g,gsg,x)
842
843!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
844! !
845! the outward integration is now complete !
846! !
847! we first decide if the kinetic energy is sufficiently non !
848! negative to permit use of the numerov eq at rmax. if !
849! it is not, then zero-value boundary condition is used !
850! !
851!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
852
853
854 yn = m_zero
855 if (.not. (n < nmax .or. abs(dr) > 1.e3_real64)) then
856 call bcrmax(e,dr,rmax,h,s,n,yn,a,b)
857 end if
858
859
860!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
861! !
862! bcrmax computes yn, which embodies the boundary condition !
863! satisfied by the radial wavefunction at rmax !
864! !
865!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
866
867 ! n should be larger than 1
868 n = max(n, 2)
869 knk = max(knk, 1)
870 call numin(e,h,s,y,n,nndin,yn,gin,gsgin,xin,knk)
871
872
873!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
874! !
875! numin performs the inward integration by the Numerov method !
876! !
877! the energy increment is now evaluated from the kink in psi !
878! !
879!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
881
882 ratio = g / gin
883 xin = xin * ratio
884 gsg = gsg + gsgin * ratio * ratio
885 t = h(knk) - e * s(knk)
886 de = g * (x + xin + t * g) / gsg
887 nnode = nnode + nndin
888 if (de < m_zero) nnode = nnode + 1
889
890 do i = knk, n
891 y(i) = y(i) * ratio
892 end do
893
894 end subroutine yofe
895
896
897 subroutine nrmlzg(namespace, g, s, n)
898
899
900!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
919 type(namespace_t), intent(in) :: namespace
920 real(real64), intent(inout) :: g(:)
921 real(real64), intent(in) :: s(:)
922 integer, intent(in) :: n
923
924 integer :: nm1,nm2,i
925 real(real64) :: norm, srnrm
926
927 push_sub(nrmlzg)
928
929!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
930! determine the norm of g(i) using Simpson`s rule !
931!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
932
933
934 if (mod(n,2) /= 1) then
935 write(message(1),'(a,i6)') ' nrmlzg: n should be odd. n =', n
936 call messages_warning(1, namespace=namespace)
937 end if
938
939 norm = m_zero
940 nm1 = n - 1
941 do i = 2,nm1,2
942 norm=norm + g(i)*s(i)*g(i)
943 end do
944 norm = norm * m_two
945 nm2 = n - 2
946 do i = 3,nm2,2
947 norm=norm + g(i)*s(i)*g(i)
948 end do
949 norm = norm * m_two
950 nm2 = n - 2
951 do i = 1, n, nm1
952 norm = norm + g(i) * s(i) * g(i)
953 end do
954 norm = norm / m_three
955 srnrm = sqrt(norm)
956 do i = 1, n
957 g(i) = g(i) / srnrm
958 end do
959
960 pop_sub(nrmlzg)
961 end subroutine nrmlzg
962
963
964 subroutine bcorgn(e,h,s,l,zdr,y2)
965 real(real64), intent(in) :: e
966 real(real64), intent(in) :: h(:), s(:)
967 integer, intent(in) :: l
968 real(real64), intent(in) :: zdr
969 real(real64), intent(out) :: y2
970
971 real(real64) :: t2, t3, d2, c0, c1, c2
972
973 ! no PUSH SUB, called too often
974
975!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
976! the quantity called d(i) in the program is actually the inverse !
977! of the diagonal of the tri-diagonal Numerov matrix !
978!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
979
980 t2 = h(2) - e * s(2)
981 d2 = -((24._real64 + 10.0_real64 * t2) / (12._real64 - t2))
982
983!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
984! The following section deals with the fact that the independent !
985! variable "y" in the Numerov equation is not zero at the origin !
986! for l less than 2. !
987! The l=0 solution g vanishes, but the first and second !
988! derivatives are finite, making the Numerov variable y finite. !
989! The l=1 solution g vanishes, and gprime also vanishes, but !
990! The second derivative gpp is finite making y finite. For l > 1, !
991! g and its first two derivatives vanish, making y zero. !
992!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
993
994 if (l < 2) then
995 if (l <= 0) then
996 c0=zdr/6.0_real64
997 c0=c0/(m_one-0.75_real64*zdr)
998 else
999 c0=m_one/12._real64
1000 c0=(-c0)*8.0_real64/m_three
1001 end if
1002
1003 c1=c0*(12._real64+13._real64*t2)/(12._real64-t2)
1004 t3=h(3)-e*s(3)
1005 c2=(-m_half)*c0*(24._real64-t3)/(12._real64-t3)
1006 d2=(d2-c1)/(m_one-c2)
1007 end if
1008 y2=(-m_one)/d2
1009
1010 if (l < 2) then
1011 if (l <= 0) then
1012 c0=zdr/6.0_real64
1013 c0=c0/(m_one-0.75_real64*zdr)
1014 else
1015 c0=m_one/12._real64
1016 c0=(-c0)*8.0_real64/m_three
1017 end if
1018 c1=c0*(12._real64+13._real64*t2)/(12._real64-t2)
1019 t3=h(3)-e*s(3)
1020 c2=(-m_half)*c0*(24._real64-t3)/(12._real64-t3)
1021 d2=(d2-c1)/(m_one-c2)
1022 end if
1023 y2=(-m_one)/d2
1024
1025 end subroutine bcorgn
1026
1027
1028!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1029! 22.7.85 !
1030!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1031 subroutine bcrmax(e, dr, rmax, h, s, n, yn, a, b)
1032 real(real64), intent(in) :: e, dr, rmax
1033 real(real64), intent(in) :: h(:), s(:)
1034 integer, intent(in) :: n
1035 real(real64), intent(out) :: yn
1036 real(real64), intent(in) :: a, b
1037
1038 real(real64) :: tnm1, tn, tnp1, beta, dg, c1, c2, c3, dn
1039
1040 push_sub(bcrmax)
1041
1042 tnm1=h(n-1)-e*s(n-1)
1043 tn =h(n )-e*s(n )
1044 tnp1=h(n+1)-e*s(n+1)
1045 beta=m_one+b/rmax
1046 dg=a*beta*(dr+m_one-m_half/beta)
1047
1048 c2=24._real64*dg/(12._real64-tn)
1049 dn=-((24._real64+10._real64*tn)/(12._real64-tn))
1050
1051 c1= (m_one-tnm1/6.0_real64)/(m_one-tnm1/12._real64)
1052 c3=-((m_one-tnp1/6.0_real64)/(m_one-tnp1/12._real64))
1053 yn=-((m_one-c1/c3)/(dn-c2/c3))
1054
1055 pop_sub(bcrmax)
1056 end subroutine bcrmax
1057
1058
1059 subroutine numin(e,h,s,y,n,nnode,yn,g,gsg,x,knk)
1060 real(real64), intent(in) :: e, h(:), s(:)
1061 real(real64), intent(inout) :: y(:)
1062 integer, intent(in) :: n
1063 integer, intent(out) :: nnode
1064 real(real64), intent(in) :: yn
1065 real(real64), intent(out) :: g, gsg, x
1066 integer, intent(in) :: knk
1067
1068 integer :: i
1069 real(real64) :: t
1070
1071 ! no PUSH SUB, called too often
1072
1073 y(n)=yn
1074 t=h(n)-e*s(n)
1075 g=y(n)/(m_one-t/12._real64)
1076 gsg=g*s(n)*g
1077 i=n-1
1078 y(i)=m_one
1079 t=h(i)-e*s(i)
1080 g=y(i)/(m_one-t/12._real64)
1081 gsg=gsg+g*s(i)*g
1082 x=y(i)-y(n)
1083 nnode=0
1084
1085
1086!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1087! !
1088! begin the inward integration by the Numerov method !
1089! !
1090!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1091
1092 assert(knk > 0)
1093
1094 do i = n - 2, knk, -1
1095 x = x + t * g
1096 y(i) = y(i+1) + x
1097 ! y(i) * y(i+1) can lead to overflows. This is now prevented
1098 ! by comparing only the signs
1099 if (sign(m_one,y(i)) * sign(m_one,y(i+1)) < m_zero) nnode = nnode + 1
1100 t = h(i) - e * s(i)
1101 g = y(i) / (m_one - t / 12._real64)
1102 gsg = gsg + g * s(i) * g
1103 end do
1104
1105
1106!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1107! !
1108! The kink radius is defined as the point where !
1109! psi first turns downward. this usually means at the outermost !
1110! maximum !
1111! !
1112! the inward integration is now complete !
1113! !
1114!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1115
1116 end subroutine numin
1117
1118
1119 subroutine numout(e, h, s, y, ncor, knk, nnode, y2, g, gsg, x)
1120 real(real64), intent(in) :: e, h(:), s(:)
1121 real(real64), intent(inout) :: y(:)
1122 integer, intent(in) :: ncor
1123 integer, intent(inout) :: knk
1124 integer, intent(out) :: nnode
1125 real(real64), intent(in) :: y2
1126 real(real64), intent(out) :: g, gsg, x
1127
1128 integer :: i, nm4
1129 real(real64) :: t, xl
1130
1131 ! no PUSH SUB, called too often
1132
1133 y(1)=m_zero
1134 y(2)=y2
1135 t=h(2)-e*s(2)
1136 g=y(2)/(m_one-t/12._real64)
1137 gsg=g*s(2)*g
1138 y(3)=m_one
1139 t=h(3)-e*s(3)
1140 g=y(3)/(m_one-t/12._real64)
1141 gsg=gsg+g*s(3)*g
1142 x=y(3)-y(2)
1143 nnode=0
1144
1145!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1146! !
1147! begin the outward integration by the Numerov method !
1148! !
1149!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1150
1151 nm4 = knk - 4
1152 knk = nm4
1153 do i = 4, nm4
1154 xl = x
1155 x = x + t * g
1156 y(i) = y(i-1) + x
1157 ! y(i) * y(i-1) can lead to overflows. This is now prevented
1158 ! by comparing only the signs
1159 if (sign(m_one,y(i)) * sign(m_one,y(i-1)) < m_zero) nnode = nnode + 1
1160 t = h(i) - e * s(i)
1161 g = y(i) / (m_one - t / 12._real64)
1162 gsg = gsg + g * s(i) * g
1163 if (.not. (nnode < ncor .or. xl*x > m_zero)) then
1164 knk = i
1165 exit
1166 end if
1167 end do
1168
1169!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1170! !
1171! the outward integration is now complete !
1172! !
1173!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1174 end subroutine numout
1175
1176end module atomic_oct_m
1177
1178!! Local Variables:
1179!! mode: f90
1180!! coding: utf-8
1181!! End:
This is the common interface to a sorting routine. It performs the shell algorithm,...
Definition: sort.F90:156
double sqrt(double __x) __attribute__((__nothrow__
subroutine bcrmax(e, dr, rmax, h, s, n, yn, a, b)
Definition: atomic.F90:1127
subroutine, public atomhxc(namespace, functl, g, nspin, dens, v, extra)
Definition: atomic.F90:229
subroutine nrmlzg(namespace, g, s, n)
Definition: atomic.F90:993
subroutine yofe(e, de, dr, rmax, h, s, y, nmax, l, ncor, nnode, z, a, b)
Definition: atomic.F90:881
subroutine, public atomxc(namespace, FUNCTL, AUTHOR, NR, MAXR, RMESH, NSPIN, DENS, EX, EC, DX, DC, V_XC)
Finds total exchange-correlation energy and potential for a ! spherical electron density distribution...
Definition: atomic.F90:338
subroutine, public valconf_unpolarized_to_polarized(conf)
Definition: atomic.F90:209
integer, parameter, public valconf_string_length
Following stuff is for the "val_conf_t" data type, made to handle atomic configurations.
Definition: atomic.F90:146
subroutine, public valconf_copy(cout, cin)
Definition: atomic.F90:165
subroutine, public valconf_sort_by_l(c, lmax)
Definition: atomic.F90:183
subroutine numin(e, h, s, y, n, nnode, yn, g, gsg, x, knk)
Definition: atomic.F90:1155
subroutine numout(e, h, s, y, ncor, knk, nnode, y2, g, gsg, x)
Definition: atomic.F90:1215
subroutine bcorgn(e, h, s, l, zdr, y2)
Definition: atomic.F90:1060
subroutine, public vhrtre(rho, v, r, drdi, srdrdi, nr, a)
VHRTRE CONSTRUCTS THE ELECTROSTATIC POTENTIAL DUE TO A SUPPLIED ! ELECTRON DENSITY....
Definition: atomic.F90:603
subroutine, public egofv(namespace, h, s, n, e, g, y, l, z, a, b, rmax, nprin, nnode, dr, ierr)
egofv determines the eigenenergy and wavefunction corresponding ! to a particular l,...
Definition: atomic.F90:741
real(real64), parameter, public m_two
Definition: global.F90:193
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_epsilon
Definition: global.F90:207
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
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:524
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
This module is intended to contain "only mathematical" functions and procedures.
Definition: sort.F90:119
int true(void)