Octopus
cube.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2011 M. Marques, A. Castro, A. Rubio,
2!! G. Bertsch, M. Oliveira, J. Alberdi-Rodriguez
3!!
4!! This program is free software; you can redistribute it and/or modify
5!! it under the terms of the GNU General Public License as published by
6!! the Free Software Foundation; either version 2, or (at your option)
7!! any later version.
8!!
9!! This program is distributed in the hope that it will be useful,
10!! but WITHOUT ANY WARRANTY; without even the implied warranty of
11!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12!! GNU General Public License for more details.
13!!
14!! You should have received a copy of the GNU General Public License
15!! along with this program; if not, write to the Free Software
16!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
17!! 02110-1301, USA.
18!!
19
20#include "global.h"
21
22module cube_oct_m
23 use accel_oct_m
26 use debug_oct_m
27 use fft_oct_m
28 use global_oct_m
29 use io_oct_m
30 use, intrinsic :: iso_fortran_env
33 use mesh_oct_m
35 use mpi_oct_m
37#ifdef HAVE_NFFT
38 use nfft_oct_m
39#endif
40 use parser_oct_m
41 use pfft_oct_m
43 use space_oct_m
44
45 implicit none
46 private
47 public :: &
48 cube_t, &
50 cube_init, &
57
58 type cube_t
59 ! Components are public by default
60 logical :: parallel_in_domains
61 type(mpi_grp_t) :: mpi_grp
62
63 integer :: rs_n_global(1:3)
64 integer :: fs_n_global(1:3)
65 integer :: rs_n(1:3)
66 integer :: fs_n(1:3)
67 integer :: rs_istart(1:3)
68 integer :: fs_istart(1:3)
69 integer :: center(1:3)
70
71 real(real64), allocatable :: Lrs(:,:)
72 real(real64), allocatable :: Lfs(:,:)
73
74 integer, allocatable :: np_local(:)
75 integer, allocatable :: xlocal(:)
76 integer, allocatable :: local(:,:)
77 integer, allocatable :: np_local_fs(:)
78 integer, allocatable :: xlocal_fs(:)
79 integer, allocatable :: local_fs(:,:)
80
81
82 type(fft_t), allocatable :: fft
83 logical, private :: has_cube_mapping = .false.
85
86 real(real64) :: spacing(3)
87 ! latt is declared as allocatable as a work-around for a bug in gfortran when invoking the finalizer of latt.
88 type(lattice_vectors_t), allocatable :: latt
89
90 type(mesh_cube_map_t) :: cube_map
91 logical :: cube_map_present = .false.
92 end type cube_t
93
99 type dimensions_t
100 integer :: start_xyz(1:3)
101 integer :: end_xyz(1:3)
102 end type dimensions_t
103
104contains
105
106 ! ---------------------------------------------------------
107 subroutine cube_init(cube, nn, namespace, space, spacing, coord_system, fft_type, fft_library, dont_optimize, nn_out, &
108 mpi_grp, need_partition, tp_enlarge, blocksize)
109 type(cube_t), intent(out) :: cube
110 integer, intent(in) :: nn(:)
111 type(namespace_t), intent(in) :: namespace
112 class(space_t), intent(in) :: space
113 real(real64), intent(in) :: spacing(:)
114 class(coordinate_system_t), intent(in) :: coord_system
115 integer, optional, intent(in) :: fft_type
116 integer, optional, intent(in) :: fft_library
117 logical, optional, intent(in) :: dont_optimize
118 integer, optional, intent(out) :: nn_out(3)
120 type(mpi_grp_t), optional, intent(in) :: mpi_grp
121 logical, optional, intent(in) :: need_partition
122 real(real64), optional, intent(in) :: tp_enlarge(3)
125 integer, optional, intent(in) :: blocksize
127
128 type(MPI_Comm) :: comm
129 integer :: tmp_n(3), fft_type_, optimize_parity(3), fft_library_, nn3d(3)
130 integer :: effdim_fft, my_n(3), idir, idir2
131 logical :: optimize(3)
132 type(mpi_grp_t) :: mpi_grp_
133 real(real64) :: tp_enlarge_(3), lattice_vectors(3, 3)
134 type(space_t) :: cube_space
135
136 push_sub(cube_init)
137
138 assert(all(nn > 0))
139 assert(space%dim <= 3)
140
141 nn3d(1:space%dim) = nn(1:space%dim)
142 nn3d(space%dim+1:3) = 1
143
144 cube%spacing(1:space%dim) = spacing(1:space%dim)
145 cube%spacing(space%dim+1:3) = -m_one
146
147 fft_type_ = optional_default(fft_type, fft_none)
148 tp_enlarge_(:) = (/m_one, m_one, m_one/)
149 if (present(tp_enlarge)) tp_enlarge_(:)=tp_enlarge(:)
150
151 effdim_fft = min(3, space%dim)
152
153 mpi_grp_ = mpi_world
154 if (present(mpi_grp)) mpi_grp_ = mpi_grp
155
156 if (fft_type_ /= fft_none) then
158 if (present(fft_library)) then
159 fft_library_ = fft_library
160 else
161 fft_library_ = fft_default_lib
162 end if
163
164#ifndef HAVE_PFFT
165 if (fft_library_ == fftlib_pfft) then
166 write(message(1),'(a)')'You have selected the PFFT for FFT, but it was not linked.'
167 call messages_fatal(1, namespace=namespace)
168 end if
169#endif
171 else
172 fft_library_ = fftlib_none
173 end if
174
175 ! Note: later we set parallel_in_domains if blocksize is given, too
176 cube%parallel_in_domains = (fft_library_ == fftlib_pfft .or. fft_library_ == fftlib_pnfft)
177 if (present(blocksize)) then
178 assert(present(need_partition).and.need_partition)
179 assert(fft_library_ == fftlib_none)
180 ! For all the different FFT libraries there are strange (?)
181 ! rules about how the decomposition is chosen. What we want
182 ! (for libvdwxc) is a cube parallelized according to the simple
183 ! but contrary rule "just do what I say". Hence the blocksize
184 ! parameter. (Later to be expanded to allow 2D distributions.)
185 cube%rs_n_global = nn3d
186 cube%fs_n_global = nn3d ! not to be used
187 cube%fs_n = cube%fs_n_global ! not to be used
188 cube%fs_istart = 1 ! not to be used
189
190 comm = mpi_grp_%comm
191 cube%parallel_in_domains = (mpi_grp_%size > 1) ! XXX whether comm size > 1
192 call cube_set_blocksize(cube%rs_n_global, blocksize, mpi_grp_%rank, cube%rs_n, cube%rs_istart)
193 else if (fft_library_ == fftlib_none) then
194 cube%rs_n_global = nn3d
195 cube%fs_n_global = nn3d
196 cube%rs_n = cube%rs_n_global
197 cube%fs_n = cube%fs_n_global
198 cube%rs_istart = 1
199 cube%fs_istart = 1
201 if (present(nn_out)) nn_out(1:3) = nn(1:3)
202 else
203 safe_allocate(cube%fft)
204 tmp_n = nn3d
205
206 optimize(1:3) = .false.
207 optimize_parity(1:3) = 0
208 optimize(space%periodic_dim + 1:effdim_fft) = .true.
209 optimize_parity(space%periodic_dim + 1:effdim_fft) = 1
210
211 if (present(dont_optimize)) then
212 if (dont_optimize) optimize = .false.
213 end if
214
215 if (present(tp_enlarge)) call cube_tp_fft_defaults(cube, fft_library_)
216
217 call fft_init(cube%fft, tmp_n, space%dim, fft_type_, fft_library_, optimize, optimize_parity, &
218 comm=comm, mpi_grp = mpi_grp_, use_aligned=.true.)
219 if (present(nn_out)) nn_out(1:3) = tmp_n(1:3)
220
221 call fft_get_dims(cube%fft, cube%rs_n_global, cube%fs_n_global, cube%rs_n, cube%fs_n, &
222 cube%rs_istart, cube%fs_istart)
223
224 if (present(tp_enlarge)) then
225 call cube_init_coords(cube, tp_enlarge_, cube%spacing, fft_library_)
226 end if
227
228 if (fft_library_ == fftlib_nfft .or. fft_library_ == fftlib_pnfft) then
229 call fft_init_stage1(cube%fft, namespace, cube%Lrs, cube%rs_n_global)
230 !set local dimensions after stage1 - needed for PNFFT
231 call fft_get_dims(cube%fft, cube%rs_n_global, cube%fs_n_global, cube%rs_n, cube%fs_n, &
232 cube%rs_istart, cube%fs_istart)
233 end if
234
235 end if
236
237 if (.not. allocated(cube%Lrs)) then
238 call cube_init_coords(cube, tp_enlarge_, cube%spacing, fft_library_)
239 end if
240
241
242 cube%center(1:3) = cube%rs_n_global(1:3)/2 + 1
243
244
245 call mpi_grp_init(cube%mpi_grp, comm)
246
247 ! Initialize mapping only if needed
248 if (present(need_partition) .and. cube%parallel_in_domains) then
249 cube%has_cube_mapping = need_partition
250 else
251 cube%has_cube_mapping = .false.
252 end if
253 if (cube%has_cube_mapping) then
254 call cube_do_mapping(cube, fs = fft_library_ == fftlib_pnfft)
255 end if
256
257 if (cube%parallel_in_domains) call cube_partition_messages_debug(cube, namespace)
258
259 select type (coord_system)
260 class is (affine_coordinates_t)
261 ! We are constructing a lattice vector for the cube
262 ! This differs from the actual lattice vectors if it is not an integer multiple of the spacing
263 ! mesh%idx%ll is "general" in aperiodic directions,
264 ! but "periodic" in periodic directions.
265 my_n(1:space%periodic_dim) = cube%rs_n_global(1:space%periodic_dim) + 1
266 my_n(space%periodic_dim + 1:space%dim) = cube%rs_n_global(space%periodic_dim + 1:space%dim)
267
268 lattice_vectors = m_zero
269 do idir = 1, space%dim
270 do idir2 = 1, space%dim
271 lattice_vectors(idir2, idir) = cube%spacing(idir) * (my_n(idir) - 1) * coord_system%basis%vectors(idir2, idir)
272 end do
273 end do
274 do idir = space%dim + 1, 3
275 lattice_vectors(idir, idir) = m_one
276 end do
277
278 cube_space%dim = 3
279 cube_space%periodic_dim = space%periodic_dim
280 safe_allocate(cube%latt)
281 cube%latt = lattice_vectors_t(namespace, cube_space, lattice_vectors)
282 class default
283 message(1) = "The cube only support affine coordinate systems."
284 call messages_fatal(1, namespace=namespace)
285 end select
286
287 pop_sub(cube_init)
288 end subroutine cube_init
289
290 ! ---------------------------------------------------------
291 subroutine cube_end(cube)
292 type(cube_t), intent(inout) :: cube
293
294 push_sub(cube_end)
295
296 if (allocated(cube%fft)) then
297 call fft_end(cube%fft)
298 safe_deallocate_a(cube%fft)
299 end if
300
301 if (cube%has_cube_mapping) then
302 safe_deallocate_a(cube%np_local)
303 safe_deallocate_a(cube%xlocal)
304 safe_deallocate_a(cube%local)
305
306 safe_deallocate_a(cube%np_local_fs)
307 safe_deallocate_a(cube%xlocal_fs)
308 safe_deallocate_a(cube%local_fs)
309 end if
310
311 if (cube%cube_map_present) then
312 call mesh_cube_map_end(cube%cube_map)
313 end if
314
315 safe_deallocate_a(cube%Lrs)
316 safe_deallocate_a(cube%Lfs)
317
318 safe_deallocate_a(cube%latt)
319
320 pop_sub(cube_end)
321 end subroutine cube_end
322
323
324 ! ---------------------------------------------------------
325 subroutine cube_tp_fft_defaults(cube, fft_library)
326 type(cube_t), intent(inout) :: cube
327 integer, intent(in) :: fft_library
328
329 push_sub(cube_tp_fft_defaults)
330 select case (fft_library)
331 case (fftlib_nfft)
332#ifdef HAVE_NFFT
333 !Set NFFT defaults to values that gives good performance for two-point enlargement
334 !These values are overridden by the NFFT options in the input file
335 cube%fft%nfft%set_defaults = .true.
336 cube%fft%nfft%guru = .true.
337 cube%fft%nfft%mm = 2
338 cube%fft%nfft%sigma = 1.1_real64
339 cube%fft%nfft%precompute = nfft_pre_psi
340#endif
341 case (fftlib_pnfft)
342 cube%fft%pnfft%set_defaults = .true.
343 cube%fft%pnfft%m = 2
344 cube%fft%pnfft%sigma = 1.1_real64
345
346 case default
347 !do nothing
348 end select
349
350 pop_sub(cube_tp_fft_defaults)
351 end subroutine cube_tp_fft_defaults
352
353
354 ! ---------------------------------------------------------
355 subroutine cube_init_coords(cube, tp_enlarge, spacing, fft_library)
356 type(cube_t), intent(inout) :: cube
357 real(real64), intent(in) :: tp_enlarge(3)
358 real(real64), intent(in) :: spacing(3)
359 integer, intent(in) :: fft_library
360
361 real(real64) :: temp
362 integer :: ii, nn(3), maxn, idim
363
364 push_sub(cube_init_coords)
365
366
367 nn(1:3) = cube%fs_n_global(1:3)
368
369 maxn = maxval(nn)
370 safe_allocate(cube%Lrs(1:maxn, 1:3))
371 cube%Lrs(:,:) = m_zero
372
373 !! Real space coordinates
374 do idim = 1,3
375 if (tp_enlarge(idim) > m_one) then
376 do ii = 2, nn(idim) - 1
377 cube%Lrs(ii, idim) = (ii - int(nn(idim)/2) -1) * spacing(idim)
378 end do
379 cube%Lrs(1, idim) = (-int(nn(idim)/2)) * spacing(idim) * tp_enlarge(idim)
380 cube%Lrs(nn(idim), idim) = (int(nn(idim)/2)) * spacing(idim) * tp_enlarge(idim)
381 else
382 do ii = 1, nn(idim)
383 cube%Lrs(ii, idim) = (ii - int(nn(idim)/2) -1) * spacing(idim)
384 end do
385 end if
386 end do
387
388
389 !! Fourier space coordinates
390 if (fft_library /= fftlib_none) then
391
392 safe_allocate(cube%Lfs(1:maxn, 1:3))
393 cube%Lfs(:,:) = m_zero
394
395 do idim = 1,3
396 temp = m_two * m_pi / (nn(idim) * spacing(idim))
397!temp = M_PI / (nn * spacing(1))
398 do ii = 1, nn(idim)
399 if (fft_library == fftlib_nfft .or. fft_library == fftlib_pnfft) then
400 !The Fourier space is shrunk by the tp_enlarge factor
401 !cube%Lfs(ii, 1:3) = (ii - nn/2 - 1)*temp/tp_enlarge
402!HH NOTE:
403!not sure this is the right general factor
404 cube%Lfs(ii, idim) = (ii - nn(idim)/2 - 1)*temp/tp_enlarge(idim)
405 else
406 cube%Lfs(ii, idim) = pad_feq(ii,nn(idim), .true.) * temp
407 end if
408 end do
409 end do
410 end if
412 pop_sub(cube_init_coords)
413 end subroutine cube_init_coords
414
415
416 ! ---------------------------------------------------------
419 logical function cube_global2local(cube, ixyz, lxyz) result(is_here)
420 type(cube_t), intent(in) :: cube
421 integer, intent(in) :: ixyz(3)
422 integer, intent(out) :: lxyz(3)
423
424 lxyz(1) = ixyz(1) - cube%rs_istart(1) + 1
425 lxyz(2) = ixyz(2) - cube%rs_istart(2) + 1
426 lxyz(3) = ixyz(3) - cube%rs_istart(3) + 1
427 is_here = lxyz(1) >= 1 .and. lxyz(1) <= cube%rs_n(1) .and. &
428 lxyz(2) >= 1 .and. lxyz(2) <= cube%rs_n(2) .and. &
429 lxyz(3) >= 1 .and. lxyz(3) <= cube%rs_n(3)
430
431 end function cube_global2local
432
433
434 ! ---------------------------------------------------------
439 integer function cube_getfftlibrary(cube) result(fft_library)
440 type(cube_t), intent(in) :: cube
442 if (allocated(cube%fft)) then
443 fft_library = cube%fft%library
444 else
445 fft_library = fftlib_none
446 end if
447 end function cube_getfftlibrary
448
449 ! ---------------------------------------------------------
451 subroutine cube_do_mapping(cube, fs)
452 type(cube_t), intent(inout) :: cube
453 logical, intent(in) :: fs
454
455 integer :: tmp_local(6), position, process, ix, iy, iz, index
456 integer, allocatable :: local_sizes(:)
457 integer(int64) :: number_points
458
459 push_sub(cube_do_mapping)
460
461 !!BEGIN:gather the local information into a unique vector.
462 !!do a gather in 3d of all the box, into a loop
463 tmp_local(1) = cube%rs_istart(1)
464 tmp_local(2) = cube%rs_istart(2)
465 tmp_local(3) = cube%rs_istart(3)
466 tmp_local(4) = cube%rs_n(1)
467 tmp_local(5) = cube%rs_n(2)
468 tmp_local(6) = cube%rs_n(3)
469
470 if (cube%parallel_in_domains) then
471 safe_allocate(local_sizes(1:6*cube%mpi_grp%size))
472 call profiling_in("CUBE_GAT")
473 call cube%mpi_grp%allgather(tmp_local, 6, mpi_integer, local_sizes, 6, mpi_integer)
474 call profiling_out("CUBE_GAT")
475 else
476 safe_allocate(local_sizes(1:6))
477 local_sizes = tmp_local
478 end if
479
480 call profiling_in("CUBE_MAP")
481
482 safe_allocate(cube%xlocal(1:cube%mpi_grp%size))
483 safe_allocate(cube%np_local(1:cube%mpi_grp%size))
484 ! make sure we do not run into integer overflow here
485 number_points = cube%rs_n_global(1) * cube%rs_n_global(2)
486 number_points = number_points * cube%rs_n_global(3)
487 if (number_points >= huge(0)) then
488 message(1) = "Error: too many points for the normal cube. Please try to use a distributed FFT."
489 call messages_fatal(1)
490 end if
491 safe_allocate(cube%local(1:cube%rs_n_global(1)*cube%rs_n_global(2)*cube%rs_n_global(3), 1:3))
492
493 index = 1
494 do process = 1, cube%mpi_grp%size
495 position = ((process-1)*6)+1
496 if (position == 1) then
497 cube%xlocal(1) = 1
498 cube%np_local(1) = local_sizes(4)*local_sizes(5)*local_sizes(6)
499 else
500 ! calculate the begin index and size of each process
501 cube%xlocal(process) = cube%xlocal(process-1) + cube%np_local(process-1)
502 cube%np_local(process) = local_sizes(position+3)*local_sizes(position+4)*local_sizes(position+5)
503 end if
504
505 ! save the mapping between the global x,y,z and the global index
506 ! and determine which partition the point belongs to
507 do iz = local_sizes(position+2), local_sizes(position+2)+local_sizes(position+5)-1
508 do iy = local_sizes(position+1), local_sizes(position+1)+local_sizes(position+4)-1
509 do ix = local_sizes(position), local_sizes(position)+local_sizes(position+3)-1
510 cube%local(index, 1) = ix
511 cube%local(index, 2) = iy
512 cube%local(index, 3) = iz
513 index = index + 1
514 end do
515 end do
516 end do
517 end do
518
519 call profiling_out("CUBE_MAP")
520
521 if (optional_default(fs,.false.)) then
522
523 tmp_local(1) = cube%fs_istart(1)
524 tmp_local(2) = cube%fs_istart(2)
525 tmp_local(3) = cube%fs_istart(3)
526 tmp_local(4) = cube%fs_n(1)
527 tmp_local(5) = cube%fs_n(2)
528 tmp_local(6) = cube%fs_n(3)
529
530 local_sizes = 0
531 if (cube%parallel_in_domains) then
532 call profiling_in("CUBE_GAT_FS")
533 call cube%mpi_grp%allgather(tmp_local, 6, mpi_integer, local_sizes, 6, mpi_integer)
534 call profiling_out("CUBE_GAT_FS")
535 else
536 local_sizes = tmp_local
537 end if
538
539 call profiling_in("CUBE_MAP_FS")
540
541 safe_allocate(cube%xlocal_fs(1:cube%mpi_grp%size))
542 safe_allocate(cube%np_local_fs(1:cube%mpi_grp%size))
543 ! make sure we do not run into integer overflow here
544 number_points = cube%fs_n_global(1) * cube%fs_n_global(2)
545 number_points = number_points * cube%fs_n_global(3)
546 if (number_points >= huge(0)) then
547 message(1) = "Error: too many points for the normal cube. Please try to use a distributed FFT."
548 call messages_fatal(1)
549 end if
550 safe_allocate(cube%local_fs(1:cube%fs_n_global(1)*cube%fs_n_global(2)*cube%fs_n_global(3), 1:3))
551
552 index = 1
553 do process = 1, cube%mpi_grp%size
554 position = ((process-1)*6)+1
555 if (position == 1) then
556 cube%xlocal_fs(1) = 1
557 cube%np_local_fs(1) = local_sizes(4)*local_sizes(5)*local_sizes(6)
558 else
559 ! calculate the begin index and size of each process
560 cube%xlocal_fs(process) = cube%xlocal_fs(process-1) + cube%np_local_fs(process-1)
561 cube%np_local_fs(process) = local_sizes(position+3)*local_sizes(position+4)*local_sizes(position+5)
562 end if
563
564 ! save the mapping between the global x,y,z and the global index
565 ! and determine which partition the point belongs to
566 do iz = local_sizes(position+2), local_sizes(position+2)+local_sizes(position+5)-1
567 do iy = local_sizes(position+1), local_sizes(position+1)+local_sizes(position+4)-1
568 do ix = local_sizes(position), local_sizes(position)+local_sizes(position+3)-1
569 cube%local_fs(index, 1) = ix
570 cube%local_fs(index, 2) = iy
571 cube%local_fs(index, 3) = iz
572 index = index + 1
573 end do
574 end do
575 end do
576 end do
577
578 call profiling_out("CUBE_MAP_FS")
579
580 end if
581
582
583
584 safe_deallocate_a(local_sizes)
585
586 pop_sub(cube_do_mapping)
587 end subroutine cube_do_mapping
588
589 !!> Given a x, y, z point of the cube, it returns the corresponding process
590 !!
591 !! last_found is used to speed-up the search
592 integer pure function cube_point_to_process(xyz, part) result(process)
593 integer, intent(in) :: xyz(1:3)
594 type(dimensions_t), intent(in) :: part(:)
595
596 integer :: proc
597 logical :: found
598
599 ! No PUSH/POP because it is a PURE function
600
601 found = .false.
602 do proc = 1, mpi_world%size
603 !Compare XYZ index
604 if (all(xyz >= part(proc)%start_xyz) .and. all(xyz <= part(proc)%end_xyz)) then
605 process = proc
606 found = .true.
607 exit
608 end if
609 end do
610
611 ! An error message should be raised, if this point is reached
612 if (.not. found) then
613 process = -1
614 end if
615
616 end function cube_point_to_process
617
618 ! Sets a 1D decomposition with fixed-size blocks over the last (least-contiguous) axis.
619 ! Each core will have <blocksize> slices except the last one which will typically have
620 ! less. (In some cases, there can be multiple trailing cores without any slices.)
621 subroutine cube_set_blocksize(rs_n_global, blocksize, rank, rs_n, rs_istart)
622 integer, intent(in) :: rs_n_global(1:3)
623 integer, intent(in) :: blocksize
624 integer, intent(in) :: rank
625 integer, intent(out) :: rs_n(1:3)
626 integer, intent(out) :: rs_istart(1:3)
627
628 integer :: imin, imax
629
630 rs_n = rs_n_global
631 rs_istart = 1
632
633 imin = min(blocksize * rank, rs_n_global(3))
634 imax = min(imin + blocksize, rs_n_global(3))
635 rs_istart(3) = 1 + imin
636 rs_n(3) = imax - imin
637 end subroutine cube_set_blocksize
638
639 ! ---------------------------------------------------------
640 subroutine cube_partition(cube, part)
641 type(cube_t), intent(in) :: cube
642 type(dimensions_t), intent(out) :: part(:)
643
644 integer :: tmp_local(6), position, process
645 integer, allocatable :: local_sizes(:)
646
647 push_sub(cube_partition)
648
649 !!gather the local information into a unique vector.
650 tmp_local(1) = cube%rs_istart(1)
651 tmp_local(2) = cube%rs_istart(2)
652 tmp_local(3) = cube%rs_istart(3)
653 tmp_local(4) = cube%rs_n(1)
654 tmp_local(5) = cube%rs_n(2)
655 tmp_local(6) = cube%rs_n(3)
656
657 if (cube%parallel_in_domains) then
658 safe_allocate(local_sizes(1:6*cube%mpi_grp%size))
659 call cube%mpi_grp%allgather(tmp_local, 6, mpi_integer, local_sizes, 6, mpi_integer)
660 else
661 safe_allocate(local_sizes(1:6))
662 local_sizes = tmp_local
663 end if
664
665 do process = 1, cube%mpi_grp%size
666 position = ((process-1)*6)+1
667
668 part(process)%start_xyz(1) = local_sizes(position)
669 part(process)%start_xyz(2) = local_sizes(position+1)
670 part(process)%start_xyz(3) = local_sizes(position+2)
671 part(process)%end_xyz(1) = local_sizes(position)+local_sizes(position+3)-1
672 part(process)%end_xyz(2) = local_sizes(position+1)+local_sizes(position+4)-1
673 part(process)%end_xyz(3) = local_sizes(position+2)+local_sizes(position+5)-1
674
675 end do
676
677 pop_sub(cube_partition)
678 end subroutine cube_partition
679
680 ! ---------------------------------------------------------
681 subroutine cube_partition_messages_debug(cube, namespace)
682 type(cube_t), intent(in) :: cube
683 type(namespace_t), intent(in) :: namespace
684
685 integer :: nn, ii, jj, kk ! Counters.
686 integer :: ixyz(3) ! Current value of xyz
687 integer :: npart
688 integer :: iunit ! For debug output to files.
689 character(len=3) :: filenum
690 type(dimensions_t), allocatable :: part(:)
691
693
694 if (debug%info) then
695 safe_allocate(part(1:cube%mpi_grp%size))
696 call cube_partition(cube, part)
697
698 if (mpi_grp_is_root(mpi_world)) then
699 call io_mkdir('debug/cube_partition', namespace)
700 npart = cube%mpi_grp%size
701
702 ! Debug output. Write points of each partition in a different file.
703 do nn = 1, npart
704
705 write(filenum, '(i3.3)') nn
706
707 iunit = io_open('debug/cube_partition/cube_partition.'//filenum, &
708 namespace, action='write')
709 do kk = 1, cube%rs_n_global(3)
710 do jj = 1, cube%rs_n_global(2)
711 do ii = 1, cube%rs_n_global(1)
712 ixyz(1) = ii
713 ixyz(2) = jj
714 ixyz(3) = kk
715 if (cube_point_to_process(ixyz, part) == nn) then
716 write(iunit, '(3i8)') ii, jj, kk
717 end if
718 end do
719 end do
720 end do
721 call io_close(iunit)
722 end do
723
724
725 end if
727 safe_deallocate_a(part)
728 end if
729
730 call mpi_world%barrier()
731
733 end subroutine cube_partition_messages_debug
734
735 ! ---------------------------------------------------------
736 subroutine cube_init_cube_map(cube, mesh)
737 type(cube_t), intent(inout) :: cube
738 class(mesh_t), intent(in) :: mesh
739
740 push_sub(cube_init_cube_map)
741
742 call mesh_cube_map_init(cube%cube_map, mesh, mesh%np)
743 cube%cube_map_present = .true.
744
745 pop_sub(cube_init_cube_map)
746 end subroutine cube_init_cube_map
747end module cube_oct_m
748
749
750!! Local Variables:
751!! mode: f90
752!! coding: utf-8
753!! End:
subroutine cube_set_blocksize(rs_n_global, blocksize, rank, rs_n, rs_istart)
Definition: cube.F90:708
subroutine, public cube_init(cube, nn, namespace, space, spacing, coord_system, fft_type, fft_library, dont_optimize, nn_out, mpi_grp, need_partition, tp_enlarge, blocksize)
Definition: cube.F90:202
subroutine cube_do_mapping(cube, fs)
do the mapping between global and local points of the cube
Definition: cube.F90:538
subroutine, public cube_end(cube)
Definition: cube.F90:378
logical function, public cube_global2local(cube, ixyz, lxyz)
True if global coordinates belong to this process. On output lxyz contains the local coordinates.
Definition: cube.F90:506
subroutine cube_tp_fft_defaults(cube, fft_library)
Definition: cube.F90:412
integer function, public cube_getfftlibrary(cube)
Returns the FFT library of the cube. Possible values are FFTLIB_NONE, FFTLIB_FFTW,...
Definition: cube.F90:526
subroutine cube_init_coords(cube, tp_enlarge, spacing, fft_library)
Definition: cube.F90:442
subroutine, public cube_partition(cube, part)
Definition: cube.F90:727
integer pure function, public cube_point_to_process(xyz, part)
Definition: cube.F90:679
subroutine, public cube_init_cube_map(cube, mesh)
Definition: cube.F90:823
subroutine cube_partition_messages_debug(cube, namespace)
Definition: cube.F90:768
Fast Fourier Transform module. This module provides a single interface that works with different FFT ...
Definition: fft.F90:118
subroutine, public fft_init(this, nn, dim, type, library, optimize, optimize_parity, comm, mpi_grp, use_aligned)
Definition: fft.F90:417
integer, parameter, public fft_none
global constants
Definition: fft.F90:177
subroutine, public fft_end(this)
Definition: fft.F90:995
integer, public fft_default_lib
Definition: fft.F90:258
pure integer function, public pad_feq(ii, nn, mode)
convert between array index and G-vector
Definition: fft.F90:1105
subroutine, public fft_get_dims(fft, rs_n_global, fs_n_global, rs_n, fs_n, rs_istart, fs_istart)
Definition: fft.F90:1079
integer, parameter, public fftlib_nfft
Definition: fft.F90:182
integer, parameter, public fftlib_none
Definition: fft.F90:182
integer, parameter, public fftlib_pnfft
Definition: fft.F90:182
integer, parameter, public fftlib_pfft
Definition: fft.F90:182
subroutine, public fft_init_stage1(this, namespace, XX, nn)
Some fft-libraries (only NFFT for the moment) need an additional precomputation stage that depends on...
Definition: fft.F90:953
real(real64), parameter, public m_two
Definition: global.F90:189
real(real64), parameter, public m_zero
Definition: global.F90:187
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:185
real(real64), parameter, public m_one
Definition: global.F90:188
Definition: io.F90:114
subroutine, public mesh_cube_map_end(this)
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:420
type(mpi_comm), parameter, public mpi_comm_undefined
used to indicate a communicator has not been initialized
Definition: mpi.F90:136
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:266
subroutine mpi_grp_init(grp, comm)
Initialize MPI group instance.
Definition: mpi.F90:346
integer, parameter, public nfft_pre_psi
Definition: nfft.F90:148
The low level module to work with the PFFT library. http:
Definition: pfft.F90:164
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
It is intended to be used within a vector.
Definition: cube.F90:192
int true(void)