Skip to content

Commit 95434d8

Browse files
authored
Fix possible truncation of particle optical properties values (#491)
1 parent 9c70fe7 commit 95434d8

2 files changed

Lines changed: 30 additions & 55 deletions

File tree

src/aero_particle.hpp

Lines changed: 10 additions & 40 deletions
Original file line numberDiff line numberDiff line change
@@ -207,36 +207,18 @@ struct AeroParticle {
207207
}
208208

209209
static auto scatter_cross_sect(const AeroParticle &self) {
210-
int len = n_swbands;
211-
double val;
212-
f_aero_particle_scatter_cross_sect(
213-
self.ptr.f_arg(),
214-
&val,
215-
&len
216-
);
217-
return val;
210+
auto fn = f_aero_particle_scatter_cross_sect;
211+
return pypartmc::get_array_values_set_len(self, fn, n_swbands);
218212
}
219213

220214
static auto absorb_cross_sect(const AeroParticle &self) {
221-
int len = n_swbands;
222-
double val;
223-
f_aero_particle_absorb_cross_sect(
224-
self.ptr.f_arg(),
225-
&val,
226-
&len
227-
);
228-
return val;
215+
auto fn = f_aero_particle_absorb_cross_sect;
216+
return pypartmc::get_array_values_set_len(self, fn, n_swbands);
229217
}
230218

231219
static auto asymmetry(const AeroParticle &self) {
232-
int len = n_swbands;
233-
double val;
234-
f_aero_particle_asymmetry(
235-
self.ptr.f_arg(),
236-
&val,
237-
&len
238-
);
239-
return val;
220+
auto fn = f_aero_particle_asymmetry;
221+
return pypartmc::get_array_values_set_len(self, fn, n_swbands);
240222
}
241223

242224
static auto sources(const AeroParticle &self) {
@@ -262,25 +244,13 @@ struct AeroParticle {
262244
}
263245

264246
static auto refract_shell(const AeroParticle &self) {
265-
int len = n_swbands;
266-
std::complex<double> refract_shell;
267-
f_aero_particle_refract_shell(
268-
self.ptr.f_arg(),
269-
&refract_shell,
270-
&len
271-
);
272-
return refract_shell;
247+
auto fn = f_aero_particle_refract_shell;
248+
return pypartmc::get_array_values_set_len<std::complex<double>>(self, fn, n_swbands);
273249
}
274250

275251
static auto refract_core(const AeroParticle &self) {
276-
int len = n_swbands;
277-
std::complex<double> refract_core;
278-
f_aero_particle_refract_core(
279-
self.ptr.f_arg(),
280-
&refract_core,
281-
&len
282-
);
283-
return refract_core;
252+
auto fn = f_aero_particle_refract_core;
253+
return pypartmc::get_array_values_set_len<std::complex<double>>(self, fn, n_swbands);
284254
}
285255

286256
static auto get_weight_class(const AeroParticle &self) {

tests/test_aero_particle.py

Lines changed: 20 additions & 15 deletions
Original file line numberDiff line numberDiff line change
@@ -441,59 +441,64 @@ def test_absorb_cross_sect():
441441
sut = ppmc.AeroParticle(ppmc.AeroData(AERO_DATA_CTOR_ARG_MINIMAL), [44])
442442

443443
# act
444-
value = sut.absorb_cross_sect
444+
value = list(sut.absorb_cross_sect)
445445

446446
# assert
447-
assert value == 0
448-
assert isinstance(value, float)
447+
assert len(value) >= 1
448+
assert all(v == 0 for v in value)
449+
assert isinstance(value[0], float)
449450

450451
@staticmethod
451452
def test_scatter_cross_sect():
452453
# arrange
453454
sut = ppmc.AeroParticle(ppmc.AeroData(AERO_DATA_CTOR_ARG_MINIMAL), [44])
454455

455456
# act
456-
value = sut.scatter_cross_sect
457+
value = list(sut.scatter_cross_sect)
457458

458459
# assert
459-
assert value == 0
460-
assert isinstance(value, float)
460+
assert len(value) >= 1
461+
assert all(v == 0 for v in value)
462+
assert isinstance(value[0], float)
461463

462464
@staticmethod
463465
def test_asymmetry():
464466
# arrange
465467
sut = ppmc.AeroParticle(ppmc.AeroData(AERO_DATA_CTOR_ARG_MINIMAL), [44])
466468

467469
# act
468-
value = sut.asymmetry
470+
value = list(sut.asymmetry)
469471

470472
# assert
471-
assert value == 0
472-
assert isinstance(value, float)
473+
assert len(value) >= 1
474+
assert all(v == 0 for v in value)
475+
assert isinstance(value[0], float)
473476

474477
@staticmethod
475478
def test_refract_shell():
476479
# arrange
477480
sut = ppmc.AeroParticle(ppmc.AeroData(AERO_DATA_CTOR_ARG_MINIMAL), [44])
478481

479482
# act
480-
value = sut.refract_shell
483+
value = list(sut.refract_shell)
481484

482485
# assert
483-
assert value == 0 + 0j
484-
assert isinstance(value, complex)
486+
assert len(value) >= 1
487+
assert all(v == 0 + 0j for v in value)
488+
assert isinstance(value[0], complex)
485489

486490
@staticmethod
487491
def test_refract_core():
488492
# arrange
489493
sut = ppmc.AeroParticle(ppmc.AeroData(AERO_DATA_CTOR_ARG_MINIMAL), [44])
490494

491495
# act
492-
value = sut.refract_core
496+
value = list(sut.refract_core)
493497

494498
# assert
495-
assert value == 0 + 0j
496-
assert isinstance(value, complex)
499+
assert len(value) >= 1
500+
assert all(v == 0 + 0j for v in value)
501+
assert isinstance(value[0], complex)
497502

498503
@staticmethod
499504
def test_sources():

0 commit comments

Comments
 (0)