source: sasmodels/sasmodels/generate.py @ b3f6bc3

core_shell_microgelscostrafo411magnetic_modelrelease_v0.94release_v0.95ticket-1257-vesicle-productticket_1156ticket_1265_superballticket_822_more_unit_tests
Last change on this file since b3f6bc3 was b3f6bc3, checked in by pkienzle, 9 years ago

support pure python model Iq/Iqxy?

  • Property mode set to 100644
File size: 26.3 KB
Line 
1"""
2SAS model constructor.
3
4Small angle scattering models are defined by a set of kernel functions:
5
6    *Iq(q, p1, p2, ...)* returns the scattering at q for a form with
7    particular dimensions averaged over all orientations.
8
9    *Iqxy(qx, qy, p1, p2, ...)* returns the scattering at qx,qy for a form
10    with particular dimensions for a single orientation.
11
12    *Imagnetic(qx, qy, result[], p1, p2, ...)* returns the scattering for the
13    polarized neutron spin states (up-up, up-down, down-up, down-down) for
14    a form with particular dimensions for a single orientation.
15
16    *form_volume(p1, p2, ...)* returns the volume of the form with particular
17    dimension.
18
19    *ER(p1, p2, ...)* returns the effective radius of the form with
20    particular dimensions.
21
22    *VR(p1, p2, ...)* returns the volume ratio for core-shell style forms.
23
24These functions are defined in a kernel module .py script and an associated
25set of .c files.  The model constructor will use them to create models with
26polydispersity across volume and orientation parameters, and provide
27scale and background parameters for each model.
28
29*Iq*, *Iqxy*, *Imagnetic* and *form_volume* should be stylized C-99
30functions written for OpenCL.  All functions need prototype declarations
31even if the are defined before they are used.  OpenCL does not support
32*#include* preprocessor directives, so instead the list of includes needs
33to be given as part of the metadata in the kernel module definition.
34The included files should be listed using a path relative to the kernel
35module, or if using "lib/file.c" if it is one of the standard includes
36provided with the sasmodels source.  The includes need to be listed in
37order so that functions are defined before they are used.
38
39Floating point values should be declared as *double*.  For single precision
40calculations, *double* will be replaced by *float*.  The single precision
41conversion will also tag floating point constants with "f" to make them
42single precision constants.  When using integral values in floating point
43expressions, they should be expressed as floating point values by including
44a decimal point.  This includes 0., 1. and 2.
45
46OpenCL has a *sincos* function which can improve performance when both
47the *sin* and *cos* values are needed for a particular argument.  Since
48this function does not exist in C99, all use of *sincos* should be
49replaced by the macro *SINCOS(value,sn,cn)* where *sn* and *cn* are
50previously declared *double* variables.  When compiled for systems without
51OpenCL, *SINCOS* will be replaced by *sin* and *cos* calls.   If *value* is
52an expression, it will appear twice in this case; whether or not it will be
53evaluated twice depends on the quality of the compiler.
54
55If the input parameters are invalid, the scattering calculator should
56return a negative number. Particularly with polydispersity, there are
57some sets of shape parameters which lead to nonsensical forms, such
58as a capped cylinder where the cap radius is smaller than the
59cylinder radius.  The polydispersity calculation will ignore these points,
60effectively chopping the parameter weight distributions at the boundary
61of the infeasible region.  The resulting scattering will be set to
62background.  This will work correctly even when polydispersity is off.
63
64*ER* and *VR* are python functions which operate on parameter vectors.
65The constructor code will generate the necessary vectors for computing
66them with the desired polydispersity.
67
68The available kernel parameters are defined as a list, with each parameter
69defined as a sublist with the following elements:
70
71    *name* is the name that will be used in the call to the kernel
72    function and the name that will be displayed to the user.  Names
73    should be lower case, with words separated by underscore.  If
74    acronyms are used, the whole acronym should be upper case.
75
76    *units* should be one of *degrees* for angles, *Ang* for lengths,
77    *1e-6/Ang^2* for SLDs.
78
79    *default value* will be the initial value for  the model when it
80    is selected, or when an initial value is not otherwise specified.
81
82    [*lb*, *ub*] are the hard limits on the parameter value, used to limit
83    the polydispersity density function.  In the fit, the parameter limits
84    given to the fit are the limits  on the central value of the parameter.
85    If there is polydispersity, it will evaluate parameter values outside
86    the fit limits, but not outside the hard limits specified in the model.
87    If there are no limits, use +/-inf imported from numpy.
88
89    *type* indicates how the parameter will be used.  "volume" parameters
90    will be used in all functions.  "orientation" parameters will be used
91    in *Iqxy* and *Imagnetic*.  "magnetic* parameters will be used in
92    *Imagnetic* only.  If *type* is the empty string, the parameter will
93    be used in all of *Iq*, *Iqxy* and *Imagnetic*.
94
95    *description* is a short description of the parameter.  This will
96    be displayed in the parameter table and used as a tool tip for the
97    parameter value in the user interface.
98
99The kernel module must set variables defining the kernel meta data:
100
101    *name* is the model name
102
103    *title* is a short description of the model, suitable for a tool tip,
104    or a one line model summary in a table of models.
105
106    *description* is an extended description of the model to be displayed
107    while the model parameters are being edited.
108
109    *parameters* is the list of parameters.  Parameters in the kernel
110    functions must appear in the same order as they appear in the
111    parameters list.  Two additional parameters, *scale* and *background*
112    are added to the beginning of the parameter list.  They will show up
113    in the documentation as model parameters, but they are never sent to
114    the kernel functions.
115
116    *source* is the list of C-99 source files that must be joined to
117    create the OpenCL kernel functions.  The files defining the functions
118    need to be listed before the files which use the functions.
119
120    *ER* is a python function defining the effective radius.  If it is
121    not present, the effective radius is 0.
122
123    *VR* is a python function defining the volume ratio.  If it is not
124    present, the volume ratio is 1.
125
126    *form_volume*, *Iq*, *Iqxy*, *Imagnetic* are strings containing the
127    C source code for the body of the volume, Iq, and Iqxy functions
128    respectively.  These can also be defined in the last source file.
129
130An *info* dictionary is constructed from the kernel meta data and
131returned to the caller.  It includes the additional fields
132
133
134The model evaluator, function call sequence consists of q inputs and the return vector,
135followed by the loop value/weight vector, followed by the values for
136the non-polydisperse parameters, followed by the lengths of the
137polydispersity loops.  To construct the call for 1D models, the
138categories *fixed-1d* and *pd-1d* list the names of the parameters
139of the non-polydisperse and the polydisperse parameters respectively.
140Similarly, *fixed-2d* and *pd-2d* provide parameter names for 2D models.
141The *pd-rel* category is a set of those parameters which give
142polydispersitiy as a portion of the value (so a 10% length dispersity
143would use a polydispersity value of 0.1) rather than absolute
144dispersity such as an angle plus or minus 15 degrees.
145
146The *volume* category lists the volume parameters in order for calls
147to volume within the kernel (used for volume normalization) and for
148calls to ER and VR for effective radius and volume ratio respectively.
149
150The *orientation* and *magnetic* categories list the orientation and
151magnetic parameters.  These are used by the sasview interface.  The
152blank category is for parameters such as scale which don't have any
153other marking.
154
155The doc string at the start of the kernel module will be used to
156construct the model documentation web pages.  Embedded figures should
157appear in the subdirectory "img" beside the model definition, and tagged
158with the kernel module name to avoid collision with other models.  Some
159file systems are case-sensitive, so only use lower case characters for
160file names and extensions.
161
162
163The function :func:`make` loads the metadata from the module and returns
164the kernel source.  The function :func:`doc` extracts the doc string
165and adds the parameter table to the top.  The function :func:`sources`
166returns a list of files required by the model.
167"""
168
169# TODO: identify model files which have changed since loading and reload them.
170
171__all__ = ["make, doc", "sources", "use_single"]
172
173import sys
174import os
175import os.path
176import re
177
178import numpy as np
179
180F64 = np.dtype('float64')
181F32 = np.dtype('float32')
182
183# Scale and background, which are parameters common to every form factor
184COMMON_PARAMETERS = [
185    [ "scale", "", 1, [0, np.inf], "", "Source intensity" ],
186    [ "background", "1/cm", 0, [0, np.inf], "", "Source background" ],
187    ]
188
189
190# Conversion from units defined in the parameter table for each model
191# to units displayed in the sphinx documentation.
192RST_UNITS = {
193    "Ang": "|Ang|",
194    "1/Ang^2": "|Ang^-2|",
195    "1e-6/Ang^2": "|1e-6Ang^-2|",
196    "degrees": "degree",
197    "1/cm": "|cm^-1|",
198    "": "None",
199    }
200
201# Headers for the parameters tables in th sphinx documentation
202PARTABLE_HEADERS = [
203    "Parameter",
204    "Description",
205    "Units",
206    "Default value",
207    ]
208
209# Minimum width for a default value (this is shorter than the column header
210# width, so will be ignored).
211PARTABLE_VALUE_WIDTH = 10
212
213# Header included before every kernel.
214# This makes sure that the appropriate math constants are defined, and does
215# whatever is required to make the kernel compile as pure C rather than
216# as an OpenCL kernel.
217KERNEL_HEADER = """\
218// GENERATED CODE --- DO NOT EDIT ---
219// Code is produced by sasmodels.gen from sasmodels/models/MODEL.c
220
221#ifdef __OPENCL_VERSION__
222# define USE_OPENCL
223#endif
224
225// If opencl is not available, then we are compiling a C function
226// Note: if using a C++ compiler, then define kernel as extern "C"
227#ifndef USE_OPENCL
228#  ifdef __cplusplus
229     #include <cmath>
230     #if defined(_MSC_VER)
231     #define kernel extern "C" __declspec( dllexport )
232     #else
233     #define kernel extern "C"
234     #endif
235     using namespace std;
236     inline void SINCOS(double angle, double &svar, double &cvar)
237       { svar=sin(angle); cvar=cos(angle); }
238#  else
239     #include <math.h>
240     #if defined(_MSC_VER)
241     #define kernel __declspec( dllexport )
242     #else
243     #define kernel
244     #endif
245     #define SINCOS(angle,svar,cvar) do {svar=sin(angle);cvar=cos(angle);} while (0)
246#  endif
247#  define global
248#  define local
249#  define constant const
250#  define powr(a,b) pow(a,b)
251#else
252#  ifdef USE_SINCOS
253#    define SINCOS(angle,svar,cvar) svar=sincos(angle,&cvar)
254#  else
255#    define SINCOS(angle,svar,cvar) do {svar=sin(angle);cvar=cos(angle);} while (0)
256#  endif
257#endif
258
259// Standard mathematical constants:
260//   M_E, M_LOG2E, M_LOG10E, M_LN2, M_LN10, M_PI, M_PI_2=pi/2, M_PI_4=pi/4,
261//   M_1_PI=1/pi, M_2_PI=2/pi, M_2_SQRTPI=2/sqrt(pi), SQRT2, SQRT1_2=sqrt(1/2)
262// OpenCL defines M_constant_F for float constants, and nothing if double
263// is not enabled on the card, which is why these constants may be missing
264#ifndef M_PI
265#  define M_PI 3.141592653589793
266#endif
267#ifndef M_PI_2
268#  define M_PI_2 1.570796326794897
269#endif
270#ifndef M_PI_4
271#  define M_PI_4 0.7853981633974483
272#endif
273
274// Non-standard pi/180, used for converting between degrees and radians
275#ifndef M_PI_180
276#  define M_PI_180 0.017453292519943295
277#endif
278"""
279
280
281# The I(q) kernel and the I(qx, qy) kernel have one and two q parameters
282# respectively, so the template builder will need to do extra work to
283# declare, initialize and pass the q parameters.
284KERNEL_1D = {
285    'fn': "Iq",
286    'q_par_decl': "global const double *q,",
287    'qinit': "const double qi = q[i];",
288    'qcall': "qi",
289    'qwork': ["q"],
290    }
291
292KERNEL_2D = {
293    'fn': "Iqxy",
294    'q_par_decl': "global const double *qx,\n    global const double *qy,",
295    'qinit': "const double qxi = qx[i];\n    const double qyi = qy[i];",
296    'qcall': "qxi, qyi",
297    'qwork': ["qx", "qy"],
298    }
299
300# Generic kernel template for the polydispersity loop.
301# This defines the opencl kernel that is available to the host.  The same
302# structure is used for Iq and Iqxy kernels, so extra flexibility is needed
303# for q parameters.  The polydispersity loop is built elsewhere and
304# substituted into this template.
305KERNEL_TEMPLATE = """\
306kernel void %(name)s(
307    %(q_par_decl)s
308    global double *result,
309#ifdef USE_OPENCL
310    global double *loops_g,
311#else
312    const int Nq,
313#endif
314    local double *loops,
315    const double cutoff,
316    %(par_decl)s
317    )
318{
319#ifdef USE_OPENCL
320  // copy loops info to local memory
321  event_t e = async_work_group_copy(loops, loops_g, (%(pd_length)s)*2, 0);
322  wait_group_events(1, &e);
323
324  int i = get_global_id(0);
325  int Nq = get_global_size(0);
326#endif
327
328#ifdef USE_OPENCL
329  if (i < Nq)
330#else
331  #pragma omp parallel for
332  for (int i=0; i < Nq; i++)
333#endif
334  {
335    %(qinit)s
336    double ret=0.0, norm=0.0;
337    double vol=0.0, norm_vol=0.0;
338%(loops)s
339    if (vol*norm_vol != 0.0) {
340      ret *= norm_vol/vol;
341    }
342    result[i] = scale*ret/norm+background;
343  }
344}
345"""
346
347# Polydispersity loop level.
348# This pulls the parameter value and weight from the looping vector in order
349# in preperation for a nested loop.
350LOOP_OPEN="""\
351for (int %(name)s_i=0; %(name)s_i < N%(name)s; %(name)s_i++) {
352  const double %(name)s = loops[2*(%(name)s_i%(offset)s)];
353  const double %(name)s_w = loops[2*(%(name)s_i%(offset)s)+1];\
354"""
355
356
357
358##########################################################
359#                                                        #
360#   !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!   #
361#   !!                                              !!   #
362#   !!  KEEP THIS CODE CONSISTENT WITH PYKERNEL.PY  !!   #
363#   !!                                              !!   #
364#   !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!   #
365#                                                        #
366##########################################################
367
368
369# Polydispersity loop body.
370# This computes the weight, and if it is sufficient, calls the scattering
371# function and adds it to the total.  If there is a volume normalization,
372# it will also be added here.
373LOOP_BODY="""\
374const double weight = %(weight_product)s;
375if (weight > cutoff) {
376  const double I = %(fn)s(%(qcall)s, %(pcall)s);
377  if (I>=0.0) { // scattering cannot be negative
378    ret += weight*I%(sasview_spherical)s;
379    norm += weight;
380    %(volume_norm)s
381  }
382  //else { printf("exclude qx,qy,I:%%g,%%g,%%g\\n",%(qcall)s,I); }
383}
384//else { printf("exclude weight:%%g\\n",weight); }\
385"""
386
387# Use this when integrating over orientation
388SPHERICAL_CORRECTION="""\
389// Correction factor for spherical integration p(theta) I(q) sin(theta) dtheta
390double spherical_correction = (Ntheta>1 ? fabs(sin(M_PI_180*theta)) : 1.0);\
391"""
392# Use this to reproduce sasview behaviour
393SASVIEW_SPHERICAL_CORRECTION="""\
394// Correction factor for spherical integration p(theta) I(q) sin(theta) dtheta
395double spherical_correction = (Ntheta>1 ? fabs(cos(M_PI_180*theta))*M_PI_2 : 1.0);\
396"""
397
398# Volume normalization.
399# If there are "volume" polydispersity parameters, then these will be used
400# to call the form_volume function from the user supplied kernel, and accumulate
401# a normalized weight.
402VOLUME_NORM="""const double vol_weight = %(vol_weight)s;
403    vol += vol_weight*form_volume(%(vol_pars)s);
404    norm_vol += vol_weight;\
405"""
406
407# functions defined as strings in the .py module
408WORK_FUNCTION="""\
409double %(name)s(%(pars)s);
410double %(name)s(%(pars)s)
411{
412%(body)s
413}\
414"""
415
416# Documentation header for the module, giving the model name, its short
417# description and its parameter table.  The remainder of the doc comes
418# from the module docstring.
419DOC_HEADER=""".. _%(name)s:
420
421%(label)s
422=======================================================
423
424%(title)s
425
426%(parameters)s
427
428The returned value is scaled to units of |cm^-1|.
429
430%(docs)s
431"""
432
433def indent(s, depth):
434    """
435    Indent a string of text with *depth* additional spaces on each line.
436    """
437    spaces = " "*depth
438    sep = "\n"+spaces
439    return spaces + sep.join(s.split("\n"))
440
441
442def kernel_name(info, is_2D):
443    """
444    Name of the exported kernel symbol.
445    """
446    return info['name'] + "_" + ("Iqxy" if is_2D else "Iq")
447
448
449def use_single(source):
450    """
451    Convert code from double precision to single precision.
452    """
453    # Convert double keyword to float.  Accept an 'n' parameter for vector
454    # values, where n is 2, 4, 8 or 16. Assume complex numbers are represented
455    # as cdouble which is typedef'd to double2.
456    source = re.sub(r'(^|[^a-zA-Z0-9_]c?)double(([248]|16)?($|[^a-zA-Z0-9_]))',
457                    r'\1float\2', source)
458    # Convert floating point constants to single by adding 'f' to the end.
459    # OS/X driver complains if you don't do this.
460    source = re.sub(r'[^a-zA-Z_](\d*[.]\d+|\d+[.]\d*)([eE][+-]?\d+)?',
461                    r'\g<0>f', source)
462    return source
463
464
465def make_kernel(info, is_2D):
466    """
467    Build a kernel call from metadata supplied by the user.
468
469    *info* is the json object defined in the kernel file.
470
471    *form* is either "Iq" or "Iqxy".
472
473    This does not create a complete OpenCL kernel source, only the top
474    level kernel call with polydispersity and a call to the appropriate
475    Iq or Iqxy function.
476    """
477
478    # If we are building the Iqxy kernel, we need to propagate qx,qy
479    # parameters, otherwise we can
480    dim = "2d" if is_2D else "1d"
481    fixed_pars = info['partype']['fixed-'+dim]
482    pd_pars = info['partype']['pd-'+dim]
483    vol_pars = info['partype']['volume']
484    q_pars = KERNEL_2D if is_2D else KERNEL_1D
485    fn = q_pars['fn']
486
487    # Build polydispersity loops
488    depth = 4
489    offset = ""
490    loop_head = []
491    loop_end = []
492    for name in pd_pars:
493        subst = { 'name': name, 'offset': offset }
494        loop_head.append(indent(LOOP_OPEN%subst, depth))
495        loop_end.insert(0, (" "*depth) + "}")
496        offset += '+N'+name
497        depth += 2
498
499    # The volume parameters in the inner loop are used to call the volume()
500    # function in the kernel, with the parameters defined in vol_pars and the
501    # weight product defined in weight.  If there are no volume parameters,
502    # then there will be no volume normalization.
503    if vol_pars:
504        subst = {
505            'vol_weight': "*".join(p+"_w" for p in vol_pars),
506            'vol_pars': ", ".join(vol_pars),
507            }
508        volume_norm = VOLUME_NORM%subst
509    else:
510        volume_norm = ""
511
512    # Define the inner loop function call
513    # The parameters to the f(q,p1,p2...) call should occur in the same
514    # order as given in the parameter info structure.  This may be different
515    # from the parameter order in the call to the kernel since the kernel
516    # call places all fixed parameters before all polydisperse parameters.
517    fq_pars = [p[0] for p in info['parameters'][len(COMMON_PARAMETERS):]
518               if p[0] in set(fixed_pars+pd_pars)]
519    if False and "theta" in pd_pars:
520        spherical_correction = [indent(SPHERICAL_CORRECTION, depth)]
521        weights = [p+"_w" for p in pd_pars]+['spherical_correction']
522        sasview_spherical = ""
523    elif True and "theta" in pd_pars:
524        spherical_correction = [indent(SASVIEW_SPHERICAL_CORRECTION,depth)]
525        weights = [p+"_w" for p in pd_pars]
526        sasview_spherical = "*spherical_correction"
527    else:
528        spherical_correction = []
529        weights = [p+"_w" for p in pd_pars]
530        sasview_spherical = ""
531    subst = {
532        'weight_product': "*".join(weights),
533        'volume_norm': volume_norm,
534        'fn': fn,
535        'qcall': q_pars['qcall'],
536        'pcall': ", ".join(fq_pars), # skip scale and background
537        'sasview_spherical': sasview_spherical,
538        }
539    loop_body = [indent(LOOP_BODY%subst, depth)]
540    loops = "\n".join(loop_head+spherical_correction+loop_body+loop_end)
541
542    # declarations for non-pd followed by pd pars
543    # e.g.,
544    #     const double sld,
545    #     const int Nradius
546    fixed_par_decl = ",\n    ".join("const double %s"%p for p in fixed_pars)
547    pd_par_decl = ",\n    ".join("const int N%s"%p for p in pd_pars)
548    if fixed_par_decl and pd_par_decl:
549        par_decl = ",\n    ".join((fixed_par_decl, pd_par_decl))
550    elif fixed_par_decl:
551        par_decl = fixed_par_decl
552    else:
553        par_decl = pd_par_decl
554
555    # Finally, put the pieces together in the kernel.
556    subst = {
557        # kernel name is, e.g., cylinder_Iq
558        'name': kernel_name(info, is_2D),
559        # to declare, e.g., global double q[],
560        'q_par_decl': q_pars['q_par_decl'],
561        # to declare, e.g., double sld, int Nradius, int Nlength
562        'par_decl': par_decl,
563        # to copy global to local pd pars we need, e.g., Nradius+Nlength
564        'pd_length': "+".join('N'+p for p in pd_pars),
565        # the q initializers, e.g., double qi = q[i];
566        'qinit': q_pars['qinit'],
567        # the actual polydispersity loop
568        'loops': loops,
569        }
570    kernel = KERNEL_TEMPLATE%subst
571
572    # If the working function is defined in the kernel metadata as a
573    # string, translate the string to an actual function definition
574    # and put it before the kernel.
575    if info[fn]:
576        subst = {
577            'name': fn,
578            'pars': ", ".join("double "+p for p in q_pars['qwork']+fq_pars),
579            'body': info[fn],
580            }
581        kernel = "\n".join((WORK_FUNCTION%subst, kernel))
582    return kernel
583
584def make_partable(pars):
585    """
586    Generate the parameter table to include in the sphinx documentation.
587    """
588    pars = COMMON_PARAMETERS + pars
589    column_widths = [
590        max(len(p[0]) for p in pars),
591        max(len(p[-1]) for p in pars),
592        max(len(RST_UNITS[p[1]]) for p in pars),
593        PARTABLE_VALUE_WIDTH,
594        ]
595    column_widths = [max(w, len(h))
596                     for w,h in zip(column_widths, PARTABLE_HEADERS)]
597
598    sep = " ".join("="*w for w in column_widths)
599    lines = [
600        sep,
601        " ".join("%-*s"%(w,h) for w,h in zip(column_widths, PARTABLE_HEADERS)),
602        sep,
603        ]
604    for p in pars:
605        lines.append(" ".join([
606            "%-*s"%(column_widths[0],p[0]),
607            "%-*s"%(column_widths[1],p[-1]),
608            "%-*s"%(column_widths[2],RST_UNITS[p[1]]),
609            "%*g"%(column_widths[3],p[2]),
610            ]))
611    lines.append(sep)
612    return "\n".join(lines)
613
614def _search(search_path, filename):
615    """
616    Find *filename* in *search_path*.
617
618    Raises ValueError if file does not exist.
619    """
620    for path in search_path:
621        target = os.path.join(path, filename)
622        if os.path.exists(target):
623            return target
624    raise ValueError("%r not found in %s"%(filename, search_path))
625
626def sources(info):
627    """
628    Return a list of the sources file paths for the module.
629    """
630    from os.path import abspath, dirname, join as joinpath
631    search_path = [ dirname(info['filename']),
632                    abspath(joinpath(dirname(__file__),'models')) ]
633    return [_search(search_path, f) for f in info['source']]
634
635def make_model(info):
636    """
637    Generate the code for the kernel defined by info, using source files
638    found in the given search path.
639    """
640    source = [open(f).read() for f in sources(info)]
641    # If the form volume is defined as a string, then wrap it in a
642    # function definition and place it after the external sources but
643    # before the kernel functions.  If the kernel functions are strings,
644    # they will be translated in the make_kernel call.
645    if info['form_volume']:
646        subst = {
647            'name': "form_volume",
648            'pars': ", ".join("double "+p for p in info['partype']['volume']),
649            'body': info['form_volume'],
650            }
651        source.append(WORK_FUNCTION%subst)
652    kernel_Iq = make_kernel(info, is_2D=False) if not callable(info['Iq']) else ""
653    kernel_Iqxy = make_kernel(info, is_2D=True) if not callable(info['Iqxy']) else ""
654    kernel = "\n\n".join([KERNEL_HEADER]+source+[kernel_Iq, kernel_Iqxy])
655    return kernel
656
657def categorize_parameters(pars):
658    """
659    Build parameter categories out of the the parameter definitions.
660
661    Returns a dictionary of categories.
662    """
663    partype = {
664        'volume': [], 'orientation': [], 'magnetic': [], '': [],
665        'fixed-1d': [], 'fixed-2d': [], 'pd-1d': [], 'pd-2d': [],
666        'pd-rel': set(),
667    }
668
669    for p in pars:
670        name,ptype = p[0],p[4]
671        if ptype == 'volume':
672            partype['pd-1d'].append(name)
673            partype['pd-2d'].append(name)
674            partype['pd-rel'].add(name)
675        elif ptype == 'magnetic':
676            partype['fixed-2d'].append(name)
677        elif ptype == 'orientation':
678            partype['pd-2d'].append(name)
679        elif ptype == '':
680            partype['fixed-1d'].append(name)
681            partype['fixed-2d'].append(name)
682        else:
683            raise ValueError("unknown parameter type %r"%ptype)
684        partype[ptype].append(name)
685
686    return partype
687
688def make(kernel_module):
689    """
690    Build an OpenCL/ctypes function from the definition in *kernel_module*.
691
692    The module can be loaded with a normal python import statement if you
693    know which module you need, or with __import__('sasmodels.model.'+name)
694    if the name is in a string.
695    """
696    # TODO: allow Iq and Iqxy to be defined in python
697    from os.path import abspath
698    #print kernelfile
699    info = dict(
700        filename = abspath(kernel_module.__file__),
701        name = kernel_module.name,
702        title = kernel_module.title,
703        description = kernel_module.description,
704        parameters = COMMON_PARAMETERS + kernel_module.parameters,
705        source = getattr(kernel_module, 'source', []),
706        )
707    # Fill in attributes which default to None
708    info.update((k,getattr(kernel_module, k, None))
709                for k in ('ER', 'VR', 'form_volume', 'Iq', 'Iqxy'))
710    # Fill in the derived attributes
711    info['limits'] = dict((p[0],p[3]) for p in info['parameters'])
712    info['partype'] = categorize_parameters(info['parameters'])
713
714    source = make_model(info)
715
716    return source, info
717
718def doc(kernel_module):
719    """
720    Return the documentation for the model.
721    """
722    subst = dict(name=kernel_module.name.replace('_','-'),
723                 label=" ".join(kernel_module.name.split('_')).capitalize(),
724                 title=kernel_module.title,
725                 parameters=make_partable(kernel_module.parameters),
726                 docs=kernel_module.__doc__)
727    return DOC_HEADER%subst
728
729
730
731def demo_time():
732    import datetime
733    from .models import cylinder
734    toc = lambda: (datetime.datetime.now()-tic).total_seconds()
735    tic = datetime.datetime.now()
736    source, info = make(cylinder)
737    print "time:",toc()
738
739def demo():
740    from .models import cylinder
741    source, info = make(cylinder)
742    #print doc(cylinder)
743    print source
744
745if __name__ == "__main__":
746    demo()
Note: See TracBrowser for help on using the repository browser.