rbmatlab  1.13.10
 All Classes Namespaces Files Functions Variables Groups Pages
advection_fv_output_vconst.m
Go to the documentation of this file.
1 function advection_fv_output_vconst(step)
2 %function advection_fv_output_vconst(step)
3 % small script implementing a simple advection example
4 % for producing matrices to be used in the RB-DS framework
5 % Discretization with FV Functions. The steps generate the
6 % results for the at-Automatisierungstechnik 2010 Paper
7 %
8 % mu_names = {'cone_weight','vx_weight','vy_weight'};
9 % mu_ranges = [0,1]^3
10 %
11 % possible steps:
12 %
13 % step = 1; % initialization of model and plot of model data
14 % step = 2; % detailed simulation of lin_evol
15 % step = 2.5; % conversion to DS primal model, lin_ds detailed simulation
16 % step = 3; % conversion to DS primal model, DS detailed simulation
17 % and comparison with lin_evol
18 % step = 4 % compute several trajectories.
19 % step = 4.5 % compute reduced basis based on trajectories.
20 % step = 4.75 % comparison of detailed and reduced simulation
21 % step = 5 % figures of detailed simulations based on step 4
22 % step = 6 % figures of reduced simulations based on step 4
23 % step = 7 % time measurements based on step 4
24 % step = 8; % compare original and accelerated ds_model
25 % step = 9; % generate greedy basis on subset of parameter domain
26 % step = 10; % reduced simulations and error plots
27 % step = 11; % test of matlab ode solver for resulting sparse system
28 % step = 12; % further plots for at-Automatisierungstechnik paper
29 
30 % Bernard Haasdonk 25.8.2009
31 
32 % structured triangular grid
33 
34 if nargin < 1
35  step = 1;
36 end;
37 
38 params = [];
39 model = [];
40 %params.coarse_factor = 2; % fact 2 => ndofs = 16384, nt = 1024, 15s detailed
41 %params.coarse_factor = 1; % fact 1 => ndofs = 32768, nt = 2048, 60s detailed
42 params.coarse_factor = 1;
43 %verbose(1);
44 trajectory_fnbase = 'advection_fv_output_vconst_trajectory';
45 detailed_data_fname = 'advection_fv_output_vconst_detailed_data';
46 rb_fname = 'advection_fv_output_vconst_rb';
47 
48 switch step
49 
50  %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
51  %step = 1; % initialization of model and plot of model data
52  case 1
53  %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
54 
55 % verbose(1,'initializaton of lin_evol_model and plot model_data');
56  disp('initializaton of lin_evol_model and plot model_data');
57  model = advection_fv_output_vconst_model(params);
58  model_data = gen_model_data(model);
59  params.axis_equal = 1;
60  params.axis_tight = 1;
61  plot(model_data.grid,params);
62  % inspect(model_data.grid)
63 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
64  %step = 2; % detailed simulation of lin_evol
65 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
66  case 2
67 % verbose(1,'detailed simulation of lin_evol_model and plot');
68  disp('detailed simulation of lin_evol_model and plot');
69 % params.coarse_factor = 8;
70  model = advection_fv_output_vconst_model(params);
71  model_data = gen_model_data(model);
72  % model.debug = 1;
73  % disp('model debugging turned on.')
74 % model.cone_weight = 0.5;
75 % model.vx_weight = 0.0;
76 % model.vy_weight = 0.0;
77 % model.vx_weight = 1.0;
78 % model.vy_weight = 1.0;
79 % model.cone_amplitude = 100;
80 % model.vx_weight = 1;
81 % model.vy_weight = 0;
82  sim_data = detailed_simulation(model,model_data);
83  plot_sim_data(model,model_data,sim_data,[]);
84 
85 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
86 %step = 2.5; % conversion to DS primal model, DS detailed simulation
87 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
88  case 2.5
89  disp('initialization of lin_ds model from lin_evol model');
90  % verbose(1,'initialization of lin_ds model from lin_evol model and ...
91  % plot');
92 
93  lin_evol_model = advection_fv_output_vconst_model(params);
94  %lin_evol_model.debug = 1;
95  lin_evol_model_data = gen_model_data(lin_evol_model);
96 
97  % set some additional fields in model which will be copied into
98  % ds model depending on whether constants are known or not:
99  estimate_constants = 0;
100  if estimate_constants
101  %%% if constant estimaton has to be performed:
102  lin_evol_model.estimate_lin_ds_nmu = 4;
103  lin_evol_model.estimate_lin_ds_nX = 10;
104  lin_evol_model.estimate_lin_ds_nt = 4;
105  else
106  %%% otherwise:
107  % the following values are rough bounds, so could be
108  % non-rigorous, then larger search must be performed.
109  lin_evol_model.state_bound_constant_C1 = 1;
110  lin_evol_model.output_bound_constant_C2 = 1.2916;
111  end;
112 
113  ds_model = lin_ds_from_lin_evol_model(lin_evol_model);
114  %ds_model.u_function_ptr = @(params) 0.1; % define new
115 % ds_model.theta = 1;
116  ds_model.theta = 0;
117  ds_model_data = gen_model_data(ds_model);
118  mu = [0.5,1,1];
119  % set mu in model and base_model!!!
120  ds_model = ds_model.set_mu(ds_model,mu);
121 
122  ds_sim_data = detailed_simulation(ds_model,ds_model_data);
123 
124  params.no_lines = 1;
125  params.clim = [0,1];
126  params.plot_title='detailed implicit lin_ds';
127  lin_evol_from_ds_sim_data.U = ds_sim_data.X;
128  lin_evol_from_ds_sim_data.y = ds_sim_data.Y;
129  plot_sim_data(lin_evol_model,lin_evol_model_data,lin_evol_from_ds_sim_data,params);
130 
131  %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
132 %step = 3; % conversion to DS primal model, DS detailed simulation
133  % and comparison with lin_evol
134 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
135 
136  case 3 % initialization of DS model and lin_ds detailed simulation
137  disp('initialization of lin_ds model from lin_evol model and plot');
138 % verbose(1,'initialization of lin_ds model from lin_evol model and ...
139 % plot');
140 
141  lin_evol_model = advection_fv_output_vconst_model(params);
142  %lin_evol_model.debug = 1;
143  lin_evol_model_data = gen_model_data(lin_evol_model);
144 
145  % set some additional fields in model which will be copied into
146  % ds model depending on whether constants are known or not:
147  estimate_constants = 0;
148  if estimate_constants
149  %%% if constant estimaton has to be performed:
150  lin_evol_model.estimate_lin_ds_nmu = 4;
151  lin_evol_model.estimate_lin_ds_nX = 10;
152  lin_evol_model.estimate_lin_ds_nt = 4;
153  else
154  %%% otherwise:
155  % the following values are rough bounds, so could be
156  % non-rigorous, then larger search must be performed.
157  lin_evol_model.state_bound_constant_C1 = 1;
158  lin_evol_model.output_bound_constant_C2 = 1.2916;
159  end;
160 
161  ds_model = lin_ds_from_lin_evol_model(lin_evol_model);
162 % ds_model.theta = 0.0;
163  %ds_model.u_function_ptr = @(params) 0.1; % define new
164  ds_model_data = gen_model_data(ds_model);
165  mu = [0.5,1,1];
166  % set mu in model and base_model!!!
167  ds_model = ds_model.set_mu(ds_model,mu);
168 
169  ds_sim_data = detailed_simulation(ds_model,ds_model_data);
170  lin_evol_model = lin_evol_model.set_mu(lin_evol_model,mu);
171  lin_evol_sim_data = detailed_simulation(lin_evol_model,...
172  lin_evol_model_data);
173  err = lin_evol_sim_data.U - ds_sim_data.X;
174  disp(['max error in DOFs of lin_evol and ds simulation: ',...
175  num2str(max(abs(err(:))))]);
176 
177  params.plot_title='detailed lin_evol';
178  params.no_lines = 1;
179  plot_sim_data(lin_evol_model,lin_evol_model_data,lin_evol_sim_data,params);
180 
181  params.plot_title='detailed lin_ds';
182 % lin_evol_from_ds_sim_data.U = ds_sim_data.X;
183 % lin_evol_from_ds_sim_data.y = ds_sim_data.Y;
184 % plot_sim_data(lin_evol_model,lin_evol_model_data,lin_evol_from_ds_sim_data,params);
185 % lin_ds_from_lin_evol_plot_sim_data(ds_model,ds_model_data,ds_sim_data,params);
186  plot_sim_data(ds_model,ds_model_data,ds_sim_data,params);
187 
188  params.plot_title='error explicit lin_ds and lin_evol';
189  error_sim_data = lin_evol_sim_data;
190  error_sim_data.U = error_sim_data.U-ds_sim_data.X;
191  plot_sim_data(lin_evol_model,lin_evol_model_data,error_sim_data,params);
192  keyboard;
193 
194  %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
195  %step = 4 % reduced Basis generation by single trajectories and
196  %comparison of detailed and reduced simulation
197  %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
198 
199  case 4 % generate some trajectories
200 
201  ds_model = fast_model(params);
202 % lin_evol_model = advection_fv_output_vconst_model(params);
203  %lin_evol_model.debug = 1;
204 % lin_evol_model_data = gen_model_data(lin_evol_model);
205 
206  % set some additional fields in model which will be copied into
207  % ds model depending on whether constants are known or not:
208 % estimate_constants = 0;
209 % if estimate_constants
210 % %%% if constant estimaton has to be performed:
211 % lin_evol_model.estimate_lin_ds_nmu = 4;
212 % lin_evol_model.estimate_lin_ds_nX = 10;
213 % lin_evol_model.estimate_lin_ds_nt = 4;
214 % else
215 % %%% otherwise:
216 % % the following values are rough bounds, so could be
217 % % non-rigorous, then larger search must be performed.
218 % lin_evol_model.state_bound_constant_C1 = 1;
219 % lin_evol_model.output_bound_constant_C2 = 1.2916;
220 % end;
221 %
222 % ds_model = lin_ds_from_lin_evol_model(lin_evol_model);
223  %ds_model.u_function_ptr = @(params) 0.1; % define new
224  ds_model_data = gen_model_data(ds_model);
225 % mu = [1,1,0.5];
226 % mu = [0.5,0.5,0.5];
227 
228  mus = {[0,0,0],[1,1,1],[1,1,0.5],[0,0.5,1]};
229  for m = 1:length(mus)
230  disp(['computing trajectory number ',num2str(m)]);
231  mu = mus{m};
232  ds_model = ds_model.set_mu(ds_model,mu);
233  %ds_model.theta = 0;
234  %ds_model.theta = 1;
235  %disp('set theta to 1!!!');
236  tic;
237  ds_sim_data = detailed_simulation(ds_model,ds_model_data);
238  t_det = toc;
239  disp(['time for full simulation: ',num2str(t_det)]);
240  save(fullfile(rbmatlabresult,[trajectory_fnbase,num2str(m)]),...
241  'ds_model','ds_model_data','ds_sim_data');
242  end;
243 
244  case 4.5 % compute reduced basis.
245  % make sure that the previous trajectories are computed!!
246 
247  disp('PCA of snapshots takes a while...');
248  load(fullfile(rbmatlabresult,[trajectory_fnbase,'1']),...
249  'ds_model','ds_model_data','ds_sim_data');
250  % smaller epsilon does result in more than 1 vector for mu=0,0,0
251  epsilon = 1e-4;
252  RB = PCA_fixspace(ds_sim_data.X,[],ds_model_data.G,[],[],epsilon);
253  %RB = orthonormalize_qr(ds_sim_data.X,ds_model_data.G);
254  % select correct orthogonal set
255  K = RB' * ds_model_data.G * RB;
256  i = find(abs(diag(K)-1)>1e-2);
257  i = min(i);
258  disp(['found ',num2str(size(RB,2)),' basis vectors after traj. 1']);
259  if isempty(i)
260  i = size(RB,2)+1;
261  end;
262 % keyboard;
263  RB = RB(:,1:(i-1));
264  disp(['reduced to ',num2str(size(RB,2)),' basis vectors.']);
265 % keyboard;
266 
267  load(fullfile(rbmatlabresult,[trajectory_fnbase,'2']),...
268  'ds_model','ds_model_data','ds_sim_data');
269  RB2 = PCA_fixspace(ds_sim_data.X,RB,ds_model_data.G,[],[],epsilon);
270  %RB = orthonormalize_qr(ds_sim_data.X,ds_model_data.G);
271  % select correct orthogonal set
272  K = RB2' * ds_model_data.G * RB2;
273  i = find(abs(diag(K)-1)>1e-2);
274  i = min(i);
275  if isempty(i)
276  i = size(RB2,2)+1;
277  end;
278  disp(['found ',num2str(size(RB2,2)),' basis vectors after traj. 2']);
279 % keyboard;
280  RB2 = RB2(:,1:(i-1));
281  disp(['reduced to ',num2str(size(RB2,2)),' basis vectors after traj. 2']);
282  RB = [RB, RB2];
283 % RB = orthonormalize(RB);
284 
285  V = RB;
286  W = ds_model_data.G * V;
287  save(fullfile(rbmatlabresult,rb_fname),'RB');
288 
289  ds_model.RB_basis_filename = rb_fname;
290  ds_model.RB_generation_mode = 'file_load';
291 
292  detailed_data = gen_detailed_data(ds_model,ds_model_data);
293 
294  model = ds_model;
295  save(fullfile(rbmatlabresult,detailed_data_fname),...
296  'model','detailed_data')
297 
298  case 4.75 % compare reduced and detailed simulation
299 
300  load(fullfile(rbmatlabresult,detailed_data_fname),...
301  'model','detailed_data')
302  load(fullfile(rbmatlabresult,[trajectory_fnbase,'2']),...
303  'ds_sim_data');
304 
305  ds_model = model;
306 
307  %perform single reduced simulation for given parameter
308  reduced_data = gen_reduced_data(ds_model,detailed_data);
309  ds_model.N = reduced_data.N;
310  tic,
311 % ds_model.nt = ds_model.nt/4;
312 % disp('set nt to quarter!!');
313  rb_sim_data = rb_simulation(ds_model,reduced_data);
314  t_red = toc;
315  rb_sim_data = rb_reconstruction(ds_model,detailed_data,rb_sim_data)
316  disp(['time for reduced simulation: ',num2str(t_red)]);
317 
318  % plot of reduced trajectory
319  lin_evol_model = ds_model.base_model;
320  lin_evol_model_data = gen_model_data(lin_evol_model);
321 
322  params.no_lines = 1;
323  params.show_colorbar = 1;
324  params.axis_equal = 1;
325  params.plot = lin_evol_model.plot;
326  % plot single snapshot
327  %lin_evol_model.plot(rb_sim_data.X(:,1),lin_evol_model_data.grid, ...
328  % params);
329 
330  % plot sequence:
331  params.title = 'reduced solution';
332  plot_sequence(rb_sim_data.X,lin_evol_model_data.grid,params)
333 
334  % plot sequence:
335  params.title = 'error';
336  plot_sequence((rb_sim_data.X-ds_sim_data.X),...
337  lin_evol_model_data.grid,params)
338 
339  % plot of reduced output
340  figure,plot(rb_sim_data.time,[rb_sim_data.Y;ds_sim_data.Y]);
341  legend('reduced output','detailed output');
342 
343  keyboard;
344 
345 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
346  case 5 % generate plots of detailed data for paper:
347  % step = 5 % figures of detailed simuations based on step 4
348 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
349 
350  load(fullfile(rbmatlabresult,detailed_data_fname),...
351  'model','detailed_data')
352 
353 % model = fast_model(params);
354 % rb_fname = 'advection_fv_output_basis';
355 % load(fullfile(rbmatlabresult,rb_fname),'ds_model','ds_model_data');
356 
357  % make sure, that these trajectories exist.
358  trajectories = [2,3,4];
359 
360  for tr = 1:length(trajectories);
361  trname = [trajectory_fnbase,num2str(trajectories(tr))];
362  % plots of trajectory
363  load(fullfile(rbmatlabresult,trname));
364  lin_evol_model = ds_model.base_model;
365  lin_evol_model_data = gen_model_data(lin_evol_model);
366 
367  params.no_lines = 1;
368  params.show_colorbar = 0;
369  params.axis_equal = 1;
370  params.plot = lin_evol_model.plot;
371  lin_evol_model.plot(ds_sim_data.X(:,1),lin_evol_model_data.grid, ...
372  params);
373  set(gca,'Xtick',[])
374  set(gca,'Ytick',[])
375  saveas(gcf,fullfile(rbmatlabresult,[trname,'_t1.eps']),'epsc');
376  close(gcf);
377 
378  lin_evol_model.plot(ds_sim_data.X(:,64),lin_evol_model_data.grid, ...
379  params);
380  set(gca,'Xtick',[])
381  set(gca,'Ytick',[])
382  saveas(gcf,fullfile(rbmatlabresult,[trname,'_tmiddle.eps']),'epsc');
383  close(gcf);
384 
385  lin_evol_model.plot(ds_sim_data.X(:,end),lin_evol_model_data.grid, ...
386  params);
387  set(gca,'Xtick',[])
388  set(gca,'Ytick',[])
389  saveas(gcf,fullfile(rbmatlabresult,[trname,'_tend.eps']),'epsc');
390  close(gcf);
391 
392  p = ds_model.plot_sim_data(ds_model, ds_model_data, ds_sim_data,params);
393  % Ylim = get(gca,'Ylim');
394  % Ylim = [Ylim(1),Ylim(2)*1.1];
395  set(gca,'Ylim',[0,0.09]);
396  set(p(2),'Color',[0,0.8,0]);
397  saveas(gcf,fullfile(rbmatlabresult,[trname,'_output.eps']),'epsc');
398  close(gcf);
399  close(gcf);
400  end;
401 
402 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
403 % step = 6 % figures of reduced simuations based on step 4
404 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
405  case 6 % generate plots of reduced model for poster/paper
406 
407  disp('generate plots of reduced simulations output estimation');
408 
409  load(fullfile(rbmatlabresult,detailed_data_fname),...
410  'model','detailed_data')
411  ds_model = model;
412 
413  % as reference
414 % load(fullfile(rbmatlabresult,'trajectory1'));
415  load(fullfile(rbmatlabresult,[trajectory_fnbase,'2']));
416 
417 % ds_model.error_estimation = 1;
418 % ds_model.RB_generation_mode = 'file_load'
419 % ds_model.RB_basis_filename = fullfile(...
420 % rbmatlabresult,'advection_fv_output_basis');
421 % detailed_data = gen_detailed_data(ds_model,ds_model_data);
422  reduced_data = gen_reduced_data(ds_model,detailed_data);
423 
424  Ns = [150:-10:10];
425 
426  Ns = Ns(find(reduced_data.N>=Ns));
427 
428  for nind = 1:length(Ns)
429  % perform reduced simulation
430  ds_model.N = Ns(nind);
431  rb_sim_data = rb_simulation(ds_model,reduced_data);
432  rb_sim_data = rb_reconstruction(ds_model,detailed_data,rb_sim_data)
433 
434  params.no_lines = 1;
435  params.show_colorbar = 1;
436  params.axis_equal = 1;
437 
438  % plot of reduced output
439  figure;
440  p = plot(rb_sim_data.time,...
441  [rb_sim_data.Y;...
442  rb_sim_data.Y+rb_sim_data.DeltaY; ...
443  rb_sim_data.Y-rb_sim_data.DeltaY; ...
444  ds_sim_data.Y]);
445  colors = {[1,0,0],[1,0,0],[1,0,0],[0,1,0]};
446  linestyles = {'-','-.','-.','-'};
447  widths = [2,2,2,2];
448  for i = 1:length(p)
449  set(p(i),'Color',colors{i},'Linestyle',linestyles{i},...
450  'Linewidth',widths(i));
451  end;
452  % set marker on rb_sim_data:
453  Xdata = get(p(1),'XData');
454  Ydata = get(p(1),'YData');
455  Zdata = get(p(1),'ZData');
456  set(p(1),'XData',Xdata(1:4:end));
457  set(p(1),'YData',Ydata(1:4:end));
458  set(p(1),'ZData',Zdata(1:4:end));
459  set(p(1),'Marker','o');
460 
461  l = legend(['y\^(t): reduced output'],...
462  'y\^(t)+\Delta_y(t)',...
463  'y\^(t)-\Delta_y(t)',...
464  ['y(t): detailed output'],'Location','NorthWest');
465  set(l,'Fontsize',12);
466  title(['k = ',num2str(ds_model.N)],'Fontsize',20);
467  xlabel('time t','Fontsize',15);
468  ylabel('output','Fontsize',15);
469  % keyboard;
470  saveas(gcf,...
471  fullfile(rbmatlabresult,...
472  ['advection_fv_output_vconst_output_estimation',num2str(ds_model.N),'.eps']),'epsc');
473  close(gcf);
474  end;
475 
476 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
477 % step = 7 % time measurements
478 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
479  case 7
480 
481 % clear all;
482  nreps = 10;
483  t_detailed = zeros(1,nreps);
484 
485  % initialize detailed model and simulate
486 
487  load(fullfile(rbmatlabresult,detailed_data_fname),...
488  'model','detailed_data')
489  ds_model = model;
490 
491 % lin_evol_model = advection_fv_output_vconst_model(params);
492  %lin_evol_model.debug = 1;
493  %lin_evol_model_data = gen_model_data(lin_evol_model);
494 
495  % set some additional fields in model which will be copied into
496  % ds model depending on whether constants are known or not:
497 % estimate_constants = 0;
498 % if estimate_constants
499 % %%% if constant estimaton has to be performed:
500 % lin_evol_model.estimate_lin_ds_nmu = 4;
501 % lin_evol_model.estimate_lin_ds_nX = 10;
502 % lin_evol_model.estimate_lin_ds_nt = 4;
503 % else
504 % %%% otherwise:
505 % % the following values are rough bounds, so could be
506 % % non-rigorous, then larger search must be performed.
507 % lin_evol_model.state_bound_constant_C1 = 1;
508 % lin_evol_model.output_bound_constant_C2 = 1.2916;
509 % end;
510 %
511 % ds_model = lin_ds_from_lin_evol_model(lin_evol_model);
512 % % set accelerated data functions from below
513 % ds_model.A_function_ptr = @(model,model_data) ...
514 % eval_affine_decomp(@A_components,...
515 % @A_coefficients,...
516 % model,model_data);
517 %
518 % ds_model.B_function_ptr = @(model,model_data) ...
519 % eval_affine_decomp(@B_components,...
520 % @B_coefficients,...
521 % model,model_data);
522 % %ds_model.u_function_ptr = @(params) 0.1; % define new
523  ds_model_data = gen_model_data(ds_model);
524 
525  disp('skipping detailed simulations...')
526 % if 0
527 
528  for r = 1:nreps
529  disp(['repetition = ',num2str(r)])
530  mu = [1,1,0.5]-rand(1,3)*0.01;
531  ds_model = ds_model.set_mu(ds_model,mu);
532  tic,
533  ds_sim_data = detailed_simulation(ds_model,ds_model_data);
534  t_detailed(1,r) = toc;
535  end;
536 
537 % end;
538 
539  disp('--------------------------------------------------')
540  disp('detailed simulation')
541  disp(t_detailed);
542  disp(['mean t_detailed=',num2str(mean(t_detailed))]);
543  disp('--------------------------------------------------')
544 
545  % reduced simulation with and without error estimation
546 
547  clear('ds_sim_data');
548  ds_model.RB_generation_mode = 'file_load';
549  ds_model.RB_basis_filename = ...
550  fullfile(rbmatlabresult,rb_fname);
551  detailed_data = gen_detailed_data(ds_model,ds_model_data);
552  reduced_data = gen_reduced_data(ds_model,detailed_data);
553 
554  Ns = [150:-10:10]; %,30,20,10];
555  Ns = Ns(find(reduced_data.N>=Ns));
556 
557  t_reduced_with_err_est = zeros(length(Ns),nreps);
558  t_reduced_without_err_est = zeros(length(Ns),nreps);
559 
560  % reduced simulation with and without error estimation
561 
562  for err_est = 0:1;
563  ds_model.error_estimation = err_est;
564  disp(['err_est= ',num2str(err_est)]);
565  for n_ind = 1:length(Ns)
566  disp(['M = ',num2str(Ns(n_ind))]);
567  ds_model.N = Ns(n_ind);
568  for r = 1:nreps
569  disp(['repetition = ',num2str(r)])
570  mu = [1,1,1]-rand(1,3)*0.01;
571  ds_model = ds_model.set_mu(ds_model,mu);
572  tic;
573  rb_sim_data = rb_simulation(ds_model,reduced_data);
574  if size(rb_sim_data.Xr,1)~=Ns(n_ind)
575  disp('N dimension not set correctly!');
576  keyboard;
577  end;
578  if err_est == 1
579  t_reduced_with_err_est(n_ind,r) = toc;
580  else
581  t_reduced_without_err_est(n_ind,r) = toc;
582  end;
583  end;
584  end;
585  end;
586 
587  % output measurements
588 
589  disp('--------------------------------------------------')
590  disp('detailed simulation')
591  disp(t_detailed);
592  disp(['mean t_detailed=',num2str(mean(t_detailed))]);
593 
594  disp('--------------------------------------------------')
595  disp('reduced simulation with error estimation')
596  disp(t_reduced_with_err_est);
597  disp(['mean t_reduced_with_err_est=',...
598  num2str(mean(t_reduced_with_err_est,2)')]);
599  disp(['Ns =',num2str(Ns)]);
600 
601  disp('--------------------------------------------------')
602  disp('reduced simulation without error estimation')
603  disp(t_reduced_without_err_est);
604  disp(['mean t_reduced_without_err_est=',...
605  num2str(mean(t_reduced_without_err_est,2)')]);
606  disp(['Ns =',num2str(Ns)]);
607  disp('--------------------------------------------------')
608 
609  save(fullfile(rbmatlabresult,'time_measurements'));
610 
611 % dim_x = 65536, nt = 2048
612 %
613 %mean t_detailed=64.3806
614 %--------------------------------------------------
615 %reduced simulation with error estimation
616 % 3.4542 3.5969 3.5659 3.4216 3.6170
617 % 3.3648 3.3816 3.4218 3.4391 3.4552
618 % 3.1643 3.1312 3.1244 3.1433 3.1720
619 % 3.0594 3.0496 3.0628 3.0486 3.0387
620 % 2.9993 2.9953 3.0078 2.9413 2.9603
621 % 2.9555 2.9093 2.9091 2.9510 2.9256%%%
622 %
623 %mean t_reduced_with_err_est=3.5311 3.4125 3.147 3.0518 2.9808 2.9301
624 %Ns =60 50 40 30 20 10
625 %--------------------------------------------------
626 %reduced simulation without error estimation
627 % 2.3975 2.4220 2.2924 2.4490 2.3302
628 % 2.3322 2.3198 2.3425 2.3056 2.3512
629 % 2.2227 2.2734 2.2210 2.2582 2.2504
630 % 2.2531 2.2180 2.2429 2.2212 2.2612
631 % 2.2283 2.1979 2.2581 2.2241 2.1957
632 % 2.2274 2.2014 2.2577 2.1856 2.2268
633 %
634 %mean t_reduced_without_err_est=2.3782 2.3303 2.2451 2.2393 2.2208 2.2198
635 %Ns =60 50 40 30 20 10
636 %--------------------------------------------------
637 
638  case 8 % compare original and accelerated detailed model
639 
640  disp('compare original and accelerated lin_ds models');
641 
642  % clear all;
643  % initialize detailed model and simulate
644  lin_evol_model = advection_fv_output_vconst_model(params);
645  estimate_constants = 0;
646  if estimate_constants
647  %%% if constant estimaton has to be performed:
648  lin_evol_model.estimate_lin_ds_nmu = 4;
649  lin_evol_model.estimate_lin_ds_nX = 10;
650  lin_evol_model.estimate_lin_ds_nt = 4;
651  else
652  %%% otherwise:
653  % the following values are rough bounds, so could be
654  % non-rigorous, then larger search must be performed.
655  lin_evol_model.state_bound_constant_C1 = 1;
656  lin_evol_model.output_bound_constant_C2 = 1.2916;
657  end;
658 
659  ds_model = lin_ds_from_lin_evol_model(lin_evol_model);
660  %ds_model.u_function_ptr = @(params) 0.1; % define new
661  ds_model_data = gen_model_data(ds_model);
662 
663  % set accelerated data functions
664  ds_model_fast = lin_ds_from_lin_evol_model(lin_evol_model);
665  ds_model_fast.A_function_ptr = @(model,model_data) ...
666  eval_affine_decomp(@A_components,...
667  @A_coefficients,...
668  model,model_data);
669 
670  ds_model_fast.B_function_ptr = @(model,model_data) ...
671  eval_affine_decomp(@B_components,...
672  @B_coefficients,...
673  model,model_data);
674 
675  % ds_model_fast.C_function_ptr = @(model,model_data) ...
676  % eval_affine_decomp(@C_components,...
677  % @C_coefficients,...
678  % model,model_data);
679 
680  % disp('skipping detailed simulations...')
681  % if 0
682 
683  mu = [1,1,0.5]-rand(1,3)*0.01;
684  ds_model = ds_model.set_mu(ds_model,mu);
685  ds_model_fast = ds_model_fast.set_mu(ds_model_fast,mu);
686  tic,
687  ds_sim_data = detailed_simulation(ds_model,ds_model_data);
688  t_detailed = toc;
689  tic,
690  ds_sim_data_fast = detailed_simulation(ds_model_fast,ds_model_data);
691  t_detailed_fast = toc;
692 
693  disp(['original t = ',num2str(t_detailed),...
694  ', accelerated t = ',num2str(t_detailed_fast)]);
695  if ~isequal(ds_sim_data,ds_sim_data_fast)
696  disp('simulation results of original and accelerated model differ!')
697  else
698  disp('simulation results of original and accelerated model equal!')
699  end;
700 
701  keyboard;
702 
703  case 9 % generate basis from POD-Greedy
704 
705  filecache_clear;
706  ds_model = fast_model(params);
707  ds_model.cone_range = [0.4,0.6];
708  ds_model.base_model.cone_range = ds_model.cone_range;
709  % ds_model.
710  % shrink range to 'single point'
711 % r = [0.4,0.6]; % => 2 -> 0.8 error nach 10 basisvektoren
712 % r = [0.0,1.0]; % => 4.6 -> 1.9 error nach 10 basisvektoren
713 % r = [0.45,0.55]; % => 2.38 -> 0.51 error nach 10 basisvektoren
714  ds_model.mu_ranges = {[0,1],[0.4,0.6],[0.4,0.6]};
715  % ds_model.rb_basis_generation_ptr = ...
716  % @basisgen_POD_greedy_uniform_fixed;
717 
718  % ds_model.RB_generation_mode = 'uniform_fixed';
719  % ds_model.RB_num_intervals = [2,2,2];
720 
721  ds_model_data = gen_model_data(ds_model);
722 
723  disp('generating basis from POD-Greedy');
724 
725  % for basis generation, set further fields:
726 
727  ds_model.verbose = 5; % => only dots in greedy loop
728  ds_model.RB_generation_mode = 'greedy_uniform_fixed';
729  ds_model.RB_numintervals = [1,1,1]; % => 0.5 in cone_pos obtained
730 % ds_model.RB_numintervals = [4,4,4];
731  ds_model.RB_stop_epsilon = 1e-5;
732  ds_model.RB_stop_timeout = 3600;
733  ds_model.RB_stop_Nmax =150;
734 % ds_model.RB_stop_timeout = 600;
735 % ds_model.RB_stop_Nmax =10;
736  ds_model.RB_extension_algorithm = @RB_extension_PCA_fixspace;
737  ds_model.RB_error_indicator = 'estimator';
738  % default:
739  % ds_model.get_estimator_from_sim_data = @(sim_data) sim_data.DeltaX(end);
740 
741  detailed_data = gen_detailed_data(ds_model,ds_model_data);
742 
743  model = ds_model;
744  dfname = fullfile(rbmatlabresult,[detailed_data_fname,'_greedy']);
745  save(dfname,'model','detailed_data');
746 
747  disp('successfully generated and saved detailed_data');
748 
749  params.plot_title = 'reduced basis';
750  params.plot = model.base_model.plot;
751  lin_evol_model_data = gen_model_data(model.base_model);
752  plot_sequence(detailed_data.V,lin_evol_model_data.grid,params)
753 
754  keyboard;
755 
756  case 10 % reduced simulations and error plots
757 
758  %basis_from_step = 9; % step 9
759  basis_from_step = 4; % step 4
760 
761  if basis_from_step == 9
762  dfname = fullfile(rbmatlabresult,[detailed_data_fname,'_greedy']);
763  load(dfname);
764  disp('loaded basis from step 9')
765  else
766  load(fullfile(rbmatlabresult,detailed_data_fname),...
767  'model','detailed_data')
768 
769 % model = fast_model(params);
770  model_data = gen_model_data(model);
771  % dfname = fullfile(rbmatlabresult,'advection_fv_output_detailed_data_146');
772 % model.RB_basis_filename = ...
773 % fullfile(rbmatlabresult,'advection_fv_output_basis');
774 % model.RB_generation_mode = 'file_load';
775 % detailed_data = gen_detailed_data(model,model_data);
776  % dfname = fullfile(rbmatlabresult,'advection_fv_output_detailed');
777  % load(dfname);
778  tmp = load(fullfile(rbmatlabresult,[trajectory_fnbase,'2']));
779  ds_sim_data = tmp.ds_sim_data;
780  mu_from_file = get_mu(tmp.ds_model);
781  disp('loaded basis from step 4')
782  model = model.set_mu(model,mu_from_file);
783  end;
784  % trajectory2 above is made with this
785  % model = model.set_mu(model,[1,1,1]);
786  % model = model.set_mu(model,[1,1,1]);
787  % trajectory3 above is made with this:
788  %model = model.set_mu(model,[1,1,0.5]);
789  %model = model.set_mu(model,[0.5,0.6,0.6]);
790  %model = model.set_mu(model,[0.5,0.4,0.4]);
791  %model = model.set_mu(model,[0.5,0.5,0.5]);
792  model.N = size(detailed_data.V,2);
793  %model.N = 150;
794  reduced_data = gen_reduced_data(model,detailed_data);
795  rb_sim_data = rb_simulation(model,reduced_data);
796 
797  params.no_lines = 1;
798  params.show_colorbar = 1;
799  params.axis_equal = 1;
800 
801  % plot of reduced output
802  figure;
803  p = plot(rb_sim_data.time,...
804  [rb_sim_data.Y;...
805  rb_sim_data.Y+rb_sim_data.DeltaY; ...
806  rb_sim_data.Y-rb_sim_data.DeltaY]);
807  colors = {[1,0,0],[1,0,0],[1,0,0],[0,1,0]};
808  linestyles = {'-','-.','-.','-'};
809  widths = [2,2,2,2];
810  for i = 1:length(p)
811  set(p(i),'Color',colors{i},'Linestyle',linestyles{i},...
812  'Linewidth',widths(i));
813  end;
814  % set marker on rb_sim_data:
815  Xdata = get(p(1),'XData');
816  Ydata = get(p(1),'YData');
817  Zdata = get(p(1),'ZData');
818  set(p(1),'XData',Xdata(1:4:end));
819  set(p(1),'YData',Ydata(1:4:end));
820  set(p(1),'ZData',Zdata(1:4:end));
821  set(p(1),'Marker','o');
822 
823  legend(['y\^(t): reduced output'],...
824  'y\^(t)+\Delta_y(t)',...
825  'y\^(t)-\Delta_y(t)','Location','NorthWest');
826  title(['k = ',num2str(model.N)],'Fontsize',15);
827 
828  % table of error, error estimator and effectivity
829 
830  reduced_data = gen_reduced_data(model,detailed_data);
831 
832  Ns = 150:-10:10;
833  Ns = Ns(find(reduced_data.N>=Ns));
834 
835  err = zeros(1,length(Ns));
836  err_estim = zeros(1,length(Ns));
837  effectivity = zeros(1,length(Ns));
838 
839  for ni = 1:length(Ns)
840  model.N = Ns(ni);
841  rb_sim_data = rb_simulation(model, reduced_data);
842  err(ni) = abs(ds_sim_data.Y(end) - rb_sim_data.Y(end));
843  err_estim(ni) = rb_sim_data.DeltaY(end);
844  effectivity(ni) = err_estim(ni)/err(ni);
845  end;
846  Ns
847  err
848  err_estim
849  effectivity
850  disp('finished computing effectivities')
851  keyboard;
852 
853  case 11 % test of matlab ode solver for resulting sparse system
854 
855  params.coarse_factor = 8; % nake model smaller
856  model = fast_model(params);
857  model.RB_generation_mode = 'from_detailed_data';
858 
859  % detailed simulation:
860  model_data = gen_model_data(model);
861  model_data.RB = zeros(size(model_data.G,1),0);
862  detailed_data = gen_detailed_data(model,model_data);
863 
864  model.decomp_mode = 0; % complete
865  model = model.set_time(model,0)
866  model = model.set_mu(model,[0.5,1,1]);
867 
868  A = @(model,model_data,t) ...
869  model.A_function_ptr(model.set_time(model,t),model_data);
870 
871  B = @(model,model_data,t) ...
872  model.B_function_ptr(model.set_time(model,t),model_data);
873 
874  x0 = model.x0_function_ptr(model,model_data);
875  tic;
876  [t,x] = ode15s(@(t,x) A(model,model_data,t)*x + B(model,model_data,t), ...
877  [0 1], x0);
878  time_ode15 = toc;
879  disp(['ode15 takes ',num2str(time_ode15),' sec']);
880 
881  params.plot_title = 'ode15 solution';
882  params.plot = model.base_model.plot;
883  params.axis_equal = 1;
884  lin_evol_model_data = gen_model_data(model.base_model);
885  plot_sequence((x'),lin_evol_model_data.grid,params)
886 
887  tic;
888  sim_data = detailed_simulation(model,model_data);
889  time_rbmatlab = toc;
890  disp(['rbmatlab takes ',num2str(time_rbmatlab),' sec']);
891  params.plot_title = 'rbmatlab solution';
892  params.plot = model.base_model.plot;
893  params.axis_equal = 1;
894  lin_evol_model_data = gen_model_data(model.base_model);
895  plot_sequence(sim_data.X,lin_evol_model_data.grid,params)
896 
897  % further plots for MCMCS paper
898  case 12
899 
900  load(fullfile(rbmatlabresult,detailed_data_fname),...
901  'model','detailed_data')
902 
903  % plot of error estimator over parameter variation
904 % model = fast_model(params);
905 % model_data = gen_model_data(model);
906  % dfname = fullfile(rbmatlabresult,'advection_fv_output_detailed_data_146');
907 % model.RB_basis_filename = ...
908 % fullfile(rbmatlabresult,'advection_fv_output_basis');
909 % model.RB_generation_mode = 'file_load';
910 % detailed_data = gen_detailed_data(model,model_data);
911 % save(fullfile(rbmatlabresult,...
912 % 'advection_fv_output_detailed_2trajectories.mat'),...
913 % 'model','detailed_data')
914 % keyboard;
915  % dfname = fullfile(rbmatlabresult,'advection_fv_output_detailed');
916  % load(dfname);
917  tmp = load(fullfile(rbmatlabresult,[trajectory_fnbase,'2']));
918  ds_sim_data = tmp.ds_sim_data;
919  disp('loaded basis from step 4')
920 
921  % factors = 0:0.05:1;
922  reduced_data = gen_reduced_data(model,detailed_data);
923  Ns = [1,8,16,32,64];
924  Ns = Ns(find(reduced_data.N>=Ns));
925  factors = 0:0.05:1;
926  err_estim = zeros(length(factors),length(Ns));
927  for ni = 1:length(Ns)
928  model.N = Ns(ni);
929  for i =1:length(factors)
930  model = model.set_mu(model, factors(i)*[1,1,1]);
931  rb_sim_data = rb_simulation(model, reduced_data);
932  % err(ni) = abs(ds_sim_data.Y(end) - rb_sim_data.Y(end));
933  err_estim(i,ni) = rb_sim_data.DeltaY(end);
934  % effectivity(ni) = err_estim(ni)/err(ni);
935  end;
936  end;
937 % keyboard;
938  p = plot(factors,err_estim');
939  colors = {[1,0,0],[0,0.6,0],[0,0,1],[0,0,0]};
940  linestyles = {'-','-.','-','-.'};
941  for i = 1:4;
942  set(p(i),'Linewidth',2,'Color',colors{i},'Linestyle',linestyles{i});
943  end;
944 
945  xlabel('parameter \mu_1=\mu_2=\mu_3');
946  ylabel('\Delta_y(\mu)');
947  title('error estimator from parameter sweep')
948  legend({'k=1','k=8','k=16','k=32','k=64'},'Location','NorthWest');
949  saveas(gcf,fullfile(rbmatlabresult,'parameter_sweep.eps'),'epsc');
950  keyboard;
951  otherwise
952  disp('step number unknown');
953 end;
954 
955 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
956 % generation of ds_model from lin_evol_model and replace
957 % by fast vectors
958 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
959 
960 function ds_model = fast_model(params);
961  % fast_model
962 lin_evol_model = advection_fv_output_vconst_model(params);
963 %lin_evol_model.debug = 1;
964 %lin_evol_model_data = gen_model_data(lin_evol_model);
965 
966 % set some additional fields in model which will be copied into
967 % ds model depending on whether constants are known or not:
968 estimate_constants = 0;
969 if estimate_constants
970  %%% if constant estimaton has to be performed:
971  lin_evol_model.estimate_lin_ds_nmu = 4;
972  lin_evol_model.estimate_lin_ds_nX = 10;
973  lin_evol_model.estimate_lin_ds_nt = 4;
974 else
975  %%% otherwise:
976  % the following values are rough bounds, so could be
977  % non-rigorous, then larger search must be performed.
978  lin_evol_model.state_bound_constant_C1 = 1;
979  lin_evol_model.output_bound_constant_C2 = 1.2916;
980 end;
981 
982 ds_model = lin_ds_from_lin_evol_model(lin_evol_model);
983 % set accelerated data functions from below
984 ds_model.A_function_ptr = @(model,model_data) ...
985  eval_affine_decomp(@A_components,...
986  @A_coefficients,...
987  model,model_data);
988 
989 ds_model.B_function_ptr = @(model,model_data) ...
990  eval_affine_decomp(@B_components,...
991  @B_coefficients,...
992  model,model_data);
993 %ds_model.u_function_ptr = @(params) 0.1; % define new
994 
995 %ds_model.theta = 1;
996 %disp('set theta to 1!!!');
997 %keyboard;
998 
999 
1000 
1001 
1002 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1003 % auxialiary coefficient functions for acceleration of
1004 % advection model: explicit implementation of coefficients
1005 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1006 
1007 function Acomp = A_components(model,model_data)
1008  % A_components
1009 Acomp = model_data.A_components;
1010 
1011 function Acoeff = A_coefficients(model);
1012  % A_coefficients
1013 %model.base_model.decomp_mode = 2;
1014 %model.base_model.t = model.t;
1015 %mu = get_mu(model);
1016 %model.base_model = model.set_mu(model.base_model,mu);
1017 %[L_I_coeff, L_E_coeff, b_coeff] = ...
1018 % model.base_model.operators_ptr(model.base_model, );
1019 %% L_E = Id + Delta t * A
1020 %Acoeff = L_E_coeff(2:end)/model.base_model.dt;
1021 %keyboard;
1022 %Acoeff = -model.base_model.dt * ...
1023 % model.base_model.velocity_coefficients_ptr(model.base_model);
1024 %if model.t>0
1025  Acoeff = - [model.vx_weight, model.vy_weight]*0.5; %*(1-model.t);
1026 % keyboard;
1027 %end;
1028 
1029 function Bcomp = B_components(model,model_data)
1030  % B_components
1031 Bcomp = model_data.B_components;
1032 
1033 function Bcoeff = B_coefficients(model);
1034  % B_coeffiecients
1035 % hmmm the coefficients are now called twice in A and B
1036 % should somehow be cached later...
1037 %model.base_model.decomp_mode = 2;
1038 %model.base_model.t = model.t;
1039 %mu = get_mu(model);
1040 %model.base_model = model.set_mu(model.base_model,mu);
1041 %[L_I_coeff, L_E_coeff, b_coeff] = ...
1042 % model.base_model.operators_ptr(model.base_model, );
1043 %% L_E = Id + Delta t * A
1044 %Bcoeff2 = b_coeff/model.base_model.dt;
1045 vcoeff = - [model.vx_weight, model.vy_weight]*0.5;%*(1-model.t);
1046 Q_0 = model.base_model.cone_number;
1047 dcoeff = zeros(1,Q_0);
1048 max_pos = model.base_model.cone_weight * (Q_0-1)+1;
1049 t = model.t;
1050 for q = 1:Q_0
1051  dcoeff(q) = (1-(max_pos-q))*(1-t) * ((max_pos>=q) && (max_pos < q+1)) ...
1052  + (1+(max_pos-q))*(1-t) * ((max_pos>=q-1) && (max_pos < q));
1053 end;
1054 v = -(vcoeff(:)*(dcoeff(:)'));
1055 Bcoeff = v(:)*model.cone_amplitude;
1056 %if model.t>0
1057 %keyboard;
1058 %end;
1059 
1060 %function Ccomp = C_components(model,model_data)
1061 %Ccomp = model_data.C_components;
1062 %
1063 %function Ccoeff = C_coefficients(model);
1064 %model.base_model.decomp_mode = 2;
1065 %model.base_model.t = model.t;
1066 %mu = get_mu(model);
1067 %model.base_model = model.set_mu(model.base_model,mu);
1068 %Ccoeff = ...
1069 % model.base_model.operators_output(model.base_model);
1070 
1071