Octopus
grid.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
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
23
24module grid_oct_m
26 use box_oct_m
31 use cube_oct_m
36 use global_oct_m
38 use mesh_oct_m
41 use mpi_oct_m
45 use parser_oct_m
47 use space_oct_m
52 use unit_oct_m
55
56 implicit none
57
58 private
59 public :: &
60 grid_t, &
64 grid_end, &
73
75 type, extends(mesh_t) :: grid_t
76 ! Components are public by default
77 type(derivatives_t) :: der
78 type(stencil_t) :: stencil
79 type(symmetries_t) :: symm
80 type(symmetrizer_t) :: symmetrizer
81 end type grid_t
82
83 integer, parameter :: &
84 CURV_AFFINE = 1, &
85 curv_gygi = 2, &
86 curv_briggs = 3, &
87 curv_modine = 4
88
89contains
90
91 !-------------------------------------------------------------------
100 subroutine grid_init_stage_1(gr, namespace, space, symm, latt, n_sites, site_position)
101 type(grid_t), intent(inout) :: gr
102 type(namespace_t), intent(in) :: namespace
103 class(space_t), intent(in) :: space
104 type(symmetries_t), optional, intent(in) :: symm
105 type(lattice_vectors_t), optional, intent(in) :: latt
106 integer, optional, intent(in) :: n_sites
107 real(real64), optional, intent(in) :: site_position(:,:)
108
109 type(stencil_t) :: cube
110 integer :: enlarge(space%dim)
111 type(block_t) :: blk
112 integer :: idir
113 real(real64) :: grid_spacing(space%dim)
114
115 push_sub(grid_init_stage_1)
116
117 gr%box => box_factory_create(namespace, space, latt=latt, n_sites=n_sites, site_position=site_position)
118
119 if (present(symm)) then
120 gr%symm = symm
121 end if
122
123 !%Variable Spacing
124 !%Type float
125 !%Section Mesh
126 !%Description
127 !% The spacing between the points in the mesh. This controls the
128 !% quality of the discretization: smaller spacing gives more
129 !% precise results but increased computational cost.
130 !%
131 !% When using curvilinear coordinates, this is a canonical spacing
132 !% that will be changed locally by the transformation. In periodic
133 !% directions, your spacing may be slightly different than what
134 !% you request here, since the box size must be an integer
135 !% multiple of the spacing.
136 !%
137 !% The default value is defined by the image resolution if <tt>BoxShape =
138 !% box_image</tt>. Othewise here is no default otherwise.
139 !%
140 !% It is possible to have a different spacing in each one of the Cartesian directions
141 !% if we define <tt>Spacing</tt> as block of the form
142 !%
143 !% <tt>%Spacing
144 !% <br>&nbsp;&nbsp;spacing_x | spacing_y | spacing_z
145 !% <br>%</tt>
146 !%End
147
148 if(parse_block(namespace, 'Spacing', blk) == 0) then
149 if(parse_block_cols(blk,0) < space%dim) call messages_input_error(namespace, 'Spacing')
150 do idir = 1, space%dim
151 call parse_block_float(blk, 0, idir - 1, grid_spacing(idir), units_inp%length)
152 end do
153 call parse_block_end(blk)
154 else
155 call parse_variable(namespace, 'Spacing', -m_one, grid_spacing(1), units_inp%length)
156 grid_spacing = grid_spacing(1)
157 end if
158
159 if (associated(gr%box)) then
160 select type (box => gr%box)
161 type is (box_image_t)
162 do idir = 1, space%dim
163 ! default grid_spacing is determined from the pixel size such that one grid point = one pixel.
164 if (grid_spacing(idir) < m_zero) then
165 grid_spacing(idir) = box%pixel_size(idir)
166 end if
167 end do
168 end select
169 end if
171 if (any(grid_spacing(1:space%dim) < m_epsilon)) then
172 message(1) = " Your input for 'Spacing' is negative or zero."
173 call messages_fatal(1, namespace=namespace)
174 end if
175
176 !%Variable PeriodicBoundaryMask
177 !%Type block
178 !%Section Mesh
179 !%Description
180 !% (Experimental) Defines a mask for which periodic boundaries are replaced by zero boundary conditions.
181 !%End
182 if (parse_block(namespace, 'PeriodicBoundaryMask', blk) < 0) then
183 gr%masked_periodic_boundaries = .false.
184 else
185 gr%masked_periodic_boundaries = .true.
186 call parse_block_string(blk, 0, 0, gr%periodic_boundary_mask)
187 call messages_experimental('PeriodicBoundaryMask')
188 call parse_block_end(blk)
189 end if
190
191 ! Initialize coordinate system
192 call initialize_coordinate_system(gr, namespace, space, latt, n_sites, site_position)
194 ! initialize derivatives
195 call derivatives_init(gr%der, namespace, space, gr%coord_system)
196 ! the stencil used to generate the grid is a union of a cube (for
197 ! multigrid) and the Laplacian
198 call stencil_cube_get_lapl(cube, space%dim, order = 2)
199 call stencil_union(cube, gr%der%lapl%stencil, gr%stencil)
200 call stencil_end(cube)
201
202 enlarge = 2
203 enlarge = max(enlarge, gr%der%n_ghost)
204
205 call mesh_init_stage_1(gr, namespace, space, gr%box, gr%coord_system, grid_spacing, enlarge)
206 call mesh_init_stage_2(gr, namespace, space, gr%box, gr%stencil)
207
208 pop_sub(grid_init_stage_1)
209 end subroutine grid_init_stage_1
210
211 !--------------------------------------------------------------------
213 subroutine grid_init_from_grid_stage_1(gr, original_gr, namespace, space, box, spacing_prefactor, latt, n_sites, site_position)
214 type(grid_t), intent(inout) :: gr
215 type(grid_t), intent(in) :: original_gr
216 type(namespace_t), intent(in) :: namespace
217 class(space_t), intent(in) :: space
218 class(box_t), target, optional, intent(in) :: box
219 real(real64), optional, intent(in) :: spacing_prefactor(:)
220 type(lattice_vectors_t), optional, intent(in) :: latt
221 integer, optional, intent(in) :: n_sites
222 real(real64), optional, intent(in) :: site_position(:,:)
223
224 type(stencil_t) :: cube
225 integer :: enlarge(space%dim)
226
228
229 ! First of all, create and associate the box that contains the grid
230 if (.not. present(box)) then
231 gr%box => box_factory_create(namespace, space)
232 else
233 gr%box => box
234 end if
235
236 !!! Copy all data !!!
237 ! Symmetries
238 gr%symm = original_gr%symm
239 ! Spacing
240 if (present(spacing_prefactor)) then
241 gr%spacing = spacing_prefactor * original_gr%spacing
242 if (any(abs(gr%spacing - m_one) > m_epsilon)) then
243 write(message(1), '(a)') "The two grids have different spacings, it will not be possible to establish a map between them"
244 call messages_warning(1)
245 end if
246 else
247 gr%spacing = original_gr%spacing
248 end if
249 ! Periodic boundaries
250 gr%masked_periodic_boundaries = original_gr%masked_periodic_boundaries
251 if (gr%masked_periodic_boundaries) gr%periodic_boundary_mask = original_gr%periodic_boundary_mask
252 ! Coordinate System
253 call initialize_coordinate_system(gr, namespace, space, latt, n_sites, site_position)
254
255 ! initialize derivatives
256 call derivatives_init(gr%der, namespace, space, gr%coord_system)
257 ! the stencil used to generate the grid is a union of a cube (for
258 ! multigrid) and the Laplacian
259 call stencil_cube_get_lapl(cube, space%dim, order = 2)
260 call stencil_union(cube, gr%der%lapl%stencil, gr%stencil)
261 call stencil_end(cube)
262
263 enlarge = 2
264 enlarge = max(enlarge, gr%der%n_ghost)
265
266 call mesh_init_stage_1(gr, namespace, space, gr%box, gr%coord_system, gr%spacing, enlarge)
267 call mesh_init_stage_2(gr, namespace, space, gr%box, gr%stencil)
268
270 end subroutine grid_init_from_grid_stage_1
271
272
273 !--------------------------------------------------------------------
277 subroutine initialize_coordinate_system(gr, namespace, space, latt, n_sites, site_position)
278 type(grid_t), intent(inout) :: gr
279 type(namespace_t), intent(in) :: namespace
280 class(space_t), intent(in) :: space
281 type(lattice_vectors_t), optional, intent(in) :: latt
282 integer, optional, intent(in) :: n_sites
283 real(real64), optional, intent(in) :: site_position(:,:)
284
285 integer :: cv_method
286
288
289 !%Variable CurvMethod
290 !%Type integer
291 !%Default curv_uniform
292 !%Section Mesh::Curvilinear
293 !%Description
294 !% The relevant functions in octopus are represented on a mesh in real space.
295 !% This mesh may be an evenly spaced regular rectangular grid (standard mode),
296 !% or else an adaptive or curvilinear grid. We have implemented
297 !% three kinds of adaptive meshes, although only one is currently working,
298 !% the one invented by F. Gygi (<tt>curv_gygi</tt>). The code will stop if any of
299 !% the other two is invoked. All are experimental with domain parallelization.
300 !%Option curv_affine 1
301 !% Regular, uniform rectangular grid.
302 !%Option curv_gygi 2
303 !% The deformation of the grid is done according to the scheme described by
304 !% F. Gygi [F. Gygi and G. Galli, <i>Phys. Rev. B</i> <b>52</b>, R2229 (1995)].
305 !%Option curv_briggs 3
306 !% The deformation of the grid is done according to the scheme described by
307 !% Briggs [E.L. Briggs, D.J. Sullivan, and J. Bernholc, <i>Phys. Rev. B</i> <b>54</b> 14362 (1996)]
308 !% (NOT WORKING).
309 !%Option curv_modine 4
310 !% The deformation of the grid is done according to the scheme described by
311 !% Modine [N.A. Modine, G. Zumbach and E. Kaxiras, <i>Phys. Rev. B</i> <b>55</b>, 10289 (1997)]
312 !% (NOT WORKING).
313 !%End
314 call parse_variable(namespace, 'CurvMethod', curv_affine, cv_method)
315 if (.not. varinfo_valid_option('CurvMethod', cv_method)) call messages_input_error(namespace, 'CurvMethod')
316 if (present(latt)) then
317 if (cv_method /= curv_affine .and. latt%nonorthogonal) then
318 call messages_input_error(namespace, 'CurvMethod', 'Curvilinear coordinates with non-orthogonal cells are not implemented')
319 end if
320 end if
321 call messages_print_var_option("CurvMethod", cv_method, namespace=namespace)
322
323 ! FIXME: The other two methods are apparently not working
324 if (cv_method > curv_gygi) then
325 call messages_experimental('Selected curvilinear coordinates method')
326 end if
327
328 select case (cv_method)
329 case (curv_briggs)
330 gr%coord_system => curv_briggs_t(namespace, space%dim, &
331 gr%box%bounding_box_l(1:space%dim), gr%spacing(1:space%dim))
332 case (curv_gygi)
333 if (present(site_position) .and. present(n_sites)) then
334 gr%coord_system => curv_gygi_t(namespace, space%dim, n_sites, site_position)
335 else
336 message(1) = "Option CurvMethod = curv_gygi is not currently implemented without ions present."
337 call messages_fatal(1, namespace=namespace)
338 end if
339 case (curv_modine)
340 if (present(site_position) .and. present(n_sites)) then
341 gr%coord_system => curv_modine_t(namespace, space%dim, n_sites, site_position, &
342 gr%box%bounding_box_l(1:space%dim), gr%spacing(1:space%dim))
343 else
344 message(1) = "Option CurvMethod = curv_modine is not currently implemented without ions present."
345 call messages_fatal(1, namespace=namespace)
346 end if
347 case (curv_affine)
348 if (present(latt)) then
349 if (latt%nonorthogonal) then
350 gr%coord_system => affine_coordinates_t(namespace, space%dim, latt%rlattice_primitive)
351 else
352 gr%coord_system => cartesian_t(namespace, space%dim)
353 end if
354 else
355 gr%coord_system => cartesian_t(namespace, space%dim)
356 end if
357 end select
358
360 end subroutine initialize_coordinate_system
361
362
363 !-------------------------------------------------------------------
371 subroutine grid_init_stage_2(gr, namespace, space, mc, qvector)
372 type(grid_t), target, intent(inout) :: gr
373 type(namespace_t), intent(in) :: namespace
374 class(space_t), intent(in) :: space
375 type(multicomm_t), intent(in) :: mc
376 real(real64), optional, intent(in) :: qvector(:)
377
378 push_sub(grid_init_stage_2)
379
380 call mesh_init_stage_3(gr, namespace, space, gr%stencil, mc)
381
382 call nl_operator_global_init(namespace)
383 call derivatives_build(gr%der, namespace, space, gr, qvector)
384
385 ! print info concerning the grid
386 call grid_write_info(gr, namespace=namespace)
387
388 if(space%dim == 3) then
389 call symmetrizer_init(gr%symmetrizer, gr, gr%symm)
390 end if
391
392 pop_sub(grid_init_stage_2)
393 end subroutine grid_init_stage_2
394
395
396 !-------------------------------------------------------------------
399 subroutine grid_end(gr)
400 type(grid_t), intent(inout) :: gr
401
402 class(coordinate_system_t), pointer :: coord_system
403 class(box_t), pointer :: box
404
405 push_sub(grid_end)
406
408
409 call derivatives_end(gr%der)
410
411 ! We need to take a pointer here, otherwise we run into a gfortran bug.
412 coord_system => gr%coord_system
413 safe_deallocate_p(coord_system)
414 box => gr%box
415 safe_deallocate_p(box)
416
417 call mesh_end(gr)
418
419 call stencil_end(gr%stencil)
420
421 call symmetrizer_end(gr%symmetrizer)
422
423 pop_sub(grid_end)
424 end subroutine grid_end
425
426
427 !-------------------------------------------------------------------
428 subroutine grid_write_info(gr, iunit, namespace)
429 type(grid_t), intent(in) :: gr
430 integer, optional, intent(in) :: iunit
431 type(namespace_t), optional, intent(in) :: namespace
432
433 push_sub(grid_write_info)
434
435 call messages_print_with_emphasis(msg="Grid", iunit=iunit, namespace=namespace)
436
437 message(1) = 'Simulation Box:'
438 call messages_info(1, iunit=iunit, namespace=namespace)
439 call gr%box%write_info(iunit, namespace)
440
441 message(1) = "Main mesh:"
442 call messages_info(1, iunit=iunit, namespace=namespace)
443 call mesh_write_info(gr, iunit, namespace)
444
445 if (gr%use_curvilinear) then
446 call gr%coord_system%write_info(iunit, namespace)
447 end if
448
449 call messages_print_with_emphasis(iunit=iunit, namespace=namespace)
450
451 pop_sub(grid_write_info)
452 end subroutine grid_write_info
453
454 !-------------------------------------------------------------------
458 subroutine grid_lattice_vectors_update(gr, space, namespace, mc, new_latt)
459 type(grid_t), intent(inout) :: gr
460 type(space_t), intent(in) :: space
461 type(namespace_t), intent(in) :: namespace
462 type(multicomm_t), intent(in) :: mc
463 type(lattice_vectors_t), intent(in) :: new_latt
465 real(real64) :: new_box_bounds(2, space%dim)
466
468
469 call new_latt%write_info(namespace)
470
471 ! Regenerate the coordinate system
472 select type (coord_system=>gr%coord_system)
473 type is (affine_coordinates_t)
474 if (new_latt%nonorthogonal) then
475 deallocate(gr%coord_system)
476 gr%coord_system => affine_coordinates_t(namespace, space%dim, new_latt%rlattice_primitive)
477 end if
478 end select
479
480 ! Regenerate the grid keeping the same number of points
481 select type (coord_system=>gr%coord_system)
482 class is (affine_coordinates_t)
483 new_box_bounds = gr%box%bounds(coord_system%basis)
484 gr%spacing = (new_box_bounds(2,:)-new_box_bounds(1,:))/gr%idx%ll(:)
485 write(message(1), '(a)') trim(gr%box%short_info(unit_angstrom))
486 call messages_info(1, namespace=namespace)
487 call messages_new_line()
488 end select
489
490 safe_deallocate_a(gr%x)
491 safe_deallocate_a(gr%vol_pp)
492 call mesh_init_stage_3(gr, namespace, space, gr%stencil, mc, regenerate=.true.)
493
494 call derivatives_build(gr%der, namespace, space, gr, regenerate=.true., verbose=.false.)
495
497 end subroutine grid_lattice_vectors_update
498
499#include "undef.F90"
500#include "real.F90"
501#include "grid_inc.F90"
502
503#include "undef.F90"
504#include "complex.F90"
505#include "grid_inc.F90"
506
507end module grid_oct_m
508
509!! Local Variables:
510!! mode: f90
511!! coding: utf-8
512!! End:
Module, implementing a factory for boxes.
class(box_t) function, pointer, public box_factory_create(namespace, space, latt, n_sites, site_position)
initialize a box of any type
This module implements the curvilinear coordinates given in E.L. Briggs, D.J. Sullivan,...
This module implements the curvilinear coordinates given in F. Gygi and G. Galli, PRB 52 R2229 (1996)...
Definition: curv_gygi.F90:117
This module implements the curvilinear coordinates given in N. A. Modine, G. Zumbach,...
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
subroutine, public derivatives_build(der, namespace, space, mesh, qvector, regenerate, verbose)
build the derivatives object:
subroutine, public derivatives_init(der, namespace, space, coord_system, order)
subroutine, public derivatives_end(der)
real(real64), parameter, public m_zero
Definition: global.F90:187
real(real64), parameter, public m_epsilon
Definition: global.F90:203
real(real64), parameter, public m_one
Definition: global.F90:188
This module implements the underlying real-space grid.
Definition: grid.F90:117
subroutine, public grid_init_stage_1(gr, namespace, space, symm, latt, n_sites, site_position)
First stage of the grid initialization.
Definition: grid.F90:194
integer, parameter curv_gygi
Definition: grid.F90:176
integer, parameter curv_briggs
Definition: grid.F90:176
subroutine, public zgrid_symmetrize_scalar_field(gr, field, suppress_warning)
Definition: grid.F90:805
subroutine initialize_coordinate_system(gr, namespace, space, latt, n_sites, site_position)
this subroutine initializes the coordinate system
Definition: grid.F90:371
subroutine, public grid_init_from_grid_stage_1(gr, original_gr, namespace, space, box, spacing_prefactor, latt, n_sites, site_position)
this subroutine allows to create a grid from an existing grid
Definition: grid.F90:307
subroutine, public zgrid_symmetrize_vector_field(gr, field, suppress_warning)
Definition: grid.F90:831
subroutine, public grid_init_stage_2(gr, namespace, space, mc, qvector)
Second stage of the grid initialization.
Definition: grid.F90:465
subroutine, public dgrid_symmetrize_vector_field(gr, field, suppress_warning)
Definition: grid.F90:687
subroutine, public grid_write_info(gr, iunit, namespace)
Definition: grid.F90:522
subroutine, public dgrid_symmetrize_single(gr, iop, field, symm_field, suppress_warning)
Definition: grid.F90:715
subroutine, public grid_lattice_vectors_update(gr, space, namespace, mc, new_latt)
Regenerate the grid information after update of the lattice vectors.
Definition: grid.F90:552
subroutine, public zgrid_symmetrize_single(gr, iop, field, symm_field, suppress_warning)
Definition: grid.F90:859
subroutine, public grid_end(gr)
finalize a grid object
Definition: grid.F90:493
subroutine, public dgrid_symmetrize_scalar_field(gr, field, suppress_warning)
Definition: grid.F90:661
integer, parameter curv_modine
Definition: grid.F90:176
This module contains subroutines, related to the initialization of the mesh.
Definition: mesh_init.F90:117
subroutine, public mesh_init_stage_3(mesh, namespace, space, stencil, mc, parent, regenerate)
When running parallel in domains, stencil and np_stencil are needed to compute the ghost points....
Definition: mesh_init.F90:538
subroutine, public mesh_init_stage_1(mesh, namespace, space, box, coord_system, spacing, enlarge)
First stage mesh initialization.
Definition: mesh_init.F90:170
subroutine, public mesh_init_stage_2(mesh, namespace, space, box, stencil, regenerate)
This subroutine creates the global array of spatial indices and the inverse mapping.
Definition: mesh_init.F90:290
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public mesh_write_info(this, iunit, namespace)
Definition: mesh.F90:304
recursive subroutine, public mesh_end(this)
Definition: mesh.F90:680
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:930
character(len=512), private msg
Definition: messages.F90:165
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:543
subroutine, public messages_new_line()
Definition: messages.F90:1146
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:420
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:723
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1097
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:624
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:145
This module defines non-local operators.
subroutine, public nl_operator_global_init(namespace)
initialize global settings for non-local operators
subroutine, public nl_operator_global_end()
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618
This module defines routines, generating operators for a cubic stencil.
subroutine, public stencil_cube_get_lapl(this, dim, order)
This module defines stencils used in Octopus.
Definition: stencil.F90:135
subroutine, public stencil_end(this)
Definition: stencil.F90:214
subroutine, public stencil_union(st1, st2, stu)
Definition: stencil.F90:226
subroutine, public symmetrizer_end(this)
subroutine, public symmetrizer_init(this, mesh, symm)
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
This module defines the unit system, used for input and output.
type(unit_t), public unit_angstrom
For XYZ files.
type(unit_system_t), public units_inp
the units systems for reading and writing
Class implementing a box generated from a 2D image.
Definition: box_image.F90:133
class to tell whether a point is inside or outside
Definition: box.F90:141
abstract class to describe coordinate systems
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:168
Describes mesh distribution to nodes.
Definition: mesh.F90:186
Stores all communicators and groups.
Definition: multicomm.F90:206
The class representing the stencil, which is used for non-local mesh operations.
Definition: stencil.F90:163
int true(void)