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