source: trunk/GSASIIobj.py @ 2184

Last change on this file since 2184 was 2062, checked in by vondreele, 10 years ago

exclude non variable SS parms from constraint lists
define ZigZag? & Block position modulations; eliminate Sawtooth (it's a variant of ZigZag?)
fix SS names in constraint lists
implement ZigZag? & Block position wave plots
implement ZigZag? & Block atom motion in structure plots
add movie making option (hidden - no file output for it yet)
fix LS I/O for ZigZag? & Block waves

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