Line data Source code
1 : ! Include shortname defintions, so that the F77 code does not have to be modified to 2 : ! reference the CARMA structure. 3 : #include "carma_globaer.h" 4 : 5 : !! This routine calculates the total amount of condensate associated with each gas. 6 : !! 7 : !! @author Chuck Bardeen 8 : !! @version Nov-2009 9 2837406876 : subroutine totalcondensate(carma, cstate, iz, total_ice, total_liquid, rc) 10 : 11 : ! types 12 : use carma_precision_mod 13 : use carma_enums_mod 14 : use carma_constants_mod 15 : use carma_types_mod 16 : use carmastate_mod 17 : use carma_mod 18 : 19 : implicit none 20 : 21 : type(carma_type), intent(in) :: carma !! the carma object 22 : type(carmastate_type), intent(inout) :: cstate !! the carma state object 23 : integer, intent(in) :: iz !! z index 24 : real(kind=f), intent(out) :: total_ice(NGAS) !! total ice at the start of substep 25 : real(kind=f), intent(out) :: total_liquid(NGAS) !! total liquid at the start of substep 26 : integer, intent(inout) :: rc !! return code, negative indicates failure 27 : 28 : ! Local declarations 29 : integer :: igroup ! group index 30 : integer :: icore ! core index 31 : integer :: igas ! gas index 32 : integer :: ibin ! bin index 33 : integer :: ielem ! element index 34 : integer :: i 35 : real(kind=f) :: coremass 36 : real(kind=f) :: volatilemass 37 : 38 : 39 : ! Initialize local variables for keeping track of gas changes due 40 : ! to nucleation and growth in each particle group. 41 8512220628 : total_ice(:) = 0._f 42 8512220628 : total_liquid(:) = 0._f 43 : 44 : ! Iterate over each particle type and total up that ones that interact 45 : ! with the gases. 46 : ! 47 : ! This code assumes that all changes in condensate are associated with 48 : ! growth in a particular gas. This doesn't handle all possible changes 49 : ! associated with nucleation, if the group do not also participate in 50 : ! growth. 51 8512220628 : do igroup = 1,NGROUP 52 : 53 5674813752 : ielem = ienconc(igroup) ! element of particle number concentration 54 : 55 5674813752 : igas = igrowgas(ielem) ! condensing gas 56 : 57 8512220628 : if ((itype(ielem) == I_VOLATILE) .and. (igas /= 0)) then 58 : 59 >11917*10^7 : do ibin = 1, NBIN 60 : 61 : ! If this group has core masses, then determine the involatile component. 62 >11349*10^7 : coremass = 0._f 63 : 64 >39723*10^7 : do i = 1, ncore(igroup) 65 >28374*10^7 : icore = icorelem(i, igroup) 66 >39723*10^7 : coremass = coremass + pc(iz, ibin, icore) 67 : end do 68 : 69 >11349*10^7 : volatilemass = (pc(iz, ibin, ielem) * rmass(ibin, igroup)) - coremass 70 : 71 : ! There seem to be times when the coremass becomes larger than the total 72 : ! mass. This shouldn't happen, but check for it here. 73 : ! 74 : ! NOTE: This can be caused by advection in the parent model or sedimentation 75 : ! in this model. 76 >11917*10^7 : if (volatilemass > 0._f) then 77 >11257*10^7 : if (is_grp_ice(igroup)) then 78 0 : total_ice(igas) = total_ice(igas) + volatilemass 79 : else 80 >11257*10^7 : total_liquid(igas) = total_liquid(igas) + volatilemass 81 : end if 82 : end if 83 : end do 84 : end if 85 : end do 86 : 87 2837406876 : return 88 2837406876 : end