/************************************************************************* * * * Open Dynamics Engine, Copyright (C) 2001,2002 Russell L. Smith. * * All rights reserved. Email: russ@q12.org Web: www.q12.org * * * * This library is free software; you can redistribute it and/or * * modify it under the terms of EITHER: * * (1) The GNU Lesser General Public License as published by the Free * * Software Foundation; either version 2.1 of the License, or (at * * your option) any later version. The text of the GNU Lesser * * General Public License is included with this library in the * * file LICENSE.TXT. * * (2) The BSD-style license that is included with this library in * * the file LICENSE-BSD.TXT. * * * * This library is distributed in the hope that it will be useful, * * but WITHOUT ANY WARRANTY; without even the implied warranty of * * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the files * * LICENSE.TXT and LICENSE-BSD.TXT for more details. * * * *************************************************************************/ #include "objects.h" #include "joint.h" #include #include #include #include #include #include #include "lcp.h" #include "util.h" //**************************************************************************** // misc defines #define FAST_FACTOR //#define TIMING // memory allocation system #ifdef dUSE_MALLOC_FOR_ALLOCA unsigned int dMemoryFlag; #define REPORT_OUT_OF_MEMORY fprintf(stderr, "Insufficient memory to complete rigid body simulation. Results will not be accurate.\n") #define ALLOCA(t,v,s) t* v=(t*)malloc(s) #define UNALLOCA(t) free(t) #else #define ALLOCA(t,v,s) t* v=(t*)dALLOCA16(s) #define UNALLOCA(t) /* nothing */ #endif //**************************************************************************** // debugging - comparison of various vectors and matrices produced by the // slow and fast versions of the stepper. //#define COMPARE_METHODS #ifdef COMPARE_METHODS #include "testing.h" dMatrixComparison comparator; #endif // undef to use the fast decomposition #define DIRECT_CHOLESKY #undef REPORT_ERROR //**************************************************************************** // special matrix multipliers // this assumes the 4th and 8th rows of B and C are zero. static void Multiply2_p8r (dReal *A, dReal *B, dReal *C, int p, int r, int Askip) { int i,j; dReal sum,*bb,*cc; dIASSERT (p>0 && r>0 && A && B && C); bb = B; for (i=p; i; i--) { cc = C; for (j=r; j; j--) { sum = bb[0]*cc[0]; sum += bb[1]*cc[1]; sum += bb[2]*cc[2]; sum += bb[4]*cc[4]; sum += bb[5]*cc[5]; sum += bb[6]*cc[6]; *(A++) = sum; cc += 8; } A += Askip - r; bb += 8; } } // this assumes the 4th and 8th rows of B and C are zero. static void MultiplyAdd2_p8r (dReal *A, dReal *B, dReal *C, int p, int r, int Askip) { int i,j; dReal sum,*bb,*cc; dIASSERT (p>0 && r>0 && A && B && C); bb = B; for (i=p; i; i--) { cc = C; for (j=r; j; j--) { sum = bb[0]*cc[0]; sum += bb[1]*cc[1]; sum += bb[2]*cc[2]; sum += bb[4]*cc[4]; sum += bb[5]*cc[5]; sum += bb[6]*cc[6]; *(A++) += sum; cc += 8; } A += Askip - r; bb += 8; } } // this assumes the 4th and 8th rows of B are zero. static void Multiply0_p81 (dReal *A, dReal *B, dReal *C, int p) { int i; dIASSERT (p>0 && A && B && C); dReal sum; for (i=p; i; i--) { sum = B[0]*C[0]; sum += B[1]*C[1]; sum += B[2]*C[2]; sum += B[4]*C[4]; sum += B[5]*C[5]; sum += B[6]*C[6]; *(A++) = sum; B += 8; } } // this assumes the 4th and 8th rows of B are zero. static void MultiplyAdd0_p81 (dReal *A, dReal *B, dReal *C, int p) { int i; dIASSERT (p>0 && A && B && C); dReal sum; for (i=p; i; i--) { sum = B[0]*C[0]; sum += B[1]*C[1]; sum += B[2]*C[2]; sum += B[4]*C[4]; sum += B[5]*C[5]; sum += B[6]*C[6]; *(A++) += sum; B += 8; } } // this assumes the 4th and 8th rows of B are zero. static void MultiplyAdd1_8q1 (dReal *A, dReal *B, dReal *C, int q) { int k; dReal sum; dIASSERT (q>0 && A && B && C); sum = 0; for (k=0; k0 && A && B && C); sum = 0; for (k=0; ktag = i; // make a local copy of the joint array, because we might want to modify it. // (the "dxJoint *const*" declaration says we're allowed to modify the joints // but not the joint array, because the caller might need it unchanged). ALLOCA(dxJoint*,joint,nj*sizeof(dxJoint*)); #ifdef dUSE_MALLOC_FOR_ALLOCA if (joint == NULL) { dMemoryFlag = d_MEMORY_OUT_OF_MEMORY; return; } #endif memcpy (joint,_joint,nj * sizeof(dxJoint*)); // for all bodies, compute the inertia tensor and its inverse in the global // frame, and compute the rotational force and add it to the torque // accumulator. // @@@ check computation of rotational force. ALLOCA(dReal,I,3*nb*4*sizeof(dReal)); #ifdef dUSE_MALLOC_FOR_ALLOCA if (I == NULL) { UNALLOCA(joint); dMemoryFlag = d_MEMORY_OUT_OF_MEMORY; return; } #endif ALLOCA(dReal,invI,3*nb*4*sizeof(dReal)); #ifdef dUSE_MALLOC_FOR_ALLOCA if (invI == NULL) { UNALLOCA(I); UNALLOCA(joint); dMemoryFlag = d_MEMORY_OUT_OF_MEMORY; return; } #endif //dSetZero (I,3*nb*4); //dSetZero (invI,3*nb*4); for (i=0; imass.I,body[i]->posr.R); dMULTIPLY0_333 (I+i*12,body[i]->posr.R,tmp); // compute inverse inertia tensor in global frame dMULTIPLY2_333 (tmp,body[i]->invI,body[i]->posr.R); dMULTIPLY0_333 (invI+i*12,body[i]->posr.R,tmp); // compute rotational force dMULTIPLY0_331 (tmp,I+i*12,body[i]->avel); dCROSS (body[i]->tacc,-=,body[i]->avel,tmp); } // add the gravity force to all bodies for (i=0; iflags & dxBodyNoGravity)==0) { body[i]->facc[0] += body[i]->mass.mass * world->gravity[0]; body[i]->facc[1] += body[i]->mass.mass * world->gravity[1]; body[i]->facc[2] += body[i]->mass.mass * world->gravity[2]; } } // get m = total constraint dimension, nub = number of unbounded variables. // create constraint offset array and number-of-rows array for all joints. // the constraints are re-ordered as follows: the purely unbounded // constraints, the mixed unbounded + LCP constraints, and last the purely // LCP constraints. // // joints with m=0 are inactive and are removed from the joints array // entirely, so that the code that follows does not consider them. int m = 0; ALLOCA(dxJoint::Info1,info,nj*sizeof(dxJoint::Info1)); #ifdef dUSE_MALLOC_FOR_ALLOCA if (info == NULL) { UNALLOCA(invI); UNALLOCA(I); UNALLOCA(joint); dMemoryFlag = d_MEMORY_OUT_OF_MEMORY; return; } #endif ALLOCA(int,ofs,nj*sizeof(int)); #ifdef dUSE_MALLOC_FOR_ALLOCA if (ofs == NULL) { UNALLOCA(info); UNALLOCA(invI); UNALLOCA(I); UNALLOCA(joint); dMemoryFlag = d_MEMORY_OUT_OF_MEMORY; return; } #endif for (i=0, j=0; jvtable->getInfo1 (joint[j],info+i); dIASSERT (info[i].m >= 0 && info[i].m <= 6 && info[i].nub >= 0 && info[i].nub <= info[i].m); if (info[i].m > 0) { joint[i] = joint[j]; i++; } } nj = i; // the purely unbounded constraints for (i=0; i 0 && info[i].nub < info[i].m) { ofs[i] = m; m += info[i].m; } // the purely LCP constraints for (i=0; iinvMass; MM[nskip+1] = body[i]->invMass; MM[2*nskip+2] = body[i]->invMass; MM += 3*nskip+3; for (j=0; j<3; j++) for (k=0; k<3; k++) { MM[j*nskip+k] = invI[i*12+j*4+k]; } } // assemble some body vectors: fe = external forces, v = velocities ALLOCA(dReal,fe,n6*sizeof(dReal)); #ifdef dUSE_MALLOC_FOR_ALLOCA if (fe == NULL) { UNALLOCA(invM); UNALLOCA(ofs); UNALLOCA(info); UNALLOCA(invI); UNALLOCA(I); UNALLOCA(joint); dMemoryFlag = d_MEMORY_OUT_OF_MEMORY; return; } #endif ALLOCA(dReal,v,n6*sizeof(dReal)); #ifdef dUSE_MALLOC_FOR_ALLOCA if (v == NULL) { UNALLOCA(fe); UNALLOCA(invM); UNALLOCA(ofs); UNALLOCA(info); UNALLOCA(invI); UNALLOCA(I); UNALLOCA(joint); dMemoryFlag = d_MEMORY_OUT_OF_MEMORY; return; } #endif //dSetZero (fe,n6); //dSetZero (v,n6); for (i=0; ifacc[j]; for (j=0; j<3; j++) fe[i*6+3+j] = body[i]->tacc[j]; for (j=0; j<3; j++) v[i*6+j] = body[i]->lvel[j]; for (j=0; j<3; j++) v[i*6+3+j] = body[i]->avel[j]; } // this will be set to the velocity update ALLOCA(dReal,vnew,n6*sizeof(dReal)); #ifdef dUSE_MALLOC_FOR_ALLOCA if (vnew == NULL) { UNALLOCA(v); UNALLOCA(fe); UNALLOCA(invM); UNALLOCA(ofs); UNALLOCA(info); UNALLOCA(invI); UNALLOCA(I); UNALLOCA(joint); dMemoryFlag = d_MEMORY_OUT_OF_MEMORY; return; } #endif dSetZero (vnew,n6); // if there are constraints, compute cforce if (m > 0) { // create a constraint equation right hand side vector `c', a constraint // force mixing vector `cfm', and LCP low and high bound vectors, and an // 'findex' vector. ALLOCA(dReal,c,m*sizeof(dReal)); #ifdef dUSE_MALLOC_FOR_ALLOCA if (c == NULL) { UNALLOCA(vnew); UNALLOCA(v); UNALLOCA(fe); UNALLOCA(invM); UNALLOCA(ofs); UNALLOCA(info); UNALLOCA(invI); UNALLOCA(I); UNALLOCA(joint); dMemoryFlag = d_MEMORY_OUT_OF_MEMORY; return; } #endif ALLOCA(dReal,cfm,m*sizeof(dReal)); #ifdef dUSE_MALLOC_FOR_ALLOCA if (cfm == NULL) { UNALLOCA(c); UNALLOCA(vnew); UNALLOCA(v); UNALLOCA(fe); UNALLOCA(invM); UNALLOCA(ofs); UNALLOCA(info); UNALLOCA(invI); UNALLOCA(I); UNALLOCA(joint); dMemoryFlag = d_MEMORY_OUT_OF_MEMORY; return; } #endif ALLOCA(dReal,lo,m*sizeof(dReal)); #ifdef dUSE_MALLOC_FOR_ALLOCA if (lo == NULL) { UNALLOCA(cfm); UNALLOCA(c); UNALLOCA(vnew); UNALLOCA(v); UNALLOCA(fe); UNALLOCA(invM); UNALLOCA(ofs); UNALLOCA(info); UNALLOCA(invI); UNALLOCA(I); UNALLOCA(joint); dMemoryFlag = d_MEMORY_OUT_OF_MEMORY; return; } #endif ALLOCA(dReal,hi,m*sizeof(dReal)); #ifdef dUSE_MALLOC_FOR_ALLOCA if (hi == NULL) { UNALLOCA(lo); UNALLOCA(cfm); UNALLOCA(c); UNALLOCA(vnew); UNALLOCA(v); UNALLOCA(fe); UNALLOCA(invM); UNALLOCA(ofs); UNALLOCA(info); UNALLOCA(invI); UNALLOCA(I); UNALLOCA(joint); dMemoryFlag = d_MEMORY_OUT_OF_MEMORY; return; } #endif ALLOCA(int,findex,m*sizeof(int)); #ifdef dUSE_MALLOC_FOR_ALLOCA if (findex == NULL) { UNALLOCA(hi); UNALLOCA(lo); UNALLOCA(cfm); UNALLOCA(c); UNALLOCA(vnew); UNALLOCA(v); UNALLOCA(fe); UNALLOCA(invM); UNALLOCA(ofs); UNALLOCA(info); UNALLOCA(invI); UNALLOCA(I); UNALLOCA(joint); dMemoryFlag = d_MEMORY_OUT_OF_MEMORY; return; } #endif dSetZero (c,m); dSetValue (cfm,m,world->global_cfm); dSetValue (lo,m,-dInfinity); dSetValue (hi,m, dInfinity); for (i=0; iglobal_erp; for (i=0; inode[0].body->tag; Jinfo.J1a = Jinfo.J1l + 3; if (joint[i]->node[1].body) { Jinfo.J2l = J + nskip*ofs[i] + 6*joint[i]->node[1].body->tag; Jinfo.J2a = Jinfo.J2l + 3; } else { Jinfo.J2l = 0; Jinfo.J2a = 0; } Jinfo.c = c + ofs[i]; Jinfo.cfm = cfm + ofs[i]; Jinfo.lo = lo + ofs[i]; Jinfo.hi = hi + ofs[i]; Jinfo.findex = findex + ofs[i]; joint[i]->vtable->getInfo2 (joint[i],&Jinfo); // adjust returned findex values for global index numbering for (j=0; j= 0) findex[ofs[i] + j] += ofs[i]; } } // compute A = J*invM*J' # ifdef TIMING dTimerNow ("compute A"); # endif ALLOCA(dReal,JinvM,m*nskip*sizeof(dReal)); #ifdef dUSE_MALLOC_FOR_ALLOCA if (JinvM == NULL) { UNALLOCA(J); UNALLOCA(findex); UNALLOCA(hi); UNALLOCA(lo); UNALLOCA(cfm); UNALLOCA(c); UNALLOCA(vnew); UNALLOCA(v); UNALLOCA(fe); UNALLOCA(invM); UNALLOCA(ofs); UNALLOCA(info); UNALLOCA(invI); UNALLOCA(I); UNALLOCA(joint); dMemoryFlag = d_MEMORY_OUT_OF_MEMORY; return; } #endif //dSetZero (JinvM,m*nskip); dMultiply0 (JinvM,J,invM,m,n6,n6); int mskip = dPAD(m); ALLOCA(dReal,A,m*mskip*sizeof(dReal)); #ifdef dUSE_MALLOC_FOR_ALLOCA if (A == NULL) { UNALLOCA(JinvM); UNALLOCA(J); UNALLOCA(findex); UNALLOCA(hi); UNALLOCA(lo); UNALLOCA(cfm); UNALLOCA(c); UNALLOCA(vnew); UNALLOCA(v); UNALLOCA(fe); UNALLOCA(invM); UNALLOCA(ofs); UNALLOCA(info); UNALLOCA(invI); UNALLOCA(I); UNALLOCA(joint); dMemoryFlag = d_MEMORY_OUT_OF_MEMORY; return; } #endif //dSetZero (A,m*mskip); dMultiply2 (A,JinvM,J,m,n6,m); // add cfm to the diagonal of A for (i=0; ilvel[j] = vnew[i*6+j]; for (j=0; j<3; j++) body[i]->avel[j] = vnew[i*6+3+j]; } // update the position and orientation from the new linear/angular velocity // (over the given timestep) #ifdef TIMING dTimerNow ("update position"); #endif for (i=0; ifacc[0] = 0; body[i]->facc[1] = 0; body[i]->facc[2] = 0; body[i]->facc[3] = 0; body[i]->tacc[0] = 0; body[i]->tacc[1] = 0; body[i]->tacc[2] = 0; body[i]->tacc[3] = 0; } #ifdef TIMING dTimerEnd(); if (m > 0) dTimerReport (stdout,1); #endif UNALLOCA(joint); UNALLOCA(I); UNALLOCA(invI); UNALLOCA(info); UNALLOCA(ofs); UNALLOCA(invM); UNALLOCA(fe); UNALLOCA(v); UNALLOCA(vnew); } //**************************************************************************** // an optimized version of dInternalStepIsland1() void dInternalStepIsland_x2 (dxWorld *world, dxBody * const *body, int nb, dxJoint * const *_joint, int nj, dReal stepsize) { int i,j,k; #ifdef TIMING dTimerStart("preprocessing"); #endif dReal stepsize1 = dRecip(stepsize); // number all bodies in the body list - set their tag values for (i=0; itag = i; // make a local copy of the joint array, because we might want to modify it. // (the "dxJoint *const*" declaration says we're allowed to modify the joints // but not the joint array, because the caller might need it unchanged). ALLOCA(dxJoint*,joint,nj*sizeof(dxJoint*)); #ifdef dUSE_MALLOC_FOR_ALLOCA if (joint == NULL) { dMemoryFlag = d_MEMORY_OUT_OF_MEMORY; return; } #endif memcpy (joint,_joint,nj * sizeof(dxJoint*)); // for all bodies, compute the inertia tensor and its inverse in the global // frame, and compute the rotational force and add it to the torque // accumulator. I and invI are vertically stacked 3x4 matrices, one per body. // @@@ check computation of rotational force. #ifdef dGYROSCOPIC ALLOCA(dReal,I,3*nb*4*sizeof(dReal)); # ifdef dUSE_MALLOC_FOR_ALLOCA if (I == NULL) { UNALLOCA(joint); dMemoryFlag = d_MEMORY_OUT_OF_MEMORY; return; } # endif #else dReal *I = NULL; #endif // for dGYROSCOPIC ALLOCA(dReal,invI,3*nb*4*sizeof(dReal)); #ifdef dUSE_MALLOC_FOR_ALLOCA if (invI == NULL) { UNALLOCA(I); UNALLOCA(joint); dMemoryFlag = d_MEMORY_OUT_OF_MEMORY; return; } #endif //dSetZero (I,3*nb*4); //dSetZero (invI,3*nb*4); for (i=0; iinvI,body[i]->posr.R); dMULTIPLY0_333 (invI+i*12,body[i]->posr.R,tmp); #ifdef dGYROSCOPIC // compute inertia tensor in global frame dMULTIPLY2_333 (tmp,body[i]->mass.I,body[i]->posr.R); dMULTIPLY0_333 (I+i*12,body[i]->posr.R,tmp); // compute rotational force dMULTIPLY0_331 (tmp,I+i*12,body[i]->avel); dCROSS (body[i]->tacc,-=,body[i]->avel,tmp); #endif } // add the gravity force to all bodies for (i=0; iflags & dxBodyNoGravity)==0) { body[i]->facc[0] += body[i]->mass.mass * world->gravity[0]; body[i]->facc[1] += body[i]->mass.mass * world->gravity[1]; body[i]->facc[2] += body[i]->mass.mass * world->gravity[2]; } } // get m = total constraint dimension, nub = number of unbounded variables. // create constraint offset array and number-of-rows array for all joints. // the constraints are re-ordered as follows: the purely unbounded // constraints, the mixed unbounded + LCP constraints, and last the purely // LCP constraints. this assists the LCP solver to put all unbounded // variables at the start for a quick factorization. // // joints with m=0 are inactive and are removed from the joints array // entirely, so that the code that follows does not consider them. // also number all active joints in the joint list (set their tag values). // inactive joints receive a tag value of -1. int m = 0; ALLOCA(dxJoint::Info1,info,nj*sizeof(dxJoint::Info1)); #ifdef dUSE_MALLOC_FOR_ALLOCA if (info == NULL) { UNALLOCA(invI); UNALLOCA(I); UNALLOCA(joint); dMemoryFlag = d_MEMORY_OUT_OF_MEMORY; return; } #endif ALLOCA(int,ofs,nj*sizeof(int)); #ifdef dUSE_MALLOC_FOR_ALLOCA if (ofs == NULL) { UNALLOCA(info); UNALLOCA(invI); UNALLOCA(I); UNALLOCA(joint); dMemoryFlag = d_MEMORY_OUT_OF_MEMORY; return; } #endif for (i=0, j=0; jvtable->getInfo1 (joint[j],info+i); dIASSERT (info[i].m >= 0 && info[i].m <= 6 && info[i].nub >= 0 && info[i].nub <= info[i].m); if (info[i].m > 0) { joint[i] = joint[j]; joint[i]->tag = i; i++; } else { joint[j]->tag = -1; } } nj = i; // the purely unbounded constraints for (i=0; i 0 && info[i].nub < info[i].m) { ofs[i] = m; m += info[i].m; } // the purely LCP constraints for (i=0; i 0) { // create a constraint equation right hand side vector `c', a constraint // force mixing vector `cfm', and LCP low and high bound vectors, and an // 'findex' vector. ALLOCA(dReal,c,m*sizeof(dReal)); #ifdef dUSE_MALLOC_FOR_ALLOCA if (c == NULL) { UNALLOCA(cforce); UNALLOCA(ofs); UNALLOCA(info); UNALLOCA(invI); UNALLOCA(I); UNALLOCA(joint); dMemoryFlag = d_MEMORY_OUT_OF_MEMORY; return; } #endif ALLOCA(dReal,cfm,m*sizeof(dReal)); #ifdef dUSE_MALLOC_FOR_ALLOCA if (cfm == NULL) { UNALLOCA(c); UNALLOCA(cforce); UNALLOCA(ofs); UNALLOCA(info); UNALLOCA(invI); UNALLOCA(I); UNALLOCA(joint); dMemoryFlag = d_MEMORY_OUT_OF_MEMORY; return; } #endif ALLOCA(dReal,lo,m*sizeof(dReal)); #ifdef dUSE_MALLOC_FOR_ALLOCA if (lo == NULL) { UNALLOCA(cfm); UNALLOCA(c); UNALLOCA(cforce); UNALLOCA(ofs); UNALLOCA(info); UNALLOCA(invI); UNALLOCA(I); UNALLOCA(joint); dMemoryFlag = d_MEMORY_OUT_OF_MEMORY; return; } #endif ALLOCA(dReal,hi,m*sizeof(dReal)); #ifdef dUSE_MALLOC_FOR_ALLOCA if (hi == NULL) { UNALLOCA(lo); UNALLOCA(cfm); UNALLOCA(c); UNALLOCA(cforce); UNALLOCA(ofs); UNALLOCA(info); UNALLOCA(invI); UNALLOCA(I); UNALLOCA(joint); dMemoryFlag = d_MEMORY_OUT_OF_MEMORY; return; } #endif ALLOCA(int,findex,m*sizeof(int)); #ifdef dUSE_MALLOC_FOR_ALLOCA if (findex == NULL) { UNALLOCA(hi); UNALLOCA(lo); UNALLOCA(cfm); UNALLOCA(c); UNALLOCA(cforce); UNALLOCA(ofs); UNALLOCA(info); UNALLOCA(invI); UNALLOCA(I); UNALLOCA(joint); dMemoryFlag = d_MEMORY_OUT_OF_MEMORY; return; } #endif dSetZero (c,m); dSetValue (cfm,m,world->global_cfm); dSetValue (lo,m,-dInfinity); dSetValue (hi,m, dInfinity); for (i=0; iglobal_erp; for (i=0; ivtable->getInfo2 (joint[i],&Jinfo); // adjust returned findex values for global index numbering for (j=0; j= 0) findex[ofs[i] + j] += ofs[i]; } } // compute A = J*invM*J'. first compute JinvM = J*invM. this has the same // format as J so we just go through the constraints in J multiplying by // the appropriate scalars and matrices. # ifdef TIMING dTimerNow ("compute A"); # endif ALLOCA(dReal,JinvM,2*m*8*sizeof(dReal)); #ifdef dUSE_MALLOC_FOR_ALLOCA if (JinvM == NULL) { UNALLOCA(J); UNALLOCA(findex); UNALLOCA(hi); UNALLOCA(lo); UNALLOCA(cfm); UNALLOCA(c); UNALLOCA(cforce); UNALLOCA(ofs); UNALLOCA(info); UNALLOCA(invI); UNALLOCA(I); UNALLOCA(joint); dMemoryFlag = d_MEMORY_OUT_OF_MEMORY; return; } #endif dSetZero (JinvM,2*m*8); for (i=0; inode[0].body->tag; dReal body_invMass = body[b]->invMass; dReal *body_invI = invI + b*12; dReal *Jsrc = J + 2*8*ofs[i]; dReal *Jdst = JinvM + 2*8*ofs[i]; for (j=info[i].m-1; j>=0; j--) { for (k=0; k<3; k++) Jdst[k] = Jsrc[k] * body_invMass; dMULTIPLY0_133 (Jdst+4,Jsrc+4,body_invI); Jsrc += 8; Jdst += 8; } if (joint[i]->node[1].body) { b = joint[i]->node[1].body->tag; body_invMass = body[b]->invMass; body_invI = invI + b*12; for (j=info[i].m-1; j>=0; j--) { for (k=0; k<3; k++) Jdst[k] = Jsrc[k] * body_invMass; dMULTIPLY0_133 (Jdst+4,Jsrc+4,body_invI); Jsrc += 8; Jdst += 8; } } } // now compute A = JinvM * J'. A's rows and columns are grouped by joint, // i.e. in the same way as the rows of J. block (i,j) of A is only nonzero // if joints i and j have at least one body in common. this fact suggests // the algorithm used to fill A: // // for b = all bodies // n = number of joints attached to body b // for i = 1..n // for j = i+1..n // ii = actual joint number for i // jj = actual joint number for j // // (ii,jj) will be set to all pairs of joints around body b // compute blockwise: A(ii,jj) += JinvM(ii) * J(jj)' // // this algorithm catches all pairs of joints that have at least one body // in common. it does not compute the diagonal blocks of A however - // another similar algorithm does that. int mskip = dPAD(m); ALLOCA(dReal,A,m*mskip*sizeof(dReal)); #ifdef dUSE_MALLOC_FOR_ALLOCA if (A == NULL) { UNALLOCA(JinvM); UNALLOCA(J); UNALLOCA(findex); UNALLOCA(hi); UNALLOCA(lo); UNALLOCA(cfm); UNALLOCA(c); UNALLOCA(cforce); UNALLOCA(ofs); UNALLOCA(info); UNALLOCA(invI); UNALLOCA(I); UNALLOCA(joint); dMemoryFlag = d_MEMORY_OUT_OF_MEMORY; return; } #endif dSetZero (A,m*mskip); for (i=0; ifirstjoint; n1; n1=n1->next) { for (dxJointNode *n2=n1->next; n2; n2=n2->next) { // get joint numbers and ensure ofs[j1] >= ofs[j2] int j1 = n1->joint->tag; int j2 = n2->joint->tag; if (ofs[j1] < ofs[j2]) { int tmp = j1; j1 = j2; j2 = tmp; } // if either joint was tagged as -1 then it is an inactive (m=0) // joint that should not be considered if (j1==-1 || j2==-1) continue; // determine if body i is the 1st or 2nd body of joints j1 and j2 int jb1 = (joint[j1]->node[1].body == body[i]); int jb2 = (joint[j2]->node[1].body == body[i]); // jb1/jb2 must be 0 for joints with only one body dIASSERT(joint[j1]->node[1].body || jb1==0); dIASSERT(joint[j2]->node[1].body || jb2==0); // set block of A MultiplyAdd2_p8r (A + ofs[j1]*mskip + ofs[j2], JinvM + 2*8*ofs[j1] + jb1*8*info[j1].m, J + 2*8*ofs[j2] + jb2*8*info[j2].m, info[j1].m,info[j2].m, mskip); } } } // compute diagonal blocks of A for (i=0; inode[1].body) { MultiplyAdd2_p8r (A + ofs[i]*(mskip+1), JinvM + 2*8*ofs[i] + 8*info[i].m, J + 2*8*ofs[i] + 8*info[i].m, info[i].m,info[i].m, mskip); } } // add cfm to the diagonal of A for (i=0; iinvMass; dReal *body_invI = invI + i*12; for (j=0; j<3; j++) tmp1[i*8+j] = body[i]->facc[j] * body_invMass + body[i]->lvel[j] * stepsize1; dMULTIPLY0_331 (tmp1 + i*8 + 4,body_invI,body[i]->tacc); for (j=0; j<3; j++) tmp1[i*8+4+j] += body[i]->avel[j] * stepsize1; } // put J*tmp1 into rhs ALLOCA(dReal,rhs,m*sizeof(dReal)); #ifdef dUSE_MALLOC_FOR_ALLOCA if (rhs == NULL) { UNALLOCA(tmp1); UNALLOCA(A); UNALLOCA(JinvM); UNALLOCA(J); UNALLOCA(findex); UNALLOCA(hi); UNALLOCA(lo); UNALLOCA(cfm); UNALLOCA(c); UNALLOCA(cforce); UNALLOCA(ofs); UNALLOCA(info); UNALLOCA(invI); UNALLOCA(I); UNALLOCA(joint); dMemoryFlag = d_MEMORY_OUT_OF_MEMORY; return; } #endif //dSetZero (rhs,m); for (i=0; inode[0].body->tag, info[i].m); if (joint[i]->node[1].body) { MultiplyAdd0_p81 (rhs+ofs[i],JJ + 8*info[i].m, tmp1 + 8*joint[i]->node[1].body->tag, info[i].m); } } // complete rhs for (i=0; inode[0].body; dxBody* b2 = joint[i]->node[1].body; dJointFeedback *fb = joint[i]->feedback; if (fb) { // the user has requested feedback on the amount of force that this // joint is applying to the bodies. we use a slightly slower // computation that splits out the force components and puts them // in the feedback structure. dReal data1[8],data2[8]; Multiply1_8q1 (data1, JJ, lambda+ofs[i], info[i].m); dReal *cf1 = cforce + 8*b1->tag; cf1[0] += (fb->f1[0] = data1[0]); cf1[1] += (fb->f1[1] = data1[1]); cf1[2] += (fb->f1[2] = data1[2]); cf1[4] += (fb->t1[0] = data1[4]); cf1[5] += (fb->t1[1] = data1[5]); cf1[6] += (fb->t1[2] = data1[6]); if (b2){ Multiply1_8q1 (data2, JJ + 8*info[i].m, lambda+ofs[i], info[i].m); dReal *cf2 = cforce + 8*b2->tag; cf2[0] += (fb->f2[0] = data2[0]); cf2[1] += (fb->f2[1] = data2[1]); cf2[2] += (fb->f2[2] = data2[2]); cf2[4] += (fb->t2[0] = data2[4]); cf2[5] += (fb->t2[1] = data2[5]); cf2[6] += (fb->t2[2] = data2[6]); } } else { // no feedback is required, let's compute cforce the faster way MultiplyAdd1_8q1 (cforce + 8*b1->tag,JJ, lambda+ofs[i], info[i].m); if (b2) { MultiplyAdd1_8q1 (cforce + 8*b2->tag, JJ + 8*info[i].m, lambda+ofs[i], info[i].m); } } } UNALLOCA(c); UNALLOCA(cfm); UNALLOCA(lo); UNALLOCA(hi); UNALLOCA(findex); UNALLOCA(J); UNALLOCA(JinvM); UNALLOCA(A); UNALLOCA(tmp1); UNALLOCA(rhs); UNALLOCA(lambda); UNALLOCA(residual); } // compute the velocity update #ifdef TIMING dTimerNow ("compute velocity update"); #endif // add fe to cforce for (i=0; ifacc[j]; for (j=0; j<3; j++) cforce[i*8+4+j] += body[i]->tacc[j]; } // multiply cforce by stepsize for (i=0; i < nb*8; i++) cforce[i] *= stepsize; // add invM * cforce to the body velocity for (i=0; iinvMass; dReal *body_invI = invI + i*12; for (j=0; j<3; j++) body[i]->lvel[j] += body_invMass * cforce[i*8+j]; dMULTIPLYADD0_331 (body[i]->avel,body_invI,cforce+i*8+4); } // update the position and orientation from the new linear/angular velocity // (over the given timestep) # ifdef TIMING dTimerNow ("update position"); # endif for (i=0; ilvel[j]; for (j=0; j<3; j++) tmp_vnew[i*6+3+j] = body[i]->avel[j]; } comparator.nextMatrix (tmp_vnew,nb*6,1,0,"vnew"); UNALLOCA(tmp); #endif #ifdef TIMING dTimerNow ("tidy up"); #endif // zero all force accumulators for (i=0; ifacc[0] = 0; body[i]->facc[1] = 0; body[i]->facc[2] = 0; body[i]->facc[3] = 0; body[i]->tacc[0] = 0; body[i]->tacc[1] = 0; body[i]->tacc[2] = 0; body[i]->tacc[3] = 0; } #ifdef TIMING dTimerEnd(); if (m > 0) dTimerReport (stdout,1); #endif UNALLOCA(joint); UNALLOCA(I); UNALLOCA(invI); UNALLOCA(info); UNALLOCA(ofs); UNALLOCA(cforce); } //**************************************************************************** void dInternalStepIsland (dxWorld *world, dxBody * const *body, int nb, dxJoint * const *joint, int nj, dReal stepsize) { #ifdef dUSE_MALLOC_FOR_ALLOCA dMemoryFlag = d_MEMORY_OK; #endif #ifndef COMPARE_METHODS dInternalStepIsland_x2 (world,body,nb,joint,nj,stepsize); #ifdef dUSE_MALLOC_FOR_ALLOCA if (dMemoryFlag == d_MEMORY_OUT_OF_MEMORY) { REPORT_OUT_OF_MEMORY; return; } #endif #endif #ifdef COMPARE_METHODS int i; // save body state ALLOCA(dxBody,state,nb*sizeof(dxBody)); #ifdef dUSE_MALLOC_FOR_ALLOCA if (state == NULL) { dMemoryFlag = d_MEMORY_OUT_OF_MEMORY; REPORT_OUT_OF_MEMORY; return; } #endif for (i=0; i