1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
|
#include <cassert>
#include <iostream>
#include <cuda.h>
#include <cuda_runtime.h>
#include "DataFormats/Common/interface/DataFrame.h"
#include "DataFormats/DetId/interface/DetId.h"
#include "DataFormats/HcalDetId/interface/HcalDetId.h"
#include "DataFormats/HcalDigi/interface/HBHEDataFrame.h"
#include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
#include "DataFormats/HcalDigi/interface/QIE10DataFrame.h"
#include "DataFormats/HcalDigi/interface/QIE11DataFrame.h"
#include "HeterogeneousCore/CUDAUtilities/interface/requireDevices.h"
__global__ void kernel_test_hcal_qiesample(HcalQIESample *sample, uint16_t value) {
printf("kernel: testing hcal qie sampel\n");
printf("%f %f %f\n", nominal_adc2fc[0], nominal_adc2fc[1], nominal_adc2fc[2]);
HcalQIESample tmp{value};
*sample = tmp;
}
__global__ void kernel_test_hcal_qie8_hbhedf(HBHEDataFrame *df) {
printf("kernel: testing hcal hbhe dataframe\n");
df->setSize(10);
for (auto i = 0; i < 10; i++)
df->setSample(i, HcalQIESample(100));
df->setReadoutIds(HcalElectronicsId(100));
}
void test_hcal_qiesample() {
HcalQIESample h_sample, h_test_sample0{100}, h_test_sample1;
HcalQIESample *d_sample;
cudaMalloc((void **)&d_sample, sizeof(HcalQIESample));
cudaMemcpy(d_sample, &h_sample, sizeof(HcalQIESample), cudaMemcpyHostToDevice);
kernel_test_hcal_qiesample<<<1, 1>>>(d_sample, 100);
cudaMemcpy(&h_sample, d_sample, sizeof(HcalQIESample), cudaMemcpyDeviceToHost);
assert(h_sample() == h_test_sample0());
assert(h_sample() != h_test_sample1());
}
template <typename TDF>
__global__ void kernel_test_hcal_qie8_digis(TDF *pdfs, uint32_t *out) {
int id = threadIdx.x;
uint32_t sum = 0;
for (auto i = 0; i < 10; i++)
sum += pdfs[id].sample(i).raw();
out[id] = sum;
}
template <typename TDF>
__global__ void kernel_test_hcal_qie1011_digis(uint16_t *pdfs, uint32_t *out, int samples) {
printf("kernel: testing hcal qie1011 df\n");
int id = threadIdx.x;
uint32_t sum = 0;
int nwords = TDF::WORDS_PER_SAMPLE * samples + TDF::HEADER_WORDS + TDF::FLAG_WORDS;
TDF df(edm::DataFrame(0, pdfs + id * nwords, nwords));
for (auto i = 0; i < df.samples(); i++) {
sum += df[i].adc();
}
out[id] = sum;
}
template <typename TDF>
void test_hcal_qie1011_digis() {
constexpr int size = 10;
constexpr int samples = 10;
constexpr int detid = 2;
HcalDataFrameContainer<TDF> coll{samples, detid};
uint16_t *d_data;
uint32_t *d_out;
uint32_t h_out[size], h_test_out[size];
for (auto i = 0; i < size; i++) {
// #words per single TDF
uint16_t tmp[TDF::WORDS_PER_SAMPLE * samples + TDF::HEADER_WORDS + TDF::FLAG_WORDS];
h_test_out[i] = 0;
for (auto j = TDF::HEADER_WORDS; j < TDF::WORDS_PER_SAMPLE * samples + TDF::HEADER_WORDS; j++) {
tmp[j] = 100;
}
TDF df(edm::DataFrame(0, tmp, TDF::WORDS_PER_SAMPLE * samples + TDF::HEADER_WORDS + TDF::FLAG_WORDS));
for (auto j = 0; j < df.samples(); j++)
h_test_out[i] += df[j].adc();
coll.addDataFrame(DetId{(uint32_t)i}, (uint16_t *)&tmp);
}
cudaMalloc((void **)&d_data,
size * (TDF::WORDS_PER_SAMPLE * samples + TDF::HEADER_WORDS + TDF::FLAG_WORDS) * sizeof(uint16_t));
cudaMalloc((void **)&d_out, size * sizeof(uint32_t));
cudaMemcpy(d_data,
coll.frame(0),
size * (TDF::WORDS_PER_SAMPLE * samples + TDF::HEADER_WORDS + TDF::FLAG_WORDS) * sizeof(uint16_t),
cudaMemcpyHostToDevice);
kernel_test_hcal_qie1011_digis<TDF><<<1, size>>>(d_data, d_out, samples);
cudaDeviceSynchronize();
auto code = cudaGetLastError();
if (code != cudaSuccess)
std::cout << cudaGetErrorString(code);
cudaMemcpy(&h_out, d_out, size * sizeof(uint32_t), cudaMemcpyDeviceToHost);
// comparison
for (auto i = 0; i < size; i++) {
std::cout << h_out[i] << " == " << h_test_out[i] << std::endl;
assert(h_out[i] == h_test_out[i]);
}
}
template <typename TDF>
void test_hcal_qie8_digis() {
constexpr int n = 10;
edm::SortedCollection<TDF> coll{n};
TDF *d_dfs;
uint32_t *d_out;
uint32_t h_out[n], h_test_out[n];
for (auto i = 0; i < n; i++) {
TDF &df = coll[i];
df.setSize(10);
h_test_out[i] = 0;
uint32_t test = 0;
for (auto j = 0; j < 10; j++) {
df.setSample(j, HcalQIESample(100));
h_test_out[i] += df.sample(j).raw();
test += df.sample(j).raw();
}
}
cudaMalloc((void **)&d_dfs, n * sizeof(TDF));
cudaMalloc((void **)&d_out, n * sizeof(uint32_t));
cudaMemcpy(d_dfs, &(*coll.begin()), n * sizeof(TDF), cudaMemcpyHostToDevice);
kernel_test_hcal_qie8_digis<<<1, n>>>(d_dfs, d_out);
cudaMemcpy(&h_out, d_out, n * sizeof(uint32_t), cudaMemcpyDeviceToHost);
std::cout << "collection size = " << coll.size() << std::endl;
// comparison
for (auto i = 0; i < n; i++) {
std::cout << h_out[i] << " == " << h_test_out[i] << std::endl;
assert(h_out[i] == h_test_out[i]);
}
}
void test_hcal_qie8_hbhedf() {
HBHEDataFrame h_df, h_test_df;
HBHEDataFrame *d_df;
h_test_df.setSize(10);
for (auto i = 0; i < 10; i++)
h_test_df.setSample(i, HcalQIESample(100));
h_test_df.setReadoutIds(HcalElectronicsId(100));
cudaMalloc((void **)&d_df, sizeof(HBHEDataFrame));
cudaMemcpy(d_df, &h_df, sizeof(HBHEDataFrame), cudaMemcpyHostToDevice);
kernel_test_hcal_qie8_hbhedf<<<1, 1>>>(d_df);
cudaMemcpy(&h_df, d_df, sizeof(HBHEDataFrame), cudaMemcpyDeviceToHost);
assert(h_df.size() == h_test_df.size());
assert(h_df.elecId() == h_test_df.elecId());
for (auto i = 0; i < 10; i++)
assert(h_df[i].raw() == h_test_df[i].raw());
}
int main(int argc, char **argv) {
cms::cudatest::requireDevices();
// qie8
test_hcal_qiesample();
test_hcal_qie8_hbhedf();
test_hcal_qie8_digis<HBHEDataFrame>();
test_hcal_qie8_digis<HFDataFrame>();
test_hcal_qie8_digis<HODataFrame>();
// qie1011
test_hcal_qie1011_digis<QIE10DataFrame>();
test_hcal_qie1011_digis<QIE11DataFrame>();
return 0;
}
|