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