C Integration

This example demonstrates how to use fluids from C.

Source Code

  1#define PY_SSIZE_T_CLEAN
  2#include <Python.h>
  3#include <stdio.h>
  4#include <time.h>
  5#include <assert.h>
  6#include <math.h>
  7
  8// Helper function to initialize Python and import fluids
  9PyObject* init_fluids(void) {
 10    Py_Initialize();
 11    PyObject* module = PyImport_ImportModule("fluids");
 12    if (!module) {
 13        PyErr_Print();
 14        fprintf(stderr, "Failed to import fluids module\n");
 15        return NULL;
 16    }
 17    return module;
 18}
 19
 20// Helper to measure time in microseconds
 21long long microseconds_now(void) {
 22    struct timespec ts;
 23    clock_gettime(CLOCK_MONOTONIC, &ts);
 24    return (long long)ts.tv_sec * 1000000LL + ts.tv_nsec / 1000LL;
 25}
 26
 27void test_atmosphere(PyObject* fluids) {
 28    printf("\nTesting atmosphere at 5000m elevation:\n");
 29    
 30    // Create ATMOSPHERE_1976 instance
 31    PyObject* atm_class = PyObject_GetAttrString(fluids, "ATMOSPHERE_1976");
 32    if (!atm_class) {
 33        PyErr_Print();
 34        return;
 35    }
 36    
 37    PyObject* args = PyTuple_New(0);
 38    PyObject* kwargs = PyDict_New();
 39    PyDict_SetItemString(kwargs, "Z", PyFloat_FromDouble(5000));
 40    
 41    PyObject* atm = PyObject_Call(atm_class, args, kwargs);
 42    if (!atm) {
 43        PyErr_Print();
 44        Py_DECREF(args);
 45        Py_DECREF(kwargs);
 46        Py_DECREF(atm_class);
 47        return;
 48    }
 49    
 50    // Test various properties
 51    PyObject* temp = PyObject_GetAttrString(atm, "T");
 52    PyObject* pressure = PyObject_GetAttrString(atm, "P");
 53    PyObject* density = PyObject_GetAttrString(atm, "rho");
 54    PyObject* gravity = PyObject_GetAttrString(atm, "g");
 55    PyObject* viscosity = PyObject_GetAttrString(atm, "mu");
 56    PyObject* conductivity = PyObject_GetAttrString(atm, "k");
 57    PyObject* sonic = PyObject_GetAttrString(atm, "v_sonic");
 58    
 59    printf("✓ Temperature: %.4f\n", PyFloat_AsDouble(temp));
 60    printf("✓ Pressure: %.4f\n", PyFloat_AsDouble(pressure));
 61    printf("✓ Density: %.6f\n", PyFloat_AsDouble(density));
 62    printf("✓ Gravity: %.6f\n", PyFloat_AsDouble(gravity));
 63    printf("✓ Viscosity: %.6e\n", PyFloat_AsDouble(viscosity));
 64    printf("✓ Thermal conductivity: %.4f\n", PyFloat_AsDouble(conductivity));
 65    printf("✓ Sonic velocity: %.4f\n", PyFloat_AsDouble(sonic));
 66    
 67    // Test static gravity method
 68    PyObject* gravity_method = PyObject_GetAttrString(atm_class, "gravity");
 69    PyObject* gravity_args = PyTuple_New(0);
 70    PyObject* gravity_kwargs = PyDict_New();
 71    PyDict_SetItemString(gravity_kwargs, "Z", PyFloat_FromDouble(1E5));
 72    
 73    PyObject* high_gravity = PyObject_Call(gravity_method, gravity_args, gravity_kwargs);
 74    printf("✓ High altitude gravity: %.6f\n", PyFloat_AsDouble(high_gravity));
 75    
 76    // Cleanup
 77    Py_DECREF(temp);
 78    Py_DECREF(pressure);
 79    Py_DECREF(density);
 80    Py_DECREF(gravity);
 81    Py_DECREF(viscosity);
 82    Py_DECREF(conductivity);
 83    Py_DECREF(sonic);
 84    Py_DECREF(gravity_method);
 85    Py_DECREF(gravity_args);
 86    Py_DECREF(gravity_kwargs);
 87    Py_DECREF(high_gravity);
 88    Py_DECREF(atm);
 89    Py_DECREF(args);
 90    Py_DECREF(kwargs);
 91    Py_DECREF(atm_class);
 92}
 93
 94void test_expanded_reynolds(PyObject* fluids) {
 95    printf("\nTesting Reynolds number calculations:\n");
 96    
 97    // Get Reynolds function
 98    PyObject* reynolds_func = PyObject_GetAttrString(fluids, "Reynolds");
 99    if (!reynolds_func) {
100        PyErr_Print();
101        return;
102    }
103    
104    // Test with density and viscosity
105    PyObject* args1 = PyTuple_New(0);
106    PyObject* kwargs1 = PyDict_New();
107    PyDict_SetItemString(kwargs1, "V", PyFloat_FromDouble(2.5));
108    PyDict_SetItemString(kwargs1, "D", PyFloat_FromDouble(0.25));
109    PyDict_SetItemString(kwargs1, "rho", PyFloat_FromDouble(1.1613));
110    PyDict_SetItemString(kwargs1, "mu", PyFloat_FromDouble(1.9E-5));
111    
112    PyObject* re1 = PyObject_Call(reynolds_func, args1, kwargs1);
113    double re1_val = PyFloat_AsDouble(re1);
114    printf("✓ Re (with rho, mu): %.4f\n", re1_val);
115    assert(fabs(re1_val - 38200.6579) < 0.1);
116    
117    // Test with kinematic viscosity
118    PyObject* args2 = PyTuple_New(0);
119    PyObject* kwargs2 = PyDict_New();
120    PyDict_SetItemString(kwargs2, "V", PyFloat_FromDouble(2.5));
121    PyDict_SetItemString(kwargs2, "D", PyFloat_FromDouble(0.25));
122    PyDict_SetItemString(kwargs2, "nu", PyFloat_FromDouble(1.636e-05));
123    
124    PyObject* re2 = PyObject_Call(reynolds_func, args2, kwargs2);
125    double re2_val = PyFloat_AsDouble(re2);
126    printf("✓ Re (with nu): %.4f\n", re2_val);
127    assert(fabs(re2_val - 38202.934) < 0.1);
128    
129    // Cleanup
130    Py_DECREF(reynolds_func);
131    Py_DECREF(args1);
132    Py_DECREF(kwargs1);
133    Py_DECREF(re1);
134    Py_DECREF(args2);
135    Py_DECREF(kwargs2);
136    Py_DECREF(re2);
137}
138
139void test_psd(PyObject* fluids) {
140    printf("\nTesting particle size distribution functionality:\n");
141    
142    // Get particle_size_distribution module
143    PyObject* psd_module = PyObject_GetAttrString(fluids, "particle_size_distribution");
144    if (!psd_module) {
145        PyErr_Print();
146        return;
147    }
148    
149    // Create discrete PSD
150    PyObject* psd_class = PyObject_GetAttrString(psd_module, "ParticleSizeDistribution");
151    if (!psd_class) {
152        PyErr_Print();
153        Py_DECREF(psd_module);
154        return;
155    }
156    
157    // Create lists for ds and numbers
158    PyObject* ds_list = PyList_New(15);
159    double ds[] = {240, 360, 450, 562.5, 703, 878, 1097, 1371, 1713, 
160                   2141, 2676, 3345, 4181, 5226, 6532};
161    for (int i = 0; i < 15; i++) {
162        PyList_SetItem(ds_list, i, PyFloat_FromDouble(ds[i]));
163    }
164    
165    PyObject* numbers_list = PyList_New(14);
166    double numbers[] = {65, 119, 232, 410, 629, 849, 990, 981, 825, 
167                       579, 297, 111, 21, 1};
168    for (int i = 0; i < 14; i++) {
169        PyList_SetItem(numbers_list, i, PyFloat_FromDouble(numbers[i]));
170    }
171    
172    PyObject* args = PyTuple_New(0);
173    PyObject* kwargs = PyDict_New();
174    PyDict_SetItemString(kwargs, "ds", ds_list);
175    PyDict_SetItemString(kwargs, "fractions", numbers_list);
176    PyDict_SetItemString(kwargs, "order", PyLong_FromLong(0));
177    
178    PyObject* psd = PyObject_Call(psd_class, args, kwargs);
179    printf("✓ Created discrete PSD\n");
180    
181    // Test mean sizes
182    PyObject* mean_size_method = PyObject_GetAttrString(psd, "mean_size");
183    PyObject* mean_args = PyTuple_Pack(2, PyLong_FromLong(2), PyLong_FromLong(1));
184    PyObject* d21 = PyObject_Call(mean_size_method, mean_args, NULL);
185    printf("✓ Size-weighted mean diameter: %.4f\n", PyFloat_AsDouble(d21));
186    
187    PyObject* mean_args2 = PyTuple_Pack(2, PyLong_FromLong(1), PyLong_FromLong(0));
188    PyObject* d10 = PyObject_Call(mean_size_method, mean_args2, NULL);
189    printf("✓ Arithmetic mean diameter: %.4f\n", PyFloat_AsDouble(d10));
190    
191    // Test percentile calculations
192    PyObject* dn_method = PyObject_GetAttrString(psd, "dn");
193    PyObject* d10_args = PyTuple_Pack(1, PyFloat_FromDouble(0.1));
194    PyObject* d90_args = PyTuple_Pack(1, PyFloat_FromDouble(0.9));
195    
196    PyObject* d10_percentile = PyObject_Call(dn_method, d10_args, NULL);
197    PyObject* d90_percentile = PyObject_Call(dn_method, d90_args, NULL);
198    
199    printf("✓ D10: %.4f\n", PyFloat_AsDouble(d10_percentile));
200    printf("✓ D90: %.4f\n", PyFloat_AsDouble(d90_percentile));
201    
202    // Cleanup
203    Py_DECREF(psd_module);
204    Py_DECREF(psd_class);
205    Py_DECREF(ds_list);
206    Py_DECREF(numbers_list);
207    Py_DECREF(args);
208    Py_DECREF(kwargs);
209    Py_DECREF(psd);
210    Py_DECREF(mean_size_method);
211    Py_DECREF(mean_args);
212    Py_DECREF(d21);
213    Py_DECREF(mean_args2);
214    Py_DECREF(d10);
215    Py_DECREF(dn_method);
216    Py_DECREF(d10_args);
217    Py_DECREF(d90_args);
218    Py_DECREF(d10_percentile);
219    Py_DECREF(d90_percentile);
220}
221
222void test_fluids(PyObject* fluids) {
223    printf("Running fluids tests from C...\n");
224
225    // Get version
226    PyObject* version = PyObject_GetAttrString(fluids, "__version__");
227    if (version) {
228        const char* version_str = PyUnicode_AsUTF8(version);
229        printf("✓ Successfully imported fluids\n");
230        printf("✓ Fluids version: %s\n", version_str);
231        Py_DECREF(version);
232    }
233
234    // Test Reynolds number calculation
235    PyObject* reynolds_func = PyObject_GetAttrString(fluids, "Reynolds");
236    if (reynolds_func) {
237        PyObject* args = PyTuple_New(0);
238        PyObject* kwargs = PyDict_New();
239        PyDict_SetItemString(kwargs, "V", PyFloat_FromDouble(2.5));
240        PyDict_SetItemString(kwargs, "D", PyFloat_FromDouble(0.1));
241        PyDict_SetItemString(kwargs, "rho", PyFloat_FromDouble(1000));
242        PyDict_SetItemString(kwargs, "mu", PyFloat_FromDouble(0.001));
243        
244        PyObject* result = PyObject_Call(reynolds_func, args, kwargs);
245        if (result) {
246            double re = PyFloat_AsDouble(result);
247            printf("✓ Reynolds number calculation successful: %f\n", re);
248            assert(re > 0);
249            Py_DECREF(result);
250        }
251        
252        Py_DECREF(args);
253        Py_DECREF(kwargs);
254        Py_DECREF(reynolds_func);
255    }
256}
257
258void benchmark_fluids(PyObject* fluids) {
259    printf("\nRunning benchmarks:\n");
260    
261    // Get friction_factor function
262    PyObject* friction_func = PyObject_GetAttrString(fluids, "friction_factor");
263    if (!friction_func) {
264        PyErr_Print();
265        return;
266    }
267
268    // Get TANK class
269    PyObject* tank_class = PyObject_GetAttrString(fluids, "TANK");
270    if (!tank_class) {
271        PyErr_Print();
272        Py_DECREF(friction_func);
273        return;
274    }
275
276    // Benchmark friction_factor
277    printf("\nBenchmarking friction_factor:\n");
278    long long start = microseconds_now();
279    
280    for(int i = 0; i < 1000000; i++) {
281        PyObject* args = PyTuple_New(0);
282        PyObject* kwargs = PyDict_New();
283        PyDict_SetItemString(kwargs, "Re", PyFloat_FromDouble(1e5));
284        PyDict_SetItemString(kwargs, "eD", PyFloat_FromDouble(0.0001));
285        
286        PyObject* result = PyObject_Call(friction_func, args, kwargs);
287        
288        Py_DECREF(result);
289        Py_DECREF(args);
290        Py_DECREF(kwargs);
291    }
292    
293    long long duration = microseconds_now() - start;
294    printf("Time for 1e6 friction_factor calls: %lld microseconds\n", duration);
295    printf("Average time per call: %.6f microseconds\n", duration / 1000000.0);
296
297    // Benchmark TANK creation
298    printf("\nBenchmarking TANK creation:\n");
299    start = microseconds_now();
300    
301    for(int i = 0; i < 1000; i++) {
302        PyObject* args = PyTuple_New(0);
303        PyObject* kwargs = PyDict_New();
304        PyDict_SetItemString(kwargs, "L", PyFloat_FromDouble(3));
305        PyDict_SetItemString(kwargs, "D", PyFloat_FromDouble(5));
306        PyDict_SetItemString(kwargs, "horizontal", Py_False);
307        PyDict_SetItemString(kwargs, "sideA", PyUnicode_FromString("torispherical"));
308        PyDict_SetItemString(kwargs, "sideB", PyUnicode_FromString("torispherical"));
309        PyDict_SetItemString(kwargs, "sideA_f", PyFloat_FromDouble(1));
310        PyDict_SetItemString(kwargs, "sideA_k", PyFloat_FromDouble(0.1));
311        PyDict_SetItemString(kwargs, "sideB_f", PyFloat_FromDouble(1));
312        PyDict_SetItemString(kwargs, "sideB_k", PyFloat_FromDouble(0.1));
313        
314        PyObject* tank = PyObject_Call(tank_class, args, kwargs);
315        
316        Py_DECREF(tank);
317        Py_DECREF(args);
318        Py_DECREF(kwargs);
319    }
320    
321    duration = microseconds_now() - start;
322    printf("Average time per tank creation: %.6f microseconds\n", duration / 1000.0);
323
324    Py_DECREF(friction_func);
325    Py_DECREF(tank_class);
326}
327
328void test_tank(PyObject* fluids) {
329    printf("\nTesting tank calculations:\n");
330    
331    // Get TANK class
332    PyObject* tank_class = PyObject_GetAttrString(fluids, "TANK");
333    if (!tank_class) {
334        PyErr_Print();
335        return;
336    }
337    
338    // Test basic tank creation
339    PyObject* args1 = PyTuple_New(0);
340    PyObject* kwargs1 = PyDict_New();
341    PyDict_SetItemString(kwargs1, "V", PyFloat_FromDouble(10));
342    PyDict_SetItemString(kwargs1, "L_over_D", PyFloat_FromDouble(0.7));
343    PyDict_SetItemString(kwargs1, "sideB", PyUnicode_FromString("conical"));
344    PyDict_SetItemString(kwargs1, "horizontal", Py_False);
345    
346    PyObject* T1 = PyObject_Call(tank_class, args1, kwargs1);
347    if (!T1) {
348        PyErr_Print();
349        goto cleanup1;
350    }
351    
352    PyObject* length = PyObject_GetAttrString(T1, "L");
353    PyObject* diameter = PyObject_GetAttrString(T1, "D");
354    printf("✓ Tank length: %.6f\n", PyFloat_AsDouble(length));
355    printf("✓ Tank diameter: %.6f\n", PyFloat_AsDouble(diameter));
356    
357    // Test torispherical tank
358    PyObject* args2 = PyTuple_New(0);
359    PyObject* kwargs2 = PyDict_New();
360    PyDict_SetItemString(kwargs2, "L", PyFloat_FromDouble(3));
361    PyDict_SetItemString(kwargs2, "D", PyFloat_FromDouble(5));
362    PyDict_SetItemString(kwargs2, "horizontal", Py_False);
363    PyDict_SetItemString(kwargs2, "sideA", PyUnicode_FromString("torispherical"));
364    PyDict_SetItemString(kwargs2, "sideB", PyUnicode_FromString("torispherical"));
365    PyDict_SetItemString(kwargs2, "sideA_f", PyFloat_FromDouble(1));
366    PyDict_SetItemString(kwargs2, "sideA_k", PyFloat_FromDouble(0.1));
367    PyDict_SetItemString(kwargs2, "sideB_f", PyFloat_FromDouble(1));
368    PyDict_SetItemString(kwargs2, "sideB_k", PyFloat_FromDouble(0.1));
369    
370    PyObject* DIN = PyObject_Call(tank_class, args2, kwargs2);
371    if (!DIN) {
372        PyErr_Print();
373        goto cleanup2;
374    }
375    
376    // Get max height
377    PyObject* h_max = PyObject_GetAttrString(DIN, "h_max");
378    printf("✓ Tank max height: %.6f\n", PyFloat_AsDouble(h_max));
379    
380    // Test h_from_V method
381    PyObject* h_from_V = PyObject_GetAttrString(DIN, "h_from_V");
382    PyObject* h_args = PyTuple_Pack(1, PyFloat_FromDouble(40));
383    PyObject* h_result = PyObject_CallObject(h_from_V, h_args);
384    printf("✓ Height at V=40: %.6f\n", PyFloat_AsDouble(h_result));
385    
386    // Test V_from_h method
387    PyObject* V_from_h = PyObject_GetAttrString(DIN, "V_from_h");
388    PyObject* v_args = PyTuple_Pack(1, PyFloat_FromDouble(4.1));
389    PyObject* v_result = PyObject_CallObject(V_from_h, v_args);
390    printf("✓ Volume at h=4.1: %.5f\n", PyFloat_AsDouble(v_result));
391    
392    // Cleanup
393    Py_DECREF(v_result);
394    Py_DECREF(v_args);
395    Py_DECREF(V_from_h);
396    Py_DECREF(h_result);
397    Py_DECREF(h_args);
398    Py_DECREF(h_from_V);
399    Py_DECREF(h_max);
400    Py_DECREF(DIN);
401cleanup2:
402    Py_DECREF(args2);
403    Py_DECREF(kwargs2);
404    Py_DECREF(diameter);
405    Py_DECREF(length);
406    Py_DECREF(T1);
407cleanup1:
408    Py_DECREF(args1);
409    Py_DECREF(kwargs1);
410    Py_DECREF(tank_class);
411}
412
413int main() {
414    PyObject* fluids = init_fluids();
415    if (!fluids) {
416        return 1;
417    }
418
419    test_fluids(fluids);
420    test_atmosphere(fluids);
421    test_tank(fluids);
422    test_expanded_reynolds(fluids);
423    test_psd(fluids);
424    benchmark_fluids(fluids);
425
426    Py_DECREF(fluids);
427    Py_Finalize();
428    printf("\nAll tests completed!\n");
429    return 0;
430}        

Requirements

  • Python with fluids installed

  • Development python libraries

Usage Notes

  • The example demonstrates basic integration with fluids

  • No overhead (1.3 microsecond friction factor, 7 microsecond tank creation)