Drake
System< T > Class Template Referenceabstract

Base class for all System functionality that is dependent on the templatized scalar type T for input, state, parameters, and outputs. More...

#include <drake/systems/framework/witness_function.h>

Public Member Functions

 ~System () override=default
 
void GetWitnessFunctions (const Context< T > &context, std::vector< const WitnessFunction< T > *> *w) const
 Gets the witness functions active for the given state. More...
 
CalcWitnessValue (const Context< T > &context, const WitnessFunction< T > &witness_func) const
 Evaluates a witness function at the given context. More...
 
virtual void AddTriggeredWitnessFunctionToCompositeEventCollection (Event< T > *event, CompositeEventCollection< T > *events) const =0
 Add event to events due to a witness function triggering. More...
 
Does not allow copy, move, or assignment
 System (const System &)=delete
 
Systemoperator= (const System &)=delete
 
 System (System &&)=delete
 
Systemoperator= (System &&)=delete
 
Resource allocation and initialization

These methods are used to allocate and initialize Context resources.

std::unique_ptr< Context< T > > AllocateContext () const
 Returns a Context<T> suitable for use with this System<T>. More...
 
virtual std::unique_ptr< CompositeEventCollection< T > > AllocateCompositeEventCollection () const =0
 Allocates a CompositeEventCollection for this system. More...
 
std::unique_ptr< BasicVector< T > > AllocateInputVector (const InputPort< T > &input_port) const
 Given an input port, allocates the vector storage. More...
 
std::unique_ptr< AbstractValueAllocateInputAbstract (const InputPort< T > &input_port) const
 Given an input port, allocates the abstract storage. More...
 
std::unique_ptr< SystemOutput< T > > AllocateOutput () const
 Returns a container that can hold the values of all of this System's output ports. More...
 
virtual std::unique_ptr< ContinuousState< T > > AllocateTimeDerivatives () const
 Returns a ContinuousState of the same size as the continuous_state allocated in CreateDefaultContext. More...
 
virtual std::unique_ptr< DiscreteValues< T > > AllocateDiscreteVariables () const
 Returns a DiscreteState of the same dimensions as the discrete_state allocated in CreateDefaultContext. More...
 
std::unique_ptr< Context< T > > CreateDefaultContext () const
 This convenience method allocates a context using AllocateContext() and sets its default values using SetDefaultContext(). More...
 
virtual void SetDefaultState (const Context< T > &context, State< T > *state) const =0
 Assigns default values to all elements of the state. More...
 
virtual void SetDefaultParameters (const Context< T > &context, Parameters< T > *parameters) const =0
 Assigns default values to all parameters. More...
 
void SetDefaultContext (Context< T > *context) const
 
virtual void SetRandomState (const Context< T > &context, State< T > *state, RandomGenerator *generator) const
 Assigns random values to all elements of the state. More...
 
virtual void SetRandomParameters (const Context< T > &context, Parameters< T > *parameters, RandomGenerator *generator) const
 Assigns random values to all parameters. More...
 
void SetRandomContext (Context< T > *context, RandomGenerator *generator) const
 
void AllocateFixedInputs (Context< T > *context) const
 For each input port, allocates a fixed input of the concrete type that this System requires, and binds it to the port, disconnecting any prior input. More...
 
virtual std::multimap< int, intGetDirectFeedthroughs () const =0
 Reports all direct feedthroughs from input ports to output ports. More...
 
bool HasAnyDirectFeedthrough () const
 Returns true if any of the inputs to the system might be directly fed through to any of its outputs and false otherwise. More...
 
bool HasDirectFeedthrough (int output_port) const
 Returns true if there might be direct-feedthrough from any input port to the given output_port, and false otherwise. More...
 
bool HasDirectFeedthrough (int input_port, int output_port) const
 Returns true if there might be direct-feedthrough from the given input_port to the given output_port, and false otherwise. More...
 
Publishing

Publishing is the primary mechanism for a System to communicate with the world outside the System abstraction during a simulation.

Publishing occurs at user-specified times or events and can generate side-effect results such as terminal output, visualization, logging, plotting, and network messages. Other than computational cost, publishing has no effect on the progress of a simulation.

void Publish (const Context< T > &context, const EventCollection< PublishEvent< T >> &events) const
 This method is the public entry point for dispatching all publish event handlers. More...
 
void Publish (const Context< T > &context) const
 Forces a publish on the system, given a context. More...
 
Cached evaluations

Given the values in a Context, a Drake System must be able to provide the results of particular computations needed for analysis and simulation of the System.

These results are maintained in a mutable cache within the Context so that a result need be computed only once, the first time it is requested after a change to one of its prerequisite values.

The Eval methods in this group return a reference to the already-computed result in the given Context's cache. If the current value is out of date, they first update the cache entry using the corresponding Calc method from the "Calculations" group. Evaluations of input ports instead delegate to the containing Diagram, which arranges to have the appropriate subsystem evaluate the source output port.

Methods in this group that specify preconditions operate as follows: The preconditions will be checked in Debug builds but some or all might not be checked in Release builds for performance reasons. If we do check and a precondition is violated, an std::logic_error will be thrown with a helpful message.

const ContinuousState< T > & EvalTimeDerivatives (const Context< T > &context) const
 Returns a reference to the cached value of the continuous state variable time derivatives, evaluating first if necessary using CalcTimeDerivatives(). More...
 
const T & EvalPotentialEnergy (const Context< T > &context) const
 Returns a reference to the cached value of the potential energy (PE), evaluating first if necessary using CalcPotentialEnergy(). More...
 
const T & EvalKineticEnergy (const Context< T > &context) const
 Returns a reference to the cached value of the kinetic energy (KE), evaluating first if necessary using CalcKineticEnergy(). More...
 
const T & EvalConservativePower (const Context< T > &context) const
 Returns a reference to the cached value of the conservative power (Pc), evaluating first if necessary using CalcConservativePower(). More...
 
const T & EvalNonConservativePower (const Context< T > &context) const
 Returns a reference to the cached value of the non-conservative power (Pnc), evaluating first if necessary using CalcNonConservativePower(). More...
 
template<template< typename > class Vec = BasicVector>
const Vec< T > * EvalVectorInput (const Context< T > &context, int port_index) const
 Returns the value of the vector-valued input port with the given port_index as a BasicVector or a specific subclass Vec derived from BasicVector. More...
 
Eigen::VectorBlock< const VectorX< T > > EvalEigenVectorInput (const Context< T > &context, int port_index) const
 Returns the value of the vector-valued input port with the given port_index as an Eigen vector. More...
 
Constraint-related functions
int get_num_constraint_equations (const Context< T > &context) const
 Gets the number of constraint equations for this system using the given context (useful in case the number of constraints is dependent upon the current state (as might be the case with a system modeled using piecewise differential algebraic equations). More...
 
Eigen::VectorXd EvalConstraintEquations (const Context< T > &context) const
 Evaluates the constraint equations for the system at the generalized coordinates and generalized velocity specified by the context. More...
 
Eigen::VectorXd EvalConstraintEquationsDot (const Context< T > &context) const
 Computes the time derivative of each constraint equation, evaluated at the generalized coordinates and generalized velocity specified by the context. More...
 
Eigen::VectorXd CalcVelocityChangeFromConstraintImpulses (const Context< T > &context, const Eigen::MatrixXd &J, const Eigen::VectorXd &lambda) const
 Computes the change in velocity from applying the given constraint forces to the system at the given context. More...
 
double CalcConstraintErrorNorm (const Context< T > &context, const Eigen::VectorXd &error) const
 Computes the norm on constraint error (used as a metric for comparing errors between the outputs of algebraic equations applied to two different state variable instances). More...
 
Calculations

A Drake System defines a set of common computations that are understood by the framework.

Most of these are embodied in a Calc method that unconditionally performs the calculation into an output argument of the appropriate type, using only values from the given Context. These are paired with an Eval method that returns a reference to an already-calculated result residing in the cache; if needed that result is first obtained using the Calc method. See the "Evaluations" group for more information.

This group also includes additional System-specific operations that depend on both Context and additional input arguments.

void CalcTimeDerivatives (const Context< T > &context, ContinuousState< T > *derivatives) const
 Calculates the time derivatives xcdot of the continuous state xc into a given output argument. More...
 
void CalcDiscreteVariableUpdates (const Context< T > &context, const EventCollection< DiscreteUpdateEvent< T >> &events, DiscreteValues< T > *discrete_state) const
 This method is the public entry point for dispatching all discrete variable update event handlers. More...
 
void CalcDiscreteVariableUpdates (const Context< T > &context, DiscreteValues< T > *discrete_state) const
 This method forces a discrete update on the system given a context, and the updated discrete state is stored in discrete_state. More...
 
void CalcUnrestrictedUpdate (const Context< T > &context, const EventCollection< UnrestrictedUpdateEvent< T >> &events, State< T > *state) const
 This method is the public entry point for dispatching all unrestricted update event handlers. More...
 
void CalcUnrestrictedUpdate (const Context< T > &context, State< T > *state) const
 This method forces an unrestricted update on the system given a context, and the updated state is stored in state. More...
 
CalcNextUpdateTime (const Context< T > &context, CompositeEventCollection< T > *events) const
 This method is called by a Simulator during its calculation of the size of the next continuous step to attempt. More...
 
void GetPerStepEvents (const Context< T > &context, CompositeEventCollection< T > *events) const
 This method is called by Simulator::Initialize() to gather all update and publish events that are to be handled in StepTo() at the point before Simulator integrates continuous state. More...
 
void GetInitializationEvents (const Context< T > &context, CompositeEventCollection< T > *events) const
 This method is called by Simulator::Initialize() to gather all update and publish events that need to be handled at initialization before the simulator starts integration. More...
 
optional< PeriodicEventDataGetUniquePeriodicDiscreteUpdateAttribute () const
 Gets whether there exists a unique periodic attribute that triggers one or more discrete update events (and, if so, returns that unique periodic attribute). More...
 
std::map< PeriodicEventData, std::vector< const Event< T > * >, PeriodicEventDataComparatorGetPeriodicEvents () const
 Gets all periodic triggered events for a system. More...
 
void CalcOutput (const Context< T > &context, SystemOutput< T > *outputs) const
 Utility method that computes for every output port i the value y(i) that should result from the current contents of the given Context. More...
 
CalcPotentialEnergy (const Context< T > &context) const
 Calculates and returns the potential energy represented by the current configuration provided in context. More...
 
CalcKineticEnergy (const Context< T > &context) const
 Calculates and returns the kinetic energy represented by the current configuration and velocity provided in context. More...
 
CalcConservativePower (const Context< T > &context) const
 Calculates and returns the conservative power represented by the current contents of the given context. More...
 
CalcNonConservativePower (const Context< T > &context) const
 Calculates and returns the non-conservative power represented by the current contents of the given context. More...
 
void MapVelocityToQDot (const Context< T > &context, const VectorBase< T > &generalized_velocity, VectorBase< T > *qdot) const
 Transforms a given generalized velocity v to the time derivative qdot of the generalized configuration q taken from the supplied Context. More...
 
void MapVelocityToQDot (const Context< T > &context, const Eigen::Ref< const VectorX< T >> &generalized_velocity, VectorBase< T > *qdot) const
 Transforms the given generalized velocity to the time derivative of generalized configuration. More...
 
void MapQDotToVelocity (const Context< T > &context, const VectorBase< T > &qdot, VectorBase< T > *generalized_velocity) const
 Transforms the time derivative qdot of the generalized configuration q to generalized velocities v. More...
 
void MapQDotToVelocity (const Context< T > &context, const Eigen::Ref< const VectorX< T >> &qdot, VectorBase< T > *generalized_velocity) const
 Transforms the given time derivative qdot of generalized configuration q to generalized velocity v. More...
 
Utility methods
std::string GetMemoryObjectName () const
 Returns a name for this System based on a stringification of its type name and memory address. More...
 
const InputPort< T > & get_input_port (int port_index) const
 Returns the typed input port at index port_index. More...
 
const InputPort< T > & GetInputPort (const std::string &port_name) const
 Returns the typed input port with the unique name port_name. More...
 
const OutputPort< T > & get_output_port (int port_index) const
 Returns the typed output port at index port_index. More...
 
const OutputPort< T > & GetOutputPort (const std::string &port_name) const
 Returns the typed output port with the unique name port_name. More...
 
int get_num_constraints () const
 Returns the number of constraints specified for the system. More...
 
virtual int get_num_continuous_states () const =0
 Returns the dimension of the continuous state vector that has been declared until now. More...
 
const SystemConstraint< T > & get_constraint (SystemConstraintIndex constraint_index) const
 Returns the constraint at index constraint_index. More...
 
boolean< T > CheckSystemConstraintsSatisfied (const Context< T > &context, double tol) const
 Returns true if context satisfies all of the registered SystemConstraints with tolerance tol. More...
 
void CheckValidOutput (const SystemOutput< T > *output) const
 Checks that output is consistent with the number and size of output ports declared by the system. More...
 
template<typename T1 = T>
void CheckValidContextT (const Context< T1 > &context) const
 Checks that context is consistent for this System template. More...
 
VectorX< T > CopyContinuousStateVector (const Context< T > &context) const
 Returns a copy of the continuous state vector xc into an Eigen vector. More...
 
Graphviz methods
std::string GetGraphvizString (int max_depth=std::numeric_limits< int >::max()) const
 Returns a Graphviz string describing this System. More...
 
virtual void GetGraphvizFragment (int max_depth, std::stringstream *dot) const
 Appends a Graphviz fragment to the dot stream. More...
 
virtual void GetGraphvizInputPortToken (const InputPort< T > &port, int max_depth, std::stringstream *dot) const
 Appends a fragment to the dot stream identifying the graphviz node representing port. More...
 
virtual void GetGraphvizOutputPortToken (const OutputPort< T > &port, int max_depth, std::stringstream *dot) const
 Appends a fragment to the dot stream identifying the graphviz node representing port. More...
 
int64_t GetGraphvizId () const
 Returns an opaque integer that uniquely identifies this system in the Graphviz output. More...
 
Scalar type conversion utilities
void FixInputPortsFrom (const System< double > &other_system, const Context< double > &other_context, Context< T > *target_context) const
 Fixes all of the input ports in target_context to their current values in other_context, as evaluated by other_system. More...
 
const SystemScalarConverterget_system_scalar_converter () const
 (Advanced) Returns the SystemScalarConverter for this object. More...
 
- Public Member Functions inherited from SystemBase
 ~SystemBase () override
 
void set_name (const std::string &name)
 Sets the name of the system. More...
 
const std::string & get_name () const
 Returns the name last supplied to set_name(), if any. More...
 
const std::string & GetSystemName () const final
 Returns a human-readable name for this system, for use in messages and logging. More...
 
std::string GetSystemPathname () const final
 Generates and returns a human-readable full path name of this subsystem, for use in messages and logging. More...
 
std::string GetSystemType () const final
 Returns the most-derived type of this concrete System object as a human-readable string suitable for use in error messages. More...
 
void ThrowIfContextNotCompatible (const ContextBase &context) const final
 Throws an exception with an appropriate message if the given context is not compatible with this System. More...
 
std::unique_ptr< ContextBaseAllocateContext () const
 Returns a Context suitable for use with this System. More...
 
int get_num_input_ports () const
 Returns the number of input ports currently allocated in this System. More...
 
int get_num_output_ports () const
 Returns the number of output ports currently allocated in this System. More...
 
const InputPortBaseget_input_port_base (InputPortIndex port_index) const
 Returns a reference to an InputPort given its port_index. More...
 
const OutputPortBaseget_output_port_base (OutputPortIndex port_index) const
 Returns a reference to an OutputPort given its port_index. More...
 
int get_num_total_inputs () const
 Returns the total dimension of all of the vector-valued input ports (as if they were muxed). More...
 
int get_num_total_outputs () const
 Returns the total dimension of all of the vector-valued output ports (as if they were muxed). More...
 
int num_cache_entries () const
 Returns the number nc of cache entries currently allocated in this System. More...
 
const CacheEntryget_cache_entry (CacheIndex index) const
 Return a reference to a CacheEntry given its index. More...
 
void CheckValidContext (const ContextBase &context) const
 Checks whether the given context is valid for this System and throws an exception with a helpful message if not. More...
 
 SystemBase (const SystemBase &)=delete
 
SystemBaseoperator= (const SystemBase &)=delete
 
 SystemBase (SystemBase &&)=delete
 
SystemBaseoperator= (SystemBase &&)=delete
 
const AbstractValueEvalAbstractInput (const ContextBase &context, int port_index) const
 Returns the value of the input port with the given port_index as an AbstractValue, which is permitted for ports of any type. More...
 
template<typename V >
const V * EvalInputValue (const ContextBase &context, int port_index) const
 Returns the value of an abstract-valued input port with the given port_index as a value of known type V. More...
 
const CacheEntryDeclareCacheEntry (std::string description, CacheEntry::AllocCallback alloc_function, CacheEntry::CalcCallback calc_function, std::set< DependencyTicket > prerequisites_of_calc={ all_sources_ticket()})
 Declares a new CacheEntry in this System using the least-restrictive definitions for the associated functions. More...
 
template<class MySystem , class MyContext , typename ValueType >
const CacheEntryDeclareCacheEntry (std::string description, ValueType(MySystem::*make)() const, void(MySystem::*calc)(const MyContext &, ValueType *) const, std::set< DependencyTicket > prerequisites_of_calc={ all_sources_ticket()})
 Declares a cache entry by specifying member functions to use both for the allocator and calculator. More...
 
template<class MySystem , class MyContext , typename ValueType >
const CacheEntryDeclareCacheEntry (std::string description, const ValueType &model_value, void(MySystem::*calc)(const MyContext &, ValueType *) const, std::set< DependencyTicket > prerequisites_of_calc={ all_sources_ticket()})
 Declares a cache entry by specifying a model value of concrete type ValueType and a calculator function that is a class member function (method) with signature: More...
 
template<class MySystem , class MyContext , typename ValueType >
const CacheEntryDeclareCacheEntry (std::string description, const ValueType &model_value, ValueType(MySystem::*calc)(const MyContext &) const, std::set< DependencyTicket > prerequisites_of_calc={ all_sources_ticket()})
 Declares a cache entry by specifying a model value of concrete type ValueType and a calculator function that is a class member function (method) with signature: More...
 
template<class MySystem , class MyContext , typename ValueType >
const CacheEntryDeclareCacheEntry (std::string description, void(MySystem::*calc)(const MyContext &, ValueType *) const, std::set< DependencyTicket > prerequisites_of_calc={ all_sources_ticket()})
 Declares a cache entry by specifying only a calculator function that is a class member function (method) with signature: More...
 
template<class MySystem , class MyContext , typename ValueType >
const CacheEntryDeclareCacheEntry (std::string description, ValueType(MySystem::*calc)(const MyContext &) const, std::set< DependencyTicket > prerequisites_of_calc={ all_sources_ticket()})
 Declares a cache entry by specifying only a calculator function that is a class member function (method) with signature: More...
 
DependencyTicket discrete_state_ticket (DiscreteStateIndex index) const
 Returns a ticket indicating dependence on a particular discrete state variable xdᵢ (may be a vector). More...
 
DependencyTicket abstract_state_ticket (AbstractStateIndex index) const
 Returns a ticket indicating dependence on a particular abstract state variable xaᵢ. More...
 
DependencyTicket numeric_parameter_ticket (NumericParameterIndex index) const
 Returns a ticket indicating dependence on a particular numeric parameter pnᵢ (may be a vector). More...
 
DependencyTicket abstract_parameter_ticket (AbstractParameterIndex index) const
 Returns a ticket indicating dependence on a particular abstract parameter paᵢ. More...
 
DependencyTicket input_port_ticket (InputPortIndex index)
 Returns a ticket indicating dependence on input port uᵢ indicated by index. More...
 
DependencyTicket cache_entry_ticket (CacheIndex index)
 Returns a ticket indicating dependence on the cache entry indicated by index. More...
 
DependencyTicket output_port_ticket (OutputPortIndex index)
 (Internal use only) Returns a ticket indicating dependence on the output port indicated by index. More...
 
int num_discrete_state_groups () const
 Returns the number of declared discrete state groups (each group is a vector-valued discrete state variable). More...
 
int num_abstract_states () const
 Returns the number of declared abstract state variables. More...
 
int num_numeric_parameter_groups () const
 Returns the number of declared numeric parameters (each of these is a vector-valued parameter). More...
 
 DRAKE_DEPRECATED ("Use num_numeric_parameter_groups(). This method will be" " removed after 2/15/19.") int num_numeric_parameters() const
 
int num_abstract_parameters () const
 Returns the number of declared abstract parameters. More...
 

Protected Member Functions

virtual T DoCalcWitnessValue (const Context< T > &context, const WitnessFunction< T > &witness_func) const =0
 Derived classes will implement this method to evaluate a witness function at the given context. More...
 
virtual void DoGetWitnessFunctions (const Context< T > &, std::vector< const WitnessFunction< T > *> *) const
 Derived classes can override this method to provide witness functions active for the given state. More...
 
SystemConstraintIndex AddConstraint (std::unique_ptr< SystemConstraint< T >> constraint)
 Adds an already-created constraint to the list of constraints for this System. More...
 
const EventCollection< PublishEvent< T > > & get_forced_publish_events () const
 
const EventCollection< DiscreteUpdateEvent< T > > & get_forced_discrete_update_events () const
 
const EventCollection< UnrestrictedUpdateEvent< T > > & get_forced_unrestricted_update_events () const
 
void set_forced_publish_events (std::unique_ptr< EventCollection< PublishEvent< T >>> forced)
 
void set_forced_discrete_update_events (std::unique_ptr< EventCollection< DiscreteUpdateEvent< T >>> forced)
 
void set_forced_unrestricted_update_events (std::unique_ptr< EventCollection< UnrestrictedUpdateEvent< T >>> forced)
 
Event handler dispatch mechanism

For a LeafSystem (or user implemented equivalent classes), these functions need to call the appropriate LeafSystem::DoX event handler.

E.g. LeafSystem::DispatchPublishHandler() calls LeafSystem::DoPublish(). User supplied custom event callbacks embedded in each individual event need to be further dispatched in the LeafSystem::DoX handlers if desired. For a LeafSystem, the pseudo code of the complete default publish event handler dispatching is roughly:

  leaf_sys.Publish(context, event_collection)
  -> leaf_sys.DispatchPublishHandler(context, event_collection)
     -> leaf_sys.DoPublish(context, event_collection.get_events())
        -> for (event : event_collection_events):
             if (event.has_handler)
               event.handler(context)

Discrete update events and unrestricted update events are dispatched similarly for a LeafSystem.

For a Diagram (or user implemented equivalent classes), these functions must iterate through all subsystems, extract their corresponding subcontext and subevent collections from context and events, and pass those to the subsystems' public non-virtual event handlers if the subevent collection is nonempty (e.g. System::Publish() for publish events).

All of these functions are only called from their corresponding public non-virtual event dispatchers, where context is error checked. The derived implementations can assume that context is valid. See, e.g., LeafSystem::DispatchPublishHandler() and Diagram::DispatchPublishHandler() for more details.

virtual void DispatchPublishHandler (const Context< T > &context, const EventCollection< PublishEvent< T >> &events) const =0
 This function dispatches all publish events to the appropriate handlers. More...
 
virtual void DispatchDiscreteVariableUpdateHandler (const Context< T > &context, const EventCollection< DiscreteUpdateEvent< T >> &events, DiscreteValues< T > *discrete_state) const =0
 This function dispatches all discrete update events to the appropriate handlers. More...
 
virtual void DispatchUnrestrictedUpdateHandler (const Context< T > &context, const EventCollection< UnrestrictedUpdateEvent< T >> &events, State< T > *state) const =0
 This function dispatches all unrestricted update events to the appropriate handlers. More...
 
System construction

Authors of derived Systems can use these methods in the constructor for those Systems.

 System (SystemScalarConverter converter)
 Constructs an empty System base class object and allocates base class resources, possibly supporting scalar-type conversion support (AutoDiff, etc.) using converter. More...
 
const InputPort< T > & DeclareInputPort (variant< std::string, UseDefaultName > name, PortDataType type, int size, optional< RandomDistribution > random_type=nullopt)
 Adds a port with the specified type and size to the input topology. More...
 
To-be-deprecated declarations

Methods in this section leave out the port name parameter and are the same as invoking the corresponding method with kUseDefaultName as the name.

We intend to make specifying the name required and will deprecate these soon. Don't use them.

const InputPort< T > & DeclareInputPort (PortDataType type, int size, optional< RandomDistribution > random_type=nullopt)
 See the nearly identical signature with an additional (first) argument specifying the port name. More...
 
Virtual methods for calculations

These virtuals allow concrete systems to implement the calculations defined by the Calc methods in the public interface.

Most have default implementations that are usable for simple systems, but you are likely to need to override some or all of these in your concrete system to produce meaningful calculations.

These methods are invoked by the corresponding method in the public interface that has the same name with Do removed. The public method performs error checking on the arguments so you do not need to do so in your implementation. Users cannot invoke these directly since they are protected. You should place your overrides in the protected or private sections of your concrete class.

virtual void DoCalcTimeDerivatives (const Context< T > &context, ContinuousState< T > *derivatives) const
 Override this if you have any continuous state variables xc in your concrete System to calculate their time derivatives. More...
 
virtual void DoCalcNextUpdateTime (const Context< T > &context, CompositeEventCollection< T > *events, T *time) const
 Computes the next time at which this System must perform a discrete action. More...
 
virtual std::map< PeriodicEventData, std::vector< const Event< T > * >, PeriodicEventDataComparatorDoGetPeriodicEvents () const =0
 Implement this method to return all periodic triggered events. More...
 
virtual void DoGetPerStepEvents (const Context< T > &context, CompositeEventCollection< T > *events) const
 Implement this method to return any events to be handled before the simulator integrates the system's continuous state at each time step. More...
 
virtual void DoGetInitializationEvents (const Context< T > &context, CompositeEventCollection< T > *events) const
 Implement this method to return any events to be handled at the simulator's initialization step. More...
 
virtual T DoCalcPotentialEnergy (const Context< T > &context) const
 Override this method for physical systems to calculate the potential energy PE currently stored in the configuration provided in the given Context. More...
 
virtual T DoCalcKineticEnergy (const Context< T > &context) const
 Override this method for physical systems to calculate the kinetic energy KE currently present in the motion provided in the given Context. More...
 
virtual T DoCalcConservativePower (const Context< T > &context) const
 Override this method to return the rate Pc at which mechanical energy is being converted from potential energy to kinetic energy by this system in the given Context. More...
 
virtual T DoCalcNonConservativePower (const Context< T > &context) const
 Override this method to return the rate Pnc at which work W is done on the system by non-conservative forces. More...
 
virtual void DoMapQDotToVelocity (const Context< T > &context, const Eigen::Ref< const VectorX< T >> &qdot, VectorBase< T > *generalized_velocity) const
 Provides the substantive implementation of MapQDotToVelocity(). More...
 
virtual void DoMapVelocityToQDot (const Context< T > &context, const Eigen::Ref< const VectorX< T >> &generalized_velocity, VectorBase< T > *qdot) const
 Provides the substantive implementation of MapVelocityToQDot(). More...
 
Constraint-related functions (protected).
virtual int do_get_num_constraint_equations (const Context< T > &context) const
 Gets the number of constraint equations for this system from the given context. More...
 
virtual Eigen::VectorXd DoEvalConstraintEquations (const Context< T > &context) const
 Evaluates the constraint equations for the system at the generalized coordinates and generalized velocity specified by the context. More...
 
virtual Eigen::VectorXd DoEvalConstraintEquationsDot (const Context< T > &context) const
 Computes the time derivative of each constraint equation, evaluated at the generalized coordinates and generalized velocity specified by the context. More...
 
virtual Eigen::VectorXd DoCalcVelocityChangeFromConstraintImpulses (const Context< T > &context, const Eigen::MatrixXd &J, const Eigen::VectorXd &lambda) const
 Computes the change in velocity from applying the given constraint forces to the system at the given context. More...
 
virtual double DoCalcConstraintErrorNorm (const Context< T > &context, const Eigen::VectorXd &error) const
 Computes the norm of the constraint error. More...
 
Utility methods (protected)
Eigen::VectorBlock< VectorX< T > > GetMutableOutputVector (SystemOutput< T > *output, int port_index) const
 Returns a mutable Eigen expression for a vector valued output port with index port_index in this system. More...
 
- Protected Member Functions inherited from SystemBase
 SystemBase ()=default
 (Internal use only) Default constructor. More...
 
void AddInputPort (std::unique_ptr< InputPortBase > port)
 (Internal use only) Adds an already-constructed input port to this System. More...
 
void AddOutputPort (std::unique_ptr< OutputPortBase > port)
 (Internal use only) Adds an already-constructed output port to this System. More...
 
std::string NextInputPortName (variant< std::string, UseDefaultName > given_name) const
 (Internal use only) Returns a name for the next input port, using the given name if it isn't kUseDefaultName, otherwise making up a name like "u3" from the next available input port index. More...
 
std::string NextOutputPortName (variant< std::string, UseDefaultName > given_name) const
 (Internal use only) Returns a name for the next output port, using the given name if it isn't kUseDefaultName, otherwise making up a name like "y3" from the next available output port index. More...
 
void AddDiscreteStateGroup (DiscreteStateIndex index)
 (Internal use only) Assigns a ticket to a new discrete variable group with the given index. More...
 
void AddAbstractState (AbstractStateIndex index)
 (Internal use only) Assigns a ticket to a new abstract state variable with the given index. More...
 
void AddNumericParameter (NumericParameterIndex index)
 (Internal use only) Assigns a ticket to a new numeric parameter with the given index. More...
 
void AddAbstractParameter (AbstractParameterIndex index)
 (Internal use only) Assigns a ticket to a new abstract parameter with the given index. More...
 
const CacheEntryDeclareCacheEntryWithKnownTicket (DependencyTicket known_ticket, std::string description, CacheEntry::AllocCallback alloc_function, CacheEntry::CalcCallback calc_function, std::set< DependencyTicket > prerequisites_of_calc={ all_sources_ticket()})
 (Internal use only) This is for cache entries associated with pre-defined tickets, for example the cache entry for time derivatives. More...
 
const internal::SystemParentServiceInterface * get_parent_service () const
 Returns a pointer to the service interface of the immediately enclosing Diagram if one has been set, otherwise nullptr. More...
 
DependencyTicket assign_next_dependency_ticket ()
 (Internal use only) Assigns the next unused dependency ticket number, unique only within a particular system. More...
 
const AbstractValueEvalAbstractInputImpl (const char *func, const ContextBase &context, InputPortIndex port_index) const
 (Internal use only) Shared code for updating an input port and returning a pointer to its abstract value, or nullptr if the port is not connected. More...
 
void ThrowNegativePortIndex (const char *func, int port_index) const
 Throws std::out_of_range to report a negative port_index that was passed to API method func. More...
 
void ThrowInputPortIndexOutOfRange (const char *func, InputPortIndex port_index) const
 Throws std::out_of_range to report bad input port_index that was passed to API method func. More...
 
void ThrowOutputPortIndexOutOfRange (const char *func, OutputPortIndex port_index) const
 Throws std::out_of_range to report bad output port_index that was passed to API method func. More...
 
void ThrowNotAVectorInputPort (const char *func, InputPortIndex port_index) const
 Throws std::logic_error because someone misused API method func, that is only allowed for declared-vector input ports, on an abstract port whose index is given here. More...
 
void ThrowInputPortHasWrongType (const char *func, InputPortIndex port_index, const std::string &expected_type, const std::string &actual_type) const
 Throws std::logic_error because someone called API method func claiming the input port had some value type that was wrong. More...
 
void ThrowCantEvaluateInputPort (const char *func, InputPortIndex port_index) const
 Throws std::logic_error because someone called API method func, that requires this input port to be evaluatable, but the port was neither fixed nor connected. More...
 
const InputPortBaseGetInputPortBaseOrThrow (const char *func, int port_index) const
 (Internal use only) Returns the InputPortBase at index port_index, throwing std::out_of_range we don't like the port index. More...
 
const OutputPortBaseGetOutputPortBaseOrThrow (const char *func, int port_index) const
 (Internal use only) Returns the OutputPortBase at index port_index, throwing std::out_of_range we don't like the port index. More...
 
void InitializeContextBase (ContextBase *context) const
 This method must be invoked from within derived class DoAllocateContext() implementations right after the concrete Context object has been allocated. More...
 
virtual std::unique_ptr< ContextBaseDoAllocateContext () const =0
 Derived class implementations should allocate a suitable concrete Context type, then invoke the above InitializeContextBase() method. More...
 

Friends

class SystemImpl
 

Automatic differentiation

From a System templatized by double, you can obtain an identical system templatized by an automatic differentation scalar providing machine-precision computation of partial derivatives of any numerical result of the System with respect to any of the numerical values that can be contained in a Context (time, inputs, parameters, and state).

std::unique_ptr< System< AutoDiffXd > > ToAutoDiffXd () const
 Creates a deep copy of this System, transmogrified to use the autodiff scalar type, with a dynamic-sized vector of partial derivatives. More...
 
std::unique_ptr< System< AutoDiffXd > > ToAutoDiffXdMaybe () const
 Creates a deep copy of this system exactly like ToAutoDiffXd(), but returns nullptr if this System does not support autodiff, instead of throwing an exception. More...
 
template<template< typename > class S = ::drake::systems::System>
static std::unique_ptr< S< AutoDiffXd > > ToAutoDiffXd (const S< T > &from)
 Creates a deep copy of from, transmogrified to use the autodiff scalar type, with a dynamic-sized vector of partial derivatives. More...
 

Symbolics

From a System templatized by double, you can obtain an identical system templatized by a symbolic expression scalar.

std::unique_ptr< System< symbolic::Expression > > ToSymbolic () const
 Creates a deep copy of this System, transmogrified to use the symbolic scalar type. More...
 
std::unique_ptr< System< symbolic::Expression > > ToSymbolicMaybe () const
 Creates a deep copy of this system exactly like ToSymbolic(), but returns nullptr if this System does not support symbolic, instead of throwing an exception. More...
 
template<template< typename > class S = ::drake::systems::System>
static std::unique_ptr< S< symbolic::Expression > > ToSymbolic (const S< T > &from)
 Creates a deep copy of from, transmogrified to use the symbolic scalar type. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from SystemBase
static DependencyTicket nothing_ticket ()
 Returns a ticket indicating that a computation does not depend on any source value; that is, it is a constant. More...
 
static DependencyTicket time_ticket ()
 Returns a ticket indicating dependence on time. More...
 
static DependencyTicket accuracy_ticket ()
 Returns a ticket indicating dependence on the accuracy setting in the Context. More...
 
static DependencyTicket q_ticket ()
 Returns a ticket indicating that a computation depends on configuration state variables q. More...
 
static DependencyTicket v_ticket ()
 Returns a ticket indicating dependence on velocity state variables v. More...
 
static DependencyTicket z_ticket ()
 Returns a ticket indicating dependence on any or all of the miscellaneous continuous state variables z. More...
 
static DependencyTicket xc_ticket ()
 Returns a ticket indicating dependence on all of the continuous state variables q, v, or z. More...
 
static DependencyTicket xd_ticket ()
 Returns a ticket indicating dependence on all of the numerical discrete state variables, in any discrete variable group. More...
 
static DependencyTicket xa_ticket ()
 Returns a ticket indicating dependence on all of the abstract state variables in the current Context. More...
 
static DependencyTicket all_state_ticket ()
 Returns a ticket indicating dependence on all state variables x in this system, including continuous variables xc, discrete (numeric) variables xd, and abstract state variables xa. More...
 
static DependencyTicket pn_ticket ()
 Returns a ticket indicating dependence on all of the numerical parameters in the current Context. More...
 
static DependencyTicket pa_ticket ()
 Returns a ticket indicating dependence on all of the abstract parameters pa in the current Context. More...
 
static DependencyTicket all_parameters_ticket ()
 Returns a ticket indicating dependence on all parameters p in this system, including numeric parameters pn, and abstract parameters pa. More...
 
static DependencyTicket all_input_ports_ticket ()
 Returns a ticket indicating dependence on all input ports u of this system. More...
 
static DependencyTicket all_sources_ticket ()
 Returns a ticket indicating dependence on every possible independent source value, including time, accuracy, state, input ports, and parameters (but not cache entries). More...
 
static DependencyTicket configuration_ticket ()
 Returns a ticket indicating dependence on all source values that may affect configuration-dependent computations. More...
 
static DependencyTicket kinematics_ticket ()
 Returns a ticket indicating dependence on all source values that may affect configuration- or velocity-dependent computations. More...
 
static DependencyTicket xcdot_ticket ()
 Returns a ticket for the cache entry that holds time derivatives of the continuous variables. More...
 
static DependencyTicket pe_ticket ()
 Returns a ticket for the cache entry that holds the potential energy calculation. More...
 
static DependencyTicket ke_ticket ()
 Returns a ticket for the cache entry that holds the kinetic energy calculation. More...
 
static DependencyTicket pc_ticket ()
 Returns a ticket for the cache entry that holds the conservative power calculation. More...
 
static DependencyTicket pnc_ticket ()
 Returns a ticket for the cache entry that holds the non-conservative power calculation. More...
 
- Static Protected Member Functions inherited from SystemBase
static void set_parent_service (SystemBase *child, const internal::SystemParentServiceInterface *parent_service)
 (Internal use only) Declares that parent_service is the service interface of the Diagram that owns this subsystem. More...
 
static void ThrowInputPortHasWrongType (const char *func, const std::string &system_pathname, InputPortIndex, const std::string &expected_type, const std::string &actual_type)
 Throws std::logic_error because someone called API method func claiming the input port had some value type that was wrong. More...
 

Detailed Description

template<typename T>
class drake::systems::System< T >

Base class for all System functionality that is dependent on the templatized scalar type T for input, state, parameters, and outputs.

Template Parameters
TThe vector element type, which must be a valid Eigen scalar.

Constructor & Destructor Documentation

◆ System() [1/3]

System ( const System< T > &  )
delete

◆ System() [2/3]

System ( System< T > &&  )
delete

◆ ~System()

~System ( )
overridedefault

◆ System() [3/3]

System ( SystemScalarConverter  converter)
inlineexplicitprotected

Constructs an empty System base class object and allocates base class resources, possibly supporting scalar-type conversion support (AutoDiff, etc.) using converter.

See System Scalar Conversion for detailed background and examples related to scalar-type conversion support.

Member Function Documentation

◆ AddConstraint()

SystemConstraintIndex AddConstraint ( std::unique_ptr< SystemConstraint< T >>  constraint)
inlineprotected

Adds an already-created constraint to the list of constraints for this System.

Ownership of the SystemConstraint is transferred to this system.

◆ AddTriggeredWitnessFunctionToCompositeEventCollection()

virtual void AddTriggeredWitnessFunctionToCompositeEventCollection ( Event< T > *  event,
CompositeEventCollection< T > *  events 
) const
pure virtual

Add event to events due to a witness function triggering.

events should be allocated with this system's AllocateCompositeEventCollection. Neither event nor events can be nullptr. Additionally, event must contain event data (event->get_event_data() must not be nullptr) and the type of that data must be WitnessTriggeredEventData.

Implemented in Diagram< T >, Diagram< double >, LeafSystem< T >, LeafSystem< double >, and LeafSystem< AutoDiffXd >.

◆ AllocateCompositeEventCollection()

virtual std::unique_ptr<CompositeEventCollection<T> > AllocateCompositeEventCollection ( ) const
pure virtual

Allocates a CompositeEventCollection for this system.

The allocated instance is used for registering events; for example, Simulator passes this object to System::CalcNextUpdateTime() to allow the system to register upcoming events.

Implemented in LeafSystem< T >, LeafSystem< double >, LeafSystem< AutoDiffXd >, Diagram< T >, and Diagram< double >.

◆ AllocateContext()

std::unique_ptr<Context<T> > AllocateContext ( ) const
inline

Returns a Context<T> suitable for use with this System<T>.

◆ AllocateDiscreteVariables()

virtual std::unique_ptr<DiscreteValues<T> > AllocateDiscreteVariables ( ) const
inlinevirtual

Returns a DiscreteState of the same dimensions as the discrete_state allocated in CreateDefaultContext.

The simulator will provide this state as the output argument to Update. By default, allocates nothing. Systems with discrete state variables should override.

Reimplemented in Diagram< T >, Diagram< double >, LeafSystem< T >, LeafSystem< double >, and LeafSystem< AutoDiffXd >.

◆ AllocateFixedInputs()

void AllocateFixedInputs ( Context< T > *  context) const
inline

For each input port, allocates a fixed input of the concrete type that this System requires, and binds it to the port, disconnecting any prior input.

Does not assign any values to the fixed inputs.

◆ AllocateInputAbstract()

std::unique_ptr<AbstractValue> AllocateInputAbstract ( const InputPort< T > &  input_port) const
inline

Given an input port, allocates the abstract storage.

The input_port must match a port declared via DeclareInputPort.

◆ AllocateInputVector()

std::unique_ptr<BasicVector<T> > AllocateInputVector ( const InputPort< T > &  input_port) const
inline

Given an input port, allocates the vector storage.

The input_port must match a port declared via DeclareInputPort.

◆ AllocateOutput()

std::unique_ptr<SystemOutput<T> > AllocateOutput ( ) const
inline

Returns a container that can hold the values of all of this System's output ports.

It is sized with the number of output ports and uses each output port's allocation method to provide an object of the right type for that port.

◆ AllocateTimeDerivatives()

virtual std::unique_ptr<ContinuousState<T> > AllocateTimeDerivatives ( ) const
inlinevirtual

Returns a ContinuousState of the same size as the continuous_state allocated in CreateDefaultContext.

The simulator will provide this state as the output argument to EvalTimeDerivatives.

By default, allocates no derivatives. Systems with continuous state variables should override.

Reimplemented in Diagram< T >, Diagram< double >, LeafSystem< T >, LeafSystem< double >, and LeafSystem< AutoDiffXd >.

◆ CalcConservativePower()

T CalcConservativePower ( const Context< T > &  context) const
inline

Calculates and returns the conservative power represented by the current contents of the given context.

Prefer EvalConservativePower() to avoid unnecessary recalculation.

See also
EvalConservativePower() for more information.

◆ CalcConstraintErrorNorm()

double CalcConstraintErrorNorm ( const Context< T > &  context,
const Eigen::VectorXd &  error 
) const
inline

Computes the norm on constraint error (used as a metric for comparing errors between the outputs of algebraic equations applied to two different state variable instances).

This norm need be neither continuous nor differentiable.

Exceptions
std::logic_errorif the dimension of err is not equivalent to the output of get_num_constraint_equations().

◆ CalcDiscreteVariableUpdates() [1/2]

void CalcDiscreteVariableUpdates ( const Context< T > &  context,
const EventCollection< DiscreteUpdateEvent< T >> &  events,
DiscreteValues< T > *  discrete_state 
) const
inline

This method is the public entry point for dispatching all discrete variable update event handlers.

Using all the discrete update handlers in events, the method calculates the update xd(n+1) to discrete variables xd(n) in context and outputs the results to discrete_state. See documentation for DispatchDiscreteVariableUpdateHandler() for more details.

◆ CalcDiscreteVariableUpdates() [2/2]

void CalcDiscreteVariableUpdates ( const Context< T > &  context,
DiscreteValues< T > *  discrete_state 
) const
inline

This method forces a discrete update on the system given a context, and the updated discrete state is stored in discrete_state.

The discrete update event will have a trigger type of kForced, with no attribute or custom callback.

◆ CalcKineticEnergy()

T CalcKineticEnergy ( const Context< T > &  context) const
inline

Calculates and returns the kinetic energy represented by the current configuration and velocity provided in context.

Prefer EvalKineticEnergy() to avoid unnecessary recalculation.

See also
EvalKineticEnergy() for more information.

◆ CalcNextUpdateTime()

T CalcNextUpdateTime ( const Context< T > &  context,
CompositeEventCollection< T > *  events 
) const
inline

This method is called by a Simulator during its calculation of the size of the next continuous step to attempt.

The System returns the next time at which some discrete action must be taken, and records what those actions ought to be in events. Upon reaching that time, the simulator will merge events with the other CompositeEventCollection instances scheduled through mechanisms (e.g. GetPerStepEvents()), and the merged CompositeEventCollection will be passed to all event handling mechanisms.

events cannot be null. events will be cleared on entry.

◆ CalcNonConservativePower()

T CalcNonConservativePower ( const Context< T > &  context) const
inline

Calculates and returns the non-conservative power represented by the current contents of the given context.

Prefer EvalNonConservativePower() to avoid unnecessary recalculation.

See also
EvalNonConservativePower() for more information.

◆ CalcOutput()

void CalcOutput ( const Context< T > &  context,
SystemOutput< T > *  outputs 
) const
inline

Utility method that computes for every output port i the value y(i) that should result from the current contents of the given Context.

Note that individual output port values can be calculated using get_output_port(i).Calc(); this method invokes that for each output port in index order. The result may depend on time and the current values of input ports, parameters, and state variables. The result is written to outputs which must already have been allocated to have the right number of entries of the right types.

◆ CalcPotentialEnergy()

T CalcPotentialEnergy ( const Context< T > &  context) const
inline

Calculates and returns the potential energy represented by the current configuration provided in context.

Prefer EvalPotentialEnergy() to avoid unnecessary recalculation.

See also
EvalPotentialEnergy() for more information.

◆ CalcTimeDerivatives()

void CalcTimeDerivatives ( const Context< T > &  context,
ContinuousState< T > *  derivatives 
) const
inline

Calculates the time derivatives xcdot of the continuous state xc into a given output argument.

Prefer EvalTimeDerivatives() instead to avoid unnecessary recomputation.

See also
EvalTimeDerivatives() for more information.
Parameters
contextThe Context whose contents will be used to evaluate the derivatives.
derivativesThe time derivatives xcdot. Must be the same size as the continuous state vector in context.

◆ CalcUnrestrictedUpdate() [1/2]

void CalcUnrestrictedUpdate ( const Context< T > &  context,
const EventCollection< UnrestrictedUpdateEvent< T >> &  events,
State< T > *  state 
) const
inline

This method is the public entry point for dispatching all unrestricted update event handlers.

Using all the unrestricted update handers in events, it updates any state variables in the context, and outputs the results to state. It does not allow the dimensionality of the state variables to change. See the documentation for DispatchUnrestrictedUpdateHandler() for more details.

Exceptions
std::logic_errorif the dimensionality of the state variables changes in the callback.

◆ CalcUnrestrictedUpdate() [2/2]

void CalcUnrestrictedUpdate ( const Context< T > &  context,
State< T > *  state 
) const
inline

This method forces an unrestricted update on the system given a context, and the updated state is stored in state.

The unrestricted update event will have a trigger type of kForced, with no additional data, attribute or custom callback.

See also
CalcUnrestrictedUpdate(const Context<T>&, const EventCollection<UnrestrictedUpdateEvent<T>>*, State<T>* state) for more information.

◆ CalcVelocityChangeFromConstraintImpulses()

Eigen::VectorXd CalcVelocityChangeFromConstraintImpulses ( const Context< T > &  context,
const Eigen::MatrixXd &  J,
const Eigen::VectorXd &  lambda 
) const
inline

Computes the change in velocity from applying the given constraint forces to the system at the given context.

Parameters
contextthe current system state, provision of which also yields the ability of the constraints to be dependent upon the current system state (as might be the case with a piecewise differential algebraic equation).
Ja m × n constraint Jacobian matrix of the m constraint equations g() differentiated with respect to the n configuration variables q (i.e., J should be ∂g/∂q). If the time derivatives of the generalized coordinates of the system are not identical to the generalized velocity (in general they need not be, e.g., if generalized coordinates use unit unit quaternions to represent 3D orientation), J should instead be defined as ∂g/∂q⋅N, where N ≡ ∂q/∂ꝗ is the Jacobian matrix (dependent on q) of the generalized coordinates with respect to the quasi-coordinates (ꝗ, pronounced "qbar", where dꝗ/dt are the generalized velocities).
lambdathe vector of constraint forces (of same dimension as the number of rows in the Jacobian matrix, J)
Returns
a n dimensional vector, where n is the dimension of the quasi-coordinates.

◆ CalcWitnessValue()

T CalcWitnessValue ( const Context< T > &  context,
const WitnessFunction< T > &  witness_func 
) const
inline

Evaluates a witness function at the given context.

◆ CheckSystemConstraintsSatisfied()

boolean<T> CheckSystemConstraintsSatisfied ( const Context< T > &  context,
double  tol 
) const
inline

Returns true if context satisfies all of the registered SystemConstraints with tolerance tol.

See also
SystemConstraint::CheckSatisfied.

◆ CheckValidContextT()

void CheckValidContextT ( const Context< T1 > &  context) const
inline

Checks that context is consistent for this System template.

Supports any scalar type, but expects T by default.

Exceptions
std::exceptionunless context is valid for this system.
Template Parameters
T1the scalar type of the Context to check.

◆ CheckValidOutput()

void CheckValidOutput ( const SystemOutput< T > *  output) const
inline

Checks that output is consistent with the number and size of output ports declared by the system.

Exceptions
std::exceptionunless output is non-null and valid for this system.

◆ CopyContinuousStateVector()

VectorX<T> CopyContinuousStateVector ( const Context< T > &  context) const
inline

Returns a copy of the continuous state vector xc into an Eigen vector.

◆ CreateDefaultContext()

std::unique_ptr<Context<T> > CreateDefaultContext ( ) const
inline

This convenience method allocates a context using AllocateContext() and sets its default values using SetDefaultContext().

◆ DeclareInputPort() [1/2]

const InputPort<T>& DeclareInputPort ( variant< std::string, UseDefaultName name,
PortDataType  type,
int  size,
optional< RandomDistribution random_type = nullopt 
)
inlineprotected

Adds a port with the specified type and size to the input topology.

Input port names must be unique for this system (passing in a duplicate name will throw std::logic_error). If name is given as kUseDefaultName, then a default value of e.g. "u2", where 2 is the input number will be provided. An empty name is not permitted.

If the port is intended to model a random noise or disturbance input, random_type can (optionally) be used to label it as such; doing so enables algorithms for design and analysis (e.g. state estimation) to reason explicitly about randomness at the system level. All random input ports are assumed to be statistically independent.

Precondition
name must not be empty.
Exceptions
std::logic_errorfor a duplicate port name.
Returns
the declared port.

◆ DeclareInputPort() [2/2]

const InputPort<T>& DeclareInputPort ( PortDataType  type,
int  size,
optional< RandomDistribution random_type = nullopt 
)
inlineprotected

See the nearly identical signature with an additional (first) argument specifying the port name.

This version will be deprecated as discussed in #9447.

◆ DispatchDiscreteVariableUpdateHandler()

virtual void DispatchDiscreteVariableUpdateHandler ( const Context< T > &  context,
const EventCollection< DiscreteUpdateEvent< T >> &  events,
DiscreteValues< T > *  discrete_state 
) const
protectedpure virtual

This function dispatches all discrete update events to the appropriate handlers.

discrete_state cannot be null.

◆ DispatchPublishHandler()

virtual void DispatchPublishHandler ( const Context< T > &  context,
const EventCollection< PublishEvent< T >> &  events 
) const
protectedpure virtual

This function dispatches all publish events to the appropriate handlers.

◆ DispatchUnrestrictedUpdateHandler()

virtual void DispatchUnrestrictedUpdateHandler ( const Context< T > &  context,
const EventCollection< UnrestrictedUpdateEvent< T >> &  events,
State< T > *  state 
) const
protectedpure virtual

This function dispatches all unrestricted update events to the appropriate handlers.

state cannot be null.

◆ do_get_num_constraint_equations()

virtual int do_get_num_constraint_equations ( const Context< T > &  context) const
inlineprotectedvirtual

Gets the number of constraint equations for this system from the given context.

The context is supplied in case the number of constraints is dependent upon the current state (as might be the case with a piecewise differential algebraic equation). Derived classes can override this function, which is called by get_num_constraint_equations().

See also
get_num_constraint_equations() for parameter documentation.
Returns
zero by default

◆ DoCalcConservativePower()

virtual T DoCalcConservativePower ( const Context< T > &  context) const
inlineprotectedvirtual

Override this method to return the rate Pc at which mechanical energy is being converted from potential energy to kinetic energy by this system in the given Context.

By default, returns zero. Physical systems should override. You may assume that context has already been validated before it is passed to you here.

See EvalConservativePower() for details on what you must compute here. In particular, this quantity must be positive when potential energy is decreasing, and your conservative power method must not depend explicitly on time or any input port values.

Reimplemented in RigidBodyPlant< T >, RigidBodyPlant< double >, SpringMassSystem< T >, and SpringMassSystem< double >.

◆ DoCalcConstraintErrorNorm()

virtual double DoCalcConstraintErrorNorm ( const Context< T > &  context,
const Eigen::VectorXd &  error 
) const
inlineprotectedvirtual

Computes the norm of the constraint error.

This default implementation computes a Euclidean norm of the error. Derived classes can override this function, which is called by CalcConstraintErrorNorm(). This norm need be neither continuous nor differentiable.

See also
CalcConstraintErrorNorm() for parameter documentation.

◆ DoCalcKineticEnergy()

virtual T DoCalcKineticEnergy ( const Context< T > &  context) const
inlineprotectedvirtual

Override this method for physical systems to calculate the kinetic energy KE currently present in the motion provided in the given Context.

The default implementation returns 0 which is correct for non-physical systems. You may assume that context has already been validated before it is passed to you here.

See EvalKineticEnergy() for details on what you must compute here. In particular, your kinetic energy method must not depend explicitly on time or any input port values.

Reimplemented in RigidBodyPlant< T >, RigidBodyPlant< double >, SpringMassSystem< T >, and SpringMassSystem< double >.

◆ DoCalcNextUpdateTime()

virtual void DoCalcNextUpdateTime ( const Context< T > &  context,
CompositeEventCollection< T > *  events,
T *  time 
) const
inlineprotectedvirtual

Computes the next time at which this System must perform a discrete action.

Override this method if your System has any discrete actions which must interrupt the continuous simulation. This method is called only from the public non-virtual CalcNextUpdateTime() which will already have error-checked the parameters so you don't have to. You may assume that context has already been validated and events pointer is not null.

The default implementation returns with the next sample time being Infinity and no events added to events.

Reimplemented in Diagram< T >, Diagram< double >, LeafSystem< T >, LeafSystem< double >, LeafSystem< AutoDiffXd >, LcmSubscriberSystem, and LcmLogPlaybackSystem.

◆ DoCalcNonConservativePower()

virtual T DoCalcNonConservativePower ( const Context< T > &  context) const
inlineprotectedvirtual

Override this method to return the rate Pnc at which work W is done on the system by non-conservative forces.

By default, returns zero. Physical systems should override. You may assume that context has already been validated before it is passed to you here.

See EvalNonConservativePower() for details on what you must compute here. In particular, this quantity must be negative if the non-conservative forces are dissipative, positive otherwise. Your non-conservative power method can depend on anything you find in the given Context, including time and input ports.

Reimplemented in RigidBodyPlant< T >, RigidBodyPlant< double >, SpringMassSystem< T >, and SpringMassSystem< double >.

◆ DoCalcPotentialEnergy()

virtual T DoCalcPotentialEnergy ( const Context< T > &  context) const
inlineprotectedvirtual

Override this method for physical systems to calculate the potential energy PE currently stored in the configuration provided in the given Context.

The default implementation returns 0 which is correct for non-physical systems. You may assume that context has already been validated before it is passed to you here.

See EvalPotentialEnergy() for details on what you must compute here. In particular, your potential energy method must not depend explicitly on time, velocities, or any input port values.

Reimplemented in RigidBodyPlant< T >, RigidBodyPlant< double >, SpringMassSystem< T >, and SpringMassSystem< double >.

◆ DoCalcTimeDerivatives()

virtual void DoCalcTimeDerivatives ( const Context< T > &  context,
ContinuousState< T > *  derivatives 
) const
inlineprotectedvirtual

Override this if you have any continuous state variables xc in your concrete System to calculate their time derivatives.

The derivatives vector will correspond elementwise with the state vector Context.state.continuous_state.get_state(). Thus, if the state in the Context has second-order structure xc=[q,v,z], that same structure applies to the derivatives.

This method is called only from the public non-virtual CalcTimeDerivatives() which will already have error-checked the parameters so you don't have to. In particular, implementations may assume that the given Context is valid for this System; that the derivatives pointer is non-null, and that the referenced object has the same constituent structure as was produced by AllocateTimeDerivatives().

The default implementation does nothing if the derivatives vector is size zero and aborts otherwise.

Reimplemented in RigidBodyPlant< T >, RigidBodyPlant< double >, Diagram< T >, Diagram< double >, SpringMassSystem< T >, SpringMassSystem< double >, PidController< T >, VectorSystem< T >, VectorSystem< double >, SpringMassDamperSystem< T >, TimeVaryingAffineSystem< T >, TimeVaryingAffineSystem< double >, StiffDoubleMassSpringSystem< T >, DiscontinuousSpringMassDamperSystem< T >, ODESystem< T >, and RobertsonSystem< T >.

◆ DoCalcVelocityChangeFromConstraintImpulses()

virtual Eigen::VectorXd DoCalcVelocityChangeFromConstraintImpulses ( const Context< T > &  context,
const Eigen::MatrixXd &  J,
const Eigen::VectorXd &  lambda 
) const
inlineprotectedvirtual

Computes the change in velocity from applying the given constraint forces to the system at the given context.

Derived classes can override this function, which is called by CalcVelocityChangeFromConstraintImpulses().

Returns
the zero vector of dimension of the dimension of the quasi-coordinates, by default.
See also
CalcVelocityChangeFromConstraintImpulses() for parameter documentation.

◆ DoCalcWitnessValue()

virtual T DoCalcWitnessValue ( const Context< T > &  context,
const WitnessFunction< T > &  witness_func 
) const
protectedpure virtual

Derived classes will implement this method to evaluate a witness function at the given context.

Implemented in Diagram< T >, Diagram< double >, LeafSystem< T >, LeafSystem< double >, and LeafSystem< AutoDiffXd >.

◆ DoEvalConstraintEquations()

virtual Eigen::VectorXd DoEvalConstraintEquations ( const Context< T > &  context) const
inlineprotectedvirtual

Evaluates the constraint equations for the system at the generalized coordinates and generalized velocity specified by the context.

The context allows the set of constraints to be dependent upon the current system state (as might be the case with a piecewise differential algebraic equation). The default implementation of this function returns a zero-dimensional vector. Derived classes can override this function, which is called by EvalConstraintEquations().

See also
EvalConstraintEquations() for parameter documentation.
Returns
a vector of dimension get_num_constraint_equations(); the zero vector indicates that the algebraic constraints are all satisfied.

◆ DoEvalConstraintEquationsDot()

virtual Eigen::VectorXd DoEvalConstraintEquationsDot ( const Context< T > &  context) const
inlineprotectedvirtual

Computes the time derivative of each constraint equation, evaluated at the generalized coordinates and generalized velocity specified by the context.

The context allows the set of constraints to be dependent upon the current system state (as might be the case with a piecewise differential algebraic equation). The default implementation of this function returns a zero-dimensional vector. Derived classes can override this function, which is called by EvalConstraintEquationsDot().

Returns
a vector of dimension get_num_constraint_equations().
See also
EvalConstraintEquationsDot() for parameter documentation.

◆ DoGetInitializationEvents()

virtual void DoGetInitializationEvents ( const Context< T > &  context,
CompositeEventCollection< T > *  events 
) const
inlineprotectedvirtual

Implement this method to return any events to be handled at the simulator's initialization step.

events is cleared in the public non-virtual GetInitializationEvents(). You may assume that context has already been validated and that events is not null. events can be changed freely by the overriding implementation.

The default implementation returns without changing events.

See also
GetInitializationEvents()

◆ DoGetPeriodicEvents()

virtual std::map<PeriodicEventData, std::vector<const Event<T>*>, PeriodicEventDataComparator> DoGetPeriodicEvents ( ) const
protectedpure virtual

Implement this method to return all periodic triggered events.

See also
GetPeriodicEvents() for a detailed description of the returned variable.
Note
The default implementation returns an empty map.

◆ DoGetPerStepEvents()

virtual void DoGetPerStepEvents ( const Context< T > &  context,
CompositeEventCollection< T > *  events 
) const
inlineprotectedvirtual

Implement this method to return any events to be handled before the simulator integrates the system's continuous state at each time step.

events is cleared in the public non-virtual GetPerStepEvents() before that method calls this function. You may assume that context has already been validated and that events is not null. events can be changed freely by the overriding implementation.

The default implementation returns without changing events.

See also
GetPerStepEvents()

◆ DoGetWitnessFunctions()

virtual void DoGetWitnessFunctions ( const Context< T > &  ,
std::vector< const WitnessFunction< T > *> *   
) const
inlineprotectedvirtual

Derived classes can override this method to provide witness functions active for the given state.

The default implementation does nothing. On entry to this function, the context will have already been validated and the vector of witness functions will have been validated to be both empty and non-null.

Reimplemented in Diagram< T >, Diagram< double >, and StatelessSystem< T >.

◆ DoMapQDotToVelocity()

virtual void DoMapQDotToVelocity ( const Context< T > &  context,
const Eigen::Ref< const VectorX< T >> &  qdot,
VectorBase< T > *  generalized_velocity 
) const
inlineprotectedvirtual

Provides the substantive implementation of MapQDotToVelocity().

The default implementation uses the identity mapping, and correctly does nothing if the System does not have second-order state variables. It throws std::runtime_error if the generalized_velocity and qdot are not the same size, but that is not enough to guarantee that the default implementation is adequate. Child classes must override this function if qdot != v (even if they are the same size). This occurs, for example, if a joint uses roll-pitch-yaw rotation angles for orientation but angular velocity for rotational rate rather than rotation angle derivatives.

If you implement this method you are required to use no more than O(nq) time where nq is the size of qdot, so that the System can meet the performance guarantee made for the public interface, and you must also implement DoMapVelocityToQDot(). Implementations may assume that qdot has already been validated to be the same size as q in the given Context, and that generalized_velocity is non-null.

Reimplemented in Diagram< T >, Diagram< double >, RigidBodyPlant< T >, and RigidBodyPlant< double >.

◆ DoMapVelocityToQDot()

virtual void DoMapVelocityToQDot ( const Context< T > &  context,
const Eigen::Ref< const VectorX< T >> &  generalized_velocity,
VectorBase< T > *  qdot 
) const
inlineprotectedvirtual

Provides the substantive implementation of MapVelocityToQDot().

The default implementation uses the identity mapping, and correctly does nothing if the System does not have second-order state variables. It throws std::runtime_error if the generalized_velocity (v) and qdot are not the same size, but that is not enough to guarantee that the default implementation is adequate. Child classes must override this function if qdot != v (even if they are the same size). This occurs, for example, if a joint uses roll-pitch-yaw rotation angles for orientation but angular velocity for rotational rate rather than rotation angle derivatives.

If you implement this method you are required to use no more than O(nq) time where nq is the size of qdot, so that the System can meet the performance guarantee made for the public interface, and you must also implement DoMapQDotToVelocity(). Implementations may assume that generalized_velocity has already been validated to be the same size as v in the given Context, and that qdot is non-null.

Reimplemented in Diagram< T >, Diagram< double >, RigidBodyPlant< T >, and RigidBodyPlant< double >.

◆ EvalConservativePower()

const T& EvalConservativePower ( const Context< T > &  context) const
inline

Returns a reference to the cached value of the conservative power (Pc), evaluating first if necessary using CalcConservativePower().

The returned Pc represents the rate at which mechanical energy is being converted from potential energy (PE) to kinetic energy (KE) by this system in the given Context. This quantity will be positive when PE is decreasing. By definition here, conservative power may depend only on quantities that explicitly contribute to PE and KE. See EvalPotentialEnergy() and EvalKineticEnergy() for details.

Power due to non-conservative forces (e.g. dampers) can contribute to the rate of change of KE. Therefore this method alone cannot be used to determine whether KE is increasing or decreasing, only whether the conservative power is adding or removing kinetic energy. EvalNonConservativePower() can be used in conjunction with this method to find the total rate of change of KE.

Non-physical systems where Pc is not meaningful will return Pc = 0.

Parameters
contextThe Context whose contents may be used to evaluate conservative power.
Return values
PcThe conservative power in watts (W or J/s) represented by the contents of the given context.
See also
CalcConservativePower(), EvalNonConservativePower(), EvalPotentialEnergy(), EvalKineticEnergy()

◆ EvalConstraintEquations()

Eigen::VectorXd EvalConstraintEquations ( const Context< T > &  context) const
inline

Evaluates the constraint equations for the system at the generalized coordinates and generalized velocity specified by the context.

The context allows the set of constraints to be dependent upon the current system state (as might be the case with a system modeled using piecewise differential algebraic equations).

Returns
a vector of dimension get_num_constraint_equations(); the zero vector indicates that the algebraic constraints are all satisfied.

◆ EvalConstraintEquationsDot()

Eigen::VectorXd EvalConstraintEquationsDot ( const Context< T > &  context) const
inline

Computes the time derivative of each constraint equation, evaluated at the generalized coordinates and generalized velocity specified by the context.

The context allows the set of constraints to be dependent upon the current system state (as might be the case with a system modeled using piecewise differential algebraic equations).

Returns
a vector of dimension get_num_constraint_equations().

◆ EvalEigenVectorInput()

Eigen::VectorBlock<const VectorX<T> > EvalEigenVectorInput ( const Context< T > &  context,
int  port_index 
) const
inline

Returns the value of the vector-valued input port with the given port_index as an Eigen vector.

Causes the value to become up to date first if necessary. See EvalAbstractInput() for more information.

Precondition
port_index selects an existing input port of this System.
the port must have been declared to be vector-valued.
the port must be evaluable (connected or fixed).
See also
EvalVectorInput()

◆ EvalKineticEnergy()

const T& EvalKineticEnergy ( const Context< T > &  context) const
inline

Returns a reference to the cached value of the kinetic energy (KE), evaluating first if necessary using CalcKineticEnergy().

By definition here, kinetic energy depends only on "configuration" and "velocity" (e.g. angular and translational velocity) of moving masses which includes a subset of the state variables, and parameters that affect configuration, velocities, or mass properties. The calculated value may also be affected by the accuracy value supplied in the Context. KE cannot depend explicitly on time (∂KE/∂t = 0) or input port values (∂KE/∂u = 0).

Non-physical systems where KE is not meaningful will return KE = 0.

Parameters
contextThe Context whose configuration and velocity variables may be used to evaluate kinetic energy.
Return values
KEThe kinetic energy in joules (J) represented by the configuration and velocity given in context.
See also
CalcKineticEnergy()

◆ EvalNonConservativePower()

const T& EvalNonConservativePower ( const Context< T > &  context) const
inline

Returns a reference to the cached value of the non-conservative power (Pnc), evaluating first if necessary using CalcNonConservativePower().

The returned Pnc represents the rate at which work W is done on the system by non-conservative forces. Pnc is negative if the non-conservative forces are dissipative, positive otherwise. Time integration of Pnc yields work W, and the total mechanical energy E = PE + KE − W should be conserved by any physically-correct model, to within integration accuracy of W. Power is in watts (J/s). (Watts are abbreviated W but not to be confused with work!) Any values in the supplied Context (including time and input ports) may contribute to the computation of non-conservative power.

Non-physical systems where Pnc is not meaningful will return Pnc = 0.

Parameters
contextThe Context whose contents may be used to evaluate non-conservative power.
Return values
PncThe non-conservative power in watts (W or J/s) represented by the contents of the given context.
See also
CalcNonConservativePower(), EvalConservativePower()

◆ EvalPotentialEnergy()

const T& EvalPotentialEnergy ( const Context< T > &  context) const
inline

Returns a reference to the cached value of the potential energy (PE), evaluating first if necessary using CalcPotentialEnergy().

By definition here, potential energy depends only on "configuration" (e.g. orientation and position), which includes a subset of the state variables, and parameters that affect configuration or conservative forces (such as lengths and masses). The calculated value may also be affected by the accuracy value supplied in the Context. PE cannot depend explicitly on time (∂PE/∂t = 0), velocities (∂PE/∂v = 0), or input port values (∂PE/∂u = 0).

Non-physical systems where PE is not meaningful will return PE = 0.

Parameters
contextThe Context whose configuration variables may be used to evaluate potential energy.
Return values
PEThe potential energy in joules (J) represented by the configuration given in context.
See also
CalcPotentialEnergy()

◆ EvalTimeDerivatives()

const ContinuousState<T>& EvalTimeDerivatives ( const Context< T > &  context) const
inline

Returns a reference to the cached value of the continuous state variable time derivatives, evaluating first if necessary using CalcTimeDerivatives().

This method returns the time derivatives xcdot of the continuous state xc. The referenced return object will correspond elementwise with the continuous state in the given Context. Thus, if the state in the Context has second-order structure xc=[q v z], that same structure applies to the derivatives so we will have xcdot=[qdot vdot zdot].

Parameters
contextThe Context whose time, input port, parameter, state, and accuracy values may be used to evaluate the derivatives.
Return values
xcdotThe time derivatives of xc returned as a reference to an object of the same type and size as this Context's continuous state.
See also
CalcTimeDerivatives()

◆ EvalVectorInput()

const Vec<T>* EvalVectorInput ( const Context< T > &  context,
int  port_index 
) const
inline

Returns the value of the vector-valued input port with the given port_index as a BasicVector or a specific subclass Vec derived from BasicVector.

Causes the value to become up to date first if necessary. See EvalAbstractInput() for more information.

The result is returned as a pointer to the input port's value of type Vec<T> or nullptr if the port is not connected.

Precondition
port_index selects an existing input port of this System.
the port must have been declared to be vector-valued.
the port's value must be of type Vec<T>.
Template Parameters
VecThe template type of the input vector, which must be a subclass of BasicVector.

◆ FixInputPortsFrom()

void FixInputPortsFrom ( const System< double > &  other_system,
const Context< double > &  other_context,
Context< T > *  target_context 
) const
inline

Fixes all of the input ports in target_context to their current values in other_context, as evaluated by other_system.

Exceptions
std::exceptionunless other_context and target_context both have the same shape as this System, and the other_system. Ignores disconnected inputs.

◆ get_constraint()

const SystemConstraint<T>& get_constraint ( SystemConstraintIndex  constraint_index) const
inline

Returns the constraint at index constraint_index.

Exceptions
std::out_of_rangefor an invalid constraint_index.

◆ get_forced_discrete_update_events()

const EventCollection<DiscreteUpdateEvent<T> >& get_forced_discrete_update_events ( ) const
inlineprotected

◆ get_forced_publish_events()

const EventCollection<PublishEvent<T> >& get_forced_publish_events ( ) const
inlineprotected

◆ get_forced_unrestricted_update_events()

const EventCollection<UnrestrictedUpdateEvent<T> >& get_forced_unrestricted_update_events ( ) const
inlineprotected

◆ get_input_port()

const InputPort<T>& get_input_port ( int  port_index) const
inline

Returns the typed input port at index port_index.

◆ get_num_constraint_equations()

int get_num_constraint_equations ( const Context< T > &  context) const
inline

Gets the number of constraint equations for this system using the given context (useful in case the number of constraints is dependent upon the current state (as might be the case with a system modeled using piecewise differential algebraic equations).

◆ get_num_constraints()

int get_num_constraints ( ) const
inline

Returns the number of constraints specified for the system.

◆ get_num_continuous_states()

virtual int get_num_continuous_states ( ) const
pure virtual

Returns the dimension of the continuous state vector that has been declared until now.

Implemented in Diagram< T >, Diagram< double >, LeafSystem< T >, LeafSystem< double >, and LeafSystem< AutoDiffXd >.

◆ get_output_port()

const OutputPort<T>& get_output_port ( int  port_index) const
inline

Returns the typed output port at index port_index.

◆ get_system_scalar_converter()

const SystemScalarConverter& get_system_scalar_converter ( ) const
inline

(Advanced) Returns the SystemScalarConverter for this object.

This is an expert-level API intended for framework authors. Most users should prefer the convenience helpers such as System::ToAutoDiffXd.

◆ GetDirectFeedthroughs()

virtual std::multimap<int, int> GetDirectFeedthroughs ( ) const
pure virtual

Reports all direct feedthroughs from input ports to output ports.

For a system with m input ports: I = i₀, i₁, ..., iₘ₋₁, and n output ports, O = o₀, o₁, ..., oₙ₋₁, the return map will contain pairs (u, v) such that

  • 0 ≤ u < m,
  • 0 ≤ v < n,
  • and there might be a direct feedthrough from input iᵤ to each output oᵥ.

Implemented in LeafSystem< T >, LeafSystem< double >, LeafSystem< AutoDiffXd >, Diagram< T >, and Diagram< double >.

◆ GetGraphvizFragment()

virtual void GetGraphvizFragment ( int  max_depth,
std::stringstream *  dot 
) const
inlinevirtual

Appends a Graphviz fragment to the dot stream.

The fragment must be valid Graphviz when wrapped in a digraph or subgraph stanza. Does nothing by default.

Parameters
max_depthSets a limit to the depth of nested diagrams to

Reimplemented in LeafSystem< T >, LeafSystem< double >, LeafSystem< AutoDiffXd >, Diagram< T >, Diagram< double >, and PidController< T >.

◆ GetGraphvizId()

int64_t GetGraphvizId ( ) const
inline

Returns an opaque integer that uniquely identifies this system in the Graphviz output.

◆ GetGraphvizInputPortToken()

virtual void GetGraphvizInputPortToken ( const InputPort< T > &  port,
int  max_depth,
std::stringstream *  dot 
) const
inlinevirtual

Appends a fragment to the dot stream identifying the graphviz node representing port.

Does nothing by default.

Reimplemented in Diagram< T >, Diagram< double >, LeafSystem< T >, LeafSystem< double >, and LeafSystem< AutoDiffXd >.

◆ GetGraphvizOutputPortToken()

virtual void GetGraphvizOutputPortToken ( const OutputPort< T > &  port,
int  max_depth,
std::stringstream *  dot 
) const
inlinevirtual

Appends a fragment to the dot stream identifying the graphviz node representing port.

Does nothing by default.

Reimplemented in Diagram< T >, Diagram< double >, LeafSystem< T >, LeafSystem< double >, and LeafSystem< AutoDiffXd >.

◆ GetGraphvizString()

std::string GetGraphvizString ( int  max_depth = std::numeric_limits<int>::max()) const
inline

Returns a Graphviz string describing this System.

To render the string, use the Graphviz tool, dot. http://www.graphviz.org/Documentation/dotguide.pdf

Parameters
max_depthSets a limit to the depth of nested diagrams to

◆ GetInitializationEvents()

void GetInitializationEvents ( const Context< T > &  context,
CompositeEventCollection< T > *  events 
) const
inline

This method is called by Simulator::Initialize() to gather all update and publish events that need to be handled at initialization before the simulator starts integration.

events cannot be null. events will be cleared on entry.

◆ GetInputPort()

const InputPort<T>& GetInputPort ( const std::string &  port_name) const
inline

Returns the typed input port with the unique name port_name.

The current implementation performs a linear search over strings; prefer get_input_port() when performance is a concern.

Exceptions
std::logic_errorif port_name is not found.

◆ GetMemoryObjectName()

std::string GetMemoryObjectName ( ) const
inline

Returns a name for this System based on a stringification of its type name and memory address.

This is intended for use in diagnostic output and should not be used for behavioral logic, because the stringification of the type name may produce differing results across platforms and because the address can vary from run to run.

◆ GetMutableOutputVector()

Eigen::VectorBlock<VectorX<T> > GetMutableOutputVector ( SystemOutput< T > *  output,
int  port_index 
) const
inlineprotected

Returns a mutable Eigen expression for a vector valued output port with index port_index in this system.

All input ports that directly depend on this output port will be notified that upstream data has changed, and may invalidate cache entries as a result.

◆ GetOutputPort()

const OutputPort<T>& GetOutputPort ( const std::string &  port_name) const
inline

Returns the typed output port with the unique name port_name.

The current implementation performs a linear search over strings; prefer get_output_port() when performance is a concern.

Exceptions
std::logic_errorif port_name is not found.

◆ GetPeriodicEvents()

std::map<PeriodicEventData, std::vector<const Event<T>*>, PeriodicEventDataComparator> GetPeriodicEvents ( ) const
inline

Gets all periodic triggered events for a system.

Each periodic attribute (offset and period, in seconds) is mapped to one or more update events that are to be triggered at the proper times.

◆ GetPerStepEvents()

void GetPerStepEvents ( const Context< T > &  context,
CompositeEventCollection< T > *  events 
) const
inline

This method is called by Simulator::Initialize() to gather all update and publish events that are to be handled in StepTo() at the point before Simulator integrates continuous state.

It is assumed that these events remain constant throughout the simulation. The "step" here refers to the major time step taken by the Simulator. During every simulation step, the simulator will merge events with the other CompositeEventCollection instances generated by other types of event triggering mechanism (e.g., CalcNextUpdateTime()), and the merged CompositeEventCollection objects will be passed to the appropriate handlers before Simulator integrates the continuous state.

events cannot be null. events will be cleared on entry.

◆ GetUniquePeriodicDiscreteUpdateAttribute()

optional<PeriodicEventData> GetUniquePeriodicDiscreteUpdateAttribute ( ) const
inline

Gets whether there exists a unique periodic attribute that triggers one or more discrete update events (and, if so, returns that unique periodic attribute).

Thus, this method can be used (1) as a test to determine whether a system's dynamics are at least partially governed by difference equations and (2) to obtain the difference equation update times.

Returns
optional<PeriodicEventData> Contains the periodic trigger attributes if the unique periodic attribute exists, otherwise nullopt.

◆ GetWitnessFunctions()

void GetWitnessFunctions ( const Context< T > &  context,
std::vector< const WitnessFunction< T > *> *  w 
) const
inline

Gets the witness functions active for the given state.

DoGetWitnessFunctions() does the actual work. The vector of active witness functions are expected to change only upon an unrestricted update.

Parameters
contexta valid context for the System (aborts if not true).
[out]wa valid pointer to an empty vector that will store pointers to the witness functions active for the current state. The method aborts if witnesses is null or non-empty.

◆ HasAnyDirectFeedthrough()

bool HasAnyDirectFeedthrough ( ) const
inline

Returns true if any of the inputs to the system might be directly fed through to any of its outputs and false otherwise.

◆ HasDirectFeedthrough() [1/2]

bool HasDirectFeedthrough ( int  output_port) const
inline

Returns true if there might be direct-feedthrough from any input port to the given output_port, and false otherwise.

◆ HasDirectFeedthrough() [2/2]

bool HasDirectFeedthrough ( int  input_port,
int  output_port 
) const
inline

Returns true if there might be direct-feedthrough from the given input_port to the given output_port, and false otherwise.

◆ MapQDotToVelocity() [1/2]

void MapQDotToVelocity ( const Context< T > &  context,
const VectorBase< T > &  qdot,
VectorBase< T > *  generalized_velocity 
) const
inline

Transforms the time derivative qdot of the generalized configuration q to generalized velocities v.

v and qdot are related linearly by qdot = N(q) * v, where N is a block diagonal matrix. For example, in a multibody system there will be one block of N per tree joint. Although N is not necessarily square, its left pseudo-inverse N+ can be used to invert that relationship without residual error, provided that qdot is in the range space of N (that is, if it could have been produced as qdot=N*v for some v). Using the configuration q from the given Context this method calculates v = N+ * qdot (where N+=N+(q)) for a given qdot. This computation requires only O(nq) time where nq is the size of qdot. Note that this method does not take qdot from the Context.

See the alternate signature if you already have qdot in an Eigen VectorX object; this signature will copy the VectorBase into an Eigen object before performing the computation.

See also
MapVelocityToQDot()

◆ MapQDotToVelocity() [2/2]

void MapQDotToVelocity ( const Context< T > &  context,
const Eigen::Ref< const VectorX< T >> &  qdot,
VectorBase< T > *  generalized_velocity 
) const
inline

Transforms the given time derivative qdot of generalized configuration q to generalized velocity v.

This signature takes qdot as an Eigen VectorX object for faster speed. See the other signature of MapQDotToVelocity() for additional information.

◆ MapVelocityToQDot() [1/2]

void MapVelocityToQDot ( const Context< T > &  context,
const VectorBase< T > &  generalized_velocity,
VectorBase< T > *  qdot 
) const
inline

Transforms a given generalized velocity v to the time derivative qdot of the generalized configuration q taken from the supplied Context.

v and qdot are related linearly by qdot = N(q) * v, where N is a block diagonal matrix. For example, in a multibody system there will be one block of N per tree joint. This computation requires only O(nq) time where nq is the size of qdot. Note that v is not taken from the Context; it is given as an argument here.

See the alternate signature if you already have the generalized velocity in an Eigen VectorX object; this signature will copy the VectorBase into an Eigen object before performing the computation.

See also
MapQDotToVelocity()

◆ MapVelocityToQDot() [2/2]

void MapVelocityToQDot ( const Context< T > &  context,
const Eigen::Ref< const VectorX< T >> &  generalized_velocity,
VectorBase< T > *  qdot 
) const
inline

Transforms the given generalized velocity to the time derivative of generalized configuration.

See the other signature of MapVelocityToQDot() for more information.

◆ operator=() [1/2]

System& operator= ( const System< T > &  )
delete

◆ operator=() [2/2]

System& operator= ( System< T > &&  )
delete

◆ Publish() [1/2]

void Publish ( const Context< T > &  context,
const EventCollection< PublishEvent< T >> &  events 
) const
inline

This method is the public entry point for dispatching all publish event handlers.

It checks the validity of context, and directly calls DispatchPublishHandler. events is a homogeneous collection of publish events, which is typically the publish portion of the heterogeneous event collection generated by CalcNextUpdateTime or GetPerStepEvents.

Note
When publishing is scheduled at particular times, those times likely will not coincide with integrator step times. A Simulator may interpolate to generate a suitable Context, or it may adjust the integrator step size so that a step begins exactly at the next publication time. In the latter case the change in step size may affect the numerical result somewhat since a smaller integrator step produces a more accurate solution.

◆ Publish() [2/2]

void Publish ( const Context< T > &  context) const
inline

Forces a publish on the system, given a context.

The publish event will have a trigger type of kForced, with no additional data, attribute or custom callback. The Simulator can be configured to call this in Simulator::Initialize() and at the start of each continuous integration step. See the Simulator API for more details.

◆ set_forced_discrete_update_events()

void set_forced_discrete_update_events ( std::unique_ptr< EventCollection< DiscreteUpdateEvent< T >>>  forced)
inlineprotected

◆ set_forced_publish_events()

void set_forced_publish_events ( std::unique_ptr< EventCollection< PublishEvent< T >>>  forced)
inlineprotected

◆ set_forced_unrestricted_update_events()

void set_forced_unrestricted_update_events ( std::unique_ptr< EventCollection< UnrestrictedUpdateEvent< T >>>  forced)
inlineprotected

◆ SetDefaultContext()

void SetDefaultContext ( Context< T > *  context) const
inline

◆ SetDefaultParameters()

virtual void SetDefaultParameters ( const Context< T > &  context,
Parameters< T > *  parameters 
) const
pure virtual

Assigns default values to all parameters.

Overrides must not change the number of parameters.

Implemented in LeafSystem< T >, LeafSystem< double >, LeafSystem< AutoDiffXd >, Diagram< T >, and Diagram< double >.

◆ SetDefaultState()

virtual void SetDefaultState ( const Context< T > &  context,
State< T > *  state 
) const
pure virtual

◆ SetRandomContext()

void SetRandomContext ( Context< T > *  context,
RandomGenerator generator 
) const
inline

◆ SetRandomParameters()

virtual void SetRandomParameters ( const Context< T > &  context,
Parameters< T > *  parameters,
RandomGenerator generator 
) const
inlinevirtual

Assigns random values to all parameters.

This default implementation calls SetDefaultParameters; override this method to provide random parameters using the stdc++ random library, e.g.:

std::uniform_real_distribution<T> uniform();
parameters->get_mutable_numeric_parameter(0)
->SetAtIndex(0, uniform(*generator));

Overrides must not change the number of state variables.

See also
Stochastic Systems

Reimplemented in Diagram< T >, and Diagram< double >.

◆ SetRandomState()

virtual void SetRandomState ( const Context< T > &  context,
State< T > *  state,
RandomGenerator generator 
) const
inlinevirtual

Assigns random values to all elements of the state.

This default implementation calls SetDefaultState; override this method to provide random initial conditions using the stdc++ random library, e.g.:

std::normal_distribution<T> gaussian();
state->get_mutable_continuous_state()->get_mutable_vector()
->SetAtIndex(0, gaussian(*generator));

Overrides must not change the number of state variables.

See also
Stochastic Systems

Reimplemented in Diagram< T >, and Diagram< double >.

◆ ToAutoDiffXd() [1/2]

std::unique_ptr<System<AutoDiffXd> > ToAutoDiffXd ( ) const
inline

Creates a deep copy of this System, transmogrified to use the autodiff scalar type, with a dynamic-sized vector of partial derivatives.

The result is never nullptr.

Exceptions
std::exceptionif this System does not support autodiff

See System Scalar Conversion for detailed background and examples related to scalar-type conversion support.

◆ ToAutoDiffXd() [2/2]

static std::unique_ptr<S<AutoDiffXd> > ToAutoDiffXd ( const S< T > &  from)
inlinestatic

Creates a deep copy of from, transmogrified to use the autodiff scalar type, with a dynamic-sized vector of partial derivatives.

The result is never nullptr.

Exceptions
std::exceptionif from does not support autodiff

Usage:

MySystem<double> plant;
std::unique_ptr<MySystem<AutoDiffXd>> ad_plant =
Template Parameters
SThe specific System type to accept and return.

See System Scalar Conversion for detailed background and examples related to scalar-type conversion support.

◆ ToAutoDiffXdMaybe()

std::unique_ptr<System<AutoDiffXd> > ToAutoDiffXdMaybe ( ) const
inline

Creates a deep copy of this system exactly like ToAutoDiffXd(), but returns nullptr if this System does not support autodiff, instead of throwing an exception.

◆ ToSymbolic() [1/2]

std::unique_ptr<System<symbolic::Expression> > ToSymbolic ( ) const
inline

Creates a deep copy of this System, transmogrified to use the symbolic scalar type.

The result is never nullptr.

Exceptions
std::exceptionif this System does not support symbolic

See System Scalar Conversion for detailed background and examples related to scalar-type conversion support.

◆ ToSymbolic() [2/2]

static std::unique_ptr<S<symbolic::Expression> > ToSymbolic ( const S< T > &  from)
inlinestatic

Creates a deep copy of from, transmogrified to use the symbolic scalar type.

The result is never nullptr.

Exceptions
std::exceptionif this System does not support symbolic

Usage:

MySystem<double> plant;
std::unique_ptr<MySystem<symbolic::Expression>> sym_plant =
Template Parameters
SThe specific System pointer type to return.

See System Scalar Conversion for detailed background and examples related to scalar-type conversion support.

◆ ToSymbolicMaybe()

std::unique_ptr<System<symbolic::Expression> > ToSymbolicMaybe ( ) const
inline

Creates a deep copy of this system exactly like ToSymbolic(), but returns nullptr if this System does not support symbolic, instead of throwing an exception.

Friends And Related Function Documentation

◆ SystemImpl

friend class SystemImpl
friend

The documentation for this class was generated from the following files: