/[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.3 - (hide annotations) (download)
Fri Mar 20 23:56:20 2015 UTC (10 years, 4 months ago) by gforget
Branch: MAIN
CVS Tags: checkpoint65r, checkpoint65p, checkpoint65q
Changes since 1.2: +1 -0 lines
- warning messages: report help gcmfaces_demo

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 gforget 1.3 help gcmfaces_demo;
41 gforget 1.1 warning(['skipping example_transports (missing ' myenv.nctilesdir ')']);
42     return;
43     end;
44    
45     % warning('skipping grid_load\n');
46     grid_load(dirGrid,nF,fileFormat);
47    
48     if myenv.verbose>0;
49     gcmfaces_msg('* call gcmfaces_lines_zonal : determine grid lines that closely follow');
50     gcmfaces_msg('parallel lines and will be used in zonal mean and overturning computations',' ');
51     end;
52 gforget 1.2 gcmfaces_lines_zonal;
53    
54 gforget 1.1 if myenv.verbose>0;
55     gcmfaces_msg('* call gcmfaces_lines_transp : determine grid lines that closely follow');
56     gcmfaces_msg('great circles and will be used to compute transsects transports',' ');
57     end;
58     % warning('skipping gcmfaces_lines_transp\n');
59 gforget 1.2 [lonPairs,latPairs,names]=line_greatC_TUV_MASKS_v4;
60 gforget 1.1 gcmfaces_lines_transp(lonPairs,latPairs,names);
61    
62     %%%%%%%%%%%%%%%%%
63     %do computations:
64     %%%%%%%%%%%%%%%%%
65    
66     % [listTimes]=diags_list_times;
67     diags.listTimes=1;
68    
69     if myenv.verbose>0; gcmfaces_msg('* call rdmds2gcmfaces : load velocity fields');end;
70     listVars={'UVELMASS','VVELMASS'};
71     listDiags={'fldBAR','gloOV','fldTRANSPORTS','gloMT_FW'};
72     for vvv=1:length(listVars);
73     vv=listVars{vvv};
74     tmp1=read_nctiles([myenv.nctilesdir vv '/' vv],vv);
75     tmp1=mean(tmp1,4);
76     tmp1(mygrid.mskC==0)=NaN;
77     eval([vv '=tmp1;']);
78     end;
79    
80     UVELMASS=UVELMASS.*mygrid.mskW;
81     VVELMASS=VVELMASS.*mygrid.mskS;
82    
83     if myenv.verbose>0; gcmfaces_msg('* call calc_barostream : comp. barotropic stream function');end;
84     [fldBAR]=calc_barostream(UVELMASS,VVELMASS);
85     if myenv.verbose>0; gcmfaces_msg('* call calc_overturn : comp. overturning stream function');end;
86     [gloOV]=calc_overturn(UVELMASS,VVELMASS);
87     if myenv.verbose>0; gcmfaces_msg('* call calc_transports : comp. transects transports');end;
88     [fldTRANSPORTS]=1e-6*calc_transports(UVELMASS,VVELMASS,mygrid.LINES_MASKS,{'dh','dz'});
89     if myenv.verbose>0; gcmfaces_msg('* call calc_MeridionalTransport : comp. meridional seawater transport');end;
90     [gloMT_FW]=1e-6*calc_MeridionalTransport(UVELMASS,VVELMASS,1);
91    
92     if myenv.verbose>0; gcmfaces_msg('* load tracer and transports fields');end;
93     listVars={'THETA','SALT','ADVx_TH','ADVy_TH','ADVx_SLT','ADVy_SLT'};
94     listVars={listVars{:},'DFxE_TH','DFyE_TH','DFxE_SLT','DFyE_SLT'};
95     listDiags={listDiags{:},'fldTzonmean','fldSzonmean','gloMT_H','gloMT_SLT'};
96     for vvv=1:length(listVars);
97     vv=listVars{vvv};
98     tmp1=read_nctiles([myenv.nctilesdir vv '/' vv],vv);
99     tmp1=mean(tmp1,4);
100     tmp1(mygrid.mskC==0)=NaN;
101     eval([vv '=tmp1;']);
102     end;
103    
104     if myenv.verbose>0; gcmfaces_msg('* call calc_zonmean_T : comp. zonal mean temperature');end;
105     [fldTzonmean]=calc_zonmean_T(THETA);
106     if myenv.verbose>0; gcmfaces_msg('* call calc_zonmean_T : comp. zonal mean salinity');end;
107     [fldSzonmean]=calc_zonmean_T(SALT);
108    
109     if myenv.verbose>0; gcmfaces_msg('* call calc_MeridionalTransport : comp. meridional heat transport');end;
110     tmpU=(ADVx_TH+DFxE_TH); tmpV=(ADVy_TH+DFyE_TH);
111     [gloMT_H]=1e-15*4e6*calc_MeridionalTransport(tmpU,tmpV,0);
112     if myenv.verbose>0; gcmfaces_msg('* call calc_MeridionalTransport : comp. meridional salt transport');end;
113     tmpU=(ADVx_SLT+DFxE_SLT); tmpV=(ADVy_SLT+DFyE_SLT);
114     [gloMT_SLT]=1e-6*calc_MeridionalTransport(tmpU,tmpV,0);
115    
116     %format output:
117     for ddd=1:length(listDiags);
118     dd=listDiags{ddd};
119     eval(['diags.' dd '=' dd ';']);
120     end;
121    
122     if myenv.verbose>0;
123     gcmfaces_msg('*** leaving example_transports');
124     gcmfaces_msg('===============================================','');
125     end;
126    

  ViewVC Help
Powered by ViewVC 1.1.22