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
61 void copy_atom(t_symtab
*tab
,t_atoms
*a1
,int i1
,t_atoms
*a2
,int i2
,
62 rvec xin
[],rvec xout
[],rvec vin
[],rvec vout
[])
64 a2
->atom
[i2
] = a1
->atom
[i1
];
65 a2
->atomname
[i2
] = put_symtab(tab
,*a1
->atomname
[i1
]);
66 a2
->resname
[a2
->atom
[i2
].resnr
] =
67 put_symtab(tab
,*a1
->resname
[a1
->atom
[i1
].resnr
]);
68 copy_rvec(xin
[i1
],xout
[i2
]);
69 copy_rvec(vin
[i1
],vout
[i2
]);
72 static void rotate_x(int natom
,rvec xin
[],real angle
,rvec xout
[],
73 bool bZ
,bool bUpsideDown
,real dz
)
81 mat
[XX
][XX
] = cos(angle
);
82 mat
[XX
][YY
] = sin(angle
);
83 mat
[YY
][XX
] = -sin(angle
);
84 mat
[YY
][YY
] = cos(angle
);
89 mat
[YY
][YY
] = cos(angle
);
90 mat
[YY
][ZZ
] = sin(angle
);
91 mat
[ZZ
][YY
] = -sin(angle
);
92 mat
[ZZ
][ZZ
] = cos(angle
);
95 for(i
=0; (i
<natom
); i
++) {
96 mvmul(mat
,xin
[i
],xout
[i
]);
103 static void prep_x(int natom
,rvec x
[],real rDist
,real rAngleZ
,real rAngleX
)
109 /* Center on Z-axis */
111 for(i
=0; (i
<natom
); i
++) {
119 for(i
=0; (i
<natom
); i
++) {
126 rotate_x(natom
,x
,rAngleZ
,xx
,TRUE
,FALSE
,0);
127 for(i
=0; (i
<natom
); i
++)
128 copy_rvec(xx
[i
],x
[i
]);
133 rotate_x(natom
,x
,rAngleX
,xx
,FALSE
,FALSE
,0);
134 for(i
=0; (i
<natom
); i
++)
135 copy_rvec(xx
[i
],x
[i
]);
139 for(i
=0; (i
<natom
); i
++)
144 int main(int argc
, char *argv
[])
147 static char *desc
[] = {
148 "hexamer takes a single input coordinate file and makes five symmetry",
151 #define NPA asize(pa)
153 { efSTX
, "-f", NULL
, ffREAD
},
154 { efPDB
, "-o", NULL
, ffWRITE
}
156 #define NFILE asize(fnm)
157 bool bCenter
= FALSE
;
158 bool bTrimer
= FALSE
;
159 bool bAlternate
= FALSE
;
160 real rDist
= 0,rAngleZ
= 0,rAngleX
= 0, alterz
= 0;
162 { "-center", FALSE
, etBOOL
, {&bCenter
},
163 "Center molecule on Z-axis first" },
164 { "-trimer", FALSE
, etBOOL
, {&bTrimer
},
165 "Make trimer rather than hexamer" },
166 { "-alternate",FALSE
, etBOOL
, {&bAlternate
},
167 "Turn every other molecule upside down" },
168 { "-alterz", FALSE
, etREAL
, {&alterz
},
169 "Add this amount to Z-coordinate in every other molecule" },
170 { "-radius", FALSE
, etREAL
, {&rDist
},
171 "Distance of protein axis from Z-axis (implies -center)" },
172 { "-anglez", FALSE
, etREAL
, {&rAngleZ
},
173 "Initial angle of rotation around Z-axis of protein" },
174 { "-anglex", FALSE
, etREAL
, {&rAngleX
},
175 "Initial angle of rotation around X-axis of protein" }
177 #define NPA asize(pa)
179 int i
,iout
,now
,natom
;
180 rvec
*xin
,*vin
,*xout
;
183 char *infile
,*outfile
,title
[256],buf
[32];
185 CopyRight(stderr
,argv
[0]);
186 parse_common_args(&argc
,argv
,0,NFILE
,fnm
,NPA
,pa
,
187 asize(desc
),desc
,0,NULL
);
188 bCenter
= bCenter
|| (rDist
> 0) || bAlternate
;
190 infile
= ftp2fn(efSTX
,NFILE
,fnm
);
191 outfile
= ftp2fn(efPDB
,NFILE
,fnm
);
193 get_stx_coordnum(infile
,&natom
);
194 init_t_atoms(&atoms
,natom
,TRUE
);
198 read_stx_conf(infile
,title
,&atoms
,xin
,vin
,box
);
199 printf("Read %d atoms\n",atoms
.nr
);
202 prep_x(atoms
.nr
,xin
,rDist
,rAngleZ
,rAngleX
);
204 fp
= ffopen(outfile
,"w");
205 for(i
=0; (i
<(bTrimer
? 3 : 6)); i
++) {
206 rotate_x(atoms
.nr
,xin
,i
*(bTrimer
? 120.0 : 60.0),xout
,TRUE
,
207 bAlternate
&& ((i
% 2) != 0),alterz
*(((i
% 2) == 0) ? 0 : 1));
208 sprintf(buf
,"Rotated %d degrees",i
*(bTrimer
? 120 : 60));
209 write_pdbfile(fp
,buf
,&atoms
,xout
,box
,'A'+i
,1+i
);