KerMor  0.9
Model order reduction for nonlinear dynamical systems and nonlinear approximation
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
AMuscleConfig.m
Go to the documentation of this file.
1 namespace models{
2 namespace muscle{
3 
4 
5 /* (Autoinserted by mtoc++)
6  * This source code has been filtered by the mtoc++ executable,
7  * which generates code that can be processed by the doxygen documentation tool.
8  *
9  * On the other hand, it can neither be interpreted by MATLAB, nor can it be compiled with a C++ compiler.
10  * Except for the comments, the function bodies of your M-file functions are untouched.
11  * Consequently, the FILTER_SOURCE_FILES doxygen switch (default in our Doxyfile.template) will produce
12  * attached source files that are highly readable by humans.
13  *
14  * Additionally, links in the doxygen generated documentation to the source code of functions and class members refer to
15  * the correct locations in the source code browser.
16  * However, the line numbers most likely do not correspond to the line numbers in the original MATLAB source files.
17  */
18 
20  :public fem.AFEMConfig {
27  public:
28 
30 
31 
32  public:
33 
35 
36 
38 
39 
61  public:
62 
93  private:
94 
95  iP;
96 
97  optArgs = {""};
98 
99 
100  public:
101 
103  this = this@fem.AFEMConfig(varargin[:]);
104  /* Force-length function type */
105  this.addOption(" FL ",1);
106  }
107 
108 
109  function m = createModel() {
110  m = models.muscle.Model(this);
111  }
119  function configureModelFinal() {
120  configureModelFinal@fem.AFEMConfig(this);
121  /* Automatically set the reducable dimensions for u,v
122  * separately. */
123  m = this.Model;
124  s = m.System;
125  /* targetd = s.NumStateDofs + (1:s.NumDerivativeDofs); */
126  targetd = 1:s.NumStateDofs;
127  /* Dont reduce those dofs subject to velocity BCs! */
128  targetd(s.DerivativeDirichletPosInStateDofs) = [];
129  m.SpaceReducer.TargetDimensions= targetd;
130  }
131 
132 
133  function u = getInputs() {
134  u = [];
135  }
155  function x0 = getX0(x0) {
156  }
165  function setForceLengthFun(f) {
166 
167  /* Steps to produce derivative of exponential version
168  * lw = sym('lw'); rw = sym('rw'); lexpo = sym('lexpo'); rexpo = sym('rexpo'); ratio = sym('ratio');
169  * f(ratio) = exp(-((1-ratio)/lw).^lexpo);
170  * df = diff(f)
171  * f2(ratio) = exp(-((ratio-1)/rw).^rexpo);
172  * df2 = diff(f2) */
173 
174  switch this.Options.FL
175  case 1
176  /* Linear (Gordon 66)
177  * Exported here to separate class for tidyness */
178  g = models.muscle.functions.Gordon66SarcoForceLength(f.mu(14));
179  [fun, dfun] = g.getFunction;
180  case 2
181  /* Exponential (Schmitt) */
182  lexpo = 4; /* 4 */
183 
184  rexpo = 3; /* 3 */
185 
186  lw = f.mu(14); /* orig .57 */
187 
188  rw = f.mu(14)*1.3; /* orig .14 */
189 
190  fun = @(ratio)(ratio<=1).*exp(-((1-ratio)/lw).^lexpo) ...
191  + (ratio>1).*exp(-((ratio-1)/rw).^rexpo);
192  dfun = @(ratio)(ratio<=1).*((lexpo*exp(-(-(ratio - 1)/lw)^lexpo)*(-(ratio - 1)/lw)^(lexpo - 1))/lw) ...
193  + (ratio > 1) .* (-(rexpo*exp(-((ratio - 1)/rw)^rexpo)*((ratio - 1)/rw)^(rexpo - 1))/rw);
194  case 3
195  /* Quadratic Polynomial (Heidlauf) */
196  fun = @(ratio)(-6.25*ratio.*ratio + 12.5*ratio - 5.25) .* (ratio >= .6) .* (ratio <= 1.4);
197  dfun = @(ratio)(12.5*ratio.*(1-ratio)) .* (ratio >= .6) .* (ratio <= 1.4);
198  end
199  f.ForceLengthFun= fun;
200  f.ForceLengthFunDeriv= dfun;
201  }
224  function alpha = getAlphaRamp(ramptime,double alphamax,double starttime) {
225  if nargin < 4
226  starttime = this.ActivationRampOffset;
227  if nargin < 3
228  alphamax = this.ActivationRampMax;
229  end
230  end
231  ramp = general.functions.Ramp(ramptime, alphamax, starttime);
232  alpha = ramp.getFunction;
233  }
252  function tmr = getTendonMuscleRatio(unused1) {
253  tmr = []; /* zeros(1,size(x,2)); */
254 
255  }
273  protected:
274 
275 
276  function init() {
277  init@fem.AFEMConfig(this);
278 
279  /* % Get the geometry */
280  geo = this.Geometry;
281  if isa(geo," fem.geometry.Cube8Node ")
282  press_geo = geo;
283  elseif isa(geo," fem.geometry.Cube20Node ") || isa(geo," fem.geometry.Cube27Node ")
284  press_geo = geo.toCube8Node;
285  else
286  error(" Scenario not yet implemented for geometry class '%s' ", class(geo));
287  end
288  this.PressureFEM= fem.HexahedronTrilinear(press_geo);
289  /* this.PressureFEM = fem.HexahedronSerendipity(press_geo.toCube20Node);
290  *this.PressureFEM = fem.HexahedronTriquadratic(press_geo.toCube27Node); */
291  }
292 
293 
294  function anull = seta0(anull) {
295  }
304  function ftw = getFibreTypeWeights() {
305  fe = this.FEM;
306  geo = fe.Geometry;
307  ftw = zeros(fe.GaussPointsPerElem,length(this.FibreTypes),geo.NumElements);
308  }
319  public: /* ( Sealed ) */
320 
321  function anull = geta0() {
322  fe = this.FEM;
323  g = this.Geometry;
324  anull = zeros(3,fe.GaussPointsPerElem,g.NumElements);
325  anull = this.seta0(anull);
326  /* Normalize anull vectors */
327  for m = 1:g.NumElements
328  anull(:,:,m) = anull(:,:,m) ./ ([1;1;1]*Norm.L2(anull(:,:,m)));
329  end
330  }
331 
332 
333 
334 };
335 }
336 }
337 
338 
339 
function setForceLengthFun(f)
Provided here only for convenient outside access.
function ftw = getFibreTypeWeights()
This is a lazy pre-implementation as fullmodels.muscle.Models always have fibre types and thus weight...
function addOption(name, default, varargin)
Definition: AFEMConfig.m:219
function alpha = getAlphaRamp(ramptime,double alphamax,double starttime)
Creates a linearly increasing scalar function starting at starttime milliseconds ranging from zero to...
A double value.
function x0 = getX0(x0)
do nothing
function anull = geta0()
function m = createModel()
Convenience method.
A variable number of input arguments.
function anull = seta0(anull)
do nothing!
AModelConfig.
Definition: AFEMConfig.m:18
double ActivationRampOffset
Determines the default number of milliseconds to wait before activation is started.
Definition: AMuscleConfig.m:78
double ActivationRampMax
Determines the default value for maximum activation in activation ramps.
Definition: AMuscleConfig.m:63
static function rowvec< double > n = L2(matrix< double > x)
Returns the discrete norm for each column vector in x.
Definition: Norm.m:39
char a0CoordinateSystem
The coordinate system in which to interpret the a0 vectors of fibre directions.
Definition: AMuscleConfig.m:40
Norm: Static class for commonly used norms on sets of vectors.
Definition: Norm.m:17
A MatLab character array.
function u = getInputs()
Returns the inputs of the model.
function tmr = getTendonMuscleRatio(unused1)
Returns the [0,1] ratio between tendon and muscle at all gauss points of all elements.