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

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.2release_4.0.1ticket-1009ticket-1094-headlessticket-1242-2d-resolutionticket-1243ticket-1249ticket885unittest-saveload
Last change on this file since 3d164b9 was 3d164b9, checked in by smk78, 8 years ago

Under review - minor changes

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