1 |
function HT = calcHeatTransDirect(varargin) |
2 |
|
3 |
% HT = calcHeatTransDirect(d,g,dE,GM,time,flu,blkFile,[mask,phi_ref]); |
4 |
% |
5 |
% Input arguements: |
6 |
% The incoming field data (d) and grid data (g) must be in a structured |
7 |
% array format (which is the format that comes from rdmnc): |
8 |
% d [Structure] Various inputs depending on oceanic or |
9 |
% atmospheric computation, rstar or not, |
10 |
% and GM form/ouputs |
11 |
% g [Structure] drF,dxG,dyG,dxC,dyC,HFacW,HFacS,rA |
12 |
% dE [Structure] More field data |
13 |
% GM [Structure] Gent anf McWilliams data |
14 |
% Other input parameters: |
15 |
% time [vector] Time levels to analyze ([] for all) |
16 |
% flu [string] 'O' or 'A' for ocean or atmosphere |
17 |
% blkFile [string] Broken line file (eg 'isoLat_cs32_59.mat') |
18 |
% Optional parameters: |
19 |
% mask [Structure] W and S mask |
20 |
% phi_ref [vector] Reference geopotential height |
21 |
% |
22 |
% Output: |
23 |
% HT_Out is a structured array with the following fields: |
24 |
% time ([nt,1]) Time axis |
25 |
% ylatHT ([ny,1]) Heat transport latitude axis |
26 |
% ylatSF ([ny-1,1]) Implied heating latitude axis |
27 |
% AreaZon ([ny,1]) Area in zonal bands |
28 |
% SenHT ([ny,nBas,nt,nFld]) Sensible heat transport |
29 |
% SenSF ([ny-1,nBas,nt,nFld]) Implied heating above |
30 |
% LatHT ([ny,nBas,nt,nFld]) Latent heat transport (atm only) |
31 |
% LatSF ([ny-1,nBas,nt,nFld]) Implied heating from aboce (atm only) |
32 |
% Currently, the routine is only configured to handle the global basin, |
33 |
% so nBas = 1 for the output. ny is defined by the broken line file used |
34 |
% for the cube calculation. nFld is the heat transport component: |
35 |
% nFld = 1 = Eulerian circulation HT |
36 |
% nFld = 2 = HT by time mean circulation |
37 |
% nFld = 3 = Residual [3=1-2] |
38 |
% nFld = 4 = HT by zonal mean circulation |
39 |
% nFld = 5 = Residual [5=2-4] |
40 |
% Ocn only: |
41 |
% If GM advective form: |
42 |
% nFld = 6 = total (GM+Redi) eddy transport: K*dT/dy |
43 |
% WARNING: This assumes K_gm=K_iso to get |
44 |
% the cancellation of the off-diagonal term |
45 |
% (K_iso-K_gm)*S_y*dT/dz |
46 |
% nFld = 7 = GM (advective) eddy transport |
47 |
% If Skew flux form: |
48 |
% nFld = 6 = total (GM+Redi) eddy transport |
49 |
% nFld = 7 = GM (advective) eddy transport |
50 |
% |
51 |
% Description: |
52 |
% Calculation heat transport, and to degree possible, decomposition. |
53 |
% Heat transport is given in PW and the implied surface heating/flux |
54 |
% in W/m^2. The incoming data arrays are all assumed to be of the |
55 |
% dimensions [6*nc,nc,nr,nt]. |
56 |
% |
57 |
% Original Author: Jean-Michel Campin |
58 |
% Modifications: Daniel Enderton |
59 |
|
60 |
% Default constants |
61 |
LhVap = 2501; |
62 |
grav = 9.81; |
63 |
CpO = 3994; |
64 |
RhoO = 1030; |
65 |
CpA = 1004; |
66 |
kappa = 2/7; |
67 |
masking=0; |
68 |
|
69 |
PHIref = grav * [ 431.199 2105.998 5536.327 10440.915 18190.889]; |
70 |
%PHIref = grav * [ 173.888 529.598 896.724 1276.210 1669.052 ... |
71 |
% 2076.380 2499.554 2940.137 3399.937 3881.061 ... |
72 |
% 4385.993 4917.689 5480.084 6078.962 6722.700 ... |
73 |
% 7420.744 8183.151 9023.181 9959.131 11019.028... |
74 |
% 12253.003 13723.604 15599.238 18412.089 24539.570]; |
75 |
|
76 |
% Read input parameters. |
77 |
d = varargin{1}; |
78 |
g = varargin{2}; |
79 |
dE= varargin{3}; |
80 |
GM= varargin{4}; |
81 |
nt= varargin{5}; |
82 |
flu = varargin{6}; |
83 |
blkFile = varargin{7}; |
84 |
if length(varargin) >= 8 & ~isempty(varargin{8}) |
85 |
mask = varargin{8}; |
86 |
masking = 1; |
87 |
disp('Masking is on...') |
88 |
end |
89 |
if length(varargin) == 9 |
90 |
PHIref= varargin{9}; |
91 |
end |
92 |
|
93 |
nBas = 0; |
94 |
if isequal(flu,'A'), nout = 5; end |
95 |
if isequal(flu,'O'), nout = 7; end |
96 |
flag_tot=0; |
97 |
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
98 |
% Prepare / reform incoming data % |
99 |
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
100 |
|
101 |
% Determine time indecies. |
102 |
%if isempty(time), time = d.iters_read_from_file; i_time = 1:length(time); |
103 |
%else [dump,i_time] = ismember(time,d.iters_read_from_file); end |
104 |
|
105 |
nc = size(g.HFacC,2); |
106 |
nr = size(g.HFacC,3); |
107 |
i_time=1:nt; |
108 |
|
109 |
if nr~=5 |
110 |
LhVap = 2500e3; |
111 |
grav = 9.79764; |
112 |
RDGAS = 287.04; |
113 |
CpA = RDGAS/kappa; |
114 |
end |
115 |
|
116 |
ac = reshape(g.rA ,[6*nc*nc, 1]); |
117 |
hw = reshape(g.HFacW(1:6*nc,1:nc,1:nr),[6*nc*nc,nr]); |
118 |
hs = reshape(g.HFacS(1:6*nc,1:nc,1:nr),[6*nc*nc,nr]); |
119 |
dxc = reshape(g.dxC(1:6*nc,1:nc) ,[6*nc*nc, 1]); |
120 |
dyc = reshape(g.dyC(1:6*nc,1:nc) ,[6*nc*nc, 1]); |
121 |
dxg = reshape(g.dxG(1:6*nc,1:nc) ,[6*nc*nc, 1]); |
122 |
dyg = reshape(g.dyG(1:6*nc,1:nc) ,[6*nc*nc, 1]); |
123 |
drf = reshape(g.drF,[1,length(g.drF)]); |
124 |
Z = reshape(g.Z,[1,length(g.Z)]); |
125 |
|
126 |
%u = reshape(d.UVEL(1:6*nc,1:nc,:,i_time),[6*nc*nc,nr,nt]); |
127 |
%v = reshape(d.VVEL(1:6*nc,1:nc,:,i_time),[6*nc*nc,nr,nt]); |
128 |
ut = reshape(dE.UTHMASS(1:6*nc,1:nc,:,i_time) ,[6*nc*nc,nr,nt]); |
129 |
vt = reshape(dE.VTHMASS(1:6*nc,1:nc,:,i_time) ,[6*nc*nc,nr,nt]); |
130 |
hu = reshape(dE.UVELMASS(1:6*nc,1:nc,:,i_time) ,[6*nc*nc,nr,nt]); |
131 |
hv = reshape(dE.VVELMASS(1:6*nc,1:nc,:,i_time) ,[6*nc*nc,nr,nt]); |
132 |
t = d.THETA(1:6*nc,1:nc,:,i_time); |
133 |
if isequal(flu,'O') |
134 |
Kux = reshape(GM.GM_Kux(1:6*nc,1:nc,:,i_time),[6*nc*nc,nr,nt]); |
135 |
Kvy = reshape(GM.GM_Kvy(1:6*nc,1:nc,:,i_time),[6*nc*nc,nr,nt]); |
136 |
|
137 |
if isfield(GM,'GM_ubT') %%% if advective form of GM |
138 |
UbTgm = reshape(GM.GM_ubT(1:6*nc,1:nc,:,i_time),[6*nc*nc,nr,nt]); |
139 |
VbTgm = reshape(GM.GM_vbT(1:6*nc,1:nc,:,i_time),[6*nc*nc,nr,nt]); |
140 |
elseif isfield(GM,'GM_KuzTz') %%% if skew-flux form of GM |
141 |
UbTgm = reshape(GM.GM_KuzTz(1:6*nc,1:nc,:,i_time),[6*nc*nc,nr,nt]); |
142 |
VbTgm = reshape(GM.GM_KvzTz(1:6*nc,1:nc,:,i_time),[6*nc*nc,nr,nt]); |
143 |
elseif isfield(dE,'ADVx_TH') %%% if total advective field |
144 |
flag_tot = 1; |
145 |
UbTgm = reshape(dE.ADVx_TH(1:6*nc,1:nc,:,i_time),[6*nc*nc,nr,nt]); |
146 |
VbTgm = reshape(dE.ADVy_TH(1:6*nc,1:nc,:,i_time),[6*nc*nc,nr,nt]); |
147 |
else |
148 |
UbTgm = zeros(6*nc*nc,nr,nt); |
149 |
VbTgm = zeros(6*nc*nc,nr,nt); |
150 |
end |
151 |
end |
152 |
|
153 |
if isequal(flu,'A') |
154 |
up = reshape(dE.UVELPHI(1:6*nc,1:nc,:,i_time) ,[6*nc*nc,nr,nt]); |
155 |
vp = reshape(dE.VVELPHI(1:6*nc,1:nc,:,i_time) ,[6*nc*nc,nr,nt]); |
156 |
uq = reshape(dE.USLTMASS(1:6*nc,1:nc,:,i_time) ,[6*nc*nc,nr,nt]); |
157 |
vq = reshape(dE.VSLTMASS(1:6*nc,1:nc,:,i_time) ,[6*nc*nc,nr,nt]); |
158 |
q = d.SALT(1:6*nc,1:nc,:,i_time); |
159 |
p = d.PHIHYD(1:6*nc,1:nc,:,i_time); |
160 |
ut = ut .* ( repmat(Z,[6*nc*nc 1 nt]) / 1e5 ).^kappa; |
161 |
vt = vt .* ( repmat(Z,[6*nc*nc 1 nt]) / 1e5 ).^kappa; |
162 |
t = t .* ( repmat(reshape(Z,[1 1 nr]),[6*nc nc 1 nt]) / 1e5 ).^kappa; |
163 |
end |
164 |
|
165 |
if masking == 1 |
166 |
hu = repmat(reshape(mask.maskW(:,:,1),6*nc*nc,1),[1 nr nt]) .* hu; |
167 |
hv = repmat(reshape(mask.maskS(:,:,1),6*nc*nc,1),[1 nr nt]) .* hv; |
168 |
ut = repmat(reshape(mask.maskW(:,:,1),6*nc*nc,1),[1 nr nt]) .* ut; |
169 |
vt = repmat(reshape(mask.maskS(:,:,1),6*nc*nc,1),[1 nr nt]) .* vt; |
170 |
if isequal(flu,'O') |
171 |
Kux = repmat(reshape(mask.maskW(:,:,1),6*nc*nc,1),[1 nr nt]) .* Kux; |
172 |
Kvy = repmat(reshape(mask.maskS(:,:,1),6*nc*nc,1),[1 nr nt]) .* Kvy; |
173 |
UbTgm = repmat(reshape(mask.maskW(:,:,1),6*nc*nc,1),[1 nr nt]) .* UbTgm; |
174 |
VbTgm = repmat(reshape(mask.maskS(:,:,1),6*nc*nc,1),[1 nr nt]) .* VbTgm; |
175 |
end |
176 |
end |
177 |
|
178 |
|
179 |
% Load broken line information. Compute (tracer point) cell area between |
180 |
% broken lines for each basin. There are nbkl broken lines and nbkl+1 |
181 |
% bands between broken lines. The variable "bkl_Zon" gives the zone number |
182 |
% (nbkl+1 total) for a given index between 0 and nbkl, that is, nbkl+1 |
183 |
% total zones. Comments block is for eventual addition of multiple basin |
184 |
% calculations. |
185 |
load(blkFile); |
186 |
ydim=length(bkl_Ylat); |
187 |
AreaZon=zeros(ydim+1,1+nBas); |
188 |
for j = 1:ydim+1 |
189 |
izon = find(bkl_Zon == j-1); |
190 |
AreaZon(j,1) = sum(ac(izon)); |
191 |
% for b = 1:nBas, |
192 |
% tmp = ac.*mskBc(:,b); |
193 |
% AreaZon(j,1+b) = sum(tmp(izon)); |
194 |
% end |
195 |
end |
196 |
|
197 |
% Latitute plotting information. Average latitude of a broken line |
198 |
% (ylatHF) is calculated from a mean value of the y value of all the edges. |
199 |
% The latitude at the surface flux points is a mean of the broken line mean |
200 |
% values. |
201 |
YlatAv=sum(bkl_Ysg,1)./(1+bkl_Npts'); %' |
202 |
ylatHT = [-90,YlatAv,90]; |
203 |
ylatSF = ( ylatHT(2:ydim+2) + ylatHT(1:ydim+1) )./2; |
204 |
|
205 |
% The variable "bkl_Flg" is -1/1 if edge (on a given broken) has a u point |
206 |
% and -2/2 if it has a v point. Positive/negative values contribute |
207 |
% positively/negatively to northward heat transport (this depends on the |
208 |
% orientation of the cell). A zero value indicates an end of edges that |
209 |
% contribute to a broken line. The u and v information is parced into two |
210 |
% seperate fields, ufac and vfac (-2/2 are reduced to -1/1 for vfac). |
211 |
ufac = zeros([size(bkl_Flg),1+nBas]); |
212 |
vfac = zeros([size(bkl_Flg),1+nBas]); |
213 |
ufac(:,:,1) = rem(bkl_Flg,2); |
214 |
vfac(:,:,1) = fix(bkl_Flg/2); |
215 |
% for jl=1:ydim, |
216 |
% ie=bkl_Npts(jl); |
217 |
% for b=1:nBas, |
218 |
% ufac(1:ie,jl,1+b)=mskBw(bkl_IJuv(1:ie,jl),b).*ufac(1:ie,jl,1); |
219 |
% vfac(1:ie,jl,1+b)=mskBs(bkl_IJuv(1:ie,jl),b).*vfac(1:ie,jl,1); |
220 |
% end |
221 |
% end |
222 |
ufacabs = abs(ufac); |
223 |
vfacabs = abs(vfac); |
224 |
|
225 |
% Prepare mask(s). |
226 |
% ??? I temporarily took out the code to configure the masks beyond this |
227 |
% global one. Does this need to account for a ridge if present? |
228 |
mskG=ones(ydim+1,1+nBas); |
229 |
|
230 |
% Area factors. "ArW_Dif" and "ArS_Dif" the areas for the western and |
231 |
% southern edge of cells, respectively The "_Dif" suffix indicates the |
232 |
% areas used for the diffusivity because there is an extra "dxc" or "dyc" |
233 |
% from the gradient (here for computational efficiency reasons). The |
234 |
% division by "dxc" and "dyc" is associates with gradient of temperature. |
235 |
% |
236 |
% ??? Why is the land mask (hw and hs) only in the diffusive area term? |
237 |
ArW = dyg*reshape(drf,[1,length(drf)]); |
238 |
ArS = dxg*reshape(drf,[1,length(drf)]); |
239 |
ArW_Dif=hw.*((dyg./dxc)*reshape(drf,[1,length(drf)])); |
240 |
ArS_Dif=hs.*((dxg./dyc)*reshape(drf,[1,length(drf)])); |
241 |
|
242 |
% Compute the temperature and its gradient and the velocity points: |
243 |
% tbi/tdi = temperature between/difference i points (at u points) |
244 |
% tbj/tdj = temperature between/difference j points (at v points) |
245 |
% The cube is first split and the extra points are added to the files. |
246 |
% Then the means and differences are taken. Note that the division of dxc |
247 |
% and dyc for the gradient is not applied until later for computational |
248 |
% efficiency purposes. The arrays are then croped and reshaped to the |
249 |
% format for the other field variables, [6*nc*nc,nr]. Note that the |
250 |
% split_C_cub function adds a row/column of tracer points in front of the |
251 |
% first row/column of the tile matries. This, when the differences and |
252 |
% gradients are computed and cropped, the off indecies are selected from |
253 |
% [2:nc+1] rather than [1:nc]. (This was a bit mystifying to me). |
254 |
t6bi=zeros(nc,nc+1,nr,nt,6); t6di=t6bi; q6bi=t6bi; |
255 |
t6bj=zeros(nc+1,nc,nr,nt,6); t6dj=t6bj; q6bj=t6bj; |
256 |
t6t=split_C_cub(t); |
257 |
t6bi([1:nc],:,:,:,:) = ( t6t([1:nc],:,:,:,:) + t6t([2:nc+1],:,:,:,:) )./2; |
258 |
t6bj(:,[1:nc],:,:,:) = ( t6t(:,[1:nc],:,:,:) + t6t(:,[2:nc+1],:,:,:) )./2; |
259 |
t6di([1:nc],:,:,:,:) = ( t6t([1:nc],:,:,:,:) - t6t([2:nc+1],:,:,:,:) ); |
260 |
t6dj(:,[1:nc],:,:,:) = ( t6t(:,[1:nc],:,:,:) - t6t(:,[2:nc+1],:,:,:) ); |
261 |
tbi = t6bi([1:nc],[2:nc+1],:,:,:); |
262 |
tbj = t6bj([2:nc+1],[1:nc],:,:,:); |
263 |
tdi = t6di([1:nc],[2:nc+1],:,:,:); |
264 |
tdj = t6dj([2:nc+1],[1:nc],:,:,:); |
265 |
tbi=reshape(permute(tbi,[1,5,2,3,4]),[6*nc*nc,nr,nt]); |
266 |
tbj=reshape(permute(tbj,[1,5,2,3,4]),[6*nc*nc,nr,nt]); |
267 |
tdi=reshape(permute(tdi,[1,5,2,3,4]),[6*nc*nc,nr,nt]); |
268 |
tdj=reshape(permute(tdj,[1,5,2,3,4]),[6*nc*nc,nr,nt]); |
269 |
if isequal(flu,'A') |
270 |
q6t=split_C_cub(q); |
271 |
q6bi([1:nc],:,:,:,:) = (q6t([1:nc],:,:,:,:)+q6t([2:nc+1],:,:,:,:))./2; |
272 |
q6bj(:,[1:nc],:,:,:) = (q6t(:,[1:nc],:,:,:)+q6t(:,[2:nc+1],:,:,:))./2; |
273 |
qbi = q6bi([1:nc],[2:nc+1],:,:,:); |
274 |
qbj = q6bj([2:nc+1],[1:nc],:,:,:); |
275 |
qbi=reshape(permute(qbi,[1,5,2,3,4]),[6*nc*nc,nr,nt]); |
276 |
qbj=reshape(permute(qbj,[1,5,2,3,4]),[6*nc*nc,nr,nt]); |
277 |
|
278 |
p6t=split_C_cub(p); |
279 |
p6bi([1:nc],:,:,:,:) = (p6t([1:nc],:,:,:,:)+p6t([2:nc+1],:,:,:,:))./2; |
280 |
p6bj(:,[1:nc],:,:,:) = (p6t(:,[1:nc],:,:,:)+p6t(:,[2:nc+1],:,:,:))./2; |
281 |
pbi = p6bi([1:nc],[2:nc+1],:,:,:); |
282 |
pbj = p6bj([2:nc+1],[1:nc],:,:,:); |
283 |
pbi=reshape(permute(pbi,[1,5,2,3,4]),[6*nc*nc,nr,nt]); |
284 |
pbj=reshape(permute(pbj,[1,5,2,3,4]),[6*nc*nc,nr,nt]); |
285 |
end |
286 |
|
287 |
% Prepare output arrays. "nout" is the number of transport output fields. |
288 |
% It is currently hard-coded, but could eventually be an input parameters |
289 |
% to set which output fields are desired if some of then become |
290 |
% computationally expensive. |
291 |
% SenHT = Sensible heat transport |
292 |
% SenSF = Sensible implied surface flux |
293 |
% IntV = Integrated volume transport |
294 |
% IntT = Integrated temperature |
295 |
SenHT = zeros(ydim+2,1+nBas,nt,nout); |
296 |
SenSF = zeros(ydim+1,1+nBas,nt,nout); |
297 |
IntV = zeros(ydim,nr,1+nBas,nt); |
298 |
IntT = zeros(ydim,nr,1+nBas,nt); |
299 |
IntM = zeros(ydim,nr,1+nBas,nt); |
300 |
if isequal(flu,'A') |
301 |
LatHT = zeros(ydim+2,1+nBas,nt,nout); |
302 |
PotHT = zeros(ydim+2,1+nBas,nt,nout); |
303 |
LatSF = zeros(ydim+1,1+nBas,nt,nout); |
304 |
IntQ = zeros(ydim,nr,1+nBas,nt); |
305 |
IntP = zeros(ydim,nr,1+nBas,nt); |
306 |
end |
307 |
Psi = zeros(ydim+2,1+nBas,nt); |
308 |
|
309 |
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
310 |
% Make heat transport calculations % |
311 |
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
312 |
|
313 |
% Preparation for calculation of zonal average temperature. The |
314 |
% tempereature multiplied by the appropriate length scale ("tbi_temp", |
315 |
% "tbj_temp") is summed up ("IntT" in the next section) and |
316 |
% divided by the total length ("IntM", composed from summing "hw_temp", |
317 |
% "hs_temp"). |
318 |
hw_temp = zeros(size(hw)); |
319 |
hs_temp = zeros(size(hs)); |
320 |
for k=1:nr, |
321 |
hw_temp(:,k) = dyg.*hw(:,k); |
322 |
hs_temp(:,k) = dxg.*hs(:,k); |
323 |
end |
324 |
|
325 |
for it = 1:length(i_time) |
326 |
|
327 |
% uz / vz = Volume transport though cell faces (velocity times area). |
328 |
% Used for zonal mean volume transport (4). |
329 |
% utz1/vtz1 = Eulerian sensible heat transport through cell faces (1). |
330 |
% utz2/vtz2 = Sensible heat transport through cell faces by Eulerian |
331 |
% mean circulations (2). |
332 |
% dtx1/dty1 = Temperatude gradient at cell face times the area (when |
333 |
% multiplied by the diffusion will be the horizontal |
334 |
% diffusion heat transport) (6). |
335 |
uz = ArW.*hu(:,:,it); |
336 |
vz = ArS.*hv(:,:,it); |
337 |
utz1 = sum(ArW.*ut(:,:,it),2); |
338 |
vtz1 = sum(ArS.*vt(:,:,it),2); |
339 |
utz2 = sum(ArW.*hu(:,:,it).*tbi(:,:,it),2); |
340 |
vtz2 = sum(ArS.*hv(:,:,it).*tbj(:,:,it),2); |
341 |
if isequal(flu,'A') |
342 |
uqz1 = sum(ArW.*uq(:,:,it),2); |
343 |
vqz1 = sum(ArS.*vq(:,:,it),2); |
344 |
uqz2 = sum(ArW.*hu(:,:,it).*qbi(:,:,it),2); |
345 |
vqz2 = sum(ArS.*hv(:,:,it).*qbj(:,:,it),2); |
346 |
upz1 = sum(ArW.*(up(:,:,it)+hu(:,:,it).*repmat(PHIref,[6*nc*nc 1])),2); |
347 |
vpz1 = sum(ArS.*(vp(:,:,it)+hv(:,:,it).*repmat(PHIref,[6*nc*nc 1])),2); |
348 |
upz2 = sum(ArW.*hu(:,:,it).*(pbi(:,:,it)+repmat(PHIref,[6*nc*nc 1])),2); |
349 |
vpz2 = sum(ArS.*hv(:,:,it).*(pbj(:,:,it)+repmat(PHIref,[6*nc*nc 1])),2); |
350 |
end |
351 |
if isequal(flu,'O') |
352 |
dtx1 = sum(ArW_Dif.*Kux(:,:,it).*tdi(:,:,it),2); |
353 |
dty1 = sum(ArS_Dif.*Kvy(:,:,it).*tdj(:,:,it),2); |
354 |
Kuztz1 = sum(UbTgm(:,:,it),2); |
355 |
Kvztz1 = sum(VbTgm(:,:,it),2); |
356 |
end |
357 |
|
358 |
% Preparation for calculation of zonal average temperature. The |
359 |
% temperature multiplied by the appropriate length scale ("tbi_temp", |
360 |
% "tbj_temp") is summed up ("IntT" in the next section) and |
361 |
% divided by the total length ("IntM", composed from summing "hw_temp", |
362 |
% "hs_temp"). |
363 |
tbi_temp = hw_temp.*tbi(:,:,it); |
364 |
tbj_temp = hs_temp.*tbj(:,:,it); |
365 |
if isequal(flu,'A') |
366 |
qbi_temp = hw_temp.*qbi(:,:,it); |
367 |
qbj_temp = hs_temp.*qbj(:,:,it); |
368 |
pbi_temp = hw_temp.*(pbi(:,:,it)+repmat(PHIref,[6*nc*nc 1])); |
369 |
pbj_temp = hs_temp.*(pbj(:,:,it)+repmat(PHIref,[6*nc*nc 1])); |
370 |
end |
371 |
|
372 |
% Block 1: |
373 |
% With the vertical integral of heat transport calculated across cell |
374 |
% edges, the zonal integral (along the broken line) is computed to |
375 |
% determine the total northward heat transport. The first for loop is |
376 |
% over the individual broken lines (determining the northward HT at the |
377 |
% representative latitude). The second loop is over the basins. The |
378 |
% third loop is over the individual edges along a specific broken line. |
379 |
% Note that an individual cell can have two edges (u and v) that have a |
380 |
% HT contributions. Hence the variable containing indecies of cells |
381 |
% with edges along the broken lines, "bkl_IJuv", has some repeats. |
382 |
% Note that the variable "bkl_Npts" is the number of edges along a |
383 |
% given broken line. Note also that the latitude axis starts at 2 |
384 |
% because heat transport at extremes (latitute = -90/90) is zero by |
385 |
% definition. Recall that index (1) is total Eulerian transport, (2) |
386 |
% is from the mean circulations, and (3) is from the horizontal |
387 |
% diffusion. |
388 |
% |
389 |
% Block 2: |
390 |
% Here zonal average circulation and temperature / moisture is |
391 |
% calculated. The zonal average volume transport v (IntV) and t |
392 |
% (IntT/IntM) are computed first and are multiplied together at the |
393 |
% end. |
394 |
for jl=1:ydim |
395 |
ie=bkl_Npts(jl); |
396 |
for b=1:1+nBas |
397 |
ij=bkl_IJuv(1:ie,jl); |
398 |
% Block 1: |
399 |
SenHT(1+jl,b,it,1) = SenHT(1+jl,b,it,1) + ... |
400 |
sum(ufac(1:ie,jl,b).*utz1(ij) + ... |
401 |
vfac(1:ie,jl,b).*vtz1(ij)); |
402 |
SenHT(1+jl,b,it,2) = SenHT(1+jl,b,it,2) + ... |
403 |
sum(ufac(1:ie,jl,b).*utz2(ij) + ... |
404 |
vfac(1:ie,jl,b).*vtz2(ij)); |
405 |
if isequal(flu,'A') |
406 |
LatHT(1+jl,b,it,1) = LatHT(1+jl,b,it,1) + ... |
407 |
sum(ufac(1:ie,jl,b).*uqz1(ij) + ... |
408 |
vfac(1:ie,jl,b).*vqz1(ij)); |
409 |
LatHT(1+jl,b,it,2) = LatHT(1+jl,b,it,2) + ... |
410 |
sum(ufac(1:ie,jl,b).*uqz2(ij) + ... |
411 |
vfac(1:ie,jl,b).*vqz2(ij)); |
412 |
PotHT(1+jl,b,it,1) = PotHT(1+jl,b,it,1) + ... |
413 |
sum(ufac(1:ie,jl,b).*upz1(ij) + ... |
414 |
vfac(1:ie,jl,b).*vpz1(ij)); |
415 |
PotHT(1+jl,b,it,2) = PotHT(1+jl,b,it,2) + ... |
416 |
sum(ufac(1:ie,jl,b).*upz2(ij) + ... |
417 |
vfac(1:ie,jl,b).*vpz2(ij)); |
418 |
end |
419 |
if isequal(flu,'O') |
420 |
SenHT(1+jl,b,it,6) = SenHT(1+jl,b,it,6) + ... |
421 |
sum(ufac(1:ie,jl,b).*dtx1(ij) + ... |
422 |
vfac(1:ie,jl,b).*dty1(ij)); |
423 |
SenHT(1+jl,b,it,7) = SenHT(1+jl,b,it,7) + ... |
424 |
sum(ufac(1:ie,jl,b).*Kuztz1(ij) + ... |
425 |
vfac(1:ie,jl,b).*Kvztz1(ij)); |
426 |
end |
427 |
% Block 2: |
428 |
IntV(jl,:,b,it) = IntV(jl,:,b,it) ... |
429 |
+ ufac(1:ie,jl,b)'*uz(ij,:) ... |
430 |
+ vfac(1:ie,jl,b)'*vz(ij,:) ; |
431 |
IntT(jl,:,b,it) = IntT(jl,:,b,it) ... |
432 |
+ ufacabs(1:ie,jl,b)'*tbi_temp(ij,:) ... |
433 |
+ vfacabs(1:ie,jl,b)'*tbj_temp(ij,:); |
434 |
IntM(jl,:,b,it) = IntM(jl,:,b,it) ... |
435 |
+ ufacabs(1:ie,jl,b)'*hw_temp(ij,:) ... |
436 |
+ vfacabs(1:ie,jl,b)'*hs_temp(ij,:); |
437 |
if isequal(flu,'A') |
438 |
IntQ(jl,:,b,it) = IntQ(jl,:,b,it) ... |
439 |
+ ufacabs(1:ie,jl,b)'*qbi_temp(ij,:) ... |
440 |
+ vfacabs(1:ie,jl,b)'*qbj_temp(ij,:); |
441 |
IntP(jl,:,b,it) = IntP(jl,:,b,it) ... |
442 |
+ ufacabs(1:ie,jl,b)'*pbi_temp(ij,:) ... |
443 |
+ vfacabs(1:ie,jl,b)'*pbj_temp(ij,:); |
444 |
end |
445 |
end |
446 |
end |
447 |
|
448 |
% Prepare HT output: Calculate HT by zonal flows and tabulate |
449 |
% residuals. Also, the multiplicative constants (Cp,rho,grav,LhVap) are |
450 |
% applied here to put moisture and potential temperature fluxes in |
451 |
% terms of heat transports. |
452 |
tmp = IntV(:,:,:,it) .* IntT(:,:,:,it) ... |
453 |
./ change(IntM(:,:,:,it),'==',0,NaN); |
454 |
SenHT(2:ydim+1,:,it,4) = sum(change(tmp,'==',NaN,0),2); |
455 |
% SenHT(2:ydim+1,:,it,4) = sum( IntV(:,:,:,it) ... |
456 |
% .* IntT(:,:,:,it) ... |
457 |
% ./ IntM(:,:,:,it),2); |
458 |
SenHT(:,:,it,3) = SenHT(:,:,it,1) - SenHT(:,:,it,2); |
459 |
SenHT(:,:,it,5) = SenHT(:,:,it,2) - SenHT(:,:,it,4); |
460 |
if isequal(flu,'O') |
461 |
% SenHT(:,:,it,6) = DiffKh.*SenHT(:,:,it,6); |
462 |
SenHT(:,:,it,:) = (CpO*RhoO)*SenHT(:,:,it,:); |
463 |
if flag_tot==1 |
464 |
SenHT(:,:,it,7)=SenHT(:,:,it,7)-SenHT(:,:,it,1); |
465 |
end |
466 |
elseif isequal(flu,'A') |
467 |
SenHT(:,:,it,:) = (CpA./grav).*SenHT(:,:,it,:); |
468 |
LatHT(2:ydim+1,:,it,4) = sum( IntV(:,:,:,it) ... |
469 |
.* IntQ(:,:,:,it) ... |
470 |
./ IntM(:,:,:,it),2); |
471 |
LatHT(:,:,it,3) = LatHT(:,:,it,1) - LatHT(:,:,it,2); |
472 |
LatHT(:,:,it,5) = LatHT(:,:,it,2) - LatHT(:,:,it,4); |
473 |
LatHT(:,:,it,:) = (LhVap./grav).*LatHT(:,:,it,:); |
474 |
|
475 |
PotHT(2:ydim+1,:,it,4) = sum( IntV(:,:,:,it) ... |
476 |
.* IntP(:,:,:,it) ... |
477 |
./ IntM(:,:,:,it),2); |
478 |
PotHT(:,:,it,3) = PotHT(:,:,it,1) - PotHT(:,:,it,2); |
479 |
PotHT(:,:,it,5) = PotHT(:,:,it,2) - PotHT(:,:,it,4); |
480 |
PotHT(:,:,it,:) = 1/grav * PotHT(:,:,it,:); |
481 |
end |
482 |
Psi(2:ydim+1,:,it) = - sum( IntV(:,:,:,it) ,2 ); |
483 |
|
484 |
% Implied surface heat flux from heat transports (implied heating). |
485 |
% Tabulated as the difference in heat transports between two broken |
486 |
% lines divided by the zonal band area. |
487 |
mskG = reshape(mskG,(ydim+1)*(1+nBas),1); |
488 |
I = find(mskG==0); |
489 |
mskG = reshape(mskG,ydim+1,1+nBas); |
490 |
var = zeros(ydim+1,1+nBas); |
491 |
for n=1:min(nout,6), |
492 |
varT = SenHT([2:ydim+2],:,it,n) ... |
493 |
- SenHT([1:ydim+1],:,it,n); |
494 |
varT = reshape(varT,(ydim+1)*(1+nBas),1); varT(I)=NaN; |
495 |
varT = reshape(varT,ydim+1,1+nBas); |
496 |
SenSF([1:ydim+1],:,it,n) = varT./AreaZon; |
497 |
if isequal(flu,'A') |
498 |
for n=1:min(nout,6) |
499 |
varQ = LatHT([2:ydim+2],:,it,n) ... |
500 |
- LatHT([1:ydim+1],:,it,n); |
501 |
varQ = reshape(varQ,(ydim+1)*(1+nBas),1); varQ(I)=NaN; |
502 |
varQ = reshape(varQ,ydim+1,1+nBas); |
503 |
LatSF([1:ydim+1],:,it,n) = varQ./AreaZon; |
504 |
end |
505 |
end |
506 |
end |
507 |
end |
508 |
|
509 |
|
510 |
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
511 |
% Assign outputs, put in units of PW % |
512 |
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
513 |
|
514 |
SenHT = SenHT*1e-15; |
515 |
|
516 |
HT.time = i_time; |
517 |
HT.SenHT = reshape(SenHT,ydim+2,nt,nout); |
518 |
HT.SenSF = reshape(SenSF,ydim+1,nt,nout); |
519 |
HT.ylatHT = ylatHT; |
520 |
HT.ylatSF = ylatSF; |
521 |
HT.AreaZon = AreaZon; |
522 |
HT.PsiSurf = 1e-6*reshape(Psi,ydim+2,nt); |
523 |
|
524 |
if isequal(flu,'A') |
525 |
LatHT = LatHT*1e-15; |
526 |
PotHT = PotHT*1e-15; |
527 |
HT.LatHT = reshape(LatHT,ydim+2,nt,nout); |
528 |
HT.LatSF = reshape(LatSF,ydim+1,nt,nout); |
529 |
HT.PotHT = reshape(PotHT,ydim+2,nt,nout); |
530 |
end |
531 |
|