source: sasview/sansmodels/src/c_models/ellipticalcylinder.cpp @ cd2f75d

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 cd2f75d was 82c11d3, checked in by Mathieu Doucet <doucetm@…>, 13 years ago

refactored bunch of models

  • Property mode set to 100644
File size: 13.2 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>
25#include <stdlib.h>
26using namespace std;
27#include "elliptical_cylinder.h"
28
29extern "C" {
30#include "libCylinder.h"
31#include "libStructureFactor.h"
32}
33
34typedef struct {
35  double scale;
36  double r_minor;
37  double r_ratio;
38  double length;
39  double sldCyl;
40  double sldSolv;
41  double background;
42  double cyl_theta;
43  double cyl_phi;
44  double cyl_psi;
45} EllipticalCylinderParameters;
46
47
48static double elliptical_cylinder_kernel(EllipticalCylinderParameters *pars, double q, double alpha, double nu) {
49  double qr;
50  double qL;
51  double Be,Si;
52  double r_major;
53  double kernel;
54
55  r_major = pars->r_ratio * pars->r_minor;
56
57  qr = q*sin(alpha)*sqrt( r_major*r_major*sin(nu)*sin(nu) + pars->r_minor*pars->r_minor*cos(nu)*cos(nu) );
58  qL = q*pars->length*cos(alpha)/2.0;
59
60  if (qr==0){
61    Be = 0.5;
62  }else{
63    Be = NR_BessJ1(qr)/qr;
64  }
65  if (qL==0){
66    Si = 1.0;
67  }else{
68    Si = sin(qL)/qL;
69  }
70
71
72  kernel = 2.0*Be * Si;
73  return kernel*kernel;
74}
75
76/**
77 * Function to evaluate 2D scattering function
78 * @param pars: parameters of the cylinder
79 * @param q: q-value
80 * @param q_x: q_x / q
81 * @param q_y: q_y / q
82 * @return: function value
83 */
84static double elliptical_cylinder_analytical_2D_scaled(EllipticalCylinderParameters *pars, double q, double q_x, double q_y) {
85  double cyl_x, cyl_y, cyl_z;
86  double ell_x, ell_y;
87  double q_z;
88  double alpha, vol, cos_val;
89  double nu, cos_nu;
90  double answer;
91  //convert angle degree to radian
92  double pi = 4.0*atan(1.0);
93  double theta = pars->cyl_theta * pi/180.0;
94  double phi = pars->cyl_phi * pi/180.0;
95  double psi = pars->cyl_psi * pi/180.0;
96
97  //Cylinder orientation
98  cyl_x = sin(theta) * cos(phi);
99  cyl_y = sin(theta) * sin(phi);
100  cyl_z = cos(theta);
101
102  // q vector
103  q_z = 0;
104
105  // Compute the angle btw vector q and the
106  // axis of the cylinder
107  cos_val = cyl_x*q_x + cyl_y*q_y + cyl_z*q_z;
108
109  // The following test should always pass
110  if (fabs(cos_val)>1.0) {
111    printf("cyl_ana_2D: Unexpected error: cos(alpha)>1\n");
112    return 0;
113  }
114
115  // Note: cos(alpha) = 0 and 1 will get an
116  // undefined value from CylKernel
117  alpha = acos( cos_val );
118
119  //ellipse orientation:
120  // the elliptical corss section was transformed and projected
121  // into the detector plane already through sin(alpha)and furthermore psi remains as same
122  // on the detector plane.
123  // So, all we need is to calculate the angle (nu) of the minor axis of the ellipse wrt
124  // the wave vector q.
125
126  //x- y- component on the detector plane.
127  ell_x =  cos(psi);
128  ell_y =  sin(psi);
129
130  // calculate the axis of the ellipse wrt q-coord.
131  cos_nu = ell_x*q_x + ell_y*q_y;
132  nu = acos(cos_nu);
133
134  // The following test should always pass
135  if (fabs(cos_nu)>1.0) {
136    printf("cyl_ana_2D: Unexpected error: cos(nu)>1\n");
137    return 0;
138  }
139
140  answer = elliptical_cylinder_kernel(pars, q, alpha,nu);
141
142  // Multiply by contrast^2
143  answer *= (pars->sldCyl - pars->sldSolv) * (pars->sldCyl - pars->sldSolv);
144
145  //normalize by cylinder volume
146  //NOTE that for this (Fournet) definition of the integral, one must MULTIPLY by Vcyl
147  vol = acos(-1.0) * pars->r_minor * pars->r_minor * pars->r_ratio * pars->length;
148  answer *= vol;
149
150  //convert to [cm-1]
151  answer *= 1.0e8;
152
153  //Scale
154  answer *= pars->scale;
155
156  // add in the background
157  answer += pars->background;
158
159  return answer;
160}
161
162
163/**
164 * Function to evaluate 2D scattering function
165 * @param pars: parameters of the cylinder
166 * @param q: q-value
167 * @return: function value
168 */
169static double elliptical_cylinder_analytical_2DXY(EllipticalCylinderParameters *pars, double qx, double qy) {
170  double q;
171  q = sqrt(qx*qx+qy*qy);
172  return elliptical_cylinder_analytical_2D_scaled(pars, q, qx/q, qy/q);
173}
174
175EllipticalCylinderModel :: EllipticalCylinderModel() {
176  scale      = Parameter(1.0);
177  r_minor    = Parameter(20.0, true);
178  r_minor.set_min(0.0);
179  r_ratio    = Parameter(1.5, true);
180  r_ratio.set_min(0.0);
181  length     = Parameter(400.0, true);
182  length.set_min(0.0);
183  sldCyl   = Parameter(4.e-6);
184  sldSolv   = Parameter(1.e-6);
185  background = Parameter(0.0);
186  cyl_theta  = Parameter(57.325, true);
187  cyl_phi    = Parameter(0.0, true);
188  cyl_psi    = Parameter(0.0, true);
189}
190
191/**
192 * Function to evaluate 1D scattering function
193 * The NIST IGOR library is used for the actual calculation.
194 * @param q: q-value
195 * @return: function value
196 */
197double EllipticalCylinderModel :: operator()(double q) {
198  double dp[7];
199
200  dp[0] = scale();
201  dp[1] = r_minor();
202  dp[2] = r_ratio();
203  dp[3] = length();
204  dp[4] = sldCyl();
205  dp[5] = sldSolv();
206  dp[6] = 0.0;
207
208  // Get the dispersion points for the r_minor
209  vector<WeightPoint> weights_rad;
210  r_minor.get_weights(weights_rad);
211
212  // Get the dispersion points for the r_ratio
213  vector<WeightPoint> weights_rat;
214  r_ratio.get_weights(weights_rat);
215
216  // Get the dispersion points for the length
217  vector<WeightPoint> weights_len;
218  length.get_weights(weights_len);
219
220  // Perform the computation, with all weight points
221  double sum = 0.0;
222  double norm = 0.0;
223  double vol = 0.0;
224
225  // Loop over r_minor weight points
226  for(size_t i=0; i<weights_rad.size(); i++) {
227    dp[1] = weights_rad[i].value;
228
229    // Loop over r_ratio weight points
230    for(size_t j=0; j<weights_rat.size(); j++) {
231      dp[2] = weights_rat[j].value;
232
233      // Loop over length weight points
234      for(size_t k=0; k<weights_len.size(); k++) {
235        dp[3] = weights_len[k].value;
236        //Un-normalize  by volume
237        sum += weights_rad[i].weight
238            * weights_len[k].weight
239            * weights_rat[j].weight
240            * EllipCyl20(dp, q)
241        * pow(weights_rad[i].value,2) * weights_rat[j].value
242        * weights_len[k].value;
243        //Find average volume
244        vol += weights_rad[i].weight
245            * weights_len[k].weight
246            * weights_rat[j].weight
247            * pow(weights_rad[i].value,2) * weights_rat[j].value
248            * weights_len[k].value;
249        norm += weights_rad[i].weight
250            * weights_len[k].weight
251            * weights_rat[j].weight;
252      }
253    }
254  }
255
256  if (vol != 0.0 && norm != 0.0) {
257    //Re-normalize by avg volume
258    sum = sum/(vol/norm);}
259
260  return sum/norm + background();
261}
262
263/**
264 * Function to evaluate 2D scattering function
265 * @param q_x: value of Q along x
266 * @param q_y: value of Q along y
267 * @return: function value
268 */
269double EllipticalCylinderModel :: operator()(double qx, double qy) {
270  EllipticalCylinderParameters dp;
271  // Fill parameter array
272  dp.scale      = scale();
273  dp.r_minor    = r_minor();
274  dp.r_ratio    = r_ratio();
275  dp.length     = length();
276  dp.sldCyl   = sldCyl();
277  dp.sldSolv   = sldSolv();
278  dp.background = 0.0;
279  dp.cyl_theta  = cyl_theta();
280  dp.cyl_phi    = cyl_phi();
281  dp.cyl_psi    = cyl_psi();
282
283  // Get the dispersion points for the r_minor
284  vector<WeightPoint> weights_rad;
285  r_minor.get_weights(weights_rad);
286
287  // Get the dispersion points for the r_ratio
288  vector<WeightPoint> weights_rat;
289  r_ratio.get_weights(weights_rat);
290
291  // Get the dispersion points for the length
292  vector<WeightPoint> weights_len;
293  length.get_weights(weights_len);
294
295  // Get angular averaging for theta
296  vector<WeightPoint> weights_theta;
297  cyl_theta.get_weights(weights_theta);
298
299  // Get angular averaging for phi
300  vector<WeightPoint> weights_phi;
301  cyl_phi.get_weights(weights_phi);
302
303  // Get angular averaging for psi
304  vector<WeightPoint> weights_psi;
305  cyl_psi.get_weights(weights_psi);
306
307  // Perform the computation, with all weight points
308  double sum = 0.0;
309  double norm = 0.0;
310  double norm_vol = 0.0;
311  double vol = 0.0;
312  double pi = 4.0*atan(1.0);
313  // Loop over minor radius weight points
314  for(size_t i=0; i<weights_rad.size(); i++) {
315    dp.r_minor = weights_rad[i].value;
316
317
318    // Loop over length weight points
319    for(size_t j=0; j<weights_len.size(); j++) {
320      dp.length = weights_len[j].value;
321
322      // Loop over r_ration weight points
323      for(size_t m=0; m<weights_rat.size(); m++) {
324        dp.r_ratio = weights_rat[m].value;
325
326        // Average over theta distribution
327        for(size_t k=0; k<weights_theta.size(); k++) {
328          dp.cyl_theta = weights_theta[k].value;
329
330          // Average over phi distribution
331          for(size_t l=0; l<weights_phi.size(); l++) {
332            dp.cyl_phi = weights_phi[l].value;
333
334            // Average over phi distribution
335            for(size_t o=0; o<weights_psi.size(); o++) {
336              dp.cyl_psi = weights_psi[o].value;
337              //Un-normalize by volume
338              double _ptvalue = weights_rad[i].weight
339                  * weights_len[j].weight
340                  * weights_rat[m].weight
341                  * weights_theta[k].weight
342                  * weights_phi[l].weight
343                  * weights_psi[o].weight
344                  * elliptical_cylinder_analytical_2DXY(&dp, qx, qy)
345              * pow(weights_rad[i].value,2)
346              * weights_len[j].value
347              * weights_rat[m].value;
348              if (weights_theta.size()>1) {
349                _ptvalue *= fabs(sin(weights_theta[k].value*pi/180.0));
350              }
351              sum += _ptvalue;
352              //Find average volume
353              vol += weights_rad[i].weight
354                  * weights_len[j].weight
355                  * weights_rat[m].weight
356                  * pow(weights_rad[i].value,2)
357              * weights_len[j].value
358              * weights_rat[m].value;
359              //Find norm for volume
360              norm_vol += weights_rad[i].weight
361                  * weights_len[j].weight
362                  * weights_rat[m].weight;
363
364              norm += weights_rad[i].weight
365                  * weights_len[j].weight
366                  * weights_rat[m].weight
367                  * weights_theta[k].weight
368                  * weights_phi[l].weight
369                  * weights_psi[o].weight;
370
371            }
372          }
373        }
374      }
375    }
376  }
377  // Averaging in theta needs an extra normalization
378  // factor to account for the sin(theta) term in the
379  // integration (see documentation).
380  if (weights_theta.size()>1) norm = norm / asin(1.0);
381
382  if (vol != 0.0 && norm_vol != 0.0) {
383    //Re-normalize by avg volume
384    sum = sum/(vol/norm_vol);}
385
386  return sum/norm + background();
387
388}
389
390/**
391 * Function to evaluate 2D scattering function
392 * @param pars: parameters of the cylinder
393 * @param q: q-value
394 * @param phi: angle phi
395 * @return: function value
396 */
397double EllipticalCylinderModel :: evaluate_rphi(double q, double phi) {
398  double qx = q*cos(phi);
399  double qy = q*sin(phi);
400  return (*this).operator()(qx, qy);
401}
402/**
403 * Function to calculate effective radius
404 * @return: effective radius value
405 */
406double EllipticalCylinderModel :: calculate_ER() {
407  EllipticalCylinderParameters dp;
408  dp.r_minor    = r_minor();
409  dp.r_ratio    = r_ratio();
410  dp.length     = length();
411  double rad_out = 0.0;
412  double suf_rad = sqrt(dp.r_minor*dp.r_minor*dp.r_ratio);
413
414  // Perform the computation, with all weight points
415  double sum = 0.0;
416  double norm = 0.0;
417
418  // Get the dispersion points for the r_minor
419  vector<WeightPoint> weights_rad;
420  r_minor.get_weights(weights_rad);
421
422  // Get the dispersion points for the r_ratio
423  vector<WeightPoint> weights_rat;
424  r_ratio.get_weights(weights_rat);
425
426  // Get the dispersion points for the length
427  vector<WeightPoint> weights_len;
428  length.get_weights(weights_len);
429
430  // Loop over minor radius weight points
431  for(size_t i=0; i<weights_rad.size(); i++) {
432    dp.r_minor = weights_rad[i].value;
433
434    // Loop over length weight points
435    for(size_t j=0; j<weights_len.size(); j++) {
436      dp.length = weights_len[j].value;
437
438      // Loop over r_ration weight points
439      for(size_t m=0; m<weights_rat.size(); m++) {
440        dp.r_ratio = weights_rat[m].value;
441        //Calculate surface averaged radius
442        suf_rad = sqrt(dp.r_minor * dp.r_minor * dp.r_ratio);
443
444        //Note: output of "DiamCyl(dp.length,dp.radius)" is DIAMETER.
445        sum +=weights_rad[i].weight *weights_len[j].weight
446            * weights_rat[m].weight*DiamCyl(dp.length, suf_rad)/2.0;
447        norm += weights_rad[i].weight *weights_len[j].weight* weights_rat[m].weight;
448      }
449    }
450  }
451  if (norm != 0){
452    //return the averaged value
453    rad_out =  sum/norm;}
454  else{
455    //return normal value
456    //Note: output of "DiamCyl(dp.length,dp.radius)" is DIAMETER.
457    rad_out = DiamCyl(dp.length, suf_rad)/2.0;}
458
459  return rad_out;
460}
Note: See TracBrowser for help on using the repository browser.