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