-
Notifications
You must be signed in to change notification settings - Fork 21
Expand file tree
/
Copy pathparticles.hpp
More file actions
314 lines (280 loc) · 13 KB
/
particles.hpp
File metadata and controls
314 lines (280 loc) · 13 KB
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
#pragma once
#include "extincl.hpp"
#include "opts.hpp"
#include "../common/output.hpp"
#include "opts_init.hpp"
#include "arrinfo.hpp"
#include "backend.hpp"
namespace libcloudphxx
{
namespace lgrngn
{
// to allow storing instances for multiple backends in one container/pointer
template <typename real_t>
struct particles_proto_t
{
// initialisation
virtual void init(
const arrinfo_t<real_t> th,
const arrinfo_t<real_t> rv,
const arrinfo_t<real_t> rhod,
const arrinfo_t<real_t> p = arrinfo_t<real_t>(),
const arrinfo_t<real_t> courant_x = arrinfo_t<real_t>(),
const arrinfo_t<real_t> courant_y = arrinfo_t<real_t>(),
const arrinfo_t<real_t> courant_z = arrinfo_t<real_t>(),
const std::map<enum common::chem::chem_species_t, const arrinfo_t<real_t> > ambient_chem = std::map<enum common::chem::chem_species_t, const arrinfo_t<real_t> >()
) {
assert(false);
}
// stuff that requires Eulerian component to wait
virtual void step_sync(
const opts_t<real_t> &,
arrinfo_t<real_t> th,
arrinfo_t<real_t> rv,
const arrinfo_t<real_t> rhod = arrinfo_t<real_t>(),
const arrinfo_t<real_t> courant_x = arrinfo_t<real_t>(),
const arrinfo_t<real_t> courant_y = arrinfo_t<real_t>(),
const arrinfo_t<real_t> courant_z = arrinfo_t<real_t>(),
const arrinfo_t<real_t> diss_rate = arrinfo_t<real_t>(), // TKE dissipation rate (epsilon)
std::map<enum common::chem::chem_species_t, arrinfo_t<real_t> > ambient_chem = std::map<enum common::chem::chem_species_t, arrinfo_t<real_t> >()
) {
assert(false);
}
virtual void sync_in(
arrinfo_t<real_t> th,
arrinfo_t<real_t> rv,
const arrinfo_t<real_t> rhod = arrinfo_t<real_t>(),
const arrinfo_t<real_t> courant_x = arrinfo_t<real_t>(),
const arrinfo_t<real_t> courant_y = arrinfo_t<real_t>(),
const arrinfo_t<real_t> courant_z = arrinfo_t<real_t>(),
const arrinfo_t<real_t> diss_rate = arrinfo_t<real_t>(),
std::map<enum common::chem::chem_species_t, arrinfo_t<real_t> > ambient_chem = std::map<enum common::chem::chem_species_t, arrinfo_t<real_t> >()
) {
assert(false);
}
virtual void step_cond(
const opts_t<real_t> &,
arrinfo_t<real_t> th,
arrinfo_t<real_t> rv,
std::map<enum common::chem::chem_species_t, arrinfo_t<real_t> > ambient_chem = std::map<enum common::chem::chem_species_t, arrinfo_t<real_t> >()
) {
assert(false);
}
// returns accumulated rainfall
virtual void step_async(
const opts_t<real_t> &
) {
assert(false);
}
// method for accessing super-droplet statistics
virtual void diag_sd_conc() { assert(false); }
virtual void diag_pressure() { assert(false); }
virtual void diag_temperature() { assert(false); }
virtual void diag_RH() { assert(false); }
virtual void diag_all() { assert(false); }
virtual void diag_rw_ge_rc() { assert(false); }
virtual void diag_RH_ge_Sc() { assert(false); }
virtual void diag_dry_rng(const real_t&, const real_t&) { assert(false); }
virtual void diag_wet_rng(const real_t&, const real_t&) { assert(false); }
virtual void diag_kappa_rng(const real_t&, const real_t&) { assert(false); }
// The 3 following functions are for consecutive selection of SDs.
// It allows the user to select SDs based on multiple characteristics, e.g. wet radius (0.5, 1) and kappa (0.1, 0.2):
// diag_wet_rng(0.5, 1); diag_kappa_rng_cons(0.1, 0.2);
// NOTE: the call with "cons" needs to be right after the previous call to diag_X_rng!
// otherwise some other function used in between can overwrite n_filtered used for selection of moments
// NOTE2: We cannot implement this as an argument to diag_X_rng, because we would like it to default to null
// and Boost Python does not work well with virtual member functions that have default arguments
virtual void diag_dry_rng_cons(const real_t&, const real_t&) { assert(false); }
virtual void diag_wet_rng_cons(const real_t&, const real_t&) { assert(false); }
virtual void diag_kappa_rng_cons(const real_t&, const real_t&) { assert(false); }
virtual void diag_dry_mom(const int&) { assert(false); }
virtual void diag_wet_mom(const int&) { assert(false); }
virtual void diag_wet_mass_dens(const real_t&, const real_t&) { assert(false); }
virtual void diag_chem(const enum common::chem::chem_species_t&) { assert(false); }
virtual void diag_precip_rate() { assert(false); }
virtual void diag_kappa_mom(const int&) { assert(false); }
virtual void diag_up_mom(const int&) { assert(false); }
virtual void diag_vp_mom(const int&) { assert(false); }
virtual void diag_wp_mom(const int&) { assert(false); }
virtual void diag_incloud_time_mom(const int&) { assert(false); } // requires opts_init.diag_incloud_time==true
virtual void diag_max_rw() { assert(false); }
virtual void diag_vel_div() { assert(false); }
virtual void diag_coal_tele_mass_flux() { assert(false); }
virtual std::map<libcloudphxx::common::output_t, real_t> diag_puddle() { assert(false); return std::map<libcloudphxx::common::output_t, real_t>(); }
virtual real_t *outbuf() { assert(false); return NULL; }
// storing a pointer to opts_init (e.g. for interrogatin about
// dimensions in Python bindings)
opts_init_t<real_t> *opts_init;
// virtual destructor
virtual ~particles_proto_t() {};
};
// prototype of what's implemented in the .tpp file
template <typename real_t, backend_t backend>
struct particles_t: particles_proto_t<real_t>
{
// initialisation
void init(
const arrinfo_t<real_t> th,
const arrinfo_t<real_t> rv,
const arrinfo_t<real_t> rhod,
const arrinfo_t<real_t> p,
const arrinfo_t<real_t> courant_x,
const arrinfo_t<real_t> courant_y,
const arrinfo_t<real_t> courant_z,
const std::map<enum common::chem::chem_species_t, const arrinfo_t<real_t> > ambient_chem
);
// time-stepping methods
void step_sync(
const opts_t<real_t> &,
arrinfo_t<real_t> th,
arrinfo_t<real_t> rv,
const arrinfo_t<real_t> rhod,
const arrinfo_t<real_t> courant_x,
const arrinfo_t<real_t> courant_y,
const arrinfo_t<real_t> courant_z,
const arrinfo_t<real_t> diss_rate,
std::map<enum common::chem::chem_species_t, arrinfo_t<real_t> > ambient_chem
);
void sync_in(
arrinfo_t<real_t> th,
arrinfo_t<real_t> rv,
const arrinfo_t<real_t> rhod,
const arrinfo_t<real_t> courant_x,
const arrinfo_t<real_t> courant_y,
const arrinfo_t<real_t> courant_z,
const arrinfo_t<real_t> diss_rate,
std::map<enum common::chem::chem_species_t, arrinfo_t<real_t> > ambient_chem
);
void step_cond(
const opts_t<real_t> &,
arrinfo_t<real_t> th,
arrinfo_t<real_t> rv,
std::map<enum common::chem::chem_species_t, arrinfo_t<real_t> > ambient_chem = std::map<enum common::chem::chem_species_t, arrinfo_t<real_t> >()
);
void step_async(
const opts_t<real_t> &
);
// diagnostic methods
void diag_sd_conc();
void diag_pressure();
void diag_temperature();
void diag_RH();
void diag_dry_rng(const real_t &r_mi, const real_t &r_mx);
void diag_wet_rng(const real_t &r_mi, const real_t &r_mx);
void diag_kappa_rng(const real_t &r_mi, const real_t &r_mx);
void diag_dry_rng_cons(const real_t &r_mi, const real_t &r_mx);
void diag_wet_rng_cons(const real_t &r_mi, const real_t &r_mx);
void diag_kappa_rng_cons(const real_t &r_mi, const real_t &r_mx);
void diag_dry_mom(const int &k);
void diag_wet_mom(const int &k);
void diag_kappa_mom(const int &k);
void diag_up_mom(const int&);
void diag_vp_mom(const int&);
void diag_wp_mom(const int&);
void diag_incloud_time_mom(const int &k);
void diag_wet_mass_dens(const real_t&, const real_t&);
void diag_chem(const enum common::chem::chem_species_t&);
void diag_rw_ge_rc();
void diag_RH_ge_Sc();
void diag_all();
void diag_precip_rate();
void diag_max_rw();
void diag_vel_div();
void diag_coal_tele_mass_flux();
std::map<libcloudphxx::common::output_t, real_t> diag_puddle();
real_t *outbuf();
struct impl;
std::unique_ptr<impl> pimpl;
// constructor
particles_t(opts_init_t<real_t> opts_init, int n_x_tot = 0); // NOTE: python bindings fail if more than one default value...
// declare destructor to delay it's definition until impl is defined
~particles_t();
// helper typedef
typedef particles_proto_t<real_t> parent_t;
};
// specialization for the multi_GPU backend
// the interface is the same as for other backends (above)
template <typename real_t>
struct particles_t<real_t, multi_CUDA>: particles_proto_t<real_t>
{
// initialisation
void init(
const arrinfo_t<real_t> th,
const arrinfo_t<real_t> rv,
const arrinfo_t<real_t> rhod,
const arrinfo_t<real_t> p,
const arrinfo_t<real_t> courant_x,
const arrinfo_t<real_t> courant_y,
const arrinfo_t<real_t> courant_z,
const std::map<enum common::chem::chem_species_t, const arrinfo_t<real_t> > ambient_chem
);
// time-stepping methods
void step_sync(
const opts_t<real_t> &,
arrinfo_t<real_t> th,
arrinfo_t<real_t> rv,
const arrinfo_t<real_t> rhod ,
const arrinfo_t<real_t> courant_x,
const arrinfo_t<real_t> courant_y,
const arrinfo_t<real_t> courant_z,
const arrinfo_t<real_t> diss_rate,
std::map<enum common::chem::chem_species_t, arrinfo_t<real_t> > ambient_chem
);
void sync_in(
arrinfo_t<real_t> th,
arrinfo_t<real_t> rv,
const arrinfo_t<real_t> rhod ,
const arrinfo_t<real_t> courant_x,
const arrinfo_t<real_t> courant_y,
const arrinfo_t<real_t> courant_z,
const arrinfo_t<real_t> diss_rate,
std::map<enum common::chem::chem_species_t, arrinfo_t<real_t> > ambient_chem
);
void step_cond(
const opts_t<real_t> &,
arrinfo_t<real_t> th,
arrinfo_t<real_t> rv,
std::map<enum common::chem::chem_species_t, arrinfo_t<real_t> > ambient_chem = std::map<enum common::chem::chem_species_t, arrinfo_t<real_t> >()
);
void step_async(
const opts_t<real_t> &
);
// diagnostic methods
void diag_sd_conc();
void diag_pressure();
void diag_temperature();
void diag_RH();
void diag_dry_rng(const real_t &r_mi, const real_t &r_mx);
void diag_wet_rng(const real_t &r_mi, const real_t &r_mx);
void diag_kappa_rng(const real_t &r_mi, const real_t &r_mx);
void diag_dry_rng_cons(const real_t &r_mi, const real_t &r_mx);
void diag_wet_rng_cons(const real_t &r_mi, const real_t &r_mx);
void diag_kappa_rng_cons(const real_t &r_mi, const real_t &r_mx);
void diag_dry_mom(const int &k);
void diag_wet_mom(const int &k);
void diag_kappa_mom(const int&);
void diag_up_mom(const int&);
void diag_vp_mom(const int&);
void diag_wp_mom(const int&);
void diag_incloud_time_mom(const int&);
void diag_wet_mass_dens(const real_t&, const real_t&);
real_t *outbuf();
void diag_chem(const enum common::chem::chem_species_t&);
void diag_rw_ge_rc();
void diag_RH_ge_Sc();
void diag_all();
void diag_precip_rate();
void diag_max_rw();
void diag_vel_div();
std::map<libcloudphxx::common::output_t, real_t> diag_puddle();
struct impl;
std::unique_ptr<impl> pimpl;
// constructors
particles_t(opts_init_t<real_t> opts_init);
// dtor
~particles_t();
// helper typedef
typedef particles_proto_t<real_t> parent_t;
};
};
};