364bdc1578d9454bab556ccbffe7aa34016ae8e5
[alexxy/gromacs.git] / src / kernel / gctio.c
1 /*
2  * 
3  *                This source code is part of
4  * 
5  *                 G   R   O   M   A   C   S
6  * 
7  *          GROningen MAchine for Chemical Simulations
8  * 
9  *                        VERSION 3.2.0
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  * And Hey:
33  * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include "typedefs.h"
40 #include "xmdrun.h"
41 #include "futil.h"
42 #include "xvgr.h"
43 #include "macros.h"
44 #include "physics.h"
45 #include "network.h"
46 #include "smalloc.h"
47 #include "string2.h"
48 #include "readinp.h"
49 #include "readir.h"
50 #include "filenm.h"
51 #include "names.h"
52 #include "gmxfio.h"
53
54 const char *eoNames[eoNR] = { 
55   "Pres", "Epot", "Vir", "Dist", "Mu", "Force", "Fx", "Fy", "Fz",
56   "Px", "Py", "Pz",
57   "Polarizability", "Dipole", "Memory", "UseEinter", "UseVirial",
58   "CombinationRule"
59 };
60
61 static int Name2eo(char *s)
62 {
63   int i,res;
64   
65   res=-1;
66   
67   for(i=0; (i<eoNR); i++)
68     if (gmx_strcasecmp(s,eoNames[i]) == 0) {
69       res=i;
70       fprintf(stderr,"Coupling to observable %d (%s)\n",res,eoNames[res]);
71       break;
72     }
73   
74   return res;
75 }
76
77 #define  block_bc(cr,   d) gmx_bcast(     sizeof(d),     &(d),(cr))
78 #define nblock_bc(cr,nr,d) gmx_bcast((nr)*sizeof((d)[0]), (d),(cr))
79 #define   snew_bc(cr,d,nr) { if (!MASTER(cr)) snew((d),(nr)); }
80
81 static void low_comm_tcr(t_commrec *cr,t_coupl_rec *tcr)
82 {
83   nblock_bc(cr,eoObsNR,tcr->ref_value);
84   
85   block_bc(cr,tcr->nLJ);
86   snew_bc(cr,tcr->tcLJ,tcr->nLJ);
87   nblock_bc(cr,tcr->nLJ,tcr->tcLJ);
88   
89   block_bc(cr,tcr->nBU);
90   snew_bc(cr,tcr->tcBU,tcr->nBU);
91   nblock_bc(cr,tcr->nBU,tcr->tcBU);
92   
93   block_bc(cr,tcr->nQ);
94   snew_bc(cr,tcr->tcQ,tcr->nQ);
95   nblock_bc(cr,tcr->nQ,tcr->tcQ);
96   
97   block_bc(cr,tcr->nmemory);
98   block_bc(cr,tcr->bInter);
99   block_bc(cr,tcr->bVirial);
100   block_bc(cr,tcr->combrule);
101 }
102
103 void comm_tcr(FILE *log,t_commrec *cr,t_coupl_rec **tcr)
104 {
105   if (!MASTER(cr))
106     snew(*tcr,1);
107   
108   low_comm_tcr(cr,*tcr);
109
110
111 static void clear_lj(t_coupl_LJ *tc)
112 {
113   tc->at_i   = 0;
114   tc->at_j   = 0;
115   tc->eObs   = -1;
116   tc->bPrint = TRUE;
117   tc->c6     = 0.0;
118   tc->c12    = 0.0;
119   tc->xi_6   = 0.0;
120   tc->xi_12  = 0.0;
121 }
122
123 static void clear_bu(t_coupl_BU *tc)
124 {
125   tc->at_i   = 0;
126   tc->at_j   = 0;
127   tc->eObs   = -1;
128   tc->bPrint = TRUE;
129   tc->a      = 0.0;
130   tc->b      = 0.0;
131   tc->c      = 0.0;
132   tc->xi_a   = 0.0;
133   tc->xi_b   = 0.0;
134   tc->xi_c   = 0.0;
135 }
136
137 static void clear_q(t_coupl_Q *tc)
138 {
139   tc->at_i   = 0;
140   tc->eObs   = -1;
141   tc->bPrint = TRUE;
142   tc->Q      = 0.0;
143   tc->xi_Q   = 0.0;
144 }
145
146 void copy_ff(t_coupl_rec *tcr,t_forcerec *fr,t_mdatoms *md,t_idef *idef)
147 {
148   int        i,j,ati,atj,type;
149   t_coupl_LJ *tclj;
150   t_coupl_BU *tcbu;
151   t_coupl_Q  *tcq;
152   
153   /* Save values for printing */
154   for(i=0; (i<tcr->nLJ); i++) {
155     tclj = &(tcr->tcLJ[i]);
156     
157     ati  = tclj->at_i;
158     atj  = tclj->at_j;
159     if (atj == -1)
160       atj = ati;
161     tclj->c6  = C6(fr->nbfp,fr->ntype,ati,atj);
162     tclj->c12 = C12(fr->nbfp,fr->ntype,ati,atj);
163   }
164   
165   for(i=0; (i<tcr->nBU); i++) {
166     tcbu = &(tcr->tcBU[i]);
167     
168     ati  = tcbu->at_i;
169     atj  = tcbu->at_j;
170     if (atj == -1)
171       atj = ati;
172     tcbu->a = BHAMA(fr->nbfp,fr->ntype,ati,atj);
173     tcbu->b = BHAMB(fr->nbfp,fr->ntype,ati,atj);
174     tcbu->c = BHAMC(fr->nbfp,fr->ntype,ati,atj);
175   }
176   
177   for(i=0; (i<tcr->nQ); i++) {
178     tcq = &(tcr->tcQ[i]);
179     for(j=0; (j<md->nr); j++) {
180       if (md->typeA[j] == tcq->at_i) {
181         tcr->tcQ[i].Q = md->chargeA[j];
182         break;
183       }
184     }
185   }
186   for(i=0; (i<tcr->nIP); i++) {
187     /* Let's just copy the whole struct !*/
188     type = tcr->tIP[i].type;
189     tcr->tIP[i].iprint=idef->iparams[type];
190   }
191 }
192
193 void write_gct(const char *fn,t_coupl_rec *tcr,t_idef *idef)
194 {
195   FILE *fp;
196   int  i,ftype;
197   
198   fp=gmx_fio_fopen(fn,"w");
199   nice_header(fp,fn);
200   fprintf(fp,"%-15s = %12g  ; Reference pressure for coupling\n",
201           eoNames[eoPres],tcr->ref_value[eoPres]);
202   fprintf(fp,"%-15s = %12g  ; Reference potential energy\n",
203           eoNames[eoEpot],tcr->ref_value[eoEpot]);
204   fprintf(fp,"%-15s = %12g  ; Reference distance\n",
205           eoNames[eoDist],tcr->ref_value[eoDist]);
206   fprintf(fp,"%-15s = %12g  ; Reference dipole\n",
207           eoNames[eoMu],tcr->ref_value[eoMu]);
208   fprintf(fp,"%-15s = %12g  ; Reference force\n",
209           eoNames[eoForce],tcr->ref_value[eoForce]);
210   fprintf(fp,"%-15s = %12g  ; Reference force in X dir\n",
211           eoNames[eoFx],tcr->ref_value[eoFx]);
212   fprintf(fp,"%-15s = %12g  ; Reference force in Y dir\n",
213           eoNames[eoFy],tcr->ref_value[eoFy]);
214   fprintf(fp,"%-15s = %12g  ; Reference force in Z dir\n",
215           eoNames[eoFz],tcr->ref_value[eoFz]);
216   fprintf(fp,"%-15s = %12g  ; Reference pres in X dir\n",
217           eoNames[eoPx],tcr->ref_value[eoPx]);
218   fprintf(fp,"%-15s = %12g  ; Reference pres in Y dir\n",
219           eoNames[eoPy],tcr->ref_value[eoPy]);
220   fprintf(fp,"%-15s = %12g  ; Reference pres in Z dir\n",
221           eoNames[eoPz],tcr->ref_value[eoPz]);
222   fprintf(fp,"%-15s = %12g  ; Polarizability used for the Epot correction\n",
223           eoNames[eoPolarizability],tcr->ref_value[eoPolarizability]);
224   fprintf(fp,"%-15s = %12g  ; Gas phase dipole moment used for Epot correction\n", 
225           eoNames[eoDipole],tcr->ref_value[eoDipole]);
226   fprintf(fp,"%-15s = %12d  ; Memory for coupling. Makes it converge faster.\n",
227           eoNames[eoMemory],tcr->nmemory);
228   fprintf(fp,"%-15s = %12s  ; Use intermolecular Epot only (LJ+Coul)\n",
229           eoNames[eoInter],yesno_names[tcr->bInter]);
230   fprintf(fp,"%-15s = %12s  ; Use virial iso pressure\n",
231           eoNames[eoUseVirial],yesno_names[tcr->bVirial]);
232   fprintf(fp,"%-15s = %12d  ; Combination rule, same coding as in grompp.\n",
233           eoNames[eoCombRule],tcr->combrule);
234   
235   fprintf(fp,"\n; Q-Coupling   %6s  %12s\n","type","xi");
236   for(i=0; (i<tcr->nQ); i++) {
237     fprintf(fp,"%-8s = %8s  %6d  %12g\n",
238             "Q",eoNames[tcr->tcQ[i].eObs],tcr->tcQ[i].at_i,tcr->tcQ[i].xi_Q);
239   }
240   
241   fprintf(fp,"\n; %8s %8s  %6s  %6s  %12s  %12s\n","Couple","To",
242           "i-type","j-type","xi-c6","xi-c12");
243   fprintf(fp,"; j-type == -1 means mixing rules will be applied!\n");
244   for(i=0; (i<tcr->nLJ); i++) {
245     fprintf(fp,"%-8s = %8s  %6d  %6d  %12g  %12g\n",
246             "LJ",eoNames[tcr->tcLJ[i].eObs],
247             tcr->tcLJ[i].at_i,tcr->tcLJ[i].at_j,
248             tcr->tcLJ[i].xi_6,tcr->tcLJ[i].xi_12);
249   }
250   
251   fprintf(fp,"\n; %8s %8s  %6s  %6s  %12s  %12s  %12s\n","Couple","To",
252           "i-type","j-type","xi-A","xi-B","xi-C");
253   fprintf(fp,"; j-type == -1 means mixing rules will be applied!\n");
254   for(i=0; (i<tcr->nBU); i++) {
255     fprintf(fp,"%-8s = %8s  %6d  %6d  %12g  %12g  %12g\n",
256             "BU",eoNames[tcr->tcBU[i].eObs],
257             tcr->tcBU[i].at_i,tcr->tcBU[i].at_j,
258             tcr->tcBU[i].xi_a,tcr->tcBU[i].xi_b,tcr->tcBU[i].xi_c);
259   }
260   
261   fprintf(fp,"\n; More Coupling\n");
262   for(i=0; (i<tcr->nIP); i++) {
263     ftype=idef->functype[tcr->tIP[i].type];
264     switch (ftype) {
265     case F_BONDS:
266       fprintf(fp,"%-15s = %-8s  %4d  %12g  %12g\n",
267               "Bonds",eoNames[tcr->tIP[i].eObs],tcr->tIP[i].type,
268               tcr->tIP[i].xi.harmonic.krA,
269               tcr->tIP[i].xi.harmonic.rA);
270       break;
271     default:
272       fprintf(stderr,"ftype %s not supported (yet)\n",
273               interaction_function[ftype].longname);
274     }
275   }
276   gmx_fio_fclose(fp);
277 }
278
279 static bool add_lj(int *nLJ,t_coupl_LJ **tcLJ,char *s,bool bObsUsed[])
280 {
281   int       j,ati,atj,eo;
282   char      buf[256];
283   double    xi6,xi12;
284   
285   if (sscanf(s,"%s%d%d%lf%lf",buf,&ati,&atj,&xi6,&xi12) != 5) 
286     return TRUE;
287   if ((eo=Name2eo(buf)) == -1)
288     gmx_fatal(FARGS,"Invalid observable for LJ coupling: %s",buf);
289   
290   for(j=0; (j<*nLJ); j++) {
291     if ((((*tcLJ)[j].at_i == ati) && ((*tcLJ)[j].at_j == atj)) &&
292         ((*tcLJ)[j].xi_6 || (*tcLJ)[j].xi_12) &&
293         ((*tcLJ)[j].eObs == eo))
294       break;
295   }
296   if (j == *nLJ) {
297     ++(*nLJ);
298     srenew((*tcLJ),*nLJ);
299   }
300   else
301     fprintf(stderr,"\n*** WARNING: overwriting entry for LJ coupling '%s'\n",s);
302   
303   clear_lj(&((*tcLJ)[j]));
304   if (((*tcLJ)[j].eObs = eo) == -1) {
305     gmx_fatal(FARGS,"Invalid observable for LJ coupling: %s",buf);
306   }
307   (*tcLJ)[j].at_i   = ati;
308   (*tcLJ)[j].at_j   = atj;
309   (*tcLJ)[j].xi_6   = xi6;
310   (*tcLJ)[j].xi_12  = xi12;
311   bObsUsed[eo] = TRUE;
312   
313   return FALSE;
314 }
315
316 static bool add_bu(int *nBU,t_coupl_BU **tcBU,char *s,bool bObsUsed[])
317 {
318   int       j,ati,atj,eo;
319   char      buf[256];
320   double    xia,xib,xic;
321   
322   if (sscanf(s,"%s%d%d%lf%lf%lf",buf,&ati,&atj,&xia,&xib,&xic) != 6) 
323     return TRUE;
324   if ((eo=Name2eo(buf)) == -1)
325     gmx_fatal(FARGS,"Invalid observable for BU coupling: %s",buf);
326   
327   for(j=0; (j<*nBU); j++) {
328     if ((((*tcBU)[j].at_i == ati) && ((*tcBU)[j].at_j == atj)) &&
329         ((*tcBU)[j].xi_a || (*tcBU)[j].xi_b || (*tcBU)[j].xi_c ) &&
330         ((*tcBU)[j].eObs == eo))
331       break;
332   }
333   if (j == *nBU) {
334     ++(*nBU);
335     srenew((*tcBU),*nBU);
336   }
337   else
338     fprintf(stderr,"\n*** WARNING: overwriting entry for BU coupling '%s'\n",s);
339   
340   clear_bu(&((*tcBU)[j]));
341   if (((*tcBU)[j].eObs = eo) == -1) {
342     gmx_fatal(FARGS,"Invalid observable for BU coupling: %s",buf);
343   }
344   (*tcBU)[j].at_i   = ati;
345   (*tcBU)[j].at_j   = atj;
346   (*tcBU)[j].xi_a   = xia;
347   (*tcBU)[j].xi_b   = xib;
348   (*tcBU)[j].xi_c   = xic;
349   bObsUsed[eo] = TRUE;
350
351   return FALSE;
352 }
353
354 static bool add_ip(int *nIP,t_coupl_iparams **tIP,char *s,int ftype,bool bObsUsed[])
355 {
356   int    i,eo,type;
357   char   buf[256];
358   double kb,b0;
359   
360   switch (ftype) {
361   case F_BONDS:
362     /* Pick out the type */
363     if (sscanf(s,"%s%d",buf,&type) != 2)
364       return TRUE;
365     if ((eo=Name2eo(buf)) == -1)
366       gmx_fatal(FARGS,"Invalid observable for IP coupling: %s",buf);
367       
368     /* Check whether this entry is there already */
369     for(i=0; (i<*nIP); i++) {
370       if ((*tIP)[i].type == type)
371         break;
372     }
373     if (i < *nIP) {
374       fprintf(stderr,"*** WARNING: overwriting entry for type %d\n",type);
375     }
376     else {
377       i=*nIP;
378       srenew((*tIP),i+1);
379       (*nIP)++;
380     }
381     if (sscanf(s,"%s%d%lf%lf",buf,&type,&kb,&b0) != 4)
382       return TRUE;
383     (*tIP)[i].type=type;
384     (*tIP)[i].eObs=eo;
385     (*tIP)[i].xi.harmonic.krA = kb;
386     (*tIP)[i].xi.harmonic.rA  = b0;
387     bObsUsed[eo] = TRUE;
388     break;
389   default:
390     fprintf(stderr,"ftype %s not supported (yet)\n",
391             interaction_function[ftype].longname);
392     return TRUE;
393   }
394   return FALSE;
395 }
396
397 static bool add_q(int *nQ,t_coupl_Q **tcQ,char *s,bool bObsUsed[])
398 {
399   int       j,ati,eo;
400   char      buf[256];
401   double    xiQ;
402   
403   if (sscanf(s,"%s%d%lf",buf,&ati,&xiQ) != 3) 
404     return TRUE;
405   
406   for(j=0; (j<*nQ); j++) {
407     if ((*tcQ)[j].at_i == ati)
408       break;
409   }
410   if (j == *nQ) {
411     ++(*nQ);
412     srenew((*tcQ),*nQ);
413   }
414   else
415     fprintf(stderr,"\n*** WARNING: overwriting entry for Q coupling '%s'\n",s);
416   
417   clear_q(&((*tcQ)[j]));
418   eo = (*tcQ)[j].eObs = Name2eo(buf);
419   if ((*tcQ)[j].eObs == -1) {
420     gmx_fatal(FARGS,"Invalid observable for Q coupling: %s",buf);
421   }
422   (*tcQ)[j].at_i   = ati;
423   (*tcQ)[j].xi_Q  = xiQ;
424   bObsUsed[eo] = TRUE;
425   
426   return FALSE;
427 }
428
429 void read_gct(const char *fn,t_coupl_rec *tcr)
430 {
431   warninp_t wi;
432   t_inpfile *inp;
433   int       i,j,ninp,nQ,nLJ,nBU,nIP;
434   bool      bWrong;
435   
436   wi = init_warning(FALSE,0);
437
438   inp=read_inpfile(fn,&ninp,NULL,wi);
439
440   for(i=0; (i<eoObsNR); i++) {
441     tcr->bObsUsed[i] = FALSE;
442     RTYPE (eoNames[i],  tcr->ref_value[i],      0.0);
443   }
444   ITYPE (eoNames[eoMemory],     tcr->nmemory,   1);
445   ETYPE (eoNames[eoInter],      tcr->bInter,    yesno_names);
446   ETYPE (eoNames[eoUseVirial],  tcr->bVirial,   yesno_names);
447   ITYPE (eoNames[eoCombRule],   tcr->combrule,  1);
448   tcr->tcLJ=NULL;
449   tcr->tcBU=NULL;
450   tcr->tcQ=NULL;
451   tcr->tIP=NULL;
452   nQ=nLJ=nBU=nIP=0;
453   
454   for(i=0; (i<ninp); i++) {
455     bWrong=FALSE;
456     if (gmx_strcasecmp(inp[i].name,"LJ") == 0) 
457       bWrong=add_lj(&nLJ,&(tcr->tcLJ),inp[i].value,tcr->bObsUsed);
458     else if (gmx_strcasecmp(inp[i].name,"BU") == 0) 
459       bWrong=add_bu(&nBU,&(tcr->tcBU),inp[i].value,tcr->bObsUsed);
460     else if (gmx_strcasecmp(inp[i].name,"Q") == 0) 
461       bWrong=add_q(&nQ,&(tcr->tcQ),inp[i].value,tcr->bObsUsed);
462     else if (gmx_strcasecmp(inp[i].name,"Bonds") == 0)
463       bWrong=add_ip(&nIP,&(tcr->tIP),inp[i].value,F_BONDS,tcr->bObsUsed);
464       
465     if (bWrong)
466       fprintf(stderr,"Wrong line in %s: '%s = %s'\n",
467               fn,inp[i].name,inp[i].value);
468     /*sfree(inp[i].name);
469       sfree(inp[i].value);*/
470   }
471   /* Check which ones have to be printed */
472   for(i=1; (i<nQ); i++)
473     for(j=0; (j<i); j++) {
474       if (tcr->tcQ[i].at_i == tcr->tcQ[j].at_i)
475         tcr->tcQ[j].bPrint=FALSE;
476     }
477   for(i=1; (i<nLJ); i++)
478     for(j=0; (j<i); j++) {
479       if (((tcr->tcLJ[i].at_i == tcr->tcLJ[j].at_i) &&
480            (tcr->tcLJ[i].at_j == tcr->tcLJ[j].at_j)) ||
481           ((tcr->tcLJ[i].at_i == tcr->tcLJ[j].at_j) &&
482            (tcr->tcLJ[i].at_j == tcr->tcLJ[j].at_i))) 
483         tcr->tcLJ[j].bPrint=FALSE;
484     }
485   
486   for(i=1; (i<nBU); i++)
487     for(j=0; (j<i); j++) {
488       if (((tcr->tcBU[i].at_i == tcr->tcBU[j].at_i) &&
489            (tcr->tcBU[i].at_j == tcr->tcBU[j].at_j)) ||
490           ((tcr->tcBU[i].at_i == tcr->tcBU[j].at_j) &&
491            (tcr->tcBU[i].at_j == tcr->tcBU[j].at_i))) 
492         tcr->tcBU[j].bPrint=FALSE;
493     }
494   
495   tcr->nQ  = nQ;
496   tcr->nLJ = nLJ;
497   tcr->nBU = nBU;
498   tcr->nIP = nIP;
499   
500   sfree(inp);
501
502   done_warning(wi,FARGS);
503 }
504