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.
bio_medusa_diag.F90 in branches/NERC/dev_r5518_GO6_split_trcbiomedusa/NEMOGCM/NEMO/TOP_SRC/MEDUSA – NEMO

source: branches/NERC/dev_r5518_GO6_split_trcbiomedusa/NEMOGCM/NEMO/TOP_SRC/MEDUSA/bio_medusa_diag.F90 @ 8395

Last change on this file since 8395 was 8395, checked in by jpalmier, 7 years ago

JPALM -- GMED #339 - split trcbio_medusa only

File size: 9.4 KB
RevLine 
[8395]1MODULE bio_medusa_diag_mod
2   !!======================================================================
3   !!                         ***  MODULE bio_medusa_diag_mod  ***
4   !! Calculates diagnostics
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   bio_medusa_diag        ! 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 bio_medusa_diag( jk )
28      !!-------------------------------------------------------------------
29      !!                     ***  ROUTINE bio_medusa_diag  ***
30      !! This called from TRC_BIO_MEDUSA and calculates diagnostics
31      !!-------------------------------------------------------------------
32      USE bio_med_diag_iomput_mod,  ONLY: bio_med_diag_iomput
33      USE bio_med_diag_trc_mod,     ONLY: bio_med_diag_trc
34      USE bio_medusa_mod
35      USE dom_oce,                  ONLY: e3t_0, e3t_n,                  &
36                                          gdepw_0, gdepw_n, tmask
37      USE in_out_manager,           ONLY: lwp, numout
38# if defined key_iomput
39      USE iom,                      ONLY: lk_iomput
40# endif
41      USE par_oce,                  ONLY: jpim1, jpjm1
42      USE sms_medusa,               ONLY: xrfn, xthetapd, xthetapn,      &
43                                          xthetazme, xthetazmi
44      USE trc,                      ONLY: ln_diatrc, med_diag 
45# if defined key_roam
46      USE trcoxy_medusa,            ONLY: oxy_sato
47# endif
48
49   !!* Substitution
50#  include "domzgr_substitute.h90"
51
52      !! level
53      INTEGER, INTENT( in ) :: jk
54
55      !! Loop avariables
56      INTEGER :: ji, jj, jn
57
58# if defined key_trc_diabio
59      !!==========================================================
60      !! LOCAL GRID CELL DIAGNOSTICS
61      !!==========================================================
62      !!
63      !!----------------------------------------------------------
64      !! Full diagnostics key_trc_diabio:
65      !! LOBSTER and PISCES support full diagnistics option
66      !! key_trc_diabio which gives an option of FULL output of
67      !! biological sourses and sinks. I cannot see any reason
68      !! for doing this. If needed, it can be done as shown
69      !! below.
70      !!----------------------------------------------------------
71      !!
72      IF(lwp) WRITE(numout,*) ' MEDUSA does not support key_trc_diabio'
73# endif
74
75      !! The section below, down to calculation of zo2min, was moved
76      !! from before the call to AIR_SEA in trcbio_medusa.F90 - marc 9/5/17
77      IF( lk_iomput ) THEN
78         DO jj = 2,jpjm1
79            DO ji = 2,jpim1
80               if (tmask(ji,jj,jk) == 1) then
81                  !! sum tracers for inventory checks
82                  IF ( med_diag%INVTN%dgsave )   THEN
83                     ftot_n(ji,jj)  = ftot_n(ji,jj) +                        &
84                        (fse3t(ji,jj,jk) * (zphn(ji,jj) + zphd(ji,jj) +      &
85                                            zzmi(ji,jj) + zzme(ji,jj) +      &
86                                            zdet(ji,jj) + zdin(ji,jj)))
87                  ENDIF
88                  IF ( med_diag%INVTSI%dgsave )  THEN
89                     ftot_si(ji,jj) = ftot_si(ji,jj) +                       & 
90                       (fse3t(ji,jj,jk) * (zpds(ji,jj) + zsil(ji,jj)))
91                  ENDIF
92                  IF ( med_diag%INVTFE%dgsave )  THEN
93                     ftot_fe(ji,jj) = ftot_fe(ji,jj) +                       & 
94                        (fse3t(ji,jj,jk) * (xrfn *                           &
95                                            (zphn(ji,jj) + zphd(ji,jj) +     &
96                                             zzmi(ji,jj) + zzme(ji,jj) +     &
97                                             zdet(ji,jj)) +                  &
98                                            zfer(ji,jj)))
99                  ENDIF
100               ENDIF
101            ENDDO
102         ENDDO
103
104# if defined key_roam
105         DO jj = 2,jpjm1
106            DO ji = 2,jpim1
107               if (tmask(ji,jj,jk) == 1) then
108                  IF ( med_diag%INVTC%dgsave )  THEN
109                     ftot_c(ji,jj)  = ftot_c(ji,jj) +                        & 
110                        (fse3t(ji,jj,jk) * ((xthetapn * zphn(ji,jj)) +       &
111                                            (xthetapd * zphd(ji,jj)) +       &
112                                            (xthetazmi * zzmi(ji,jj)) +      &
113                                            (xthetazme * zzme(ji,jj)) +      &
114                                            zdtc(ji,jj) + zdic(ji,jj)))
115                  ENDIF
116                  IF ( med_diag%INVTALK%dgsave ) THEN
117                     ftot_a(ji,jj)  = ftot_a(ji,jj) + (fse3t(ji,jj,jk) *     &
118                                                       zalk(ji,jj))
119                  ENDIF
120                  IF ( med_diag%INVTO2%dgsave )  THEN
121                     ftot_o2(ji,jj) = ftot_o2(ji,jj) + (fse3t(ji,jj,jk) *    &
122                                                        zoxy(ji,jj))
123                  ENDIF
124               ENDIF
125            ENDDO
126         ENDDO
127
128         DO jj = 2,jpjm1
129            DO ji = 2,jpim1
130               if (tmask(ji,jj,jk) == 1) then
131                  IF ( med_diag%INVTC%dgsave )  THEN
132                     !!
133                     !! AXY (10/11/16): CMIP6 diagnostics
134                     IF ( med_diag%INTDISSIC%dgsave ) THEN
135                        intdissic(ji,jj) = intdissic(ji,jj) +                &
136                                           (fse3t(ji,jj,jk) * zdic(ji,jj))
137                     ENDIF
138                     IF ( med_diag%INTDISSIN%dgsave ) THEN
139                        intdissin(ji,jj) = intdissin(ji,jj) +                &
140                                           (fse3t(ji,jj,jk) * zdin(ji,jj))
141                     ENDIF
142                     IF ( med_diag%INTDISSISI%dgsave ) THEN
143                        intdissisi(ji,jj) = intdissisi(ji,jj) +              &
144                                            (fse3t(ji,jj,jk) * zsil(ji,jj))
145                     ENDIF
146                     IF ( med_diag%INTTALK%dgsave ) THEN
147                        inttalk(ji,jj) = inttalk(ji,jj) +                    &
148                                         (fse3t(ji,jj,jk) * zalk(ji,jj))
149                     ENDIF
150                  ENDIF
151               ENDIF
152            ENDDO
153         ENDDO
154
155         DO jj = 2,jpjm1
156            DO ji = 2,jpim1
157               if (tmask(ji,jj,jk) == 1) then
158                  IF ( med_diag%O2min%dgsave ) THEN
159                     if ( zoxy(ji,jj) < o2min(ji,jj) ) then
160                        o2min(ji,jj)  = zoxy(ji,jj)
161                        IF ( med_diag%ZO2min%dgsave ) THEN
162                           !! layer midpoint
163                           zo2min(ji,jj) = (fsdepw(ji,jj,jk) +               &
164                                            fdep1(ji,jj)) / 2.0
165                        ENDIF
166                     endif
167                  ENDIF
168               ENDIF
169            ENDDO
170         ENDDO
171# endif
172      ENDIF
173
174# if defined key_roam
175      !! This section is moved from just below CALL to AIR_SEA in
176      !! trcbio_medusa.F90 - marc 9/5/17
177      DO jj = 2,jpjm1
178         DO ji = 2,jpim1
179            !! OPEN wet point IF..THEN loop
180            if (tmask(ji,jj,jk) == 1) then
181
182               !! AXY (11/11/16): CMIP6 oxygen saturation 3D diagnostic
183               IF ( med_diag%O2SAT3%dgsave ) THEN
184! Remove f_o2sat3 - marc 9/5/17
185!                  call oxy_sato( ztmp(ji,jj), zsal(ji,jj), f_o2sat3 )
186!                  o2sat3(ji, jj, jk) = f_o2sat3
187                  call oxy_sato( ztmp(ji,jj), zsal(ji,jj),                   &
188                                 o2sat3(ji,jj,jk) )
189               ENDIF
190            ENDIF
191         ENDDO
192      ENDDO
193# endif
194
195      IF( lk_iomput  .AND.  .NOT.  ln_diatrc  ) THEN
196
197         !!---------------------------------------------------------------
198         !! Calculates the diagnostics used with iom_put
199         !!---------------------------------------------------------------
200         CALL bio_med_diag_iomput( jk )
201
202      ELSE IF( ln_diatrc ) THEN
203
204         !!---------------------------------------------------------------
205         !! The diagnostics without using iom_use
206         !!---------------------------------------------------------------
207         CALL bio_med_diag_trc( jk )
208
209      ENDIF
210
211   END SUBROUTINE bio_medusa_diag
212
213#else
214   !!======================================================================
215   !!  Dummy module :                                   No MEDUSA bio-model
216   !!======================================================================
217CONTAINS
218   SUBROUTINE bio_medusa_diag( )                    ! Empty routine
219      WRITE(*,*) 'bio_medusa_diag: You should not have seen this print! error?'
220   END SUBROUTINE bio_medusa_diag
221#endif 
222
223   !!======================================================================
224END MODULE bio_medusa_diag_mod
Note: See TracBrowser for help on using the repository browser.