initialize_profiles Subroutine

public subroutine initialize_profiles(params, P, F)

KORC can run using either analytical and pre-computed plasma profiles. Pre-computed plasma profiles, as in the case of pre-computed electric or magnetic fields, are interpolated to electrons' position in korc_profiles.

There are two types of analytical plasma profiles that can be used in KORC: 3rd degree polynomial radial plasma profiles,

and radial plasma profiles with a dependency:

,

where is the radial coordinate in toroidal coordinates, is a given plasma parameter at the magnetic axis, and is the plasma radius as measured from the magnetic axis to the last closed flux surface. Notice that the larger is, the more uniform the radial profiles are.

Arguments

Type IntentOptional AttributesName
type(KORC_PARAMS), intent(in) :: params

Core KORC simulation parameters.

type(PROFILES), intent(out) :: P

An instance of KORC's derived type PROFILES containing all the information about the plasma profiles used in the simulation. See korc_types and korc_profiles.

type(FIELDS), intent(in) :: F

String containing the type of electron density profile to be used in the simulation.

String containing the type of electron temperature profile to be used in the simulation.

String containing the type of profile to be used in the simulation.

Full path to the HDF5 file containing the pre-computed plasma profiles.

Plasma radius as measured from the magnetic axis.

Electron density at the magnetic axis .

Electron temperature at the magnetic axis .

at the magnetic axis .

Exponent used in of the electron density profile.

Exponent used in of the electron temperature profile.

Exponent used in of the profile.

Coefficients of the polynomial electron density profile. See detailed description above, a_ne=(,,,).

Coefficients of the polynomial electron temperature profile. See detailed description above, a_Te=(,,,).

Coefficients of the profile. See detailed description above, a_Zeff=(,,,).

Flag to indicate if the plasma profiles are axisymmetric.


Calls

proc~~initialize_profiles~~CallsGraph proc~initialize_profiles initialize_profiles proc~allocate_2d_profiles_arrays ALLOCATE_2D_PROFILES_ARRAYS proc~initialize_profiles->proc~allocate_2d_profiles_arrays proc~load_profiles_data_from_hdf5 load_profiles_data_from_hdf5 proc~initialize_profiles->proc~load_profiles_data_from_hdf5 proc~load_profiles_data_from_hdf5->proc~allocate_2d_profiles_arrays h5fclose_f h5fclose_f proc~load_profiles_data_from_hdf5->h5fclose_f interface~load_array_from_hdf5 load_array_from_hdf5 proc~load_profiles_data_from_hdf5->interface~load_array_from_hdf5 proc~allocate_3d_profiles_arrays ALLOCATE_3D_PROFILES_ARRAYS proc~load_profiles_data_from_hdf5->proc~allocate_3d_profiles_arrays h5fopen_f h5fopen_f proc~load_profiles_data_from_hdf5->h5fopen_f interface~load_from_hdf5 load_from_hdf5 proc~load_profiles_data_from_hdf5->interface~load_from_hdf5 proc~rload_1d_array_from_hdf5 rload_1d_array_from_hdf5 interface~load_array_from_hdf5->proc~rload_1d_array_from_hdf5 proc~rload_2d_array_from_hdf5 rload_2d_array_from_hdf5 interface~load_array_from_hdf5->proc~rload_2d_array_from_hdf5 proc~rload_3d_array_from_hdf5 rload_3d_array_from_hdf5 interface~load_array_from_hdf5->proc~rload_3d_array_from_hdf5 proc~iload_from_hdf5 iload_from_hdf5 interface~load_from_hdf5->proc~iload_from_hdf5 proc~rload_from_hdf5 rload_from_hdf5 interface~load_from_hdf5->proc~rload_from_hdf5 proc~korc_abort korc_abort proc~iload_from_hdf5->proc~korc_abort h5dclose_f h5dclose_f proc~iload_from_hdf5->h5dclose_f h5dread_f h5dread_f proc~iload_from_hdf5->h5dread_f proc~rload_from_hdf5->proc~korc_abort proc~rload_from_hdf5->h5dclose_f proc~rload_from_hdf5->h5dread_f proc~rload_1d_array_from_hdf5->proc~korc_abort proc~rload_1d_array_from_hdf5->h5dclose_f proc~rload_1d_array_from_hdf5->h5dread_f proc~rload_2d_array_from_hdf5->proc~korc_abort proc~rload_2d_array_from_hdf5->h5dclose_f proc~rload_2d_array_from_hdf5->h5dread_f proc~rload_3d_array_from_hdf5->proc~korc_abort proc~rload_3d_array_from_hdf5->h5dclose_f proc~rload_3d_array_from_hdf5->h5dread_f mpi_abort mpi_abort proc~korc_abort->mpi_abort

Called by

proc~~initialize_profiles~~CalledByGraph proc~initialize_profiles initialize_profiles program~main main program~main->proc~initialize_profiles

Contents

Source Code


Source Code

  subroutine initialize_profiles(params,P,F)
    !! @note Subroutine that initializes the parameters of analytical 
    !! or pre-computed plasma profiles for being used in the
    !! simulation. @endnote
    !! KORC can run using either analytical and pre-computed plasma
    !! profiles. Pre-computed plasma profiles, as in the case of
    !! pre-computed electric or magnetic fields, are interpolated to
    !! electrons' position in [[korc_profiles]].
    !!
    !! There are two types of analytical plasma profiles that can be
    !! used in KORC: 3rd degree polynomial radial plasma profiles,
    !!
    !! $$f(r) = a_3r^3 + a_2r^2 +a_1r + a_0,$$
    !!
    !! and radial plasma profiles with a \(\tanh(r)\) dependency:
    !!
    !! $$f(r) = f_0\left[1 - \tanh^n\left(\frac{2r}{a}\right)\right]$$,
    !!
    !! where \(r\) is the radial coordinate in toroidal coordinates,
    !! \(f_0\) is a given plasma parameter at the magnetic axis,
    !! and \(a\) is the plasma radius as measured from the magnetic
    !! axis to the last closed flux surface. Notice that the larger
    !! \(n\) is, the more uniform the radial profiles are.

    TYPE(KORC_PARAMS), INTENT(IN)   :: params
    !! Core KORC simulation parameters.
    TYPE(PROFILES), INTENT(OUT)     :: P
    !! An instance of KORC's derived type PROFILES containing all 
    !! the information about the plasma profiles used in the
    !! simulation. See [[korc_types]] and [[korc_profiles]].
    TYPE(FIELDS), INTENT(IN)     :: F
    !CHARACTER(MAX_STRING_LENGTH)    :: ne_profile
    !! String containing the type of electron density profile 
    !! to be used in the simulation.
    !CHARACTER(MAX_STRING_LENGTH)    :: Te_profile
    !! String containing the type of electron temperature profile 
    !! to be used in the simulation.
    !CHARACTER(MAX_STRING_LENGTH)    :: Zeff_profile  
    !! String containing the type of \(Z_{eff}\) profile to be used
    !! in the simulation.
    !CHARACTER(MAX_STRING_LENGTH)    :: filename
    !! Full path to the HDF5 file containing the pre-computed
    !! plasma profiles.
    !REAL(rp)                        :: radius_profile
    !! Plasma radius \(a\) as measured from the magnetic axis.
    !REAL(rp)                        :: neo
    !! Electron density at the magnetic axis \(f_0 = n_{e,0}\).
    !REAL(rp)                        :: Teo
    !! Electron temperature at the magnetic axis \(f_0 = T_{e,0}\).
    !REAL(rp)                        :: Zeffo
    !! \(Z_{eff}\) at the magnetic axis \(f_0 = Z_{eff,0}\).
    !REAL(rp)                        :: n_ne
    !! Exponent \(n\) used in \(\tanh^n(r)\) of the electron
    !! density profile.
    !REAL(rp)                        :: n_Te
    !! Exponent \(n\) used in \(\tanh^n(r)\) of the electron
    !! temperature profile.
    !REAL(rp)                        :: n_Zeff
    !! Exponent \(n\) used in \(\tanh^n(r)\) of the \(Z_{eff}\) profile.
    !REAL(rp), DIMENSION(4)          :: a_ne
    !! Coefficients of the polynomial electron density profile. 
    !! See detailed description above,
    !! a_ne=(\(a_{0}\),\(a_{2}\),\(a_{3}\),\(a_{4}\)).
    !REAL(rp), DIMENSION(4)          :: a_Te
    !! Coefficients of the polynomial electron temperature profile. 
    !! See detailed description above,
    !! a_Te=(\(a_{0}\),\(a_{2}\),\(a_{3}\),\(a_{4}\)).
    !REAL(rp), DIMENSION(4)          :: a_Zeff
    !! Coefficients of the \(Z_{eff}\) profile. See detailed
    !! description above, a_Zeff=(\(a_{0}\),\(a_{2}\),\(a_{3}\),\(a_{4}\)).
    !LOGICAL                         :: axisymmetric
    !! Flag to indicate if the plasma profiles are axisymmetric.
    INTEGER  :: ii,kk
    !REAL(rp)  ::  n_REr0
    !REAL(rp)  ::  n_tauion
    !REAL(rp)  ::  n_lamfront
    !REAL(rp)  ::  n_lamback,n_lamshelf,n_shelfdelay,n_tauin,n_tauout,n_shelf
    REAL(rp)  ::  rm,r_a!,psiN_0

    !NAMELIST /plasmaProfiles/ radius_profile,ne_profile,neo,n_ne,a_ne, &
    !     Te_profile,Teo,n_Te,a_Te,n_REr0,n_tauion,n_lamfront,n_lamback, &
    !     Zeff_profile,Zeffo,n_Zeff,a_Zeff,filename,axisymmetric, &
    !     n_lamshelf,n_shelfdelay,n_tauin,n_tauout,n_shelf,psiN_0

    if (params%mpi_params%rank .EQ. 0) then
       write(output_unit_write,'("* * * * * * * * INITIALIZING PROFILES * * * * * * * *")')
    end if
    
    if (params%profile_model(1:10).eq.'ANALYTICAL') then
       !open(unit=default_unit_open,file=TRIM(params%path_to_inputs), &
       !     status='OLD',form='formatted')
       !read(default_unit_open,nml=plasmaProfiles)
       !close(default_unit_open)

       P%a = radius_profile
       P%R0=F%Ro
       P%Z0=F%Zo 

       
       P%ne_profile = TRIM(ne_profile)
       P%neo = neo
       P%n_ne = n_ne
       P%a_ne = a_ne

       P%n_REr0=n_REr0
       P%n_tauion=n_tauion
       P%n_tauin=n_tauin
       P%n_tauout=n_tauout
       P%n_shelfdelay=n_shelfdelay
       P%n_lamfront=n_lamfront
       P%n_lamback=n_lamback
       P%n_lamshelf=n_lamshelf
       P%n_shelf=n_shelf
       P%psiN_0=psiN_0
        
       P%Te_profile = TRIM(Te_profile)
       P%Teo = Teo*C_E ! Converted to Joules
       P%n_Te = n_Te
       P%a_Te = a_Te

       P%Zeff_profile = TRIM(Zeff_profile)
       P%Zeffo = Zeffo
       P%n_Zeff = n_Zeff
       P%a_Zeff = a_Zeff


       if (params%mpi_params%rank .EQ. 0) then
          write(output_unit_write,'("ANALYTICAL")')
          write(output_unit_write,'("ne profile: ",A20)') P%ne_profile
          write(output_unit_write,'("Te profile: ",A20)') P%Te_profile
          write(output_unit_write,'("Zeff profile: ",A20)') P%Zeff_profile
       end if
       
       if (params%field_eval.eq.'interp') then

          P%axisymmetric = axisymmetric

          P%dims(1) = F%dims(1)
          P%dims(3) = F%dims(3)

          call ALLOCATE_2D_PROFILES_ARRAYS(params,P)

          P%X%R=F%X%R
          P%X%Z=F%X%Z

          do ii=1_idef,P%dims(1)
             do kk=1_idef,P%dims(3)

                rm=sqrt((P%X%R(ii)-P%R0)**2+(P%X%Z(kk)-P%Z0)**2)
                r_a=rm/P%a

                SELECT CASE (TRIM(P%ne_profile))
                CASE('FLAT')
                   P%ne_2D(ii,kk) = P%neo
                CASE('SPONG')
                   P%ne_2D(ii,kk) = P%neo*(1._rp-0.2*r_a**8)+P%n_ne
                CASE('RE-EVO')                   
                   !flat profile placeholder, updates every timestep
                   P%ne_2D(ii,kk) = P%neo
                CASE('RE-EVO1')                   
                   !flat profile placeholder, updates every timestep
                   P%ne_2D(ii,kk) = P%neo
                CASE('RE-EVO-PSI')                   
                   !flat profile placeholder, updates every timestep
                   P%ne_2D(ii,kk) = P%neo
                CASE('RE-EVO-PSIN-SG')                   
                   !flat profile placeholder, updates every timestep
                   P%ne_2D(ii,kk) = P%neo
                CASE('RE-EVO-PSIP-G')                   
                   !flat profile placeholder, updates every timestep
                   P%ne_2D(ii,kk) = P%neo
                CASE('RE-EVO-PSIP-G1')                   
                   !flat profile placeholder, updates every timestep
                   P%ne_2D(ii,kk) = P%neo
                CASE('MST_FSA')
                   P%ne_2D(ii,kk) = P%neo*(1._rp-r_a**4._rp)**4._rp
                CASE DEFAULT
                   P%ne_2D(ii,kk) = P%neo
                END SELECT

                SELECT CASE (TRIM(P%Te_profile))
                CASE('FLAT')       
                   P%Te_2D(ii,kk) = P%Teo
                CASE('SPONG')
                   P%Te_2D(ii,kk) = (P%Teo-P%n_Te)*(1._rp-0.6*r_a**2)**2+P%n_Te
                CASE('MST_FSA')
                   P%Te_2D(ii,kk) = P%Teo*(1._rp-r_a**8._rp)**4._rp
                CASE DEFAULT
                   P%Te_2D(ii,kk) = P%Teo
                END SELECT

                SELECT CASE (TRIM(P%Zeff_profile))
                CASE('FLAT') 
                   P%Zeff_2D(ii,kk) = P%Zeffo
                CASE('SPONG')
                   P%Zeff_2D(ii,kk) = P%Zeffo
                CASE DEFAULT
                   P%Zeff_2D(ii,kk) = P%Zeffo
                END SELECT

                
             end do
          end do


          
       end if

    else if (params%profile_model(1:8).eq.'EXTERNAL') then
       !open(unit=default_unit_open,file=TRIM(params%path_to_inputs), &
       !     status='OLD',form='formatted')
       !read(default_unit_open,nml=plasmaProfiles)
       !close(default_unit_open)

       P%a = radius_profile
       P%R0=F%Ro
       P%Z0=F%Zo 
       
       P%ne_profile = TRIM(ne_profile)
       P%neo = neo
       P%Te_profile = TRIM(Te_profile)
       P%Teo = Teo*C_E ! Converted to Joules
       P%Zeff_profile = TRIM(Zeff_profile)
       P%Zeffo = Zeffo

       if (params%mpi_params%rank .EQ. 0) then
          write(output_unit_write,'("EXTERNAL")')
       end if
       
       P%filename = TRIM(filename)
       P%axisymmetric = axisymmetric

       call load_profiles_data_from_hdf5(params,P)
    else if (params%profile_model.eq.'UNIFORM') then
       !open(unit=default_unit_open,file=TRIM(params%path_to_inputs), &
       !     status='OLD',form='formatted')
       !read(default_unit_open,nml=plasmaProfiles)
       !close(default_unit_open)

       P%a = radius_profile
       P%R0=F%Ro
       P%Z0=F%Zo 
       
       P%ne_profile = TRIM(ne_profile)
       P%neo = neo
       P%n_ne = 0.0_rp
       P%a_ne = (/0.0_rp,0.0_rp,0.0_rp,0.0_rp/)

       P%Te_profile = TRIM(Te_profile)
       P%Teo = Teo*C_E ! Converted to Joules
       P%n_Te = 0.0_rp
       P%a_Te = (/0.0_rp,0.0_rp,0.0_rp,0.0_rp/)

       P%Zeff_profile = TRIM(Zeff_profile)
       P%Zeffo = Zeffo
       P%n_Zeff = 0.0_rp
       P%a_Zeff = (/0.0_rp,0.0_rp,0.0_rp,0.0_rp/)

       if (params%mpi_params%rank .EQ. 0) then
          write(output_unit_write,'("UNIFORM")')
       end if
       
    endif

    if (params%mpi_params%rank .EQ. 0) then
       write(output_unit_write,'("* * * * * * * * * * * * * * * * * * * * * * * * *")')
    end if
    
  end subroutine initialize_profiles