pest_washp.f90 Source File


This file depends on

sourcefile~~pest_washp.f90~~EfferentGraph sourcefile~pest_washp.f90 pest_washp.f90 sourcefile~constituent_mass_module.f90 constituent_mass_module.f90 sourcefile~pest_washp.f90->sourcefile~constituent_mass_module.f90 sourcefile~hru_module.f90 hru_module.f90 sourcefile~pest_washp.f90->sourcefile~hru_module.f90 sourcefile~output_ls_pesticide_module.f90 output_ls_pesticide_module.f90 sourcefile~pest_washp.f90->sourcefile~output_ls_pesticide_module.f90 sourcefile~pesticide_data_module.f90 pesticide_data_module.f90 sourcefile~pest_washp.f90->sourcefile~pesticide_data_module.f90 sourcefile~plant_module.f90 plant_module.f90 sourcefile~pest_washp.f90->sourcefile~plant_module.f90 sourcefile~soil_module.f90 soil_module.f90 sourcefile~pest_washp.f90->sourcefile~soil_module.f90 sourcefile~carbon_module.f90 carbon_module.f90 sourcefile~soil_module.f90->sourcefile~carbon_module.f90

Source Code

      subroutine pest_washp

!!    ~ ~ ~ PURPOSE ~ ~ ~
!!    this subroutine calculates the amount of pesticide washed off the plant
!!    foliage and onto the soil

      use pesticide_data_module
      use output_ls_pesticide_module
      use hru_module, only :  ihru
      use soil_module
      use constituent_mass_module
      use plant_module
      
      implicit none       
      
      integer :: j = 0    !none          |HRU number
      integer :: k = 0    !none          |pesticide number
      integer :: ipl = 0  !none          |plant number
      integer :: ipest_db = 0!none          |pesticide number from pest.dat
      real :: pest_soil = 0.!kg/ha         |amount of pesticide in soil   

      j = ihru

      if (cs_db%num_pests == 0) return

      do k = 1, cs_db%num_pests
        ipest_db = cs_db%pest_num(k)
        !! adjust foliar pesticide for wash off
        do ipl = 1, pcom(j)%npl
          if (cs_pl(j)%pl_on(ipl)%pest(k) >= 0.0001) then
            if (ipest_db > 0) then
              pest_soil = pestdb(ipest_db)%washoff * cs_pl(j)%pl_on(ipl)%pest(k)
              if (pest_soil > cs_pl(j)%pl_on(ipl)%pest(k)) pest_soil = cs_pl(j)%pl_on(ipl)%pest(k)
              cs_soil(j)%ly(1)%pest(k) = cs_soil(j)%ly(1)%pest(k) + pest_soil
              cs_pl(j)%pl_on(ipl)%pest(k) = cs_pl(j)%pl_on(ipl)%pest(k) - pest_soil
              hpestb_d(j)%pest(k)%wash = pest_soil
            end if
          end if 
        end do
      end do

      return
      end subroutine pest_washp