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 | |
---|
24 | IMPLICIT NONE |
---|
25 | PRIVATE |
---|
26 | |
---|
27 | PUBLIC dia_ar5 ! routine called in step.F90 module |
---|
28 | PUBLIC dia_ar5_init ! routine called in opa.F90 module |
---|
29 | PUBLIC dia_ar5_alloc ! routine called in nemogcm.F90 module |
---|
30 | |
---|
31 | LOGICAL, PUBLIC, PARAMETER :: lk_diaar5 = .TRUE. ! coupled flag |
---|
32 | |
---|
33 | REAL(wp) :: vol0 ! ocean volume (interior domain) |
---|
34 | REAL(wp) :: area_tot ! total ocean surface (interior domain) |
---|
35 | REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,: ) :: area ! cell surface (interior domain) |
---|
36 | REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,: ) :: thick0 ! ocean thickness (interior domain) |
---|
37 | REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: sn0 ! initial salinity |
---|
38 | |
---|
39 | !! * Substitutions |
---|
40 | # include "domzgr_substitute.h90" |
---|
41 | !!---------------------------------------------------------------------- |
---|
42 | !! NEMO/OPA 3.3 , NEMO Consortium (2010) |
---|
43 | !! $Id$ |
---|
44 | !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) |
---|
45 | !!---------------------------------------------------------------------- |
---|
46 | CONTAINS |
---|
47 | |
---|
48 | FUNCTION dia_ar5_alloc() |
---|
49 | !!---------------------------------------------------------------------- |
---|
50 | !! *** ROUTINE dia_ar5_alloc *** |
---|
51 | !!---------------------------------------------------------------------- |
---|
52 | INTEGER :: dia_ar5_alloc |
---|
53 | !!---------------------------------------------------------------------- |
---|
54 | ! |
---|
55 | ALLOCATE( area(jpi,jpj), thick0(jpi,jpj) , sn0(jpi,jpj,jpk) , STAT=dia_ar5_alloc ) |
---|
56 | ! |
---|
57 | IF( lk_mpp ) CALL mpp_sum ( dia_ar5_alloc ) |
---|
58 | IF( dia_ar5_alloc /= 0 ) CALL ctl_warn('dia_ar5_alloc: failed to allocate arrays') |
---|
59 | ! |
---|
60 | END FUNCTION dia_ar5_alloc |
---|
61 | |
---|
62 | |
---|
63 | SUBROUTINE dia_ar5( kt ) |
---|
64 | !!---------------------------------------------------------------------- |
---|
65 | !! *** ROUTINE dia_ar5 *** |
---|
66 | !! |
---|
67 | !! ** Purpose : compute and output some AR5 diagnostics |
---|
68 | !!---------------------------------------------------------------------- |
---|
69 | ! |
---|
70 | INTEGER, INTENT( in ) :: kt ! ocean time-step index |
---|
71 | ! |
---|
72 | INTEGER :: ji, jj, jk ! dummy loop arguments |
---|
73 | REAL(wp) :: zvolssh, zvol, zssh_steric, zztmp, zarho, ztemp, zsal, zmass |
---|
74 | ! |
---|
75 | REAL(wp), POINTER, DIMENSION(:,:) :: zarea_ssh , zbotpres ! 2D workspace |
---|
76 | REAL(wp), POINTER, DIMENSION(:,:,:) :: zrhd , zrhop ! 3D workspace |
---|
77 | REAL(wp), POINTER, DIMENSION(:,:,:,:) :: ztsn ! 4D workspace |
---|
78 | !!-------------------------------------------------------------------- |
---|
79 | IF( nn_timing == 1 ) CALL timing_start('dia_ar5') |
---|
80 | |
---|
81 | CALL wrk_alloc( jpi , jpj , zarea_ssh , zbotpres ) |
---|
82 | CALL wrk_alloc( jpi , jpj , jpk , zrhd , zrhop ) |
---|
83 | CALL wrk_alloc( jpi , jpj , jpk , jpts , ztsn ) |
---|
84 | |
---|
85 | zarea_ssh(:,:) = area(:,:) * sshn(:,:) |
---|
86 | |
---|
87 | ! ! total volume of liquid seawater |
---|
88 | zvolssh = SUM( zarea_ssh(:,:) ) |
---|
89 | IF( lk_mpp ) CALL mpp_sum( zvolssh ) |
---|
90 | zvol = vol0 + zvolssh |
---|
91 | |
---|
92 | CALL iom_put( 'voltot', zvol ) |
---|
93 | CALL iom_put( 'sshtot', zvolssh / area_tot ) |
---|
94 | |
---|
95 | ! |
---|
96 | ztsn(:,:,:,jp_tem) = tsn(:,:,:,jp_tem) ! thermosteric ssh |
---|
97 | ztsn(:,:,:,jp_sal) = sn0(:,:,:) |
---|
98 | CALL eos( ztsn, zrhd, fsdept_n(:,:,:) ) ! now in situ density using initial salinity |
---|
99 | ! |
---|
100 | zbotpres(:,:) = 0._wp ! no atmospheric surface pressure, levitating sea-ice |
---|
101 | DO jk = 1, jpkm1 |
---|
102 | zbotpres(:,:) = zbotpres(:,:) + fse3t(:,:,jk) * zrhd(:,:,jk) |
---|
103 | END DO |
---|
104 | IF( .NOT.lk_vvl ) THEN |
---|
105 | IF ( ln_isfcav ) THEN |
---|
106 | DO ji=1,jpi |
---|
107 | DO jj=1,jpj |
---|
108 | zbotpres(ji,jj) = zbotpres(ji,jj) + sshn(ji,jj) * zrhd(ji,jj,mikt(ji,jj)) + riceload(ji,jj) |
---|
109 | END DO |
---|
110 | END DO |
---|
111 | ELSE |
---|
112 | zbotpres(:,:) = zbotpres(:,:) + sshn(:,:) * zrhd(:,:,1) |
---|
113 | END IF |
---|
114 | END IF |
---|
115 | ! |
---|
116 | zarho = SUM( area(:,:) * zbotpres(:,:) ) |
---|
117 | IF( lk_mpp ) CALL mpp_sum( zarho ) |
---|
118 | zssh_steric = - zarho / area_tot |
---|
119 | CALL iom_put( 'sshthster', zssh_steric ) |
---|
120 | |
---|
121 | ! ! steric sea surface height |
---|
122 | CALL eos( tsn, zrhd, zrhop, fsdept_n(:,:,:) ) ! now in situ and potential density |
---|
123 | zrhop(:,:,jpk) = 0._wp |
---|
124 | CALL iom_put( 'rhop', zrhop ) |
---|
125 | ! |
---|
126 | zbotpres(:,:) = 0._wp ! no atmospheric surface pressure, levitating sea-ice |
---|
127 | DO jk = 1, jpkm1 |
---|
128 | zbotpres(:,:) = zbotpres(:,:) + fse3t(:,:,jk) * zrhd(:,:,jk) |
---|
129 | END DO |
---|
130 | IF( .NOT.lk_vvl ) THEN |
---|
131 | IF ( ln_isfcav ) THEN |
---|
132 | DO ji=1,jpi |
---|
133 | DO jj=1,jpj |
---|
134 | zbotpres(ji,jj) = zbotpres(ji,jj) + sshn(ji,jj) * zrhd(ji,jj,mikt(ji,jj)) + riceload(ji,jj) |
---|
135 | END DO |
---|
136 | END DO |
---|
137 | ELSE |
---|
138 | zbotpres(:,:) = zbotpres(:,:) + sshn(:,:) * zrhd(:,:,1) |
---|
139 | END IF |
---|
140 | END IF |
---|
141 | ! |
---|
142 | zarho = SUM( area(:,:) * zbotpres(:,:) ) |
---|
143 | IF( lk_mpp ) CALL mpp_sum( zarho ) |
---|
144 | zssh_steric = - zarho / area_tot |
---|
145 | CALL iom_put( 'sshsteric', zssh_steric ) |
---|
146 | |
---|
147 | ! ! ocean bottom pressure |
---|
148 | zztmp = rau0 * grav * 1.e-4_wp ! recover pressure from pressure anomaly and cover to dbar = 1.e4 Pa |
---|
149 | zbotpres(:,:) = zztmp * ( zbotpres(:,:) + sshn(:,:) + thick0(:,:) ) |
---|
150 | CALL iom_put( 'botpres', zbotpres ) |
---|
151 | |
---|
152 | ! ! Mean density anomalie, temperature and salinity |
---|
153 | ztemp = 0._wp |
---|
154 | zsal = 0._wp |
---|
155 | DO jk = 1, jpkm1 |
---|
156 | DO jj = 1, jpj |
---|
157 | DO ji = 1, jpi |
---|
158 | zztmp = area(ji,jj) * fse3t(ji,jj,jk) |
---|
159 | ztemp = ztemp + zztmp * tsn(ji,jj,jk,jp_tem) |
---|
160 | zsal = zsal + zztmp * tsn(ji,jj,jk,jp_sal) |
---|
161 | END DO |
---|
162 | END DO |
---|
163 | END DO |
---|
164 | IF( .NOT.lk_vvl ) THEN |
---|
165 | IF ( ln_isfcav ) THEN |
---|
166 | DO ji=1,jpi |
---|
167 | DO jj=1,jpj |
---|
168 | ztemp = ztemp + zarea_ssh(ji,jj) * tsn(ji,jj,mikt(ji,jj),jp_tem) |
---|
169 | zsal = zsal + zarea_ssh(ji,jj) * tsn(ji,jj,mikt(ji,jj),jp_sal) |
---|
170 | END DO |
---|
171 | END DO |
---|
172 | ELSE |
---|
173 | ztemp = ztemp + SUM( zarea_ssh(:,:) * tsn(:,:,1,jp_tem) ) |
---|
174 | zsal = zsal + SUM( zarea_ssh(:,:) * tsn(:,:,1,jp_sal) ) |
---|
175 | END IF |
---|
176 | ENDIF |
---|
177 | IF( lk_mpp ) THEN |
---|
178 | CALL mpp_sum( ztemp ) |
---|
179 | CALL mpp_sum( zsal ) |
---|
180 | END IF |
---|
181 | ! |
---|
182 | zmass = rau0 * ( zarho + zvol ) ! total mass of liquid seawater |
---|
183 | ztemp = ztemp / zvol ! potential temperature in liquid seawater |
---|
184 | zsal = zsal / zvol ! Salinity of liquid seawater |
---|
185 | ! |
---|
186 | CALL iom_put( 'masstot', zmass ) |
---|
187 | CALL iom_put( 'temptot', ztemp ) |
---|
188 | CALL iom_put( 'saltot' , zsal ) |
---|
189 | ! |
---|
190 | CALL wrk_dealloc( jpi , jpj , zarea_ssh , zbotpres ) |
---|
191 | CALL wrk_dealloc( jpi , jpj , jpk , zrhd , zrhop ) |
---|
192 | CALL wrk_dealloc( jpi , jpj , jpk , jpts , ztsn ) |
---|
193 | ! |
---|
194 | IF( nn_timing == 1 ) CALL timing_stop('dia_ar5') |
---|
195 | ! |
---|
196 | END SUBROUTINE dia_ar5 |
---|
197 | |
---|
198 | |
---|
199 | SUBROUTINE dia_ar5_init |
---|
200 | !!---------------------------------------------------------------------- |
---|
201 | !! *** ROUTINE dia_ar5_init *** |
---|
202 | !! |
---|
203 | !! ** Purpose : initialization for AR5 diagnostic computation |
---|
204 | !!---------------------------------------------------------------------- |
---|
205 | INTEGER :: inum |
---|
206 | INTEGER :: ik |
---|
207 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
208 | REAL(wp) :: zztmp |
---|
209 | REAL(wp), POINTER, DIMENSION(:,:,:,:) :: zsaldta ! Jan/Dec levitus salinity |
---|
210 | !!---------------------------------------------------------------------- |
---|
211 | ! |
---|
212 | IF( nn_timing == 1 ) CALL timing_start('dia_ar5_init') |
---|
213 | ! |
---|
214 | CALL wrk_alloc( jpi , jpj , jpk, jpts, zsaldta ) |
---|
215 | ! ! allocate dia_ar5 arrays |
---|
216 | IF( dia_ar5_alloc() /= 0 ) CALL ctl_stop( 'STOP', 'dia_ar5_init : unable to allocate arrays' ) |
---|
217 | |
---|
218 | area(:,:) = e1t(:,:) * e2t(:,:) * tmask_i(:,:) |
---|
219 | |
---|
220 | area_tot = SUM( area(:,:) ) ; IF( lk_mpp ) CALL mpp_sum( area_tot ) |
---|
221 | |
---|
222 | vol0 = 0._wp |
---|
223 | thick0(:,:) = 0._wp |
---|
224 | DO jk = 1, jpkm1 |
---|
225 | vol0 = vol0 + SUM( area (:,:) * tmask(:,:,jk) * e3t_0(:,:,jk) ) |
---|
226 | thick0(:,:) = thick0(:,:) + tmask_i(:,:) * tmask(:,:,jk) * e3t_0(:,:,jk) |
---|
227 | END DO |
---|
228 | IF( lk_mpp ) CALL mpp_sum( vol0 ) |
---|
229 | |
---|
230 | CALL iom_open ( 'data_1m_salinity_nomask', inum ) |
---|
231 | CALL iom_get ( inum, jpdom_data, 'vosaline', zsaldta(:,:,:,1), 1 ) |
---|
232 | CALL iom_get ( inum, jpdom_data, 'vosaline', zsaldta(:,:,:,2), 12 ) |
---|
233 | CALL iom_close( inum ) |
---|
234 | sn0(:,:,:) = 0.5_wp * ( zsaldta(:,:,:,1) + zsaldta(:,:,:,2) ) |
---|
235 | sn0(:,:,:) = sn0(:,:,:) * tmask(:,:,:) |
---|
236 | IF( ln_zps ) THEN ! z-coord. partial steps |
---|
237 | DO jj = 1, jpj ! interpolation of salinity at the last ocean level (i.e. the partial step) |
---|
238 | DO ji = 1, jpi |
---|
239 | ik = mbkt(ji,jj) |
---|
240 | IF( ik > 1 ) THEN |
---|
241 | zztmp = ( gdept_1d(ik) - gdept_0(ji,jj,ik) ) / ( gdept_1d(ik) - gdept_1d(ik-1) ) |
---|
242 | sn0(ji,jj,ik) = ( 1._wp - zztmp ) * sn0(ji,jj,ik) + zztmp * sn0(ji,jj,ik-1) |
---|
243 | ENDIF |
---|
244 | END DO |
---|
245 | END DO |
---|
246 | ENDIF |
---|
247 | ! |
---|
248 | CALL wrk_dealloc( jpi , jpj , jpk, jpts, zsaldta ) |
---|
249 | ! |
---|
250 | IF( nn_timing == 1 ) CALL timing_stop('dia_ar5_init') |
---|
251 | ! |
---|
252 | END SUBROUTINE dia_ar5_init |
---|
253 | |
---|
254 | #else |
---|
255 | !!---------------------------------------------------------------------- |
---|
256 | !! Default option : NO diaar5 |
---|
257 | !!---------------------------------------------------------------------- |
---|
258 | LOGICAL, PUBLIC, PARAMETER :: lk_diaar5 = .FALSE. ! coupled flag |
---|
259 | CONTAINS |
---|
260 | SUBROUTINE dia_ar5_init ! Dummy routine |
---|
261 | END SUBROUTINE dia_ar5_init |
---|
262 | SUBROUTINE dia_ar5( kt ) ! Empty routine |
---|
263 | INTEGER :: kt |
---|
264 | WRITE(*,*) 'dia_ar5: You should not have seen this print! error?', kt |
---|
265 | END SUBROUTINE dia_ar5 |
---|
266 | #endif |
---|
267 | |
---|
268 | !!====================================================================== |
---|
269 | END MODULE diaar5 |
---|