Mantid
Loading...
Searching...
No Matches
RefinePowderInstrumentParameters.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 +
8
12
17#include "MantidAPI/IFunction.h"
20#include "MantidAPI/TableRow.h"
21#include "MantidAPI/TextAxis.h"
23
29
31
32#include <boost/algorithm/string.hpp>
33#include <boost/algorithm/string/split.hpp>
34
35#include <fstream>
36#include <iomanip>
37#include <utility>
38
39#include <gsl/gsl_sf_erf.h>
40
41using namespace Mantid;
42using namespace Mantid::API;
43using namespace Mantid::Kernel;
44using namespace Mantid::DataObjects;
48using namespace Mantid::HistogramData;
49
50using namespace std;
51
53
55
56//----------------------------------------------------------------------------------------------
60 : m_BestGSLChi2(0.0), m_MinSigma(0.0), m_MinNumFittedPeaks(0), m_MaxNumberStoredParameters(0) {
61 this->useAlgorithm("RefinePowderInstrumentParameters", 3);
62}
63
64//----------------------------------------------------------------------------------------------
68 // Input/output peaks table workspace
69 declareProperty(std::make_unique<API::WorkspaceProperty<DataObjects::TableWorkspace>>("BraggPeakParameterWorkspace",
70 "Anonymous", Direction::Input),
71 "TableWorkspace containg all peaks' parameters.");
72
73 // Input and output instrument parameters table workspace
75 "InstrumentParameterWorkspace", "AnonymousInstrument", Direction::InOut),
76 "TableWorkspace containg instrument's parameters.");
77
78 // Output workspace
79 declareProperty(std::make_unique<API::WorkspaceProperty<DataObjects::Workspace2D>>("OutputWorkspace", "AnonymousOut",
81 "Output Workspace2D for the d-TOF curves. ");
82
83 // Workspace to output fitted peak parameters
84 declareProperty(std::make_unique<WorkspaceProperty<TableWorkspace>>("OutputInstrumentParameterWorkspace",
85 "AnonymousOut2", Direction::Output),
86 "Output TableWorkspace for the fitted peak parameters for each peak.");
87
88 // Workspace to output N best MC parameters
89 declareProperty(std::make_unique<WorkspaceProperty<TableWorkspace>>("OutputBestResultsWorkspace", "",
91 "Output TableWorkspace for the N best MC fitting results. ");
92
93 // Lower limit on number of peaks for fitting
94 declareProperty("MinNumberFittedPeaks", 5, "Minimum number of fitted peaks for refining instrument parameters.");
95
96 // Refinement algorithm
97 vector<string> algoptions{"DirectFit", "MonteCarlo"};
98 auto validator = std::make_shared<Kernel::StringListValidator>(algoptions);
99 declareProperty("RefinementAlgorithm", "MonteCarlo", validator, "Algorithm to refine the instrument parameters.");
100
101 declareProperty("RandomWalkSteps", 10000, "Number of Monte Carlo random walk steps. ");
102
103 // Parameters to fit
104 declareProperty(std::make_unique<Kernel::ArrayProperty<std::string>>("ParametersToFit"),
105 "Names of the parameters to fit. ");
106
107 // Mininum allowed peak's sigma (avoid wrong fitting peak with very narrow
108 // width)
109 declareProperty("MinSigma", 1.0, "Minimum allowed value for Sigma of a peak.");
110
111 // Method to calcualte the standard error of peaks
112 vector<string> stdoptions{"ConstantValue", "InvertedPeakHeight"};
113 auto listvalidator = std::make_shared<Kernel::StringListValidator>(stdoptions);
114 declareProperty("StandardError", "ConstantValue", listvalidator,
115 "Algorithm to calculate the standard error of peak positions.");
116
117 declareProperty("NumberBestFitRecorded", 1, "Number of best fits (Monte Carlo) recorded and output. ");
118
119 declareProperty("MonteCarloRandomSeed", 0, "Random seed for Monte Carlo simulation. ");
120}
121
122//----------------------------------------------------------------------------------------------
126 enum { DirectFit, MonteCarlo } refinealgorithm;
127
128 // 1. Get input
129 DataObjects::TableWorkspace_sptr peakWS = this->getProperty("BraggPeakParameterWorkspace");
130 DataObjects::TableWorkspace_sptr parameterWS = this->getProperty("InstrumentParameterWorkspace");
131
132 m_MinSigma = getProperty("MinSigma");
133
134 int tempint = getProperty("MinNumberFittedPeaks");
135 if (tempint <= 1) {
136 g_log.error() << "Input MinNumberFittedPeaks = " << tempint << " is too small. \n";
137 throw std::invalid_argument("Input MinNumberFittedPeaks is too small.");
138 }
139 m_MinNumFittedPeaks = static_cast<size_t>(tempint);
140
141 tempint = getProperty("NumberBestFitRecorded");
142 if (tempint <= 0)
143 throw runtime_error("Input NumberBestFitRecorded cannot be less and equal to 0. ");
144 m_MaxNumberStoredParameters = static_cast<size_t>(tempint);
145
146 string algoption = getProperty("RefinementAlgorithm");
147 if (algoption == "DirectFit")
148 refinealgorithm = DirectFit;
149 else if (algoption == "MonteCarlo")
150 refinealgorithm = MonteCarlo;
151 else
152 throw runtime_error("RefinementAlgorithm other than DirectFit and "
153 "MonteCarlo are not supported.");
154
155 // 2. Parse input table workspace
156 genPeaksFromTable(peakWS);
159
160 // 3. Generate a cener workspace as function of d-spacing.
161 bool usemc = false;
162 if (refinealgorithm == MonteCarlo)
163 usemc = true;
165
166 // 4. Fit instrument geometry function
167 stringstream errss;
168 TableWorkspace_sptr mcresultws;
169
170 switch (refinealgorithm) {
171 case DirectFit:
172 // a) Simple (directly) fit all parameters
174 break;
175
176 case MonteCarlo:
177 // b) Use Monte Carlo/Annealing method to search global minimum
178 refineInstrumentParametersMC(parameterWS, true);
179 mcresultws = genMCResultTable();
180 setProperty("OutputBestResultsWorkspace", mcresultws);
181 break;
182
183 default:
184 // c) Unsupported
185 errss << "Refinement algorithm " << algoption << " is not supported. Quit!";
186 g_log.error(errss.str());
187 throw invalid_argument(errss.str());
188 break;
189 }
190
191 // 5. Set output workspace
192 this->setProperty("OutputWorkspace", m_dataWS);
193
194 // 6. Output new instrument parameters
196 this->setProperty("OutputInstrumentParameterWorkspace", fitparamws);
197}
198
199//------- Related to Fitting Instrument Geometry Function -------------------
200
204 g_log.debug() << "=========== Method [FitInstrumentParameters] ===============\n";
205
206 // 1. Initialize the fitting function
207 m_Function = std::make_shared<ThermalNeutronDtoTOFFunction>();
208 m_Function->initialize();
209
210 API::FunctionDomain1DVector domain(m_dataWS->x(1).rawData());
211 API::FunctionValues values(domain);
212 const auto &rawY = m_dataWS->y(0);
213 const auto &rawE = m_dataWS->e(0);
214
215 // 2. Set up parameters values
216 std::vector<std::string> funparamnames = m_Function->getParameterNames();
217
218 std::vector<std::string> paramtofit = getProperty("ParametersToFit");
219 std::sort(paramtofit.begin(), paramtofit.end());
220
221 stringstream msgss;
222 msgss << "Set Instrument Function Parameter : \n";
223 std::map<std::string, double>::iterator paramiter;
224 for (const auto &parname : funparamnames) {
225 paramiter = m_FuncParameters.find(parname);
226 if (paramiter == m_FuncParameters.end()) {
227 // Not found and thus skip
228 continue;
229 }
230
231 double parvalue = paramiter->second;
232 m_Function->setParameter(parname, parvalue);
233 msgss << setw(10) << parname << " = " << parvalue << '\n';
234 }
235 g_log.debug() << msgss.str();
236
237 // 2b. Calculate the statistic of the starting values
238 double gslchi2 = calculateFunctionStatistic(m_Function, m_dataWS, 0);
239 double homchi2 = calculateD2TOFFunction(m_Function, domain, values, rawY, rawE);
240 g_log.debug() << "Fit Starting Value: Chi^2 (GSL) = " << gslchi2 << ", Chi2^2 (Home) = " << homchi2 << '\n';
241
242 // 3. Fix parameters that are not listed in parameter-to-fit. Unfix the rest
243 size_t numparams = funparamnames.size();
244 for (size_t i = 0; i < numparams; ++i) {
245 string parname = funparamnames[i];
246 vector<string>::iterator vsiter;
247 vsiter = std::find(paramtofit.begin(), paramtofit.end(), parname);
248
249 if (vsiter == paramtofit.end())
250 m_Function->fix(i);
251 else
252 m_Function->unfix(i);
253 }
254
255 // 4. Select minimizer. Use Simplex for more than 1 parameters to fit.
256 // Levenberg-MarquardtMD otherwise
257 string minimizer("Levenberg-MarquardtMD");
258 if (paramtofit.size() > 1) {
259 minimizer = "Simplex";
260 }
261 g_log.information() << "Fit use minizer: " << minimizer << '\n';
262
263 // 5. Create and setup fit algorithm
264 g_log.information() << "Fit instrument geometry: " << m_Function->asString() << '\n';
265
266 stringstream outss;
267 for (size_t i = 0; i < m_dataWS->x(0).size(); ++i)
268 outss << m_dataWS->x(0)[i] << "\t\t" << m_dataWS->y(0)[i] << "\t\t" << m_dataWS->e(0)[i] << '\n';
269 g_log.debug() << "Input Peak Position Workspace To Fit: \n" << outss.str() << '\n';
270
271 auto fitalg = createChildAlgorithm("Fit", 0.0, 0.2, true);
272 fitalg->initialize();
273
274 fitalg->setProperty("Function", std::dynamic_pointer_cast<API::IFunction>(m_Function));
275 fitalg->setProperty("InputWorkspace", m_dataWS);
276 fitalg->setProperty("WorkspaceIndex", 0);
277 fitalg->setProperty("Minimizer", minimizer);
278 fitalg->setProperty("CostFunction", "Least squares");
279 fitalg->setProperty("MaxIterations", 1000);
280
281 bool successfulfit = fitalg->execute();
282 if (!fitalg->isExecuted() || !successfulfit) {
283 // Early return due to bad fit
284 g_log.error() << "Fitting to instrument geometry function failed. \n";
285 throw std::runtime_error("Fitting failed.");
286 }
287
288 double chi2 = fitalg->getProperty("OutputChi2overDoF");
289 std::string fitstatus = fitalg->getProperty("OutputStatus");
290
291 g_log.debug() << "Fit Result (GSL): Chi^2 = " << chi2 << "; Fit Status = " << fitstatus << '\n';
292
293 API::IFunction_sptr fitfunc = fitalg->getProperty("Function");
294
295 // 4. Set the output data (model and diff)
296 m_Function->function(domain, values);
297
298 for (size_t i = 0; i < domain.size(); ++i) {
299 m_dataWS->mutableY(1)[i] = values[i];
300 m_dataWS->mutableY(2)[i] = m_dataWS->y(0)[i] - values[i];
301 }
302
303 double selfchi2 = calculateD2TOFFunction(m_Function, domain, values, rawY, rawE);
304 g_log.debug() << "Homemade Chi^2 = " << selfchi2 << '\n';
305
306 // 5. Update fitted parameters
307 for (const auto &parname : funparamnames) {
308 double parvalue = fitfunc->getParameter(parname);
309 m_FuncParameters[parname] = parvalue;
310 }
311
312 // 6. Pretty screen output
313 stringstream dbss;
314 dbss << "************ Fit Parameter Result *************\n";
315 for (paramiter = m_FuncParameters.begin(); paramiter != m_FuncParameters.end(); ++paramiter) {
316 std::string parname = paramiter->first;
317 double inpparvalue = m_OrigParameters[parname];
318 double parvalue = paramiter->second;
319 dbss << setw(20) << parname << " = " << setw(15) << setprecision(6) << parvalue << "\t\tFrom " << setw(15)
320 << setprecision(6) << inpparvalue << "\t\tDiff = " << inpparvalue - parvalue << '\n';
321 }
322 dbss << "*********************************************\n";
323 g_log.debug() << dbss.str();
324
325 // 7. Play with Zscore: template<typename TYPE>
326 // std::vector<double> getZscore(const std::vector<TYPE>& data, const bool
327 // sorted=false);
328 vector<double> z0 = Kernel::getZscore(m_dataWS->y(0).rawData());
329 vector<double> z1 = Kernel::getZscore(m_dataWS->y(1).rawData());
330 vector<double> z2 = Kernel::getZscore(m_dataWS->y(2).rawData());
331 stringstream zss;
332 zss << setw(20) << "d_h" << setw(20) << "Z DataY" << setw(20) << "Z ModelY" << setw(20) << "Z DiffY" << setw(20)
333 << "DiffY\n";
334 const auto &X = m_dataWS->x(0);
335 const auto &Y = m_dataWS->y(2);
336 for (size_t i = 0; i < z0.size(); ++i) {
337 double d_h = X[i];
338 double zdatay = z0[i];
339 double zmodely = z1[i];
340 double zdiffy = z2[i];
341 double diffy = Y[i];
342 zss << setw(20) << d_h << setw(20) << zdatay << setw(20) << zmodely << setw(20) << zdiffy << setw(20) << diffy
343 << '\n';
344 }
345 g_log.debug() << "Zscore Survey: \b" << zss.str();
346}
347
351 auto fitalg = createChildAlgorithm("Fit", 0.0, 0.2, true);
352 fitalg->initialize();
353
354 fitalg->setProperty("Function", std::dynamic_pointer_cast<API::IFunction>(func));
355 fitalg->setProperty("InputWorkspace", m_dataWS);
356 fitalg->setProperty("WorkspaceIndex", 0);
357 fitalg->setProperty("Minimizer", "Simplex");
358 fitalg->setProperty("CostFunction", "Least squares");
359 fitalg->setProperty("MaxIterations", 1000);
360
361 bool successfulfit = fitalg->execute();
362 if (!fitalg->isExecuted() || !successfulfit) {
363 // Early return due to bad fit
364 g_log.error() << "Fitting to instrument geometry function failed. \n";
365 throw std::runtime_error("Fitting failed.");
366 }
367
368 gslchi2 = fitalg->getProperty("OutputChi2overDoF");
369 std::string fitstatus = fitalg->getProperty("OutputStatus");
370
371 g_log.debug() << "Function Fit: Chi^2 = " << gslchi2 << "; Fit Status = " << fitstatus << '\n';
372
373 bool fitgood = (fitstatus == "success");
374
375 return fitgood;
376}
377
381 const MatrixWorkspace_sptr &dataws,
382 size_t workspaceindex) {
383 // 1. Fix all parameters of the function
384 vector<string> funcparameters = func->getParameterNames();
385 size_t numparams = funcparameters.size();
386 for (size_t i = 0; i < numparams; ++i) {
387 func->fix(i);
388 }
389
390 // 2. Call a non fit refine
391 auto fitalg = createChildAlgorithm("Fit", 0.0, 0.2, true);
392 fitalg->initialize();
393
394 fitalg->setProperty("Function", std::dynamic_pointer_cast<API::IFunction>(func));
395 fitalg->setProperty("InputWorkspace", dataws);
396 fitalg->setProperty("WorkspaceIndex", static_cast<int>(workspaceindex));
397 fitalg->setProperty("Minimizer", "Levenberg-MarquardtMD");
398 fitalg->setProperty("CostFunction", "Least squares");
399 fitalg->setProperty("MaxIterations", 2);
400
401 bool successfulfit = fitalg->execute();
402 if (!fitalg->isExecuted() || !successfulfit) {
403 // Early return due to bad fit
404 g_log.error() << "Fitting to instrument geometry function failed. \n";
405 throw std::runtime_error("Fitting failed.");
406 }
407
408 double chi2 = fitalg->getProperty("OutputChi2overDoF");
409 std::string fitstatus = fitalg->getProperty("OutputStatus");
410
411 g_log.debug() << "Function calculation [L.M]: Chi^2 = " << chi2 << "; Fit Status = " << fitstatus << '\n';
412
413 return chi2;
414}
415
419 // 1. Get function's parameter names
421
422 // 2. Parse parameter (table) workspace
423 vector<double> stepsizes, lowerbounds, upperbounds;
424 importMonteCarloParametersFromTable(parameterWS, m_PeakFunctionParameterNames, stepsizes, lowerbounds, upperbounds);
425
426 stringstream dbss;
427 for (size_t i = 0; i < m_PeakFunctionParameterNames.size(); ++i) {
428 dbss << setw(20) << m_PeakFunctionParameterNames[i] << ": Min = " << setw(15) << setprecision(6) << lowerbounds[i]
429 << ", Max = " << setw(15) << setprecision(6) << upperbounds[i] << ", Step Size = " << setw(15)
430 << setprecision(6) << stepsizes[i] << '\n';
431 }
432 g_log.notice() << "Monte Carlo Parameters: \n" << dbss.str();
433
434 // 3. Maximum step size
435 size_t maxsteps;
436 int tempint = getProperty("RandomWalkSteps");
437 if (tempint > 0)
438 maxsteps = static_cast<size_t>(tempint);
439 else
440 throw runtime_error("RandomwWalkSteps cannot be less than or equal to 0. ");
441
442 // 4. Random seed and step size rescale factor
443 int randomseed = getProperty("MonteCarloRandomSeed");
444 srand(randomseed);
445
446 double stepsizescalefactor = 1.1;
447
448 // 5. Monte Carlo simulation
449 doParameterSpaceRandomWalk(m_PeakFunctionParameterNames, lowerbounds, upperbounds, stepsizes, maxsteps,
450 stepsizescalefactor, fit2);
451
452 // 6. Record the result
453 const auto &X = m_dataWS->x(0);
454 const auto &Y = m_dataWS->y(0);
455 const auto &E = m_dataWS->e(0);
456 FunctionDomain1DVector domain(X.rawData());
457 FunctionValues values(domain);
458 for (size_t i = 0; i < m_BestFitParameters.size(); ++i) {
459 // a. Set the function with the
460 for (size_t j = 0; j < m_PeakFunctionParameterNames.size(); ++j) {
461 m_Function->setParameter(m_PeakFunctionParameterNames[j], m_BestFitParameters[i].second[j]);
462 }
463
464 // b. Calculate
465 calculateD2TOFFunction(m_Function, domain, values, Y, E);
466
467 vector<double> vec_n;
469
470 // c. Put the data to output workspace
471 auto &newY = m_dataWS->mutableY(3 * i + 1);
472 auto &newD = m_dataWS->mutableY(3 * i + 2);
473 auto &newN = m_dataWS->mutableY(3 * i + 3);
474 for (size_t j = 0; j < newY.size(); ++j) {
475 newY[j] = values[j];
476 newD[j] = Y[j] - values[j];
477 newN[j] = vec_n[j];
478 }
479 }
480}
481
487 const vector<double> &lowerbounds,
488 const vector<double> &upperbounds,
489 vector<double> &stepsizes, size_t maxsteps,
490 double stepsizescalefactor, bool fit2) {
491
492 // 1. Set up starting values, esp. to m_Function
493 size_t numparameters = parnames.size();
494 vector<double> paramvalues;
495 for (size_t i = 0; i < numparameters; ++i) {
496 string parname = parnames[i];
497 double parvalue = m_FuncParameters[parname];
498 paramvalues.emplace_back(parvalue);
499 m_Function->setParameter(parname, parvalue);
500 }
501
502 // Calcualte the function's initial statistic
504 g_log.debug() << "Function with starting values has Chi2 = " << m_BestGSLChi2 << " (GSL L.M) \n";
505
506 const auto &X = m_dataWS->x(0);
507 const auto &rawY = m_dataWS->y(0);
508 const auto &rawE = m_dataWS->e(0);
509 FunctionDomain1DVector domain(X.rawData());
510 FunctionValues values(domain);
511
512 // 2. Determine the parameters to fit
513 vector<string> paramstofit = getProperty("ParametersToFit");
514 set<string> paramstofitset;
515 bool refineallparams;
516 if (paramstofit.empty()) {
517 // Default case to refine all parameters
518 refineallparams = true;
519 } else {
520 // Refine part of the parameters
521 refineallparams = false;
522 vector<string>::iterator vsiter;
523 for (vsiter = paramstofit.begin(); vsiter != paramstofit.end(); ++vsiter) {
524 paramstofitset.insert(*vsiter);
525 }
526 }
527
528 stringstream dbss;
529 set<string>::iterator setiter;
530 for (setiter = paramstofitset.begin(); setiter != paramstofitset.end(); ++setiter) {
531 string paramName = *setiter;
532 dbss << setw(20) << paramName;
533 }
534 g_log.notice() << "Parameters to refine: " << dbss.str() << '\n';
535
536 // 3. Create a local function for fit and set the parameters unfixed
538 std::shared_ptr<ThermalNeutronDtoTOFFunction>(new ThermalNeutronDtoTOFFunction());
539 func4fit->initialize();
540 for (size_t i = 0; i < numparameters; ++i) {
541 string parname = parnames[i];
542 // Fit or fix
543 if (paramstofitset.count(parname) > 0)
544 func4fit->unfix(i);
545 else
546 func4fit->fix(i);
547 // Constraint
548 double lowerb = lowerbounds[i];
549 double upperb = upperbounds[i];
550 auto newconstraint = std::make_unique<BoundaryConstraint>(func4fit.get(), parname, lowerb, upperb);
551 func4fit->addConstraint(std::move(newconstraint));
552 }
553 g_log.debug() << "Function for fitting in MC: " << func4fit->asString() << '\n';
554
555 // 4. Do MC loops
556 double curchi2 = calculateD2TOFFunction(m_Function, domain, values, rawY, rawE);
557
558 g_log.notice() << "Monte Carlo Random Walk Starting Chi^2 = " << curchi2 << '\n';
559
560 size_t paramindex = 0;
561 size_t numacceptance = 0;
562 for (size_t istep = 0; istep < maxsteps; ++istep) {
563 // a. Determine whether to refine this parameter
564 if (!refineallparams) {
565 if (paramstofitset.count(parnames[paramindex]) == 0) {
566 ++paramindex;
567 if (paramindex >= parnames.size())
568 paramindex = 0;
569 continue;
570 }
571 }
572
573 // b. Propose for a new value
574 double randomnumber = static_cast<double>(rand()) / static_cast<double>(RAND_MAX);
575 double newvalue = paramvalues[paramindex] + (randomnumber - 1.0) * stepsizes[paramindex];
576 if (newvalue > upperbounds[paramindex]) {
577 newvalue = lowerbounds[paramindex] + (newvalue - upperbounds[paramindex]);
578 } else if (newvalue < lowerbounds[paramindex]) {
579 newvalue = upperbounds[paramindex] - (lowerbounds[paramindex] - newvalue);
580 }
581
582 try {
583 m_Function->setParameter(parnames[paramindex], newvalue);
584 } catch (runtime_error &) {
585 stringstream errss;
586 errss << "New Value = " << newvalue << ", Random Number = " << randomnumber
587 << "Step Size = " << stepsizes[paramindex] << ", Step size rescale factor = " << stepsizescalefactor;
588 g_log.error(errss.str());
589 throw;
590 }
591
592 // b. Calcualte the new
593 double newchi2 = calculateD2TOFFunction(m_Function, domain, values, rawY, rawE);
594
595 // Optionally fit
596 if (fit2) {
597 // i. Copy the parameters
598 for (size_t i = 0; i < numparameters; ++i) {
599 double parvalue = m_Function->getParameter(i);
600 func4fit->setParameter(i, parvalue);
601 }
602
603 // ii. Fit function
604 double gslchi2;
605 bool fitgood = fitFunction(func4fit, gslchi2);
606
607 if (fitgood) {
608 // iii. Caculate
609 double homchi2 = calculateD2TOFFunction(func4fit, domain, values, rawY, rawE);
610
611 if (gslchi2 < m_BestGSLChi2)
612 m_BestGSLChi2 = gslchi2;
613
614 // iv. Archive
615 vector<double> newparvalues;
616 newparvalues.reserve(numparameters);
617 for (size_t i = 0; i < numparameters; ++i) {
618 double parvalue = func4fit->getParameter(i);
619 newparvalues.emplace_back(parvalue);
620 }
621 m_BestFitParameters.emplace_back(homchi2, newparvalues);
622 m_BestFitChi2s.emplace_back(homchi2, gslchi2);
623
624 // v. Sort and keep in size
625 sort(m_BestFitParameters.begin(), m_BestFitParameters.end());
626 sort(m_BestFitChi2s.begin(), m_BestFitChi2s.end());
628 m_BestFitParameters.pop_back();
629 m_BestFitChi2s.pop_back();
630 }
631 }
632 }
633
634 // c. Accept?
635 bool accept;
636 double prob = exp(-(newchi2 - curchi2) / curchi2);
637 double randnumber = static_cast<double>(rand()) / static_cast<double>(RAND_MAX);
638 accept = randnumber < prob;
639
640 // d. Adjust step size
641 if (false) {
642 if (newchi2 < curchi2) {
643 // Decrease step size if it approaches a minima
644 stepsizes[paramindex] = stepsizes[paramindex] / stepsizescalefactor;
645 } else {
646 // Increase step size if it diverges
647 double newstepsize = stepsizes[paramindex] * stepsizescalefactor;
648 double maxstepsize = upperbounds[paramindex] - lowerbounds[paramindex];
649 if (newstepsize >= maxstepsize) {
650 newstepsize = maxstepsize;
651 }
652 stepsizes[paramindex] = newstepsize;
653 }
654 }
655
656 // e. Record the solution
657 if (accept) {
658 // i. Accept the new value
659 paramvalues[paramindex] = newvalue;
660
661 // ii. Add the new values to vector
662 m_BestMCParameters.emplace_back(newchi2, paramvalues);
663
664 // iii. Sort and delete the last if necessary
665 sort(m_BestMCParameters.begin(), m_BestMCParameters.end());
667 m_BestMCParameters.pop_back();
668
669 // iv. Update chi2 and ...
670 curchi2 = newchi2;
671 ++numacceptance;
672 }
673
674 // z. Update the parameter index for next movement
675 ++paramindex;
676 if (paramindex >= numparameters)
677 paramindex = 0;
678
679 } // FOREACH MC Step
680
681 // 3. Debug output
682 stringstream mcresult;
683 mcresult << "Monte Carlo Result for " << m_BestMCParameters.size() << " Best Results\n";
684 mcresult << "Number of acceptance = " << numacceptance << ", out of " << maxsteps << " MC steps."
685 << "Accept ratio = " << static_cast<double>(numacceptance) / static_cast<double>(maxsteps) << '\n';
686 mcresult << "Best " << m_BestMCParameters.size() << " Monte Carlo (no fit) results: \n";
687 for (size_t i = 0; i < m_BestMCParameters.size(); ++i) {
688 mcresult << setw(3) << i << ": Chi^2 = " << m_BestMCParameters[i].first << '\n';
689 }
690 mcresult << "Best " << m_BestMCParameters.size() << " fitting results. Best Chi^2 = " << m_BestGSLChi2 << '\n';
691 for (size_t i = 0; i < m_BestFitParameters.size(); ++i) {
692 mcresult << setw(3) << i << ": Chi^2 = " << m_BestFitParameters[i].first
693 << ", GSL Chi^2 = " << m_BestFitChi2s[i].second << '\n';
694 }
695 g_log.notice() << mcresult.str();
696}
697
701 // 1. Clear output
702 parnames.clear();
703
704 // 2. Get the parameter names from function
705 m_Function = std::make_shared<ThermalNeutronDtoTOFFunction>();
706 std::vector<std::string> funparamnames = m_Function->getParameterNames();
707
708 // 3. Copy
709 parnames = funparamnames;
710}
711
715 const API::FunctionDomain1DVector &domain,
716 API::FunctionValues &values,
717 const Mantid::HistogramData::HistogramY &rawY,
718 const Mantid::HistogramData::HistogramE &rawE) {
719 // 1. Check validity
720 if (!func) {
721 throw std::runtime_error("m_Function has not been initialized!");
722 }
723
724 if (domain.size() != values.size() || domain.size() != rawY.size() || rawY.size() != rawE.size()) {
725 throw std::runtime_error("Input domain, values and raw data have different sizes.");
726 }
727
728 // 2. Calculate vlaues
729 func->function(domain, values);
730
731 // 3. Calculate the difference
732 double chi2 = 0;
733 for (size_t i = 0; i < domain.size(); ++i) {
734 double temp = (values[i] - rawY[i]) / rawE[i];
735 chi2 += temp * temp;
736 // cout << "Peak " << i << ": Model = " << values[i] << ", Data = " <<
737 // rawY[i] << ". Standard Error = " << rawE[i] << '\n';
738 }
739
740 return chi2;
741}
742
743//------------------------------- Processing Inputs
744//----------------------------------------
749 // 1. Check and clear input and output
750 if (!peakparamws) {
751 g_log.error() << "Input tableworkspace for peak parameters is invalid!\n";
752 throw std::invalid_argument("Invalid input table workspace for peak parameters");
753 }
754
755 m_Peaks.clear();
756
757 // 2. Parse table workspace rows to generate peaks
758 vector<string> colnames = peakparamws->getColumnNames();
759 size_t numrows = peakparamws->rowCount();
760
761 for (size_t ir = 0; ir < numrows; ++ir) {
762 // a) Generate peak
763 BackToBackExponential_sptr newpeakptr = std::make_shared<BackToBackExponential>();
764 newpeakptr->initialize();
765
766 // b) Parse parameters
767 int h, k, l;
768 double alpha, beta, tof_h, sigma, sigma2, chi2, height, dbtemp;
769 string strtemp;
770 sigma2 = -1;
771
772 API::TableRow row = peakparamws->getRow(ir);
773 for (const auto &colname : colnames) {
774 if (colname == "H")
775 row >> h;
776 else if (colname == "K")
777 row >> k;
778 else if (colname == "L")
779 row >> l;
780 else if (colname == "Alpha")
781 row >> alpha;
782 else if (colname == "Beta")
783 row >> beta;
784 else if (colname == "Sigma2")
785 row >> sigma2;
786 else if (colname == "Sigma")
787 row >> sigma;
788 else if (colname == "Chi2")
789 row >> chi2;
790 else if (colname == "Height")
791 row >> height;
792 else if (colname == "TOF_h")
793 row >> tof_h;
794 else {
795 try {
796 row >> dbtemp;
797 } catch (runtime_error &) {
798 row >> strtemp;
799 }
800 }
801 }
802
803 if (sigma2 > 0)
804 sigma = sqrt(sigma2);
805
806 // c) Set peak parameters and etc.
807 newpeakptr->setParameter("A", alpha);
808 newpeakptr->setParameter("B", beta);
809 newpeakptr->setParameter("S", sigma);
810 newpeakptr->setParameter("X0", tof_h);
811 newpeakptr->setParameter("I", height);
812
813 std::vector<int> hkl;
814 hkl.emplace_back(h);
815 hkl.emplace_back(k);
816 hkl.emplace_back(l);
817
818 m_Peaks.emplace(hkl, newpeakptr);
819
820 m_PeakErrors.emplace(hkl, chi2);
821
822 g_log.information() << "[Generatem_Peaks] Peak " << ir << " HKL = [" << hkl[0] << ", " << hkl[1] << ", " << hkl[2]
823 << "], Input Center = " << setw(10) << setprecision(6) << newpeakptr->centre() << '\n';
824
825 } // ENDFOR Each potential peak
826}
827
832 std::map<std::string, double> &parameters) {
833 // 1. Check column orders
834 std::vector<std::string> colnames = parameterWS->getColumnNames();
835 if (colnames.size() < 2) {
836 g_log.error() << "Input parameter table workspace does not have enough "
837 "number of columns. "
838 << " Number of columns = " << colnames.size() << " < 3 as required. \n";
839 throw std::runtime_error("Input parameter workspace is wrong. ");
840 }
841
842 if (colnames[0] != "Name" || colnames[1] != "Value") {
843 g_log.error() << "Input parameter table workspace does not have the "
844 "columns in order. "
845 << " It must be Name, Value, FitOrTie.\n";
846 throw std::runtime_error("Input parameter workspace is wrong. ");
847 }
848
849 // 2. Import data to maps
850 std::string parname;
851 double value;
852
853 size_t numrows = parameterWS->rowCount();
854
855 for (size_t ir = 0; ir < numrows; ++ir) {
856 try {
857 API::TableRow trow = parameterWS->getRow(ir);
858 trow >> parname >> value;
859 parameters.emplace(parname, value);
860 } catch (runtime_error &) {
861 g_log.error() << "Import table workspace " << parameterWS->getName() << " error in line " << ir << ". "
862 << " Requires [string, double] in the first 2 columns.\n";
863 throw;
864 }
865 }
866}
867
872 const vector<string> &parameternames,
873 vector<double> &stepsizes,
874 vector<double> &lowerbounds,
875 vector<double> &upperbounds) {
876 // 1. Get column information
877 vector<string> colnames = tablews->getColumnNames();
878 size_t imax, imin, istep;
879 imax = colnames.size() + 1;
880 imin = imax;
881 istep = imax;
882
883 for (size_t i = 0; i < colnames.size(); ++i) {
884 if (colnames[i] == "Max")
885 imax = i;
886 else if (colnames[i] == "Min")
887 imin = i;
888 else if (colnames[i] == "StepSize")
889 istep = i;
890 }
891
892 if (imax > colnames.size() || imin > colnames.size() || istep > colnames.size()) {
893 stringstream errss;
894 errss << "Input parameter workspace misses information for Monte Carlo "
895 "minimizer. "
896 << "One or more of the following columns are missing (Max, Min, "
897 "StepSize).";
898 g_log.error(errss.str());
899 throw runtime_error(errss.str());
900 }
901
902 // 2. Parse input to a map
903 stepsizes.clear();
904 lowerbounds.clear();
905 upperbounds.clear();
906
907 map<string, vector<double>> mcparameters;
908 size_t numrows = tablews->rowCount();
909 for (size_t ir = 0; ir < numrows; ++ir) {
910 TableRow row = tablews->getRow(ir);
911 string parname;
912 double tmax = 0, tmin = 0, tstepsize = 0;
913 row >> parname;
914 for (size_t ic = 1; ic < colnames.size(); ++ic) {
915 double tmpdbl = std::numeric_limits<float>::quiet_NaN();
916 try {
917 row >> tmpdbl;
918 } catch (runtime_error &) {
919 g_log.error() << "Import MC parameter " << colnames[ic] << " error in row " << ir << " of workspace "
920 << tablews->getName() << '\n';
921 string tmpstr;
922 row >> tmpstr;
923 g_log.error() << "Should be " << tmpstr << '\n';
924 }
925
926 if (ic == imax)
927 tmax = tmpdbl;
928 else if (ic == imin)
929 tmin = tmpdbl;
930 else if (ic == istep)
931 tstepsize = tmpdbl;
932 }
933 vector<double> tmpvec;
934 tmpvec.emplace_back(tmin);
935 tmpvec.emplace_back(tmax);
936 tmpvec.emplace_back(tstepsize);
937 mcparameters.emplace(parname, tmpvec);
938 }
939
940 // 3. Retrieve the information for geometry parameters
941 for (const auto &parname : parameternames) {
942 // a) Get on hold of the MC parameter vector
943 auto mit = mcparameters.find(parname);
944 if (mit == mcparameters.end()) {
945 // Not found the parameter. raise error!
946 stringstream errss;
947 errss << "Input instrument parameter workspace does not have parameter " << parname
948 << ". Information is incomplete for Monte Carlo simulation.\n";
949 g_log.error(errss.str());
950 throw runtime_error(errss.str());
951 }
952 vector<double> mcparvalues = mit->second;
953
954 // b) Build for the output
955 lowerbounds.emplace_back(mcparvalues[0]);
956 upperbounds.emplace_back(mcparvalues[1]);
957 stepsizes.emplace_back(mcparvalues[2]);
958 }
959}
960
979 const HistogramX &xVals, vector<double> &vec_n) {
980 if (m_Function->name() != "ThermalNeutronDtoTOFFunction") {
981 g_log.warning() << "Function (" << m_Function->name()
982 << " is not ThermalNeutronDtoTOFFunction. And it is not "
983 "required to calculate n.\n";
984 for (size_t i = 0; i < xVals.size(); ++i)
985 vec_n.emplace_back(0);
986 }
987
988 double width = m_Function->getParameter("Width");
989 double tcross = m_Function->getParameter("Tcross");
990
991 for (double dh : xVals) {
992 double n = 0.5 * gsl_sf_erfc(width * (tcross - 1 / dh));
993 vec_n.emplace_back(n);
994 }
995}
996
997//------- Related to Algorith's Output
998//-------------------------------------------------------------------
999
1004void RefinePowderInstrumentParameters::genPeakCentersWorkspace(bool montecarlo, size_t numbestfit) {
1005 // 1. Collect values in a vector for sorting
1006 double lattice = m_FuncParameters["LatticeConstant"];
1007 if (lattice < 1.0E-5) {
1008 std::stringstream errmsg;
1009 errmsg << "Input Lattice constant = " << lattice << " is wrong or not set up right. ";
1010 throw std::invalid_argument(errmsg.str());
1011 }
1012
1013 string stdoption = getProperty("StandardError");
1014 enum { ConstantValue, InvertedPeakHeight } stderroroption;
1015 if (stdoption == "ConstantValue") {
1016 stderroroption = ConstantValue;
1017 } else if (stdoption == "InvertedPeakHeight") {
1018 stderroroption = InvertedPeakHeight;
1019 } else {
1020 stringstream errss;
1021 errss << "Input StandardError (" << stdoption << ") is not supported. ";
1022 g_log.error(errss.str());
1023 throw runtime_error(errss.str());
1024 }
1025
1026 std::map<std::vector<int>, BackToBackExponential_sptr>::iterator peakiter;
1027 std::vector<std::pair<double, std::pair<double, double>>> peakcenters; // d_h [TOF_h, CHI2]
1028
1029 Geometry::UnitCell unitcell(lattice, lattice, lattice, 90.0, 90.0, 90.0);
1030
1031 for (peakiter = m_Peaks.begin(); peakiter != m_Peaks.end(); ++peakiter) {
1032 vector<int> hkl = peakiter->first;
1033 BackToBackExponential_sptr peak = peakiter->second;
1034
1035 double sigma = peak->getParameter("S");
1036 if (sigma < m_MinSigma) {
1037 g_log.information() << "Peak (" << hkl[0] << ", " << hkl[1] << ", " << hkl[2]
1038 << ") has unphysically small Sigma = " << sigma << ". "
1039 << "It is thus excluded. \n";
1040 continue;
1041 }
1042
1043 /* Replaced by UnitCell
1044 double dh = calculateDspaceValue(hkl, lattice);
1045 */
1046 double dh = unitcell.d(hkl[0], hkl[1], hkl[2]);
1047 double center = peak->centre();
1048 double height = peak->height();
1049 double chi2;
1050 if (stderroroption == ConstantValue) {
1051 chi2 = 1.0;
1052 } else if (stderroroption == InvertedPeakHeight) {
1053 chi2 = sqrt(1.0 / height);
1054 } else {
1055 throw runtime_error("Standard error option is not supported. ");
1056 }
1057
1058 peakcenters.emplace_back(dh, make_pair(center, chi2));
1059 }
1060
1061 // 2. Sort by d-spacing value
1062 std::sort(peakcenters.begin(), peakcenters.end());
1063
1064 // 3. Create output workspace
1065 size_t size = peakcenters.size();
1066 size_t nspec;
1067 if (montecarlo) {
1068 // Monte Carlo, 1 + 2*N spectra: raw, N x (refined, diff, n)
1069 nspec = 1 + 3 * numbestfit;
1070 } else {
1071 // Regular fit, 3 spectra: raw, refined, diff, n
1072 nspec = 1 + 3;
1073 }
1074
1075 m_dataWS = std::dynamic_pointer_cast<DataObjects::Workspace2D>(
1076 API::WorkspaceFactory::Instance().create("Workspace2D", nspec, size, size));
1077 m_dataWS->getAxis(0)->setUnit("dSpacing");
1078
1079 // 4. Put data to output workspace
1080 for (size_t i = 0; i < peakcenters.size(); ++i) {
1081 for (size_t j = 0; j < nspec; ++j) {
1082 m_dataWS->mutableX(j)[i] = peakcenters[i].first;
1083 }
1084 m_dataWS->mutableY(0)[i] = peakcenters[i].second.first;
1085 m_dataWS->mutableE(0)[i] = peakcenters[i].second.second;
1086 }
1087}
1088
1093 // 1. Create table workspace
1094 DataObjects::TableWorkspace_sptr tablews = std::make_shared<TableWorkspace>();
1095
1096 tablews->addColumn("double", "Chi2");
1097 tablews->addColumn("double", "GSLChi2");
1098 for (auto &peakFunctionParameterName : m_PeakFunctionParameterNames) {
1099 tablews->addColumn("double", peakFunctionParameterName);
1100 }
1101
1102 // 2. Put values in
1103 for (size_t ib = 0; ib < m_BestFitParameters.size(); ++ib) {
1104 TableRow newrow = tablews->appendRow();
1105 double chi2 = m_BestFitParameters[ib].first;
1106 double gslchi2 = m_BestFitChi2s[ib].second;
1107 newrow << chi2 << gslchi2;
1108 for (size_t ip = 0; ip < m_PeakFunctionParameterNames.size(); ++ip) {
1109 double tempdbl = m_BestFitParameters[ib].second[ip];
1110 newrow << tempdbl;
1111 }
1112 } // ENDFOR 1 Best Answer
1113
1114 return tablews;
1115}
1116
1123 // TableWorkspace is not copyable (default CC is incorrect and no point in
1124 // writing a non-default one)
1126 std::shared_ptr<DataObjects::TableWorkspace>(new DataObjects::TableWorkspace());
1127 newtablews->addColumn("str", "Name");
1128 newtablews->addColumn("double", "Value");
1129 newtablews->addColumn("str", "FitOrTie");
1130 newtablews->addColumn("double", "Min");
1131 newtablews->addColumn("double", "Max");
1132 newtablews->addColumn("double", "StepSize");
1133
1134 std::map<std::string, double>::iterator pariter;
1135
1136 for (pariter = m_FuncParameters.begin(); pariter != m_FuncParameters.end(); ++pariter) {
1137 API::TableRow newrow = newtablews->appendRow();
1138 std::string parname = pariter->first;
1139 double parvalue = pariter->second;
1140 newrow << parname << parvalue;
1141 }
1142
1143 return newtablews;
1144}
1145
1146} // namespace Mantid::CurveFitting::Algorithms
#define DECLARE_ALGORITHM(classname)
Definition Algorithm.h:538
double value
The value of the point.
Definition FitMW.cpp:51
double sigma
Definition GetAllEi.cpp:156
double height
Definition GetAllEi.cpp:155
void declareProperty(std::unique_ptr< Kernel::Property > p, const std::string &doc="") override
Add a property to the list of managed properties.
TypedValue getProperty(const std::string &name) const override
Get the value of a property.
virtual std::shared_ptr< Algorithm > createChildAlgorithm(const std::string &name, const double startProgress=-1., const double endProgress=-1., const bool enableLogging=true, const int &version=-1)
Create a Child Algorithm.
Kernel::Logger & g_log
Definition Algorithm.h:422
Implements FunctionDomain1D with its own storage in form of a std::vector.
size_t size() const override
Return the number of arguments in the domain.
A class to store values calculated by a function.
size_t size() const
Return the number of values.
TableRow represents a row in a TableWorkspace.
Definition TableRow.h:39
A property class for workspaces.
RefinePowderInstrumentParameters : Algorithm to refine instrument geometry parameters only.
void genPeaksFromTable(const DataObjects::TableWorkspace_sptr &peakparamws)
Generate peaks from table (workspace)
void importMonteCarloParametersFromTable(const DataObjects::TableWorkspace_sptr &tablews, const std::vector< std::string > &parameternames, std::vector< double > &stepsizes, std::vector< double > &lowerbounds, std::vector< double > &upperbounds)
Import the Monte Carlo related parameters from table.
std::map< std::string, double > m_OrigParameters
Map to store the original (input) parameters.
bool fitFunction(const API::IFunction_sptr &func, double &gslchi2)
Fit function to data.
std::vector< std::string > m_PeakFunctionParameterNames
Peak function parameter names.
double calculateD2TOFFunction(const API::IFunction_sptr &func, const API::FunctionDomain1DVector &domain, API::FunctionValues &values, const Mantid::HistogramData::HistogramY &rawY, const Mantid::HistogramData::HistogramE &rawE)
Calculate the value and chi2.
void importParametersFromTable(const DataObjects::TableWorkspace_sptr &parameterWS, std::map< std::string, double > &parameters)
Import instrument parameter from table (workspace)
void fitInstrumentParameters()
Fit instrument geometry parameters by ThermalNeutronDtoTOFFunction.
Functions::ThermalNeutronDtoTOFFunction_sptr m_Function
Modelling function.
std::map< std::string, double > m_FuncParameters
Map for function (instrument parameter)
std::vector< std::pair< double, std::vector< double > > > m_BestFitParameters
N sets of the peak parameter values for the best N chi2 for MC.
std::map< std::vector< int >, double > m_PeakErrors
Map for all peaks' error (fitted vs. experimental): [HKL]: Chi^2.
std::map< std::vector< int >, Functions::BackToBackExponential_sptr > m_Peaks
Map for all peaks to fit individually.
double calculateFunctionStatistic(const API::IFunction_sptr &func, const API::MatrixWorkspace_sptr &dataws, size_t workspaceindex)
Calculate function's statistic.
DataObjects::Workspace2D_sptr m_dataWS
Output Workspace containing the dspacing ~ TOF peak positions.
DataObjects::TableWorkspace_sptr genMCResultTable()
Generate an output table workspace for N best fitting result.
void doParameterSpaceRandomWalk(const std::vector< std::string > &parnames, const std::vector< double > &lowerbounds, const std::vector< double > &upperbounds, std::vector< double > &stepsizes, size_t maxsteps, double stepsizescalefactor, bool fit2)
Core Monte Carlo random walk on parameter-space.
void calculateThermalNeutronSpecial(const API::IFunction_sptr &m_Function, const Mantid::HistogramData::HistogramX &xVals, std::vector< double > &vec_n)
Calculate d-space value from peak's miller index for thermal neutron.
void getD2TOFFuncParamNames(std::vector< std::string > &parnames)
Get the names of the parameters of D-TOF conversion function.
std::vector< std::pair< double, std::vector< double > > > m_BestMCParameters
N sets of the peak parameter values for the best N chi2 for MC.
void genPeakCentersWorkspace(bool montecarlo, size_t numbestfit)
Generate (output) workspace of peak centers.
std::vector< std::pair< double, double > > m_BestFitChi2s
N sets of the homemade chi2 and gsl chi2.
DataObjects::TableWorkspace_sptr genOutputInstrumentParameterTable()
Generate (output) table worksspace for instrument parameters.
void refineInstrumentParametersMC(const DataObjects::TableWorkspace_sptr &parameterWS, bool fit2=false)
Set up and run a monte carlo simulation to refine the peak parameters.
TableWorkspace is an implementation of Workspace in which the data are organised in columns of same s...
Class to implement unit cell of crystals.
Definition UnitCell.h:44
double d(double h, double k, double l) const
Return d-spacing ( ) for a given h,k,l coordinate.
Definition UnitCell.cpp:700
Support for a property that holds an array of values.
IPropertyManager * setProperty(const std::string &name, const T &value)
Templated method to set the value of a PropertyWithValue.
void debug(const std::string &msg)
Logs at debug level.
Definition Logger.cpp:145
void notice(const std::string &msg)
Logs at notice level.
Definition Logger.cpp:126
void error(const std::string &msg)
Logs at error level.
Definition Logger.cpp:108
void warning(const std::string &msg)
Logs at warning level.
Definition Logger.cpp:117
void information(const std::string &msg)
Logs at information level.
Definition Logger.cpp:136
static T & Instance()
Return a reference to the Singleton instance, creating it if it does not already exist Creation is do...
std::shared_ptr< IFunction > IFunction_sptr
shared pointer to the function base class
Definition IFunction.h:743
std::shared_ptr< MatrixWorkspace > MatrixWorkspace_sptr
shared pointer to the matrix workspace base class
std::shared_ptr< ThermalNeutronDtoTOFFunction > ThermalNeutronDtoTOFFunction_sptr
std::shared_ptr< BackToBackExponential > BackToBackExponential_sptr
std::shared_ptr< TableWorkspace > TableWorkspace_sptr
shared pointer to Mantid::DataObjects::TableWorkspace
std::unique_ptr< T > create(const P &parent, const IndexArg &indexArg, const HistArg &histArg)
This is the create() method that all the other create() methods call.
std::vector< double > getZscore(const std::vector< TYPE > &data)
Return the Z score values for a dataset.
Helper class which provides the Collimation Length for SANS instruments.
STL namespace.
@ InOut
Both an input & output workspace.
Definition Property.h:55
@ Input
An input workspace.
Definition Property.h:53
@ Output
An output workspace.
Definition Property.h:54