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

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

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


Revision 1.1 - (hide annotations) (download)
Sat Nov 8 18:21:38 2014 UTC (10 years, 8 months ago) by gforget
Branch: MAIN
CVS Tags: checkpoint65r, checkpoint65p, checkpoint65q
- simple display on the sphere, with option for
  overlaid grid lines (without additional toolbox).

1 gforget 1.1 function []=gcmfaces_sphere(fld,cc,nn,vw,co);
2     %object: display field and grid mesh on the sphere
3     %inputs: fld is the field to be plotted
4     % cc is the color range(s)
5     % nn subsampling rate for grid lines ([] for no grid lines)
6     % vw is either [],{az,elev},'N','S'
7     % co color scheme option (1,2, or 3)
8    
9     %mask=mygrid.mskC(:,:,1);
10     %bathy=-mygrid.Depth;
11     %fld=-bathy/1e4.*mask);
12     %cc=[[0:0.05:0.5] [0.6 0.75 1 1.25]];
13    
14     gcmfaces_global;
15    
16     if isempty(who('nn')); nn=[]; end;
17     if isempty(who('vw')); vw=[]; end;
18    
19     if iscell(vw);
20     az=vw{1}; elev=vw{2};
21     elseif isempty(vw);
22     az=0; elev=20;
23     elseif strcmp(vw,'N');
24     az=0; elev=90;
25     elseif strcmp(vw,'S');
26     az=0; elev=-90;
27     end;
28    
29     %overlay the mesh:
30     if ~isempty(nn);
31     myfac=1.001;
32     elev=myfac*elev;
33     end;
34     if co==1;
35     listCol='kkkkkk';
36     elseif co==2;
37     listCol='bbbbbb';
38     elseif co==3;
39     listCol='gmbcrk'; if mygrid.nFaces==4; listCol='rgmc'; end;
40     %%listCol='rkbcmg'; if mygrid.nFaces==4; listCol='mrkc'; end;
41     end;
42    
43     figureL;
44    
45     %underlay plain sphere:
46     [xs,ys,zs]=sphere(100); aa=0.99;
47     aa=surf(aa*xs,aa*ys,aa*zs); hold on;
48     set(aa,'LineStyle','none','FaceColor',[1 1 1]*0.7);
49    
50     %plot bathy:
51     test1=1; test2=0;
52     XC=exch_T_N(mygrid.XC); YC=exch_T_N(mygrid.YC);
53     FLD=exch_T_N(fld);
54     for ii=1:mygrid.nFaces;
55     x=sin(pi/2-YC{ii}*pi/180).*cos(XC{ii}*pi/180);
56     y=sin(pi/2-YC{ii}*pi/180).*sin(XC{ii}*pi/180);
57     z=cos(pi/2-YC{ii}*pi/180);
58     c=FLD{ii};
59     %subsampling
60     k0=1; i0=1:k0:size(x,1); j0=1:k0:size(x,2);
61     surf(x(i0,j0),y(i0,j0),z(i0,j0),c(i0,j0),'LineStyle','none'); hold on;
62     test1=nanmin([test1;c(:)]); test2=nanmax([test2;c(:)]);
63     end;%for ii=1:mygrid.nFaces;
64     axis equal;
65     view(az,elev);
66     %
67     if co==1;
68     cb=gcmfaces_cmap_cbar(cc);
69     elseif co==2;
70     cb=gcmfaces_cmap_cbar(cc,{'myBW',-2},{'myCmap','hot'});
71     elseif co==3;
72     cb=gcmfaces_cmap_cbar(cc,{'myBW',-2},{'myCmap','gray'});
73     end;
74     %
75     set(cb,'Position',[0.78 0.2 0.02 0.6]);
76     %delete(cb);
77    
78     if ~isempty(nn);
79     XG=exch_Z(mygrid.XG); YG=exch_Z(mygrid.YG);
80     %now do the plot:
81     for ii=1:mygrid.nFaces;
82     x=myfac*sin(pi/2-YG{ii}*pi/180).*cos(XG{ii}*pi/180);
83     y=myfac*sin(pi/2-YG{ii}*pi/180).*sin(XG{ii}*pi/180);
84     z=myfac*cos(pi/2-YG{ii}*pi/180);
85     [n1,n2]=size(x); i1=[1:nn:n1]; i2=[1:nn:n2];
86     mesh(x(i1,i2),y(i1,i2),z(i1,i2),0*x(i1,i2),'FaceColor','none','EdgeColor',listCol(ii)); hold on;
87     end;
88     end;
89    
90     %find tune plot:
91     axis equal; axis off;
92     %set(gca,'XTick',[],'YTick',[],'ZTick',[]);
93     %set(gca,'Color','none');
94     %set(gcf,'Renderer','zbuffer');
95     %box on;
96    
97    

  ViewVC Help
Powered by ViewVC 1.1.22