source: trunk/GSASIIobj.py @ 1604

Last change on this file since 1604 was 1604, checked in by vondreele, 9 years ago

make modulation wave maps
fix atom index bugs
begin modulated structure imput to least squares

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