source: sasview/sansmodels/src/c_models/parallelepiped.cpp @ df88829

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 df88829 was 8343e18, checked in by Mathieu Doucet <doucetm@…>, 12 years ago

refactored parallelepiped

  • Property mode set to 100644
File size: 13.6 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 "models.hh"
25#include "parameters.hh"
26#include <stdio.h>
27using namespace std;
28
29extern "C" {
30        #include "libCylinder.h"
31        #include "libStructureFactor.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} ParallelepipedParameters;
48
49
50static double pkernel(double a, double b,double c, double ala, double alb, double alc){
51    // mu passed in is really mu*sqrt(1-sig^2)
52    double argA,argB,argC,tmp1,tmp2,tmp3;     //local variables
53
54    //handle arg=0 separately, as sin(t)/t -> 1 as t->0
55    argA = a*ala/2.0;
56    argB = b*alb/2.0;
57    argC = c*alc/2.0;
58    if(argA==0.0) {
59    tmp1 = 1.0;
60  } else {
61    tmp1 = sin(argA)*sin(argA)/argA/argA;
62    }
63    if (argB==0.0) {
64    tmp2 = 1.0;
65  } else {
66    tmp2 = sin(argB)*sin(argB)/argB/argB;
67    }
68
69    if (argC==0.0) {
70    tmp3 = 1.0;
71  } else {
72    tmp3 = sin(argC)*sin(argC)/argC/argC;
73    }
74
75    return (tmp1*tmp2*tmp3);
76
77}
78
79/**
80 * Function to evaluate 2D scattering function
81 * @param pars: parameters of the parallelepiped
82 * @param q: q-value
83 * @param q_x: q_x / q
84 * @param q_y: q_y / q
85 * @return: function value
86 */
87static double parallelepiped_analytical_2D_scaled(ParallelepipedParameters *pars, double q, double q_x, double q_y) {
88  double cparallel_x, cparallel_y, cparallel_z, bparallel_x, bparallel_y, parallel_x, parallel_y;
89  double q_z;
90  double alpha, vol, cos_val_c, cos_val_b, cos_val_a, edgeA, edgeB, edgeC;
91
92  double answer;
93  double pi = 4.0*atan(1.0);
94
95  //convert angle degree to radian
96  double theta = pars->parallel_theta * pi/180.0;
97  double phi = pars->parallel_phi * pi/180.0;
98  double psi = pars->parallel_psi * pi/180.0;
99
100  edgeA = pars->short_a;
101  edgeB = pars->short_b;
102  edgeC = pars->long_c;
103
104
105    // parallelepiped c axis orientation
106    cparallel_x = sin(theta) * cos(phi);
107    cparallel_y = sin(theta) * sin(phi);
108    cparallel_z = cos(theta);
109
110    // q vector
111    q_z = 0.0;
112
113    // Compute the angle btw vector q and the
114    // axis of the parallelepiped
115    cos_val_c = cparallel_x*q_x + cparallel_y*q_y + cparallel_z*q_z;
116    alpha = acos(cos_val_c);
117
118    // parallelepiped a axis orientation
119    parallel_x = sin(psi);//cos(pars->parallel_theta) * sin(pars->parallel_phi)*sin(pars->parallel_psi);
120    parallel_y = cos(psi);//cos(pars->parallel_theta) * cos(pars->parallel_phi)*cos(pars->parallel_psi);
121
122    cos_val_a = parallel_x*q_x + parallel_y*q_y;
123
124
125
126    // parallelepiped b axis orientation
127    bparallel_x = sqrt(1.0-sin(theta)*cos(phi))*cos(psi);//cos(pars->parallel_theta) * cos(pars->parallel_phi)* cos(pars->parallel_psi);
128    bparallel_y = sqrt(1.0-sin(theta)*cos(phi))*sin(psi);//cos(pars->parallel_theta) * sin(pars->parallel_phi)* sin(pars->parallel_psi);
129    // axis of the parallelepiped
130    cos_val_b = sin(acos(cos_val_a)) ;
131
132
133
134    // The following test should always pass
135    if (fabs(cos_val_c)>1.0) {
136      printf("parallel_ana_2D: Unexpected error: cos(alpha)>1\n");
137      return 0;
138    }
139
140  // Call the IGOR library function to get the kernel
141  answer = pkernel( q*edgeA, q*edgeB, q*edgeC, sin(alpha)*cos_val_a,sin(alpha)*cos_val_b,cos_val_c);
142
143  // Multiply by contrast^2
144  answer *= (pars->sldPipe - pars->sldSolv) * (pars->sldPipe - pars->sldSolv);
145
146  //normalize by cylinder volume
147  //NOTE that for this (Fournet) definition of the integral, one must MULTIPLY by Vparallel
148    vol = edgeA* edgeB * edgeC;
149  answer *= vol;
150
151  //convert to [cm-1]
152  answer *= 1.0e8;
153
154  //Scale
155  answer *= pars->scale;
156
157  // add in the background
158  answer += pars->background;
159
160  return answer;
161}
162
163/**
164 * Function to evaluate 2D scattering function
165 * @param pars: parameters of the parallelepiped
166 * @param q: q-value
167 * @return: function value
168 */
169static double parallelepiped_analytical_2DXY(ParallelepipedParameters *pars, double qx, double qy) {
170  double q;
171  q = sqrt(qx*qx+qy*qy);
172    return parallelepiped_analytical_2D_scaled(pars, q, qx/q, qy/q);
173}
174
175ParallelepipedModel :: ParallelepipedModel() {
176        scale      = Parameter(1.0);
177        short_a     = Parameter(35.0, true);
178        short_a.set_min(1.0);
179        short_b     = Parameter(75.0, true);
180        short_b.set_min(1.0);
181        long_c     = Parameter(400.0, true);
182        long_c.set_min(1.0);
183        sldPipe   = Parameter(6.3e-6);
184        sldSolv   = Parameter(1.0e-6);
185        background = Parameter(0.0);
186        parallel_theta  = Parameter(0.0, true);
187        parallel_phi    = Parameter(0.0, true);
188        parallel_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 ParallelepipedModel :: operator()(double q) {
198        double dp[7];
199
200        // Fill parameter array for IGOR library
201        // Add the background after averaging
202        dp[0] = scale();
203        dp[1] = short_a();
204        dp[2] = short_b();
205        dp[3] = long_c();
206        dp[4] = sldPipe();
207        dp[5] = sldSolv();
208        dp[6] = 0.0;
209
210        // Get the dispersion points for the short_edgeA
211        vector<WeightPoint> weights_short_a;
212        short_a.get_weights(weights_short_a);
213
214        // Get the dispersion points for the longer_edgeB
215        vector<WeightPoint> weights_short_b;
216        short_b.get_weights(weights_short_b);
217
218        // Get the dispersion points for the longuest_edgeC
219        vector<WeightPoint> weights_long_c;
220        long_c.get_weights(weights_long_c);
221
222
223
224        // Perform the computation, with all weight points
225        double sum = 0.0;
226        double norm = 0.0;
227        double vol = 0.0;
228
229        // Loop over short_edgeA weight points
230        for(int i=0; i< (int)weights_short_a.size(); i++) {
231                dp[1] = weights_short_a[i].value;
232
233                // Loop over longer_edgeB weight points
234                for(int j=0; j< (int)weights_short_b.size(); j++) {
235                        dp[2] = weights_short_b[j].value;
236
237                        // Loop over longuest_edgeC weight points
238                        for(int k=0; k< (int)weights_long_c.size(); k++) {
239                                dp[3] = weights_long_c[k].value;
240                                //Un-normalize  by volume
241                                sum += weights_short_a[i].weight * weights_short_b[j].weight
242                                        * weights_long_c[k].weight * Parallelepiped(dp, q)
243                                        * weights_short_a[i].value*weights_short_b[j].value
244                                        * weights_long_c[k].value;
245                                //Find average volume
246                                vol += weights_short_a[i].weight * weights_short_b[j].weight
247                                        * weights_long_c[k].weight
248                                        * weights_short_a[i].value * weights_short_b[j].value
249                                        * weights_long_c[k].value;
250
251                                norm += weights_short_a[i].weight
252                                         * weights_short_b[j].weight * weights_long_c[k].weight;
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 * 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 ParallelepipedModel :: operator()(double qx, double qy) {
269        ParallelepipedParameters dp;
270        // Fill parameter array
271        dp.scale      = scale();
272        dp.short_a   = short_a();
273        dp.short_b   = short_b();
274        dp.long_c  = long_c();
275        dp.sldPipe   = sldPipe();
276        dp.sldSolv   = sldSolv();
277        dp.background = 0.0;
278        //dp.background = background();
279        dp.parallel_theta  = parallel_theta();
280        dp.parallel_phi    = parallel_phi();
281        dp.parallel_psi    = parallel_psi();
282
283
284        // Get the dispersion points for the short_edgeA
285        vector<WeightPoint> weights_short_a;
286        short_a.get_weights(weights_short_a);
287
288        // Get the dispersion points for the longer_edgeB
289        vector<WeightPoint> weights_short_b;
290        short_b.get_weights(weights_short_b);
291
292        // Get angular averaging for the longuest_edgeC
293        vector<WeightPoint> weights_long_c;
294        long_c.get_weights(weights_long_c);
295
296        // Get angular averaging for theta
297        vector<WeightPoint> weights_parallel_theta;
298        parallel_theta.get_weights(weights_parallel_theta);
299
300        // Get angular averaging for phi
301        vector<WeightPoint> weights_parallel_phi;
302        parallel_phi.get_weights(weights_parallel_phi);
303
304        // Get angular averaging for psi
305        vector<WeightPoint> weights_parallel_psi;
306        parallel_psi.get_weights(weights_parallel_psi);
307
308        // Perform the computation, with all weight points
309        double sum = 0.0;
310        double norm = 0.0;
311        double norm_vol = 0.0;
312        double vol = 0.0;
313        double pi = 4.0*atan(1.0);
314        // Loop over radius weight points
315        for(int i=0; i< (int)weights_short_a.size(); i++) {
316                dp.short_a = weights_short_a[i].value;
317
318                // Loop over longer_edgeB weight points
319                for(int j=0; j< (int)weights_short_b.size(); j++) {
320                        dp.short_b = weights_short_b[j].value;
321
322                        // Average over longuest_edgeC distribution
323                        for(int k=0; k< (int)weights_long_c.size(); k++) {
324                                dp.long_c = weights_long_c[k].value;
325
326                                // Average over theta distribution
327                                for(int l=0; l< (int)weights_parallel_theta.size(); l++) {
328                                dp.parallel_theta = weights_parallel_theta[l].value;
329
330                                        // Average over phi distribution
331                                        for(int m=0; m< (int)weights_parallel_phi.size(); m++) {
332                                                dp.parallel_phi = weights_parallel_phi[m].value;
333
334                                                // Average over phi distribution
335                                                for(int n=0; n< (int)weights_parallel_psi.size(); n++) {
336                                                        dp.parallel_psi = weights_parallel_psi[n].value;
337                                                        //Un-normalize by volume
338                                                        double _ptvalue = weights_short_a[i].weight
339                                                                * weights_short_b[j].weight
340                                                                * weights_long_c[k].weight
341                                                                * weights_parallel_theta[l].weight
342                                                                * weights_parallel_phi[m].weight
343                                                                * weights_parallel_psi[n].weight
344                                                                * parallelepiped_analytical_2DXY(&dp, qx, qy)
345                                                                * weights_short_a[i].value*weights_short_b[j].value
346                                                                * weights_long_c[k].value;
347
348                                                        if (weights_parallel_theta.size()>1) {
349                                                                _ptvalue *= fabs(sin(weights_parallel_theta[l].value*pi/180.0));
350                                                        }
351                                                        sum += _ptvalue;
352                                                        //Find average volume
353                                                        vol += weights_short_a[i].weight
354                                                                * weights_short_b[j].weight
355                                                                * weights_long_c[k].weight
356                                                                * weights_short_a[i].value*weights_short_b[j].value
357                                                                * weights_long_c[k].value;
358                                                        //Find norm for volume
359                                                        norm_vol += weights_short_a[i].weight
360                                                                * weights_short_b[j].weight
361                                                                * weights_long_c[k].weight;
362
363                                                        norm += weights_short_a[i].weight
364                                                                * weights_short_b[j].weight
365                                                                * weights_long_c[k].weight
366                                                                * weights_parallel_theta[l].weight
367                                                                * weights_parallel_phi[m].weight
368                                                                * weights_parallel_psi[n].weight;
369                                                }
370                                        }
371
372                                }
373                        }
374                }
375        }
376        // Averaging in theta needs an extra normalization
377        // factor to account for the sin(theta) term in the
378        // integration (see documentation).
379        if (weights_parallel_theta.size()>1) norm = norm / asin(1.0);
380
381        if (vol != 0.0 && norm_vol != 0.0) {
382                //Re-normalize by avg volume
383                sum = sum/(vol/norm_vol);}
384
385        return sum/norm + background();
386}
387
388
389/**
390 * Function to evaluate 2D scattering function
391 * @param pars: parameters of the cylinder
392 * @param q: q-value
393 * @param phi: angle phi
394 * @return: function value
395 */
396double ParallelepipedModel :: evaluate_rphi(double q, double phi) {
397        double qx = q*cos(phi);
398        double qy = q*sin(phi);
399        return (*this).operator()(qx, qy);
400}
401/**
402 * Function to calculate effective radius
403 * @return: effective radius value
404 */
405double ParallelepipedModel :: calculate_ER() {
406        ParallelepipedParameters dp;
407        dp.short_a   = short_a();
408        dp.short_b   = short_b();
409        dp.long_c  = long_c();
410        double rad_out = 0.0;
411        double pi = 4.0*atan(1.0);
412        double suf_rad = sqrt(dp.short_a*dp.short_b/pi);
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 short_edgeA
419        vector<WeightPoint> weights_short_a;
420        short_a.get_weights(weights_short_a);
421
422        // Get the dispersion points for the longer_edgeB
423        vector<WeightPoint> weights_short_b;
424        short_b.get_weights(weights_short_b);
425
426        // Get angular averaging for the longuest_edgeC
427        vector<WeightPoint> weights_long_c;
428        long_c.get_weights(weights_long_c);
429
430        // Loop over radius weight points
431        for(int i=0; i< (int)weights_short_a.size(); i++) {
432                dp.short_a = weights_short_a[i].value;
433
434                // Loop over longer_edgeB weight points
435                for(int j=0; j< (int)weights_short_b.size(); j++) {
436                        dp.short_b = weights_short_b[j].value;
437
438                        // Average over longuest_edgeC distribution
439                        for(int k=0; k< (int)weights_long_c.size(); k++) {
440                                dp.long_c = weights_long_c[k].value;
441                                //Calculate surface averaged radius
442                                //This is rough approximation.
443                                suf_rad = sqrt(dp.short_a*dp.short_b/pi);
444
445                                //Note: output of "DiamCyl(dp.length,dp.radius)" is DIAMETER.
446                                sum +=weights_short_a[i].weight* weights_short_b[j].weight
447                                        * weights_long_c[k].weight*DiamCyl(dp.long_c, suf_rad)/2.0;
448                                norm += weights_short_a[i].weight* weights_short_b[j].weight*weights_long_c[k].weight;
449                        }
450                }
451        }
452
453        if (norm != 0){
454                //return the averaged value
455                rad_out =  sum/norm;}
456        else{
457                //return normal value
458                //Note: output of "DiamCyl(length,radius)" is DIAMETER.
459                rad_out = DiamCyl(dp.long_c, suf_rad)/2.0;}
460        return rad_out;
461
462}
Note: See TracBrowser for help on using the repository browser.