LCOV - code coverage report
Current view: top level - physics/clubb/src/CLUBB_core - fill_holes.F90 (source / functions) Hit Total Coverage
Test: coverage.info Lines: 64 268 23.9 %
Date: 2025-03-13 18:42:46 Functions: 1 7 14.3 %

          Line data    Source code
       1             : !-----------------------------------------------------------------------
       2             : ! $Id: fill_holes.F90 8738 2018-07-19 19:58:53Z bmg2@uwm.edu $
       3             : !===============================================================================
       4             : module fill_holes
       5             : 
       6             :   implicit none
       7             : 
       8             :   public :: fill_holes_driver, &
       9             :             fill_holes_vertical, &
      10             :             hole_filling_hm_one_lev, &
      11             :             fill_holes_hydromet, &
      12             :             fill_holes_wv, &
      13             :             clip_hydromet_conc_mvr, &
      14             :             setup_stats_indices
      15             : 
      16             :   private ! Set Default Scope
      17             : 
      18             :   contains
      19             : 
      20             :   !=============================================================================
      21    10588320 :   subroutine fill_holes_vertical( nz, ngrdcol, num_draw_pts, threshold, upper_hf_level, &
      22    10588320 :                                   dz, rho_ds, &
      23    10588320 :                                   field )
      24             : 
      25             :     ! Description:
      26             :     !   This subroutine clips values of 'field' that are below 'threshold' as much
      27             :     !   as possible (i.e. "fills holes"), but conserves the total integrated mass
      28             :     !   of 'field'.  This prevents clipping from acting as a spurious source.
      29             :     !
      30             :     !   Mass is conserved by reducing the clipped field everywhere by a constant
      31             :     !   multiplicative coefficient.
      32             :     !
      33             :     !   This subroutine does not guarantee that the clipped field will exceed
      34             :     !   threshold everywhere; blunt clipping is needed for that.
      35             :     !
      36             :     !   The lowest level (k=1) should not be included, as the hole-filling scheme
      37             :     !   should not alter the set value of 'field' at the surface (for momentum
      38             :     !   level variables), or consider the value of 'field' at a level below the
      39             :     !   surface (for thermodynamic level variables).  
      40             :     !
      41             :     !   For momentum level variables only, the hole-filling scheme should not 
      42             :     !   alter the set value of 'field' at the upper boundary level (k=nz).
      43             :     !   So for momemtum level variables, call with upper_hf_level=nz-1, and
      44             :     !   for thermodynamic level variables, call with upper_hf_level=nz.
      45             :     !
      46             :     ! References:
      47             :     !   ``Numerical Methods for Wave Equations in Geophysical Fluid
      48             :     !     Dynamics'', Durran (1999), p. 292.
      49             :     !-----------------------------------------------------------------------
      50             : 
      51             :     use grid_class, only: & 
      52             :         grid ! Type
      53             : 
      54             :     use clubb_precision, only: &
      55             :         core_rknd ! Variable(s)
      56             : 
      57             :     use constants_clubb, only: &
      58             :         eps, &
      59             :         one
      60             : 
      61             :     implicit none
      62             :     
      63             :     ! --------------------- Input variables ---------------------
      64             :     integer, intent(in) :: &
      65             :       nz, &
      66             :       ngrdcol
      67             :     
      68             :     real( kind = core_rknd ), dimension(ngrdcol,nz), intent(in) :: &
      69             :       dz      ! Spacing between thermodynamic grid levels; centered over
      70             :               ! momentum grid levels
      71             :               ! OR
      72             :               ! Spcaing between momentum grid levels; centered over
      73             :               ! thermodynamic grid levels
      74             :                   
      75             :     integer, intent(in) :: & 
      76             :       num_draw_pts, & ! The number of points on either side of the hole;
      77             :                       ! Mass is drawn from these points to fill the hole.  []
      78             :       upper_hf_level  ! Upper grid level of global hole-filling range      []
      79             : 
      80             :     real( kind = core_rknd ), intent(in) :: & 
      81             :       threshold  ! A threshold (e.g. w_tol*w_tol) below which field must not
      82             :                  ! fall                           [Units vary; same as field]
      83             : 
      84             :     real( kind = core_rknd ), dimension(ngrdcol,nz), intent(in) ::  & 
      85             :       rho_ds       ! Dry, static density on thermodynamic or momentum levels    [kg/m^3]
      86             : 
      87             :     ! --------------------- Input/Output variable ---------------------
      88             :     real( kind = core_rknd ), dimension(ngrdcol,nz), intent(inout) :: & 
      89             :       field  ! The field (e.g. wp2) that contains holes [Units same as threshold]
      90             :  
      91             :     ! --------------------- Local Variables ---------------------
      92             :     integer :: & 
      93             :       i,             & ! Loop index for column                           []
      94             :       k,             & ! Loop index for absolute grid level              []
      95             :       j, &
      96             :       k_start,     & ! Lower grid level of local hole-filling range    []
      97             :       k_end          ! Upper grid level of local hole-filling range    []
      98             : 
      99             :     real( kind = core_rknd ), dimension(ngrdcol,nz)  ::  & 
     100    21176640 :       rho_ds_dz,            & ! rho_ds * dz
     101    21176640 :       invrs_denom_integral, & ! Inverse of the integral in the denominator (see description)
     102    21176640 :       field_clipped           ! The raw field (e.g. wp2) that contains no holes
     103             :                               !                          [Units same as threshold]
     104             : 
     105             :     real( kind = core_rknd ) ::  & 
     106             :       field_avg,          & ! Vertical average of field [Units of field]
     107             :       field_clipped_avg,  & ! Vertical average of clipped field [Units of field]
     108             :       mass_fraction         ! Coefficient that multiplies clipped field
     109             :                             ! in order to conserve mass.                      []
     110             : 
     111             :     real( kind = core_rknd ), dimension(ngrdcol) ::  & 
     112    21176640 :       denom_integral_global,  & ! Integral in the denominator for global filling
     113    21176640 :       numer_integral_global,  & ! Integral in the numerator for global filling
     114    21176640 :       field_avg_global,       & ! Vertical average of field [Units of field]
     115    21176640 :       mass_fraction_global      ! Coefficient that multiplies clipped field
     116             :                                 ! in order to conserve mass.                      []
     117             : 
     118             :     logical :: &
     119             :       l_field_below_threshold
     120             : 
     121             :     ! --------------------- Begin Code --------------------- 
     122             : 
     123             :     !$acc enter data create( invrs_denom_integral, field_clipped, denom_integral_global, rho_ds_dz, &
     124             :     !$acc                    numer_integral_global, field_avg_global, mass_fraction_global )
     125             : 
     126    10588320 :     l_field_below_threshold = .false.
     127             : 
     128             :     !$acc parallel loop gang vector collapse(2) default(present) &
     129             :     !$acc          reduction(.or.:l_field_below_threshold)
     130   910595520 :     do k = 1, nz
     131 15038615520 :       do i = 1, ngrdcol
     132 15028027200 :         if ( field(i,k) < threshold ) then
     133   355220211 :           l_field_below_threshold = .true.
     134             :         end if
     135             :       end do
     136             :     end do
     137             :     !$acc end parallel loop
     138             : 
     139             :     ! If all field values are above the specified threshold, no hole filling is required
     140    10588320 :     if ( .not. l_field_below_threshold ) then
     141             :       !$acc exit data delete( invrs_denom_integral, field_clipped, denom_integral_global, rho_ds_dz, &
     142             :       !$acc                   numer_integral_global, field_avg_global, mass_fraction_global )
     143             :       return
     144             :     end if
     145             : 
     146             :     !$acc parallel loop gang vector collapse(2) default(present)
     147    95935838 :     do k = 1, nz
     148  1584830263 :       do i = 1, ngrdcol
     149  1583714730 :         rho_ds_dz(i,k) = rho_ds(i,k) * dz(i,k)
     150             :       end do  
     151             :     end do
     152             :     !$acc end parallel loop
     153             :     
     154             : 
     155             :     ! denom_integral does not change throughout the hole filling algorithm
     156             :     ! so we can calculate it before hand. This results in unneccesary computations,
     157             :     ! but is parallelizable and reduces the cost of the serial k loop
     158             :     !$acc parallel loop gang vector collapse(2) default(present)
     159    18631938 :     do i = 1, ngrdcol
     160  1402645370 :       do k = 2+num_draw_pts, upper_hf_level-num_draw_pts
     161  1384013432 :         k_start = k - num_draw_pts
     162  1384013432 :         k_end   = k + num_draw_pts
     163  8321596997 :         invrs_denom_integral(i,k) = one / sum( rho_ds_dz(i,k_start:k_end) )
     164             :       end do  
     165             :     end do
     166             :     !$acc end parallel loop
     167             :     
     168             :     !$acc parallel loop gang vector default(present)
     169    18631938 :     do i = 1, ngrdcol
     170             : 
     171             :       ! Make one pass up the profile, filling holes as much as we can using
     172             :       ! nearby mass, ignoring the first level.
     173             :       !
     174             :       ! This loop can only be done in serial due to the field values required for the next
     175             :       ! iteration potentially changing. We could in theory expose more parallelism in cases 
     176             :       ! where there are large enough gaps between vertical levels which need hole-filling,
     177             :       ! but levels which require hole-filling are often close or consecutive.
     178  1402645370 :       do k = 2+num_draw_pts, upper_hf_level-num_draw_pts
     179             : 
     180  1384013432 :         k_start = k - num_draw_pts
     181  1384013432 :         k_end   = k + num_draw_pts
     182             : 
     183  6759057310 :         if ( any( field(i,k_start:k_end) < threshold ) ) then
     184             : 
     185             :           ! Compute the field's vertical average cenetered at k, which we must conserve,
     186             :           ! see description of the vertical_avg function in advance_helper_module
     187             :           field_avg = sum( rho_ds_dz(i,k_start:k_end) * field(i,k_start:k_end) ) &
     188  2212098792 :                       * invrs_denom_integral(i,k)
     189             : 
     190             :           ! Clip small or negative values from field.
     191   368683132 :           if ( field_avg >= threshold ) then
     192             :             ! We know we can fill in holes completely
     193   293027670 :             field_clipped(i,k_start:k_end) = max( threshold, field(i,k_start:k_end) )
     194             :           else
     195             :             ! We can only fill in holes partly;
     196             :             ! to do so, we remove all mass above threshold.
     197  1919071122 :             field_clipped(i,k_start:k_end) = min( threshold, field(i,k_start:k_end) )
     198             :           endif
     199             : 
     200             :           ! Compute the clipped field's vertical integral.
     201             :           ! clipped_total_mass >= original_total_mass,
     202             :           ! see description of the vertical_avg function in advance_helper_module
     203             :           field_clipped_avg = sum( rho_ds_dz(i,k_start:k_end) * field_clipped(i,k_start:k_end) ) &
     204  2212098792 :                               * invrs_denom_integral(i,k)
     205             : 
     206             :           ! Avoid divide by zero issues by doing nothing if field_clipped_avg ~= threshold
     207   368683132 :           if ( abs(field_clipped_avg-threshold) > abs(field_clipped_avg+threshold)*eps/2) then
     208             :             ! Compute coefficient that makes the clipped field have the same mass as the
     209             :             ! original field.  We should always have mass_fraction > 0.
     210             :             mass_fraction = ( field_avg - threshold ) &
     211   368683132 :                             / ( field_clipped_avg - threshold )
     212             : 
     213             :             ! Calculate normalized, filled field
     214             :             field(i,k_start:k_end) = threshold &
     215  2212098792 :                         + mass_fraction * ( field_clipped(i,k_start:k_end) - threshold )
     216             :           endif
     217             : 
     218             :         endif
     219             : 
     220             :       end do
     221             : 
     222             :     end do
     223             :     !$acc end parallel loop
     224             : 
     225             :     l_field_below_threshold = .false.
     226             : 
     227             :     !$acc parallel loop gang vector collapse(2) default(present) reduction(.or.:l_field_below_threshold)
     228    95935838 :     do k = 1, nz
     229  1584830263 :       do i = 1, ngrdcol
     230  1583714730 :         if ( field(i,k) < threshold ) then
     231   301028512 :           l_field_below_threshold = .true.
     232             :         end if
     233             :       end do
     234             :     end do
     235             :     !$acc end parallel loop
     236             : 
     237             :     ! If all field values are above the threshold, no further hole filling is required
     238     1115533 :     if ( .not. l_field_below_threshold ) then
     239             :       !$acc exit data delete( invrs_denom_integral, field_clipped, denom_integral_global, rho_ds_dz, &
     240             :       !$acc                   numer_integral_global, field_avg_global, mass_fraction_global )
     241             :       return
     242             :     end if
     243             : 
     244             : 
     245             :     ! Now we fill holes globally to maximize the chance that all holes are filled.
     246             :     ! To improve parallelism we assume that global hole filling needs to be done 
     247             :     ! for each grid column, perform all calculations required, and only check
     248             :     ! if any holes need filling before the final step of updating the field. 
     249             : 
     250             :     ! Compute the numerator and denominator integrals
     251             :     !$acc parallel loop gang vector default(present)
     252     7477634 :     do i = 1, ngrdcol
     253     7030014 :       numer_integral_global(i) = 0.0_core_rknd
     254     7477634 :       denom_integral_global(i) = 0.0_core_rknd
     255             :     end do
     256             :     !$acc end parallel loop
     257             : 
     258             :     !$acc parallel loop gang vector default(present)
     259     7477634 :     do i = 1, ngrdcol
     260   590977601 :       do k = 2, upper_hf_level
     261   583499967 :         numer_integral_global(i) = numer_integral_global(i) + rho_ds_dz(i,k) * field(i,k)
     262             : 
     263   590529981 :         denom_integral_global(i) = denom_integral_global(i) + rho_ds_dz(i,k)
     264             :       end do  
     265             :     end do
     266             :     !$acc end parallel loop
     267             : 
     268             :     
     269             :     !$acc parallel loop gang vector default(present)
     270     7477634 :     do i = 1, ngrdcol
     271             : 
     272             :       ! Find the vertical average of field, using the precomputed numerator and denominator,
     273             :       ! see description of the vertical_avg function in advance_helper_module
     274     7030014 :       field_avg_global(i) = numer_integral_global(i) / denom_integral_global(i)
     275             : 
     276             :       ! Clip small or negative values from field.
     277     7477634 :       if ( field_avg_global(i) >= threshold ) then
     278             :         ! We know we can fill in holes completely
     279   586943097 :         field_clipped(i,2:upper_hf_level) = max( threshold, field(i,2:upper_hf_level) )
     280             :       else
     281             :         ! We can only fill in holes partly;
     282             :         ! to do so, we remove all mass above threshold.
     283     3586884 :         field_clipped(i,2:upper_hf_level) = min( threshold, field(i,2:upper_hf_level) )
     284             :       endif
     285             : 
     286             :     end do
     287             :     !$acc end parallel loop
     288             : 
     289             :     ! To compute the clipped field's vertical integral we only need to recompute the numerator
     290             :     !$acc parallel loop gang vector default(present)
     291     7477634 :     do i = 1, ngrdcol
     292     7030014 :       numer_integral_global(i) = 0.0_core_rknd
     293   590977601 :       do k = 2, upper_hf_level
     294   590529981 :         numer_integral_global(i) = numer_integral_global(i) + rho_ds_dz(i,k) * field_clipped(i,k)
     295             :       end do  
     296             :     end do
     297             :     !$acc end parallel loop
     298             : 
     299             :     !$acc parallel loop gang vector default(present)
     300     7477634 :     do i = 1, ngrdcol
     301             : 
     302             :       ! Do not complete calculations or update field values for this 
     303             :       ! column if there are no holes that need filling
     304   232751184 :       if ( any( field(i,2:upper_hf_level) < threshold ) ) then
     305             : 
     306             :         ! Compute the clipped field's vertical integral,
     307             :         ! see description of the vertical_avg function in advance_helper_module
     308     5761933 :         field_clipped_avg = numer_integral_global(i) / denom_integral_global(i)
     309             : 
     310             :         ! Avoid divide by zero issues by doing nothing if field_clipped_avg ~= threshold
     311     5761933 :         if ( abs(field_clipped_avg-threshold) > abs(field_clipped_avg+threshold)*eps/2) then
     312             : 
     313             :           ! Compute coefficient that makes the clipped field have the same mass as the
     314             :           ! original field.  We should always have mass_fraction > 0.
     315             :           mass_fraction_global(i) = ( field_avg_global(i) - threshold ) &
     316     5761933 :                                     / ( field_clipped_avg - threshold )
     317             : 
     318             :           ! Calculate normalized, filled field
     319             :           field(i,2:upper_hf_level) = threshold + mass_fraction_global(i) &
     320             :                                                   * ( field_clipped(i,2:upper_hf_level) &
     321   484002372 :                                                       - threshold )
     322             :         end if
     323             :       end if
     324             : 
     325             :     end do
     326             :     !$acc end parallel loop
     327             :     
     328             :     !$acc exit data delete( invrs_denom_integral, field_clipped, denom_integral_global, rho_ds_dz, &
     329             :     !$acc                   numer_integral_global, field_avg_global, mass_fraction_global )
     330             : 
     331             :     return
     332             : 
     333             :   end subroutine fill_holes_vertical
     334             : 
     335             :   !===============================================================================
     336           0 :   subroutine hole_filling_hm_one_lev( num_hm_fill, hm_one_lev, & ! Intent(in)
     337           0 :                                    hm_one_lev_filled ) ! Intent(out)
     338             : 
     339             :   ! Description:
     340             :   ! Fills holes between same-phase (i.e. either liquid or frozen) hydrometeors for
     341             :   ! one height level.
     342             :   !
     343             :   ! Warning: Do not input hydrometeors of different phases, e.g. liquid and frozen.
     344             :   ! Otherwise heat will not be conserved.
     345             :   !
     346             :   ! References:
     347             :   !
     348             :   ! None
     349             :   !-----------------------------------------------------------------------
     350             : 
     351             :     use constants_clubb, only: &
     352             :         one, & ! Variable(s)
     353             :         zero, &
     354             :         eps
     355             : 
     356             :     use clubb_precision, only: &
     357             :         core_rknd ! Variable(s)
     358             : 
     359             :     use error_code, only: &
     360             :         clubb_at_least_debug_level  ! Procedure
     361             : 
     362             :     implicit none
     363             : 
     364             :     ! Input Variables
     365             :     integer, intent(in) :: num_hm_fill ! number of hydrometeors involved
     366             : 
     367             :     real(kind = core_rknd), dimension(num_hm_fill), intent(in) :: hm_one_lev
     368             : 
     369             :     ! Output Variables
     370             :     real(kind = core_rknd), dimension(num_hm_fill), intent(out) :: hm_one_lev_filled
     371             : 
     372             :     ! Local Variables
     373             :     integer :: num_neg_hm ! number of holes
     374             : 
     375             :     real(kind = core_rknd) :: &
     376             :       total_hole, & ! Size of the hole ( missing mass, less than 0 )
     377             :       total_mass    ! Total mass to fill the hole
     378             :       ! total mass of water substance = total_mass + total_hole
     379             : 
     380             :     integer :: i ! loop iterator
     381             : 
     382             :   !-----------------------------------------------------------------------
     383             : 
     384             :     !----- Begin Code -----
     385             : 
     386             :     ! Initialization
     387           0 :     hm_one_lev_filled = 0._core_rknd
     388             :     total_hole = 0._core_rknd
     389             :     total_mass = 0._core_rknd
     390             :     num_neg_hm = 0
     391             : 
     392             :     ! Determine the total size of the hole and the number of neg. hydrometeors
     393             :     ! and the total mass of hole filling material
     394           0 :     do i=1, num_hm_fill
     395             : !       print *, "hm_one_lev(",i,") = ", hm_one_lev(i)
     396           0 :        if ( hm_one_lev(i) < zero ) then
     397           0 :           total_hole = total_hole + hm_one_lev(i) ! less than zero
     398             :           num_neg_hm = num_neg_hm + 1
     399             :        else
     400           0 :           total_mass = total_mass + hm_one_lev(i)
     401             :        endif
     402             : 
     403             :     enddo
     404             : 
     405             : !    print *, "total_hole = ", total_hole
     406             : !    print *, "total_mass = ", total_mass
     407             : !    print *, "num_neg_hm = ", num_neg_hm
     408             : 
     409             :     ! There is no water substance at all to fill the hole
     410           0 :     if ( abs(total_mass) < eps ) then
     411             : 
     412           0 :        if ( clubb_at_least_debug_level( 2 ) ) then
     413           0 :           print *, "Warning: One-level hole filling was not successful! total_mass ~= 0"
     414             :        endif
     415             : 
     416           0 :        hm_one_lev_filled = hm_one_lev
     417             : 
     418             :        return
     419             :     endif
     420             : 
     421             :     ! Fill the holes and adjust the remaining quantities:
     422             :     ! hm_filled(i) = 0, if hm(i) < 0
     423             :     ! or
     424             :     ! hm_filled(i) = (1 + total_hole/total_mass)*hm(i), if hm(i) > 0
     425           0 :     do i=1, num_hm_fill
     426             : 
     427             :        ! if there is not enough material, fill the holes partially with all the material available
     428           0 :        if ( abs(total_hole) > total_mass ) then
     429             : 
     430           0 :           if ( clubb_at_least_debug_level( 2 ) ) then
     431             :              print *, "Warning: One-level hole filling was not able to fill holes completely!" // &
     432           0 :                       " The holes were filled partially. |total_hole| > total_mass"
     433             :           endif
     434             : 
     435           0 :           hm_one_lev_filled(i) = min(hm_one_lev(i), zero) * ( one + total_mass / total_hole )
     436             : 
     437             :        else ! fill holes completely
     438           0 :           hm_one_lev_filled(i) = max(hm_one_lev(i), zero) * ( one + total_hole / total_mass )
     439             : 
     440             :        endif
     441             : 
     442             :     enddo
     443             : 
     444             :     ! Assertion checks (water substance conservation, non-negativity)
     445           0 :     if ( clubb_at_least_debug_level( 2 ) ) then
     446             : 
     447           0 :        if ( abs(sum( hm_one_lev ) - sum(hm_one_lev_filled)) > &
     448             :             abs(sum( hm_one_lev ) + sum(hm_one_lev_filled)) * eps/2 ) then
     449           0 :           print *, "Warning: Hole filling was not conservative!"
     450             :        endif
     451             : 
     452           0 :        if ( any( hm_one_lev_filled < zero ) ) then
     453           0 :           print *, "Warning: Hole filling failed! A hole could not be filled."
     454             :        endif
     455             : 
     456             :     endif
     457             : 
     458             :     return
     459             : 
     460             :   end subroutine hole_filling_hm_one_lev
     461             :   !-----------------------------------------------------------------------
     462             : 
     463             :   !-----------------------------------------------------------------------
     464           0 :   subroutine fill_holes_hydromet( nz, hydromet_dim, hydromet, & ! Intent(in)
     465           0 :                                   hydromet_filled ) ! Intent(out)
     466             : 
     467             :   ! Description:
     468             :   ! Fills holes between same-phase hydrometeors(i.e. for frozen hydrometeors).
     469             :   ! The hole filling conserves water substance between all same-phase (frozen or liquid)
     470             :   ! hydrometeors at each height level.
     471             :   !
     472             :   ! Attention: The hole filling for the liquid phase hydrometeors is not yet implemented
     473             :   !
     474             :   ! Attention: l_frozen_hm and l_mix_rat_hm need to be set up before this subroutine is called!
     475             :   !
     476             :   ! References:
     477             :   !
     478             :   ! None
     479             :   !-----------------------------------------------------------------------
     480             : 
     481             :     use clubb_precision, only: &
     482             :         core_rknd
     483             : 
     484             :     use array_index, only: &
     485             :         l_frozen_hm, & ! Variable(s)
     486             :         l_mix_rat_hm
     487             : 
     488             :     use constants_clubb, only: &
     489             :         zero
     490             : 
     491             :     implicit none
     492             : 
     493             :     ! Input Variables
     494             :     integer, intent(in) :: hydromet_dim, nz
     495             : 
     496             :     real( kind = core_rknd ), dimension(nz,hydromet_dim), intent(in) :: &
     497             :       hydromet
     498             : 
     499             :     ! Output Variables
     500             :     real( kind = core_rknd ), dimension(nz,hydromet_dim), intent(out) :: &
     501             :       hydromet_filled
     502             : 
     503             :     ! Local Variables
     504             :     integer :: i,j ! Loop iterators
     505             : 
     506             :     integer :: num_frozen_hm ! Number of frozen hydrometeor mixing ratios
     507             : 
     508             :     real( kind = core_rknd ), dimension(:,:), allocatable :: &
     509           0 :       hydromet_frozen,       & ! Frozen hydrometeor mixing ratios
     510           0 :       hydromet_frozen_filled   ! Frozen hydrometeor mixing ratios after hole filling
     511             : 
     512             :   !-----------------------------------------------------------------------
     513             : 
     514             :     !----- Begin Code -----
     515             : 
     516             :     ! Determine the number of frozen hydrometeor mixing ratios
     517           0 :     num_frozen_hm = 0
     518           0 :     do i=1,hydromet_dim
     519           0 :        if ( l_frozen_hm(i) .and. l_mix_rat_hm(i) ) then
     520           0 :           num_frozen_hm = num_frozen_hm + 1
     521             :        endif
     522             :     enddo
     523             : 
     524             :     ! Allocation
     525           0 :     allocate( hydromet_frozen(nz,num_frozen_hm) )
     526           0 :     allocate( hydromet_frozen_filled(nz,num_frozen_hm) )
     527             : 
     528             :     ! Determine frozen hydrometeor mixing ratios
     529           0 :     j = 1
     530           0 :     do i = 1,hydromet_dim
     531           0 :        if ( l_frozen_hm(i) .and. l_mix_rat_hm(i) ) then
     532           0 :           hydromet_frozen(:,j) = hydromet(:,i)
     533           0 :           j = j+1
     534             :        endif
     535             :     enddo
     536             : 
     537             :     ! Fill holes for the frozen hydrometeors
     538           0 :     do i=1,nz
     539           0 :        if ( any( hydromet_frozen(i,:) < zero ) ) then
     540             :           call hole_filling_hm_one_lev( num_frozen_hm, hydromet_frozen(i,:), & ! Intent(in)
     541           0 :                                      hydromet_frozen_filled(i,:) ) ! Intent(out)
     542             :        else
     543           0 :           hydromet_frozen_filled(i,:) = hydromet_frozen(i,:)
     544             :        endif
     545             :     enddo
     546             : 
     547             :     ! Setup the filled hydromet array
     548             :     j = 1
     549           0 :     do i=1, hydromet_dim
     550           0 :        if ( l_frozen_hm(i) .and. l_mix_rat_hm(i) ) then
     551           0 :           hydromet_filled(:,i) = hydromet_frozen_filled(:,j)
     552           0 :           j = j+1
     553             :        else
     554           0 :           hydromet_filled(:,i) = hydromet(:,i)
     555             :        endif
     556             :     enddo
     557             : 
     558             :     !!! Here we could do the same hole filling for all the liquid phase hydrometeors
     559             : 
     560           0 :     return
     561           0 :   end subroutine fill_holes_hydromet
     562             :   !-----------------------------------------------------------------------
     563             : 
     564             :     !-----------------------------------------------------------------------
     565           0 :   subroutine fill_holes_wv( nz, dt, exner, hydromet_name, & ! Intent(in)
     566           0 :                             rvm_mc, thlm_mc, hydromet )! Intent(inout)
     567             : 
     568             :   ! Description:
     569             :   ! Fills holes using the cloud water mixing ratio from the current height level.
     570             :   !
     571             :   ! References:
     572             :   !
     573             :   ! None
     574             :   !-----------------------------------------------------------------------
     575             : 
     576             :     use clubb_precision, only: &
     577             :         core_rknd
     578             : 
     579             :     use constants_clubb, only: &
     580             :         zero_threshold, &
     581             :         Lv, &
     582             :         Ls, &
     583             :         Cp
     584             : 
     585             :     implicit none
     586             : 
     587             :     ! Input Variables
     588             :     integer, intent(in) :: nz
     589             : 
     590             :     real( kind = core_rknd ), intent(in) ::  &
     591             :       dt           ! Timestep         [s]
     592             : 
     593             :     character(len=10), intent(in) :: hydromet_name
     594             : 
     595             :     real( kind = core_rknd ), dimension(nz), intent(in) :: &
     596             :       exner        ! Exner function                            [-]
     597             : 
     598             :     ! Input/Output Variables
     599             :     real( kind = core_rknd ), dimension(nz), intent(inout) :: &
     600             :       hydromet, &  ! Hydrometeor array                         [units vary]
     601             :       rvm_mc, &
     602             :       thlm_mc
     603             : 
     604             :     ! Local Variables
     605             :     integer :: k ! Loop iterator
     606             : 
     607             :     real( kind = core_rknd ) :: rvm_clip_tndcy
     608             :   !-----------------------------------------------------------------------
     609             : 
     610             :     !----- Begin Code -----
     611             : 
     612           0 :     do k = 2, nz, 1
     613             : 
     614           0 :        if ( hydromet(k) < zero_threshold ) then
     615             : 
     616             :           ! Set rvm_clip_tndcy to the time tendency applied to vapor and removed
     617             :           ! from the hydrometeor.
     618           0 :           rvm_clip_tndcy = hydromet(k) / dt
     619             : 
     620             :           ! Adjust the tendency rvm_mc accordingly
     621           0 :           rvm_mc(k) = rvm_mc(k) + rvm_clip_tndcy
     622             : 
     623             :           ! Adjust the tendency of thlm_mc according to whether the
     624             :           ! effect is an evaporation or sublimation tendency.
     625           0 :           select case ( trim( hydromet_name ) )
     626             :           case( "rrm" )
     627           0 :              thlm_mc(k) = thlm_mc(k) - rvm_clip_tndcy * ( Lv / ( Cp*exner(k) ) )
     628             :           case( "rim", "rsm", "rgm" )
     629           0 :              thlm_mc(k) = thlm_mc(k) - rvm_clip_tndcy * ( Ls / ( Cp*exner(k) ) )
     630             :           case default
     631           0 :              error stop "Fatal error in microphys_driver"
     632             :           end select
     633             : 
     634             :           ! Set the mixing ratio to 0
     635           0 :           hydromet(k) = zero_threshold
     636             : 
     637             :        endif ! hydromet(k,i) < 0
     638             : 
     639             :     enddo ! k = 2..gr%nz
     640             : 
     641           0 :     return
     642             :   end subroutine fill_holes_wv
     643             :   !-----------------------------------------------------------------------
     644             : 
     645             :   !-----------------------------------------------------------------------
     646           0 :   subroutine fill_holes_driver( gr, nz, dt, hydromet_dim,        & ! Intent(in)
     647             :                                 l_fill_holes_hm,             & ! Intent(in)
     648           0 :                                 rho_ds_zm, rho_ds_zt, exner, & ! Intent(in)
     649             :                                 stats_metadata,              & ! Intent(in)
     650             :                                 stats_zt,                    & ! intent(inout)
     651           0 :                                 thlm_mc, rvm_mc, hydromet )    ! Intent(inout)
     652             : 
     653             :   ! Description:
     654             :   ! Fills holes between same-phase hydrometeors(i.e. for frozen hydrometeors).
     655             :   ! The hole filling conserves water substance between all same-phase (frozen or liquid)
     656             :   ! hydrometeors at each height level.
     657             :   !
     658             :   ! Attention: The hole filling for the liquid phase hydrometeors is not yet implemented
     659             :   !
     660             :   ! Attention: l_frozen_hm and l_mix_rat_hm need to be set up before this subroutine is called!
     661             :   !
     662             :   ! References:
     663             :   !
     664             :   ! None
     665             :   !-----------------------------------------------------------------------
     666             : 
     667             :     use grid_class, only: &
     668             :         grid ! Type
     669             : 
     670             :     use clubb_precision, only: &
     671             :         core_rknd   ! Variable(s)
     672             : 
     673             :     use constants_clubb, only: &
     674             :         zero,            &
     675             :         zero_threshold,  &
     676             :         Lv,              &
     677             :         Ls,              &
     678             :         Cp,              &
     679             :         fstderr,         &
     680             :         num_hf_draw_points
     681             : 
     682             :     use array_index, only: &
     683             :         hydromet_list, & ! Names of the hydrometeor species
     684             :         hydromet_tol
     685             : 
     686             :     use array_index, only: &
     687             :         l_mix_rat_hm, & ! Variable(s)
     688             :         l_frozen_hm
     689             : 
     690             :     use index_mapping, only: &
     691             :         Nx2rx_hm_idx, & ! Procedure(s)
     692             :         mvr_hm_max
     693             : 
     694             :     use stats_type_utilities, only: &
     695             :         stat_begin_update, & ! Subroutines
     696             :         stat_end_update
     697             : 
     698             :     use stats_variables, only: &
     699             :         stats_metadata_type
     700             : 
     701             :     use error_code, only: &
     702             :         clubb_at_least_debug_level  ! Procedure
     703             : 
     704             :     use stats_type, only: stats ! Type
     705             : 
     706             :     implicit none
     707             : 
     708             :     !----------------------- Input Variables -----------------------
     709             :     type (grid), target, intent(in) :: gr
     710             : 
     711             :     integer, intent(in) :: hydromet_dim, nz
     712             : 
     713             :     logical, intent(in) :: l_fill_holes_hm
     714             : 
     715             :     real( kind = core_rknd ), intent(in) ::  &
     716             :       dt           ! Timestep         [s]
     717             : 
     718             :     real( kind = core_rknd ), dimension(nz), intent(in) :: &
     719             :       rho_ds_zm, & ! Dry, static density on momentum levels   [kg/m^3]
     720             :       rho_ds_zt    ! Dry, static density on thermo. levels    [kg/m^3]
     721             : 
     722             :     real( kind = core_rknd ), dimension(nz), intent(in) :: &
     723             :       exner  ! Exner function                                       [-]
     724             : 
     725             :     type (stats_metadata_type), intent(in) :: &
     726             :       stats_metadata
     727             : 
     728             :     !----------------------- Input/Output Variables -----------------------
     729             :     type (stats), target, intent(inout) :: &
     730             :       stats_zt
     731             : 
     732             :     real( kind = core_rknd ), dimension(nz, hydromet_dim), intent(inout) :: &
     733             :       hydromet    ! Mean of hydrometeor fields    [units vary]
     734             : 
     735             :     real( kind = core_rknd ), dimension(nz), intent(inout) :: &
     736             :       rvm_mc,  & ! Microphysics contributions to vapor water           [kg/kg/s]
     737             :       thlm_mc    ! Microphysics contributions to liquid potential temp [K/s]
     738             : 
     739             :     !----------------------- Local Variables -----------------------
     740             :     integer :: i, k ! Loop iterators
     741             : 
     742             :     real( kind = core_rknd ), dimension(nz, hydromet_dim) :: &
     743           0 :       hydromet_filled,  & ! Frozen hydrometeor mixing ratios after hole filling
     744           0 :       hydromet_clipped    ! Clipped mean of hydrometeor fields    [units vary]
     745             : 
     746             :     character( len = 10 ) :: hydromet_name
     747             : 
     748             :     real( kind = core_rknd ) :: &
     749             :       max_velocity    ! Maximum sedimentation velocity                     [m/s]
     750             : 
     751             :     integer :: ixrm_hf, ixrm_wvhf, ixrm_cl, &
     752             :                ixrm_bt, ixrm_mc
     753             : 
     754             :     logical :: l_hole_fill = .true.
     755             : 
     756             :     !----------------------- Begin Code -----------------------
     757             : 
     758             :     ! Start stats output for the _hf variables (changes in the hydromet array
     759             :     ! due to fill_holes_hydromet and fill_holes_vertical)
     760           0 :     if ( stats_metadata%l_stats_samp ) then
     761             : 
     762           0 :        do i = 1, hydromet_dim
     763             : 
     764             :           ! Set up the stats indices for hydrometeor at index i
     765             :           call setup_stats_indices( i, stats_metadata,           & ! Intent(in)
     766             :                                     ixrm_bt, ixrm_hf, ixrm_wvhf, & ! Intent(inout)
     767             :                                     ixrm_cl, ixrm_mc,            & ! Intent(inout)
     768           0 :                                     max_velocity )                 ! Intent(inout)
     769             : 
     770           0 :           call stat_begin_update( gr%nz, ixrm_hf, hydromet(:,i) / dt, & ! intent(in)
     771           0 :                                   stats_zt ) ! intent(inout)
     772             : 
     773             :        enddo ! i = 1, hydromet_dim
     774             : 
     775             :     endif ! stats_metadata%l_stats_samp
     776             : 
     777             :     ! If we're dealing with negative hydrometeors, we first try to fill the
     778             :     ! holes proportionally from other same-phase hydrometeors at each height
     779             :     ! level.
     780           0 :     if ( any( hydromet < zero_threshold ) .and. l_fill_holes_hm ) then
     781             : 
     782             :        call fill_holes_hydromet( nz, hydromet_dim, hydromet, & ! Intent(in)
     783           0 :                                  hydromet_filled ) ! Intent(out)
     784             : 
     785           0 :        hydromet = hydromet_filled
     786             : 
     787             :     endif ! any( hydromet < zero ) .and. l_fill_holes_hm
     788             : 
     789           0 :     hydromet_filled = zero
     790             : 
     791           0 :     do i = 1, hydromet_dim
     792             : 
     793             :       ! Set up the stats indices for hydrometeor at index i
     794             :       call setup_stats_indices( i, stats_metadata,           & ! Intent(in)
     795             :                                 ixrm_bt, ixrm_hf, ixrm_wvhf, & ! Intent(inout)
     796             :                                 ixrm_cl, ixrm_mc,            & ! Intent(inout)
     797           0 :                                 max_velocity )                 ! Intent(inout)
     798             : 
     799           0 :       hydromet_name = hydromet_list(i)
     800             : 
     801             :       ! Print warning message if any hydrometeor species has a value < 0.
     802           0 :       if ( clubb_at_least_debug_level( 1 ) ) then
     803           0 :          if ( any( hydromet(:,i) < zero_threshold ) ) then
     804             : 
     805           0 :             do k = 1, nz
     806           0 :                if ( hydromet(k,i) < zero_threshold ) then
     807           0 :                   write(fstderr,*) trim( hydromet_name ) //" < ", &
     808           0 :                                    zero_threshold, &
     809           0 :                                    " in fill_holes_driver at k= ", k
     810             :                endif ! hydromet(k,i) < 0
     811             :             enddo ! k = 1, nz
     812             : 
     813             :          endif ! hydromet(:,i) < 0       
     814             :       endif ! clubb_at_least_debug_level( 1 )
     815             : 
     816             : 
     817             :       ! Store the previous value of the hydrometeor for the effect of the
     818             :       ! hole-filling scheme.
     819             : !      if ( stats_metadata%l_stats_samp ) then
     820             : !         call stat_begin_update( ixrm_hf, hydromet(:,i) &
     821             : !                                          / dt, stats_zt )
     822             : !      endif
     823             : 
     824             :       ! If we're dealing with a mixing ratio and hole filling is enabled,
     825             :       ! then we apply the hole filling algorithm
     826           0 :       if ( any( hydromet(:,i) < zero_threshold ) ) then
     827             : 
     828           0 :          if ( hydromet_name(1:1) == "r" .and. l_hole_fill ) then
     829             : 
     830             :             !$acc data copyin( gr, gr%dzt, rho_ds_zt ) &
     831             :             !$acc        copy( hydromet(:,i) )
     832             : 
     833             :             ! Apply the hole filling algorithm
     834             :             ! upper_hf_level = nz since we are filling the zt levels
     835             :             call fill_holes_vertical( gr%nz, 1, num_hf_draw_points, zero_threshold, gr%nz, & ! In
     836             :                                       gr%dzt, rho_ds_zt,                                   & ! In
     837           0 :                                       hydromet(:,i) )                                      ! InOut
     838             : 
     839             :             !$acc end data
     840             : 
     841             :          endif ! Variable is a mixing ratio and l_hole_fill is true
     842             : 
     843             :       endif ! hydromet(:,i) < 0
     844             : 
     845             :       ! Enter the new value of the hydrometeor for the effect of the
     846             :       ! hole-filling scheme.
     847           0 :       if ( stats_metadata%l_stats_samp ) then
     848           0 :          call stat_end_update( gr%nz, ixrm_hf, hydromet(:,i) / dt, & ! intent(in)
     849           0 :                                stats_zt ) ! intent(inout)
     850             :       endif
     851             : 
     852             :       ! Store the previous value of the hydrometeor for the effect of the water
     853             :       ! vapor hole-filling scheme.
     854           0 :       if ( stats_metadata%l_stats_samp ) then
     855           0 :          call stat_begin_update( gr%nz, ixrm_wvhf, hydromet(:,i) / dt, & ! intent(in)
     856           0 :                                  stats_zt ) ! intent(inout)
     857             :       endif
     858             : 
     859           0 :       if ( any( hydromet(:,i) < zero_threshold ) ) then
     860             : 
     861           0 :          if ( hydromet_name(1:1) == "r" .and. l_hole_fill ) then
     862             : 
     863             :             ! If the hole filling algorithm failed, then we attempt to fill
     864             :             ! the missing mass with water vapor mixing ratio.
     865             :             ! We noticed this is needed for ASEX A209, particularly if Latin
     866             :             ! hypercube sampling is enabled.  -dschanen 11 Nov 2010
     867             :             call fill_holes_wv( nz, dt, exner, hydromet_name, & ! Intent(in)
     868           0 :                                 rvm_mc, thlm_mc, hydromet(:,i) )   ! Intent(out)
     869             : 
     870             :          endif ! Variable is a mixing ratio and l_hole_fill is true
     871             : 
     872             :       endif ! hydromet(:,i) < 0
     873             : 
     874             :       ! Enter the new value of the hydrometeor for the effect of the water vapor
     875             :       ! hole-filling scheme.
     876           0 :       if ( stats_metadata%l_stats_samp ) then
     877           0 :          call stat_end_update( gr%nz, ixrm_wvhf, hydromet(:,i) / dt, & ! intent(in)
     878           0 :                                stats_zt ) ! intent(inout)
     879             :       endif
     880             : 
     881             :       ! Clipping for hydrometeor mixing ratios.
     882           0 :       if ( l_mix_rat_hm(i) ) then
     883             : 
     884             :          ! Store the previous value of the hydrometeor for the effect of
     885             :          ! clipping.
     886           0 :          if ( stats_metadata%l_stats_samp ) then
     887           0 :             call stat_begin_update( gr%nz, ixrm_cl, hydromet(:,i) / dt, & ! intent(in)
     888           0 :                                     stats_zt ) ! intent(inout)
     889             :          endif
     890             : 
     891           0 :          if ( any( hydromet(:,i) < zero_threshold ) ) then
     892             : 
     893             :             ! Clip any remaining negative values of precipitating hydrometeor
     894             :             ! mixing ratios to 0.
     895           0 :             where ( hydromet(:,i) < zero_threshold )
     896             :                hydromet(:,i) = zero_threshold
     897             :             end where
     898             : 
     899             :          endif ! hydromet(:,i) < 0
     900             : 
     901             :          ! Eliminate very small values of mean precipitating hydrometeor mixing
     902             :          ! ratios by setting them to 0.
     903           0 :          do k = 2, gr%nz, 1
     904             : 
     905           0 :             if ( hydromet(k,i) <= hydromet_tol(i) ) then
     906             : 
     907             :                rvm_mc(k) &
     908             :                = rvm_mc(k) &
     909           0 :                  + ( hydromet(k,i) / dt )
     910             : 
     911           0 :                if ( .not. l_frozen_hm(i) ) then
     912             : 
     913             :                   ! Rain water mixing ratio
     914             :    
     915             :                   thlm_mc(k) &
     916             :                   = thlm_mc(k) &
     917             :                     - ( Lv / ( Cp * exner(k) ) ) &
     918           0 :                       * ( hydromet(k,i) / dt )
     919             : 
     920             :                else ! Frozen hydrometeor mixing ratio
     921             : 
     922             :                   thlm_mc(k) &
     923             :                   = thlm_mc(k) &
     924             :                     - ( Ls / ( Cp * exner(k) ) ) &
     925           0 :                       * ( hydromet(k,i) / dt )
     926             : 
     927             :                endif ! l_frozen_hm(i)
     928             : 
     929           0 :                hydromet(k,i) = zero
     930             : 
     931             :             endif ! hydromet(k,i) <= hydromet_tol(i)
     932             : 
     933             :          enddo ! k = 2, gr%nz, 1
     934             : 
     935             : 
     936             :          ! Enter the new value of the hydrometeor for the effect of clipping.
     937           0 :          if ( stats_metadata%l_stats_samp ) then
     938             :             call stat_end_update( gr%nz, ixrm_cl, hydromet(:,i) / dt, & ! intent(in)
     939           0 :                                   stats_zt ) ! intent(inout)
     940             :          endif
     941             : 
     942             :       endif ! l_mix_rat_hm(i)
     943             : 
     944             :     enddo ! i = 1, hydromet_dim, 1
     945             : 
     946             :     ! Calculate clipping for hydrometeor concentrations.
     947             :     call clip_hydromet_conc_mvr( gr%nz, hydromet_dim, hydromet, & ! Intent(in)
     948           0 :                                  hydromet_clipped )        ! Intent(out)
     949             : 
     950             :     ! Clip hydrometeor concentrations and output stats.
     951           0 :     do i = 1, hydromet_dim
     952             : 
     953           0 :        if ( .not. l_mix_rat_hm(i) ) then
     954             : 
     955           0 :           if ( stats_metadata%l_stats_samp ) then
     956             : 
     957             :              ! Set up the stats indices for hydrometeor at index i
     958             :              call setup_stats_indices( i, stats_metadata,           & ! In
     959             :                                        ixrm_bt, ixrm_hf, ixrm_wvhf, & ! In/Out
     960             :                                        ixrm_cl, ixrm_mc,            & ! In/Out
     961           0 :                                        max_velocity )                 ! In/Out
     962             : 
     963             :              ! Store the previous value of the hydrometeor for the effect of
     964             :              ! clipping.
     965           0 :              call stat_begin_update( gr%nz, ixrm_cl, hydromet(:,i) / dt, & ! intent(in)
     966           0 :                                      stats_zt ) ! intent(inout)
     967             : 
     968             :           endif ! stats_metadata%l_stats_samp
     969             : 
     970             :           ! Apply clipping of hydrometeor concentrations.
     971           0 :           hydromet(:,i) = hydromet_clipped(:,i)
     972             : 
     973             :           ! Enter the new value of the hydrometeor for the effect of clipping.
     974           0 :           if ( stats_metadata%l_stats_samp ) then
     975             :              call stat_end_update( gr%nz, ixrm_cl, hydromet(:,i) / dt, & ! intent(in)
     976           0 :                                    stats_zt ) ! intent(inout)
     977             :           endif
     978             : 
     979             :        endif ! .not. l_mix_rat_hm(i)
     980             : 
     981             :     enddo ! i = 1, hydromet_dim, 1
     982             : 
     983             : 
     984           0 :     return
     985             : 
     986             :   end subroutine fill_holes_driver
     987             : 
     988             :   !=============================================================================
     989           0 :   subroutine clip_hydromet_conc_mvr( nz, hydromet_dim, hydromet, & ! Intent(in)
     990           0 :                                      hydromet_clipped )        ! Intent(out)
     991             : 
     992             :     ! Description:
     993             :     ! Increases the value of a hydrometeor concentration when it is too small,
     994             :     ! according to the hydrometeor mixing ratio (which remains unchanged) and
     995             :     ! the maximum acceptable drop or particle mean volume radius for that
     996             :     ! hydrometeor species.
     997             : 
     998             :     ! References:
     999             :     !-----------------------------------------------------------------------
    1000             : 
    1001             :     use grid_class, only: &
    1002             :         grid ! Type
    1003             : 
    1004             :     use constants_clubb, only: &
    1005             :         pi,          & ! Variable(s)
    1006             :         four_thirds, &
    1007             :         one,         &
    1008             :         zero,        &
    1009             :         rho_lw,      &
    1010             :         rho_ice
    1011             : 
    1012             :     use array_index, only: &
    1013             :         l_mix_rat_hm, & ! Variable(s)
    1014             :         l_frozen_hm
    1015             : 
    1016             :     use index_mapping, only: &
    1017             :         Nx2rx_hm_idx, & ! Procedure(s)
    1018             :         mvr_hm_max
    1019             : 
    1020             :     use clubb_precision, only: &
    1021             :         core_rknd    ! Variable(s)
    1022             : 
    1023             :     implicit none
    1024             : 
    1025             :     ! Input Variables
    1026             :     integer, intent(in) :: &
    1027             :       nz
    1028             :     
    1029             :     integer, intent(in) :: &
    1030             :       hydromet_dim    ! Number of hydrometeor fields
    1031             : 
    1032             :     real( kind = core_rknd ), dimension(nz,hydromet_dim), intent(in) :: &
    1033             :       hydromet    ! Mean of hydrometeor fields    [units vary]
    1034             : 
    1035             :     ! Output Variable
    1036             :     real( kind = core_rknd ), dimension(nz,hydromet_dim), intent(out) :: &
    1037             :       hydromet_clipped    ! Clipped mean of hydrometeor fields    [units vary]
    1038             : 
    1039             :     ! Local Variables
    1040             :     real( kind = core_rknd ) :: &
    1041             :       Nxm_min_coef    ! Coefficient for min mean value of a concentration [1/kg]
    1042             : 
    1043             :     integer :: &
    1044             :       k,   & ! Vertical grid index
    1045             :       idx    ! Hydrometeor species index
    1046             : 
    1047             : 
    1048             :     ! Clipping for hydrometeor concentrations.
    1049           0 :     do idx = 1, hydromet_dim
    1050             : 
    1051           0 :        if ( .not. l_mix_rat_hm(idx) ) then
    1052             : 
    1053             :           ! The field is a hydrometeor concentration.
    1054             : 
    1055           0 :           if ( .not. l_frozen_hm(idx) ) then
    1056             : 
    1057             :              ! Clipping for mean rain drop concentration, <Nr>.
    1058             :              !
    1059             :              ! When mean rain water mixing ratio, <rr>, is found at a grid
    1060             :              ! level, mean rain drop concentration must be at least a minimum
    1061             :              ! value so that average rain drop mean volume radius stays within
    1062             :              ! an upper bound.  Otherwise, mean rain drop concentration is 0.
    1063             :              ! This can also be applied to values at a sample point, rather than
    1064             :              ! a grid-box mean.
    1065             : 
    1066             :              ! The minimum mean rain drop concentration is given by:
    1067             :              !
    1068             :              ! <Nr> = <rr> / ( (4/3) * pi * rho_lw * mvr_rain_max^3 );
    1069             :              !
    1070             :              ! where mvr_rain_max is the maximum acceptable average rain drop
    1071             :              ! mean volume radius.
    1072             : 
    1073             :              Nxm_min_coef &
    1074           0 :              = one / ( four_thirds * pi * rho_lw * mvr_hm_max(idx)**3 )
    1075             : 
    1076             :           else ! l_frozen_hm(idx)
    1077             : 
    1078             :              ! Clipping for mean frozen hydrometeor concentration, <Nx>, where x
    1079             :              ! stands for any frozen hydrometeor species.
    1080             :              !
    1081             :              ! When mean frozen hydrometeor mixing ratio, <rx>, is found at a
    1082             :              ! grid level, mean frozen hydrometeor concentration must be at
    1083             :              ! least a minimum value so that average frozen hydrometeor mean
    1084             :              ! volume radius stays within an upper bound.  Otherwise, mean
    1085             :              ! frozen hydrometeor concentration is 0.  This can also be applied
    1086             :              ! to values at a sample point, rather than a grid-box mean.
    1087             : 
    1088             :              ! The minimum mean frozen hydrometeor concentration is given by:
    1089             :              !
    1090             :              ! <Nx> = <rx> / ( (4/3) * pi * rho_ice * mvr_x_max^3 );
    1091             :              !
    1092             :              ! where mvr_x_max is the maximum acceptable average frozen
    1093             :              ! hydrometeor mean volume radius for frozen hydrometeor species, x.
    1094             : 
    1095             :              Nxm_min_coef &
    1096           0 :              = one / ( four_thirds * pi * rho_ice * mvr_hm_max(idx)**3 )
    1097             : 
    1098             :           endif ! .not. l_frozen_hm(idx)
    1099             : 
    1100             :           ! Loop over vertical levels and increase hydrometeor concentrations
    1101             :           ! when necessary.
    1102           0 :           do k = 2, nz, 1
    1103             : 
    1104           0 :              if ( hydromet(k,Nx2rx_hm_idx(idx)) > zero ) then
    1105             : 
    1106             :                 ! Hydrometeor mixing ratio, <rx>, is found at the grid level.
    1107           0 :                 hydromet_clipped(k,idx) &
    1108             :                 = max( hydromet(k,idx), &
    1109           0 :                        Nxm_min_coef * hydromet(k,Nx2rx_hm_idx(idx)) )
    1110             : 
    1111             :              else ! <rx> = 0
    1112             : 
    1113           0 :                 hydromet_clipped(k,idx) = zero
    1114             : 
    1115             :              endif ! hydromet(k,Nx2rx_hm_idx(idx)) > 0
    1116             : 
    1117             :           enddo ! k = 2, nz, 1
    1118             : 
    1119             :           ! The lowest thermodynamic level is below the model's lower boundary.
    1120           0 :           hydromet_clipped(1,idx) = hydromet(1,idx)
    1121             : 
    1122             :        else ! l_mix_rat_hm(idx)
    1123             : 
    1124             :           ! The field is a hydrometeor mixing ratio.
    1125           0 :           hydromet_clipped(:,idx) = hydromet(:,idx)
    1126             : 
    1127             :        endif ! .not. l_mix_rat_hm(idx)
    1128             : 
    1129             :     enddo ! idx = 1, hydromet_dim, 1
    1130             : 
    1131             : 
    1132           0 :     return
    1133             : 
    1134             :   end subroutine clip_hydromet_conc_mvr
    1135             : 
    1136             :   !-----------------------------------------------------------------------
    1137           0 :   subroutine setup_stats_indices( ihm, stats_metadata,         & ! Intent(in)
    1138             :                                   ixrm_bt, ixrm_hf, ixrm_wvhf, & ! Intent(inout)
    1139             :                                   ixrm_cl, ixrm_mc,            & ! Intent(inout)
    1140             :                                   max_velocity )                 ! Intent(inout)
    1141             : 
    1142             :   ! Description:
    1143             :   !
    1144             :   ! Determines the stats output indices depending on the hydrometeor.
    1145             : 
    1146             :   ! Attention: hydromet_list needs to be set up before this routine is called.
    1147             :   !
    1148             :   ! Bogus example
    1149             :   ! References:
    1150             :   !
    1151             :   ! None
    1152             :   !-----------------------------------------------------------------------
    1153             : 
    1154             : 
    1155             :     use array_index, only: &
    1156             :         hydromet_list  ! Names of the hydrometeor species
    1157             : 
    1158             :     use clubb_precision, only: &
    1159             :         core_rknd
    1160             : 
    1161             :     use constants_clubb, only: &
    1162             :         zero
    1163             : 
    1164             :     use stats_variables, only: & 
    1165             :         stats_metadata_type
    1166             : 
    1167             :     implicit none
    1168             : 
    1169             :     ! Input Variables
    1170             :     integer, intent(in) :: ihm
    1171             : 
    1172             :     type (stats_metadata_type), intent(in) :: &
    1173             :       stats_metadata
    1174             : 
    1175             :     ! Input/Output Variables
    1176             :     real( kind = core_rknd ), intent(inout) :: &
    1177             :       max_velocity ! Maximum sedimentation velocity [m/s]
    1178             : 
    1179             :     integer, intent(inout) :: ixrm_hf, ixrm_wvhf, ixrm_cl, &
    1180             :                               ixrm_bt, ixrm_mc
    1181             : 
    1182             :   !-----------------------------------------------------------------------
    1183             : 
    1184             :     !----- Begin Code -----
    1185             : 
    1186             :     ! Initializing max_velocity in order to avoid a compiler warning.
    1187             :     ! Regardless of the case, it will be reset in the 'select case'
    1188             :     ! statement immediately below.
    1189           0 :     max_velocity = zero
    1190             : 
    1191           0 :     select case ( trim( hydromet_list(ihm) ) )
    1192             :       case ( "rrm" )
    1193           0 :         ixrm_bt   = stats_metadata%irrm_bt
    1194           0 :         ixrm_hf   = stats_metadata%irrm_hf
    1195           0 :         ixrm_wvhf = stats_metadata%irrm_wvhf
    1196           0 :         ixrm_cl   = stats_metadata%irrm_cl
    1197           0 :         ixrm_mc   = stats_metadata%irrm_mc
    1198             : 
    1199           0 :         max_velocity = -9.1_core_rknd ! m/s
    1200             : 
    1201             :       case ( "rim" )
    1202           0 :         ixrm_bt   = stats_metadata%irim_bt
    1203           0 :         ixrm_hf   = stats_metadata%irim_hf
    1204           0 :         ixrm_wvhf = stats_metadata%irim_wvhf
    1205           0 :         ixrm_cl   = stats_metadata%irim_cl
    1206           0 :         ixrm_mc   = stats_metadata%irim_mc
    1207             : 
    1208           0 :         max_velocity = -1.2_core_rknd ! m/s
    1209             : 
    1210             :       case ( "rsm" )
    1211           0 :         ixrm_bt   = stats_metadata%irsm_bt
    1212           0 :         ixrm_hf   = stats_metadata%irsm_hf
    1213           0 :         ixrm_wvhf = stats_metadata%irsm_wvhf
    1214           0 :         ixrm_cl   = stats_metadata%irsm_cl
    1215           0 :         ixrm_mc   = stats_metadata%irsm_mc
    1216             : 
    1217             :         ! Morrison limit
    1218             : !         max_velocity = -1.2_core_rknd ! m/s
    1219             :         ! Made up limit.  The literature suggests that it is quite possible
    1220             :         ! that snow flake might achieve a terminal velocity of 2 m/s, and this
    1221             :         ! happens in the COAMPS microphysics -dschanen 29 Sept 2009
    1222           0 :         max_velocity = -2.0_core_rknd ! m/s
    1223             : 
    1224             :       case ( "rgm" )
    1225           0 :         ixrm_bt   = stats_metadata%irgm_bt
    1226           0 :         ixrm_hf   = stats_metadata%irgm_hf
    1227           0 :         ixrm_wvhf = stats_metadata%irgm_wvhf
    1228           0 :         ixrm_cl   = stats_metadata%irgm_cl
    1229           0 :         ixrm_mc   = stats_metadata%irgm_mc
    1230             : 
    1231           0 :         max_velocity = -20._core_rknd ! m/s
    1232             : 
    1233             :       case ( "Nrm" )
    1234           0 :         ixrm_bt   = stats_metadata%iNrm_bt
    1235           0 :         ixrm_hf   = 0
    1236           0 :         ixrm_wvhf = 0
    1237           0 :         ixrm_cl   = stats_metadata%iNrm_cl
    1238           0 :         ixrm_mc   = stats_metadata%iNrm_mc
    1239             : 
    1240           0 :         max_velocity = -9.1_core_rknd ! m/s
    1241             : 
    1242             :       case ( "Nim" )
    1243           0 :         ixrm_bt   = stats_metadata%iNim_bt
    1244           0 :         ixrm_hf   = 0
    1245           0 :         ixrm_wvhf = 0
    1246           0 :         ixrm_cl   = stats_metadata%iNim_cl
    1247           0 :         ixrm_mc   = stats_metadata%iNim_mc
    1248             : 
    1249           0 :         max_velocity = -1.2_core_rknd ! m/s
    1250             : 
    1251             :       case ( "Nsm" )
    1252           0 :         ixrm_bt   = stats_metadata%iNsm_bt
    1253           0 :         ixrm_hf   = 0
    1254           0 :         ixrm_wvhf = 0
    1255           0 :         ixrm_cl   = stats_metadata%iNsm_cl
    1256           0 :         ixrm_mc   = stats_metadata%iNsm_mc
    1257             : 
    1258             :         ! Morrison limit
    1259             : !         max_velocity = -1.2_core_rknd ! m/s
    1260             :         ! Made up limit.  The literature suggests that it is quite possible
    1261             :         ! that snow flake might achieve a terminal velocity of 2 m/s, and this
    1262             :         ! happens in the COAMPS microphysics -dschanen 29 Sept 2009
    1263           0 :         max_velocity = -2.0_core_rknd ! m/s
    1264             : 
    1265             :       case ( "Ngm" )
    1266           0 :         ixrm_bt   = stats_metadata%iNgm_bt
    1267           0 :         ixrm_hf   = 0
    1268           0 :         ixrm_wvhf = 0
    1269           0 :         ixrm_cl   = stats_metadata%iNgm_cl
    1270           0 :         ixrm_mc   = stats_metadata%iNgm_mc
    1271             : 
    1272           0 :         max_velocity = -20._core_rknd ! m/s
    1273             : 
    1274             :       case ( "Ncm" )
    1275           0 :         ixrm_bt   = stats_metadata%iNcm_bt
    1276           0 :         ixrm_hf   = 0
    1277           0 :         ixrm_wvhf = 0
    1278           0 :         ixrm_cl   = stats_metadata%iNcm_cl
    1279           0 :         ixrm_mc   = stats_metadata%iNcm_mc
    1280             : 
    1281             :         ! Use the rain water limit, since Morrison has no explicit limit on
    1282             :         ! cloud water.  Presumably these numbers are never large.
    1283             :         ! -dschanen 28 Sept 2009
    1284           0 :         max_velocity = -9.1_core_rknd ! m/s
    1285             : 
    1286             :       case default
    1287           0 :         ixrm_bt   = 0
    1288           0 :         ixrm_hf   = 0
    1289           0 :         ixrm_wvhf = 0
    1290           0 :         ixrm_cl   = 0
    1291           0 :         ixrm_mc   = 0
    1292             : 
    1293           0 :         max_velocity = -9.1_core_rknd ! m/s
    1294             : 
    1295             :     end select
    1296             : 
    1297             : 
    1298           0 :     return
    1299             : 
    1300             :   end subroutine setup_stats_indices
    1301             :   !-----------------------------------------------------------------------
    1302             : 
    1303             : end module fill_holes

Generated by: LCOV version 1.14