source: sasview/src/sas/sasgui/perspectives/fitting/media/plugin.rst @ f485ba0

ESS_GUIESS_GUI_DocsESS_GUI_batch_fittingESS_GUI_bumps_abstractionESS_GUI_iss1116ESS_GUI_iss879ESS_GUI_iss959ESS_GUI_openclESS_GUI_orderingESS_GUI_sync_sascalccostrafo411magnetic_scattrelease-4.2.2ticket-1009ticket-1094-headlessticket-1242-2d-resolutionticket-1243ticket-1249ticket885unittest-saveload
Last change on this file since f485ba0 was f485ba0, checked in by smk78, 7 years ago

Updated fitting_help.rst amd plugin.rst in light of work on #767 & #861

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