30 use,
intrinsic :: iso_fortran_env
84 real(real64),
allocatable :: center(:)
85 real(real64) :: radius =
m_zero
86 class(box_t),
pointer :: box => null()
88 integer,
allocatable :: map(:)
89 integer :: num_regions
90 integer,
allocatable :: regions(:)
91 type(accel_mem_t) :: buff_map
92 real(real64),
allocatable :: rel_x(:,:)
93 real(real64),
allocatable :: r(:)
94 type(mesh_t),
pointer :: mesh => null()
97 integer :: np_global = -1
98 real(real64),
allocatable :: rel_x_global(:,:)
99 integer,
allocatable :: part_v(:)
100 integer,
allocatable :: global2local(:)
102 type(mesh_cube_map_t) :: cube_map
119 recursive real(real64) function f_n(dims)
result(fn)
124 else if (dims == 1)
then
134 type(submesh_t),
intent(inout) :: this
135 class(space_t),
intent(in) :: space
136 class(mesh_t),
target,
intent(in) :: mesh
137 class(box_t),
target,
intent(in) :: box
138 real(real64),
intent(in) :: center(:)
140 integer :: original_mesh_index, submesh_index
141 logical,
allocatable :: submesh_points(:)
142 real(real64),
allocatable :: xtmp(:, :), rtmp(:)
146 assert(.not. space%is_periodic())
154 safe_allocate(submesh_points(1:mesh%np))
155 submesh_points = box%contains_points(mesh%np, mesh%x)
156 this%np = count(submesh_points)
158 assert(this%np <= mesh%np)
160 safe_allocate(this%map(1:this%np))
161 safe_allocate(xtmp(1:space%dim, 1:this%np))
162 safe_allocate(rtmp(1:this%np))
168 do original_mesh_index = 1, mesh%np
170 if (submesh_points(original_mesh_index))
then
171 submesh_index = submesh_index + 1
172 this%map(submesh_index) = original_mesh_index
174 xtmp(:, submesh_index) = mesh%x(original_mesh_index, :) - center(:)
175 rtmp(submesh_index) = norm2(xtmp(:, submesh_index))
181 safe_deallocate_a(xtmp)
182 safe_deallocate_a(rtmp)
183 safe_deallocate_a(submesh_points)
189 subroutine submesh_init(this, space, mesh, latt, center, rc)
192 class(
mesh_t),
target,
intent(in) :: mesh
194 real(real64),
intent(in) :: center(1:space%dim)
195 real(real64),
intent(in) :: rc
197 real(real64) :: r2, rc2, xx(space%dim), rc_norm_n
198 real(real64),
allocatable :: center_copies(:,:), xtmp(:, :), rtmp(:), x(:,:)
199 integer :: icell, is, ip, ix, iy, iz
200 integer(int64) :: max_elements_count
201 type(lattice_iterator_t) :: latt_iter
202 integer,
allocatable :: map_inv(:), map_temp(:)
203 integer :: nmax(3), nmin(3)
204 real(real64),
parameter :: tol = 1e-13_real64
210 assert(space%dim <= 3)
214 safe_allocate(this%center(1:space%dim))
215 this%center(:) = center(:)
223 if (.not. space%is_periodic())
then
226 safe_allocate(map_inv(0:this%mesh%np))
227 map_inv(0:this%mesh%np) = 0
233 nmin(1:space%dim) = int((center(1:space%dim) - abs(rc))/mesh%spacing(1:space%dim)) - 1
234 nmax(1:space%dim) = int((center(1:space%dim) + abs(rc))/mesh%spacing(1:space%dim)) + 1
238 nmin(1:space%dim) = max(mesh%idx%nr(1, 1:space%dim), nmin(1:space%dim))
239 nmax(1:space%dim) = min(mesh%idx%nr(2, 1:space%dim), nmax(1:space%dim))
243 do iz = nmin(3), nmax(3)
244 do iy = nmin(2), nmax(2)
245 do ix = nmin(1), nmax(1)
247 if (ip == 0 .or. ip > mesh%np) cycle
248 r2 = sum((mesh%x(ip, :) - center)**2)
260 safe_allocate(this%map(1:this%np))
261 safe_allocate(xtmp(1:space%dim, 1:this%np))
262 safe_allocate(rtmp(1:this%np))
265 do iz = nmin(3), nmax(3)
266 do iy = nmin(2), nmax(2)
267 do ix = nmin(1), nmax(1)
269 if (ip == 0 .or. ip > mesh%np) cycle
273 xtmp(:, is) = mesh%x(ip, :) - center
274 rtmp(is) = norm2(xtmp(:,is))
279 safe_deallocate_a(map_inv)
292 safe_allocate(center_copies(1:space%dim, 1:latt_iter%n_cells))
293 do icell = 1, latt_iter%n_cells
294 center_copies(:, icell) = center + latt_iter%get(icell)
299 rc_norm_n = product(ceiling(rc / mesh%spacing(1:space%dim), int64) +
m_one)
300 if (mesh%use_curvilinear) rc_norm_n = rc_norm_n / mesh%coord_system%min_mesh_scaling_product
301 max_elements_count = 3**space%dim * int(
m_pi**
floor(0.5 * space%dim) * rc_norm_n *
f_n(space%dim), int64)
303 safe_allocate(x(1:space%dim,1:mesh%np))
305 x(:,ip) = mesh%x(ip,:)
308 safe_allocate(map_temp(1:max_elements_count))
309 safe_allocate(xtmp(1:space%dim, 1:max_elements_count))
310 safe_allocate(rtmp(1:max_elements_count))
314 do icell = 1, latt_iter%n_cells
315 xx = x(:,ip) - center_copies(:, icell)
316 if(any(abs(xx)>rc+tol)) cycle
326 assert(is < huge(is))
329 safe_allocate(this%map(1:this%np))
330 this%map(1:this%np) = map_temp(1:this%np)
333 safe_deallocate_a(map_temp)
334 safe_deallocate_a(center_copies)
342 safe_deallocate_a(xtmp)
343 safe_deallocate_a(rtmp)
350 class(
space_t),
intent(in) :: space
351 real(real64),
intent(in) :: xtmp(:, :), rtmp(:)
353 integer :: ip, i_region, offset
354 integer,
allocatable :: order(:), order_new(:)
355 integer,
allocatable :: map_new(:)
356 integer,
allocatable :: np_region(:), tmp_array(:)
363 safe_allocate(order(1:this%np))
364 safe_allocate(this%rel_x(1:space%dim, 1:this%np))
365 safe_allocate(this%r(1:this%np))
372 call sort(this%map, order)
375 this%overlap = .false.
376 do ip = 1, this%np - 1
377 if (this%map(ip) == this%map(ip + 1))
then
379 this%overlap = .
true.
387 if(this%overlap)
then
391 safe_allocate(tmp_array(1:this%np))
392 safe_allocate(order_new(1:this%np))
393 safe_allocate(np_region(1:this%np))
394 safe_allocate(map_new( 1:this%np))
401 if (this%map(ip) == this%map(ip - 1))
then
402 i_region = i_region + 1
403 if (i_region > this%num_regions)
then
404 this%num_regions = i_region
405 np_region(i_region) = 0
410 tmp_array(ip) = i_region
411 np_region(i_region) = np_region(i_region) + 1
414 assert( .not.
allocated(this%regions))
417 safe_allocate(this%regions(1:this%num_regions+1))
421 if(this%num_regions > 1)
then
422 do i_region = 1, this%num_regions
423 this%regions(i_region + 1) = this%regions(i_region) + np_region(i_region)
426 this%regions(2) = this%np + 1
429 np_region(1:this%np) = 0
430 order_new(1:this%np) = -1
431 map_new(1:this%np) = -1
436 i_region = tmp_array(ip)
437 np_region(i_region) = np_region(i_region) + 1
438 offset = this%regions(i_region) - 1
439 map_new( offset + np_region(i_region) ) = this%map(ip)
440 order_new( offset + np_region(i_region) ) = order(ip)
443 order(1:this%np) = order_new(1:this%np)
444 this%map(1:this%np) = map_new(1:this%np)
446 safe_deallocate_a(tmp_array)
447 safe_deallocate_a(order_new)
448 safe_deallocate_a(np_region)
449 safe_deallocate_a(map_new)
454 safe_allocate(this%regions(1:2))
456 this%regions(2) = this%np + 1
461 this%rel_x(:, ip) = xtmp(:, order(ip))
462 this%r(ip) = rtmp(order(ip))
465 safe_deallocate_a(order)
476 class(
space_t),
intent(in) :: space
477 class(
mesh_t),
target,
intent(in) :: mesh
480 real(real64),
optional,
intent(in) :: shift(:)
484 real(real64) :: xx(space%dim), diff_centers(space%dim)
491 safe_allocate(this%center(1:space%dim))
492 this%center(:) = sm1%center(:)
493 this%radius = sm1%radius
497 this%overlap = sm1%overlap .or. sm2%overlap
499 diff_centers = sm1%center - sm2%center
500 if (
present(shift)) diff_centers = diff_centers - shift
507 xx = sm2%rel_x(:, ip) - diff_centers
509 if (sum(xx**2) > sm1%radius**2) is = is + 1
514 safe_allocate(this%map(1:this%np))
515 safe_allocate(this%rel_x(1:space%dim, 1:this%np))
516 safe_allocate(this%r(1:this%np))
517 this%map(1:sm1%np) = sm1%map(1:sm1%np)
518 this%rel_x(:, 1:sm1%np) = sm1%rel_x(:, 1:sm1%np)
519 this%r(1:sm1%np) = sm1%r(1:sm1%np)
524 xx = sm2%rel_x(:, ip) - diff_centers
526 if (r2 > sm1%radius**2)
then
528 this%map(is) = sm2%map(ip)
529 this%r(is) =
sqrt(r2)
530 this%rel_x(:, is) = xx
542 class(
space_t),
intent(in) :: space
543 real(real64),
intent(in) :: newcenter(:)
545 real(real64) :: xx(space%dim), diff_centers(space%dim), oldcenter(space%dim)
551 oldcenter = this%center
552 this%center(:) = newcenter(:)
554 diff_centers = newcenter - oldcenter
557 xx = this%rel_x(:, ip) - diff_centers
558 this%r(ip) = norm2(xx)
559 this%rel_x(:, ip) = xx
569 class(
space_t),
intent(in) :: space
570 type(
mesh_t),
target,
intent(in) :: mesh
571 real(real64),
intent(in) :: center(1:space%dim)
572 real(real64),
intent(in) :: radius
573 integer,
intent(in) :: root
576 integer :: nparray(1:3)
581 if (root /= mpi_grp%rank)
then
583 safe_allocate(this%center(1:space%dim))
584 this%center(:) = center(:)
588 if (mpi_grp%size > 1)
then
590 if (root == mpi_grp%rank)
then
592 nparray(2) = this%num_regions
593 if (this%overlap)
then
600 call mpi_grp%bcast(nparray, 3, mpi_integer, root)
602 this%num_regions = nparray(2)
603 this%overlap = (nparray(3) == 1)
605 if (root /= mpi_grp%rank)
then
606 safe_allocate(this%map(1:this%np))
607 safe_allocate(this%rel_x(1:space%dim, 1:this%np))
608 safe_allocate(this%r(1:this%np))
609 safe_allocate(this%regions(1:this%num_regions+1))
612 call mpi_grp%bcast(this%regions(1), this%num_regions+1, mpi_integer, root)
614 if (this%np > 0)
then
615 call mpi_grp%bcast(this%map(1), this%np, mpi_integer, root)
616 call mpi_grp%bcast(this%rel_x(1, 1), this%np*space%dim, mpi_double_precision, root)
617 call mpi_grp%bcast(this%r(1), this%np, mpi_double_precision, root)
629 real(real64),
intent(in) :: radius
630 real(real64),
intent(in) :: center(:)
631 real(real64),
intent(in) :: dx
634 if (
allocated(this%center))
then
637 if (radius <= this%radius+dx*1e-6_real64 .and. all(abs(this%center - center) < 0.25*dx))
then
652 if (this%np /= -1)
then
654 safe_deallocate_a(this%center)
655 safe_deallocate_a(this%map)
656 safe_deallocate_a(this%rel_x)
657 safe_deallocate_a(this%r)
658 safe_deallocate_a(this%regions)
672 integer,
intent(out) :: map_inv(:)
678 map_inv(1:this%mesh%np_part) = 0
680 map_inv(this%map(is)) = is
690 class(
space_t),
intent(in) :: space
692 integer :: ii, jj, dd
697 if (.not. space%is_periodic())
then
699 distance = sum((sm1%center - sm2%center)**2)
700 overlap =
distance <= (1.5_real64*(sm1%radius + sm2%radius))**2
703 if (.not. overlap)
return
711 do while(ii <= sm1%np .and. jj <= sm2%np)
712 dd = sm1%map(ii) - sm2%map(jj)
715 else if (dd > 0)
then
723 if (sm1%mesh%parallel_in_domains)
then
724 call sm1%mesh%mpi_grp%allreduce_inplace(overlap, 1, mpi_logical, mpi_lor)
732 class(
space_t),
intent(in) :: space
734 integer,
allocatable :: part_np(:)
735 integer :: ipart, ind, ip
739 if (.not. this%mesh%parallel_in_domains)
then
740 this%np_global = this%np
745 safe_allocate(part_np(this%mesh%pv%npart))
747 part_np(this%mesh%pv%partno) = this%np
749 call this%mesh%allreduce(part_np)
750 this%np_global = sum(part_np)
752 safe_allocate(this%rel_x_global(1:space%dim, 1:this%np_global))
753 safe_allocate(this%part_v(1:this%np_global))
754 safe_allocate(this%global2local(1:this%np_global))
755 this%rel_x_global(1:space%dim, 1:this%np_global) =
m_zero
756 this%part_v(1:this%np_global) = 0
757 this%global2local(1:this%np_global) = 0
760 do ipart = 1, this%mesh%pv%npart
761 if (ipart == this%mesh%pv%partno)
then
763 this%rel_x_global(:, ind + ip) = this%rel_x(:, ip)
764 this%part_v(ind + ip) = this%mesh%pv%partno
765 this%global2local(ind + ip) = ip
768 ind = ind + part_np(ipart)
771 call this%mesh%allreduce(this%rel_x_global)
772 call this%mesh%allreduce(this%part_v)
773 call this%mesh%allreduce(this%global2local)
775 safe_deallocate_a(part_np)
786 safe_deallocate_a(this%rel_x_global)
788 safe_deallocate_a(this%part_v)
789 safe_deallocate_a(this%global2local)
798 complex(real64),
intent(in) :: sphi(:)
799 complex(real64),
intent(inout) :: phi(:)
800 complex(real64),
optional,
intent(in) :: factor
806 if (
present(factor))
then
810 phi(this%map(ip)) = phi(this%map(ip)) + factor*sphi(ip)
812 if (this%np.ge.4)
then
813 do ip = m+1, this%np, 4
814 phi(this%map(ip)) = phi(this%map(ip)) + factor*sphi(ip)
815 phi(this%map(ip+1)) = phi(this%map(ip+1)) + factor*sphi(ip+1)
816 phi(this%map(ip+2)) = phi(this%map(ip+2)) + factor*sphi(ip+2)
817 phi(this%map(ip+3)) = phi(this%map(ip+3)) + factor*sphi(ip+3)
823 phi(this%map(ip)) = phi(this%map(ip)) + sphi(ip)
825 if (this%np.ge.4)
then
826 do ip = m+1, this%np, 4
827 phi(this%map(ip)) = phi(this%map(ip)) + sphi(ip)
828 phi(this%map(ip+1)) = phi(this%map(ip+1)) + sphi(ip+1)
829 phi(this%map(ip+2)) = phi(this%map(ip+2)) + sphi(ip+2)
830 phi(this%map(ip+3)) = phi(this%map(ip+3)) + sphi(ip+3)
841 complex(real64),
intent(in) :: sphi(:)
842 complex(real64),
intent(in) :: phi(:)
843 logical,
optional,
intent(in) :: reduce
851 if (this%mesh%use_curvilinear)
then
853 dotp = dotp + this%mesh%vol_pp(this%map(is))*phi(this%map(is))*conjg(sphi(is))
858 dotp = dotp + phi(this%map(ip))*conjg(sphi(ip))
860 if (this%np.ge.4)
then
861 do ip = m+1, this%np, 4
862 dotp = dotp + phi(this%map(ip))*conjg(sphi(ip)) &
863 + phi(this%map(ip+1))*conjg(sphi(ip+1)) &
864 + phi(this%map(ip+2))*conjg(sphi(ip+2)) &
865 + phi(this%map(ip+3))*conjg(sphi(ip+3))
868 dotp = dotp*this%mesh%vol_pp(1)
873 call this%mesh%allreduce(dotp)
883 type(
submesh_t),
target,
intent(in) :: sm
884 class(
space_t),
intent(in) :: space
885 integer,
intent(out) :: db(1:space%dim)
888 real(real64) :: chi(space%dim), max_chi(space%dim)
889 integer :: db_red(1:space%dim)
890 real(real64),
parameter :: tol=1.0e-10_real64
896 chi = sm%mesh%coord_system%from_cartesian(sm%rel_x(:, ip))
897 do idir = 1, space%dim
898 max_chi(idir) = max(max_chi(idir), abs(chi(idir))/sm%mesh%spacing(idir))
902 do idir = 1, space%dim
903 db(idir) = nint(max_chi(idir)-tol)
906 if(sm%mesh%parallel_in_domains)
then
907 call sm%mesh%mpi_grp%allreduce(db(1), db_red(1), space%dim, mpi_integer, mpi_max)
908 db(1:space%dim) = db_red(1:space%dim)
918 type(
submesh_t),
target,
intent(inout) :: sm
919 class(
space_t),
intent(in) :: space
923 real(real64) :: chi(space%dim), shift(space%dim)
927 sm%cube_map%nmap = sm%np
929 safe_allocate(sm%cube_map%map(1:space%dim, 1:sm%cube_map%nmap))
933 chi = sm%mesh%coord_system%from_cartesian(sm%center)
934 do idir = 1, space%dim
935 shift(idir) = nint(chi(idir)/sm%mesh%spacing(idir))*sm%mesh%spacing(idir)
937 shift = sm%mesh%coord_system%to_cartesian(shift)
938 shift = shift - sm%center
940 do ip = 1, sm%cube_map%nmap
941 chi = sm%mesh%coord_system%from_cartesian(sm%rel_x(:,ip) - shift)
942 do idir = 1, space%dim
943 sm%cube_map%map(idir, ip) = nint(chi(idir)/sm%mesh%spacing(idir))
949 call accel_write_buffer(sm%cube_map%map_buffer, sm%cube_map%nmap*space%dim, sm%cube_map%map)
974 class(
batch_t),
intent(in) :: ss
977 integer :: ist, idim, jdim, is
981 assert(.not. mm%is_packed())
984 assert(ss%nst_linear == mm%nst_linear)
985 assert(ss%status() == mm%status())
986 assert(ss%dim == mm%dim)
988 assert(mm%nst == ss%nst)
993 jdim = min(idim, ss%dim)
996 mm%zff(this%map(is), idim, ist) = &
997 mm%zff(this%map(is), idim, ist) + ss%dff(is, jdim, ist)
1010#include "submesh_inc.F90"
1013#include "complex.F90"
1014#include "submesh_inc.F90"
This is the common interface to a sorting routine. It performs the shell algorithm,...
double sqrt(double __x) __attribute__((__nothrow__
double floor(double __x) __attribute__((__nothrow__
double fn(const gsl_vector *v, void *params)
subroutine, public accel_release_buffer(this)
pure logical function, public accel_is_enabled()
integer, parameter, public accel_mem_read_only
This module implements batches of mesh functions.
Module implementing boundary conditions in Octopus.
real(real64), parameter, public m_two
real(real64), parameter, public m_zero
real(real64), parameter, public m_pi
some mathematical constants
complex(real64), parameter, public m_z0
real(real64), parameter, public m_one
This module implements the index, used for the mesh points.
subroutine, public mesh_cube_map_end(this)
This module defines the meshes, which are used in Octopus.
integer function, public mesh_local_index_from_coords(mesh, ix)
This function returns the local index of the point for a given vector of integer coordinates.
Some general things and nomenclature:
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.
This module is intended to contain "only mathematical" functions and procedures.
real(real64) function, public dsm_integrate_frommesh(mesh, sm, ff, reduce)
subroutine zdsubmesh_add_to_mesh(this, sphi, phi, factor)
subroutine, public submesh_init_cube_map(sm, space)
logical function, public submesh_overlap(sm1, sm2, space)
subroutine, public zsubmesh_batch_dotp_matrix(this, mm, ss, dot, reduce)
real(real64) function, public dsm_nrm2(sm, ff, reduce)
this function returns the the norm of a vector
subroutine zzsubmesh_add_to_mesh(this, sphi, phi, factor)
subroutine, public submesh_end_global(this)
subroutine, public zsubmesh_copy_from_mesh(this, phi, sphi, conjugate)
subroutine submesh_reorder_points(this, space, xtmp, rtmp)
complex(real64) function, public zsm_integrate(mesh, sm, ff, reduce)
subroutine, public submesh_shift_center(this, space, newcenter)
subroutine, public submesh_init_box(this, space, mesh, box, center)
This subroutine creates a submesh which has the same box of the one that the user passes.
complex(real64) function zzsubmesh_to_mesh_dotp(this, sphi, phi, reduce)
real(real64) function, public dsm_integrate(mesh, sm, ff, reduce)
subroutine, public dsubmesh_batch_add_matrix(this, factor, ss, mm)
The following functions takes a batch of functions defined in submesh (ss) and adds all of them to ea...
recursive real(real64) function f_n(dims)
subroutine, public dsubmesh_copy_from_mesh(this, phi, sphi, conjugate)
real(real64) function ddsubmesh_to_mesh_dotp(this, sphi, phi, reduce)
subroutine, public submesh_merge(this, space, mesh, sm1, sm2, shift)
subroutine, public dzsubmesh_batch_add(this, ss, mm)
The following function takes a batch of functions defined in submesh (ss) and adds one of them to eac...
subroutine, public zsubmesh_batch_add(this, ss, mm)
The following function takes a batch of functions defined in submesh (ss) and adds one of them to eac...
real(real64) function, public zsm_nrm2(sm, ff, reduce)
this function returns the the norm of a vector
subroutine, public zsubmesh_batch_add_matrix(this, factor, ss, mm)
The following functions takes a batch of functions defined in submesh (ss) and adds all of them to ea...
subroutine, public submesh_end_cube_map(sm)
subroutine, public submesh_end(this)
subroutine ddsubmesh_add_to_mesh(this, sphi, phi, factor)
subroutine, public submesh_get_cube_dim(sm, space, db)
finds the dimension of a box containing the submesh
subroutine, public dsubmesh_batch_dotp_matrix(this, mm, ss, dot, reduce)
complex(real64) function, public zsm_integrate_frommesh(mesh, sm, ff, reduce)
logical function, public submesh_compatible(this, radius, center, dx)
subroutine, public submesh_build_global(this, space)
subroutine, public submesh_broadcast(this, space, mesh, center, radius, root, mpi_grp)
subroutine, public submesh_get_inv(this, map_inv)
subroutine, public dsubmesh_batch_add(this, ss, mm)
The following function takes a batch of functions defined in submesh (ss) and adds one of them to eac...
subroutine, public submesh_init(this, space, mesh, latt, center, rc)
complex(real64) function zdsubmesh_to_mesh_dotp(this, sphi, phi, reduce)
type(type_t), parameter, public type_cmplx
type(type_t), parameter, public type_integer
type(type_t), parameter, public type_float
Class defining batches of mesh functions.
The following class implements a lattice iterator. It allows one to loop over all cells that are with...
Describes mesh distribution to nodes.
This is defined even when running serial.
A submesh is a type of mesh, used for the projectors in the pseudopotentials It contains points on a ...
void distance(const int iatom, const int jatom, const double coordinates[], double *rr, double *rr2, double *rr6, double *rr7)