analytic_divertor_equilibrium_m.f90 Source File


This file depends on

sourcefile~~analytic_divertor_equilibrium_m.f90~~EfferentGraph sourcefile~analytic_divertor_equilibrium_m.f90 analytic_divertor_equilibrium_m.f90 sourcefile~divertor_equilibrium_m.f90 divertor_equilibrium_m.f90 sourcefile~analytic_divertor_equilibrium_m.f90->sourcefile~divertor_equilibrium_m.f90 sourcefile~precision_m.f90 precision_m.f90 sourcefile~analytic_divertor_equilibrium_m.f90->sourcefile~precision_m.f90 sourcefile~divertor_equilibrium_m.f90->sourcefile~precision_m.f90 sourcefile~descriptors_m.f90 descriptors_m.f90 sourcefile~divertor_equilibrium_m.f90->sourcefile~descriptors_m.f90 sourcefile~equilibrium_m.f90 equilibrium_m.f90 sourcefile~divertor_equilibrium_m.f90->sourcefile~equilibrium_m.f90 sourcefile~error_handling_m.f90 error_handling_m.f90 sourcefile~divertor_equilibrium_m.f90->sourcefile~error_handling_m.f90 sourcefile~polygon_m.f90 polygon_m.f90 sourcefile~divertor_equilibrium_m.f90->sourcefile~polygon_m.f90 sourcefile~screen_io_m.f90 screen_io_m.f90 sourcefile~divertor_equilibrium_m.f90->sourcefile~screen_io_m.f90 sourcefile~status_codes_m.f90 status_codes_m.f90 sourcefile~divertor_equilibrium_m.f90->sourcefile~status_codes_m.f90 sourcefile~descriptors_m.f90->sourcefile~screen_io_m.f90 sourcefile~equilibrium_m.f90->sourcefile~precision_m.f90 sourcefile~error_handling_m.f90->sourcefile~precision_m.f90 sourcefile~error_handling_m.f90->sourcefile~screen_io_m.f90 sourcefile~error_handling_m.f90->sourcefile~status_codes_m.f90 sourcefile~comm_handling_m.f90 comm_handling_m.f90 sourcefile~error_handling_m.f90->sourcefile~comm_handling_m.f90 sourcefile~polygon_m.f90->sourcefile~precision_m.f90 sourcefile~polygon_m.f90->sourcefile~descriptors_m.f90 sourcefile~polygon_m.f90->sourcefile~screen_io_m.f90 sourcefile~polygon_m.f90->sourcefile~comm_handling_m.f90 sourcefile~screen_io_m.f90->sourcefile~precision_m.f90

Files dependent on this one

sourcefile~~analytic_divertor_equilibrium_m.f90~~AfferentGraph sourcefile~analytic_divertor_equilibrium_m.f90 analytic_divertor_equilibrium_m.f90 sourcefile~carthy_equilibrium_m.f90 carthy_equilibrium_m.f90 sourcefile~carthy_equilibrium_m.f90->sourcefile~analytic_divertor_equilibrium_m.f90 sourcefile~cerfons_equilibrium_m.f90 cerfons_equilibrium_m.f90 sourcefile~cerfons_equilibrium_m.f90->sourcefile~analytic_divertor_equilibrium_m.f90 sourcefile~equilibrium_factory_m.f90 equilibrium_factory_m.f90 sourcefile~equilibrium_factory_m.f90->sourcefile~carthy_equilibrium_m.f90 sourcefile~equilibrium_factory_m.f90->sourcefile~cerfons_equilibrium_m.f90 sourcefile~benchmark_helmholtz_solvers.f90 benchmark_helmholtz_solvers.f90 sourcefile~benchmark_helmholtz_solvers.f90->sourcefile~equilibrium_factory_m.f90 sourcefile~diagnose_poincare.f90 diagnose_poincare.f90 sourcefile~diagnose_poincare.f90->sourcefile~equilibrium_factory_m.f90 sourcefile~test_diffusion.f90 test_diffusion.f90 sourcefile~test_diffusion.f90->sourcefile~equilibrium_factory_m.f90

Source Code

module analytic_divertor_equilibrium_m
    use divertor_equilibrium_m, only: divertor_equilibrium_t, max_polygon_N_pts
    use precision_m, only : FP

    implicit none

    type, abstract, extends(divertor_equilibrium_t) :: &
        analytic_divertor_equilibrium_t
        real(FP) :: R0
        !! Non-normalised magentic axis radius
        real(FP) :: Z0
        !! Non-normalised magentic axis vertical position
        real(FP) :: RX
        !! Normalised primary x-point radius
        real(FP) :: ZX
        !! Normalised primary x-point vertical position
        real(FP) :: rhomin_privflux
    contains
        procedure, public :: check_privflux_regions
    end type analytic_divertor_equilibrium_t

contains

    subroutine check_privflux_regions(self, x, y, phi, local_rhomin, &
                        local_rhomin_exists, local_rhomax, local_rhomax_exists)
        !! For simple X-point geometries, can set any points below the X-point
        !! as private flux if not in the SOL
        class(analytic_divertor_equilibrium_t), intent(in) :: self
        real(FP), intent(in) :: x, y, phi
        real(FP), intent(inout) :: local_rhomax, local_rhomin
        logical, intent(out) :: local_rhomin_exists, local_rhomax_exists

        local_rhomin_exists = .false.
        local_rhomax_exists = .false.

        if ((y <= self%ZX) .neqv. self%flip_Z) then
            local_rhomin = max(self%rhomin_privflux, local_rhomin)
            local_rhomin_exists = .true.
        endif

    end subroutine check_privflux_regions

end module analytic_divertor_equilibrium_m