uEMEP_set_pollutant_loop Subroutine

public subroutine uEMEP_set_pollutant_loop()

Uses

  • proc~~uemep_set_pollutant_loop~~UsesGraph proc~uemep_set_pollutant_loop uEMEP_set_pollutant_loop module~uemep_definitions uEMEP_definitions proc~uemep_set_pollutant_loop->module~uemep_definitions

Arguments

None

Called by

proc~~uemep_set_pollutant_loop~~CalledByGraph proc~uemep_set_pollutant_loop uEMEP_set_pollutant_loop program~uemep uEMEP program~uemep->proc~uemep_set_pollutant_loop

Source Code

    subroutine uEMEP_set_pollutant_loop

        use uEMEP_definitions

        implicit none

        integer p_loop

        !Set the pollutant index loops after reading in pollutant_index
        !Remove the sand and salt PM2.5, not necessary. Fixed ratio if needed n_pollutant_loop=6
        if (pollutant_index.eq.all_sand_salt_nc_index) then
            n_emep_pollutant_loop=3
            !if (use_GNFR19_emissions_from_EMEP_flag) n_emep_pollutant_loop=4 !Include exhaust
            n_pollutant_loop=6
            pollutant_loop_index(1)=nox_nc_index
            pollutant_loop_index(2)=pm25_nc_index
            pollutant_loop_index(3)=pm10_nc_index
            pollutant_loop_index(4)=pmex_nc_index
            pollutant_loop_index(5)=pm10_sand_nc_index
            pollutant_loop_index(6)=pm10_salt_nc_index
            pollutant_loop_index(7)=pm25_sand_nc_index
            pollutant_loop_index(8)=pm25_salt_nc_index
            pollutant_loop_back_index(nox_nc_index)=1
            pollutant_loop_back_index(pm25_nc_index)=2
            pollutant_loop_back_index(pm10_nc_index)=3
            pollutant_loop_back_index(pmex_nc_index)=4
            pollutant_loop_back_index(pm10_sand_nc_index)=5
            pollutant_loop_back_index(pm10_salt_nc_index)=6
            pollutant_loop_back_index(pm25_sand_nc_index)=7
            pollutant_loop_back_index(pm25_salt_nc_index)=8
        elseif (pollutant_index.eq.all_salt_nc_index) then
            n_emep_pollutant_loop=3
            !if (use_GNFR19_emissions_from_EMEP_flag) n_emep_pollutant_loop=4 !Include exhaust
            n_pollutant_loop=5
            pollutant_loop_index(1)=nox_nc_index
            pollutant_loop_index(2)=pm25_nc_index
            pollutant_loop_index(3)=pm10_nc_index
            pollutant_loop_index(4)=pmex_nc_index
            pollutant_loop_index(5)=pm10_salt_nc_index
            pollutant_loop_index(6)=pm25_salt_nc_index
            pollutant_loop_back_index(nox_nc_index)=1
            pollutant_loop_back_index(pm25_nc_index)=2
            pollutant_loop_back_index(pm10_nc_index)=3
            pollutant_loop_back_index(pmex_nc_index)=4
            pollutant_loop_back_index(pm10_salt_nc_index)=5
            pollutant_loop_back_index(pm25_salt_nc_index)=6
        elseif (pollutant_index.eq.all_sand_nc_index) then
            n_emep_pollutant_loop=3
            !if (use_GNFR19_emissions_from_EMEP_flag) n_emep_pollutant_loop=4 !Include exhaust
            n_pollutant_loop=5
            pollutant_loop_index(1)=nox_nc_index
            pollutant_loop_index(2)=pm25_nc_index
            pollutant_loop_index(3)=pm10_nc_index
            pollutant_loop_index(4)=pmex_nc_index
            pollutant_loop_index(5)=pm10_sand_nc_index
            pollutant_loop_index(6)=pm25_sand_nc_index
            pollutant_loop_back_index(nox_nc_index)=1
            pollutant_loop_back_index(pm25_nc_index)=2
            pollutant_loop_back_index(pm10_nc_index)=3
            pollutant_loop_back_index(pmex_nc_index)=4
            pollutant_loop_back_index(pm10_sand_nc_index)=5
            pollutant_loop_back_index(pm25_sand_nc_index)=6
        elseif (pollutant_index.eq.all_nc_index) then
            n_emep_pollutant_loop=3
            !if (use_GNFR19_emissions_from_EMEP_flag) n_emep_pollutant_loop=4 !Include exhaust
            n_pollutant_loop=4
            pollutant_loop_index(1)=nox_nc_index
            pollutant_loop_index(2)=pm25_nc_index
            pollutant_loop_index(3)=pm10_nc_index
            pollutant_loop_index(4)=pmex_nc_index
            pollutant_loop_back_index(nox_nc_index)=1
            pollutant_loop_back_index(pm25_nc_index)=2
            pollutant_loop_back_index(pm10_nc_index)=3
            pollutant_loop_back_index(pmex_nc_index)=4
        elseif (pollutant_index.eq.all_totals_nc_index) then
            n_emep_pollutant_loop=3
            n_pollutant_loop=3
            pollutant_loop_index(1)=nox_nc_index
            pollutant_loop_index(2)=pm25_nc_index
            pollutant_loop_index(3)=pm10_nc_index
            pollutant_loop_back_index(nox_nc_index)=1
            pollutant_loop_back_index(pm25_nc_index)=2
            pollutant_loop_back_index(pm10_nc_index)=3
        elseif (pollutant_index.eq.aaqd_totals_nc_index) then
            n_emep_pollutant_loop=6
            n_pollutant_loop=6
            pollutant_loop_index(1)=nox_nc_index
            pollutant_loop_index(2)=pm25_nc_index
            pollutant_loop_index(3)=pm10_nc_index
            pollutant_loop_index(4)=co_nc_index
            pollutant_loop_index(5)=bap_nc_index
            pollutant_loop_index(6)=c6h6_nc_index
            pollutant_loop_back_index(nox_nc_index)=1
            pollutant_loop_back_index(pm25_nc_index)=2
            pollutant_loop_back_index(pm10_nc_index)=3
            pollutant_loop_back_index(co_nc_index)=4
            pollutant_loop_back_index(bap_nc_index)=5
            pollutant_loop_back_index(c6h6_nc_index)=6
            extract_benzene_from_voc_emissions=.true.
        elseif (pollutant_index.eq.gp_totals_nc_index) then
            n_emep_pollutant_loop=4
            n_pollutant_loop=4
            pollutant_loop_index(1)=nox_nc_index
            pollutant_loop_index(2)=pm25_nc_index
            pollutant_loop_index(3)=pm10_nc_index
            pollutant_loop_index(4)=co_nc_index
            pollutant_loop_back_index(nox_nc_index)=1
            pollutant_loop_back_index(pm25_nc_index)=2
            pollutant_loop_back_index(pm10_nc_index)=3
            pollutant_loop_back_index(co_nc_index)=4
        elseif (pollutant_index.eq.op_totals_nc_index) then
            n_emep_pollutant_loop=3
            n_pollutant_loop=3
            pollutant_loop_index(1)=nox_nc_index
            pollutant_loop_index(2)=pm25_nc_index
            pollutant_loop_index(3)=pm10_nc_index
            pollutant_loop_back_index(nox_nc_index)=1
            pollutant_loop_back_index(pm25_nc_index)=2
            pollutant_loop_back_index(pm10_nc_index)=3
        elseif (pollutant_index.eq.pm_nc_index) then
            n_emep_pollutant_loop=2
            n_pollutant_loop=3
            pollutant_loop_index(1)=pm25_nc_index
            pollutant_loop_index(2)=pm10_nc_index
            pollutant_loop_index(3)=pmex_nc_index
            pollutant_loop_back_index(pm25_nc_index)=1
            pollutant_loop_back_index(pm10_nc_index)=2
            pollutant_loop_back_index(pmex_nc_index)=3
        elseif (pollutant_index.eq.nh3_nc_index) then
            n_emep_pollutant_loop=1
            n_pollutant_loop=1
            pollutant_loop_index(1)=nh3_nc_index
            pollutant_loop_index(2)=nh4_nc_index
            pollutant_loop_back_index(nh3_nc_index)=1
            pollutant_loop_back_index(nh4_nc_index)=2
        else
            n_emep_pollutant_loop=1
            n_pollutant_loop=1
            pollutant_loop_index(1)=pollutant_index
            pollutant_loop_back_index(pollutant_index)=1
            !write(*,*) pollutant_loop_index(1),pollutant_index

        endif

        !Set indexing for additional compounds. Only used when reading in EMEP data
        do p_loop=1,n_pollutant_loop
            if (pollutant_loop_index(p_loop).eq.nox_nc_index) then
                n_pollutant_compound_loop(p_loop)=3
                pollutant_compound_loop_index(p_loop,1)=nox_nc_index
                pollutant_compound_loop_index(p_loop,2)=no2_nc_index
                pollutant_compound_loop_index(p_loop,3)=o3_nc_index
                !Add addition values to be read and saved
                if (save_EMEP_somo35) then
                    n_pollutant_compound_loop(p_loop)=n_pollutant_compound_loop(p_loop)+1
                    pollutant_compound_loop_index(p_loop,n_pollutant_compound_loop(p_loop))=somo35_nc_index
                endif
                if (save_EMEP_o3max) then
                    n_pollutant_compound_loop(p_loop)=n_pollutant_compound_loop(p_loop)+1
                    pollutant_compound_loop_index(p_loop,n_pollutant_compound_loop(p_loop))=o3max_nc_index
                endif
                if (save_EMEP_o3_26th) then
                    n_pollutant_compound_loop(p_loop)=n_pollutant_compound_loop(p_loop)+1
                    pollutant_compound_loop_index(p_loop,n_pollutant_compound_loop(p_loop))=o3_26th_nc_index
                endif
                if (save_EMEP_so2) then
                    n_pollutant_compound_loop(p_loop)=n_pollutant_compound_loop(p_loop)+1
                    pollutant_compound_loop_index(p_loop,n_pollutant_compound_loop(p_loop))=so2_nc_index
                endif
            elseif (pollutant_loop_index(p_loop).eq.nh3_nc_index) then
                n_pollutant_compound_loop(p_loop)=2
                pollutant_compound_loop_index(p_loop,1)=nh3_nc_index
                pollutant_compound_loop_index(p_loop,2)=nh4_nc_index
            elseif (pollutant_loop_index(p_loop).eq.co_nc_index) then
                n_pollutant_compound_loop(p_loop)=1
                pollutant_compound_loop_index(p_loop,1)=co_nc_index
                if (save_EMEP_comax) then
                    n_pollutant_compound_loop(p_loop)=n_pollutant_compound_loop(p_loop)+1
                    pollutant_compound_loop_index(p_loop,n_pollutant_compound_loop(p_loop))=comax_nc_index
                endif
            else
                n_pollutant_compound_loop(p_loop)=1
                pollutant_compound_loop_index(p_loop,1)=pollutant_loop_index(p_loop)
            endif
        enddo



        write(unit_logfile,'(a,i)') 'Number of pollutants=',n_pollutant_loop
        write(unit_logfile,'(a,i)') 'Number of EMEP pollutants=',n_emep_pollutant_loop

    end subroutine uEMEP_set_pollutant_loop