added Verlet scheme and NxN non-bonded functionality
[gromacs.git] / src / mdlib / shellfc.c
blob86ee3c05aa58c64890beeda8b72fda638d5afe71
1 /*
2 *
3 * This source code is part of
4 *
5 * G R O M A C S
6 *
7 * GROningen MAchine for Chemical Simulations
8 *
9 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
10 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
11 * Copyright (c) 2001-2008, The GROMACS development team,
12 * check out http://www.gromacs.org for more information.
14 * This program is free software; you can redistribute it and/or
15 * modify it under the terms of the GNU General Public License
16 * as published by the Free Software Foundation; either version 2
17 * of the License, or (at your option) any later version.
19 * If you want to redistribute modifications, please consider that
20 * scientific software is very special. Version control is crucial -
21 * bugs must be traceable. We will be happy to consider code for
22 * inclusion in the official distribution, but derived work must not
23 * be called official GROMACS. Details are found in the README & COPYING
24 * files - if they are missing, get the official version at www.gromacs.org.
26 * To help us fund GROMACS development, we humbly ask that you cite
27 * the papers on the package - you can find them in the top README file.
29 * For more info, check our website at http://www.gromacs.org
31 * And Hey:
32 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
34 #ifdef HAVE_CONFIG_H
35 #include <config.h>
36 #endif
38 #include <string.h>
39 #include "typedefs.h"
40 #include "smalloc.h"
41 #include "gmx_fatal.h"
42 #include "vec.h"
43 #include "txtdump.h"
44 #include "mdrun.h"
45 #include "partdec.h"
46 #include "mdatoms.h"
47 #include "vsite.h"
48 #include "network.h"
49 #include "names.h"
50 #include "constr.h"
51 #include "domdec.h"
52 #include "partdec.h"
53 #include "physics.h"
54 #include "copyrite.h"
55 #include "shellfc.h"
56 #include "mtop_util.h"
57 #include "chargegroup.h"
60 typedef struct {
61 int nnucl;
62 atom_id shell; /* The shell id */
63 atom_id nucl1,nucl2,nucl3; /* The nuclei connected to the shell */
64 /* gmx_bool bInterCG; */ /* Coupled to nuclei outside cg? */
65 real k; /* force constant */
66 real k_1; /* 1 over force constant */
67 rvec xold;
68 rvec fold;
69 rvec step;
70 } t_shell;
72 typedef struct gmx_shellfc {
73 int nshell_gl; /* The number of shells in the system */
74 t_shell *shell_gl; /* All the shells (for DD only) */
75 int *shell_index_gl; /* Global shell index (for DD only) */
76 gmx_bool bInterCG; /* Are there inter charge-group shells? */
77 int nshell; /* The number of local shells */
78 t_shell *shell; /* The local shells */
79 int shell_nalloc; /* The allocation size of shell */
80 gmx_bool bPredict; /* Predict shell positions */
81 gmx_bool bRequireInit; /* Require initialization of shell positions */
82 int nflexcon; /* The number of flexible constraints */
83 rvec *x[2]; /* Array for iterative minimization */
84 rvec *f[2]; /* Array for iterative minimization */
85 int x_nalloc; /* The allocation size of x and f */
86 rvec *acc_dir; /* Acceleration direction for flexcon */
87 rvec *x_old; /* Old coordinates for flexcon */
88 int flex_nalloc; /* The allocation size of acc_dir and x_old */
89 rvec *adir_xnold; /* Work space for init_adir */
90 rvec *adir_xnew; /* Work space for init_adir */
91 int adir_nalloc; /* Work space for init_adir */
92 } t_gmx_shellfc;
95 static void pr_shell(FILE *fplog,int ns,t_shell s[])
97 int i;
99 fprintf(fplog,"SHELL DATA\n");
100 fprintf(fplog,"%5s %8s %5s %5s %5s\n",
101 "Shell","Force k","Nucl1","Nucl2","Nucl3");
102 for(i=0; (i<ns); i++) {
103 fprintf(fplog,"%5d %8.3f %5d",s[i].shell,1.0/s[i].k_1,s[i].nucl1);
104 if (s[i].nnucl == 2)
105 fprintf(fplog," %5d\n",s[i].nucl2);
106 else if (s[i].nnucl == 3)
107 fprintf(fplog," %5d %5d\n",s[i].nucl2,s[i].nucl3);
108 else
109 fprintf(fplog,"\n");
113 static void predict_shells(FILE *fplog,rvec x[],rvec v[],real dt,
114 int ns,t_shell s[],
115 real mass[],gmx_mtop_t *mtop,gmx_bool bInit)
117 int i,m,s1,n1,n2,n3;
118 real dt_1,dt_2,dt_3,fudge,tm,m1,m2,m3;
119 rvec *ptr;
120 gmx_mtop_atomlookup_t alook=NULL;
121 t_atom *atom;
123 if (mass == NULL) {
124 alook = gmx_mtop_atomlookup_init(mtop);
127 /* We introduce a fudge factor for performance reasons: with this choice
128 * the initial force on the shells is about a factor of two lower than
129 * without
131 fudge = 1.0;
133 if (bInit) {
134 if (fplog)
135 fprintf(fplog,"RELAX: Using prediction for initial shell placement\n");
136 ptr = x;
137 dt_1 = 1;
139 else {
140 ptr = v;
141 dt_1 = fudge*dt;
144 for(i=0; (i<ns); i++) {
145 s1 = s[i].shell;
146 if (bInit)
147 clear_rvec(x[s1]);
148 switch (s[i].nnucl) {
149 case 1:
150 n1 = s[i].nucl1;
151 for(m=0; (m<DIM); m++)
152 x[s1][m]+=ptr[n1][m]*dt_1;
153 break;
154 case 2:
155 n1 = s[i].nucl1;
156 n2 = s[i].nucl2;
157 if (mass) {
158 m1 = mass[n1];
159 m2 = mass[n2];
160 } else {
161 /* Not the correct masses with FE, but it is just a prediction... */
162 m1 = atom[n1].m;
163 m2 = atom[n2].m;
165 tm = dt_1/(m1+m2);
166 for(m=0; (m<DIM); m++)
167 x[s1][m]+=(m1*ptr[n1][m]+m2*ptr[n2][m])*tm;
168 break;
169 case 3:
170 n1 = s[i].nucl1;
171 n2 = s[i].nucl2;
172 n3 = s[i].nucl3;
173 if (mass) {
174 m1 = mass[n1];
175 m2 = mass[n2];
176 m3 = mass[n3];
177 } else {
178 /* Not the correct masses with FE, but it is just a prediction... */
179 gmx_mtop_atomnr_to_atom(alook,n1,&atom);
180 m1 = atom->m;
181 gmx_mtop_atomnr_to_atom(alook,n2,&atom);
182 m2 = atom->m;
183 gmx_mtop_atomnr_to_atom(alook,n3,&atom);
184 m3 = atom->m;
186 tm = dt_1/(m1+m2+m3);
187 for(m=0; (m<DIM); m++)
188 x[s1][m]+=(m1*ptr[n1][m]+m2*ptr[n2][m]+m3*ptr[n3][m])*tm;
189 break;
190 default:
191 gmx_fatal(FARGS,"Shell %d has %d nuclei!",i,s[i].nnucl);
195 if (mass == NULL) {
196 gmx_mtop_atomlookup_destroy(alook);
200 gmx_shellfc_t init_shell_flexcon(FILE *fplog,
201 gmx_mtop_t *mtop,int nflexcon,
202 rvec *x)
204 struct gmx_shellfc *shfc;
205 t_shell *shell;
206 int *shell_index=NULL,*at2cg;
207 t_atom *atom;
208 int n[eptNR],ns,nshell,nsi;
209 int i,j,nmol,type,mb,mt,a_offset,cg,mol,ftype,nra;
210 real qS,alpha;
211 int aS,aN=0; /* Shell and nucleus */
212 int bondtypes[] = { F_BONDS, F_HARMONIC, F_CUBICBONDS, F_POLARIZATION, F_ANHARM_POL, F_WATER_POL };
213 #define NBT asize(bondtypes)
214 t_iatom *ia;
215 gmx_mtop_atomloop_block_t aloopb;
216 gmx_mtop_atomloop_all_t aloop;
217 gmx_ffparams_t *ffparams;
218 gmx_molblock_t *molb;
219 gmx_moltype_t *molt;
220 t_block *cgs;
222 /* Count number of shells, and find their indices */
223 for(i=0; (i<eptNR); i++) {
224 n[i] = 0;
227 aloopb = gmx_mtop_atomloop_block_init(mtop);
228 while (gmx_mtop_atomloop_block_next(aloopb,&atom,&nmol)) {
229 n[atom->ptype] += nmol;
232 if (fplog) {
233 /* Print the number of each particle type */
234 for(i=0; (i<eptNR); i++) {
235 if (n[i] != 0) {
236 fprintf(fplog,"There are: %d %ss\n",n[i],ptype_str[i]);
241 nshell = n[eptShell];
243 if (nshell == 0 && nflexcon == 0) {
244 return NULL;
247 snew(shfc,1);
248 shfc->nflexcon = nflexcon;
250 if (nshell == 0) {
251 return shfc;
254 /* We have shells: fill the shell data structure */
256 /* Global system sized array, this should be avoided */
257 snew(shell_index,mtop->natoms);
259 aloop = gmx_mtop_atomloop_all_init(mtop);
260 nshell = 0;
261 while (gmx_mtop_atomloop_all_next(aloop,&i,&atom)) {
262 if (atom->ptype == eptShell) {
263 shell_index[i] = nshell++;
267 snew(shell,nshell);
269 /* Initiate the shell structures */
270 for(i=0; (i<nshell); i++) {
271 shell[i].shell = NO_ATID;
272 shell[i].nnucl = 0;
273 shell[i].nucl1 = NO_ATID;
274 shell[i].nucl2 = NO_ATID;
275 shell[i].nucl3 = NO_ATID;
276 /* shell[i].bInterCG=FALSE; */
277 shell[i].k_1 = 0;
278 shell[i].k = 0;
281 ffparams = &mtop->ffparams;
283 /* Now fill the structures */
284 shfc->bInterCG = FALSE;
285 ns = 0;
286 a_offset = 0;
287 for(mb=0; mb<mtop->nmolblock; mb++) {
288 molb = &mtop->molblock[mb];
289 molt = &mtop->moltype[molb->type];
291 cgs = &molt->cgs;
292 snew(at2cg,molt->atoms.nr);
293 for(cg=0; cg<cgs->nr; cg++) {
294 for(i=cgs->index[cg]; i<cgs->index[cg+1]; i++) {
295 at2cg[i] = cg;
299 atom = molt->atoms.atom;
300 for(mol=0; mol<molb->nmol; mol++) {
301 for(j=0; (j<NBT); j++) {
302 ia = molt->ilist[bondtypes[j]].iatoms;
303 for(i=0; (i<molt->ilist[bondtypes[j]].nr); ) {
304 type = ia[0];
305 ftype = ffparams->functype[type];
306 nra = interaction_function[ftype].nratoms;
308 /* Check whether we have a bond with a shell */
309 aS = NO_ATID;
311 switch (bondtypes[j]) {
312 case F_BONDS:
313 case F_HARMONIC:
314 case F_CUBICBONDS:
315 case F_POLARIZATION:
316 case F_ANHARM_POL:
317 if (atom[ia[1]].ptype == eptShell) {
318 aS = ia[1];
319 aN = ia[2];
321 else if (atom[ia[2]].ptype == eptShell) {
322 aS = ia[2];
323 aN = ia[1];
325 break;
326 case F_WATER_POL:
327 aN = ia[4]; /* Dummy */
328 aS = ia[5]; /* Shell */
329 break;
330 default:
331 gmx_fatal(FARGS,"Death Horror: %s, %d",__FILE__,__LINE__);
334 if (aS != NO_ATID) {
335 qS = atom[aS].q;
337 /* Check whether one of the particles is a shell... */
338 nsi = shell_index[a_offset+aS];
339 if ((nsi < 0) || (nsi >= nshell))
340 gmx_fatal(FARGS,"nsi is %d should be within 0 - %d. aS = %d",
341 nsi,nshell,aS);
342 if (shell[nsi].shell == NO_ATID) {
343 shell[nsi].shell = a_offset + aS;
344 ns ++;
346 else if (shell[nsi].shell != a_offset+aS)
347 gmx_fatal(FARGS,"Weird stuff in %s, %d",__FILE__,__LINE__);
349 if (shell[nsi].nucl1 == NO_ATID) {
350 shell[nsi].nucl1 = a_offset + aN;
351 } else if (shell[nsi].nucl2 == NO_ATID) {
352 shell[nsi].nucl2 = a_offset + aN;
353 } else if (shell[nsi].nucl3 == NO_ATID) {
354 shell[nsi].nucl3 = a_offset + aN;
355 } else {
356 if (fplog)
357 pr_shell(fplog,ns,shell);
358 gmx_fatal(FARGS,"Can not handle more than three bonds per shell\n");
360 if (at2cg[aS] != at2cg[aN]) {
361 /* shell[nsi].bInterCG = TRUE; */
362 shfc->bInterCG = TRUE;
365 switch (bondtypes[j]) {
366 case F_BONDS:
367 case F_HARMONIC:
368 shell[nsi].k += ffparams->iparams[type].harmonic.krA;
369 break;
370 case F_CUBICBONDS:
371 shell[nsi].k += ffparams->iparams[type].cubic.kb;
372 break;
373 case F_POLARIZATION:
374 case F_ANHARM_POL:
375 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS*10))
376 gmx_fatal(FARGS,"polarize can not be used with qA(%e) != qB(%e) for atom %d of molecule block %d", qS, atom[aS].qB, aS+1, mb+1);
377 shell[nsi].k += sqr(qS)*ONE_4PI_EPS0/
378 ffparams->iparams[type].polarize.alpha;
379 break;
380 case F_WATER_POL:
381 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS*10))
382 gmx_fatal(FARGS,"water_pol can not be used with qA(%e) != qB(%e) for atom %d of molecule block %d", qS, atom[aS].qB, aS+1, mb+1);
383 alpha = (ffparams->iparams[type].wpol.al_x+
384 ffparams->iparams[type].wpol.al_y+
385 ffparams->iparams[type].wpol.al_z)/3.0;
386 shell[nsi].k += sqr(qS)*ONE_4PI_EPS0/alpha;
387 break;
388 default:
389 gmx_fatal(FARGS,"Death Horror: %s, %d",__FILE__,__LINE__);
391 shell[nsi].nnucl++;
393 ia += nra+1;
394 i += nra+1;
397 a_offset += molt->atoms.nr;
399 /* Done with this molecule type */
400 sfree(at2cg);
403 /* Verify whether it's all correct */
404 if (ns != nshell)
405 gmx_fatal(FARGS,"Something weird with shells. They may not be bonded to something");
407 for(i=0; (i<ns); i++)
408 shell[i].k_1 = 1.0/shell[i].k;
410 if (debug)
411 pr_shell(debug,ns,shell);
414 shfc->nshell_gl = ns;
415 shfc->shell_gl = shell;
416 shfc->shell_index_gl = shell_index;
418 shfc->bPredict = (getenv("GMX_NOPREDICT") == NULL);
419 shfc->bRequireInit = FALSE;
420 if (!shfc->bPredict) {
421 if (fplog)
422 fprintf(fplog,"\nWill never predict shell positions\n");
423 } else {
424 shfc->bRequireInit = (getenv("GMX_REQUIRE_SHELL_INIT") != NULL);
425 if (shfc->bRequireInit && fplog)
426 fprintf(fplog,"\nWill always initiate shell positions\n");
429 if (shfc->bPredict) {
430 if (x) {
431 predict_shells(fplog,x,NULL,0,shfc->nshell_gl,shfc->shell_gl,
432 NULL,mtop,TRUE);
435 if (shfc->bInterCG) {
436 if (fplog)
437 fprintf(fplog,"\nNOTE: there all shells that are connected to particles outside thier own charge group, will not predict shells positions during the run\n\n");
438 shfc->bPredict = FALSE;
442 return shfc;
445 void make_local_shells(t_commrec *cr,t_mdatoms *md,
446 struct gmx_shellfc *shfc)
448 t_shell *shell;
449 int a0,a1,*ind,nshell,i;
450 gmx_domdec_t *dd=NULL;
452 if (PAR(cr)) {
453 if (DOMAINDECOMP(cr)) {
454 dd = cr->dd;
455 a0 = 0;
456 a1 = dd->nat_home;
457 } else {
458 pd_at_range(cr,&a0,&a1);
460 } else {
461 /* Single node: we need all shells, just copy the pointer */
462 shfc->nshell = shfc->nshell_gl;
463 shfc->shell = shfc->shell_gl;
465 return;
468 ind = shfc->shell_index_gl;
470 nshell = 0;
471 shell = shfc->shell;
472 for(i=a0; i<a1; i++) {
473 if (md->ptype[i] == eptShell) {
474 if (nshell+1 > shfc->shell_nalloc) {
475 shfc->shell_nalloc = over_alloc_dd(nshell+1);
476 srenew(shell,shfc->shell_nalloc);
478 if (dd) {
479 shell[nshell] = shfc->shell_gl[ind[dd->gatindex[i]]];
480 } else {
481 shell[nshell] = shfc->shell_gl[ind[i]];
483 /* With inter-cg shells we can no do shell prediction,
484 * so we do not need the nuclei numbers.
486 if (!shfc->bInterCG) {
487 shell[nshell].nucl1 = i + shell[nshell].nucl1 - shell[nshell].shell;
488 if (shell[nshell].nnucl > 1)
489 shell[nshell].nucl2 = i + shell[nshell].nucl2 - shell[nshell].shell;
490 if (shell[nshell].nnucl > 2)
491 shell[nshell].nucl3 = i + shell[nshell].nucl3 - shell[nshell].shell;
493 shell[nshell].shell = i;
494 nshell++;
498 shfc->nshell = nshell;
499 shfc->shell = shell;
502 static void do_1pos(rvec xnew,rvec xold,rvec f,real step)
504 real xo,yo,zo;
505 real dx,dy,dz;
507 xo=xold[XX];
508 yo=xold[YY];
509 zo=xold[ZZ];
511 dx=f[XX]*step;
512 dy=f[YY]*step;
513 dz=f[ZZ]*step;
515 xnew[XX]=xo+dx;
516 xnew[YY]=yo+dy;
517 xnew[ZZ]=zo+dz;
520 static void do_1pos3(rvec xnew,rvec xold,rvec f,rvec step)
522 real xo,yo,zo;
523 real dx,dy,dz;
525 xo=xold[XX];
526 yo=xold[YY];
527 zo=xold[ZZ];
529 dx=f[XX]*step[XX];
530 dy=f[YY]*step[YY];
531 dz=f[ZZ]*step[ZZ];
533 xnew[XX]=xo+dx;
534 xnew[YY]=yo+dy;
535 xnew[ZZ]=zo+dz;
538 static void directional_sd(FILE *log,rvec xold[],rvec xnew[],rvec acc_dir[],
539 int start,int homenr,real step)
541 int i;
543 for(i=start; i<homenr; i++)
544 do_1pos(xnew[i],xold[i],acc_dir[i],step);
547 static void shell_pos_sd(FILE *log,rvec xcur[],rvec xnew[],rvec f[],
548 int ns,t_shell s[],int count)
550 const real step_scale_min = 0.8,
551 step_scale_increment = 0.2,
552 step_scale_max = 1.2,
553 step_scale_multiple = (step_scale_max - step_scale_min) / step_scale_increment;
554 int i,shell,d;
555 real dx,df,k_est;
556 #ifdef PRINT_STEP
557 real step_min,step_max;
559 step_min = 1e30;
560 step_max = 0;
561 #endif
562 for(i=0; (i<ns); i++) {
563 shell = s[i].shell;
564 if (count == 1) {
565 for(d=0; d<DIM; d++) {
566 s[i].step[d] = s[i].k_1;
567 #ifdef PRINT_STEP
568 step_min = min(step_min,s[i].step[d]);
569 step_max = max(step_max,s[i].step[d]);
570 #endif
572 } else {
573 for(d=0; d<DIM; d++) {
574 dx = xcur[shell][d] - s[i].xold[d];
575 df = f[shell][d] - s[i].fold[d];
576 /* -dx/df gets used to generate an interpolated value, but would
577 * cause a NaN if df were binary-equal to zero. Values close to
578 * zero won't cause problems (because of the min() and max()), so
579 * just testing for binary inequality is OK. */
580 if (0.0 != df)
582 k_est = -dx/df;
583 /* Scale the step size by a factor interpolated from
584 * step_scale_min to step_scale_max, as k_est goes from 0 to
585 * step_scale_multiple * s[i].step[d] */
586 s[i].step[d] =
587 step_scale_min * s[i].step[d] +
588 step_scale_increment * min(step_scale_multiple * s[i].step[d], max(k_est, 0));
590 else
592 /* Here 0 == df */
593 if (gmx_numzero(dx)) /* 0 == dx */
595 /* Likely this will never happen, but if it does just
596 * don't scale the step. */
598 else /* 0 != dx */
600 s[i].step[d] *= step_scale_max;
603 #ifdef PRINT_STEP
604 step_min = min(step_min,s[i].step[d]);
605 step_max = max(step_max,s[i].step[d]);
606 #endif
609 copy_rvec(xcur[shell],s[i].xold);
610 copy_rvec(f[shell], s[i].fold);
612 do_1pos3(xnew[shell],xcur[shell],f[shell],s[i].step);
614 if (gmx_debug_at) {
615 fprintf(debug,"shell[%d] = %d\n",i,shell);
616 pr_rvec(debug,0,"fshell",f[shell],DIM,TRUE);
617 pr_rvec(debug,0,"xold",xcur[shell],DIM,TRUE);
618 pr_rvec(debug,0,"step",s[i].step,DIM,TRUE);
619 pr_rvec(debug,0,"xnew",xnew[shell],DIM,TRUE);
622 #ifdef PRINT_STEP
623 printf("step %.3e %.3e\n",step_min,step_max);
624 #endif
627 static void decrease_step_size(int nshell,t_shell s[])
629 int i;
631 for(i=0; i<nshell; i++)
632 svmul(0.8,s[i].step,s[i].step);
635 static void print_epot(FILE *fp,gmx_large_int_t mdstep,int count,real epot,real df,
636 int ndir,real sf_dir)
638 char buf[22];
640 fprintf(fp,"MDStep=%5s/%2d EPot: %12.8e, rmsF: %6.2e",
641 gmx_step_str(mdstep,buf),count,epot,df);
642 if (ndir)
643 fprintf(fp,", dir. rmsF: %6.2e\n",sqrt(sf_dir/ndir));
644 else
645 fprintf(fp,"\n");
649 static real rms_force(t_commrec *cr,rvec f[],int ns,t_shell s[],
650 int ndir,real *sf_dir,real *Epot)
652 int i,shell,ntot;
653 double buf[4];
655 buf[0] = *sf_dir;
656 for(i=0; i<ns; i++) {
657 shell = s[i].shell;
658 buf[0] += norm2(f[shell]);
660 ntot = ns;
662 if (PAR(cr)) {
663 buf[1] = ntot;
664 buf[2] = *sf_dir;
665 buf[3] = *Epot;
666 gmx_sumd(4,buf,cr);
667 ntot = (int)(buf[1] + 0.5);
668 *sf_dir = buf[2];
669 *Epot = buf[3];
671 ntot += ndir;
673 return (ntot ? sqrt(buf[0]/ntot) : 0);
676 static void check_pbc(FILE *fp,rvec x[],int shell)
678 int m,now;
680 now = shell-4;
681 for(m=0; (m<DIM); m++)
682 if (fabs(x[shell][m]-x[now][m]) > 0.3) {
683 pr_rvecs(fp,0,"SHELL-X",x+now,5);
684 break;
688 static void dump_shells(FILE *fp,rvec x[],rvec f[],real ftol,int ns,t_shell s[])
690 int i,shell;
691 real ft2,ff2;
693 ft2 = sqr(ftol);
695 for(i=0; (i<ns); i++) {
696 shell = s[i].shell;
697 ff2 = iprod(f[shell],f[shell]);
698 if (ff2 > ft2)
699 fprintf(fp,"SHELL %5d, force %10.5f %10.5f %10.5f, |f| %10.5f\n",
700 shell,f[shell][XX],f[shell][YY],f[shell][ZZ],sqrt(ff2));
701 check_pbc(fp,x,shell);
705 static void init_adir(FILE *log,gmx_shellfc_t shfc,
706 gmx_constr_t constr,t_idef *idef,t_inputrec *ir,
707 t_commrec *cr,int dd_ac1,
708 gmx_large_int_t step,t_mdatoms *md,int start,int end,
709 rvec *x_old,rvec *x_init,rvec *x,
710 rvec *f,rvec *acc_dir,
711 gmx_bool bMolPBC,matrix box,
712 real *lambda,real *dvdlambda,t_nrnb *nrnb)
714 rvec *xnold,*xnew;
715 double w_dt;
716 int gf,ga,gt;
717 real dt,scale;
718 int n,d;
719 unsigned short *ptype;
720 rvec p,dx;
722 if (DOMAINDECOMP(cr))
723 n = dd_ac1;
724 else
725 n = end - start;
726 if (n > shfc->adir_nalloc) {
727 shfc->adir_nalloc = over_alloc_dd(n);
728 srenew(shfc->adir_xnold,shfc->adir_nalloc);
729 srenew(shfc->adir_xnew ,shfc->adir_nalloc);
731 xnold = shfc->adir_xnold;
732 xnew = shfc->adir_xnew;
734 ptype = md->ptype;
736 dt = ir->delta_t;
738 /* Does NOT work with freeze or acceleration groups (yet) */
739 for (n=start; n<end; n++) {
740 w_dt = md->invmass[n]*dt;
742 for (d=0; d<DIM; d++) {
743 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell)) {
744 xnold[n-start][d] = x[n][d] - (x_init[n][d] - x_old[n][d]);
745 xnew[n-start][d] = 2*x[n][d] - x_old[n][d] + f[n][d]*w_dt*dt;
746 } else {
747 xnold[n-start][d] = x[n][d];
748 xnew[n-start][d] = x[n][d];
752 constrain(log,FALSE,FALSE,constr,idef,ir,NULL,cr,step,0,md,
753 x,xnold-start,NULL,bMolPBC,box,
754 lambda[efptBONDED],&(dvdlambda[efptBONDED]),
755 NULL,NULL,nrnb,econqCoord,FALSE,0,0);
756 constrain(log,FALSE,FALSE,constr,idef,ir,NULL,cr,step,0,md,
757 x,xnew-start,NULL,bMolPBC,box,
758 lambda[efptBONDED],&(dvdlambda[efptBONDED]),
759 NULL,NULL,nrnb,econqCoord,FALSE,0,0);
761 for (n=start; n<end; n++) {
762 for(d=0; d<DIM; d++)
763 xnew[n-start][d] =
764 -(2*x[n][d]-xnold[n-start][d]-xnew[n-start][d])/sqr(dt)
765 - f[n][d]*md->invmass[n];
766 clear_rvec(acc_dir[n]);
769 /* Project the acceleration on the old bond directions */
770 constrain(log,FALSE,FALSE,constr,idef,ir,NULL,cr,step,0,md,
771 x_old,xnew-start,acc_dir,bMolPBC,box,
772 lambda[efptBONDED],&(dvdlambda[efptBONDED]),
773 NULL,NULL,nrnb,econqDeriv_FlexCon,FALSE,0,0);
776 int relax_shell_flexcon(FILE *fplog,t_commrec *cr,gmx_bool bVerbose,
777 gmx_large_int_t mdstep,t_inputrec *inputrec,
778 gmx_bool bDoNS,int force_flags,
779 gmx_bool bStopCM,
780 gmx_localtop_t *top,
781 gmx_mtop_t* mtop,
782 gmx_constr_t constr,
783 gmx_enerdata_t *enerd,t_fcdata *fcd,
784 t_state *state,rvec f[],
785 tensor force_vir,
786 t_mdatoms *md,
787 t_nrnb *nrnb,gmx_wallcycle_t wcycle,
788 t_graph *graph,
789 gmx_groups_t *groups,
790 struct gmx_shellfc *shfc,
791 t_forcerec *fr,
792 gmx_bool bBornRadii,
793 double t,rvec mu_tot,
794 int natoms,gmx_bool *bConverged,
795 gmx_vsite_t *vsite,
796 FILE *fp_field)
798 int nshell;
799 t_shell *shell;
800 t_idef *idef;
801 rvec *pos[2],*force[2],*acc_dir=NULL,*x_old=NULL;
802 real Epot[2],df[2];
803 rvec dx;
804 real sf_dir,invdt;
805 real ftol,xiH,xiS,dum=0;
806 char sbuf[22];
807 gmx_bool bCont,bInit;
808 int nat,dd_ac0,dd_ac1=0,i;
809 int start=md->start,homenr=md->homenr,end=start+homenr,cg0,cg1;
810 int nflexcon,g,number_steps,d,Min=0,count=0;
811 #define Try (1-Min) /* At start Try = 1 */
813 bCont = (mdstep == inputrec->init_step) && inputrec->bContinuation;
814 bInit = (mdstep == inputrec->init_step) || shfc->bRequireInit;
815 ftol = inputrec->em_tol;
816 number_steps = inputrec->niter;
817 nshell = shfc->nshell;
818 shell = shfc->shell;
819 nflexcon = shfc->nflexcon;
821 idef = &top->idef;
823 if (DOMAINDECOMP(cr)) {
824 nat = dd_natoms_vsite(cr->dd);
825 if (nflexcon > 0) {
826 dd_get_constraint_range(cr->dd,&dd_ac0,&dd_ac1);
827 nat = max(nat,dd_ac1);
829 } else {
830 nat = state->natoms;
833 if (nat > shfc->x_nalloc) {
834 /* Allocate local arrays */
835 shfc->x_nalloc = over_alloc_dd(nat);
836 for(i=0; (i<2); i++) {
837 srenew(shfc->x[i],shfc->x_nalloc);
838 srenew(shfc->f[i],shfc->x_nalloc);
841 for(i=0; (i<2); i++) {
842 pos[i] = shfc->x[i];
843 force[i] = shfc->f[i];
846 /* With particle decomposition this code only works
847 * when all particles involved with each shell are in the same cg.
850 if (bDoNS && inputrec->ePBC != epbcNONE && !DOMAINDECOMP(cr)) {
851 /* This is the only time where the coordinates are used
852 * before do_force is called, which normally puts all
853 * charge groups in the box.
855 if (PARTDECOMP(cr)) {
856 pd_cg_range(cr,&cg0,&cg1);
857 } else {
858 cg0 = 0;
859 cg1 = top->cgs.nr;
861 put_charge_groups_in_box(fplog,cg0,cg1,fr->ePBC,state->box,
862 &(top->cgs),state->x,fr->cg_cm);
863 if (graph)
864 mk_mshift(fplog,graph,fr->ePBC,state->box,state->x);
867 /* After this all coordinate arrays will contain whole molecules */
868 if (graph)
869 shift_self(graph,state->box,state->x);
871 if (nflexcon) {
872 if (nat > shfc->flex_nalloc) {
873 shfc->flex_nalloc = over_alloc_dd(nat);
874 srenew(shfc->acc_dir,shfc->flex_nalloc);
875 srenew(shfc->x_old,shfc->flex_nalloc);
877 acc_dir = shfc->acc_dir;
878 x_old = shfc->x_old;
879 for(i=0; i<homenr; i++) {
880 for(d=0; d<DIM; d++)
881 shfc->x_old[i][d] =
882 state->x[start+i][d] - state->v[start+i][d]*inputrec->delta_t;
886 /* Do a prediction of the shell positions */
887 if (shfc->bPredict && !bCont) {
888 predict_shells(fplog,state->x,state->v,inputrec->delta_t,nshell,shell,
889 md->massT,NULL,bInit);
892 /* do_force expected the charge groups to be in the box */
893 if (graph)
894 unshift_self(graph,state->box,state->x);
896 /* Calculate the forces first time around */
897 if (gmx_debug_at) {
898 pr_rvecs(debug,0,"x b4 do_force",state->x + start,homenr);
900 do_force(fplog,cr,inputrec,mdstep,nrnb,wcycle,top,mtop,groups,
901 state->box,state->x,&state->hist,
902 force[Min],force_vir,md,enerd,fcd,
903 state->lambda,graph,
904 fr,vsite,mu_tot,t,fp_field,NULL,bBornRadii,
905 (bDoNS ? GMX_FORCE_NS : 0) | force_flags);
907 sf_dir = 0;
908 if (nflexcon) {
909 init_adir(fplog,shfc,
910 constr,idef,inputrec,cr,dd_ac1,mdstep,md,start,end,
911 shfc->x_old-start,state->x,state->x,force[Min],
912 shfc->acc_dir-start,
913 fr->bMolPBC,state->box,state->lambda,&dum,nrnb);
915 for(i=start; i<end; i++)
916 sf_dir += md->massT[i]*norm2(shfc->acc_dir[i-start]);
919 Epot[Min] = enerd->term[F_EPOT];
921 df[Min]=rms_force(cr,shfc->f[Min],nshell,shell,nflexcon,&sf_dir,&Epot[Min]);
922 df[Try]=0;
923 if (debug) {
924 fprintf(debug,"df = %g %g\n",df[Min],df[Try]);
927 if (gmx_debug_at) {
928 pr_rvecs(debug,0,"force0",force[Min],md->nr);
931 if (nshell+nflexcon > 0) {
932 /* Copy x to pos[Min] & pos[Try]: during minimization only the
933 * shell positions are updated, therefore the other particles must
934 * be set here.
936 memcpy(pos[Min],state->x,nat*sizeof(state->x[0]));
937 memcpy(pos[Try],state->x,nat*sizeof(state->x[0]));
940 if (bVerbose && MASTER(cr))
941 print_epot(stdout,mdstep,0,Epot[Min],df[Min],nflexcon,sf_dir);
943 if (debug) {
944 fprintf(debug,"%17s: %14.10e\n",
945 interaction_function[F_EKIN].longname,enerd->term[F_EKIN]);
946 fprintf(debug,"%17s: %14.10e\n",
947 interaction_function[F_EPOT].longname,enerd->term[F_EPOT]);
948 fprintf(debug,"%17s: %14.10e\n",
949 interaction_function[F_ETOT].longname,enerd->term[F_ETOT]);
950 fprintf(debug,"SHELLSTEP %s\n",gmx_step_str(mdstep,sbuf));
953 /* First check whether we should do shells, or whether the force is
954 * low enough even without minimization.
956 *bConverged = (df[Min] < ftol);
958 for(count=1; (!(*bConverged) && (count < number_steps)); count++) {
959 if (vsite)
960 construct_vsites(fplog,vsite,pos[Min],nrnb,inputrec->delta_t,state->v,
961 idef->iparams,idef->il,
962 fr->ePBC,fr->bMolPBC,graph,cr,state->box);
964 if (nflexcon) {
965 init_adir(fplog,shfc,
966 constr,idef,inputrec,cr,dd_ac1,mdstep,md,start,end,
967 x_old-start,state->x,pos[Min],force[Min],acc_dir-start,
968 fr->bMolPBC,state->box,state->lambda,&dum,nrnb);
970 directional_sd(fplog,pos[Min],pos[Try],acc_dir-start,start,end,
971 fr->fc_stepsize);
974 /* New positions, Steepest descent */
975 shell_pos_sd(fplog,pos[Min],pos[Try],force[Min],nshell,shell,count);
977 /* do_force expected the charge groups to be in the box */
978 if (graph)
979 unshift_self(graph,state->box,pos[Try]);
981 if (gmx_debug_at) {
982 pr_rvecs(debug,0,"RELAX: pos[Min] ",pos[Min] + start,homenr);
983 pr_rvecs(debug,0,"RELAX: pos[Try] ",pos[Try] + start,homenr);
985 /* Try the new positions */
986 do_force(fplog,cr,inputrec,1,nrnb,wcycle,
987 top,mtop,groups,state->box,pos[Try],&state->hist,
988 force[Try],force_vir,
989 md,enerd,fcd,state->lambda,graph,
990 fr,vsite,mu_tot,t,fp_field,NULL,bBornRadii,
991 force_flags);
993 if (gmx_debug_at) {
994 pr_rvecs(debug,0,"RELAX: force[Min]",force[Min] + start,homenr);
995 pr_rvecs(debug,0,"RELAX: force[Try]",force[Try] + start,homenr);
997 sf_dir = 0;
998 if (nflexcon) {
999 init_adir(fplog,shfc,
1000 constr,idef,inputrec,cr,dd_ac1,mdstep,md,start,end,
1001 x_old-start,state->x,pos[Try],force[Try],acc_dir-start,
1002 fr->bMolPBC,state->box,state->lambda,&dum,nrnb);
1004 for(i=start; i<end; i++)
1005 sf_dir += md->massT[i]*norm2(acc_dir[i-start]);
1008 Epot[Try] = enerd->term[F_EPOT];
1010 df[Try]=rms_force(cr,force[Try],nshell,shell,nflexcon,&sf_dir,&Epot[Try]);
1012 if (debug)
1013 fprintf(debug,"df = %g %g\n",df[Min],df[Try]);
1015 if (debug) {
1016 if (gmx_debug_at)
1017 pr_rvecs(debug,0,"F na do_force",force[Try] + start,homenr);
1018 if (gmx_debug_at) {
1019 fprintf(debug,"SHELL ITER %d\n",count);
1020 dump_shells(debug,pos[Try],force[Try],ftol,nshell,shell);
1024 if (bVerbose && MASTER(cr))
1025 print_epot(stdout,mdstep,count,Epot[Try],df[Try],nflexcon,sf_dir);
1027 *bConverged = (df[Try] < ftol);
1029 if ((df[Try] < df[Min])) {
1030 if (debug)
1031 fprintf(debug,"Swapping Min and Try\n");
1032 if (nflexcon) {
1033 /* Correct the velocities for the flexible constraints */
1034 invdt = 1/inputrec->delta_t;
1035 for(i=start; i<end; i++) {
1036 for(d=0; d<DIM; d++)
1037 state->v[i][d] += (pos[Try][i][d] - pos[Min][i][d])*invdt;
1040 Min = Try;
1041 } else {
1042 decrease_step_size(nshell,shell);
1045 if (MASTER(cr) && !(*bConverged)) {
1046 /* Note that the energies and virial are incorrect when not converged */
1047 if (fplog)
1048 fprintf(fplog,
1049 "step %s: EM did not converge in %d iterations, RMS force %.3f\n",
1050 gmx_step_str(mdstep,sbuf),number_steps,df[Min]);
1051 fprintf(stderr,
1052 "step %s: EM did not converge in %d iterations, RMS force %.3f\n",
1053 gmx_step_str(mdstep,sbuf),number_steps,df[Min]);
1056 /* Copy back the coordinates and the forces */
1057 memcpy(state->x,pos[Min],nat*sizeof(state->x[0]));
1058 memcpy(f,force[Min],nat*sizeof(f[0]));
1060 return count;