3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
9 * VERSION 3.3.99_development_20071104
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-2006, 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
33 * Groningen Machine for Chemical Simulation
42 #include "gromacs/fileio/pdbio.h"
43 #include "gromacs/fileio/confio.h"
44 #include "gromacs/utility/smalloc.h"
45 #include "gromacs/topology/symtab.h"
48 #include "gromacs/commandline/pargs.h"
49 #include "gromacs/math/vec.h"
52 #include "gromacs/math/units.h"
54 void copy_atom(t_symtab
*tab
,t_atoms
*a1
,int i1
,t_atoms
*a2
,int i2
,
55 rvec xin
[],rvec xout
[],rvec vin
[],rvec vout
[])
57 a2
->atom
[i2
] = a1
->atom
[i1
];
58 a2
->atomname
[i2
] = put_symtab(tab
,*a1
->atomname
[i1
]);
59 a2
->resname
[a2
->atom
[i2
].resnr
] =
60 put_symtab(tab
,*a1
->resname
[a1
->atom
[i1
].resnr
]);
61 copy_rvec(xin
[i1
],xout
[i2
]);
62 copy_rvec(vin
[i1
],vout
[i2
]);
65 static void rotate_x(int natom
,rvec xin
[],real angle
,rvec xout
[],
66 gmx_bool bZ
,gmx_bool bUpsideDown
,real dz
)
74 mat
[XX
][XX
] = cos(angle
);
75 mat
[XX
][YY
] = sin(angle
);
76 mat
[YY
][XX
] = -sin(angle
);
77 mat
[YY
][YY
] = cos(angle
);
82 mat
[YY
][YY
] = cos(angle
);
83 mat
[YY
][ZZ
] = sin(angle
);
84 mat
[ZZ
][YY
] = -sin(angle
);
85 mat
[ZZ
][ZZ
] = cos(angle
);
88 for(i
=0; (i
<natom
); i
++) {
89 mvmul(mat
,xin
[i
],xout
[i
]);
96 static void prep_x(int natom
,rvec x
[],real rDist
,real rAngleZ
,real rAngleX
)
102 /* Center on Z-axis */
104 for(i
=0; (i
<natom
); i
++) {
112 for(i
=0; (i
<natom
); i
++) {
119 rotate_x(natom
,x
,rAngleZ
,xx
,TRUE
,FALSE
,0);
120 for(i
=0; (i
<natom
); i
++)
121 copy_rvec(xx
[i
],x
[i
]);
126 rotate_x(natom
,x
,rAngleX
,xx
,FALSE
,FALSE
,0);
127 for(i
=0; (i
<natom
); i
++)
128 copy_rvec(xx
[i
],x
[i
]);
132 for(i
=0; (i
<natom
); i
++)
137 int main(int argc
, char *argv
[])
140 static char *desc
[] = {
141 "[TT]hexamer[tt] takes a single input coordinate file and makes five symmetry",
144 #define NPA asize(pa)
146 { efSTX
, "-f", NULL
, ffREAD
},
147 { efPDB
, "-o", NULL
, ffWRITE
}
149 #define NFILE asize(fnm)
150 gmx_bool bCenter
= FALSE
;
151 gmx_bool bTrimer
= FALSE
;
152 gmx_bool bAlternate
= FALSE
;
153 real rDist
= 0,rAngleZ
= 0,rAngleX
= 0, alterz
= 0;
155 { "-center", FALSE
, etBOOL
, {&bCenter
},
156 "Center molecule on Z-axis first" },
157 { "-trimer", FALSE
, etBOOL
, {&bTrimer
},
158 "Make trimer rather than hexamer" },
159 { "-alternate",FALSE
, etBOOL
, {&bAlternate
},
160 "Turn every other molecule upside down" },
161 { "-alterz", FALSE
, etREAL
, {&alterz
},
162 "Add this amount to Z-coordinate in every other molecule" },
163 { "-radius", FALSE
, etREAL
, {&rDist
},
164 "Distance of protein axis from Z-axis (implies [TT]-center[tt])" },
165 { "-anglez", FALSE
, etREAL
, {&rAngleZ
},
166 "Initial angle of rotation around Z-axis of protein" },
167 { "-anglex", FALSE
, etREAL
, {&rAngleX
},
168 "Initial angle of rotation around X-axis of protein" }
170 #define NPA asize(pa)
172 int i
,iout
,now
,natom
;
173 rvec
*xin
,*vin
,*xout
;
176 char *infile
,*outfile
,title
[256],buf
[32];
178 CopyRight(stderr
,argv
[0]);
179 parse_common_args(&argc
,argv
,0,NFILE
,fnm
,NPA
,pa
,
180 asize(desc
),desc
,0,NULL
);
181 bCenter
= bCenter
|| (rDist
> 0) || bAlternate
;
183 infile
= ftp2fn(efSTX
,NFILE
,fnm
);
184 outfile
= ftp2fn(efPDB
,NFILE
,fnm
);
186 get_stx_coordnum(infile
,&natom
);
187 init_t_atoms(&atoms
,natom
,TRUE
);
191 read_stx_conf(infile
,title
,&atoms
,xin
,vin
,box
);
192 printf("Read %d atoms\n",atoms
.nr
);
195 prep_x(atoms
.nr
,xin
,rDist
,rAngleZ
,rAngleX
);
197 fp
= gmx_ffopen(outfile
,"w");
198 for(i
=0; (i
<(bTrimer
? 3 : 6)); i
++) {
199 rotate_x(atoms
.nr
,xin
,i
*(bTrimer
? 120.0 : 60.0),xout
,TRUE
,
200 bAlternate
&& ((i
% 2) != 0),alterz
*(((i
% 2) == 0) ? 0 : 1));
201 sprintf(buf
,"Rotated %d degrees",i
*(bTrimer
? 120 : 60));
202 write_pdbfile(fp
,buf
,&atoms
,xout
,box
,'A'+i
,1+i
);