28 use,
intrinsic :: iso_fortran_env
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(:)
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(:, :)
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
99 assert(spin_channels == 1 .or. spin_channels == 2)
104 select type (species)
116 do isp = 1, spin_channels
117 do iorb = 1, species%get_niwfs()
118 call species%get_iwf_ilm(iorb, isp, ii, ll, mm)
120 call species%get_iwf_n(iorb, isp, nn)
122 radius = species%get_iwf_radius(nn, isp)
124 radius = max(radius,
m_two*maxval(mesh%spacing))
126 call submesh_init(sphere, space, mesh, latt, pos, radius)
127 safe_allocate(dorbital(1:sphere%np))
134 dorbital(ip) = species%conf%occ(ii, isp)/real(2*ll+1, real64) *dorbital(ip)*dorbital(ip)
137 safe_deallocate_a(dorbital)
151 rmax = latt%max_length()
154 do icell = 1, latt_iter%n_cells
155 yy = latt_iter%get(icell)
157 call mesh_r(mesh, ip, rr, origin = pos, coords = xx)
163 rho(ip, 1) = rho(ip, 1) + rerho
169 if (spin_channels > 1)
then
171 rho(:, 2) = rho(:, 1)
175 do isp = 1, spin_channels
177 rho(1:mesh%np, isp) = x * rho(1:mesh%np, isp)
185 rmax = latt%max_length()
188 do icell = 1, latt_iter%n_cells
189 yy = latt_iter%get(icell)
191 call mesh_r(mesh, ip, rr, origin = pos, coords = xx)
197 rho(ip, 1) = rho(ip, 1) + rerho
201 if (spin_channels > 1)
then
202 rho(:, 1) =
m_half*rho(:, 1)
203 rho(:, 2) = rho(:, 1)
207 do isp = 1, spin_channels
209 rho(1:mesh%np, isp) = x * rho(1:mesh%np, isp)
216 call mesh_r(mesh, ip, rr, origin = pos)
217 if (rr <= species%radius())
then
218 in_points = in_points + 1
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
227 if (in_points > 0)
then
230 if (mesh%use_curvilinear)
then
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) )
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) )
250 density = species%get_density(mesh%box%bounding_box_l) / spin_channels
253 rr = abs(mesh%x(3, ip) - pos(3))
254 if (rr <= species%thickness() /
m_two)
then
255 rho(ip, 1:spin_channels) = density
266 assert(
allocated(ps%density))
269 do isp = 1, spin_channels
270 rmax = max(rmax, ps%density(isp)%x_threshold)
274 do icell = 1, latt_iter%n_cells
275 pos_pc = pos + latt_iter%get(icell)
277 call mesh_r(mesh, ip, rr, origin = pos_pc)
280 do isp = 1, spin_channels
282 rho(ip, isp) = rho(ip, isp) +
spline_eval(ps%density(isp), rr)
293 do icell = 1, latt_iter%n_cells
294 pos_pc = pos + latt_iter%get(icell)
296 call mesh_r(mesh, ip, rr, origin = pos_pc)
301 do isp = 1, spin_channels
310 do isp = 1, spin_channels
314 rho(1:mesh%np, 1:spin_channels) = rho(1:mesh%np, 1:spin_channels)*species%get_zval()/nrm
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)
337 class(
species_t),
target,
intent(in) :: species
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(:, :)
345 real(real64) :: rr, nrm
346 type(
ps_t),
pointer :: ps
360 assert(
allocated(ps%density))
363 do isp = 1, spin_channels
366 call mesh_r(mesh, ip, rr, origin = pos)
369 rho(ip, isp) = rho(ip, isp) +
spline_eval(ps%density(isp), rr)
380 call mesh_r(mesh, ip, rr, origin = pos)
385 do isp = 1, spin_channels
393 do isp = 1, spin_channels
397 rho(1:mesh%np, 1:spin_channels) = rho(1:mesh%np, 1:spin_channels)*species%get_zval()/nrm
413 class(
species_t),
target,
intent(in) :: species
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(:, :)
422 type(
ps_t),
pointer :: ps
434 do isp = 1, spin_channels
437 call mesh_r(mesh, ip, rr, origin = pos)
440 drho(ip, isp) = drho(ip, isp) +
spline_eval(ps%density_der(isp), rr)
465 class(
species_t),
target,
intent(in) :: species
467 class(
space_t),
intent(in) :: space
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(:, :, :)
474 integer :: isp, ip, icell, idir
475 real(real64) :: rr, pos_pc(space%dim), range, spline
476 type(
ps_t),
pointer :: ps
481 assert(spin_channels == 1 .or. spin_channels == 2)
493 range = ps%density_der(1)%x_threshold
494 if (spin_channels == 2) range = max(range, ps%density_der(2)%x_threshold)
497 do icell = 1, latt_iter%n_cells
498 pos_pc = pos + latt_iter%get(icell)
501 call mesh_r(mesh, ip, rr, origin = pos_pc)
504 do isp = 1, spin_channels
508 if(abs(spline) < 1e-150_real64) cycle
510 do idir = 1, space%dim
511 drho(ip, isp, idir) = drho(ip, isp, idir) - spline*(mesh%x(idir, ip) - pos_pc(idir))/rr
535 class(
species_t),
target,
intent(in) :: species
537 class(
space_t),
intent(in) :: space
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(:,:)
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
552 type(
ps_t),
pointer :: ps
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
565 if(
present(nlr_x))
then
566 assert(species%is_ps())
569 select type (species)
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
580 call submesh_init(sphere_local, space, mesh, latt, pos, radius)
581 sphere => sphere_local
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)
590 if(sphere%r(ip) <= radius_nlr)
then
602 norm_factor = abs(species%get_zval()/
dsm_integrate(mesh, sphere, rho_sphere))
604 rho(sphere%map(ip)) = rho(sphere%map(ip)) + norm_factor*rho_sphere(ip)
607 if (
present(nlr_x))
then
608 do idir = 1, space%dim
610 nlr_x(sphere%map(ip), idir) = nlr_x(sphere%map(ip), idir) + norm_factor*rho_sphere(ip)*sphere%rel_x(idir, ip)
615 safe_deallocate_a(rho_sphere)
617 if ( .not.
present(sphere_inout) )
then
628 if (mesh%mpi_grp%rank /= rankmin) have_point = .false.
631 if (mesh%use_curvilinear)
then
632 rho(ipos) = -species%get_z()/mesh%vol_pp(ipos)
634 rho(ipos) = -species%get_z()/mesh%vol_pp(1)
638 write(
message(1),
'(3a,f5.2,3a)') &
639 "Info: species_full_delta species ", trim(species%get_label()), &
647 if (space%is_periodic())
then
653 if (mesh%use_curvilinear)
then
663 safe_allocate(rho_p(1:mesh%np))
664 safe_allocate(grho_p(1:mesh%np, 1:space%dim))
670 delta = mesh%spacing(1)
671 alpha =
sqrt(
m_two)*species%get_sigma()*delta
675 startval(1:space%dim) = pos
683 write(
message(1),
'(a)')
'Root finding in species_get_density did not converge.'
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),
')'
693 rho = -species%get_z()*rho_p
697 safe_deallocate_a(grho_p)
698 safe_deallocate_a(rho_p)
709 range = latt%max_length()
713 do icell = 1, latt_iter%n_cells
714 yy = latt_iter%get(icell)
716 call mesh_r(mesh, ip, rr, origin = pos, coords = xx)
722 rho(ip) = rho(ip) - rerho
730 range = latt%max_length()
734 do icell = 1, latt_iter%n_cells
735 yy = latt_iter%get(icell)
737 call mesh_r(mesh, ip, rr, origin = pos, coords = xx)
743 rho(ip) = rho(ip) - rerho
748 rho(1:mesh%np) = rr * rho(1:mesh%np)
760 subroutine func(xin, ff, jacobian)
761 real(real64),
intent(in) :: xin(:)
762 real(real64),
intent(out) :: ff(:), jacobian(:,:)
764 real(real64),
allocatable :: xrho(:)
765 integer :: idir, jdir, dim, ip
772 safe_allocate(xrho(1:mesh_p%np))
778 xrho(ip) = rho_p(ip) * mesh_p%x_t(ip, idir)
788 xrho(ip) = grho_p(ip, jdir) * mesh_p%x_t(ip, idir)
794 safe_deallocate_a(xrho)
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
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
808 real(real64) :: center(space%dim), rr
811 type(
ps_t),
pointer :: ps
822 do icell = 1, latt_iter%n_cells
823 center = pos + latt_iter%get(icell)
825 rr = norm2(mesh%x(1:space%dim, ip) - center)
827 rho_core(ip) = rho_core(ip) +
spline_eval(ps%core, rr)
841 subroutine getrho(dim, xin)
842 integer,
intent(in) :: dim
843 real(real64),
intent(in) :: xin(1:dim)
846 real(real64) :: r2, chi(dim), norm, threshold
852 threshold = -
log(0.001_real64)*alpha2_p
856 chi(1:dim) = mesh_p%x(1:dim, ip)
857 r2 = sum((chi - xin(1:dim))**2)
859 if (r2 < threshold)
then
860 rho_p(ip) =
exp(-r2/alpha2_p)
866 grho_p(ip, idir) = (chi(idir) - xin(idir)) * rho_p(ip)
881 class(
species_t),
target,
intent(in) :: species
883 class(
space_t),
intent(in) :: space
885 real(real64),
intent(in) :: pos(1:space%dim)
886 type(
mesh_t),
intent(in) :: mesh
887 real(real64),
intent(out) :: vl(:)
889 real(real64) :: a1, a2, Rb2, range, density
890 real(real64) :: xx(space%dim), pos_pc(space%dim), r, r2, threshold
891 integer :: ip, err, icell
892 complex(real64) :: zpot
894 real(real64) :: aa, bb
902 call parse_variable(namespace,
'SpeciesProjectorSphereThreshold', 0.001_real64, threshold)
908 do icell = 1, latt_iter%n_cells
909 pos_pc = pos + latt_iter%get(icell)
911 call mesh_r(mesh, ip, r, origin = pos_pc)
913 vl(ip) = vl(ip) -species%get_zval()/
sqrt(r2+species%get_softening2())
919 range = 5.0_real64 * latt%max_length()
922 do icell = 1, latt_iter%n_cells
923 pos_pc = pos + latt_iter%get(icell)
925 call mesh_r(mesh, ip, r, origin = pos_pc, coords = xx)
927 zpot = species%user_pot(space%dim, xx, r)
928 vl(ip) = vl(ip) + real(zpot, real64)
934 call dio_function_input(trim(species%get_filename()), namespace, space, mesh, vl, err)
936 write(
message(1),
'(a)')
'Error loading file '//trim(species%get_filename())//
'.'
937 write(
message(2),
'(a,i4)')
'Error code returned = ', err
943 assert(.not. space%is_periodic())
945 a1 = species%get_z()/(
m_two*species%radius()**3)
946 a2 = species%get_z()/species%radius()
947 rb2= species%radius()**2
951 xx = mesh%x(:, ip) - pos(1:space%dim)
954 if (r <= species%radius())
then
955 vl(ip) = (a1*(r*r - rb2) - a2)
957 vl(ip) = -species%get_z()/r
966 density = species%get_density(mesh%box%bounding_box_l)
971 r = abs(mesh%x(3, ip) - pos(3))
973 if (r <= species%thickness()/
m_two)
then
974 vl(ip) = a1 * (r * r / species%thickness() + species%thickness() /
m_four)
983 assert(.not. space%is_periodic())
987 r = norm2(mesh%x(:, ip) - pos)
993 if (space%is_periodic())
then
1003 r2 = sum((mesh%x(:, ip) - pos)**2)*(species%get_z()*aa)**2
1014 vl(ip) = vl(ip) * (species%get_z())**2
Copies a vector x, to a vector y.
scales a vector by a constant
Both the filling of the function, and the retrieval of the values may be done using single- or double...
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
real(real64), parameter, public m_two
real(real64), parameter, public r_small
real(real64), parameter, public m_zero
real(real64), parameter, public m_four
real(real64), parameter, public m_pi
some mathematical constants
real(real64), parameter, public m_half
real(real64), parameter, public m_one
real(real64), parameter, public m_three
This module implements the index, used for the mesh points.
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.
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.
pure subroutine, public mesh_r(mesh, ip, rr, origin, coords)
return the distance to the origin for a given grid point
subroutine, public messages_not_implemented(feature, namespace)
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
subroutine, public messages_experimental(name, namespace)
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
pure logical function, public ps_has_density(ps)
real(real64) pure function, public long_range_potential(r, sigma, z_val)
Evaluate the long-range potential at a given distance.
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.
real(real64) function, public spline_eval(spl, x)
real(real64) pure function, public spline_range_max(this)
real(real64) function, public dsm_integrate(mesh, sm, ff, reduce)
subroutine, public submesh_end(this)
subroutine, public submesh_init(this, space, mesh, latt, center, rc)
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
character(len=20) pure function, public units_abbrev(this)
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)
logical function, public volume_in_volume(space, vol, xx)
subroutine, public volume_end(vol)
subroutine, public volume_init(vol)
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.
A type storing the information and data about a pseudopotential.
An abstract class for species. Derived classes include jellium, all electron, and pseudopotential spe...
A submesh is a type of mesh, used for the projectors in the pseudopotentials It contains points on a ...