efit_read_and_set Subroutine

private subroutine efit_read_and_set(self, inputs)

This subroutine reads an EFIT output file containing the axisymmetric magnetic field geometry on a rectangular domain defined by the coordinates (R,Z). It reads and stores the following quantities (among others).

  • eq_R: R grid
  • eq_Z: Z grid
  • fp: F on psibar grid
  • eqpsi_2d: psi on (R,Z) grid
  • R_mag: major radius of the magnetic axis
  • Z_mag: elevation of the magnetic axis
  • rwid: total width of the domain
  • rleft: position of leftmost point of domain
  • zhei: total height of domain

Type Bound

eeq_type

Arguments

Type IntentOptional Attributes Name
class(eeq_type), intent(inout) :: self
type(geo_input_type), intent(in) :: inputs

Contents

Source Code


Source Code

  subroutine efit_read_and_set(self, inputs)
    implicit none
    class(eeq_type), intent(in out) :: self
    type(geo_input_type), intent(in) :: inputs
    real :: xdum, rwid, rleft, zhei, rcentr, bcentr, zmid
    real, dimension(:), allocatable :: dummy, spsi_bar, sefit_R, sefit_Z, rbbbs, zbbbs, f, p, q
    real, dimension(:, :), allocatable :: sefit_psi
    character(10) ::  char
    integer :: i, j, ndum, nw_in, nh_in, funit, nbbbs
    self%type_name = 'EFIT'
    self%filename = trim(adjustl(inputs%eqfile))
    open(newunit = funit,  file = self%filename, status='old',form='formatted')
    read(funit, '(5(a10),i2,i4,i4)') char, char, char, char, char, i, nw_in, nh_in
    allocate(spsi_bar(nw_in), sefit_R(nw_in), sefit_Z(nh_in), dummy(nw_in))
    allocate(sefit_psi(nw_in, nh_in), f(nw_in), p(nw_in), q(nw_in))

    read(funit, 2020) rwid, zhei, rcentr, rleft, zmid
    read(funit, 2020) self%R_mag, self%Z_mag, self%psi_0, self%psi_a, bcentr
    read(funit, 2020) xdum,xdum,xdum,xdum,xdum
    read(funit, 2020) xdum,xdum,xdum,xdum,xdum

    read(funit, 2020) (f(j), j = 1, nw_in)
    read(funit, 2020) (p(j), j = 1, nw_in)
    read(funit, 2020) (dummy(j), j = 1, nw_in)
    read(funit, 2020) (dummy(j), j = 1 ,nw_in)
    read(funit, 2020) ((sefit_psi(i,j) , i = 1, nw_in) , j = 1, nh_in)
    read(funit, 2020) (q(j), j = 1, nw_in)

    read(funit, '(2i5)') nbbbs, ndum
    allocate(rbbbs(nbbbs), zbbbs(nbbbs))
    read(funit, 2020) (rbbbs(i), zbbbs(i) , i = 1, nbbbs)
    close(funit)

    ! pbar is defined by pbar = (psi-psi_0)/(psi_a-psi_0) fp and q are
    ! functions of pbar
    do i = 1, nw_in
       spsi_bar(i) = float(i - 1) / float(nw_in - 1)
       sefit_R(i) = rleft + rwid * float(i - 1) / float(nw_in - 1)
    end do

    do j = 1, nh_in
       sefit_Z(j) = ((float(j - 1) / float(nh_in - 1) ) - 0.5) * zhei !Should we add zmid?
    end do

    call self%shared_setup(nw_in, nh_in, spsi_bar, sefit_R, sefit_Z, f, p, q, sefit_psi, &
         nbbbs, rbbbs, zbbbs, bcentr, rleft, rwid, -0.5*zhei, zhei, inputs%big, .true.) !Should -0.5*zhei be zmid
2020 format (5e16.9)
  end subroutine efit_read_and_set