00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020 #ifndef TBASKETSQL_CXX
00021 #define TBASKETSQL_CXX
00022
00023 #include "TBasket.h"
00024 #include "TTree.h"
00025 #include "TBranch.h"
00026 #include "TFile.h"
00027 #include "TMath.h"
00028 #include "TBasketSQL.h"
00029 #include <Riostream.h>
00030 #include <vector>
00031 #include "TTreeSQL.h"
00032 #include "TBufferSQL.h"
00033
00034 ClassImp(TBasketSQL)
00035
00036 namespace std {} using namespace std;
00037
00038
00039 TBasketSQL::TBasketSQL() : TBasket()
00040 {
00041
00042 }
00043
00044
00045 TBasketSQL::TBasketSQL(const char *name, const char *title, TBranch *branch,
00046 TSQLResult ** rs, TString *insert_query,
00047 vector<Int_t> *vc, TSQLRow **r) :
00048 fResultPtr(rs),fRowPtr(r)
00049 {
00050
00051
00052 SetName(name);
00053 SetTitle(title);
00054 fClassName = "TBasketSQL";
00055 fBufferSize = branch->GetBasketSize();
00056 fNevBufSize = branch->GetEntryOffsetLen();
00057 fNevBuf = 0;
00058 fEntryOffset = 0;
00059 fDisplacement= 0;
00060 fBuffer = 0;
00061 fInsertQuery = insert_query;
00062
00063 if (vc==0) {
00064 fBufferRef = 0;
00065 } else {
00066 fBufferRef = new TBufferSQL(TBuffer::kWrite, fBufferSize, vc, fInsertQuery, fRowPtr);
00067 }
00068 fHeaderOnly = kTRUE;
00069 fLast = 0;
00070
00071 fBuffer = 0;
00072 fBranch = branch;
00073 fHeaderOnly = kFALSE;
00074 if (fNevBufSize) {
00075 fEntryOffset = new Int_t[fNevBufSize];
00076 for (Int_t i=0;i<fNevBufSize;i++) fEntryOffset[i] = 0;
00077 }
00078 branch->GetTree()->IncrementTotalBuffers(fBufferSize);
00079 }
00080
00081
00082 TBasketSQL::~TBasketSQL()
00083 {
00084
00085 }
00086
00087
00088 void TBasketSQL::CreateBuffer(const char *name, TString title,
00089 vector<Int_t> *vc,
00090 TBranch *branch, TSQLResult ** rs)
00091 {
00092
00093
00094 fResultPtr = rs;
00095 SetName(name);
00096 SetTitle(title);
00097 fClassName = "TBasketSQL";
00098 fBufferSize = branch->GetBasketSize();
00099 fNevBufSize = branch->GetEntryOffsetLen();
00100 fNevBuf = 0;
00101 fEntryOffset = 0;
00102 fDisplacement= 0;
00103 fBuffer = 0;
00104
00105 if(vc==0) {
00106 fBufferRef = 0;
00107 Error("CreateBuffer","Need a vector of columns\n");
00108 } else {
00109 fBufferRef = new TBufferSQL(TBuffer::kWrite, fBufferSize, vc, fInsertQuery, fRowPtr);
00110 }
00111 fHeaderOnly = kTRUE;
00112 fLast = 0;
00113
00114 fBuffer = 0;
00115 fBranch = branch;
00116 fHeaderOnly = kFALSE;
00117 if (fNevBufSize) {
00118 fEntryOffset = new Int_t[fNevBufSize];
00119 for (Int_t i=0;i<fNevBufSize;i++) fEntryOffset[i] = 0;
00120 }
00121 branch->GetTree()->IncrementTotalBuffers(fBufferSize);
00122 }
00123
00124
00125 void TBasketSQL::PrepareBasket(Long64_t entry)
00126 {
00127
00128
00129 ((TBufferSQL*)fBufferRef)->ResetOffset();
00130 ((TTreeSQL*)fBranch->GetTree())->PrepEntry(entry);
00131 fBufferRef->Reset();
00132 }
00133
00134
00135 Int_t TBasketSQL::ReadBasketBytes(Long64_t , TFile *)
00136 {
00137
00138
00139 Error("ReadBasketBytes","This member function should not be called!");
00140 return 0;
00141 }
00142
00143
00144 Int_t TBasketSQL::ReadBasketBuffers(Long64_t , Int_t, TFile *)
00145 {
00146
00147
00148 Error("ReadBasketBuffers","This member function should not be called!");
00149 return 0;
00150 }
00151
00152
00153 void TBasketSQL::Reset()
00154 {
00155
00156
00157 TBasket::Reset();
00158 }
00159
00160
00161
00162 void TBasketSQL::Update(Int_t, Int_t)
00163 {
00164
00165
00166 ((TBufferSQL*)fBufferRef)->ResetOffset();
00167 fNevBuf++;
00168 }
00169
00170
00171 #endif