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
44 #include "gromacs/utility/smalloc.h"
45 #include "gromacs/commandline/pargs.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"
52 #include "gromacs/fileio/tpxio.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/fileio/confio.h"
63 static void process_multiprot_output(const char *fn
, real
*rmsd
, int *nres
, rvec rotangles
,
64 rvec translation
, bool bCountres
, t_countres
*countres
)
73 mpoutput
=gmx_ffopen (fn
,"r");
77 fgets(line
, 256, mpoutput
);
78 } while (strstr(line
,"Match List") == NULL
);
79 fgets(line
, 256, mpoutput
);
81 string
= strtok (line
,".");
82 while (i
<2 && string
!= NULL
) {
83 string
= strtok (NULL
,". ");
89 while (countres
[j
].resnr
!=res
)
93 fgets(line
, 256, mpoutput
);
94 } while (strstr(line
,"End of Match List") == NULL
);
98 fgets(line
, 256, mpoutput
);
99 } while (strstr(line
,"Trans : ") == NULL
);
101 string
= strtok (line
," :");
102 string
= strtok (NULL
," ");
103 while (i
<3 && string
!= NULL
) {
104 string
= strtok (NULL
," ");
105 rotangles
[i
]=atof(string
);
109 while (i
<3 && string
!= NULL
) {
110 string
= strtok (NULL
," ");
111 translation
[i
]=atof(string
)/10;
115 gmx_warning("Not enough values for rotation and translation vectors in the output of multiprot");
118 rotangles
[YY
]=rotangles
[YY
]*(-1);
121 fgets(line
, 256, mpoutput
);
122 if (strstr(line
,"RMSD : ") != NULL
) {
123 string
= strtok (line
,":");
124 string
= strtok (NULL
,":");
125 (*rmsd
)=atof(string
)/10;
129 fgets(line
,256, mpoutput
);
130 if (strstr(line
,"Match List") != NULL
) {
131 string
= strtok (line
,":");
132 string
= strtok (NULL
,":");
133 (*nres
) = atoi(string
);
136 gmx_ffclose(mpoutput
);
139 int main(int argc
,char *argv
[])
141 const char *desc
[] = {
142 "[TT]do_multiprot[tt] ",
143 "reads a trajectory file and aligns it to a reference structure ",
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.",
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]"
169 static bool bVerbose
;
171 { "-v", FALSE
, etBOOL
, {&bVerbose
},
172 "HIDDENGenerate miles of useless information" }
175 const char *bugs
[] = {
176 "The program is very slow, since multiprot is run externally"
180 t_trxstatus
*trxout
=NULL
;
181 FILE *tapein
,*fo
,*frc
,*tmpf
,*out
=NULL
,*fres
=NULL
;
183 const char *fn
="2_sol.res";
186 t_atoms
*atoms
,ratoms
,useatoms
;
191 int i
,j
,natoms
,nratoms
,nframe
=0,model_nr
=-1;
192 int cur_res
,prev_res
;
194 t_countres
*countres
=NULL
;
200 char pdbfile
[32],refpdb
[256],title
[256],rtitle
[256],filemode
[5];
202 char multiprot
[256],*mptr
;
206 bool bTrjout
,bCountres
;
207 const char *TrjoutFile
=NULL
;
209 static rvec translation
={0,0,0},rotangles
={0,0,0};
210 gmx_rmpbc_t gpbc
=NULL
;
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
}
222 #define NFILE asize(fnm)
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
229 fnRef
=opt2fn("-r",NFILE
,fnm
);
230 bTrjout
= opt2bSet("-ox",NFILE
,fnm
);
231 bCountres
= opt2bSet("-cr",NFILE
,fnm
);
234 TrjoutFile
= opt2fn_null("-ox",NFILE
,fnm
);
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
);
243 get_stx_coordnum(fnRef
,&nratoms
);
244 init_t_atoms(&ratoms
,nratoms
,TRUE
);
246 read_stx_conf(fnRef
,rtitle
,&ratoms
,xr
,NULL
,&ePBC
,rbox
);
249 fprintf(stderr
,"Read %d atoms\n",atoms
->nr
);
250 fprintf(stderr
,"Read %d reference atoms\n",ratoms
.nr
);
253 snew(countres
,ratoms
.nres
);
256 for (i
=0;i
<ratoms
.nr
;i
++) {
258 cur_res
=ratoms
.atom
[i
].resind
;
259 if (cur_res
!= prev_res
) {
260 countres
[j
].resnr
=cur_res
;
266 get_index(atoms
,ftp2fn_null(efNDX
,NFILE
,fnm
),1,&gnx
,&index
,&grpnm
);
269 for(i
=0; (i
<gnx
); i
++) {
270 if (atoms
->atom
[index
[i
]].resind
!= nr0
) {
271 nr0
=atoms
->atom
[index
[i
]].resind
;
275 fprintf(stderr
,"There are %d residues in your selected group\n",nres
);
277 strcpy(pdbfile
,"ddXXXXXX");
279 if ((tmpf
= fopen(pdbfile
,"w")) == NULL
) {
280 sprintf(pdbfile
,"%ctmp%cfilterXXXXXX",DIR_SEPARATOR
,DIR_SEPARATOR
);
282 if ((tmpf
= fopen(pdbfile
,"w")) == NULL
) {
283 gmx_fatal(FARGS
,"Can not open tmp file %s",pdbfile
);
291 strcpy(refpdb
,"ddXXXXXX");
293 strcat(refpdb
,".pdb");
294 write_sto_conf(refpdb
,rtitle
,&ratoms
,xr
,NULL
,ePBC
,rbox
);
297 strcpy(refpdb
,fnRef
);
300 if ((mptr
=getenv("MULTIPROT")) == NULL
) {
301 mptr
="/usr/local/bin/multiprot";
303 if (!gmx_fexist(mptr
)) {
304 gmx_fatal(FARGS
,"MULTIPROT executable (%s) does not exist (use setenv MULTIPROT)",
307 sprintf (multiprot
,"%s %s %s > /dev/null %s",
308 mptr
, refpdb
, pdbfile
, "2> /dev/null");
311 fprintf(stderr
,"multiprot cmd='%s'\n",multiprot
);
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
));
320 outftp
=fn2ftp(TrjoutFile
);
322 fprintf(stderr
,"Will write %s: %s\n",ftp2ext(ftp
),ftp2desc(outftp
));
323 strcpy(filemode
,"w");
330 trxout
= open_trx(TrjoutFile
,filemode
);
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);
346 out
=gmx_ffopen(TrjoutFile
,filemode
);
351 if (natoms
> atoms
->nr
) {
352 gmx_fatal(FARGS
,"\nTrajectory does not match topology!");
355 gmx_fatal(FARGS
,"\nTrajectory does not match selected group!");
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
);
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
);
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
);
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
);
384 write_trxframe(trxout
,&fr
,NULL
);
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
));
393 write_hconf_p(out
,out_title
,&useatoms
,prec2ndec(fr
.prec
),
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
) {
408 write_pdbfile(out
,out_title
,&useatoms
,fr
.x
,ePBC
,fr
.box
,' ',model_nr
,NULL
,TRUE
);
411 fr
.title
= out_title
;
412 fr
.bTitle
= (nframe
== 0);
416 write_g96_conf(out
,&fr
,-1,NULL
);
422 } while(read_next_frame(oenv
,status
,&fr
));
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
);
432 fprintf(stderr
,"\n");
434 if (trxout
!= NULL
) {
437 else if (out
!= NULL
) {
440 view_all(oenv
,NFILE
, fnm
);
445 free_t_atoms(&ratoms
,TRUE
);
447 if (outftp
==efPDB
|| outftp
==efGRO
|| outftp
==efG96
) {
448 free_t_atoms(&useatoms
,TRUE
);