source: sasview/realSpaceModeling/pointsmodelpy/libpointsmodelpy/points_model.cc @ b6b9d76

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 b6b9d76 was bbfad0a, checked in by Mathieu Doucet <doucetm@…>, 17 years ago

Updated realSpaceModeling to be thread-friendly

  • Property mode set to 100644
File size: 11.6 KB
Line 
1/** \file points_model.cc */
2
3#include <vector>
4#include <algorithm>
5#include <fstream>
6#include <stdio.h>
7//#include <exception>
8#include <stdexcept>
9#include "points_model.h"
10#include "Point3D.h"
11
12PointsModel::PointsModel()
13{
14  r_grids_num_ = 2000;
15  rmax_ = 0;
16  cormax_ = 0;
17  rstep_ = 0;
18}
19
20void PointsModel::CalculateIQ(IQ *iq)
21{
22  //fourier transform of the returned Array2D<double> from ddFunction()
23  int nIpoints = iq->GetNumI();
24  double qstep = (iq->GetQmax()) / (nIpoints-1);
25  vector<double> fint(nIpoints, 0);
26
27  //I(0) is calculated seperately
28  int num_rstep = pr_.dim1();
29
30  for (int k = 1; k<nIpoints; k++){
31
32    double q = k * qstep;
33
34    double r =0;
35    double debeye = 0;
36    double fadd =0;
37
38
39    for (int i = 1; i < num_rstep; ++i){
40      r = i*rstep_;  //r should start from 1* rstep
41      double qr = q*r;
42      debeye = sin(qr)/qr;
43      fadd = pr_[i][1]*debeye;
44      fint[k] = fint[k] + fadd;
45    }
46  }
47
48  //I(0)
49  double Izero = 0;
50  for (int i = 0; i < num_rstep; ++i)
51    Izero += pr_[i][1];
52  fint[0] = Izero;
53
54  //assign I(Q) with normalization
55  for(int j = 0; j < nIpoints; ++j){
56    (*iq).iq_data[j][0] = j * qstep;
57    (*iq).iq_data[j][1] = fint[j];
58    // remove normalization Izero;
59  }
60}
61
62//return I with a single q value
63double PointsModel::CalculateIQ(double q)
64{
65  //fourier transform of the returned Array2D<double> from ddFunction()
66  int num_rstep = pr_.dim1();
67
68  double r =0;
69  double debeye = 0;
70  double fadd = 0;
71  double Irelative = 0;
72
73  //I(0) is calculated seperately
74  if (q == 0){
75    //I(0)
76    double Izero = 0;
77    for (int i = 0; i < num_rstep; ++i)
78      Izero += pr_[i][1];
79    Irelative = Izero;
80  }
81  else {
82    for (int i = 1; i < num_rstep; ++i){
83      r = i*rstep_;  //r should start from 1* rstep
84      double qr = q*r;
85      debeye = sin(qr)/qr;
86      fadd = pr_[i][1]*debeye;
87      Irelative = Irelative + fadd;
88    }
89  }
90  return Irelative;
91}
92
93double PointsModel::CalculateIQError(double q)
94{
95  //fourier transform of the returned Array2D<double> from ddFunction()
96  int num_rstep = pr_.dim1();
97
98  double r =0;
99  double debeye = 0;
100  double fadd = 0;
101  double Irelative = 0;
102
103  //I(0) is calculated seperately
104  for (int i = 1; i < num_rstep; ++i){
105          r = i*rstep_;  //r should start from 1* rstep
106          double qr = q*r;
107          debeye = sin(qr)/qr;
108          fadd = fabs(pr_[i][2])*debeye*debeye
109                + rstep_*rstep_/4.0/r/r*(cos(qr)*cos(qr) + debeye*debeye);
110          Irelative = Irelative + fadd;
111  }
112  return sqrt(Irelative);
113}
114
115//pass in a vector of points, and calculate the P(r)
116double PointsModel::DistDistribution(const vector<Point3D> &vp)
117{
118  //get r axis:0,rstep,2rstep,3rstep......d_bound
119  int sizeofpr = r_grids_num_ + 1; //+1 just for overflow prevention
120
121  double d_bound = GetDimBound();
122  rstep_ = CalculateRstep(r_grids_num_,d_bound);
123
124  Array2D<double> pr(sizeofpr, 3); //third column is left for error for the future
125  pr = 0;
126
127  for (int i = 1; i != sizeofpr; ++i)
128    pr[i][0] = pr[i-1][0] + rstep_ ;  //column 1:  distance
129
130  int size = vp.size();
131
132  for (int i1 = 0; i1 < size - 1; ++i1) {
133    for (int i2 = i1 + 1; i2 < size; ++i2) {
134      //dist_.push_back(vp[i1].distanceToPoint(vp[i2]));
135      //product_sld_.push_back(vp[i1].getSLD() * vp[i2].getSLD());
136      double a_dist = vp[i1].distanceToPoint(vp[i2]);
137      double its_sld = vp[i1].getSLD() * vp[i2].getSLD();
138
139      //save maximum distance
140      if (a_dist>rmax_) {
141        rmax_ = a_dist;
142      }
143      //insert into pr array
144      int l = int(floor(a_dist/rstep_));
145
146      //cout << "i1,i2,l,a_dist"<<vp[i1]<<" "<<vp[i2]<<" "<<l<<" "<<a_dist<<endl;     
147      //overflow check
148      if (l >= sizeofpr) {
149        cerr << "one distance is out of range: " << l <<endl;
150        //throw "Out of range";
151      }
152      else {
153        pr[l][1] += its_sld;  //column 2 intermediate: sum of SLD of the points with specific distance     
154        // Estimate uncertainty (squared)
155        pr[l][2] += its_sld*its_sld; 
156        //keep maxium Pr absolute number, in order to normalize
157        //if (pr[l][1] > cormax_)         cormax_ = pr[l][1];
158      }
159    }
160  }
161 
162  //normalize Pr
163  for (int j = 0; j != sizeofpr; ++j){           //final column2 for P(r)
164    //pr[j][1] = pr[j][1]/cormax_;
165
166    // 'Size' is the number of space points, without double counting (excluding
167    // overlapping regions between shapes). The volume of the combined shape
168    // is given by  V = size * (sum of all sub-volumes) / (Total number of points)
169    //              V = size / (lores_density)
170   
171    // - To transform the integral to a sum, we need to give a weight
172    // to each entry equal to the average space volume of a point (w = V/N = 1/lores_density).
173    // The final output, I(q), should therefore be multiplied by V*V/N*N.
174    // Since we will be interested in P(r)/V, we only need to multiply by 1/N*(V/N) = 1/N/lores_density.
175    // We don't have access to lores_density from this class; we will therefore apply
176    // this correction externally.
177    //
178    // - Since the loop goes through half the points, multiply by 2.
179    // TODO: have access to lores_density from this class.
180    //
181    pr[j][1] = 2.0*pr[j][1]/size;
182    pr[j][2] = 4.0*pr[j][2]/size/size;
183  }
184  pr_ = pr;
185
186  return rmax_;
187}
188
189Array2D<double> PointsModel::GetPr()
190{
191  return pr_;
192}
193
194
195double PointsModel::CalculateRstep(int num_grids, double rmax)
196{
197  assert(num_grids > 0);
198
199  double rstep;
200  rstep = rmax / num_grids;
201
202  return rstep;
203}
204
205void PointsModel::OutputPR(const string &fpr){
206  ofstream outfile(fpr.c_str());
207  if (!outfile) {
208    cerr << "error: unable to open output file: "
209         << outfile << endl;
210    exit(1);
211  }
212  for (int i = 0;  i < pr_.dim1(); ++i){
213    outfile << pr_[i][0] << "       " << pr_[i][1] << endl;
214  }
215}
216
217void PointsModel::OutputPDB(const vector<Point3D> &vp,const char *fpr){
218  FILE *outfile=NULL;
219  outfile = fopen(fpr,"w+");
220  if (!outfile) {
221    cerr << "error: unable to open output file: "
222         << outfile << endl;
223    exit(1);
224  }
225  int size = vp.size();
226  int index = 0;
227  for (int i = 0;  i < size; ++i){
228    ++index;
229    fprintf(outfile,"ATOM%7d  C%24.3lf%8.3lf%8.3lf%6.3lf\n", \
230           index,vp[i].getX(),vp[i].getY(),vp[i].getZ(),vp[i].getSLD());
231  }
232  fclose(outfile);
233}
234
235PointsModel::~PointsModel()
236{
237}
238
239void PointsModel::DistDistributionXY(const vector<Point3D> &vp)
240{
241  //the max box get from 3D should be more than enough for 2D,but doesn't hurt
242  double d_bound = GetDimBound();
243
244  //using 1A for rstep, so the total bins is the max distance for the object
245  int sizeofpr = ceil(d_bound) + 1; //+1 just for overflow prevention
246  rstep_ = 1;
247
248  Array2D<double> pr_xy(sizeofpr,sizeofpr); //2D histogram
249
250  //the max frequency in the correlation histogram
251  double cormax_xy_ = 0;
252
253  //initialization
254  pr_xy = 0;
255
256  for (int i = 1; i != sizeofpr; ++i){
257    pr_xy[i][0] = pr_xy[i-1][0] + rstep_ ;  //column 1:  distance
258  }
259
260  int size = vp.size();
261
262  for (int i1 = 0; i1 < size - 1; ++i1) {
263    for (int i2 = i1 + 1; i2 < size; ++i2) {
264      int jx = int(floor(fabs(vp[i1].getX()-vp[i2].getX())/rstep_));
265      int jy = int(floor(fabs(vp[i1].getY()-vp[i2].getY())/rstep_));
266      //the sld for the pair of points
267      double its_sld = vp[i1].getSLD()*vp[i2].getSLD();
268
269      //overflow check
270      if ((jx >= sizeofpr) || (jy >= sizeofpr))
271      {
272        cerr << "one distance is out of range: " <<endl;
273        //throw "Out of range";
274      }
275      else{
276        pr_xy[jx][jy] += its_sld;
277        if (pr_xy[jx][jy] > cormax_xy_ ) cormax_xy_ = pr_xy[jx][jy];
278      }
279    }
280  }
281
282  //normalize Pr_xy
283  for (int m = 0; m != sizeofpr; ++m){           //final column2 for P(r)
284    for (int n = 0; n != sizeofpr; ++n){
285      pr_xy[m][n] = pr_xy[m][n]/cormax_xy_;
286      //cout << "m n:"<<m<<" "<<n<<" "<<pr_xy[m][n]<<endl;
287    }
288  }
289
290  pr_xy_ = pr_xy;
291}
292
293void PointsModel::OutputPR_XY(const std::string &fpr)
294{
295  ofstream outfile(fpr.c_str());
296  if (!outfile) {
297    cerr << "error: unable to open output file: "
298         << outfile << endl;
299    exit(1);
300  }
301  int size = pr_xy_.dim1();
302  //pr_xy_ is a N x N array
303  for (int i = 0;  i != size; ++i){
304    for (int j = 0; j != size; ++j)
305    {
306      outfile << i << "    " << j <<"    "<< pr_xy_[i][j] << endl;
307    }
308  }
309}
310
311void PointsModel::CalculateIQ_2D(IQ *iq,double phi)
312{
313  int nIpoints = iq->GetNumI();
314  double qstep = (iq->GetQmax()) / (nIpoints-1);
315  vector<double> fint(nIpoints, 0);
316  double Izero = 0;
317
318  //number of bins on x and y axis
319  int size_r = pr_xy_.dim1();
320  //rstep is set to one, otherwise should be cos(phi)*rstep
321  double cosphi = cos(phi);
322  double sinphi = sin(phi);
323
324  for(int k = 1; k != nIpoints; ++k){
325    double q = k * qstep;
326    double tmp = cos(q*(cosphi+sinphi));
327
328    for(int i=0; i!=size_r; ++i){
329      for(int j = 0; j!=size_r; ++j){
330        fint[k] += pr_xy_[i][j]*tmp;
331      }
332    }
333  }
334
335  for(int i=0; i!=size_r; ++i){
336    for(int j = 0; j!=size_r; ++j){
337      Izero += pr_xy_[i][j];
338    }
339  }
340  fint[0] = Izero;
341
342  //assign I(Q) with normalization
343  for(int j = 0; j < nIpoints; ++j){
344    (*iq).iq_data[j][0] = j * qstep;
345    (*iq).iq_data[j][1] = fint[j] / Izero;
346  }
347}
348
349vector<double> PointsModel::GetCenter()
350{
351  vector<double> vp(3,0);
352  return vp;
353}
354
355double PointsModel::CalculateIQ_2D(double qx, double qy)
356{
357  //for each (Qx,Qy) on 2D detector, calculate I
358  double q = sqrt(qx*qx+qy*qy);
359  double I = 0;
360
361  double cosphi = qx/q;
362  double sinphi = qy/q;
363  double tmp = cos(q*(cosphi+sinphi));
364
365  //loop through P(r) on xy plane
366  int size_r = pr_xy_.dim1();
367  for(int i=-size_r+1; i!=size_r; ++i){
368    for(int j = -size_r+1; j!=size_r; ++j){
369      //rstep is set to one, left out from calculation
370      I += pr_xy_[abs(i)][abs(j)]*cos(q*(cosphi*i+sinphi*j));
371    }
372  }
373
374  //return I, without normalization
375  return I;
376}
377
378/*
379 * 2D simulation for oriented systems
380 * The beam direction is assumed to be in the z direction.
381 *
382 * @param points: vector of space points
383 * @param qx: qx [A-1]
384 * @param qy: qy [A-1]
385 * @return: I(qx, qy) for the system described by the space points [cm-1]
386 *
387 */
388double PointsModel::CalculateIQ_2D(const vector<Point3D>&points, double qx, double qy){
389        /*
390         * TODO: the vector of points should really be part of the class
391         *           This is a design flaw inherited from the original programmer.
392         */
393       
394        int size = points.size();
395
396        double cos_term = 0;
397        double sin_term = 0;
398        for (int i = 0; i < size; i++) {
399                //the sld for the pair of points
400       
401                double phase = qx*points[i].getX() + qy*points[i].getY();
402                               
403                cos_term += cos(phase) * points[i].getSLD();
404                sin_term += sin(phase) * points[i].getSLD();
405       
406        }                       
407
408        // P(q) = 1/V I(q) = (V/N)^2 (1/V) (cos_term^2 + sin_term^2)
409        // We divide by N here and we will multiply by the density later.
410
411        return (cos_term*cos_term + sin_term*sin_term)/size;   
412}
413
414double PointsModel::CalculateIQ_2D_Error(const vector<Point3D>&points, double qx, double qy){
415       
416        int size = points.size();
417
418        double delta_x, delta_y;
419        double q_t2 = qx*qx + qy*qy;
420        double cos_term = 0;
421        double sin_term = 0;
422        double cos_err = 0;
423        double sin_err = 0;
424       
425        // Estimate the error on the position of each point
426        // in x or y as V^(1/3)/N
427       
428        for (int i = 0; i < size; i++) {
429               
430               
431                //the sld for the pair of points
432       
433                double phase = qx*points[i].getX() + qy*points[i].getY();
434                double sld_fac = points[i].getSLD() * points[i].getSLD();
435               
436                cos_term += cos(phase) * points[i].getSLD();
437                sin_term += sin(phase) * points[i].getSLD();
438
439                sin_err += cos(phase) * cos(phase) * sld_fac;
440                cos_err += sin(phase) * sin(phase) * sld_fac;
441       
442        }                       
443
444        // P(q) = 1/V I(q) = (V/N)^2 (1/V) (cos_term^2 + sin_term^2)
445        // We divide by N here and we will multiply by the density later.
446
447        // We will need to multiply this error by V^(1/3)/N.
448        // We don't have access to V from within this class.
449        return 2*sqrt(cos_term*cos_term*cos_err*cos_err + sin_term*sin_term*sin_err*sin_err)/size;     
450}
451
Note: See TracBrowser for help on using the repository browser.