42 use,
intrinsic :: iso_fortran_env
107 integer,
parameter :: &
108 DER_BC_ZERO_F = 0, & !< function is zero at the boundaries
112 integer,
parameter,
public :: &
119 integer,
parameter :: &
127 type(boundaries_t) :: boundaries
128 type(mesh_t),
pointer :: mesh => null()
130 integer,
private :: periodic_dim = 0
132 integer :: stencil_type = 0
135 logical,
private :: stencil_primitive_coordinates
137 class(coordinate_system_t),
private,
pointer :: coord_system
138 logical :: remove_zero_weight_points
140 real(real64),
allocatable :: masses(:)
144 real(real64),
private :: lapl_cutoff =
m_zero
146 type(nl_operator_t),
allocatable,
private :: op(:)
148 type(nl_operator_t),
pointer :: lapl => null()
149 type(nl_operator_t),
pointer :: grad(:) => null()
151 integer,
allocatable :: n_ghost(:)
153 integer,
private :: comm_method = 0
155 type(derivatives_t),
pointer :: finer => null()
156 type(derivatives_t),
pointer :: coarser => null()
157 type(transfer_table_t),
pointer :: to_finer => null()
158 type(transfer_table_t),
pointer :: to_coarser => null()
165 type(par_vec_handle_batch_t) :: pv_h
167 type(derivatives_t),
pointer :: der
168 type(nl_operator_t),
pointer :: op
169 type(batch_t),
pointer :: ff
170 type(batch_t),
pointer :: opff
171 logical :: ghost_update
172 logical :: factor_present
173 real(real64) :: factor
176 type(accel_kernel_t) :: kernel_uvw_xyz, kernel_dcurl, kernel_zcurl
182 type(derivatives_t),
target,
intent(inout) :: der
183 type(namespace_t),
intent(in) :: namespace
184 class(space_t),
intent(in) :: space
185 class(coordinate_system_t),
target,
intent(in) :: coord_system
186 integer,
optional,
intent(in) :: order
189 integer :: default_stencil
190 character(len=40) :: flags
191 logical :: stencil_primitive_coordinates
196 assert(.not. coord_system%local_basis .or. coord_system%orthogonal)
200 der%periodic_dim = space%periodic_dim
201 der%coord_system => coord_system
226 default_stencil = der_star
230 call parse_variable(namespace,
'DerivativesStencil', default_stencil, der%stencil_type)
270 if (
present(order))
then
282 call parse_variable(namespace,
'DerivativesRemoveZeroWeightPoints', .
true., der%remove_zero_weight_points)
323 stencil_primitive_coordinates = .not. coord_system%local_basis
324 call parse_variable(namespace,
'StencilPrimitiveCoordinates', stencil_primitive_coordinates, der%stencil_primitive_coordinates)
327 safe_allocate(der%masses(1:space%dim))
331 safe_allocate(der%op(1:der%dim + 1))
333 der%lapl => der%op(der%dim + 1)
339 safe_allocate(der%n_ghost(1:der%dim))
342 der%n_ghost(idir) = maxval(abs(der%lapl%stencil%points(idir, :)))
347 nullify(der%to_coarser)
348 nullify(der%to_finer)
352 select type (coord_system)
356 if (der%dim > 3)
then
357 message(1) =
"Calculation of derivatives on the GPU with dimension > 3 are only implemented for Cartesian coordinates."
360 write(flags,
'(A,I1.1)')
' -DDIMENSION=', der%dim
380 assert(
allocated(der%op))
382 do idim = 1, der%dim+1
386 safe_deallocate_a(der%masses)
388 safe_deallocate_a(der%n_ghost)
390 safe_deallocate_a(der%op)
391 nullify(der%lapl, der%grad)
395 nullify(der%to_coarser)
396 nullify(der%to_finer)
408 class(
space_t),
intent(in) :: space
410 character(len=80),
optional,
intent(in) :: name
411 integer,
optional,
intent(in) :: order
413 character(len=80) :: name_
425 select case (der%stencil_type)
445 real(real64),
intent(out) :: lapl(:)
449 assert(ubound(lapl, dim=1) >= der%mesh%np)
465 character :: dir_label
469 assert(
associated(der%grad))
479 select case (der%stencil_type)
542 subroutine derivatives_build(der, namespace, space, mesh, qvector, regenerate, verbose)
545 class(
space_t),
intent(in) :: space
546 class(
mesh_t),
target,
intent(in) :: mesh
547 real(real64),
optional,
intent(in) :: qvector(:)
549 logical,
optional,
intent(in) :: regenerate
550 logical,
optional,
intent(in) :: verbose
554 integer :: np_zero_bc
562 assert(
allocated(der%op))
563 assert(der%stencil_type >= der_star .and. der%stencil_type <=
der_stargeneral)
564 assert(.not.(der%stencil_type ==
der_variational .and. mesh%use_curvilinear))
571 if (mesh%use_curvilinear) const_w_ = .false.
578 safe_deallocate_a(der%op(i)%w)
586 call nl_operator_build(space, mesh, der%op(i), der%mesh%np, const_w = const_w_, &
587 regenerate=regenerate)
591 assert(np_zero_bc > mesh%np .and. np_zero_bc <= mesh%np_part)
593 select case (der%stencil_type)
596 do i = 1, der%dim + 1
619 subroutine stencil_pol_grad(stencil_type, dim, direction, order, polynomials)
620 integer,
intent(in) :: stencil_type
621 integer,
intent(in) :: dim
622 integer,
intent(in) :: direction
623 integer,
intent(in) :: order
624 integer,
intent(inout) :: polynomials(:, :)
626 select case (stencil_type)
639 integer,
intent(in) :: stencil_type
641 integer,
intent(in) :: dim
642 integer,
intent(in) :: order
643 integer,
intent(inout) :: polynomials(:, :)
645 select case (stencil_type)
658 subroutine get_rhs_lapl(coord_system, polynomials, chi, rhs, stencil_primitive_coordinates)
660 integer,
intent(in) :: polynomials(:,:)
661 real(real64),
intent(in) :: chi(:)
662 real(real64),
intent(out) :: rhs(:)
663 logical,
intent(in) :: stencil_primitive_coordinates
665 integer :: i, j, k, dim
666 real(real64) :: metric_inverse(1:size(polynomials, dim=1), 1:size(polynomials, dim=1))
667 real(real64) :: hessian_trace(1:size(polynomials, dim=1))
668 integer :: powers(0:2)
672 dim =
size(polynomials, dim=1)
674 if (stencil_primitive_coordinates)
then
676 metric_inverse = coord_system%metric_inverse(chi)
677 hessian_trace = coord_system%trace_hessian(chi)
682 metric_inverse(i, i) =
m_one
689 do j = 1,
size(polynomials, dim=2)
693 if (polynomials(i, j) <= 2)
then
694 powers(polynomials(i, j)) = powers(polynomials(i, j)) + 1
700 if (powers(2) == 1 .and. powers(0) == dim - 1)
then
702 if (polynomials(i, j) == 2)
then
703 rhs(j) =
m_two*metric_inverse(i, i)
709 if (powers(1) == 2 .and. powers(0) == dim - 2)
then
711 if (polynomials(i, j) == 1)
then
713 if (polynomials(k, j) == 1)
then
714 rhs(j) = metric_inverse(i, k) + metric_inverse(k, i)
723 if (powers(1) == 1 .and. powers(0) == dim - 1)
then
725 if (polynomials(i, j) == 1)
then
726 rhs(j) = hessian_trace(i)
737 integer,
intent(in) :: polynomials(:,:)
738 integer,
intent(in) :: dir
739 real(real64),
intent(out) :: rhs(:)
746 dim =
size(polynomials, dim=1)
750 do j = 1,
size(polynomials, dim=2)
753 if (k == dir .and. polynomials(k, j) /= 1) this_one = .false.
754 if (k /= dir .and. polynomials(k, j) /= 0) this_one = .false.
756 if (this_one) rhs(j) =
m_one
765 space, name, verbose, order)
768 integer,
intent(in) :: nderiv
769 integer,
intent(in) :: ideriv
771 type(
space_t),
intent(in) :: space
772 character(len=80),
optional,
intent(in) :: name
773 logical,
optional,
intent(in) :: verbose
774 integer,
optional,
intent(in) :: order
776 integer :: p, p_max, i, j, k, pow_max, order_
777 real(real64) :: x(der%dim)
778 real(real64),
allocatable :: mat(:,:), sol(:,:), powers(:,:)
779 integer,
allocatable :: pol(:,:)
780 real(real64),
allocatable :: rhs(:,:)
781 character(len=80) :: name_
787 safe_allocate(rhs(1:op(1)%stencil%size, nderiv))
788 safe_allocate(pol(1:der%dim, 1:op(1)%stencil%size))
790 if (nderiv == 1)
then
791 if (ideriv <= der%dim)
then
795 call stencil_pol_lapl(der%stencil_type, op(1)%stencil, der%dim, order_, pol)
798 else if (nderiv == der%dim)
then
805 message(1) =
"Error: derivatives_make_discretization can only be called with nderiv = 1 or der%dim."
809 safe_allocate(mat(1:op(1)%stencil%size, 1:op(1)%stencil%size))
810 safe_allocate(sol(1:op(1)%stencil%size, 1:nderiv))
814 message(1) =
'Info: Generating weights for finite-difference discretization of ' // trim(name_)
819 pow_max = maxval(pol)
820 safe_allocate(powers(1:der%dim, 0:pow_max))
825 if (op(1)%const_w) p_max = 1
830 if (nderiv == 1)
then
831 if (ideriv <= der%dim)
then
835 call get_rhs_lapl(der%mesh%coord_system, pol, der%mesh%chi(:, p), rhs(:,1), der%stencil_primitive_coordinates)
837 else if (nderiv == der%dim)
then
848 do i = 1, op(1)%stencil%size
849 if (ideriv <= der%dim)
then
851 x = real(op(1)%stencil%points(:, i), real64)*der%mesh%spacing
854 if (der%stencil_primitive_coordinates)
then
855 x = real(op(1)%stencil%points(:, i), real64)*der%mesh%spacing
857 x = der%mesh%x(:, p + op(1)%ri(i, op(1)%rimap(p))) - der%mesh%x(:, p)
862 x = x*
sqrt(der%masses)
867 powers(:, k) = x*powers(:, k-1)
872 do j = 2, op(1)%stencil%size
873 mat(j, i) = powers(1, pol(1, j))
875 mat(j, i) = mat(j, i)*powers(k, pol(k, j))
887 op(i)%w(:, p) = sol(:, i)
893 if (der%remove_zero_weight_points)
then
897 if (op(i)%const_w)
then
910 safe_deallocate_a(mat)
911 safe_deallocate_a(sol)
912 safe_deallocate_a(powers)
913 safe_deallocate_a(pol)
914 safe_deallocate_a(rhs)
921 logical function derivatives_overlap(this)
result(overlap)
924 push_sub(derivatives_overlap)
926 overlap = this%comm_method /= blocking
928 pop_sub(derivatives_overlap)
929 end function derivatives_overlap
937 class(
space_t),
intent(in) :: space
938 character(len=80),
intent(in) :: name
939 integer,
intent(in) :: order
944 call nl_operator_build(space, this%mesh, op(1), this%mesh%np, const_w = .not. this%mesh%use_curvilinear)
966 logical :: mask(1:this%mesh%np)
967 integer :: ip, is, index
971 do ip = 1, this%mesh%np
973 do is = 1, this%lapl%stencil%size
977 if (index > this%mesh%np + this%mesh%pv%np_ghost)
then
991 real(real64) function derivatives_lapl_get_max_eigenvalue(this)
996 push_sub(derivatives_lapl_get_max_eigenvalue)
998 derivatives_lapl_get_max_eigenvalue =
m_zero
999 if ((this%stencil_type == der_star .or. this%stencil_type ==
der_stargeneral) &
1000 .and. .not. this%mesh%use_curvilinear)
then
1002 do i = 1, this%lapl%stencil%size
1003 derivatives_lapl_get_max_eigenvalue = derivatives_lapl_get_max_eigenvalue + &
1004 (-1)**maxval(abs(this%lapl%stencil%points(:, i)))*this%lapl%w(i, 1)
1006 derivatives_lapl_get_max_eigenvalue = abs(derivatives_lapl_get_max_eigenvalue)
1010 derivatives_lapl_get_max_eigenvalue = derivatives_lapl_get_max_eigenvalue + &
1011 m_pi**2/this%mesh%spacing(i)**2
1015 pop_sub(derivatives_lapl_get_max_eigenvalue)
1020 class(coordinate_system_t),
target,
intent(in) :: coord_system
1022 der%coord_system => coord_system
1027#include "derivatives_inc.F90"
1030#include "complex.F90"
1031#include "derivatives_inc.F90"
double sqrt(double __x) __attribute__((__nothrow__
subroutine, public accel_kernel_build(this, file_name, kernel_name, flags)
pure logical function, public accel_is_enabled()
This module implements batches of mesh functions.
This module implements common operations on batches of mesh functions.
Module implementing boundary conditions in Octopus.
subroutine, public boundaries_end(this)
subroutine, public boundaries_init(this, namespace, space, mesh, qvector)
initialize the boundary contitions
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 zderivatives_curl(der, ff, op_ff, ghost_update, set_bc)
apply the curl operator to a vector of mesh functions
subroutine, public dderivatives_perform(op, der, ff, op_ff, ghost_update, set_bc, factor)
apply a nl_operator to a mesh function
subroutine, public derivatives_lapl_diag(der, lapl)
Returns the diagonal elements of the Laplacian, needed for preconditioning.
integer, parameter, public der_cube
subroutine, public dderivatives_batch_curl_from_gradient(der, ffb, gradb)
calculate the curl from a batch and its gradient
subroutine derivatives_make_discretization(namespace, der, nderiv, ideriv, op, space, name, verbose, order)
subroutine, public derivatives_build(der, namespace, space, mesh, qvector, regenerate, verbose)
build the derivatives object:
subroutine, public dderivatives_batch_grad(der, ffb, opffb, ghost_update, set_bc, to_cartesian, factor)
apply the gradient to a batch of mesh functions
subroutine stencil_pol_grad(stencil_type, dim, direction, order, polynomials)
subroutine stencil_pol_lapl(stencil_type, stencil, dim, order, polynomials)
subroutine, public zderivatives_batch_curl_from_gradient(der, ffb, gradb)
calculate the curl from a batch and its gradient
subroutine, public zderivatives_batch_div(der, ffb, opffb, ghost_update, set_bc, to_cartesian)
calculate the divergence of a vector of batches
subroutine, public dderivatives_partial(der, ff, op_ff, dir, ghost_update, set_bc)
apply the partial derivative along dir to a mesh function
subroutine get_rhs_lapl(coord_system, polynomials, chi, rhs, stencil_primitive_coordinates)
subroutine, public dderivatives_batch_perform(op, der, ff, opff, ghost_update, set_bc, factor, async)
apply an operator to a bach of mesh functions
subroutine, public zderivatives_div(der, ff, op_ff, ghost_update, set_bc, to_cartesian)
apply the divergence operator to a vector of mesh functions
integer, parameter der_bc_period
boundary is periodic
subroutine, public dderivatives_batch_finish(handle, async)
apply a non-local operator to a batch (2nd part)
subroutine, public dderivatives_curl(der, ff, op_ff, ghost_update, set_bc)
apply the curl operator to a vector of mesh functions
subroutine derivatives_get_stencil_grad(der)
subroutine, public dderivatives_grad(der, ff, op_ff, ghost_update, set_bc, to_cartesian)
apply the gradient to a mesh function
subroutine derivatives_get_stencil_lapl(der, lapl, space, coord_system, name, order)
subroutine, public derivatives_init(der, namespace, space, coord_system, order)
subroutine, public derivatives_end(der)
integer, parameter der_bc_zero_df
first derivative of the function is zero
subroutine, public zderivatives_batch_grad(der, ffb, opffb, ghost_update, set_bc, to_cartesian, factor)
apply the gradient to a batch of mesh functions
subroutine, public dderivatives_batch_div(der, ffb, opffb, ghost_update, set_bc, to_cartesian)
calculate the divergence of a vector of batches
subroutine, public derivatives_get_lapl(this, namespace, op, space, name, order)
subroutine, public zderivatives_batch_perform(op, der, ff, opff, ghost_update, set_bc, factor, async)
apply an operator to a bach of mesh functions
subroutine, public dderivatives_batch_start(op, der, ff, opff, handle, ghost_update, set_bc, factor)
apply a non-local operator to a batch (1st part)
subroutine get_rhs_grad(polynomials, dir, rhs)
subroutine, public dderivatives_lapl(der, ff, op_ff, ghost_update, set_bc, factor)
apply the Laplacian to a mesh function
real(real64) function, public derivatives_lapl_get_max_eigenvalue(this)
Get maximum eigenvalue of discrete Laplacian. For the star and star_general stencils,...
subroutine, public dderivatives_batch_curl(der, ffb, ghost_update, set_bc)
apply the curl to a batch of mesh functions
subroutine, public zderivatives_grad(der, ff, op_ff, ghost_update, set_bc, to_cartesian)
apply the gradient to a mesh function
subroutine, public zderivatives_partial(der, ff, op_ff, dir, ghost_update, set_bc)
apply the partial derivative along dir to a mesh function
subroutine, public zderivatives_batch_start(op, der, ff, opff, handle, ghost_update, set_bc, factor)
apply a non-local operator to a batch (1st part)
logical function, dimension(1:this%mesh%np), public derivatives_get_inner_boundary_mask(this)
This function tells whether a point in the grid is contained in a layer of the width of the stencil b...
integer, parameter, public der_starplus
subroutine, public zderivatives_perform(op, der, ff, op_ff, ghost_update, set_bc, factor)
apply a nl_operator to a mesh function
subroutine, public zderivatives_batch_curl(der, ffb, ghost_update, set_bc)
apply the curl to a batch of mesh functions
integer, parameter, public der_variational
subroutine, public zderivatives_lapl(der, ff, op_ff, ghost_update, set_bc, factor)
apply the Laplacian to a mesh function
subroutine, public derivates_set_coordinates_system(der, coord_system)
subroutine, public zderivatives_batch_finish(handle, async)
apply a non-local operator to a batch (2nd part)
subroutine, public dderivatives_div(der, ff, op_ff, ghost_update, set_bc, to_cartesian)
apply the divergence operator to a vector of mesh functions
integer, parameter non_blocking
integer, parameter, public der_stargeneral
real(real64), parameter, public m_two
real(real64), parameter, public m_zero
real(real64), parameter, public m_pi
some mathematical constants
real(real64), parameter, public m_one
This module is intended to contain "only mathematical" functions and procedures.
This module defines the meshes, which are used in Octopus.
subroutine, public messages_obsolete_variable(namespace, name, rep)
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_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
This module defines non-local operators.
subroutine, public dnl_operator_operate_diag(op, fo)
subroutine, public nl_operator_init(op, label, symm)
initialize an instance of a non-local operator by setting the label
subroutine, public nl_operator_build_symmetric_weights(op, max_size)
Builds (or rebuild) the necessary arrays for symmetric and antisymmetric stencils.
subroutine, public nl_operator_update_gpu_buffers(op)
subroutine, public nl_operator_output_weights(this)
subroutine, public nl_operator_end(op)
subroutine, public nl_operator_remove_zero_weight_points(op, space, mesh)
Removes the zero-weight points for constant weight stencils.
integer, parameter, public op_symmetric
subroutine, public nl_operator_build(space, mesh, op, np, const_w, regenerate)
Creates the nonlocal operators for the stencils used for finite differences.
integer pure function, public nl_operator_np_zero_bc(op)
integer, parameter, public op_antisymmetric
integer pure function, public nl_operator_get_index(op, is, ip)
subroutine, public nl_operator_allocate_gpu_buffers(op)
Some general things and nomenclature:
This module defines routines, generating operators for a cubic stencil.
subroutine, public stencil_cube_get_lapl(this, dim, order)
subroutine, public stencil_cube_polynomials_lapl(dim, order, pol)
subroutine, public stencil_cube_get_grad(this, dim, order)
This module defines stencils used in Octopus.
This module defines routines, generating operators for a star stencil.
subroutine, public stencil_star_get_grad(this, dim, dir, order)
subroutine, public stencil_star_get_lapl(this, dim, order)
subroutine, public stencil_star_polynomials_grad(dir, order, pol)
subroutine, public stencil_star_polynomials_lapl(dim, order, pol)
This module defines routines, generating operators for a generalized star stencil.
subroutine, public stencil_stargeneral_get_arms(this, dim, coord_system)
Finds the direction of the arms of the star-general stencil as described in Natan et al....
subroutine, public stencil_stargeneral_pol_lapl(this, dim, order, pol)
subroutine, public stencil_stargeneral_get_lapl(this, dim, order)
This module defines routines, generating operators for a stencil consisting of a star and a cross....
subroutine, public stencil_starplus_pol_grad(dim, dir, order, pol)
subroutine, public stencil_starplus_pol_lapl(dim, order, pol)
subroutine, public stencil_starplus_get_lapl(this, dim, order)
subroutine, public stencil_starplus_get_grad(this, dim, dir, order)
Implements the variational discretization of the Laplacian as proposed by P. Maragakis,...
subroutine, public stencil_variational_coeff_lapl(dim, order, h, lapl, alpha)
This module is intended to contain simple general-purpose utility functions and procedures.
character pure function, public index2axis(idir)
abstract class to describe coordinate systems
handle to transfer data from the start() to finish() calls.
class representing derivatives
Describes mesh distribution to nodes.
data type for non local operators
The class representing the stencil, which is used for non-local mesh operations.