source: issm/trunk-jpl/src/m/classes/model.m@ 21827

Last change on this file since 21827 was 21827, checked in by schlegel, 8 years ago

CHG: add a python test for 443, regional output

File size: 66.6 KB
Line 
1%MODEL class definition
2%
3% Usage:
4% md = model(varargin)
5
6classdef model
7 properties (SetAccess=public) %Model fields
8 % {{{
9 %Careful here: no other class should be used as default value this is a bug of matlab
10 mesh = 0;
11 mask = 0;
12
13 geometry = 0;
14 constants = 0;
15 smb = 0;
16 basalforcings = 0;
17 materials = 0;
18 damage = 0;
19 friction = 0;
20 flowequation = 0;
21 timestepping = 0;
22 initialization = 0;
23 rifts = 0;
24 slr = 0;
25
26 debug = 0;
27 verbose = 0;
28 settings = 0;
29 toolkits = 0;
30 cluster = 0;
31
32 balancethickness = 0;
33 stressbalance = 0;
34 groundingline = 0;
35 hydrology = 0;
36 masstransport = 0;
37 thermal = 0;
38 steadystate = 0;
39 transient = 0;
40 levelset = 0;
41 calving = 0;
42 gia = 0;
43 esa = 0;
44
45 autodiff = 0;
46 inversion = 0;
47 qmu = 0;
48 amr = 0;
49
50 results = 0;
51 outputdefinition = 0;
52 radaroverlay = 0;
53 miscellaneous = 0;
54 private = 0;
55
56 %}}}
57 end
58 methods (Static)
59 function md = loadobj(md) % {{{
60 % This function is directly called by matlab when a model object is
61 % loaded. If the input is a struct it is an old version of model and
62 % old fields must be recovered (make sure they are in the deprecated
63 % model properties)
64
65 if verLessThan('matlab','7.9'),
66 disp('Warning: your matlab version is old and there is a risk that load does not work correctly');
67 disp(' if the model is not loaded correctly, rename temporarily loadobj so that matlab does not use it');
68
69 % This is a Matlab bug: all the fields of md have their default value
70 % Example of error message:
71 % Warning: Error loading an object of class 'model':
72 % Undefined function or method 'exist' for input arguments of type 'cell'
73 %
74 % This has been fixed in MATLAB 7.9 (R2009b) and later versions
75 end
76
77 if isstruct(md)
78 disp('Recovering model object from a previous version');
79 md = structtomodel(model,md);
80 end
81
82 %2012 August 4th
83 if isa(md.materials,'materials'),
84 disp('Recovering old materials');
85 if numel(md.materials.rheology_Z)==1 & isnan(md.materials.rheology_Z),
86 md.materials=matice(md.materials);
87 else
88 md.materials=matdamageice(md.materials);
89 end
90 end
91 %2013 April 12
92 if numel(md.stressbalance.loadingforce==1)
93 md.stressbalance.loadingforce=0*ones(md.mesh.numberofvertices,3);
94 end
95 %2013 April 17
96 if isa(md.hydrology,'hydrology'),
97 disp('Recovering old hydrology class');
98 md.hydrology=hydrologyshreve(md.materials);
99 end
100 %2013 October 9
101 if ~isa(md.damage,'damage'),
102 md.damage=damage();
103 md.damage.D=zeros(md.mesh.numberofvertices,1);
104 md.damage.spcdamage=NaN*ones(md.mesh.numberofvertices,1);
105 end
106 %2013 November 18
107 if ~isa(md.outputdefinition,'outputdefinition'),
108 md.outputdefinition=outputdefinition();
109 end
110 %2014 March 26th
111 if isa(md.mesh,'mesh'),
112 disp('Recovering old mesh class');
113 if isprop(md.mesh,'dimension'),
114 if md.mesh.dimension==2,
115 md.mesh=mesh2d(md.mesh);
116 else
117 md.mesh=mesh3dprisms(md.mesh);
118 end
119 else
120 md.mesh=mesh2dvertical(md.mesh);
121 end
122 end
123 %2014 November 12
124 if isa(md.calving,'double'); md.calving=calving(); end
125 %2016 February 3
126 if isa(md.slr,'double'); md.slr=slr(); end
127 %2016 October 11
128 if isa(md.esa,'double'); md.esa=esa(); end
129 %2017 February 10th
130 if md.settings.solver_residue_threshold==0,
131 md.settings.solver_residue_threshold = 1e-6;
132 end
133 %2017 April 10th
134 if isa(md.gia,'gia'), md.gia=giaivins(); end
135 %2017 May 4th
136 if isa(md.amr,'double'); md.amr=amr(); end
137
138 end% }}}
139 end
140 methods
141 function md = model(varargin) % {{{
142
143 switch nargin
144 case 0
145 md=setdefaultparameters(md);
146 case 1
147 error('model constructor not supported yet');
148
149 otherwise
150 error('model constructor error message: 0 of 1 argument only in input.');
151 end
152 end
153 %}}}
154 function md = checkmessage(md,string) % {{{
155 if(nargout~=1) error('wrong usage, model must be an output'); end
156 disp(['model not consistent: ' string]);
157 md.private.isconsistent=false;
158 end
159 %}}}
160 function md = collapse(md)% {{{
161 %COLLAPSE - collapses a 3d mesh into a 2d mesh
162 %
163 % This routine collapses a 3d model into a 2d model
164 % and collapses all the fileds of the 3d model by
165 % taking their depth-averaged values
166 %
167 % Usage:
168 % md=collapse(md)
169 %
170 % See also: EXTRUDE, MODELEXTRACT
171
172 %Check that the model is really a 3d model
173 if ~strcmp(md.mesh.elementtype(),'Penta'),
174 error('collapse error message: only 3d mesh can be collapsed')
175 end
176
177 %Start with changing all the fields from the 3d mesh
178
179 %dealing with the friction law
180 %drag is limited to nodes that are on the bedrock.
181 if isa(md.friction,'friction'),
182 md.friction.coefficient=project2d(md,md.friction.coefficient,1);
183 md.friction.p=project2d(md,md.friction.p,1);
184 md.friction.q=project2d(md,md.friction.q,1);
185 elseif isa(md.friction,'frictioncoulomb'),
186 md.friction.coefficient=project2d(md,md.friction.coefficient,1);
187 md.friction.coefficientcoulomb=project2d(md,md.friction.coefficientcoulomb,1);
188 md.friction.p=project2d(md,md.friction.p,1);
189 md.friction.q=project2d(md,md.friction.q,1);
190 elseif isa(md.friction,'frictionhydro'),
191 md.friction.q=project2d(md,md.friction.q,1);
192 md.friction.C=project2d(md,md.friction.C,1);
193 md.friction.As=project2d(md,md.friction.As,1);
194 md.friction.effective_pressure=project2d(md,md.friction.effective_pressure,1);
195 elseif isa(md.friction,'frictionwaterlayer'),
196 md.friction.coefficient=project2d(md,md.friction.coefficient,1);
197 md.friction.p=project2d(md,md.friction.p,1);
198 md.friction.q=project2d(md,md.friction.q,1);
199 md.friction.water_layer=project2d(md,md.friction.water_layer,1);
200 elseif isa(md.friction,'frictionweertman'),
201 md.friction.C=project2d(md,md.friction.C,1);
202 md.friction.m=project2d(md,md.friction.m,1);
203 elseif isa(md.friction,'frictionweertmantemp'),
204 md.friction.C=project2d(md,md.friction.C,1);
205 md.friction.m=project2d(md,md.friction.m,1);
206 else
207 disp('friction type not supported');
208 end
209
210 %observations
211 if ~isnan(md.inversion.vx_obs), md.inversion.vx_obs=project2d(md,md.inversion.vx_obs,md.mesh.numberoflayers); end;
212 if ~isnan(md.inversion.vy_obs), md.inversion.vy_obs=project2d(md,md.inversion.vy_obs,md.mesh.numberoflayers); end;
213 if ~isnan(md.inversion.vel_obs), md.inversion.vel_obs=project2d(md,md.inversion.vel_obs,md.mesh.numberoflayers); end;
214 if ~isnan(md.inversion.cost_functions_coefficients), md.inversion.cost_functions_coefficients=project2d(md,md.inversion.cost_functions_coefficients,md.mesh.numberoflayers); end;
215 if numel(md.inversion.min_parameters)>1, md.inversion.min_parameters=project2d(md,md.inversion.min_parameters,md.mesh.numberoflayers); end;
216 if numel(md.inversion.max_parameters)>1, md.inversion.max_parameters=project2d(md,md.inversion.max_parameters,md.mesh.numberoflayers); end;
217 if isa(md.smb,'SMBforcing') & ~isnan(md.smb.mass_balance),
218 md.smb.mass_balance=project2d(md,md.smb.mass_balance,md.mesh.numberoflayers);
219 elseif isa(md.smb,'SMBhenning') & ~isnan(md.smb.smbref),
220 md.smb.smbref=project2d(md,md.smb.smbref,md.mesh.numberoflayers);
221 end;
222
223 %results
224 if ~isnan(md.initialization.vx),md.initialization.vx=DepthAverage(md,md.initialization.vx);end;
225 if ~isnan(md.initialization.vy),md.initialization.vy=DepthAverage(md,md.initialization.vy);end;
226 if ~isnan(md.initialization.vz),md.initialization.vz=DepthAverage(md,md.initialization.vz);end;
227 if ~isnan(md.initialization.vel),md.initialization.vel=DepthAverage(md,md.initialization.vel);end;
228 if ~isnan(md.initialization.temperature),md.initialization.temperature=DepthAverage(md,md.initialization.temperature);end;
229 if ~isnan(md.initialization.pressure),md.initialization.pressure=project2d(md,md.initialization.pressure,1);end;
230 if ~isnan(md.initialization.sediment_head),md.initialization.sediment_head=project2d(md,md.initialization.sediment_head,1);end;
231 if ~isnan(md.initialization.epl_head),md.initialization.epl_head=project2d(md,md.initialization.epl_head,1);end;
232 if ~isnan(md.initialization.epl_thickness),md.initialization.epl_thickness=project2d(md,md.initialization.epl_thickness,1);end;
233 if ~isnan(md.initialization.waterfraction),md.initialization.waterfraction=project2d(md,md.initialization.waterfraction,1);end;
234 if ~isnan(md.initialization.watercolumn),md.initialization.watercolumn=project2d(md,md.initialization.watercolumn,1);end;
235 %giaivins
236 if ~isnan(md.gia.mantle_viscosity), md.gia.mantle_viscosity=project2d(md,md.gia.mantle_viscosity,1); end
237 if ~isnan(md.gia.lithosphere_thickness), md.gia.lithosphere_thickness=project2d(md,md.gia.lithosphere_thickness,1); end
238
239 %elementstype
240 if ~isnan(md.flowequation.element_equation)
241 md.flowequation.element_equation=project2d(md,md.flowequation.element_equation,1);
242 md.flowequation.vertex_equation=project2d(md,md.flowequation.vertex_equation,1);
243 md.flowequation.borderSSA=project2d(md,md.flowequation.borderSSA,1);
244 md.flowequation.borderHO=project2d(md,md.flowequation.borderHO,1);
245 md.flowequation.borderFS=project2d(md,md.flowequation.borderFS,1);
246 end
247
248 %boundary conditions
249 md.stressbalance.spcvx=project2d(md,md.stressbalance.spcvx,md.mesh.numberoflayers);
250 md.stressbalance.spcvy=project2d(md,md.stressbalance.spcvy,md.mesh.numberoflayers);
251 md.stressbalance.spcvz=project2d(md,md.stressbalance.spcvz,md.mesh.numberoflayers);
252 md.stressbalance.referential=project2d(md,md.stressbalance.referential,md.mesh.numberoflayers);
253 md.stressbalance.loadingforce=project2d(md,md.stressbalance.loadingforce,md.mesh.numberoflayers);
254 md.masstransport.spcthickness=project2d(md,md.masstransport.spcthickness,md.mesh.numberoflayers);
255 if numel(md.damage.spcdamage)>1, md.damage.spcdamage=project2d(md,md.damage.spcdamage,md.mesh.numberoflayers); end
256 if numel(md.levelset.spclevelset)>1, md.levelset.spclevelset=project2d(md,md.levelset.spclevelset,md.mesh.numberoflayers); end
257 md.thermal.spctemperature=project2d(md,md.thermal.spctemperature,md.mesh.numberoflayers);
258
259 % Hydrologydc variables
260 if isa(md.hydrology,'hydrologydc');
261 md.hydrology.spcsediment_head=project2d(md,md.hydrology.spcsediment_head,1);
262 md.hydrology.mask_eplactive_node=project2d(md,md.hydrology.mask_eplactive_node,1);
263 md.hydrology.sediment_transmitivity=project2d(md,md.hydrology.sediment_transmitivity,1);
264 md.hydrology.basal_moulin_input=project2d(md,md.hydrology.basal_moulin_input,1);
265 if(md.hydrology.isefficientlayer==1)
266 md.hydrology.spcepl_head=project2d(md,md.hydrology.spcepl_head,1);
267 end
268 end
269
270 %materials
271 md.materials.rheology_B=DepthAverage(md,md.materials.rheology_B);
272 md.materials.rheology_n=project2d(md,md.materials.rheology_n,1);
273
274 %damage:
275 if md.damage.isdamage,
276 md.damage.D=DepthAverage(md,md.damage.D);
277 end
278
279 %special for thermal modeling:
280 if ~isnan(md.basalforcings.groundedice_melting_rate),
281 md.basalforcings.groundedice_melting_rate=project2d(md,md.basalforcings.groundedice_melting_rate,1);
282 end
283 if isprop(md.basalforcings,'floatingice_melting_rate') & ~isnan(md.basalforcings.floatingice_melting_rate),
284 md.basalforcings.floatingice_melting_rate=project2d(md,md.basalforcings.floatingice_melting_rate,1);
285 end
286 md.basalforcings.geothermalflux=project2d(md,md.basalforcings.geothermalflux,1); %bedrock only gets geothermal flux
287
288 if isprop(md.calving,'coeff') & ~isnan(md.calving.coeff),
289 md.calving.coeff=project2d(md,md.calving.coeff,1);
290 end
291 if isprop(md.calving,'meltingrate') & ~isnan(md.calving.meltingrate),
292 md.calving.meltingrate=project2d(md,md.calving.meltingrate,1);
293 end
294
295 %update of connectivity matrix
296 md.mesh.average_vertex_connectivity=25;
297
298 %Collapse the mesh
299 nodes2d=md.mesh.numberofvertices2d;
300 elements2d=md.mesh.numberofelements2d;
301
302 %parameters
303 md.geometry.surface=project2d(md,md.geometry.surface,1);
304 md.geometry.thickness=project2d(md,md.geometry.thickness,1);
305 md.geometry.base=project2d(md,md.geometry.base,1);
306 if ~isnan(md.geometry.bed),
307 md.geometry.bed=project2d(md,md.geometry.bed,1);
308 end
309
310 if ~isnan(md.mask.groundedice_levelset),
311 md.mask.groundedice_levelset=project2d(md,md.mask.groundedice_levelset,1);
312 end
313 if ~isnan(md.mask.ice_levelset),
314 md.mask.ice_levelset=project2d(md,md.mask.ice_levelset,1);
315 end
316
317 %outputdefinitions
318 for i=1:length(md.outputdefinition.definitions)
319 if isobject(md.outputdefinition.definitions{i})
320 %get subfields
321 solutionsubfields=fields(md.outputdefinition.definitions{i});
322 for j=1:length(solutionsubfields),
323 field=md.outputdefinition.definitions{i}.(solutionsubfields{j});
324 if length(field)==md.mesh.numberofvertices | length(field)==md.mesh.numberofelements,
325 md.outputdefinition.definitions{i}.(solutionsubfields{j})=project2d(md,md.outputdefinition.definitions{i}.(solutionsubfields{j}),1);
326 end
327 end
328 end
329 end
330
331 %Initialize with the 2d mesh
332 mesh=mesh2d();
333 mesh.x=md.mesh.x2d;
334 mesh.y=md.mesh.y2d;
335 mesh.numberofvertices=md.mesh.numberofvertices2d;
336 mesh.numberofelements=md.mesh.numberofelements2d;
337 mesh.elements=md.mesh.elements2d;
338 if numel(md.mesh.lat) ==md.mesh.numberofvertices, mesh.lat=project2d(md,md.mesh.lat,1); end
339 if numel(md.mesh.long)==md.mesh.numberofvertices, mesh.long=project2d(md,md.mesh.long,1); end
340 mesh.epsg=md.mesh.epsg;
341 if ~isnan(md.mesh.vertexonboundary), mesh.vertexonboundary=project2d(md,md.mesh.vertexonboundary,1); end
342 if ~isnan(md.mesh.elementconnectivity), mesh.elementconnectivity=project2d(md,md.mesh.elementconnectivity,1); end
343 md.mesh=mesh;
344 md.mesh.vertexconnectivity=NodeConnectivity(md.mesh.elements,md.mesh.numberofvertices);
345 md.mesh.elementconnectivity=ElementConnectivity(md.mesh.elements,md.mesh.vertexconnectivity);
346 md.mesh.segments=contourenvelope(md.mesh);
347
348 end % }}}
349 function md2 = extract(md,area) % {{{
350 %extract - extract a model according to an Argus contour or flag list
351 %
352 % This routine extracts a submodel from a bigger model with respect to a given contour
353 % md must be followed by the corresponding exp file or flags list
354 % It can either be a domain file (argus type, .exp extension), or an array of element flags.
355 % If user wants every element outside the domain to be
356 % extract2d, add '~' to the name of the domain file (ex: '~HO.exp');
357 % an empty string '' will be considered as an empty domain
358 % a string 'all' will be considered as the entire domain
359 %
360 % Usage:
361 % md2=extract(md,area);
362 %
363 % Examples:
364 % md2=extract(md,'Domain.exp');
365 %
366 % See also: EXTRUDE, COLLAPSE
367
368 %copy model
369 md1=md;
370
371 %some checks
372 if ((nargin~=2) | (nargout~=1)),
373 help extract
374 error('extract error message: bad usage');
375 end
376
377 %get elements that are inside area
378 flag_elem=FlagElements(md1,area);
379 if ~any(flag_elem),
380 error('extracted model is empty');
381 end
382
383 %kick out all elements with 3 dirichlets
384 spc_elem=find(~flag_elem);
385 spc_node=sort(unique(md1.mesh.elements(spc_elem,:)));
386 flag=ones(md1.mesh.numberofvertices,1);
387 flag(spc_node)=0;
388 pos=find(sum(flag(md1.mesh.elements),2)==0);
389 flag_elem(pos)=0;
390
391 %extracted elements and nodes lists
392 pos_elem=find(flag_elem);
393 pos_node=sort(unique(md1.mesh.elements(pos_elem,:)));
394
395 %keep track of some fields
396 numberofvertices1=md1.mesh.numberofvertices;
397 numberofelements1=md1.mesh.numberofelements;
398 numberofvertices2=length(pos_node);
399 numberofelements2=length(pos_elem);
400 flag_node=zeros(numberofvertices1,1);
401 flag_node(pos_node)=1;
402
403 %Create Pelem and Pnode (transform old nodes in new nodes and same thing for the elements)
404 Pelem=zeros(numberofelements1,1);
405 Pelem(pos_elem)=[1:numberofelements2]';
406 Pnode=zeros(numberofvertices1,1);
407 Pnode(pos_node)=[1:numberofvertices2]';
408
409 %renumber the elements (some nodes won't exist anymore)
410 elements_1=md1.mesh.elements;
411 elements_2=elements_1(pos_elem,:);
412 elements_2(:,1)=Pnode(elements_2(:,1));
413 elements_2(:,2)=Pnode(elements_2(:,2));
414 elements_2(:,3)=Pnode(elements_2(:,3));
415 if isa(md1.mesh,'mesh3dprisms'),
416 elements_2(:,4)=Pnode(elements_2(:,4));
417 elements_2(:,5)=Pnode(elements_2(:,5));
418 elements_2(:,6)=Pnode(elements_2(:,6));
419 end
420
421 %OK, now create the new model!
422
423 %take every field from model
424 md2=md1;
425
426 %automatically modify fields
427
428 %loop over model fields
429 model_fields=fields(md1);
430 for i=1:length(model_fields),
431 %get field
432 field=md1.(model_fields{i});
433 fieldsize=size(field);
434 if isobject(field), %recursive call
435 object_fields=fields(md1.(model_fields{i}));
436 for j=1:length(object_fields),
437 %get field
438 field=md1.(model_fields{i}).(object_fields{j});
439 fieldsize=size(field);
440 %size = number of nodes * n
441 if fieldsize(1)==numberofvertices1
442 md2.(model_fields{i}).(object_fields{j})=field(pos_node,:);
443 elseif (fieldsize(1)==numberofvertices1+1)
444 md2.(model_fields{i}).(object_fields{j})=[field(pos_node,:); field(end,:)];
445 %size = number of elements * n
446 elseif fieldsize(1)==numberofelements1
447 md2.(model_fields{i}).(object_fields{j})=field(pos_elem,:);
448 elseif (fieldsize(1)==numberofelements1+1)
449 md2.(model_fields{i}).(object_fields{j})=[field(pos_elem,:); field(end,:)];
450 end
451 end
452 else
453 %size = number of nodes * n
454 if fieldsize(1)==numberofvertices1
455 md2.(model_fields{i})=field(pos_node,:);
456 elseif (fieldsize(1)==numberofvertices1+1)
457 md2.(model_fields{i})=[field(pos_node,:); field(end,:)];
458 %size = number of elements * n
459 elseif fieldsize(1)==numberofelements1
460 md2.(model_fields{i})=field(pos_elem,:);
461 elseif (fieldsize(1)==numberofelements1+1)
462 md2.(model_fields{i})=[field(pos_elem,:); field(end,:)];
463 end
464 end
465 end
466
467 %modify some specific fields
468
469 %Mesh
470 md2.mesh.numberofelements=numberofelements2;
471 md2.mesh.numberofvertices=numberofvertices2;
472 md2.mesh.elements=elements_2;
473
474 %mesh.uppervertex mesh.lowervertex
475 if isa(md1.mesh,'mesh3dprisms'),
476 md2.mesh.uppervertex=md1.mesh.uppervertex(pos_node);
477 pos=find(~isnan(md2.mesh.uppervertex));
478 md2.mesh.uppervertex(pos)=Pnode(md2.mesh.uppervertex(pos));
479
480 md2.mesh.lowervertex=md1.mesh.lowervertex(pos_node);
481 pos=find(~isnan(md2.mesh.lowervertex));
482 md2.mesh.lowervertex(pos)=Pnode(md2.mesh.lowervertex(pos));
483
484 md2.mesh.upperelements=md1.mesh.upperelements(pos_elem);
485 pos=find(~isnan(md2.mesh.upperelements));
486 md2.mesh.upperelements(pos)=Pelem(md2.mesh.upperelements(pos));
487
488 md2.mesh.lowerelements=md1.mesh.lowerelements(pos_elem);
489 pos=find(~isnan(md2.mesh.lowerelements));
490 md2.mesh.lowerelements(pos)=Pelem(md2.mesh.lowerelements(pos));
491 end
492
493 %Initial 2d mesh
494 if isa(md1.mesh,'mesh3dprisms'),
495 flag_elem_2d=flag_elem(1:md1.mesh.numberofelements2d);
496 pos_elem_2d=find(flag_elem_2d);
497 flag_node_2d=flag_node(1:md1.mesh.numberofvertices2d);
498 pos_node_2d=find(flag_node_2d);
499
500 md2.mesh.numberofelements2d=length(pos_elem_2d);
501 md2.mesh.numberofvertices2d=length(pos_node_2d);
502 md2.mesh.elements2d=md1.mesh.elements2d(pos_elem_2d,:);
503 md2.mesh.elements2d(:,1)=Pnode(md2.mesh.elements2d(:,1));
504 md2.mesh.elements2d(:,2)=Pnode(md2.mesh.elements2d(:,2));
505 md2.mesh.elements2d(:,3)=Pnode(md2.mesh.elements2d(:,3));
506
507 md2.mesh.x2d=md1.mesh.x(pos_node_2d);
508 md2.mesh.y2d=md1.mesh.y(pos_node_2d);
509 end
510
511 %Edges
512 if(dimension(md.mesh)==2),
513 if size(md2.mesh.edges,2)>1, %do not use ~isnan because there are some NaNs...
514 %renumber first two columns
515 pos=find(md2.mesh.edges(:,4)~=-1);
516 md2.mesh.edges(: ,1)=Pnode(md2.mesh.edges(:,1));
517 md2.mesh.edges(: ,2)=Pnode(md2.mesh.edges(:,2));
518 md2.mesh.edges(: ,3)=Pelem(md2.mesh.edges(:,3));
519 md2.mesh.edges(pos,4)=Pelem(md2.mesh.edges(pos,4));
520 %remove edges when the 2 vertices are not in the domain.
521 md2.mesh.edges=md2.mesh.edges(find(md2.mesh.edges(:,1) & md2.mesh.edges(:,2)),:);
522 %Replace all zeros by -1 in the last two columns
523 pos=find(md2.mesh.edges(:,3)==0);
524 md2.mesh.edges(pos,3)=-1;
525 pos=find(md2.mesh.edges(:,4)==0);
526 md2.mesh.edges(pos,4)=-1;
527 %Invert -1 on the third column with last column (Also invert first two columns!!)
528 pos=find(md2.mesh.edges(:,3)==-1);
529 md2.mesh.edges(pos,3)=md2.mesh.edges(pos,4);
530 md2.mesh.edges(pos,4)=-1;
531 values=md2.mesh.edges(pos,2);
532 md2.mesh.edges(pos,2)=md2.mesh.edges(pos,1);
533 md2.mesh.edges(pos,1)=values;
534 %Finally remove edges that do not belong to any element
535 pos=find(md2.mesh.edges(:,3)==-1 & md2.mesh.edges(:,4)==-1);
536 md2.mesh.edges(pos,:)=[];
537 end
538 end
539
540 %Penalties
541 if ~isnan(md2.stressbalance.vertex_pairing),
542 for i=1:size(md1.stressbalance.vertex_pairing,1);
543 md2.stressbalance.vertex_pairing(i,:)=Pnode(md1.stressbalance.vertex_pairing(i,:));
544 end
545 md2.stressbalance.vertex_pairing=md2.stressbalance.vertex_pairing(find(md2.stressbalance.vertex_pairing(:,1)),:);
546 end
547 if ~isnan(md2.masstransport.vertex_pairing),
548 for i=1:size(md1.masstransport.vertex_pairing,1);
549 md2.masstransport.vertex_pairing(i,:)=Pnode(md1.masstransport.vertex_pairing(i,:));
550 end
551 md2.masstransport.vertex_pairing=md2.masstransport.vertex_pairing(find(md2.masstransport.vertex_pairing(:,1)),:);
552 end
553
554 %recreate segments
555 if isa(md1.mesh,'mesh2d') | isa(md1.mesh','mesh3dsurface'),
556 md2.mesh.vertexconnectivity=NodeConnectivity(md2.mesh.elements,md2.mesh.numberofvertices);
557 md2.mesh.elementconnectivity=ElementConnectivity(md2.mesh.elements,md2.mesh.vertexconnectivity);
558 md2.mesh.segments=contourenvelope(md2.mesh);
559 md2.mesh.vertexonboundary=zeros(numberofvertices2,1); md2.mesh.vertexonboundary(md2.mesh.segments(:,1:2))=1;
560 else
561 %First do the connectivity for the contourenvelope in 2d
562 md2.mesh.vertexconnectivity=NodeConnectivity(md2.mesh.elements2d,md2.mesh.numberofvertices2d);
563 md2.mesh.elementconnectivity=ElementConnectivity(md2.mesh.elements2d,md2.mesh.vertexconnectivity);
564 segments=contourenvelope(md2.mesh);
565 md2.mesh.vertexonboundary=zeros(numberofvertices2/md2.mesh.numberoflayers,1); md2.mesh.vertexonboundary(segments(:,1:2))=1;
566 md2.mesh.vertexonboundary=repmat(md2.mesh.vertexonboundary,md2.mesh.numberoflayers,1);
567 %Then do it for 3d as usual
568 md2.mesh.vertexconnectivity=NodeConnectivity(md2.mesh.elements,md2.mesh.numberofvertices);
569 md2.mesh.elementconnectivity=ElementConnectivity(md2.mesh.elements,md2.mesh.vertexconnectivity);
570 end
571
572 %Boundary conditions: Dirichlets on new boundary
573 %Catch the elements that have not been extracted
574 orphans_elem=find(~flag_elem);
575 orphans_node=unique(md1.mesh.elements(orphans_elem,:))';
576 %Figure out which node are on the boundary between md2 and md1
577 nodestoflag1=intersect(orphans_node,pos_node);
578 nodestoflag2=Pnode(nodestoflag1);
579 if numel(md1.stressbalance.spcvx)>1 & numel(md1.stressbalance.spcvy)>2 & numel(md1.stressbalance.spcvz)>2,
580 if numel(md1.inversion.vx_obs)>1 & numel(md1.inversion.vy_obs)>1
581 md2.stressbalance.spcvx(nodestoflag2)=md2.inversion.vx_obs(nodestoflag2);
582 md2.stressbalance.spcvy(nodestoflag2)=md2.inversion.vy_obs(nodestoflag2);
583 else
584 md2.stressbalance.spcvx(nodestoflag2)=NaN;
585 md2.stressbalance.spcvy(nodestoflag2)=NaN;
586 disp(' ')
587 disp('!! extract warning: spc values should be checked !!')
588 disp(' ')
589 end
590 %put 0 for vz
591 md2.stressbalance.spcvz(nodestoflag2)=0;
592 end
593 if ~isnan(md1.thermal.spctemperature),
594 md2.thermal.spctemperature(nodestoflag2,1)=1;
595 end
596
597 %Results fields
598 if isstruct(md1.results),
599 md2.results=struct();
600 solutionfields=fields(md1.results);
601 for i=1:length(solutionfields),
602 if isstruct(md1.results.(solutionfields{i}))
603 %get subfields
604 solutionsubfields=fields(md1.results.(solutionfields{i}));
605 for j=1:length(solutionsubfields),
606 field=md1.results.(solutionfields{i}).(solutionsubfields{j});
607 if length(field)==numberofvertices1,
608 md2.results.(solutionfields{i}).(solutionsubfields{j})=field(pos_node);
609 elseif length(field)==numberofelements1,
610 md2.results.(solutionfields{i}).(solutionsubfields{j})=field(pos_elem);
611 else
612 md2.results.(solutionfields{i}).(solutionsubfields{j})=field;
613 end
614 end
615 else
616 field=md1.results.(solutionfields{i});
617 if length(field)==numberofvertices1,
618 md2.results.(solutionfields{i})=field(pos_node);
619 elseif length(field)==numberofelements1,
620 md2.results.(solutionfields{i})=field(pos_elem);
621 else
622 md2.results.(solutionfields{i})=field;
623 end
624 end
625 end
626 end
627
628 %OutputDefinitions fields
629 for i=1:length(md1.outputdefinition.definitions),
630 if isobject(md1.outputdefinition.definitions{i})
631 %get subfields
632 solutionsubfields=fields(md1.outputdefinition.definitions{i});
633 for j=1:length(solutionsubfields),
634 field=md1.outputdefinition.definitions{i}.(solutionsubfields{j});
635 if length(field)==numberofvertices1,
636 md2.outputdefinition.definitions{i}.(solutionsubfields{j})=field(pos_node);
637 elseif length(field)==numberofelements1,
638 md2.outputdefinition.definitions{i}.(solutionsubfields{j})=field(pos_elem);
639 end
640 end
641 end
642 end
643
644 %Keep track of pos_node and pos_elem
645 md2.mesh.extractedvertices=pos_node;
646 md2.mesh.extractedelements=pos_elem;
647 end % }}}
648 function md = extrude(md,varargin) % {{{
649 %EXTRUDE - vertically extrude a 2d mesh
650 %
651 % vertically extrude a 2d mesh and create corresponding 3d mesh.
652 % The vertical distribution can:
653 % - follow a polynomial law
654 % - follow two polynomial laws, one for the lower part and one for the upper part of the mesh
655 % - be discribed by a list of coefficients (between 0 and 1)
656 %
657 %
658 % Usage:
659 % md=extrude(md,numlayers,extrusionexponent);
660 % md=extrude(md,numlayers,lowerexponent,upperexponent);
661 % md=extrude(md,listofcoefficients);
662 %
663 % Example:
664 % md=extrude(md,15,1.3);
665 % md=extrude(md,15,1.3,1.2);
666 % md=extrude(md,[0 0.2 0.5 0.7 0.9 0.95 1]);
667 %
668 % See also: MODELEXTRACT, COLLAPSE
669
670 %some checks on list of arguments
671 if ((nargin>4) | (nargin<2) | (nargout~=1)),
672 help extrude;
673 error('extrude error message');
674 end
675
676 %Extrude the mesh
677 if nargin==2, %list of coefficients
678 clist=varargin{1};
679 if any(clist<0) | any(clist>1),
680 error('extrusioncoefficients must be between 0 and 1');
681 end
682 extrusionlist=sort(unique([clist(:);0;1]));
683 numlayers=length(extrusionlist);
684 elseif nargin==3, %one polynomial law
685 if varargin{2}<=0,
686 help extrude;
687 error('extrusionexponent must be >=0');
688 end
689 numlayers=varargin{1};
690 extrusionlist=((0:1:numlayers-1)/(numlayers-1)).^varargin{2};
691 elseif nargin==4, %two polynomial laws
692 numlayers=varargin{1};
693 lowerexp=varargin{2};
694 upperexp=varargin{3};
695
696 if varargin{2}<=0 | varargin{3}<=0,
697 help extrude;
698 error('lower and upper extrusionexponents must be >=0');
699 end
700
701 lowerextrusionlist=[(0:2/(numlayers-1):1).^lowerexp]/2;
702 upperextrusionlist=[(0:2/(numlayers-1):1).^upperexp]/2;
703 extrusionlist=sort(unique([lowerextrusionlist 1-upperextrusionlist]));
704
705 end
706
707 if numlayers<2,
708 error('number of layers should be at least 2');
709 end
710 if strcmp(md.mesh.domaintype(),'3D')
711 error('Cannot extrude a 3d mesh (extrude cannot be called more than once)');
712 end
713
714 %Initialize with the 2d mesh
715 mesh2d = md.mesh;
716 md.mesh=mesh3dprisms();
717 md.mesh.x = mesh2d.x;
718 md.mesh.y = mesh2d.y;
719 md.mesh.elements = mesh2d.elements;
720 md.mesh.numberofelements = mesh2d.numberofelements;
721 md.mesh.numberofvertices = mesh2d.numberofvertices;
722
723 md.mesh.lat = mesh2d.lat;
724 md.mesh.long = mesh2d.long;
725 md.mesh.epsg = mesh2d.epsg;
726
727 md.mesh.vertexonboundary = mesh2d.vertexonboundary;
728 md.mesh.vertexconnectivity = mesh2d.vertexconnectivity;
729 md.mesh.elementconnectivity = mesh2d.elementconnectivity;
730 md.mesh.average_vertex_connectivity = mesh2d.average_vertex_connectivity;
731
732 md.mesh.extractedvertices = mesh2d.extractedvertices;
733 md.mesh.extractedelements = mesh2d.extractedelements;
734
735 x3d=[];
736 y3d=[];
737 z3d=[]; %the lower node is on the bed
738 thickness3d=md.geometry.thickness; %thickness and bed for these nodes
739 bed3d=md.geometry.base;
740
741 %Create the new layers
742 for i=1:numlayers,
743 x3d=[x3d; md.mesh.x];
744 y3d=[y3d; md.mesh.y];
745 %nodes are distributed between bed and surface accordingly to the given exponent
746 z3d=[z3d; bed3d+thickness3d*extrusionlist(i)];
747 end
748 number_nodes3d=size(x3d,1); %number of 3d nodes for the non extruded part of the mesh
749
750 %Extrude elements
751 elements3d=[];
752 for i=1:numlayers-1,
753 elements3d=[elements3d;[md.mesh.elements+(i-1)*md.mesh.numberofvertices md.mesh.elements+i*md.mesh.numberofvertices]]; %Create the elements of the 3d mesh for the non extruded part
754 end
755 number_el3d=size(elements3d,1); %number of 3d nodes for the non extruded part of the mesh
756
757 %Keep a trace of lower and upper nodes
758 lowervertex=NaN*ones(number_nodes3d,1);
759 uppervertex=NaN*ones(number_nodes3d,1);
760 lowervertex(md.mesh.numberofvertices+1:end)=1:(numlayers-1)*md.mesh.numberofvertices;
761 uppervertex(1:(numlayers-1)*md.mesh.numberofvertices)=md.mesh.numberofvertices+1:number_nodes3d;
762 md.mesh.lowervertex=lowervertex;
763 md.mesh.uppervertex=uppervertex;
764
765 %same for lower and upper elements
766 lowerelements=NaN*ones(number_el3d,1);
767 upperelements=NaN*ones(number_el3d,1);
768 lowerelements(md.mesh.numberofelements+1:end)=1:(numlayers-2)*md.mesh.numberofelements;
769 upperelements(1:(numlayers-2)*md.mesh.numberofelements)=md.mesh.numberofelements+1:(numlayers-1)*md.mesh.numberofelements;
770 md.mesh.lowerelements=lowerelements;
771 md.mesh.upperelements=upperelements;
772
773 %Save old mesh
774 md.mesh.x2d=md.mesh.x;
775 md.mesh.y2d=md.mesh.y;
776 md.mesh.elements2d=md.mesh.elements;
777 md.mesh.numberofelements2d=md.mesh.numberofelements;
778 md.mesh.numberofvertices2d=md.mesh.numberofvertices;
779
780 %Build global 3d mesh
781 md.mesh.elements=elements3d;
782 md.mesh.x=x3d;
783 md.mesh.y=y3d;
784 md.mesh.z=z3d;
785 md.mesh.numberofelements=number_el3d;
786 md.mesh.numberofvertices=number_nodes3d;
787 md.mesh.numberoflayers=numlayers;
788
789 %Ok, now deal with the other fields from the 2d mesh:
790
791 %bedinfo and surface info
792 md.mesh.vertexonbase=project3d(md,'vector',ones(md.mesh.numberofvertices2d,1),'type','node','layer',1);
793 md.mesh.vertexonsurface=project3d(md,'vector',ones(md.mesh.numberofvertices2d,1),'type','node','layer',md.mesh.numberoflayers);
794 md.mesh.vertexonboundary=project3d(md,'vector',md.mesh.vertexonboundary,'type','node');
795
796 %lat long
797 md.mesh.lat=project3d(md,'vector',md.mesh.lat,'type','node');
798 md.mesh.long=project3d(md,'vector',md.mesh.long,'type','node');
799
800 md.geometry=extrude(md.geometry,md);
801 md.friction = extrude(md.friction,md);
802 md.inversion = extrude(md.inversion,md);
803 md.smb = extrude(md.smb,md);
804 md.initialization = extrude(md.initialization,md);
805
806 md.flowequation=md.flowequation.extrude(md);
807 md.stressbalance=extrude(md.stressbalance,md);
808 md.thermal=md.thermal.extrude(md);
809 md.masstransport=md.masstransport.extrude(md);
810 md.levelset=extrude(md.levelset,md);
811 md.calving=extrude(md.calving,md);
812 md.hydrology = extrude(md.hydrology,md);
813
814 %connectivity
815 if ~isnan(md.mesh.elementconnectivity)
816 md.mesh.elementconnectivity=repmat(md.mesh.elementconnectivity,numlayers-1,1);
817 md.mesh.elementconnectivity(find(md.mesh.elementconnectivity==0))=NaN;
818 for i=2:numlayers-1,
819 md.mesh.elementconnectivity((i-1)*md.mesh.numberofelements2d+1:(i)*md.mesh.numberofelements2d,:)...
820 =md.mesh.elementconnectivity((i-1)*md.mesh.numberofelements2d+1:(i)*md.mesh.numberofelements2d,:)+md.mesh.numberofelements2d;
821 end
822 md.mesh.elementconnectivity(find(isnan(md.mesh.elementconnectivity)))=0;
823 end
824
825 md.materials=extrude(md.materials,md);
826 md.damage=extrude(md.damage,md);
827 md.mask=extrude(md.mask,md);
828 md.qmu=extrude(md.qmu,md);
829 md.basalforcings=extrude(md.basalforcings,md);
830 md.outputdefinition=extrude(md.outputdefinition,md);
831
832 %increase connectivity if less than 25:
833 if md.mesh.average_vertex_connectivity<=25,
834 md.mesh.average_vertex_connectivity=100;
835 end
836 end % }}}
837 function md = structtomodel(md,structmd) % {{{
838
839 if ~isstruct(structmd) error('input model is not a structure'); end
840
841 %loaded model is a struct, initialize output and recover all fields
842 md = structtoobj(model,structmd);
843
844 %Old field now classes
845 if (isfield(structmd,'timestepping') & isnumeric(md.timestepping)), md.timestepping=timestepping(); end
846 if (isfield(structmd,'mask') & isnumeric(md.mask)),md.mask=mask(); end
847
848 %Field name change
849 if isfield(structmd,'drag'), md.friction.coefficient=structmd.drag; end
850 if isfield(structmd,'p'), md.friction.p=structmd.p; end
851 if isfield(structmd,'q'), md.friction.q=structmd.p; end
852 if isfield(structmd,'melting'), md.basalforcings.floatingice_melting_rate=structmd.melting; end
853 if isfield(structmd,'melting_rate'), md.basalforcings.floatingice_melting_rate=structmd.melting_rate; end
854 if isfield(structmd,'melting_rate'), md.basalforcings.groundedice_melting_rate=structmd.melting_rate; end
855 if isfield(structmd,'accumulation'), md.smb.mass_balance=structmd.accumulation; end
856 if isfield(structmd,'numberofgrids'), md.mesh.numberofvertices=structmd.numberofgrids; end
857 if isfield(structmd,'numberofgrids2d'), md.mesh.numberofvertices2d=structmd.numberofgrids2d; end
858 if isfield(structmd,'uppergrids'), md.mesh.uppervertex=structmd.uppergrids; end
859 if isfield(structmd,'lowergrids'), md.mesh.lowervertex=structmd.lowergrids; end
860 if isfield(structmd,'gridonbase'), md.mesh.vertexonbase=structmd.gridonbase; end
861 if isfield(structmd,'gridonsurface'), md.mesh.vertexonsurface=structmd.gridonsurface; end
862 if isfield(structmd,'extractedgrids'), md.mesh.extractedvertices=structmd.extractedgrids; end
863 if isfield(structmd,'gridonboundary'), md.mesh.vertexonboundary=structmd.gridonboundary; end
864 if isfield(structmd,'petscoptions') & ~isempty(structmd.petscoptions), md.toolkits=structmd.petscoptions; end
865 if isfield(structmd,'g'), md.constants.g=structmd.g; end
866 if isfield(structmd,'yts'), md.constants.yts=structmd.yts; end
867 if isfield(structmd,'surface_mass_balance'), md.smb.mass_balance=structmd.surface_mass_balance; end
868 if isfield(structmd,'basal_melting_rate'), md.basalforcings.floatingice_melting_rate=structmd.basal_melting_rate; end
869 if isfield(structmd,'geothermalflux'), md.basalforcings.geothermalflux=structmd.geothermalflux; end
870 if isfield(structmd,'drag'), md.friction.coefficient=structmd.drag; end
871 if isfield(structmd,'drag_coefficient'), md.friction.coefficient=structmd.drag_coefficient; end
872 if isfield(structmd,'drag_p'), md.friction.p=structmd.drag_p; end
873 if isfield(structmd,'drag_q'), md.friction.q=structmd.drag_q; end
874 if isfield(structmd,'riftproperties'), %old implementation
875 md.rifts=rifts();
876 md.rifts.riftproperties=structmd.riftproperties;
877 md.rifts.riftstruct=structmd.rifts;
878 md.rifts.riftproperties=structmd.riftinfo;
879 end
880 if isfield(structmd,'bamg'), md.private.bamg=structmd.bamg; end
881 if isfield(structmd,'lowmem'), md.settings.lowmem=structmd.lowmem; end
882 if isfield(structmd,'io_gather'), md.settings.io_gather=structmd.io_gather; end
883 if isfield(structmd,'spcwatercolumn'), md.hydrology.spcwatercolumn=structmd.spcwatercolumn; end
884 if isfield(structmd,'hydro_n'), md.hydrology.n=structmd.hydro_n; end
885 if isfield(structmd,'hydro_p'), md.hydrology.p=structmd.hydro_p; end
886 if isfield(structmd,'hydro_q'), md.hydrology.q=structmd.hydro_q; end
887 if isfield(structmd,'hydro_CR'), md.hydrology.CR=structmd.hydro_CR; end
888 if isfield(structmd,'hydro_kn'), md.hydrology.kn=structmd.hydro_kn; end
889 if isfield(structmd,'spctemperature'), md.thermal.spctemperature=structmd.spctemperature; end
890 if isfield(structmd,'min_thermal_constraints'), md.thermal.penalty_threshold=structmd.min_thermal_constraints; end
891 if isfield(structmd,'artificial_diffusivity'), md.thermal.stabilization=structmd.artificial_diffusivity; end
892 if isfield(structmd,'max_nonlinear_iterations'), md.thermal.maxiter=structmd.max_nonlinear_iterations; end
893 if isfield(structmd,'stabilize_constraints'), md.thermal.penalty_lock=structmd.stabilize_constraints; end
894 if isfield(structmd,'penalty_offset'), md.thermal.penalty_factor=structmd.penalty_offset; end
895 if isfield(structmd,'name'), md.miscellaneous.name=structmd.name; end
896 if isfield(structmd,'notes'), md.miscellaneous.notes=structmd.notes; end
897 if isfield(structmd,'dummy'), md.miscellaneous.dummy=structmd.dummy; end
898 if isfield(structmd,'dt'), md.timestepping.time_step=structmd.dt; end
899 if isfield(structmd,'ndt'), md.timestepping.final_time=structmd.ndt; end
900 if isfield(structmd,'time_adapt'), md.timestepping.time_adapt=structmd.time_adapt; end
901 if isfield(structmd,'cfl_coefficient'), md.timestepping.cfl_coefficient=structmd.cfl_coefficient; end
902 if isfield(structmd,'spcthickness'), md.masstransport.spcthickness=structmd.spcthickness; end
903 if isfield(structmd,'artificial_diffusivity'), md.masstransport.stabilization=structmd.artificial_diffusivity; end
904 if isfield(structmd,'hydrostatic_adjustment'), md.masstransport.hydrostatic_adjustment=structmd.hydrostatic_adjustment; end
905 if isfield(structmd,'penalties'), md.masstransport.vertex_pairing=structmd.penalties; end
906 if isfield(structmd,'penalty_offset'), md.masstransport.penalty_factor=structmd.penalty_offset; end
907 if isfield(structmd,'B'), md.materials.rheology_B=structmd.B; end
908 if isfield(structmd,'n'), md.materials.rheology_n=structmd.n; end
909 if isfield(structmd,'rheology_B'), md.materials.rheology_B=structmd.rheology_B; end
910 if isfield(structmd,'rheology_n'), md.materials.rheology_n=structmd.rheology_n; end
911 if isfield(structmd,'rheology_Z'), md.damage.D=(1-structmd.rheology_Z); end
912 if isfield(structmd,'spcthickness'), md.balancethickness.spcthickness=structmd.spcthickness; end
913 if isfield(structmd,'artificial_diffusivity'), md.balancethickness.stabilization=structmd.artificial_diffusivity; end
914 if isfield(structmd,'dhdt'), md.balancethickness.thickening_rate=structmd.dhdt; end
915 if isfield(structmd,'isSIA'), md.flowequation.isSIA=structmd.isSIA; end
916 if isfield(structmd,'isFS'), md.flowequation.isFS=structmd.isFS; end
917 if isfield(structmd,'elements_type'), md.flowequation.element_equation=structmd.elements_type; end
918 if isfield(structmd,'vertices_type'), md.flowequation.vertex_equation=structmd.vertices_type; end
919 if isfield(structmd,'eps_rel'), md.steadystate.reltol=structmd.eps_rel; end
920 if isfield(structmd,'max_steadystate_iterations'), md.steadystate.maxiter=structmd.max_steadystate_iterations; end
921 if isfield(structmd,'isdiagnostic'), md.transient.isstressbalance=structmd.isdiagnostic; end
922 if isfield(structmd,'isprognostic'), md.transient.ismasstransport=structmd.isprognostic; end
923 if isfield(structmd,'isthermal'), md.transient.isthermal=structmd.isthermal; end
924 if isfield(structmd,'control_analysis'), md.inversion.iscontrol=structmd.control_analysis; end
925 if isfield(structmd,'weights'), md.inversion.cost_functions_coefficients=structmd.weights; end
926 if isfield(structmd,'nsteps'), md.inversion.nsteps=structmd.nsteps; end
927 if isfield(structmd,'maxiter_per_step'), md.inversion.maxiter_per_step=structmd.maxiter_per_step; end
928 if isfield(structmd,'cm_min'), md.inversion.min_parameters=structmd.cm_min; end
929 if isfield(structmd,'cm_max'), md.inversion.max_parameters=structmd.cm_max; end
930 if isfield(structmd,'vx_obs'), md.inversion.vx_obs=structmd.vx_obs; end
931 if isfield(structmd,'vy_obs'), md.inversion.vy_obs=structmd.vy_obs; end
932 if isfield(structmd,'vel_obs'), md.inversion.vel_obs=structmd.vel_obs; end
933 if isfield(structmd,'thickness_obs'), md.inversion.thickness_obs=structmd.thickness_obs; end
934 if isfield(structmd,'vx'), md.initialization.vx=structmd.vx; end
935 if isfield(structmd,'vy'), md.initialization.vy=structmd.vy; end
936 if isfield(structmd,'vz'), md.initialization.vz=structmd.vz; end
937 if isfield(structmd,'vel'), md.initialization.vel=structmd.vel; end
938 if isfield(structmd,'pressure'), md.initialization.pressure=structmd.pressure; end
939 if isfield(structmd,'temperature'), md.initialization.temperature=structmd.temperature; end
940 if isfield(structmd,'waterfraction'), md.initialization.waterfraction=structmd.waterfraction; end
941 if isfield(structmd,'watercolumn'), md.initialization.watercolumn=structmd.watercolumn; end
942 if isfield(structmd,'surface'), md.geometry.surface=structmd.surface; end
943 if isfield(structmd,'bed'), md.geometry.base=structmd.bed; end
944 if isfield(structmd,'thickness'), md.geometry.thickness=structmd.thickness; end
945 if isfield(structmd,'bathymetry'), md.geometry.bed=structmd.bathymetry; end
946 if isfield(structmd,'thickness_coeff'), md.geometry.hydrostatic_ratio=structmd.thickness_coeff; end
947 if isfield(structmd,'connectivity'), md.mesh.average_vertex_connectivity=structmd.connectivity; end
948 if isfield(structmd,'extractednodes'), md.mesh.extractedvertices=structmd.extractednodes; end
949 if isfield(structmd,'extractedelements'), md.mesh.extractedelements=structmd.extractedelements; end
950 if isfield(structmd,'nodeonboundary'), md.mesh.vertexonboundary=structmd.nodeonboundary; end
951 if isfield(structmd,'lat'), md.mesh.lat=structmd.lat; end
952 if isfield(structmd,'long'), md.mesh.long=structmd.long; end
953 if isfield(structmd,'segments'), md.mesh.segments=structmd.segments; end
954 if isfield(structmd,'segmentmarkers'), md.mesh.segmentmarkers=structmd.segmentmarkers; end
955 if isfield(structmd,'numlayers'), md.mesh.numberoflayers=structmd.numlayers; end
956 if isfield(structmd,'numberofelements'), md.mesh.numberofelements=structmd.numberofelements; end
957 if isfield(structmd,'numberofvertices'), md.mesh.numberofvertices=structmd.numberofvertices; end
958 if isfield(structmd,'numberofnodes'), md.mesh.numberofvertices=structmd.numberofnodes; end
959 if isfield(structmd,'numberofedges'), md.mesh.numberofedges=structmd.numberofedges; end
960 if isfield(structmd,'numberofelements2d'), md.mesh.numberofelements2d=structmd.numberofelements2d; end
961 if isfield(structmd,'numberofnodes2d'), md.mesh.numberofvertices2d=structmd.numberofnodes2d; end
962 if isfield(structmd,'nodeconnectivity'), md.mesh.vertexconnectivity=structmd.nodeconnectivity; end
963 if isfield(structmd,'elementconnectivity'), md.mesh.elementconnectivity=structmd.elementconnectivity; end
964 if isfield(structmd,'uppernodes'), md.mesh.uppervertex=structmd.uppernodes; end
965 if isfield(structmd,'lowernodes'), md.mesh.lowervertex=structmd.lowernodes; end
966 if isfield(structmd,'upperelements'), md.mesh.upperelements=structmd.upperelements; end
967 if isfield(structmd,'lowerelements'), md.mesh.lowerelements=structmd.lowerelements; end
968 if isfield(structmd,'nodeonsurface'), md.mesh.vertexonsurface=structmd.nodeonsurface; end
969 if isfield(structmd,'nodeonbase'), md.mesh.vertexonbase=structmd.nodeonbase; end
970 if isfield(structmd,'elements2d'), md.mesh.elements2d=structmd.elements2d; end
971 if isfield(structmd,'y2d'), md.mesh.y2d=structmd.y2d; end
972 if isfield(structmd,'x2d'), md.mesh.x2d=structmd.x2d; end
973 if isfield(structmd,'elements'), md.mesh.elements=structmd.elements; end
974 if isfield(structmd,'edges'),
975 md.mesh.edges=structmd.edges;
976 md.mesh.edges(isnan(md.mesh.edges))=-1;
977 end
978 if isfield(structmd,'y'), md.mesh.y=structmd.y; end
979 if isfield(structmd,'x'), md.mesh.x=structmd.x; end
980 if isfield(structmd,'z'), md.mesh.z=structmd.z; end
981 if isfield(structmd,'diagnostic_ref'), md.stressbalance.referential=structmd.diagnostic_ref; end
982 if isfield(structmd,'npart'); md.qmu.numberofpartitions=structmd.npart; end
983 if isfield(structmd,'part'); md.qmu.partition=structmd.part; end
984
985 if isnumeric(md.verbose),
986 md.verbose=verbose;
987 end
988
989 if isfield(structmd,'spcvelocity'),
990 md.stressbalance.spcvx=NaN*ones(md.mesh.numberofvertices,1);
991 md.stressbalance.spcvy=NaN*ones(md.mesh.numberofvertices,1);
992 md.stressbalance.spcvz=NaN*ones(md.mesh.numberofvertices,1);
993 pos=find(structmd.spcvelocity(:,1)); md.stressbalance.spcvx(pos)=structmd.spcvelocity(pos,4);
994 pos=find(structmd.spcvelocity(:,2)); md.stressbalance.spcvy(pos)=structmd.spcvelocity(pos,5);
995 pos=find(structmd.spcvelocity(:,3)); md.stressbalance.spcvz(pos)=structmd.spcvelocity(pos,6);
996 end
997 if isfield(structmd,'spcvx'),
998 md.stressbalance.spcvx=NaN*ones(md.mesh.numberofvertices,1);
999 pos=find(~isnan(structmd.spcvx)); md.stressbalance.spcvx(pos)=structmd.spcvx(pos);
1000 end
1001 if isfield(structmd,'spcvy'),
1002 md.stressbalance.spcvy=NaN*ones(md.mesh.numberofvertices,1);
1003 pos=find(~isnan(structmd.spcvy)); md.stressbalance.spcvy(pos)=structmd.spcvy(pos);
1004 end
1005 if isfield(structmd,'spcvz'),
1006 md.stressbalance.spcvz=NaN*ones(md.mesh.numberofvertices,1);
1007 pos=find(~isnan(structmd.spcvz)); md.stressbalance.spcvz(pos)=structmd.spcvz(pos);
1008 end
1009 if isfield(structmd,'pressureload'),
1010 if ~isempty(structmd.pressureload) & ismember(structmd.pressureload(end,end),[118 119 120]),
1011 pos=find(structmd.pressureload(:,end)==120); md.stressbalance.icefront(pos,end)=0;
1012 pos=find(structmd.pressureload(:,end)==118); md.stressbalance.icefront(pos,end)=1;
1013 pos=find(structmd.pressureload(:,end)==119); md.stressbalance.icefront(pos,end)=2;
1014 end
1015 end
1016 if isfield(structmd,'elements_type') & structmd.elements_type(end,end)>50,
1017 pos=find(structmd.elements_type==59); md.flowequation.element_equation(pos,end)=0;
1018 pos=find(structmd.elements_type==55); md.flowequation.element_equation(pos,end)=1;
1019 pos=find(structmd.elements_type==56); md.flowequation.element_equation(pos,end)=2;
1020 pos=find(structmd.elements_type==60); md.flowequation.element_equation(pos,end)=3;
1021 pos=find(structmd.elements_type==62); md.flowequation.element_equation(pos,end)=4;
1022 pos=find(structmd.elements_type==57); md.flowequation.element_equation(pos,end)=5;
1023 pos=find(structmd.elements_type==58); md.flowequation.element_equation(pos,end)=6;
1024 pos=find(structmd.elements_type==61); md.flowequation.element_equation(pos,end)=7;
1025 end
1026 if isfield(structmd,'vertices_type') & structmd.vertices_type(end,end)>50,
1027 pos=find(structmd.vertices_type==59); md.flowequation.vertex_equation(pos,end)=0;
1028 pos=find(structmd.vertices_type==55); md.flowequation.vertex_equation(pos,end)=1;
1029 pos=find(structmd.vertices_type==56); md.flowequation.vertex_equation(pos,end)=2;
1030 pos=find(structmd.vertices_type==60); md.flowequation.vertex_equation(pos,end)=3;
1031 pos=find(structmd.vertices_type==62); md.flowequation.vertex_equation(pos,end)=4;
1032 pos=find(structmd.vertices_type==57); md.flowequation.vertex_equation(pos,end)=5;
1033 pos=find(structmd.vertices_type==58); md.flowequation.vertex_equation(pos,end)=6;
1034 pos=find(structmd.vertices_type==61); md.flowequation.vertex_equation(pos,end)=7;
1035 end
1036 if isfield(structmd,'rheology_law') & isnumeric(structmd.rheology_law),
1037 if (structmd.rheology_law==272), md.materials.rheology_law='None'; end
1038 if (structmd.rheology_law==368), md.materials.rheology_law='Paterson'; end
1039 if (structmd.rheology_law==369), md.materials.rheology_law='Arrhenius'; end
1040 end
1041 if isfield(structmd,'groundingline_migration') & isnumeric(structmd.groundingline_migration),
1042 if (structmd.groundingline_migration==272), md.groundingline.migration='None'; end
1043 if (structmd.groundingline_migration==273), md.groundingline.migration='AggressiveMigration'; end
1044 if (structmd.groundingline_migration==274), md.groundingline.migration='SoftMigration'; end
1045 end
1046 if isfield(structmd,'control_type') & isnumeric(structmd.control_type),
1047 if (structmd.control_type==143), md.inversion.control_parameters={'FrictionCoefficient'}; end
1048 if (structmd.control_type==190), md.inversion.control_parameters={'RheologyBbar'}; end
1049 if (structmd.control_type==147), md.inversion.control_parameters={'Thickeningrate'}; end
1050 end
1051 if isfield(structmd,'cm_responses') & ismember(structmd.cm_responses(end,end),[165:170 383 388 389]),
1052 pos=find(structmd.cm_responses==166); md.inversion.cost_functions(pos)=101;
1053 pos=find(structmd.cm_responses==167); md.inversion.cost_functions(pos)=102;
1054 pos=find(structmd.cm_responses==168); md.inversion.cost_functions(pos)=103;
1055 pos=find(structmd.cm_responses==169); md.inversion.cost_functions(pos)=104;
1056 pos=find(structmd.cm_responses==170); md.inversion.cost_functions(pos)=105;
1057 pos=find(structmd.cm_responses==165); md.inversion.cost_functions(pos)=201;
1058 pos=find(structmd.cm_responses==389); md.inversion.cost_functions(pos)=501;
1059 pos=find(structmd.cm_responses==388); md.inversion.cost_functions(pos)=502;
1060 pos=find(structmd.cm_responses==382); md.inversion.cost_functions(pos)=503;
1061 end
1062
1063 if isfield(structmd,'artificial_diffusivity') & structmd.artificial_diffusivity==2,
1064 md.thermal.stabilization=2;
1065 md.masstransport.stabilization=1;
1066 md.balancethickness.stabilization=1;
1067 end
1068 if isnumeric(md.masstransport.hydrostatic_adjustment)
1069 if md.masstransport.hydrostatic_adjustment==269,
1070 md.masstransport.hydrostatic_adjustment='Incremental';
1071 else
1072 md.masstransport.hydrostatic_adjustment='Absolute';
1073 end
1074 end
1075
1076 %New fields
1077 if ~isfield(structmd,'upperelements') & isa(md.mesh,'mesh3dprisms')
1078 md.mesh.upperelements=transpose(1:md.mesh.numberofelements)+md.mesh.numberofelements2d;
1079 md.mesh.upperelements(end-md.mesh.numberofelements2d+1:end)=NaN;
1080 end
1081 if ~isfield(structmd,'lowerelements') & isa(md.mesh,'mesh3dprisms')
1082 md.mesh.lowerelements=transpose(1:md.mesh.numberofelements)-md.mesh.numberofelements2d;
1083 md.mesh.lowerelements(1:md.mesh.numberofelements2d)=NaN;
1084 end
1085 if ~isfield(structmd,'diagnostic_ref');
1086 md.stressbalance.referential=NaN*ones(md.mesh.numberofvertices,6);
1087 end
1088 if ~isfield(structmd,'loadingforce');
1089 md.stressbalance.loadingforce=0*ones(md.mesh.numberofvertices,3);
1090 end
1091
1092 %2013 August 9
1093 if isfield(structmd,'prognostic') & isa(structmd.prognostic,'prognostic'),
1094 disp('Recovering old prognostic class');
1095 md.masstransport=masstransport(structmd.prognostic);
1096 end
1097 %2013 August 9
1098 if isfield(structmd,'diagnostic') & (isa(structmd.diagnostic,'diagnostic') || isa(structmd.diagnostic,'stressbalance')),
1099 disp('Recovering old diagnostic class');
1100 md.stressbalance=stressbalance(structmd.diagnostic);
1101 end
1102 %2014 January 9th
1103 if isfield(structmd,'surfaceforcings') & isa(md.smb,'surfaceforcings'),
1104 disp('Recovering old surfaceforcings class');
1105 mass_balance=structmd.surfaceforcings.mass_balance;
1106 md.smb=SMB();
1107 md.smb.mass_balance=mass_balance;
1108 end
1109 %2015 September 10
1110 if isfield(structmd,'surfaceforcings') & isa(structmd.surfaceforcings,'SMB'),
1111 disp('Recovering old SMB class');
1112 md.smb=SMBforcing(structmd.surfaceforcings);
1113 end
1114 if isfield(structmd,'surfaceforcings') & isa(structmd.surfaceforcings,'SMBhenning'),
1115 disp('Recovering old SMBhenning class');
1116 md.smb=SMBhenning(structmd.surfaceforcings);
1117 end
1118
1119 end% }}}
1120 function md = setdefaultparameters(md) % {{{
1121
1122 %initialize subclasses
1123 md.mesh = mesh2d();
1124 md.mask = mask();
1125 md.constants = constants();
1126 md.geometry = geometry();
1127 md.initialization = initialization();
1128 md.smb = SMBforcing();
1129 md.basalforcings = basalforcings();
1130 md.friction = friction();
1131 md.rifts = rifts();
1132 md.slr = slr();
1133 md.timestepping = timestepping();
1134 md.groundingline = groundingline();
1135 md.materials = matice();
1136 md.damage = damage();
1137 md.flowequation = flowequation();
1138 md.debug = debug();
1139 md.verbose = verbose();
1140 md.settings = settings();
1141 md.toolkits = toolkits();
1142 md.cluster = generic();
1143 md.balancethickness = balancethickness();
1144 md.stressbalance = stressbalance();
1145 md.hydrology = hydrologyshreve();
1146 md.masstransport = masstransport();
1147 md.thermal = thermal();
1148 md.steadystate = steadystate();
1149 md.transient = transient();
1150 md.levelset = levelset();
1151 md.calving = calving();
1152 md.gia = giaivins();
1153 md.esa = esa();
1154 md.autodiff = autodiff();
1155 md.inversion = inversion();
1156 md.qmu = qmu();
1157 md.amr = amr();
1158 md.radaroverlay = radaroverlay();
1159 md.results = struct();
1160 md.outputdefinition = outputdefinition();
1161 md.miscellaneous = miscellaneous();
1162 md.private = private();
1163 end
1164 %}}}
1165 function md = tetras(md,varargin) % {{{
1166 %TETRAS - split 3d prismatic mesh into 3 tetrahedrons
1167 %
1168 % Usage:
1169 % md=tetra(md)
1170
1171 if ~isa(md.mesh,'mesh3dprisms')
1172 error('mesh is not a 3d prismatic mesh');
1173 end
1174
1175 %Initialize tetra mesh
1176 md.mesh=mesh3dtetras(md.mesh);
1177
1178 %Subdivision from Philipp Furnstahl (http://studierstube.icg.tugraz.at/thesis/fuernstahl_thesis.pdf)
1179 steiner = 0;
1180 nbv = md.mesh.numberofvertices;
1181 nbt = 3*md.mesh.numberofelements;
1182 elements = zeros(nbt,4);
1183 for i=1:md.mesh.numberofelements
1184 v1=md.mesh.elements(i,1); v2=md.mesh.elements(i,2); v3=md.mesh.elements(i,3);
1185 v4=md.mesh.elements(i,4); v5=md.mesh.elements(i,5); v6=md.mesh.elements(i,6);
1186 if(min(v2,v4)<min(v1,v5) & min(v1,v6)<min(v3,v4) & min(v3,v5)<min(v2,v6)),
1187 steiner = steiner+1; nbv = nbv+1; nbt = nbt+5; v7 = nbv;
1188 md.mesh.x=[md.mesh.x; mean(md.mesh.x(md.mesh.elements(i,:)))];
1189 md.mesh.y=[md.mesh.y; mean(md.mesh.y(md.mesh.elements(i,:)))];
1190 md.mesh.z=[md.mesh.z; mean(md.mesh.z(md.mesh.elements(i,:)))];
1191 elements(3*(i-1)+1,:) = [v1 v2 v3 v7];
1192 elements(3*(i-1)+2,:) = [v1 v2 v4 v7];
1193 elements(3*(i-1)+3,:) = [v2 v4 v5 v7];
1194 elements(end+1,:) = [v2 v3 v5 v7];
1195 elements(end+1,:) = [v3 v5 v6 v7];
1196 elements(end+1,:) = [v1 v3 v6 v7];
1197 elements(end+1,:) = [v1 v4 v6 v7];
1198 elements(end+1,:) = [v4 v5 v6 v7];
1199 elseif(min(v2,v4)<min(v1,v5) & min(v1,v6)<min(v3,v4) & min(v3,v5)>min(v2,v6)),
1200 elements(3*(i-1)+1,:) = [v1 v2 v4 v6];
1201 elements(3*(i-1)+2,:) = [v2 v4 v5 v6];
1202 elements(3*(i-1)+3,:) = [v1 v2 v3 v6];
1203 elseif(min(v2,v4)<min(v1,v5) & min(v1,v6)>min(v3,v4) & min(v3,v5)<min(v2,v6)),
1204 elements(3*(i-1)+1,:) = [v1 v2 v3 v4];
1205 elements(3*(i-1)+2,:) = [v2 v3 v4 v5];
1206 elements(3*(i-1)+3,:) = [v3 v4 v5 v6];
1207 elseif(min(v2,v4)<min(v1,v5) & min(v1,v6)>min(v3,v4) & min(v3,v5)>min(v2,v6)),
1208 elements(3*(i-1)+1,:) = [v1 v2 v3 v4];
1209 elements(3*(i-1)+2,:) = [v2 v4 v5 v6];
1210 elements(3*(i-1)+3,:) = [v2 v3 v4 v6];
1211 elseif(min(v2,v4)>min(v1,v5) & min(v1,v6)<min(v3,v4) & min(v3,v5)<min(v2,v6)),
1212 elements(3*(i-1)+1,:) = [v1 v4 v5 v6];
1213 elements(3*(i-1)+2,:) = [v1 v2 v3 v5];
1214 elements(3*(i-1)+3,:) = [v1 v3 v5 v6];
1215 elseif(min(v2,v4)>min(v1,v5) & min(v1,v6)<min(v3,v4) & min(v3,v5)>min(v2,v6)),
1216 elements(3*(i-1)+1,:) = [v1 v4 v5 v6];
1217 elements(3*(i-1)+2,:) = [v1 v2 v5 v6];
1218 elements(3*(i-1)+3,:) = [v1 v2 v3 v6];
1219 elseif(min(v2,v4)>min(v1,v5) & min(v1,v6)>min(v3,v4) & min(v3,v5)<min(v2,v6)),
1220 elements(3*(i-1)+1,:) = [v1 v3 v4 v5];
1221 elements(3*(i-1)+2,:) = [v1 v2 v3 v5];
1222 elements(3*(i-1)+3,:) = [v3 v4 v5 v6];
1223 elseif(min(v2,v4)>min(v1,v5) & min(v1,v6)<min(v3,v4) & min(v3,v5)<min(v2,v6)),
1224 elements(3*(i-1)+1,:) = [v1 v5 v6 v4];
1225 elements(3*(i-1)+2,:) = [v1 v2 v3 v5];
1226 elements(3*(i-1)+3,:) = [v5 v6 v3 v1];
1227 elseif(min(v2,v4)>min(v1,v5) & min(v1,v6)>min(v3,v4) & min(v3,v5)>min(v2,v6)),
1228 steiner = steiner+1; nbv = nbv+1; nbt = nbt+5; v7 = nbv;
1229 md.mesh.x=[md.mesh.x; mean(md.mesh.x(md.mesh.elements(i,:)))];
1230 md.mesh.y=[md.mesh.y; mean(md.mesh.y(md.mesh.elements(i,:)))];
1231 md.mesh.z=[md.mesh.z; mean(md.mesh.z(md.mesh.elements(i,:)))];
1232 elements(3*(i-1)+1,:) = [v1 v2 v3 v7];
1233 elements(3*(i-1)+2,:) = [v1 v4 v5 v7];
1234 elements(3*(i-1)+3,:) = [v1 v2 v5 v7];
1235 elements(end+1,:) = [v2 v5 v6 v7];
1236 elements(end+1,:) = [v2 v3 v6 v7];
1237 elements(end+1,:) = [v3 v4 v6 v7];
1238 elements(end+1,:) = [v1 v3 v4 v7];
1239 elements(end+1,:) = [v4 v5 v6 v7];
1240 else
1241 error('Case not supported'); %not supposed to happen!
1242 end
1243 %Reorder elements to make sure they are direct
1244 for j=1:3
1245 element = elements(3*(i-1)+j,:);
1246 matrix = [md.mesh.x(element), md.mesh.y(element), md.mesh.z(element), ones(4,1)];
1247 if det(matrix)>0,
1248 elements(3*(i-1)+j,1)=element(2);
1249 elements(3*(i-1)+j,2)=element(1);
1250 end
1251 end
1252 end
1253 %%Split in 3 tetras
1254 %subelement1 = [1 2 3 5];
1255 %subelement2 = [4 6 5 1];
1256 %subelement3 = [5 6 3 1];
1257 %elements=[md.mesh.elements(:,subelement1);md.mesh.elements(:,subelement2);md.mesh.elements(:,subelement3)];
1258 if steiner==0,
1259 disp('No Steiner point required to split prismatic mesh into tets');
1260 else
1261 disp([num2str(steiner) ' Steiner points had to be included'])
1262 error('Steiner point not supported yet');
1263 end
1264
1265 pos_elements = repmat([1:md.mesh.numberofelements]',3,1);
1266
1267 md.mesh.elements=elements;
1268 md.mesh.numberofelements=size(elements,1);
1269
1270 %p and q (same deal, except for element that are on the bedrock: )
1271 if ~isnan(md.friction.p),
1272 md.friction.p=md.friction.p(pos_elements);
1273 md.friction.q=md.friction.q(pos_elements);
1274 end
1275
1276 %elementstype
1277 if ~isnan(md.flowequation.element_equation)
1278 oldelements_type=md.flowequation.element_equation;
1279 md.flowequation.element_equation=md.flowequation.element_equation(pos_elements);
1280 end
1281
1282 %connectivity
1283 md.mesh.elementconnectivity=NaN;
1284
1285 %materials
1286 if ~isnan(md.materials.rheology_n),
1287 md.materials.rheology_n=md.materials.rheology_n(pos_elements);
1288 end
1289
1290 %increase connectivity if less than 25:
1291 if md.mesh.average_vertex_connectivity<=25,
1292 md.mesh.average_vertex_connectivity=100;
1293 end
1294 end % }}}
1295 function disp(self) % {{{
1296 disp(sprintf('%19s: %-22s -- %s','mesh' ,['[1x1 ' class(self.mesh) ']'],'mesh properties'));
1297 disp(sprintf('%19s: %-22s -- %s','mask' ,['[1x1 ' class(self.mask) ']'],'defines grounded and floating elements'));
1298 disp(sprintf('%19s: %-22s -- %s','geometry' ,['[1x1 ' class(self.geometry) ']'],'surface elevation, bedrock topography, ice thickness,...'));
1299 disp(sprintf('%19s: %-22s -- %s','constants' ,['[1x1 ' class(self.constants) ']'],'physical constants'));
1300 disp(sprintf('%19s: %-22s -- %s','smb' ,['[1x1 ' class(self.smb) ']'],'surface mass balance'));
1301 disp(sprintf('%19s: %-22s -- %s','basalforcings' ,['[1x1 ' class(self.basalforcings) ']'],'bed forcings'));
1302 disp(sprintf('%19s: %-22s -- %s','materials' ,['[1x1 ' class(self.materials) ']'],'material properties'));
1303 disp(sprintf('%19s: %-22s -- %s','damage' ,['[1x1 ' class(self.damage) ']'],'parameters for damage evolution solution'));
1304 disp(sprintf('%19s: %-22s -- %s','friction' ,['[1x1 ' class(self.friction) ']'],'basal friction/drag properties'));
1305 disp(sprintf('%19s: %-22s -- %s','flowequation' ,['[1x1 ' class(self.flowequation) ']'],'flow equations'));
1306 disp(sprintf('%19s: %-22s -- %s','timestepping' ,['[1x1 ' class(self.timestepping) ']'],'time stepping for transient models'));
1307 disp(sprintf('%19s: %-22s -- %s','initialization' ,['[1x1 ' class(self.initialization) ']'],'initial guess/state'));
1308 disp(sprintf('%19s: %-22s -- %s','rifts' ,['[1x1 ' class(self.rifts) ']'],'rifts properties'));
1309 disp(sprintf('%19s: %-22s -- %s','slr' ,['[1x1 ' class(self.slr) ']'],'slr forcings'));
1310 disp(sprintf('%19s: %-22s -- %s','debug' ,['[1x1 ' class(self.debug) ']'],'debugging tools (valgrind, gprof)'));
1311 disp(sprintf('%19s: %-22s -- %s','verbose' ,['[1x1 ' class(self.verbose) ']'],'verbosity level in solve'));
1312 disp(sprintf('%19s: %-22s -- %s','settings' ,['[1x1 ' class(self.settings) ']'],'settings properties'));
1313 disp(sprintf('%19s: %-22s -- %s','toolkits' ,['[1x1 ' class(self.toolkits) ']'],'PETSc options for each solution'));
1314 disp(sprintf('%19s: %-22s -- %s','cluster' ,['[1x1 ' class(self.cluster) ']'],'cluster parameters (number of cpus...)'));
1315 disp(sprintf('%19s: %-22s -- %s','balancethickness',['[1x1 ' class(self.balancethickness) ']'],'parameters for balancethickness solution'));
1316 disp(sprintf('%19s: %-22s -- %s','stressbalance' ,['[1x1 ' class(self.stressbalance) ']'],'parameters for stressbalance solution'));
1317 disp(sprintf('%19s: %-22s -- %s','groundingline' ,['[1x1 ' class(self.groundingline) ']'],'parameters for groundingline solution'));
1318 disp(sprintf('%19s: %-22s -- %s','hydrology' ,['[1x1 ' class(self.hydrology) ']'],'parameters for hydrology solution'));
1319 disp(sprintf('%19s: %-22s -- %s','masstransport' ,['[1x1 ' class(self.masstransport) ']'],'parameters for masstransport solution'));
1320 disp(sprintf('%19s: %-22s -- %s','thermal' ,['[1x1 ' class(self.thermal) ']'],'parameters for thermal solution'));
1321 disp(sprintf('%19s: %-22s -- %s','steadystate' ,['[1x1 ' class(self.steadystate) ']'],'parameters for steadystate solution'));
1322 disp(sprintf('%19s: %-22s -- %s','transient' ,['[1x1 ' class(self.transient) ']'],'parameters for transient solution'));
1323 disp(sprintf('%19s: %-22s -- %s','levelset' ,['[1x1 ' class(self.levelset) ']'],'parameters for moving boundaries (level-set method)'));
1324 disp(sprintf('%19s: %-22s -- %s','calving' ,['[1x1 ' class(self.calving) ']'],'parameters for calving'));
1325 disp(sprintf('%19s: %-22s -- %s','gia' ,['[1x1 ' class(self.gia) ']'],'parameters for gia solution'));
1326 disp(sprintf('%19s: %-22s -- %s','esa' ,['[1x1 ' class(self.esa) ']'],'parameters for elastic adjustment solution'));
1327 disp(sprintf('%19s: %-22s -- %s','autodiff' ,['[1x1 ' class(self.autodiff) ']'],'automatic differentiation parameters'));
1328 disp(sprintf('%19s: %-22s -- %s','inversion' ,['[1x1 ' class(self.inversion) ']'],'parameters for inverse methods'));
1329 disp(sprintf('%19s: %-22s -- %s','qmu' ,['[1x1 ' class(self.qmu) ']'],'dakota properties'));
1330 disp(sprintf('%19s: %-22s -- %s','amr' ,['[1x1 ' class(self.amr) ']'],'adaptive mesh refinement properties'));
1331 disp(sprintf('%19s: %-22s -- %s','outputdefinition',['[1x1 ' class(self.outputdefinition) ']'],'output definition'));
1332 disp(sprintf('%19s: %-22s -- %s','results' ,['[1x1 ' class(self.results) ']'],'model results'));
1333 disp(sprintf('%19s: %-22s -- %s','radaroverlay' ,['[1x1 ' class(self.radaroverlay) ']'],'radar image for plot overlay'));
1334 disp(sprintf('%19s: %-22s -- %s','miscellaneous' ,['[1x1 ' class(self.miscellaneous) ']'],'miscellaneous fields'));
1335 end % }}}
1336 function memory(self) % {{{
1337
1338 disp(sprintf('\nMemory imprint:\n'));
1339
1340 fields=properties('model');
1341 mem=0;
1342
1343 for i=1:length(fields),
1344 field=self.(fields{i});
1345 s=whos('field');
1346 mem=mem+s.bytes/1e6;
1347 disp(sprintf('%19s: %6.2f Mb',fields{i},s.bytes/1e6));
1348 end
1349 disp(sprintf('%19s--%10s','--------------','--------------'));
1350 disp(sprintf('%19s: %g Mb','Total',mem));
1351 end % }}}
1352 function netcdf(self,filename) % {{{
1353 %NETCDF - save model as netcdf
1354 %
1355 % Usage:
1356 % netcdf(md,filename)
1357 %
1358 % Example:
1359 % netcdf(md,'model.nc');
1360
1361 disp('Saving model as NetCDF');
1362 %1. Create NetCDF file
1363 ncid=netcdf.create(filename,'CLOBBER');
1364 netcdf.putAtt(ncid,netcdf.getConstant('NC_GLOBAL'),'Conventions','CF-1.4');
1365 netcdf.putAtt(ncid,netcdf.getConstant('NC_GLOBAL'),'Title',['ISSM model (' self.miscellaneous.name ')']);
1366 netcdf.putAtt(ncid,netcdf.getConstant('NC_GLOBAL'),'Author',getenv('USER'));
1367 netcdf.putAtt(ncid,netcdf.getConstant('NC_GLOBAL'),'Date',datestr(now));
1368
1369 %Preallocate variable id, needed to write variables in netcdf file
1370 var_id=zeros(1000,1);%preallocate
1371
1372 for step=1:2,
1373 counter=0;
1374 [var_id,counter]=structtonc(ncid,'md',self,0,var_id,counter,step);
1375 if step==1, netcdf.endDef(ncid); end
1376 end
1377
1378 if counter>1000,
1379 warning(['preallocation of var_id need to be updated from ' num2str(1000) ' to ' num2str(counter)]);
1380 end
1381
1382 netcdf.close(ncid)
1383 end % }}}
1384 function xylim(self) % {{{
1385
1386 xlim([min(self.mesh.x) max(self.mesh.x)]);
1387 ylim([min(self.mesh.y) max(self.mesh.y)])
1388 end % }}}
1389 function md=upload(md) % {{{
1390 %the goal of this routine is to upload the model onto a server, and to empty it.
1391 %So first, save the model with a unique name and upload the file to the server:
1392 random_part=fix(rand(1)*10000);
1393 id=[md.miscellaneous.name '-' regexprep(datestr(now),'[^\w'']','') '-' num2str(random_part) '-' getenv('USER') '-' oshostname() '.upload'];
1394 eval(['save ' id ' md']);
1395
1396 %Now, upload the file:
1397 issmscpout(md.settings.upload_server,md.settings.upload_path,md.settings.upload_login,md.settings.upload_port,{id},1);
1398
1399 %Now, empty this model of everything except settings, and record name of file we just uploaded!
1400 settings_back=md.settings;
1401 md=model();
1402 md.settings=settings_back;
1403 md.settings.upload_filename=id;
1404
1405 %get locally rid of file that was uploaded
1406 eval(['delete ' id]);
1407
1408 end % }}}
1409 function md=download(md) % {{{
1410
1411 %the goal of this routine is to download the internals of the current model from a server, because
1412 %this model is empty, except for the settings which tell us where to go and find this model!
1413
1414 %Download the file:
1415 issmscpin(md.settings.upload_server, md.settings.upload_login, md.settings.upload_port, md.settings.upload_path, {md.settings.upload_filename});
1416
1417 name=md.settings.upload_filename;
1418
1419 %Now, load this model:
1420 md=loadmodel(md.settings.upload_filename);
1421
1422 %get locally rid of file that was downloaded
1423 eval(['delete ' name]);
1424
1425 end % }}}
1426 function savemodeljs(md,modelname,websiteroot) % {{{
1427
1428 %the goal of this routine is to save the model as a javascript array that can be included in any html
1429 %file:
1430
1431 %disp:
1432 disp(['saving model ''' modelname ''' in file ' websiteroot '/js/' modelname '.js']);
1433
1434 %open file for writing and declare the model:
1435 fid=fopen([websiteroot '/js/' modelname '.js'],'w');
1436 fprintf(fid,'var %s=new model();\n',modelname);
1437
1438 %now go through all the classes and fwrite all the corresponding fields:
1439
1440 fields=properties('model');
1441 for i=1:length(fields),
1442 field=fields{i};
1443
1444 %Some properties do not need to be saved
1445 if ismember(field,{'results','cluster' }),
1446 continue;
1447 end
1448
1449 %Check that current field is an object
1450 if ~isobject(md.(field))
1451 error(['field ''' char(field) ''' is not an object']);
1452 end
1453
1454 %savemodeljs for current object
1455 %disp(['javascript saving ' field '...']);
1456 savemodeljs(md.(field),fid,modelname);
1457 end
1458
1459 %done, close file:
1460 fclose(fid);
1461 end
1462 end
1463 end
Note: See TracBrowser for help on using the repository browser.