source: trunk/GSASIIobj.py @ 2675

Last change on this file since 2675 was 2675, checked in by toby, 8 years ago

Make new PDF names unique; plot G(r); minor fixes

  • Property svn:eol-style set to native
  • Property svn:keywords set to Date Author Revision URL Id
File size: 101.6 KB
Line 
1# -*- coding: utf-8 -*-
2#GSASIIobj - data objects for GSAS-II
3########### SVN repository information ###################
4# $Date: 2017-01-30 20:38:50 +0000 (Mon, 30 Jan 2017) $
5# $Author: toby $
6# $Revision: 2675 $
7# $URL: trunk/GSASIIobj.py $
8# $Id: GSASIIobj.py 2675 2017-01-30 20:38:50Z 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: 2675 $")
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,'HatomFix':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'''
885def StripUnicode(string,subs='.'):
886    '''Strip non-ASCII characters from strings
887   
888    :param str string: string to strip Unicode characters from
889    :param str subs: character(s) to place into string in place of each
890      Unicode character. Defaults to '.'
891
892    :returns: a new string with only ASCII characters
893    '''
894    s = ''
895    for c in string:
896        if ord(c) < 128:
897            s += c
898        else:
899            s += subs
900    return s.encode('ascii','replace')
901
902def MakeUniqueLabel(lbl,labellist):
903    '''Make sure that every a label is unique against a list by adding
904    digits at the end until it is not found in list.
905
906    :param str lbl: the input label
907    :param list labellist: the labels that have already been encountered
908    :returns: lbl if not found in labellist or lbl with ``_1-9`` (or
909      ``_10-99``, etc.) appended at the end
910    '''
911    lbl = StripUnicode(lbl.strip(),'_')
912    if not lbl: # deal with a blank label
913        lbl = '_1'
914    if lbl not in labellist:
915        labellist.append(lbl)
916        return lbl
917    i = 1
918    prefix = lbl
919    if '_' in lbl:
920        prefix = lbl[:lbl.rfind('_')]
921        suffix = lbl[lbl.rfind('_')+1:]
922        try:
923            i = int(suffix)+1
924        except: # suffix could not be parsed
925            i = 1
926            prefix = lbl
927    while prefix+'_'+str(i) in labellist:
928        i += 1
929    else:
930        lbl = prefix+'_'+str(i)
931        labellist.append(lbl)
932    return lbl
933
934PhaseIdLookup = {}
935'''dict listing phase name and random Id keyed by sequential phase index as a str;
936best to access this using :func:`LookupPhaseName`
937'''
938PhaseRanIdLookup = {}
939'''dict listing phase sequential index keyed by phase random Id;
940best to access this using :func:`LookupPhaseId`
941'''
942HistIdLookup = {}
943'''dict listing histogram name and random Id, keyed by sequential histogram index as a str;
944best to access this using :func:`LookupHistName`
945'''
946HistRanIdLookup = {}
947'''dict listing histogram sequential index keyed by histogram random Id;
948best to access this using :func:`LookupHistId`
949'''
950AtomIdLookup = {}
951'''dict listing for each phase index as a str, the atom label and atom random Id,
952keyed by atom sequential index as a str;
953best to access this using :func:`LookupAtomLabel`
954'''
955AtomRanIdLookup = {}
956'''dict listing for each phase the atom sequential index keyed by atom random Id;
957best to access this using :func:`LookupAtomId`
958'''
959ShortPhaseNames = {}
960'''a dict containing a possibly shortened and when non-unique numbered
961version of the phase name. Keyed by the phase sequential index.
962'''
963ShortHistNames = {}
964'''a dict containing a possibly shortened and when non-unique numbered
965version of the histogram name. Keyed by the histogram sequential index.
966'''
967
968VarDesc = {}
969''' This dictionary lists descriptions for GSAS-II variables,
970as set in :func:`CompileVarDesc`. See that function for a description
971for how keys and values are written.
972'''
973
974reVarDesc = {}
975''' This dictionary lists descriptions for GSAS-II variables with
976the same values as :attr:`VarDesc` except that keys have been compiled as
977regular expressions. Initialized in :func:`CompileVarDesc`.
978'''
979
980def IndexAllIds(Histograms,Phases):
981    '''Scan through the used phases & histograms and create an index
982    to the random numbers of phases, histograms and atoms. While doing this,
983    confirm that assigned random numbers are unique -- just in case lightning
984    strikes twice in the same place.
985
986    Note: this code assumes that the atom random Id (ranId) is the last
987    element each atom record.
988
989    This is called in three places (only): :func:`GSASIIstrIO.GetUsedHistogramsAndPhases`
990    (which loads the histograms and phases from a GPX file),
991    :meth:`~GSASII.GSASII.GetUsedHistogramsAndPhasesfromTree`
992    (which loads the histograms and phases from the data tree.) and
993    :meth:`GSASIIconstrGUI.UpdateConstraints`
994    (which displays & edits the constraints in a GUI)
995
996    TODO: do we need a lookup for rigid body variables?
997    '''
998    # process phases and atoms
999    PhaseIdLookup.clear()
1000    PhaseRanIdLookup.clear()   
1001    AtomIdLookup.clear()
1002    AtomRanIdLookup.clear()
1003    ShortPhaseNames.clear()
1004    for ph in Phases:
1005        cx,ct,cs,cia = Phases[ph]['General']['AtomPtrs']
1006        ranId = Phases[ph]['ranId'] 
1007        while ranId in PhaseRanIdLookup:
1008            # Found duplicate random Id! note and reassign
1009            print ("\n\n*** Phase "+str(ph)+" has repeated ranId. Fixing.\n")
1010            Phases[ph]['ranId'] = ranId = ran.randint(0,sys.maxint)
1011        pId = str(Phases[ph]['pId'])
1012        PhaseIdLookup[pId] = (ph,ranId)
1013        PhaseRanIdLookup[ranId] = pId
1014        shortname = ph  #[:10]
1015        while shortname in ShortPhaseNames.values():
1016            shortname = ph[:8] + ' ('+ pId + ')'
1017        ShortPhaseNames[pId] = shortname
1018        AtomIdLookup[pId] = {}
1019        AtomRanIdLookup[pId] = {}
1020        for iatm,at in enumerate(Phases[ph]['Atoms']):
1021            ranId = at[cia+8]
1022            while ranId in AtomRanIdLookup[pId]: # check for dups
1023                print ("\n\n*** Phase "+str(ph)+" atom "+str(iatm)+" has repeated ranId. Fixing.\n")
1024                at[cia+8] = ranId = ran.randint(0,sys.maxint)
1025            AtomRanIdLookup[pId][ranId] = str(iatm)
1026            if Phases[ph]['General']['Type'] == 'macromolecular':
1027                label = '%s_%s_%s_%s'%(at[ct-1],at[ct-3],at[ct-4],at[ct-2])
1028            else:
1029                label = at[ct-1]
1030            AtomIdLookup[pId][str(iatm)] = (label,ranId)
1031    # process histograms
1032    HistIdLookup.clear()
1033    HistRanIdLookup.clear()
1034    ShortHistNames.clear()
1035    for hist in Histograms:
1036        ranId = Histograms[hist]['ranId']
1037        while ranId in HistRanIdLookup:
1038            # Found duplicate random Id! note and reassign
1039            print ("\n\n*** Histogram "+str(hist)+" has repeated ranId. Fixing.\n")
1040            Histograms[hist]['ranId'] = ranId = ran.randint(0,sys.maxint)
1041        hId = str(Histograms[hist]['hId'])
1042        HistIdLookup[hId] = (hist,ranId)
1043        HistRanIdLookup[ranId] = hId
1044        shortname = hist[:15]
1045        while shortname in ShortHistNames.values():
1046            shortname = hist[:11] + ' ('+ hId + ')'
1047        ShortHistNames[hId] = shortname
1048
1049def LookupAtomId(pId,ranId):
1050    '''Get the atom number from a phase and atom random Id
1051
1052    :param int/str pId: the sequential number of the phase
1053    :param int ranId: the random Id assigned to an atom
1054
1055    :returns: the index number of the atom (str)
1056    '''
1057    if not AtomRanIdLookup:
1058        raise Exception,'Error: LookupAtomId called before IndexAllIds was run'
1059    if pId is None or pId == '':
1060        raise KeyError,'Error: phase is invalid (None or blank)'
1061    pId = str(pId)
1062    if pId not in AtomRanIdLookup:
1063        raise KeyError,'Error: LookupAtomId does not have phase '+pId
1064    if ranId not in AtomRanIdLookup[pId]:
1065        raise KeyError,'Error: LookupAtomId, ranId '+str(ranId)+' not in AtomRanIdLookup['+pId+']'
1066    return AtomRanIdLookup[pId][ranId]
1067
1068def LookupAtomLabel(pId,index):
1069    '''Get the atom label from a phase and atom index number
1070
1071    :param int/str pId: the sequential number of the phase
1072    :param int index: the index of the atom in the list of atoms
1073
1074    :returns: the label for the atom (str) and the random Id of the atom (int)
1075    '''
1076    if not AtomIdLookup:
1077        raise Exception,'Error: LookupAtomLabel called before IndexAllIds was run'
1078    if pId is None or pId == '':
1079        raise KeyError,'Error: phase is invalid (None or blank)'
1080    pId = str(pId)
1081    if pId not in AtomIdLookup:
1082        raise KeyError,'Error: LookupAtomLabel does not have phase '+pId
1083    if index not in AtomIdLookup[pId]:
1084        raise KeyError,'Error: LookupAtomLabel, ranId '+str(index)+' not in AtomRanIdLookup['+pId+']'
1085    return AtomIdLookup[pId][index]
1086
1087def LookupPhaseId(ranId):
1088    '''Get the phase number and name from a phase random Id
1089
1090    :param int ranId: the random Id assigned to a phase
1091    :returns: the sequential Id (pId) number for the phase (str)
1092    '''
1093    if not PhaseRanIdLookup:
1094        raise Exception,'Error: LookupPhaseId called before IndexAllIds was run'
1095    if ranId not in PhaseRanIdLookup:
1096        raise KeyError,'Error: LookupPhaseId does not have ranId '+str(ranId)
1097    return PhaseRanIdLookup[ranId]
1098
1099def LookupPhaseName(pId):
1100    '''Get the phase number and name from a phase Id
1101
1102    :param int/str pId: the sequential assigned to a phase
1103    :returns:  (phase,ranId) where phase is the name of the phase (str)
1104      and ranId is the random # id for the phase (int)
1105    '''
1106    if not PhaseIdLookup:
1107        raise Exception,'Error: LookupPhaseName called before IndexAllIds was run'
1108    if pId is None or pId == '':
1109        raise KeyError,'Error: phase is invalid (None or blank)'
1110    pId = str(pId)
1111    if pId not in PhaseIdLookup:
1112        raise KeyError,'Error: LookupPhaseName does not have index '+pId
1113    return PhaseIdLookup[pId]
1114
1115def LookupHistId(ranId):
1116    '''Get the histogram number and name from a histogram random Id
1117
1118    :param int ranId: the random Id assigned to a histogram
1119    :returns: the sequential Id (hId) number for the histogram (str)
1120    '''
1121    if not HistRanIdLookup:
1122        raise Exception,'Error: LookupHistId called before IndexAllIds was run'
1123    if ranId not in HistRanIdLookup:
1124        raise KeyError,'Error: LookupHistId does not have ranId '+str(ranId)
1125    return HistRanIdLookup[ranId]
1126
1127def LookupHistName(hId):
1128    '''Get the histogram number and name from a histogram Id
1129
1130    :param int/str hId: the sequential assigned to a histogram
1131    :returns:  (hist,ranId) where hist is the name of the histogram (str)
1132      and ranId is the random # id for the histogram (int)
1133    '''
1134    if not HistIdLookup:
1135        raise Exception,'Error: LookupHistName called before IndexAllIds was run'
1136    if hId is None or hId == '':
1137        raise KeyError,'Error: histogram is invalid (None or blank)'
1138    hId = str(hId)
1139    if hId not in HistIdLookup:
1140        raise KeyError,'Error: LookupHistName does not have index '+hId
1141    return HistIdLookup[hId]
1142
1143def fmtVarDescr(varname):
1144    '''Return a string with a more complete description for a GSAS-II variable
1145
1146    :param str varname: A full G2 variable name with 2 or 3 or 4
1147       colons (<p>:<h>:name[:<a>] or <p>::RBname:<r>:<t>])
1148       
1149    :returns: a string with the description
1150    '''
1151    s,l = VarDescr(varname)
1152    return s+": "+l
1153
1154def VarDescr(varname):
1155    '''Return two strings with a more complete description for a GSAS-II variable
1156
1157    :param str name: A full G2 variable name with 2 or 3 or 4
1158       colons (<p>:<h>:name[:<a>] or <p>::RBname:<r>:<t>])
1159       
1160    :returns: (loc,meaning) where loc describes what item the variable is mapped
1161      (phase, histogram, etc.) and meaning describes what the variable does.
1162    '''
1163   
1164    # special handling for parameter names without a colons
1165    # for now, assume self-defining
1166    if varname.find(':') == -1:
1167        return "Global",varname
1168       
1169    l = getVarDescr(varname)
1170    if not l:
1171        return ("invalid variable name ("+str(varname)+")!"),""
1172#        return "invalid variable name!",""
1173
1174    if not l[-1]:
1175        l[-1] = "(variable needs a definition! Set it in CompileVarDesc)"
1176
1177    if len(l) == 3:         #SASD variable name!
1178        s = 'component:'+l[1]
1179        return s,l[-1]
1180    s = ""
1181    if l[0] is not None and l[1] is not None: # HAP: keep short
1182        if l[2] == "Scale": # fix up ambigous name
1183            l[5] = "Phase fraction"
1184        if l[0] == '*':
1185            lbl = 'Seq. ref.'
1186        else:
1187            lbl = ShortPhaseNames.get(l[0],'? #'+str(l[0]))
1188        if l[1] == '*':
1189            hlbl = 'Seq. ref.'
1190        else:
1191            hlbl = ShortHistNames.get(l[1],'? #'+str(l[1]))
1192        if hlbl[:4] == 'HKLF':
1193            hlbl = 'Xtl='+hlbl[5:]
1194        elif hlbl[:4] == 'PWDR':
1195            hlbl = 'Pwd='+hlbl[5:]
1196        else:
1197            hlbl = 'Hist='+hlbl
1198        s = "Ph="+str(lbl)+" * "+str(hlbl)
1199    else:
1200        if l[2] == "Scale": # fix up ambigous name: must be scale factor, since not HAP
1201            l[5] = "Scale factor"
1202        if l[2] == 'Back': # background parameters are "special", alas
1203            s = 'Hist='+ShortHistNames.get(l[1],'? #'+str(l[1]))
1204            l[-1] += ' #'+str(l[3])
1205        elif l[4] is not None: # rigid body parameter or modulation parm
1206            lbl = ShortPhaseNames.get(l[0],'phase?')
1207            if 'RB' in l[2]:    #rigid body parm
1208                s = "Res #"+str(l[3])+" body #"+str(l[4])+" in "+str(lbl)
1209            else: #modulation parm
1210                s = 'Atom %s wave %s in %s'%(LookupAtomLabel(l[0],l[3])[0],l[4],lbl)
1211        elif l[3] is not None: # atom parameter,
1212            lbl = ShortPhaseNames.get(l[0],'phase?')
1213            try:
1214                albl = LookupAtomLabel(l[0],l[3])[0]
1215            except KeyError:
1216                albl = 'Atom?'
1217            s = "Atom "+str(albl)+" in "+str(lbl)
1218        elif l[0] == '*':
1219            s = "All phases "
1220        elif l[0] is not None:
1221            lbl = ShortPhaseNames.get(l[0],'phase?')
1222            s = "Phase "+str(lbl)
1223        elif l[1] == '*':
1224            s = 'All hists'
1225        elif l[1] is not None:
1226            hlbl = ShortHistNames.get(l[1],'? #'+str(l[1]))
1227            if hlbl[:4] == 'HKLF':
1228                hlbl = 'Xtl='+hlbl[5:]
1229            elif hlbl[:4] == 'PWDR':
1230                hlbl = 'Pwd='+hlbl[5:]
1231            else:
1232                hlbl = 'Hist='+hlbl
1233            s = str(hlbl)
1234    if not s:
1235        s = 'Global'
1236    return s,l[-1]
1237
1238def getVarDescr(varname):
1239    '''Return a short description for a GSAS-II variable
1240
1241    :param str name: A full G2 variable name with 2 or 3 or 4
1242       colons (<p>:<h>:name[:<a1>][:<a2>])
1243     
1244    :returns: a six element list as [`p`,`h`,`name`,`a1`,`a2`,`description`],
1245      where `p`, `h`, `a1`, `a2` are str values or `None`, for the phase number,
1246      the histogram number and the atom number; `name` will always be
1247      a str; and `description` is str or `None`.
1248      If the variable name is incorrectly formed (for example, wrong
1249      number of colons), `None` is returned instead of a list.
1250    '''
1251    l = varname.split(':')
1252    if len(l) == 2:     #SASD parameter name
1253        return varname,l[0],getDescr(l[1])
1254    if len(l) == 3:
1255        l += [None,None]
1256    elif len(l) == 4:
1257        l += [None]
1258    elif len(l) != 5:
1259        return None
1260    for i in (0,1,3,4):
1261        if l[i] == "":
1262            l[i] = None
1263    l += [getDescr(l[2])]
1264    return l
1265   
1266def CompileVarDesc():
1267    '''Set the values in the variable description lookup table (:attr:`VarDesc`)
1268    into :attr:`reVarDesc`. This is called in :func:`getDescr` so the initialization
1269    is always done before use.
1270
1271    Note that keys may contain regular expressions, where '[xyz]'
1272    matches 'x' 'y' or 'z' (equivalently '[x-z]' describes this as range of values).
1273    '.*' matches any string. For example::
1274
1275    'AUiso':'Atomic isotropic displacement parameter',
1276
1277    will match variable ``'p::AUiso:a'``.
1278    If parentheses are used in the key, the contents of those parentheses can be
1279    used in the value, such as::
1280
1281    'AU([123][123])':'Atomic anisotropic displacement parameter U\\1',
1282
1283    will match ``AU11``, ``AU23``,.. and `U11`, `U23` etc will be displayed
1284    in the value when used.
1285   
1286    '''
1287    if reVarDesc: return # already done
1288    for key,value in {
1289        # derived or other sequential vars
1290        '([abc])$' : 'Lattice parameter, \\1, from Ai and Djk', # N.B. '$' prevents match if any characters follow
1291        u'\u03B1' : u'Lattice parameter, \u03B1, from Ai and Djk',
1292        u'\u03B2' : u'Lattice parameter, \u03B2, from Ai and Djk',
1293        u'\u03B3' : u'Lattice parameter, \u03B3, from Ai and Djk',
1294        # ambiguous, alas:
1295        'Scale' : 'Phase or Histogram scale factor',
1296        # Phase vars (p::<var>)
1297        'A([0-5])' : 'Reciprocal metric tensor component \\1',
1298        '[vV]ol' : 'Unit cell volume', # probably an error that both upper and lower case are used
1299        # Atom vars (p::<var>:a)
1300        'dA([xyz])$' : 'change to atomic coordinate, \\1',
1301        'A([xyz])$' : '\\1 fractional atomic coordinate',
1302        'AUiso':'Atomic isotropic displacement parameter',
1303        'AU([123][123])':'Atomic anisotropic displacement parameter U\\1',
1304        'Afrac': 'Atomic occupancy parameter',
1305        'Amul': 'Atomic site multiplicity value',
1306        'AM([xyz])$' : 'Atomic magnetic moment parameter, \\1',
1307        # Hist & Phase (HAP) vars (p:h:<var>)
1308        'Back': 'Background term',
1309        'BkPkint;(.*)':'Background peak #\\1 intensity',
1310        'BkPkpos;(.*)':'Background peak #\\1 position',
1311        'BkPksig;(.*)':'Background peak #\\1 Gaussian width',
1312        'BkPkgam;(.*)':'Background peak #\\1 Cauchy width',
1313        'Bab([AU])': 'Babinet solvent scattering coef. \\1',
1314        'D([123][123])' : 'Anisotropic strain coef. \\1',
1315        'Extinction' : 'Extinction coef.',
1316        'MD' : 'March-Dollase coef.',
1317        'Mustrain;.*' : 'Microstrain coef.',
1318        'Size;.*' : 'Crystallite size value',
1319        'eA$' : 'Cubic mustrain value',
1320        'Ep$' : 'Primary extinction',
1321        'Es$' : 'Secondary type II extinction',
1322        'Eg$' : 'Secondary type I extinction',
1323        'Flack' : 'Flack parameter',
1324        'TwinFr' : 'Twin fraction',
1325        #Histogram vars (:h:<var>)
1326        'Absorption' : 'Absorption coef.',
1327        'Displace([XY])' : 'Debye-Scherrer sample displacement \\1',
1328        'Lam' : 'Wavelength',
1329        'Polariz\.' : 'Polarization correction',
1330        'SH/L' : 'FCJ peak asymmetry correction',
1331        '([UVW])$' : 'Gaussian instrument broadening \\1',
1332        '([XY])$' : 'Cauchy instrument broadening \\1',
1333        'Zero' : 'Debye-Scherrer zero correction',
1334        'nDebye' : 'Debye model background corr. terms',
1335        'nPeaks' : 'Fixed peak background corr. terms',
1336        'RBV.*' : 'Vector rigid body parameter',
1337        'RBR.*' : 'Residue rigid body parameter',
1338        'RBRO([aijk])' : 'Residue rigid body orientation parameter',
1339        'RBRP([xyz])' : 'Residue rigid body position parameter',
1340        'RBRTr;.*' : 'Residue rigid body torsion parameter',
1341        'RBR([TLS])([123AB][123AB])' : 'Residue rigid body group disp. param.',
1342        'constr([0-9]*)' : 'Parameter from constraint',
1343        # supersymmetry parameters  p::<var>:a:o 'Flen','Fcent'?
1344        'mV([0-2])$' : 'Modulation vector component \\1',
1345        'Fsin'  :   'Sin site fraction modulation',
1346        'Fcos'  :   'Cos site fraction modulation',
1347        'Fzero'  :   'Crenel function offset',      #may go away
1348        'Fwid'   :   'Crenel function width',
1349        'Tmin'   :   'ZigZag/Block min location',
1350        'Tmax'   :   'ZigZag/Block max location',
1351        '([XYZ])max': 'ZigZag/Block max value for \\1',
1352        '([XYZ])sin'  : 'Sin position wave for \\1',
1353        '([XYZ])cos'  : 'Cos position wave for \\1',
1354        'U([123][123])sin$' :  'Sin thermal wave for U\\1',
1355        'U([123][123])cos$' :  'Cos thermal wave for U\\1',
1356        'M([XYZ])sin$' :  'Sin mag. moment wave for \\1',
1357        'M([XYZ])cos$' :  'Cos mag. moment wave for \\1',
1358        # SASD vars (l:<var>;l = component)
1359        'Aspect ratio' : 'Particle aspect ratio',
1360        'Length' : 'Cylinder length',
1361        'Diameter' : 'Cylinder/disk diameter',
1362        'Thickness' : 'Disk thickness',
1363        'Shell thickness' : 'Multiplier to get inner(<1) or outer(>1) sphere radius',
1364        'Dist' : 'Interparticle distance',
1365        'VolFr' : 'Dense scatterer volume fraction',
1366        'epis' : 'Sticky sphere epsilon',
1367        'Sticky' : 'Stickyness',
1368        'Depth' : 'Well depth',
1369        'Width' : 'Well width',
1370        'Volume' : 'Particle volume',
1371        'Radius' : 'Sphere/cylinder/disk radius',
1372        'Mean' : 'Particle mean radius',
1373        'StdDev' : 'Standard deviation in Mean',
1374        'G$': 'Guinier prefactor',
1375        'Rg$': 'Guinier radius of gyration',
1376        'B$': 'Porod prefactor',
1377        'P$': 'Porod power',
1378        'Cutoff': 'Porod cutoff',
1379        'PkInt': 'Bragg peak intensity',
1380        'PkPos': 'Bragg peak position',
1381        'PkSig': 'Bragg peak sigma',
1382        'PkGam': 'Bragg peak gamma',
1383        'e([12][12])' : 'strain tensor e\1',   # strain vars e11, e22, e12
1384        'Dcalc': 'Calc. d-spacing',
1385        'Back$': 'background parameter',
1386        'pos$': 'peak position',
1387        'int$': 'peak intensity',
1388        'WgtFrac':'phase weight fraction',
1389        }.items():
1390        VarDesc[key] = value
1391        reVarDesc[re.compile(key)] = value
1392
1393def getDescr(name):
1394    '''Return a short description for a GSAS-II variable
1395
1396    :param str name: The descriptive part of the variable name without colons (:)
1397     
1398    :returns: a short description or None if not found
1399    '''
1400
1401    CompileVarDesc() # compile the regular expressions, if needed
1402    for key in reVarDesc:
1403        m = key.match(name)
1404        if m:
1405            reVarDesc[key]
1406            return m.expand(reVarDesc[key])
1407    return None
1408
1409def GenWildCard(varlist):
1410    '''Generate wildcard versions of G2 variables. These introduce '*'
1411    for a phase, histogram or atom number (but only for one of these
1412    fields) but only when there is more than one matching variable in the
1413    input variable list. So if the input is this::
1414   
1415      varlist = ['0::AUiso:0', '0::AUiso:1', '1::AUiso:0']
1416
1417    then the output will be this::
1418   
1419       wildList = ['*::AUiso:0', '0::AUiso:*']
1420
1421    :param list varlist: an input list of GSAS-II variable names
1422      (such as 0::AUiso:0)
1423
1424    :returns: wildList, the generated list of wild card variable names.
1425    '''
1426    wild = []
1427    for i in (0,1,3):
1428        currentL = varlist[:]
1429        while currentL:
1430            item1 = currentL.pop(0)
1431            i1splt = item1.split(':')
1432            if i >= len(i1splt): continue
1433            if i1splt[i]:
1434                nextL = []
1435                i1splt[i] = '[0-9]+'
1436                rexp = re.compile(':'.join(i1splt))
1437                matchlist = [item1]
1438                for nxtitem in currentL:
1439                    if rexp.match(nxtitem):
1440                        matchlist += [nxtitem]
1441                    else:
1442                        nextL.append(nxtitem)
1443                if len(matchlist) > 1:
1444                    i1splt[i] = '*'
1445                    wild.append(':'.join(i1splt))
1446                currentL = nextL
1447    return wild
1448
1449def LookupWildCard(varname,varlist):
1450    '''returns a list of variable names from list varname
1451    that match wildcard name in varname
1452   
1453    :param str varname: a G2 variable name containing a wildcard
1454      (such as \*::var)
1455    :param list varlist: the list of all variable names used in
1456      the current project
1457    :returns: a list of matching GSAS-II variables (may be empty) 
1458    '''
1459    rexp = re.compile(varname.replace('*','[0-9]+'))
1460    return sorted([var for var in varlist if rexp.match(var)])
1461
1462
1463def _lookup(dic,key):
1464    '''Lookup a key in a dictionary, where None returns an empty string
1465    but an unmatched key returns a question mark. Used in :class:`G2VarObj`
1466    '''
1467    if key is None:
1468        return ""
1469    elif key == "*":
1470        return "*"
1471    else:
1472        return dic.get(key,'?')
1473
1474class G2VarObj(object):
1475    '''Defines a GSAS-II variable either using the phase/atom/histogram
1476    unique Id numbers or using a character string that specifies
1477    variables by phase/atom/histogram number (which can change).
1478    Note that :func:`LoadID` should be used to (re)load the current Ids
1479    before creating or later using the G2VarObj object.
1480
1481    This can store rigid body variables, but does not translate the residue # and
1482    body # to/from random Ids
1483
1484    A :class:`G2VarObj` object can be created with a single parameter:
1485   
1486    :param str/tuple varname: a single value can be used to create a :class:`G2VarObj`
1487      object. If a string, it must be of form "p:h:var" or "p:h:var:a", where
1488
1489     * p is the phase number (which may be left blank or may be '*' to indicate all phases);
1490     * h is the histogram number (which may be left blank or may be '*' to indicate all histograms);
1491     * a is the atom number (which may be left blank in which case the third colon is omitted).
1492       The atom number can be specified as '*' if a phase number is specified (not as '*').
1493       For rigid body variables, specify a will be a string of form "residue:body#"
1494
1495      Alternately a single tuple of form (Phase,Histogram,VarName,AtomID) can be used, where
1496      Phase, Histogram, and AtomID are None or are ranId values (or one can be '*')
1497      and VarName is a string. Note that if Phase is '*' then the AtomID is an atom number.
1498      For a rigid body variables, AtomID is a string of form "residue:body#".
1499
1500    If four positional arguments are supplied, they are:
1501
1502    :param str/int phasenum: The number for the phase (or None or '*')
1503    :param str/int histnum: The number for the histogram (or None or '*')
1504    :param str varname: a single value can be used to create a :class:`G2VarObj`
1505    :param str/int atomnum: The number for the atom (or None or '*')
1506   
1507    '''
1508    IDdict = {}
1509    IDdict['phases'] = {}
1510    IDdict['hists'] = {}
1511    IDdict['atoms'] = {}
1512    def __init__(self,*args):
1513        self.phase = None
1514        self.histogram = None
1515        self.name = ''
1516        self.atom = None
1517        if len(args) == 1 and (type(args[0]) is list or type(args[0]) is tuple) and len(args[0]) == 4:
1518            # single arg with 4 values
1519            self.phase,self.histogram,self.name,self.atom = args[0]
1520        elif len(args) == 1 and ':' in args[0]:
1521            #parse a string
1522            lst = args[0].split(':')
1523            if lst[0] == '*':
1524                self.phase = '*'
1525                if len(lst) > 3:
1526                    self.atom = lst[3]
1527                self.histogram = HistIdLookup.get(lst[1],[None,None])[1]
1528            elif lst[1] == '*':           
1529                self.histogram = '*'
1530                self.phase = PhaseIdLookup.get(lst[0],[None,None])[1]
1531            else:
1532                self.histogram = HistIdLookup.get(lst[1],[None,None])[1]
1533                self.phase = PhaseIdLookup.get(lst[0],[None,None])[1]
1534                if len(lst) == 4:
1535                    if lst[3] == '*':
1536                        self.atom = '*'
1537                    else:
1538                        self.atom = AtomIdLookup[lst[0]].get(lst[3],[None,None])[1]
1539                elif len(lst) == 5:
1540                    self.atom = lst[3]+":"+lst[4]
1541                elif len(lst) == 3:
1542                    pass
1543                else:
1544                    raise Exception,"Too many colons in var name "+str(args[0])
1545            self.name = lst[2]
1546        elif len(args) == 4:
1547            if args[0] == '*':
1548                self.phase = '*'
1549                self.atom = args[3]
1550            else:
1551                self.phase = PhaseIdLookup.get(str(args[0]),[None,None])[1]
1552                if args[3] == '*':
1553                    self.atom = '*'
1554                elif args[0] is not None:
1555                    self.atom = AtomIdLookup[args[0]].get(str(args[3]),[None,None])[1]
1556            if args[1] == '*':
1557                self.histogram = '*'
1558            else:
1559                self.histogram = HistIdLookup.get(str(args[1]),[None,None])[1]
1560            self.name = args[2]
1561        else:
1562            raise Exception,"Incorrectly called GSAS-II parameter name"
1563
1564        #print "DEBUG: created ",self.phase,self.histogram,self.name,self.atom
1565
1566    def __str__(self):
1567        return self.varname()
1568
1569    def varname(self):
1570        '''Formats the GSAS-II variable name as a "traditional" GSAS-II variable
1571        string (p:h:<var>:a) or (p:h:<var>)
1572
1573        :returns: the variable name as a str
1574        '''
1575        a = ""
1576        if self.phase == "*":
1577            ph = "*"
1578            if self.atom:
1579                a = ":" + str(self.atom)
1580        else:
1581            ph = _lookup(PhaseRanIdLookup,self.phase)
1582            if self.atom == '*':
1583                a = ':*'
1584            elif self.atom:
1585                if ":" in str(self.atom):
1586                    a = ":" + str(self.atom)
1587                elif ph in AtomRanIdLookup:
1588                    a = ":" + AtomRanIdLookup[ph].get(self.atom,'?')
1589                else:
1590                    a = ":?"
1591        if self.histogram == "*":
1592            hist = "*"
1593        else:
1594            hist = _lookup(HistRanIdLookup,self.histogram)
1595        s = (ph + ":" + hist + ":" + str(self.name)) + a
1596        return s
1597   
1598    def __repr__(self):
1599        '''Return the detailed contents of the object
1600        '''
1601        s = "<"
1602        if self.phase == '*':
1603            s += "Phases: all; "
1604            if self.atom is not None:
1605                if ":" in str(self.atom):
1606                    s += "Rigid body" + str(self.atom) + "; "
1607                else:
1608                    s += "Atom #" + str(self.atom) + "; "
1609        elif self.phase is not None:
1610            ph =  _lookup(PhaseRanIdLookup,self.phase)
1611            s += "Phase: rId=" + str(self.phase) + " (#"+ ph + "); "
1612            if self.atom == '*':
1613                s += "Atoms: all; "
1614            elif ":" in str(self.atom):
1615                s += "Rigid body" + str(self.atom) + "; "
1616            elif self.atom is not None:
1617                s += "Atom rId=" + str(self.atom)
1618                if ph in AtomRanIdLookup:
1619                    s += " (#" + AtomRanIdLookup[ph].get(self.atom,'?') + "); "
1620                else:
1621                    s += " (#? -- not found!); "
1622        if self.histogram == '*':
1623            s += "Histograms: all; "
1624        elif self.histogram is not None:
1625            hist = _lookup(HistRanIdLookup,self.histogram)
1626            s += "Histogram: rId=" + str(self.histogram) + " (#"+ hist + "); "
1627        s += 'Variable name="' + str(self.name) + '">'
1628        return s+" ("+self.varname()+")"
1629
1630    def __eq__(self, other):
1631        if type(other) is type(self):
1632            return (self.phase == other.phase and
1633                    self.histogram == other.histogram and
1634                    self.name == other.name and
1635                    self.atom == other.atom)
1636        return False
1637
1638    def _show(self):
1639        'For testing, shows the current lookup table'
1640        print 'phases', self.IDdict['phases']
1641        print 'hists', self.IDdict['hists']
1642        print 'atomDict', self.IDdict['atoms']
1643
1644#==========================================================================
1645# shortcut routines
1646exp = np.exp
1647sind = sin = s = lambda x: np.sin(x*np.pi/180.)
1648cosd = cos = c = lambda x: np.cos(x*np.pi/180.)
1649tand = tan = t = lambda x: np.tan(x*np.pi/180.)
1650sqrt = sq = lambda x: np.sqrt(x)
1651pi = lambda: np.pi
1652class ExpressionObj(object):
1653    '''Defines an object with a user-defined expression, to be used for
1654    secondary fits or restraints. Object is created null, but is changed
1655    using :meth:`LoadExpression`. This contains only the minimum
1656    information that needs to be stored to save and load the expression
1657    and how it is mapped to GSAS-II variables.
1658    '''
1659    def __init__(self):
1660        self.expression = ''
1661        'The expression as a text string'
1662        self.assgnVars = {}
1663        '''A dict where keys are label names in the expression mapping to a GSAS-II
1664        variable. The value a G2 variable name.
1665        Note that the G2 variable name may contain a wild-card and correspond to
1666        multiple values.
1667        '''
1668        self.freeVars = {}
1669        '''A dict where keys are label names in the expression mapping to a free
1670        parameter. The value is a list with:
1671
1672         * a name assigned to the parameter
1673         * a value for to the parameter and
1674         * a flag to determine if the variable is refined.
1675        ''' 
1676        self.depVar = None
1677
1678        self.lastError = ('','')
1679        '''Shows last encountered error in processing expression
1680        (list of 1-3 str values)'''
1681
1682        self.distance_dict  = None  # to be used for defining atom phase/symmetry info
1683        self.distance_atoms = None  # to be used for defining atom distances
1684
1685    def LoadExpression(self,expr,exprVarLst,varSelect,varName,varValue,varRefflag):
1686        '''Load the expression and associated settings into the object. Raises
1687        an exception if the expression is not parsed, if not all functions
1688        are defined or if not all needed parameter labels in the expression
1689        are defined.
1690
1691        This will not test if the variable referenced in these definitions
1692        are actually in the parameter dictionary. This is checked when the
1693        computation for the expression is done in :meth:`SetupCalc`.
1694       
1695        :param str expr: the expression
1696        :param list exprVarLst: parameter labels found in the expression
1697        :param dict varSelect: this will be 0 for Free parameters
1698          and non-zero for expression labels linked to G2 variables.
1699        :param dict varName: Defines a name (str) associated with each free parameter
1700        :param dict varValue: Defines a value (float) associated with each free parameter
1701        :param dict varRefflag: Defines a refinement flag (bool)
1702          associated with each free parameter
1703        '''
1704        self.expression = expr
1705        self.compiledExpr = None
1706        self.freeVars = {}
1707        self.assgnVars = {}
1708        for v in exprVarLst:
1709            if varSelect[v] == 0:
1710                self.freeVars[v] = [
1711                    varName.get(v),
1712                    varValue.get(v),
1713                    varRefflag.get(v),
1714                    ]
1715            else:
1716                self.assgnVars[v] = varName[v]
1717        self.CheckVars()
1718
1719    def EditExpression(self,exprVarLst,varSelect,varName,varValue,varRefflag):
1720        '''Load the expression and associated settings from the object into
1721        arrays used for editing.
1722
1723        :param list exprVarLst: parameter labels found in the expression
1724        :param dict varSelect: this will be 0 for Free parameters
1725          and non-zero for expression labels linked to G2 variables.
1726        :param dict varName: Defines a name (str) associated with each free parameter
1727        :param dict varValue: Defines a value (float) associated with each free parameter
1728        :param dict varRefflag: Defines a refinement flag (bool)
1729          associated with each free parameter
1730
1731        :returns: the expression as a str
1732        '''
1733        for v in self.freeVars:
1734            varSelect[v] = 0
1735            varName[v] = self.freeVars[v][0]
1736            varValue[v] = self.freeVars[v][1]
1737            varRefflag[v] = self.freeVars[v][2]
1738        for v in self.assgnVars:
1739            varSelect[v] = 1
1740            varName[v] = self.assgnVars[v]
1741        return self.expression
1742
1743    def GetVaried(self):
1744        'Returns the names of the free parameters that will be refined'
1745        return ["::"+self.freeVars[v][0] for v in self.freeVars if self.freeVars[v][2]]
1746
1747    def GetVariedVarVal(self):
1748        'Returns the names and values of the free parameters that will be refined'
1749        return [("::"+self.freeVars[v][0],self.freeVars[v][1]) for v in self.freeVars if self.freeVars[v][2]]
1750
1751    def UpdateVariedVars(self,varyList,values):
1752        'Updates values for the free parameters (after a refinement); only updates refined vars'
1753        for v in self.freeVars:
1754            if not self.freeVars[v][2]: continue
1755            if "::"+self.freeVars[v][0] not in varyList: continue
1756            indx = varyList.index("::"+self.freeVars[v][0])
1757            self.freeVars[v][1] = values[indx]
1758
1759    def GetIndependentVars(self):
1760        'Returns the names of the required independent parameters used in expression'
1761        return [self.assgnVars[v] for v in self.assgnVars]
1762
1763    def CheckVars(self):
1764        '''Check that the expression can be parsed, all functions are
1765        defined and that input loaded into the object is internally
1766        consistent. If not an Exception is raised.
1767
1768        :returns: a dict with references to packages needed to
1769          find functions referenced in the expression.
1770        '''
1771        ret = self.ParseExpression(self.expression)
1772        if not ret:
1773            raise Exception("Expression parse error")
1774        exprLblList,fxnpkgdict = ret
1775        # check each var used in expression is defined
1776        defined = self.assgnVars.keys() + self.freeVars.keys()
1777        notfound = []
1778        for var in exprLblList:
1779            if var not in defined:
1780                notfound.append(var)
1781        if notfound:
1782            msg = 'Not all variables defined'
1783            msg1 = 'The following variables were not defined: '
1784            msg2 = ''
1785            for var in notfound:
1786                if msg: msg += ', '
1787                msg += var
1788            self.lastError = (msg1,'  '+msg2)
1789            raise Exception(msg)
1790        return fxnpkgdict
1791
1792    def ParseExpression(self,expr):
1793        '''Parse an expression and return a dict of called functions and
1794        the variables used in the expression. Returns None in case an error
1795        is encountered. If packages are referenced in functions, they are loaded
1796        and the functions are looked up into the modules global
1797        workspace.
1798       
1799        Note that no changes are made to the object other than
1800        saving an error message, so that this can be used for testing prior
1801        to the save.
1802
1803        :returns: a list of used variables
1804        '''
1805        self.lastError = ('','')
1806        import ast
1807        def FindFunction(f):
1808            '''Find the object corresponding to function f
1809            :param str f: a function name such as 'numpy.exp'
1810            :returns: (pkgdict,pkgobj) where pkgdict contains a dict
1811              that defines the package location(s) and where pkgobj
1812              defines the object associated with the function.
1813              If the function is not found, pkgobj is None.
1814            '''
1815            df = f.split('.')
1816            pkgdict = {}
1817            # no listed package, try in current namespace
1818            if len(df) == 1: 
1819                try:
1820                    fxnobj = eval(f)
1821                    return pkgdict,fxnobj
1822                except (AttributeError, NameError):
1823                    return None,None
1824            else:
1825                try:
1826                    fxnobj = eval(f)
1827                    pkgdict[df[0]] = eval(df[0])
1828                    return pkgdict,fxnobj
1829                except (AttributeError, NameError):
1830                    pass
1831            # includes a package, lets try to load the packages
1832            pkgname = ''
1833            path = sys.path
1834            for pkg in f.split('.')[:-1]: # if needed, descend down the tree
1835                if pkgname:
1836                    pkgname += '.' + pkg
1837                else:
1838                    pkgname = pkg
1839                fp = None
1840                try:
1841                    fp, fppath,desc = imp.find_module(pkg,path)
1842                    pkgobj = imp.load_module(pkg,fp,fppath,desc)
1843                    pkgdict[pkgname] = pkgobj
1844                    path = [fppath]
1845                except Exception as msg:
1846                    print('load of '+pkgname+' failed with error='+str(msg))
1847                    return {},None
1848                finally:
1849                    if fp: fp.close()
1850                try:
1851                    #print 'before',pkgdict.keys()
1852                    fxnobj = eval(f,globals(),pkgdict)
1853                    #print 'after 1',pkgdict.keys()
1854                    #fxnobj = eval(f,pkgdict)
1855                    #print 'after 2',pkgdict.keys()
1856                    return pkgdict,fxnobj
1857                except:
1858                    continue
1859            return None # not found
1860        def ASTtransverse(node,fxn=False):
1861            '''Transverse a AST-parsed expresson, compiling a list of variables
1862            referenced in the expression. This routine is used recursively.
1863
1864            :returns: varlist,fxnlist where
1865              varlist is a list of referenced variable names and
1866              fxnlist is a list of used functions
1867            '''
1868            varlist = []
1869            fxnlist = []
1870            if isinstance(node, list):
1871                for b in node:
1872                    v,f = ASTtransverse(b,fxn)
1873                    varlist += v
1874                    fxnlist += f
1875            elif isinstance(node, ast.AST):
1876                for a, b in ast.iter_fields(node):
1877                    if isinstance(b, ast.AST):
1878                        if a == 'func': 
1879                            fxnlist += ['.'.join(ASTtransverse(b,True)[0])]
1880                            continue
1881                        v,f = ASTtransverse(b,fxn)
1882                        varlist += v
1883                        fxnlist += f
1884                    elif isinstance(b, list):
1885                        v,f = ASTtransverse(b,fxn)
1886                        varlist += v
1887                        fxnlist += f
1888                    elif node.__class__.__name__ == "Name":
1889                        varlist += [b]
1890                    elif fxn and node.__class__.__name__ == "Attribute":
1891                        varlist += [b]
1892            return varlist,fxnlist
1893        try:
1894            exprast = ast.parse(expr)
1895        except SyntaxError:
1896            s = ''
1897            import traceback
1898            for i in traceback.format_exc().splitlines()[-3:-1]:
1899                if s: s += "\n"
1900                s += str(i)
1901            self.lastError = ("Error parsing expression:",s)
1902            return
1903        # find the variables & functions
1904        v,f = ASTtransverse(exprast)
1905        varlist = sorted(list(set(v)))
1906        fxnlist = list(set(f))
1907        pkgdict = {}
1908        # check the functions are defined
1909        for fxn in fxnlist:
1910            fxndict,fxnobj = FindFunction(fxn)
1911            if not fxnobj:
1912                self.lastError = ("Error: Invalid function",fxn,
1913                                  "is not defined")
1914                return
1915            if not hasattr(fxnobj,'__call__'):
1916                self.lastError = ("Error: Not a function.",fxn,
1917                                  "cannot be called as a function")
1918                return
1919            pkgdict.update(fxndict)
1920        return varlist,pkgdict
1921
1922    def GetDepVar(self):
1923        'return the dependent variable, or None'
1924        return self.depVar
1925
1926    def SetDepVar(self,var):
1927        'Set the dependent variable, if used'
1928        self.depVar = var
1929#==========================================================================
1930class ExpressionCalcObj(object):
1931    '''An object used to evaluate an expression from a :class:`ExpressionObj`
1932    object.
1933   
1934    :param ExpressionObj exprObj: a :class:`~ExpressionObj` expression object with
1935      an expression string and mappings for the parameter labels in that object.
1936    '''
1937    def __init__(self,exprObj):
1938        self.eObj = exprObj
1939        'The expression and mappings; a :class:`ExpressionObj` object'
1940        self.compiledExpr = None
1941        'The expression as compiled byte-code'
1942        self.exprDict = {}
1943        '''dict that defines values for labels used in expression and packages
1944        referenced by functions
1945        '''
1946        self.lblLookup = {}
1947        '''Lookup table that specifies the expression label name that is
1948        tied to a particular GSAS-II parameters in the parmDict.
1949        '''
1950        self.fxnpkgdict = {}
1951        '''a dict with references to packages needed to
1952        find functions referenced in the expression.
1953        '''
1954        self.varLookup = {}
1955        '''Lookup table that specifies the GSAS-II variable(s)
1956        indexed by the expression label name. (Used for only for diagnostics
1957        not evaluation of expression.)
1958        '''
1959        self.su = None
1960        '''Standard error evaluation where supplied by the evaluator
1961        '''
1962        # Patch: for old-style expressions with a (now removed step size)
1963        for v in self.eObj.assgnVars:
1964            if not isinstance(self.eObj.assgnVars[v], basestring):
1965                self.eObj.assgnVars[v] = self.eObj.assgnVars[v][0]
1966        self.parmDict = {}
1967        '''A copy of the parameter dictionary, for distance and angle computation
1968        '''
1969
1970    def SetupCalc(self,parmDict):
1971        '''Do all preparations to use the expression for computation.
1972        Adds the free parameter values to the parameter dict (parmDict).
1973        '''
1974        if self.eObj.expression.startswith('Dist') or self.eObj.expression.startswith('Angle'):
1975            return
1976        self.fxnpkgdict = self.eObj.CheckVars()
1977        # all is OK, compile the expression
1978        self.compiledExpr = compile(self.eObj.expression,'','eval')
1979
1980        # look at first value in parmDict to determine its type
1981        parmsInList = True
1982        for key in parmDict:
1983            val = parmDict[key]
1984            if isinstance(val, basestring):
1985                parmsInList = False
1986                break
1987            try: # check if values are in lists
1988                val = parmDict[key][0]
1989            except (TypeError,IndexError):
1990                parmsInList = False
1991            break
1992           
1993        # set up the dicts needed to speed computations
1994        self.exprDict = {}
1995        self.lblLookup = {}
1996        self.varLookup = {}
1997        for v in self.eObj.freeVars:
1998            varname = self.eObj.freeVars[v][0]
1999            varname = "::" + varname.lstrip(':').replace(' ','_').replace(':',';')
2000            self.lblLookup[varname] = v
2001            self.varLookup[v] = varname
2002            if parmsInList:
2003                parmDict[varname] = [self.eObj.freeVars[v][1],self.eObj.freeVars[v][2]]
2004            else:
2005                parmDict[varname] = self.eObj.freeVars[v][1]
2006            self.exprDict[v] = self.eObj.freeVars[v][1]
2007        for v in self.eObj.assgnVars:
2008            varname = self.eObj.assgnVars[v]
2009            if '*' in varname:
2010                varlist = LookupWildCard(varname,parmDict.keys())
2011                if len(varlist) == 0:
2012                    raise Exception,"No variables match "+str(v)
2013                for var in varlist:
2014                    self.lblLookup[var] = v
2015                if parmsInList:
2016                    self.exprDict[v] = np.array([parmDict[var][0] for var in varlist])
2017                else:
2018                    self.exprDict[v] = np.array([parmDict[var] for var in varlist])
2019                self.varLookup[v] = [var for var in varlist]
2020            elif varname in parmDict:
2021                self.lblLookup[varname] = v
2022                self.varLookup[v] = varname
2023                if parmsInList:
2024                    self.exprDict[v] = parmDict[varname][0]
2025                else:
2026                    self.exprDict[v] = parmDict[varname]
2027            else:
2028                raise Exception,"No value for variable "+str(v)
2029        self.exprDict.update(self.fxnpkgdict)
2030
2031    def UpdateVars(self,varList,valList):
2032        '''Update the dict for the expression with a set of values
2033        :param list varList: a list of variable names
2034        :param list valList: a list of corresponding values
2035        '''
2036        for var,val in zip(varList,valList):
2037            self.exprDict[self.lblLookup.get(var,'undefined: '+var)] = val
2038
2039    def UpdateDict(self,parmDict):
2040        '''Update the dict for the expression with values in a dict
2041        :param list parmDict: a dict of values some of which may be in use here
2042        '''
2043        if self.eObj.expression.startswith('Dist') or self.eObj.expression.startswith('Angle'):
2044            self.parmDict = parmDict
2045            return
2046        for var in parmDict:
2047            if var in self.lblLookup:
2048                self.exprDict[self.lblLookup[var]] = parmDict[var]
2049           
2050    def EvalExpression(self):
2051        '''Evaluate an expression. Note that the expression
2052        and mapping are taken from the :class:`ExpressionObj` expression object
2053        and the parameter values were specified in :meth:`SetupCalc`.
2054        :returns: a single value for the expression. If parameter
2055        values are arrays (for example, from wild-carded variable names),
2056        the sum of the resulting expression is returned.
2057
2058        For example, if the expression is ``'A*B'``,
2059        where A is 2.0 and B maps to ``'1::Afrac:*'``, which evaluates to::
2060
2061        [0.5, 1, 0.5]
2062
2063        then the result will be ``4.0``.
2064        '''
2065        self.su = None
2066        if self.eObj.expression.startswith('Dist'):
2067#            GSASIIpath.IPyBreak()
2068            dist = G2mth.CalcDist(self.eObj.distance_dict, self.eObj.distance_atoms, self.parmDict)
2069            return dist
2070        elif self.eObj.expression.startswith('Angle'):
2071            angle = G2mth.CalcAngle(self.eObj.angle_dict, self.eObj.angle_atoms, self.parmDict)
2072            return angle
2073        if self.compiledExpr is None:
2074            raise Exception,"EvalExpression called before SetupCalc"
2075        val = eval(self.compiledExpr,globals(),self.exprDict)
2076        if not np.isscalar(val):
2077            val = np.sum(val)
2078        return val
2079       
2080class G2Exception(Exception):
2081    def __init__(self,msg):
2082        self.msg = msg
2083    def __str__(self):
2084        return repr(self.msg)
2085
2086def CreatePDFitems(G2frame,PWDRtree,ElList,Qlimits,PDFnames=[]):
2087    '''Create and initialize a new set of PDF tree entries
2088
2089    :param Frame G2frame: main GSAS-II tree frame object
2090    :param str PWDRtree: name of PWDR to be used to create PDF item
2091    :param dict ElList: data structure with composition
2092    :param list Qlimits: Q limits to be used for computing the PDF
2093    :returns: the Id of the newly created PDF entry
2094    '''
2095    PDFname = 'PDF '+PWDRtree[4:] # this places two spaces after PDF which is needed is some places
2096    PDFname = MakeUniqueLabel(PDFname,PDFnames)
2097    Id = G2frame.PatternTree.AppendItem(parent=G2frame.root,text=PDFname)
2098    Data = {
2099        'Sample':{'Name':PWDRtree,'Mult':1.0},
2100        'Sample Bkg.':{'Name':'','Mult':-1.0},
2101        'Container':{'Name':'','Mult':-1.0},
2102        'Container Bkg.':{'Name':'','Mult':-1.0,'Add':0.0},'ElList':ElList,
2103        'Geometry':'Cylinder','Diam':1.0,'Pack':0.50,'Form Vol':10.0,
2104        'DetType':'Image plate','ObliqCoeff':0.2,'Ruland':0.025,'QScaleLim':Qlimits,
2105        'Lorch':False,'BackRatio':0.0,'Rmax':100.,'noRing':False,'IofQmin':1.0,
2106        'I(Q)':[],'S(Q)':[],'F(Q)':[],'G(R)':[]}
2107    G2frame.PatternTree.SetItemPyData(G2frame.PatternTree.AppendItem(Id,text='PDF Controls'),Data)
2108    G2frame.PatternTree.SetItemPyData(G2frame.PatternTree.AppendItem(Id,text='PDF Peaks'),
2109        {'Limits':[1.,5.],'Background':[2,[0.,-0.2*np.pi],False],'Peaks':[]})
2110    return Id       
2111
2112
2113if __name__ == "__main__":
2114    # test equation evaluation
2115    def showEQ(calcobj):
2116        print 50*'='
2117        print calcobj.eObj.expression,'=',calcobj.EvalExpression()
2118        for v in sorted(calcobj.varLookup):
2119            print "  ",v,'=',calcobj.exprDict[v],'=',calcobj.varLookup[v]
2120        # print '  Derivatives'
2121        # for v in calcobj.derivStep.keys():
2122        #     print '    d(Expr)/d('+v+') =',calcobj.EvalDeriv(v)
2123
2124    obj = ExpressionObj()
2125
2126    obj.expression = "A*np.exp(B)"
2127    obj.assgnVars =  {'B': '0::Afrac:1'}
2128    obj.freeVars =  {'A': [u'A', 0.5, True]}
2129    #obj.CheckVars()
2130    calcobj = ExpressionCalcObj(obj)
2131
2132    obj1 = ExpressionObj()
2133    obj1.expression = "A*np.exp(B)"
2134    obj1.assgnVars =  {'B': '0::Afrac:*'}
2135    obj1.freeVars =  {'A': [u'Free Prm A', 0.5, True]}
2136    #obj.CheckVars()
2137    calcobj1 = ExpressionCalcObj(obj1)
2138
2139    obj2 = ExpressionObj()
2140    obj2.distance_stuff = np.array([[0,1],[1,-1]])
2141    obj2.expression = "Dist(1,2)"
2142    GSASIIpath.InvokeDebugOpts()
2143    parmDict2 = {'0::Afrac:0':[0.0,True], '0::Afrac:1': [1.0,False]}
2144    calcobj2 = ExpressionCalcObj(obj2)
2145    calcobj2.SetupCalc(parmDict2)
2146    showEQ(calcobj2)
2147   
2148    parmDict1 = {'0::Afrac:0':1.0, '0::Afrac:1': 1.0}
2149    print '\nDict = ',parmDict1
2150    calcobj.SetupCalc(parmDict1)
2151    showEQ(calcobj)
2152    calcobj1.SetupCalc(parmDict1)
2153    showEQ(calcobj1)
2154
2155    parmDict2 = {'0::Afrac:0':[0.0,True], '0::Afrac:1': [1.0,False]}
2156    print 'Dict = ',parmDict2
2157    calcobj.SetupCalc(parmDict2)
2158    showEQ(calcobj)
2159    calcobj1.SetupCalc(parmDict2)
2160    showEQ(calcobj1)
2161    calcobj2.SetupCalc(parmDict2)
2162    showEQ(calcobj2)
Note: See TracBrowser for help on using the repository browser.