Skip to content

Commit c869837

Browse files
committed
adding diam2vol to aero_data with tests
1 parent 6defa9b commit c869837

File tree

4 files changed

+58
-3
lines changed

4 files changed

+58
-3
lines changed

src/aero_data.F90

+12
Original file line numberDiff line numberDiff line change
@@ -149,4 +149,16 @@ subroutine f_aero_data_vol2rad(ptr_c, vol, radius) bind(C)
149149

150150
end subroutine
151151

152+
subroutine f_aero_data_diam2vol(ptr_c, diam, vol) bind(C)
153+
type(aero_data_t), pointer :: ptr_f => null()
154+
type(c_ptr), intent(in) :: ptr_c
155+
real(c_double), intent(in) :: diam
156+
real(c_double), intent(out) :: vol
157+
158+
call c_f_pointer(ptr_c, ptr_f)
159+
160+
vol = aero_data_diam2vol(ptr_f, diam)
161+
162+
end subroutine
163+
152164
end module

src/aero_data.hpp

+7
Original file line numberDiff line numberDiff line change
@@ -22,6 +22,7 @@ extern "C" void f_aero_data_set_prime_radius(void *ptr, const double*) noexcept;
2222
extern "C" void f_aero_data_get_prime_radius(const void *ptr, double*) noexcept;
2323
extern "C" void f_aero_data_rad2vol(const void *ptr, const double*, double*) noexcept;
2424
extern "C" void f_aero_data_vol2rad(const void *ptr, const double*, double*) noexcept;
25+
extern "C" void f_aero_data_diam2vol(const void *ptr, const double*, double*) noexcept;
2526

2627
struct AeroData {
2728
PMCResource ptr;
@@ -90,5 +91,11 @@ struct AeroData {
9091
return radius;
9192
}
9293

94+
static double diam2vol(const AeroData &self, const double diam) {
95+
double vol;
96+
f_aero_data_diam2vol(&self.ptr, &diam, &vol);
97+
return vol;
98+
}
99+
93100
};
94101

src/pypartmc.cpp

+4-2
Original file line numberDiff line numberDiff line change
@@ -72,9 +72,11 @@ PYBIND11_MODULE(_PyPartMC, m) {
7272
.def_property("vol_fill_factor", &AeroData::get_vol_fill_factor, &AeroData::set_vol_fill_factor)
7373
.def_property("prime_radius", &AeroData::get_prime_radius, &AeroData::set_prime_radius)
7474
.def("rad2vol", AeroData::rad2vol,
75-
"Convert geometric radius to mass-equivalent volume")
75+
"Convert geometric radius (m) to mass-equivalent volume (m^3).")
7676
.def("vol2rad", AeroData::vol2rad,
77-
"Convert mass-equivalent volume to geometric radius")
77+
"Convert mass-equivalent volume (m^3) to geometric radius (m)")
78+
.def("diam2vol", AeroData::diam2vol,
79+
"Convert geometric diameter (m) to mass-equivalent volume (m^3).")
7880
;
7981

8082
py::class_<AeroParticle>(m, "AeroParticle",

tests/test_aero_data.py

+35-1
Original file line numberDiff line numberDiff line change
@@ -150,7 +150,7 @@ def test_vol2rad_sphere():
150150
{"frac_dim": 2.5, "vol_fill_factor": 1.1, "prime_radius": 1e-8},
151151
{"frac_dim": 2.2, "vol_fill_factor": 1.3, 'prime_radius': 1e-6}))
152152
def test_vol2rad_fractal(aero_data_params:dict):
153-
#arrange
153+
# arrange
154154
sut = ppmc.AeroData(AERO_DATA_CTOR_ARG_MINIMAL)
155155
vol = 4.19e-18
156156
for key,value in aero_data_params.items():
@@ -163,3 +163,37 @@ def test_vol2rad_fractal(aero_data_params:dict):
163163
np.testing.assert_almost_equal(value, sut.prime_radius *
164164
(((3*vol/4/np.pi)**(1/3)/sut.prime_radius)**3 * sut.vol_fill_factor)**
165165
(1/sut.frac_dim))
166+
167+
@staticmethod
168+
def test_diam2vol_sphere():
169+
# arrange
170+
sut = ppmc.AeroData(AERO_DATA_CTOR_ARG_MINIMAL)
171+
diam = 2e-6
172+
sut.frac_dim = 3.0
173+
sut.vol_fill_factor = 1.0
174+
sut.prime_radius = 1e-8
175+
176+
# act
177+
vol = sut.diam2vol(diam)
178+
179+
# assert
180+
np.testing.assert_almost_equal(vol, (np.pi/6)*diam**3)
181+
182+
@staticmethod
183+
@pytest.mark.parametrize("aero_data_params", (
184+
{"frac_dim": 2.4, "vol_fill_factor": 1.2, 'prime_radius': 1e-7},
185+
{"frac_dim": 2.5, "vol_fill_factor": 1.1, "prime_radius": 1e-8},
186+
{"frac_dim": 2.2, "vol_fill_factor": 1.3, 'prime_radius': 1e-6}))
187+
def test_diam2vol_fractal(aero_data_params:dict):
188+
# arrange
189+
sut = ppmc.AeroData(AERO_DATA_CTOR_ARG_MINIMAL)
190+
diam = 2e-6
191+
for key,value in aero_data_params.items():
192+
setattr(sut, key, value)
193+
194+
# act
195+
value = sut.diam2vol(diam)
196+
197+
# assert
198+
np.testing.assert_almost_equal(value, (4/3)*np.pi*(sut.prime_radius)**3 *
199+
(1e-6/sut.prime_radius)**sut.frac_dim / sut.vol_fill_factor)

0 commit comments

Comments
 (0)