source: sasmodels/sasmodels/generate.py @ 5ff1b03

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

working kerneldll

  • Property mode set to 100644
File size: 33.2 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.
71The kernel module must set variables defining the kernel meta data:
72
73    *id* is an implicit variable formed from the filename.  It will be
74    a valid python identifier, and will be used as the reference into
75    the html documentation, with '_' replaced by '-'.
76
77    *name* is the model name as displayed to the user.  If it is missing,
78    it will be constructed from the id.
79
80    *title* is a short description of the model, suitable for a tool tip,
81    or a one line model summary in a table of models.
82
83    *description* is an extended description of the model to be displayed
84    while the model parameters are being edited.
85
86    *parameters* is the list of parameters.  Parameters in the kernel
87    functions must appear in the same order as they appear in the
88    parameters list.  Two additional parameters, *scale* and *background*
89    are added to the beginning of the parameter list.  They will show up
90    in the documentation as model parameters, but they are never sent to
91    the kernel functions.  Note that *effect_radius* and *volfraction*
92    must occur first in structure factor calculations.
93
94    *category* is the default category for the model.  The category is
95    two level structure, with the form "group:section", indicating where
96    in the manual the model will be located.  Models are alphabetical
97    within their section.
98
99    *source* is the list of C-99 source files that must be joined to
100    create the OpenCL kernel functions.  The files defining the functions
101    need to be listed before the files which use the functions.
102
103    *ER* is a python function defining the effective radius.  If it is
104    not present, the effective radius is 0.
105
106    *VR* is a python function defining the volume ratio.  If it is not
107    present, the volume ratio is 1.
108
109    *form_volume*, *Iq*, *Iqxy*, *Imagnetic* are strings containing the
110    C source code for the body of the volume, Iq, and Iqxy functions
111    respectively.  These can also be defined in the last source file.
112
113    *Iq* and *Iqxy* also be instead be python functions defining the
114    kernel.  If they are marked as *Iq.vectorized = True* then the
115    kernel is passed the entire *q* vector at once, otherwise it is
116    passed values one *q* at a time.  The performance improvement of
117    this step is significant.
118
119    *demo* is a dictionary of parameter=value defining a set of
120    parameters to use by default when *compare* is called.  Any
121    parameter not set in *demo* gets the initial value from the
122    parameter list.  *demo* is mostly needed to set the default
123    polydispersity values for tests.
124
125    *oldname* is the name of the model in sasview before sasmodels
126    was split into its own package, and *oldpars* is a dictionary
127    of *parameter: old_parameter* pairs defining the new names for
128    the parameters.  This is used by *compare* to check the values
129    of the new model against the values of the old model before
130    you are ready to add the new model to sasmodels.
131
132
133An *model_info* dictionary is constructed from the kernel meta data and
134returned to the caller.
135
136The model evaluator, function call sequence consists of q inputs and the return vector,
137followed by the loop value/weight vector, followed by the values for
138the non-polydisperse parameters, followed by the lengths of the
139polydispersity loops.  To construct the call for 1D models, the
140categories *fixed-1d* and *pd-1d* list the names of the parameters
141of the non-polydisperse and the polydisperse parameters respectively.
142Similarly, *fixed-2d* and *pd-2d* provide parameter names for 2D models.
143The *pd-rel* category is a set of those parameters which give
144polydispersitiy as a portion of the value (so a 10% length dispersity
145would use a polydispersity value of 0.1) rather than absolute
146dispersity such as an angle plus or minus 15 degrees.
147
148The *volume* category lists the volume parameters in order for calls
149to volume within the kernel (used for volume normalization) and for
150calls to ER and VR for effective radius and volume ratio respectively.
151
152The *orientation* and *magnetic* categories list the orientation and
153magnetic parameters.  These are used by the sasview interface.  The
154blank category is for parameters such as scale which don't have any
155other marking.
156
157The doc string at the start of the kernel module will be used to
158construct the model documentation web pages.  Embedded figures should
159appear in the subdirectory "img" beside the model definition, and tagged
160with the kernel module name to avoid collision with other models.  Some
161file systems are case-sensitive, so only use lower case characters for
162file names and extensions.
163
164
165The function :func:`make` loads the metadata from the module and returns
166the kernel source.  The function :func:`make_doc` extracts the doc string
167and adds the parameter table to the top.  The function :func:`model_sources`
168returns a list of files required by the model.
169
170Code follows the C99 standard with the following extensions and conditions::
171
172    M_PI_180 = pi/180
173    M_4PI_3 = 4pi/3
174    square(x) = x*x
175    cube(x) = x*x*x
176    sinc(x) = sin(x)/x, with sin(0)/0 -> 1
177    all double precision constants must include the decimal point
178    all double declarations may be converted to half, float, or long double
179    FLOAT_SIZE is the number of bytes in the converted variables
180"""
181from __future__ import print_function
182
183#TODO: determine which functions are useful outside of generate
184#__all__ = ["model_info", "make_doc", "make_source", "convert_type"]
185
186import sys
187from os.path import abspath, dirname, join as joinpath, exists, basename, \
188    splitext, getmtime
189import re
190import string
191import warnings
192
193import numpy as np
194
195from .modelinfo import ModelInfo, Parameter, make_parameter_table
196
197# TODO: identify model files which have changed since loading and reload them.
198
199TEMPLATE_ROOT = dirname(__file__)
200
201MAX_PD = 4
202
203F16 = np.dtype('float16')
204F32 = np.dtype('float32')
205F64 = np.dtype('float64')
206try:  # CRUFT: older numpy does not support float128
207    F128 = np.dtype('float128')
208except TypeError:
209    F128 = None
210
211# Conversion from units defined in the parameter table for each model
212# to units displayed in the sphinx documentation.
213RST_UNITS = {
214    "Ang": "|Ang|",
215    "1/Ang": "|Ang^-1|",
216    "1/Ang^2": "|Ang^-2|",
217    "1e-6/Ang^2": "|1e-6Ang^-2|",
218    "degrees": "degree",
219    "1/cm": "|cm^-1|",
220    "Ang/cm": "|Ang*cm^-1|",
221    "g/cm3": "|g/cm^3|",
222    "mg/m2": "|mg/m^2|",
223    "": "None",
224    }
225
226# Headers for the parameters tables in th sphinx documentation
227PARTABLE_HEADERS = [
228    "Parameter",
229    "Description",
230    "Units",
231    "Default value",
232    ]
233
234# Minimum width for a default value (this is shorter than the column header
235# width, so will be ignored).
236PARTABLE_VALUE_WIDTH = 10
237
238# Documentation header for the module, giving the model name, its short
239# description and its parameter table.  The remainder of the doc comes
240# from the module docstring.
241DOC_HEADER = """.. _%(id)s:
242
243%(name)s
244=======================================================
245
246%(title)s
247
248%(parameters)s
249
250%(returns)s
251
252%(docs)s
253"""
254
255def format_units(units):
256    """
257    Convert units into ReStructured Text format.
258    """
259    return "string" if isinstance(units, list) else RST_UNITS.get(units, units)
260
261def make_partable(pars):
262    """
263    Generate the parameter table to include in the sphinx documentation.
264    """
265    column_widths = [
266        max(len(p.name) for p in pars),
267        max(len(p.description) for p in pars),
268        max(len(format_units(p.units)) for p in pars),
269        PARTABLE_VALUE_WIDTH,
270        ]
271    column_widths = [max(w, len(h))
272                     for w, h in zip(column_widths, PARTABLE_HEADERS)]
273
274    sep = " ".join("="*w for w in column_widths)
275    lines = [
276        sep,
277        " ".join("%-*s" % (w, h)
278                 for w, h in zip(column_widths, PARTABLE_HEADERS)),
279        sep,
280        ]
281    for p in pars:
282        lines.append(" ".join([
283            "%-*s" % (column_widths[0], p.name),
284            "%-*s" % (column_widths[1], p.description),
285            "%-*s" % (column_widths[2], format_units(p.units)),
286            "%*g" % (column_widths[3], p.default),
287            ]))
288    lines.append(sep)
289    return "\n".join(lines)
290
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:
298        target = joinpath(path, filename)
299        if exists(target):
300            return target
301    raise ValueError("%r not found in %s" % (filename, search_path))
302
303
304def model_sources(model_info):
305    """
306    Return a list of the sources file paths for the module.
307    """
308    search_path = [dirname(model_info['filename']),
309                   abspath(joinpath(dirname(__file__), 'models'))]
310    return [_search(search_path, f) for f in model_info['source']]
311
312def timestamp(model_info):
313    """
314    Return a timestamp for the model corresponding to the most recently
315    changed file or dependency.
316    """
317    source_files = (model_sources(model_info)
318                    + model_templates()
319                    + [model_info['filename']])
320    newest = max(getmtime(f) for f in source_files)
321    return newest
322
323def convert_type(source, dtype):
324    """
325    Convert code from double precision to the desired type.
326
327    Floating point constants are tagged with 'f' for single precision or 'L'
328    for long double precision.
329    """
330    if dtype == F16:
331        fbytes = 2
332        source = _convert_type(source, "float", "f")
333    elif dtype == F32:
334        fbytes = 4
335        source = _convert_type(source, "float", "f")
336    elif dtype == F64:
337        fbytes = 8
338        # no need to convert if it is already double
339    elif dtype == F128:
340        fbytes = 16
341        source = _convert_type(source, "long double", "L")
342    else:
343        raise ValueError("Unexpected dtype in source conversion: %s"%dtype)
344    return ("#define FLOAT_SIZE %d\n"%fbytes)+source
345
346
347def _convert_type(source, type_name, constant_flag):
348    """
349    Replace 'double' with *type_name* in *source*, tagging floating point
350    constants with *constant_flag*.
351    """
352    # Convert double keyword to float/long double/half.
353    # Accept an 'n' # parameter for vector # values, where n is 2, 4, 8 or 16.
354    # Assume complex numbers are represented as cdouble which is typedef'd
355    # to double2.
356    source = re.sub(r'(^|[^a-zA-Z0-9_]c?)double(([248]|16)?($|[^a-zA-Z0-9_]))',
357                    r'\1%s\2'%type_name, source)
358    # Convert floating point constants to single by adding 'f' to the end,
359    # or long double with an 'L' suffix.  OS/X complains if you don't do this.
360    source = re.sub(r'[^a-zA-Z_](\d*[.]\d+|\d+[.]\d*)([eE][+-]?\d+)?',
361                    r'\g<0>%s'%constant_flag, source)
362    return source
363
364
365def kernel_name(model_info, is_2d):
366    """
367    Name of the exported kernel symbol.
368    """
369    return model_info['name'] + "_" + ("Iqxy" if is_2d else "Iq")
370
371
372def indent(s, depth):
373    """
374    Indent a string of text with *depth* additional spaces on each line.
375    """
376    spaces = " "*depth
377    sep = "\n" + spaces
378    return spaces + sep.join(s.split("\n"))
379
380
381_template_cache = {}
382def load_template(filename):
383    path = joinpath(TEMPLATE_ROOT, filename)
384    mtime = getmtime(path)
385    if filename not in _template_cache or mtime > _template_cache[filename][0]:
386        with open(path) as fid:
387            _template_cache[filename] = (mtime, fid.read(), path)
388    return _template_cache[filename][1]
389
390def model_templates():
391    # TODO: fails DRY; templates are listed in two places.
392    # should instead have model_info contain a list of paths
393    return [joinpath(TEMPLATE_ROOT, filename)
394            for filename in ('kernel_header.c', 'kernel_iq.c')]
395
396
397_FN_TEMPLATE = """\
398double %(name)s(%(pars)s);
399double %(name)s(%(pars)s) {
400    %(body)s
401}
402
403
404"""
405
406def _gen_fn(name, pars, body):
407    """
408    Generate a function given pars and body.
409
410    Returns the following string::
411
412         double fn(double a, double b, ...);
413         double fn(double a, double b, ...) {
414             ....
415         }
416    """
417    par_decl = ', '.join(p.as_function_argument() for p in pars) if pars else 'void'
418    return _FN_TEMPLATE % {'name': name, 'body': body, 'pars': par_decl}
419
420def _call_pars(prefix, pars):
421    """
422    Return a list of *prefix.parameter* from parameter items.
423    """
424    return [p.as_call_reference(prefix) for p in pars]
425
426_IQXY_PATTERN = re.compile("^((inline|static) )? *(double )? *Iqxy *([(]|$)",
427                           flags=re.MULTILINE)
428def _have_Iqxy(sources):
429    """
430    Return true if any file defines Iqxy.
431
432    Note this is not a C parser, and so can be easily confused by
433    non-standard syntax.  Also, it will incorrectly identify the following
434    as having Iqxy::
435
436        /*
437        double Iqxy(qx, qy, ...) { ... fill this in later ... }
438        */
439
440    If you want to comment out an Iqxy function, use // on the front of the
441    line instead.
442    """
443    for code in sources:
444        if _IQXY_PATTERN.search(code):
445            return True
446    else:
447        return False
448
449def make_source(model_info):
450    """
451    Generate the OpenCL/ctypes kernel from the module info.
452
453    Uses source files found in the given search path.
454    """
455    if callable(model_info['Iq']):
456        return None
457
458    # TODO: need something other than volume to indicate dispersion parameters
459    # No volume normalization despite having a volume parameter.
460    # Thickness is labelled a volume in order to trigger polydispersity.
461    # May want a separate dispersion flag, or perhaps a separate category for
462    # disperse, but not volume.  Volume parameters also use relative values
463    # for the distribution rather than the absolute values used by angular
464    # dispersion.  Need to be careful that necessary parameters are available
465    # for computing volume even if we allow non-disperse volume parameters.
466
467    partable = model_info['parameters']
468
469    # Identify parameters for Iq, Iqxy, Iq_magnetic and form_volume.
470    # Note that scale and volume are not possible types.
471
472    # Load templates and user code
473    kernel_header = load_template('kernel_header.c')
474    kernel_code = load_template('kernel_iq.c')
475    user_code = [open(f).read() for f in model_sources(model_info)]
476
477    # Build initial sources
478    source = [kernel_header] + user_code
479
480    vol_parameters = partable.kernel_pars('volume')
481    iq_parameters = partable.kernel_pars('1d')
482    iqxy_parameters = partable.kernel_pars('2d')
483
484    # Make parameters for q, qx, qy so that we can use them in declarations
485    q, qx, qy = [Parameter(name=v) for v in ('q', 'qx', 'qy')]
486    # Generate form_volume function, etc. from body only
487    if model_info['form_volume'] is not None:
488        pars = vol_parameters
489        source.append(_gen_fn('form_volume', pars, model_info['form_volume']))
490    if model_info['Iq'] is not None:
491        pars = [q] + iq_parameters
492        source.append(_gen_fn('Iq', pars, model_info['Iq']))
493    if model_info['Iqxy'] is not None:
494        pars = [qx, qy] + iqxy_parameters
495        source.append(_gen_fn('Iqxy', pars, model_info['Iqxy']))
496
497    # Define the parameter table
498    source.append("#define PARAMETER_TABLE \\")
499    source.append("\\\n".join(p.as_definition()
500                                  for p in model_info['parameters'][2:]))
501
502    # Define the function calls
503    if vol_parameters:
504        refs = _call_pars("v.", vol_parameters)
505        call_volume = "#define CALL_VOLUME(v) form_volume(%s)" % (",".join(refs))
506    else:
507        # Model doesn't have volume.  We could make the kernel run a little
508        # faster by not using/transferring the volume normalizations, but
509        # the ifdef's reduce readability more than is worthwhile.
510        call_volume = "#define CALL_VOLUME(v) 1.0"
511    source.append(call_volume)
512
513    refs = ["q[i]"] + _call_pars("v.", iq_parameters)
514    call_iq = "#define CALL_IQ(q,i,v) Iq(%s)" % (",".join(refs))
515    if _have_Iqxy(user_code):
516        # Call 2D model
517        refs = ["q[2*i]", "q[2*i+1]"] + _call_pars("v.", iqxy_parameters)
518        call_iqxy = "#define CALL_IQ(q,i,v) Iqxy(%s)" % (",".join(refs))
519    else:
520        # Call 1D model with sqrt(qx^2 + qy^2)
521        warnings.warn("Creating Iqxy = Iq(sqrt(qx^2 + qy^2))")
522        # still defined:: refs = ["q[i]"] + _call_pars("v", iq_parameters)
523        pars_sqrt = ["sqrt(q[2*i]*q[2*i]+q[2*i+1]*q[2*i+1])"] + refs[1:]
524        call_iqxy = "#define CALL_IQ(q,i,v) Iq(%s)" % (",".join(pars_sqrt))
525
526    # Fill in definitions for numbers of parameters
527    source.append("#define MAX_PD %s"%model_info['max_pd'])
528    source.append("#define NPARS %d"%(len(partable.kernel_pars())))
529
530    # TODO: allow mixed python/opencl kernels?
531
532    # define the Iq kernel
533    source.append("#define KERNEL_NAME %s_Iq"%model_info['name'])
534    source.append(call_iq)
535    source.append(kernel_code)
536    source.append("#undef CALL_IQ")
537    source.append("#undef KERNEL_NAME")
538
539    # define the Iqxy kernel from the same source with different #defines
540    source.append("#define KERNEL_NAME %s_Iqxy"%model_info['name'])
541    source.append(call_iqxy)
542    source.append(kernel_code)
543    source.append("#undef CALL_IQ")
544    source.append("#undef KERNEL_NAME")
545
546    return '\n'.join(source)
547
548def categorize_parameters(pars):
549    """
550    Categorize the parameters by use:
551
552    * *pd* list of polydisperse parameters in order; gui should test whether
553      they are in *2d* or *magnetic* as appropriate for the data
554    * *1d* set of parameters that are used to compute 1D patterns
555    * *2d* set of parameters that are used to compute 2D patterns (which
556      includes all 1D parameters)
557    * *magnetic* set of parameters that are used to compute magnetic
558      patterns (which includes all 1D and 2D parameters)
559    * *pd_relative* is the set of parameters with relative distribution
560      width (e.g., radius +/- 10%) rather than absolute distribution
561      width (e.g., theta +/- 6 degrees).
562    * *theta_par* is the index of the polar angle polydispersion parameter
563      or -1 if no such parameter exists
564    """
565    par_set = {}
566
567def process_parameters(model_info):
568    """
569    Process parameter block, precalculating parameter details.
570    """
571    partable = model_info['parameters']
572    if model_info.get('demo', None) is None:
573        model_info['demo'] = partable.defaults
574
575    # Don't use more polydisperse parameters than are available in the model
576    # Note: we can do polydispersity on arbitrary parameters, so it is not
577    # clear that this is a good idea; it does however make the poly_details
578    # code easier to write, so we will leave it in for now.
579    model_info['max_pd'] = min(partable.num_pd, MAX_PD)
580
581class CoordinationDetails(object):
582    def __init__(self, model_info):
583        max_pd = model_info['max_pd']
584        npars = len(model_info['parameters'].kernel_pars())
585        par_offset = 4*max_pd
586        self.details = np.zeros(par_offset + 3*npars + 4, 'i4')
587
588        # generate views on different parts of the array
589        self._pd_par     = self.details[0*max_pd:1*max_pd]
590        self._pd_length  = self.details[1*max_pd:2*max_pd]
591        self._pd_offset  = self.details[2*max_pd:3*max_pd]
592        self._pd_stride  = self.details[3*max_pd:4*max_pd]
593        self._par_offset = self.details[par_offset+0*npars:par_offset+1*npars]
594        self._par_coord  = self.details[par_offset+1*npars:par_offset+2*npars]
595        self._pd_coord   = self.details[par_offset+2*npars:par_offset+3*npars]
596
597        # theta_par is fixed
598        self.details[-1] = model_info['parameters'].theta_par
599
600    @property
601    def ctypes(self): return self.details.ctypes
602    @property
603    def pd_par(self): return self._pd_par
604    @property
605    def pd_length(self): return self._pd_length
606    @property
607    def pd_offset(self): return self._pd_offset
608    @property
609    def pd_stride(self): return self._pd_stride
610    @property
611    def pd_coord(self): return self._pd_coord
612    @property
613    def par_coord(self): return self._par_coord
614    @property
615    def par_offset(self): return self._par_offset
616    @property
617    def num_coord(self): return self.details[-2]
618    @num_coord.setter
619    def num_coord(self, v): self.details[-2] = v
620    @property
621    def total_pd(self): return self.details[-3]
622    @total_pd.setter
623    def total_pd(self, v): self.details[-3] = v
624    @property
625    def num_active(self): return self.details[-4]
626    @num_active.setter
627    def num_active(self, v): self.details[-4] = v
628
629    def show(self):
630        print("total_pd", self.total_pd)
631        print("num_active", self.num_active)
632        print("pd_par", self.pd_par)
633        print("pd_length", self.pd_length)
634        print("pd_offset", self.pd_offset)
635        print("pd_stride", self.pd_stride)
636        print("par_offsets", self.par_offset)
637        print("num_coord", self.num_coord)
638        print("par_coord", self.par_coord)
639        print("pd_coord", self.pd_coord)
640        print("theta par", self.details[-1])
641
642def mono_details(model_info):
643    details = CoordinationDetails(model_info)
644    # The zero defaults for monodisperse systems are mostly fine
645    details.par_offset[:] = np.arange(2, len(details.par_offset)+2)
646    return details
647
648def poly_details(model_info, weights):
649    weights = weights[2:]
650    max_pd = model_info['max_pd']
651
652    # Decreasing list of polydispersity lengths
653    # Note: the reversing view, x[::-1], does not require a copy
654    pd_length = np.array([len(w) for w in weights])
655    num_active = np.sum(pd_length>1)
656    if num_active > max_pd:
657        raise ValueError("Too many polydisperse parameters")
658
659    pd_offset = np.cumsum(np.hstack((0, pd_length)))
660    idx = np.argsort(pd_length)[::-1][:num_active]
661    par_length = np.array([max(len(w),1) for w in weights])
662    pd_stride = np.cumprod(np.hstack((1, par_length[idx])))
663    par_offsets = np.cumsum(np.hstack((2, par_length)))
664
665    details = CoordinationDetails(model_info)
666    details.pd_par[:num_active] = idx
667    details.pd_length[:num_active] = pd_length[idx]
668    details.pd_offset[:num_active] = pd_offset[idx]
669    details.pd_stride[:num_active] = pd_stride[:-1]
670    details.par_offset[:] = par_offsets[:-1]
671    details.total_pd = pd_stride[-1]
672    details.num_active = num_active
673    # Without constraints coordinated parameters are just the pd parameters
674    details.par_coord[:num_active] = idx
675    details.pd_coord[:num_active] = 2**np.arange(num_active)
676    details.num_coord = num_active
677    #details.show()
678    return details
679
680def constrained_poly_details(model_info, weights, constraints):
681    # Need to find the independently varying pars and sort them
682    # Need to build a coordination list for the dependent variables
683    # Need to generate a constraints function which takes values
684    # and weights, returning par blocks
685    raise NotImplementedError("Can't handle constraints yet")
686
687
688def create_default_functions(model_info):
689    """
690    Autogenerate missing functions, such as Iqxy from Iq.
691
692    This only works for Iqxy when Iq is written in python. :func:`make_source`
693    performs a similar role for Iq written in C.
694    """
695    if callable(model_info['Iq']) and model_info['Iqxy'] is None:
696        partable = model_info['parameters']
697        if partable.type['1d'] != partable.type['2d']:
698            raise ValueError("Iqxy model is missing")
699        Iq = model_info['Iq']
700        def Iqxy(qx, qy, **kw):
701            return Iq(np.sqrt(qx**2 + qy**2), **kw)
702        model_info['Iqxy'] = Iqxy
703
704
705def make_model_info(kernel_module):
706    """
707    Interpret the model definition file, categorizing the parameters.
708
709    The module can be loaded with a normal python import statement if you
710    know which module you need, or with __import__('sasmodels.model.'+name)
711    if the name is in a string.
712
713    The *model_info* structure contains the following fields:
714
715    * *id* is the id of the kernel
716    * *name* is the display name of the kernel
717    * *filename* is the full path to the module defining the file (if any)
718    * *title* is a short description of the kernel
719    * *description* is a long description of the kernel (this doesn't seem
720      very useful since the Help button on the model page brings you directly
721      to the documentation page)
722    * *docs* is the docstring from the module.  Use :func:`make_doc` to
723    * *category* specifies the model location in the docs
724    * *parameters* is the model parameter table
725    * *single* is True if the model allows single precision
726    * *structure_factor* is True if the model is useable in a product
727    * *variant_info* contains the information required to select between
728      model variants (e.g., the list of cases) or is None if there are no
729      model variants
730    * *par_type* categorizes the model parameters. See
731      :func:`categorize_parameters` for details.
732    * *demo* contains the *{parameter: value}* map used in compare (and maybe
733      for the demo plot, if plots aren't set up to use the default values).
734      If *demo* is not given in the file, then the default values will be used.
735    * *tests* is a set of tests that must pass
736    * *source* is the list of library files to include in the C model build
737    * *Iq*, *Iqxy*, *form_volume*, *ER*, *VR* and *sesans* are python functions
738      implementing the kernel for the module, or None if they are not
739      defined in python
740    * *oldname* is the model name in pre-4.0 Sasview
741    * *oldpars* is the *{new: old}* parameter translation table
742      from pre-4.0 Sasview
743    * *composition* is None if the model is independent, otherwise it is a
744      tuple with composition type ('product' or 'mixture') and a list of
745      *model_info* blocks for the composition objects.  This allows us to
746      build complete product and mixture models from just the info.
747    * *max_pd* is the max polydispersity dimension.  This is constant and
748      should not be reset.  You may be able to change it when the program
749      starts by setting *sasmodels.generate.MAX_PD*.
750
751    """
752    # TODO: maybe turn model_info into a class ModelDefinition
753    #print("make parameter table", kernel_module.parameters)
754    parameters = make_parameter_table(kernel_module.parameters)
755    filename = abspath(kernel_module.__file__)
756    kernel_id = splitext(basename(filename))[0]
757    name = getattr(kernel_module, 'name', None)
758    if name is None:
759        name = " ".join(w.capitalize() for w in kernel_id.split('_'))
760    model_info = dict(
761        id=kernel_id,  # string used to load the kernel
762        filename=abspath(kernel_module.__file__),
763        name=name,
764        title=getattr(kernel_module, 'title', name+" model"),
765        description=getattr(kernel_module, 'description', 'no description'),
766        parameters=parameters,
767        composition=None,
768        docs=kernel_module.__doc__,
769        category=getattr(kernel_module, 'category', None),
770        single=getattr(kernel_module, 'single', True),
771        structure_factor=getattr(kernel_module, 'structure_factor', False),
772        variant_info=getattr(kernel_module, 'invariant_info', None),
773        demo=getattr(kernel_module, 'demo', None),
774        source=getattr(kernel_module, 'source', []),
775        oldname=getattr(kernel_module, 'oldname', None),
776        oldpars=getattr(kernel_module, 'oldpars', {}),
777        tests=getattr(kernel_module, 'tests', []),
778        )
779    process_parameters(model_info)
780    # Check for optional functions
781    functions = "ER VR form_volume Iq Iqxy shape sesans".split()
782    model_info.update((k, getattr(kernel_module, k, None)) for k in functions)
783    create_default_functions(model_info)
784    # Precalculate the monodisperse parameters
785    # TODO: make this a lazy evaluator
786    # make_model_info is called for every model on sasview startup
787    model_info['mono_details'] = mono_details(model_info)
788    return model_info
789
790section_marker = re.compile(r'\A(?P<first>[%s])(?P=first)*\Z'
791                            %re.escape(string.punctuation))
792def _convert_section_titles_to_boldface(lines):
793    """
794    Do the actual work of identifying and converting section headings.
795    """
796    prior = None
797    for line in lines:
798        if prior is None:
799            prior = line
800        elif section_marker.match(line):
801            if len(line) >= len(prior):
802                yield "".join(("**", prior, "**"))
803                prior = None
804            else:
805                yield prior
806                prior = line
807        else:
808            yield prior
809            prior = line
810    if prior is not None:
811        yield prior
812
813def convert_section_titles_to_boldface(s):
814    """
815    Use explicit bold-face rather than section headings so that the table of
816    contents is not polluted with section names from the model documentation.
817
818    Sections are identified as the title line followed by a line of punctuation
819    at least as long as the title line.
820    """
821    return "\n".join(_convert_section_titles_to_boldface(s.split('\n')))
822
823def make_doc(model_info):
824    """
825    Return the documentation for the model.
826    """
827    Iq_units = "The returned value is scaled to units of |cm^-1| |sr^-1|, absolute scale."
828    Sq_units = "The returned value is a dimensionless structure factor, $S(q)$."
829    docs = convert_section_titles_to_boldface(model_info['docs'])
830    subst = dict(id=model_info['id'].replace('_', '-'),
831                 name=model_info['name'],
832                 title=model_info['title'],
833                 parameters=make_partable(model_info['parameters']),
834                 returns=Sq_units if model_info['structure_factor'] else Iq_units,
835                 docs=docs)
836    return DOC_HEADER % subst
837
838
839
840def demo_time():
841    """
842    Show how long it takes to process a model.
843    """
844    from .models import cylinder
845    import datetime
846    tic = datetime.datetime.now()
847    make_source(make_model_info(cylinder))
848    toc = (datetime.datetime.now() - tic).total_seconds()
849    print("time: %g"%toc)
850
851def main():
852    """
853    Program which prints the source produced by the model.
854    """
855    if len(sys.argv) <= 1:
856        print("usage: python -m sasmodels.generate modelname")
857    else:
858        name = sys.argv[1]
859        import sasmodels.models
860        __import__('sasmodels.models.' + name)
861        model = getattr(sasmodels.models, name)
862        model_info = make_model_info(model)
863        source = make_source(model_info)
864        print(source)
865
866if __name__ == "__main__":
867    main()
Note: See TracBrowser for help on using the repository browser.