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, &
58 dsm_nrm2, &
59 zsm_nrm2, &
77
80 type submesh_t
81 ! Components are public by default
82 real(real64), allocatable :: center(:)
83 real(real64) :: radius = m_zero
84 class(box_t), pointer :: box => null()
85 integer :: np = -1
86 integer, allocatable :: map(:)
87 integer :: num_regions
88 integer, allocatable :: regions(:)
89 type(accel_mem_t) :: buff_map
90 real(real64), allocatable :: rel_x(:,:)
91 real(real64), allocatable :: r(:)
92 type(mesh_t), pointer :: mesh => null()
93 logical :: overlap
95 integer :: np_global = -1
96 real(real64), allocatable :: rel_x_global(:,:)
97 integer, allocatable :: part_v(:)
98 integer, allocatable :: global2local(:)
99
100 type(mesh_cube_map_t) :: cube_map
101 end type submesh_t
102
103 interface submesh_add_to_mesh
105 end interface submesh_add_to_mesh
106
107 interface submesh_to_mesh_dotp
109 end interface submesh_to_mesh_dotp
110
111contains
112
113 ! -------------------------------------------------------------
114 ! Multipliers for recursive formulation of n-ellipsoid volume
115 ! simplifying the Gamma function
116 ! f(n) = 2f(n-2)/n, f(0)=1, f(1)=2
117 recursive real(real64) function f_n(dims) result(fn)
118 integer :: dims
119
120 if (dims == 0) then
121 fn = m_one
122 else if (dims == 1) then
123 fn = m_two
124 else
125 fn = m_two * f_n(dims - 2) / dims
126 end if
127
128 end function f_n
129
130 ! -------------------------------------------------------------
131 subroutine submesh_init(this, space, mesh, latt, center, rc)
132 type(submesh_t), intent(inout) :: this
133 class(space_t), intent(in) :: space
134 class(mesh_t), target, intent(in) :: mesh
135 type(lattice_vectors_t), intent(in) :: latt
136 real(real64), intent(in) :: center(1:space%dim)
137 real(real64), intent(in) :: rc
138
139 real(real64) :: r2, rc2, xx(space%dim), rc_norm_n
140 real(real64), allocatable :: center_copies(:,:), xtmp(:, :), rtmp(:)
141 integer :: icell, is, ip, ix, iy, iz
142 integer(int64) :: max_elements_count
143 type(lattice_iterator_t) :: latt_iter
144 integer, allocatable :: map_inv(:), map_temp(:)
145 integer :: nmax(3), nmin(3)
146 real(real64), parameter :: tol = 1e-13_real64
147
148
149 push_sub(submesh_init)
150 call profiling_in("SUBMESH_INIT")
151
152 assert(space%dim <= 3)
153
154 this%mesh => mesh
155
156 safe_allocate(this%center(1:space%dim))
157 this%center(:) = center(:)
158
159 this%radius = rc
160 ! We add a small number of avoid instabilities due to rounding errors
161 rc2 = rc**2 + tol
162
163 ! The spheres are generated differently for periodic coordinates,
164 ! mainly for performance reasons.
165 if (.not. space%is_periodic()) then
166
167 call profiling_in("SUBMESH_INIT_MAP_INV")
168 safe_allocate(map_inv(0:this%mesh%np))
169 map_inv(0:this%mesh%np) = 0
170
171 nmin = 0
172 nmax = 0
173
174 ! get a cube of points that contains the sphere
175 nmin(1:space%dim) = int((center(1:space%dim) - abs(rc))/mesh%spacing(1:space%dim)) - 1
176 nmax(1:space%dim) = int((center(1:space%dim) + abs(rc))/mesh%spacing(1:space%dim)) + 1
178 ! make sure that the cube is inside the grid
179 ! parts of the cube which would fall outside the simulation box are chopped off.
180 nmin(1:space%dim) = max(mesh%idx%nr(1, 1:space%dim), nmin(1:space%dim))
181 nmax(1:space%dim) = min(mesh%idx%nr(2, 1:space%dim), nmax(1:space%dim))
183 ! Get the total number of points inside the sphere
184 is = 0 ! this index counts inner points
185 do iz = nmin(3), nmax(3)
186 do iy = nmin(2), nmax(2)
187 do ix = nmin(1), nmax(1)
188 ip = mesh_local_index_from_coords(mesh, [ix, iy, iz])
189 if (ip == 0 .or. ip > mesh%np) cycle
190 r2 = sum((mesh%x(:, ip) - center)**2)
191 if (r2 <= rc2) then
192 is = is + 1
193 map_inv(ip) = is
194 end if
195 end do
196 end do
197 end do
198 this%np = is
199 call profiling_out("SUBMESH_INIT_MAP_INV")
200
201 call profiling_in("SUBMESH_INIT_RTMP")
202 safe_allocate(this%map(1:this%np))
203 safe_allocate(xtmp(1:space%dim, 1:this%np))
204 safe_allocate(rtmp(1:this%np))
205
206 ! Generate the table and the positions
207 do iz = nmin(3), nmax(3)
208 do iy = nmin(2), nmax(2)
209 do ix = nmin(1), nmax(1)
210 ip = mesh_local_index_from_coords(mesh, [ix, iy, iz])
211 if (ip == 0 .or. ip > mesh%np) cycle
212 is = map_inv(ip)
213 if (is == 0) cycle
214 this%map(is) = ip
215 xtmp(:, is) = mesh%x(:, ip) - center
216 rtmp(is) = norm2(xtmp(:,is))
217 end do
218 end do
219 end do
220
221 safe_deallocate_a(map_inv)
222 call profiling_out("SUBMESH_INIT_RTMP")
223
224 ! This is the case for a periodic system
225 else
227 ! Get the total number of points inside the sphere considering
228 ! replicas along PBCs
229
230 ! this requires some optimization
231
232 latt_iter = lattice_iterator_t(latt, rc)
233
234 safe_allocate(center_copies(1:space%dim, 1:latt_iter%n_cells))
235 do icell = 1, latt_iter%n_cells
236 center_copies(:, icell) = center + latt_iter%get(icell)
237 end do
238
239 !Recursive formulation for the volume of n-ellipsoid
240 !Garry Tee, NZ J. Mathematics Vol. 34 (2005) p. 165 eqs. 53,55
241 rc_norm_n = product(ceiling(rc / mesh%spacing(1:space%dim), int64) + m_one)
242 if (mesh%use_curvilinear) rc_norm_n = rc_norm_n / mesh%coord_system%min_mesh_scaling_product
243 max_elements_count = 3**space%dim * int(m_pi**floor(0.5 * space%dim) * rc_norm_n * f_n(space%dim), int64)
244
245 call profiling_in("SUBMESH_INIT_PERIODIC_R")
246 safe_allocate(map_temp(1:max_elements_count))
247 safe_allocate(xtmp(1:space%dim, 1:max_elements_count))
248 safe_allocate(rtmp(1:max_elements_count))
249
250 is = 0
251 do ip = 1, mesh%np
252 do icell = 1, latt_iter%n_cells
253 xx = mesh%x(:,ip) - center_copies(:, icell)
254 if(any(abs(xx)>rc+tol)) cycle
255 r2 = sum(xx**2)
256 if (r2 > rc2) cycle
257 is = is + 1
258 map_temp(is) = ip
259 rtmp(is) = sqrt(r2)
260 xtmp(:, is) = xx
261 ! Note that xx can be outside the unit cell
262 end do
263 end do
264 assert(is < huge(is))
265 this%np = is
266
267 safe_allocate(this%map(1:this%np))
268 this%map(1:this%np) = map_temp(1:this%np)
269 call profiling_out("SUBMESH_INIT_PERIODIC_R")
270
271 safe_deallocate_a(map_temp)
272 safe_deallocate_a(center_copies)
273
274 end if
275
276 call submesh_reorder_points(this, space, xtmp, rtmp)
277 call profiling_out("SUBMESH_INIT")
278
279 safe_deallocate_a(xtmp)
280 safe_deallocate_a(rtmp)
281
282 pop_sub(submesh_init)
283 end subroutine submesh_init
284
285 subroutine submesh_reorder_points(this, space, xtmp, rtmp)
286 type(submesh_t), intent(inout) :: this
287 class(space_t), intent(in) :: space
288 real(real64), intent(in) :: xtmp(:, :), rtmp(:)
289
290 integer :: ip, i_region, offset
291 integer, allocatable :: order(:), order_new(:)
292 integer, allocatable :: map_new(:)
293 integer, allocatable :: np_region(:), tmp_array(:)
294
295 push_sub(submesh_reorder_points)
296
297 ! now order points for better locality
298
299 call profiling_in("SUBMESH_INIT_ORDER")
300 safe_allocate(order(1:this%np))
301 safe_allocate(this%rel_x(1:space%dim, 1:this%np))
302 safe_allocate(this%r(1:this%np))
303
304 do ip = 1, this%np
305 order(ip) = ip
306 end do
307
308 ! First we just reorder in order to determine overlap:
309 call sort(this%map, order)
310
311 !check whether points overlap (i.e. whether a submesh contains the same point more than once)
312 this%overlap = .false.
313 do ip = 1, this%np - 1
314 if (this%map(ip) == this%map(ip + 1)) then
315 ! this simplified test works, as the points are ordered.
316 this%overlap = .true.
317 exit
318 end if
319 end do
320
321 this%num_regions = 1
322 call profiling_out("SUBMESH_INIT_ORDER")
323
324 if(this%overlap) then
325 call profiling_in("SUBMESH_INIT_OVERLAP")
326 !disentangle the map into injective regions
327
328 safe_allocate(tmp_array(1:this%np))
329 safe_allocate(order_new(1:this%np))
330 safe_allocate(np_region(1:this%np))
331 safe_allocate(map_new( 1:this%np))
332
333 np_region(1) = 1
334 tmp_array(1) = 1
335 i_region = 1
336
337 do ip = 2, this%np
338 if (this%map(ip) == this%map(ip - 1)) then
339 i_region = i_region + 1
340 if (i_region > this%num_regions) then
341 this%num_regions = i_region
342 np_region(i_region) = 0
343 end if
344 else
345 i_region = 1
346 end if
347 tmp_array(ip) = i_region ! which region does ip belong to
348 np_region(i_region) = np_region(i_region) + 1 ! increase number of points in i_region
349 end do
350
351 assert( .not. allocated(this%regions))
352
353 ! construct array of offsets
354 safe_allocate(this%regions(1:this%num_regions+1))
355
356 this%regions(1) = 1
357
358 if(this%num_regions > 1) then
359 do i_region = 1, this%num_regions
360 this%regions(i_region + 1) = this%regions(i_region) + np_region(i_region)
361 end do
362 else
363 this%regions(2) = this%np + 1
364 end if
365
366 np_region(1:this%np) = 0
367 order_new(1:this%np) = -1
368 map_new(1:this%np) = -1
369
370
371 !reassemble regions into global map array
372 do ip = 1, this%np
373 i_region = tmp_array(ip)
374 np_region(i_region) = np_region(i_region) + 1
375 offset = this%regions(i_region) - 1
376 map_new( offset + np_region(i_region) ) = this%map(ip)
377 order_new( offset + np_region(i_region) ) = order(ip)
378 end do
379
380 order(1:this%np) = order_new(1:this%np)
381 this%map(1:this%np) = map_new(1:this%np)
382
383 safe_deallocate_a(tmp_array)
384 safe_deallocate_a(order_new)
385 safe_deallocate_a(np_region)
386 safe_deallocate_a(map_new)
387 call profiling_out("SUBMESH_INIT_OVERLAP")
388
389 else
390 this%num_regions = 1
391 safe_allocate(this%regions(1:2))
392 this%regions(1) = 1
393 this%regions(2) = this%np + 1
394 end if
395
396 ! Lastly, reorder the points according to the new scheme
397 do ip = 1, this%np
398 this%rel_x(:, ip) = xtmp(:, order(ip))
399 this%r(ip) = rtmp(order(ip))
400 end do
401
402 safe_deallocate_a(order)
403
405 end subroutine submesh_reorder_points
406
407
408 ! --------------------------------------------------------------
409 !This routine takes two submeshes and merge them into a bigger submesh
410 !The grid is centered on the first center
411 subroutine submesh_merge(this, space, mesh, sm1, sm2, shift)
412 type(submesh_t), intent(inout) :: this
413 class(space_t), intent(in) :: space
414 class(mesh_t), target, intent(in) :: mesh
415 type(submesh_t), intent(in) :: sm1
416 type(submesh_t), intent(in) :: sm2
417 real(real64), optional, intent(in) :: shift(:)
418
419 real(real64) :: r2
420 integer :: ip, is
421 real(real64) :: xx(space%dim), diff_centers(space%dim)
422
423 push_sub(submesh_merge)
424 call profiling_in("SUBMESH_MERGE")
425
426 this%mesh => mesh
427
428 safe_allocate(this%center(1:space%dim))
429 this%center(:) = sm1%center(:)
430 this%radius = sm1%radius
431
432 ! This is a quick fix to prevent uninitialized variables. To properly check the self-overlap,
433 ! a similar approach as in submesh_init should be taken with respect to the merged map.
434 this%overlap = sm1%overlap .or. sm2%overlap
435
436 diff_centers = sm1%center - sm2%center
437 if (present(shift)) diff_centers = diff_centers - shift
438
439 !As we take the union of the two submeshes, we know that we have all the points from the first one included.
440 !The extra points from the second submesh are those which are not included in the first one
441 is = sm1%np
442 do ip = 1, sm2%np
443 !sm2%x contains points coordinates defined with respect to sm2%center
444 xx = sm2%rel_x(:, ip) - diff_centers
445 !If the point is not in sm1, we add it
446 if (sum(xx**2) > sm1%radius**2) is = is + 1
447 end do
448
449 this%np = is
450
451 safe_allocate(this%map(1:this%np))
452 safe_allocate(this%rel_x(1:space%dim, 1:this%np))
453 safe_allocate(this%r(1:this%np))
454 this%map(1:sm1%np) = sm1%map(1:sm1%np)
455 this%rel_x(:, 1:sm1%np) = sm1%rel_x(:, 1:sm1%np)
456 this%r(1:sm1%np) = sm1%r(1:sm1%np)
457
458 !iterate again to fill the tables
459 is = sm1%np
460 do ip = 1, sm2%np
461 xx = sm2%rel_x(:, ip) - diff_centers
462 r2 = sum(xx**2)
463 if (r2 > sm1%radius**2) then
464 is = is + 1
465 this%map(is) = sm2%map(ip)
466 this%r(is) = sqrt(r2)
467 this%rel_x(:, is) = xx
468 end if
469 end do
470
471 call profiling_out("SUBMESH_MERGE")
472 pop_sub(submesh_merge)
473 end subroutine submesh_merge
474
475 ! --------------------------------------------------------------
476 !This routine shifts the center of a submesh, without changing the grid points
477 subroutine submesh_shift_center(this, space, newcenter)
478 type(submesh_t), intent(inout) :: this
479 class(space_t), intent(in) :: space
480 real(real64), intent(in) :: newcenter(:)
481
482 real(real64) :: xx(space%dim), diff_centers(space%dim), oldcenter(space%dim)
483 integer :: ip
484
485 push_sub(submesh_shift_center)
486 call profiling_in("SUBMESH_SHIFT")
487
488 oldcenter = this%center
489 this%center(:) = newcenter(:)
490
491 diff_centers = newcenter - oldcenter
492
493 do ip = 1, this%np
494 xx = this%rel_x(:, ip) - diff_centers
495 this%r(ip) = norm2(xx)
496 this%rel_x(:, ip) = xx
497 end do
498
499 call profiling_out("SUBMESH_SHIFT")
500 pop_sub(submesh_shift_center)
501 end subroutine submesh_shift_center
502
503 ! --------------------------------------------------------------
504 subroutine submesh_broadcast(this, space, mesh, center, radius, root, mpi_grp)
505 type(submesh_t), intent(inout) :: this
506 class(space_t), intent(in) :: space
507 type(mesh_t), target, intent(in) :: mesh
508 real(real64), intent(in) :: center(1:space%dim)
509 real(real64), intent(in) :: radius
510 integer, intent(in) :: root
511 type(mpi_grp_t), intent(in) :: mpi_grp
512
513 integer :: nparray(1:3)
514
515 push_sub(submesh_broadcast)
516 call profiling_in('SUBMESH_BCAST')
517
518 if (root /= mpi_grp%rank) then
519 this%mesh => mesh
520 safe_allocate(this%center(1:space%dim))
521 this%center(:) = center(:)
522 this%radius = radius
523 end if
524
525 if (mpi_grp%size > 1) then
526
527 if (root == mpi_grp%rank) then
528 nparray(1) = this%np
529 nparray(2) = this%num_regions
530 if (this%overlap) then
531 nparray(3) = 1
532 else
533 nparray(3) = 0
534 end if
535 end if
536
537 call mpi_grp%bcast(nparray, 3, mpi_integer, root)
538 this%np = nparray(1)
539 this%num_regions = nparray(2)
540 this%overlap = (nparray(3) == 1)
541
542 if (root /= mpi_grp%rank) then
543 safe_allocate(this%map(1:this%np))
544 safe_allocate(this%rel_x(1:space%dim, 1:this%np))
545 safe_allocate(this%r(1:this%np))
546 safe_allocate(this%regions(1:this%num_regions+1))
547 end if
548
549 call mpi_grp%bcast(this%regions(1), this%num_regions+1, mpi_integer, root)
550
551 if (this%np > 0) then
552 call mpi_grp%bcast(this%map(1), this%np, mpi_integer, root)
553 call mpi_grp%bcast(this%rel_x(1, 1), this%np*space%dim, mpi_double_precision, root)
554 call mpi_grp%bcast(this%r(1), this%np, mpi_double_precision, root)
555 end if
556
557 end if
558
559 call profiling_out('SUBMESH_BCAST')
560 pop_sub(submesh_broadcast)
561 end subroutine submesh_broadcast
562
563 ! --------------------------------------------------------------
564 logical function submesh_compatible(this, radius, center, dx) result(compatible)
565 type(submesh_t), intent(in) :: this
566 real(real64), intent(in) :: radius
567 real(real64), intent(in) :: center(:)
568 real(real64), intent(in) :: dx
569
570 compatible =.false.
571 if (allocated(this%center)) then
574 if (radius <= this%radius+dx*1e-6_real64 .and. all(abs(this%center - center) < 0.25*dx)) then
575 compatible = .true.
576 end if
577 end if
578
579 end function submesh_compatible
580
581 ! --------------------------------------------------------------
582 subroutine submesh_end(this)
583 type(submesh_t), intent(inout) :: this
584
585 push_sub(submesh_end)
586
587 nullify(this%mesh)
588 nullify(this%box)
589 if (this%np /= -1) then
590 this%np = -1
591 safe_deallocate_a(this%center)
592 safe_deallocate_a(this%map)
593 safe_deallocate_a(this%rel_x)
594 safe_deallocate_a(this%r)
595 safe_deallocate_a(this%regions)
596 end if
597
598 if (accel_is_enabled()) then
599 call accel_free_buffer(this%buff_map)
600 end if
601
602 pop_sub(submesh_end)
603 end subroutine submesh_end
604
605
606 ! --------------------------------------------------------------
607 logical function submesh_overlap(sm1, sm2, space) result(overlap)
608 type(submesh_t), intent(in) :: sm1
609 type(submesh_t), intent(in) :: sm2
610 class(space_t), intent(in) :: space
611
612 integer :: ii, jj, dd
613 real(real64) :: distance
614
615 !no PUSH_SUB, called too often
616
617 if (.not. space%is_periodic()) then
618 !first check the distance
619 distance = sum((sm1%center - sm2%center)**2)
620 overlap = distance <= (1.5_real64*(sm1%radius + sm2%radius))**2
621
622 ! if they are very far, no need to check in detail
623 if (.not. overlap) return
624 end if
625
626 ! Otherwise check whether they have the some point in common. We
627 ! can make the comparison faster using that the arrays are sorted.
628 overlap = .false.
629 ii = 1
630 jj = 1
631 do while(ii <= sm1%np .and. jj <= sm2%np)
632 dd = sm1%map(ii) - sm2%map(jj)
633 if (dd < 0) then
634 ii = ii + 1
635 else if (dd > 0) then
636 jj = jj + 1
637 else
638 overlap = .true.
639 exit
640 end if
641 end do
642
643 if (sm1%mesh%parallel_in_domains) then
644 call sm1%mesh%mpi_grp%allreduce_inplace(overlap, 1, mpi_logical, mpi_lor)
645 end if
646
647 end function submesh_overlap
648
649 ! -------------------------------------------------------------
650 subroutine submesh_build_global(this, space)
651 type(submesh_t), intent(inout) :: this
652 class(space_t), intent(in) :: space
653
654 integer, allocatable :: part_np(:)
655 integer :: ipart, ind, ip
656
657 push_sub(submesh_build_global)
658
659 if (.not. this%mesh%parallel_in_domains) then
660 this%np_global = this%np
661 pop_sub(submesh_build_global)
662 return
663 end if
664
665 safe_allocate(part_np(this%mesh%pv%npart))
666 part_np = 0
667 part_np(this%mesh%pv%partno) = this%np
668
669 call this%mesh%allreduce(part_np)
670 this%np_global = sum(part_np)
671
672 safe_allocate(this%rel_x_global(1:space%dim, 1:this%np_global))
673 safe_allocate(this%part_v(1:this%np_global))
674 safe_allocate(this%global2local(1:this%np_global))
675 this%rel_x_global(1:space%dim, 1:this%np_global) = m_zero
676 this%part_v(1:this%np_global) = 0
677 this%global2local(1:this%np_global) = 0
678
679 ind = 0
680 do ipart = 1, this%mesh%pv%npart
681 if (ipart == this%mesh%pv%partno) then
682 do ip = 1, this%np
683 this%rel_x_global(:, ind + ip) = this%rel_x(:, ip)
684 this%part_v(ind + ip) = this%mesh%pv%partno
685 this%global2local(ind + ip) = ip
686 end do
687 end if
688 ind = ind + part_np(ipart)
689 end do
690
691 call this%mesh%allreduce(this%rel_x_global)
692 call this%mesh%allreduce(this%part_v)
693 call this%mesh%allreduce(this%global2local)
694
695 safe_deallocate_a(part_np)
696
697 pop_sub(submesh_build_global)
698 end subroutine submesh_build_global
699
700 ! -----------------------------------------------------------
701 subroutine submesh_end_global(this)
702 type(submesh_t), intent(inout) :: this
703
704 push_sub(submesh_end_global)
705
706 safe_deallocate_a(this%rel_x_global)
707 this%np_global = -1
708 safe_deallocate_a(this%part_v)
709 safe_deallocate_a(this%global2local)
710
711 pop_sub(submesh_end_global)
712 end subroutine submesh_end_global
713
714
715 ! -----------------------------------------------------------
716 subroutine zzsubmesh_add_to_mesh(this, sphi, phi, factor)
717 type(submesh_t), intent(in) :: this
718 complex(real64), intent(in) :: sphi(:)
719 complex(real64), intent(inout) :: phi(:)
720 complex(real64), optional, intent(in) :: factor
721
722 integer :: ip, m
723
724 push_sub(zzsubmesh_add_to_mesh)
725
726 if (present(factor)) then
727 !Loop unrolling inspired by BLAS axpy routine
728 m = mod(this%np,4)
729 do ip = 1, m
730 phi(this%map(ip)) = phi(this%map(ip)) + factor*sphi(ip)
731 end do
732 if (this%np.ge.4) then
733 do ip = m+1, this%np, 4
734 phi(this%map(ip)) = phi(this%map(ip)) + factor*sphi(ip)
735 phi(this%map(ip+1)) = phi(this%map(ip+1)) + factor*sphi(ip+1)
736 phi(this%map(ip+2)) = phi(this%map(ip+2)) + factor*sphi(ip+2)
737 phi(this%map(ip+3)) = phi(this%map(ip+3)) + factor*sphi(ip+3)
738 end do
739 end if
740 else
741 m = mod(this%np,4)
742 do ip = 1, m
743 phi(this%map(ip)) = phi(this%map(ip)) + sphi(ip)
744 end do
745 if (this%np.ge.4) then
746 do ip = m+1, this%np, 4
747 phi(this%map(ip)) = phi(this%map(ip)) + sphi(ip)
748 phi(this%map(ip+1)) = phi(this%map(ip+1)) + sphi(ip+1)
749 phi(this%map(ip+2)) = phi(this%map(ip+2)) + sphi(ip+2)
750 phi(this%map(ip+3)) = phi(this%map(ip+3)) + sphi(ip+3)
751 end do
752 end if
753 end if
754
755 pop_sub(zzsubmesh_add_to_mesh)
756 end subroutine zzsubmesh_add_to_mesh
757
758 !------------------------------------------------------------
759 complex(real64) function zzsubmesh_to_mesh_dotp(this, sphi, phi, reduce) result(dotp)
760 type(submesh_t), intent(in) :: this
761 complex(real64), intent(in) :: sphi(:)
762 complex(real64), intent(in) :: phi(:)
763 logical, optional, intent(in) :: reduce
764
765 integer :: is, m, ip
766
767 push_sub(zzsubmesh_to_mesh_dotp)
768
769 dotp = m_z0
770
771 if (this%mesh%use_curvilinear) then
772 do is = 1, this%np
773 dotp = dotp + this%mesh%vol_pp(this%map(is))*phi(this%map(is))*conjg(sphi(is))
774 end do
775 else
776 m = mod(this%np,4)
777 do ip = 1, m
778 dotp = dotp + phi(this%map(ip))*conjg(sphi(ip))
779 end do
780 if (this%np.ge.4) then
781 do ip = m+1, this%np, 4
782 dotp = dotp + phi(this%map(ip))*conjg(sphi(ip)) &
783 + phi(this%map(ip+1))*conjg(sphi(ip+1)) &
784 + phi(this%map(ip+2))*conjg(sphi(ip+2)) &
785 + phi(this%map(ip+3))*conjg(sphi(ip+3))
786 end do
787 end if
788 dotp = dotp*this%mesh%vol_pp(1)
789 end if
790
791 if (optional_default(reduce, .true.)) then
792 call profiling_in("SM_REDUCE_DOTP")
793 call this%mesh%allreduce(dotp)
794 call profiling_out("SM_REDUCE_DOTP")
795 end if
798 end function zzsubmesh_to_mesh_dotp
799
800 !------------------------------------------------------------
802 subroutine submesh_get_cube_dim(sm, space, db)
803 type(submesh_t), target, intent(in) :: sm
804 class(space_t), intent(in) :: space
805 integer, intent(out) :: db(1:space%dim)
806
807 integer :: ip, idir
808 real(real64) :: chi(space%dim), max_chi(space%dim)
809 integer :: db_red(1:space%dim)
810 real(real64), parameter :: tol=1.0e-10_real64
812 push_sub(submesh_get_cube_dim)
813
814 max_chi = m_zero
815 do ip = 1, sm%np
816 chi = sm%mesh%coord_system%from_cartesian(sm%rel_x(:, ip))
817 do idir = 1, space%dim
818 max_chi(idir) = max(max_chi(idir), abs(chi(idir))/sm%mesh%spacing(idir))
819 end do
820 end do
821
822 do idir = 1, space%dim
823 db(idir) = nint(max_chi(idir)-tol)
824 end do
825
826 if(sm%mesh%parallel_in_domains) then
827 call sm%mesh%mpi_grp%allreduce(db(1), db_red(1), space%dim, mpi_integer, mpi_max)
828 db(1:space%dim) = db_red(1:space%dim)
829 end if
830
831 db = 2 * db + 1
832
833 pop_sub(submesh_get_cube_dim)
834 end subroutine submesh_get_cube_dim
835
836 !------------------------------------------------------------
837 subroutine submesh_init_cube_map(sm, space)
838 type(submesh_t), target, intent(inout) :: sm
839 class(space_t), intent(in) :: space
840
841 integer(int64) :: ip
842 integer ::idir
843 real(real64) :: chi(space%dim), shift(space%dim)
844
845 push_sub(submesh_init_cube_map)
846
847 sm%cube_map%nmap = sm%np
848
849 safe_allocate(sm%cube_map%map(1:space%dim, 1:sm%cube_map%nmap))
850
851 !The center of the submesh does not belong to the mesh
852 !So we first need to find the closest grid point, and center the cube to it
853 chi = sm%mesh%coord_system%from_cartesian(sm%center)
854 do idir = 1, space%dim
855 shift(idir) = nint(chi(idir)/sm%mesh%spacing(idir))*sm%mesh%spacing(idir)
856 end do
857 shift = sm%mesh%coord_system%to_cartesian(shift)
858 shift = shift - sm%center
859
860 do ip = 1, sm%cube_map%nmap
861 chi = sm%mesh%coord_system%from_cartesian(sm%rel_x(:,ip) - shift)
862 do idir = 1, space%dim
863 sm%cube_map%map(idir, ip) = nint(chi(idir)/sm%mesh%spacing(idir))
864 end do
865 end do
866
867 if (accel_is_enabled()) then
868 call accel_create_buffer(sm%cube_map%map_buffer, accel_mem_read_only, type_integer, sm%cube_map%nmap*space%dim)
869 call accel_write_buffer(sm%cube_map%map_buffer, space%dim, sm%cube_map%nmap, sm%cube_map%map)
870 end if
871
872 pop_sub(submesh_init_cube_map)
873 end subroutine submesh_init_cube_map
874
875 !------------------------------------------------------------
876 subroutine submesh_end_cube_map(sm)
877 type(submesh_t), intent(inout) :: sm
878
879 push_sub(submesh_end_cube_map)
880
881 call mesh_cube_map_end(sm%cube_map)
882
883 pop_sub(submesh_end_cube_map)
884 end subroutine submesh_end_cube_map
885
892 subroutine dzsubmesh_batch_add(this, ss, mm)
893 type(submesh_t), intent(in) :: this
894 class(batch_t), intent(in) :: ss
895 class(batch_t), intent(inout) :: mm
896
897 integer :: ist, idim, jdim, is
898
899 push_sub(dzsubmesh_batch_add)
900
901 assert(.not. mm%is_packed())
902 assert(ss%type() == type_float)
903 assert(mm%type() == type_cmplx)
904 assert(ss%nst_linear == mm%nst_linear)
905 assert(ss%status() == mm%status())
906 assert(ss%dim == mm%dim)
907
908 assert(mm%nst == ss%nst)
909
910 !$omp parallel do private(ist, idim, jdim, is) if(.not. this%overlap)
911 do ist = 1, mm%nst
912 do idim = 1, mm%dim
913 jdim = min(idim, ss%dim)
914
915 do is = 1, this%np
916 mm%zff(this%map(is), idim, ist) = &
917 mm%zff(this%map(is), idim, ist) + ss%dff(is, jdim, ist)
918 end do
919
920 end do
921 end do
922 !$omp end parallel do
923
924 pop_sub(dzsubmesh_batch_add)
925 end subroutine dzsubmesh_batch_add
926
927
928#include "undef.F90"
929#include "real.F90"
930#include "submesh_inc.F90"
931
932#include "undef.F90"
933#include "complex.F90"
934#include "submesh_inc.F90"
935
936end module submesh_oct_m
937
938!! Local Variables:
939!! mode: f90
940!! coding: utf-8
941!! 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__
double floor(double __x) __attribute__((__nothrow__
double fn(const gsl_vector *v, void *params)
subroutine, public accel_free_buffer(this, async)
Definition: accel.F90:1005
pure logical function, public accel_is_enabled()
Definition: accel.F90:402
integer, parameter, public accel_mem_read_only
Definition: accel.F90:185
This module implements batches of mesh functions.
Definition: batch.F90:135
Module implementing boundary conditions in Octopus.
Definition: boundaries.F90:124
real(real64), parameter, public m_two
Definition: global.F90:202
real(real64), parameter, public m_zero
Definition: global.F90:200
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:198
complex(real64), parameter, public m_z0
Definition: global.F90:210
real(real64), parameter, public m_one
Definition: global.F90:201
This module implements the index, used for the mesh points.
Definition: index.F90:124
subroutine, public mesh_cube_map_end(this)
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
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:938
Some general things and nomenclature:
Definition: par_vec.F90:173
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:631
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:554
This module is intended to contain "only mathematical" functions and procedures.
Definition: sort.F90:119
real(real64) function, public dsm_integrate_frommesh(mesh, sm, ff, reduce)
Definition: submesh.F90:1132
subroutine zdsubmesh_add_to_mesh(this, sphi, phi, factor)
Definition: submesh.F90:1702
subroutine, public submesh_init_cube_map(sm, space)
Definition: submesh.F90:933
logical function, public submesh_overlap(sm1, sm2, space)
Definition: submesh.F90:703
subroutine, public zsubmesh_batch_dotp_matrix(this, mm, ss, dot, reduce)
Definition: submesh.F90:2007
real(real64) function, public dsm_nrm2(sm, ff, reduce)
this function returns the the norm of a vector
Definition: submesh.F90:1243
subroutine zzsubmesh_add_to_mesh(this, sphi, phi, factor)
Definition: submesh.F90:812
subroutine, public submesh_end_global(this)
Definition: submesh.F90:797
subroutine, public zsubmesh_copy_from_mesh(this, phi, sphi, conjugate)
Definition: submesh.F90:1750
subroutine submesh_reorder_points(this, space, xtmp, rtmp)
Definition: submesh.F90:381
complex(real64) function, public zsm_integrate(mesh, sm, ff, reduce)
Definition: submesh.F90:1629
subroutine, public submesh_shift_center(this, space, newcenter)
Definition: submesh.F90:573
complex(real64) function zzsubmesh_to_mesh_dotp(this, sphi, phi, reduce)
Definition: submesh.F90:855
real(real64) function, public dsm_integrate(mesh, sm, ff, reduce)
Definition: submesh.F90:1092
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:1325
recursive real(real64) function f_n(dims)
Definition: submesh.F90:213
subroutine, public dsubmesh_copy_from_mesh(this, phi, sphi, conjugate)
Definition: submesh.F90:1213
real(real64) function ddsubmesh_to_mesh_dotp(this, sphi, phi, reduce)
Definition: submesh.F90:1280
subroutine, public submesh_merge(this, space, mesh, sm1, sm2, shift)
Definition: submesh.F90:507
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:988
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:1955
real(real64) function, public zsm_nrm2(sm, ff, reduce)
this function returns the the norm of a vector
Definition: submesh.F90:1780
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:1862
subroutine, public submesh_end_cube_map(sm)
Definition: submesh.F90:972
subroutine, public submesh_end(this)
Definition: submesh.F90:678
subroutine ddsubmesh_add_to_mesh(this, sphi, phi, factor)
Definition: submesh.F90:1165
subroutine, public submesh_get_cube_dim(sm, space, db)
finds the dimension of a box containing the submesh
Definition: submesh.F90:898
subroutine, public dsubmesh_batch_dotp_matrix(this, mm, ss, dot, reduce)
Definition: submesh.F90:1470
complex(real64) function, public zsm_integrate_frommesh(mesh, sm, ff, reduce)
Definition: submesh.F90:1669
logical function, public submesh_compatible(this, radius, center, dx)
Definition: submesh.F90:660
subroutine, public submesh_build_global(this, space)
Definition: submesh.F90:746
subroutine, public submesh_broadcast(this, space, mesh, center, radius, root, mpi_grp)
Definition: submesh.F90:600
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:1418
subroutine, public submesh_init(this, space, mesh, latt, center, rc)
Definition: submesh.F90:227
complex(real64) function zdsubmesh_to_mesh_dotp(this, sphi, phi, reduce)
Definition: submesh.F90:1817
type(type_t), public type_float
Definition: types.F90:135
type(type_t), public type_cmplx
Definition: types.F90:136
type(type_t), public type_integer
Definition: types.F90:137
Class defining batches of mesh functions.
Definition: batch.F90:161
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:187
This is defined even when running serial.
Definition: mpi.F90:144
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