nc_loop Subroutine

public subroutine nc_loop(file_id, nout, time, phi2, apar2, bpar2, igomega, h, hk, omega, omegaavg, woutunits, tunits, phitot, ql_metric, write_omega, write_heating, write_ql_metric)

Write various fields, heating, and frequency diagnostics to netCDF

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: file_id

NetCDF ID of the file to write to

integer, intent(in) :: nout

Current output timestep

real, intent(in) :: time

Current simulation time

real, intent(out) :: phi2

Volume average of fields squared

real, intent(out) :: apar2

Volume average of fields squared

real, intent(out) :: bpar2

Volume average of fields squared

integer, intent(in) :: igomega

Theta-index to measure fields at

type(heating_diagnostics), intent(in) :: h

Heating diagnostics

type(heating_diagnostics), intent(in), dimension(:,:) :: hk

Heating per-mode

complex, intent(in), dimension(:,:) :: omega

Fields, frequency and time-averaged frequency

complex, intent(in), dimension(:,:) :: omegaavg

Fields, frequency and time-averaged frequency

real, intent(in), dimension (:) :: woutunits

Normalisation factor for frequency

real, intent(in), dimension (:) :: tunits

Normalisation factor for time

real, intent(in), dimension(:,:) :: phitot

FIXME: add documentation

real, intent(in), dimension(:,:) :: ql_metric

FIXME: add documentation

logical, intent(in) :: write_omega
logical, intent(in) :: write_heating
logical, intent(in) :: write_ql_metric

Contents

Source Code


Source Code

  subroutine nc_loop (file_id, nout, time, &
       phi2, apar2, bpar2, igomega, &
       h, hk, omega, omegaavg, woutunits, tunits, phitot, ql_metric, &
       write_omega, write_heating, write_ql_metric)
    use gs2_heating, only: heating_diagnostics
    use run_parameters, only: has_phi, has_apar, has_bpar
    use fields_arrays, only: phinew, aparnew, bparnew
# ifdef NETCDF
    use run_parameters, only: wstar_units
    use neasyf, only: neasyf_write
# else
    use theta_grid, only: ntgrid
    use kt_grids, only: naky, ntheta0
# endif
    implicit none
    !> NetCDF ID of the file to write to
    integer, intent(in) :: file_id
    !> Current output timestep
    integer, intent (in) :: nout
    !> Current simulation time
    real, intent (in) :: time
    !> Volume average of fields squared
    real, intent(out) :: phi2, apar2, bpar2
    !> Theta-index to measure fields at
    integer, intent(in) :: igomega
    !> Normalisation factor for frequency
    real, dimension (:), intent (in) :: woutunits
    !> Normalisation factor for time
    real, dimension (:), intent (in) :: tunits
    !> Fields, frequency and time-averaged frequency
    complex, dimension(:,:), intent (in) :: omega, omegaavg
    !> FIXME: add documentation
    real, dimension(:,:), intent (in) :: phitot, ql_metric
    !> Heating diagnostics
    type(heating_diagnostics), intent (in) :: h
    !> Heating per-mode
    type(heating_diagnostics), dimension(:,:), intent (in) :: hk
    logical, intent(in) :: write_omega, write_heating, write_ql_metric
# ifndef NETCDF
    real, dimension(ntheta0, naky) :: field2_by_mode
    real, dimension(ntheta0) :: field2_by_kx
    real, dimension(naky) :: field2_by_ky
# endif
    ! Default set the intent(out) args
    phi2 = 0. ; apar2 = 0. ; bpar2 = 0.
# ifdef NETCDF

    call neasyf_write(file_id, time_dim, time, dim_names=[time_dim], start=[nout])

    if (wstar_units) then
       call neasyf_write(file_id, "t_wstar", time * tunits, dim_names=[ky_dim, time_dim], &
            start=[1, nout], long_name="Time (wstar)", units="L/vt")
    end if

    if (has_phi) then
       call nc_write_field(file_id, nout, phinew, "phi", "Electrostatic potential", &
            igomega, write_phi_t, phi2)
    end if

    if (has_apar) then
       call nc_write_field(file_id, nout, aparnew, "apar", "Parallel magnetic potential", &
            igomega, write_apar_t, apar2)
    end if

    if (has_bpar) then
       call nc_write_field(file_id, nout, bparnew, "bpar", "Parallel magnetic field", &
            igomega, write_bpar_t, bpar2)
    end if

    if (write_heating) call nc_write_heating(file_id, nout, h, hk)

    if (write_omega) call nc_write_omega(file_id, nout, omega, omegaavg, woutunits)

    if (write_ql_metric) call nc_write_ql_metric(file_id, nout, ql_metric)

    call neasyf_write(file_id, "phtot", phitot, dim_names=mode_dims, start=starts(3, nout), &
         long_name="Roughly the field line averaged summed field amplitudes")
# else
    ! Make sure we set the intent(out) variables
    if (has_phi) call calculate_field_sums(phinew, field2_by_mode, field2_by_kx, &
         field2_by_ky, phi2)
    if (has_apar) call calculate_field_sums(aparnew, field2_by_mode, field2_by_kx, &
         field2_by_ky, apar2)
    if (has_bpar) call calculate_field_sums(bparnew, field2_by_mode, field2_by_kx, &
         field2_by_ky, bpar2)
    UNUSED_DUMMY(file_id); UNUSED_DUMMY(nout); UNUSED_DUMMY(time); UNUSED_DUMMY(igomega)
    UNUSED_DUMMY(h); UNUSED_DUMMY(hk); UNUSED_DUMMY(omega); UNUSED_DUMMY(omegaavg)
    UNUSED_DUMMY(woutunits); UNUSED_DUMMY(tunits); UNUSED_DUMMY(phitot)
    UNUSED_DUMMY(ql_metric); UNUSED_DUMMY(write_omega); UNUSED_DUMMY(write_heating)
    UNUSED_DUMMY(write_ql_metric)
# endif
  end subroutine nc_loop