# coding: utf-8 # vim: set et sw=4 ts=4 sts=4 nu ai cc=+0 fdm=indent mouse=a: """ Module calculator ================= This module contains different classes used to define a new calculator for specific spectroscopies understood by MsSpec. These spectroscopies are listed :ref:`here `. There is one *calculator* class for each spectroscopy. The class name is based on the spectroscopy name. For instance, the class for PhotoElectron Diffraction is called :py:class:`_PED`. The helper function :py:func:`calculator.MSSPEC` is used to create objects from these classes by passing the kind of spectroscopy as a keyword argument. For more information on MsSpec, follow this `link `__ """ import os import sys import re import inspect from subprocess import Popen, PIPE from shutil import copyfile, rmtree from datetime import datetime import time from io import StringIO from collections import OrderedDict from ase.calculators.calculator import Calculator import ase.data import ase.atom import ase.atoms import numpy as np from msspec import iodata from msspec.data import electron_be from msspec.config import Config from msspec.misc import (UREG, LOGGER, get_call_info, get_level_from_electron_configuration, XRaySource, set_log_output, log_process_output) from msspec.utils import get_atom_index from msspec.parameters import (PhagenParameters, PhagenMallocParameters, SpecParameters, SpecMallocParameters, GlobalParameters, MuffintinParameters, TMatrixParameters, SourceParameters, DetectorParameters, ScanParameters, CalculationParameters, PEDParameters, EIGParameters) from msspec.calcio import PhagenIO, SpecIO from msspec.phagen.libphagen import main as do_phagen from msspec.spec.libspec import run as do_spec from terminaltables.ascii_table import AsciiTable try: MSSPEC_ROOT = os.environ['MSSPEC_ROOT'] except KeyError: cfg = Config() MSSPEC_ROOT = cfg.get('path') if MSSPEC_ROOT == str(None): raise NameError('No path to the MsSpec distribution found !!') def init_msspec(): LOGGER.debug('Initialization of the msspec module') ase.atom.names['mt_radius'] = ('mt_radii', 0.) ase.atom.names['mt_radius_scale'] = ('mt_radii_scale', 1.) ase.atom.names['proto_index'] = ('proto_indices', 1) ase.atom.names['mean_square_vibration'] = ('mean_square_vibrations', 0.) ase.atom.names['forward_angle'] = ('forward_angles', 20.) ase.atom.names['backward_angle'] = ('backward_angles', 20.) ase.atom.names['RA_cut_off'] = ('RA_cuts_off', 1) ase.atoms.Atoms.absorber = None init_msspec() class _MSCALCULATOR(Calculator): """ This class defines an ASE calculator for doing Multiple scattering calculations. """ implemented_properties = ['', ] __data = {} def __init__(self, spectroscopy='PED', algorithm='expansion', polarization=None, dichroism=None, spinpol=False, folder='./calc', txt='-', **kwargs): stdout = sys.stdout if isinstance(txt, str) and txt != '-': stdout = open(txt, 'w') #elif isinstance(txt, buffer): # stdout = txt elif txt == None: stdout = open('/dev/null', 'a') #set_log_output(stdout) ######################################################################## LOGGER.debug('Initialization of %s', self.__class__.__name__) LOGGER.debug(get_call_info(inspect.currentframe())) ######################################################################## # Init the upper class Calculator.__init__(self, **kwargs) ######################################################################## LOGGER.debug(' create low level parameters') ######################################################################## self.phagen_parameters = PhagenParameters() self.phagen_malloc_parameters = PhagenMallocParameters() self.spec_parameters = SpecParameters() self.spec_malloc_parameters = SpecMallocParameters() ######################################################################## LOGGER.debug(' create higher level parameters') ######################################################################## self.tmatrix_parameters = TMatrixParameters(self.phagen_parameters) self.muffintin_parameters = MuffintinParameters(self.phagen_parameters, self.spec_parameters) self.global_parameters = GlobalParameters(self.phagen_parameters, self.spec_parameters) if spectroscopy == 'PED': self.spectroscopy_parameters = PEDParameters(self.phagen_parameters, self.spec_parameters) elif spectroscopy == 'EIG': self.spectroscopy_parameters = EIGParameters(self.phagen_parameters, self.spec_parameters) #pass else: raise NameError('No such spectrosopy') self.source_parameters = SourceParameters(self.global_parameters, self.phagen_parameters, self.spec_parameters) self.detector_parameters = DetectorParameters(self.global_parameters, self.phagen_parameters, self.spec_parameters) self.scan_parameters = ScanParameters(self.global_parameters, self.phagen_parameters, self.spec_parameters) self.calculation_parameters = CalculationParameters( self.global_parameters, self.phagen_parameters, self.spec_parameters) # updated global parameters with provided keywords self.global_parameters.spectroscopy = spectroscopy self.global_parameters.algorithm = algorithm self.global_parameters.polarization = polarization self.global_parameters.dichroism = dichroism self.global_parameters.spinpol = spinpol self.global_parameters.folder = folder self.phagenio = PhagenIO(self.phagen_parameters, self.phagen_malloc_parameters) self.specio = SpecIO(self.spec_parameters, self.spec_malloc_parameters, self.phagenio) ######################################################################## LOGGER.debug(' create a space dedicated to the calculation') ######################################################################## self.init_folder = os.getcwd() self.msspec_folder = os.path.join(MSSPEC_ROOT) self.tmp_folder = os.path.abspath(folder) LOGGER.debug(' folder: \'%s\'', self.tmp_folder) if not os.path.exists(self.tmp_folder): os.makedirs(self.tmp_folder) os.makedirs(os.path.join(self.tmp_folder, 'input')) os.makedirs(os.path.join(self.tmp_folder, 'output')) #copyfile(os.path.join(self.msspec_folder, 'ase', 'Makefile'), # os.path.join(self.tmp_folder, 'Makefile')) os.chdir(self.tmp_folder) inv = cor = 'NO' if algorithm == 'expansion': pass elif algorithm == 'inversion': inv = 'YES' elif algorithm == 'correlation': cor = 'YES' # spin orbit resolved (not yet) sorb = 'NO' # spin resolved dichro_spinpol = False if dichroism in ('sum_over_spin', 'spin_resolved'): dichro_spinpol = True spin = 'NO' if spinpol or dichro_spinpol: spin = 'YES' if spin == 'YES': LOGGER.error('Option not implemented!') raise NotImplementedError( 'Spin polarization is not implemeted yet!') calctype_spectro = self.spec_parameters.get_parameter('calctype_spectro') calctype_spectro = calctype_spectro.value self._make_opts = (MSSPEC_ROOT, calctype_spectro, inv, cor, spin, sorb, self.tmp_folder) # Initialize the working environment #self._make('init') self.modified = False self.resources = {} ######################################################################## LOGGER.debug(' initialization done.\n') ######################################################################## def _make(self, target): LOGGER.debug(get_call_info(inspect.currentframe())) os.chdir(self.tmp_folder) cmd = ("make__SPACE__ROOT_FOLDER=\"{}\"__SPACE__SPEC=\"{}\"__SPACE__INV=\"{}\"__SPACE__COR=\"{" "}\"__SPACE__" "SPIN=\"{}\"__SPACE__SO=\"{}\"__SPACE__CALC_FOLDER=\"{}\"__SPACE__{}").format(*(self._make_opts + ( target,))).split('__SPACE__') #cmd = cmd.replace(' ', '\ ') #cmd = cmd.split('__SPACE__') LOGGER.debug(' the full command is: %s', cmd) child = Popen(cmd,stdout=PIPE, stderr=PIPE) logger_name = LOGGER.name if target == 'tmatrix': logger_name = 'Phagen' elif target == 'compute': logger_name = 'Spec' log_process_output(child, logger=logger_name) os.chdir(self.init_folder) if child.returncode != 0: LOGGER.error("Unable to successfully run the target: {}".format(target)) sys.exit(1) def _guess_ke(self, level): """ Try to guess the kinetic energy based on the level and the source energy. If the kinetic energy cannot be infered because the level is not reported in the database, the returned value is None. """ try: state = get_level_from_electron_configuration(level) absorber_atomic_number = self.atoms[self.atoms.absorber].number lines = electron_be[absorber_atomic_number] binding_energy = lines[state] except KeyError: # unable to find a binding energy in the database return None # let's assume work function energy (in eV) wf = 4.5 source_energy = self.source_parameters.get_parameter('energy').value ke = source_energy - binding_energy - wf #return np.array(ke, dtype=np.float).flatten() return ke def run_phagen(self): #include_fname = os.path.join(self.tmp_folder, 'src/msxas3.inc') input_fname = os.path.join(self.tmp_folder, 'input/input.ms') #mod0 = self.phagenio.write_include_file(filename=include_fname) mod1 = self.phagenio.write_input_file(filename=input_fname) self.modified = self.modified or mod1 #or mod0 or mod1 if self.modified: # run phagen #self._make('tmatrix') os.chdir(os.path.join(self.tmp_folder, 'output')) do_phagen() # rename some output files to be more explicit os.rename('fort.10', 'cluster.clu') os.rename('fort.35', 'tmatrix.tl') os.rename('fort.55', 'tmatrix.rad') def run_spec(self): def get_li(level): orbitals = 'spdfghi' m = re.match(r'\d(?P[%s])(\d/2)?' % orbitals, level) return orbitals.index(m.group('l')) #include_fname = os.path.join(self.tmp_folder, 'src/spec.inc') input_fname = os.path.join(self.tmp_folder, 'input/spec.dat') kdirs_fname = os.path.join(self.tmp_folder, 'input/kdirs.dat') mod0 = self.specio.write_input_file(filename=input_fname) #mod1 = self.specio.write_include_file(filename=include_fname) mod2 = self.specio.write_kdirs_file(filename=kdirs_fname) #self.modified = self.modified or mod0 or mod1 or mod2 self.modified = self.modified or mod0 or mod2 #self._make('tmatrix') #self._make('bin/spec') #t0 = time.time() #self._make('compute') #t1 = time.time() #self.resources['spec_time'] = t1 - t0 if self.modified: #self.get_tmatrix() t0 = time.time() os.chdir(os.path.join(self.tmp_folder, 'output')) # set/get the dimension values requirements = OrderedDict({ 'NATP_M' : self.phagenio.nat, 'NATCLU_M' : len(self.atoms), 'NAT_EQ_M' : self.phagenio.nateqm, 'N_CL_L_M' : 1, 'NE_M' : self.phagenio.ne, 'NL_M' : self.phagenio.nlmax + 1, 'LI_M' : get_li(self.spec_parameters.extra_level), 'NEMET_M' : 1, 'NO_ST_M' : self.spec_parameters.calc_no, 'NDIF_M' : 10, 'NSO_M' : 2, 'NTEMP_M' : 1, 'NODES_EX_M' : 3, 'NSPIN_M' : 1, # to change for spin dependent 'NTH_M' : 2000, 'NPH_M' : 2000, 'NDIM_M' : 100000, 'N_TILT_M' : 11, # to change see extdir.f 'N_ORD_M' : 200, 'NPATH_M' : 500, 'NGR_M' : 10,}) for key, value in requirements.items(): setattr(self.spec_malloc_parameters, key, value) do_spec(*requirements.values()) t1 = time.time() self.resources['spec_time'] = t1 - t0 def get_tmatrix(self): LOGGER.info("Getting the TMatrix...") LOGGER.debug(get_call_info(inspect.currentframe())) self.run_phagen() filename = os.path.join(self.tmp_folder, 'output/tmatrix.tl') tl = self.phagenio.load_tl_file(filename) filename = os.path.join(self.tmp_folder, 'output/cluster.clu') self.phagenio.load_cluster_file(filename) tl_threshold = self.tmatrix_parameters.get_parameter('tl_threshold') if tl_threshold.value != None: LOGGER.debug(" applying tl_threshold to %s...", tl_threshold.value) go_on = True while go_on: go_on = False for ia in range(self.phagenio.nat): for ie in range(self.phagenio.ne): last_tl = tl[ia][ie][-1, -2:] # convert to complex last_tl = last_tl[0] + last_tl[1] * 1j if np.abs(last_tl) < tl_threshold.value: # remove last line of tl tl[ia][ie] = tl[ia][ie][:-1, :] go_on = True max_tl = self.tmatrix_parameters.get_parameter('max_tl').value cluster = self.phagen_parameters.get_parameter('atoms').value proto_indices = cluster.get_array('proto_indices') if max_tl != None: LOGGER.debug(" applying max_tl: %s", max_tl) for ia in range(self.phagenio.nat): for ie in range(self.phagenio.ne): try: # for each set of tl: # 1. get the symbol of the prototipical atom j = np.where(proto_indices == ia+1) symbol = cluster[j][0].symbol # 2. get the number of max tl allowed ntl = max_tl[symbol] # 3. reshape the tl set accordingly tl[ia][ie] = tl[ia][ie][:ntl, :] except KeyError: pass self.phagenio.write_tl_file( os.path.join(self.tmp_folder, 'output/tmatrix.tl')) # update spec extra parameters here self.spec_parameters.set_parameter('extra_nat', self.phagenio.nat) self.spec_parameters.set_parameter('extra_nlmax', self.phagenio.nlmax) def set_atoms(self, atoms): """Defines the cluster on which the calculator will work. :param atoms: The cluster to attach the calculator to. :type atoms: :py:class:`ase.Atoms` """ if atoms.absorber == None: LOGGER.error("You must define the absorber before setting the atoms to the" "calculator.") self.atoms = atoms self.phagen_parameters.set_parameter('atoms', atoms) self.spec_parameters.set_parameter('extra_atoms', atoms) def get_parameters(self): """Get all the defined parameters in the calculator. :return: A list of all parameters objects. :rtype: List of :py:class:`parameters.Parameter` """ _ = [] for section in ('global', 'muffintin', 'tmatrix', 'spectroscopy', 'source', 'detector', 'scan', 'calculation'): parameters = getattr(self, section + '_parameters') for p in parameters: _.append(p) return _ def shutdown(self): """Removes the temporary folder and all its content. The user may whish to keep the calculation folder (see :ref:`globalparameters-folder`) so it is not removed at the end of the calculation. The calculation folder contains raw results from *Phagen* and *Spec* programs as well as their input files and configuration. It allows the program to save some time by not repeating some tasks (such as the Fortran code generation, building the binaries, computing things that did not changed between runs...). Calling this function at the end of the calculation will erase this calculation folder. .. warning:: Calling this function will erase the produced data without prompting you for confirmation, so take care of explicitly saving your results in your script, by using the :py:func:`iodata.Data.save` method for example. """ LOGGER.info('Deleting temporary files...') rmtree(self.tmp_folder) class _PED(_MSCALCULATOR): """This class creates a calculator object for PhotoElectron DIffraction spectroscopy. :param algorithm: The algorithm to use for the computation. See :ref:`globalparameters-algorithm` for more details about the allowed values and the type. :param polarization: The incoming light polarization (see :ref:`globalparameters-polarization`) :param dichroism: Wether to enable or not the dichroism (see :ref:`globalparameters-dichroism`) :param spinpol: Enable or disable the spin polarization in the calculation (see :ref:`globalparameters-spinpol`) :param folder: The path to the temporary folder for the calculations. See :ref:`globalparameters-folder` :param txt: The name of a file where to redirect standard output. The string '-' will redirect the standard output to the screen (default). :type txt: str .. note:: This class constructor is not meant to be called directly by the user. Use the :py:func:`MSSPEC` to instanciate any calculator. """ def __init__(self, algorithm='expansion', polarization=None, dichroism=None, spinpol=False, folder='./calc', txt='-'): _MSCALCULATOR.__init__(self, spectroscopy='PED', algorithm=algorithm, polarization=polarization, dichroism=dichroism, spinpol=spinpol, folder=folder, txt=txt) self.iodata = iodata.Data('PED Simulation') def _get_scan(self, scan_type='theta', phi=0, theta=np.linspace(-70, 70, 141), level=None, kinetic_energy=None, data=None): LOGGER.info("Computting the %s scan...", scan_type) if data: self.iodata = data if kinetic_energy is None: # try to guess the kinetic energy kinetic_energy = self._guess_ke(level) # if still None... if kinetic_energy is None: LOGGER.error('Unable to guess the kinetic energy!') raise ValueError('You must define a kinetic_energy value.') # update the parameters self.scan_parameters.set_parameter('kinetic_energy', kinetic_energy) all_ke = self.scan_parameters.get_parameter('ke_array') if np.any(all_ke.value < 0): LOGGER.error('Source energy is not high enough or level too deep!') raise ValueError('Kinetic energy is < 0! ({})'.format( kinetic_energy)) self.scan_parameters.set_parameter('type', scan_type) # make sure there is only one energy point in scatf scan if scan_type == 'scatf': assert len(all_ke) == 1, ('kinetic_energy should not be an array ' 'in scatf scan') if scan_type != 'scatf': self.scan_parameters.set_parameter('phi', phi) self.scan_parameters.set_parameter('theta', theta) self.spectroscopy_parameters.set_parameter('level', level) self.get_tmatrix() self.run_spec() # Now load the data ndset = len(self.iodata) dset = self.iodata.add_dset('{} scan [{:d}]'.format(scan_type, ndset)) for p in self.get_parameters(): bundle = {'group': str(p.group), 'name': str(p.name), 'value': str(p.value), 'unit': '' if p.unit is None else str(p.unit)} dset.add_parameter(**bundle) if scan_type in ('theta', 'phi', 'energy'): results_fname = os.path.join(self.tmp_folder, 'output/results.dat') data = self.specio.load_results(results_fname) for _plane, _theta, _phi, _energy, _dirsig, _cs in data.T: if _plane == -1: dset.add_row(theta=_theta, phi=_phi, energy=_energy, cross_section=_cs, direct_signal=_dirsig) elif scan_type in ('scatf',): results_fname = os.path.join(self.tmp_folder, 'output/facdif1.dat') data = self.specio.load_facdif(results_fname) data = data[:, [1, 4, 5, 6, 8]].T _proto, _sf_real, _sf_imag, _theta, _energy = data _sf = _sf_real + _sf_imag * 1j dset.add_columns(proto_index=_proto, sf_real=np.real(_sf), sf_imag=np.imag(_sf), sf_module=np.abs(_sf), theta=_theta, energy=_energy) elif scan_type in ('theta_phi',): results_fname = os.path.join(self.tmp_folder, 'output/results.dat') data = self.specio.load_results(results_fname) #theta_c, phi_c = data[[2, 3], :] #xsec_c = data[-1, :] #dirsig_c = data[-2, :] #dset.add_columns(theta=theta_c) #dset.add_columns(phi=phi_c) #dset.add_columns(cross_section=xsec_c) #dset.add_columns(direct_signal=dirsig_c) for _plane, _theta, _phi, _energy, _dirsig, _cs in data.T: if _plane == -1: dset.add_row(theta=_theta, phi=_phi, energy=_energy, cross_section=_cs, direct_signal=_dirsig) # create a view title = '' for ke in all_ke.value: if scan_type == 'theta': absorber_symbol = self.atoms[self.atoms.absorber].symbol title = 'Polar scan of {}({}) at {:.2f} eV'.format( absorber_symbol, level, ke) xlabel = r'Angle $\theta$($\degree$)' ylabel = r'Signal (a. u.)' view = dset.add_view("E = {:.2f} eV".format(ke), title=title, xlabel=xlabel, ylabel=ylabel) for angle_phi in self.scan_parameters.get_parameter( 'phi').value: where = ("energy=={:.2f} and phi=={:.2f}" "").format(ke, angle_phi) legend = r'$\phi$ = {:.1f} $\degree$'.format(angle_phi) view.select('theta', 'cross_section', where=where, legend=legend) if scan_type == 'phi': absorber_symbol = self.atoms[self.atoms.absorber].symbol title = 'Azimuthal scan of {}({}) at {:.2f} eV'.format( absorber_symbol, level, ke) xlabel = r'Angle $\phi$($\degree$)' ylabel = r'Signal (a. u.)' view = dset.add_view("E = {:.2f} eV".format(ke), title=title, xlabel=xlabel, ylabel=ylabel) for angle_theta in self.scan_parameters.get_parameter( 'theta').value: where = ("energy=={:.2f} and theta=={:.2f}" "").format(ke, angle_theta) legend = r'$\theta$ = {:.1f} $\degree$'.format(angle_theta) view.select('phi', 'cross_section', where=where, legend=legend) if scan_type == 'theta_phi': absorber_symbol = self.atoms[self.atoms.absorber].symbol title = ('Stereographic projection of {}({}) at {:.2f} eV' '').format(absorber_symbol, level, ke) xlabel = r'Angle $\phi$($\degree$)' ylabel = r'Signal (a. u.)' view = dset.add_view("E = {:.2f} eV".format(ke), title=title, xlabel=xlabel, ylabel=ylabel, projection='stereo', colorbar=True, autoscale=True) view.select('theta', 'phi', 'cross_section') if scan_type == 'scatf': for i in range(self.phagenio.nat): proto_index = i+1 title = 'Scattering factor at {:.3f} eV'.format(kinetic_energy) view = dset.add_view("Proto. atom #{:d}".format(proto_index), title=title, projection='polar') where = "proto_index=={:d}".format(proto_index) view.select('theta', 'sf_module', where=where, legend=r'$|f(\theta)|$') view.select('theta', 'sf_real', where=where, legend=r'$\Im(f(\theta))$') view.select('theta', 'sf_imag', where=where, legend=r'$\Re(f(\theta))$') # save the cluster clusbuf = StringIO() self.atoms.info['absorber'] = self.atoms.absorber self.atoms.write(clusbuf, format='xyz') dset.add_parameter(group='Cluster', name='cluster', value=clusbuf.getvalue(), hidden="True") LOGGER.info('%s scan computing done!', scan_type) return self.iodata def get_potential(self, atom_index=None, data=None, units={'energy': 'eV', 'space': 'angstrom'}): """Computes the coulombic part of the atomic potential. :param atom_index: The atom indices to get the potential of, either as a list or as a single integer :param data: The data object to store the results to :param units: The units to be used. A dictionary with the keys 'energy' and 'space' :return: A Data object """ LOGGER.info("Getting the Potential...") LOGGER.debug(get_call_info(inspect.currentframe())) _units = {'energy': 'eV', 'space': 'angstrom'} _units.update(units) if data: self.iodata = data self.run_phagen() filename = os.path.join(self.tmp_folder, 'output/tmatrix.tl') tl = self.phagenio.load_tl_file(filename) filename = os.path.join(self.tmp_folder, 'output/cluster.clu') self.phagenio.load_cluster_file(filename) filename = os.path.join(self.tmp_folder, 'bin/plot/plot_vc.dat') pot_data = self.phagenio.load_potential_file(filename) cluster = self.phagen_parameters.get_parameter('atoms').value dset = self.iodata.add_dset('Potential [{:d}]'.format(len(self.iodata))) r = [] v = [] index = np.empty((0,1), dtype=int) absorber_position = cluster[cluster.absorber].position for _pot_data in pot_data: # find the proto index of these data at_position = (_pot_data['coord'] * UREG.bohr_radius).to('angstrom').magnitude + absorber_position at_index = get_atom_index(cluster, *at_position) at_proto_index = cluster[at_index].get('proto_index') #values = np.asarray(_pot_data['values']) values = _pot_data['values'] index = np.append(index, np.ones(values.shape[0], dtype=int) * at_proto_index) r = np.append(r, (values[:, 0] * UREG.bohr_radius).to(_units['space']).magnitude) v = np.append(v, (values[:, 1] * UREG.rydberg).to(_units['energy']).magnitude) dset.add_columns(distance=r, potential=v, index=index) view = dset.add_view('potential data', title='Potential energy of atoms', xlabel='distance from atomic center [{:s}]'.format(_units['space']), ylabel='energy [{:s}]'.format(_units['energy']), scale='linear', autoscale=True) if atom_index == None: for i in range(pot_data[len(pot_data) - 1]['index']): view.select('distance', 'potential', where="index=={:d}".format(i), legend="Atom index #{:d}".format(i + 1)) else: for i in atom_index: view.select('distance', 'potential', where="index=={:d}".format(cluster[i].get('proto_index') - 1), legend="Atom index #{:d}".format(i)) return self.iodata def get_scattering_factors(self, level='1s', kinetic_energy=None, data=None): """Computes the scattering factors of all prototypical atoms in the cluster. This function computes the real and imaginery parts of the scattering factor as well as its modulus for each non symetrically equivalent atom in the cluster. The results are stored in the *data* object if provided as a parameter. :param level: The electronic level. See :ref:`pedparameters-level`. :param kinetic_energy: see :ref:`scanparameters-kinetic_energy`. :param data: a :py:class:`iodata.Data` object to append the results to or None. :returns: The modified :py:class:`iodata.Data` object passed as an argument or a new :py:class:`iodata.Data` object. """ data = self._get_scan(scan_type='scatf', level=level, data=data, kinetic_energy=kinetic_energy) return data def get_theta_scan(self, phi=0, theta=np.linspace(-70, 70, 141), level=None, kinetic_energy=None, data=None): """Computes a polar scan of the emitted photoelectrons. :param phi: The azimuthal angle in degrees. See :ref:`scanparameters-phi`. :param theta: All the values of the polar angle to be computed. See :ref:`scanparameters-theta`. :param level: The electronic level. See :ref:`pedparameters-level`. :param kinetic_energy: see :ref:`scanparameters-kinetic_energy`. :param data: a :py:class:`iodata.Data` object to append the results to or None. :returns: The modified :py:class:`iodata.Data` object passed as an argument or a new :py:class:`iodata.Data` object. """ data = self._get_scan(scan_type='theta', level=level, theta=theta, phi=phi, kinetic_energy=kinetic_energy, data=data) return data def get_phi_scan(self, phi=np.linspace(0, 359, 359), theta=0, level=None, kinetic_energy=None, data=None): """Computes an azimuthal scan of the emitted photoelectrons. :param phi: All the values of the azimuthal angle to be computed. See :ref:`scanparameters-phi`. :param theta: The polar angle in degrees. See :ref:`scanparameters-theta`. :param level: The electronic level. See :ref:`pedparameters-level`. :param kinetic_energy: see :ref:`scanparameters-kinetic_energy`. :param data: a :py:class:`iodata.Data` object to append the results to or None. :returns: The modified :py:class:`iodata.Data` object passed as an argument or a new :py:class:`iodata.Data` object. """ data = self._get_scan(scan_type='phi', level=level, theta=theta, phi=phi, kinetic_energy=kinetic_energy, data=data) return data def get_theta_phi_scan(self, phi=np.linspace(0, 360), theta=np.linspace(0, 90, 45), level=None, kinetic_energy=None, data=None): """Computes a stereographic scan of the emitted photoelectrons. The azimuth ranges from 0 to 360° and the polar angle ranges from 0 to 90°. :param level: The electronic level. See :ref:`pedparameters-level`. :param kinetic_energy: see :ref:`scanparameters-kinetic_energy`. :param data: a :py:class:`iodata.Data` object to append the results to or None. :returns: The modified :py:class:`iodata.Data` object passed as an argument or a new :py:class:`iodata.Data` object. """ self.spec_malloc_parameters.NPH_M = 8000 data = self._get_scan(scan_type='theta_phi', level=level, theta=theta, phi=phi, kinetic_energy=kinetic_energy, data=data) return data class _EIG(_MSCALCULATOR): """ .. note:: This class constructor is not meant to be called directly by the user. Use the :py:func:`MSSPEC` to instanciate any calculator. """ def __init__(self, algorithm='inversion', polarization=None, dichroism=None, spinpol=False, folder='./calc', txt='-'): _MSCALCULATOR.__init__(self, spectroscopy='EIG', algorithm=algorithm, polarization=polarization, dichroism=dichroism, spinpol=spinpol, folder=folder, txt=txt) if algorithm not in ('inversion', 'power'): LOGGER.error("Only the 'inversion' or the 'power' algorithms " "are supported in EIG spectroscopy mode") exit(1) self.iodata = iodata.Data('EIG Simulation') def get_eigen_values(self, level=None, kinetic_energy=None, data=None): LOGGER.info("Computting the eigen values...") if data: self.iodata = data if kinetic_energy is None: # try to guess the kinetic energy kinetic_energy = self._guess_ke(level) # if still None... if kinetic_energy is None: LOGGER.error('Unable to guess the kinetic energy!') raise ValueError('You must define a kinetic_energy value.') # update the parameters self.scan_parameters.set_parameter('kinetic_energy', kinetic_energy) all_ke = self.scan_parameters.get_parameter('ke_array') if np.any(all_ke < 0): LOGGER.error('Source energy is not high enough or level too deep!') raise ValueError('Kinetic energy is < 0! ({})'.format( kinetic_energy)) self.spectroscopy_parameters.set_parameter('level', level) self.get_tmatrix() self.run_spec() # Now load the data ndset = len(self.iodata) dset = self.iodata.add_dset('Eigen values calculation [{:d}]'.format(ndset)) for p in self.get_parameters(): bundle = {'group': str(p.group), 'name': str(p.name), 'value': str(p.value), 'unit': '' if p.unit is None else str(p.unit)} dset.add_parameter(**bundle) results_fname = os.path.join(self.tmp_folder, 'output/results.dat') data = self.specio.load_results(results_fname) dset.add_columns(energy=data[:,0], eigen_value=data[:,1]) return self.iodata def MSSPEC(spectroscopy='PED', **kwargs): """ The MsSpec calculator constructor. Instanciates a calculator for the given spectroscopy. :param spectroscopy: See :ref:`globalparameters-spectroscopy`. :param kwargs: Other keywords are passed to the spectroscopy-specific constructor. :returns: A :py:class:`calculator._MSCALCULATOR` object. """ module = sys.modules[__name__] cls = getattr(module, '_' + spectroscopy) return cls(**kwargs) if __name__ == "__main__": pass