Revised wording in pdb2gmx.c, hopefully clearer now.
[gromacs/rigid-bodies.git] / src / gmxlib / orires.c
blobd57d9ac7e85c3f72bb21549a300c1cfe077d74ec
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
3 *
4 * This source code is part of
5 *
6 * G R O M A C S
7 *
8 * GROningen MAchine for Chemical Simulations
9 *
10 * VERSION 3.2.0
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
33 * And Hey:
34 * GROningen Mixture of Alchemy and Childrens' Stories
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
40 #include "typedefs.h"
41 #include "smalloc.h"
42 #include "vec.h"
43 #include "nrjac.h"
44 #include "network.h"
45 #include "orires.h"
46 #include "do_fit.h"
47 #include "main.h"
48 #include "copyrite.h"
49 #include "pbc.h"
50 #include "mtop_util.h"
52 void init_orires(FILE *fplog,const gmx_mtop_t *mtop,
53 rvec xref[],
54 const t_inputrec *ir,
55 const gmx_multisim_t *ms,t_oriresdata *od,
56 t_state *state)
58 int i,j,d,ex,nmol,nr,*nr_ex;
59 double mtot;
60 rvec com;
61 gmx_mtop_ilistloop_t iloop;
62 t_ilist *il;
63 gmx_mtop_atomloop_all_t aloop;
64 t_atom *atom;
66 od->fc = ir->orires_fc;
67 od->nex = 0;
68 od->S = NULL;
70 od->M=NULL;
71 od->eig=NULL;
72 od->v=NULL;
74 od->nr = gmx_mtop_ftype_count(mtop,F_ORIRES);
75 if (od->nr == 0)
77 return;
80 nr_ex = NULL;
82 iloop = gmx_mtop_ilistloop_init(mtop);
83 while (gmx_mtop_ilistloop_next(iloop,&il,&nmol))
85 for(i=0; i<il[F_ORIRES].nr; i+=3)
87 ex = mtop->ffparams.iparams[il[F_ORIRES].iatoms[i]].orires.ex;
88 if (ex >= od->nex)
90 srenew(nr_ex,ex+1);
91 for(j=od->nex; j<ex+1; j++)
93 nr_ex[j] = 0;
95 od->nex = ex+1;
97 nr_ex[ex]++;
100 snew(od->S,od->nex);
101 /* When not doing time averaging, the instaneous and time averaged data
102 * are indentical and the pointers can point to the same memory.
104 snew(od->Dinsl,od->nr);
105 if (ms)
107 snew(od->Dins,od->nr);
109 else
111 od->Dins = od->Dinsl;
114 if (ir->orires_tau == 0)
116 od->Dtav = od->Dins;
117 od->edt = 0.0;
118 od->edt_1= 1.0;
120 else
122 snew(od->Dtav,od->nr);
123 od->edt = exp(-ir->delta_t/ir->orires_tau);
124 od->edt_1= 1.0 - od->edt;
126 /* Extend the state with the orires history */
127 state->flags |= (1<<estORIRE_INITF);
128 state->hist.orire_initf = 1;
129 state->flags |= (1<<estORIRE_DTAV);
130 state->hist.norire_Dtav = od->nr*5;
131 snew(state->hist.orire_Dtav,state->hist.norire_Dtav);
134 snew(od->oinsl,od->nr);
135 if (ms)
137 snew(od->oins,od->nr);
139 else
141 od->oins = od->oinsl;
143 if (ir->orires_tau == 0) {
144 od->otav = od->oins;
146 else
148 snew(od->otav,od->nr);
150 snew(od->tmp,od->nex);
151 snew(od->TMP,od->nex);
152 for(ex=0; ex<od->nex; ex++)
154 snew(od->TMP[ex],5);
155 for(i=0; i<5; i++)
157 snew(od->TMP[ex][i],5);
161 od->nref = 0;
162 for(i=0; i<mtop->natoms; i++)
164 if (ggrpnr(&mtop->groups,egcORFIT,i) == 0)
166 od->nref++;
169 snew(od->mref,od->nref);
170 snew(od->xref,od->nref);
171 snew(od->xtmp,od->nref);
173 snew(od->eig,od->nex*12);
175 /* Determine the reference structure on the master node.
176 * Copy it to the other nodes after checking multi compatibility,
177 * so we are sure the subsystems match before copying.
179 clear_rvec(com);
180 mtot = 0.0;
181 j = 0;
182 aloop = gmx_mtop_atomloop_all_init(mtop);
183 while(gmx_mtop_atomloop_all_next(aloop,&i,&atom))
185 if (mtop->groups.grpnr[egcORFIT] == NULL ||
186 mtop->groups.grpnr[egcORFIT][i] == 0)
188 /* Not correct for free-energy with changing masses */
189 od->mref[j] = atom->m;
190 if (ms==NULL || MASTERSIM(ms))
192 copy_rvec(xref[i],od->xref[j]);
193 for(d=0; d<DIM; d++)
195 com[d] += od->mref[j]*xref[i][d];
198 mtot += od->mref[j];
199 j++;
202 svmul(1.0/mtot,com,com);
203 if (ms==NULL || MASTERSIM(ms))
205 for(j=0; j<od->nref; j++)
207 rvec_dec(od->xref[j],com);
211 fprintf(fplog,"Found %d orientation experiments\n",od->nex);
212 for(i=0; i<od->nex; i++)
214 fprintf(fplog," experiment %d has %d restraints\n",i+1,nr_ex[i]);
217 sfree(nr_ex);
219 fprintf(fplog," the fit group consists of %d atoms and has total mass %g\n",
220 od->nref,mtot);
222 if (ms)
224 fprintf(fplog," the orientation restraints are ensemble averaged over %d systems\n",ms->nsim);
226 check_multi_int(fplog,ms,od->nr,
227 "the number of orientation restraints");
228 check_multi_int(fplog,ms,od->nref,
229 "the number of fit atoms for orientation restraining");
230 check_multi_int(fplog,ms,ir->nsteps,"nsteps");
231 /* Copy the reference coordinates from the master to the other nodes */
232 gmx_sum_sim(DIM*od->nref,od->xref[0],ms);
235 please_cite(fplog,"Hess2003");
238 void diagonalize_orires_tensors(t_oriresdata *od)
240 int ex,i,j,nrot,ord[DIM],t;
241 matrix S,TMP;
243 if (od->M == NULL)
245 snew(od->M,DIM);
246 for(i=0; i<DIM; i++)
248 snew(od->M[i],DIM);
250 snew(od->eig_diag,DIM);
251 snew(od->v,DIM);
252 for(i=0; i<DIM; i++)
254 snew(od->v[i],DIM);
258 for(ex=0; ex<od->nex; ex++)
260 /* Rotate the S tensor back to the reference frame */
261 mmul(od->R,od->S[ex],TMP);
262 mtmul(TMP,od->R,S);
263 for(i=0; i<DIM; i++)
265 for(j=0; j<DIM; j++)
267 od->M[i][j] = S[i][j];
271 jacobi(od->M,DIM,od->eig_diag,od->v,&nrot);
273 for(i=0; i<DIM; i++)
275 ord[i] = i;
277 for(i=0; i<DIM; i++)
279 for(j=i+1; j<DIM; j++)
281 if (sqr(od->eig_diag[ord[j]]) > sqr(od->eig_diag[ord[i]]))
283 t = ord[i];
284 ord[i] = ord[j];
285 ord[j] = t;
290 for(i=0; i<DIM; i++)
292 od->eig[ex*12 + i] = od->eig_diag[ord[i]];
294 for(i=0; i<DIM; i++)
296 for(j=0; j<DIM; j++)
298 od->eig[ex*12 + 3 + 3*i + j] = od->v[j][ord[i]];
304 void print_orires_log(FILE *log,t_oriresdata *od)
306 int ex,i;
307 real *eig;
309 diagonalize_orires_tensors(od);
311 for(ex=0; ex<od->nex; ex++)
313 eig = od->eig + ex*12;
314 fprintf(log," Orientation experiment %d:\n",ex+1);
315 fprintf(log," order parameter: %g\n",eig[0]);
316 for(i=0; i<DIM; i++)
318 fprintf(log," eig: %6.3f %6.3f %6.3f %6.3f\n",
319 (eig[0] != 0) ? eig[i]/eig[0] : eig[i],
320 eig[DIM+i*DIM+XX],
321 eig[DIM+i*DIM+YY],
322 eig[DIM+i*DIM+ZZ]);
324 fprintf(log,"\n");
328 real calc_orires_dev(const gmx_multisim_t *ms,
329 int nfa,const t_iatom forceatoms[],const t_iparams ip[],
330 const t_mdatoms *md,const rvec x[],const t_pbc *pbc,
331 t_fcdata *fcd,history_t *hist)
333 int fa,d,i,j,type,ex,nref;
334 real edt,edt_1,invn,pfac,r2,invr,corrfac,weight,wsv2,sw,dev;
335 tensor *S,R,TMP;
336 rvec5 *Dinsl,*Dins,*Dtav,*rhs;
337 real *mref,***T;
338 double mtot;
339 rvec *xref,*xtmp,com,r_unrot,r;
340 t_oriresdata *od;
341 gmx_bool bTAV;
342 const real two_thr=2.0/3.0;
344 od = &(fcd->orires);
346 if (od->nr == 0)
348 /* This means that this is not the master node */
349 gmx_fatal(FARGS,"Orientation restraints are only supported on the master node, use less processors");
352 bTAV = (od->edt != 0);
353 edt = od->edt;
354 edt_1= od->edt_1;
355 S = od->S;
356 Dinsl= od->Dinsl;
357 Dins = od->Dins;
358 Dtav = od->Dtav;
359 T = od->TMP;
360 rhs = od->tmp;
361 nref = od->nref;
362 mref = od->mref;
363 xref = od->xref;
364 xtmp = od->xtmp;
366 if (bTAV)
368 od->exp_min_t_tau = hist->orire_initf*edt;
370 /* Correction factor to correct for the lack of history
371 * at short times.
373 corrfac = 1.0/(1.0 - od->exp_min_t_tau);
375 else
377 corrfac = 1.0;
380 if (ms)
382 invn = 1.0/ms->nsim;
384 else
386 invn = 1.0;
389 clear_rvec(com);
390 mtot = 0;
391 j=0;
392 for(i=0; i<md->nr; i++)
394 if (md->cORF[i] == 0)
396 copy_rvec(x[i],xtmp[j]);
397 mref[j] = md->massT[i];
398 for(d=0; d<DIM; d++)
400 com[d] += mref[j]*xref[j][d];
402 mtot += mref[j];
403 j++;
406 svmul(1.0/mtot,com,com);
407 for(j=0; j<nref; j++)
409 rvec_dec(xtmp[j],com);
411 /* Calculate the rotation matrix to rotate x to the reference orientation */
412 calc_fit_R(DIM,nref,mref,xref,xtmp,R);
413 copy_mat(R,od->R);
415 d = 0;
416 for(fa=0; fa<nfa; fa+=3)
418 type = forceatoms[fa];
419 if (pbc)
421 pbc_dx_aiuc(pbc,x[forceatoms[fa+1]],x[forceatoms[fa+2]],r_unrot);
423 else
425 rvec_sub(x[forceatoms[fa+1]],x[forceatoms[fa+2]],r_unrot);
427 mvmul(R,r_unrot,r);
428 r2 = norm2(r);
429 invr = gmx_invsqrt(r2);
430 /* Calculate the prefactor for the D tensor, this includes the factor 3! */
431 pfac = ip[type].orires.c*invr*invr*3;
432 for(i=0; i<ip[type].orires.power; i++)
434 pfac *= invr;
436 Dinsl[d][0] = pfac*(2*r[0]*r[0] + r[1]*r[1] - r2);
437 Dinsl[d][1] = pfac*(2*r[0]*r[1]);
438 Dinsl[d][2] = pfac*(2*r[0]*r[2]);
439 Dinsl[d][3] = pfac*(2*r[1]*r[1] + r[0]*r[0] - r2);
440 Dinsl[d][4] = pfac*(2*r[1]*r[2]);
442 if (ms)
444 for(i=0; i<5; i++)
446 Dins[d][i] = Dinsl[d][i]*invn;
450 d++;
453 if (ms)
455 gmx_sum_sim(5*od->nr,Dins[0],ms);
458 /* Calculate the order tensor S for each experiment via optimization */
459 for(ex=0; ex<od->nex; ex++)
461 for(i=0; i<5; i++)
463 rhs[ex][i] = 0;
464 for(j=0; j<=i; j++)
466 T[ex][i][j] = 0;
470 d = 0;
471 for(fa=0; fa<nfa; fa+=3)
473 if (bTAV)
475 /* Here we update Dtav in t_fcdata using the data in history_t.
476 * Thus the results stay correct when this routine
477 * is called multiple times.
479 for(i=0; i<5; i++)
481 Dtav[d][i] = edt*hist->orire_Dtav[d*5+i] + edt_1*Dins[d][i];
485 type = forceatoms[fa];
486 ex = ip[type].orires.ex;
487 weight = ip[type].orires.kfac;
488 /* Calculate the vector rhs and half the matrix T for the 5 equations */
489 for(i=0; i<5; i++) {
490 rhs[ex][i] += Dtav[d][i]*ip[type].orires.obs*weight;
491 for(j=0; j<=i; j++)
493 T[ex][i][j] += Dtav[d][i]*Dtav[d][j]*weight;
496 d++;
498 /* Now we have all the data we can calculate S */
499 for(ex=0; ex<od->nex; ex++)
501 /* Correct corrfac and copy one half of T to the other half */
502 for(i=0; i<5; i++)
504 rhs[ex][i] *= corrfac;
505 T[ex][i][i] *= sqr(corrfac);
506 for(j=0; j<i; j++)
508 T[ex][i][j] *= sqr(corrfac);
509 T[ex][j][i] = T[ex][i][j];
512 m_inv_gen(T[ex],5,T[ex]);
513 /* Calculate the orientation tensor S for this experiment */
514 S[ex][0][0] = 0;
515 S[ex][0][1] = 0;
516 S[ex][0][2] = 0;
517 S[ex][1][1] = 0;
518 S[ex][1][2] = 0;
519 for(i=0; i<5; i++)
521 S[ex][0][0] += 1.5*T[ex][0][i]*rhs[ex][i];
522 S[ex][0][1] += 1.5*T[ex][1][i]*rhs[ex][i];
523 S[ex][0][2] += 1.5*T[ex][2][i]*rhs[ex][i];
524 S[ex][1][1] += 1.5*T[ex][3][i]*rhs[ex][i];
525 S[ex][1][2] += 1.5*T[ex][4][i]*rhs[ex][i];
527 S[ex][1][0] = S[ex][0][1];
528 S[ex][2][0] = S[ex][0][2];
529 S[ex][2][1] = S[ex][1][2];
530 S[ex][2][2] = -S[ex][0][0] - S[ex][1][1];
533 wsv2 = 0;
534 sw = 0;
536 d = 0;
537 for(fa=0; fa<nfa; fa+=3)
539 type = forceatoms[fa];
540 ex = ip[type].orires.ex;
542 od->otav[d] = two_thr*
543 corrfac*(S[ex][0][0]*Dtav[d][0] + S[ex][0][1]*Dtav[d][1] +
544 S[ex][0][2]*Dtav[d][2] + S[ex][1][1]*Dtav[d][3] +
545 S[ex][1][2]*Dtav[d][4]);
546 if (bTAV)
548 od->oins[d] = two_thr*(S[ex][0][0]*Dins[d][0] + S[ex][0][1]*Dins[d][1] +
549 S[ex][0][2]*Dins[d][2] + S[ex][1][1]*Dins[d][3] +
550 S[ex][1][2]*Dins[d][4]);
552 if (ms)
554 /* When ensemble averaging is used recalculate the local orientation
555 * for output to the energy file.
557 od->oinsl[d] = two_thr*
558 (S[ex][0][0]*Dinsl[d][0] + S[ex][0][1]*Dinsl[d][1] +
559 S[ex][0][2]*Dinsl[d][2] + S[ex][1][1]*Dinsl[d][3] +
560 S[ex][1][2]*Dinsl[d][4]);
563 dev = od->otav[d] - ip[type].orires.obs;
565 wsv2 += ip[type].orires.kfac*sqr(dev);
566 sw += ip[type].orires.kfac;
568 d++;
570 od->rmsdev = sqrt(wsv2/sw);
572 /* Rotate the S matrices back, so we get the correct grad(tr(S D)) */
573 for(ex=0; ex<od->nex; ex++)
575 tmmul(R,S[ex],TMP);
576 mmul(TMP,R,S[ex]);
579 return od->rmsdev;
581 /* Approx. 120*nfa/3 flops */
584 real orires(int nfa,const t_iatom forceatoms[],const t_iparams ip[],
585 const rvec x[],rvec f[],rvec fshift[],
586 const t_pbc *pbc,const t_graph *g,
587 real lambda,real *dvdlambda,
588 const t_mdatoms *md,t_fcdata *fcd,
589 int *global_atom_index)
591 atom_id ai,aj;
592 int fa,d,i,type,ex,power,ki=CENTRAL;
593 ivec dt;
594 real r2,invr,invr2,fc,smooth_fc,dev,devins,pfac;
595 rvec r,Sr,fij;
596 real vtot;
597 const t_oriresdata *od;
598 gmx_bool bTAV;
600 vtot = 0;
601 od = &(fcd->orires);
603 if (od->fc != 0)
605 bTAV = (od->edt != 0);
607 smooth_fc = od->fc;
608 if (bTAV)
610 /* Smoothly switch on the restraining when time averaging is used */
611 smooth_fc *= (1.0 - od->exp_min_t_tau);
614 d = 0;
615 for(fa=0; fa<nfa; fa+=3)
617 type = forceatoms[fa];
618 ai = forceatoms[fa+1];
619 aj = forceatoms[fa+2];
620 if (pbc)
622 ki = pbc_dx_aiuc(pbc,x[ai],x[aj],r);
624 else
626 rvec_sub(x[ai],x[aj],r);
628 r2 = norm2(r);
629 invr = gmx_invsqrt(r2);
630 invr2 = invr*invr;
631 ex = ip[type].orires.ex;
632 power = ip[type].orires.power;
633 fc = smooth_fc*ip[type].orires.kfac;
634 dev = od->otav[d] - ip[type].orires.obs;
636 /* NOTE:
637 * there is no real potential when time averaging is applied
639 vtot += 0.5*fc*sqr(dev);
641 if (bTAV)
643 /* Calculate the force as the sqrt of tav times instantaneous */
644 devins = od->oins[d] - ip[type].orires.obs;
645 if (dev*devins <= 0)
647 dev = 0;
649 else
651 dev = sqrt(dev*devins);
652 if (devins < 0)
654 dev = -dev;
659 pfac = fc*ip[type].orires.c*invr2;
660 for(i=0; i<power; i++)
662 pfac *= invr;
664 mvmul(od->S[ex],r,Sr);
665 for(i=0; i<DIM; i++)
667 fij[i] =
668 -pfac*dev*(4*Sr[i] - 2*(2+power)*invr2*iprod(Sr,r)*r[i]);
671 if (g)
673 ivec_sub(SHIFT_IVEC(g,ai),SHIFT_IVEC(g,aj),dt);
674 ki=IVEC2IS(dt);
677 for(i=0; i<DIM; i++)
679 f[ai][i] += fij[i];
680 f[aj][i] -= fij[i];
681 fshift[ki][i] += fij[i];
682 fshift[CENTRAL][i] -= fij[i];
684 d++;
688 return vtot;
690 /* Approx. 80*nfa/3 flops */
693 void update_orires_history(t_fcdata *fcd,history_t *hist)
695 t_oriresdata *od;
696 int pair,i;
698 od = &(fcd->orires);
699 if (od->edt != 0)
701 /* Copy the new time averages that have been calculated
702 * in calc_orires_dev.
704 hist->orire_initf = od->exp_min_t_tau;
705 for(pair=0; pair<od->nr; pair++)
707 for(i=0; i<5; i++)
709 hist->orire_Dtav[pair*5+i] = od->Dtav[pair][i];