changed reading hint
[gromacs/adressmacs.git] / src / kernel / add_par.c
blob758a59e04bbffef2550129d7179fcbcde921ecef
1 /*
2 * $Id$
3 *
4 * This source code is part of
5 *
6 * G R O M A C S
7 *
8 * GROningen MAchine for Chemical Simulations
9 *
10 * VERSION 2.0
12 * Copyright (c) 1991-1999
13 * BIOSON Research Institute, Dept. of Biophysical Chemistry
14 * University of Groningen, The Netherlands
16 * Please refer to:
17 * GROMACS: A message-passing parallel molecular dynamics implementation
18 * H.J.C. Berendsen, D. van der Spoel and R. van Drunen
19 * Comp. Phys. Comm. 91, 43-56 (1995)
21 * Also check out our WWW page:
22 * http://md.chem.rug.nl/~gmx
23 * or e-mail to:
24 * gromacs@chem.rug.nl
26 * And Hey:
27 * GRowing Old MAkes el Chrono Sweat
29 static char *SRCID_add_par_c = "$Id$";
31 #include <string.h>
32 #include "typedefs.h"
33 #include "smalloc.h"
34 #include "grompp.h"
35 #include "toputil.h"
36 #include "hackblock.h"
38 static void clear_atom_list(int i0, atom_id a[])
40 int i;
42 for (i=i0; i < MAXATOMLIST; i++)
43 a[i]=-1;
46 static void clear_force_param(int i0, real c[])
48 int i;
50 for (i=i0; i < MAXFORCEPARAM; i++)
51 c[i]=NOTSET;
54 void add_param(t_params *ps,int ai,int aj, real *c, char *s)
56 int i;
58 if ((ai < 0) || (aj < 0))
59 fatal_error(0,"Trying to add impossible atoms: ai=%d, aj=%d",ai,aj);
60 pr_alloc(1,ps);
61 ps->param[ps->nr].AI=ai;
62 ps->param[ps->nr].AJ=aj;
63 clear_atom_list(2, ps->param[ps->nr].a);
64 if (c==NULL)
65 clear_force_param(0, ps->param[ps->nr].c);
66 else
67 for(i=0; (i < MAXFORCEPARAM); i++)
68 ps->param[ps->nr].c[i]=c[i];
69 if (s)
70 ps->param[ps->nr].s = strdup(s);
71 else
72 ps->param[ps->nr].s = strdup("");
73 ps->nr++;
76 void add_imp_param(t_params *ps,int ai,int aj,int ak,int al,real c0, real c1,
77 char *s)
79 pr_alloc(1,ps);
80 ps->param[ps->nr].AI=ai;
81 ps->param[ps->nr].AJ=aj;
82 ps->param[ps->nr].AK=ak;
83 ps->param[ps->nr].AL=al;
84 clear_atom_list (4, ps->param[ps->nr].a);
85 ps->param[ps->nr].C0=c0;
86 ps->param[ps->nr].C1=c1;
87 clear_force_param(2, ps->param[ps->nr].c);
88 if (s)
89 ps->param[ps->nr].s = strdup(s);
90 else
91 ps->param[ps->nr].s = strdup("");
92 ps->nr++;
95 void add_dih_param(t_params *ps,int ai,int aj,int ak,int al,real c0, real c1,
96 real c2,char *s)
98 pr_alloc(1,ps);
99 ps->param[ps->nr].AI=ai;
100 ps->param[ps->nr].AJ=aj;
101 ps->param[ps->nr].AK=ak;
102 ps->param[ps->nr].AL=al;
103 clear_atom_list (4, ps->param[ps->nr].a);
104 ps->param[ps->nr].C0=c0;
105 ps->param[ps->nr].C1=c1;
106 ps->param[ps->nr].C2=c2;
107 clear_force_param(3, ps->param[ps->nr].c);
108 if (s)
109 ps->param[ps->nr].s = strdup(s);
110 else
111 ps->param[ps->nr].s = strdup("");
112 ps->nr++;
115 void add_dum2_atoms(t_params *ps,int ai,int aj,int ak)
117 pr_alloc(1,ps);
118 ps->param[ps->nr].AI=ai;
119 ps->param[ps->nr].AJ=aj;
120 ps->param[ps->nr].AK=ak;
121 clear_atom_list (3, ps->param[ps->nr].a);
122 clear_force_param(0, ps->param[ps->nr].c);
123 ps->param[ps->nr].s = strdup("");
124 ps->nr++;
127 void add_dum3_param(t_params *ps,int ai,int aj,int ak,int al,
128 real c0, real c1)
130 pr_alloc(1,ps);
131 ps->param[ps->nr].AI=ai;
132 ps->param[ps->nr].AJ=aj;
133 ps->param[ps->nr].AK=ak;
134 ps->param[ps->nr].AL=al;
135 clear_atom_list (4, ps->param[ps->nr].a);
136 ps->param[ps->nr].C0=c0;
137 ps->param[ps->nr].C1=c1;
138 clear_force_param(2, ps->param[ps->nr].c);
139 ps->param[ps->nr].s = strdup("");
140 ps->nr++;
143 void add_dum3_atoms(t_params *ps,int ai,int aj,int ak,int al, bool bSwapParity)
145 pr_alloc(1,ps);
146 ps->param[ps->nr].AI=ai;
147 ps->param[ps->nr].AJ=aj;
148 ps->param[ps->nr].AK=ak;
149 ps->param[ps->nr].AL=al;
150 clear_atom_list (4, ps->param[ps->nr].a);
151 clear_force_param(0, ps->param[ps->nr].c);
152 if (bSwapParity)
153 ps->param[ps->nr].C1=-1;
154 ps->param[ps->nr].s = strdup("");
155 ps->nr++;
158 void add_dum4_atoms(t_params *ps,int ai,int aj,int ak,int al,int am)
160 pr_alloc(1,ps);
161 ps->param[ps->nr].AI=ai;
162 ps->param[ps->nr].AJ=aj;
163 ps->param[ps->nr].AK=ak;
164 ps->param[ps->nr].AL=al;
165 ps->param[ps->nr].AM=am;
166 clear_atom_list (5, ps->param[ps->nr].a);
167 clear_force_param(0, ps->param[ps->nr].c);
168 ps->param[ps->nr].s = strdup("");
169 ps->nr++;
172 int search_jtype(t_restp *rtp,char *name,bool bNterm)
174 int j,k,kmax,jmax;
175 char *rtpname,searchname[12];
177 strcpy(searchname,name);
179 /* Do a best match comparison */
180 /* for protein N-terminus, rename H1, H2 and H3 to H */
181 if ( bNterm && (strlen(searchname)==2) && (searchname[0] == 'H') &&
182 ( (searchname[1] == '1') || (searchname[1] == '2') ||
183 (searchname[1] == '3') ) )
184 searchname[1]='\0';
185 kmax=0;
186 jmax=-1;
187 for(j=0; (j<rtp->natom); j++) {
188 rtpname=*(rtp->atomname[j]);
189 if (strcasecmp(searchname,rtpname) == 0) {
190 jmax=j;
191 kmax=strlen(searchname);
192 break;
194 for(k=0; k < min(strlen(searchname), strlen(rtpname)); k++)
195 if (searchname[k] != rtpname[k])
196 break;
197 if (k > kmax) {
198 kmax=k;
199 jmax=j;
202 if (jmax == -1)
203 fatal_error(0,"Atom %s not found in database in residue %s",
204 searchname,rtp->resname);
205 if (kmax != strlen(searchname))
206 fatal_error(0,"Atom %s not found in database in residue %s, "
207 "it looks a bit like %s",
208 searchname,rtp->resname,*(rtp->atomname[jmax]));
209 return jmax;