• MIT的J-P. Péraud and N. G. Hadjiconstantinou用Matlab写的一维动力学MC教程 https://github.com/jeanphilippeperaud/1D_KMC_matlab
  • C++的二维系综MC https://github.com/jeanphilippeperaud/Phonon-Code
  • 印度人写的Matlab写的系综MC https://github.com/abhipath90/MCBTE
  • 东京大学Nomura实验室的Roman Anufriev 的基于Python的跟踪MC https://github.com/anufrievroman/freepaths
  • 全带的almaBTE https://almabte.bitbucket.io
  • 东京大学Shiomi和邵成博士的pTrans https://github.com/TEEL-UTokyo/P-TRANS/tree/main

这些项目都很好,但很多规模比较小,很多项目发布之后就没人再维护了,自己写程序的时候也总是觉得内力太差,因为的确需要在程序开发、物理本身等等上都要很熟悉,像北大的陈默涵老师或者中科院的李新亮老师一样...

之前的文章里也介绍过Python的一维系综和动力学MC,这篇文章是来自于并行计算的一个课程作业,把之前一维系综MC的程序用C++写了一遍,记录一下自己的学习历程,回顾一下本科时学习的C++知识...

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
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
#include <iostream>
#include <vector>
#include <random>
#include <cmath>
#include <algorithm>
#include <iomanip>
#include <fstream>
// #include <omp.h>
#include <chrono>

double hbar = 1.05457e-034;
double k_B = 1.380649e-23 ;

// 声子结构体
struct Phonon {
double x; // 当前位置
double next_x; // 下一位置
double mu; // 方向
int property_index; // 频率索引
int sign; // 方向标志
bool remove; // 是否删除
};

class PhononMonteCarlo {
public:
double Teq; // 参考温度
double Tl; // 左侧高温边界
double Tr; // 右侧低温边界
double L; // 薄膜长度
int Nx; // 空间网格数目
double * x_grid_lines; // 空间网格线
double dx; // 空间步长
double dt; // 时间步长
int Nt; // 时间步数
double Eeff; // 等效声子团能量
double energy_hot_emitted_per_dt; // 每个时间步高温边界发射的能量数
double energy_cold_emitted_per_dt; // 每个时间步低温边界发射的能量数
int phonon_hot_emitted_per_dt; // 每个时间步高温边界发射的声子数
int phonon_cold_emitted_per_dt; // 每个时间步低温边界发射的声子数

std::chrono::nanoseconds delete_time{0};
std::chrono::nanoseconds drift_time{0};
std::chrono::nanoseconds boundary_time{0};
std::chrono::nanoseconds temperature_time{0};
std::chrono::nanoseconds scatter_time{0};
std::chrono::nanoseconds delete_scatter_time{0};


// 通过文件读取以及后续计算得到的Phonon性质

std::vector<double> frequency; // 频率
std::vector<double> DoS; // 态密度
std::vector<double> velocity; // 速度
std::vector<double> frequency_width; // 频率宽度
std::vector<double> relaxation_time; // 弛豫时间
std::vector<int> branch; // 分支
std::vector<double> heat_capacity; // 声子比热容
double total_heat_capacity; // 总比热容

// 抽样矩阵信息
std::vector<double> cumsum_scatter; // 累计散射抽样矩阵
std::vector<double> cumsum_emit; // 累计边界发射抽样矩阵

// 储存温度场、热流等信息

std::vector<Phonon> phonons; // 声子数组
std::vector<std::vector<int>> scatterd_phonon_index; // 发生散射的声子索引
double ** temperature; // 温度场
std::vector<double> drift_projection; // 热流投影

// 随机数发生器
std::uniform_real_distribution<double> rand_dist;

PhononMonteCarlo(double Teq, double Tl, double Tr, double L, int Nx, double dt, int Nt, double Eeff = 1E-5) : Teq(Teq), Tl(Tl), Tr(Tr), L(L), Nx(Nx), dt(dt), Nt(Nt), Eeff(Eeff), rand_dist(0, 1), scatterd_phonon_index(Nx){

// 初始化温度场
this->temperature = new double*[Nt];
for (int i = 0; i < Nt; ++i) {
this->temperature[i] = new double[Nx];
for (int j = 0; j < Nx; ++j) {
temperature[i][j] = Teq;
}
}

// 划分网格
this->dx = L / Nx;
this->x_grid_lines = new double[Nx + 1];
for (int i = 0; i < Nx + 1; ++i) {
x_grid_lines[i] = i * this->dx;
}

}

~PhononMonteCarlo() {
// 释放动态分配的内存
for (int i = 0; i < Nt; ++i) {
delete[] temperature[i];
}
delete[] temperature;
delete[] x_grid_lines;
}

void loadDispersionData(const std::string& filename) {
std::ifstream file(filename);
if (!file) {
std::cerr << "Failed to open file: " << filename << std::endl;
return;
}

double omega, dos, v, freq_width, relax_time;
int p_b;
while (file >> omega >> dos >> v >> freq_width >> relax_time >> p_b) {
this->frequency.push_back(omega);
this->DoS.push_back(dos);
this->velocity.push_back(v);
this->frequency_width.push_back(freq_width);
this->relaxation_time.push_back(relax_time);
this->branch.push_back(p_b);
this->heat_capacity.push_back(this->calculateEnergyTemperatureDerivative(omega) * freq_width * dos);
}
this->total_heat_capacity = std::accumulate(heat_capacity.begin(), heat_capacity.end(), 0);
file.close();
}

void initialize() {

// 生成散射过程的概率抽样矩阵
std::vector<double> scatter_term(this->frequency.size());
this->cumsum_scatter = std::vector<double>(this->frequency.size());

for (int i = 0; i < this->frequency.size(); ++i) {
scatter_term[i] = this->DoS[i] * this->calculateEnergyTemperatureDerivative(this->frequency[i]) * this->frequency_width[i] / this->relaxation_time[i];
}

std::partial_sum(scatter_term.begin(), scatter_term.end(), this->cumsum_scatter.begin());
double total_sum = this->cumsum_scatter.back();
std::transform(this->cumsum_scatter.begin(), this->cumsum_scatter.end(), this->cumsum_scatter.begin(), [total_sum](double value) { return value / total_sum; });

// 生成边界发射过程的概率抽样矩阵
std::vector<double> emit_term(this->frequency.size());
this->cumsum_emit = std::vector<double>(this->frequency.size());

for (int i = 0; i < this->frequency.size(); ++i) {
emit_term[i] = scatter_term[i] * this->relaxation_time[i] * this->velocity[i];
}

std::partial_sum(emit_term.begin(), emit_term.end(), this->cumsum_emit.begin());
total_sum = this->cumsum_emit.back();
std::transform(this->cumsum_emit.begin(), this->cumsum_emit.end(), this->cumsum_emit.begin(), [total_sum](double value) { return value / total_sum;});

// 判断边界发射的声子数
this->energy_hot_emitted_per_dt = total_sum * this->dt * (this->Tl - this->Teq) / 4;
this->energy_cold_emitted_per_dt = total_sum * this->dt * (this->Teq - this->Tr) / 4;
this->phonon_hot_emitted_per_dt = static_cast<int>(this->energy_hot_emitted_per_dt / this->Eeff);
this->phonon_cold_emitted_per_dt = static_cast<int>(this->energy_cold_emitted_per_dt / this->Eeff);

// 初始化drift-热流数组
this->drift_projection = std::vector<double>(this->frequency.size());
}

int checkSign(double num) {
if (num > 0) {
return 1;
} else if (num < 0) {
return -1;
} else {
return 0;
}
}

// 累积分布函数抽样函数
int sampleFromCDF(const std::vector<double>& cdf, double randomValue) {

for (std::size_t i = 0; i < cdf.size(); ++i) {
if (randomValue <= cdf[i]) {
return i;
}
}

return -1; // 未找到合适的索引
}

// 在已排序的 double 数组中搜索给定值的函数
int findPhononPosition(double x) {
auto it = std::lower_bound(this->x_grid_lines, this->x_grid_lines+Nx, x);
return std::distance(x_grid_lines, it) - 1;
}

double boseEinsteinDistribution(double omega, double T) {
return 1.0 / (std::exp(hbar * omega / (k_B * T)) - 1.0);
}

double calculateEnergyTemperatureDerivative(double omega) {
return std::pow(hbar * omega, 2) / k_B / (this->Teq * this->Teq) * this->boseEinsteinDistribution(omega, this->Teq) * (this->boseEinsteinDistribution(omega, this->Teq) +1);
}

void runSimulation() {

std::random_device rd;
std::mt19937 gen(rd());

std::vector<int> removed_phonon_index; // 将被移除的声子索引
std::vector<int> selected_phonon_index;
std::vector<int> unselected_phonon_index;

for (int t = 0; t < Nt; ++t) {

auto start = std::chrono::high_resolution_clock::now();

std::fill(drift_projection.begin(), drift_projection.end(), 0.0);
removed_phonon_index.clear();
// Phonon drift
// #pragma omp parallel for
for (int phonon_index = 0; phonon_index < phonons.size(); ++phonon_index) {
Phonon& phonon = phonons[phonon_index];
phonon.next_x = phonon.x + velocity[phonon.property_index] * phonon.mu * dt;
if (phonon.next_x < 0) {
phonon.next_x = 0;
phonon.remove = true;
}
else if (phonon.next_x > L) {
phonon.next_x = L;
phonon.remove = true;
}
drift_projection[phonon.property_index] += phonon.sign * (phonon.next_x - phonon.x);
}

auto end = std::chrono::high_resolution_clock::now();

auto duration = std::chrono::duration_cast<std::chrono::nanoseconds>(end - start);

drift_time += duration;

start = std::chrono::high_resolution_clock::now();

// 从后往前遍历 phonons 列表
for (auto it = phonons.rbegin(); it != phonons.rend(); ++it) {
// 检查 remove 的值
if (it->remove) {
// 删除对应的对象
phonons.erase(std::next(it).base());
}
}

end = std::chrono::high_resolution_clock::now();

duration = std::chrono::duration_cast<std::chrono::nanoseconds>(end - start);

delete_time += duration;

// Boundary emission

start = std::chrono::high_resolution_clock::now();

// hot boundary emission
std::vector<Phonon> private_phonons; // 私有的phonon数组
std::vector<double> private_drift_projection(drift_projection.size(), 0.0); // 私有的漂移投影

std::fill(private_drift_projection.begin(), private_drift_projection.end(), 0.0);
// #pragma omp parallel private(private_phonons)
{
private_phonons.clear(); // 清空临时私有变量

// #pragma omp for
for (int i = 0; i < phonon_hot_emitted_per_dt; ++i) {
double mu = std::sqrt(rand_dist(gen));
int phonon_index = sampleFromCDF(this->cumsum_emit, rand_dist(gen));
Phonon new_phonon = { 0, 0, mu, phonon_index, 1, false };
new_phonon.next_x = new_phonon.x + velocity[new_phonon.property_index] * new_phonon.mu * dt * rand_dist(gen);
if (new_phonon.next_x >= L) {
new_phonon.next_x = L;
private_drift_projection[new_phonon.property_index] += new_phonon.sign * (new_phonon.next_x - new_phonon.x);
}
else {
private_drift_projection[new_phonon.property_index] += new_phonon.sign * (new_phonon.next_x - new_phonon.x);
private_phonons.push_back(new_phonon);
}
}

// #pragma omp critical
{
for (size_t i = 0; i < drift_projection.size(); ++i) {
drift_projection[i] += private_drift_projection[i];
}
phonons.insert(phonons.end(), private_phonons.begin(), private_phonons.end());
}
}

// cold boundary emission
std::fill(private_drift_projection.begin(), private_drift_projection.end(), 0.0);

// #pragma omp parallel private(private_phonons)
{
private_phonons.clear(); // 清空临时私有变量

// #pragma omp for
for (int i = 0; i < phonon_cold_emitted_per_dt; ++i) {
double mu = - std::sqrt(rand_dist(gen));
int phonon_index = sampleFromCDF(this->cumsum_emit, rand_dist(gen));
Phonon new_phonon = { L, 0, mu, phonon_index, -1, false};
new_phonon.next_x = new_phonon.x + velocity[new_phonon.property_index] * new_phonon.mu * dt * rand_dist(gen);
if (new_phonon.next_x <= 0) {
new_phonon.next_x = 0;
private_drift_projection[new_phonon.property_index] += new_phonon.sign * (new_phonon.next_x - new_phonon.x);
}
else {
private_drift_projection[new_phonon.property_index] += new_phonon.sign * (new_phonon.next_x - new_phonon.x);
private_phonons.push_back(new_phonon);
}
}

// #pragma omp critical
{
for (size_t i = 0; i < drift_projection.size(); ++i) {
drift_projection[i] += private_drift_projection[i];
}
phonons.insert(phonons.end(), private_phonons.begin(), private_phonons.end());
}
}

end = std::chrono::high_resolution_clock::now();

duration = std::chrono::duration_cast<std::chrono::nanoseconds>(end - start);

boundary_time += duration;


// Count heat flux
double heat_flux = std::accumulate(drift_projection.begin(), drift_projection.end(), 0.0) / dt * Eeff / L;

// std::cout << heat_flux / (2.0 / L) << std::endl;


// Move phonons
// #pragma omp parallel for
for (int i = 0; i < phonons.size(); ++i) {
phonons[i].x = phonons[i].next_x;
}

// #pragma omp parallel for
for (auto& array : scatterd_phonon_index) {
array.clear();
}

start = std::chrono::high_resolution_clock::now();

// Calculate temperature distribution and found the scattered phonons
for (int i = 0; i < phonons.size(); ++i) {
int position_index = findPhononPosition(phonons[i].x);
temperature[t][position_index] += Eeff * phonons[i].sign / dx / total_heat_capacity;
if (rand_dist(gen) < 1 - std::exp(- dt / relaxation_time[phonons[i].property_index])) {
scatterd_phonon_index[position_index].push_back(i);
}
}

end = std::chrono::high_resolution_clock::now();

duration = std::chrono::duration_cast<std::chrono::nanoseconds>(end - start);

temperature_time += duration;

start = std::chrono::high_resolution_clock::now();

// Scatter process
unselected_phonon_index.clear();

for (int position_index = 0; position_index < Nx; ++position_index) {

int net_scatterd_phonon_number = 0;

for (auto& phonon_index : scatterd_phonon_index[position_index]) {
net_scatterd_phonon_number += phonons[phonon_index].sign;
}
net_scatterd_phonon_number = std::abs(net_scatterd_phonon_number);
int newly_generated_phonon_sign = checkSign(temperature[t][position_index] - Teq);
// 发生散射的phonon,在当前的网格里正负抵消,只剩下剩下的净能量的phonon,这些净能量的phonon数量大大小于之前的发生散射的这些phonon,在发生散射的phonon里随机挑一些对应这些净能量phonon,然后重新抽样phonon性质。
selected_phonon_index.clear();

if (net_scatterd_phonon_number > 0 && !scatterd_phonon_index[position_index].empty()) {
std::sample(scatterd_phonon_index[position_index].begin(), scatterd_phonon_index[position_index].end(), std::back_inserter(selected_phonon_index),
net_scatterd_phonon_number, gen);
}

// 把这些被选中的phonon,重新抽样
for (auto& phonon_index : selected_phonon_index) {
phonons[phonon_index].mu = 2 * rand_dist(gen) - 1;
phonons[phonon_index].sign = newly_generated_phonon_sign;
phonons[phonon_index].property_index = sampleFromCDF(this->cumsum_scatter, rand_dist(gen));
}

std::sort(scatterd_phonon_index[position_index].begin(), scatterd_phonon_index[position_index].end());
std::sort(selected_phonon_index.begin(), selected_phonon_index.end());
std::set_difference(scatterd_phonon_index[position_index].begin(), scatterd_phonon_index[position_index].end(), selected_phonon_index.begin(), selected_phonon_index.end(), std::back_inserter(unselected_phonon_index));
}

end = std::chrono::high_resolution_clock::now();

duration = std::chrono::duration_cast<std::chrono::nanoseconds>(end - start);

scatter_time += duration;

start = std::chrono::high_resolution_clock::now();

std::sort(unselected_phonon_index.begin(), unselected_phonon_index.end());
for (auto it = unselected_phonon_index.rbegin(); it != unselected_phonon_index.rend(); ++it) {
phonons.erase(phonons.begin() + *it);
}

end = std::chrono::high_resolution_clock::now();

duration = std::chrono::duration_cast<std::chrono::nanoseconds>(end - start);

delete_scatter_time += duration;
}
}

void exportTemperatureDistribution(const std::string& filename) {

// 创建新文件并写入数据
std::remove(filename.c_str());
std::ofstream file(filename);
for (int t = 0; t < Nt; ++t) {
for (int i = 0; i < Nx; ++i) {
file << temperature[t][i] << "\n";
}
}

file.close();
}

};

int main() {

double Teq = 300.0;
double Tl = 301.0;
double Tr = 299.0;
double L = 1e-7;
int Nx = 20;
double dt = 5e-12;
int Nt = 50;

PhononMonteCarlo simulation(Teq, Tl, Tr, L, Nx, dt, Nt, 1e-6);
simulation.loadDispersionData("dataSi.txt");
simulation.initialize();

auto start = std::chrono::high_resolution_clock::now();
simulation.runSimulation();
auto end = std::chrono::high_resolution_clock::now();

// 计算函数执行时间
auto duration = std::chrono::duration_cast<std::chrono::nanoseconds>(end - start);

// 输出时间
std::cout << "删除运行时间: " << simulation.delete_time.count() << " 毫秒" << std::endl;
std::cout << "漂移运行时间: " << simulation.drift_time.count() << " 毫秒" << std::endl;
std::cout << "边界发射运行时间: " << simulation.boundary_time.count() << " 毫秒" << std::endl;
std::cout << "计算温度: " << simulation.temperature_time.count() << " 毫秒" << std::endl;
std::cout << "散射运行时间: " << simulation.scatter_time.count() << " 毫秒" << std::endl;
std::cout << "删除散射运行时间: " << simulation.delete_scatter_time.count() << " 毫秒" << std::endl;
std::cout << "函数运行时间: " << duration.count() << " 毫秒" << std::endl;

simulation.exportTemperatureDistribution("temperature-ensemble.txt");

return 0;
}

需要用C++17以上去编译,另外在源文件目录下下载好硅的色散文件(https://github.com/jeanphilippeperaud/1D_KMC_matlab)。模拟的体系和参数设置为左侧为高温热沉,温度为301K,右侧为低温热沉,温度为299K,参考温度选取为300K。模拟的薄膜长度为100nm,均匀离散为了20个网格,模拟的时间步长为1ps,模拟的时间步数为100步,其中每个声子能量团代表的能量设置为5E-6J。

image-20230718153637732

loadDispersionData用于从文件中加载声子的色散关系数据,并依据给定的性质计算声子的比热容等属性。

initialize依照给定的声子性质,计算散射过程的概率抽样矩阵以及边界发射过程的概率抽样矩阵。加载的色散关系数据将声子的性质按照频率离散成了数组,计算得到的概率抽样数组和声子性质数组是一一对应的。

checkSign用于检查给定数值的正负,主要用于判断在某一个点格点发生散射后新生成的声子能量团的符号。

sampleFromCDF用于根据累积分布函数抽样索引,用于确定边界发射以及散射后的声子性质。

findPhononPosition在已排序的网格线数组中找到给定位置对应的网格索引,主要用于确定当前声子能量团所在的格点。

boseEinsteinDistribution根据给定的频率和温度,计算玻色-爱因斯坦分布函数的值并返回。

calculateEnergyTemperatureDerivative根据给定的频率,计算能量对温度的导数的值并返回。

exportTemperatureDistribution将温度分布导出到文件。

runSimulation运行声子蒙特卡洛模拟,整个模拟沿着时间步进行推进,首先对所有声子进行漂移过程的更新,根据漂移后的位置判断是否需要删除该声子。之后边界发射新的声子并更新声子的位置,之后对于每个声子,根据其位置更新温度场并记录发生散射的声子索引。在散射索引记录以后,遍历整个空间网格,统计每个网格内发生散射的声子数目,处理网格内的散射过程。假如发生散射的包括\(N^+\)个正能量声子团,\(N^-\)个负能量声子团,那么最后需要在这些声子团中随机选出\(|N^+ - N^-|\)个声子,重新抽样他们的性质,剩下\(|N^+ + N^-| - |N^+ - N^-|\)的声子被删除。处理结束后,进入下一时间步,直到模拟结束。

image-20230718152056835

在整个模拟过程中,各个部分所占用的时间比例大约为:

• 漂移过程: 0.35%

• 删除抵达边界的声子: 28.84%

• 边界发射: 15.99%

• 计算温度: 4.94%

• 遍历空间网格进行散射: 18.68%

• 删除散射过程中要删去的声子: 31.00%

删除声子的过程占据了绝大部分的程序运行时间,这是由系综蒙特卡洛模拟的特性决定的。在整个模拟过程中,由于存在着边界发射的声子诞生过程,以及散射和抵达边界的声子湮灭过程,系统中总的声子数目是不确定的,这个时候程序中必须维护着一个动态数组,在当前的程序中,phonons数组储存了体系中所有的声子。频繁地对数组进行插入和删除操作是非常耗时的,而且由于删除操作会改变整个数组的索引,因此这部分程序也不容易并行化。而容易实现并行化的过程,如漂移和计算温度,在耗时上又只占用了很小的比例。因此系综蒙特卡洛程序本身的并行性并不好,可以优化的部分就是边界发射和遍历空间网格进行散射的过程。虽然可以用MPI_Bcast将声子分发到各个进程中,再使用MPI_Gather等操作将结果归约到0进程中,但整体来看优化效果其实十分有限。在这一点上,声子跟踪蒙特卡洛模拟有很大的优势。

所用到的C++知识

没想到什么太好的方式去组织这篇文章,因为本质上就是零碎的知识..

模版类

动态数组

VSCODE中C++项目配置

VSCODE只是一个文本编辑器,在VSCODE中写C++的程序时,

这篇文章分成这几块内容,

  • Vscode里C++的项目配置;

  • STL模板类的使用;

  • C++科学计算库的管理;

今天先简单看一些模板类和Vector的一些东西,装上vcpkg,理解它的工作逻辑。

c_cpp_properties配置对编译过程不影响,他只是一个提示插件的配置。只有在Tasks中的-I配置才会对实际建构过程产生影响,因此即使在前一个文件中设置了include外部库的路径,还是一样会报错。Ctrl+Shift+P会build task的任务,如果定义了多个Task,会接下来再进行选择。

模板类是先用尖括号给定参数生成一个类,然后后面就相当于用一个类去实例化一个对象了。

哦这个就像python的标准库math一样,但是不同Python库的命名空间都是不一样的,而C++的标准库统一放在了std这个命名空间之中。

在C++的标准库中,数学函数(如expsincos等)被定义在<cmath>头文件中,并位于std命名空间中。这是因为C++标准库中的函数和类型都位于std命名空间中,以避免潜在的命名冲突。

为了使用这些数学函数,我们需要在代码中显式指定命名空间std,例如std::exp。这样做是为了明确指定我们使用的是标准库中的函数,而不是其他可能存在的同名函数(可能来自其他库或用户自定义的函数)。

如果我们不使用std命名空间,而是省略前缀std::,编译器可能无法识别所使用的函数,或者会将其解释为其他非预期的符号。因此,为了编写符合标准的C++代码,我们应该使用std命名空间来访问标准库中的函数和类型。

总的来说,借助 GDB 调试器可以实现以下几个功能:

  1. 程序启动时,可以按照我们自定义的要求运行程序,例如设置参数和环境变量;
  2. 可使被调试程序在指定代码处暂停运行,并查看当前程序的运行状态(例如当前变量的值,函数的执行结果等),即支持断点调试;
  3. 程序执行过程中,可以改变某个变量的值,还可以改变代码的执行顺序,从而尝试修改程序中出现的逻辑错误。

Vscode中,继续运行就是继续运行直到遇到下一个断点,但是不会进入到里面的调用中。

  1. 单步调试(Step Into):单步调试允许你逐行执行代码,进入函数调用并跟踪函数内部的执行过程。当你在单步调试模式下,遇到函数调用时,你可以选择单步进入该函数并跟踪其内部执行。这对于逐步调试代码以了解每个语句的执行情况非常有用。

  2. 单步跳过(Step Over):单步跳过允许你跳过当前行的执行,并直接执行下一行。如果当前行包含函数调用,单步跳过将执行该函数调用并直接跳到函数调用的下一行。这在你对当前行不感兴趣,希望快速执行到下一行时很有用。

区别就是这一行还是会执行,但是不再跟踪这一行的调用内部的情况。

单步跳出主要是涉及嵌套调用的层之间的关系。

下面的属性是每个launch配置文件必须具备的:

  • type - 调试器的类型。
  • request - 请求类型,目前支持launch和attach
  • name - 在Debug启动配置下拉菜单中显示的名字

在C++中想要对STL动态数组做和Numpy一样的操作的时候,可以用STL提供的一些函数,比如transform,可以把一个数组中的每一个元素逐个输入到函数中,然后把返回的值赋值给另一个数组。accumulate则是为每个函数都执行一遍这个函数的操作,然后把函数的返回值分别累加起来。std::partial_sum还可以用来生成累积分布函数,

感觉有了Copilot和AI以后,我要做的就是把变量和函数的名字写的清楚一点。

和Python一样,STL有些操作是原位操作,不会返回值,而是直接修改传入到函数中的对象的数值。

C++的函数参数的配置和Python是一样的,支持默认参数,但是默认参数必须放到最后面。