read_parameters Subroutine

private subroutine read_parameters(theta_grid_parameters_config_in)

FIXME : Add documentation

Arguments

Type IntentOptional Attributes Name
type(theta_grid_parameters_config_type), intent(in), optional :: theta_grid_parameters_config_in

Contents

Source Code


Source Code

  subroutine read_parameters(theta_grid_parameters_config_in)
    use unit_tests, only: debug_message
    use geo_utils, only: geo_type_miller, geo_type_global, &
         geo_type_generalized_elongation, geo_type_fourier_series, &
         geo_type_miller_extended_harmonic, mxh_max_moments
    use warning_helpers, only: not_exactly_equal, is_not_zero
    implicit none
    type(theta_grid_parameters_config_type), intent(in), optional :: theta_grid_parameters_config_in
    integer, parameter :: verb=3
    real :: c0_mxh
    call debug_message(verb, "theta_grid_params::read_parameters start")

    if (present(theta_grid_parameters_config_in)) theta_grid_parameters_config = theta_grid_parameters_config_in

    call theta_grid_parameters_config%init(name = 'theta_grid_parameters', requires_index = .false.)

    ! Copy out internal values into module level parameters
    associate(self => theta_grid_parameters_config)
#include "theta_grid_parameters_copy_out_auto_gen.inc"
    end associate

    if (kp > 0.) pk = 2.*kp
    ! Make sure ntheta is even
    ntheta = 2 * (ntheta / 2)

    ! Print warning if the user has specified non-default values for the
    ! geometrical input parameters of other geometry types ! JRB
    select case (geoType)
       case (geo_type_miller)
          if (mMode/=2) write (*,*) "WARNING: ignoring value of mMode, value not needed"
          if (nMode/=3) write (*,*) "WARNING: ignoring value of nMode, value not needed"
          if (not_exactly_equal(deltam, 1.0)) write (*,*) "WARNING: ignoring value of deltam, did you mean akappa?"
          if (not_exactly_equal(deltan, 1.0)) write (*,*) "WARNING: ignoring value of deltan, did you mean tri?"
          if (not_exactly_equal(delta2, 1.0)) write (*,*) "WARNING: ignoring value of delta2, did you mean akappa?"
          if (not_exactly_equal(delta3, 1.0)) write (*,*) "WARNING: ignoring value of delta3, did you mean tri?"
          if (is_not_zero(deltampri)) write (*,*) "WARNING: ignoring value of deltampri, did you mean akappri?"
          if (is_not_zero(deltanpri)) write (*,*) "WARNING: ignoring value of deltanpri, did you mean tripri?"
          if (is_not_zero(thetam)) write (*,*) "WARNING: ignoring value of thetam, did you mean thetak?"
          if (is_not_zero(thetan)) write (*,*) "WARNING: ignoring value of thetan, did you mean thetad?"
          if (is_not_zero(theta2)) write (*,*) "WARNING: ignoring value of theta2, did you mean thetak?"
          if (is_not_zero(theta3)) write (*,*) "WARNING: ignoring value of theta3, did you mean thetad?"
          if (not_exactly_equal(aSurf, 1.0)) write (*,*) "WARNING: ignoring value of aSurf, value not needed"
          if (not_exactly_equal(raxis, 3.0)) write (*,*) "WARNING: ignoring value of raxis, did you mean shift?"
          if (is_not_zero(zaxis)) write (*,*) "WARNING: ignoring value of zaxis, did you mean shiftVert?"
          deltam = akappa
          deltan = tri
          deltampri = akappri
          deltanpri = tripri
          thetam = thetak
          thetan = thetad
       case (geo_type_global)
          if (mMode/=2) write (*,*) "WARNING: ignoring value of mMode, value not needed"
          if (nMode/=3) write (*,*) "WARNING: ignoring value of nMode, value not needed"
          if (not_exactly_equal(deltam, 1.0)) write (*,*) "WARNING: ignoring value of deltam, did you mean delta2?"
          if (not_exactly_equal(deltan, 1.0)) write (*,*) "WARNING: ignoring value of deltan, did you mean delta3?"
          if (not_exactly_equal(akappa, 1.0)) write (*,*) "WARNING: ignoring value of akappa, did you mean delta2?"
          if (is_not_zero(tri)) write (*,*) "WARNING: ignoring value of tri, did you mean delta3?"
          if (is_not_zero(deltampri)) write (*,*) "WARNING: ignoring value of deltampri, value not needed"
          if (is_not_zero(deltanpri)) write (*,*) "WARNING: ignoring value of deltanpri, value not needed"
          if (is_not_zero(akappri)) write (*,*) "WARNING: ignoring value of akappri, value not needed"
          if (is_not_zero(tripri)) write (*,*) "WARNING: ignoring value of tripri, value not needed"
          if (is_not_zero(thetam)) write (*,*) "WARNING: ignoring value of thetam, did you mean theta2?"
          if (is_not_zero(thetan)) write (*,*) "WARNING: ignoring value of thetan, did you mean theta3?"
          if (is_not_zero(thetak)) write (*,*) "WARNING: ignoring value of thetak, did you mean theta2?"
          if (is_not_zero(thetad)) write (*,*) "WARNING: ignoring value of thetad, did you mean theta3?"
          if (is_not_zero(shift)) write (*,*) "WARNING: ignoring value of shift, did you mean Raxis?"
          if (is_not_zero(shiftVert)) write (*,*) "WARNING: ignoring value of shiftVert, did you mean Zaxis?"

          if (raxis>=Rmaj+aSurf) write (*,*) "WARNING: Raxis>=Rmag+aSurf, may specify non-nested flux surfaces"
          if (raxis<=Rmaj-aSurf) write (*,*) "WARNING: Raxis<=Rmag-aSurf, may specify non-nested flux surfaces"
          if (zaxis>=aSurf) write (*,*) "WARNING: Zaxis>=aSurf, may specify non-nested flux surfaces"
          if (zaxis<=-aSurf) write (*,*) "WARNING: Zaxis<=-aSurf, may specify non-nested flux surfaces"

          shift = raxis
          shiftVert = zaxis
          deltam = delta2
          deltan = delta3
          thetam = theta2
          thetan = theta3
       case (geo_type_generalized_elongation)
          if (not_exactly_equal(delta2, 1.0)) write (*,*) "WARNING: ignoring value of delta2, did you mean deltam?"
          if (not_exactly_equal(delta3, 1.0)) write (*,*) "WARNING: ignoring value of delta3, did you mean deltan?"
          if (not_exactly_equal(akappa, 1.0)) write (*,*) "WARNING: ignoring value of akappa, did you mean deltam?"
          if (is_not_zero(tri)) write (*,*) "WARNING: ignoring value of tri, did you mean deltan?"
          if (is_not_zero(akappri)) write (*,*) "WARNING: ignoring value of akappri, did you mean deltampri?"
          if (is_not_zero(tripri)) write (*,*) "WARNING: ignoring value of tripri, did you mean deltanpri?"
          if (is_not_zero(theta2)) write (*,*) "WARNING: ignoring value of theta2, did you mean thetam?"
          if (is_not_zero(theta3)) write (*,*) "WARNING: ignoring value of theta3, did you mean thetan?"
          if (is_not_zero(thetak)) write (*,*) "WARNING: ignoring value of thetak, did you mean thetam?"
          if (is_not_zero(thetad)) write (*,*) "WARNING: ignoring value of thetad, did you mean thetan?"
          if (not_exactly_equal(aSurf, 1.0)) write (*,*) "WARNING: ignoring value of aSurf, value not needed"
          if (not_exactly_equal(raxis, 3.0)) write (*,*) "WARNING: ignoring value of raxis, did you mean shift?"
          if (is_not_zero(zaxis)) write (*,*) "WARNING: ignoring value of zaxis, did you mean shiftVert?"

       case (geo_type_fourier_series)
          if (not_exactly_equal(delta2, 1.0)) write (*,*) "WARNING: ignoring value of delta2, did you mean deltam?"
          if (not_exactly_equal(delta3, 1.0)) write (*,*) "WARNING: ignoring value of delta3, did you mean deltan?"
          if (not_exactly_equal(akappa, 1.0)) write (*,*) "WARNING: ignoring value of akappa, did you mean deltam?"
          if (is_not_zero(tri)) write (*,*) "WARNING: ignoring value of tri, did you mean deltan?"
          if (is_not_zero(akappri)) write (*,*) "WARNING: ignoring value of akappri, did you mean deltampri?"
          if (is_not_zero(tripri)) write (*,*) "WARNING: ignoring value of tripri, did you mean deltanpri?"
          if (is_not_zero(theta2)) write (*,*) "WARNING: ignoring value of theta2, did you mean thetam?"
          if (is_not_zero(theta3)) write (*,*) "WARNING: ignoring value of theta3, did you mean thetan?"
          if (is_not_zero(thetak)) write (*,*) "WARNING: ignoring value of thetak, did you mean thetam?"
          if (is_not_zero(thetad)) write (*,*) "WARNING: ignoring value of thetad, did you mean thetan?"
          if (not_exactly_equal(aSurf, 1.0)) write (*,*) "WARNING: ignoring value of aSurf, value not needed"
          if (not_exactly_equal(raxis, 3.0)) write (*,*) "WARNING: ignoring value of raxis, did you mean shift?"
          if (is_not_zero(zaxis)) write (*,*) "WARNING: ignoring value of zaxis, did you mean shiftVert?"

       case (geo_type_miller_extended_harmonic)
          if (not_exactly_equal(deltam, 1.0)) write (*,*) "WARNING: ignoring value of deltam, did you mean akappa?"
          if (is_not_zero(deltampri)) write (*,*) "WARNING: ignoring value of deltampri, did you mean akappri?"

          if (mMode/=2) write (*,*) "WARNING: ignoring value of mMode, value not needed"
          if (nMode/=3) write (*,*) "WARNING: ignoring value of nMode, value not needed"
          if (not_exactly_equal(deltan, 1.0)) write (*,*) "WARNING: ignoring value of deltan, value not needed"
          if (not_exactly_equal(delta2, 1.0)) write (*,*) "WARNING: ignoring value of delta2, value not needed"
          if (not_exactly_equal(delta3, 1.0)) write (*,*) "WARNING: ignoring value of delta3, value not needed"
          if (is_not_zero(deltampri)) write (*,*) "WARNING: ignoring value of deltampri, value not needed"
          if (is_not_zero(deltanpri)) write (*,*) "WARNING: ignoring value of deltanpri, value not needed"
          if (is_not_zero(thetam)) write (*,*) "WARNING: ignoring value of thetam, value not needed"
          if (is_not_zero(thetan)) write (*,*) "WARNING: ignoring value of thetan, value not needed"
          if (is_not_zero(theta2)) write (*,*) "WARNING: ignoring value of theta2, value not needed"
          if (is_not_zero(theta3)) write (*,*) "WARNING: ignoring value of theta3, value not needed"
          if (not_exactly_equal(aSurf, 1.0)) write (*,*) "WARNING: ignoring value of aSurf, value not needed"
          if (not_exactly_equal(raxis, 3.0)) write (*,*) "WARNING: ignoring value of raxis, did you mean shift?"
          if (is_not_zero(zaxis)) write (*,*) "WARNING: ignoring value of zaxis, did you mean shiftVert?"

          deltam = akappa
          deltampri = akappri

          if (theta_grid_parameters_config%n_mxh > mxh_max_moments) then
             write(*, '(A, I0, A, I0, A, I0)') &
                  "WARNING: ''n_mxh'' (", &
                  theta_grid_parameters_config%n_mxh, &
                  ") is larger than maximum number of moments for MXH geometry (", &
                  mxh_max_moments, "), truncating to ", mxh_max_moments
          end if

          if (is_not_zero(theta_grid_parameters_config%c0_mxh)) then
             write(*,'("Warning: c0_mxh input is deprecated, set via c_mxh now.")')
          end if

       case default ! We should abort here, but we don't know if local_eq is true so maybe this is fine.
          write (*,*) "ERROR: invalid analytic geometry specification"
    end select

  end subroutine read_parameters