/[MITgcm]/MITgcm_contrib/high_res_cube/eddy_flux/c20a/plot_c20a.m
ViewVC logotype

Annotation of /MITgcm_contrib/high_res_cube/eddy_flux/c20a/plot_c20a.m

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


Revision 1.2 - (hide annotations) (download)
Fri Feb 25 19:21:30 2005 UTC (20 years, 5 months ago) by edhill
Branch: MAIN
Changes since 1.1: +119 -17 lines
 o a simple and reasonably fast cube --> lat/lon converter

1 edhill 1.1 %=======================================================
2     %
3 edhill 1.2 % $Header: /u/gcmpack/MITgcm_contrib/high_res_cube/eddy_flux/c20a/plot_c20a.m,v 1.1 2005/02/24 15:38:35 edhill Exp $
4 edhill 1.1 %
5     % Ed Hill
6     %
7    
8     % The following are the MatLAB commands used to create the various
9     % plots related to eddy fluxes using average velocities and densities
10     % (called bouyancy or "b" in many of the variables) from Dimitris'
11     % "cube_20a" or "c20a" integration.
12    
13     % Groups of commands contained within the following "UNUSED" comments
14     % we not used in this analyssis. They are "left over" from the
15     % previous temperature-based calculations and have been kept for
16     % reference purposes only.
17    
18     %------- UNUSED -----------------------
19     % ...Commands...
20     %------- UNUSED -----------------------
21    
22    
23     % ssh eddy
24     % cd /r/r0/edhill/eddy_stats/eddy_flux_c20a
25    
26     % matlab -nojvm
27     % matlab -nojvm -nodisplay
28    
29     clear all
30     close all
31    
32    
33 edhill 1.2 %==================================================================
34     % Read the tile00?.mitgrid files
35 edhill 1.1 gvars = { 'XC','YC','DXF','DYF','RA','XG','YG','DXV', ...
36     'DYU','RAZ','DXC','DYC','RAW','RAS','DXG','DYG' };
37    
38 edhill 1.2 ne = 510;
39     nep1 = ne + 1;
40 edhill 1.1 iface = 1;
41     for iface = 1:6
42     fname = sprintf('grid/tile%03d.mitgrid', iface);
43     gid = fopen(fname, 'r', 'ieee-be');
44 edhill 1.2 tmp = reshape(fread(gid,inf,'real*8',0,'ieee-be'),[nep1,nep1,16]);
45 edhill 1.1 fclose(gid);
46     % surf(tmp(:,:,1)), view(2), shading interp
47     % for jj = 1:length(gvars)
48     for jj = 1:7
49     comm = sprintf('%s(:,:,%d) = tmp(:,:,%d);', ...
50     [gvars{jj}], iface, jj);
51     eval(comm);
52     end
53     end
54     % surf(XC(:,:,1)), view(2), shading interp
55     % subplot(2,1,1), a = [1:10]; surf(XC(a,a,1)), view(2)
56 edhill 1.2 % subplot(2,1,2), a = [(nep1-10):nep1]; surf(XC(a,a,1)), view(2)
57 edhill 1.1 % surf(YC(:,:,1)), view(2), shading interp
58     % surf(XG(:,:,1)), view(2), shading interp
59     % surf(YG(:,:,1)), view(2), shading interp
60 edhill 1.2 is = [1:ne];
61     vs = { 'XC','YC','DXF','DYF','RA' };
62     for i = 1:length(vs)
63     eval(sprintf('%s = %s(is,is,:);',vs{i},vs{i}));
64 edhill 1.1 end
65    
66 edhill 1.2 delR = [ ...
67     10.00, 10.00, 10.00, 10.00, 10.00, 10.00, 10.00, 10.01, ...
68     10.03, 10.11, 10.32, 10.80, 11.76, 13.42, 16.04 , 19.82, 24.85, ...
69     31.10, 38.42, 46.50, 55.00, 63.50, 71.58, 78.90, 85.15, 90.18, ...
70     93.96, 96.58, 98.25, 99.25,100.01,101.33,104.56,111.33,122.83, ...
71     139.09,158.94,180.83,203.55,226.50,249.50,272.50,295.50,318.50, ...
72     341.50,364.50,387.50,410.50,433.50,456.50 ];
73     R = cumsum(delR) - 0.5*delR;
74    
75     %==================================================================
76     % Project fields to lower-res 1-degree Lat-Lon and write
77     % as NetCDF for viewing with Ingrid
78     %
79 edhill 1.1 % !echo `ls -1 ave__92_99 | sed -e 's|\.ave||g' \
80     % | sed -e "s|^|'|g" | sed -e "s|$|',|g"`
81 edhill 1.2 fields_3d = { ...
82     'DRHODR', 'RHOANOSQ', 'RHOAnoma', 'SALT', ...
83     'SALTSQ', ...
84     'THETA', 'THETASQ', 'URHOMASS', ...
85 edhill 1.1 'USLTMASS', 'UTHMASS', 'UVEL', 'UVELMASS', 'UVELSQ', ...
86     'UV_VEL_Z', 'VRHOMASS', 'VSLTMASS', 'VTHMASS', ...
87     'VVEL', 'VVELMASS', 'VVELSQ', 'WRHOMASS', 'WSLTMASS', ...
88     'WTHMASS', 'WU_VEL', 'WVELMASS', 'WVELSQ', 'WV_VEL' };
89 edhill 1.2 fields_2d = { ...
90     'ETANSQ', 'SFLUX', 'SRELAX', 'TAUX', 'TAUY', 'TFLUX', ...
91     'TICE', 'TRELAX' };
92     % fid = fopen('ave__92_99/DRHODR.ave','r','ieee-be');
93     % tmp = fread(gid,inf,'real*4',0,'ieee-be');
94     % fclose(fid);
95    
96     ne = 510;
97     nf = 6;
98     nz = 50;
99     nslab = ne*ne*nf;
100     adir = 'ave__92_99';
101     lat = [-90:90];
102     lon = [0:360];
103    
104     % ! rm -f cube_20a_*.nc
105    
106     ifld = 5;
107     fields = union(fields_2d, fields_3d);
108     for ifld = 1:length(fields)
109    
110     ir = [ 1 2 3 5 10 15 20 35 30 35 40 50 ];
111     if ismember(fields{ifld},fields_2d)
112     ir = [ 1 ];
113     end
114    
115     disp([ ' ' fields{ifld} ' :' ]);
116    
117     nc = netcdf(['cube_20a_' fields{ifld} '.nc'], 'clobber');
118     nc.reference = [ 'Results from Dimitris Menemenlis''' ...
119     '"cube 20a" integrations' ];
120     nc.author = 'Ed Hill <eh3@mit.edu>';
121     nc.date = 'Feb 24, 2005';
122     nc('X') = length(lon);
123     nc('Y') = length(lat);
124     nc('Z') = length(ir);
125     nc{'X'} = 'X';
126     nc{'Y'} = 'Y';
127     nc{'Z'} = 'Z';
128     nc{'X'}.uniquename = 'X';
129     nc{'X'}.long_name = 'longitude';
130     nc{'X'}.gridtype = ncint(1);
131     nc{'X'}.units = 'degree_east';
132     nc{'Y'}.uniquename = 'Y';
133     nc{'Y'}.long_name = 'latitude';
134     nc{'Y'}.gridtype = ncint(0);
135     nc{'Y'}.units = 'degree_north';
136     nc{'Z'}.uniquename = 'Z';
137     nc{'Z'}.long_name = 'depth';
138     nc{'Z'}.gridtype = ncint(0);
139     nc{'Z'}.units = 'm';
140     nc{'X'}(:) = lon;
141     nc{'Y'}(:) = lat;
142     nc{'Z'}(:) = R(ir);
143    
144     llval = zeros(length(lon),length(lat),50);
145     fname = sprintf('%s/%s.ave',adir,fields{ifld});
146     fid = fopen(fname,'r','ieee-be');
147    
148     id = fields{ifld};
149     nc{ id } = { 'Z' 'Y' 'X' };
150     nc{ id }.missing_value = ncdouble(NaN);
151     nc{ id }.FillValue_ = ncdouble(0.0);
152    
153     ii = 1;
154     for ii = 1:length(ir)
155    
156     iz = ir(ii);
157     disp(sprintf(' iz = %g',iz));
158    
159     fseek(fid,nslab*4*(iz-1),'bof');
160     tmp = fread(fid,nslab,'real*4',0,'ieee-be');
161     tr = permute(reshape(tmp,[ 510 6 510 ]),[1 3 2]);
162     % surf(tr(:,:,1)), view(2), shading interp
163     xc360 = XC + 180;
164     trn = tr;
165     trn(find(tr == 0.0)) = NaN;
166     clear tmp tr
167     % v = sdac_regrid(xc360,YC,trn,lonm,latm);
168     v = ll_regrid(xc360,YC,trn,lon,lat);
169     % surf(lon,lat,v'), caxis([25 40]), view(2), shading interp, colorbar
170    
171     nc{ id }(ii,:,:) = permute(v,[2 1]);
172    
173     end
174    
175     fclose(fid);
176     nc = close(nc);
177 edhill 1.1
178 edhill 1.2 end
179 edhill 1.1
180 edhill 1.2 % !ncdump cube_20a.nc | more
181 edhill 1.1
182 edhill 1.2 % !scp Zonal_ave_92-03.nc ingrid.mit.edu:INGRID_PEOPLE/EH3/eddy_flux/cube_5/
183 edhill 1.1
184    
185    
186    
187     %------- UNUSED -----------------------
188    
189    
190     %=======================================================
191     % Compute u'b' and v'b'
192    
193     clear all; close all
194    
195     bfid = fopen( 'b-ave-1992-2003', 'r', 'ieee-be');
196     ufid = fopen( 'U_ave_1992-2003', 'r', 'ieee-be');
197     vfid = fopen( 'V_ave_1992-2003', 'r', 'ieee-be');
198     ubfid = fopen('ub-ave-1992-2003', 'r', 'ieee-be');
199     vbfid = fopen('vb-ave-1992-2003', 'r', 'ieee-be');
200    
201     afid = [ ufid vfid ];
202     av = [ 'ua'; 'va' ];
203    
204     nztot = 50;
205     ne = 510;
206     nps = 6 * ne * ne;
207     iz = 1;
208     for iz = 1:nztot
209    
210     disp(sprintf('iz = %d', iz));
211    
212     offset = (iz - 1)*nps*4;
213     id = 1;
214     for id = 1:length(afid)
215     fseek(afid(id), offset, 'bof');
216     % Convert the horrible JPL mdsio layout to six saner faces
217     ne = 510; tx = 85; ty = 85; nt = 216;
218     phi = unmangleJPL1( reshape(fread(afid(id), nps, 'real*4'), ...
219     tx*nt,ty), ne, tx );
220     comm = [ av(id,:) ' = phi;' ];
221     % disp(comm);
222     eval(comm);
223     end
224     clear phi
225     % surf(ua(:,:,1)), view(2), shading interp, axis equal
226     % surf(va(:,:,1)), view(2), shading interp, axis equal
227    
228     fseek(bfid, offset, 'bof');
229     b = reshape(fread(bfid, nps, 'real*4'),ne,ne,6);
230     % surf(b(:,:,1)), view(2), shading interp, axis equal
231     fseek(ubfid, offset, 'bof');
232     ub = reshape(fread(ubfid, nps, 'real*4'),ne,ne,6);
233     % surf(ub(:,:,1)), view(2), shading interp, axis equal
234     fseek(vbfid, offset, 'bof');
235     vb = reshape(fread(vbfid, nps, 'real*4'),ne,ne,6);
236     % surf(vb(:,:,1)), view(2), shading interp, axis equal
237    
238     % Fill the (ne+1)*(ne+1) grid with B values
239     ni = ne; nip1 = ni + 1;
240     nj = ne; njp1 = nj + 1;
241     bcubep1 = zeros( [ nip1 njp1 6 ] );
242     for icf = 1:6
243     bcubep1(2:nip1,2:njp1,icf) = b(:,:,icf);
244     end
245    
246     % Do the upwind-edge B exchanges
247     bcubep1(1,2:njp1,1) = b(ni:-1:1,nj,5); % -
248     bcubep1(2:nip1,1,1) = b(1:ni,nj,6); % -
249     bcubep1(1,2:njp1,2) = b(ni,1:nj,1); % -
250     bcubep1(2:nip1,1,2) = b(ni,nj:-1:1,6); % -
251     bcubep1(1,2:njp1,3) = b(ni:-1:1,nj,1); % -
252     bcubep1(2:nip1,1,3) = b(1:ni,nj,2); % -
253     bcubep1(1,2:njp1,4) = b(ni,1:nj,3); % -
254     bcubep1(2:nip1,1,4) = b(ni,nj:-1:1,2); % -
255     bcubep1(1,2:njp1,5) = b(ni:-1:1,nj,3); % -
256     bcubep1(2:nip1,1,5) = b(1:ni,nj,4); % -
257     bcubep1(1,2:njp1,6) = b(ni,1:nj,5); % -
258     bcubep1(2:nip1,1,6) = b(ni,nj:-1:1,4); % -
259    
260     % Get B values on the U grid points
261     masku = 1.0 - abs(diff(double(bcubep1(2:nip1,:,:) == 0.0),1,2));
262     diffu = 0.5*diff(bcubep1(2:nip1,:,:),1,2);
263     bonu = b(:,:,:) + masku.*diffu;
264    
265     % Get B values on the V grid points
266     maskv = 1.0 - abs(diff(double(bcubep1(:,2:nip1,:) == 0.0),1,1));
267     diffv = 0.5*diff(bcubep1(:,2:nip1,:),1,1);
268     bonv = b(:,:,:) + maskv.*diffv;
269    
270     % Compute u'B' and v'B'
271     upbp = ub - ua .* bonu;
272     vpbp = vb - va .* bonv;
273    
274     % Write the results
275     ubid = fopen('upbp-1992-2003', 'a', 'ieee-be');
276     vbid = fopen('vpbp-1992-2003', 'a', 'ieee-be');
277     fwrite(ubid, upbp, 'real*4');
278     fwrite(vbid, vpbp, 'real*4');
279     fclose(ubid);
280     fclose(vbid);
281    
282     end
283    
284     %=======================================================
285     % Plot u'B' and v'B' on a Lat-Lon grid
286    
287     clear all, close all
288     tx = 85;
289     ty = 85;
290     nt = 216;
291     cx = 510;
292     cy = 510;
293     nz = 1;
294     ne = 510;
295     nps = ne * ne * 6;
296    
297     % XCS YCS XGS YGS
298     fnam = [ 'XC' ; 'YC'; 'XG'; 'YG' ];
299     for in = 1:4
300     uid = fopen([fnam(in,:) '.data'], 'r', 'ieee-be');
301     phi = unmangleJPL1( reshape(fread(uid, nps, 'real*4'), ...
302     tx*nt,ty), ne, tx );
303     fclose(uid);
304     eval([lower(fnam(in,:)) ' = phi;']);
305     a = zeros(6*510,510);
306     for i = 1:6
307     xb = (i-1)*510 + 1; xe = xb + 510 - 1;
308     yb = 1; ye = 510;
309     a(xb:xe,yb:ye) = phi(:,:,i);
310     end
311     eval([lower(fnam(in,:)) 's = a;'])
312     end
313     xcs = xcs + -360.0*(xcs > 180.0);
314     xgs = xgs + -360.0*(xgs > 180.0);
315     clear phi a
316    
317     %=======================================================
318     % Calculation rotations for Lat-Lon vector alignment
319    
320     ne = 510;
321     n1 = ne - 1;
322     dux = zeros(size(xg));
323     duy = zeros(size(xg));
324     dvx = zeros(size(xg));
325     dvy = zeros(size(xg));
326     dux(1:n1,:,:) = diff(xg,1,1); dux(ne,:,:) = dux(n1,:,:);
327     dvx(:,1:n1,:) = diff(xg,1,2); dvx(:,ne,:) = dvx(:,n1,:);
328     duy(1:n1,:,:) = diff(yg,1,1); duy(ne,:,:) = duy(n1,:,:);
329     dvy(:,1:n1,:) = diff(yg,1,2); dvy(:,ne,:) = dvy(:,n1,:);
330     dux = dux + 360*double(dux < 180);
331     dux = dux - 360*double(dux > 180); % [ min(min(dux)) max(max(dux)) ]
332     duy = duy + 360*double(duy < 180);
333     duy = duy - 360*double(duy > 180); % [ min(min(duy)) max(max(duy)) ]
334     dvx = dvx + 360*double(dvx < 180);
335     dvx = dvx - 360*double(dvx > 180); % [ min(min(dvx)) max(max(dvx)) ]
336     dvy = dvy + 360*double(dvy < 180);
337     dvy = dvy - 360*double(dvy > 180); % [ min(min(dvy)) max(max(dvy)) ]
338     llux = dux ./ sqrt(dux.^2 + duy.^2);
339     lluy = duy ./ sqrt(dux.^2 + duy.^2);
340     llvx = dvx ./ sqrt(dvx.^2 + dvy.^2);
341     llvy = dvy ./ sqrt(dvx.^2 + dvy.^2);
342    
343     % surf(xg(:,:,2)), view(2), axis equal, shading interp, colorbar
344     % surf(dux(:,:,2)), view(2), axis equal, shading interp, colorbar
345     % surf(duy(:,:,2)), view(2), axis equal, shading interp, colorbar
346    
347     delR = [
348     10.00, 10.00, 10.00, 10.00, 10.00, 10.00, 10.00, 10.01, ...
349     10.03, 10.11, 10.32, 10.80, 11.76, 13.42, 16.04 , 19.82, 24.85, ...
350     31.10, 38.42, 46.50, 55.00, 63.50, 71.58, 78.90, 85.15, 90.18, ...
351     93.96, 96.58, 98.25, 99.25,100.01,101.33,104.56,111.33,122.83, ...
352     139.09,158.94,180.83,203.55,226.50,249.50,272.50,295.50,318.50, ...
353     341.50,364.50,387.50,410.50,433.50,456.50 ];
354     R = cumsum(delR) - 0.5*delR;
355     Rmid = R(1:(length(R)-1)) + 0.5*diff(R);
356    
357     ne = 510;
358     nps = ne * ne * 6;
359    
360     % u'T' and v'T'
361     upbpid = fopen('upbp-1992-2003', 'r', 'ieee-be');
362     vpbpid = fopen('vpbp-1992-2003', 'r', 'ieee-be');
363    
364     clow = [ -1.0 ]; % -0.2;
365     chigh = [ 1.0 ]; % 0.2;
366    
367     iz = 1;
368     iz = 15;
369     iz = 25;
370     for iz = [ 1:50 ]
371     disp(sprintf('iz = %d',iz));
372     % iz = 21;
373     offset = (iz - 1)*nps*4;
374     fseek(upbpid, offset, 'bof');
375     fseek(vpbpid, offset, 'bof');
376     upbp = reshape(fread(upbpid, nps, 'real*4'),ne,ne,6);
377     vpbp = reshape(fread(vpbpid, nps, 'real*4'),ne,ne,6);
378    
379    
380     % Look at u'b' and v'b' on ll coords
381     % surf(xg(:,:,2), yg(:,:,2), upbp(:,:,2)), view(2), shading interp, colorbar
382     % surf(xg(:,:,3), yg(:,:,3), upbp(:,:,3)), view(2), shading interp, colorbar
383     % surf(xg(:,:,4), yg(:,:,4), upbp(:,:,4)), view(2), shading interp, colorbar
384     llupbp = upbp .* llux + vpbp .* llvx;
385     llvpbp = upbp .* lluy + vpbp .* llvy;
386     % surf(xg(:,:,2), yg(:,:,2), lluptp(:,:,2)), view(2), shading interp, colorbar
387     % surf(xg(:,:,2), yg(:,:,2), llvptp(:,:,2)), view(2), shading
388     % interp, colorbar
389    
390     % save for later use
391     upbpllid = fopen('upbpll-1992-2003', 'a', 'ieee-be');
392     vpbpllid = fopen('vpbpll-1992-2003', 'a', 'ieee-be');
393     fwrite(upbpllid, llupbp, 'real*4');
394     fwrite(vpbpllid, llvpbp, 'real*4');
395     fclose(upbpllid);
396     fclose(vpbpllid);
397    
398     % Look at lat-lon projections of u'T' and v'T'
399     pllupbp = zeros(6*510,510);
400     pllvpbp = zeros(6*510,510);
401     for i = 1:6
402     xb = (i-1)*510 + 1; xe = xb + 510 - 1;
403     yb = 1; ye = 510;
404     pllupbp(xb:xe,yb:ye) = llupbp(:,:,i);
405     pllvpbp(xb:xe,yb:ye) = llvpbp(:,:,i);
406     end
407     shift=-1;
408     grph_CS(sq(pllupbp),xcs,ycs,xgs,ygs,[clow],[chigh],shift);
409     title([ 'u''B'' at ' ...
410     sprintf('%g',R(iz)) ...
411     'm depth on the 510x510x6 cubesphere for 1992--2003 ["cube5"]']);
412     % print('-painters', '-dpng', '-r650', ...
413     % ['upTp_' sprintf('%02d_%07.1f',iz,R(iz)) 'm_650.png'])
414     print('-painters', '-dpng', '-r150', ...
415     ['upBp_' sprintf('%02d_%07.1f',iz,R(iz)) 'm_150.png'])
416    
417     grph_CS(sq(pllvpbp),xcs,ycs,xgs,ygs,[clow],[chigh],shift);
418     title([ 'v''B'' at ' ...
419     sprintf('%g',R(iz)) ...
420     'm depth on the 510x510x6 cubesphere for 1992--2003 ["cube5"]']);
421     % print('-painters', '-dpng', '-r650', ...
422     % ['vpTp_' sprintf('%02d_%07.1f',iz,R(iz)) 'm_650.png'])
423     print('-painters', '-dpng', '-r150', ...
424     ['vpBp_' sprintf('%02d_%07.1f',iz,R(iz)) 'm_150.png'])
425    
426     end
427    
428    
429     %=======================================================
430     % Calculate : (v'B')/(dB/dz)
431    
432     bid = fopen( 'b-ave-1992-2003', 'r', 'ieee-be');
433     bm1id = fopen('bm1-ave-1992-2003', 'r', 'ieee-be');
434    
435     iz = 1;
436     offset = (iz - 1)*nps*4;
437     fseek(upbpid, offset, 'bof');
438     fseek(vpbpid, offset, 'bof');
439     upbp = reshape(fread(upbpid, nps, 'real*4'),ne,ne,6);
440     vpbp = reshape(fread(vpbpid, nps, 'real*4'),ne,ne,6);
441     llvpbp = upbp .* lluy + vpbp .* llvy;
442     llvpbpm1 = llvpbp;
443     ne = 510; tx = 85; ty = 85; nt = 216;
444    
445     iz = 2;
446     iz = 15;
447     iz = 16;
448     for iz = 2:50,
449     izm1 = iz - 1;
450     disp(sprintf('iz = %d',iz));
451     % iz = 21;
452     offset = (iz - 1)*nps*4;
453     fseek(upbpid, offset, 'bof');
454     fseek(vpbpid, offset, 'bof');
455     upbp = reshape(fread(upbpid, nps, 'real*4'),ne,ne,6);
456     vpbp = reshape(fread(vpbpid, nps, 'real*4'),ne,ne,6);
457    
458     % v'B' on ll coords
459     llvpbp = upbp .* lluy + vpbp .* llvy;
460    
461     % B and Bm1
462     fseek(bid, offset, 'bof');
463     b = reshape(fread(bid, nps, 'real*4'),ne,ne,6);
464     offm1 = (iz - 2)*nps*4;
465     fseek(bm1id, offm1, 'bof');
466     bm1 =reshape(fread(bm1id, nps, 'real*4'),ne,ne,6);
467    
468     % ( v'B' )/( dB/dz )
469     dbdz = (b - bm1)/(R(iz) - R(izm1));
470     ind0 = find(dbdz==0.0);
471     dbdz(ind0) = 1.0;
472     rdbdz = 1.0 / (dbdz);
473     rdbdz(ind0) = 0.0;
474     vpbpdbdz = 0.5*(llvpbp + llvpbpm1) .* rdbdz;
475    
476     % Write the results
477     mid = fopen('vpbpdbdz-1992-2003', 'a', 'ieee-be');
478     fwrite(mid, vpbpdbdz, 'real*4');
479     fclose(mid);
480    
481     % Plot results
482     clow = [ -10 ]; % -20;
483     chigh = [ 10 ]; % 20;
484     ll = zeros(6*510,510);
485     for i = 1:6
486     xb = (i-1)*510 + 1; xe = xb + 510 - 1;
487     yb = 1; ye = 510;
488     ll(xb:xe,yb:ye) = vpbpdbdz(:,:,i);
489     end
490     shift=-1;
491     grph_CS(sq(ll),xcs,ycs,xgs,ygs,[clow],[chigh],shift);
492     title([ '(v''B'')/(dB/dz) at ' ...
493     sprintf('%g',Rmid(iz)) ...
494     'm depth on the 510x510x6 cubesphere for 1992--2003 ["cube5"]']);
495     % print('-painters', '-dpng', '-r650', ...
496     % ['vpTpdTdz_' sprintf('%02d_%07.1f',iz,Rmid(iz)) 'm_650.png'])
497     print('-painters', '-dpng', '-r150', ...
498     ['vpBpdBdz_' sprintf('%02d_%07.1f',iz,Rmid(iz)) 'm_150.png'])
499    
500    
501     % Stress: calc, write, and plot
502     stress = 1000 * (2*pi/(24*3600)*sin(pi*yc/180));
503     stress = stress .* vpbpdbdz;
504     sid = fopen('stress-b-1992-2003', 'a', 'ieee-be');
505     fwrite(sid, stress, 'real*4');
506     fclose(sid);
507     clow = []; % [ -1 ];
508     chigh = []; % [ 1 ];
509     ll = zeros(6*510,510);
510     for i = 1:6
511     xb = (i-1)*510 + 1; xe = xb + 510 - 1;
512     yb = 1; ye = 510;
513     ll(xb:xe,yb:ye) = stress(:,:,i);
514     end
515     shift=-1;
516     grph_CS(sq(ll),xcs,ycs,xgs,ygs,[clow],[chigh],shift);
517     title([ 'Stress at ' ...
518     sprintf('%g',Rmid(iz)) ...
519     'm depth on the 510x510x6 cubesphere for 1994--2003 ["cube5"]']);
520     % print('-painters', '-dpng', '-r650', ...
521     % ['stress_' sprintf('%02d_%07.1f',iz,Rmid(iz)) 'm_650.png'])
522     print('-painters', '-dpng', '-r150', ...
523     ['stress_' sprintf('%02d_%07.1f',iz,Rmid(iz)) 'm_150.png'])
524    
525    
526     % Next level
527     llvpbpm1 = llvpbp;
528    
529     end
530    
531    
532     %=======================================================
533     % Zonally average vpbpdbdz and stress
534    
535     clear all ; close all
536     tx = 85;
537     ty = 85;
538     nt = 216;
539     cx = 510;
540     cy = 510;
541     nz = 1;
542     ne = 510;
543     nps = ne * ne * 6;
544    
545     % XCS YCS XGS YGS
546     fnam = [ 'XC' ; 'YC'; 'XG'; 'YG' ];
547     for in = 1:4
548     uid = fopen([fnam(in,:) '.data'], 'r', 'ieee-be');
549     phi = unmangleJPL1( reshape(fread(uid, nps, 'real*4'), ...
550     tx*nt,ty), ne, tx );
551     fclose(uid);
552     eval([lower(fnam(in,:)) ' = phi;']);
553     a = zeros(6*510,510);
554     for i = 1:6
555     xb = (i-1)*510 + 1; xe = xb + 510 - 1;
556     yb = 1; ye = 510;
557     a(xb:xe,yb:ye) = phi(:,:,i);
558     end
559     eval([lower(fnam(in,:)) 's = a;'])
560     end
561     xcs = xcs + -360.0*(xcs > 180.0);
562     xgs = xgs + -360.0*(xgs > 180.0);
563     clear phi a
564    
565     nz = 1;
566     nr = 50;
567     nrm1 = nr - 1;
568     nlat = 181; nlatm1 = nlat - 1;
569    
570     hvals = linspace(-90,90,nlat);
571     % save indicies for zonal averages
572     i = 2;
573     for i = 2:nlat
574     inds = find(hvals(i-1)<yg & yg<hvals(i));
575     comm = sprintf('inds%04d = uint32(inds);',i-1);
576     eval(comm);
577     end
578    
579     % Zonally average vpbp
580     acc = zeros(nlatm1, nrm1);
581     num = zeros(size(acc));
582     ne = 510;
583     nps = ne * ne * 6;
584     zid = fopen('vpbpll-1992-2003', 'r', 'ieee-be');
585     iz = 1;
586     for iz = 1:50,
587     disp(sprintf('iz = %d',iz));
588     offset = (iz - 1)*nps*4;
589     fseek(zid, offset, 'bof');
590     vpbp = reshape(fread(zid, nps, 'real*4'),ne,ne,6);
591    
592     for jj = 1:nlatm1
593     eval( sprintf('clear inds; inds = inds%04d;',jj) );
594     tmp = vpbp(inds);
595     nzinds = find(tmp ~= 0.0);
596     num(jj,iz) = length(nzinds);
597     acc(jj,iz) = sum(tmp(nzinds));
598     end
599     end
600     fclose(zid);
601     zvpbp = acc ./ num;
602     % save zvpbp zvpbp
603    
604     % zonally average vpbpdbdz
605     acc = zeros(nlatm1, nrm1);
606     num = zeros(size(acc));
607     ne = 510;
608     nps = ne * ne * 6;
609     zid = fopen('vpbpdbdz-1992-2003', 'r', 'ieee-be');
610     iz = 1;
611     for iz = 1:49,
612     disp(sprintf('iz = %d',iz));
613     offset = (iz - 1)*nps*4;
614     fseek(zid, offset, 'bof');
615     vpbpdbdz = reshape(fread(zid, nps, 'real*4'),ne,ne,6);
616    
617     for jj = 1:nlatm1
618     eval( sprintf('clear inds; inds = inds%04d;',jj) );
619     tmp = vpbpdbdz(inds);
620     nzinds = find(tmp ~= 0.0);
621     num(jj,iz) = length(nzinds);
622     acc(jj,iz) = sum(tmp(nzinds));
623     end
624     end
625     fclose(zid);
626     zvpbpdbdz = acc ./ num;
627     % save zvpbpdbdz zvpbpdbdz
628    
629     % zonally average stress
630     acc = zeros(nlatm1, nrm1);
631     num = zeros(size(acc));
632     ne = 510;
633     nps = ne * ne * 6;
634     zid = fopen('stress-b-1992-2003', 'r', 'ieee-be');
635     iz = 1;
636     for iz = 1:49,
637     disp(sprintf('iz = %d',iz));
638     offset = (iz - 1)*nps*4;
639     fseek(zid, offset, 'bof');
640     stress = reshape(fread(zid, nps, 'real*4'),ne,ne,6);
641     jj = 22;
642     for jj = 1:nlatm1
643     eval( sprintf('clear inds; inds = inds%04d;',jj) );
644     tmp = stress(inds);
645     nzinds = find(tmp ~= 0.0);
646     num(jj,iz) = length(nzinds);
647     acc(jj,iz) = sum(tmp(nzinds));
648     end
649     end
650     fclose(zid);
651     stress = acc ./ num;
652     % save stress stress
653    
654     %------- UNUSED -----------------------
655     % zonally average u
656     acc = zeros(nlatm1, 50);
657     num = zeros(size(acc));
658     ne = 510;
659     nps = ne * ne * 6;
660     zid = fopen('Ull_ave_1994--2003', 'r', 'ieee-be');
661     iz = 1;
662     for iz = 1:50,
663     disp(sprintf('iz = %d',iz));
664     offset = (iz - 1)*nps*4;
665     fseek(zid, offset, 'bof');
666     U = reshape(fread(zid, nps, 'real*4'),ne,ne,6);
667     jj = 22;
668     for jj = 1:nlatm1
669     eval( sprintf('clear inds; inds = inds%04d;',jj) );
670     tmp = U(inds);
671     nzinds = find(tmp ~= 0.0);
672     num(jj,iz) = length(nzinds);
673     acc(jj,iz) = sum(tmp(nzinds));
674     end
675     end
676     fclose(zid);
677     zonalu = acc ./ num;
678     % save zonalu zonalu
679     %------- UNUSED -----------------------
680    
681     % zonally average B
682     acc = zeros(nlatm1, 50);
683     num = zeros(size(acc));
684     ne = 510;
685     nps = ne * ne * 6;
686     zid = fopen('b-ave-1992-2003', 'r', 'ieee-be');
687     iz = 1;
688     for iz = 1:50,
689     disp(sprintf('iz = %d',iz));
690     offset = (iz - 1)*nps*4;
691     fseek(zid, offset, 'bof');
692     B = reshape(fread(zid, nps, 'real*4'),ne,ne,6);
693     jj = 22;
694     for jj = 1:nlatm1
695     eval( sprintf('clear inds; inds = inds%04d;',jj) );
696     tmp = B(inds);
697     nzinds = find(tmp ~= 0.0);
698     num(jj,iz) = length(nzinds);
699     acc(jj,iz) = sum(tmp(nzinds));
700     end
701     end
702     fclose(zid);
703     zonalB = acc ./ num;
704     % save zonalB zonalB
705    
706     %------- UNUSED -----------------------
707     % Vertical gradient of zonally averaged B (dB/dz)
708     load zonalT
709     load allR
710     nz = size(zonalT,1);
711     dzaTdz = zeros(nz,size(zonalT,2)-1);
712     for iz = 2:50
713     dzaTdz(:,iz-1) = (zonalT(:,iz) - zonalT(:,iz-1)) ./ (R(iz) - R(iz-1));
714     end
715     % save dzaTdz dzaTdz
716     %------- UNUSED -----------------------
717    
718    
719     delR = [
720     10.00, 10.00, 10.00, 10.00, 10.00, 10.00, 10.00, 10.01, ...
721     10.03, 10.11, 10.32, 10.80, 11.76, 13.42, 16.04 , 19.82, 24.85, ...
722     31.10, 38.42, 46.50, 55.00, 63.50, 71.58, 78.90, 85.15, 90.18, ...
723     93.96, 96.58, 98.25, 99.25,100.01,101.33,104.56,111.33,122.83, ...
724     139.09,158.94,180.83,203.55,226.50,249.50,272.50,295.50,318.50, ...
725     341.50,364.50,387.50,410.50,433.50,456.50 ];
726     R = cumsum(delR) - 0.5*delR;
727     Rmid = R(1:(length(R)-1)) + 0.5*diff(R);
728     hmid = [ -89.5:1:89.5 ];
729     % save allR delR R Rmid hmid
730    
731     %------- UNUSED -----------------------
732     surf(hmid,-Rmid,num'), view(2), shading interp, colorbar
733     caxis([ -1 1 ])
734     caxis('manual')
735     surf(hmid,-Rmid,mu'), view(2), shading interp, colorbar
736     title('Zonally Averaged Stress for 1994--2003 [cube5]')
737     xlabel('Latitude [deg]')
738     zlabel('Depth [m]')
739     axis([ -90 90 -5500 200 ])
740    
741     % print -dps -painters t_001.ps
742     print('-painters', '-dpng', '-r150', 'za_str_94--03_r150.png')
743     print('-painters', '-dpng', '-r650', 'za_str_94--03_r650.png')
744    
745     % plot and print bar(U)
746     load zonalu
747     surf(hmid,-R,zonalu'), view(2), shading interp, colorbar
748     contour(hmid,-R,mu'), colorbar
749     ac = 25;
750     ac = [-0.3:.005:.3];
751     [c,h] = contourf(hmid,-R,mu',ac); colorbar
752     hold on, contour(hmid,-R,mu',ac), hold off
753     title('Zonally Averaged U for 1994--2003 [cube5]')
754     xlabel('Latitude [deg]')
755     zlabel('Depth [m]')
756     axis([ -90 90 -6000 200 ])
757     print('-painters', '-dpng', '-r150', 'Ull_94--03_r150.png')
758     print('-painters', '-dpng', '-r650', 'Ull_94--03_r650.png')
759     %------- UNUSED -----------------------
760    
761    
762     %=======================================================
763     % Put it all in NetCDF
764    
765     clear all ; close all
766     load allR
767     % load zonalu
768     load zvpbpdbdz
769     load stress
770     load zonalB
771     % load dzaTdz
772     load zvpbp
773    
774     % id0 = 'U';
775     id1 = 'zvpbpdbdz';
776     id2 = 'stress';
777     id3 = 'B';
778     % id4 = 'dzaTdz';
779     id5 = 'zvpbp';
780     % units0 = 'm/s';
781     units1 = 'm^2/s';
782     units2 = 'N/m^2';
783     units3 = 'kg/m^3';
784     % units4 = 'deg C/m';
785     units5 = 'm*kg/m^3';
786    
787     !rm -f Zonal_ave_92-03.nc
788     nc = netcdf('Zonal_ave_92-03.nc', 'clobber');
789     nc.reference = 'Zonal averages from Dimitris "cube 5" integrations';
790     nc.author = 'Ed Hill <eh3@mit.edu>';
791     nc.date = 'Aug 11, 2004';
792    
793     nc('X') = length(hmid);
794     nc('Y') = length(R);
795     nc('Ym1') = length(R) - 1;
796     nc{'X'} = 'X';
797     nc{'Y'} = 'Y';
798     nc{'Ym1'} = 'Ym1';
799     % nc{ id0 } = { 'Y', 'X' };
800     nc{ id1 } = { 'Ym1', 'X' };
801     nc{ id2 } = { 'Ym1', 'X' };
802     nc{ id3 } = { 'Y', 'X' };
803     % nc{ id4 } = { 'Ym1', 'X' };
804     nc{ id5 } = { 'Y', 'X' };
805     nc{'X'}.uniquename = 'X';
806     nc{'X'}.long_name = 'latitude';
807     nc{'X'}.gridtype = ncint(0);
808     nc{'X'}.units = 'degree_north';
809     nc{'Y'}.uniquename = 'Y';
810     nc{'Y'}.long_name = 'depth';
811     nc{'Y'}.gridtype = ncint(1);
812     nc{'Y'}.units = 'm';
813     nc{'Ym1'}.uniquename = 'Ym1';
814     nc{'Ym1'}.long_name = 'depth';
815     nc{'Ym1'}.gridtype = ncint(1);
816     nc{'Ym1'}.units = 'm';
817     % nc{ id0 }.units = units0;
818     % nc{ id0 }.long_name = id0;
819     % nc{ id0 }.missing_value = ncdouble(NaN);
820     % nc{ id0 }.FillValue_ = ncdouble(0.);
821     nc{ id1 }.units = units1;
822     nc{ id1 }.long_name = id1;
823     nc{ id1 }.missing_value = ncdouble(NaN);
824     nc{ id1 }.FillValue_ = ncdouble(0.);
825     nc{ id2 }.units = units2;
826     nc{ id2 }.long_name = id2;
827     nc{ id2 }.missing_value = ncdouble(NaN);
828     nc{ id2 }.FillValue_ = ncdouble(0.);
829     nc{ id3 }.units = units3;
830     nc{ id3 }.long_name = id3;
831     nc{ id3 }.missing_value = ncdouble(NaN);
832     nc{ id3 }.FillValue_ = ncdouble(0.);
833     % nc{ id4 }.units = units4;
834     % nc{ id4 }.long_name = id4;
835     % nc{ id4 }.missing_value = ncdouble(NaN);
836     % nc{ id4 }.FillValue_ = ncdouble(0.);
837     nc{ id5 }.units = units5;
838     nc{ id5 }.long_name = id5;
839     nc{ id5 }.missing_value = ncdouble(NaN);
840     nc{ id5 }.FillValue_ = ncdouble(0.);
841     nc{'X'}(:) = hmid;
842     nc{'Y'}(:) = -R;
843     nc{'Ym1'}(:) = -Rmid;
844     % nc{ id0 }(:) = permute(zonalu,[ 2 1 ]);
845     nc{ id1 }(:) = permute(zvpbpdbdz,[ 2 1 ]);
846     nc{ id2 }(:) = permute(stress,[ 2 1 ]);
847     nc{ id3 }(:) = permute(zonalB,[ 2 1 ]);
848     % nc{ id4 }(:) = permute(dzaTdz,[ 2 1 ]);
849     nc{ id5 }(:) = permute(zvpbp,[ 2 1 ]);
850     nc = close(nc);
851    
852     % !scp Zonal_ave_92-03.nc ingrid.mit.edu:INGRID_PEOPLE/EH3/eddy_flux/cube_5/

  ViewVC Help
Powered by ViewVC 1.1.22