source: sasview/sansmodels/src/c_models/sc.cpp @ e08bd5b

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 e08bd5b was e08bd5b, checked in by Jae Cho <jhjcho@…>, 12 years ago

c models fix: scale fix for P*S

  • Property mode set to 100644
File size: 9.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 */
21
22#include <math.h>
23#include "parameters.hh"
24#include <stdio.h>
25using namespace std;
26#include "sc.h"
27
28extern "C" {
29#include "libSphere.h"
30}
31
32// Convenience structure
33typedef struct {
34  double scale;
35  double dnn;
36  double d_factor;
37  double radius;
38  double sldSph;
39  double sldSolv;
40  double background;
41  double theta;
42  double phi;
43  double psi;
44} SCParameters;
45
46/**
47 * Function to evaluate 2D scattering function
48 * @param pars: parameters of the SCCrystalModel
49 * @param q: q-value
50 * @param q_x: q_x / q
51 * @param q_y: q_y / q
52 * @return: function value
53 */
54static double sc_analytical_2D_scaled(SCParameters *pars, double q, double q_x, double q_y) {
55  double a3_x, a3_y, a3_z, a2_x, a2_y, a1_x, a1_y;
56  double q_z;
57  double alpha, cos_val_a3, cos_val_a2, cos_val_a1;
58  double a1_dot_q, a2_dot_q,a3_dot_q;
59  double answer;
60  double Pi = 4.0*atan(1.0);
61  double aa, Da, qDa_2, latticeScale, Zq;
62
63  double dp[5];
64  //convert angle degree to radian
65  double theta = pars->theta * Pi/180.0;
66  double phi = pars->phi * Pi/180.0;
67  double psi = pars->psi * Pi/180.0;
68  dp[0] = 1.0;
69  dp[1] = pars->radius;
70  dp[2] = pars->sldSph;
71  dp[3] = pars->sldSolv;
72  dp[4] = 0.0;
73
74
75  aa = pars->dnn;
76  Da = pars->d_factor*aa;
77  qDa_2 = pow(q*Da,2.0);
78
79  latticeScale = (4.0/3.0)*Pi*(dp[1]*dp[1]*dp[1])/pow(aa,3.0);
80  /// Angles here are respect to detector coordinate instead of against q coordinate(PRB 36, 3, 1754)
81  // a3 axis orientation
82  a3_x = sin(theta) * cos(phi);//negative sign here???
83  a3_y = sin(theta) * sin(phi);
84  a3_z = cos(theta);
85
86  // q vector
87  q_z = 0.0;
88
89  // Compute the angle btw vector q and the a3 axis
90  cos_val_a3 = a3_x*q_x + a3_y*q_y + a3_z*q_z;
91  alpha = acos(cos_val_a3);
92  //alpha = acos(cos_val_a3);
93  a3_dot_q = aa*q*cos_val_a3;
94  // a1 axis orientation
95  a1_x = sin(psi);
96  a1_y = cos(psi);
97
98  cos_val_a1 = a1_x*q_x + a1_y*q_y;
99  a1_dot_q = aa*q*cos_val_a1*sin(alpha);
100
101  // a2 axis orientation
102  a2_x = sqrt(1.0-sin(theta)*cos(phi))*cos(psi);
103  a2_y = sqrt(1.0-sin(theta)*cos(phi))*sin(psi);
104  // a2 axis
105  cos_val_a2 =  sin(acos(cos_val_a1));//a2_x*q_x + a2_y*q_y;
106  a2_dot_q = aa*q*cos_val_a2*sin(alpha);
107
108  // The following test should always pass
109  if (fabs(cos_val_a3)>1.0) {
110    printf("parallel_ana_2D: Unexpected error: cos(alpha)>1\n");
111    return 0;
112  }
113  // Call Zq=Z1*Z2*Z3
114  Zq = (1.0-exp(-qDa_2))/(1.0-2.0*exp(-0.5*qDa_2)*cos(a1_dot_q)+exp(-qDa_2));
115  Zq = Zq * (1.0-exp(-qDa_2))/(1.0-2.0*exp(-0.5*qDa_2)*cos(a2_dot_q)+exp(-qDa_2));
116  Zq = Zq * (1.0-exp(-qDa_2))/(1.0-2.0*exp(-0.5*qDa_2)*cos(a3_dot_q)+exp(-qDa_2));
117
118  // Use SphereForm directly from libigor
119  answer = SphereForm(dp,q)*Zq;
120
121  //consider scales
122  answer *= latticeScale * pars->scale;
123
124  // This FIXES a singualrity the kernel in libigor.
125  if ( answer == INFINITY || answer == NAN){
126    answer = 0.0;
127  }
128
129  // add background
130  answer += pars->background;
131
132  return answer;
133}
134
135/**
136 * Function to evaluate 2D scattering function
137 * @param pars: parameters of the SC_ParaCrystal
138 * @param q: q-value
139 * @return: function value
140 */
141static double sc_analytical_2DXY(SCParameters *pars, double qx, double qy){
142  double q;
143  q = sqrt(qx*qx+qy*qy);
144  return sc_analytical_2D_scaled(pars, q, qx/q, qy/q);
145}
146
147SCCrystalModel :: SCCrystalModel() {
148  scale      = Parameter(1.0);
149  dnn           = Parameter(220.0);
150  d_factor = Parameter(0.06);
151  radius     = Parameter(40.0, true);
152  radius.set_min(0.0);
153  sldSph   = Parameter(3.0e-6);
154  sldSolv   = Parameter(6.3e-6);
155  background = Parameter(0.0);
156  theta  = Parameter(0.0, true);
157  phi    = Parameter(0.0, true);
158  psi    = Parameter(0.0, true);
159}
160
161/**
162 * Function to evaluate 1D scattering function
163 * The NIST IGOR library is used for the actual calculation.
164 * @param q: q-value
165 * @return: function value
166 */
167double SCCrystalModel :: operator()(double q) {
168  double dp[7];
169
170  // Fill parameter array for IGOR library
171  // Add the background after averaging
172  dp[0] = scale();
173  dp[1] = dnn();
174  dp[2] = d_factor();
175  dp[3] = radius();
176  dp[4] = sldSph();
177  dp[5] = sldSolv();
178  dp[6] = 0.0;
179
180  // Get the dispersion points for the radius
181  vector<WeightPoint> weights_rad;
182  radius.get_weights(weights_rad);
183
184  // Perform the computation, with all weight points
185  double sum = 0.0;
186  double norm = 0.0;
187  double vol = 0.0;
188  double result;
189
190  // Loop over radius weight points
191  for(size_t i=0; i<weights_rad.size(); i++) {
192    dp[3] = weights_rad[i].value;
193
194    //Un-normalize SphereForm by volume
195    result = SC_ParaCrystal(dp, q) * pow(weights_rad[i].value,3);
196    // This FIXES a singualrity the kernel in libigor.
197    if ( result == INFINITY || result == NAN){
198      result = 0.0;
199    }
200    sum += weights_rad[i].weight
201        * result;
202    //Find average volume
203    vol += weights_rad[i].weight
204        * pow(weights_rad[i].value,3);
205
206    norm += weights_rad[i].weight;
207  }
208
209  if (vol != 0.0 && norm != 0.0) {
210    //Re-normalize by avg volume
211    sum = sum/(vol/norm);}
212  return sum/norm + background();
213}
214
215/**
216 * Function to evaluate 2D scattering function
217 * @param q_x: value of Q along x
218 * @param q_y: value of Q along y
219 * @return: function value
220 */
221double SCCrystalModel :: operator()(double qx, double qy) {
222  SCParameters dp;
223  dp.scale      = scale();
224  dp.dnn   = dnn();
225  dp.d_factor   = d_factor();
226  dp.radius  = radius();
227  dp.sldSph   = sldSph();
228  dp.sldSolv   = sldSolv();
229  dp.background = 0.0;
230  dp.theta  = theta();
231  dp.phi    = phi();
232  dp.psi    = psi();
233
234  // Get the dispersion points for the radius
235  vector<WeightPoint> weights_rad;
236  radius.get_weights(weights_rad);
237
238  // Get angular averaging for theta
239  vector<WeightPoint> weights_theta;
240  theta.get_weights(weights_theta);
241
242  // Get angular averaging for phi
243  vector<WeightPoint> weights_phi;
244  phi.get_weights(weights_phi);
245
246  // Get angular averaging for psi
247  vector<WeightPoint> weights_psi;
248  psi.get_weights(weights_psi);
249
250  // Perform the computation, with all weight points
251  double sum = 0.0;
252  double norm = 0.0;
253  double norm_vol = 0.0;
254  double vol = 0.0;
255  double pi = 4.0*atan(1.0);
256  // Loop over radius weight points
257  for(size_t i=0; i<weights_rad.size(); i++) {
258    dp.radius = weights_rad[i].value;
259    // Average over theta distribution
260    for(size_t j=0; j< weights_theta.size(); j++) {
261      dp.theta = weights_theta[j].value;
262      // Average over phi distribution
263      for(size_t k=0; k< weights_phi.size(); k++) {
264        dp.phi = weights_phi[k].value;
265        // Average over phi distribution
266        for(size_t l=0; l< weights_psi.size(); l++) {
267          dp.psi = weights_psi[l].value;
268          //Un-normalize SphereForm by volume
269          double _ptvalue = weights_rad[i].weight
270              * weights_theta[j].weight
271              * weights_phi[k].weight
272              * weights_psi[l].weight
273              * sc_analytical_2DXY(&dp, qx, qy);
274          //* pow(weights_rad[i].value,3.0);
275          // Consider when there is infinte or nan.
276          if ( _ptvalue == INFINITY || _ptvalue == NAN){
277            _ptvalue = 0.0;
278          }
279          if (weights_theta.size()>1) {
280            _ptvalue *= fabs(sin(weights_theta[j].value*pi/180.0));
281          }
282          sum += _ptvalue;
283          // This model dose not need the volume of spheres correction!!!
284          //Find average volume
285          //vol += weights_rad[i].weight
286          //    * pow(weights_rad[i].value,3);
287          //Find norm for volume
288          //norm_vol += weights_rad[i].weight;
289
290          norm += weights_rad[i].weight
291              * weights_theta[j].weight
292              * weights_phi[k].weight
293              * weights_psi[l].weight;
294        }
295      }
296    }
297  }
298  // Averaging in theta needs an extra normalization
299  // factor to account for the sin(theta) term in the
300  // integration (see documentation).
301  if (weights_theta.size()>1) norm = norm / asin(1.0);
302
303  if (vol != 0.0 && norm_vol != 0.0) {
304    //Re-normalize by avg volume
305    sum = sum/(vol/norm_vol);}
306
307  return sum/norm + background();
308}
309
310/**
311 * Function to evaluate 2D scattering function
312 * @param pars: parameters of the SCCrystal
313 * @param q: q-value
314 * @param phi: angle phi
315 * @return: function value
316 */
317double SCCrystalModel :: evaluate_rphi(double q, double phi) {
318  return (*this).operator()(q);
319}
320
321/**
322 * Function to calculate effective radius
323 * @return: effective radius value
324 */
325double SCCrystalModel :: calculate_ER() {
326  //NOT implemented yet!!!
327  return 0.0;
328}
329double SCCrystalModel :: calculate_VR() {
330  return 1.0;
331}
Note: See TracBrowser for help on using the repository browser.