/[MITgcm]/MITgcm_contrib/dfer/matlab_stuff/calcBolusVelCube.m
ViewVC logotype

Contents of /MITgcm_contrib/dfer/matlab_stuff/calcBolusVelCube.m

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.2 - (show annotations) (download)
Wed Mar 7 22:01:18 2018 UTC (7 years, 4 months ago) by dfer
Branch: MAIN
Changes since 1.1: +48 -10 lines
Update with various little adjustments.

1 function [ub,vb,wb]=calcBolusVelCube(d,g,GMform);
2
3 % [ub,vb] = calcBolusVelCube(d,g,GMform);
4 %
5 % Input arguments:
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 [Field data] Kwx,Kwy
9 % g [Grid data ] drF,rA,dxC,dyC,dxG,dyG,HFacW,HFacS
10 % GMform [string] GM form, 'Skew' or 'Advc'
11 %
12 % Output arguments:
13 % ub, vb: GM-Bolus mass-weigthed velocity (i.e include
14 % implicitly hFac factor)
15 %
16 % Comments:
17 % For Skew-Flux form: uses Kwx & Kwy divided by 2
18 % compute Volume Stream function psiX,psiY above uVel.vVel
19 % (at interface between 2 levels), units=m^3/s :
20 % psiX=(rAc*kwx)_i / dXc ; psiY=(rAc*kwy)_j / dYc ;
21 % and then the bolus velocity (m/s):
22 % ub = d_k(psiX)/rAw/drF ; vb = d_k(psiY)/rAs/drF ;
23 %
24 %---------------------------------------------------------------------
25
26
27 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
28 % Prepare / reform incoming data %
29 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
30
31 nc = size(g.XC,2);
32 nr = length(g.drF);
33
34 switch GMform
35 case 'Skew'
36 nt = size(d.GM_Kwx,4);
37 case 'Advc'
38 nt = size(d.GM_PsiX,4);
39 end
40
41 dr = g.drF;
42 hw = reshape(g.HFacW(1:6*nc,1:nc,1:nr),[6*nc*nc,nr]);
43 hs = reshape(g.HFacS(1:6*nc,1:nc,1:nr),[6*nc*nc,nr]);
44 ra = reshape(g.rA(1:6*nc,1:nc) ,[6*nc*nc,1]);
45 rAs = reshape(g.rAs(1:6*nc,1:nc) ,[6*nc*nc,1]);
46 rAw = reshape(g.rAw(1:6*nc,1:nc) ,[6*nc*nc,1]);
47 dxc = reshape(g.dxC(1:6*nc,1:nc),[6*nc*nc,1]);
48 dyc = reshape(g.dyC(1:6*nc,1:nc),[6*nc*nc,1]);
49 dxg = reshape(g.dxG(1:6*nc,1:nc),[6*nc*nc,1]);
50 dyg = reshape(g.dyG(1:6*nc,1:nc),[6*nc*nc,1]);
51
52 %rAw=dxc.*dyg;
53 %rAs=dyc.*dxg;
54
55 %--- recip_hFac & mask :
56 mw=ceil(hw); mw=min(1,mw);
57 ms=ceil(hs); ms=min(1,ms);
58
59 hw(find(hw==0))=Inf;
60 hs(find(hs==0))=Inf;
61 hw_recip=1./hw; %hw_recip(find(hw==0))=0;
62 hs_recip=1./hs; %hs_recip(find(hs==0))=0;
63
64 ub_all = zeros(6*nc,nc,nr,nt);
65 vb_all = zeros(6*nc,nc,nr,nt);
66 wb_all = zeros(6*nc,nc,nr,nt);
67
68 switch GMform
69
70 %%%%%%% Skew-flux form case
71 case 'Skew'
72 kwx_all = reshape(d.GM_Kwx,[6*nc*nc,nr,nt]);
73 kwy_all = reshape(d.GM_Kwy,[6*nc*nc,nr,nt]);
74
75 kwx_all = 0.5*kwx_all;
76 kwy_all = 0.5*kwy_all;
77
78 for it = 1:nt
79 kwx = kwx_all(:,:,it);
80 kwy = kwy_all(:,:,it);
81
82 %-- K*ra + add 1 overlap :
83 kwx = (ra*ones(1,nr)).*kwx;
84 kwy = (ra*ones(1,nr)).*kwy;
85 kwx = reshape(kwx,[6*nc,nc,nr]);
86 kwy = reshape(kwy,[6*nc,nc,nr]);
87 v6X = split_C_cub(kwx,1);
88 v6Y = split_C_cub(kwy,1);
89 k6x = v6X(:,[2:nc+1],:,:);
90 k6y = v6Y([2:nc+1],:,:,:);
91
92 %-----------------
93 v6X = zeros(nc,nc,nr,6);
94 v6Y = zeros(nc,nc,nr,6);
95
96 v6X([1:nc],:,:,:) = k6x([2:nc+1],:,:,:) + k6x([1:nc],:,:,:);
97 v6Y(:,[1:nc],:,:) = k6y(:,[2:nc+1],:,:) + k6y(:,[1:nc],:,:);
98
99 v6X = v6X/2;
100 v6Y = v6Y/2;
101
102 psiX = zeros(6*nc,nc,nr+1);
103 psiY = zeros(6*nc,nc,nr+1);
104
105 for n = 1:6
106 is = 1+nc*(n-1);ie=nc*n;
107 psiX([is:ie],[1:nc],[1:nr]) = v6X([1:nc],[1:nc],[1:nr],n);
108 psiY([is:ie],[1:nc],[1:nr]) = v6Y([1:nc],[1:nc],[1:nr],n);
109 end
110
111 psiX = reshape(psiX,6*nc*nc,nr+1);
112 psiY = reshape(psiY,6*nc*nc,nr+1);
113
114 psiX(:,[1:nr]) = mw.*psiX(:,[1:nr]);
115 psiY(:,[1:nr]) = ms.*psiY(:,[1:nr]);
116 ub = psiX(:,[2:nr+1]) - psiX(:,[1:nr]);
117 vb = psiY(:,[2:nr+1]) - psiY(:,[1:nr]);
118
119 dr = reshape(dr,[1,length(dr)]);
120 % ub = reshape(hw_recip.*ub./(rAw*dr),[6*nc,nc,nr]);
121 ub = reshape(ub./(rAw*dr),[6*nc,nc,nr]);
122 % vb = reshape(hs_recip.*vb./(rAs*dr),[6*nc,nc,nr]);
123 vb = reshape(vb./(rAs*dr),[6*nc,nc,nr]);
124
125 ub_all(:,:,:,it) = ub;
126 vb_all(:,:,:,it) = vb;
127
128 %%%%%%%%%%%%%
129 [u6t,v6t] =split_UV_cub(ub,vb,0,1);
130 [hw6t,hs6t]=split_UV_cub(reshape(hw,6*nc,nc,nr),reshape(hs,6*nc,nc,nr),0,1);
131 [dy6t,dx6t]=split_UV_cub(reshape(dyg,6*nc,nc),reshape(dxg,6*nc,nc),0,1);
132
133 %F6tX = u6t.*hw6t.*permute(repmat(dy6t,[1 1 1 nr]),[1 2 4 3]);
134 %F6tY = v6t.*hs6t.*permute(repmat(dx6t,[1 1 1 nr]),[1 2 4 3]);
135 F6tX = u6t.*permute(repmat(dy6t,[1 1 1 nr]),[1 2 4 3]);
136 F6tY = v6t.*permute(repmat(dx6t,[1 1 1 nr]),[1 2 4 3]);
137
138 Hdiv = zeros(nc,nc,nr,6);
139 Hdiv = F6tX([2:nc+1],:,:,:) - F6tX([1:nc],:,:,:) ...
140 + F6tY(:,[2:nc+1],:,:) - F6tY(:,[1:nc],:,:);
141 for k=1:nr
142 Hdiv(:,:,k,:) = -Hdiv(:,:,k,:)*dr(k);
143 end
144
145 %psiX = zeros(6*nc,nc,nr);
146 %for n = 1:6
147 % is = 1+nc*(n-1);ie=nc*n;
148 % psiX([is:ie],[1:nc],[1:nr]) = Hdiv([1:nc],[1:nc],[1:nr],n);
149 %end
150 psiX = reshape(permute(Hdiv,[1 4 2 3]),6*nc,nc,nr);
151 %wb = psiX ./ repmat(reshape(ra,6*nc,nc),[1 1 nr]);
152
153 wb = zeros(6*nc,nc,nr+1);
154 for k=nr:-1:1
155 wb(:,:,k)= psiX(:,:,k) ./ reshape(ra,6*nc,nc) + wb(:,:,k+1);
156 end
157
158 wb_all(:,:,:,it) = wb(:,:,1:30);
159 %%%%%%%%%%%%
160 end
161
162 %%%%%%% Advective form case
163 case 'Advc'
164
165 PsiX_all = reshape(d.GM_PsiX(1:6*nc,1:nc,1:nr,:),[6*nc*nc,nr,nt]);
166 PsiY_all = reshape(d.GM_PsiY(1:6*nc,1:nc,1:nr,:),[6*nc*nc,nr,nt]);
167
168 dr3d = ones(6*nc*nc,1)*reshape(dr,[1,length(dr)]);
169
170 for it = 1:nt
171
172 psiX = zeros(6*nc*nc,nr+1);
173 psiY = zeros(6*nc*nc,nr+1);
174
175 psiX(:,1:nr) = mw.*PsiX_all(:,:,it);
176 psiY(:,1:nr) = ms.*PsiY_all(:,:,it);
177
178 % psiX(:,[1:nr]) = mw.*psiX(:,[1:nr]);
179 % psiY(:,[1:nr]) = ms.*psiY(:,[1:nr]);
180 ub = psiX(:,[2:nr+1]) - psiX(:,[1:nr]);
181 vb = psiY(:,[2:nr+1]) - psiY(:,[1:nr]);
182
183 % ub = reshape(hw_recip.*ub./(rAw*dr),[6*nc,nc,nr]);
184 % ub = reshape(ub./(rAw*dr),[6*nc,nc,nr]);
185 ub = reshape(ub./dr3d,[6*nc,nc,nr]);
186 % vb = reshape(hs_recip.*vb./(rAs*dr),[6*nc,nc,nr]);
187 % vb = reshape(vb./(rAs*dr),[6*nc,nc,nr]);
188 vb = reshape(vb./dr3d,[6*nc,nc,nr]);
189
190 ub_all(:,:,:,it) = ub;
191 vb_all(:,:,:,it) = vb;
192 end
193
194 otherwise
195 disp('Ce portnawak')
196 end
197
198 ub = ub_all;
199 vb = vb_all;
200 wb = wb_all;

  ViewVC Help
Powered by ViewVC 1.1.22