Fixed bugs and added a few comments after first testing.

This commit is contained in:
Bastian M. Wojek 2008-05-24 15:05:18 +00:00
parent 728ec00816
commit a8046f80d0
12 changed files with 642 additions and 175 deletions

View File

@ -5,7 +5,7 @@
Author: Bastian M. Wojek Author: Bastian M. Wojek
e-mail: bastian.wojek@psi.ch e-mail: bastian.wojek@psi.ch
2008/05/23 2008/05/24
***************************************************************************/ ***************************************************************************/
@ -14,6 +14,9 @@
#include <iostream> #include <iostream>
#include <algorithm> #include <algorithm>
//------------------
// Method returning the field minimum of the model
//------------------
double TBofZCalc::GetBmin() const { double TBofZCalc::GetBmin() const {
vector<double>::const_iterator iter; vector<double>::const_iterator iter;
@ -22,6 +25,10 @@ double TBofZCalc::GetBmin() const {
return *iter; return *iter;
} }
//------------------
// Method returning the field maximum of the model
//------------------
double TBofZCalc::GetBmax() const { double TBofZCalc::GetBmax() const {
vector<double>::const_iterator iter; vector<double>::const_iterator iter;
iter = max_element(fBZ.begin(),fBZ.end()); iter = max_element(fBZ.begin(),fBZ.end());
@ -29,6 +36,10 @@ double TBofZCalc::GetBmax() const {
return *iter; return *iter;
} }
//------------------
// Method returning the field B(z) at a given z according to the model
//------------------
double TBofZCalc::GetBofZ(double zz) const { double TBofZCalc::GetBofZ(double zz) const {
bool found = false; bool found = false;
@ -50,10 +61,13 @@ double TBofZCalc::GetBofZ(double zz) const {
} }
TLondon1D_1L::TLondon1D_1L(const vector<double> &param) { //------------------
// Constructor of the TLondon1D_1L class
// 1D-London screening in a thin superconducting film
// Parameters: Bext[G], deadlayer[nm], thickness[nm], lambda[nm]
//------------------
// Parameters for 1D-London screening in a thin superconducting film TLondon1D_1L::TLondon1D_1L(const vector<double> &param) {
// Bext[G], deadlayer[nm], thickness[nm], lambda[nm]
unsigned int n(2000); // number of steps for the calculation unsigned int n(2000); // number of steps for the calculation
@ -71,12 +85,15 @@ TLondon1D_1L::TLondon1D_1L(const vector<double> &param) {
} }
//------------------
// Constructor of the TLondon1D_2L class
// 1D-London screening in a thin superconducting film, two layers, two lambda
// Parameters: Bext[G], deadlayer[nm], thickness1[nm], thickness2[nm], lambda1[nm], lambda2[nm]
//------------------
TLondon1D_2L::TLondon1D_2L(const vector<double> &param) { TLondon1D_2L::TLondon1D_2L(const vector<double> &param) {
// Parameters for 1D-London screening in a thin superconducting film, two layers, two lambda unsigned int n(4000); // number of steps for the calculation
// 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 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); double N2(4.0*N1);
@ -97,12 +114,15 @@ TLondon1D_2L::TLondon1D_2L(const vector<double> &param) {
} }
//------------------
// Constructor of the TLondon1D_3L class
// 1D-London screening in a thin superconducting film, three layers, two lambdas (top and bottom layer: lambda1)
// Parameters: Bext[G], deadlayer[nm], thickness1[nm], thickness2[nm], thickness3[nm], lambda1[nm], lambda2[nm]
//------------------
TLondon1D_3L::TLondon1D_3L(const vector<double> &param) { TLondon1D_3L::TLondon1D_3L(const vector<double> &param) {
// Parameters for 1D-London screening in a thin superconducting film, three layers, two lambdas (top and bottom layer: lambda1) unsigned int n(4000); // number of steps for the calculation
// 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 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])));

View File

@ -5,7 +5,7 @@
Author: Bastian M. Wojek Author: Bastian M. Wojek
e-mail: bastian.wojek@psi.ch e-mail: bastian.wojek@psi.ch
2008/05/23 2008/05/24
***************************************************************************/ ***************************************************************************/
@ -13,6 +13,11 @@
#include <iostream> #include <iostream>
using namespace std; using namespace std;
//------------------
// Constructor of the TFitPofB class -- reading available implantation profiles and
// creates (a pointer to) the TPofTCalc object (with the FFT plan)
//------------------
TFitPofB::TFitPofB(const vector<unsigned int> &parNo, const vector<double> &par) : fCalcNeeded(true) { TFitPofB::TFitPofB(const vector<unsigned int> &parNo, const vector<double> &par) : fCalcNeeded(true) {
for(unsigned int i(0); i<parNo.size(); i++) { for(unsigned int i(0); i<parNo.size(); i++) {
@ -41,27 +46,42 @@ TFitPofB::TFitPofB(const vector<unsigned int> &parNo, const vector<double> &par)
} }
//------------------
// Destructor of the TFitPofB class -- cleaning up
//------------------
TFitPofB::~TFitPofB() { TFitPofB::~TFitPofB() {
fPar.clear(); fPar.clear();
fImpProfile = 0;
delete fImpProfile; delete fImpProfile;
fPofT = 0; fImpProfile = 0;
delete fPofT; delete fPofT;
fPofT = 0;
} }
//------------------
// Method that calls the procedures to create B(z), p(B) and P(t)
// It finally returns P(t) for a given t.
// Parameters: all the parameters for the function to be fitted through TFitPofB
//------------------
double TFitPofB::Eval(double t, const vector<double> &par) const { double TFitPofB::Eval(double t, const vector<double> &par) const {
// check if any parameter has changed // check if any parameter has changed
bool par_changed(false); bool par_changed(false);
bool only_phase_changed(false);
for (unsigned int i(0); i<fPar.size(); i++) { for (unsigned int i(0); i<fPar.size(); i++) {
if( fPar[i]-par[i] ) { if( fPar[i]-par[i] ) {
fPar[i] = par[i]; fPar[i] = par[i];
par_changed = true; par_changed = true;
if(i == 2 || i == 3) { if(i == 0 || i == 2 || i == 3 || i == 4) {
cout << "You are varying dt or dB! These parameters have to be fixed! Quitting..." << endl; cout << "You are varying the model-flag, dt, dB or E! These parameters have to be fixed! Quitting..." << endl;
exit(0); exit(0);
} else if (i == 1) {
only_phase_changed = true;
} else {
only_phase_changed = false;
} }
} }
} }
@ -69,49 +89,81 @@ double TFitPofB::Eval(double t, const vector<double> &par) const {
if (par_changed) if (par_changed)
fCalcNeeded = true; fCalcNeeded = true;
// if parameters have changed, recalculate B(z), P(B) and P(t) // if model parameters have changed, recalculate B(z), P(B) and P(t)
/* DEBUGGING CODE COMMENTED -- quite a mess... sorry*/
if (fCalcNeeded) { if (fCalcNeeded) {
vector<double> par_for_BofZ;
vector<double> par_for_PofB;
vector<double> par_for_PofT; vector<double> par_for_PofT;
for (unsigned int i(5); i<par.size(); i++) // cout << "par_for_PofT: ";
par_for_BofZ.push_back(par[i]);
for (unsigned int i(2); i<5; i++) for (unsigned int i(1); i<4; i++) {
par_for_PofB.push_back(par[i]);
for (unsigned int i(1); i<4; i++)
par_for_PofT.push_back(par[i]); par_for_PofT.push_back(par[i]);
// cout << par[i] << " ";
}
// cout << endl;
switch (int(par[0])) { if(!only_phase_changed) {
case 1:
{ cout << " Parameters have changed, (re-)calculating p(B) and P(t) now..." << endl;
TLondon1D_1L BofZ1(par_for_BofZ);
TPofBCalc PofB1(BofZ1, *fImpProfile, par_for_PofB); vector<double> par_for_BofZ;
fPofT->DoFFT(PofB1, par_for_PofT); vector<double> par_for_PofB;
// cout << "par_for_BofZ: ";
for (unsigned int i(5); i<par.size(); i++) {
par_for_BofZ.push_back(par[i]);
// cout << par[i] << " ";
} }
break; // cout << endl;
case 2:
{ // cout << "par_for_PofB: ";
TLondon1D_2L BofZ2(par_for_BofZ);
TPofBCalc PofB2(BofZ2, *fImpProfile, par_for_PofB); for (unsigned int i(2); i<5; i++) {
fPofT->DoFFT(PofB2, par_for_PofT); par_for_PofB.push_back(par[i]);
// cout << par[i] << " ";
} }
break; // cout << endl;
case 3:
{ switch (int(par[0])) {
TLondon1D_3L BofZ3(par_for_BofZ); case 1:
TPofBCalc PofB3(BofZ3, *fImpProfile, par_for_PofB); {
fPofT->DoFFT(PofB3, par_for_PofT); // cout << "Found the 1D-London model." << endl;
} TLondon1D_1L BofZ1(par_for_BofZ);
break; TPofBCalc PofB1(BofZ1, *fImpProfile, par_for_PofB);
default: fPofT->DoFFT(PofB1);
cout << "The user function you specified, does not exist. Quitting..." << endl; }
exit(0); break;
case 2:
{
// cout << "Found the 1D-London model.2L" << endl;
TLondon1D_2L BofZ2(par_for_BofZ);
TPofBCalc PofB2(BofZ2, *fImpProfile, par_for_PofB);
fPofT->DoFFT(PofB2);
}
break;
case 3:
{
// cout << "Found the 1D-London model.3L" << endl;
TLondon1D_3L BofZ3(par_for_BofZ);
TPofBCalc PofB3(BofZ3, *fImpProfile, par_for_PofB);
fPofT->DoFFT(PofB3);
}
break;
default:
cout << "The user function you specified with the first parameter, does not exist. Quitting..." << endl;
exit(0);
} }
} else {
cout << "Only the phase parameter has changed, (re-)calculating P(t) now..." << endl;
}
fPofT->CalcPol(par_for_PofT);
fCalcNeeded = false;
} }

View File

@ -5,24 +5,32 @@
Author: Bastian M. Wojek Author: Bastian M. Wojek
e-mail: bastian.wojek@psi.ch e-mail: bastian.wojek@psi.ch
2008/05/23 2008/05/24
***************************************************************************/ ***************************************************************************/
#include "TPofBCalc.h" #include "TPofBCalc.h"
#include <cmath> #include <cmath>
#include <iostream> #include <iostream>
#include <fstream>
//#include <cstdio>
//-----------
// Constructor that does the P(B) calculation for given B(z) and n(z)
// Parameters: dt[us], dB[G], Energy[keV]
//-----------
TPofBCalc::TPofBCalc( const TBofZCalc &BofZ, const TTrimSPData &dataTrimSP, const vector<double> &para ) { TPofBCalc::TPofBCalc( const TBofZCalc &BofZ, const TTrimSPData &dataTrimSP, const vector<double> &para ) {
// Parameters: dt[us], dB[G], Energy[keV] fBmin = BofZ.GetBmin();
fBmax = BofZ.GetBmax();
double BB, BBnext; double BB, BBnext;
double zm, zp, zNextm, zNextp, dz; double zm, zp, zNextm, zNextp, dz;
// fill not used Bs before Bmin with 0.0 // fill not used Bs before Bmin with 0.0
for ( BB = 0.0 ; BB < BofZ.GetBmin() ; BB += para[1] ) { for ( BB = 0.0 ; BB < fBmin ; BB += para[1] ) {
fB.push_back(BB); fB.push_back(BB);
fPB.push_back(0.0); fPB.push_back(0.0);
} }
@ -34,10 +42,25 @@ TPofBCalc::TPofBCalc( const TBofZCalc &BofZ, const TTrimSPData &dataTrimSP, cons
vector<double> bofzZ(BofZ.DataZ()); vector<double> bofzZ(BofZ.DataZ());
vector<double> bofzBZ(BofZ.DataBZ()); vector<double> bofzBZ(BofZ.DataBZ());
double ddZ(BofZ.GetdZ()); double ddZ(BofZ.GetdZ());
double Bmax(BofZ.GetBmax());
/* USED FOR DEBUGGING-----------------------------------
cout << "Bmin = " << fBmin << ", Bmax = " << fBmax << endl;
char debugfile[50];
int n = sprintf (debugfile, "test_Bz_%f.dat", fBmin);
if (n > 0) {
ofstream of(debugfile);
for (unsigned int i(0); i<bofzZ.size(); i++) {
of << bofzZ[i] << " " << bofzBZ[i] << endl;
}
of.close();
}
---------------------------------------------------------*/
double nn; double nn;
for ( ; BB <= Bmax ; BB += para[1]) { for ( ; BB <= fBmax ; BB += para[1]) {
BBnext = BB + para[1]; BBnext = BB + para[1];
fB.push_back(BB); fB.push_back(BB);
fPB.push_back(0.0); fPB.push_back(0.0);
@ -89,23 +112,23 @@ TPofBCalc::TPofBCalc( const TBofZCalc &BofZ, const TTrimSPData &dataTrimSP, cons
double BmaxFFT(1.0/gBar/para[0]); double BmaxFFT(1.0/gBar/para[0]);
//// cout << "N = " << int(BmaxFFT/para[1]+1.0) << endl; // cout << "N = " << int(BmaxFFT/para[1]+1.0) << endl;
for ( ; BB <= BmaxFFT ; BB += para[1] ) { for ( ; BB <= BmaxFFT ; BB += para[1] ) {
fB.push_back(BB); fB.push_back(BB);
fPB.push_back(0.0); fPB.push_back(0.0);
} }
// make sure that we have an even number of elements in p(B) for FFT // make sure that we have an even number of elements in p(B) for FFT, so we do not have to care later
if (fB.size() % 2) { if (fB.size() % 2) {
fB.push_back(BB); fB.push_back(BB);
fPB.push_back(0.0); fPB.push_back(0.0);
} }
cout << "size of B = " << fB.size() << ", size of p(B) = " << fPB.size() << endl; // cout << "size of B = " << fB.size() << ", size of p(B) = " << fPB.size() << endl;
// normalize pB // normalize p(B)
double pBsum = 0.0; double pBsum = 0.0;
for (unsigned int i(firstZerosEnd); i<lastZerosStart; i++) for (unsigned int i(firstZerosEnd); i<lastZerosStart; i++)

View File

@ -5,7 +5,7 @@
Author: Bastian M. Wojek Author: Bastian M. Wojek
e-mail: bastian.wojek@psi.ch e-mail: bastian.wojek@psi.ch
2008/05/23 2008/05/24
***************************************************************************/ ***************************************************************************/
@ -13,14 +13,18 @@
#include "fftw3.h" #include "fftw3.h"
#include <cmath> #include <cmath>
#include <cstdio> #include <cstdio>
#include <fstream>
#include <iostream> #include <iostream>
//------------------
// Constructor of the TPofTCalc class - it creates the FFT plan
// Parameters: phase, dt, dB
//------------------
TPofTCalc::TPofTCalc (const vector<double> &par) { 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); 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]; 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); fFFTin = (double *)malloc(sizeof(double) * fNFFT);
fFFTout = (fftw_complex *)fftw_malloc(sizeof(fftw_complex) * (fNFFT/2+1)); fFFTout = (fftw_complex *)fftw_malloc(sizeof(fftw_complex) * (fNFFT/2+1));
@ -48,10 +52,31 @@ TPofTCalc::TPofTCalc (const vector<double> &par) {
} }
void TPofTCalc::DoFFT(const TPofBCalc &PofB, const vector<double> &par) { //--------------
// Method that does the FFT of a given p(B)
//--------------
void TPofTCalc::DoFFT(const TPofBCalc &PofB) {
vector<double> pB(PofB.DataPB()); vector<double> pB(PofB.DataPB());
/* USED FOR DEBUGGING -----------------------
vector<double> B(PofB.DataB());
double Bmin(PofB.GetBmin());
char debugfile[50];
int n = sprintf (debugfile, "test_PB_%f.dat", Bmin);
if (n > 0) {
ofstream of(debugfile);
for (unsigned int i(0); i<B.size(); i++) {
of << B[i] << " " << pB[i] << endl;
}
of.close();
}
--------------------------------------------*/
for (unsigned int i(0); i<fNFFT; i++) { for (unsigned int i(0); i<fNFFT; i++) {
fFFTin[i] = pB[i]; fFFTin[i] = pB[i];
} }
@ -64,20 +89,32 @@ void TPofTCalc::DoFFT(const TPofBCalc &PofB, const vector<double> &par) {
fftw_execute(fFFTplan); fftw_execute(fFFTplan);
// Calculate polarisation }
//---------------------
// Method for calculating the muon spin polarization P(t) from the Fourier transformed p(B)
// Parameters: phase, dt, dB
//---------------------
void TPofTCalc::CalcPol(const vector<double> &par) {
double sinph(sin(par[0]*PI/180.0)), cosph(cos(par[0]*PI/180.0)); double sinph(sin(par[0]*PI/180.0)), cosph(cos(par[0]*PI/180.0));
double polTemp(0.0); double polTemp(0.0);
fT.clear();
fPT.clear();
for (unsigned int i(0); i<fNFFT/2+1; i++){ for (unsigned int i(0); i<fNFFT/2+1; i++){
fT.push_back(double(i)*fTBin); fT.push_back(double(i)*fTBin);
polTemp = cosph*fFFTout[i][0]*par[2] + sinph*fFFTout[i][1]*par[2]; polTemp = cosph*fFFTout[i][0]*par[2] + sinph*fFFTout[i][1]*par[2];
fPT.push_back(polTemp); fPT.push_back(polTemp);
} }
} }
//---------------------
// Method for evaluating P(t) at a given t
//---------------------
double TPofTCalc::Eval(double t) const { double TPofTCalc::Eval(double t) const {
for (unsigned int i(0); i<fT.size()-1; i++) { for (unsigned int i(0); i<fT.size()-1; i++) {
if (t < fT[i+1]) { if (t < fT[i+1]) {
@ -85,12 +122,16 @@ double TPofTCalc::Eval(double t) const {
} }
} }
cout << "No data for the time " << t << " us available! Quitting..." << endl; cout << "No data for the time " << t << " us available! Returning -999.0 ..." << endl;
exit(0); return -999.0;
} }
//---------------------
// Destructor of the TPofTCalc class - it saves the FFT plan and cleans up
//---------------------
TPofTCalc::~TPofTCalc() { TPofTCalc::~TPofTCalc() {
// export wisdom // export wisdom so it has not to be checked for the FFT-plan next time
FILE *wordsOfWisdomW; FILE *wordsOfWisdomW;
wordsOfWisdomW = fopen("WordsOfWisdom.dat", "w"); wordsOfWisdomW = fopen("WordsOfWisdom.dat", "w");
@ -110,4 +151,3 @@ TPofTCalc::~TPofTCalc() {
} }

View File

@ -5,7 +5,7 @@
Author: Bastian M. Wojek Author: Bastian M. Wojek
e-mail: bastian.wojek@psi.ch e-mail: bastian.wojek@psi.ch
2008/05/23 2008/05/24
***************************************************************************/ ***************************************************************************/
@ -17,20 +17,24 @@
using namespace std; using namespace std;
// TrimSPData constructor -- reading all available trim.SP-rge-files into std::vectors //--------------------
// Constructor of the TrimSPData class -- reading all available trim.SP-rge-files with a given name into std::vectors
// The rge-file names have to contain the Implantation energy just before the rge-extension in the format %02u_%1u
// Example: string path("/home/user/TrimSP/SomeSample-");
// string energyArr[] = {"02_1", "02_5", "03_5", "05_0", "07_5", "10_0", "12_5"};
// vector<string> energyVec(energyArr, energyArr+(sizeof(energyArr)/sizeof(energyArr[0])));
//
// This will read the files "/home/user/TrimSP/SomeSample-02_1.rge", "/home/user/TrimSP/SomeSample-02_5.rge" and so on.
//--------------------
TTrimSPData::TTrimSPData(const string &path, vector<string> &energyVec) { 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); double zz(0.0), nzz(0.0);
vector<double> vzz, vnzz; vector<double> vzz, vnzz;
string word, energyStr; string word, energyStr;
for(unsigned int i(0); i<energyVec.size(); i++){ 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"; energyStr = path + energyVec[i] + ".rge";
fEnergy.push_back(atof(energyVec[i].replace(2,1,".").c_str())); fEnergy.push_back(atof(energyVec[i].replace(2,1,".").c_str()));
@ -63,7 +67,9 @@ TTrimSPData::TTrimSPData(const string &path, vector<string> &energyVec) {
} }
} }
// vector<double> DataZ(double) -- returning z-vector calculated by trim.SP for energy given energy //---------------------
// Method returning z-vector calculated by trim.SP for given energy[keV]
//---------------------
vector<double> TTrimSPData::DataZ(double e) const { vector<double> TTrimSPData::DataZ(double e) const {
@ -73,13 +79,15 @@ vector<double> TTrimSPData::DataZ(double e) const {
return fDataZ[i]; return fDataZ[i];
} }
} }
// default // default
cout << "No implantation profile available for the specified energy... You got back the first one." << endl; cout << "No implantation profile available for the specified energy... You get back the first one." << endl;
return fDataZ[0]; return fDataZ[0];
} }
// vector<double> DataNZ(double) -- returning n(z)-vector calculated by trim.SP for given energy //---------------------
// Method returning n(z)-vector calculated by trim.SP for given energy[keV]
//---------------------
vector<double> TTrimSPData::DataNZ(double e) const { vector<double> TTrimSPData::DataNZ(double e) const {
@ -89,12 +97,14 @@ vector<double> TTrimSPData::DataNZ(double e) const {
} }
} }
// default // default
cout << "No implantation profile available for the specified energy... You got back the first one." << endl; cout << "No implantation profile available for the specified energy... You get back the first one." << endl;
return fDataNZ[0]; return fDataNZ[0];
} }
// double GetNofZ(double, double) -- return NofZ for given z in nanometers //---------------------
// Method returning n(z) for given z[nm] and energy[keV]
//---------------------
double TTrimSPData::GetNofZ(double zz, double e) const { double TTrimSPData::GetNofZ(double zz, double e) const {

View File

@ -5,7 +5,7 @@
Author: Bastian M. Wojek Author: Bastian M. Wojek
e-mail: bastian.wojek@psi.ch e-mail: bastian.wojek@psi.ch
2008/05/23 2008/05/24
***************************************************************************/ ***************************************************************************/
@ -15,7 +15,10 @@
#include <vector> #include <vector>
using namespace std; using namespace std;
//--------------------
// Base class for any kind of theory function B(z) // Base class for any kind of theory function B(z)
// In principle only constructors for the different models have to be implemented
//--------------------
class TBofZCalc { class TBofZCalc {
@ -41,7 +44,9 @@ protected:
double fDZ; double fDZ;
}; };
//--------------------
// Class "for Meissner screening" in a thin superconducting film // Class "for Meissner screening" in a thin superconducting film
//--------------------
class TLondon1D_1L : public TBofZCalc { class TLondon1D_1L : public TBofZCalc {
@ -51,8 +56,9 @@ public:
}; };
//--------------------
// Class "for Meissner screening" in a thin superconducting film - bilayer with two different lambdas // Class "for Meissner screening" in a thin superconducting film - bilayer with two different lambdas
//--------------------
class TLondon1D_2L : public TBofZCalc { class TLondon1D_2L : public TBofZCalc {
@ -62,7 +68,9 @@ public:
}; };
//--------------------
// Class "for Meissner screening" in a thin superconducting film - tri-layer with two different lambdas // Class "for Meissner screening" in a thin superconducting film - tri-layer with two different lambdas
//--------------------
class TLondon1D_3L : public TBofZCalc { class TLondon1D_3L : public TBofZCalc {

View File

@ -5,7 +5,7 @@
Author: Bastian M. Wojek Author: Bastian M. Wojek
e-mail: bastian.wojek@psi.ch e-mail: bastian.wojek@psi.ch
2008/05/23 2008/05/24
***************************************************************************/ ***************************************************************************/

View File

@ -5,7 +5,7 @@
Author: Bastian M. Wojek Author: Bastian M. Wojek
e-mail: bastian.wojek@psi.ch e-mail: bastian.wojek@psi.ch
2008/05/23 2008/05/24
***************************************************************************/ ***************************************************************************/
@ -27,10 +27,14 @@ public:
vector<double> DataB() const {return fB;} vector<double> DataB() const {return fB;}
vector<double> DataPB() const {return fPB;} vector<double> DataPB() const {return fPB;}
double GetBmin() const {return fBmin;}
double GetBmax() const {return fBmax;}
private: private:
vector<double> fB; vector<double> fB;
vector<double> fPB; vector<double> fPB;
double fBmin;
double fBmax;
static const double gBar = 0.0135538817; static const double gBar = 0.0135538817;
}; };

View File

@ -5,7 +5,7 @@
Author: Bastian M. Wojek Author: Bastian M. Wojek
e-mail: bastian.wojek@psi.ch e-mail: bastian.wojek@psi.ch
2008/05/23 2008/05/24
***************************************************************************/ ***************************************************************************/
@ -23,7 +23,8 @@ public:
vector<double> DataT() const {return fT;} vector<double> DataT() const {return fT;}
vector<double> DataPT() const {return fPT;} vector<double> DataPT() const {return fPT;}
void DoFFT(const TPofBCalc&, const vector<double>&); void DoFFT(const TPofBCalc&);
void CalcPol(const vector<double>&);
double Eval(double) const; double Eval(double) const;
private: private:

View File

@ -5,7 +5,7 @@
Author: Bastian M. Wojek Author: Bastian M. Wojek
e-mail: bastian.wojek@psi.ch e-mail: bastian.wojek@psi.ch
2008/05/23 2008/05/24
***************************************************************************/ ***************************************************************************/

View File

@ -0,0 +1,33 @@
CXX = g++
CXXFLAGS = -g -Wall
LOCALINCLUDE = ../include/
INCLUDES = -I$(LOCALINCLUDE)
LD = g++
LDFLAGS = -g -L../classes/ -lTFitPofB -lfftw3 -lm
EXEC = test
# some definitions: headers, sources, objects,...
OBJS =
OBJS += $(EXEC).o
# make the executable:
#
all: $(EXEC)
$(EXEC): $(OBJS)
@echo "---> Building $(EXEC) ..."
$(LD) $(LDFLAGS) $(OBJS) -o $(EXEC)
@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)
@echo "---> removing $(OBJS)"
#
$(OBJS): %.o: %.cpp
$(CXX) $(INCLUDES) $(CXXFLAGS) -c $<

276
src/external/TFitPofB-lib/test/test.cpp vendored Normal file
View File

@ -0,0 +1,276 @@
#include "TFitPofB.h"
#include <iostream>
#include <fstream>
using namespace std;
int main(){
/* 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])));
TTrimSPData calcData(rge_path, energy_vec);
vector<double> energies(calcData.Energy());
for (unsigned int i(0); i<energies.size(); i++)
cout << energies[i] << endl;
vector<double> z(calcData.DataZ(2.5));
vector<double> nz(calcData.DataNZ(2.5));
vector<double> z2(calcData.DataZ(25.0));
vector<double> nz2(calcData.DataNZ(25.0));
ofstream of("test_out1.dat");
for (unsigned int i(0); i<z.size(); i++) {
of << z[i] << " " << nz[i] << endl;
}
of.close();
ofstream of2("test_out2.dat");
for (unsigned int i(0); i<z2.size(); i++) {
of2 << z2[i] << " " << nz2[i] << endl;
}
of2.close();
ofstream of3("test_out4.dat");
for (unsigned int i(0); i<z.size(); i++) {
of3 << i << " " << calcData.GetNofZ(double(i), 25.0) << endl;
}
of3.close();
double param[8] = {100.0, 5.0, 50.0, 100.0, 40.0, 0.01, 0.1, 15.0};
vector<double> parameter(param,param+8);
vector<double> param_for_BofZ;
vector<double> param_for_PofB;
vector<double> param_for_PofT;
for (unsigned int i(0); i<4; i++)
param_for_BofZ.push_back(parameter[i]);
for (unsigned int i(5); i<8; i++)
param_for_PofB.push_back(parameter[i]);
for (unsigned int i(4); i<7; i++)
param_for_PofT.push_back(parameter[i]);
TLondon1D_1L bofz(param_for_BofZ);
cout << "Bmin = " << bofz.GetBmin() << endl;
cout << "Bmax = " << bofz.GetBmax() << endl;
cout << "dZ = " << bofz.GetdZ() << endl;
ofstream of5("test_Bz.dat");
for (double i(0); i<50.; i+=0.1) {
of5 << i << " " << bofz.GetBofZ(i) << endl;
}
of5.close();
ofstream of6("test_zBz.dat");
for (unsigned int i(0); i<bofz.DataZ().size(); i++) {
of6 << bofz.DataZ()[i] << " " << bofz.DataBZ()[i] << endl;
}
of6.close();
TPofBCalc pofb(bofz, calcData, param_for_PofB);
cout << "Output to file now..." << endl;
vector<double> hurgaB(pofb.DataB());
vector<double> hurgaPB(pofb.DataPB());
ofstream of7("test_BpB.dat");
for (unsigned int i(0); i<hurgaB.size(); i++) {
of7 << hurgaB[i] << " " << hurgaPB[i] << endl;
}
of7.close();
TPofTCalc poft(param_for_PofT);
poft.DoFFT(pofb, param_for_PofT);
ofstream of8("test_tpt4.dat");
for (double i(0.); i<12.0; i+=0.003) {
of8 << i << " " << poft(i) << endl;
}
of8.close();
*/
unsigned int parNo_arr[] = {1, 3, 5, 7, 9, 11, 12, 13, 14, 15, 16, 17};
double par_arr[] = {3.0, 999.0, 0.0, 999.0, 0.01, 999.0, 0.05, 999.0, 20.0, 999.0, 100.0, 10.0, 80.0, 50.0, 80.0, 180.0, 500.0};
vector<unsigned int> parNo_vec(parNo_arr, parNo_arr+(sizeof(parNo_arr)/sizeof(parNo_arr[0])));
vector<double> par_vec(par_arr, par_arr+(sizeof(par_arr)/sizeof(par_arr[0])));
vector<double> par_vec_sub;
for(unsigned int i(0); i<parNo_vec.size(); i++) {
par_vec_sub.push_back(par_vec[parNo_vec[i]-1]);
}
TFitPofB fitter(parNo_vec, par_vec);
ofstream of01("test_fitter01.dat");
ofstream of02("test_fitter02.dat");
ofstream of03("test_fitter03.dat");
ofstream of04("test_fitter04.dat");
ofstream of05("test_fitter05.dat");
ofstream of06("test_fitter06.dat");
ofstream of07("test_fitter07.dat");
ofstream of08("test_fitter08.dat");
ofstream of09("test_fitter09.dat");
ofstream of10("test_fitter10.dat");
for (double i(0.); i<12.0; i+=0.003) {
of01 << i << " " << fitter.Eval(i, par_vec_sub) << endl;
}
of01.close();
par_vec_sub[1] += 10.0;
par_vec_sub[11] -= 20.0;
for (double i(0.); i<12.0; i+=0.003) {
of02 << i << " " << fitter.Eval(i, par_vec_sub) << endl;
}
of02.close();
par_vec_sub[1] += 10.0;
par_vec_sub[11] -= 20.0;
for (double i(0.); i<12.0; i+=0.003) {
of03 << i << " " << fitter.Eval(i, par_vec_sub) << endl;
}
of03.close();
par_vec_sub[1] += 10.0;
for (double i(0.); i<12.0; i+=0.003) {
of04 << i << " " << fitter.Eval(i, par_vec_sub) << endl;
}
of04.close();
par_vec_sub[1] += 10.0;
par_vec_sub[11] -= 20.0;
for (double i(0.); i<12.0; i+=0.003) {
of05 << i << " " << fitter.Eval(i, par_vec_sub) << endl;
}
of05.close();
par_vec_sub[1] += 10.0;
par_vec_sub[11] -= 20.0;
for (double i(0.); i<12.0; i+=0.003) {
of06 << i << " " << fitter.Eval(i, par_vec_sub) << endl;
}
of06.close();
par_vec_sub[1] += 10.0;
for (double i(0.); i<12.0; i+=0.003) {
of07 << i << " " << fitter.Eval(i, par_vec_sub) << endl;
}
of07.close();
par_vec_sub[1] += 10.0;
par_vec_sub[7] = 75.0;
for (double i(0.); i<12.0; i+=0.003) {
of08 << i << " " << fitter.Eval(i, par_vec_sub) << endl;
}
of08.close();
par_vec_sub[1] += 10.0;
par_vec_sub[11] -= 20.0;
for (double i(0.); i<12.0; i+=0.003) {
of09 << i << " " << fitter.Eval(i, par_vec_sub) << endl;
}
of09.close();
par_vec_sub[1] += 10.0;
for (double i(0.); i<12.0; i+=0.003) {
of10 << i << " " << fitter.Eval(i, par_vec_sub) << endl;
}
of10.close();
/* vector<double> data01, data02, data03, data04, data05, data06, data07, data08, data09, data10;
for (double i(0.); i<12.0; i+=0.003)
data01.push_back(fitter.Eval(i, par_vec_sub));
par_vec_sub[1] += 10.0;
par_vec_sub[8] -= 10.0;
for (double i(0.); i<12.0; i+=0.003)
data02.push_back(fitter.Eval(i, par_vec_sub));
par_vec_sub[1] += 10.0;
par_vec_sub[8] -= 10.0;
for (double i(0.); i<12.0; i+=0.003)
data03.push_back(fitter.Eval(i, par_vec_sub));
par_vec_sub[1] += 10.0;
for (double i(0.); i<12.0; i+=0.003)
data04.push_back(fitter.Eval(i, par_vec_sub));
par_vec_sub[1] += 10.0;
par_vec_sub[8] -= 10.0;
for (double i(0.); i<12.0; i+=0.003)
data05.push_back(fitter.Eval(i, par_vec_sub));
par_vec_sub[1] += 10.0;
par_vec_sub[8] -= 10.0;
for (double i(0.); i<12.0; i+=0.003)
data06.push_back(fitter.Eval(i, par_vec_sub));
par_vec_sub[1] += 10.0;
for (double i(0.); i<12.0; i+=0.003)
data07.push_back(fitter.Eval(i, par_vec_sub));
par_vec_sub[1] += 10.0;
par_vec_sub[7] = 190.0;
for (double i(0.); i<12.0; i+=0.003)
data08.push_back(fitter.Eval(i, par_vec_sub));
par_vec_sub[1] += 10.0;
par_vec_sub[8] -= 10.0;
for (double i(0.); i<12.0; i+=0.003)
data09.push_back(fitter.Eval(i, par_vec_sub));
par_vec_sub[1] += 10.0;
for (double i(0.); i<12.0; i+=0.003)
data10.push_back(fitter.Eval(i, par_vec_sub));
*/
parNo_vec.clear();
par_vec.clear();
par_vec_sub.clear();
return EXIT_SUCCESS;
}