1 | MODULE diaar5 |
---|
2 | !!====================================================================== |
---|
3 | !! *** MODULE diaar5 *** |
---|
4 | !! AR5 diagnostics |
---|
5 | !!====================================================================== |
---|
6 | !! History : 3.2 ! 2009-11 (S. Masson) Original code |
---|
7 | !! 3.3 ! 2010-10 (C. Ethe, G. Madec) reorganisation of initialisation phase + merge TRC-TRA |
---|
8 | !!---------------------------------------------------------------------- |
---|
9 | #if defined key_diaar5 || defined key_esopa |
---|
10 | !!---------------------------------------------------------------------- |
---|
11 | !! 'key_diaar5' : activate ar5 diagnotics |
---|
12 | !!---------------------------------------------------------------------- |
---|
13 | !! dia_ar5 : AR5 diagnostics |
---|
14 | !! dia_ar5_init : initialisation of AR5 diagnostics |
---|
15 | !!---------------------------------------------------------------------- |
---|
16 | USE oce ! ocean dynamics and active tracers |
---|
17 | USE dom_oce ! ocean space and time domain |
---|
18 | USE eosbn2 ! equation of state (eos_bn2 routine) |
---|
19 | USE lib_mpp ! distribued memory computing library |
---|
20 | USE iom ! I/O manager library |
---|
21 | USE timing ! preformance summary |
---|
22 | USE wrk_nemo ! working arrays |
---|
23 | USE fldread ! type FLD_N |
---|
24 | USE phycst ! physical constant |
---|
25 | USE in_out_manager ! I/O manager |
---|
26 | |
---|
27 | IMPLICIT NONE |
---|
28 | PRIVATE |
---|
29 | |
---|
30 | PUBLIC dia_ar5 ! routine called in step.F90 module |
---|
31 | PUBLIC dia_ar5_init ! routine called in opa.F90 module |
---|
32 | PUBLIC dia_ar5_alloc ! routine called in nemogcm.F90 module |
---|
33 | |
---|
34 | LOGICAL, PUBLIC, PARAMETER :: lk_diaar5 = .TRUE. ! coupled flag |
---|
35 | |
---|
36 | REAL(wp) :: vol0 ! ocean volume (interior domain) |
---|
37 | REAL(wp) :: area_tot ! total ocean surface (interior domain) |
---|
38 | REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,: ) :: area ! cell surface (interior domain) |
---|
39 | REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,: ) :: thick0 ! ocean thickness (interior domain) |
---|
40 | REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: sn0 ! initial salinity |
---|
41 | |
---|
42 | !! * Substitutions |
---|
43 | # include "domzgr_substitute.h90" |
---|
44 | !!---------------------------------------------------------------------- |
---|
45 | !! NEMO/OPA 3.3 , NEMO Consortium (2010) |
---|
46 | !! $Id$ |
---|
47 | !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) |
---|
48 | !!---------------------------------------------------------------------- |
---|
49 | CONTAINS |
---|
50 | |
---|
51 | FUNCTION dia_ar5_alloc() |
---|
52 | !!---------------------------------------------------------------------- |
---|
53 | !! *** ROUTINE dia_ar5_alloc *** |
---|
54 | !!---------------------------------------------------------------------- |
---|
55 | INTEGER :: dia_ar5_alloc |
---|
56 | !!---------------------------------------------------------------------- |
---|
57 | ! |
---|
58 | ALLOCATE( area(jpi,jpj), thick0(jpi,jpj) , sn0(jpi,jpj,jpk) , STAT=dia_ar5_alloc ) |
---|
59 | ! |
---|
60 | IF( lk_mpp ) CALL mpp_sum ( dia_ar5_alloc ) |
---|
61 | IF( dia_ar5_alloc /= 0 ) CALL ctl_warn('dia_ar5_alloc: failed to allocate arrays') |
---|
62 | ! |
---|
63 | END FUNCTION dia_ar5_alloc |
---|
64 | |
---|
65 | |
---|
66 | SUBROUTINE dia_ar5( kt ) |
---|
67 | !!---------------------------------------------------------------------- |
---|
68 | !! *** ROUTINE dia_ar5 *** |
---|
69 | !! |
---|
70 | !! ** Purpose : compute and output some AR5 diagnostics |
---|
71 | !!---------------------------------------------------------------------- |
---|
72 | ! |
---|
73 | INTEGER, INTENT( in ) :: kt ! ocean time-step index |
---|
74 | ! |
---|
75 | INTEGER :: ji, jj, jk ! dummy loop arguments |
---|
76 | REAL(wp) :: zvolssh, zvol, zssh_steric, zztmp, zarho, ztemp, zsal, zmass |
---|
77 | ! |
---|
78 | REAL(wp), POINTER, DIMENSION(:,:) :: zarea_ssh , zbotpres ! 2D workspace |
---|
79 | REAL(wp), POINTER, DIMENSION(:,:,:) :: zrhd , zrhop ! 3D workspace |
---|
80 | REAL(wp), POINTER, DIMENSION(:,:,:,:) :: ztsn ! 4D workspace |
---|
81 | !!-------------------------------------------------------------------- |
---|
82 | IF( nn_timing == 1 ) CALL timing_start('dia_ar5') |
---|
83 | |
---|
84 | CALL wrk_alloc( jpi , jpj , zarea_ssh , zbotpres ) |
---|
85 | CALL wrk_alloc( jpi , jpj , jpk , zrhd , zrhop ) |
---|
86 | CALL wrk_alloc( jpi , jpj , jpk , jpts , ztsn ) |
---|
87 | |
---|
88 | zarea_ssh(:,:) = area(:,:) * sshn(:,:) |
---|
89 | |
---|
90 | ! ! total volume of liquid seawater |
---|
91 | zvolssh = SUM( zarea_ssh(:,:) ) |
---|
92 | IF( lk_mpp ) CALL mpp_sum( zvolssh ) |
---|
93 | zvol = vol0 + zvolssh |
---|
94 | |
---|
95 | CALL iom_put( 'voltot', zvol ) |
---|
96 | CALL iom_put( 'sshtot', zvolssh / area_tot ) |
---|
97 | |
---|
98 | ! |
---|
99 | ztsn(:,:,:,jp_tem) = tsn(:,:,:,jp_tem) ! thermosteric ssh |
---|
100 | ztsn(:,:,:,jp_sal) = sn0(:,:,:) |
---|
101 | CALL eos( ztsn, zrhd, fsdept_n(:,:,:) ) ! now in situ density using initial salinity |
---|
102 | ! |
---|
103 | zbotpres(:,:) = 0._wp ! no atmospheric surface pressure, levitating sea-ice |
---|
104 | DO jk = 1, jpkm1 |
---|
105 | zbotpres(:,:) = zbotpres(:,:) + fse3t(:,:,jk) * zrhd(:,:,jk) |
---|
106 | END DO |
---|
107 | IF( .NOT.lk_vvl ) THEN |
---|
108 | IF ( ln_isfcav ) THEN |
---|
109 | DO ji=1,jpi |
---|
110 | DO jj=1,jpj |
---|
111 | zbotpres(ji,jj) = zbotpres(ji,jj) + sshn(ji,jj) * zrhd(ji,jj,mikt(ji,jj)) + riceload(ji,jj) |
---|
112 | END DO |
---|
113 | END DO |
---|
114 | ELSE |
---|
115 | zbotpres(:,:) = zbotpres(:,:) + sshn(:,:) * zrhd(:,:,1) |
---|
116 | END IF |
---|
117 | END IF |
---|
118 | ! |
---|
119 | zarho = SUM( area(:,:) * zbotpres(:,:) ) |
---|
120 | IF( lk_mpp ) CALL mpp_sum( zarho ) |
---|
121 | zssh_steric = - zarho / area_tot |
---|
122 | CALL iom_put( 'sshthster', zssh_steric ) |
---|
123 | |
---|
124 | ! ! steric sea surface height |
---|
125 | CALL eos( tsn, zrhd, zrhop, fsdept_n(:,:,:) ) ! now in situ and potential density |
---|
126 | zrhop(:,:,jpk) = 0._wp |
---|
127 | CALL iom_put( 'rhop', zrhop ) |
---|
128 | ! |
---|
129 | zbotpres(:,:) = 0._wp ! no atmospheric surface pressure, levitating sea-ice |
---|
130 | DO jk = 1, jpkm1 |
---|
131 | zbotpres(:,:) = zbotpres(:,:) + fse3t(:,:,jk) * zrhd(:,:,jk) |
---|
132 | END DO |
---|
133 | IF( .NOT.lk_vvl ) THEN |
---|
134 | IF ( ln_isfcav ) THEN |
---|
135 | DO ji=1,jpi |
---|
136 | DO jj=1,jpj |
---|
137 | zbotpres(ji,jj) = zbotpres(ji,jj) + sshn(ji,jj) * zrhd(ji,jj,mikt(ji,jj)) + riceload(ji,jj) |
---|
138 | END DO |
---|
139 | END DO |
---|
140 | ELSE |
---|
141 | zbotpres(:,:) = zbotpres(:,:) + sshn(:,:) * zrhd(:,:,1) |
---|
142 | END IF |
---|
143 | END IF |
---|
144 | ! |
---|
145 | zarho = SUM( area(:,:) * zbotpres(:,:) ) |
---|
146 | IF( lk_mpp ) CALL mpp_sum( zarho ) |
---|
147 | zssh_steric = - zarho / area_tot |
---|
148 | CALL iom_put( 'sshsteric', zssh_steric ) |
---|
149 | |
---|
150 | ! ! ocean bottom pressure |
---|
151 | zztmp = rau0 * grav * 1.e-4_wp ! recover pressure from pressure anomaly and cover to dbar = 1.e4 Pa |
---|
152 | zbotpres(:,:) = zztmp * ( zbotpres(:,:) + sshn(:,:) + thick0(:,:) ) |
---|
153 | CALL iom_put( 'botpres', zbotpres ) |
---|
154 | |
---|
155 | ! ! Mean density anomalie, temperature and salinity |
---|
156 | ztemp = 0._wp |
---|
157 | zsal = 0._wp |
---|
158 | DO jk = 1, jpkm1 |
---|
159 | DO jj = 1, jpj |
---|
160 | DO ji = 1, jpi |
---|
161 | zztmp = area(ji,jj) * fse3t(ji,jj,jk) |
---|
162 | ztemp = ztemp + zztmp * tsn(ji,jj,jk,jp_tem) |
---|
163 | zsal = zsal + zztmp * tsn(ji,jj,jk,jp_sal) |
---|
164 | END DO |
---|
165 | END DO |
---|
166 | END DO |
---|
167 | IF( .NOT.lk_vvl ) THEN |
---|
168 | IF ( ln_isfcav ) THEN |
---|
169 | DO ji=1,jpi |
---|
170 | DO jj=1,jpj |
---|
171 | ztemp = ztemp + zarea_ssh(ji,jj) * tsn(ji,jj,mikt(ji,jj),jp_tem) |
---|
172 | zsal = zsal + zarea_ssh(ji,jj) * tsn(ji,jj,mikt(ji,jj),jp_sal) |
---|
173 | END DO |
---|
174 | END DO |
---|
175 | ELSE |
---|
176 | ztemp = ztemp + SUM( zarea_ssh(:,:) * tsn(:,:,1,jp_tem) ) |
---|
177 | zsal = zsal + SUM( zarea_ssh(:,:) * tsn(:,:,1,jp_sal) ) |
---|
178 | END IF |
---|
179 | ENDIF |
---|
180 | IF( lk_mpp ) THEN |
---|
181 | CALL mpp_sum( ztemp ) |
---|
182 | CALL mpp_sum( zsal ) |
---|
183 | END IF |
---|
184 | ! |
---|
185 | zmass = rau0 * ( zarho + zvol ) ! total mass of liquid seawater |
---|
186 | ztemp = ztemp / zvol ! potential temperature in liquid seawater |
---|
187 | zsal = zsal / zvol ! Salinity of liquid seawater |
---|
188 | ! |
---|
189 | CALL iom_put( 'masstot', zmass ) |
---|
190 | CALL iom_put( 'temptot', ztemp ) |
---|
191 | CALL iom_put( 'saltot' , zsal ) |
---|
192 | ! |
---|
193 | CALL wrk_dealloc( jpi , jpj , zarea_ssh , zbotpres ) |
---|
194 | CALL wrk_dealloc( jpi , jpj , jpk , zrhd , zrhop ) |
---|
195 | CALL wrk_dealloc( jpi , jpj , jpk , jpts , ztsn ) |
---|
196 | ! |
---|
197 | IF( nn_timing == 1 ) CALL timing_stop('dia_ar5') |
---|
198 | ! |
---|
199 | END SUBROUTINE dia_ar5 |
---|
200 | |
---|
201 | |
---|
202 | SUBROUTINE dia_ar5_init |
---|
203 | !!---------------------------------------------------------------------- |
---|
204 | !! *** ROUTINE dia_ar5_init *** |
---|
205 | !! |
---|
206 | !! ** Purpose : initialization for AR5 diagnostic computation |
---|
207 | !!---------------------------------------------------------------------- |
---|
208 | INTEGER :: inum |
---|
209 | INTEGER :: ik |
---|
210 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
211 | REAL(wp) :: zztmp |
---|
212 | REAL(wp), POINTER, DIMENSION(:,:,:,:) :: zsaldta ! Jan/Dec levitus salinity |
---|
213 | ! reading initial file |
---|
214 | LOGICAL :: ln_tsd_init !: T & S data flag |
---|
215 | LOGICAL :: ln_tsd_tradmp !: internal damping toward input data flag |
---|
216 | CHARACTER(len=100) :: cn_dir |
---|
217 | TYPE(FLD_N) :: sn_tem,sn_sal |
---|
218 | INTEGER :: ios=0 |
---|
219 | |
---|
220 | NAMELIST/namtsd/ ln_tsd_init,ln_tsd_tradmp,cn_dir,sn_tem,sn_sal |
---|
221 | ! |
---|
222 | |
---|
223 | REWIND( numnam_ref ) ! Namelist namtsd in reference namelist : |
---|
224 | READ ( numnam_ref, namtsd, IOSTAT = ios, ERR = 901) |
---|
225 | 901 IF( ios /= 0 ) CALL ctl_nam ( ios , ' namtsd in reference namelist for dia_ar5', lwp ) |
---|
226 | REWIND( numnam_cfg ) ! Namelist namtsd in configuration namelist : Parameters of the run |
---|
227 | READ ( numnam_cfg, namtsd, IOSTAT = ios, ERR = 902 ) |
---|
228 | 902 IF( ios /= 0 ) CALL ctl_nam ( ios , ' namtsd in configuration namelist for dia_ar5', lwp ) |
---|
229 | IF(lwm) WRITE ( numond, namtsd ) |
---|
230 | ! |
---|
231 | !!---------------------------------------------------------------------- |
---|
232 | ! |
---|
233 | IF( nn_timing == 1 ) CALL timing_start('dia_ar5_init') |
---|
234 | ! |
---|
235 | CALL wrk_alloc( jpi , jpj , jpk, jpts, zsaldta ) |
---|
236 | ! ! allocate dia_ar5 arrays |
---|
237 | IF( dia_ar5_alloc() /= 0 ) CALL ctl_stop( 'STOP', 'dia_ar5_init : unable to allocate arrays' ) |
---|
238 | |
---|
239 | area(:,:) = e1t(:,:) * e2t(:,:) * tmask_i(:,:) |
---|
240 | |
---|
241 | area_tot = SUM( area(:,:) ) ; IF( lk_mpp ) CALL mpp_sum( area_tot ) |
---|
242 | |
---|
243 | vol0 = 0._wp |
---|
244 | thick0(:,:) = 0._wp |
---|
245 | DO jk = 1, jpkm1 |
---|
246 | vol0 = vol0 + SUM( area (:,:) * tmask(:,:,jk) * e3t_0(:,:,jk) ) |
---|
247 | thick0(:,:) = thick0(:,:) + tmask_i(:,:) * tmask(:,:,jk) * e3t_0(:,:,jk) |
---|
248 | END DO |
---|
249 | IF( lk_mpp ) CALL mpp_sum( vol0 ) |
---|
250 | |
---|
251 | CALL iom_open ( TRIM( cn_dir )//TRIM(sn_sal%clname), inum ) |
---|
252 | CALL iom_get ( inum, jpdom_data, TRIM(sn_sal%clvar), zsaldta(:,:,:,1), 1 ) |
---|
253 | CALL iom_get ( inum, jpdom_data, TRIM(sn_sal%clvar), zsaldta(:,:,:,2), 12 ) |
---|
254 | CALL iom_close( inum ) |
---|
255 | sn0(:,:,:) = 0.5_wp * ( zsaldta(:,:,:,1) + zsaldta(:,:,:,2) ) |
---|
256 | sn0(:,:,:) = sn0(:,:,:) * tmask(:,:,:) |
---|
257 | IF( ln_zps ) THEN ! z-coord. partial steps |
---|
258 | DO jj = 1, jpj ! interpolation of salinity at the last ocean level (i.e. the partial step) |
---|
259 | DO ji = 1, jpi |
---|
260 | ik = mbkt(ji,jj) |
---|
261 | IF( ik > 1 ) THEN |
---|
262 | zztmp = ( gdept_1d(ik) - gdept_0(ji,jj,ik) ) / ( gdept_1d(ik) - gdept_1d(ik-1) ) |
---|
263 | sn0(ji,jj,ik) = ( 1._wp - zztmp ) * sn0(ji,jj,ik) + zztmp * sn0(ji,jj,ik-1) |
---|
264 | ENDIF |
---|
265 | END DO |
---|
266 | END DO |
---|
267 | ENDIF |
---|
268 | ! |
---|
269 | CALL wrk_dealloc( jpi , jpj , jpk, jpts, zsaldta ) |
---|
270 | ! |
---|
271 | IF( nn_timing == 1 ) CALL timing_stop('dia_ar5_init') |
---|
272 | ! |
---|
273 | END SUBROUTINE dia_ar5_init |
---|
274 | |
---|
275 | #else |
---|
276 | !!---------------------------------------------------------------------- |
---|
277 | !! Default option : NO diaar5 |
---|
278 | !!---------------------------------------------------------------------- |
---|
279 | LOGICAL, PUBLIC, PARAMETER :: lk_diaar5 = .FALSE. ! coupled flag |
---|
280 | CONTAINS |
---|
281 | SUBROUTINE dia_ar5_init ! Dummy routine |
---|
282 | END SUBROUTINE dia_ar5_init |
---|
283 | SUBROUTINE dia_ar5( kt ) ! Empty routine |
---|
284 | INTEGER :: kt |
---|
285 | WRITE(*,*) 'dia_ar5: You should not have seen this print! error?', kt |
---|
286 | END SUBROUTINE dia_ar5 |
---|
287 | #endif |
---|
288 | |
---|
289 | !!====================================================================== |
---|
290 | END MODULE diaar5 |
---|