• Skip to content
  • Skip to link menu
KDE 4.2 API Reference
  • KDE API Reference
  • kdeedu
  • Sitemap
  • Contact Us
 

step/stepcore

world.cc

Go to the documentation of this file.
00001 /* This file is part of StepCore library.
00002    Copyright (C) 2007 Vladimir Kuznetsov <ks.vladimir@gmail.com>
00003 
00004    StepCore library is free software; you can redistribute it and/or modify
00005    it under the terms of the GNU General Public License as published by
00006    the Free Software Foundation; either version 2 of the License, or
00007    (at your option) any later version.
00008 
00009    StepCore library is distributed in the hope that it will be useful,
00010    but WITHOUT ANY WARRANTY; without even the implied warranty of
00011    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00012    GNU General Public License for more details.
00013 
00014    You should have received a copy of the GNU General Public License
00015    along with StepCore; if not, write to the Free Software
00016    Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA  02110-1301  USA
00017 */
00018 
00019 #include "world.h"
00020 #include "solver.h"
00021 #include "collisionsolver.h"
00022 #include "constraintsolver.h"
00023 
00024 #include <algorithm>
00025 #include <cmath>
00026 
00027 namespace StepCore
00028 {
00029 
00030 STEPCORE_META_OBJECT(Item, "Item", MetaObject::ABSTRACT, STEPCORE_SUPER_CLASS(Object),
00031         STEPCORE_PROPERTY_RW(StepCore::Color, color, STEPCORE_UNITS_NULL, "Item color", color, setColor))
00032 STEPCORE_META_OBJECT(Body, "Body", MetaObject::ABSTRACT,,)
00033 STEPCORE_META_OBJECT(Force, "Force", MetaObject::ABSTRACT,,)
00034 STEPCORE_META_OBJECT(Joint, "Joint", MetaObject::ABSTRACT,,)
00035 STEPCORE_META_OBJECT(Tool, "Tool", MetaObject::ABSTRACT,,)
00036 
00037 STEPCORE_META_OBJECT(ObjectErrors, "ObjectErrors", MetaObject::ABSTRACT, STEPCORE_SUPER_CLASS(Object),)
00038 
00039 STEPCORE_META_OBJECT(ItemGroup, "ItemGroup", 0, STEPCORE_SUPER_CLASS(Item),)
00040 
00041 STEPCORE_META_OBJECT(World, "World", 0, STEPCORE_SUPER_CLASS(ItemGroup),
00042         STEPCORE_PROPERTY_RW_D(double, time, "s", "Current time", time, setTime)
00043         STEPCORE_PROPERTY_RW  (double, timeScale, STEPCORE_UNITS_1, "Simulation speed scale", timeScale, setTimeScale)
00044         STEPCORE_PROPERTY_RW  (bool, errorsCalculation, STEPCORE_UNITS_NULL,
00045                         "Enable global errors calculation", errorsCalculation, setErrorsCalculation))
00046 
00047 Item& Item::operator=(const Item& item)
00048 {
00049     Object::operator=(item);
00050 
00051     _world = item._world;
00052     _group = item._group;
00053 
00054     if(item._objectErrors) {
00055         _objectErrors = static_cast<ObjectErrors*>(
00056             item._objectErrors->metaObject()->cloneObject(*item._objectErrors) );
00057         _objectErrors->setOwner(this);
00058     } else {
00059         _objectErrors = NULL;
00060     }
00061 
00062     _color = item._color;
00063 
00064     return *this;
00065 }
00066 
00067 ObjectErrors* Item::objectErrors()
00068 {
00069     if(!_objectErrors) _objectErrors = createObjectErrors();
00070     return _objectErrors;
00071 }
00072 
00073 ItemGroup::ItemGroup(const ItemGroup& group)
00074     : Item()
00075 {
00076     *this = group;
00077 }
00078 
00079 void ItemGroup::setWorld(World* world)
00080 {
00081     ItemList::const_iterator end = _items.end();
00082     for(ItemList::const_iterator it = _items.begin(); it != end; ++it)
00083         (*it)->setWorld(world);
00084     Item::setWorld(world);
00085 }
00086 
00087 void ItemGroup::worldItemRemoved(Item* item)
00088 {
00089     ItemList::const_iterator end = _items.end();
00090     for(ItemList::const_iterator it = _items.begin(); it != end; ++it)
00091         (*it)->worldItemRemoved(item);
00092 }
00093 
00094 void ItemGroup::addItem(Item* item)
00095 {
00096     _items.push_back(item);
00097 
00098     if(world()) world()->worldItemAdded(item);
00099 
00100     item->setGroup(this);
00101     item->setWorld(this->world());
00102 }
00103 
00104 void ItemGroup::removeItem(Item* item)
00105 {
00106     item->setWorld(NULL);
00107     item->setGroup(NULL);
00108 
00109     if(world()) world()->worldItemRemoved(item);
00110 
00111     ItemList::iterator i = std::find(_items.begin(), _items.end(), item);
00112     STEPCORE_ASSERT_NOABORT(i != _items.end());
00113     _items.erase(i);
00114 }
00115 
00116 void ItemGroup::clear()
00117 {
00118     ItemList::const_iterator end = _items.end();
00119     for(ItemList::const_iterator it = _items.begin(); it != end; ++it) {
00120         (*it)->setWorld(NULL);
00121         (*it)->setGroup(NULL);
00122         if(world()) world()->worldItemRemoved(*it);
00123     }
00124 
00125     for(ItemList::const_iterator it = _items.begin(); it != end; ++it) {
00126         delete *it;
00127     }
00128 
00129     _items.clear();
00130 }
00131 
00132 ItemGroup::~ItemGroup()
00133 {
00134     clear();
00135 }
00136 
00137 ItemGroup& ItemGroup::operator=(const ItemGroup& group)
00138 {
00139 
00140     /*
00141     item->setGroup(this);
00142     item->setWorld(this->world());
00143     _items.push_back(item);
00144 
00145     if(world()) {
00146         //world()->worldItemAdded(item);
00147         ItemGroup* gr = dynamic_cast<ItemGroup*>(item);
00148         if(gr) gr->groupItemsAdded();
00149     }
00150     */
00151 
00152     if(this == &group) return *this;
00153 
00154     clear();
00155 
00156     _items.reserve(group._items.size());
00157 
00158     const ItemList::const_iterator gr_end = group._items.end();
00159     for(ItemList::const_iterator it = group._items.begin(); it != gr_end; ++it) {
00160         StepCore::Item* item = static_cast<Item*>( (*it)->metaObject()->cloneObject(*(*it)) );
00161         _items.push_back(item);
00162     }
00163 
00164     const ItemList::const_iterator end = _items.end();
00165     for(ItemList::const_iterator it = _items.begin(); it != end; ++it) {
00166         (*it)->setGroup(this);
00167     }
00168 
00169     Item::operator=(group);
00170 
00171     // NOTE: We don't change world() here
00172 
00173     return *this;
00174 }
00175 
00176 int ItemGroup::childItemIndex(const Item* item) const
00177 {
00178     ItemList::const_iterator o = std::find(_items.begin(), _items.end(), item);
00179     STEPCORE_ASSERT_NOABORT(o != _items.end());
00180     return std::distance(_items.begin(), o);
00181 }
00182 
00183 Item* ItemGroup::childItem(const QString& name) const
00184 {
00185     ItemList::const_iterator end = _items.end();
00186     for(ItemList::const_iterator it = _items.begin(); it != end; ++it)
00187         if((*it)->name() == name) return *it;
00188     return NULL;
00189 }
00190 
00191 Item* ItemGroup::item(const QString& name) const
00192 {
00193     if(name.isEmpty()) return NULL;
00194     ItemList::const_iterator end = _items.end();
00195     for(ItemList::const_iterator it = _items.begin(); it != end; ++it) {
00196         if((*it)->name() == name) return *it;
00197         if((*it)->metaObject()->inherits<ItemGroup>()) {
00198             Item* ret = static_cast<ItemGroup*>(*it)->item(name);
00199             if(ret) return ret;
00200         }
00201     }
00202     return NULL;
00203 }
00204 
00205 void ItemGroup::allItems(ItemList* items) const
00206 {
00207     items->reserve(_items.size());
00208     ItemList::const_iterator end = _items.end();
00209     for(ItemList::const_iterator it = _items.begin(); it != end; ++it) {
00210         items->push_back(*it);
00211         if((*it)->metaObject()->inherits<ItemGroup>())
00212             static_cast<ItemGroup*>(*it)->allItems(items);
00213     }
00214 }
00215 
00216 void ConstraintsInfo::setDimension(int newVariablesCount, int newConstraintsCount)
00217 {
00218     constraintsCount += contactsCount;
00219     contactsCount = 0;
00220 
00221     if(variablesCount != newConstraintsCount || constraintsCount != newConstraintsCount) {
00222         jacobian.resize(newConstraintsCount, newVariablesCount);
00223         jacobianDerivative.resize(newConstraintsCount, newVariablesCount);
00224     }
00225 
00226     if(variablesCount != newVariablesCount) {
00227         variablesCount = newVariablesCount;
00228         inverseMass.resize(variablesCount, variablesCount);
00229         force.resize(variablesCount);
00230     }
00231 
00232     if(constraintsCount != newConstraintsCount) {
00233         constraintsCount = newConstraintsCount;
00234         value.resize(constraintsCount);
00235         derivative.resize(constraintsCount);
00236         forceMin.resize(constraintsCount);
00237         forceMax.resize(constraintsCount);
00238     }
00239 }
00240 
00241 int ConstraintsInfo::addContact()
00242 {
00243     int n = constraintsCount + contactsCount; ++contactsCount;
00244     int s = n+1;
00245 
00246     jacobian.resize(s, variablesCount);
00247     jacobianDerivative.resize(s, variablesCount);
00248 
00249     value.resize(s); value[n] = 0;
00250     derivative.resize(s); derivative[n] = 0;
00251     forceMin.resize(s); forceMin[n] = -HUGE_VAL;
00252     forceMax.resize(s); forceMax[n] = HUGE_VAL;
00253 
00254     return n;
00255 }
00256 
00257 void ConstraintsInfo::clearContacts()
00258 {
00259     if(contactsCount == 0) return;
00260     contactsCount = 0;
00261     jacobian.resize(constraintsCount, variablesCount);
00262     jacobianDerivative.resize(constraintsCount, variablesCount);
00263     value.resize(constraintsCount);
00264     derivative.resize(constraintsCount);
00265     forceMin.resize(constraintsCount);
00266     forceMax.resize(constraintsCount);
00267 }
00268 
00269 void ConstraintsInfo::clear()
00270 {
00271     clearContacts();
00272 
00273     gmm::clear(jacobian);
00274     gmm::clear(jacobianDerivative);
00275     gmm::clear(inverseMass);
00276 
00277     std::fill(forceMin.begin(), forceMin.end(), -HUGE_VAL);
00278     std::fill(forceMax.begin(), forceMax.end(),  HUGE_VAL);
00279 
00280     collisionFlag = false;
00281 }
00282 
00283 World::World()
00284     : _time(0), _timeScale(1), _errorsCalculation(false),
00285       _solver(NULL), _collisionSolver(NULL), _constraintSolver(NULL),
00286       _variablesCount(0)
00287 {
00288     setColor(0xffffffff);
00289     setWorld(this);
00290     clear();
00291 }
00292 
00293 World::World(const World& world)
00294     : ItemGroup(), _time(0), _timeScale(1), _errorsCalculation(false),
00295       _solver(NULL), _collisionSolver(NULL), _constraintSolver(NULL),
00296       _variablesCount(0)
00297 {
00298     *this = world;
00299 }
00300 
00301 World::~World()
00302 {
00303     clear();
00304 }
00305 
00306 void World::clear()
00307 {
00308     // Avoid erasing each element individually in the cache
00309     if(_collisionSolver) _collisionSolver->resetCaches();
00310 
00311     // clear _items
00312     ItemGroup::clear();
00313 
00314     STEPCORE_ASSERT_NOABORT(_bodies.empty());
00315     STEPCORE_ASSERT_NOABORT(_forces.empty());
00316     STEPCORE_ASSERT_NOABORT(_joints.empty());
00317     //_bodies.clear();
00318     //_forces.clear();
00319 
00320     delete _solver; _solver = NULL;
00321     delete _collisionSolver; _collisionSolver = NULL;
00322     delete _constraintSolver; _constraintSolver = NULL;
00323 
00324     _variablesCount = 0;
00325     _variables.resize(0);
00326     _variances.resize(0);
00327     _tempArray.resize(0);
00328 
00329     _constraintsInfo.setDimension(0, 0);
00330 
00331     setColor(0xffffffff);
00332     deleteObjectErrors();
00333 
00334     _time = 0;
00335     _timeScale = 1;
00336     _errorsCalculation = false;
00337 
00338     _stopOnCollision = false;
00339     _stopOnIntersection = false;
00340     _evolveAbort = false;
00341 
00342 #ifdef STEPCORE_WITH_QT
00343     setName(QString());
00344 #endif
00345 }
00346 
00347 World& World::operator=(const World& world)
00348 {
00349     if(this == &world) return *this;
00350 
00351     clear();
00352     ItemGroup::operator=(world);
00353 
00354     if(world._solver) {
00355         setSolver(static_cast<Solver*>(
00356                 world._solver->metaObject()->cloneObject(*(world._solver))));
00357     } else setSolver(0);
00358 
00359     if(world._collisionSolver) {
00360         setCollisionSolver(static_cast<CollisionSolver*>(
00361                world._collisionSolver->metaObject()->cloneObject(*(world._collisionSolver))));
00362     } else setCollisionSolver(0);
00363 
00364     if(world._constraintSolver) setConstraintSolver(static_cast<ConstraintSolver*>(
00365                world._constraintSolver->metaObject()->cloneObject(*(world._constraintSolver))));
00366     else setConstraintSolver(0);
00367 
00368     _time = world._time;
00369     _timeScale = world._timeScale;
00370     _errorsCalculation = world._errorsCalculation;
00371 
00372     _stopOnCollision = world._stopOnCollision;
00373     _stopOnIntersection = world._stopOnIntersection;
00374     _evolveAbort = world._evolveAbort;
00375 
00376     // Fix links
00377     QHash<const Object*, Object*> copyMap;
00378     copyMap.insert(NULL, NULL);
00379     copyMap.insert(&world, this);
00380     if(_solver) copyMap.insert(world._solver, _solver);
00381     if(_collisionSolver) copyMap.insert(world._collisionSolver, _collisionSolver);
00382     if(_constraintSolver) copyMap.insert(world._constraintSolver, _constraintSolver);
00383     fillCopyMap(&copyMap, &world, this);
00384 
00385     applyCopyMap(&copyMap, this);
00386     if(_solver) applyCopyMap(&copyMap, _solver);
00387     if(_collisionSolver) applyCopyMap(&copyMap, _collisionSolver);
00388     if(_constraintSolver) applyCopyMap(&copyMap, _constraintSolver);
00389 
00390     const ItemList::const_iterator end = items().end();
00391     for(ItemList::const_iterator it = items().begin(); it != end; ++it) {
00392         worldItemCopied(&copyMap, *it);
00393     }
00394 
00395     setWorld(this);
00396 
00397     checkVariablesCount();
00398 
00399     return *this;
00400 }
00401 
00402 void World::fillCopyMap(QHash<const Object*, Object*>* map,
00403                         const ItemGroup* g1, ItemGroup* g2)
00404 {
00405     const ItemList::const_iterator end = g1->items().end();
00406     for(ItemList::const_iterator it1 = g1->items().begin(),
00407                                  it2 = g2->items().begin();
00408                                  it1 != end; ++it1, ++it2) {
00409         map->insert(*it1, *it2);
00410         if((*it1)->metaObject()->inherits<StepCore::ItemGroup>())
00411             fillCopyMap(map, static_cast<ItemGroup*>(*it1), static_cast<ItemGroup*>(*it2));
00412     }
00413 }
00414 
00415 void World::applyCopyMap(QHash<const Object*, Object*>* map, Object* obj)
00416 {
00417     const StepCore::MetaObject* mobj = obj->metaObject();
00418     for(int i=0; i<mobj->propertyCount(); ++i) {
00419         const StepCore::MetaProperty* pr = mobj->property(i);
00420         if(pr->userTypeId() == qMetaTypeId<Object*>()) {
00421             QVariant v = pr->readVariant(obj);
00422             v = QVariant::fromValue(map->value(v.value<Object*>(), NULL));
00423             pr->writeVariant(obj, v);
00424         }
00425     }
00426 }
00427 
00428 void World::worldItemCopied(QHash<const Object*, Object*>* map, Item* item)
00429 {
00430     applyCopyMap(map, item);
00431 
00432     if(item->metaObject()->inherits<Force>())
00433         _forces.push_back(dynamic_cast<Force*>(item));
00434     if(item->metaObject()->inherits<Joint>())
00435         _joints.push_back(dynamic_cast<Joint*>(item));
00436     if(item->metaObject()->inherits<Body>())
00437         _bodies.push_back(dynamic_cast<Body*>(item));
00438 
00439     if(item->metaObject()->inherits<ItemGroup>()) {
00440         ItemGroup* group = static_cast<ItemGroup*>(item);
00441         ItemList::const_iterator end = group->items().end();
00442         for(ItemList::const_iterator it = group->items().begin(); it != end; ++it) {
00443             worldItemCopied(map, *it);
00444         }
00445     }
00446 }
00447 
00448 void World::worldItemAdded(Item* item)
00449 {
00450     if(item->metaObject()->inherits<Force>())
00451         _forces.push_back(dynamic_cast<Force*>(item));
00452 
00453     if(item->metaObject()->inherits<Joint>())
00454         _joints.push_back(dynamic_cast<Joint*>(item));
00455 
00456     if(item->metaObject()->inherits<Body>()) {
00457         Body* body = dynamic_cast<Body*>(item);
00458         _bodies.push_back(body);
00459         if(_collisionSolver) _collisionSolver->bodyAdded(_bodies, body);
00460     }
00461 
00462     if(item->metaObject()->inherits<ItemGroup>()) {
00463         ItemGroup* group = static_cast<ItemGroup*>(item);
00464         ItemList::const_iterator end = group->items().end();
00465         for(ItemList::const_iterator it = group->items().begin(); it != end; ++it) {
00466             worldItemAdded(*it);
00467         }
00468     }
00469 
00470     checkVariablesCount();
00471 }
00472 
00473 void World::worldItemRemoved(Item* item)
00474 {
00475     if(item->metaObject()->inherits<ItemGroup>()) {
00476         ItemGroup* group = static_cast<ItemGroup*>(item);
00477         ItemList::const_iterator end = group->items().end();
00478         for(ItemList::const_iterator it = group->items().begin(); it != end; ++it) {
00479             worldItemRemoved(*it);
00480         }
00481     }
00482 
00483     const ItemList::const_iterator end = items().end();
00484     for(ItemList::const_iterator it = items().begin(); it != end; ++it) {
00485         (*it)->worldItemRemoved(item);
00486     }
00487 
00488     if(item->metaObject()->inherits<Body>()) {
00489         Body* body = dynamic_cast<Body*>(item);
00490         if(_collisionSolver) _collisionSolver->bodyRemoved(_bodies, body);
00491         BodyList::iterator b = std::find(_bodies.begin(), _bodies.end(), body);
00492         STEPCORE_ASSERT_NOABORT(b != _bodies.end());
00493         _bodies.erase(b);
00494     }
00495 
00496     if(item->metaObject()->inherits<Joint>()) {
00497         JointList::iterator j = std::find(_joints.begin(), _joints.end(),
00498                                             dynamic_cast<Joint*>(item));
00499         STEPCORE_ASSERT_NOABORT(j != _joints.end());
00500         _joints.erase(j);
00501     }
00502 
00503     if(item->metaObject()->inherits<Force>()) {
00504         ForceList::iterator f = std::find(_forces.begin(), _forces.end(),
00505                                             dynamic_cast<Force*>(item));
00506         STEPCORE_ASSERT_NOABORT(f != _forces.end());
00507         _forces.erase(f);
00508     }
00509 
00510     // XXX: on ItemGroup::clear this will be called on each object !
00511     checkVariablesCount();
00512 }
00513 
00514 /*
00515 void World::addItem(Item* item)
00516 {
00517     _items.push_back(item);
00518     item->setWorld(this);
00519     Force* force = dynamic_cast<Force*>(item);
00520     if(force) _forces.push_back(force);
00521     Body* body = dynamic_cast<Body*>(item);
00522     if(body) _bodies.push_back(body);
00523 }
00524 */
00525 
00526 /*void World::removeItem(Item* item)
00527 {
00528     const ItemList::const_iterator it_end = _items.end();
00529     for(ItemList::iterator it = _items.begin(); it != it_end; ++it)
00530         (*it)->worldItemRemoved(item);
00531 
00532     item->setWorld(NULL);
00533 
00534     ItemList::iterator i = std::find(_items.begin(), _items.end(), item);
00535     STEPCORE_ASSERT_NOABORT(i != _items.end());
00536     _items.erase(i);
00537 
00538     Force* force = dynamic_cast<Force*>(item);
00539     if(force) {
00540         ForceList::iterator f = std::find(_forces.begin(), _forces.end(), force);
00541         STEPCORE_ASSERT_NOABORT(f != _forces.end());
00542         _forces.erase(f);
00543     }
00544 
00545     Body* body = dynamic_cast<Body*>(item);
00546     if(body) {
00547         BodyList::iterator b = std::find(_bodies.begin(), _bodies.end(), body);
00548         STEPCORE_ASSERT_NOABORT(b != _bodies.end());
00549         _bodies.erase(b);
00550     }
00551 }
00552 */
00553 
00554 /*
00555 int World::itemIndex(const Item* item) const
00556 {
00557     ItemList::const_iterator o = std::find(_items.begin(), _items.end(), item);
00558     STEPCORE_ASSERT_NOABORT(o != _items.end());
00559     return std::distance(_items.begin(), o);
00560 }
00561 */
00562 
00563 /*
00564 Item* World::item(const QString& name) const
00565 {
00566     for(ItemList::const_iterator o = _items.begin(); o != _items.end(); ++o) {
00567         if((*o)->name() == name) return *o;
00568     }
00569     return NULL;
00570 }
00571 */
00572 
00573 Object* World::object(const QString& name)
00574 {
00575     if(name.isEmpty()) return NULL;
00576     if(this->name() == name) return this;
00577     else if(_solver && _solver->name() == name) return _solver;
00578     else if(_collisionSolver && _collisionSolver->name() == name) return _collisionSolver;
00579     else if(_constraintSolver && _constraintSolver->name() == name) return _constraintSolver;
00580     else return item(name);
00581 }
00582 
00583 void World::setSolver(Solver* solver)
00584 {
00585     delete _solver;
00586     _solver = solver;
00587     if(_solver != 0) {
00588         _solver->setDimension(_variablesCount*2);
00589         _solver->setFunction(solverFunction);
00590         _solver->setParams(this);
00591     }
00592 }
00593 
00594 Solver* World::removeSolver()
00595 {
00596     Solver* solver = _solver;
00597     _solver = NULL;
00598     return solver;
00599 }
00600 
00601 void World::setCollisionSolver(CollisionSolver* collisionSolver)
00602 {
00603     delete _collisionSolver;
00604     _collisionSolver = collisionSolver;
00605 }
00606 
00607 CollisionSolver* World::removeCollisionSolver()
00608 {
00609     CollisionSolver* collisionSolver = _collisionSolver;
00610     _collisionSolver = NULL;
00611     return collisionSolver;
00612 }
00613 
00614 void World::setConstraintSolver(ConstraintSolver* constraintSolver)
00615 {
00616     delete _constraintSolver;
00617     _constraintSolver = constraintSolver;
00618 }
00619 
00620 ConstraintSolver* World::removeConstraintSolver()
00621 {
00622     ConstraintSolver* constraintSolver = _constraintSolver;
00623     _constraintSolver = NULL;
00624     return constraintSolver;
00625 }
00626 
00627 void World::checkVariablesCount()
00628 {
00629     int variablesCount = 0;
00630     for(BodyList::iterator b = _bodies.begin(); b != _bodies.end(); ++b) {
00631         (*b)->setVariablesOffset(variablesCount);
00632         variablesCount += (*b)->variablesCount();
00633     }
00634 
00635     int constraintsCount = 0;
00636     for(JointList::iterator j = _joints.begin(); j != _joints.end(); ++j) {
00637         constraintsCount += (*j)->constraintsCount();
00638     }
00639 
00640     _constraintsInfo.setDimension(variablesCount, constraintsCount);
00641 
00642     if(variablesCount != _variablesCount) {
00643         _variablesCount = variablesCount;
00644         _variables.resize(_variablesCount*2);
00645         _variances.resize(_variablesCount*2);
00646         _tempArray.resize(_variablesCount*2);
00647         if(_solver) _solver->setDimension(_variablesCount*2);
00648     }
00649 }
00650 
00651 void World::gatherAccelerations(double* acceleration, double* accelerationVariance)
00652 {
00653     if(accelerationVariance)
00654         memset(accelerationVariance, 0, _variablesCount*sizeof(*accelerationVariance));
00655 
00656     int index = 0;
00657     const BodyList::const_iterator it_end = _bodies.end();
00658     for(BodyList::iterator b = _bodies.begin(); b != it_end; ++b) {
00659         (*b)->getAccelerations(acceleration + index, accelerationVariance ? accelerationVariance + index : NULL);
00660         index += (*b)->variablesCount();
00661     }
00662 }
00663 
00664 void World::gatherVariables(double* variables, double* variances)
00665 {
00666     if(variances) memset(variances, 0, _variablesCount*2*sizeof(*variances));
00667 
00668     int index = 0;
00669     const BodyList::const_iterator it_end = _bodies.end();
00670     for(BodyList::iterator b = _bodies.begin(); b != it_end; ++b) {
00671         (*b)->getVariables(variables + index, variables + _variablesCount + index,
00672                            variances ? variances + index : NULL,
00673                            variances ? variances + _variablesCount + index : NULL);
00674         index += (*b)->variablesCount();
00675     }
00676 }
00677 
00678 void World::scatterVariables(const double* variables, const double* variances)
00679 {
00680     int index = 0;
00681     const BodyList::const_iterator it_end = _bodies.end();
00682     for(BodyList::iterator b = _bodies.begin(); b != it_end; ++b) {
00683         (*b)->setVariables(variables + index,  variables + _variablesCount + index,
00684                            variances ? variances + index : NULL,
00685                            variances ? variances + _variablesCount + index : NULL);
00686         index += (*b)->variablesCount();
00687     }
00688 }
00689 
00690 void World::gatherJointsInfo(ConstraintsInfo* info)
00691 {
00692     info->clear();
00693 
00694     int offset = 0;
00695     const BodyList::const_iterator b_end = _bodies.end();
00696     for(BodyList::iterator body = _bodies.begin(); body != b_end; ++body) {
00697         (*body)->getInverseMass(&(info->inverseMass), NULL, offset);
00698         offset += (*body)->variablesCount();
00699     }
00700 
00701     offset = 0;
00702     const JointList::const_iterator j_end = _joints.end();
00703     for(JointList::iterator joint = _joints.begin(); joint != j_end; ++joint) {
00704         (*joint)->getConstraintsInfo(info, offset);
00705         offset += (*joint)->constraintsCount();
00706     }
00707 }
00708 
00709 int World::doCalcFn()
00710 {
00711     STEPCORE_ASSERT_NOABORT(_solver != NULL);
00712 
00713     //if(_collisionSolver) _collisionSolver->resetCaches();
00714 
00715     _stopOnCollision = false;
00716     _stopOnIntersection = false;
00717     checkVariablesCount();
00718     double* variances = _errorsCalculation ? &_variances[0] : NULL;
00719     gatherVariables(&_variables[0], variances);
00720     return _solver->doCalcFn(&_time, &_variables[0], variances, NULL, variances);
00721 }
00722 
00723 int World::doEvolve(double delta)
00724 {
00725     STEPCORE_ASSERT_NOABORT(_solver != NULL);
00726 
00727     checkVariablesCount();
00728     gatherVariables(&_variables[0], _errorsCalculation ? &_variances[0] : NULL);
00729 
00730     int ret = Solver::OK;
00731     double targetTime = _time + delta*_timeScale;
00732     
00733     if(_collisionSolver) {
00734         //_collisionSolver->resetCaches();
00735         if(Contact::Intersected == _collisionSolver->checkContacts(_bodies, NULL))
00736             return Solver::IntersectionDetected;
00737     }
00738 
00739     while(_time < targetTime) {
00740         STEPCORE_ASSERT_NOABORT( targetTime - _time > _solver->stepSize() / 1000 );
00741         if( !(   targetTime - _time > _solver->stepSize() / 1000 ) ) {
00742                     qDebug("* %e %e %e", targetTime, _time, _solver->stepSize());
00743         }
00744         double time = _time;
00745 
00746         _stopOnCollision = true;
00747         _stopOnIntersection = true;
00748         ret = _solver->doEvolve(&time, targetTime, &_variables[0],
00749                             _errorsCalculation ? &_variances[0] : NULL);
00750         _time = time;
00751 
00752         if(ret == Solver::CollisionDetected ||
00753            ret == Solver::IntersectionDetected) {
00754             // If we have stopped on collision
00755             // 1. Decrease timestep to stop before collision
00756             // 2. Proceed with decresed timestep until
00757             //    - we have meet collision again: go to 1
00758             //    - we pass collision point: it means that we have come close enough
00759             //      to collision point and CollisionSolver have resolved collision
00760             // We can't simply change Solver::stepSize since adaptive solvers can
00761             // abuse our settings so we have to step manually
00762             //STEPCORE_ASSERT_NOABORT(_collisionTime <= targetTime);
00763             //STEPCORE_ASSERT_NOABORT(_collisionTime > _time);
00764             double stepSize = fmin(_solver->stepSize() / 2, targetTime - _time);
00765             double collisionEndTime = targetTime - _time > stepSize*3.01 ? _time + stepSize*3 : targetTime;
00766 
00767             _stopOnCollision = false;
00768 
00769             do {
00770                 double endTime = collisionEndTime - time > stepSize*1.01 ? time+stepSize : collisionEndTime;
00771 
00772                 ret = _solver->doEvolve(&time, endTime, &_variables[0],
00773                             _errorsCalculation ? &_variances[0] : NULL);
00774                 _time = time;
00775 
00776                 if(ret == Solver::IntersectionDetected || ret == Solver::CollisionDetected) {
00777                     //STEPCORE_ASSERT_NOABORT(_collisionTime > _time);
00778                     //STEPCORE_ASSERT_NOABORT(_collisionTime < _collisionExpectedTime);
00779                     stepSize = fmin(stepSize/2, targetTime - _time);
00780                     collisionEndTime = targetTime - _time > stepSize*3.01 ? _time + stepSize*3 : targetTime;
00781 
00782                     //STEPCORE_ASSERT_NOABORT(_time + stepSize != _time);
00783                     // XXX: what to do if stepSize becomes too small ?
00784 
00785                 } else if(ret == Solver::OK) {
00786                     // We can be close to the collision point.
00787                     //
00788                     // Now we will calculate impulses required to fix collisions.
00789                     // All joints should be taken into account, but since we are
00790                     // calculating impulses we should "froze" the jacobian i.e.
00791                     // ignore it derivative
00792                     scatterVariables(&_variables[0], _errorsCalculation ? &_variances[0] : NULL);
00793 
00794                     std::fill(_tempArray.begin(), _tempArray.end(), 0);
00795                     _constraintsInfo.position = GmmArrayVector(const_cast<double*>(&_variables[0]), _variablesCount);
00796                     _constraintsInfo.velocity = GmmArrayVector(const_cast<double*>(&_variables[0]+_variablesCount), _variablesCount);
00797                     _constraintsInfo.acceleration = GmmArrayVector(const_cast<double*>(&_tempArray[0]), _variablesCount);
00798 
00799                     gatherJointsInfo(&_constraintsInfo);
00800                     gmm::clear(_constraintsInfo.jacobianDerivative);
00801 
00802                     int status = _collisionSolver->checkContacts(_bodies, &_constraintsInfo, true);
00803                     if(status >= CollisionSolver::InternalError) { ret = status; goto out; }
00804 
00805                         //GmmArrayVector cforce(const_cast<double*>(&_constraintsTotalForce[0]), _variablesCount);
00806                     if(_constraintsInfo.collisionFlag) {
00807                         ret = _constraintSolver->solve(&_constraintsInfo);
00808                         if(ret != Solver::OK) goto out;
00809 
00810                         // XXX: variances
00811                         gmm::mult_add(_constraintsInfo.inverseMass, _constraintsInfo.force,
00812                                                     _constraintsInfo.velocity);
00813                     }
00814 
00815                 } else goto out;
00816 
00817             } while(_time + stepSize/100 <= collisionEndTime); // XXX
00818         } else if(ret != Solver::OK) goto out;
00819     }
00820 
00821 out:
00822     scatterVariables(&_variables[0], _errorsCalculation ? &_variances[0] : NULL);
00823     return ret;
00824 }
00825 
00826 inline int World::solverFunction(double t, const double* y,
00827                     const double* yvar, double* f, double* fvar)
00828 {
00829     if(_evolveAbort) return Solver::Aborted;
00830 
00831     _time = t;
00832     scatterVariables(y, yvar); // this will reset force
00833 
00834     // 1. Forces
00835     bool calcVariances = (fvar != NULL);
00836     const ForceList::const_iterator f_end = _forces.end();
00837     for(ForceList::iterator force = _forces.begin(); force != f_end; ++force) {
00838         (*force)->calcForce(calcVariances);
00839     }
00840 
00841     std::memcpy(f, y+_variablesCount, _variablesCount*sizeof(*f));
00842     if(fvar) std::memcpy(fvar, y+_variablesCount, _variablesCount*sizeof(*fvar));
00843     gatherAccelerations(f+_variablesCount, fvar ? fvar+_variablesCount : NULL);
00844 
00845     // 2. Joints
00846     if(_constraintSolver) {
00847         _constraintsInfo.position = GmmArrayVector(const_cast<double*>(y), _variablesCount);
00848         _constraintsInfo.velocity = GmmArrayVector(const_cast<double*>(y+_variablesCount), _variablesCount);
00849         _constraintsInfo.acceleration = GmmArrayVector(const_cast<double*>(f+_variablesCount), _variablesCount);
00850 
00851         gatherJointsInfo(&_constraintsInfo);
00852     }
00853 
00854     // 3. Collisions (TODO: variances for collisions)
00855     if(_collisionSolver && _constraintSolver) {
00856         int state = _collisionSolver->checkContacts(_bodies, &_constraintsInfo, false);
00857         if(state == Contact::Intersected) {
00858             if(_stopOnIntersection) return Solver::IntersectionDetected;
00859         } else if(state == Contact::Colliding) {
00860             if(_stopOnCollision) return Solver::CollisionDetected;
00861             // XXX: We are not stopping on colliding contact
00862             // and resolving them only at the end of timestep
00863             // XXX: is it right solution ? Shouldn't we try to find
00864             // contact point more exactly for example using binary search ?
00865             //_collisionTime = t;
00866             //_collisionTime = t;
00867             //if(t < _collisionExpectedTime)
00868             //    return DantzigLCPCollisionSolver::CollisionDetected;
00869         } else if(state >= CollisionSolver::InternalError) {
00870             return state;
00871         }
00872     }
00873 
00874     // 4. Solve constraints
00875     if(_constraintSolver &&
00876             _constraintsInfo.constraintsCount + _constraintsInfo.contactsCount > 0) {
00877 
00878         int state = _constraintSolver->solve(&_constraintsInfo);
00879         if(state != Solver::OK) return state;
00880 
00881         int offset = 0;
00882         const BodyList::const_iterator b_end = _bodies.end();
00883         for(BodyList::iterator body = _bodies.begin(); body != b_end; ++body) {
00884             (*body)->addForce(&_constraintsInfo.force[offset], NULL);
00885             (*body)->getAccelerations(f+_variablesCount+offset, NULL);
00886             offset += (*body)->variablesCount();
00887         }
00888     }
00889 
00890     return 0;
00891 }
00892 
00893 int World::solverFunction(double t, const double* y,
00894                 const double* yvar, double* f, double* fvar, void* params)
00895 {
00896     return static_cast<World*>(params)->solverFunction(t, y, yvar, f, fvar);
00897 }
00898 
00899 } // namespace StepCore
00900 

step/stepcore

Skip menu "step/stepcore"
  • Main Page
  • Modules
  • Namespace List
  • Class Hierarchy
  • Alphabetical List
  • Class List
  • File List
  • Namespace Members
  • Class Members
  • Related Pages

kdeedu

Skip menu "kdeedu"
  • kalzium
  • kanagram
  • kig
  •   lib
  • klettres
  • kstars
  • libkdeedu
  •   keduvocdocument
  •   docs
  •   src
  • parley
  •   stepcore
Generated for kdeedu by doxygen 1.5.4
This website is maintained by Adriaan de Groot and Allen Winter.
KDE® and the K Desktop Environment® logo are registered trademarks of KDE e.V. | Legal