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), allocatable :: x_t(:,:)
115 real(real64), allocatable :: chi(:,:)
116 real(real64) :: volume_element
117 real(real64), allocatable :: vol_pp(:)
118 real(real64), allocatable :: jacobian_inverse(:,:,:)
119
120 logical :: masked_periodic_boundaries
121 character(len=256) :: periodic_boundary_mask
123 contains
124 procedure :: end => mesh_end
125 procedure :: init => mesh_init
126 procedure :: write_info => mesh_write_info
127 procedure :: dmesh_allreduce_0, zmesh_allreduce_0, imesh_allreduce_0, lmesh_allreduce_0
128 procedure :: dmesh_allreduce_1, zmesh_allreduce_1, imesh_allreduce_1, lmesh_allreduce_1
129 procedure :: dmesh_allreduce_2, zmesh_allreduce_2, imesh_allreduce_2, lmesh_allreduce_2
130 procedure :: dmesh_allreduce_3, zmesh_allreduce_3, imesh_allreduce_3, lmesh_allreduce_3
131 procedure :: dmesh_allreduce_4, zmesh_allreduce_4, imesh_allreduce_4, lmesh_allreduce_4
132 procedure :: dmesh_allreduce_5, zmesh_allreduce_5, imesh_allreduce_5, lmesh_allreduce_5
133 procedure :: dmesh_allreduce_6, zmesh_allreduce_6, imesh_allreduce_6, lmesh_allreduce_6
134 generic :: allreduce => dmesh_allreduce_0, zmesh_allreduce_0, imesh_allreduce_0, lmesh_allreduce_0
135 generic :: allreduce => dmesh_allreduce_1, zmesh_allreduce_1, imesh_allreduce_1, lmesh_allreduce_1
136 generic :: allreduce => dmesh_allreduce_2, zmesh_allreduce_2, imesh_allreduce_2, lmesh_allreduce_2
137 generic :: allreduce => dmesh_allreduce_3, zmesh_allreduce_3, imesh_allreduce_3, lmesh_allreduce_3
138 generic :: allreduce => dmesh_allreduce_4, zmesh_allreduce_4, imesh_allreduce_4, lmesh_allreduce_4
139 generic :: allreduce => dmesh_allreduce_5, zmesh_allreduce_5, imesh_allreduce_5, lmesh_allreduce_5
140 generic :: allreduce => dmesh_allreduce_6, zmesh_allreduce_6, imesh_allreduce_6, lmesh_allreduce_6
141 end type mesh_t
142
154 type mesh_plane_t
155 ! Components are public by default
156 real(real64) :: n(3)
157 real(real64) :: u(3), v(3)
158 real(real64) :: origin(3)
159 real(real64) :: spacing
160 integer :: nu, mu, nv, mv
161 end type mesh_plane_t
162
166 type mesh_line_t
167 ! Components are public by default
168 real(real64) :: n(2)
169 real(real64) :: u(2)
170 real(real64) :: origin(2)
171 real(real64) :: spacing
172 integer :: nu, mu
173 end type mesh_line_t
174
175contains
176
177 subroutine mesh_init(this)
178 class(mesh_t), intent(inout) :: this
179
180 push_sub(mesh_init)
181
182 call this%set_time_dependent(.false.)
183
184 pop_sub(mesh_init)
185 end subroutine mesh_init
186
187! ---------------------------------------------------------
189 subroutine mesh_double_box(space, mesh, alpha, db)
190 class(space_t), intent(in) :: space
191 type(mesh_t), intent(in) :: mesh
192 real(real64), intent(in) :: alpha
193 integer, intent(out) :: db(:)
195 integer :: idir
196
197 push_sub(mesh_double_box)
198
199 db = 1
201 ! double mesh with 2n points
202 do idir = 1, space%periodic_dim
203 db(idir) = mesh%idx%ll(idir)
204 end do
205 do idir = space%periodic_dim + 1, space%dim
206 db(idir) = nint(alpha * (mesh%idx%ll(idir) - 1)) + 1
207 end do
210 end subroutine mesh_double_box
213 ! ---------------------------------------------------------
214 subroutine mesh_write_info(this, iunit, namespace)
215 class(mesh_t), intent(in) :: this
216 integer, optional, intent(in) :: iunit
217 type(namespace_t), optional, intent(in) :: namespace
218
219 integer :: ii
220 real(real64) :: cutoff
224 write(message(1),'(3a)') ' Spacing [', trim(units_abbrev(units_out%length)), '] = ('
225 do ii = 1, this%box%dim
226 if (ii > 1) write(message(1), '(2a)') trim(message(1)), ','
227 write(message(1), '(a,f6.3)') trim(message(1)), units_from_atomic(units_out%length, this%spacing(ii))
228 end do
229 write(message(1), '(5a,f12.5)') trim(message(1)), ') ', &
230 ' volume/point [', trim(units_abbrev(units_out%length**this%box%dim)), '] = ', &
231 units_from_atomic(units_out%length**this%box%dim, this%vol_pp(1))
233 write(message(2),'(a, i10)') ' # inner mesh = ', this%np_global
234 write(message(3),'(a, i10)') ' # total mesh = ', this%np_part_global
236 cutoff = mesh_gcutoff(this)**2 / m_two
237 write(message(4),'(3a,f12.6,a,f12.6)') ' Grid Cutoff [', trim(units_abbrev(units_out%energy)),'] = ', &
238 units_from_atomic(units_out%energy, cutoff), ' Grid Cutoff [Ry] = ', cutoff * m_two
239 call messages_info(4, iunit=iunit, namespace=namespace)
240
241 pop_sub(mesh_write_info)
242 end subroutine mesh_write_info
243
244
246 pure subroutine mesh_r(mesh, ip, rr, origin, coords)
247 class(mesh_t), intent(in) :: mesh
248 integer, intent(in) :: ip
249 real(real64), intent(out) :: rr
250 real(real64), intent(in), optional :: origin(:)
251 real(real64), intent(out), optional :: coords(:)
253 real(real64) :: xx(1:mesh%box%dim)
255 xx = mesh%x(:, ip)
256 if (present(origin)) xx = xx - origin
257 rr = norm2(xx)
258
259 if (present(coords)) then
260 coords = xx
261 end if
262
263 end subroutine mesh_r
266 subroutine mesh_r_global(mesh, ipg, rr, origin, coords)
267 class(mesh_t), intent(in) :: mesh
268 integer(int64),intent(in) :: ipg
269 real(real64), intent(out) :: rr
270 real(real64), intent(in), optional :: origin(:)
271 real(real64), intent(out), optional :: coords(:)
273 real(real64) :: xx(1:mesh%box%dim)
274
275 xx = mesh_x_global(mesh, ipg)
276 if (present(origin)) xx = xx - origin
277 rr = norm2(xx)
278
279 if (present(coords)) then
280 coords = xx
281 end if
282
283 end subroutine mesh_r_global
285
290 integer function mesh_nearest_point(mesh, pos, dmin, rankmin) result(ind)
291 class(mesh_t),intent(in) :: mesh
292 real(real64), intent(in) :: pos(:)
293 real(real64), intent(out) :: dmin
294 integer, intent(out) :: rankmin
295 ! This returns 0 in the absence of domain decomposition
296
297 real(real64) :: dd
298 integer :: imin, ip
299
300 push_sub(mesh_nearest_point)
301
302 ! find the point of the grid that is closer to the atom
303 dmin = m_zero
304 do ip = 1, mesh%np
305 dd = sum((pos - mesh%x(:, ip))**2)
306 if ((dd < dmin) .or. (ip == 1)) then
307 imin = ip
308 dmin = dd
309 end if
310 end do
311
312 call mesh_minmaxloc(mesh, dmin, rankmin, mpi_minloc)
313 call mesh%mpi_grp%bcast(imin, 1, mpi_integer, rankmin)
314
315 ind = imin
316 pop_sub(mesh_nearest_point)
317 end function mesh_nearest_point
318
320 subroutine mesh_discretize_values_to_mesh(mesh, values)
321 class(mesh_t), intent(in ) :: mesh
322 real(real64), intent(inout) :: values(:, :)
323 ! Out: Values discretized to mesh points
324 integer :: i, ip, ndim
325 integer :: process
326 real(real64) :: dummy
327
329
330 if (mesh%parallel_in_domains) then
331 ndim = size(values, 1)
332 do i = 1, size(values, 2)
333 ip = mesh_nearest_point(mesh, values(:, i), dummy, process)
334 ! ip is defined for `process`, only
335 if (mesh%mpi_grp%rank == process) values(:, i) = mesh%x(:, ip)
336 call mesh%mpi_grp%bcast(values(:, i), ndim, mpi_double_precision, process)
337 enddo
338 else
339 do i = 1, size(values, 2)
340 ip = mesh_nearest_point(mesh, values(:, i), dummy, process)
341 values(:, i) = mesh%x(:, ip)
342 enddo
343 endif
344
346
347 end subroutine mesh_discretize_values_to_mesh
348
349 ! --------------------------------------------------------------
353 ! --------------------------------------------------------------
354 real(real64) function mesh_gcutoff(mesh) result(gmax)
355 class(mesh_t), intent(in) :: mesh
356
357 push_sub(mesh_gcutoff)
358 gmax = m_pi / (maxval(mesh%spacing))
359
360 pop_sub(mesh_gcutoff)
361 end function mesh_gcutoff
362
363 ! --------------------------------------------------------------
364 subroutine mesh_write_fingerprint(mesh, dir, filename, mpi_grp, namespace, ierr)
365 type(mesh_t), intent(in) :: mesh
366 character(len=*), intent(in) :: dir
367 character(len=*), intent(in) :: filename
368 type(mpi_grp_t), intent(in) :: mpi_grp
369 type(namespace_t),intent(in) :: namespace
370 integer, intent(out) :: ierr
371
372 integer :: iunit, ii
373
374 push_sub(mesh_write_fingerprint)
375
376 ierr = 0
377
378 iunit = io_open(trim(dir)//"/"//trim(filename), namespace, action='write', &
379 die=.false., grp=mpi_grp)
380 if (iunit == -1) then
381 message(1) = "Unable to open file '"//trim(dir)//"/"//trim(filename)//"'."
382 call messages_warning(1, namespace=namespace)
383 ierr = ierr + 1
384 else
385 if (mpi_grp%is_root()) then
386 write(iunit, '(a20,i21)') 'np_part_global= ', mesh%np_part_global
387 write(iunit, '(a20,i21)') 'np_global= ', mesh%np_global
388 write(iunit, '(a20,i21)') 'algorithm= ', 1
389 write(iunit, '(a20,i21)') 'checksum= ', mesh%idx%checksum
390 write(iunit, '(a20,i21)') 'bits= ', mesh%idx%bits
391 write(iunit, '(a20,i21)') 'dim= ', mesh%idx%dim
392 write(iunit, '(a20,i21)') 'type= ', mesh%idx%type
393 do ii = 1, mesh%idx%dim
394 write(iunit, '(a7,i2,a11,i21)') 'offset(',ii,')= ', mesh%idx%offset(ii)
395 end do
396 do ii = 1, mesh%idx%dim
397 write(iunit, '(a7,i2,a11,i21)') 'nn(',ii,')= ', mesh%idx%nr(2, ii) - mesh%idx%nr(1, ii) + 1
398 end do
399 end if
400 call io_close(iunit, grp=mpi_grp)
401 end if
402
404 end subroutine mesh_write_fingerprint
405
406
407 ! -----------------------------------------------------------------------
412 subroutine mesh_read_fingerprint(mesh, dir, filename, mpi_grp, namespace, &
413 read_np_part, read_np, bits, type, offset, nn, ierr)
414 type(mesh_t), intent(in) :: mesh
415 character(len=*), intent(in) :: dir
416 character(len=*), intent(in) :: filename
417 type(mpi_grp_t), intent(in) :: mpi_grp
418 type(namespace_t),intent(in) :: namespace
419 integer(int64), intent(out) :: read_np_part
420 integer(int64), intent(out) :: read_np
421 integer, intent(out) :: bits
422 integer, intent(out) :: type
423 integer, intent(out) :: offset(1:mesh%idx%dim)
424 integer, intent(out) :: nn(1:mesh%idx%dim)
425 integer, intent(out) :: ierr
426
427 character(len=20) :: str
428 character(len=100) :: lines(7)
429 integer :: iunit, algorithm, dim, err, ii
430 integer(int64) :: checksum
431
432 push_sub(mesh_read_fingerprint)
433
434 ierr = 0
435
436 read_np_part = 0_int64
437 read_np = 0_int64
438
439 iunit = io_open(trim(dir)//"/"//trim(filename), namespace, action='read', &
440 status='old', die=.false., grp=mpi_grp)
441 if (iunit == -1) then
442 ierr = ierr + 1
443 message(1) = "Unable to open file '"//trim(dir)//"/"//trim(filename)//"'."
444 call messages_warning(1, namespace=namespace)
445 else
446 call iopar_read(mpi_grp, iunit, lines, 7, err)
447 if (err /= 0) then
448 ierr = ierr + 4
449 else
450 read(lines(1), '(a20,i21)') str, read_np_part
451 read(lines(2), '(a20,i21)') str, read_np
452 read(lines(3), '(a20,i21)') str, algorithm
453 read(lines(4), '(a20,i21)') str, checksum
454 read(lines(5), '(a20,i21)') str, bits
455 read(lines(6), '(a20,i21)') str, dim
456 read(lines(7), '(a20,i21)') str, type
457 ! only allow restarting simulations with the same dimensions
458 if (dim /= mesh%idx%dim) then
459 ierr = ierr + 8
460 else
461 ! read offset, has dim lines
462 call iopar_read(mpi_grp, iunit, lines, dim, err)
463 if (err /= 0) then
464 ierr = ierr + 4
465 else
466 do ii = 1, dim
467 read(lines(ii), '(a20,i21)') str, offset(ii)
468 end do
469 end if
470
471 ! read nn, has dim lines
472 call iopar_read(mpi_grp, iunit, lines, dim, err)
473 if (err /= 0) then
474 ierr = ierr + 4
475 else
476 do ii = 1, dim
477 read(lines(ii), '(a20,i21)') str, nn(ii)
478 end do
479 end if
480 end if
481
482 assert(read_np_part >= read_np)
483
484 if (read_np_part == mesh%np_part_global &
485 .and. read_np == mesh%np_global &
486 .and. algorithm == 1 &
487 .and. checksum == mesh%idx%checksum) then
488 read_np_part = 0
489 read_np = 0
490 end if
491 end if
492
493 call io_close(iunit, grp=mpi_grp)
494 end if
495
496 pop_sub(mesh_read_fingerprint)
497 end subroutine mesh_read_fingerprint
498
499 ! ---------------------------------------------------------
500 subroutine mesh_check_dump_compatibility(mesh, dir, filename, namespace, mpi_grp, grid_changed, grid_reordered, map, ierr)
501 type(mesh_t), intent(in) :: mesh
502 character(len=*), intent(in) :: dir
503 character(len=*), intent(in) :: filename
504 type(namespace_t), intent(in) :: namespace
505 type(mpi_grp_t), intent(in) :: mpi_grp
506 logical, intent(out) :: grid_changed
507 logical, intent(out) :: grid_reordered
508 integer(int64), allocatable, intent(out) :: map(:)
509 integer, intent(out) :: ierr
510
511 integer(int64) :: ipg, ipg_new, read_np_part, read_np
512 integer :: err, idir
513 integer :: bits, type, offset(mesh%idx%dim), point(mesh%idx%dim), nn(mesh%idx%dim)
514 integer(int64), allocatable :: read_indices(:)
515 type(index_t) :: idx_old
516
518
519 ierr = 0
520
521 grid_changed = .false.
522 grid_reordered = .false.
523
524 ! Read the mesh fingerprint
525 call mesh_read_fingerprint(mesh, dir, filename, mpi_grp, namespace, read_np_part, read_np, &
526 bits, type, offset, nn, err)
527 if (err /= 0) then
528 ierr = 1
529 message(1) = "Unable to read mesh fingerprint from '"//trim(dir)//"/"//trim(filename)//"'."
530 call messages_warning(1, namespace=namespace)
531
532 else if (read_np > 0) then
533 if (.not. associated(mesh%box)) then
534 ! We can only check the compatibility of two meshes that have different fingerprints if we also
535 ! have the simulation box. In the case we do not, we will assume that the fingerprint is enough.
536 ierr = ierr + 2
537 else
538 grid_changed = .true.
539
540 ! perhaps only the order of the points changed, this can only
541 ! happen if the number of points is the same and no points maps
542 ! to zero (this is checked below)
543 grid_reordered = (read_np == mesh%np_global)
544
545 ! the grid is different, so we read the coordinates.
546 safe_allocate(read_indices(1:read_np_part))
547 call io_binary_read(trim(io_workpath(dir, namespace))//"/indices.obf", read_np_part, &
548 read_indices, err)
549 if (err /= 0) then
550 ierr = ierr + 4
551 message(1) = "Unable to read index map from '"//trim(dir)//"'."
552 call messages_warning(1, namespace=namespace)
553 else
554 ! dummy index object
555 call index_init(idx_old, mesh%idx%dim)
556 idx_old%type = type
557 idx_old%bits = bits
558 idx_old%nr(1, :) = -offset
559 idx_old%nr(2, :) = -offset + nn - 1
560 idx_old%offset = offset
561 idx_old%stride(1) = 1
562 do idir = 2, mesh%idx%dim
563 idx_old%stride(idir) = idx_old%stride(idir-1) * nn(idir-1)
564 end do
565 ! generate the map
566 safe_allocate(map(1:read_np))
567 do ipg = 1, read_np
568 ! get nd-index from old 1d index
569 call index_spatial_to_point(idx_old, mesh%idx%dim, read_indices(ipg), point)
570 ! get new global index
571 ipg_new = mesh_global_index_from_coords(mesh, point)
572 map(ipg) = ipg_new
573 ! ignore boundary points
574 if (map(ipg) > mesh%np_global) map(ipg) = 0
575 ! if the map is zero for one point, it is not a simple reordering
576 if (map(ipg) == 0) grid_reordered = .false.
577 end do
578 call index_end(idx_old)
579 end if
580
581 safe_deallocate_a(read_indices)
582 end if
583 end if
584
586 end subroutine mesh_check_dump_compatibility
587
588
589 ! --------------------------------------------------------------
590 recursive subroutine mesh_end(this)
591 class(mesh_t), intent(inout) :: this
592
593 push_sub(mesh_end)
594
595#ifdef HAVE_MPI
596 call lmpi_destroy_shared_memory_window(this%idx%window_grid_to_spatial)
597 call lmpi_destroy_shared_memory_window(this%idx%window_spatial_to_grid)
598 nullify(this%idx%grid_to_spatial_global)
599 nullify(this%idx%spatial_to_grid_global)
600#else
601 safe_deallocate_p(this%idx%grid_to_spatial_global)
602 safe_deallocate_p(this%idx%spatial_to_grid_global)
603#endif
604
605 safe_deallocate_a(this%x)
606 safe_deallocate_a(this%chi)
607 safe_deallocate_a(this%vol_pp)
608 safe_deallocate_a(this%jacobian_inverse)
609
610 if (this%parallel_in_domains) then
611 call par_vec_end(this%pv)
612 call partition_end(this%partition)
613 end if
614
615 call index_end(this%idx)
616 safe_deallocate_a(this%spacing)
617
618 pop_sub(mesh_end)
619 end subroutine mesh_end
620
621
627 ! ---------------------------------------------------------
628 integer(int64) function mesh_periodic_point(mesh, space, ip) result(ipg)
629 class(mesh_t), intent(in) :: mesh
630 class(space_t),intent(in) :: space
631 integer, intent(in) :: ip
632
633 integer :: ix(space%dim), nr(2, space%dim), idim
634 real(real64) :: xx(space%dim), rr, ufn_re, ufn_im
635
636 ! no push_sub, called too frequently
637
638 call mesh_local_index_to_coords(mesh, ip, ix)
639 nr(1, :) = mesh%idx%nr(1, :) + mesh%idx%enlarge(:)
640 nr(2, :) = mesh%idx%nr(2, :) - mesh%idx%enlarge(:)
641
642 do idim = 1, space%periodic_dim
643 do while (ix(idim) < nr(1, idim))
644 ix(idim) = ix(idim) + mesh%idx%ll(idim)
645 end do
646 do while (ix(idim) > nr(2, idim))
647 ix(idim) = ix(idim) - mesh%idx%ll(idim)
648 end do
649 end do
650
651 ipg = mesh_global_index_from_coords(mesh, ix)
652 assert(ipg > 0)
653
654 if (mesh%masked_periodic_boundaries) then
655 call mesh_r(mesh, ip, rr, coords = xx)
656 call parse_expression(ufn_re, ufn_im, space%dim, xx, rr, m_zero, mesh%periodic_boundary_mask)
657 if (int(ufn_re) == 0) ipg = mesh_local2global(mesh, ip) ! Nothing will be done
658 end if
659
660 end function mesh_periodic_point
661
662 integer(int64) function mesh_periodic_point_global(mesh, space, ipg_in) result(ipg)
663 class(mesh_t), intent(in) :: mesh
664 class(space_t),intent(in) :: space
665 integer(int64),intent(in) :: ipg_in
666
667 integer :: ix(space%dim), nr(2, space%dim), idim
668 real(real64) :: xx(space%dim), rr, ufn_re, ufn_im
669
670 ! no push_sub, called too frequently
671
672 call mesh_global_index_to_coords(mesh, ipg_in, ix)
673 nr(1, :) = mesh%idx%nr(1, :) + mesh%idx%enlarge(:)
674 nr(2, :) = mesh%idx%nr(2, :) - mesh%idx%enlarge(:)
675
676 do idim = 1, space%periodic_dim
677 if (ix(idim) < nr(1, idim)) ix(idim) = ix(idim) + mesh%idx%ll(idim)
678 if (ix(idim) > nr(2, idim)) ix(idim) = ix(idim) - mesh%idx%ll(idim)
679 end do
680
681 ipg = mesh_global_index_from_coords(mesh, ix)
682 assert(ipg > 0)
683
684 if (mesh%masked_periodic_boundaries) then
685 call mesh_r_global(mesh, ipg_in, rr, coords = xx)
686 call parse_expression(ufn_re, ufn_im, space%dim, xx, rr, m_zero, mesh%periodic_boundary_mask)
687 if (int(ufn_re) == 0) ipg = ipg_in ! Nothing will be done
688 end if
689
690 end function mesh_periodic_point_global
691
692
693
694 ! ---------------------------------------------------------
695 real(real64) pure function mesh_global_memory(mesh) result(memory)
696 use iso_c_binding, only: c_sizeof, c_long_long
697 class(mesh_t), intent(in) :: mesh
698
699 ! 2 global index arrays
700 memory = c_sizeof(c_long_long) * real(mesh%np_part_global, real64) * 2
701
702 end function mesh_global_memory
703
704
705 ! ---------------------------------------------------------
706 real(real64) pure function mesh_local_memory(mesh) result(memory)
707 use iso_c_binding, only: c_sizeof, c_long_long
708 class(mesh_t), intent(in) :: mesh
709
710 memory = m_zero
711
712 ! x
713 memory = memory + real64 * real(mesh%np_part, real64) * mesh%idx%dim
714 ! local index arrays
715 memory = memory + c_sizeof(c_long_long) * real(mesh%np_part, real64) * 2
716 end function mesh_local_memory
717
718
721 function mesh_x_global(mesh, ipg) result(xx)
722 class(mesh_t), intent(in) :: mesh
723 integer(int64), intent(in) :: ipg
724 real(real64) :: xx(1:mesh%box%dim)
725
726 real(real64) :: chi(1:mesh%box%dim)
727 integer :: ix(1:mesh%box%dim)
728
729 ! no push_sub because function is called too frequently
730
731 call mesh_global_index_to_coords(mesh, ipg, ix)
732 chi = ix * mesh%spacing
733 xx = mesh%coord_system%to_cartesian(chi)
734
735 end function mesh_x_global
736
737
738 ! ---------------------------------------------------------
739 subroutine mesh_check_symmetries(mesh, symm, periodic_dim)
740 class(mesh_t), intent(in) :: mesh
741 type(symmetries_t), intent(in) :: symm
742 integer, intent(in) :: periodic_dim
743
744 integer :: iop, ip, idim, nops, ix(1:3)
745 real(real64) :: destpoint(1:3), srcpoint(1:3), lsize(1:3), offset(1:3)
746 real(real64), parameter :: tol_spacing = 1e-12_real64
747
748 !If all the axis have the same spacing and the same length
749 !the grid is by obviously symmetric
750 !Indeed, reduced coordinates are proportional to the point index
751 !and the reduced rotation are integer matrices
752 !The result of the product is also proportional to an integer
753 !and therefore belong to the grid.
754 if (mesh%idx%ll(1) == mesh%idx%ll(2) .and. &
755 mesh%idx%ll(2) == mesh%idx%ll(3) .and. &
756 abs(mesh%spacing(1) - mesh%spacing(2)) < tol_spacing .and. &
757 abs(mesh%spacing(2) - mesh%spacing(3)) < tol_spacing) return
758
759 push_sub(mesh_check_symmetries)
760
761 message(1) = "Checking if the real-space grid is symmetric"
762 call messages_info(1)
763
764 lsize(1:3) = real(mesh%idx%ll(1:3), real64)
765 offset(1:3) = real(mesh%idx%nr(1, 1:3) + mesh%idx%enlarge(1:3), real64)
766
767 nops = symmetries_number(symm)
768
769 do ip = 1, mesh%np
770 !We use floating point coordinates to check if the symmetric point
771 !belong to the grid.
772 !If yes, it should have integer reduced coordinates
773 call mesh_local_index_to_coords(mesh, ip, ix)
774 destpoint(1:3) = real(ix(1:3), real64) - offset(1:3)
775 ! offset moves corner of cell to origin, in integer mesh coordinates
776 assert(all(destpoint >= 0))
777 assert(all(destpoint < lsize))
778
779 ! move to center of cell in real coordinates
780 destpoint = destpoint - real(int(lsize)/2, real64)
781
782 !convert to proper reduced coordinates
783 do idim = 1, 3
784 destpoint(idim) = destpoint(idim)/lsize(idim)
785 end do
786
787 ! iterate over all points that go to this point by a symmetry operation
788 do iop = 1, nops
789 srcpoint = symm_op_apply_red(symm%ops(iop), destpoint)
791 !We now come back to what should be an integer, if the symmetric point beloings to the grid
792 do idim = 1, 3
793 srcpoint(idim) = srcpoint(idim)*lsize(idim)
794 end do
795
796 ! move back to reference to origin at corner of cell
797 srcpoint = srcpoint + real(int(lsize)/2, real64)
798
799 ! apply periodic boundary conditions in periodic directions
800 do idim = 1, periodic_dim
801 if (nint(srcpoint(idim)) < 0 .or. nint(srcpoint(idim)) >= lsize(idim)) then
802 srcpoint(idim) = modulo(srcpoint(idim)+m_half*symprec, lsize(idim))
803 end if
804 end do
805 assert(all(srcpoint >= -symprec))
806 assert(all(srcpoint < lsize))
807
808 srcpoint(1:3) = srcpoint(1:3) + offset(1:3)
809
810 if (any(srcpoint-anint(srcpoint)> symprec*m_two)) then
811 message(1) = "The real-space grid breaks at least one of the symmetries of the system."
812 message(2) = "Change your spacing or use SymmetrizeDensity=no."
813 call messages_fatal(2)
814 end if
815 end do
816 end do
817
818 pop_sub(mesh_check_symmetries)
819 end subroutine
820
823 integer(int64) function mesh_global_index_from_coords(mesh, ix) result(index)
824 class(mesh_t), intent(in) :: mesh
825 integer, intent(in) :: ix(:)
826
827 index = index_from_coords(mesh%idx, ix)
829
832 subroutine mesh_global_index_to_coords(mesh, ipg, ix)
833 type(mesh_t), intent(in) :: mesh
834 integer(int64), intent(in) :: ipg
835 integer, intent(out) :: ix(:)
836
837 call index_to_coords(mesh%idx, ipg, ix)
838 end subroutine mesh_global_index_to_coords
839
842 integer function mesh_local_index_from_coords(mesh, ix) result(ip)
843 type(mesh_t), intent(in) :: mesh
844 integer, intent(in) :: ix(:)
845
846 integer(int64) :: ipg
847
848 ipg = index_from_coords(mesh%idx, ix)
849 ip = mesh_global2local(mesh, ipg)
851
854 subroutine mesh_local_index_to_coords(mesh, ip, ix)
855 type(mesh_t), intent(in) :: mesh
856 integer, intent(in) :: ip
857 integer, intent(out) :: ix(:)
858
859 integer(int64) :: ipg
860
861 ipg = mesh_local2global(mesh, ip)
862 call index_to_coords(mesh%idx, ipg, ix)
863 end subroutine mesh_local_index_to_coords
864
866 integer(int64) function mesh_local2global(mesh, ip) result(ipg)
867 type(mesh_t), intent(in) :: mesh
868 integer, intent(in) :: ip
869
870 ipg = par_vec_local2global(mesh%pv, ip)
871 end function mesh_local2global
872
876 integer function mesh_global2local(mesh, ipg) result(ip)
877 type(mesh_t), intent(in) :: mesh
878 integer(int64), intent(in) :: ipg
879
880 ip = par_vec_global2local(mesh%pv, ipg)
881 end function mesh_global2local
882
883 !-----------------------------------------------------------------------------
885 subroutine mesh_minmaxloc(mesh, min_or_max, rank_min_or_max, op)
886 type(mesh_t), intent(in) :: mesh
887 real(real64), intent(inout) :: min_or_max
888 integer, intent(out) :: rank_min_or_max
889 type(mpi_op), intent(in) :: op
890
891 real(real64) :: loc_in(2), loc_out(2)
892
893 push_sub(mesh_minmaxloc)
894
895 assert(op == mpi_minloc .or. op == mpi_maxloc)
896
897 rank_min_or_max = 0
898 if (mesh%parallel_in_domains) then
899 loc_in(1) = min_or_max
900 loc_in(2) = mesh%mpi_grp%rank
901 call mesh%mpi_grp%allreduce(loc_in, loc_out, 1, mpi_2double_precision, op)
902 min_or_max = loc_out(1)
903 rank_min_or_max = nint(loc_out(2))
904 end if
905
906 pop_sub(mesh_minmaxloc)
907 end subroutine mesh_minmaxloc
908
909
910#include "undef.F90"
911#include "real.F90"
912#include "mesh_inc.F90"
913
914#include "undef.F90"
915#include "complex.F90"
916#include "mesh_inc.F90"
917
918#include "undef.F90"
919#include "integer.F90"
920#include "mesh_inc.F90"
921
922#include "undef.F90"
923#include "integer8.F90"
924#include "mesh_inc.F90"
925
926end module mesh_oct_m
928
929!! Local Variables:
930!! mode: f90
931!! coding: utf-8
932!! End:
real(real64), parameter, public m_two
Definition: global.F90:193
real(real64), parameter, public m_zero
Definition: global.F90:191
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:189
This module implements a simple hash table for non-negative integer keys and integer values.
Definition: iihash.F90:127
This module implements the index, used for the mesh points.
Definition: index.F90:124
Definition: io.F90:116
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
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:724
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:928
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
subroutine, public mesh_double_box(space, mesh, alpha, db)
finds the dimension of a box doubled in the non-periodic dimensions
Definition: mesh.F90:285
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:919
subroutine, public mesh_check_symmetries(mesh, symm, periodic_dim)
Definition: mesh.F90:835
subroutine, public mesh_write_info(this, iunit, namespace)
Definition: mesh.F90:310
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:950
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:416
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:386
integer(int64) function, public mesh_periodic_point_global(mesh, space, ipg_in)
Definition: mesh.F90:758
real(real64) pure function, public mesh_global_memory(mesh)
Definition: mesh.F90:791
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:509
pure subroutine, public mesh_r(mesh, ip, rr, origin, coords)
return the distance to the origin for a given grid point
Definition: mesh.F90:342
subroutine, public mesh_check_dump_compatibility(mesh, dir, filename, namespace, mpi_grp, grid_changed, grid_reordered, map, ierr)
Definition: mesh.F90:596
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:450
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:981
integer function, public mesh_global2local(mesh, ipg)
This function returns the local mesh index for a given global index.
Definition: mesh.F90:972
recursive subroutine, public mesh_end(this)
Definition: mesh.F90:686
real(real64) pure function, public mesh_local_memory(mesh)
Definition: mesh.F90:802
integer(int64) function, public mesh_local2global(mesh, ip)
This function returns the global mesh index for a given local index.
Definition: mesh.F90:962
real(real64) function, dimension(1:mesh%box%dim), public mesh_x_global(mesh, ipg)
Given a global point index, this function returns the coordinates of the point.
Definition: mesh.F90:817
subroutine mesh_r_global(mesh, ipg, rr, origin, coords)
return the distance to the origin for a given grid point
Definition: mesh.F90:362
subroutine, public mesh_write_fingerprint(mesh, dir, filename, mpi_grp, namespace, ierr)
Definition: mesh.F90:460
subroutine mesh_init(this)
Definition: mesh.F90:273
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:594
This module contains some common usage patterns of MPI routines.
Definition: mpi_lib.F90:117
Some general things and nomenclature:
Definition: par_vec.F90:173
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:134
character(len=20) pure function, public units_abbrev(this)
Definition: unit.F90:225
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:261
define a grid on a plane.
Definition: mesh.F90:249
Describes mesh distribution to nodes.
Definition: mesh.F90:187
int true(void)
real(real64) function values(xx)
Definition: test.F90:2031