Octopus
helmholtz_decomposition_m Module Reference

The Helmholtz decomposition is intended to contain "only mathematical" functions and procedures to compute the Helmholtz decomposition of a generic field. More...

Detailed Description

The Helmholtz decomposition is intended to contain "only mathematical" functions and procedures to compute the Helmholtz decomposition of a generic field.

Given a generic vector field F, it is possibile to define a curl-free component div(phi) and a divergence free component curl(A) such that \( \mathbf{F} = -div(phi) + curl(A) \), where

\[ \phi = \frac{1}{4\pi} \left( \int \frac{\nabla \cdot \mathbf{F}}{|\mathbf{r} - \mathbf{r_1}|} \,dr_1 - \oint \mathbf{n} \cdot \frac{\mathbf{F}}{|\mathbf{r} - \mathbf{r_1}|} \,dS_1 \right) \]

\[\mathbf{A} = \frac{1}{4\pi} \left( \int \frac{\nabla \times \mathbf{F}}{|\mathbf{r} - \mathbf{r_1}|} \,dr_1 - \oint \mathbf{n} \times \frac{\mathbf{F}}{|\mathbf{r} - \mathbf{r_1}|} \,dS_1 \right) \]

We start from a field defined on the system grid, i.e. the grid used for the Maxwell simulation.

VECTOR and SCALAR POTENTIAL

If the user asks to compute the vector potential A or the scalar potential phi, then they are computed on the system grid. All the points up to system_gridnp will have physical meaning

TRANSVERSE and LONGITUDINAL FIELDS

If the user asks to compute the transverse or the longitudinal field, then they are computed on the system grid. However, the points the points contained in a layer of the width of the stencil between the innermost layer of points and last row of points in the grid will be set to zero. This is a visual representation for stencil = 2: 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 0 0 0 1 1 1 1 0 0 0 1 1 1 1 0 0 0 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 To understand this, one has to remember that when we get the scalar or vector potential, we know the solution only up to system_gridnp. When we take the final divergence or curl, also the points from system_gridnp+1 to system_gridnp_part matter (due to the way in which the finite differences are computed). Since potentials are set to zero in this non-physical region, there might be spikes at the border. Therefore, we set to zero the points of the fields in the aforementioned mask, so that we do not have spikes

 SURFACE CORRECTION

In the equations above, there are two integrals, the volume one (always computed) and the surface one. The surface integral is important when the total field is significantly different from zero at the boundaries of the box

Data Types

type  helmholtz_decomposition_t
 

Functions/Subroutines

subroutine helmholtz_decomposition_init (this, namespace, sys_grid, system_mc, space)
 Initialize Helmholtz decomposition object. More...
 
subroutine get_indices_from_mask (np, mask, indices)
 
subroutine helmholtz_finalize (this)
 
subroutine helmholtz_visualize_boxes (this, namespace, space)
 Visualise boxes for use in Helmholtz Decomposition. More...
 
subroutine zapply_inner_stencil_mask (this, field)
 
subroutine zsubtract_average_from_inner_stencil (this, field)
 
subroutine zget_vector_potential (this, namespace, vector_potential, total_field, trans_field, apply_boundary)
 
subroutine zget_trans_field (this, namespace, transverse_field, total_field, vector_potential, apply_boundary)
 
subroutine zcompute_surface_correction_trans_field (this, namespace, field, surface_correction)
 
subroutine zget_scalar_potential (this, namespace, scalar_potential, total_field, apply_boundary)
 
subroutine zget_long_field (this, namespace, longitudinal_field, total_field, scalar_potential, apply_boundary)
 
subroutine zcompute_surface_correction_long_field (this, namespace, field, surface_correction)
 
subroutine dapply_inner_stencil_mask (this, field)
 
subroutine dsubtract_average_from_inner_stencil (this, field)
 
subroutine dget_vector_potential (this, namespace, vector_potential, total_field, trans_field, apply_boundary)
 
subroutine dget_trans_field (this, namespace, transverse_field, total_field, vector_potential, apply_boundary)
 
subroutine dcompute_surface_correction_trans_field (this, namespace, field, surface_correction)
 
subroutine dget_scalar_potential (this, namespace, scalar_potential, total_field, apply_boundary)
 
subroutine dget_long_field (this, namespace, longitudinal_field, total_field, scalar_potential, apply_boundary)
 
subroutine dcompute_surface_correction_long_field (this, namespace, field, surface_correction)
 

Function/Subroutine Documentation

◆ helmholtz_decomposition_init()

subroutine helmholtz_decomposition_m::helmholtz_decomposition_init ( class(helmholtz_decomposition_t), intent(inout)  this,
type(namespace_t), intent(in)  namespace,
type(grid_t), intent(in), target  sys_grid,
type(multicomm_t), intent(in)  system_mc,
class(space_t), intent(in)  space 
)
private

Initialize Helmholtz decomposition object.

Parameters
[in,out]thisHelmholtz decomposition type to initialize
[in]sys_gridGrid of the system calling Helmholtz
[in]system_mcMultiCommunicator of the system calling Helmholtz
[in]spaceSpace of the system calling Helmholtz

Definition at line 225 of file helmholtz_decomposition.F90.

◆ get_indices_from_mask()

subroutine helmholtz_decomposition_m::get_indices_from_mask ( integer, intent(in)  np,
logical, dimension(:), intent(in)  mask,
integer, dimension(:), intent(out)  indices 
)
private

Definition at line 335 of file helmholtz_decomposition.F90.

◆ helmholtz_finalize()

subroutine helmholtz_decomposition_m::helmholtz_finalize ( type(helmholtz_decomposition_t), intent(inout)  this)
private

Definition at line 354 of file helmholtz_decomposition.F90.

◆ helmholtz_visualize_boxes()

subroutine helmholtz_decomposition_m::helmholtz_visualize_boxes ( type(helmholtz_decomposition_t), intent(inout)  this,
type(namespace_t), intent(in)  namespace,
class(space_t), intent(in)  space 
)
private

Visualise boxes for use in Helmholtz Decomposition.

Definition at line 367 of file helmholtz_decomposition.F90.

◆ zapply_inner_stencil_mask()

subroutine helmholtz_decomposition_m::zapply_inner_stencil_mask ( class(helmholtz_decomposition_t), intent(inout)  this,
complex(real64), dimension(:,:), intent(out)  field 
)
private

Definition at line 480 of file helmholtz_decomposition.F90.

◆ zsubtract_average_from_inner_stencil()

subroutine helmholtz_decomposition_m::zsubtract_average_from_inner_stencil ( class(helmholtz_decomposition_t), intent(inout)  this,
complex(real64), dimension(:,:), intent(inout)  field 
)
private

Definition at line 497 of file helmholtz_decomposition.F90.

◆ zget_vector_potential()

subroutine helmholtz_decomposition_m::zget_vector_potential ( class(helmholtz_decomposition_t), intent(inout)  this,
type(namespace_t), intent(in)  namespace,
complex(real64), dimension(:,:), intent(out), contiguous  vector_potential,
complex(real64), dimension(:,:), intent(inout), optional, contiguous  total_field,
complex(real64), dimension(:,:), intent(inout), optional, contiguous  trans_field,
logical, intent(in), optional  apply_boundary 
)
private
Parameters
[in,out]thisHelmholtz object
[out]vector_potentialVector potential obtained from the total_field. The field will be defined on the grid of the Helmholtz type (thisgrid), UP UNTIL THISGRIDNP vector_potential(1:thisgridnp_part, 1:thisgridboxdim)
[in,out]total_fieldTotal field to decompose. To ensure the correct behaviour this field must be correctly defined on the grid of the system that calls Helmholtz (sys_gr), UP UNTIL SYS_GRNP_PART total_field(1:sys_grnp_part, 1:sys_grboxdim)
[in,out]trans_fieldTotal field to decompose. To ensure the correct behaviour this field must be correctly defined on the grid of the system that calls Helmholtz (sys_gr), UP UNTIL SYS_GRNP_PART total_field(1:sys_grnp_part, 1:sys_grboxdim)
[in]apply_boundaryShould the curl apply boundary conditions?

Definition at line 526 of file helmholtz_decomposition.F90.

◆ zget_trans_field()

subroutine helmholtz_decomposition_m::zget_trans_field ( class(helmholtz_decomposition_t), intent(inout)  this,
type(namespace_t), intent(in)  namespace,
complex(real64), dimension(:,:), intent(out), contiguous  transverse_field,
complex(real64), dimension(:,:), intent(inout), optional, contiguous  total_field,
complex(real64), dimension(:,:), intent(inout), optional, contiguous  vector_potential,
logical, intent(in), optional  apply_boundary 
)
private
Parameters
[in,out]thisHelmholtz object
[out]transverse_fieldTransverse field obtained from the total_field. It will be defined on the inner_grid of the Helmholtz type (thisinner_grid), UP UNTIL INNER_GRNP transverse_field(1:thisinner_gridnp_part, 1:thisinner_gridboxdim)
[in,out]total_fieldTotal field to decompose. To ensure the correct behaviour this field must be correctly defined on the grid of the system that calls Helmholtz (sys_gr), UP UNTIL SYS_GRNP_PART total_field(1:sys_grnp_part, 1:sys_grboxdim)
[in,out]vector_potentialVector potential from which trans field whould be computed
[in]apply_boundaryShould the curl apply boundary conditions?

Definition at line 690 of file helmholtz_decomposition.F90.

◆ zcompute_surface_correction_trans_field()

subroutine helmholtz_decomposition_m::zcompute_surface_correction_trans_field ( class(helmholtz_decomposition_t), intent(in)  this,
type(namespace_t), intent(in)  namespace,
complex(real64), dimension(:,:), intent(in)  field,
complex(real64), dimension(:,:), intent(out)  surface_correction 
)
private
Parameters
[in]thisHelmholtz object

Definition at line 735 of file helmholtz_decomposition.F90.

◆ zget_scalar_potential()

subroutine helmholtz_decomposition_m::zget_scalar_potential ( class(helmholtz_decomposition_t), intent(inout)  this,
type(namespace_t), intent(in)  namespace,
complex(real64), dimension(:), intent(out), contiguous  scalar_potential,
complex(real64), dimension(:,:), intent(inout), contiguous  total_field,
logical, intent(in), optional  apply_boundary 
)
private
Parameters
[in,out]thisHelmholtz object
[out]scalar_potentialVector potential obtained from the total_field. It will be defined on the grid of the Helmholtz type (thisgrid), UP UNTIL THISGRIDNP: scalar_potential(1:thisgridnp_part)
[in,out]total_fieldTotal field to decompose. To ensure the correct behaviour this field must be correctly defined on the grid of the system that calls Helmholtz (sys_gr), UP UNTIL SYS_GRNP_PART: total_field(1:sys_grnp_part, 1:sys_grboxdim)
[in]apply_boundaryShould the divergence apply boundary conditions?

Definition at line 785 of file helmholtz_decomposition.F90.

◆ zget_long_field()

subroutine helmholtz_decomposition_m::zget_long_field ( class(helmholtz_decomposition_t), intent(inout)  this,
type(namespace_t), intent(in)  namespace,
complex(real64), dimension(:,:), intent(out), contiguous  longitudinal_field,
complex(real64), dimension(:,:), intent(inout), optional, contiguous  total_field,
complex(real64), dimension(:), intent(inout), optional, contiguous  scalar_potential,
logical, intent(in), optional  apply_boundary 
)
private
Parameters
[in,out]thisHelmholtz object
[out]longitudinal_fieldTransverse field obtained from the total_field. It will be defined on the inner_grid of the Helmholtz type (thisinner_grid), UP UNTIL INNER_GRNP longitudinal_field(1:thisinner_gridnp_part, 1:thisinner_gridboxdim)
[in,out]total_fieldTotal field to decompose. To ensure the correct behaviour this field must be correctly defined on the grid of the system that calls Helmholtz (sys_gr) UP UNTIL SYS_GRNP_PART total_field(1:sys_grnp_part, 1:sys_grboxdim)
[in,out]scalar_potentialScalar potential from which long field whould be computed
[in]apply_boundaryShould the divergence apply boundary conditions?

Definition at line 833 of file helmholtz_decomposition.F90.

◆ zcompute_surface_correction_long_field()

subroutine helmholtz_decomposition_m::zcompute_surface_correction_long_field ( class(helmholtz_decomposition_t), intent(in)  this,
type(namespace_t), intent(in)  namespace,
complex(real64), dimension(:,:), intent(in)  field,
complex(real64), dimension(:), intent(out)  surface_correction 
)
private
Parameters
[in]thisHelmholtz object
[in]field(1:thissys_gridnp, 1:thissys_gridboxdim)
[out]surface_correction(1:thissys_gridnp)

Definition at line 880 of file helmholtz_decomposition.F90.

◆ dapply_inner_stencil_mask()

subroutine helmholtz_decomposition_m::dapply_inner_stencil_mask ( class(helmholtz_decomposition_t), intent(inout)  this,
real(real64), dimension(:,:), intent(out)  field 
)
private

Definition at line 1004 of file helmholtz_decomposition.F90.

◆ dsubtract_average_from_inner_stencil()

subroutine helmholtz_decomposition_m::dsubtract_average_from_inner_stencil ( class(helmholtz_decomposition_t), intent(inout)  this,
real(real64), dimension(:,:), intent(inout)  field 
)
private

Definition at line 1021 of file helmholtz_decomposition.F90.

◆ dget_vector_potential()

subroutine helmholtz_decomposition_m::dget_vector_potential ( class(helmholtz_decomposition_t), intent(inout)  this,
type(namespace_t), intent(in)  namespace,
real(real64), dimension(:,:), intent(out), contiguous  vector_potential,
real(real64), dimension(:,:), intent(inout), optional, contiguous  total_field,
real(real64), dimension(:,:), intent(inout), optional, contiguous  trans_field,
logical, intent(in), optional  apply_boundary 
)
private
Parameters
[in,out]thisHelmholtz object
[out]vector_potentialVector potential obtained from the total_field. The field will be defined on the grid of the Helmholtz type (thisgrid), UP UNTIL THISGRIDNP vector_potential(1:thisgridnp_part, 1:thisgridboxdim)
[in,out]total_fieldTotal field to decompose. To ensure the correct behaviour this field must be correctly defined on the grid of the system that calls Helmholtz (sys_gr), UP UNTIL SYS_GRNP_PART total_field(1:sys_grnp_part, 1:sys_grboxdim)
[in,out]trans_fieldTotal field to decompose. To ensure the correct behaviour this field must be correctly defined on the grid of the system that calls Helmholtz (sys_gr), UP UNTIL SYS_GRNP_PART total_field(1:sys_grnp_part, 1:sys_grboxdim)
[in]apply_boundaryShould the curl apply boundary conditions?

Definition at line 1050 of file helmholtz_decomposition.F90.

◆ dget_trans_field()

subroutine helmholtz_decomposition_m::dget_trans_field ( class(helmholtz_decomposition_t), intent(inout)  this,
type(namespace_t), intent(in)  namespace,
real(real64), dimension(:,:), intent(out), contiguous  transverse_field,
real(real64), dimension(:,:), intent(inout), optional, contiguous  total_field,
real(real64), dimension(:,:), intent(inout), optional, contiguous  vector_potential,
logical, intent(in), optional  apply_boundary 
)
private
Parameters
[in,out]thisHelmholtz object
[out]transverse_fieldTransverse field obtained from the total_field. It will be defined on the inner_grid of the Helmholtz type (thisinner_grid), UP UNTIL INNER_GRNP transverse_field(1:thisinner_gridnp_part, 1:thisinner_gridboxdim)
[in,out]total_fieldTotal field to decompose. To ensure the correct behaviour this field must be correctly defined on the grid of the system that calls Helmholtz (sys_gr), UP UNTIL SYS_GRNP_PART total_field(1:sys_grnp_part, 1:sys_grboxdim)
[in,out]vector_potentialVector potential from which trans field whould be computed
[in]apply_boundaryShould the curl apply boundary conditions?

Definition at line 1214 of file helmholtz_decomposition.F90.

◆ dcompute_surface_correction_trans_field()

subroutine helmholtz_decomposition_m::dcompute_surface_correction_trans_field ( class(helmholtz_decomposition_t), intent(in)  this,
type(namespace_t), intent(in)  namespace,
real(real64), dimension(:,:), intent(in)  field,
real(real64), dimension(:,:), intent(out)  surface_correction 
)
private
Parameters
[in]thisHelmholtz object

Definition at line 1259 of file helmholtz_decomposition.F90.

◆ dget_scalar_potential()

subroutine helmholtz_decomposition_m::dget_scalar_potential ( class(helmholtz_decomposition_t), intent(inout)  this,
type(namespace_t), intent(in)  namespace,
real(real64), dimension(:), intent(out), contiguous  scalar_potential,
real(real64), dimension(:,:), intent(inout), contiguous  total_field,
logical, intent(in), optional  apply_boundary 
)
private
Parameters
[in,out]thisHelmholtz object
[out]scalar_potentialVector potential obtained from the total_field. It will be defined on the grid of the Helmholtz type (thisgrid), UP UNTIL THISGRIDNP: scalar_potential(1:thisgridnp_part)
[in,out]total_fieldTotal field to decompose. To ensure the correct behaviour this field must be correctly defined on the grid of the system that calls Helmholtz (sys_gr), UP UNTIL SYS_GRNP_PART: total_field(1:sys_grnp_part, 1:sys_grboxdim)
[in]apply_boundaryShould the divergence apply boundary conditions?

Definition at line 1309 of file helmholtz_decomposition.F90.

◆ dget_long_field()

subroutine helmholtz_decomposition_m::dget_long_field ( class(helmholtz_decomposition_t), intent(inout)  this,
type(namespace_t), intent(in)  namespace,
real(real64), dimension(:,:), intent(out), contiguous  longitudinal_field,
real(real64), dimension(:,:), intent(inout), optional, contiguous  total_field,
real(real64), dimension(:), intent(inout), optional, contiguous  scalar_potential,
logical, intent(in), optional  apply_boundary 
)
private
Parameters
[in,out]thisHelmholtz object
[out]longitudinal_fieldTransverse field obtained from the total_field. It will be defined on the inner_grid of the Helmholtz type (thisinner_grid), UP UNTIL INNER_GRNP longitudinal_field(1:thisinner_gridnp_part, 1:thisinner_gridboxdim)
[in,out]total_fieldTotal field to decompose. To ensure the correct behaviour this field must be correctly defined on the grid of the system that calls Helmholtz (sys_gr) UP UNTIL SYS_GRNP_PART total_field(1:sys_grnp_part, 1:sys_grboxdim)
[in,out]scalar_potentialScalar potential from which long field whould be computed
[in]apply_boundaryShould the divergence apply boundary conditions?

Definition at line 1357 of file helmholtz_decomposition.F90.

◆ dcompute_surface_correction_long_field()

subroutine helmholtz_decomposition_m::dcompute_surface_correction_long_field ( class(helmholtz_decomposition_t), intent(in)  this,
type(namespace_t), intent(in)  namespace,
real(real64), dimension(:,:), intent(in)  field,
real(real64), dimension(:), intent(out)  surface_correction 
)
private
Parameters
[in]thisHelmholtz object
[in]field(1:thissys_gridnp, 1:thissys_gridboxdim)
[out]surface_correction(1:thissys_gridnp)

Definition at line 1404 of file helmholtz_decomposition.F90.