! source file: /Users/jfidler/work/UVic_ESCM/2.9/updates/idealage_tx/source/mom/tracer.F subroutine tracer (joff, js, je, is, ie) !======================================================================= ! compute tracers at "tau+1" for rows js through je in the MW. ! input: ! joff = offset relating "j" in the MW to latitude "jrow" ! js = starting row in the MW ! je = ending row in the MW ! is = starting longitude index in the MW ! ie = ending longitude index in the MW !======================================================================= implicit none character(120) :: fname, new_file_name integer istrt, iend, i, k, j, ip, kr, jq, n, jp, jrow, iou, js integer je, limit, joff, is, ie, kmx, m, kb, idiag, index integer it(10), iu(10), ib(10), ic(10), nfnpzd, mfnpzd, mxfnpzd integer id_time, id_xt, id_yt, id_zt,fe_jlo,fe_m,fe_k,fe_n parameter (fe_n = 14) logical inqvardef, exists real rctheta, declin, gl, impo, expo, npp, time real remi, excr, graz, morp, morpt, morz, temp, swr, dayfrac real graz_Det, graz_Z, avej, avej_D, gmax, no3P, po4P, po4_D real phin, dz, prca, dprca, nud, bct, tap, fo2, so2, ai, hi, hs real npp_D, graz_D, morp_D, no3flag, deni, nfix, felimit real t_i, t_j, dz_t2r, dz_tr, dz_wtr, dx_t2r, dx_tr, dy_t2r real dy_tr, adv_tx, adv_ty, adv_tz, adv_txiso, adv_tyiso real adv_tziso, diff_tx, diff_ty, diff_tz, zmax, cont, drho real drhom1, wt, ahbi_cstr, ahbi_csu_dyur, gamma, rrstd, fy, fyz real fe_dy, fe_conc, fe_x(fe_n), fe_y(fe_n), bctz, felimit_D include "size.h" include "param.h" include "pconst.h" include "stdunits.h" include "accel.h" include "coord.h" include "cregin.h" include "csbc.h" include "emode.h" include "grdvar.h" include "hmixc.h" include "levind.h" include "mw.h" include "scalar.h" include "switch.h" include "timeavgs.h" include "tmngr.h" include "vmixc.h" include "diaga.h" include "ice.h" include "atm.h" include "npzd.h" parameter (istrt=2, iend=imt-1) real twodt(km) real snpzd(ntnpzd), tnpzd(ntnpzd) real src(imt,km,jsmw:jemw,nsrc) include "isopyc.h" include "fdift.h" !----------------------------------------------------------------------- ! bail out if starting row exceeds ending row !----------------------------------------------------------------------- if (js .gt. je) return !----------------------------------------------------------------------- ! limit the longitude indices based on those from the argument list ! Note: this is currently bypassed. istrt and iend are set as ! parameters to optimize performance !----------------------------------------------------------------------- ! istrt = max(2,is) ! iend = min(imt-1,ie) !----------------------------------------------------------------------- ! build coefficients to minimize advection and diffusion computation !----------------------------------------------------------------------- limit = min(je+1+joff,jmt) - joff do j=js,limit jrow = j + joff do i=istrt-1,iend cstdxtr(i,j) = cstr(jrow)*dxtr(i) cstdxt2r(i,j) = cstr(jrow)*dxtr(i)*p5 cstdxur(i,j) = cstr(jrow)*dxur(i) ah_cstdxur(i,j) = diff_cet*cstr(jrow)*dxur(i) enddo enddo !----------------------------------------------------------------------- ! calculation of biological interactions !----------------------------------------------------------------------- declin = sin((mod(relyr,1.) - 0.22)*2.*pi)*0.4 ! declination do k=1,km twodt(k) = c2dtts*dtxcel(k) nbio(k) = twodt(k)/dtnpzd dtbio(k) = twodt(k)/nbio(k) rdtts(k) = 1./twodt(k) rnbio(k) = 1./nbio(k) enddo tap = 2.*alpha*par do j=js,je jrow = j + joff do i=is,ie if (kmt(i,jrow) .gt. 0) then ai = aice(i,jrow,2) hi = hice(i,jrow,2) hs = hsno(i,jrow,2) ! calculate day fraction and incoming solar ! angle of incidence = lat - declin, refraction index = 1.33 rctheta = max(-1.5, min(1.5, tlat(i,jrow)/radian - declin)) rctheta = kw/sqrt(1. - (1. - cos(rctheta)**2.)/1.33**2.) dayfrac = min( 1., -tan(tlat(i,jrow)/radian)*tan(declin)) dayfrac = max(1e-12, acos(max(-1., dayfrac))/pi) swr = dnswr(i,jrow)*1e-3*(1. + ai*(exp(-ki*(hi + hs)) - 1.)) expo = 0.0 impo = 0.0 phin = 0.0 ! integrated phytoplankton prca = 0.0 ! integrated production of calcite kmx = min(kmt(i,jrow), kpzd) do k=1,kmx !----------------------------------------------------------------------- ! limit tracers to positive values !----------------------------------------------------------------------- tnpzd(1) = max(t(i,k,j,ipo4,taum1), trcmin) tnpzd(2) = max(t(i,k,j,iphyt,taum1), trcmin) tnpzd(3) = max(t(i,k,j,izoop,taum1), trcmin) tnpzd(4) = max(t(i,k,j,idetr,taum1), trcmin) tnpzd(5) = max(t(i,k,j,ino3,taum1), trcmin) tnpzd(6) = max(t(i,k,j,idiaz,taum1), trcmin) swr = swr*exp(-kc*phin) phin = phin + (tnpzd(6) + tnpzd(2))*dzt(k) gl = tap*swr*exp(ztt(k)*rctheta) impo = expo*dztr(k) bct = bbio**(cbio*t(i,k,j,itemp,taum1)) if (t(i,k,j,itemp,taum1).gt.20) then bctz = (0.5*(tanh(t(i,k,j,io2,taum1)*1000. - 8.)+1)) & *bbio**(cbio*20) else bctz = (0.5*(tanh(t(i,k,j,io2,taum1)*1000. - 8.)+1)) & *bct end if if (k.le.3) then ! create x (time) and y (data) arrays for fe interpolation ! first zero them out just in case fe_y(:) = 0 fe_x(:) = 0 ! write in values for days 0 and 365 for interp. bounds fe_y(1) = fe_dissolved(i,j,k,1) fe_x(1) = 0 fe_y(14) = fe_dissolved(i,j,k,12) fe_x(14) = 365 ! now create the rest of the array based on the BLING ! data which is from day 16 of each month do m=2,13 fe_y(m) = fe_dissolved(i,j,k,m-1) if (m.eq.2) then fe_x(2) = 16 else if (m.eq.3) then fe_x(3) = 47 else if (m.eq.4) then fe_x(4) = 75 else if (m.eq.5) then fe_x(5) = 106 else if (m.eq.6) then fe_x(6) = 136 else if (m.eq.7) then fe_x(7) = 167 else if (m.eq.8) then fe_x(8) = 197 else if (m.eq.9) then fe_x(9) = 228 else if (m.eq.10) then fe_x(10) = 259 else if (m.eq.11) then fe_x(11) = 289 else if (m.eq.12) then fe_x(12) = 320 else if (m.eq.13) then fe_x(13) = 350 endif enddo fe_jlo = 2 fe_m = 4 ! find fe data day index call hunt (fe_x,fe_n,dayoyr,fe_jlo) ! initialize the fe data array at the right day fe_k = min(max(fe_jlo-(fe_m-1)/2,1),fe_n+1-fe_m) ! interpolate the fe data, note this does not use the whole array call polint (fe_x(fe_k),fe_y(fe_k) &, fe_m,dayoyr,fe_conc,fe_dy) ! calculate the fe limitation term felimit = fe_conc/(kfe + fe_conc) felimit_D = fe_conc/(kfe_D + fe_conc) else felimit = 1 felimit_D = 1 end if ! print*,'felimit=',felimit ! decrease remineralisation rate in oxygen minimum zone nud = nud0*(0.65+0.35*tanh(t(i,k,j,io2,taum1)*1000.-6.)) !----------------------------------------------------------------------- ! call the npzd model !----------------------------------------------------------------------- call npzd_src (tnpzd, nbio(k), dtbio(k), gl, bct, impo &, dzt(k), dayfrac, wd(k), rkwz(k), nud &, snpzd, expo, graz, morp, morz, graz_Det &, graz_Z &, npp, morpt, remi, excr &, npp_D, graz_D, morp_D, nfix &, felimit, felimit_D &, bctz & ) ! These are source/sink terms snpzd(1:4) = snpzd(1:4)*rdtts(k) snpzd(5:6) = snpzd(5:6)*rdtts(k) expo = expo*rnbio(k) rexpo(i,k,j) = expo rgraz(i,k,j) = graz*rnbio(k) rgraz_Det(i,k,j) = graz_Det*rnbio(k) rgraz_Z(i,k,j) = graz_Z*rnbio(k) rmorp(i,k,j) = morp*rnbio(k) rmorz(i,k,j) = morz*rnbio(k) rnpp(i,k,j) = npp*rnbio(k) rmorpt(i,k,j) = morpt*rnbio(k) rremi(i,k,j) = remi*rnbio(k) rexcr(i,k,j) = excr*rnbio(k) rnpp_D(i,k,j) = npp_D*rnbio(k) rgraz_D(i,k,j) = graz_D*rnbio(k) rmorp_D(i,k,j) = morp_D*rnbio(k) rnfix(i,k,j) = nfix*rnbio(k) !----------------------------------------------------------------------- ! calculate detritus at the bottom and remineralize !----------------------------------------------------------------------- if (k .eq. kmt(i,jrow)) then rremi(i,k,j) = rremi(i,k,j) + expo snpzd(1) = snpzd(1) + redptn*expo snpzd(5) = snpzd(5) + expo endif !----------------------------------------------------------------------- ! set source/sink terms !----------------------------------------------------------------------- src(i,k,j,ispo4) = snpzd(1) src(i,k,j,isidealage) = 1 / 3.6e3 / 24 ! JG 15.1.14 src(i,k,j,istracer1) = 0.0 ! JG 12.6.14 src(i,k,j,istracer2) = 0.0 src(i,k,j,istracer3) = 0.0 ! if (k==1 .AND. jrow .lt. 14) then ! JG 22.1.14 ! src(i,k,j,istracer1) = 1.0 ! else ! src(i,k,j,istracer1) = 0.0 ! endif ! if (k==1 .AND. jrow > 14 .AND. jrow < 23) then ! src(i,k,j,istracer2) = 1.0 ! else ! src(i,k,j,istracer2) = 0.0 ! endif ! if (k==1 .AND. jrow > 22 .AND. jrow < 29) then ! src(i,k,j,istracer3) = 1.0 ! else ! src(i,k,j,istracer3) = 0.0 ! endif src(i,k,j,isphyt) = snpzd(2) src(i,k,j,iszoop) = snpzd(3) src(i,k,j,isdetr) = snpzd(4) src(i,k,j,isno3) = snpzd(5) src(i,k,j,isdiaz) = snpzd(6) ! production of calcite dprca = (morp+morz+(graz+graz_Z)*(1.-gamma1)) & *capr*redctn*rnbio(k) prca = prca + dprca*dzt(k) ! These are sources and sinks of DIC (i.e. remin - pp) ! all are based on po4 uptake and remineralization ! dprca is a correction term src(i,k,j,isdic) = (snpzd(1)*redctp - dprca) src(i,k,j,isalk) = (-snpzd(1)*redntp*1.e-3 - 2.*dprca) !----------------------------------------------------------------------- ! accumulate time averages !----------------------------------------------------------------------- if (timavgperts .and. .not. euler2) then ta_rnpp(i,k,jrow) = ta_rnpp(i,k,jrow) + rnpp(i,k,j) ta_rgraz(i,k,jrow) = ta_rgraz(i,k,jrow) + rgraz(i,k,j) ta_rgraz_Z(i,k,jrow) = ta_rgraz_Z(i,k,jrow) & + rgraz_Z(i,k,j) ta_rgraz_Det(i,k,jrow) = ta_rgraz_Det(i,k,jrow) & + rgraz_Det(i,k,j) ta_rmorp(i,k,jrow) = ta_rmorp(i,k,jrow) + rmorp(i,k,j) ta_rmorpt(i,k,jrow)= ta_rmorpt(i,k,jrow) + rmorpt(i,k,j) ta_rmorz(i,k,jrow) = ta_rmorz(i,k,jrow) + rmorz(i,k,j) ta_rexcr(i,k,jrow) = ta_rexcr(i,k,jrow) + rexcr(i,k,j) ta_rnpp_D(i,k,jrow) = ta_rnpp_D(i,k,jrow) & + rnpp_D(i,k,j) ta_rgraz_D(i,k,jrow) = ta_rgraz_D(i,k,jrow) & + rgraz_D(i,k,j) ta_rmorp_D(i,k,jrow) = ta_rmorp_D(i,k,jrow) & + rmorp_D(i,k,j) ta_rnfix(i,k,jrow) = ta_rnfix(i,k,jrow) + rnfix(i,k,j) endif ! calculate total export to get total import for next layer expo = expo*dzt(k) enddo kmx = kmt(i,jrow) do k=1,kmx ! limit oxygen consumption below concentrations of ! 5umol/kg as recommended in OCMIP fo2 = 0.5*tanh(t(i,k,j,io2,taum1)*1000. - 5.) ! sink of oxygen so2 = src(i,k,j,ispo4)*redotp src(i,k,j,iso2) = -so2*(0.5 + fo2) ! add denitrification as source term for NO3 no3flag = 0.5+sign(0.5,t(i,k,j,ino3,taum1)-trcmin) ! 800 = 0.8*1000 = (elec/mol O2)/(elec/mol NO3)*(mmol/mol) deni = 800.*no3flag*so2*(0.5 - fo2) src(i,k,j,isno3) = src(i,k,j,isno3) - deni rdeni(i,k,jrow) = deni enddo !----------------------------------------------------------------------- ! remineralize calcite !----------------------------------------------------------------------- kmx = kmt(i,jrow) do k=1,kmx-1 src(i,k,j,isdic) = src(i,k,j,isdic) + prca*rcak(k) src(i,k,j,isalk) = src(i,k,j,isalk) + 2.*prca*rcak(k) enddo src(i,kmx,j,isdic) = src(i,kmx,j,isdic) + prca*rcab(kmx) src(i,kmx,j,isalk) = src(i,kmx,j,isalk) + 2.*prca*rcab(kmx) !----------------------------------------------------------------------- ! accumulate time averages for full depth variables !----------------------------------------------------------------------- if (timavgperts .and. .not. euler2) then kmx = kmt(i,jrow) expo = prca ta_rprocal(i,jrow) = ta_rprocal(i,jrow) + prca do k=1,kmx expo = expo*dztr(k) ta_rremi(i,k,jrow) = ta_rremi(i,k,jrow) + rremi(i,k,j) ta_rexpo(i,k,jrow) = ta_rexpo(i,k,jrow) + rexpo(i,k,j) expo = expo - prca*rcak(k) ta_rexpocal(i,k,jrow) = ta_rexpocal(i,k,jrow) + expo expo = expo*dzt(k) ta_rdeni(i,k,jrow) = ta_rdeni(i,k,jrow) + rdeni(i,k,j) enddo endif endif enddo enddo !----------------------------------------------------------------------- ! set source for c14 !----------------------------------------------------------------------- do j=js,je jrow = j + joff do i=is,ie if (kmt(i,jrow) .gt. 0) then do k=1,kmt(i,jrow) src(i,k,j,isc14) = src(i,k,j,isdic)*rstd & - 3.836e-12*t(i,k,j,ic14,taum1) enddo endif enddo enddo !----------------------------------------------------------------------- ! solve for one tracer at a time ! n = 1 => temperature ! n = 2 => salinity ! n > 2 => other tracers (if applicable) !----------------------------------------------------------------------- do n=1,nt !----------------------------------------------------------------------- ! calculate advective tracer flux !----------------------------------------------------------------------- call adv_flux (joff, js, je, is, ie, n) !----------------------------------------------------------------------- ! calculate diffusive flux across eastern and northern faces ! of "T" cells due to various parameterizations for diffusion. !----------------------------------------------------------------------- ! diffusive flux on eastern face of "T" cells do j=js,je do k=1,km do i=istrt-1,iend diff_fe(i,k,j) = & ah_cstdxur(i,j)* & (t(i+1,k,j,n,taum1) - t(i,k,j,n,taum1)) enddo enddo enddo ! diffusive flux on northern face of "T" cells ! (background for isopycnal mixing) do j=js-1,je jrow = j + joff do k=1,km do i=istrt,iend diff_fn(i,k,j) = & diff_cnt* & csu_dyur(jrow)*(t(i,k,j+1,n,taum1) - t(i,k,j,n,taum1)) enddo enddo enddo !----------------------------------------------------------------------- ! calculate diffusive flux across bottom face of "T" cells !----------------------------------------------------------------------- do j=js,je do k=1,km-1 do i=istrt,iend diff_fb(i,k,j) = diff_cbt(i,k,j)*dzwr(k)* & (t(i,k,j,n,taum1) - t(i,k+1,j,n,taum1)) enddo enddo enddo !----------------------------------------------------------------------- ! compute isopycnal diffusive flux through east, north, ! and bottom faces of T cells. !----------------------------------------------------------------------- call isoflux (joff, js, je, is, ie, n) !----------------------------------------------------------------------- ! set surface and bottom vert b.c. on "T" cells for diffusion ! and advection. for isopycnal diffusion, set adiabatic boundary ! conditions. ! note: the b.c. at adv_fb(i,k=bottom,j) is set by the above code. ! However, it is not set when k=km so it is set below. ! adv_fb(i,km,j) is always zero (to within roundoff). !----------------------------------------------------------------------- do j=js,je jrow = j + joff do i=istrt,iend kb = kmt(i,jrow) diff_fb(i,0,j) = stf(i,j,n) diff_fb(i,kb,j) = btf(i,j,n) adv_fb(i,0,j) = adv_vbt(i,0,j)*(t(i,1,j,n,tau) + & t(i,1,j,n,tau)) adv_fb(i,km,j) = adv_vbt(i,km,j)*t(i,km,j,n,tau) enddo enddo !----------------------------------------------------------------------- ! set source term for "T" cells !----------------------------------------------------------------------- source(:,:,:) = c0 if (itrc(n) .ne. 0) then do j=js,je do k=1,km do i=istrt,iend source(i,k,j) = src(i,k,j,itrc(n)) enddo enddo enddo endif !----------------------------------------------------------------------- ! solve for "tau+1" tracer using statement functions to represent ! each component of the calculation !----------------------------------------------------------------------- ! 1st: solve using all components which are treated explicitly do j=js,je jrow = j + joff do k=1,km twodt(k) = c2dtts*dtxcel(k) do i=istrt,iend t(i,k,j,n,taup1) = t(i,k,j,n,taum1) + twodt(k)*( & DIFF_Tx(i,k,j) + DIFF_Ty(i,k,j,jrow,n) + DIFF_Tz(i,k,j) & - ADV_Tx(i,k,j) - ADV_Ty(i,k,j,jrow,n) - ADV_Tz(i,k,j) & + source(i,k,j) & )*tmask(i,k,j) enddo enddo enddo ! 2nd: add in portion of vertical diffusion handled implicitly call ivdift (joff, js, je, istrt, iend, n, twodt) do j=js,je call setbcx (t(1,1,j,n,taup1), imt, km) enddo !----------------------------------------------------------------------- ! construct diagnostics associated with tracer "n" !----------------------------------------------------------------------- call diagt1 (joff, js, je, istrt, iend, n, twodt) !----------------------------------------------------------------------- ! end of tracer component "n" loop !----------------------------------------------------------------------- enddo !----------------------------------------------------------------------- ! explicit convection: adjust column if gravitationally unstable !----------------------------------------------------------------------- call convct2 (t(1,1,1,1,taup1), joff, js, je, is, ie, kmt) do j=js,je do n=1,nt call setbcx (t(1,1,j,n,taup1), imt, km) enddo enddo if (timavgperts .and. eots) then if (joff .eq. 0) nta_conv = nta_conv + 1 do j=js,je jrow = j + joff do i=istrt,iend ta_totalk(i,jrow) = ta_totalk(i,jrow) + totalk(i,j) ta_vdepth(i,jrow) = ta_vdepth(i,jrow) + vdepth(i,j) ta_pe(i,jrow) = ta_pe(i,jrow) + pe(i,j) enddo enddo endif !----------------------------------------------------------------------- ! construct diagnostics after convection !----------------------------------------------------------------------- idiag = 10 call diagt2 (joff, js, je, istrt, iend, idiag) !----------------------------------------------------------------------- ! filter tracers at high latitudes !----------------------------------------------------------------------- if (istrt .eq. 2 .and. iend .eq. imt-1) then call filt (joff, js, je) else write (stdout,'(a)') & 'Error: filtering requires is=2 and ie=imt-1 in tracer' stop '=>tracer' endif do n=1,nt do j=js,je call setbcx (t(1,1,j,n,taup1), imt, km) enddo enddo ! JG 15.1.14 do j=js,je jrow = j + joff do i=is,ie t(i,1,j,iidealage,taup1)=0.0 if (jrow .lt. 14) then ! JG 12.1.14 t(i,1,j,itracer1,taup1) = 1.0 else t(i,1,j,itracer1,taup1) = 0.0 endif if (jrow > 14 .AND. jrow < 23) then t(i,1,j,itracer2,taup1) = 1.0 else t(i,1,j,itracer2,taup1) = 0.0 endif if (jrow > 22 .AND. jrow < 29) then t(i,1,j,itracer3,taup1) = 1.0 else t(i,1,j,itracer3,taup1) = 0.0 endif enddo enddo !----------------------------------------------------------------------- ! construct diagnostics after filtering (for total dT/dt) !----------------------------------------------------------------------- idiag = 1 call diagt2 (joff, js, je, istrt, iend, idiag) !----------------------------------------------------------------------- ! if needed, construct the Atmos S.B.C.(surface boundary conditions) ! averaged over this segment ! eg: SST and possibly SSS !----------------------------------------------------------------------- call asbct (joff, js, je, istrt, iend, isst, itemp) call asbct (joff, js, je, istrt, iend, isss, isalt) call asbct (joff, js, je, istrt, iend, issdic, idic) call asbct (joff, js, je, istrt, iend, issc14, ic14) call asbct (joff, js, je, istrt, iend, issalk, ialk) call asbct (joff, js, je, istrt, iend, isso2, io2) call asbct (joff, js, je, istrt, iend, isspo4, ipo4) call asbct (joff, js, je, istrt, iend, issno3, ino3) !----------------------------------------------------------------------- ! calculate diagnostic delta carbon 14 !----------------------------------------------------------------------- if (tsiperts .or. timavgperts) then rrstd = 1000./rstd do j=js,je jrow = j + joff do k=1,km do i=istrt,iend dc14(i,k,j) = (rrstd*t(i,k,j,ic14,taup1) & /(t(i,k,j,idic,taup1) + epsln) - 1000.) & *tmask(i,k,j) enddo enddo enddo endif if (tsiperts .and. eots) then if (js+joff .eq. 2) dc14bar = 0. do j=js,je jrow = j + joff fy = cst(jrow)*dyt(jrow) do k=1,km fyz = fy*dzt(k) do i=istrt,iend dc14bar = dc14bar + dc14(i,k,j)*dxt(i)*fyz*tmask(i,k,j) enddo enddo enddo endif if (timavgperts .and. .not. euler2) then do j=js,je jrow = j + joff do k=1,km do i=istrt,iend ta_dc14(i,k,jrow) = ta_dc14(i,k,jrow) + dc14(i,k,j) enddo enddo enddo endif return end subroutine diagt1 (joff, js, je, is, ie, n, twodt) !----------------------------------------------------------------------- ! construct diagnostics associated with tracer component "n" ! input: ! joff = offset relating "j" in the MW to latitude "jrow" ! js = starting row in the MW ! je = ending row in the MW ! is = starting longitude index in the MW ! ie = ending longitude index in the MW ! n = (1,2) = (u,v) velocity component ! twodt = (2*dtts,dtts) on (leapfrog,mixing) time steps !----------------------------------------------------------------------- implicit none integer i, k, j, ip, kr, jq, n, jp, jrow, js, je, joff, is, ie integer mask, m real t_i, t_j, dz_t2r, dz_tr, dz_wtr, dx_t2r, dx_tr, dy_t2r real dy_tr, adv_tx, adv_ty, adv_tz, adv_txiso, adv_tyiso real adv_tziso, diff_tx, diff_ty, diff_tz, dtdx, dtdy, dtdz real r2dt, cosdyt, fx, darea, boxar, rtwodt, sumdx, delx real sumdxr, dxdy, dxdydz include "size.h" include "param.h" include "pconst.h" include "stdunits.h" include "accel.h" include "coord.h" include "cregin.h" include "csbc.h" include "ctavg.h" include "diag.h" include "diaga.h" include "emode.h" include "grdvar.h" include "hmixc.h" include "levind.h" include "mw.h" include "scalar.h" include "switch.h" include "vmixc.h" real temp1(imt,km), temp2(imt,km), temp3(imt,km) real twodt(km) include "isopyc.h" include "fdift.h" !----------------------------------------------------------------------- ! diagnostic: integrate |d(tracer)/dt| and tracer variance on "tau" ! globally !----------------------------------------------------------------------- if (tsiperts .and. eots) then do j=js,je jrow = j + joff r2dt = c1/c2dtts cosdyt = cst(jrow)*dyt(jrow) do k=1,km fx = r2dt/dtxcel(k) do i=is,ie darea = dzt(k)*dxt(i)*cosdyt*tmask(i,k,j) temp3(i,k) = t(i,k,j,n,tau)*darea temp1(i,k) = t(i,k,j,n,tau)**2*darea temp2(i,k) = abs(t(i,k,j,n,taup1)-t(i,k,j,n,taum1))* & darea*fx enddo do i=is,ie tbar(k,n,jrow) = tbar(k,n,jrow) + temp3(i,k) travar(k,n,jrow) = travar(k,n,jrow) + temp1(i,k) dtabs(k,n,jrow) = dtabs(k,n,jrow) + temp2(i,k) enddo enddo enddo endif !----------------------------------------------------------------------- ! diagnostic: accumulate tracers for averages under horizontal ! regions (use units of meters, rather than cm) !----------------------------------------------------------------------- if (tavgts .and. eots) then do j=js,je jrow = j + joff do i=is,ie mask = mskhr(i,jrow) if (mask .ne. 0) then boxar = cst(jrow)*dxt(i)*dyt(jrow)*tmask(i,1,j)*0.0001 sumbf(mask,n) = sumbf(mask,n) + stf(i,j,n)*boxar do k=1,km sumbk(mask,k,n) = sumbk(mask,k,n) + t(i,k,j,n,tau) & *boxar*dzt(k)*tmask(i,k,j)*0.01 enddo endif enddo enddo endif !----------------------------------------------------------------------- ! diagnostic: compute the northward transport components of ! each tracer !----------------------------------------------------------------------- if (gyrets .and. eots) call gyre (joff, js, je, is, ie, n) !----------------------------------------------------------------------- ! diagnostic: integrate r.h.s. terms in the tracer equations ! over specified regional volumes. !----------------------------------------------------------------------- if (trmbts .and. eots) call ttb1 (joff, js, je, is, ie, n) return end subroutine diagt2 (joff, js, je, is, ie, idiag) !----------------------------------------------------------------------- ! construct d(tracer)/dt diagnostics ! input: ! joff = offset relating "j" in the MW to latitude "jrow" ! js = starting row in the MW ! je = ending row in the MW ! is = starting longitude index in the MW ! ie = ending longitude index in the MW ! idiag = 1 => total tracer change ! idiag = 10 => change of tracer due to filtering(also convection) !----------------------------------------------------------------------- implicit none integer idiag, j, js, je, k, i, joff, iocv, jrow, is, ie real rdt, reltim, period include "size.h" include "param.h" include "pconst.h" include "stdunits.h" include "coord.h" include "diaga.h" include "iounit.h" include "mw.h" include "scalar.h" include "switch.h" include "tmngr.h" include "timeavgs.h" !----------------------------------------------------------------------- ! diagnostic: integrate d/dt(tracer) over specified regional volumes ! after convection and filtering !----------------------------------------------------------------------- if (trmbts .and. eots) call ttb2 (joff, js, je, is, ie, idiag) return end subroutine asbct (joff, js, je, is, ie, isbc, itr) !----------------------------------------------------------------------- ! construct the Atmos S.B.C. (surface boundary conditions) ! input: ! joff = offset relating "j" in the MW to latitude "jrow" ! js = starting row in the MW ! je = ending row in the MW ! is = starting longitude index in the MW ! ie = ending longitude index in the MW ! isbc = index for sbc ! itr = index for tracer !----------------------------------------------------------------------- implicit none integer isbc, itr, j, js, je, jrow, joff, i, is, ie real rts include "size.h" include "param.h" include "pconst.h" include "stdunits.h" include "csbc.h" include "levind.h" include "mw.h" include "scalar.h" include "switch.h" ! initialize the Atmos S.B.C. at the start of each ocean segment ! (do not alter values in land) if (isbc .le. 0 .or. itr .le. 0) return if (eots .and. osegs) then do j=js,je jrow = j + joff do i=is,ie if (kmt(i,jrow) .ne. 0) sbc(i,jrow,isbc) = c0 enddo enddo endif ! accumulate surface tracers for the Atmos S.B.C. every time step if (eots) then do j=js,je jrow = j + joff do i=is,ie sbc(i,jrow,isbc) = sbc(i,jrow,isbc)+t(i,1,j,itr,taup1) enddo enddo endif ! average the surface tracers for the Atmos S.B.C. at the end of ! each ocean segment. (do not alter values in land) if (eots .and. osege) then rts = c1/ntspos do j=js,je jrow = j + joff do i=is,ie if (kmt(i,jrow) .ne. 0) & sbc(i,jrow,isbc) = rts*sbc(i,jrow,isbc) enddo enddo endif return end subroutine ivdift (joff, js, je, is, ie, n, twodt) !----------------------------------------------------------------------- ! solve vertical diffusion of tracers implicitly ! input: ! joff = offset relating "j" in the MW to latitude "jrow" ! js = starting row in the MW ! je = ending row in the MW ! is = starting longitude index in the MW ! ie = ending longitude index in the MW ! n = tracer component ! twodt = (2*dtts, dtts) on (leapfrog, mixing) time steps !----------------------------------------------------------------------- implicit none integer j, js, je, k, i, is, ie, n, joff real rc2dt include "size.h" include "param.h" include "pconst.h" include "stdunits.h" include "levind.h" include "mw.h" include "switch.h" include "vmixc.h" real twodt(km) ! store terms to compute implicit vertical mixing on ! diagnostic time steps if (trmbts .and. eots) then do j=js,je do k=1,km do i=is,ie zzi(i,k,j) = t(i,k,j,n,taup1) enddo enddo enddo endif call invtri (t(1,1,1,n,taup1), stf(1,1,n), btf(1,1,n) &, diff_cbt(1,1,jsmw), twodt, kmt, tmask(1,1,1), is, ie &, joff, js, je) ! compute residual implicit vertical mixing if (trmbts .and. eots) then do j=js,je do k=1,km rc2dt = c1/twodt(k) do i=is,ie zzi(i,k,j) = rc2dt*(t(i,k,j,n,taup1) - zzi(i,k,j)) enddo enddo enddo endif return end subroutine swflux0 (joff, js, je, is, ie, source) return end ! subroutine polint (xa,ya,n,x,y,dy) ! implicit none ! integer n,i,m,ns, nmax parameter (nmax = 10) ! largest anticipated value of n real dy,x,y,xa(n),ya(n) real den,dif,dift,ho,hp,w,c(nmax),d(nmax) ! ! Given arrays xa and ya, each of length n, and a given value x, this routine ! returns a value y, and an error estimate dy. If P(x) is the polynomial of ! degree N-1 such that P(xai) = yai, i=1,....,n then the returned value y=P(x) ! ! ns=1 dif=abs(x-xa(1)) ! do i=1,n dift = abs(x-xa(i)) IF (dift.lt.dif) THEN ns = i dif = dift END IF ! c(i) = ya(i) d(i) = ya(i) ! enddo ! y = ya (ns) ! this is the initial approximation to y ns = ns-1 ! do m=1,n-1 do i=1,n-m ho = xa(i) - x hp = xa (i+m) -x w = c(i+1) - d(i) den = ho - hp IF (den.eq.0.) PAUSE 'failure in polint' den = w/den d(i) = hp*den c(i) = ho*den enddo IF (2*ns.lt.n-m) THEN dy = c(ns+1) ELSE dy = d(ns) ns = ns - 1 END IF y = y + dy enddo return end ! subroutine hunt (xa,n,x,jlo) ! integer jlo,n,inc,jhi,jm real x, xa(n) logical ascnd ! ! Given an array xa(1:n), and given a value x, returns a value jlo ! such that x is between xa(jlo) and xa(jlo +1). xa(1:n) must be ! monotonic, either increasing or decreasing. jlo=0 or jlo=n is ! returned to indicate that x is out of range. jlo on input is taken ! as the initial guess for jlo on output ! ascnd = xa(n).ge.xa(1) ! True if ascending order of table, false otherwise IF (jlo.le.0.or.jlo.gt.n) THEN ! Input guess not useful. Go to bisection jlo = 0 jhi = n+1 goto 3 END IF inc = 1 ! Set the hunting increment IF (x.ge.xa(jlo).eqv.ascnd) THEN ! Hunting up: 1 jhi = jlo + inc IF (jhi.gt.n) THEN ! Done hunting, since off end of table jhi = n + 1 ELSE IF (x.ge.xa(jhi).eqv.ascnd) THEN ! Not done hunting jlo = jhi inc = inc + inc ! so double the increment goto 1 ! and try again END IF ! Done hunting, value bracketed. ELSE ! Hunt down: jhi = jlo 2 jlo = jhi-inc IF (jlo.lt.1) THEN !Done hunting, since off end of table jlo = 0 ELSE IF (x.lt.xa(jlo).eqv.ascnd) THEN ! Not done hunting jhi = jlo inc = inc + inc ! so double the increment goto 2 ! and try again END IF ! Done hunting, value bracketed END IF ! Hunt is done, so begin the final bisection phase: ! 3 IF (jhi-jlo.eq.1) THEN IF (x.eq.xa(n)) jlo = n-1 IF (x.eq.xa(1)) jlo = 1 RETURN END IF jm = (jhi + jlo)/2 IF (x.ge.xa(jm).eqv.ascnd) THEN jlo = jm ELSE jhi = jm END IF ! GOTO 3 ! end