1 |
% m-file: mit_loadglobal.m |
2 |
% plot various diagnostics of the global 4x4 degree runs |
3 |
|
4 |
% uses |
5 |
% mit_loadgrid |
6 |
% mit_getparm |
7 |
% mit_oceanmasks |
8 |
% mit_timesteps |
9 |
% mit_loadglobaldata |
10 |
% mit_globalmovie |
11 |
% mit_plotstreamfunctions |
12 |
% mit_plotzonalvelocity |
13 |
% mit_globalmean |
14 |
% mit_zonalmean |
15 |
% mit_meridflux |
16 |
% mit_plotmeandrift |
17 |
|
18 |
% $Header: /u/gcmpack/MITgcm/verification/tutorial_global_oce_latlon/diags_matlab/mit_loadglobal.m,v 1.4 2016/09/15 22:12:09 heimbach Exp $ |
19 |
% $Name: $ |
20 |
|
21 |
% read in all grid files, etc. This has to be done at the very beginning! |
22 |
grd = mit_loadgrid('.'); |
23 |
dname = grd.dname; |
24 |
% add some masks |
25 |
grd = mit_oceanmasks(grd); |
26 |
% |
27 |
if ~exist('meanfields','var') |
28 |
meanfields = 1; % read in mean fields instead of instantaneous ones, |
29 |
% include some more diagnostics |
30 |
% meanfields = 1 is the default |
31 |
end |
32 |
% |
33 |
mit_timesteps |
34 |
|
35 |
if meanfields |
36 |
disp('reading averaged fields (*tave.*)') |
37 |
% |
38 |
mit_loadglobaldata |
39 |
% |
40 |
u=rdmds('uVeltave',timesteps(2:nt)'); |
41 |
v=rdmds('vVeltave',timesteps(2:nt)'); |
42 |
w=rdmds('wVeltave',timesteps(2:nt)'); |
43 |
t=rdmds('Ttave',timesteps(2:nt)'); |
44 |
s=rdmds('Stave',timesteps(2:nt)'); |
45 |
eta=rdmds('ETAtave',timesteps(2:nt)').*repmat(squeeze(grd.cmask(:,:,1)),[1 1 nt-1]); |
46 |
k=1; |
47 |
u=cat(4,zeros([nx ny nz 1]),u); |
48 |
v=cat(4,zeros([nx ny nz 1]),v); |
49 |
w=cat(4,zeros([nx ny nz 1]),w); |
50 |
% hack |
51 |
t=cat(4,rdmds('T',timesteps(1)'),t); |
52 |
s=cat(4,rdmds('S',timesteps(1)'),s); |
53 |
eta=cat(3,0.*grd.cmask(:,:,1),eta); |
54 |
else |
55 |
% load snap shots |
56 |
disp('reading snap shot fields') |
57 |
u=rdmds('U',timesteps'); |
58 |
v=rdmds('V',timesteps'); |
59 |
w=rdmds('W',timesteps'); |
60 |
t=rdmds('T',timesteps'); |
61 |
s=rdmds('S',timesteps'); |
62 |
eta=rdmds('Eta',timesteps').*repmat(squeeze(grd.cmask(:,:,1)),[1 1 nt]); |
63 |
end |
64 |
% orientation and masking |
65 |
if ~strcmp(grd.buoyancy,'OCEANIC'); |
66 |
u = u(:,:,end:-1:1,:); |
67 |
v = v(:,:,end:-1:1,:); |
68 |
w = w(:,:,end:-1:1,:).*repmat(grd.cmask,[1 1 1 nt]); |
69 |
t = t(:,:,end:-1:1,:).*repmat(grd.cmask,[1 1 1 nt]); |
70 |
s = s(:,:,end:-1:1,:).*repmat(grd.cmask,[1 1 1 nt]); |
71 |
else |
72 |
% $$$ u = u; |
73 |
% $$$ v = v; |
74 |
w = w.*repmat(grd.cmask,[1 1 1 nt]); |
75 |
t = t.*repmat(grd.cmask,[1 1 1 nt]); |
76 |
s = s.*repmat(grd.cmask,[1 1 1 nt]); |
77 |
end |
78 |
|
79 |
% transport through Drake Passage |
80 |
TDP = NaN*ones(nt,1); |
81 |
kx = 73; |
82 |
kyg = [5 6 7]; |
83 |
%kyg = [4 5 6 7]; |
84 |
da = grd.dz*grd.dyg(kx,kyg); |
85 |
for k=kt; |
86 |
TDP(k) = sum(nansum(squeeze(u(kx,kyg,:,k))'.*da))*1e-6; % in Sv |
87 |
TDPz(:,k) = nansum(squeeze(u(kx,kyg,:,k)).*da')'*1e-6; % in Sv |
88 |
end |
89 |
|
90 |
% set some depth |
91 |
if ~meanfields |
92 |
mit_globalmovie |
93 |
end |
94 |
|
95 |
% $$$ [x y z]=meshgrid(grd.lonc,grd.grd.latc,-grd.zc); |
96 |
% $$$ tk = permute(squeeze(t(:,:,:,k)),[2 1 3]); |
97 |
% $$$ xplane = [0:30:360]; |
98 |
% $$$ yplane = 30; |
99 |
% $$$ zplane = -[50 500:500:6000]; |
100 |
% $$$ slice(x,y,z,tk,xplane,yplane,zplane); |
101 |
% $$$ shading flat, colorbar('h') |
102 |
|
103 |
mit_plotstreamfunctions |
104 |
|
105 |
mit_plotzonalvelocity |
106 |
|
107 |
if meanfields |
108 |
% plot some mean fields |
109 |
figure('PaperPosition',[0.31 0.25 10.5 7.88],'PaperOrientation','landscape') |
110 |
k=kmax; |
111 |
mySST = t(:,:,1,k); |
112 |
mySSS = s(:,:,1,k); |
113 |
myETA = eta(:,:,k); |
114 |
myVEL = sqrt(u(:,:,1,k).*u(:,:,1,k)+v(:,:,1,k).*v(:,:,1,k)); |
115 |
% |
116 |
clear sh; |
117 |
sh(1) = subplot(2,2,1);imagesc(grd.lonc,grd.latc,mySST'.*grd.hfacc(:,:,1)'); |
118 |
title(['Sea Surface Temperature [degC]']) |
119 |
sh(2) = subplot(2,2,2);imagesc(grd.lonc,grd.latc,mySSS'.*grd.hfacc(:,:,1)'); |
120 |
title(['Sea Surface Salinity [psu]']) |
121 |
sh(3) = subplot(2,2,3);imagesc(grd.lonc,grd.latc,myETA'); |
122 |
title(['Sea Surface Height [m]']) |
123 |
sh(4) = subplot(2,2,4);imagesc(grd.lonc,grd.latc,myVEL'.*grd.hfacc(:,:,1)'); |
124 |
title(['Sea Surface Speed [m/s]']) |
125 |
% |
126 |
axis(sh,'image'); axis(sh,'xy') |
127 |
set(gcf,'currentAxes',sh(1));colorbar('h') |
128 |
set(gcf,'currentAxes',sh(2));colorbar('h') |
129 |
set(gcf,'currentAxes',sh(3));colorbar('h') |
130 |
set(gcf,'currentAxes',sh(4));colorbar('h') |
131 |
end |
132 |
|
133 |
if meanfields |
134 |
% compute actual heat and E-P fluxes from time averages |
135 |
% this only makes sense when you load the time averaged fields |
136 |
tauTheta = mit_getparm('data','tauThetaClimRelax'); |
137 |
tauSalt = mit_getparm('data','tauSaltClimRelax'); |
138 |
rhonil = mit_getparm('data','rhonil'); |
139 |
if isempty(rhonil) |
140 |
rhonil = 1035; |
141 |
end |
142 |
Cp = 3994; |
143 |
if isempty(tauTheta) |
144 |
tauTheta = 0; |
145 |
end |
146 |
if tauTheta==0 |
147 |
recip_tauTheta = 0.; |
148 |
else |
149 |
recip_tauTheta = 1./tauTheta; |
150 |
end |
151 |
if isempty(tauSalt) |
152 |
tauSalt = 0; |
153 |
end |
154 |
if tauSalt==0 |
155 |
recip_tauSalt = 0.; |
156 |
else |
157 |
recip_tauSalt = 1./tauSalt; |
158 |
end |
159 |
figure('PaperPosition',[0.31 0.25 10.5 7.88],'PaperOrientation','landscape') |
160 |
k=kmax; |
161 |
qnet_diag = (t(:,:,1,k)-tdmean)*grd.dz(1).*grd.hfacc(:,:,1)*Cp*rhonil*recip_tauTheta; |
162 |
empr_diag = - (s(:,:,1,k)-sdmean)*grd.dz(1).*grd.hfacc(:,:,1)*recip_tauSalt/35; |
163 |
clear sh; |
164 |
sh(1) = subplot(2,2,1);imagesc(grd.lonc,grd.latc,qnet_diag'); |
165 |
title(['Q_{diag} from T-T_{lev} at the surface']) |
166 |
sh(2) =subplot(2,2,2);imagesc(grd.lonc,grd.latc,qnet'.*grd.cmask(:,:,1)'); |
167 |
title('annual mean Q_{net} (data)') |
168 |
sh(3) =subplot(2,2,3);imagesc(grd.lonc,grd.latc,empr_diag'); |
169 |
title(['(E-P)_{diag} from S-S_{lev} at the surface']) |
170 |
sh(4) =subplot(2,2,4);imagesc(grd.lonc,grd.latc,empr'.*grd.cmask(:,:,1)'); |
171 |
title('annual mean E-P (data)') |
172 |
axis(sh,'image'); axis(sh,'xy') |
173 |
suptitle(['experiment ' dname ', timestep = ' num2str(timesteps(k)) ... |
174 |
', ' tuname ' = ' num2str(tim(k))]) |
175 |
set(sh(1:2),'clim',[-1 1]*200); |
176 |
set(gcf,'currentAxes',sh(1));colorbar('h') |
177 |
set(gcf,'currentAxes',sh(2));colorbar('h') |
178 |
set(sh(3:4),'clim',[-1 1]*1.e-7); |
179 |
set(gcf,'currentAxes',sh(3));colorbar('h') |
180 |
set(gcf,'currentAxes',sh(4));colorbar('h') |
181 |
|
182 |
clear sh |
183 |
|
184 |
if length(tim) > 1 |
185 |
rac3d = repmat(grd.rac,[1 1 nz]).*grd.hfacc; |
186 |
% horizontal averages of t and s in a hovmoeller-type diagramme |
187 |
mdt = repmat(NaN,grd.nz,length(kt)); |
188 |
mdtt= mdt; |
189 |
mdtstd= mdt; |
190 |
mdt_alt = mdt; |
191 |
mdtt_alt= mdt; |
192 |
mdtstd_alt = mdt; |
193 |
mdt_pac = mdt; |
194 |
mdtt_pac= mdt; |
195 |
mdtstd_pac = mdt; |
196 |
mdt_sou = mdt; |
197 |
mdtt_sou= mdt; |
198 |
mdtstd_sou = mdt; |
199 |
mds = mdt; |
200 |
mdss= mdt; |
201 |
mdsstd = mdt; |
202 |
mds_alt = mdt; |
203 |
mdss_alt= mdt; |
204 |
mdsstd_alt = mdt; |
205 |
mds_pac = mdt; |
206 |
mdss_pac= mdt; |
207 |
mdsstd_pac = mdt; |
208 |
mds_sou = mdt; |
209 |
mdss_sou= mdt; |
210 |
mdsstd_sou = mdt; |
211 |
indopac_hfacc = ... |
212 |
change(change(grd.pacific_hfacc,'==',NaN,0) ... |
213 |
+change(grd.indic_hfacc,'==',NaN,0),'==',0,NaN); |
214 |
indopac_hfacc(:,1:12,:) = NaN; |
215 |
southern_hfacc = grd.hfacc; |
216 |
southern_hfacc(:,13:end,:) = NaN; |
217 |
atl_hfacc = grd.atlantic_hfacc; |
218 |
atl_hfacc(:,1:12,:) = NaN; |
219 |
for k=kt; |
220 |
dt = t(:,:,:,k)-tdatamean; % anomalies |
221 |
ds = s(:,:,:,k)-sdatamean; % anomalies |
222 |
[mdt(:,k) mdtt(:,k) mdtstd(:,k)] = mit_globalmean(dt,rac3d); |
223 |
[mds(:,k) mdss(:,k) mdsstd(:,k)] = mit_globalmean(ds,rac3d); |
224 |
[mdt_atl(:,k) mdtt_atl(:,k) mdtstd_atl(:,k)] = mit_globalmean(dt,rac3d.*atl_hfacc); |
225 |
[mds_atl(:,k) mdss_atl(:,k) mdsstd_atl(:,k)] = mit_globalmean(ds,rac3d.*atl_hfacc); |
226 |
[mdt_pac(:,k) mdtt_pac(:,k) mdtstd_pac(:,k)] = mit_globalmean(dt,rac3d.*indopac_hfacc); |
227 |
[mds_pac(:,k) mdss_pac(:,k) mdsstd_pac(:,k)] = mit_globalmean(ds,rac3d.*indopac_hfacc); |
228 |
[mdt_sou(:,k) mdtt_sou(:,k) mdtstd_sou(:,k)] = mit_globalmean(dt,rac3d.*southern_hfacc); |
229 |
[mds_sou(:,k) mdss_sou(:,k) mdsstd_sou(:,k)] = mit_globalmean(ds,rac3d.*southern_hfacc); |
230 |
end |
231 |
|
232 |
figure('PaperPosition',[0.31 0.25 10.5 7.88],'PaperOrientation','landscape') |
233 |
if length(tim) > 2 |
234 |
clear sh |
235 |
sh(1) = subplot(2,2,1); |
236 |
contourf(tim(2:end),-grd.zc,mdt(:,2:end),20); |
237 |
caxis([-1 1]*max(abs(caxis))); shading flat; colorbar('v') |
238 |
title('T-T_{lev} horizontally averaged') |
239 |
xlabel(['Time [' timeunit ']']); ylabel('Depth [m]') |
240 |
sh(2) = subplot(2,2,2); |
241 |
contourf(tim(2:end),-grd.zc,mds(:,2:end),20); |
242 |
caxis([-1 1]*max(abs(caxis))); shading flat; colorbar('v') |
243 |
title('S-S_{lev} horizontally averaged') |
244 |
xlabel(['Time [' timeunit ']']); ylabel('Depth [m]') |
245 |
% |
246 |
sh(3) = subplot(2,2,3); |
247 |
contourf(tim(2:end),-grd.zc,mdtt(:,2:end),20); |
248 |
shading flat; colorbar('v') |
249 |
title('|T-T_{lev}| horizontally averaged'); |
250 |
xlabel(['Time [' timeunit ']']); ylabel('Depth [m]') |
251 |
sh(4) = subplot(2,2,4); |
252 |
contourf(tim(2:end),-grd.zc,mdss(:,2:end),20); |
253 |
shading flat; colorbar('v') |
254 |
title('|S-S_{lev}| horizontally averaged') |
255 |
xlabel(['Time [' timeunit ']']); ylabel('Depth [m]') |
256 |
set(sh,'layer','top'); |
257 |
suptitle(['experiment ' dname]) |
258 |
end %if length(tim) > 2 |
259 |
|
260 |
k=kmax; |
261 |
figure('PaperOrientation','landscape','PaperPosition',[0.25 0.3125 10.5 7.875]) |
262 |
clear sh |
263 |
sh(1) = subplot(2,4,1); |
264 |
plot(mdt(:,k),-grd.zc,'b-',mdt(:,k)+mdtstd(:,k),-grd.zc,'r--',mdt(:,k)-mdtstd(:,k),-grd.zc,'r--'); |
265 |
xlabel('T-T_{lev} [degC]'); ylabel('depth [m]'); title('global') |
266 |
legend('horiz. mean','std',3); |
267 |
sh(2) = subplot(2,4,2); |
268 |
plot(mdt_atl(:,k),-grd.zc,'b-',mdt_atl(:,k)+mdtstd_atl(:,k),-grd.zc,'r--',mdt_atl(:,k)-mdtstd_atl(:,k),-grd.zc,'r--'); |
269 |
xlabel('T-T_{lev} [degC]'); %ylabel('depth [m]'); |
270 |
title('atlantic') |
271 |
sh(3) = subplot(2,4,3); |
272 |
plot(mdt_pac(:,k),-grd.zc,'b-',mdt_pac(:,k)+mdtstd_pac(:,k),-grd.zc,'r--',mdt_pac(:,k)-mdtstd_pac(:,k),-grd.zc,'r--'); |
273 |
xlabel('T-T_{lev} [degC]'); %ylabel('depth [m]'); |
274 |
title('indo-pacific') |
275 |
sh(4) = subplot(2,4,4); |
276 |
plot(mdt_sou(:,k),-grd.zc,'b-',mdt_sou(:,k)+mdtstd_sou(:,k),-grd.zc,'r--',mdt_sou(:,k)-mdtstd_sou(:,k),-grd.zc,'r--'); |
277 |
xlabel('T-T_{lev} [degC]'); ylabel('depth [m]'); title('southern') |
278 |
% |
279 |
sh(5) = subplot(2,4,5); |
280 |
plot(mds(:,k),-grd.zc,'b-',mds(:,k)+mdsstd(:,k),-grd.zc,'r--',mds(:,k)-mdsstd(:,k),-grd.zc,'r--'); |
281 |
xlabel('S-S_{lev} [PSU]'); ylabel('depth [m]'); |
282 |
sh(6) = subplot(2,4,6); |
283 |
plot(mds_atl(:,k),-grd.zc,'b-',mds_atl(:,k)+mdsstd_atl(:,k),-grd.zc,'r--',mds_atl(:,k)-mdsstd_atl(:,k),-grd.zc,'r--'); |
284 |
xlabel('S-S_{lev} [PSU]'); %ylabel('depth [m]'); |
285 |
sh(7) = subplot(2,4,7); |
286 |
plot(mds_pac(:,k),-grd.zc,'b-',mds_pac(:,k)+mdsstd_pac(:,k),-grd.zc,'r--',mds_pac(:,k)-mdsstd_pac(:,k),-grd.zc,'r--'); |
287 |
xlabel('S-S_{lev} [PSU]'); %ylabel('depth [m]'); |
288 |
sh(8) = subplot(2,4,8); |
289 |
plot(mds_sou(:,k),-grd.zc,'b-',mds_sou(:,k)+mdsstd_sou(:,k),-grd.zc,'r--',mds_sou(:,k)-mdsstd_sou(:,k),-grd.zc,'r--'); |
290 |
xlabel('S-S_{lev} [PSU]'); ylabel('depth [m]'); |
291 |
set(sh,'xgrid','on','ygrid','on') |
292 |
set([sh(4); sh(8)],'YAxisLocation','right') |
293 |
suptitle(['experiment ' dname ', timestep = ' num2str(timesteps(k)) ... |
294 |
', ' tuname ' = ' num2str(tim(k))]) |
295 |
end %if length(tim) > 1 |
296 |
|
297 |
% zonal mean of temperature and salinity for different ocean basins |
298 |
k=kmax; |
299 |
tdzm_glo = mit_zonalmean(tdatamean,grd.hfacc,grd.dxc); |
300 |
tdzm_atl = mit_zonalmean(tdatamean,grd.atlantic_hfacc,grd.dxc); |
301 |
tdzm_pac = mit_zonalmean(tdatamean,grd.pacific_hfacc,grd.dxc); |
302 |
tdzm_ind = mit_zonalmean(tdatamean,grd.indic_hfacc,grd.dxc); |
303 |
sdzm_glo = mit_zonalmean(sdatamean,grd.hfacc,grd.dxc); |
304 |
sdzm_atl = mit_zonalmean(sdatamean,grd.atlantic_hfacc,grd.dxc); |
305 |
sdzm_pac = mit_zonalmean(sdatamean,grd.pacific_hfacc,grd.dxc); |
306 |
sdzm_ind = mit_zonalmean(sdatamean,grd.indic_hfacc,grd.dxc); |
307 |
|
308 |
figure('PaperPosition',[0.25 0.368552 8 10.2629]); |
309 |
tlev = -5:.5:5; |
310 |
slev = -1:.1:1; |
311 |
clear sh clh |
312 |
clh = cell(8,1); |
313 |
delay = 1; |
314 |
for k = kmax |
315 |
tzm_glo = mit_zonalmean(t(:,:,:,k),grd.hfacc,grd.dxc); |
316 |
tzm_atl = mit_zonalmean(t(:,:,:,k),grd.atlantic_hfacc,grd.dxc); |
317 |
tzm_pac = mit_zonalmean(t(:,:,:,k),grd.pacific_hfacc,grd.dxc); |
318 |
tzm_ind = mit_zonalmean(t(:,:,:,k),grd.indic_hfacc,grd.dxc); |
319 |
szm_glo = mit_zonalmean(s(:,:,:,k),grd.hfacc,grd.dxc); |
320 |
szm_atl = mit_zonalmean(s(:,:,:,k),grd.atlantic_hfacc,grd.dxc); |
321 |
szm_pac = mit_zonalmean(s(:,:,:,k),grd.pacific_hfacc,grd.dxc); |
322 |
szm_ind = mit_zonalmean(s(:,:,:,k),grd.indic_hfacc,grd.dxc); |
323 |
caxt = [min(tzm_glo(:)-tdzm_glo(:)) max(tzm_glo(:)-tdzm_glo(:))]; |
324 |
tlev = -5:.5:5; |
325 |
if max(abs(caxt)) < 1; |
326 |
tlev = -1:.1:1; |
327 |
end |
328 |
if max(abs(caxt)) < .1; |
329 |
tlev = .1*tlev; |
330 |
end |
331 |
caxs = [min(szm_glo(:)-sdzm_glo(:)) max(szm_glo(:)-sdzm_glo(:))]; |
332 |
slev = -1.:.05:1; |
333 |
if max(abs(caxs)) < 0.2; |
334 |
slev = -.20:.01:.20; |
335 |
end |
336 |
if max(abs(caxs)) < .02; |
337 |
slev = .1*slev; |
338 |
end |
339 |
sh(1) = subplot(4,2,1); |
340 |
[cs h] = contourf(grd.latc,-grd.zc,(tzm_glo-tdzm_glo)',tlev); |
341 |
if ~isempty(h); |
342 |
clh{1} = clabel(cs,h,[0 0]); |
343 |
end |
344 |
title('\theta-\theta_{lev} [degC]: global ocean') |
345 |
sh(3) = subplot(4,2,3); |
346 |
[cs h] = contourf(grd.latc,-grd.zc,(tzm_atl-tdzm_atl)',tlev); |
347 |
if ~isempty(h); |
348 |
clh{3} = clabel(cs,h,[0 0]); |
349 |
end |
350 |
title('atlantic ocean') |
351 |
sh(5) = subplot(4,2,5); |
352 |
[cs h] = contourf(grd.latc,-grd.zc,(tzm_pac-tdzm_pac)',tlev); |
353 |
if ~isempty(h); |
354 |
clh{5} = clabel(cs,h,[0 0]); |
355 |
end |
356 |
title('pacific ocean') |
357 |
sh(7) = subplot(4,2,7); |
358 |
[cs h] = contourf(grd.latc,-grd.zc,(tzm_ind-tdzm_ind)',tlev); |
359 |
if ~isempty(h); |
360 |
clh{7} = clabel(cs,h,[0 0]); |
361 |
end |
362 |
title('indian ocean') |
363 |
set(sh(1:2:end),'clim',[tlev(1) tlev(end)]) |
364 |
colorbar('h') |
365 |
sh(2) = subplot(4,2,2); |
366 |
[cs h] = contourf(grd.latc,-grd.zc,(szm_glo-sdzm_glo)',slev); |
367 |
if ~isempty(h); |
368 |
clh{2} = clabel(cs,h,[0 0]); |
369 |
end |
370 |
title('S-S_{lev} [PSU]: global ocean') |
371 |
sh(4) = subplot(4,2,4); |
372 |
[cs h] = contourf(grd.latc,-grd.zc,(szm_atl-sdzm_atl)',slev); |
373 |
if ~isempty(h); |
374 |
clh{4} = clabel(cs,h,[0 0]); |
375 |
end |
376 |
title('atlantic ocean') |
377 |
sh(6) = subplot(4,2,6); |
378 |
[cs h] = contourf(grd.latc,-grd.zc,(szm_pac-sdzm_pac)',slev); |
379 |
if ~isempty(h); |
380 |
clh{6} = clabel(cs,h,[0 0]); |
381 |
end |
382 |
title('pacific ocean') |
383 |
sh(8) = subplot(4,2,8); |
384 |
[cs h] = contourf(grd.latc,-grd.zc,(szm_ind-sdzm_ind)',slev); |
385 |
if ~isempty(h); |
386 |
clh{8} = clabel(cs,h,[0 0]); |
387 |
end |
388 |
title('indian ocean') |
389 |
set(sh(2:2:end),'clim',[slev(1) slev(end)]) |
390 |
colorbar('h') |
391 |
set(sh,'layer','top') |
392 |
if ~isempty(cat(1,clh{:})) |
393 |
set(cat(1,clh{:}),'fontsize',8); |
394 |
end |
395 |
|
396 |
suptitle(['experiment ' dname ', timestep = ' num2str(timesteps(k)) ... |
397 |
', ' tuname ' = ' num2str(tim(k)) ', zonal averages']) |
398 |
drawnow; |
399 |
pause(delay) |
400 |
end |
401 |
clear clh sh |
402 |
|
403 |
% meridional transport of heat (internal energy) and fresh water |
404 |
% diagnosed from surface fluxes |
405 |
petawatts = 1e-15; |
406 |
sverdrups = 1e-6; |
407 |
heat_data = -mit_meridflux(qnet,grd.dxc,grd.dyc)*petawatts; |
408 |
heat_diag = -mit_meridflux(qnet_diag,grd.dxc,grd.dyc)*petawatts; |
409 |
freshwater_data = -mit_meridflux(empr,grd.dxc,grd.dyc)*sverdrups; |
410 |
freshwater_diag = -mit_meridflux(empr_diag,grd.dxc,grd.dyc)*sverdrups; |
411 |
figure |
412 |
latgf = [2*grd.latc-grd.latg]; |
413 |
clear sh |
414 |
sh(1) = subplot(2,1,1); |
415 |
hh = plot(latgf,heat_data,'-',latgf,heat_diag,'--', ... |
416 |
latgf,heat_data+heat_diag,'-.'); |
417 |
title('heat (internal energy) flux [PW]') |
418 |
suptitle(['experiment ' dname ', timestep = ' num2str(timesteps(k)) ... |
419 |
', ' tuname ' = ' num2str(tim(k))]) |
420 |
legend('surface flux','restoring','net',2) |
421 |
sh(2) = subplot(2,1,2); |
422 |
fwh = plot(latgf,freshwater_data,'-',latgf,freshwater_diag,'--', ... |
423 |
latgf,freshwater_data+freshwater_diag,'-.'); |
424 |
title('freshwater flux [Sv]') |
425 |
legend('surface flux','restoring','net',3) |
426 |
set(sh,'xlim',[latgf(1) latgf(end)],'xgrid','on','ygrid','on'); |
427 |
|
428 |
end % meanfields |
429 |
|
430 |
mit_plotmeandrift |
431 |
|
432 |
fname = [grd.dname '_' sprintf('%u',tim(kmax)) timeunit]; |
433 |
in = findstr(fname,'\'); |
434 |
fname(in) = []; |
435 |
fn = gcf; |
436 |
if meanfields |
437 |
fname = [fname '.ps']; |
438 |
% change required by newer Matlab versions |
439 |
% f0 = fn-7; |
440 |
f0=fn.Number-7; |
441 |
else |
442 |
fname = [fname '_snap.ps']; |
443 |
% change required by newer Matlab versions |
444 |
% f0 = fn-4; |
445 |
f0=fn.Number-4; |
446 |
end |
447 |
print(f0,'-dpsc',fname); |
448 |
for k = f0+1:fn |
449 |
print(k,'-dpsc',fname,'-append') |
450 |
end |