/[MITgcm]/MITgcm_contrib/gael/matlab_class/gcmfaces_diags/diags_set_A.m
ViewVC logotype

Annotation of /MITgcm_contrib/gael/matlab_class/gcmfaces_diags/diags_set_A.m

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.11 - (hide annotations) (download)
Sun Aug 3 19:51:49 2014 UTC (10 years, 11 months ago) by gforget
Branch: MAIN
Changes since 1.10: +1 -1 lines
- use trsp_2d_set1 if available, trsp_3d_set2 otherwise.

1 gforget 1.1
2     if userStep==1;%diags to be computed
3     listDiags='fldBAR fldTRANSPORTS gloOV gloOVbolus gloOVres gloMT_H gloMT_FW gloMT_SLT';
4 gforget 1.3 if sum([90 1170]~=mygrid.ioSize)==0;
5     listDiags=[listDiags ' atlOV atlOVbolus atlOVres atlMT_H atlMT_FW atlMT_SLT'];
6     listDiags=[listDiags ' pacindOV pacindOVbolus pacindOVres pacindMT_H pacindMT_FW pacindMT_SLT'];
7     end;
8 gforget 1.1 elseif userStep==2;%input files and variables
9     listFlds={'UVELMASS','VVELMASS','GM_PsiX','GM_PsiY'};
10     listFlds={listFlds{:},'ADVx_TH ','ADVy_TH ','DFxE_TH ','DFyE_TH '};
11     listFlds={listFlds{:},'ADVx_SLT','ADVy_SLT','DFxE_SLT','DFyE_SLT'};
12     listFldsNames=deblank(listFlds);
13 gforget 1.11 listFiles={'trsp_3d_set1','trsp_2d_set1','trsp_3d_set2'};
14 gforget 1.9 listSubdirs={[dirModel 'diags/TRSP/']};
15 gforget 1.1 elseif userStep==3;%computational part;
16     %mask fields:
17     fldU=UVELMASS.*mygrid.mskW; fldV=VVELMASS.*mygrid.mskS;
18     if ~isempty(whos('ADVx_TH'));%for backward compatibility
19     if size(ADVx_TH{1},3)>1; %assume full 3D fields
20     mskW=mygrid.mskW; mskS=mygrid.mskS;
21     else; %assume vertically integrated (2D fields)
22     mskW=mygrid.mskW(:,:,1); mskS=mygrid.mskS(:,:,1);
23     end;
24     ADVx_TH=ADVx_TH.*mskW; ADVy_TH=ADVy_TH.*mskS;
25     ADVx_SLT=ADVx_SLT.*mskW; ADVy_SLT=ADVy_SLT.*mskS;
26     DFxE_TH=DFxE_TH.*mskW; DFyE_TH=DFyE_TH.*mskS;
27     DFxE_SLT=DFxE_SLT.*mskW; DFyE_SLT=DFyE_SLT.*mskS;
28     end;
29     if ~isempty(whos('GM_PsiX'));%for backward compatibility
30     [fldUbolus,fldVbolus,fldWbolus]=calc_bolus(GM_PsiX,GM_PsiY);
31     fldUbolus=fldUbolus.*mygrid.mskW; fldVbolus=fldVbolus.*mygrid.mskS;
32     fldUres=fldU+fldUbolus; fldVres=fldV+fldVbolus;
33     end;
34    
35     %compute barotropic stream function:
36     [fldBAR]=calc_barostream(fldU,fldV);
37     %compute transports along transects:
38 gforget 1.2 [fldTRANSPORTS]=1e-6*calc_transports(fldU,fldV,mygrid.LINES_MASKS,{'dh','dz'});
39 gforget 1.3
40 gforget 1.1 %compute overturning stream functions:
41 gforget 1.3 if sum([90 1170]~=mygrid.ioSize)==0;
42     listBasins=[1:3];
43     else;
44     listBasins=1;
45     end;
46     for bb=listBasins;
47 gforget 1.1 %mask : global, atlantic or Pac+Ind
48     if bb==1; mskC=mygrid.mskC; mskC(:)=1;
49     elseif bb==2; mskC=v4_basin({'atlExt'}); mskC=mk3D(mskC,fldU);
50     elseif bb==3; mskC=v4_basin({'pacExt','indExt'}); mskC=mk3D(mskC,fldU);
51     end;
52     %note: while mskC is a basin mask for tracer points, it can be applied to U/V below
53     %compute overturning: eulerian contribution
54     [fldOV]=calc_overturn(fldU.*mskC,fldV.*mskC);
55     if ~isempty(whos('GM_PsiX'));%for backward compatibility
56     %compute overturning: eddy contribution
57     [fldOVbolus]=calc_overturn(fldUbolus.*mskC,fldVbolus.*mskC);
58     %compute overturning: residual overturn
59     [fldOVres]=calc_overturn(fldUres.*mskC,fldVres.*mskC);
60     else;
61     fldOVbolus=NaN*fldOV; fldOVres=NaN*fldOV;
62     end;
63    
64     if ~isempty(whos('ADVx_TH'));%for backward compatibility
65     %compute meridional heat transports:
66     tmpU=(ADVx_TH+DFxE_TH).*mskC(:,:,1:size(ADVx_TH{1},3));
67     tmpV=(ADVy_TH+DFyE_TH).*mskC(:,:,1:size(ADVx_TH{1},3));
68     [fldMT_H]=1e-15*4e6*calc_MeridionalTransport(tmpU,tmpV,0);
69     %compute meridional fresh water transports:
70     %... using the virtual salt flux formula:
71     %[fldMT_FW]=1e-6/35*calc_MeridionalTransport(ADVx_SLT+DFxE_SLT,ADVy_SLT+DFyE_SLT,0);
72     %[fldMT_FW]=1e-6/35*calc_MeridionalTransport(ADVx_SLT,ADVy_SLT,0);
73     %... using the real freshwater flux formula:
74     [fldMT_FW]=1e-6*calc_MeridionalTransport(fldU.*mskC,fldV.*mskC,1);
75     %compute meridional salt transports:
76     tmpU=(ADVx_SLT+DFxE_SLT).*mskC(:,:,1:size(ADVx_TH{1},3));
77     tmpV=(ADVy_SLT+DFyE_SLT).*mskC(:,:,1:size(ADVx_TH{1},3));
78     [fldMT_SLT]=1e-6*calc_MeridionalTransport(tmpU,tmpV,0);
79     else;
80     fldMT_H=NaN*mygrid.LATS; fldMT_FW=NaN*mygrid.LATS; fldMT_SLT=NaN*mygrid.LATS;
81     end;
82    
83     %store to global, atlantic or Pac+Ind arrays:
84     if bb==1;
85     gloMT_H=fldMT_H; gloMT_FW=fldMT_FW; gloMT_SLT=fldMT_SLT;
86     gloOV=fldOV; gloOVbolus=fldOVbolus; gloOVres=fldOVres;
87     elseif bb==2;
88 gforget 1.5 kk=find(mygrid.LATS<-35|mygrid.LATS>70);
89 gforget 1.1 fldMT_H(kk,:)=NaN; fldMT_FW(kk,:)=NaN; fldMT_SLT(kk,:)=NaN;
90     fldOV(kk,:)=NaN; fldOVbolus(kk,:)=NaN; fldOVres(kk,:)=NaN;
91     atlMT_H=fldMT_H; atlMT_FW=fldMT_FW; atlMT_SLT=fldMT_SLT;
92     atlOV=fldOV; atlOVbolus=fldOVbolus; atlOVres=fldOVres;
93     elseif bb==3;
94 gforget 1.5 kk=find(mygrid.LATS<-35|mygrid.LATS>65);
95 gforget 1.1 fldMT_H(kk,:)=NaN; fldMT_FW(kk,:)=NaN; fldMT_SLT(kk,:)=NaN;
96     fldOV(kk,:)=NaN; fldOVbolus(kk,:)=NaN; fldOVres(kk,:)=NaN;
97     pacindMT_H=fldMT_H; pacindMT_FW=fldMT_FW; pacindMT_SLT=fldMT_SLT;
98     pacindOV=fldOV; pacindOVbolus=fldOVbolus; pacindOVres=fldOVres;
99     end;
100     end;
101 gforget 1.4
102     %===================== COMPUTATIONAL SEQUENCE ENDS =========================%
103     %===================== PLOTTING SEQUENCE BEGINS =========================%
104    
105     elseif userStep==-1;%plotting
106    
107     if isempty(setDiagsParams);
108     choicePlot={'all'};
109     else;
110     choicePlot=setDiagsParams;
111     end;
112    
113     if sum(strcmp(choicePlot,'all'))|sum(strcmp(choicePlot,'bf'));
114    
115 gforget 1.6 if addToTex; write2tex(fileTex,1,'barotropic streamfunction',2); end;
116    
117 gforget 1.4 %barotropic streamfunction:
118     fld=mean(alldiag.fldBAR(:,:,tt),3);
119     cc=[[-80:20:-40] [-25 -15:5:15 25] [40:40:200]]; title0='Horizontal Stream Function';
120     if doAnomalies; cc=scaleAnom*[-5:0.5:5]; end;
121     figureL; m_map_gcmfaces(fld,0,{'myCaxis',cc},{'do_m_coast',1},{'myTitle',title0});
122     myCaption={myYmeanTxt,'mean -- barotropic streamfunction (Sv)'};
123     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
124     %and the corresponding standard deviation:
125     if multiTimes;
126     fld=std(alldiag.fldBAR(:,:,tt),[],3);
127     cc=[0:0.5:3 4 5 7 10 15:5:25 35 50]; title0='Horizontal Stream Function';
128     if doAnomalies; cc=scaleAnom*[0:0.25:2]; end;
129     figureL; m_map_gcmfaces(fld,0,{'myCaxis',cc},{'do_m_coast',1},{'myTitle',title0});
130     myCaption={myYmeanTxt,' standard deviation -- barotropic streamfunction (Sv)'};
131     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
132     end;
133    
134     end;
135    
136     if sum(strcmp(choicePlot,'all'))|sum(strcmp(choicePlot,'ov'));
137    
138 gforget 1.6 if addToTex; write2tex(fileTex,1,'meridional streamfunction',2); end;
139    
140 gforget 1.4 %meridional streamfunction (Eulerian) :
141     fld=mean(alldiag.gloOV(:,:,tt),3); fld(fld==0)=NaN;
142     X=mygrid.LATS*ones(1,length(mygrid.RF)); Y=ones(length(mygrid.LATS),1)*(mygrid.RF');
143     cc=[[-50:10:-30] [-24:3:24] [30:10:50]]; title0='Meridional Stream Function';
144     if doAnomalies; cc=scaleAnom*[-1:0.1:1]; end;
145     figureL; set(gcf,'Renderer','zbuffer'); %set(gcf,'Units','Normalized','Position',[0.05 0.1 0.4 0.8]);
146     depthStretchPlot('pcolor',{X,Y,fld}); shading interp; cbar=gcmfaces_cmap_cbar(cc); title(title0);
147     myCaption={myYmeanTxt,'mean -- overturning streamfunction (Sv)'};
148     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
149    
150     %meridional streamfunction (residual) :
151     fld=mean(alldiag.gloOVres(:,:,tt),3); fld(fld==0)=NaN;
152     X=mygrid.LATS*ones(1,length(mygrid.RF)); Y=ones(length(mygrid.LATS),1)*(mygrid.RF');
153     cc=[[-50:10:-30] [-24:3:24] [30:10:50]]; title0='Meridional Stream Function (incl. GM)';
154     if doAnomalies; cc=scaleAnom*[-1:0.1:1]; end;
155     figureL; set(gcf,'Renderer','zbuffer'); %set(gcf,'Units','Normalized','Position',[0.05 0.1 0.4 0.8]);
156     depthStretchPlot('pcolor',{X,Y,fld}); shading interp; cbar=gcmfaces_cmap_cbar(cc); title(title0);
157     myCaption={myYmeanTxt,'mean -- overturning streamfunction incl. GM (Sv)'};
158     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
159    
160     if multiBasins;
161     %meridional streamfunction (Eulerian):
162     fld=mean(alldiag.atlOV(:,:,tt),3); fld(fld==0)=NaN;
163     X=mygrid.LATS*ones(1,length(mygrid.RF)); Y=ones(length(mygrid.LATS),1)*(mygrid.RF');
164     cc=[[-50:10:-30] [-24:3:24] [30:10:50]]; title0='Atlantic Meridional Stream Function';
165     if doAnomalies; cc=scaleAnom*[-1:0.1:1]; end;
166     figureL; set(gcf,'Renderer','zbuffer'); %set(gcf,'Units','Normalized','Position',[0.05 0.1 0.4 0.8]);
167     depthStretchPlot('pcolor',{X,Y,fld}); shading interp; cbar=gcmfaces_cmap_cbar(cc); title(title0);
168     myCaption={myYmeanTxt,'mean -- Atlantic overturning streamfunction (Sv)'};
169     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
170     %meridional streamfunction (residual):
171     fld=mean(alldiag.pacindOV(:,:,tt),3); fld(fld==0)=NaN;
172     X=mygrid.LATS*ones(1,length(mygrid.RF)); Y=ones(length(mygrid.LATS),1)*(mygrid.RF');
173     cc=[[-50:10:-30] [-24:3:24] [30:10:50]]; title0='Pac+Ind Meridional Stream Function';
174     if doAnomalies; cc=scaleAnom*[-1:0.1:1]; end;
175     figureL; set(gcf,'Renderer','zbuffer'); %set(gcf,'Units','Normalized','Position',[0.05 0.1 0.4 0.8]);
176     depthStretchPlot('pcolor',{X,Y,fld}); shading interp; cbar=gcmfaces_cmap_cbar(cc); title(title0);
177     myCaption={myYmeanTxt,'mean -- Pac+Ind overturning streamfunction (Sv)'};
178     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
179     end;
180    
181     if multiTimes;
182     %and the corresponding standard deviation:
183     fld=std(alldiag.gloOV(:,:,tt),[],3); fld(fld==0)=NaN;
184     X=mygrid.LATS*ones(1,length(mygrid.RF)); Y=ones(length(mygrid.LATS),1)*(mygrid.RF');
185     cc=[0:0.5:3 4 5 7 10 15:5:25 35 50]; title0='Meridional Stream Function';
186     if doAnomalies; cc=scaleAnom*[0:0.1:1]; end;
187     figureL; set(gcf,'Renderer','zbuffer'); %set(gcf,'Units','Normalized','Position',[0.05 0.1 0.4 0.8]);
188     depthStretchPlot('pcolor',{X,Y,fld}); shading interp; cbar=gcmfaces_cmap_cbar(cc); title(title0);
189     myCaption={myYmeanTxt,' standard deviation -- overturning streamfunction (Sv)'};
190     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
191     end;
192    
193     if multiTimes&multiBasins;
194     %and the corresponding standard deviation:
195     fld=std(alldiag.atlOV(:,:,tt),[],3); fld(fld==0)=NaN;
196     X=mygrid.LATS*ones(1,length(mygrid.RF)); Y=ones(length(mygrid.LATS),1)*(mygrid.RF');
197     cc=[0:0.5:3 4 5 7 10 15:5:25 35 50]; title0='Atlantic Meridional Stream Function';
198     if doAnomalies; cc=scaleAnom*[0:0.1:1]; end;
199     figureL; set(gcf,'Renderer','zbuffer'); %set(gcf,'Units','Normalized','Position',[0.05 0.1 0.4 0.8]);
200     depthStretchPlot('pcolor',{X,Y,fld}); shading interp; cbar=gcmfaces_cmap_cbar(cc); title(title0);
201     myCaption={myYmeanTxt,' standard deviation -- Atlantic overturning streamfunction (Sv)'};
202     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
203     end;
204    
205     end;
206    
207     if sum(strcmp(choicePlot,'all'))|sum(strcmp(choicePlot,'ovtime'));
208    
209 gforget 1.6 if addToTex; write2tex(fileTex,1,'meridional streamfunction (time series)',2); end;
210    
211 gforget 1.4 %time series
212     if multiTimes;
213     figureL;
214     tmp1=abs(-mygrid.RC-1000); kk=find(tmp1==min(tmp1)); kk=kk(1);
215     gloOV1000=squeeze(alldiag.gloOV(:,kk,:)); gloOV1000=runmean(gloOV1000,myNmean,2);
216     plot(TT,gloOV1000(90+25,:),'LineWidth',2); hold on; plot(TT,gloOV1000(90+35,:),'k','LineWidth',2);
217     plot(TT,gloOV1000(90+45,:),'r','LineWidth',2); plot(TT,gloOV1000(90+55,:),'g','LineWidth',2);
218 heimbach 1.8 title0='annual global overturning at \approx 1000m depth (Sv)'; aa=axis; aa()
219 gforget 1.4 legend('25N','35N','45N','55N'); title(title0);
220     aa(3:4)=[0 20]; axis(aa); grid on;
221     if doAnomalies; aa(3:4)=scaleAnom*[-1 1]; axis(aa); end;
222     myCaption={'annual global overturning at select latitudes at $\\approx$ 1000m depth'};
223     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
224     end;
225    
226     if multiTimes&multiBasins;
227     figureL;
228     tmp1=abs(-mygrid.RC-1000); kk=find(tmp1==min(tmp1)); kk=kk(1);
229     atlOV1000=squeeze(alldiag.atlOV(:,kk,:)); atlOV1000=runmean(atlOV1000,myNmean,2);
230     plot(TT,atlOV1000(90+25,:),'LineWidth',2); hold on; plot(TT,atlOV1000(90+35,:),'k','LineWidth',2);
231     plot(TT,atlOV1000(90+45,:),'r','LineWidth',2); plot(TT,atlOV1000(90+55,:),'g','LineWidth',2);
232     title0='annual atlantic overturning at \approx 1000m depth';
233     legend('25N','35N','45N','55N'); title(title0);
234     aa(3:4)=[0 20]; axis(aa); grid on;
235     if doAnomalies; aa(3:4)=scaleAnom*[-1 1]*2; axis(aa); end;
236 heimbach 1.8 myCaption={'annual Atlantic overturning at select latitudes at $\\approx$ 1000m depth (Sv)'};
237 gforget 1.4 if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
238     end;
239    
240     end;
241    
242     if sum(strcmp(choicePlot,'all'))|sum(strcmp(choicePlot,'mht'));
243    
244 gforget 1.6 if addToTex; write2tex(fileTex,1,'meridional heat transport',2); end;
245    
246 gforget 1.4 %meridional heat transport
247     figureL;
248     fld=mean(alldiag.gloMT_H(:,tt),2);
249     plot(mygrid.LATS,fld,'LineWidth',2);
250     if multiBasins;
251     atl=mean(alldiag.atlMT_H(:,tt),2); pacind=mean(alldiag.pacindMT_H(:,tt),2);
252     hold on; plot(mygrid.LATS,atl,'r','LineWidth',2);
253     plot(mygrid.LATS,pacind,'g','LineWidth',2);
254     legend('global','Atlantic','Pacific+Indian','Location','SouthEast');
255     end;
256     set(gca,'FontSize',14); grid on; axis([-90 90 -2 2]);
257     if doAnomalies; aa=axis; aa(3:4)=scaleAnom*[-1 1]*0.05; axis(aa); end;
258 heimbach 1.7 title('Meridional Heat Transport (in PW)');
259     myCaption={myYmeanTxt,'mean -- meridional heat transport (PW)'};
260 gforget 1.4 if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
261    
262     %and the corresponding standard deviation:
263     if multiTimes;
264     figureL;
265     fld=std(alldiag.gloMT_H(:,tt),[],2);
266     plot(mygrid.LATS,fld,'LineWidth',2);
267     if multiBasins;
268     atl=std(alldiag.atlMT_H(:,tt),[],2);
269     pacind=std(alldiag.pacindMT_H(:,tt),[],2);
270     hold on; plot(mygrid.LATS,atl,'r','LineWidth',2);
271     plot(mygrid.LATS,pacind,'g','LineWidth',2);
272     legend('global','Atlantic','Pacific+Indian');
273     end;
274 gforget 1.10 set(gca,'FontSize',14); grid on; axis([-90 90 0 4]);
275 gforget 1.4 if doAnomalies; aa=axis; aa(3:4)=scaleAnom*[0 1]*0.1; axis(aa); end;
276 heimbach 1.7 title('Meridional Heat Transport (in PW)');
277     myCaption={myYmeanTxt,' standard deviation -- meridional heat transport (PW)'};
278 gforget 1.4 if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
279     end;
280    
281     end;
282    
283     if sum(strcmp(choicePlot,'all'))|sum(strcmp(choicePlot,'mfwt'));
284    
285 heimbach 1.8 if addToTex; write2tex(fileTex,1,'meridional freshwater transport',2); end;
286 gforget 1.6
287 gforget 1.4 %meridional freshwater transport
288     figureL;
289     fld=mean(alldiag.gloMT_FW(:,tt),2);
290     plot(mygrid.LATS,fld,'LineWidth',2);
291     if multiBasins;
292     atl=mean(alldiag.atlMT_FW(:,tt),2); pacind=mean(alldiag.pacindMT_FW(:,tt),2);
293     hold on; plot(mygrid.LATS,atl,'r','LineWidth',2);
294     plot(mygrid.LATS,pacind,'g','LineWidth',2);
295     legend('global','Atlantic','Pacific+Indian');
296     end;
297 gforget 1.10 set(gca,'FontSize',14); grid on; axis([-90 90 -1.5 2.0]);
298 gforget 1.4 if doAnomalies; aa=axis; aa(3:4)=scaleAnom*[-1 1]*0.1; axis(aa); end;
299 heimbach 1.8 title('Meridional Freshwater Transport (in Sv)');
300 gforget 1.4 myCaption={myYmeanTxt,'mean -- meridional freshwater transport (Sv)'};
301     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
302    
303     %and the corresponding standard deviation:
304     if multiTimes;
305     figureL;
306     fld=std(alldiag.gloMT_FW(:,tt),[],2);
307     plot(mygrid.LATS,fld,'LineWidth',2);
308     if multiBasins;
309     atl=std(alldiag.atlMT_FW(:,tt),[],2); pacind=std(alldiag.pacindMT_FW(:,tt),[],2);
310     hold on; plot(mygrid.LATS,atl,'r','LineWidth',2);
311     plot(mygrid.LATS,pacind,'g','LineWidth',2);
312     legend('global','Atlantic','Pacific+Indian');
313     end;
314     set(gca,'FontSize',14); grid on; axis([-90 90 0 2]);
315     if doAnomalies; aa=axis; aa(3:4)=scaleAnom*[0 1]*0.2; axis(aa); end;
316 heimbach 1.8 title('Meridional Freshwater Transport (in Sv)');
317 gforget 1.4 myCaption={myYmeanTxt,' standard deviation -- meridional freshwater transport (Sv)'};
318     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
319     end;
320    
321     end;
322    
323     if sum(strcmp(choicePlot,'all'))|sum(strcmp(choicePlot,'mslt'));
324    
325 gforget 1.6 if addToTex; write2tex(fileTex,1,'meridional salt transport',2); end;
326    
327 gforget 1.4 %meridional salt transport
328     figureL;
329     fld=mean(alldiag.gloMT_SLT(:,tt),2);
330     plot(mygrid.LATS,fld,'LineWidth',2);
331     if sum([90 1170]~=mygrid.ioSize)==0;
332     atl=mean(alldiag.atlMT_SLT(:,tt),2); pacind=mean(alldiag.pacindMT_SLT(:,tt),2);
333     hold on; plot(mygrid.LATS,atl,'r','LineWidth',2);
334     plot(mygrid.LATS,pacind,'g','LineWidth',2);
335     legend('global','Atlantic','Pacific+Indian');
336     end;
337     set(gca,'FontSize',14); grid on; axis([-90 90 -50 50]);
338     if doAnomalies; aa=axis; aa(3:4)=scaleAnom*[-1 1]*2; axis(aa); end;
339     title('Meridional Salt Transport (in psu.Sv)');
340     myCaption={myYmeanTxt,'mean -- meridional salt transport (psu.Sv)'};
341     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
342     %and the corresponding standard deviation:
343    
344     if multiTimes;
345     figureL;
346     fld=std(alldiag.gloMT_SLT(:,tt),[],2);
347     plot(mygrid.LATS,fld,'LineWidth',2);
348     if sum([90 1170]~=mygrid.ioSize)==0;
349     atl=std(alldiag.atlMT_SLT(:,tt),[],2); pacind=std(alldiag.pacindMT_SLT(:,tt),[],2);
350     hold on; plot(mygrid.LATS,atl,'r','LineWidth',2);
351     plot(mygrid.LATS,pacind,'g','LineWidth',2);
352     legend('global','Atlantic','Pacific+Indian');
353     end;
354 gforget 1.10 set(gca,'FontSize',14); grid on; axis([-90 90 0 60]);
355 gforget 1.4 if doAnomalies; aa=axis; aa(3:4)=scaleAnom*[0 1]*2; axis(aa); end;
356     title('Meridional Salt Transport (in psu.Sv)');
357     myCaption={myYmeanTxt,' standard deviation -- meridional salt transport (psu.Sv)'};
358     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
359     end;
360    
361     end;
362    
363     if multiTimes&(sum(strcmp(choicePlot,'all'))|sum(strcmp(choicePlot,'mttime')));
364    
365 gforget 1.6 if addToTex; write2tex(fileTex,1,'meridional transports (time series)',2); end;
366    
367 gforget 1.4 %meridional heat transport
368     fld=squeeze(alldiag.gloMT_H(:,tt))';
369     x=TT*ones(1,length(mygrid.LATS)); y=ones(nt,1)*mygrid.LATS';
370     fld=runmean(fld,myNmean,1);
371     cc=[[-250:50:-100] [-75 -50] [-35:10:35] [50 75] [100:50:250]]/50;
372     if doAnomalies; cc=scaleAnom*[-1:0.1:1]*0.05; end;
373     figureL; pcolor(x,y,fld); shading flat; axis([TT(1) TT(end) -90 90]);
374 heimbach 1.7 gcmfaces_cmap_cbar(cc); title('Meridional Heat Transport (in PW)');
375     myCaption={'meridional heat transport (PW, annual mean)'};
376 gforget 1.4 if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
377    
378     %meridional freshwater transport
379     fld=squeeze(alldiag.gloMT_FW(:,tt))';
380     x=TT*ones(1,length(mygrid.LATS)); y=ones(nt,1)*mygrid.LATS';
381     fld=runmean(fld,myNmean,1);
382     cc=[[-250:50:-100] [-75 -50] [-35:10:35] [50 75] [100:50:250]]/100;
383     if doAnomalies; cc=scaleAnom*[-1:0.1:1]*0.05; end;
384     figureL; pcolor(x,y,fld); shading flat; axis([TT(1) TT(end) -90 90]);
385 heimbach 1.8 gcmfaces_cmap_cbar(cc); title('Meridional Freshwater Transport (in Sv)');
386 gforget 1.4 myCaption={'meridional freshwater transport (Sv, annual mean)'};
387     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
388    
389     %meridional salt transport
390     fld=squeeze(alldiag.gloMT_SLT(:,tt))';
391     x=TT*ones(1,length(mygrid.LATS)); y=ones(nt,1)*mygrid.LATS';
392     fld=runmean(fld,myNmean,1);
393     cc=[[-250:50:-100] [-75 -50] [-35:10:35] [50 75] [100:50:250]]/10;
394     if doAnomalies; cc=scaleAnom*[-1:0.1:1]*1; end;
395     figureL; pcolor(x,y,fld); shading flat; axis([TT(1) TT(end) -90 90]);
396     gcmfaces_cmap_cbar(cc); title('Meridional Salt Transport (in psu.Sv)');
397     myCaption={'meridional salt transport (psu.Sv, annual mean)'};
398     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
399    
400     end;
401    
402    
403     if (sum(strcmp(choicePlot,'all'))|sum(strcmp(choicePlot,'sectime')));
404    
405 gforget 1.6 if addToTex; write2tex(fileTex,1,'transects transport',2); end;
406    
407 gforget 1.4 %Bering Strait and Arctic/Atlantic exchanges:
408     if multiTimes; figureL; end;
409     iiList=[1 8:12]; rrList=[[-1 3];[-3 1];[-7 -2];[0 9];[-4 4];[-0.5 0.5]];
410     for iii=1:length(iiList);
411     ii=iiList(iii);
412     if multiTimes; subplot(3,2,iii); end;
413     trsp=squeeze(alldiag.fldTRANSPORTS(ii,:,:));
414     txt=[mygrid.LINES_MASKS(ii).name ' (>0 to Arctic)'];
415     ylim=rrList(iii,:);
416     if doAnomalies; ylim=scaleAnom*[-1 1]*0.4; end;
417     disp_transport(trsp,TT,txt,{'ylim',ylim},{'nmean',myNmean});
418     end;
419     myCaption={'volume transports entering the Arctic (Sv, annual mean)'};
420 gforget 1.5 if addToTex&multiTimes; write2tex(fileTex,2,myCaption,gcf); end;
421 gforget 1.4
422     %Florida Strait:
423     if multiTimes; figureL; end;
424     iiList=[3 4 6 7]; rrList=[[20 40];[20 40];[-1 3];[-6 2]];
425     for iii=1:length(iiList);
426     ii=iiList(iii);
427     if multiTimes; subplot(2,2,iii); end;
428     trsp=squeeze(alldiag.fldTRANSPORTS(ii,:,:));
429     txt=[mygrid.LINES_MASKS(ii).name ' (>0 to Atlantic)'];
430     ylim=rrList(iii,:);
431     if doAnomalies; ylim=scaleAnom*[-1 1]*2; end;
432     disp_transport(trsp,TT,txt,{'ylim',ylim},{'nmean',myNmean});
433     end;
434     myCaption={'volume transports entering the Atlantic (Sv, annual mean)'};
435 gforget 1.5 if addToTex&multiTimes; write2tex(fileTex,2,myCaption,gcf); end;
436 gforget 1.4
437     %Gibraltar special case:
438     if multiTimes; figureL; end;
439     trsp=squeeze(alldiag.fldTRANSPORTS(2,:,:));
440     txt='Gibraltar Overturn (upper ocean transport towards Med.)';
441 gforget 1.10 ylim=[0.4 1.2];
442 gforget 1.4 if doAnomalies; ylim=scaleAnom*[-1 1]*0.05; end;
443     disp_transport(trsp,TT,txt,{'ylim',ylim},{'nmean',myNmean},{'choicePlot',2});
444     myCaption={'Gibraltar Overturn (Sv, annual mean)'};
445 gforget 1.5 if addToTex&multiTimes; write2tex(fileTex,2,myCaption,gcf); end;
446 gforget 1.4
447     %Drake, ACC etc:
448     if multiTimes; figureL; end;
449     iiList=[13 20 19 18]; rrList=[[120 180];[140 200];[-40 -10];[140 200]];
450     for iii=1:length(iiList);
451     ii=iiList(iii);
452     if multiTimes; subplot(2,2,iii); end;
453     trsp=squeeze(alldiag.fldTRANSPORTS(ii,:,:));
454     txt=[mygrid.LINES_MASKS(ii).name ' (>0 westward)'];
455     ylim=rrList(iii,:);
456     if doAnomalies; ylim=scaleAnom*[-1 1]*5; end;
457     disp_transport(trsp,TT,txt,{'ylim',ylim},{'nmean',myNmean});
458     end;
459     myCaption={'ACC volume transports (Sv, annual mean)'};
460 gforget 1.5 if addToTex&multiTimes; write2tex(fileTex,2,myCaption,gcf); end;
461 gforget 1.4
462     %Indonesian Throughflow special case:
463     if multiTimes; figureL; end;
464     %trsp=-squeeze(nansum(alldiag.fldTRANSPORTS([15:17],:,:),1));
465     trsp=-squeeze(nansum(alldiag.fldTRANSPORTS([14:17],:,:),1));%needed to fix sign for the (small) No.14 transport
466 heimbach 1.8 txt='Indonesian Throughflow (>0 toward Indian Ocean)';
467 gforget 1.4 ylim=[5 25];
468     if doAnomalies; ylim=scaleAnom*[-1 1]*0.5; end;
469     disp_transport(trsp,TT,txt,{'ylim',ylim},{'nmean',myNmean});
470     myCaption={'Indonesian Throughflow (Sv, annual mean)'};
471 gforget 1.5 if addToTex&multiTimes; write2tex(fileTex,2,myCaption,gcf); end;
472 gforget 1.4
473     end;
474    
475    
476 gforget 1.1 end;

  ViewVC Help
Powered by ViewVC 1.1.22