fields_arrays.f90 Source File


Contents

Source Code


Source Code

!> Provides storage for the EM potentials
!! and some related data
module fields_arrays
  use constants, only: run_name_size
  implicit none

  private

  public :: phi, apar, bpar, phinew, aparnew, bparnew
  public :: gf_phi, gf_apar, gf_bpar, gf_phinew, gf_aparnew, gf_bparnew
  public :: apar_ext, time_field, response_file, time_field_mpi
  public :: time_field_invert, time_field_invert_mpi
  public :: time_dump_response, time_read_response
  public :: get_specific_response_file_name

  !Main fields
  complex, dimension (:,:,:), allocatable :: phi,    apar,    bpar
  complex, dimension (:,:,:), allocatable :: phinew, aparnew, bparnew
  complex, dimension (:,:,:), allocatable :: gf_phi, gf_apar, gf_bpar
  complex, dimension (:,:,:), allocatable :: gf_phinew, gf_aparnew, gf_bparnew
  !For antenna
  complex, dimension (:,:,:), allocatable :: apar_ext
  ! (-ntgrid:ntgrid,ntheta0,naky) replicated

  !Timing data
  real :: time_field(2) = 0., time_field_mpi = 0.
  real :: time_field_invert(2) = 0., time_field_invert_mpi = 0.
  real :: time_dump_response(2) = 0.
  real :: time_read_response(2) = 0.

  !For response data
  character(run_name_size) :: response_file = 'NOT_SET'

contains

  !> This function returns the response file name when given
  !> characteristic data.
  function get_specific_response_file_name(ik, is, dt, suffix) result(file_name)
    use optionals, only: get_option_with_default
    implicit none
    integer, intent(in) :: ik, is
    real, intent(in) :: dt
    character(len = *), optional, intent(in) :: suffix
    character(len = 256) :: file_name

    character(len = 64) :: suffix_local
    character(len = 64), parameter :: suffix_default='.response'
    character(len = 14) :: dt_tmp

    !Set file suffix
    suffix_local = get_option_with_default(suffix, suffix_default)

    !First write the time step into a temporary so that we can strip
    !whitespace.
    write(dt_tmp,'(E14.6E3)') dt

    write(file_name,'(A,"_ik_",I0,"_is_",I0,"_dt_",A,A)') trim(response_file),ik,is,trim(adjustl(dt_tmp)),trim(suffix_local)
  end function get_specific_response_file_name

end module fields_arrays