New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
air_sea.F90 in branches/UKMO/dev_r5518_GO6_package_asm_surf_bgc_v2/NEMOGCM/NEMO/TOP_SRC/MEDUSA – NEMO

source: branches/UKMO/dev_r5518_GO6_package_asm_surf_bgc_v2/NEMOGCM/NEMO/TOP_SRC/MEDUSA/air_sea.F90 @ 8495

Last change on this file since 8495 was 8495, checked in by dford, 7 years ago

Merge in changes from dev_r5518_GO6_package_asm_surf_bgc, and adapt to the updated MEDUSA structure.

File size: 27.1 KB
Line 
1MODULE air_sea_mod
2   !!======================================================================
3   !!                         ***  MODULE air_sea_mod  ***
4   !! Calculate the carbon chemistry for the whole ocean
5   !!======================================================================
6   !! History :
7   !!   -   ! 2017-04 (M. Stringer)        Code taken from trcbio_medusa.F90
8   !!----------------------------------------------------------------------
9#if defined key_medusa
10   !!----------------------------------------------------------------------
11   !!                                                   MEDUSA bio-model
12   !!----------------------------------------------------------------------
13
14   IMPLICIT NONE
15   PRIVATE
16     
17   PUBLIC   air_sea        ! Called in trcbio_medusa.F90
18
19   !!----------------------------------------------------------------------
20   !! NEMO/TOP 2.0 , LOCEAN-IPSL (2007)
21   !! $Id$
22   !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt)
23   !!----------------------------------------------------------------------
24
25CONTAINS
26
27   SUBROUTINE air_sea( kt )
28      !!---------------------------------------------------------------------
29      !!                     ***  ROUTINE air_sea  ***
30      !! This called from TRC_BIO_MEDUSA and
31      !!  - calculate air-sea gas exchange
32      !!  - river inputs
33      !!----------------------------------------------------------------------
34      USE bio_medusa_mod,    ONLY: f_riv_alk, f_riv_c, f_riv_n,           &
35                                   f_riv_si, f_runoff,                    & 
36                                   fgco2, zphn, zphd,                     &
37# if defined key_roam
38                                   dms_andr, dms_andr2d, dms_aran,        &
39                                   dms_aran2d, dms_hall, dms_hall2d,      &
40                                   dms_simo, dms_simo2d, dms_surf,        &
41                                   dms_surf2d, dms_andm, dms_andm2d,      &
42                                   dms_nlim, dms_wtkn,                    &
43                                   f_co2flux, f_co2flux2d,                &
44                                   f_co2starair_2d, f_co3,                &
45                                   f_dcf, f_fco2a_2d, f_fco2w_2d,         &
46                                   f_h2co3, f_hco3, f_henry,              &
47                                   f_kw660, f_kw6602d,                    &
48                                   f_o2flux, f_o2flux2d, f_o2sat,         &
49                                   f_o2sat2d, f_ocndpco2_2d,              &
50                                   f_ocnk0_2d, f_ocnkwco2_2d,             &
51                                   f_ocnrhosw_2d, f_ocnschco2_2d,         &
52                                   f_omarg, f_omcal,                      &
53                                   f_pco2a2d, f_pco2atm, f_pco2w,         &
54                                   f_pco2w2d, f_ph, f_pp0, f_pp02d,       &
55                                   f_TALK, f_TALK2d, f_TDIC, f_TDIC2d,    &
56                                   f_xco2a, f_xco2a_2d,                   &
57                                   zalk, zdic, zoxy, zsal, ztmp,          &
58# endif
59# if defined key_mocsy
60                                   zpho,                                  &
61# endif
62                                   zchd, zchn, zdin, zsil
63      USE dom_oce,           ONLY: e3t_0, e3t_n, gphit, tmask
64# if defined key_iomput
65      USE iom,               ONLY: lk_iomput
66# endif
67      USE in_out_manager,    ONLY: lwp, numout
68      USE oce,               ONLY: PCO2a_in_cpl
69      USE par_kind,          ONLY: wp
70      USE par_oce,           ONLY: jpi, jpim1, jpj, jpjm1
71      USE sbc_oce,           ONLY: fr_i, lk_oasis, qsr, wndm
72      USE sms_medusa,        ONLY: jdms, jdms_input, jdms_model,          &
73                                   jriver_alk, jriver_c,                  &
74                                   jriver_n, jriver_si,                   &
75                                   riv_alk, riv_c, riv_n, riv_si,         &
76                                   zn_dms_chd, zn_dms_chn, zn_dms_din,    &
77                                   zn_dms_mld, zn_dms_qsr,                &
78# if defined key_foam_medusa
79                                   f2_pco2w, f2_fco2w,                    &
80# endif
81                                   xnln, xnld 
82      USE trc,               ONLY: med_diag
83      USE zdfmxl,            ONLY: hmld
84
85# if defined key_roam
86      USE gastransfer,       ONLY: gas_transfer
87#  if defined key_mocsy
88      USE mocsy_wrapper,     ONLY: mocsy_interface
89#  else
90      USE trcco2_medusa,     ONLY: trc_co2_medusa
91#   if defined key_foam_medusa
92      USE mocsy_mainmod,     ONLY: p2fCO2
93#   endif
94#  endif
95      USE trcdms_medusa,     ONLY: trc_dms_medusa
96      USE trcoxy_medusa,     ONLY: trc_oxy_medusa
97# endif
98
99   !!* Substitution
100#  include "domzgr_substitute.h90"
101
102      !! time (integer timestep)
103      INTEGER, INTENT( in ) :: kt
104
105      !! Loop variables
106      INTEGER :: ji, jj
107
108# if defined key_roam
109      !! jpalm 14-07-2016: convert CO2flux diag from mmol/m2/d to kg/m2/s
110      REAL, PARAMETER :: weight_CO2_mol = 44.0095  !! g / mol
111      REAL, PARAMETER :: secs_in_day    = 86400.0  !! s / d
112      REAL, PARAMETER :: CO2flux_conv   = (1.e-6 * weight_CO2_mol) / secs_in_day
113
114      INTEGER :: iters
115
116      !! AXY (23/06/15): additional diagnostics for MOCSY and oxygen
117      REAL(wp), DIMENSION(jpi,jpj) :: f_fco2w, f_rhosw
118      REAL(wp), DIMENSION(jpi,jpj) :: f_fco2atm
119      REAL(wp), DIMENSION(jpi,jpj) :: f_schmidtco2, f_kwco2, f_K0
120      REAL(wp), DIMENSION(jpi,jpj) :: f_co2starair, f_dpco2
121      !! Output arguments from mocsy_interface, which aren't used
122      REAL(wp) :: f_BetaD_dum, f_opres_dum
123      REAL(wp) :: f_insitut_dum
124      REAL(wp) :: f_kwo2_dum
125# endif
126
127
128# if defined key_roam
129      !!-----------------------------------------------------------
130      !! Air-sea gas exchange
131      !!-----------------------------------------------------------
132
133#   if defined key_debug_medusa
134               IF (lwp) write (numout,*)                     & 
135               'air-sea: gas_transfer kt = ', kt
136               CALL flush(numout)
137#   endif
138      DO jj = 2,jpjm1
139         DO ji = 2,jpim1
140            !! OPEN wet point IF..THEN loop
141            if (tmask(ji,jj,1) == 1) then
142               IF (lk_oasis) THEN
143                  !! use 2D atm xCO2 from atm coupling
144                  f_xco2a(ji,jj) = PCO2a_in_cpl(ji,jj)
145               ENDIF
146               !!
147               !! AXY (23/06/15): as part of an effort to update the
148               !!                 carbonate chemistry in MEDUSA, the gas
149               !!                 transfer velocity used in the carbon
150               !!                 and oxygen cycles has been harmonised
151               !!                 and is calculated by the same function
152               !!                 here; this harmonisation includes
153               !!                 changes to the PML carbonate chemistry
154               !!                 scheme so that it too makes use of the
155               !!                 same gas transfer velocity; the
156               !!                 preferred parameterisation of this is
157               !!                 Wanninkhof (2014), option 7
158               !!
159               CALL gas_transfer( wndm(ji,jj), 1, 7,         &  ! inputs
160                                  f_kw660(ji,jj) )              ! outputs
161            ENDIF
162         ENDDO
163      ENDDO
164
165#   if defined key_debug_medusa
166               IF (lwp) write (numout,*)                     &
167               'air-sea: carb-chem kt = ', kt
168               CALL flush(numout)
169#   endif
170      DO jj = 2,jpjm1
171         DO ji = 2,jpim1
172            if (tmask(ji,jj,1) == 1) then
173               !! air pressure (atm); ultimately this will use air
174               !! pressure at the base of the UKESM1 atmosphere
175               !!                                     
176               f_pp0(ji,jj)   = 1.0
177               !!
178               !! IF(lwp) WRITE(numout,*) ' MEDUSA ztmp    =', ztmp(ji,jj)
179               !! IF(lwp) WRITE(numout,*) ' MEDUSA wndm    =', wndm(ji,jj)
180               !! IF(lwp) WRITE(numout,*) ' MEDUSA fr_i    =', fr_i(ji,jj)
181               !!
182#  if defined key_axy_carbchem
183#   if defined key_mocsy
184               !!
185               !! AXY (22/06/15): use Orr & Epitalon (2015) MOCSY-2 carbonate
186               !!                 chemistry package; note that depth is set to
187               !!                 zero in this call
188               CALL mocsy_interface(ztmp(ji,jj),zsal(ji,jj),zalk(ji,jj),     &
189                                    zdic(ji,jj),zsil(ji,jj),zpho(ji,jj),     &
190                                    f_pp0(ji,jj),0.0,                        &
191                                    gphit(ji,jj),f_kw660(ji,jj),             &
192                                    f_xco2a(ji,jj),1,f_ph(ji,jj),            &
193                                    f_pco2w(ji,jj),f_fco2w(ji,jj),           &
194                                    f_h2co3(ji,jj),f_hco3(ji,jj),            &
195                                    f_co3(ji,jj),f_omarg(ji,jj),             &
196                                    f_omcal(ji,jj),f_BetaD_dum,              &
197                                    f_rhosw(ji,jj),f_opres_dum,              &
198                                    f_insitut_dum,f_pco2atm(ji,jj),          &
199                                    f_fco2atm(ji,jj),f_schmidtco2(ji,jj),    &
200                                    f_kwco2(ji,jj),f_K0(ji,jj),              &
201                                    f_co2starair(ji,jj),f_co2flux(ji,jj),    &
202                                    f_dpco2(ji,jj))
203               !! mmol / m3 -> umol / kg
204               f_TDIC(ji,jj) = (zdic(ji,jj) / f_rhosw(ji,jj)) * 1000.
205               !! meq / m3 ->  ueq / kg
206               f_TALK(ji,jj) = (zalk(ji,jj) / f_rhosw(ji,jj)) * 1000.
207               f_dcf(ji,jj)  = f_rhosw(ji,jj)
208            ENDIF
209         ENDDO
210      ENDDO
211
212#   else   
213
214      DO jj = 2,jpjm1
215         DO ji = 2,jpim1
216            if (tmask(ji,jj,1) == 1) then     
217               iters = 0
218               !!
219               !! carbon dioxide (CO2); Jerry Blackford code (ostensibly
220               !! OCMIP-2, but not)
221               CALL trc_co2_medusa(ztmp(ji,jj),zsal(ji,jj),zdic(ji,jj),      &
222                                   zalk(ji,jj),0.0,                          &
223                                   f_kw660(ji,jj),f_xco2a(ji,jj),            &
224                                   f_ph(ji,jj),                              &
225                                   f_pco2w(ji,jj),f_h2co3(ji,jj),            &
226                                   f_hco3(ji,jj),f_co3(ji,jj),               &
227                                   f_omcal(ji,jj),f_omarg(ji,jj),            &
228                                   f_co2flux(ji,jj),f_TDIC(ji,jj),           &
229                                   f_TALK(ji,jj),f_dcf(ji,jj),               &
230                                   f_henry(ji,jj),iters)
231               !!
232               !! AXY (09/01/14): removed iteration and NaN checks; these have
233               !!                 been moved to trc_co2_medusa together with a
234               !!                 fudge that amends erroneous values (this is
235               !!                 intended to be a temporary fudge!); the
236               !!                 output warnings are retained here so that
237               !!                 failure position can be determined
238               if (iters .eq. 25) then
239                  IF(lwp) WRITE(numout,*) 'air-sea: ITERS WARNING, ',       &
240                     iters, ' AT (', ji, ', ', jj, ', 1) AT ', kt
241               endif
242#    if defined key_foam_medusa
243               !! DAF (Aug 2017): calculate fCO2 for observation operator
244               CALL p2fCO2( f_pco2w, ztmp, f_pp0, 0.0, 1, f_fco2w )
245#    endif
246            ENDIF
247         ENDDO
248      ENDDO
249
250#   endif
251#  else
252
253      DO jj = 2,jpjm1
254         DO ji = 2,jpim1
255            if (tmask(ji,jj,1) == 1) then
256               !! AXY (18/04/13): switch off carbonate chemistry
257               !!                 calculations; provide quasi-sensible
258               !!                 alternatives
259               f_ph(ji,jj)           = 8.1
260               f_pco2w(ji,jj)        = f_xco2a(ji,jj)
261               f_h2co3(ji,jj)        = 0.005 * zdic(ji,jj)
262               f_hco3(ji,jj)         = 0.865 * zdic(ji,jj)
263               f_co3(ji,jj)          = 0.130 * zdic(ji,jj)
264               f_omcal(ji,jj) = 4.
265               f_omarg(ji,jj) = 2.
266               f_co2flux(ji,jj)      = 0.
267               f_TDIC(ji,jj)         = zdic(ji,jj)
268               f_TALK(ji,jj)         = zalk(ji,jj)
269               f_dcf(ji,jj)          = 1.026
270               f_henry(ji,jj)        = 1.
271               !! AXY (23/06/15): add in some extra MOCSY diagnostics
272               f_fco2w(ji,jj)        = f_xco2a(ji,jj)
273! This doesn't seem to be used - marc 16/5/17
274!               f_BetaD(ji,jj)        = 1.
275               f_rhosw(ji,jj)        = 1.026
276! This doesn't seem to be used - marc 16/5/17
277!               f_opres(ji,jj)        = 0.
278!               f_insitut(ji,jj)      = ztmp(ji,jj)
279               f_pco2atm(ji,jj)      = f_xco2a(ji,jj)
280               f_fco2atm(ji,jj)      = f_xco2a(ji,jj)
281               f_schmidtco2(ji,jj)   = 660.
282               f_kwco2(ji,jj)        = 0.
283               f_K0(ji,jj)           = 0.
284               f_co2starair(ji,jj)   = f_xco2a(ji,jj)
285               f_dpco2(ji,jj)        = 0.
286            ENDIF
287         ENDDO
288      ENDDO
289#  endif
290
291      DO jj = 2,jpjm1
292         DO ji = 2,jpim1
293            if (tmask(ji,jj,1) == 1) then
294               !! mmol/m2/s -> mmol/m3/d; correct for sea-ice; divide
295               !! through by layer thickness
296               f_co2flux(ji,jj) = (1. - fr_i(ji,jj)) * f_co2flux(ji,jj) *    &
297                                  86400. / fse3t(ji,jj,1)
298               !!
299               !! oxygen (O2); OCMIP-2 code
300               !! AXY (23/06/15): amend input list for oxygen to account
301               !!                 for common gas transfer velocity
302               CALL trc_oxy_medusa(ztmp(ji,jj),zsal(ji,jj),f_kw660(ji,jj),   &
303                                   f_pp0(ji,jj),zoxy(ji,jj),                 &
304                                   f_kwo2_dum,f_o2flux(ji,jj),               &
305                                   f_o2sat(ji,jj))
306               !!
307               !! mmol/m2/s -> mol/m3/d; correct for sea-ice; divide
308               !! through by layer thickness
309               f_o2flux(ji,jj)  = (1. - fr_i(ji,jj)) * f_o2flux(ji,jj) *     &
310                                  86400. / fse3t(ji,jj,1)
311            ENDIF
312         ENDDO
313      ENDDO
314
315      !! Jpalm (08-2014)
316      !! DMS surface concentration calculation
317      !! initialy added for UKESM1 model.
318      !! using MET-OFFICE subroutine.
319      !! DMS module only needs Chl concentration and MLD
320      !! to get an aproximate value of DMS concentration.
321      !! air-sea fluxes are calculated by atmospheric chemitry model
322      !! from atm and oc-surface concentrations.
323      !!
324      !! AXY (13/03/15): this is amended to calculate all of the DMS
325      !!                 estimates examined during UKESM1 (see
326      !!                 comments in trcdms_medusa.F90)
327      !!
328      !! AXY (25/05/17): amended to additionally pass DIN limitation as well as [DIN];
329      !!                 accounts for differences in nutrient half-saturations; changes
330      !!                 also made in trc_dms_medusa; this permits an additional DMS
331      !!                 calculation while retaining the existing Anderson one
332      !!
333      IF (jdms == 1) THEN
334         DO jj = 2,jpjm1
335            DO ji = 2,jpim1
336               if (tmask(ji,jj,1) == 1) then
337                  !! calculate weighted half-saturation for DIN uptake
338                  dms_wtkn(ji,jj) = ((zphn(ji,jj) * xnln) +            &
339                                     (zphd(ji,jj) * xnld)) /           &
340                                     (zphn(ji,jj) + zphd(ji,jj))       
341                  !!
342                  !! feed in correct inputs
343                  if (jdms_input .eq. 0) then
344                     !! use instantaneous inputs
345                     dms_nlim(ji,jj) = zdin(ji,jj) / (zdin(ji,jj) + dms_wtkn(ji,jj))
346                     !!
347                     CALL trc_dms_medusa(zchn(ji,jj),zchd(ji,jj),             &
348                                         hmld(ji,jj),qsr(ji,jj),              &
349                                         zdin(ji,jj), dms_nlim(ji,jj),        &
350                                         dms_andr,dms_simo,dms_aran,dms_hall, & 
351                                         dms_andm)
352                  else
353                     !! use diel-average inputs
354                     dms_nlim(ji,jj) = zn_dms_din(ji,jj) /                    &
355                                      (zn_dms_din(ji,jj) + dms_wtkn(ji,jj))
356                     !!
357                     CALL trc_dms_medusa(zn_dms_chn(ji,jj),zn_dms_chd(ji,jj), &
358                                         zn_dms_mld(ji,jj),zn_dms_qsr(ji,jj), &
359                                         zn_dms_din(ji,jj),dms_nlim(ji,jj),   &
360                                         dms_andr,dms_simo,dms_aran,dms_hall, & 
361                                         dms_andm)
362                  endif
363                  !!
364                  !! assign correct output to variable passed to atmosphere
365                  if (jdms_model .eq. 1) then
366                     dms_surf = dms_andr
367                  elseif (jdms_model .eq. 2) then
368                     dms_surf = dms_simo
369                  elseif (jdms_model .eq. 3) then
370                     dms_surf = dms_aran
371                  elseif (jdms_model .eq. 4) then
372                     dms_surf = dms_hall
373                  elseif (jdms_model .eq. 5) then
374                     dms_surf = dms_andm
375                  endif
376                  !!
377                  !! 2D diag through iom_use
378                  IF( med_diag%DMS_SURF%dgsave ) THEN
379                     dms_surf2d(ji,jj) = dms_surf
380                  ENDIF
381                  IF( med_diag%DMS_ANDR%dgsave ) THEN
382                     dms_andr2d(ji,jj) = dms_andr
383                  ENDIF
384                  IF( med_diag%DMS_SIMO%dgsave ) THEN
385                     dms_simo2d(ji,jj) = dms_simo
386                  ENDIF
387                  IF( med_diag%DMS_ARAN%dgsave ) THEN
388                     dms_aran2d(ji,jj) = dms_aran
389                  ENDIF
390                  IF( med_diag%DMS_HALL%dgsave ) THEN
391                     dms_hall2d(ji,jj) = dms_hall
392                  ENDIF
393                  IF( med_diag%DMS_ANDM%dgsave ) THEN
394                     dms_andm2d(ji,jj) = dms_andm
395                  ENDIF
396               ENDIF
397            ENDDO
398         ENDDO
399#   if defined key_debug_medusa
400         IF (lwp) write (numout,*)                                &
401            'air-sea: finish calculating dms kt = ',kt
402            CALL flush(numout)
403#   endif
404      ENDIF                  !! End IF (jdms == 1)
405
406      !!
407      !! store 2D outputs
408      !!
409      !! JPALM -- 17-11-16 -- put fgco2 out of diag request
410      !!       is needed for coupling; pass through restart
411      DO jj = 2,jpjm1
412         DO ji = 2,jpim1
413            if (tmask(ji,jj,1) == 1) then
414               !! IF( med_diag%FGCO2%dgsave ) THEN
415               !! convert from  mol/m2/day to kg/m2/s
416               !! mmol-C/m3/d -> kg-CO2/m2/s
417               fgco2(ji,jj) = f_co2flux(ji,jj) * fse3t(ji,jj,1) *             &
418                              CO2flux_conv
419               !! ENDIF
420#   if defined key_foam_medusa
421               !! DAF (Aug 2017): Save pCO2 and fCO2 for observation operator
422               f2_pco2w(ji,jj) = f_pco2w(ji,jj)
423               f2_fco2w(ji,jj) = f_pco2w(ji,jj)
424#   endif
425               IF ( lk_iomput ) THEN
426                  IF( med_diag%ATM_PCO2%dgsave ) THEN
427                     f_pco2a2d(ji,jj) = f_pco2atm(ji,jj)
428                  ENDIF
429                  IF( med_diag%OCN_PCO2%dgsave ) THEN
430                     f_pco2w2d(ji,jj) = f_pco2w(ji,jj)
431                  ENDIF
432                  IF( med_diag%CO2FLUX%dgsave ) THEN
433                     !! mmol/m3/d -> mmol/m2/d
434                     f_co2flux2d(ji,jj) = f_co2flux(ji,jj) * fse3t(ji,jj,1)
435                  ENDIF
436                  IF( med_diag%TCO2%dgsave ) THEN
437                     f_TDIC2d(ji,jj) = f_TDIC(ji,jj)
438                  ENDIF
439                  IF( med_diag%TALK%dgsave ) THEN
440                     f_TALK2d(ji,jj) = f_TALK(ji,jj)
441                  ENDIF
442                  IF( med_diag%KW660%dgsave ) THEN
443                      f_kw6602d(ji,jj) = f_kw660(ji,jj)
444                  ENDIF
445                  IF( med_diag%ATM_PP0%dgsave ) THEN
446                      f_pp02d(ji,jj) = f_pp0(ji,jj)
447                  ENDIF
448                  IF( med_diag%O2FLUX%dgsave ) THEN
449                     f_o2flux2d(ji,jj) = f_o2flux(ji,jj)
450                  ENDIF
451                  IF( med_diag%O2SAT%dgsave ) THEN
452                     f_o2sat2d(ji,jj) = f_o2sat(ji,jj)
453                  ENDIF
454                  !! AXY (24/11/16): add in extra MOCSY diagnostics
455                  IF( med_diag%ATM_XCO2%dgsave ) THEN
456                     f_xco2a_2d(ji,jj) = f_xco2a(ji,jj)
457                  ENDIF
458                  IF( med_diag%OCN_FCO2%dgsave ) THEN
459                     f_fco2w_2d(ji,jj) = f_fco2w(ji,jj)
460                  ENDIF
461                  IF( med_diag%ATM_FCO2%dgsave ) THEN
462                     f_fco2a_2d(ji,jj) = f_fco2atm(ji,jj)
463                  ENDIF
464                  IF( med_diag%OCN_RHOSW%dgsave ) THEN
465                     f_ocnrhosw_2d(ji,jj) = f_rhosw(ji,jj)
466                  ENDIF
467                  IF( med_diag%OCN_SCHCO2%dgsave ) THEN
468                     f_ocnschco2_2d(ji,jj) = f_schmidtco2(ji,jj)
469                  ENDIF
470                  IF( med_diag%OCN_KWCO2%dgsave ) THEN
471                     f_ocnkwco2_2d(ji,jj) = f_kwco2(ji,jj)
472                  ENDIF
473                  IF( med_diag%OCN_K0%dgsave ) THEN
474                     f_ocnk0_2d(ji,jj) = f_K0(ji,jj)
475                  ENDIF
476                  IF( med_diag%CO2STARAIR%dgsave ) THEN
477                     f_co2starair_2d(ji,jj) = f_co2starair(ji,jj)
478                  ENDIF
479                  IF( med_diag%OCN_DPCO2%dgsave ) THEN
480                     f_ocndpco2_2d(ji,jj) = f_dpco2(ji,jj)
481                  ENDIF
482               ENDIF
483            ENDIF
484         ENDDO
485      ENDDO
486
487# endif
488
489      !!-----------------------------------------------------------------
490      !! River inputs
491      !!-----------------------------------------------------------------
492      DO jj = 2,jpjm1
493         DO ji = 2,jpim1
494            !! OPEN wet point IF..THEN loop
495            if (tmask(ji,jj,1) == 1) then
496               !!
497               !! runoff comes in as        kg / m2 / s
498               !! used and written out as   m3 / m2 / d (= m / d)
499               !! where                     1000 kg / m2 / d =
500               !!                             1 m3 / m2 / d = 1 m / d
501               !!
502               !! AXY (17/07/14): the compiler doesn't like this line for
503               !!                 some reason; as MEDUSA doesn't even use
504               !!                 runoff for riverine inputs, a temporary
505               !!                 solution is to switch off runoff entirely
506               !!                 here; again, this change is one of several
507               !!                 that will need revisiting once MEDUSA has
508               !!                 bedded down in UKESM1; particularly so if
509               !!                 the land scheme provides information
510               !!                 concerning nutrient fluxes
511               !!
512               !! f_runoff(ji,jj) = sf_rnf(1)%fnow(ji,jj,1) / 1000. * 60. *   &
513               !!                   60. * 24.
514               f_runoff(ji,jj) = 0.0
515               !!
516               !! nutrients are added via rivers to the model in one of
517               !! two ways:
518               !!   1. via river concentration; i.e. the average nutrient
519               !!      concentration of a river water is described by a
520               !!      spatial file, and this is multiplied by runoff to
521               !!      give a nutrient flux
522               !!   2. via direct river flux; i.e. the average nutrient
523               !!      flux due to rivers is described by a spatial file,
524               !!      and this is simply applied as a direct nutrient
525               !!      flux (i.e. it does not relate or respond to model
526               !!      runoff) nutrient fields are derived from the
527               !!      GlobalNEWS 2 database; carbon and alkalinity are
528               !!      derived from continent-scale DIC estimates (Huang et
529               !!      al., 2012) and some Arctic river alkalinity
530               !!      estimates (Katya?)
531               !!
532               !! as of 19/07/12, riverine nutrients can now be spread
533               !! vertically across several grid cells rather than just
534               !! poured into the surface box; this block of code is still
535               !! executed, however, to set up the total amounts of
536               !! nutrient entering via rivers
537               !!
538               !! nitrogen
539               if (jriver_n .eq. 1) then
540                  !! river concentration specified; use runoff to
541                  !! calculate input
542                  f_riv_n(ji,jj) = f_runoff(ji,jj) * riv_n(ji,jj)
543               elseif (jriver_n .eq. 2) then
544                  !! river flux specified; independent of runoff
545                  f_riv_n(ji,jj) = riv_n(ji,jj)
546               endif
547               !!
548               !! silicon
549               if (jriver_si .eq. 1) then
550                  !! river concentration specified; use runoff to
551                  !! calculate input
552                  f_riv_si(ji,jj) = f_runoff(ji,jj) * riv_si(ji,jj)
553               elseif (jriver_si .eq. 2) then
554                  !! river flux specified; independent of runoff
555                  f_riv_si(ji,jj) = riv_si(ji,jj)
556               endif
557               !!
558               !! carbon
559               if (jriver_c .eq. 1) then
560                  !! river concentration specified; use runoff to
561                  !! calculate input
562                  f_riv_c(ji,jj) = f_runoff(ji,jj) * riv_c(ji,jj)
563               elseif (jriver_c .eq. 2) then
564                  !! river flux specified; independent of runoff
565                  f_riv_c(ji,jj) = riv_c(ji,jj)
566               endif
567               !!
568               !! alkalinity
569               if (jriver_alk .eq. 1) then
570                  !! river concentration specified; use runoff to
571                  !! calculate input
572                  f_riv_alk(ji,jj) = f_runoff(ji,jj) * riv_alk(ji,jj)
573               elseif (jriver_alk .eq. 2) then
574                  !! river flux specified; independent of runoff
575                  f_riv_alk(ji,jj) = riv_alk(ji,jj)
576               endif
577            ENDIF
578         ENDDO
579      ENDDO
580
581   END SUBROUTINE air_sea
582
583#else
584   !!======================================================================
585   !!  Dummy module :                                   No MEDUSA bio-model
586   !!======================================================================
587CONTAINS
588   SUBROUTINE air_sea( )                    ! Empty routine
589      WRITE(*,*) 'air_sea: You should not have seen this print! error?'
590   END SUBROUTINE air_sea
591#endif 
592
593   !!======================================================================
594END MODULE air_sea_mod
Note: See TracBrowser for help on using the repository browser.