source: sasmodels/doc/guide/plugin.rst @ b15f768

Last change on this file since b15f768 was 0a9fcab, checked in by Paul Kienzle <pkienzle@…>, 7 years ago

simplify handling of gauss.z in py2c

  • Property mode set to 100644
File size: 41.2 KB
RevLine 
[990d8df]1.. _Writing_a_Plugin:
2
3Writing a Plugin Model
4======================
5
6Overview
7^^^^^^^^
8
9In addition to the models provided with the sasmodels package, you are free to
10create your own models.
11
12Models can be of three types:
13
14- A pure python model : Example -
15  `broadpeak.py <https://github.com/SasView/sasmodels/blob/master/sasmodels/models/broad_peak.py>`_
16
17- A python model with embedded C : Example -
18  `sphere.py <https://github.com/SasView/sasmodels/blob/master/sasmodels/models/sphere.py>`_
19
20- A python wrapper with separate C code : Example -
21  `cylinder.py <https://github.com/SasView/sasmodels/blob/master/sasmodels/models/cylinder.py>`_,
22  `cylinder.c <https://github.com/SasView/sasmodels/blob/master/sasmodels/models/cylinder.c>`_
23
24When using SasView, plugin models should be saved to the SasView
25*plugin_models* folder *C:\\Users\\{username}\\.sasview\\plugin_models*
26(on Windows) or */Users/{username}/.sasview\\plugin_models* (on Mac).
27The next time SasView is started it will compile the plugin and add
28it to the list of *Plugin Models* in a FitPage.  Scripts can load
29the models from anywhere.
30
31The built-in modules are available in the *models* subdirectory
32of the sasmodels package.  For SasView on Windows, these will
33be found in *C:\\Program Files (x86)\\SasView\\sasmodels-data\\models*.
34On Mac OSX, these will be within the application bundle as
35*/Applications/SasView 4.0.app/Contents/Resources/sasmodels-data/models*.
36
37Other models are available for download from the
38`Model Marketplace <http://marketplace.sasview.org/>`_. You can contribute your
39own models to the Marketplace as well.
40
41Create New Model Files
42^^^^^^^^^^^^^^^^^^^^^^
43
44Copy the appropriate files to your plugin models directory (we recommend
45using the examples above as templates) as mymodel.py (and mymodel.c, etc)
46as required, where "mymodel" is the name for the model you are creating.
47
48*Please follow these naming rules:*
49
50- No capitalization and thus no CamelCase
51- If necessary use underscore to separate words (i.e. barbell not BarBell or
52  broad_peak not BroadPeak)
53- Do not include "model" in the name (i.e. barbell not BarBellModel)
54
55
56Edit New Model Files
57^^^^^^^^^^^^^^^^^^^^
58
59Model Contents
60..............
61
62The model interface definition is in the .py file.  This file contains:
63
64- a **model name**:
65   - this is the **name** string in the *.py* file
66   - titles should be:
67
68    - all in *lower* case
69    - without spaces (use underscores to separate words instead)
70    - without any capitalization or CamelCase
71    - without incorporating the word "model"
72    - examples: *barbell* **not** *BarBell*; *broad_peak* **not** *BroadPeak*;
73      *barbell* **not** *BarBellModel*
74
75- a **model title**:
76   - this is the **title** string in the *.py* file
77   - this is a one or two line description of the model, which will appear
78     at the start of the model documentation and as a tooltip in the SasView GUI
79
[3048ec6]80- a **short description**:
[990d8df]81   - this is the **description** string in the *.py* file
82   - this is a medium length description which appears when you click
83     *Description* on the model FitPage
84
85- a **parameter table**:
86   - this will be auto-generated from the *parameters* in the *.py* file
87
88- a **long description**:
89   - this is ReStructuredText enclosed between the r""" and """ delimiters
90     at the top of the *.py* file
91   - what you write here is abstracted into the SasView help documentation
92   - this is what other users will refer to when they want to know what
93     your model does; so please be helpful!
94
95- a **definition** of the model:
96   - as part of the **long description**
97
98- a **formula** defining the function the model calculates:
99   - as part of the **long description**
100
101- an **explanation of the parameters**:
102   - as part of the **long description**
103   - explaining how the symbols in the formula map to the model parameters
104
105- a **plot** of the function, with a **figure caption**:
106   - this is automatically generated from your default parameters
107
108- at least one **reference**:
109   - as part of the **long description**
110   - specifying where the reader can obtain more information about the model
111
112- the **name of the author**
113   - as part of the **long description**
114   - the *.py* file should also contain a comment identifying *who*
115     converted/created the model file
116
117Models that do not conform to these requirements will *never* be incorporated
118into the built-in library.
119
120
121Model Documentation
122...................
123
124The *.py* file starts with an r (for raw) and three sets of quotes
125to start the doc string and ends with a second set of three quotes.
126For example::
127
128    r"""
129    Definition
130    ----------
131
132    The 1D scattering intensity of the sphere is calculated in the following
133    way (Guinier, 1955)
134
135    .. math::
136
137        I(q) = \frac{\text{scale}}{V} \cdot \left[
138            3V(\Delta\rho) \cdot \frac{\sin(qr) - qr\cos(qr))}{(qr)^3}
139            \right]^2 + \text{background}
140
141    where *scale* is a volume fraction, $V$ is the volume of the scatterer,
142    $r$ is the radius of the sphere and *background* is the background level.
143    *sld* and *sld_solvent* are the scattering length densities (SLDs) of the
144    scatterer and the solvent respectively, whose difference is $\Delta\rho$.
145
146    You can included figures in your documentation, as in the following
147    figure for the cylinder model.
148
149    .. figure:: img/cylinder_angle_definition.jpg
150
151        Definition of the angles for oriented cylinders.
152
153    References
154    ----------
155
156    A Guinier, G Fournet, *Small-Angle Scattering of X-Rays*,
157    John Wiley and Sons, New York, (1955)
158    """
159
160This is where the FULL documentation for the model goes (to be picked up by
161the automatic documentation system).  Although it feels odd, you
162should start the documentation immediately with the **definition**---the model
163name, a brief description and the parameter table are automatically inserted
164above the definition, and the a plot of the model is automatically inserted
165before the **reference**.
166
167Figures can be included using the *figure* command, with the name
168of the *.png* file containing the figure and a caption to appear below the
169figure.  Figure numbers will be added automatically.
170
171See this `Sphinx cheat sheet <http://matplotlib.org/sampledoc/cheatsheet.html>`_
172for a quick guide to the documentation layout commands, or the
173`Sphinx Documentation <http://www.sphinx-doc.org/en/stable/>`_ for
174complete details.
175
176The model should include a **formula** written using LaTeX markup.
177The example above uses the *math* command to make a displayed equation.  You
178can also use *\$formula\$* for an inline formula. This is handy for defining
179the relationship between the model parameters and formula variables, such
180as the phrase "\$r\$ is the radius" used above.  The live demo MathJax
181page `<http://www.mathjax.org/>`_ is handy for checking that the equations
182will look like you intend.
183
184Math layout uses the `amsmath <http://www.ams.org/publications/authors/tex/amslatex>`_
185package for aligning equations (see amsldoc.pdf on that page for complete
186documentation). You will automatically be in an aligned environment, with
187blank lines separating the lines of the equation.  Place an ampersand before
188the operator on which to align.  For example::
189
190    .. math::
191
192      x + y &= 1 \\
193      y &= x - 1
194
195produces
196
197.. math::
198
199      x + y &= 1 \\
200      y &= x - 1
201
202If you need more control, use::
203
204    .. math::
205        :nowrap:
206
207
208Model Definition
209................
210
211Following the documentation string, there are a series of definitions::
212
213    name = "sphere"  # optional: defaults to the filename without .py
214
215    title = "Spheres with uniform scattering length density"
216
217    description = """\
218    P(q)=(scale/V)*[3V(sld-sld_solvent)*(sin(qr)-qr cos(qr))
219                    /(qr)^3]^2 + background
220        r: radius of sphere
221        V: The volume of the scatter
222        sld: the SLD of the sphere
223        sld_solvent: the SLD of the solvent
224    """
225
226    category = "shape:sphere"
227
228    single = True   # optional: defaults to True
229
230    opencl = False  # optional: defaults to False
231
232    structure_factor = False  # optional: defaults to False
233
234**name = "mymodel"** defines the name of the model that is shown to the user.
[3048ec6]235If it is not provided it will use the name of the model file. The name must
236be a valid variable name, starting with a letter and contains only letters,
237numbers or underscore.  Spaces, dashes, and other symbols are not permitted.
[990d8df]238
239**title = "short description"** is short description of the model which
240is included after the model name in the automatically generated documentation.
241The title can also be used for a tooltip.
242
243**description = """doc string"""** is a longer description of the model. It
244shows up when you press the "Description" button of the SasView FitPage.
245It should give a brief description of the equation and the parameters
246without the need to read the entire model documentation. The triple quotes
247allow you to write the description over multiple lines. Keep the lines
248short since the GUI will wrap each one separately if they are too long.
249**Make sure the parameter names in the description match the model definition!**
250
251**category = "shape:sphere"** defines where the model will appear in the
252model documentation.  In this example, the model will appear alphabetically
253in the list of spheroid models in the *Shape* category.
254
255**single = True** indicates that the model can be run using single
256precision floating point values.  Set it to False if the numerical
257calculation for the model is unstable, which is the case for about 20 of
258the built in models.  It is worthwhile modifying the calculation to support
259single precision, allowing models to run up to 10 times faster.  The
260section `Test_Your_New_Model`_  describes how to compare model values for
261single vs. double precision so you can decide if you need to set
262single to False.
263
264**opencl = False** indicates that the model should not be run using OpenCL.
265This may be because the model definition includes code that cannot be
266compiled for the GPU (for example, goto statements).  It can also be used
267for large models which can't run on most GPUs.  This flag has not been
268used on any of the built in models; models which were failing were
269streamlined so this flag was not necessary.
270
271**structure_factor = True** indicates that the model can be used as a
272structure factor to account for interactions between particles.  See
273`Form_Factors`_ for more details.
274
275Model Parameters
276................
277
278Next comes the parameter table.  For example::
279
280    # pylint: disable=bad-whitespace, line-too-long
281    #   ["name",        "units", default, [min, max], "type",    "description"],
282    parameters = [
283        ["sld",         "1e-6/Ang^2",  1, [-inf, inf], "sld",    "Layer scattering length density"],
284        ["sld_solvent", "1e-6/Ang^2",  6, [-inf, inf], "sld",    "Solvent scattering length density"],
285        ["radius",      "Ang",        50, [0, inf],    "volume", "Sphere radius"],
286    ]
287    # pylint: enable=bad-whitespace, line-too-long
288
289**parameters = [["name", "units", default, [min,max], "type", "tooltip"],...]**
290defines the parameters that form the model.
291
292**Note: The order of the parameters in the definition will be the order of the
293parameters in the user interface and the order of the parameters in Iq(),
294Iqxy() and form_volume(). And** *scale* **and** *background* **parameters are
295implicit to all models, so they do not need to be included in the parameter table.**
296
297- **"name"** is the name of the parameter shown on the FitPage.
298
[3048ec6]299  - the name must be a valid variable name, starting with a letter and
300    containing only letters, numbers and underscore.
301
[990d8df]302  - parameter names should follow the mathematical convention; e.g.,
303    *radius_core* not *core_radius*, or *sld_solvent* not *solvent_sld*.
304
305  - model parameter names should be consistent between different models,
306    so *sld_solvent*, for example, should have exactly the same name
307    in every model.
308
309  - to see all the parameter names currently in use, type the following in the
310    python shell/editor under the Tools menu::
311
312       import sasmodels.list_pars
313       sasmodels.list_pars.list_pars()
314
315    *re-use* as many as possible!!!
316
317  - use "name[n]" for multiplicity parameters, where *n* is the name of
318    the parameter defining the number of shells/layers/segments, etc.
319
320- **"units"** are displayed along with the parameter name
321
322  - every parameter should have units; use "None" if there are no units.
323
324  - **sld's should be given in units of 1e-6/Ang^2, and not simply
325    1/Ang^2 to be consistent with the builtin models.  Adjust your formulas
326    appropriately.**
327
328  - fancy units markup is available for some units, including::
329
330        Ang, 1/Ang, 1/Ang^2, 1e-6/Ang^2, degrees, 1/cm, Ang/cm, g/cm^3, mg/m^2
331
332  - the list of units is defined in the variable *RST_UNITS* within
333    `sasmodels/generate.py <https://github.com/SasView/sasmodels/tree/master/sasmodels/generate.py>`_
334
335    - new units can be added using the macros defined in *doc/rst_prolog*
336      in the sasmodels source.
337    - units should be properly formatted using sub-/super-scripts
338      and using negative exponents instead of the / operator, though
339      the unit name should use the / operator for consistency.
340    - please post a message to the SasView developers mailing list with your changes.
341
342- **default** is the initial value for the parameter.
343
344  - **the parameter default values are used to auto-generate a plot of
345    the model function in the documentation.**
346
347- **[min, max]** are the lower and upper limits on the parameter.
348
349  - lower and upper limits can be any number, or *-inf* or *inf*.
350
351  - the limits will show up as the default limits for the fit making it easy,
352    for example, to force the radius to always be greater than zero.
353
354  - these are hard limits defining the valid range of parameter values;
355    polydisperity distributions will be truncated at the limits.
356
357- **"type"** can be one of: "", "sld", "volume", or "orientation".
358
359  - "sld" parameters can have magnetic moments when fitting magnetic models;
360    depending on the spin polarization of the beam and the $q$ value being
361    examined, the effective sld for that material will be used to compute the
362    scattered intensity.
363
364  - "volume" parameters are passed to Iq(), Iqxy(), and form_volume(), and
365    have polydispersity loops generated automatically.
366
367  - "orientation" parameters are only passed to Iqxy(), and have angular
368    dispersion.
369
[9844c3a]370Some models will have integer parameters, such as number of pearls in the
371pearl necklace model, or number of shells in the multi-layer vesicle model.
372The optimizers in BUMPS treat all parameters as floating point numbers which
373can take arbitrary values, even for integer parameters, so your model should
374round the incoming parameter value to the nearest integer inside your model
375you should round to the nearest integer.  In C code, you can do this using::
376
377    static double
378    Iq(double q, ..., double fp_n, ...)
379    {
380        int n = (int)(fp_n + 0.5);
381        ...
382    }
383
384in python::
385
386    def Iq(q, ..., n, ...):
387        n = int(n+0.5)
388        ...
389
[3048ec6]390Derivative based optimizers such as Levenberg-Marquardt will not work
[9844c3a]391for integer parameters since the partial derivative is always zero, but
392the remaining optimizers (DREAM, differential evolution, Nelder-Mead simplex)
393will still function.
394
[990d8df]395Model Computation
396.................
397
398Models can be defined as pure python models, or they can be a mixture of
399python and C models.  C models are run on the GPU if it is available,
400otherwise they are compiled and run on the CPU.
401
402Models are defined by the scattering kernel, which takes a set of parameter
403values defining the shape, orientation and material, and returns the
404expected scattering. Polydispersity and angular dispersion are defined
405by the computational infrastructure.  Any parameters defined as "volume"
406parameters are polydisperse, with polydispersity defined in proportion
407to their value.  "orientation" parameters use angular dispersion defined
408in degrees, and are not relative to the current angle.
409
410Based on a weighting function $G(x)$ and a number of points $n$, the
411computed value is
412
413.. math::
414
415     \hat I(q)
416     = \frac{\int G(x) I(q, x)\,dx}{\int G(x) V(x)\,dx}
417     \approx \frac{\sum_{i=1}^n G(x_i) I(q,x_i)}{\sum_{i=1}^n G(x_i) V(x_i)}
418
[3048ec6]419That is, the individual models do not need to include polydispersity
[990d8df]420calculations, but instead rely on numerical integration to compute the
421appropriately smeared pattern.   Angular dispersion values over polar angle
422$\theta$ requires an additional $\cos \theta$ weighting due to decreased
423arc length for the equatorial angle $\phi$ with increasing latitude.
424
425Python Models
426.............
427
428For pure python models, define the *Iq* function::
429
430      import numpy as np
431      from numpy import cos, sin, ...
432
433      def Iq(q, par1, par2, ...):
434          return I(q, par1, par2, ...)
435      Iq.vectorized = True
436
437The parameters *par1, par2, ...* are the list of non-orientation parameters
438to the model in the order that they appear in the parameter table.
[3048ec6]439**Note that the auto-generated model file uses** *x* **rather than** *q*.
[990d8df]440
441The *.py* file should import trigonometric and exponential functions from
442numpy rather than from math.  This lets us evaluate the model for the whole
443range of $q$ values at once rather than looping over each $q$ separately in
444python.  With $q$ as a vector, you cannot use if statements, but must instead
445do tricks like
446
447::
448
449     a = x*q*(q>0) + y*q*(q<=0)
450
451or
452
453::
454
455     a = np.empty_like(q)
456     index = q>0
457     a[index] = x*q[index]
458     a[~index] = y*q[~index]
459
460which sets $a$ to $q \cdot x$ if $q$ is positive or $q \cdot y$ if $q$
461is zero or negative. If you have not converted your function to use $q$
462vectors, you can set the following and it will only receive one $q$
463value at a time::
464
465    Iq.vectorized = False
466
467Return np.NaN if the parameters are not valid (e.g., cap_radius < radius in
468barbell).  If I(q; pars) is NaN for any $q$, then those parameters will be
469ignored, and not included in the calculation of the weighted polydispersity.
470
471Similar to *Iq*, you can define *Iqxy(qx, qy, par1, par2, ...)* where the
472parameter list includes any orientation parameters.  If *Iqxy* is not defined,
473then it will default to *Iqxy = Iq(sqrt(qx**2+qy**2), par1, par2, ...)*.
474
475Models should define *form_volume(par1, par2, ...)* where the parameter
476list includes the *volume* parameters in order.  This is used for a weighted
477volume normalization so that scattering is on an absolute scale.  If
478*form_volume* is not defined, then the default *form_volume = 1.0* will be
479used.
480
481Embedded C Models
482.................
483
484Like pure python models, inline C models need to define an *Iq* function::
485
486    Iq = """
487        return I(q, par1, par2, ...);
488    """
489
490This expands into the equivalent C code::
491
492    #include <math.h>
493    double Iq(double q, double par1, double par2, ...);
494    double Iq(double q, double par1, double par2, ...)
495    {
496        return I(q, par1, par2, ...);
497    }
498
499*Iqxy* is similar to *Iq*, except it uses parameters *qx, qy* instead of *q*,
500and it includes orientation parameters.
501
502*form_volume* defines the volume of the shape. As in python models, it
503includes only the volume parameters.
504
505*Iqxy* will default to *Iq(sqrt(qx**2 + qy**2), par1, ...)* and
506*form_volume* will default to 1.0.
507
508**source=['fn.c', ...]** includes the listed C source files in the
509program before *Iq* and *Iqxy* are defined. This allows you to extend the
510library of C functions available to your model.
511
512Models are defined using double precision declarations for the
513parameters and return values.  When a model is run using single
514precision or long double precision, each variable is converted
515to the target type, depending on the precision requested.
516
517**Floating point constants must include the decimal point.**  This allows us
518to convert values such as 1.0 (double precision) to 1.0f (single precision)
519so that expressions that use these values are not promoted to double precision
520expressions.  Some graphics card drivers are confused when functions
521that expect floating point values are passed integers, such as 4*atan(1); it
522is safest to not use integers in floating point expressions.  Even better,
523use the builtin constant M_PI rather than 4*atan(1); it is faster and smaller!
524
525The C model operates on a single $q$ value at a time.  The code will be
526run in parallel across different $q$ values, either on the graphics card
527or the processor.
528
529Rather than returning NAN from Iq, you must define the *INVALID(v)*.  The
530*v* parameter lets you access all the parameters in the model using
531*v.par1*, *v.par2*, etc. For example::
532
533    #define INVALID(v) (v.bell_radius < v.radius)
534
535Special Functions
536.................
537
538The C code follows the C99 standard, with the usual math functions,
539as defined in
540`OpenCL <https://www.khronos.org/registry/cl/sdk/1.1/docs/man/xhtml/mathFunctions.html>`_.
541This includes the following:
542
543    M_PI, M_PI_2, M_PI_4, M_SQRT1_2, M_E:
544        $\pi$, $\pi/2$, $\pi/4$, $1/\sqrt{2}$ and Euler's constant $e$
[0a9fcab]545    exp, log, pow(x,y), expm1, log1p, sqrt, cbrt:
546        Power functions $e^x$, $\ln x$, $x^y$, $e^x - 1$, $\ln 1 + x$,
547        $\sqrt{x}$, $\sqrt[3]{x}$. The functions expm1(x) and log1p(x)
548        are accurate across all $x$, including $x$ very close to zero.
[990d8df]549    sin, cos, tan, asin, acos, atan:
550        Trigonometry functions and inverses, operating on radians.
551    sinh, cosh, tanh, asinh, acosh, atanh:
552        Hyperbolic trigonometry functions.
553    atan2(y,x):
554        Angle from the $x$\ -axis to the point $(x,y)$, which is equal to
555        $\tan^{-1}(y/x)$ corrected for quadrant.  That is, if $x$ and $y$ are
556        both negative, then atan2(y,x) returns a value in quadrant III where
557        atan(y/x) would return a value in quadrant I. Similarly for
558        quadrants II and IV when $x$ and $y$ have opposite sign.
[0a9fcab]559    fabs(x), fmin(x,y), fmax(x,y), trunc, rint:
[990d8df]560        Floating point functions.  rint(x) returns the nearest integer.
561    NAN:
562        NaN, Not a Number, $0/0$.  Use isnan(x) to test for NaN.  Note that
563        you cannot use :code:`x == NAN` to test for NaN values since that
[0a9fcab]564        will always return false.  NAN does not equal NAN!  The alternative,
565        :code:`x != x` may fail if the compiler optimizes the test away.
[990d8df]566    INFINITY:
567        $\infty, 1/0$.  Use isinf(x) to test for infinity, or isfinite(x)
568        to test for finite and not NaN.
569    erf, erfc, tgamma, lgamma:  **do not use**
570        Special functions that should be part of the standard, but are missing
571        or inaccurate on some platforms. Use sas_erf, sas_erfc and sas_gamma
572        instead (see below). Note: lgamma(x) has not yet been tested.
573
574Some non-standard constants and functions are also provided:
575
576    M_PI_180, M_4PI_3:
577        $\frac{\pi}{180}$, $\frac{4\pi}{3}$
578    SINCOS(x, s, c):
579        Macro which sets s=sin(x) and c=cos(x). The variables *c* and *s*
580        must be declared first.
581    square(x):
582        $x^2$
583    cube(x):
584        $x^3$
585    sas_sinx_x(x):
586        $\sin(x)/x$, with limit $\sin(0)/0 = 1$.
587    powr(x, y):
588        $x^y$ for $x \ge 0$; this is faster than general $x^y$ on some GPUs.
589    pown(x, n):
590        $x^n$ for $n$ integer; this is faster than general $x^n$ on some GPUs.
591    FLOAT_SIZE:
592        The number of bytes in a floating point value.  Even though all
593        variables are declared double, they may be converted to single
594        precision float before running. If your algorithm depends on
595        precision (which is not uncommon for numerical algorithms), use
596        the following::
597
598            #if FLOAT_SIZE>4
599            ... code for double precision ...
600            #else
601            ... code for single precision ...
602            #endif
603    SAS_DOUBLE:
604        A replacement for :code:`double` so that the declared variable will
605        stay double precision; this should generally not be used since some
606        graphics cards do not support double precision.  There is no provision
607        for forcing a constant to stay double precision.
608
609The following special functions and scattering calculations are defined in
610`sasmodels/models/lib <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib>`_.
611These functions have been tuned to be fast and numerically stable down
612to $q=0$ even in single precision.  In some cases they work around bugs
613which appear on some platforms but not others, so use them where needed.
614Add the files listed in :code:`source = ["lib/file.c", ...]` to your *model.py*
615file in the order given, otherwise these functions will not be available.
616
617    polevl(x, c, n):
618        Polynomial evaluation $p(x) = \sum_{i=0}^n c_i x^i$ using Horner's
619        method so it is faster and more accurate.
620
621        $c = \{c_n, c_{n-1}, \ldots, c_0 \}$ is the table of coefficients,
622        sorted from highest to lowest.
623
624        :code:`source = ["lib/polevl.c", ...]` (`link to code <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/polevl.c>`_)
625
626    p1evl(x, c, n):
627        Evaluation of normalized polynomial $p(x) = x^n + \sum_{i=0}^{n-1} c_i x^i$
628        using Horner's method so it is faster and more accurate.
629
630        $c = \{c_{n-1}, c_{n-2} \ldots, c_0 \}$ is the table of coefficients,
631        sorted from highest to lowest.
632
633        :code:`source = ["lib/polevl.c", ...]`
[870a2f4]634        (`polevl.c <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/polevl.c>`_)
[990d8df]635
636    sas_gamma(x):
[30b60d2]637        Gamma function sas_gamma\ $(x) = \Gamma(x)$.
[990d8df]638
639        The standard math function, tgamma(x) is unstable for $x < 1$
640        on some platforms.
641
[870a2f4]642        :code:`source = ["lib/sas_gamma.c", ...]`
643        (`sas_gamma.c <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/sas_gamma.c>`_)
[990d8df]644
645    sas_erf(x), sas_erfc(x):
646        Error function
[30b60d2]647        sas_erf\ $(x) = \frac{2}{\sqrt\pi}\int_0^x e^{-t^2}\,dt$
[990d8df]648        and complementary error function
[30b60d2]649        sas_erfc\ $(x) = \frac{2}{\sqrt\pi}\int_x^{\infty} e^{-t^2}\,dt$.
[990d8df]650
651        The standard math functions erf(x) and erfc(x) are slower and broken
652        on some platforms.
653
654        :code:`source = ["lib/polevl.c", "lib/sas_erf.c", ...]`
[870a2f4]655        (`sas_erf.c <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/sas_erf.c>`_)
[990d8df]656
657    sas_J0(x):
[30b60d2]658        Bessel function of the first kind sas_J0\ $(x)=J_0(x)$ where
[990d8df]659        $J_0(x) = \frac{1}{\pi}\int_0^\pi \cos(x\sin(\tau))\,d\tau$.
660
661        The standard math function j0(x) is not available on all platforms.
662
663        :code:`source = ["lib/polevl.c", "lib/sas_J0.c", ...]`
[870a2f4]664        (`sas_J0.c <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/sas_J0.c>`_)
[990d8df]665
666    sas_J1(x):
[30b60d2]667        Bessel function of the first kind  sas_J1\ $(x)=J_1(x)$ where
[990d8df]668        $J_1(x) = \frac{1}{\pi}\int_0^\pi \cos(\tau - x\sin(\tau))\,d\tau$.
669
670        The standard math function j1(x) is not available on all platforms.
671
672        :code:`source = ["lib/polevl.c", "lib/sas_J1.c", ...]`
[870a2f4]673        (`sas_J1.c <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/sas_J1.c>`_)
[990d8df]674
675    sas_JN(n, x):
[30b60d2]676        Bessel function of the first kind and integer order $n$,
677        sas_JN\ $(n, x) =J_n(x)$ where
[990d8df]678        $J_n(x) = \frac{1}{\pi}\int_0^\pi \cos(n\tau - x\sin(\tau))\,d\tau$.
[30b60d2]679        If $n$ = 0 or 1, it uses sas_J0($x$) or sas_J1($x$), respectively.
[990d8df]680
681        The standard math function jn(n, x) is not available on all platforms.
682
683        :code:`source = ["lib/polevl.c", "lib/sas_J0.c", "lib/sas_J1.c", "lib/sas_JN.c", ...]`
[870a2f4]684        (`sas_JN.c <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/sas_JN.c>`_)
[990d8df]685
686    sas_Si(x):
[30b60d2]687        Sine integral Si\ $(x) = \int_0^x \tfrac{\sin t}{t}\,dt$.
[990d8df]688
689        This function uses Taylor series for small and large arguments:
690
691        For large arguments,
692
693        .. math::
694
695             \text{Si}(x) \sim \frac{\pi}{2}
696             - \frac{\cos(x)}{x}\left(1 - \frac{2!}{x^2} + \frac{4!}{x^4} - \frac{6!}{x^6} \right)
697             - \frac{\sin(x)}{x}\left(\frac{1}{x} - \frac{3!}{x^3} + \frac{5!}{x^5} - \frac{7!}{x^7}\right)
698
699        For small arguments,
700
701        .. math::
702
703           \text{Si}(x) \sim x
704           - \frac{x^3}{3\times 3!} + \frac{x^5}{5 \times 5!} - \frac{x^7}{7 \times 7!}
705           + \frac{x^9}{9\times 9!} - \frac{x^{11}}{11\times 11!}
706
707        :code:`source = ["lib/Si.c", ...]`
[870a2f4]708        (`Si.c <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/Si.c>`_)
[990d8df]709
710    sas_3j1x_x(x):
711        Spherical Bessel form
[30b60d2]712        sph_j1c\ $(x) = 3 j_1(x)/x = 3 (\sin(x) - x \cos(x))/x^3$,
[990d8df]713        with a limiting value of 1 at $x=0$, where $j_1(x)$ is the spherical
714        Bessel function of the first kind and first order.
715
716        This function uses a Taylor series for small $x$ for numerical accuracy.
717
718        :code:`source = ["lib/sas_3j1x_x.c", ...]`
[870a2f4]719        (`sas_3j1x_x.c <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/sas_3j1x_x.c>`_)
[990d8df]720
721
722    sas_2J1x_x(x):
[30b60d2]723        Bessel form sas_J1c\ $(x) = 2 J_1(x)/x$, with a limiting value
[990d8df]724        of 1 at $x=0$, where $J_1(x)$ is the Bessel function of first kind
725        and first order.
726
727        :code:`source = ["lib/polevl.c", "lib/sas_J1.c", ...]`
[870a2f4]728        (`sas_J1.c <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/sas_J1.c>`_)
[990d8df]729
730
731    Gauss76Z[i], Gauss76Wt[i]:
732        Points $z_i$ and weights $w_i$ for 76-point Gaussian quadrature, respectively,
733        computing $\int_{-1}^1 f(z)\,dz \approx \sum_{i=1}^{76} w_i\,f(z_i)$.
734
735        Similar arrays are available in :code:`gauss20.c` for 20-point
736        quadrature and in :code:`gauss150.c` for 150-point quadrature.
[0a9fcab]737        The macros :code:`GAUSS_N`, :code:`GAUSS_Z` and :code:`GAUSS_W` are
738        defined so that you can change the order of the integration by
739        selecting an different source without touching the C code.
[990d8df]740
741        :code:`source = ["lib/gauss76.c", ...]`
[870a2f4]742        (`gauss76.c <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/gauss76.c>`_)
[990d8df]743
744
745
746Problems with C models
747......................
748
749The graphics processor (GPU) in your computer is a specialized computer tuned
750for certain kinds of problems.  This leads to strange restrictions that you
751need to be aware of.  Your code may work fine on some platforms or for some
752models, but then return bad values on other platforms.  Some examples of
753particular problems:
754
755  **(1) Code is too complex, or uses too much memory.** GPU devices only
756  have a limited amount of memory available for each processor. If you run
757  programs which take too much memory, then rather than running multiple
758  values in parallel as it usually does, the GPU may only run a single
759  version of the code at a time, making it slower than running on the CPU.
760  It may fail to run on some platforms, or worse, cause the screen to go
761  blank or the system to reboot.
762
763  **(2) Code takes too long.** Because GPU devices are used for the computer
764  display, the OpenCL drivers are very careful about the amount of time they
765  will allow any code to run. For example, on OS X, the model will stop
766  running after 5 seconds regardless of whether the computation is complete.
767  You may end up with only some of your 2D array defined, with the rest
768  containing random data. Or it may cause the screen to go blank or the
769  system to reboot.
770
771  **(3) Memory is not aligned**. The GPU hardware is specialized to operate
772  on multiple values simultaneously. To keep the GPU simple the values in
773  memory must be aligned with the different GPU compute engines. Not
774  following these rules can lead to unexpected values being loaded into
775  memory, and wrong answers computed. The conclusion from a very long and
776  strange debugging session was that any arrays that you declare in your
777  model should be a multiple of four. For example::
778
779      double Iq(q, p1, p2, ...)
780      {
781          double vector[8];  // Only going to use seven slots, but declare 8
782          ...
783      }
784
785The first step when your model is behaving strangely is to set
786**single=False**. This automatically restricts the model to only run on the
787CPU, or on high-end GPU cards. There can still be problems even on high-end
788cards, so you can force the model off the GPU by setting **opencl=False**.
789This runs the model as a normal C program without any GPU restrictions so
790you know that strange results are probably from your code rather than the
791environment. Once the code is debugged, you can compare your output to the
792output on the GPU.
793
794Although it can be difficult to get your model to work on the GPU, the reward
795can be a model that runs 1000x faster on a good card.  Even your laptop may
796show a 50x improvement or more over the equivalent pure python model.
797
798External C Models
799.................
800
801External C models are very much like embedded C models, except that
802*Iq*, *Iqxy* and *form_volume* are defined in an external source file
803loaded using the *source=[...]* statement. You need to supply the function
804declarations for each of these that you need instead of building them
805automatically from the parameter table.
806
807
808.. _Form_Factors:
809
810Form Factors
811............
812
813Away from the dilute limit you can estimate scattering including
814particle-particle interactions using $I(q) = P(q)*S(q)$ where $P(q)$
815is the form factor and $S(q)$ is the structure factor.  The simplest
816structure factor is the *hardsphere* interaction, which
817uses the effective radius of the form factor as an input to the structure
818factor model.  The effective radius is the average radius of the
819form averaged over all the polydispersity values.
820
821::
822
823    def ER(radius, thickness):
824        """Effective radius of a core-shell sphere."""
825        return radius + thickness
826
827Now consider the *core_shell_sphere*, which has a simple effective radius
828equal to the radius of the core plus the thickness of the shell, as
829shown above. Given polydispersity over *(r1, r2, ..., rm)* in radius and
830*(t1, t2, ..., tn)* in thickness, *ER* is called with a mesh
831grid covering all possible combinations of radius and thickness.
832That is, *radius* is *(r1, r2, ..., rm, r1, r2, ..., rm, ...)*
833and *thickness* is *(t1, t1, ... t1, t2, t2, ..., t2, ...)*.
834The *ER* function returns one effective radius for each combination.
835The effective radius calculator weights each of these according to
836the polydispersity distributions and calls the structure factor
837with the average *ER*.
838
839::
840
841    def VR(radius, thickness):
842        """Sphere and shell volumes for a core-shell sphere."""
843        whole = 4.0/3.0 * pi * (radius + thickness)**3
844        core = 4.0/3.0 * pi * radius**3
845        return whole, whole - core
846
847Core-shell type models have an additional volume ratio which scales
848the structure factor.  The *VR* function returns the volume of
849the whole sphere and the volume of the shell. Like *ER*, there is
850one return value for each point in the mesh grid.
851
852*NOTE: we may be removing or modifying this feature soon. As of the
853time of writing, core-shell sphere returns (1., 1.) for VR, giving a volume
854ratio of 1.0.*
855
856Unit Tests
857..........
858
859THESE ARE VERY IMPORTANT. Include at least one test for each model and
860PLEASE make sure that the answer value is correct (i.e. not a random number).
861
862::
863
864    tests = [
865        [{}, 0.2, 0.726362],
866        [{"scale": 1., "background": 0., "sld": 6., "sld_solvent": 1.,
867          "radius": 120., "radius_pd": 0.2, "radius_pd_n":45},
868         0.2, 0.228843],
869        [{"radius": 120., "radius_pd": 0.2, "radius_pd_n":45}, "ER", 120.],
870        [{"radius": 120., "radius_pd": 0.2, "radius_pd_n":45}, "VR", 1.],
871    ]
872
873
874**tests=[[{parameters}, q, result], ...]** is a list of lists.
875Each list is one test and contains, in order:
876
877- a dictionary of parameter values. This can be *{}* using the default
878  parameters, or filled with some parameters that will be different from the
879  default, such as *{"radius":10.0, "sld":4}*. Unlisted parameters will
880  be given the default values.
881- the input $q$ value or tuple of $(q_x, q_y)$ values.
882- the output $I(q)$ or $I(q_x,q_y)$ expected of the model for the parameters
883  and input value given.
884- input and output values can themselves be lists if you have several
885  $q$ values to test for the same model parameters.
886- for testing *ER* and *VR*, give the inputs as "ER" and "VR" respectively;
887  the output for *VR* should be the sphere/shell ratio, not the individual
888  sphere and shell values.
889
890.. _Test_Your_New_Model:
891
892Test Your New Model
893^^^^^^^^^^^^^^^^^^^
894
895Minimal Testing
896...............
897
898From SasView either open the Python shell (*Tools* > *Python Shell/Editor*)
899or the plugin editor (*Fitting* > *Plugin Model Operations* > *Advanced
900Plugin Editor*), load your model, and then select *Run > Check Model* from
901the menu bar. An *Info* box will appear with the results of the compilation
902and a check that the model runs.
903
904If you are not using sasmodels from SasView, skip this step.
905
906Recommended Testing
907...................
908
909If the model compiles and runs, you can next run the unit tests that
910you have added using the **test =** values.
911
912From SasView, switch to the *Shell* tab and type the following::
913
914    from sasmodels.model_test import run_one
915    run_one("~/.sasview/plugin_models/model.py")
916
917This should print::
918
919    test_model_python (sasmodels.model_test.ModelTestCase) ... ok
920
921To check whether single precision is good enough, type the following::
922
923    from sasmodels.compare import main as compare
924    compare("~/.sasview/plugin_models/model.py")
925
926This will pop up a plot showing the difference between single precision
927and double precision on a range of $q$ values.
928
929::
930
931  demo = dict(scale=1, background=0,
932              sld=6, sld_solvent=1,
933              radius=120,
934              radius_pd=.2, radius_pd_n=45)
935
936**demo={'par': value, ...}** in the model file sets the default values for
937the comparison. You can include polydispersity parameters such as
938*radius_pd=0.2, radius_pd_n=45* which would otherwise be zero.
939
940These commands can also be run directly in the python interpreter:
941
942    $ python -m sasmodels.model_test -v ~/.sasview/plugin_models/model.py
943    $ python -m sasmodels.compare ~/.sasview/plugin_models/model.py
944
945The options to compare are quite extensive; type the following for help::
946
947    compare()
948
949Options will need to be passed as separate strings.
950For example to run your model with a random set of parameters::
951
952    compare("-random", "-pars", "~/.sasview/plugin_models/model.py")
953
954For the random models,
955
956- *sld* will be in the range (-0.5,10.5),
957- angles (*theta, phi, psi*) will be in the range (-180,180),
958- angular dispersion will be in the range (0,45),
959- polydispersity will be in the range (0,1)
960- other values will be in the range (0, 2\ *v*), where *v* is the value
961  of the parameter in demo.
962
963Dispersion parameters *n*\, *sigma* and *type* will be unchanged from
964demo so that run times are more predictable (polydispersity calculated
965across multiple parameters can be very slow).
966
[3048ec6]967If your model has 2D orientation calculation, then you should also
[990d8df]968test with::
969
970    compare("-2d", "~/.sasview/plugin_models/model.py")
971
972Check The Docs
973^^^^^^^^^^^^^^
974
975You can get a rough idea of how the documentation will look using the
976following::
977
978    compare("-help", "~/.sasview/plugin_models/model.py")
979
980This does not use the same styling as the rest of the docs, but it will
981allow you to check that your ReStructuredText and LaTeX formatting.
982Here are some tools to help with the inevitable syntax errors:
983
984- `Sphinx cheat sheet <http://matplotlib.org/sampledoc/cheatsheet.html>`_
985- `Sphinx Documentation <http://www.sphinx-doc.org/en/stable/>`_
986- `MathJax <http://www.mathjax.org/>`_
987- `amsmath <http://www.ams.org/publications/authors/tex/amslatex>`_
988
989There is also a neat online WYSIWYG ReStructuredText editor at
990http://rst.ninjs.org\ .
991
992
993Clean Lint - (Developer Version Only)
994^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
995
996**NB: For now we are not providing pylint with the installer version
997of SasView; so unless you have a SasView build environment available,
998you can ignore this section!**
999
1000Run the lint check with::
1001
1002    python -m pylint --rcfile=extra/pylint.rc ~/.sasview/plugin_models/model.py
1003
1004We are not aiming for zero lint just yet, only keeping it to a minimum.
1005For now, don't worry too much about *invalid-name*. If you really want a
1006variable name *Rg* for example because $R_g$ is the right name for the model
1007parameter then ignore the lint errors.  Also, ignore *missing-docstring*
1008for standard model functions *Iq*, *Iqxy*, etc.
1009
1010We will have delinting sessions at the SasView Code Camps, where we can
1011decide on standards for model files, parameter names, etc.
1012
1013For now, you can tell pylint to ignore things.  For example, to align your
1014parameters in blocks::
1015
1016    # pylint: disable=bad-whitespace,line-too-long
1017    #   ["name",                  "units", default, [lower, upper], "type", "description"],
1018    parameters = [
1019        ["contrast_factor",       "barns",    10.0,  [-inf, inf], "", "Contrast factor of the polymer"],
1020        ["bjerrum_length",        "Ang",       7.1,  [0, inf],    "", "Bjerrum length"],
1021        ["virial_param",          "1/Ang^2",  12.0,  [-inf, inf], "", "Virial parameter"],
1022        ["monomer_length",        "Ang",      10.0,  [0, inf],    "", "Monomer length"],
1023        ["salt_concentration",    "mol/L",     0.0,  [-inf, inf], "", "Concentration of monovalent salt"],
1024        ["ionization_degree",     "",          0.05, [0, inf],    "", "Degree of ionization"],
1025        ["polymer_concentration", "mol/L",     0.7,  [0, inf],    "", "Polymer molar concentration"],
1026        ]
1027    # pylint: enable=bad-whitespace,line-too-long
1028
1029Don't put in too many pylint statements, though, since they make the code ugly.
1030
1031Share Your Model!
1032^^^^^^^^^^^^^^^^^
1033
1034Once compare and the unit test(s) pass properly and everything is done,
1035consider adding your model to the
1036`Model Marketplace <http://marketplace.sasview.org/>`_ so that others may use it!
1037
1038.. ZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZ
1039
1040*Document History*
1041
1042| 2016-10-25 Steve King
1043| 2017-05-07 Paul Kienzle - Moved from sasview to sasmodels docs
Note: See TracBrowser for help on using the repository browser.