Subroutine that loads the simulation parameters from the file specified in params\%path_to_inputs
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(KORC_PARAMS), | intent(inout) | :: | params | Core KORC simulation parameters. |
subroutine load_korc_params(params)
!! @note Subroutine that loads the simulation parameters from the
!! file specified in params\%path_to_inputs @endnote
TYPE (KORC_PARAMS), INTENT(INOUT) :: params
!! Core KORC simulation parameters.
LOGICAL :: restart
!! Flag to indicate if the simulations restarts (restart=T) or not
!! (restart=F). Restart simulation that exited before simulation_time
!! reached.
LOGICAL :: proceed
!! Flag to indicate if the simulations proceeds (proceed=T) or not
!! (proceed=F). Append simulation results after previous simulation_time
!! reached.
LOGICAL :: reinit
!! Flag to begin a new simulation, reinitializing from restart file state
REAL(rp) :: simulation_time
!! Total simulation time in seconds.
REAL(rp) :: snapshot_frequency
!! Time between snapshots in time of the simulation.
REAL(rp) :: restart_overwrite_frequency
!! Time between overwrites of restart file in time of the simulation.
REAL(rp) :: dt
!! Time step in the simulation as a fraction of the relativistic
!! electron gyro-period @f$\tau_e = 2\pi\gamma m_e/eB_0@f$
REAL(rp) :: minimum_particle_energy
!! Minimum allowed relativistic factor @f$\gamma@f$ of simulated electrons.
LOGICAL :: radiation
!! Flag to indicate if synchrotron radiation losses are included
!! (radiation=T) or not (radiation=F).
LOGICAL :: collisions
!! Flag to indicate if collisionsare included (collisions=T) or not
!! (collisions=F).
CHARACTER(MAX_STRING_LENGTH) :: GC_rad_model
CHARACTER(MAX_STRING_LENGTH) :: collisions_model
!! String with the name of the collisions model to be used in the simulation.
CHARACTER(MAX_STRING_LENGTH) :: bound_electron_model
CHARACTER(MAX_STRING_LENGTH) :: profile_model
!! String with the name of the model for the plasma profiles.
CHARACTER(MAX_STRING_LENGTH) :: field_model
!! String with the name of the model for the field profiles.
CHARACTER(MAX_STRING_LENGTH) :: magnetic_field_filename
!! String with the name of the model for the fields and plasma profiles.
CHARACTER(MAX_STRING_LENGTH) :: outputs_list
!! List of electron variables to include in the outputs.
INTEGER :: num_species
!! Number of different populations of simulated relativistic electrons
!! in KORC.
INTEGER :: imax
!! Auxiliary variable used to parse the output_list
INTEGER :: imin
!! Auxiliary variable used to parse the output_list
INTEGER :: ii
!! Iterator used to parse the output_list
INTEGER :: jj
!! Iterator used to parse the output_list
INTEGER :: num_outputs
!! Auxiliary variable used to parse the output_list
INTEGER, DIMENSION(2) :: indices
!! Auxiliary variable used to parse the output_list
LOGICAL :: HDF5_error_handling
!! Flag for HDF5 error handling
LOGICAL :: FO_GC_compare
CHARACTER(MAX_STRING_LENGTH) :: orbit_model
!! String with the name of the orbit model ('FO' or 'GC').
CHARACTER(MAX_STRING_LENGTH) :: field_eval
!! String with the name of the field evaluation method for
!! analytical fields ('interp' or 'eqn')
LOGICAL :: FokPlan
!! Flag to decouple spatial-dependence of evolution
LOGICAL :: SameRandSeed
LOGICAL :: SC_E
LOGICAL :: SC_E_add
NAMELIST /input_parameters/ restart,field_model,magnetic_field_filename, &
simulation_time,snapshot_frequency,dt,num_species,radiation, &
collisions,collisions_model,outputs_list,minimum_particle_energy, &
HDF5_error_handling,orbit_model,field_eval,proceed,profile_model, &
restart_overwrite_frequency,FokPlan,GC_rad_model,bound_electron_model, &
FO_GC_compare,SameRandSeed,SC_E,reinit,SC_E_add
open(unit=default_unit_open,file=TRIM(params%path_to_inputs), &
status='OLD',form='formatted')
read(default_unit_open,nml=input_parameters)
close(default_unit_open)
params%restart = restart
params%proceed = proceed
params%reinit = reinit
params%simulation_time = simulation_time
params%snapshot_frequency = snapshot_frequency
params%restart_overwrite_frequency=restart_overwrite_frequency
params%dt = dt
params%num_species = num_species
params%profile_model = TRIM(profile_model)
params%field_model = TRIM(field_model)
params%magnetic_field_filename = TRIM(magnetic_field_filename)
params%minimum_particle_energy = minimum_particle_energy*C_E
params%minimum_particle_g = 1.0_rp + params%minimum_particle_energy/ &
(C_ME*C_C**2) ! Minimum value of relativistic gamma factor
params%radiation = radiation
params%collisions = collisions
params%collisions_model = TRIM(collisions_model)
params%bound_electron_model = TRIM(bound_electron_model)
params%GC_rad_model = TRIM(GC_rad_model)
if (HDF5_error_handling) then
params%HDF5_error_handling = 1_idef
else
params%HDF5_error_handling = 0_idef
end if
params%orbit_model = orbit_model
params%FO_GC_compare = FO_GC_compare
params%field_eval = field_eval
params%GC_coords=.FALSE.
params%FokPlan=FokPlan
params%SameRandSeed = SameRandSeed
params%SC_E=SC_E
params%SC_E_add=SC_E_add
! Loading list of output parameters (parsing)
imin = SCAN(outputs_list,'{')
imax = SCAN(outputs_list,'}')
ii = 1_idef
jj = 1_idef
num_outputs = 1_idef
do while (ii.NE.0)
ii = SCAN(outputs_list(jj:),",")
if (ii.NE.0) then
jj = jj + ii
num_outputs = num_outputs + 1_idef
end if
end do
ALLOCATE(params%outputs_list(num_outputs))
if (num_outputs.GT.1_idef) then
indices = 0_idef
indices(2) = SCAN(outputs_list,",")
params%outputs_list(1) = TRIM(outputs_list(imin+1_idef:indices(2)-1_idef))
indices(1) = indices(1) + indices(2) + 1_idef
do ii=2_idef,num_outputs
indices(2) = SCAN(outputs_list(indices(1):),",")
if (indices(2).EQ.0_idef) then
params%outputs_list(ii) = TRIM(outputs_list(indices(1):imax-1_idef))
else
params%outputs_list(ii) = TRIM(outputs_list(indices(1):indices(1)+indices(2)-2_idef))
indices(1) = indices(1) + indices(2)
end if
end do
else
params%outputs_list(1) = TRIM(outputs_list(imin+1_idef:imax-1_idef))
end if
if (params%mpi_params%rank .EQ. 0) then
write(6,'(/,"* * * * * SIMULATION PARAMETERS * * * * *")')
write(6,'("Restarting simulation: ",L1)') params%restart
write(6,'("Continuing simulation: ",L1)') params%proceed
write(6,'("Number of electron populations: ",I16)') params%num_species
write(6,'("Orbit model: ",A50)') TRIM(params%orbit_model)
write(6,'("Magnetic field model: ",A50)') TRIM(params%field_model)
write(6,'("Magnetic field evaluation: ",A50)') TRIM(params%field_eval)
if (TRIM(params%field_model).EQ.'EXTERNAL') then
write(6,'("Magnetic field file: ",A100)') TRIM(params%magnetic_field_filename)
end if
write(6,'("Radiation losses included: ",L1)') params%radiation
if (params%radiation.and.(params%orbit_model(1:2).eq.'GC')) then
write(6,'("Radiation model: ",A50)') TRIM(params%GC_rad_model)
end if
write(6,'("Collisions losses included: ",L1)') params%collisions
if (params%collisions) then
write(6,'("Collisions model: ",A50)') TRIM(params%collisions_model)
write(6,'("Bound electron model: ",A50)') &
TRIM(params%bound_electron_model)
end if
write(6,'("Self-consistent E included: ",L1)') params%SC_E
write(6,'("* * * * * * * * * * * * * * * * * * * * *",/)')
end if
end subroutine load_korc_params