static void write_constr_pdb(const char *fn, const char *title, gmx_mtop_t *mtop, int start, int homenr, t_commrec *cr, rvec x[], matrix box) { char fname[STRLEN], format[STRLEN]; FILE *out; int dd_ac0 = 0, dd_ac1 = 0, i, ii, resnr; gmx_domdec_t *dd; char *anm, *resnm; dd = NULL; if (DOMAINDECOMP(cr)) { dd = cr->dd; dd_get_constraint_range(dd, &dd_ac0, &dd_ac1); start = 0; homenr = dd_ac1; } if (PAR(cr)) { sprintf(fname, "%s_n%d.pdb", fn, cr->sim_nodeid); } else { sprintf(fname, "%s.pdb", fn); } sprintf(format, "%s\n", get_pdbformat()); out = gmx_fio_fopen(fname, "w"); fprintf(out, "TITLE %s\n", title); gmx_write_pdb_box(out, -1, box); for (i = start; i < start+homenr; i++) { if (dd != NULL) { if (i >= dd->nat_home && i < dd_ac0) { continue; } ii = dd->gatindex[i]; } else { ii = i; } gmx_mtop_atominfo_global(mtop, ii, &anm, &resnr, &resnm); fprintf(out, format, "ATOM", (ii+1)%100000, anm, resnm, ' ', resnr%10000, ' ', 10*x[i][XX], 10*x[i][YY], 10*x[i][ZZ]); } fprintf(out, "TER\n"); gmx_fio_fclose(out); }
static void dump_conf(gmx_domdec_t *dd,struct gmx_lincsdata *li, t_blocka *at2con, char *name,gmx_bool bAll,rvec *x,matrix box) { char str[STRLEN]; FILE *fp; int ac0,ac1,i; dd_get_constraint_range(dd,&ac0,&ac1); sprintf(str,"%s_%d_%d_%d.pdb",name,dd->ci[XX],dd->ci[YY],dd->ci[ZZ]); fp = gmx_fio_fopen(str,"w"); fprintf(fp,"CRYST1%9.3f%9.3f%9.3f%7.2f%7.2f%7.2f P 1 1\n", 10*norm(box[XX]),10*norm(box[YY]),10*norm(box[ZZ]), 90.0,90.0,90.0); for(i=0; i<ac1; i++) { if (i < dd->nat_home || (bAll && i >= ac0 && i < ac1)) { fprintf(fp,"%-6s%5u %-4.4s%3.3s %c%4d %8.3f%8.3f%8.3f%6.2f%6.2f\n", "ATOM",ddglatnr(dd,i),"C","ALA",' ',i+1, 10*x[i][XX],10*x[i][YY],10*x[i][ZZ], 1.0,i<dd->nat_tot ? 0.0 : 1.0); } } if (bAll) { for(i=0; i<li->nc; i++) { fprintf(fp,"CONECT%5d%5d\n", ddglatnr(dd,li->bla[2*i]), ddglatnr(dd,li->bla[2*i+1])); } } gmx_fio_fclose(fp); }
int relax_shell_flexcon(FILE *fplog, t_commrec *cr, gmx_bool bVerbose, gmx_int64_t mdstep, t_inputrec *inputrec, gmx_bool bDoNS, int force_flags, gmx_localtop_t *top, gmx_constr_t constr, gmx_enerdata_t *enerd, t_fcdata *fcd, t_state *state, rvec f[], tensor force_vir, t_mdatoms *md, t_nrnb *nrnb, gmx_wallcycle_t wcycle, t_graph *graph, gmx_groups_t *groups, struct gmx_shellfc *shfc, t_forcerec *fr, gmx_bool bBornRadii, double t, rvec mu_tot, gmx_bool *bConverged, gmx_vsite_t *vsite, FILE *fp_field) { int nshell; t_shell *shell; t_idef *idef; rvec *pos[2], *force[2], *acc_dir = NULL, *x_old = NULL; real Epot[2], df[2]; rvec dx; real sf_dir, invdt; real ftol, xiH, xiS, dum = 0; char sbuf[22]; gmx_bool bCont, bInit; int nat, dd_ac0, dd_ac1 = 0, i; int start = 0, homenr = md->homenr, end = start+homenr, cg0, cg1; int nflexcon, g, number_steps, d, Min = 0, count = 0; #define Try (1-Min) /* At start Try = 1 */ bCont = (mdstep == inputrec->init_step) && inputrec->bContinuation; bInit = (mdstep == inputrec->init_step) || shfc->bRequireInit; ftol = inputrec->em_tol; number_steps = inputrec->niter; nshell = shfc->nshell; shell = shfc->shell; nflexcon = shfc->nflexcon; idef = &top->idef; if (DOMAINDECOMP(cr)) { nat = dd_natoms_vsite(cr->dd); if (nflexcon > 0) { dd_get_constraint_range(cr->dd, &dd_ac0, &dd_ac1); nat = max(nat, dd_ac1); } } else { nat = state->natoms; } if (nat > shfc->x_nalloc) { /* Allocate local arrays */ shfc->x_nalloc = over_alloc_dd(nat); for (i = 0; (i < 2); i++) { srenew(shfc->x[i], shfc->x_nalloc); srenew(shfc->f[i], shfc->x_nalloc); } } for (i = 0; (i < 2); i++) { pos[i] = shfc->x[i]; force[i] = shfc->f[i]; } /* When we had particle decomposition, this code only worked with * PD when all particles involved with each shell were in the same * charge group. Not sure if this is still relevant. */ if (bDoNS && inputrec->ePBC != epbcNONE && !DOMAINDECOMP(cr)) { /* This is the only time where the coordinates are used * before do_force is called, which normally puts all * charge groups in the box. */ cg0 = 0; cg1 = top->cgs.nr; put_charge_groups_in_box(fplog, cg0, cg1, fr->ePBC, state->box, &(top->cgs), state->x, fr->cg_cm); if (graph) { mk_mshift(fplog, graph, fr->ePBC, state->box, state->x); } } /* After this all coordinate arrays will contain whole molecules */ if (graph) { shift_self(graph, state->box, state->x); } if (nflexcon) { if (nat > shfc->flex_nalloc) { shfc->flex_nalloc = over_alloc_dd(nat); srenew(shfc->acc_dir, shfc->flex_nalloc); srenew(shfc->x_old, shfc->flex_nalloc); } acc_dir = shfc->acc_dir; x_old = shfc->x_old; for (i = 0; i < homenr; i++) { for (d = 0; d < DIM; d++) { shfc->x_old[i][d] = state->x[start+i][d] - state->v[start+i][d]*inputrec->delta_t; } } } /* Do a prediction of the shell positions */ if (shfc->bPredict && !bCont) { predict_shells(fplog, state->x, state->v, inputrec->delta_t, nshell, shell, md->massT, NULL, bInit); } /* do_force expected the charge groups to be in the box */ if (graph) { unshift_self(graph, state->box, state->x); } /* Calculate the forces first time around */ if (gmx_debug_at) { pr_rvecs(debug, 0, "x b4 do_force", state->x + start, homenr); } do_force(fplog, cr, inputrec, mdstep, nrnb, wcycle, top, groups, state->box, state->x, &state->hist, force[Min], force_vir, md, enerd, fcd, state->lambda, graph, fr, vsite, mu_tot, t, fp_field, NULL, bBornRadii, (bDoNS ? GMX_FORCE_NS : 0) | force_flags); sf_dir = 0; if (nflexcon) { init_adir(fplog, shfc, constr, idef, inputrec, cr, dd_ac1, mdstep, md, start, end, shfc->x_old-start, state->x, state->x, force[Min], shfc->acc_dir-start, fr->bMolPBC, state->box, state->lambda, &dum, nrnb); for (i = start; i < end; i++) { sf_dir += md->massT[i]*norm2(shfc->acc_dir[i-start]); } } Epot[Min] = enerd->term[F_EPOT]; df[Min] = rms_force(cr, shfc->f[Min], nshell, shell, nflexcon, &sf_dir, &Epot[Min]); df[Try] = 0; if (debug) { fprintf(debug, "df = %g %g\n", df[Min], df[Try]); } if (gmx_debug_at) { pr_rvecs(debug, 0, "force0", force[Min], md->nr); } if (nshell+nflexcon > 0) { /* Copy x to pos[Min] & pos[Try]: during minimization only the * shell positions are updated, therefore the other particles must * be set here. */ memcpy(pos[Min], state->x, nat*sizeof(state->x[0])); memcpy(pos[Try], state->x, nat*sizeof(state->x[0])); } if (bVerbose && MASTER(cr)) { print_epot(stdout, mdstep, 0, Epot[Min], df[Min], nflexcon, sf_dir); } if (debug) { fprintf(debug, "%17s: %14.10e\n", interaction_function[F_EKIN].longname, enerd->term[F_EKIN]); fprintf(debug, "%17s: %14.10e\n", interaction_function[F_EPOT].longname, enerd->term[F_EPOT]); fprintf(debug, "%17s: %14.10e\n", interaction_function[F_ETOT].longname, enerd->term[F_ETOT]); fprintf(debug, "SHELLSTEP %s\n", gmx_step_str(mdstep, sbuf)); } /* First check whether we should do shells, or whether the force is * low enough even without minimization. */ *bConverged = (df[Min] < ftol); for (count = 1; (!(*bConverged) && (count < number_steps)); count++) { if (vsite) { construct_vsites(vsite, pos[Min], inputrec->delta_t, state->v, idef->iparams, idef->il, fr->ePBC, fr->bMolPBC, cr, state->box); } if (nflexcon) { init_adir(fplog, shfc, constr, idef, inputrec, cr, dd_ac1, mdstep, md, start, end, x_old-start, state->x, pos[Min], force[Min], acc_dir-start, fr->bMolPBC, state->box, state->lambda, &dum, nrnb); directional_sd(pos[Min], pos[Try], acc_dir-start, start, end, fr->fc_stepsize); } /* New positions, Steepest descent */ shell_pos_sd(pos[Min], pos[Try], force[Min], nshell, shell, count); /* do_force expected the charge groups to be in the box */ if (graph) { unshift_self(graph, state->box, pos[Try]); } if (gmx_debug_at) { pr_rvecs(debug, 0, "RELAX: pos[Min] ", pos[Min] + start, homenr); pr_rvecs(debug, 0, "RELAX: pos[Try] ", pos[Try] + start, homenr); } /* Try the new positions */ do_force(fplog, cr, inputrec, 1, nrnb, wcycle, top, groups, state->box, pos[Try], &state->hist, force[Try], force_vir, md, enerd, fcd, state->lambda, graph, fr, vsite, mu_tot, t, fp_field, NULL, bBornRadii, force_flags); if (gmx_debug_at) { pr_rvecs(debug, 0, "RELAX: force[Min]", force[Min] + start, homenr); pr_rvecs(debug, 0, "RELAX: force[Try]", force[Try] + start, homenr); } sf_dir = 0; if (nflexcon) { init_adir(fplog, shfc, constr, idef, inputrec, cr, dd_ac1, mdstep, md, start, end, x_old-start, state->x, pos[Try], force[Try], acc_dir-start, fr->bMolPBC, state->box, state->lambda, &dum, nrnb); for (i = start; i < end; i++) { sf_dir += md->massT[i]*norm2(acc_dir[i-start]); } } Epot[Try] = enerd->term[F_EPOT]; df[Try] = rms_force(cr, force[Try], nshell, shell, nflexcon, &sf_dir, &Epot[Try]); if (debug) { fprintf(debug, "df = %g %g\n", df[Min], df[Try]); } if (debug) { if (gmx_debug_at) { pr_rvecs(debug, 0, "F na do_force", force[Try] + start, homenr); } if (gmx_debug_at) { fprintf(debug, "SHELL ITER %d\n", count); dump_shells(debug, pos[Try], force[Try], ftol, nshell, shell); } } if (bVerbose && MASTER(cr)) { print_epot(stdout, mdstep, count, Epot[Try], df[Try], nflexcon, sf_dir); } *bConverged = (df[Try] < ftol); if ((df[Try] < df[Min])) { if (debug) { fprintf(debug, "Swapping Min and Try\n"); } if (nflexcon) { /* Correct the velocities for the flexible constraints */ invdt = 1/inputrec->delta_t; for (i = start; i < end; i++) { for (d = 0; d < DIM; d++) { state->v[i][d] += (pos[Try][i][d] - pos[Min][i][d])*invdt; } } } Min = Try; } else { decrease_step_size(nshell, shell); } } if (MASTER(cr) && !(*bConverged)) { /* Note that the energies and virial are incorrect when not converged */ if (fplog) { fprintf(fplog, "step %s: EM did not converge in %d iterations, RMS force %.3f\n", gmx_step_str(mdstep, sbuf), number_steps, df[Min]); } fprintf(stderr, "step %s: EM did not converge in %d iterations, RMS force %.3f\n", gmx_step_str(mdstep, sbuf), number_steps, df[Min]); } /* Copy back the coordinates and the forces */ memcpy(state->x, pos[Min], nat*sizeof(state->x[0])); memcpy(f, force[Min], nat*sizeof(f[0])); return count; }
void set_lincs(t_idef *idef,t_mdatoms *md, gmx_bool bDynamics,t_commrec *cr, struct gmx_lincsdata *li) { int start,natoms,nflexcon; t_blocka at2con; t_iatom *iatom; int i,k,ncc_alloc,ni,con,nconnect,concon; int type,a1,a2; real lenA=0,lenB; gmx_bool bLocal; li->nc = 0; li->ncc = 0; /* This is the local topology, so there are only F_CONSTR constraints */ if (idef->il[F_CONSTR].nr == 0) { /* There are no constraints, * we do not need to fill any data structures. */ return; } if (debug) { fprintf(debug,"Building the LINCS connectivity\n"); } if (DOMAINDECOMP(cr)) { if (cr->dd->constraints) { dd_get_constraint_range(cr->dd,&start,&natoms); } else { natoms = cr->dd->nat_home; } start = 0; } else if(PARTDECOMP(cr)) { pd_get_constraint_range(cr->pd,&start,&natoms); } else { start = md->start; natoms = md->homenr; } at2con = make_at2con(start,natoms,idef->il,idef->iparams,bDynamics, &nflexcon); if (idef->il[F_CONSTR].nr/3 > li->nc_alloc || li->nc_alloc == 0) { li->nc_alloc = over_alloc_dd(idef->il[F_CONSTR].nr/3); srenew(li->bllen0,li->nc_alloc); srenew(li->ddist,li->nc_alloc); srenew(li->bla,2*li->nc_alloc); srenew(li->blc,li->nc_alloc); srenew(li->blc1,li->nc_alloc); srenew(li->blnr,li->nc_alloc+1); srenew(li->bllen,li->nc_alloc); srenew(li->tmpv,li->nc_alloc); srenew(li->tmp1,li->nc_alloc); srenew(li->tmp2,li->nc_alloc); srenew(li->tmp3,li->nc_alloc); srenew(li->lambda,li->nc_alloc); if (li->ncg_triangle > 0) { /* This is allocating too much, but it is difficult to improve */ srenew(li->triangle,li->nc_alloc); srenew(li->tri_bits,li->nc_alloc); } } iatom = idef->il[F_CONSTR].iatoms; ncc_alloc = li->ncc_alloc; li->blnr[0] = 0; ni = idef->il[F_CONSTR].nr/3; con = 0; nconnect = 0; li->blnr[con] = nconnect; for(i=0; i<ni; i++) { bLocal = TRUE; type = iatom[3*i]; a1 = iatom[3*i+1]; a2 = iatom[3*i+2]; lenA = idef->iparams[type].constr.dA; lenB = idef->iparams[type].constr.dB; /* Skip the flexible constraints when not doing dynamics */ if (bDynamics || lenA!=0 || lenB!=0) { li->bllen0[con] = lenA; li->ddist[con] = lenB - lenA; /* Set the length to the topology A length */ li->bllen[con] = li->bllen0[con]; li->bla[2*con] = a1; li->bla[2*con+1] = a2; /* Construct the constraint connection matrix blbnb */ for(k=at2con.index[a1-start]; k<at2con.index[a1-start+1]; k++) { concon = at2con.a[k]; if (concon != i) { if (nconnect >= ncc_alloc) { ncc_alloc = over_alloc_small(nconnect+1); srenew(li->blbnb,ncc_alloc); } li->blbnb[nconnect++] = concon; } } for(k=at2con.index[a2-start]; k<at2con.index[a2-start+1]; k++) { concon = at2con.a[k]; if (concon != i) { if (nconnect+1 > ncc_alloc) { ncc_alloc = over_alloc_small(nconnect+1); srenew(li->blbnb,ncc_alloc); } li->blbnb[nconnect++] = concon; } } li->blnr[con+1] = nconnect; if (cr->dd == NULL) { /* Order the blbnb matrix to optimize memory access */ qsort(&(li->blbnb[li->blnr[con]]),li->blnr[con+1]-li->blnr[con], sizeof(li->blbnb[0]),int_comp); } /* Increase the constraint count */ con++; } } done_blocka(&at2con); /* This is the real number of constraints, * without dynamics the flexible constraints are not present. */ li->nc = con; li->ncc = li->blnr[con]; if (cr->dd == NULL) { /* Since the matrix is static, we can free some memory */ ncc_alloc = li->ncc; srenew(li->blbnb,ncc_alloc); } if (ncc_alloc > li->ncc_alloc) { li->ncc_alloc = ncc_alloc; srenew(li->blmf,li->ncc_alloc); srenew(li->blmf1,li->ncc_alloc); srenew(li->tmpncc,li->ncc_alloc); } if (debug) { fprintf(debug,"Number of constraints is %d, couplings %d\n", li->nc,li->ncc); } set_lincs_matrix(li,md->invmass,md->lambda); }