1 /* -------------------------------------------------------------------------- *
2  *                            OpenSim:  Model.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): Frank C. Anderson, Peter Loan, Ayman Habib, Ajay Seth,          *
12  *            Michael Sherman                                                 *
13  *                                                                            *
14  * Licensed under the Apache License, Version 2.0 (the "License"); you may    *
15  * not use this file except in compliance with the License. You may obtain a  *
16  * copy of the License at http://www.apache.org/licenses/LICENSE-2.0.         *
17  *                                                                            *
18  * Unless required by applicable law or agreed to in writing, software        *
19  * distributed under the License is distributed on an "AS IS" BASIS,          *
20  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.   *
21  * See the License for the specific language governing permissions and        *
22  * limitations under the License.                                             *
23  * -------------------------------------------------------------------------- */
24 
25 //=============================================================================
26 // INCLUDES
27 //=============================================================================
28 
29 #include <OpenSim/Common/IO.h>
30 #include <OpenSim/Common/XMLDocument.h>
31 #include <OpenSim/Common/ScaleSet.h>
32 #include <OpenSim/Common/Storage.h>
33 #include <OpenSim/Simulation/SimbodyEngine/FreeJoint.h>
34 #include <OpenSim/Simulation/SimbodyEngine/SimbodyEngine.h>
35 #include <OpenSim/Simulation/SimbodyEngine/WeldConstraint.h>
36 #include <OpenSim/Simulation/SimbodyEngine/PointConstraint.h>
37 #include <OpenSim/Common/Constant.h>
38 #include <OpenSim/Simulation/CoordinateReference.h>
39 
40 #include "SimTKcommon/internal/SystemGuts.h"
41 
42 #include "Model.h"
43 
44 #include "Actuator.h"
45 #include "AnalysisSet.h"
46 #include "BodySet.h"
47 #include "ComponentSet.h"
48 #include "ContactGeometrySet.h"
49 #include "ControllerSet.h"
50 #include "CoordinateSet.h"
51 #include "ForceSet.h"
52 #include "Ligament.h"
53 #include "MarkerSet.h"
54 #include "ProbeSet.h"
55 
56 #include <iostream>
57 #include <string>
58 
59 #include <OpenSim/Simulation/AssemblySolver.h>
60 
61 using namespace std;
62 using namespace OpenSim;
63 using namespace SimTK;
64 
65 
66 //=============================================================================
67 // STATICS
68 //=============================================================================
69 
70 
71 //=============================================================================
72 // CONSTRUCTOR(S) AND DESTRUCTOR
73 //=============================================================================
74 //_____________________________________________________________________________
75 /*
76  * Default constructor.
77  */
Model()78 Model::Model() : ModelComponent(),
79     _fileName("Unassigned"),
80     _analysisSet(AnalysisSet()),
81     _coordinateSet(CoordinateSet()),
82     _workingState(),
83     _useVisualizer(false),
84     _allControllersEnabled(true)
85 {
86     constructProperties();
87     setNull();
88     finalizeFromProperties();
89 }
90 //_____________________________________________________________________________
91 /*
92  * Constructor from an XML file
93  */
Model(const string & aFileName)94 Model::Model(const string &aFileName) :
95     ModelComponent(aFileName, false),
96     _fileName("Unassigned"),
97     _analysisSet(AnalysisSet()),
98     _coordinateSet(CoordinateSet()),
99     _workingState(),
100     _useVisualizer(false),
101     _allControllersEnabled(true)
102 {
103     constructProperties();
104     setNull();
105     updateFromXMLDocument();
106     // Check is done below because only Model files have migration issues, version is not available until
107     // updateFromXMLDocument is called. Fixes core issue #2395
108     OPENSIM_THROW_IF(getDocument()->getDocumentVersion() < 10901,
109         Exception,
110         "Model file " + aFileName + " is using unsupported file format"
111         ". Please open model and save it in OpenSim version 3.3 to upgrade.");
112 
113     _fileName = aFileName;
114     cout << "Loaded model " << getName() << " from file " << getInputFileName() << endl;
115 
116     try {
117         finalizeFromProperties();
118     }
119     catch(const InvalidPropertyValue& err) {
120         cout << "WARNING: Model was unable to finalizeFromProperties.\n" <<
121             "Update the model file and reload OR update the property and call "
122             "finalizeFromProperties() on the model.\n" <<
123             "(details: " << err.what() << ")." << endl;
124     }
125 }
126 
clone() const127 Model* Model::clone() const
128 {
129     // Invoke default copy constructor.
130     Model* clone = new Model(*this);
131 
132     try {
133         clone->finalizeFromProperties();
134     }
135     catch (const InvalidPropertyValue& err) {
136         cout << "WARNING: clone() was unable to finalizeFromProperties.\n" <<
137             "Update the model and call clone() again OR update the clone's "
138             "property and call finalizeFromProperties() on it.\n"
139             "(details: " << err.what() << ")." << endl;
140     }
141 
142     return clone;
143 }
144 
145 //_____________________________________________________________________________
146 /*
147  * Override default implementation by object to intercept and fix the XML node
148  * underneath the model to match current version
149  */
150 /*virtual*/
updateFromXMLNode(SimTK::Xml::Element & aNode,int versionNumber)151 void Model::updateFromXMLNode(SimTK::Xml::Element& aNode, int versionNumber)
152 {
153     if (versionNumber < XMLDocument::getLatestVersion()){
154         cout << "Updating Model file from " << versionNumber
155             << " to latest format..." << endl;
156         // Version has to be 1.6 or later, otherwise assert
157         if (versionNumber == 10600){
158             // Get node for DynamicsEngine
159             SimTK::Xml::element_iterator engIter =
160                 aNode.element_begin("DynamicsEngine");
161             //Get node for SimbodyEngine
162             if (engIter != aNode.element_end()){
163                 SimTK::Xml::element_iterator simbodyEngIter =
164                     engIter->element_begin("SimbodyEngine");
165                 // Move all Children of simbodyEngineNode to be children of _node
166                 // we'll keep inserting before enginesNode then remove it;
167                 SimTK::Array_<SimTK::Xml::Element> elts =
168                     simbodyEngIter->getAllElements();
169                 while (elts.size() != 0){
170                     // get first child and move it to Model
171                     aNode.insertNodeAfter(aNode.element_end(),
172                         simbodyEngIter->removeNode(simbodyEngIter->element_begin()));
173                     elts = simbodyEngIter->getAllElements();
174                 }
175                 engIter->eraseNode(simbodyEngIter);
176             }
177             // Now handling the rename of ActuatorSet to ForceSet
178             XMLDocument::renameChildNode(aNode, "ActuatorSet", "ForceSet");
179         }
180         if (versionNumber < 30501) {
181             // Create JointSet node after BodySet under <OpenSimDocument>/<Model>
182             SimTK::Xml::Element jointSetElement("JointSet");
183             SimTK::Xml::Element jointObjects("objects");
184             jointSetElement.insertNodeBefore(jointSetElement.element_begin(), jointObjects);
185             SimTK::Xml::element_iterator bodySetNode = aNode.element_begin("BodySet");
186             aNode.insertNodeAfter(bodySetNode, jointSetElement);
187             // Cycle through Bodies and move their Joint nodes under the Model's JointSet
188             SimTK::Xml::element_iterator  objects_node =
189                 bodySetNode->element_begin("objects");
190             SimTK::Xml::element_iterator bodyIter =
191                 objects_node->element_begin("Body");
192             for (; bodyIter != objects_node->element_end(); ++bodyIter) {
193                 std::string body_name = bodyIter->getOptionalAttributeValue("name");
194                 if (body_name == "ground") {
195                     SimTK::Xml::Element newGroundElement("Ground");
196                     newGroundElement.setAttributeValue("name", "ground");
197 
198                     SimTK::Xml::element_iterator geometryIter =
199                         bodyIter->element_begin("geometry");
200                     if (geometryIter != bodyIter->element_end()) {
201                         SimTK::Xml::Element cloneOfGeomety = geometryIter->clone();
202                         newGroundElement.insertNodeAfter(newGroundElement.node_end(),
203                             cloneOfGeomety);
204                     }
205 
206                     SimTK::Xml::element_iterator visObjIter =
207                         bodyIter->element_begin("VisibleObject");
208                     if (visObjIter != bodyIter->element_end()) {
209                         SimTK::Xml::Element cloneOfVisObj = visObjIter->clone();
210                         newGroundElement.insertNodeAfter(newGroundElement.node_end(),
211                             cloneOfVisObj);
212                     }
213 
214                     SimTK::Xml::element_iterator wrapSetIter =
215                         bodyIter->element_begin("WrapObjectSet");
216                     if (wrapSetIter != bodyIter->element_end()) {
217                         SimTK::Xml::Element cloneOfWrapSet = wrapSetIter->clone();
218                         newGroundElement.insertNodeAfter(newGroundElement.node_end(),
219                             cloneOfWrapSet);
220                     }
221 
222                     String test;
223                     newGroundElement.writeToString(test);
224 
225                     objects_node->eraseNode(bodyIter);
226 
227                     aNode.insertNodeBefore(bodySetNode, newGroundElement);
228                     break;
229                 }
230             }
231             bodyIter = objects_node->element_begin("Body");
232             for (; bodyIter != objects_node->element_end(); ++bodyIter) {
233                 std::string body_name = bodyIter->getOptionalAttributeValue("name");
234                 SimTK::Xml::element_iterator  joint_node =
235                     bodyIter->element_begin("Joint");
236                 if (joint_node->element_begin() != joint_node->element_end()){
237                     SimTK::Xml::Element detach_joint_node = joint_node->clone();
238                     SimTK::Xml::element_iterator concreteJointNode =
239                         detach_joint_node.element_begin();
240                     detach_joint_node.removeNode(concreteJointNode);
241                     SimTK::Xml::element_iterator parentBodyElement =
242                         concreteJointNode->element_begin("parent_body");
243                     SimTK::String parent_name = "ground";
244                     parentBodyElement->getValueAs<SimTK::String>(parent_name);
245 
246                     // In version 30501, this Joint is 1 level deep (in the
247                     // model's JointSet), and the Bodies are necessarily 1
248                     // level deep (in the model's BodySet). Prepend "../" to
249                     // get the correct relative path.
250                     std::string parent_frame = parent_name;
251                     if (!parent_frame.empty())
252                         parent_frame = "../" + parent_frame;
253                     XMLDocument::addConnector(*concreteJointNode,
254                         "Connector_PhysicalFrame_", "parent_frame",
255                         parent_frame);
256                     std::string child_frame = body_name;
257                     if (!child_frame.empty())
258                         child_frame = "../" + child_frame;
259                     XMLDocument::addConnector(*concreteJointNode,
260                         "Connector_PhysicalFrame_", "child_frame",
261                         child_frame);
262                     concreteJointNode->eraseNode(parentBodyElement);
263                     jointObjects.insertNodeAfter(jointObjects.node_end(),
264                         *concreteJointNode);
265                     detach_joint_node.clearOrphan();
266                 }
267                 bodyIter->eraseNode(joint_node);
268             }
269         }
270         if (versionNumber < 30512) {
271             // FrameSet was removed from Model as of 30512 and this update
272             // is responsible for moving the Frames in the FrameSet to
273             // the Model's list of components.
274             SimTK::Xml::element_iterator componentsNode =
275                 aNode.element_begin("components");
276             SimTK::Xml::element_iterator framesNode =
277                 aNode.element_begin("FrameSet");
278 
279             // If no FrameSet nothing to be done
280             if (framesNode != aNode.element_end()) {
281                 if (componentsNode == aNode.element_end()) {
282                     // if no 'components' list element, create one and insert it
283                     SimTK::Xml::Element componentsElement("components");
284                     aNode.insertNodeBefore(framesNode, componentsElement);
285                 }
286 
287                 componentsNode = aNode.element_begin("components");
288 
289                 SimTK::Xml::element_iterator  objects_node =
290                     framesNode->element_begin("objects");
291 
292                 SimTK::Xml::element_iterator frameIter =
293                     objects_node->element_begin();
294                 for (; frameIter != objects_node->element_end(); ++frameIter) {
295                     SimTK::Xml::Element cloneOfFrame = frameIter->clone();
296                     componentsNode->insertNodeAfter(componentsNode->node_end(),
297                         cloneOfFrame);
298                 }
299                 // now delete the FrameSet
300                 framesNode->getParentElement().eraseNode(framesNode);
301             }
302         }
303     }
304      // Call base class now assuming _node has been corrected for current version
305      Super::updateFromXMLNode(aNode, versionNumber);
306      setDefaultProperties();
307 }
308 
309 
310 //=============================================================================
311 // CONSTRUCTION METHODS
312 //=============================================================================
313 //_____________________________________________________________________________
314 /**
315  * Set the values of all data members to an appropriate "null" value.
316  */
setNull()317 void Model::setNull()
318 {
319     _useVisualizer = false;
320     _allControllersEnabled = true;
321 
322     _validationLog="";
323 
324     _analysisSet.setMemoryOwner(false);
325 }
326 
constructProperties()327 void Model::constructProperties()
328 {
329     constructProperty_assembly_accuracy(1e-9);
330 
331     constructProperty_ground(Ground());
332 
333     constructProperty_gravity(SimTK::Vec3(0.0, -9.80665, 0.0));
334 
335     constructProperty_credits("Frank Anderson, Peter Loan, Ayman Habib, Ajay Seth, Michael Sherman");
336 
337     constructProperty_publications("List of publications related to model...");
338 
339     constructProperty_length_units("meters");
340     _lengthUnits = Units::Meters;
341 
342     constructProperty_force_units("N");
343     _forceUnits = Units::Newtons;
344 
345     BodySet bodies;
346     bodies.setName(IO::Lowercase(bodies.getConcreteClassName()));
347     constructProperty_BodySet(bodies);
348 
349     JointSet joints;
350     joints.setName(IO::Lowercase(joints.getConcreteClassName()));
351     constructProperty_JointSet(joints);
352 
353     ControllerSet controllers;
354     controllers.setName(IO::Lowercase(controllers.getConcreteClassName()));
355     constructProperty_ControllerSet(controllers);
356 
357     ConstraintSet constraints;
358     constraints.setName(IO::Lowercase(constraints.getConcreteClassName()));
359     constructProperty_ConstraintSet(constraints);
360 
361     ForceSet forces;
362     forces.setName(IO::Lowercase(forces.getConcreteClassName()));
363     constructProperty_ForceSet(forces);
364 
365     MarkerSet markers;
366     markers.setName(IO::Lowercase(markers.getConcreteClassName()));
367     constructProperty_MarkerSet(markers);
368 
369     ContactGeometrySet contacts;
370     contacts.setName(IO::Lowercase(contacts.getConcreteClassName()));
371     constructProperty_ContactGeometrySet(contacts);
372 
373     ProbeSet probes;
374     probes.setName(IO::Lowercase(probes.getConcreteClassName()));
375     constructProperty_ProbeSet(probes);
376 
377     ComponentSet miscComponents;
378     miscComponents.setName(IO::Lowercase(miscComponents.getConcreteClassName()));
379     constructProperty_ComponentSet(miscComponents);
380 
381     ModelVisualPreferences mvps;
382     mvps.setName(IO::Lowercase(mvps.getConcreteClassName()));
383     constructProperty_ModelVisualPreferences(mvps);
384 }
385 
386 // Append to the Model's validation log
appendToValidationLog(const std::string & note)387 void Model::appendToValidationLog(const std::string& note) {
388     _validationLog.append(note);
389 }
390 
391 //------------------------------------------------------------------------------
392 //                                BUILD SYSTEM
393 //------------------------------------------------------------------------------
394 // Perform some final checks on the Model, wire up all its components, and then
395 // build a computational System for it.
buildSystem()396 void Model::buildSystem() {
397     // Finish connecting up the Model.
398     setup();
399 
400     // Create the computational System representing this Model.
401     createMultibodySystem();
402 
403     // Create a Visualizer for this Model if one has been requested. This adds
404     // necessary elements to the System. Doesn't initialize geometry yet.
405     if (getUseVisualizer())
406         _modelViz.reset(new ModelVisualizer(*this));
407 }
408 
409 
410 //------------------------------------------------------------------------------
411 //                            INITIALIZE STATE
412 //------------------------------------------------------------------------------
413 // Requires that buildSystem() has already been called.
initializeState()414 SimTK::State& Model::initializeState() {
415     if (!hasSystem())
416         throw Exception("Model::initializeState(): call buildSystem() first.");
417 
418     // This tells Simbody to finalize the System.
419     getMultibodySystem().invalidateSystemTopologyCache();
420     getMultibodySystem().realizeTopology();
421 
422     // Set the model's operating state (internal member variable) to the
423     // default state that is stored inside the System.
424     _workingState = getMultibodySystem().getDefaultState();
425 
426     // Set the Simbody modeling option that tells any joints that use
427     // quaternions to use Euler angles instead.
428     _matter->setUseEulerAngles(_workingState, true);
429 
430     // Process the modified modeling option.
431     getMultibodySystem().realizeModel(_workingState);
432 
433     // Invoke the ModelComponent interface for initializing the state.
434     initStateFromProperties(_workingState);
435 
436     // Realize instance variables that may have been set above. This
437     // means floating point parameters such as mass properties and
438     // geometry placements are frozen.
439     getMultibodySystem().realize(_workingState, Stage::Instance);
440 
441     // Realize the initial configuration in preparation for assembly. This
442     // initial configuration does not necessarily satisfy constraints.
443     getMultibodySystem().realize(_workingState, Stage::Position);
444 
445     // Reset (initialize) all underlying Probe SimTK::Measures
446     for (int i=0; i<getProbeSet().getSize(); ++i)
447         getProbeSet().get(i).reset(_workingState);
448 
449     // Do the assembly
450     createAssemblySolver(_workingState);
451     assemble(_workingState);
452     // We can now collect up all the fixed geometry, which needs full configuration.
453     if (getUseVisualizer())
454         _modelViz->collectFixedGeometry(_workingState);
455 
456     return _workingState;
457 }
458 
459 
updWorkingState()460 SimTK::State& Model::updWorkingState()
461 {
462     if (!isValidSystem())
463         throw Exception("Model::updWorkingState(): call initializeState() first.");
464 
465     return _workingState;
466 }
467 
468 
getWorkingState() const469 const SimTK::State& Model::getWorkingState() const
470 {
471     if (!isValidSystem())
472         throw Exception("Model::getWorkingState(): call initializeState() first.");
473 
474     return _workingState;
475 }
476 
477 
assemble(SimTK::State & s,const Coordinate * coord,double weight)478 void Model::assemble(SimTK::State& s, const Coordinate *coord, double weight)
479 {
480 
481     bool constrained = false;
482     const CoordinateSet &coords = getCoordinateSet();
483     for(int i=0; i<coords.getSize(); ++i){
484         constrained = constrained || coords[i].isConstrained(s);
485     }
486 
487     // Don't bother assembling if the model has no constraints
488     if(get_ConstraintSet().getSize()< 1){
489         // just realize the current state to position
490         getMultibodySystem().realize(s, Stage::Position);
491 
492         // if a coordinate is locked or prescribed, then project will suffice
493         if(constrained){
494             // correct position constraint violations due to prescribed motion
495             getMultibodySystem().projectQ(s, 1e-10);
496 
497             // Have a new working configuration so should realize to velocity
498             getMultibodySystem().realize(s, Stage::Velocity);
499             // correct velocity constraint violations due to prescribed motion
500             getMultibodySystem().projectU(s, 1e-10);
501         }
502         return;
503     }
504 
505     if (!_assemblySolver){
506         createAssemblySolver(s);
507     }
508     const Array_<CoordinateReference>& coordRefs = _assemblySolver->getCoordinateReferences();
509 
510     for(unsigned int i=0; i<coordRefs.size(); i++){
511         const string &coordName = coordRefs[i].getName();
512         Coordinate& c = _coordinateSet.get(coordName);
513         _assemblySolver->updateCoordinateReference(coordName, c.getValue(s));
514     }
515 
516     if(coord) // use specified weighting for coordinate being set
517         _assemblySolver->updateCoordinateReference(coord->getName(), coord->getValue(s), weight);
518 
519 
520     try{
521         // Try to track first with model satisfying the constraints exactly.
522         _assemblySolver->track(s);
523     }
524     catch (const std::exception&)    {
525         try{
526             // Otherwise try to do a full-blown assemble
527             _assemblySolver->assemble(s);
528         }
529         catch (const std::exception& ex){
530             // Constraints are probably infeasible so try again relaxing constraints
531             cout << "Model unable to assemble: " << ex.what() << endl;
532             cout << "Model relaxing constraints and trying again." << endl;
533 
534             try{
535                 // Try to satisfy with constraints as errors weighted heavily.
536                 _assemblySolver->setConstraintWeight(20.0);
537                 _assemblySolver->assemble(s);
538             }
539             catch (const std::exception& ex){
540                 cout << "Model unable to assemble with relaxed constraints: " << ex.what() << endl;
541             }
542         }
543     }
544 
545     // Have a new working configuration so should realize to velocity
546     getMultibodySystem().realize(s, Stage::Velocity);
547 
548 }
549 
invalidateSystem()550 void Model::invalidateSystem()
551 {
552     if (_system)
553         _system->getSystemGuts().invalidateSystemTopologyCache();
554 }
555 
isValidSystem() const556 bool Model::isValidSystem() const
557 {
558     if (_system)
559         return _system->systemTopologyHasBeenRealized();
560     else
561         return false;
562 }
563 
564 
565 //_____________________________________________________________________________
566 /**
567  * Create the multibody system.
568  *
569  */
createMultibodySystem()570 void Model::createMultibodySystem()
571 {
572     // We must reset these unique_ptr's before deleting the System (through
573     // reset()), since deleting the System puts a null handle pointer inside
574     // the subsystems (since System deletes the subsystems).
575     _matter.reset();
576     _forceSubsystem.reset();
577     _contactSubsystem.reset();
578     // create system
579     _system.reset(new SimTK::MultibodySystem);
580     _matter.reset(new SimTK::SimbodyMatterSubsystem(*_system));
581     _forceSubsystem.reset(new SimTK::GeneralForceSubsystem(*_system));
582     _contactSubsystem.reset(new SimTK::GeneralContactSubsystem(*_system));
583 
584     // create gravity force, a direction is needed even if magnitude=0 for
585     // PotentialEnergy purposes.
586     double magnitude = get_gravity().norm();
587     SimTK::UnitVec3 direction = magnitude==0 ? SimTK::UnitVec3(0,-1,0) : SimTK::UnitVec3(get_gravity()/magnitude);
588     _gravityForce.reset(new SimTK::Force::Gravity(*_forceSubsystem, *_matter,
589                 direction, magnitude));
590 
591     addToSystem(*_system);
592 }
593 
594 
595 std::vector<SimTK::ReferencePtr<const Coordinate>>
getCoordinatesInMultibodyTreeOrder() const596     Model::getCoordinatesInMultibodyTreeOrder() const
597 {
598     OPENSIM_THROW_IF_FRMOBJ(!isValidSystem(), Exception,
599         "Cannot order Coordinates without a valid MultibodySystem.");
600 
601     int nc = getNumCoordinates();
602     auto coordinates = getComponentList<Coordinate>();
603 
604     std::vector<SimTK::ReferencePtr<const Coordinate>>
605         coordinatesInTreeOrder(nc,
606             SimTK::ReferencePtr<const Coordinate>(*coordinates.begin()));
607 
608     // We have a valid MultibodySystem underlying the Coordinates
609     const SimTK::State& s = getWorkingState();
610     SimTK_ASSERT_ALWAYS(nc <= s.getNQ(),
611         "Number of Coordinates exceeds the number of generalized coordinates in "
612         "the underlying MultibodySystem.");
613 
614     auto& matter = getSystem().getMatterSubsystem();
615 
616     int cnt = 0;
617 
618     for (auto& coord : coordinates) {
619         auto mbix = coord.getBodyIndex();
620         auto mqix = coord.getMobilizerQIndex();
621 
622         int cix = matter.getMobilizedBody(mbix).getFirstUIndex(s) + mqix;
623 
624         SimTK_ASSERT_ALWAYS(cix < nc, "Index exceeds the number of Coordinates "
625             "in this Model.");
626 
627         coordinatesInTreeOrder.at(cix).reset(&coord);
628         cnt++;
629     }
630 
631     SimTK_ASSERT_ALWAYS(cnt == nc,
632         "The number of ordered Coordinates does not correspond to the number of "
633         "Coordinates in the Model's CoordinateSet.");
634 
635     return coordinatesInTreeOrder;
636 }
637 
getWarningMesssageForMotionTypeInconsistency() const638 std::string Model::getWarningMesssageForMotionTypeInconsistency() const
639 {
640     std::string message;
641 
642     auto enumToString = [](Coordinate::MotionType mt)->std::string {
643         switch (mt) {
644         case Coordinate::MotionType::Rotational :
645             return "Rotational";
646         case Coordinate::MotionType::Translational :
647             return "Translational";
648         case Coordinate::MotionType::Coupled :
649             return "Coupled";
650         default:
651             return "Undefined";
652         }
653     };
654 
655     auto coordinates = getComponentList<Coordinate>();
656     for (auto& coord : coordinates) {
657         const Coordinate::MotionType oldMotionType =
658             coord.getUserSpecifiedMotionTypePriorTo40();
659         const Coordinate::MotionType motionType = coord.getMotionType();
660 
661         if( (oldMotionType != Coordinate::MotionType::Undefined ) &&
662             (oldMotionType != motionType) ){
663             message += "Coordinate '" + coord.getName() +
664                 "' was labeled as '" + enumToString(oldMotionType) +
665                 "' but was found to be '" + enumToString(motionType) + "' based on the joint definition.\n";
666         }
667     }
668 
669     // We have a reason to provide a warning. Add more details about the model
670     // and how to resolve future issues.
671     if (message.size()) {
672         message = "\nModel '" + getName() + "' has inconsistencies:\n" + message;
673         message +=
674             "You must update any motion files you generated in versions prior to 4.0. You can:\n"
675             "  (1) Run the updatePre40KinematicsFilesFor40MotionType() utility (in the scripting shell) OR\n"
676             "  (2) Re-run the Inverse Kinematics Tool with the updated model in 4.0.\n"
677             "In versions prior to 4.0, we allowed some Coupled coordinates to be incorrectly\n"
678             "labeled as Rotational. This leads to incorrect motion when playing back a pre-4.0\n"
679             "motion file (.mot or .sto in degrees) and incorrect inverse dynamics and\n"
680             "static optimization results.";
681     }
682 
683     return message;
684 }
685 
extendFinalizeFromProperties()686 void Model::extendFinalizeFromProperties()
687 {
688     Super::extendFinalizeFromProperties();
689 
690     // wipe-out the existing System
691     _matter.reset();
692     _forceSubsystem.reset();
693     _contactSubsystem.reset();
694     _system.reset();
695 
696     if(getForceSet().getSize()>0)
697     {
698         ForceSet &fs = updForceSet();
699         // Update internal subsets of the ForceSet
700         fs.updActuators();
701         fs.updMuscles();
702     }
703 
704     std::string warn = getWarningMesssageForMotionTypeInconsistency();
705     appendToValidationLog(warn);
706 
707     updCoordinateSet().populate(*this);
708 }
709 
createMultibodyTree()710 void Model::createMultibodyTree()
711 {
712     // building the system for the first time, need to tell
713     // multibodyTree builder what joints are available
714     _multibodyTree.clearGraph();
715     _multibodyTree.setWeldJointTypeName("WeldJoint");
716     _multibodyTree.setFreeJointTypeName("FreeJoint");
717 
718     ArrayPtrs<OpenSim::Joint> availablJointTypes;
719     Object::getRegisteredObjectsOfGivenType<OpenSim::Joint>(availablJointTypes);
720     for (int i = 0; i< availablJointTypes.getSize(); i++) {
721         OpenSim::Joint* jt = availablJointTypes[i];
722         if ((jt->getConcreteClassName() == "WeldJoint") ||
723             (jt->getConcreteClassName() == "FreeJoint")) {
724             continue;
725         }
726         else {
727             _multibodyTree.addJointType(
728                 availablJointTypes[i]->getConcreteClassName(),
729                 availablJointTypes[i]->numCoordinates(), false);
730         }
731     }
732 
733     Ground& ground = updGround();
734 
735     // assemble a multibody tree according to the PhysicalFrames in the
736     // OpenSim model, which include Ground and Bodies
737     _multibodyTree.addBody(ground.getAbsolutePathString(),
738                            0, false, &ground);
739 
740     auto bodies = getComponentList<Body>();
741     for (auto& body : bodies) {
742         _multibodyTree.addBody( body.getAbsolutePathString(),
743                                 body.getMass(), false,
744                                 const_cast<Body*>(&body) );
745     }
746 
747     std::vector<SimTK::ReferencePtr<Joint>> joints;
748     for (auto& joint : updComponentList<Joint>()) {
749         joints.push_back(SimTK::ReferencePtr<Joint>(&joint));
750     }
751 
752     // Complete multibody tree description by indicating how (mobilized)
753     // "bodies" are connected by joints.
754     for (auto& joint : joints) {
755         std::string name = joint->getAbsolutePathString();
756         IO::TrimWhitespace(name);
757 
758         // Currently we need to take a first pass at connecting the joints
759         // in order to ask the joint for the frames that they attach to and
760         // to determine their underlying base (physical) frames.
761         joint->finalizeConnections(*this);
762 
763         // Verify that the underlying frames are also connected so we can
764         // traverse to the base frame and get its name. This allows the
765         // (offset) frames to satisfy the sockets of Joint to be added
766         // to a Body, for example, and not just joint itself.
767         const PhysicalFrame& parent = joint->getParentFrame();
768         const PhysicalFrame& child = joint->getChildFrame();
769         const_cast<PhysicalFrame&>(parent).finalizeConnections(*this);
770         const_cast<PhysicalFrame&>(child).finalizeConnections(*this);
771 
772         OPENSIM_THROW_IF(&parent.findBaseFrame() == &child.findBaseFrame(),
773             JointFramesHaveSameBaseFrame, getName(),
774             parent.getName(), child.getName(), child.findBaseFrame().getName());
775 
776         // Use joints to define the underlying multibody tree
777         _multibodyTree.addJoint(name,
778             joint->getConcreteClassName(),
779             // Multibody tree builder only cares about bodies not intermediate
780             // frames that joints actually connect to.
781             parent.findBaseFrame().getAbsolutePathString(),
782             child.findBaseFrame().getAbsolutePathString(),
783             false,
784             joint.get());
785     }
786 }
787 
extendConnectToModel(Model & model)788 void Model::extendConnectToModel(Model &model)
789 {
790     Super::extendConnectToModel(model);
791 
792     if (&model != this){
793         cout << "Model::" << getName() <<
794             " is being connected to model " <<
795             model.getName() << "." << endl;
796         // if part of another Model, that Model is in charge
797         // of creating a valid Multibody tree that includes
798         // Components of this Model.
799         return;
800     }
801 
802     // Create the Multibody tree according to the components that
803     // form this model.
804     createMultibodyTree();
805 
806     // generate the graph of the Multibody tree to determine the order in
807     // which subcomponents will be added to the MultibodySystem (in addToSystem)
808     _multibodyTree.generateGraph();
809     //_multibodyTree.dumpGraph(cout);
810     //cout << endl;
811 
812     // Ground is the first component of any MultibodySystem
813     Ground& ground = updGround();
814     setNextSubcomponentInSystem(ground);
815 
816     // The JointSet of the Model is only being manipulated for consistency with
817     // Tool expectations. TODO fix Tools and remove
818     JointSet& joints = upd_JointSet();
819 
820     bool isMemoryOwner = joints.getMemoryOwner();
821     //Temporarily set owner ship to false so we
822     //can swap to rearrange order of the joints
823     joints.setMemoryOwner(false);
824 
825     // Run through all the mobilizers in the multibody tree, adding
826     // a joint in the correct sequence. Also add massless bodies,
827     // loop closure constraints, etc... to form the valid tree.
828     for (int m = 0; m < _multibodyTree.getNumMobilizers(); ++m) {
829         // Get a mobilizer from the tree, then extract its corresponding
830         // joint and bodies. Note that these should have equivalents in OpenSim.
831         const MultibodyGraphMaker::Mobilizer& mob
832             = _multibodyTree.getMobilizer(m);
833 
834         if (mob.isSlaveMobilizer()){
835             // add the slave body and joint
836             Body* outbMaster = static_cast<Body*>(mob.getOutboardMasterBodyRef());
837             Joint* useJoint = static_cast<Joint*>(mob.getJointRef());
838             Body* outb = static_cast<Body*>(mob.getOutboardBodyRef());
839 
840             // the joint must be added to the system next
841             setNextSubcomponentInSystem(*useJoint);
842 
843             if (!outb) {
844                 outb = outbMaster->addSlave();
845                 useJoint->setSlaveBodyForChild(*outb);
846                 SimTK::Transform o(SimTK::Vec3(0));
847                 //Now add the constraints that weld the slave to the master at the
848                 // body origin
849                 std::string pathName = outb->getAbsolutePathString();
850                 WeldConstraint* weld = new WeldConstraint(outb->getName()+"_weld",
851                                                           *outbMaster, o, *outb, o);
852 
853                 // include within adopted list of owned components
854                 adoptSubcomponent(weld);
855                 setNextSubcomponentInSystem(*weld);
856             }
857         }
858 
859         if (mob.isAddedBaseMobilizer()){
860             // create and add the base joint to enable these dofs
861             Body* child = static_cast<Body*>(mob.getOutboardBodyRef());
862             cout << "Body '" << child->getName() << "' not connected by a Joint.\n"
863                 << "A FreeJoint will be added to connect it to ground." << endl;
864             Ground* ground = static_cast<Ground*>(mob.getInboardBodyRef());
865 
866             std::string jname = "free_" + child->getName();
867             SimTK::Vec3 zeroVec(0.0);
868             Joint* free = new FreeJoint(jname, *ground, *child);
869             free->isReversed = mob.isReversedFromJoint();
870             // TODO: Joints are currently required to be in the JointSet
871             // When the reordering of Joints is eliminated (see following else block)
872             // this limitation can be removed and the free joint adopted as in
873             // internal subcomponent (similar to the weld constraint above)
874             addJoint(free);
875             setNextSubcomponentInSystem(*free);
876         }
877         else{
878             // Update the directionality of the joint according to tree's
879             // preferential direction
880             static_cast<Joint*>(mob.getJointRef())->isReversed =
881                 mob.isReversedFromJoint();
882 
883             // order the joint components in the order of the multibody tree
884             Joint* joint = static_cast<Joint*>(mob.getJointRef());
885             setNextSubcomponentInSystem(*joint);
886 
887 
888         }
889     }
890     joints.setMemoryOwner(isMemoryOwner);
891 
892     // Add the loop joints if any.
893     for (int lcx = 0; lcx < _multibodyTree.getNumLoopConstraints(); ++lcx) {
894         const MultibodyGraphMaker::LoopConstraint& loop =
895             _multibodyTree.getLoopConstraint(lcx);
896 
897         Joint& joint = *(Joint*)loop.getJointRef();
898         Body&  parent = *(Body*)loop.getParentBodyRef();
899         Body&  child = *(Body*)loop.getChildBodyRef();
900 
901         if (joint.getConcreteClassName() == "WeldJoint") {
902             WeldConstraint* weld = new WeldConstraint( joint.getName()+"_Loop",
903                 parent, joint.getParentFrame().findTransformInBaseFrame(),
904                 child, joint.getChildFrame().findTransformInBaseFrame());
905             adoptSubcomponent(weld);
906             setNextSubcomponentInSystem(*weld);
907         }
908         else if (joint.getConcreteClassName() == "BallJoint") {
909             PointConstraint* point = new PointConstraint(
910                 parent, joint.getParentFrame().findTransformInBaseFrame().p(),
911                 child, joint.getChildFrame().findTransformInBaseFrame().p());
912             point->setName(joint.getName() + "_Loop");
913             adoptSubcomponent(point);
914             setNextSubcomponentInSystem(*point);
915         }
916         else if (joint.getConcreteClassName() == "FreeJoint") {
917             // A "free" loop constraint is no constraint at all so we can
918             // just ignore it. It might be more convenient if there were
919             // a 0-constraint Constraint::Free, just as there is a 0-mobility
920             // MobilizedBody::Weld.
921         }
922         else
923             throw std::runtime_error(
924             "Unrecognized loop constraint type '" + joint.getConcreteClassName() + "'.");
925     }
926 
927     // Now include all remaining Components beginning with PhysicalOffsetFrames
928     // since Constraints, Forces and other components can only be applied to
929     // PhyicalFrames, which include PhysicalOffsetFrames.
930     // PhysicalOffsetFrames require that their parent frame (a PhysicalFrame)
931     // be added to the System first. So for each PhysicalOffsetFrame locate its
932     // parent and verify its presence in the _orderedList otherwise add it first.
933     auto poFrames = updComponentList<PhysicalOffsetFrame>();
934     for (auto& pof : poFrames) {
935         // Ground and Body type PhysicalFrames are included in the Multibody graph
936         // PhysicalOffsetFrame can be listed in any order and may be attached
937         // to any other PhysicalOffsetFrame, so we need to find their parent(s)
938         // in the tree and add them first.
939         pof.finalizeConnections(*this);
940         const PhysicalOffsetFrame* parentPof =
941             dynamic_cast<const PhysicalOffsetFrame*>(&pof.getParentFrame());
942         std::vector<const PhysicalOffsetFrame*> parentPofs;
943         while (parentPof) {
944             const auto found =
945                 std::find(parentPofs.begin(), parentPofs.end(), parentPof);
946             OPENSIM_THROW_IF_FRMOBJ(found != parentPofs.end(),
947                 PhysicalOffsetFramesFormLoop, (*found)->getName());
948             parentPofs.push_back(parentPof);
949             // Given a chain of offsets, the most proximal must have Ground or
950             // Body as its parent. When that happens we can stop.
951             parentPof =
952                 dynamic_cast<const PhysicalOffsetFrame*>(&parentPof->getParentFrame());
953         }
954         while (parentPofs.size()) {
955             setNextSubcomponentInSystem(*parentPofs.back());
956             parentPofs.pop_back();
957         }
958         setNextSubcomponentInSystem(pof);
959     }
960 
961     // and everything else including Forces, Controllers, etc...
962     // belonging to the model
963     auto components = getImmediateSubcomponents();
964     for (const auto& comp : components) {
965         setNextSubcomponentInSystem(comp.getRef());
966     }
967 
968     // Populate the model's list of Coordinates as references into
969     // the individual Joints.
970     updCoordinateSet().populate(*this);
971     updForceSet().setupGroups();
972     updControllerSet().setActuators(updActuators());
973 
974     // TODO: Get rid of the SimbodyEngine
975     updSimbodyEngine().connectSimbodyEngineToModel(*this);
976 }
977 
978 
979 // ModelComponent interface enables this model to be a subcomponent of another
980 // model. In that case, it adds itself to the parent model's system.
extendAddToSystem(SimTK::MultibodySystem & system) const981 void Model::extendAddToSystem(SimTK::MultibodySystem& system) const
982 {
983     Super::extendAddToSystem(system);
984 
985     Model *mutableThis = const_cast<Model *>(this);
986 
987     //Analyses are not Components so add them after legit
988     //Components have been wired-up correctly.
989     mutableThis->updAnalysisSet().setModel(*mutableThis);
990 
991     // Reset the vector of all controls' defaults
992     mutableThis->_defaultControls.resize(0);
993 
994     // Create the shared cache that will hold all model controls
995     // This must be created before Actuator.extendAddToSystem() since Actuator
996     // will append its "slots" and retain its index by accessing this cached Vector.
997     // Value depends on velocity and invalidates dynamics BUT should not trigger
998     // re-computation of the controls which are necessary for dynamics
999     Measure_<Vector>::Result modelControls(system.updDefaultSubsystem(),
1000         Stage::Velocity, Stage::Acceleration);
1001 
1002     mutableThis->_modelControlsIndex = modelControls.getSubsystemMeasureIndex();
1003 }
1004 
1005 
1006 // Add any Component derived from ModelComponent to the Model
addModelComponent(ModelComponent * component)1007 void Model::addModelComponent(ModelComponent* component)
1008 {
1009     if(component){
1010         upd_ComponentSet().adoptAndAppend(component);
1011         finalizeFromProperties();
1012         prependComponentPathToConnecteePath(*component);
1013     }
1014 }
1015 
1016 // Add a body to the Model
addBody(OpenSim::Body * body)1017 void Model::addBody(OpenSim::Body* body)
1018 {
1019     if (body){
1020         updBodySet().adoptAndAppend(body);
1021         finalizeFromProperties();
1022         prependComponentPathToConnecteePath(*body);
1023     }
1024 }
1025 
1026 // Add a Marker to the Model
addMarker(OpenSim::Marker * marker)1027 void Model::addMarker(OpenSim::Marker* marker)
1028 {
1029     if (marker){
1030         updMarkerSet().adoptAndAppend(marker);
1031         finalizeFromProperties();
1032         prependComponentPathToConnecteePath(*marker);
1033     }
1034 }
1035 
1036 // Add a joint to the Model
addJoint(Joint * joint)1037 void Model::addJoint(Joint* joint)
1038 {
1039     if (joint){
1040         updJointSet().adoptAndAppend(joint);
1041         finalizeFromProperties();
1042         updCoordinateSet().populate(*this);
1043         prependComponentPathToConnecteePath(*joint);
1044     }
1045 }
1046 
1047 // Add a constraint to the Model
addConstraint(OpenSim::Constraint * constraint)1048 void Model::addConstraint(OpenSim::Constraint *constraint)
1049 {
1050     if(constraint){
1051         updConstraintSet().adoptAndAppend(constraint);
1052         finalizeFromProperties();
1053         prependComponentPathToConnecteePath(*constraint);
1054     }
1055 }
1056 
1057 // Add a force to the Model
addForce(OpenSim::Force * force)1058 void Model::addForce(OpenSim::Force *force)
1059 {
1060     if(force){
1061         updForceSet().adoptAndAppend(force);
1062         finalizeFromProperties();
1063         prependComponentPathToConnecteePath(*force);
1064     }
1065 }
1066 
1067 // Add a probe to the Model
addProbe(OpenSim::Probe * probe)1068 void Model::addProbe(OpenSim::Probe *probe)
1069 {
1070     if(probe){
1071         updProbeSet().adoptAndAppend(probe);
1072         finalizeFromProperties();
1073         prependComponentPathToConnecteePath(*probe);
1074     }
1075 }
1076 
1077 // Remove a probe from the Model. Probe will be deleted since model owns it
removeProbe(OpenSim::Probe * aProbe)1078 void Model::removeProbe(OpenSim::Probe *aProbe)
1079 {
1080     clearConnections();
1081     updProbeSet().remove(aProbe);
1082     finalizeFromProperties();
1083 }
1084 
1085 // Add a contact geometry to the Model
addContactGeometry(OpenSim::ContactGeometry * contactGeometry)1086 void Model::addContactGeometry(OpenSim::ContactGeometry *contactGeometry)
1087 {
1088     if (contactGeometry) {
1089         updContactGeometrySet().adoptAndAppend(contactGeometry);
1090         finalizeFromProperties();
1091         prependComponentPathToConnecteePath(*contactGeometry);
1092     }
1093 }
1094 
1095 // Add a controller to the Model
addController(Controller * controller)1096 void Model::addController(Controller *controller)
1097 {
1098     if (controller) {
1099         updControllerSet().adoptAndAppend(controller);
1100         finalizeFromProperties();
1101         prependComponentPathToConnecteePath(*controller);
1102     }
1103 }
1104 //_____________________________________________________________________________
1105 /* Perform some setup functions that happen after the
1106  * object has been deserialized. This method is
1107  * not yet designed to be called after a model has been
1108  * copied.
1109  */
setup()1110 void Model::setup()
1111 {
1112     // finalize the model and its subcomponents from its properties
1113     // automatically marks properties that are Components as subcomponents
1114     finalizeFromProperties();
1115     //now connect the Model and all its subcomponents all up
1116     finalizeConnections();
1117 }
1118 
1119 //_____________________________________________________________________________
1120 /* Perform some clean up functions that are normally done from the destructor
1121  * however this gives the GUI a way to proactively do the cleaning without
1122  * waiting for garbage collection to do the actual cleanup.
1123  */
cleanup()1124 void Model::cleanup()
1125 {
1126     clearConnections();
1127     upd_ForceSet().setSize(0);
1128 }
1129 
setDefaultProperties()1130 void Model::setDefaultProperties()
1131 {
1132     // Initialize the length and force units from the strings specified in the model file.
1133     // If they were not specified, use meters and Newtons.
1134     _lengthUnits = Units(get_length_units());
1135     _forceUnits = Units(get_force_units());
1136 }
1137 
extendInitStateFromProperties(SimTK::State & state) const1138 void Model::extendInitStateFromProperties(SimTK::State& state) const
1139 {
1140     Super::extendInitStateFromProperties(state);
1141     // Allocate the size and default values for controls
1142     // Actuators will have a const view into the cache
1143     Measure_<Vector>::Result controlsCache =
1144         Measure_<Vector>::Result::getAs(updSystem().updDefaultSubsystem().getMeasure(_modelControlsIndex));
1145     controlsCache.updValue(state).resize(_defaultControls.size());
1146     controlsCache.updValue(state) = _defaultControls;
1147 }
1148 
extendSetPropertiesFromState(const SimTK::State & state)1149 void Model::extendSetPropertiesFromState(const SimTK::State& state)
1150 {
1151     Super::extendSetPropertiesFromState(state);
1152 }
1153 
generateDecorations(bool fixed,const ModelDisplayHints & hints,const SimTK::State & state,SimTK::Array_<SimTK::DecorativeGeometry> & appendToThis) const1154 void Model::generateDecorations
1155        (bool                                        fixed,
1156         const ModelDisplayHints&                    hints,
1157         const SimTK::State&                         state,
1158         SimTK::Array_<SimTK::DecorativeGeometry>&   appendToThis) const
1159 {
1160     ComponentList<const Component> allComps = getComponentList();
1161     ComponentList<Component>::const_iterator iter = allComps.begin();
1162     while (iter != allComps.end()){
1163         //std::string cn = iter->getConcreteClassName();
1164         //std::cout << cn << ":" << iter->getName() << std::endl;
1165         iter->generateDecorations(fixed, hints, state, appendToThis);
1166         iter++;
1167     }
1168 }
1169 
equilibrateMuscles(SimTK::State & state)1170 void Model::equilibrateMuscles(SimTK::State& state)
1171 {
1172     getMultibodySystem().realize(state, Stage::Velocity);
1173 
1174     bool failed = false;
1175     string errorMsg = "";
1176 
1177     auto muscles = getComponentList<Muscle>();
1178 
1179     for (auto& muscle : muscles) {
1180         if (muscle.appliesForce(state)){
1181             try{
1182                 muscle.computeEquilibrium(state);
1183             }
1184             catch (const std::exception& e) {
1185                 if(!failed){ // haven't failed to equilibrate other muscles yet
1186                     errorMsg = e.what();
1187                     failed = true;
1188                 }
1189                 // just because one muscle failed to equilibrate doesn't mean
1190                 // it isn't still useful to have remaining muscles equilibrate
1191                 // in an analysis, for example, we might not be reporting about
1192                 // all muscles, so continue with the rest.
1193                 continue;
1194             }
1195         }
1196     }
1197 
1198     if(failed) // Notify the caller of the failure to equilibrate
1199         throw Exception("Model::equilibrateMuscles() "+errorMsg, __FILE__, __LINE__);
1200 }
1201 
1202 //=============================================================================
1203 // GRAVITY
1204 //=============================================================================
1205 //_____________________________________________________________________________
1206 /**
1207  * Get the gravity vector in the global frame.
1208  *
1209  * @return the XYZ gravity vector in the global frame is returned here.
1210  */
getGravity() const1211 SimTK::Vec3 Model::getGravity() const
1212 {
1213     if(_gravityForce)
1214         return _gravityForce->getDefaultGravityVector();
1215 
1216     return get_gravity();
1217 }
1218 //_____________________________________________________________________________
1219 /**
1220  * Set the gravity vector in the global frame.
1221  *
1222  * @param aGrav the XYZ gravity vector
1223  * @return Whether or not the gravity vector was successfully set.
1224  */
setGravity(const SimTK::Vec3 & aGrav)1225 bool Model::setGravity(const SimTK::Vec3& aGrav)
1226 {
1227     upd_gravity() = aGrav;
1228 
1229     if(_gravityForce)
1230         _gravityForce->setDefaultGravityVector(aGrav);
1231 
1232     return true;
1233 }
1234 
1235 
1236 //=============================================================================
1237 // NUMBERS
1238 //=============================================================================
1239 //_____________________________________________________________________________
1240 /**
1241  * Get the number of markers in the model.
1242  *
1243  * @return Number of markers.
1244  */
getNumMarkers() const1245 int Model::getNumMarkers() const
1246 {
1247     return get_MarkerSet().getSize();
1248 }
1249 /**
1250  * Get the number of ContactGeometry objects in the model.
1251  *
1252  * @return Number of ContactGeometries.
1253  */
getNumContactGeometries() const1254 int Model::getNumContactGeometries() const
1255 {
1256     return get_ContactGeometrySet().getSize();
1257 }
1258 
1259 /**
1260  * Get the number of Muscle state variables in the model.
1261  *
1262  * @return Number of MuscleStates.
1263  */
1264 
getNumMuscleStates() const1265 int Model::getNumMuscleStates() const {
1266 
1267     int n = 0;
1268     for(int i=0;i<get_ForceSet().getSize();i++){
1269         Muscle *mus = dynamic_cast<Muscle*>( &get_ForceSet().get(i) );
1270         if(mus!=NULL) {
1271             n += mus->getNumStateVariables();
1272         }
1273     }
1274     return(n);
1275 }
1276 
1277 /**
1278  * Get the number of Probe state variables in the model.
1279  *
1280  * @return Number of MuscleStates.
1281  */
1282 
getNumProbeStates() const1283 int Model::getNumProbeStates() const {
1284 
1285     int n = 0;
1286     for(int i=0;i<get_ProbeSet().getSize();i++){
1287         Probe *p = dynamic_cast<Probe*>( &get_ProbeSet().get(i) );
1288         if(p!=NULL) {
1289             n += p->getNumInternalMeasureStates();
1290         }
1291     }
1292     return(n);
1293 }
1294 
1295 //_____________________________________________________________________________
1296 /**
1297  * Get the total number of bodies in the model.
1298  *
1299  * @return Number of bodies.
1300  */
getNumBodies() const1301 int Model::getNumBodies() const
1302 {
1303     return  getBodySet().getSize();
1304 }
1305 
1306 //_____________________________________________________________________________
1307 /**
1308  * Get the total number of joints in the model.
1309  *
1310  * @return Number of joints.
1311  */
getNumJoints() const1312 int Model::getNumJoints() const
1313 {
1314     return  getJointSet().getSize();
1315 }
1316 //_____________________________________________________________________________
1317 /**
1318  * Get the total number of coordinates in the model.
1319  *
1320  * @return Number of coordinates.
1321  */
getNumCoordinates() const1322 int Model::getNumCoordinates() const
1323 {
1324     return _coordinateSet.getSize();
1325 }
1326 
1327 /**
1328  * Get the total number of coordinates = number of speeds in the model.
1329  *
1330  * @return Number of coordinates.
1331  */
getNumSpeeds() const1332 int Model::getNumSpeeds() const
1333 {
1334     return _coordinateSet.getSize();
1335 }
1336 
1337 /**
1338 * Get the total number of constraints in the model.
1339 * @return Number of constraints.
1340 */
getNumConstraints() const1341 int Model::getNumConstraints() const
1342 {
1343     return get_ConstraintSet().getSize();
1344 }
1345 
1346 
1347 /**
1348  * Get the total number of probes in the model.
1349  *
1350  * @return Number of probes.
1351  */
getNumProbes() const1352 int Model::getNumProbes() const
1353 {
1354     return get_ProbeSet().getSize();
1355 }
1356 
1357 //_____________________________________________________________________________
1358 /**
1359  * Get the subset of Forces in the model which are actuators
1360  *
1361  * @return The set of Actuators
1362  */
getActuators() const1363 const Set<Actuator>& Model::getActuators() const
1364 {
1365     return get_ForceSet().getActuators();
1366 }
updActuators()1367 Set<Actuator>& Model::updActuators()
1368 {
1369     return upd_ForceSet().updActuators();
1370 }
1371 
1372 //_____________________________________________________________________________
1373 /**
1374  * Get the subset of Forces in the model which are muscles
1375  *
1376  * @return The set of Muscles
1377  */
getMuscles() const1378 const Set<Muscle>& Model::getMuscles() const
1379 {
1380     return get_ForceSet().getMuscles();
1381 }
updMuscles()1382 Set<Muscle>& Model::updMuscles()
1383 {
1384     return upd_ForceSet().updMuscles();
1385 }
1386 
1387 //_____________________________________________________________________________
1388 /**
1389  * Get the number of analyses in the model.
1390  *
1391  * @return The number of analyses
1392  */
getNumAnalyses() const1393 int Model::getNumAnalyses() const
1394 {
1395     return _analysisSet.getSize();
1396 }
1397 
1398 //_____________________________________________________________________________
1399 
1400 //=============================================================================
1401 // TIME NORMALIZATION CONSTANT
1402 //=============================================================================
1403 //_____________________________________________________________________________
1404 
1405 
1406 //=============================================================================
1407 // STATES
1408 //=============================================================================
1409 //_____________________________________________________________________________
1410 
1411 
1412 
1413 
1414 
1415 //_____________________________________________________________________________
1416 /**
1417  * Add an analysis to the model.
1418  *
1419  * @param aAnalysis pointer to the analysis to add
1420  */
addAnalysis(Analysis * aAnalysis)1421 void Model::addAnalysis(Analysis *aAnalysis)
1422 {
1423     if (aAnalysis )
1424     {
1425 //      aAnalysis->setModel(this);
1426         _analysisSet.adoptAndAppend(aAnalysis);
1427     }
1428 }
1429 //_____________________________________________________________________________
1430 /**
1431  * Remove an analysis from the model
1432  *
1433  * @param aAnalysis Pointer to the analysis to remove.
1434  * If deleteIt is true (default) the Analysis object itself is destroyed
1435  * else only removed from the list which is the desired behavior when the Analysis
1436  * is created from the GUI.
1437  */
removeAnalysis(Analysis * aAnalysis,bool deleteIt)1438 void Model::removeAnalysis(Analysis *aAnalysis, bool deleteIt)
1439 {
1440     // CHECK FOR NULL
1441     if(aAnalysis==NULL) {
1442         cout << "Model.removeAnalysis:  ERROR- NULL analysis.\n" << endl;
1443     }
1444     if (!deleteIt){
1445         bool saveStatus = _analysisSet.getMemoryOwner();
1446         _analysisSet.setMemoryOwner(false);
1447         _analysisSet.remove(aAnalysis);
1448         _analysisSet.setMemoryOwner(saveStatus);
1449     }
1450     else
1451         _analysisSet.remove(aAnalysis);
1452 }
1453 
1454 //_____________________________________________________________________________
1455 /**
1456  * Remove a controller from the model
1457  *
1458  * @param aController Pointer to the controller to remove.
1459  */
removeController(Controller * aController)1460 void Model::removeController(Controller *aController)
1461 {
1462     // CHECK FOR NULL
1463     if(aController==NULL) {
1464         cout << "Model.removeController:  ERROR- NULL controller.\n" << endl;
1465     }
1466 
1467     upd_ControllerSet().remove(aController);
1468 }
1469 
1470 
1471 
1472 //==============================================================================
1473 // OPERATIONS
1474 //==============================================================================
1475 
1476 //------------------------------------------------------------------------------
1477 // SCALE
1478 //------------------------------------------------------------------------------
scale(SimTK::State & s,const ScaleSet & scaleSet,bool preserveMassDist,double finalMass)1479 bool Model::scale(SimTK::State& s, const ScaleSet& scaleSet,
1480                   bool preserveMassDist, double finalMass)
1481 {
1482     // Save the model's current pose.
1483     SimTK::Vector savedConfiguration = s.getY();
1484 
1485     // Put the model in a default pose so that GeometryPath lengths can be
1486     // computed and stored. These lengths will be required for adjusting
1487     // properties after the rest of the model has been scaled.
1488     applyDefaultConfiguration(s);
1489 
1490     // Call preScale() on each ModelComponent owned by the model to store
1491     // GeometryPath lengths (and perform any other necessary computations).
1492     for (ModelComponent& comp : updComponentList<ModelComponent>())
1493         comp.preScale(s, scaleSet);
1494 
1495     // Call scale() on each ModelComponent owned by the model. Each
1496     // ModelComponent is responsible for scaling itself. All scaling operations
1497     // are performed here except scaling inertial properties of bodies, which is
1498     // done below.
1499     for (ModelComponent& comp : updComponentList<ModelComponent>())
1500         comp.scale(s, scaleSet);
1501 
1502     // Scale the inertial properties of bodies. If "preserve mass distribution"
1503     // is true, then the masses are not scaled (but inertias are still updated).
1504     for (Body& body : updComponentList<Body>())
1505         body.scaleInertialProperties(scaleSet, !preserveMassDist);
1506 
1507     // When bodies are scaled, the properties of the model are changed. The
1508     // general rule is that you MUST recreate and initialize the system when
1509     // properties of the model change. We must do that here or we will be
1510     // querying a stale system (e.g., wrong body properties!).
1511     s = initSystem();
1512 
1513     // Now that the masses of the individual bodies have been scaled (if
1514     // preserveMassDist == false), get the total mass and compare it to
1515     // finalMass in order to determine how much to scale the body masses again,
1516     // so that the total model mass comes out to finalMass.
1517     if (finalMass > 0.0)
1518     {
1519         const double mass = getTotalMass(s);
1520         if (mass > 0.0)
1521         {
1522             const double factor = finalMass / mass;
1523             for (Body& body : updComponentList<Body>())
1524                 body.scaleMass(factor);
1525 
1526             // Recreate the system and update the state after updating masses.
1527             s = initSystem();
1528 
1529             // Ensure the final model mass is correct.
1530             const double newMass = getTotalMass(s);
1531             const double normDiffMass = abs(finalMass - newMass) / finalMass;
1532             if (normDiffMass > SimTK::SignificantReal) {
1533                 throw Exception("Model::scale() scaled model mass does not match specified subject mass.");
1534             }
1535         }
1536     }
1537 
1538     // Call postScale() on all ModelComponents owned by the model so that
1539     // components like muscles, ligaments, and path springs can update their
1540     // properties based on their new path length.
1541     for (ModelComponent& comp : updComponentList<ModelComponent>())
1542         comp.postScale(s, scaleSet);
1543 
1544     // Changed the model after scaling path actuators. Have to recreate system!
1545     s = initSystem();
1546 
1547     // Put the model back in its original pose.
1548     s.updY() = savedConfiguration;
1549     getMultibodySystem().realize( s, SimTK::Stage::Velocity );
1550 
1551     return true;
1552 }
1553 
1554 
1555 //=============================================================================
1556 // PRINT
1557 //=============================================================================
printBasicInfo(std::ostream & aOStream) const1558 void Model::printBasicInfo(std::ostream& aOStream) const
1559 {
1560     OPENSIM_THROW_IF_FRMOBJ(!isObjectUpToDateWithProperties(), Exception,
1561         "Model::finalizeFromProperties() must be called first.");
1562 
1563     aOStream
1564         << "\n               MODEL: " << getName()
1565         << "\n         coordinates: " << countNumComponents<Coordinate>()
1566         << "\n              forces: " << countNumComponents<Force>()
1567         << "\n           actuators: " << countNumComponents<Actuator>()
1568         << "\n             muscles: " << countNumComponents<Muscle>()
1569         << "\n            analyses: " << getNumAnalyses()
1570         << "\n              probes: " << countNumComponents<Probe>()
1571         << "\n              bodies: " << countNumComponents<Body>()
1572         << "\n              joints: " << countNumComponents<Joint>()
1573         << "\n         constraints: " << countNumComponents<Constraint>()
1574         << "\n             markers: " << countNumComponents<Marker>()
1575         << "\n         controllers: " << countNumComponents<Controller>()
1576         << "\n  contact geometries: " << countNumComponents<ContactGeometry>()
1577         << "\nmisc modelcomponents: " << getMiscModelComponentSet().getSize()
1578         << std::endl;
1579 }
1580 
printDetailedInfo(const SimTK::State & s,std::ostream & aOStream) const1581 void Model::printDetailedInfo(const SimTK::State& s, std::ostream& aOStream) const
1582 {
1583     aOStream << "MODEL: " << getName() << std::endl;
1584 
1585     aOStream
1586         << "\nnumStates = "      << s.getNY()
1587         << "\nnumCoordinates = " << countNumComponents<Coordinate>()
1588         << "\nnumSpeeds = "      << countNumComponents<Coordinate>()
1589         << "\nnumActuators = "   << countNumComponents<Actuator>()
1590         << "\nnumBodies = "      << countNumComponents<Body>()
1591         << "\nnumConstraints = " << countNumComponents<Constraint>()
1592         << "\nnumProbes = "      << countNumComponents<Probe>()
1593         << std::endl;
1594 
1595     aOStream << "\nANALYSES (total: " << getNumAnalyses() << ")" << std::endl;
1596     for (int i = 0; i < _analysisSet.getSize(); ++i)
1597         aOStream << "analysis[" << i << "] = " << _analysisSet.get(i).getName()
1598                  << std::endl;
1599 
1600     aOStream << "\nBODIES (total: " << countNumComponents<Body>() << ")"
1601              << std::endl;
1602     unsigned bodyNum = 0u;
1603     for (const auto& body : getComponentList<Body>()) {
1604         const auto inertia = body.getInertia();
1605         aOStream << "body[" << bodyNum << "] = " << body.getName() << ". "
1606             << "mass: " << body.get_mass()
1607             << "\n              moments of inertia:  " << inertia.getMoments()
1608             << "\n              products of inertia: " << inertia.getProducts()
1609             << std::endl;
1610         ++bodyNum;
1611     }
1612 
1613     aOStream << "\nJOINTS (total: " << countNumComponents<Joint>() << ")"
1614              << std::endl;
1615     unsigned jointNum = 0u;
1616     for (const auto& joint : getComponentList<Joint>()) {
1617         aOStream << "joint[" << jointNum << "] = " << joint.getName() << "."
1618                  << " parent: " << joint.getParentFrame().getName()
1619                  << ", child: " << joint.getChildFrame().getName() << std::endl;
1620         ++jointNum;
1621     }
1622 
1623     aOStream << "\nACTUATORS (total: " << countNumComponents<Actuator>() << ")"
1624              << std::endl;
1625     unsigned actuatorNum = 0u;
1626     for (const auto& actuator : getComponentList<Actuator>()) {
1627         aOStream << "actuator[" << actuatorNum << "] = " << actuator.getName()
1628                  << std::endl;
1629         ++actuatorNum;
1630     }
1631 
1632     /*
1633     int n;
1634     aOStream<<"MODEL: "<<getName()<<std::endl;
1635 
1636     n = getNumBodies();
1637     aOStream<<"\nBODIES ("<<n<<")" << std::endl;
1638     for(i=0;i<n;i++) aOStream<<"body["<<i<<"] = "<<getFrameName(i)<<std::endl;
1639 
1640     n = getNQ();
1641     aOStream<<"\nGENERALIZED COORDINATES ("<<n<<")" << std::endl;
1642     for(i=0;i<n;i++) aOStream<<"q["<<i<<"] = "<<getCoordinateName(i)<<std::endl;
1643 
1644     n = getNU();
1645     aOStream<<"\nGENERALIZED SPEEDS ("<<n<<")" << std::endl;
1646     for(i=0;i<n;i++) aOStream<<"u["<<i<<"] = "<<getSpeedName(i)<<std::endl;
1647 
1648     n = getNA();
1649     aOStream<<"\nACTUATORS ("<<n<<")" << std::endl;
1650     for(i=0;i<n;i++) aOStream<<"actuator["<<i<<"] = "<<getActuatorName(i)<<std::endl;
1651 
1652     n = getNP();
1653     aOStream<<"\nCONTACTS ("<<n<<")" << std::endl;
1654     */
1655 
1656     Array<string> stateNames = getStateVariableNames();
1657     aOStream << "\nSTATES (total: " << stateNames.getSize() << ")" << std::endl;
1658     for (int i = 0; i < stateNames.getSize(); ++i)
1659         aOStream << "y[" << i << "] = " << stateNames[i] << std::endl;
1660 }
1661 
1662 //--------------------------------------------------------------------------
1663 // CONFIGURATION
1664 //--------------------------------------------------------------------------
1665 
1666 //_____________________________________________________________________________
1667 /**
1668  * Apply the default configuration to the model.  This means setting the
1669  * generalized coordinates and speeds to their default values.
1670  */
applyDefaultConfiguration(SimTK::State & s)1671 void Model::applyDefaultConfiguration(SimTK::State& s)
1672 {
1673     int i;
1674 
1675     // Coordinates
1676     int ncoords = getCoordinateSet().getSize();
1677 
1678     for(i=0; i<ncoords; i++) {
1679         Coordinate& coord = getCoordinateSet().get(i);
1680         coord.setValue(s, coord.getDefaultValue(), false);
1681         coord.setSpeedValue(s, coord.getDefaultSpeedValue());
1682     }
1683 
1684     // Satisfy the constraints.
1685     assemble(s);
1686 }
1687 
1688 //_____________________________________________________________________________
1689 /**
1690  * createAssemblySolver anew, old solver is deleted first. This should be invoked whenever changes in locks and/or constraint status
1691  * that has the potential to affect model assembly.
1692  */
createAssemblySolver(const SimTK::State & s)1693 void Model::createAssemblySolver(const SimTK::State& s)
1694 {
1695     // Allocate on stack and pass to AssemblySolver to make working copy.
1696     SimTK::Array_<CoordinateReference> coordsToTrack;
1697 
1698     for(int i=0; i<getNumCoordinates(); ++i){
1699         // Iff a coordinate is dependent on other coordinates for its value,
1700         // do not give it a reference/goal.
1701         if(!_coordinateSet[i].isDependent(s)){
1702             Constant reference(_coordinateSet[i].getValue(s));
1703             CoordinateReference coordRef(_coordinateSet[i].getName(), reference);
1704             coordsToTrack.push_back(coordRef);
1705         }
1706     }
1707 
1708     // Use the assembler to generate the initial pose from Coordinate defaults
1709     // that also satisfies the constraints. AssemblySolver makes copy of
1710     // coordsToTrack
1711     _assemblySolver.reset(new AssemblySolver(*this, coordsToTrack));
1712     _assemblySolver->setConstraintWeight(SimTK::Infinity);
1713     _assemblySolver->setAccuracy(get_assembly_accuracy());
1714 }
1715 
updateAssemblyConditions(SimTK::State & s)1716 void Model::updateAssemblyConditions(SimTK::State& s)
1717 {
1718     createAssemblySolver(s);
1719 }
1720 //--------------------------------------------------------------------------
1721 // MARKERS
1722 //--------------------------------------------------------------------------
1723 //_____________________________________________________________________________
1724 /**
1725  * Write an XML file of all the markers in the model.
1726  *
1727  * @param aFileName the name of the file to create
1728  */
writeMarkerFile(const string & aFileName)1729 void Model::writeMarkerFile(const string& aFileName)
1730 {
1731     upd_MarkerSet().print(aFileName);
1732 }
1733 
1734 
1735 //_____________________________________________________________________________
updateMarkerSet(MarkerSet & newMarkerSet)1736 void Model::updateMarkerSet(MarkerSet& newMarkerSet)
1737 {
1738     for (int i = 0; i < newMarkerSet.getSize(); i++) {
1739         Marker& updatingMarker = newMarkerSet.get(i);
1740 
1741         /* If there is already a marker in the model with that name,
1742          * replace it with the updating marker.
1743          */
1744         if (updMarkerSet().contains(updatingMarker.getName())) {
1745             Marker& modelMarker = updMarkerSet().get(updatingMarker.getName());
1746             // Delete the marker from the model and add the updating one
1747             upd_MarkerSet().remove(&modelMarker);
1748         }
1749         // append the marker to the model's Set
1750         addMarker(updatingMarker.clone());
1751     }
1752 
1753     cout << "Updated markers in model " << getName() << endl;
1754 }
1755 
1756 //_____________________________________________________________________________
1757 /**
1758  * Remove all markers from the model that are not in the passed-in list.
1759  *
1760  * @param aMarkerNames array of marker names not to be deleted
1761  * @return Number of markers deleted
1762  *
1763  * @Todo_AYMAN make sure visuals adjust as well
1764  */
deleteUnusedMarkers(const OpenSim::Array<string> & aMarkerNames)1765 int Model::deleteUnusedMarkers(const OpenSim::Array<string>& aMarkerNames)
1766 {
1767     int i, numDeleted = 0;
1768 
1769     for (i = 0; i < get_MarkerSet().getSize(); )
1770     {
1771         int index = aMarkerNames.findIndex(get_MarkerSet().get(i).getName());
1772         if (index < 0)
1773         {
1774             // Delete the marker, but don't increment i or else you'll
1775             // skip over the marker right after the deleted one.
1776             upd_MarkerSet().remove(i);
1777             numDeleted++;
1778         }
1779         else
1780         {
1781             i++;
1782         }
1783     }
1784 
1785     cout << "Deleted " << numDeleted << " unused markers from model " << getName() << endl;
1786 
1787     return numDeleted;
1788 }
1789 
1790 /**
1791  ** Get a flat list of Joints contained in the model
1792  **
1793  **/
updJointSet()1794 JointSet& Model::updJointSet()
1795 {
1796     return upd_JointSet();
1797 }
1798 
getJointSet() const1799 const JointSet& Model::getJointSet() const
1800 {
1801     return get_JointSet();
1802 }
1803 
getGround() const1804 const Ground& Model::getGround() const
1805 {
1806     return get_ground();
1807 }
1808 
updGround()1809 Ground& Model::updGround()
1810 {
1811     return upd_ground();
1812 }
1813 
1814 //--------------------------------------------------------------------------
1815 // CONTROLS
1816 //--------------------------------------------------------------------------
1817 /** Get the number of controls for this the model.
1818  * Throws an exception if called before Model::initSystem()  */
getNumControls() const1819 int Model::getNumControls() const
1820 {
1821     if(!_system){
1822         throw Exception("Model::getNumControls() requires an initialized Model./n"
1823             "Prior Model::initSystem() required.");
1824     }
1825 
1826     return _defaultControls.size();
1827 }
1828 
1829 /** Update the controls for this the model at a given state.
1830  * Throws an exception if called before Model::initSystem() */
updControls(const SimTK::State & s) const1831 Vector& Model::updControls(const SimTK::State &s) const
1832 {
1833     if( (!_system) || (!_modelControlsIndex.isValid()) ){
1834         throw Exception("Model::updControls() requires an initialized Model./n"
1835             "Prior call to Model::initSystem() is required.");
1836     }
1837 
1838     // direct the system shared cache
1839     Measure_<Vector>::Result controlsCache =
1840         Measure_<Vector>::Result::getAs(_system->updDefaultSubsystem()
1841             .getMeasure(_modelControlsIndex));
1842     return controlsCache.updValue(s);
1843 }
1844 
markControlsAsValid(const SimTK::State & s) const1845 void Model::markControlsAsValid(const SimTK::State& s) const
1846 {
1847     if( (!_system) || (!_modelControlsIndex.isValid()) ){
1848         throw Exception("Model::markControlsAsValid() requires an initialized Model./n"
1849             "Prior call to Model::initSystem() is required.");
1850     }
1851 
1852     Measure_<Vector>::Result controlsCache =
1853         Measure_<Vector>::Result::getAs(_system->updDefaultSubsystem()
1854             .getMeasure(_modelControlsIndex));
1855     controlsCache.markAsValid(s);
1856 }
1857 
setControls(const SimTK::State & s,const SimTK::Vector & controls) const1858 void Model::setControls(const SimTK::State& s, const SimTK::Vector& controls) const
1859 {
1860     if( (!_system) || (!_modelControlsIndex.isValid()) ){
1861         throw Exception("Model::setControls() requires an initialized Model./n"
1862             "Prior call to Model::initSystem() is required.");
1863     }
1864 
1865     // direct the system shared cache
1866     Measure_<Vector>::Result controlsCache =
1867         Measure_<Vector>::Result::getAs(_system->updDefaultSubsystem()
1868         .getMeasure(_modelControlsIndex));
1869     controlsCache.setValue(s, controls);
1870 
1871     // Make sure to re-realize dynamics to make sure controls can affect forces
1872     // and not just derivatives
1873     if(s.getSystemStage() == Stage::Dynamics)
1874         s.invalidateAllCacheAtOrAbove(Stage::Dynamics);
1875 }
1876 
1877 /** Const access to controls does not invalidate dynamics */
getControls(const SimTK::State & s) const1878 const Vector& Model::getControls(const SimTK::State &s) const
1879 {
1880     if( (!_system) || (!_modelControlsIndex.isValid()) ){
1881         throw Exception("Model::getControls() requires an initialized Model./n"
1882             "Prior call to Model::initSystem() is required.");
1883     }
1884 
1885     // direct the system shared cache
1886     Measure_<Vector>::Result controlsCache = Measure_<Vector>::Result::getAs(_system->updDefaultSubsystem().getMeasure(_modelControlsIndex));
1887 
1888     if(!controlsCache.isValid(s)){
1889         // Always reset controls to their default values before computing controls
1890         // since default behavior is for controllers to "addInControls" so there should be valid
1891         // values to begin with.
1892         controlsCache.updValue(s) = _defaultControls;
1893         computeControls(s, controlsCache.updValue(s));
1894         controlsCache.markAsValid(s);
1895     }
1896 
1897     return controlsCache.getValue(s);
1898 }
1899 
1900 
1901 /** Compute the controls the model */
computeControls(const SimTK::State & s,SimTK::Vector & controls) const1902 void Model::computeControls(const SimTK::State& s, SimTK::Vector &controls) const
1903 {
1904     for (auto& controller : getComponentList<Controller>()) {
1905         if (controller.isEnabled()) {
1906             controller.computeControls(s, controls);
1907         }
1908     }
1909 }
1910 
1911 
1912 /** Get a flag indicating if the model needs controls to operate its actuators */
isControlled() const1913 bool Model::isControlled() const
1914 {
1915     bool isControlled = getActuators().getSize() > 0;
1916 
1917     return isControlled;
1918 }
1919 
getControllerSet() const1920 const ControllerSet& Model::getControllerSet() const{
1921     return(get_ControllerSet());
1922 }
updControllerSet()1923 ControllerSet& Model::updControllerSet() {
1924     return(upd_ControllerSet());
1925 }
storeControls(const SimTK::State & s,int step)1926 void Model::storeControls( const SimTK::State& s, int step ) {
1927     upd_ControllerSet().storeControls(s, step);
1928     return;
1929 }
printControlStorage(const string & fileName) const1930 void Model::printControlStorage(const string& fileName ) const {
1931     get_ControllerSet().printControlStorage(fileName);
1932 }
1933 
getControlsTable() const1934 TimeSeriesTable Model::getControlsTable() const {
1935     return get_ControllerSet().getControlTable();
1936 }
1937 
getAllControllersEnabled() const1938 bool Model::getAllControllersEnabled() const{
1939   return( _allControllersEnabled );
1940 }
setAllControllersEnabled(bool enabled)1941 void Model::setAllControllersEnabled( bool enabled ) {
1942     _allControllersEnabled = enabled;
1943 }
1944 
formStateStorage(const Storage & originalStorage,Storage & statesStorage,bool warnUnspecifiedStates) const1945 void Model::formStateStorage(const Storage& originalStorage,
1946                              Storage& statesStorage,
1947                              bool warnUnspecifiedStates) const
1948 {
1949     Array<string> rStateNames = getStateVariableNames();
1950     int numStates = getNumStateVariables();
1951     // make sure same size, otherwise warn
1952     if (originalStorage.getSmallestNumberOfStates() != rStateNames.getSize() && warnUnspecifiedStates){
1953         cout << "Number of columns does not match in formStateStorage. Found "
1954             << originalStorage.getSmallestNumberOfStates() << " Expected  " << rStateNames.getSize() << "." << endl;
1955     }
1956 
1957     OPENSIM_THROW_IF_FRMOBJ(originalStorage.isInDegrees(), Exception,
1958         "Input Storage must have values for Coordinates in meters or radians (not degrees).\n"
1959         "Please convert input Storage to meters and/or degrees first.");
1960 
1961     // when the state value is not found in the storage use its default value in the State
1962     SimTK::Vector defaultStateValues = getStateVariableValues(getWorkingState());
1963 
1964     // Create a list with entry for each desiredName telling which column in originalStorage has the data
1965     Array<int> mapColumns(-1, rStateNames.getSize());
1966     for(int i=0; i< rStateNames.getSize(); i++){
1967         // the index is -1 if not found, >=1 otherwise since time has index 0 by defn.
1968         int fix = originalStorage.getStateIndex(rStateNames[i]);
1969         mapColumns[i] = fix;
1970         if (fix==-1 && warnUnspecifiedStates){
1971             cout << "Column "<< rStateNames[i] <<
1972                 " not found by Model::formStateStorage(). "
1973                 "Assuming its default value of "
1974                 << defaultStateValues[i] << endl;
1975         }
1976     }
1977     // Now cycle through each state (row of Storage) and form the Model consistent
1978     // order for the state values
1979     double assignedValue = SimTK::NaN;
1980     for (int row =0; row< originalStorage.getSize(); row++){
1981         StateVector* originalVec = originalStorage.getStateVector(row);
1982         StateVector stateVec{originalVec->getTime()};
1983         stateVec.getData().setSize(numStates);
1984         for(int column=0; column< numStates; column++) {
1985             if (mapColumns[column] != -1)
1986                 originalVec->getDataValue(mapColumns[column], assignedValue);
1987             else
1988                 assignedValue = defaultStateValues[column];
1989 
1990             stateVec.setDataValue(column, assignedValue);
1991         }
1992         statesStorage.append(stateVec);
1993     }
1994     rStateNames.insert(0, "time");
1995     statesStorage.setColumnLabels(rStateNames);
1996 }
1997 
1998 /**
1999  * Model::formStateStorage is intended to take any storage and populate qStorage.
2000  * stateStorage is supposed to be a Storage with labels identical to those obtained by calling
2001  * Model::getStateNames(). Columns/entries found in the "originalStorage" are copied to the
2002  * output qStorage. Entries not found are populated with 0s.
2003  */
formQStorage(const Storage & originalStorage,Storage & qStorage)2004 void Model::formQStorage(const Storage& originalStorage, Storage& qStorage) {
2005 
2006     int nq = getWorkingState().getNQ();
2007     Array<string> qNames;
2008     getCoordinateSet().getNames(qNames);
2009 
2010 
2011     int* mapColumns = new int[qNames.getSize()];
2012     for(int i=0; i< nq; i++){
2013         // the index is -1 if not found, >=1 otherwise since time has index 0 by defn.
2014         mapColumns[i] = originalStorage.getColumnLabels().findIndex(qNames[i]);
2015         if (mapColumns[i]==-1)
2016             cout << "\n Column "<< qNames[i] << " not found in formQStorage, assuming 0.\n" << endl;
2017     }
2018 
2019 
2020     // Now cycle through and shuffle each
2021     for (int row =0; row< originalStorage.getSize(); row++){
2022         StateVector* originalVec = originalStorage.getStateVector(row);
2023         StateVector* stateVec = new StateVector(originalVec->getTime());
2024         stateVec->getData().setSize(nq);  // default value 0f 0.
2025         for(int column=0; column< nq; column++){
2026             double valueInOriginalStorage=0.0;
2027             if (mapColumns[column]!=-1)
2028                 originalVec->getDataValue(mapColumns[column]-1, valueInOriginalStorage);
2029 
2030             stateVec->setDataValue(column, valueInOriginalStorage);
2031         }
2032         qStorage.append(*stateVec);
2033     }
2034     qNames.insert(0, "time");
2035 
2036     qStorage.setColumnLabels(qNames);
2037     // Since we're copying data from a Storage file, keep the inDegrees flag consistent
2038     qStorage.setInDegrees(originalStorage.isInDegrees());
2039 }
disownAllComponents()2040 void Model::disownAllComponents()
2041 {
2042     updMiscModelComponentSet().setMemoryOwner(false);
2043     updBodySet().setMemoryOwner(false);
2044     updJointSet().setMemoryOwner(false);
2045     updConstraintSet().setMemoryOwner(false);
2046     updForceSet().setMemoryOwner(false);
2047     updContactGeometrySet().setMemoryOwner(false);
2048     updControllerSet().setMemoryOwner(false);
2049     updAnalysisSet().setMemoryOwner(false);
2050     updMarkerSet().setMemoryOwner(false);
2051     updProbeSet().setMemoryOwner(false);
2052 }
2053 
overrideAllActuators(SimTK::State & s,bool flag)2054 void Model::overrideAllActuators( SimTK::State& s, bool flag) {
2055      Set<Actuator>& as = updActuators();
2056 
2057      for(int i=0;i<as.getSize();i++) {
2058          ScalarActuator* act = dynamic_cast<ScalarActuator*>(&as[i]);
2059          act->overrideActuation(s, flag);
2060      }
2061 
2062 }
2063 
getObjectByTypeAndName(const std::string & typeString,const std::string & nameString)2064 const Object& Model::getObjectByTypeAndName(const std::string& typeString, const std::string& nameString) {
2065     if (typeString=="Body")
2066         return getBodySet().get(nameString);
2067     else if (typeString == "Joint")
2068         return getJointSet().get(nameString);
2069     else if (typeString=="Force")
2070         return getForceSet().get(nameString);
2071     else if (typeString=="Constraint")
2072         return getConstraintSet().get(nameString);
2073     else if (typeString=="Coordinate")
2074         return getCoordinateSet().get(nameString);
2075     else if (typeString=="Marker")
2076         return getMarkerSet().get(nameString);
2077     else if (typeString=="Controller")
2078         return getControllerSet().get(nameString);
2079     else if (typeString == "Probe")
2080         return getProbeSet().get(nameString);
2081     throw Exception("Model::getObjectByTypeAndName: no object of type " + typeString +
2082         " and name "+nameString+" was found in the model.");
2083 
2084 }
2085 
2086 //------------------------------------------------------------------------------
2087 //          REALIZE THE SYSTEM TO THE REQUIRED COMPUTATIONAL STAGE
2088 //------------------------------------------------------------------------------
realizeTime(const SimTK::State & state) const2089 void Model::realizeTime(const SimTK::State& state) const
2090 {
2091     getSystem().realize(state, Stage::Time);
2092 }
2093 
realizePosition(const SimTK::State & state) const2094 void Model::realizePosition(const SimTK::State& state) const
2095 {
2096     getSystem().realize(state, Stage::Position);
2097 }
2098 
realizeVelocity(const SimTK::State & state) const2099 void Model::realizeVelocity(const SimTK::State& state) const
2100 {
2101     getSystem().realize(state, Stage::Velocity);
2102 }
2103 
realizeDynamics(const SimTK::State & state) const2104 void Model::realizeDynamics(const SimTK::State& state) const
2105 {
2106     getSystem().realize(state, Stage::Dynamics);
2107 }
2108 
realizeAcceleration(const SimTK::State & state) const2109 void Model::realizeAcceleration(const SimTK::State& state) const
2110 {
2111     getSystem().realize(state, Stage::Acceleration);
2112 }
2113 
realizeReport(const SimTK::State & state) const2114 void Model::realizeReport(const SimTK::State& state) const
2115 {
2116     getSystem().realize(state, Stage::Report);
2117 }
2118 
2119 
2120 /**
2121  * Compute the derivatives of the generalized coordinates and speeds.
2122  */
computeStateVariableDerivatives(const SimTK::State & s) const2123 void Model::computeStateVariableDerivatives(const SimTK::State &s) const
2124 {
2125     realizeAcceleration(s);
2126 }
2127 
2128 /**
2129  * Get the total mass of the model
2130  *
2131  * @return the mass of the model
2132  */
getTotalMass(const SimTK::State & s) const2133 double Model::getTotalMass(const SimTK::State &s) const
2134 {
2135     return getMatterSubsystem().calcSystemMass(s);
2136 }
2137 /**
2138  * getInertiaAboutMassCenter
2139  *
2140  */
getInertiaAboutMassCenter(const SimTK::State & s) const2141 SimTK::Inertia Model::getInertiaAboutMassCenter(const SimTK::State &s) const
2142 {
2143     SimTK::Inertia inertia = getMatterSubsystem().calcSystemCentralInertiaInGround(s);
2144 
2145     return inertia;
2146 }
2147 /**
2148  * Return the position vector of the system mass center, measured from the Ground origin, and expressed in Ground.
2149  *
2150  */
calcMassCenterPosition(const SimTK::State & s) const2151 SimTK::Vec3 Model::calcMassCenterPosition(const SimTK::State &s) const
2152 {
2153     getMultibodySystem().realize(s, Stage::Position);
2154     return getMatterSubsystem().calcSystemMassCenterLocationInGround(s);
2155 }
2156 /**
2157  * Return the velocity vector of the system mass center, measured from the Ground origin, and expressed in Ground.
2158  *
2159  */
calcMassCenterVelocity(const SimTK::State & s) const2160 SimTK::Vec3 Model::calcMassCenterVelocity(const SimTK::State &s) const
2161 {
2162     getMultibodySystem().realize(s, Stage::Velocity);
2163     return getMatterSubsystem().calcSystemMassCenterVelocityInGround(s);
2164 }
2165 /**
2166  * Return the acceleration vector of the system mass center, measured from the Ground origin, and expressed in Ground.
2167  *
2168  */
calcMassCenterAcceleration(const SimTK::State & s) const2169 SimTK::Vec3 Model::calcMassCenterAcceleration(const SimTK::State &s) const
2170 {
2171     getMultibodySystem().realize(s, Stage::Acceleration);
2172     return getMatterSubsystem().calcSystemMassCenterAccelerationInGround(s);
2173 }
2174 
2175 /**
2176 * Construct outputs
2177 *
2178 **/
2179 
2180