source: trunk/GSASIIobj.py @ 1338

Last change on this file since 1338 was 1338, checked in by vondreele, 8 years ago

put in SASD parameter descriptions in CompileVarDesc?; required mods for VarDescr? & getVarDescr

  • Property svn:eol-style set to native
  • Property svn:keywords set to Date Author Revision URL Id
File size: 95.0 KB
Line 
1# -*- coding: utf-8 -*-
2#GSASIIobj - data objects for GSAS-II
3########### SVN repository information ###################
4# $Date: 2014-05-08 15:45:45 +0000 (Thu, 08 May 2014) $
5# $Author: vondreele $
6# $Revision: 1338 $
7# $URL: trunk/GSASIIobj.py $
8# $Id: GSASIIobj.py 1338 2014-05-08 15:45:45Z vondreele $
9########### SVN repository information ###################
10
11'''
12*GSASIIobj: Data objects*
13=========================
14
15This module defines and/or documents the data structures used in GSAS-II, as well
16as provides misc. support routines.
17
18Constraints Tree Item
19----------------------
20
21.. _Constraints_table:
22
23.. index::
24   single: Constraints object description
25   single: Data object descriptions; Constraints
26
27Constraints are stored in a dict, separated into groups.
28Note that parameter are named in the following pattern,
29p:h:<var>:n, where p is the phase number, h is the histogram number
30<var> is a variable name and n is the parameter number.
31If a parameter does not depend on a histogram or phase or is unnumbered, that
32number is omitted.
33Note that the contents of each dict item is a List where each element in the
34list is a :ref:`constraint definition objects <Constraint_definitions_table>`.
35The constraints in this form are converted in
36:func:`GSASIIstrIO.ProcessConstraints` to the form used in :mod:`GSASIImapvars`
37
38The keys in the Constraints dict are:
39
40.. tabularcolumns:: |l|p{4.5in}|
41
42==========  ====================================================
43  key         explanation
44==========  ====================================================
45Hist        This specifies a list of constraints on
46            histogram-related parameters,
47            which will be of form :h:<var>:n.
48HAP         This specifies a list of constraints on parameters
49            that are defined for every histogram in each phase
50            and are of form p:h:<var>:n.           
51Phase       This specifies a list of constraints on phase
52            parameters,
53            which will be of form p::<var>:n.
54Global      This specifies a list of constraints on parameters
55            that are not tied to a histogram or phase and
56            are of form ::<var>:n
57==========  ====================================================
58
59.. _Constraint_definitions_table:
60
61.. index::
62   single: Constraint definition object description
63   single: Data object descriptions; Constraint Definition
64
65Each constraint is defined as an item in a list. Each constraint is of form::
66
67[[<mult1>, <var1>], [<mult2>, <var2>],..., <fixedval>, <varyflag>, <constype>]
68
69Where the variable pair list item containing two values [<mult>, <var>], where:
70
71  * <mult> is a multiplier for the constraint (float)
72  * <var> a :class:`G2VarObj` object (previously a str variable name of form
73      'p:h:name[:at]')
74
75Note that the last three items in the list play a special role:
76
77 * <fixedval> is the fixed value for a `constant equation` (``constype=c``)
78   constraint or is None. For a `New variable` (``constype=f``) constraint,
79   a variable name can be specified as a str (used for externally
80   generated constraints)
81 * <varyflag> is True or False for `New variable` (``constype=f``) constraints
82   or is None. This will be implemented in the future to indicate if these variables
83   should be refined.
84 * <constype> is one of four letters, 'e', 'c', 'h', 'f' that determines the type of constraint:
85
86    * 'e' defines a set of equivalent variables. Only the first variable is refined (if the
87      appropriate refine flag is set) and and all other equivalent variables in the list
88      are generated from that variable, using the appropriate multipliers.
89    * 'c' defines a constraint equation of form,
90      :math:`m_1 \\times var_1 + m_2 \\times var_2 + ... = c`
91    * 'h' defines a variable to hold (not vary). Any variable on this list is not varied,
92      even if its refinement flag is set. Only one [mult,var] pair is allowed in a hold
93      constraint and the mult value is ignored.
94      This is of particular value when needing to hold one or more variables where a
95      single flag controls a set of variables such as, coordinates,
96      the reciprocal metric tensor or anisotropic displacement parameter.
97    * 'f' defines a new variable (function) according to relationship
98      :math:`newvar = m_1 \\times var_1 + m_2 \\times var_2 + ...`
99
100Covariance Tree Item
101--------------------
102
103.. _Covariance_table:
104
105.. index::
106   single: Covariance description
107   single: Data object descriptions; Covariance
108
109The Covariance tree item has results from the last least-squares run. They
110are stored in a dict with these keys:
111
112.. tabularcolumns:: |l|l|p{4in}|
113
114=============  ===============  ====================================================
115  key            sub-key        explanation
116=============  ===============  ====================================================
117newCellDict    \                dict with lattice parameters computed by
118                                :func:`GSASIIstrMath.GetNewCellParms` (dict)
119title          \                Name of gpx file(?) (str)
120variables      \                Values for all N refined variables
121                                (list of float values, length N,
122                                ordered to match varyList)
123sig            \                Uncertainty values for all N refined variables
124                                (list of float values, length N,
125                                ordered to match varyList)
126varyList       \                List of directly refined variables
127                                (list of str values, length N)
128newAtomDict    \                dict with atom position values computed in
129                                :func:`GSASIIstrMath.ApplyXYZshifts` (dict)
130Rvals          \                R-factors, GOF, Marquardt value for last
131                                refinement cycle (dict)
132\              Nobs             Number of observed data points (int)
133\              Rwp              overall weighted profile R-factor (%, float)
134\              chisq            sum[w*(Iobs-Icalc)**2] for all data
135                                note this is not the reduced chi squared (float)
136\              lamMax           Marquardt value applied to Hessian diagonal
137                                (float)
138\              GOF              The goodness-of-fit, aka square root of
139                                the reduced chi squared. (float)
140covMatrix      \                The (NxN) covVariance matrix (np.array)
141=============  ===============  ====================================================
142
143Phase Tree Items
144----------------
145
146.. _Phase_table:
147
148.. index::
149   single: Phase object description
150   single: Data object descriptions; Phase
151
152Phase information is stored in the GSAS-II data tree as children of the
153Phases item in a dict with keys:
154
155.. tabularcolumns:: |l|l|p{4in}|
156
157==========  ===============  ====================================================
158  key         sub-key        explanation
159==========  ===============  ====================================================
160General         \            Overall information for the phase (dict)
161  \         AtomPtrs         list of four locations to use to pull info
162                             from the atom records (list)
163  \         F000X            x-ray F(000) intensity (float)
164  \         F000N            neutron F(000) intensity (float)
165  \         Mydir            directory of current .gpx file (str)
166  \         MCSA controls    Monte Carlo-Simulated Annealing controls (dict)
167  \         Cell             List with 8 items: cell refinement flag (bool)
168                             a, b, c, (Angstrom, float)
169                             alpha, beta & gamma (degrees, float)
170                             volume (A^3, float)
171  \         Type             'nuclear' or 'macromolecular' for now (str)
172  \         Map              dict of map parameters
173  \         SH Texture       dict of spherical harmonic preferred orientation
174                             parameters
175  \         Isotope          dict of isotopes for each atom type
176  \         Isotopes         dict of scattering lengths for each isotope
177                             combination for each element in phase 
178  \         Name             phase name (str)
179  \         SGData           Space group details as a :ref:`space group (SGData) object <SGData_table>`
180                             as defined in :func:`GSASIIspc.SpcGroup`.
181  \         Pawley neg wt    Restraint value for negative Pawley intensities
182                             (float)
183  \         Flip             dict of Charge flip controls
184  \         Data plot type   data plot type ('Mustrain', 'Size' or
185                             'Preferred orientation') for powder data (str)
186  \         Mass             Mass of unit cell contents in g/mol
187  \         POhkl            March-Dollase preferred orientation direction
188  \         Z                dict of atomic numbers for each atom type
189  \         vdWRadii         dict of van der Waals radii for each atom type
190  \         Color            Colors for atoms (list of (r,b,g) triplets)
191  \         AtomTypes        List of atom types
192  \         AtomMass         List of masses for atoms
193  \         doPawley         Flag for Pawley intensity extraction (bool)
194  \         NoAtoms          Number of atoms per unit cell of each type (dict)
195  \         Pawley dmin      maximum Q (as d-space) to use for Pawley
196                             extraction (float)
197  \         BondRadii        Default radius for each atom used to compute
198                             interatomic distances (list of floats)
199  \         AngleRadii       Default radius for each atom used to compute
200                             interatomic angles (list of floats)
201  \         DisAglCtls       Dict with distance/angle search controls,
202                             which has keys 'Name', 'AtomTypes',
203                             'BondRadii', 'AngleRadii' which are as above
204                             except are possibly edited. Also contains
205                             'Factors', which is a 2 element list with
206                             a multiplier for bond and angle search range
207                             [typically (0.85,0.85)].
208ranId           \            unique random number Id for phase (int)
209pId             \            Phase Id number for current project (int).
210Atoms           \            Atoms in phase as a list of lists. The outer list
211                             is for each atom, the inner list contains varying
212                             items depending on the type of phase, see
213                             the :ref:`Atom Records <Atoms_table>` description.
214                             (list of lists)
215Drawing         \            Display parameters (dict)
216\           ballScale        Size of spheres in ball-and-stick display (float)
217\           bondList         dict with bonds
218\           contourLevel     map contour level in e/A^3 (float)
219\           showABC          Flag to show view point triplet (bool). True=show.
220\           viewDir          cartesian viewing direction (np.array with three
221                             elements)
222\           Zclip            clipping distance in A (float)
223\           backColor        background for plot as and R,G,B triplet
224                             (default = [0, 0, 0], black).
225                             (list with three atoms)
226\           selectedAtoms    List of selected atoms (list of int values)
227\           showRigidBodies  Flag to highlight rigid body placement
228\           sizeH            Size ratio for H atoms (float)
229\           bondRadius       Size of binds in A (float)
230\           atomPtrs         positions of x, type, site sym, ADP flag in Draw Atoms (list)
231\           viewPoint        list of lists. First item in list is [x,y,z]
232                             in fractional coordinates for the center of
233                             the plot. Second item list of previous & current
234                             atom number viewed (may be [0,0])
235\           showHydrogen     Flag to control plotting of H atoms.
236\           unitCellBox      Flag to control display of the unit cell.
237\           ellipseProb      Probability limit for display of thermal
238                             ellipsoids in % (float).
239\           vdwScale         Multiplier of van der Waals radius for
240                             display of vdW spheres.
241\           Atoms            A list of lists with an entry for each atom
242                             that is plotted.
243\           Zstep            Step to de/increase Z-clip (float)
244\           Quaternion       Viewing quaternion (4 element np.array)
245\           radiusFactor     Distance ratio for searching for bonds. ? Bonds
246                             are located that are within r(Ra+Rb) and (Ra+Rb)/r
247                             where Ra and Rb are the atomic radii.
248\           oldxy            previous view point (list with two floats)
249\           cameraPos        Viewing position in A for plot (float)
250\           depthFog         True if use depthFog on plot - set currently as False (bool)
251RBModels        \            Rigid body assignments (note Rigid body definitions
252                             are stored in their own main top-level tree entry.)
253Pawley ref      \            Pawley reflections
254Histograms      \            A dict of dicts. The key for the outer dict is
255                             the histograms tied to this phase. The inner
256                             dict contains the combined phase/histogram
257                             parameters for items such as scale factors,
258                             size and strain parameters. (dict)
259MCSA            \            Monte-Carlo simulated annealing parameters (dict)
260\           
261==========  ===============  ====================================================
262
263Rigid Body Objects
264------------------
265
266.. _RBData_table:
267
268.. index::
269   single: Rigid Body Data description
270   single: Data object descriptions; Rigid Body Data
271   
272Rigid body descriptions are available for two types of rigid bodies: 'Vector'
273and 'Residue'. Vector rigid bodies are developed by a sequence of translations each
274with a refinable magnitude and Residue rigid bodies are described as Cartesian coordinates
275with defined refinable torsion angles.
276
277.. tabularcolumns:: |l|l|p{4in}|
278
279==========  ===============  ====================================================
280  key         sub-key        explanation
281==========  ===============  ====================================================
282Vector      RBId             vector rigid bodies (dict of dict)
283\           AtInfo           Drad, Color: atom drawing radius & color for each atom type (dict)
284\           RBname           Name assigned by user to rigid body (str)
285\           VectMag          vector magnitudes in A (list)
286\           rbXYZ            Cartesian coordinates for Vector rigid body (list of 3 float)
287\           rbRef            3 assigned reference atom nos. in rigid body for origin
288                             definition, use center of atoms flag (list of 3 int & 1 bool)
289\           VectRef          refinement flags for VectMag values (list of bool)
290\           rbTypes          Atom types for each atom in rigid body (list of str)
291\           rbVect           Cartesian vectors for each translation used to build rigid body (list of lists)
292\           useCount         Number of times rigid body is used in any structure (int)
293Residue     RBId             residue rigid bodies (dict of dict)
294\           AtInfo           Drad, Color: atom drawing radius & color for each atom type(dict)
295\           RBname           Name assigned by user to rigid body (str)
296\           rbXYZ            Cartesian coordinates for Residue rigid body (list of 3 float)
297\           rbTypes          Atom types for each atom in rigid body (list of str)
298\           atNames          Names of each atom in rigid body (e.g. C1,N2...) (list of str)
299\           rbRef            3 assigned reference atom nos. in rigid body for origin
300                             definition, use center of atoms flag (list of 3 int & 1 bool)
301\           rbSeq            Orig,Piv,angle,Riding (list): definition of internal rigid body
302                             torsion; origin atom (int), pivot atom (int), torsion angle (float),
303                             riding atoms (list of int)
304\           SelSeq           [int,int] used by SeqSizer to identify objects
305\           useCount         Number of times rigid body is used in any structure (int)
306RBIds           \            unique Ids generated upon creation of each rigid body (dict)
307\           Vector           Ids for each Vector rigid body (list)
308\           Residue          Ids for each Residue rigid body (list)
309==========  ===============  ====================================================
310
311Space Group Objects
312-------------------
313
314.. _SGData_table:
315
316.. index::
317   single: Space Group Data description
318   single: Data object descriptions; Space Group Data
319
320Space groups are interpreted by :func:`GSASIIspc.SpcGroup`
321and the information is placed in a SGdata object
322which is a dict with these keys:
323
324.. tabularcolumns:: |l|p{4.5in}|
325
326==========  ====================================================
327  key         explanation
328==========  ====================================================
329SpGrp       space group symbol (str)
330Laue        one of the following 14 Laue classes:
331            -1, 2/m, mmm, 4/m, 4/mmm, 3R,
332            3mR, 3, 3m1, 31m, 6/m, 6/mmm, m3, m3m (str)
333SGInv       True if centrosymmetric, False if not (bool)
334SGLatt      Lattice centering type. Will be one of
335            P, A, B, C, I, F, R (str)
336SGUniq      unique axis if monoclinic. Will be
337            a, b, or c for monoclinic space groups.
338            Will be blank for non-monoclinic. (str)
339SGCen       Symmetry cell centering vectors. A (n,3) np.array
340            of centers. Will always have at least one row:
341            ``np.array([[0, 0, 0]])``
342SGOps       symmetry operations as a list of form
343            ``[[M1,T1], [M2,T2],...]``
344            where :math:`M_n` is a 3x3 np.array
345            and :math:`T_n` is a length 3 np.array.
346            Atom coordinates are transformed where the
347            Asymmetric unit coordinates [X is (x,y,z)]
348            are transformed using
349            :math:`X^\prime = M_n*X+T_n`
350SGSys       symmetry unit cell: type one of
351            'triclinic', 'monoclinic', 'orthorhombic',
352            'tetragonal', 'rhombohedral', 'trigonal',
353            'hexagonal', 'cubic' (str)
354SGPolax     Axes for space group polarity. Will be one of
355            '', 'x', 'y', 'x y', 'z', 'x z', 'y z',
356            'xyz'. In the case where axes are arbitrary
357            '111' is used (P 1, and ?).
358==========  ====================================================
359
360Atom Records
361------------
362
363.. _Atoms_table:
364
365.. index::
366   single: Atoms record description
367   single: Data object descriptions; Atoms record
368
369
370If ``phasedict`` points to the phase information in the data tree, then
371atoms are contained in a list of atom records (list) in
372``phasedict['Atoms']``. Also needed to read atom information
373are four pointers, ``cx,ct,cs,cia = phasedict['General']['atomPtrs']``,
374which define locations in the atom record, as shown below. Items shown are
375always present; additional ones for macromolecular phases are marked 'mm'
376
377.. tabularcolumns:: |l|p{4.5in}|
378
379==============   ====================================================
380location         explanation
381==============   ====================================================
382ct-4              mm - residue number (str)
383ct-3              mm - residue name (e.g. ALA) (str)
384ct-2              mm - chain label (str)
385ct-1              atom label (str)
386ct                atom type (str)
387ct+1              refinement flags; combination of 'F', 'X', 'U' (str)
388cx,cx+1,cx+2      the x,y and z coordinates (3 floats)
389cs                site symmetry (str)
390cs+1              site multiplicity (int)
391cia               ADP flag: Isotropic ('I') or Anisotropic ('A')
392cia+1             Uiso (float)
393cia+2...cia+6     U11, U22, U33, U12, U13, U23 (6 floats)
394atom[-1]                unique atom identifier (int)
395==============   ====================================================
396
397Drawing Atom Records
398--------------------
399
400.. _Drawing_atoms_table:
401
402.. index::
403   single: Drawing atoms record description
404   single: Data object descriptions; Drawing atoms record
405
406
407If ``phasedict`` points to the phase information in the data tree, then
408drawing atoms are contained in a list of drawing atom records (list) in
409``phasedict['Drawing']['Atoms']``. Also needed to read atom information
410are four pointers, ``cx,ct,cs,ci = phasedict['Drawing']['AtomPtrs']``,
411which define locations in the atom record, as shown below. Items shown are
412always present; additional ones for macromolecular phases are marked 'mm'
413
414.. tabularcolumns:: |l|p{4.5in}|
415
416==============   ====================================================
417location         explanation
418==============   ====================================================
419ct-4              mm - residue number (str)
420ct-3              mm - residue name (e.g. ALA) (str)
421ct-2              mm - chain label (str)
422ct-1              atom label (str)
423ct                atom type (str)
424cx,cx+1,cx+2      the x,y and z coordinates (3 floats)
425cs-1              Sym Op symbol; sym. op number + unit cell id (e.g. '1,0,-1') (str)
426cs                atom drawing style; e.g. 'balls & sticks' (str)
427cs+1              atom label style (e.g. 'name') (str)
428cs+2              atom color (RBG triplet) (int)
429cs+3              ADP flag: Isotropic ('I') or Anisotropic ('A')
430cs+4              Uiso (float)
431cs+5...cs+11      U11, U22, U33, U12, U13, U23 (6 floats)
432ci                unique atom identifier; matches source atom Id in Atom Records (int)
433==============   ====================================================
434
435Powder Diffraction Tree Items
436-----------------------------
437
438.. _Powder_table:
439
440.. index::
441   single: Powder data object description
442   single: Data object descriptions; Powder Data
443
444Every powder diffraction histogram is stored in the GSAS-II data tree
445with a top-level entry named beginning with the string "PWDR ". The
446diffraction data for that information are directly associated with
447that tree item and there are a series of children to that item. The
448routines :func:`GSASII.GSASII.GetUsedHistogramsAndPhasesfromTree`
449and :func:`GSASIIstrIO.GetUsedHistogramsAndPhases` will
450load this information into a dictionary where the child tree name is
451used as a key, and the information in the main entry is assigned
452a key of ``Data``, as outlined below.
453
454.. tabularcolumns:: |l|l|p{4in}|
455
456======================  ===============  ====================================================
457  key                      sub-key        explanation
458======================  ===============  ====================================================
459Comments                      \           Text strings extracted from the original powder
460                                          data header. These cannot be changed by the user;
461                                          it may be empty.
462Limits                       \            A list of two two element lists, as [[Ld,Hd],[L,H]]
463                                          where L and Ld are the current and default lowest
464                                          two-theta value to be used and
465                                          where H and Hd are the current and default highest
466                                          two-theta value to be used.
467Reflection Lists              \           A dict with an entry for each phase in the
468                                          histogram. The contents of each dict item
469                                          is a dict containing reflections, as described in
470                                          the :ref:`Powder Reflections <PowderRefl_table>`
471                                          description.
472Instrument Parameters         \           A list containing two dicts where the possible
473                                          keys in each dict are listed below. The value
474                                          for each item is a list containing three values:
475                                          the initial value, the current value and a
476                                          refinement flag which can have a value of
477                                          True, False or 0 where 0 indicates a value that
478                                          cannot be refined. The first and second
479                                          values are floats unless otherwise noted.
480                                          Items in the first dict are noted as [1]
481\                         Lam             Specifies a wavelength in Angstroms [1]
482\                         Lam1            Specifies the primary wavelength in
483                                          Angstrom, when an alpha1, alpha2
484                                          source is used [1]
485\                         Lam2            Specifies the secondary wavelength in
486                                          Angstrom, when an alpha1, alpha2
487                                          source is used [1]
488                          I(L2)/I(L1)     Ratio of Lam2 to Lam1 [1]           
489\                         Type            Histogram type (str) [1]:
490                                           * 'PXC' for constant wavelength x-ray
491                                           * 'PNC' for constant wavelength neutron
492                                           * 'PNT' for time of flight neutron
493\                         Zero            Two-theta zero correction in *degrees* [1]
494\                         Azimuth         Azimuthal setting angle for data recorded
495                                          with differing setting angles [1]
496\                         U, V, W         Cagliotti profile coefficients
497                                          for Gaussian instrumental broadening, where the
498                                          FWHM goes as
499                                          :math:`U \\tan^2\\theta + V \\tan\\theta + W` [1]
500\                         X, Y            Cauchy (Lorentzian) instrumental broadening
501                                          coefficients [1]
502\                         SH/L            Variant of the Finger-Cox-Jephcoat asymmetric
503                                          peak broadening ratio. Note that this is the
504                                          average between S/L and H/L where S is
505                                          sample height, H is the slit height and
506                                          L is the goniometer diameter. [1]
507\                         Polariz.        Polarization coefficient. [1]
508wtFactor                      \           A weighting factor to increase or decrease
509                                          the leverage of data in the histogram (float).
510                                          A value of 1.0 weights the data with their
511                                          standard uncertainties and a larger value
512                                          increases the weighting of the data (equivalent
513                                          to decreasing the uncertainties).
514Sample Parameters             \           Specifies a dict with parameters that describe how
515                                          the data were collected, as listed
516                                          below. Refinable parameters are a list containing
517                                          a float and a bool, where the second value
518                                          specifies if the value is refined, otherwise
519                                          the value is a float unless otherwise noted.
520\                         Scale           The histogram scale factor (refinable)
521\                         Absorption      The sample absorption coefficient as
522                                          :math:`\\mu r` where r is the radius
523                                          (refinable). Only valid for Debye-Scherrer geometry.
524\                         SurfaceRoughA   Surface roughness parameter A as defined by
525                                          Surotti,J. Appl. Cryst, 5,325-331, 1972.(refinable -
526                                          only valid for Bragg-Brentano geometry)                                         
527\                         SurfaceRoughB   Surface roughness parameter B (refinable -
528                                          only valid for Bragg-Brentano geometry)                                         
529\                         DisplaceX,      Sample displacement from goniometer center
530                          DisplaceY       where Y is along the beam direction and
531                                          X is perpendicular. Units are :math:`\\mu m`
532                                          (refinable).
533\                         Phi, Chi,       Goniometer sample setting angles, in degrees.
534                          Omega
535\                         Gonio. radius   Radius of the diffractometer in mm
536\                         InstrName       A name for the instrument, used in preparing
537                                          a CIF (str).
538\                         Force,          Variables that describe how the measurement
539                          Temperature,    was performed. Not used directly in
540                          Humidity,       any computations.
541                          Pressure,
542                          Voltage
543\                         ranId           The random-number Id for the histogram
544                                          (same value as where top-level key is ranId)
545\                         Type            Type of diffraction data, may be 'Debye-Scherrer'
546                                          or 'Bragg-Brentano' (str).
547\                         Diffuse         not in use?
548hId                           \           The number assigned to the histogram when
549                                          the project is loaded or edited (can change)
550ranId                         \           A random number id for the histogram
551                                          that does not change
552Background                    \           The background is stored as a list with where
553                                          the first item in the list is list and the second
554                                          item is a dict. The list contains the background
555                                          function and its coefficients; the dict contains
556                                          Debye diffuse terms and background peaks.
557                                          (TODO: this needs to be expanded.)
558Data                          \           The data consist of a list of 6 np.arrays
559                                          containing in order:
560
561                                           0. the x-postions (two-theta in degrees),
562                                           1. the intensity values (Yobs),
563                                           2. the weights for each Yobs value
564                                           3. the computed intensity values (Ycalc)
565                                           4. the background values
566                                           5. Yobs-Ycalc
567======================  ===============  ====================================================
568
569Powder Reflection Data Structure
570--------------------------------
571
572.. _PowderRefl_table:
573
574.. index::
575   single: Powder reflection object description
576   single: Data object descriptions; Powder Reflections
577   
578For every phase in a histogram, the ``Reflection Lists`` value is a dict
579one element of which is `'RefList'`, which is a np.array containing
580reflections. The columns in that array are documented below.
581
582==========  ====================================================
583  index         explanation
584==========  ====================================================
585 0,1,2       h,k,l (float)
586 3           multiplicity
587 4           d-space, Angstrom
588 5           pos, two-theta
589 6           sig, Gaussian width
590 7           gam, Lorenzian width
591 8           :math:`F_{obs}^2`
592 9           :math:`F_{calc}^2`
593 10          reflection phase, in degrees
594 11          intensity correction for reflection, this times
595             :math:`F_{obs}^2` or :math:`F_{calc}^2` gives Iobs or Icalc
596==========  ====================================================
597
598Single Crystal Tree Items
599-------------------------
600
601.. _Xtal_table:
602
603.. index::
604   single: Single Crystal data object description
605   single: Data object descriptions; Single crystal data
606
607Every single crystal diffraction histogram is stored in the GSAS-II data tree
608with a top-level entry named beginning with the string "HKLF ". The
609diffraction data for that information are directly associated with
610that tree item and there are a series of children to that item. The
611routines :func:`GSASII.GSASII.GetUsedHistogramsAndPhasesfromTree`
612and :func:`GSASIIstrIO.GetUsedHistogramsAndPhases` will
613load this information into a dictionary where the child tree name is
614used as a key, and the information in the main entry is assigned
615a key of ``Data``, as outlined below.
616
617.. tabularcolumns:: |l|l|p{4in}|
618
619======================  ===============  ====================================================
620  key                      sub-key        explanation
621======================  ===============  ====================================================
622Data                          \           A dict that contains the
623                                          reflection table,
624                                          as described in the
625                                          :ref:`Single Crystal Reflections
626                                          <XtalRefl_table>`
627                                          description.
628
629Instrument Parameters         \           A list containing two dicts where the possible
630                                          keys in each dict are listed below. The value
631                                          for most items is a list containing two values:
632                                          the initial value, the current value.
633                                          The first and second
634                                          values are floats unless otherwise noted.
635\                         Lam             Specifies a wavelength in Angstroms (two floats)
636\                         Type            Histogram type (two str values):
637                                           * 'SXC' for constant wavelength x-ray
638                                           * 'SNC' for constant wavelength neutron
639                                           * 'SNT' for time of flight neutron
640\                         InstrName       A name for the instrument, used in preparing
641                                          a CIF (str).
642
643wtFactor                      \           A weighting factor to increase or decrease
644                                          the leverage of data in the histogram (float).
645                                          A value of 1.0 weights the data with their
646                                          standard uncertainties and a larger value
647                                          increases the weighting of the data (equivalent
648                                          to decreasing the uncertainties).
649
650hId                           \           The number assigned to the histogram when
651                                          the project is loaded or edited (can change)
652ranId                         \           A random number id for the histogram
653                                          that does not change
654======================  ===============  ====================================================
655
656Single Crystal Reflection Data Structure
657----------------------------------------
658
659.. _XtalRefl_table:
660
661.. index::
662   single: Single Crystal reflection object description
663   single: Data object descriptions; Single Crystal Reflections
664   
665For every single crystal a histogram, the ``'Data'`` item contains
666the structure factors as an np.array in item `'RefList'`.
667The columns in that array are documented below.
668
669==========  ====================================================
670  index         explanation
671==========  ====================================================
672 0,1,2       h,k,l (float)
673 3           multiplicity
674 4           d-space, Angstrom
675 5           :math:`F_{obs}^2`
676 6           :math:`\sigma(F_{obs}^2)`
677 7           :math:`F_{calc}^2`
678 8           :math:`F_{obs}^2T`
679 9           :math:`F_{calc}^2T`
680 10          reflection phase, in degrees
681 11          intensity correction for reflection, this times
682             :math:`F_{obs}^2` or :math:`F_{calc}^2`
683             gives Iobs or Icalc
684==========  ====================================================
685
686Image Data Structure
687--------------------
688
689.. _Image_table:
690
691.. index::
692   image: Image data object description
693   image: Image object descriptions
694   
695Every 2-dimensional image is stored in the GSAS-II data tree
696with a top-level entry named beginning with the string "IMG ". The
697image data are directly associated with that tree item and there
698are a series of children to that item. The routines :func:`GSASII.GSASII.GetUsedHistogramsAndPhasesfromTree`
699and :func:`GSASIIstrIO.GetUsedHistogramsAndPhases` will
700load this information into a dictionary where the child tree name is
701used as a key, and the information in the main entry is assigned
702a key of ``Data``, as outlined below.
703
704.. tabularcolumns:: |l|l|p{4in}|
705
706======================  ======================  ====================================================
707  key                      sub-key              explanation
708======================  ======================  ====================================================
709Comments                       \                Text strings extracted from the original image data
710                                                header or a metafile. These cannot be changed by 
711                                                the user; it may be empty.                                               
712Image Controls              azmthOff            (float) The offset to be applied to an azimuthal
713                                                value. Accomodates
714                                                detector orientations other than with the detector
715                                                X-axis
716                                                horizontal.
717\                           background image    (list:str,float) The name of a tree item ("IMG ...") that is to be subtracted
718                                                during image integration multiplied by value. It must have the same size/shape as
719                                                the integrated image. NB: value < 0 for subtraction.
720\                           calibrant           (str) The material used for determining the position/orientation
721                                                of the image. The data is obtained from :func:`ImageCalibrants`
722                                                and UserCalibrants.py (supplied by user).
723\                           calibdmin           (float) The minimum d-spacing used during the last calibration run.
724\                           calibskip           (int) The number of expected diffraction lines skipped during the last
725                                                calibration run.
726\                           center              (list:floats) The [X,Y] point in detector coordinates (mm) where the direct beam
727                                                strikes the detector plane as determined by calibration. This point
728                                                does not have to be within the limits of the detector boundaries.
729\                           centerAzm           (bool) If True then the azimuth reported for the integrated slice
730                                                of the image is at the center line otherwise it is at the leading edge.
731\                           color               (str) The name of the colormap used to display the image. Default = 'Paired'.
732\                           cutoff              (float) The minimum value of I/Ib for a point selected in a diffraction ring for
733                                                calibration calculations. See pixLimit for details as how point is found.           
734\                           DetDepth            (float) Coefficient for penetration correction to distance; accounts for diffraction
735                                                ring offset at higher angles. Optionally determined by calibration.
736\                           DetDepthRef         (bool) If True then refine DetDepth during calibration/recalibration calculation.
737\                           distance            (float) The distance (mm) from sample to detector plane.
738\                           ellipses            (list:lists) Each object in ellipses is a list [center,phi,radii,color] where
739                                                center (list) is location (mm) of the ellipse center on the detector plane, phi is the
740                                                rotation of the ellipse minor axis from the x-axis, and radii are the minor & major
741                                                radii of the ellipse. If radii[0] is negative then parameters describe a hyperbola. Color
742                                                is the selected drawing color (one of 'b', 'g' ,'r') for the ellipse/hyperbola.
743\                           edgemin             (float) Not used;  parameter in EdgeFinder code.
744\                           fullIntegrate       (bool) If True then integrate over full 360 deg azimuthal range.
745\                           GonioAngles         (list:floats) The 'Omega','Chi','Phi' goniometer angles used for this image.
746                                                Required for texture calculations.
747\                           invert_x            (bool) If True display the image with the x-axis inverted.
748\                           invert_y            (bool) If True display the image with the y-axis inverted.
749\                           IOtth               (list:floats) The minimum and maximum 2-theta values to be used for integration.
750\                           LRazimuth           (list:floats) The minimum and maximum azimuth values to be used for integration.
751\                           Oblique             (list:float,bool) If True apply a detector absorption correction using the value to the
752                                                intensities obtained during integration.
753\                           outAzimuths         (int) The number of azimuth pie slices.
754\                           outChannels         (int) The number of 2-theta steps.
755\                           pixelSize           (list:ints) The X,Y dimensions (microns) of each pixel.
756\                           pixLimit            (int) A box in the image with 2*pixLimit+1 edges is searched to find the maximum.
757                                                This value (I) along with the minimum (Ib) in the box is reported by :func:`GSASIIimage.ImageLocalMax`
758                                                and subject to cutoff in :func:`GSASIIimage.makeRing`.
759                                                Locations are used to construct rings of points for calibration calcualtions.
760\                           PolaVal             (list:float,bool) If type='SASD' and if True, apply polarization correction to intensities from
761                                                integration using value.
762\                           rings               (list:lists) Each entry is [X,Y,dsp] where X & Y are lists of x,y coordinates around a
763                                                diffraction ring with the same d-spacing (dsp)
764\                           ring                (list) The x,y coordinates of the >5 points on an inner ring
765                                                selected by the user,
766\                           Range               (list) The minimum & maximum values of the image
767\                           rotation            (float) The angle between the x-axis and the vector about which the
768                                                detector is tilted. Constrained to -180 to 180 deg.     
769\                           SampleShape         (str) Currently only 'Cylinder'. Sample shape for Debye-Scherrer experiments; used for absorption
770                                                calculations.
771\                           SampleAbs           (list: float,bool) Value of absorption coefficient for Debye-Scherrer experimnents, flag if True
772                                                to cause correction to be applied.
773\                           setDefault          (bool) If True the use the image controls values for all new images to be read. (might be removed)
774\                           setRings            (bool) If True then display all the selected x,y ring positions (vida supra rings) used in the calibration.           
775\                           showLines           (bool) If True then isplay the integration limits to be used.
776\                           size                (list:int) The number of pixels on the image x & y axes
777\                           type                (str) One of 'PWDR', 'SASD' or 'REFL' for powder, small angle or reflectometry data, respectively.
778\                           tilt                (float) The angle the detector normal makes with the incident beam; range -90 to 90.
779\                           wavelength          (float) Tha radiation wavelength (Angstroms) as entered by the user (or someday obtained from the image header).
780                                               
781Masks                       Arcs                (list: lists) Each entry [2-theta,[azimuth[0],azimuth[1]],thickness] describes an arc mask
782                                                to be excluded from integration
783\                           Frames              (list:lists) Each entry describes the x,y points (3 or more - mm) that describe a frame outside
784                                                of which is excluded from recalibration and integration. Only one frame is allowed.
785\                           Points              (list:lists) Each entry [x,y,radius] (mm) describes an excluded spot on the image to be excluded
786                                                from integration.
787\                           Polygons            (list:lists) Each entry is a list of 3+ [x,y] points (mm) that describe a polygon on the image
788                                                to be excluded from integration.
789\                           Rings               (list: lists) Each entry [2-theta,thickness] describes a ring mask
790                                                to be excluded from integration.
791\                           Thresholds          (list:[tuple,list]) [(Imin,Imax),[Imin,Imax]] This gives lower and upper limits for points on the image to be included
792                                                in integrsation. The tuple is the image intensity limits and the list are those set by the user.   
793                                               
794Stress/Strain               Sample phi          (float) Sample rotation about vertical axis.
795\                           Sample z            (float) Sample translation from the calibration sample position (for Sample phi = 0)
796                                                These will be restricted by space group symmetry; result of strain fit refinement.
797\                           Type                (str) 'True' or 'Conventional': The strain model used for the calculation.
798\                           d-zero              (list:dict) Each item is for a diffraction ring on the image; all items are from the same phase
799                                                and are used to determine the strain tensor.
800                                                The dictionary items are:
801                                                'Dset': (float) True d-spacing for the diffraction ring; entered by the user.
802                                                'Dcalc': (float) Average calculated d-spacing determined from strain coeff.
803                                                'Emat': (list: float) The strain tensor elements e11, e12 & e22 (e21=e12, rest are 0)
804                                                'Esig': (list: float) Esds for Emat from fitting.
805                                                'pixLimit': (int) Search range to find highest point on ring for each data point
806                                                'cutoff': (float) I/Ib cutoff for searching.
807                                                'ImxyObs': (list: lists) [[X],[Y]] observed points to be used for strain calculations.
808                                                'ImtaObs': (list: lists) [[d],[azm]] transformed via detector calibration from ImxyObs.
809                                                'ImtaCalc': (list: lists [[d],[azm]] calculated d-spacing & azimuth from fit.
810                                               
811======================  ======================  ====================================================
812
813Parameter Dictionary
814-------------------------
815
816.. _parmDict_table:
817
818.. index::
819   single: Parameter dictionary
820
821The parameter dictionary contains all of the variable parameters for the refinement.
822The dictionary keys are the name of the parameter (<phase>:<hist>:<name>:<atom>).
823It is prepared in two ways. When loaded from the tree
824(in :meth:`GSASII.GSASII.MakeLSParmDict` and
825:meth:`GSASIIIO.ExportBaseclass.loadParmDict`),
826the values are lists with two elements: ``[value, refine flag]``
827
828When loaded from the GPX file (in
829:func:`GSASIIstrMain.Refine` and :func:`GSASIIstrMain.SeqRefine`), the value in the
830dict is the actual parameter value (usually a float, but sometimes a
831letter or string flag value (such as I or A for iso/anisotropic).
832
833
834*Classes and routines*
835----------------------
836
837'''
838import re
839import imp
840import random as ran
841import sys
842import GSASIIpath
843import GSASIImath as G2mth
844import numpy as np
845
846GSASIIpath.SetVersionNumber("$Revision: 1338 $")
847
848DefaultControls = {
849    'deriv type':'analytic Hessian',    #default controls
850    'min dM/M':0.0001,'shift factor':1.,'max cyc':3,'F**2':True,
851    'minF/sig':0,
852    'Author':'no name',
853    'FreeVar1':'Sample humidity (%)',
854    'FreeVar2':'Sample voltage (V)',
855    'FreeVar3':'Applied load (MN)',
856    }
857'''Values to be used as defaults for the initial contents of the ``Controls``
858data tree item.
859'''
860
861def MakeUniqueLabel(lbl,labellist):
862    '''Make sure that every a label is unique against a list by adding
863    digits at the end until it is not found in list.
864
865    :param str lbl: the input label
866    :param list labellist: the labels that have already been encountered
867    :returns: lbl if not found in labellist or lbl with ``_1-9`` (or
868      ``_10-99``, etc.) appended at the end
869    '''
870    lbl = lbl.strip()
871    if not lbl: # deal with a blank label
872        lbl = '_1'
873    if lbl not in labellist:
874        labellist.append(lbl)
875        return lbl
876    i = 1
877    prefix = lbl
878    if '_' in lbl:
879        prefix = lbl[:lbl.rfind('_')]
880        suffix = lbl[lbl.rfind('_')+1:]
881        try:
882            i = int(suffix)+1
883        except: # suffix could not be parsed
884            i = 1
885            prefix = lbl
886    while prefix+'_'+str(i) in labellist:
887        i += 1
888    else:
889        lbl = prefix+'_'+str(i)
890        labellist.append(lbl)
891    return lbl
892
893PhaseIdLookup = {}
894'''dict listing phase name and random Id keyed by sequential phase index as a str;
895best to access this using :func:`LookupPhaseName`
896'''
897PhaseRanIdLookup = {}
898'''dict listing phase sequential index keyed by phase random Id;
899best to access this using :func:`LookupPhaseId`
900'''
901HistIdLookup = {}
902'''dict listing histogram name and random Id, keyed by sequential histogram index as a str;
903best to access this using :func:`LookupHistName`
904'''
905HistRanIdLookup = {}
906'''dict listing histogram sequential index keyed by histogram random Id;
907best to access this using :func:`LookupHistId`
908'''
909AtomIdLookup = {}
910'''dict listing for each phase index as a str, the atom label and atom random Id,
911keyed by atom sequential index as a str;
912best to access this using :func:`LookupAtomLabel`
913'''
914AtomRanIdLookup = {}
915'''dict listing for each phase the atom sequential index keyed by atom random Id;
916best to access this using :func:`LookupAtomId`
917'''
918ShortPhaseNames = {}
919'''a dict containing a possibly shortened and when non-unique numbered
920version of the phase name. Keyed by the phase sequential index.
921'''
922ShortHistNames = {}
923'''a dict containing a possibly shortened and when non-unique numbered
924version of the histogram name. Keyed by the histogram sequential index.
925'''
926
927VarDesc = {}
928''' This dictionary lists descriptions for GSAS-II variables,
929as set in :func:`CompileVarDesc`. See that function for a description
930for how keys and values are written.
931'''
932
933reVarDesc = {}
934''' This dictionary lists descriptions for GSAS-II variables with
935the same values as :attr:`VarDesc` except that keys have been compiled as
936regular expressions. Initialized in :func:`CompileVarDesc`.
937'''
938
939def IndexAllIds(Histograms,Phases):
940    '''Scan through the used phases & histograms and create an index
941    to the random numbers of phases, histograms and atoms. While doing this,
942    confirm that assigned random numbers are unique -- just in case lightning
943    strikes twice in the same place.
944
945    Note: this code assumes that the atom random Id (ranId) is the last
946    element each atom record.
947
948    This is called in two places (only) :func:`GSASIIstrIO.GetUsedHistogramsAndPhases`
949    (which loads the histograms and phases from a GPX file) and
950    :meth:`GSASII.GSASII.GetUsedHistogramsAndPhasesfromTree`
951    (which loads the histograms and phases from the data tree.)
952
953    TODO: do we need a lookup for rigid body variables?
954    '''
955    # process phases and atoms
956    PhaseIdLookup.clear()
957    PhaseRanIdLookup.clear()   
958    AtomIdLookup.clear()
959    AtomRanIdLookup.clear()
960    ShortPhaseNames.clear()
961    for ph in Phases:
962        cx,ct,cs,cia = Phases[ph]['General']['AtomPtrs']
963        ranId = Phases[ph]['ranId'] 
964        while ranId in PhaseRanIdLookup:
965            # Found duplicate random Id! note and reassign
966            print ("\n\n*** Phase "+str(ph)+" has repeated ranId. Fixing.\n")
967            Phases[ph]['ranId'] = ranId = ran.randint(0,sys.maxint)
968        pId = str(Phases[ph]['pId'])
969        PhaseIdLookup[pId] = (ph,ranId)
970        PhaseRanIdLookup[ranId] = pId
971        shortname = ph[:10]
972        while shortname in ShortPhaseNames.values():
973            shortname = ph[:8] + ' ('+ pId + ')'
974        ShortPhaseNames[pId] = shortname
975        AtomIdLookup[pId] = {}
976        AtomRanIdLookup[pId] = {}
977        for iatm,at in enumerate(Phases[ph]['Atoms']):
978            ranId = at[-1]
979            while ranId in AtomRanIdLookup[pId]: # check for dups
980                print ("\n\n*** Phase "+str(ph)+" atom "+str(iatm)+" has repeated ranId. Fixing.\n")
981                at[-1] = ranId = ran.randint(0,sys.maxint)
982            AtomRanIdLookup[pId][ranId] = str(iatm)
983            if Phases[ph]['General']['Type'] == 'macromolecular':
984                label = '%s_%s_%s_%s'%(at[ct-1],at[ct-3],at[ct-4],at[ct-2])
985            else:
986                label = at[ct-1]
987            AtomIdLookup[pId][str(iatm)] = (label,ranId)
988    # process histograms
989    HistIdLookup.clear()
990    HistRanIdLookup.clear()
991    ShortHistNames.clear()
992    for hist in Histograms:
993        ranId = Histograms[hist]['ranId']
994        while ranId in HistRanIdLookup:
995            # Found duplicate random Id! note and reassign
996            print ("\n\n*** Histogram "+str(hist)+" has repeated ranId. Fixing.\n")
997            Histograms[hist]['ranId'] = ranId = ran.randint(0,sys.maxint)
998        hId = str(Histograms[hist]['hId'])
999        HistIdLookup[hId] = (hist,ranId)
1000        HistRanIdLookup[ranId] = hId
1001        shortname = hist[:15]
1002        while shortname in ShortHistNames.values():
1003            shortname = hist[:11] + ' ('+ hId + ')'
1004        ShortHistNames[hId] = shortname
1005
1006def LookupAtomId(pId,ranId):
1007    '''Get the atom number from a phase and atom random Id
1008
1009    :param int/str pId: the sequential number of the phase
1010    :param int ranId: the random Id assigned to an atom
1011
1012    :returns: the index number of the atom (str)
1013    '''
1014    if not AtomRanIdLookup:
1015        raise Exception,'Error: LookupAtomId called before IndexAllIds was run'
1016    if pId is None or pId == '':
1017        raise KeyError,'Error: phase is invalid (None or blank)'
1018    pId = str(pId)
1019    if pId not in AtomRanIdLookup:
1020        raise KeyError,'Error: LookupAtomId does not have phase '+pId
1021    if ranId not in AtomRanIdLookup[pId]:
1022        raise KeyError,'Error: LookupAtomId, ranId '+str(ranId)+' not in AtomRanIdLookup['+pId+']'
1023    return AtomRanIdLookup[pId][ranId]
1024
1025def LookupAtomLabel(pId,index):
1026    '''Get the atom label from a phase and atom index number
1027
1028    :param int/str pId: the sequential number of the phase
1029    :param int index: the index of the atom in the list of atoms
1030
1031    :returns: the label for the atom (str) and the random Id of the atom (int)
1032    '''
1033    if not AtomIdLookup:
1034        raise Exception,'Error: LookupAtomLabel called before IndexAllIds was run'
1035    if pId is None or pId == '':
1036        raise KeyError,'Error: phase is invalid (None or blank)'
1037    pId = str(pId)
1038    if pId not in AtomIdLookup:
1039        raise KeyError,'Error: LookupAtomLabel does not have phase '+pId
1040    if index not in AtomIdLookup[pId]:
1041        raise KeyError,'Error: LookupAtomLabel, ranId '+str(index)+' not in AtomRanIdLookup['+pId+']'
1042    return AtomIdLookup[pId][index]
1043
1044def LookupPhaseId(ranId):
1045    '''Get the phase number and name from a phase random Id
1046
1047    :param int ranId: the random Id assigned to a phase
1048    :returns: the sequential Id (pId) number for the phase (str)
1049    '''
1050    if not PhaseRanIdLookup:
1051        raise Exception,'Error: LookupPhaseId called before IndexAllIds was run'
1052    if ranId not in PhaseRanIdLookup:
1053        raise KeyError,'Error: LookupPhaseId does not have ranId '+str(ranId)
1054    return PhaseRanIdLookup[ranId]
1055
1056def LookupPhaseName(pId):
1057    '''Get the phase number and name from a phase Id
1058
1059    :param int/str pId: the sequential assigned to a phase
1060    :returns:  (phase,ranId) where phase is the name of the phase (str)
1061      and ranId is the random # id for the phase (int)
1062    '''
1063    if not PhaseIdLookup:
1064        raise Exception,'Error: LookupPhaseName called before IndexAllIds was run'
1065    if pId is None or pId == '':
1066        raise KeyError,'Error: phase is invalid (None or blank)'
1067    pId = str(pId)
1068    if pId not in PhaseIdLookup:
1069        raise KeyError,'Error: LookupPhaseName does not have index '+pId
1070    return PhaseIdLookup[pId]
1071
1072def LookupHistId(ranId):
1073    '''Get the histogram number and name from a histogram random Id
1074
1075    :param int ranId: the random Id assigned to a histogram
1076    :returns: the sequential Id (hId) number for the histogram (str)
1077    '''
1078    if not HistRanIdLookup:
1079        raise Exception,'Error: LookupHistId called before IndexAllIds was run'
1080    if ranId not in HistRanIdLookup:
1081        raise KeyError,'Error: LookupHistId does not have ranId '+str(ranId)
1082    return HistRanIdLookup[ranId]
1083
1084def LookupHistName(hId):
1085    '''Get the histogram number and name from a histogram Id
1086
1087    :param int/str hId: the sequential assigned to a histogram
1088    :returns:  (hist,ranId) where hist is the name of the histogram (str)
1089      and ranId is the random # id for the histogram (int)
1090    '''
1091    if not HistIdLookup:
1092        raise Exception,'Error: LookupHistName called before IndexAllIds was run'
1093    if hId is None or hId == '':
1094        raise KeyError,'Error: histogram is invalid (None or blank)'
1095    hId = str(hId)
1096    if hId not in HistIdLookup:
1097        raise KeyError,'Error: LookupHistName does not have index '+hId
1098    return HistIdLookup[hId]
1099
1100def fmtVarDescr(varname):
1101    '''Return a string with a more complete description for a GSAS-II variable
1102
1103    :param str varname: A full G2 variable name with 2 or 3 or 4
1104       colons (<p>:<h>:name[:<a>] or <p>::RBname:<r>:<t>])
1105       
1106    :returns: a string with the description
1107    '''
1108    s,l = VarDescr(varname)
1109    return s+": "+l
1110
1111def VarDescr(varname):
1112    '''Return two strings with a more complete description for a GSAS-II variable
1113
1114    :param str name: A full G2 variable name with 2 or 3 or 4
1115       colons (<p>:<h>:name[:<a>] or <p>::RBname:<r>:<t>])
1116       
1117    :returns: (loc,meaning) where loc describes what item the variable is mapped
1118      (phase, histogram, etc.) and meaning describes what the variable does.
1119    '''
1120   
1121    # special handling for parameter names without a colons
1122    # for now, assume self-defining
1123    if varname.find(':') == -1:
1124        return "Global",varname
1125       
1126    l = getVarDescr(varname)
1127    if not l:
1128        return ("invalid variable name ("+str(varname)+")!"),""
1129#        return "invalid variable name!",""
1130
1131    if not l[-1]:
1132        l[-1] = "(variable needs a definition!)"
1133
1134    if len(l) == 3:         #SASD variable name!
1135        s = 'component:'+l[1]
1136        return s,l[-1]
1137    s = ""
1138    if l[0] is not None and l[1] is not None: # HAP: keep short
1139        if l[0] == '*':
1140            lbl = 'all'
1141        else:
1142            lbl = ShortPhaseNames.get(l[0],'? #'+str(l[0]))
1143        if l[1] == '*':
1144            hlbl = 'all'
1145        else:
1146            hlbl = ShortHistNames.get(l[1],'? #'+str(l[1]))
1147        if hlbl[:4] == 'HKLF':
1148            hlbl = 'Xtl='+hlbl[5:]
1149        elif hlbl[:4] == 'PWDR':
1150            hlbl = 'Pwd='+hlbl[5:]
1151        else:
1152            hlbl = 'Hist='+hlbl
1153        s = "Ph="+str(lbl)+" * "+str(hlbl)
1154    elif l[2] == 'Back': # background parameters are "special", alas
1155        s = 'Hist='+ShortHistNames.get(l[1],'? #'+str(l[1]))
1156        l[-1] += ' #'+str(l[3])
1157    elif l[4] is not None: # rigid body parameter
1158        lbl = ShortPhaseNames.get(l[0],'phase?')
1159        s = "Res #"+str(l[3])+" body #"+str(l[4])+" in "+str(lbl)
1160    elif l[3] is not None: # atom parameter,
1161        lbl = ShortPhaseNames.get(l[0],'phase?')
1162        try:
1163            albl = LookupAtomLabel(l[0],l[3])[0]
1164        except KeyError:
1165            albl = 'Atom?'
1166        s = "Atom "+str(albl)+" in "+str(lbl)
1167    elif l[0] == '*':
1168        s = "All phases "
1169    elif l[0] is not None:
1170        lbl = ShortPhaseNames.get(l[0],'phase?')
1171        s = "Phase "+str(lbl)
1172    elif l[1] == '*':
1173        s = 'All hists'
1174    elif l[1] is not None:
1175        hlbl = ShortHistNames.get(l[1],'? #'+str(l[1]))
1176        if hlbl[:4] == 'HKLF':
1177            hlbl = 'Xtl='+hlbl[5:]
1178        elif hlbl[:4] == 'PWDR':
1179            hlbl = 'Pwd='+hlbl[5:]
1180        else:
1181            hlbl = 'Hist='+hlbl
1182        s = str(hlbl)
1183    if not s:
1184        s = 'Global'
1185    return s,l[-1]
1186
1187def getVarDescr(varname):
1188    '''Return a short description for a GSAS-II variable
1189
1190    :param str name: A full G2 variable name with 2 or 3 or 4
1191       colons (<p>:<h>:name[:<a1>][:<a2>])
1192     
1193    :returns: a six element list as [`p`,`h`,`name`,`a1`,`a2`,`description`],
1194      where `p`, `h`, `a1`, `a2` are str values or `None`, for the phase number,
1195      the histogram number and the atom number; `name` will always be
1196      a str; and `description` is str or `None`.
1197      If the variable name is incorrectly formed (for example, wrong
1198      number of colons), `None` is returned instead of a list.
1199    '''
1200    l = varname.split(':')
1201    if len(l) == 2:     #SASD parameter name
1202        return varname,l[0],getDescr(l[1])
1203    if len(l) == 3:
1204        l += [None,None]
1205    elif len(l) == 4:
1206        l += [None]
1207    elif len(l) != 5:
1208        return None
1209    for i in (0,1,3,4):
1210        if l[i] == "":
1211            l[i] = None
1212    l += [getDescr(l[2])]
1213    return l
1214   
1215def CompileVarDesc():
1216    '''Set the values in the variable description lookup table (:attr:`VarDesc`)
1217    into :attr:`reVarDesc`. This is called in :func:`getDescr` so the initialization
1218    is always done before use.
1219
1220    Note that keys may contain regular expressions, where '[xyz]'
1221    matches 'x' 'y' or 'z' (equivalently '[x-z]' describes this as range of values).
1222    '.*' matches any string. For example::
1223
1224    'AUiso':'Atomic isotropic displacement parameter',
1225
1226    will match variable ``'p::AUiso:a'``.
1227    If parentheses are used in the key, the contents of those parentheses can be
1228    used in the value, such as::
1229
1230    'AU([123][123])':'Atomic anisotropic displacement parameter U\\1',
1231
1232    will match ``AU11``, ``AU23``,.. and `U11`, `U23` etc will be displayed
1233    in the value when used.
1234   
1235    '''
1236    if reVarDesc: return # already done
1237    for key,value in {
1238        # derived or other sequential vars
1239        '([abc])$' : 'Lattice parameter, \\1, from Ai and Djk', # N.B. '$' prevents match if any characters follow
1240        u'\u03B1' : u'Lattice parameter, \u03B1, from Ai and Djk',
1241        u'\u03B2' : u'Lattice parameter, \u03B2, from Ai and Djk',
1242        u'\u03B3' : u'Lattice parameter, \u03B3, from Ai and Djk',
1243        # ambiguous, alas:
1244        'Scale' : 'Phase or Histogram scale factor',
1245        # Phase vars (p::<var>)
1246        'A([0-5])' : 'Reciprocal metric tensor component \\1',
1247        'Vol' : 'Unit cell volume',
1248        # Atom vars (p::<var>:a)
1249        'dA([xyz])' : 'change to atomic coordinate, \\1',
1250        'A([xyz])$' : '\\1 fractional atomic coordinate',
1251        'AUiso':'Atomic isotropic displacement parameter',
1252        'AU([123][123])':'Atomic anisotropic displacement parameter U\\1',
1253        'Afrac': 'Atomic occupancy parameter',
1254        # Hist & Phase (HAP) vars (p:h:<var>)
1255        'Back': 'Background term',
1256        'BkPkint;(.*)':'Background peak #\\1 intensity',
1257        'BkPkpos;(.*)':'Background peak #\\1 position',
1258        'BkPksig;(.*)':'Background peak #\\1 Gaussian width',
1259        'BkPkgam;(.*)':'Background peak #\\1 Cauchy width',
1260        'Bab([AU])': 'Babinet solvent scattering coef. \\1',
1261        'D([123][123])' : 'Anisotropic strain coef. \\1',
1262        'Extinction' : 'Extinction coef.',
1263        'MD' : 'March-Dollase coef.',
1264        'Mustrain;.*' : 'Microstrain coef.',
1265        'Size;.*' : 'Crystallite size value',
1266        'eA' : '?',
1267        #Histogram vars (:h:<var>)
1268        'Absorption' : 'Absorption coef.',
1269        'Displace([XY])' : 'Debye-Scherrer sample displacement \\1',
1270        'Lam' : 'Wavelength',
1271        'Polariz\.' : 'Polarization correction',
1272        'SH/L' : 'FCJ peak asymmetry correction',
1273        '([UVW])' : 'Gaussian instrument broadening \\1',
1274        '([XY])' : 'Cauchy instrument broadening \\1',
1275        'Zero' : 'Debye-Scherrer zero correction',
1276        'nDebye' : 'Debye model background corr. terms',
1277        'nPeaks' : 'Fixed peak background corr. terms',
1278        'RBV.*' : 'Vector rigid body parameter',
1279        'RBR.*' : 'Residue rigid body parameter',
1280        'RBRO([aijk])' : 'Residue rigid body orientation parameter',
1281        'RBRP([xyz])' : 'Residue rigid body position parameter',
1282        'RBRTr;.*' : 'Residue rigid body torsion parameter',
1283        'RBR([TLS])([123AB][123AB])' : 'Residue rigid body group disp. param.',
1284        # Global vars (::<var>)
1285        'constr([0-9]*)' : 'Parameter from constraint',
1286        # SASD vars (l:<var>;l = component)
1287        'Aspect ratio' : 'Particle aspect ratio',
1288        'Length' : 'Cylinder length',
1289        'Diameter' : 'Cylinder/disk diameter',
1290        'Thickness' : 'Disk thickness',
1291        'Dist' : 'Interparticle distance',
1292        'VolFr' : 'Dense scatterer volume fraction',
1293        'epis' : 'Sticky sphere epsilon',
1294        'Sticky' : 'Stickyness',
1295        'Depth' : 'Well depth',
1296        'Width' : 'Well width',
1297        'Volume' : 'Particle volume',
1298        'Radius' : 'Sphere/cylinder/disk radius',
1299        'Mean' : 'Particle mean radius',
1300        'StdDev' : 'Standard deviation in Mean',
1301        'G': 'Guinier prefactor',
1302        'Rg': 'Guinier radius of gyration',
1303        'B': 'Porod prefactor',
1304        'P': 'Porod power',
1305        'Cutoff': 'Porod cutoff',
1306        'PkInt': 'Bragg peak intensity',
1307        'PkPos': 'Bragg peak position',
1308        'PkSig': 'Bragg peak sigma',
1309        'PkGam': 'Bragg peak gamma',   
1310        }.items():
1311        VarDesc[key] = value
1312        reVarDesc[re.compile(key)] = value
1313
1314def getDescr(name):
1315    '''Return a short description for a GSAS-II variable
1316
1317    :param str name: The descriptive part of the variable name without colons (:)
1318     
1319    :returns: a short description or None if not found
1320    '''
1321
1322    CompileVarDesc() # compile the regular expressions, if needed
1323    for key in reVarDesc:
1324        m = key.match(name)
1325        if m:
1326            return m.expand(reVarDesc[key])
1327    return None
1328
1329def GenWildCard(varlist):
1330    '''Generate wildcard versions of G2 variables. These introduce '*'
1331    for a phase, histogram or atom number (but only for one of these
1332    fields) but only when there is more than one matching variable in the
1333    input variable list. So if the input is this::
1334   
1335      varlist = ['0::AUiso:0', '0::AUiso:1', '1::AUiso:0']
1336
1337    then the output will be this::
1338   
1339       wildList = ['*::AUiso:0', '0::AUiso:*']
1340
1341    :param list varlist: an input list of GSAS-II variable names
1342      (such as 0::AUiso:0)
1343
1344    :returns: wildList, the generated list of wild card variable names.
1345    '''
1346    wild = []
1347    for i in (0,1,3):
1348        currentL = varlist[:]
1349        while currentL:
1350            item1 = currentL.pop(0)
1351            i1splt = item1.split(':')
1352            if i >= len(i1splt): continue
1353            if i1splt[i]:
1354                nextL = []
1355                i1splt[i] = '[0-9]+'
1356                rexp = re.compile(':'.join(i1splt))
1357                matchlist = [item1]
1358                for nxtitem in currentL:
1359                    if rexp.match(nxtitem):
1360                        matchlist += [nxtitem]
1361                    else:
1362                        nextL.append(nxtitem)
1363                if len(matchlist) > 1:
1364                    i1splt[i] = '*'
1365                    wild.append(':'.join(i1splt))
1366                currentL = nextL
1367    return wild
1368
1369def LookupWildCard(varname,varlist):
1370    '''returns a list of variable names from list varname
1371    that match wildcard name in varname
1372   
1373    :param str varname: a G2 variable name containing a wildcard
1374      (such as \*::var)
1375    :param list varlist: the list of all variable names used in
1376      the current project
1377    :returns: a list of matching GSAS-II variables (may be empty) 
1378    '''
1379    rexp = re.compile(varname.replace('*','[0-9]+'))
1380    return sorted([var for var in varlist if rexp.match(var)])
1381
1382
1383def _lookup(dic,key):
1384    '''Lookup a key in a dictionary, where None returns an empty string
1385    but an unmatched key returns a question mark. Used in :class:`G2VarObj`
1386    '''
1387    if key is None:
1388        return ""
1389    elif key == "*":
1390        return "*"
1391    else:
1392        return dic.get(key,'?')
1393
1394class G2VarObj(object):
1395    '''Defines a GSAS-II variable either using the phase/atom/histogram
1396    unique Id numbers or using a character string that specifies
1397    variables by phase/atom/histogram number (which can change).
1398    Note that :func:`LoadID` should be used to (re)load the current Ids
1399    before creating or later using the G2VarObj object.
1400
1401    This can store rigid body variables, but does not translate the residue # and
1402    body # to/from random Ids
1403
1404    A :class:`G2VarObj` object can be created with a single parameter:
1405   
1406    :param str/tuple varname: a single value can be used to create a :class:`G2VarObj`
1407      object. If a string, it must be of form "p:h:var" or "p:h:var:a", where
1408
1409     * p is the phase number (which may be left blank or may be '*' to indicate all phases);
1410     * h is the histogram number (which may be left blank or may be '*' to indicate all histograms);
1411     * a is the atom number (which may be left blank in which case the third colon is omitted).
1412       The atom number can be specified as '*' if a phase number is specified (not as '*').
1413       For rigid body variables, specify a will be a string of form "residue:body#"
1414
1415      Alternately a single tuple of form (Phase,Histogram,VarName,AtomID) can be used, where
1416      Phase, Histogram, and AtomID are None or are ranId values (or one can be '*')
1417      and VarName is a string. Note that if Phase is '*' then the AtomID is an atom number.
1418      For a rigid body variables, AtomID is a string of form "residue:body#".
1419
1420    If four positional arguments are supplied, they are:
1421
1422    :param str/int phasenum: The number for the phase (or None or '*')
1423    :param str/int histnum: The number for the histogram (or None or '*')
1424    :param str varname: a single value can be used to create a :class:`G2VarObj`
1425    :param str/int atomnum: The number for the atom (or None or '*')
1426   
1427    '''
1428    IDdict = {}
1429    IDdict['phases'] = {}
1430    IDdict['hists'] = {}
1431    IDdict['atoms'] = {}
1432    def __init__(self,*args):
1433        self.phase = None
1434        self.histogram = None
1435        self.name = ''
1436        self.atom = None
1437        if len(args) == 1 and (type(args[0]) is list or type(args[0]) is tuple) and len(args[0]) == 4:
1438            # single arg with 4 values
1439            self.phase,self.histogram,self.name,self.atom = args[0]
1440        elif len(args) == 1 and ':' in args[0]:
1441            #parse a string
1442            lst = args[0].split(':')
1443            if lst[0] == '*':
1444                self.phase = '*'
1445                if len(lst) > 3:
1446                    self.atom = lst[3]
1447                self.histogram = HistIdLookup.get(lst[1],[None,None])[1]
1448            elif lst[1] == '*':           
1449                self.histogram = '*'
1450                self.phase = PhaseIdLookup.get(lst[0],[None,None])[1]
1451            else:
1452                self.histogram = HistIdLookup.get(lst[1],[None,None])[1]
1453                self.phase = PhaseIdLookup.get(lst[0],[None,None])[1]
1454                if len(lst) == 4:
1455                    if lst[3] == '*':
1456                        self.atom = '*'
1457                    else:
1458                        self.atom = AtomIdLookup[lst[0]].get(lst[3],[None,None])[1]
1459                elif len(lst) == 5:
1460                    self.atom = lst[3]+":"+lst[4]
1461                elif len(lst) == 3:
1462                    pass
1463                else:
1464                    raise Exception,"Too many colons in var name "+str(args[0])
1465            self.name = lst[2]
1466        elif len(args) == 4:
1467            if args[0] == '*':
1468                self.phase = '*'
1469                self.atom = args[3]
1470            else:
1471                self.phase = PhaseIdLookup.get(str(args[0]),[None,None])[1]
1472                if args[3] == '*':
1473                    self.atom = '*'
1474                elif args[0] is not None:
1475                    self.atom = AtomIdLookup[args[0]].get(str(args[3]),[None,None])[1]
1476            if args[1] == '*':
1477                self.histogram = '*'
1478            else:
1479                self.histogram = HistIdLookup.get(str(args[1]),[None,None])[1]
1480            self.name = args[2]
1481        else:
1482            raise Exception,"Incorrectly called GSAS-II parameter name"
1483
1484        #print "DEBUG: created ",self.phase,self.histogram,self.name,self.atom
1485
1486    def __str__(self):
1487        return self.varname()
1488
1489    def varname(self):
1490        '''Formats the GSAS-II variable name as a "traditional" GSAS-II variable
1491        string (p:h:<var>:a) or (p:h:<var>)
1492
1493        :returns: the variable name as a str
1494        '''
1495        a = ""
1496        if self.phase == "*":
1497            ph = "*"
1498            if self.atom:
1499                a = ":" + str(self.atom)
1500        else:
1501            ph = _lookup(PhaseRanIdLookup,self.phase)
1502            if self.atom == '*':
1503                a = ':*'
1504            elif self.atom:
1505                if ":" in str(self.atom):
1506                    a = ":" + str(self.atom)
1507                elif ph in AtomRanIdLookup:
1508                    a = ":" + AtomRanIdLookup[ph].get(self.atom,'?')
1509                else:
1510                    a = ":?"
1511        if self.histogram == "*":
1512            hist = "*"
1513        else:
1514            hist = _lookup(HistRanIdLookup,self.histogram)
1515        s = (ph + ":" + hist + ":" + str(self.name)) + a
1516        return s
1517   
1518    def __repr__(self):
1519        '''Return the detailed contents of the object
1520        '''
1521        s = "<"
1522        if self.phase == '*':
1523            s += "Phases: all; "
1524            if self.atom is not None:
1525                if ":" in str(self.atom):
1526                    s += "Rigid body" + str(self.atom) + "; "
1527                else:
1528                    s += "Atom #" + str(self.atom) + "; "
1529        elif self.phase is not None:
1530            ph =  _lookup(PhaseRanIdLookup,self.phase)
1531            s += "Phase: rId=" + str(self.phase) + " (#"+ ph + "); "
1532            if self.atom == '*':
1533                s += "Atoms: all; "
1534            elif ":" in self(self.atom):
1535                s += "Rigid body" + str(self.atom) + "; "
1536            elif self.atom is not None:
1537                s += "Atom rId=" + str(self.atom)
1538                if ph in AtomRanIdLookup:
1539                    s += " (#" + AtomRanIdLookup[ph].get(self.atom,'?') + "); "
1540                else:
1541                    s += " (#? -- not found!); "
1542        if self.histogram == '*':
1543            s += "Histograms: all; "
1544        elif self.histogram is not None:
1545            hist = _lookup(HistRanIdLookup,self.histogram)
1546            s += "Histogram: rId=" + str(self.histogram) + " (#"+ hist + "); "
1547        s += 'Variable name="' + str(self.name) + '">'
1548        return s+" ("+self.varname()+")"
1549
1550    def __eq__(self, other):
1551        if type(other) is type(self):
1552            return (self.phase == other.phase and
1553                    self.histogram == other.histogram and
1554                    self.name == other.name and
1555                    self.atom == other.atom)
1556        return False
1557
1558    def _show(self):
1559        'For testing, shows the current lookup table'
1560        print 'phases', self.IDdict['phases']
1561        print 'hists', self.IDdict['hists']
1562        print 'atomDict', self.IDdict['atoms']
1563
1564#==========================================================================
1565# shortcut routines
1566exp = np.exp
1567sind = lambda x: np.sin(x*np.pi/180.)
1568tand = lambda x: np.tan(x*np.pi/180.)
1569cosd = lambda x: np.cos(x*np.pi/180.)
1570sind = sin = s = lambda x: np.sin(x*np.pi/180.)
1571cosd = cos = c = lambda x: np.cos(x*np.pi/180.)
1572tand = tan = t = lambda x: np.tan(x*np.pi/180.)
1573sqrt = sq = lambda x: np.sqrt(x)
1574pi = lambda: np.pi
1575class ExpressionObj(object):
1576    '''Defines an object with a user-defined expression, to be used for
1577    secondary fits or restraints. Object is created null, but is changed
1578    using :meth:`LoadExpression`. This contains only the minimum
1579    information that needs to be stored to save and load the expression
1580    and how it is mapped to GSAS-II variables.
1581    '''
1582    def __init__(self):
1583        self.expression = ''
1584        'The expression as a text string'
1585        self.assgnVars = {}
1586        '''A dict where keys are label names in the expression mapping to a GSAS-II
1587        variable. The value is a list with a G2 variable name and derivative step size.
1588        Note that the G2 variable name may contain a wild-card and correspond to
1589        multiple values.
1590        '''
1591        self.freeVars = {}
1592        '''A dict where keys are label names in the expression mapping to a free
1593        parameter. The value is a list with:
1594
1595         * a name assigned to the parameter
1596         * a value for to the parameter and
1597         * a flag to determine if the variable is refined.
1598        ''' 
1599        self.depVar = None
1600
1601        self.lastError = ('','')
1602        '''Shows last encountered error in processing expression
1603        (list of 1-3 str values)'''
1604
1605    def LoadExpression(self,expr,exprVarLst,varSelect,varName,varValue,varRefflag):
1606        '''Load the expression and associated settings into the object. Raises
1607        an exception if the expression is not parsed, if not all functions
1608        are defined or if not all needed parameter labels in the expression
1609        are defined.
1610
1611        This will not test if the variable referenced in these definitions
1612        are actually in the parameter dictionary. This is checked when the
1613        computation for the expression is done in :meth:`SetupCalc`.
1614       
1615        :param str expr: the expression
1616        :param list exprVarLst: parameter labels found in the expression
1617        :param dict varSelect: this will be 0 for Free parameters
1618          and non-zero for expression labels linked to G2 variables.
1619        :param dict varName: Defines a name (str) associated with each free parameter
1620        :param dict varValue: Defines a value (float) associated with each free parameter
1621        :param dict varRefflag: Defines a refinement flag (bool)
1622          associated with each free parameter
1623        '''
1624        self.expression = expr
1625        self.compiledExpr = None
1626        self.freeVars = {}
1627        self.assgnVars = {}
1628        for v in exprVarLst:
1629            if varSelect[v] == 0:
1630                self.freeVars[v] = [
1631                    varName.get(v),
1632                    varValue.get(v),
1633                    varRefflag.get(v),
1634                    ]
1635            else:
1636                self.assgnVars[v] = varName[v]
1637        self.CheckVars()
1638
1639    def EditExpression(self,exprVarLst,varSelect,varName,varValue,varRefflag):
1640        '''Load the expression and associated settings from the object into
1641        arrays used for editing.
1642
1643        :param list exprVarLst: parameter labels found in the expression
1644        :param dict varSelect: this will be 0 for Free parameters
1645          and non-zero for expression labels linked to G2 variables.
1646        :param dict varName: Defines a name (str) associated with each free parameter
1647        :param dict varValue: Defines a value (float) associated with each free parameter
1648        :param dict varRefflag: Defines a refinement flag (bool)
1649          associated with each free parameter
1650
1651        :returns: the expression as a str
1652        '''
1653        for v in self.freeVars:
1654            varSelect[v] = 0
1655            varName[v] = self.freeVars[v][0]
1656            varValue[v] = self.freeVars[v][1]
1657            varRefflag[v] = self.freeVars[v][2]
1658        for v in self.assgnVars:
1659            varSelect[v] = 1
1660            varName[v] = self.assgnVars[v]
1661        return self.expression
1662
1663    def GetVaried(self):
1664        'Returns the names of the free parameters that will be refined'
1665        return ["::"+self.freeVars[v][0] for v in self.freeVars if self.freeVars[v][2]]
1666
1667    def GetVariedVarVal(self):
1668        'Returns the names and values of the free parameters that will be refined'
1669        return [("::"+self.freeVars[v][0],self.freeVars[v][1]) for v in self.freeVars if self.freeVars[v][2]]
1670
1671    def UpdateVariedVars(self,varyList,values):
1672        'Updates values for the free parameters (after a refinement); only updates refined vars'
1673        for v in self.freeVars:
1674            if not self.freeVars[v][2]: continue
1675            if "::"+self.freeVars[v][0] not in varyList: continue
1676            indx = varyList.index("::"+self.freeVars[v][0])
1677            self.freeVars[v][1] = values[indx]
1678
1679    def GetIndependentVars(self):
1680        'Returns the names of the required independent parameters used in expression'
1681        return [self.assgnVars[v] for v in self.assgnVars]
1682
1683    def CheckVars(self):
1684        '''Check that the expression can be parsed, all functions are
1685        defined and that input loaded into the object is internally
1686        consistent. If not an Exception is raised.
1687
1688        :returns: a dict with references to packages needed to
1689          find functions referenced in the expression.
1690        '''
1691        ret = self.ParseExpression(self.expression)
1692        if not ret:
1693            raise Exception("Expression parse error")
1694        exprLblList,fxnpkgdict = ret
1695        # check each var used in expression is defined
1696        defined = self.assgnVars.keys() + self.freeVars.keys()
1697        notfound = []
1698        for var in exprLblList:
1699            if var not in defined:
1700                notfound.append(var)
1701        if notfound:
1702            msg = 'Not all variables defined'
1703            msg1 = 'The following variables were not defined: '
1704            msg2 = ''
1705            for var in notfound:
1706                if msg: msg += ', '
1707                msg += var
1708            self.lastError = (msg1,'  '+msg2)
1709            raise Exception(msg)
1710        return fxnpkgdict
1711
1712    def ParseExpression(self,expr):
1713        '''Parse an expression and return a dict of called functions and
1714        the variables used in the expression. Returns None in case an error
1715        is encountered. If packages are referenced in functions, they are loaded
1716        and the functions are looked up into the modules global
1717        workspace.
1718       
1719        Note that no changes are made to the object other than
1720        saving an error message, so that this can be used for testing prior
1721        to the save.
1722
1723        :returns: a list of used variables
1724        '''
1725        self.lastError = ('','')
1726        import ast
1727        def FindFunction(f):
1728            '''Find the object corresponding to function f
1729            :param str f: a function name such as 'numpy.exp'
1730            :returns: (pkgdict,pkgobj) where pkgdict contains a dict
1731              that defines the package location(s) and where pkgobj
1732              defines the object associated with the function.
1733              If the function is not found, pkgobj is None.
1734            '''
1735            df = f.split('.')
1736            pkgdict = {}
1737            # no listed package, try in current namespace
1738            if len(df) == 1: 
1739                try:
1740                    fxnobj = eval(f)
1741                    return pkgdict,fxnobj
1742                except (AttributeError, NameError):
1743                    return None,None
1744            else:
1745                try:
1746                    fxnobj = eval(f)
1747                    pkgdict[df[0]] = eval(df[0])
1748                    return pkgdict,fxnobj
1749                except (AttributeError, NameError):
1750                    pass
1751            # includes a package, lets try to load the packages
1752            pkgname = ''
1753            path = sys.path
1754            for pkg in f.split('.')[:-1]: # if needed, descend down the tree
1755                if pkgname:
1756                    pkgname += '.' + pkg
1757                else:
1758                    pkgname = pkg
1759                fp = None
1760                try:
1761                    fp, fppath,desc = imp.find_module(pkg,path)
1762                    pkgobj = imp.load_module(pkg,fp,fppath,desc)
1763                    pkgdict[pkgname] = pkgobj
1764                    path = [fppath]
1765                except Exception as msg:
1766                    print('load of '+pkgname+' failed with error='+str(msg))
1767                    return {},None
1768                finally:
1769                    if fp: fp.close()
1770                try:
1771                    #print 'before',pkgdict.keys()
1772                    fxnobj = eval(f,globals(),pkgdict)
1773                    #print 'after 1',pkgdict.keys()
1774                    #fxnobj = eval(f,pkgdict)
1775                    #print 'after 2',pkgdict.keys()
1776                    return pkgdict,fxnobj
1777                except:
1778                    continue
1779            return None # not found
1780        def ASTtransverse(node,fxn=False):
1781            '''Transverse a AST-parsed expresson, compiling a list of variables
1782            referenced in the expression. This routine is used recursively.
1783
1784            :returns: varlist,fxnlist where
1785              varlist is a list of referenced variable names and
1786              fxnlist is a list of used functions
1787            '''
1788            varlist = []
1789            fxnlist = []
1790            if isinstance(node, list):
1791                for b in node:
1792                    v,f = ASTtransverse(b,fxn)
1793                    varlist += v
1794                    fxnlist += f
1795            elif isinstance(node, ast.AST):
1796                for a, b in ast.iter_fields(node):
1797                    if isinstance(b, ast.AST):
1798                        if a == 'func': 
1799                            fxnlist += ['.'.join(ASTtransverse(b,True)[0])]
1800                            continue
1801                        v,f = ASTtransverse(b,fxn)
1802                        varlist += v
1803                        fxnlist += f
1804                    elif isinstance(b, list):
1805                        v,f = ASTtransverse(b,fxn)
1806                        varlist += v
1807                        fxnlist += f
1808                    elif node.__class__.__name__ == "Name":
1809                        varlist += [b]
1810                    elif fxn and node.__class__.__name__ == "Attribute":
1811                        varlist += [b]
1812            return varlist,fxnlist
1813        try:
1814            exprast = ast.parse(expr)
1815        except SyntaxError as err:
1816            s = ''
1817            import traceback
1818            for i in traceback.format_exc().splitlines()[-3:-1]:
1819                if s: s += "\n"
1820                s += str(i)
1821            self.lastError = ("Error parsing expression:",s)
1822            return
1823        # find the variables & functions
1824        v,f = ASTtransverse(exprast)
1825        varlist = sorted(list(set(v)))
1826        fxnlist = list(set(f))
1827        pkgdict = {}
1828        # check the functions are defined
1829        for fxn in fxnlist:
1830            fxndict,fxnobj = FindFunction(fxn)
1831            if not fxnobj:
1832                self.lastError = ("Error: Invalid function",fxn,
1833                                  "is not defined")
1834                return
1835            if not hasattr(fxnobj,'__call__'):
1836                self.lastError = ("Error: Not a function.",fxn,
1837                                  "cannot be called as a function")
1838                return
1839            pkgdict.update(fxndict)
1840        return varlist,pkgdict
1841
1842    def GetDepVar(self):
1843        'return the dependent variable, or None'
1844        return self.depVar
1845
1846    def SetDepVar(self,var):
1847        'Set the dependent variable, if used'
1848        self.depVar = var
1849#==========================================================================
1850class ExpressionCalcObj(object):
1851    '''An object used to evaluate an expression from a :class:`ExpressionObj`
1852    object.
1853   
1854    :param ExpressionObj exprObj: a :class:`~ExpressionObj` expression object with
1855      an expression string and mappings for the parameter labels in that object.
1856    '''
1857    def __init__(self,exprObj):
1858        self.eObj = exprObj
1859        'The expression and mappings; a :class:`ExpressionObj` object'
1860        self.compiledExpr = None
1861        'The expression as compiled byte-code'
1862        self.exprDict = {}
1863        '''dict that defines values for labels used in expression and packages
1864        referenced by functions
1865        '''
1866        self.lblLookup = {}
1867        '''Lookup table that specifies the expression label name that is
1868        tied to a particular GSAS-II parameters in the parmDict.
1869        '''
1870        self.fxnpkgdict = {}
1871        '''a dict with references to packages needed to
1872        find functions referenced in the expression.
1873        '''
1874        self.varLookup = {}
1875        '''Lookup table that specifies the GSAS-II variable(s)
1876        indexed by the expression label name. (Used for only for diagnostics
1877        not evaluation of expression.)
1878        '''
1879        # Patch: for old-style expressions with a (now removed step size)
1880        for v in self.eObj.assgnVars:
1881            if not isinstance(self.eObj.assgnVars[v], basestring):
1882                self.eObj.assgnVars[v] = self.eObj.assgnVars[v][0]
1883
1884    def SetupCalc(self,parmDict):
1885        '''Do all preparations to use the expression for computation.
1886        Adds the free parameter values to the parameter dict (parmDict).
1887        '''
1888        self.fxnpkgdict = self.eObj.CheckVars()
1889        # all is OK, compile the expression
1890        self.compiledExpr = compile(self.eObj.expression,'','eval')
1891
1892        # look at first value in parmDict to determine its type
1893        parmsInList = True
1894        for key in parmDict:
1895            val = parmDict[key]
1896            if isinstance(val, basestring):
1897                parmsInList = False
1898                break
1899            try: # check if values are in lists
1900                val = parmDict[key][0]
1901            except (TypeError,IndexError):
1902                parmsInList = False
1903            break
1904           
1905        # set up the dicts needed to speed computations
1906        self.exprDict = {}
1907        self.lblLookup = {}
1908        self.varLookup = {}
1909        for v in self.eObj.freeVars:
1910            varname = self.eObj.freeVars[v][0]
1911            varname = "::" + varname.lstrip(':').replace(' ','_').replace(':',';')
1912            self.lblLookup[varname] = v
1913            self.varLookup[v] = varname
1914            if parmsInList:
1915                parmDict[varname] = [self.eObj.freeVars[v][1],self.eObj.freeVars[v][2]]
1916            else:
1917                parmDict[varname] = self.eObj.freeVars[v][1]
1918            self.exprDict[v] = self.eObj.freeVars[v][1]
1919        for v in self.eObj.assgnVars:
1920            varname = self.eObj.assgnVars[v]
1921            if '*' in varname:
1922                varlist = LookupWildCard(varname,parmDict.keys())
1923                if len(varlist) == 0:
1924                    raise Exception,"No variables match "+str(v)
1925                for var in varlist:
1926                    self.lblLookup[var] = v
1927                if parmsInList:
1928                    self.exprDict[v] = np.array([parmDict[var][0] for var in varlist])
1929                else:
1930                    self.exprDict[v] = np.array([parmDict[var] for var in varlist])
1931                self.varLookup[v] = [var for var in varlist]
1932            elif varname in parmDict:
1933                self.lblLookup[varname] = v
1934                self.varLookup[v] = varname
1935                if parmsInList:
1936                    self.exprDict[v] = parmDict[varname][0]
1937                else:
1938                    self.exprDict[v] = parmDict[varname]
1939            else:
1940                raise Exception,"No value for variable "+str(v)
1941        self.exprDict.update(self.fxnpkgdict)
1942
1943    def UpdateVars(self,varList,valList):
1944        '''Update the dict for the expression with a set of values
1945        :param list varList: a list of variable names
1946        :param list valList: a list of corresponding values
1947        '''
1948        for var,val in zip(varList,valList):
1949            self.exprDict[self.lblLookup.get(var,'undefined: '+var)] = val
1950
1951    def UpdateDict(self,parmDict):
1952        '''Update the dict for the expression with values in a dict
1953        :param list parmDict: a dict of values some of which may be in use here
1954        '''
1955        for var in parmDict:
1956            if var in self.lblLookup:
1957                self.exprDict[self.lblLookup[var]] = parmDict[var]
1958           
1959    def EvalExpression(self):
1960        '''Evaluate an expression. Note that the expression
1961        and mapping are taken from the :class:`ExpressionObj` expression object
1962        and the parameter values were specified in :meth:`SetupCalc`.
1963        :returns: a single value for the expression. If parameter
1964        values are arrays (for example, from wild-carded variable names),
1965        the sum of the resulting expression is returned.
1966
1967        For example, if the expression is ``'A*B'``,
1968        where A is 2.0 and B maps to ``'1::Afrac:*'``, which evaluates to::
1969
1970        [0.5, 1, 0.5]
1971
1972        then the result will be ``4.0``.
1973        '''
1974        if self.compiledExpr is None:
1975            raise Exception,"EvalExpression called before SetupCalc"
1976        val = eval(self.compiledExpr,globals(),self.exprDict)
1977        if not np.isscalar(val):
1978            val = np.sum(val)
1979        return val
1980
1981
1982if __name__ == "__main__":
1983    # test equation evaluation
1984    def showEQ(calcobj):
1985        print 50*'='
1986        print calcobj.eObj.expression,'=',calcobj.EvalExpression()
1987        for v in sorted(calcobj.varLookup):
1988            print "  ",v,'=',calcobj.exprDict[v],'=',calcobj.varLookup[v]
1989        # print '  Derivatives'
1990        # for v in calcobj.derivStep.keys():
1991        #     print '    d(Expr)/d('+v+') =',calcobj.EvalDeriv(v)
1992
1993    obj = ExpressionObj()
1994
1995    obj.expression = "A*np.exp(B)"
1996    obj.assgnVars =  {'B': '0::Afrac:1'}
1997    obj.freeVars =  {'A': [u'A', 0.5, True]}
1998    #obj.CheckVars()
1999    parmDict2 = {'0::Afrac:0':[0.0,True], '0::Afrac:1': [1.0,False]}
2000    calcobj = ExpressionCalcObj(obj)
2001    calcobj.SetupCalc(parmDict2)
2002    showEQ(calcobj)
2003
2004    obj.expression = "A*np.exp(B)"
2005    obj.assgnVars =  {'B': '0::Afrac:*'}
2006    obj.freeVars =  {'A': [u'Free Prm A', 0.5, True]}
2007    #obj.CheckVars()
2008    parmDict1 = {'0::Afrac:0':1.0, '0::Afrac:1': 1.0}
2009    calcobj = ExpressionCalcObj(obj)
2010    calcobj.SetupCalc(parmDict1)
2011    showEQ(calcobj)
2012
2013    calcobj.SetupCalc(parmDict2)
2014    showEQ(calcobj)
Note: See TracBrowser for help on using the repository browser.