source: sasmodels/sasmodels/generate.py @ d5ac45f

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

refactoring generate/kernel_template in process…

  • Property mode set to 100644
File size: 30.7 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
24    #define INVALID(v) (expr)  returns False if v.parameter is invalid
25    for some parameter or other (e.g., v.bell_radius < v.radius).  If
26    necessary, the expression can call a function.
27
28These functions are defined in a kernel module .py script and an associated
29set of .c files.  The model constructor will use them to create models with
30polydispersity across volume and orientation parameters, and provide
31scale and background parameters for each model.
32
33*Iq*, *Iqxy*, *Imagnetic* and *form_volume* should be stylized C-99
34functions written for OpenCL.  All functions need prototype declarations
35even if the are defined before they are used.  OpenCL does not support
36*#include* preprocessor directives, so instead the list of includes needs
37to be given as part of the metadata in the kernel module definition.
38The included files should be listed using a path relative to the kernel
39module, or if using "lib/file.c" if it is one of the standard includes
40provided with the sasmodels source.  The includes need to be listed in
41order so that functions are defined before they are used.
42
43Floating point values should be declared as *double*.  For single precision
44calculations, *double* will be replaced by *float*.  The single precision
45conversion will also tag floating point constants with "f" to make them
46single precision constants.  When using integral values in floating point
47expressions, they should be expressed as floating point values by including
48a decimal point.  This includes 0., 1. and 2.
49
50OpenCL has a *sincos* function which can improve performance when both
51the *sin* and *cos* values are needed for a particular argument.  Since
52this function does not exist in C99, all use of *sincos* should be
53replaced by the macro *SINCOS(value, sn, cn)* where *sn* and *cn* are
54previously declared *double* variables.  When compiled for systems without
55OpenCL, *SINCOS* will be replaced by *sin* and *cos* calls.   If *value* is
56an expression, it will appear twice in this case; whether or not it will be
57evaluated twice depends on the quality of the compiler.
58
59If the input parameters are invalid, the scattering calculator should
60return a negative number. Particularly with polydispersity, there are
61some sets of shape parameters which lead to nonsensical forms, such
62as a capped cylinder where the cap radius is smaller than the
63cylinder radius.  The polydispersity calculation will ignore these points,
64effectively chopping the parameter weight distributions at the boundary
65of the infeasible region.  The resulting scattering will be set to
66background.  This will work correctly even when polydispersity is off.
67
68*ER* and *VR* are python functions which operate on parameter vectors.
69The constructor code will generate the necessary vectors for computing
70them with the desired polydispersity.
71
72The available kernel parameters are defined as a list, with each parameter
73defined as a sublist with the following elements:
74
75    *name* is the name that will be used in the call to the kernel
76    function and the name that will be displayed to the user.  Names
77    should be lower case, with words separated by underscore.  If
78    acronyms are used, the whole acronym should be upper case.
79
80    *units* should be one of *degrees* for angles, *Ang* for lengths,
81    *1e-6/Ang^2* for SLDs.
82
83    *default value* will be the initial value for  the model when it
84    is selected, or when an initial value is not otherwise specified.
85
86    *limits = [lb, ub]* are the hard limits on the parameter value, used to
87    limit the polydispersity density function.  In the fit, the parameter limits
88    given to the fit are the limits  on the central value of the parameter.
89    If there is polydispersity, it will evaluate parameter values outside
90    the fit limits, but not outside the hard limits specified in the model.
91    If there are no limits, use +/-inf imported from numpy.
92
93    *type* indicates how the parameter will be used.  "volume" parameters
94    will be used in all functions.  "orientation" parameters will be used
95    in *Iqxy* and *Imagnetic*.  "magnetic* parameters will be used in
96    *Imagnetic* only.  If *type* is the empty string, the parameter will
97    be used in all of *Iq*, *Iqxy* and *Imagnetic*.  "sld" parameters
98    can automatically be promoted to magnetic parameters, each of which
99    will have a magnitude and a direction, which may be different from
100    other sld parameters.
101
102    *description* is a short description of the parameter.  This will
103    be displayed in the parameter table and used as a tool tip for the
104    parameter value in the user interface.
105
106The kernel module must set variables defining the kernel meta data:
107
108    *id* is an implicit variable formed from the filename.  It will be
109    a valid python identifier, and will be used as the reference into
110    the html documentation, with '_' replaced by '-'.
111
112    *name* is the model name as displayed to the user.  If it is missing,
113    it will be constructed from the id.
114
115    *title* is a short description of the model, suitable for a tool tip,
116    or a one line model summary in a table of models.
117
118    *description* is an extended description of the model to be displayed
119    while the model parameters are being edited.
120
121    *parameters* is the list of parameters.  Parameters in the kernel
122    functions must appear in the same order as they appear in the
123    parameters list.  Two additional parameters, *scale* and *background*
124    are added to the beginning of the parameter list.  They will show up
125    in the documentation as model parameters, but they are never sent to
126    the kernel functions.  Note that *effect_radius* and *volfraction*
127    must occur first in structure factor calculations.
128
129    *category* is the default category for the model.  The category is
130    two level structure, with the form "group:section", indicating where
131    in the manual the model will be located.  Models are alphabetical
132    within their section.
133
134    *source* is the list of C-99 source files that must be joined to
135    create the OpenCL kernel functions.  The files defining the functions
136    need to be listed before the files which use the functions.
137
138    *ER* is a python function defining the effective radius.  If it is
139    not present, the effective radius is 0.
140
141    *VR* is a python function defining the volume ratio.  If it is not
142    present, the volume ratio is 1.
143
144    *form_volume*, *Iq*, *Iqxy*, *Imagnetic* are strings containing the
145    C source code for the body of the volume, Iq, and Iqxy functions
146    respectively.  These can also be defined in the last source file.
147
148    *Iq* and *Iqxy* also be instead be python functions defining the
149    kernel.  If they are marked as *Iq.vectorized = True* then the
150    kernel is passed the entire *q* vector at once, otherwise it is
151    passed values one *q* at a time.  The performance improvement of
152    this step is significant.
153
154    *demo* is a dictionary of parameter=value defining a set of
155    parameters to use by default when *compare* is called.  Any
156    parameter not set in *demo* gets the initial value from the
157    parameter list.  *demo* is mostly needed to set the default
158    polydispersity values for tests.
159
160    *oldname* is the name of the model in sasview before sasmodels
161    was split into its own package, and *oldpars* is a dictionary
162    of *parameter: old_parameter* pairs defining the new names for
163    the parameters.  This is used by *compare* to check the values
164    of the new model against the values of the old model before
165    you are ready to add the new model to sasmodels.
166
167
168An *model_info* dictionary is constructed from the kernel meta data and
169returned to the caller.
170
171The model evaluator, function call sequence consists of q inputs and the return vector,
172followed by the loop value/weight vector, followed by the values for
173the non-polydisperse parameters, followed by the lengths of the
174polydispersity loops.  To construct the call for 1D models, the
175categories *fixed-1d* and *pd-1d* list the names of the parameters
176of the non-polydisperse and the polydisperse parameters respectively.
177Similarly, *fixed-2d* and *pd-2d* provide parameter names for 2D models.
178The *pd-rel* category is a set of those parameters which give
179polydispersitiy as a portion of the value (so a 10% length dispersity
180would use a polydispersity value of 0.1) rather than absolute
181dispersity such as an angle plus or minus 15 degrees.
182
183The *volume* category lists the volume parameters in order for calls
184to volume within the kernel (used for volume normalization) and for
185calls to ER and VR for effective radius and volume ratio respectively.
186
187The *orientation* and *magnetic* categories list the orientation and
188magnetic parameters.  These are used by the sasview interface.  The
189blank category is for parameters such as scale which don't have any
190other marking.
191
192The doc string at the start of the kernel module will be used to
193construct the model documentation web pages.  Embedded figures should
194appear in the subdirectory "img" beside the model definition, and tagged
195with the kernel module name to avoid collision with other models.  Some
196file systems are case-sensitive, so only use lower case characters for
197file names and extensions.
198
199
200The function :func:`make` loads the metadata from the module and returns
201the kernel source.  The function :func:`make_doc` extracts the doc string
202and adds the parameter table to the top.  The function :func:`model_sources`
203returns a list of files required by the model.
204
205Code follows the C99 standard with the following extensions and conditions::
206
207    M_PI_180 = pi/180
208    M_4PI_3 = 4pi/3
209    square(x) = x*x
210    cube(x) = x*x*x
211    sinc(x) = sin(x)/x, with sin(0)/0 -> 1
212    all double precision constants must include the decimal point
213    all double declarations may be converted to half, float, or long double
214    FLOAT_SIZE is the number of bytes in the converted variables
215"""
216from __future__ import print_function
217
218# TODO: identify model files which have changed since loading and reload them.
219
220import sys
221from os.path import abspath, dirname, join as joinpath, exists, basename, \
222    splitext, getmtime
223import re
224import string
225import warnings
226from collections import namedtuple
227
228import numpy as np
229
230# TODO: promote Parameter and model_info to classes
231PARAMETER_FIELDS = ['name', 'units', 'default', 'limits', 'type', 'description']
232Parameter = namedtuple('Parameter', PARAMETER_FIELDS)
233
234#TODO: determine which functions are useful outside of generate
235#__all__ = ["model_info", "make_doc", "make_source", "convert_type"]
236
237TEMPLATE_ROOT = dirname(__file__)
238
239F16 = np.dtype('float16')
240F32 = np.dtype('float32')
241F64 = np.dtype('float64')
242try:  # CRUFT: older numpy does not support float128
243    F128 = np.dtype('float128')
244except TypeError:
245    F128 = None
246
247# Scale and background, which are parameters common to every form factor
248COMMON_PARAMETERS = [
249    ["scale", "", 1, [0, np.inf], "", "Source intensity"],
250    ["background", "1/cm", 1e-3, [0, np.inf], "", "Source background"],
251    ]
252
253# Conversion from units defined in the parameter table for each model
254# to units displayed in the sphinx documentation.
255RST_UNITS = {
256    "Ang": "|Ang|",
257    "1/Ang": "|Ang^-1|",
258    "1/Ang^2": "|Ang^-2|",
259    "1e-6/Ang^2": "|1e-6Ang^-2|",
260    "degrees": "degree",
261    "1/cm": "|cm^-1|",
262    "Ang/cm": "|Ang*cm^-1|",
263    "g/cm3": "|g/cm^3|",
264    "mg/m2": "|mg/m^2|",
265    "": "None",
266    }
267
268# Headers for the parameters tables in th sphinx documentation
269PARTABLE_HEADERS = [
270    "Parameter",
271    "Description",
272    "Units",
273    "Default value",
274    ]
275
276# Minimum width for a default value (this is shorter than the column header
277# width, so will be ignored).
278PARTABLE_VALUE_WIDTH = 10
279
280# Documentation header for the module, giving the model name, its short
281# description and its parameter table.  The remainder of the doc comes
282# from the module docstring.
283DOC_HEADER = """.. _%(id)s:
284
285%(name)s
286=======================================================
287
288%(title)s
289
290%(parameters)s
291
292%(returns)s
293
294%(docs)s
295"""
296
297def format_units(units):
298    """
299    Convert units into ReStructured Text format.
300    """
301    return "string" if isinstance(units, list) else RST_UNITS.get(units, units)
302
303def make_partable(pars):
304    """
305    Generate the parameter table to include in the sphinx documentation.
306    """
307    column_widths = [
308        max(len(p.name) for p in pars),
309        max(len(p.description) for p in pars),
310        max(len(format_units(p.units)) for p in pars),
311        PARTABLE_VALUE_WIDTH,
312        ]
313    column_widths = [max(w, len(h))
314                     for w, h in zip(column_widths, PARTABLE_HEADERS)]
315
316    sep = " ".join("="*w for w in column_widths)
317    lines = [
318        sep,
319        " ".join("%-*s" % (w, h)
320                 for w, h in zip(column_widths, PARTABLE_HEADERS)),
321        sep,
322        ]
323    for p in pars:
324        lines.append(" ".join([
325            "%-*s" % (column_widths[0], p.name),
326            "%-*s" % (column_widths[1], p.description),
327            "%-*s" % (column_widths[2], format_units(p.units)),
328            "%*g" % (column_widths[3], p.default),
329            ]))
330    lines.append(sep)
331    return "\n".join(lines)
332
333def _search(search_path, filename):
334    """
335    Find *filename* in *search_path*.
336
337    Raises ValueError if file does not exist.
338    """
339    for path in search_path:
340        target = joinpath(path, filename)
341        if exists(target):
342            return target
343    raise ValueError("%r not found in %s" % (filename, search_path))
344
345
346def model_sources(model_info):
347    """
348    Return a list of the sources file paths for the module.
349    """
350    search_path = [dirname(model_info['filename']),
351                   abspath(joinpath(dirname(__file__), 'models'))]
352    return [_search(search_path, f) for f in model_info['source']]
353
354
355def convert_type(source, dtype):
356    """
357    Convert code from double precision to the desired type.
358
359    Floating point constants are tagged with 'f' for single precision or 'L'
360    for long double precision.
361    """
362    if dtype == F16:
363        fbytes = 2
364        source = _convert_type(source, "float", "f")
365    elif dtype == F32:
366        fbytes = 4
367        source = _convert_type(source, "float", "f")
368    elif dtype == F64:
369        fbytes = 8
370        # no need to convert if it is already double
371    elif dtype == F128:
372        fbytes = 16
373        source = _convert_type(source, "long double", "L")
374    else:
375        raise ValueError("Unexpected dtype in source conversion: %s"%dtype)
376    return ("#define FLOAT_SIZE %d\n"%fbytes)+source
377
378
379def _convert_type(source, type_name, constant_flag):
380    """
381    Replace 'double' with *type_name* in *source*, tagging floating point
382    constants with *constant_flag*.
383    """
384    # Convert double keyword to float/long double/half.
385    # Accept an 'n' # parameter for vector # values, where n is 2, 4, 8 or 16.
386    # Assume complex numbers are represented as cdouble which is typedef'd
387    # to double2.
388    source = re.sub(r'(^|[^a-zA-Z0-9_]c?)double(([248]|16)?($|[^a-zA-Z0-9_]))',
389                    r'\1%s\2'%type_name, source)
390    # Convert floating point constants to single by adding 'f' to the end,
391    # or long double with an 'L' suffix.  OS/X complains if you don't do this.
392    source = re.sub(r'[^a-zA-Z_](\d*[.]\d+|\d+[.]\d*)([eE][+-]?\d+)?',
393                    r'\g<0>%s'%constant_flag, source)
394    return source
395
396
397def kernel_name(model_info, is_2d):
398    """
399    Name of the exported kernel symbol.
400    """
401    return model_info['name'] + "_" + ("Iqxy" if is_2d else "Iq")
402
403
404def indent(s, depth):
405    """
406    Indent a string of text with *depth* additional spaces on each line.
407    """
408    spaces = " "*depth
409    sep = "\n" + spaces
410    return spaces + sep.join(s.split("\n"))
411
412
413_template_cache = {}
414def load_template(filename):
415    path = joinpath(TEMPLATE_ROOT, filename)
416    mtime = getmtime(path)
417    if filename not in _template_cache or mtime > _template_cache[filename][0]:
418        with open(path) as fid:
419            _template_cache[filename] = (mtime, fid.read())
420    return _template_cache[filename][1]
421
422def _gen_fn(name, pars, body):
423    """
424    Generate a function given pars and body.
425
426    Returns the following string::
427
428         double fn(double a, double b, ...);
429         double fn(double a, double b, ...) {
430             ....
431         }
432    """
433    template = """\
434double %(name)s(%(pars)s);
435double %(name)s(%(pars)s) {
436    %(body)s
437}
438
439
440"""
441    par_decl = ', '.join('double ' + p for p in pars) if pars else 'void'
442    return template % {'name': name, 'body': body, 'pars': par_decl}
443
444def _gen_call_pars(name, pars):
445    name += "."
446    return ",".join(name+p for p in pars)
447
448def make_source(model_info):
449    """
450    Generate the OpenCL/ctypes kernel from the module info.
451
452    Uses source files found in the given search path.
453    """
454    if callable(model_info['Iq']):
455        return None
456
457    # TODO: need something other than volume to indicate dispersion parameters
458    # No volume normalization despite having a volume parameter.
459    # Thickness is labelled a volume in order to trigger polydispersity.
460    # May want a separate dispersion flag, or perhaps a separate category for
461    # disperse, but not volume.  Volume parameters also use relative values
462    # for the distribution rather than the absolute values used by angular
463    # dispersion.  Need to be careful that necessary parameters are available
464    # for computing volume even if we allow non-disperse volume parameters.
465
466    # Load template
467    source = [load_template('kernel_header.c')]
468
469    # Load additional sources
470    source += [open(f).read() for f in model_sources(model_info)]
471
472    # Prepare defines
473    defines = []
474
475    iq_parameters = [p.name
476                     for p in model_info['parameters'][2:]  # skip scale, background
477                     if p.name in model_info['par_set']['1d']]
478    iqxy_parameters = [p.name
479                       for p in model_info['parameters'][2:]  # skip scale, background
480                       if p.name in model_info['par_set']['2d']]
481    volume_parameters = model_info['par_type']['volume']
482
483    # Generate form_volume function, etc. from body only
484    if model_info['form_volume'] is not None:
485        pnames = [p.name for p in volume_parameters]
486        source.append(_gen_fn('form_volume', pnames, model_info['form_volume']))
487    if model_info['Iq'] is not None:
488        pnames = ['q'] + [p.name for p in iq_parameters]
489        source.append(_gen_fn('Iq', pnames, model_info['Iq']))
490    if model_info['Iqxy'] is not None:
491        pnames = ['qx', 'qy'] + [p.name for p in iqxy_parameters]
492        source.append(_gen_fn('Iqxy', pnames, model_info['Iqxy']))
493
494    # Fill in definitions for volume parameters
495    if volume_parameters:
496        deref_vol = ",".join("v."+p.name for p in volume_parameters)
497        defines.append(('CALL_VOLUME(v)', 'form_volume(%s)\n'%deref_vol))
498    else:
499        # Model doesn't have volume.  We could make the kernel run a little
500        # faster by not using/transferring the volume normalizations, but
501        # the ifdef's reduce readability more than is worthwhile.
502        defines.append(('CALL_VOLUME(v)', '0.0'))
503
504    # Fill in definitions for Iq parameters
505    defines.append(('KERNEL_NAME', model_info['name']))
506    defines.append(('IQ_PARAMETERS', ', '.join(iq_parameters)))
507    if fixed_1d:
508        defines.append(('IQ_FIXED_PARAMETER_DECLARATIONS',
509                        ', \\\n    '.join('const double %s' % p for p in fixed_1d)))
510    # Fill in definitions for Iqxy parameters
511    defines.append(('IQXY_KERNEL_NAME', model_info['name'] + '_Iqxy'))
512    defines.append(('IQXY_PARAMETERS', ', '.join(iqxy_parameters)))
513    if fixed_2d:
514        defines.append(('IQXY_FIXED_PARAMETER_DECLARATIONS',
515                        ', \\\n    '.join('const double %s' % p for p in fixed_2d)))
516    if pd_2d:
517        defines.append(('IQXY_WEIGHT_PRODUCT',
518                        '*'.join(p + '_w' for p in pd_2d)))
519        defines.append(('IQXY_DISPERSION_LENGTH_DECLARATIONS',
520                        ', \\\n    '.join('const int N%s' % p for p in pd_2d)))
521        defines.append(('IQXY_DISPERSION_LENGTH_SUM',
522                        '+'.join('N' + p for p in pd_2d)))
523        open_loops, close_loops = build_polydispersity_loops(pd_2d)
524        defines.append(('IQXY_OPEN_LOOPS',
525                        open_loops.replace('\n', ' \\\n')))
526        defines.append(('IQXY_CLOSE_LOOPS',
527                        close_loops.replace('\n', ' \\\n')))
528    # Need to know if we have a theta parameter for Iqxy; it is not there
529    # for the magnetic sphere model, for example, which has a magnetic
530    # orientation but no shape orientation.
531    if 'theta' in pd_2d:
532        defines.append(('IQXY_HAS_THETA', '1'))
533
534    #for d in defines: print(d)
535    defines = '\n'.join('#define %s %s' % (k, v) for k, v in defines)
536    sources = '\n\n'.join(source)
537    return C_KERNEL_TEMPLATE % {
538        'DEFINES': defines,
539        'SOURCES': sources,
540        }
541
542def categorize_parameters(pars):
543    """
544    Categorize the parameters by use:
545
546    * *pd* list of polydisperse parameters in order; gui should test whether
547      they are in *2d* or *magnetic* as appropriate for the data
548    * *1d* set of parameters that are used to compute 1D patterns
549    * *2d* set of parameters that are used to compute 2D patterns (which
550      includes all 1D parameters)
551    * *magnetic* set of parameters that are used to compute magnetic
552      patterns (which includes all 1D and 2D parameters)
553    * *sesans* set of parameters that are used to compute sesans patterns
554     (which is just 1D without background)
555    * *pd-relative* is the set of parameters with relative distribution
556      width (e.g., radius +/- 10%) rather than absolute distribution
557      width (e.g., theta +/- 6 degrees).
558    """
559    par_set = {}
560    par_set['1d'] = [p for p in pars if p.type not in ('orientation', 'magnetic')]
561    par_set['2d'] = [p for p in pars if p.type != 'magnetic']
562    par_set['magnetic'] = [p for p in pars]
563    par_set['pd'] = [p for p in pars if p.type in ('volume', 'orientation')]
564    par_set['pd_relative'] = [p for p in pars if p.type == 'volume']
565    return par_set
566
567def collect_types(pars):
568    """
569    Build parameter categories out of the the parameter definitions.
570
571    Returns a dictionary of categories.
572
573    Note: these categories are subject to change, depending on the needs of
574    the UI and the needs of the kernel calling function.
575
576    The categories are as follows:
577
578    * *volume* list of volume parameter names
579    * *orientation* list of orientation parameters
580    * *magnetic* list of magnetic parameters
581    * *sld* list of parameters that have no type info
582    * *other* list of parameters that have no type info
583
584    Each parameter is in one and only one category.
585    """
586    par_type = {
587        'volume': [], 'orientation': [], 'magnetic': [], 'sld': [], 'other': [],
588    }
589    for p in pars:
590        par_type[p.type if p.type else 'other'].append(p.name)
591    return  par_type
592
593
594def process_parameters(model_info):
595    """
596    Process parameter block, precalculating parameter details.
597    """
598    # convert parameters into named tuples
599    for p in model_info['parameters']:
600        if p[4] == '' and (p[0].startswith('sld') or p[0].endswith('sld')):
601            p[4] = 'sld'
602            # TODO: make sure all models explicitly label their sld parameters
603            #raise ValueError("%s.%s needs to be explicitly set to type 'sld'" %(model_info['id'], p[0]))
604
605    pars = [Parameter(*p) for p in model_info['parameters']]
606    # Fill in the derived attributes
607    model_info['parameters'] = pars
608    partype = categorize_parameters(pars)
609    model_info['limits'] = dict((p.name, p.limits) for p in pars)
610    model_info['par_type'] = collect_types(pars)
611    model_info['par_set'] = categorize_parameters(pars)
612    model_info['defaults'] = dict((p.name, p.default) for p in pars)
613    if model_info.get('demo', None) is None:
614        model_info['demo'] = model_info['defaults']
615    model_info['has_2d'] = partype['orientation'] or partype['magnetic']
616
617def make_model_info(kernel_module):
618    """
619    Interpret the model definition file, categorizing the parameters.
620
621    The module can be loaded with a normal python import statement if you
622    know which module you need, or with __import__('sasmodels.model.'+name)
623    if the name is in a string.
624
625    The *model_info* structure contains the following fields:
626
627    * *id* is the id of the kernel
628    * *name* is the display name of the kernel
629    * *filename* is the full path to the module defining the file (if any)
630    * *title* is a short description of the kernel
631    * *description* is a long description of the kernel (this doesn't seem
632      very useful since the Help button on the model page brings you directly
633      to the documentation page)
634    * *docs* is the docstring from the module.  Use :func:`make_doc` to
635    * *category* specifies the model location in the docs
636    * *parameters* is the model parameter table
637    * *single* is True if the model allows single precision
638    * *structure_factor* is True if the model is useable in a product
639    * *variant_info* contains the information required to select between
640      model variants (e.g., the list of cases) or is None if there are no
641      model variants
642    * *defaults* is the *{parameter: value}* table built from the parameter
643      description table.
644    * *limits* is the *{parameter: [min, max]}* table built from the
645      parameter description table.
646    * *partypes* categorizes the model parameters. See
647      :func:`categorize_parameters` for details.
648    * *demo* contains the *{parameter: value}* map used in compare (and maybe
649      for the demo plot, if plots aren't set up to use the default values).
650      If *demo* is not given in the file, then the default values will be used.
651    * *tests* is a set of tests that must pass
652    * *source* is the list of library files to include in the C model build
653    * *Iq*, *Iqxy*, *form_volume*, *ER*, *VR* and *sesans* are python functions
654      implementing the kernel for the module, or None if they are not
655      defined in python
656    * *oldname* is the model name in pre-4.0 Sasview
657    * *oldpars* is the *{new: old}* parameter translation table
658      from pre-4.0 Sasview
659    * *composition* is None if the model is independent, otherwise it is a
660      tuple with composition type ('product' or 'mixture') and a list of
661      *model_info* blocks for the composition objects.  This allows us to
662      build complete product and mixture models from just the info.
663
664    """
665    # TODO: maybe turn model_info into a class ModelDefinition
666    parameters = COMMON_PARAMETERS + kernel_module.parameters
667    filename = abspath(kernel_module.__file__)
668    kernel_id = splitext(basename(filename))[0]
669    name = getattr(kernel_module, 'name', None)
670    if name is None:
671        name = " ".join(w.capitalize() for w in kernel_id.split('_'))
672    model_info = dict(
673        id=kernel_id,  # string used to load the kernel
674        filename=abspath(kernel_module.__file__),
675        name=name,
676        title=kernel_module.title,
677        description=kernel_module.description,
678        parameters=parameters,
679        composition=None,
680        docs=kernel_module.__doc__,
681        category=getattr(kernel_module, 'category', None),
682        single=getattr(kernel_module, 'single', True),
683        structure_factor=getattr(kernel_module, 'structure_factor', False),
684        variant_info=getattr(kernel_module, 'invariant_info', None),
685        demo=getattr(kernel_module, 'demo', None),
686        source=getattr(kernel_module, 'source', []),
687        oldname=getattr(kernel_module, 'oldname', None),
688        oldpars=getattr(kernel_module, 'oldpars', {}),
689        tests=getattr(kernel_module, 'tests', []),
690        )
691    process_parameters(model_info)
692    # Check for optional functions
693    functions = "ER VR form_volume Iq Iqxy shape sesans".split()
694    model_info.update((k, getattr(kernel_module, k, None)) for k in functions)
695    return model_info
696
697section_marker = re.compile(r'\A(?P<first>[%s])(?P=first)*\Z'
698                            %re.escape(string.punctuation))
699def _convert_section_titles_to_boldface(lines):
700    """
701    Do the actual work of identifying and converting section headings.
702    """
703    prior = None
704    for line in lines:
705        if prior is None:
706            prior = line
707        elif section_marker.match(line):
708            if len(line) >= len(prior):
709                yield "".join(("**", prior, "**"))
710                prior = None
711            else:
712                yield prior
713                prior = line
714        else:
715            yield prior
716            prior = line
717    if prior is not None:
718        yield prior
719
720def convert_section_titles_to_boldface(s):
721    """
722    Use explicit bold-face rather than section headings so that the table of
723    contents is not polluted with section names from the model documentation.
724
725    Sections are identified as the title line followed by a line of punctuation
726    at least as long as the title line.
727    """
728    return "\n".join(_convert_section_titles_to_boldface(s.split('\n')))
729
730def make_doc(model_info):
731    """
732    Return the documentation for the model.
733    """
734    Iq_units = "The returned value is scaled to units of |cm^-1| |sr^-1|, absolute scale."
735    Sq_units = "The returned value is a dimensionless structure factor, $S(q)$."
736    docs = convert_section_titles_to_boldface(model_info['docs'])
737    subst = dict(id=model_info['id'].replace('_', '-'),
738                 name=model_info['name'],
739                 title=model_info['title'],
740                 parameters=make_partable(model_info['parameters']),
741                 returns=Sq_units if model_info['structure_factor'] else Iq_units,
742                 docs=docs)
743    return DOC_HEADER % subst
744
745
746
747def demo_time():
748    """
749    Show how long it takes to process a model.
750    """
751    from .models import cylinder
752    import datetime
753    tic = datetime.datetime.now()
754    make_source(make_model_info(cylinder))
755    toc = (datetime.datetime.now() - tic).total_seconds()
756    print("time: %g"%toc)
757
758def main():
759    """
760    Program which prints the source produced by the model.
761    """
762    if len(sys.argv) <= 1:
763        print("usage: python -m sasmodels.generate modelname")
764    else:
765        name = sys.argv[1]
766        import sasmodels.models
767        __import__('sasmodels.models.' + name)
768        model = getattr(sasmodels.models, name)
769        model_info = make_model_info(model)
770        source = make_source(model_info)
771        print(source)
772
773if __name__ == "__main__":
774    main()
Note: See TracBrowser for help on using the repository browser.