/[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.6 - (hide annotations) (download)
Sun Jan 24 16:25:51 2016 UTC (9 years, 6 months ago) by gforget
Branch: MAIN
Changes since 1.5: +2 -5 lines
- gcmfaces_demo.m: improve help section, add comments wrt gcmfaces_global etc
- example*.m: homogeneize treatment of mygrid (don't reload if not needed)
  and path (use just gcmfaces_global)

1 gforget 1.2 function [diags]=example_transports();
2 gforget 1.4 %EXAMPLE_TRANSPORTS computes transports (and zonal averages)
3 gforget 1.5 %
4     % Note: this demonstration routine includes its own call to grid_load
5     %
6 gforget 1.6 % Example:
7     % addpath gcmfaces/; gcmfaces_global; myenv.verbose=1;
8 gforget 1.5 % [diags]=example_transports();
9     % example_transports_disp(diags);
10 gforget 1.4
11 gforget 1.1 gcmfaces_global;
12     if myenv.verbose>0;
13     gcmfaces_msg('===============================================');
14     gcmfaces_msg(['*** entering example_transports ' ...
15     'that will load a set of variables form file ' ...
16     'and compute a series of derived diagnostics '],'');
17     % gcmfaces_msg('*** entering example_transports','');
18     % gcmfaces_msg('that will load a set of variables form file',' ');
19     % gcmfaces_msg('and compute a series of derived diagnostics',' ');
20     end;
21    
22     %%%%%%%%%%%%%%%%%
23 gforget 1.4 %load grid:
24 gforget 1.1 %%%%%%%%%%%%%%%%%
25    
26 gforget 1.4 if isempty(mygrid);
27     grid_load;
28 gforget 1.1 end;
29 gforget 1.2
30 gforget 1.4 myenv.nctilesdir=fullfile('release1',filesep,'nctiles_climatology',filesep);
31 gforget 1.1
32     if ~isdir(myenv.nctilesdir);
33     diags=[];
34 gforget 1.3 help gcmfaces_demo;
35 gforget 1.1 warning(['skipping example_transports (missing ' myenv.nctilesdir ')']);
36     return;
37     end;
38    
39     if myenv.verbose>0;
40     gcmfaces_msg('* call gcmfaces_lines_zonal : determine grid lines that closely follow');
41     gcmfaces_msg('parallel lines and will be used in zonal mean and overturning computations',' ');
42     end;
43 gforget 1.2 gcmfaces_lines_zonal;
44    
45 gforget 1.1 if myenv.verbose>0;
46     gcmfaces_msg('* call gcmfaces_lines_transp : determine grid lines that closely follow');
47     gcmfaces_msg('great circles and will be used to compute transsects transports',' ');
48     end;
49     % warning('skipping gcmfaces_lines_transp\n');
50 gforget 1.5 [lonPairs,latPairs,names]=gcmfaces_lines_pairs;
51 gforget 1.1 gcmfaces_lines_transp(lonPairs,latPairs,names);
52    
53     %%%%%%%%%%%%%%%%%
54     %do computations:
55     %%%%%%%%%%%%%%%%%
56    
57 gforget 1.4 listDiags={};
58    
59 gforget 1.1 % [listTimes]=diags_list_times;
60     diags.listTimes=1;
61    
62 gforget 1.4 %part 1:
63    
64     listVars={'UVELMASS','VVELMASS'};
65     missingVars={};
66     for vv=1:length(listVars);
67     tmp1=[myenv.nctilesdir listVars{vv} filesep listVars{vv} '*nc'];
68     if isempty(dir(tmp1)); missingVars={missingVars{:},tmp1}; end;
69     end;
70    
71     if ~isempty(missingVars);
72     fprintf('\n example_transports could not find the following files ---> skipping related computation!\n');
73     disp(missingVars');
74     else;
75    
76 gforget 1.1 if myenv.verbose>0; gcmfaces_msg('* call rdmds2gcmfaces : load velocity fields');end;
77 gforget 1.4
78 gforget 1.1 for vvv=1:length(listVars);
79     vv=listVars{vvv};
80     tmp1=read_nctiles([myenv.nctilesdir vv '/' vv],vv);
81     tmp1=mean(tmp1,4);
82     tmp1(mygrid.mskC==0)=NaN;
83     eval([vv '=tmp1;']);
84     end;
85    
86     UVELMASS=UVELMASS.*mygrid.mskW;
87     VVELMASS=VVELMASS.*mygrid.mskS;
88    
89 gforget 1.4 listDiags={listDiags{:},'fldBAR','gloOV','fldTRANSPORTS','gloMT_FW'};
90 gforget 1.1 if myenv.verbose>0; gcmfaces_msg('* call calc_barostream : comp. barotropic stream function');end;
91     [fldBAR]=calc_barostream(UVELMASS,VVELMASS);
92     if myenv.verbose>0; gcmfaces_msg('* call calc_overturn : comp. overturning stream function');end;
93     [gloOV]=calc_overturn(UVELMASS,VVELMASS);
94     if myenv.verbose>0; gcmfaces_msg('* call calc_transports : comp. transects transports');end;
95     [fldTRANSPORTS]=1e-6*calc_transports(UVELMASS,VVELMASS,mygrid.LINES_MASKS,{'dh','dz'});
96     if myenv.verbose>0; gcmfaces_msg('* call calc_MeridionalTransport : comp. meridional seawater transport');end;
97     [gloMT_FW]=1e-6*calc_MeridionalTransport(UVELMASS,VVELMASS,1);
98    
99 gforget 1.4 end;%if ~isempty(missingVars);
100    
101     %part 2:
102    
103 gforget 1.1 listVars={'THETA','SALT','ADVx_TH','ADVy_TH','ADVx_SLT','ADVy_SLT'};
104     listVars={listVars{:},'DFxE_TH','DFyE_TH','DFxE_SLT','DFyE_SLT'};
105 gforget 1.4
106     missingVars={};
107     for vv=1:length(listVars);
108     tmp1=[myenv.nctilesdir listVars{vv} filesep listVars{vv} '*nc'];
109     if isempty(dir(tmp1)); missingVars={missingVars{:},tmp1}; end;
110     end;
111    
112     if ~isempty(missingVars);
113     fprintf('\n example_transports could not find the following files ---> skipping related computation!\n');
114     disp(missingVars');
115     else;
116    
117     if myenv.verbose>0; gcmfaces_msg('* load tracer and transports fields');end;
118 gforget 1.1 listDiags={listDiags{:},'fldTzonmean','fldSzonmean','gloMT_H','gloMT_SLT'};
119     for vvv=1:length(listVars);
120     vv=listVars{vvv};
121     tmp1=read_nctiles([myenv.nctilesdir vv '/' vv],vv);
122     tmp1=mean(tmp1,4);
123     tmp1(mygrid.mskC==0)=NaN;
124     eval([vv '=tmp1;']);
125     end;
126    
127     if myenv.verbose>0; gcmfaces_msg('* call calc_zonmean_T : comp. zonal mean temperature');end;
128     [fldTzonmean]=calc_zonmean_T(THETA);
129     if myenv.verbose>0; gcmfaces_msg('* call calc_zonmean_T : comp. zonal mean salinity');end;
130     [fldSzonmean]=calc_zonmean_T(SALT);
131    
132     if myenv.verbose>0; gcmfaces_msg('* call calc_MeridionalTransport : comp. meridional heat transport');end;
133     tmpU=(ADVx_TH+DFxE_TH); tmpV=(ADVy_TH+DFyE_TH);
134     [gloMT_H]=1e-15*4e6*calc_MeridionalTransport(tmpU,tmpV,0);
135     if myenv.verbose>0; gcmfaces_msg('* call calc_MeridionalTransport : comp. meridional salt transport');end;
136     tmpU=(ADVx_SLT+DFxE_SLT); tmpV=(ADVy_SLT+DFyE_SLT);
137     [gloMT_SLT]=1e-6*calc_MeridionalTransport(tmpU,tmpV,0);
138    
139 gforget 1.4 end;%if ~isempty(missingVars);
140    
141     %part 3: format output
142    
143 gforget 1.1 for ddd=1:length(listDiags);
144     dd=listDiags{ddd};
145     eval(['diags.' dd '=' dd ';']);
146     end;
147    
148     if myenv.verbose>0;
149     gcmfaces_msg('*** leaving example_transports');
150     gcmfaces_msg('===============================================','');
151     end;
152    

  ViewVC Help
Powered by ViewVC 1.1.22