specific dirdyn function for the integrators More...
#include <math.h>#include "integrator.h"#include "mbs_define.h"#include "mbs_dirdyn_struct.h"#include "mbs_matrix.h"#include "mbs_1D_array.h"#include "mbs_rk4.h"#include "mbs_dopri5.h"#include "mbs_euler_explicit.h"#include "mbs_euler_implicit.h"#include "mbs_rosenbrock.h"#include "mbs_custom.h"#include "mbs_bader.h"#include "mbs_w_methods.h"#include "mbs_message.h"Functions | |
| int | mbs_estim_jac_acc (double x, double htry, double y[], double dydx[], int compute_dfdx, double dfdx[], double **dfdy, int n, int(*derivs)(double, double[], double[], MbsData *, MbsDirdyn *), MbsData *s, MbsDirdyn *dd) |
| Evaluate the Jacobian of a function accelerations using finite difference. More... | |
| int | mbs_freeze_jac (int freeze_index, int *next_freeze_index, double x, double h, double y[], double dydx[], int compute_dfdx, double dfdx[], double **dfdy, double dydx_freeze[], double dfdx_freeze[], double **dfdy_freeze, int n, int(*derivs)(double, double[], double[], MbsData *, MbsDirdyn *), MbsData *s, MbsDirdyn *mbs_dd) |
| Freeze the Jacobian of an integrator structure until mbs_dd->options->n_freeze. More... | |
| void | set_integrator (MbsDirdyn *mbs_dd) |
| Set the function pointer in MbsDirdyn. More... | |
| void | print_warnings_integrator (MbsData *mbs_data, MbsDirdyn *mbs_dd, int type_of_integrator) |
| Check the options set by user to warn him when he modified an unused option for the integrator. More... | |
| void | print_warnings_constant_step_integrator (MbsDirdyn *mbs_dd, char *integrator_name) |
| Print the warning message, called by print_warnings_integrator for constant step size dt integrators. More... | |
| void | print_warnings_explicit_integrator (MbsDirdyn *mbs_dd, char *integrator_name) |
| Print the warning message, called by print_warnings_integrator for explicit integrators. More... | |
specific dirdyn function for the integrators
Creation date: May 2015
Modification date: June 2018 \modified by Sebastien Timmermans
(c) Universite catholique de Louvain
| int mbs_estim_jac_acc | ( | double | x, |
| double | htry, | ||
| double | y[], | ||
| double | dydx[], | ||
| int | compute_dfdx, | ||
| double | dfdx[], | ||
| double ** | dfdy, | ||
| int | n, | ||
| int(*)(double, double[], double[], MbsData *, MbsDirdyn *) | derivs, | ||
| MbsData * | s, | ||
| MbsDirdyn * | dd | ||
| ) |
Evaluate the Jacobian of a function accelerations using finite difference.
Modify and do not restore dydx. WARNING: All vector/matrix indexes start at 0! Update the value of dfdx and dfdy.
| x | current time value, at which the function needs to be computed |
| htry | time step |
| y | state vector of variable |
| compute_dfdx | flag to ON to compute the dfdx terms, default = 1 ; if compute_dfdx = 0 => dfdx is filled with zeroes |
| dydx | derivative of y by x |
| dfdx | derivative of f by x |
| dfdy | derivative of f by y |
| n | number of variables |
| derivs | the function computing f' |
| s | the MbsData structure of the model on which dirdyn analysis is computed. |
| dd | the MbsDirdyn structure related to the integration. |
| int mbs_freeze_jac | ( | int | freeze_index, |
| int * | next_freeze_index, | ||
| double | x, | ||
| double | h, | ||
| double | y[], | ||
| double | dydx[], | ||
| int | compute_dfdx, | ||
| double | dfdx[], | ||
| double ** | dfdy, | ||
| double | dydx_freeze[], | ||
| double | dfdx_freeze[], | ||
| double ** | dfdy_freeze, | ||
| int | n, | ||
| int(*)(double, double[], double[], MbsData *, MbsDirdyn *) | derivs, | ||
| MbsData * | s, | ||
| MbsDirdyn * | mbs_dd | ||
| ) |
Freeze the Jacobian of an integrator structure until mbs_dd->options->n_freeze.
WARNING: All vector/matrix indexes start at 0!
| freeze_index | current count of refreezing steps (if 0, recompute Jac) |
| next_freeze_index | [out] pointor to next refreezing step |
| x | current time value, at which the function needs to be computed |
| h | time step |
| y | state vector of variable |
| dydx | derivative of y by x |
| compute_dfdx | flag to ON to compute the dfdx terms, default = 1 ; if compute_dfdx = 0 => dfdx is filled with zeroes |
| dfdx | derivative of f by x |
| dfdy | derivative of f by y |
| dydx_freeze | pointor to freezed derivative of y by x |
| dfdx_freeze | pointor to freezed derivative of f by x |
| dfdy_freeze | pointor to freezed derivative of f by y |
| n | number of variables |
| derivs | the function computing f' |
| s | the MbsData structure of the model on which dirdyn analysis is computed. |
| mbs_dd | the MbsDirdyn structure related to the integration. |
| void print_warnings_constant_step_integrator | ( | MbsDirdyn * | mbs_dd, |
| char * | integrator_name | ||
| ) |
Print the warning message, called by print_warnings_integrator for constant step size dt integrators.
| mbs_dd | MbsDirdyn structure |
| integrator_name | string to print in the warning message, corresponding to the integrator (RK4, EulerEx, ...) |
| void print_warnings_explicit_integrator | ( | MbsDirdyn * | mbs_dd, |
| char * | integrator_name | ||
| ) |
Print the warning message, called by print_warnings_integrator for explicit integrators.
| mbs_dd | MbsDirdyn structure |
| integrator_name | string to print in the warning message, corresponding to the integrator (RK4, EulerEx, ...) |
1.8.17