/[MITgcm]/MITgcm_contrib/gael/matlab_class/sample_analysis/example_transports.m
ViewVC logotype

Annotation of /MITgcm_contrib/gael/matlab_class/sample_analysis/example_transports.m

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


Revision 1.2 - (hide annotations) (download)
Fri Mar 20 15:41:27 2015 UTC (10 years, 4 months ago) by gforget
Branch: MAIN
Changes since 1.1: +10 -23 lines
- remove deprecated choiceGrid (or choiceV3orV4) switches and arguments

1 gforget 1.2 function [diags]=example_transports();
2 gforget 1.1 %object: illustrates a series of standard computations
3     % (streamfunctions, transports, zonal means, etc.)
4    
5     gcmfaces_global;
6     if myenv.verbose>0;
7     gcmfaces_msg('===============================================');
8     gcmfaces_msg(['*** entering example_transports ' ...
9     'that will load a set of variables form file ' ...
10     'and compute a series of derived diagnostics '],'');
11     % gcmfaces_msg('*** entering example_transports','');
12     % gcmfaces_msg('that will load a set of variables form file',' ');
13     % gcmfaces_msg('and compute a series of derived diagnostics',' ');
14     end;
15    
16     %%%%%%%%%%%%%%%%%
17     %load parameters:
18     %%%%%%%%%%%%%%%%%
19    
20     if myenv.verbose>0;
21     gcmfaces_msg('* set grid files path and format, and number of faces for this grid');
22     end;
23     myenv.nctiles=1;
24 gforget 1.2 nF=5;
25     fileFormat='compact';
26     dirGrid=fullfile(myenv.gcmfaces_dir,'..',filesep,'GRID',filesep);
27     myenv.nctilesdir=fullfile(myenv.gcmfaces_dir,'..',filesep,...
28     'release1',filesep,'nctiles_climatology',filesep);
29    
30 gforget 1.1 if myenv.verbose>0;
31     gcmfaces_msg('* call grid_load : load grid to memory (mygrid) according to');
32     gcmfaces_msg(['dirGrid = ' dirGrid],' ');
33     gcmfaces_msg(['nFaces = ' num2str(nF)],' ');
34     gcmfaces_msg(['fileFormat = ' fileFormat],' ');
35     % fprintf([' > dirGrid = ' dirGrid '\n']);
36     end;
37    
38     if ~isdir(myenv.nctilesdir);
39     diags=[];
40     warning(['skipping example_transports (missing ' myenv.nctilesdir ')']);
41     return;
42     end;
43    
44     % warning('skipping grid_load\n');
45     grid_load(dirGrid,nF,fileFormat);
46    
47     if myenv.verbose>0;
48     gcmfaces_msg('* call gcmfaces_lines_zonal : determine grid lines that closely follow');
49     gcmfaces_msg('parallel lines and will be used in zonal mean and overturning computations',' ');
50     end;
51 gforget 1.2 gcmfaces_lines_zonal;
52    
53 gforget 1.1 if myenv.verbose>0;
54     gcmfaces_msg('* call gcmfaces_lines_transp : determine grid lines that closely follow');
55     gcmfaces_msg('great circles and will be used to compute transsects transports',' ');
56     end;
57     % warning('skipping gcmfaces_lines_transp\n');
58 gforget 1.2 [lonPairs,latPairs,names]=line_greatC_TUV_MASKS_v4;
59 gforget 1.1 gcmfaces_lines_transp(lonPairs,latPairs,names);
60    
61     %%%%%%%%%%%%%%%%%
62     %do computations:
63     %%%%%%%%%%%%%%%%%
64    
65     % [listTimes]=diags_list_times;
66     diags.listTimes=1;
67    
68     if myenv.verbose>0; gcmfaces_msg('* call rdmds2gcmfaces : load velocity fields');end;
69     listVars={'UVELMASS','VVELMASS'};
70     listDiags={'fldBAR','gloOV','fldTRANSPORTS','gloMT_FW'};
71     for vvv=1:length(listVars);
72     vv=listVars{vvv};
73     tmp1=read_nctiles([myenv.nctilesdir vv '/' vv],vv);
74     tmp1=mean(tmp1,4);
75     tmp1(mygrid.mskC==0)=NaN;
76     eval([vv '=tmp1;']);
77     end;
78    
79     UVELMASS=UVELMASS.*mygrid.mskW;
80     VVELMASS=VVELMASS.*mygrid.mskS;
81    
82     if myenv.verbose>0; gcmfaces_msg('* call calc_barostream : comp. barotropic stream function');end;
83     [fldBAR]=calc_barostream(UVELMASS,VVELMASS);
84     if myenv.verbose>0; gcmfaces_msg('* call calc_overturn : comp. overturning stream function');end;
85     [gloOV]=calc_overturn(UVELMASS,VVELMASS);
86     if myenv.verbose>0; gcmfaces_msg('* call calc_transports : comp. transects transports');end;
87     [fldTRANSPORTS]=1e-6*calc_transports(UVELMASS,VVELMASS,mygrid.LINES_MASKS,{'dh','dz'});
88     if myenv.verbose>0; gcmfaces_msg('* call calc_MeridionalTransport : comp. meridional seawater transport');end;
89     [gloMT_FW]=1e-6*calc_MeridionalTransport(UVELMASS,VVELMASS,1);
90    
91     if myenv.verbose>0; gcmfaces_msg('* load tracer and transports fields');end;
92     listVars={'THETA','SALT','ADVx_TH','ADVy_TH','ADVx_SLT','ADVy_SLT'};
93     listVars={listVars{:},'DFxE_TH','DFyE_TH','DFxE_SLT','DFyE_SLT'};
94     listDiags={listDiags{:},'fldTzonmean','fldSzonmean','gloMT_H','gloMT_SLT'};
95     for vvv=1:length(listVars);
96     vv=listVars{vvv};
97     tmp1=read_nctiles([myenv.nctilesdir vv '/' vv],vv);
98     tmp1=mean(tmp1,4);
99     tmp1(mygrid.mskC==0)=NaN;
100     eval([vv '=tmp1;']);
101     end;
102    
103     if myenv.verbose>0; gcmfaces_msg('* call calc_zonmean_T : comp. zonal mean temperature');end;
104     [fldTzonmean]=calc_zonmean_T(THETA);
105     if myenv.verbose>0; gcmfaces_msg('* call calc_zonmean_T : comp. zonal mean salinity');end;
106     [fldSzonmean]=calc_zonmean_T(SALT);
107    
108     if myenv.verbose>0; gcmfaces_msg('* call calc_MeridionalTransport : comp. meridional heat transport');end;
109     tmpU=(ADVx_TH+DFxE_TH); tmpV=(ADVy_TH+DFyE_TH);
110     [gloMT_H]=1e-15*4e6*calc_MeridionalTransport(tmpU,tmpV,0);
111     if myenv.verbose>0; gcmfaces_msg('* call calc_MeridionalTransport : comp. meridional salt transport');end;
112     tmpU=(ADVx_SLT+DFxE_SLT); tmpV=(ADVy_SLT+DFyE_SLT);
113     [gloMT_SLT]=1e-6*calc_MeridionalTransport(tmpU,tmpV,0);
114    
115     %format output:
116     for ddd=1:length(listDiags);
117     dd=listDiags{ddd};
118     eval(['diags.' dd '=' dd ';']);
119     end;
120    
121     if myenv.verbose>0;
122     gcmfaces_msg('*** leaving example_transports');
123     gcmfaces_msg('===============================================','');
124     end;
125    

  ViewVC Help
Powered by ViewVC 1.1.22