<p><b>laura@ucar.edu</b> 2012-03-01 16:20:34 -0700 (Thu, 01 Mar 2012)</p><p>added calculation of accumulated longwave and shortwave radiation diagnostics for regional climate simulations<br>
</p><hr noshade><pre><font color="gray">Modified: branches/atmos_physics/src/core_atmos_physics/Makefile
===================================================================
--- branches/atmos_physics/src/core_atmos_physics/Makefile        2012-03-01 23:17:44 UTC (rev 1568)
+++ branches/atmos_physics/src/core_atmos_physics/Makefile        2012-03-01 23:20:34 UTC (rev 1569)
@@ -194,23 +194,28 @@
         mpas_atmphys_vars.o
 
 mpas_atmphys_driver_microphysics.o: \
-        ./physics_wrf/module_mp_kessler.o     \
-        ./physics_wrf/module_mp_thompson.o    \
-        ./physics_wrf/module_mp_wsm6.o        \
+        ./physics_wrf/module_mp_kessler.o   \
+        ./physics_wrf/module_mp_thompson.o  \
+        ./physics_wrf/module_mp_wsm6.o      \
         mpas_atmphys_constants.o            \
         mpas_atmphys_interface_nhyd.o       \
         mpas_atmphys_vars.o
 
 mpas_atmphys_driver.o: \
-        mpas_atmphys_driver_convection_deep.o       \
-        mpas_atmphys_driver_pbl.o                   \
-        mpas_atmphys_driver_radiation_lw.o          \
-        mpas_atmphys_driver_radiation_sw.o          \
-        mpas_atmphys_driver_sfclayer.o              \
-        mpas_atmphys_constants.o            \
-        mpas_atmphys_interface_nhyd.o       \
+        mpas_atmphys_driver_convection_deep.o \
+        mpas_atmphys_driver_pbl.o             \
+        mpas_atmphys_driver_radiation_lw.o    \
+        mpas_atmphys_driver_radiation_sw.o    \
+        mpas_atmphys_driver_sfclayer.o        \
+        mpas_atmphys_constants.o              \
+        mpas_atmphys_interface_nhyd.o         \
+        mpas_atmphys_update.o                 \
         mpas_atmphys_vars.o
 
+mpas_atmphys_update.o: \
+        mpas_atmphys_driver_convection_deep.o \
+        mpas_atmphys_vars.o
+
 endif
 
 clean:

Modified: branches/atmos_physics/src/core_atmos_physics/mpas_atmphys_driver.F
===================================================================
--- branches/atmos_physics/src/core_atmos_physics/mpas_atmphys_driver.F        2012-03-01 23:17:44 UTC (rev 1568)
+++ branches/atmos_physics/src/core_atmos_physics/mpas_atmphys_driver.F        2012-03-01 23:20:34 UTC (rev 1569)
@@ -11,6 +11,7 @@
  use mpas_atmphys_driver_radiation_lw
  use mpas_atmphys_driver_sfclayer
  use mpas_atmphys_constants
+ use mpas_atmphys_update
  use mpas_atmphys_vars
 #ifdef non_hydrostatic_core
  use mpas_atmphys_interface_nhyd
@@ -84,8 +85,12 @@
                               block%diag_physics,block%sfc_input,block%tend_physics)
     endif
     if(l_camlw .and. config_radt_lw_scheme .eq. 'cam_lw') &amp;
-                              call radiation_camlw_to_MPAS(block%diag_physics) 
+                              call radiation_camlw_to_MPAS(block%diag_physics)
 
+    !call to accumulate long- and short-wave diagnostics if needed:
+     if(config_bucket_update /= 'none' .and. config_bucket_radt .gt. 0._RKIND) &amp;
+        call update_radiation_diagnostics(config_bucket_radt,block%mesh,block%diag_physics)
+
     !deallocate all radiation arrays:
     if(config_radt_sw_scheme.ne.'off' .or. config_radt_lw_scheme.ne.'off') &amp;
        call deallocate_cloudiness

Modified: branches/atmos_physics/src/core_atmos_physics/mpas_atmphys_driver_convection_deep.F
===================================================================
--- branches/atmos_physics/src/core_atmos_physics/mpas_atmphys_driver_convection_deep.F        2012-03-01 23:17:44 UTC (rev 1568)
+++ branches/atmos_physics/src/core_atmos_physics/mpas_atmphys_driver_convection_deep.F        2012-03-01 23:20:34 UTC (rev 1569)
@@ -413,12 +413,12 @@
  end subroutine convection_to_MPAS
 
 !=============================================================================================
- subroutine update_convection_deep(rainc_bucket,mesh,diag_physics)
+ subroutine update_convection_deep(bucket_rainc,mesh,diag_physics)
 !=============================================================================================
 
 !input arguments:
  type(mesh_type),intent(in):: mesh
- real(kind=RKIND),intent(in):: rainc_bucket
+ real(kind=RKIND),intent(in):: bucket_rainc
 
 !inout arguments:
  type(diag_physics_type),intent(inout):: diag_physics
@@ -433,11 +433,11 @@
     diag_physics % rainc % array(iCell) = diag_physics % rainc % array(iCell) &amp;
                                     + diag_physics % cuprec % array(iCell) * dt_dyn
 
-    if(l_rain .and. rainc_bucket.gt.0._RKIND .and. &amp;
-       diag_physics%rainc%array(iCell).gt.rainc_bucket) then
+    if(l_acrain .and. bucket_rainc.gt.0._RKIND .and. &amp;
+       diag_physics%rainc%array(iCell).gt.bucket_rainc) then
        diag_physics % i_rainc % array(iCell) = diag_physics % i_rainc % array(iCell) + 1
        diag_physics % rainc % array(iCell) = diag_physics % rainc % array(iCell) &amp;
-                                           - rainc_bucket
+                                           - bucket_rainc
     endif
 
  enddo

Modified: branches/atmos_physics/src/core_atmos_physics/mpas_atmphys_driver_microphysics.F
===================================================================
--- branches/atmos_physics/src/core_atmos_physics/mpas_atmphys_driver_microphysics.F        2012-03-01 23:17:44 UTC (rev 1568)
+++ branches/atmos_physics/src/core_atmos_physics/mpas_atmphys_driver_microphysics.F        2012-03-01 23:20:34 UTC (rev 1569)
@@ -374,11 +374,11 @@
  end subroutine precip_from_MPAS
 
 !=============================================================================================
- subroutine precip_to_MPAS(rainnc_bucket,diag_physics)
+ subroutine precip_to_MPAS(bucket_rainnc,diag_physics)
 !=============================================================================================
 
 !output variables:
- real(kind=RKIND),intent(in):: rainnc_bucket
+ real(kind=RKIND),intent(in):: bucket_rainnc
  type(diag_physics_type),intent(inout):: diag_physics
 
 !local variables:
@@ -398,10 +398,10 @@
     diag_physics % rainnc % array(i) = diag_physics % rainnc % array(i) &amp;
                                      + diag_physics % rainncv % array(i)
 
-   if(l_rain .and. rainnc_bucket.gt.0._RKIND .and. &amp;
-      diag_physics%rainnc%array(i).gt.rainnc_bucket) then
+   if(l_acrain .and. bucket_rainnc.gt.0._RKIND .and. &amp;
+      diag_physics%rainnc%array(i).gt.bucket_rainnc) then
       diag_physics % i_rainnc % array(i) = diag_physics % i_rainnc % array(i) + 1
-      diag_physics % rainnc % array(i) = diag_physics % rainnc % array(i) - rainnc_bucket
+      diag_physics % rainnc % array(i) = diag_physics % rainnc % array(i) - bucket_rainnc
    endif
  
  enddo

Modified: branches/atmos_physics/src/core_atmos_physics/mpas_atmphys_update.F
===================================================================
--- branches/atmos_physics/src/core_atmos_physics/mpas_atmphys_update.F        2012-03-01 23:17:44 UTC (rev 1568)
+++ branches/atmos_physics/src/core_atmos_physics/mpas_atmphys_update.F        2012-03-01 23:20:34 UTC (rev 1569)
@@ -4,10 +4,12 @@
  use mpas_grid_types
 
  use mpas_atmphys_driver_convection_deep
+ use mpas_atmphys_vars
 
  implicit none
  private
- public:: physics_update
+ public:: physics_update, &amp;
+          update_radiation_diagnostics
 
  contains
  
@@ -40,5 +42,117 @@
  end subroutine physics_update
 
 !=============================================================================================
+ subroutine update_radiation_diagnostics(bucket_radt,mesh,diag)
+!=============================================================================================
+
+!input arguments:
+ real(kind=RKIND),intent(in):: bucket_radt
+ type(mesh_type),intent(in):: mesh
+
+!inout arguments:
+ type(diag_physics_type),intent(inout):: diag
+
+!local variables:
+ integer:: iCell
+
+!--------------------------------------------------------------------------------------------
+
+ do iCell = 1, mesh%nCellsSolve
+    !short-wave radiation:
+    diag%acswdnb %array(iCell) = diag%acswdnb %array(iCell) + diag%swdnb %array(iCell)*dt_dyn
+    diag%acswdnbc%array(iCell) = diag%acswdnbc%array(iCell) + diag%swdnbc%array(iCell)*dt_dyn
+    diag%acswdnt %array(iCell) = diag%acswdnt %array(iCell) + diag%swdnt %array(iCell)*dt_dyn
+    diag%acswdntc%array(iCell) = diag%acswdntc%array(iCell) + diag%swdntc%array(iCell)*dt_dyn
+    diag%acswupb %array(iCell) = diag%acswupb %array(iCell) + diag%swupb %array(iCell)*dt_dyn
+    diag%acswupbc%array(iCell) = diag%acswupbc%array(iCell) + diag%swupbc%array(iCell)*dt_dyn
+    diag%acswupt %array(iCell) = diag%acswupt %array(iCell) + diag%swupt %array(iCell)*dt_dyn
+    diag%acswuptc%array(iCell) = diag%acswuptc%array(iCell) + diag%swuptc%array(iCell)*dt_dyn
+    !long-wave radiation:
+    diag%aclwdnb %array(iCell) = diag%aclwdnb %array(iCell) + diag%lwdnb %array(iCell)*dt_dyn
+    diag%aclwdnbc%array(iCell) = diag%aclwdnbc%array(iCell) + diag%lwdnbc%array(iCell)*dt_dyn
+    diag%aclwdnt %array(iCell) = diag%aclwdnt %array(iCell) + diag%lwdnt %array(iCell)*dt_dyn
+    diag%aclwdntc%array(iCell) = diag%aclwdntc%array(iCell) + diag%lwdntc%array(iCell)*dt_dyn
+    diag%aclwupb %array(iCell) = diag%aclwupb %array(iCell) + diag%lwupb %array(iCell)*dt_dyn
+    diag%aclwupbc%array(iCell) = diag%aclwupbc%array(iCell) + diag%lwupbc%array(iCell)*dt_dyn
+    diag%aclwupt %array(iCell) = diag%aclwupt %array(iCell) + diag%lwupt %array(iCell)*dt_dyn
+    diag%aclwuptc%array(iCell) = diag%aclwuptc%array(iCell) + diag%lwuptc%array(iCell)*dt_dyn
+ enddo
+
+ if(l_acradt .and. bucket_radt.gt.0._RKIND) then
+
+    do iCell = 1, mesh%nCellsSolve
+       !short-wave radiation:
+       if(diag%acswdnb%array(iCell) .gt. bucket_radt) then
+          diag%i_acswdnb%array(iCell) = diag%i_acswdnb%array(iCell) + 1
+          diag%acswdnb%array(iCell) = diag%acswdnb%array(iCell) - bucket_radt
+       endif   
+       if(diag%acswdnbc%array(iCell) .gt. bucket_radt) then
+          diag%i_acswdnbc%array(iCell) = diag%i_acswdnbc%array(iCell) + 1
+          diag%acswdnbc%array(iCell) = diag%acswdnbc%array(iCell) - bucket_radt
+       endif   
+       if(diag%acswdnt%array(iCell) .gt. bucket_radt) then
+          diag%i_acswdnt%array(iCell) = diag%i_acswdnt%array(iCell) + 1
+          diag%acswdnt%array(iCell) = diag%acswdnt%array(iCell) - bucket_radt
+       endif   
+       if(diag%acswdntc%array(iCell) .gt. bucket_radt) then
+          diag%i_acswdntc%array(iCell) = diag%i_acswdntc%array(iCell) + 1
+          diag%acswdntc%array(iCell) = diag%acswdntc%array(iCell) - bucket_radt
+       endif
+       if(diag%acswupb%array(iCell) .gt. bucket_radt) then
+          diag%i_acswupb%array(iCell) = diag%i_acswupb%array(iCell) + 1
+          diag%acswupb%array(iCell) = diag%acswupb%array(iCell) - bucket_radt
+       endif   
+       if(diag%acswupbc%array(iCell) .gt. bucket_radt) then
+          diag%i_acswupbc%array(iCell) = diag%i_acswupbc%array(iCell) + 1
+          diag%acswupbc%array(iCell) = diag%acswupbc%array(iCell) - bucket_radt
+       endif   
+       if(diag%acswupt%array(iCell) .gt. bucket_radt) then
+          diag%i_acswupt%array(iCell) = diag%i_acswupt%array(iCell) + 1
+          diag%acswupt%array(iCell) = diag%acswupt%array(iCell) - bucket_radt
+       endif   
+       if(diag%acswuptc%array(iCell) .gt. bucket_radt) then
+          diag%i_acswuptc%array(iCell) = diag%i_acswuptc%array(iCell) + 1
+          diag%acswuptc%array(iCell) = diag%acswuptc%array(iCell) - bucket_radt
+       endif
+       !long-wave radiation:
+       if(diag%aclwdnb%array(iCell) .gt. bucket_radt) then
+          diag%i_aclwdnb%array(iCell) = diag%i_aclwdnb%array(iCell) + 1
+          diag%aclwdnb%array(iCell) = diag%aclwdnb%array(iCell) - bucket_radt
+       endif   
+       if(diag%aclwdnbc%array(iCell) .gt. bucket_radt) then
+          diag%i_aclwdnbc%array(iCell) = diag%i_aclwdnbc%array(iCell) + 1
+          diag%aclwdnbc%array(iCell) = diag%aclwdnbc%array(iCell) - bucket_radt
+       endif   
+       if(diag%aclwdnt%array(iCell) .gt. bucket_radt) then
+          diag%i_aclwdnt%array(iCell) = diag%i_aclwdnt%array(iCell) + 1
+          diag%aclwdnt%array(iCell) = diag%aclwdnt%array(iCell) - bucket_radt
+       endif   
+       if(diag%aclwdntc%array(iCell) .gt. bucket_radt) then
+          diag%i_aclwdntc%array(iCell) = diag%i_aclwdntc%array(iCell) + 1
+          diag%aclwdntc%array(iCell) = diag%aclwdntc%array(iCell) - bucket_radt
+       endif
+       if(diag%aclwupb%array(iCell) .gt. bucket_radt) then
+          diag%i_aclwupb%array(iCell) = diag%i_aclwupb%array(iCell) + 1
+          diag%aclwupb%array(iCell) = diag%aclwupb%array(iCell) - bucket_radt
+       endif   
+       if(diag%aclwupbc%array(iCell) .gt. bucket_radt) then
+          diag%i_aclwupbc%array(iCell) = diag%i_aclwupbc%array(iCell) + 1
+          diag%aclwupbc%array(iCell) = diag%aclwupbc%array(iCell) - bucket_radt
+       endif   
+       if(diag%aclwupt%array(iCell) .gt. bucket_radt) then
+          diag%i_aclwupt%array(iCell) = diag%i_aclwupt%array(iCell) + 1
+          diag%aclwupt%array(iCell) = diag%aclwupt%array(iCell) - bucket_radt
+       endif   
+       if(diag%aclwuptc%array(iCell) .gt. bucket_radt) then
+          diag%i_aclwuptc%array(iCell) = diag%i_aclwuptc%array(iCell) + 1
+          diag%aclwuptc%array(iCell) = diag%aclwuptc%array(iCell) - bucket_radt
+       endif
+    enddo
+
+ endif
+
+ end subroutine update_radiation_diagnostics
+
+!=============================================================================================
  end module mpas_atmphys_update
 !=============================================================================================

Modified: branches/atmos_physics/src/core_atmos_physics/mpas_atmphys_vars.F
===================================================================
--- branches/atmos_physics/src/core_atmos_physics/mpas_atmphys_vars.F        2012-03-01 23:17:44 UTC (rev 1568)
+++ branches/atmos_physics/src/core_atmos_physics/mpas_atmphys_vars.F        2012-03-01 23:20:34 UTC (rev 1569)
@@ -29,7 +29,8 @@
  logical:: l_radtlw                   !controls call to longwave radiation parameterization.
  logical:: l_radtsw                   !controls call to shortwave radiationparameterization.
  logical:: l_camlw                    !controls when to save local CAM LW abs and ems arrays.
- logical:: l_rain                     !when .true., limit to accumulated rain is applied.
+ logical:: l_acrain                   !when .true., limit to accumulated rain is applied.
+ logical:: l_acradt                   !when .true., limit to lw and sw radiation is applied.
 
  integer,public:: ids,ide,jds,jde,kds,kde
  integer,public:: ims,ime,jms,jme,kms,kme

</font>
</pre>