/[MITgcm]/MITgcm_contrib/gael/matlab_class/gcmfaces_maps/m_map_gcmfaces.m
ViewVC logotype

Annotation of /MITgcm_contrib/gael/matlab_class/gcmfaces_maps/m_map_gcmfaces.m

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


Revision 1.10 - (hide annotations) (download)
Wed Feb 12 05:05:40 2014 UTC (11 years, 5 months ago) by gforget
Branch: MAIN
Changes since 1.9: +1 -1 lines
- set myShading='interp'; as default.

1 gforget 1.1 function []=m_map_gcmfaces(fld,varargin);
2 gforget 1.3 %object: gcmfaces front end to m_map
3 gforget 1.7 %inputs: fld is the 2D field to be mapped, or a cell (see below).
4 gforget 1.4 %optional: proj is either the index (integer; 0 by default) of pre-defined
5     % projection(s) or parameters to pass to m_proj (cell)
6 gforget 1.7 %more so: other optional paramaters can be provided ONLY AFTER proj,
7 gforget 1.4 % and must take the following form {'name',param1,param2,...}
8 gforget 1.6 % those that are currently active are
9 gforget 1.4 % {'myCaxis',myCaxis} is the color axis ('auto' by default)
10     % {'do_m_coast',do_m_coast} adds call to m_coast (1; default) or not (0).
11 gforget 1.6 % {'myCmap',myCmap} is the colormap name ('jet' by default)
12 gforget 1.7 % {'doHold',1} indicates to 'hold on' and superimpose e.g. a contour
13 gforget 1.8 % {'doCbar',1} indicates to include the colorbar (1 by default)
14     % {'doLabl',1} indicates to include the axes labels (1 by default)
15     % {'doFit',1} indicates to exclude white space (0 by default)
16 gforget 1.9 % {'myShading','flat'} indicates to use flat shading instead of interp
17 gforget 1.3 %
18 gforget 1.4 %notes: - for proj==0 (i.e. the default) three panels will be plotted :
19     % a mercator projection over mid-latitude, and two polar
20 gforget 1.7 % stereographic projections. The other predefined proj are
21 gforget 1.3 % -1 mollweide or cylindrical projection
22     % 1 mercator projection only
23     % 2 arctic stereographic projection only
24     % 3 antarctic stereographic projection only
25 gforget 1.7 % 1.1 and 1.2 are atlantic mercator projections
26 gforget 1.4 % - myTitle is currently not used; it will take the form
27     % {'myTitle',myTitle} is the title (none by default)
28 gforget 1.7 % - myCaxis can be specified with more than 2 values, in which
29 gforget 1.8 % case gcmfaces_cb will be used.
30 gforget 1.7 % - if fld is a 2D field we use pcolor to display it
31     % - if fld is a cell one can specify the plotting tool.
32     % The cell must then have the form {plotType,FLD,varargin}
33 gforget 1.8 % where plotType is e.g. 'pcolor' or 'contour', FLD is
34     % the 2D field to be plotted, and varargin are options
35 gforget 1.7 % to pass over to the contour command.
36 gforget 1.8 % - Hence .e.g. fld={'contour',mygrid.Depth,'r'} will draw
37 gforget 1.7 % red contours, while fld=mygrid.Depth will color shade.
38 gforget 1.1
39 gforget 1.2 %check that m_map is in the path
40     aa=which('m_proj'); if isempty(aa); error('this function requires m_map that is missing'); end;
41    
42 gforget 1.1 global mygrid;
43    
44 gforget 1.4 %get optional parameters
45     if nargin>1; proj=varargin{1}; else; proj=0; end;
46 gforget 1.7 if iscell(proj);
47     error('not yet implemented');
48 gforget 1.4 else;
49     choicePlot=proj;
50     end;
51 gforget 1.7 %determine the type of plot
52     if iscell(fld); myPlot=fld{1}; else; myPlot='pcolor'; fld={'pcolor',fld}; end;
53 gforget 1.4 %set more optional paramaters to default values
54 gforget 1.10 myCaxis=[]; myTitle=''; myShading='interp'; myCmap='jet';
55 gforget 1.9 do_m_coast=1; doHold=0; doCbar=1; doLabl=1; doFit=0;
56 gforget 1.4 %set more optional paramaters to user defined values
57     for ii=2:nargin-1;
58     if ~iscell(varargin{ii});
59     warning('inputCheck:m_map_gcmfaces_1',...
60     ['As of june 2011, m_map_gcmfaces expects \n'...
61     ' its optional parameters as cell arrays. \n'...
62     ' Argument no. ' num2str(ii+1) ' was ignored \n'...
63     ' Type ''help m_map_gcmfaces'' for details.']);
64     elseif ~ischar(varargin{ii}{1});
65     warning('inputCheck:m_map_gcmfaces_2',...
66     ['As of june 2011, m_map_gcmfaces expects \n'...
67     ' its optional parameters cell arrays \n'...
68     ' to start with character string. \n'...
69     ' Argument no. ' num2str(ii+1) ' was ignored \n'...
70     ' Type ''help m_map_gcmfaces'' for details.']);
71     else;
72     if strcmp(varargin{ii}{1},'myCaxis')|...
73 gforget 1.6 strcmp(varargin{ii}{1},'myCmap')|...
74 gforget 1.9 strcmp(varargin{ii}{1},'myShading')|...
75     strcmp(varargin{ii}{1},'myTitle')|...
76 gforget 1.7 strcmp(varargin{ii}{1},'doHold')|...
77 gforget 1.8 strcmp(varargin{ii}{1},'doCbar')|...
78     strcmp(varargin{ii}{1},'doLabl')|...
79     strcmp(varargin{ii}{1},'doFit')|...
80 gforget 1.9 strcmp(varargin{ii}{1},'do_m_coast');
81 gforget 1.4 eval([varargin{ii}{1} '=varargin{ii}{2};']);
82     else;
83     warning('inputCheck:m_map_gcmfaces_3',...
84     ['unknown option ''' varargin{ii}{1} ''' was ignored']);
85     end;
86     end;
87     end;
88 gforget 1.7
89 gforget 1.4 %make parameter inferences
90 gforget 1.7 if length(myCaxis)==0;
91     plotCBAR=0;
92     elseif length(myCaxis)==2;
93     plotCBAR=1;
94     else;
95     plotCBAR=2;
96     end;
97     %
98     if choicePlot==0&~doHold;
99     clf;
100     elseif ~doHold;
101     cla;
102     else;
103     hold on;
104     end;
105    
106     %re-group param:
107     param.plotCBAR=plotCBAR;
108     param.myCmap=myCmap;
109     param.myCaxis=myCaxis;
110     param.do_m_coast=do_m_coast;
111     param.doHold=doHold;
112 gforget 1.8 param.doCbar=doCbar;
113     param.doLabl=doLabl;
114     param.doFit=doFit;
115 gforget 1.7 param.myPlot=myPlot;
116    
117     %do the plotting:
118     if (choicePlot~=0&choicePlot~=1&choicePlot~=2&choicePlot~=3);
119 gforget 1.9 do_my_plot(fld,param,choicePlot,myShading);
120 gforget 1.7 end;%if choicePlot==0|choicePlot==1;
121 gforget 1.1
122     if choicePlot==0; subplot(2,1,1); end;
123 gforget 1.7 if choicePlot==0|choicePlot==1;
124 gforget 1.8 if mygrid.nFaces~=5&mygrid.nFaces~=6;
125 gforget 1.9 do_my_plot(fld,param,1,myShading);
126 gforget 1.8 else;
127 gforget 1.9 do_my_plot(fld,param,1.1,myShading);
128 gforget 1.8 end;
129 gforget 1.7 end;%if choicePlot==0|choicePlot==1;
130 gforget 1.3
131 gforget 1.1 if choicePlot==0; subplot(2,2,3); end;
132 gforget 1.7 if choicePlot==0|choicePlot==2;
133 gforget 1.9 do_my_plot(fld,param,2,myShading);
134 gforget 1.7 end;%if choicePlot==0|choicePlot==1;
135 gforget 1.1
136     if choicePlot==0; subplot(2,2,4); end;
137 gforget 1.7 if choicePlot==0|choicePlot==3;
138 gforget 1.9 do_my_plot(fld,param,3,myShading);
139 gforget 1.7 end;%if choicePlot==0|choicePlot==1;
140 gforget 1.1
141 gforget 1.8 if plotCBAR==2&strcmp(myPlot,'pcolor')&doCbar;
142 gforget 1.6 cbar=gcmfaces_cmap_cbar(myCaxis,{'myCmap',myCmap});
143 gforget 1.7 if choicePlot==0;
144 gforget 1.3 set(cbar,'Position',[0.92 0.15 0.02 0.75]);
145     elseif choicePlot==-1;
146     set(cbar,'Position',[0.92 0.35 0.02 0.3]);
147 gforget 1.7 elseif choicePlot==1;
148 gforget 1.3 set(cbar,'Position',[0.92 0.34 0.02 0.35]);
149     else;
150     set(cbar,'Position',[0.92 0.23 0.02 0.55]);
151 gforget 1.7 end;
152     end;
153    
154 gforget 1.8 if doFit;
155     if doLabl&~doCbar;
156     set(gca,'LooseInset',[0.05 0.02 0.03 0.03]);
157     else;
158     set(gca,'LooseInset',[0.01 0 0.03 0.03]);
159     end;
160     tmp1=get(gca,'PlotBoxAspectRatio');
161     set(gcf,'PaperPosition',[0 4 8 8/tmp1(1)]);
162     end;
163 gforget 1.7
164 gforget 1.9 function []=do_my_plot(fld,param,proj,shad);
165 gforget 1.7
166     gcmfaces_global;
167    
168 gforget 1.8 %default m_grid params:
169     xloc='bottom'; xtic=[-180:60:180]; xticlab=1;
170     yloc='left'; ytic=[-60:30:60]; yticlab=1;
171    
172     %choice of projection:
173 gforget 1.7 if proj==-1;
174     %%m_proj('Miller Cylindrical','lat',[-90 90]);
175     %m_proj('Equidistant cylindrical','lat',[-90 90]);
176     %m_proj('mollweide','lon',[-180 180],'lat',[-80 80]);
177     m_proj('mollweide','lon',[-180 180],'lat',[-88 88]);
178 gforget 1.8 myConv='pcol'; xticlab=0; yticlab=0;
179 gforget 1.7 elseif proj==1;
180     m_proj('Mercator','lat',[-70 70]);
181 gforget 1.8 myConv='pcol';
182     elseif proj==1.1;
183     m_proj('Mercator','lat',[-70 70],'lon',[0 360]+20);
184     myConv='pcol';
185     xtic=[-360:60:360]; ytic=[-60:30:60];
186 gforget 1.7 elseif proj==2;
187     m_proj('Stereographic','lon',0,'lat',90,'rad',40);
188 gforget 1.8 myConv='convert2arctic';
189     yloc='bottom'; ytic=[50:10:90];
190     elseif proj==2.1;
191     m_proj('Stereographic','lon',0,'lat',90,'rad',60);
192     myConv='convert2arctic';
193     yloc='bottom'; ytic=[30:10:90];
194 gforget 1.7 elseif proj==3;
195     m_proj('Stereographic','lon',0,'lat',-90,'rad',40);
196 gforget 1.8 myConv='convert2southern';
197     xloc='top'; ytic=[-90:10:-50];
198     elseif proj==3.1;
199     m_proj('Stereographic','lon',0,'lat',-90,'rad',60);
200     myConv='convert2southern';
201     xloc='top'; ytic=[-90:10:-30];
202     elseif proj==4.1;
203     m_proj('mollweide','lat',[25 75],'lon',[-100 30]);
204     myConv='pcol';
205     xtic=[-100:20:30]; ytic=[30:10:70];
206     elseif proj==4.2;
207     m_proj('mollweide','lat',[-30 30],'lon',[-65 20]);
208     myConv='pcol';
209     xtic=[-60:10:20]; ytic=[-30:10:30];
210     elseif proj==4.3;
211     m_proj('mollweide','lat',[-75 -25],'lon',[-75 25]);
212     myConv='pcol';
213     xtic=[-70:20:20]; ytic=[-70:10:-30]; xloc='top';
214     elseif proj==4.4;
215     m_proj('mollweide','lat',[25 75],'lon',[-240 -110]);
216     myConv='pcol';
217     xtic=[-240:20:-120]; ytic=[30:10:60];
218     elseif proj==4.5;
219     m_proj('mollweide','lat',[-30 30],'lon',[-240 -70]);
220     myConv='pcol';
221     xtic=[-240:20:-70]; ytic=[-30:10:30];
222     elseif proj==4.6;
223     m_proj('mollweide','lat',[-75 -25],'lon',[-215 -60]);
224     myConv='pcol';
225     xtic=[-210:20:-70]; ytic=[-70:10:-30]; xloc='top';
226     elseif proj==4.7;
227     m_proj('mollweide','lat',[-30 30],'lon',[15 155]);
228     myConv='pcol';
229     xtic=[10:20:160]; ytic=[-30:10:30]; xloc='top';
230     elseif proj==4.8;
231     m_proj('mollweide','lat',[-75 -25],'lon',[15 155]);
232     myConv='pcol';
233     xtic=[10:20:160]; ytic=[-70:10:-30]; xloc='top';
234     else;
235     error('undefined projection');
236 gforget 1.3 end;
237 gforget 1.1
238 gforget 1.8 if ~param.doLabl; xticlab=0; yticlab=0;end; %omit labels
239    
240     m_grid_opt=['''XaxisLocation'',xloc,''YaxisLocation'',yloc,''xtick'',xtic,''ytick'',ytic'];
241     if xticlab==0; m_grid_opt=[m_grid_opt ',''xticklabel'',[]']; end;
242     if yticlab==0; m_grid_opt=[m_grid_opt ',''yticklabel'',[]']; end;
243    
244 gforget 1.7 if strcmp(param.myPlot,'pcolor')|strcmp(param.myPlot,'contour');
245 gforget 1.8 x=mygrid.XC;
246     %mask out the XC padded zeros
247     if isfield(mygrid,'xtrct');
248     pt1=mygrid.xtrct.pt1face; pt2=mygrid.xtrct.pt2face;
249     if pt1~=pt2;
250     for iF=1:6; if iF~=pt1&iF~=pt2; x{iF}(:)=NaN; end; end;
251     end;
252     x(x<0)=x(x<0)+360;
253     end;
254     if strcmp(myConv,'pcol');
255     [xx,yy,z]=convert2pcol(x,mygrid.YC,fld{2});
256 gforget 1.7 else;
257 gforget 1.8 eval(['xx=' myConv '(x);']);
258     eval(['yy=' myConv '(mygrid.YC);']);
259     eval(['z=' myConv '(fld{2});']);
260 gforget 1.7 end;
261     [x,y]=m_ll2xy(xx,yy);
262     if strcmp(param.myPlot,'pcolor');
263 gforget 1.9 if sum(~isnan(x(:)))>0; pcolor(x,y,z); eval(['shading ' shad ';']); end;
264 gforget 1.8 if param.plotCBAR==0;
265     colormap(param.myCmap); if param.doCbar; colorbar; end;
266     elseif param.plotCBAR==1;
267     caxis(param.myCaxis); colormap(param.myCmap); if param.doCbar; colorbar; end;
268     else;
269     cbar=gcmfaces_cmap_cbar(param.myCaxis,{'myCmap',param.myCmap}); delete(cbar);
270 gforget 1.7 end;
271     if param.do_m_coast; m_coast('patch',[1 1 1]*.7,'edgecolor','none'); end;
272 gforget 1.8 eval(['m_grid(' m_grid_opt ');']);
273 gforget 1.7 elseif strcmp(param.myPlot,'contour');
274     if ~param.doHold;
275     if param.do_m_coast; m_coast('patch',[1 1 1]*.7,'edgecolor','none'); end;
276 gforget 1.8 eval(['m_grid(' m_grid_opt ');']);
277 gforget 1.7 end;
278     if length(fld)==2; fld{3}='k'; end;
279     hold on; contour(x,y,z,fld{3:end});
280     end;
281     elseif strcmp(param.myPlot,'plot');
282     if ~param.doHold;
283     if param.do_m_coast; m_coast('patch',[1 1 1]*.7,'edgecolor','none'); end;
284 gforget 1.8 eval(['m_grid(' m_grid_opt ');']);
285 gforget 1.7 end;
286     [x,y]=m_ll2xy(fld{2},fld{3});
287     hold on; plot(x,y,fld{4:end});
288 gforget 1.8 elseif strcmp(param.myPlot,'text');
289     if ~param.doHold;
290     if param.do_m_coast; m_coast('patch',[1 1 1]*.7,'edgecolor','none'); end;
291     eval(['m_grid(' m_grid_opt ');']);
292     end;
293     [x,y]=m_ll2xy(fld{2},fld{3});
294     if length(fld)>4; cc=fld{5}; else; cc='k'; end;
295     hold on; hold on; text(x,y,fld{4},'Color',cc,'FontSize',16,'FontWeight','bold');
296 gforget 1.7 end;
297 gforget 1.1
298 gforget 1.8 %add tag spec. to map & proj generated with this routine
299     set(gca,'Tag',['gfMap' num2str(proj)]);
300    

  ViewVC Help
Powered by ViewVC 1.1.22