source: sasmodels/sasmodels/generate.py @ cd3dba0

core_shell_microgelscostrafo411magnetic_modelrelease_v0.94release_v0.95ticket-1257-vesicle-productticket_1156ticket_1265_superballticket_822_more_unit_tests
Last change on this file since cd3dba0 was cd3dba0, checked in by Paul Kienzle <pkienzle@…>, 8 years ago

improve compare.py so that parameters can be constrained to valid values

  • Property mode set to 100644
File size: 23.8 KB
RevLine 
[a7684e5]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
[5d4777d]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
[994d77f]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.
[5d4777d]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
[994d77f]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.
[5d4777d]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.
[a7684e5]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
[ff7119b]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
[a7684e5]99The kernel module must set variables defining the kernel meta data:
100
[cd3dba0]101    *id* is an implicit variable formed from the filename.  It will be
102    a valid python identifier, and will be used as the reference into
103    the html documentation, with '_' replaced by '-'.
104
105    *name* is the model name as displayed to the user.  If it is missing,
106    it will be constructed from the id.
[a7684e5]107
108    *title* is a short description of the model, suitable for a tool tip,
109    or a one line model summary in a table of models.
110
111    *description* is an extended description of the model to be displayed
112    while the model parameters are being edited.
113
[ff7119b]114    *parameters* is the list of parameters.  Parameters in the kernel
115    functions must appear in the same order as they appear in the
116    parameters list.  Two additional parameters, *scale* and *background*
117    are added to the beginning of the parameter list.  They will show up
118    in the documentation as model parameters, but they are never sent to
119    the kernel functions.
[a7684e5]120
[cd3dba0]121    *category* is the default category for the model.  Models in the
122    *structure-factor* category do not have *scale* and *background*
123    added.
124
[a7684e5]125    *source* is the list of C-99 source files that must be joined to
126    create the OpenCL kernel functions.  The files defining the functions
127    need to be listed before the files which use the functions.
128
129    *ER* is a python function defining the effective radius.  If it is
130    not present, the effective radius is 0.
131
132    *VR* is a python function defining the volume ratio.  If it is not
133    present, the volume ratio is 1.
134
[5d4777d]135    *form_volume*, *Iq*, *Iqxy*, *Imagnetic* are strings containing the
136    C source code for the body of the volume, Iq, and Iqxy functions
137    respectively.  These can also be defined in the last source file.
138
[a503bfd]139    *Iq* and *Iqxy* also be instead be python functions defining the
140    kernel.  If they are marked as *Iq.vectorized = True* then the
141    kernel is passed the entire *q* vector at once, otherwise it is
142    passed values one *q* at a time.  The performance improvement of
143    this step is significant.
144
145    *demo* is a dictionary of parameter=value defining a set of
[cd3dba0]146    parameters to use by default when *compare* is called.  Any
147    parameter not set in *demo* gets the initial value from the
148    parameter list.  *demo* is mostly needed to set the default
149    polydispersity values for tests.
[ff7119b]150
[a503bfd]151    *oldname* is the name of the model in sasview before sasmodels
152    was split into its own package, and *oldpars* is a dictionary
153    of *parameter: old_parameter* pairs defining the new names for
154    the parameters.  This is used by *compare* to check the values
155    of the new model against the values of the old model before
156    you are ready to add the new model to sasmodels.
[ff7119b]157
[cd3dba0]158
159An *info* dictionary is constructed from the kernel meta data and
160returned to the caller.
161
[ff7119b]162The model evaluator, function call sequence consists of q inputs and the return vector,
163followed by the loop value/weight vector, followed by the values for
164the non-polydisperse parameters, followed by the lengths of the
165polydispersity loops.  To construct the call for 1D models, the
166categories *fixed-1d* and *pd-1d* list the names of the parameters
167of the non-polydisperse and the polydisperse parameters respectively.
168Similarly, *fixed-2d* and *pd-2d* provide parameter names for 2D models.
169The *pd-rel* category is a set of those parameters which give
170polydispersitiy as a portion of the value (so a 10% length dispersity
171would use a polydispersity value of 0.1) rather than absolute
172dispersity such as an angle plus or minus 15 degrees.
173
174The *volume* category lists the volume parameters in order for calls
175to volume within the kernel (used for volume normalization) and for
176calls to ER and VR for effective radius and volume ratio respectively.
177
178The *orientation* and *magnetic* categories list the orientation and
179magnetic parameters.  These are used by the sasview interface.  The
180blank category is for parameters such as scale which don't have any
181other marking.
182
[a7684e5]183The doc string at the start of the kernel module will be used to
184construct the model documentation web pages.  Embedded figures should
185appear in the subdirectory "img" beside the model definition, and tagged
186with the kernel module name to avoid collision with other models.  Some
187file systems are case-sensitive, so only use lower case characters for
188file names and extensions.
189
190
191The function :func:`make` loads the metadata from the module and returns
[ff7119b]192the kernel source.  The function :func:`doc` extracts the doc string
193and adds the parameter table to the top.  The function :func:`sources`
194returns a list of files required by the model.
[a7684e5]195"""
196
197# TODO: identify model files which have changed since loading and reload them.
198
[99e6860]199__all__ = ["make", "doc", "sources", "use_single"]
[14de349]200
[5d4777d]201import sys
[cd3dba0]202from os.path import abspath, dirname, join as joinpath, exists, basename
[994d77f]203import re
[14de349]204
205import numpy as np
[f734e7d]206C_KERNEL_TEMPLATE_PATH = joinpath(dirname(__file__), 'kernel_template.c')
[14de349]207
208F64 = np.dtype('float64')
209F32 = np.dtype('float32')
210
[ce27e21]211# Scale and background, which are parameters common to every form factor
212COMMON_PARAMETERS = [
[33e91b1]213    ["scale", "", 1, [0, np.inf], "", "Source intensity"],
214    ["background", "1/cm", 0, [0, np.inf], "", "Source background"],
[ce27e21]215    ]
216
217
[14de349]218# Conversion from units defined in the parameter table for each model
219# to units displayed in the sphinx documentation.
220RST_UNITS = {
221    "Ang": "|Ang|",
[a5d0d00]222    "1/Ang": "|Ang^-1|",
[14de349]223    "1/Ang^2": "|Ang^-2|",
224    "1e-6/Ang^2": "|1e-6Ang^-2|",
225    "degrees": "degree",
226    "1/cm": "|cm^-1|",
227    "": "None",
228    }
229
230# Headers for the parameters tables in th sphinx documentation
231PARTABLE_HEADERS = [
[19dcb933]232    "Parameter",
233    "Description",
[14de349]234    "Units",
235    "Default value",
236    ]
237
[ff7119b]238# Minimum width for a default value (this is shorter than the column header
239# width, so will be ignored).
[14de349]240PARTABLE_VALUE_WIDTH = 10
241
[a7684e5]242# Documentation header for the module, giving the model name, its short
243# description and its parameter table.  The remainder of the doc comes
244# from the module docstring.
[cd3dba0]245DOC_HEADER = """.. _%(id)s:
[a7684e5]246
[cd3dba0]247%(name)s
[a7684e5]248=======================================================
249
250%(title)s
251
252%(parameters)s
253
254The returned value is scaled to units of |cm^-1|.
255
256%(docs)s
257"""
[ce27e21]258
[4eac427]259def format_units(par):
260    return RST_UNITS.get(par, par)
261
[19dcb933]262def make_partable(pars):
[ff7119b]263    """
264    Generate the parameter table to include in the sphinx documentation.
265    """
[14de349]266    column_widths = [
267        max(len(p[0]) for p in pars),
[19dcb933]268        max(len(p[-1]) for p in pars),
[4eac427]269        max(len(format_units(p[1])) for p in pars),
[14de349]270        PARTABLE_VALUE_WIDTH,
271        ]
272    column_widths = [max(w, len(h))
[33e91b1]273                     for w, h in zip(column_widths, PARTABLE_HEADERS)]
[14de349]274
275    sep = " ".join("="*w for w in column_widths)
276    lines = [
277        sep,
[33e91b1]278        " ".join("%-*s" % (w, h) for w, h in zip(column_widths, PARTABLE_HEADERS)),
[14de349]279        sep,
280        ]
281    for p in pars:
282        lines.append(" ".join([
[33e91b1]283            "%-*s" % (column_widths[0], p[0]),
284            "%-*s" % (column_widths[1], p[-1]),
[4eac427]285            "%-*s" % (column_widths[2], format_units(p[1])),
[33e91b1]286            "%*g" % (column_widths[3], p[2]),
[14de349]287            ]))
288    lines.append(sep)
289    return "\n".join(lines)
290
[32c160a]291def _search(search_path, filename):
292    """
293    Find *filename* in *search_path*.
294
295    Raises ValueError if file does not exist.
296    """
297    for path in search_path:
[f734e7d]298        target = joinpath(path, filename)
299        if exists(target):
[32c160a]300            return target
[33e91b1]301    raise ValueError("%r not found in %s" % (filename, search_path))
[14de349]302
[ff7119b]303def sources(info):
304    """
305    Return a list of the sources file paths for the module.
306    """
[33e91b1]307    search_path = [dirname(info['filename']),
308                   abspath(joinpath(dirname(__file__), 'models'))]
[5d4777d]309    return [_search(search_path, f) for f in info['source']]
[ff7119b]310
[f734e7d]311def use_single(source):
[ff7119b]312    """
[f734e7d]313    Convert code from double precision to single precision.
[ff7119b]314    """
[f734e7d]315    # Convert double keyword to float.  Accept an 'n' parameter for vector
316    # values, where n is 2, 4, 8 or 16. Assume complex numbers are represented
317    # as cdouble which is typedef'd to double2.
318    source = re.sub(r'(^|[^a-zA-Z0-9_]c?)double(([248]|16)?($|[^a-zA-Z0-9_]))',
319                    r'\1float\2', source)
320    # Convert floating point constants to single by adding 'f' to the end.
321    # OS/X driver complains if you don't do this.
322    source = re.sub(r'[^a-zA-Z_](\d*[.]\d+|\d+[.]\d*)([eE][+-]?\d+)?',
323                    r'\g<0>f', source)
324    return source
325
326
327def kernel_name(info, is_2D):
328    """
329    Name of the exported kernel symbol.
330    """
331    return info['name'] + "_" + ("Iqxy" if is_2D else "Iq")
332
[14de349]333
[ce27e21]334def categorize_parameters(pars):
[14de349]335    """
[ce27e21]336    Build parameter categories out of the the parameter definitions.
337
338    Returns a dictionary of categories.
[14de349]339    """
[ce27e21]340    partype = {
341        'volume': [], 'orientation': [], 'magnetic': [], '': [],
342        'fixed-1d': [], 'fixed-2d': [], 'pd-1d': [], 'pd-2d': [],
343        'pd-rel': set(),
[14de349]344    }
345
[ce27e21]346    for p in pars:
[33e91b1]347        name, ptype = p[0], p[4]
[ce27e21]348        if ptype == 'volume':
349            partype['pd-1d'].append(name)
350            partype['pd-2d'].append(name)
351            partype['pd-rel'].add(name)
352        elif ptype == 'magnetic':
353            partype['fixed-2d'].append(name)
354        elif ptype == 'orientation':
355            partype['pd-2d'].append(name)
356        elif ptype == '':
357            partype['fixed-1d'].append(name)
358            partype['fixed-2d'].append(name)
359        else:
[33e91b1]360            raise ValueError("unknown parameter type %r" % ptype)
[ce27e21]361        partype[ptype].append(name)
[14de349]362
[ce27e21]363    return partype
[14de349]364
[f734e7d]365def indent(s, depth):
366    """
367    Indent a string of text with *depth* additional spaces on each line.
368    """
369    spaces = " "*depth
[33e91b1]370    sep = "\n" + spaces
[f734e7d]371    return spaces + sep.join(s.split("\n"))
372
373
374def build_polydispersity_loops(pd_pars):
375    """
376    Build polydispersity loops
377
378    Returns loop opening and loop closing
379    """
[33e91b1]380    LOOP_OPEN = """\
[f734e7d]381for (int %(name)s_i=0; %(name)s_i < N%(name)s; %(name)s_i++) {
382  const double %(name)s = loops[2*(%(name)s_i%(offset)s)];
383  const double %(name)s_w = loops[2*(%(name)s_i%(offset)s)+1];\
384"""
385    depth = 4
386    offset = ""
387    loop_head = []
388    loop_end = []
389    for name in pd_pars:
[33e91b1]390        subst = {'name': name, 'offset': offset}
391        loop_head.append(indent(LOOP_OPEN % subst, depth))
[f734e7d]392        loop_end.insert(0, (" "*depth) + "}")
[33e91b1]393        offset += '+N' + name
[f734e7d]394        depth += 2
395    return "\n".join(loop_head), "\n".join(loop_end)
396
[33e91b1]397C_KERNEL_TEMPLATE = None
[f734e7d]398def make_model(info):
399    """
400    Generate the code for the kernel defined by info, using source files
401    found in the given search path.
402    """
403    # TODO: need something other than volume to indicate dispersion parameters
404    # No volume normalization despite having a volume parameter.
405    # Thickness is labelled a volume in order to trigger polydispersity.
406    # May want a separate dispersion flag, or perhaps a separate category for
407    # disperse, but not volume.  Volume parameters also use relative values
408    # for the distribution rather than the absolute values used by angular
409    # dispersion.  Need to be careful that necessary parameters are available
410    # for computing volume even if we allow non-disperse volume parameters.
411
412    # Load template
413    global C_KERNEL_TEMPLATE
414    if C_KERNEL_TEMPLATE is None:
415        with open(C_KERNEL_TEMPLATE_PATH) as fid:
416            C_KERNEL_TEMPLATE = fid.read()
417
418    # Load additional sources
419    source = [open(f).read() for f in sources(info)]
420
421    # Prepare defines
422    defines = []
423    partype = info['partype']
424    pd_1d = partype['pd-1d']
425    pd_2d = partype['pd-2d']
426    fixed_1d = partype['fixed-1d']
427    fixed_2d = partype['fixed-1d']
428
429    iq_parameters = [p[0]
[33e91b1]430                     for p in info['parameters'][2:] # skip scale, background
431                     if p[0] in set(fixed_1d + pd_1d)]
[f734e7d]432    iqxy_parameters = [p[0]
[33e91b1]433                       for p in info['parameters'][2:] # skip scale, background
434                       if p[0] in set(fixed_2d + pd_2d)]
[f734e7d]435    volume_parameters = [p[0]
[33e91b1]436                         for p in info['parameters']
437                         if p[4] == 'volume']
[f734e7d]438
439    # Fill in defintions for volume parameters
440    if volume_parameters:
441        defines.append(('VOLUME_PARAMETERS',
442                        ','.join(volume_parameters)))
443        defines.append(('VOLUME_WEIGHT_PRODUCT',
[33e91b1]444                        '*'.join(p + '_w' for p in volume_parameters)))
[f734e7d]445
446    # Generate form_volume function from body only
447    if info['form_volume'] is not None:
[6137124]448        if volume_parameters:
[33e91b1]449            vol_par_decl = ', '.join('double ' + p for p in volume_parameters)
[6137124]450        else:
451            vol_par_decl = 'void'
[f734e7d]452        defines.append(('VOLUME_PARAMETER_DECLARATIONS',
[6137124]453                        vol_par_decl))
[f734e7d]454        fn = """\
455double form_volume(VOLUME_PARAMETER_DECLARATIONS);
456double form_volume(VOLUME_PARAMETER_DECLARATIONS) {
457    %(body)s
458}
[33e91b1]459""" % {'body':info['form_volume']}
[f734e7d]460        source.append(fn)
461
462    # Fill in definitions for Iq parameters
[33e91b1]463    defines.append(('IQ_KERNEL_NAME', info['name'] + '_Iq'))
[f734e7d]464    defines.append(('IQ_PARAMETERS', ', '.join(iq_parameters)))
465    if fixed_1d:
466        defines.append(('IQ_FIXED_PARAMETER_DECLARATIONS',
[33e91b1]467                        ', \\\n    '.join('const double %s' % p for p in fixed_1d)))
[f734e7d]468    if pd_1d:
469        defines.append(('IQ_WEIGHT_PRODUCT',
[33e91b1]470                        '*'.join(p + '_w' for p in pd_1d)))
[f734e7d]471        defines.append(('IQ_DISPERSION_LENGTH_DECLARATIONS',
[33e91b1]472                        ', \\\n    '.join('const int N%s' % p for p in pd_1d)))
[f734e7d]473        defines.append(('IQ_DISPERSION_LENGTH_SUM',
[33e91b1]474                        '+'.join('N' + p for p in pd_1d)))
[f734e7d]475        open_loops, close_loops = build_polydispersity_loops(pd_1d)
476        defines.append(('IQ_OPEN_LOOPS',
[33e91b1]477                        open_loops.replace('\n', ' \\\n')))
[f734e7d]478        defines.append(('IQ_CLOSE_LOOPS',
[33e91b1]479                        close_loops.replace('\n', ' \\\n')))
[f734e7d]480    if info['Iq'] is not None:
481        defines.append(('IQ_PARAMETER_DECLARATIONS',
[33e91b1]482                        ', '.join('double ' + p for p in iq_parameters)))
[f734e7d]483        fn = """\
484double Iq(double q, IQ_PARAMETER_DECLARATIONS);
485double Iq(double q, IQ_PARAMETER_DECLARATIONS) {
486    %(body)s
487}
[33e91b1]488""" % {'body':info['Iq']}
[f734e7d]489        source.append(fn)
490
491    # Fill in definitions for Iqxy parameters
[33e91b1]492    defines.append(('IQXY_KERNEL_NAME', info['name'] + '_Iqxy'))
[f734e7d]493    defines.append(('IQXY_PARAMETERS', ', '.join(iqxy_parameters)))
494    if fixed_2d:
495        defines.append(('IQXY_FIXED_PARAMETER_DECLARATIONS',
[33e91b1]496                        ', \\\n    '.join('const double %s' % p for p in fixed_2d)))
[f734e7d]497    if pd_2d:
498        defines.append(('IQXY_WEIGHT_PRODUCT',
[33e91b1]499                        '*'.join(p + '_w' for p in pd_2d)))
[f734e7d]500        defines.append(('IQXY_DISPERSION_LENGTH_DECLARATIONS',
[33e91b1]501                        ', \\\n    '.join('const int N%s' % p for p in pd_2d)))
[f734e7d]502        defines.append(('IQXY_DISPERSION_LENGTH_SUM',
[33e91b1]503                        '+'.join('N' + p for p in pd_2d)))
[f734e7d]504        open_loops, close_loops = build_polydispersity_loops(pd_2d)
505        defines.append(('IQXY_OPEN_LOOPS',
[33e91b1]506                        open_loops.replace('\n', ' \\\n')))
[f734e7d]507        defines.append(('IQXY_CLOSE_LOOPS',
[33e91b1]508                        close_loops.replace('\n', ' \\\n')))
[f734e7d]509    if info['Iqxy'] is not None:
510        defines.append(('IQXY_PARAMETER_DECLARATIONS',
[33e91b1]511                        ', '.join('double ' + p for p in iqxy_parameters)))
[f734e7d]512        fn = """\
513double Iqxy(double qx, double qy, IQXY_PARAMETER_DECLARATIONS);
514double Iqxy(double qx, double qy, IQXY_PARAMETER_DECLARATIONS) {
515    %(body)s
516}
[33e91b1]517""" % {'body':info['Iqxy']}
[f734e7d]518        source.append(fn)
519
520    # Need to know if we have a theta parameter for Iqxy; it is not there
521    # for the magnetic sphere model, for example, which has a magnetic
522    # orientation but no shape orientation.
523    if 'theta' in pd_2d:
524        defines.append(('IQXY_HAS_THETA', '1'))
525
526    #for d in defines: print d
[33e91b1]527    DEFINES = '\n'.join('#define %s %s' % (k, v) for k, v in defines)
528    SOURCES = '\n\n'.join(source)
529    return C_KERNEL_TEMPLATE % {
[f734e7d]530        'DEFINES':DEFINES,
531        'SOURCES':SOURCES,
532        }
533
[cd3dba0]534def make_info(kernel_module):
[14de349]535    """
[cd3dba0]536    Interpret the model definition file, categorizing the parameters.
[14de349]537    """
[ce27e21]538    #print kernelfile
[cd3dba0]539    category = getattr(kernel_module, 'category', None)
540    parameters = COMMON_PARAMETERS + kernel_module.parameters
541    # Default the demo parameters to the starting values for the individual
542    # parameters if an explicit demo parameter set has not been specified.
543    demo_parameters = getattr(kernel_module, 'demo', None)
544    if demo_parameters is None:
545        demo_parameters = dict((p[0],p[2]) for p in parameters)
546    filename = abspath(kernel_module.__file__)
547    kernel_id = basename(filename)[:-3]
548    name = getattr(kernel_module, 'name', None)
549    if name is None:
550        name = " ".join(w.capitalize() for w in kernel_id.split('_'))
[32c160a]551    info = dict(
[cd3dba0]552        id = kernel_id,  # string used to load the kernel
[33e91b1]553        filename=abspath(kernel_module.__file__),
[cd3dba0]554        name=name,
[33e91b1]555        title=kernel_module.title,
556        description=kernel_module.description,
[cd3dba0]557        category=category,
558        parameters=parameters,
559        demo=demo_parameters,
[33e91b1]560        source=getattr(kernel_module, 'source', []),
561        oldname=kernel_module.oldname,
562        oldpars=kernel_module.oldpars,
[32c160a]563        )
[5d4777d]564    # Fill in attributes which default to None
[33e91b1]565    info.update((k, getattr(kernel_module, k, None))
[5d4777d]566                for k in ('ER', 'VR', 'form_volume', 'Iq', 'Iqxy'))
567    # Fill in the derived attributes
[33e91b1]568    info['limits'] = dict((p[0], p[3]) for p in info['parameters'])
[32c160a]569    info['partype'] = categorize_parameters(info['parameters'])
[33e91b1]570    info['defaults'] = dict((p[0], p[2]) for p in info['parameters'])
[cd3dba0]571    return info
572
573def make(kernel_module):
574    """
575    Build an OpenCL/ctypes function from the definition in *kernel_module*.
[32c160a]576
[cd3dba0]577    The module can be loaded with a normal python import statement if you
578    know which module you need, or with __import__('sasmodels.model.'+name)
579    if the name is in a string.
580    """
581    info = make_info(kernel_module)
[f734e7d]582    # Assume if one part of the kernel is python then all parts are.
583    source = make_model(info) if not callable(info['Iq']) else None
[32c160a]584    return source, info
[14de349]585
[a7684e5]586def doc(kernel_module):
587    """
588    Return the documentation for the model.
589    """
[cd3dba0]590    info = make_info(kernel_module)
591    subst = dict(id=info['id'].replace('_', '-'),
592                 name=info['name'],
593                 title=info['title'],
594                 parameters=make_partable(info['parameters']),
[19dcb933]595                 docs=kernel_module.__doc__)
[33e91b1]596    return DOC_HEADER % subst
[a7684e5]597
[ff7119b]598
[14de349]599
600def demo_time():
[b3f6bc3]601    from .models import cylinder
[3c56da87]602    import datetime
[b3f6bc3]603    tic = datetime.datetime.now()
[33e91b1]604    make(cylinder)
[3c56da87]605    toc = (datetime.datetime.now() - tic).total_seconds()
606    print "time:", toc
[14de349]607
[ff58110]608def main():
609    if len(sys.argv) <= 1:
610        print "usage: python -m sasmodels.generate modelname"
611    else:
612        name = sys.argv[1]
613        import sasmodels.models
[33e91b1]614        __import__('sasmodels.models.' + name)
[ff58110]615        model = getattr(sasmodels.models, name)
[33e91b1]616        source, _ = make(model)
617        print source
[14de349]618
619if __name__ == "__main__":
[ff58110]620    main()
Note: See TracBrowser for help on using the repository browser.