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