/[MITgcm]/MITgcm_contrib/high_res_cube/matlab/mk_psiLine_CS.m
ViewVC logotype

Annotation of /MITgcm_contrib/high_res_cube/matlab/mk_psiLine_CS.m

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


Revision 1.1 - (hide annotations) (download)
Wed Mar 28 18:07:24 2007 UTC (18 years, 3 months ago) by dimitri
Branch: MAIN
CVS Tags: HEAD
Files contributed by JM for computation of barotropic streamfunction on cs510.

1 dimitri 1.1 kgr=0; nRgr=3; kwr=1;
2     dbug=1;
3    
4     % $Header: /u/gcmpack/MITgcm/utils/matlab/cs_grid/bk_line/mk_psiLine_CS.m,v 1.1 2005/09/15 16:46:28 jmc Exp $
5     % $Name: $
6    
7     %rac='/home/jmc/grid_cs32/';
8     rac='grid_files/';
9     load_cs;
10     nc=size(xcs,2); nPt=6*nc*nc ; ncx=6*nc; nPt2=nPt+2; ncp=nc+1;
11    
12     xx2=split_Z_cub(xcg);
13     yy2=split_Z_cub(ycg);
14     xx3=reshape(xx2,ncp*ncp*6,1); yy3=reshape(yy2,ncp*ncp*6,1);
15    
16     %- make a "small list" of Z point (= with no double point)
17     xx1=zeros(nPt2,1); yy1=zeros(nPt2,1);
18     xx1(1:nPt,1)=reshape(xcg,nPt,1);
19     yy1(1:nPt,1)=reshape(ycg,nPt,1);
20    
21     %- add to the list the 2 missing corners:
22     ij1=1+nc*ncp; ij2=ncp+ncp*ncp;
23     xx1(nPt+1,1)=xx2(1,ncp,1); xx1(nPt+2,1)=xx2(ncp,1,2);
24     yy1(nPt+1,1)=yy2(1,ncp,1); yy1(nPt+2,1)=yy2(ncp,1,2);
25     fprintf(' add 1rst corner: %i %8.3f %8.3f\n',nPt+1,xx1(nPt+1,1),yy1(nPt+1,1));
26     fprintf(' add 2nd corner: %i %8.3f %8.3f\n',nPt+2,xx1(nPt+2,1),yy1(nPt+2,1));
27    
28     %-- start from the N.Pole:
29     list_C=zeros(ncx,ncx); list_P=zeros(ncx,ncx);
30     listUV=zeros(ncx,ncx); list_T=zeros(ncx,ncx);
31     list_N=zeros(1,ncx); tagC=zeros(nPt2,1);
32     [IJ]=find(yy1==max(yy1));
33     nRd=1; list_C(1,nRd)=IJ(1); tagC(IJ(1),1)=nRd;
34     nbp=length(IJ); list_N(1,nRd)=nbp;
35    
36     doRd=1;
37     fprintf('nRd= %i ; nbp= %i \n',nRd,nbp);
38     while doRd > 0,
39     %- next round :
40     n2p=0;
41     list_c=zeros(ncx,1); list_p=zeros(ncx,1);
42     listuv=zeros(ncx,1); list_t=zeros(ncx,1);
43    
44     %- for all point of the previous list, find the neighbour point not yet tagged:
45     for n=1:nbp,
46     %- find the 4 neighbours:
47     ij=list_C(n,nRd);
48     %i=1+rem(ij-1,nc); j=fix((ij-1)/nc); k=1+fix(j/nc); j=1+rem(j,nc);
49     i=rem(ij-1,ncx); j=1+fix((ij-1)/ncx); k=1+fix(i/nc); i=1+rem(i,nc);
50     if j == ncp | i*j==1,
51     %-- Corner point:
52     nbv=3;
53     [I]=find( xx3==xx1(ij) & yy3==yy1(ij) ) ;
54     if dbug>0, fprintf('Corner: i,j,k= %2i %2i %i ',i,j,k); end
55     if length(I) ~= 3, fprintf('Pb C: stop\n'); return; end
56     for l=1:length(I),
57     i2=I(l);
58     j2=fix((i2-1)/ncp); k2=1+fix(j2/ncp); j2=1+rem(j2,ncp);
59     i2=1+rem(i2-1,ncp);
60     if dbug>0, fprintf('; l=%i: %2i %2i %i ',l,i2,j2,k2); end
61     if ij==nPt+2,
62     %- 3 times (ncp,1) :
63     lc(l)=nc+(k2-1)*nc+(j2-1)*ncx; luv(l)=lc(l); typ(l)=-2;
64     elseif ij==nPt+1,
65     %- 3 times (1,ncp) :
66     lc(l)=i2+(k2-1)*nc+(nc-1)*ncx; luv(l)=lc(l); typ(l)=1;
67     elseif i2==ncp & j2==1,
68     lc(2)=nc+(k2-1)*nc+(j2-1)*ncx; luv(2)=lc(2); typ(2)=-2;
69     elseif i2==1 & j2==ncp,
70     lc(2)=i2+(k2-1)*nc+(nc-1)*ncx; luv(2)=lc(2); typ(2)=1;
71     end
72     end
73     if dbug>0, fprintf('\n'); end
74     if j < ncp,
75     lc(1)=ij+1; luv(1)=ij; typ(1)=2;
76     lc(3)=ij+ncx; luv(3)=ij; typ(3)=-1;
77     end
78     % fprintf('Pb C: stop\n'); return;
79     else
80     %-- Not a corner point:
81     nbv=4;
82     %--- set 1rst neighbour i+1 :
83     lc(1)=ij+1; luv(1)=ij; typ(1)=2;
84     if i==nc,
85     %- edge point:
86     [I]=find( xx1==xx2(i+1,j,k) & yy1==yy2(i+1,j,k) ) ;
87     if length(I) ~= 1, fprintf('Pb A-i: stop\n'); return; end
88     lc(1)=I(1);
89     end
90     %--- set 2nd neighbour i-1 :
91     if i==1,
92     %- edge point:
93     [I]=find( xx3==xx1(ij) & yy3==yy1(ij) ) ;
94     if dbug>1, fprintf('Edge: i,j,k= %2i %2i %i ',i,j,k); end
95     if length(I) ~= 2, fprintf('Pb E-i: stop\n'); return; end
96     for l=1:length(I),
97     i2=I(l);
98     j2=fix((i2-1)/ncp); k2=1+fix(j2/ncp); j2=1+rem(j2,ncp);
99     i2=1+rem(i2-1,ncp);
100     if dbug>1, fprintf('; l= %i : i,j,k_2= %2i %2i %i ',l,i2,j2,k2); end
101     if k2 ~= k & i2==ncp,
102     lc(2)=nc+(k2-1)*nc+(j2-1)*ncx; luv(2)=lc(2); typ(2)=-2;
103     elseif k2 ~= k & j2==ncp,
104     lc(2)=i2+(k2-1)*nc+(nc-1)*ncx; luv(2)=lc(2); typ(2)=1;
105     end
106     end
107     if dbug>1, fprintf('\n'); end
108     else
109     lc(2)=ij-1; luv(2)=ij-1; typ(2)=-2;
110     end
111     %--- set 3rd neighbour j+1 :
112     lc(3)=ij+ncx; luv(3)=ij; typ(3)=-1;
113     if j==nc,
114     %- edge point:
115     [I]=find( xx1==xx2(i,j+1,k) & yy1==yy2(i,j+1,k) ) ;
116     if length(I) ~= 1, fprintf('Pb A-j: stop\n'); return; end
117     lc(3)=I(1);
118     end
119     %--- set 4th neighbour j-1 :
120     if j==1,
121     %- edge point:
122     [I]=find( xx3==xx1(ij) & yy3==yy1(ij) ) ;
123     if dbug>1, fprintf('Edge: i,j,k= %2i %2i %i ',i,j,k); end
124     if length(I) ~= 2, fprintf('Pb E-j: stop\n'); return; end
125     for l=1:length(I),
126     i2=I(l);
127     j2=fix((i2-1)/ncp); k2=1+fix(j2/ncp); j2=1+rem(j2,ncp);
128     i2=1+rem(i2-1,ncp);
129     if dbug>1, fprintf('; l= %i : i,j,k_2= %2i %2i %i ',l,i2,j2,k2); end
130     if k2 ~= k & j2==ncp,
131     lc(4)=i2+(k2-1)*nc+(nc-1)*ncx; luv(4)=lc(4); typ(4)=1;
132     elseif k2 ~= k & i2==ncp,
133     lc(4)=nc+(k2-1)*nc+(j2-1)*ncx; luv(4)=lc(4); typ(4)=-2;
134     end
135     end
136     if dbug>1, fprintf('\n'); end
137     else
138     lc(4)=ij-ncx; luv(4)=ij-ncx; typ(4)=1;
139     end
140     end
141     %- keep the untagged points in a temp list:
142     for l=1:nbv, if tagC(lc(l))==0,
143     n2p=n2p+1;
144     list_p(n2p)=ij; list_c(n2p)=lc(l);
145     listuv(n2p)=luv(l); list_t(n2p)=typ(l);
146     end; end;
147     end
148     if kgr == 1 & nRd==nRgr,
149     for n=1:n2p,
150     xloc=[xx1(list_p(n)) xx1(list_c(n))];
151     yloc=[yy1(list_p(n)) yy1(list_c(n))];
152     plot(xloc,yloc,'b-');
153     if n==1, hold on ; end
154     end
155     hold off ; grid
156     elseif kgr == 2 & nRd==nRgr,
157     xloc=xx1(list_c(1:n2p)); yloc=yy1(list_c(1:n2p));
158     plot(xloc,yloc,'ro'); hold on ;
159     xloc=xx1(list_p(1:n2p)); yloc=yy1(list_p(1:n2p));
160     plot(xloc,yloc,'bo');
161     hold off ; grid ;
162     end
163     if n2p == 0, doRd=0; else nRd=nRd+1 ; end
164     %- deal with non-single points:
165     nbp=0;
166     for n=1:n2p,
167     ij=list_c(n); dbl=0;
168     if n > 1, [I]=find(list_C(1:nbp,nRd)==ij); dbl=length(I); end
169     if dbl==0,
170     nbp=nbp+1;
171     list_C(nbp,nRd)=list_c(n);
172     list_P(nbp,nRd)=list_p(n);
173     listUV(nbp,nRd)=listuv(n);
174     list_T(nbp,nRd)=list_t(n);
175     tagC(ij,1)=nRd;
176     elseif dbl==1,
177     %- choose between 2 origin points: select the smaller Delta-X
178     xx=xx1(ij);
179     dx1=xx1(list_P(I,nRd)); dx1=rem(dx1-xx+540,360)-180;
180     dx2=xx1(list_p(n)); dx2=rem(dx2-xx+540,360)-180;
181     if abs(dx2) < abs(dx1),
182     %- replace I with n :
183     list_C(I,nRd)=list_c(n);
184     list_P(I,nRd)=list_p(n);
185     listUV(I,nRd)=listuv(n);
186     list_T(I,nRd)=list_t(n);
187     else
188     %- remove point n :
189     list_c(n)=0;
190     end
191     else
192     fprintf('Pb D: stop\n'); return;
193     end
194     end
195     if n2p > 0,
196     list_N(1,nRd)=nbp;
197     fprintf('nRd= %i ; nbp= %i \n',nRd,nbp);
198     end
199     %if nRd > 64, doRd=0; end
200     end
201    
202     if kgr == 3,
203     for n=2:nRd,
204     if n == nRd, linC='r-o'; else linC='b-'; end
205     for i=1:list_N(1,n),
206     xloc=[xx1(list_P(i,n)) xx1(list_C(i,n))];
207     yloc=[yy1(list_P(i,n)) yy1(list_C(i,n))];
208     plot(xloc,yloc,linC);
209     if i==1 & n==2, hold on ; end
210     end
211     end
212     hold off ; AA=axis;
213     AA(1)=max(-180,AA(1)); AA(2)=min(183,AA(2));
214     AA(3)=max(-90,AA(3)); AA(4)=min(90,AA(4));
215     axis(AA); grid
216     elseif kgr == 4,
217     for n=1:nRd,
218     if n == nRd, linC='ro'; else linC='ko'; end
219     nbp=list_N(1,n);
220     xloc=xx1(list_C(1:nbp,n)); yloc=yy1(list_C(1:nbp,n));
221     plot(xloc,yloc,linC);
222     if n==1, hold on ; end
223     end
224     xloc=xx1(list_P(1:nbp,nRd)); yloc=yy1(list_P(1:nbp,nRd));
225     plot(xloc,yloc,'bx');
226     hold off ; grid ;
227     end
228    
229     if kwr==1,
230     %-- copy to the right sixe arrays:
231     psiNx=max(list_N);
232     psiNy=nRd;
233     psi_N=zeros(1,psiNy);
234     psi_C=zeros(psiNx,psiNy); psi_P=zeros(psiNx,psiNy);
235     psiUV=zeros(psiNx,psiNy); psi_T=zeros(psiNx,psiNy);
236     psi_N=list_N(1,1:psiNy);
237     psi_C=list_C(1:psiNx,1:psiNy);
238     psi_P=list_P(1:psiNx,1:psiNy);
239     psiUV=listUV(1:psiNx,1:psiNy);
240     psi_T=list_T(1:psiNx,1:psiNy);
241     %- write to a file :
242     namf=['psiLine_N2S_cs',int2str(nc)];
243     fprintf(['write on file : ',namf,'.mat ']);
244     save(namf,'psi_N','psi_C','psi_P','psiUV','psi_T');
245     fprintf(' <== O.K. \n');
246     end

  ViewVC Help
Powered by ViewVC 1.1.22