Merge release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / gmxlib / sfactor.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  * GROningen Mixture of Alchemy and Childrens' Stories
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include <ctype.h>
40 #include "sysstuff.h"
41 #include "smalloc.h"
42 #include "string2.h"
43 #include "futil.h"
44 #include "maths.h"
45 #include "gmx_fatal.h"
46 #include "vec.h"
47 #include "macros.h"
48 #include "index.h"
49 #include "strdb.h"
50 #include "copyrite.h"
51 #include "tpxio.h"
52 #include "typedefs.h"
53 #include "statutil.h"
54 #include "oenv.h"
55 #include "gmxfio.h"
56 #include "xvgr.h"
57 #include "matio.h"
58 #include "gmx_ana.h"
59 #include "names.h"
60 #include "sfactor.h"
61
62
63 typedef struct gmx_structurefactors {
64     int nratoms;
65     int *p;        /* proton number */
66     int *n;        /* neutron number */
67     /* Parameters for the Cromer Mann fit */
68     real **a;    /* parameter a */
69     real **b;    /* parameter b */
70     real *c;       /* parameter c */
71     char **atomnm;  /* atomname */
72
73 } gmx_structurefactors;
74
75 typedef struct reduced_atom{
76     rvec x;
77     int  t;
78 } reduced_atom;
79
80
81 typedef struct structure_factor
82 {
83   int     n_angles;
84   int     n_groups;
85   double  lambda;
86   double  energy;
87   double  momentum;
88   double  ref_k;
89   double  **F;
90   int     nSteps;
91   int     total_n_atoms;
92 } structure_factor;
93
94
95 extern int * create_indexed_atom_type (reduced_atom_t * atm, int size)
96 {
97 /*
98  * create an index of the atom types found in a  group
99  * i.e.: for water index_atp[0]=type_number_of_O and
100  *                 index_atp[1]=type_number_of_H
101  *
102  * the last element is set to 0
103  */
104     int *index_atp, i, i_tmp, j;
105
106     reduced_atom *att=(reduced_atom *)atm;
107
108     snew (index_atp, 1);
109     i_tmp = 1;
110     index_atp[0] = att[0].t;
111     for (i = 1; i < size; i++) {
112         for (j = 0; j < i_tmp; j++)
113             if (att[i].t == index_atp[j])
114                 break;
115         if (j == i_tmp) {       /* i.e. no indexed atom type is  == to atm[i].t */
116             i_tmp++;
117             srenew (index_atp, i_tmp * sizeof (int));
118             index_atp[i_tmp - 1] = att[i].t;
119         }
120     }
121     i_tmp++;
122     srenew (index_atp, i_tmp * sizeof (int));
123     index_atp[i_tmp - 1] = 0;
124     return index_atp;
125 }
126
127
128
129 extern t_complex *** rc_tensor_allocation(int x, int y, int z)
130 {
131   t_complex ***t;
132   int i,j;
133
134   t = (t_complex ***)calloc(x,sizeof(t_complex**));
135   if(!t) exit(fprintf(stderr,"\nallocation error"));
136   t[0] = (t_complex **)calloc(x*y,sizeof(t_complex*));
137   if(!t[0]) exit(fprintf(stderr,"\nallocation error"));
138   t[0][0] = (t_complex *)calloc(x*y*z,sizeof(t_complex));
139   if(!t[0][0]) exit(fprintf(stderr,"\nallocation error"));
140
141   for( j = 1 ; j < y ; j++)
142     t[0][j] = t[0][j-1] + z;
143   for( i = 1 ; i < x ; i++) {
144     t[i] = t[i-1] + y;
145     t[i][0] = t[i-1][0] + y*z;
146     for( j = 1 ; j < y ; j++)
147       t[i][j] = t[i][j-1] + z;
148   }
149   return t;
150 }
151
152
153 extern void compute_structure_factor (structure_factor_t * sft, matrix box,
154                                reduced_atom_t * red, int isize, real start_q,
155                                real end_q, int group,real **sf_table)
156 {
157     structure_factor *sf=(structure_factor *)sft;
158     reduced_atom *redt=(reduced_atom *)red;
159
160     t_complex ***tmpSF;
161     rvec k_factor;
162     real kdotx, asf, kx, ky, kz, krr;
163     int kr, maxkx, maxky, maxkz, i, j, k, p, *counter;
164
165
166     k_factor[XX] = 2 * M_PI / box[XX][XX];
167     k_factor[YY] = 2 * M_PI / box[YY][YY];
168     k_factor[ZZ] = 2 * M_PI / box[ZZ][ZZ];
169
170     maxkx = (int) (end_q / k_factor[XX] + 0.5);
171     maxky = (int) (end_q / k_factor[YY] + 0.5);
172     maxkz = (int) (end_q / k_factor[ZZ] + 0.5);
173
174     snew (counter, sf->n_angles);
175
176     tmpSF = rc_tensor_allocation(maxkx,maxky,maxkz);
177 /*
178  * The big loop...
179  * compute real and imaginary part of the structure factor for every
180  * (kx,ky,kz))
181  */
182     fprintf(stderr,"\n");
183     for (i = 0; i < maxkx; i++) {
184         fprintf (stderr,"\rdone %3.1f%%     ", (double)(100.0*(i+1))/maxkx);
185         kx = i * k_factor[XX];
186         for (j = 0; j < maxky; j++) {
187             ky = j * k_factor[YY];
188             for (k = 0; k < maxkz; k++)
189                 if (i != 0 || j != 0 || k != 0) {
190                     kz = k * k_factor[ZZ];
191                     krr = sqrt (sqr (kx) + sqr (ky) + sqr (kz));
192                     if (krr >= start_q && krr <= end_q) {
193                         kr = (int) (krr/sf->ref_k + 0.5);
194                         if (kr < sf->n_angles) {
195                             counter[kr]++;  /* will be used for the copmutation
196                                                of the average*/
197                             for (p = 0; p < isize; p++) {
198                                 asf = sf_table[redt[p].t][kr];
199
200                                 kdotx = kx * redt[p].x[XX] +
201                                     ky * redt[p].x[YY] + kz * redt[p].x[ZZ];
202
203                                 tmpSF[i][j][k].re += cos (kdotx) * asf;
204                                 tmpSF[i][j][k].im += sin (kdotx) * asf;
205                             }
206                         }
207                     }
208                 }
209         }
210     }                           /* end loop on i */
211 /*
212  *  compute the square modulus of the structure factor, averaging on the surface
213  *  kx*kx + ky*ky + kz*kz = krr*krr
214  *  note that this is correct only for a (on the macroscopic scale)
215  *  isotropic system.
216  */
217     for (i = 0; i < maxkx; i++) {
218         kx = i * k_factor[XX]; for (j = 0; j < maxky; j++) {
219             ky = j * k_factor[YY]; for (k = 0; k < maxkz; k++) {
220                 kz = k * k_factor[ZZ]; krr = sqrt (sqr (kx) + sqr (ky)
221                 + sqr (kz)); if (krr >= start_q && krr <= end_q) {
222                     kr = (int) (krr / sf->ref_k + 0.5);
223                         if (kr < sf->n_angles && counter[kr] != 0)
224                                 sf->F[group][kr] +=
225                             (sqr (tmpSF[i][j][k].re) +
226                              sqr (tmpSF[i][j][k].im))/ counter[kr];
227                 }
228             }
229         }
230     } sfree (counter); free(tmpSF[0][0]); free(tmpSF[0]); free(tmpSF);
231 }
232
233
234 extern gmx_structurefactors_t *gmx_structurefactors_init(const char *datfn) {
235     
236         /* Read the database for the structure factor of the different atoms */
237
238     FILE *fp;
239     char line[STRLEN];
240     gmx_structurefactors *gsf;
241     double a1,a2,a3,a4,b1,b2,b3,b4,c;
242     int p;
243     int i;
244     int nralloc=10;
245     int line_no;
246     char atomn[32];
247     fp=libopen(datfn);
248     line_no = 0;
249     snew(gsf,1);
250
251     snew(gsf->atomnm,nralloc);
252     snew(gsf->a,nralloc);
253     snew(gsf->b,nralloc);
254     snew(gsf->c,nralloc);
255     snew(gsf->p,nralloc);
256     gsf->n=NULL;
257     gsf->nratoms=line_no;
258     while(get_a_line(fp,line,STRLEN)) {
259         i=line_no;
260         if (sscanf(line,"%s %d %lf %lf %lf %lf %lf %lf %lf %lf %lf",
261                    atomn,&p,&a1,&a2,&a3,&a4,&b1,&b2,&b3,&b4,&c) == 11) {
262             gsf->atomnm[i]=strdup(atomn);
263             gsf->p[i]=p;
264             snew(gsf->a[i],4);
265             snew(gsf->b[i],4);
266             gsf->a[i][0]=a1;
267             gsf->a[i][1]=a2;
268             gsf->a[i][2]=a3;
269             gsf->a[i][3]=a4;
270             gsf->b[i][0]=b1;
271             gsf->b[i][1]=b2;
272             gsf->b[i][2]=b3;
273             gsf->b[i][3]=b4;
274             gsf->c[i]=c;
275             line_no++;            
276             gsf->nratoms=line_no;
277             if (line_no==nralloc){
278                 nralloc+=10;
279                 srenew(gsf->atomnm,nralloc);
280                 srenew(gsf->a,nralloc);
281                 srenew(gsf->b,nralloc);
282                 srenew(gsf->c,nralloc);
283                 srenew(gsf->p,nralloc);
284             }
285         }
286         else 
287             fprintf(stderr,"WARNING: Error in file %s at line %d ignored\n",
288                     datfn,line_no);
289     }
290
291     srenew(gsf->atomnm,gsf->nratoms);
292     srenew(gsf->a,gsf->nratoms);
293     srenew(gsf->b,gsf->nratoms);
294     srenew(gsf->c,gsf->nratoms);
295     srenew(gsf->p,gsf->nratoms);
296
297     fclose(fp);
298
299     return (gmx_structurefactors_t *) gsf;
300
301 }
302
303
304 extern void rearrange_atoms (reduced_atom_t * positions, t_trxframe *fr, atom_id * index,
305                       int isize, t_topology * top, gmx_bool flag,gmx_structurefactors_t *gsf)
306 /* given the group's index, return the (continuous) array of atoms */
307 {
308   int i;
309
310   reduced_atom *pos=(reduced_atom *)positions;
311
312   if (flag)
313     for (i = 0; i < isize; i++)
314       pos[i].t =
315         return_atom_type (*(top->atoms.atomname[index[i]]),gsf);
316   for (i = 0; i < isize; i++)
317     copy_rvec (fr->x[index[i]], pos[i].x);
318
319    positions=(reduced_atom_t *)pos;
320 }
321
322
323 extern int return_atom_type (const char *name,gmx_structurefactors_t *gsf)
324 {
325   typedef struct {
326     const char *name;
327     int  nh;
328   } t_united_h;
329   t_united_h uh[] = {
330     { "CH1", 1 }, { "CH2", 2 }, { "CH3", 3 },
331     { "CS1", 1 }, { "CS2", 2 }, { "CS3", 3 },
332     { "CP1", 1 }, { "CP2", 2 }, { "CP3", 3 }
333   };
334   int i,cnt=0;
335   int *tndx;
336   int nrc;
337   int fndx=0;
338   int NCMT;
339
340   gmx_structurefactors *gsft=(gmx_structurefactors *)gsf;
341
342   NCMT=gsft->nratoms;
343
344   snew(tndx,NCMT);
345
346   for(i=0; (i<asize(uh)); i++)
347     if (strcmp(name,uh[i].name) == 0)
348       return NCMT-1+uh[i].nh;
349
350   for(i=0; (i<NCMT); i++){
351     if (strncmp (name, gsft->atomnm[i],strlen(gsft->atomnm[i])) == 0){
352       tndx[cnt]=i;
353       cnt++;
354     }
355   }
356
357   if (cnt==0)
358   gmx_fatal(FARGS,"\nError: atom (%s) not in list (%d types checked)!\n",
359             name,i);
360   else{
361           nrc=0;
362           for(i=0;i<cnt;i++){
363                   if(strlen(gsft->atomnm[tndx[i]])>(size_t)nrc){
364                           nrc=strlen(gsft->atomnm[tndx[i]]);
365                       fndx=tndx[i];
366                   }
367           }
368      
369           return fndx;
370   }
371
372   return 0;
373 }
374
375 extern int gmx_structurefactors_get_sf(gmx_structurefactors_t *gsf, int elem, real a[4], real b[4], real *c){
376
377         int success;
378         int i;
379     gmx_structurefactors *gsft=(gmx_structurefactors *)gsf;
380         success=0;
381
382         for(i=0;i<4;i++){
383             a[i]=gsft->a[elem][i];
384             b[i]=gsft->b[elem][i];
385             *c=gsft->c[elem];
386         }
387
388         success+=1;
389         return success;
390 }
391
392 extern int do_scattering_intensity (const char* fnTPS, const char* fnNDX,
393                              const char* fnXVG, const char *fnTRX,
394                              const char* fnDAT,
395                              real start_q,real end_q,
396                              real energy,int ng,const output_env_t oenv)
397 {
398     int i,*isize,flags = TRX_READ_X,**index_atp;
399     t_trxstatus *status;
400     char **grpname,title[STRLEN];
401     atom_id **index;
402     t_topology top;
403     int ePBC;
404     t_trxframe fr;
405     reduced_atom_t **red;
406     structure_factor *sf;
407     rvec *xtop;
408     real **sf_table;
409     int nsftable;
410     matrix box;
411     double r_tmp;
412
413     gmx_structurefactors_t *gmx_sf;
414     real *a,*b,c;
415     int success;
416
417     snew(a,4);
418     snew(b,4);
419
420
421     gmx_sf=gmx_structurefactors_init(fnDAT);
422
423     success=gmx_structurefactors_get_sf(gmx_sf,0, a, b, &c);
424
425     snew (sf, 1);
426     sf->energy = energy;
427
428     /* Read the topology informations */
429     read_tps_conf (fnTPS, title, &top, &ePBC, &xtop, NULL, box, TRUE);
430     sfree (xtop);
431
432     /* groups stuff... */
433     snew (isize, ng);
434     snew (index, ng);
435     snew (grpname, ng);
436
437     fprintf (stderr, "\nSelect %d group%s\n", ng,
438              ng == 1 ? "" : "s");
439     if (fnTPS)
440         get_index (&top.atoms, fnNDX, ng, isize, index, grpname);
441     else
442         rd_index (fnNDX, ng, isize, index, grpname);
443
444     /* The first time we read data is a little special */
445     read_first_frame (oenv,&status, fnTRX, &fr, flags);
446
447     sf->total_n_atoms = fr.natoms;
448
449     snew (red, ng);
450     snew (index_atp, ng);
451
452     r_tmp = max (box[XX][XX], box[YY][YY]);
453     r_tmp = (double) max (box[ZZ][ZZ], r_tmp);
454
455     sf->ref_k = (2.0 * M_PI) / (r_tmp);
456     /* ref_k will be the reference momentum unit */
457     sf->n_angles = (int) (end_q / sf->ref_k + 0.5);
458
459     snew (sf->F, ng);
460     for (i = 0; i < ng; i++)
461         snew (sf->F[i], sf->n_angles);
462     for (i = 0; i < ng; i++) {
463         snew (red[i], isize[i]);
464         rearrange_atoms (red[i], &fr, index[i], isize[i], &top, TRUE,gmx_sf);
465         index_atp[i] = create_indexed_atom_type (red[i], isize[i]);
466     }
467
468     sf_table = compute_scattering_factor_table (gmx_sf,(structure_factor_t *)sf,&nsftable);
469
470
471     /* This is the main loop over frames */
472
473     do {
474         sf->nSteps++;
475         for (i = 0; i < ng; i++) {
476             rearrange_atoms (red[i], &fr, index[i], isize[i], &top,FALSE,gmx_sf);
477
478             compute_structure_factor ((structure_factor_t *)sf, box, red[i], isize[i],
479                                       start_q, end_q, i, sf_table);
480         }
481     }
482
483     while (read_next_frame (oenv,status, &fr));
484
485     save_data ((structure_factor_t *)sf, fnXVG, ng, start_q, end_q,oenv);
486
487
488     sfree(a);
489     sfree(b);
490
491     gmx_structurefactors_done(gmx_sf);
492
493     return 0;
494 }
495
496
497 extern void save_data (structure_factor_t *sft, const char *file, int ngrps,
498                 real start_q, real end_q, const output_env_t oenv)
499 {
500
501     FILE *fp;
502     int i, g = 0;
503     double *tmp, polarization_factor, A;
504
505     structure_factor *sf=(structure_factor *)sft;
506     
507     fp = xvgropen (file, "Scattering Intensity", "q (1/nm)",
508                    "Intensity (a.u.)",oenv);
509
510     snew (tmp, ngrps);
511
512     for (g = 0; g < ngrps; g++)
513         for (i = 0; i < sf->n_angles; i++) {
514             
515 /*
516  *          theta is half the angle between incoming and scattered vectors.
517  *
518  *          polar. fact. = 0.5*(1+cos^2(2*theta)) = 1 - 0.5 * sin^2(2*theta)
519  *
520  *          sin(theta) = q/(2k) := A  ->  sin^2(theta) = 4*A^2 (1-A^2) ->
521  *          -> 0.5*(1+cos^2(2*theta)) = 1 - 2 A^2 (1-A^2)
522  */
523             A = (double) (i * sf->ref_k) / (2.0 * sf->momentum);
524             polarization_factor = 1 - 2.0 * sqr (A) * (1 - sqr (A));
525             sf->F[g][i] *= polarization_factor;
526         }
527     for (i = 0; i < sf->n_angles; i++) {
528         if (i * sf->ref_k >= start_q && i * sf->ref_k <= end_q) {
529             fprintf (fp, "%10.5f  ", i * sf->ref_k);
530             for (g = 0; g < ngrps; g++)
531                fprintf (fp, "  %10.5f ", (sf->F[g][i]) /( sf->total_n_atoms*
532                                                           sf->nSteps));
533             fprintf (fp, "\n");
534         }
535     }
536
537     ffclose (fp);
538 }
539
540
541 extern double CMSF (gmx_structurefactors_t *gsf,int type,int nh,double lambda, double sin_theta)
542 /*
543  * return Cromer-Mann fit for the atomic scattering factor:
544  * sin_theta is the sine of half the angle between incoming and scattered
545  * vectors. See g_sq.h for a short description of CM fit.
546  */
547 {
548   int i,success;
549   double tmp = 0.0, k2;
550   real *a,*b;
551   real c;
552
553   snew(a,4);
554   snew(b,4);
555
556  /*
557   *
558   * f0[k] = c + [SUM a_i*EXP(-b_i*(k^2)) ]
559   *             i=1,4
560   */
561
562   /*
563    *  united atoms case
564    *  CH2 / CH3 groups
565    */
566   if (nh > 0) {
567     tmp = (CMSF (gsf,return_atom_type ("C",gsf),0,lambda, sin_theta) +
568            nh*CMSF (gsf,return_atom_type ("H",gsf),0,lambda, sin_theta));
569   }
570   /* all atom case */
571   else {
572     k2 = (sqr (sin_theta) / sqr (10.0 * lambda));
573     success=gmx_structurefactors_get_sf(gsf,type,a,b,&c);
574     tmp = c;
575     for (i = 0; (i < 4); i++)
576       tmp += a[i] * exp (-b[i] * k2);
577   }
578   return tmp;
579 }
580
581
582
583 extern real **gmx_structurefactors_table(gmx_structurefactors_t *gsf,real momentum, real ref_k, real lambda, int n_angles){
584
585         int NCMT;
586         int nsftable;
587         int i,j;
588         double q,sin_theta;
589     real **sf_table;
590     gmx_structurefactors *gsft=(gmx_structurefactors *)gsf;
591
592     NCMT=gsft->nratoms;
593         nsftable = NCMT+3;
594
595     snew (sf_table,nsftable);
596     for (i = 0; (i < nsftable); i++) {
597         snew (sf_table[i], n_angles);
598         for (j = 0; j < n_angles; j++) {
599                 q = ((double) j * ref_k);
600                 /* theta is half the angle between incoming
601                            and scattered wavevectors. */
602                 sin_theta = q / (2.0 * momentum);
603                 if (i < NCMT){
604                     sf_table[i][j] = CMSF (gsf,i,0,lambda, sin_theta);
605                 }
606                 else
607                         sf_table[i][j] = CMSF (gsf,i,i-NCMT+1,lambda, sin_theta);
608         }
609     }
610     return sf_table;
611 }
612
613 extern void gmx_structurefactors_done(gmx_structurefactors_t *gsf){
614
615         int i;
616         gmx_structurefactors *sf;
617         sf=(gmx_structurefactors *) gsf;
618
619         for(i=0;i<sf->nratoms;i++){
620         sfree(sf->a[i]);
621                 sfree(sf->b[i]);
622                 sfree(sf->atomnm[i]);
623         }
624
625         sfree(sf->a);
626         sfree(sf->b);
627         sfree(sf->atomnm);
628         sfree(sf->p);
629         sfree(sf->c);
630
631         sfree(sf);
632
633 }
634
635 extern real **compute_scattering_factor_table (gmx_structurefactors_t *gsf,structure_factor_t *sft,int *nsftable)
636 {
637 /*
638  *  this function build up a table of scattering factors for every atom
639  *  type and for every scattering angle.
640  */
641    
642     double hc=1239.842;
643     real ** sf_table;
644
645     structure_factor *sf=(structure_factor *)sft;
646
647
648     /* \hbar \omega \lambda = hc = 1239.842 eV * nm */
649     sf->momentum = ((double) (2. * 1000.0 * M_PI * sf->energy) / hc);
650     sf->lambda = hc / (1000.0 * sf->energy);
651     fprintf (stderr, "\nwavelenght = %f nm\n", sf->lambda);
652
653     sf_table=gmx_structurefactors_table(gsf,sf->momentum,sf->ref_k,sf->lambda,sf->n_angles);
654
655     return sf_table;
656 }
657
658