77 type(derivatives_t) :: der
78 type(stencil_t) :: stencil
79 type(symmetries_t) :: symm
80 type(symmetrizer_t) :: symmetrizer
83 integer,
parameter :: &
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(:,:)
109 type(stencil_t) :: cube
110 integer :: enlarge(space%dim)
113 real(real64) :: grid_spacing(space%dim)
117 gr%box =>
box_factory_create(namespace, space, latt=latt, n_sites=n_sites, site_position=site_position)
119 if (
present(symm))
then
148 if(
parse_block(namespace,
'Spacing', blk) == 0)
then
150 do idir = 1, space%dim
156 grid_spacing = grid_spacing(1)
159 if (
associated(gr%box))
then
160 select type (box => gr%box)
162 do idir = 1, space%dim
164 if (grid_spacing(idir) <
m_zero)
then
165 grid_spacing(idir) = box%pixel_size(idir)
172 message(1) =
" Your input for 'Spacing' is negative or zero."
182 if (
parse_block(namespace,
'PeriodicBoundaryMask', blk) < 0)
then
183 gr%masked_periodic_boundaries = .false.
185 gr%masked_periodic_boundaries = .
true.
203 enlarge = max(enlarge, gr%der%n_ghost)
205 call mesh_init_stage_1(gr, namespace, space, gr%box, gr%coord_system, grid_spacing, enlarge)
214 type(
grid_t),
intent(inout) :: gr
215 type(
grid_t),
intent(in) :: original_gr
217 class(
space_t),
intent(in) :: space
218 class(
box_t),
target,
optional,
intent(in) :: box
219 real(real64),
optional,
intent(in) :: spacing_prefactor(:)
221 integer,
optional,
intent(in) :: n_sites
222 real(real64),
optional,
intent(in) :: site_position(:,:)
225 integer :: enlarge(space%dim)
230 if (.not.
present(box))
then
238 gr%symm = original_gr%symm
240 if (
present(spacing_prefactor))
then
241 gr%spacing = spacing_prefactor * original_gr%spacing
243 write(
message(1),
'(a)')
"The two grids have different spacings, it will not be possible to establish a map between them"
247 gr%spacing = original_gr%spacing
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
264 enlarge = max(enlarge, gr%der%n_ghost)
266 call mesh_init_stage_1(gr, namespace, space, gr%box, gr%coord_system, gr%spacing, enlarge)
278 type(
grid_t),
intent(inout) :: gr
280 class(
space_t),
intent(in) :: space
282 integer,
optional,
intent(in) :: n_sites
283 real(real64),
optional,
intent(in) :: site_position(:,:)
314 call parse_variable(namespace,
'CurvMethod', curv_affine, cv_method)
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')
328 select case (cv_method)
331 gr%box%bounding_box_l(1:space%dim), gr%spacing(1:space%dim))
333 if (
present(site_position) .and.
present(n_sites))
then
334 gr%coord_system =>
curv_gygi_t(namespace, space%dim, n_sites, site_position)
336 message(1) =
"Option CurvMethod = curv_gygi is not currently implemented without ions present."
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))
344 message(1) =
"Option CurvMethod = curv_modine is not currently implemented without ions present."
348 if (
present(latt))
then
349 if (latt%nonorthogonal)
then
352 gr%coord_system =>
cartesian_t(namespace, space%dim)
355 gr%coord_system =>
cartesian_t(namespace, space%dim)
372 type(
grid_t),
target,
intent(inout) :: gr
374 class(
space_t),
intent(in) :: space
376 real(real64),
optional,
intent(in) :: qvector(:)
388 if(space%dim == 3)
then
400 type(
grid_t),
intent(inout) :: gr
403 class(
box_t),
pointer :: box
412 coord_system => gr%coord_system
413 safe_deallocate_p(coord_system)
415 safe_deallocate_p(box)
429 type(
grid_t),
intent(in) :: gr
430 integer,
optional,
intent(in) :: iunit
431 type(
namespace_t),
optional,
intent(in) :: namespace
439 call gr%box%write_info(iunit, namespace)
445 if (gr%use_curvilinear)
then
446 call gr%coord_system%write_info(iunit, namespace)
459 type(
grid_t),
intent(inout) :: gr
460 type(
space_t),
intent(in) :: space
465 real(real64) :: new_box_bounds(2, space%dim)
469 call new_latt%write_info(namespace)
472 select type (coord_system=>gr%coord_system)
474 if (new_latt%nonorthogonal)
then
475 deallocate(gr%coord_system)
481 select type (coord_system=>gr%coord_system)
483 new_box_bounds = gr%box%bounds(coord_system%basis)
484 gr%spacing = (new_box_bounds(2,:)-new_box_bounds(1,:))/gr%idx%ll(:)
490 safe_deallocate_a(gr%x)
491 safe_deallocate_a(gr%vol_pp)
501#include "grid_inc.F90"
504#include "complex.F90"
505#include "grid_inc.F90"
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)...
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
real(real64), parameter, public m_epsilon
real(real64), parameter, public m_one
This module implements the underlying real-space grid.
subroutine, public grid_init_stage_1(gr, namespace, space, symm, latt, n_sites, site_position)
First stage of the grid initialization.
integer, parameter curv_gygi
integer, parameter curv_briggs
subroutine, public zgrid_symmetrize_scalar_field(gr, field, suppress_warning)
subroutine initialize_coordinate_system(gr, namespace, space, latt, n_sites, site_position)
this subroutine initializes the coordinate system
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
subroutine, public zgrid_symmetrize_vector_field(gr, field, suppress_warning)
subroutine, public grid_init_stage_2(gr, namespace, space, mc, qvector)
Second stage of the grid initialization.
subroutine, public dgrid_symmetrize_vector_field(gr, field, suppress_warning)
subroutine, public grid_write_info(gr, iunit, namespace)
subroutine, public dgrid_symmetrize_single(gr, iop, field, symm_field, suppress_warning)
subroutine, public grid_lattice_vectors_update(gr, space, namespace, mc, new_latt)
Regenerate the grid information after update of the lattice vectors.
subroutine, public zgrid_symmetrize_single(gr, iop, field, symm_field, suppress_warning)
subroutine, public grid_end(gr)
finalize a grid object
subroutine, public dgrid_symmetrize_scalar_field(gr, field, suppress_warning)
integer, parameter curv_modine
This module contains subroutines, related to the initialization of the mesh.
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....
subroutine, public mesh_init_stage_1(mesh, namespace, space, box, coord_system, spacing, enlarge)
First stage mesh initialization.
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.
This module defines the meshes, which are used in Octopus.
subroutine, public mesh_write_info(this, iunit, namespace)
recursive subroutine, public mesh_end(this)
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
character(len=512), private msg
subroutine, public messages_warning(no_lines, all_nodes, namespace)
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
subroutine, public messages_new_line()
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
subroutine, public messages_input_error(namespace, var, details, row, column)
subroutine, public messages_experimental(name, namespace)
This module handles the communicators for the various parallelization strategies.
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_)
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.
subroutine, public stencil_end(this)
subroutine, public stencil_union(st1, st2, stu)
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.
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.
class to tell whether a point is inside or outside
abstract class to describe coordinate systems
Description of the grid, containing information on derivatives, stencil, and symmetries.
Describes mesh distribution to nodes.
Stores all communicators and groups.
The class representing the stencil, which is used for non-local mesh operations.