public release 2.2.0 - see README.md and CHANGES.md for details

This commit is contained in:
2020-09-04 16:22:42 +02:00
parent fbd2d4fa8c
commit 7c61eb1b41
67 changed files with 2934 additions and 682 deletions

View File

@ -42,7 +42,7 @@ class Calculator(object):
or <code>output_file + '.etpai'</code> depending on scan mode.
all other intermediate files are deleted unless keep_temp_files is True.
@param params: a pmsco.project.Params object with all necessary values except cluster and output files set.
@param params: a pmsco.project.CalculatorParams object with all necessary values except cluster and output files set.
@param cluster: a pmsco.cluster.Cluster(format=FMT_EDAC) object with all atom positions set.

View File

@ -49,7 +49,7 @@ class EdacCalculator(calculator.Calculator):
if alpha is defined, theta is implicitly set to normal emission! (to be generalized)
@param params: a pmsco.project.Params object with all necessary values except cluster and output files set.
@param params: a pmsco.project.CalculatorParams object with all necessary values except cluster and output files set.
@param scan: a pmsco.project.Scan() object describing the experimental scanning scheme.
@ -173,7 +173,7 @@ class EdacCalculator(calculator.Calculator):
f.write(" ".join(format(order, "d") for order in params.orders) + "\n")
f.write("emission angle window {0:F}\n".format(params.angular_resolution / 2.0))
# scattering factor output (see project.Params.phase_output_classes)
# scattering factor output (see project.CalculatorParams.phase_output_classes)
if params.phase_output_classes is not None:
fn = "{0}.clu".format(params.output_file)
f.write("cluster output l(A) {fn}\n".format(fn=fn))
@ -197,7 +197,7 @@ class EdacCalculator(calculator.Calculator):
"""
run EDAC with the given parameters and cluster.
@param params: a pmsco.project.Params object with all necessary values except cluster and output files set.
@param params: a pmsco.project.CalculatorParams object with all necessary values except cluster and output files set.
@param cluster: a pmsco.cluster.Cluster(format=FMT_EDAC) object with all atom positions set.

View File

@ -62,7 +62,7 @@ class MscCalculator(calculator.Calculator):
"""
run the MSC program with the given parameters and cluster.
@param params: a project.Params() object with all necessary values except cluster and output files set.
@param params: a project.CalculatorParams() object with all necessary values except cluster and output files set.
@param cluster: a cluster.Cluster(format=FMT_MSC) object with all atom positions set.

View File

@ -28,7 +28,7 @@ PYTHON_EXT_SUFFIX ?= $(shell ${PYTHON_CONFIG} --extension-suffix)
all: phagen
phagen: phagen.exe phagen$(EXT_SUFFIX)
phagen: phagen.exe phagen$(PYTHON_EXT_SUFFIX)
phagen.exe: phagen_scf.f msxas3.inc msxasc3.inc
$(FC) $(FCOPTS) -o phagen.exe phagen_scf.f
@ -36,7 +36,7 @@ phagen.exe: phagen_scf.f msxas3.inc msxasc3.inc
phagen.pyf: | phagen_scf.f
$(F2PY) -h phagen.pyf -m phagen phagen_scf.f only: libmain
phagen$(EXT_SUFFIX): phagen_scf.f phagen.pyf msxas3.inc msxasc3.inc
phagen$(PYTHON_EXT_SUFFIX): phagen_scf.f phagen.pyf msxas3.inc msxasc3.inc
$(F2PY) -c $(F2PYOPTS) -m phagen phagen.pyf phagen_scf.f
clean:

View File

@ -61,7 +61,7 @@ class PhagenCalculator(AtomicCalculator):
because PHAGEN generates a lot of files with hard-coded names,
the function creates a temporary directory for PHAGEN and deletes it before returning.
@param params: pmsco.project.Params object.
@param params: pmsco.project.CalculatorParams object.
the phase_files attribute is updated with the paths of the scattering files.
@param cluster: pmsco.cluster.Cluster object.
@ -76,6 +76,8 @@ class PhagenCalculator(AtomicCalculator):
@return (None, dict) where dict is a list of output files with their category.
the category is "atomic" for all output files.
"""
assert cluster.get_emitter_count() == 1, "PHAGEN cannot handle more than one emitter at a time"
transl = Translator()
transl.params.set_params(params)
transl.params.set_cluster(cluster)
@ -132,6 +134,14 @@ class PhagenCalculator(AtomicCalculator):
except IOError:
logger.error("error loading phagen cluster file {fi}".format(fi=clufile))
try:
listfile = outfile + ".list"
report_listfile = os.path.join(prev_wd, output_file + ".phagen.list")
shutil.copy(listfile, report_listfile)
files[report_listfile] = "log"
except IOError:
logger.error("error loading phagen list file {fi}".format(fi=listfile))
finally:
os.chdir(prev_wd)

View File

@ -19,6 +19,7 @@ from __future__ import print_function
import numpy as np
from pmsco.cluster import Cluster
from pmsco.compat import open
## rydberg energy in electron volts
@ -72,7 +73,7 @@ class TranslationParams(object):
"""
set the translation parameters.
@param params: a pmsco.project.Params object or
@param params: a pmsco.project.CalculatorParams object or
a dictionary containing some or all public fields of this class.
@return: None
"""
@ -125,6 +126,44 @@ class Translator(object):
6. call write_edac_scattering to produce the EDAC scattering matrix files.
7. call write_edac_emission to produce the EDAC emission matrix file.
"""
## @var params
#
# project parameters needed for translation.
#
# fill the attributes of this object before using any translator methods.
## @var scattering
#
# t-matrix storage
#
# the t-matrix is stored in a flat, one-dimensional numpy structured array consisting of the following fields:
# @arg e (float) energy (eV)
# @arg a (int) atom index (1-based)
# @arg l (int) angular momentum quantum number l
# @arg t (complex) scattering matrix element, t = exp(-i * delta) * sin delta
#
# @note PHAGEN uses the convention t = exp(-i * delta) * sin delta,
# whereas EDAC uses t = exp(i * delta) * sin delta (complex conjugate).
# this object stores the t-matrix according to the PHAGEN convention.
# the conversion to the EDAC convention occurs in write_edac_scattering_file().
## @var emission
#
# radial matrix element storage
#
# the radial matrix elemnts are stored in a flat, one-dimensional numpy structured array
# consisting of the following fields:
# @arg e (float) energy (eV)
# @arg dw (complex) matrix element for the transition to l-1
# @arg up (complex) matrix element for the transition to l+1
## @var cluster
#
# cluster object for PHAGEN
#
# this object is created by translate_cluster().
def __init__(self):
"""
initialize the object instance.
@ -134,18 +173,33 @@ class Translator(object):
self.scattering = np.empty(0, dtype=dt)
dt = [('e', 'f4'), ('dw', 'c16'), ('up', 'c16')]
self.emission = np.empty(0, dtype=dt)
self.cluster = None
def translate_cluster(self):
"""
translate the cluster into a form suitable for PHAGEN.
specifically, move the (first and hopefully only) emitter to the first atom position.
the method copies the cluster from self.params into a new object
and stores it under self.cluster.
@return: None
"""
self.cluster = Cluster()
self.cluster.copy_from(self.params.cluster)
ems = self.cluster.get_emitters(['i'])
self.cluster.move_to_first(idx=ems[0][0]-1)
def write_cluster(self, f):
"""
write the cluster section of the PHAGEN input file.
requires a valid pmsco.cluster.Cluster in self.params.cluster.
@param f: file or output stream (an object with a write method)
@return: None
"""
for atom in self.params.cluster.data:
for atom in self.cluster.data:
d = {k: atom[k] for k in atom.dtype.names}
f.write("{s} {t} {x} {y} {z}\n".format(**d))
f.write("-1 -1 0. 0. 0.\n")
@ -163,7 +217,7 @@ class Translator(object):
@return: None
"""
data = self.params.cluster.data
data = self.cluster.data
elements = np.unique(data['t'])
for element in elements:
idx = np.where(data['t'] == element)
@ -181,29 +235,34 @@ class Translator(object):
@return: None
"""
phagen_params = {}
self.translate_cluster()
phagen_params['absorber'] = 1
phagen_params['emin'] = self.params.kinetic_energies.min() / ERYDBERG
phagen_params['emax'] = self.params.kinetic_energies.max() / ERYDBERG
phagen_params['delta'] = (phagen_params['emax'] - phagen_params['emin']) / \
(self.params.kinetic_energies.shape[0] - 1)
if phagen_params['delta'] < 0.0001:
if self.params.kinetic_energies.shape[0] > 1:
phagen_params['delta'] = (phagen_params['emax'] - phagen_params['emin']) / \
(self.params.kinetic_energies.shape[0] - 1)
else:
phagen_params['delta'] = 0.1
phagen_params['edge'] = state_to_edge(self.params.initial_state) # possibly not used
phagen_params['edge'] = state_to_edge(self.params.initial_state)
phagen_params['edge1'] = 'm4' # auger not supported
phagen_params['edge2'] = 'm4' # auger not supported
phagen_params['cip'] = self.params.binding_energy / ERYDBERG
if phagen_params['cip'] < 0.001:
raise ValueError("binding energy parameter is zero.")
if np.sum(np.abs(self.params.cluster.data['q']) >= 0.001) > 0:
if np.sum(np.abs(self.cluster.data['q'])) > 0.:
phagen_params['ionzst'] = 'ionic'
else:
phagen_params['ionzst'] = 'neutral'
if hasattr(f, "write"):
if hasattr(f, "write") and callable(f.write):
f.write("&job\n")
f.write("calctype='xpd',\n")
f.write("coor='angs',\n")
f.write("cip={cip},\n".format(**phagen_params))
f.write("absorber={absorber},\n".format(**phagen_params))
f.write("edge='{edge}',\n".format(**phagen_params))
f.write("edge1='{edge1}',\n".format(**phagen_params))
f.write("edge2='{edge1}',\n".format(**phagen_params))
@ -254,13 +313,18 @@ class Translator(object):
@arg l angular momentum quantum number l
@arg t complex scattering matrix element
@note PHAGEN uses the convention t = exp(-i * delta) * sin delta,
whereas EDAC uses t = exp(i * delta) * sin delta (complex conjugate).
this class stores the t-matrix according to the PHAGEN convention.
the conversion to the EDAC convention occurs in write_edac_scattering_file().
@param f: file or path (any file-like or path-like object that can be passed to numpy.genfromtxt).
@return: None
"""
dt = [('e', 'f4'), ('x1', 'f4'), ('x2', 'f4'), ('na', 'i4'), ('nl', 'i4'),
('tr', 'f8'), ('ti', 'f8'), ('ph', 'f4')]
data = np.genfromtxt(f, dtype=dt)
data = np.atleast_1d(np.genfromtxt(f, dtype=dt))
self.scattering = np.resize(self.scattering, data.shape)
scat = self.scattering
@ -308,7 +372,7 @@ class Translator(object):
@return: None
"""
if hasattr(f, "write"):
if hasattr(f, "write") and callable(f.write):
energies = np.unique(scat['e'])
ne = energies.shape[0]
lmax = scat['l'].max()
@ -323,7 +387,7 @@ class Translator(object):
if ne > 1:
f.write("{0:.3f} ".format(energy))
for item in energy_scat:
f.write(" {0:.6f} {1:.6f}".format(item['t'].real, item['t'].imag))
f.write(" {0:.6f} {1:.6f}".format(item['t'].real, -item['t'].imag))
for i in range(len(energy_scat), lmax + 1):
f.write(" 0 0")
f.write("\n")
@ -341,7 +405,7 @@ class Translator(object):
@return: None
"""
if hasattr(f, "write"):
if hasattr(f, "write") and callable(f.write):
energies = np.unique(scat['e'])
ne = energies.shape[0]
lmax = scat['l'].max()
@ -356,7 +420,8 @@ class Translator(object):
if ne > 1:
f.write("{0:.3f} ".format(energy))
for item in energy_scat:
f.write(" {0:.6f}".format(np.angle(item['t'])))
pha = np.sign(item['t'].real) * np.arcsin(np.sqrt(np.abs(item['t'].imag)))
f.write(" {0:.6f}".format(pha))
for i in range(len(energy_scat), lmax + 1):
f.write(" 0")
f.write("\n")
@ -373,7 +438,7 @@ class Translator(object):
@return: None
"""
dt = [('ar', 'f8'), ('ai', 'f8'), ('br', 'f8'), ('bi', 'f8')]
data = np.genfromtxt(f, dtype=dt)
data = np.atleast_1d(np.genfromtxt(f, dtype=dt))
self.emission = np.resize(self.emission, data.shape)
emission = self.emission
@ -390,7 +455,7 @@ class Translator(object):
@return: None
"""
if hasattr(f, "write"):
if hasattr(f, "write") and callable(f.write):
l0 = self.params.l_init
energies = self.params.kinetic_energies
emission = self.emission