source: sasmodels/sasmodels/gen.py @ 8fff00e

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

merge

  • Property mode set to 100644
File size: 25.8 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# Polydispersity loop body.
357# This computes the weight, and if it is sufficient, calls the scattering
358# function and adds it to the total.  If there is a volume normalization,
359# it will also be added here.
360LOOP_BODY="""\
361const double weight = %(weight_product)s;
362if (weight > cutoff) {
363  const double I = %(fn)s(%(qcall)s, %(pcall)s);
364  if (I>=0.0) { // scattering cannot be negative
365    ret += weight*I%(sasview_spherical)s;
366    norm += weight;
367    %(volume_norm)s
368  }
369  //else { printf("exclude qx,qy,I:%%g,%%g,%%g\\n",%(qcall)s,I); }
370}
371//else { printf("exclude weight:%%g\\n",weight); }\
372"""
373
374# Use this when integrating over orientation
375SPHERICAL_CORRECTION="""\
376// Correction factor for spherical integration p(theta) I(q) sin(theta) dtheta
377double spherical_correction = (Ntheta>1 ? fabs(sin(M_PI_180*theta)) : 1.0);\
378"""
379# Use this to reproduce sasview behaviour
380SASVIEW_SPHERICAL_CORRECTION="""\
381// Correction factor for spherical integration p(theta) I(q) sin(theta) dtheta
382double spherical_correction = (Ntheta>1 ? fabs(cos(M_PI_180*theta))*M_PI_2 : 1.0);\
383"""
384
385# Volume normalization.
386# If there are "volume" polydispersity parameters, then these will be used
387# to call the form_volume function from the user supplied kernel, and accumulate
388# a normalized weight.
389VOLUME_NORM="""const double vol_weight = %(weight)s;
390    vol += vol_weight*form_volume(%(pars)s);
391    norm_vol += vol_weight;\
392"""
393
394# functions defined as strings in the .py module
395WORK_FUNCTION="""\
396double %(name)s(%(pars)s);
397double %(name)s(%(pars)s)
398{
399%(body)s
400}\
401"""
402
403# Documentation header for the module, giving the model name, its short
404# description and its parameter table.  The remainder of the doc comes
405# from the module docstring.
406DOC_HEADER=""".. _%(name)s:
407
408%(label)s
409=======================================================
410
411%(title)s
412
413%(parameters)s
414
415The returned value is scaled to units of |cm^-1|.
416
417%(docs)s
418"""
419
420def indent(s, depth):
421    """
422    Indent a string of text with *depth* additional spaces on each line.
423    """
424    spaces = " "*depth
425    sep = "\n"+spaces
426    return spaces + sep.join(s.split("\n"))
427
428
429def kernel_name(info, is_2D):
430    """
431    Name of the exported kernel symbol.
432    """
433    return info['name'] + "_" + ("Iqxy" if is_2D else "Iq")
434
435
436def use_single(source):
437    """
438    Convert code from double precision to single precision.
439    """
440    # Convert double keyword to float.  Accept an 'n' parameter for vector
441    # values, where n is 2, 4, 8 or 16. Assume complex numbers are represented
442    # as cdouble which is typedef'd to double2.
443    source = re.sub(r'(^|[^a-zA-Z0-9_]c?)double(([248]|16)?($|[^a-zA-Z0-9_]))',
444                    r'\1float\2', source)
445    # Convert floating point constants to single by adding 'f' to the end.
446    # OS/X driver complains if you don't do this.
447    source = re.sub(r'[^a-zA-Z_](\d*[.]\d+|\d+[.]\d*)([eE][+-]?\d+)?',
448                    r'\g<0>f', source)
449    return source
450
451
452def make_kernel(info, is_2D):
453    """
454    Build a kernel call from metadata supplied by the user.
455
456    *info* is the json object defined in the kernel file.
457
458    *form* is either "Iq" or "Iqxy".
459
460    This does not create a complete OpenCL kernel source, only the top
461    level kernel call with polydispersity and a call to the appropriate
462    Iq or Iqxy function.
463    """
464
465    # If we are building the Iqxy kernel, we need to propagate qx,qy
466    # parameters, otherwise we can
467    dim = "2d" if is_2D else "1d"
468    fixed_pars = info['partype']['fixed-'+dim]
469    pd_pars = info['partype']['pd-'+dim]
470    vol_pars = info['partype']['volume']
471    q_pars = KERNEL_2D if is_2D else KERNEL_1D
472    fn = q_pars['fn']
473
474    # Build polydispersity loops
475    depth = 4
476    offset = ""
477    loop_head = []
478    loop_end = []
479    for name in pd_pars:
480        subst = { 'name': name, 'offset': offset }
481        loop_head.append(indent(LOOP_OPEN%subst, depth))
482        loop_end.insert(0, (" "*depth) + "}")
483        offset += '+N'+name
484        depth += 2
485
486    # The volume parameters in the inner loop are used to call the volume()
487    # function in the kernel, with the parameters defined in vol_pars and the
488    # weight product defined in weight.  If there are no volume parameters,
489    # then there will be no volume normalization.
490    if vol_pars:
491        subst = {
492            'weight': "*".join(p+"_w" for p in vol_pars),
493            'pars': ", ".join(vol_pars),
494            }
495        volume_norm = VOLUME_NORM%subst
496    else:
497        volume_norm = ""
498
499    # Define the inner loop function call
500    # The parameters to the f(q,p1,p2...) call should occur in the same
501    # order as given in the parameter info structure.  This may be different
502    # from the parameter order in the call to the kernel since the kernel
503    # call places all fixed parameters before all polydisperse parameters.
504    fq_pars = [p[0] for p in info['parameters'][len(COMMON_PARAMETERS):]
505               if p[0] in set(fixed_pars+pd_pars)]
506    if False and "theta" in pd_pars:
507        spherical_correction = [indent(SPHERICAL_CORRECTION, depth)]
508        weights = [p+"_w" for p in pd_pars]+['spherical_correction']
509        sasview_spherical = ""
510    elif True and "theta" in pd_pars:
511        spherical_correction = [indent(SASVIEW_SPHERICAL_CORRECTION,depth)]
512        weights = [p+"_w" for p in pd_pars]
513        sasview_spherical = "*spherical_correction"
514    else:
515        spherical_correction = []
516        weights = [p+"_w" for p in pd_pars]
517        sasview_spherical = ""
518    subst = {
519        'weight_product': "*".join(weights),
520        'volume_norm': volume_norm,
521        'fn': fn,
522        'qcall': q_pars['qcall'],
523        'pcall': ", ".join(fq_pars), # skip scale and background
524        'sasview_spherical': sasview_spherical,
525        }
526    loop_body = [indent(LOOP_BODY%subst, depth)]
527    loops = "\n".join(loop_head+spherical_correction+loop_body+loop_end)
528
529    # declarations for non-pd followed by pd pars
530    # e.g.,
531    #     const double sld,
532    #     const int Nradius
533    fixed_par_decl = ",\n    ".join("const double %s"%p for p in fixed_pars)
534    pd_par_decl = ",\n    ".join("const int N%s"%p for p in pd_pars)
535    if fixed_par_decl and pd_par_decl:
536        par_decl = ",\n    ".join((fixed_par_decl, pd_par_decl))
537    elif fixed_par_decl:
538        par_decl = fixed_par_decl
539    else:
540        par_decl = pd_par_decl
541
542    # Finally, put the pieces together in the kernel.
543    subst = {
544        # kernel name is, e.g., cylinder_Iq
545        'name': kernel_name(info, is_2D),
546        # to declare, e.g., global double q[],
547        'q_par_decl': q_pars['q_par_decl'],
548        # to declare, e.g., double sld, int Nradius, int Nlength
549        'par_decl': par_decl,
550        # to copy global to local pd pars we need, e.g., Nradius+Nlength
551        'pd_length': "+".join('N'+p for p in pd_pars),
552        # the q initializers, e.g., double qi = q[i];
553        'qinit': q_pars['qinit'],
554        # the actual polydispersity loop
555        'loops': loops,
556        }
557    kernel = KERNEL_TEMPLATE%subst
558
559    # If the working function is defined in the kernel metadata as a
560    # string, translate the string to an actual function definition
561    # and put it before the kernel.
562    if info[fn]:
563        subst = {
564            'name': fn,
565            'pars': ", ".join("double "+p for p in q_pars['qwork']+fq_pars),
566            'body': info[fn],
567            }
568        kernel = "\n".join((WORK_FUNCTION%subst, kernel))
569    return kernel
570
571def make_partable(pars):
572    """
573    Generate the parameter table to include in the sphinx documentation.
574    """
575    pars = COMMON_PARAMETERS + pars
576    column_widths = [
577        max(len(p[0]) for p in pars),
578        max(len(p[-1]) for p in pars),
579        max(len(RST_UNITS[p[1]]) for p in pars),
580        PARTABLE_VALUE_WIDTH,
581        ]
582    column_widths = [max(w, len(h))
583                     for w,h in zip(column_widths, PARTABLE_HEADERS)]
584
585    sep = " ".join("="*w for w in column_widths)
586    lines = [
587        sep,
588        " ".join("%-*s"%(w,h) for w,h in zip(column_widths, PARTABLE_HEADERS)),
589        sep,
590        ]
591    for p in pars:
592        lines.append(" ".join([
593            "%-*s"%(column_widths[0],p[0]),
594            "%-*s"%(column_widths[1],p[-1]),
595            "%-*s"%(column_widths[2],RST_UNITS[p[1]]),
596            "%*g"%(column_widths[3],p[2]),
597            ]))
598    lines.append(sep)
599    return "\n".join(lines)
600
601def _search(search_path, filename):
602    """
603    Find *filename* in *search_path*.
604
605    Raises ValueError if file does not exist.
606    """
607    for path in search_path:
608        target = os.path.join(path, filename)
609        if os.path.exists(target):
610            return target
611    raise ValueError("%r not found in %s"%(filename, search_path))
612
613def sources(info):
614    """
615    Return a list of the sources file paths for the module.
616    """
617    from os.path import abspath, dirname, join as joinpath
618    search_path = [ dirname(info['filename']),
619                    abspath(joinpath(dirname(__file__),'models')) ]
620    return [_search(search_path, f) for f in info['source']]
621
622def make_model(info):
623    """
624    Generate the code for the kernel defined by info, using source files
625    found in the given search path.
626    """
627    source = [open(f).read() for f in sources(info)]
628    # If the form volume is defined as a string, then wrap it in a
629    # function definition and place it after the external sources but
630    # before the kernel functions.  If the kernel functions are strings,
631    # they will be translated in the make_kernel call.
632    if info['form_volume']:
633        subst = {
634            'name': "form_volume",
635            'pars': ", ".join("double "+p for p in info['partype']['volume']),
636            'body': info['form_volume'],
637            }
638        source.append(WORK_FUNCTION%subst)
639    kernel_Iq = make_kernel(info, is_2D=False)
640    kernel_Iqxy = make_kernel(info, is_2D=True)
641    kernel = "\n\n".join([KERNEL_HEADER]+source+[kernel_Iq, kernel_Iqxy])
642    return kernel
643
644def categorize_parameters(pars):
645    """
646    Build parameter categories out of the the parameter definitions.
647
648    Returns a dictionary of categories.
649    """
650    partype = {
651        'volume': [], 'orientation': [], 'magnetic': [], '': [],
652        'fixed-1d': [], 'fixed-2d': [], 'pd-1d': [], 'pd-2d': [],
653        'pd-rel': set(),
654    }
655
656    for p in pars:
657        name,ptype = p[0],p[4]
658        if ptype == 'volume':
659            partype['pd-1d'].append(name)
660            partype['pd-2d'].append(name)
661            partype['pd-rel'].add(name)
662        elif ptype == 'magnetic':
663            partype['fixed-2d'].append(name)
664        elif ptype == 'orientation':
665            partype['pd-2d'].append(name)
666        elif ptype == '':
667            partype['fixed-1d'].append(name)
668            partype['fixed-2d'].append(name)
669        else:
670            raise ValueError("unknown parameter type %r"%ptype)
671        partype[ptype].append(name)
672
673    return partype
674
675def make(kernel_module):
676    """
677    Build an OpenCL/ctypes function from the definition in *kernel_module*.
678
679    The module can be loaded with a normal python import statement if you
680    know which module you need, or with __import__('sasmodels.model.'+name)
681    if the name is in a string.
682    """
683    # TODO: allow Iq and Iqxy to be defined in python
684    from os.path import abspath
685    #print kernelfile
686    info = dict(
687        filename = abspath(kernel_module.__file__),
688        name = kernel_module.name,
689        title = kernel_module.title,
690        description = kernel_module.description,
691        parameters = COMMON_PARAMETERS + kernel_module.parameters,
692        source = getattr(kernel_module, 'source', []),
693        )
694    # Fill in attributes which default to None
695    info.update((k,getattr(kernel_module, k, None))
696                for k in ('ER', 'VR', 'form_volume', 'Iq', 'Iqxy'))
697    # Fill in the derived attributes
698    info['limits'] = dict((p[0],p[3]) for p in info['parameters'])
699    info['partype'] = categorize_parameters(info['parameters'])
700
701    source = make_model(info)
702
703    return source, info
704
705def doc(kernel_module):
706    """
707    Return the documentation for the model.
708    """
709    subst = dict(name=kernel_module.name.replace('_','-'),
710                 label=" ".join(kernel_module.name.split('_')).capitalize(),
711                 title=kernel_module.title,
712                 parameters=make_partable(kernel_module.parameters),
713                 docs=kernel_module.__doc__)
714    return DOC_HEADER%subst
715
716
717
718def demo_time():
719    import datetime
720    tic = datetime.datetime.now()
721    toc = lambda: (datetime.datetime.now()-tic).total_seconds()
722    path = os.path.dirname("__file__")
723    doc, c = make_model(os.path.join(path, "models", "cylinder.c"))
724    print "time:",toc()
725
726def demo():
727    from os.path import join as joinpath, dirname
728    c, info, doc = make_model(joinpath(dirname(__file__), "models", "cylinder.c"))
729    #print doc
730    #print c
731
732if __name__ == "__main__":
733    demo()
Note: See TracBrowser for help on using the repository browser.