report Subroutine

public subroutine report(header)

FIXME : Add documentation

Arguments

Type IntentOptional Attributes Name
type(standard_header_type), intent(in), optional :: header

Header for files with build and run information


Contents

Source Code


Source Code

  subroutine report(header)
    use antenna, only: check_antenna
    use collisions, only: check_collisions
    use dist_fn, only: check_dist_fn
    use fields, only: check_fields
    use gs2_diagnostics, only: check_gs2_diagnostics, nwrite, save_for_restart
    use gs2_diagnostics, only: nsave, make_movie, nmovie, exit_when_converged, omegatol
    use gs2_reinit, only: delt_adj
    use gs2_time, only: code_dt, user_dt
    use hyper, only: check_hyper
    use init_g, only: check_init_g
    use kt_grids, only: check_kt_grids, is_box, naky, ntheta0, nx, ny
    use le_grids, only: negrid, nlambda
    use nonlinear_terms, only: nonlin, check_nonlinear_terms
    use run_parameters, only: use_old_diagnostics, check_run_parameters
    use run_parameters, only: beta, tite, nstep, wstar_units
    use species, only: check_species, nspec, has_electron_species
    use theta_grid, only: check_theta_grid, ntgrid
    use collisions, only: collision_model_switch
    use collisions, only: collision_model_full, collision_model_ediffuse
    use collisions, only: collision_model_lorentz, collision_model_lorentz_test
    use collisions, only: use_le_layout
    use standard_header, only: standard_header_type
    use file_utils, only: run_name, open_output_file, close_output_file
    implicit none
    !> 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
    real :: alne, dbetadrho_spec, local_omegatol
    logical :: local_exit_when_converged, local_save_for_restart, local_make_movie
    integer :: nmesh, local_nwrite, local_nmovie, local_nsave, report_unit, i
    integer, dimension(4) :: pfacs

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

    call open_output_file (report_unit, ".report")
    write (*, *) "Writing report to '" // trim(run_name) // ".report'"

    write (report_unit, *)

    write (report_unit, fmt='(a)') local_header%to_string(file_description="Input file report from ingen")
    call write_separator(report_unit)

    if (nonlin) then
       if (.not. is_box) then
          write (report_unit, *) 
          write (report_unit, fmt="('################# WARNING #######################')")
          write (report_unit, fmt="('Nonlinear runs must be carried out in a box.')") 
          write (report_unit, fmt="('Set grid_option to box in the kt_grids_knobs namelist')") 
          write (report_unit, fmt="('or set nonlinear_mode to off in the nonlinear_knobs namelist.')") 
          write (report_unit, fmt="('THIS IS AN ERROR.')") 
          write (report_unit, fmt="('################# WARNING #######################')")
          write (report_unit, *) 
       end if

       write (report_unit, fmt="('nmesh=(2*ntgrid+1)*2*nlambda*negrid*nx*ny*nspec')")
       nmesh = (2*ntgrid+1)*2*nlambda*negrid*nx*ny*nspec
       write (report_unit, fmt="('Number of meshpoints:    ',i16)") nmesh

       !
       ! check that nx, ny have no large prime factors
       !
       pfacs = 0
       call pfactors (ny, pfacs)
       if (pfacs(4) /= 0) then
          write (report_unit, *) 
          write (report_unit, fmt="('################# WARNING #######################')")
          write (report_unit, fmt="('ny is a multiple of ',i4)") pfacs(4)
          write (report_unit, fmt="('THIS IS PROBABLY AN ERROR.')")
          write (report_unit, fmt="('ny should have only 2, 3, and/or 5 as prime factors.')")
          write (report_unit, fmt="('################# WARNING #######################')")
          write (report_unit, *) 
       end if
       i = 1
       if (pfacs(1) > 0) i=2**pfacs(1)
       if (pfacs(2) > 0) i=3**pfacs(2)*i
       if (pfacs(3) > 0) i=5**pfacs(3)*i
       if (i /= ny) then
          write (report_unit, *) 
          write (report_unit, fmt="('################# WARNING #######################')")
          write (report_unit, fmt="('ny = ',i3)") ny
          write (report_unit, fmt="('THIS IS PROBABLY AN ERROR.')")
          write (report_unit, fmt="('ny should have only 2, 3, and/or 5 as prime factors.')")
          write (report_unit, fmt="('################# WARNING #######################')")
          write (report_unit, *) 
       end if

       pfacs = 0
       call pfactors (nx, pfacs)
       if (pfacs(4) /= 0) then
          write (report_unit, *) 
          write (report_unit, fmt="('################# WARNING #######################')")
          write (report_unit, fmt="('nx is a multiple of ',i3)") pfacs(4)
          write (report_unit, fmt="('THIS IS PROBABLY AN ERROR.')")
          write (report_unit, fmt="('nx should have only 2, 3, and/or 5 as prime factors.')")
          write (report_unit, fmt="('################# WARNING #######################')")
          write (report_unit, *) 
       end if
       i = 1
       if (pfacs(1) > 0) i=2**pfacs(1)
       if (pfacs(2) > 0) i=3**pfacs(2)*i
       if (pfacs(3) > 0) i=5**pfacs(3)*i
       if (i /= nx) then
          write (report_unit, *) 
          write (report_unit, fmt="('################# WARNING #######################')")
          write (report_unit, fmt="('nx = ',i3)") nx
          write (report_unit, fmt="('THIS IS PROBABLY AN ERROR.')")
          write (report_unit, fmt="('nx should have only 2, 3, and/or 5 as prime factors.')")
          write (report_unit, fmt="('################# WARNING #######################')")
          write (report_unit, *) 
       end if
    else
       write (report_unit, fmt="('nmesh=(2*ntgrid+1)*2*nlambda*negrid*ntheta0*naky*nspec')")
       nmesh = (2*ntgrid+1)*2*nlambda*negrid*ntheta0*naky*nspec
       write (report_unit, fmt="('Number of meshpoints:    ',i14)") nmesh
    end if
    
    write (report_unit, fmt="(T12,' ntgrid=',i12)") ntgrid
    write (report_unit, fmt="(T9,'2*ntgrid+1=',i12)") 2*ntgrid+1
    write (report_unit, fmt="(T12,'nlambda=',i12)") nlambda
    write (report_unit, fmt="(T12,' negrid=',i12)") negrid
    write (report_unit, fmt="(T12,'ntheta0=',i12)") ntheta0
    write (report_unit, fmt="(T12,'     nx=',i12)") nx
    write (report_unit, fmt="(T12,'   naky=',i12)") naky
    write (report_unit, fmt="(T12,'     ny=',i12)") ny
    write (report_unit, fmt="(T12,'  nspec=',i12,/)") nspec

    call nprocs (nmesh, report_unit)
    if (nonlin) then
       write (report_unit, fmt="(/'Nonlinear run => consider #proc sweetspots for xxf+yxf objects!')") 
       call nprocs_xxf(report_unit)
       call nprocs_yxf(report_unit)
    endif

    if(use_le_layout) then
       write (report_unit, fmt="(/'Collisions using le_lo')") 
       call nprocs_le(report_unit)
    else
       select case (collision_model_switch)
       case (collision_model_full)
          write (report_unit, fmt="(/'Collisions using lz_lo')") 
          call nprocs_lz(report_unit)
          write (report_unit, fmt="(/'Collisions using e_lo')") 
          call nprocs_e(report_unit)
       case(collision_model_lorentz,collision_model_lorentz_test)
          write (report_unit, fmt="(/'Collisions using lz_lo')") 
          call nprocs_lz(report_unit)
       case (collision_model_ediffuse)
          write (report_unit, fmt="(/'Collisions using e_lo')") 
          call nprocs_e(report_unit)
       end select
    end if

    call write_separator(report_unit)
    call check_species(report_unit,beta,tite,alne,dbetadrho_spec)
    call write_separator(report_unit)
    call check_theta_grid(report_unit,alne,dbetadrho_spec)
    call write_separator(report_unit)

    if (coll_on) then 
       call check_collisions(report_unit) 
    else
       write (report_unit, fmt="('All collisionality parameters (vnewk) are zero.')")
       write (report_unit, fmt="('No collision operator will be used.')")
    end if

    call write_separator(report_unit)
    call check_fields(report_unit) 
    call write_separator(report_unit)
    call check_init_g(report_unit)
    call write_separator(report_unit)
    call check_nonlinear_terms(report_unit,delt_adj)

    if (.not. nonlin) then
       write (report_unit, *)
       write (report_unit, fmt="('This is a linear calculation.')")
       write (report_unit, *)
       write (report_unit, fmt="('The time step (code_dt) = ',e11.4)") code_dt
       write (report_unit, fmt="('The time step (user_dt) = ',e11.4)") user_dt
       write (report_unit, fmt="('The maximum number of time steps is ',i7)") nstep
       if (use_old_diagnostics) then
          local_exit_when_converged = exit_when_converged
          local_omegatol = omegatol
       else
          local_exit_when_converged = gnostics%exit_when_converged
          local_omegatol = gnostics%omegatol
       end if

       if (local_exit_when_converged) then
          write (report_unit, fmt="('When the frequencies for each k have converged, the run will stop.')")
          write (report_unit, fmt="('The convergence has to be better than one part in ',e11.4)") 1./local_omegatol
       end if

       if (wstar_units) then
          write (report_unit, *)
          write (report_unit, fmt="('The timestep for each ky is scaled by a factor of 1/ky.')")
       end if
    end if

    if (use_old_diagnostics) then
       local_nwrite = nwrite
       local_make_movie = make_movie
       local_nmovie = nmovie
       local_save_for_restart = save_for_restart
       local_nsave = nsave
    else
       local_nwrite = gnostics%nwrite
       local_make_movie = gnostics%make_movie
       local_nmovie = gnostics%nmovie
       local_save_for_restart = gnostics%save_for_restart
       local_nsave = gnostics%nsave
    end if

    write (report_unit, *)
    write (report_unit, fmt="('Data will be written to ',a,' every ',i4,' timesteps.')") trim(run_name)//'.out.nc', local_nwrite
    write (report_unit, *)

    if (local_make_movie) then
       write (report_unit, *)
       write (report_unit, fmt="('Movie data will be written to runname.movie.nc every ',i4,' timesteps.')") local_nmovie
       write (report_unit, *)
    end if

    if (local_save_for_restart .and. local_nsave > 0) then
       write (report_unit, *)
       write (report_unit, fmt="('Restart data will be written every ',i4,' timesteps.')") local_nsave
       write (report_unit, *)
    end if

    call write_separator(report_unit)
    call check_kt_grids(report_unit)
    call write_separator(report_unit)
    call check_dist_fn(report_unit)
    call write_separator(report_unit)
    call check_antenna(report_unit)
    call check_hyper(report_unit)
    call write_separator(report_unit)
    call check_run_parameters(report_unit)

    ! diagnostic controls:
    if (use_old_diagnostics) call check_gs2_diagnostics(report_unit)
    call close_output_file (report_unit)
  end subroutine report