source: issm/trunk/src/m/solutions/macayeal/control.m@ 408

Last change on this file since 408 was 408, checked in by Mathieu Morlighem, 16 years ago

minor

File size: 12.1 KB
Line 
1function md=control(md)
2%CONTROL - launch a control method using MacAyeal solution
3%
4% the routine is used for a control method. It determines the most adapted viscosity
5% field so that the calculated velocity field is as close as possible to an observed velocity field
6%
7% Usage:
8% md=control(md)
9
10%First check we do have the correct argument number
11if ((nargin~=1) || (nargout~=1)),
12 velfinderusage();
13 error('macayealcontrol error message: incorrect number of input and output arguments');
14end
15
16if ~isa(md,'model'),
17 macayealcontrolusage();
18 error('macayealcontrol error message: input argument is not a @model object');
19end
20
21%Check that control is done on flow law, and not drag (not supported yet):
22%if ~strcmp(md.control_type,'B')
23% error('macayealcontrol error message: only ''B'' inversion supported yet');
24%end
25
26%Transfer model fields into matlab variables
27x=md.x;
28y=md.y;
29index=md.elements;
30index=sort(index,2); %necessary
31nods=md.numberofgrids;
32nel=md.numberofelements;
33z_thick=md.thickness;
34z_surf=md.surface;
35z_bed=md.bed;
36z_thick_bar=(z_thick(index(:,1))+z_thick(index(:,2))+z_thick(index(:,3)))/3;
37rho_ice=md.rho_ice;
38rho_water=md.rho_water;
39g=md.g;
40index_icefront=md.segmentonneumann_diag; index_icefront=index_icefront(:,1:2); %we strip the last column, which holds the element number for the boundary segment
41nodes_on_boundary=md.gridonboundary;
42nodes_on_dirichlet=md.gridondirichlet_diag;
43nodes_on_icefront=zeros(nods,1); nodes_on_icefront(index_icefront)=1;
44nodes_on_iceshelf=md.gridoniceshelf;
45
46criterion=md.eps_rel;
47yts=md.yts;
48tolx=md.tolx;
49maxiter=md.maxiter(1);
50if strcmp(md.control_type,'B')
51 B_ini=md.B;
52 B_to_p=10^-3*yts^(-1/3);
53else
54 drag_coeff_ini=md.drag;
55end
56glen_coeff=md.n;
57vx_obs=md.vx_obs/md.yts; %From m/a to m/s
58vy_obs=md.vy_obs/md.yts;
59nsteps=md.nsteps;
60
61%Build length_icefront and normal_icefront:
62[length_icefront,normal_icefront]=buildicefrontnormal(x,y,index_icefront);
63
64%Building shape functions and derivative operators
65[alpha, beta, gamma, area]=shape(index,x,y,nel,nods);
66
67[matrix_bar, matrix_xbar, matrix_ybar]=...
68 bar_maker(nel,nods,index,alpha,beta);
69
70%initialize some data
71create_el2nod_matrices
72old_direction=zeros(nods,1);
73
74%setup some fake distributions.
75z_thick_bar=matrix_bar*z_thick;
76z_surf_bar=matrix_bar*z_surf;
77z_surf_xbar=matrix_xbar*z_surf;
78z_surf_ybar=matrix_ybar*z_surf;
79z_surf_x=el2nod\(el2nodRhs*z_surf_xbar);
80z_surf_y=el2nod\(el2nodRhs*z_surf_ybar);
81z_bed_bar=matrix_bar*z_bed;
82z_bed_xbar=matrix_xbar*z_bed;
83z_bed_ybar=matrix_ybar*z_bed;
84
85%data_elements and data_nodes are the elements and nodes where
86%we have observations
87data_elements=[1:nel]';
88data_nodes=ones(nods,1);
89
90%initialize misfit between model velocity and observed velocity
91J=zeros(nsteps,1);
92
93%Weighting: one areas where we don't want the control method to optimize
94%the misfit, we can set weighting to 0.
95weighting=ones(nods,1);
96
97%optimization parameters for the matlab functoin fminbnd
98options=optimset('fminbnd');
99options=optimset(options,'Display','iter');
100options=optimset(options,'MaxFunEvals',maxiter);
101options=optimset(options,'MaxIter',100);
102options=optimset(options,'TolX',tolx);
103
104%build useful matrices.
105setup_control
106
107%setup initial drag paramter distribution to startup the optimization.
108drag_type=md.drag_type;
109qcoeff=md.q;
110pcoeff=md.p;
111
112if strcmp(md.control_type,'B')
113 %setup initial flow law paramter distribution to startup the optimization.
114 B_bar_ini=(B_ini(index(:,1))+B_ini(index(:,2))+B_ini(index(:,3)))/3;
115 B=B_ini; %variables used in optimziation are B and B_bar.
116 B_bar=B_bar_ini;
117else
118 B_bar=(md.B(index(:,1))+md.B(index(:,2))+md.B(index(:,3)))/3;
119 if (drag_type~=2), error('md.drag_type must be 2 for control methods'); end
120 drag_coeff_bar_ini=(drag_coeff_ini(index(:,1))+drag_coeff_ini(index(:,2))+drag_coeff_ini(index(:,3)))/3;
121 drag_coeff=drag_coeff_ini; %variables used in optimziation are drag_coeff and drag_coeff_bar.
122 drag_coeff_bar=drag_coeff_bar_ini;
123end
124
125
126
127%check that the model is not a pure ice shelf (no friction on ice shelves)
128if strcmp(md.control_type,'drag') & (length(find(drag_coeff))==0),
129 disp(sprintf('\n No drag pure for ice shelves => STOP'))
130 return
131end
132
133%nu_bar=10^14*ones(nel,1); %initial element viscosity distribution.
134nu_bar=viscosity(index,nel,alpha,beta,[],[],B_bar,glen_coeff);
135
136%%%AK velfinder; %forward model that determines first velocity field used
137%to start optimization.
138disp('calculating the velocity with the initial parameters')
139c_velfinder;
140
141if strcmp(md.control_type,'B'),
142 for niteration=1:nsteps,
143
144 disp([' Step #' num2str(niteration) ' on ' num2str(nsteps)]);
145
146 % Compute search direction -
147 [u,v,adjointu,adjointv,direction]= ...
148 grad_J_flow(Rhs,S,F,P,P0,area,specified_velocity,nods,vx_obs,vy_obs,...
149 index,x,y,nel,rho_ice,g,weighting,alpha,beta,z_thick_bar,B_bar);
150
151 if md.plot==1,
152 if niteration==1
153 scrsz = get(0,'ScreenSize');
154 figure('Position',[scrsz(3)/3 scrsz(4)*1.8/3 scrsz(3)/3 scrsz(4)/3])
155 end
156 plotmodel(md,'data',sqrt(u.^2+v.^2)*yts,'title','Modeled velocity',...
157 'data',sqrt(adjointu.^2+adjointv.^2)*yts,'title','Adjoint vectors','colorbar#all', 'on',...
158 'data',sqrt(vx_obs.^2+vy_obs.^2)*yts,'title','Observed velocity',...
159 'data',100*abs(sqrt(vx_obs.^2+vy_obs.^2)-sqrt(u.^2+v.^2))./sqrt(vx_obs.^2+vy_obs.^2),'title','Relative misfit','caxis#4',[0 100],'figure',1);drawnow;
160 end
161
162 %Keep track of u and v to use in objectivefunction_C:
163 u_objective=u;
164 v_objective=v;
165
166 %Orthogonalize direction
167 direction=real(direction);
168 direction=direction/sqrt(direction'*direction);
169
170
171 % rough orthagonalization
172 direction=direction-(direction'*old_direction)*old_direction;
173 old_direction=direction;
174
175 %visualize direction
176 if md.plot==1,
177 if niteration==1
178 scrsz = get(0,'ScreenSize');
179 figure('Position',[10 scrsz(4)*1/3 scrsz(3)/3 scrsz(4)/3])
180 end
181 plotmodel(md,'data',direction,'title','Orthogonal direction','figure',2);drawnow;
182 end
183
184
185 % normalize direction to 10^7, so that when variations on B are computed
186 %they will be significant.
187 if abs(max(direction))>abs(min(direction))
188 direction=10^7*direction/abs(max(direction));
189 else
190 direction=10^7*direction/abs(min(direction));
191 end
192
193 %during optimization, bounds on B variations can vary. Here, they are less
194 %strict in the first iteration.
195 if niteration<=2,
196 upperbound=20;
197 lowerbound=0;
198 else
199 upperbound=10;
200 lowerbound=0;
201 end
202
203 %search the multiplicative constant to the direction, that will
204 %be used to modify B.
205
206 search_constant=fminbnd('objectivefunction_C_flow',lowerbound,upperbound, ...
207 options,B,glen_coeff,direction,Rhs,S,F,P,specified_velocity,nods,nel,vx_obs, ...
208 vy_obs,u_objective, v_objective, index,alpha, beta, area,weighting,rowD,colD,rowDshort,colDshort,valueuu,valuevv,valueuv,matrix_bar,criterion);
209
210
211 %update value of B
212 B_old=B;
213 B_new=B+search_constant*direction;
214 B=B_new;
215 pos=find(B<0);B(pos)=-B(pos);
216
217 %Average B over elements:
218 B_bar=(B(index(:,1))+B(index(:,2))+B(index(:,3)))/3;
219
220 %visualize new distribution.
221 if md.plot==1,
222 if niteration==1
223 scrsz = get(0,'ScreenSize');
224 figure('Position',[scrsz(3)*1.97/3 scrsz(4)*1/3 scrsz(3)/3 scrsz(4)/3])
225 end
226 plotmodel(md,'data',B*B_to_p,'title',['B at iteration ' num2str(niteration)],'caxis',[200 900],'colorbar','on','figure',3);drawnow;
227 end
228
229 %evaluate new misfit.
230 %@@@AK load F_file %this file was created in objectivefunction_C, and is reloaded
231 %here to win some computation time.
232 load F_file
233
234 J(niteration)=objectivefunction_C_flow(0,B,glen_coeff,direction,Rhs,S,F,P,specified_velocity,nods,nel,vx_obs, ...
235 vy_obs,u_objective, v_objective, index,alpha, beta, area,weighting,rowD,colD,rowDshort,colDshort,valueuu,valuevv,valueuv,matrix_bar,criterion);
236 disp(J(niteration));
237
238 %do a backup every 5 iterations
239 if(mod(niteration,5)==0),
240 save temporary_control_results_flow B B_bar u v J direction
241 end
242 end
243
244 %Load results onto md:
245 md.cont_J=J;
246 md.cont_parameter=B;
247
248elseif strcmp(md.control_type,'drag'),
249 for niteration=1:nsteps,
250
251 disp([' Step #' num2str(niteration) ' on ' num2str(nsteps)]);
252
253 % Compute search direction -
254 [u,v,adjointu,adjointv,direction]= ...
255 grad_J_drag(Rhs,S,F,P,P0,area,specified_velocity,nods,vx_obs,vx_obs,...
256 index,x,y,nel,rho_ice,g,weighting,alpha,beta,z_thick_bar,drag_coeff,drag_coeff_bar);
257
258 if md.plot==1,
259 if niteration==1
260 scrsz = get(0,'ScreenSize');
261 figure('Position',[scrsz(3)/3 scrsz(4)*1.8/3 scrsz(3)/3 scrsz(4)/3])
262 end
263 plotmodel(md,'data',sqrt(u.^2+v.^2),'title','Modeled velocity',...
264 'data',sqrt(adjointu.^2+adjointv.^2),'title','Adjoint vectors','colorbar', 'all',...
265 'data',sqrt(vx_obs.^2+vy_obs.^2),'title','Observed velocity',...
266 'data',100*abs(sqrt(vx_obs.^2+vy_obs.^2)-sqrt(u.^2+v.^2))./sqrt(vx_obs.^2+vy_obs.^2),'title','Relative misfit','caxis#3',[0 100],'figure',1);drawnow;
267 end
268
269 %Keep track of u and v to use in objectivefunction_C:
270 u_objective=u;
271 v_objective=v;
272
273 %Orthogonalize direction
274 direction=real(direction);
275 direction=direction/sqrt(direction'*direction);
276 pos=find(nodes_on_iceshelf);
277 direction(pos)=0;
278
279
280 % rough orthagonalization
281 direction=direction-(direction'*old_direction)*old_direction;
282 old_direction=direction;
283
284 %visualize direction
285 if md.plot==1,
286 if niteration==1
287 scrsz = get(0,'ScreenSize');
288 figure('Position',[10 scrsz(4)*1/3 scrsz(3)/3 scrsz(4)/3])
289 end
290 plotmodel(md,'data',direction,'title','Orthogonal direction','figure',2);drawnow;
291 end
292
293
294 % normalize direction to 50, so that when variations on drag are computed
295 %they will be significant.
296 direction=50*direction/max(abs(direction));
297
298 %during optimization, bounds on drag variations can vary. Here, they are less
299 %strict in the first iteration.
300 if niteration<=2,
301 upperbound=2;
302 lowerbound=0;
303 else
304 upperbound=1;
305 lowerbound=0;
306 end
307
308 %search the multiplicative constant to the direction, that will
309 %be used to modify drag.
310
311 search_constant=fminbnd('objectivefunction_C_drag',lowerbound,upperbound, ...
312 options,B,glen_coeff,drag_coeff,direction,Rhs,S,F,P,specified_velocity,nods,nel,vx_obs, ...
313 vy_obs,u_objective, v_objective, index,alpha, beta, area,weighting,rowD,colD,rowDshort,colDshort,valueuu,valuevv,valueuv,matrix_bar,criterion);
314
315 %update value of drag
316 drag_coeff_old=drag_coeff;
317 drag_coeff_new=drag_coeff+search_constant*direction;
318 drag_coeff=drag_coeff_new;
319
320 if length(find(drag_coeff<0))~=0,
321 disp(sprintf('\n Some basal drag coefficient negative => STOP'))
322 break
323 end
324
325 %Average drag over elements:
326 if niteration==1
327 scrsz = get(0,'ScreenSize');
328 figure('Position',[scrsz(3)*1.97/3 scrsz(4)*1/3 scrsz(3)/3 scrsz(4)/3])
329 end
330 drag_coeff_bar=(drag_coeff(index(:,1))+drag_coeff(index(:,2))+drag_coeff(index(:,3)))/3;
331
332 %visualize new distribution.
333 if md.plot==1,
334 plotmodel(md,'data',drag_coeff,'title',['Drag at iteration ' num2str(niteration)],'caxis',[0 1000],'colorbar','on','figure',3);drawnow;
335 end
336
337 %evaluate new misfit.
338 %@@@AK load F_file %this file was created in objectivefunction_C, and is reloaded
339 %here to win some computation time.
340 load F_file
341
342 J(niteration)=objectivefunction_C_drag(0,B,glen_coeff,drag_coeff,direction,Rhs,S,F,P,specified_velocity,nods,nel,vx_obs, ...
343 vy_obs,u_objective, v_objective, index,alpha, beta, area,weighting,rowD,colD,rowDshort,colDshort,valueuu,valuevv,valueuv,matrix_bar,criterion);
344 disp(J(niteration));
345
346 %do a backup every 5 iterations
347 if(mod(niteration,5)==0),
348 save temporary_control_results_drag drag_coeff drag_coeff_bar u v J direction
349 end
350 end
351end
352
353function macayealcontrolusage();
354disp('md=macayealcontrol(md)');
355disp(' where md is a structure of class @model');
Note: See TracBrowser for help on using the repository browser.