/[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.5 - (hide annotations) (download)
Mon Jan 7 22:13:01 2013 UTC (12 years, 6 months ago) by gforget
Branch: MAIN
Changes since 1.4: +7 -7 lines
- set latitude limits to -35,+70 for atlExt, and -35,+65 for pacExt
  consistent with Patricks extended atl mask.

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     listFiles={'trsp_3d_set1','trsp_3d_set2'};
14     listSubdirs={[dirModel 'diags/TRSP/']};
15     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     %barotropic streamfunction:
116     fld=mean(alldiag.fldBAR(:,:,tt),3);
117     cc=[[-80:20:-40] [-25 -15:5:15 25] [40:40:200]]; title0='Horizontal Stream Function';
118     if doAnomalies; cc=scaleAnom*[-5:0.5:5]; end;
119     figureL; m_map_gcmfaces(fld,0,{'myCaxis',cc},{'do_m_coast',1},{'myTitle',title0});
120     myCaption={myYmeanTxt,'mean -- barotropic streamfunction (Sv)'};
121     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
122     %and the corresponding standard deviation:
123     if multiTimes;
124     fld=std(alldiag.fldBAR(:,:,tt),[],3);
125     cc=[0:0.5:3 4 5 7 10 15:5:25 35 50]; title0='Horizontal Stream Function';
126     if doAnomalies; cc=scaleAnom*[0:0.25:2]; end;
127     figureL; m_map_gcmfaces(fld,0,{'myCaxis',cc},{'do_m_coast',1},{'myTitle',title0});
128     myCaption={myYmeanTxt,' standard deviation -- barotropic streamfunction (Sv)'};
129     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
130     end;
131    
132     end;
133    
134     if sum(strcmp(choicePlot,'all'))|sum(strcmp(choicePlot,'ov'));
135    
136     %meridional streamfunction (Eulerian) :
137     fld=mean(alldiag.gloOV(:,:,tt),3); fld(fld==0)=NaN;
138     X=mygrid.LATS*ones(1,length(mygrid.RF)); Y=ones(length(mygrid.LATS),1)*(mygrid.RF');
139     cc=[[-50:10:-30] [-24:3:24] [30:10:50]]; title0='Meridional Stream Function';
140     if doAnomalies; cc=scaleAnom*[-1:0.1:1]; end;
141     figureL; set(gcf,'Renderer','zbuffer'); %set(gcf,'Units','Normalized','Position',[0.05 0.1 0.4 0.8]);
142     depthStretchPlot('pcolor',{X,Y,fld}); shading interp; cbar=gcmfaces_cmap_cbar(cc); title(title0);
143     myCaption={myYmeanTxt,'mean -- overturning streamfunction (Sv)'};
144     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
145    
146     %meridional streamfunction (residual) :
147     fld=mean(alldiag.gloOVres(:,:,tt),3); fld(fld==0)=NaN;
148     X=mygrid.LATS*ones(1,length(mygrid.RF)); Y=ones(length(mygrid.LATS),1)*(mygrid.RF');
149     cc=[[-50:10:-30] [-24:3:24] [30:10:50]]; title0='Meridional Stream Function (incl. GM)';
150     if doAnomalies; cc=scaleAnom*[-1:0.1:1]; end;
151     figureL; set(gcf,'Renderer','zbuffer'); %set(gcf,'Units','Normalized','Position',[0.05 0.1 0.4 0.8]);
152     depthStretchPlot('pcolor',{X,Y,fld}); shading interp; cbar=gcmfaces_cmap_cbar(cc); title(title0);
153     myCaption={myYmeanTxt,'mean -- overturning streamfunction incl. GM (Sv)'};
154     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
155    
156     if multiBasins;
157     %meridional streamfunction (Eulerian):
158     fld=mean(alldiag.atlOV(:,:,tt),3); fld(fld==0)=NaN;
159     X=mygrid.LATS*ones(1,length(mygrid.RF)); Y=ones(length(mygrid.LATS),1)*(mygrid.RF');
160     cc=[[-50:10:-30] [-24:3:24] [30:10:50]]; title0='Atlantic Meridional Stream Function';
161     if doAnomalies; cc=scaleAnom*[-1:0.1:1]; end;
162     figureL; set(gcf,'Renderer','zbuffer'); %set(gcf,'Units','Normalized','Position',[0.05 0.1 0.4 0.8]);
163     depthStretchPlot('pcolor',{X,Y,fld}); shading interp; cbar=gcmfaces_cmap_cbar(cc); title(title0);
164     myCaption={myYmeanTxt,'mean -- Atlantic overturning streamfunction (Sv)'};
165     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
166     %meridional streamfunction (residual):
167     fld=mean(alldiag.pacindOV(:,:,tt),3); fld(fld==0)=NaN;
168     X=mygrid.LATS*ones(1,length(mygrid.RF)); Y=ones(length(mygrid.LATS),1)*(mygrid.RF');
169     cc=[[-50:10:-30] [-24:3:24] [30:10:50]]; title0='Pac+Ind Meridional Stream Function';
170     if doAnomalies; cc=scaleAnom*[-1:0.1:1]; end;
171     figureL; set(gcf,'Renderer','zbuffer'); %set(gcf,'Units','Normalized','Position',[0.05 0.1 0.4 0.8]);
172     depthStretchPlot('pcolor',{X,Y,fld}); shading interp; cbar=gcmfaces_cmap_cbar(cc); title(title0);
173     myCaption={myYmeanTxt,'mean -- Pac+Ind overturning streamfunction (Sv)'};
174     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
175     end;
176    
177     if multiTimes;
178     %and the corresponding standard deviation:
179     fld=std(alldiag.gloOV(:,:,tt),[],3); fld(fld==0)=NaN;
180     X=mygrid.LATS*ones(1,length(mygrid.RF)); Y=ones(length(mygrid.LATS),1)*(mygrid.RF');
181     cc=[0:0.5:3 4 5 7 10 15:5:25 35 50]; title0='Meridional Stream Function';
182     if doAnomalies; cc=scaleAnom*[0:0.1:1]; end;
183     figureL; set(gcf,'Renderer','zbuffer'); %set(gcf,'Units','Normalized','Position',[0.05 0.1 0.4 0.8]);
184     depthStretchPlot('pcolor',{X,Y,fld}); shading interp; cbar=gcmfaces_cmap_cbar(cc); title(title0);
185     myCaption={myYmeanTxt,' standard deviation -- overturning streamfunction (Sv)'};
186     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
187     end;
188    
189     if multiTimes&multiBasins;
190     %and the corresponding standard deviation:
191     fld=std(alldiag.atlOV(:,:,tt),[],3); fld(fld==0)=NaN;
192     X=mygrid.LATS*ones(1,length(mygrid.RF)); Y=ones(length(mygrid.LATS),1)*(mygrid.RF');
193     cc=[0:0.5:3 4 5 7 10 15:5:25 35 50]; title0='Atlantic Meridional Stream Function';
194     if doAnomalies; cc=scaleAnom*[0:0.1:1]; end;
195     figureL; set(gcf,'Renderer','zbuffer'); %set(gcf,'Units','Normalized','Position',[0.05 0.1 0.4 0.8]);
196     depthStretchPlot('pcolor',{X,Y,fld}); shading interp; cbar=gcmfaces_cmap_cbar(cc); title(title0);
197     myCaption={myYmeanTxt,' standard deviation -- Atlantic overturning streamfunction (Sv)'};
198     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
199     end;
200    
201     end;
202    
203     if sum(strcmp(choicePlot,'all'))|sum(strcmp(choicePlot,'ovtime'));
204    
205     %time series
206     if multiTimes;
207     figureL;
208     tmp1=abs(-mygrid.RC-1000); kk=find(tmp1==min(tmp1)); kk=kk(1);
209     gloOV1000=squeeze(alldiag.gloOV(:,kk,:)); gloOV1000=runmean(gloOV1000,myNmean,2);
210     plot(TT,gloOV1000(90+25,:),'LineWidth',2); hold on; plot(TT,gloOV1000(90+35,:),'k','LineWidth',2);
211     plot(TT,gloOV1000(90+45,:),'r','LineWidth',2); plot(TT,gloOV1000(90+55,:),'g','LineWidth',2);
212     title0='annual global overturning at \approx 1000m depth'; aa=axis; aa()
213     legend('25N','35N','45N','55N'); title(title0);
214     aa(3:4)=[0 20]; axis(aa); grid on;
215     if doAnomalies; aa(3:4)=scaleAnom*[-1 1]; axis(aa); end;
216     myCaption={'annual global overturning at select latitudes at $\\approx$ 1000m depth'};
217     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
218     end;
219    
220     if multiTimes&multiBasins;
221     figureL;
222     tmp1=abs(-mygrid.RC-1000); kk=find(tmp1==min(tmp1)); kk=kk(1);
223     atlOV1000=squeeze(alldiag.atlOV(:,kk,:)); atlOV1000=runmean(atlOV1000,myNmean,2);
224     plot(TT,atlOV1000(90+25,:),'LineWidth',2); hold on; plot(TT,atlOV1000(90+35,:),'k','LineWidth',2);
225     plot(TT,atlOV1000(90+45,:),'r','LineWidth',2); plot(TT,atlOV1000(90+55,:),'g','LineWidth',2);
226     title0='annual atlantic overturning at \approx 1000m depth';
227     legend('25N','35N','45N','55N'); title(title0);
228     aa(3:4)=[0 20]; axis(aa); grid on;
229     if doAnomalies; aa(3:4)=scaleAnom*[-1 1]*2; axis(aa); end;
230     myCaption={'annual atlantic overturning at select latitudes at $\\approx$ 1000m depth'};
231     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
232     end;
233    
234     end;
235    
236     if sum(strcmp(choicePlot,'all'))|sum(strcmp(choicePlot,'mht'));
237    
238     %meridional heat transport
239     figureL;
240     fld=mean(alldiag.gloMT_H(:,tt),2);
241     plot(mygrid.LATS,fld,'LineWidth',2);
242     if multiBasins;
243     atl=mean(alldiag.atlMT_H(:,tt),2); pacind=mean(alldiag.pacindMT_H(:,tt),2);
244     hold on; plot(mygrid.LATS,atl,'r','LineWidth',2);
245     plot(mygrid.LATS,pacind,'g','LineWidth',2);
246     legend('global','Atlantic','Pacific+Indian','Location','SouthEast');
247     end;
248     set(gca,'FontSize',14); grid on; axis([-90 90 -2 2]);
249     if doAnomalies; aa=axis; aa(3:4)=scaleAnom*[-1 1]*0.05; axis(aa); end;
250     title('Meridional Heat Transport (in pW)');
251     myCaption={myYmeanTxt,'mean -- meridional heat transport (pW)'};
252     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
253    
254     %and the corresponding standard deviation:
255     if multiTimes;
256     figureL;
257     fld=std(alldiag.gloMT_H(:,tt),[],2);
258     plot(mygrid.LATS,fld,'LineWidth',2);
259     if multiBasins;
260     atl=std(alldiag.atlMT_H(:,tt),[],2);
261     pacind=std(alldiag.pacindMT_H(:,tt),[],2);
262     hold on; plot(mygrid.LATS,atl,'r','LineWidth',2);
263     plot(mygrid.LATS,pacind,'g','LineWidth',2);
264     legend('global','Atlantic','Pacific+Indian');
265     end;
266     set(gca,'FontSize',14); grid on; axis([-90 90 0 3]);
267     if doAnomalies; aa=axis; aa(3:4)=scaleAnom*[0 1]*0.1; axis(aa); end;
268     title('Meridional Heat Transport (in pW)');
269     myCaption={myYmeanTxt,' standard deviation -- meridional heat transport (pW)'};
270     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
271     end;
272    
273     end;
274    
275     if sum(strcmp(choicePlot,'all'))|sum(strcmp(choicePlot,'mfwt'));
276    
277     %meridional freshwater transport
278     figureL;
279     fld=mean(alldiag.gloMT_FW(:,tt),2);
280     plot(mygrid.LATS,fld,'LineWidth',2);
281     if multiBasins;
282     atl=mean(alldiag.atlMT_FW(:,tt),2); pacind=mean(alldiag.pacindMT_FW(:,tt),2);
283     hold on; plot(mygrid.LATS,atl,'r','LineWidth',2);
284     plot(mygrid.LATS,pacind,'g','LineWidth',2);
285     legend('global','Atlantic','Pacific+Indian');
286     end;
287     set(gca,'FontSize',14); grid on; axis([-90 90 -1.5 1.5]);
288     if doAnomalies; aa=axis; aa(3:4)=scaleAnom*[-1 1]*0.1; axis(aa); end;
289     title('Meridional Fresh Water Transport (in Sv)');
290     myCaption={myYmeanTxt,'mean -- meridional freshwater transport (Sv)'};
291     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
292    
293     %and the corresponding standard deviation:
294     if multiTimes;
295     figureL;
296     fld=std(alldiag.gloMT_FW(:,tt),[],2);
297     plot(mygrid.LATS,fld,'LineWidth',2);
298     if multiBasins;
299     atl=std(alldiag.atlMT_FW(:,tt),[],2); pacind=std(alldiag.pacindMT_FW(:,tt),[],2);
300     hold on; plot(mygrid.LATS,atl,'r','LineWidth',2);
301     plot(mygrid.LATS,pacind,'g','LineWidth',2);
302     legend('global','Atlantic','Pacific+Indian');
303     end;
304     set(gca,'FontSize',14); grid on; axis([-90 90 0 2]);
305     if doAnomalies; aa=axis; aa(3:4)=scaleAnom*[0 1]*0.2; axis(aa); end;
306     title('Meridional Fresh Water Transport (in Sv)');
307     myCaption={myYmeanTxt,' standard deviation -- meridional freshwater transport (Sv)'};
308     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
309     end;
310    
311     end;
312    
313     if sum(strcmp(choicePlot,'all'))|sum(strcmp(choicePlot,'mslt'));
314    
315     %meridional salt transport
316     figureL;
317     fld=mean(alldiag.gloMT_SLT(:,tt),2);
318     plot(mygrid.LATS,fld,'LineWidth',2);
319     if sum([90 1170]~=mygrid.ioSize)==0;
320     atl=mean(alldiag.atlMT_SLT(:,tt),2); pacind=mean(alldiag.pacindMT_SLT(:,tt),2);
321     hold on; plot(mygrid.LATS,atl,'r','LineWidth',2);
322     plot(mygrid.LATS,pacind,'g','LineWidth',2);
323     legend('global','Atlantic','Pacific+Indian');
324     end;
325     set(gca,'FontSize',14); grid on; axis([-90 90 -50 50]);
326     if doAnomalies; aa=axis; aa(3:4)=scaleAnom*[-1 1]*2; axis(aa); end;
327     title('Meridional Salt Transport (in psu.Sv)');
328     myCaption={myYmeanTxt,'mean -- meridional salt transport (psu.Sv)'};
329     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
330     %and the corresponding standard deviation:
331    
332     if multiTimes;
333     figureL;
334     fld=std(alldiag.gloMT_SLT(:,tt),[],2);
335     plot(mygrid.LATS,fld,'LineWidth',2);
336     if sum([90 1170]~=mygrid.ioSize)==0;
337     atl=std(alldiag.atlMT_SLT(:,tt),[],2); pacind=std(alldiag.pacindMT_SLT(:,tt),[],2);
338     hold on; plot(mygrid.LATS,atl,'r','LineWidth',2);
339     plot(mygrid.LATS,pacind,'g','LineWidth',2);
340     legend('global','Atlantic','Pacific+Indian');
341     end;
342     set(gca,'FontSize',14); grid on; axis([-90 90 0 40]);
343     if doAnomalies; aa=axis; aa(3:4)=scaleAnom*[0 1]*2; axis(aa); end;
344     title('Meridional Salt Transport (in psu.Sv)');
345     myCaption={myYmeanTxt,' standard deviation -- meridional salt transport (psu.Sv)'};
346     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
347     end;
348    
349     end;
350    
351     if multiTimes&(sum(strcmp(choicePlot,'all'))|sum(strcmp(choicePlot,'mttime')));
352    
353     %meridional heat transport
354     fld=squeeze(alldiag.gloMT_H(:,tt))';
355     x=TT*ones(1,length(mygrid.LATS)); y=ones(nt,1)*mygrid.LATS';
356     fld=runmean(fld,myNmean,1);
357     cc=[[-250:50:-100] [-75 -50] [-35:10:35] [50 75] [100:50:250]]/50;
358     if doAnomalies; cc=scaleAnom*[-1:0.1:1]*0.05; end;
359     figureL; pcolor(x,y,fld); shading flat; axis([TT(1) TT(end) -90 90]);
360     gcmfaces_cmap_cbar(cc); title('Meridional Heat Transport (in pW)');
361     myCaption={'meridional heat transport (pW, annual mean)'};
362     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
363    
364     %meridional freshwater transport
365     fld=squeeze(alldiag.gloMT_FW(:,tt))';
366     x=TT*ones(1,length(mygrid.LATS)); y=ones(nt,1)*mygrid.LATS';
367     fld=runmean(fld,myNmean,1);
368     cc=[[-250:50:-100] [-75 -50] [-35:10:35] [50 75] [100:50:250]]/100;
369     if doAnomalies; cc=scaleAnom*[-1:0.1:1]*0.05; end;
370     figureL; pcolor(x,y,fld); shading flat; axis([TT(1) TT(end) -90 90]);
371     gcmfaces_cmap_cbar(cc); title('Meridional Fresh Water Transport (in Sv)');
372     myCaption={'meridional freshwater transport (Sv, annual mean)'};
373     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
374    
375     %meridional salt transport
376     fld=squeeze(alldiag.gloMT_SLT(:,tt))';
377     x=TT*ones(1,length(mygrid.LATS)); y=ones(nt,1)*mygrid.LATS';
378     fld=runmean(fld,myNmean,1);
379     cc=[[-250:50:-100] [-75 -50] [-35:10:35] [50 75] [100:50:250]]/10;
380     if doAnomalies; cc=scaleAnom*[-1:0.1:1]*1; end;
381     figureL; pcolor(x,y,fld); shading flat; axis([TT(1) TT(end) -90 90]);
382     gcmfaces_cmap_cbar(cc); title('Meridional Salt Transport (in psu.Sv)');
383     myCaption={'meridional salt transport (psu.Sv, annual mean)'};
384     if addToTex; write2tex(fileTex,2,myCaption,gcf); end;
385    
386     end;
387    
388    
389     if (sum(strcmp(choicePlot,'all'))|sum(strcmp(choicePlot,'sectime')));
390    
391     %Bering Strait and Arctic/Atlantic exchanges:
392     if multiTimes; figureL; end;
393     iiList=[1 8:12]; rrList=[[-1 3];[-3 1];[-7 -2];[0 9];[-4 4];[-0.5 0.5]];
394     for iii=1:length(iiList);
395     ii=iiList(iii);
396     if multiTimes; subplot(3,2,iii); end;
397     trsp=squeeze(alldiag.fldTRANSPORTS(ii,:,:));
398     txt=[mygrid.LINES_MASKS(ii).name ' (>0 to Arctic)'];
399     ylim=rrList(iii,:);
400     if doAnomalies; ylim=scaleAnom*[-1 1]*0.4; end;
401     disp_transport(trsp,TT,txt,{'ylim',ylim},{'nmean',myNmean});
402     end;
403     myCaption={'volume transports entering the Arctic (Sv, annual mean)'};
404 gforget 1.5 if addToTex&multiTimes; write2tex(fileTex,2,myCaption,gcf); end;
405 gforget 1.4
406     %Florida Strait:
407     if multiTimes; figureL; end;
408     iiList=[3 4 6 7]; rrList=[[20 40];[20 40];[-1 3];[-6 2]];
409     for iii=1:length(iiList);
410     ii=iiList(iii);
411     if multiTimes; subplot(2,2,iii); end;
412     trsp=squeeze(alldiag.fldTRANSPORTS(ii,:,:));
413     txt=[mygrid.LINES_MASKS(ii).name ' (>0 to Atlantic)'];
414     ylim=rrList(iii,:);
415     if doAnomalies; ylim=scaleAnom*[-1 1]*2; end;
416     disp_transport(trsp,TT,txt,{'ylim',ylim},{'nmean',myNmean});
417     end;
418     myCaption={'volume transports entering the Atlantic (Sv, annual mean)'};
419 gforget 1.5 if addToTex&multiTimes; write2tex(fileTex,2,myCaption,gcf); end;
420 gforget 1.4
421     %Gibraltar special case:
422     if multiTimes; figureL; end;
423     trsp=squeeze(alldiag.fldTRANSPORTS(2,:,:));
424     txt='Gibraltar Overturn (upper ocean transport towards Med.)';
425     ylim=[0 2];
426     if doAnomalies; ylim=scaleAnom*[-1 1]*0.05; end;
427     disp_transport(trsp,TT,txt,{'ylim',ylim},{'nmean',myNmean},{'choicePlot',2});
428     myCaption={'Gibraltar Overturn (Sv, annual mean)'};
429 gforget 1.5 if addToTex&multiTimes; write2tex(fileTex,2,myCaption,gcf); end;
430 gforget 1.4
431     %Drake, ACC etc:
432     if multiTimes; figureL; end;
433     iiList=[13 20 19 18]; rrList=[[120 180];[140 200];[-40 -10];[140 200]];
434     for iii=1:length(iiList);
435     ii=iiList(iii);
436     if multiTimes; subplot(2,2,iii); end;
437     trsp=squeeze(alldiag.fldTRANSPORTS(ii,:,:));
438     txt=[mygrid.LINES_MASKS(ii).name ' (>0 westward)'];
439     ylim=rrList(iii,:);
440     if doAnomalies; ylim=scaleAnom*[-1 1]*5; end;
441     disp_transport(trsp,TT,txt,{'ylim',ylim},{'nmean',myNmean});
442     end;
443     myCaption={'ACC volume transports (Sv, annual mean)'};
444 gforget 1.5 if addToTex&multiTimes; write2tex(fileTex,2,myCaption,gcf); end;
445 gforget 1.4
446     %Indonesian Throughflow special case:
447     if multiTimes; figureL; end;
448     %trsp=-squeeze(nansum(alldiag.fldTRANSPORTS([15:17],:,:),1));
449     trsp=-squeeze(nansum(alldiag.fldTRANSPORTS([14:17],:,:),1));%needed to fix sign for the (small) No.14 transport
450     txt='Indonesian Throughflow (>0 toward Indian)';
451     ylim=[5 25];
452     if doAnomalies; ylim=scaleAnom*[-1 1]*0.5; end;
453     disp_transport(trsp,TT,txt,{'ylim',ylim},{'nmean',myNmean});
454     myCaption={'Indonesian Throughflow (Sv, annual mean)'};
455 gforget 1.5 if addToTex&multiTimes; write2tex(fileTex,2,myCaption,gcf); end;
456 gforget 1.4
457     end;
458    
459    
460 gforget 1.1 end;

  ViewVC Help
Powered by ViewVC 1.1.22