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

core_shell_microgelscostrafo411magnetic_modelticket-1257-vesicle-productticket_1156ticket_1265_superballticket_822_more_unit_tests
Last change on this file since ba2ed54 was 870a2f4, checked in by Paul Kienzle <pkienzle@…>, 7 years ago

fix doc builds

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