source: sasview/park-1.2.1/park/lib/resolution.c @ d7b46e71

ESS_GUIESS_GUI_DocsESS_GUI_batch_fittingESS_GUI_bumps_abstractionESS_GUI_iss1116ESS_GUI_iss879ESS_GUI_iss959ESS_GUI_openclESS_GUI_orderingESS_GUI_sync_sascalccostrafo411magnetic_scattrelease-4.1.1release-4.1.2release-4.2.2release_4.0.1ticket-1009ticket-1094-headlessticket-1242-2d-resolutionticket-1243ticket-1249ticket885unittest-saveload
Last change on this file since d7b46e71 was d36facf, checked in by Mathieu Doucet <doucetm@…>, 13 years ago

Re #5 fixing park compilation on MSVC (turning code into proper C)

  • Property mode set to 100644
File size: 9.4 KB
Line 
1/* This program is public domain. */
2
3#include <stdio.h>
4#include <assert.h>
5#include <string.h>
6#include <math.h>
7#ifdef SGI
8#include <ieeefp.h>
9#endif
10#if defined(_MSC_VER)
11double erf(double x)
12{
13    // constants
14    double a1 =  0.254829592;
15    double a2 = -0.284496736;
16    double a3 =  1.421413741;
17    double a4 = -1.453152027;
18    double a5 =  1.061405429;
19    double p  =  0.3275911;
20
21    // Save the sign of x
22    int sign = 1;
23
24    double t = 1.0;
25    double y = 1.0;
26    if (x < 0)
27        sign = -1;
28    x = fabs(x);
29
30    // A&S formula 7.1.26
31    t = 1.0/(1.0 + p*x);
32    y = 1.0 - (((((a5*t + a4)*t) + a3)*t + a2)*t + a1)*t*exp(-x*x);
33
34    return sign*y;
35}
36#endif
37
38/* What to do at the endpoints --- USE_TRUNCATED_NORMALIZATION will
39 * avoid the assumption that the data is zero where it hasn't been
40 * measured.
41 */
42#define USE_TRUNCATED_NORMALIZATION
43
44/* Computed using extended precision with Octave's symbolic toolbox. */
45#define PI4          12.56637061435917295385
46#define PI_180        0.01745329251994329576
47#define LN256         5.54517744447956247533
48#define SQRT2         1.41421356237309504880
49#define SQRT2PI       2.50662827463100050241
50
51/* Choose the resolution limit based on relative contribution to resolution
52 * rather than absolute contribution.  So if we ignore everything below,
53 * e.g. 0.1% of the peak, that limit occurs when G(x)/G(0) = 0.001 for
54 * gaussian G of width sima, or when x = sqrt(-2 sigma^2 log(0.001)). */
55#define LOG_RESLIMIT -6.90775527898213703123
56
57/** \file
58The resolution function returns the convolution of the reflectometry
59curve with a Q-dependent gaussian.
60
61We provide the following function:
62   resolution(Nin, Qin, Rin, N, Q, dQ, R)  returns convolution
63
64where
65   Nin is the number of theory points
66   Qin,Rin are the computed theory points
67   N is the number of Q points to calculate
68   Q are the locations of the measured data points
69   dQ are the width (sigma) of the convolution at each measured point
70   R is the returned convolution.
71
72Note that FWHM = sqrt(8 ln 2) dQ, so scale dQ appropriately.
73
74The contribution of Q to a resolution of width dQo at point Qo is:
75
76   p(Q) = 1/sqrt(2 pi dQo^2) exp ( (Q-Qo)^2/(2 dQo^2) )
77
78We are approximating the convolution at Qo using a numerical
79approximation to the integral over the measured points.  For
80efficiency, the integral is limited to p(Q_i)/p(0)>=0.001. 
81
82Note that the function we are convoluting is falling off as Q^4.
83That means the correct convolution should uniformly sample across
84the entire width of the Gaussian.  This is not possible at the
85end points unless you calculate the reflectivity beyond what is
86strictly needed for the data. For a given dQ and step size,
87you need enough steps that the following is true:
88
89    (n*step)^2 > -2 dQ^2 * ln 0.001
90
91The choice of sampling density is particularly important near the
92critical edge.  This is where the resolution calculation has the
93largest effect on the reflectivity curve. In one particular model,
94calculating every 0.001 rather than every 0.02 changed one value
95above the critical edge by 15%.  This is likely to be a problem for
96any system with a well defined critical edge.  The solution is to
97compute the theory function over a finer mesh where the derivative
98is changing rapidly.  For the critical edge, I have found a sampling
99density of 0.005 to be good enough.
100
101For systems involving thick layers, the theory function oscillates
102rapidly around the measured points.  This is a problem when the
103period of the oscillation, 2 pi/d for total sample depth d, is on
104the order of the width of the resolution function. This is true even
105for gradually changing profiles in materials with very high roughness
106values.  In these systems, the theory function should be oversampled
107around the measured points Q.  With a single thick layer, oversampling
108can be limited to just one period.  With multiple thick layers,
109oscillations will show interference patterns and it will be necessary
110to oversample uniformly between the measured points.  When oversampled
111spacing is less than about 2 pi/7 d, it is possible to see aliasing
112effects. 
113
114FIXME is it better to use random sampling or strictly
115regular spacing when you are undersampling?
116
117===============================================================
118*/
119
120#undef USE_TRAPEZOID_RULE
121#ifdef USE_TRAPEZOID_RULE
122#warning This code does strange things with small sigma and large spacing
123/* FIXME trapezoid performs very badly for large spacing unless
124   we normalize to the unit width.  For very small sigma, the gaussian
125   is a spike, but we are approximating it by a triangle so it is
126   not surprising it works so poorly.  A slightly better solution is
127   to use the inner limits rather than the outer limits, but this will
128   still break down if the Q spacing is approximately equal to limit.
129   Best is to not use trapezoid.
130*/
131
132/* Trapezoid rule for numerical integration of convolution */
133double 
134convolve_point(const double Qin[], const double Rin[], int k, int n,
135        double Qo, double limit, double sigma)
136{
137  const double two_sigma_sq = 2. * sigma * sigma;
138  double z, Glo, RGlo, R, norm;
139 
140  z = Qo - Qin[k];
141  Glo = exp(-z*z/two_sigma_sq);
142  RGlo = Rin[k]*Glo;
143  norm = R = 0.;
144  while (++k < n) {
145    /* Compute the next endpoint */
146    const double zhi = Qo - Qin[k];
147    const double Ghi = exp(-zhi*zhi/two_sigma_sq);
148    const double RGhi = Rin[k] * Ghi;
149    const double halfstep = 0.5*(Qin[k] - Qin[k-1]);
150   
151    /* Add the trapezoidal area. */
152    norm += halfstep * (Ghi + Glo);
153    R += halfstep * (RGhi + RGlo);
154   
155    /* Save the endpoint for next trapezoid. */
156    Glo = Ghi;
157    RGlo = RGhi;
158   
159    /* Check if we've calculated far enough */
160    if (Qin[k] >= Qo+limit) break;
161  }
162
163  /* Scale to area of the linear spline distribution we actually used. */
164  return R / norm;
165
166  /* Scale to gaussian of unit area */
167  /* Fails badly for small sigma or large Q steps---do not use. */
168  /* return R / sigma*SQRT2PI; */
169}
170
171#else /* !USE_TRAPEZOID_RULE */
172
173/* Analytic convolution of gaussian with linear spline */
174/* More expensive but more reliable */
175double 
176convolve_point(const double Qin[], const double Rin[], int k, int n,
177               double Qo, double limit, double sigma)
178{
179  const double two_sigma_sq = 2. * sigma * sigma;
180  double z, Glo, erflo, erfmin, R;
181  double zhi, Ghi, erfhi, m, b;
182 
183  z = Qo - Qin[k];
184  Glo = exp(-z*z/two_sigma_sq);
185  erfmin = erflo = erf(-z/(SQRT2*sigma));
186  R = 0.;
187  /* printf("%5.3f: (%5.3f,%11.5g)",Qo,Qin[k],Rin[k]); */
188  while (++k < n) {
189        /* No additional contribution from duplicate points. */
190        if (Qin[k] == Qin[k-1]) continue;
191 
192    /* Compute the next endpoint */
193    zhi = Qo - Qin[k];
194    Ghi = exp(-zhi*zhi/two_sigma_sq);
195    erfhi = erf(-zhi/(SQRT2*sigma));
196    m = (Rin[k]-Rin[k-1])/(Qin[k]-Qin[k-1]);
197    b = Rin[k] - m * Qin[k];
198
199    /* Add the integrals. */
200    R += 0.5*(m*Qo+b)*(erfhi-erflo) - sigma/SQRT2PI*m*(Ghi-Glo);
201
202    /* Debug computation failures. */
203    /*
204    if isnan(R) {
205        printf("NaN from %d: zhi=%g, Ghi=%g, erfhi=%g, m=%g, b=%g\n",
206               k,zhi,Ghi,erfhi,m,b);
207    }
208    */
209   
210    /* Save the endpoint for next trapezoid. */
211    Glo = Ghi;
212    erflo = erfhi;
213   
214    /* Check if we've calculated far enough */
215    if (Qin[k] >= Qo+limit) break;
216  }
217  /* printf(" (%5.3f,%11.5g)",Qin[k<n?k:n-1],Rin[k<n?k:n-1]); */
218
219#ifdef USE_TRUNCATED_NORMALIZATION
220  /* Normalize by the area of the truncated gaussian */
221  /* At this point erflo = erfmax */
222  /* printf ("---> %11.5g\n",2*R/(erflo-erfmin)); */
223  return 2 * R / (erflo - erfmin);
224#else
225  /* Return unnormalized (we used a gaussian of unit area) */
226  /* printf ("---> %11.5g\n",R); */
227  return R;
228#endif
229}
230
231#endif /* !USE_TRAPEZOID_RULE */
232
233void
234resolution(int Nin, const double Qin[], const double Rin[],
235           int N, const double Q[], const double dQ[], double R[])
236{
237  int lo,out;
238  double sigma, Qo, limit;
239
240  /* FIXME fails if Qin are not sorted; slow if Q not sorted */
241  assert(Nin>1);
242
243  /* Scan through all Q values to be calculated */
244  lo = 0;
245  for (out=0; out < N; out++) {
246    /* width of resolution window for Q is w = 2 dQ^2. */
247    sigma = dQ[out];
248    Qo = Q[out];
249    limit = sqrt(-2.*sigma*sigma* LOG_RESLIMIT);
250
251    /* if (out%20==0) printf("%d: Q,dQ = %g,%g\n",out,Qo,sigma); */
252
253    /* Line up the left edge of the convolution window */
254    /* It is probably forward from the current position, */
255    /* but if the next dQ is a lot higher than the current */
256    /* dQ or if the Q are not sorted, then it may be before */
257    /* the current position. */
258    /* FIXME verify that the convolution window is just right */
259    while (lo < Nin-1 && Qin[lo] < Qo-limit) lo++;
260    while (lo > 0 && Qin[lo] > Qo-limit) lo--;
261
262    /* Special handling to avoid 0/0 for w=0. */
263    if (sigma > 0.) {
264      R[out] = convolve_point(Qin,Rin,lo,Nin,Qo,limit,sigma);
265    } else if (lo < Nin-1) {
266      /* Linear interpolation */
267      double m = (Rin[lo+1]-Rin[lo])/(Qin[lo+1]-Qin[lo]);
268      double b = Rin[lo] - m*Qin[lo];
269      R[out] = m*Qo + b;
270    } else if (lo > 0) {
271      /* Linear extrapolation */
272      double m = (Rin[lo]-Rin[lo-1])/(Qin[lo]-Qin[lo-1]);
273      double b = Rin[lo] - m*Qin[lo];
274      R[out] = m*Qo + b;
275    } else {
276      /* Can't happen because there is more than one point in Qin. */
277      assert(Nin>1);
278    }
279  }
280
281}
Note: See TracBrowser for help on using the repository browser.