MBDyn-1.7.3
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups
InverseDynamicsStepSolver Class Reference

#include <stepsol.h>

Inheritance diagram for InverseDynamicsStepSolver:
Collaboration diagram for InverseDynamicsStepSolver:

Public Member Functions

 InverseDynamicsStepSolver (const integer MaxIt, const doublereal dT, const doublereal dSolutionTol, const integer stp, const integer sts, const bool bmod_res_test)
 
 ~InverseDynamicsStepSolver (void)
 
virtual void EvalProd (doublereal Tau, const VectorHandler &f0, const VectorHandler &w, VectorHandler &z) const
 
virtual doublereal TestScale (const NonlinearSolverTest *pTest, doublereal &dCoef) const
 
virtual doublereal Advance (Solver *pS, const doublereal TStep, const doublereal dAlph, const StepChange StType, std::deque< MyVectorHandler * > &qX, std::deque< MyVectorHandler * > &qXPrime, MyVectorHandler *const pX, MyVectorHandler *const pXPrime, integer &EffIter, doublereal &Err, doublereal &SolErr)
 
virtual doublereal Advance (InverseSolver *pS, const doublereal TStep, const StepChange StType, MyVectorHandler *const pX, MyVectorHandler *const pXPrime, MyVectorHandler *const pXPrimePrime, MyVectorHandler *const pLambda, integer &EffIter, doublereal &Err, doublereal &SolErr)
 
void Residual (VectorHandler *pRes) const
 
void Jacobian (MatrixHandler *pJac) const
 
void Update (const VectorHandler *pSol) const
 
void SetOrder (InverseDynamics::Order iOrder)
 
InverseDynamics::Order GetOrder (void) const
 
bool bJacobian (void) const
 
- Public Member Functions inherited from StepIntegrator
 StepIntegrator (const integer MaxIt, const doublereal dT, const doublereal dSolutionTol, const integer stp, const integer sts)
 
virtual ~StepIntegrator (void)
 
void SetDataManager (DataManager *pDatMan)
 
virtual integer GetIntegratorNumPreviousStates (void) const
 
virtual integer GetIntegratorNumUnknownStates (void) const
 
virtual integer GetIntegratorMaxIters (void) const
 
virtual doublereal GetIntegratorDTol (void) const
 
virtual doublereal GetIntegratorDSolTol (void) const
 
virtual void OutputTypes (const bool fpred)
 
virtual void SetDriveHandler (const DriveHandler *pDH)
 
- Public Member Functions inherited from NonlinearProblem
virtual ~NonlinearProblem (void)
 

Protected Attributes

VectorHandlerpXCurr
 
VectorHandlerpXPrimeCurr
 
VectorHandlerpXPrimePrimeCurr
 
VectorHandlerpLambdaCurr
 
bool bModResTest
 
- Protected Attributes inherited from StepIntegrator
DataManagerpDM
 
const DataManager::DofVecTypepDofs
 
bool outputPred
 
integer MaxIters
 
doublereal dTol
 
doublereal dSolTol
 
integer steps
 
integer unkstates
 

Private Attributes

MyVectorHandler XTau
 
MyVectorHandler SavedState
 
MyVectorHandler SavedDerState
 
bool bEvalProdCalledFirstTime
 
InverseDynamics::Order iOrder
 
bool m_bJacobian
 

Additional Inherited Members

- Public Types inherited from StepIntegrator
enum  { DIFFERENTIAL = 0, ALGEBRAIC = 1 }
 
enum  StepChange { NEWSTEP, REPEATSTEP }
 
- Protected Member Functions inherited from StepIntegrator
template<class T >
void UpdateLoop (const T *const t, void(T::*pUpd)(const int DCount, const DofOrder::Order Order, const VectorHandler *const pSol) const, const VectorHandler *const pSol=0) const
 

Detailed Description

Definition at line 666 of file stepsol.h.

Constructor & Destructor Documentation

InverseDynamicsStepSolver::InverseDynamicsStepSolver ( const integer  MaxIt,
const doublereal  dT,
const doublereal  dSolutionTol,
const integer  stp,
const integer  sts,
const bool  bmod_res_test 
)

Definition at line 1378 of file stepsol.cc.

1384 : StepIntegrator(MaxIt, dT, dSolutionTol, stp, sts),
1385 XTau(0),
1386 SavedState(0),
1387 SavedDerState(0),
1390 m_bJacobian(true),
1391 pXCurr(0),
1392 pXPrimeCurr(0),
1393 pXPrimePrimeCurr(0),
1394 pLambdaCurr(0),
1395 bModResTest(bmod_res_test)
1396 {
1397  return;
1398 }
StepIntegrator(const integer MaxIt, const doublereal dT, const doublereal dSolutionTol, const integer stp, const integer sts)
Definition: stepsol.cc:51
VectorHandler * pXPrimeCurr
Definition: stepsol.h:684
MyVectorHandler XTau
Definition: stepsol.h:672
InverseDynamics::Order iOrder
Definition: stepsol.h:679
VectorHandler * pXPrimePrimeCurr
Definition: stepsol.h:685
MyVectorHandler SavedState
Definition: stepsol.h:673
MyVectorHandler SavedDerState
Definition: stepsol.h:676
VectorHandler * pXCurr
Definition: stepsol.h:683
VectorHandler * pLambdaCurr
Definition: stepsol.h:686
InverseDynamicsStepSolver::~InverseDynamicsStepSolver ( void  )

Definition at line 1400 of file stepsol.cc.

References NO_OP.

1401 {
1402  NO_OP;
1403 }
#define NO_OP
Definition: myassert.h:74

Member Function Documentation

virtual doublereal InverseDynamicsStepSolver::Advance ( Solver pS,
const doublereal  TStep,
const doublereal  dAlph,
const StepChange  StType,
std::deque< MyVectorHandler * > &  qX,
std::deque< MyVectorHandler * > &  qXPrime,
MyVectorHandler *const  pX,
MyVectorHandler *const  pXPrime,
integer EffIter,
doublereal Err,
doublereal SolErr 
)
inlinevirtual

Implements StepIntegrator.

Definition at line 708 of file stepsol.h.

References MBDYN_EXCEPT_ARGS.

719  {
720  silent_cerr("InverseDynamicsStepSolver::Advance()");
722  };
#define MBDYN_EXCEPT_ARGS
Definition: except.h:63
doublereal InverseDynamicsStepSolver::Advance ( InverseSolver pS,
const doublereal  TStep,
const StepChange  StType,
MyVectorHandler *const  pX,
MyVectorHandler *const  pXPrime,
MyVectorHandler *const  pXPrimePrime,
MyVectorHandler *const  pLambda,
integer EffIter,
doublereal Err,
doublereal SolErr 
)
virtual

Definition at line 1528 of file stepsol.cc.

References InverseDynamics::ACCELERATION, ASSERT, bJacobian(), StepIntegrator::dSolTol, StepIntegrator::dTol, InverseDynamics::FULLY_ACTUATED_COLLOCATED, DataManager::GetSolver(), DataManager::IDAfterConvergence(), InverseDynamics::INVERSE_DYNAMICS, Jacobian(), DataManager::LinkToSolution(), SolutionManager::MatrReset(), StepIntegrator::MaxIters, SolverDiagnostics::outputJac(), SolverDiagnostics::outputRes(), SolverDiagnostics::outputSol(), StepIntegrator::pDM, Solver::pGetNonlinearSolver(), Solver::pGetSolutionManager(), pLambdaCurr, SolutionManager::pMatHdl(), InverseDynamics::POSITION, SolutionManager::pResHdl(), Solver::PrintResidual(), Solver::PrintSolution(), SolutionManager::pSolHdl(), pXCurr, pXPrimeCurr, pXPrimePrimeCurr, VectorHandler::Reset(), Residual(), SetOrder(), SolutionManager::Solve(), NonlinearSolver::Solve(), SolutionManager::SolveT(), Update(), and InverseDynamics::VELOCITY.

1538 {
1539  ASSERT(pDM != NULL);
1540  pXCurr = pX;
1541  pXPrimeCurr = pXPrime;
1542  pXPrimePrimeCurr = pXPrimePrime;
1543  pLambdaCurr = pLambda;
1544 
1546 
1547  Err = 0.;
1548  NonlinearSolver *pNLSolver = pS->pGetNonlinearSolver();
1549 
1550  /* Position */
1552 
1553  /* Setting iOrder = 0, the residual is called only
1554  * for constraints, on positions */
1555  pNLSolver->Solve(this, pS, MaxIters, dTol,
1556  EffIter, Err, dSolTol, SolErr);
1557 
1558  SolutionManager *pSM = pS->pGetSolutionManager();
1559  MatrixHandler *pMat = pSM->pMatHdl();
1560  VectorHandler *pRes = pSM->pResHdl();
1561  VectorHandler *pSol = pSM->pSolHdl();
1562 
1563  /* use position */
1564  // Update(pSol);
1565 
1566  /* Velocity */
1568 
1569  pRes->Reset();
1570  pSol->Reset();
1571  /* there's no need to check changes in
1572  * equation structure... it is already
1573  * performed by NonlinearSolver->Solve() */
1574  Residual(pRes);
1575 
1576  if (pS->outputRes()) {
1577  silent_cout("Residual(velocity):" << std::endl);
1578  pS->PrintResidual(*pRes, 0);
1579  }
1580 
1581  if (bJacobian()) {
1582  pSM->MatrReset();
1583  Jacobian(pMat);
1584 
1585  if (pS->outputJac()) {
1586  silent_cout("Jacobian(velocity):" << std::endl << *pMat);
1587  }
1588 
1589  EffIter++;
1590  }
1591 
1592  pSM->Solve();
1593 
1594  if (pS->outputSol()) {
1595  silent_cout("Solution(velocity):" << std::endl);
1596  pS->PrintSolution(*pSol, 0);
1597  }
1598 
1599  /* use velocity */
1600  Update(pSol);
1601 
1602  // TODO: if UNDERDETERMINED_UNDERACTUATED_COLLOCATED,
1603  // InverseDynamics::ACCELERATION and InverseDynamics::INVERSE_DYNAMICS
1604  // are solved together
1605 
1606  /* Acceleration */
1608 
1609  pRes->Reset();
1610  pSol->Reset();
1611  Residual(pRes);
1612 
1613  if (pS->outputRes()) {
1614  silent_cout("Residual(acceleration):" << std::endl);
1615  pS->PrintResidual(*pRes, 0);
1616  }
1617 
1618  if (bJacobian()) {
1619  pSM->MatrReset();
1620  Jacobian(pMat);
1621 
1622  if (pS->outputJac()) {
1623  silent_cout("Jacobian(acceleration):" << std::endl << *pMat);
1624  }
1625 
1626  EffIter++;
1627  }
1628 
1629  pSM->Solve();
1630 
1631  if (pS->outputSol()) {
1632  silent_cout("Solution(acceleration):" << std::endl);
1633  pS->PrintSolution(*pSol, 0);
1634  }
1635 
1636  Update(pSol);
1637 
1638  /* Forces */
1640 
1641  pRes->Reset();
1642  pSol->Reset();
1643  Residual(pRes);
1644 
1645  if (pS->outputRes()) {
1646  silent_cout("Residual(inverseDynamics):" << std::endl);
1647  pS->PrintResidual(*pRes, 0);
1648  }
1649 
1650  if (bJacobian()) {
1651  pSM->MatrReset();
1652  Jacobian(pMat);
1653 
1654  if (pS->outputJac()) {
1655  silent_cout("Jacobian(inverseDynamics):" << std::endl << *pMat);
1656  }
1657 
1658  EffIter++;
1659  }
1660 
1661  switch (dynamic_cast<const InverseSolver *>(pDM->GetSolver())->GetProblemType()) {
1663  pSM->SolveT();
1664  break;
1665 
1666  default:
1667  pSM->Solve();
1668  }
1669 
1670  if (pS->outputSol()) {
1671  silent_cout("Solution(inverseDynamics):" << std::endl);
1672  pS->PrintSolution(*pSol, 0);
1673  }
1674 
1675  Update(pSol);
1676 
1678 
1679  return Err;
1680 }
virtual void Reset(void)=0
virtual VectorHandler * pResHdl(void) const =0
virtual SolutionManager * pGetSolutionManager(void) const
Definition: solver.h:398
VectorHandler * pXPrimeCurr
Definition: stepsol.h:684
void Jacobian(MatrixHandler *pJac) const
Definition: stepsol.cc:1698
virtual void PrintResidual(const VectorHandler &Res, integer iIterCnt) const
Definition: solver.cc:5194
doublereal dSolTol
Definition: stepsol.h:99
void SetOrder(InverseDynamics::Order iOrder)
Definition: stepsol.cc:1510
bool outputRes(void) const
bool outputJac(void) const
bool bJacobian(void) const
Definition: stepsol.cc:1522
virtual void IDAfterConvergence(void) const
Definition: invdataman.cc:834
bool outputSol(void) const
void LinkToSolution(VectorHandler &XCurr, VectorHandler &XPrimeCurr)
Definition: dataman2.cc:172
virtual void SolveT(void)
Definition: solman.cc:78
virtual MatrixHandler * pMatHdl(void) const =0
DataManager * pDM
Definition: stepsol.h:93
VectorHandler * pXPrimePrimeCurr
Definition: stepsol.h:685
virtual NonlinearSolver * pGetNonlinearSolver(void) const
Definition: solver.h:404
integer MaxIters
Definition: stepsol.h:98
virtual void MatrReset(void)=0
virtual void Solve(void)=0
#define ASSERT(expression)
Definition: colamd.c:977
void Residual(VectorHandler *pRes) const
Definition: stepsol.cc:1683
VectorHandler * pXCurr
Definition: stepsol.h:683
void Update(const VectorHandler *pSol) const
Definition: stepsol.cc:1705
virtual void PrintSolution(const VectorHandler &Sol, integer iIterCnt) const
Definition: solver.cc:5200
const Solver * GetSolver(void) const
Definition: dataman.h:343
virtual VectorHandler * pSolHdl(void) const =0
virtual void Solve(const NonlinearProblem *pNLP, Solver *pS, const integer iMaxIter, const doublereal &Tol, integer &iIterCnt, doublereal &dErr, const doublereal &SolTol, doublereal &dSolErr)=0
VectorHandler * pLambdaCurr
Definition: stepsol.h:686
doublereal dTol
Definition: stepsol.h:99

Here is the call graph for this function:

bool InverseDynamicsStepSolver::bJacobian ( void  ) const

Definition at line 1522 of file stepsol.cc.

References m_bJacobian.

Referenced by Advance().

1523 {
1524  return m_bJacobian;
1525 }
void InverseDynamicsStepSolver::EvalProd ( doublereal  Tau,
const VectorHandler f0,
const VectorHandler w,
VectorHandler z 
) const
virtual

Implements NonlinearProblem.

Definition at line 1406 of file stepsol.cc.

References ASSERT, bEvalProdCalledFirstTime, copysign(), f0, grad::fabs(), VectorHandler::iGetSize(), VectorHandler::InnerProd(), VectorHandler::Norm(), StepIntegrator::pDM, pXCurr, pXPrimeCurr, VectorHandler::Reset(), MyVectorHandler::Reset(), Residual(), MyVectorHandler::Resize(), SavedDerState, SavedState, VectorHandler::ScalarMul(), MyVectorHandler::ScalarMul(), DataManager::Update(), Update(), and XTau.

1408 {
1409  /* matrix-free product
1410  *
1411  * J(XCurr) * w = -||w|| * (Res(XCurr + sigma * Tau * w/||w||) - f0) / (sigma * Tau)
1412  *
1413  */
1415  XTau.Resize(w.iGetSize());
1416  SavedState.Resize(w.iGetSize());
1418  bEvalProdCalledFirstTime = false;
1419  }
1420 
1421  SavedState = *pXCurr;
1423 
1424  /* if w = 0; J * w = 0 */
1425  ASSERT(pDM != NULL);
1426 
1427  doublereal nw = w.Norm();
1428  if (nw < std::numeric_limits<doublereal>::epsilon()) {
1429  z.Reset();
1430  return;
1431  }
1432  doublereal sigma = pXCurr->InnerProd(w);
1433  sigma /= nw;
1434  if (fabs(sigma) > std::numeric_limits<doublereal>::epsilon()) {
1435  doublereal xx = (fabs( sigma) <= 1.) ? 1. : fabs(sigma);
1436  Tau = copysign(Tau*xx, sigma);
1437  }
1438  Tau /= nw;
1439 #ifdef DEBUG_ITERATIVE
1440  std::cout << "Tau " << Tau << std::endl;
1441 #endif /* DEBUG_ITERATIVE */
1442 
1443  XTau.Reset();
1444  z.Reset();
1445  XTau.ScalarMul(w, Tau);
1446  Update(&XTau);
1447 #ifdef USE_EXTERNAL
1448  External::SendFreeze();
1449 #endif /* USE_EXTERNAL */
1450  /* deal with throwing elements: do not honor their requests while perfoming matrix free update */
1451  try {
1452  Residual(&z);
1453  }
1455  }
1456  XTau.ScalarMul(XTau, -1.);
1457 
1458  /* riporta tutto nelle condizioni inziali */
1459 #if 0
1460  Update(&XTau);
1461 #endif
1462  *pXCurr = SavedState;
1464  pDM->Update();
1465  z -= f0;
1466  z.ScalarMul(z, -1./Tau);
1467 }
virtual void Reset(void)=0
VectorHandler * pXPrimeCurr
Definition: stepsol.h:684
MyVectorHandler XTau
Definition: stepsol.h:672
static double * f0
virtual doublereal InnerProd(const VectorHandler &VH) const
Definition: vh.cc:269
virtual integer iGetSize(void) const =0
virtual doublereal Norm(void) const
Definition: vh.cc:262
GradientExpression< UnaryExpr< FuncFabs, Expr > > fabs(const GradientExpression< Expr > &u)
Definition: gradient.h:2973
virtual VectorHandler & ScalarMul(const VectorHandler &VH, const doublereal &d)
Definition: vh.cc:519
virtual void Resize(integer iSize)
Definition: vh.cc:347
DataManager * pDM
Definition: stepsol.h:93
doublereal copysign(doublereal x, doublereal y)
Definition: gradient.h:97
virtual VectorHandler & ScalarMul(const VectorHandler &VH, const doublereal &d)
Definition: vh.cc:145
#define ASSERT(expression)
Definition: colamd.c:977
virtual void Reset(void)
Definition: vh.cc:459
MyVectorHandler SavedState
Definition: stepsol.h:673
void Residual(VectorHandler *pRes) const
Definition: stepsol.cc:1683
MyVectorHandler SavedDerState
Definition: stepsol.h:676
VectorHandler * pXCurr
Definition: stepsol.h:683
void Update(const VectorHandler *pSol) const
Definition: stepsol.cc:1705
virtual void Update(void) const
Definition: dataman2.cc:2511
double doublereal
Definition: colamd.c:52

Here is the call graph for this function:

InverseDynamics::Order InverseDynamicsStepSolver::GetOrder ( void  ) const

Definition at line 1516 of file stepsol.cc.

References iOrder.

Referenced by DataManager::AssConstrJac().

1517 {
1518  return iOrder;
1519 }
InverseDynamics::Order iOrder
Definition: stepsol.h:679
void InverseDynamicsStepSolver::Jacobian ( MatrixHandler pJac) const
virtual

Implements NonlinearProblem.

Definition at line 1698 of file stepsol.cc.

References DataManager::AssConstrJac(), ASSERT, and StepIntegrator::pDM.

Referenced by Advance().

1699 {
1700  ASSERT(pDM != NULL);
1701  pDM->AssConstrJac(*pJac);
1702 }
DataManager * pDM
Definition: stepsol.h:93
#define ASSERT(expression)
Definition: colamd.c:977
virtual void AssConstrJac(MatrixHandler &JacHdl)
Definition: invdataman.cc:94

Here is the call graph for this function:

void InverseDynamicsStepSolver::Residual ( VectorHandler pRes) const
virtual

Implements NonlinearProblem.

Definition at line 1683 of file stepsol.cc.

References DataManager::AssConstrRes(), ASSERT, DataManager::AssRes(), InverseDynamics::INVERSE_DYNAMICS, iOrder, and StepIntegrator::pDM.

Referenced by Advance(), and EvalProd().

1684 {
1685  ASSERT(pDM != NULL);
1686  switch (iOrder) {
1688  pDM->AssRes(*pRes);
1689  break;
1690 
1691  default:
1692  pDM->AssConstrRes(*pRes, iOrder);
1693  break;
1694  }
1695 }
virtual void AssRes(VectorHandler &ResHdl, doublereal dCoef)
Definition: elman.cc:498
InverseDynamics::Order iOrder
Definition: stepsol.h:679
DataManager * pDM
Definition: stepsol.h:93
#define ASSERT(expression)
Definition: colamd.c:977
virtual void AssConstrRes(VectorHandler &ResHdl, InverseDynamics::Order iOrder)
Definition: invdataman.cc:303

Here is the call graph for this function:

void InverseDynamicsStepSolver::SetOrder ( InverseDynamics::Order  iOrder)

Definition at line 1510 of file stepsol.cc.

References iOrder.

Referenced by Advance().

1511 {
1512  this->iOrder = iOrder;
1513 }
InverseDynamics::Order iOrder
Definition: stepsol.h:679
doublereal InverseDynamicsStepSolver::TestScale ( const NonlinearSolverTest pTest,
doublereal dCoef 
) const
virtual

Implements NonlinearProblem.

Definition at line 1471 of file stepsol.cc.

References ASSERT, bModResTest, DofOrder::DIFFERENTIAL, NonlinearSolverTest::dScaleCoef(), VectorHandler::iGetSize(), StepIntegrator::pDofs, and pXPrimeCurr.

1472 {
1473  dCoef = 1.;
1474 
1475  if (bModResTest) {
1476 #ifdef USE_MPI
1477 #warning "InverseDynamicsStepSolver::TestScale() not available with Schur solution"
1478 #endif /* USE_MPI */
1479 
1480  doublereal dXPr = 0.;
1481 
1482  DataManager::DofIterator_const CurrDof = pDofs->begin();
1483 
1484  for (int iCntp1 = 1; iCntp1 <= pXPrimeCurr->iGetSize();
1485  iCntp1++, ++CurrDof)
1486  {
1487  ASSERT(CurrDof != pDofs->end());
1488 
1489  if (CurrDof->Order == DofOrder::DIFFERENTIAL) {
1490  doublereal d = pXPrimeCurr->operator()(iCntp1);
1491  doublereal d2 = d*d;
1492 
1493  doublereal ds = pTest->dScaleCoef(iCntp1);
1494  doublereal ds2 = ds*ds;
1495  d2 *= ds2;
1496 
1497  dXPr += d2;
1498  }
1499  /* else if ALGEBRAIC: non aggiunge nulla */
1500  }
1501 
1502  return 1./(1. + dXPr);
1503 
1504  } else {
1505  return 1.;
1506  }
1507 }
VectorHandler * pXPrimeCurr
Definition: stepsol.h:684
virtual integer iGetSize(void) const =0
const DataManager::DofVecType * pDofs
Definition: stepsol.h:94
DofVecType::const_iterator DofIterator_const
Definition: dataman.h:802
#define ASSERT(expression)
Definition: colamd.c:977
virtual const doublereal & dScaleCoef(const integer &iIndex) const
Definition: nonlin.cc:167
double doublereal
Definition: colamd.c:52

Here is the call graph for this function:

void InverseDynamicsStepSolver::Update ( const VectorHandler pSol) const
virtual

Implements NonlinearProblem.

Definition at line 1705 of file stepsol.cc.

References InverseDynamics::ACCELERATION, ASSERT, DEBUGCOUTFNAME, InverseDynamics::FULLY_ACTUATED_COLLOCATED, InverseDynamics::FULLY_ACTUATED_NON_COLLOCATED, DataManager::GetSolver(), InverseDynamics::INVERSE_DYNAMICS, iOrder, m_bJacobian, MBDYN_EXCEPT_ARGS, StepIntegrator::pDM, pLambdaCurr, InverseDynamics::POSITION, pXCurr, pXPrimeCurr, pXPrimePrimeCurr, InverseDynamics::UNDERDETERMINED_FULLY_ACTUATED, InverseDynamics::UNDERDETERMINED_UNDERACTUATED_COLLOCATED, DataManager::Update(), and InverseDynamics::VELOCITY.

Referenced by Advance(), and EvalProd().

1706 {
1707  DEBUGCOUTFNAME("InverseDynamicsStepSolver::Update");
1708  ASSERT(pDM != NULL);
1709 
1710  switch (iOrder) {
1712  *pXCurr += *pSol;
1713  break;
1714 
1716  *pXPrimeCurr = *pSol;
1717  break;
1718 
1720  *pXPrimePrimeCurr = *pSol;
1721  break;
1722 
1724  *pLambdaCurr = *pSol;
1725  break;
1726 
1727  default:
1728  ASSERT(0);
1730  }
1731 
1732  pDM->Update(iOrder);
1733 
1734  // prepare m_bJacobian for next phase
1735  switch (dynamic_cast<const InverseSolver *>(pDM->GetSolver())->GetProblemType()) {
1737  switch (iOrder) {
1740  m_bJacobian = true;
1741  break;
1742 
1745  m_bJacobian = false;
1746  break;
1747 
1748  default:
1749  break;
1750  }
1751  break;
1752 
1754  switch (iOrder) {
1758  m_bJacobian = true;
1759  break;
1760 
1762  m_bJacobian = false;
1763  break;
1764 
1765  default:
1766  break;
1767  }
1768  break;
1769 
1771  // TODO
1773 
1775  m_bJacobian = true;
1776  break;
1777 
1778  default:
1779  break;
1780  }
1781 }
#define MBDYN_EXCEPT_ARGS
Definition: except.h:63
VectorHandler * pXPrimeCurr
Definition: stepsol.h:684
#define DEBUGCOUTFNAME(fname)
Definition: myassert.h:256
InverseDynamics::Order iOrder
Definition: stepsol.h:679
DataManager * pDM
Definition: stepsol.h:93
VectorHandler * pXPrimePrimeCurr
Definition: stepsol.h:685
#define ASSERT(expression)
Definition: colamd.c:977
VectorHandler * pXCurr
Definition: stepsol.h:683
const Solver * GetSolver(void) const
Definition: dataman.h:343
virtual void Update(void) const
Definition: dataman2.cc:2511
VectorHandler * pLambdaCurr
Definition: stepsol.h:686

Here is the call graph for this function:

Member Data Documentation

bool InverseDynamicsStepSolver::bEvalProdCalledFirstTime
mutableprivate

Definition at line 677 of file stepsol.h.

Referenced by EvalProd().

bool InverseDynamicsStepSolver::bModResTest
protected

Definition at line 687 of file stepsol.h.

Referenced by TestScale().

InverseDynamics::Order InverseDynamicsStepSolver::iOrder
private

Definition at line 679 of file stepsol.h.

Referenced by GetOrder(), Residual(), SetOrder(), and Update().

bool InverseDynamicsStepSolver::m_bJacobian
mutableprivate

Definition at line 680 of file stepsol.h.

Referenced by bJacobian(), and Update().

VectorHandler* InverseDynamicsStepSolver::pLambdaCurr
protected

Definition at line 686 of file stepsol.h.

Referenced by Advance(), and Update().

VectorHandler* InverseDynamicsStepSolver::pXCurr
protected

Definition at line 683 of file stepsol.h.

Referenced by Advance(), EvalProd(), and Update().

VectorHandler* InverseDynamicsStepSolver::pXPrimeCurr
protected

Definition at line 684 of file stepsol.h.

Referenced by Advance(), EvalProd(), TestScale(), and Update().

VectorHandler* InverseDynamicsStepSolver::pXPrimePrimeCurr
protected

Definition at line 685 of file stepsol.h.

Referenced by Advance(), and Update().

MyVectorHandler InverseDynamicsStepSolver::SavedDerState
mutableprivate

Definition at line 676 of file stepsol.h.

Referenced by EvalProd().

MyVectorHandler InverseDynamicsStepSolver::SavedState
mutableprivate

Definition at line 673 of file stepsol.h.

Referenced by EvalProd().

MyVectorHandler InverseDynamicsStepSolver::XTau
mutableprivate

Definition at line 672 of file stepsol.h.

Referenced by EvalProd().


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