source: sasview/src/sas/models/c_extension/c_models/csparallelepiped.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: 16.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 2010, 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#include <math.h>
22#include "parameters.hh"
23#include <stdio.h>
24using namespace std;
25
26extern "C" {
27#include "libCylinder.h"
28#include "libStructureFactor.h"
29}
30#include "csparallelepiped.h"
31
32// Convenience parameter structure
33typedef struct {
34  double scale;
35  double shortA;
36  double midB;
37  double longC;
38  double rimA;
39  double rimB;
40  double rimC;
41  double sld_rimA;
42  double sld_rimB;
43  double sld_rimC;
44  double sld_pcore;
45  double sld_solv;
46  double background;
47  double parallel_theta;
48  double parallel_phi;
49  double parallel_psi;
50} CSParallelepipedParameters;
51
52static double cspkernel(double dp[],double q, double ala, double alb, double alc){
53  // mu passed in is really mu*sqrt(1-sig^2)
54  double argA,argB,argC,argtA,argtB,argtC,tmp1,tmp2,tmp3,tmpt1,tmpt2,tmpt3;     //local variables
55
56  double aa,bb,cc, ta,tb,tc;
57  double Vin,Vot,V1,V2,V3;
58  double rhoA,rhoB,rhoC, rhoP, rhosolv;
59  double dr0, drA,drB, drC;
60  double retVal;
61
62  aa = dp[1];
63  bb = dp[2];
64  cc = dp[3];
65  ta = dp[4];
66  tb = dp[5];
67  tc = dp[6];
68  rhoA=dp[7];
69  rhoB=dp[8];
70  rhoC=dp[9];
71  rhoP=dp[10];
72  rhosolv=dp[11];
73  dr0=rhoP-rhosolv;
74  drA=rhoA-rhosolv;
75  drB=rhoB-rhosolv;
76  drC=rhoC-rhosolv;
77  Vin=(aa*bb*cc);
78  Vot=(aa*bb*cc+2.0*ta*bb*cc+2.0*aa*tb*cc+2.0*aa*bb*tc);
79  V1=(2.0*ta*bb*cc);   //  incorrect V1 (aa*bb*cc+2*ta*bb*cc)
80  V2=(2.0*aa*tb*cc);  // incorrect V2(aa*bb*cc+2*aa*tb*cc)
81  V3=(2.0*aa*bb*tc);
82  //aa = aa/bb;
83  ta=(aa+2.0*ta);///bb;
84  tb=(aa+2.0*tb);///bb;
85  tc=(aa+2.0*tc);
86  //handle arg=0 separately, as sin(t)/t -> 1 as t->0
87  argA = q*aa*ala/2.0;
88  argB = q*bb*alb/2.0;
89  argC = q*cc*alc/2.0;
90  argtA = q*ta*ala/2.0;
91  argtB = q*tb*alb/2.0;
92  argtC = q*tc*alc/2.0;
93
94  if(argA==0.0) {
95    tmp1 = 1.0;
96  } else {
97    tmp1 = sin(argA)/argA;
98  }
99  if (argB==0.0) {
100    tmp2 = 1.0;
101  } else {
102    tmp2 = sin(argB)/argB;
103  }
104
105  if (argC==0.0) {
106    tmp3 = 1.0;
107  } else {
108    tmp3 = sin(argC)/argC;
109  }
110  if(argtA==0.0) {
111    tmpt1 = 1.0;
112  } else {
113    tmpt1 = sin(argtA)/argtA;
114  }
115  if (argtB==0.0) {
116    tmpt2 = 1.0;
117  } else {
118    tmpt2 = sin(argtB)/argtB;
119  }
120  if (argtC==0.0) {
121    tmpt3 = 1.0;
122  } else {
123    tmpt3 = sin(argtC)*sin(argtC)/argtC/argtC;
124  }
125  // This expression is different from NIST/IGOR package (I strongly believe the IGOR is wrong!!!). 10/15/2010.
126  retVal =( dr0*tmp1*tmp2*tmp3*Vin + drA*(tmpt1-tmp1)*tmp2*tmp3*V1+ drB*tmp1*(tmpt2-tmp2)*tmp3*V2 + drC*tmp1*tmp2*(tmpt3-tmp3)*V3)*
127      ( dr0*tmp1*tmp2*tmp3*Vin + drA*(tmpt1-tmp1)*tmp2*tmp3*V1+ drB*tmp1*(tmpt2-tmp2)*tmp3*V2 + drC*tmp1*tmp2*(tmpt3-tmp3)*V3);   //  correct FF : square of sum of phase factors
128  //retVal *= (tmp3*tmp3);
129  retVal /= Vot;
130
131  return (retVal);
132
133}//Function cspkernel()
134
135/**
136 * Function to evaluate 2D scattering function
137 * @param pars: parameters of the CSparallelepiped
138 * @param q: q-value
139 * @param q_x: q_x / q
140 * @param q_y: q_y / q
141 * @return: function value
142 */
143static double csparallelepiped_analytical_2D_scaled(CSParallelepipedParameters *pars, double q, double q_x, double q_y) {
144  double dp[13];
145  double cparallel_x, cparallel_y, bparallel_x, bparallel_y, parallel_x, parallel_y;
146  //double q_z;
147  double cos_val_c, cos_val_b, cos_val_a, edgeA, edgeB, edgeC;
148
149  double answer;
150  //convert angle degree to radian
151  double pi = 4.0*atan(1.0);
152  double theta = pars->parallel_theta * pi/180.0;
153  double phi = pars->parallel_phi * pi/180.0;
154  double psi = pars->parallel_psi* pi/180.0;
155
156  // Fill paramater array
157  dp[0] = 1.0;
158  dp[1] = pars->shortA;
159  dp[2] = pars->midB;
160  dp[3] = pars->longC;
161  dp[4] = pars->rimA;
162  dp[5] = pars->rimB;
163  dp[6] = pars->rimC;
164  dp[7] = pars->sld_rimA;
165  dp[8] = pars->sld_rimB;
166  dp[9] = pars->sld_rimC;
167  dp[10] = pars->sld_pcore;
168  dp[11] = pars->sld_solv;
169  dp[12] = 0.0;
170
171
172  edgeA = pars->shortA;
173  edgeB = pars->midB;
174  edgeC = pars->longC;
175
176
177  // parallelepiped c axis orientation
178  cparallel_x = cos(theta) * cos(phi);
179  cparallel_y = sin(theta);
180  //cparallel_z = -cos(theta) * sin(phi);
181
182  // q vector
183  //q_z = 0.0;
184
185  // Compute the angle btw vector q and the
186  // axis of the parallelepiped
187  cos_val_c = cparallel_x*q_x + cparallel_y*q_y;// + cparallel_z*q_z;
188  //alpha = acos(cos_val_c);
189
190  // parallelepiped a axis orientation
191  parallel_x = -cos(phi)*sin(psi) * sin(theta)+sin(phi)*cos(psi);
192  parallel_y = sin(psi)*cos(theta);
193
194  cos_val_a = parallel_x*q_x + parallel_y*q_y;
195
196
197
198  // parallelepiped b axis orientation
199  bparallel_x = -sin(theta)*cos(psi)*cos(phi)-sin(psi)*sin(phi);
200  bparallel_y = cos(theta)*cos(psi);
201  // axis of the parallelepiped
202  cos_val_b = bparallel_x*q_x + bparallel_y*q_y;
203
204  // The following test should always pass
205  if (fabs(cos_val_c)>1.0) {
206    //printf("parallel_ana_2D: Unexpected error: cos(alpha)>1\n");
207    cos_val_c = 1.0;
208  }
209  if (fabs(cos_val_a)>1.0) {
210    //printf("parallel_ana_2D: Unexpected error: cos(alpha)>1\n");
211    cos_val_a = 1.0;
212  }
213  if (fabs(cos_val_b)>1.0) {
214    //printf("parallel_ana_2D: Unexpected error: cos(alpha)>1\n");
215    cos_val_b = 1.0;
216  }
217  // Call the IGOR library function to get the kernel
218  answer = cspkernel( dp, q, cos_val_a, cos_val_b, cos_val_c);
219
220  //convert to [cm-1]
221  answer *= 1.0e8;
222
223  //Scale
224  answer *= pars->scale;
225
226  // add in the background
227  answer += pars->background;
228
229  return answer;
230}
231
232/**
233 * Function to evaluate 2D scattering function
234 * @param pars: parameters of the CSparallelepiped
235 * @param q: q-value
236 * @return: function value
237 */
238static double csparallelepiped_analytical_2DXY(CSParallelepipedParameters *pars, double qx, double qy) {
239  double q;
240  q = sqrt(qx*qx+qy*qy);
241  return csparallelepiped_analytical_2D_scaled(pars, q, qx/q, qy/q);
242}
243
244
245
246
247CSParallelepipedModel :: CSParallelepipedModel() {
248  scale      = Parameter(1.0);
249  shortA     = Parameter(35.0, true);
250  shortA.set_min(1.0);
251  midB     = Parameter(75.0, true);
252  midB.set_min(1.0);
253  longC    = Parameter(400.0, true);
254  longC.set_min(1.0);
255  rimA     = Parameter(10.0, true);
256  rimB     = Parameter(10.0, true);
257  rimC     = Parameter(10.0, true);
258  sld_rimA     = Parameter(2.0e-6, true);
259  sld_rimB     = Parameter(4.0e-6, true);
260  sld_rimC    = Parameter(2.0e-6, true);
261  sld_pcore   = Parameter(1.0e-6);
262  sld_solv   = Parameter(6.0e-6);
263  background = Parameter(0.06);
264  parallel_theta  = Parameter(0.0, true);
265  parallel_phi    = Parameter(0.0, true);
266  parallel_psi    = Parameter(0.0, true);
267}
268
269/**
270 * Function to evaluate 1D scattering function
271 * The NIST IGOR library is used for the actual calculation.
272 * @param q: q-value
273 * @return: function value
274 */
275double CSParallelepipedModel :: operator()(double q) {
276  double dp[13];
277
278  // Fill parameter array for IGOR library
279  // Add the background after averaging
280  dp[0] = scale();
281  dp[1] = shortA();
282  dp[2] = midB();
283  dp[3] = longC();
284  dp[4] = rimA();
285  dp[5] = rimB();
286  dp[6] = rimC();
287  dp[7] = sld_rimA();
288  dp[8] = sld_rimB();
289  dp[9] = sld_rimC();
290  dp[10] = sld_pcore();
291  dp[11] = sld_solv();
292  dp[12] = 0.0;
293
294  // Get the dispersion points for the short_edgeA
295  vector<WeightPoint> weights_shortA;
296  shortA.get_weights(weights_shortA);
297
298  // Get the dispersion points for the longer_edgeB
299  vector<WeightPoint> weights_midB;
300  midB.get_weights(weights_midB);
301
302  // Get the dispersion points for the longuest_edgeC
303  vector<WeightPoint> weights_longC;
304  longC.get_weights(weights_longC);
305
306
307
308  // Perform the computation, with all weight points
309  double sum = 0.0;
310  double norm = 0.0;
311  double vol = 0.0;
312
313  // Loop over short_edgeA weight points
314  for(int i=0; i< (int)weights_shortA.size(); i++) {
315    dp[1] = weights_shortA[i].value;
316
317    // Loop over longer_edgeB weight points
318    for(int j=0; j< (int)weights_midB.size(); j++) {
319      dp[2] = weights_midB[j].value;
320
321      // Loop over longuest_edgeC weight points
322      for(int k=0; k< (int)weights_longC.size(); k++) {
323        dp[3] = weights_longC[k].value;
324        //Un-normalize  by volume
325        sum += weights_shortA[i].weight * weights_midB[j].weight
326            * weights_longC[k].weight * CSParallelepiped(dp, q)
327        * weights_shortA[i].value*weights_midB[j].value
328        * weights_longC[k].value;
329        //Find average volume
330        vol += weights_shortA[i].weight * weights_midB[j].weight
331            * weights_longC[k].weight
332            * weights_shortA[i].value * weights_midB[j].value
333            * weights_longC[k].value;
334
335        norm += weights_shortA[i].weight
336            * weights_midB[j].weight * weights_longC[k].weight;
337      }
338    }
339  }
340  if (vol != 0.0 && norm != 0.0) {
341    //Re-normalize by avg volume
342    sum = sum/(vol/norm);}
343
344  return sum/norm + background();
345}
346/**
347 * Function to evaluate 2D scattering function
348 * @param q_x: value of Q along x
349 * @param q_y: value of Q along y
350 * @return: function value
351 */
352double CSParallelepipedModel :: operator()(double qx, double qy) {
353  CSParallelepipedParameters dp;
354  // Fill parameter array
355  dp.scale      = scale();
356  dp.shortA   = shortA();
357  dp.midB   = midB();
358  dp.longC  = longC();
359  dp.rimA   = rimA();
360  dp.rimB   = rimB();
361  dp.rimC  = rimC();
362  dp.sld_rimA   = sld_rimA();
363  dp.sld_rimB   = sld_rimB();
364  dp.sld_rimC  = sld_rimC();
365  dp.sld_pcore   = sld_pcore();
366  dp.sld_solv   = sld_solv();
367  dp.background = 0.0;
368  //dp.background = background();
369  dp.parallel_theta  = parallel_theta();
370  dp.parallel_phi    = parallel_phi();
371  dp.parallel_psi    = parallel_psi();
372
373
374
375  // Get the dispersion points for the short_edgeA
376  vector<WeightPoint> weights_shortA;
377  shortA.get_weights(weights_shortA);
378
379  // Get the dispersion points for the longer_edgeB
380  vector<WeightPoint> weights_midB;
381  midB.get_weights(weights_midB);
382
383  // Get the dispersion points for the longuest_edgeC
384  vector<WeightPoint> weights_longC;
385  longC.get_weights(weights_longC);
386
387  // Get angular averaging for theta
388  vector<WeightPoint> weights_parallel_theta;
389  parallel_theta.get_weights(weights_parallel_theta);
390
391  // Get angular averaging for phi
392  vector<WeightPoint> weights_parallel_phi;
393  parallel_phi.get_weights(weights_parallel_phi);
394
395  // Get angular averaging for psi
396  vector<WeightPoint> weights_parallel_psi;
397  parallel_psi.get_weights(weights_parallel_psi);
398
399  // Perform the computation, with all weight points
400  double sum = 0.0;
401  double norm = 0.0;
402  double norm_vol = 0.0;
403  double vol = 0.0;
404  double pi = 4.0*atan(1.0);
405
406  // Loop over radius weight points
407  for(int i=0; i< (int)weights_shortA.size(); i++) {
408    dp.shortA = weights_shortA[i].value;
409
410    // Loop over longer_edgeB weight points
411    for(int j=0; j< (int)weights_midB.size(); j++) {
412      dp.midB = weights_midB[j].value;
413
414      // Average over longuest_edgeC distribution
415      for(int k=0; k< (int)weights_longC.size(); k++) {
416        dp.longC = weights_longC[k].value;
417
418        // Average over theta distribution
419        for(int l=0; l< (int)weights_parallel_theta.size(); l++) {
420          dp.parallel_theta = weights_parallel_theta[l].value;
421
422          // Average over phi distribution
423          for(int m=0; m< (int)weights_parallel_phi.size(); m++) {
424            dp.parallel_phi = weights_parallel_phi[m].value;
425
426            // Average over phi distribution
427            for(int n=0; n< (int)weights_parallel_psi.size(); n++) {
428              dp.parallel_psi = weights_parallel_psi[n].value;
429              //Un-normalize by volume
430              double _ptvalue = weights_shortA[i].weight
431                  * weights_midB[j].weight
432                  * weights_longC[k].weight
433                  * weights_parallel_theta[l].weight
434                  * weights_parallel_phi[m].weight
435                  * weights_parallel_psi[n].weight
436                  * csparallelepiped_analytical_2DXY(&dp, qx, qy)
437              * weights_shortA[i].value*weights_midB[j].value
438              * weights_longC[k].value;
439
440              if (weights_parallel_theta.size()>1) {
441                _ptvalue *= fabs(cos(weights_parallel_theta[l].value*pi/180.0));
442              }
443              sum += _ptvalue;
444              //Find average volume
445              vol += weights_shortA[i].weight
446                  * weights_midB[j].weight
447                  * weights_longC[k].weight
448                  * weights_shortA[i].value*weights_midB[j].value
449                  * weights_longC[k].value;
450              //Find norm for volume
451              norm_vol += weights_shortA[i].weight
452                  * weights_midB[j].weight
453                  * weights_longC[k].weight;
454
455              norm += weights_shortA[i].weight
456                  * weights_midB[j].weight
457                  * weights_longC[k].weight
458                  * weights_parallel_theta[l].weight
459                  * weights_parallel_phi[m].weight
460                  * weights_parallel_psi[n].weight;
461            }
462          }
463
464        }
465      }
466    }
467  }
468  // Averaging in theta needs an extra normalization
469  // factor to account for the sin(theta) term in the
470  // integration (see documentation).
471  if (weights_parallel_theta.size()>1) norm = norm / asin(1.0);
472
473  if (vol != 0.0 && norm_vol != 0.0) {
474    //Re-normalize by avg volume
475    sum = sum/(vol/norm_vol);}
476
477  return sum/norm + background();
478}
479
480
481/**
482 * Function to evaluate 2D scattering function
483 * @param pars: parameters of the cylinder
484 * @param q: q-value
485 * @param phi: angle phi
486 * @return: function value
487 */
488double CSParallelepipedModel :: evaluate_rphi(double q, double phi) {
489  double qx = q*cos(phi);
490  double qy = q*sin(phi);
491  return (*this).operator()(qx, qy);
492}
493/**
494 * Function to calculate effective radius
495 * @return: effective radius value
496 */
497double CSParallelepipedModel :: calculate_ER() {
498  CSParallelepipedParameters dp;
499  dp.shortA   = shortA();
500  dp.midB   = midB();
501  dp.longC  = longC();
502  dp.rimA   = rimA();
503  dp.rimB   = rimB();
504  dp.rimC  = rimC();
505
506  double rad_out = 0.0;
507  double pi = 4.0*atan(1.0);
508  double suf_rad = sqrt((dp.shortA*dp.midB+2.0*dp.rimA*dp.midB+2.0*dp.rimA*dp.shortA)/pi);
509  double height =(dp.longC + 2.0*dp.rimC);
510  // Perform the computation, with all weight points
511  double sum = 0.0;
512  double norm = 0.0;
513
514  // Get the dispersion points for the short_edgeA
515  vector<WeightPoint> weights_shortA;
516  shortA.get_weights(weights_shortA);
517
518  // Get the dispersion points for the longer_edgeB
519  vector<WeightPoint> weights_midB;
520  midB.get_weights(weights_midB);
521
522  // Get angular averaging for the longuest_edgeC
523  vector<WeightPoint> weights_longC;
524  longC.get_weights(weights_longC);
525
526  // Loop over radius weight points
527  for(int i=0; i< (int)weights_shortA.size(); i++) {
528    dp.shortA = weights_shortA[i].value;
529
530    // Loop over longer_edgeB weight points
531    for(int j=0; j< (int)weights_midB.size(); j++) {
532      dp.midB = weights_midB[j].value;
533
534      // Average over longuest_edgeC distribution
535      for(int k=0; k< (int)weights_longC.size(); k++) {
536        dp.longC = weights_longC[k].value;
537        //Calculate surface averaged radius
538        //This is rough approximation.
539        suf_rad = sqrt((dp.shortA*dp.midB+2.0*dp.rimA*dp.midB+2.0*dp.rimA*dp.shortA)/pi);
540        height =(dp.longC + 2.0*dp.rimC);
541        //Note: output of "DiamCyl(dp.length,dp.radius)" is a DIAMETER.
542        sum +=weights_shortA[i].weight* weights_midB[j].weight
543            * weights_longC[k].weight*DiamCyl(height, suf_rad)/2.0;
544        norm += weights_shortA[i].weight* weights_midB[j].weight*weights_longC[k].weight;
545      }
546    }
547  }
548
549  if (norm != 0){
550    //return the averaged value
551    rad_out =  sum/norm;}
552  else{
553    //return normal value
554    //Note: output of "DiamCyl(length,radius)" is DIAMETER.
555    rad_out = DiamCyl(dp.longC, suf_rad)/2.0;}
556  return rad_out;
557
558}
559double CSParallelepipedModel :: calculate_VR() {
560  return 1.0;
561}
Note: See TracBrowser for help on using the repository browser.