RFOF - RF Library for Orbit Following Codes
RFOF_kick.F90 File Reference
#include "config.h"
+ Include dependency graph for RFOF_kick.F90:

Go to the source code of this file.

Data Types

module  RFOF_kick
 Module handling wave acceleartion physics, i.e. the Monte Carlo kicks. More...
type  RFOF_kick::coeff_for_guiding_centre_kick
 Coefficients defining the direction of the RF acceleration (for a single mode) in the configuration space $(\mu, \, v_\|, \, \psi, \, \theta)$. More...

Functions/Subroutines

subroutine wrap_RF_diffusion_and_move_marker (Iperp, state, diffusion, errorFlag)
 Wrapper for the subroutines.
real(8) function quasilinear_RF_diffusion_coeff (marker, RFlocal, mem)
 Quasilinear diffusion coefficient for resonant interactions between RF wave field and charged particle in an inhomogeous magnetic field. This is a local operator given by the local magnetic fields and wave quantities and their time derivatives in the frame of the particles.
real(8) function max_RFkick_in_single_pass (marker, wave, mem)
 Maximum possible kick in $ v_\perp $ during one crossing of the resonance.
real(8) function tau_RF_phase_integral (time, omega_res, Nelements)
 Phase integral describing the effective time tau which the particle is in resonance with the wave. [T. Johnson et al 2006 Nucl. Fusion 46 S433, http://iopscience.iop.org/0029-5515/46/7/S06].
type(coeff_for_guiding_centre_kick)
function 
coeff_RF_characteristic (marker, Blocal, RFlocal, nharm)
 Monte Carlo time stepping scheme for quasilinear resonant interactions between RF wave field and charged particle in an inhomogeous magnetic field. This is a local operator given by the local magnetic fields and wave quantities and their time derivatives in the frame of the particles.
subroutine wrap_move_marker_on_characteristic (Iperp, state, errorFlag)
 Move the marker along the RF characteristic, coeff, to $I_\perp$=Iperp.
subroutine move_marker_on_characteristic (marker, dIperp, Blocal, coeff)
 Move the marker a distance dIperp along the RF characteristic coeff.
subroutine get_dvpar_on_characteristic (marker, Blocal, RFlocal, dvperp2, dvpar)
real(8) function timestep_factor (x_a, x_c, state_a)
 compares the one step relative error with prescribed values and return a factor<1 , if the eror is too large, >1 if the error is close to machine precision, and exactly 1 if OK
logical function close_to_boundary (x, state)
 return .TRUE. of the distance to the boundary is less then prescribed value Neumann boundary conditions

Function/Subroutine Documentation

subroutine quasilinear_RF_kick_steinbrecher_integrator::wrap_RF_diffusion_and_move_marker ( real(8), intent(in)  Iperp,
type(RFOF_state), intent(inout)  state,
real(8), intent(out)  diffusion,
integer, intent(out)  errorFlag 
)

Wrapper for the subroutines.

  • quasilinear_RF_diffusion_coeff
  • wrap_move_marker_on_characteristic Returns the quasilinear RF diffusion coefficient. Before evaluating this diffusion coefficient, the state should be updated from the input state in the following way: starting from the input state, establish the characteristic lines (in phase space) of the RF wave-particle interactions; then move the marker (statemarker) from the present position, along the characteristic to the point where Iperp = stateIperp. Once the state has been updated (using wrap_move_marker_on_characteristic) the diffusion coefficient can be evaluated using quasilinear_RF_diffusion_coeff.
Todo:
Optimize, the present version is a quick wrapper including already wrapped routines.

INPUT

Parameters:
IperpNew value of $I_\perp$ to which the state should be moved

INPUT/OUTPUT

Parameters:
stateThe state of a collection of RFOF structures.

OUTPUT

Parameters:
diffusion(out) Quasilinear scattering coefficient : $ D = <dI_\perp dI_\perp> $ during one crossing of the resonance.
errorFlagFlag to be raised if the routine failed. errorFlag=0 marker successfully moved errorFlag>0 failed to move marker errorFlag=1 new state does not exists errorFlag=2 step too long, could not perform move
[in]IperpMagnetic moment normalised for RF
[in,out]stateState of RFOF (Note: may not be latest state - in case stateIperp is different from Iperp, then state needed updating)
[out]errorFlagError flag
[out]diffusiondiffusion coefficient

Definition at line 244 of file RFOF_kick.F90.

References quasilinear_RF_diffusion_coeff(), and wrap_move_marker_on_characteristic().

Referenced by RFOF_kick::quasilinear_RF_kick_steinbrecher_integrator().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

real(8) function quasilinear_RF_kick_steinbrecher_integrator::quasilinear_RF_diffusion_coeff ( type(particle), intent(in)  marker,
type(rf_wave_local), intent(in)  RFlocal,
type(resonance_memory), intent(in)  mem 
)

Quasilinear diffusion coefficient for resonant interactions between RF wave field and charged particle in an inhomogeous magnetic field. This is a local operator given by the local magnetic fields and wave quantities and their time derivatives in the frame of the particles.

The diffusion coefficient in $ I_\perp $ is:

\[ D = \frac{1}{2} \left| \tau Ze v_\perp E_{eff} \right|^2 \]

See [ L.-G. Eriksson, Nuclear Fusion, 2005 ]

Attention:
dIperp is in units [MeV].

INPUT

Parameters:
marker(in) Properties of the marker (weight, mass, charge, velocity...)
RFlocal(in) Local properties of the RF wave field ( $|B|, RB_\phi,\dots$)
mem(in) Short term memory of the resonance condition along the orbit

OUTPUT

Returns:
diffusion (out) Quasilinear scattering coefficient : $ D = <dI_\perp dI_\perp> $ during one crossing of the resonance.

Definition at line 285 of file RFOF_kick.F90.

References max_RFkick_in_single_pass().

Referenced by wrap_RF_diffusion_and_move_marker().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

real(8) function quasilinear_RF_kick_steinbrecher_integrator::max_RFkick_in_single_pass ( type(particle), intent(in)  marker,
type(rf_wave_local), intent(in)  wave,
type(resonance_memory), intent(in)  mem 
)

Maximum possible kick in $ v_\perp $ during one crossing of the resonance.

INPUT

Parameters:
markerProperties of the marker
waveLocal properties of the RF wave field
memShort term memory of the resonance condition along the orbit

OUTPUT

Returns:
max_kick Maximum possible kick in $ v_\perp $ during one crossing of the resonance

Definition at line 317 of file RFOF_kick.F90.

References tau_RF_phase_integral().

Referenced by quasilinear_RF_diffusion_coeff().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

real(8) function quasilinear_RF_kick_steinbrecher_integrator::tau_RF_phase_integral ( real(8), dimension(nelements), intent(in)  time,
real(8), dimension(nelements), intent(in)  omega_res,
integer, intent(in)  Nelements 
)

Phase integral describing the effective time tau which the particle is in resonance with the wave. [T. Johnson et al 2006 Nucl. Fusion 46 S433, http://iopscience.iop.org/0029-5515/46/7/S06].

INPUT

Parameters:
timeTime vector (used for time-derivatives of the resonance function)
omega_resVector of the resonance function $ n\Omega+\mathbf{k}\cdot\mathbf{v}-\omega $
NelementsNumber of elements in the vectors time and omega_res

OUTPUT

Returns:
tau The effective time tau which the particle is in resonance with the wave; given by the phase integral.

Definition at line 377 of file RFOF_kick.F90.

Referenced by max_RFkick_in_single_pass().

+ Here is the caller graph for this function:

type(coeff_for_guiding_centre_kick) function quasilinear_RF_kick_steinbrecher_integrator::coeff_RF_characteristic ( type(particle), intent(in)  marker,
type(magnetic_field_local), intent(in)  Blocal,
type(rf_wave_local), intent(in)  RFlocal,
integer, intent(in)  nharm 
)

Monte Carlo time stepping scheme for quasilinear resonant interactions between RF wave field and charged particle in an inhomogeous magnetic field. This is a local operator given by the local magnetic fields and wave quantities and their time derivatives in the frame of the particles.

Attention:
dIperp is in units [MeV].

INPUT

Parameters:
markerProperties of the marker
BlocalLocal properties of the magnetic field
RFlocalLocal properties of the RF wave field
nharmHarmonic number of the resonance

OUTPUT

Returns:
coeff Coefficients describing the characteristic of RF acceleration

Definition at line 440 of file RFOF_kick.F90.

Referenced by wrap_move_marker_on_characteristic().

+ Here is the caller graph for this function:

subroutine quasilinear_RF_kick_steinbrecher_integrator::wrap_move_marker_on_characteristic ( real(8), intent(in)  Iperp,
type(RFOF_state), intent(inout)  state,
integer, intent(out)  errorFlag 
)

Move the marker along the RF characteristic, coeff, to $I_\perp$=Iperp.

Attention:
Iperp is in units [MeV].
Todo:
Make this into routine which is not wrapped

INPUT

Parameters:
IperpNew value of $I_\perp$ to which the state should be moved

INPUT/OUTPUT

Parameters:
stateThe state of a collection of RFOF structures.

OUTPUT

Parameters:
errorFlagFlag to be raised if the routine failed. errorFlag=0 marker successfully moved errorFlag>0 failed to move marker errorFlag=1 new state does not exists errorFlag=2 step too long, could not perform move

Definition at line 505 of file RFOF_kick.F90.

References coeff_RF_characteristic(), and move_marker_on_characteristic().

Referenced by RFOF_kick::quasilinear_RF_kick_steinbrecher_integrator(), and wrap_RF_diffusion_and_move_marker().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

subroutine quasilinear_RF_kick_steinbrecher_integrator::move_marker_on_characteristic ( type(particle), intent(inout)  marker,
real(8)  dIperp,
type(magnetic_field_local), intent(in)  Blocal,
type(coeff_for_guiding_centre_kick), intent(in)  coeff 
)

Move the marker a distance dIperp along the RF characteristic coeff.

Attention:
dIperp is in units [MeV].
Todo:
Change the kick in vpar to stay exactly on characteristics.

INPUT

Parameters:
dIperpPerturbation of $I_\perp$
BlocalLocal properties of the magnetic field
coeffCoefficients describing the characteristic of RF acceleration

INPUT/OUTPUT

Parameters:
markerProperties of the marker

Definition at line 571 of file RFOF_kick.F90.

Referenced by wrap_move_marker_on_characteristic().

+ Here is the caller graph for this function:

subroutine quasilinear_RF_kick_steinbrecher_integrator::get_dvpar_on_characteristic ( type(particle), intent(in)  marker,
type(magnetic_field_local), intent(in)  Blocal,
type(rf_wave_local), intent(in)  RFlocal,
real(8), intent(in)  dvperp2,
real(8), intent(out)  dvpar 
)

Definition at line 621 of file RFOF_kick.F90.

real(8) function quasilinear_RF_kick_steinbrecher_integrator::timestep_factor ( real(8), intent(in)  x_a,
real(8), intent(in)  x_c,
type(RFOF_state), intent(in)  state_a 
)

compares the one step relative error with prescribed values and return a factor<1 , if the eror is too large, >1 if the error is close to machine precision, and exactly 1 if OK

Definition at line 667 of file RFOF_kick.F90.

Referenced by MC_kick_steinbrecher_integrator().

+ Here is the caller graph for this function:

logical function quasilinear_RF_kick_steinbrecher_integrator::close_to_boundary ( real(8), intent(in)  x,
type(RFOF_state), intent(in)  state 
)

return .TRUE. of the distance to the boundary is less then prescribed value Neumann boundary conditions

Definition at line 682 of file RFOF_kick.F90.

Referenced by MC_kick_steinbrecher_integrator().

+ Here is the caller graph for this function: