ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/mdtools/md_code/Symplectic.cpp
Revision: 189
Committed: Tue Nov 26 21:04:43 2002 UTC (21 years, 7 months ago) by mmeineke
File size: 15047 byte(s)
Log Message:
*** empty log message ***

File Contents

# User Rev Content
1 mmeineke 10 #include <iostream>
2     #include <cstdlib>
3    
4     #include "Integrator.hpp"
5     #include "Thermo.hpp"
6     #include "ReadWrite.hpp"
7 mmeineke 184 #include "simError.h"
8 mmeineke 10
9     extern "C"{
10    
11     void v_constrain_a_( double &dt, int &n_atoms, double* mass,
12     double* Rx, double* Ry, double* Rz,
13     double* Vx, double* Vy, double* Vz,
14     double* Fx, double* Fy, double* Fz,
15     int &n_constrained, double *constr_sqr,
16     int* constr_i, int* constr_j,
17     double &box_x, double &box_y, double &box_z );
18    
19     void v_constrain_b_( double &dt, int &n_atoms, double* mass,
20     double* Rx, double* Ry, double* Rz,
21     double* Vx, double* Vy, double* Vz,
22     double* Fx, double* Fy, double* Fz,
23     double &Kinetic,
24     int &n_constrained, double *constr_sqr,
25     int* constr_i, int* constr_j,
26     double &box_x, double &box_y, double &box_z );
27     }
28    
29    
30    
31    
32     Symplectic::Symplectic( SimInfo* the_entry_plug ){
33     entry_plug = the_entry_plug;
34     isFirst = 1;
35    
36     srInteractions = entry_plug->sr_interactions;
37     longRange = entry_plug->longRange;
38     nSRI = entry_plug->n_SRI;
39    
40     // give a little love back to the SimInfo object
41    
42     if( entry_plug->the_integrator != NULL ) delete entry_plug->the_integrator;
43     entry_plug->the_integrator = this;
44    
45     // grab the masses
46    
47     mass = new double[entry_plug->n_atoms];
48     for(int i = 0; i < entry_plug->n_atoms; i++){
49     mass[i] = entry_plug->atoms[i]->getMass();
50     }
51    
52    
53     // check for constraints
54    
55     is_constrained = 0;
56    
57     Constraint *temp_con;
58     Constraint *dummy_plug;
59     temp_con = new Constraint[nSRI];
60     n_constrained = 0;
61     int constrained = 0;
62    
63     for(int i = 0; i < nSRI; i++){
64    
65     constrained = srInteractions[i]->is_constrained();
66    
67     if(constrained){
68    
69     dummy_plug = srInteractions[i]->get_constraint();
70     temp_con[n_constrained].set_a( dummy_plug->get_a() );
71     temp_con[n_constrained].set_b( dummy_plug->get_b() );
72     temp_con[n_constrained].set_dsqr( dummy_plug->get_dsqr() );
73    
74     n_constrained++;
75     constrained = 0;
76     }
77     }
78    
79     if(n_constrained > 0){
80    
81     is_constrained = 1;
82     constrained_i = new int[n_constrained];
83     constrained_j = new int[n_constrained];
84     constrained_dsqr = new double[n_constrained];
85    
86     for( int i = 0; i < n_constrained; i++){
87    
88     /* add 1 to the index for the fortran arrays. */
89    
90     constrained_i[i] = temp_con[i].get_a() + 1;
91     constrained_j[i] = temp_con[i].get_b() + 1;
92     constrained_dsqr[i] = temp_con[i].get_dsqr();
93     }
94     }
95    
96     delete[] temp_con;
97     }
98    
99     Symplectic::~Symplectic() {
100    
101     if( n_constrained ){
102     delete[] constrained_i;
103     delete[] constrained_j;
104     delete[] constrained_dsqr;
105     }
106    
107     }
108    
109    
110     void Symplectic::integrate( void ){
111    
112     const double e_convert = 4.184e-4; // converts kcal/mol -> amu*A^2/fs^2
113    
114     int i, j; // loop counters
115     int nAtoms = entry_plug->n_atoms; // the number of atoms
116     double kE = 0.0; // the kinetic energy
117     double rot_kE;
118     double trans_kE;
119     int tl; // the time loop conter
120     double dt2; // half the dt
121    
122     double vx, vy, vz; // the velocities
123 chuckv 134 // double vx2, vy2, vz2; // the square of the velocities
124 mmeineke 10 double rx, ry, rz; // the postitions
125    
126     double ji[3]; // the body frame angular momentum
127     double jx2, jy2, jz2; // the square of the angular momentums
128     double Tb[3]; // torque in the body frame
129     double angle; // the angle through which to rotate the rotation matrix
130     double A[3][3]; // the rotation matrix
131    
132 mmeineke 25 int time;
133 mmeineke 10
134     double dt = entry_plug->dt;
135     double runTime = entry_plug->run_time;
136     double sampleTime = entry_plug->sampleTime;
137     double statusTime = entry_plug->statusTime;
138     double thermalTime = entry_plug->thermalTime;
139    
140     int n_loops = (int)( runTime / dt );
141     int sample_n = (int)( sampleTime / dt );
142     int status_n = (int)( statusTime / dt );
143     int vel_n = (int)( thermalTime / dt );
144    
145     Thermo *tStats = new Thermo( entry_plug );
146    
147     StatWriter* e_out = new StatWriter( entry_plug );
148     DumpWriter* dump_out = new DumpWriter( entry_plug );
149    
150     Atom** atoms = entry_plug->atoms;
151     DirectionalAtom* dAtom;
152     dt2 = 0.5 * dt;
153    
154     // initialize the forces the before the first step
155    
156    
157     for(i = 0; i < nAtoms; i++){
158     atoms[i]->zeroForces();
159     }
160    
161     for(i = 0; i < nSRI; i++){
162    
163     srInteractions[i]->calc_forces();
164     }
165    
166     longRange->calc_forces();
167    
168     if( entry_plug->setTemp ){
169    
170     tStats->velocitize();
171     }
172    
173 mmeineke 25 dump_out->writeDump( 0.0 );
174     e_out->writeStat( 0.0 );
175    
176 mmeineke 10 if( n_constrained ){
177    
178     double *Rx = new double[nAtoms];
179     double *Ry = new double[nAtoms];
180     double *Rz = new double[nAtoms];
181    
182     double *Vx = new double[nAtoms];
183     double *Vy = new double[nAtoms];
184     double *Vz = new double[nAtoms];
185    
186     double *Fx = new double[nAtoms];
187     double *Fy = new double[nAtoms];
188     double *Fz = new double[nAtoms];
189    
190    
191     for( tl=0; tl < n_loops; tl++ ){
192    
193     for( j=0; j<nAtoms; j++ ){
194    
195     Rx[j] = atoms[j]->getX();
196     Ry[j] = atoms[j]->getY();
197     Rz[j] = atoms[j]->getZ();
198    
199     Vx[j] = atoms[j]->get_vx();
200     Vy[j] = atoms[j]->get_vy();
201     Vz[j] = atoms[j]->get_vz();
202    
203     Fx[j] = atoms[j]->getFx();
204     Fy[j] = atoms[j]->getFy();
205     Fz[j] = atoms[j]->getFz();
206    
207     }
208    
209     v_constrain_a_( dt, nAtoms, mass, Rx, Ry, Rz, Vx, Vy, Vz,
210     Fx, Fy, Fz,
211     n_constrained, constrained_dsqr,
212     constrained_i, constrained_j,
213     entry_plug->box_x,
214     entry_plug->box_y,
215     entry_plug->box_z );
216    
217     for( j=0; j<nAtoms; j++ ){
218    
219     atoms[j]->setX(Rx[j]);
220     atoms[j]->setY(Ry[j]);
221     atoms[j]->setZ(Rz[j]);
222    
223     atoms[j]->set_vx(Vx[j]);
224     atoms[j]->set_vy(Vy[j]);
225     atoms[j]->set_vz(Vz[j]);
226     }
227    
228    
229     for( i=0; i<nAtoms; i++ ){
230     if( atoms[i]->isDirectional() ){
231    
232     dAtom = (DirectionalAtom *)atoms[i];
233    
234     // get and convert the torque to body frame
235    
236     Tb[0] = dAtom->getTx();
237     Tb[1] = dAtom->getTy();
238     Tb[2] = dAtom->getTz();
239    
240     dAtom->lab2Body( Tb );
241    
242     // get the angular momentum, and propagate a half step
243    
244     ji[0] = dAtom->getJx() + ( dt2 * Tb[0] ) * e_convert;
245     ji[1] = dAtom->getJy() + ( dt2 * Tb[1] ) * e_convert;
246     ji[2] = dAtom->getJz() + ( dt2 * Tb[2] ) * e_convert;
247    
248     // get the atom's rotation matrix
249    
250     A[0][0] = dAtom->getAxx();
251     A[0][1] = dAtom->getAxy();
252     A[0][2] = dAtom->getAxz();
253    
254     A[1][0] = dAtom->getAyx();
255     A[1][1] = dAtom->getAyy();
256     A[1][2] = dAtom->getAyz();
257    
258     A[2][0] = dAtom->getAzx();
259     A[2][1] = dAtom->getAzy();
260     A[2][2] = dAtom->getAzz();
261    
262    
263     // use the angular velocities to propagate the rotation matrix a
264     // full time step
265    
266    
267     angle = dt2 * ji[0] / dAtom->getIxx();
268     this->rotate( 1, 2, angle, ji, A ); // rotate about the x-axis
269    
270     angle = dt2 * ji[1] / dAtom->getIyy();
271     this->rotate( 2, 0, angle, ji, A ); // rotate about the y-axis
272    
273     angle = dt * ji[2] / dAtom->getIzz();
274     this->rotate( 0, 1, angle, ji, A ); // rotate about the z-axis
275    
276     angle = dt2 * ji[1] / dAtom->getIyy();
277     this->rotate( 2, 0, angle, ji, A ); // rotate about the y-axis
278    
279     angle = dt2 * ji[0] / dAtom->getIxx();
280     this->rotate( 1, 2, angle, ji, A ); // rotate about the x-axis
281    
282    
283     dAtom->setA( A );
284     dAtom->setJx( ji[0] );
285     dAtom->setJy( ji[1] );
286     dAtom->setJz( ji[2] );
287     }
288     }
289    
290     // calculate the forces
291    
292     for(j = 0; j < nAtoms; j++){
293     atoms[j]->zeroForces();
294     }
295    
296     for(j = 0; j < nSRI; j++){
297     srInteractions[j]->calc_forces();
298     }
299    
300     longRange->calc_forces();
301    
302     // move b
303    
304     for( j=0; j<nAtoms; j++ ){
305    
306     Rx[j] = atoms[j]->getX();
307     Ry[j] = atoms[j]->getY();
308     Rz[j] = atoms[j]->getZ();
309    
310     Vx[j] = atoms[j]->get_vx();
311     Vy[j] = atoms[j]->get_vy();
312     Vz[j] = atoms[j]->get_vz();
313    
314     Fx[j] = atoms[j]->getFx();
315     Fy[j] = atoms[j]->getFy();
316     Fz[j] = atoms[j]->getFz();
317     }
318    
319     v_constrain_b_( dt, nAtoms, mass, Rx, Ry, Rz, Vx, Vy, Vz,
320     Fx, Fy, Fz,
321     kE, n_constrained, constrained_dsqr,
322     constrained_i, constrained_j,
323     entry_plug->box_x,
324     entry_plug->box_y,
325     entry_plug->box_z );
326    
327     for( j=0; j<nAtoms; j++ ){
328    
329     atoms[j]->setX(Rx[j]);
330     atoms[j]->setY(Ry[j]);
331     atoms[j]->setZ(Rz[j]);
332    
333     atoms[j]->set_vx(Vx[j]);
334     atoms[j]->set_vy(Vy[j]);
335     atoms[j]->set_vz(Vz[j]);
336     }
337    
338     for( i=0; i< nAtoms; i++ ){
339    
340     if( atoms[i]->isDirectional() ){
341    
342     dAtom = (DirectionalAtom *)atoms[i];
343    
344     // get and convert the torque to body frame
345    
346     Tb[0] = dAtom->getTx();
347     Tb[1] = dAtom->getTy();
348     Tb[2] = dAtom->getTz();
349    
350     dAtom->lab2Body( Tb );
351    
352     // get the angular momentum, and complete the angular momentum
353     // half step
354    
355     ji[0] = dAtom->getJx() + ( dt2 * Tb[0] ) * e_convert;
356     ji[1] = dAtom->getJy() + ( dt2 * Tb[1] ) * e_convert;
357     ji[2] = dAtom->getJz() + ( dt2 * Tb[2] ) * e_convert;
358    
359     dAtom->setJx( ji[0] );
360     dAtom->setJy( ji[1] );
361     dAtom->setJz( ji[2] );
362     }
363     }
364 mmeineke 25
365     time = tl + 1;
366    
367 mmeineke 10 if( entry_plug->setTemp ){
368 mmeineke 25 if( !(time % vel_n) ) tStats->velocitize();
369 mmeineke 10 }
370 mmeineke 25 if( !(time % sample_n) ) dump_out->writeDump( time * dt );
371     if( !(time % status_n) ) e_out->writeStat( time * dt );
372 mmeineke 10 }
373     }
374     else{
375    
376     for( tl=0; tl<n_loops; tl++ ){
377    
378     kE = 0.0;
379     rot_kE= 0.0;
380     trans_kE = 0.0;
381    
382     for( i=0; i<nAtoms; i++ ){
383    
384     // velocity half step
385    
386     vx = atoms[i]->get_vx() +
387     ( dt2 * atoms[i]->getFx() / atoms[i]->getMass() ) * e_convert;
388     vy = atoms[i]->get_vy() +
389     ( dt2 * atoms[i]->getFy() / atoms[i]->getMass() ) * e_convert;
390     vz = atoms[i]->get_vz() +
391     ( dt2 * atoms[i]->getFz() / atoms[i]->getMass() ) * e_convert;
392    
393     // position whole step
394    
395     rx = atoms[i]->getX() + dt * vx;
396     ry = atoms[i]->getY() + dt * vy;
397     rz = atoms[i]->getZ() + dt * vz;
398    
399     atoms[i]->setX( rx );
400     atoms[i]->setY( ry );
401     atoms[i]->setZ( rz );
402    
403     atoms[i]->set_vx( vx );
404     atoms[i]->set_vy( vy );
405     atoms[i]->set_vz( vz );
406    
407     if( atoms[i]->isDirectional() ){
408    
409     dAtom = (DirectionalAtom *)atoms[i];
410    
411     // get and convert the torque to body frame
412    
413     Tb[0] = dAtom->getTx();
414     Tb[1] = dAtom->getTy();
415     Tb[2] = dAtom->getTz();
416    
417     dAtom->lab2Body( Tb );
418    
419     // get the angular momentum, and propagate a half step
420    
421     ji[0] = dAtom->getJx() + ( dt2 * Tb[0] ) * e_convert;
422     ji[1] = dAtom->getJy() + ( dt2 * Tb[1] ) * e_convert;
423     ji[2] = dAtom->getJz() + ( dt2 * Tb[2] ) * e_convert;
424    
425     // get the atom's rotation matrix
426    
427     A[0][0] = dAtom->getAxx();
428     A[0][1] = dAtom->getAxy();
429     A[0][2] = dAtom->getAxz();
430    
431     A[1][0] = dAtom->getAyx();
432     A[1][1] = dAtom->getAyy();
433     A[1][2] = dAtom->getAyz();
434    
435     A[2][0] = dAtom->getAzx();
436     A[2][1] = dAtom->getAzy();
437     A[2][2] = dAtom->getAzz();
438    
439    
440     // use the angular velocities to propagate the rotation matrix a
441     // full time step
442    
443    
444     angle = dt2 * ji[0] / dAtom->getIxx();
445     this->rotate( 1, 2, angle, ji, A ); // rotate about the x-axis
446    
447     angle = dt2 * ji[1] / dAtom->getIyy();
448     this->rotate( 2, 0, angle, ji, A ); // rotate about the y-axis
449    
450     angle = dt * ji[2] / dAtom->getIzz();
451     this->rotate( 0, 1, angle, ji, A ); // rotate about the z-axis
452    
453     angle = dt2 * ji[1] / dAtom->getIyy();
454     this->rotate( 2, 0, angle, ji, A ); // rotate about the y-axis
455    
456     angle = dt2 * ji[0] / dAtom->getIxx();
457     this->rotate( 1, 2, angle, ji, A ); // rotate about the x-axis
458    
459    
460     dAtom->setA( A );
461     dAtom->setJx( ji[0] );
462     dAtom->setJy( ji[1] );
463     dAtom->setJz( ji[2] );
464     }
465     }
466    
467     // calculate the forces
468    
469     for(j = 0; j < nAtoms; j++){
470     atoms[j]->zeroForces();
471     }
472    
473     for(j = 0; j < nSRI; j++){
474     srInteractions[j]->calc_forces();
475     }
476    
477     longRange->calc_forces();
478    
479     for( i=0; i< nAtoms; i++ ){
480    
481     // complete the velocity half step
482    
483     vx = atoms[i]->get_vx() +
484     ( dt2 * atoms[i]->getFx() / atoms[i]->getMass() ) * e_convert;
485     vy = atoms[i]->get_vy() +
486     ( dt2 * atoms[i]->getFy() / atoms[i]->getMass() ) * e_convert;
487     vz = atoms[i]->get_vz() +
488     ( dt2 * atoms[i]->getFz() / atoms[i]->getMass() ) * e_convert;
489    
490     atoms[i]->set_vx( vx );
491     atoms[i]->set_vy( vy );
492     atoms[i]->set_vz( vz );
493    
494 chuckv 134 // vx2 = vx * vx;
495     // vy2 = vy * vy;
496     // vz2 = vz * vz;
497 mmeineke 10
498     if( atoms[i]->isDirectional() ){
499    
500     dAtom = (DirectionalAtom *)atoms[i];
501    
502     // get and convert the torque to body frame
503    
504     Tb[0] = dAtom->getTx();
505     Tb[1] = dAtom->getTy();
506     Tb[2] = dAtom->getTz();
507    
508     dAtom->lab2Body( Tb );
509    
510     // get the angular momentum, and complete the angular momentum
511     // half step
512    
513     ji[0] = dAtom->getJx() + ( dt2 * Tb[0] ) * e_convert;
514     ji[1] = dAtom->getJy() + ( dt2 * Tb[1] ) * e_convert;
515     ji[2] = dAtom->getJz() + ( dt2 * Tb[2] ) * e_convert;
516    
517     jx2 = ji[0] * ji[0];
518     jy2 = ji[1] * ji[1];
519     jz2 = ji[2] * ji[2];
520    
521     rot_kE += (jx2 / dAtom->getIxx()) + (jy2 / dAtom->getIyy())
522     + (jz2 / dAtom->getIzz());
523    
524     dAtom->setJx( ji[0] );
525     dAtom->setJy( ji[1] );
526     dAtom->setJz( ji[2] );
527     }
528     }
529 mmeineke 25
530     time = tl + 1;
531    
532 mmeineke 10 if( entry_plug->setTemp ){
533 mmeineke 25 if( !(time % vel_n) ) tStats->velocitize();
534 mmeineke 10 }
535 mmeineke 25 if( !(time % sample_n) ) dump_out->writeDump( time * dt );
536     if( !(time % status_n) ) e_out->writeStat( time * dt );
537 mmeineke 10 }
538     }
539    
540     dump_out->writeFinal();
541    
542     delete dump_out;
543     delete e_out;
544     }
545    
546     void Symplectic::rotate( int axes1, int axes2, double angle, double ji[3],
547     double A[3][3] ){
548    
549     int i,j,k;
550     double sinAngle;
551     double cosAngle;
552     double angleSqr;
553     double angleSqrOver4;
554     double top, bottom;
555     double rot[3][3];
556     double tempA[3][3];
557     double tempJ[3];
558    
559     // initialize the tempA
560    
561     for(i=0; i<3; i++){
562     for(j=0; j<3; j++){
563     tempA[i][j] = A[i][j];
564     }
565     }
566    
567     // initialize the tempJ
568    
569     for( i=0; i<3; i++) tempJ[i] = ji[i];
570    
571     // initalize rot as a unit matrix
572    
573     rot[0][0] = 1.0;
574     rot[0][1] = 0.0;
575     rot[0][2] = 0.0;
576    
577     rot[1][0] = 0.0;
578     rot[1][1] = 1.0;
579     rot[1][2] = 0.0;
580    
581     rot[2][0] = 0.0;
582     rot[2][1] = 0.0;
583     rot[2][2] = 1.0;
584    
585     // use a small angle aproximation for sin and cosine
586    
587     angleSqr = angle * angle;
588     angleSqrOver4 = angleSqr / 4.0;
589     top = 1.0 - angleSqrOver4;
590     bottom = 1.0 + angleSqrOver4;
591    
592     cosAngle = top / bottom;
593     sinAngle = angle / bottom;
594    
595     rot[axes1][axes1] = cosAngle;
596     rot[axes2][axes2] = cosAngle;
597    
598     rot[axes1][axes2] = sinAngle;
599     rot[axes2][axes1] = -sinAngle;
600    
601     // rotate the momentum acoording to: ji[] = rot[][] * ji[]
602    
603     for(i=0; i<3; i++){
604     ji[i] = 0.0;
605     for(k=0; k<3; k++){
606     ji[i] += rot[i][k] * tempJ[k];
607     }
608     }
609    
610     // rotate the Rotation matrix acording to:
611     // A[][] = A[][] * transpose(rot[][])
612    
613    
614     // NOte for as yet unknown reason, we are setting the performing the
615     // calculation as:
616     // transpose(A[][]) = transpose(A[][]) * transpose(rot[][])
617    
618     for(i=0; i<3; i++){
619     for(j=0; j<3; j++){
620     A[j][i] = 0.0;
621     for(k=0; k<3; k++){
622     A[j][i] += tempA[k][i] * rot[j][k];
623     }
624     }
625     }
626     }