ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/Integrator.cpp
(Generate patch)

Comparing trunk/OOPSE/libmdtools/Integrator.cpp (file contents):
Revision 572 by mmeineke, Wed Jul 2 21:26:55 2003 UTC vs.
Revision 677 by tim, Mon Aug 11 19:41:36 2003 UTC

# Line 11 | Line 11 | Integrator::Integrator( SimInfo *theInfo, ForceFields*
11   #include "simError.h"
12  
13  
14 < Integrator::Integrator( SimInfo *theInfo, ForceFields* the_ff ){
14 > template<typename T> Integrator<T>::Integrator( SimInfo *theInfo, ForceFields* the_ff ) {
15    
16    info = theInfo;
17    myFF = the_ff;
# Line 41 | Line 41 | Integrator::~Integrator() {
41    checkConstraints();
42   }
43  
44 < Integrator::~Integrator() {
44 > template<typename T> Integrator<T>::~Integrator() {
45    
46    if( nConstrained ){
47      delete[] constrainedA;
# Line 54 | Line 54 | void Integrator::checkConstraints( void ){
54    
55   }
56  
57 < void Integrator::checkConstraints( void ){
57 > template<typename T> void Integrator<T>::checkConstraints( void ){
58  
59  
60    isConstrained = 0;
# Line 72 | Line 72 | void Integrator::checkConstraints( void ){
72      for(int j=0; j<molecules[i].getNBonds(); j++){
73        
74        constrained = theArray[j]->is_constrained();
75 <      
75 >
76        if(constrained){
77 <        
77 >
78          dummy_plug = theArray[j]->get_constraint();
79          temp_con[nConstrained].set_a( dummy_plug->get_a() );
80          temp_con[nConstrained].set_b( dummy_plug->get_b() );
# Line 82 | Line 82 | void Integrator::checkConstraints( void ){
82          
83          nConstrained++;
84          constrained = 0;
85 <      }
85 >      }
86      }
87  
88      theArray = (SRI**) molecules[i].getMyBends();
# Line 155 | Line 155 | void Integrator::integrate( void ){
155   }
156  
157  
158 < void Integrator::integrate( void ){
158 > template<typename T> void Integrator<T>::integrate( void ){
159  
160    int i, j;                         // loop counters
161  
# Line 167 | Line 167 | void Integrator::integrate( void ){
167    double currSample;
168    double currThermal;
169    double currStatus;
170  double currTime;
170  
171    int calcPot, calcStress;
172    int isError;
173  
175
176
174    tStats   = new Thermo( info );
175    statOut  = new StatWriter( info );
176    dumpOut  = new DumpWriter( info );
# Line 186 | Line 183 | void Integrator::integrate( void ){
183  
184    // initialize the forces before the first step
185  
186 <  myFF->doForces(1,1);
186 >  calcForce(1, 1);
187    
188    if( info->setTemp ){
189      
190 <    tStats->velocitize();
190 >    thermalize();
191    }
192    
196  dumpOut->writeDump( 0.0 );
197  statOut->writeStat( 0.0 );
198  
193    calcPot     = 0;
194    calcStress  = 0;
195    currSample  = sampleTime;
196    currThermal = thermalTime;
197    currStatus  = statusTime;
204  currTime    = 0.0;;
198  
199 +  dumpOut->writeDump( info->getTime() );
200 +  statOut->writeStat( info->getTime() );
201  
202    readyCheck();
203  
# Line 212 | Line 207 | void Integrator::integrate( void ){
207    MPIcheckPoint();
208   #endif // is_mpi
209  
210 +  while( info->getTime() < runTime ){
211  
212 <  pos  = Atom::getPosArray();
217 <  vel  = Atom::getVelArray();
218 <  frc  = Atom::getFrcArray();
219 <  trq  = Atom::getTrqArray();
220 <  Amat = Atom::getAmatArray();
221 <
222 <  while( currTime < runTime ){
223 <
224 <    if( (currTime+dt) >= currStatus ){
212 >    if( (info->getTime()+dt) >= currStatus ){
213        calcPot = 1;
214        calcStress = 1;
215      }
216  
217      integrateStep( calcPot, calcStress );
218        
219 <    currTime += dt;
219 >    info->incrTime(dt);
220  
221      if( info->setTemp ){
222 <      if( currTime >= currThermal ){
223 <        tStats->velocitize();
222 >      if( info->getTime() >= currThermal ){
223 >        thermalize();
224          currThermal += thermalTime;
225        }
226      }
227  
228 <    if( currTime >= currSample ){
229 <      dumpOut->writeDump( currTime );
228 >    if( info->getTime() >= currSample ){
229 >      dumpOut->writeDump( info->getTime() );
230        currSample += sampleTime;
231      }
232  
233 <    if( currTime >= currStatus ){
234 <      statOut->writeStat( currTime );
233 >    if( info->getTime() >= currStatus ){
234 >      statOut->writeStat( info->getTime() );
235        calcPot = 0;
236        calcStress = 0;
237        currStatus += statusTime;
# Line 257 | Line 245 | void Integrator::integrate( void ){
245  
246    }
247  
248 <  dumpOut->writeFinal(currTime);
248 >  dumpOut->writeFinal(info->getTime());
249  
250    delete dumpOut;
251    delete statOut;
252   }
253  
254 < void Integrator::integrateStep( int calcPot, int calcStress ){
254 > template<typename T> void Integrator<T>::integrateStep( int calcPot, int calcStress ){
255  
256  
257        
# Line 272 | Line 260 | void Integrator::integrateStep( int calcPot, int calcS
260    preMove();
261    moveA();
262    if( nConstrained ) constrainA();
263 +
264 +  
265 + #ifdef IS_MPI
266 +  strcpy( checkPointMsg, "Succesful moveA\n" );
267 +  MPIcheckPoint();
268 + #endif // is_mpi
269 +  
270  
271    // calc forces
272  
273 <  myFF->doForces(calcPot,calcStress);
273 >  calcForce(calcPot,calcStress);
274  
275 + #ifdef IS_MPI
276 +  strcpy( checkPointMsg, "Succesful doForces\n" );
277 +  MPIcheckPoint();
278 + #endif // is_mpi
279 +  
280 +
281    // finish the velocity  half step
282    
283    moveB();
284    if( nConstrained ) constrainB();
285 <  
285 >  
286 > #ifdef IS_MPI
287 >  strcpy( checkPointMsg, "Succesful moveB\n" );
288 >  MPIcheckPoint();
289 > #endif // is_mpi
290 >  
291 >
292   }
293  
294  
295 < void Integrator::moveA( void ){
295 > template<typename T> void Integrator<T>::moveA( void ){
296    
297 <  int i,j,k;
291 <  int atomIndex, aMatIndex;
297 >  int i, j;
298    DirectionalAtom* dAtom;
299 <  double Tb[3];
300 <  double ji[3];
299 >  double Tb[3], ji[3];
300 >  double A[3][3], I[3][3];
301    double angle;
302 +  double vel[3], pos[3], frc[3];
303 +  double mass;
304  
305 +  for( i=0; i<nAtoms; i++ ){
306  
307 +    atoms[i]->getVel( vel );
308 +    atoms[i]->getPos( pos );
309 +    atoms[i]->getFrc( frc );
310  
311 <  for( i=0; i<nAtoms; i++ ){
300 <    atomIndex = i * 3;
301 <    aMatIndex = i * 9;
311 >    mass = atoms[i]->getMass();
312  
313 <    // velocity half step
314 <    for( j=atomIndex; j<(atomIndex+3); j++ )
315 <      vel[j] += ( dt2 * frc[j] / atoms[i]->getMass() ) * eConvert;
313 >    for (j=0; j < 3; j++) {
314 >      // velocity half step
315 >      vel[j] += ( dt2 * frc[j] / mass ) * eConvert;
316 >      // position whole step
317 >      pos[j] += dt * vel[j];
318 >    }
319  
320 <    // position whole step    
321 <    for( j=atomIndex; j<(atomIndex+3); j++ ) pos[j] += dt * vel[j];
322 <    
320 >    atoms[i]->setVel( vel );
321 >    atoms[i]->setPos( pos );
322 >
323      if( atoms[i]->isDirectional() ){
324  
325        dAtom = (DirectionalAtom *)atoms[i];
326            
327        // get and convert the torque to body frame
328        
329 <      Tb[0] = dAtom->getTx();
317 <      Tb[1] = dAtom->getTy();
318 <      Tb[2] = dAtom->getTz();
319 <      
329 >      dAtom->getTrq( Tb );
330        dAtom->lab2Body( Tb );
331 <      
331 >
332        // get the angular momentum, and propagate a half step
333 +
334 +      dAtom->getJ( ji );
335 +
336 +      for (j=0; j < 3; j++)
337 +        ji[j] += (dt2 * Tb[j]) * eConvert;
338        
324      ji[0] = dAtom->getJx() + ( dt2 * Tb[0] ) * eConvert;
325      ji[1] = dAtom->getJy() + ( dt2 * Tb[1] ) * eConvert;
326      ji[2] = dAtom->getJz() + ( dt2 * Tb[2] ) * eConvert;
327      
339        // use the angular velocities to propagate the rotation matrix a
340        // full time step
341 <      
341 >
342 >      dAtom->getA(A);
343 >      dAtom->getI(I);
344 >    
345        // rotate about the x-axis      
346 <      angle = dt2 * ji[0] / dAtom->getIxx();
347 <      this->rotate( 1, 2, angle, ji, &Amat[aMatIndex] );
348 <      
346 >      angle = dt2 * ji[0] / I[0][0];
347 >      this->rotate( 1, 2, angle, ji, A );
348 >
349        // rotate about the y-axis
350 <      angle = dt2 * ji[1] / dAtom->getIyy();
351 <      this->rotate( 2, 0, angle, ji, &Amat[aMatIndex] );
350 >      angle = dt2 * ji[1] / I[1][1];
351 >      this->rotate( 2, 0, angle, ji, A );
352        
353        // rotate about the z-axis
354 <      angle = dt * ji[2] / dAtom->getIzz();
355 <      this->rotate( 0, 1, angle, ji, &Amat[aMatIndex] );
354 >      angle = dt * ji[2] / I[2][2];
355 >      this->rotate( 0, 1, angle, ji, A);
356        
357        // rotate about the y-axis
358 <      angle = dt2 * ji[1] / dAtom->getIyy();
359 <      this->rotate( 2, 0, angle, ji, &Amat[aMatIndex] );
358 >      angle = dt2 * ji[1] / I[1][1];
359 >      this->rotate( 2, 0, angle, ji, A );
360        
361         // rotate about the x-axis
362 <      angle = dt2 * ji[0] / dAtom->getIxx();
363 <      this->rotate( 1, 2, angle, ji, &Amat[aMatIndex] );
362 >      angle = dt2 * ji[0] / I[0][0];
363 >      this->rotate( 1, 2, angle, ji, A );
364        
365 <      dAtom->setJx( ji[0] );
366 <      dAtom->setJy( ji[1] );
367 <      dAtom->setJz( ji[2] );
368 <    }
369 <    
365 >
366 >      dAtom->setJ( ji );
367 >      dAtom->setA( A  );
368 >          
369 >    }    
370    }
371   }
372  
373  
374 < void Integrator::moveB( void ){
375 <  int i,j,k;
362 <  int atomIndex;
374 > template<typename T> void Integrator<T>::moveB( void ){
375 >  int i, j;
376    DirectionalAtom* dAtom;
377 <  double Tb[3];
378 <  double ji[3];
377 >  double Tb[3], ji[3];
378 >  double vel[3], frc[3];
379 >  double mass;
380  
381    for( i=0; i<nAtoms; i++ ){
382 <    atomIndex = i * 3;
382 >
383 >    atoms[i]->getVel( vel );
384 >    atoms[i]->getFrc( frc );
385  
386 <    // velocity half step
371 <    for( j=atomIndex; j<(atomIndex+3); j++ )
372 <      vel[j] += ( dt2 * frc[j] / atoms[i]->getMass() ) * eConvert;
386 >    mass = atoms[i]->getMass();
387  
388 +    // velocity half step
389 +    for (j=0; j < 3; j++)
390 +      vel[j] += ( dt2 * frc[j] / mass ) * eConvert;
391 +    
392 +    atoms[i]->setVel( vel );
393 +
394      if( atoms[i]->isDirectional() ){
395 <      
395 >
396        dAtom = (DirectionalAtom *)atoms[i];
397 <      
398 <      // get and convert the torque to body frame
399 <      
400 <      Tb[0] = dAtom->getTx();
381 <      Tb[1] = dAtom->getTy();
382 <      Tb[2] = dAtom->getTz();
383 <      
397 >
398 >      // get and convert the torque to body frame      
399 >
400 >      dAtom->getTrq( Tb );
401        dAtom->lab2Body( Tb );
402 +
403 +      // get the angular momentum, and propagate a half step
404 +
405 +      dAtom->getJ( ji );
406 +
407 +      for (j=0; j < 3; j++)
408 +        ji[j] += (dt2 * Tb[j]) * eConvert;
409        
410 <      // get the angular momentum, and complete the angular momentum
411 <      // half step
388 <      
389 <      ji[0] = dAtom->getJx() + ( dt2 * Tb[0] ) * eConvert;
390 <      ji[1] = dAtom->getJy() + ( dt2 * Tb[1] ) * eConvert;
391 <      ji[2] = dAtom->getJz() + ( dt2 * Tb[2] ) * eConvert;
392 <      
393 <      dAtom->setJx( ji[0] );
394 <      dAtom->setJy( ji[1] );
395 <      dAtom->setJz( ji[2] );
410 >
411 >      dAtom->setJ( ji );
412      }
413    }
398
414   }
415  
416 < void Integrator::preMove( void ){
417 <  int i;
416 > template<typename T> void Integrator<T>::preMove( void ){
417 >  int i, j;
418 >  double pos[3];
419  
420    if( nConstrained ){
421  
422 <    for(i=0; i<(nAtoms*3); i++) oldPos[i] = pos[i];
423 <  }
424 < }  
422 >    for(i=0; i < nAtoms; i++) {
423 >
424 >      atoms[i]->getPos( pos );
425  
426 < void Integrator::constrainA(){
426 >      for (j = 0; j < 3; j++) {        
427 >        oldPos[3*i + j] = pos[j];
428 >      }
429  
430 +    }
431 +  }  
432 + }
433 +
434 + template<typename T> void Integrator<T>::constrainA(){
435 +
436    int i,j,k;
437    int done;
438 +  double posA[3], posB[3];
439 +  double velA[3], velB[3];
440    double pab[3];
441    double rab[3];
442    int a, b, ax, ay, az, bx, by, bz;
# Line 422 | Line 448 | void Integrator::constrainA(){
448    double gab;
449    int iteration;
450  
451 <
426 <  
427 <  for( i=0; i<nAtoms; i++){
428 <    
451 >  for( i=0; i<nAtoms; i++){    
452      moving[i] = 0;
453      moved[i]  = 1;
454    }
# Line 449 | Line 472 | void Integrator::constrainA(){
472        bz = (b*3) + 2;
473  
474        if( moved[a] || moved[b] ){
475 <        
476 <        pab[0] = pos[ax] - pos[bx];
477 <        pab[1] = pos[ay] - pos[by];
478 <        pab[2] = pos[az] - pos[bz];
479 <
475 >        
476 >        atoms[a]->getPos( posA );
477 >        atoms[b]->getPos( posB );
478 >        
479 >        for (j = 0; j < 3; j++ )
480 >          pab[j] = posA[j] - posB[j];
481 >        
482          //periodic boundary condition
483  
484          info->wrapVector( pab );
# Line 498 | Line 523 | void Integrator::constrainA(){
523            dy = rab[1] * gab;
524            dz = rab[2] * gab;
525  
526 <          pos[ax] += rma * dx;
527 <          pos[ay] += rma * dy;
528 <          pos[az] += rma * dz;
526 >          posA[0] += rma * dx;
527 >          posA[1] += rma * dy;
528 >          posA[2] += rma * dz;
529  
530 <          pos[bx] -= rmb * dx;
506 <          pos[by] -= rmb * dy;
507 <          pos[bz] -= rmb * dz;
530 >          atoms[a]->setPos( posA );
531  
532 +          posB[0] -= rmb * dx;
533 +          posB[1] -= rmb * dy;
534 +          posB[2] -= rmb * dz;
535 +
536 +          atoms[b]->setPos( posB );
537 +
538            dx = dx / dt;
539            dy = dy / dt;
540            dz = dz / dt;
541  
542 <          vel[ax] += rma * dx;
514 <          vel[ay] += rma * dy;
515 <          vel[az] += rma * dz;
542 >          atoms[a]->getVel( velA );
543  
544 <          vel[bx] -= rmb * dx;
545 <          vel[by] -= rmb * dy;
546 <          vel[bz] -= rmb * dz;
544 >          velA[0] += rma * dx;
545 >          velA[1] += rma * dy;
546 >          velA[2] += rma * dz;
547  
548 +          atoms[a]->setVel( velA );
549 +
550 +          atoms[b]->getVel( velB );
551 +
552 +          velB[0] -= rmb * dx;
553 +          velB[1] -= rmb * dy;
554 +          velB[2] -= rmb * dz;
555 +
556 +          atoms[b]->setVel( velB );
557 +
558            moving[a] = 1;
559            moving[b] = 1;
560            done = 0;
# Line 545 | Line 582 | void Integrator::constrainB( void ){
582  
583   }
584  
585 < void Integrator::constrainB( void ){
585 > template<typename T> void Integrator<T>::constrainB( void ){
586    
587    int i,j,k;
588    int done;
589 +  double posA[3], posB[3];
590 +  double velA[3], velB[3];
591    double vxab, vyab, vzab;
592    double rab[3];
593    int a, b, ax, ay, az, bx, by, bz;
# Line 584 | Line 623 | void Integrator::constrainB( void ){
623        bz = (b*3) + 2;
624  
625        if( moved[a] || moved[b] ){
587        
588        vxab = vel[ax] - vel[bx];
589        vyab = vel[ay] - vel[by];
590        vzab = vel[az] - vel[bz];
626  
627 <        rab[0] = pos[ax] - pos[bx];
628 <        rab[1] = pos[ay] - pos[by];
629 <        rab[2] = pos[az] - pos[bz];
630 <        
627 >        atoms[a]->getVel( velA );
628 >        atoms[b]->getVel( velB );
629 >          
630 >        vxab = velA[0] - velB[0];
631 >        vyab = velA[1] - velB[1];
632 >        vzab = velA[2] - velB[2];
633 >
634 >        atoms[a]->getPos( posA );
635 >        atoms[b]->getPos( posB );
636 >
637 >        for (j = 0; j < 3; j++)
638 >          rab[j] = posA[j] - posB[j];
639 >          
640          info->wrapVector( rab );
641          
642          rma = 1.0 / atoms[a]->getMass();
# Line 607 | Line 651 | void Integrator::constrainB( void ){
651            dx = rab[0] * gab;
652            dy = rab[1] * gab;
653            dz = rab[2] * gab;
654 <          
655 <          vel[ax] += rma * dx;
656 <          vel[ay] += rma * dy;
657 <          vel[az] += rma * dz;
654 >        
655 >          velA[0] += rma * dx;
656 >          velA[1] += rma * dy;
657 >          velA[2] += rma * dz;
658  
659 <          vel[bx] -= rmb * dx;
660 <          vel[by] -= rmb * dy;
661 <          vel[bz] -= rmb * dz;
659 >          atoms[a]->setVel( velA );
660 >
661 >          velB[0] -= rmb * dx;
662 >          velB[1] -= rmb * dy;
663 >          velB[2] -= rmb * dz;
664 >
665 >          atoms[b]->setVel( velB );
666            
667            moving[a] = 1;
668            moving[b] = 1;
# Line 630 | Line 678 | void Integrator::constrainB( void ){
678      
679      iteration++;
680    }
681 <
681 >  
682    if( !done ){
683  
684    
# Line 643 | Line 691 | void Integrator::constrainB( void ){
691  
692   }
693  
694 + template<typename T> void Integrator<T>::rotate( int axes1, int axes2, double angle, double ji[3],
695 +                         double A[3][3] ){
696  
647
648
649
650
651
652 void Integrator::rotate( int axes1, int axes2, double angle, double ji[3],
653                         double A[9] ){
654
697    int i,j,k;
698    double sinAngle;
699    double cosAngle;
# Line 666 | Line 708 | void Integrator::rotate( int axes1, int axes2, double
708  
709    for(i=0; i<3; i++){
710      for(j=0; j<3; j++){
711 <      tempA[j][i] = A[3*i + j];
711 >      tempA[j][i] = A[i][j];
712      }
713    }
714  
# Line 723 | Line 765 | void Integrator::rotate( int axes1, int axes2, double
765  
766    for(i=0; i<3; i++){
767      for(j=0; j<3; j++){
768 <      A[3*j + i] = 0.0;
768 >      A[j][i] = 0.0;
769        for(k=0; k<3; k++){
770 <        A[3*j + i] += tempA[i][k] * rot[j][k];
770 >        A[j][i] += tempA[i][k] * rot[j][k];
771        }
772      }
773    }
774   }
775 +
776 + template<typename T> void Integrator<T>::calcForce( int calcPot, int calcStress ){
777 +   myFF->doForces(calcPot,calcStress);
778 +  
779 + }
780 +
781 + template<typename T> void Integrator<T>::thermalize(){
782 +  tStats->velocitize();  
783 + }

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines