Added a first version of code for a library to fit data to any p(B) -- at the moment it is untested and subject to change...
This commit is contained in:
parent
badcc28f34
commit
728ec00816
39
src/external/TFitPofB-lib/classes/Makefile.TFitPofB
vendored
Normal file
39
src/external/TFitPofB-lib/classes/Makefile.TFitPofB
vendored
Normal file
@ -0,0 +1,39 @@
|
||||
CXX = g++
|
||||
CXXFLAGS = -g -Wall -Wno-trigraphs -fPIC
|
||||
LOCALINCLUDE = ../include/
|
||||
INCLUDES = -I$(LOCALINCLUDE)
|
||||
LD = g++
|
||||
LDFLAGS = -g
|
||||
SOFLAGS = -O -shared
|
||||
|
||||
# some definitions: headers (used to generate *Dict* stuff), sources, objects,...
|
||||
OBJS =
|
||||
OBJS += TTrimSPDataHandler.o
|
||||
OBJS += TBofZCalc.o
|
||||
OBJS += TPofBCalc.o
|
||||
OBJS += TPofTCalc.o
|
||||
OBJS += TFitPofB.o
|
||||
|
||||
SHLIB = libTFitPofB.so
|
||||
|
||||
# make the shared lib:
|
||||
#
|
||||
all: $(SHLIB)
|
||||
|
||||
$(SHLIB): $(OBJS)
|
||||
@echo "---> Building shared library $(SHLIB) ..."
|
||||
/bin/rm -f $(SHLIB)
|
||||
$(LD) $(OBJS) $(EXTOBJS) $(SOFLAGS) -o $(SHLIB)
|
||||
@echo "done"
|
||||
|
||||
# clean up: remove all object file (and core files)
|
||||
# semicolon needed to tell make there is no source
|
||||
# for this target!
|
||||
#
|
||||
clean:; @rm -f $(OBJS) *Dict* core*
|
||||
@echo "---> removing $(OBJS)"
|
||||
|
||||
#
|
||||
$(OBJS): %.o: %.cpp
|
||||
$(CXX) $(INCLUDES) $(CXXFLAGS) -c $<
|
||||
|
129
src/external/TFitPofB-lib/classes/TBofZCalc.cpp
vendored
Normal file
129
src/external/TFitPofB-lib/classes/TBofZCalc.cpp
vendored
Normal file
@ -0,0 +1,129 @@
|
||||
/***************************************************************************
|
||||
|
||||
TBofZCalc.cpp
|
||||
|
||||
Author: Bastian M. Wojek
|
||||
e-mail: bastian.wojek@psi.ch
|
||||
|
||||
2008/05/23
|
||||
|
||||
***************************************************************************/
|
||||
|
||||
#include "TBofZCalc.h"
|
||||
#include <cmath>
|
||||
#include <iostream>
|
||||
#include <algorithm>
|
||||
|
||||
|
||||
double TBofZCalc::GetBmin() const {
|
||||
vector<double>::const_iterator iter;
|
||||
iter = min_element(fBZ.begin(),fBZ.end());
|
||||
|
||||
return *iter;
|
||||
}
|
||||
|
||||
double TBofZCalc::GetBmax() const {
|
||||
vector<double>::const_iterator iter;
|
||||
iter = max_element(fBZ.begin(),fBZ.end());
|
||||
|
||||
return *iter;
|
||||
}
|
||||
|
||||
double TBofZCalc::GetBofZ(double zz) const {
|
||||
|
||||
bool found = false;
|
||||
unsigned int i;
|
||||
for (i=0; i<fZ.size(); i++) {
|
||||
if (fZ[i] > zz) {
|
||||
found = true;
|
||||
break;
|
||||
}
|
||||
}
|
||||
|
||||
if (!found || i == 0) {
|
||||
cout << "B(z) cannot be calculated for z = " << zz << " !" << endl;
|
||||
cout << "Check your theory function!" << endl;
|
||||
return -1.0;
|
||||
}
|
||||
|
||||
return fBZ[i-1]+(fBZ[i]-fBZ[i-1])*(zz-fZ[i-1])/(fZ[i]-fZ[i-1]);
|
||||
|
||||
}
|
||||
|
||||
TLondon1D_1L::TLondon1D_1L(const vector<double> ¶m) {
|
||||
|
||||
// Parameters for 1D-London screening in a thin superconducting film
|
||||
// Bext[G], deadlayer[nm], thickness[nm], lambda[nm]
|
||||
|
||||
unsigned int n(2000); // number of steps for the calculation
|
||||
|
||||
double N(cosh(param[2]/2.0/param[3]));
|
||||
|
||||
fDZ = param[2]/double(n);
|
||||
double ZZ, BBz;
|
||||
|
||||
for (unsigned int j(0); j<n; j++) {
|
||||
ZZ = param[1] + (double)j*fDZ;
|
||||
fZ.push_back(ZZ);
|
||||
BBz = param[0]*cosh((param[2]/2.0-(ZZ-param[1]))/param[3])/N;
|
||||
fBZ.push_back(BBz);
|
||||
}
|
||||
|
||||
}
|
||||
|
||||
TLondon1D_2L::TLondon1D_2L(const vector<double> ¶m) {
|
||||
|
||||
// Parameters for 1D-London screening in a thin superconducting film, two layers, two lambda
|
||||
// Bext[G], deadlayer[nm], thickness1[nm], thickness2[nm], lambda1[nm], lambda2[nm]
|
||||
|
||||
unsigned int n(2000); // number of steps for the calculation
|
||||
|
||||
double N1(param[5]*cosh(param[3]/param[5])*sinh(param[2]/param[4]) + param[4]*cosh(param[2]/param[4])*sinh(param[3]/param[5]));
|
||||
double N2(4.0*N1);
|
||||
|
||||
fDZ = (param[2]+param[3])/double(n);
|
||||
double ZZ, BBz;
|
||||
|
||||
for (unsigned int j(0); j<n; j++) {
|
||||
ZZ = param[1] + (double)j*fDZ;
|
||||
fZ.push_back(ZZ);
|
||||
if (ZZ < param[1]+param[2]) {
|
||||
BBz = param[0]*(param[4]*cosh((param[2]+param[1]-ZZ)/param[4])*sinh(param[3]/param[5]) - param[5]*sinh((param[1]-ZZ)/param[4]) + param[5]*cosh(param[3]/param[5])*sinh((param[2]+param[1]-ZZ)/param[4]))/N1;
|
||||
} else {
|
||||
BBz = param[0]*(exp(-param[2]/param[4]-(param[2]+param[3]+param[1]+ZZ)/param[5]) * (exp((param[3]+2.0*param[1])/param[5])*(exp(2.0*param[2]*(param[4]+param[5])/param[4]/param[5]) * (param[5]-param[4]) - exp(2.0*param[2]/param[5]) * (param[4]-2.0*param[4]*exp(param[2]/param[4]+param[3]/param[5])+param[5])) + exp(2.0*ZZ/param[5])*((param[4]-param[5]+(param[4]+param[5]) * exp(2.0*param[2]/param[4]))*exp(param[3]/param[5])-2.0*param[4]*exp(param[2]/param[4]))))/N2;
|
||||
}
|
||||
fBZ.push_back(BBz);
|
||||
}
|
||||
|
||||
}
|
||||
|
||||
TLondon1D_3L::TLondon1D_3L(const vector<double> ¶m) {
|
||||
|
||||
// Parameters for 1D-London screening in a thin superconducting film, three layers, two lambdas (top and bottom layer: lambda1)
|
||||
// Bext[G], deadlayer[nm], thickness1[nm], thickness2[nm], thickness3[nm], lambda1[nm], lambda2[nm]
|
||||
|
||||
unsigned int n(2000); // number of steps for the calculation
|
||||
|
||||
double N1(8.0*(param[5]*param[6]*cosh(param[3]/param[6])*sinh((param[2]+param[4])/param[5]) + ((param[5]*param[5]*cosh(param[2]/param[5])*cosh(param[4]/param[5])) + (param[6]*param[6]*sinh(param[2]/param[5])*sinh(param[4]/param[5])))*sinh(param[3]/param[6])));
|
||||
|
||||
double N2(2.0*param[5]*param[6]*cosh(param[3]/param[6])*sinh((param[2]+param[4])/param[5]) + 2.0*(param[5]*param[5]*cosh(param[2]/param[5])*cosh(param[4]/param[5]) + param[6]*param[6]*sinh(param[2]/param[5])*sinh(param[4]/param[5]))*sinh(param[3]/param[6]));
|
||||
|
||||
double N3(8.0*(param[5]*param[6]*cosh(param[3]/param[6])*sinh((param[2]+param[4])/param[5]) + (param[5]*param[5]*cosh(param[2]/param[5])*cosh(param[4]/param[5]) + param[6]*param[6]*sinh(param[2]/param[5])*sinh(param[4]/param[5]))*sinh(param[3]/param[6])));
|
||||
|
||||
fDZ = (param[2]+param[3]+param[4])/double(n);
|
||||
double ZZ, BBz;
|
||||
|
||||
for (unsigned int j(0); j<n; j++) {
|
||||
ZZ = param[1] + (double)j*fDZ;
|
||||
fZ.push_back(ZZ);
|
||||
if (ZZ < param[1]+param[2]) {
|
||||
BBz = (param[0]*exp(-1.0*((param[2]+param[4]+param[1]+ZZ)/param[5]+(param[3]/param[6]))) *
|
||||
(-4.0*exp((param[2]+param[4])/param[5]+(param[3]/param[6]))*(exp((2.0*param[1])/param[5]) - exp((2.0*ZZ)/param[5]))*param[5]*param[6]+exp((2.0*param[3])/param[6])*(param[5]-param[6] + (exp((2.0*param[4])/param[5])*(param[5]+param[6])))*((exp((2.0*ZZ)/param[5])*(param[5]-param[6])) + (exp((2.0*(param[2]+param[1]))/param[5])*(param[5]+param[6]))) - 4.0*((param[5]*cosh(param[4]/param[5]))-(param[6]*sinh(param[4]/param[5])))*((param[5]*cosh((param[2]+param[1]-ZZ)/param[5])) - (param[6]*sinh((param[2]+param[1]-ZZ)/param[5])))*(cosh((param[2]+param[4]+param[1]+ZZ)/param[5]) + sinh((param[2]+param[4]+param[1]+ZZ)/param[5]))))/N1;
|
||||
} else if (ZZ < param[1]+param[2]+param[3]) {
|
||||
BBz = (param[0]*param[5]*(2.0*param[6]*cosh((param[2]+param[3]+param[1]-ZZ)/param[6])*sinh(param[4]/param[5]) + (param[5]+param[6])*sinh(param[2]/param[5]-(param[2]+param[1]-ZZ)/param[6]) - (param[5]-param[6])*sinh(param[2]/param[5]+(param[2]+param[1]-ZZ)/param[6]) + 2.0*param[5]*cosh(param[4]/param[5])*sinh((param[2]+param[3]+param[1]-ZZ)/param[6])))/N2;
|
||||
} else {
|
||||
BBz = (param[0]*exp(-((2.0*param[2]+param[3]+param[4]+param[1]+ZZ)/param[5])-param[3]/param[6]) * (4.0*exp(param[2]/param[5]+param[3]/param[6])*param[5]*param[6]*((-1.0)*exp(2.0*ZZ/param[5]) + cosh(2.0*(param[2]+param[3]+param[4]+param[1])/param[5])+sinh(2.0*(param[2]+param[3]+param[4]+param[1])/param[5])) + 4.0*exp(((2.0*param[2]+param[3]+param[4]+param[1]+ZZ)/param[5]) + ((2.0*param[3])/param[6]))*(param[5]*cosh(param[2]/param[5]) + param[6]*sinh(param[2]/param[5]))*(param[5]*cosh((param[2]+param[3]+param[1]-ZZ)/param[5]) - param[6]*sinh((param[2]+param[3]+param[1]-ZZ)/param[5])) - 4.0*(param[5]*cosh(param[2]/param[5])-param[6]*sinh(param[2]/param[5])) * (param[5]*cosh((param[2]+param[3]+param[1]-ZZ)/param[5]) + param[6]*sinh((param[2]+param[3]+param[1]-ZZ)/param[5]))*(cosh((2.0*param[2]+param[3]+param[4]+param[1]+ZZ)/param[5]) + sinh((2.0*param[2]+param[3]+param[4]+param[1]+ZZ)/param[5]))))/N3;
|
||||
}
|
||||
fBZ.push_back(BBz);
|
||||
}
|
||||
}
|
121
src/external/TFitPofB-lib/classes/TFitPofB.cpp
vendored
Normal file
121
src/external/TFitPofB-lib/classes/TFitPofB.cpp
vendored
Normal file
@ -0,0 +1,121 @@
|
||||
/***************************************************************************
|
||||
|
||||
TFitPofB.cpp
|
||||
|
||||
Author: Bastian M. Wojek
|
||||
e-mail: bastian.wojek@psi.ch
|
||||
|
||||
2008/05/23
|
||||
|
||||
***************************************************************************/
|
||||
|
||||
#include "TFitPofB.h"
|
||||
#include <iostream>
|
||||
using namespace std;
|
||||
|
||||
TFitPofB::TFitPofB(const vector<unsigned int> &parNo, const vector<double> &par) : fCalcNeeded(true) {
|
||||
|
||||
for(unsigned int i(0); i<parNo.size(); i++) {
|
||||
fPar.push_back(par[parNo[i]-1]);
|
||||
}
|
||||
|
||||
string rge_path("/home/l_wojek/nt/wojek/g/Bastian/ImplantationDepth/YBCO_PBCO-");
|
||||
string energy_arr[] = {"02_1", "02_5", "03_5", "05_0", "07_5", "10_0", "12_5", "15_0", "17_5", "19_0", "20_0", "22_5", "25_0"};
|
||||
|
||||
vector<string> energy_vec(energy_arr, energy_arr+(sizeof(energy_arr)/sizeof(energy_arr[0])));
|
||||
|
||||
vector<double> par_for_PofT;
|
||||
|
||||
for (unsigned int i(1); i<4; i++)
|
||||
par_for_PofT.push_back(fPar[i]);
|
||||
|
||||
TTrimSPData *x = new TTrimSPData(rge_path, energy_vec);
|
||||
fImpProfile = x;
|
||||
x = 0;
|
||||
delete x;
|
||||
|
||||
TPofTCalc *y = new TPofTCalc(par_for_PofT);
|
||||
fPofT = y;
|
||||
y = 0;
|
||||
delete y;
|
||||
|
||||
}
|
||||
|
||||
TFitPofB::~TFitPofB() {
|
||||
fPar.clear();
|
||||
fImpProfile = 0;
|
||||
delete fImpProfile;
|
||||
fPofT = 0;
|
||||
delete fPofT;
|
||||
}
|
||||
|
||||
double TFitPofB::Eval(double t, const vector<double> &par) const {
|
||||
|
||||
// check if any parameter has changed
|
||||
|
||||
bool par_changed(false);
|
||||
|
||||
for (unsigned int i(0); i<fPar.size(); i++) {
|
||||
if( fPar[i]-par[i] ) {
|
||||
fPar[i] = par[i];
|
||||
par_changed = true;
|
||||
if(i == 2 || i == 3) {
|
||||
cout << "You are varying dt or dB! These parameters have to be fixed! Quitting..." << endl;
|
||||
exit(0);
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
if (par_changed)
|
||||
fCalcNeeded = true;
|
||||
|
||||
// if parameters have changed, recalculate B(z), P(B) and P(t)
|
||||
|
||||
if (fCalcNeeded) {
|
||||
|
||||
vector<double> par_for_BofZ;
|
||||
vector<double> par_for_PofB;
|
||||
vector<double> par_for_PofT;
|
||||
|
||||
for (unsigned int i(5); i<par.size(); i++)
|
||||
par_for_BofZ.push_back(par[i]);
|
||||
|
||||
for (unsigned int i(2); i<5; i++)
|
||||
par_for_PofB.push_back(par[i]);
|
||||
|
||||
for (unsigned int i(1); i<4; i++)
|
||||
par_for_PofT.push_back(par[i]);
|
||||
|
||||
switch (int(par[0])) {
|
||||
case 1:
|
||||
{
|
||||
TLondon1D_1L BofZ1(par_for_BofZ);
|
||||
TPofBCalc PofB1(BofZ1, *fImpProfile, par_for_PofB);
|
||||
fPofT->DoFFT(PofB1, par_for_PofT);
|
||||
}
|
||||
break;
|
||||
case 2:
|
||||
{
|
||||
TLondon1D_2L BofZ2(par_for_BofZ);
|
||||
TPofBCalc PofB2(BofZ2, *fImpProfile, par_for_PofB);
|
||||
fPofT->DoFFT(PofB2, par_for_PofT);
|
||||
}
|
||||
break;
|
||||
case 3:
|
||||
{
|
||||
TLondon1D_3L BofZ3(par_for_BofZ);
|
||||
TPofBCalc PofB3(BofZ3, *fImpProfile, par_for_PofB);
|
||||
fPofT->DoFFT(PofB3, par_for_PofT);
|
||||
}
|
||||
break;
|
||||
default:
|
||||
cout << "The user function you specified, does not exist. Quitting..." << endl;
|
||||
exit(0);
|
||||
}
|
||||
|
||||
}
|
||||
|
||||
return fPofT->Eval(t);
|
||||
|
||||
}
|
||||
|
117
src/external/TFitPofB-lib/classes/TPofBCalc.cpp
vendored
Normal file
117
src/external/TFitPofB-lib/classes/TPofBCalc.cpp
vendored
Normal file
@ -0,0 +1,117 @@
|
||||
/***************************************************************************
|
||||
|
||||
TPofBCalc.cpp
|
||||
|
||||
Author: Bastian M. Wojek
|
||||
e-mail: bastian.wojek@psi.ch
|
||||
|
||||
2008/05/23
|
||||
|
||||
***************************************************************************/
|
||||
|
||||
#include "TPofBCalc.h"
|
||||
#include <cmath>
|
||||
#include <iostream>
|
||||
|
||||
TPofBCalc::TPofBCalc( const TBofZCalc &BofZ, const TTrimSPData &dataTrimSP, const vector<double> ¶ ) {
|
||||
|
||||
// Parameters: dt[us], dB[G], Energy[keV]
|
||||
|
||||
double BB, BBnext;
|
||||
double zm, zp, zNextm, zNextp, dz;
|
||||
|
||||
// fill not used Bs before Bmin with 0.0
|
||||
|
||||
for ( BB = 0.0 ; BB < BofZ.GetBmin() ; BB += para[1] ) {
|
||||
fB.push_back(BB);
|
||||
fPB.push_back(0.0);
|
||||
}
|
||||
|
||||
unsigned int firstZerosEnd(fB.size());
|
||||
|
||||
// calculate p(B) from B(z)
|
||||
|
||||
vector<double> bofzZ(BofZ.DataZ());
|
||||
vector<double> bofzBZ(BofZ.DataBZ());
|
||||
double ddZ(BofZ.GetdZ());
|
||||
double Bmax(BofZ.GetBmax());
|
||||
|
||||
double nn;
|
||||
for ( ; BB <= Bmax ; BB += para[1]) {
|
||||
BBnext = BB + para[1];
|
||||
fB.push_back(BB);
|
||||
fPB.push_back(0.0);
|
||||
|
||||
for ( unsigned int j(0); j < bofzZ.size() - 1; j++ ) {
|
||||
if ( bofzBZ[j] >= BB && bofzBZ[j+1] <= BB ) {
|
||||
zm = (BB-bofzBZ[j])*ddZ/(bofzBZ[j+1]-bofzBZ[j]) + bofzZ[j];
|
||||
|
||||
for (unsigned int k(0); k < j; k++) {
|
||||
if ( ( bofzBZ[j-k] <= BBnext && bofzBZ[j-k-1] >= BBnext ) ) {
|
||||
// cout << "1 " << j << " " << k << endl;
|
||||
zNextm = (BBnext-bofzBZ[j-k-1])*ddZ/(bofzBZ[j-k]-bofzBZ[j-k-1]) + bofzZ[j-k-1];
|
||||
break;
|
||||
}
|
||||
}
|
||||
|
||||
dz = zNextm-zm;
|
||||
nn = dataTrimSP.GetNofZ(zm, para[2]);
|
||||
if (nn != -1.0) {
|
||||
// cout << "zNext = " << zNextm << ", zm = " << zm << ", dz = " << dz << endl;
|
||||
*(fPB.end()-1) += nn*fabs(dz/para[1]);
|
||||
}
|
||||
|
||||
} else if (bofzBZ[j] <= BB && bofzBZ[j+1] >= BB ) {
|
||||
zp = (BB-bofzBZ[j])*ddZ/(bofzBZ[j+1]-bofzBZ[j]) + bofzZ[j];
|
||||
|
||||
for (unsigned int k(0); k < bofzZ.size() - j - 1; k++) {
|
||||
if ( ( bofzBZ[j+k] <= BBnext && bofzBZ[j+k+1] >= BBnext ) ) {
|
||||
// cout << "2 " << j << " " << k << endl;
|
||||
zNextp = (BBnext-bofzBZ[j+k])*ddZ/(bofzBZ[j+k+1]-bofzBZ[j+k]) + bofzZ[j+k];
|
||||
break;
|
||||
}
|
||||
}
|
||||
|
||||
dz = zNextp-zp;
|
||||
nn = dataTrimSP.GetNofZ(zp, para[2]);
|
||||
if (nn != -1.0) {
|
||||
// cout << "zNext = " << zNextp << ", zp = " << zp << ", dz = " << dz << endl;
|
||||
*(fPB.end()-1) += nn*fabs(dz/para[1]);
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
}
|
||||
|
||||
unsigned int lastZerosStart(fB.size());
|
||||
|
||||
// fill not used Bs after Bext with 0.0
|
||||
|
||||
double BmaxFFT(1.0/gBar/para[0]);
|
||||
|
||||
//// cout << "N = " << int(BmaxFFT/para[1]+1.0) << endl;
|
||||
|
||||
for ( ; BB <= BmaxFFT ; BB += para[1] ) {
|
||||
fB.push_back(BB);
|
||||
fPB.push_back(0.0);
|
||||
}
|
||||
|
||||
// make sure that we have an even number of elements in p(B) for FFT
|
||||
|
||||
if (fB.size() % 2) {
|
||||
fB.push_back(BB);
|
||||
fPB.push_back(0.0);
|
||||
}
|
||||
|
||||
cout << "size of B = " << fB.size() << ", size of p(B) = " << fPB.size() << endl;
|
||||
|
||||
// normalize pB
|
||||
|
||||
double pBsum = 0.0;
|
||||
for (unsigned int i(firstZerosEnd); i<lastZerosStart; i++)
|
||||
pBsum += fPB[i];
|
||||
pBsum *= para[1];
|
||||
for (unsigned int i(firstZerosEnd); i<lastZerosStart; i++)
|
||||
fPB[i] /= pBsum;
|
||||
|
||||
}
|
113
src/external/TFitPofB-lib/classes/TPofTCalc.cpp
vendored
Normal file
113
src/external/TFitPofB-lib/classes/TPofTCalc.cpp
vendored
Normal file
@ -0,0 +1,113 @@
|
||||
/***************************************************************************
|
||||
|
||||
TPofTCalc.cpp
|
||||
|
||||
Author: Bastian M. Wojek
|
||||
e-mail: bastian.wojek@psi.ch
|
||||
|
||||
2008/05/23
|
||||
|
||||
***************************************************************************/
|
||||
|
||||
#include "TPofTCalc.h"
|
||||
#include "fftw3.h"
|
||||
#include <cmath>
|
||||
#include <cstdio>
|
||||
#include <iostream>
|
||||
|
||||
TPofTCalc::TPofTCalc (const vector<double> &par) {
|
||||
|
||||
fNFFT = ( int(1.0/gBar/par[1]/par[2]+1.0) % 2 ) ? int(1.0/gBar/par[1]/par[2]+2.0) : int(1.0/gBar/par[1]/par[2]+1.0);
|
||||
fTBin = 1.0/gBar/double(fNFFT-1)/par[2];
|
||||
// tBin = 1.0/gBar/(*(PofB.DataB().end()-1));
|
||||
// nFFT = PofB.DataB().size();
|
||||
|
||||
fFFTin = (double *)malloc(sizeof(double) * fNFFT);
|
||||
fFFTout = (fftw_complex *)fftw_malloc(sizeof(fftw_complex) * (fNFFT/2+1));
|
||||
|
||||
cout << "Check for the FFT plan..." << endl;
|
||||
|
||||
// Load wisdom from file
|
||||
|
||||
int wisdomLoaded(0);
|
||||
|
||||
FILE *wordsOfWisdomR;
|
||||
wordsOfWisdomR = fopen("WordsOfWisdom.dat", "r");
|
||||
if (wordsOfWisdomR == NULL) {
|
||||
cout << "Couldn't open wisdom file ..." << endl;
|
||||
} else {
|
||||
wisdomLoaded = fftw_import_wisdom_from_file(wordsOfWisdomR);
|
||||
fclose(wordsOfWisdomR);
|
||||
}
|
||||
|
||||
if (!wisdomLoaded) {
|
||||
cout << "No wisdom is imported..." << endl;
|
||||
}
|
||||
|
||||
fFFTplan = fftw_plan_dft_r2c_1d(fNFFT, fFFTin, fFFTout, FFTW_EXHAUSTIVE);
|
||||
|
||||
}
|
||||
|
||||
void TPofTCalc::DoFFT(const TPofBCalc &PofB, const vector<double> &par) {
|
||||
|
||||
vector<double> pB(PofB.DataPB());
|
||||
|
||||
for (unsigned int i(0); i<fNFFT; i++) {
|
||||
fFFTin[i] = pB[i];
|
||||
}
|
||||
for (unsigned int i(0); i<fNFFT/2+1; i++) {
|
||||
fFFTout[i][0] = 0.0;
|
||||
fFFTout[i][1] = 0.0;
|
||||
}
|
||||
|
||||
cout << "perform the Fourier transform..." << endl;
|
||||
|
||||
fftw_execute(fFFTplan);
|
||||
|
||||
// Calculate polarisation
|
||||
|
||||
double sinph(sin(par[0]*PI/180.0)), cosph(cos(par[0]*PI/180.0));
|
||||
|
||||
double polTemp(0.0);
|
||||
|
||||
for (unsigned int i(0); i<fNFFT/2+1; i++){
|
||||
fT.push_back(double(i)*fTBin);
|
||||
polTemp = cosph*fFFTout[i][0]*par[2] + sinph*fFFTout[i][1]*par[2];
|
||||
fPT.push_back(polTemp);
|
||||
}
|
||||
|
||||
}
|
||||
|
||||
double TPofTCalc::Eval(double t) const {
|
||||
for (unsigned int i(0); i<fT.size()-1; i++) {
|
||||
if (t < fT[i+1]) {
|
||||
return fPT[i]+(fPT[i+1]-fPT[i])/(fT[i+1]-fT[i])*(t-fT[i]);
|
||||
}
|
||||
}
|
||||
|
||||
cout << "No data for the time " << t << " us available! Quitting..." << endl;
|
||||
exit(0);
|
||||
}
|
||||
|
||||
TPofTCalc::~TPofTCalc() {
|
||||
// export wisdom
|
||||
|
||||
FILE *wordsOfWisdomW;
|
||||
wordsOfWisdomW = fopen("WordsOfWisdom.dat", "w");
|
||||
if (wordsOfWisdomW == NULL) {
|
||||
cout << "couldn't open file ... No wisdom is exported..." << endl;
|
||||
}
|
||||
|
||||
fftw_export_wisdom_to_file(wordsOfWisdomW);
|
||||
|
||||
fclose(wordsOfWisdomW);
|
||||
|
||||
fftw_destroy_plan(fFFTplan);
|
||||
free(fFFTin);
|
||||
fftw_free(fFFTout);
|
||||
fT.clear();
|
||||
fPT.clear();
|
||||
|
||||
}
|
||||
|
||||
|
127
src/external/TFitPofB-lib/classes/TTrimSPDataHandler.cpp
vendored
Normal file
127
src/external/TFitPofB-lib/classes/TTrimSPDataHandler.cpp
vendored
Normal file
@ -0,0 +1,127 @@
|
||||
/***************************************************************************
|
||||
|
||||
TTrimSPDataHandler.cpp
|
||||
|
||||
Author: Bastian M. Wojek
|
||||
e-mail: bastian.wojek@psi.ch
|
||||
|
||||
2008/05/23
|
||||
|
||||
***************************************************************************/
|
||||
|
||||
#include "TTrimSPDataHandler.h"
|
||||
#include <iostream>
|
||||
#include <fstream>
|
||||
#include <string>
|
||||
#include <cmath>
|
||||
|
||||
using namespace std;
|
||||
|
||||
// TrimSPData constructor -- reading all available trim.SP-rge-files into std::vectors
|
||||
|
||||
TTrimSPData::TTrimSPData(const string &path, vector<string> &energyVec) {
|
||||
|
||||
//string energyArr[] = {"02_1", "02_5", "03_5", "05_0", "07_5", "10_0", "12_5", "15_0", "17_5", "19_0", "20_0", "22_5", "25_0"};
|
||||
|
||||
double zz(0.0), nzz(0.0);
|
||||
vector<double> vzz, vnzz;
|
||||
string word, energyStr;
|
||||
|
||||
for(unsigned int i(0); i<energyVec.size(); i++){
|
||||
|
||||
// energyStr = "/home/l_wojek/nt/wojek/g/Bastian/ImplantationDepth/YBCO_PBCO-" + energyVec[i] + ".rge";
|
||||
|
||||
energyStr = path + energyVec[i] + ".rge";
|
||||
|
||||
fEnergy.push_back(atof(energyVec[i].replace(2,1,".").c_str()));
|
||||
|
||||
ifstream *rgeFile = new ifstream(energyStr.c_str());
|
||||
if(! *rgeFile) {
|
||||
cout << "rge-file not found! Exit now." << endl;
|
||||
exit(-1);
|
||||
} else {
|
||||
|
||||
while(*rgeFile >> word)
|
||||
if(word == "PARTICLES") break;
|
||||
|
||||
while(!rgeFile->eof()) {
|
||||
*rgeFile >> zz >> nzz;
|
||||
vzz.push_back(zz);
|
||||
vnzz.push_back(nzz);
|
||||
}
|
||||
|
||||
fDataZ.push_back(vzz);
|
||||
fDataNZ.push_back(vnzz);
|
||||
|
||||
rgeFile->close();
|
||||
delete rgeFile;
|
||||
rgeFile = 0;
|
||||
vzz.clear();
|
||||
vnzz.clear();
|
||||
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
// vector<double> DataZ(double) -- returning z-vector calculated by trim.SP for energy given energy
|
||||
|
||||
vector<double> TTrimSPData::DataZ(double e) const {
|
||||
|
||||
for(unsigned int i(0); i<fEnergy.size(); i++) {
|
||||
// cout << tEnergy[i] << " " << e << " " << tEnergy[i] - e << endl;
|
||||
if(!(fEnergy[i] - e)) {
|
||||
return fDataZ[i];
|
||||
}
|
||||
}
|
||||
// default
|
||||
cout << "No implantation profile available for the specified energy... You got back the first one." << endl;
|
||||
return fDataZ[0];
|
||||
|
||||
}
|
||||
|
||||
// vector<double> DataNZ(double) -- returning n(z)-vector calculated by trim.SP for given energy
|
||||
|
||||
vector<double> TTrimSPData::DataNZ(double e) const {
|
||||
|
||||
for(unsigned int i(0); i<fEnergy.size(); i++) {
|
||||
if(!(fEnergy[i] - e)) {
|
||||
return fDataNZ[i];
|
||||
}
|
||||
}
|
||||
// default
|
||||
cout << "No implantation profile available for the specified energy... You got back the first one." << endl;
|
||||
return fDataNZ[0];
|
||||
|
||||
}
|
||||
|
||||
// double GetNofZ(double, double) -- return NofZ for given z in nanometers
|
||||
|
||||
double TTrimSPData::GetNofZ(double zz, double e) const {
|
||||
|
||||
vector<double> z, nz;
|
||||
|
||||
for(unsigned int i(0); i<fEnergy.size(); i++) {
|
||||
if(!(fEnergy[i] - e)) {
|
||||
z = DataZ(e);
|
||||
nz = DataNZ(e);
|
||||
break;
|
||||
}
|
||||
}
|
||||
|
||||
bool found = false;
|
||||
unsigned int i;
|
||||
for (i=0; i<z.size(); i++) {
|
||||
if (z[i]/10.0 >= zz) {
|
||||
found = true;
|
||||
break;
|
||||
}
|
||||
}
|
||||
|
||||
if (!found)
|
||||
return -1.0;
|
||||
|
||||
if (i == 0)
|
||||
return nz[0]*10.0*zz/z[0];
|
||||
|
||||
return fabs(nz[i-1]+(nz[i]-nz[i-1])*(10.0*zz-z[i-1])/(z[i]-z[i-1]));
|
||||
}
|
75
src/external/TFitPofB-lib/include/TBofZCalc.h
vendored
Normal file
75
src/external/TFitPofB-lib/include/TBofZCalc.h
vendored
Normal file
@ -0,0 +1,75 @@
|
||||
/***************************************************************************
|
||||
|
||||
TBofZCalc.h
|
||||
|
||||
Author: Bastian M. Wojek
|
||||
e-mail: bastian.wojek@psi.ch
|
||||
|
||||
2008/05/23
|
||||
|
||||
***************************************************************************/
|
||||
|
||||
#ifndef _TBofZCalc_H_
|
||||
#define _TBofZCalc_H_
|
||||
|
||||
#include <vector>
|
||||
using namespace std;
|
||||
|
||||
// Base class for any kind of theory function B(z)
|
||||
|
||||
class TBofZCalc {
|
||||
|
||||
public:
|
||||
|
||||
TBofZCalc() {}
|
||||
|
||||
~TBofZCalc() {
|
||||
fZ.clear();
|
||||
fBZ.clear();
|
||||
}
|
||||
|
||||
vector<double> DataZ() const {return fZ;}
|
||||
vector<double> DataBZ() const {return fBZ;}
|
||||
double GetBmin() const;
|
||||
double GetBmax() const;
|
||||
double GetBofZ(double) const;
|
||||
double GetdZ() const {return fDZ;}
|
||||
|
||||
protected:
|
||||
vector<double> fZ;
|
||||
vector<double> fBZ;
|
||||
double fDZ;
|
||||
};
|
||||
|
||||
// Class "for Meissner screening" in a thin superconducting film
|
||||
|
||||
class TLondon1D_1L : public TBofZCalc {
|
||||
|
||||
public:
|
||||
|
||||
TLondon1D_1L( const vector<double>& );
|
||||
|
||||
};
|
||||
|
||||
|
||||
// Class "for Meissner screening" in a thin superconducting film - bilayer with two different lambdas
|
||||
|
||||
class TLondon1D_2L : public TBofZCalc {
|
||||
|
||||
public:
|
||||
|
||||
TLondon1D_2L( const vector<double>& );
|
||||
|
||||
};
|
||||
|
||||
// Class "for Meissner screening" in a thin superconducting film - tri-layer with two different lambdas
|
||||
|
||||
class TLondon1D_3L : public TBofZCalc {
|
||||
|
||||
public:
|
||||
|
||||
TLondon1D_3L( const vector<double>& );
|
||||
|
||||
};
|
||||
|
||||
#endif // _BofZCalc_H_
|
33
src/external/TFitPofB-lib/include/TFitPofB.h
vendored
Normal file
33
src/external/TFitPofB-lib/include/TFitPofB.h
vendored
Normal file
@ -0,0 +1,33 @@
|
||||
/***************************************************************************
|
||||
|
||||
TFitPofB.h
|
||||
|
||||
Author: Bastian M. Wojek
|
||||
e-mail: bastian.wojek@psi.ch
|
||||
|
||||
2008/05/23
|
||||
|
||||
***************************************************************************/
|
||||
|
||||
#ifndef _TFitPofB_H_
|
||||
#define _TFitPofB_H_
|
||||
|
||||
#include "TPofTCalc.h"
|
||||
|
||||
class TFitPofB {
|
||||
|
||||
public:
|
||||
TFitPofB(const vector<unsigned int>& , const vector<double>&);
|
||||
~TFitPofB();
|
||||
|
||||
double Eval(double, const vector<double>&) const;
|
||||
|
||||
private:
|
||||
mutable vector<double> fPar;
|
||||
TTrimSPData *fImpProfile;
|
||||
TPofTCalc *fPofT;
|
||||
mutable bool fCalcNeeded;
|
||||
|
||||
};
|
||||
|
||||
#endif //_TFitPofB_H_
|
38
src/external/TFitPofB-lib/include/TPofBCalc.h
vendored
Normal file
38
src/external/TFitPofB-lib/include/TPofBCalc.h
vendored
Normal file
@ -0,0 +1,38 @@
|
||||
/***************************************************************************
|
||||
|
||||
TPofBCalc.h
|
||||
|
||||
Author: Bastian M. Wojek
|
||||
e-mail: bastian.wojek@psi.ch
|
||||
|
||||
2008/05/23
|
||||
|
||||
***************************************************************************/
|
||||
|
||||
#ifndef _TPofBCalc_H_
|
||||
#define _TPofBCalc_H_
|
||||
|
||||
#include "TBofZCalc.h"
|
||||
#include "TTrimSPDataHandler.h"
|
||||
|
||||
class TPofBCalc {
|
||||
|
||||
public:
|
||||
|
||||
TPofBCalc( const TBofZCalc&, const TTrimSPData&, const vector<double>& );
|
||||
~TPofBCalc() {
|
||||
fB.clear();
|
||||
fPB.clear();
|
||||
}
|
||||
|
||||
vector<double> DataB() const {return fB;}
|
||||
vector<double> DataPB() const {return fPB;}
|
||||
|
||||
private:
|
||||
vector<double> fB;
|
||||
vector<double> fPB;
|
||||
static const double gBar = 0.0135538817;
|
||||
|
||||
};
|
||||
|
||||
#endif // _TPofBCalc_H_
|
42
src/external/TFitPofB-lib/include/TPofTCalc.h
vendored
Normal file
42
src/external/TFitPofB-lib/include/TPofTCalc.h
vendored
Normal file
@ -0,0 +1,42 @@
|
||||
/***************************************************************************
|
||||
|
||||
TPofTCalc.h
|
||||
|
||||
Author: Bastian M. Wojek
|
||||
e-mail: bastian.wojek@psi.ch
|
||||
|
||||
2008/05/23
|
||||
|
||||
***************************************************************************/
|
||||
|
||||
#ifndef _TPofTCalc_H_
|
||||
#define _TPofTCalc_H_
|
||||
|
||||
#include "TPofBCalc.h"
|
||||
#include "fftw3.h"
|
||||
|
||||
class TPofTCalc {
|
||||
|
||||
public:
|
||||
TPofTCalc(const vector<double>&);
|
||||
~TPofTCalc();
|
||||
|
||||
vector<double> DataT() const {return fT;}
|
||||
vector<double> DataPT() const {return fPT;}
|
||||
void DoFFT(const TPofBCalc&, const vector<double>&);
|
||||
double Eval(double) const;
|
||||
|
||||
private:
|
||||
fftw_plan fFFTplan;
|
||||
fftw_complex *fFFTout;
|
||||
double *fFFTin;
|
||||
vector<double> fT;
|
||||
vector<double> fPT;
|
||||
double fTBin;
|
||||
unsigned int fNFFT;
|
||||
static const double PI = 3.14159265358979323846;
|
||||
static const double gBar = 0.0135538817;
|
||||
|
||||
};
|
||||
|
||||
#endif // _TPofTCalc_H_
|
41
src/external/TFitPofB-lib/include/TTrimSPDataHandler.h
vendored
Normal file
41
src/external/TFitPofB-lib/include/TTrimSPDataHandler.h
vendored
Normal file
@ -0,0 +1,41 @@
|
||||
/***************************************************************************
|
||||
|
||||
TTrimSPDataHandler.h
|
||||
|
||||
Author: Bastian M. Wojek
|
||||
e-mail: bastian.wojek@psi.ch
|
||||
|
||||
2008/05/23
|
||||
|
||||
***************************************************************************/
|
||||
|
||||
#ifndef _TTrimSPDataHandler_H_
|
||||
#define _TTrimSPDataHandler_H_
|
||||
|
||||
#include <vector>
|
||||
#include <string>
|
||||
using namespace std;
|
||||
|
||||
class TTrimSPData {
|
||||
|
||||
public:
|
||||
|
||||
TTrimSPData(const string&, vector<string>&);
|
||||
|
||||
~TTrimSPData() {
|
||||
fDataZ.clear();
|
||||
fDataNZ.clear();
|
||||
}
|
||||
|
||||
vector<double> Energy() const {return fEnergy;}
|
||||
vector<double> DataZ(double) const;
|
||||
vector<double> DataNZ(double) const;
|
||||
double GetNofZ(double, double) const;
|
||||
|
||||
private:
|
||||
vector<double> fEnergy;
|
||||
vector< vector<double> > fDataZ;
|
||||
vector< vector<double> > fDataNZ;
|
||||
};
|
||||
|
||||
#endif // _TTrimSPDataHandler_H_
|
5
src/external/scripts/mlog2db
vendored
5
src/external/scripts/mlog2db
vendored
@ -223,12 +223,9 @@ awk -v parFOUR=$4 -v parFIVE=$5 -v parSIX=$6 '{
|
||||
dataArray[i+3] = runNumber[j]
|
||||
break
|
||||
}
|
||||
else dataArray[i+3] = "0000"
|
||||
}
|
||||
|
||||
# if(runNumber[2] == "his") {dataArray[i+3] = substr(runNumber[3],1,4) }
|
||||
# else if(runNumber[2] == "pta") {dataArray[i+3] = substr(runNumber[4],1,4)}
|
||||
# else {dataArray[i+3] = runNumber[2] }
|
||||
|
||||
negErrArray[i] = ""
|
||||
negErrArray[i+1] = ""
|
||||
negErrArray[i+2] = ""
|
||||
|
Loading…
x
Reference in New Issue
Block a user