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

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.1.1release-4.1.2release-4.2.2ticket-1009ticket-1094-headlessticket-1242-2d-resolutionticket-1243ticket-1249ticket885unittest-saveload
Last change on this file since 2510b9b was ca6cbc1c, checked in by wojciech, 8 years ago

Updated documentation for lib functions in sasmodels. Refers to ticket 804

  • Property mode set to 100644
File size: 40.3 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
[7f23423]21Overview
22^^^^^^^^
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
[907186d]29it to the list of *Customized Models* in a FitPage.
[05829fb]30
[3e1c9e5]31SasView models can be of three types:
[05829fb]32
[3e1c9e5]33- A pure python model : Example -
[05829fb]34  `broadpeak.py <https://github.com/SasView/sasmodels/blob/master/sasmodels/models/broad_peak.py>`_
[3e1c9e5]35- A python model with embedded C : Example -
[05829fb]36  `sphere.py <https://github.com/SasView/sasmodels/blob/master/sasmodels/models/sphere.py>`_
[3e1c9e5]37- A python wrapper with separate C code : Example -
[05829fb]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
[05829fb]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
[05829fb]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)
[05829fb]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
[05829fb]62  broad_peak not BroadPeak)
[907186d]63- Do not include "model" in the name (i.e. barbell not BarBellModel)
[05829fb]64
65
66Edit New Model Files
67^^^^^^^^^^^^^^^^^^^^
68
[7f23423]69Model Contents
70..............
71
[05829fb]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"
[05829fb]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
[05829fb]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
[05829fb]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!
[05829fb]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
[05829fb]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
[05829fb]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
[05829fb]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.
[05829fb]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
[05829fb]230    title = "Spheres with uniform scattering length density"
[7f23423]231
[05829fb]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
[05829fb]241    category = "shape:sphere"
[7f23423]242
[05829fb]243    single = True   # optional: defaults to True
[7f23423]244
[05829fb]245    opencl = False  # optional: defaults to False
[7f23423]246
[05829fb]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.
[05829fb]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.
[05829fb]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!**
[05829fb]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.
[05829fb]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
[05829fb]304
305**parameters = [["name", "units", default, [min,max], "type", "tooltip"],...]**
[7f23423]306defines the parameters that form the model.
[05829fb]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.**
[05829fb]312
[7f23423]313- **"name"** is the name of the parameter shown on the FitPage.
[05829fb]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
[05829fb]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
[05829fb]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
[05829fb]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
[05829fb]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
[05829fb]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.
[05829fb]354
[7f23423]355- **default** is the initial value for the parameter.
[05829fb]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*.
[05829fb]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
[7f23423]367- **"type"** can be one of: "", "sld", "volume", or "orientation".
[05829fb]368
369  - "sld" parameters can have magnetic moments when fitting magnetic models;
370    depending on the spin polarization of the beam and the $q$ value being
371    examined, the effective sld for that material will be used to compute the
[7f23423]372    scattered intensity.
373
[05829fb]374  - "volume" parameters are passed to Iq(), Iqxy(), and form_volume(), and
375    have polydispersity loops generated automatically.
[7f23423]376
[05829fb]377  - "orientation" parameters are only passed to Iqxy(), and have angular
378    dispersion.
379
380
381Model Computation
382.................
383
384Models can be defined as pure python models, or they can be a mixture of
385python and C models.  C models are run on the GPU if it is available,
386otherwise they are compiled and run on the CPU.
387
388Models are defined by the scattering kernel, which takes a set of parameter
389values defining the shape, orientation and material, and returns the
390expected scattering. Polydispersity and angular dispersion are defined
391by the computational infrastructure.  Any parameters defined as "volume"
392parameters are polydisperse, with polydispersity defined in proportion
393to their value.  "orientation" parameters use angular dispersion defined
394in degrees, and are not relative to the current angle.
395
396Based on a weighting function $G(x)$ and a number of points $n$, the
397computed value is
398
399.. math::
400
401     \hat I(q)
402     = \frac{\int G(x) I(q, x)\,dx}{\int G(x) V(x)\,dx}
403     \approx \frac{\sum_{i=1}^n G(x_i) I(q,x_i)}{\sum_{i=1}^n G(x_i) V(x_i)}
404
405That is, the indivdual models do not need to include polydispersity
406calculations, but instead rely on numerical integration to compute the
407appropriately smeared pattern.   Angular dispersion values over polar angle
408$\theta$ requires an additional $\cos \theta$ weighting due to decreased
409arc length for the equatorial angle $\phi$ with increasing latitude.
410
411Python Models
412.............
413
[7f23423]414For pure python models, define the *Iq* function::
[05829fb]415
416      import numpy as np
417      from numpy import cos, sin, ...
418
419      def Iq(q, par1, par2, ...):
420          return I(q, par1, par2, ...)
421      Iq.vectorized = True
422
423The parameters *par1, par2, ...* are the list of non-orientation parameters
424to the model in the order that they appear in the parameter table.
[7f23423]425**Note that the autogenerated model file uses** *x* **rather than** *q*.
[05829fb]426
427The *.py* file should import trigonometric and exponential functions from
[7f23423]428numpy rather than from math.  This lets us evaluate the model for the whole
[05829fb]429range of $q$ values at once rather than looping over each $q$ separately in
430python.  With $q$ as a vector, you cannot use if statements, but must instead
431do tricks like
432
433::
434
435     a = x*q*(q>0) + y*q*(q<=0)
436
437or
438
439::
440
441     a = np.empty_like(q)
442     index = q>0
443     a[index] = x*q[index]
444     a[~index] = y*q[~index]
445
446which sets $a$ to $q \cdot x$ if $q$ is positive or $q \cdot y$ if $q$
447is zero or negative. If you have not converted your function to use $q$
448vectors, you can set the following and it will only receive one $q$
449value at a time::
450
451    Iq.vectorized = False
452
453Return np.NaN if the parameters are not valid (e.g., cap_radius < radius in
454barbell).  If I(q; pars) is NaN for any $q$, then those parameters will be
455ignored, and not included in the calculation of the weighted polydispersity.
456
457Similar to *Iq*, you can define *Iqxy(qx, qy, par1, par2, ...)* where the
458parameter list includes any orientation parameters.  If *Iqxy* is not defined,
459then it will default to *Iqxy = Iq(sqrt(qx**2+qy**2), par1, par2, ...)*.
460
461Models should define *form_volume(par1, par2, ...)* where the parameter
462list includes the *volume* parameters in order.  This is used for a weighted
463volume normalization so that scattering is on an absolute scale.  If
[7f23423]464*form_volume* is not defined, then the default *form_volume = 1.0* will be
[05829fb]465used.
466
467Embedded C Models
468.................
469
[7f23423]470Like pure python models, inline C models need to define an *Iq* function::
[05829fb]471
472    Iq = """
473        return I(q, par1, par2, ...);
474    """
475
476This expands into the equivalent C code::
477
478    #include <math.h>
479    double Iq(double q, double par1, double par2, ...);
480    double Iq(double q, double par1, double par2, ...)
481    {
482        return I(q, par1, par2, ...);
483    }
484
[af6de50]485*Iqxy* is similar to *Iq*, except it uses parameters *qx, qy* instead of *q*,
486and it includes orientation parameters.
487
[907186d]488*form_volume* defines the volume of the shape. As in python models, it
[af6de50]489includes only the volume parameters.
490
491*Iqxy* will default to *Iq(sqrt(qx**2 + qy**2), par1, ...)* and
492*form_volume* will default to 1.0.
493
494**source=['fn.c', ...]** includes the listed C source files in the
495program before *Iq* and *Iqxy* are defined. This allows you to extend the
496library of C functions available to your model.
497
498Models are defined using double precision declarations for the
499parameters and return values.  When a model is run using single
500precision or long double precision, each variable is converted
501to the target type, depending on the precision requested.
502
503**Floating point constants must include the decimal point.**  This allows us
504to convert values such as 1.0 (double precision) to 1.0f (single precision)
505so that expressions that use these values are not promoted to double precision
506expressions.  Some graphics card drivers are confused when functions
507that expect floating point values are passed integers, such as 4*atan(1); it
508is safest to not use integers in floating point expressions.  Even better,
509use the builtin constant M_PI rather than 4*atan(1); it is faster and smaller!
510
[05829fb]511The C model operates on a single $q$ value at a time.  The code will be
512run in parallel across different $q$ values, either on the graphics card
513or the processor.
514
515Rather than returning NAN from Iq, you must define the *INVALID(v)*.  The
516*v* parameter lets you access all the parameters in the model using
517*v.par1*, *v.par2*, etc. For example::
518
519    #define INVALID(v) (v.bell_radius < v.radius)
520
[af6de50]521Special Functions
522.................
[05829fb]523
[af6de50]524The C code follows the C99 standard, with the usual math functions,
[05829fb]525as defined in
526`OpenCL <https://www.khronos.org/registry/cl/sdk/1.1/docs/man/xhtml/mathFunctions.html>`_.
[af6de50]527This includes the following:
528
529    M_PI, M_PI_2, M_PI_4, M_SQRT1_2, M_E:
530        $\pi$, $\pi/2$, $\pi/4$, $1/\sqrt{2}$ and Euler's constant $e$
531    exp, log, pow(x,y), expm1, sqrt:
532        Power functions $e^x$, $\ln x$, $x^y$, $e^x - 1$, $\sqrt{x}$.
533        The function expm1(x) is accurate across all $x$, including $x$
534        very close to zero.
535    sin, cos, tan, asin, acos, atan:
536        Trigonometry functions and inverses, operating on radians.
537    sinh, cos, tanh, asinh, acosh, atanh:
538        Hyperbolic trigonometry functions.
539    atan2(y,x):
540        Angle from the $x$\ -axis to the point $(x,y)$, which is equal to
541        $\tan^{-1}(y/x)$ corrected for quadrant.  That is, if $x$ and $y$ are
542        both negative, then atan2(y,x) returns a value in quadrant III where
543        atan(y/x) would return a value in quadrant I. Similarly for
544        quadrants II and IV when $x$ and $y$ have opposite sign.
545    fmin(x,y), fmax(x,y), trunc, rint:
546        Floating point functions.  rint(x) returns the nearest integer.
547    NAN:
548        NaN, Not a Number, $0/0$.  Use isnan(x) to test for NaN.  Note that
549        you cannot use :code:`x == NAN` to test for NaN values since that
550        will always return false.  NAN does not equal NAN!
551    INFINITY:
552        $\infty, 1/0$.  Use isinf(x) to test for infinity, or isfinite(x)
553        to test for finite and not NaN.
554    erf, erfc, tgamma, lgamma:  **do not use**
555        Special functions that should be part of the standard, but are missing
556        or inaccurate on some platforms. Use sas_erf, sas_erfc and sas_gamma
557        instead (see below). Note: lgamma(x) has not yet been tested.
558
559Some non-standard constants and functions are also provided:
560
561    M_PI_180, M_4PI_3:
[ca1eaeb]562        $\frac{\pi}{180}$, $\frac{4\pi}{3}$
[af6de50]563    SINCOS(x, s, c):
564        Macro which sets s=sin(x) and c=cos(x). The variables *c* and *s*
565        must be declared first.
566    square(x):
567        $x^2$
568    cube(x):
569        $x^3$
[ca6cbc1c]570    sas_sinx_x(x):
[af6de50]571        $\sin(x)/x$, with limit $\sin(0)/0 = 1$.
572    powr(x, y):
573        $x^y$ for $x \ge 0$; this is faster than general $x^y$ on some GPUs.
574    pown(x, n):
575        $x^n$ for $n$ integer; this is faster than general $x^n$ on some GPUs.
576    FLOAT_SIZE:
577        The number of bytes in a floating point value.  Even though all
578        variables are declared double, they may be converted to single
579        precision float before running. If your algorithm depends on
580        precision (which is not uncommon for numerical algorithms), use
581        the following::
582
583            #if FLOAT_SIZE>4
584            ... code for double precision ...
585            #else
586            ... code for single precision ...
587            #endif
588    SAS_DOUBLE:
589        A replacement for :code:`double` so that the declared variable will
590        stay double precision; this should generally not be used since some
591        graphics cards do not support double precision.  There is no provision
592        for forcing a constant to stay double precision.
593
594The following special functions and scattering calculations are defined in
595`sasmodels/models/lib <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib>`_.
[05829fb]596These functions have been tuned to be fast and numerically stable down
597to $q=0$ even in single precision.  In some cases they work around bugs
[ca1eaeb]598which appear on some platforms but not others, so use them where needed.
599Add the files listed in :code:`source = ["lib/file.c", ...]` to your *model.py*
600file in the order given, otherwise these functions will not be available.
[05829fb]601
[af6de50]602    polevl(x, c, n):
[ca1eaeb]603        Polynomial evaluation $p(x) = \sum_{i=0}^n c_i x^i$ using Horner's
[af6de50]604        method so it is faster and more accurate.
[7f23423]605
[ca1eaeb]606        $c = \{c_n, c_{n-1}, \ldots, c_0 \}$ is the table of coefficients,
607        sorted from highest to lowest.
[48cd5b3]608
[ca1eaeb]609        :code:`source = ["lib/polevl.c", ...]` (`link to code <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/polevl.c>`_)
[48cd5b3]610
611    p1evl(x, c, n):
[ca1eaeb]612        Evaluation of normalized polynomial $p(x) = x^n + \sum_{i=0}^{n-1} c_i x^i$
613        using Horner's method so it is faster and more accurate.
[48cd5b3]614
[ca1eaeb]615        $c = \{c_{n-1}, c_{n-2} \ldots, c_0 \}$ is the table of coefficients,
616        sorted from highest to lowest.
[48cd5b3]617
[af6de50]618        :code:`source = ["lib/polevl.c", ...]`
[48cd5b3]619        (`link to code <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/polevl.c>`_)
[af6de50]620
[48cd5b3]621    sas_gamma(x):
622        Gamma function $\text{sas_gamma}(x) = \Gamma(x)$.
623
[ca1eaeb]624        The standard math function, tgamma(x) is unstable for $x < 1$
625        on some platforms.
[af6de50]626
627        :code:`source = ["lib/sasgamma.c", ...]`
[48cd5b3]628        (`link to code <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/sas_gamma.c>`_)
[af6de50]629
[48cd5b3]630    sas_erf(x), sas_erfc(x):
[af6de50]631        Error function
[ca1eaeb]632        $\text{sas_erf}(x) = \frac{2}{\sqrt\pi}\int_0^x e^{-t^2}\,dt$
[af6de50]633        and complementary error function
[ca1eaeb]634        $\text{sas_erfc}(x) = \frac{2}{\sqrt\pi}\int_x^{\infty} e^{-t^2}\,dt$.
[48cd5b3]635
[ca1eaeb]636        The standard math functions erf(x) and erfc(x) are slower and broken
[af6de50]637        on some platforms.
638
639        :code:`source = ["lib/polevl.c", "lib/sas_erf.c", ...]`
[48cd5b3]640        (`link to error functions' code <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/sas_erf.c>`_)
[af6de50]641
[48cd5b3]642    sas_J0(x):
643        Bessel function of the first kind $\text{sas_J0}(x)=J_0(x)$ where
[af6de50]644        $J_0(x) = \frac{1}{\pi}\int_0^\pi \cos(x\sin(\tau))\,d\tau$.
645
[ca1eaeb]646        The standard math function j0(x) is not available on all platforms.
647
[af6de50]648        :code:`source = ["lib/polevl.c", "lib/sas_J0.c", ...]`
[48cd5b3]649        (`link to Bessel function's code <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/sas_J0.c>`_)
[af6de50]650
[48cd5b3]651    sas_J1(x):
652        Bessel function of the first kind  $\text{sas_J1}(x)=J_1(x)$ where
[af6de50]653        $J_1(x) = \frac{1}{\pi}\int_0^\pi \cos(\tau - x\sin(\tau))\,d\tau$.
654
[ca1eaeb]655        The standard math function j1(x) is not available on all platforms.
656
[af6de50]657        :code:`source = ["lib/polevl.c", "lib/sas_J1.c", ...]`
[48cd5b3]658        (`link to Bessel function's code <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/sas_J1.c>`_)
[af6de50]659
[48cd5b3]660    sas_JN(n, x):
[ca1eaeb]661        Bessel function of the first kind and integer order $n$:
662        $\text{sas_JN}(n, x)=J_n(x)$ where
[af6de50]663        $J_n(x) = \frac{1}{\pi}\int_0^\pi \cos(n\tau - x\sin(\tau))\,d\tau$.
[48cd5b3]664        If $n$ = 0 or 1, it uses sas_J0(x) or sas_J1(x), respectively.
665
[ca1eaeb]666        The standard math function jn(n, x) is not available on all platforms.
667
[af6de50]668        :code:`source = ["lib/polevl.c", "lib/sas_J0.c", "lib/sas_J1.c", "lib/sas_JN.c", ...]`
[48cd5b3]669        (`link to Bessel function's code <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/sas_JN.c>`_)
[af6de50]670
[ca6cbc1c]671    sas_Si(x):
[af6de50]672        Sine integral $\text{Si}(x) = \int_0^x \tfrac{\sin t}{t}\,dt$.
673
[48cd5b3]674        This function uses Taylor series for small and large arguments:
675
676        For large arguments,
677
678        .. math::
[af6de50]679
[ca1eaeb]680             \text{Si}(x) \sim \frac{\pi}{2}
681             - \frac{\cos(x)}{x}\left(1 - \frac{2!}{x^2} + \frac{4!}{x^4} - \frac{6!}{x^6} \right)
682             - \frac{\sin(x)}{x}\left(\frac{1}{x} - \frac{3!}{x^3} + \frac{5!}{x^5} - \frac{7!}{x^7}\right)
[48cd5b3]683
684        For small arguments,
685
686        .. math::
687
[ca1eaeb]688           \text{Si}(x) \sim x
689           - \frac{x^3}{3\times 3!} + \frac{x^5}{5 \times 5!} - \frac{x^7}{7 \times 7!}
[48cd5b3]690           + \frac{x^9}{9\times 9!} - \frac{x^{11}}{11\times 11!}
691
692        :code:`source = ["lib/Si.c", ...]`
693        (`link to code <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/Si.c>`_)
694
[ca6cbc1c]695    sas_3j1x_x(x):
[af6de50]696        Spherical Bessel form
[ca1eaeb]697        $\text{sph_j1c}(x) = 3 j_1(x)/x = 3 (\sin(x) - x \cos(x))/x^3$,
698        with a limiting value of 1 at $x=0$, where $j_1(x)$ is the spherical
699        Bessel function of the first kind and first order.
[48cd5b3]700
701        This function uses a Taylor series for small $x$ for numerical accuracy.
[af6de50]702
[ca6cbc1c]703        :code:`source = ["lib/sas_3j1x_x.c", ...]`
704        (`link to code <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/sas_3j1x_x.c>`_)
[af6de50]705
[48cd5b3]706
[ca6cbc1c]707    sas_2J1x_x(x):
[ca1eaeb]708        Bessel form $\text{sas_J1c}(x) = 2 J_1(x)/x$, with a limiting value
709        of 1 at $x=0$, where $J_1(x)$ is the Bessel function of first kind
710        and first order.
[af6de50]711
[20cfa23]712        :code:`source = ["lib/polevl.c", "lib/sas_J1.c", ...]`
[48cd5b3]713        (`link to Bessel form's code <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/sas_J1.c>`_)
714
[af6de50]715
[48cd5b3]716    Gauss76Z[i], Gauss76Wt[i]:
717        Points $z_i$ and weights $w_i$ for 76-point Gaussian quadrature, respectively,
[ca1eaeb]718        computing $\int_{-1}^1 f(z)\,dz \approx \sum_{i=1}^{76} w_i\,f(z_i)$.
[af6de50]719
[48cd5b3]720        Similar arrays are available in :code:`gauss20.c` for 20-point
721        quadrature and in :code:`gauss150.c` for 150-point quadrature.
722
723        :code:`source = ["lib/gauss76.c", ...]`
724        (`link to code <https://github.com/SasView/sasmodels/tree/master/sasmodels/models/lib/gauss76.c>`_)
725
726
[af6de50]727
728Problems with C models
729......................
730
731The graphics processor (GPU) in your computer is a specialized computer tuned
732for certain kinds of problems.  This leads to strange restrictions that you
733need to be aware of.  Your code may work fine on some platforms or for some
734models, but then return bad values on other platforms.  Some examples of
735particular problems:
736
[907186d]737  **(1) Code is too complex, or uses too much memory.**  GPU devices only have a
[af6de50]738  limited amount of memory available for each processor.  If you run programs
739  which take too much memory, then rather than running multiple values in parallel
740  as it usually does, the GPU may only run a single version of the code at a
741  time, making it slower than running on the CPU.  It may fail to run on
742  some platforms, or worse, cause the screen to go blank or the system to reboot.
[05829fb]743
[907186d]744  **(2) Code takes too long.**  Because GPU devices are used for the computer
[af6de50]745  display, the OpenCL drivers are very careful about the amount of time they
746  will allow any code to run.  For example, on OS X, the model will stop running
[907186d]747  after 5 seconds regardless of whether the computation is complete.  You may end up
748  with only some of your 2D array defined, with the rest containing random
[af6de50]749  data. Or it may cause the screen to go blank or the system to reboot.
[05829fb]750
[907186d]751  **(3) Memory is not aligned**.  The GPU hardware is specialized to operate on
752  multiple values simultaneously.  To keep the GPU simple the values in memory
[af6de50]753  must be aligned with the different GPU compute engines.  Not following these
754  rules can lead to unexpected values being loaded into memory, and wrong answers
755  computed.  The conclusion from a very long and strange debugging session was
756  that any arrays that you declare in your model should be a multiple of four.
[907186d]757  For example::
[05829fb]758
[af6de50]759      double Iq(q, p1, p2, ...)
760      {
761          double vector[8];  // Only going to use seven slots, but declare 8
762          ...
763      }
[05829fb]764
[af6de50]765The first step when your model is behaving strangely is to set **single=False**.
[907186d]766This automatically restricts the model to only run on the CPU, or on high-end
767GPU cards.  There can still be problems even on high-end cards, so you can force
[af6de50]768the model off the GPU by setting **opencl=False**.  This runs the model
769as a normal C program without any GPU restrictions so you know that
770strange results are probably from your code rather than the environment.  Once
771the code is debugged, you can compare your output to the output on the GPU.
772
773Although it can be difficult to get your model to work on the GPU, the reward
774can be a model that runs 1000x faster on a good card.  Even your laptop may
775show a 50x improvement or more over the equivalent pure python model.
[05829fb]776
777External C Models
778.................
779
780External C models are very much like embedded C models, except that
781*Iq*, *Iqxy* and *form_volume* are defined in an external source file
[af6de50]782loaded using the *source=[...]* statement. You need to supply the function
[05829fb]783declarations for each of these that you need instead of building them
784automatically from the parameter table.
785
786
787.. _Form_Factors:
788
789Form Factors
790............
791
792Away from the dilute limit you can estimate scattering including
793particle-particle interactions using $I(q) = P(q)*S(q)$ where $P(q)$
794is the form factor and $S(q)$ is the structure factor.  The simplest
795structure factor is the *hardsphere* interaction, which
796uses the effective radius of the form factor as an input to the structure
797factor model.  The effective radius is the average radius of the
798form averaged over all the polydispersity values.
799
[31d7803]800::
801
802    def ER(radius, thickness):
803        """Effective radius of a core-shell sphere."""
804        return radius + thickness
805
806Now consider the *core_shell_sphere*, which has a simple effective radius
[05829fb]807equal to the radius of the core plus the thickness of the shell, as
808shown above. Given polydispersity over *(r1, r2, ..., rm)* in radius and
809*(t1, t2, ..., tn)* in thickness, *ER* is called with a mesh
810grid covering all possible combinations of radius and thickness.
811That is, *radius* is *(r1, r2, ..., rm, r1, r2, ..., rm, ...)*
812and *thickness* is *(t1, t1, ... t1, t2, t2, ..., t2, ...)*.
813The *ER* function returns one effective radius for each combination.
814The effective radius calculator weights each of these according to
815the polydispersity distributions and calls the structure factor
816with the average *ER*.
817
818::
819
820    def VR(radius, thickness):
821        """Sphere and shell volumes for a core-shell sphere."""
822        whole = 4.0/3.0 * pi * (radius + thickness)**3
823        core = 4.0/3.0 * pi * radius**3
824        return whole, whole - core
825
826Core-shell type models have an additional volume ratio which scales
827the structure factor.  The *VR* function returns the volume of
828the whole sphere and the volume of the shell. Like *ER*, there is
829one return value for each point in the mesh grid.
830
[31d7803]831*NOTE: we may be removing or modifying this feature soon. As of the
832time of writing, core-shell sphere returns (1., 1.) for VR, giving a volume
833ratio of 1.0.*
[05829fb]834
835Unit Tests
836..........
837
838THESE ARE VERY IMPORTANT. Include at least one test for each model and
839PLEASE make sure that the answer value is correct (i.e. not a random number).
840
841::
842
843    tests = [
844        [{}, 0.2, 0.726362],
845        [{"scale": 1., "background": 0., "sld": 6., "sld_solvent": 1.,
846          "radius": 120., "radius_pd": 0.2, "radius_pd_n":45},
847         0.2, 0.228843],
848        [{"radius": 120., "radius_pd": 0.2, "radius_pd_n":45}, "ER", 120.],
849        [{"radius": 120., "radius_pd": 0.2, "radius_pd_n":45}, "VR", 1.],
850    ]
851
852
853**tests=[[{parameters}, q, result], ...]** is a list of lists.
854Each list is one test and contains, in order:
855
856- a dictionary of parameter values. This can be {} using the default
857  parameters, or filled with some parameters that will be different
[cbbb6a4]858  from the default, such as {‘radius’:10.0, ‘sld’:4}. Unlisted parameters
[05829fb]859  will be given the default values.
860- the input $q$ value or tuple of $(q_x, q_y)$ values.
861- the output $I(q)$ or $I(q_x,q_y)$ expected of the model for the parameters
862  and input value given.
863- input and output values can themselves be lists if you have several
864  $q$ values to test for the same model parameters.
865- for testing *ER* and *VR*, give the inputs as "ER" and "VR" respectively;
866  the output for *VR* should be the sphere/shell ratio, not the individual
867  sphere and shell values.
868
869.. _Test_Your_New_Model:
870
871Test Your New Model
872^^^^^^^^^^^^^^^^^^^
873
[3e1c9e5]874Minimal Testing
875...............
[e925f61]876
[3e1c9e5]877Either open the :ref:`Python_shell` (*Tools* > *Python Shell/Editor*) or the :ref:`Advanced_Plugin_Editor` (*Fitting* > *Plugin Model Operations* > *Advanced
878Plugin Editor*), load your model, and then select *Run > Check Model* from the
879menu bar.
[05829fb]880
[3e1c9e5]881An *Info* box will appear with the results of the compilation and a check that
882the model runs.
[e925f61]883
[3e1c9e5]884Recommended Testing
885...................
[e925f61]886
[05829fb]887If the model compiles and runs, you can next run the unit tests that
[31d7803]888you have added using the **test =** values. Switch to the *Shell* tab
[05829fb]889and type the following::
890
891    from sasmodels.model_test import run_one
892    run_one("~/.sasview/plugin_models/model.py")
893
894This should print::
895
896    test_model_python (sasmodels.model_test.ModelTestCase) ... ok
897
898To check whether single precision is good enough, type the following::
899
900    from sasmodels.compare import main
901    main("~/.sasview/plugin_models/model.py")
902
903This will pop up a plot showing the difference between single precision
904and double precision on a range of $q$ values.
905
906::
907
908  demo = dict(scale=1, background=0,
909              sld=6, sld_solvent=1,
910              radius=120,
911              radius_pd=.2, radius_pd_n=45)
912
913**demo={'par': value, ...}** in the model file sets the default values for
914the comparison. You can include polydispersity parameters such as
915*radius_pd=0.2, radius_pd_n=45* which would otherwise be zero.
916
917The options to compare are quite extensive; type the following for help::
918
919    main()
920
921Options will need to be passed as separate strings.
922For example to run your model with a random set of parameters::
923
924    main("-random", "-pars", "~/.sasview/plugin_models/model.py")
925
926For the random models,
927
[31d7803]928- *sld* will be in the range (-0.5,10.5),
929- angles (*theta, phi, psi*) will be in the range (-180,180),
930- angular dispersion will be in the range (0,45),
931- polydispersity will be in the range (0,1)
932- other values will be in the range (0, 2\ *v*), where *v* is the value of the parameter in demo.
[05829fb]933
[31d7803]934Dispersion parameters *n*\, *sigma* and *type* will be unchanged from demo so that
[05829fb]935run times are predictable.
936
937If your model has 2D orientational calculation, then you should also
938test with::
939
940    main("-2d", "~/.sasview/plugin_models/model.py")
941
942
[e925f61]943Clean Lint - (Developer Version Only)
944^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
[05829fb]945
[e925f61]946**NB: For now we are not providing pylint with the installer version of SasView;
947so unless you have a SasView build environment available, you can ignore this section!**
[05829fb]948
949Run the lint check with::
950
951    python -m pylint --rcfile=extra/pylint.rc ~/.sasview/plugin_models/model.py
952
953We are not aiming for zero lint just yet, only keeping it to a minimum.
954For now, don't worry too much about *invalid-name*. If you really want a
955variable name *Rg* for example because $R_g$ is the right name for the model
956parameter then ignore the lint errors.  Also, ignore *missing-docstring*
957for standard model functions *Iq*, *Iqxy*, etc.
958
[31d7803]959We will have delinting sessions at the SasView Code Camps, where we can
[05829fb]960decide on standards for model files, parameter names, etc.
961
[31d7803]962For now, you can tell pylint to ignore things.  For example, to align your
[05829fb]963parameters in blocks::
964
965    # pylint: disable=bad-whitespace,line-too-long
966    #   ["name",                  "units", default, [lower, upper], "type", "description"],
967    parameters = [
968        ["contrast_factor",       "barns",    10.0,  [-inf, inf], "", "Contrast factor of the polymer"],
969        ["bjerrum_length",        "Ang",       7.1,  [0, inf],    "", "Bjerrum length"],
970        ["virial_param",          "1/Ang^2",  12.0,  [-inf, inf], "", "Virial parameter"],
971        ["monomer_length",        "Ang",      10.0,  [0, inf],    "", "Monomer length"],
972        ["salt_concentration",    "mol/L",     0.0,  [-inf, inf], "", "Concentration of monovalent salt"],
973        ["ionization_degree",     "",          0.05, [0, inf],    "", "Degree of ionization"],
974        ["polymer_concentration", "mol/L",     0.7,  [0, inf],    "", "Polymer molar concentration"],
975        ]
976    # pylint: enable=bad-whitespace,line-too-long
977
978Don't put in too many pylint statements, though, since they make the code ugly.
979
[e925f61]980Check The Docs - (Developer Version Only)
981^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
[05829fb]982
983You can get a rough idea of how the documentation will look using the
984following::
985
986    from sasmodels.generate import view_html
987    view_html('~/.sasview/plugin_models/model.py')
988
989This does not use the same styling as the SasView docs, but it will allow
990you to check that your ReStructuredText and LaTeX formatting.  Here are
991some tools to help with the inevitable syntax errors:
992
993- `Sphinx cheat sheet <http://matplotlib.org/sampledoc/cheatsheet.html>`_
994- `Sphinx Documentation <http://www.sphinx-doc.org/en/stable/>`_
995- `MathJax <http://www.mathjax.org/>`_
996- `amsmath <http://www.ams.org/publications/authors/tex/amslatex>`_
997
[31d7803]998There is also a neat online WYSIWYG ReStructuredText editor at http://rst.ninjs.org\ .
999
[e925f61]1000Share Your Model!
1001^^^^^^^^^^^^^^^^^
[05829fb]1002
1003Once compare and the unit test(s) pass properly and everything is done,
1004consider adding your model to the
[e925f61]1005`Model Marketplace <http://marketplace.sasview.org/>`_ so that others may use it!
[3e1c9e5]1006
1007.. ZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZZ
1008
[907186d]1009.. note::  This help document was last changed by Steve King, 25Oct2016
Note: See TracBrowser for help on using the repository browser.