1 | # -*- coding: utf-8 -*- |
---|
2 | #GSASIIobj - data objects for GSAS-II |
---|
3 | ########### SVN repository information ################### |
---|
4 | # $Date: 2023-09-29 20:47:55 +0000 (Fri, 29 Sep 2023) $ |
---|
5 | # $Author: vondreele $ |
---|
6 | # $Revision: 5663 $ |
---|
7 | # $URL: trunk/GSASIIobj.py $ |
---|
8 | # $Id: GSASIIobj.py 5663 2023-09-29 20:47:55Z vondreele $ |
---|
9 | ########### SVN repository information ################### |
---|
10 | |
---|
11 | ''' |
---|
12 | Classes and routines defined in :mod:`GSASIIobj` follow. |
---|
13 | ''' |
---|
14 | # Note that documentation for GSASIIobj.py has been moved |
---|
15 | # to file docs/source/GSASIIobj.rst |
---|
16 | |
---|
17 | from __future__ import division, print_function |
---|
18 | import platform |
---|
19 | import re |
---|
20 | import random as ran |
---|
21 | import sys |
---|
22 | import os.path |
---|
23 | if '2' in platform.python_version_tuple()[0]: |
---|
24 | import cPickle |
---|
25 | else: |
---|
26 | import pickle as cPickle |
---|
27 | import GSASIIpath |
---|
28 | import GSASIImath as G2mth |
---|
29 | import GSASIIspc as G2spc |
---|
30 | import numpy as np |
---|
31 | |
---|
32 | GSASIIpath.SetVersionNumber("$Revision: 5663 $") |
---|
33 | |
---|
34 | DefaultControls = { |
---|
35 | 'deriv type':'analytic Hessian', |
---|
36 | 'min dM/M':0.001,'shift factor':1.,'max cyc':3,'F**2':False,'SVDtol':1.e-6, |
---|
37 | 'UsrReject':{'minF/sig':0.,'MinExt':0.01,'MaxDF/F':100.,'MaxD':500.,'MinD':0.05}, |
---|
38 | 'Copy2Next':False,'Reverse Seq':False,'HatomFix':False, |
---|
39 | 'Author':'no name','newLeBail':False, |
---|
40 | 'FreePrm1':'Sample humidity (%)', |
---|
41 | 'FreePrm2':'Sample voltage (V)', |
---|
42 | 'FreePrm3':'Applied load (MN)', |
---|
43 | 'ShowCell':False, |
---|
44 | } |
---|
45 | '''Values to be used as defaults for the initial contents of the ``Controls`` |
---|
46 | data tree item. |
---|
47 | ''' |
---|
48 | def StripUnicode(string,subs='.'): |
---|
49 | '''Strip non-ASCII characters from strings |
---|
50 | |
---|
51 | :param str string: string to strip Unicode characters from |
---|
52 | :param str subs: character(s) to place into string in place of each |
---|
53 | Unicode character. Defaults to '.' |
---|
54 | |
---|
55 | :returns: a new string with only ASCII characters |
---|
56 | ''' |
---|
57 | s = '' |
---|
58 | for c in string: |
---|
59 | if ord(c) < 128: |
---|
60 | s += c |
---|
61 | else: |
---|
62 | s += subs |
---|
63 | return s |
---|
64 | # return s.encode('ascii','replace') |
---|
65 | |
---|
66 | def MakeUniqueLabel(lbl,labellist): |
---|
67 | '''Make sure that every a label is unique against a list by adding |
---|
68 | digits at the end until it is not found in list. |
---|
69 | |
---|
70 | :param str lbl: the input label |
---|
71 | :param list labellist: the labels that have already been encountered |
---|
72 | :returns: lbl if not found in labellist or lbl with ``_1-9`` (or |
---|
73 | ``_10-99``, etc.) appended at the end |
---|
74 | ''' |
---|
75 | lbl = StripUnicode(lbl.strip(),'_') |
---|
76 | if not lbl: # deal with a blank label |
---|
77 | lbl = '_1' |
---|
78 | if lbl not in labellist: |
---|
79 | labellist.append(lbl) |
---|
80 | return lbl |
---|
81 | i = 1 |
---|
82 | prefix = lbl |
---|
83 | if '_' in lbl: |
---|
84 | prefix = lbl[:lbl.rfind('_')] |
---|
85 | suffix = lbl[lbl.rfind('_')+1:] |
---|
86 | try: |
---|
87 | i = int(suffix)+1 |
---|
88 | except: # suffix could not be parsed |
---|
89 | i = 1 |
---|
90 | prefix = lbl |
---|
91 | while prefix+'_'+str(i) in labellist: |
---|
92 | i += 1 |
---|
93 | else: |
---|
94 | lbl = prefix+'_'+str(i) |
---|
95 | labellist.append(lbl) |
---|
96 | return lbl |
---|
97 | |
---|
98 | PhaseIdLookup = {} |
---|
99 | '''dict listing phase name and random Id keyed by sequential phase index as a str; |
---|
100 | best to access this using :func:`LookupPhaseName` |
---|
101 | ''' |
---|
102 | PhaseRanIdLookup = {} |
---|
103 | '''dict listing phase sequential index keyed by phase random Id; |
---|
104 | best to access this using :func:`LookupPhaseId` |
---|
105 | ''' |
---|
106 | HistIdLookup = {} |
---|
107 | '''dict listing histogram name and random Id, keyed by sequential histogram index as a str; |
---|
108 | best to access this using :func:`LookupHistName` |
---|
109 | ''' |
---|
110 | HistRanIdLookup = {} |
---|
111 | '''dict listing histogram sequential index keyed by histogram random Id; |
---|
112 | best to access this using :func:`LookupHistId` |
---|
113 | ''' |
---|
114 | AtomIdLookup = {} |
---|
115 | '''dict listing for each phase index as a str, the atom label and atom random Id, |
---|
116 | keyed by atom sequential index as a str; |
---|
117 | best to access this using :func:`LookupAtomLabel` |
---|
118 | ''' |
---|
119 | AtomRanIdLookup = {} |
---|
120 | '''dict listing for each phase the atom sequential index keyed by atom random Id; |
---|
121 | best to access this using :func:`LookupAtomId` |
---|
122 | ''' |
---|
123 | ShortPhaseNames = {} |
---|
124 | '''a dict containing a possibly shortened and when non-unique numbered |
---|
125 | version of the phase name. Keyed by the phase sequential index. |
---|
126 | ''' |
---|
127 | ShortHistNames = {} |
---|
128 | '''a dict containing a possibly shortened and when non-unique numbered |
---|
129 | version of the histogram name. Keyed by the histogram sequential index. |
---|
130 | ''' |
---|
131 | |
---|
132 | #VarDesc = {} # removed 1/30/19 BHT as no longer needed (I think) |
---|
133 | #''' This dictionary lists descriptions for GSAS-II variables, |
---|
134 | #as set in :func:`CompileVarDesc`. See that function for a description |
---|
135 | #for how keys and values are written. |
---|
136 | #''' |
---|
137 | |
---|
138 | reVarDesc = {} |
---|
139 | ''' This dictionary lists descriptions for GSAS-II variables where |
---|
140 | keys are compiled regular expressions that will match the name portion |
---|
141 | of a parameter name. Initialized in :func:`CompileVarDesc`. |
---|
142 | ''' |
---|
143 | |
---|
144 | reVarStep = {} |
---|
145 | ''' This dictionary lists the preferred step size for numerical |
---|
146 | derivative computation w/r to a GSAS-II variable. Keys are compiled |
---|
147 | regular expressions and values are the step size for that parameter. |
---|
148 | Initialized in :func:`CompileVarDesc`. |
---|
149 | ''' |
---|
150 | # create a default space group object for P1; N.B. fails when building documentation |
---|
151 | try: |
---|
152 | P1SGData = G2spc.SpcGroup('P 1')[1] # data structure for default space group |
---|
153 | except: |
---|
154 | pass |
---|
155 | |
---|
156 | def GetPhaseNames(fl): |
---|
157 | ''' Returns a list of phase names found under 'Phases' in GSASII gpx file |
---|
158 | NB: there is another one of these in GSASIIstrIO.py that uses the gpx filename |
---|
159 | |
---|
160 | :param file fl: opened .gpx file |
---|
161 | :return: list of phase names |
---|
162 | ''' |
---|
163 | PhaseNames = [] |
---|
164 | while True: |
---|
165 | try: |
---|
166 | data = cPickle.load(fl) |
---|
167 | except EOFError: |
---|
168 | break |
---|
169 | datum = data[0] |
---|
170 | if 'Phases' == datum[0]: |
---|
171 | for datus in data[1:]: |
---|
172 | PhaseNames.append(datus[0]) |
---|
173 | fl.seek(0) #reposition file |
---|
174 | return PhaseNames |
---|
175 | |
---|
176 | def SetNewPhase(Name='New Phase',SGData=None,cell=None,Super=None): |
---|
177 | '''Create a new phase dict with default values for various parameters |
---|
178 | |
---|
179 | :param str Name: Name for new Phase |
---|
180 | |
---|
181 | :param dict SGData: space group data from :func:`GSASIIspc:SpcGroup`; |
---|
182 | defaults to data for P 1 |
---|
183 | |
---|
184 | :param list cell: unit cell parameter list; defaults to |
---|
185 | [1.0,1.0,1.0,90.,90,90.,1.] |
---|
186 | |
---|
187 | ''' |
---|
188 | if SGData is None: SGData = P1SGData |
---|
189 | if cell is None: cell=[1.0,1.0,1.0,90.,90.,90.,1.] |
---|
190 | phaseData = { |
---|
191 | 'ranId':ran.randint(0,sys.maxsize), |
---|
192 | 'General':{ |
---|
193 | 'Name':Name, |
---|
194 | 'Type':'nuclear', |
---|
195 | 'Modulated':False, |
---|
196 | 'AtomPtrs':[3,1,7,9], |
---|
197 | 'SGData':SGData, |
---|
198 | 'Cell':[False,]+cell, |
---|
199 | 'Pawley dmin':1.0, |
---|
200 | 'Data plot type':'None', |
---|
201 | 'SH Texture':{ |
---|
202 | 'Order':0, |
---|
203 | 'Model':'cylindrical', |
---|
204 | 'Sample omega':[False,0.0], |
---|
205 | 'Sample chi':[False,0.0], |
---|
206 | 'Sample phi':[False,0.0], |
---|
207 | 'SH Coeff':[False,{}], |
---|
208 | 'SHShow':False, |
---|
209 | 'PFhkl':[0,0,1], |
---|
210 | 'PFxyz':[0,0,1], |
---|
211 | 'PlotType':'Pole figure', |
---|
212 | 'Penalty':[['',],0.1,False,1.0]}}, |
---|
213 | 'Atoms':[], |
---|
214 | 'Drawing':{}, |
---|
215 | 'Histograms':{}, |
---|
216 | 'Pawley ref':[], |
---|
217 | 'RBModels':{}, |
---|
218 | } |
---|
219 | if Super and Super.get('Use',False): |
---|
220 | phaseData['General'].update({'Modulated':True,'Super':True,'SuperSg':Super['ssSymb']}) |
---|
221 | phaseData['General']['SSGData'] = G2spc.SSpcGroup(SGData,Super['ssSymb'])[1] |
---|
222 | phaseData['General']['SuperVec'] = [Super['ModVec'],False,Super['maxH']] |
---|
223 | |
---|
224 | return phaseData |
---|
225 | |
---|
226 | def ReadCIF(URLorFile): |
---|
227 | '''Open a CIF, which may be specified as a file name or as a URL using PyCifRW |
---|
228 | (from James Hester). |
---|
229 | The open routine gets confused with DOS names that begin with a letter and colon |
---|
230 | "C:\\dir\" so this routine will try to open the passed name as a file and if that |
---|
231 | fails, try it as a URL |
---|
232 | |
---|
233 | :param str URLorFile: string containing a URL or a file name. Code will try first |
---|
234 | to open it as a file and then as a URL. |
---|
235 | |
---|
236 | :returns: a PyCifRW CIF object. |
---|
237 | ''' |
---|
238 | import CifFile as cif # PyCifRW from James Hester |
---|
239 | |
---|
240 | # alternate approach: |
---|
241 | #import urllib |
---|
242 | #ciffile = 'file:'+urllib.pathname2url(filename) |
---|
243 | |
---|
244 | try: |
---|
245 | fp = open(URLorFile,'r') |
---|
246 | cf = cif.ReadCif(fp) |
---|
247 | fp.close() |
---|
248 | return cf |
---|
249 | except IOError: |
---|
250 | return cif.ReadCif(URLorFile) |
---|
251 | |
---|
252 | def TestIndexAll(): |
---|
253 | '''Test if :func:`IndexAllIds` has been called to index all phases and |
---|
254 | histograms (this is needed before :func:`G2VarObj` can be used. |
---|
255 | |
---|
256 | :returns: Returns True if indexing is needed. |
---|
257 | ''' |
---|
258 | if PhaseIdLookup or AtomIdLookup or HistIdLookup: |
---|
259 | return False |
---|
260 | return True |
---|
261 | |
---|
262 | def IndexAllIds(Histograms,Phases): |
---|
263 | '''Scan through the used phases & histograms and create an index |
---|
264 | to the random numbers of phases, histograms and atoms. While doing this, |
---|
265 | confirm that assigned random numbers are unique -- just in case lightning |
---|
266 | strikes twice in the same place. |
---|
267 | |
---|
268 | Note: this code assumes that the atom random Id (ranId) is the last |
---|
269 | element each atom record. |
---|
270 | |
---|
271 | This is called when phases & histograms are looked up |
---|
272 | in these places (only): |
---|
273 | |
---|
274 | * :func:`GSASIIstrIO.GetUsedHistogramsAndPhases` (which loads the histograms and phases from a GPX file), |
---|
275 | * :meth:`~GSASIIdataGUI.GSASII.GetUsedHistogramsAndPhasesfromTree` (which does the same thing but from the data tree.) |
---|
276 | * :meth:`~GSASIIdataGUI.GSASII.OnFileClose` (clears out an old project) |
---|
277 | |
---|
278 | Note that globals :data:`PhaseIdLookup` and :data:`PhaseRanIdLookup` are |
---|
279 | also set in :func:`AddPhase2Index` to temporarily assign a phase number |
---|
280 | as a phase is being imported. |
---|
281 | |
---|
282 | TODO: do we need a lookup for rigid body variables? |
---|
283 | ''' |
---|
284 | # process phases and atoms |
---|
285 | PhaseIdLookup.clear() |
---|
286 | PhaseRanIdLookup.clear() |
---|
287 | AtomIdLookup.clear() |
---|
288 | AtomRanIdLookup.clear() |
---|
289 | ShortPhaseNames.clear() |
---|
290 | for ph in Phases: |
---|
291 | cx,ct,cs,cia = Phases[ph]['General']['AtomPtrs'] |
---|
292 | ranId = Phases[ph]['ranId'] |
---|
293 | while ranId in PhaseRanIdLookup: |
---|
294 | # Found duplicate random Id! note and reassign |
---|
295 | print ("\n\n*** Phase "+str(ph)+" has repeated ranId. Fixing.\n") |
---|
296 | Phases[ph]['ranId'] = ranId = ran.randint(0,sys.maxsize) |
---|
297 | pId = str(Phases[ph]['pId']) |
---|
298 | PhaseIdLookup[pId] = (ph,ranId) |
---|
299 | PhaseRanIdLookup[ranId] = pId |
---|
300 | shortname = ph #[:10] |
---|
301 | while shortname in ShortPhaseNames.values(): |
---|
302 | shortname = ph[:8] + ' ('+ pId + ')' |
---|
303 | ShortPhaseNames[pId] = shortname |
---|
304 | AtomIdLookup[pId] = {} |
---|
305 | AtomRanIdLookup[pId] = {} |
---|
306 | for iatm,at in enumerate(Phases[ph]['Atoms']): |
---|
307 | ranId = at[cia+8] |
---|
308 | while ranId in AtomRanIdLookup[pId]: # check for dups |
---|
309 | print ("\n\n*** Phase "+str(ph)+" atom "+str(iatm)+" has repeated ranId. Fixing.\n") |
---|
310 | at[cia+8] = ranId = ran.randint(0,sys.maxsize) |
---|
311 | AtomRanIdLookup[pId][ranId] = str(iatm) |
---|
312 | if Phases[ph]['General']['Type'] == 'macromolecular': |
---|
313 | label = '%s_%s_%s_%s'%(at[ct-1],at[ct-3],at[ct-4],at[ct-2]) |
---|
314 | else: |
---|
315 | label = at[ct-1] |
---|
316 | AtomIdLookup[pId][str(iatm)] = (label,ranId) |
---|
317 | # process histograms |
---|
318 | HistIdLookup.clear() |
---|
319 | HistRanIdLookup.clear() |
---|
320 | ShortHistNames.clear() |
---|
321 | for hist in Histograms: |
---|
322 | ranId = Histograms[hist]['ranId'] |
---|
323 | while ranId in HistRanIdLookup: |
---|
324 | # Found duplicate random Id! note and reassign |
---|
325 | print ("\n\n*** Histogram "+str(hist)+" has repeated ranId. Fixing.\n") |
---|
326 | Histograms[hist]['ranId'] = ranId = ran.randint(0,sys.maxsize) |
---|
327 | hId = str(Histograms[hist]['hId']) |
---|
328 | HistIdLookup[hId] = (hist,ranId) |
---|
329 | HistRanIdLookup[ranId] = hId |
---|
330 | shortname = hist[:15] |
---|
331 | while shortname in ShortHistNames.values(): |
---|
332 | shortname = hist[:11] + ' ('+ hId + ')' |
---|
333 | ShortHistNames[hId] = shortname |
---|
334 | |
---|
335 | def AddPhase2Index(rdObj,filename): |
---|
336 | '''Add a phase to the index during reading |
---|
337 | Used where constraints are generated during import (ISODISTORT CIFs) |
---|
338 | ''' |
---|
339 | ranId = rdObj.Phase['ranId'] |
---|
340 | ph = 'from '+filename # phase is not named yet |
---|
341 | if ranId in PhaseRanIdLookup: return |
---|
342 | maxph = -1 |
---|
343 | for r in PhaseRanIdLookup: |
---|
344 | maxph = max(maxph,int(PhaseRanIdLookup[r])) |
---|
345 | PhaseRanIdLookup[ranId] = pId = str(maxph + 1) |
---|
346 | PhaseIdLookup[pId] = (ph,ranId) |
---|
347 | shortname = 'from '+ os.path.splitext((os.path.split(filename))[1])[0] |
---|
348 | while shortname in ShortPhaseNames.values(): |
---|
349 | shortname = ph[:8] + ' ('+ pId + ')' |
---|
350 | ShortPhaseNames[pId] = shortname |
---|
351 | AtomIdLookup[pId] = {} |
---|
352 | AtomRanIdLookup[pId] = {} |
---|
353 | for iatm,at in enumerate(rdObj.Phase['Atoms']): |
---|
354 | ranId = at[-1] |
---|
355 | while ranId in AtomRanIdLookup[pId]: # check for dups |
---|
356 | print ("\n\n*** Phase "+str(ph)+" atom "+str(iatm)+" has repeated ranId. Fixing.\n") |
---|
357 | at[-1] = ranId = ran.randint(0,sys.maxsize) |
---|
358 | AtomRanIdLookup[pId][ranId] = str(iatm) |
---|
359 | #if Phases[ph]['General']['Type'] == 'macromolecular': |
---|
360 | # label = '%s_%s_%s_%s'%(at[ct-1],at[ct-3],at[ct-4],at[ct-2]) |
---|
361 | #else: |
---|
362 | # label = at[ct-1] |
---|
363 | label = at[0] |
---|
364 | AtomIdLookup[pId][str(iatm)] = (label,ranId) |
---|
365 | |
---|
366 | def LookupAtomId(pId,ranId): |
---|
367 | '''Get the atom number from a phase and atom random Id |
---|
368 | |
---|
369 | :param int/str pId: the sequential number of the phase |
---|
370 | :param int ranId: the random Id assigned to an atom |
---|
371 | |
---|
372 | :returns: the index number of the atom (str) |
---|
373 | ''' |
---|
374 | if not AtomRanIdLookup: |
---|
375 | raise Exception('Error: LookupAtomId called before IndexAllIds was run') |
---|
376 | if pId is None or pId == '': |
---|
377 | raise KeyError('Error: phase is invalid (None or blank)') |
---|
378 | pId = str(pId) |
---|
379 | if pId not in AtomRanIdLookup: |
---|
380 | raise KeyError('Error: LookupAtomId does not have phase '+pId) |
---|
381 | if ranId not in AtomRanIdLookup[pId]: |
---|
382 | raise KeyError('Error: LookupAtomId, ranId '+str(ranId)+' not in AtomRanIdLookup['+pId+']') |
---|
383 | return AtomRanIdLookup[pId][ranId] |
---|
384 | |
---|
385 | def LookupAtomLabel(pId,index): |
---|
386 | '''Get the atom label from a phase and atom index number |
---|
387 | |
---|
388 | :param int/str pId: the sequential number of the phase |
---|
389 | :param int index: the index of the atom in the list of atoms |
---|
390 | |
---|
391 | :returns: the label for the atom (str) and the random Id of the atom (int) |
---|
392 | ''' |
---|
393 | if not AtomIdLookup: |
---|
394 | raise Exception('Error: LookupAtomLabel called before IndexAllIds was run') |
---|
395 | if pId is None or pId == '': |
---|
396 | raise KeyError('Error: phase is invalid (None or blank)') |
---|
397 | pId = str(pId) |
---|
398 | if pId not in AtomIdLookup: |
---|
399 | raise KeyError('Error: LookupAtomLabel does not have phase '+pId) |
---|
400 | if index not in AtomIdLookup[pId]: |
---|
401 | raise KeyError('Error: LookupAtomLabel, ranId '+str(index)+' not in AtomRanIdLookup['+pId+']') |
---|
402 | return AtomIdLookup[pId][index] |
---|
403 | |
---|
404 | def LookupPhaseId(ranId): |
---|
405 | '''Get the phase number and name from a phase random Id |
---|
406 | |
---|
407 | :param int ranId: the random Id assigned to a phase |
---|
408 | :returns: the sequential Id (pId) number for the phase (str) |
---|
409 | ''' |
---|
410 | if not PhaseRanIdLookup: |
---|
411 | raise Exception('Error: LookupPhaseId called before IndexAllIds was run') |
---|
412 | if ranId not in PhaseRanIdLookup: |
---|
413 | raise KeyError('Error: LookupPhaseId does not have ranId '+str(ranId)) |
---|
414 | return PhaseRanIdLookup[ranId] |
---|
415 | |
---|
416 | def LookupPhaseName(pId): |
---|
417 | '''Get the phase number and name from a phase Id |
---|
418 | |
---|
419 | :param int/str pId: the sequential assigned to a phase |
---|
420 | :returns: (phase,ranId) where phase is the name of the phase (str) |
---|
421 | and ranId is the random # id for the phase (int) |
---|
422 | ''' |
---|
423 | if not PhaseIdLookup: |
---|
424 | raise Exception('Error: LookupPhaseName called before IndexAllIds was run') |
---|
425 | if pId is None or pId == '': |
---|
426 | raise KeyError('Error: phase is invalid (None or blank)') |
---|
427 | pId = str(pId) |
---|
428 | if pId not in PhaseIdLookup: |
---|
429 | raise KeyError('Error: LookupPhaseName does not have index '+pId) |
---|
430 | return PhaseIdLookup[pId] |
---|
431 | |
---|
432 | def LookupHistId(ranId): |
---|
433 | '''Get the histogram number and name from a histogram random Id |
---|
434 | |
---|
435 | :param int ranId: the random Id assigned to a histogram |
---|
436 | :returns: the sequential Id (hId) number for the histogram (str) |
---|
437 | ''' |
---|
438 | if not HistRanIdLookup: |
---|
439 | raise Exception('Error: LookupHistId called before IndexAllIds was run') |
---|
440 | if ranId not in HistRanIdLookup: |
---|
441 | raise KeyError('Error: LookupHistId does not have ranId '+str(ranId)) |
---|
442 | return HistRanIdLookup[ranId] |
---|
443 | |
---|
444 | def LookupHistName(hId): |
---|
445 | '''Get the histogram number and name from a histogram Id |
---|
446 | |
---|
447 | :param int/str hId: the sequential assigned to a histogram |
---|
448 | :returns: (hist,ranId) where hist is the name of the histogram (str) |
---|
449 | and ranId is the random # id for the histogram (int) |
---|
450 | ''' |
---|
451 | if not HistIdLookup: |
---|
452 | raise Exception('Error: LookupHistName called before IndexAllIds was run') |
---|
453 | if hId is None or hId == '': |
---|
454 | raise KeyError('Error: histogram is invalid (None or blank)') |
---|
455 | hId = str(hId) |
---|
456 | if hId not in HistIdLookup: |
---|
457 | raise KeyError('Error: LookupHistName does not have index '+hId) |
---|
458 | return HistIdLookup[hId] |
---|
459 | |
---|
460 | def fmtVarDescr(varname): |
---|
461 | '''Return a string with a more complete description for a GSAS-II variable |
---|
462 | |
---|
463 | :param str varname: A full G2 variable name with 2 or 3 or 4 |
---|
464 | colons (<p>:<h>:name[:<a>] or <p>::RBname:<r>:<t>]) |
---|
465 | |
---|
466 | :returns: a string with the description |
---|
467 | ''' |
---|
468 | s,l = VarDescr(varname) |
---|
469 | return s+": "+l |
---|
470 | |
---|
471 | def VarDescr(varname): |
---|
472 | '''Return two strings with a more complete description for a GSAS-II variable |
---|
473 | |
---|
474 | :param str name: A full G2 variable name with 2 or 3 or 4 |
---|
475 | colons (<p>:<h>:name[:<a>] or <p>::RBname:<r>:<t>]) |
---|
476 | |
---|
477 | :returns: (loc,meaning) where loc describes what item the variable is mapped |
---|
478 | (phase, histogram, etc.) and meaning describes what the variable does. |
---|
479 | ''' |
---|
480 | |
---|
481 | # special handling for parameter names without a colon |
---|
482 | # for now, assume self-defining |
---|
483 | if varname.find(':') == -1: |
---|
484 | return "Global",varname |
---|
485 | |
---|
486 | l = getVarDescr(varname) |
---|
487 | if not l: |
---|
488 | return ("invalid variable name ("+str(varname)+")!"),"" |
---|
489 | # return "invalid variable name!","" |
---|
490 | |
---|
491 | if not l[-1]: |
---|
492 | l[-1] = "(variable needs a definition! Set it in CompileVarDesc)" |
---|
493 | |
---|
494 | if len(l) == 3: #SASD variable name! |
---|
495 | s = 'component:'+l[1] |
---|
496 | return s,l[-1] |
---|
497 | s = "" |
---|
498 | if l[0] is not None and l[1] is not None: # HAP: keep short |
---|
499 | if l[2] == "Scale": # fix up ambigous name |
---|
500 | l[5] = "Phase fraction" |
---|
501 | if l[0] == '*': |
---|
502 | lbl = 'Seq. ref.' |
---|
503 | else: |
---|
504 | lbl = ShortPhaseNames.get(l[0],'? #'+str(l[0])) |
---|
505 | if l[1] == '*': |
---|
506 | hlbl = 'Seq. ref.' |
---|
507 | else: |
---|
508 | hlbl = ShortHistNames.get(l[1],'? #'+str(l[1])) |
---|
509 | if hlbl[:4] == 'HKLF': |
---|
510 | hlbl = 'Xtl='+hlbl[5:] |
---|
511 | elif hlbl[:4] == 'PWDR': |
---|
512 | hlbl = 'Pwd='+hlbl[5:] |
---|
513 | else: |
---|
514 | hlbl = 'Hist='+hlbl |
---|
515 | s = "Ph="+str(lbl)+" * "+str(hlbl) |
---|
516 | else: |
---|
517 | if l[2] == "Scale": # fix up ambigous name: must be scale factor, since not HAP |
---|
518 | l[5] = "Scale factor" |
---|
519 | if l[2] == 'Back': # background parameters are "special", alas |
---|
520 | s = 'Hist='+ShortHistNames.get(l[1],'? #'+str(l[1])) |
---|
521 | l[-1] += ' #'+str(l[3]) |
---|
522 | elif l[4] is not None: # rigid body parameter or modulation parm |
---|
523 | lbl = ShortPhaseNames.get(l[0],'phase?') |
---|
524 | if 'RB' in l[2]: #rigid body parm |
---|
525 | s = "RB body #"+str(l[3])+" (type "+str(l[4])+") in "+str(lbl) + ',' |
---|
526 | else: #modulation parm |
---|
527 | s = 'Atom %s wave %s in %s'%(LookupAtomLabel(l[0],l[3])[0],l[4],lbl) |
---|
528 | elif l[3] is not None: # atom parameter, |
---|
529 | lbl = ShortPhaseNames.get(l[0],'phase?') |
---|
530 | try: |
---|
531 | albl = LookupAtomLabel(l[0],l[3])[0] |
---|
532 | except KeyError: |
---|
533 | albl = 'Atom?' |
---|
534 | s = "Atom "+str(albl)+" in "+str(lbl) |
---|
535 | elif l[0] == '*': |
---|
536 | s = "All phases " |
---|
537 | elif l[0] is not None: |
---|
538 | lbl = ShortPhaseNames.get(l[0],'phase?') |
---|
539 | s = "Phase "+str(lbl) |
---|
540 | elif l[1] == '*': |
---|
541 | s = 'All hists' |
---|
542 | elif l[1] is not None: |
---|
543 | hlbl = ShortHistNames.get(l[1],'? #'+str(l[1])) |
---|
544 | if hlbl[:4] == 'HKLF': |
---|
545 | hlbl = 'Xtl='+hlbl[5:] |
---|
546 | elif hlbl[:4] == 'PWDR': |
---|
547 | hlbl = 'Pwd='+hlbl[5:] |
---|
548 | else: |
---|
549 | hlbl = 'Hist='+hlbl |
---|
550 | s = str(hlbl) |
---|
551 | if not s: |
---|
552 | s = 'Global' |
---|
553 | return s,l[-1] |
---|
554 | |
---|
555 | def getVarDescr(varname): |
---|
556 | '''Return a short description for a GSAS-II variable |
---|
557 | |
---|
558 | :param str name: A full G2 variable name with 2 or 3 or 4 |
---|
559 | colons (<p>:<h>:name[:<a1>][:<a2>]) |
---|
560 | |
---|
561 | :returns: a six element list as [`p`,`h`,`name`,`a1`,`a2`,`description`], |
---|
562 | where `p`, `h`, `a1`, `a2` are str values or `None`, for the phase number, |
---|
563 | the histogram number and the atom number; `name` will always be |
---|
564 | a str; and `description` is str or `None`. |
---|
565 | If the variable name is incorrectly formed (for example, wrong |
---|
566 | number of colons), `None` is returned instead of a list. |
---|
567 | ''' |
---|
568 | l = varname.split(':') |
---|
569 | if len(l) == 2: #SASD parameter name |
---|
570 | return varname,l[0],getDescr(l[1]) |
---|
571 | if len(l) == 3: |
---|
572 | l += [None,None] |
---|
573 | elif len(l) == 4: |
---|
574 | l += [None] |
---|
575 | elif len(l) != 5: |
---|
576 | return None |
---|
577 | for i in (0,1,3,4): |
---|
578 | if l[i] == "": |
---|
579 | l[i] = None |
---|
580 | l += [getDescr(l[2])] |
---|
581 | return l |
---|
582 | |
---|
583 | def CompileVarDesc(): |
---|
584 | '''Set the values in the variable lookup tables |
---|
585 | (:attr:`reVarDesc` and :attr:`reVarStep`). |
---|
586 | This is called in :func:`getDescr` and :func:`getVarStep` so this |
---|
587 | initialization is always done before use. These variables are |
---|
588 | also used in script `makeVarTbl.py` which creates the table in section 3.2 |
---|
589 | of the Sphinx docs (:ref:`VarNames_table`). |
---|
590 | |
---|
591 | Note that keys may contain regular expressions, where '[xyz]' |
---|
592 | matches 'x' 'y' or 'z' (equivalently '[x-z]' describes this as range |
---|
593 | of values). '.*' matches any string. For example:: |
---|
594 | |
---|
595 | 'AUiso':'Atomic isotropic displacement parameter', |
---|
596 | |
---|
597 | will match variable ``'p::AUiso:a'``. |
---|
598 | If parentheses are used in the key, the contents of those parentheses can be |
---|
599 | used in the value, such as:: |
---|
600 | |
---|
601 | 'AU([123][123])':'Atomic anisotropic displacement parameter U\\1', |
---|
602 | |
---|
603 | will match ``AU11``, ``AU23``,... and `U11`, `U23` etc will be displayed |
---|
604 | in the value when used. |
---|
605 | |
---|
606 | ''' |
---|
607 | if reVarDesc: return # already done |
---|
608 | for key,value in { |
---|
609 | # derived or other sequential vars |
---|
610 | '([abc])$' : 'Lattice parameter, \\1, from Ai and Djk', # N.B. '$' prevents match if any characters follow |
---|
611 | u'\u03B1' : u'Lattice parameter, \u03B1, computed from both Ai and Djk', |
---|
612 | u'\u03B2' : u'Lattice parameter, \u03B2, computed from both Ai and Djk', |
---|
613 | u'\u03B3' : u'Lattice parameter, \u03B3, computed from both Ai and Djk', |
---|
614 | # ambiguous, alas: |
---|
615 | 'Scale' : 'Phase fraction (as p:h:Scale) or Histogram scale factor (as :h:Scale)', |
---|
616 | # Phase vars (p::<var>) |
---|
617 | 'A([0-5])' : ('Reciprocal metric tensor component \\1',1e-5), |
---|
618 | '([vV]ol)' : 'Unit cell volume', # probably an error that both upper and lower case are used |
---|
619 | # Atom vars (p::<var>:a) |
---|
620 | 'dA([xyz])$' : ('Refined change to atomic coordinate, \\1',1e-6), |
---|
621 | 'A([xyz])$' : 'Fractional atomic coordinate, \\1', |
---|
622 | 'AUiso':('Atomic isotropic displacement parameter',1e-4), |
---|
623 | 'AU([123][123])':('Atomic anisotropic displacement parameter U\\1',1e-4), |
---|
624 | 'Afrac': ('Atomic site fraction parameter',1e-5), |
---|
625 | 'Amul': 'Atomic site multiplicity value', |
---|
626 | 'AM([xyz])$' : 'Atomic magnetic moment parameter, \\1', |
---|
627 | # Atom deformation parameters |
---|
628 | 'Akappa([0-6])' : ' Atomic orbital softness for orbital, \\1', |
---|
629 | 'ANe([01])' : ' Atomic <j0> orbital population for orbital, \\1', |
---|
630 | 'AD\\([0-6],[0-6]\\)([0-6])' : ' Atomic sp. harm. coeff for orbital, \\1', |
---|
631 | 'AD\\([0-6],-[0-6]\\)([0-6])' : ' Atomic sp. harm. coeff for orbital, \\1', #need both! |
---|
632 | # Hist (:h:<var>) & Phase (HAP) vars (p:h:<var>) |
---|
633 | 'Back(.*)': 'Background term #\\1', |
---|
634 | 'BkPkint;(.*)':'Background peak #\\1 intensity', |
---|
635 | 'BkPkpos;(.*)':'Background peak #\\1 position', |
---|
636 | 'BkPksig;(.*)':'Background peak #\\1 Gaussian width', |
---|
637 | 'BkPkgam;(.*)':'Background peak #\\1 Cauchy width', |
---|
638 | # 'Back File' : 'Background file name', |
---|
639 | 'BF mult' : 'Background file multiplier', |
---|
640 | 'Bab([AU])': 'Babinet solvent scattering coef. \\1', |
---|
641 | 'D([123][123])' : 'Anisotropic strain coef. \\1', |
---|
642 | 'Extinction' : 'Extinction coef.', |
---|
643 | 'MD' : 'March-Dollase coef.', |
---|
644 | 'Mustrain;.*' : 'Microstrain coefficient (delta Q/Q x 10**6)', |
---|
645 | 'Size;.*' : 'Crystallite size value (in microns)', |
---|
646 | 'eA$' : 'Cubic mustrain value', |
---|
647 | 'Ep$' : 'Primary extinction', |
---|
648 | 'Es$' : 'Secondary type II extinction', |
---|
649 | 'Eg$' : 'Secondary type I extinction', |
---|
650 | 'Flack' : 'Flack parameter', |
---|
651 | 'TwinFr' : 'Twin fraction', |
---|
652 | 'Layer Disp' : 'Layer displacement along beam', |
---|
653 | #Histogram vars (:h:<var>) |
---|
654 | 'Absorption' : 'Absorption coef.', |
---|
655 | 'LayerDisp' : 'Bragg-Brentano Layer displacement', |
---|
656 | 'Displace([XY])' : ('Debye-Scherrer sample displacement \\1',0.1), |
---|
657 | 'Lam' : ('Wavelength',1e-6), |
---|
658 | 'I\\(L2\\)\\/I\\(L1\\)' : ('Ka2/Ka1 intensity ratio',0.001), |
---|
659 | 'Polariz.' : ('Polarization correction',1e-3), |
---|
660 | 'SH/L' : ('FCJ peak asymmetry correction',1e-4), |
---|
661 | '([UVW])$' : ('Gaussian instrument broadening \\1',1e-5), |
---|
662 | '([XYZ])$' : ('Cauchy instrument broadening \\1',1e-5), |
---|
663 | 'Zero' : 'Debye-Scherrer zero correction', |
---|
664 | 'Shift' : 'Bragg-Brentano sample displ.', |
---|
665 | 'SurfRoughA' : 'Bragg-Brenano surface roughness A', |
---|
666 | 'SurfRoughB' : 'Bragg-Brenano surface roughness B', |
---|
667 | 'Transparency' : 'Bragg-Brentano sample tranparency', |
---|
668 | 'DebyeA' : 'Debye model amplitude', |
---|
669 | 'DebyeR' : 'Debye model radius', |
---|
670 | 'DebyeU' : 'Debye model Uiso', |
---|
671 | 'RBV.*' : 'Vector rigid body parameter', |
---|
672 | 'RBVO([aijk])' : 'Vector rigid body orientation parameter \\1', |
---|
673 | 'RBVP([xyz])' : 'Vector rigid body \\1 position parameter', |
---|
674 | 'RBVf' : 'Vector rigid body site fraction', |
---|
675 | 'RBV([TLS])([123AB][123AB])' : 'Residue rigid body group disp. param.', |
---|
676 | 'RBVU' : 'Residue rigid body group Uiso param.', |
---|
677 | 'RBRO([aijk])' : 'Residue rigid body orientation parameter \\1', |
---|
678 | 'RBRP([xyz])' : 'Residue rigid body \\1 position parameter', |
---|
679 | 'RBRTr;.*' : 'Residue rigid body torsion parameter', |
---|
680 | 'RBRf' : 'Residue rigid body site fraction', |
---|
681 | 'RBR([TLS])([123AB][123AB])' : 'Residue rigid body group disp. param.', |
---|
682 | 'RBRU' : 'Residue rigid body group Uiso param.', |
---|
683 | 'RBSAtNo' : 'Atom number for spinning rigid body', |
---|
684 | 'RBSO([aijk])' : 'Spinning rigid body orientation parameter \\1', |
---|
685 | 'RBSP([xyz])' : 'Spinning rigid body \\1 position parameter', |
---|
686 | 'RBSShRadius' : 'Spinning rigid body shell radius', |
---|
687 | 'RBSShC([1-20,1-20])' : 'Spinning rigid body sph. harmonics term', |
---|
688 | 'constr([0-9]*)' : 'Generated degree of freedom from constraint', |
---|
689 | 'nv-(.+)' : 'New variable assignment with name \\1', |
---|
690 | # supersymmetry parameters p::<var>:a:o 'Flen','Fcent'? |
---|
691 | 'mV([0-2])$' : 'Modulation vector component \\1', |
---|
692 | 'Fsin' : 'Sin site fraction modulation', |
---|
693 | 'Fcos' : 'Cos site fraction modulation', |
---|
694 | 'Fzero' : 'Crenel function offset', #may go away |
---|
695 | 'Fwid' : 'Crenel function width', |
---|
696 | 'Tmin' : 'ZigZag/Block min location', |
---|
697 | 'Tmax' : 'ZigZag/Block max location', |
---|
698 | '([XYZ])max': 'ZigZag/Block max value for \\1', |
---|
699 | '([XYZ])sin' : 'Sin position wave for \\1', |
---|
700 | '([XYZ])cos' : 'Cos position wave for \\1', |
---|
701 | 'U([123][123])sin$' : 'Sin thermal wave for U\\1', |
---|
702 | 'U([123][123])cos$' : 'Cos thermal wave for U\\1', |
---|
703 | 'M([XYZ])sin$' : 'Sin mag. moment wave for \\1', |
---|
704 | 'M([XYZ])cos$' : 'Cos mag. moment wave for \\1', |
---|
705 | # PDF peak parms (l:<var>;l = peak no.) |
---|
706 | 'PDFpos' : 'PDF peak position', |
---|
707 | 'PDFmag' : 'PDF peak magnitude', |
---|
708 | 'PDFsig' : 'PDF peak std. dev.', |
---|
709 | # SASD vars (l:<var>;l = component) |
---|
710 | 'Aspect ratio' : 'Particle aspect ratio', |
---|
711 | 'Length' : 'Cylinder length', |
---|
712 | 'Diameter' : 'Cylinder/disk diameter', |
---|
713 | 'Thickness' : 'Disk thickness', |
---|
714 | 'Shell thickness' : 'Multiplier to get inner(<1) or outer(>1) sphere radius', |
---|
715 | 'Dist' : 'Interparticle distance', |
---|
716 | 'VolFr' : 'Dense scatterer volume fraction', |
---|
717 | 'epis' : 'Sticky sphere epsilon', |
---|
718 | 'Sticky' : 'Stickyness', |
---|
719 | 'Depth' : 'Well depth', |
---|
720 | 'Width' : 'Well width', |
---|
721 | 'Volume' : 'Particle volume', |
---|
722 | 'Radius' : 'Sphere/cylinder/disk radius', |
---|
723 | 'Mean' : 'Particle mean radius', |
---|
724 | 'StdDev' : 'Standard deviation in Mean', |
---|
725 | 'G$': 'Guinier prefactor', |
---|
726 | 'Rg$': 'Guinier radius of gyration', |
---|
727 | 'B$': 'Porod prefactor', |
---|
728 | 'P$': 'Porod power', |
---|
729 | 'Cutoff': 'Porod cutoff', |
---|
730 | 'PkInt': 'Bragg peak intensity', |
---|
731 | 'PkPos': 'Bragg peak position', |
---|
732 | 'PkSig': 'Bragg peak sigma', |
---|
733 | 'PkGam': 'Bragg peak gamma', |
---|
734 | 'e([12][12])' : 'strain tensor e\\1', # strain vars e11, e22, e12 |
---|
735 | 'Dcalc': 'Calc. d-spacing', |
---|
736 | 'Back$': 'background parameter', |
---|
737 | 'pos$': 'peak position', |
---|
738 | 'int$': 'peak intensity', |
---|
739 | 'WgtFrac':'phase weight fraction', |
---|
740 | 'alpha':'TOF profile term', |
---|
741 | 'alpha-([01])':'Pink profile term', |
---|
742 | 'beta-([01q])':'TOF/Pink profile term', |
---|
743 | 'sig-([012q])':'TOF profile term', |
---|
744 | 'dif([ABC])':'TOF to d-space calibration', |
---|
745 | 'C\\([0-9]*,[0-9]*\\)' : 'spherical harmonics preferred orientation coef.', |
---|
746 | 'Pressure': 'Pressure level for measurement in MPa', |
---|
747 | 'Temperature': 'T value for measurement, K', |
---|
748 | 'FreePrm([123])': 'User defined measurement parameter \\1', |
---|
749 | 'Gonio. radius': 'Distance from sample to detector, mm', |
---|
750 | }.items(): |
---|
751 | # Needs documentation: HAP: LeBail, newLeBail |
---|
752 | # hist: Azimuth, Chi, Omega, Phi, Bank, nDebye, nPeaks |
---|
753 | |
---|
754 | if len(value) == 2: |
---|
755 | #VarDesc[key] = value[0] |
---|
756 | reVarDesc[re.compile(key)] = value[0] |
---|
757 | reVarStep[re.compile(key)] = value[1] |
---|
758 | else: |
---|
759 | #VarDesc[key] = value |
---|
760 | reVarDesc[re.compile(key)] = value |
---|
761 | |
---|
762 | def removeNonRefined(parmList): |
---|
763 | '''Remove items from variable list that are not refined and should not |
---|
764 | appear as options for constraints |
---|
765 | |
---|
766 | :param list parmList: a list of strings of form "p:h:VAR:a" where |
---|
767 | VAR is the variable name |
---|
768 | |
---|
769 | :returns: a list after removing variables where VAR matches a |
---|
770 | entry in local variable NonRefinedList |
---|
771 | ''' |
---|
772 | NonRefinedList = ['Omega','Type','Chi','Phi', 'Azimuth','Gonio. radius', |
---|
773 | 'Lam1','Lam2','Back','Temperature','Pressure', |
---|
774 | 'FreePrm1','FreePrm2','FreePrm3', |
---|
775 | 'Source','nPeaks','LeBail','newLeBail','Bank', |
---|
776 | 'nDebye', #'', |
---|
777 | ] |
---|
778 | return [prm for prm in parmList if prm.split(':')[2] not in NonRefinedList] |
---|
779 | |
---|
780 | def getDescr(name): |
---|
781 | '''Return a short description for a GSAS-II variable |
---|
782 | |
---|
783 | :param str name: The descriptive part of the variable name without colons (:) |
---|
784 | |
---|
785 | :returns: a short description or None if not found |
---|
786 | ''' |
---|
787 | |
---|
788 | CompileVarDesc() # compile the regular expressions, if needed |
---|
789 | for key in reVarDesc: |
---|
790 | m = key.match(name) |
---|
791 | if m: |
---|
792 | reVarDesc[key] |
---|
793 | try: |
---|
794 | return m.expand(reVarDesc[key]) |
---|
795 | except: |
---|
796 | print('Error in key: %s'%key) |
---|
797 | return None |
---|
798 | |
---|
799 | def getVarStep(name,parmDict=None): |
---|
800 | '''Return a step size for computing the derivative of a GSAS-II variable |
---|
801 | |
---|
802 | :param str name: A complete variable name (with colons, :) |
---|
803 | :param dict parmDict: A dict with parameter values or None (default) |
---|
804 | |
---|
805 | :returns: a float that should be an appropriate step size, either from |
---|
806 | the value supplied in :func:`CompileVarDesc` or based on the value for |
---|
807 | name in parmDict, if supplied. If not found or the value is zero, |
---|
808 | a default value of 1e-5 is used. If parmDict is None (default) and |
---|
809 | no value is provided in :func:`CompileVarDesc`, then None is returned. |
---|
810 | ''' |
---|
811 | CompileVarDesc() # compile the regular expressions, if needed |
---|
812 | for key in reVarStep: |
---|
813 | m = key.match(name) |
---|
814 | if m: |
---|
815 | return reVarStep[key] |
---|
816 | if parmDict is None: return None |
---|
817 | val = parmDict.get(key,0.0) |
---|
818 | if abs(val) > 0.05: |
---|
819 | return abs(val)/1000. |
---|
820 | else: |
---|
821 | return 1e-5 |
---|
822 | |
---|
823 | def GenWildCard(varlist): |
---|
824 | '''Generate wildcard versions of G2 variables. These introduce '*' |
---|
825 | for a phase, histogram or atom number (but only for one of these |
---|
826 | fields) but only when there is more than one matching variable in the |
---|
827 | input variable list. So if the input is this:: |
---|
828 | |
---|
829 | varlist = ['0::AUiso:0', '0::AUiso:1', '1::AUiso:0'] |
---|
830 | |
---|
831 | then the output will be this:: |
---|
832 | |
---|
833 | wildList = ['*::AUiso:0', '0::AUiso:*'] |
---|
834 | |
---|
835 | :param list varlist: an input list of GSAS-II variable names |
---|
836 | (such as 0::AUiso:0) |
---|
837 | |
---|
838 | :returns: wildList, the generated list of wild card variable names. |
---|
839 | ''' |
---|
840 | wild = [] |
---|
841 | for i in (0,1,3): |
---|
842 | currentL = varlist[:] |
---|
843 | while currentL: |
---|
844 | item1 = currentL.pop(0) |
---|
845 | i1splt = item1.split(':') |
---|
846 | if i >= len(i1splt): continue |
---|
847 | if i1splt[i]: |
---|
848 | nextL = [] |
---|
849 | i1splt[i] = '[0-9]+' |
---|
850 | rexp = re.compile(':'.join(i1splt)) |
---|
851 | matchlist = [item1] |
---|
852 | for nxtitem in currentL: |
---|
853 | if rexp.match(nxtitem): |
---|
854 | matchlist += [nxtitem] |
---|
855 | else: |
---|
856 | nextL.append(nxtitem) |
---|
857 | if len(matchlist) > 1: |
---|
858 | i1splt[i] = '*' |
---|
859 | wild.append(':'.join(i1splt)) |
---|
860 | currentL = nextL |
---|
861 | return wild |
---|
862 | |
---|
863 | def LookupWildCard(varname,varlist): |
---|
864 | '''returns a list of variable names from list varname |
---|
865 | that match wildcard name in varname |
---|
866 | |
---|
867 | :param str varname: a G2 variable name containing a wildcard |
---|
868 | (such as \\*::var) |
---|
869 | :param list varlist: the list of all variable names used in |
---|
870 | the current project |
---|
871 | :returns: a list of matching GSAS-II variables (may be empty) |
---|
872 | ''' |
---|
873 | rexp = re.compile(varname.replace('*','[0-9]+')) |
---|
874 | return sorted([var for var in varlist if rexp.match(var)]) |
---|
875 | |
---|
876 | def prmLookup(name,prmDict): |
---|
877 | '''Looks for a parameter in a min/max dictionary, optionally |
---|
878 | considering a wild card for histogram or atom number (use of |
---|
879 | both will never occur at the same time). |
---|
880 | |
---|
881 | :param name: a GSAS-II parameter name (str, see :func:`getVarDescr` |
---|
882 | and :func:`CompileVarDesc`) or a :class:`G2VarObj` object. |
---|
883 | :param dict prmDict: a min/max dictionary, (parmMinDict |
---|
884 | or parmMaxDict in Controls) where keys are :class:`G2VarObj` |
---|
885 | objects. |
---|
886 | :returns: Two values, (**matchname**, **value**), are returned where: |
---|
887 | |
---|
888 | * **matchname** *(str)* is the :class:`G2VarObj` object |
---|
889 | corresponding to the actual matched name, |
---|
890 | which could contain a wildcard even if **name** does not; and |
---|
891 | * **value** *(float)* which contains the parameter limit. |
---|
892 | ''' |
---|
893 | for key,value in prmDict.items(): |
---|
894 | if str(key) == str(name): return key,value |
---|
895 | if key == name: return key,value |
---|
896 | return None,None |
---|
897 | |
---|
898 | |
---|
899 | def _lookup(dic,key): |
---|
900 | '''Lookup a key in a dictionary, where None returns an empty string |
---|
901 | but an unmatched key returns a question mark. Used in :class:`G2VarObj` |
---|
902 | ''' |
---|
903 | if key is None: |
---|
904 | return "" |
---|
905 | elif key == "*": |
---|
906 | return "*" |
---|
907 | else: |
---|
908 | return dic.get(key,'?') |
---|
909 | |
---|
910 | def SortVariables(varlist): |
---|
911 | '''Sorts variable names in a sensible manner |
---|
912 | ''' |
---|
913 | def cvnnums(var): |
---|
914 | v = [] |
---|
915 | for i in var.split(':'): |
---|
916 | # if i == '' or i == '*': |
---|
917 | # v.append(-1) |
---|
918 | # continue |
---|
919 | try: |
---|
920 | v.append(int(i)) |
---|
921 | except: |
---|
922 | v.append(-1) |
---|
923 | return v |
---|
924 | return sorted(varlist,key=cvnnums) |
---|
925 | |
---|
926 | class G2VarObj(object): |
---|
927 | '''Defines a GSAS-II variable either using the phase/atom/histogram |
---|
928 | unique Id numbers or using a character string that specifies |
---|
929 | variables by phase/atom/histogram number (which can change). |
---|
930 | Note that :func:`GSASIIstrIO.GetUsedHistogramsAndPhases`, |
---|
931 | which calls :func:`IndexAllIds` (or |
---|
932 | :func:`GSASIIscriptable.G2Project.index_ids`) should be used to |
---|
933 | (re)load the current Ids |
---|
934 | before creating or later using the G2VarObj object. |
---|
935 | |
---|
936 | This can store rigid body variables, but does not translate the residue # and |
---|
937 | body # to/from random Ids |
---|
938 | |
---|
939 | A :class:`G2VarObj` object can be created with a single parameter: |
---|
940 | |
---|
941 | :param str/tuple varname: a single value can be used to create a :class:`G2VarObj` |
---|
942 | object. If a string, it must be of form "p:h:var" or "p:h:var:a", where |
---|
943 | |
---|
944 | * p is the phase number (which may be left blank or may be '*' to indicate all phases); |
---|
945 | * h is the histogram number (which may be left blank or may be '*' to indicate all histograms); |
---|
946 | * a is the atom number (which may be left blank in which case the third colon is omitted). |
---|
947 | The atom number can be specified as '*' if a phase number is specified (not as '*'). |
---|
948 | For rigid body variables, specify a will be a string of form "residue:body#" |
---|
949 | |
---|
950 | Alternately a single tuple of form (Phase,Histogram,VarName,AtomID) can be used, where |
---|
951 | Phase, Histogram, and AtomID are None or are ranId values (or one can be '*') |
---|
952 | and VarName is a string. Note that if Phase is '*' then the AtomID is an atom number. |
---|
953 | For a rigid body variables, AtomID is a string of form "residue:body#". |
---|
954 | |
---|
955 | If four positional arguments are supplied, they are: |
---|
956 | |
---|
957 | :param str/int phasenum: The number for the phase (or None or '*') |
---|
958 | :param str/int histnum: The number for the histogram (or None or '*') |
---|
959 | :param str varname: a single value can be used to create a :class:`G2VarObj` |
---|
960 | :param str/int atomnum: The number for the atom (or None or '*') |
---|
961 | |
---|
962 | ''' |
---|
963 | IDdict = {} |
---|
964 | IDdict['phases'] = {} |
---|
965 | IDdict['hists'] = {} |
---|
966 | IDdict['atoms'] = {} |
---|
967 | def __init__(self,*args): |
---|
968 | self.phase = None |
---|
969 | self.histogram = None |
---|
970 | self.name = '' |
---|
971 | self.atom = None |
---|
972 | if len(args) == 1 and (type(args[0]) is list or type(args[0]) is tuple) and len(args[0]) == 4: |
---|
973 | # single arg with 4 values |
---|
974 | self.phase,self.histogram,self.name,self.atom = args[0] |
---|
975 | elif len(args) == 1 and ':' in args[0]: |
---|
976 | #parse a string |
---|
977 | lst = args[0].split(':') |
---|
978 | if lst[0] == '*': |
---|
979 | self.phase = '*' |
---|
980 | if len(lst) > 3: |
---|
981 | self.atom = lst[3] |
---|
982 | self.histogram = HistIdLookup.get(lst[1],[None,None])[1] |
---|
983 | elif lst[1] == '*': |
---|
984 | self.histogram = '*' |
---|
985 | self.phase = PhaseIdLookup.get(lst[0],[None,None])[1] |
---|
986 | else: |
---|
987 | self.histogram = HistIdLookup.get(lst[1],[None,None])[1] |
---|
988 | self.phase = PhaseIdLookup.get(lst[0],[None,None])[1] |
---|
989 | if len(lst) == 4: |
---|
990 | if lst[3] == '*': |
---|
991 | self.atom = '*' |
---|
992 | else: |
---|
993 | self.atom = AtomIdLookup[lst[0]].get(lst[3],[None,None])[1] |
---|
994 | elif len(lst) == 5: |
---|
995 | self.atom = lst[3]+":"+lst[4] |
---|
996 | elif len(lst) == 3: |
---|
997 | pass |
---|
998 | else: |
---|
999 | raise Exception("Incorrect number of colons in var name "+str(args[0])) |
---|
1000 | self.name = lst[2] |
---|
1001 | elif len(args) == 4: |
---|
1002 | if args[0] == '*': |
---|
1003 | self.phase = '*' |
---|
1004 | self.atom = args[3] |
---|
1005 | else: |
---|
1006 | self.phase = PhaseIdLookup.get(str(args[0]),[None,None])[1] |
---|
1007 | if args[3] == '*': |
---|
1008 | self.atom = '*' |
---|
1009 | elif args[0] is not None: |
---|
1010 | self.atom = AtomIdLookup[args[0]].get(str(args[3]),[None,None])[1] |
---|
1011 | if args[1] == '*': |
---|
1012 | self.histogram = '*' |
---|
1013 | else: |
---|
1014 | self.histogram = HistIdLookup.get(str(args[1]),[None,None])[1] |
---|
1015 | self.name = args[2] |
---|
1016 | else: |
---|
1017 | raise Exception("Incorrectly called GSAS-II parameter name") |
---|
1018 | |
---|
1019 | #print "DEBUG: created ",self.phase,self.histogram,self.name,self.atom |
---|
1020 | |
---|
1021 | def __str__(self): |
---|
1022 | return self.varname() |
---|
1023 | |
---|
1024 | def __hash__(self): |
---|
1025 | 'Allow G2VarObj to be a dict key by implementing hashing' |
---|
1026 | return hash(self.varname()) |
---|
1027 | |
---|
1028 | def varname(self,hist=None): |
---|
1029 | '''Formats the GSAS-II variable name as a "traditional" GSAS-II variable |
---|
1030 | string (p:h:<var>:a) or (p:h:<var>) |
---|
1031 | |
---|
1032 | :param str/int hist: if specified, overrides the histogram number |
---|
1033 | with the specified value |
---|
1034 | :returns: the variable name as a str |
---|
1035 | ''' |
---|
1036 | a = "" |
---|
1037 | if self.phase == "*": |
---|
1038 | ph = "*" |
---|
1039 | if self.atom: |
---|
1040 | a = ":" + str(self.atom) |
---|
1041 | else: |
---|
1042 | ph = _lookup(PhaseRanIdLookup,self.phase) |
---|
1043 | if self.atom == '*': |
---|
1044 | a = ':*' |
---|
1045 | elif self.atom: |
---|
1046 | if ":" in str(self.atom): |
---|
1047 | a = ":" + str(self.atom) |
---|
1048 | elif ph in AtomRanIdLookup: |
---|
1049 | a = ":" + AtomRanIdLookup[ph].get(self.atom,'?') |
---|
1050 | else: |
---|
1051 | a = ":?" |
---|
1052 | if hist is not None and self.histogram: |
---|
1053 | hist = str(hist) |
---|
1054 | elif self.histogram == "*": |
---|
1055 | hist = "*" |
---|
1056 | else: |
---|
1057 | hist = _lookup(HistRanIdLookup,self.histogram) |
---|
1058 | s = (ph + ":" + hist + ":" + str(self.name)) + a |
---|
1059 | return s |
---|
1060 | |
---|
1061 | def __repr__(self): |
---|
1062 | '''Return the detailed contents of the object |
---|
1063 | ''' |
---|
1064 | s = "<" |
---|
1065 | if self.phase == '*': |
---|
1066 | s += "Phases: all; " |
---|
1067 | if self.atom is not None: |
---|
1068 | if ":" in str(self.atom): |
---|
1069 | s += "Rigid body" + str(self.atom) + "; " |
---|
1070 | else: |
---|
1071 | s += "Atom #" + str(self.atom) + "; " |
---|
1072 | elif self.phase is not None: |
---|
1073 | ph = _lookup(PhaseRanIdLookup,self.phase) |
---|
1074 | s += "Phase: rId=" + str(self.phase) + " (#"+ ph + "); " |
---|
1075 | if self.atom == '*': |
---|
1076 | s += "Atoms: all; " |
---|
1077 | elif ":" in str(self.atom): |
---|
1078 | s += "Rigid body" + str(self.atom) + "; " |
---|
1079 | elif self.atom is not None: |
---|
1080 | s += "Atom rId=" + str(self.atom) |
---|
1081 | if ph in AtomRanIdLookup: |
---|
1082 | s += " (#" + AtomRanIdLookup[ph].get(self.atom,'?') + "); " |
---|
1083 | else: |
---|
1084 | s += " (#? -- not found!); " |
---|
1085 | if self.histogram == '*': |
---|
1086 | s += "Histograms: all; " |
---|
1087 | elif self.histogram is not None: |
---|
1088 | hist = _lookup(HistRanIdLookup,self.histogram) |
---|
1089 | s += "Histogram: rId=" + str(self.histogram) + " (#"+ hist + "); " |
---|
1090 | s += 'Variable name="' + str(self.name) + '">' |
---|
1091 | return s+" ("+self.varname()+")" |
---|
1092 | |
---|
1093 | def __eq__(self, other): |
---|
1094 | '''Allow comparison of G2VarObj to other G2VarObj objects or strings. |
---|
1095 | If any field is a wildcard ('*') that field matches. |
---|
1096 | ''' |
---|
1097 | if type(other) is str: |
---|
1098 | other = G2VarObj(other) |
---|
1099 | elif type(other) is not G2VarObj: |
---|
1100 | raise Exception("Invalid type ({}) for G2VarObj comparison with {}" |
---|
1101 | .format(type(other),other)) |
---|
1102 | if self.phase != other.phase and self.phase != '*' and other.phase != '*': |
---|
1103 | return False |
---|
1104 | if self.histogram != other.histogram and self.histogram != '*' and other.histogram != '*': |
---|
1105 | return False |
---|
1106 | if self.atom != other.atom and self.atom != '*' and other.atom != '*': |
---|
1107 | return False |
---|
1108 | if self.name != other.name: |
---|
1109 | return False |
---|
1110 | return True |
---|
1111 | |
---|
1112 | def fmtVarByMode(self, seqmode, note, warnmsg): |
---|
1113 | '''Format a parameter object for display. Note that these changes |
---|
1114 | are only temporary and are only shown only when the Constraints |
---|
1115 | data tree is selected. |
---|
1116 | |
---|
1117 | * In a non-sequential refinement or where the mode is 'use-all', the |
---|
1118 | name is converted unchanged to a str |
---|
1119 | * In a sequential refinement when the mode is 'wildcards-only' the |
---|
1120 | name is converted unchanged to a str but a warning is added |
---|
1121 | for non-wildcarded HAP or Histogram parameters |
---|
1122 | * In a sequential refinement or where the mode is 'auto-wildcard', |
---|
1123 | a histogram number is converted to a wildcard (*) and then |
---|
1124 | converted to str |
---|
1125 | |
---|
1126 | :param str mode: the sequential mode (see above) |
---|
1127 | :param str note: value displayed on the line of the constraint/equiv. |
---|
1128 | :param str warnmsg: a message saying the constraint is not used |
---|
1129 | |
---|
1130 | :returns: varname, explain, note, warnmsg (all str values) where: |
---|
1131 | |
---|
1132 | * varname is the parameter expressed as a string, |
---|
1133 | * explain is blank unless there is a warning explanation about |
---|
1134 | the parameter or blank |
---|
1135 | * note is the previous value unless overridden |
---|
1136 | * warnmsg is the previous value unless overridden |
---|
1137 | ''' |
---|
1138 | explain = '' |
---|
1139 | s = self.varname() |
---|
1140 | if seqmode == 'auto-wildcard': |
---|
1141 | if self.histogram: s = self.varname('*') |
---|
1142 | elif seqmode == 'wildcards-only' and self.histogram: |
---|
1143 | if self.histogram != '*': |
---|
1144 | warnmsg = 'Ignored due to use of a non-wildcarded histogram number' |
---|
1145 | note = 'Ignored' |
---|
1146 | explain = '\nIgnoring: '+self.varname()+' does not contain a wildcard.\n' |
---|
1147 | elif seqmode != 'use-all' and seqmode != 'wildcards-only': |
---|
1148 | print('Unexpected mode',seqmode,' in fmtVarByMode') |
---|
1149 | return s,explain,note,warnmsg |
---|
1150 | |
---|
1151 | def _show(self): |
---|
1152 | 'For testing, shows the current lookup table' |
---|
1153 | print ('phases'+ self.IDdict['phases']) |
---|
1154 | print ('hists'+ self.IDdict['hists']) |
---|
1155 | print ('atomDict'+ self.IDdict['atoms']) |
---|
1156 | |
---|
1157 | #========================================================================== |
---|
1158 | def SetDefaultSample(): |
---|
1159 | 'Fills in default items for the Sample dictionary for Debye-Scherrer & SASD' |
---|
1160 | return { |
---|
1161 | 'InstrName':'', |
---|
1162 | 'ranId':ran.randint(0,sys.maxsize), |
---|
1163 | 'Scale':[1.0,True],'Type':'Debye-Scherrer','Absorption':[0.0,False], |
---|
1164 | 'DisplaceX':[0.0,False],'DisplaceY':[0.0,False], |
---|
1165 | 'Temperature':300.,'Pressure':0.1,'Time':0.0, |
---|
1166 | 'FreePrm1':0.,'FreePrm2':0.,'FreePrm3':0., |
---|
1167 | 'Gonio. radius':200.0, |
---|
1168 | 'Omega':0.0,'Chi':0.0,'Phi':0.0,'Azimuth':0.0, |
---|
1169 | #SASD items |
---|
1170 | 'Materials':[{'Name':'vacuum','VolFrac':1.0,},{'Name':'vacuum','VolFrac':0.0,}], |
---|
1171 | 'Thick':1.0,'Contrast':[0.0,0.0], #contrast & anomalous contrast |
---|
1172 | 'Trans':1.0, #measured transmission |
---|
1173 | 'SlitLen':0.0, #Slit length - in Q(A-1) |
---|
1174 | } |
---|
1175 | ###################################################################### |
---|
1176 | class ImportBaseclass(object): |
---|
1177 | '''Defines a base class for the reading of input files (diffraction |
---|
1178 | data, coordinates,...). See :ref:`Writing a Import Routine<import_routines>` |
---|
1179 | for an explanation on how to use a subclass of this class. |
---|
1180 | ''' |
---|
1181 | class ImportException(Exception): |
---|
1182 | '''Defines an Exception that is used when an import routine hits an expected error, |
---|
1183 | usually in .Reader. |
---|
1184 | |
---|
1185 | Good practice is that the Reader should define a value in self.errors that |
---|
1186 | tells the user some information about what is wrong with their file. |
---|
1187 | ''' |
---|
1188 | pass |
---|
1189 | |
---|
1190 | UseReader = True # in __init__ set value of self.UseReader to False to skip use of current importer |
---|
1191 | def __init__(self,formatName,longFormatName=None, |
---|
1192 | extensionlist=[],strictExtension=False,): |
---|
1193 | self.formatName = formatName # short string naming file type |
---|
1194 | if longFormatName: # longer string naming file type |
---|
1195 | self.longFormatName = longFormatName |
---|
1196 | else: |
---|
1197 | self.longFormatName = formatName |
---|
1198 | # define extensions that are allowed for the file type |
---|
1199 | # for windows, remove any extensions that are duplicate, as case is ignored |
---|
1200 | if sys.platform == 'windows' and extensionlist: |
---|
1201 | extensionlist = list(set([s.lower() for s in extensionlist])) |
---|
1202 | self.extensionlist = extensionlist |
---|
1203 | # If strictExtension is True, the file will not be read, unless |
---|
1204 | # the extension matches one in the extensionlist |
---|
1205 | self.strictExtension = strictExtension |
---|
1206 | self.errors = '' |
---|
1207 | self.warnings = '' |
---|
1208 | self.SciPy = False #image reader needed scipy |
---|
1209 | # used for readers that will use multiple passes to read |
---|
1210 | # more than one data block |
---|
1211 | self.repeat = False |
---|
1212 | self.selections = [] |
---|
1213 | self.repeatcount = 0 |
---|
1214 | self.readfilename = '?' |
---|
1215 | self.scriptable = False |
---|
1216 | #print 'created',self.__class__ |
---|
1217 | |
---|
1218 | def ReInitialize(self): |
---|
1219 | 'Reinitialize the Reader to initial settings' |
---|
1220 | self.errors = '' |
---|
1221 | self.warnings = '' |
---|
1222 | self.SciPy = False #image reader needed scipy |
---|
1223 | self.repeat = False |
---|
1224 | self.repeatcount = 0 |
---|
1225 | self.readfilename = '?' |
---|
1226 | |
---|
1227 | |
---|
1228 | # def Reader(self, filename, filepointer, ParentFrame=None, **unused): |
---|
1229 | # '''This method must be supplied in the child class to read the file. |
---|
1230 | # if the read fails either return False or raise an Exception |
---|
1231 | # preferably of type ImportException. |
---|
1232 | # ''' |
---|
1233 | # #start reading |
---|
1234 | # raise ImportException("Error occurred while...") |
---|
1235 | # self.errors += "Hint for user on why the error occur |
---|
1236 | # return False # if an error occurs |
---|
1237 | # return True # if read OK |
---|
1238 | |
---|
1239 | def ExtensionValidator(self, filename): |
---|
1240 | '''This methods checks if the file has the correct extension |
---|
1241 | |
---|
1242 | :returns: |
---|
1243 | |
---|
1244 | * False if this filename will not be supported by this reader (only |
---|
1245 | when strictExtension is True) |
---|
1246 | * True if the extension matches the list supplied by the reader |
---|
1247 | * None if the reader allows un-registered extensions |
---|
1248 | |
---|
1249 | ''' |
---|
1250 | if filename: |
---|
1251 | ext = os.path.splitext(filename)[1] |
---|
1252 | if not ext and self.strictExtension: return False |
---|
1253 | for ext in self.extensionlist: |
---|
1254 | if sys.platform == 'windows': |
---|
1255 | if filename.lower().endswith(ext): return True |
---|
1256 | else: |
---|
1257 | if filename.endswith(ext): return True |
---|
1258 | if self.strictExtension: |
---|
1259 | return False |
---|
1260 | else: |
---|
1261 | return None |
---|
1262 | |
---|
1263 | def ContentsValidator(self, filename): |
---|
1264 | '''This routine will attempt to determine if the file can be read |
---|
1265 | with the current format. |
---|
1266 | This will typically be overridden with a method that |
---|
1267 | takes a quick scan of [some of] |
---|
1268 | the file contents to do a "sanity" check if the file |
---|
1269 | appears to match the selected format. |
---|
1270 | the file must be opened here with the correct format (binary/text) |
---|
1271 | ''' |
---|
1272 | #filepointer.seek(0) # rewind the file pointer |
---|
1273 | return True |
---|
1274 | |
---|
1275 | def CIFValidator(self, filepointer): |
---|
1276 | '''A :meth:`ContentsValidator` for use to validate CIF files. |
---|
1277 | ''' |
---|
1278 | filepointer.seek(0) |
---|
1279 | for i,l in enumerate(filepointer): |
---|
1280 | if i >= 1000: return True |
---|
1281 | '''Encountered only blank lines or comments in first 1000 |
---|
1282 | lines. This is unlikely, but assume it is CIF anyway, since we are |
---|
1283 | even less likely to find a file with nothing but hashes and |
---|
1284 | blank lines''' |
---|
1285 | line = l.strip() |
---|
1286 | if len(line) == 0: # ignore blank lines |
---|
1287 | continue |
---|
1288 | elif line.startswith('#'): # ignore comments |
---|
1289 | continue |
---|
1290 | elif line.startswith('data_'): # on the right track, accept this file |
---|
1291 | return True |
---|
1292 | else: # found something invalid |
---|
1293 | self.errors = 'line '+str(i+1)+' contains unexpected data:\n' |
---|
1294 | if all([ord(c) < 128 and ord(c) != 0 for c in str(l)]): # show only if ASCII |
---|
1295 | self.errors += ' '+str(l) |
---|
1296 | else: |
---|
1297 | self.errors += ' (binary)' |
---|
1298 | self.errors += '\n Note: a CIF should only have blank lines or comments before' |
---|
1299 | self.errors += '\n a data_ statement begins a block.' |
---|
1300 | return False |
---|
1301 | |
---|
1302 | ###################################################################### |
---|
1303 | class ImportPhase(ImportBaseclass): |
---|
1304 | '''Defines a base class for the reading of files with coordinates |
---|
1305 | |
---|
1306 | Objects constructed that subclass this (in import/G2phase_*.py etc.) will be used |
---|
1307 | in :meth:`GSASIIdataGUI.GSASII.OnImportPhase` and in |
---|
1308 | :func:`GSASIIscriptable.import_generic`. |
---|
1309 | See :ref:`Writing a Import Routine<import_routines>` |
---|
1310 | for an explanation on how to use this class. |
---|
1311 | |
---|
1312 | ''' |
---|
1313 | def __init__(self,formatName,longFormatName=None,extensionlist=[], |
---|
1314 | strictExtension=False,): |
---|
1315 | # call parent __init__ |
---|
1316 | ImportBaseclass.__init__(self,formatName,longFormatName, |
---|
1317 | extensionlist,strictExtension) |
---|
1318 | self.Phase = None # a phase must be created with G2IO.SetNewPhase in the Reader |
---|
1319 | self.SymOps = {} # specified when symmetry ops are in file (e.g. CIF) |
---|
1320 | self.Constraints = None |
---|
1321 | |
---|
1322 | ###################################################################### |
---|
1323 | class ImportStructFactor(ImportBaseclass): |
---|
1324 | '''Defines a base class for the reading of files with tables |
---|
1325 | of structure factors. |
---|
1326 | |
---|
1327 | Structure factors are read with a call to :meth:`GSASIIdataGUI.GSASII.OnImportSfact` |
---|
1328 | which in turn calls :meth:`GSASIIdataGUI.GSASII.OnImportGeneric`, which calls |
---|
1329 | methods :meth:`ExtensionValidator`, :meth:`ContentsValidator` and |
---|
1330 | :meth:`Reader`. |
---|
1331 | |
---|
1332 | See :ref:`Writing a Import Routine<import_routines>` |
---|
1333 | for an explanation on how to use import classes in general. The specifics |
---|
1334 | for reading a structure factor histogram require that |
---|
1335 | the ``Reader()`` routine in the import |
---|
1336 | class need to do only a few things: It |
---|
1337 | should load :attr:`RefDict` item ``'RefList'`` with the reflection list, |
---|
1338 | and set :attr:`Parameters` with the instrument parameters |
---|
1339 | (initialized with :meth:`InitParameters` and set with :meth:`UpdateParameters`). |
---|
1340 | ''' |
---|
1341 | def __init__(self,formatName,longFormatName=None,extensionlist=[], |
---|
1342 | strictExtension=False,): |
---|
1343 | ImportBaseclass.__init__(self,formatName,longFormatName, |
---|
1344 | extensionlist,strictExtension) |
---|
1345 | |
---|
1346 | # define contents of Structure Factor entry |
---|
1347 | self.Parameters = [] |
---|
1348 | 'self.Parameters is a list with two dicts for data parameter settings' |
---|
1349 | self.InitParameters() |
---|
1350 | self.RefDict = {'RefList':[],'FF':{},'Super':0} |
---|
1351 | self.Banks = [] #for multi bank data (usually TOF) |
---|
1352 | '''self.RefDict is a dict containing the reflection information, as read from the file. |
---|
1353 | Item 'RefList' contains the reflection information. See the |
---|
1354 | :ref:`Single Crystal Reflection Data Structure<XtalRefl_table>` |
---|
1355 | for the contents of each row. Dict element 'FF' |
---|
1356 | contains the form factor values for each element type; if this entry |
---|
1357 | is left as initialized (an empty list) it will be initialized as needed later. |
---|
1358 | ''' |
---|
1359 | def ReInitialize(self): |
---|
1360 | 'Reinitialize the Reader to initial settings' |
---|
1361 | ImportBaseclass.ReInitialize(self) |
---|
1362 | self.InitParameters() |
---|
1363 | self.Banks = [] #for multi bank data (usually TOF) |
---|
1364 | self.RefDict = {'RefList':[],'FF':{},'Super':0} |
---|
1365 | |
---|
1366 | def InitParameters(self): |
---|
1367 | 'initialize the instrument parameters structure' |
---|
1368 | Lambda = 0.70926 |
---|
1369 | HistType = 'SXC' |
---|
1370 | self.Parameters = [{'Type':[HistType,HistType], # create the structure |
---|
1371 | 'Lam':[Lambda,Lambda] |
---|
1372 | }, {}] |
---|
1373 | 'Parameters is a list with two dicts for data parameter settings' |
---|
1374 | |
---|
1375 | def UpdateParameters(self,Type=None,Wave=None): |
---|
1376 | 'Revise the instrument parameters' |
---|
1377 | if Type is not None: |
---|
1378 | self.Parameters[0]['Type'] = [Type,Type] |
---|
1379 | if Wave is not None: |
---|
1380 | self.Parameters[0]['Lam'] = [Wave,Wave] |
---|
1381 | |
---|
1382 | ###################################################################### |
---|
1383 | class ImportPowderData(ImportBaseclass): |
---|
1384 | '''Defines a base class for the reading of files with powder data. |
---|
1385 | |
---|
1386 | Objects constructed that subclass this (in import/G2pwd_*.py etc.) will be used |
---|
1387 | in :meth:`GSASIIdataGUI.GSASII.OnImportPowder` and in |
---|
1388 | :func:`GSASIIscriptable.import_generic`. |
---|
1389 | See :ref:`Writing a Import Routine<import_routines>` |
---|
1390 | for an explanation on how to use this class. |
---|
1391 | ''' |
---|
1392 | def __init__(self,formatName,longFormatName=None, |
---|
1393 | extensionlist=[],strictExtension=False,): |
---|
1394 | ImportBaseclass.__init__(self,formatName,longFormatName, |
---|
1395 | extensionlist,strictExtension) |
---|
1396 | self.clockWd = None # used in TOF |
---|
1397 | self.ReInitialize() |
---|
1398 | |
---|
1399 | def ReInitialize(self): |
---|
1400 | 'Reinitialize the Reader to initial settings' |
---|
1401 | ImportBaseclass.ReInitialize(self) |
---|
1402 | self.powderentry = ['',None,None] # (filename,Pos,Bank) |
---|
1403 | self.powderdata = [] # Powder dataset |
---|
1404 | '''A powder data set is a list with items [x,y,w,yc,yb,yd]: |
---|
1405 | np.array(x), # x-axis values |
---|
1406 | np.array(y), # powder pattern intensities |
---|
1407 | np.array(w), # 1/sig(intensity)^2 values (weights) |
---|
1408 | np.array(yc), # calc. intensities (zero) |
---|
1409 | np.array(yb), # calc. background (zero) |
---|
1410 | np.array(yd), # obs-calc profiles |
---|
1411 | ''' |
---|
1412 | self.comments = [] |
---|
1413 | self.idstring = '' |
---|
1414 | self.Sample = SetDefaultSample() # default sample parameters |
---|
1415 | self.Controls = {} # items to be placed in top-level Controls |
---|
1416 | self.GSAS = None # used in TOF |
---|
1417 | self.repeat_instparm = True # Should a parm file be |
---|
1418 | # used for multiple histograms? |
---|
1419 | self.instparm = None # name hint from file of instparm to use |
---|
1420 | self.instfile = '' # full path name to instrument parameter file |
---|
1421 | self.instbank = '' # inst parm bank number |
---|
1422 | self.instmsg = '' # a label that gets printed to show |
---|
1423 | # where instrument parameters are from |
---|
1424 | self.numbanks = 1 |
---|
1425 | self.instdict = {} # place items here that will be transferred to the instrument parameters |
---|
1426 | self.pwdparms = {} # place parameters that are transferred directly to the tree |
---|
1427 | # here (typically from an existing GPX file) |
---|
1428 | ###################################################################### |
---|
1429 | class ImportSmallAngleData(ImportBaseclass): |
---|
1430 | '''Defines a base class for the reading of files with small angle data. |
---|
1431 | See :ref:`Writing a Import Routine<import_routines>` |
---|
1432 | for an explanation on how to use this class. |
---|
1433 | ''' |
---|
1434 | def __init__(self,formatName,longFormatName=None,extensionlist=[], |
---|
1435 | strictExtension=False,): |
---|
1436 | |
---|
1437 | ImportBaseclass.__init__(self,formatName,longFormatName,extensionlist, |
---|
1438 | strictExtension) |
---|
1439 | self.ReInitialize() |
---|
1440 | |
---|
1441 | def ReInitialize(self): |
---|
1442 | 'Reinitialize the Reader to initial settings' |
---|
1443 | ImportBaseclass.ReInitialize(self) |
---|
1444 | self.smallangleentry = ['',None,None] # (filename,Pos,Bank) |
---|
1445 | self.smallangledata = [] # SASD dataset |
---|
1446 | '''A small angle data set is a list with items [x,y,w,yc,yd]: |
---|
1447 | np.array(x), # x-axis values |
---|
1448 | np.array(y), # powder pattern intensities |
---|
1449 | np.array(w), # 1/sig(intensity)^2 values (weights) |
---|
1450 | np.array(yc), # calc. intensities (zero) |
---|
1451 | np.array(yd), # obs-calc profiles |
---|
1452 | np.array(yb), # preset bkg |
---|
1453 | ''' |
---|
1454 | self.comments = [] |
---|
1455 | self.idstring = '' |
---|
1456 | self.Sample = SetDefaultSample() |
---|
1457 | self.GSAS = None # used in TOF |
---|
1458 | self.clockWd = None # used in TOF |
---|
1459 | self.numbanks = 1 |
---|
1460 | self.instdict = {} # place items here that will be transferred to the instrument parameters |
---|
1461 | |
---|
1462 | ###################################################################### |
---|
1463 | class ImportReflectometryData(ImportBaseclass): |
---|
1464 | '''Defines a base class for the reading of files with reflectometry data. |
---|
1465 | See :ref:`Writing a Import Routine<import_routines>` |
---|
1466 | for an explanation on how to use this class. |
---|
1467 | ''' |
---|
1468 | def __init__(self,formatName,longFormatName=None,extensionlist=[], |
---|
1469 | strictExtension=False,): |
---|
1470 | |
---|
1471 | ImportBaseclass.__init__(self,formatName,longFormatName,extensionlist, |
---|
1472 | strictExtension) |
---|
1473 | self.ReInitialize() |
---|
1474 | |
---|
1475 | def ReInitialize(self): |
---|
1476 | 'Reinitialize the Reader to initial settings' |
---|
1477 | ImportBaseclass.ReInitialize(self) |
---|
1478 | self.reflectometryentry = ['',None,None] # (filename,Pos,Bank) |
---|
1479 | self.reflectometrydata = [] # SASD dataset |
---|
1480 | '''A small angle data set is a list with items [x,y,w,yc,yd]: |
---|
1481 | np.array(x), # x-axis values |
---|
1482 | np.array(y), # powder pattern intensities |
---|
1483 | np.array(w), # 1/sig(intensity)^2 values (weights) |
---|
1484 | np.array(yc), # calc. intensities (zero) |
---|
1485 | np.array(yd), # obs-calc profiles |
---|
1486 | np.array(yb), # preset bkg |
---|
1487 | ''' |
---|
1488 | self.comments = [] |
---|
1489 | self.idstring = '' |
---|
1490 | self.Sample = SetDefaultSample() |
---|
1491 | self.GSAS = None # used in TOF |
---|
1492 | self.clockWd = None # used in TOF |
---|
1493 | self.numbanks = 1 |
---|
1494 | self.instdict = {} # place items here that will be transferred to the instrument parameters |
---|
1495 | |
---|
1496 | ###################################################################### |
---|
1497 | class ImportPDFData(ImportBaseclass): |
---|
1498 | '''Defines a base class for the reading of files with PDF G(R) data. |
---|
1499 | See :ref:`Writing a Import Routine<import_routines>` |
---|
1500 | for an explanation on how to use this class. |
---|
1501 | ''' |
---|
1502 | def __init__(self,formatName,longFormatName=None,extensionlist=[], |
---|
1503 | strictExtension=False,): |
---|
1504 | |
---|
1505 | ImportBaseclass.__init__(self,formatName,longFormatName,extensionlist, |
---|
1506 | strictExtension) |
---|
1507 | self.ReInitialize() |
---|
1508 | |
---|
1509 | def ReInitialize(self): |
---|
1510 | 'Reinitialize the Reader to initial settings' |
---|
1511 | ImportBaseclass.ReInitialize(self) |
---|
1512 | self.pdfentry = ['',None,None] # (filename,Pos,Bank) |
---|
1513 | self.pdfdata = [] # PDF G(R) dataset |
---|
1514 | '''A pdf g(r) data set is a list with items [x,y]: |
---|
1515 | np.array(x), # r-axis values |
---|
1516 | np.array(y), # pdf g(r) |
---|
1517 | ''' |
---|
1518 | self.comments = [] |
---|
1519 | self.idstring = '' |
---|
1520 | self.numbanks = 1 |
---|
1521 | |
---|
1522 | ###################################################################### |
---|
1523 | class ImportImage(ImportBaseclass): |
---|
1524 | '''Defines a base class for the reading of images |
---|
1525 | |
---|
1526 | Images are read in only these places: |
---|
1527 | |
---|
1528 | * Initial reading is typically done from a menu item |
---|
1529 | with a call to :meth:`GSASIIdataGUI.GSASII.OnImportImage` |
---|
1530 | which in turn calls :meth:`GSASIIdataGUI.GSASII.OnImportGeneric`. That calls |
---|
1531 | methods :meth:`ExtensionValidator`, :meth:`ContentsValidator` and |
---|
1532 | :meth:`Reader`. This returns a list of reader objects for each read image. |
---|
1533 | Also used in :func:`GSASIIscriptable.import_generic`. |
---|
1534 | |
---|
1535 | * Images are read alternatively in :func:`GSASIIIO.ReadImages`, which puts image info |
---|
1536 | directly into the data tree. |
---|
1537 | |
---|
1538 | * Images are reloaded with :func:`GSASIIIO.GetImageData`. |
---|
1539 | |
---|
1540 | When reading an image, the ``Reader()`` routine in the ImportImage class |
---|
1541 | should set: |
---|
1542 | |
---|
1543 | * :attr:`Comments`: a list of strings (str), |
---|
1544 | * :attr:`Npix`: the number of pixels in the image (int), |
---|
1545 | * :attr:`Image`: the actual image as a numpy array (np.array) |
---|
1546 | * :attr:`Data`: a dict defining image parameters (dict). Within this dict the following |
---|
1547 | data items are needed: |
---|
1548 | |
---|
1549 | * 'pixelSize': size of each pixel in microns (such as ``[200.,200.]``. |
---|
1550 | * 'wavelength': wavelength in :math:`\\AA`. |
---|
1551 | * 'distance': distance of detector from sample in cm. |
---|
1552 | * 'center': uncalibrated center of beam on detector (such as ``[204.8,204.8]``. |
---|
1553 | * 'size': size of image (such as ``[2048,2048]``). |
---|
1554 | * 'ImageTag': image number or other keyword used to retrieve image from |
---|
1555 | a multi-image data file (defaults to ``1`` if not specified). |
---|
1556 | * 'sumfile': holds sum image file name if a sum was produced from a multi image file |
---|
1557 | |
---|
1558 | optional data items: |
---|
1559 | |
---|
1560 | * :attr:`repeat`: set to True if there are additional images to |
---|
1561 | read in the file, False otherwise |
---|
1562 | * :attr:`repeatcount`: set to the number of the image. |
---|
1563 | |
---|
1564 | Note that the above is initialized with :meth:`InitParameters`. |
---|
1565 | (Also see :ref:`Writing a Import Routine<import_routines>` |
---|
1566 | for an explanation on how to use import classes in general.) |
---|
1567 | ''' |
---|
1568 | def __init__(self,formatName,longFormatName=None,extensionlist=[], |
---|
1569 | strictExtension=False,): |
---|
1570 | ImportBaseclass.__init__(self,formatName,longFormatName, |
---|
1571 | extensionlist,strictExtension) |
---|
1572 | self.InitParameters() |
---|
1573 | |
---|
1574 | def ReInitialize(self): |
---|
1575 | 'Reinitialize the Reader to initial settings -- not used at present' |
---|
1576 | ImportBaseclass.ReInitialize(self) |
---|
1577 | self.InitParameters() |
---|
1578 | |
---|
1579 | def InitParameters(self): |
---|
1580 | 'initialize the instrument parameters structure' |
---|
1581 | self.Comments = ['No comments'] |
---|
1582 | self.Data = {'samplechangerpos':0.0,'det2theta':0.0,'Gain map':''} |
---|
1583 | self.Npix = 0 |
---|
1584 | self.Image = None |
---|
1585 | self.repeat = False |
---|
1586 | self.repeatcount = 1 |
---|
1587 | self.sumfile = '' |
---|
1588 | |
---|
1589 | def LoadImage(self,ParentFrame,imagefile,imagetag=None): |
---|
1590 | '''Optionally, call this after reading in an image to load it into the tree. |
---|
1591 | This saves time by preventing a reread of the same information. |
---|
1592 | ''' |
---|
1593 | if ParentFrame: |
---|
1594 | ParentFrame.ImageZ = self.Image # store the image for plotting |
---|
1595 | ParentFrame.oldImagefile = imagefile # save the name of the last image file read |
---|
1596 | ParentFrame.oldImageTag = imagetag # save the tag of the last image file read |
---|
1597 | |
---|
1598 | ################################################################################################# |
---|
1599 | # shortcut routines |
---|
1600 | exp = np.exp |
---|
1601 | sind = sin = s = lambda x: np.sin(x*np.pi/180.) |
---|
1602 | cosd = cos = c = lambda x: np.cos(x*np.pi/180.) |
---|
1603 | tand = tan = t = lambda x: np.tan(x*np.pi/180.) |
---|
1604 | sqrt = sq = lambda x: np.sqrt(x) |
---|
1605 | pi = lambda: np.pi |
---|
1606 | |
---|
1607 | def FindFunction(f): |
---|
1608 | '''Find the object corresponding to function f |
---|
1609 | |
---|
1610 | :param str f: a function name such as 'numpy.exp' |
---|
1611 | :returns: (pkgdict,pkgobj) where pkgdict contains a dict |
---|
1612 | that defines the package location(s) and where pkgobj |
---|
1613 | defines the object associated with the function. |
---|
1614 | If the function is not found, pkgobj is None. |
---|
1615 | ''' |
---|
1616 | df = f.split('.') |
---|
1617 | pkgdict = {} |
---|
1618 | # no listed module name, try in current namespace |
---|
1619 | if len(df) == 1: |
---|
1620 | try: |
---|
1621 | fxnobj = eval(f) |
---|
1622 | return pkgdict,fxnobj |
---|
1623 | except (AttributeError, NameError): |
---|
1624 | return None,None |
---|
1625 | |
---|
1626 | # includes a package, see if package is already imported |
---|
1627 | pkgnam = '.'.join(df[:-1]) |
---|
1628 | try: |
---|
1629 | fxnobj = eval(f) |
---|
1630 | pkgdict[pkgnam] = eval(pkgnam) |
---|
1631 | return pkgdict,fxnobj |
---|
1632 | except (AttributeError, NameError): |
---|
1633 | pass |
---|
1634 | # package not yet imported, so let's try |
---|
1635 | if '.' not in sys.path: sys.path.append('.') |
---|
1636 | pkgnam = '.'.join(df[:-1]) |
---|
1637 | #for pkg in f.split('.')[:-1]: # if needed, descend down the tree |
---|
1638 | # if pkgname: |
---|
1639 | # pkgname += '.' + pkg |
---|
1640 | # else: |
---|
1641 | # pkgname = pkg |
---|
1642 | try: |
---|
1643 | exec('import '+pkgnam) |
---|
1644 | pkgdict[pkgnam] = eval(pkgnam) |
---|
1645 | fxnobj = eval(f) |
---|
1646 | except Exception as msg: |
---|
1647 | print('load of '+pkgnam+' failed with error='+str(msg)) |
---|
1648 | return {},None |
---|
1649 | # can we access the function? I am not exactly sure what |
---|
1650 | # I intended this to test originally (BHT) |
---|
1651 | try: |
---|
1652 | fxnobj = eval(f,globals(),pkgdict) |
---|
1653 | return pkgdict,fxnobj |
---|
1654 | except Exception as msg: |
---|
1655 | print('call to',f,' failed with error=',str(msg)) |
---|
1656 | return None,None # not found |
---|
1657 | |
---|
1658 | class ExpressionObj(object): |
---|
1659 | '''Defines an object with a user-defined expression, to be used for |
---|
1660 | secondary fits or restraints. Object is created null, but is changed |
---|
1661 | using :meth:`LoadExpression`. This contains only the minimum |
---|
1662 | information that needs to be stored to save and load the expression |
---|
1663 | and how it is mapped to GSAS-II variables. |
---|
1664 | ''' |
---|
1665 | def __init__(self): |
---|
1666 | self.expression = '' |
---|
1667 | 'The expression as a text string' |
---|
1668 | self.assgnVars = {} |
---|
1669 | '''A dict where keys are label names in the expression mapping to a GSAS-II |
---|
1670 | variable. The value a G2 variable name. |
---|
1671 | Note that the G2 variable name may contain a wild-card and correspond to |
---|
1672 | multiple values. |
---|
1673 | ''' |
---|
1674 | self.freeVars = {} |
---|
1675 | '''A dict where keys are label names in the expression mapping to a free |
---|
1676 | parameter. The value is a list with: |
---|
1677 | |
---|
1678 | * a name assigned to the parameter |
---|
1679 | * a value for to the parameter and |
---|
1680 | * a flag to determine if the variable is refined. |
---|
1681 | ''' |
---|
1682 | self.depVar = None |
---|
1683 | |
---|
1684 | self.lastError = ('','') |
---|
1685 | '''Shows last encountered error in processing expression |
---|
1686 | (list of 1-3 str values)''' |
---|
1687 | |
---|
1688 | self.distance_dict = None # to be used for defining atom phase/symmetry info |
---|
1689 | self.distance_atoms = None # to be used for defining atom distances |
---|
1690 | |
---|
1691 | def LoadExpression(self,expr,exprVarLst,varSelect,varName,varValue,varRefflag): |
---|
1692 | '''Load the expression and associated settings into the object. Raises |
---|
1693 | an exception if the expression is not parsed, if not all functions |
---|
1694 | are defined or if not all needed parameter labels in the expression |
---|
1695 | are defined. |
---|
1696 | |
---|
1697 | This will not test if the variable referenced in these definitions |
---|
1698 | are actually in the parameter dictionary. This is checked when the |
---|
1699 | computation for the expression is done in :meth:`SetupCalc`. |
---|
1700 | |
---|
1701 | :param str expr: the expression |
---|
1702 | :param list exprVarLst: parameter labels found in the expression |
---|
1703 | :param dict varSelect: this will be 0 for Free parameters |
---|
1704 | and non-zero for expression labels linked to G2 variables. |
---|
1705 | :param dict varName: Defines a name (str) associated with each free parameter |
---|
1706 | :param dict varValue: Defines a value (float) associated with each free parameter |
---|
1707 | :param dict varRefflag: Defines a refinement flag (bool) |
---|
1708 | associated with each free parameter |
---|
1709 | ''' |
---|
1710 | self.expression = expr |
---|
1711 | self.compiledExpr = None |
---|
1712 | self.freeVars = {} |
---|
1713 | self.assgnVars = {} |
---|
1714 | for v in exprVarLst: |
---|
1715 | if varSelect[v] == 0: |
---|
1716 | self.freeVars[v] = [ |
---|
1717 | varName.get(v), |
---|
1718 | varValue.get(v), |
---|
1719 | varRefflag.get(v), |
---|
1720 | ] |
---|
1721 | else: |
---|
1722 | self.assgnVars[v] = varName[v] |
---|
1723 | self.CheckVars() |
---|
1724 | |
---|
1725 | def EditExpression(self,exprVarLst,varSelect,varName,varValue,varRefflag): |
---|
1726 | '''Load the expression and associated settings from the object into |
---|
1727 | arrays used for editing. |
---|
1728 | |
---|
1729 | :param list exprVarLst: parameter labels found in the expression |
---|
1730 | :param dict varSelect: this will be 0 for Free parameters |
---|
1731 | and non-zero for expression labels linked to G2 variables. |
---|
1732 | :param dict varName: Defines a name (str) associated with each free parameter |
---|
1733 | :param dict varValue: Defines a value (float) associated with each free parameter |
---|
1734 | :param dict varRefflag: Defines a refinement flag (bool) |
---|
1735 | associated with each free parameter |
---|
1736 | |
---|
1737 | :returns: the expression as a str |
---|
1738 | ''' |
---|
1739 | for v in self.freeVars: |
---|
1740 | varSelect[v] = 0 |
---|
1741 | varName[v] = self.freeVars[v][0] |
---|
1742 | varValue[v] = self.freeVars[v][1] |
---|
1743 | varRefflag[v] = self.freeVars[v][2] |
---|
1744 | for v in self.assgnVars: |
---|
1745 | varSelect[v] = 1 |
---|
1746 | varName[v] = self.assgnVars[v] |
---|
1747 | return self.expression |
---|
1748 | |
---|
1749 | def GetVaried(self): |
---|
1750 | 'Returns the names of the free parameters that will be refined' |
---|
1751 | return ["::"+self.freeVars[v][0] for v in self.freeVars if self.freeVars[v][2]] |
---|
1752 | |
---|
1753 | def GetVariedVarVal(self): |
---|
1754 | 'Returns the names and values of the free parameters that will be refined' |
---|
1755 | return [("::"+self.freeVars[v][0],self.freeVars[v][1]) for v in self.freeVars if self.freeVars[v][2]] |
---|
1756 | |
---|
1757 | def UpdateVariedVars(self,varyList,values): |
---|
1758 | 'Updates values for the free parameters (after a refinement); only updates refined vars' |
---|
1759 | for v in self.freeVars: |
---|
1760 | if not self.freeVars[v][2]: continue |
---|
1761 | if "::"+self.freeVars[v][0] not in varyList: continue |
---|
1762 | indx = list(varyList).index("::"+self.freeVars[v][0]) |
---|
1763 | self.freeVars[v][1] = values[indx] |
---|
1764 | |
---|
1765 | def GetIndependentVars(self): |
---|
1766 | 'Returns the names of the required independent parameters used in expression' |
---|
1767 | return [self.assgnVars[v] for v in self.assgnVars] |
---|
1768 | |
---|
1769 | def CheckVars(self): |
---|
1770 | '''Check that the expression can be parsed, all functions are |
---|
1771 | defined and that input loaded into the object is internally |
---|
1772 | consistent. If not an Exception is raised. |
---|
1773 | |
---|
1774 | :returns: a dict with references to packages needed to |
---|
1775 | find functions referenced in the expression. |
---|
1776 | ''' |
---|
1777 | ret = self.ParseExpression(self.expression) |
---|
1778 | if not ret: |
---|
1779 | raise Exception("Expression parse error") |
---|
1780 | exprLblList,fxnpkgdict = ret |
---|
1781 | # check each var used in expression is defined |
---|
1782 | defined = list(self.assgnVars.keys()) + list(self.freeVars.keys()) |
---|
1783 | notfound = [] |
---|
1784 | for var in exprLblList: |
---|
1785 | if var not in defined: |
---|
1786 | notfound.append(var) |
---|
1787 | if notfound: |
---|
1788 | msg = 'Not all variables defined' |
---|
1789 | msg1 = 'The following variables were not defined: ' |
---|
1790 | msg2 = '' |
---|
1791 | for var in notfound: |
---|
1792 | if msg: msg += ', ' |
---|
1793 | msg += var |
---|
1794 | self.lastError = (msg1,' '+msg2) |
---|
1795 | raise Exception(msg) |
---|
1796 | return fxnpkgdict |
---|
1797 | |
---|
1798 | def ParseExpression(self,expr): |
---|
1799 | '''Parse an expression and return a dict of called functions and |
---|
1800 | the variables used in the expression. Returns None in case an error |
---|
1801 | is encountered. If packages are referenced in functions, they are loaded |
---|
1802 | and the functions are looked up into the modules global |
---|
1803 | workspace. |
---|
1804 | |
---|
1805 | Note that no changes are made to the object other than |
---|
1806 | saving an error message, so that this can be used for testing prior |
---|
1807 | to the save. |
---|
1808 | |
---|
1809 | :returns: a list of used variables |
---|
1810 | ''' |
---|
1811 | self.lastError = ('','') |
---|
1812 | import ast |
---|
1813 | def ASTtransverse(node,fxn=False): |
---|
1814 | '''Transverse a AST-parsed expresson, compiling a list of variables |
---|
1815 | referenced in the expression. This routine is used recursively. |
---|
1816 | |
---|
1817 | :returns: varlist,fxnlist where |
---|
1818 | varlist is a list of referenced variable names and |
---|
1819 | fxnlist is a list of used functions |
---|
1820 | ''' |
---|
1821 | varlist = [] |
---|
1822 | fxnlist = [] |
---|
1823 | if isinstance(node, list): |
---|
1824 | for b in node: |
---|
1825 | v,f = ASTtransverse(b,fxn) |
---|
1826 | varlist += v |
---|
1827 | fxnlist += f |
---|
1828 | elif isinstance(node, ast.AST): |
---|
1829 | for a, b in ast.iter_fields(node): |
---|
1830 | if isinstance(b, ast.AST): |
---|
1831 | if a == 'func': |
---|
1832 | fxnlist += ['.'.join(ASTtransverse(b,True)[0])] |
---|
1833 | continue |
---|
1834 | v,f = ASTtransverse(b,fxn) |
---|
1835 | varlist += v |
---|
1836 | fxnlist += f |
---|
1837 | elif isinstance(b, list): |
---|
1838 | v,f = ASTtransverse(b,fxn) |
---|
1839 | varlist += v |
---|
1840 | fxnlist += f |
---|
1841 | elif node.__class__.__name__ == "Name": |
---|
1842 | varlist += [b] |
---|
1843 | elif fxn and node.__class__.__name__ == "Attribute": |
---|
1844 | varlist += [b] |
---|
1845 | return varlist,fxnlist |
---|
1846 | try: |
---|
1847 | exprast = ast.parse(expr) |
---|
1848 | except SyntaxError: |
---|
1849 | s = '' |
---|
1850 | import traceback |
---|
1851 | for i in traceback.format_exc().splitlines()[-3:-1]: |
---|
1852 | if s: s += "\n" |
---|
1853 | s += str(i) |
---|
1854 | self.lastError = ("Error parsing expression:",s) |
---|
1855 | return |
---|
1856 | # find the variables & functions |
---|
1857 | v,f = ASTtransverse(exprast) |
---|
1858 | varlist = sorted(list(set(v))) |
---|
1859 | fxnlist = list(set(f)) |
---|
1860 | pkgdict = {} |
---|
1861 | # check the functions are defined |
---|
1862 | for fxn in fxnlist: |
---|
1863 | fxndict,fxnobj = FindFunction(fxn) |
---|
1864 | if not fxnobj: |
---|
1865 | self.lastError = ("Error: Invalid function",fxn, |
---|
1866 | "is not defined") |
---|
1867 | return |
---|
1868 | if not hasattr(fxnobj,'__call__'): |
---|
1869 | self.lastError = ("Error: Not a function.",fxn, |
---|
1870 | "cannot be called as a function") |
---|
1871 | return |
---|
1872 | pkgdict.update(fxndict) |
---|
1873 | return varlist,pkgdict |
---|
1874 | |
---|
1875 | def GetDepVar(self): |
---|
1876 | 'return the dependent variable, or None' |
---|
1877 | return self.depVar |
---|
1878 | |
---|
1879 | def SetDepVar(self,var): |
---|
1880 | 'Set the dependent variable, if used' |
---|
1881 | self.depVar = var |
---|
1882 | #========================================================================== |
---|
1883 | class ExpressionCalcObj(object): |
---|
1884 | '''An object used to evaluate an expression from a :class:`ExpressionObj` |
---|
1885 | object. |
---|
1886 | |
---|
1887 | :param ExpressionObj exprObj: a :class:`~ExpressionObj` expression object with |
---|
1888 | an expression string and mappings for the parameter labels in that object. |
---|
1889 | ''' |
---|
1890 | def __init__(self,exprObj): |
---|
1891 | self.eObj = exprObj |
---|
1892 | 'The expression and mappings; a :class:`ExpressionObj` object' |
---|
1893 | self.compiledExpr = None |
---|
1894 | 'The expression as compiled byte-code' |
---|
1895 | self.exprDict = {} |
---|
1896 | '''dict that defines values for labels used in expression and packages |
---|
1897 | referenced by functions |
---|
1898 | ''' |
---|
1899 | self.lblLookup = {} |
---|
1900 | '''Lookup table that specifies the expression label name that is |
---|
1901 | tied to a particular GSAS-II parameters in the parmDict. |
---|
1902 | ''' |
---|
1903 | self.fxnpkgdict = {} |
---|
1904 | '''a dict with references to packages needed to |
---|
1905 | find functions referenced in the expression. |
---|
1906 | ''' |
---|
1907 | self.varLookup = {} |
---|
1908 | '''Lookup table that specifies the GSAS-II variable(s) |
---|
1909 | indexed by the expression label name. (Used for only for diagnostics |
---|
1910 | not evaluation of expression.) |
---|
1911 | ''' |
---|
1912 | self.su = None |
---|
1913 | '''Standard error evaluation where supplied by the evaluator |
---|
1914 | ''' |
---|
1915 | # Patch: for old-style expressions with a (now removed step size) |
---|
1916 | if '2' in platform.python_version_tuple()[0]: |
---|
1917 | basestr = basestring |
---|
1918 | else: |
---|
1919 | basestr = str |
---|
1920 | for v in self.eObj.assgnVars: |
---|
1921 | if not isinstance(self.eObj.assgnVars[v], basestr): |
---|
1922 | self.eObj.assgnVars[v] = self.eObj.assgnVars[v][0] |
---|
1923 | self.parmDict = {} |
---|
1924 | '''A copy of the parameter dictionary, for distance and angle computation |
---|
1925 | ''' |
---|
1926 | |
---|
1927 | def SetupCalc(self,parmDict): |
---|
1928 | '''Do all preparations to use the expression for computation. |
---|
1929 | Adds the free parameter values to the parameter dict (parmDict). |
---|
1930 | ''' |
---|
1931 | if self.eObj.expression.startswith('Dist') or self.eObj.expression.startswith('Angle'): |
---|
1932 | return |
---|
1933 | self.fxnpkgdict = self.eObj.CheckVars() |
---|
1934 | # all is OK, compile the expression |
---|
1935 | self.compiledExpr = compile(self.eObj.expression,'','eval') |
---|
1936 | |
---|
1937 | # look at first value in parmDict to determine its type |
---|
1938 | parmsInList = True |
---|
1939 | if '2' in platform.python_version_tuple()[0]: |
---|
1940 | basestr = basestring |
---|
1941 | else: |
---|
1942 | basestr = str |
---|
1943 | for key in parmDict: |
---|
1944 | val = parmDict[key] |
---|
1945 | if isinstance(val, basestr): |
---|
1946 | parmsInList = False |
---|
1947 | break |
---|
1948 | try: # check if values are in lists |
---|
1949 | val = parmDict[key][0] |
---|
1950 | except (TypeError,IndexError): |
---|
1951 | parmsInList = False |
---|
1952 | break |
---|
1953 | |
---|
1954 | # set up the dicts needed to speed computations |
---|
1955 | self.exprDict = {} |
---|
1956 | self.lblLookup = {} |
---|
1957 | self.varLookup = {} |
---|
1958 | for v in self.eObj.freeVars: |
---|
1959 | varname = self.eObj.freeVars[v][0] |
---|
1960 | varname = "::" + varname.lstrip(':').replace(' ','_').replace(':',';') |
---|
1961 | self.lblLookup[varname] = v |
---|
1962 | self.varLookup[v] = varname |
---|
1963 | if parmsInList: |
---|
1964 | parmDict[varname] = [self.eObj.freeVars[v][1],self.eObj.freeVars[v][2]] |
---|
1965 | else: |
---|
1966 | parmDict[varname] = self.eObj.freeVars[v][1] |
---|
1967 | self.exprDict[v] = self.eObj.freeVars[v][1] |
---|
1968 | for v in self.eObj.assgnVars: |
---|
1969 | varname = self.eObj.assgnVars[v] |
---|
1970 | if varname in parmDict: |
---|
1971 | self.lblLookup[varname] = v |
---|
1972 | self.varLookup[v] = varname |
---|
1973 | if parmsInList: |
---|
1974 | self.exprDict[v] = parmDict[varname][0] |
---|
1975 | else: |
---|
1976 | self.exprDict[v] = parmDict[varname] |
---|
1977 | elif '*' in varname: |
---|
1978 | varlist = LookupWildCard(varname,list(parmDict.keys())) |
---|
1979 | if len(varlist) == 0: |
---|
1980 | self.exprDict[v] = None |
---|
1981 | self.lblLookup[v] = varname # needed? |
---|
1982 | self.exprDict.update(self.fxnpkgdict) # needed? |
---|
1983 | return |
---|
1984 | for var in varlist: |
---|
1985 | self.lblLookup[var] = v |
---|
1986 | if parmsInList: |
---|
1987 | self.exprDict[v] = np.array([parmDict[var][0] for var in varlist]) |
---|
1988 | else: |
---|
1989 | self.exprDict[v] = np.array([parmDict[var] for var in varlist]) |
---|
1990 | self.varLookup[v] = [var for var in varlist] |
---|
1991 | else: |
---|
1992 | self.exprDict[v] = None |
---|
1993 | # raise Exception,"No value for variable "+str(v) |
---|
1994 | self.exprDict.update(self.fxnpkgdict) |
---|
1995 | |
---|
1996 | def UpdateVars(self,varList,valList): |
---|
1997 | '''Update the dict for the expression with a set of values |
---|
1998 | :param list varList: a list of variable names |
---|
1999 | :param list valList: a list of corresponding values |
---|
2000 | ''' |
---|
2001 | for var,val in zip(varList,valList): |
---|
2002 | self.exprDict[self.lblLookup.get(var,'undefined: '+var)] = val |
---|
2003 | |
---|
2004 | def UpdateDict(self,parmDict): |
---|
2005 | '''Update the dict for the expression with values in a dict |
---|
2006 | :param dict parmDict: a dict of values, items not in use are ignored |
---|
2007 | ''' |
---|
2008 | if self.eObj.expression.startswith('Dist') or self.eObj.expression.startswith('Angle'): |
---|
2009 | self.parmDict = parmDict |
---|
2010 | return |
---|
2011 | for var in parmDict: |
---|
2012 | if var in self.lblLookup: |
---|
2013 | self.exprDict[self.lblLookup[var]] = parmDict[var] |
---|
2014 | |
---|
2015 | def EvalExpression(self): |
---|
2016 | '''Evaluate an expression. Note that the expression |
---|
2017 | and mapping are taken from the :class:`ExpressionObj` expression object |
---|
2018 | and the parameter values were specified in :meth:`SetupCalc`. |
---|
2019 | :returns: a single value for the expression. If parameter |
---|
2020 | values are arrays (for example, from wild-carded variable names), |
---|
2021 | the sum of the resulting expression is returned. |
---|
2022 | |
---|
2023 | For example, if the expression is ``'A*B'``, |
---|
2024 | where A is 2.0 and B maps to ``'1::Afrac:*'``, which evaluates to:: |
---|
2025 | |
---|
2026 | [0.5, 1, 0.5] |
---|
2027 | |
---|
2028 | then the result will be ``4.0``. |
---|
2029 | ''' |
---|
2030 | self.su = None |
---|
2031 | if self.eObj.expression.startswith('Dist'): |
---|
2032 | # GSASIIpath.IPyBreak() |
---|
2033 | dist = G2mth.CalcDist(self.eObj.distance_dict, self.eObj.distance_atoms, self.parmDict) |
---|
2034 | return dist |
---|
2035 | elif self.eObj.expression.startswith('Angle'): |
---|
2036 | angle = G2mth.CalcAngle(self.eObj.angle_dict, self.eObj.angle_atoms, self.parmDict) |
---|
2037 | return angle |
---|
2038 | if self.compiledExpr is None: |
---|
2039 | raise Exception("EvalExpression called before SetupCalc") |
---|
2040 | try: |
---|
2041 | val = eval(self.compiledExpr,globals(),self.exprDict) |
---|
2042 | except TypeError: |
---|
2043 | val = None |
---|
2044 | except NameError: |
---|
2045 | val = None |
---|
2046 | if not np.isscalar(val): |
---|
2047 | val = np.sum(val) |
---|
2048 | return val |
---|
2049 | |
---|
2050 | def makeAngleObj(Phase,Oatom,Tatoms): |
---|
2051 | General = Phase['General'] |
---|
2052 | cx,ct = General['AtomPtrs'][:2] |
---|
2053 | pId = Phase['pId'] |
---|
2054 | SGData = General['SGData'] |
---|
2055 | Atoms = Phase['Atoms'] |
---|
2056 | aNames = [atom[ct-1] for atom in Atoms] |
---|
2057 | tIds = [] |
---|
2058 | symNos = [] |
---|
2059 | cellNos = [] |
---|
2060 | oId = aNames.index(Oatom) |
---|
2061 | for Tatom in Tatoms.split(';'): |
---|
2062 | sB = Tatom.find('(')+1 |
---|
2063 | symNo = 0 |
---|
2064 | if sB: |
---|
2065 | sF = Tatom.find(')') |
---|
2066 | symNo = int(Tatom[sB:sF]) |
---|
2067 | symNos.append(symNo) |
---|
2068 | cellNo = [0,0,0] |
---|
2069 | cB = Tatom.find('[') |
---|
2070 | if cB>0: |
---|
2071 | cF = Tatom.find(']')+1 |
---|
2072 | cellNo = eval(Tatom[cB:cF]) |
---|
2073 | cellNos.append(cellNo) |
---|
2074 | tIds.append(aNames.index(Tatom.split('+')[0])) |
---|
2075 | # create an expression object |
---|
2076 | obj = ExpressionObj() |
---|
2077 | obj.expression = 'Angle(%s,%s,\n%s)'%(Tatoms[0],Oatom,Tatoms[1]) |
---|
2078 | obj.angle_dict = {'pId':pId,'SGData':SGData,'symNo':symNos,'cellNo':cellNos} |
---|
2079 | obj.angle_atoms = [oId,tIds] |
---|
2080 | return obj |
---|
2081 | |
---|
2082 | class G2Exception(Exception): |
---|
2083 | 'A generic GSAS-II exception class' |
---|
2084 | def __init__(self,msg): |
---|
2085 | self.msg = msg |
---|
2086 | def __str__(self): |
---|
2087 | return repr(self.msg) |
---|
2088 | |
---|
2089 | class G2RefineCancel(Exception): |
---|
2090 | 'Raised when Cancel is pressed in a refinement dialog' |
---|
2091 | def __init__(self,msg): |
---|
2092 | self.msg = msg |
---|
2093 | def __str__(self): |
---|
2094 | return repr(self.msg) |
---|
2095 | |
---|
2096 | def HowDidIgetHere(wherecalledonly=False): |
---|
2097 | '''Show a traceback with calls that brought us to the current location. |
---|
2098 | Used for debugging. |
---|
2099 | ''' |
---|
2100 | import traceback |
---|
2101 | if wherecalledonly: |
---|
2102 | i = traceback.format_list(traceback.extract_stack()[:-1])[-2] |
---|
2103 | print(i.strip().rstrip()) |
---|
2104 | else: |
---|
2105 | print (70*'*') |
---|
2106 | for i in traceback.format_list(traceback.extract_stack()[:-1]): print(i.strip().rstrip()) |
---|
2107 | print (70*'*') |
---|
2108 | |
---|
2109 | # Note that this is GUI code and should be moved at somepoint |
---|
2110 | def CreatePDFitems(G2frame,PWDRtree,ElList,Qlimits,numAtm=1,FltBkg=0,PDFnames=[]): |
---|
2111 | '''Create and initialize a new set of PDF tree entries |
---|
2112 | |
---|
2113 | :param Frame G2frame: main GSAS-II tree frame object |
---|
2114 | :param str PWDRtree: name of PWDR to be used to create PDF item |
---|
2115 | :param dict ElList: data structure with composition |
---|
2116 | :param list Qlimits: Q limits to be used for computing the PDF |
---|
2117 | :param float numAtm: no. atom in chemical formula |
---|
2118 | :param float FltBkg: flat background value |
---|
2119 | :param list PDFnames: previously used PDF names |
---|
2120 | |
---|
2121 | :returns: the Id of the newly created PDF entry |
---|
2122 | ''' |
---|
2123 | PDFname = 'PDF '+PWDRtree[4:] # this places two spaces after PDF, which is needed is some places |
---|
2124 | if PDFname in PDFnames: |
---|
2125 | print('Skipping, entry already exists: '+PDFname) |
---|
2126 | return None |
---|
2127 | #PDFname = MakeUniqueLabel(PDFname,PDFnames) |
---|
2128 | Id = G2frame.GPXtree.AppendItem(parent=G2frame.root,text=PDFname) |
---|
2129 | Data = { |
---|
2130 | 'Sample':{'Name':PWDRtree,'Mult':1.0}, |
---|
2131 | 'Sample Bkg.':{'Name':'','Mult':-1.0,'Refine':False}, |
---|
2132 | 'Container':{'Name':'','Mult':-1.0,'Refine':False}, |
---|
2133 | 'Container Bkg.':{'Name':'','Mult':-1.0},'ElList':ElList, |
---|
2134 | 'Geometry':'Cylinder','Diam':1.0,'Pack':0.50,'Form Vol':10.0*numAtm,'Flat Bkg':FltBkg, |
---|
2135 | 'DetType':'Area detector','ObliqCoeff':0.3,'Ruland':0.025,'QScaleLim':Qlimits, |
---|
2136 | 'Lorch':False,'BackRatio':0.0,'Rmax':100.,'noRing':False,'IofQmin':1.0,'Rmin':1.0, |
---|
2137 | 'I(Q)':[],'S(Q)':[],'F(Q)':[],'G(R)':[], |
---|
2138 | #items for sequential PDFfit |
---|
2139 | 'Datarange':[0.,30.],'Fitrange':[0.,30.],'qdamp':[0.03,False],'qbroad':[0,False],'Temp':300} |
---|
2140 | G2frame.GPXtree.SetItemPyData(G2frame.GPXtree.AppendItem(Id,text='PDF Controls'),Data) |
---|
2141 | G2frame.GPXtree.SetItemPyData(G2frame.GPXtree.AppendItem(Id,text='PDF Peaks'), |
---|
2142 | {'Limits':[1.,5.],'Background':[2,[0.,-0.2*np.pi],False],'Peaks':[]}) |
---|
2143 | return Id |
---|
2144 | |
---|
2145 | class ShowTiming(object): |
---|
2146 | '''An object to use for timing repeated sections of code. |
---|
2147 | |
---|
2148 | Create the object with:: |
---|
2149 | tim0 = ShowTiming() |
---|
2150 | |
---|
2151 | Tag sections of code to be timed with:: |
---|
2152 | tim0.start('start') |
---|
2153 | tim0.start('in section 1') |
---|
2154 | tim0.start('in section 2') |
---|
2155 | |
---|
2156 | etc. (Note that each section should have a unique label.) |
---|
2157 | |
---|
2158 | After the last section, end timing with:: |
---|
2159 | tim0.end() |
---|
2160 | |
---|
2161 | Show timing results with:: |
---|
2162 | tim0.show() |
---|
2163 | |
---|
2164 | ''' |
---|
2165 | def __init__(self): |
---|
2166 | self.timeSum = [] |
---|
2167 | self.timeStart = [] |
---|
2168 | self.label = [] |
---|
2169 | self.prev = None |
---|
2170 | def start(self,label): |
---|
2171 | import time |
---|
2172 | if label in self.label: |
---|
2173 | i = self.label.index(label) |
---|
2174 | self.timeStart[i] = time.time() |
---|
2175 | else: |
---|
2176 | i = len(self.label) |
---|
2177 | self.timeSum.append(0.0) |
---|
2178 | self.timeStart.append(time.time()) |
---|
2179 | self.label.append(label) |
---|
2180 | if self.prev is not None: |
---|
2181 | self.timeSum[self.prev] += self.timeStart[i] - self.timeStart[self.prev] |
---|
2182 | self.prev = i |
---|
2183 | def end(self): |
---|
2184 | import time |
---|
2185 | if self.prev is not None: |
---|
2186 | self.timeSum[self.prev] += time.time() - self.timeStart[self.prev] |
---|
2187 | self.prev = None |
---|
2188 | def show(self): |
---|
2189 | sumT = sum(self.timeSum) |
---|
2190 | print('Timing results (total={:.2f} sec)'.format(sumT)) |
---|
2191 | for i,(lbl,val) in enumerate(zip(self.label,self.timeSum)): |
---|
2192 | print('{} {:20} {:8.2f} ms {:5.2f}%'.format(i,lbl,1000.*val,100*val/sumT)) |
---|
2193 | |
---|
2194 | def validateAtomDrawType(typ,generalData={}): |
---|
2195 | '''Confirm that the selected Atom drawing type is valid for the current |
---|
2196 | phase. If not, use 'vdW balls'. This is currently used only for setting a |
---|
2197 | default when atoms are added to the atoms draw list. |
---|
2198 | ''' |
---|
2199 | if typ in ('lines','vdW balls','sticks','balls & sticks','ellipsoids'): |
---|
2200 | return typ |
---|
2201 | # elif generalData.get('Type','') == 'macromolecular': |
---|
2202 | # if typ in ('backbone',): |
---|
2203 | # return typ |
---|
2204 | return 'vdW balls' |
---|
2205 | |
---|
2206 | if __name__ == "__main__": |
---|
2207 | # test variable descriptions |
---|
2208 | for var in '0::Afrac:*',':1:Scale','1::dAx:0','::undefined': |
---|
2209 | v = var.split(':')[2] |
---|
2210 | print(var+':\t', getDescr(v),getVarStep(v)) |
---|
2211 | import sys; sys.exit() |
---|
2212 | # test equation evaluation |
---|
2213 | def showEQ(calcobj): |
---|
2214 | print (50*'=') |
---|
2215 | print (calcobj.eObj.expression+'='+calcobj.EvalExpression()) |
---|
2216 | for v in sorted(calcobj.varLookup): |
---|
2217 | print (" "+v+'='+calcobj.exprDict[v]+'='+calcobj.varLookup[v]) |
---|
2218 | # print ' Derivatives' |
---|
2219 | # for v in calcobj.derivStep.keys(): |
---|
2220 | # print ' d(Expr)/d('+v+') =',calcobj.EvalDeriv(v) |
---|
2221 | |
---|
2222 | obj = ExpressionObj() |
---|
2223 | |
---|
2224 | obj.expression = "A*np.exp(B)" |
---|
2225 | obj.assgnVars = {'B': '0::Afrac:1'} |
---|
2226 | obj.freeVars = {'A': [u'A', 0.5, True]} |
---|
2227 | #obj.CheckVars() |
---|
2228 | calcobj = ExpressionCalcObj(obj) |
---|
2229 | |
---|
2230 | obj1 = ExpressionObj() |
---|
2231 | obj1.expression = "A*np.exp(B)" |
---|
2232 | obj1.assgnVars = {'B': '0::Afrac:*'} |
---|
2233 | obj1.freeVars = {'A': [u'Free Prm A', 0.5, True]} |
---|
2234 | #obj.CheckVars() |
---|
2235 | calcobj1 = ExpressionCalcObj(obj1) |
---|
2236 | |
---|
2237 | obj2 = ExpressionObj() |
---|
2238 | obj2.distance_stuff = np.array([[0,1],[1,-1]]) |
---|
2239 | obj2.expression = "Dist(1,2)" |
---|
2240 | GSASIIpath.InvokeDebugOpts() |
---|
2241 | parmDict2 = {'0::Afrac:0':[0.0,True], '0::Afrac:1': [1.0,False]} |
---|
2242 | calcobj2 = ExpressionCalcObj(obj2) |
---|
2243 | calcobj2.SetupCalc(parmDict2) |
---|
2244 | showEQ(calcobj2) |
---|
2245 | |
---|
2246 | parmDict1 = {'0::Afrac:0':1.0, '0::Afrac:1': 1.0} |
---|
2247 | print ('\nDict = '+parmDict1) |
---|
2248 | calcobj.SetupCalc(parmDict1) |
---|
2249 | showEQ(calcobj) |
---|
2250 | calcobj1.SetupCalc(parmDict1) |
---|
2251 | showEQ(calcobj1) |
---|
2252 | |
---|
2253 | parmDict2 = {'0::Afrac:0':[0.0,True], '0::Afrac:1': [1.0,False]} |
---|
2254 | print ('Dict = '+parmDict2) |
---|
2255 | calcobj.SetupCalc(parmDict2) |
---|
2256 | showEQ(calcobj) |
---|
2257 | calcobj1.SetupCalc(parmDict2) |
---|
2258 | showEQ(calcobj1) |
---|
2259 | calcobj2.SetupCalc(parmDict2) |
---|
2260 | showEQ(calcobj2) |
---|