REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
src/TRestDetectorGas.cxx
1/*************************************************************************
2 * This file is part of the REST software framework. *
3 * *
4 * Copyright (C) 2016 GIFNA/TREX (University of Zaragoza) *
5 * For more information see http://gifna.unizar.es/trex *
6 * *
7 * REST is free software: you can redistribute it and/or modify *
8 * it under the terms of the GNU General Public License as published by *
9 * the Free Software Foundation, either version 3 of the License, or *
10 * (at your option) any later version. *
11 * *
12 * REST is distributed in the hope that it will be useful, *
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
15 * GNU General Public License for more details. *
16 * *
17 * You should have a copy of the GNU General Public License along with *
18 * REST in $REST_PATH/LICENSE. *
19 * If not, see http://www.gnu.org/licenses/. *
20 * For the list of contributors see $REST_PATH/CREDITS. *
21 *************************************************************************/
22
215
216#include "TRestDetectorGas.h"
217
218#include <TRestDataBase.h>
219
220#include <algorithm>
221
222using namespace std;
223
224ClassImp(TRestDetectorGas);
225
226// const char* defaultServer = "https://rest-for-physics.github.io/gas-files/";
227
232 Initialize();
233
234 fGasGeneration = false;
235}
236
249TRestDetectorGas::TRestDetectorGas(const char* configFilename, const string& name, bool gasGeneration,
250 bool test)
252 Initialize();
253 fGasGeneration = gasGeneration;
254
255 fTest = test;
256
257 if (strcmp(configFilename, "server") == 0) {
258 LoadConfigFromElement(StringToElement("<TRestDetectorGas name=\"" + name + R"(" file="server"/>)"),
259 nullptr);
260 } else {
261 fConfigFileName = configFilename;
263 }
264
265 // if ( fStatus == RESTGAS_CFG_LOADED ) LoadGasFile( );
266}
267
271TRestDetectorGas::~TRestDetectorGas() {
272 RESTDebug << "Entering ... TRestDetectorGas() destructor." << RESTendl;
273
274#if defined USE_Garfield
275 delete fGasMedium;
276#endif
277}
278
284 RESTDebug << "TRestDetectorGas. Entering ... Initialize()." << RESTendl;
285
286 SetSectionName(this->ClassName());
287 SetLibraryVersion(LIBRARY_VERSION);
288
289 fPressureInAtm = 1;
290 fTemperatureInK = 300;
291
292 fNofGases = 0;
293
294 fGasComponentName.clear();
295 fGasComponentFraction.clear();
296
297 fStatus = RESTGAS_INTITIALIZED;
298
299 fGasFilename = "";
300 fGasFileContent = "";
301
302#if defined USE_Garfield
303 fGasMedium = new Garfield::MediumMagboltz();
304#else
305 fGasMedium = nullptr;
306#endif
307
309 // This must be commented. If not when we specify gasGeneration=true on the
310 // constructor, it will be overriden inside LoadConfigFromFile
311 //
312 // fGasGeneration = false;
314
315 fEmin = 10;
316 fEmax = 1000;
317 fEnodes = 20;
318}
319
328void TRestDetectorGas::LoadGasFile() {
329 RESTDebug << "Entering ... TRestDetectorGas::LoadGasFile()." << RESTendl;
330
331#if defined USE_Garfield
332 RESTDebug << "fGasFilename = " << fGasFilename << RESTendl;
333 if (!TRestTools::fileExists((string)(fGasFilename))) {
334 RESTError << __PRETTY_FUNCTION__ << RESTendl;
335 RESTError << "The gas file does not exist. (name:" << fGasFilename << ")" << RESTendl;
336 fStatus = RESTGAS_ERROR;
337 return;
338 }
339
340 fGasMedium->LoadGasFile((string)(fGasFilename));
341
342 fEFields.clear();
343 fBFields.clear();
344 fAngles.clear();
345 fGasMedium->GetFieldGrid(fEFields, fBFields, fAngles);
346
347 fStatus = RESTGAS_GASFILE_LOADED;
348 RESTInfo << "TRestDetectorGas. Gas file loaded!" << RESTendl;
349
350 for (unsigned int i = 0; i < fEFields.size(); i++)
351 RESTDebug << "node " << i << " Field : " << fEFields[i] << " V/cm" << RESTendl;
352
353 if (fGasMedium && fGasMedium->GetW() == 0.) {
354 fGasMedium->SetW(GetWvalue());
355 } // as it is probably not computed by Magboltz
356#else
357 cout << "This REST is not complied with garfield, it cannot load any gas file!" << endl;
358#endif
359}
360
361void TRestDetectorGas::CalcGarField(double Emin, double Emax, int n) {
362 RESTDebug << "Entering ... TRestDetectorGas::CalcGarField( Emin=" << Emin << " , Emax=" << Emax << " )"
363 << RESTendl;
364
365#if defined USE_Garfield
366 if (fEnodes <= 0) {
367 cout << "REST ERROR : The number of nodes is not a positive number!!. Gas "
368 "file generation cancelled."
369 << endl;
370 fStatus = RESTGAS_ERROR;
371 return;
372 }
373 if (fEmin >= fEmax) {
374 cout << "REST ERROR : The Electric field grid boundaries are not properly "
375 "defined."
376 << endl;
377 fStatus = RESTGAS_ERROR;
378 return;
379 }
380
381 string gasStr[3];
382 for (int i = 0; i < fNofGases; i++) {
383 gasStr[i] = (string)fGasComponentName[i];
384 if (i == 2) break;
385 }
386
387 if (fNofGases == 1) fGasMedium->SetComposition(gasStr[0], fGasComponentFraction[0] * 100.);
388
389 if (fNofGases == 2)
390 fGasMedium->SetComposition(gasStr[0], fGasComponentFraction[0] * 100., gasStr[1],
391 fGasComponentFraction[1] * 100.);
392
393 if (fNofGases == 3)
394 fGasMedium->SetComposition(gasStr[0], fGasComponentFraction[0] * 100., gasStr[1],
395 fGasComponentFraction[1] * 100., gasStr[2],
396 fGasComponentFraction[2] * 100.);
397
398 if (fNofGases > 3) {
399 cout << "REST ERROR : Number of gas components higher than 3 not allowed" << endl;
400 fStatus = RESTGAS_ERROR;
401 return;
402 }
403
404 fGasMedium->SetTemperature(fTemperatureInK);
405
406 if (fPressureInAtm != 1)
407 RESTWarning << "-- Warning : The gas will be generated for gas pressure = 1atm" << RESTendl;
408
409 fGasMedium->SetPressure(1 * REST_Units::torr);
410
411 fGasMedium->SetFieldGrid(Emin, Emax, n, n > 1);
412
413 fGasMedium->SetMaxElectronEnergy(fMaxElectronEnergy);
414
415 cout << "Garfield: calculating..." << endl;
416
417 if (fVerboseLevel >= TRestStringOutput::REST_Verbose_Level::REST_Info) fGasMedium->EnableDebugging();
418 fGasMedium->Initialise();
419 if (fVerboseLevel >= TRestStringOutput::REST_Verbose_Level::REST_Info) fGasMedium->DisableDebugging();
420
421 fGasMedium->GenerateGasTable(fNCollisions, true);
422 if (fPressureInAtm != 1) {
423 RESTWarning << "-- Warning : Restoring the gas pressure" << RESTendl;
424 fGasMedium->SetPressure(fPressureInAtm * 760.);
425 }
426#else
427 cout << "This REST is not complied with garfield, it cannot calculate "
428 "garfield!"
429 << endl;
430#endif
431}
432
442void TRestDetectorGas::AddGasComponent(const string& gasName, Double_t fraction) {
443 RESTDebug << "Entering ... TRestDetectorGas::AddGasComponent( gasName=" << gasName
444 << " , fraction=" << fraction << " )" << RESTendl;
445
446 fGasComponentName.emplace_back(gasName);
447 fGasComponentFraction.push_back(fraction);
448 fNofGases++;
449}
450
451/*
452// This was just a test to try to Get the calculated W for the gas definition.
453// However, I tested with Xe+TMA and I got an error message that TMA
454// photoncrossection database is not available
455
456void TRestDetectorGas::GetGasWorkFunction() {
457#if defined USE_Garfield
458 RESTEssential << __PRETTY_FUNCTION__ << RESTendl;
459 RESTEssential << "This method has never been validated to operate properly" << RESTendl;
460 RESTEssential << "If we manage to make it work we could use this method to "
461 "obtain the calculated W of the gas"
462 << RESTendl;
463
464 // Gas gap [cm].
465 const double width = 1.;
466 SolidBox* box = new SolidBox(width / 2., 0., 0., width / 2., 10., 10.);
467 GeometrySimple* geo = new GeometrySimple();
468 geo->AddSolid(box, fGasMedium);
469
470 ComponentConstant* cmp = new ComponentConstant();
471 cmp->SetGeometry(geo);
472 cmp->SetElectricField(100., 0., 0.);
473
474 Sensor* sensor = new Sensor();
475 sensor->AddComponent(cmp);
476
477 TrackHeed* heed = new TrackHeed();
478 heed->SetSensor(sensor);
479 // Set the particle type.
480 heed->SetParticle("pi");
481 // Set the particle momentum (in eV/c).
482 heed->SetMomentum(120.e9);
483
484 // Switch on debugging to print out some information (stopping power, W value,
485 // ...)
486 heed->EnableDebugging();
487 // Initial position
488 double x0 = 0., y0 = 0., z0 = 0., t0 = 0.;
489 // Direction of the track (perpendicular incidence)
490 double dx0 = 1., dy0 = 0., dz0 = 0.;
491 heed->NewTrack(x0, y0, z0, t0, dx0, dy0, dz0);
492 cout << "W : " << heed->GetW() << endl;
493#else
494 cout << "This REST is not complied with garfield, it cannot calculate "
495 "garfield!"
496 << endl;
497#endif
498}
499*/
500
501// Get the fano factor from Garfield::MediumMagboltz
502// User need to have installed the last version of
503// Garfield to this to work
504
505Double_t TRestDetectorGas::GetGasFanoFactor() const {
506#if defined USE_Garfield
507 if (fStatus != RESTGAS_GASFILE_LOADED) {
508 RESTDebug << "-- Error : " << __PRETTY_FUNCTION__ << RESTendl;
509 RESTDebug << "-- Error : Gas file was not loaded!" << RESTendl;
510 return 0;
511 }
512
513 RESTInfo << "Calling Garfield directly to fetch Fano factor" << RESTendl;
514 const auto fanoFactor = fGasMedium->GetFanoFactor();
515
516 if (fanoFactor == 0.) {
517 runtime_error("Fano Factor is 0! This REST is not compiled with the last version of Garfield");
518 }
519 return fanoFactor;
520#else
521 cerr << "This REST is not compiled with garfield, Do not use Fano "
522 "Factor from TRestDetectorGas!"
523 << endl
524 << "Please define the Fano factor in each process!" << endl;
525 throw runtime_error("This REST is not compiled with garfield, cannot retrieve the Fano Factor");
526#endif
527}
528
530#if defined USE_Garfield
531 if (fStatus != RESTGAS_GASFILE_LOADED) {
532 RESTDebug << "-- Error : " << __PRETTY_FUNCTION__ << RESTendl;
533 RESTDebug << "-- Error : Gas file was not loaded!" << RESTendl;
534 return 0;
535 }
536
537 RESTInfo << "Calling Garfield directly to fetch the work function" << RESTendl;
538 const auto workFunction = fGasMedium->GetW();
539
540 if (workFunction == 0.) {
541 throw runtime_error("Work Function is 0! This should never happen");
542 }
543 return workFunction;
544#else
545 throw runtime_error(
546 "This REST is not compiled with garfield, cannot retrieve the work function from the gas");
547#endif
548}
549
559 // if (GetVerboseLevel() <= REST_Info) fVerboseLevel = REST_Info;
560
561 RESTDebug << "Entering ... TRestDetectorGas::InitFromConfigFile()" << RESTendl;
562
563 // read config parameters of base class
565 if (fW == -1) {
566 fW = 21.9;
567 cout << "Setting default W-value : " << fW << endl;
568 }
569
570 // match the database, id=0(any), type="GAS_SERVER"
571 string _gasServer = gDataBase->query_data(DBEntry(0, "GAS_SERVER")).value;
572 if (_gasServer == "") _gasServer = "none";
573 fGasServer = GetParameter("gasServer", _gasServer);
574
575 // add gas component
576 TiXmlElement* gasComponentDefinition = GetElement("gasComponent");
577 while (gasComponentDefinition != nullptr) {
578 string gasName = GetFieldValue("name", gasComponentDefinition);
579 Double_t gasFraction = StringToDouble(GetFieldValue("fraction", gasComponentDefinition));
580 AddGasComponent(gasName, gasFraction);
581 gasComponentDefinition = GetNextElement(gasComponentDefinition);
582 }
583 if (fNofGases == 0 && !fMaterial.empty()) {
584 vector<string> componentsdef = Split(fMaterial, " ");
585 for (const auto& componentdef : componentsdef) {
586 vector<string> componentdefpair = Split(componentdef, ":");
587 if (componentdefpair.size() != 2) {
588 continue;
589 }
590 AddGasComponent(componentdefpair[0], StringToDouble(componentdefpair[1]));
591 }
592 }
593 if (fNofGases == 0) {
594 RESTError << "TRestDetectorGas: No gas components added!" << RESTendl;
595 }
596 double sum = 0;
597 for (int i = 0; i < fNofGases; i++) sum += GetGasComponentFraction(i);
598 if (sum - 1 < 1.e12)
599 fStatus = RESTGAS_CFG_LOADED;
600 else {
601 RESTWarning << "TRestDetectorGas : The total gas fractions is NOT 1." << RESTendl;
602 fStatus = RESTGAS_ERROR;
603 return;
604 }
605
606 // setup e-field calculation range and gas file generation parameters
607 TiXmlElement* eFieldDefinition = GetElement("eField");
608 fEmax = StringToDouble(GetFieldValue("Emax", eFieldDefinition));
609 fEmin = StringToDouble(GetFieldValue("Emin", eFieldDefinition));
610 fEnodes = StringToInteger(GetFieldValue("nodes", eFieldDefinition));
611 fNCollisions = StringToInteger(GetParameter("nCollisions"));
612 fMaxElectronEnergy = StringToDouble(GetParameter("maxElectronEnergy", "40"));
613 if (ToUpper(GetParameter("generate")) == "ON" || ToUpper(GetParameter("generate")) == "TRUE")
614 fGasGeneration = true;
615 fGasOutputPath = GetParameter("gasOutputPath", "./");
617 RESTWarning << "-- Warning : The specified gasOutputPath is not writable!" << RESTendl;
618 RESTWarning << "-- Warning : The output path will be changed to ./" << RESTendl;
619 fGasOutputPath = "./";
620 }
621 fGDMLMaterialRef = GetParameter("GDMLMaterialRef", "");
622
623 // construct the gas file name and try to find it, either locally or from gas server
624 fGasFilename = ConstructFilename();
625 RESTDebug << "TRestDetectorGas::InitFromConfigFile. ConstructFilename. fGasFilename = " << fGasFilename
626 << RESTendl;
627 fGasFilename = FindGasFile((string)fGasFilename);
628 RESTDebug << "TRestDetectorGas::InitFromConfigFile. FindGasFile. fGasFilename = " << fGasFilename
629 << RESTendl;
630
631 // If we found the gasFile then obviously we disable the gas generation
632 if (fGasGeneration && TRestTools::fileExists((string)fGasFilename)) {
633 fGasGeneration = false;
634
635 RESTWarning << "TRestDetectorGas gasFile generation is enabled, but the "
636 "gasFile already exists!!"
637 << RESTendl;
638 RESTWarning << "Filename : " << fGasFilename << RESTendl;
639 RESTWarning << "fGasGeneration should be disabled to remove this "
640 "warning."
641 << RESTendl;
642 RESTWarning << "If you really want to re-generate the gas file you "
643 "will need to disable the gasServer."
644 << RESTendl;
645 RESTWarning << "And/or remove any local copies that are found by "
646 "SearchPath."
647 << RESTendl;
648 }
649 fStatus = RESTGAS_CFG_LOADED;
650
651#if defined USE_Garfield
652 // calling garfield, either to generate gas file or load existing gas file
653 if (fGasGeneration) {
654 RESTEssential << "Starting gas generation" << RESTendl;
655
656 CalcGarField(fEmin, fEmax, fEnodes);
657 GenerateGasFile();
658 fStatus = RESTGAS_GASFILE_LOADED;
659 } else {
660 LoadGasFile();
661 fGasMedium->SetPressure(fPressureInAtm * REST_Units::torr);
662 }
663 if (fGasMedium && fGasMedium->GetW() == 0.)
664 fGasMedium->SetW(fW); // as it is probably not computed by Magboltz
665#endif
666
668}
669
671 RESTDebug << "Entering ... TRestDetectorGas::InitFromRootFile()" << RESTendl;
673 if (fGasFileContent != "") // use gas file content by default
674 {
675 fGasFilename = "/tmp/restGasFile_" + REST_USER + ".gas";
676 ofstream outf;
677 outf.open(fGasFilename, ios::ate);
678 outf << fGasFileContent << endl;
679 outf.close();
680 LoadGasFile();
681 int z = system("rm " + fGasFilename);
682 if (z != 0) RESTError << "Problem removing gas file: " << fGasFilename << RESTendl;
683 } else {
684 fGasFilename = FindGasFile((string)fGasFilename);
685 if (fGasFilename != "") {
686 LoadGasFile();
687 }
688 }
689}
690
691void TRestDetectorGas::UploadGasToServer(string absoluteGasFilename) {
692 RESTEssential << "uploading gas file and gas definition rmls" << RESTendl;
693
694 if (!fTest && (fMaxElectronEnergy < 400 || fNCollisions < 10 || fEnodes < 20)) {
695 RESTWarning << "-- Warning : The gas file does not fulfill the requirements "
696 "for being uploaded to the gasServer"
697 << RESTendl;
698 RESTWarning << "-- Warning : maxElectronEnergy >= 400. Number of collisions >= "
699 "10. Number of E nodes >= 20."
700 << RESTendl;
701 RESTWarning << "-- Warning : The generated file will NOT be uploaded to the "
702 "server but preserved locally."
703 << RESTendl;
704 return;
705 }
706
708 // We add the gas definition we used to generate the gas file and prepare it
709 // to upload/update in the gasServer
710 string cmd;
711 int a;
712 // We download (probably again) the original version
713 string url = gDataBase->query_data(DBEntry(0, "META_RML", "TRestDetectorGas")).value;
714 string fname = TRestTools::DownloadRemoteFile(url);
715
716// We remove the last line. I.e. the enclosing </gases> in the original file
717#ifdef __APPLE__
718 cmd = "sed -i '' -e '$ d' " + fname;
719#else
720 cmd = "sed -i '$ d' " + fname;
721#endif
722
723 a = system(cmd.c_str());
724
725 if (a != 0) {
726 RESTError << "-- Error : " << __PRETTY_FUNCTION__ << RESTendl;
727 RESTError << "-- Error : problem removing last line from " << fname << RESTendl;
728 return;
729 }
730
731 // We add some header before the gas definition. We might add also date an
732 // other information essential to identify the gasFile submission
733 ofstream outf;
734 outf.open(fname, ios::app);
735 outf << endl;
736 outf << "//------- User : " << REST_USER << " ---- REST version : " << REST_RELEASE
737 << " ---------------------------" << endl;
738 outf.close();
739
740 // We write the TRestDetectorGas section
741 this->WriteConfigBuffer(fname);
742
743 // We re-write the enclosing </gases> tag
744 outf.open(fname, ios::app);
745 outf << "\n" << endl;
746 outf << "</gases>" << endl;
747 outf.close();
749
750 // We transfer the new gas definitions to the gasServer
751 TRestTools::UploadToServer(fname, (string)fGasServer, "ssh://gasUser@:22");
752
753 // We transfer the gasFile to the gasServer
754 TRestTools::UploadToServer(absoluteGasFilename, (string)fGasServer, "ssh://gasUser@:22");
755
756 // We remove the local file (afterwards, the remote copy will be used)
757 // cmd = "rm " + _name;
758 // a = system(cmd.c_str());
759
760 // if (a != 0) {
761 // ferr << "-- Error : " << __PRETTY_FUNCTION__ << endl;
762 // ferr << "-- Error : problem removing the locally generated gas file" << endl;
763 // ferr << "-- Error : Please report this problem at "
764 // "http://gifna.unizar.es/rest-forum/"
765 // << endl;
766 // return;
767 //}
768
769 RESTSuccess << "-- Success : Gasfile server database was updated successfully!" << RESTendl;
770}
771
786string TRestDetectorGas::FindGasFile(string name) {
787 RESTDebug << "Entering ... TRestDetectorGas::FindGasFile( name=" << name << " )" << RESTendl;
788
789 string absoluteName = "";
790
791 if (!fGasGeneration && fGasServer != "none") {
792 absoluteName = TRestTools::DownloadRemoteFile((string)fGasServer + "/" + name, true);
793 }
794
795 if (absoluteName.empty()) {
796 RESTInfo << "Trying to find the gasFile locally" << RESTendl;
797 absoluteName = SearchFile(name);
798 if (absoluteName.empty()) {
799 RESTWarning << "-- Warning : No sucess finding local gas file definition." << RESTendl;
800 RESTWarning << "-- Warning : Gas file definition does not exist." << RESTendl;
801 RESTInfo << "To generate a new gasFile enable gas generation in TRestDetectorGas "
802 "constructor"
803 << RESTendl;
804 RESTInfo << R"(TRestDetectorGas ( "gasDefinition.rml", "gas Name", true );)" << RESTendl;
805 RESTInfo << "Further details can be found at TRestDetectorGas class definition and "
806 "tutorial."
807 << RESTendl;
808
809 absoluteName = name;
810 }
811 }
812
813 return absoluteName;
814}
815
819TString TRestDetectorGas::GetGasMixture() {
820 RESTDebug << "Entering ... TRestDetectorGas::GetGasMixture( )" << RESTendl;
821
822 TString gasMixture;
823 char tmpStr[64];
824 for (int n = 0; n < fNofGases; n++) {
825 if (n > 0) gasMixture += "-";
826 gasMixture += GetGasComponentName(n) + "_";
827 sprintf(tmpStr, "%03.1lf", GetGasComponentFraction(n) * 100.);
828 gasMixture += (TString)tmpStr;
829 }
830 return gasMixture;
831}
832
839string TRestDetectorGas::ConstructFilename() {
840 RESTDebug << "Entering ... TRestDetectorGas::ConstructFilename( )" << RESTendl;
841
842 string name;
843 char tmpStr[256];
844 for (int n = 0; n < fNofGases; n++) {
845 if (n > 0) name += "-";
846 name += GetGasComponentName(n) + "_";
847 if (GetGasComponentFraction(n) >= 0.001)
848 sprintf(tmpStr, "%03.1lf", GetGasComponentFraction(n) * 100.);
849 else
850 sprintf(tmpStr, "%03.1lfppm", GetGasComponentFraction(n) * 1.e6);
851
852 name += (TString)tmpStr;
853 }
854
855 // The filename is constructed always at 1 atm pressure.
856 // We keep E_vs_P to remind the field calculation range will
857 // depend on pressure.
858
859 name += "-E_vs_P_";
860 sprintf(tmpStr, "%03.1lf", fEmin);
861 name += (TString)tmpStr;
862
863 name += "_";
864 sprintf(tmpStr, "%03.1lf", fEmax);
865 name += (TString)tmpStr;
866
867 name += "_nodes_";
868 sprintf(tmpStr, "%02d", fEnodes);
869 name += (TString)tmpStr;
870
871 name += "-nCol_";
872 sprintf(tmpStr, "%02d", fNCollisions);
873 name += (TString)tmpStr;
874
875 name += "-maxE_";
876 sprintf(tmpStr, "%03d", (Int_t)fMaxElectronEnergy);
877 name += (TString)tmpStr;
878
879 name += ".gas";
880
881 RESTDebug << "Constructed filename : " << name << RESTendl;
882 return name;
883}
884
887void TRestDetectorGas::GenerateGasFile() {
888 RESTDebug << "Entering ... TRestDetectorGas::GenerateGasFile( )" << RESTendl;
889
890#if defined USE_Garfield
891
892 fGasFilename = ConstructFilename();
893 RESTDebug << " TRestDetectorGas::GenerateGasFile. fGasFilename = " << fGasFilename << RESTendl;
894
896 cout << endl;
897 RESTWarning << "-- Warning: REST ERROR. TRestDetectorGas. Path is not writtable." << RESTendl;
898 RESTWarning << "-- Warning: Path : " << fGasOutputPath << RESTendl;
899 RESTWarning << "-- Warning: Make sure the final data path is writtable before "
900 "proceed to gas generation."
901 << RESTendl;
902 RESTWarning << "-- Warning: or change the gas data path ... " << RESTendl;
903 RESTWarning << RESTendl;
904 GetChar();
905 return;
906 }
907
908 cout << "Writing gas file : " << endl;
909 cout << "-----------------" << endl;
910 cout << "Path : " << fGasOutputPath << endl;
911 cout << "Filename : " << fGasFilename << endl;
912
913 fGasMedium->WriteGasFile((string)(fGasOutputPath + "/" + fGasFilename));
914
915 if (fGasServer != "none") UploadGasToServer((string)(fGasOutputPath + "/" + fGasFilename));
916
917#else
918 cout << "This REST is not complied with garfield, it cannot save any gas file!" << endl;
919#endif
920}
921
931void TRestDetectorGas::SetPressure(Double_t pressure) {
932 RESTDebug << "Entering ... TRestDetectorGas::SetPressure( pressure=" << pressure << " )" << RESTendl;
933
934 fPressureInAtm = pressure;
935#if defined USE_Garfield
936 fGasMedium->SetPressure(fPressureInAtm * REST_Units::torr);
937#endif
938}
939
942void TRestDetectorGas::SetTemperature(Double_t temperature) {
943 RESTDebug << "Entering ... TRestDetectorGas::SetPressure( temperature=" << temperature << " )"
944 << RESTendl;
945
946 fTemperatureInK = temperature;
947#if defined USE_Garfield
948 fGasMedium->SetTemperature(temperature);
949#endif
950}
951
960void TRestDetectorGas::PlotDriftVelocity(Double_t eMin, Double_t eMax, Int_t nSteps) {
961 RESTDebug << "Entering ... TRestDetectorGas::PlotDriftVelocity( eMin=" << eMin << " , eMax=" << eMax
962 << ", nSteps=" << nSteps << " )" << RESTendl;
963
964 vector<Double_t> eField(nSteps), driftVel(nSteps);
965
966 for (int i = 0; i < nSteps; i++) {
967 eField[i] = (eMin + (double)i * (eMax - eMin) / nSteps);
968
969 this->SetElectricField(eField[i] / units("V/cm"));
970 driftVel[i] = GetDriftVelocity() * units("cm/us");
971 }
972
973 auto c = new TCanvas("Drift velocity", " ");
974 auto fDriftVel = new TGraph(nSteps, &eField[0], &driftVel[0]);
975 TString str;
976 str.Form("Drift Velocity for %s (Pressure: %3.1lf bar)", GetName(), this->GetPressure());
977 fDriftVel->SetTitle(str);
978 fDriftVel->GetXaxis()->SetTitle("E [V/cm]");
979 fDriftVel->GetYaxis()->SetTitle("Drift velocity [cm/#mus]");
980 fDriftVel->GetYaxis()->SetTitleOffset(2);
981 fDriftVel->Draw("");
982 c->Update();
983}
984
993void TRestDetectorGas::PlotLongitudinalDiffusion(Double_t eMin, Double_t eMax, Int_t nSteps) {
994 RESTDebug << "Entering ... TRestDetectorGas::PlotLongitudinalDiffusion( eMin=" << eMin
995 << " , eMax=" << eMax << ", nSteps=" << nSteps << " )" << RESTendl;
996
997 vector<Double_t> eField(nSteps), longDiff(nSteps);
998
999 for (int i = 0; i < nSteps; i++) {
1000 eField[i] = eMin + (double)i * (eMax - eMin) / nSteps;
1001
1002 this->SetElectricField(eField[i] / units("V/cm"));
1003 longDiff[i] = 10. * GetLongitudinalDiffusion(); // to express it in mm/sqrt(cm)
1004 }
1005
1006 auto c = new TCanvas("Longitudinal diffusion", " ");
1007 auto fLongDiff = new TGraph(nSteps, &eField[0], &longDiff[0]);
1008 TString str;
1009 str.Form("Longitudinal diffusion for %s (Pressure: %3.1lf bar)", GetName(), this->GetPressure());
1010 fLongDiff->SetTitle(str);
1011 fLongDiff->GetXaxis()->SetTitle("E [V/cm]");
1012 fLongDiff->GetYaxis()->SetTitle("Longitudinal diffusion [mm/#sqrt{cm}]");
1013 fLongDiff->GetYaxis()->SetTitleOffset(2);
1014 fLongDiff->Draw("");
1015 c->Update();
1016}
1017
1026void TRestDetectorGas::PlotTransversalDiffusion(Double_t eMin, Double_t eMax, Int_t nSteps) {
1027 RESTDebug << "Entering ... TRestDetectorGas::PlotTransversalDiffusion( eMin=" << eMin
1028 << " , eMax=" << eMax << ", nSteps=" << nSteps << " )" << RESTendl;
1029
1030 vector<Double_t> eField(nSteps), transDiff(nSteps);
1031
1032 for (int i = 0; i < nSteps; i++) {
1033 eField[i] = eMin + (double)i * (eMax - eMin) / nSteps;
1034
1035 this->SetElectricField(eField[i] / units("V/cm"));
1036 transDiff[i] = 10. * GetTransversalDiffusion(); // to express it in mm/sqrt(cm)
1037 }
1038
1039 auto c = new TCanvas("Transitudinal diffusion", " ");
1040 auto fTransDiff = new TGraph(nSteps, &eField[0], &transDiff[0]);
1041 TString str;
1042 str.Form("Transversal diffusion for %s (Pressure: %3.1lf bar)", GetName(), this->GetPressure());
1043 fTransDiff->SetTitle(str);
1044 fTransDiff->GetXaxis()->SetTitle("E [V/cm]");
1045 fTransDiff->GetYaxis()->SetTitle("Transversal diffusion [mm/#sqrt{cm}]");
1046 fTransDiff->GetYaxis()->SetTitleOffset(2);
1047 fTransDiff->Draw("");
1048 c->Update();
1049}
1050
1059void TRestDetectorGas::PlotTownsendCoefficient(Double_t eMin, Double_t eMax, Int_t nSteps) {
1060 RESTDebug << "Entering ... TRestDetectorGas::PlotTownsendCoefficient( eMin=" << eMin << " , eMax=" << eMax
1061 << ", nSteps=" << nSteps << " )" << RESTendl;
1062
1063 vector<Double_t> eField(nSteps), townsendCoeff(nSteps);
1064
1065 for (int i = 0; i < nSteps; i++) {
1066 eField[i] = eMin + (double)i * (eMax - eMin) / nSteps;
1067
1068 townsendCoeff[i] = GetTownsendCoefficient(eField[i]);
1069 }
1070
1071 auto c = new TCanvas("Townsend coefficient", " ");
1072 auto fTownsend = new TGraph(nSteps, &eField[0], &townsendCoeff[0]);
1073 TString str;
1074 str.Form("Townsend coefficient for %s", GetName());
1075 fTownsend->SetTitle(str);
1076 fTownsend->GetXaxis()->SetTitle("E [V/cm]");
1077 fTownsend->GetYaxis()->SetTitle("Townsend coefficient [1/cm]");
1078 fTownsend->GetYaxis()->SetTitleOffset(2);
1079 fTownsend->Draw("");
1080 c->Update();
1081}
1082
1087Double_t TRestDetectorGas::GetDriftVelocity(Double_t E) {
1088 RESTDebug << "Entering ... TRestDetectorGas::GetDriftVelocity( E=" << E << " )" << RESTendl;
1089
1090 SetPressure(fPressureInAtm);
1091
1092#if defined USE_Garfield
1093 if (fStatus != RESTGAS_GASFILE_LOADED) {
1094 RESTDebug << "-- Error : " << __PRETTY_FUNCTION__ << RESTendl;
1095 RESTDebug << "-- Error : Gas file was not loaded!" << RESTendl;
1096 return 0;
1097 }
1098
1099 RESTInfo << "Calling Garfield directly. Please be aware that the unit is different "
1100 << "from REST standard unit. E is V/cm. The return is cm/us" << RESTendl;
1101
1102 Double_t vx, vy, vz;
1103 fGasMedium->ElectronVelocity(0., 0, -E, 0, 0, 0, vx, vy, vz);
1104 return vz * 1000.;
1105#else
1106 std::cout << "This REST is not complied with garfield, Do not use Drift Velocity "
1107 "from TRestDetectorGas!"
1108 << std::endl;
1109 std::cout << "Please define the Drift Velocity in each process!" << std::endl;
1110 return 0.001;
1111#endif
1112}
1113
1118Double_t TRestDetectorGas::GetLongitudinalDiffusion(Double_t E) {
1119 RESTDebug << "Entering ... TRestDetectorGas::GetLongitudinalDiffusion( E=" << E << " )" << RESTendl;
1120
1121 SetPressure(fPressureInAtm);
1122
1123#if defined USE_Garfield
1124 if (fStatus != RESTGAS_GASFILE_LOADED) {
1125 RESTDebug << "-- Error : " << __PRETTY_FUNCTION__ << RESTendl;
1126 RESTDebug << "-- Error : Gas file was not loaded!" << RESTendl;
1127 return 0;
1128 }
1129
1130 RESTInfo << "Calling Garfield directly. Please be aware that the unit is different "
1131 << "from REST standard unit. E is V/cm. The return is cm^1/2" << RESTendl;
1132
1133 Double_t dl, dt;
1134 fGasMedium->ElectronDiffusion(0., 0, -E, 0, 0, 0, dl, dt);
1135 return dl;
1136#else
1137 std::cout << "This REST is not compiled with garfield, Do not use Longitudinal "
1138 "Diffusion from TRestDetectorGas!"
1139 << std::endl;
1140 std::cout << "Please define the Longitudinal Diffusion in each process!" << std::endl;
1141 return 0;
1142#endif
1143}
1144
1149Double_t TRestDetectorGas::GetTransversalDiffusion(Double_t E) {
1150 RESTDebug << "Entering ... TRestDetectorGas::GetTransversalDiffusion( E=" << E << " )" << RESTendl;
1151
1152 SetPressure(fPressureInAtm);
1153#if defined USE_Garfield
1154 if (fStatus != RESTGAS_GASFILE_LOADED) {
1155 RESTDebug << "-- Error : " << __PRETTY_FUNCTION__ << RESTendl;
1156 RESTDebug << "-- Error : Gas file was not loaded!" << RESTendl;
1157 return 0;
1158 }
1159
1160 RESTInfo << "Calling Garfield directly. Please be aware that the unit is different "
1161 << "from REST standard unit. E is V/cm. The return is cm^1/2" << RESTendl;
1162
1163 Double_t dl, dt;
1164 fGasMedium->ElectronDiffusion(0., 0, -E, 0, 0, 0, dl, dt);
1165 return dt;
1166#else
1167 std::cout << "This REST is not complied with garfield, Do not use Transversal "
1168 "Diffusion from TRestDetectorGas!"
1169 << std::endl;
1170 std::cout << "Please define the Transversal Diffusion in each process!" << std::endl;
1171 return 0;
1172#endif
1173}
1174
1179Double_t TRestDetectorGas::GetTownsendCoefficient(Double_t E) {
1180 RESTDebug << "Entering ... TRestDetectorGas::GetTownsendCoefficient( E=" << E << " )" << RESTendl;
1181
1182 SetPressure(fPressureInAtm);
1183#if defined USE_Garfield
1184 if (fStatus != RESTGAS_GASFILE_LOADED) {
1185 RESTDebug << "-- Error : " << __PRETTY_FUNCTION__ << RESTendl;
1186 RESTDebug << "-- Error : Gas file was not loaded!" << RESTendl;
1187 return 0;
1188 }
1189
1190 RESTInfo << "Calling Garfield directly. Please be aware that the unit is different "
1191 << "from REST standard unit. E is V/cm. The return is V/cm" << RESTendl;
1192
1193 Double_t alpha;
1194 fGasMedium->ElectronTownsend(0., 0, -E, 0, 0, 0, alpha);
1195 return alpha;
1196#else
1197 std::cout << "This REST is not complied with garfield, Do not use Townsend "
1198 "Coefficient from TRestDetectorGas!"
1199 << std::endl;
1200 std::cout << "Please define the Townsend Coefficient in each process!" << std::endl;
1201 return 0;
1202#endif
1203}
1204
1209Double_t TRestDetectorGas::GetAttachmentCoefficient(Double_t E) {
1210 RESTDebug << "Entering ... TRestDetectorGas::GetAttachmentCoefficient( E=" << E << " )" << RESTendl;
1211
1212 SetPressure(fPressureInAtm);
1213#if defined USE_Garfield
1214 if (fStatus != RESTGAS_GASFILE_LOADED) {
1215 RESTDebug << "-- Error : " << __PRETTY_FUNCTION__ << RESTendl;
1216 RESTDebug << "-- Error : Gas file was not loaded!" << RESTendl;
1217 return 0;
1218 }
1219
1220 RESTInfo << "Calling Garfield directly. Please be aware that the unit is different "
1221 << "from REST standard unit. E is V/cm. The return is V/cm" << RESTendl;
1222
1223 Double_t eta;
1224 fGasMedium->ElectronAttachment(0., 0, -E, 0, 0, 0, eta);
1225 return eta;
1226#else
1227 std::cout << "This REST is not complied with garfield, Do not use Attachment "
1228 "Coefficient from TRestDetectorGas!"
1229 << std::endl;
1230 std::cout << "Please define the Attachment Coefficient in each process!" << std::endl;
1231 return 0;
1232#endif
1233}
1234
1238void TRestDetectorGas::PrintGasInfo() {
1239 RESTDebug << "Entering ... TRestDetectorGas::PrintGasInfo( )" << RESTendl;
1240
1242
1243 RESTMetadata << "Status : ";
1244 if (fStatus == RESTGAS_INTITIALIZED) RESTMetadata << "Initialized";
1245 if (fStatus == RESTGAS_CFG_LOADED) RESTMetadata << "Configuration loaded";
1246 if (fStatus == RESTGAS_GASFILE_LOADED) RESTMetadata << "Gasfile loaded";
1247 if (fStatus == RESTGAS_ERROR) RESTMetadata << "Error";
1248 RESTMetadata << RESTendl;
1249
1250 RESTMetadata << "Gas filename : " << TRestTools::GetPureFileName((string)fGasFilename) << RESTendl;
1251 RESTMetadata << "Pressure : " << fPressureInAtm << " atm" << RESTendl;
1252 RESTMetadata << "Temperature : " << fTemperatureInK << " K" << RESTendl;
1253 RESTMetadata << "Electric Field : " << fElectricField * units("V/cm") << " V/cm " << RESTendl;
1254 RESTMetadata << "W-value : " << fW << " eV" << RESTendl;
1255 RESTMetadata << "Drift velocity : " << GetDriftVelocity(fElectricField * units("V/cm")) / units("cm/us")
1256 << " mm/us " << RESTendl;
1257 RESTMetadata << "Max. Electron energy : " << fMaxElectronEnergy << " eV" << RESTendl;
1258 RESTMetadata << "Field grid nodes : " << fEnodes << RESTendl;
1259 RESTMetadata << "Efield range : ( " << fEmin << " , " << fEmax << " ) V/cm " << RESTendl;
1260 RESTMetadata << "Number of Gases : " << fNofGases << RESTendl;
1261 for (int i = 0; i < fNofGases; i++)
1262 RESTMetadata << "Gas id : " << i << ", Name : " << fGasComponentName[i]
1263 << ", Fraction : " << fGasComponentFraction[i] << RESTendl;
1264 RESTMetadata << "******************************************" << RESTendl;
1265 RESTMetadata << RESTendl;
1266 RESTMetadata << RESTendl;
1267}
1268
1269Int_t TRestDetectorGas::Write(const char* name, Int_t option, Int_t bufsize) {
1270 RESTDebug << "Entering ... TRestDetectorGas::Write( name=" << name << " option=" << option
1271 << " bufsize=" << bufsize << " )" << RESTendl;
1272
1273 if (fGasFileContent == "" && GasFileLoaded()) {
1274 ifstream infile;
1275 infile.open(fGasFilename);
1276 if (!infile) {
1277 cout << "TRestDetectorGas: error reading gas file, gas file content won't be "
1278 "saved!"
1279 << endl;
1280 } else {
1281 string str;
1282 while (getline(infile, str)) {
1283 fGasFileContent += str + "\n";
1284 }
1285 }
1286 }
1287 return TRestMetadata::Write();
1288}
virtual DBEntry query_data(DBEntry info)
void InitFromConfigFile() override
To make settings from rml file. This method must be implemented in the derived class.
virtual void SetElectricField(double value)
Sets the electric field. Must be given in V/mm.
Double_t GetGasFanoFactor() const
Returns the gas fano factor.
void InitFromConfigFile() override
To make settings from rml file. This method must be implemented in the derived class.
void InitFromRootFile() override
Method called after the object is retrieved from root file.
TString GetGasComponentName(Int_t n)
Returns the gas component n.
bool GasFileLoaded() const
Returns true if the gas file has been properly loaded. False otherwise.
Double_t GetTransversalDiffusion() override
Returns the transversal diffusion in (cm)^1/2.
TRestDetectorGas()
TRestDetectorGas default constructor.
Double_t GetDriftVelocity() override
Returns the drift velocity in mm/us.
Int_t Write(const char *name=nullptr, Int_t option=0, Int_t bufsize=0) override
overwriting the write() method with fStore considered
void Initialize() override
Making default settings.
Double_t GetGasWorkFunction() const
Returns the gas work function in eV.
Double_t GetGasComponentFraction(Int_t n)
Returns the gas fraction in volume for component n.
Double_t GetLongitudinalDiffusion() override
Returns the longitudinal diffusion in (cm)^1/2.
virtual void PrintMetadata()
Implemented it in the derived metadata class to print out specific metadata information.
endl_t RESTendl
Termination flag object for TRestStringOutput.
TiXmlElement * GetElement(std::string eleDeclare, TiXmlElement *e=nullptr)
Get an xml element from a given parent element, according to its declaration.
Int_t LoadConfigFromFile(const std::string &configFilename, const std::string &sectionName="")
Give the file name, find out the corresponding section. Then call the main starter.
TRestStringOutput::REST_Verbose_Level fVerboseLevel
Verbose level used to print debug info.
TiXmlElement * StringToElement(std::string definition)
Parsing a string into TiXmlElement object.
void SetLibraryVersion(TString version)
Set the library version of this metadata class.
TRestStringOutput::REST_Verbose_Level GetVerboseLevel()
returns the verboselevel in type of REST_Verbose_Level enumerator
std::string GetFieldValue(std::string parName, TiXmlElement *e)
Returns the field value of an xml element which has the specified name.
virtual void InitFromRootFile()
Method called after the object is retrieved from root file.
Int_t LoadConfigFromElement(TiXmlElement *eSectional, TiXmlElement *eGlobal, std::map< std::string, std::string > envs={})
Main starter method.
void SetSectionName(std::string sName)
set the section name, clear the section content
std::string SearchFile(std::string filename)
Search files in current directory and directories specified in "searchPath" section.
void WriteConfigBuffer(std::string fName)
Writes the config buffer to a file in append mode.
std::string fConfigFileName
Full name of the rml file.
virtual Int_t Write(const char *name=nullptr, Int_t option=0, Int_t bufsize=0)
overwriting the write() method with fStore considered
TiXmlElement * GetNextElement(TiXmlElement *e)
Get the next sibling xml element of this element, with same eleDeclare.
std::string GetParameter(std::string parName, TiXmlElement *e, TString defaultValue=PARAMETER_NOT_FOUND_STR)
Returns the value for the parameter named parName in the given section.
@ REST_Info
+show most of the information for each steps
static std::string GetPureFileName(const std::string &fullPathFileName)
Removes all directories in the full path filename description given in the argument.
Definition: TRestTools.cxx:863
static bool fileExists(const std::string &filename)
Returns true if the file (or directory) with path filename exists.
Definition: TRestTools.cxx:728
static std::string DownloadRemoteFile(const std::string &remoteFile, bool pidPrefix=false)
download the remote file automatically, returns the downloaded file name.
static bool isPathWritable(const std::string &path)
Returns true if the path given by argument is writable.
Definition: TRestTools.cxx:778
static int UploadToServer(std::string localFile, std::string remoteFile, std::string methodUrl="")
Int_t GetChar(std::string hint="Press a KEY to continue ...")
Helps to pause the program, printing a message before pausing.
std::vector< std::string > Split(std::string in, std::string separator, bool allowBlankString=false, bool removeWhiteSpaces=false, int startPos=-1)
Split the input string according to the given separator. Returning a vector of fragments.
Double_t StringToDouble(std::string in)
Gets a double from a string.
std::string ToUpper(std::string in)
Convert string to its upper case. Alternative of TString::ToUpper.
Int_t StringToInteger(std::string in)
Gets an integer from a string.