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