1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
3 * This source code is NOT REALLY part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
32 * Author: do_multiprot was written by Ran Friedman <r.friedman@bioc.uzh.ch>
35 * Green Red Orange Magenta Azure Cyan Skyblue
45 #include "gromacs/utility/smalloc.h"
47 #include "gromacs/commandline/pargs.h"
49 #include "gromacs/fileio/pdbio.h"
50 #include "gmx_fatal.h"
52 #include "gromacs/fileio/matio.h"
55 #include "gromacs/fileio/tpxio.h"
59 #include "gromacs/fileio/confio.h"
60 #include "gromacs/fileio/gmxfio.h"
67 static void process_multiprot_output(const char *fn, real *rmsd, int *nres, rvec rotangles,
68 rvec translation, bool bCountres, t_countres *countres)
77 mpoutput=gmx_ffopen (fn,"r");
81 fgets(line, 256, mpoutput);
82 } while (strstr(line,"Match List") == NULL);
83 fgets(line, 256, mpoutput);
85 string = strtok (line,".");
86 while (i<2 && string != NULL) {
87 string = strtok (NULL,". ");
93 while (countres[j].resnr!=res)
97 fgets(line, 256, mpoutput);
98 } while (strstr(line,"End of Match List") == NULL);
102 fgets(line, 256, mpoutput);
103 } while (strstr(line,"Trans : ") == NULL);
105 string = strtok (line," :");
106 string = strtok (NULL," ");
107 while (i<3 && string != NULL) {
108 string = strtok (NULL," ");
109 rotangles[i]=atof(string);
113 while (i<3 && string != NULL) {
114 string = strtok (NULL," ");
115 translation[i]=atof(string)/10;
119 gmx_warning("Not enough values for rotation and translation vectors in the output of multiprot");
122 rotangles[YY]=rotangles[YY]*(-1);
125 fgets(line, 256, mpoutput);
126 if (strstr(line,"RMSD : ") != NULL) {
127 string = strtok (line,":");
128 string = strtok (NULL,":");
129 (*rmsd)=atof(string)/10;
133 fgets(line,256, mpoutput);
134 if (strstr(line,"Match List") != NULL) {
135 string = strtok (line,":");
136 string = strtok (NULL,":");
137 (*nres) = atoi(string);
140 gmx_ffclose(mpoutput);
143 int main(int argc,char *argv[])
145 const char *desc[] = {
146 "[TT]do_multiprot[tt] ",
147 "reads a trajectory file and aligns it to a reference structure ",
149 "calling the multiprot program. This allows you to use a reference",
150 "structure whose sequence is different than that of the protein in the ",
151 "trajectory, since the alignment is based on the geometry, not sequence.",
152 "The output of [TT]do_multiprot[tt] includes the rmsd and the number of residues",
153 "on which it was calculated.",
155 "An aligned trajectory file is generated with the [TT]-ox[tt] option.[PAR]",
156 "With the [TT]-cr[tt] option, the number of hits in the alignment is given",
157 "per residue. This number can be between 0 and the number of frames, and",
158 "indicates the structural conservation of this residue.[PAR]",
159 "If you do not have the [TT]multiprot[tt] program, get it. [TT]do_multiprot[tt] assumes",
160 "that the [TT]multiprot[tt] executable is [TT]/usr/local/bin/multiprot[tt]. If this is ",
161 "not the case, then you should set an environment variable [BB]MULTIPROT[bb]",
162 "pointing to the [TT]multiprot[tt] executable, e.g.: [PAR]",
163 "[TT]setenv MULTIPROT /usr/MultiProtInstall/multiprot.Linux[tt][PAR]",
164 "Note that at the current implementation only binary alignment (your",
165 "molecule to a reference) is supported. In addition, note that the ",
166 "by default [TT]multiprot[tt] aligns the two proteins on their C-alpha carbons.",
167 "and that this depends on the [TT]multiprot[tt] parameters which are not dealt ",
168 "with here. Thus, the C-alpha carbons is expected to give the same "
169 "results as choosing the whole protein and will be slightly faster.[PAR]",
170 "For information about [TT]multiprot[tt], see:",
171 "http://bioinfo3d.cs.tau.ac.il/MultiProt/.[PAR]"
173 static bool bVerbose;
175 { "-v", FALSE, etBOOL, {&bVerbose},
176 "HIDDENGenerate miles of useless information" }
179 const char *bugs[] = {
180 "The program is very slow, since multiprot is run externally"
184 t_trxstatus *trxout=NULL;
185 FILE *tapein,*fo,*frc,*tmpf,*out=NULL,*fres=NULL;
187 const char *fn="2_sol.res";
190 t_atoms *atoms,ratoms,useatoms;
195 int i,j,natoms,nratoms,nframe=0,model_nr=-1;
196 int cur_res,prev_res;
198 t_countres *countres=NULL;
204 char pdbfile[32],refpdb[256],title[256],rtitle[256],filemode[5];
206 char multiprot[256],*mptr;
210 bool bTrjout,bCountres;
211 const char *TrjoutFile=NULL;
213 static rvec translation={0,0,0},rotangles={0,0,0};
214 gmx_rmpbc_t gpbc=NULL;
217 { efTRX, "-f", NULL, ffREAD },
218 { efTPS, NULL, NULL, ffREAD },
219 { efNDX, NULL, NULL, ffOPTRD },
220 { efSTX, "-r", NULL , ffREAD },
221 { efXVG, "-o", "rmss", ffWRITE },
222 { efXVG, "-rc", "rescount", ffWRITE},
223 { efXVG, "-cr", "countres", ffOPTWR},
224 { efTRX, "-ox", "aligned", ffOPTWR }
226 #define NFILE asize(fnm)
228 CopyRight(stderr,argv[0]);
229 parse_common_args(&argc,argv,PCA_CAN_TIME | PCA_CAN_VIEW | PCA_TIME_UNIT | PCA_BE_NICE ,
230 NFILE,fnm, asize(pa),pa, asize(desc),desc,
231 asize(bugs),bugs,&oenv
233 fnRef=opt2fn("-r",NFILE,fnm);
234 bTrjout = opt2bSet("-ox",NFILE,fnm);
235 bCountres= opt2bSet("-cr",NFILE,fnm);
238 TrjoutFile = opt2fn_null("-ox",NFILE,fnm);
241 read_tps_conf(ftp2fn(efTPS,NFILE,fnm),title,&top,&ePBC,&xp,NULL,box,FALSE);
242 gpbc = gmx_rmpbc_init(&top.idef,ePBC,top.atoms.nr,box);
247 get_stx_coordnum(fnRef,&nratoms);
248 init_t_atoms(&ratoms,nratoms,TRUE);
250 read_stx_conf(fnRef,rtitle,&ratoms,xr,NULL,&ePBC,rbox);
253 fprintf(stderr,"Read %d atoms\n",atoms->nr);
254 fprintf(stderr,"Read %d reference atoms\n",ratoms.nr);
257 snew(countres,ratoms.nres);
260 for (i=0;i<ratoms.nr;i++) {
262 cur_res=ratoms.atom[i].resind;
263 if (cur_res != prev_res) {
264 countres[j].resnr=cur_res;
270 get_index(atoms,ftp2fn_null(efNDX,NFILE,fnm),1,&gnx,&index,&grpnm);
273 for(i=0; (i<gnx); i++) {
274 if (atoms->atom[index[i]].resind != nr0) {
275 nr0=atoms->atom[index[i]].resind;
279 fprintf(stderr,"There are %d residues in your selected group\n",nres);
281 strcpy(pdbfile,"ddXXXXXX");
283 if ((tmpf = fopen(pdbfile,"w")) == NULL) {
284 sprintf(pdbfile,"%ctmp%cfilterXXXXXX",DIR_SEPARATOR,DIR_SEPARATOR);
286 if ((tmpf = fopen(pdbfile,"w")) == NULL) {
287 gmx_fatal(FARGS,"Can not open tmp file %s",pdbfile);
295 strcpy(refpdb,"ddXXXXXX");
297 strcat(refpdb,".pdb");
298 write_sto_conf(refpdb,rtitle,&ratoms,xr,NULL,ePBC,rbox);
301 strcpy(refpdb,fnRef);
304 if ((mptr=getenv("MULTIPROT")) == NULL) {
305 mptr="/usr/local/bin/multiprot";
307 if (!gmx_fexist(mptr)) {
308 gmx_fatal(FARGS,"MULTIPROT executable (%s) does not exist (use setenv MULTIPROT)",
311 sprintf (multiprot,"%s %s %s > /dev/null %s",
312 mptr, refpdb, pdbfile, "2> /dev/null");
315 fprintf(stderr,"multiprot cmd='%s'\n",multiprot);
317 if (!read_first_frame(oenv,&status,ftp2fn(efTRX,NFILE,fnm),&fr,TRX_READ_X))
318 gmx_fatal(FARGS,"Could not read a frame from %s",ftp2fn(efTRX,NFILE,fnm));
324 outftp=fn2ftp(TrjoutFile);
326 fprintf(stderr,"Will write %s: %s\n",ftp2ext(ftp),ftp2desc(outftp));
327 strcpy(filemode,"w");
334 trxout = open_trx(TrjoutFile,filemode);
339 /* Make atoms struct for output in GRO or PDB files */
340 /* get memory for stuff to go in pdb file */
341 init_t_atoms(&useatoms,nout,FALSE);
342 sfree(useatoms.resinfo);
343 useatoms.resinfo=atoms->resinfo;
344 for(i=0;(i<nout);i++) {
345 useatoms.atomname[i]=atoms->atomname[i];
346 useatoms.atom[i]=atoms->atom[i];
347 useatoms.nres=max(useatoms.nres,useatoms.atom[i].resind+1);
350 out=gmx_ffopen(TrjoutFile,filemode);
354 fprintf(gmx_fio_getfp(trx_get_fileio(trxout)),"Generated by %s. #atoms=%d, a BOX is"
355 " stored in this file.\n",ShortProgram(),nout);
358 if (natoms > atoms->nr) {
359 gmx_fatal(FARGS,"\nTrajectory does not match topology!");
362 gmx_fatal(FARGS,"\nTrajectory does not match selected group!");
365 fo = xvgropen(opt2fn("-o",NFILE,fnm),"RMSD","Time (ps)","RMSD (nm)",oenv);
366 frc = xvgropen(opt2fn("-rc",NFILE,fnm),"Number of Residues in the alignment","Time (ps)","Residues",oenv);
369 t = output_env_conv_time(oenv,fr.time);
370 gmx_rmpbc(gpbc,natoms,fr.box,fr.x);
371 tapein=gmx_ffopen(pdbfile,"w");
372 write_pdbfile_indexed(tapein,NULL,atoms,fr.x,ePBC,fr.box,' ',-1,gnx,index,NULL,TRUE);
376 process_multiprot_output(fn, &rmsd, &nres2,rotangles,translation,bCountres,countres);
377 fprintf(fo,"%12.7f",t);
378 fprintf(fo," %12.7f\n",rmsd);
379 fprintf(frc,"%12.7f",t);
380 fprintf(frc,"%12d\n",nres2);
382 rotate_conf(natoms,fr.x,NULL,rotangles[XX],rotangles[YY],rotangles[ZZ]);
383 for(i=0; i<natoms; i++) {
384 rvec_inc(fr.x[i],translation);
391 write_trxframe(trxout,&fr,NULL);
396 sprintf(out_title,"Generated by do_multiprot : %s t= %g %s",
397 title,output_env_conv_time(oenv,fr.time),output_env_get_time_unit(oenv));
400 write_hconf_p(out,out_title,&useatoms,prec2ndec(fr.prec),
404 fprintf(out,"REMARK GENERATED BY DO_MULTIPROT\n");
405 sprintf(out_title,"%s t= %g %s",title,output_env_conv_time(oenv,fr.time),output_env_get_time_unit(oenv));
406 /* if reading from pdb, we want to keep the original
407 model numbering else we write the output frame
408 number plus one, because model 0 is not allowed in pdb */
409 if (ftp==efPDB && fr.step > model_nr) {
415 write_pdbfile(out,out_title,&useatoms,fr.x,ePBC,fr.box,' ',model_nr,NULL,TRUE);
418 fr.title = out_title;
419 fr.bTitle = (nframe == 0);
423 write_g96_conf(out,&fr,-1,NULL);
429 } while(read_next_frame(oenv,status,&fr));
431 fres= xvgropen(opt2fn("-cr",NFILE,fnm),"Number of frames in which the residues are aligned to","Residue","Number",oenv);
432 for (i=0;i<ratoms.nres;i++) {
433 fprintf(fres,"%10d %12d\n",countres[i].resnr,countres[i].count);
439 fprintf(stderr,"\n");
441 if (trxout != NULL) {
444 else if (out != NULL) {
447 view_all(oenv,NFILE, fnm);
452 free_t_atoms(&ratoms,TRUE);
454 if (outftp==efPDB || outftp==efGRO || outftp==efG96) {
455 free_t_atoms(&useatoms,TRUE);