Octopus
mesh.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2011 M. Marques, A. Castro, A. Rubio, G. Bertsch, M. Oliveira
2!! Copyright (C) 2021 S. Ohlmann
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
24
25module mesh_oct_m
27 use box_oct_m
28 use comm_oct_m
30 use debug_oct_m
31 use global_oct_m
32 use iihash_oct_m
33 use index_oct_m
34 use io_oct_m
37 use mpi_oct_m
42 use parser_oct_m
44 use space_oct_m
47 use unit_oct_m
49
50 implicit none
51
52 private
53 public :: &
54 mesh_t, &
58 mesh_end, &
60 mesh_r, &
80
92 type, extends(basis_set_abst_t) :: mesh_t
93 ! Components are public by default
94 class(box_t), pointer :: box
95 class(coordinate_system_t), pointer :: coord_system
96 type(index_t) :: idx
97 logical :: use_curvilinear
98
99 real(real64), allocatable :: spacing(:)
100
101 ! When running serially, the local number of points is
102 ! equal to the global number of points.
103 ! Otherwise, the next two are different on each node.
104 integer :: np
105 integer :: np_part
106 integer(int64) :: np_global
107 integer(int64) :: np_part_global
108 logical :: parallel_in_domains
109 type(mpi_grp_t) :: mpi_grp
110 type(par_vec_t) :: pv
111 type(partition_t) :: partition
112
113 real(real64), allocatable :: x(:,:)
114 real(real64) :: volume_element
115 real(real64), allocatable :: vol_pp(:)
116
117 logical :: masked_periodic_boundaries
118 character(len=256) :: periodic_boundary_mask
120 contains
121 procedure :: end => mesh_end
122 procedure :: init => mesh_init
123 procedure :: write_info => mesh_write_info
124 procedure :: dmesh_allreduce_0, zmesh_allreduce_0, imesh_allreduce_0
125 procedure :: dmesh_allreduce_1, zmesh_allreduce_1, imesh_allreduce_1
126 procedure :: dmesh_allreduce_2, zmesh_allreduce_2, imesh_allreduce_2
127 procedure :: dmesh_allreduce_3, zmesh_allreduce_3, imesh_allreduce_3
128 procedure :: dmesh_allreduce_4, zmesh_allreduce_4, imesh_allreduce_4
129 procedure :: dmesh_allreduce_5, zmesh_allreduce_5, imesh_allreduce_5
130 generic :: allreduce => dmesh_allreduce_0, zmesh_allreduce_0, imesh_allreduce_0
131 generic :: allreduce => dmesh_allreduce_1, zmesh_allreduce_1, imesh_allreduce_1
132 generic :: allreduce => dmesh_allreduce_2, zmesh_allreduce_2, imesh_allreduce_2
133 generic :: allreduce => dmesh_allreduce_3, zmesh_allreduce_3, imesh_allreduce_3
134 generic :: allreduce => dmesh_allreduce_4, zmesh_allreduce_4, imesh_allreduce_4
135 generic :: allreduce => dmesh_allreduce_5, zmesh_allreduce_5, imesh_allreduce_5
136 end type mesh_t
137
149 type mesh_plane_t
150 ! Components are public by default
151 real(real64) :: n(3)
152 real(real64) :: u(3), v(3)
153 real(real64) :: origin(3)
154 real(real64) :: spacing
155 integer :: nu, mu, nv, mv
156 end type mesh_plane_t
157
161 type mesh_line_t
162 ! Components are public by default
163 real(real64) :: n(2)
164 real(real64) :: u(2)
165 real(real64) :: origin(2)
166 real(real64) :: spacing
167 integer :: nu, mu
168 end type mesh_line_t
169
170contains
171
172 subroutine mesh_init(this)
173 class(mesh_t), intent(inout) :: this
174
175 push_sub(mesh_init)
176
177 call this%set_time_dependent(.false.)
178
179 pop_sub(mesh_init)
180 end subroutine mesh_init
181
182! ---------------------------------------------------------
184 subroutine mesh_double_box(space, mesh, alpha, db)
185 class(space_t), intent(in) :: space
186 type(mesh_t), intent(in) :: mesh
187 real(real64), intent(in) :: alpha
188 integer, intent(out) :: db(:)
190 integer :: idir
191
193
194 db = 1
195
196 ! double mesh with 2n points
197 do idir = 1, space%periodic_dim
198 db(idir) = mesh%idx%ll(idir)
199 end do
200 do idir = space%periodic_dim + 1, space%dim
201 db(idir) = nint(alpha * (mesh%idx%ll(idir) - 1)) + 1
202 end do
205 end subroutine mesh_double_box
208 ! ---------------------------------------------------------
209 subroutine mesh_write_info(this, iunit, namespace)
210 class(mesh_t), intent(in) :: this
211 integer, optional, intent(in) :: iunit
212 type(namespace_t), optional, intent(in) :: namespace
213
214 integer :: ii
215 real(real64) :: cutoff
219 write(message(1),'(3a)') ' Spacing [', trim(units_abbrev(units_out%length)), '] = ('
220 do ii = 1, this%box%dim
221 if (ii > 1) write(message(1), '(2a)') trim(message(1)), ','
222 write(message(1), '(a,f6.3)') trim(message(1)), units_from_atomic(units_out%length, this%spacing(ii))
223 end do
224 write(message(1), '(5a,f12.5)') trim(message(1)), ') ', &
225 ' volume/point [', trim(units_abbrev(units_out%length**this%box%dim)), '] = ', &
226 units_from_atomic(units_out%length**this%box%dim, this%vol_pp(1))
228 write(message(2),'(a, i10)') ' # inner mesh = ', this%np_global
229 write(message(3),'(a, i10)') ' # total mesh = ', this%np_part_global
230
231 cutoff = mesh_gcutoff(this)**2 / m_two
232 write(message(4),'(3a,f12.6,a,f12.6)') ' Grid Cutoff [', trim(units_abbrev(units_out%energy)),'] = ', &
233 units_from_atomic(units_out%energy, cutoff), ' Grid Cutoff [Ry] = ', cutoff * m_two
234 call messages_info(4, iunit=iunit, namespace=namespace)
235
236 pop_sub(mesh_write_info)
237 end subroutine mesh_write_info
238
239
241 pure subroutine mesh_r(mesh, ip, rr, origin, coords)
242 class(mesh_t), intent(in) :: mesh
243 integer, intent(in) :: ip
244 real(real64), intent(out) :: rr
245 real(real64), intent(in), optional :: origin(:)
246 real(real64), intent(out), optional :: coords(:)
248 real(real64) :: xx(1:mesh%box%dim)
249
250 xx = mesh%x(ip, :)
251 if (present(origin)) xx = xx - origin
252 rr = norm2(xx)
253
254 if (present(coords)) then
255 coords = xx
256 end if
258 end subroutine mesh_r
265 integer function mesh_nearest_point(mesh, pos, dmin, rankmin) result(ind)
266 class(mesh_t),intent(in) :: mesh
267 real(real64), intent(in) :: pos(:)
268 real(real64), intent(out) :: dmin
269 integer, intent(out) :: rankmin
270 ! This returns 0 in the absence of domain decomposition
271
272 real(real64) :: dd
273 integer :: imin, ip
274
275 push_sub(mesh_nearest_point)
276
277 ! find the point of the grid that is closer to the atom
278 dmin = m_zero
279 do ip = 1, mesh%np
280 dd = sum((pos - mesh%x(ip, :))**2)
281 if ((dd < dmin) .or. (ip == 1)) then
282 imin = ip
283 dmin = dd
284 end if
285 end do
286
287 call mesh_minmaxloc(mesh, dmin, rankmin, mpi_minloc)
288 call mesh%mpi_grp%bcast(imin, 1, mpi_integer, rankmin)
289
290 ind = imin
291 pop_sub(mesh_nearest_point)
292 end function mesh_nearest_point
293
295 subroutine mesh_discretize_values_to_mesh(mesh, values)
296 class(mesh_t), intent(in ) :: mesh
297 real(real64), intent(inout) :: values(:, :)
298 ! Out: Values discretized to mesh points
299 integer :: i, ip, ndim
300 integer :: process
301 real(real64) :: dummy
304
305 if (mesh%parallel_in_domains) then
306 ndim = size(values, 1)
307 do i = 1, size(values, 2)
308 ip = mesh_nearest_point(mesh, values(:, i), dummy, process)
309 ! ip is defined for `process`, only
310 if (mesh%mpi_grp%rank == process) values(:, i) = mesh%x(ip, :)
311 call mesh%mpi_grp%bcast(values(:, i), ndim, mpi_double_precision, process)
312 enddo
313 else
314 do i = 1, size(values, 2)
315 ip = mesh_nearest_point(mesh, values(:, i), dummy, process)
316 values(:, i) = mesh%x(ip, :)
317 enddo
318 endif
319
321
322 end subroutine mesh_discretize_values_to_mesh
323
324 ! --------------------------------------------------------------
328 ! --------------------------------------------------------------
329 real(real64) function mesh_gcutoff(mesh) result(gmax)
330 class(mesh_t), intent(in) :: mesh
331
332 push_sub(mesh_gcutoff)
333 gmax = m_pi / (maxval(mesh%spacing))
335 pop_sub(mesh_gcutoff)
336 end function mesh_gcutoff
337
338 ! --------------------------------------------------------------
339 subroutine mesh_write_fingerprint(mesh, dir, filename, mpi_grp, namespace, ierr)
340 type(mesh_t), intent(in) :: mesh
341 character(len=*), intent(in) :: dir
342 character(len=*), intent(in) :: filename
343 type(mpi_grp_t), intent(in) :: mpi_grp
344 type(namespace_t),intent(in) :: namespace
345 integer, intent(out) :: ierr
346
347 integer :: iunit, ii
348
349 push_sub(mesh_write_fingerprint)
350
351 ierr = 0
352
353 iunit = io_open(trim(dir)//"/"//trim(filename), namespace, action='write', &
354 die=.false., grp=mpi_grp)
355 if (iunit <= 0) then
356 message(1) = "Unable to open file '"//trim(dir)//"/"//trim(filename)//"'."
357 call messages_warning(1, namespace=namespace)
358 ierr = ierr + 1
359 else
360 if (mpi_grp_is_root(mpi_grp)) then
361 write(iunit, '(a20,i21)') 'np_part_global= ', mesh%np_part_global
362 write(iunit, '(a20,i21)') 'np_global= ', mesh%np_global
363 write(iunit, '(a20,i21)') 'algorithm= ', 1
364 write(iunit, '(a20,i21)') 'checksum= ', mesh%idx%checksum
365 write(iunit, '(a20,i21)') 'bits= ', mesh%idx%bits
366 write(iunit, '(a20,i21)') 'dim= ', mesh%idx%dim
367 write(iunit, '(a20,i21)') 'type= ', mesh%idx%type
368 do ii = 1, mesh%idx%dim
369 write(iunit, '(a7,i2,a11,i21)') 'offset(',ii,')= ', mesh%idx%offset(ii)
370 end do
371 do ii = 1, mesh%idx%dim
372 write(iunit, '(a7,i2,a11,i21)') 'nn(',ii,')= ', mesh%idx%nr(2, ii) - mesh%idx%nr(1, ii) + 1
373 end do
374 end if
375 call io_close(iunit, grp=mpi_grp)
376 end if
377
379 end subroutine mesh_write_fingerprint
380
381
382 ! -----------------------------------------------------------------------
387 subroutine mesh_read_fingerprint(mesh, dir, filename, mpi_grp, namespace, &
388 read_np_part, read_np, bits, type, offset, nn, ierr)
389 type(mesh_t), intent(in) :: mesh
390 character(len=*), intent(in) :: dir
391 character(len=*), intent(in) :: filename
392 type(mpi_grp_t), intent(in) :: mpi_grp
393 type(namespace_t),intent(in) :: namespace
394 integer(int64), intent(out) :: read_np_part
395 integer(int64), intent(out) :: read_np
396 integer, intent(out) :: bits
397 integer, intent(out) :: type
398 integer, intent(out) :: offset(1:mesh%idx%dim)
399 integer, intent(out) :: nn(1:mesh%idx%dim)
400 integer, intent(out) :: ierr
401
402 character(len=20) :: str
403 character(len=100) :: lines(7)
404 integer :: iunit, algorithm, dim, err, ii
405 integer(int64) :: checksum
406
407 push_sub(mesh_read_fingerprint)
408
409 ierr = 0
410
411 read_np_part = 0_int64
412 read_np = 0_int64
413
414 iunit = io_open(trim(dir)//"/"//trim(filename), namespace, action='read', &
415 status='old', die=.false., grp=mpi_grp)
416 if (iunit <= 0) then
417 ierr = ierr + 1
418 message(1) = "Unable to open file '"//trim(dir)//"/"//trim(filename)//"'."
419 call messages_warning(1, namespace=namespace)
420 else
421 call iopar_read(mpi_grp, iunit, lines, 7, err)
422 if (err /= 0) then
423 ierr = ierr + 4
424 else
425 read(lines(1), '(a20,i21)') str, read_np_part
426 read(lines(2), '(a20,i21)') str, read_np
427 read(lines(3), '(a20,i21)') str, algorithm
428 read(lines(4), '(a20,i21)') str, checksum
429 read(lines(5), '(a20,i21)') str, bits
430 read(lines(6), '(a20,i21)') str, dim
431 read(lines(7), '(a20,i21)') str, type
432 ! only allow restarting simulations with the same dimensions
433 if (dim /= mesh%idx%dim) then
434 ierr = ierr + 8
435 else
436 ! read offset, has dim lines
437 call iopar_read(mpi_grp, iunit, lines, dim, err)
438 if (err /= 0) then
439 ierr = ierr + 4
440 else
441 do ii = 1, dim
442 read(lines(ii), '(a20,i21)') str, offset(ii)
443 end do
444 end if
445
446 ! read nn, has dim lines
447 call iopar_read(mpi_grp, iunit, lines, dim, err)
448 if (err /= 0) then
449 ierr = ierr + 4
450 else
451 do ii = 1, dim
452 read(lines(ii), '(a20,i21)') str, nn(ii)
453 end do
454 end if
455 end if
456
457 assert(read_np_part >= read_np)
458
459 if (read_np_part == mesh%np_part_global &
460 .and. read_np == mesh%np_global &
461 .and. algorithm == 1 &
462 .and. checksum == mesh%idx%checksum) then
463 read_np_part = 0
464 read_np = 0
465 end if
466 end if
467
468 call io_close(iunit, grp=mpi_grp)
469 end if
470
471 pop_sub(mesh_read_fingerprint)
472 end subroutine mesh_read_fingerprint
473
474 ! ---------------------------------------------------------
475 subroutine mesh_check_dump_compatibility(mesh, dir, filename, namespace, mpi_grp, grid_changed, grid_reordered, map, ierr)
476 type(mesh_t), intent(in) :: mesh
477 character(len=*), intent(in) :: dir
478 character(len=*), intent(in) :: filename
479 type(namespace_t), intent(in) :: namespace
480 type(mpi_grp_t), intent(in) :: mpi_grp
481 logical, intent(out) :: grid_changed
482 logical, intent(out) :: grid_reordered
483 integer(int64), allocatable, intent(out) :: map(:)
484 integer, intent(out) :: ierr
485
486 integer(int64) :: ipg, ipg_new, read_np_part, read_np
487 integer :: err, idir
488 integer :: bits, type, offset(mesh%idx%dim), point(mesh%idx%dim), nn(mesh%idx%dim)
489 integer(int64), allocatable :: read_indices(:)
490 type(index_t) :: idx_old
491
493
494 ierr = 0
495
496 grid_changed = .false.
497 grid_reordered = .false.
498
499 ! Read the mesh fingerprint
500 call mesh_read_fingerprint(mesh, dir, filename, mpi_grp, namespace, read_np_part, read_np, &
501 bits, type, offset, nn, err)
502 if (err /= 0) then
503 ierr = ierr + 1
504 message(1) = "Unable to read mesh fingerprint from '"//trim(dir)//"/"//trim(filename)//"'."
505 call messages_warning(1, namespace=namespace)
506
507 else if (read_np > 0) then
508 if (.not. associated(mesh%box)) then
509 ! We can only check the compatibility of two meshes that have different fingerprints if we also
510 ! have the simulation box. In the case we do not, we will assume that the fingerprint is enough.
511 ierr = ierr + 2
512 else
513 grid_changed = .true.
514
515 ! perhaps only the order of the points changed, this can only
516 ! happen if the number of points is the same and no points maps
517 ! to zero (this is checked below)
518 grid_reordered = (read_np == mesh%np_global)
519
520 ! the grid is different, so we read the coordinates.
521 safe_allocate(read_indices(1:read_np_part))
522 call io_binary_read(trim(io_workpath(dir, namespace))//"/indices.obf", read_np_part, &
523 read_indices, err)
524 if (err /= 0) then
525 ierr = ierr + 4
526 message(1) = "Unable to read index map from '"//trim(dir)//"'."
527 call messages_warning(1, namespace=namespace)
528 else
529 ! dummy index object
530 call index_init(idx_old, mesh%idx%dim)
531 idx_old%type = type
532 idx_old%bits = bits
533 idx_old%nr(1, :) = -offset
534 idx_old%nr(2, :) = -offset + nn - 1
535 idx_old%offset = offset
536 idx_old%stride(1) = 1
537 do idir = 2, mesh%idx%dim
538 idx_old%stride(idir) = idx_old%stride(idir-1) * nn(idir-1)
539 end do
540 ! generate the map
541 safe_allocate(map(1:read_np))
542 do ipg = 1, read_np
543 ! get nd-index from old 1d index
544 call index_spatial_to_point(idx_old, mesh%idx%dim, read_indices(ipg), point)
545 ! get new global index
546 ipg_new = mesh_global_index_from_coords(mesh, point)
547 map(ipg) = ipg_new
548 ! ignore boundary points
549 if (map(ipg) > mesh%np_global) map(ipg) = 0
550 ! if the map is zero for one point, it is not a simple reordering
551 if (map(ipg) == 0) grid_reordered = .false.
552 end do
553 call index_end(idx_old)
554 end if
555
556 safe_deallocate_a(read_indices)
557 end if
558 end if
559
561 end subroutine mesh_check_dump_compatibility
562
563
564 ! --------------------------------------------------------------
565 recursive subroutine mesh_end(this)
566 class(mesh_t), intent(inout) :: this
567
568 push_sub(mesh_end)
569
570#ifdef HAVE_MPI
571 call lmpi_destroy_shared_memory_window(this%idx%window_grid_to_spatial)
572 call lmpi_destroy_shared_memory_window(this%idx%window_spatial_to_grid)
573 nullify(this%idx%grid_to_spatial_global)
574 nullify(this%idx%spatial_to_grid_global)
575#else
576 safe_deallocate_p(this%idx%grid_to_spatial_global)
577 safe_deallocate_p(this%idx%spatial_to_grid_global)
578#endif
579
580 safe_deallocate_a(this%x)
581 safe_deallocate_a(this%vol_pp)
582
583 if (this%parallel_in_domains) then
584 call par_vec_end(this%pv)
585 call partition_end(this%partition)
586 end if
587
588 call index_end(this%idx)
589 safe_deallocate_a(this%spacing)
590
591 pop_sub(mesh_end)
592 end subroutine mesh_end
593
594
600 ! ---------------------------------------------------------
601 integer(int64) function mesh_periodic_point(mesh, space, ip) result(ipg)
602 class(mesh_t), intent(in) :: mesh
603 class(space_t),intent(in) :: space
604 integer, intent(in) :: ip
605
606 integer :: ix(space%dim), nr(2, space%dim), idim
607 real(real64) :: xx(space%dim), rr, ufn_re, ufn_im
608
609 ! no push_sub, called too frequently
610
611 call mesh_local_index_to_coords(mesh, ip, ix)
612 nr(1, :) = mesh%idx%nr(1, :) + mesh%idx%enlarge(:)
613 nr(2, :) = mesh%idx%nr(2, :) - mesh%idx%enlarge(:)
614
615 do idim = 1, space%periodic_dim
616 if (ix(idim) < nr(1, idim)) ix(idim) = ix(idim) + mesh%idx%ll(idim)
617 if (ix(idim) > nr(2, idim)) ix(idim) = ix(idim) - mesh%idx%ll(idim)
618 end do
619
620 ipg = mesh_global_index_from_coords(mesh, ix)
621 assert(ipg > 0)
622
623 if (mesh%masked_periodic_boundaries) then
624 call mesh_r(mesh, ip, rr, coords = xx)
625 call parse_expression(ufn_re, ufn_im, space%dim, xx, rr, m_zero, mesh%periodic_boundary_mask)
626 if (int(ufn_re) == 0) ipg = mesh_local2global(mesh, ip) ! Nothing will be done
627 end if
628
629 end function mesh_periodic_point
630
631
632 ! ---------------------------------------------------------
633 real(real64) pure function mesh_global_memory(mesh) result(memory)
634 use iso_c_binding, only: c_sizeof, c_long_long
635 class(mesh_t), intent(in) :: mesh
636
637 ! 2 global index arrays
638 memory = c_sizeof(c_long_long) * real(mesh%np_part_global, real64) * 2
639
640 end function mesh_global_memory
641
642
643 ! ---------------------------------------------------------
644 real(real64) pure function mesh_local_memory(mesh) result(memory)
645 use iso_c_binding, only: c_sizeof, c_long_long
646 class(mesh_t), intent(in) :: mesh
647
648 memory = m_zero
649
650 ! x
651 memory = memory + real64 * real(mesh%np_part, real64) * mesh%idx%dim
652 ! local index arrays
653 memory = memory + c_sizeof(c_long_long) * real(mesh%np_part, real64) * 2
654 end function mesh_local_memory
655
656
657 ! ---------------------------------------------------------
658 function mesh_x_global(mesh, ipg) result(xx)
659 class(mesh_t), intent(in) :: mesh
660 integer(int64), intent(in) :: ipg
661 real(real64) :: xx(1:mesh%box%dim)
662
663 real(real64) :: chi(1:mesh%box%dim)
664 integer :: ix(1:mesh%box%dim)
665
666! no push_sub because function is called too frequently
667
668 call mesh_global_index_to_coords(mesh, ipg, ix)
669 chi = ix * mesh%spacing
670 xx = mesh%coord_system%to_cartesian(chi)
671
672 end function mesh_x_global
673
674
675 ! ---------------------------------------------------------
676 logical pure function mesh_compact_boundaries(mesh) result(cb)
677 type(mesh_t), intent(in) :: mesh
678
679 cb = .not. mesh%use_curvilinear .and. &
680 .not. mesh%parallel_in_domains
681
682 end function mesh_compact_boundaries
683
684
685 ! ---------------------------------------------------------
686 subroutine mesh_check_symmetries(mesh, symm, periodic_dim)
687 class(mesh_t), intent(in) :: mesh
688 type(symmetries_t), intent(in) :: symm
689 integer, intent(in) :: periodic_dim
690
691 integer :: iop, ip, idim, nops, ix(1:3)
692 real(real64) :: destpoint(1:3), srcpoint(1:3), lsize(1:3), offset(1:3)
693 real(real64), parameter :: tol_spacing = 1e-12_real64
695 !If all the axis have the same spacing and the same length
696 !the grid is by obviously symmetric
697 !Indeed, reduced coordinates are proportional to the point index
698 !and the reduced rotation are integer matrices
699 !The result of the product is also proportional to an integer
700 !and therefore belong to the grid.
701 if (mesh%idx%ll(1) == mesh%idx%ll(2) .and. &
702 mesh%idx%ll(2) == mesh%idx%ll(3) .and. &
703 abs(mesh%spacing(1) - mesh%spacing(2)) < tol_spacing .and. &
704 abs(mesh%spacing(2) - mesh%spacing(3)) < tol_spacing) return
705
706 push_sub(mesh_check_symmetries)
707
708 message(1) = "Checking if the real-space grid is symmetric"
709 call messages_info(1)
710
711 lsize(1:3) = real(mesh%idx%ll(1:3), real64)
712 offset(1:3) = real(mesh%idx%nr(1, 1:3) + mesh%idx%enlarge(1:3), real64)
713
714 nops = symmetries_number(symm)
715
716 do ip = 1, mesh%np
717 !We use floating point coordinates to check if the symmetric point
718 !belong to the grid.
719 !If yes, it should have integer reduced coordinates
720 call mesh_local_index_to_coords(mesh, ip, ix)
721 destpoint(1:3) = real(ix(1:3), real64) - offset(1:3)
722 ! offset moves corner of cell to origin, in integer mesh coordinates
723 assert(all(destpoint >= 0))
724 assert(all(destpoint < lsize))
725
726 ! move to center of cell in real coordinates
727 destpoint = destpoint - real(int(lsize)/2, real64)
728
729 !convert to proper reduced coordinates
730 do idim = 1, 3
731 destpoint(idim) = destpoint(idim)/lsize(idim)
732 end do
733
734 ! iterate over all points that go to this point by a symmetry operation
735 do iop = 1, nops
736 srcpoint = symm_op_apply_red(symm%ops(iop), destpoint)
738 !We now come back to what should be an integer, if the symmetric point beloings to the grid
739 do idim = 1, 3
740 srcpoint(idim) = srcpoint(idim)*lsize(idim)
741 end do
742
743 ! move back to reference to origin at corner of cell
744 srcpoint = srcpoint + real(int(lsize)/2, real64)
745
746 ! apply periodic boundary conditions in periodic directions
747 do idim = 1, periodic_dim
748 if (nint(srcpoint(idim)) < 0 .or. nint(srcpoint(idim)) >= lsize(idim)) then
749 srcpoint(idim) = modulo(srcpoint(idim)+m_half*symprec, lsize(idim))
750 end if
751 end do
752 assert(all(srcpoint >= -symprec))
753 assert(all(srcpoint < lsize))
754
755 srcpoint(1:3) = srcpoint(1:3) + offset(1:3)
756
757 if (any(srcpoint-anint(srcpoint)> symprec*m_two)) then
758 message(1) = "The real-space grid breaks at least one of the symmetries of the system."
759 message(2) = "Change your spacing or use SymmetrizeDensity=no."
760 call messages_fatal(2)
761 end if
762 end do
763 end do
764
765 pop_sub(mesh_check_symmetries)
766 end subroutine
767
770 integer(int64) function mesh_global_index_from_coords(mesh, ix) result(index)
771 class(mesh_t), intent(in) :: mesh
772 integer, intent(in) :: ix(:)
773
774 index = index_from_coords(mesh%idx, ix)
776
779 subroutine mesh_global_index_to_coords(mesh, ipg, ix)
780 type(mesh_t), intent(in) :: mesh
781 integer(int64), intent(in) :: ipg
782 integer, intent(out) :: ix(:)
783
784 call index_to_coords(mesh%idx, ipg, ix)
785 end subroutine mesh_global_index_to_coords
786
789 integer function mesh_local_index_from_coords(mesh, ix) result(ip)
790 type(mesh_t), intent(in) :: mesh
791 integer, intent(in) :: ix(:)
792
793 integer(int64) :: ipg
794
795 ipg = index_from_coords(mesh%idx, ix)
796 ip = mesh_global2local(mesh, ipg)
798
801 subroutine mesh_local_index_to_coords(mesh, ip, ix)
802 type(mesh_t), intent(in) :: mesh
803 integer, intent(in) :: ip
804 integer, intent(out) :: ix(:)
805
806 integer(int64) :: ipg
807
808 ipg = mesh_local2global(mesh, ip)
809 call index_to_coords(mesh%idx, ipg, ix)
810 end subroutine mesh_local_index_to_coords
811
813 integer(int64) function mesh_local2global(mesh, ip) result(ipg)
814 type(mesh_t), intent(in) :: mesh
815 integer, intent(in) :: ip
816
817 ipg = par_vec_local2global(mesh%pv, ip)
818 end function mesh_local2global
819
821 integer function mesh_global2local(mesh, ipg) result(ip)
822 type(mesh_t), intent(in) :: mesh
823 integer(int64), intent(in) :: ipg
824
825 ip = par_vec_global2local(mesh%pv, ipg)
826 end function mesh_global2local
827
828 !-----------------------------------------------------------------------------
830 subroutine mesh_minmaxloc(mesh, min_or_max, rank_min_or_max, op)
831 type(mesh_t), intent(in) :: mesh
832 real(real64), intent(inout) :: min_or_max
833 integer, intent(out) :: rank_min_or_max
834 type(mpi_op), intent(in) :: op
835
836 real(real64) :: loc_in(2), loc_out(2)
837
838 push_sub(mesh_minmaxloc)
839
840 assert(op == mpi_minloc .or. op == mpi_maxloc)
841
842 rank_min_or_max = 0
843 if (mesh%parallel_in_domains) then
844 loc_in(1) = min_or_max
845 loc_in(2) = mesh%mpi_grp%rank
846 call mesh%mpi_grp%allreduce(loc_in, loc_out, 1, mpi_2double_precision, op)
847 min_or_max = loc_out(1)
848 rank_min_or_max = nint(loc_out(2))
849 end if
850
851 pop_sub(mesh_minmaxloc)
852 end subroutine mesh_minmaxloc
853
854
855#include "undef.F90"
856#include "real.F90"
857#include "mesh_inc.F90"
858
859#include "undef.F90"
860#include "complex.F90"
861#include "mesh_inc.F90"
862
863#include "undef.F90"
864#include "integer.F90"
865#include "mesh_inc.F90"
866
867end module mesh_oct_m
868
869
870!! Local Variables:
871!! mode: f90
872!! coding: utf-8
873!! End:
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
This module implements a simple hash table for non-negative integer keys and integer values.
Definition: iihash.F90:125
This module implements the index, used for the mesh points.
Definition: index.F90:122
Definition: io.F90:114
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
integer(int64) function, public mesh_periodic_point(mesh, space, ip)
This function returns the point inside the grid corresponding to a boundary point when PBCs are used....
Definition: mesh.F90:695
subroutine, public mesh_global_index_to_coords(mesh, ipg, ix)
Given a global point index, this function returns the set of integer coordinates of the point.
Definition: mesh.F90:873
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:883
subroutine, public mesh_double_box(space, mesh, alpha, db)
finds the dimension of a box doubled in the non-periodic dimensions
Definition: mesh.F90:278
integer(int64) function, public mesh_global_index_from_coords(mesh, ix)
This function returns the true global index of the point for a given vector of integer coordinates.
Definition: mesh.F90:864
subroutine, public mesh_check_symmetries(mesh, symm, periodic_dim)
Definition: mesh.F90:780
subroutine, public mesh_write_info(this, iunit, namespace)
Definition: mesh.F90:303
subroutine, public mesh_local_index_to_coords(mesh, ip, ix)
Given a local point index, this function returns the set of integer coordinates of the point.
Definition: mesh.F90:895
subroutine, public mesh_discretize_values_to_mesh(mesh, values)
Assign a set of values to their nearest discrete points on the mesh.
Definition: mesh.F90:389
integer function, public mesh_nearest_point(mesh, pos, dmin, rankmin)
Returns the index of the point which is nearest to a given vector position pos.
Definition: mesh.F90:359
real(real64) pure function, public mesh_global_memory(mesh)
Definition: mesh.F90:727
subroutine, public mesh_read_fingerprint(mesh, dir, filename, mpi_grp, namespace, read_np_part, read_np, bits, type, offset, nn, ierr)
This function reads the fingerprint of a mesh written in filename. If the meshes are equal (same fing...
Definition: mesh.F90:482
pure subroutine, public mesh_r(mesh, ip, rr, origin, coords)
return the distance to the origin for a given grid point
Definition: mesh.F90:335
subroutine, public mesh_check_dump_compatibility(mesh, dir, filename, namespace, mpi_grp, grid_changed, grid_reordered, map, ierr)
Definition: mesh.F90:569
real(real64) function, public mesh_gcutoff(mesh)
mesh_gcutoff returns the "natural" band limitation of the grid mesh, in terms of the maximum G vector...
Definition: mesh.F90:423
subroutine, public mesh_minmaxloc(mesh, min_or_max, rank_min_or_max, op)
Given a local min/max this returns the global min/max and the rank where this is located.
Definition: mesh.F90:924
integer function, public mesh_global2local(mesh, ipg)
This function returns the local mesh index for a given global index.
Definition: mesh.F90:915
recursive subroutine, public mesh_end(this)
Definition: mesh.F90:659
logical pure function, public mesh_compact_boundaries(mesh)
Definition: mesh.F90:770
real(real64) pure function, public mesh_local_memory(mesh)
Definition: mesh.F90:738
integer(int64) function, public mesh_local2global(mesh, ip)
This function returns the global mesh index for a given local index.
Definition: mesh.F90:907
real(real64) function, dimension(1:mesh%box%dim), public mesh_x_global(mesh, ipg)
Definition: mesh.F90:752
subroutine, public mesh_write_fingerprint(mesh, dir, filename, mpi_grp, namespace, ierr)
Definition: mesh.F90:433
subroutine mesh_init(this)
Definition: mesh.F90:266
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
Definition: messages.F90:624
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
This module contains some common usage patterns of MPI routines.
Definition: mpi_lib.F90:115
Some general things and nomenclature:
Definition: par_vec.F90:171
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
character(len=20) pure function, public units_abbrev(this)
Definition: unit.F90:223
This module defines the unit system, used for input and output.
type(unit_system_t), public units_out
abstract class for basis sets
This data type defines a line, and a regular grid defined on this line (or rather,...
Definition: mesh.F90:254
define a grid on a plane.
Definition: mesh.F90:242
Describes mesh distribution to nodes.
Definition: mesh.F90:185
int true(void)
real(real64) function values(xx)
Definition: test.F90:1862