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