00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025 #include <stdio.h>
00026 #include <stdlib.h>
00027 #include <math.h>
00028 #include "Measure.h"
00029 #include "AtomSel.h"
00030 #include "utilities.h"
00031 #include "ResizeArray.h"
00032 #include "MoleculeList.h"
00033 #include "Inform.h"
00034 #include "Timestep.h"
00035 #include "VMDApp.h"
00036 #include "WKFThreads.h"
00037 #include "WKFUtils.h"
00038 #include "SpatialSearch.h"
00039
00040 class AtomSelThr : public AtomSel
00041 {
00042 public:
00043 AtomSelThr(AtomSel *osel, wkf_mutex_t *olock)
00044 : AtomSel(NULL, osel->molid()),
00045 sel(osel), lock(olock) {
00046 if (sel) {
00047 selected=sel->selected;
00048 num_atoms=sel->num_atoms;
00049 which_frame=sel->which_frame;
00050 if (sel->on) {
00051 on = new int[num_atoms];
00052 memcpy(on,sel->on,num_atoms*sizeof(int));
00053 }
00054 } else {
00055 selected=-1;
00056 num_atoms=-1;
00057 which_frame=-1;
00058 }
00059 }
00060
00061 ~AtomSelThr() {
00062 sel=NULL;
00063 }
00064
00065
00066 private:
00067 AtomSelThr() : AtomSel(NULL,-1) {};
00068 AtomSelThr& operator=(const AtomSelThr &) { return *this; };
00069 AtomSelThr(AtomSelThr &) : AtomSel(NULL,-1) {};
00070 int change(const char *newcmd, DrawMolecule *mol) { return NO_PARSE; }
00071
00072 public:
00073
00074
00075 void update( DrawMolecule *mol, const int frame) {
00076 if (!sel) return;
00077
00078 wkf_mutex_lock(lock);
00079
00080 sel->which_frame=frame;
00081 which_frame=frame;
00082
00083 if (sel->change(NULL, mol) != AtomSel::PARSE_SUCCESS)
00084 msgErr << "AtomSelThr::update(): failed to evaluate atom selection update";
00085
00086 num_atoms=sel->num_atoms;
00087 selected=sel->selected;
00088 if (!on) on = new int[num_atoms];
00089 memcpy(on,sel->on,num_atoms*sizeof(int));
00090
00091 wkf_mutex_unlock(lock);
00092 }
00093
00094 protected:
00095 AtomSel *sel;
00096 wkf_mutex_t *lock;
00097 };
00098
00099
00100
00101
00102
00103
00104
00105
00106
00107
00108
00109
00110
00111
00112
00113
00114
00115
00116
00117
00118
00119
00120
00121
00122
00123
00124
00125
00126
00127
00128
00129
00130
00131
00132
00133
00134
00135
00136
00137
00138
00139
00140
00141
00142
00143
00144
00145
00146
00147
00148
00149
00150
00151
00152
00153
00154
00155
00156
00157
00158
00159
00160
00161
00162
00163
00164
00165
00166
00167
00168
00169
00170
00171
00172
00173
00174
00175 #ifdef __cplusplus
00176 extern "C" {
00177 #endif
00178
00179 #ifdef R_SMALL
00180 #undef R_SMALL
00181 #endif
00182 #define R_SMALL 0.000000001
00183
00184 static void normalize3d(double *v) {
00185 double vlen;
00186 vlen = sqrt(v[0]*v[0] + v[1]*v[1] + v[2]*v[2]);
00187 if (vlen > R_SMALL) {
00188 v[0] /= vlen;
00189 v[1] /= vlen;
00190 v[2] /= vlen;
00191 } else {
00192 v[0] = 0;
00193 v[1] = 0;
00194 v[2] = 0;
00195 }
00196 }
00197
00198
00199 static float MyMatrixFitRMS(int n, float *v1, float *v2, const float *wt, const double tol)
00200 {
00201
00202
00203
00204
00205
00206
00207
00208
00209 float *vv1,*vv2;
00210 double m[3][3],aa[3][3];
00211 double sumwt, sig, gam;
00212 double sg, bb, cc, tmp, err, etmp;
00213 int a, b, c, maxiter, iters, ix, iy, iz;
00214 double t1[3],t2[3];
00215 double aatmp[9];
00216
00217
00218
00219 for(a=0;a<3;a++) {
00220 for(b=0;b<3;b++) {
00221 m[a][b] = 0.0;
00222 aa[a][b] = 0.0;
00223 aatmp[3*a+b] = 0;
00224 }
00225 m[a][a] = 1.0;
00226 t1[a]=0.0;
00227 t2[a]=0.0;
00228 }
00229
00230
00231 maxiter = 1000;
00232
00233
00234
00235 vv1=v1;
00236 vv2=v2;
00237 sumwt = 0.0;
00238
00239 for(c=0;c<n;c++) {
00240 double w = wt ? wt[c] : 1;
00241 t1[0] += w * vv1[0];
00242 t1[1] += w * vv1[1];
00243 t1[2] += w * vv1[2];
00244 t2[0] += w * vv2[0];
00245 t2[1] += w * vv2[1];
00246 t2[2] += w * vv2[2];
00247 sumwt += w;
00248 vv1+=3;
00249 vv2+=3;
00250 }
00251 for(a=0;a<3;a++) {
00252 t1[a] /= sumwt;
00253 t2[a] /= sumwt;
00254 }
00255
00256
00257 vv1=v1;
00258 vv2=v2;
00259 for(c=0;c<n;c++) {
00260 double w = wt ? wt[c] : 1;
00261 double x1 = w * (vv1[0] - t1[0]);
00262 double y1 = w * (vv1[1] - t1[1]);
00263 double z1 = w * (vv1[2] - t1[2]);
00264
00265
00266 double x2 = (vv2[0] - t2[0]);
00267 double y2 = (vv2[1] - t2[1]);
00268 double z2 = (vv2[2] - t2[2]);
00269 aatmp[0] += x2 * x1;
00270 aatmp[1] += x2 * y1;
00271 aatmp[2] += x2 * z1;
00272 aatmp[3] += y2 * x1;
00273 aatmp[4] += y2 * y1;
00274 aatmp[5] += y2 * z1;
00275 aatmp[6] += z2 * x1;
00276 aatmp[7] += z2 * y1;
00277 aatmp[8] += z2 * z1;
00278 vv1+=3;
00279 vv2+=3;
00280 }
00281
00282 for (a=0; a<3; a++)
00283 for (b=0; b<3; b++)
00284 aa[a][b] = aatmp[3*a+b];
00285
00286 if(n>1) {
00287
00288 iters = 0;
00289 while(1) {
00290
00291 iz = (iters+1) % 3;
00292 iy = (iz+1) % 3;
00293 ix = (iy+1) % 3;
00294 sig = aa[iz][iy] - aa[iy][iz];
00295 gam = aa[iy][iy] + aa[iz][iz];
00296
00297 if(iters>=maxiter) {
00298 fprintf(stderr,
00299 " Matrix: Warning: no convergence (%.15f>%.15f after %d iterations).\n",
00300 fabs(sig),tol*fabs(gam),iters);
00301 break;
00302 }
00303
00304
00305
00306 tmp = sig*sig + gam*gam;
00307 sg = sqrt(tmp);
00308 if( (sg > 0.0) && (fabs(sig)>(tol*fabs(gam))) ) {
00309 sg = 1.0 / sg;
00310 for(a=0;a<3;a++) {
00311 bb = gam*aa[iy][a] + sig*aa[iz][a];
00312 cc = gam*aa[iz][a] - sig*aa[iy][a];
00313 aa[iy][a] = bb*sg;
00314 aa[iz][a] = cc*sg;
00315
00316 bb = gam*m[iy][a] + sig*m[iz][a];
00317 cc = gam*m[iz][a] - sig*m[iy][a];
00318 m[iy][a] = bb*sg;
00319 m[iz][a] = cc*sg;
00320 }
00321 } else {
00322 break;
00323 }
00324 iters++;
00325 }
00326 }
00327
00328
00329 err=0.0;
00330 vv1=v1;
00331 vv2=v2;
00332
00333 normalize3d(m[0]);
00334 normalize3d(m[1]);
00335 normalize3d(m[2]);
00336
00337 for(c=0;c<n;c++) {
00338 etmp = 0.0;
00339 for(a=0;a<3;a++) {
00340 tmp = m[a][0]*(vv2[0]-t2[0])
00341 + m[a][1]*(vv2[1]-t2[1])
00342 + m[a][2]*(vv2[2]-t2[2]);
00343 tmp = (vv1[a]-t1[a])-tmp;
00344 etmp += tmp*tmp;
00345 }
00346
00347 if(wt)
00348 err += wt[c] * etmp;
00349 else
00350 err += etmp;
00351
00352 vv1+=3;
00353 vv2+=3;
00354 }
00355
00356 err=err/sumwt;
00357 err=sqrt(err);
00358 return (float)err;
00359 }
00360
00361 #ifdef __cplusplus
00362 }
00363 #endif
00364
00365
00366
00367
00368
00369
00370 static float cluster_get_rmsd(const float *Frame1Pos, const float *Frame2Pos,
00371 AtomSel *sel, float *weights) {
00372 float distance = 0.0f;
00373 measure_rmsd(sel, sel, sel->num_atoms, Frame1Pos, Frame2Pos, weights, &distance);
00374 return distance;
00375 }
00376
00377
00378
00379
00380 static float cluster_get_rgyrd(const float *Frame1Pos, const float *Frame2Pos,
00381 AtomSel *sel, float *weights) {
00382
00383 float distance = 10000000.0f;
00384
00385
00386 float com1[3], com2[3];
00387 int ret_val;
00388
00389 ret_val = measure_center(sel, Frame1Pos, weights, com1);
00390 if (ret_val < 0)
00391 return distance;
00392
00393 ret_val = measure_center(sel, Frame2Pos, weights, com2);
00394 if (ret_val < 0)
00395 return distance;
00396
00397
00398 int i, j;
00399 float total_w, w, sum1, sum2;
00400 total_w=sum1=sum2=0.0f;
00401 for (j=0,i=sel->firstsel; i<=sel->lastsel; i++) {
00402 if (sel->on[i]) {
00403 w = weights[j];
00404 total_w += w;
00405 sum1 += w * distance2(Frame1Pos + 3*i, com1);
00406 sum2 += w * distance2(Frame2Pos + 3*i, com2);
00407 j++;
00408 }
00409 }
00410
00411 if (total_w == 0.0f)
00412 return distance;
00413
00414
00415 distance = sqrtf(sum1/total_w) - sqrtf(sum2/total_w);
00416 return fabsf(distance);
00417 }
00418
00419
00420
00421
00422
00423
00424 static float cluster_get_fitrmsd(const float *Frame1Pos, const float *Frame2Pos,
00425 AtomSel *sel, float *weights, const double tol) {
00426
00427 int num = sel->selected;
00428
00429
00430 float distance = 10000000.0f;
00431
00432
00433
00434
00435
00436
00437
00438 if (sel->selected < 4)
00439 return distance;
00440
00441 int i, j, k;
00442 float *v1, *v2, *wt;
00443 v1 = new float[3*num];
00444 v2 = new float[3*num];
00445 wt = new float[num];
00446 for (j=0,k=0,i=sel->firstsel; i<=sel->lastsel; i++) {
00447 if (sel->on[i]) {
00448 int ind = 3 * i;
00449 wt[j] = weights[i];
00450 ++j;
00451 v1[k] = Frame1Pos[ind];
00452 v2[k] = Frame2Pos[ind];
00453 v1[k+1] = Frame1Pos[ind+1];
00454 v2[k+1] = Frame2Pos[ind+1];
00455 v1[k+2] = Frame1Pos[ind+2];
00456 v2[k+2] = Frame2Pos[ind+2];
00457 k+=3;
00458 }
00459 }
00460 distance = MyMatrixFitRMS(num, v1, v2, wt, tol);
00461
00462 delete [] v1;
00463 delete [] v2;
00464 delete [] wt;
00465
00466 return distance;
00467 }
00468
00469
00470
00471
00472 typedef struct {
00473 int threadid;
00474 int threadcount;
00475
00476 int max_cluster_size;
00477 const int *skip_list;
00478 int *new_skip_list;
00479 int *max_cluster;
00480
00481 int istart;
00482 int iend;
00483 int *frames_list;
00484 int numframes;
00485
00486 AtomSelThr *sel;
00487 Molecule *mol;
00488 int selupdate;
00489 float cutoff;
00490 int likeness;
00491 float *weights;
00492 } clusterparms_t;
00493
00494
00495
00496 extern "C" void * find_cluster_thr(void *voidparms)
00497 {
00498
00499 clusterparms_t *parms = (clusterparms_t *)voidparms;
00500 const int istart = parms->istart;
00501 const int iend = parms->iend;
00502 int *framesList = parms->frames_list;
00503 const int numframes = parms->numframes;
00504
00505 const int selupdate = parms->selupdate;
00506 const int likeness = parms->likeness;
00507 float cutoff = parms->cutoff;
00508 float *weights = parms->weights;
00509
00510 AtomSelThr *sel = parms->sel;
00511 Molecule *mymol = parms->mol;
00512 const int *skipList = parms->skip_list;
00513
00514 int *maxCluster = parms->max_cluster;
00515 memset(maxCluster, 0, numframes*sizeof(int));
00516 int *newSkipList = parms->new_skip_list;
00517 memset(newSkipList, 0, numframes*sizeof(int));
00518
00519 int maxClusterSize = 0, tempClusterSize = 0;
00520 int *tempCluster = new int[numframes];
00521 int *tempSkipList = new int[numframes];
00522
00523
00524
00525
00526 double tol = 1e-15;
00527 const char *TOL = getenv( "VMDFITRMSTOLERANCE" );
00528 if (TOL)
00529 tol = atof(TOL);
00530
00531
00532 int i,j;
00533 for (i = istart; i < iend; i++) {
00534 memset(tempSkipList, 0, numframes*sizeof(int));
00535 memset(tempCluster, 0, numframes*sizeof(int));
00536
00537 if (skipList[i]==0) {
00538 if (selupdate)
00539 sel->update(mymol,framesList[i]);
00540
00541 const Timestep *tsMain = mymol->get_frame(framesList[i]);
00542 const float *framePos = tsMain->pos;
00543
00544 tempCluster[0] = i;
00545 tempSkipList[i] = 1;
00546 tempClusterSize = 1;
00547
00548
00549 for (j = 0; j < numframes; j++) {
00550 if (skipList[j]==0 && j != i) {
00551 const Timestep *ts2;
00552 ts2 = mymol->get_frame(framesList[j]);
00553 float distance;
00554
00555
00556 switch(likeness) {
00557 case MEASURE_DIST_RMSD:
00558 distance = cluster_get_rmsd(framePos, ts2->pos, sel, weights);
00559 break;
00560
00561 case MEASURE_DIST_FITRMSD:
00562 distance = cluster_get_fitrmsd(framePos, ts2->pos, sel, weights, tol);
00563 break;
00564
00565 case MEASURE_DIST_RGYRD:
00566 distance = cluster_get_rgyrd(framePos, ts2->pos, sel, weights);
00567 break;
00568
00569 default:
00570 distance = 10000000.0f;
00571 }
00572
00573 if (distance <= cutoff) {
00574 tempCluster[tempClusterSize] = j;
00575 ++tempClusterSize;
00576 tempSkipList[j] = 1;
00577 }
00578 }
00579 }
00580
00581
00582 if (tempClusterSize > maxClusterSize) {
00583 int *temp;
00584 maxClusterSize = tempClusterSize;
00585
00586 temp = maxCluster;
00587 maxCluster = tempCluster;
00588 tempCluster = temp;
00589
00590 temp = newSkipList;
00591 newSkipList = tempSkipList;
00592 tempSkipList = temp;
00593 }
00594 }
00595 }
00596
00597
00598 parms->max_cluster_size = maxClusterSize;
00599 parms->max_cluster = maxCluster;
00600 parms->new_skip_list = newSkipList;
00601
00602
00603 delete[] tempCluster;
00604 delete[] tempSkipList;
00605
00606 return MEASURE_NOERR;
00607 }
00608
00609
00611 static int *find_next_cluster(Molecule *mymol, int *framesList, const int numframes,
00612 const int remframes, const int *skipList, int **newSkipList,
00613 const int likeness, AtomSel *sel, const int selupdate,
00614 const double cutoff, float *weights)
00615 {
00616 int i,j;
00617
00618
00619 wkf_thread_t *threads;
00620 clusterparms_t *parms;
00621
00622 #if defined(VMDTHREADS)
00623 int numprocs = wkf_thread_numprocessors();
00624 #else
00625 int numprocs = 1;
00626 #endif
00627
00628 int delta = remframes / numprocs;
00629 int istart = 0;
00630 int iend = 0;
00631
00632
00633 if (delta < 1) {
00634 numprocs=1;
00635 delta=numframes;
00636 }
00637
00638 threads = new wkf_thread_t[numprocs];
00639 memset(threads, 0, numprocs * sizeof(wkf_thread_t));
00640 wkf_mutex_t *atomsel_lock = new wkf_mutex_t;
00641 wkf_mutex_init(atomsel_lock);
00642
00643
00644 parms = new clusterparms_t[numprocs];
00645 for (i=0; i<numprocs; ++i) {
00646 parms[i].threadid = i;
00647 parms[i].threadcount = numprocs;
00648
00649 parms[i].max_cluster_size = 1;
00650 parms[i].skip_list = skipList;
00651 parms[i].new_skip_list = new int[numframes];
00652 parms[i].max_cluster = new int[numframes];
00653
00654
00655
00656
00657
00658
00659 parms[i].sel = new AtomSelThr(sel,atomsel_lock);
00660 parms[i].mol = mymol;
00661
00662
00663
00664 parms[i].istart = istart;
00665 int nframe=0;
00666 for (j=istart; (j < numframes) && (nframe < delta); ++j) {
00667 if (skipList[framesList[j]]==0)
00668 ++nframe;
00669 iend=j;
00670 }
00671 parms[i].iend = iend;
00672 istart=iend;
00673
00674 parms[i].frames_list = framesList;
00675 parms[i].numframes = numframes;
00676 parms[i].selupdate = selupdate;
00677 parms[i].likeness = likeness;
00678 parms[i].cutoff = (float) cutoff;
00679 parms[i].weights = weights;
00680 }
00681 parms[numprocs-1].iend=numframes;
00682
00683
00684 #if defined(VMDTHREADS)
00685 if (numprocs > 1) {
00686 for (i=0; i<numprocs; ++i) {
00687 wkf_thread_create(&threads[i], find_cluster_thr, &parms[i]);
00688 }
00689 for (i=0; i<numprocs; ++i) {
00690 wkf_thread_join(threads[i], NULL);
00691 }
00692 } else
00693 #endif
00694 find_cluster_thr(&parms[0]);
00695
00696 int maxClusterSize = parms[0].max_cluster_size;
00697 int *maxCluster = parms[0].max_cluster;
00698 delete[] *newSkipList;
00699 *newSkipList= parms[0].new_skip_list;
00700
00701
00702
00703 if (numprocs > 1) {
00704 for (i = 1; i < numprocs; i++) {
00705 if (parms[i].max_cluster_size > maxClusterSize) {
00706 maxClusterSize = parms[i].max_cluster_size;
00707 delete[] maxCluster;
00708 maxCluster = parms[i].max_cluster;
00709 delete[] *newSkipList;
00710 *newSkipList = parms[i].new_skip_list;
00711 } else {
00712 delete[] parms[i].max_cluster;
00713 delete[] parms[i].new_skip_list;
00714 }
00715 }
00716 }
00717
00718
00719 for (i = 0; i < numframes; i++) {
00720 maxCluster[i] = framesList[maxCluster[i]];
00721 }
00722
00723
00724 wkf_mutex_destroy(atomsel_lock);
00725 delete atomsel_lock;
00726
00727 if (selupdate) {
00728 for (i=0; i<numprocs; ++i)
00729 delete parms[i].sel;
00730 }
00731 delete[] threads;
00732 delete[] parms;
00733
00734 return maxCluster;
00735 }
00736
00737 int measure_cluster(AtomSel *sel, MoleculeList *mlist,
00738 const int numcluster, const int algorithm,
00739 const int likeness, const double cutoff,
00740 int *clustersize, int **clusterlist,
00741 int first, int last, int step, int selupdate,
00742 float *weights)
00743 {
00744 Molecule *mymol = mlist->mol_from_id(sel->molid());
00745 int maxframe = mymol->numframes()-1;
00746
00747 if (last == -1) last = maxframe;
00748
00749 if ((last < first) || (last < 0) || (step <=0) || (first < 0)
00750 || (last > maxframe)) {
00751 msgErr << "measure cluster: bad frame range given."
00752 << " max. allowed frame#: " << maxframe << sendmsg;
00753 return MEASURE_ERR_BADFRAMERANGE;
00754 }
00755
00756 int numframes = (last-first+1)/step;
00757 int remframes = numframes;
00758
00759
00760 int *framesList = new int[numframes];
00761 int frame_count = 0;
00762 int n;
00763
00764 for(n = first; n <= last; n += step)
00765 framesList[frame_count++] = n;
00766
00767
00768 int *skipList = new int[numframes];
00769
00770 int *newSkipList = new int[numframes];
00771
00772 memset(skipList, 0, numframes*sizeof(int));
00773
00774 wkfmsgtimer *msgtp = wkf_msg_timer_create(5);
00775
00776
00777 for(n = 0; n < numcluster; ++n){
00778
00779
00780 memset(newSkipList, 0, numframes*sizeof(int));
00781 clusterlist[n] = find_next_cluster(mymol, framesList, numframes, remframes,
00782 skipList, &newSkipList, likeness,
00783 sel, selupdate, cutoff, weights);
00784 int n_cluster=0;
00785 for(int i = 0; i < numframes; ++i){
00786 if (newSkipList[i] == 1) {
00787 skipList[i] = 1;
00788 n_cluster++;
00789 }
00790 }
00791 clustersize[n]=n_cluster;
00792 remframes -= n_cluster;
00793
00794
00795 if (msgtp && wkf_msg_timer_timeout(msgtp)) {
00796 char tmpbuf[1024];
00797 sprintf(tmpbuf, "cluster %d of %d: (%6.2f%% complete). %d frames of %d left.",
00798 n+1, numcluster, 100.0f*(n+1)/((float) numcluster), remframes, numframes);
00799 msgInfo << "measure cluster: " << tmpbuf << sendmsg;
00800 }
00801 }
00802
00803
00804 int *unclustered = new int[numframes];
00805 int numunclustered = 0;
00806 for (n = 0; n < numframes; ++n) {
00807 if (skipList[n] == 0) {
00808 unclustered[numunclustered] = framesList[n];
00809 ++numunclustered;
00810 }
00811 }
00812
00813
00814 clusterlist[numcluster] = unclustered;
00815 clustersize[numcluster] = numunclustered;
00816
00817
00818 delete[] newSkipList;
00819 delete[] skipList;
00820 wkf_msg_timer_destroy(msgtp);
00821
00822 return MEASURE_NOERR;
00823 }
00824
00825
00826
00827
00828 typedef ResizeArray<int> intlist;
00829
00830
00831
00832 static void assemble_cluster(intlist &cluster_list, intlist &candidate_list,
00833 intlist **neighbor_grid, int atom, int numshared, int *idxmap) {
00834 int idx, nn, i,j;
00835
00836
00837 candidate_list.clear();
00838 cluster_list.clear();
00839
00840 idx = idxmap[atom];
00841 nn = neighbor_grid[idx]->num();
00842 for (i = 0; i < nn; i++) {
00843 int bn = (*neighbor_grid[idx])[i];
00844 if (neighbor_grid[idxmap[bn]]) {
00845 candidate_list.append(atom);
00846 candidate_list.append(bn);
00847 }
00848 }
00849
00850
00851 int curidx=0;
00852
00853 while (curidx < candidate_list.num()) {
00854
00855
00856
00857
00858 int count = 0;
00859
00860 if (numshared > 0) {
00861 int ida = idxmap[candidate_list[curidx]];
00862 int idb = idxmap[candidate_list[curidx+1]];
00863 int nna = neighbor_grid[ida]->num();
00864 int nnb = neighbor_grid[idb]->num();
00865
00866 for (i = 0; i < nna; i++) {
00867 if (neighbor_grid[ida]) {
00868 for (j = 0; j < nnb; j++) {
00869 if (neighbor_grid[idb]) {
00870 if ( (*neighbor_grid[ida])[i] == (*neighbor_grid[idb])[j] ) {
00871 ++count;
00872 if (count == numshared)
00873 goto exit;
00874 }
00875 }
00876 }
00877 }
00878 }
00879 }
00880 exit:
00881
00882 if (count == numshared) {
00883
00884
00885
00886
00887 int atma = candidate_list[curidx];
00888 if (cluster_list.find(atma) < 0)
00889 cluster_list.append(atma);
00890
00891
00892
00893
00894 int atmb = candidate_list[curidx+1];
00895 idx = idxmap[atmb];
00896
00897 if (cluster_list.find(atmb) < 0) {
00898 cluster_list.append(atmb);
00899
00900 int nnb = neighbor_grid[idx]->num();
00901 for (i = 0; i < nnb; i++) {
00902 int bn = (*neighbor_grid[idx])[i];
00903 if ((neighbor_grid[idxmap[bn]]) && (cluster_list.find(bn) < 0)) {
00904 candidate_list.append(atmb);
00905 candidate_list.append(bn);
00906 }
00907 }
00908 }
00909 }
00910
00911 ++curidx;++curidx;
00912 }
00913
00914 return;
00915 }
00916
00917
00918 int measure_clustsize(const AtomSel *sel, MoleculeList *mlist,
00919 const double cutoff, int *clustersize,
00920 int *clusternum, int *clusteridx,
00921 int minsize, int numshared, int usepbc) {
00922 int i,j;
00923
00924 const float *framepos = sel->coordinates(mlist);
00925 const int num_selected = sel->selected;
00926 const int num_atoms = sel->num_atoms;
00927 const int *selected = sel->on;
00928
00929
00930
00931 int *idxmap = new int[num_atoms];
00932 int *idxrev = new int[num_selected];
00933
00934 for (i=0; i<sel->firstsel; i++)
00935 idxmap[i]=-1;
00936
00937 for(j=0,i=sel->firstsel; i<=sel->lastsel; i++) {
00938 if (sel->on[i]) {
00939 idxrev[j]=i;
00940 idxmap[i]=j++;
00941 } else {
00942 idxmap[i]=-1;
00943 }
00944 }
00945
00946 for (i=sel->lastsel+1; i<sel->num_atoms; i++)
00947 idxmap[i]=-1;
00948
00949
00950 intlist **neighbor_grid;
00951 neighbor_grid = new intlist *[num_selected];
00952 for (i = 0; i < num_selected; i++)
00953 neighbor_grid[i] = new intlist;
00954
00955
00956 GridSearchPair *pairlist, *currentPair, *nextPair;
00957 pairlist = vmd_gridsearch1(framepos, num_atoms, selected, (float)cutoff, 0, -1);
00958
00959
00960 for (currentPair = pairlist; currentPair != NULL; currentPair = nextPair) {
00961 neighbor_grid[idxmap[currentPair->ind1]]->append(currentPair->ind2);
00962 neighbor_grid[idxmap[currentPair->ind2]]->append(currentPair->ind1);
00963 nextPair = currentPair->next;
00964 free(currentPair);
00965 }
00966
00967
00968 int currentClusterNum = 0;
00969 int currentPosition = 0;
00970 intlist cluster_list(64);
00971 intlist candidate_list(128);
00972
00973 while (currentPosition < num_selected) {
00974 if (neighbor_grid[currentPosition]) {
00975
00976
00977
00978
00979
00980
00981 if (neighbor_grid[currentPosition]->num() > numshared) {
00982 assemble_cluster(cluster_list, candidate_list, neighbor_grid,
00983 idxrev[currentPosition], numshared, idxmap);
00984
00985 if (minsize <= cluster_list.num()) {
00986
00987 for (i = 0; i < cluster_list.num(); i++) {
00988 int idx = idxmap[cluster_list[i]];
00989 delete neighbor_grid[idx];
00990 neighbor_grid[idx] = 0;
00991 }
00992
00993
00994
00995 for (i = 0; i < cluster_list.num(); i++) {
00996 int atom = idxmap[cluster_list[i]];
00997 clusteridx[atom] = cluster_list[i];
00998 clusternum[atom] = currentClusterNum;
00999 clustersize[atom] = cluster_list.num();
01000 }
01001 currentClusterNum++;
01002 }
01003 }
01004 }
01005 ++currentPosition;
01006 }
01007
01008 for(i=0; i < num_selected; ++i) {
01009 if (neighbor_grid[i])
01010 delete neighbor_grid[i];
01011 }
01012 delete[] neighbor_grid;
01013
01014 return MEASURE_NOERR;
01015 }