rbmatlab  1.16.09
 All Classes Namespaces Files Functions Variables Modules Pages
newton_model.m
1 function model = newton_model(params)
2 % function model = newton_model(params)
3 %
4 % optional fields of params:
5 % model_size: - 'big' simulations are done on a larger grid with very
6 % small time steps
7 % - 'small' simulations are done on a smaller grid with
8 % bigger time steps
9 % model_type: - 'buckley_leverett' -- A Buckley-Leverett problem
10 % - 'linear_diffusion' -- A linear parabolic problem with
11 % constant diffusion
12 % `` \partial_{t} u - \nabla \cdot ( \mu_1 \nabla u ) = 0 ``
13 % - 'exponential_diffusion' -- A non-linear parabolic problem
14 % with exponential function as diffusion coefficient
15 % `` \partial_{t} u
16 % - \nabla \cdot ( k_0 + \mu_1 u^{\mu_2} \nabla u ) ``
17 % - 'eoc' -- A linear diffusion problem with known exact
18 % solution for validation of numerical scheme
19 % - 'eoc_nonlinear' -- A nonlinear diffusion problem with known
20 % exact solution for validation of numerical scheme
21 %
22  model = nonlin_evol_model_default;
23  model.name = 'newton_model';
24  model.rb_problem_type = 'newton';
25  %% data that is definitely used outside of the detailed simulations
26  % model.mu_names = {'diff_m', 'diff_p', 'hill_height'};
27  % model.mu_ranges = {[0,0.5],[0.01,100],[0,0.3]};
28  model.mu_names = {'diff_k0'};
29  model.mu_ranges = {[0.01, 0.5]};
30 
31 
32  %% data that might be of interest outside of the detailed simulations
33  model.T = 0.2;
34  model.nt = 80;
35 
36  %% finite volume settings
37 % model.diffusivity_ptr = @diffusivity_buckley_leverett;
38 % model.diffusivity_derivative_ptr = @diffusivity_buckley_leverett_derivative;
39  model.diffusivity_ptr = @diffusivity_exponential;
40  model.diffusivity_derivative_ptr = @diffusivity_exponential_derivative;
41  model.diff_k0 = 0.00;
42  model.diff_m = 1.0;
43  model.diff_p = 2.0;
44 
45  model.num_conv_flux_ptr = @fv_num_conv_flux_engquist_osher;
46  model.num_diff_flux_ptr = @fv_num_diff_flux_gradient;
47  model.flux_linear = 0;
48  model.name_flux = 'exponential';
49  model.geometry_transformation = 'none';
50 
51  model.newton_solver = 1;
52  model.newton_steps = 60;
53 
54 
55  model.operators_diff_implicit = @fv_operators_diff_implicit_gradient;
56  model.operators_conv_implicit = @fv_operators_zero;
57  model.operators_neumann_implicit = @fv_operators_neumann_explicit;
58  model.fv_expl_conv_weight = 0.0;
59  model.fv_impl_diff_weight = 1.0;
60  model.implicit_gradient_operators_algorithm = ...
62 
63  % lxf_lambda = 1.25;
64  % lxf_lambda = fv_search_max_lxf_lambda([],
65 
66 
67  %% define the grid
68  model = unitcube(model);
69  model.xnumintervals = 100;
70  model.ynumintervals = 100;
71 
72  % pointer to function in rbmatlab/datafunc/dirichlet_values
73 % model.dirichlet_values_ptr = @dirichlet_values_uplow;
74 % model.c_dir = 0.2;
75 % model.dir_box_xrange = [-1 2];
76 % model.dir_box_yrange = 0.5 + [1/320 3/320];
77 % model.c_dir_left = 1;
78 % model.c_dir_right = 1;
79 % model.c_dir_up = 0;
80 % model.c_dir_low = 0.5;
81 % model.dir_middle = 0.7;
82  model.dirichlet_values_ptr = @dirichlet_values_leftright;
83  model.c_dir_left = 0.1;
84  model.c_dir_right = 0.6;
85  model.dir_middle = 0.5;
86 % model.c_dir = 0.1;
87 
88  % pointer to function in rbmatlab/datafunc/neumann_values
89  model.neumann_values_ptr = @neumann_values_homogeneous;
90 % model.neumann_values_ptr = @neumann_values_outflow;
91 % model.neumann_values_ptr = @neumann_values_rightflow;
92  model.c_neu = 0.00;
93 
94  % name of function in rbmatlab/datafunc/init_values/
95 % model.init_values_ptr = @init_values_transformed_blobs;
96  model.init_values_ptr = @init_values_bars;
97 % model.init_values_ptr = @init_values_as_dirichlet;
98  % parameters for data functions
99  model.blob_height = 0.4;
100  model.c_init = 0.1;
101 
102  model.filecache_velocity_matrixfile_extract = 0;
103 
104  model.geometry_spline_type = 'linear';
105  model.hill_height = 0.0;
106  model.geometry_transformation_spline_x = [ 0 0.5 1 ];
107  model.geometry_transformation_spline_y = [ 0 -0.033 0 ];
108 
109  model.k = 0.005;
110  model.gravity = 0;
111 
112 
113  % settings for CRB generation
114  model.CRB_generation_mode = 'param-time-space-grid';
115  model.RB_generation_mode = 'greedy_uniform_fixed';
116  model.Mmax = 250;
117  model.sMmax = {50, 50};
118  model.M = 40;
119  model.Mstrich = 5;
120  model.ei_stop_on_Mmax = 1;
121  % build ei-basis with WHOLE trajectory
122  model.ei_target_error = 'interpol';
123  model.ei_numintervals = [8];
124 
125  model.RB_stop_Nmax = 50;
126  model.RB_stop_epsilon = 1e-9;
127  model.RB_stop_max_val_train_ratio = inf;
128  model.RB_numintervals = [8];
129 
130  % target accuracy epsilon:
131 
132  % decide, whether estimator or true error is error-indicator for greedy
133  model.RB_error_indicator = 'error'; % true error
134  % RB_error_indicator = 'estimator'; % Delta from rb_simulation
135  % RB_error_indicator = 'ei_estimator_test'; % Delta from rb_simulation testet
136  % against true error
137 
138  model.divclean_mode = false;
139  model.flux_quad_degree = 1;
140 
141  model.conv_flux_ptr = @conv_flux_linear;
142 
143  model.velocity_ptr = @velocity_transport;
144  model.transport_x = 0.1;
145  model.transport_y = 0.1;
146 
147  model.conv_flux_derivative_ptr = @(glob, U, params) ...
148  params.velocity_ptr(glob, params);
149 
150 
151  model.data_const_in_time = 0;
152 
153  % set all to dirichlet-boundary by specifying "rectangles", all
154  % boundary within is set to boundary type by bnd_rect_index
155  model.bnd_rect_corner1 = [-0.999999, 0; ...
156  0, -0.9999999 ];
157  model.bnd_rect_corner2 = [ 2, 0.999999; ...
158  0.999999, 2 ];
159  % -1 means dirichlet, -2 means neumann
160  model.bnd_rect_index = [ -2, -2 ];
161 
162  if isfield(params, 'separate_CRBs') && params.separate_CRBs
163  model.separate_CRBs = params.separate_CRBs;
164  savepath_infix = '';
165  else
166  savepath_infix = '_one_CRB';
167  end
168 
169  model.ei_detailed_savepath = [model.name,'_',params.model_type,...
170  savepath_infix,'_ei_data_interpol'];
171  model.ei_operator_savepath = [model.name,'_',params.model_type,...
172  savepath_infix,'_ei_operators_interpol'];
173 
174  model.RB_detailed_train_savepath = model.ei_detailed_savepath;
175 
176  if isfield(params, 'model_size')
177  if strcmp(params.model_size, 'big') == 1
178  model.xnumintervals = 200;
179  model.ynumintervals = 200;
180  end
181  if strcmp(params.model_size, 'small') == 1
182  model.xnumintervals = 50;
183  model.ynumintervals = 50;
184  model.T = 2;
185  model.nt = 20;
186  end
187  end
188 
189  if isfield(params, 'model_type')
190  if strcmp(params.model_type, 'linear_diffusion')
191  model.mu_names = {'diff_k0'};
192  model.mu_ranges = {[0.01, 0.5]};
193  model.ei_numintervals = [8];
194  model.RB_numintervals = [8];
195  elseif strcmp(params.model_type, 'buckley_leverett')
196  model.bl_M = 0.50;
197  model.bl_k = -1.0;
198 % model.xnumintervals = 100;
199 % model.ynumintervals = 100;
200 % model.nt = 100;
201  model.xnumintervals = 50;
202  model.ynumintervals = 50;
203  model.nt = 70;
204  model.T = 0.6;
205  model.conv_flux_ptr = @conv_flux_buckley_leverett;
206  model.conv_flux_derivative_ptr = @conv_flux_buckley_leverett_derivative;
207 % model.num_diff_flux_ptr = @fv_num_conv_flux_lax_friedrichs;
208  model.num_conv_flux_ptr = @fv_num_conv_flux_engquist_osher;
209  model.num_diff_flux_ptr = @fv_num_diff_flux_gradient;
210  model.operators_diff_implicit = @fv_operators_diff_implicit_gradient;
211  model.implicit_gradient_operators_algorithm = ...
212  @fv_operators_conv_explicit_engquist_osher;
214  model.operators_conv_implicit = ...
215  @fv_operators_conv_explicit_engquist_osher;
216 % model.operators_conv_explicit = ...
217 % @fv_operators_conv_explicit_engquist_osher;
219  model.fv_impl_conv_weight = 1.0;
220 % model.fv_expl_conv_weight = 1.0;
221  model.fv_impl_diff_weight = 1.0;
222  model.neumann_values_ptr = @neumann_values_rightflow;
223  model.lxf_lambda = 0.3;
224  model.init_values_ptr = @init_values;
225  model.blob_height = 0.5;
226  model.init_values_ptr = @init_values_as_dirichlet;
227  model.dirichlet_values_ptr = @init_values;
228  model.diffusivity_ptr = @diffusivity_homogeneous;
229  model.k = 0.0000;
230 % model.c_dir = -0.1;
231  model.c_init = 0.1;
232  model.c_dir_left = 0.2;
233  model.c_dir_right = 0.2;
234  model.yrange = [0.2 1];
235  model.ei_time_splits = 5;
236  model.Mmax = 300;
237  model.RB_stop_timeout = inf;
238 % model.dirichlet_values_ptr = @dirichlet_values_quarter_of_five;
239 % model.c_dir_right = 0.0;
240 % model.c_dir_left = 0.6;
241 % model.bnd_rect_corner1 = [-0.1, -0.1; ...
242 % 0.9, 0.9; ...
243 % -0.1, 0.099; ...
244 % 0.099, -0.1 ]';
245 % model.bnd_rect_corner2 = [ 0.1, 0.1; ...
246 % 1.1, 1.1; ...
247 % 0.901, 1.1; ...
248 % 1.1, 0.901 ]';
249 % model.bnd_rect_index = [-1, -1, -2, -2];
250  model.conv_a = 1.0;
251  model.bnd_rect_index = [-1, -1];
252 % model.mu_names = {'k','blob_height','c_init','conv_a','bl_M'};
253 % model.mu_ranges = { [0 1e-3], [0.1 0.5], [-0.1 0.1], [0 1], [0.2, 0.5] };
254  model.mu_names = {'blob_height','conv_a','bl_k'};
255  model.mu_ranges = { [0.2 0.5], [0.2 0.5], [-1, -0.5]};
256  model.ei_numintervals = [6, 6, 6];
257  model.RB_numintervals = [6, 6, 6];
258  % purely implicit scheme: only interpolate the implicit operator therefore!
259  model.ei_space_operators = { model.L_I_local_ptr };
260  elseif strcmp(params.model_type, 'exponential_diffusion')
261  model.ei_numintervals = [6,10];
262  model.RB_numintervals = [6,10];
263  model.diff_k0 = 0.0000;
264  model.diff_m = 0.01;
265  model.c_init = 0.0001;
266  model.newton_regularisation = 0;
267  model.diff_p = 2.0;
268  model.T = 1;
269  model.mu_names = {'diff_m','diff_p','diff_k0'};
270  model.mu_ranges = {[0.0,0.01],[1,5],[0,0.0001]};
271  model.ei_numintervals = [4,6,3];
272  model.RB_numintervals = [4,6,3];
273  model.nt = 50;
274  model.RB_stop_Nmax = 50;
275  model.RB_stop_epsilon = 1e-6;
276  model.RB_stop_timeout = inf;
277  model.bnd_rect_index = [-2,-2];
278  model.separate_CRBs = false;
279  model.ei_space_operators = { model.L_I_local_ptr };
280 
281  if isfield(params, 'use_laplacian') && params.use_laplacian
282  model.laplacian_derivative_ptr = @(glob, U, model) ...
283  real((model.diff_m*U.^(model.diff_p-1)));
284  model.laplacian_ptr = @(glob, U, model) ...
285  real((1/model.diff_p * model.diff_m * U.^(model.diff_p)));
286  model.diffusivity_ptr = @(glob, U, model) ...
287  struct('K',1,'epsilon',0);
288  model.diffusivity_derivative = @(glob, U, model) ...
289  struct('K',1,'epsilon',0);
290  model.implicit_gradient_operators_algorithm = ...
292  end
293  elseif strcmp(params.model_type, 'eoc')
294  model.xnumintervals = 10;
295  model.ynumintervals = 10;
296  model.diff_k0 = 0.05;
297  model.diff_m = 0;
298  model.nt = 50;
299  model.init_values_ptr = @exact_function_heat_equation;
300  model.dirichlet_values_ptr = @exact_function_heat_equation;
301  model.bnd_rect_index = [ -1, -1 ];
302  elseif strcmp(params.model_type, 'eoc_nonlinear')
303  model.xnumintervals = 10;
304  model.ynumintervals = 10;
305  model.newton_regularisation = 0;
306  model.diff_k0 = 0.0;
307  model.diff_m = 1;
308  model.diff_p = 1.0;
309  model.nt = 150;
310  model.T = 1.0;
311  model.init_values_ptr = @exact_function_plaplace;
312  model.dirichlet_values_ptr = @exact_function_plaplace;
313  model.data_const_in_time = false;
314  model.bnd_rect_index = [ -1, -2 ];
315  elseif strcmp(params.model_type, 'buckley_leverett')
316  error('buckley leverett problem type is not implemented yet.');
317  end
318  end
319 
320  model.model_type = params.model_type;
321 
322  model.implicit_nonlinear = true;
323 
324  model.data_const_in_time = false;
325 
326  model = model_default(model);
327 
328  if ~isfield(model, 'ei_space_operators')
329  model.ei_space_operators = { model.L_E_local_ptr, model.L_I_local_ptr };
330  end
331 
332  if isfield(params, 'verbose')
333  model.verbose = params.verbose;
334  else
335  model.verbose = 0;
336  end
337  model.debug = 0;
338 
339 end
340 
341 function U0=init_values(glob, params)
342  if ~isfield(params, 'blob_radius')
343  params.blob_radius = 0.2;
344  end
345  params.blob_height = params.blob_height - params.c_init;
346  if ~isfield(params, 'blob_height2')
347  params.blob_height2 = params.blob_height;
348  end
349  params.blob_height2 = params.blob_height2 - params.blob_height;
350 
351  decomp_mode = params.decomp_mode;
352  if decomp_mode == 2
353  U0 = [params.c_init params.blob_height params.blob_height2];
354  else
355  X = glob(:,1);
356  Y = glob(:,2);
357  B1 = sqrt((X(:)-0.75).^2+(Y(:)-0.75).^2);
358  if decomp_mode == 0
359  U0 = params.c_init * ones(size(X));
360  U0 = U0 + params.blob_height * (B1 < params.blob_radius);
361  U0 = U0 + params.blob_height2 * (B1 < 1/2*params.blob_radius);
362  elseif decomp_mode == 1
363  U0 = cell(3,1);
364  U0{1} = ones(size(X(:)));
365  U0{2} = (B1 < params.blob_radius);
366  U0{3} = (B1 < 1/2*params.blob_radius);
367  end
368  end
369 end
370 
371 % vim: set et sw=2:
372 %| \docupdate
function [ L_I_diff , bdir_I_diff ] = fv_operators_diff_implicit_gradient(model, model_data, U, NU_ind)
computes diffusion contributions to finite volume time evolution matrices, or their Frechet derivati...
function Udirichlet = dirichlet_values(model, X, Y)
UDIRICHLET = DIRICHLET_VALUES([X],[Y], MODEL) Examples dirichlet_values([0,1,2],[1,1,1],struct(name_dirichlet_values, homogeneous, ... c_dir, 1)) dirichlet_values([0:0.1:1],[0],struct(name_dirichlet_values, xstripes, ... c_dir, [0 1 2], ... dir_borders, [0.3 0.6])) dirichlet_values([0:0.1:1],[0],struct(name_dirichlet_values, box, ... c_dir, 1, ... dir_box_xrange, [0.3 0.6], ... dir_box_yrange, [-0.1 0.1]))
function num_flux_mat = fv_num_diff_flux_gradient(model, model_data, U, NU_ind)
computes a numerical diffusive flux for a diffusion problem
function [ L_E_neu , b_E_neu ] = fv_operators_neumann_explicit(model, model_data, U, NU_ind)
computes a neumann contribution matrix for finite volume time evolution operators, or their Frechet derivative
function model = model_default(model, T, nt)
model = model_default(model)
Definition: model_default.m:17
function model = unitcube(model)
function adding fields to model for generating a 2D rectgrid with 100 x 100 elements on the unit-squa...
Definition: unitcube.m:17
function [ L_I_diff_jac , bdir_I_diff_jac ] = fv_frechet_operators_diff_implicit_gradient(model, model_data, U, NU_ind)
computes a jacobian of implicit non-linear diffusion contributions to time evolution matrices at a po...
function num_flux = fv_num_conv_flux_lax_friedrichs(model, model_data, U, NU_ind)
Function computing a numerical convective Lax-Friedrichs flux matrix.
function buckley_leverett()
small script demonstrating a buckley leverett problem discretization and RB model reduction ...
function [ L_E_conv , bdir_E_conv ] = fv_operators_conv_explicit_lax_friedrichs(model, model_data, U, NU_ind)
computes convection contribution to finite volume time evolution matrices, or their Frechet derivati...
function num_flux = fv_num_conv_flux_engquist_osher(model, model_data, U, NU_ind)
Function computing a numerical convective Engquist-Osher flux matrix.
function [ flux , lambda ] = conv_flux_linear(glob, U, params)
function computing the convective flux of a convection problem.