--- trunk/OOPSE/libmdtools/Integrator.cpp 2003/07/14 21:28:54 597 +++ trunk/OOPSE/libmdtools/Integrator.cpp 2003/08/11 19:41:36 677 @@ -11,7 +11,7 @@ Integrator::Integrator( SimInfo *theInfo, ForceFields* #include "simError.h" -Integrator::Integrator( SimInfo *theInfo, ForceFields* the_ff ){ +template Integrator::Integrator( SimInfo *theInfo, ForceFields* the_ff ) { info = theInfo; myFF = the_ff; @@ -27,8 +27,6 @@ Integrator::Integrator( SimInfo *theInfo, ForceFields* nAtoms = info->n_atoms; - std::cerr << "integ nAtoms = " << nAtoms << "\n"; - // check for constraints constrainedA = NULL; @@ -43,7 +41,7 @@ Integrator::~Integrator() { checkConstraints(); } -Integrator::~Integrator() { +template Integrator::~Integrator() { if( nConstrained ){ delete[] constrainedA; @@ -56,7 +54,7 @@ void Integrator::checkConstraints( void ){ } -void Integrator::checkConstraints( void ){ +template void Integrator::checkConstraints( void ){ isConstrained = 0; @@ -75,12 +73,7 @@ void Integrator::checkConstraints( void ){ constrained = theArray[j]->is_constrained(); - std::cerr << "Is the folowing bond constrained \n"; - theArray[j]->printMe(); - if(constrained){ - - std::cerr << "Yes\n"; dummy_plug = theArray[j]->get_constraint(); temp_con[nConstrained].set_a( dummy_plug->get_a() ); @@ -90,7 +83,6 @@ void Integrator::checkConstraints( void ){ nConstrained++; constrained = 0; } - else std::cerr << "No.\n"; } theArray = (SRI**) molecules[i].getMyBends(); @@ -163,7 +155,7 @@ void Integrator::integrate( void ){ } -void Integrator::integrate( void ){ +template void Integrator::integrate( void ){ int i, j; // loop counters @@ -175,13 +167,10 @@ void Integrator::integrate( void ){ double currSample; double currThermal; double currStatus; - double currTime; int calcPot, calcStress; int isError; - - tStats = new Thermo( info ); statOut = new StatWriter( info ); dumpOut = new DumpWriter( info ); @@ -194,23 +183,21 @@ void Integrator::integrate( void ){ // initialize the forces before the first step - myFF->doForces(1,1); + calcForce(1, 1); if( info->setTemp ){ - tStats->velocitize(); + thermalize(); } - dumpOut->writeDump( 0.0 ); - statOut->writeStat( 0.0 ); - calcPot = 0; calcStress = 0; currSample = sampleTime; currThermal = thermalTime; currStatus = statusTime; - currTime = 0.0;; + dumpOut->writeDump( info->getTime() ); + statOut->writeStat( info->getTime() ); readyCheck(); @@ -220,38 +207,31 @@ void Integrator::integrate( void ){ MPIcheckPoint(); #endif // is_mpi + while( info->getTime() < runTime ){ - pos = Atom::getPosArray(); - vel = Atom::getVelArray(); - frc = Atom::getFrcArray(); - - while( currTime < runTime ){ - - if( (currTime+dt) >= currStatus ){ + if( (info->getTime()+dt) >= currStatus ){ calcPot = 1; calcStress = 1; } - std::cerr << currTime << "\n"; - integrateStep( calcPot, calcStress ); - currTime += dt; + info->incrTime(dt); if( info->setTemp ){ - if( currTime >= currThermal ){ - tStats->velocitize(); + if( info->getTime() >= currThermal ){ + thermalize(); currThermal += thermalTime; } } - if( currTime >= currSample ){ - dumpOut->writeDump( currTime ); + if( info->getTime() >= currSample ){ + dumpOut->writeDump( info->getTime() ); currSample += sampleTime; } - if( currTime >= currStatus ){ - statOut->writeStat( currTime ); + if( info->getTime() >= currStatus ){ + statOut->writeStat( info->getTime() ); calcPot = 0; calcStress = 0; currStatus += statusTime; @@ -265,13 +245,13 @@ void Integrator::integrate( void ){ } - dumpOut->writeFinal(currTime); + dumpOut->writeFinal(info->getTime()); delete dumpOut; delete statOut; } -void Integrator::integrateStep( int calcPot, int calcStress ){ +template void Integrator::integrateStep( int calcPot, int calcStress ){ @@ -279,173 +259,184 @@ void Integrator::integrateStep( int calcPot, int calcS preMove(); moveA(); - //if( nConstrained ) constrainA(); + if( nConstrained ) constrainA(); + + +#ifdef IS_MPI + strcpy( checkPointMsg, "Succesful moveA\n" ); + MPIcheckPoint(); +#endif // is_mpi + // calc forces - myFF->doForces(calcPot,calcStress); + calcForce(calcPot,calcStress); +#ifdef IS_MPI + strcpy( checkPointMsg, "Succesful doForces\n" ); + MPIcheckPoint(); +#endif // is_mpi + + // finish the velocity half step moveB(); if( nConstrained ) constrainB(); - + +#ifdef IS_MPI + strcpy( checkPointMsg, "Succesful moveB\n" ); + MPIcheckPoint(); +#endif // is_mpi + + } -void Integrator::moveA( void ){ +template void Integrator::moveA( void ){ - int i,j,k; - int atomIndex, aMatIndex; + int i, j; DirectionalAtom* dAtom; - double Tb[3]; - double ji[3]; + double Tb[3], ji[3]; + double A[3][3], I[3][3]; double angle; - double A[3][3], At[3][3]; + double vel[3], pos[3], frc[3]; + double mass; - for( i=0; igetMass() ) * eConvert; + atoms[i]->getVel( vel ); + atoms[i]->getPos( pos ); + atoms[i]->getFrc( frc ); + mass = atoms[i]->getMass(); - // position whole step - for( j=atomIndex; j<(atomIndex+3); j++ ) pos[j] += dt * vel[j]; - + for (j=0; j < 3; j++) { + // velocity half step + vel[j] += ( dt2 * frc[j] / mass ) * eConvert; + // position whole step + pos[j] += dt * vel[j]; + } + atoms[i]->setVel( vel ); + atoms[i]->setPos( pos ); + if( atoms[i]->isDirectional() ){ dAtom = (DirectionalAtom *)atoms[i]; // get and convert the torque to body frame - Tb[0] = dAtom->getTx(); - Tb[1] = dAtom->getTy(); - Tb[2] = dAtom->getTz(); - + dAtom->getTrq( Tb ); dAtom->lab2Body( Tb ); // get the angular momentum, and propagate a half step + + dAtom->getJ( ji ); + + for (j=0; j < 3; j++) + ji[j] += (dt2 * Tb[j]) * eConvert; - ji[0] = dAtom->getJx() + ( dt2 * Tb[0] ) * eConvert; - ji[1] = dAtom->getJy() + ( dt2 * Tb[1] ) * eConvert; - ji[2] = dAtom->getJz() + ( dt2 * Tb[2] ) * eConvert; - // use the angular velocities to propagate the rotation matrix a // full time step - + + dAtom->getA(A); + dAtom->getI(I); + // rotate about the x-axis - angle = dt2 * ji[0] / dAtom->getIxx(); - this->rotate( 1, 2, angle, ji, &Amat[aMatIndex] ); + angle = dt2 * ji[0] / I[0][0]; + this->rotate( 1, 2, angle, ji, A ); // rotate about the y-axis - angle = dt2 * ji[1] / dAtom->getIyy(); - this->rotate( 2, 0, angle, ji, &Amat[aMatIndex] ); + angle = dt2 * ji[1] / I[1][1]; + this->rotate( 2, 0, angle, ji, A ); // rotate about the z-axis - angle = dt * ji[2] / dAtom->getIzz(); - this->rotate( 0, 1, angle, ji, &Amat[aMatIndex] ); + angle = dt * ji[2] / I[2][2]; + this->rotate( 0, 1, angle, ji, A); // rotate about the y-axis - angle = dt2 * ji[1] / dAtom->getIyy(); - this->rotate( 2, 0, angle, ji, &Amat[aMatIndex] ); + angle = dt2 * ji[1] / I[1][1]; + this->rotate( 2, 0, angle, ji, A ); // rotate about the x-axis - angle = dt2 * ji[0] / dAtom->getIxx(); - this->rotate( 1, 2, angle, ji, &Amat[aMatIndex] ); + angle = dt2 * ji[0] / I[0][0]; + this->rotate( 1, 2, angle, ji, A ); - dAtom->setJx( ji[0] ); - dAtom->setJy( ji[1] ); - dAtom->setJz( ji[2] ); - std::cerr << "Amat[" << i << "]\n"; - info->printMat9( &Amat[aMatIndex] ); + dAtom->setJ( ji ); + dAtom->setA( A ); - std::cerr << "ji[" << i << "]\t" - << ji[0] << "\t" - << ji[1] << "\t" - << ji[2] << "\n"; - - } - + } } } -void Integrator::moveB( void ){ - int i,j,k; - int atomIndex, aMatIndex; +template void Integrator::moveB( void ){ + int i, j; DirectionalAtom* dAtom; - double Tb[3]; - double ji[3]; + double Tb[3], ji[3]; + double vel[3], frc[3]; + double mass; for( i=0; igetVel( vel ); + atoms[i]->getFrc( frc ); - // velocity half step - for( j=atomIndex; j<(atomIndex+3); j++ ) - vel[j] += ( dt2 * frc[j] / atoms[i]->getMass() ) * eConvert; + mass = atoms[i]->getMass(); + // velocity half step + for (j=0; j < 3; j++) + vel[j] += ( dt2 * frc[j] / mass ) * eConvert; + + atoms[i]->setVel( vel ); if( atoms[i]->isDirectional() ){ - + dAtom = (DirectionalAtom *)atoms[i]; - - // get and convert the torque to body frame - - Tb[0] = dAtom->getTx(); - Tb[1] = dAtom->getTy(); - Tb[2] = dAtom->getTz(); - - std::cerr << "TrqB[" << i << "]\t" - << Tb[0] << "\t" - << Tb[1] << "\t" - << Tb[2] << "\n"; + + // get and convert the torque to body frame + dAtom->getTrq( Tb ); dAtom->lab2Body( Tb ); - - // get the angular momentum, and complete the angular momentum - // half step - - ji[0] = dAtom->getJx() + ( dt2 * Tb[0] ) * eConvert; - ji[1] = dAtom->getJy() + ( dt2 * Tb[1] ) * eConvert; - ji[2] = dAtom->getJz() + ( dt2 * Tb[2] ) * eConvert; - - dAtom->setJx( ji[0] ); - dAtom->setJy( ji[1] ); - dAtom->setJz( ji[2] ); + // get the angular momentum, and propagate a half step - std::cerr << "Amat[" << i << "]\n"; - info->printMat9( &Amat[aMatIndex] ); - - std::cerr << "ji[" << i << "]\t" - << ji[0] << "\t" - << ji[1] << "\t" - << ji[2] << "\n"; + dAtom->getJ( ji ); + + for (j=0; j < 3; j++) + ji[j] += (dt2 * Tb[j]) * eConvert; + + + dAtom->setJ( ji ); } } - } -void Integrator::preMove( void ){ - int i; +template void Integrator::preMove( void ){ + int i, j; + double pos[3]; if( nConstrained ){ - for(i=0; i<(nAtoms*3); i++) oldPos[i] = pos[i]; - } -} + for(i=0; i < nAtoms; i++) { + + atoms[i]->getPos( pos ); -void Integrator::constrainA(){ + for (j = 0; j < 3; j++) { + oldPos[3*i + j] = pos[j]; + } + } + } +} + +template void Integrator::constrainA(){ + int i,j,k; int done; + double posA[3], posB[3]; + double velA[3], velB[3]; double pab[3]; double rab[3]; int a, b, ax, ay, az, bx, by, bz; @@ -457,8 +448,7 @@ void Integrator::constrainA(){ double gab; int iteration; - for( i=0; igetPos( posA ); + atoms[b]->getPos( posB ); + + for (j = 0; j < 3; j++ ) + pab[j] = posA[j] - posB[j]; + //periodic boundary condition info->wrapVector( pab ); @@ -531,26 +523,38 @@ void Integrator::constrainA(){ dy = rab[1] * gab; dz = rab[2] * gab; - pos[ax] += rma * dx; - pos[ay] += rma * dy; - pos[az] += rma * dz; + posA[0] += rma * dx; + posA[1] += rma * dy; + posA[2] += rma * dz; - pos[bx] -= rmb * dx; - pos[by] -= rmb * dy; - pos[bz] -= rmb * dz; + atoms[a]->setPos( posA ); + posB[0] -= rmb * dx; + posB[1] -= rmb * dy; + posB[2] -= rmb * dz; + + atoms[b]->setPos( posB ); + dx = dx / dt; dy = dy / dt; dz = dz / dt; - vel[ax] += rma * dx; - vel[ay] += rma * dy; - vel[az] += rma * dz; + atoms[a]->getVel( velA ); - vel[bx] -= rmb * dx; - vel[by] -= rmb * dy; - vel[bz] -= rmb * dz; + velA[0] += rma * dx; + velA[1] += rma * dy; + velA[2] += rma * dz; + atoms[a]->setVel( velA ); + + atoms[b]->getVel( velB ); + + velB[0] -= rmb * dx; + velB[1] -= rmb * dy; + velB[2] -= rmb * dz; + + atoms[b]->setVel( velB ); + moving[a] = 1; moving[b] = 1; done = 0; @@ -578,10 +582,12 @@ void Integrator::constrainB( void ){ } -void Integrator::constrainB( void ){ +template void Integrator::constrainB( void ){ int i,j,k; int done; + double posA[3], posB[3]; + double velA[3], velB[3]; double vxab, vyab, vzab; double rab[3]; int a, b, ax, ay, az, bx, by, bz; @@ -617,15 +623,20 @@ void Integrator::constrainB( void ){ bz = (b*3) + 2; if( moved[a] || moved[b] ){ - - vxab = vel[ax] - vel[bx]; - vyab = vel[ay] - vel[by]; - vzab = vel[az] - vel[bz]; - rab[0] = pos[ax] - pos[bx]; - rab[1] = pos[ay] - pos[by]; - rab[2] = pos[az] - pos[bz]; - + atoms[a]->getVel( velA ); + atoms[b]->getVel( velB ); + + vxab = velA[0] - velB[0]; + vyab = velA[1] - velB[1]; + vzab = velA[2] - velB[2]; + + atoms[a]->getPos( posA ); + atoms[b]->getPos( posB ); + + for (j = 0; j < 3; j++) + rab[j] = posA[j] - posB[j]; + info->wrapVector( rab ); rma = 1.0 / atoms[a]->getMass(); @@ -640,14 +651,18 @@ void Integrator::constrainB( void ){ dx = rab[0] * gab; dy = rab[1] * gab; dz = rab[2] * gab; - - vel[ax] += rma * dx; - vel[ay] += rma * dy; - vel[az] += rma * dz; + + velA[0] += rma * dx; + velA[1] += rma * dy; + velA[2] += rma * dz; - vel[bx] -= rmb * dx; - vel[by] -= rmb * dy; - vel[bz] -= rmb * dz; + atoms[a]->setVel( velA ); + + velB[0] -= rmb * dx; + velB[1] -= rmb * dy; + velB[2] -= rmb * dz; + + atoms[b]->setVel( velB ); moving[a] = 1; moving[b] = 1; @@ -663,7 +678,7 @@ void Integrator::constrainB( void ){ iteration++; } - + if( !done ){ @@ -676,15 +691,9 @@ void Integrator::constrainB( void ){ } +template void Integrator::rotate( int axes1, int axes2, double angle, double ji[3], + double A[3][3] ){ - - - - - -void Integrator::rotate( int axes1, int axes2, double angle, double ji[3], - double A[9] ){ - int i,j,k; double sinAngle; double cosAngle; @@ -695,12 +704,11 @@ void Integrator::rotate( int axes1, int axes2, double double tempA[3][3]; double tempJ[3]; - // initialize the tempA for(i=0; i<3; i++){ for(j=0; j<3; j++){ - tempA[j][i] = A[3*i+j]; + tempA[j][i] = A[i][j]; } } @@ -757,10 +765,19 @@ void Integrator::rotate( int axes1, int axes2, double for(i=0; i<3; i++){ for(j=0; j<3; j++){ - A[3*j+i] = 0.0; + A[j][i] = 0.0; for(k=0; k<3; k++){ - A[3*j+i] += tempA[i][k] * rot[j][k]; + A[j][i] += tempA[i][k] * rot[j][k]; } } } } + +template void Integrator::calcForce( int calcPot, int calcStress ){ + myFF->doForces(calcPot,calcStress); + +} + +template void Integrator::thermalize(){ + tStats->velocitize(); +}