source: sasmodels/sasmodels/generate.py @ f734e7d

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

restructure c code generation for maintainability; extend test harness to allow opencl and ctypes tests

  • Property mode set to 100644
File size: 22.3 KB
Line 
1"""
2SAS model constructor.
3
4Small angle scattering models are defined by a set of kernel functions:
5
6    *Iq(q, p1, p2, ...)* returns the scattering at q for a form with
7    particular dimensions averaged over all orientations.
8
9    *Iqxy(qx, qy, p1, p2, ...)* returns the scattering at qx,qy for a form
10    with particular dimensions for a single orientation.
11
12    *Imagnetic(qx, qy, result[], p1, p2, ...)* returns the scattering for the
13    polarized neutron spin states (up-up, up-down, down-up, down-down) for
14    a form with particular dimensions for a single orientation.
15
16    *form_volume(p1, p2, ...)* returns the volume of the form with particular
17    dimension.
18
19    *ER(p1, p2, ...)* returns the effective radius of the form with
20    particular dimensions.
21
22    *VR(p1, p2, ...)* returns the volume ratio for core-shell style forms.
23
24These functions are defined in a kernel module .py script and an associated
25set of .c files.  The model constructor will use them to create models with
26polydispersity across volume and orientation parameters, and provide
27scale and background parameters for each model.
28
29*Iq*, *Iqxy*, *Imagnetic* and *form_volume* should be stylized C-99
30functions written for OpenCL.  All functions need prototype declarations
31even if the are defined before they are used.  OpenCL does not support
32*#include* preprocessor directives, so instead the list of includes needs
33to be given as part of the metadata in the kernel module definition.
34The included files should be listed using a path relative to the kernel
35module, or if using "lib/file.c" if it is one of the standard includes
36provided with the sasmodels source.  The includes need to be listed in
37order so that functions are defined before they are used.
38
39Floating point values should be declared as *double*.  For single precision
40calculations, *double* will be replaced by *float*.  The single precision
41conversion will also tag floating point constants with "f" to make them
42single precision constants.  When using integral values in floating point
43expressions, they should be expressed as floating point values by including
44a decimal point.  This includes 0., 1. and 2.
45
46OpenCL has a *sincos* function which can improve performance when both
47the *sin* and *cos* values are needed for a particular argument.  Since
48this function does not exist in C99, all use of *sincos* should be
49replaced by the macro *SINCOS(value,sn,cn)* where *sn* and *cn* are
50previously declared *double* variables.  When compiled for systems without
51OpenCL, *SINCOS* will be replaced by *sin* and *cos* calls.   If *value* is
52an expression, it will appear twice in this case; whether or not it will be
53evaluated twice depends on the quality of the compiler.
54
55If the input parameters are invalid, the scattering calculator should
56return a negative number. Particularly with polydispersity, there are
57some sets of shape parameters which lead to nonsensical forms, such
58as a capped cylinder where the cap radius is smaller than the
59cylinder radius.  The polydispersity calculation will ignore these points,
60effectively chopping the parameter weight distributions at the boundary
61of the infeasible region.  The resulting scattering will be set to
62background.  This will work correctly even when polydispersity is off.
63
64*ER* and *VR* are python functions which operate on parameter vectors.
65The constructor code will generate the necessary vectors for computing
66them with the desired polydispersity.
67
68The available kernel parameters are defined as a list, with each parameter
69defined as a sublist with the following elements:
70
71    *name* is the name that will be used in the call to the kernel
72    function and the name that will be displayed to the user.  Names
73    should be lower case, with words separated by underscore.  If
74    acronyms are used, the whole acronym should be upper case.
75
76    *units* should be one of *degrees* for angles, *Ang* for lengths,
77    *1e-6/Ang^2* for SLDs.
78
79    *default value* will be the initial value for  the model when it
80    is selected, or when an initial value is not otherwise specified.
81
82    [*lb*, *ub*] are the hard limits on the parameter value, used to limit
83    the polydispersity density function.  In the fit, the parameter limits
84    given to the fit are the limits  on the central value of the parameter.
85    If there is polydispersity, it will evaluate parameter values outside
86    the fit limits, but not outside the hard limits specified in the model.
87    If there are no limits, use +/-inf imported from numpy.
88
89    *type* indicates how the parameter will be used.  "volume" parameters
90    will be used in all functions.  "orientation" parameters will be used
91    in *Iqxy* and *Imagnetic*.  "magnetic* parameters will be used in
92    *Imagnetic* only.  If *type* is the empty string, the parameter will
93    be used in all of *Iq*, *Iqxy* and *Imagnetic*.
94
95    *description* is a short description of the parameter.  This will
96    be displayed in the parameter table and used as a tool tip for the
97    parameter value in the user interface.
98
99The kernel module must set variables defining the kernel meta data:
100
101    *name* is the model name
102
103    *title* is a short description of the model, suitable for a tool tip,
104    or a one line model summary in a table of models.
105
106    *description* is an extended description of the model to be displayed
107    while the model parameters are being edited.
108
109    *parameters* is the list of parameters.  Parameters in the kernel
110    functions must appear in the same order as they appear in the
111    parameters list.  Two additional parameters, *scale* and *background*
112    are added to the beginning of the parameter list.  They will show up
113    in the documentation as model parameters, but they are never sent to
114    the kernel functions.
115
116    *source* is the list of C-99 source files that must be joined to
117    create the OpenCL kernel functions.  The files defining the functions
118    need to be listed before the files which use the functions.
119
120    *ER* is a python function defining the effective radius.  If it is
121    not present, the effective radius is 0.
122
123    *VR* is a python function defining the volume ratio.  If it is not
124    present, the volume ratio is 1.
125
126    *form_volume*, *Iq*, *Iqxy*, *Imagnetic* are strings containing the
127    C source code for the body of the volume, Iq, and Iqxy functions
128    respectively.  These can also be defined in the last source file.
129
130    *Iq* and *Iqxy* also be instead be python functions defining the
131    kernel.  If they are marked as *Iq.vectorized = True* then the
132    kernel is passed the entire *q* vector at once, otherwise it is
133    passed values one *q* at a time.  The performance improvement of
134    this step is significant.
135
136An *info* dictionary is constructed from the kernel meta data and
137returned to the caller.
138
139Additional fields can be defined in the kernel definition file that
140are not needed for sas modelling.
141
142    *demo* is a dictionary of parameter=value defining a set of
143    parameters to use by default when *compare* is called.
144
145    *oldname* is the name of the model in sasview before sasmodels
146    was split into its own package, and *oldpars* is a dictionary
147    of *parameter: old_parameter* pairs defining the new names for
148    the parameters.  This is used by *compare* to check the values
149    of the new model against the values of the old model before
150    you are ready to add the new model to sasmodels.
151
152The model evaluator, function call sequence consists of q inputs and the return vector,
153followed by the loop value/weight vector, followed by the values for
154the non-polydisperse parameters, followed by the lengths of the
155polydispersity loops.  To construct the call for 1D models, the
156categories *fixed-1d* and *pd-1d* list the names of the parameters
157of the non-polydisperse and the polydisperse parameters respectively.
158Similarly, *fixed-2d* and *pd-2d* provide parameter names for 2D models.
159The *pd-rel* category is a set of those parameters which give
160polydispersitiy as a portion of the value (so a 10% length dispersity
161would use a polydispersity value of 0.1) rather than absolute
162dispersity such as an angle plus or minus 15 degrees.
163
164The *volume* category lists the volume parameters in order for calls
165to volume within the kernel (used for volume normalization) and for
166calls to ER and VR for effective radius and volume ratio respectively.
167
168The *orientation* and *magnetic* categories list the orientation and
169magnetic parameters.  These are used by the sasview interface.  The
170blank category is for parameters such as scale which don't have any
171other marking.
172
173The doc string at the start of the kernel module will be used to
174construct the model documentation web pages.  Embedded figures should
175appear in the subdirectory "img" beside the model definition, and tagged
176with the kernel module name to avoid collision with other models.  Some
177file systems are case-sensitive, so only use lower case characters for
178file names and extensions.
179
180
181The function :func:`make` loads the metadata from the module and returns
182the kernel source.  The function :func:`doc` extracts the doc string
183and adds the parameter table to the top.  The function :func:`sources`
184returns a list of files required by the model.
185"""
186
187# TODO: identify model files which have changed since loading and reload them.
188
189__all__ = ["make", "doc", "sources", "use_single"]
190
191import sys
192from os.path import abspath, dirname, join as joinpath, exists
193import re
194
195import numpy as np
196C_KERNEL_TEMPLATE_PATH = joinpath(dirname(__file__), 'kernel_template.c')
197
198F64 = np.dtype('float64')
199F32 = np.dtype('float32')
200
201# Scale and background, which are parameters common to every form factor
202COMMON_PARAMETERS = [
203    [ "scale", "", 1, [0, np.inf], "", "Source intensity" ],
204    [ "background", "1/cm", 0, [0, np.inf], "", "Source background" ],
205    ]
206
207
208# Conversion from units defined in the parameter table for each model
209# to units displayed in the sphinx documentation.
210RST_UNITS = {
211    "Ang": "|Ang|",
212    "1/Ang^2": "|Ang^-2|",
213    "1e-6/Ang^2": "|1e-6Ang^-2|",
214    "degrees": "degree",
215    "1/cm": "|cm^-1|",
216    "": "None",
217    }
218
219# Headers for the parameters tables in th sphinx documentation
220PARTABLE_HEADERS = [
221    "Parameter",
222    "Description",
223    "Units",
224    "Default value",
225    ]
226
227# Minimum width for a default value (this is shorter than the column header
228# width, so will be ignored).
229PARTABLE_VALUE_WIDTH = 10
230
231# Documentation header for the module, giving the model name, its short
232# description and its parameter table.  The remainder of the doc comes
233# from the module docstring.
234DOC_HEADER=""".. _%(name)s:
235
236%(label)s
237=======================================================
238
239%(title)s
240
241%(parameters)s
242
243The returned value is scaled to units of |cm^-1|.
244
245%(docs)s
246"""
247
248def make_partable(pars):
249    """
250    Generate the parameter table to include in the sphinx documentation.
251    """
252    pars = COMMON_PARAMETERS + pars
253    column_widths = [
254        max(len(p[0]) for p in pars),
255        max(len(p[-1]) for p in pars),
256        max(len(RST_UNITS[p[1]]) for p in pars),
257        PARTABLE_VALUE_WIDTH,
258        ]
259    column_widths = [max(w, len(h))
260                     for w,h in zip(column_widths, PARTABLE_HEADERS)]
261
262    sep = " ".join("="*w for w in column_widths)
263    lines = [
264        sep,
265        " ".join("%-*s"%(w,h) for w,h in zip(column_widths, PARTABLE_HEADERS)),
266        sep,
267        ]
268    for p in pars:
269        lines.append(" ".join([
270            "%-*s"%(column_widths[0],p[0]),
271            "%-*s"%(column_widths[1],p[-1]),
272            "%-*s"%(column_widths[2],RST_UNITS[p[1]]),
273            "%*g"%(column_widths[3],p[2]),
274            ]))
275    lines.append(sep)
276    return "\n".join(lines)
277
278def _search(search_path, filename):
279    """
280    Find *filename* in *search_path*.
281
282    Raises ValueError if file does not exist.
283    """
284    for path in search_path:
285        target = joinpath(path, filename)
286        if exists(target):
287            return target
288    raise ValueError("%r not found in %s"%(filename, search_path))
289
290def sources(info):
291    """
292    Return a list of the sources file paths for the module.
293    """
294    search_path = [ dirname(info['filename']),
295                    abspath(joinpath(dirname(__file__),'models')) ]
296    return [_search(search_path, f) for f in info['source']]
297
298def use_single(source):
299    """
300    Convert code from double precision to single precision.
301    """
302    # Convert double keyword to float.  Accept an 'n' parameter for vector
303    # values, where n is 2, 4, 8 or 16. Assume complex numbers are represented
304    # as cdouble which is typedef'd to double2.
305    source = re.sub(r'(^|[^a-zA-Z0-9_]c?)double(([248]|16)?($|[^a-zA-Z0-9_]))',
306                    r'\1float\2', source)
307    # Convert floating point constants to single by adding 'f' to the end.
308    # OS/X driver complains if you don't do this.
309    source = re.sub(r'[^a-zA-Z_](\d*[.]\d+|\d+[.]\d*)([eE][+-]?\d+)?',
310                    r'\g<0>f', source)
311    return source
312
313
314def kernel_name(info, is_2D):
315    """
316    Name of the exported kernel symbol.
317    """
318    return info['name'] + "_" + ("Iqxy" if is_2D else "Iq")
319
320
321def categorize_parameters(pars):
322    """
323    Build parameter categories out of the the parameter definitions.
324
325    Returns a dictionary of categories.
326    """
327    partype = {
328        'volume': [], 'orientation': [], 'magnetic': [], '': [],
329        'fixed-1d': [], 'fixed-2d': [], 'pd-1d': [], 'pd-2d': [],
330        'pd-rel': set(),
331    }
332
333    for p in pars:
334        name,ptype = p[0],p[4]
335        if ptype == 'volume':
336            partype['pd-1d'].append(name)
337            partype['pd-2d'].append(name)
338            partype['pd-rel'].add(name)
339        elif ptype == 'magnetic':
340            partype['fixed-2d'].append(name)
341        elif ptype == 'orientation':
342            partype['pd-2d'].append(name)
343        elif ptype == '':
344            partype['fixed-1d'].append(name)
345            partype['fixed-2d'].append(name)
346        else:
347            raise ValueError("unknown parameter type %r"%ptype)
348        partype[ptype].append(name)
349
350    return partype
351
352def indent(s, depth):
353    """
354    Indent a string of text with *depth* additional spaces on each line.
355    """
356    spaces = " "*depth
357    sep = "\n"+spaces
358    return spaces + sep.join(s.split("\n"))
359
360
361def build_polydispersity_loops(pd_pars):
362    """
363    Build polydispersity loops
364
365    Returns loop opening and loop closing
366    """
367    LOOP_OPEN="""\
368for (int %(name)s_i=0; %(name)s_i < N%(name)s; %(name)s_i++) {
369  const double %(name)s = loops[2*(%(name)s_i%(offset)s)];
370  const double %(name)s_w = loops[2*(%(name)s_i%(offset)s)+1];\
371"""
372    depth = 4
373    offset = ""
374    loop_head = []
375    loop_end = []
376    for name in pd_pars:
377        subst = { 'name': name, 'offset': offset }
378        loop_head.append(indent(LOOP_OPEN%subst, depth))
379        loop_end.insert(0, (" "*depth) + "}")
380        offset += '+N'+name
381        depth += 2
382    return "\n".join(loop_head), "\n".join(loop_end)
383
384C_KERNEL_TEMPLATE=None
385def make_model(info):
386    """
387    Generate the code for the kernel defined by info, using source files
388    found in the given search path.
389    """
390    # TODO: need something other than volume to indicate dispersion parameters
391    # No volume normalization despite having a volume parameter.
392    # Thickness is labelled a volume in order to trigger polydispersity.
393    # May want a separate dispersion flag, or perhaps a separate category for
394    # disperse, but not volume.  Volume parameters also use relative values
395    # for the distribution rather than the absolute values used by angular
396    # dispersion.  Need to be careful that necessary parameters are available
397    # for computing volume even if we allow non-disperse volume parameters.
398
399    # Load template
400    global C_KERNEL_TEMPLATE
401    if C_KERNEL_TEMPLATE is None:
402        with open(C_KERNEL_TEMPLATE_PATH) as fid:
403            C_KERNEL_TEMPLATE = fid.read()
404
405    # Load additional sources
406    source = [open(f).read() for f in sources(info)]
407
408    # Prepare defines
409    defines = []
410    partype = info['partype']
411    pd_1d = partype['pd-1d']
412    pd_2d = partype['pd-2d']
413    fixed_1d = partype['fixed-1d']
414    fixed_2d = partype['fixed-1d']
415
416    iq_parameters = [p[0]
417        for p in info['parameters'][2:] # skip scale, background
418        if p[0] in set(fixed_1d+pd_1d)]
419    iqxy_parameters = [p[0]
420        for p in info['parameters'][2:] # skip scale, background
421        if p[0] in set(fixed_2d+pd_2d)]
422    volume_parameters = [p[0]
423        for p in info['parameters']
424        if p[4]=='volume']
425
426    # Fill in defintions for volume parameters
427    if volume_parameters:
428        defines.append(('VOLUME_PARAMETERS',
429                        ','.join(volume_parameters)))
430        defines.append(('VOLUME_WEIGHT_PRODUCT',
431                        '*'.join(p+'_w' for p in volume_parameters)))
432
433    # Generate form_volume function from body only
434    if info['form_volume'] is not None:
435        defines.append(('VOLUME_PARAMETER_DECLARATIONS',
436                        ', '.join('double '+p for p in volume_parameters)))
437        fn = """\
438double form_volume(VOLUME_PARAMETER_DECLARATIONS);
439double form_volume(VOLUME_PARAMETER_DECLARATIONS) {
440    %(body)s
441}
442"""%{'body':info['form_volume']}
443        source.append(fn)
444
445    # Fill in definitions for Iq parameters
446    defines.append(('IQ_KERNEL_NAME', info['name']+'_Iq'))
447    defines.append(('IQ_PARAMETERS', ', '.join(iq_parameters)))
448    if fixed_1d:
449        defines.append(('IQ_FIXED_PARAMETER_DECLARATIONS',
450                        ', \\\n    '.join('const double %s'%p for p in fixed_1d)))
451    if pd_1d:
452        defines.append(('IQ_WEIGHT_PRODUCT',
453                        '*'.join(p+'_w' for p in pd_1d)))
454        defines.append(('IQ_DISPERSION_LENGTH_DECLARATIONS',
455                        ', \\\n    '.join('const int N%s'%p for p in pd_1d)))
456        defines.append(('IQ_DISPERSION_LENGTH_SUM',
457                        '+'.join('N'+p for p in pd_1d)))
458        open_loops, close_loops = build_polydispersity_loops(pd_1d)
459        defines.append(('IQ_OPEN_LOOPS',
460                        open_loops.replace('\n',' \\\n')))
461        defines.append(('IQ_CLOSE_LOOPS',
462                        close_loops.replace('\n',' \\\n')))
463    if info['Iq'] is not None:
464        defines.append(('IQ_PARAMETER_DECLARATIONS',
465                       ', '.join('double '+p for p in iq_parameters)))
466        fn = """\
467double Iq(double q, IQ_PARAMETER_DECLARATIONS);
468double Iq(double q, IQ_PARAMETER_DECLARATIONS) {
469    %(body)s
470}
471"""%{'body':info['Iq']}
472        source.append(fn)
473
474    # Fill in definitions for Iqxy parameters
475    defines.append(('IQXY_KERNEL_NAME', info['name']+'_Iqxy'))
476    defines.append(('IQXY_PARAMETERS', ', '.join(iqxy_parameters)))
477    if fixed_2d:
478        defines.append(('IQXY_FIXED_PARAMETER_DECLARATIONS',
479                        ', \\\n    '.join('const double %s'%p for p in fixed_2d)))
480    if pd_2d:
481        defines.append(('IQXY_WEIGHT_PRODUCT',
482                        '*'.join(p+'_w' for p in pd_2d)))
483        defines.append(('IQXY_DISPERSION_LENGTH_DECLARATIONS',
484                        ', \\\n    '.join('const int N%s'%p for p in pd_2d)))
485        defines.append(('IQXY_DISPERSION_LENGTH_SUM',
486                        '+'.join('N'+p for p in pd_2d)))
487        open_loops, close_loops = build_polydispersity_loops(pd_2d)
488        defines.append(('IQXY_OPEN_LOOPS',
489                        open_loops.replace('\n',' \\\n')))
490        defines.append(('IQXY_CLOSE_LOOPS',
491                        close_loops.replace('\n',' \\\n')))
492    if info['Iqxy'] is not None:
493        defines.append(('IQXY_PARAMETER_DECLARATIONS',
494                       ', '.join('double '+p for p in iqxy_parameters)))
495        fn = """\
496double Iqxy(double qx, double qy, IQXY_PARAMETER_DECLARATIONS);
497double Iqxy(double qx, double qy, IQXY_PARAMETER_DECLARATIONS) {
498    %(body)s
499}
500"""%{'body':info['Iqxy']}
501        source.append(fn)
502
503    # Need to know if we have a theta parameter for Iqxy; it is not there
504    # for the magnetic sphere model, for example, which has a magnetic
505    # orientation but no shape orientation.
506    if 'theta' in pd_2d:
507        defines.append(('IQXY_HAS_THETA', '1'))
508
509    #for d in defines: print d
510    DEFINES='\n'.join('#define %s %s'%(k,v) for k,v in defines)
511    SOURCES='\n\n'.join(source)
512    return C_KERNEL_TEMPLATE%{
513        'DEFINES':DEFINES,
514        'SOURCES':SOURCES,
515        }
516
517def make(kernel_module):
518    """
519    Build an OpenCL/ctypes function from the definition in *kernel_module*.
520
521    The module can be loaded with a normal python import statement if you
522    know which module you need, or with __import__('sasmodels.model.'+name)
523    if the name is in a string.
524    """
525    # TODO: allow Iq and Iqxy to be defined in python
526    #print kernelfile
527    info = dict(
528        filename = abspath(kernel_module.__file__),
529        name = kernel_module.name,
530        title = kernel_module.title,
531        description = kernel_module.description,
532        parameters = COMMON_PARAMETERS + kernel_module.parameters,
533        source = getattr(kernel_module, 'source', []),
534        oldname = kernel_module.oldname,
535        oldpars = kernel_module.oldpars,
536        )
537    # Fill in attributes which default to None
538    info.update((k,getattr(kernel_module, k, None))
539                for k in ('ER', 'VR', 'form_volume', 'Iq', 'Iqxy'))
540    # Fill in the derived attributes
541    info['limits'] = dict((p[0],p[3]) for p in info['parameters'])
542    info['partype'] = categorize_parameters(info['parameters'])
543    info['defaults'] = dict((p[0],p[2]) for p in info['parameters'])
544
545    # Assume if one part of the kernel is python then all parts are.
546    source = make_model(info) if not callable(info['Iq']) else None
547    return source, info
548
549def doc(kernel_module):
550    """
551    Return the documentation for the model.
552    """
553    subst = dict(name=kernel_module.name.replace('_','-'),
554                 label=" ".join(kernel_module.name.split('_')).capitalize(),
555                 title=kernel_module.title,
556                 parameters=make_partable(kernel_module.parameters),
557                 docs=kernel_module.__doc__)
558    return DOC_HEADER%subst
559
560
561
562def demo_time():
563    import datetime
564    from .models import cylinder
565    toc = lambda: (datetime.datetime.now()-tic).total_seconds()
566    tic = datetime.datetime.now()
567    source, info = make(cylinder)
568    print "time:",toc()
569
570def main():
571    if len(sys.argv) <= 1:
572        print "usage: python -m sasmodels.generate modelname"
573    else:
574        name = sys.argv[1]
575        import sasmodels.models
576        __import__('sasmodels.models.'+name)
577        model = getattr(sasmodels.models, name)
578        source, info = make(model); print source
579        #print doc(model)
580
581if __name__ == "__main__":
582    main()
Note: See TracBrowser for help on using the repository browser.