@@ -187,6 +187,240 @@ std::tuple<HighsStatus, HighsRanging> highs_getRanging(Highs* h) {
187187 return std::make_tuple (status, ranging);
188188}
189189
190+ std::tuple<HighsStatus, dense_array_t <HighsInt>>
191+ highs_getBasicVariables (Highs* h) {
192+ HighsInt num_row = h->getNumRow ();
193+
194+ HighsStatus status = HighsStatus::kOk ;
195+ std::vector<HighsInt> basic_variables (num_row);
196+ HighsInt* basic_variables_ptr = static_cast <HighsInt*>(basic_variables.data ());
197+ if (num_row > 0 ) status = h->getBasicVariables (basic_variables_ptr);
198+ return std::make_tuple (status, py::cast (basic_variables));
199+ }
200+
201+ std::tuple<HighsStatus, dense_array_t <double >>
202+ highs_getBasisInverseRow (Highs* h, HighsInt row) {
203+ HighsInt num_row = h->getNumRow ();
204+
205+ HighsStatus status = HighsStatus::kOk ;
206+ std::vector<double > solution_vector (num_row);
207+ double * solution_vector_ptr = static_cast <double *>(solution_vector.data ());
208+
209+ if (num_row > 0 ) status = h->getBasisInverseRow (row, solution_vector_ptr);
210+ return std::make_tuple (status, py::cast (solution_vector));
211+ }
212+
213+ std::tuple<HighsStatus, dense_array_t <double >, HighsInt, dense_array_t <HighsInt>>
214+ highs_getBasisInverseRowSparse (Highs* h, HighsInt row) {
215+ HighsInt num_row = h->getNumRow ();
216+
217+ HighsStatus status = HighsStatus::kOk ;
218+ HighsInt solution_num_nz = 0 ;
219+ std::vector<double > solution_vector (num_row);
220+ std::vector<HighsInt> solution_index (num_row);
221+ double * solution_vector_ptr = static_cast <double *>(solution_vector.data ());
222+ HighsInt* solution_index_ptr = static_cast <HighsInt*>(solution_index.data ());
223+
224+ if (num_row > 0 ) status = h->getBasisInverseRow (row, solution_vector_ptr, &solution_num_nz, solution_index_ptr);
225+ return std::make_tuple (status, py::cast (solution_vector), solution_num_nz, py::cast (solution_index));
226+ }
227+
228+ std::tuple<HighsStatus, dense_array_t <double >>
229+ highs_getBasisInverseCol (Highs* h, HighsInt col) {
230+ HighsInt num_row = h->getNumRow ();
231+
232+ HighsStatus status = HighsStatus::kOk ;
233+ std::vector<double > solution_vector (num_row);
234+ double * solution_vector_ptr = static_cast <double *>(solution_vector.data ());
235+
236+ if (num_row > 0 ) status = h->getBasisInverseCol (col, solution_vector_ptr);
237+ return std::make_tuple (status, py::cast (solution_vector));
238+ }
239+
240+ std::tuple<HighsStatus, dense_array_t <double >, HighsInt, dense_array_t <HighsInt>>
241+ highs_getBasisInverseColSparse (Highs* h, HighsInt col) {
242+ HighsInt num_row = h->getNumRow ();
243+
244+ HighsStatus status = HighsStatus::kOk ;
245+ HighsInt solution_num_nz = 0 ;
246+ std::vector<double > solution_vector (num_row);
247+ std::vector<HighsInt> solution_index (num_row);
248+ double * solution_vector_ptr = static_cast <double *>(solution_vector.data ());
249+ HighsInt* solution_index_ptr = static_cast <HighsInt*>(solution_index.data ());
250+
251+ if (num_row > 0 ) status = h->getBasisInverseCol (col, solution_vector_ptr, &solution_num_nz, solution_index_ptr);
252+ return std::make_tuple (status, py::cast (solution_vector), solution_num_nz, py::cast (solution_index));
253+ }
254+
255+ std::tuple<HighsStatus, dense_array_t <double >>
256+ highs_getBasisSolve (Highs* h, dense_array_t <double > rhs) {
257+ HighsInt num_row = h->getNumRow ();
258+
259+ py::buffer_info rhs_info = rhs.request ();
260+ double * rhs_ptr = static_cast <double *>(rhs_info.ptr );
261+
262+ HighsStatus status = HighsStatus::kOk ;
263+ std::vector<double > solution_vector (num_row);
264+ double * solution_vector_ptr = static_cast <double *>(solution_vector.data ());
265+
266+ if (num_row > 0 ) status = h->getBasisSolve (rhs_ptr, solution_vector_ptr);
267+ return std::make_tuple (status, py::cast (solution_vector));
268+ }
269+
270+ std::tuple<HighsStatus, dense_array_t <double >, HighsInt, dense_array_t <HighsInt>>
271+ highs_getBasisSolveSparse (Highs* h, dense_array_t <double > rhs) {
272+ HighsInt num_row = h->getNumRow ();
273+
274+ py::buffer_info rhs_info = rhs.request ();
275+ double * rhs_ptr = static_cast <double *>(rhs_info.ptr );
276+
277+ HighsStatus status = HighsStatus::kOk ;
278+ HighsInt solution_num_nz = 0 ;
279+ std::vector<double > solution_vector (num_row);
280+ std::vector<HighsInt> solution_index (num_row);
281+ double * solution_vector_ptr = static_cast <double *>(solution_vector.data ());
282+ HighsInt* solution_index_ptr = static_cast <HighsInt*>(solution_index.data ());
283+
284+ if (num_row > 0 ) status = h->getBasisSolve (rhs_ptr, solution_vector_ptr, &solution_num_nz, solution_index_ptr);
285+ return std::make_tuple (status, py::cast (solution_vector), solution_num_nz, py::cast (solution_index));
286+ }
287+
288+ std::tuple<HighsStatus, dense_array_t <double >> highs_getBasisTransposeSolve (Highs* h, dense_array_t <double > rhs) {
289+ HighsInt num_row = h->getNumRow ();
290+
291+ py::buffer_info rhs_info = rhs.request ();
292+ double * rhs_ptr = static_cast <double *>(rhs_info.ptr );
293+
294+ HighsStatus status = HighsStatus::kOk ;
295+ std::vector<double > solution_vector (num_row);
296+ double * solution_vector_ptr = static_cast <double *>(solution_vector.data ());
297+
298+ if (num_row > 0 ) status = h->getBasisTransposeSolve (rhs_ptr, solution_vector_ptr);
299+ return std::make_tuple (status, py::cast (solution_vector));
300+ }
301+
302+ std::tuple<HighsStatus, dense_array_t <double >, HighsInt, dense_array_t <HighsInt>>
303+ highs_getBasisTransposeSolveSparse (Highs* h, dense_array_t <double > rhs) {
304+ HighsInt num_row = h->getNumRow ();
305+
306+ py::buffer_info rhs_info = rhs.request ();
307+ double * rhs_ptr = static_cast <double *>(rhs_info.ptr );
308+
309+ HighsStatus status = HighsStatus::kOk ;
310+ HighsInt solution_num_nz = 0 ;
311+ std::vector<double > solution_vector (num_row);
312+ std::vector<HighsInt> solution_index (num_row);
313+ double * solution_vector_ptr = static_cast <double *>(solution_vector.data ());
314+ HighsInt* solution_index_ptr = static_cast <HighsInt*>(solution_index.data ());
315+
316+ if (num_row > 0 ) status = h->getBasisTransposeSolve (rhs_ptr, solution_vector_ptr, &solution_num_nz, solution_index_ptr);
317+ return std::make_tuple (status, py::cast (solution_vector), solution_num_nz, py::cast (solution_index));
318+ }
319+
320+ std::tuple<HighsStatus, dense_array_t <double >>
321+ highs_getReducedRow (Highs* h, HighsInt row) {
322+ HighsInt num_col = h->getNumCol ();
323+ HighsInt num_row = h->getNumRow ();
324+
325+ HighsStatus status = HighsStatus::kOk ;
326+ std::vector<double > solution_vector (num_col);
327+ double * solution_vector_ptr = static_cast <double *>(solution_vector.data ());
328+
329+ if (num_row > 0 ) status = h->getReducedRow (row, solution_vector_ptr);
330+ return std::make_tuple (status, py::cast (solution_vector));
331+ }
332+
333+ std::tuple<HighsStatus, dense_array_t <double >, HighsInt, dense_array_t <HighsInt>>
334+ highs_getReducedRowSparse (Highs* h, HighsInt row) {
335+ HighsInt num_col = h->getNumCol ();
336+ HighsInt num_row = h->getNumRow ();
337+
338+ HighsStatus status = HighsStatus::kOk ;
339+ HighsInt solution_num_nz = 0 ;
340+ std::vector<double > solution_vector (num_row);
341+ std::vector<HighsInt> solution_index (num_row);
342+ double * solution_vector_ptr = static_cast <double *>(solution_vector.data ());
343+ HighsInt* solution_index_ptr = static_cast <HighsInt*>(solution_index.data ());
344+
345+ if (num_row > 0 ) status = h->getReducedRow (row, solution_vector_ptr, &solution_num_nz, solution_index_ptr);
346+ return std::make_tuple (status, py::cast (solution_vector), solution_num_nz, py::cast (solution_index));
347+ }
348+
349+ std::tuple<HighsStatus, dense_array_t <double >>
350+ highs_getReducedColumn (Highs* h, HighsInt col) {
351+ HighsInt num_row = h->getNumRow ();
352+
353+ HighsStatus status = HighsStatus::kOk ;
354+ std::vector<double > solution_vector (num_row);
355+ double * solution_vector_ptr = static_cast <double *>(solution_vector.data ());
356+
357+ if (num_row > 0 ) status = h->getReducedColumn (col, solution_vector_ptr);
358+ return std::make_tuple (status, py::cast (solution_vector));
359+ }
360+
361+ std::tuple<HighsStatus, dense_array_t <double >, HighsInt, dense_array_t <HighsInt>>
362+ highs_getReducedColumnSparse (Highs* h, HighsInt col) {
363+ HighsInt num_row = h->getNumRow ();
364+
365+ HighsStatus status = HighsStatus::kOk ;
366+ HighsInt solution_num_nz = 0 ;
367+ std::vector<double > solution_vector (num_row);
368+ std::vector<HighsInt> solution_index (num_row);
369+ double * solution_vector_ptr = static_cast <double *>(solution_vector.data ());
370+ HighsInt* solution_index_ptr = static_cast <HighsInt*>(solution_index.data ());
371+
372+ if (num_row > 0 ) status = h->getReducedColumn (col, solution_vector_ptr, &solution_num_nz, solution_index_ptr);
373+ return std::make_tuple (status, py::cast (solution_vector), solution_num_nz, py::cast (solution_index));
374+ }
375+
376+ std::tuple<HighsStatus, bool > highs_getDualRayExist (Highs* h) {
377+ bool has_dual_ray;
378+ HighsStatus status = h->getDualRay (has_dual_ray);
379+ return std::make_tuple (status, has_dual_ray);
380+ }
381+
382+ std::tuple<HighsStatus, bool , dense_array_t <double >> highs_getDualRay (Highs* h) {
383+ HighsInt num_row = h->getNumRow ();
384+ bool has_dual_ray;
385+ HighsStatus status = HighsStatus::kOk ;
386+ std::vector<double > value (num_row);
387+ double * value_ptr = static_cast <double *>(value.data ());
388+ if (num_row > 0 ) status = h->getDualRay (has_dual_ray, value_ptr);
389+ return std::make_tuple (status, has_dual_ray, py::cast (value));
390+ }
391+
392+ std::tuple<HighsStatus, bool > highs_getDualUnboundednessDirectionExist (Highs* h) {
393+ bool has_dual_unboundedness_direction;
394+ HighsStatus status = h->getDualUnboundednessDirection (has_dual_unboundedness_direction);
395+ return std::make_tuple (status, has_dual_unboundedness_direction);
396+ }
397+
398+ std::tuple<HighsStatus, bool , dense_array_t <double >> highs_getDualUnboundednessDirection (Highs* h) {
399+ HighsInt num_col = h->getNumCol ();
400+ bool has_dual_unboundedness_direction;
401+ HighsStatus status = HighsStatus::kOk ;
402+ std::vector<double > value (num_col);
403+ double * value_ptr = static_cast <double *>(value.data ());
404+ if (num_col > 0 ) status = h->getDualUnboundednessDirection (has_dual_unboundedness_direction, value_ptr);
405+ return std::make_tuple (status, has_dual_unboundedness_direction, py::cast (value));
406+ }
407+
408+ std::tuple<HighsStatus, bool > highs_getPrimalRayExist (Highs* h) {
409+ bool has_primal_ray;
410+ HighsStatus status = h->getPrimalRay (has_primal_ray);
411+ return std::make_tuple (status, has_primal_ray);
412+ }
413+
414+ std::tuple<HighsStatus, bool , dense_array_t <double >> highs_getPrimalRay (Highs* h) {
415+ HighsInt num_col = h->getNumCol ();
416+ bool has_primal_ray;
417+ HighsStatus status = HighsStatus::kOk ;
418+ std::vector<double > value (num_col);
419+ double * value_ptr = static_cast <double *>(value.data ());
420+ if (num_col > 0 ) status = h->getPrimalRay (has_primal_ray, value_ptr);
421+ return std::make_tuple (status, has_primal_ray, py::cast (value));
422+ }
423+
190424HighsStatus highs_addRow (Highs* h, double lower, double upper,
191425 HighsInt num_new_nz, dense_array_t <HighsInt> indices,
192426 dense_array_t <double > values) {
@@ -1036,6 +1270,26 @@ PYBIND11_MODULE(_core, m, py::mod_gil_not_used()) {
10361270 .def (" getModelPresolveStatus" , &Highs::getModelPresolveStatus)
10371271 .def (" getRanging" , &highs_getRanging)
10381272 .def (" getObjectiveValue" , &Highs::getObjectiveValue)
1273+ .def (" getDualObjectiveValue" , &Highs::getDualObjectiveValue)
1274+ .def (" getBasicVariables" , &highs_getBasicVariables)
1275+ .def (" getBasisInverseRow" , &highs_getBasisInverseRow)
1276+ .def (" getBasisInverseRowSparse" , &highs_getBasisInverseRowSparse)
1277+ .def (" getBasisInverseCol" , &highs_getBasisInverseCol)
1278+ .def (" getBasisInverseColSparse" , &highs_getBasisInverseColSparse)
1279+ .def (" getBasisSolve" , &highs_getBasisSolve)
1280+ .def (" getBasisSolveSparse" , &highs_getBasisSolveSparse)
1281+ .def (" getBasisTransposeSolve" , &highs_getBasisTransposeSolve)
1282+ .def (" getBasisTransposeSolveSparse" , &highs_getBasisTransposeSolveSparse)
1283+ .def (" getReducedRow" , &highs_getReducedRow)
1284+ .def (" getReducedRowSparse" , &highs_getReducedRowSparse)
1285+ .def (" getReducedColumn" , &highs_getReducedColumn)
1286+ .def (" getReducedColumnSparse" , &highs_getReducedColumnSparse)
1287+ .def (" getDualRayExist" , &highs_getDualRayExist)
1288+ .def (" getDualRay" , &highs_getDualRay)
1289+ .def (" getDualUnboundednessDirectionExist" , &highs_getDualUnboundednessDirectionExist)
1290+ .def (" getDualUnboundednessDirection" , &highs_getDualUnboundednessDirection)
1291+ .def (" getPrimalRayExist" , &highs_getPrimalRayExist)
1292+ .def (" getPrimalRay" , &highs_getPrimalRay)
10391293 .def (" getNumCol" , &Highs::getNumCol)
10401294 .def (" getNumRow" , &Highs::getNumRow)
10411295 .def (" getNumNz" , &Highs::getNumNz)
0 commit comments