do_write_verr Subroutine

private subroutine do_write_verr()

Write some error estimates.

Error estimates are calculated by: - comparing standard integral with less-accurate integral - monitoring amplitudes of legendre polynomial coefficients

Only writes to text file, requires lpc_unit to be open

Arguments

None

Contents

Source Code


Source Code

  subroutine do_write_verr
    use mp, only: proc0
    use le_grids, only: grid_has_trapped_particles
    use fields_arrays, only: phinew, bparnew
    use gs2_time, only: user_time
    use collisions, only: vnmult
    use species, only: spec
    use diagnostics_velocity_space, only: get_verr, get_gtran
    implicit none
    real, dimension(5,2) :: errest
    integer, dimension (5,3) :: erridx
    real :: geavg, glavg, gtavg

    errest = 0.0; erridx = 0

    ! error estimate obtained by comparing standard integral with less-accurate integral
    call get_verr (errest, erridx, phinew, bparnew)

    ! error estimate based on monitoring amplitudes of legendre polynomial coefficients
    call get_gtran (geavg, glavg, gtavg, phinew, bparnew)

    if (proc0) then
       ! write error estimates to .nc file          
       !          call nc_loop_vres (nout, errest_by_mode, lpcoef_by_mode)

       ! write error estimates for ion dist. fn. at outboard midplane with ik=it=1 to ascii files
       if (write_ascii) then
          if (grid_has_trapped_particles()) then
             write(lpc_unit,"(4(1x,e13.6))") user_time, geavg, glavg, gtavg
          else
             write(lpc_unit,"(3(1x,e13.6))") user_time, geavg, glavg
          end if
          write(res_unit,"(8(1x,e13.6))") user_time, errest(1,2), errest(2,2), errest(3,2), &
               errest(4,2), errest(5,2), vnmult(1)*spec(1)%vnewk, vnmult(2)*spec(1)%vnewk
          if (write_max_verr) then
             write(res_unit2,"(3(i8),(1x,e13.6),3(i8),(1x,e13.6),3(i8),(1x,e13.6),3(i8),(1x,e13.6),3(i8),(1x,e13.6))") &
                  erridx(1,1), erridx(1,2), erridx(1,3), errest(1,1), &
                  erridx(2,1), erridx(2,2), erridx(2,3), errest(2,1), &
                  erridx(3,1), erridx(3,2), erridx(3,3), errest(3,1), &
                  erridx(4,1), erridx(4,2), erridx(4,3), errest(4,1), &
                  erridx(5,1), erridx(5,2), erridx(5,3), errest(5,1)
          end if
       end if
    end if
  end subroutine do_write_verr