Octopus
submesh.F90
Go to the documentation of this file.
1!! Copyright (C) 2007 X. Andrade
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 submesh_oct_m
22 use accel_oct_m
23 use batch_oct_m
25 use box_oct_m
26 use comm_oct_m
27 use debug_oct_m
28 use global_oct_m
29 use index_oct_m
30 use, intrinsic :: iso_fortran_env
34 use sort_oct_m
35 use mesh_oct_m
37 use mpi_oct_m
40 use space_oct_m
41 use types_oct_m
42
43 implicit none
44 private
45
46 public :: &
47 submesh_t, &
60 dsm_nrm2, &
61 zsm_nrm2, &
79
82 type submesh_t
83 ! Components are public by default
84 real(real64), allocatable :: center(:)
85 real(real64) :: radius = m_zero
86 class(box_t), pointer :: box => null()
87 integer :: np = -1
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()
95 logical :: overlap
97 integer :: np_global = -1
98 real(real64), allocatable :: rel_x_global(:,:)
99 integer, allocatable :: part_v(:)
100 integer, allocatable :: global2local(:)
101
102 type(mesh_cube_map_t) :: cube_map
103 end type submesh_t
104
105 interface submesh_add_to_mesh
107 end interface submesh_add_to_mesh
108
109 interface submesh_to_mesh_dotp
111 end interface submesh_to_mesh_dotp
112
113contains
115 ! -------------------------------------------------------------
116 ! Multipliers for recursive formulation of n-ellipsoid volume
117 ! simplifying the Gamma function
118 ! f(n) = 2f(n-2)/n, f(0)=1, f(1)=2
119 recursive real(real64) function f_n(dims) result(fn)
120 integer :: dims
121
122 if (dims == 0) then
123 fn = m_one
124 else if (dims == 1) then
125 fn = m_two
126 else
127 fn = m_two * f_n(dims - 2) / dims
128 end if
129
130 end function f_n
131
133 subroutine submesh_init_box(this, space, mesh, box, center)
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(:)
139
140 integer :: original_mesh_index, submesh_index
141 logical, allocatable :: submesh_points(:)
142 real(real64), allocatable :: xtmp(:, :), rtmp(:)
143
144 push_sub(submesh_init_box)
145
146 assert(.not. space%is_periodic())
147
148 ! Assign the mesh and the box
149 this%mesh => mesh
150 this%box => box
151 this%center = center
152
153 ! Get the total number of points inside the submesh. First of all get the number of points inside the submesh box
154 safe_allocate(submesh_points(1:mesh%np))
155 submesh_points = box%contains_points(mesh%np, mesh%x)
156 this%np = count(submesh_points)
157 ! Verify that the box of the submesh is contained into the box of the mesh
158 assert(this%np <= mesh%np)
159
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))
163 xtmp = m_zero
164 rtmp = m_zero
165
166 ! Now we generate a map between the points of the submesh and the points of the original mesh
167 submesh_index = 0
168 do original_mesh_index = 1, mesh%np
169 ! check that the point is contained in the new box
170 if (submesh_points(original_mesh_index)) then
171 submesh_index = submesh_index + 1
172 this%map(submesh_index) = original_mesh_index
173 ! Assign the coordinates to this new point
174 xtmp(:, submesh_index) = mesh%x(original_mesh_index, :) - center(:)
175 rtmp(submesh_index) = norm2(xtmp(:, submesh_index))
176 end if
177 end do
179 call submesh_reorder_points(this, space, xtmp, rtmp)
181 safe_deallocate_a(xtmp)
182 safe_deallocate_a(rtmp)
183 safe_deallocate_a(submesh_points)
186 end subroutine submesh_init_box
188 ! -------------------------------------------------------------
189 subroutine submesh_init(this, space, mesh, latt, center, rc)
190 type(submesh_t), intent(inout) :: this
191 class(space_t), intent(in) :: space
192 class(mesh_t), target, intent(in) :: mesh
193 type(lattice_vectors_t), intent(in) :: latt
194 real(real64), intent(in) :: center(1:space%dim)
195 real(real64), intent(in) :: rc
196
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
205
206
207 push_sub(submesh_init)
208 call profiling_in("SUBMESH_INIT")
209
210 assert(space%dim <= 3)
211
212 this%mesh => mesh
213
214 safe_allocate(this%center(1:space%dim))
215 this%center(:) = center(:)
216
217 this%radius = rc
218 ! We add a small number of avoid instabilities due to rounding errors
219 rc2 = rc**2 + tol
220
221 ! The spheres are generated differently for periodic coordinates,
222 ! mainly for performance reasons.
223 if (.not. space%is_periodic()) then
224
225 call profiling_in("SUBMESH_INIT_MAP_INV")
226 safe_allocate(map_inv(0:this%mesh%np))
227 map_inv(0:this%mesh%np) = 0
228
229 nmin = 0
230 nmax = 0
231
232 ! get a cube of points that contains the sphere
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
235
236 ! make sure that the cube is inside the grid
237 ! parts of the cube which would fall outside the simulation box are chopped off.
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))
240
241 ! Get the total number of points inside the sphere
242 is = 0 ! this index counts inner points
243 do iz = nmin(3), nmax(3)
244 do iy = nmin(2), nmax(2)
245 do ix = nmin(1), nmax(1)
246 ip = mesh_local_index_from_coords(mesh, [ix, iy, iz])
247 if (ip == 0 .or. ip > mesh%np) cycle
248 r2 = sum((mesh%x(ip, :) - center)**2)
249 if (r2 <= rc2) then
250 is = is + 1
251 map_inv(ip) = is
252 end if
253 end do
254 end do
255 end do
256 this%np = is
257 call profiling_out("SUBMESH_INIT_MAP_INV")
258
259 call profiling_in("SUBMESH_INIT_RTMP")
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))
263
264 ! Generate the table and the positions
265 do iz = nmin(3), nmax(3)
266 do iy = nmin(2), nmax(2)
267 do ix = nmin(1), nmax(1)
268 ip = mesh_local_index_from_coords(mesh, [ix, iy, iz])
269 if (ip == 0 .or. ip > mesh%np) cycle
270 is = map_inv(ip)
271 if (is == 0) cycle
272 this%map(is) = ip
273 xtmp(:, is) = mesh%x(ip, :) - center
274 rtmp(is) = norm2(xtmp(:,is))
275 end do
276 end do
277 end do
278
279 safe_deallocate_a(map_inv)
280 call profiling_out("SUBMESH_INIT_RTMP")
281
282 ! This is the case for a periodic system
283 else
284
285 ! Get the total number of points inside the sphere considering
286 ! replicas along PBCs
287
288 ! this requires some optimization
289
290 latt_iter = lattice_iterator_t(latt, rc)
291
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)
295 end do
296
297 !Recursive formulation for the volume of n-ellipsoid
298 !Garry Tee, NZ J. Mathematics Vol. 34 (2005) p. 165 eqs. 53,55
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)
302
303 safe_allocate(x(1:space%dim,1:mesh%np))
304 do ip = 1, mesh%np
305 x(:,ip) = mesh%x(ip,:)
306 end do
307 call profiling_in("SUBMESH_INIT_PERIODIC_R")
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))
311
312 is = 0
313 do ip = 1, mesh%np
314 do icell = 1, latt_iter%n_cells
315 xx = x(:,ip) - center_copies(:, icell)
316 if(any(abs(xx)>rc+tol)) cycle
317 r2 = sum(xx**2)
318 if (r2 > rc2) cycle
319 is = is + 1
320 map_temp(is) = ip
321 rtmp(is) = sqrt(r2)
322 xtmp(:, is) = xx
323 ! Note that xx can be outside the unit cell
324 end do
325 end do
326 assert(is < huge(is))
327 this%np = is
328
329 safe_allocate(this%map(1:this%np))
330 this%map(1:this%np) = map_temp(1:this%np)
331 call profiling_out("SUBMESH_INIT_PERIODIC_R")
332
333 safe_deallocate_a(map_temp)
334 safe_deallocate_a(center_copies)
335 safe_deallocate_a(x)
336
337 end if
338
339 call submesh_reorder_points(this, space, xtmp, rtmp)
340 call profiling_out("SUBMESH_INIT")
341
342 safe_deallocate_a(xtmp)
343 safe_deallocate_a(rtmp)
344
345 pop_sub(submesh_init)
346 end subroutine submesh_init
347
348 subroutine submesh_reorder_points(this, space, xtmp, rtmp)
349 type(submesh_t), intent(inout) :: this
350 class(space_t), intent(in) :: space
351 real(real64), intent(in) :: xtmp(:, :), rtmp(:)
352
353 integer :: ip, i_region, offset
354 integer, allocatable :: order(:), order_new(:)
355 integer, allocatable :: map_new(:)
356 integer, allocatable :: np_region(:), tmp_array(:)
357
358 push_sub(submesh_reorder_points)
359
360 ! now order points for better locality
361
362 call profiling_in("SUBMESH_INIT_ORDER")
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))
366
367 do ip = 1, this%np
368 order(ip) = ip
369 end do
370
371 ! First we just reorder in order to determine overlap:
372 call sort(this%map, order)
373
374 !check whether points overlap (i.e. whether a submesh contains the same point more than once)
375 this%overlap = .false.
376 do ip = 1, this%np - 1
377 if (this%map(ip) == this%map(ip + 1)) then
378 ! this simplified test works, as the points are ordered.
379 this%overlap = .true.
380 exit
381 end if
382 end do
383
384 this%num_regions = 1
385 call profiling_out("SUBMESH_INIT_ORDER")
386
387 if(this%overlap) then
388 call profiling_in("SUBMESH_INIT_OVERLAP")
389 !disentangle the map into injective regions
390
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))
395
396 np_region(1) = 1
397 tmp_array(1) = 1
398 i_region = 1
399
400 do ip = 2, 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
406 end if
407 else
408 i_region = 1
409 end if
410 tmp_array(ip) = i_region ! which region does ip belong to
411 np_region(i_region) = np_region(i_region) + 1 ! increase number of points in i_region
412 end do
413
414 assert( .not. allocated(this%regions))
415
416 ! construct array of offsets
417 safe_allocate(this%regions(1:this%num_regions+1))
418
419 this%regions(1) = 1
420
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)
424 end do
425 else
426 this%regions(2) = this%np + 1
427 end if
428
429 np_region(1:this%np) = 0
430 order_new(1:this%np) = -1
431 map_new(1:this%np) = -1
432
433
434 !reassemble regions into global map array
435 do ip = 1, this%np
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)
441 end do
442
443 order(1:this%np) = order_new(1:this%np)
444 this%map(1:this%np) = map_new(1:this%np)
445
446 safe_deallocate_a(tmp_array)
447 safe_deallocate_a(order_new)
448 safe_deallocate_a(np_region)
449 safe_deallocate_a(map_new)
450 call profiling_out("SUBMESH_INIT_OVERLAP")
451
452 else
453 this%num_regions = 1
454 safe_allocate(this%regions(1:2))
455 this%regions(1) = 1
456 this%regions(2) = this%np + 1
457 end if
458
459 ! Lastly, reorder the points according to the new scheme
460 do ip = 1, this%np
461 this%rel_x(:, ip) = xtmp(:, order(ip))
462 this%r(ip) = rtmp(order(ip))
463 end do
464
465 safe_deallocate_a(order)
466
468 end subroutine submesh_reorder_points
469
470
471 ! --------------------------------------------------------------
472 !This routine takes two submeshes and merge them into a bigger submesh
473 !The grid is centered on the first center
474 subroutine submesh_merge(this, space, mesh, sm1, sm2, shift)
475 type(submesh_t), intent(inout) :: this
476 class(space_t), intent(in) :: space
477 class(mesh_t), target, intent(in) :: mesh
478 type(submesh_t), intent(in) :: sm1
479 type(submesh_t), intent(in) :: sm2
480 real(real64), optional, intent(in) :: shift(:)
481
482 real(real64) :: r2
483 integer :: ip, is
484 real(real64) :: xx(space%dim), diff_centers(space%dim)
485
486 push_sub(submesh_merge)
487 call profiling_in("SUBMESH_MERGE")
488
489 this%mesh => mesh
490
491 safe_allocate(this%center(1:space%dim))
492 this%center(:) = sm1%center(:)
493 this%radius = sm1%radius
494
495 ! This is a quick fix to prevent uninitialized variables. To properly check the self-overlap,
496 ! a similar approach as in submesh_init should be taken with respect to the merged map.
497 this%overlap = sm1%overlap .or. sm2%overlap
498
499 diff_centers = sm1%center - sm2%center
500 if (present(shift)) diff_centers = diff_centers - shift
501
502 !As we take the union of the two submeshes, we know that we have all the points from the first one included.
503 !The extra points from the second submesh are those which are not included in the first one
504 is = sm1%np
505 do ip = 1, sm2%np
506 !sm2%x contains points coordinates defined with respect to sm2%center
507 xx = sm2%rel_x(:, ip) - diff_centers
508 !If the point is not in sm1, we add it
509 if (sum(xx**2) > sm1%radius**2) is = is + 1
510 end do
511
512 this%np = is
513
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)
520
521 !iterate again to fill the tables
522 is = sm1%np
523 do ip = 1, sm2%np
524 xx = sm2%rel_x(:, ip) - diff_centers
525 r2 = sum(xx**2)
526 if (r2 > sm1%radius**2) then
527 is = is + 1
528 this%map(is) = sm2%map(ip)
529 this%r(is) = sqrt(r2)
530 this%rel_x(:, is) = xx
531 end if
532 end do
533
534 call profiling_out("SUBMESH_MERGE")
535 pop_sub(submesh_merge)
536 end subroutine submesh_merge
537
538 ! --------------------------------------------------------------
539 !This routine shifts the center of a submesh, without changing the grid points
540 subroutine submesh_shift_center(this, space, newcenter)
541 type(submesh_t), intent(inout) :: this
542 class(space_t), intent(in) :: space
543 real(real64), intent(in) :: newcenter(:)
544
545 real(real64) :: xx(space%dim), diff_centers(space%dim), oldcenter(space%dim)
546 integer :: ip
547
548 push_sub(submesh_shift_center)
549 call profiling_in("SUBMESH_SHIFT")
550
551 oldcenter = this%center
552 this%center(:) = newcenter(:)
553
554 diff_centers = newcenter - oldcenter
555
556 do ip = 1, this%np
557 xx = this%rel_x(:, ip) - diff_centers
558 this%r(ip) = norm2(xx)
559 this%rel_x(:, ip) = xx
560 end do
561
562 call profiling_out("SUBMESH_SHIFT")
563 pop_sub(submesh_shift_center)
564 end subroutine submesh_shift_center
565
566 ! --------------------------------------------------------------
567 subroutine submesh_broadcast(this, space, mesh, center, radius, root, mpi_grp)
568 type(submesh_t), intent(inout) :: this
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
574 type(mpi_grp_t), intent(in) :: mpi_grp
575
576 integer :: nparray(1:3)
577
578 push_sub(submesh_broadcast)
579 call profiling_in('SUBMESH_BCAST')
580
581 if (root /= mpi_grp%rank) then
582 this%mesh => mesh
583 safe_allocate(this%center(1:space%dim))
584 this%center(:) = center(:)
585 this%radius = radius
586 end if
587
588 if (mpi_grp%size > 1) then
589
590 if (root == mpi_grp%rank) then
591 nparray(1) = this%np
592 nparray(2) = this%num_regions
593 if (this%overlap) then
594 nparray(3) = 1
595 else
596 nparray(3) = 0
597 end if
598 end if
599
600 call mpi_grp%bcast(nparray, 3, mpi_integer, root)
601 this%np = nparray(1)
602 this%num_regions = nparray(2)
603 this%overlap = (nparray(3) == 1)
604
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))
610 end if
611
612 call mpi_grp%bcast(this%regions(1), this%num_regions+1, mpi_integer, root)
613
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)
618 end if
619
620 end if
621
622 call profiling_out('SUBMESH_BCAST')
623 pop_sub(submesh_broadcast)
624 end subroutine submesh_broadcast
625
626 ! --------------------------------------------------------------
627 logical function submesh_compatible(this, radius, center, dx) result(compatible)
628 type(submesh_t), intent(in) :: this
629 real(real64), intent(in) :: radius
630 real(real64), intent(in) :: center(:)
631 real(real64), intent(in) :: dx
632
633 compatible =.false.
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
638 compatible = .true.
639 end if
640 end if
641
642 end function submesh_compatible
643
644 ! --------------------------------------------------------------
645 subroutine submesh_end(this)
646 type(submesh_t), intent(inout) :: this
647
648 push_sub(submesh_end)
649
650 nullify(this%mesh)
651 nullify(this%box)
652 if (this%np /= -1) then
653 this%np = -1
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)
659 end if
661 if (accel_is_enabled()) then
662 call accel_release_buffer(this%buff_map)
663 end if
664
665 pop_sub(submesh_end)
666 end subroutine submesh_end
667
668 ! --------------------------------------------------------------
669
670 subroutine submesh_get_inv(this, map_inv)
671 type(submesh_t), intent(in) :: this
672 integer, intent(out) :: map_inv(:)
673
674 integer :: is
675
676 push_sub(submesh_get_inv)
677
678 map_inv(1:this%mesh%np_part) = 0
679 do is = 1, this%np
680 map_inv(this%map(is)) = is
681 end do
682
683 pop_sub(submesh_get_inv)
684 end subroutine submesh_get_inv
685
686 ! --------------------------------------------------------------
687 logical function submesh_overlap(sm1, sm2, space) result(overlap)
688 type(submesh_t), intent(in) :: sm1
689 type(submesh_t), intent(in) :: sm2
690 class(space_t), intent(in) :: space
691
692 integer :: ii, jj, dd
693 real(real64) :: distance
694
695 !no PUSH_SUB, called too often
696
697 if (.not. space%is_periodic()) then
698 !first check the distance
699 distance = sum((sm1%center - sm2%center)**2)
700 overlap = distance <= (1.5_real64*(sm1%radius + sm2%radius))**2
701
702 ! if they are very far, no need to check in detail
703 if (.not. overlap) return
704 end if
705
706 ! Otherwise check whether they have the some point in common. We
707 ! can make the comparison faster using that the arrays are sorted.
708 overlap = .false.
709 ii = 1
710 jj = 1
711 do while(ii <= sm1%np .and. jj <= sm2%np)
712 dd = sm1%map(ii) - sm2%map(jj)
713 if (dd < 0) then
714 ii = ii + 1
715 else if (dd > 0) then
716 jj = jj + 1
717 else
718 overlap = .true.
719 exit
720 end if
721 end do
722
723 if (sm1%mesh%parallel_in_domains) then
724 call sm1%mesh%mpi_grp%allreduce_inplace(overlap, 1, mpi_logical, mpi_lor)
725 end if
726
727 end function submesh_overlap
728
729 ! -------------------------------------------------------------
730 subroutine submesh_build_global(this, space)
731 type(submesh_t), intent(inout) :: this
732 class(space_t), intent(in) :: space
733
734 integer, allocatable :: part_np(:)
735 integer :: ipart, ind, ip
736
737 push_sub(submesh_build_global)
739 if (.not. this%mesh%parallel_in_domains) then
740 this%np_global = this%np
741 pop_sub(submesh_build_global)
742 return
743 end if
744
745 safe_allocate(part_np(this%mesh%pv%npart))
746 part_np = 0
747 part_np(this%mesh%pv%partno) = this%np
748
749 call this%mesh%allreduce(part_np)
750 this%np_global = sum(part_np)
751
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
758
759 ind = 0
760 do ipart = 1, this%mesh%pv%npart
761 if (ipart == this%mesh%pv%partno) then
762 do ip = 1, this%np
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
766 end do
767 end if
768 ind = ind + part_np(ipart)
769 end do
770
771 call this%mesh%allreduce(this%rel_x_global)
772 call this%mesh%allreduce(this%part_v)
773 call this%mesh%allreduce(this%global2local)
774
775 safe_deallocate_a(part_np)
776
777 pop_sub(submesh_build_global)
778 end subroutine submesh_build_global
779
780 ! -----------------------------------------------------------
781 subroutine submesh_end_global(this)
782 type(submesh_t), intent(inout) :: this
783
784 push_sub(submesh_end_global)
785
786 safe_deallocate_a(this%rel_x_global)
787 this%np_global = -1
788 safe_deallocate_a(this%part_v)
789 safe_deallocate_a(this%global2local)
790
791 pop_sub(submesh_end_global)
792 end subroutine submesh_end_global
793
794
795 ! -----------------------------------------------------------
796 subroutine zzsubmesh_add_to_mesh(this, sphi, phi, factor)
797 type(submesh_t), intent(in) :: this
798 complex(real64), intent(in) :: sphi(:)
799 complex(real64), intent(inout) :: phi(:)
800 complex(real64), optional, intent(in) :: factor
801
802 integer :: ip, m
803
804 push_sub(zzsubmesh_add_to_mesh)
805
806 if (present(factor)) then
807 !Loop unrolling inspired by BLAS axpy routine
808 m = mod(this%np,4)
809 do ip = 1, m
810 phi(this%map(ip)) = phi(this%map(ip)) + factor*sphi(ip)
811 end do
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)
818 end do
819 end if
820 else
821 m = mod(this%np,4)
822 do ip = 1, m
823 phi(this%map(ip)) = phi(this%map(ip)) + sphi(ip)
824 end do
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)
831 end do
832 end if
833 end if
834
835 pop_sub(zzsubmesh_add_to_mesh)
836 end subroutine zzsubmesh_add_to_mesh
837
838 !------------------------------------------------------------
839 complex(real64) function zzsubmesh_to_mesh_dotp(this, sphi, phi, reduce) result(dotp)
840 type(submesh_t), intent(in) :: this
841 complex(real64), intent(in) :: sphi(:)
842 complex(real64), intent(in) :: phi(:)
843 logical, optional, intent(in) :: reduce
844
845 integer :: is, m, ip
846
847 push_sub(zzsubmesh_to_mesh_dotp)
848
849 dotp = m_z0
850
851 if (this%mesh%use_curvilinear) then
852 do is = 1, this%np
853 dotp = dotp + this%mesh%vol_pp(this%map(is))*phi(this%map(is))*conjg(sphi(is))
854 end do
855 else
856 m = mod(this%np,4)
857 do ip = 1, m
858 dotp = dotp + phi(this%map(ip))*conjg(sphi(ip))
859 end do
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))
866 end do
867 end if
868 dotp = dotp*this%mesh%vol_pp(1)
869 end if
870
871 if (optional_default(reduce, .true.)) then
872 call profiling_in("SM_REDUCE_DOTP")
873 call this%mesh%allreduce(dotp)
874 call profiling_out("SM_REDUCE_DOTP")
875 end if
876
878 end function zzsubmesh_to_mesh_dotp
879
880 !------------------------------------------------------------
882 subroutine submesh_get_cube_dim(sm, space, db)
883 type(submesh_t), target, intent(in) :: sm
884 class(space_t), intent(in) :: space
885 integer, intent(out) :: db(1:space%dim)
886
887 integer :: ip, idir
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
891
892 push_sub(submesh_get_cube_dim)
893
894 max_chi = m_zero
895 do ip = 1, sm%np
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))
899 end do
900 end do
901
902 do idir = 1, space%dim
903 db(idir) = nint(max_chi(idir)-tol)
904 end do
905
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)
909 end if
910
911 db = 2 * db + 1
912
913 pop_sub(submesh_get_cube_dim)
914 end subroutine submesh_get_cube_dim
915
916 !------------------------------------------------------------
917 subroutine submesh_init_cube_map(sm, space)
918 type(submesh_t), target, intent(inout) :: sm
919 class(space_t), intent(in) :: space
920
921 integer(int64) :: ip
922 integer ::idir
923 real(real64) :: chi(space%dim), shift(space%dim)
924
925 push_sub(submesh_init_cube_map)
926
927 sm%cube_map%nmap = sm%np
928
929 safe_allocate(sm%cube_map%map(1:space%dim, 1:sm%cube_map%nmap))
930
931 !The center of the submesh does not belong to the mesh
932 !So we first need to find the closest grid point, and center the cube to it
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)
936 end do
937 shift = sm%mesh%coord_system%to_cartesian(shift)
938 shift = shift - sm%center
939
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))
944 end do
945 end do
946
947 if (accel_is_enabled()) then
948 call accel_create_buffer(sm%cube_map%map_buffer, accel_mem_read_only, type_integer, sm%cube_map%nmap*space%dim)
949 call accel_write_buffer(sm%cube_map%map_buffer, sm%cube_map%nmap*space%dim, sm%cube_map%map)
950 end if
951
952 pop_sub(submesh_init_cube_map)
953 end subroutine submesh_init_cube_map
954
955 !------------------------------------------------------------
956 subroutine submesh_end_cube_map(sm)
957 type(submesh_t), intent(inout) :: sm
958
959 push_sub(submesh_end_cube_map)
960
961 call mesh_cube_map_end(sm%cube_map)
962
963 pop_sub(submesh_end_cube_map)
964 end subroutine submesh_end_cube_map
965
972 subroutine dzsubmesh_batch_add(this, ss, mm)
973 type(submesh_t), intent(in) :: this
974 class(batch_t), intent(in) :: ss
975 class(batch_t), intent(inout) :: mm
976
977 integer :: ist, idim, jdim, is
978
979 push_sub(dzsubmesh_batch_add)
980
981 assert(.not. mm%is_packed())
982 assert(ss%type() == type_float)
983 assert(mm%type() == type_cmplx)
984 assert(ss%nst_linear == mm%nst_linear)
985 assert(ss%status() == mm%status())
986 assert(ss%dim == mm%dim)
987
988 assert(mm%nst == ss%nst)
989
990 !$omp parallel do private(ist, idim, jdim, is) if(.not. this%overlap)
991 do ist = 1, mm%nst
992 do idim = 1, mm%dim
993 jdim = min(idim, ss%dim)
994
995 do is = 1, this%np
996 mm%zff(this%map(is), idim, ist) = &
997 mm%zff(this%map(is), idim, ist) + ss%dff(is, jdim, ist)
998 end do
999
1000 end do
1001 end do
1002 !$omp end parallel do
1003
1004 pop_sub(dzsubmesh_batch_add)
1005 end subroutine dzsubmesh_batch_add
1006
1007
1008#include "undef.F90"
1009#include "real.F90"
1010#include "submesh_inc.F90"
1011
1012#include "undef.F90"
1013#include "complex.F90"
1014#include "submesh_inc.F90"
1015
1016end module submesh_oct_m
1017
1018!! Local Variables:
1019!! mode: f90
1020!! coding: utf-8
1021!! End:
This is the common interface to a sorting routine. It performs the shell algorithm,...
Definition: sort.F90:149
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)
Definition: accel.F90:1248
pure logical function, public accel_is_enabled()
Definition: accel.F90:401
integer, parameter, public accel_mem_read_only
Definition: accel.F90:183
This module implements batches of mesh functions.
Definition: batch.F90:133
Module implementing boundary conditions in Octopus.
Definition: boundaries.F90:122
real(real64), parameter, public m_two
Definition: global.F90:190
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:186
complex(real64), parameter, public m_z0
Definition: global.F90:198
real(real64), parameter, public m_one
Definition: global.F90:189
This module implements the index, used for the mesh points.
Definition: index.F90:122
subroutine, public mesh_cube_map_end(this)
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
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.
Definition: mesh.F90:939
Some general things and nomenclature:
Definition: par_vec.F90:171
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
This module is intended to contain "only mathematical" functions and procedures.
Definition: sort.F90:117
real(real64) function, public dsm_integrate_frommesh(mesh, sm, ff, reduce)
Definition: submesh.F90:1210
subroutine zdsubmesh_add_to_mesh(this, sphi, phi, factor)
Definition: submesh.F90:1780
subroutine, public submesh_init_cube_map(sm, space)
Definition: submesh.F90:1011
logical function, public submesh_overlap(sm1, sm2, space)
Definition: submesh.F90:781
subroutine, public zsubmesh_batch_dotp_matrix(this, mm, ss, dot, reduce)
Definition: submesh.F90:2085
real(real64) function, public dsm_nrm2(sm, ff, reduce)
this function returns the the norm of a vector
Definition: submesh.F90:1321
subroutine zzsubmesh_add_to_mesh(this, sphi, phi, factor)
Definition: submesh.F90:890
subroutine, public submesh_end_global(this)
Definition: submesh.F90:875
subroutine, public zsubmesh_copy_from_mesh(this, phi, sphi, conjugate)
Definition: submesh.F90:1828
subroutine submesh_reorder_points(this, space, xtmp, rtmp)
Definition: submesh.F90:442
complex(real64) function, public zsm_integrate(mesh, sm, ff, reduce)
Definition: submesh.F90:1707
subroutine, public submesh_shift_center(this, space, newcenter)
Definition: submesh.F90:634
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.
Definition: submesh.F90:227
complex(real64) function zzsubmesh_to_mesh_dotp(this, sphi, phi, reduce)
Definition: submesh.F90:933
real(real64) function, public dsm_integrate(mesh, sm, ff, reduce)
Definition: submesh.F90:1170
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...
Definition: submesh.F90:1403
recursive real(real64) function f_n(dims)
Definition: submesh.F90:213
subroutine, public dsubmesh_copy_from_mesh(this, phi, sphi, conjugate)
Definition: submesh.F90:1291
real(real64) function ddsubmesh_to_mesh_dotp(this, sphi, phi, reduce)
Definition: submesh.F90:1358
subroutine, public submesh_merge(this, space, mesh, sm1, sm2, shift)
Definition: submesh.F90:568
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...
Definition: submesh.F90:1066
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...
Definition: submesh.F90:2033
real(real64) function, public zsm_nrm2(sm, ff, reduce)
this function returns the the norm of a vector
Definition: submesh.F90:1858
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...
Definition: submesh.F90:1940
subroutine, public submesh_end_cube_map(sm)
Definition: submesh.F90:1050
subroutine, public submesh_end(this)
Definition: submesh.F90:739
subroutine ddsubmesh_add_to_mesh(this, sphi, phi, factor)
Definition: submesh.F90:1243
subroutine, public submesh_get_cube_dim(sm, space, db)
finds the dimension of a box containing the submesh
Definition: submesh.F90:976
subroutine, public dsubmesh_batch_dotp_matrix(this, mm, ss, dot, reduce)
Definition: submesh.F90:1548
complex(real64) function, public zsm_integrate_frommesh(mesh, sm, ff, reduce)
Definition: submesh.F90:1747
logical function, public submesh_compatible(this, radius, center, dx)
Definition: submesh.F90:721
subroutine, public submesh_build_global(this, space)
Definition: submesh.F90:824
subroutine, public submesh_broadcast(this, space, mesh, center, radius, root, mpi_grp)
Definition: submesh.F90:661
subroutine, public submesh_get_inv(this, map_inv)
Definition: submesh.F90:764
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...
Definition: submesh.F90:1496
subroutine, public submesh_init(this, space, mesh, latt, center, rc)
Definition: submesh.F90:283
complex(real64) function zdsubmesh_to_mesh_dotp(this, sphi, phi, reduce)
Definition: submesh.F90:1895
type(type_t), parameter, public type_cmplx
Definition: types.F90:134
type(type_t), parameter, public type_integer
Definition: types.F90:135
type(type_t), parameter, public type_float
Definition: types.F90:133
Class defining batches of mesh functions.
Definition: batch.F90:159
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
This is defined even when running serial.
Definition: mpi.F90:142
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)
void distance(const int iatom, const int jatom, const double coordinates[], double *rr, double *rr2, double *rr6, double *rr7)
Definition: vdw_ts_low.c:2059