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