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, &
81
93 type, extends(basis_set_abst_t) :: mesh_t
94 ! Components are public by default
95 class(box_t), pointer :: box
96 class(coordinate_system_t), pointer :: coord_system
97 type(index_t) :: idx
98 logical :: use_curvilinear
99
100 real(real64), allocatable :: spacing(:)
101
102 ! When running serially, the local number of points is
103 ! equal to the global number of points.
104 ! Otherwise, the next two are different on each node.
105 integer :: np
106 integer :: np_part
107 integer(int64) :: np_global
108 integer(int64) :: np_part_global
109 logical :: parallel_in_domains
110 type(mpi_grp_t) :: mpi_grp
111 type(par_vec_t) :: pv
112 type(partition_t) :: partition
113
114 real(real64), allocatable :: x(:,:)
115 real(real64) :: volume_element
116 real(real64), allocatable :: vol_pp(:)
117
118 logical :: masked_periodic_boundaries
119 character(len=256) :: periodic_boundary_mask
121 contains
122 procedure :: end => mesh_end
123 procedure :: init => mesh_init
124 procedure :: write_info => mesh_write_info
125 procedure :: dmesh_allreduce_0, zmesh_allreduce_0, imesh_allreduce_0, lmesh_allreduce_0
126 procedure :: dmesh_allreduce_1, zmesh_allreduce_1, imesh_allreduce_1, lmesh_allreduce_1
127 procedure :: dmesh_allreduce_2, zmesh_allreduce_2, imesh_allreduce_2, lmesh_allreduce_2
128 procedure :: dmesh_allreduce_3, zmesh_allreduce_3, imesh_allreduce_3, lmesh_allreduce_3
129 procedure :: dmesh_allreduce_4, zmesh_allreduce_4, imesh_allreduce_4, lmesh_allreduce_4
130 procedure :: dmesh_allreduce_5, zmesh_allreduce_5, imesh_allreduce_5, lmesh_allreduce_5
131 generic :: allreduce => dmesh_allreduce_0, zmesh_allreduce_0, imesh_allreduce_0, lmesh_allreduce_0
132 generic :: allreduce => dmesh_allreduce_1, zmesh_allreduce_1, imesh_allreduce_1, lmesh_allreduce_1
133 generic :: allreduce => dmesh_allreduce_2, zmesh_allreduce_2, imesh_allreduce_2, lmesh_allreduce_2
134 generic :: allreduce => dmesh_allreduce_3, zmesh_allreduce_3, imesh_allreduce_3, lmesh_allreduce_3
135 generic :: allreduce => dmesh_allreduce_4, zmesh_allreduce_4, imesh_allreduce_4, lmesh_allreduce_4
136 generic :: allreduce => dmesh_allreduce_5, zmesh_allreduce_5, imesh_allreduce_5, lmesh_allreduce_5
137 end type mesh_t
138
150 type mesh_plane_t
151 ! Components are public by default
152 real(real64) :: n(3)
153 real(real64) :: u(3), v(3)
154 real(real64) :: origin(3)
155 real(real64) :: spacing
156 integer :: nu, mu, nv, mv
157 end type mesh_plane_t
158
162 type mesh_line_t
163 ! Components are public by default
164 real(real64) :: n(2)
165 real(real64) :: u(2)
166 real(real64) :: origin(2)
167 real(real64) :: spacing
168 integer :: nu, mu
169 end type mesh_line_t
170
171contains
172
173 subroutine mesh_init(this)
174 class(mesh_t), intent(inout) :: this
175
176 push_sub(mesh_init)
177
178 call this%set_time_dependent(.false.)
179
180 pop_sub(mesh_init)
181 end subroutine mesh_init
182
183! ---------------------------------------------------------
185 subroutine mesh_double_box(space, mesh, alpha, db)
186 class(space_t), intent(in) :: space
187 type(mesh_t), intent(in) :: mesh
188 real(real64), intent(in) :: alpha
189 integer, intent(out) :: db(:)
191 integer :: idir
192
194
195 db = 1
196
197 ! double mesh with 2n points
198 do idir = 1, space%periodic_dim
199 db(idir) = mesh%idx%ll(idir)
200 end do
201 do idir = space%periodic_dim + 1, space%dim
202 db(idir) = nint(alpha * (mesh%idx%ll(idir) - 1)) + 1
203 end do
206 end subroutine mesh_double_box
209 ! ---------------------------------------------------------
210 subroutine mesh_write_info(this, iunit, namespace)
211 class(mesh_t), intent(in) :: this
212 integer, optional, intent(in) :: iunit
213 type(namespace_t), optional, intent(in) :: namespace
214
215 integer :: ii
216 real(real64) :: cutoff
220 write(message(1),'(3a)') ' Spacing [', trim(units_abbrev(units_out%length)), '] = ('
221 do ii = 1, this%box%dim
222 if (ii > 1) write(message(1), '(2a)') trim(message(1)), ','
223 write(message(1), '(a,f6.3)') trim(message(1)), units_from_atomic(units_out%length, this%spacing(ii))
224 end do
225 write(message(1), '(5a,f12.5)') trim(message(1)), ') ', &
226 ' volume/point [', trim(units_abbrev(units_out%length**this%box%dim)), '] = ', &
227 units_from_atomic(units_out%length**this%box%dim, this%vol_pp(1))
229 write(message(2),'(a, i10)') ' # inner mesh = ', this%np_global
230 write(message(3),'(a, i10)') ' # total mesh = ', this%np_part_global
231
232 cutoff = mesh_gcutoff(this)**2 / m_two
233 write(message(4),'(3a,f12.6,a,f12.6)') ' Grid Cutoff [', trim(units_abbrev(units_out%energy)),'] = ', &
234 units_from_atomic(units_out%energy, cutoff), ' Grid Cutoff [Ry] = ', cutoff * m_two
235 call messages_info(4, iunit=iunit, namespace=namespace)
236
237 pop_sub(mesh_write_info)
238 end subroutine mesh_write_info
239
240
242 pure subroutine mesh_r(mesh, ip, rr, origin, coords)
243 class(mesh_t), intent(in) :: mesh
244 integer, intent(in) :: ip
245 real(real64), intent(out) :: rr
246 real(real64), intent(in), optional :: origin(:)
247 real(real64), intent(out), optional :: coords(:)
249 real(real64) :: xx(1:mesh%box%dim)
250
251 xx = mesh%x(ip, :)
252 if (present(origin)) xx = xx - origin
253 rr = norm2(xx)
254
255 if (present(coords)) then
256 coords = xx
257 end if
259 end subroutine mesh_r
262 subroutine mesh_r_global(mesh, ipg, rr, origin, coords)
263 class(mesh_t), intent(in) :: mesh
264 integer(int64),intent(in) :: ipg
265 real(real64), intent(out) :: rr
266 real(real64), intent(in), optional :: origin(:)
267 real(real64), intent(out), optional :: coords(:)
268
269 real(real64) :: xx(1:mesh%box%dim)
270
271 xx = mesh_x_global(mesh, ipg)
272 if (present(origin)) xx = xx - origin
273 rr = norm2(xx)
274
275 if (present(coords)) then
276 coords = xx
277 end if
279 end subroutine mesh_r_global
280
281
286 integer function mesh_nearest_point(mesh, pos, dmin, rankmin) result(ind)
287 class(mesh_t),intent(in) :: mesh
288 real(real64), intent(in) :: pos(:)
289 real(real64), intent(out) :: dmin
290 integer, intent(out) :: rankmin
291 ! This returns 0 in the absence of domain decomposition
292
293 real(real64) :: dd
294 integer :: imin, ip
295
296 push_sub(mesh_nearest_point)
297
298 ! find the point of the grid that is closer to the atom
299 dmin = m_zero
300 do ip = 1, mesh%np
301 dd = sum((pos - mesh%x(ip, :))**2)
302 if ((dd < dmin) .or. (ip == 1)) then
303 imin = ip
304 dmin = dd
305 end if
306 end do
307
308 call mesh_minmaxloc(mesh, dmin, rankmin, mpi_minloc)
309 call mesh%mpi_grp%bcast(imin, 1, mpi_integer, rankmin)
310
311 ind = imin
312 pop_sub(mesh_nearest_point)
313 end function mesh_nearest_point
314
316 subroutine mesh_discretize_values_to_mesh(mesh, values)
317 class(mesh_t), intent(in ) :: mesh
318 real(real64), intent(inout) :: values(:, :)
319 ! Out: Values discretized to mesh points
320 integer :: i, ip, ndim
321 integer :: process
322 real(real64) :: dummy
323
325
326 if (mesh%parallel_in_domains) then
327 ndim = size(values, 1)
328 do i = 1, size(values, 2)
329 ip = mesh_nearest_point(mesh, values(:, i), dummy, process)
330 ! ip is defined for `process`, only
331 if (mesh%mpi_grp%rank == process) values(:, i) = mesh%x(ip, :)
332 call mesh%mpi_grp%bcast(values(:, i), ndim, mpi_double_precision, process)
333 enddo
334 else
335 do i = 1, size(values, 2)
336 ip = mesh_nearest_point(mesh, values(:, i), dummy, process)
337 values(:, i) = mesh%x(ip, :)
338 enddo
339 endif
340
342
343 end subroutine mesh_discretize_values_to_mesh
344
345 ! --------------------------------------------------------------
349 ! --------------------------------------------------------------
350 real(real64) function mesh_gcutoff(mesh) result(gmax)
351 class(mesh_t), intent(in) :: mesh
352
353 push_sub(mesh_gcutoff)
354 gmax = m_pi / (maxval(mesh%spacing))
356 pop_sub(mesh_gcutoff)
357 end function mesh_gcutoff
358
359 ! --------------------------------------------------------------
360 subroutine mesh_write_fingerprint(mesh, dir, filename, mpi_grp, namespace, ierr)
361 type(mesh_t), intent(in) :: mesh
362 character(len=*), intent(in) :: dir
363 character(len=*), intent(in) :: filename
364 type(mpi_grp_t), intent(in) :: mpi_grp
365 type(namespace_t),intent(in) :: namespace
366 integer, intent(out) :: ierr
367
368 integer :: iunit, ii
369
370 push_sub(mesh_write_fingerprint)
371
372 ierr = 0
373
374 iunit = io_open(trim(dir)//"/"//trim(filename), namespace, action='write', &
375 die=.false., grp=mpi_grp)
376 if (iunit <= 0) then
377 message(1) = "Unable to open file '"//trim(dir)//"/"//trim(filename)//"'."
378 call messages_warning(1, namespace=namespace)
379 ierr = ierr + 1
380 else
381 if (mpi_grp_is_root(mpi_grp)) then
382 write(iunit, '(a20,i21)') 'np_part_global= ', mesh%np_part_global
383 write(iunit, '(a20,i21)') 'np_global= ', mesh%np_global
384 write(iunit, '(a20,i21)') 'algorithm= ', 1
385 write(iunit, '(a20,i21)') 'checksum= ', mesh%idx%checksum
386 write(iunit, '(a20,i21)') 'bits= ', mesh%idx%bits
387 write(iunit, '(a20,i21)') 'dim= ', mesh%idx%dim
388 write(iunit, '(a20,i21)') 'type= ', mesh%idx%type
389 do ii = 1, mesh%idx%dim
390 write(iunit, '(a7,i2,a11,i21)') 'offset(',ii,')= ', mesh%idx%offset(ii)
391 end do
392 do ii = 1, mesh%idx%dim
393 write(iunit, '(a7,i2,a11,i21)') 'nn(',ii,')= ', mesh%idx%nr(2, ii) - mesh%idx%nr(1, ii) + 1
394 end do
395 end if
396 call io_close(iunit, grp=mpi_grp)
397 end if
398
400 end subroutine mesh_write_fingerprint
401
402
403 ! -----------------------------------------------------------------------
408 subroutine mesh_read_fingerprint(mesh, dir, filename, mpi_grp, namespace, &
409 read_np_part, read_np, bits, type, offset, nn, ierr)
410 type(mesh_t), intent(in) :: mesh
411 character(len=*), intent(in) :: dir
412 character(len=*), intent(in) :: filename
413 type(mpi_grp_t), intent(in) :: mpi_grp
414 type(namespace_t),intent(in) :: namespace
415 integer(int64), intent(out) :: read_np_part
416 integer(int64), intent(out) :: read_np
417 integer, intent(out) :: bits
418 integer, intent(out) :: type
419 integer, intent(out) :: offset(1:mesh%idx%dim)
420 integer, intent(out) :: nn(1:mesh%idx%dim)
421 integer, intent(out) :: ierr
422
423 character(len=20) :: str
424 character(len=100) :: lines(7)
425 integer :: iunit, algorithm, dim, err, ii
426 integer(int64) :: checksum
427
428 push_sub(mesh_read_fingerprint)
429
430 ierr = 0
431
432 read_np_part = 0_int64
433 read_np = 0_int64
434
435 iunit = io_open(trim(dir)//"/"//trim(filename), namespace, action='read', &
436 status='old', die=.false., grp=mpi_grp)
437 if (iunit <= 0) then
438 ierr = ierr + 1
439 message(1) = "Unable to open file '"//trim(dir)//"/"//trim(filename)//"'."
440 call messages_warning(1, namespace=namespace)
441 else
442 call iopar_read(mpi_grp, iunit, lines, 7, err)
443 if (err /= 0) then
444 ierr = ierr + 4
445 else
446 read(lines(1), '(a20,i21)') str, read_np_part
447 read(lines(2), '(a20,i21)') str, read_np
448 read(lines(3), '(a20,i21)') str, algorithm
449 read(lines(4), '(a20,i21)') str, checksum
450 read(lines(5), '(a20,i21)') str, bits
451 read(lines(6), '(a20,i21)') str, dim
452 read(lines(7), '(a20,i21)') str, type
453 ! only allow restarting simulations with the same dimensions
454 if (dim /= mesh%idx%dim) then
455 ierr = ierr + 8
456 else
457 ! read offset, has dim lines
458 call iopar_read(mpi_grp, iunit, lines, dim, err)
459 if (err /= 0) then
460 ierr = ierr + 4
461 else
462 do ii = 1, dim
463 read(lines(ii), '(a20,i21)') str, offset(ii)
464 end do
465 end if
466
467 ! read nn, has dim lines
468 call iopar_read(mpi_grp, iunit, lines, dim, err)
469 if (err /= 0) then
470 ierr = ierr + 4
471 else
472 do ii = 1, dim
473 read(lines(ii), '(a20,i21)') str, nn(ii)
474 end do
475 end if
476 end if
477
478 assert(read_np_part >= read_np)
479
480 if (read_np_part == mesh%np_part_global &
481 .and. read_np == mesh%np_global &
482 .and. algorithm == 1 &
483 .and. checksum == mesh%idx%checksum) then
484 read_np_part = 0
485 read_np = 0
486 end if
487 end if
488
489 call io_close(iunit, grp=mpi_grp)
490 end if
491
492 pop_sub(mesh_read_fingerprint)
493 end subroutine mesh_read_fingerprint
494
495 ! ---------------------------------------------------------
496 subroutine mesh_check_dump_compatibility(mesh, dir, filename, namespace, mpi_grp, grid_changed, grid_reordered, map, ierr)
497 type(mesh_t), intent(in) :: mesh
498 character(len=*), intent(in) :: dir
499 character(len=*), intent(in) :: filename
500 type(namespace_t), intent(in) :: namespace
501 type(mpi_grp_t), intent(in) :: mpi_grp
502 logical, intent(out) :: grid_changed
503 logical, intent(out) :: grid_reordered
504 integer(int64), allocatable, intent(out) :: map(:)
505 integer, intent(out) :: ierr
506
507 integer(int64) :: ipg, ipg_new, read_np_part, read_np
508 integer :: err, idir
509 integer :: bits, type, offset(mesh%idx%dim), point(mesh%idx%dim), nn(mesh%idx%dim)
510 integer(int64), allocatable :: read_indices(:)
511 type(index_t) :: idx_old
512
514
515 ierr = 0
516
517 grid_changed = .false.
518 grid_reordered = .false.
519
520 ! Read the mesh fingerprint
521 call mesh_read_fingerprint(mesh, dir, filename, mpi_grp, namespace, read_np_part, read_np, &
522 bits, type, offset, nn, err)
523 if (err /= 0) then
524 ierr = ierr + 1
525 message(1) = "Unable to read mesh fingerprint from '"//trim(dir)//"/"//trim(filename)//"'."
526 call messages_warning(1, namespace=namespace)
527
528 else if (read_np > 0) then
529 if (.not. associated(mesh%box)) then
530 ! We can only check the compatibility of two meshes that have different fingerprints if we also
531 ! have the simulation box. In the case we do not, we will assume that the fingerprint is enough.
532 ierr = ierr + 2
533 else
534 grid_changed = .true.
535
536 ! perhaps only the order of the points changed, this can only
537 ! happen if the number of points is the same and no points maps
538 ! to zero (this is checked below)
539 grid_reordered = (read_np == mesh%np_global)
540
541 ! the grid is different, so we read the coordinates.
542 safe_allocate(read_indices(1:read_np_part))
543 call io_binary_read(trim(io_workpath(dir, namespace))//"/indices.obf", read_np_part, &
544 read_indices, err)
545 if (err /= 0) then
546 ierr = ierr + 4
547 message(1) = "Unable to read index map from '"//trim(dir)//"'."
548 call messages_warning(1, namespace=namespace)
549 else
550 ! dummy index object
551 call index_init(idx_old, mesh%idx%dim)
552 idx_old%type = type
553 idx_old%bits = bits
554 idx_old%nr(1, :) = -offset
555 idx_old%nr(2, :) = -offset + nn - 1
556 idx_old%offset = offset
557 idx_old%stride(1) = 1
558 do idir = 2, mesh%idx%dim
559 idx_old%stride(idir) = idx_old%stride(idir-1) * nn(idir-1)
560 end do
561 ! generate the map
562 safe_allocate(map(1:read_np))
563 do ipg = 1, read_np
564 ! get nd-index from old 1d index
565 call index_spatial_to_point(idx_old, mesh%idx%dim, read_indices(ipg), point)
566 ! get new global index
567 ipg_new = mesh_global_index_from_coords(mesh, point)
568 map(ipg) = ipg_new
569 ! ignore boundary points
570 if (map(ipg) > mesh%np_global) map(ipg) = 0
571 ! if the map is zero for one point, it is not a simple reordering
572 if (map(ipg) == 0) grid_reordered = .false.
573 end do
574 call index_end(idx_old)
575 end if
576
577 safe_deallocate_a(read_indices)
578 end if
579 end if
580
582 end subroutine mesh_check_dump_compatibility
583
584
585 ! --------------------------------------------------------------
586 recursive subroutine mesh_end(this)
587 class(mesh_t), intent(inout) :: this
588
589 push_sub(mesh_end)
590
591#ifdef HAVE_MPI
592 call lmpi_destroy_shared_memory_window(this%idx%window_grid_to_spatial)
593 call lmpi_destroy_shared_memory_window(this%idx%window_spatial_to_grid)
594 nullify(this%idx%grid_to_spatial_global)
595 nullify(this%idx%spatial_to_grid_global)
596#else
597 safe_deallocate_p(this%idx%grid_to_spatial_global)
598 safe_deallocate_p(this%idx%spatial_to_grid_global)
599#endif
600
601 safe_deallocate_a(this%x)
602 safe_deallocate_a(this%vol_pp)
603
604 if (this%parallel_in_domains) then
605 call par_vec_end(this%pv)
606 call partition_end(this%partition)
607 end if
608
609 call index_end(this%idx)
610 safe_deallocate_a(this%spacing)
611
612 pop_sub(mesh_end)
613 end subroutine mesh_end
614
615
621 ! ---------------------------------------------------------
622 integer(int64) function mesh_periodic_point(mesh, space, ip) result(ipg)
623 class(mesh_t), intent(in) :: mesh
624 class(space_t),intent(in) :: space
625 integer, intent(in) :: ip
626
627 integer :: ix(space%dim), nr(2, space%dim), idim
628 real(real64) :: xx(space%dim), rr, ufn_re, ufn_im
629
630 ! no push_sub, called too frequently
631
632 call mesh_local_index_to_coords(mesh, ip, ix)
633 nr(1, :) = mesh%idx%nr(1, :) + mesh%idx%enlarge(:)
634 nr(2, :) = mesh%idx%nr(2, :) - mesh%idx%enlarge(:)
635
636 do idim = 1, space%periodic_dim
637 if (ix(idim) < nr(1, idim)) ix(idim) = ix(idim) + mesh%idx%ll(idim)
638 if (ix(idim) > nr(2, idim)) ix(idim) = ix(idim) - mesh%idx%ll(idim)
639 end do
640
641 ipg = mesh_global_index_from_coords(mesh, ix)
642 assert(ipg > 0)
643
644 if (mesh%masked_periodic_boundaries) then
645 call mesh_r(mesh, ip, rr, coords = xx)
646 call parse_expression(ufn_re, ufn_im, space%dim, xx, rr, m_zero, mesh%periodic_boundary_mask)
647 if (int(ufn_re) == 0) ipg = mesh_local2global(mesh, ip) ! Nothing will be done
648 end if
649
650 end function mesh_periodic_point
651
652 integer(int64) function mesh_periodic_point_global(mesh, space, ipg_in) result(ipg)
653 class(mesh_t), intent(in) :: mesh
654 class(space_t),intent(in) :: space
655 integer(int64),intent(in) :: ipg_in
656
657 integer :: ix(space%dim), nr(2, space%dim), idim
658 real(real64) :: xx(space%dim), rr, ufn_re, ufn_im
659
660 ! no push_sub, called too frequently
661
662 call mesh_global_index_to_coords(mesh, ipg_in, ix)
663 nr(1, :) = mesh%idx%nr(1, :) + mesh%idx%enlarge(:)
664 nr(2, :) = mesh%idx%nr(2, :) - mesh%idx%enlarge(:)
665
666 do idim = 1, space%periodic_dim
667 if (ix(idim) < nr(1, idim)) ix(idim) = ix(idim) + mesh%idx%ll(idim)
668 if (ix(idim) > nr(2, idim)) ix(idim) = ix(idim) - mesh%idx%ll(idim)
669 end do
670
671 ipg = mesh_global_index_from_coords(mesh, ix)
672 assert(ipg > 0)
673
674 if (mesh%masked_periodic_boundaries) then
675 call mesh_r_global(mesh, ipg_in, rr, coords = xx)
676 call parse_expression(ufn_re, ufn_im, space%dim, xx, rr, m_zero, mesh%periodic_boundary_mask)
677 if (int(ufn_re) == 0) ipg = ipg_in ! Nothing will be done
678 end if
680 end function mesh_periodic_point_global
681
682
683
684 ! ---------------------------------------------------------
685 real(real64) pure function mesh_global_memory(mesh) result(memory)
686 use iso_c_binding, only: c_sizeof, c_long_long
687 class(mesh_t), intent(in) :: mesh
688
689 ! 2 global index arrays
690 memory = c_sizeof(c_long_long) * real(mesh%np_part_global, real64) * 2
691
692 end function mesh_global_memory
693
694
695 ! ---------------------------------------------------------
696 real(real64) pure function mesh_local_memory(mesh) result(memory)
697 use iso_c_binding, only: c_sizeof, c_long_long
698 class(mesh_t), intent(in) :: mesh
699
700 memory = m_zero
701
702 ! x
703 memory = memory + real64 * real(mesh%np_part, real64) * mesh%idx%dim
704 ! local index arrays
705 memory = memory + c_sizeof(c_long_long) * real(mesh%np_part, real64) * 2
706 end function mesh_local_memory
707
708
709 ! ---------------------------------------------------------
710 function mesh_x_global(mesh, ipg) result(xx)
711 class(mesh_t), intent(in) :: mesh
712 integer(int64), intent(in) :: ipg
713 real(real64) :: xx(1:mesh%box%dim)
714
715 real(real64) :: chi(1:mesh%box%dim)
716 integer :: ix(1:mesh%box%dim)
717
718! no push_sub because function is called too frequently
719
720 call mesh_global_index_to_coords(mesh, ipg, ix)
721 chi = ix * mesh%spacing
722 xx = mesh%coord_system%to_cartesian(chi)
723
724 end function mesh_x_global
725
726
727 ! ---------------------------------------------------------
728 logical pure function mesh_compact_boundaries(mesh) result(cb)
729 type(mesh_t), intent(in) :: mesh
730
731 cb = .not. mesh%use_curvilinear .and. &
732 .not. mesh%parallel_in_domains
733
734 end function mesh_compact_boundaries
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
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:189
real(real64), parameter, public m_zero
Definition: global.F90:187
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:185
This module implements a simple hash table for non-negative integer keys and integer values.
Definition: iihash.F90:125
This module implements the index, used for the mesh points.
Definition: index.F90:122
Definition: io.F90:114
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
integer(int64) function, public mesh_periodic_point(mesh, space, ip)
This function returns the point inside the grid corresponding to a boundary point when PBCs are used....
Definition: mesh.F90:716
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:279
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:304
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:410
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:380
integer(int64) function, public mesh_periodic_point_global(mesh, space, ipg_in)
Definition: mesh.F90:746
real(real64) pure function, public mesh_global_memory(mesh)
Definition: mesh.F90:779
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:503
pure subroutine, public mesh_r(mesh, ip, rr, origin, coords)
return the distance to the origin for a given grid point
Definition: mesh.F90:336
subroutine, public mesh_check_dump_compatibility(mesh, dir, filename, namespace, mpi_grp, grid_changed, grid_reordered, map, ierr)
Definition: mesh.F90:590
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:444
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:976
integer function, public mesh_global2local(mesh, ipg)
This function returns the local mesh index for a given global index.
Definition: mesh.F90:967
recursive subroutine, public mesh_end(this)
Definition: mesh.F90:680
logical pure function, public mesh_compact_boundaries(mesh)
Definition: mesh.F90:822
real(real64) pure function, public mesh_local_memory(mesh)
Definition: mesh.F90:790
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)
Definition: mesh.F90:804
subroutine mesh_r_global(mesh, ipg, rr, origin, coords)
return the distance to the origin for a given grid point
Definition: mesh.F90:356
subroutine, public mesh_write_fingerprint(mesh, dir, filename, mpi_grp, namespace, ierr)
Definition: mesh.F90:454
subroutine mesh_init(this)
Definition: mesh.F90:267
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:624
This module contains some common usage patterns of MPI routines.
Definition: mpi_lib.F90:115
Some general things and nomenclature:
Definition: par_vec.F90:171
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
character(len=20) pure function, public units_abbrev(this)
Definition: unit.F90:223
This module defines the unit system, used for input and output.
type(unit_system_t), public units_out
abstract class for basis sets
This data type defines a line, and a regular grid defined on this line (or rather,...
Definition: mesh.F90:255
define a grid on a plane.
Definition: mesh.F90:243
Describes mesh distribution to nodes.
Definition: mesh.F90:186
int true(void)
real(real64) function values(xx)
Definition: test.F90:1867