** Warning **

Issuing rollback() due to DESTROY without explicit disconnect() of DBD::mysql::db handle dbname=MITgcm at /usr/local/share/lxr/lib/LXR/Common.pm line 1224.

Last-Modified: Wed, 21 May 2024 05:11:33 GMT Content-Type: text/html; charset=utf-8 MITgcm/MITgcm/pkg/atm_phys/atm_phys_driver.F
Back to home page

MITgcm

 
 

    


File indexing completed on 2019-01-25 06:10:04 UTC

view on githubraw file Latest commit 88391fb6 on 2019-01-24 19:38:27 UTC
b2ea1d2979 Jean*0001 #include "ATM_PHYS_OPTIONS.h"
                0002 
                0003 CBOP
                0004 C !ROUTINE: ATM_PHYS_DRIVER
                0005 
                0006 C !INTERFACE: ==========================================================
                0007       SUBROUTINE ATM_PHYS_DRIVER(
                0008      I                    myTime, myIter, myThid )
                0009 
                0010 C !DESCRIPTION:
                0011 C Calculate custom tendency terms outside k-loop in DO_OCEANIC_PHYS
                0012 
                0013 C !USES: ===============================================================
                0014       use radiation_mod
                0015       use lscale_cond_mod
                0016       use dargan_bettsmiller_mod
                0017       use surface_flux_mod
                0018       use vert_turb_driver_mod
                0019       use vert_diff_mod, only: gcm_vert_diff_down,
                0020      &                         gcm_vert_diff_up,
                0021      &                         surf_diff_type
                0022       use mixed_layer_mod, only: mixed_layer
6c614dc1a3 Jean*0023       use constants_mod, only:  HLv
b2ea1d2979 Jean*0024 
                0025       IMPLICIT NONE
                0026 #include "SIZE.h"
                0027 #include "EEPARAMS.h"
                0028 #include "PARAMS.h"
                0029 #include "GRID.h"
                0030 #include "SURFACE.h"
6c614dc1a3 Jean*0031 #include "FFIELDS.h"
b2ea1d2979 Jean*0032 #include "ATM_PHYS_PARAMS.h"
                0033 #include "ATM_PHYS_VARS.h"
                0034 
                0035 C !INPUT PARAMETERS: ===================================================
e259fb67de Jean*0036 C  myTime   :: Current time in simulation
                0037 C  myIter   :: Current time-step number
                0038 C  myThid   :: my Thread Id number
b2ea1d2979 Jean*0039       _RL     myTime
                0040       INTEGER myIter, myThid
                0041 
                0042 C !OUTPUT PARAMETERS: ==================================================
                0043 
                0044 C !LOCAL VARIABLES: ====================================================
e259fb67de Jean*0045 C  bi,bj    :: Tile indices
                0046 C  lat2d    :: latitude of grid-cell center          [rad]
                0047 C pHalf3d   :: pressure at interface between 2 levels [Pa]
                0048 C pFull3d   :: pressure at level center               [Pa]
                0049 C zHalf3d   :: height of interface between 2 levels   [m]
                0050 C zFull3d   :: height of level center                 [m]
                0051 C  t3d      :: absolute temperature                   [K]
                0052 C  q3d      :: specific humidity                    [kg/kg]
                0053 C  u3d      :: wind speed, 1rst component (X-dir)    [m/s]
                0054 C  v3d      :: wind speed, 2nd  component (Y-dir)    [m/s]
b2ea1d2979 Jean*0055       INTEGER bi, bj
                0056       _RL lat2d   (sNx,sNy)
                0057       _RL pHalf3d (sNx,sNy,Nr+1)
                0058       _RL pFull3d (sNx,sNy,Nr)
                0059       _RL zHalf3d (sNx,sNy,Nr+1)
                0060       _RL zFull3d (sNx,sNy,Nr)
                0061       _RL t3d     (sNx,sNy,Nr)
                0062       _RL q3d     (sNx,sNy,Nr)
                0063       _RL u3d     (sNx,sNy,Nr)
                0064       _RL v3d     (sNx,sNy,Nr)
                0065       _RL tdt3d   (sNx,sNy,Nr)
                0066       _RL qdt3d   (sNx,sNy,Nr)
                0067       _RL udt3d   (sNx,sNy,Nr)
                0068       _RL vdt3d   (sNx,sNy,Nr)
e259fb67de Jean*0069 CEOP
b2ea1d2979 Jean*0070       _RL s_sw_dwn(sNx,sNy)
                0071       _RL s_lw_dwn(sNx,sNy)
                0072       _RL t_surf  (sNx,sNy)
                0073 C-- radiation fields:
                0074       _RL albedo_2d (sNx,sNy)
88391fb671 jm-c 0075       _RL ozone_3d  (sNx,sNy,Nr)
                0076       _RL ozoneColmn(sNx,sNy,Nr+1)
b2ea1d2979 Jean*0077       _RL dtrans_3d (sNx,sNy,Nr)
a690051f17 Jean*0078       _RL dtrans_win(sNx,sNy,Nr)
b2ea1d2979 Jean*0079       _RL b_3d      (sNx,sNy,Nr)
a690051f17 Jean*0080       _RL b_win     (sNx,sNy,Nr)
b2ea1d2979 Jean*0081       _RL lw_down_3d(sNx,sNy,Nr+1)
                0082       _RL sw_down_3d(sNx,sNy,Nr+1)
a690051f17 Jean*0083       _RL sw_net_3d (sNx,sNy,Nr+1)
                0084       _RL lw_net_3d (sNx,sNy,Nr+1)
                0085       _RL adj_lw_up (sNx,sNy)
                0086       _RL rad_dt_tg (sNx,sNy,Nr)
                0087 
b2ea1d2979 Jean*0088       LOGICAL coldT (sNx,sNy)
                0089 C-- output from convection & LSC :
                0090       _RL t3d_tmp    (sNx,sNy,Nr)
                0091       _RL q3d_tmp    (sNx,sNy,Nr)
                0092       _RL cond_dt_tg (sNx,sNy,Nr)
                0093       _RL cond_dt_qg (sNx,sNy,Nr)
                0094       _RL rain2d     (sNx,sNy)
                0095       _RL snow2d     (sNx,sNy)
                0096       _RL q_ref      (sNx,sNy,Nr)
                0097       _RL t_ref      (sNx,sNy,Nr)
                0098       _RL bmflag     (sNx,sNy)
                0099       _RL klzbs      (sNx,sNy)
                0100       _RL cape       (sNx,sNy)
                0101       _RL cin        (sNx,sNy)
                0102       _RL invtau_bm_t(sNx,sNy)
                0103       _RL invtau_bm_q(sNx,sNy)
                0104       _RL capeflag   (sNx,sNy)
                0105 C-- Input/Output for surface flux:
                0106       _RL  q_surf(sNx,sNy)
                0107       _RL  u_surf(sNx,sNy), v_surf(sNx,sNy)
                0108       _RL  rough_mom(sNx,sNy), rough_heat(sNx,sNy)
                0109       _RL  rough_moist(sNx,sNy), gust(sNx,sNy)
                0110       _RL  flux_t(sNx,sNy), flux_q(sNx,sNy), flux_r(sNx,sNy)
                0111       _RL  flux_u(sNx,sNy), flux_v(sNx,sNy)
                0112       _RL  drag_m(sNx,sNy), drag_t(sNx,sNy), drag_q(sNx,sNy)
                0113       _RL  w_atm(sNx,sNy)
                0114       _RL  ustar(sNx,sNy), bstar(sNx,sNy), qstar(sNx,sNy)
                0115       _RL  dhdt_surf(sNx,sNy), dedt_surf(sNx,sNy), dedq_surf(sNx,sNy)
                0116       _RL  drdt_surf(sNx,sNy)
                0117       _RL  dhdt_atm(sNx,sNy), dedq_atm(sNx,sNy), dtaudv_atm(sNx,sNy)
                0118       LOGICAL land(sNx,sNy), avail(sNx,sNy)
                0119 C-- Input for turb:
                0120       _RL fracland(sNx,sNy)
                0121       _RL rough(sNx,sNy)
                0122 C-- Output from turbulence driver:
                0123       _RL  diff_t(sNx,sNy,Nr), diff_m(sNx,sNy,Nr)
                0124       _RL diff_dt_tg (sNx,sNy,Nr)
                0125       _RL diff_dt_qg (sNx,sNy,Nr)
                0126       _RL diss_heat  (sNx,sNy,Nr)
                0127 c     TYPE(surf_diff_type) :: tri_surf ! used by gcm_vert_diff
                0128       _RL tri_surf_dtmass(sNx,sNy)
                0129       _RL tri_surf_dflux_t(sNx,sNy), tri_surf_dflux_q(sNx,sNy)
                0130       _RL tri_surf_delta_t(sNx,sNy), tri_surf_delta_q(sNx,sNy)
                0131       _RL e_global(sNx,sNy,Nr-1)
                0132       _RL f_t_global(sNx,sNy,Nr-1), f_q_global(sNx,sNy,Nr-1)
                0133 C-- Mixed Layer fields:
                0134       _RL ocean_qflux(sNx,sNy)
c9694dc201 Jean*0135       _RL mixLayDepth(sNx,sNy)
b2ea1d2979 Jean*0136       _RL delta_t_surf(sNx,sNy)
                0137 C--
55a26a1b95 Jean*0138       _RL dpFac(sNx,sNy)
e259fb67de Jean*0139       _RL conv_T2theta
                0140       INTEGER k, kc
b2ea1d2979 Jean*0141 c     INTEGER ioUnit
                0142 c     _RS     dummyRS(1)
6c614dc1a3 Jean*0143 #ifdef COMPONENT_MODULE
31f8711b60 Jean*0144       INTEGER i, j
6c614dc1a3 Jean*0145       _RL taux(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0146       _RL tauy(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0147 #endif /* COMPONENT_MODULE */
b2ea1d2979 Jean*0148 
                0149 C--   Loops on tile indices bi,bj
                0150       DO bj=myByLo(myThid),myByHi(myThid)
                0151        DO bi=myBxLo(myThid),myBxHi(myThid)
                0152 
e259fb67de Jean*0153 C--   Get surface temperature
                0154         t_surf(:,:) = atmPhys_SST(1:sNx,1:sNy,bi,bj)
6c614dc1a3 Jean*0155 #ifdef COMPONENT_MODULE
e259fb67de Jean*0156         IF ( useCoupler ) THEN
6c614dc1a3 Jean*0157 C--   take surface data from the ocean component
                0158 C     to replace MxL fields (if use sea-ice) or directly atm_phys SST
e259fb67de Jean*0159           CALL ATM_APPLY_IMPORT(
                0160      I             maskInC,
                0161      U             t_surf, taux,
                0162 c    U             sst1(1,myThid), oice1,
                0163      I             myTime, myIter, bi, bj, myThid )
                0164         ENDIF
6c614dc1a3 Jean*0165 #endif /* COMPONENT_MODULE */
                0166 #ifdef ALLOW_DIAGNOSTICS
e259fb67de Jean*0167         IF ( useDiagnostics ) THEN
6c614dc1a3 Jean*0168 C--   fill-up Atm_Phys diagnostics for state variables
                0169          CALL DIAGNOSTICS_FILL( t_surf, 'AtPh_SST',
                0170      &                          0, 1, 3, bi, bj, myThid )
e259fb67de Jean*0171         ENDIF
6c614dc1a3 Jean*0172 #endif /* ALLOW_DIAGNOSTICS */
e259fb67de Jean*0173 c       IF ( myIter.EQ.nIter0 ) THEN
                0174 c        ioUnit = 0
                0175 c        CALL MDS_WRITEVEC_LOC(
                0176 c    I                     'SST.check', writeBinaryPrec, ioUnit,
                0177 c    I                     'RL', sNx*sNy, t_surf, dummyRS,
                0178 c    I                     bi, bj, 1, myIter, myThid )
                0179 c       ENDIF
                0180         ocean_qflux(:,:) = atmPhys_Qflx(1:sNx,1:sNy,bi,bj)
c9694dc201 Jean*0181         mixLayDepth(:,:) = atmPhys_MxLD(1:sNx,1:sNy,bi,bj)
                0182         albedo_2d(:,:) = atmPhys_Albedo(1:sNx,1:sNy,bi,bj)
88391fb671 jm-c 0183 C--   copy atmPhys_Ozone to ozone_3d for rad scheme
                0184         DO k=1,Nr
                0185           kc = Nr-k+1
                0186           ozone_3d(:,:,k) = atmPhys_Ozone(1:sNx,1:sNy,kc,bi,bj)
                0187         ENDDO
e259fb67de Jean*0188 
                0189 C--   Get grid and dynamical fields from main model common blocks
                0190         CALL ATM_PHYS_DYN2PHYS(
                0191      O                    lat2d, pHalf3d, pFull3d,
                0192      O                    zHalf3d, zFull3d,
                0193      O                    t3d, q3d, u3d, v3d,
                0194      I                    bi, bj, myTime, myIter, myThid )
                0195 
                0196 C--   initialise
                0197         coldT(:,:) = .FALSE.
                0198         EmPmR(:,:,bi,bj) = 0.
                0199         rain2d(:,:) = 0.
                0200         snow2d(:,:) = 0.
b2ea1d2979 Jean*0201 
                0202 C--   initialise tendency
e259fb67de Jean*0203         tdt3d = 0.
                0204         qdt3d = 0.
                0205         udt3d = 0.
                0206         vdt3d = 0.
                0207         cond_dt_tg = 0.
                0208         cond_dt_qg = 0.
                0209 
                0210         IF (lwet_convection) THEN
                0211           CALL DARGAN_BETTSMILLER(
b2ea1d2979 Jean*0212      I                   deltaT, t3d, q3d, pFull3d, pHalf3d, coldT,
                0213      O                   rain2d, snow2d, cond_dt_tg, cond_dt_qg,
                0214      O                   q_ref,      bmflag,
                0215      O                   klzbs,        cape,
                0216      O                   cin,         t_ref,
                0217      O                   invtau_bm_t,     invtau_bm_q,
                0218      O                   capeflag,
                0219      I                   bi,bj,myIter,myThid )
                0220 
                0221 C-    store updated state to pass to LSC (addition from POG)
e259fb67de Jean*0222           t3d_tmp = t3d + cond_dt_tg
                0223           q3d_tmp = q3d + cond_dt_qg
                0224 
                0225           cond_dt_tg = cond_dt_tg / deltaT
                0226           cond_dt_qg = cond_dt_qg / deltaT
                0227           rain2d = rain2d / deltaT
                0228           EmPmR(1:sNx,1:sNy,bi,bj) = -rain2d(:,:)
                0229 
                0230           tdt3d = tdt3d + cond_dt_tg
                0231           qdt3d = qdt3d + cond_dt_qg
                0232 
                0233 #ifdef ALLOW_DIAGNOSTICS
                0234           IF ( useDiagnostics ) THEN
                0235            CALL DIAGNOSTICS_FILL( rain2d , 'AtPhCnvP',
                0236      &                            0, 1, 3, bi, bj, myThid )
                0237            CALL DIAGNOSTICS_FILL( cape   , 'AtPhCAPE',
                0238      &                            0, 1, 3, bi, bj, myThid )
                0239            CALL DIAGNOSTICS_FILL( cin    , 'AtPhCnIn',
                0240      &                            0, 1, 3, bi, bj, myThid )
                0241            CALL DIAGNOSTICS_FILL( klzbs  , 'AtPhKlzb',
                0242      &                            0, 1, 3, bi, bj, myThid )
                0243            CALL DIAGNOSTICS_FILL( bmflag , 'AtPhConv',
                0244      &                            0, 1, 3, bi, bj, myThid )
                0245            CALL DIAGNOSTICS_FILL( invtau_bm_t, 'AtPhRlxT',
                0246      &                            0, 1, 3, bi, bj, myThid )
                0247            CALL DIAGNOSTICS_FILL( invtau_bm_q, 'AtPhRlxQ',
                0248      &                            0, 1, 3, bi, bj, myThid )
                0249            CALL DIAGNOSTICS_FILL( t_ref  , 'AtPh_Trf',
                0250      &                           -1,Nr, 3, bi, bj, myThid )
                0251            CALL DIAGNOSTICS_FILL( q_ref  , 'AtPh_Qrf',
                0252      &                           -1,Nr, 3, bi, bj, myThid )
a690051f17 Jean*0253           CALL DIAGNOSTICS_FILL( cond_dt_tg, 'AtPhdTcv',
                0254      &                           -1,Nr, 3, bi, bj, myThid )
e259fb67de Jean*0255 c          snow2d = snow2d / deltaT
                0256 c          CALL DIAGNOSTICS_FILL( snow2d , 'SDIAG1  ',
                0257 c    &                            0, 1, 3, bi, bj, myThid )
                0258           ENDIF
                0259 #endif /* ALLOW_DIAGNOSTICS */
                0260         ELSE
                0261           t3d_tmp = t3d
                0262           q3d_tmp = q3d
                0263         ENDIF
b2ea1d2979 Jean*0264 
e259fb67de Jean*0265         cond_dt_tg = 0.
                0266         cond_dt_qg = 0.
                0267         rain2d(:,:) = 0.
                0268         CALL LSCALE_COND(
                0269      I              t3d_tmp, q3d_tmp, pFull3d, pHalf3d, coldT,
                0270      O              rain2d, snow2d, cond_dt_tg, cond_dt_qg, q_ref,
                0271      I              myThid )
b2ea1d2979 Jean*0272         cond_dt_tg = cond_dt_tg / deltaT
                0273         cond_dt_qg = cond_dt_qg / deltaT
                0274         rain2d = rain2d / deltaT
e259fb67de Jean*0275         EmPmR(1:sNx,1:sNy,bi,bj) = EmPmR(1:sNx,1:sNy,bi,bj)
                0276      &                           - rain2d(:,:)
b2ea1d2979 Jean*0277 
                0278         tdt3d = tdt3d + cond_dt_tg
                0279         qdt3d = qdt3d + cond_dt_qg
                0280 
                0281 #ifdef ALLOW_DIAGNOSTICS
                0282         IF ( useDiagnostics ) THEN
e259fb67de Jean*0283            CALL DIAGNOSTICS_FILL( rain2d , 'AtPhLscP',
                0284      &                            0, 1, 3, bi, bj, myThid )
b2ea1d2979 Jean*0285 C     Re-use "q_ref" array to calculate Relative Humidity (in %)
88391fb671 jm-c 0286            q_ref = 100. _d 0 * q3d_tmp/MAX( q_ref, 1. _d -12 )
e259fb67de Jean*0287            CALL DIAGNOSTICS_FILL( q_ref , 'RELHUM  ',
                0288      &                           -1, Nr, 3, bi, bj, myThid )
                0289         ENDIF
b2ea1d2979 Jean*0290 #endif /* ALLOW_DIAGNOSTICS */
                0291 
e259fb67de Jean*0292         IF ( two_stream ) THEN
                0293           CALL RADIATION_DOWN(
                0294      I                   sNx,sNy, myTime, lat2d, pHalf3d, t3d, q3d,
88391fb671 jm-c 0295      I                   albedo_2d, ozone_3d,
                0296      O                   ozoneColmn, s_sw_dwn, s_lw_dwn,
a690051f17 Jean*0297      O                   dtrans_3d, dtrans_win, b_3d, b_win,
e259fb67de Jean*0298      O                   lw_down_3d, sw_down_3d, myThid )
b2ea1d2979 Jean*0299 
e259fb67de Jean*0300 c         write(errorMessageUnit,'(A,I5,I3,4F9.3)')
                0301 c    &      'it,bj,sw,lw:',myIter,bj,
b2ea1d2979 Jean*0302 c    &      sw_down_3d(1,1,15), sw_down_3d(1,sNy,15),
                0303 c    &      lw_down_3d(1,1,15), lw_down_3d(1,sNy,15)
e259fb67de Jean*0304         ENDIF
b2ea1d2979 Jean*0305 
e259fb67de Jean*0306         IF (.TRUE.) THEN
                0307          land = .false.
                0308          avail = .true.
                0309          rough_mom = roughness_mom
                0310          rough_heat = roughness_heat
                0311          rough_moist = roughness_moist
                0312          gust = 1.0
                0313          u_surf = 0.
                0314          v_surf = 0.
                0315          CALL SURFACE_FLUX(
b2ea1d2979 Jean*0316      I            t3d(:,:,Nr), q3d(:,:,Nr), u3d(:,:,Nr), v3d(:,:,Nr),
                0317      I            pFull3d(:,:,Nr), zFull3d(:,:,Nr), pHalf3d(:,:,Nr+1),
                0318      I            t_surf, t_surf,
                0319      U            q_surf,
                0320      I            u_surf, v_surf,
                0321      I            rough_mom, rough_heat, rough_moist, gust,
                0322      O            flux_t, flux_q, flux_r, flux_u, flux_v,
                0323      O            drag_m, drag_t, drag_q, w_atm,
                0324      O            ustar, bstar, qstar,
                0325      O            dhdt_surf, dedt_surf, dedq_surf, drdt_surf,
                0326      O            dhdt_atm, dedq_atm, dtaudv_atm,
                0327      I            deltaT, land(:,:), avail(:,:), myThid  )
e259fb67de Jean*0328         ENDIF
b2ea1d2979 Jean*0329 
e259fb67de Jean*0330         IF ( two_stream ) THEN
a690051f17 Jean*0331           rad_dt_tg = tdt3d
e259fb67de Jean*0332           CALL RADIATION_UP(
                0333      I                   sNx,sNy, myTime, lat2d, pHalf3d, t_surf, t3d,
a690051f17 Jean*0334      U                   tdt3d, lw_net_3d, sw_net_3d,
88391fb671 jm-c 0335      I                   albedo_2d, ozoneColmn, dtrans_3d, dtrans_win,
a690051f17 Jean*0336      I                   b_3d, b_win, lw_down_3d, sw_down_3d, myThid )
b2ea1d2979 Jean*0337 
e259fb67de Jean*0338 c         write(errorMessageUnit,'(A,I5,I3,4F9.3)')
                0339 c    &      'it,bj,t,tdt:',myIter,bj,
b2ea1d2979 Jean*0340 c    &      t3d(1,1,15), t3d(1,sNy,15),
                0341 c    &      tdt3d(1,1,15)*86400., tdt3d(1,sNy,15)*86400.
a690051f17 Jean*0342           rad_dt_tg = tdt3d - rad_dt_tg
                0343         ELSE
                0344           rad_dt_tg = 0.
e259fb67de Jean*0345         ENDIF
b2ea1d2979 Jean*0346 
e259fb67de Jean*0347         IF (turb) THEN
                0348          fracland = 0.0
                0349          rough = roughness_mom
                0350          CALL VERT_TURB_DRIVER(  1, 1, myTime, myTime+deltaT, deltaT,
                0351      I             fracland(:,:), pHalf3d, pFull3d, zHalf3d, zFull3d,
                0352      I             ustar, bstar, rough,
                0353      I             u3d, v3d, t3d, q3d,
                0354      O             diff_t(:,:,:), diff_m(:,:,:), gust(:,:),
                0355      I             myThid )
                0356         ENDIF
b2ea1d2979 Jean*0357 
e259fb67de Jean*0358         diff_dt_tg = tdt3d
                0359         diff_dt_qg = qdt3d
b2ea1d2979 Jean*0360 
e259fb67de Jean*0361         CALL GCM_VERT_DIFF_DOWN( 1, 1, deltaT,
b2ea1d2979 Jean*0362      I           u3d, v3d, t3d, q3d,
                0363      I           diff_m(:,:,:), diff_t(:,:,:),
                0364      I           pHalf3d, pFull3d, zFull3d,
                0365      U           flux_u(:,:),  flux_v(:,:),  dtaudv_atm,
                0366      U           udt3d, vdt3d, tdt3d,
                0367      I           qdt3d,
                0368      O           diss_heat(:,:,:),
                0369      U           tri_surf_dtmass,
                0370      U           tri_surf_dflux_t, tri_surf_dflux_q,
                0371      U           tri_surf_delta_t, tri_surf_delta_q,
                0372      O           e_global, f_t_global, f_q_global,
                0373      I           myThid )
                0374 
e259fb67de Jean*0375         CALL MIXED_LAYER(
                0376      I                    myTime,
                0377      U                    t_surf(:,:),
37a6d89494 Jean*0378      U                    flux_t(:,:), flux_q(:,:), flux_r(:,:),
e259fb67de Jean*0379      I                    deltaT,
37a6d89494 Jean*0380      I                    s_sw_dwn(:,:), s_lw_dwn(:,:),
e259fb67de Jean*0381      U                    tri_surf_dtmass,
                0382      U                    tri_surf_dflux_t, tri_surf_dflux_q,
                0383      U                    tri_surf_delta_t, tri_surf_delta_q,
37a6d89494 Jean*0384      I                    dhdt_surf(:,:), dedt_surf(:,:),
                0385      I                    dedq_surf(:,:), drdt_surf(:,:),
                0386      I                    dhdt_atm(:,:),  dedq_atm(:,:),
c9694dc201 Jean*0387      I                    ocean_qflux, mixLayDepth,
e259fb67de Jean*0388      O                    delta_t_surf(:,:),
                0389      I                    myThid )
                0390 
                0391         CALL GCM_VERT_DIFF_UP ( 1, 1, deltat,
b2ea1d2979 Jean*0392      I           tri_surf_delta_t, tri_surf_delta_q,
                0393      I           e_global, f_t_global, f_q_global,
                0394      O           tdt3d, qdt3d,
                0395      I           myThid )
                0396 
e259fb67de Jean*0397         diff_dt_tg = tdt3d - diff_dt_tg
                0398         diff_dt_qg = qdt3d - diff_dt_qg
b2ea1d2979 Jean*0399 
a690051f17 Jean*0400 C-- update Upward LW (assuming adjustment of surf LW emission goes directly to space)
                0401         adj_lw_up = drdt_surf*delta_t_surf
                0402         DO k=1,Nr+1
                0403          lw_net_3d(:,:,k) = lw_net_3d(:,:,k) + adj_lw_up
                0404         ENDDO
b2ea1d2979 Jean*0405 
e259fb67de Jean*0406         DO k=1,Nr
                0407           kc = Nr-k+1
                0408           conv_T2theta = (atm_po/rC(kc))**atm_kappa
17c247df15 Jean*0409 #ifdef NONLIN_FRSURF
e259fb67de Jean*0410           IF ( select_rStar.GE.1 ) THEN
                0411            atmPhys_dT(1:sNx,1:sNy,kc,bi,bj) = tdt3d(:,:,k)
17c247df15 Jean*0412      &                   * conv_T2theta/pStarFacK(1:sNx,1:sNy,bi,bj)
e259fb67de Jean*0413           ELSE
17c247df15 Jean*0414 #else /* NONLIN_FRSURF */
e259fb67de Jean*0415           IF ( .TRUE. ) THEN
17c247df15 Jean*0416 #endif /* NONLIN_FRSURF */
e259fb67de Jean*0417            atmPhys_dT(1:sNx,1:sNy,kc,bi,bj) = tdt3d(:,:,k)
17c247df15 Jean*0418      &                   * conv_T2theta
e259fb67de Jean*0419           ENDIF
                0420           atmPhys_dQ(1:sNx,1:sNy,kc,bi,bj)  = qdt3d(:,:,k)
55a26a1b95 Jean*0421 C-    Note: multiply A-grid tend. by hFacC (<-> dpFac) here and
                0422 C     C-grid average will be divided by hFacW,S when applied to Dynamics
                0423           dpFac(:,:) = ( pHalf3d(:,:,k+1) - pHalf3d(:,:,k)
                0424      &                 )*recip_drF(kc)
                0425           atmPhys_dU(1:sNx,1:sNy,kc,bi,bj) = udt3d(:,:,k)*dpFac(:,:)
                0426           atmPhys_dV(1:sNx,1:sNy,kc,bi,bj) = vdt3d(:,:,k)*dpFac(:,:)
e259fb67de Jean*0427         ENDDO
b2ea1d2979 Jean*0428 
                0429 C--   Update SST
e259fb67de Jean*0430         IF ( atmPhys_stepSST ) THEN
                0431           atmPhys_SST(1:sNx,1:sNy,bi,bj) = t_surf(:,:)
                0432         ENDIF
b2ea1d2979 Jean*0433 
e259fb67de Jean*0434         EmPmR(1:sNx,1:sNy,bi,bj) = EmPmR(1:sNx,1:sNy,bi,bj)
                0435      &                           + flux_q(:,:)
                0436         Qnet(1:sNx,1:sNy,bi,bj) = flux_t(:,:) + flux_r(:,:)
                0437      &                          - s_lw_dwn(:,:) - s_sw_dwn(:,:)
                0438      &                          + flux_q(:,:)*HLv
                0439         Qsw (1:sNx,1:sNy,bi,bj) = -s_sw_dwn(:,:)
6c614dc1a3 Jean*0440 #ifdef COMPONENT_MODULE
                0441 C-    surface wind-stress (+ = down) at grid cell center (A-grid):
e259fb67de Jean*0442         taux(1:sNx,1:sNy,bi,bj) = -flux_u
                0443         tauy(1:sNx,1:sNy,bi,bj) = -flux_v
6c614dc1a3 Jean*0444 #endif /* COMPONENT_MODULE */
                0445 
b2ea1d2979 Jean*0446 #ifdef ALLOW_DIAGNOSTICS
e259fb67de Jean*0447         IF ( useDiagnostics ) THEN
                0448           CALL DIAGNOSTICS_FILL( atmPhys_dT , 'AtPhdTdt',
                0449      &                           0, Nr, 1, bi, bj, myThid )
                0450           CALL DIAGNOSTICS_FILL( atmPhys_dQ , 'AtPhdQdt',
                0451      &                           0, Nr, 1, bi, bj, myThid )
ea820123c4 Jean*0452 c         CALL DIAGNOSTICS_FILL( atmPhys_dU , 'AtPhdUdt',
                0453 c    &                           0, Nr, 1, bi, bj, myThid )
                0454 c         CALL DIAGNOSTICS_FILL( atmPhys_dV , 'AtPhdVdt',
                0455 c    &                           0, Nr, 1, bi, bj, myThid )
                0456           CALL DIAGNOSTICS_FILL( udt3d , 'AtPhdUdt',
                0457      &                           -1, Nr, 3, bi, bj, myThid )
                0458           CALL DIAGNOSTICS_FILL( vdt3d , 'AtPhdVdt',
                0459      &                           -1, Nr, 3, bi, bj, myThid )
e259fb67de Jean*0460           CALL DIAGNOSTICS_FILL( diff_t , 'AtPhDifT',
                0461      &                           -1, Nr, 3, bi, bj, myThid )
                0462           CALL DIAGNOSTICS_FILL( diff_m , 'AtPhDifM',
                0463      &                           -1, Nr, 3, bi, bj, myThid )
a690051f17 Jean*0464           CALL DIAGNOSTICS_FILL( rad_dt_tg , 'AtPhdTrd',
                0465      &                           -1, Nr, 3, bi, bj, myThid )
                0466           CALL DIAGNOSTICS_FILL( sw_net_3d , 'AtPhNSR ',
                0467      &                           -1, Nr+1, 3, bi, bj, myThid )
                0468           CALL DIAGNOSTICS_FILL( lw_net_3d , 'AtPhNLR ',
                0469      &                           -1, Nr+1, 3, bi, bj, myThid )
                0470           CALL DIAGNOSTICS_FILL( sw_down_3d, 'AtPhDSR ',
                0471      &                           -1, Nr+1, 3, bi, bj, myThid )
                0472           CALL DIAGNOSTICS_FILL( lw_down_3d, 'AtPhDLR ',
                0473      &                           -1, Nr+1, 3, bi, bj, myThid )
                0474           CALL DIAGNOSTICS_FILL( sw_down_3d, 'AtPhInSR',
e259fb67de Jean*0475      &                           0, 1, 3, bi, bj, myThid )
a690051f17 Jean*0476           CALL DIAGNOSTICS_FILL( sw_net_3d , 'AtPhNTSR',
e259fb67de Jean*0477      &                           0, 1, 3, bi, bj, myThid )
a690051f17 Jean*0478           CALL DIAGNOSTICS_FILL( lw_net_3d , 'AtPhOLR ',
e259fb67de Jean*0479      &                           0, 1, 3, bi, bj, myThid )
                0480           CALL DIAGNOSTICS_FILL( sw_down_3d(:,:,Nr+1),'AtPhDSSR',
                0481      &                           0, 1, 3, bi, bj, myThid )
                0482           CALL DIAGNOSTICS_FILL( s_sw_dwn,  'AtPhNSSR',
                0483      &                           0, 1, 3, bi, bj, myThid )
                0484           CALL DIAGNOSTICS_FILL( s_lw_dwn,  'AtPhDSLR',
                0485      &                           0, 1, 3, bi, bj, myThid )
                0486           CALL DIAGNOSTICS_FILL( flux_r ,   'AtPhUSLR',
                0487      &                           0, 1, 3, bi, bj, myThid )
                0488           CALL DIAGNOSTICS_FILL( flux_t ,   'AtPhSens',
                0489      &                           0, 1, 3, bi, bj, myThid )
                0490           CALL DIAGNOSTICS_FILL( flux_q ,   'AtPhEvap',
                0491      &                           0, 1, 3, bi, bj, myThid )
                0492           CALL DIAGNOSTICS_FILL( flux_u , 'AtPhTauX',
                0493      &                           0, 1, 3, bi, bj, myThid )
                0494           CALL DIAGNOSTICS_FILL( flux_v , 'AtPhTauY',
                0495      &                           0, 1, 3, bi, bj, myThid )
a690051f17 Jean*0496           CALL DIAGNOSTICS_FILL( cond_dt_tg, 'AtPhdTlc',
                0497      &                           -1, Nr, 3, bi, bj, myThid )
e259fb67de Jean*0498           CALL DIAGNOSTICS_FILL( diff_dt_tg, 'AtPhdtTg',
                0499      &                           -1, Nr, 3, bi, bj, myThid )
                0500           CALL DIAGNOSTICS_FILL( diff_dt_qg, 'AtPhdtQg',
                0501      &                           -1, Nr, 3, bi, bj, myThid )
                0502           CALL DIAGNOSTICS_FILL( diss_heat,  'AtPhDisH',
                0503      &                           -1, Nr, 3, bi, bj, myThid )
                0504         ENDIF
b2ea1d2979 Jean*0505 #endif /* ALLOW_DIAGNOSTICS */
                0506 
                0507 C--   end bi,bj loops.
                0508        ENDDO
                0509       ENDDO
                0510 
                0511       CALL EXCH_UV_AGRID_3D_RL( atmPhys_dU, atmPhys_dV,
                0512      &                          .TRUE., Nr, myThid )
6c614dc1a3 Jean*0513 #ifdef COMPONENT_MODULE
                0514 C-    surface wind-stress at grid cell center (A-grid):
                0515       IF ( useCoupler ) THEN
                0516         CALL EXCH_UV_AGRID_3D_RL( taux, tauy,
                0517      &                           .TRUE., 1, myThid )
                0518        DO bj=myByLo(myThid),myByHi(myThid)
                0519         DO bi=myBxLo(myThid),myBxHi(myThid)
                0520          DO j=2-OLy,sNy+OLy
                0521           DO i=2-OLx,sNx+OLx
                0522             fu(i,j,bi,bj) = halfRL
                0523      &                    *( taux(i-1,j,bi,bj) + taux(i,j,bi,bj) )
                0524             fv(i,j,bi,bj) = halfRL
                0525      &                    *( tauy(i,j-1,bi,bj) + tauy(i,j,bi,bj) )
                0526           ENDDO
                0527          ENDDO
                0528         ENDDO
                0529        ENDDO
                0530        CALL ATM_STORE_MY_DATA( myTime, myIter, myThid )
                0531       ENDIF
                0532 #endif /* COMPONENT_MODULE */
b2ea1d2979 Jean*0533 
                0534       RETURN
                0535       END