get_cross_phase Subroutine

private subroutine get_cross_phase(phase_tot, phase_theta)

This is a highly simplified synthetic diagnostic which calculates the cross phase between the electron density and the perpendicular electron temperature for comparisons with DIII-D.
Returns the value of the cross-phase at the outboard midplane and integrated over all v and x. We can generalize this routine to other fields at some point, but for now this is just a skeleton for a more realistic synthetic diagnostic.

Arguments

Type IntentOptional Attributes Name
real, intent(out) :: phase_tot
real, intent(out) :: phase_theta

Contents

Source Code


Source Code

  subroutine get_cross_phase (phase_tot, phase_theta)
    use species, only: nspec, spec, is_electron_species
    use kt_grids, only: ntheta0, naky
    use theta_grid, only: ntgrid
    use dist_fn, only: getemoms
    use fields_arrays, only: phinew, bparnew
    implicit none
    real, intent (out) :: phase_tot, phase_theta
    complex, dimension (:,:,:,:), allocatable :: ntot, tperp
    complex, dimension (ntheta0, naky) :: nTp_by_mode
    complex :: nTp
    !    real, dimension (ntheta0, naky) :: n2_by_mode, T2_by_mode
    !    real :: n2, T2

    integer ::  is

    allocate ( ntot(-ntgrid:ntgrid,ntheta0,naky,nspec))
    allocate (tperp(-ntgrid:ntgrid,ntheta0,naky,nspec))

    call getemoms (phinew, bparnew, ntot, tperp)

    do is = 1,nspec
       if (is_electron_species(spec(is))) then
          ! get cross_phase at outboard midplane
          call get_vol_int (ntot(0,:,:,is), tperp(0,:,:,is), nTp, nTp_by_mode)
          phase_theta = atan2(aimag(nTp),real(nTp))!/sqrt(n2*T2)
          ! get integrated cross_phase 
          call get_vol_int (ntot(:,:,:,is), tperp(:,:,:,is), nTp, nTp_by_mode)
          phase_tot = atan2(aimag(nTp),real(nTp))!/sqrt(n2*T2)
       end if
    end do

    deallocate (ntot, tperp)

  end subroutine get_cross_phase