source: sasmodels/sasmodels/kernel_iq.c @ e44432d

core_shell_microgelsmagnetic_modelticket-1257-vesicle-productticket_1156ticket_1265_superballticket_822_more_unit_tests
Last change on this file since e44432d was e44432d, checked in by Paul Kienzle <pkienzle@…>, 6 years ago

support hollow models in structure factor calculations

  • Property mode set to 100644
File size: 29.9 KB
Line 
1/*
2    ##########################################################
3    #                                                        #
4    #   !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!   #
5    #   !!                                              !!   #
6    #   !!  KEEP THIS CODE CONSISTENT WITH KERNELPY.PY  !!   #
7    #   !!                                              !!   #
8    #   !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!   #
9    #                                                        #
10    ##########################################################
11*/
12
13// NOTE: the following macros are defined in generate.py:
14//
15//  MAX_PD : the maximum number of dispersity loops allowed for this model,
16//      which will be at most modelinfo.MAX_PD.
17//  NUM_PARS : the number of parameters in the parameter table
18//  NUM_VALUES : the number of values to skip at the start of the
19//      values array before you get to the dispersity values.
20//  PARAMETER_TABLE : list of parameter declarations used to create the
21//      ParameterTable type.
22//  KERNEL_NAME : model_Iq, model_Iqxy or model_Imagnetic.  This code is
23//      included three times, once for each kernel type.
24//  MAGNETIC : defined when the magnetic kernel is being instantiated
25//  NUM_MAGNETIC : the number of magnetic parameters
26//  MAGNETIC_PARS : a comma-separated list of indices to the sld
27//      parameters in the parameter table.
28//  CALL_VOLUME(form, shell, table) : assign form and shell values
29//  CALL_EFFECTIVE_RADIUS(type, table) : call the R_eff function
30//  CALL_IQ(q, table) : call the Iq function for 1D calcs.
31//  CALL_IQ_A(q, table) : call the Iq function with |q| for 2D data.
32//  CALL_FQ(q, F1, F2, table) : call the Fq function for 1D calcs.
33//  CALL_FQ_A(q, F1, F2, table) : call the Iq function with |q| for 2D data.
34//  CALL_IQ_AC(qa, qc, table) : call the Iqxy function for symmetric shapes
35//  CALL_IQ_ABC(qa, qc, table) : call the Iqxy function for asymmetric shapes
36//  CALL_IQ_XY(qx, qy, table) : call the Iqxy function for arbitrary models
37//  INVALID(table) : test if the current point is feesible to calculate.  This
38//      will be defined in the kernel definition file.
39//  PROJECTION : equirectangular=1, sinusoidal=2
40//      see explore/jitter.py for definitions.
41
42
43#ifndef _PAR_BLOCK_ // protected block so we can include this code twice.
44#define _PAR_BLOCK_
45
46typedef struct {
47#if MAX_PD > 0
48    int32_t pd_par[MAX_PD];     // id of the nth dispersity variable
49    int32_t pd_length[MAX_PD];  // length of the nth dispersity weight vector
50    int32_t pd_offset[MAX_PD];  // offset of pd weights in the value & weight vector
51    int32_t pd_stride[MAX_PD];  // stride to move to the next index at this level
52#endif // MAX_PD > 0
53    int32_t num_eval;           // total number of voxels in hypercube
54    int32_t num_weights;        // total length of the weights vector
55    int32_t num_active;         // number of non-trivial pd loops
56    int32_t theta_par;          // id of first orientation variable
57} ProblemDetails;
58
59// Intel HD 4000 needs private arrays to be a multiple of 4 long
60typedef struct {
61    PARAMETER_TABLE
62} ParameterTable;
63typedef union {
64    ParameterTable table;
65    double vector[4*((NUM_PARS+3)/4)];
66} ParameterBlock;
67#endif // _PAR_BLOCK_
68
69#if defined(MAGNETIC) && NUM_MAGNETIC > 0
70// ===== Helper functions for magnetism =====
71
72// Return value restricted between low and high
73static double clip(double value, double low, double high)
74{
75  return (value < low ? low : (value > high ? high : value));
76}
77
78// Compute spin cross sections given in_spin and out_spin
79// To convert spin cross sections to sld b:
80//     uu * (sld - m_sigma_x);
81//     dd * (sld + m_sigma_x);
82//     ud * (m_sigma_y - 1j*m_sigma_z);
83//     du * (m_sigma_y + 1j*m_sigma_z);
84// weights for spin crosssections: dd du real, ud real, uu, du imag, ud imag
85static void set_spin_weights(double in_spin, double out_spin, double weight[6])
86{
87  in_spin = clip(in_spin, 0.0, 1.0);
88  out_spin = clip(out_spin, 0.0, 1.0);
89  // Previous version of this function took the square root of the weights,
90  // under the assumption that
91  //
92  //     w*I(q, rho1, rho2, ...) = I(q, sqrt(w)*rho1, sqrt(w)*rho2, ...)
93  //
94  // However, since the weights are applied to the final intensity and
95  // are not interned inside the I(q) function, we want the full
96  // weight and not the square root.  Any function using
97  // set_spin_weights as part of calculating an amplitude will need to
98  // manually take that square root, but there is currently no such
99  // function.
100  weight[0] = (1.0-in_spin) * (1.0-out_spin); // dd
101  weight[1] = (1.0-in_spin) * out_spin;       // du
102  weight[2] = in_spin * (1.0-out_spin);       // ud
103  weight[3] = in_spin * out_spin;             // uu
104  weight[4] = weight[1]; // du.imag
105  weight[5] = weight[2]; // ud.imag
106}
107
108// Compute the magnetic sld
109static double mag_sld(
110  const unsigned int xs, // 0=dd, 1=du.real, 2=ud.real, 3=uu, 4=du.imag, 5=ud.imag
111  const double qx, const double qy,
112  const double px, const double py,
113  const double sld,
114  const double mx, const double my, const double mz
115)
116{
117  if (xs < 4) {
118    const double perp = qy*mx - qx*my;
119    switch (xs) {
120      default: // keep compiler happy; condition ensures xs in [0,1,2,3]
121      case 0: // uu => sld - D M_perpx
122          return sld - px*perp;
123      case 1: // ud.real => -D M_perpy
124          return py*perp;
125      case 2: // du.real => -D M_perpy
126          return py*perp;
127      case 3: // dd => sld + D M_perpx
128          return sld + px*perp;
129    }
130  } else {
131    if (xs== 4) {
132      return -mz;  // du.imag => +D M_perpz
133    } else { // index == 5
134      return +mz;  // ud.imag => -D M_perpz
135    }
136  }
137}
138
139
140#endif
141
142// ===== Helper functions for orientation and jitter =====
143
144// To change the definition of the angles, run explore/angles.py, which
145// uses sympy to generate the equations.
146
147#if !defined(_QAC_SECTION) && defined(CALL_IQ_AC)
148#define _QAC_SECTION
149
150typedef struct {
151    double R31, R32;
152} QACRotation;
153
154// Fill in the rotation matrix R from the view angles (theta, phi) and the
155// jitter angles (dtheta, dphi).  This matrix can be applied to all of the
156// (qx, qy) points in the image to produce R*[qx,qy]' = [qa,qc]'
157static void
158qac_rotation(
159    QACRotation *rotation,
160    double theta, double phi,
161    double dtheta, double dphi)
162{
163    double sin_theta, cos_theta;
164    double sin_phi, cos_phi;
165
166    // reverse view matrix
167    SINCOS(theta*M_PI_180, sin_theta, cos_theta);
168    SINCOS(phi*M_PI_180, sin_phi, cos_phi);
169    const double V11 = cos_phi*cos_theta;
170    const double V12 = sin_phi*cos_theta;
171    const double V21 = -sin_phi;
172    const double V22 = cos_phi;
173    const double V31 = sin_theta*cos_phi;
174    const double V32 = sin_phi*sin_theta;
175
176    // reverse jitter matrix
177    SINCOS(dtheta*M_PI_180, sin_theta, cos_theta);
178    SINCOS(dphi*M_PI_180, sin_phi, cos_phi);
179    const double J31 = sin_theta;
180    const double J32 = -sin_phi*cos_theta;
181    const double J33 = cos_phi*cos_theta;
182
183    // reverse matrix
184    rotation->R31 = J31*V11 + J32*V21 + J33*V31;
185    rotation->R32 = J31*V12 + J32*V22 + J33*V32;
186}
187
188// Apply the rotation matrix returned from qac_rotation to the point (qx,qy),
189// returning R*[qx,qy]' = [qa,qc]'
190static void
191qac_apply(
192    QACRotation *rotation,
193    double qx, double qy,
194    double *qab_out, double *qc_out)
195{
196    // Indirect calculation of qab, from qab^2 = |q|^2 - qc^2
197    const double dqc = rotation->R31*qx + rotation->R32*qy;
198    const double dqab_sq = -dqc*dqc + qx*qx + qy*qy;
199    //*qab_out = sqrt(fabs(dqab_sq));
200    *qab_out = dqab_sq > 0.0 ? sqrt(dqab_sq) : 0.0;
201    *qc_out = dqc;
202}
203#endif // _QAC_SECTION
204
205#if !defined(_QABC_SECTION) && defined(CALL_IQ_ABC)
206#define _QABC_SECTION
207
208typedef struct {
209    double R11, R12;
210    double R21, R22;
211    double R31, R32;
212} QABCRotation;
213
214// Fill in the rotation matrix R from the view angles (theta, phi, psi) and the
215// jitter angles (dtheta, dphi, dpsi).  This matrix can be applied to all of the
216// (qx, qy) points in the image to produce R*[qx,qy]' = [qa,qb,qc]'
217static void
218qabc_rotation(
219    QABCRotation *rotation,
220    double theta, double phi, double psi,
221    double dtheta, double dphi, double dpsi)
222{
223    double sin_theta, cos_theta;
224    double sin_phi, cos_phi;
225    double sin_psi, cos_psi;
226
227    // reverse view matrix
228    SINCOS(theta*M_PI_180, sin_theta, cos_theta);
229    SINCOS(phi*M_PI_180, sin_phi, cos_phi);
230    SINCOS(psi*M_PI_180, sin_psi, cos_psi);
231    const double V11 = -sin_phi*sin_psi + cos_phi*cos_psi*cos_theta;
232    const double V12 = sin_phi*cos_psi*cos_theta + sin_psi*cos_phi;
233    const double V21 = -sin_phi*cos_psi - sin_psi*cos_phi*cos_theta;
234    const double V22 = -sin_phi*sin_psi*cos_theta + cos_phi*cos_psi;
235    const double V31 = sin_theta*cos_phi;
236    const double V32 = sin_phi*sin_theta;
237
238    // reverse jitter matrix
239    SINCOS(dtheta*M_PI_180, sin_theta, cos_theta);
240    SINCOS(dphi*M_PI_180, sin_phi, cos_phi);
241    SINCOS(dpsi*M_PI_180, sin_psi, cos_psi);
242    const double J11 = cos_psi*cos_theta;
243    const double J12 = sin_phi*sin_theta*cos_psi + sin_psi*cos_phi;
244    const double J13 = sin_phi*sin_psi - sin_theta*cos_phi*cos_psi;
245    const double J21 = -sin_psi*cos_theta;
246    const double J22 = -sin_phi*sin_psi*sin_theta + cos_phi*cos_psi;
247    const double J23 = sin_phi*cos_psi + sin_psi*sin_theta*cos_phi;
248    const double J31 = sin_theta;
249    const double J32 = -sin_phi*cos_theta;
250    const double J33 = cos_phi*cos_theta;
251
252    // reverse matrix
253    rotation->R11 = J11*V11 + J12*V21 + J13*V31;
254    rotation->R12 = J11*V12 + J12*V22 + J13*V32;
255    rotation->R21 = J21*V11 + J22*V21 + J23*V31;
256    rotation->R22 = J21*V12 + J22*V22 + J23*V32;
257    rotation->R31 = J31*V11 + J32*V21 + J33*V31;
258    rotation->R32 = J31*V12 + J32*V22 + J33*V32;
259}
260
261// Apply the rotation matrix returned from qabc_rotation to the point (qx,qy),
262// returning R*[qx,qy]' = [qa,qb,qc]'
263static void
264qabc_apply(
265    QABCRotation *rotation,
266    double qx, double qy,
267    double *qa_out, double *qb_out, double *qc_out)
268{
269    *qa_out = rotation->R11*qx + rotation->R12*qy;
270    *qb_out = rotation->R21*qx + rotation->R22*qy;
271    *qc_out = rotation->R31*qx + rotation->R32*qy;
272}
273
274#endif // _QABC_SECTION
275
276// ==================== KERNEL CODE ========================
277kernel
278void KERNEL_NAME(
279    int32_t nq,                 // number of q values
280    const int32_t pd_start,     // where we are in the dispersity loop
281    const int32_t pd_stop,      // where we are stopping in the dispersity loop
282    global const ProblemDetails *details,
283    global const double *values,
284    global const double *q, // nq q values, with padding to boundary
285    global double *result,  // nq+1 return values, again with padding
286    const double cutoff,     // cutoff in the dispersity weight product
287    int32_t effective_radius_type // which effective radius to compute
288    )
289{
290#ifdef USE_OPENCL
291  // who we are and what element we are working with
292  const int q_index = get_global_id(0);
293  if (q_index >= nq) return;
294#else
295  // Define q_index here so that debugging statements can be written to work
296  // for both OpenCL and DLL using:
297  //    if (q_index == 0) {printf(...);}
298  int q_index = 0;
299#endif
300
301  // ** Fill in the local values table **
302  // Storage for the current parameter values.
303  // These will be updated as we walk the dispersity mesh.
304  ParameterBlock local_values;
305  //   values[0] is scale
306  //   values[1] is background
307  #ifdef USE_OPENMP
308  #pragma omp parallel for
309  #endif
310  for (int i=0; i < NUM_PARS; i++) {
311    local_values.vector[i] = values[2+i];
312    //if (q_index==0) printf("p%d = %g\n",i, local_values.vector[i]);
313  }
314  //if (q_index==0) printf("NUM_VALUES:%d  NUM_PARS:%d  MAX_PD:%d\n", NUM_VALUES, NUM_PARS, MAX_PD);
315  //if (q_index==0) printf("start:%d  stop:%d\n", pd_start, pd_stop);
316
317  // ** Precompute magnatism values **
318#if defined(MAGNETIC) && NUM_MAGNETIC>0
319  // Location of the sld parameters in the parameter vector.
320  // These parameters are updated with the effective sld due to magnetism.
321  const int32_t slds[] = { MAGNETIC_PARS };
322
323  // Interpret polarization cross section.
324  //     up_frac_i = values[NUM_PARS+2];
325  //     up_frac_f = values[NUM_PARS+3];
326  //     up_angle = values[NUM_PARS+4];
327  // TODO: could precompute more magnetism parameters before calling the kernel.
328  double xs_weights[8];  // uu, ud real, du real, dd, ud imag, du imag, fill, fill
329  double cos_mspin, sin_mspin;
330  set_spin_weights(values[NUM_PARS+2], values[NUM_PARS+3], xs_weights);
331  SINCOS(-values[NUM_PARS+4]*M_PI_180, sin_mspin, cos_mspin);
332#endif // MAGNETIC
333
334  // ** Fill in the initial results **
335  // If pd_start is zero that means that we are starting a new calculation,
336  // and must initialize the result to zero.  Otherwise, we are restarting
337  // the calculation from somewhere in the middle of the dispersity mesh,
338  // and we update the value rather than reset it. Similarly for the
339  // normalization factor, which is stored as the final value in the
340  // results vector (one past the number of q values).
341  //
342  // The code differs slightly between opencl and dll since opencl is only
343  // seeing one q value (stored in the variable "this_result") while the dll
344  // version must loop over all q.
345  #ifdef USE_OPENCL
346    #if defined(CALL_FQ)
347      double weight_norm = (pd_start == 0 ? 0.0 : result[2*nq]);
348      double weighted_form = (pd_start == 0 ? 0.0 : result[2*nq+1]);
349      double weighted_shell = (pd_start == 0 ? 0.0 : result[2*nq+2]);
350      double weighted_radius = (pd_start == 0 ? 0.0 : result[2*nq+3]);
351      double this_F2 = (pd_start == 0 ? 0.0 : result[2*q_index+0]);
352      double this_F1 = (pd_start == 0 ? 0.0 : result[2*q_index+1]);
353    #else
354      double weight_norm = (pd_start == 0 ? 0.0 : result[nq]);
355      double weighted_form = (pd_start == 0 ? 0.0 : result[nq+1]);
356      double weighted_shell = (pd_start == 0 ? 0.0 : result[nq+2]);
357      double weighted_radius = (pd_start == 0 ? 0.0 : result[nq+3]);
358      double this_result = (pd_start == 0 ? 0.0 : result[q_index]);
359    #endif
360  #else // !USE_OPENCL
361    #if defined(CALL_FQ)
362      double weight_norm = (pd_start == 0 ? 0.0 : result[2*nq]);
363      double weighted_form = (pd_start == 0 ? 0.0 : result[2*nq+1]);
364      double weighted_shell = (pd_start == 0 ? 0.0 : result[2*nq+2]);
365      double weighted_radius = (pd_start == 0 ? 0.0 : result[2*nq+3]);
366    #else
367      double weight_norm = (pd_start == 0 ? 0.0 : result[nq]);
368      double weighted_form = (pd_start == 0 ? 0.0 : result[nq+1]);
369      double weighted_shell = (pd_start == 0 ? 0.0 : result[nq+2]);
370      double weighted_radius = (pd_start == 0 ? 0.0 : result[nq+3]);
371    #endif
372    if (pd_start == 0) {
373      #ifdef USE_OPENMP
374      #pragma omp parallel for
375      #endif
376      #if defined(CALL_FQ)
377          // 2*nq for F^2,F pairs
378          for (int q_index=0; q_index < 2*nq; q_index++) result[q_index] = 0.0;
379      #else
380          for (int q_index=0; q_index < nq; q_index++) result[q_index] = 0.0;
381      #endif
382    }
383    //if (q_index==0) printf("start %d %g %g\n", pd_start, weighted_shell, result[0]);
384#endif // !USE_OPENCL
385
386
387// ====== macros to set up the parts of the loop =======
388/*
389Based on the level of the loop, uses C preprocessor magic to construct
390level-specific looping variables, including these from loop level 3:
391
392  int n3 : length of loop for mesh level 3
393  int i3 : current position in the loop for level 3, which is calculated
394       from a combination of pd_start, pd_stride[3] and pd_length[3].
395  int p3 : is the index into the parameter table for mesh level 3
396  double v3[] : pointer into dispersity array to values for loop 3
397  double w3[] : pointer into dispersity array to weights for loop 3
398  double weight3 : the product of weights from levels 3 and up, computed
399       as weight5*weight4*w3[i3].  Note that we need an outermost
400       value weight5 set to 1.0 for this to work properly.
401
402After expansion, the loop struction will look like the following:
403
404  // --- PD_INIT(4) ---
405  const int n4 = pd_length[4];
406  const int p4 = pd_par[4];
407  global const double *v4 = pd_value + pd_offset[4];
408  global const double *w4 = pd_weight + pd_offset[4];
409  int i4 = (pd_start/pd_stride[4])%n4;  // position in level 4 at pd_start
410
411  // --- PD_INIT(3) ---
412  const int n3 = pd_length[3];
413  ...
414  int i3 = (pd_start/pd_stride[3])%n3;  // position in level 3 at pd_start
415
416  PD_INIT(2)
417  PD_INIT(1)
418  PD_INIT(0)
419
420  // --- PD_OUTERMOST_WEIGHT(5) ---
421  const double weight5 = 1.0;
422
423  // --- PD_OPEN(4,5) ---
424  while (i4 < n4) {
425    parameter[p4] = v4[i4];  // set the value for pd parameter 4 at this mesh point
426    const double weight4 = w4[i4] * weight5;
427
428    // from PD_OPEN(3,4)
429    while (i3 < n3) {
430      parameter[p3] = v3[i3];  // set the value for pd parameter 3 at this mesh point
431      const double weight3 = w3[i3] * weight4;
432
433      PD_OPEN(3,2)
434      PD_OPEN(2,1)
435      PD_OPEN(0,1)
436
437      // ... main loop body ...
438      APPLY_PROJECTION    // convert jitter values to spherical coords
439      BUILD_ROTATION      // construct the rotation matrix qxy => qabc
440      for each q
441          FETCH_Q         // set qx,qy from the q input vector
442          APPLY_ROTATION  // convert qx,qy to qa,qb,qc
443          CALL_KERNEL     // scattering = Iqxy(qa, qb, qc, p1, p2, ...)
444
445      ++step;  // increment counter representing position in dispersity mesh
446
447      PD_CLOSE(0)
448      PD_CLOSE(1)
449      PD_CLOSE(2)
450
451      // --- PD_CLOSE(3) ---
452      if (step >= pd_stop) break;
453      ++i3;
454    }
455    i3 = 0; // reset loop counter for next round through the loop
456
457    // --- PD_CLOSE(4) ---
458    if (step >= pd_stop) break;
459    ++i4;
460  }
461  i4 = 0; // reset loop counter even though no more rounds through the loop
462
463*/
464
465
466// ** prepare inner loops **
467
468// Depending on the shape type (radial, axial, triaxial), the variables
469// and calling parameters in the loop body will be slightly different.
470// Macros capture the differences in one spot so the rest of the code
471// is easier to read. The code below both declares variables for the
472// inner loop and defines the macros that use them.
473
474
475#if defined(CALL_FQ) // COMPUTE_F1_F2 is true
476  // unoriented 1D returning <F> and <F^2>
477  double qk;
478  double F1, F2;
479  #define FETCH_Q() do { qk = q[q_index]; } while (0)
480  #define BUILD_ROTATION() do {} while(0)
481  #define APPLY_ROTATION() do {} while(0)
482  #define CALL_KERNEL() CALL_FQ(qk,F1,F2,local_values.table)
483
484#elif defined(CALL_FQ_A)
485  // unoriented 2D return <F> and <F^2>
486  double qx, qy;
487  double F1, F2;
488  #define FETCH_Q() do { qx = q[2*q_index]; qy = q[2*q_index+1]; } while (0)
489  #define BUILD_ROTATION() do {} while(0)
490  #define APPLY_ROTATION() do {} while(0)
491  #define CALL_KERNEL() CALL_FQ_A(sqrt(qx*qx+qy*qy),F1,F2,local_values.table)
492
493#elif defined(CALL_IQ)
494  // unoriented 1D return <F^2>
495  double qk;
496  #define FETCH_Q() do { qk = q[q_index]; } while (0)
497  #define BUILD_ROTATION() do {} while(0)
498  #define APPLY_ROTATION() do {} while(0)
499  #define CALL_KERNEL() CALL_IQ(qk,local_values.table)
500
501#elif defined(CALL_IQ_A)
502  // unoriented 2D
503  double qx, qy;
504  #define FETCH_Q() do { qx = q[2*q_index]; qy = q[2*q_index+1]; } while (0)
505  #define BUILD_ROTATION() do {} while(0)
506  #define APPLY_ROTATION() do {} while(0)
507  #define CALL_KERNEL() CALL_IQ_A(sqrt(qx*qx+qy*qy), local_values.table)
508
509#elif defined(CALL_IQ_AC)
510  // oriented symmetric 2D
511  double qx, qy;
512  #define FETCH_Q() do { qx = q[2*q_index]; qy = q[2*q_index+1]; } while (0)
513  double qa, qc;
514  QACRotation rotation;
515  // theta, phi, dtheta, dphi are defined below in projection to avoid repeated code.
516  #define BUILD_ROTATION() qac_rotation(&rotation, theta, phi, dtheta, dphi);
517  #define APPLY_ROTATION() qac_apply(&rotation, qx, qy, &qa, &qc)
518  #define CALL_KERNEL() CALL_IQ_AC(qa, qc, local_values.table)
519
520#elif defined(CALL_IQ_ABC)
521  // oriented asymmetric 2D
522  double qx, qy;
523  #define FETCH_Q() do { qx = q[2*q_index]; qy = q[2*q_index+1]; } while (0)
524  double qa, qb, qc;
525  QABCRotation rotation;
526  // theta, phi, dtheta, dphi are defined below in projection to avoid repeated code.
527  // psi and dpsi are only for IQ_ABC, so they are processed here.
528  const double psi = values[details->theta_par+4];
529  local_values.table.psi = 0.;
530  #define BUILD_ROTATION() qabc_rotation(&rotation, theta, phi, psi, dtheta, dphi, local_values.table.psi)
531  #define APPLY_ROTATION() qabc_apply(&rotation, qx, qy, &qa, &qb, &qc)
532  #define CALL_KERNEL() CALL_IQ_ABC(qa, qb, qc, local_values.table)
533
534#elif defined(CALL_IQ_XY)
535  // direct call to qx,qy calculator
536  double qx, qy;
537  #define FETCH_Q() do { qx = q[2*q_index]; qy = q[2*q_index+1]; } while (0)
538  #define BUILD_ROTATION() do {} while(0)
539  #define APPLY_ROTATION() do {} while(0)
540  #define CALL_KERNEL() CALL_IQ_XY(qx, qy, local_values.table)
541#endif
542
543// Define APPLY_PROJECTION depending on model symmetries. We do this outside
544// the previous if block so that we don't need to repeat the identical
545// logic in the IQ_AC and IQ_ABC branches.  This will become more important
546// if we implement more projections, or more complicated projections.
547#if defined(CALL_IQ) || defined(CALL_IQ_A) || defined(CALL_FQ) || defined(CALL_FQ_A)
548  // no orientation
549  #define APPLY_PROJECTION() const double weight=weight0
550#elif defined(CALL_IQ_XY) // pass orientation to the model
551  // CRUFT: support oriented model which define Iqxy rather than Iqac or Iqabc
552  // Need to plug the values for the orientation angles back into parameter
553  // table in case they were overridden by the orientation offset.  This
554  // means that orientation dispersity will not work for these models, but
555  // it was broken anyway, so no matter.  Still want to provide Iqxy in case
556  // the user model wants full control of orientation/magnetism.
557  #if defined(HAVE_PSI)
558    const double theta = values[details->theta_par+2];
559    const double phi = values[details->theta_par+3];
560    const double psi = values[details->theta_par+4];
561    double weight;
562    #define APPLY_PROJECTION() do { \
563      local_values.table.theta = theta; \
564      local_values.table.phi = phi; \
565      local_values.table.psi = psi; \
566      weight=weight0; \
567    } while (0)
568  #elif defined(HAVE_THETA)
569    const double theta = values[details->theta_par+2];
570    const double phi = values[details->theta_par+3];
571    double weight;
572    #define APPLY_PROJECTION() do { \
573      local_values.table.theta = theta; \
574      local_values.table.phi = phi; \
575      weight=weight0; \
576    } while (0)
577  #else
578    #define APPLY_PROJECTION() const double weight=weight0
579  #endif
580#else // apply jitter and view before calling the model
581  // Grab the "view" angles (theta, phi, psi) from the initial parameter table.
582  const double theta = values[details->theta_par+2];
583  const double phi = values[details->theta_par+3];
584  // Make sure jitter angle defaults to zero if there is no jitter distribution
585  local_values.table.theta = 0.;
586  local_values.table.phi = 0.;
587  // The "jitter" angles (dtheta, dphi, dpsi) are stored with the
588  // dispersity values and copied to the local parameter table as
589  // we go through the mesh.
590  double dtheta, dphi, weight;
591  #if PROJECTION == 1 // equirectangular
592    #define APPLY_PROJECTION() do { \
593      dtheta = local_values.table.theta; \
594      dphi = local_values.table.phi; \
595      weight = fabs(cos(dtheta*M_PI_180)) * weight0; \
596    } while (0)
597  #elif PROJECTION == 2 // sinusoidal
598    #define APPLY_PROJECTION() do { \
599      dtheta = local_values.table.theta; \
600      dphi = local_values.table.phi; \
601      weight = weight0; \
602      if (dtheta != 90.0) dphi /= cos(dtheta*M_PI_180); \
603      else if (dphi != 0.0) weight = 0.; \
604      if (fabs(dphi) >= 180.) weight = 0.; \
605    } while (0)
606  #endif
607#endif // done defining APPLY_PROJECTION
608
609// ** define looping macros **
610
611// Define looping variables
612#define PD_INIT(_LOOP) \
613  const int n##_LOOP = details->pd_length[_LOOP]; \
614  const int p##_LOOP = details->pd_par[_LOOP]; \
615  global const double *v##_LOOP = pd_value + details->pd_offset[_LOOP]; \
616  global const double *w##_LOOP = pd_weight + details->pd_offset[_LOOP]; \
617  int i##_LOOP = (pd_start/details->pd_stride[_LOOP])%n##_LOOP;
618
619// Jump into the middle of the dispersity loop
620#define PD_OPEN(_LOOP,_OUTER) \
621  while (i##_LOOP < n##_LOOP) { \
622    local_values.vector[p##_LOOP] = v##_LOOP[i##_LOOP]; \
623    const double weight##_LOOP = w##_LOOP[i##_LOOP] * weight##_OUTER;
624
625// create the variable "weight#=1.0" where # is the outermost level+1 (=MAX_PD).
626#define _PD_OUTERMOST_WEIGHT(_n) const double weight##_n = 1.0;
627#define PD_OUTERMOST_WEIGHT(_n) _PD_OUTERMOST_WEIGHT(_n)
628
629// Close out the loop
630#define PD_CLOSE(_LOOP) \
631    if (step >= pd_stop) break; \
632    ++i##_LOOP; \
633  } \
634  i##_LOOP = 0;
635
636// ====== construct the loops =======
637
638// Pointers to the start of the dispersity and weight vectors, if needed.
639#if MAX_PD>0
640  global const double *pd_value = values + NUM_VALUES;
641  global const double *pd_weight = pd_value + details->num_weights;
642#endif
643
644// The variable "step" is the current position in the dispersity loop.
645// It will be incremented each time a new point in the mesh is accumulated,
646// and used to test whether we have reached pd_stop.
647int step = pd_start;
648
649// *** define loops for each of 0, 1, 2, ..., modelinfo.MAX_PD-1 ***
650
651// define looping variables
652#if MAX_PD>4
653  PD_INIT(4)
654#endif
655#if MAX_PD>3
656  PD_INIT(3)
657#endif
658#if MAX_PD>2
659  PD_INIT(2)
660#endif
661#if MAX_PD>1
662  PD_INIT(1)
663#endif
664#if MAX_PD>0
665  PD_INIT(0)
666#endif
667
668// open nested loops
669PD_OUTERMOST_WEIGHT(MAX_PD)
670#if MAX_PD>4
671  PD_OPEN(4,5)
672#endif
673#if MAX_PD>3
674  PD_OPEN(3,4)
675#endif
676#if MAX_PD>2
677  PD_OPEN(2,3)
678#endif
679#if MAX_PD>1
680  PD_OPEN(1,2)
681#endif
682#if MAX_PD>0
683  PD_OPEN(0,1)
684#endif
685
686//if (q_index==0) {printf("step:%d of %d, pars:",step,pd_stop); for (int i=0; i < NUM_PARS; i++) printf("p%d=%g ",i, local_values.vector[i]); printf("\n");}
687
688  // ====== loop body =======
689  #ifdef INVALID
690  if (!INVALID(local_values.table))
691  #endif
692  {
693     APPLY_PROJECTION();
694
695    // Accumulate I(q)
696    // Note: weight==0 must always be excluded
697    if (weight > cutoff) {
698      double form, shell;
699      CALL_VOLUME(form, shell, local_values.table);
700      weight_norm += weight;
701      weighted_form += weight * form;
702      weighted_shell += weight * shell;
703      if (effective_radius_type != 0) {
704        weighted_radius += weight * CALL_EFFECTIVE_RADIUS(effective_radius_type, local_values.table);
705      }
706      BUILD_ROTATION();
707
708#ifndef USE_OPENCL
709      // DLL needs to explicitly loop over the q values.
710      #ifdef USE_OPENMP
711      #pragma omp parallel for
712      #endif
713      for (q_index=0; q_index<nq; q_index++)
714#endif // !USE_OPENCL
715      {
716
717        FETCH_Q();
718        APPLY_ROTATION();
719
720        // ======= COMPUTE SCATTERING ==========
721        #if defined(MAGNETIC) && NUM_MAGNETIC > 0
722          // Compute the scattering from the magnetic cross sections.
723          double scattering = 0.0;
724          const double qsq = qx*qx + qy*qy;
725          if (qsq > 1.e-16) {
726            // TODO: what is the magnetic scattering at q=0
727            const double px = (qy*cos_mspin + qx*sin_mspin)/qsq;
728            const double py = (qy*sin_mspin - qx*cos_mspin)/qsq;
729
730            // loop over uu, ud real, du real, dd, ud imag, du imag
731            for (unsigned int xs=0; xs<6; xs++) {
732              const double xs_weight = xs_weights[xs];
733              if (xs_weight > 1.e-8) {
734                // Since the cross section weight is significant, set the slds
735                // to the effective slds for this cross section, call the
736                // kernel, and add according to weight.
737                for (int sk=0; sk<NUM_MAGNETIC; sk++) {
738                  const int32_t mag_index = NUM_PARS+5 + 3*sk;
739                  const int32_t sld_index = slds[sk];
740                  const double mx = values[mag_index];
741                  const double my = values[mag_index+1];
742                  const double mz = values[mag_index+2];
743                  local_values.vector[sld_index] =
744                    mag_sld(xs, qx, qy, px, py, values[sld_index+2], mx, my, mz);
745//if (q_index==0) printf("%d: (qx,qy)=(%g,%g) xs=%d sld%d=%g p=(%g,%g) m=(%g,%g,%g)\n",
746//  q_index, qx, qy, xs, sk, local_values.vector[sld_index], px, py, mx, my, mz);
747                }
748                scattering += xs_weight * CALL_KERNEL();
749              }
750            }
751          }
752        #else  // !MAGNETIC
753          #if defined(CALL_FQ)
754            CALL_KERNEL(); // sets F1 and F2 by reference
755          #else
756            const double scattering = CALL_KERNEL();
757          #endif
758        #endif // !MAGNETIC
759//printf("q_index:%d %g %g %g %g\n", q_index, scattering, weight0);
760
761        #ifdef USE_OPENCL
762          #if defined(CALL_FQ)
763            this_F2 += weight * F2;
764            this_F1 += weight * F1;
765          #else
766            this_result += weight * scattering;
767          #endif
768        #else // !USE_OPENCL
769          #if defined(CALL_FQ)
770            result[2*q_index+0] += weight * F2;
771            result[2*q_index+1] += weight * F1;
772          #else
773            result[q_index] += weight * scattering;
774          #endif
775        #endif // !USE_OPENCL
776      }
777    }
778  }
779// close nested loops
780++step;
781#if MAX_PD>0
782  PD_CLOSE(0)
783#endif
784#if MAX_PD>1
785  PD_CLOSE(1)
786#endif
787#if MAX_PD>2
788  PD_CLOSE(2)
789#endif
790#if MAX_PD>3
791  PD_CLOSE(3)
792#endif
793#if MAX_PD>4
794  PD_CLOSE(4)
795#endif
796
797// Remember the current result and the updated norm.
798#ifdef USE_OPENCL
799  #if defined(CALL_FQ)
800    result[2*q_index+0] = this_F2;
801    result[2*q_index+1] = this_F1;
802    if (q_index == 0) {
803      result[2*nq+0] = weight_norm;
804      result[2*nq+1] = weighted_form;
805      result[2*nq+3] = weighted_shell;
806      result[2*nq+3] = weighted_radius;
807    }
808  #else
809    result[q_index] = this_result;
810    if (q_index == 0) {
811      result[nq+0] = weight_norm;
812      result[nq+1] = weighted_form;
813      result[nq+2] = weighted_shell;
814      result[nq+3] = weighted_radius;
815    }
816  #endif
817
818//if (q_index == 0) printf("res: %g/%g\n", result[0], weighted_shell);
819#else // !USE_OPENCL
820  #if defined(CALL_FQ)
821    result[2*nq] = weight_norm;
822    result[2*nq+1] = weighted_form;
823    result[2*nq+2] = weighted_shell;
824    result[2*nq+3] = weighted_radius;
825  #else
826    result[nq] = weight_norm;
827    result[nq+1] = weighted_form;
828    result[nq+2] = weighted_shell;
829    result[nq+3] = weighted_radius;
830  #endif
831//printf("res: %g/%g\n", result[0], weighted_shell);
832#endif // !USE_OPENCL
833
834// ** clear the macros in preparation for the next kernel **
835#undef PD_INIT
836#undef PD_OPEN
837#undef PD_CLOSE
838#undef FETCH_Q
839#undef APPLY_PROJECTION
840#undef BUILD_ROTATION
841#undef APPLY_ROTATION
842#undef CALL_KERNEL
843}
Note: See TracBrowser for help on using the repository browser.