source: sasview/sansmodels/src/c_models/cylinder.cpp @ 0c2389e

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 0c2389e was 0c2389e, checked in by Mathieu Doucet <doucetm@…>, 12 years ago

refactored capped cylinder

  • Property mode set to 100644
File size: 9.3 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;
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)>1\n");
148      return 0;
149    }
150
151    // Note: cos(alpha) = 0 and 1 will get an
152    // undefined value from CylKernel
153  alpha = acos( cos_val );
154
155  // Call the IGOR library function to get the kernel
156  answer = CylKernel(q, pars->radius, pars->length/2.0, alpha) / sin(alpha);
157
158  // Multiply by contrast^2
159  answer *= (pars->sldCyl - pars->sldSolv)*(pars->sldCyl - pars->sldSolv);
160
161  //normalize by cylinder volume
162  //NOTE that for this (Fournet) definition of the integral, one must MULTIPLY by Vcyl
163    vol = acos(-1.0) * pars->radius * pars->radius * pars->length;
164  answer *= vol;
165
166  //convert to [cm-1]
167  answer *= 1.0e8;
168
169  //Scale
170  answer *= pars->scale;
171
172  // add in the background
173  answer += pars->background;
174
175  return answer;
176}
177
178/**
179 * Function to evaluate 2D scattering function
180 * @param pars: parameters of the cylinder
181 * @param q: q-value
182 * @return: function value
183 */
184static double cylinder_analytical_2DXY(CylinderParameters *pars, double qx, double qy) {
185  double q;
186  q = sqrt(qx*qx+qy*qy);
187  return cylinder_analytical_2D_scaled(pars, q, qx/q, qy/q);
188}
189
190/**
191 * Function to evaluate 2D scattering function
192 * @param q_x: value of Q along x
193 * @param q_y: value of Q along y
194 * @return: function value
195 */
196double CylinderModel :: operator()(double qx, double qy) {
197        CylinderParameters dp;
198        // Fill parameter array
199        dp.scale      = scale();
200        dp.radius     = radius();
201        dp.length     = length();
202        dp.sldCyl   = sldCyl();
203        dp.sldSolv   = sldSolv();
204        dp.background = 0.0;
205        dp.cyl_theta  = cyl_theta();
206        dp.cyl_phi    = cyl_phi();
207
208        // Get the dispersion points for the radius
209        vector<WeightPoint> weights_rad;
210        radius.get_weights(weights_rad);
211
212        // Get the dispersion points for the length
213        vector<WeightPoint> weights_len;
214        length.get_weights(weights_len);
215
216        // Get angular averaging for theta
217        vector<WeightPoint> weights_theta;
218        cyl_theta.get_weights(weights_theta);
219
220        // Get angular averaging for phi
221        vector<WeightPoint> weights_phi;
222        cyl_phi.get_weights(weights_phi);
223
224        // Perform the computation, with all weight points
225        double sum = 0.0;
226        double norm = 0.0;
227        double norm_vol = 0.0;
228        double vol = 0.0;
229        double pi = 4.0*atan(1.0);
230        // Loop over radius weight points
231        for(size_t i=0; i<weights_rad.size(); i++) {
232                dp.radius = weights_rad[i].value;
233
234
235                // Loop over length weight points
236                for(size_t j=0; j<weights_len.size(); j++) {
237                        dp.length = weights_len[j].value;
238
239                        // Average over theta distribution
240                        for(size_t k=0; k<weights_theta.size(); k++) {
241                                dp.cyl_theta = weights_theta[k].value;
242
243                                // Average over phi distribution
244                                for(size_t l=0; l<weights_phi.size(); l++) {
245                                        dp.cyl_phi = weights_phi[l].value;
246                                        //Un-normalize by volume
247                                        double _ptvalue = weights_rad[i].weight
248                                                * weights_len[j].weight
249                                                * weights_theta[k].weight
250                                                * weights_phi[l].weight
251                                                * cylinder_analytical_2DXY(&dp, qx, qy)
252                                                *pow(weights_rad[i].value,2)*weights_len[j].value;
253                                        if (weights_theta.size()>1) {
254                                                _ptvalue *= fabs(sin(weights_theta[k].value*pi/180.0));
255                                        }
256                                        sum += _ptvalue;
257                                        //Find average volume
258                                        vol += weights_rad[i].weight
259                                                        * weights_len[j].weight
260                                                        * pow(weights_rad[i].value,2)*weights_len[j].value;
261                                        //Find norm for volume
262                                        norm_vol += weights_rad[i].weight
263                                                        * weights_len[j].weight;
264
265                                        norm += weights_rad[i].weight
266                                                * weights_len[j].weight
267                                                * weights_theta[k].weight
268                                                * weights_phi[l].weight;
269
270                                }
271                        }
272                }
273        }
274        // Averaging in theta needs an extra normalization
275        // factor to account for the sin(theta) term in the
276        // integration (see documentation).
277        if (weights_theta.size()>1) norm = norm / asin(1.0);
278        if (vol != 0.0 && norm_vol != 0.0) {
279                //Re-normalize by avg volume
280                sum = sum/(vol/norm_vol);}
281
282        return sum/norm + background();
283}
284
285/**
286 * Function to evaluate 2D scattering function
287 * @param pars: parameters of the cylinder
288 * @param q: q-value
289 * @param phi: angle phi
290 * @return: function value
291 */
292double CylinderModel :: evaluate_rphi(double q, double phi) {
293        double qx = q*cos(phi);
294        double qy = q*sin(phi);
295        return (*this).operator()(qx, qy);
296}
297/**
298 * Function to calculate effective radius
299 * @return: effective radius value
300 */
301double CylinderModel :: calculate_ER() {
302        CylinderParameters dp;
303
304        dp.radius     = radius();
305        dp.length     = length();
306        double rad_out = 0.0;
307
308        // Perform the computation, with all weight points
309        double sum = 0.0;
310        double norm = 0.0;
311
312        // Get the dispersion points for the major shell
313        vector<WeightPoint> weights_length;
314        length.get_weights(weights_length);
315
316        // Get the dispersion points for the minor shell
317        vector<WeightPoint> weights_radius ;
318        radius.get_weights(weights_radius);
319
320        // Loop over major shell weight points
321        for(int i=0; i< (int)weights_length.size(); i++) {
322                dp.length = weights_length[i].value;
323                for(int k=0; k< (int)weights_radius.size(); k++) {
324                        dp.radius = weights_radius[k].value;
325                        //Note: output of "DiamCyl(dp.length,dp.radius)" is DIAMETER.
326                        sum +=weights_length[i].weight
327                                * weights_radius[k].weight*DiamCyl(dp.length,dp.radius)/2.0;
328                        norm += weights_length[i].weight* weights_radius[k].weight;
329                }
330        }
331        if (norm != 0){
332                //return the averaged value
333                rad_out =  sum/norm;}
334        else{
335                //return normal value
336                //Note: output of "DiamCyl(dp.length,dp.radius)" is DIAMETER.
337                rad_out = DiamCyl(dp.length,dp.radius)/2.0;}
338
339        return rad_out;
340}
Note: See TracBrowser for help on using the repository browser.