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/UKMO/dev_r5518_GO6_under_ice_relax_dr_hook/NEMOGCM/NEMO/TOP_SRC/MEDUSA – NEMO

source: branches/UKMO/dev_r5518_GO6_under_ice_relax_dr_hook/NEMOGCM/NEMO/TOP_SRC/MEDUSA/bio_medusa_diag.F90 @ 11738

Last change on this file since 11738 was 11738, checked in by marc, 5 years ago

The Dr Hook changes from my perl code.

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