-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathFrame.cpp
1233 lines (1156 loc) · 37.3 KB
/
Frame.cpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
#include <cmath> // sqrt
#include <cstring> // memcpy, memset
#include "Frame.h"
#include "Constants.h" // SMALL
#include "CpptrajStdio.h"
const size_t Frame::COORDSIZE_ = 3 * sizeof(double);
// ---------- CONSTRUCTION/DESTRUCTION/ASSIGNMENT ------------------------------
/// CONSTRUCTOR
Frame::Frame( ) :
natom_(0),
maxnatom_(0),
ncoord_(0),
T_(0.0),
time_(0.0),
X_(0),
V_(0),
F_(0),
memIsExternal_(false)
{}
/// DESTRUCTOR
// Defined since this class can be inherited.
Frame::~Frame( ) {
if (!memIsExternal_) {
if (X_ != 0) delete[] X_;
if (V_ != 0) delete[] V_;
if (F_ != 0) delete[] F_;
}
}
// CONSTRUCTOR
Frame::Frame(int natomIn) :
natom_(natomIn),
maxnatom_(natomIn),
ncoord_(natomIn*3),
T_(0.0),
time_(0.0),
X_(0),
V_(0),
F_(0),
Mass_(natomIn, 1.0),
memIsExternal_(false)
{
if (ncoord_ > 0)
X_ = new double[ ncoord_ ];
}
// CONSTRUCTOR
Frame::Frame(std::vector<Atom> const& atoms) :
natom_(atoms.size()),
maxnatom_(natom_),
ncoord_(natom_*3),
T_(0.0),
time_(0.0),
X_(0),
V_(0),
F_(0),
memIsExternal_(false)
{
if (ncoord_ > 0) {
X_ = new double[ ncoord_ ];
Mass_.reserve( natom_ );
for (std::vector<Atom>::const_iterator atom = atoms.begin(); atom != atoms.end(); ++atom)
Mass_.push_back( (*atom).Mass() );
}
}
// CONSTRUCTOR
Frame::Frame(Frame const& frameIn, AtomMask const& maskIn) :
natom_( maskIn.Nselected() ),
maxnatom_(natom_),
ncoord_(natom_*3),
box_(frameIn.box_),
T_( frameIn.T_ ),
time_( frameIn.time_ ),
X_(0),
V_(0),
F_(0),
remd_indices_(frameIn.remd_indices_),
memIsExternal_(false)
{
if (ncoord_ > 0) {
Mass_.reserve(natom_);
X_ = new double[ ncoord_ ];
double* newX = X_;
double* newV = NULL;
double* newF = NULL;
bool dupV = false;
bool dupF = false;
// check if we should copy velocities
if ( frameIn.V_ != 0 ){
dupV = true;
newV = new double[ ncoord_ ];
}
// check if we should copy forces
if ( frameIn.F_ != 0 ){
dupF = true;
newF = new double[ ncoord_ ];
}
// do the copying
for (AtomMask::const_iterator atom = maskIn.begin(); atom != maskIn.end(); ++atom){
int oldcrd = ((*atom) * 3);
memcpy(newX, frameIn.X_ + oldcrd, COORDSIZE_);
newX += 3;
if ( dupV ){
memcpy(newV, frameIn.V_ + oldcrd, COORDSIZE_);
newV += 3;
}
if ( dupF ){
memcpy(newF, frameIn.F_ + oldcrd, COORDSIZE_);
newF += 3;
}
Mass_.push_back( frameIn.Mass_[*atom] );
}
}
}
// CONSTRUCTOR
Frame::Frame(int natom, double* Xptr) :
natom_(natom),
maxnatom_(natom),
ncoord_(natom*3),
X_(Xptr),
V_(0),
F_(0),
Mass_(natom, 1.0),
memIsExternal_(true)
{
if (Xptr == 0) {
mprinterr("Internal Error: in Frame::Frame(int,double*) pointer is NULL.\n");
ncoord_ = maxnatom_ = natom_ = 0;
}
}
// COPY CONSTRUCTOR
Frame::Frame(const Frame& rhs) :
natom_(rhs.natom_),
maxnatom_(rhs.maxnatom_),
ncoord_(rhs.ncoord_),
box_(rhs.box_),
T_(rhs.T_),
time_(rhs.time_),
X_(0),
V_(0),
F_(0),
remd_indices_(rhs.remd_indices_),
Mass_(rhs.Mass_),
memIsExternal_(rhs.memIsExternal_)
{
if (memIsExternal_) {
X_ = rhs.X_;
V_ = rhs.V_;
F_ = rhs.F_;
} else {
// Copy coords/velo/forces; allocate for maxnatom but copy natom
int maxncoord = maxnatom_ * 3;
if (rhs.X_!=0) {
X_ = new double[ maxncoord ];
memcpy(X_, rhs.X_, natom_ * COORDSIZE_);
}
if (rhs.V_!=0) {
V_ = new double[ maxncoord ];
memcpy(V_, rhs.V_, natom_ * COORDSIZE_);
}
if (rhs.F_!=0) {
F_ = new double[ maxncoord ];
memcpy(F_, rhs.F_, natom_ * COORDSIZE_);
}
}
}
// Frame::swap()
void Frame::swap(Frame &first, Frame &second) {
using std::swap;
swap(first.natom_, second.natom_);
swap(first.maxnatom_, second.maxnatom_);
swap(first.ncoord_, second.ncoord_);
swap(first.T_, second.T_);
swap(first.time_, second.time_);
swap(first.X_, second.X_);
swap(first.V_, second.V_);
swap(first.F_, second.F_);
first.remd_indices_.swap(second.remd_indices_);
first.Mass_.swap(second.Mass_);
swap(first.memIsExternal_, second.memIsExternal_);
swap( first.box_[0], second.box_[0] );
swap( first.box_[1], second.box_[1] );
swap( first.box_[2], second.box_[2] );
swap( first.box_[3], second.box_[3] );
swap( first.box_[4], second.box_[4] );
swap( first.box_[5], second.box_[5] );
}
// Frame::operator=()
/** Assignment operator using copy/swap idiom. */
Frame &Frame::operator=(Frame rhs) {
if (memIsExternal_)
mprinterr("Internal Error: Attempting to assign to Frame with external memory.\n");
else if (rhs.memIsExternal_) {
// Do not use copy/swap here since we do not want to free external mem.
// FIXME: seems like this could be pretty inefficient. Should assignment
// not use copy/swap?
natom_ = rhs.natom_;
maxnatom_ = rhs.maxnatom_;
ncoord_ = rhs.ncoord_;
box_ = rhs.box_;
T_ = rhs.T_;
time_ = rhs.time_;
remd_indices_ = rhs.remd_indices_;
Mass_ = rhs.Mass_;
memIsExternal_ = false;
if (X_ != 0) delete[] X_;
if (V_ != 0) delete[] V_;
if (F_ != 0) delete[] F_;
F_ = V_ = X_ = 0;
if (maxnatom_ > 0) {
int maxncoord = maxnatom_ * 3;
X_ = new double[ maxncoord ];
std::copy( rhs.X_, rhs.X_ + ncoord_, X_ );
if (rhs.V_ != 0) {
V_ = new double[ maxncoord ];
std::copy( rhs.V_, rhs.V_ + ncoord_, V_ );
}
if (rhs.F_ != 0) {
F_ = new double[ maxncoord ];
std::copy( rhs.F_, rhs.F_ + ncoord_, F_ );
}
}
} else
swap(*this, rhs);
return *this;
}
// ---------- CONVERT TO/FROM CRDtype ------------------------------------------
// Frame::SetFromCRD()
void Frame::SetFromCRD(CRDtype const& farray, int numCrd, int numBoxCrd, bool hasVel) {
int f_ncoord = numCrd;
if (f_ncoord > maxnatom_*3) {
mprinterr("Error: Float array size (%i) > max #coords in frame (%i)\n",
f_ncoord, maxnatom_*3);
return;
}
ncoord_ = f_ncoord;
natom_ = ncoord_ / 3;
for (int ix = 0; ix < ncoord_; ++ix)
X_[ix] = (double)farray[ix];
if (hasVel && V_ != 0) {
for (int iv = 0; iv < ncoord_; ++iv)
V_[iv] = (double)farray[f_ncoord++];
}
for (int ib = 0; ib < numBoxCrd; ++ib)
box_[ib] = (double)farray[f_ncoord++];
}
// Frame::SetFromCRD()
void Frame::SetFromCRD(CRDtype const& crdIn, AtomMask const& mask, int numCrd,
int numBoxCrd, bool hasVel)
{
if (mask.Nselected() > maxnatom_) {
mprinterr("Internal Error: Selected # atoms in float array (%i) > max #atoms in frame (%i)\n",
mask.Nselected(), maxnatom_);
return;
}
natom_ = mask.Nselected();
ncoord_ = natom_ * 3;
unsigned int ix = 0;
unsigned int iv = 0;
for (AtomMask::const_iterator atom = mask.begin(); atom != mask.end(); ++atom) {
unsigned int xoffset = ((unsigned int)(*atom)) * 3;
X_[ix++] = (double)crdIn[xoffset ];
X_[ix++] = (double)crdIn[xoffset+1];
X_[ix++] = (double)crdIn[xoffset+2];
if (hasVel && V_ != 0) {
unsigned int voffset = numCrd + xoffset;
V_[iv++] = (double)crdIn[voffset ];
V_[iv++] = (double)crdIn[voffset+1];
V_[iv++] = (double)crdIn[voffset+2];
}
}
int f_ncoord = (int)crdIn.size() - numBoxCrd;
for (int ib = 0; ib < numBoxCrd; ++ib)
box_[ib] = (double)crdIn[f_ncoord++];
}
// Frame::ConvertToCRD()
Frame::CRDtype Frame::ConvertToCRD(int numBoxCrd, bool hasVel) const {
int nvel;
if (hasVel)
nvel = ncoord_;
else
nvel = 0;
CRDtype farray;
farray.reserve( ncoord_ + nvel + numBoxCrd );
for (int ix = 0; ix < ncoord_; ++ix)
farray.push_back( (float)X_[ix] );
for (int iv = 0; iv < nvel; ++iv )
farray.push_back( (float)V_[iv] );
for (int ib = 0; ib < numBoxCrd; ++ib)
farray.push_back( (float)box_[ib] );
return farray;
}
// ---------- ACCESS INTERNAL DATA ---------------------------------------------
// Frame::printAtomCoord()
void Frame::printAtomCoord(int atom) const {
int atmidx = atom * 3;
if (atmidx >= ncoord_) return;
mprintf("%i: %f %f %f\n",atom+1,X_[atmidx],X_[atmidx+1],X_[atmidx+2]);
}
// Frame::Info()
void Frame::Info(const char *msg) const {
if (msg!=0)
mprintf("\tFrame [%s]:",msg);
else
mprintf("\tFrame:");
mprintf("%i atoms, %i coords",natom_, ncoord_);
if (V_!=0) mprintf(", with Velocities");
if (!remd_indices_.empty()) mprintf(", with replica indices");
mprintf("\n");
}
// Frame::IncreaseX()
void Frame::IncreaseX() {
maxnatom_ += 500;
double *newX = new double[ maxnatom_ * 3 ];
if (X_!=0) {
memcpy(newX, X_, natom_ * COORDSIZE_);
if (!memIsExternal_)
delete[] X_;
else
memIsExternal_ = false;
}
X_ = newX;
}
/** Set atom/coord count to zero but do not clear memory. */
void Frame::ClearAtoms() {
natom_ = 0;
ncoord_ = 0;
}
// Frame::AddXYZ()
/** Append the given XYZ coord to this frame. */
void Frame::AddXYZ(const double *XYZin) {
if (XYZin == 0) return;
if (natom_ >= maxnatom_)
IncreaseX();
memcpy(X_ + ncoord_, XYZin, COORDSIZE_);
++natom_;
ncoord_ += 3;
}
// Frame::AddVec3()
void Frame::AddVec3(Vec3 const& vIn) {
if (natom_ >= maxnatom_)
IncreaseX();
memcpy(X_ + ncoord_, vIn.Dptr(), COORDSIZE_);
++natom_;
ncoord_ += 3;
}
void Frame::SetMass(std::vector<Atom> const& atoms) {
// No memory reallocation TODO allow atoms to be larger?
if (natom_ != (int)atoms.size()) {
mprinterr("Internal Error: Size of atoms array is %zu, Frame size is %i\n",
atoms.size(), natom_);
} else { // Assume mass has been allocated
for (unsigned int i = 0; i != atoms.size(); i++)
Mass_[i] = atoms[i].Mass();
}
}
// ---------- FRAME MEMORY ALLOCATION/REALLOCATION -----------------------------
/** \return True if reallocation of coordinate arrray must occur based on
* given number of atoms.
*/
bool Frame::ReallocateX(int natomIn) {
natom_ = natomIn;
ncoord_ = natom_ * 3;
if (natom_ > maxnatom_ || memIsExternal_) {
if (memIsExternal_)
memIsExternal_ = false;
else if (X_ != 0)
delete[] X_;
X_ = new double[ ncoord_ ];
maxnatom_ = natom_;
return true;
}
return false;
}
// Frame::SetupFrame()
int Frame::SetupFrame(int natomIn) {
ReallocateX( natomIn );
if (V_ != 0) delete[] V_;
Mass_.assign(natomIn, 1.0);
return 0;
}
// Frame::SetupFrameM()
int Frame::SetupFrameM(std::vector<Atom> const& atoms) {
return SetupFrameV( atoms, CoordinateInfo() );
}
// Frame::SetupFrameXM()
int Frame::SetupFrameXM(Darray const& Xin, Darray const& massIn) {
ReallocateX( Xin.size() / 3 );
// Copy coords
std::copy( Xin.begin(), Xin.end(), X_ );
// Copy masses, or set all to 1.0 if input masses are empty.
if (!massIn.empty())
Mass_ = massIn;
else
Mass_.assign( natom_, 1.0 );
if (V_ != 0) delete[] V_;
return 0;
}
// Frame::SetupFrameV()
int Frame::SetupFrameV(std::vector<Atom> const& atoms, CoordinateInfo const& cinfo) {
bool reallocate = ReallocateX( atoms.size() );
// Velocity
if (cinfo.HasVel()) {
if (reallocate || V_ == 0) {
if (V_ != 0) delete[] V_;
V_ = new double[ maxnatom_*3 ];
// Since velocity might not be read in, initialize it to 0.
memset(V_, 0, maxnatom_ * COORDSIZE_);
}
} else {
if (V_ != 0) delete[] V_;
V_ = 0;
}
// Force
if (cinfo.HasForce()) {
if (reallocate || F_ == 0) {
if (F_ != 0) delete[] F_;
F_ = new double [ maxnatom_*3 ];
// Since force might not be read in, initialize it to 0.
memset(F_, 0, maxnatom_ * COORDSIZE_);
}
}
// Mass
if (reallocate || Mass_.empty())
Mass_.resize(maxnatom_);
Darray::iterator mass = Mass_.begin();
for (std::vector<Atom>::const_iterator atom = atoms.begin();
atom != atoms.end(); ++atom)
*(mass++) = (*atom).Mass();
// Replica indices
remd_indices_.assign( cinfo.ReplicaDimensions().Ndims(), 0 );
return 0;
}
// Frame::SetupFrameFromMask()
/** Set up frame to hold # selected atoms in given mask. If mass
* information is passed in store the masses corresponding to
* selected atoms in the mask. No velocity info. Only reallocate
* memory if Nselected > maxnatom.
*/
int Frame::SetupFrameFromMask(AtomMask const& maskIn, std::vector<Atom> const& atoms) {
bool reallocate = ReallocateX( maskIn.Nselected() );
if (reallocate || Mass_.empty())
Mass_.resize(maxnatom_);
// Copy masses according to maskIn
Darray::iterator mass = Mass_.begin();
for (AtomMask::const_iterator atom = maskIn.begin(); atom != maskIn.end(); ++atom)
*(mass++) = atoms[ *atom ].Mass();
return 0;
}
// ---------- FRAME SETUP OF COORDINATES ---------------------------------------
// Frame::SetCoordinates()
void Frame::SetCoordinates(Frame const& frameIn, AtomMask const& maskIn) {
if (maskIn.Nselected() > maxnatom_) {
mprinterr("Error: SetCoordinates: Mask [%s] selected (%i) > max natom (%i)\n",
maskIn.MaskString(), maskIn.Nselected(), maxnatom_);
return;
}
natom_ = maskIn.Nselected();
ncoord_ = natom_ * 3;
box_ = frameIn.box_;
T_ = frameIn.T_;
time_ = frameIn.time_;
remd_indices_ = frameIn.remd_indices_;
double* newXptr = X_;
for (AtomMask::const_iterator atom = maskIn.begin(); atom != maskIn.end(); ++atom)
{
memcpy( newXptr, frameIn.X_ + ((*atom) * 3), COORDSIZE_);
newXptr += 3;
}
}
// Frame::SetCoordinates()
void Frame::SetCoordinates(Frame const& frameIn) {
if (frameIn.natom_ > maxnatom_) {
mprinterr("Error: Frame::SetCoordinates: Input frame atoms (%i) > max natom (%i)\n",
frameIn.natom_, maxnatom_);
return;
}
natom_ = frameIn.natom_;
ncoord_ = natom_ * 3;
memcpy(X_, frameIn.X_, natom_ * COORDSIZE_);
}
int Frame::SetCoordinates(int natom, double* Xptr) {
if (!memIsExternal_)
mprinterr("Internal Error: Frame memory is internal, not setting from external pointer.\n");
else if (natom != natom_)
mprinterr("Internal Error: Frame set up for %i atoms, external memory has %i atoms.\n",
natom_, natom);
else {
X_ = Xptr;
return 0;
}
return 1;
}
// Frame::SetFrame()
void Frame::SetFrame(Frame const& frameIn, AtomMask const& maskIn) {
if (maskIn.Nselected() > maxnatom_) {
mprinterr("Error: SetFrame: Mask [%s] selected (%i) > max natom (%i)\n",
maskIn.MaskString(), maskIn.Nselected(), maxnatom_);
return;
}
natom_ = maskIn.Nselected();
ncoord_ = natom_ * 3;
// Copy T/box
box_ = frameIn.box_;
T_ = frameIn.T_;
time_ = frameIn.time_;
remd_indices_ = frameIn.remd_indices_;
double* newXptr = X_;
Darray::iterator mass = Mass_.begin();
if (frameIn.V_ != 0 && V_ != 0) {
// Copy Coords/Mass/Velo
double *newVptr = V_;
for (AtomMask::const_iterator atom = maskIn.begin(); atom != maskIn.end(); ++atom)
{
int oldcrd = ((*atom) * 3);
memcpy( newXptr, frameIn.X_ + oldcrd, COORDSIZE_);
newXptr += 3;
memcpy( newVptr, frameIn.V_ + oldcrd, COORDSIZE_);
newVptr += 3;
*mass = frameIn.Mass_[*atom];
++mass;
}
} else {
// Copy coords/mass only
for (AtomMask::const_iterator atom = maskIn.begin(); atom != maskIn.end(); ++atom)
{
memcpy( newXptr, frameIn.X_ + ((*atom) * 3), COORDSIZE_);
newXptr += 3;
*mass = frameIn.Mass_[*atom];
++mass;
}
}
}
// ---------- FRAME SETUP WITH ATOM MAPPING ------------------------------------
// NOTE: Should these map routines be part of a child class used only by AtomMap?
// Frame::SetCoordinatesByMap()
/** Reorder atoms in the given target frame so that it matches the given
* atom map with format (tgtAtom = Map[refAtom]). End result is that
* this frame will be in the same order as the original reference. Assumes
* that the map is complete (i.e. no unmapped atoms).
*/
void Frame::SetCoordinatesByMap(Frame const& tgtIn, std::vector<int> const& mapIn) {
if (tgtIn.natom_ > maxnatom_) {
mprinterr("Error: SetCoordinatesByMap: # Input map frame atoms (%i) > max atoms (%i)\n",
tgtIn.natom_, maxnatom_);
return;
}
if ((int)mapIn.size() != tgtIn.natom_) {
mprinterr("Error: SetCoordinatesByMap: Input map size (%zu) != input frame natom (%i)\n",
mapIn.size(), tgtIn.natom_);
return;
}
natom_ = tgtIn.natom_;
ncoord_ = natom_ * 3;
box_ = tgtIn.box_;
T_ = tgtIn.T_;
time_ = tgtIn.time_;
remd_indices_ = tgtIn.remd_indices_;
double* newXptr = X_;
Darray::iterator newmass = Mass_.begin();
if (tgtIn.V_ != 0 && V_ != 0) {
// Copy Coords/Mass/Velo
double *newVptr = V_;
for (std::vector<int>::const_iterator refatom = mapIn.begin();
refatom != mapIn.end(); ++refatom)
{
int idx = ((*refatom) * 3);
memcpy( newXptr, tgtIn.X_ + idx, COORDSIZE_ );
newXptr += 3;
memcpy( newVptr, tgtIn.V_ + idx, COORDSIZE_ );
newVptr += 3;
*(newmass++) = tgtIn.Mass_[*refatom];
}
} else {
// Copy Coords/Mass only
for (std::vector<int>::const_iterator refatom = mapIn.begin();
refatom != mapIn.end(); ++refatom)
{
memcpy( newXptr, tgtIn.X_ + ((*refatom) * 3), COORDSIZE_ );
newXptr += 3;
*(newmass++) = tgtIn.Mass_[*refatom];
}
}
}
// Frame::StripUnmappedAtoms()
/** Set this frame to include only atoms from the given reference that are
* mapped: Map[newatom] = refatom (refatom -> this frame).
*/
void Frame::StripUnmappedAtoms(Frame const& refIn, std::vector<int> const& mapIn) {
if (refIn.natom_ > maxnatom_) {
mprinterr("Error: StripUnmappedAtoms: # Input map frame atoms (%i) > max atoms (%i)\n",
refIn.natom_, maxnatom_);
return;
}
if ((int)mapIn.size() != refIn.natom_) {
mprinterr("Error: StripUnmappedAtoms: Input map size (%zu) != input frame natom (%i)\n",
mapIn.size(), refIn.natom_);
return;
}
box_ = refIn.box_;
T_ = refIn.T_;
time_ = refIn.time_;
remd_indices_ = refIn.remd_indices_;
double* newXptr = X_;
double* refptr = refIn.X_;
for (std::vector<int>::const_iterator refatom = mapIn.begin();
refatom != mapIn.end(); ++refatom)
{
if (*refatom != -1) {
memcpy( newXptr, refptr, COORDSIZE_ );
newXptr += 3;
}
refptr += 3;
}
ncoord_ = (int)(newXptr - X_);
natom_ = ncoord_ / 3;
}
// Frame::ModifyByMap()
/** Set this frame to include only atoms from the given target frame, remapped
* according to the given atom map: Map[newatom] = oldatom (oldatom -> frameIn)
*/
void Frame::ModifyByMap(Frame const& frameIn, std::vector<int> const& mapIn) {
if ((int)mapIn.size() > maxnatom_) {
mprinterr("Error: SetTargetByMap: Input map size (%zu) > this frame max natom (%i)\n",
mapIn.size(), maxnatom_);
return;
}
box_ = frameIn.box_;
T_ = frameIn.T_;
time_ = frameIn.time_;
remd_indices_ = frameIn.remd_indices_;
double* Xptr = X_;
for (std::vector<int>::const_iterator oldatom = mapIn.begin();
oldatom != mapIn.end(); ++oldatom)
{
if (*oldatom != -1) {
memcpy( Xptr, frameIn.X_ + ((*oldatom) * 3), COORDSIZE_ );
Xptr += 3;
}
}
ncoord_ = (int)(Xptr - X_);
natom_ = ncoord_ / 3;
}
// ---------- BASIC ARITHMETIC -------------------------------------------------
// Frame::ZeroCoords()
/** Set all coords to 0.0 */
void Frame::ZeroCoords( ) {
memset(X_, 0, natom_ * COORDSIZE_);
}
// Frame::operator+=()
Frame &Frame::operator+=(const Frame& rhs) {
// For now ensure same natom
if (natom_ != rhs.natom_) {
mprinterr("Error: Frame::operator+=: Frames have different natom.\n");
return *this;
}
for (int i = 0; i < ncoord_; ++i)
X_[i] += rhs.X_[i];
return *this;
}
// Frame::operator-=()
Frame &Frame::operator-=(const Frame& rhs) {
// For now ensure same natom
if (natom_ != rhs.natom_) {
mprinterr("Error: Frame::operator-=: Frames have different natom.\n");
return *this;
}
for (int i = 0; i < ncoord_; ++i)
X_[i] -= rhs.X_[i];
return *this;
}
// Frame::operator*=()
Frame &Frame::operator*=(const Frame& rhs) {
// For now ensure same natom
if (natom_ != rhs.natom_) {
mprinterr("Error: Frame::operator*=: Frames have different natom.\n");
return *this;
}
for (int i = 0; i < ncoord_; ++i)
X_[i] *= rhs.X_[i];
return *this;
}
// Frame::operator*()
// NOTE: return const instance and not const reference to disallow
// nonsense statements like (a * b) = c
const Frame Frame::operator*(const Frame& rhs) const {
return (Frame(*this) *= rhs);
}
const Frame Frame::operator-(const Frame& rhs) const {
return (Frame(*this) -= rhs);
}
// Frame::Divide()
/** Divide all coord values of dividend by divisor and store in this frame.
*/
int Frame::Divide(Frame const& dividend, double divisor) {
if (divisor < Constants::SMALL) {
mprinterr("Error: Frame::Divide(Frame,divisor): Detected divide by 0.\n");
return 1;
}
// For now ensure same natom
if (natom_ != dividend.natom_) {
mprinterr("Error: Frame::Divide: Frames have different natom.\n");
return 1;
}
for (int i=0; i < ncoord_; ++i)
X_[i] = dividend.X_[i] / divisor;
return 0;
}
// Frame::Divide()
/** Divide all coord values of this frame by divisor.
*/
void Frame::Divide(double divisor) {
if (divisor < Constants::SMALL) {
mprinterr("Error: Frame::Divide(divisor): Detected divide by 0.\n");
return;
}
for (int i = 0; i < ncoord_; ++i)
X_[i] /= divisor;
}
void Frame::Multiply(double mult) {
for (int i = 0; i < ncoord_; i++)
X_[i] *= mult;
}
// Frame::AddByMask()
int Frame::AddByMask(Frame const& frameIn, AtomMask const& maskIn) {
if (maskIn.Nselected() > maxnatom_) {
mprinterr("Error: AddByMask: Input mask #atoms (%i) > frame #atoms (%i)\n",
maskIn.Nselected(), maxnatom_);
return 1;
}
unsigned int xidx = 0;
for (AtomMask::const_iterator atom = maskIn.begin(); atom != maskIn.end(); ++atom)
{
unsigned int fidx = (*atom) * 3;
X_[xidx ] += frameIn.X_[fidx ];
X_[xidx+1] += frameIn.X_[fidx+1];
X_[xidx+2] += frameIn.X_[fidx+2];
xidx += 3;
}
return 0;
}
// ---------- COORDINATE MANIPULATION ------------------------------------------
// Frame::Scale()
void Frame::Scale(AtomMask const& maskIn, double sx, double sy, double sz) {
for (AtomMask::const_iterator atom = maskIn.begin();
atom != maskIn.end(); atom++)
{
unsigned int xidx = *atom * 3;
X_[xidx ] *= sx;
X_[xidx+1] *= sy;
X_[xidx+2] *= sz;
}
}
// Frame::CenterOnOrigin()
/** \return translation vector from origin to original center. */
Vec3 Frame::CenterOnOrigin(bool useMassIn)
{
Vec3 center;
if (useMassIn)
center = VCenterOfMass(0, natom_); // TODO: Replace with total version?
else
center = VGeometricCenter(0, natom_);
//mprinterr(" REF FRAME CENTER: %lf %lf %lf\n",Trans[0],Trans[1],Trans[2]); //DEBUG
// center contains translation from origin -> Ref, therefore
// -center is translation from Ref -> origin.
NegTranslate(center);
return center;
}
// ---------- COORDINATE CALCULATION -------------------------------------------
// Frame::RMSD()
double Frame::RMSD( Frame& Ref, bool useMassIn ) {
Matrix_3x3 U;
Vec3 Trans, refTrans;
return RMSD( Ref, U, Trans, refTrans, useMassIn );
}
// Frame::RMSD()
/** Get the RMSD of this Frame to Ref Frame. Ref frame must contain the same
* number of atoms as this Frame - should be checked for before this routine
* is called.
* \param Ref frame to calc RMSD to. Will be translated to origin.
* \param U Will be set with best-fit rotation matrix.
* \param Trans Will be set with translation of coords to origin.
* \param refTrans Will be set with translation from origin to Ref.
* \param useMassIn If true, mass-weight everything.
* To reproduce the fit perform the first translation (Trans),
* then rotate (U), then the second translation (refTrans).
*/
double Frame::RMSD( Frame& Ref, Matrix_3x3& U, Vec3& Trans, Vec3& refTrans, bool useMassIn)
{
// Rotation will occur around geometric center/center of mass. Coords are
// shifted to common CoM first (Trans), then to original reference
// location (refTrans).
refTrans = Ref.CenterOnOrigin(useMassIn);
return RMSD_CenteredRef( Ref, U, Trans, useMassIn );
}
// Frame::RMSD_CenteredRef()
/** Calculate RMSD of this Frame to Ref Frame previously centered at origin.
*/
double Frame::RMSD_CenteredRef( Frame const& Ref, bool useMassIn ) {
Matrix_3x3 U;
Vec3 Trans;
return RMSD_CenteredRef( Ref, U, Trans, useMassIn );
}
static inline void normalize(double* vIn) {
double b = 1.0 / sqrt(vIn[0]*vIn[0] + vIn[1]*vIn[1] + vIn[2]*vIn[2]);
vIn[0] *= b;
vIn[1] *= b;
vIn[2] *= b;
}
// Frame::RMSD_CenteredRef()
/** Get the RMSD of this Frame to given Reference Frame. Ref frame must contain
* the same number of atoms as this Frame and should have already been
* translated to coordinate origin (neither is checked for in the interest
* of speed). This frame will be translated to the origin.
* \param Ref Previously-centered frame to calc RMSD to.
* \param U Will be set to the best-fit rotation matrix
* \param Trans will contain translation vector for this frame to origin.
* \param useMassIn If true, mass-weight everything.
*/
double Frame::RMSD_CenteredRef( Frame const& Ref, Matrix_3x3& U, Vec3& Trans, bool useMassIn)
{
double total_mass, cp[3], sig3, b[9];
// Rotation will occur around geometric center/center of mass
Trans.Zero();
if (useMassIn) {
Darray::iterator mass = Mass_.begin();
total_mass = 0.0;
for (int ix = 0; ix < ncoord_; ix += 3) {
total_mass += *mass;
Trans[0] += (X_[ix ] * (*mass));
Trans[1] += (X_[ix+1] * (*mass));
Trans[2] += (X_[ix+2] * (*mass));
++mass;
}
} else {
total_mass = (double)natom_;
for (int ix = 0; ix < ncoord_; ix += 3) {
Trans[0] += X_[ix ];
Trans[1] += X_[ix+1];
Trans[2] += X_[ix+2];
}
}
if (total_mass<Constants::SMALL) {
mprinterr("Error: Frame::RMSD: Divide by zero.\n");
return -1;
}
Trans[0] /= total_mass;
Trans[1] /= total_mass;
Trans[2] /= total_mass;
//mprinterr(" FRAME COM: %lf %lf %lf\n",Trans[0],Trans[1],Trans[2]); //DEBUG
// Shift to common COM
Trans.Neg();
Translate(Trans);
//for (int i = 0; i < natom_; i++) {
// mprinterr("\tSHIFTED FRAME %i: %f %f %f\n",i,X_[i*3],X_[i*3+1],X_[i*3+2]); //DEBUG
// mprinterr("\tSHIFTED REF %i : %f %f %f\n",i,Ref.X_[i*3],Ref.X_[i*3+1],Ref.X_[i*3+2]); //DEBUG
//}
// Use Kabsch algorithm to calculate optimum rotation matrix.
// U = [(RtR)^.5][R^-1]
double mwss = 0.0;
Matrix_3x3 rot(0.0);
// Calculate covariance matrix of Coords and Reference (R = Xt * Ref)
Darray::iterator mass = Mass_.begin();
double atom_mass = 1.0;
for (int i = 0; i < ncoord_; i += 3)
{
double xt = X_[i ];
double yt = X_[i+1];
double zt = X_[i+2];
double xr = Ref.X_[i ];
double yr = Ref.X_[i+1];
double zr = Ref.X_[i+2];
// Use atom_mass to hold mass for this atom if specified
if (useMassIn)
atom_mass = *(mass++);
mwss += atom_mass * ( (xt*xt)+(yt*yt)+(zt*zt)+(xr*xr)+(yr*yr)+(zr*zr) );
// Calculate the Kabsch matrix: R = (rij) = Sum(yni*xnj)
rot[0] += atom_mass*xt*xr;
rot[1] += atom_mass*xt*yr;
rot[2] += atom_mass*xt*zr;
rot[3] += atom_mass*yt*xr;
rot[4] += atom_mass*yt*yr;
rot[5] += atom_mass*yt*zr;
rot[6] += atom_mass*zt*xr;
rot[7] += atom_mass*zt*yr;
rot[8] += atom_mass*zt*zr;
}
mwss *= 0.5; // E0 = 0.5*Sum(xn^2+yn^2)
//DEBUG
//mprinterr("ROT:\n%lf %lf %lf\n%lf %lf %lf\n%lf %lf %lf\n",
// rot[0],rot[1],rot[2],rot[3],rot[4],rot[5],rot[6],rot[7],rot[8]);
//mprinterr("MWSS: %lf\n",mwss);
// calculate Kabsch matrix multiplied by its transpose: RtR
Matrix_3x3 Evector = rot.TransposeMult(rot);
// Diagonalize
Vec3 Eigenvalue;
if (Evector.Diagonalize_Sort( Eigenvalue )) return 0;
// a3 = a1 x a2
Evector[6] = (Evector[1]*Evector[5]) - (Evector[2]*Evector[4]);
Evector[7] = (Evector[2]*Evector[3]) - (Evector[0]*Evector[5]);
Evector[8] = (Evector[0]*Evector[4]) - (Evector[1]*Evector[3]);
// Evector dot transpose rot: b = R . ak
b[0] = Evector[0]*rot[0] + Evector[1]*rot[3] + Evector[2]*rot[6];
b[1] = Evector[0]*rot[1] + Evector[1]*rot[4] + Evector[2]*rot[7];
b[2] = Evector[0]*rot[2] + Evector[1]*rot[5] + Evector[2]*rot[8];
normalize(b);
b[3] = Evector[3]*rot[0] + Evector[4]*rot[3] + Evector[5]*rot[6];
b[4] = Evector[3]*rot[1] + Evector[4]*rot[4] + Evector[5]*rot[7];
b[5] = Evector[3]*rot[2] + Evector[4]*rot[5] + Evector[5]*rot[8];
normalize(b+3);
b[6] = Evector[6]*rot[0] + Evector[7]*rot[3] + Evector[8]*rot[6];
b[7] = Evector[6]*rot[1] + Evector[7]*rot[4] + Evector[8]*rot[7];
b[8] = Evector[6]*rot[2] + Evector[7]*rot[5] + Evector[8]*rot[8];
normalize(b+6);
// b3 = b1 x b2
cp[0] = (b[1]*b[5]) - (b[2]*b[4]);
cp[1] = (b[2]*b[3]) - (b[0]*b[5]);
cp[2] = (b[0]*b[4]) - (b[1]*b[3]);
if ( (cp[0]*b[6] + cp[1]*b[7] + cp[2]*b[8]) < 0.0 )
sig3 = -1.0;
else
sig3 = 1.0;
b[6] = cp[0];
b[7] = cp[1];
b[8] = cp[2];
// U has the best rotation
U[0] = (Evector[0]*b[0]) + (Evector[3]*b[3]) + (Evector[6]*b[6]);
U[1] = (Evector[1]*b[0]) + (Evector[4]*b[3]) + (Evector[7]*b[6]);
U[2] = (Evector[2]*b[0]) + (Evector[5]*b[3]) + (Evector[8]*b[6]);
U[3] = (Evector[0]*b[1]) + (Evector[3]*b[4]) + (Evector[6]*b[7]);
U[4] = (Evector[1]*b[1]) + (Evector[4]*b[4]) + (Evector[7]*b[7]);
U[5] = (Evector[2]*b[1]) + (Evector[5]*b[4]) + (Evector[8]*b[7]);
U[6] = (Evector[0]*b[2]) + (Evector[3]*b[5]) + (Evector[6]*b[8]);
U[7] = (Evector[1]*b[2]) + (Evector[4]*b[5]) + (Evector[7]*b[8]);
U[8] = (Evector[2]*b[2]) + (Evector[5]*b[5]) + (Evector[8]*b[8]);
// E=E0-sqrt(mu1)-sqrt(mu2)-sig3*sqrt(mu3)
double rms_return = mwss - sqrt(fabs(Eigenvalue[0]))
- sqrt(fabs(Eigenvalue[1]))
- (sig3*sqrt(fabs(Eigenvalue[2])));
if (rms_return<0) {
//mprinterr("RMS returned is <0 before sqrt, setting to 0 (%f)\n", rms_return);