Back to home page

MITgcm

 
 

    


Warning, /verification/tutorial_global_oce_latlon/diags_matlab/mit_plotstreamfunctions.m is written in an unsupported language. File is not indexed.

view on githubraw file Latest commit add29e06 on 2018-01-31 20:35:05 UTC
051ee7f715 Jean*0001 % m-file: mit_plotstreamfunctions.m
                0002 
                0003 % select timestep
                0004 k=kmax;
                0005 
                0006 if iscell(u)
                0007   uk = u{k};
                0008 else
                0009   uk = squeeze(u(:,:,:,k));
                0010 end
                0011 if iscell(v)
                0012   vk = v{k};
                0013 else
                0014   vk = squeeze(v(:,:,:,k));
                0015 end
                0016 
                0017 
                0018 addlayer = 1;
                0019 
                0020 clear global_psi atlantic_psi baro_psi
                0021 % global overturning
                0022 global_psi = mit_overturning(vk,grd.hfacs,grd.dxg,grd.dz,addlayer);
                0023 
                0024 % atlantic overturning
                0025 atlantic_psi = mit_overturning(vk,grd.atlantic_hfacs,grd.dxg,grd.dz,addlayer);
                0026 % pacific overturning
                0027 pacific_psi = mit_overturning(vk,grd.pacific_hfacs,grd.dxg,grd.dz,addlayer);
                0028 
                0029 clear addlayer
                0030 
                0031 % global barotropic stream function
                0032 baro_psi = mit_barostream(uk,grd.umask,grd.dyg,grd.dz);
                0033 
                0034 % plot stream functions
                0035 figure('PaperPosition',[0.31 0.25 10.5 7.88],'PaperOrientation','landscape')
                0036 clear sh
                0037 sh(1) = subplot(2,2,1);
                0038 otlev = [-60:2:60];
                0039 contourf(grd.latg,-grd.zgpsi,global_psi'*1e-6,otlev); 
                0040 hold on; 
                0041 [cs h1] = contour(grd.latg,-grd.zgpsi,global_psi'*1e-6,[0 0]); 
                0042 clh1 = clabel(cs);
                0043 hold off
                0044 caxis([-1 1]*max(abs(global_psi(:)))*1.e-6); colorbar('h')
                0045 psimin = min(min(global_psi(:,5:end)));
                0046 [iy iz] = find(abs(global_psi(:,:)-psimin)<=1e-4);
                0047 text(grd.latg(iy),-grd.zgpsi(iz), ...
                0048      ['\leftarrow ' num2str(psimin*1e-6,'%5.1f')], ...
                0049      'horizontalalignment','left')
                0050 title('global overturning streamfunction [Sv]')
                0051 sh(2) = subplot(2,2,2);
                0052 contourf(grd.latg,-grd.zgpsi,atlantic_psi'*1e-6,otlev); 
                0053 hold on; 
                0054 [cs h2] = contour(grd.latg,-grd.zgpsi,atlantic_psi'*1e-6,[0 0]); 
                0055 clh2 = clabel(cs);
                0056 hold off
                0057 caxis([-1 1]*max(abs(atlantic_psi(:)))*1.e-6); colorbar('h');
                0058 psimax = max(atlantic_psi(13,5:end));
                0059 iz = find(abs(atlantic_psi(13,:)-psimax)<=1e-4);
                0060 text(grd.latg(13),-grd.zgpsi(iz), ...
                0061      [num2str(psimax*1e-6,'%5.1f') ' \rightarrow'], ...
                0062      'horizontalalignment','right')
                0063 psimin = min(min(atlantic_psi(1:35,5:end)));
                0064 [iymin,izmin] = find(abs(atlantic_psi(:,:)-psimin)<=1e-4);
                0065 text(grd.latg(iymin),-grd.zgpsi(izmin), ...
                0066      [num2str(psimin*1e-6,'%5.1f') ' \rightarrow'], ...
                0067      'horizontalalignment','right')
                0068 title('atlantic overturning streamfunction [Sv]')
                0069 %
                0070 sh(3) = subplot(2,2,3);
                0071 contourf(grd.latg,-grd.zgpsi,pacific_psi'*1e-6,otlev); 
                0072 hold on; 
                0073 [cs h3] = contour(grd.latg,-grd.zgpsi,pacific_psi'*1e-6,[0 0]); 
                0074 clh3 = clabel(cs);
                0075 hold off
                0076 caxis([-1 1]*max(abs(pacific_psi(:)))*1.e-6); colorbar('h');
                0077 title('pacific overturning streamfunction [Sv]')
                0078 if ~isempty([h1;h2;h3])
                0079   set([h1;h2;h3],'LineWidth',2,'EdgeColor','k');
                0080 end
                0081 clh = [clh1;clh2;clh3];
                0082 if ~isempty(clh)
                0083   set(clh(2:2:end),'FontSize',8);
                0084 end
                0085 % $$$ [cs h] = contourf(grd.long,grd.latg,baro_psi'*1e-6,20); 
                0086 % $$$ if ~isempty(h);
                0087 % $$$   set(h,'edgecolor','none'); 
                0088 % $$$ end; 
                0089 % $$$ axis image; 
                0090 % $$$ caxis([-1 1]*max(abs(baro_psi(:)))*1.e-6); colorbar('h');
                0091 % $$$ title('global barotropic stream function [Sv]')
                0092 bstlev = [-200:20:200];
                0093 
                0094 sh(4) = subplot(2,2,4);
                0095 imagesc(grd.long,grd.latg,baro_psi'*1e-6); 
                0096 hold on;
                0097 [cs h ]=contour(grd.long,grd.latg,baro_psi'*1e-6,bstlev);
                0098 set(h,'edgecolor','k')
                0099 if ~isempty(h); 
                0100   clh = clabel(cs,h); 
                0101   set(clh,'Fontsize',8);
                0102 end
                0103 hold off
                0104 axis image, axis xy; 
                0105 caxis([-1 1]*max(abs(baro_psi(:)))*1.e-6); colorbar('h');
                0106 title('global barotropic stream function [Sv]')
                0107 suptitle(['experiment ' dname ', timestep = ' num2str(timesteps(k)) ...
                0108           ', ' tuname ' = ' num2str(tim(k))])
                0109 set(sh,'layer','top')
                0110 
                0111 clear addlayer