Merge branch release-5-1
[alexxy/gromacs.git] / src / contrib / do_multiprot.c
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
3  *                This source code is NOT REALLY part of
4  * 
5  *                 G   R   O   M   A   C   S
6  * 
7  *          GROningen MAchine for Chemical Simulations
8  * 
9  *                        VERSION 4.2.5
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.
14
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.
19  * 
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.
26  * 
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.
29  * 
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * Author: do_multiprot was written by Ran Friedman <r.friedman@bioc.uzh.ch>
33  * 
34  * And Hey:
35  * Green Red Orange Magenta Azure Cyan Skyblue
36  */
37
38 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 #include "typedefs.h"
43 #include "macros.h"
44 #include "gromacs/utility/smalloc.h"
45 #include "gromacs/commandline/pargs.h"
46 #include "copyrite.h"
47 #include "gromacs/fileio/pdbio.h"
48 #include "gromacs/utility/fatalerror.h"
49 #include "gromacs/fileio/xvgr.h"
50 #include "gromacs/topology/index.h"
51 #include "gstat.h"
52 #include "gromacs/fileio/tpxio.h"
53 #include "viewit.h"
54 #include "gbutil.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/fileio/confio.h"
57
58 typedef struct {
59     int resnr; 
60     int count; 
61 } t_countres;
62
63 static void process_multiprot_output(const char *fn, real *rmsd, int *nres, rvec rotangles, 
64                                      rvec translation, bool bCountres, t_countres *countres)
65 {
66     FILE       *mpoutput;
67     char       line[256];
68     char       *string;
69     int        i=0,j=0,res;
70     
71     (*rmsd)=-1;
72     (*nres)=0;
73     mpoutput=gmx_ffopen (fn,"r");
74     
75     if (bCountres) {
76         do {
77             fgets(line, 256, mpoutput);
78         } while (strstr(line,"Match List") == NULL);
79         fgets(line, 256, mpoutput);
80         do {
81             string = strtok (line,".");
82             while (i<2 && string != NULL) {
83                 string = strtok (NULL,". ");
84                 i++;
85             }
86             i=0;
87             res=atoi(string);
88             if (res > 0) {
89                 while (countres[j].resnr!=res)
90                     j++;
91                 countres[j].count++;
92             }
93             fgets(line, 256, mpoutput);
94         } while (strstr(line,"End of Match List") == NULL);
95         rewind(mpoutput);
96     }
97     do {
98         fgets(line, 256, mpoutput);
99     } while (strstr(line,"Trans : ") == NULL);
100     
101     string = strtok (line," :");
102     string = strtok (NULL," ");
103     while (i<3 && string != NULL) {
104         string = strtok (NULL," ");
105         rotangles[i]=atof(string);
106         i++;
107     }
108     i=0;
109     while (i<3 && string != NULL) {
110         string = strtok (NULL," ");
111         translation[i]=atof(string)/10;
112         i++;
113     }
114     if (i!=3) {
115         gmx_warning("Not enough values for rotation and translation vectors in the output of multiprot");
116     }
117     
118     rotangles[YY]=rotangles[YY]*(-1); 
119     
120     while ((*rmsd) <0) {
121         fgets(line, 256, mpoutput);
122         if (strstr(line,"RMSD : ") != NULL) {
123             string = strtok (line,":");
124             string = strtok (NULL,":");
125             (*rmsd)=atof(string)/10;
126         }
127     }
128     while (!(*nres)) {
129         fgets(line,256, mpoutput);
130         if (strstr(line,"Match List") != NULL) {
131             string = strtok (line,":");
132             string = strtok (NULL,":");
133             (*nres) = atoi(string);
134         }
135     }
136     gmx_ffclose(mpoutput);
137 }
138
139 int main(int argc,char *argv[])
140 {
141     const char *desc[] = {
142         "[TT]do_multiprot[tt] ", 
143         "reads a trajectory file and aligns it to a reference structure  ",
144         "each time frame",
145         "calling the multiprot program. This allows you to use a reference",
146         "structure whose sequence is different than that of the protein in the ",
147         "trajectory, since the alignment is based on the geometry, not sequence.",
148         "The output of [TT]do_multiprot[tt] includes the rmsd and the number of residues",
149         "on which it was calculated.",
150         "[PAR]",
151         "An aligned trajectory file is generated with the [TT]-ox[tt] option.[PAR]",
152         "With the [TT]-cr[tt] option, the number of hits in the alignment is given",
153         "per residue. This number can be between 0 and the number of frames, and",
154         "indicates the structural conservation of this residue.[PAR]",
155         "If you do not have the [TT]multiprot[tt] program, get it. [TT]do_multiprot[tt] assumes", 
156         "that the [TT]multiprot[tt] executable is [TT]/usr/local/bin/multiprot[tt]. If this is ",
157         "not the case, then you should set an environment variable [BB]MULTIPROT[bb]", 
158         "pointing to the [TT]multiprot[tt] executable, e.g.: [PAR]",
159         "[TT]setenv MULTIPROT /usr/MultiProtInstall/multiprot.Linux[tt][PAR]",
160         "Note that at the current implementation only binary alignment (your",
161         "molecule to a reference) is supported. In addition, note that the ",
162         "by default [TT]multiprot[tt] aligns the two proteins on their C-alpha carbons.",
163         "and that this depends on the [TT]multiprot[tt] parameters which are not dealt ",
164         "with here. Thus, the C-alpha carbons is expected to give the same "
165         "results as choosing the whole protein and will be slightly faster.[PAR]",
166         "For information about [TT]multiprot[tt], see:",
167         "http://bioinfo3d.cs.tau.ac.il/MultiProt/.[PAR]"
168     };
169     static bool bVerbose;
170     t_pargs pa[] = {
171         { "-v",  FALSE, etBOOL, {&bVerbose},
172           "HIDDENGenerate miles of useless information" }
173     };
174   
175     const char *bugs[] = { 
176         "The program is very slow, since multiprot is run externally"
177     };
178   
179     t_trxstatus *status;
180     t_trxstatus *trxout=NULL;
181     FILE        *tapein,*fo,*frc,*tmpf,*out=NULL,*fres=NULL;
182     const char  *fnRef;
183     const char  *fn="2_sol.res";
184     t_topology  top;
185     int         ePBC;
186     t_atoms     *atoms,ratoms,useatoms;
187     t_trxframe  fr;
188     t_pdbinfo   p;
189     int         nres,nres2,nr0;
190     real        t;
191     int         i,j,natoms,nratoms,nframe=0,model_nr=-1;
192     int         cur_res,prev_res;
193     int         nout;
194     t_countres  *countres=NULL;
195     matrix      box,rbox;
196     int         gnx;
197     char        *grpnm,*ss_str; 
198     atom_id     *index;
199     rvec        *xp,*x,*xr;
200     char        pdbfile[32],refpdb[256],title[256],rtitle[256],filemode[5];
201     char        out_title[256];
202     char        multiprot[256],*mptr;
203     int         ftp;
204     int         outftp=-1;
205     real        rmsd;
206     bool        bTrjout,bCountres;
207     const char  *TrjoutFile=NULL;
208     output_env_t oenv;
209     static rvec translation={0,0,0},rotangles={0,0,0};
210     gmx_rmpbc_t gpbc=NULL;
211     
212     t_filenm   fnm[] = {
213         { efTRX, "-f",   NULL,      ffREAD },
214         { efTPS, NULL,   NULL,      ffREAD },
215         { efNDX, NULL,   NULL,      ffOPTRD },
216         { efSTX, "-r",   NULL     , ffREAD },
217         { efXVG, "-o",  "rmss",     ffWRITE },
218         { efXVG, "-rc", "rescount", ffWRITE},
219         { efXVG, "-cr", "countres", ffOPTWR},
220         { efTRX, "-ox", "aligned",  ffOPTWR }
221     };
222 #define NFILE asize(fnm)
223     
224     CopyRight(stderr,argv[0]);
225     parse_common_args(&argc,argv,PCA_CAN_TIME | PCA_CAN_VIEW | PCA_TIME_UNIT,
226                       NFILE,fnm, asize(pa),pa, asize(desc),desc,
227                       asize(bugs),bugs,&oenv
228         );
229     fnRef=opt2fn("-r",NFILE,fnm);
230     bTrjout = opt2bSet("-ox",NFILE,fnm);
231     bCountres=  opt2bSet("-cr",NFILE,fnm);
232     
233     if (bTrjout) {
234         TrjoutFile = opt2fn_null("-ox",NFILE,fnm);
235     }
236     
237     read_tps_conf(ftp2fn(efTPS,NFILE,fnm),title,&top,&ePBC,&xp,NULL,box,FALSE);
238     gpbc = gmx_rmpbc_init(&top.idef,ePBC,top.atoms.nr,box);
239     atoms=&(top.atoms);
240
241     ftp=fn2ftp(fnRef);
242  
243     get_stx_coordnum(fnRef,&nratoms);
244     init_t_atoms(&ratoms,nratoms,TRUE);  
245     snew(xr,nratoms);
246     read_stx_conf(fnRef,rtitle,&ratoms,xr,NULL,&ePBC,rbox);
247     
248     if (bVerbose) {
249         fprintf(stderr,"Read %d atoms\n",atoms->nr); 
250         fprintf(stderr,"Read %d reference atoms\n",ratoms.nr); 
251     }
252     if (bCountres) {
253         snew(countres,ratoms.nres);
254         j=0;
255         cur_res=0;
256         for (i=0;i<ratoms.nr;i++) {
257             prev_res=cur_res;
258             cur_res=ratoms.atom[i].resind;
259             if (cur_res != prev_res) {
260                 countres[j].resnr=cur_res;
261                 countres[j].count=0;
262                 j++;
263             }
264         }
265     }
266     get_index(atoms,ftp2fn_null(efNDX,NFILE,fnm),1,&gnx,&index,&grpnm);
267     nres=0;
268     nr0=-1;
269     for(i=0; (i<gnx); i++) {
270         if (atoms->atom[index[i]].resind != nr0) {
271             nr0=atoms->atom[index[i]].resind;
272             nres++;
273         }
274     }
275     fprintf(stderr,"There are %d residues in your selected group\n",nres);
276     
277     strcpy(pdbfile,"ddXXXXXX");
278     gmx_tmpnam(pdbfile);
279     if ((tmpf = fopen(pdbfile,"w")) == NULL) {
280         sprintf(pdbfile,"%ctmp%cfilterXXXXXX",DIR_SEPARATOR,DIR_SEPARATOR);
281         gmx_tmpnam(pdbfile);
282         if ((tmpf = fopen(pdbfile,"w")) == NULL) {
283             gmx_fatal(FARGS,"Can not open tmp file %s",pdbfile);
284         }
285     }
286     else {
287         gmx_ffclose(tmpf);
288     }
289
290     if (ftp != efPDB) {
291         strcpy(refpdb,"ddXXXXXX");
292         gmx_tmpnam(refpdb);
293         strcat(refpdb,".pdb");
294         write_sto_conf(refpdb,rtitle,&ratoms,xr,NULL,ePBC,rbox);
295     }
296     else {
297         strcpy(refpdb,fnRef);
298     }
299
300     if ((mptr=getenv("MULTIPROT")) == NULL) {
301         mptr="/usr/local/bin/multiprot";
302     }
303     if (!gmx_fexist(mptr)) {
304         gmx_fatal(FARGS,"MULTIPROT executable (%s) does not exist (use setenv MULTIPROT)",
305                   mptr);
306     }
307     sprintf (multiprot,"%s %s %s > /dev/null %s",
308              mptr, refpdb, pdbfile, "2> /dev/null");
309     
310     if (bVerbose)
311         fprintf(stderr,"multiprot cmd='%s'\n",multiprot);
312     
313     if (!read_first_frame(oenv,&status,ftp2fn(efTRX,NFILE,fnm),&fr,TRX_READ_X)) 
314         gmx_fatal(FARGS,"Could not read a frame from %s",ftp2fn(efTRX,NFILE,fnm));
315     natoms = fr.natoms;
316
317     if (bTrjout) {
318         nout=natoms;
319         /* open file now */
320         outftp=fn2ftp(TrjoutFile);
321         if (bVerbose)
322             fprintf(stderr,"Will write %s: %s\n",ftp2ext(ftp),ftp2desc(outftp));
323         strcpy(filemode,"w");
324         switch (outftp) {
325             case efXTC:
326             case efG87:
327             case efTRR:
328             case efTRJ:
329                 out=NULL;
330                 trxout = open_trx(TrjoutFile,filemode);
331                 break;
332             case efGRO:
333             case efG96:
334             case efPDB:
335                 /* Make atoms struct for output in GRO or PDB files */
336                 /* get memory for stuff to go in pdb file */
337                 init_t_atoms(&useatoms,nout,FALSE);
338                 sfree(useatoms.resinfo);
339                 useatoms.resinfo=atoms->resinfo;
340                 for(i=0;(i<nout);i++) {
341                     useatoms.atomname[i]=atoms->atomname[i];
342                     useatoms.atom[i]=atoms->atom[i];
343                     useatoms.nres=max(useatoms.nres,useatoms.atom[i].resind+1);
344                 }
345                 useatoms.nr=nout;
346                 out=gmx_ffopen(TrjoutFile,filemode);
347                 break;
348         }
349     }
350     
351     if (natoms > atoms->nr) {
352         gmx_fatal(FARGS,"\nTrajectory does not match topology!");
353     }
354     if (gnx > natoms) {
355         gmx_fatal(FARGS,"\nTrajectory does not match selected group!");
356     }
357
358     fo = xvgropen(opt2fn("-o",NFILE,fnm),"RMSD","Time (ps)","RMSD (nm)",oenv);
359     frc = xvgropen(opt2fn("-rc",NFILE,fnm),"Number of Residues in the alignment","Time (ps)","Residues",oenv);
360     
361     do {
362         t = output_env_conv_time(oenv,fr.time);
363         gmx_rmpbc(gpbc,natoms,fr.box,fr.x);
364         tapein=gmx_ffopen(pdbfile,"w");
365         write_pdbfile_indexed(tapein,NULL,atoms,fr.x,ePBC,fr.box,' ',-1,gnx,index,NULL,TRUE); 
366         gmx_ffclose(tapein);
367         system(multiprot);
368         remove(pdbfile);
369         process_multiprot_output(fn, &rmsd, &nres2,rotangles,translation,bCountres,countres);
370         fprintf(fo,"%12.7f",t);
371         fprintf(fo," %12.7f\n",rmsd);
372         fprintf(frc,"%12.7f",t);
373         fprintf(frc,"%12d\n",nres2);
374         if (bTrjout) {
375             rotate_conf(natoms,fr.x,NULL,rotangles[XX],rotangles[YY],rotangles[ZZ]);
376             for(i=0; i<natoms; i++) {
377                 rvec_inc(fr.x[i],translation);
378             }
379             switch(outftp) {
380                 case efTRJ:
381                 case efTRR:
382                 case efG87:
383                 case efXTC:
384                     write_trxframe(trxout,&fr,NULL);
385                     break;
386                 case efGRO:
387                 case efG96:
388                 case efPDB:
389                     sprintf(out_title,"Generated by do_multiprot : %s t= %g %s",
390                             title,output_env_conv_time(oenv,fr.time),output_env_get_time_unit(oenv));
391                     switch(outftp) {
392                         case efGRO: 
393                             write_hconf_p(out,out_title,&useatoms,prec2ndec(fr.prec),
394                                           fr.x,NULL,fr.box);
395                             break;
396                         case efPDB:
397                             fprintf(out,"REMARK    GENERATED BY DO_MULTIPROT\n");
398                             sprintf(out_title,"%s t= %g %s",title,output_env_conv_time(oenv,fr.time),output_env_get_time_unit(oenv));
399                             /* if reading from pdb, we want to keep the original 
400                                model numbering else we write the output frame
401                                number plus one, because model 0 is not allowed in pdb */
402                             if (ftp==efPDB && fr.step > model_nr) {
403                                 model_nr = fr.step;
404                             }
405                             else {
406                                 model_nr++;
407                             }
408                             write_pdbfile(out,out_title,&useatoms,fr.x,ePBC,fr.box,' ',model_nr,NULL,TRUE);
409                             break;
410                         case efG96:
411                             fr.title = out_title;
412                             fr.bTitle = (nframe == 0);
413                             fr.bAtoms = FALSE;
414                             fr.bStep = TRUE;
415                             fr.bTime = TRUE;
416                             write_g96_conf(out,&fr,-1,NULL);
417                     }
418                     break;
419             }
420         }
421         nframe++;
422     } while(read_next_frame(oenv,status,&fr));
423     if (bCountres) {
424         fres=  xvgropen(opt2fn("-cr",NFILE,fnm),"Number of frames in which the residues are aligned to","Residue","Number",oenv);
425         for (i=0;i<ratoms.nres;i++) {
426             fprintf(fres,"%10d  %12d\n",countres[i].resnr,countres[i].count);
427         }
428         gmx_ffclose(fres);
429     }
430     gmx_ffclose(fo);
431     gmx_ffclose(frc);
432     fprintf(stderr,"\n");
433     close_trj(status);
434     if (trxout != NULL) {
435         close_trx(trxout);
436     }
437     else if (out != NULL) {
438         gmx_ffclose(out);
439     }
440     view_all(oenv,NFILE, fnm);
441     sfree(xr);
442     if (bCountres) {
443         sfree(countres);
444     }
445     free_t_atoms(&ratoms,TRUE);
446     if (bTrjout) {
447         if (outftp==efPDB || outftp==efGRO || outftp==efG96) {
448             free_t_atoms(&useatoms,TRUE);
449         }
450     }
451     gmx_thanx(stderr);
452     return 0;
453 }