source: sasview/sansmodels/src/c_models/cylinder.cpp @ f518eef

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 f518eef was 20efe7b, checked in by Jae Cho <jhjcho@…>, 12 years ago

block alpha=0 that breaks fitting

  • Property mode set to 100644
File size: 9.5 KB
Line 
1/**
2        This software was developed by the University of Tennessee as part of the
3        Distributed Data Analysis of Neutron Scattering Experiments (DANSE)
4        project funded by the US National Science Foundation.
5
6        If you use DANSE applications to do scientific research that leads to
7        publication, we ask that you acknowledge the use of the software with the
8        following sentence:
9
10        "This work benefited from DANSE software developed under NSF award DMR-0520547."
11
12        copyright 2008, University of Tennessee
13 */
14
15/**
16 * Scattering model classes
17 * The classes use the IGOR library found in
18 *   sansmodels/src/libigor
19 *
20 */
21
22#include <math.h>
23#include "parameters.hh"
24#include <stdio.h>
25using namespace std;
26
27extern "C" {
28        #include "libCylinder.h"
29        #include "libStructureFactor.h"
30}
31#include "cylinder.h"
32
33// Convenience parameter structure
34typedef struct {
35    double scale;
36    double radius;
37    double length;
38    double sldCyl;
39    double sldSolv;
40    double background;
41    double cyl_theta;
42    double cyl_phi;
43} CylinderParameters;
44
45CylinderModel :: CylinderModel() {
46        scale      = Parameter(1.0);
47        radius     = Parameter(20.0, true);
48        radius.set_min(0.0);
49        length     = Parameter(400.0, true);
50        length.set_min(0.0);
51        sldCyl   = Parameter(4.e-6);
52        sldSolv   = Parameter(1.e-6);
53        background = Parameter(0.0);
54        cyl_theta  = Parameter(0.0, true);
55        cyl_phi    = Parameter(0.0, true);
56}
57
58/**
59 * Function to evaluate 1D scattering function
60 * The NIST IGOR library is used for the actual calculation.
61 * @param q: q-value
62 * @return: function value
63 */
64double CylinderModel :: operator()(double q) {
65        double dp[6];
66
67        // Fill parameter array for IGOR library
68        // Add the background after averaging
69        dp[0] = scale();
70        dp[1] = radius();
71        dp[2] = length();
72        dp[3] = sldCyl();
73        dp[4] = sldSolv();
74        dp[5] = 0.0;
75
76        // Get the dispersion points for the radius
77        vector<WeightPoint> weights_rad;
78        radius.get_weights(weights_rad);
79
80        // Get the dispersion points for the length
81        vector<WeightPoint> weights_len;
82        length.get_weights(weights_len);
83
84        // Perform the computation, with all weight points
85        double sum = 0.0;
86        double norm = 0.0;
87        double vol = 0.0;
88
89        // Loop over radius weight points
90        for(size_t i=0; i<weights_rad.size(); i++) {
91                dp[1] = weights_rad[i].value;
92
93                // Loop over length weight points
94                for(size_t j=0; j<weights_len.size(); j++) {
95                        dp[2] = weights_len[j].value;
96                        //Un-normalize by volume
97                        sum += weights_rad[i].weight
98                                * weights_len[j].weight * CylinderForm(dp, q)
99                                *pow(weights_rad[i].value,2)*weights_len[j].value;
100
101                        //Find average volume
102                        vol += weights_rad[i].weight
103                                * weights_len[j].weight *pow(weights_rad[i].value,2)*weights_len[j].value;
104                        norm += weights_rad[i].weight
105                                * weights_len[j].weight;
106                }
107        }
108        if (vol != 0.0 && norm != 0.0) {
109                //Re-normalize by avg volume
110                sum = sum/(vol/norm);}
111
112        return sum/norm + background();
113}
114
115/**
116 * Function to evaluate 2D scattering function
117 * @param pars: parameters of the cylinder
118 * @param q: q-value
119 * @param q_x: q_x / q
120 * @param q_y: q_y / q
121 * @return: function value
122 */
123static double cylinder_analytical_2D_scaled(CylinderParameters *pars, double q, double q_x, double q_y) {
124  double cyl_x, cyl_y, cyl_z;
125  double q_z;
126  double alpha, vol, cos_val;
127  double answer;
128  //convert angle degree to radian
129  double pi = 4.0*atan(1.0);
130  double theta = pars->cyl_theta * pi/180.0;
131  double phi = pars->cyl_phi * pi/180.0;
132
133    // Cylinder orientation
134    cyl_x = sin(theta) * cos(phi);
135    cyl_y = sin(theta) * sin(phi);
136    cyl_z = cos(theta);
137
138    // q vector
139    q_z = 0.0;
140
141    // Compute the angle btw vector q and the
142    // axis of the cylinder
143    cos_val = cyl_x*q_x + cyl_y*q_y + cyl_z*q_z;
144
145    // The following test should always pass
146    if (fabs(cos_val)>1.0) {
147      printf("cyl_ana_2D: Unexpected error: |cos(alpha)=%g|>1\n", cos_val);
148      printf("cyl_ana_2D: at theta=%g and phi=%g.", theta, phi);
149      return 1.0;
150    }
151
152    // Note: cos(alpha) = 0 and 1 will get an
153    // undefined value from CylKernel
154  alpha = acos( cos_val );
155  if (alpha == 0.0){
156        alpha = 1.0e-26;
157        }
158  // Call the IGOR library function to get the kernel
159  answer = CylKernel(q, pars->radius, pars->length/2.0, alpha) / sin(alpha);
160
161  // Multiply by contrast^2
162  answer *= (pars->sldCyl - pars->sldSolv)*(pars->sldCyl - pars->sldSolv);
163
164  //normalize by cylinder volume
165  //NOTE that for this (Fournet) definition of the integral, one must MULTIPLY by Vcyl
166    vol = acos(-1.0) * pars->radius * pars->radius * pars->length;
167  answer *= vol;
168
169  //convert to [cm-1]
170  answer *= 1.0e8;
171
172  //Scale
173  answer *= pars->scale;
174
175  // add in the background
176  answer += pars->background;
177
178  return answer;
179}
180
181/**
182 * Function to evaluate 2D scattering function
183 * @param pars: parameters of the cylinder
184 * @param q: q-value
185 * @return: function value
186 */
187static double cylinder_analytical_2DXY(CylinderParameters *pars, double qx, double qy) {
188  double q;
189  q = sqrt(qx*qx+qy*qy);
190  return cylinder_analytical_2D_scaled(pars, q, qx/q, qy/q);
191}
192
193/**
194 * Function to evaluate 2D scattering function
195 * @param q_x: value of Q along x
196 * @param q_y: value of Q along y
197 * @return: function value
198 */
199double CylinderModel :: operator()(double qx, double qy) {
200        CylinderParameters dp;
201        // Fill parameter array
202        dp.scale      = scale();
203        dp.radius     = radius();
204        dp.length     = length();
205        dp.sldCyl   = sldCyl();
206        dp.sldSolv   = sldSolv();
207        dp.background = 0.0;
208        dp.cyl_theta  = cyl_theta();
209        dp.cyl_phi    = cyl_phi();
210
211        // Get the dispersion points for the radius
212        vector<WeightPoint> weights_rad;
213        radius.get_weights(weights_rad);
214
215        // Get the dispersion points for the length
216        vector<WeightPoint> weights_len;
217        length.get_weights(weights_len);
218
219        // Get angular averaging for theta
220        vector<WeightPoint> weights_theta;
221        cyl_theta.get_weights(weights_theta);
222
223        // Get angular averaging for phi
224        vector<WeightPoint> weights_phi;
225        cyl_phi.get_weights(weights_phi);
226
227        // Perform the computation, with all weight points
228        double sum = 0.0;
229        double norm = 0.0;
230        double norm_vol = 0.0;
231        double vol = 0.0;
232        double pi = 4.0*atan(1.0);
233        // Loop over radius weight points
234        for(size_t i=0; i<weights_rad.size(); i++) {
235                dp.radius = weights_rad[i].value;
236
237
238                // Loop over length weight points
239                for(size_t j=0; j<weights_len.size(); j++) {
240                        dp.length = weights_len[j].value;
241
242                        // Average over theta distribution
243                        for(size_t k=0; k<weights_theta.size(); k++) {
244                                dp.cyl_theta = weights_theta[k].value;
245
246                                // Average over phi distribution
247                                for(size_t l=0; l<weights_phi.size(); l++) {
248                                        dp.cyl_phi = weights_phi[l].value;
249                                        //Un-normalize by volume
250                                        double _ptvalue = weights_rad[i].weight
251                                                * weights_len[j].weight
252                                                * weights_theta[k].weight
253                                                * weights_phi[l].weight
254                                                * cylinder_analytical_2DXY(&dp, qx, qy)
255                                                *pow(weights_rad[i].value,2)*weights_len[j].value;
256                                        if (weights_theta.size()>1) {
257                                                _ptvalue *= fabs(sin(weights_theta[k].value*pi/180.0));
258                                        }
259                                        sum += _ptvalue;
260                                        //Find average volume
261                                        vol += weights_rad[i].weight
262                                                        * weights_len[j].weight
263                                                        * pow(weights_rad[i].value,2)*weights_len[j].value;
264                                        //Find norm for volume
265                                        norm_vol += weights_rad[i].weight
266                                                        * weights_len[j].weight;
267
268                                        norm += weights_rad[i].weight
269                                                * weights_len[j].weight
270                                                * weights_theta[k].weight
271                                                * weights_phi[l].weight;
272
273                                }
274                        }
275                }
276        }
277        // Averaging in theta needs an extra normalization
278        // factor to account for the sin(theta) term in the
279        // integration (see documentation).
280        if (weights_theta.size()>1) norm = norm / asin(1.0);
281        if (vol != 0.0 && norm_vol != 0.0) {
282                //Re-normalize by avg volume
283                sum = sum/(vol/norm_vol);}
284
285        return sum/norm + background();
286}
287
288/**
289 * Function to evaluate 2D scattering function
290 * @param pars: parameters of the cylinder
291 * @param q: q-value
292 * @param phi: angle phi
293 * @return: function value
294 */
295double CylinderModel :: evaluate_rphi(double q, double phi) {
296        double qx = q*cos(phi);
297        double qy = q*sin(phi);
298        return (*this).operator()(qx, qy);
299}
300/**
301 * Function to calculate effective radius
302 * @return: effective radius value
303 */
304double CylinderModel :: calculate_ER() {
305        CylinderParameters dp;
306
307        dp.radius     = radius();
308        dp.length     = length();
309        double rad_out = 0.0;
310
311        // Perform the computation, with all weight points
312        double sum = 0.0;
313        double norm = 0.0;
314
315        // Get the dispersion points for the major shell
316        vector<WeightPoint> weights_length;
317        length.get_weights(weights_length);
318
319        // Get the dispersion points for the minor shell
320        vector<WeightPoint> weights_radius ;
321        radius.get_weights(weights_radius);
322
323        // Loop over major shell weight points
324        for(int i=0; i< (int)weights_length.size(); i++) {
325                dp.length = weights_length[i].value;
326                for(int k=0; k< (int)weights_radius.size(); k++) {
327                        dp.radius = weights_radius[k].value;
328                        //Note: output of "DiamCyl(dp.length,dp.radius)" is DIAMETER.
329                        sum +=weights_length[i].weight
330                                * weights_radius[k].weight*DiamCyl(dp.length,dp.radius)/2.0;
331                        norm += weights_length[i].weight* weights_radius[k].weight;
332                }
333        }
334        if (norm != 0){
335                //return the averaged value
336                rad_out =  sum/norm;}
337        else{
338                //return normal value
339                //Note: output of "DiamCyl(dp.length,dp.radius)" is DIAMETER.
340                rad_out = DiamCyl(dp.length,dp.radius)/2.0;}
341
342        return rad_out;
343}
344double CylinderModel :: calculate_VR() {
345  return 1.0;
346}
Note: See TracBrowser for help on using the repository browser.