source: trunk/GSASIIobj.py @ 1833

Last change on this file since 1833 was 1833, checked in by toby, 7 years ago

prev. ver. also fixed linux menu problem; missed one doc correction

  • Property svn:eol-style set to native
  • Property svn:keywords set to Date Author Revision URL Id
File size: 97.5 KB
Line 
1# -*- coding: utf-8 -*-
2#GSASIIobj - data objects for GSAS-II
3########### SVN repository information ###################
4# $Date: 2015-05-02 22:32:05 +0000 (Sat, 02 May 2015) $
5# $Author: toby $
6# $Revision: 1833 $
7# $URL: trunk/GSASIIobj.py $
8# $Id: GSASIIobj.py 1833 2015-05-02 22:32:05Z toby $
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
360.. _SSGData_table:
361
362.. index::
363   single: Superspace Group Data description
364   single: Data object descriptions; Superspace Group Data
365
366Superspace groups [3+1] are interpreted by :func:`GSASIIspc.SSpcGroup`
367and the information is placed in a SSGdata object
368which is a dict with these keys:
369
370.. tabularcolumns:: |l|p{4.5in}|
371
372==========  ====================================================
373  key         explanation
374==========  ====================================================
375SSpGrp      superspace group symbol extension to space group
376            symbol, accidental spaces removed
377SSGCen      4D cell centering vectors [0,0,0,0] at least
378SSGOps      4D symmetry operations as [M,T] so that M*x+T = x'
379==========  ====================================================
380
381
382Atom Records
383------------
384
385.. _Atoms_table:
386
387.. index::
388   single: Atoms record description
389   single: Data object descriptions; Atoms record
390
391
392If ``phasedict`` points to the phase information in the data tree, then
393atoms are contained in a list of atom records (list) in
394``phasedict['Atoms']``. Also needed to read atom information
395are four pointers, ``cx,ct,cs,cia = phasedict['General']['atomPtrs']``,
396which define locations in the atom record, as shown below. Items shown are
397always present; additional ones for macromolecular phases are marked 'mm'
398
399.. tabularcolumns:: |l|p{4.5in}|
400
401==============   ====================================================
402location         explanation
403==============   ====================================================
404ct-4              mm - residue number (str)
405ct-3              mm - residue name (e.g. ALA) (str)
406ct-2              mm - chain label (str)
407ct-1              atom label (str)
408ct                atom type (str)
409ct+1              refinement flags; combination of 'F', 'X', 'U' (str)
410cx,cx+1,cx+2      the x,y and z coordinates (3 floats)
411cs                site symmetry (str)
412cs+1              site multiplicity (int)
413cia               ADP flag: Isotropic ('I') or Anisotropic ('A')
414cia+1             Uiso (float)
415cia+2...cia+7     U11, U22, U33, U12, U13, U23 (6 floats)
416atom[cia+8]       unique atom identifier (int)
417
418==============   ====================================================
419
420Drawing Atom Records
421--------------------
422
423.. _Drawing_atoms_table:
424
425.. index::
426   single: Drawing atoms record description
427   single: Data object descriptions; Drawing atoms record
428
429
430If ``phasedict`` points to the phase information in the data tree, then
431drawing atoms are contained in a list of drawing atom records (list) in
432``phasedict['Drawing']['Atoms']``. Also needed to read atom information
433are four pointers, ``cx,ct,cs,ci = phasedict['Drawing']['AtomPtrs']``,
434which define locations in the atom record, as shown below. Items shown are
435always present; additional ones for macromolecular phases are marked 'mm'
436
437.. tabularcolumns:: |l|p{4.5in}|
438
439==============   ====================================================
440location         explanation
441==============   ====================================================
442ct-4              mm - residue number (str)
443ct-3              mm - residue name (e.g. ALA) (str)
444ct-2              mm - chain label (str)
445ct-1              atom label (str)
446ct                atom type (str)
447cx,cx+1,cx+2      the x,y and z coordinates (3 floats)
448cs-1              Sym Op symbol; sym. op number + unit cell id (e.g. '1,0,-1') (str)
449cs                atom drawing style; e.g. 'balls & sticks' (str)
450cs+1              atom label style (e.g. 'name') (str)
451cs+2              atom color (RBG triplet) (int)
452cs+3              ADP flag: Isotropic ('I') or Anisotropic ('A')
453cs+4              Uiso (float)
454cs+5...cs+11      U11, U22, U33, U12, U13, U23 (6 floats)
455ci                unique atom identifier; matches source atom Id in Atom Records (int)
456==============   ====================================================
457
458Powder Diffraction Tree Items
459-----------------------------
460
461.. _Powder_table:
462
463.. index::
464   single: Powder data object description
465   single: Data object descriptions; Powder Data
466
467Every powder diffraction histogram is stored in the GSAS-II data tree
468with a top-level entry named beginning with the string "PWDR ". The
469diffraction data for that information are directly associated with
470that tree item and there are a series of children to that item. The
471routines :func:`GSASII.GSASII.GetUsedHistogramsAndPhasesfromTree`
472and :func:`GSASIIstrIO.GetUsedHistogramsAndPhases` will
473load this information into a dictionary where the child tree name is
474used as a key, and the information in the main entry is assigned
475a key of ``Data``, as outlined below.
476
477.. tabularcolumns:: |l|l|p{4in}|
478
479======================  ===============  ====================================================
480  key                      sub-key        explanation
481======================  ===============  ====================================================
482Comments                      \           Text strings extracted from the original powder
483                                          data header. These cannot be changed by the user;
484                                          it may be empty.
485Limits                       \            A list of two two element lists, as [[Ld,Hd],[L,H]]
486                                          where L and Ld are the current and default lowest
487                                          two-theta value to be used and
488                                          where H and Hd are the current and default highest
489                                          two-theta value to be used.
490Reflection Lists              \           A dict with an entry for each phase in the
491                                          histogram. The contents of each dict item
492                                          is a dict containing reflections, as described in
493                                          the :ref:`Powder Reflections <PowderRefl_table>`
494                                          description.
495Instrument Parameters         \           A list containing two dicts where the possible
496                                          keys in each dict are listed below. The value
497                                          for each item is a list containing three values:
498                                          the initial value, the current value and a
499                                          refinement flag which can have a value of
500                                          True, False or 0 where 0 indicates a value that
501                                          cannot be refined. The first and second
502                                          values are floats unless otherwise noted.
503                                          Items in the first dict are noted as [1]
504\                         Lam             Specifies a wavelength in Angstroms [1]
505\                         Lam1            Specifies the primary wavelength in
506                                          Angstrom, when an alpha1, alpha2
507                                          source is used [1]
508\                         Lam2            Specifies the secondary wavelength in
509                                          Angstrom, when an alpha1, alpha2
510                                          source is used [1]
511                          I(L2)/I(L1)     Ratio of Lam2 to Lam1 [1]           
512\                         Type            Histogram type (str) [1]:
513                                           * 'PXC' for constant wavelength x-ray
514                                           * 'PNC' for constant wavelength neutron
515                                           * 'PNT' for time of flight neutron
516\                         Zero            Two-theta zero correction in *degrees* [1]
517\                         Azimuth         Azimuthal setting angle for data recorded
518                                          with differing setting angles [1]
519\                         U, V, W         Cagliotti profile coefficients
520                                          for Gaussian instrumental broadening, where the
521                                          FWHM goes as
522                                          :math:`U \\tan^2\\theta + V \\tan\\theta + W` [1]
523\                         X, Y            Cauchy (Lorentzian) instrumental broadening
524                                          coefficients [1]
525\                         SH/L            Variant of the Finger-Cox-Jephcoat asymmetric
526                                          peak broadening ratio. Note that this is the
527                                          average between S/L and H/L where S is
528                                          sample height, H is the slit height and
529                                          L is the goniometer diameter. [1]
530\                         Polariz.        Polarization coefficient. [1]
531wtFactor                      \           A weighting factor to increase or decrease
532                                          the leverage of data in the histogram (float).
533                                          A value of 1.0 weights the data with their
534                                          standard uncertainties and a larger value
535                                          increases the weighting of the data (equivalent
536                                          to decreasing the uncertainties).
537Sample Parameters             \           Specifies a dict with parameters that describe how
538                                          the data were collected, as listed
539                                          below. Refinable parameters are a list containing
540                                          a float and a bool, where the second value
541                                          specifies if the value is refined, otherwise
542                                          the value is a float unless otherwise noted.
543\                         Scale           The histogram scale factor (refinable)
544\                         Absorption      The sample absorption coefficient as
545                                          :math:`\\mu r` where r is the radius
546                                          (refinable). Only valid for Debye-Scherrer geometry.
547\                         SurfaceRoughA   Surface roughness parameter A as defined by
548                                          Surotti,J. Appl. Cryst, 5,325-331, 1972.(refinable -
549                                          only valid for Bragg-Brentano geometry)                                         
550\                         SurfaceRoughB   Surface roughness parameter B (refinable -
551                                          only valid for Bragg-Brentano geometry)                                         
552\                         DisplaceX,      Sample displacement from goniometer center
553                          DisplaceY       where Y is along the beam direction and
554                                          X is perpendicular. Units are :math:`\\mu m`
555                                          (refinable).
556\                         Phi, Chi,       Goniometer sample setting angles, in degrees.
557                          Omega
558\                         Gonio. radius   Radius of the diffractometer in mm
559\                         InstrName       A name for the instrument, used in preparing
560                                          a CIF (str).
561\                         Force,          Variables that describe how the measurement
562                          Temperature,    was performed. Not used directly in
563                          Humidity,       any computations.
564                          Pressure,
565                          Voltage
566\                         ranId           The random-number Id for the histogram
567                                          (same value as where top-level key is ranId)
568\                         Type            Type of diffraction data, may be 'Debye-Scherrer'
569                                          or 'Bragg-Brentano' (str).
570\                         Diffuse         not in use?
571hId                           \           The number assigned to the histogram when
572                                          the project is loaded or edited (can change)
573ranId                         \           A random number id for the histogram
574                                          that does not change
575Background                    \           The background is stored as a list with where
576                                          the first item in the list is list and the second
577                                          item is a dict. The list contains the background
578                                          function and its coefficients; the dict contains
579                                          Debye diffuse terms and background peaks.
580                                          (TODO: this needs to be expanded.)
581Data                          \           The data consist of a list of 6 np.arrays
582                                          containing in order:
583
584                                           0. the x-postions (two-theta in degrees),
585                                           1. the intensity values (Yobs),
586                                           2. the weights for each Yobs value
587                                           3. the computed intensity values (Ycalc)
588                                           4. the background values
589                                           5. Yobs-Ycalc
590======================  ===============  ====================================================
591
592Powder Reflection Data Structure
593--------------------------------
594
595.. _PowderRefl_table:
596
597.. index::
598   single: Powder reflection object description
599   single: Data object descriptions; Powder Reflections
600   
601For every phase in a histogram, the ``Reflection Lists`` value is a dict
602one element of which is `'RefList'`, which is a np.array containing
603reflections. The columns in that array are documented below.
604
605==========  ====================================================
606  index         explanation
607==========  ====================================================
608 0,1,2       h,k,l (float)
609 3           multiplicity
610 4           d-space, Angstrom
611 5           pos, two-theta
612 6           sig, Gaussian width
613 7           gam, Lorenzian width
614 8           :math:`F_{obs}^2`
615 9           :math:`F_{calc}^2`
616 10          reflection phase, in degrees
617 11          intensity correction for reflection, this times
618             :math:`F_{obs}^2` or :math:`F_{calc}^2` gives Iobs or Icalc
619==========  ====================================================
620
621Single Crystal Tree Items
622-------------------------
623
624.. _Xtal_table:
625
626.. index::
627   single: Single Crystal data object description
628   single: Data object descriptions; Single crystal data
629
630Every single crystal diffraction histogram is stored in the GSAS-II data tree
631with a top-level entry named beginning with the string "HKLF ". The
632diffraction data for that information are directly associated with
633that tree item and there are a series of children to that item. The
634routines :func:`GSASII.GSASII.GetUsedHistogramsAndPhasesfromTree`
635and :func:`GSASIIstrIO.GetUsedHistogramsAndPhases` will
636load this information into a dictionary where the child tree name is
637used as a key, and the information in the main entry is assigned
638a key of ``Data``, as outlined below.
639
640.. tabularcolumns:: |l|l|p{4in}|
641
642======================  ===============  ====================================================
643  key                      sub-key        explanation
644======================  ===============  ====================================================
645Data                          \           A dict that contains the
646                                          reflection table,
647                                          as described in the
648                                          :ref:`Single Crystal Reflections
649                                          <XtalRefl_table>`
650                                          description.
651
652Instrument Parameters         \           A list containing two dicts where the possible
653                                          keys in each dict are listed below. The value
654                                          for most items is a list containing two values:
655                                          the initial value, the current value.
656                                          The first and second
657                                          values are floats unless otherwise noted.
658\                         Lam             Specifies a wavelength in Angstroms (two floats)
659\                         Type            Histogram type (two str values):
660                                           * 'SXC' for constant wavelength x-ray
661                                           * 'SNC' for constant wavelength neutron
662                                           * 'SNT' for time of flight neutron
663\                         InstrName       A name for the instrument, used in preparing
664                                          a CIF (str).
665
666wtFactor                      \           A weighting factor to increase or decrease
667                                          the leverage of data in the histogram (float).
668                                          A value of 1.0 weights the data with their
669                                          standard uncertainties and a larger value
670                                          increases the weighting of the data (equivalent
671                                          to decreasing the uncertainties).
672
673hId                           \           The number assigned to the histogram when
674                                          the project is loaded or edited (can change)
675ranId                         \           A random number id for the histogram
676                                          that does not change
677======================  ===============  ====================================================
678
679Single Crystal Reflection Data Structure
680----------------------------------------
681
682.. _XtalRefl_table:
683
684.. index::
685   single: Single Crystal reflection object description
686   single: Data object descriptions; Single Crystal Reflections
687   
688For every single crystal a histogram, the ``'Data'`` item contains
689the structure factors as an np.array in item `'RefList'`.
690The columns in that array are documented below.
691
692==========  ====================================================
693  index         explanation
694==========  ====================================================
695 0,1,2       h,k,l (float)
696 3           multiplicity
697 4           d-space, Angstrom
698 5           :math:`F_{obs}^2`
699 6           :math:`\sigma(F_{obs}^2)`
700 7           :math:`F_{calc}^2`
701 8           :math:`F_{obs}^2T`
702 9           :math:`F_{calc}^2T`
703 10          reflection phase, in degrees
704 11          intensity correction for reflection, this times
705             :math:`F_{obs}^2` or :math:`F_{calc}^2`
706             gives Iobs or Icalc
707==========  ====================================================
708
709Image Data Structure
710--------------------
711
712.. _Image_table:
713
714.. index::
715   image: Image data object description
716   image: Image object descriptions
717   
718Every 2-dimensional image is stored in the GSAS-II data tree
719with a top-level entry named beginning with the string "IMG ". The
720image data are directly associated with that tree item and there
721are a series of children to that item. The routines :func:`GSASII.GSASII.GetUsedHistogramsAndPhasesfromTree`
722and :func:`GSASIIstrIO.GetUsedHistogramsAndPhases` will
723load this information into a dictionary where the child tree name is
724used as a key, and the information in the main entry is assigned
725a key of ``Data``, as outlined below.
726
727.. tabularcolumns:: |l|l|p{4in}|
728
729======================  ======================  ====================================================
730  key                      sub-key              explanation
731======================  ======================  ====================================================
732Comments                       \                Text strings extracted from the original image data
733                                                header or a metafile. These cannot be changed by 
734                                                the user; it may be empty.                                               
735Image Controls              azmthOff            (float) The offset to be applied to an azimuthal
736                                                value. Accomodates
737                                                detector orientations other than with the detector
738                                                X-axis
739                                                horizontal.
740\                           background image    (list:str,float) The name of a tree item ("IMG ...") that is to be subtracted
741                                                during image integration multiplied by value. It must have the same size/shape as
742                                                the integrated image. NB: value < 0 for subtraction.
743\                           calibrant           (str) The material used for determining the position/orientation
744                                                of the image. The data is obtained from :func:`ImageCalibrants`
745                                                and UserCalibrants.py (supplied by user).
746\                           calibdmin           (float) The minimum d-spacing used during the last calibration run.
747\                           calibskip           (int) The number of expected diffraction lines skipped during the last
748                                                calibration run.
749\                           center              (list:floats) The [X,Y] point in detector coordinates (mm) where the direct beam
750                                                strikes the detector plane as determined by calibration. This point
751                                                does not have to be within the limits of the detector boundaries.
752\                           centerAzm           (bool) If True then the azimuth reported for the integrated slice
753                                                of the image is at the center line otherwise it is at the leading edge.
754\                           color               (str) The name of the colormap used to display the image. Default = 'Paired'.
755\                           cutoff              (float) The minimum value of I/Ib for a point selected in a diffraction ring for
756                                                calibration calculations. See pixLimit for details as how point is found.           
757\                           DetDepth            (float) Coefficient for penetration correction to distance; accounts for diffraction
758                                                ring offset at higher angles. Optionally determined by calibration.
759\                           DetDepthRef         (bool) If True then refine DetDepth during calibration/recalibration calculation.
760\                           distance            (float) The distance (mm) from sample to detector plane.
761\                           ellipses            (list:lists) Each object in ellipses is a list [center,phi,radii,color] where
762                                                center (list) is location (mm) of the ellipse center on the detector plane, phi is the
763                                                rotation of the ellipse minor axis from the x-axis, and radii are the minor & major
764                                                radii of the ellipse. If radii[0] is negative then parameters describe a hyperbola. Color
765                                                is the selected drawing color (one of 'b', 'g' ,'r') for the ellipse/hyperbola.
766\                           edgemin             (float) Not used;  parameter in EdgeFinder code.
767\                           fullIntegrate       (bool) If True then integrate over full 360 deg azimuthal range.
768\                           GonioAngles         (list:floats) The 'Omega','Chi','Phi' goniometer angles used for this image.
769                                                Required for texture calculations.
770\                           invert_x            (bool) If True display the image with the x-axis inverted.
771\                           invert_y            (bool) If True display the image with the y-axis inverted.
772\                           IOtth               (list:floats) The minimum and maximum 2-theta values to be used for integration.
773\                           LRazimuth           (list:floats) The minimum and maximum azimuth values to be used for integration.
774\                           Oblique             (list:float,bool) If True apply a detector absorption correction using the value to the
775                                                intensities obtained during integration.
776\                           outAzimuths         (int) The number of azimuth pie slices.
777\                           outChannels         (int) The number of 2-theta steps.
778\                           pixelSize           (list:ints) The X,Y dimensions (microns) of each pixel.
779\                           pixLimit            (int) A box in the image with 2*pixLimit+1 edges is searched to find the maximum.
780                                                This value (I) along with the minimum (Ib) in the box is reported by :func:`GSASIIimage.ImageLocalMax`
781                                                and subject to cutoff in :func:`GSASIIimage.makeRing`.
782                                                Locations are used to construct rings of points for calibration calcualtions.
783\                           PolaVal             (list:float,bool) If type='SASD' and if True, apply polarization correction to intensities from
784                                                integration using value.
785\                           rings               (list:lists) Each entry is [X,Y,dsp] where X & Y are lists of x,y coordinates around a
786                                                diffraction ring with the same d-spacing (dsp)
787\                           ring                (list) The x,y coordinates of the >5 points on an inner ring
788                                                selected by the user,
789\                           Range               (list) The minimum & maximum values of the image
790\                           rotation            (float) The angle between the x-axis and the vector about which the
791                                                detector is tilted. Constrained to -180 to 180 deg.     
792\                           SampleShape         (str) Currently only 'Cylinder'. Sample shape for Debye-Scherrer experiments; used for absorption
793                                                calculations.
794\                           SampleAbs           (list: float,bool) Value of absorption coefficient for Debye-Scherrer experimnents, flag if True
795                                                to cause correction to be applied.
796\                           setDefault          (bool) If True the use the image controls values for all new images to be read. (might be removed)
797\                           setRings            (bool) If True then display all the selected x,y ring positions (vida supra rings) used in the calibration.           
798\                           showLines           (bool) If True then isplay the integration limits to be used.
799\                           size                (list:int) The number of pixels on the image x & y axes
800\                           type                (str) One of 'PWDR', 'SASD' or 'REFL' for powder, small angle or reflectometry data, respectively.
801\                           tilt                (float) The angle the detector normal makes with the incident beam; range -90 to 90.
802\                           wavelength          (float) Tha radiation wavelength (Angstroms) as entered by the user (or someday obtained from the image header).
803                                               
804Masks                       Arcs                (list: lists) Each entry [2-theta,[azimuth[0],azimuth[1]],thickness] describes an arc mask
805                                                to be excluded from integration
806\                           Frames              (list:lists) Each entry describes the x,y points (3 or more - mm) that describe a frame outside
807                                                of which is excluded from recalibration and integration. Only one frame is allowed.
808\                           Points              (list:lists) Each entry [x,y,radius] (mm) describes an excluded spot on the image to be excluded
809                                                from integration.
810\                           Polygons            (list:lists) Each entry is a list of 3+ [x,y] points (mm) that describe a polygon on the image
811                                                to be excluded from integration.
812\                           Rings               (list: lists) Each entry [2-theta,thickness] describes a ring mask
813                                                to be excluded from integration.
814\                           Thresholds          (list:[tuple,list]) [(Imin,Imax),[Imin,Imax]] This gives lower and upper limits for points on the image to be included
815                                                in integrsation. The tuple is the image intensity limits and the list are those set by the user.   
816                                               
817Stress/Strain               Sample phi          (float) Sample rotation about vertical axis.
818\                           Sample z            (float) Sample translation from the calibration sample position (for Sample phi = 0)
819                                                These will be restricted by space group symmetry; result of strain fit refinement.
820\                           Type                (str) 'True' or 'Conventional': The strain model used for the calculation.
821\                           d-zero              (list:dict) Each item is for a diffraction ring on the image; all items are from the same phase
822                                                and are used to determine the strain tensor.
823                                                The dictionary items are:
824                                                'Dset': (float) True d-spacing for the diffraction ring; entered by the user.
825                                                'Dcalc': (float) Average calculated d-spacing determined from strain coeff.
826                                                'Emat': (list: float) The strain tensor elements e11, e12 & e22 (e21=e12, rest are 0)
827                                                'Esig': (list: float) Esds for Emat from fitting.
828                                                'pixLimit': (int) Search range to find highest point on ring for each data point
829                                                'cutoff': (float) I/Ib cutoff for searching.
830                                                'ImxyObs': (list: lists) [[X],[Y]] observed points to be used for strain calculations.
831                                                'ImtaObs': (list: lists) [[d],[azm]] transformed via detector calibration from ImxyObs.
832                                                'ImtaCalc': (list: lists [[d],[azm]] calculated d-spacing & azimuth from fit.
833                                               
834======================  ======================  ====================================================
835
836Parameter Dictionary
837-------------------------
838
839.. _parmDict_table:
840
841.. index::
842   single: Parameter dictionary
843
844The parameter dictionary contains all of the variable parameters for the refinement.
845The dictionary keys are the name of the parameter (<phase>:<hist>:<name>:<atom>).
846It is prepared in two ways. When loaded from the tree
847(in :meth:`GSASII.GSASII.MakeLSParmDict` and
848:meth:`GSASIIIO.ExportBaseclass.loadParmDict`),
849the values are lists with two elements: ``[value, refine flag]``
850
851When loaded from the GPX file (in
852:func:`GSASIIstrMain.Refine` and :func:`GSASIIstrMain.SeqRefine`), the value in the
853dict is the actual parameter value (usually a float, but sometimes a
854letter or string flag value (such as I or A for iso/anisotropic).
855
856
857*Classes and routines*
858----------------------
859
860'''
861import re
862import imp
863import random as ran
864import sys
865import GSASIIpath
866import GSASIImath as G2mth
867import numpy as np
868
869GSASIIpath.SetVersionNumber("$Revision: 1833 $")
870
871DefaultControls = {
872    'deriv type':'analytic Hessian',
873    'min dM/M':0.0001,'shift factor':1.,'max cyc':3,'F**2':False,
874    'UsrReject':{'minF/sig':0,'MinExt':0.01,'MaxDF/F':100.,'MaxD':500.,'MinD':0.05},
875    'Copy2Next':False,'Reverse Seq':False,
876    'Author':'no name',
877    'FreePrm1':'Sample humidity (%)',
878    'FreePrm2':'Sample voltage (V)',
879    'FreePrm3':'Applied load (MN)',
880    'SeqPseudoVars':{},'SeqParFitEqList':[],'ShowCell':False,
881    }
882'''Values to be used as defaults for the initial contents of the ``Controls``
883data tree item.
884'''
885
886def MakeUniqueLabel(lbl,labellist):
887    '''Make sure that every a label is unique against a list by adding
888    digits at the end until it is not found in list.
889
890    :param str lbl: the input label
891    :param list labellist: the labels that have already been encountered
892    :returns: lbl if not found in labellist or lbl with ``_1-9`` (or
893      ``_10-99``, etc.) appended at the end
894    '''
895    lbl = lbl.strip()
896    if not lbl: # deal with a blank label
897        lbl = '_1'
898    if lbl not in labellist:
899        labellist.append(lbl)
900        return lbl
901    i = 1
902    prefix = lbl
903    if '_' in lbl:
904        prefix = lbl[:lbl.rfind('_')]
905        suffix = lbl[lbl.rfind('_')+1:]
906        try:
907            i = int(suffix)+1
908        except: # suffix could not be parsed
909            i = 1
910            prefix = lbl
911    while prefix+'_'+str(i) in labellist:
912        i += 1
913    else:
914        lbl = prefix+'_'+str(i)
915        labellist.append(lbl)
916    return lbl
917
918PhaseIdLookup = {}
919'''dict listing phase name and random Id keyed by sequential phase index as a str;
920best to access this using :func:`LookupPhaseName`
921'''
922PhaseRanIdLookup = {}
923'''dict listing phase sequential index keyed by phase random Id;
924best to access this using :func:`LookupPhaseId`
925'''
926HistIdLookup = {}
927'''dict listing histogram name and random Id, keyed by sequential histogram index as a str;
928best to access this using :func:`LookupHistName`
929'''
930HistRanIdLookup = {}
931'''dict listing histogram sequential index keyed by histogram random Id;
932best to access this using :func:`LookupHistId`
933'''
934AtomIdLookup = {}
935'''dict listing for each phase index as a str, the atom label and atom random Id,
936keyed by atom sequential index as a str;
937best to access this using :func:`LookupAtomLabel`
938'''
939AtomRanIdLookup = {}
940'''dict listing for each phase the atom sequential index keyed by atom random Id;
941best to access this using :func:`LookupAtomId`
942'''
943ShortPhaseNames = {}
944'''a dict containing a possibly shortened and when non-unique numbered
945version of the phase name. Keyed by the phase sequential index.
946'''
947ShortHistNames = {}
948'''a dict containing a possibly shortened and when non-unique numbered
949version of the histogram name. Keyed by the histogram sequential index.
950'''
951
952VarDesc = {}
953''' This dictionary lists descriptions for GSAS-II variables,
954as set in :func:`CompileVarDesc`. See that function for a description
955for how keys and values are written.
956'''
957
958reVarDesc = {}
959''' This dictionary lists descriptions for GSAS-II variables with
960the same values as :attr:`VarDesc` except that keys have been compiled as
961regular expressions. Initialized in :func:`CompileVarDesc`.
962'''
963
964def IndexAllIds(Histograms,Phases):
965    '''Scan through the used phases & histograms and create an index
966    to the random numbers of phases, histograms and atoms. While doing this,
967    confirm that assigned random numbers are unique -- just in case lightning
968    strikes twice in the same place.
969
970    Note: this code assumes that the atom random Id (ranId) is the last
971    element each atom record.
972
973    This is called in three places (only): :func:`GSASIIstrIO.GetUsedHistogramsAndPhases`
974    (which loads the histograms and phases from a GPX file),
975    :meth:`~GSASII.GSASII.GetUsedHistogramsAndPhasesfromTree`
976    (which loads the histograms and phases from the data tree.) and
977    :meth:`GSASIIconstrGUI.UpdateConstraints`
978    (which displays & edits the constraints in a GUI)
979
980    TODO: do we need a lookup for rigid body variables?
981    '''
982    # process phases and atoms
983    PhaseIdLookup.clear()
984    PhaseRanIdLookup.clear()   
985    AtomIdLookup.clear()
986    AtomRanIdLookup.clear()
987    ShortPhaseNames.clear()
988    for ph in Phases:
989        cx,ct,cs,cia = Phases[ph]['General']['AtomPtrs']
990        ranId = Phases[ph]['ranId'] 
991        while ranId in PhaseRanIdLookup:
992            # Found duplicate random Id! note and reassign
993            print ("\n\n*** Phase "+str(ph)+" has repeated ranId. Fixing.\n")
994            Phases[ph]['ranId'] = ranId = ran.randint(0,sys.maxint)
995        pId = str(Phases[ph]['pId'])
996        PhaseIdLookup[pId] = (ph,ranId)
997        PhaseRanIdLookup[ranId] = pId
998        shortname = ph[:10]
999        while shortname in ShortPhaseNames.values():
1000            shortname = ph[:8] + ' ('+ pId + ')'
1001        ShortPhaseNames[pId] = shortname
1002        AtomIdLookup[pId] = {}
1003        AtomRanIdLookup[pId] = {}
1004        for iatm,at in enumerate(Phases[ph]['Atoms']):
1005            ranId = at[cia+8]
1006            while ranId in AtomRanIdLookup[pId]: # check for dups
1007                print ("\n\n*** Phase "+str(ph)+" atom "+str(iatm)+" has repeated ranId. Fixing.\n")
1008                at[cia+8] = ranId = ran.randint(0,sys.maxint)
1009            AtomRanIdLookup[pId][ranId] = str(iatm)
1010            if Phases[ph]['General']['Type'] == 'macromolecular':
1011                label = '%s_%s_%s_%s'%(at[ct-1],at[ct-3],at[ct-4],at[ct-2])
1012            else:
1013                label = at[ct-1]
1014            AtomIdLookup[pId][str(iatm)] = (label,ranId)
1015    # process histograms
1016    HistIdLookup.clear()
1017    HistRanIdLookup.clear()
1018    ShortHistNames.clear()
1019    for hist in Histograms:
1020        ranId = Histograms[hist]['ranId']
1021        while ranId in HistRanIdLookup:
1022            # Found duplicate random Id! note and reassign
1023            print ("\n\n*** Histogram "+str(hist)+" has repeated ranId. Fixing.\n")
1024            Histograms[hist]['ranId'] = ranId = ran.randint(0,sys.maxint)
1025        hId = str(Histograms[hist]['hId'])
1026        HistIdLookup[hId] = (hist,ranId)
1027        HistRanIdLookup[ranId] = hId
1028        shortname = hist[:15]
1029        while shortname in ShortHistNames.values():
1030            shortname = hist[:11] + ' ('+ hId + ')'
1031        ShortHistNames[hId] = shortname
1032
1033def LookupAtomId(pId,ranId):
1034    '''Get the atom number from a phase and atom random Id
1035
1036    :param int/str pId: the sequential number of the phase
1037    :param int ranId: the random Id assigned to an atom
1038
1039    :returns: the index number of the atom (str)
1040    '''
1041    if not AtomRanIdLookup:
1042        raise Exception,'Error: LookupAtomId called before IndexAllIds was run'
1043    if pId is None or pId == '':
1044        raise KeyError,'Error: phase is invalid (None or blank)'
1045    pId = str(pId)
1046    if pId not in AtomRanIdLookup:
1047        raise KeyError,'Error: LookupAtomId does not have phase '+pId
1048    if ranId not in AtomRanIdLookup[pId]:
1049        raise KeyError,'Error: LookupAtomId, ranId '+str(ranId)+' not in AtomRanIdLookup['+pId+']'
1050    return AtomRanIdLookup[pId][ranId]
1051
1052def LookupAtomLabel(pId,index):
1053    '''Get the atom label from a phase and atom index number
1054
1055    :param int/str pId: the sequential number of the phase
1056    :param int index: the index of the atom in the list of atoms
1057
1058    :returns: the label for the atom (str) and the random Id of the atom (int)
1059    '''
1060    if not AtomIdLookup:
1061        raise Exception,'Error: LookupAtomLabel called before IndexAllIds was run'
1062    if pId is None or pId == '':
1063        raise KeyError,'Error: phase is invalid (None or blank)'
1064    pId = str(pId)
1065    if pId not in AtomIdLookup:
1066        raise KeyError,'Error: LookupAtomLabel does not have phase '+pId
1067    if index not in AtomIdLookup[pId]:
1068        raise KeyError,'Error: LookupAtomLabel, ranId '+str(index)+' not in AtomRanIdLookup['+pId+']'
1069    return AtomIdLookup[pId][index]
1070
1071def LookupPhaseId(ranId):
1072    '''Get the phase number and name from a phase random Id
1073
1074    :param int ranId: the random Id assigned to a phase
1075    :returns: the sequential Id (pId) number for the phase (str)
1076    '''
1077    if not PhaseRanIdLookup:
1078        raise Exception,'Error: LookupPhaseId called before IndexAllIds was run'
1079    if ranId not in PhaseRanIdLookup:
1080        raise KeyError,'Error: LookupPhaseId does not have ranId '+str(ranId)
1081    return PhaseRanIdLookup[ranId]
1082
1083def LookupPhaseName(pId):
1084    '''Get the phase number and name from a phase Id
1085
1086    :param int/str pId: the sequential assigned to a phase
1087    :returns:  (phase,ranId) where phase is the name of the phase (str)
1088      and ranId is the random # id for the phase (int)
1089    '''
1090    if not PhaseIdLookup:
1091        raise Exception,'Error: LookupPhaseName called before IndexAllIds was run'
1092    if pId is None or pId == '':
1093        raise KeyError,'Error: phase is invalid (None or blank)'
1094    pId = str(pId)
1095    if pId not in PhaseIdLookup:
1096        raise KeyError,'Error: LookupPhaseName does not have index '+pId
1097    return PhaseIdLookup[pId]
1098
1099def LookupHistId(ranId):
1100    '''Get the histogram number and name from a histogram random Id
1101
1102    :param int ranId: the random Id assigned to a histogram
1103    :returns: the sequential Id (hId) number for the histogram (str)
1104    '''
1105    if not HistRanIdLookup:
1106        raise Exception,'Error: LookupHistId called before IndexAllIds was run'
1107    if ranId not in HistRanIdLookup:
1108        raise KeyError,'Error: LookupHistId does not have ranId '+str(ranId)
1109    return HistRanIdLookup[ranId]
1110
1111def LookupHistName(hId):
1112    '''Get the histogram number and name from a histogram Id
1113
1114    :param int/str hId: the sequential assigned to a histogram
1115    :returns:  (hist,ranId) where hist is the name of the histogram (str)
1116      and ranId is the random # id for the histogram (int)
1117    '''
1118    if not HistIdLookup:
1119        raise Exception,'Error: LookupHistName called before IndexAllIds was run'
1120    if hId is None or hId == '':
1121        raise KeyError,'Error: histogram is invalid (None or blank)'
1122    hId = str(hId)
1123    if hId not in HistIdLookup:
1124        raise KeyError,'Error: LookupHistName does not have index '+hId
1125    return HistIdLookup[hId]
1126
1127def fmtVarDescr(varname):
1128    '''Return a string with a more complete description for a GSAS-II variable
1129
1130    :param str varname: A full G2 variable name with 2 or 3 or 4
1131       colons (<p>:<h>:name[:<a>] or <p>::RBname:<r>:<t>])
1132       
1133    :returns: a string with the description
1134    '''
1135    s,l = VarDescr(varname)
1136    return s+": "+l
1137
1138def VarDescr(varname):
1139    '''Return two strings with a more complete description for a GSAS-II variable
1140
1141    :param str name: A full G2 variable name with 2 or 3 or 4
1142       colons (<p>:<h>:name[:<a>] or <p>::RBname:<r>:<t>])
1143       
1144    :returns: (loc,meaning) where loc describes what item the variable is mapped
1145      (phase, histogram, etc.) and meaning describes what the variable does.
1146    '''
1147   
1148    # special handling for parameter names without a colons
1149    # for now, assume self-defining
1150    if varname.find(':') == -1:
1151        return "Global",varname
1152       
1153    l = getVarDescr(varname)
1154    if not l:
1155        return ("invalid variable name ("+str(varname)+")!"),""
1156#        return "invalid variable name!",""
1157
1158    if not l[-1]:
1159        l[-1] = "(variable needs a definition! Set it in CompileVarDesc)"
1160
1161    if len(l) == 3:         #SASD variable name!
1162        s = 'component:'+l[1]
1163        return s,l[-1]
1164    s = ""
1165    if l[0] is not None and l[1] is not None: # HAP: keep short
1166        if l[2] == "Scale": # fix up ambigous name
1167            l[5] = "Phase fraction"
1168        if l[0] == '*':
1169            lbl = 'Seq. ref.'
1170        else:
1171            lbl = ShortPhaseNames.get(l[0],'? #'+str(l[0]))
1172        if l[1] == '*':
1173            hlbl = 'Seq. ref.'
1174        else:
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 = "Ph="+str(lbl)+" * "+str(hlbl)
1183    else:
1184        if l[2] == "Scale": # fix up ambigous name: must be scale factor, since not HAP
1185            l[5] = "Scale factor"
1186        if l[2] == 'Back': # background parameters are "special", alas
1187            s = 'Hist='+ShortHistNames.get(l[1],'? #'+str(l[1]))
1188            l[-1] += ' #'+str(l[3])
1189        elif l[4] is not None: # rigid body parameter
1190            lbl = ShortPhaseNames.get(l[0],'phase?')
1191            s = "Res #"+str(l[3])+" body #"+str(l[4])+" in "+str(lbl)
1192        elif l[3] is not None: # atom parameter,
1193            lbl = ShortPhaseNames.get(l[0],'phase?')
1194            try:
1195                albl = LookupAtomLabel(l[0],l[3])[0]
1196            except KeyError:
1197                albl = 'Atom?'
1198            s = "Atom "+str(albl)+" in "+str(lbl)
1199        elif l[0] == '*':
1200            s = "All phases "
1201        elif l[0] is not None:
1202            lbl = ShortPhaseNames.get(l[0],'phase?')
1203            s = "Phase "+str(lbl)
1204        elif l[1] == '*':
1205            s = 'All hists'
1206        elif l[1] is not None:
1207            hlbl = ShortHistNames.get(l[1],'? #'+str(l[1]))
1208            if hlbl[:4] == 'HKLF':
1209                hlbl = 'Xtl='+hlbl[5:]
1210            elif hlbl[:4] == 'PWDR':
1211                hlbl = 'Pwd='+hlbl[5:]
1212            else:
1213                hlbl = 'Hist='+hlbl
1214            s = str(hlbl)
1215    if not s:
1216        s = 'Global'
1217    return s,l[-1]
1218
1219def getVarDescr(varname):
1220    '''Return a short description for a GSAS-II variable
1221
1222    :param str name: A full G2 variable name with 2 or 3 or 4
1223       colons (<p>:<h>:name[:<a1>][:<a2>])
1224     
1225    :returns: a six element list as [`p`,`h`,`name`,`a1`,`a2`,`description`],
1226      where `p`, `h`, `a1`, `a2` are str values or `None`, for the phase number,
1227      the histogram number and the atom number; `name` will always be
1228      a str; and `description` is str or `None`.
1229      If the variable name is incorrectly formed (for example, wrong
1230      number of colons), `None` is returned instead of a list.
1231    '''
1232    l = varname.split(':')
1233    if len(l) == 2:     #SASD parameter name
1234        return varname,l[0],getDescr(l[1])
1235    if len(l) == 3:
1236        l += [None,None]
1237    elif len(l) == 4:
1238        l += [None]
1239    elif len(l) != 5:
1240        return None
1241    for i in (0,1,3,4):
1242        if l[i] == "":
1243            l[i] = None
1244    l += [getDescr(l[2])]
1245    return l
1246   
1247def CompileVarDesc():
1248    '''Set the values in the variable description lookup table (:attr:`VarDesc`)
1249    into :attr:`reVarDesc`. This is called in :func:`getDescr` so the initialization
1250    is always done before use.
1251
1252    Note that keys may contain regular expressions, where '[xyz]'
1253    matches 'x' 'y' or 'z' (equivalently '[x-z]' describes this as range of values).
1254    '.*' matches any string. For example::
1255
1256    'AUiso':'Atomic isotropic displacement parameter',
1257
1258    will match variable ``'p::AUiso:a'``.
1259    If parentheses are used in the key, the contents of those parentheses can be
1260    used in the value, such as::
1261
1262    'AU([123][123])':'Atomic anisotropic displacement parameter U\\1',
1263
1264    will match ``AU11``, ``AU23``,.. and `U11`, `U23` etc will be displayed
1265    in the value when used.
1266   
1267    '''
1268    if reVarDesc: return # already done
1269    for key,value in {
1270        # derived or other sequential vars
1271        '([abc])$' : 'Lattice parameter, \\1, from Ai and Djk', # N.B. '$' prevents match if any characters follow
1272        u'\u03B1' : u'Lattice parameter, \u03B1, from Ai and Djk',
1273        u'\u03B2' : u'Lattice parameter, \u03B2, from Ai and Djk',
1274        u'\u03B3' : u'Lattice parameter, \u03B3, from Ai and Djk',
1275        # ambiguous, alas:
1276        'Scale' : 'Phase or Histogram scale factor',
1277        # Phase vars (p::<var>)
1278        'A([0-5])' : 'Reciprocal metric tensor component \\1',
1279        'Vol' : 'Unit cell volume',
1280        # Atom vars (p::<var>:a)
1281        'dA([xyz])$' : 'change to atomic coordinate, \\1',
1282        'A([xyz])$' : '\\1 fractional atomic coordinate',
1283        'AUiso':'Atomic isotropic displacement parameter',
1284        'AU([123][123])':'Atomic anisotropic displacement parameter U\\1',
1285        'Afrac': 'Atomic occupancy parameter',
1286        # Hist & Phase (HAP) vars (p:h:<var>)
1287        'Back': 'Background term',
1288        'BkPkint;(.*)':'Background peak #\\1 intensity',
1289        'BkPkpos;(.*)':'Background peak #\\1 position',
1290        'BkPksig;(.*)':'Background peak #\\1 Gaussian width',
1291        'BkPkgam;(.*)':'Background peak #\\1 Cauchy width',
1292        'Bab([AU])': 'Babinet solvent scattering coef. \\1',
1293        'D([123][123])' : 'Anisotropic strain coef. \\1',
1294        'Extinction' : 'Extinction coef.',
1295        'MD' : 'March-Dollase coef.',
1296        'Mustrain;.*' : 'Microstrain coef.',
1297        'Size;.*' : 'Crystallite size value',
1298        'eA$' : 'Cubic mustrain value',
1299        'Ep$' : 'Primary extinction',
1300        'Es$' : 'Secondary type II extinction',
1301        'Eg$' : 'Secondary type I extinction',
1302        #Histogram vars (:h:<var>)
1303        'Absorption' : 'Absorption coef.',
1304        'Displace([XY])' : 'Debye-Scherrer sample displacement \\1',
1305        'Lam' : 'Wavelength',
1306        'Polariz\.' : 'Polarization correction',
1307        'SH/L' : 'FCJ peak asymmetry correction',
1308        '([UVW])$' : 'Gaussian instrument broadening \\1',
1309        '([XY])$' : 'Cauchy instrument broadening \\1',
1310        'Zero' : 'Debye-Scherrer zero correction',
1311        'nDebye' : 'Debye model background corr. terms',
1312        'nPeaks' : 'Fixed peak background corr. terms',
1313        'RBV.*' : 'Vector rigid body parameter',
1314        'RBR.*' : 'Residue rigid body parameter',
1315        'RBRO([aijk])' : 'Residue rigid body orientation parameter',
1316        'RBRP([xyz])' : 'Residue rigid body position parameter',
1317        'RBRTr;.*' : 'Residue rigid body torsion parameter',
1318        'RBR([TLS])([123AB][123AB])' : 'Residue rigid body group disp. param.',
1319        'constr([0-9]*)' : 'Parameter from constraint',
1320        # supersymmetry parameters  p::<var>:a:o 'Flen','Fcent'?
1321        'mV([0-2])$' : 'Modulation vector component \\1',
1322        'Fsin'  :   'Sin site fraction modulation',
1323        'Fcos'  :   'Cos site fraction modulation',
1324        'Fzero'  :   'Crenel function offset',
1325        'Fwid'   :   'Crenel function width',
1326        'Tzero'  :   'Sawtooth/ZigZag location',
1327        '([XYZ])slope': 'Sawtooth/ZigZag slope for \\1',
1328        '([XYZ])sin'  : 'Sin position wave for \\1',
1329        '([XYZ])cos'  : 'Cos position wave for \\1',
1330        'U([123][123])sin$' :  'Sin thermal wave for U\\1',
1331        'U([123][123])cos$' :  'Cos thermal wave for U\\1',
1332        'M([XYZ])sin$' :  'Sin mag. moment wave for \\1',
1333        'M([XYZ])cos$' :  'Cos mag. moment wave for \\1',
1334        # SASD vars (l:<var>;l = component)
1335        'Aspect ratio' : 'Particle aspect ratio',
1336        'Length' : 'Cylinder length',
1337        'Diameter' : 'Cylinder/disk diameter',
1338        'Thickness' : 'Disk thickness',
1339        'Dist' : 'Interparticle distance',
1340        'VolFr' : 'Dense scatterer volume fraction',
1341        'epis' : 'Sticky sphere epsilon',
1342        'Sticky' : 'Stickyness',
1343        'Depth' : 'Well depth',
1344        'Width' : 'Well width',
1345        'Volume' : 'Particle volume',
1346        'Radius' : 'Sphere/cylinder/disk radius',
1347        'Mean' : 'Particle mean radius',
1348        'StdDev' : 'Standard deviation in Mean',
1349        'G$': 'Guinier prefactor',
1350        'Rg$': 'Guinier radius of gyration',
1351        'B$': 'Porod prefactor',
1352        'P$': 'Porod power',
1353        'Cutoff': 'Porod cutoff',
1354        'PkInt': 'Bragg peak intensity',
1355        'PkPos': 'Bragg peak position',
1356        'PkSig': 'Bragg peak sigma',
1357        'PkGam': 'Bragg peak gamma',
1358        'e([12][12])' : 'strain tensor e\1',   # strain vars e11, e22, e12
1359        'Dcalc': 'Calc. d-spacing',
1360        'Back$': 'background parameter',
1361        'pos$': 'peak position',
1362        'int$': 'peak intensity',
1363        }.items():
1364        VarDesc[key] = value
1365        reVarDesc[re.compile(key)] = value
1366
1367def getDescr(name):
1368    '''Return a short description for a GSAS-II variable
1369
1370    :param str name: The descriptive part of the variable name without colons (:)
1371     
1372    :returns: a short description or None if not found
1373    '''
1374
1375    CompileVarDesc() # compile the regular expressions, if needed
1376    for key in reVarDesc:
1377        m = key.match(name)
1378        if m:
1379            reVarDesc[key]
1380            return m.expand(reVarDesc[key])
1381    return None
1382
1383def GenWildCard(varlist):
1384    '''Generate wildcard versions of G2 variables. These introduce '*'
1385    for a phase, histogram or atom number (but only for one of these
1386    fields) but only when there is more than one matching variable in the
1387    input variable list. So if the input is this::
1388   
1389      varlist = ['0::AUiso:0', '0::AUiso:1', '1::AUiso:0']
1390
1391    then the output will be this::
1392   
1393       wildList = ['*::AUiso:0', '0::AUiso:*']
1394
1395    :param list varlist: an input list of GSAS-II variable names
1396      (such as 0::AUiso:0)
1397
1398    :returns: wildList, the generated list of wild card variable names.
1399    '''
1400    wild = []
1401    for i in (0,1,3):
1402        currentL = varlist[:]
1403        while currentL:
1404            item1 = currentL.pop(0)
1405            i1splt = item1.split(':')
1406            if i >= len(i1splt): continue
1407            if i1splt[i]:
1408                nextL = []
1409                i1splt[i] = '[0-9]+'
1410                rexp = re.compile(':'.join(i1splt))
1411                matchlist = [item1]
1412                for nxtitem in currentL:
1413                    if rexp.match(nxtitem):
1414                        matchlist += [nxtitem]
1415                    else:
1416                        nextL.append(nxtitem)
1417                if len(matchlist) > 1:
1418                    i1splt[i] = '*'
1419                    wild.append(':'.join(i1splt))
1420                currentL = nextL
1421    return wild
1422
1423def LookupWildCard(varname,varlist):
1424    '''returns a list of variable names from list varname
1425    that match wildcard name in varname
1426   
1427    :param str varname: a G2 variable name containing a wildcard
1428      (such as \*::var)
1429    :param list varlist: the list of all variable names used in
1430      the current project
1431    :returns: a list of matching GSAS-II variables (may be empty) 
1432    '''
1433    rexp = re.compile(varname.replace('*','[0-9]+'))
1434    return sorted([var for var in varlist if rexp.match(var)])
1435
1436
1437def _lookup(dic,key):
1438    '''Lookup a key in a dictionary, where None returns an empty string
1439    but an unmatched key returns a question mark. Used in :class:`G2VarObj`
1440    '''
1441    if key is None:
1442        return ""
1443    elif key == "*":
1444        return "*"
1445    else:
1446        return dic.get(key,'?')
1447
1448class G2VarObj(object):
1449    '''Defines a GSAS-II variable either using the phase/atom/histogram
1450    unique Id numbers or using a character string that specifies
1451    variables by phase/atom/histogram number (which can change).
1452    Note that :func:`LoadID` should be used to (re)load the current Ids
1453    before creating or later using the G2VarObj object.
1454
1455    This can store rigid body variables, but does not translate the residue # and
1456    body # to/from random Ids
1457
1458    A :class:`G2VarObj` object can be created with a single parameter:
1459   
1460    :param str/tuple varname: a single value can be used to create a :class:`G2VarObj`
1461      object. If a string, it must be of form "p:h:var" or "p:h:var:a", where
1462
1463     * p is the phase number (which may be left blank or may be '*' to indicate all phases);
1464     * h is the histogram number (which may be left blank or may be '*' to indicate all histograms);
1465     * a is the atom number (which may be left blank in which case the third colon is omitted).
1466       The atom number can be specified as '*' if a phase number is specified (not as '*').
1467       For rigid body variables, specify a will be a string of form "residue:body#"
1468
1469      Alternately a single tuple of form (Phase,Histogram,VarName,AtomID) can be used, where
1470      Phase, Histogram, and AtomID are None or are ranId values (or one can be '*')
1471      and VarName is a string. Note that if Phase is '*' then the AtomID is an atom number.
1472      For a rigid body variables, AtomID is a string of form "residue:body#".
1473
1474    If four positional arguments are supplied, they are:
1475
1476    :param str/int phasenum: The number for the phase (or None or '*')
1477    :param str/int histnum: The number for the histogram (or None or '*')
1478    :param str varname: a single value can be used to create a :class:`G2VarObj`
1479    :param str/int atomnum: The number for the atom (or None or '*')
1480   
1481    '''
1482    IDdict = {}
1483    IDdict['phases'] = {}
1484    IDdict['hists'] = {}
1485    IDdict['atoms'] = {}
1486    def __init__(self,*args):
1487        self.phase = None
1488        self.histogram = None
1489        self.name = ''
1490        self.atom = None
1491        if len(args) == 1 and (type(args[0]) is list or type(args[0]) is tuple) and len(args[0]) == 4:
1492            # single arg with 4 values
1493            self.phase,self.histogram,self.name,self.atom = args[0]
1494        elif len(args) == 1 and ':' in args[0]:
1495            #parse a string
1496            lst = args[0].split(':')
1497            if lst[0] == '*':
1498                self.phase = '*'
1499                if len(lst) > 3:
1500                    self.atom = lst[3]
1501                self.histogram = HistIdLookup.get(lst[1],[None,None])[1]
1502            elif lst[1] == '*':           
1503                self.histogram = '*'
1504                self.phase = PhaseIdLookup.get(lst[0],[None,None])[1]
1505            else:
1506                self.histogram = HistIdLookup.get(lst[1],[None,None])[1]
1507                self.phase = PhaseIdLookup.get(lst[0],[None,None])[1]
1508                if len(lst) == 4:
1509                    if lst[3] == '*':
1510                        self.atom = '*'
1511                    else:
1512                        self.atom = AtomIdLookup[lst[0]].get(lst[3],[None,None])[1]
1513                elif len(lst) == 5:
1514                    self.atom = lst[3]+":"+lst[4]
1515                elif len(lst) == 3:
1516                    pass
1517                else:
1518                    raise Exception,"Too many colons in var name "+str(args[0])
1519            self.name = lst[2]
1520        elif len(args) == 4:
1521            if args[0] == '*':
1522                self.phase = '*'
1523                self.atom = args[3]
1524            else:
1525                self.phase = PhaseIdLookup.get(str(args[0]),[None,None])[1]
1526                if args[3] == '*':
1527                    self.atom = '*'
1528                elif args[0] is not None:
1529                    self.atom = AtomIdLookup[args[0]].get(str(args[3]),[None,None])[1]
1530            if args[1] == '*':
1531                self.histogram = '*'
1532            else:
1533                self.histogram = HistIdLookup.get(str(args[1]),[None,None])[1]
1534            self.name = args[2]
1535        else:
1536            raise Exception,"Incorrectly called GSAS-II parameter name"
1537
1538        #print "DEBUG: created ",self.phase,self.histogram,self.name,self.atom
1539
1540    def __str__(self):
1541        return self.varname()
1542
1543    def varname(self):
1544        '''Formats the GSAS-II variable name as a "traditional" GSAS-II variable
1545        string (p:h:<var>:a) or (p:h:<var>)
1546
1547        :returns: the variable name as a str
1548        '''
1549        a = ""
1550        if self.phase == "*":
1551            ph = "*"
1552            if self.atom:
1553                a = ":" + str(self.atom)
1554        else:
1555            ph = _lookup(PhaseRanIdLookup,self.phase)
1556            if self.atom == '*':
1557                a = ':*'
1558            elif self.atom:
1559                if ":" in str(self.atom):
1560                    a = ":" + str(self.atom)
1561                elif ph in AtomRanIdLookup:
1562                    a = ":" + AtomRanIdLookup[ph].get(self.atom,'?')
1563                else:
1564                    a = ":?"
1565        if self.histogram == "*":
1566            hist = "*"
1567        else:
1568            hist = _lookup(HistRanIdLookup,self.histogram)
1569        s = (ph + ":" + hist + ":" + str(self.name)) + a
1570        return s
1571   
1572    def __repr__(self):
1573        '''Return the detailed contents of the object
1574        '''
1575        s = "<"
1576        if self.phase == '*':
1577            s += "Phases: all; "
1578            if self.atom is not None:
1579                if ":" in str(self.atom):
1580                    s += "Rigid body" + str(self.atom) + "; "
1581                else:
1582                    s += "Atom #" + str(self.atom) + "; "
1583        elif self.phase is not None:
1584            ph =  _lookup(PhaseRanIdLookup,self.phase)
1585            s += "Phase: rId=" + str(self.phase) + " (#"+ ph + "); "
1586            if self.atom == '*':
1587                s += "Atoms: all; "
1588            elif ":" in self(self.atom):
1589                s += "Rigid body" + str(self.atom) + "; "
1590            elif self.atom is not None:
1591                s += "Atom rId=" + str(self.atom)
1592                if ph in AtomRanIdLookup:
1593                    s += " (#" + AtomRanIdLookup[ph].get(self.atom,'?') + "); "
1594                else:
1595                    s += " (#? -- not found!); "
1596        if self.histogram == '*':
1597            s += "Histograms: all; "
1598        elif self.histogram is not None:
1599            hist = _lookup(HistRanIdLookup,self.histogram)
1600            s += "Histogram: rId=" + str(self.histogram) + " (#"+ hist + "); "
1601        s += 'Variable name="' + str(self.name) + '">'
1602        return s+" ("+self.varname()+")"
1603
1604    def __eq__(self, other):
1605        if type(other) is type(self):
1606            return (self.phase == other.phase and
1607                    self.histogram == other.histogram and
1608                    self.name == other.name and
1609                    self.atom == other.atom)
1610        return False
1611
1612    def _show(self):
1613        'For testing, shows the current lookup table'
1614        print 'phases', self.IDdict['phases']
1615        print 'hists', self.IDdict['hists']
1616        print 'atomDict', self.IDdict['atoms']
1617
1618#==========================================================================
1619# shortcut routines
1620exp = np.exp
1621sind = sin = s = lambda x: np.sin(x*np.pi/180.)
1622cosd = cos = c = lambda x: np.cos(x*np.pi/180.)
1623tand = tan = t = lambda x: np.tan(x*np.pi/180.)
1624sqrt = sq = lambda x: np.sqrt(x)
1625pi = lambda: np.pi
1626class ExpressionObj(object):
1627    '''Defines an object with a user-defined expression, to be used for
1628    secondary fits or restraints. Object is created null, but is changed
1629    using :meth:`LoadExpression`. This contains only the minimum
1630    information that needs to be stored to save and load the expression
1631    and how it is mapped to GSAS-II variables.
1632    '''
1633    def __init__(self):
1634        self.expression = ''
1635        'The expression as a text string'
1636        self.assgnVars = {}
1637        '''A dict where keys are label names in the expression mapping to a GSAS-II
1638        variable. The value a G2 variable name.
1639        Note that the G2 variable name may contain a wild-card and correspond to
1640        multiple values.
1641        '''
1642        self.freeVars = {}
1643        '''A dict where keys are label names in the expression mapping to a free
1644        parameter. The value is a list with:
1645
1646         * a name assigned to the parameter
1647         * a value for to the parameter and
1648         * a flag to determine if the variable is refined.
1649        ''' 
1650        self.depVar = None
1651
1652        self.lastError = ('','')
1653        '''Shows last encountered error in processing expression
1654        (list of 1-3 str values)'''
1655
1656    def LoadExpression(self,expr,exprVarLst,varSelect,varName,varValue,varRefflag):
1657        '''Load the expression and associated settings into the object. Raises
1658        an exception if the expression is not parsed, if not all functions
1659        are defined or if not all needed parameter labels in the expression
1660        are defined.
1661
1662        This will not test if the variable referenced in these definitions
1663        are actually in the parameter dictionary. This is checked when the
1664        computation for the expression is done in :meth:`SetupCalc`.
1665       
1666        :param str expr: the expression
1667        :param list exprVarLst: parameter labels found in the expression
1668        :param dict varSelect: this will be 0 for Free parameters
1669          and non-zero for expression labels linked to G2 variables.
1670        :param dict varName: Defines a name (str) associated with each free parameter
1671        :param dict varValue: Defines a value (float) associated with each free parameter
1672        :param dict varRefflag: Defines a refinement flag (bool)
1673          associated with each free parameter
1674        '''
1675        self.expression = expr
1676        self.compiledExpr = None
1677        self.freeVars = {}
1678        self.assgnVars = {}
1679        for v in exprVarLst:
1680            if varSelect[v] == 0:
1681                self.freeVars[v] = [
1682                    varName.get(v),
1683                    varValue.get(v),
1684                    varRefflag.get(v),
1685                    ]
1686            else:
1687                self.assgnVars[v] = varName[v]
1688        self.CheckVars()
1689
1690    def EditExpression(self,exprVarLst,varSelect,varName,varValue,varRefflag):
1691        '''Load the expression and associated settings from the object into
1692        arrays used for editing.
1693
1694        :param list exprVarLst: parameter labels found in the expression
1695        :param dict varSelect: this will be 0 for Free parameters
1696          and non-zero for expression labels linked to G2 variables.
1697        :param dict varName: Defines a name (str) associated with each free parameter
1698        :param dict varValue: Defines a value (float) associated with each free parameter
1699        :param dict varRefflag: Defines a refinement flag (bool)
1700          associated with each free parameter
1701
1702        :returns: the expression as a str
1703        '''
1704        for v in self.freeVars:
1705            varSelect[v] = 0
1706            varName[v] = self.freeVars[v][0]
1707            varValue[v] = self.freeVars[v][1]
1708            varRefflag[v] = self.freeVars[v][2]
1709        for v in self.assgnVars:
1710            varSelect[v] = 1
1711            varName[v] = self.assgnVars[v]
1712        return self.expression
1713
1714    def GetVaried(self):
1715        'Returns the names of the free parameters that will be refined'
1716        return ["::"+self.freeVars[v][0] for v in self.freeVars if self.freeVars[v][2]]
1717
1718    def GetVariedVarVal(self):
1719        'Returns the names and values of the free parameters that will be refined'
1720        return [("::"+self.freeVars[v][0],self.freeVars[v][1]) for v in self.freeVars if self.freeVars[v][2]]
1721
1722    def UpdateVariedVars(self,varyList,values):
1723        'Updates values for the free parameters (after a refinement); only updates refined vars'
1724        for v in self.freeVars:
1725            if not self.freeVars[v][2]: continue
1726            if "::"+self.freeVars[v][0] not in varyList: continue
1727            indx = varyList.index("::"+self.freeVars[v][0])
1728            self.freeVars[v][1] = values[indx]
1729
1730    def GetIndependentVars(self):
1731        'Returns the names of the required independent parameters used in expression'
1732        return [self.assgnVars[v] for v in self.assgnVars]
1733
1734    def CheckVars(self):
1735        '''Check that the expression can be parsed, all functions are
1736        defined and that input loaded into the object is internally
1737        consistent. If not an Exception is raised.
1738
1739        :returns: a dict with references to packages needed to
1740          find functions referenced in the expression.
1741        '''
1742        ret = self.ParseExpression(self.expression)
1743        if not ret:
1744            raise Exception("Expression parse error")
1745        exprLblList,fxnpkgdict = ret
1746        # check each var used in expression is defined
1747        defined = self.assgnVars.keys() + self.freeVars.keys()
1748        notfound = []
1749        for var in exprLblList:
1750            if var not in defined:
1751                notfound.append(var)
1752        if notfound:
1753            msg = 'Not all variables defined'
1754            msg1 = 'The following variables were not defined: '
1755            msg2 = ''
1756            for var in notfound:
1757                if msg: msg += ', '
1758                msg += var
1759            self.lastError = (msg1,'  '+msg2)
1760            raise Exception(msg)
1761        return fxnpkgdict
1762
1763    def ParseExpression(self,expr):
1764        '''Parse an expression and return a dict of called functions and
1765        the variables used in the expression. Returns None in case an error
1766        is encountered. If packages are referenced in functions, they are loaded
1767        and the functions are looked up into the modules global
1768        workspace.
1769       
1770        Note that no changes are made to the object other than
1771        saving an error message, so that this can be used for testing prior
1772        to the save.
1773
1774        :returns: a list of used variables
1775        '''
1776        self.lastError = ('','')
1777        import ast
1778        def FindFunction(f):
1779            '''Find the object corresponding to function f
1780            :param str f: a function name such as 'numpy.exp'
1781            :returns: (pkgdict,pkgobj) where pkgdict contains a dict
1782              that defines the package location(s) and where pkgobj
1783              defines the object associated with the function.
1784              If the function is not found, pkgobj is None.
1785            '''
1786            df = f.split('.')
1787            pkgdict = {}
1788            # no listed package, try in current namespace
1789            if len(df) == 1: 
1790                try:
1791                    fxnobj = eval(f)
1792                    return pkgdict,fxnobj
1793                except (AttributeError, NameError):
1794                    return None,None
1795            else:
1796                try:
1797                    fxnobj = eval(f)
1798                    pkgdict[df[0]] = eval(df[0])
1799                    return pkgdict,fxnobj
1800                except (AttributeError, NameError):
1801                    pass
1802            # includes a package, lets try to load the packages
1803            pkgname = ''
1804            path = sys.path
1805            for pkg in f.split('.')[:-1]: # if needed, descend down the tree
1806                if pkgname:
1807                    pkgname += '.' + pkg
1808                else:
1809                    pkgname = pkg
1810                fp = None
1811                try:
1812                    fp, fppath,desc = imp.find_module(pkg,path)
1813                    pkgobj = imp.load_module(pkg,fp,fppath,desc)
1814                    pkgdict[pkgname] = pkgobj
1815                    path = [fppath]
1816                except Exception as msg:
1817                    print('load of '+pkgname+' failed with error='+str(msg))
1818                    return {},None
1819                finally:
1820                    if fp: fp.close()
1821                try:
1822                    #print 'before',pkgdict.keys()
1823                    fxnobj = eval(f,globals(),pkgdict)
1824                    #print 'after 1',pkgdict.keys()
1825                    #fxnobj = eval(f,pkgdict)
1826                    #print 'after 2',pkgdict.keys()
1827                    return pkgdict,fxnobj
1828                except:
1829                    continue
1830            return None # not found
1831        def ASTtransverse(node,fxn=False):
1832            '''Transverse a AST-parsed expresson, compiling a list of variables
1833            referenced in the expression. This routine is used recursively.
1834
1835            :returns: varlist,fxnlist where
1836              varlist is a list of referenced variable names and
1837              fxnlist is a list of used functions
1838            '''
1839            varlist = []
1840            fxnlist = []
1841            if isinstance(node, list):
1842                for b in node:
1843                    v,f = ASTtransverse(b,fxn)
1844                    varlist += v
1845                    fxnlist += f
1846            elif isinstance(node, ast.AST):
1847                for a, b in ast.iter_fields(node):
1848                    if isinstance(b, ast.AST):
1849                        if a == 'func': 
1850                            fxnlist += ['.'.join(ASTtransverse(b,True)[0])]
1851                            continue
1852                        v,f = ASTtransverse(b,fxn)
1853                        varlist += v
1854                        fxnlist += f
1855                    elif isinstance(b, list):
1856                        v,f = ASTtransverse(b,fxn)
1857                        varlist += v
1858                        fxnlist += f
1859                    elif node.__class__.__name__ == "Name":
1860                        varlist += [b]
1861                    elif fxn and node.__class__.__name__ == "Attribute":
1862                        varlist += [b]
1863            return varlist,fxnlist
1864        try:
1865            exprast = ast.parse(expr)
1866        except SyntaxError as err:
1867            s = ''
1868            import traceback
1869            for i in traceback.format_exc().splitlines()[-3:-1]:
1870                if s: s += "\n"
1871                s += str(i)
1872            self.lastError = ("Error parsing expression:",s)
1873            return
1874        # find the variables & functions
1875        v,f = ASTtransverse(exprast)
1876        varlist = sorted(list(set(v)))
1877        fxnlist = list(set(f))
1878        pkgdict = {}
1879        # check the functions are defined
1880        for fxn in fxnlist:
1881            fxndict,fxnobj = FindFunction(fxn)
1882            if not fxnobj:
1883                self.lastError = ("Error: Invalid function",fxn,
1884                                  "is not defined")
1885                return
1886            if not hasattr(fxnobj,'__call__'):
1887                self.lastError = ("Error: Not a function.",fxn,
1888                                  "cannot be called as a function")
1889                return
1890            pkgdict.update(fxndict)
1891        return varlist,pkgdict
1892
1893    def GetDepVar(self):
1894        'return the dependent variable, or None'
1895        return self.depVar
1896
1897    def SetDepVar(self,var):
1898        'Set the dependent variable, if used'
1899        self.depVar = var
1900#==========================================================================
1901class ExpressionCalcObj(object):
1902    '''An object used to evaluate an expression from a :class:`ExpressionObj`
1903    object.
1904   
1905    :param ExpressionObj exprObj: a :class:`~ExpressionObj` expression object with
1906      an expression string and mappings for the parameter labels in that object.
1907    '''
1908    def __init__(self,exprObj):
1909        self.eObj = exprObj
1910        'The expression and mappings; a :class:`ExpressionObj` object'
1911        self.compiledExpr = None
1912        'The expression as compiled byte-code'
1913        self.exprDict = {}
1914        '''dict that defines values for labels used in expression and packages
1915        referenced by functions
1916        '''
1917        self.lblLookup = {}
1918        '''Lookup table that specifies the expression label name that is
1919        tied to a particular GSAS-II parameters in the parmDict.
1920        '''
1921        self.fxnpkgdict = {}
1922        '''a dict with references to packages needed to
1923        find functions referenced in the expression.
1924        '''
1925        self.varLookup = {}
1926        '''Lookup table that specifies the GSAS-II variable(s)
1927        indexed by the expression label name. (Used for only for diagnostics
1928        not evaluation of expression.)
1929        '''
1930        # Patch: for old-style expressions with a (now removed step size)
1931        for v in self.eObj.assgnVars:
1932            if not isinstance(self.eObj.assgnVars[v], basestring):
1933                self.eObj.assgnVars[v] = self.eObj.assgnVars[v][0]
1934
1935    def SetupCalc(self,parmDict):
1936        '''Do all preparations to use the expression for computation.
1937        Adds the free parameter values to the parameter dict (parmDict).
1938        '''
1939        self.fxnpkgdict = self.eObj.CheckVars()
1940        # all is OK, compile the expression
1941        self.compiledExpr = compile(self.eObj.expression,'','eval')
1942
1943        # look at first value in parmDict to determine its type
1944        parmsInList = True
1945        for key in parmDict:
1946            val = parmDict[key]
1947            if isinstance(val, basestring):
1948                parmsInList = False
1949                break
1950            try: # check if values are in lists
1951                val = parmDict[key][0]
1952            except (TypeError,IndexError):
1953                parmsInList = False
1954            break
1955           
1956        # set up the dicts needed to speed computations
1957        self.exprDict = {}
1958        self.lblLookup = {}
1959        self.varLookup = {}
1960        for v in self.eObj.freeVars:
1961            varname = self.eObj.freeVars[v][0]
1962            varname = "::" + varname.lstrip(':').replace(' ','_').replace(':',';')
1963            self.lblLookup[varname] = v
1964            self.varLookup[v] = varname
1965            if parmsInList:
1966                parmDict[varname] = [self.eObj.freeVars[v][1],self.eObj.freeVars[v][2]]
1967            else:
1968                parmDict[varname] = self.eObj.freeVars[v][1]
1969            self.exprDict[v] = self.eObj.freeVars[v][1]
1970        for v in self.eObj.assgnVars:
1971            varname = self.eObj.assgnVars[v]
1972            if '*' in varname:
1973                varlist = LookupWildCard(varname,parmDict.keys())
1974                if len(varlist) == 0:
1975                    raise Exception,"No variables match "+str(v)
1976                for var in varlist:
1977                    self.lblLookup[var] = v
1978                if parmsInList:
1979                    self.exprDict[v] = np.array([parmDict[var][0] for var in varlist])
1980                else:
1981                    self.exprDict[v] = np.array([parmDict[var] for var in varlist])
1982                self.varLookup[v] = [var for var in varlist]
1983            elif varname in parmDict:
1984                self.lblLookup[varname] = v
1985                self.varLookup[v] = varname
1986                if parmsInList:
1987                    self.exprDict[v] = parmDict[varname][0]
1988                else:
1989                    self.exprDict[v] = parmDict[varname]
1990            else:
1991                raise Exception,"No value for variable "+str(v)
1992        self.exprDict.update(self.fxnpkgdict)
1993
1994    def UpdateVars(self,varList,valList):
1995        '''Update the dict for the expression with a set of values
1996        :param list varList: a list of variable names
1997        :param list valList: a list of corresponding values
1998        '''
1999        for var,val in zip(varList,valList):
2000            self.exprDict[self.lblLookup.get(var,'undefined: '+var)] = val
2001
2002    def UpdateDict(self,parmDict):
2003        '''Update the dict for the expression with values in a dict
2004        :param list parmDict: a dict of values some of which may be in use here
2005        '''
2006        for var in parmDict:
2007            if var in self.lblLookup:
2008                self.exprDict[self.lblLookup[var]] = parmDict[var]
2009           
2010    def EvalExpression(self):
2011        '''Evaluate an expression. Note that the expression
2012        and mapping are taken from the :class:`ExpressionObj` expression object
2013        and the parameter values were specified in :meth:`SetupCalc`.
2014        :returns: a single value for the expression. If parameter
2015        values are arrays (for example, from wild-carded variable names),
2016        the sum of the resulting expression is returned.
2017
2018        For example, if the expression is ``'A*B'``,
2019        where A is 2.0 and B maps to ``'1::Afrac:*'``, which evaluates to::
2020
2021        [0.5, 1, 0.5]
2022
2023        then the result will be ``4.0``.
2024        '''
2025        if self.compiledExpr is None:
2026            raise Exception,"EvalExpression called before SetupCalc"
2027        val = eval(self.compiledExpr,globals(),self.exprDict)
2028        if not np.isscalar(val):
2029            val = np.sum(val)
2030        return val
2031       
2032class G2Exception(Exception):
2033    def __init__(self,msg):
2034        self.msg = msg
2035    def __str__(self):
2036        return repr(self.msg)
2037
2038
2039if __name__ == "__main__":
2040    # test equation evaluation
2041    def showEQ(calcobj):
2042        print 50*'='
2043        print calcobj.eObj.expression,'=',calcobj.EvalExpression()
2044        for v in sorted(calcobj.varLookup):
2045            print "  ",v,'=',calcobj.exprDict[v],'=',calcobj.varLookup[v]
2046        # print '  Derivatives'
2047        # for v in calcobj.derivStep.keys():
2048        #     print '    d(Expr)/d('+v+') =',calcobj.EvalDeriv(v)
2049
2050    obj = ExpressionObj()
2051
2052    obj.expression = "A*np.exp(B)"
2053    obj.assgnVars =  {'B': '0::Afrac:1'}
2054    obj.freeVars =  {'A': [u'A', 0.5, True]}
2055    #obj.CheckVars()
2056    parmDict2 = {'0::Afrac:0':[0.0,True], '0::Afrac:1': [1.0,False]}
2057    calcobj = ExpressionCalcObj(obj)
2058    calcobj.SetupCalc(parmDict2)
2059    showEQ(calcobj)
2060
2061    obj.expression = "A*np.exp(B)"
2062    obj.assgnVars =  {'B': '0::Afrac:*'}
2063    obj.freeVars =  {'A': [u'Free Prm A', 0.5, True]}
2064    #obj.CheckVars()
2065    parmDict1 = {'0::Afrac:0':1.0, '0::Afrac:1': 1.0}
2066    calcobj = ExpressionCalcObj(obj)
2067    calcobj.SetupCalc(parmDict1)
2068    showEQ(calcobj)
2069
2070    calcobj.SetupCalc(parmDict2)
2071    showEQ(calcobj)
Note: See TracBrowser for help on using the repository browser.