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

ESS_GUIESS_GUI_DocsESS_GUI_batch_fittingESS_GUI_bumps_abstractionESS_GUI_iss1116ESS_GUI_iss879ESS_GUI_iss959ESS_GUI_openclESS_GUI_orderingESS_GUI_sync_sascalc
Last change on this file since 0215e0a was 9687d58, checked in by Piotr Rozyczko <rozyczko@…>, 8 years ago

Merge branch 'master' into ESS_GUI

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