real_init Subroutine

private subroutine real_init(is_list_run, gs2_diagnostics_config_in, header)

Read the input parameters; open the various text output files according to the relevant input flags; allocate and zero the module-level flux arrays, gs2_diagnostics_knobs and gs2_diagnostics_knobs

Arguments

Type IntentOptional Attributes Name
logical, intent(in) :: is_list_run

If true, this is a "list-mode run" and so turn off print_flux_line and print_line if set

type(diagnostics_config_type), intent(in), optional :: gs2_diagnostics_config_in

Input parameters for this module

type(standard_header_type), intent(in), optional :: header

Header for files with build and run information


Contents

Source Code


Source Code

  subroutine real_init(is_list_run, gs2_diagnostics_config_in, header)
    use run_parameters, only: has_apar
    use file_utils, only: open_output_file, get_unused_unit
    use kt_grids, only: naky, ntheta0
    use mp, only: proc0
    use standard_header, only: standard_header_type
    use diagnostics_fluxes, only: init_diagnostics_fluxes
    use diagnostics_kinetic_energy_transfer, only: init_diagnostics_kinetic_energy_transfer
    implicit none
    !> If true, this is a "list-mode run" and so turn off
    !> [[gs2_diagnostics_knobs:print_flux_line]] and
    !> [[gs2_diagnostics_knobs:print_line]] if set
    logical, intent (in) :: is_list_run
    !> Input parameters for this module
    type(diagnostics_config_type), intent(in), optional :: gs2_diagnostics_config_in
    !> Header for files with build and run information
    type(standard_header_type), intent(in), optional :: header
    ! Actual value for optional `header` input
    type(standard_header_type) :: local_header

    if (present(header)) then
      local_header = header
    else
      local_header = standard_header_type()
    end if

    call read_parameters(is_list_run, gs2_diagnostics_config_in)

    if (proc0) then
       if (write_ascii) then
          call open_output_file (out_unit, ".out")          
       end if

       if (write_cross_phase .and. write_ascii) then
          call open_output_file (phase_unit, ".phase")
       end if

       if (write_heating .and. write_ascii) then
          call open_output_file (heat_unit, ".heat")
          call open_output_file (heat_unit2, ".heat2")
       end if

       if (write_verr .and. write_ascii) then
          call open_output_file (res_unit, ".vres")
          call open_output_file (lpc_unit, ".lpc")
          if (write_max_verr) call open_output_file (res_unit2, ".vres2")
       end if

       if (write_cerr .and. write_ascii) then
          call open_output_file (cres_unit, ".cres")
       end if

       if (write_parity .and. write_ascii) then
          call open_output_file (parity_unit, ".parity")
       end if

       if (write_g .and. write_ascii) then
          call open_output_file(dist_unit, ".dist")
       end if

       if (write_gyx .and. write_ascii) then
          call open_output_file(yxdist_unit, ".yxdist")
       end if

       !GGH J_external, only if A_parallel is being calculated.
       if (write_jext .and. has_apar) then
          if (write_ascii) then
             call open_output_file (jext_unit, ".jext")
          end if
       else
          write_jext = .false.
       end if

       if (dump_check1 .and. proc0) then
          call get_unused_unit (dump_check1_unit)
          ! We should switch to open_output_file here.
          open (unit=dump_check1_unit, file="dump.check1", status="unknown")
       end if

       if (dump_check2 .and. proc0) then
          call open_output_file (dump_check2_unit, ".dc2")
       end if

       if (write_ascii) then
          write (unit=out_unit, fmt='(a)') local_header%to_string(file_description="GS2 output file")
       end if

       allocate (omegahist(0:navg-1,ntheta0,naky))
       omegahist = 0.0
       if(.not. allocated(conv_heat)) allocate(conv_heat(0:conv_nstep_av/nwrite-1))

    end if

    call init_diagnostics_fluxes()

    if (write_kinetic_energy_transfer) call init_diagnostics_kinetic_energy_transfer
    if (.not. allocated(omega)) then
       allocate(omega(ntheta0, naky))
       allocate(omegaavg(ntheta0, naky))
       omega=0
       omegaavg=0
    end if

  end subroutine real_init