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

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 e4078ed was 3e1c9e5, checked in by smk78, 8 years ago

Disinfected menu_bar.rst, fitting_help.rst and plugin.rst of 'custom
model' in favour of 'plugin model' and re-wrote sections of the
documentation to match it to the functionality of the program with
particular emphasis on creating/testing plugin models.

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