![]() |
RFOF - RF Library for Orbit Following Codes
|
Controls the workflow in RFOF. More...
Public Member Functions | |
subroutine | RFOF_master (time, time_previous_step, RFglobal, marker, resonanceMemory, diagno, MPI_nod_Id, interaction_failed_particle_overshot_resonance, time_at_resonance) |
Main routine to be called by orbit code after every time step. Note that this routine give an error flag requiring the orbit code to retake timestep. | |
subroutine | split_RFstep_in_single_mode (marker, resonanceMemory, RFglobal, nr_resonant_waves, list_resonant_waves, diagno, MPI_nod_Id) |
This routine separates each wave-particle resonance and calls the Monte Carlo routines give the marker a Monte Carlo kick. | |
subroutine | RFOF_destructor (wave, mem, diagno) |
Deallocates memory allocated by RFOF. |
Controls the workflow in RFOF.
Part of RFOF - RF operator for Orbit Following codes. Developed during 2010 by the EU-ITM, Integrated modelling project 5 under ACT4.
Public routine:
Overview of workflow: 1. After each orbit step the RF module is called to check: a) if particle if is resonating with one or more waves (near enough to resonance); a list is made of all waves which the particle is resonating with. b) If the particle has crossed at least on resonance without recieving a kick an error flag is returned together with an approximate time for the earliest missed resonance. This allows the orbit code to refine the time step and a approximate time when the resonance is expected. 4. RF operator now prepares for a Monte Carlo kick by calculating: a) the quasilinear diffucion coefficient b) the coefficients of the Monte Carlo stepping scheme 5. Monte Carlo step is performed 6. Check if step is correct; if not GOTO 4.
Definition at line 65 of file RFOF_main.F90.
subroutine RFOF_main::RFOF_master | ( | real(8), intent(in) | time, |
real(8), intent(in) | time_previous_step, | ||
type(rf_wave_global), intent(inout) | RFglobal, | ||
type(particle), intent(inout) | marker, | ||
type(resonance_memory), dimension(:,:), intent(inout) | resonanceMemory, | ||
type(RFOF_cumlative_diagnostics), intent(inout) | diagno, | ||
integer, intent(in) | MPI_nod_Id, | ||
logical, intent(out) | interaction_failed_particle_overshot_resonance, | ||
real(8), intent(out) | time_at_resonance | ||
) |
Main routine to be called by orbit code after every time step. Note that this routine give an error flag requiring the orbit code to retake timestep.
This routine controls the workflow of the RFOF library. The workflow consists of two main blocks:
Input:
time | Time in simulation. Use to provide input to the time-history stored in the "resonanceMemory" and used to evaluate time derivaties. |
time_previous_step | Time in simulation at the previous timestep for this marker. Make sure RFOF ignores resonance from before the previous timestep. |
MPI_nod_Id | Number used to identify MPI nods during parallel execusion. |
Input & Output:
marker | datatype defined in the module RFOF_markers. Includes all needed information about the particle at "present" time. |
resonanceMemory | datatype defined in the module RFOF_resonance_memory. Includes time histories of marker parameter of which we need time derivaties. |
RFglobal | datatype defined in the module RFOF_waves. Includes all needed information about the RF wave fields at "present" time. |
diagno | Structure containing time integrated diagnostics. |
Output:
interaction_failed_particle_overshot_resonance | Logial variables which states if the marker has crossed any resonance without recieving any acceleration, i.e. is has "overshot" the resonance. |
time_at_resonance | Estimate of the time at the next resonance. In case the marker overshot the resonance this time is in the past, othervise it is in the future. |
Definition at line 151 of file RFOF_main.F90.
subroutine RFOF_main::split_RFstep_in_single_mode | ( | type(particle), intent(inout) | marker, |
type(resonance_memory), dimension(:,:), intent(in) | resonanceMemory, | ||
type(rf_wave_global), intent(inout) | RFglobal, | ||
integer, intent(in) | nr_resonant_waves, | ||
real(8), dimension(:,:), intent(in), pointer | list_resonant_waves, | ||
type(RFOF_cumlative_diagnostics), intent(inout) | diagno, | ||
integer, intent(in) | MPI_nod_Id | ||
) |
This routine separates each wave-particle resonance and calls the Monte Carlo routines give the marker a Monte Carlo kick.
INPUT:
resonanceMemory | (in) Short term memory of the resonance condition along the orbit |
nr_resonant_waves | (in) Number of wave modes with which the particle is in resonance. |
list_resonant_waves | (in) List including pointers to the wave modes which the particle is resonance with. |
MPI_nod_Id | (in) Number used to identify MPI nods during parallel execusion. |
INPUT & OUTPUT:
marker | (inout) Properties of the marker |
RFglobal | (inout) Global RF wave field; describing the wave in the throughout plasma. |
diagno | (inout) Time cumulative diagnostics (to be updated as marker is accelerated by the RF wave field) |
Definition at line 248 of file RFOF_main.F90.
subroutine RFOF_main::RFOF_destructor | ( | type(rf_wave_global), intent(inout) | wave, |
type(resonance_memory), dimension(:,:), intent(out), pointer | mem, | ||
type(RFOF_cumlative_diagnostics), intent(inout) | diagno | ||
) |
Deallocates memory allocated by RFOF.
wave | (inout) Global RF wave field structure |
mem | (inout) Resonance memory structure |
diagno | (inout) Cumulative diagnostics structure |
Definition at line 301 of file RFOF_main.F90.