source: sasview/src/sas/models/c_extension/c_models/ellipticalcylinder.cpp @ 79492222

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 79492222 was 79492222, checked in by krzywon, 9 years ago

Changed the file and folder names to remove all SANS references.

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