make_passive_grid Subroutine

public subroutine make_passive_grid()

Makes [[libfida:pass_grid] from user defined inputs

Arguments

None

Calls

proc~~make_passive_grid~~CallsGraph proc~make_passive_grid make_passive_grid proc~get_plasma_extrema get_plasma_extrema proc~make_passive_grid->proc~get_plasma_extrema proc~in_plasma in_plasma proc~get_plasma_extrema->proc~in_plasma proc~cyl_to_uvw cyl_to_uvw proc~in_plasma->proc~cyl_to_uvw interface~interpol_coeff interpol_coeff proc~in_plasma->interface~interpol_coeff proc~xyz_to_uvw xyz_to_uvw proc~in_plasma->proc~xyz_to_uvw proc~cyl_interpol3d_coeff cyl_interpol3D_coeff interface~interpol_coeff->proc~cyl_interpol3d_coeff proc~interpol1d_coeff interpol1D_coeff interface~interpol_coeff->proc~interpol1d_coeff proc~interpol2d_coeff_arr interpol2D_coeff_arr interface~interpol_coeff->proc~interpol2d_coeff_arr proc~interpol2d_coeff interpol2D_coeff interface~interpol_coeff->proc~interpol2d_coeff proc~cyl_interpol3d_coeff_arr cyl_interpol3D_coeff_arr interface~interpol_coeff->proc~cyl_interpol3d_coeff_arr proc~interpol1d_coeff_arr interpol1D_coeff_arr interface~interpol_coeff->proc~interpol1d_coeff_arr proc~cyl_interpol3d_coeff->proc~interpol2d_coeff proc~interpol2d_coeff_arr->proc~interpol2d_coeff proc~cyl_interpol3d_coeff_arr->proc~cyl_interpol3d_coeff proc~cyl_interpol3d_coeff_arr->proc~interpol2d_coeff proc~interpol1d_coeff_arr->proc~interpol1d_coeff

Called by

proc~~make_passive_grid~~CalledByGraph proc~make_passive_grid make_passive_grid proc~make_diagnostic_grids make_diagnostic_grids proc~make_diagnostic_grids->proc~make_passive_grid program~fidasim fidasim program~fidasim->proc~make_diagnostic_grids

Contents

Source Code


Source Code

subroutine make_passive_grid
    !+ Makes [[libfida:pass_grid] from user defined inputs
    real(Float64), dimension(3,spec_chords%nchan+npa_chords%nchan) :: r0_arr, v0_arr
    real(Float64), dimension(2,spec_chords%nchan+npa_chords%nchan) :: xy_enter, xy_exit
    logical, dimension(8+2*(spec_chords%nchan+npa_chords%nchan))   :: yle, ygt
    logical, dimension(spec_chords%nchan+npa_chords%nchan)         :: skip
    real(Float64), dimension(:), allocatable  :: xarr, yarr
    real(Float64), dimension(3,8) :: vertices_xyz, vertices_uvw
    real(Float64), dimension(2,3) :: extrema
    real(Float64), dimension(3,3) :: xyz_axis
    real(Float64), dimension(3)   :: r0, vi
    real(Float64), dimension(8)   :: xarr_beam_grid, yarr_beam_grid
    real(Float64) :: xmin, ymin, xmax, ymax, zmin, zmax, max_length
    real(Float64) :: dlength = 3.0 !cm
    integer :: i, iin, iout, dim_le, dim_gt, dim
    logical :: inp, phi_pos

    !! Get beam grid boundaries
    ! Convert vertices_xyz to vertices_uvw
    ! Note: vertices_xyz has the following coordinate definitions:
    ! 111 = 1, 222 = 2, 112 = 3, 211 = 4, 121 = 5, 212 = 6, 122 = 7, 221 = 8
    xmin = beam_grid%xmin ; xmax = beam_grid%xmax
    ymin = beam_grid%ymin ; ymax = beam_grid%ymax
    zmin = beam_grid%zmin ; zmax = beam_grid%zmax

    !Initialize minimum and maximum vertices
    vertices_xyz(1,1) = xmin   ; vertices_xyz(2,1) = ymin   ; vertices_xyz(3,1) = zmin
    vertices_xyz(1,2) = xmax   ; vertices_xyz(2,2) = ymax   ; vertices_xyz(3,2) = zmax

    !Initialize
    vertices_xyz(:,3) = vertices_xyz(:,1)
    vertices_xyz(:,4) = vertices_xyz(:,1)
    vertices_xyz(:,5) = vertices_xyz(:,1)
    !Update
    vertices_xyz(3,3) = zmax
    vertices_xyz(1,4) = xmax
    vertices_xyz(2,5) = ymax
    !Initialize
    vertices_xyz(:,7) = vertices_xyz(:,2)
    vertices_xyz(:,8) = vertices_xyz(:,2)
    vertices_xyz(:,6) = vertices_xyz(:,2)
    !Update
    vertices_xyz(1,7) = xmin
    vertices_xyz(3,8) = zmin
    vertices_xyz(2,6) = ymin

    do i=1, 8
        call xyz_to_uvw(vertices_xyz(:,i),vertices_uvw(:,i))
        xarr_beam_grid(i) = vertices_uvw(1,i)
        yarr_beam_grid(i) = vertices_uvw(2,i)
    enddo

    !! Next consider passive diagnostic extrema relative to the plasma
    if((inputs%calc_pfida.gt.0).and.(inputs%calc_pnpa.gt.0)) then
        do i=1,(spec_chords%nchan)
            r0_arr(:,i) = spec_chords%los(i)%lens_uvw
            v0_arr(:,i) = spec_chords%los(i)%axis_uvw
        enddo
        do i=1, npa_chords%nchan
            call xyz_to_uvw(npa_chords%det(i)%detector%origin, r0_arr(:,i+spec_chords%nchan))
            xyz_axis = npa_chords%det(i)%detector%basis
            v0_arr(:,i+spec_chords%nchan) = matmul(beam_grid%basis, xyz_axis(:,3))
        enddo
    else if(inputs%calc_pfida.gt.0) then
        do i=1, spec_chords%nchan
            r0_arr(:,i) = spec_chords%los(i)%lens_uvw
            v0_arr(:,i) = spec_chords%los(i)%axis_uvw
        enddo
    else !pnpa>=1 case
        do i=1, npa_chords%nchan
            call xyz_to_uvw(npa_chords%det(i)%detector%origin, r0_arr(:,i))
            xyz_axis = npa_chords%det(i)%detector%basis
            v0_arr(:,i) = matmul(beam_grid%basis, xyz_axis(:,3))
        enddo
    endif

    call get_plasma_extrema(r0_arr,v0_arr,extrema,xarr_beam_grid,yarr_beam_grid)

    !! Store the passive neutral grid
    pass_grid%dr = inter_grid%dr
    pass_grid%dz = inter_grid%dz
    pass_grid%nr = inter_grid%nr
    pass_grid%nz = inter_grid%nz
    allocate(pass_grid%r(pass_grid%nr), pass_grid%z(pass_grid%nz))
    pass_grid%r = inter_grid%r
    pass_grid%z = inter_grid%z
    pass_grid%da = pass_grid%dr*pass_grid%dz

    pass_grid%dphi = 2*pi/100 !TODO: make this user input
    pass_grid%nphi = int(ceiling((extrema(2,3)-extrema(1,3))/pass_grid%dphi))

    allocate(pass_grid%phi(pass_grid%nphi))
    do i=1, pass_grid%nphi
        pass_grid%phi(i) = extrema(1,3) + (i-1)*pass_grid%dphi
    enddo

    pass_grid%dv = pass_grid%dr*pass_grid%dphi*pass_grid%dz
    pass_grid%dims = [pass_grid%nr, pass_grid%nz, pass_grid%nphi]

    pass_grid%ntrack = pass_grid%nr+pass_grid%nz+pass_grid%nphi
    pass_grid%ngrid  = pass_grid%nr*pass_grid%nz*pass_grid%nphi

    if(inputs%verbose.ge.1) then
        write(*,'(a)') "---- Passive grid settings ----"
        write(*,'(T2,"Nr: ", i3)') pass_grid%nr
        write(*,'(T2,"Nz: ", i3)') pass_grid%nz
        write(*,'(T2,"Nphi: ", i3)') pass_grid%nphi
        write(*,'(T2,"R  range = [",f6.2,",",f6.2,"]")') &
              pass_grid%r(1),pass_grid%r(pass_grid%nr)
        write(*,'(T2,"Z  range = [",f7.2,",",f6.2,"]")') &
              pass_grid%z(1),pass_grid%z(pass_grid%nz)
        write(*,'(T2,"Phi  range = [",f5.2,",",f5.2,"]")') &
              pass_grid%phi(1),pass_grid%phi(pass_grid%nphi)
        write(*,'(T2,"dA: ", f5.2," [cm^3]")') pass_grid%da
        write(*,*) ''
    endif

end subroutine make_passive_grid