1 /* -------------------------------------------------------------------------- *
2  *                       OpenSim:  Thelen2003Muscle.cpp                       *
3  * -------------------------------------------------------------------------- *
4  * The OpenSim API is a toolkit for musculoskeletal modeling and simulation.  *
5  * See http://opensim.stanford.edu and the NOTICE file for more information.  *
6  * OpenSim is developed at Stanford University and supported by the US        *
7  * National Institutes of Health (U54 GM072970, R24 HD065690) and by DARPA    *
8  * through the Warrior Web program.                                           *
9  *                                                                            *
10  * Copyright (c) 2005-2017 Stanford University and the Authors                *
11  * Author(s): Matthew Millard                                                 *
12  *                                                                            *
13  * Licensed under the Apache License, Version 2.0 (the "License"); you may    *
14  * not use this file except in compliance with the License. You may obtain a  *
15  * copy of the License at http://www.apache.org/licenses/LICENSE-2.0.         *
16  *                                                                            *
17  * Unless required by applicable law or agreed to in writing, software        *
18  * distributed under the License is distributed on an "AS IS" BASIS,          *
19  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.   *
20  * See the License for the specific language governing permissions and        *
21  * limitations under the License.                                             *
22  * -------------------------------------------------------------------------- */
23 //=============================================================================
24 // INCLUDES
25 //=============================================================================
26 #include <fstream>
27 #include <OpenSim/Simulation/Model/Model.h>
28 #include "Thelen2003Muscle.h"
29 
30 //=============================================================================
31 // STATICS
32 //=============================================================================
33 using namespace std;
34 using namespace OpenSim;
35 using namespace SimTK;
36 
37 //=============================================================================
38 // CONSTRUCTORS
39 //=============================================================================
40 // Uses default (compiler-generated) destructor, copy constructor, copy
41 // assignment operator.
42 
43 
44 //_____________________________________________________________________________
45 /**
46  * Default constructor.
47  */
Thelen2003Muscle()48 Thelen2003Muscle::Thelen2003Muscle()
49 {
50     setNull();
51     constructProperties();
52 }
53 
54 //_____________________________________________________________________________
55 /**
56  * Constructor.
57  */
58 Thelen2003Muscle::
Thelen2003Muscle(const std::string & aName,double aMaxIsometricForce,double aOptimalFiberLength,double aTendonSlackLength,double aPennationAngle)59 Thelen2003Muscle(const std::string& aName,  double aMaxIsometricForce,
60                   double aOptimalFiberLength,double aTendonSlackLength,
61                   double aPennationAngle)
62 {
63     setNull();
64     constructProperties();
65     setName(aName);
66 
67     setMaxIsometricForce(aMaxIsometricForce);
68     setOptimalFiberLength(aOptimalFiberLength);
69     setTendonSlackLength(aTendonSlackLength);
70     setPennationAngleAtOptimalFiberLength(aPennationAngle);
71 }
72 
73 //==============================================================================
74 // Component Interface
75 //==============================================================================
extendFinalizeFromProperties()76 void Thelen2003Muscle::extendFinalizeFromProperties()
77 {
78     Super::extendFinalizeFromProperties();
79 
80     SimTK_ERRCHK1_ALWAYS(get_FmaxTendonStrain() > 0,
81         "Thelen2003Muscle::extendFinalizeFromProperties",
82         "%s: FmaxTendonStrain must be greater than zero", getName().c_str());
83     SimTK_ERRCHK1_ALWAYS(get_FmaxMuscleStrain() > 0,
84         "Thelen2003Muscle::extendFinalizeFromProperties",
85         "%s: FmaxMuscleStrain must be greater than zero", getName().c_str());
86     SimTK_ERRCHK1_ALWAYS(get_KshapeActive() > 0,
87         "Thelen2003Muscle::extendFinalizeFromProperties",
88         "%s: KshapeActive must be greater than zero", getName().c_str());
89     SimTK_ERRCHK1_ALWAYS(get_KshapePassive() > 0,
90         "Thelen2003Muscle::extendFinalizeFromProperties",
91         "%s: KshapePassive must be greater than zero", getName().c_str());
92     SimTK_ERRCHK1_ALWAYS(get_Af() > 0,
93         "Thelen2003Muscle::extendFinalizeFromProperties",
94         "%s: Af must be greater than zero", getName().c_str());
95     SimTK_ERRCHK1_ALWAYS(get_Flen() > 1.0,
96         "Thelen2003Muscle::extendFinalizeFromProperties",
97         "%s: Flen must be greater than 1.0", getName().c_str());
98     SimTK_ERRCHK1_ALWAYS(get_fv_linear_extrap_threshold() > 1.0/get_Flen(),
99         "Thelen2003::extendFinalizeFromProperties",
100         "%s: F-v extrapolation threshold must be greater than 1.0/Flen",
101         getName().c_str());
102 
103 
104     OPENSIM_THROW_IF_FRMOBJ(get_minimum_activation() < 0.01,
105         InvalidPropertyValue, getProperty_minimum_activation().getName(),
106         "Minimum activation cannot be less than 0.01");
107 
108     OPENSIM_THROW_IF_FRMOBJ(getMinControl() < get_minimum_activation(),
109         InvalidPropertyValue, getProperty_min_control().getName(),
110         "Minimum control cannot be less than minimum activation");
111 
112     // Propagate properties down to pennation model subcomponent. If any of the
113     // new property values are invalid, restore the subcomponent's current
114     // property values (to avoid throwing again when the subcomponent's
115     // extendFinalizeFromProperties() method is called directly) and then
116     // re-throw the exception thrown by the subcomponent.
117     auto& pennMdl =
118         updMemberSubcomponent<MuscleFixedWidthPennationModel>(pennMdlIdx);
119     MuscleFixedWidthPennationModel pennMdlCopy(pennMdl);
120     pennMdl.set_optimal_fiber_length(getOptimalFiberLength());
121     pennMdl.set_pennation_angle_at_optimal(
122         getPennationAngleAtOptimalFiberLength());
123     pennMdl.set_maximum_pennation_angle(get_maximum_pennation_angle());
124     try {
125         pennMdl.finalizeFromProperties();
126     } catch (const InvalidPropertyValue&) {
127         pennMdl = pennMdlCopy;
128         throw;
129     }
130 
131     // Propagate properties down to activation dynamics model subcomponent.
132     // Handle invalid properties as above for pennation model.
133     auto& actMdl =
134         updMemberSubcomponent<MuscleFirstOrderActivationDynamicModel>(actMdlIdx);
135     MuscleFirstOrderActivationDynamicModel actMdlCopy(actMdl);
136     actMdl.set_activation_time_constant(get_activation_time_constant());
137     actMdl.set_deactivation_time_constant(get_deactivation_time_constant());
138     actMdl.set_minimum_activation(get_minimum_activation());
139     try {
140         actMdl.finalizeFromProperties();
141     } catch (const InvalidPropertyValue&) {
142         actMdl = actMdlCopy;
143         throw;
144     }
145 }
146 
147 //====================================================================
148 // Model Component Interface
149 //====================================================================
extendConnectToModel(Model & aModel)150 void Thelen2003Muscle::extendConnectToModel(Model& aModel)
151 {
152     Super::extendConnectToModel(aModel);
153 }
154 
extendInitStateFromProperties(SimTK::State & s) const155 void Thelen2003Muscle::extendInitStateFromProperties(SimTK::State& s) const
156 {
157     Super::extendInitStateFromProperties(s);
158 }
159 
160 void Thelen2003Muscle::
extendSetPropertiesFromState(const SimTK::State & s)161     extendSetPropertiesFromState(const SimTK::State& s)
162 {
163     Super::extendSetPropertiesFromState(s);
164 }
165 
166 //_____________________________________________________________________________
167 // Set the data members of this muscle to their null values.
setNull()168 void Thelen2003Muscle::setNull()
169 {
170     // no data members
171     setAuthors("Matthew Millard");
172 }
173 
174 //_____________________________________________________________________________
175 /**
176  * Populate this object's properties
177  */
constructProperties()178 void Thelen2003Muscle::constructProperties()
179 {
180     constructProperty_FmaxTendonStrain(0.04); // was 0.033
181     constructProperty_FmaxMuscleStrain(0.6);
182     constructProperty_KshapeActive(0.45);
183     constructProperty_KshapePassive(5.0);
184     constructProperty_Af(0.25);
185     constructProperty_Flen(1.4);    //was 1.8,
186     constructProperty_fv_linear_extrap_threshold(0.95);
187     //acos(0.1) = 84.26 degrees
188     constructProperty_maximum_pennation_angle(acos(0.1));
189     constructProperty_activation_time_constant(0.015);
190     constructProperty_deactivation_time_constant(0.050);
191     constructProperty_minimum_activation(0.01);
192 
193     setMinControl(get_minimum_activation());
194 }
195 
196 //=============================================================================
197 // GET
198 //=============================================================================
getActivationTimeConstant() const199 double Thelen2003Muscle::getActivationTimeConstant() const
200 {   return get_activation_time_constant(); }
201 
getDeactivationTimeConstant() const202 double Thelen2003Muscle::getDeactivationTimeConstant() const
203 {   return get_deactivation_time_constant(); }
204 
getMinimumActivation() const205 double Thelen2003Muscle::getMinimumActivation() const
206 {   return get_minimum_activation(); }
207 
208 const MuscleFirstOrderActivationDynamicModel& Thelen2003Muscle::
getActivationModel() const209 getActivationModel() const
210 {
211     return getMemberSubcomponent<MuscleFirstOrderActivationDynamicModel>(
212            actMdlIdx);
213 }
214 
215 const MuscleFixedWidthPennationModel& Thelen2003Muscle::
getPennationModel() const216 getPennationModel() const
217 {   return getMemberSubcomponent<MuscleFixedWidthPennationModel>(pennMdlIdx); }
218 
getMaximumPennationAngle() const219 double Thelen2003Muscle::getMaximumPennationAngle() const
220 {   return get_maximum_pennation_angle(); }
221 
222 //=============================================================================
223 // SET
224 //=============================================================================
setActivationTimeConstant(double actTimeConstant)225 void Thelen2003Muscle::setActivationTimeConstant(double actTimeConstant)
226 {   set_activation_time_constant(actTimeConstant); }
227 
setDeactivationTimeConstant(double deactTimeConstant)228 void Thelen2003Muscle::setDeactivationTimeConstant(double deactTimeConstant)
229 {   set_deactivation_time_constant(deactTimeConstant); }
230 
setMinimumActivation(double minimumActivation)231 void Thelen2003Muscle::setMinimumActivation(double minimumActivation)
232 {   set_minimum_activation(minimumActivation); }
233 
setMaximumPennationAngle(double maximumPennationAngle)234 void Thelen2003Muscle::setMaximumPennationAngle(double maximumPennationAngle)
235 {   set_maximum_pennation_angle(maximumPennationAngle); }
236 
237 //==============================================================================
238 //                             START OF DEPRECATED
239 //==============================================================================
240 double Thelen2003Muscle::
calcInextensibleTendonActiveFiberForce(SimTK::State & s,double aActivation) const241 calcInextensibleTendonActiveFiberForce(SimTK::State& s,
242                                        double aActivation) const
243 {
244     double inextensibleTendonActiveFiberForce = 0;
245 
246     double muscleLength = getLength(s);
247     double muscleVelocity = getLengtheningSpeed(s);
248     double tendonSlackLength = getTendonSlackLength();
249     double tendonVelocity = 0.0; //Inextensible tendon;
250 
251     double fiberLength  = getPennationModel().calcFiberLength(muscleLength,
252                                             tendonSlackLength);
253 
254     if(fiberLength > getPennationModel().getMinimumFiberLength()) {
255         double phi = getPennationModel().calcPennationAngle(fiberLength);
256 
257         double fiberVelocity = getPennationModel().calcFiberVelocity(
258                                     cos(phi),muscleVelocity,tendonVelocity);
259 
260         inextensibleTendonActiveFiberForce =
261             calcActiveFiberForceAlongTendon(    aActivation,
262                                                 fiberLength,
263                                                 fiberVelocity);
264     }
265 
266     return inextensibleTendonActiveFiberForce;
267 }
268 
269 double Thelen2003Muscle::
calcActiveFiberForceAlongTendon(double activation,double fiberLength,double fiberVelocity) const270             calcActiveFiberForceAlongTendon(double activation,
271                                             double fiberLength,
272                                             double fiberVelocity) const
273 {
274     double activeFiberForce = 0;
275     double clampedFiberLength = getPennationModel()
276                                 .clampFiberLength(fiberLength);
277 
278     //If the fiber is in a legal range, compute the force its generating
279     if(fiberLength > getPennationModel().getMinimumFiberLength())
280     {
281         //Clamp activation to a legal range
282         double clampedActivation = getActivationModel()
283                                    .clampActivation(activation);
284 
285         //Normalize fiber length and velocity
286         double normFiberLength    = clampedFiberLength/getOptimalFiberLength();
287         double normFiberVelocity  = fiberVelocity /
288                         (getOptimalFiberLength() * getMaxContractionVelocity());
289 
290 
291         //Evaluate the force active length and force velocity multipliers
292         double fal  = calcfal(normFiberLength);
293 
294         double fv = SimTK::NaN;
295         double fvTol = 1e-6;
296         int fvMaxIter= 100;
297         fv = calcfvInv( clampedActivation,
298                                             fal,
299                                             normFiberVelocity,
300                                             fvTol,
301                                             fvMaxIter);
302 
303         double fiso = getMaxIsometricForce();
304 
305         //Evaluate the pennation angle
306         double phi = getPennationModel().calcPennationAngle(fiberLength);
307 
308         //Compute the active fiber force
309         activeFiberForce = fiso * clampedActivation * fal * fv * cos(phi);
310     }
311     //Compute the active fiber force
312 
313 
314     return activeFiberForce;
315 }
316 
317 //==============================================================================
318 //                              END OF DEPRECATED
319 //==============================================================================
320 
321 //==============================================================================
322 // Muscle.h Interface
323 //==============================================================================
324 
computeActuation(const SimTK::State & s) const325 double  Thelen2003Muscle::computeActuation(const SimTK::State& s) const
326 {
327     const MuscleDynamicsInfo& mdi = getMuscleDynamicsInfo(s);
328     setActuation(s, mdi.tendonForce);
329     return( mdi.tendonForce );
330 }
331 
computeInitialFiberEquilibrium(SimTK::State & s) const332 void Thelen2003Muscle::computeInitialFiberEquilibrium(SimTK::State& s) const
333 {
334     //Initial activation and fiber length from input State, s.
335     _model->getMultibodySystem().realize(s, SimTK::Stage::Velocity);
336     double activation = getActivation(s);
337 
338     //Tolerance, in Newtons, of the desired equilibrium
339     const double tol = max( 1e-8*getMaxIsometricForce(),
340                             SimTK::SignificantReal * 10 );
341 
342     int maxIter = 20;  //Should this be user settable?
343 
344     std::pair<StatusFromInitMuscleState, ValuesFromInitMuscleState> result;
345 
346     try {
347          result = initMuscleState(s, activation, tol, maxIter);
348     }
349     catch (const std::exception& x) {
350         OPENSIM_THROW_FRMOBJ(MuscleCannotEquilibrate, x.what());
351     }
352 
353     switch(result.first) {
354 
355     case StatusFromInitMuscleState::Success_Converged:
356         setActuation(s, result.second["tendon_force"]);
357         setFiberLength(s, result.second["fiber_length"]);
358         break;
359 
360     case StatusFromInitMuscleState::Warning_FiberAtLowerBound:
361         printf("\n\nThelen2003Muscle initialization:"
362                " %s is at its minimum fiber length of %f\n",
363                getName().c_str(), result.second["fiber_length"]);
364         setActuation(s, result.second["tendon_force"]);
365         setFiberLength(s, result.second["fiber_length"]);
366         break;
367 
368     case StatusFromInitMuscleState::Failure_MaxIterationsReached:
369         // Report internal variables and throw exception.
370         std::ostringstream ss;
371         ss << "\n  Solution error " << abs(result.second["solution_error"])
372            << " exceeds tolerance of " << tol << "\n"
373            << "  Newton iterations reached limit of " << maxIter << "\n"
374            << "  Activation is " << activation << "\n"
375            << "  Fiber length is " << result.second["fiber_length"] << "\n";
376         OPENSIM_THROW_FRMOBJ(MuscleCannotEquilibrate, ss.str());
377         break;
378     }
379 }
380 
calcMuscleLengthInfo(const SimTK::State & s,MuscleLengthInfo & mli) const381 void Thelen2003Muscle::calcMuscleLengthInfo(const SimTK::State& s,
382                                             MuscleLengthInfo& mli) const
383 {
384     try{
385         double optFiberLength   = getOptimalFiberLength();
386         double mclLength        = getLength(s);
387         double tendonSlackLen   = getTendonSlackLength();
388 
389         //Clamp the minimum fiber length to its minimum physical value.
390         mli.fiberLength  = getPennationModel().clampFiberLength(
391                                 getStateVariableValue(s, STATE_FIBER_LENGTH_NAME));
392 
393         mli.normFiberLength = mli.fiberLength/optFiberLength;
394         mli.pennationAngle  = getPennationModel()
395                               .calcPennationAngle(mli.fiberLength);
396 
397         mli.cosPennationAngle = cos(mli.pennationAngle);
398         mli.sinPennationAngle = sin(mli.pennationAngle);
399 
400         mli.fiberLengthAlongTendon = mli.fiberLength*mli.cosPennationAngle;
401 
402         mli.tendonLength      = getPennationModel().calcTendonLength(
403                                     mli.cosPennationAngle,
404                                     mli.fiberLength,mclLength );
405         mli.normTendonLength  = mli.tendonLength / tendonSlackLen;
406         mli.tendonStrain      = mli.normTendonLength -  1.0;
407 
408 
409 
410         mli.fiberPassiveForceLengthMultiplier= calcfpe(mli.normFiberLength);
411         mli.fiberActiveForceLengthMultiplier = calcfal(mli.normFiberLength);
412     }catch(const std::exception &x){
413         std::string msg = "Exception caught in Thelen2003Muscle::"
414                           "calcMuscleLengthInfo\n"
415                            "of " + getName()  + "\n"
416                            + x.what();
417         throw OpenSim::Exception(msg);
418     }
419 }
420 
calcMusclePotentialEnergyInfo(const SimTK::State & s,MusclePotentialEnergyInfo & mpei) const421 void Thelen2003Muscle::calcMusclePotentialEnergyInfo(const SimTK::State& s,
422         MusclePotentialEnergyInfo& mpei) const
423 {
424     try {
425         // Get the quantities that we've already computed.
426         const MuscleLengthInfo &mli = getMuscleLengthInfo(s);
427         mpei.fiberPotentialEnergy = calcfpefisoPE(mli.fiberLength);
428         mpei.tendonPotentialEnergy= calcfsefisoPE(mli.tendonStrain);
429         mpei.musclePotentialEnergy=  mpei.fiberPotentialEnergy
430                                     + mpei.tendonPotentialEnergy;
431     }
432     catch(const std::exception &x){
433         std::string msg = "Exception caught in Thelen2003Muscle::"
434                           "calcMusclePotentialEnergyInfo\n"
435                            "of " + getName()  + "\n"
436                            + x.what();
437         throw OpenSim::Exception(msg);
438     }
439 }
440 
441 
calcFiberVelocityInfo(const SimTK::State & s,FiberVelocityInfo & fvi) const442 void Thelen2003Muscle::calcFiberVelocityInfo(const SimTK::State& s,
443                                                FiberVelocityInfo& fvi) const
444 {
445     try{
446         //Get the quantities that we've already computed
447             const MuscleLengthInfo &mli = getMuscleLengthInfo(s);
448 
449         //Get the static properties of this muscle
450             // double mclLength      = getLength(s);
451             double tendonSlackLen = getTendonSlackLength();
452             double optFiberLen    = getOptimalFiberLength();
453         //=========================================================================
454         // Compute fv by inverting the force-velocity relationship in the
455         // equilibrium equations
456         //=========================================================================
457 
458         //1. Get fiber/tendon kinematic information
459 
460         //clamp activation to a legal range
461         double a = getActivationModel().clampActivation(getStateVariableValue(s,
462                                           STATE_ACTIVATION_NAME));
463 
464 
465         double lce  = mli.fiberLength;
466         double phi  = mli.pennationAngle;
467         double cosphi=mli.cosPennationAngle;
468         double sinphi = mli.sinPennationAngle;
469 
470         //2. compute the tendon length ... because we can with this kind of model
471         double tl = mli.tendonLength;
472 
473 
474         //3. Compute force multipliers and the fiber velocity
475 
476             //This exception was causing headaches, so we're handling the special
477             //case of a fiber that is at its minimum length instead
478                 //SimTK_ERRCHK1_ALWAYS(cosphi > SimTK::Eps, fcnName.c_str(),
479                 //    "%s: Pennation angle is 90 degrees, and is causing a singularity",
480                 //    muscleName.c_str());
481 
482 
483             //6. Invert the force velocity curve to get dlce. Check for singularities
484             //   first
485 
486             //These exceptions were causing headaches, so we're handling the special
487             //case of a fiber that is at its minimum length instead
488                 //SimTK_ERRCHK1_ALWAYS(a > SimTK::Eps, fcnName.c_str(),
489                 //    "%s: Activation is 0, and is causing a singularity",
490                 //    muscleName.c_str());
491                 //SimTK_ERRCHK1_ALWAYS(fal > SimTK::Eps, fcnName.c_str(),
492                 //    "%s: The active force length curve value is 0,"
493                 //    " and is causing a singularity",
494                 //    muscleName.c_str());
495 
496             //Check for singularity conditions, and clamp output appropriately
497 
498         double dmcldt = getLengtheningSpeed(s);
499 
500         //default values that are appropriate when fiber length has been clamped
501         //to its minimum allowable value.
502 
503         double fse  = calcfse(tl/tendonSlackLen);
504         double fal  = mli.fiberActiveForceLengthMultiplier;
505         double fpe  = mli.fiberPassiveForceLengthMultiplier;
506 
507         double afalfv = ((fse/cosphi)-fpe); //we can do this without fear of
508                                               //a singularity because fiber length
509                                               //is clamped
510         double fv     = afalfv/(a*fal);
511         double dlceN  = calcdlceN(a,fal,afalfv);
512         double dlce   = dlceN*getMaxContractionVelocity()*optFiberLen;
513         double tanPhi = tan(phi);
514         double dphidt = getPennationModel().calcPennationAngularVelocity(
515                                             tanPhi,lce,dlce);
516         double dlceAT = getPennationModel().calcFiberVelocityAlongTendon(
517                             lce, dlce, sinphi, cosphi, dphidt);
518         double dtl    = getPennationModel().calcTendonVelocity(
519                             cosphi, sinphi, dphidt, lce, dlce, dmcldt);
520 
521 
522         //Switching condition: if the fiber is clamped and the tendon and the
523         //                   : fiber are out of equilibrium
524         double fiberStateClamped = 0.0;
525         if(isFiberStateClamped(s,dlceN)){
526              dlce = 0;
527              dlceAT = 0;
528              dlceN = 0;
529              dphidt = 0;
530              dtl = dmcldt;
531              fv = 1.0;
532              fiberStateClamped = 1.0;
533         }
534 
535         //Populate the struct;
536         fvi.fiberVelocity               = dlce;
537         fvi.fiberVelocityAlongTendon    = dlceAT;
538         fvi.normFiberVelocity           = dlceN;
539 
540         fvi.pennationAngularVelocity    = dphidt;
541 
542         fvi.tendonVelocity              = dtl;
543         fvi.normTendonVelocity          = dtl/getTendonSlackLength();
544 
545         fvi.fiberForceVelocityMultiplier = fv;
546 
547         fvi.userDefinedVelocityExtras.resize(2);
548         fvi.userDefinedVelocityExtras[0]=fse;
549         fvi.userDefinedVelocityExtras[1]=fiberStateClamped;
550     }
551     catch(const std::exception &x){
552         std::string msg = "Exception caught in Thelen2003Muscle::"
553                             "calcFiberVelocityInfo\n"
554                             "of " + getName()  + "\n"
555                             + x.what();
556         throw OpenSim::Exception(msg);
557     }
558 }
559 
560 
561 //=======================================
562 // computeFiberVelocityInfo helper functions
563 //=======================================
564 
calcMuscleDynamicsInfo(const SimTK::State & s,MuscleDynamicsInfo & mdi) const565 void Thelen2003Muscle::calcMuscleDynamicsInfo(const SimTK::State& s,
566                                                MuscleDynamicsInfo& mdi) const
567 {
568     try {
569         //Get the quantities that we've already computed
570         const MuscleLengthInfo &mli = getMuscleLengthInfo(s);
571         const FiberVelocityInfo &mvi = getFiberVelocityInfo(s);
572         //Get the static properties of this muscle
573         // double mclLength      = getLength(s);
574         double tendonSlackLen = getTendonSlackLength();
575         double optFiberLen    = getOptimalFiberLength();
576         double fiso           = getMaxIsometricForce();
577         double penHeight      = getPennationModel().getParallelogramHeight();
578 
579         //=========================================================================
580         // Compute required quantities
581         //=========================================================================
582         //1. Get fiber/tendon kinematic information
583         double a = getActivationModel().clampActivation(
584                        getStateVariableValue(s, STATE_ACTIVATION_NAME) );
585 
586         double lce      = mli.fiberLength;
587         double fiberStateClamped = mvi.userDefinedVelocityExtras[1];
588         double dlce     = mvi.fiberVelocity;
589         double phi      = mli.pennationAngle;
590         double cosphi   = mli.cosPennationAngle;
591         // double sinphi   = mli.sinPennationAngle;
592 
593         double tl   = mli.tendonLength;
594         double dtl  = mvi.tendonVelocity;
595         // double tlN  = mli.normTendonLength;
596 
597         //Default values appropriate when the fiber is clamped to its minimum length
598         //and is generating no force
599 
600         //These quantities should already be set to legal values from
601         //calcFiberVelocityInfo
602         double fal  = mli.fiberActiveForceLengthMultiplier;
603         double fpe  = mli.fiberPassiveForceLengthMultiplier;
604         double fv   = mvi.fiberForceVelocityMultiplier;
605         double fse  = mvi.userDefinedVelocityExtras[0];
606 
607         double aFm          = 0; //active fiber force
608         double Fm           = 0;
609         double dFm_dlce     = 0;
610         double dFmAT_dlce   = 0;
611         double dFmAT_dlceAT = 0;
612         double dFt_dtl      = 0;
613         double Ke           = 0;
614 
615         if(fiberStateClamped < 0.5){
616             aFm          = calcActiveFm(a,fal,fv,fiso);
617             Fm           = calcFm(a,fal,fv,fpe,fiso);
618             dFm_dlce     = calcDFmDlce(lce,a,fv,fiso,optFiberLen);
619             dFmAT_dlce   = calcDFmATDlce(lce,phi,cosphi,Fm,dFm_dlce,penHeight);
620 
621             //The expression below is correct only because we are using a pennation
622             //model that has a parallelogram of constant height.
623             dFmAT_dlceAT= dFmAT_dlce*cosphi;
624 
625             dFt_dtl = calcDFseDtl(tl, fiso, tendonSlackLen);
626 
627             //Compute the stiffness of the whole muscle/tendon complex
628             Ke = (dFmAT_dlceAT*dFt_dtl)/(dFmAT_dlceAT+dFt_dtl);
629         }
630 
631         mdi.activation                   = a;
632         mdi.fiberForce                   = Fm;
633         mdi.fiberForceAlongTendon        = Fm*cosphi;
634         mdi.normFiberForce               = Fm/fiso;
635         mdi.activeFiberForce             = aFm;
636         mdi.passiveFiberForce            = fpe*fiso;
637 
638         mdi.tendonForce                  = fse*fiso;
639         mdi.normTendonForce              = fse;
640 
641         mdi.fiberStiffness               = dFm_dlce;
642         mdi.fiberStiffnessAlongTendon    = dFmAT_dlceAT;
643         mdi.tendonStiffness              = dFt_dtl;
644         mdi.muscleStiffness              = Ke;
645 
646 
647 
648         //Check that the derivative of system energy less work is zero within
649         //a reasonable numerical tolerance. Throw an exception if this is not true
650         double dFibPEdt     = fpe*fiso*dlce;
651         double dTdnPEdt     = fse*fiso*dtl;
652         double dFibWdt      = -mdi.activeFiberForce*mvi.fiberVelocity;
653         double dmcldt       = getLengtheningSpeed(s);
654         double dBoundaryWdt = mdi.tendonForce * dmcldt;
655         double ddt_KEPEmW   = dFibPEdt+dTdnPEdt-dFibWdt-dBoundaryWdt;
656         SimTK::Vector userVec(1);
657         userVec(0) = ddt_KEPEmW;
658         mdi.userDefinedDynamicsExtras = userVec;
659 
660         /////////////////////////////
661         //Populate the power entries
662         /////////////////////////////
663         mdi.fiberActivePower             = dFibWdt;
664         mdi.fiberPassivePower            = -dFibPEdt;
665         mdi.tendonPower                  = -dTdnPEdt;
666         mdi.musclePower                  = -dBoundaryWdt;
667 
668     }
669     catch(const std::exception &x) {
670         std::string msg = "Exception caught in Thelen2003Muscle::"
671                             "calcMuscleDynamicsInfo\n"
672                             "of " + getName()  + "\n"
673                             + x.what();
674         throw OpenSim::Exception(msg);
675     }
676 
677 }
678 
getMinimumFiberLength() const679 double Thelen2003Muscle::getMinimumFiberLength() const
680 {
681     return getPennationModel().getMinimumFiberLength();
682 }
683 
684 
685 bool Thelen2003Muscle::
isFiberStateClamped(const SimTK::State & s,double dlceN) const686     isFiberStateClamped(const SimTK::State& s, double dlceN) const
687 {
688     bool clamped = false;
689 
690     //Is the fiber length  clamped and it is shortening, then the fiber length
691     //not valid
692     if( (getStateVariableValue(s, STATE_FIBER_LENGTH_NAME)
693             <= getMinimumFiberLength())
694         && dlceN <= 0){
695         clamped = true;
696     }
697 
698     return clamped;
699 }
700 
701 //==============================================================================
702 // ActivationFiberLength.h Interface
703 //==============================================================================
704 
705 
706 /** Get the rate change of activation */
calcActivationRate(const SimTK::State & s) const707 double Thelen2003Muscle::calcActivationRate(const SimTK::State& s) const
708 {
709     double excitation = getExcitation(s);
710     double activation = getActivation(s);
711     double dadt = getActivationModel().calcDerivative(activation,excitation);
712     return dadt;
713 }
714 
715 
716 //==============================================================================
717 // Numerical Guts: Initialization
718 //==============================================================================
719 std::pair<Thelen2003Muscle::StatusFromInitMuscleState,
720           Thelen2003Muscle::ValuesFromInitMuscleState>
initMuscleState(const SimTK::State & s,const double aActivation,const double aSolTolerance,const int aMaxIterations) const721 Thelen2003Muscle::initMuscleState(const SimTK::State& s,
722                                   const double aActivation,
723                                   const double aSolTolerance,
724                                   const int aMaxIterations) const
725 {
726     // Using short variable names to facilitate writing out long equations
727     const double ma = aActivation;
728     const double ml = getLength(s);
729     const double dml= getLengtheningSpeed(s);
730 
731     //Shorter version of the constants
732     const double tsl = getTendonSlackLength();
733     const double ofl = getOptimalFiberLength();
734     const double ophi= getPennationAngleAtOptimalFiberLength();
735     const double vol = ofl * sin(ophi);
736     const double fiso= getMaxIsometricForce();
737     const double vmax = getMaxContractionVelocity();
738 
739     //Shorter version of normalized muscle multipliers
740     double fse = 0; //Normalized tendon (series element) force
741     double fal = 0; //Normalized active force length multiplier
742     double fpe = 0; //Normalized parallel element force
743     double fv  = 0; //Normalized force-velocity multiplier
744 
745     //*******************************
746     //Position level
747     double tl  = getTendonSlackLength()*1.01;
748     double lce = getPennationModel().calcFiberLength(ml, tl);
749 
750     double phi    = 0.0;
751     double cosphi = 1.0;
752 
753     //Normalized quantities
754     double tlN  = tl/tsl;
755     double lceN = lce/ofl;
756     //Velocity level
757     double dtl      = 0;
758     double dlce = 0; // (dml - dtl) * cos(phi);
759     double dlceN = 0; // dlce / (vmax*ofl);
760     double dphi = 0; // -(dlce / lce)*tan(phi);
761     // double dlceAT   = dlce*cosphi -vol*dphi;
762 
763     //*******************************
764     //Internal variables for the loop
765     double Fm = 0;          // Muscle force
766     double FmAT=0;          // Muscle force along tendon
767     double Ft = 0;          // Tendon force
768     double ferr = SimTK::MostPositiveReal;        // Solution error
769 
770     double dFm_dlce     = 0;  // Partial derivative of muscle force w.r.t. lce
771     double dFmAT_dlce   = 0;  // Partial derivative of muscle force along
772                                // tendon w.r.t. lce
773     double dFmAT_dlceAT = 0;  // Partial derivative of muscle force along
774                                // tendon w.r.t. lce along the tendon.
775 
776     double dFt_d_lce   = 0;  // Partial derivative of tendon force w.r.t. lce
777     double dFt_d_tl   = 0;   // Partial derivative of tendon force w.r.t. tl
778 
779     double dferr_d_lce = 0;  // Partial derivative of the solution error w.r.t
780                              // lce
781     double delta_lce   = 0;  // Chance in lce
782 
783     double Ke          = 0;  // Linearized local stiffness of the muscle
784 
785     //*******************************
786     // Helper functions
787     //Update position level quantities, only if they won't go singular
788     auto positionFunc = [&] {
789         phi = getPennationModel().calcPennationAngle(lce);
790         cosphi = cos(phi);
791         tl = ml - lce*cosphi;
792         lceN = lce / ofl;
793         tlN = tl / tsl;
794     };
795 
796     // Functional to update the force multipliers
797     auto multipliersFunc = [&] {
798         fse = calcfse(tlN);
799         fal = calcfal(lceN);
800         fpe = calcfpe(lceN);
801     };
802 
803     // Functional to compute the equilibrium force error
804     auto ferrFunc = [&] {
805         Fm = (ma*fal*fv + fpe)*fiso;
806         FmAT = Fm*cosphi;
807         Ft = fse*fiso;
808         ferr = FmAT - Ft;
809     };
810 
811     // Functional to compute the partial derivative of muscle force w.r.t. lce
812     auto partialsFunc = [&] {
813         dFm_dlce = calcDFmDlce(lce, ma, fv, fiso, ofl);
814         dFmAT_dlce = calcDFmATDlce(lce, phi, cosphi, Fm, dFm_dlce, vol);
815         dFmAT_dlceAT = dFmAT_dlce*cosphi;
816 
817         dFt_d_tl = calcDFseDtl(tl, fiso, tsl);
818         dFt_d_lce = calcDFseDlce(tl, lce, phi, cosphi, fiso, tsl, vol);
819     };
820 
821     // Functional to estimate fiber velocity and force-velocity multiplier
822     // from the relative fiber and tendon stiffnesses from above partials
823     auto velocityFunc = [&] {
824         //Update velocity level quantities: share the muscle velocity
825         //between the tendon and the fiber according to their relative
826         //stiffness:
827         //
828         //Fm = Ft at equilibrium
829         //Fm = Km*lceAT
830         //Ft = Kt*xt
831         //dFm_d_xm = Km*dlceAT + dKm_d_t*lceAT (assume dKm_d_t = 0)
832         //dFt_d_xt = Kt*dtl + dKt_d_t*dtl (assume dKt_d_t = 0)
833         //
834         //This is a heuristic. The above assumptions are necessary as
835         //computing the partial derivatives of Km or Kt w.r.t. requires
836         //acceleration level knowledge, which is not available in
837         //general.
838 
839         //Stiffness of the muscle is the stiffness of the tendon and the
840         //fiber (along the tendon) in series
841 
842         //the if statement here is to handle the special case when the
843         //negative stiffness of the fiber (which could happen in this
844         // model) is equal to the positive stiffness of the tendon.
845         if (abs(dFmAT_dlceAT + dFt_d_tl) > SimTK::SignificantReal
846             && tl > getTendonSlackLength()) {
847 
848             //Ke = (dFmAT_dlceAT*dFt_d_tl) / (dFmAT_dlceAT + dFt_d_tl);
849             // resultant stiffness = k1/(k1+k2)
850             dtl = (dFmAT_dlceAT / (dFmAT_dlceAT + dFt_d_tl)) * dml;
851         }
852         else {
853             dtl = dml;
854         }
855 
856         // Update fiber velocity
857         dlce = getPennationModel().calcFiberVelocity(cosphi, dml, dtl);
858         dlceN = dlce / (vmax*ofl);
859         // Update the force-velocity multiplier
860         fv = calcfvInv(ma, fal, dlceN, aSolTolerance, 100);
861     };
862 
863     //*******************************
864     //Initialize the loop
865     int iter = 0;
866 
867     // Estimate the position level quantities (lengths, angles) of the muscle
868     positionFunc();
869 
870     // Multipliers based on initial fiber-length estimate
871     multipliersFunc();
872 
873     // Starting guess at the force-velocity multiplier
874     fv = 1.0;
875 
876     // Estimate partial derivatives of muscle forces
877     partialsFunc();
878 
879     // Update the velocity (and multiplier) estimate from partials
880     velocityFunc();
881 
882     // Compute the equilibrium error with velocity estimate
883     ferrFunc();
884 
885     // Update the partial derivatives of the force error w.r.t. lce with
886     // newly estimated fv
887     partialsFunc();
888 
889     double ferrPrev = ferr;
890     double lcePrev = lce;
891 
892     double h = 1.0;
893     while( (abs(ferr) > aSolTolerance)  && (iter < aMaxIterations)) {
894         // Compute the search direction
895         dferr_d_lce = dFmAT_dlce - dFt_d_lce;
896         h = 1.0;
897 
898         while (abs(ferr) >= abs(ferrPrev)) {
899              // Compute the Newton step
900             delta_lce = -h*ferrPrev / dferr_d_lce;
901             // Take a Newton Step if the step is nonzero
902             if (abs(delta_lce) > SimTK::SignificantReal)
903                 lce = lcePrev + delta_lce;
904             else {
905                 // We've stagnated or hit a limit; assume we are hitting local
906                 // minimum and attempt to approach from the other direction.
907                 lce = lcePrev - sign(delta_lce)*SimTK::SqrtEps;
908                 // Force a break, which will update the derivatives of
909                 // the muscle force and estimate of the fiber-velocity
910                 h = 0;
911             }
912 
913             if (lce < getMinimumFiberLength()) {
914                 lce = getMinimumFiberLength();
915             }
916 
917             // Update the muscles's position level quantities (lengths, angles)
918             positionFunc();
919 
920             // Update the muscle force multipliers
921             multipliersFunc();
922 
923             // Compute the force error assuming fiber-velocity is unchanged
924             ferrFunc();
925 
926             if (h <= SimTK::SqrtEps ) {
927                 break;
928             }
929             else
930                 h = 0.5*h;
931         }
932 
933         ferrPrev = ferr;
934         lcePrev = lce;
935 
936         // Update the partial derivative of the force error w.r.t. lce
937         partialsFunc();
938         // Update velocity estimate and velocity multiplier
939         velocityFunc();
940 
941         iter++;
942     }
943 
944     // Populate the result map.
945     ValuesFromInitMuscleState resultValues;
946 
947     if (abs(ferr) < aSolTolerance) {  // The solution converged.
948 
949         resultValues["solution_error"] = ferr;
950         resultValues["iterations"]     = (double)iter;
951         resultValues["fiber_length"]   = lce;
952         resultValues["passive_force"]  = fpe*fiso;
953         resultValues["tendon_force"]   = fse*fiso;
954 
955         return std::pair<StatusFromInitMuscleState, ValuesFromInitMuscleState>
956             (StatusFromInitMuscleState::Success_Converged, resultValues);
957     }
958 
959     // Fiber length is at or exceeds its lower bound.
960     if (lce <= getMinimumFiberLength()) {
961 
962         lce = getMinimumFiberLength();
963         phi = getPennationModel().calcPennationAngle(lce);
964         cosphi = cos(phi);
965         tl  = getPennationModel().calcTendonLength(cosphi,lce,ml);
966         lceN = lce/ofl;
967         tlN  = tl/tsl;
968         fse = calcfse(tlN);
969         fpe = calcfpe(lceN);
970 
971         resultValues["solution_error"] = ferr;
972         resultValues["iterations"]     = (double)iter;
973         resultValues["fiber_length"]   = lce;
974         resultValues["passive_force"]  = fpe*fiso;
975         resultValues["tendon_force"]   = fse*fiso;
976 
977         return std::pair<StatusFromInitMuscleState, ValuesFromInitMuscleState>
978            (StatusFromInitMuscleState::Warning_FiberAtLowerBound, resultValues);
979     }
980 
981     resultValues["solution_error"] = ferr;
982     resultValues["iterations"]     = (double)iter;
983     resultValues["fiber_length"]   = SimTK::NaN;
984     resultValues["passive_force"]  = SimTK::NaN;
985     resultValues["tendon_force"]   = SimTK::NaN;
986 
987     return std::pair<StatusFromInitMuscleState, ValuesFromInitMuscleState>
988         (StatusFromInitMuscleState::Failure_MaxIterationsReached, resultValues);
989 }
990 
991 //==============================================================================
992 //
993 // STIFFNESS RELATED FUNCTIONS
994 //
995 //==============================================================================
calcFm(double ma,double fal,double fv,double fpe,double fiso) const996 double Thelen2003Muscle::calcFm(double ma, double fal, double fv, double fpe,
997                                  double fiso) const
998 {
999     double Fm = (ma*fal*fv + fpe)*fiso;
1000     return Fm;
1001 }
1002 
calcActiveFm(double ma,double fal,double fv,double fiso) const1003 double Thelen2003Muscle::calcActiveFm(double ma, double fal, double fv,
1004                                  double fiso) const
1005 {
1006     double aFm = (ma*fal*fv)*fiso;
1007     return aFm;
1008 }
1009 
calcDFmDlce(double lce,double ma,double fv,double fiso,double ofl) const1010 double Thelen2003Muscle::calcDFmDlce(double lce, double ma, double fv,
1011                                       double fiso, double ofl) const
1012 {
1013 
1014             double lceN = lce/ofl;
1015             double dfal_d_lceN = calcDfalDlceN(lceN);
1016             double dfpe_d_lceN = calcDfpeDlceN(lceN);
1017 
1018             double dFm_d_lce = ((ma*fv)*dfal_d_lceN + dfpe_d_lceN)*fiso
1019                               *(1/ofl);
1020             return dFm_d_lce;
1021 }
1022 
calcDFmATDlce(double lce,double phi,double cosphi,double Fm,double d_Fm_d_lce,double penHeight) const1023 double Thelen2003Muscle::calcDFmATDlce(double lce, double phi, double cosphi,
1024     double Fm, double d_Fm_d_lce, double penHeight) const
1025 {
1026     //SINGULARITY: when vol*vol/(lce*lce) = 1,same as phi=pi/2
1027     double tmp1 = penHeight*penHeight;
1028     double tmp2 = lce*lce;
1029     double tmp3 = tmp2*lce;
1030     double dcosphi_d_lce = (tmp1 /(tmp3*pow((1-(tmp1/tmp2)),0.5 ) ));
1031 
1032 
1033     double d_FmAT_d_lce = d_Fm_d_lce*cosphi + Fm*dcosphi_d_lce;
1034     return d_FmAT_d_lce;
1035 }
1036 
calcDFseDlce(double tl,double lce,double phi,double cosphi,double fiso,double tsl,double vol) const1037 double Thelen2003Muscle::calcDFseDlce(double tl, double lce, double phi, double cosphi,
1038                                       double fiso, double tsl, double vol) const
1039 {
1040     double tlN = tl/tsl;
1041         //SINGULARITY: When lce = 0
1042     double tmp1 = vol/lce;
1043 
1044     //SINGULARITY: when vol/lce = 1 - equivalent to when phi = pi/2
1045     double dphi_d_lce = -vol / ( lce*lce * pow( (1-(tmp1)*(tmp1)),0.5));
1046     double dtl_d_lce  = -cos(phi) + lce*sin(phi)*dphi_d_lce;
1047 
1048     double dfse_d_tlN  = calcDfseDtlN(tlN);
1049             tmp1 = (fiso/tsl);
1050     double dFt_d_lce = dfse_d_tlN*dtl_d_lce*tmp1;
1051     return dFt_d_lce;
1052 }
1053 
calcDFseDtl(double tl,double fiso,double tsl) const1054 double Thelen2003Muscle::calcDFseDtl(double tl, double fiso, double tsl) const
1055 {
1056     double dfse_d_tlN  = calcDfseDtlN(tl/tsl);
1057     double tmp1 = (fiso/tsl);
1058     double dFt_d_tl= dfse_d_tlN*tmp1;
1059     return dFt_d_tl;
1060 }
1061 
1062 
1063 //==============================================================================
1064 //
1065 // Convenience Method
1066 //
1067 //==============================================================================
printCurveToCSVFile(const CurveType ctype,const std::string & path) const1068 void Thelen2003Muscle::printCurveToCSVFile(const CurveType ctype,
1069                                            const std::string&path) const
1070 {
1071     std::string fname = getName();
1072 
1073     switch(ctype){
1074         case FiberActiveForceLength:{
1075 
1076             fname.append("_FiberActiveForceLength.csv");
1077 
1078             SimTK::Matrix data(100,3);
1079             SimTK::Array_<std::string> colNames(data.ncol());
1080             colNames[0] = "lceN";
1081             colNames[1] = "fal";
1082             colNames[2] = "DfalDlceN";
1083 
1084             double delta = (2.0-0.0)/(data.nrow()-1.0);
1085             double lceN = 0;
1086             double lceN0 = 0;
1087             for(int i=0; i<data.nrow(); i++){
1088                 lceN = lceN0 + i*delta;
1089                 data(i,0) = lceN;
1090                 data(i,1) = calcfal(lceN);
1091                 data(i,2) = calcDfalDlceN(lceN);
1092             }
1093 
1094             printMatrixToFile(data,colNames,path, fname);
1095 
1096         }break;
1097         case FiberPassiveForceLength:{
1098             fname.append("_FiberPassiveForceLength.csv");
1099 
1100             SimTK::Matrix data(100,3);
1101             SimTK::Array_<std::string> colNames(data.ncol());
1102             colNames[0] = "lceN";
1103             colNames[1] = "fpe";
1104             colNames[2] = "DfpeDlceN";
1105 
1106             double lceNMax = get_FmaxMuscleStrain()+1.0;
1107             double lceNMin = 1.0;
1108             double lceNExtra = 0.1*(lceNMax-lceNMin);
1109 
1110             double delta = ((lceNMax + lceNExtra)-(lceNMin-lceNExtra))
1111                            /(data.nrow()-1.0);
1112             double lceN = 0;
1113             double lceN0 = lceNMin-lceNExtra;
1114 
1115             for(int i=0; i<data.nrow(); i++){
1116                 lceN = lceN0 + i*delta;
1117                 data(i,0) = lceN;
1118                 data(i,1) = calcfpe(lceN);
1119                 data(i,2) = calcDfpeDlceN(lceN);
1120             }
1121 
1122             printMatrixToFile(data,colNames,path, fname);
1123 
1124         }break;
1125         case FiberForceVelocity:{
1126 
1127             fname.append("_FiberForceVelocity.csv");
1128 
1129             SimTK::Matrix data(1000,5);
1130             SimTK::Array_<std::string> colNames(data.ncol());
1131             colNames[0] = "a";
1132             colNames[1] = "fal";
1133             colNames[2] = "dlceN";
1134             colNames[3] = "fv";
1135             colNames[4] = "DfvDdlceN";
1136 
1137             //Thelen's fv varies with activation
1138             double a = 0;
1139             double a0 = 0.1;
1140             double da = 0.1;
1141 
1142             //Thelen's fv varies with fiber velocity
1143             double dlceNMax = 1.0;
1144             double dlceNMin = -1.0;
1145             double dlceNExtra = 0.1*(dlceNMax-dlceNMin);
1146 
1147             double delta = ((dlceNMax + dlceNExtra)-(dlceNMin-dlceNExtra))
1148                            /(100.0-1.0);
1149             double dlceN = 0;
1150             double dlceN0 = dlceNMin-dlceNExtra;
1151 
1152             //Thelen's fv also varies with the value of the active force
1153             //length curve. Most experiments that establish fv do so at
1154             //a normalized fiber length of 1.0, where fal is 1. We will
1155             //also make that assumption here
1156 
1157             double fvInv;
1158             double DdlceDaFalFv = 0;
1159             double DfvDdlceN = 0;
1160 
1161             const double fal = 1.0;
1162 
1163             int idx = 0;
1164             for(int i=0; i<=9; i++){
1165                 a = a0 + da*i;
1166 
1167                 for(int j=0; j<100; j++){
1168                     dlceN = dlceN0 + j*delta;
1169                     fvInv = calcfvInv(a, fal, dlceN, 1e-6,100);
1170                     DdlceDaFalFv = calcDdlceDaFalFv(a,fal,a*fal*fvInv);
1171                     DfvDdlceN = (1/(DdlceDaFalFv*a*fal));
1172 
1173                     data(idx,0) = a;
1174                     data(idx,1) = fal;
1175                     data(idx,2) = dlceN;
1176                     data(idx,3) = fvInv;
1177                     data(idx,4) = DfvDdlceN;
1178 
1179                     idx++;
1180                 }
1181             }
1182             printMatrixToFile(data,colNames,path, fname);
1183 
1184         }break;
1185         case TendonForceLength:{
1186             fname.append("_TendonForceLength.csv");
1187 
1188             SimTK::Matrix data(100,3);
1189             SimTK::Array_<std::string> colNames(data.ncol());
1190             colNames[0] = "ltN";
1191             colNames[1] = "fse";
1192             colNames[2] = "DfseDtlN";
1193 
1194             double ltNMax = get_FmaxTendonStrain()+1.0;
1195             double ltNMin = 1.0;
1196             double ltNExtra = 0.1*(ltNMax-ltNMin);
1197 
1198             double delta = ((ltNMax + ltNExtra)-(ltNMin-ltNExtra))
1199                            /(data.nrow()-1.0);
1200             double ltN = 0;
1201             double ltN0 = ltNMin-ltNExtra;
1202 
1203             for(int i=0; i<data.nrow(); i++){
1204                 ltN = ltN0 + i*delta;
1205                 data(i,0) = ltN;
1206                 data(i,1) = calcfse(ltN);
1207                 data(i,2) = calcDfseDtlN(ltN);
1208             }
1209 
1210             printMatrixToFile(data,colNames,path, fname);
1211 
1212         }break;
1213         default:{
1214             std::string msg = "Thelen2003Muscle::printCurveToCSVFile ";
1215             msg.append(getName());
1216             msg.append(" invalid curve type");
1217             SimTK_ASSERT(false,msg.c_str());
1218         }
1219     }
1220 
1221 }
1222 
1223 /*
1224 This function will print cvs file of the column vector col0 and the matrix data
1225 
1226 @params data: A matrix of data
1227 @params filename: The name of the file to print
1228 */
1229 void Thelen2003Muscle::
printMatrixToFile(SimTK::Matrix & data,SimTK::Array_<std::string> & colNames,const std::string & path,const std::string & filename) const1230     printMatrixToFile(SimTK::Matrix& data, SimTK::Array_<std::string>& colNames,
1231     const std::string& path, const std::string& filename) const
1232 {
1233 
1234     ofstream datafile;
1235     std::string fullpath = path;
1236 
1237     if(fullpath.length() > 0)
1238         fullpath.append("/");
1239 
1240     fullpath.append(filename);
1241 
1242     datafile.open(fullpath.c_str(),std::ios::out);
1243 
1244     if(!datafile){
1245         datafile.close();
1246         string name = getName();
1247         SimTK_ERRCHK2_ALWAYS( false,
1248                 "Thelen2003Muscle::printMatrixToFile",
1249                 "%s: Failed to open the file path: %s",
1250                 name.c_str(),
1251                 fullpath.c_str());
1252     }
1253 
1254 
1255     for(int i = 0; i < (signed)colNames.size(); i++){
1256         if(i < (signed)colNames.size()-1)
1257             datafile << colNames[i] << ",";
1258         else
1259             datafile << colNames[i] << "\n";
1260     }
1261 
1262     for(int i = 0; i < data.nrow(); i++){
1263         for(int j = 0; j < data.ncol(); j++){
1264             if(j<data.ncol()-1)
1265                 datafile << data(i,j) << ",";
1266             else
1267                 datafile << data(i,j) << "\n";
1268         }
1269     }
1270     datafile.close();
1271 }
1272 
1273 //==============================================================================
1274 //
1275 // TENDON RELATED FUNCTIONS
1276 //
1277 //==============================================================================
1278 
calcfse(const double tlN) const1279 double Thelen2003Muscle::calcfse(const double tlN) const
1280 {
1281     double x = tlN-1;
1282     double e0 = get_FmaxTendonStrain();
1283 
1284     /*The paper reports etoe = 0.609e0, however, this is a severely rounded off
1285         The exact answer, to SimTK::Eps is
1286         etoe =  99*e0*e^3 / ( 166*e^3 - 67)
1287         klin =  67 /( 100*(e0 - (99*e0*e^3)/(166*e^3-67)) )
1288         See thelenINIT_20120127.mw for details
1289     */
1290     double kToe = 3.0;
1291     double Ftoe = 33.0/100.0;
1292     double t1   = exp(0.3e1);
1293     double eToe = (0.99e2*e0*t1) / (0.166e3*t1 - 0.67e2);
1294     t1 = exp(0.3e1);
1295     double klin = (0.67e2/0.100e3)
1296                 * 1.0/(e0 - (0.99e2*e0*t1) / (0.166e3*t1 - 0.67e2));
1297 
1298     //Compute tendon force
1299     double fse = 0;
1300     if (x > eToe){
1301         fse = klin*(x-eToe)+Ftoe;
1302     }else if (x>0.0){
1303         fse =(Ftoe/(exp(kToe)-1.0))*(exp(kToe*x/eToe)-1.0);
1304     }else{
1305         fse=0.;}
1306 
1307     return fse;
1308 }
1309 
1310 
calcDfseDtlN(const double tlN) const1311 double Thelen2003Muscle::calcDfseDtlN(const double tlN) const {
1312     double x = tlN-1;
1313     double e0 = get_FmaxTendonStrain();
1314 
1315     /*The paper reports etoe = 0.609e0, however, this is a severely rounded off
1316     result of the exact answer:
1317         etoe =  99*e0*e^3 / ( 166*e^3 - 67)
1318         See thelenINIT_20120127.mw for details
1319     */
1320     double kToe = 3.0;
1321     double Ftoe = 33.0/100.0;
1322 
1323     double t1   = exp(0.3e1);
1324     double eToe = (0.99e2*e0*t1) / (0.166e3*t1 - 0.67e2);
1325 
1326     /*The paper reports etoe = 0.609e0, however, this is a severely rounded off
1327     result of the exact answer:
1328         klin =  67 /( 100*(e0 - (99*e0*e^3)/(166*e^3-67)) )
1329         See thelenINIT_20120127.mw for details
1330     */
1331     t1 = exp(0.3e1);
1332     double klin = (0.67e2/0.100e3)
1333                 * 1.0/(e0 - (0.99e2*e0*t1) / (0.166e3*t1 - 0.67e2));
1334 
1335     //Compute tendon force
1336     double dfse_d_dtlN = 0;
1337     if (x > eToe){
1338         dfse_d_dtlN = klin;
1339     }else if (x>0.0){
1340         dfse_d_dtlN =(Ftoe/(exp(kToe)-1.0)) * (kToe/eToe) * (exp(kToe*x/eToe));
1341     }else{
1342         dfse_d_dtlN=0.;}
1343 
1344     return dfse_d_dtlN;
1345 }
1346 
calcfsefisoPE(double tendonStrain) const1347 double Thelen2003Muscle::calcfsefisoPE(double tendonStrain) const
1348 {
1349 
1350     double tendon_strain =  tendonStrain;
1351     double fmaxTendonStrain = get_FmaxTendonStrain();
1352 
1353     //Future optimization opportunity: precompute kToe, fToe, eToe and klin
1354     //when the muscle is initialized. Store these values rather than
1355     //computing them every time.
1356 
1357     double kToe = 3.0;
1358     double Ftoe = 33.0/100.0;
1359 
1360     double t1   = exp(0.3e1);
1361     double eToe = (0.99e2*fmaxTendonStrain*t1) / (0.166e3*t1 - 0.67e2);
1362 
1363     t1 = exp(0.3e1);
1364     double klin = (0.67e2/0.100e3)
1365                 * 1.0/(fmaxTendonStrain - (0.99e2*fmaxTendonStrain*t1)
1366                 / (0.166e3*t1 - 0.67e2));
1367 
1368     //Compute the energy stored in the tendon.
1369     //Integrals computed symbolically in muscle_kepew_20111021.mw just to check
1370     double tendonPE = 0.0;
1371     double lenR        = getTendonSlackLength();
1372     double lenTdn    = (tendon_strain+1)*lenR;
1373     double lenToe    = (eToe+1.0)*lenR;
1374     double fiso        = getMaxIsometricForce();
1375 
1376     if (tendon_strain>eToe){
1377        //compute the energy stored in the toe portion of the tendon strain curve
1378         double len = lenToe;
1379         double toePE_len = (fiso*Ftoe/(exp(kToe)-1.0))
1380                             *((lenR*eToe/kToe)
1381                             *exp(kToe*(len-lenR)/(lenR*eToe)) - len);
1382         len =  lenR;
1383         double toePE_0    =  (fiso*Ftoe/(exp(kToe)-1.0))
1384                             *((lenR*eToe/kToe)
1385                             *exp(kToe*(len-lenR)/(lenR*eToe)) - len);
1386         // double toePEtest = toePE_len-toePE_0;
1387 
1388         //compute the energy stored in the linear section of the
1389         //tendon strain curve from ..... 0 to len
1390         len = lenTdn;
1391         double linPE_len = (1.0/2.0)*(fiso*klin*(len*len)/lenR)
1392                            + fiso*len*(klin*(-1.0-eToe)+Ftoe);
1393         //ditto from 0 .... eToe
1394         len = lenToe;
1395         double linPE_eToe= (1.0/2.0)*(fiso*klin*(len*len)/lenR)
1396                             + fiso*len*(klin*(-1.0-eToe)+Ftoe);
1397 
1398         //compute the total potential energy stored in the tendon
1399          tendonPE =(toePE_len-toePE_0) + (linPE_len-linPE_eToe);
1400     }else if (tendon_strain>0.0){
1401         //PE from 0 .... len
1402         double len = lenTdn;
1403         double toePE_len = (fiso*Ftoe/(exp(kToe)-1.0)) * ((lenR*eToe/kToe)
1404             * exp(kToe*(len-lenR)/(lenR*eToe)) - len);
1405         //PE from 0 .... eToe
1406         len = lenR;
1407         double toePE_0    =  (fiso*Ftoe/(exp(kToe)-1.0)) * ((lenR*eToe/kToe)
1408             * exp(kToe*(len-lenR)/(lenR*eToe)) - len);
1409 
1410         //Compute the total PE stored in the tendon
1411         tendonPE = toePE_len-toePE_0;
1412     }else{
1413         tendonPE = 0.0;
1414     }
1415 
1416 
1417     return tendonPE;
1418 }
1419 
1420 //==============================================================================
1421 //
1422 // ACTIVE FORCE LENGTH FUNCTIONS
1423 //
1424 //==============================================================================
calcfal(const double lceN) const1425 double Thelen2003Muscle::calcfal(const double lceN) const{
1426     double kShapeActive = get_KshapeActive();
1427     double x=(lceN-1.)*(lceN-1.);
1428     double fal = exp(-x/kShapeActive);
1429     return fal;
1430 }
calcDfalDlceN(const double lceN) const1431 double Thelen2003Muscle::calcDfalDlceN(const double lceN) const {
1432     double kShapeActive = get_KshapeActive();
1433     double t1 = lceN - 0.10e1;
1434     double t2 = 0.1e1 / kShapeActive;
1435     double t4 = t1 * t1;
1436     double t6 = exp(-t4 * t2);
1437     double dfal_d_lceN = -0.2e1 * t1 * t2 * t6;
1438     return dfal_d_lceN;
1439 }
1440 
1441 //=============================================================================
1442 //
1443 // FIBER PARALLEL ELEMENT HELPER FUNCTIONS
1444 //
1445 //=============================================================================
calcfpe(const double lceN) const1446 double Thelen2003Muscle::calcfpe(const double lceN) const {
1447     double fpe = 0;
1448     double e0 = get_FmaxMuscleStrain();
1449     double kpe = get_KshapePassive();
1450 
1451     //Compute the passive force developed by the muscle
1452     if(lceN > 1.0){
1453         double t5 = exp(kpe * (lceN - 0.10e1) / e0);
1454         double t7 = exp(kpe);
1455         fpe = (t5 - 0.10e1) / (t7 - 0.10e1);
1456     }
1457     return fpe;
1458 }
1459 
calcDfpeDlceN(const double lceN) const1460 double Thelen2003Muscle::calcDfpeDlceN(const double lceN) const {
1461     double dfpe_d_lceN = 0;
1462     double e0 = get_FmaxMuscleStrain();
1463     double kpe = get_KshapePassive();
1464 
1465     if(lceN > 1.0){
1466         double t1 = 0.1e1 / e0;
1467         double t6 = exp(kpe * (lceN - 0.10e1) * t1);
1468         double t7 = exp(kpe);
1469         dfpe_d_lceN = kpe * t1 * t6 / (t7 - 0.10e1);
1470     }
1471     return dfpe_d_lceN;
1472 }
1473 
calcfpefisoPE(double lceN) const1474 double Thelen2003Muscle::calcfpefisoPE(double lceN) const
1475 {
1476     double fmaxMuscleStrain = get_FmaxMuscleStrain();
1477     double kShapePassive = get_KshapePassive();
1478 
1479     double musclePE = 0.0;
1480     //Compute the potential energy stored in the muscle
1481     if(lceN > 1.0){
1482         //Shorter variable names to make the equations readable.
1483         double lenR = getOptimalFiberLength();
1484         double fiso = getMaxIsometricForce();
1485         double len = lceN*lenR;
1486         double kpe = kShapePassive;
1487         double e0 = fmaxMuscleStrain;
1488 
1489 
1490         //PE storage at current stretch
1491         double fpePE_len = (fiso/(exp(kpe)-1))
1492             *( (lenR*e0/kpe)*exp( (kpe/e0)*( (len/lenR)-1)) - len);
1493 
1494         //PE stored between 0 and 1 for the exponential function that is
1495         //used to represent fpe for normalized fiber lengths > 1.
1496         len = lenR;
1497         double fpePE_0 = (fiso/(exp(kpe)-1))
1498             *( (lenR*e0/kpe)*exp( (kpe/e0)*( (len/lenR)-1)) - len);
1499 
1500         musclePE = fpePE_len - fpePE_0;
1501 
1502     }else{
1503         musclePE = 0.0;
1504     }
1505     return musclePE;
1506 }
1507 
1508 
1509 //=============================================================================
1510 //
1511 // FIBER FORCE - VELOCITY CURVE FUNCTIONS
1512 //
1513 //=============================================================================
1514 
calcdlceN(double act,double fal,double actFalFv) const1515 double Thelen2003Muscle::calcdlceN(double act,double fal,double actFalFv) const
1516 {
1517     //The variable names have all been switched to closely match
1518     //with the notation in Thelen 2003.
1519     double dlceN = 0.0;      //contractile element velocity
1520     double af   = get_Af();
1521 
1522     double a    = act;
1523     double afl  = a*fal; //afl = a*fl
1524     double Fm   = actFalFv;     //Fm = a*fl*fv
1525     double flen = get_Flen();
1526     // double Fmlen_afl = flen*afl;
1527 
1528     double dlcedFm = 0.0; //partial derivative of contractile element
1529                           // velocity w.r.t. Fm
1530 
1531     double b = 0;
1532     double db= 0;
1533 
1534     double Fm_asyC = 0;           //Concentric contraction asymptote
1535     double Fm_asyE = afl*flen;
1536                                 //Eccentric contraction asymptote
1537     double asyE_thresh = get_fv_linear_extrap_threshold();
1538 
1539     //If fv is in the appropriate region, use
1540     //Thelen 2003 Eqns 6 & 7 to compute dlceN
1541     if (Fm > Fm_asyC && Fm < Fm_asyE*asyE_thresh){
1542 
1543         if( Fm <= afl ){        //Muscle is concentrically contracting
1544             b = afl + Fm/af;
1545             db= 1/af;
1546         }else{                    //Muscle is eccentrically contracting
1547             b = ((2+2/af)*(afl*flen-Fm))/(flen-1);
1548             db= ((2+2/af)*(-1))/(flen-1);
1549         }
1550 
1551         dlceN = (0.25 + 0.75*a)*(Fm-afl)/b;
1552         //Scaling by VMAX is left out, and is post multiplied outside
1553         //of the function
1554 
1555 
1556     }else{  //Linear extrapolation
1557             double Fm0 = 0.0; //Last Fm value from the Thelen curve
1558 
1559             //Compute d and db/dFm from Eqn 7. of Thelen2003
1560             //for the last
1561             if(Fm <= Fm_asyC){ //Concentrically contracting
1562                 Fm0 = Fm_asyC;
1563                 b = afl + Fm0/af;
1564                 db= 1/af;
1565             }else{             //Eccentrically contracting
1566                 Fm0 = asyE_thresh*Fm_asyE;
1567                 b = ((2+2/af)*(afl*flen-Fm0))/(flen-1);
1568                 db= ((2+2/af)*(-1))/(flen-1);
1569             }
1570 
1571             //Compute the last dlceN value that falls in the region where
1572             //Thelen 2003 Eqn. 6 is valid
1573             double dlce0 = (0.25 + 0.75*a)*(Fm0-afl)/b;
1574 
1575             //Compute the dlceN/dfm of Eqn. 6 of Thelen 2003 at the last
1576             //valid point
1577             dlcedFm = (0.25 + 0.75*a)*(1)/b
1578                     - ((0.25 + 0.75*a)*(Fm0-afl)/(b*b))*db;
1579 
1580             //Linearly extrapolate Eqn. 6 from Thelen 2003 to compute
1581             //the new value for dlceN/dFm
1582             dlceN = dlce0 + dlcedFm*(Fm-Fm0);
1583         }
1584 
1585         return dlceN;
1586 }
1587 
1588 double Thelen2003Muscle::
calcfv(double aFse,double aFpe,double aFal,double aCosPhi,double aAct) const1589     calcfv(double aFse,     double aFpe, double aFal,
1590            double aCosPhi,  double aAct) const
1591 {
1592     //This only works for an equilibrium model, but its a lot less
1593     //computationally expensive (and error prone) than trying to invert the
1594     //weird function that defines the fv curve in the Thelen model
1595     double fv = ((aFse/aCosPhi) - aFpe)/(aAct*aFal);
1596     return fv;
1597 }
1598 
calcDdlceDaFalFv(double aAct,double aFal,double aFalFv) const1599 double Thelen2003Muscle::calcDdlceDaFalFv(double aAct,
1600                                           double aFal, double aFalFv) const
1601 {
1602     //The variable names have all been switched to closely match with
1603     //the notation in Thelen 2003.
1604     // double dlceN = 0.0;      //contractile element velocity
1605     double af   = get_Af();
1606 
1607     double a    = aAct;
1608     double afl  = aAct*aFal;  //afl = a*fl
1609     double Fm   = aFalFv;    //Fm = a*fl*fv
1610     double flen = get_Flen();
1611     // double Fmlen_afl = flen*aAct*aFal;
1612 
1613     double dlcedFm = 0.0; //partial derivative of contractile element
1614                           //velocity w.r.t. Fm
1615 
1616     double b = 0;
1617     double db= 0;
1618 
1619     double Fm_asyC = 0;           //Concentric contraction asymptote
1620     double Fm_asyE = aAct*aFal*flen;
1621                                 //Eccentric contraction asymptote
1622     double asyE_thresh = get_fv_linear_extrap_threshold();
1623 
1624     //If fv is in the appropriate region, use
1625     //Thelen 2003 Eqns 6 & 7 to compute dlceN
1626     if (Fm > Fm_asyC && Fm < Fm_asyE*asyE_thresh){
1627 
1628         if( Fm <= afl ){        //Muscle is concentrically contracting
1629             b = afl + Fm/af;
1630             db= 1/af;
1631         }else{                    //Muscle is eccentrically contracting
1632             b = ((2+2/af)*(afl*flen-Fm))/(flen-1);
1633             db= ((2+2/af)*(-1))/(flen-1);
1634         }
1635 
1636         //This variable may have future use outside this function
1637         dlcedFm = (0.25 + 0.75*a)*(1)/b - ((0.25 + 0.75*a)*(Fm-afl)/(b*b))*db;
1638 
1639     }else{  //Linear extrapolation
1640             double Fm0 = 0.0; //Last Fm value from the Thelen curve
1641 
1642             //Compute d and db/dFm from Eqn 7. of Thelen2003
1643             //for the last
1644             if(Fm <= Fm_asyC){ //Concentrically contracting
1645                 Fm0 = Fm_asyC;
1646                 b = afl + Fm0/af;
1647                 db= 1/af;
1648             }else{             //Eccentrically contracting
1649                 Fm0 = asyE_thresh*Fm_asyE;
1650                 b = ((2+2/af)*(afl*flen-Fm0))/(flen-1);
1651                 db= ((2+2/af)*(-1))/(flen-1);
1652             }
1653 
1654 
1655             //Compute the dlceN/dfm of Eqn. 6 of Thelen 2003 at the last
1656             //valid point
1657             dlcedFm = (0.25 + 0.75*a)*(1)/b
1658                 - ((0.25 + 0.75*a)*(Fm0-afl)/(b*b))*db;
1659 
1660         }
1661 
1662         return dlcedFm;
1663 }
1664 
1665 // Compute the force-velocity multiplier by inverting Thelen 2003' f-v
1666 // equations for fiber-velocity given the active fiber force (see calcdlceN()).
1667 // This is here because it is non-trivial to correctly invert the piece-wise
1668 // continuous force velocity curve specified by the modified Thelen2003Muscle
1669 // force velocity curve. This converges quickly and is well tested.
1670 double Thelen2003Muscle::
calcfvInv(double aAct,double aFal,double dlceN,double tolerance,int maxIterations) const1671         calcfvInv(double aAct,double aFal,double dlceN,
1672                   double tolerance, int maxIterations) const
1673 {
1674     double result = SimTK::NaN;
1675     double ferr=1;
1676     double iter= 0;
1677 
1678     double dlceN1 = 0;
1679     double dlceN1_d_Fm = 0;
1680     double fv = 1;
1681     double aFalFv = fv*aAct*aFal;
1682     double delta_aFalFv = 0;
1683 
1684     while(abs(ferr) >= tolerance && iter < maxIterations)
1685     {
1686         dlceN1 = calcdlceN(aAct,aFal, aFalFv);
1687         ferr   = dlceN1-dlceN;
1688         dlceN1_d_Fm = calcDdlceDaFalFv(aAct,aFal,aFalFv);
1689 
1690 
1691         if(abs(dlceN1_d_Fm) > SimTK::SignificantReal){
1692            delta_aFalFv = -ferr/(dlceN1_d_Fm);
1693            aFalFv = aFalFv + delta_aFalFv;
1694         }
1695         iter = iter+1;
1696     }
1697 
1698     if(abs(ferr) < tolerance){
1699         result = max(0.0, aFalFv/(aAct*aFal));
1700         return result;
1701     }
1702 
1703     OPENSIM_THROW_FRMOBJ(Exception,
1704         "Solver for force-velocity multiplier failed to converge.");
1705 }
1706