forked from open-atmos/PyPartMC
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathaero_data.hpp
298 lines (260 loc) · 8.46 KB
/
aero_data.hpp
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
/*##################################################################################################
# This file is a part of PyPartMC licensed under the GNU General Public License v3 (LICENSE file) #
# Copyright (C) 2022 University of Illinois Urbana-Champaign #
# Authors: https://github.com/open-atmos/PyPartMC/graphs/contributors #
##################################################################################################*/
#pragma once
#include "pmc_resource.hpp"
#include "json_resource.hpp"
#include "pybind11/stl.h"
#include "aero_data_parameters.hpp"
extern "C" void f_aero_data_ctor(void *ptr) noexcept;
extern "C" void f_aero_data_dtor(void *ptr) noexcept;
extern "C" void f_aero_data_from_json(const void *ptr) noexcept;
extern "C" void f_aero_data_spec_by_name(const void *ptr, int *value, const char *name_data, const int *name_size) noexcept;
extern "C" void f_aero_data_len(const void *ptr, int *len) noexcept;
extern "C" void f_aero_data_n_source(const void *ptr, int *len) noexcept;
extern "C" void f_aero_data_set_frac_dim(void *ptr, const double*) noexcept;
extern "C" void f_aero_data_get_frac_dim(const void *ptr, double*) noexcept;
extern "C" void f_aero_data_set_vol_fill_factor(void *ptr, const double*) noexcept;
extern "C" void f_aero_data_get_vol_fill_factor(const void *ptr, double*) noexcept;
extern "C" void f_aero_data_set_prime_radius(void *ptr, const double*) noexcept;
extern "C" void f_aero_data_get_prime_radius(const void *ptr, double*) noexcept;
extern "C" void f_aero_data_rad2vol(const void *ptr, const double*, double*) noexcept;
extern "C" void f_aero_data_vol2rad(const void *ptr, const double*, double*) noexcept;
extern "C" void f_aero_data_diam2vol(const void *ptr, const double*, double*) noexcept;
extern "C" void f_aero_data_vol2diam(const void *ptr, const double*, double*) noexcept;
extern "C" void f_aero_data_get_species_density(const void *ptr, const int *idx, double *val) noexcept;
extern "C" void f_aero_data_get_species_kappa(const void *ptr, const int *idx, double *val) noexcept;
extern "C" void f_aero_data_get_species_molecular_weight(const void *ptr, const int *idx, double *val) noexcept;
extern "C" void f_aero_data_source_name_by_index(const void *ptr, const int *i_source,
char *name_data) noexcept;
extern "C" void f_aero_data_spec_name_by_index(const void *ptr, const int *i_spec,
char *name_data) noexcept;
struct AeroData {
PMCResource ptr;
AeroData(const nlohmann::json &json) :
ptr(f_aero_data_ctor, f_aero_data_dtor)
{
if (!InputJSONResource::unique_keys(json))
throw std::runtime_error("Species names must be unique");
JSONResourceGuard<InputJSONResource> guard(json);
f_aero_data_from_json(this->ptr.f_arg());
}
AeroData() :
ptr(f_aero_data_ctor, f_aero_data_dtor)
{}
static auto spec_by_name(const AeroData &self, const std::string &name) {
int value;
const int name_size = name.size();
f_aero_data_spec_by_name(
self.ptr.f_arg(),
&value,
name.c_str(),
&name_size
);
if (value==0)
throw std::runtime_error("Element not found.");
return value - 1;
}
static std::size_t __len__(const AeroData &self) {
int len;
f_aero_data_len(
self.ptr.f_arg(),
&len
);
return len;
}
static void set_frac_dim(AeroData &self, const double value) {
f_aero_data_set_frac_dim(
self.ptr.f_arg_non_const(),
&value
);
}
static auto get_frac_dim(const AeroData &self) {
double value;
f_aero_data_get_frac_dim(
self.ptr.f_arg(),
&value
);
return value;
}
static void set_vol_fill_factor(AeroData &self, const double value) {
f_aero_data_set_vol_fill_factor(
self.ptr.f_arg_non_const(),
&value
);
}
static auto get_prime_radius(AeroData &self) {
double value;
f_aero_data_get_prime_radius(
self.ptr.f_arg(),
&value
);
return value;
}
static void set_prime_radius(AeroData &self, const double value) {
f_aero_data_set_prime_radius(
self.ptr.f_arg_non_const(),
&value
);
}
static auto get_vol_fill_factor(const AeroData &self) {
double value;
f_aero_data_get_vol_fill_factor(
self.ptr.f_arg(),
&value
);
return value;
}
static auto rad2vol(const AeroData &self, const double radius) {
double vol;
f_aero_data_rad2vol(
self.ptr.f_arg(),
&radius,
&vol
);
return vol;
}
static auto vol2rad(const AeroData &self, const double vol) {
double radius;
f_aero_data_vol2rad(
self.ptr.f_arg(),
&vol,
&radius
);
return radius;
}
static auto diam2vol(const AeroData &self, const double diam) {
double vol;
f_aero_data_diam2vol(
self.ptr.f_arg(),
&diam,
&vol
);
return vol;
}
static auto vol2diam(const AeroData &self, const double vol) {
double diam;
f_aero_data_vol2diam(
self.ptr.f_arg(),
&vol,
&diam
);
return diam;
}
static auto densities(const AeroData &self) {
int len;
f_aero_data_len(
self.ptr.f_arg(),
&len
);
std::valarray<double> data(len);
for (int idx = 0; idx < len; idx++) {
f_aero_data_get_species_density(
self.ptr.f_arg(),
&idx,
&data[idx]
);
}
return data;
}
static auto density(const AeroData &self, const std::string &name) {
int idx;
double data;
const int name_size = name.size();
f_aero_data_spec_by_name(
self.ptr.f_arg(),
&idx,
name.c_str(),
&name_size
);
if (idx==0)
throw std::runtime_error("Element not found.");
idx--;
f_aero_data_get_species_density(
self.ptr.f_arg(),
&idx,
&data
);
return data;
}
static std::size_t n_source(const AeroData &self) {
int len;
f_aero_data_n_source(
self.ptr.f_arg(),
&len
);
if (len == -1)
throw std::runtime_error("No sources defined.");
return len;
}
static auto sources(const AeroData &self) {
int len;
f_aero_data_n_source(
self.ptr.f_arg(),
&len
);
char name[AERO_SOURCE_NAME_LEN];
auto names = pybind11::tuple(len);
for (int idx = 0; idx < len; idx++) {
f_aero_data_source_name_by_index(
self.ptr.f_arg(),
&idx,
name
);
names[idx] = std::string(name);
}
return names;
}
static auto names(const AeroData &self) {
int len;
f_aero_data_len(
self.ptr.f_arg(),
&len
);
char name[AERO_NAME_LEN];
auto names = pybind11::tuple(len);
for (int idx = 0; idx < len; idx++) {
f_aero_data_spec_name_by_index(
self.ptr.f_arg(),
&idx,
name
);
names[idx] = std::string(name);
}
return names;
}
static auto kappa(const AeroData &self) {
int len;
f_aero_data_len(
self.ptr.f_arg(),
&len
);
std::valarray<double> data(len);
for (int idx = 0; idx < len; idx++) {
f_aero_data_get_species_kappa(
self.ptr.f_arg(),
&idx,
&data[idx]
);
}
return data;
}
static auto molecular_weights(const AeroData &self) {
int len;
f_aero_data_len(
self.ptr.f_arg(),
&len
);
std::valarray<double> data(len);
for (int idx = 0; idx < len; idx++) {
f_aero_data_get_species_molecular_weight(
self.ptr.f_arg(),
&idx,
&data[idx]
);
}
return data;
}
};