void calc_compact_unitcell_vertices(int ecenter, const matrix box, rvec vert[]) { rvec img[NTRICIMG], box_center; int n, i, j, tmp[4], d; const real oneFourth = 0.25; calc_triclinic_images(box, img); n = 0; for (i = 2; i <= 5; i += 3) { tmp[0] = i-1; if (i == 2) { tmp[1] = 8; } else { tmp[1] = 6; } tmp[2] = (i+1) % 6; tmp[3] = tmp[1]+4; for (j = 0; j < 4; j++) { for (d = 0; d < DIM; d++) { vert[n][d] = img[i][d]+img[tmp[j]][d]+img[tmp[(j+1)%4]][d]; } n++; } } for (i = 7; i <= 13; i += 6) { tmp[0] = (i-7)/2; tmp[1] = tmp[0]+1; if (i == 7) { tmp[2] = 8; } else { tmp[2] = 10; } tmp[3] = i-1; for (j = 0; j < 4; j++) { for (d = 0; d < DIM; d++) { vert[n][d] = img[i][d]+img[tmp[j]][d]+img[tmp[(j+1)%4]][d]; } n++; } } for (i = 9; i <= 11; i += 2) { if (i == 9) { tmp[0] = 3; } else { tmp[0] = 0; } tmp[1] = tmp[0]+1; if (i == 9) { tmp[2] = 6; } else { tmp[2] = 12; } tmp[3] = i-1; for (j = 0; j < 4; j++) { for (d = 0; d < DIM; d++) { vert[n][d] = img[i][d]+img[tmp[j]][d]+img[tmp[(j+1)%4]][d]; } n++; } } calc_box_center(ecenter, box, box_center); for (i = 0; i < NCUCVERT; i++) { for (d = 0; d < DIM; d++) { vert[i][d] = vert[i][d]*oneFourth+box_center[d]; } } }
int gmx_editconf(int argc, char *argv[]) { static char *desc[] = { "editconf converts generic structure format to [TT].gro[tt], [TT].g96[tt]", "or [TT].pdb[tt].", "[PAR]", "The box can be modified with options [TT]-box[tt], [TT]-d[tt] and", "[TT]-angles[tt]. Both [TT]-box[tt] and [TT]-d[tt]", "will center the system in the box, unless [TT]-noc[tt] is used.", "[PAR]", "Option [TT]-bt[tt] determines the box type: [TT]triclinic[tt] is a", "triclinic box, [TT]cubic[tt] is a rectangular box with all sides equal", "[TT]dodecahedron[tt] represents a rhombic dodecahedron and " "[TT]octahedron[tt] is a truncated octahedron.", "The last two are special cases of a triclinic box.", "The length of the three box vectors of the truncated octahedron is the", "shortest distance between two opposite hexagons.", "The volume of a dodecahedron is 0.71 and that of a truncated octahedron", "is 0.77 of that of a cubic box with the same periodic image distance.", "[PAR]", "Option [TT]-box[tt] requires only", "one value for a cubic box, dodecahedron and a truncated octahedron.", "[PAR]", "With [TT]-d[tt] and a [TT]triclinic[tt] box the size of the system in the x, y", "and z directions is used. With [TT]-d[tt] and [TT]cubic[tt],", "[TT]dodecahedron[tt] or [TT]octahedron[tt] boxes, the dimensions are set", "to the diameter of the system (largest distance between atoms) plus twice", "the specified distance.", "[PAR]", "Option [TT]-angles[tt] is only meaningful with option [TT]-box[tt] and", "a triclinic box and can not be used with option [TT]-d[tt].", "[PAR]", "When [TT]-n[tt] or [TT]-ndef[tt] is set, a group", "can be selected for calculating the size and the geometric center,", "otherwise the whole system is used.", "[PAR]", "[TT]-rotate[tt] rotates the coordinates and velocities.", "[PAR]", "[TT]-princ[tt] aligns the principal axes of the system along the", "coordinate axes, this may allow you to decrease the box volume,", "but beware that molecules can rotate significantly in a nanosecond.", "[PAR]", "Scaling is applied before any of the other operations are", "performed. Boxes and coordinates can be scaled to give a certain density (option", "[TT]-density[tt]). Note that this may be inaccurate in case a gro", "file is given as input. A special feature of the scaling option, when the", "factor -1 is given in one dimension, one obtains a mirror image,", "mirrored in one of the plains, when one uses -1 in three dimensions", "a point-mirror image is obtained.[PAR]", "Groups are selected after all operations have been applied.[PAR]", "Periodicity can be removed in a crude manner.", "It is important that the box sizes at the bottom of your input file", "are correct when the periodicity is to be removed.", "[PAR]", "When writing [TT].pdb[tt] files, B-factors can be", "added with the [TT]-bf[tt] option. B-factors are read", "from a file with with following format: first line states number of", "entries in the file, next lines state an index", "followed by a B-factor. The B-factors will be attached per residue", "unless an index is larger than the number of residues or unless the", "[TT]-atom[tt] option is set. Obviously, any type of numeric data can", "be added instead of B-factors. [TT]-legend[tt] will produce", "a row of CA atoms with B-factors ranging from the minimum to the", "maximum value found, effectively making a legend for viewing.", "[PAR]", "With the option -mead a special pdb (pqr) file for the MEAD electrostatics", "program (Poisson-Boltzmann solver) can be made. A further prerequisite", "is that the input file is a run input file.", "The B-factor field is then filled with the Van der Waals radius", "of the atoms while the occupancy field will hold the charge.", "[PAR]", "The option -grasp is similar, but it puts the charges in the B-factor", "and the radius in the occupancy.", "[PAR]", "Finally with option [TT]-label[tt] editconf can add a chain identifier", "to a pdb file, which can be useful for analysis with e.g. rasmol." "[PAR]", "To convert a truncated octrahedron file produced by a package which uses", "a cubic box with the corners cut off (such as Gromos) use:[BR]", "[TT]editconf -f <in> -rotate 0 45 35.264 -bt o -box <veclen> -o <out>[tt][BR]", "where [TT]veclen[tt] is the size of the cubic box times sqrt(3)/2." }; static char *bugs[] = { "For complex molecules, the periodicity removal routine may break down, " "in that case you can use trjconv" }; static real dist=0.0,rbox=0.0,to_diam=0.0; static bool bNDEF=FALSE,bRMPBC=FALSE,bCenter=FALSE,bReadVDW=FALSE; static bool peratom=FALSE,bLegend=FALSE,bOrient=FALSE,bMead=FALSE,bGrasp=FALSE,bSig56=FALSE; static rvec scale={1,1,1},newbox={0,0,0},newang={90,90,90}; static real rho=1000.0,rvdw=0.12; static rvec center={0,0,0},translation={0,0,0},rotangles={0,0,0}; static char *btype[]={ NULL, "triclinic", "cubic", "dodecahedron", "octahedron", NULL },*label="A"; static rvec visbox={0,0,0}; t_pargs pa[] = { { "-ndef", FALSE, etBOOL, {&bNDEF}, "Choose output from default index groups" }, { "-visbox", FALSE, etRVEC, {visbox}, "HIDDENVisualize a grid of boxes, -1 visualizes the 14 box images" }, { "-bt", FALSE, etENUM, {btype}, "Box type for -box and -d" }, { "-box", FALSE, etRVEC, {newbox}, "Box vector lengths (a,b,c)" }, { "-angles", FALSE, etRVEC, {newang}, "Angles between the box vectors (bc,ac,ab)" }, { "-d", FALSE, etREAL, {&dist}, "Distance between the solute and the box" }, { "-c", FALSE, etBOOL, {&bCenter}, "Center molecule in box (implied by -box and -d)" }, { "-center", FALSE, etRVEC, {center}, "Coordinates of geometrical center"}, { "-translate", FALSE, etRVEC, {translation}, "Translation" }, { "-rotate", FALSE, etRVEC, {rotangles}, "Rotation around the X, Y and Z axes in degrees" }, { "-princ", FALSE, etBOOL, {&bOrient}, "Orient molecule(s) along their principal axes" }, { "-scale", FALSE, etRVEC, {scale}, "Scaling factor" }, { "-density",FALSE, etREAL, {&rho}, "Density (g/l) of the output box achieved by scaling" }, { "-pbc", FALSE, etBOOL, {&bRMPBC}, "Remove the periodicity (make molecule whole again)" }, { "-grasp", FALSE, etBOOL, {&bGrasp}, "Store the charge of the atom in the B-factor field and the radius of the atom in the occupancy field" }, { "-rvdw", FALSE, etREAL, {&rvdw}, "Default Van der Waals radius (in nm) if one can not be found in the database or if no parameters are present in the topology file" }, { "-sig56", FALSE, etREAL, {&bSig56}, "Use rmin/2 (minimum in the Van der Waals potential) rather than sigma/2 " }, { "-vdwread",FALSE, etBOOL, {&bReadVDW}, "Read the Van der Waals radii from the file vdwradii.dat rather than computing the radii based on the force field" }, { "-atom", FALSE, etBOOL, {&peratom}, "Force B-factor attachment per atom" }, { "-legend", FALSE, etBOOL, {&bLegend}, "Make B-factor legend" }, { "-label", FALSE, etSTR, {&label}, "Add chain label for all residues" } }; #define NPA asize(pa) FILE *out; char *infile,*outfile,title[STRLEN]; int outftp,natom,i,j,n_bfac,itype,ntype; double *bfac=NULL,c6,c12; int *bfac_nr=NULL; t_topology *top; t_atoms atoms; char *grpname,*sgrpname; int isize,ssize,tsize; atom_id *index,*sindex,*tindex; rvec *x,*v,gc,min,max,size; int ePBC; matrix box; bool bIndex,bSetSize,bSetAng,bCubic,bDist,bSetCenter; bool bHaveV,bScale,bRho,bTranslate,bRotate,bCalcGeom,bCalcDiam; real xs,ys,zs,xcent,ycent,zcent,diam=0,mass=0,d,vdw; gmx_atomprop_t aps; t_filenm fnm[] = { { efSTX, "-f", NULL, ffREAD }, { efNDX, "-n", NULL, ffOPTRD }, { efSTO, NULL, NULL, ffOPTWR }, { efPQR, "-mead", "mead", ffOPTWR }, { efDAT, "-bf", "bfact", ffOPTRD } }; #define NFILE asize(fnm) CopyRight(stderr,argv[0]); parse_common_args(&argc,argv,PCA_CAN_VIEW,NFILE,fnm,NPA,pa, asize(desc),desc,asize(bugs),bugs); bIndex = opt2bSet("-n",NFILE,fnm) || bNDEF; bMead = opt2bSet("-mead",NFILE,fnm); bSetSize = opt2parg_bSet("-box" ,NPA,pa); bSetAng = opt2parg_bSet("-angles" ,NPA,pa); bSetCenter= opt2parg_bSet("-center" ,NPA,pa); bDist = opt2parg_bSet("-d" ,NPA,pa); /* Only automatically turn on centering without -noc */ if ((bDist || bSetSize || bSetCenter) && !opt2parg_bSet("-c",NPA,pa)) { bCenter = TRUE; } bScale = opt2parg_bSet("-scale" ,NPA,pa); bRho = opt2parg_bSet("-density",NPA,pa); bTranslate= opt2parg_bSet("-translate",NPA,pa); bRotate = opt2parg_bSet("-rotate",NPA,pa); if (bScale && bRho) fprintf(stderr,"WARNING: setting -density overrides -scale\n"); bScale = bScale || bRho; bCalcGeom = bCenter || bRotate || bOrient || bScale; bCalcDiam = btype[0][0]=='c' || btype[0][0]=='d' || btype[0][0]=='o'; infile = ftp2fn(efSTX,NFILE,fnm); if (bMead) outfile = ftp2fn(efPQR,NFILE,fnm); else outfile = ftp2fn(efSTO,NFILE,fnm); outftp = fn2ftp(outfile); aps = gmx_atomprop_init(); if (bMead && bGrasp) { printf("Incompatible options -mead and -grasp. Turning off -grasp\n"); bGrasp = FALSE; } if (bGrasp && (outftp != efPDB)) gmx_fatal(FARGS,"Output file should be a .pdb file" " when using the -grasp option\n"); if ((bMead || bGrasp) && !((fn2ftp(infile) == efTPR) || (fn2ftp(infile) == efTPA) || (fn2ftp(infile) == efTPB))) gmx_fatal(FARGS,"Input file should be a .tp[abr] file" " when using the -mead option\n"); get_stx_coordnum(infile,&natom); init_t_atoms(&atoms,natom,TRUE); snew(x,natom); snew(v,natom); read_stx_conf(infile,title,&atoms,x,v,&ePBC,box); printf("Read %d atoms\n",atoms.nr); if (ePBC != epbcNONE) { real vol = det(box); printf("Volume: %g nm^3, corresponds to roughly %d electrons\n", vol,100*((int)(vol*4.5))); } if (bMead || bGrasp) { top = read_top(infile,NULL); if (atoms.nr != top->atoms.nr) gmx_fatal(FARGS,"Atom numbers don't match (%d vs. %d)", atoms.nr,top->atoms.nr); snew(atoms.pdbinfo,top->atoms.nr); ntype = top->idef.atnr; for(i=0; (i<atoms.nr); i++) { /* Determine the Van der Waals radius from the force field */ if (bReadVDW) { if (!gmx_atomprop_query(aps,epropVDW, *top->atoms.resname[top->atoms.atom[i].resnr], *top->atoms.atomname[i],&vdw)) vdw = rvdw; } else { itype = top->atoms.atom[i].type; c12 = top->idef.iparams[itype*ntype+itype].lj.c12; c6 = top->idef.iparams[itype*ntype+itype].lj.c6; if ((c6 != 0) && (c12 != 0)) { real sig6; if (bSig56) sig6 = 2*c12/c6; else sig6 = c12/c6; vdw = 0.5*pow(sig6,1.0/6.0); } else vdw = rvdw; } /* Factor of 10 for nm -> Angstroms */ vdw *= 10; if (bMead) { atoms.pdbinfo[i].occup = top->atoms.atom[i].q; atoms.pdbinfo[i].bfac = vdw; } else { atoms.pdbinfo[i].occup = vdw; atoms.pdbinfo[i].bfac = top->atoms.atom[i].q; } } } bHaveV=FALSE; for (i=0; (i<natom) && !bHaveV; i++) for (j=0; (j<DIM) && !bHaveV; j++) bHaveV=bHaveV || (v[i][j]!=0); printf("%selocities found\n",bHaveV?"V":"No v"); if (visbox[0] > 0) { if (bIndex) gmx_fatal(FARGS,"Sorry, can not visualize box with index groups"); if (outftp != efPDB) gmx_fatal(FARGS,"Sorry, can only visualize box with a pdb file"); } else if (visbox[0] == -1) visualize_images("images.pdb",ePBC,box); /* remove pbc */ if (bRMPBC) rm_gropbc(&atoms,x,box); if (bCalcGeom) { if (bIndex) { fprintf(stderr,"\nSelect a group for determining the system size:\n"); get_index(&atoms,ftp2fn_null(efNDX,NFILE,fnm), 1,&ssize,&sindex,&sgrpname); } else { ssize = atoms.nr; sindex = NULL; } diam=calc_geom(ssize,sindex,x,gc,min,max,bCalcDiam); rvec_sub(max, min, size); printf(" system size :%7.3f%7.3f%7.3f (nm)\n", size[XX], size[YY], size[ZZ]); if (bCalcDiam) printf(" diameter :%7.3f (nm)\n",diam); printf(" center :%7.3f%7.3f%7.3f (nm)\n", gc[XX], gc[YY], gc[ZZ]); printf(" box vectors :%7.3f%7.3f%7.3f (nm)\n", norm(box[XX]), norm(box[YY]), norm(box[ZZ])); printf(" box angles :%7.2f%7.2f%7.2f (degrees)\n", norm2(box[ZZ])==0 ? 0 : RAD2DEG*acos(cos_angle_no_table(box[YY],box[ZZ])), norm2(box[ZZ])==0 ? 0 : RAD2DEG*acos(cos_angle_no_table(box[XX],box[ZZ])), norm2(box[YY])==0 ? 0 : RAD2DEG*acos(cos_angle_no_table(box[XX],box[YY]))); printf(" box volume :%7.2f (nm^3)\n",det(box)); } if (bRho || bOrient) mass = calc_mass(&atoms,!fn2bTPX(infile),aps); if (bOrient) { atom_id *index; char *grpnames; /* Get a group for principal component analysis */ fprintf(stderr,"\nSelect group for the determining the orientation\n"); get_index(&atoms,ftp2fn_null(efNDX,NFILE,fnm),1,&isize,&index,&grpnames); /* Orient the principal axes along the coordinate axes */ orient_princ(&atoms,isize,index,natom,x,bHaveV ? v : NULL, NULL); sfree(index); sfree(grpnames); } if ( bScale ) { /* scale coordinates and box */ if (bRho) { /* Compute scaling constant */ real vol,dens; vol = det(box); dens = (mass*AMU)/(vol*NANO*NANO*NANO); fprintf(stderr,"Volume of input %g (nm^3)\n",vol); fprintf(stderr,"Mass of input %g (a.m.u.)\n",mass); fprintf(stderr,"Density of input %g (g/l)\n",dens); if (vol==0 || mass==0) gmx_fatal(FARGS,"Cannot scale density with " "zero mass (%g) or volume (%g)\n",mass,vol); scale[XX] = scale[YY] = scale[ZZ] = pow(dens/rho,1.0/3.0); fprintf(stderr,"Scaling all box vectors by %g\n",scale[XX]); } scale_conf(atoms.nr,x,box,scale); } if (bTranslate) { if (bIndex) { fprintf(stderr,"\nSelect a group that you want to translate:\n"); get_index(&atoms,ftp2fn_null(efNDX,NFILE,fnm), 1,&ssize,&sindex,&sgrpname); } else { ssize = atoms.nr; sindex = NULL; } printf("Translating %d atoms (out of %d) by %g %g %g nm\n",ssize,natom, translation[XX],translation[YY],translation[ZZ]); if (sindex) { for(i=0; i<ssize; i++) rvec_inc(x[sindex[i]],translation); } else { for(i=0; i<natom; i++) rvec_inc(x[i],translation); } } if (bRotate) { /* Rotate */ printf("Rotating %g, %g, %g degrees around the X, Y and Z axis respectively\n",rotangles[XX],rotangles[YY],rotangles[ZZ]); for(i=0; i<DIM; i++) rotangles[i] *= DEG2RAD; rotate_conf(natom,x,v,rotangles[XX],rotangles[YY],rotangles[ZZ]); } if (bCalcGeom) { /* recalc geometrical center and max and min coordinates and size */ calc_geom(ssize,sindex,x,gc,min,max,FALSE); rvec_sub(max, min, size); if (bScale || bOrient || bRotate) printf("new system size : %6.3f %6.3f %6.3f\n", size[XX],size[YY],size[ZZ]); } if (bSetSize || bDist || (btype[0][0]=='t' && bSetAng)) { ePBC = epbcXYZ; if (!(bSetSize || bDist)) for (i=0; i<DIM; i++) newbox[i] = norm(box[i]); clear_mat(box); /* calculate new boxsize */ switch(btype[0][0]){ case 't': if (bDist) for(i=0; i<DIM; i++) newbox[i] = size[i]+2*dist; if (!bSetAng) { box[XX][XX] = newbox[XX]; box[YY][YY] = newbox[YY]; box[ZZ][ZZ] = newbox[ZZ]; } else { svmul(DEG2RAD,newang,newang); box[XX][XX] = newbox[XX]; box[YY][XX] = newbox[YY]*cos(newang[ZZ]); box[YY][YY] = newbox[YY]*sin(newang[ZZ]); box[ZZ][XX] = newbox[ZZ]*cos(newang[YY]); box[ZZ][YY] = newbox[ZZ] *(cos(newang[XX])-cos(newang[YY])*cos(newang[ZZ]))/sin(newang[ZZ]); box[ZZ][ZZ] = sqrt(sqr(newbox[ZZ]) -box[ZZ][XX]*box[ZZ][XX]-box[ZZ][YY]*box[ZZ][YY]); } break; case 'c': case 'd': case 'o': if (bSetSize) d = newbox[0]; else d = diam+2*dist; if (btype[0][0] == 'c') for(i=0; i<DIM; i++) box[i][i] = d; else if (btype[0][0] == 'd') { box[XX][XX] = d; box[YY][YY] = d; box[ZZ][XX] = d/2; box[ZZ][YY] = d/2; box[ZZ][ZZ] = d*sqrt(2)/2; } else { box[XX][XX] = d; box[YY][XX] = d/3; box[YY][YY] = d*sqrt(2)*2/3; box[ZZ][XX] = -d/3; box[ZZ][YY] = d*sqrt(2)/3; box[ZZ][ZZ] = d*sqrt(6)/3; } break; } } /* calculate new coords for geometrical center */ if (!bSetCenter) calc_box_center(ecenterDEF,box,center); /* center molecule on 'center' */ if (bCenter) center_conf(natom,x,center,gc); /* print some */ if (bCalcGeom) { calc_geom(ssize,sindex,x, gc, min, max, FALSE); printf("new center :%7.3f%7.3f%7.3f (nm)\n",gc[XX],gc[YY],gc[ZZ]); } if (bOrient || bScale || bDist || bSetSize) { printf("new box vectors :%7.3f%7.3f%7.3f (nm)\n", norm(box[XX]), norm(box[YY]), norm(box[ZZ])); printf("new box angles :%7.2f%7.2f%7.2f (degrees)\n", norm2(box[ZZ])==0 ? 0 : RAD2DEG*acos(cos_angle_no_table(box[YY],box[ZZ])), norm2(box[ZZ])==0 ? 0 : RAD2DEG*acos(cos_angle_no_table(box[XX],box[ZZ])), norm2(box[YY])==0 ? 0 : RAD2DEG*acos(cos_angle_no_table(box[XX],box[YY]))); printf("new box volume :%7.2f (nm^3)\n",det(box)); } if (check_box(epbcXYZ,box)) printf("\nWARNING: %s\n",check_box(epbcXYZ,box)); if (bDist && btype[0][0]=='t') { if(TRICLINIC(box)) { printf("\nWARNING: Your box is triclinic with non-orthogonal axes. In this case, the\n" "distance from the solute to a box surface along the corresponding normal\n" "vector might be somewhat smaller than your specified value %f.\n" "You can check the actual value with g_mindist -pi\n",dist); } else { printf("\nWARNING: No boxtype specified - distance condition applied in each dimension.\n" "If the molecule rotates the actual distance will be smaller. You might want\n" "to use a cubic box instead, or why not try a dodecahedron today?\n"); } } if (bIndex) { fprintf(stderr,"\nSelect a group for output:\n"); get_index(&atoms,opt2fn_null("-n",NFILE,fnm), 1,&isize,&index,&grpname); if (opt2bSet("-bf",NFILE,fnm)) gmx_fatal(FARGS,"combination not implemented: -bf -n or -bf -ndef"); else write_sto_conf_indexed(outfile,title,&atoms,x,bHaveV?v:NULL,ePBC,box, isize,index); } else { if ((outftp == efPDB) || (outftp == efPQR)) { out=ffopen(outfile,"w"); if (bMead) { set_pdb_wide_format(TRUE); fprintf(out,"REMARK " "The B-factors in this file hold atomic radii\n"); fprintf(out,"REMARK " "The occupancy in this file hold atomic charges\n"); } else if (bGrasp) { fprintf(out,"GRASP PDB FILE\nFORMAT NUMBER=1\n"); fprintf(out,"REMARK " "The B-factors in this file hold atomic charges\n"); fprintf(out,"REMARK " "The occupancy in this file hold atomic radii\n"); } else if (opt2bSet("-bf",NFILE,fnm)) { read_bfac(opt2fn("-bf",NFILE,fnm),&n_bfac,&bfac,&bfac_nr); set_pdb_conf_bfac(atoms.nr,atoms.nres,&atoms, n_bfac,bfac,bfac_nr,peratom); } if (opt2parg_bSet("-label",NPA,pa)) { for(i=0; (i<atoms.nr); i++) atoms.atom[i].chain=label[0]; } write_pdbfile(out,title,&atoms,x,ePBC,box,0,-1); if (bLegend) pdb_legend(out,atoms.nr,atoms.nres,&atoms,x); if (visbox[0] > 0) visualize_box(out,bLegend ? atoms.nr+12 : atoms.nr, bLegend? atoms.nres=12 : atoms.nres,box,visbox); fclose(out); } else write_sto_conf(outfile,title,&atoms,x,bHaveV?v:NULL,ePBC,box); } gmx_atomprop_destroy(aps); do_view(outfile,NULL); thanx(stderr); return 0; }