connection_length_m Module

Connection length and shortest distance to the target routines


Uses

  • module~~connection_length_m~~UsesGraph module~connection_length_m connection_length_m module~descriptors_m descriptors_m module~connection_length_m->module~descriptors_m module~equilibrium_m equilibrium_m module~connection_length_m->module~equilibrium_m module~error_handling_m error_handling_m module~connection_length_m->module~error_handling_m module~fieldline_tracer_m fieldline_tracer_m module~connection_length_m->module~fieldline_tracer_m module~precision_m precision_m module~connection_length_m->module~precision_m module~screen_io_m screen_io_m module~connection_length_m->module~screen_io_m module~status_codes_m status_codes_m module~connection_length_m->module~status_codes_m module~descriptors_m->module~screen_io_m module~equilibrium_m->module~precision_m module~error_handling_m->module~precision_m module~error_handling_m->module~screen_io_m module~error_handling_m->module~status_codes_m module~comm_handling_m comm_handling_m module~error_handling_m->module~comm_handling_m mpi mpi module~error_handling_m->mpi netcdf netcdf module~error_handling_m->netcdf module~fieldline_tracer_m->module~equilibrium_m module~fieldline_tracer_m->module~error_handling_m module~fieldline_tracer_m->module~precision_m module~fieldline_tracer_m->module~screen_io_m module~fieldline_tracer_m->module~status_codes_m dop853_constants dop853_constants module~fieldline_tracer_m->dop853_constants dop853_module dop853_module module~fieldline_tracer_m->dop853_module module~fieldline_tracer_m->module~comm_handling_m iso_c_binding iso_c_binding module~precision_m->iso_c_binding iso_fortran_env iso_fortran_env module~precision_m->iso_fortran_env module~precision_m->mpi module~precision_m->netcdf module~screen_io_m->module~precision_m module~screen_io_m->iso_fortran_env module~screen_io_m->netcdf module~comm_handling_m->mpi

Used by

  • module~~connection_length_m~~UsedByGraph module~connection_length_m connection_length_m module~immersed_trace_m immersed_trace_m module~immersed_trace_m->module~connection_length_m module~parbnd_taylor_m parbnd_taylor_m module~parbnd_taylor_m->module~connection_length_m module~immersed_factory_m immersed_factory_m module~immersed_factory_m->module~immersed_trace_m module~parbnd_taylor_netcdf_s parbnd_taylor_netcdf_s module~parbnd_taylor_netcdf_s->module~parbnd_taylor_m

Functions

public function connection_length(equi, x, y, phi, dphi_max, maxstep)

It computes the connection length for the field line defined by x, y, phi (the distance from boundary to boundary). If (x, y, phi) is located inside the wall, the output is negative. If (x, y, phi) is located in the closed flux surface region, the output is the largest allowable number for PARALLAX precision. If the boundary is reached along the field line, but only in one direction, then the routine will return with the arclength to the found boundary point plus the arclength up to dphi_max in the other direction.

Arguments

Type IntentOptional Attributes Name
class(equilibrium_t), intent(inout), target :: equi

Equilibrium defining the field line

real(kind=FP), intent(in) :: x

Spatial coordinate x at which the field line is located

real(kind=FP), intent(in) :: y

Spatial coordinate y at which the field line is located

real(kind=FP), intent(in) :: phi

Spatial coordinate phi at which the field line is located

real(kind=FP), intent(in) :: dphi_max

Maximum toroidal angle to be traced

real(kind=FP), intent(in), optional :: maxstep

Maximum stepsize for the DOP853 integrator

Return Value real(kind=fp)


Subroutines

public subroutine distance_to_boundary(equi, x, y, phi, dphi_max, arclength_fwd, arclength_bwd, target_reached_fwd, target_reached_bwd, maxstepsize, dphi_fwd, dphi_bwd)

Computes the distance to the target in each direction If the (x,y,phi) point is located outside of the vessel, the distances are defined negative

Arguments

Type IntentOptional Attributes Name
class(equilibrium_t), intent(inout), target :: equi

Equilibrium defining the field line

real(kind=FP), intent(in) :: x

Spatial coordinate x at which the field line is located

real(kind=FP), intent(in) :: y

Spatial coordinate y at which the field line is located

real(kind=FP), intent(in) :: phi

Spatial coordinate phi at which the field line is located

real(kind=FP), intent(in) :: dphi_max

Maximum toroidal angle to be traced

real(kind=FP), intent(out) :: arclength_fwd

For points that are inside vessel (based on equi%in_vessel): Arclength in the direction of B until the target is reached For points that are not inside vessel: Negative arclength in the direction of -B until target is reached

real(kind=FP), intent(out) :: arclength_bwd

For points that are inside vessel (based on equi%in_vessel): Arclength in the direction of -B until the target is reached For points that are not inside vessel: Negative arclength in the direction of B until target is reached

logical, intent(inout) :: target_reached_fwd

True, if the target is reached during the forward (as defined above) tracing

logical, intent(inout) :: target_reached_bwd

True, if the target is reached during the backward (as defined above) tracing

real(kind=FP), intent(in) :: maxstepsize

Maximum stepsize for the DOP853 integrator

real(kind=FP), intent(out), optional :: dphi_fwd

Toroidal angle traced in direction of B until target is reached It follows the same sign convention as arclength_fwd

real(kind=FP), intent(out), optional :: dphi_bwd

Toroidal angle traced in direction of -B until target is reached It follows the same sign convention as arclength_bwd

public subroutine shortest_dist_to_bound(equi, x, y, phi, dphi_max, shortest_dist, maxstep, dirind)

It computes the shortest distance to the boundary along B. If (x, y, phi) is located inside the wall, the output is negative. If (x, y, phi) is located in the closed flux surface region, the output is the largest allowable number for PARALLAX precision.

Arguments

Type IntentOptional Attributes Name
class(equilibrium_t), intent(inout), target :: equi

Equilibrium defining the field line

real(kind=FP), intent(in) :: x

Spatial coordinate x at which the field line is located

real(kind=FP), intent(in) :: y

Spatial coordinate y at which the field line is located

real(kind=FP), intent(in) :: phi

Spatial coordinate phi at which the field line is located

real(kind=FP), intent(in) :: dphi_max

Maximum toroidal angle to be traced

real(kind=FP), intent(out) :: shortest_dist

Shortest distance to boundary along the field line from (x,y,phi)

real(kind=FP), intent(in), optional :: maxstep

Maximum stepsize for the DOP853 integrator

integer, intent(out), optional :: dirind

Refers to the direction to the closest boundary along B: = 1, if the closest transition is in the +phi direction = -1, if the closest transition is in the -phi direction = 0, if the it could not found any boundary along the field line.