source: sasview/src/sans/models/c_extension/c_models/triaxialellipsoid.cpp @ 230f479

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 230f479 was 230f479, checked in by Mathieu Doucet <doucetm@…>, 11 years ago

Rename C source dir for models (minor updates)

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