Octopus
field_transfer_oct_m Module Reference

This module implements the field transfer. More...

Detailed Description

This module implements the field transfer.

Idea: field_transfer_t is an abstract interaction that offers regridding and interpolation and can be subclassed for concrete implementations (currently current_to_mxll_field_oct_m::current_to_mxll_field_t and mxll_e_field_to_matter_oct_m::mxll_e_field_to_matter_t).

How to get quantities from a field transfer interaction using the interpolation

In order to get correctly interpolated values also when the systems have different timesteps, at the moment, one always needs to call the interpolate function to get the field at the requested time. As an example from the Maxwell system, this is how to get the current at a certain time:

! interpolate current from interaction
call iter%start(this%interactions)
do while (iter%has_next())
select type (interaction => iter%get_next())
class is (current_to_mxll_field_t)
call interaction%interpolate(time, current_density_ext)
call lalg_axpy(this%gr%np, 3, M_ONE, current_density_ext * this%hm%current_factor, current)
end select
end do

One can loop over all interactions, call the interpolate function and add the result to get the total current in this case.

Time interpolation

The time interpolation class (time_interpolation_t) offers a interpolation in time for each point of a field for an arbitrary depth. It uses the math_oct_m::interpolate function from the math_oct_m module that is based on Lagrange polynomial interpolation. The time interpolation class implements also functionality for reading and writing restart data which is needed to be able to restart consistently with the same history of past values.

Derived interactions

Data Types

type  field_transfer_t
 class defining the field_transfer interaction More...
 

Functions/Subroutines

subroutine, public field_transfer_init (this, gr, ndim)
 the system field is allocated and initialized to 0 More...
 
subroutine field_transfer_init_from_partner (this, partner_gr, partner_space, partner_namespace)
 the partner field is allocated and initialized to 0; moreover the regridding structure is initialized More...
 
subroutine field_transfer_init_interpolation (this, depth, label, cmplx)
 the time interpolation is initialized; it needs to know the depth which is usually given by the order of the propagator; thus it can only be called after the propagator is known (this can be done in the system in init_algorithm) More...
 
subroutine field_transfer_end (this)
 
subroutine field_transfer_calculate (this)
 
subroutine field_transfer_do_mapping (this)
 perform the regridding and add the system field to the time interpolator using the time of the quantity at this point; More...
 
subroutine dfield_transfer_interpolate (this, time, field)
 return the interpolated field for a given time More...
 
subroutine zfield_transfer_interpolate (this, time, field)
 return the interpolated field for a given time More...
 
subroutine field_transfer_calculate_energy (this)
 
subroutine field_transfer_read_restart (this, mesh, space, restart, err)
 
subroutine field_transfer_write_restart (this, mesh, space, restart, err)
 

Function/Subroutine Documentation

◆ field_transfer_init()

subroutine, public field_transfer_oct_m::field_transfer_init ( class(field_transfer_t), intent(inout)  this,
type(grid_t), intent(in), target  gr,
integer, intent(in)  ndim 
)

the system field is allocated and initialized to 0

Definition at line 206 of file field_transfer.F90.

◆ field_transfer_init_from_partner()

subroutine field_transfer_oct_m::field_transfer_init_from_partner ( class(field_transfer_t), intent(inout)  this,
type(grid_t), intent(in)  partner_gr,
class(space_t), intent(in)  partner_space,
type(namespace_t), intent(in)  partner_namespace 
)
private

the partner field is allocated and initialized to 0; moreover the regridding structure is initialized

Definition at line 222 of file field_transfer.F90.

◆ field_transfer_init_interpolation()

subroutine field_transfer_oct_m::field_transfer_init_interpolation ( class(field_transfer_t), intent(inout)  this,
integer, intent(in)  depth,
character(len=*), intent(in)  label,
logical, intent(in), optional  cmplx 
)
private

the time interpolation is initialized; it needs to know the depth which is usually given by the order of the propagator; thus it can only be called after the propagator is known (this can be done in the system in init_algorithm)

Definition at line 240 of file field_transfer.F90.

◆ field_transfer_end()

subroutine field_transfer_oct_m::field_transfer_end ( class(field_transfer_t), intent(inout)  this)
private

Definition at line 255 of file field_transfer.F90.

◆ field_transfer_calculate()

subroutine field_transfer_oct_m::field_transfer_calculate ( class(field_transfer_t), intent(inout)  this)
private

Definition at line 269 of file field_transfer.F90.

◆ field_transfer_do_mapping()

subroutine field_transfer_oct_m::field_transfer_do_mapping ( class(field_transfer_t), intent(inout)  this)
private

perform the regridding and add the system field to the time interpolator using the time of the quantity at this point;

it is called by the partner system when copying its data to the interaction because the system grid does not change Note: this needs an extra function call for each of the partner systems in their copy_quantities_to_interaction function

Definition at line 290 of file field_transfer.F90.

◆ dfield_transfer_interpolate()

subroutine field_transfer_oct_m::dfield_transfer_interpolate ( class(field_transfer_t), intent(in)  this,
real(real64), intent(in)  time,
real(real64), dimension(:, :), intent(out), contiguous  field 
)
private

return the interpolated field for a given time

Definition at line 309 of file field_transfer.F90.

◆ zfield_transfer_interpolate()

subroutine field_transfer_oct_m::zfield_transfer_interpolate ( class(field_transfer_t), intent(in)  this,
real(real64), intent(in)  time,
complex(real64), dimension(:, :), intent(out), contiguous  field 
)
private

return the interpolated field for a given time

Definition at line 323 of file field_transfer.F90.

◆ field_transfer_calculate_energy()

subroutine field_transfer_oct_m::field_transfer_calculate_energy ( class(field_transfer_t), intent(inout)  this)
private

Definition at line 336 of file field_transfer.F90.

◆ field_transfer_read_restart()

subroutine field_transfer_oct_m::field_transfer_read_restart ( class(field_transfer_t), intent(inout)  this,
class(mesh_t), intent(in)  mesh,
class(space_t), intent(in)  space,
type(restart_t), intent(in)  restart,
integer, intent(out)  err 
)
private

Definition at line 348 of file field_transfer.F90.

◆ field_transfer_write_restart()

subroutine field_transfer_oct_m::field_transfer_write_restart ( class(field_transfer_t), intent(inout)  this,
class(mesh_t), intent(in)  mesh,
class(space_t), intent(in)  space,
type(restart_t), intent(in)  restart,
integer, intent(out)  err 
)
private

Definition at line 362 of file field_transfer.F90.