source: sasmodels/sasmodels/generate.py @ abc03d8

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

Bug fixes

  • Property mode set to 100644
File size: 36.1 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
239MAX_PD = 4
240
241F16 = np.dtype('float16')
242F32 = np.dtype('float32')
243F64 = np.dtype('float64')
244try:  # CRUFT: older numpy does not support float128
245    F128 = np.dtype('float128')
246except TypeError:
247    F128 = None
248
249# Scale and background, which are parameters common to every form factor
250COMMON_PARAMETERS = [
251    ["scale", "", 1, [0, np.inf], "", "Source intensity"],
252    ["background", "1/cm", 1e-3, [0, np.inf], "", "Source background"],
253    ]
254
255# Conversion from units defined in the parameter table for each model
256# to units displayed in the sphinx documentation.
257RST_UNITS = {
258    "Ang": "|Ang|",
259    "1/Ang": "|Ang^-1|",
260    "1/Ang^2": "|Ang^-2|",
261    "1e-6/Ang^2": "|1e-6Ang^-2|",
262    "degrees": "degree",
263    "1/cm": "|cm^-1|",
264    "Ang/cm": "|Ang*cm^-1|",
265    "g/cm3": "|g/cm^3|",
266    "mg/m2": "|mg/m^2|",
267    "": "None",
268    }
269
270# Headers for the parameters tables in th sphinx documentation
271PARTABLE_HEADERS = [
272    "Parameter",
273    "Description",
274    "Units",
275    "Default value",
276    ]
277
278# Minimum width for a default value (this is shorter than the column header
279# width, so will be ignored).
280PARTABLE_VALUE_WIDTH = 10
281
282# Documentation header for the module, giving the model name, its short
283# description and its parameter table.  The remainder of the doc comes
284# from the module docstring.
285DOC_HEADER = """.. _%(id)s:
286
287%(name)s
288=======================================================
289
290%(title)s
291
292%(parameters)s
293
294%(returns)s
295
296%(docs)s
297"""
298
299def format_units(units):
300    """
301    Convert units into ReStructured Text format.
302    """
303    return "string" if isinstance(units, list) else RST_UNITS.get(units, units)
304
305def make_partable(pars):
306    """
307    Generate the parameter table to include in the sphinx documentation.
308    """
309    column_widths = [
310        max(len(p.name) for p in pars),
311        max(len(p.description) for p in pars),
312        max(len(format_units(p.units)) for p in pars),
313        PARTABLE_VALUE_WIDTH,
314        ]
315    column_widths = [max(w, len(h))
316                     for w, h in zip(column_widths, PARTABLE_HEADERS)]
317
318    sep = " ".join("="*w for w in column_widths)
319    lines = [
320        sep,
321        " ".join("%-*s" % (w, h)
322                 for w, h in zip(column_widths, PARTABLE_HEADERS)),
323        sep,
324        ]
325    for p in pars:
326        lines.append(" ".join([
327            "%-*s" % (column_widths[0], p.name),
328            "%-*s" % (column_widths[1], p.description),
329            "%-*s" % (column_widths[2], format_units(p.units)),
330            "%*g" % (column_widths[3], p.default),
331            ]))
332    lines.append(sep)
333    return "\n".join(lines)
334
335def _search(search_path, filename):
336    """
337    Find *filename* in *search_path*.
338
339    Raises ValueError if file does not exist.
340    """
341    for path in search_path:
342        target = joinpath(path, filename)
343        if exists(target):
344            return target
345    raise ValueError("%r not found in %s" % (filename, search_path))
346
347
348def model_sources(model_info):
349    """
350    Return a list of the sources file paths for the module.
351    """
352    search_path = [dirname(model_info['filename']),
353                   abspath(joinpath(dirname(__file__), 'models'))]
354    return [_search(search_path, f) for f in model_info['source']]
355
356def convert_type(source, dtype):
357    """
358    Convert code from double precision to the desired type.
359
360    Floating point constants are tagged with 'f' for single precision or 'L'
361    for long double precision.
362    """
363    if dtype == F16:
364        fbytes = 2
365        source = _convert_type(source, "float", "f")
366    elif dtype == F32:
367        fbytes = 4
368        source = _convert_type(source, "float", "f")
369    elif dtype == F64:
370        fbytes = 8
371        # no need to convert if it is already double
372    elif dtype == F128:
373        fbytes = 16
374        source = _convert_type(source, "long double", "L")
375    else:
376        raise ValueError("Unexpected dtype in source conversion: %s"%dtype)
377    return ("#define FLOAT_SIZE %d\n"%fbytes)+source
378
379
380def _convert_type(source, type_name, constant_flag):
381    """
382    Replace 'double' with *type_name* in *source*, tagging floating point
383    constants with *constant_flag*.
384    """
385    # Convert double keyword to float/long double/half.
386    # Accept an 'n' # parameter for vector # values, where n is 2, 4, 8 or 16.
387    # Assume complex numbers are represented as cdouble which is typedef'd
388    # to double2.
389    source = re.sub(r'(^|[^a-zA-Z0-9_]c?)double(([248]|16)?($|[^a-zA-Z0-9_]))',
390                    r'\1%s\2'%type_name, source)
391    # Convert floating point constants to single by adding 'f' to the end,
392    # or long double with an 'L' suffix.  OS/X complains if you don't do this.
393    source = re.sub(r'[^a-zA-Z_](\d*[.]\d+|\d+[.]\d*)([eE][+-]?\d+)?',
394                    r'\g<0>%s'%constant_flag, source)
395    return source
396
397
398def kernel_name(model_info, is_2d):
399    """
400    Name of the exported kernel symbol.
401    """
402    return model_info['name'] + "_" + ("Iqxy" if is_2d else "Iq")
403
404
405def indent(s, depth):
406    """
407    Indent a string of text with *depth* additional spaces on each line.
408    """
409    spaces = " "*depth
410    sep = "\n" + spaces
411    return spaces + sep.join(s.split("\n"))
412
413
414_template_cache = {}
415def load_template(filename):
416    path = joinpath(TEMPLATE_ROOT, filename)
417    mtime = getmtime(path)
418    if filename not in _template_cache or mtime > _template_cache[filename][0]:
419        with open(path) as fid:
420            _template_cache[filename] = (mtime, fid.read(), path)
421    return _template_cache[filename][1]
422
423def model_templates():
424    # TODO: fails DRY; templates are listed in two places.
425    # should instead have model_info contain a list of paths
426    return [joinpath(TEMPLATE_ROOT, filename)
427            for filename in ('kernel_header.c', 'kernel_iq.c')]
428
429
430_FN_TEMPLATE = """\
431double %(name)s(%(pars)s);
432double %(name)s(%(pars)s) {
433    %(body)s
434}
435
436
437"""
438
439def _gen_fn(name, pars, body):
440    """
441    Generate a function given pars and body.
442
443    Returns the following string::
444
445         double fn(double a, double b, ...);
446         double fn(double a, double b, ...) {
447             ....
448         }
449    """
450    par_decl = ', '.join('double ' + p for p in pars) if pars else 'void'
451    return _FN_TEMPLATE % {'name': name, 'body': body, 'pars': par_decl}
452
453def _call_pars(prefix, pars):
454    """
455    Return a list of *prefix.parameter* from parameter items.
456    """
457    prefix += "."
458    return [prefix+p for p in pars]
459
460_IQXY_PATTERN = re.compile("^((inline|static) )? *(double )? *Iqxy *([(]|$)",
461                           flags=re.MULTILINE)
462def _have_Iqxy(sources):
463    """
464    Return true if any file defines Iqxy.
465
466    Note this is not a C parser, and so can be easily confused by
467    non-standard syntax.  Also, it will incorrectly identify the following
468    as having Iqxy::
469
470        /*
471        double Iqxy(qx, qy, ...) { ... fill this in later ... }
472        */
473
474    If you want to comment out an Iqxy function, use // on the front of the
475    line instead.
476    """
477    for code in sources:
478        if _IQXY_PATTERN.search(code):
479            return True
480    else:
481        return False
482
483def make_source(model_info):
484    """
485    Generate the OpenCL/ctypes kernel from the module info.
486
487    Uses source files found in the given search path.
488    """
489    if callable(model_info['Iq']):
490        return None
491
492    # TODO: need something other than volume to indicate dispersion parameters
493    # No volume normalization despite having a volume parameter.
494    # Thickness is labelled a volume in order to trigger polydispersity.
495    # May want a separate dispersion flag, or perhaps a separate category for
496    # disperse, but not volume.  Volume parameters also use relative values
497    # for the distribution rather than the absolute values used by angular
498    # dispersion.  Need to be careful that necessary parameters are available
499    # for computing volume even if we allow non-disperse volume parameters.
500
501    # kernel_iq assumes scale and background are the first parameters;
502    # they should be first for 1d and 2d parameter lists as well.
503    assert model_info['parameters'][0].name == 'scale'
504    assert model_info['parameters'][1].name == 'background'
505
506    # Identify parameter types
507    iq_parameters = model_info['par_type']['1d'][2:]
508    iqxy_parameters = model_info['par_type']['2d'][2:]
509    vol_parameters = model_info['par_type']['volume']
510
511    # Load templates and user code
512    kernel_header = load_template('kernel_header.c')
513    kernel_code = load_template('kernel_iq.c')
514    user_code = [open(f).read() for f in model_sources(model_info)]
515
516    # Build initial sources
517    source = [kernel_header] + user_code
518
519    # Generate form_volume function, etc. from body only
520    if model_info['form_volume'] is not None:
521        pnames = [p for p in vol_parameters]
522        source.append(_gen_fn('form_volume', pnames, model_info['form_volume']))
523    if model_info['Iq'] is not None:
524        pnames = ['q'] + [p for p in iq_parameters]
525        source.append(_gen_fn('Iq', pnames, model_info['Iq']))
526    if model_info['Iqxy'] is not None:
527        pnames = ['qx', 'qy'] + [p for p in iqxy_parameters]
528        source.append(_gen_fn('Iqxy', pnames, model_info['Iqxy']))
529
530    # Define the parameter table
531    source.append("#define PARAMETER_TABLE \\")
532    source.append("\\\n    ".join("double %s;"%p.name
533                                   for p in model_info['parameters'][2:]))
534
535    # Define the function calls
536    if vol_parameters:
537        refs = ",".join(_call_pars("v", vol_parameters))
538        call_volume = "#define CALL_VOLUME(v) form_volume(%s)"%refs
539    else:
540        # Model doesn't have volume.  We could make the kernel run a little
541        # faster by not using/transferring the volume normalizations, but
542        # the ifdef's reduce readability more than is worthwhile.
543        call_volume = "#define CALL_VOLUME(v) 0.0"
544    source.append(call_volume)
545
546    refs = ["q[i]"] + _call_pars("v", iq_parameters)
547    call_iq = "#define CALL_IQ(q,i,v) Iq(%s)" % (",".join(refs))
548    if _have_Iqxy(user_code):
549        # Call 2D model
550        refs = ["q[2*i]", "q[2*i+1]"] + _call_pars("v", iqxy_parameters)
551        call_iqxy = "#define CALL_IQ(q,i,v) Iqxy(%s)" % (",".join(refs))
552    else:
553        # Call 1D model with sqrt(qx^2 + qy^2)
554        warnings.warn("Creating Iqxy = Iq(sqrt(qx^2 + qy^2))")
555        # still defined:: refs = ["q[i]"] + _call_pars("v", iq_parameters)
556        pars_sqrt = ["sqrt(q[2*i]*q[2*i]+q[2*i+1]*q[2*i+1])"] + refs[1:]
557        call_iqxy = "#define CALL_IQ(q,i,v) Iq(%s)" % (",".join(pars_sqrt))
558
559    # Fill in definitions for numbers of parameters
560    source.append("#define MAX_PD %s"%model_info['max_pd'])
561    source.append("#define NPARS %d"%(len(model_info['parameters'])-2))
562
563    # TODO: allow mixed python/opencl kernels?
564
565    # define the Iq kernel
566    source.append("#define KERNEL_NAME %s_Iq"%model_info['name'])
567    source.append(call_iq)
568    source.append(kernel_code)
569    source.append("#undef CALL_IQ")
570    source.append("#undef KERNEL_NAME")
571
572    # define the Iqxy kernel from the same source with different #defines
573    source.append("#define KERNEL_NAME %s_Iqxy"%model_info['name'])
574    source.append(call_iqxy)
575    source.append(kernel_code)
576    source.append("#undef CALL_IQ")
577    source.append("#undef KERNEL_NAME")
578
579    return '\n'.join(source)
580
581def categorize_parameters(pars):
582    """
583    Categorize the parameters by use:
584
585    * *pd* list of polydisperse parameters in order; gui should test whether
586      they are in *2d* or *magnetic* as appropriate for the data
587    * *1d* set of parameters that are used to compute 1D patterns
588    * *2d* set of parameters that are used to compute 2D patterns (which
589      includes all 1D parameters)
590    * *magnetic* set of parameters that are used to compute magnetic
591      patterns (which includes all 1D and 2D parameters)
592    * *pd_relative* is the set of parameters with relative distribution
593      width (e.g., radius +/- 10%) rather than absolute distribution
594      width (e.g., theta +/- 6 degrees).
595    * *theta_par* is the index of the polar angle polydispersion parameter
596      or -1 if no such parameter exists
597    """
598    par_set = {}
599    par_set['1d'] = [p.name for p in pars if p.type not in ('orientation', 'magnetic')]
600    par_set['2d'] = [p.name for p in pars if p.type != 'magnetic']
601    par_set['magnetic'] = [p.name for p in pars]
602    par_set['pd'] = [p.name for p in pars if p.type in ('volume', 'orientation')]
603    par_set['pd_relative'] = [p.name for p in pars if p.type == 'volume']
604    if 'theta' in par_set['2d']:
605        # find id of theta in parameter set (or whatever variable is
606        # used for spherical normalization during polydispersity...
607        par_set['theta_par'] = [k for k,p in enumerate(pars) if p.name=='theta'][0]
608    else:
609        par_set['theta_par'] = -1
610    return par_set
611
612def collect_types(pars):
613    """
614    Build parameter categories out of the the parameter definitions.
615
616    Returns a dictionary of categories.
617
618    Note: these categories are subject to change, depending on the needs of
619    the UI and the needs of the kernel calling function.
620
621    The categories are as follows:
622
623    * *volume* list of volume parameter names
624    * *orientation* list of orientation parameters
625    * *magnetic* list of magnetic parameters
626    * *sld* list of parameters that have no type info
627    * *other* list of parameters that have no type info
628
629    Each parameter is in one and only one category.
630    """
631    par_type = {
632        'volume': [], 'orientation': [], 'magnetic': [], 'sld': [], 'other': [],
633    }
634    for p in pars:
635        par_type[p.type if p.type else 'other'].append(p.name)
636    return  par_type
637
638
639def process_parameters(model_info):
640    """
641    Process parameter block, precalculating parameter details.
642    """
643    # convert parameters into named tuples
644    for p in model_info['parameters']:
645        if p[4] == '' and (p[0].startswith('sld') or p[0].endswith('sld')):
646            p[4] = 'sld'
647            # TODO: make sure all models explicitly label their sld parameters
648            #raise ValueError("%s.%s needs to be explicitly set to type 'sld'" %(model_info['id'], p[0]))
649
650    pars = [Parameter(*p) for p in model_info['parameters']]
651    # Fill in the derived attributes
652    par_type = collect_types(pars)
653    par_type.update(categorize_parameters(pars))
654    model_info['parameters'] = pars
655    model_info['par_type'] = par_type
656    if model_info.get('demo', None) is None:
657        model_info['demo'] = dict((p.name, p.default) for p in pars)
658    model_info['has_2d'] = par_type['orientation'] or par_type['magnetic']
659    # Don't use more polydisperse parameters than are available in the model
660    # Note: we can do polydispersity on arbitrary parameters, so it is not
661    # clear that this is a good idea; it does however make the poly_details
662    # code easier to write, so we will leave it in for now.
663    model_info['max_pd'] = min(len(par_type['pd']), MAX_PD)
664
665def mono_details(model_info):
666    max_pd = model_info['max_pd']
667    npars = len(model_info['parameters']) - 2
668    p = 5*max_pd
669    c = p + 3*npars
670
671    details = np.zeros(c + 2, 'int32')
672    details[0*max_pd:1*max_pd] = range(max_pd)       # pd_par: arbitrary order; use first
673    details[1*max_pd:2*max_pd] = [1]*max_pd          # pd_length: only one element
674    details[2*max_pd:3*max_pd] = range(max_pd)       # pd_offset: consecutive 1.0 weights
675    details[3*max_pd:4*max_pd] = [1]*max_pd          # pd_stride: vectors of length 1
676    details[4*max_pd:5*max_pd] = [0]*max_pd          # pd_isvol: doens't matter if no norm
677    details[p+0*npars:p+1*npars] = range(2, npars+2) # par_offset: skip scale and background
678    details[p+1*npars:p+2*npars] = [0]*npars         # no coordination
679    #details[p+npars] = 1 # par_coord[0] is coordinated with the first par?
680    details[p+2*npars:p+3*npars] = 0 # fast coord with 0
681    details[c]   = 1     # fast_coord_count: one fast index
682    details[c+1] = -1    # theta_par: None
683    return details
684
685def poly_details(model_info, weights):
686    pars = model_info['parameters'][2:]  # skip scale and background
687    weights = weights[2:]
688    max_pd = model_info['max_pd']
689    npars = len(pars) # scale and background already removed
690    par_offset = 5*max_pd
691    constants_offset = par_offset + 3*npars
692
693    # Decreasing list of polydispersity lengths
694    # Note: the reversing view, x[::-1], does not require a copy
695    pd_length = np.array([len(w) for w in weights])
696    print (pd_length)
697    print (weights)
698    pd_offset = np.cumsum(np.hstack((0, pd_length)))
699    pd_isvol = np.array([p.type=='volume' for p in pars])
700    idx = np.argsort(pd_length)[::-1][:max_pd]
701    print (idx)
702    pd_stride = np.cumprod(np.hstack((1, pd_length[idx][:-1])))
703    par_offsets = np.cumsum(np.hstack((2, pd_length)))[:-1]
704
705    theta_par = -1
706    if 'theta_par' in model_info:
707        theta_par = model_info['theta_par']
708        if theta_par >= 0 and pd_length[theta_par] <= 1:
709            theta_par = -1
710
711    details = np.empty(constants_offset + 2, 'int32')
712    details[0*max_pd:1*max_pd] = idx             # pd_par
713    details[1*max_pd:2*max_pd] = pd_length[idx]
714    details[2*max_pd:3*max_pd] = pd_offset[idx]
715    details[3*max_pd:4*max_pd] = pd_stride
716    details[4*max_pd:5*max_pd] = pd_isvol[idx]
717    details[par_offset+0*npars:par_offset+1*npars] = par_offsets
718    details[par_offset+1*npars:par_offset+2*npars] = 0  # no coordination for most
719    details[par_offset+2*npars:par_offset+3*npars] = 0  # no fast coord with 0
720    coord_offset = par_offset+1*npars
721    for k,parameter_num in enumerate(idx):
722        details[coord_offset+parameter_num] = 2**k
723    details[constants_offset] = 1   # fast_coord_count: one fast index
724    details[constants_offset+1] = theta_par
725    print ("details",details)
726    return details
727
728def constrained_poly_details(model_info, weights, constraints):
729    # Need to find the independently varying pars and sort them
730    # Need to build a coordination list for the dependent variables
731    # Need to generate a constraints function which takes values
732    # and weights, returning par blocks
733    raise NotImplementedError("Can't handle constraints yet")
734
735
736def create_default_functions(model_info):
737    """
738    Autogenerate missing functions, such as Iqxy from Iq.
739
740    This only works for Iqxy when Iq is written in python. :func:`make_source`
741    performs a similar role for Iq written in C.
742    """
743    if model_info['Iq'] is not None and model_info['Iqxy'] is None:
744        if model_info['par_type']['1d'] != model_info['par_type']['2d']:
745            raise ValueError("Iqxy model is missing")
746        Iq = model_info['Iq']
747        def Iqxy(qx, qy, **kw):
748            return Iq(np.sqrt(qx**2 + qy**2), **kw)
749        model_info['Iqxy'] = Iqxy
750
751def make_model_info(kernel_module):
752    """
753    Interpret the model definition file, categorizing the parameters.
754
755    The module can be loaded with a normal python import statement if you
756    know which module you need, or with __import__('sasmodels.model.'+name)
757    if the name is in a string.
758
759    The *model_info* structure contains the following fields:
760
761    * *id* is the id of the kernel
762    * *name* is the display name of the kernel
763    * *filename* is the full path to the module defining the file (if any)
764    * *title* is a short description of the kernel
765    * *description* is a long description of the kernel (this doesn't seem
766      very useful since the Help button on the model page brings you directly
767      to the documentation page)
768    * *docs* is the docstring from the module.  Use :func:`make_doc` to
769    * *category* specifies the model location in the docs
770    * *parameters* is the model parameter table
771    * *single* is True if the model allows single precision
772    * *structure_factor* is True if the model is useable in a product
773    * *variant_info* contains the information required to select between
774      model variants (e.g., the list of cases) or is None if there are no
775      model variants
776    * *par_type* categorizes the model parameters. See
777      :func:`categorize_parameters` for details.
778    * *demo* contains the *{parameter: value}* map used in compare (and maybe
779      for the demo plot, if plots aren't set up to use the default values).
780      If *demo* is not given in the file, then the default values will be used.
781    * *tests* is a set of tests that must pass
782    * *source* is the list of library files to include in the C model build
783    * *Iq*, *Iqxy*, *form_volume*, *ER*, *VR* and *sesans* are python functions
784      implementing the kernel for the module, or None if they are not
785      defined in python
786    * *oldname* is the model name in pre-4.0 Sasview
787    * *oldpars* is the *{new: old}* parameter translation table
788      from pre-4.0 Sasview
789    * *composition* is None if the model is independent, otherwise it is a
790      tuple with composition type ('product' or 'mixture') and a list of
791      *model_info* blocks for the composition objects.  This allows us to
792      build complete product and mixture models from just the info.
793    * *max_pd* is the max polydispersity dimension.  This is constant and
794      should not be reset.  You may be able to change it when the program
795      starts by setting *sasmodels.generate.MAX_PD*.
796
797    """
798    # TODO: maybe turn model_info into a class ModelDefinition
799    parameters = COMMON_PARAMETERS + kernel_module.parameters
800    filename = abspath(kernel_module.__file__)
801    kernel_id = splitext(basename(filename))[0]
802    name = getattr(kernel_module, 'name', None)
803    if name is None:
804        name = " ".join(w.capitalize() for w in kernel_id.split('_'))
805    model_info = dict(
806        id=kernel_id,  # string used to load the kernel
807        filename=abspath(kernel_module.__file__),
808        name=name,
809        title=kernel_module.title,
810        description=kernel_module.description,
811        parameters=parameters,
812        composition=None,
813        docs=kernel_module.__doc__,
814        category=getattr(kernel_module, 'category', None),
815        single=getattr(kernel_module, 'single', True),
816        structure_factor=getattr(kernel_module, 'structure_factor', False),
817        variant_info=getattr(kernel_module, 'invariant_info', None),
818        demo=getattr(kernel_module, 'demo', None),
819        source=getattr(kernel_module, 'source', []),
820        oldname=getattr(kernel_module, 'oldname', None),
821        oldpars=getattr(kernel_module, 'oldpars', {}),
822        tests=getattr(kernel_module, 'tests', []),
823        )
824    process_parameters(model_info)
825    # Check for optional functions
826    functions = "ER VR form_volume Iq Iqxy shape sesans".split()
827    model_info.update((k, getattr(kernel_module, k, None)) for k in functions)
828    create_default_functions(model_info)
829    # Precalculate the monodisperse parameters
830    # TODO: make this a lazy evaluator
831    # make_model_info is called for every model on sasview startup
832    model_info['mono_details'] = mono_details(model_info)
833    return model_info
834
835section_marker = re.compile(r'\A(?P<first>[%s])(?P=first)*\Z'
836                            %re.escape(string.punctuation))
837def _convert_section_titles_to_boldface(lines):
838    """
839    Do the actual work of identifying and converting section headings.
840    """
841    prior = None
842    for line in lines:
843        if prior is None:
844            prior = line
845        elif section_marker.match(line):
846            if len(line) >= len(prior):
847                yield "".join(("**", prior, "**"))
848                prior = None
849            else:
850                yield prior
851                prior = line
852        else:
853            yield prior
854            prior = line
855    if prior is not None:
856        yield prior
857
858def convert_section_titles_to_boldface(s):
859    """
860    Use explicit bold-face rather than section headings so that the table of
861    contents is not polluted with section names from the model documentation.
862
863    Sections are identified as the title line followed by a line of punctuation
864    at least as long as the title line.
865    """
866    return "\n".join(_convert_section_titles_to_boldface(s.split('\n')))
867
868def make_doc(model_info):
869    """
870    Return the documentation for the model.
871    """
872    Iq_units = "The returned value is scaled to units of |cm^-1| |sr^-1|, absolute scale."
873    Sq_units = "The returned value is a dimensionless structure factor, $S(q)$."
874    docs = convert_section_titles_to_boldface(model_info['docs'])
875    subst = dict(id=model_info['id'].replace('_', '-'),
876                 name=model_info['name'],
877                 title=model_info['title'],
878                 parameters=make_partable(model_info['parameters']),
879                 returns=Sq_units if model_info['structure_factor'] else Iq_units,
880                 docs=docs)
881    return DOC_HEADER % subst
882
883
884
885def demo_time():
886    """
887    Show how long it takes to process a model.
888    """
889    from .models import cylinder
890    import datetime
891    tic = datetime.datetime.now()
892    make_source(make_model_info(cylinder))
893    toc = (datetime.datetime.now() - tic).total_seconds()
894    print("time: %g"%toc)
895
896def main():
897    """
898    Program which prints the source produced by the model.
899    """
900    if len(sys.argv) <= 1:
901        print("usage: python -m sasmodels.generate modelname")
902    else:
903        name = sys.argv[1]
904        import sasmodels.models
905        __import__('sasmodels.models.' + name)
906        model = getattr(sasmodels.models, name)
907        model_info = make_model_info(model)
908        source = make_source(model_info)
909        print(source)
910
911if __name__ == "__main__":
912    main()
Note: See TracBrowser for help on using the repository browser.