source: sasview/sansmodels/src/c_models/fcc.cpp @ f518eef

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