/[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.5 - (hide annotations) (download)
Wed Jan 13 16:24:11 2016 UTC (9 years, 6 months ago) by gforget
Branch: MAIN
Changes since 1.4: +8 -1 lines
- replace sample_analysis/line_greatC_TUV_MASKS_v4.m
     with gcmfaces_calc/gcmfaces_lines_pairs.m

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

  ViewVC Help
Powered by ViewVC 1.1.22