Changeset 25836 for issm/trunk/src/m/classes/hydrologydc.m
- Timestamp:
- 12/08/20 08:45:53 (4 years ago)
- Location:
- issm/trunk
- Files:
-
- 3 edited
Legend:
- Unmodified
- Added
- Removed
-
issm/trunk
- Property svn:mergeinfo changed
-
issm/trunk/src
- Property svn:mergeinfo changed
-
issm/trunk/src/m/classes/hydrologydc.m
r24313 r25836 1 1 2 %Hydrologydc class definition 2 3 % … … 12 13 rel_tol = 0; 13 14 max_iter = 0; 14 steps_per_step = 0; 15 steps_per_step = 0; 16 averaging = 0; 15 17 sedimentlimit_flag = 0; 16 18 sedimentlimit = 0; … … 19 21 leakage_factor = 0; 20 22 basal_moulin_input = NaN; 21 23 requested_outputs = {}; 22 24 23 25 spcsediment_head = NaN; 24 26 mask_thawed_node = NaN; 25 27 sediment_transmitivity = NaN; 26 28 sediment_compressibility = 0; … … 39 41 epl_conductivity = 0; 40 42 eplflip_lock = 0; 41 43 end 42 44 methods 43 % {{{ function self = extrude(self,md) 44 function self = extrude(self,md) 45 function self = extrude(self,md) % {{{ 45 46 self.spcsediment_head=project3d(md,'vector',self.spcsediment_head,'type','node','layer',1); 46 47 self.sediment_transmitivity=project3d(md,'vector',self.sediment_transmitivity,'type','node','layer',1); 47 48 self.basal_moulin_input=project3d(md,'vector',self.basal_moulin_input,'type','node','layer',1); 48 49 self.mask_thawed_node=project3d(md,'vector',self.mask_thawed_node,'type','node','layer',1); 49 50 if(self.isefficientlayer==1); 50 51 self.spcepl_head=project3d(md,'vector',self.spcepl_head,'type','node','layer',1); 51 self.mask_eplactive_node=project3d(md,'vector',self.mask_eplactive_node,'type','node','layer',1); 52 end 53 end 54 % }}} 55 % {{{ function self = hydrologydc(varargin) 56 function self = hydrologydc(varargin) 52 self.mask_eplactive_node=project3d(md,'vector',self.mask_eplactive_node,'type','node','layer',1); 53 end 54 end % }}} 55 function self = hydrologydc(varargin) % {{{ 57 56 switch nargin 58 57 case 0 … … 61 60 error('constructor not supported'); 62 61 end 63 end 64 % }}} 65 function list = defaultoutputs(self,md) % {{{ 66 list = {'SedimentHead','SedimentHeadResidual','EffectivePressure'}; 67 if self.isefficientlayer, 68 list=[list,{'EplHead','HydrologydcMaskEplactiveNode','HydrologydcMaskEplactiveElt','EplHeadSlopeX','EplHeadSlopeY','HydrologydcEplThickness'}]; 69 end 70 if self.steps_per_step>1, 71 list = [list,'EffectivePressureSubstep','SedimentHeadSubstep']; 72 if self.isefficientlayer, 73 list = [list,'EplHeadSubstep','HydrologydcEplThicknessSubstep']; 74 end 75 end 76 end % }}} 77 62 end% }}} 63 function list = defaultoutputs(self,md) % {{{ 64 list = {'SedimentHead','SedimentHeadResidual','EffectivePressure'}; 65 if self.isefficientlayer, 66 list=[list,{'EplHead','HydrologydcMaskEplactiveNode','HydrologydcMaskEplactiveElt','EplHeadSlopeX','EplHeadSlopeY','HydrologydcEplThickness'}]; 67 end 68 if self.steps_per_step>1, 69 list = [list,'EffectivePressureSubstep','SedimentHeadSubstep']; 70 if self.isefficientlayer, 71 list = [list,'EplHeadSubstep','HydrologydcEplThicknessSubstep']; 72 end 73 end 74 end % }}} 78 75 function self = initialize(self,md) % {{{ 79 76 self.epl_colapse_thickness = self.sediment_transmitivity/self.epl_conductivity; … … 84 81 85 82 end % }}} 86 % {{{ function self = setdefaultparameters(self) 87 88 function self = setdefaultparameters(self) 89 83 function self = setdefaultparameters(self)% {{{ 90 84 %Parameters from de Fleurian 2014 91 85 self.water_compressibility = 5.04e-10; … … 96 90 self.max_iter = 100; 97 91 self.steps_per_step = 1; 92 self.averaging = 0; 98 93 self.sedimentlimit_flag = 0; 99 94 self.sedimentlimit = 0; … … 101 96 self.unconfined_flag = 0; 102 97 self.leakage_factor = 1.0e-10; 103 98 self.requested_outputs = {'default'}; 104 99 105 100 … … 117 112 self.epl_max_thickness = 5.0; 118 113 self.eplflip_lock = 0; 119 end 120 121 % }}} 122 % {{{ function md = checkconsistency(self,md,solution,analyses) 123 124 function md = checkconsistency(self,md,solution,analyses) 125 %Early return 114 end % }}} 115 function md = checkconsistency(self,md,solution,analyses)% {{{ 116 %Early return 126 117 if ~ismember('HydrologyDCInefficientAnalysis',analyses) & ~ismember('HydrologyDCEfficientAnalysis',analyses), 127 118 return; … … 135 126 md = checkfield(md,'fieldname','hydrology.max_iter','>',0,'numel',1); 136 127 md = checkfield(md,'fieldname','hydrology.steps_per_step','>',0,'numel',1); 128 md = checkfield(md,'fieldname','hydrology.averaging','numel',[1],'values',[0 1 2]); 137 129 md = checkfield(md,'fieldname','hydrology.sedimentlimit_flag','numel',[1],'values',[0 1 2 3]); 138 130 md = checkfield(md,'fieldname','hydrology.transfer_flag','numel',[1],'values',[0 1]); 139 131 md = checkfield(md,'fieldname','hydrology.unconfined_flag','numel',[1],'values',[0 1]); 140 132 md = checkfield(md,'fieldname','hydrology.requested_outputs','stringrow',1); 141 133 if self.sedimentlimit_flag==1, 142 134 md = checkfield(md,'fieldname','hydrology.sedimentlimit','>',0,'numel',1); … … 152 144 md = checkfield(md,'fieldname','hydrology.sediment_thickness','>',0,'numel',1); 153 145 md = checkfield(md,'fieldname','hydrology.sediment_transmitivity','>=',0,'size',[md.mesh.numberofvertices 1]); 154 146 md = checkfield(md,'fieldname','hydrology.mask_thawed_node','size',[md.mesh.numberofvertices 1],'values',[0 1]); 155 147 156 148 if self.isefficientlayer==1, … … 169 161 end 170 162 end 171 end 172 173 % }}} 174 % {{{ function disp(self) 175 176 function disp(self) 163 end% }}} 164 function disp(self)% {{{ 177 165 disp(sprintf(' hydrology Dual Porous Continuum Equivalent parameters:')); 178 166 disp(sprintf(' - general parameters')); … … 184 172 fielddisplay(self,'max_iter','maximum number of nonlinear iteration'); 185 173 fielddisplay(self,'steps_per_step','number of hydrology steps per timestep'); 174 fielddisplay(self, 'averaging', 'averaging methods from short to long steps'); 175 disp(sprintf('%55s 0: Arithmetic (default)')); 176 disp(sprintf('%55s 0: Geometric')); 177 disp(sprintf('%55s 0: Harmonic')); 186 178 fielddisplay(self,'sedimentlimit_flag','what kind of upper limit is applied for the inefficient layer'); 187 179 disp(sprintf('%55s 0: no limit',' ')); … … 224 216 fielddisplay(self,'eplflip_lock','lock the epl activation to avoid fli-floping (default is 0, no stabilization)'); 225 217 end 226 227 end 228 229 % }}} 230 % {{{ function marshall(self,prefix,md,fid) 231 232 function marshall(self,prefix,md,fid) 218 end % }}} 219 function marshall(self,prefix,md,fid)% {{{ 233 220 WriteData(fid,prefix,'name','md.hydrology.model','data',1,'format','Integer'); 234 221 WriteData(fid,prefix,'object',self,'fieldname','water_compressibility','format','Double'); … … 239 226 WriteData(fid,prefix,'object',self,'fieldname','max_iter','format','Integer'); 240 227 WriteData(fid,prefix,'object',self,'fieldname','steps_per_step','format','Integer'); 228 WriteData(fid,prefix,'object',self,'fieldname','averaging','format','Integer'); 241 229 WriteData(fid,prefix,'object',self,'fieldname','sedimentlimit_flag','format','Integer'); 242 230 WriteData(fid,prefix,'object',self,'fieldname','transfer_flag','format','Integer'); … … 255 243 WriteData(fid,prefix,'object',self,'fieldname','sediment_thickness','format','Double'); 256 244 WriteData(fid,prefix,'object',self,'fieldname','sediment_transmitivity','format','DoubleMat','mattype',1'); 257 WriteData(fid,prefix,'object',self,'fieldname','mask_thawed_node','format','DoubleMat','mattype',1); 258 245 WriteData(fid,prefix,'object',self,'fieldname','mask_thawed_node','format','DoubleMat','mattype',1); 259 246 if self.isefficientlayer==1, 260 247 WriteData(fid,prefix,'object',self,'fieldname','spcepl_head','format','DoubleMat','mattype',1,'timeserieslength',md.mesh.numberofvertices+1,'yts',md.constants.yts); … … 269 256 WriteData(fid,prefix,'object',self,'fieldname','eplflip_lock','format','Integer'); 270 257 end 271 outputs = self.requested_outputs; 272 pos = find(ismember(outputs,'default')); 273 if ~isempty(pos), 274 outputs(pos) = []; %remove 'default' from outputs 275 outputs = [outputs defaultoutputs(self,md)]; %add defaults 276 end 277 WriteData(fid,prefix,'data',outputs,'name','md.hydrology.requested_outputs','format','StringArray'); 278 end 279 280 % }}} 258 outputs = self.requested_outputs; 259 pos = find(ismember(outputs,'default')); 260 if ~isempty(pos), 261 outputs(pos) = []; %remove 'default' from outputs 262 outputs = [outputs defaultoutputs(self,md)]; %add defaults 263 end 264 WriteData(fid,prefix,'data',outputs,'name','md.hydrology.requested_outputs','format','StringArray'); 265 end% }}} 281 266 end 282 267 end
Note:
See TracChangeset
for help on using the changeset viewer.