Mantid
Loading...
Searching...
No Matches
CrystalFieldPeaks.cpp
Go to the documentation of this file.
1// Mantid Repository : https://github.com/mantidproject/mantid
2//
3// Copyright © 2018 ISIS Rutherford Appleton Laboratory UKRI,
4// NScD Oak Ridge National Laboratory, European Spallation Source,
5// Institut Laue - Langevin & CSNS, Institute of High Energy Physics, CAS
6// SPDX - License - Identifier: GPL - 3.0 +
10
11#include <functional>
12#include <map>
13
15
16DECLARE_FUNCTION(CrystalFieldPeaks)
17
18
20
21 declareAttribute("Temperature", Attribute(1.0));
22 declareParameter("IntensityScaling", 1.0, "A scaling factor for peak intensities.");
23}
24
25std::string CrystalFieldPeaks::name() const { return "CrystalFieldPeaks"; }
26
28
30
32
34 API::FunctionValues &values) const {
35
38 int nre = 0;
39 calculateEigenSystem(en, wf, nre);
40
41 auto temperature = getAttribute("Temperature").asDouble();
42 IntFortranVector degeneration;
43 DoubleFortranVector eEnergies;
44 DoubleFortranMatrix iEnergies;
45 const double de = getAttribute("ToleranceEnergy").asDouble();
46 const double di = getAttribute("ToleranceIntensity").asDouble();
47 calculateIntensities(nre, en, wf, temperature, de, degeneration, eEnergies, iEnergies);
48
49 DoubleFortranVector eExcitations;
50 DoubleFortranVector iExcitations;
51 calculateExcitations(eEnergies, iEnergies, de, di, eExcitations, iExcitations);
52
53 size_t n = eExcitations.size();
54 if (2 * n > values.size()) {
55 values.expand(2 * n);
56 }
57
59 double scaling = getParameter("IntensityScaling");
60
61 for (size_t i = 0; i < n; ++i) {
62 values.setCalculated(i, eExcitations.get(i));
63 values.setCalculated(i + n, iExcitations.get(i) * scaling);
64 }
65}
66
67} // namespace Mantid::CurveFitting::Functions
#define DECLARE_FUNCTION(classname)
Macro for declaring a new type of function to be used with the FunctionFactory.
Mantid::API::IFunction::Attribute Attribute
Definition: IsoRotDiff.cpp:25
Represent a domain of a very general type.
A class to store values calculated by a function.
size_t size() const
Return the number of values.
void expand(size_t n)
Expand values to a new size, preserve stored values.
void setCalculated(double value)
set all calculated values to same number
IFunctionGeneral: a very general function definition.
double asDouble() const
Returns double value if attribute is a double, throws exception otherwise.
Definition: IFunction.cpp:739
virtual Attribute getAttribute(const std::string &name) const
Return a value of attribute attName.
Definition: IFunction.cpp:1394
virtual double getParameter(size_t i) const =0
Get i-th parameter.
double get(const size_t i) const
Get an element.
size_t size() const
Size of the vector.
CrystalFieldPeaks is a function that calculates crystal field peak positions and intensities.
void calculateEigenSystem(DoubleFortranVector &en, ComplexFortranMatrix &wf, ComplexFortranMatrix &ham, ComplexFortranMatrix &hz, int &nre) const
Calculate the crystal field eigensystem.
CrystalFieldPeaks is a function that calculates crystal field peak positions and intensities.
size_t getNumberValuesPerArgument() const override
Get number of values per argument in the domain.
void functionGeneral(const API::FunctionDomainGeneral &generalDomain, API::FunctionValues &values) const override
Provide a concrete function in an implementation that operates on a FunctionDomainGeneral.
size_t getDefaultDomainSize() const override
Get the default size of a domain.
size_t getNumberDomainColumns() const override
Get number of columns that the domain must have.
size_t m_defaultDomainSize
Store the default domain size after first function evaluation.
std::string name() const override
Returns the function's name.
void MANTID_CURVEFITTING_DLL calculateExcitations(const DoubleFortranVector &e_energies, const DoubleFortranMatrix &i_energies, double de, double di, DoubleFortranVector &e_excitations, DoubleFortranVector &i_excitations)
Calculate the excitations (transition energies) and their intensities.
void MANTID_CURVEFITTING_DLL calculateIntensities(int nre, const DoubleFortranVector &energies, const ComplexFortranMatrix &wavefunctions, double temperature, double de, IntFortranVector &degeneration, DoubleFortranVector &e_energies, DoubleFortranMatrix &i_energies)
Calculate the intensities of transitions.