Free cookie consent management tool by TermsFeed Policy Generator

source: branches/2925_AutoDiffForDynamicalModels/HeuristicLab.Problems.DynamicalSystemsModelling/3.3/CVODES.cs @ 16977

Last change on this file since 16977 was 16976, checked in by gkronber, 5 years ago

#2925: Update data files to have equidistant time steps. Reactivated CVODES solver.

File size: 49.2 KB
Line 
1using System;
2using System.Runtime.InteropServices;
3
4namespace HeuristicLab.Problems.DynamicalSystemsModelling {
5  public static class CVODES {
6
7    // CVODES types
8    public enum MultistepMethod : int { CV_ADAMS = 1, CV_BDF = 2 };
9    public enum NonlinearSolverIteration : int { CV_NEWTON = 1, CV_FUNCTIONAL = 2 };
10
11
12    /* itask */
13    public const int CV_NORMAL = 1;
14    public const int CV_ONE_STEP = 2;
15
16    /* ism */
17    public const int CV_SIMULTANEOUS = 1;
18    public const int CV_STAGGERED = 2;
19    public const int CV_STAGGERED1 = 3;
20
21    /* DQtype */
22    public const int CV_CENTERED = 1;
23    public const int CV_FORWARD = 2;
24
25    /* interp */
26    public const int CV_HERMITE = 1;
27    public const int CV_POLYNOMIAL = 2;
28
29    /*
30     * ----------------------------------------
31     * CVODES return flags
32     * ----------------------------------------
33     */
34
35    public const int CV_SUCCESS = 0;
36    public const int CV_TSTOP_RETURN = 1;
37    public const int CV_ROOT_RETURN = 2;
38
39    public const int CV_WARNING = 99;
40
41    public const int CV_TOO_MUCH_WORK = -1;
42    public const int CV_TOO_MUCH_ACC = -2;
43    public const int CV_ERR_FAILURE = -3;
44    public const int CV_CONV_FAILURE = -4;
45
46    public const int CV_LINIT_FAIL = -5;
47    public const int CV_LSETUP_FAIL = -6;
48    public const int CV_LSOLVE_FAIL = -7;
49    public const int CV_RHSFUNC_FAIL = -8;
50    public const int CV_FIRST_RHSFUNC_ERR = -9;
51    public const int CV_REPTD_RHSFUNC_ERR = -10;
52    public const int CV_UNREC_RHSFUNC_ERR = -11;
53    public const int CV_RTFUNC_FAIL = -12;
54    public const int CV_CONSTR_FAIL = -13;
55
56    public const int CV_MEM_FAIL = -20;
57    public const int CV_MEM_NULL = -21;
58    public const int CV_ILL_INPUT = -22;
59    public const int CV_NO_MALLOC = -23;
60    public const int CV_BAD_K = -24;
61    public const int CV_BAD_T = -25;
62    public const int CV_BAD_DKY = -26;
63    public const int CV_TOO_CLOSE = -27;
64
65    public const int CV_NO_QUAD = -30;
66    public const int CV_QRHSFUNC_FAIL = -31;
67    public const int CV_FIRST_QRHSFUNC_ERR = -32;
68    public const int CV_REPTD_QRHSFUNC_ERR = -33;
69    public const int CV_UNREC_QRHSFUNC_ERR = -34;
70
71    public const int CV_NO_SENS = -40;
72    public const int CV_SRHSFUNC_FAIL = -41;
73    public const int CV_FIRST_SRHSFUNC_ERR = -42;
74    public const int CV_REPTD_SRHSFUNC_ERR = -43;
75    public const int CV_UNREC_SRHSFUNC_ERR = -44;
76
77    public const int CV_BAD_IS = -45;
78
79    public const int CV_NO_QUADSENS = -50;
80    public const int CV_QSRHSFUNC_FAIL = -51;
81    public const int CV_FIRST_QSRHSFUNC_ERR = -52;
82    public const int CV_REPTD_QSRHSFUNC_ERR = -53;
83    public const int CV_UNREC_QSRHSFUNC_ERR = -54;
84
85    /*
86     * ----------------------------------------
87     * CVODEA return flags
88     * ----------------------------------------
89     */
90
91    public const int CV_NO_ADJ = -101;
92    public const int CV_NO_FWD = -102;
93    public const int CV_NO_BCK = -103;
94    public const int CV_BAD_TB0 = -104;
95    public const int CV_REIFWD_FAIL = -105;
96    public const int CV_FWD_FAIL = -106;
97
98
99    public const int CV_GETY_BADT = -107;
100
101    [UnmanagedFunctionPointer(CallingConvention.Cdecl)]
102    public delegate int CVRhsFunc(
103        double t, // realtype
104        IntPtr y, // N_Vector
105        IntPtr ydot, // N_Vector
106        IntPtr user_data
107      );
108
109    [UnmanagedFunctionPointer(CallingConvention.Cdecl)]
110    public delegate int CVRhsFuncB(
111        double t, // realtype
112        IntPtr y, // N_Vector
113        IntPtr yB, // N_Vector
114        IntPtr yBdot, // N_Vector
115        IntPtr user_data
116      );
117
118    [UnmanagedFunctionPointer(CallingConvention.Cdecl)]
119    public delegate int CVDlsJacFunc(
120      double t,
121      IntPtr y, // N_Vector
122      IntPtr fy, // N_Vector
123      IntPtr Jac, // SUNMatrix
124      IntPtr user_data,
125      IntPtr tmp1, // N_Vector
126      IntPtr tmp2, // N_Vector
127      IntPtr tmp3 // N_Vector
128      );
129
130    [UnmanagedFunctionPointer(CallingConvention.Cdecl)]
131    public delegate int CVDlsJacFuncB(
132      double t,
133      IntPtr y, // N_Vector
134      IntPtr yB, // N_Vector
135      IntPtr fyB, // N_Vector
136      IntPtr Jac, // SUNMatrix
137      IntPtr user_data,
138      IntPtr tmp1, // N_Vector
139      IntPtr tmp2, // N_Vector
140      IntPtr tmp3 // N_Vector
141    );
142
143    [UnmanagedFunctionPointer(CallingConvention.Cdecl)]
144    public delegate int CVQuadRhsFn(
145      double t,
146      IntPtr y, // N_Vector
147      IntPtr yQdot, // N_Vector
148      IntPtr user_data
149      );
150
151    [UnmanagedFunctionPointer(CallingConvention.Cdecl)]
152    public delegate int CVQuadRhsFnB(
153      double t,
154      IntPtr y, // N_Vector
155      IntPtr yB, // N_Vector
156      IntPtr qBdot, // N_Vector
157      IntPtr user_data
158      );
159
160
161    // to calculate sensitivities RHS for all equations at once
162    // must compute (∂f/∂y)s_i(t) + ∂f/∂p_i
163    [UnmanagedFunctionPointer(CallingConvention.Cdecl)]
164    public delegate int CVSensRhsFn(
165      int Ns,
166      double t,
167      IntPtr y, // N_Vector
168      IntPtr ydot, // N_Vector
169      IntPtr yS, // N_Vector* one vector for each parameter
170      IntPtr ySdot, // N_Vector* one vector for each parameter
171      IntPtr user_data,
172      IntPtr tmp1, // N_Vector
173      IntPtr tmp2 // N_Vector
174    );
175
176    [UnmanagedFunctionPointer(CallingConvention.Cdecl)]
177    public delegate void CVErrHandlerFn(
178     int errorCode,
179     IntPtr module, // const char *
180     IntPtr function,  // const char *
181     IntPtr msg, // char *
182     IntPtr ehdata
183    );
184
185    // returns a void* to the cvodes memory block if successful otherwise NULL
186    public static IntPtr CVodeCreate(MultistepMethod lmm, NonlinearSolverIteration iter) {
187      if (Environment.Is64BitProcess) return CVodeCreate_x64(lmm, iter);
188      else return CVodeCreate_x86(lmm, iter);
189    }
190
191    public static int CVodeInit(
192      IntPtr cvode_mem, // pointer returned by CVodeCreate
193      CVRhsFunc f,
194      double t0, // realtype, the inital value of t
195      IntPtr y0 // N_Vector the initial value of y
196    ) {
197      if (Environment.Is64BitProcess) return CVodeInit_x64(cvode_mem, f, t0, y0);
198      else return CVodeInit_x86(cvode_mem, f, t0, y0);
199    }
200
201
202    public static int CVodeSStolerances(
203      IntPtr cvode_mem,
204      double reltol,
205      double abstol
206      ) {
207      if (Environment.Is64BitProcess) return CVodeSStolerances_x64(cvode_mem, reltol, abstol);
208      else return CVodeSStolerances_x86(cvode_mem, reltol, abstol);
209    }
210
211    public static int CVDlsSetLinearSolver(
212      IntPtr cvode_mem,
213      IntPtr linearSolver, // SUNLinearSolver
214      IntPtr j // SUNMatrix
215      ) {
216      if (Environment.Is64BitProcess) return CVDlsSetLinearSolver_x64(cvode_mem, linearSolver, j);
217      else return CVDlsSetLinearSolver_x86(cvode_mem, linearSolver, j);
218    }
219
220    public static int CVDlsSetJacFn(
221      IntPtr cvode_mem,
222      CVDlsJacFunc jacFunc
223      ) {
224      if (Environment.Is64BitProcess) return CVDlsSetJacFn_x64(cvode_mem, jacFunc);
225      else return CVDlsSetJacFn_x86(cvode_mem, jacFunc);
226    }
227
228    public static int CVodeSensInit(
229      IntPtr cvode_mem,
230      int Ns, // number of parameters
231      int ism, // sensitivity solution method, CV_SIMULTANEOUS or CV_STAGGERED
232      CVSensRhsFn fS, // right hand side function which computes all sensitivity RHS at the same time
233      IntPtr yS0 // N_Vector
234      ) {
235      if (Environment.Is64BitProcess) return CVodeSensInit_x64(cvode_mem, Ns, ism, fS, yS0);
236      else return CVodeSensInit_x86(cvode_mem, Ns, ism, fS, yS0);
237    }
238
239    public static int CVodeSensSStolerances(
240      IntPtr cvode_mem,
241      double reltol,
242      double[] abstol
243    ) {
244      if (Environment.Is64BitProcess) return CVodeSensSStolerances_x64(cvode_mem, reltol, abstol);
245      else return CVodeSensSStolerances_x86(cvode_mem, reltol, abstol);
246    }
247
248    /* Call CVodeSensEEtolerances to estimate tolerances for sensitivity
249   variables based on the rolerances supplied for states variables and
250   the scaling factor pbar */
251    public static int CVodeSensEEtolerances(
252      IntPtr cvode_mem
253      ) {
254      if (Environment.Is64BitProcess) return CVodeSensEEtolerances_x64(cvode_mem);
255      else return CVodeSensEEtolerances_x86(cvode_mem);
256    }
257
258    public static int CVodeGetSens(
259      IntPtr cvode_mem,
260      ref double tret,
261      IntPtr yS //N_Vector*, one vector for each parameter
262      ) {
263      if (Environment.Is64BitProcess) return CVodeGetSens_x64(cvode_mem, ref tret, yS);
264      else return CVodeGetSens_x86(cvode_mem, ref tret, yS);
265    }
266
267    public static int CVodeQuadInit(
268      IntPtr cvode_mem,
269      CVQuadRhsFn qF,
270      IntPtr yQ0 // N_Vector, initial value of yQ
271      ) {
272      if (Environment.Is64BitProcess) return CVodeQuadInit_x64(cvode_mem, qF, yQ0);
273      else return CVodeQuadInit_x86(cvode_mem, qF, yQ0);
274    }
275
276    public static int CVodeQuadInitB(
277      IntPtr cvode_mem,
278      int indexB,
279      CVQuadRhsFnB rhsQB,
280      IntPtr yQB0 // N_Vector, initial value of yQB
281    ) {
282      if (Environment.Is64BitProcess) return CVodeQuadInitB_x64(cvode_mem, indexB, rhsQB, yQB0);
283      else return CVodeQuadInitB_x86(cvode_mem, indexB, rhsQB, yQB0);
284    }
285
286
287    public static int CVodeAdjInit(
288      IntPtr cvode_mem,
289      int nd,  // integration steps between checkpoints
290      int interpType // either CV_POLYNOMIAL or CV_HERMITE
291      ) {
292      if (Environment.Is64BitProcess) return CVodeAdjInit_x64(cvode_mem, nd, interpType);
293      else return CVodeAdjInit_x86(cvode_mem, nd, interpType);
294    }
295
296    public static int CVodeF(
297        IntPtr cvode_mem,
298        double t_out, // the next time at which a computed solution is desired
299        IntPtr y, // N_Vector, the solution vector y
300        ref double t_ret, // the time reached by the solver (output)
301        int itask, // CV_NORMAL or CV_ONE_STEP
302        ref int ncheck  // the number of internal checkpoints stored so far.
303      ) {
304      if (Environment.Is64BitProcess) return CVodeF_x64(cvode_mem, t_out, y, ref t_ret, itask, ref ncheck);
305      else return CVodeF_x86(cvode_mem, t_out, y, ref t_ret, itask, ref ncheck);
306    }
307
308    public static int CVodeGetNumSteps(
309      IntPtr cvode_mem,
310      ref long ncheck  // the number of steps taken
311    ) {
312      if (Environment.Is64BitProcess) return CVodeGetNumSteps_x64(cvode_mem, ref ncheck);
313      else return CVodeGetNumSteps_x86(cvode_mem, ref ncheck);
314    }
315
316    public static int CVodeCreateB(
317        IntPtr cvode_mem,
318        MultistepMethod lmmB,
319        NonlinearSolverIteration iterB,
320        ref int which
321      ) {
322      if (Environment.Is64BitProcess) return CVodeCreateB_x64(cvode_mem, lmmB, iterB, ref which);
323      else return CVodeCreateB_x86(cvode_mem, lmmB, iterB, ref which);
324    }
325
326    public static int CVodeInitB(
327        IntPtr cvode_mem,
328        int which,
329        CVRhsFuncB rightHandSide,
330        double tB0, // endpoint T where final conditions are provided (equal to endpoint of forward integration)
331        IntPtr yB0 // N_Vector inital value at t=tb0 of backward solution
332      ) {
333      if (Environment.Is64BitProcess) return CVodeInitB_x64(cvode_mem, which, rightHandSide, tB0, yB0);
334      else return CVodeInitB_x86(cvode_mem, which, rightHandSide, tB0, yB0);
335    }
336
337    public static int CVode(
338      IntPtr cvode_mem,
339      double tout, // next time at which a solution is desired
340      IntPtr yout, // N_Vector
341      ref double tret, // the time reached by the solver (output)
342      int itask // flag indicating the job of the solver for the next step.
343      ) {
344      if (Environment.Is64BitProcess) return CVode_x64(cvode_mem, tout, yout, ref tret, itask);
345      else return CVode_x86(cvode_mem, tout, yout, ref tret, itask);
346    }
347
348
349    public static int CVodeSStolerancesB(
350      IntPtr cvode_mem,
351      int which,
352      double reltol,
353      double abstol
354    ) {
355      if (Environment.Is64BitProcess) return CVodeSStolerancesB_x64(cvode_mem, which, reltol, abstol);
356      else return CVodeSStolerancesB_x86(cvode_mem, which, reltol, abstol);
357    }
358
359    public static int CVDlsSetLinearSolverB(
360      IntPtr cvode_mem,
361      int which,
362      IntPtr linearSolver, // SUNLinearSolver
363      IntPtr j // SUNMatrix
364    ) {
365      if (Environment.Is64BitProcess) return CVDlsSetLinearSolverB_x64(cvode_mem, which, linearSolver, j);
366      else return CVDlsSetLinearSolverB_x86(cvode_mem, which, linearSolver, j);
367    }
368
369    public static int CVDlsSetJacFnB(
370      IntPtr cvode_mem,
371      int which,
372      CVDlsJacFuncB jacFunc
373    ) {
374      if (Environment.Is64BitProcess) return CVDlsSetJacFnB_x64(cvode_mem, which, jacFunc);
375      else return CVDlsSetJacFnB_x86(cvode_mem, which, jacFunc);
376    }
377
378    public static int CVodeB(
379      IntPtr cvode_mem,
380      double tout, // next time at which a solution is desired
381      int itask // flag indicating the job of the solver for the next step.
382    ) {
383      if (Environment.Is64BitProcess) return CVodeB_x64(cvode_mem, tout, itask);
384      else return CVodeB_x86(cvode_mem, tout, itask);
385    }
386
387    public static int CVodeGetB(
388          IntPtr cvode_mem,
389          int which,
390          ref double tret,
391          IntPtr yB // the backward solution at time tret
392        ) {
393      if (Environment.Is64BitProcess) return CVodeGetB_x64(cvode_mem, which, ref tret, yB);
394      else return CVodeGetB_x86(cvode_mem, which, ref tret, yB);
395    }
396    public static int CVodeGetAdjY(
397          IntPtr cvode_mem,
398          double t,
399          IntPtr y // the forward solution y(t)
400        ) {
401      if (Environment.Is64BitProcess) return CVodeGetAdjY_x64(cvode_mem, t, y);
402      else return CVodeGetAdjY_x86(cvode_mem, t, y);
403    }
404
405    public static IntPtr CVodeGetAdjCVodeBmem(
406              IntPtr cvode_mem,
407              int which
408            ) {
409      if (Environment.Is64BitProcess) return CVodeGetAdjCVodeBmem_x64(cvode_mem, which);
410      else return CVodeGetAdjCVodeBmem_x86(cvode_mem, which);
411    }
412
413    public static int CVodeGetQuadB(
414        IntPtr cvode_mem,
415        int which,
416        ref double tret,
417        IntPtr yQB // N_Vector
418      ) {
419      if (Environment.Is64BitProcess) return CVodeGetQuadB_x64(cvode_mem, which, ref tret, yQB);
420      else return CVodeGetQuadB_x86(cvode_mem, which, ref tret, yQB);
421    }
422
423    public static int CVodeFree(ref IntPtr cvode_mem) {
424      if (Environment.Is64BitProcess) return CVodeFree_x64(ref cvode_mem);
425      else return CVodeFree_x86(ref cvode_mem);
426    }
427
428
429    internal static int CVodeSetErrHandlerFn(IntPtr cvode_mem, CVErrHandlerFn errorHandler, IntPtr data) {
430      if (Environment.Is64BitProcess) return CVodeSetErrHandlerFn_x64(cvode_mem, errorHandler, data);
431      else return CVodeSetErrHandlerFn_x86(cvode_mem, errorHandler, data);
432    }
433
434    #region matrix
435    public static IntPtr SUNDenseMatrix(long m, long n) {
436      if (Environment.Is64BitProcess) return SUNDenseMatrix_x64(m, n);
437      else return SUNDenseMatrix_x86(m, n);
438    }
439
440    public static void SUNMatDestroy(IntPtr A) {
441      if (Environment.Is64BitProcess) SUNMatDestroy_x64(A);
442      else SUNMatDestroy_x86(A);
443    }
444
445    public unsafe static double* SUNDenseMatrix_Data(IntPtr matrix) {
446      if (Environment.Is64BitProcess) return SUNDenseMatrix_Data_x64(matrix);
447      else return SUNDenseMatrix_Data_x86(matrix);
448    }
449
450    public static long SUNDenseMatrix_Cols(IntPtr matrix) {
451      if (Environment.Is64BitProcess) return SUNDenseMatrix_Cols_x64(matrix);
452      else return SUNDenseMatrix_Cols_x86(matrix);
453    }
454
455    public static long SUNDenseMatrix_Rows(IntPtr matrix) {
456      if (Environment.Is64BitProcess) return SUNDenseMatrix_Rows_x64(matrix);
457      else return SUNDenseMatrix_Rows_x86(matrix);
458    }
459
460    public unsafe static double SUNDenseMatrix_Get(IntPtr mat, long i, long j) {
461      if (Environment.Is64BitProcess) {
462        long M = SUNDenseMatrix_Rows_x64(mat);
463        // the(i, j)th element of A(with 0 <= i<M and 0 <= j<N) is given by(A->data)[j*M+i]
464        return SUNDenseMatrix_Data_x64(mat)[j * M + i];
465      } else {
466        long M = SUNDenseMatrix_Rows_x86(mat);
467        return SUNDenseMatrix_Data_x86(mat)[j * M + i];
468      }
469    }
470
471    public unsafe static void SUNDenseMatrix_Set(IntPtr mat, long i, long j, double val) {
472      if (Environment.Is64BitProcess) {
473        long M = SUNDenseMatrix_Rows_x64(mat);
474        // the(i, j)th element of A(with 0 <= i<M and 0 <= j<N) is given by(A->data)[j*M+i]
475        SUNDenseMatrix_Data_x64(mat)[j * M + i] = val;
476      } else {
477        long M = SUNDenseMatrix_Rows_x86(mat);
478        // the(i, j)th element of A(with 0 <= i<M and 0 <= j<N) is given by(A->data)[j*M+i]
479        SUNDenseMatrix_Data_x86(mat)[j * M + i] = val;
480      }
481    }
482
483    #endregion
484
485    #region linear solver
486    public static IntPtr SUNDenseLinearSolver(
487      IntPtr y, // N_Vector
488      IntPtr A // SUNMatrix
489      ) {
490      if (Environment.Is64BitProcess) return SUNDenseLinearSolver_x64(y, A);
491      else return SUNDenseLinearSolver_x86(y, A);
492    }
493
494    public static int SUNLinSolInitialize(IntPtr linearSolver) {
495      if (Environment.Is64BitProcess) return SUNLinSolInitialize_x64(linearSolver);
496      else return SUNLinSolInitialize_x86(linearSolver);
497    }
498
499    public static int SUNLinSolSetup(
500      IntPtr linearSolver,
501      IntPtr A // SUNMatrix
502      ) {
503      if (Environment.Is64BitProcess) return SUNLinSolSetup_x64(linearSolver, A);
504      else return SUNLinSolSetup_x86(linearSolver, A);
505    }
506
507    public static int SUNLinSolSolve(
508      IntPtr linearSolver,
509      IntPtr A, // SUNMatrix
510      IntPtr x, // N_Vector
511      IntPtr b, // N_Vector
512      double tol
513      ) {
514      if (Environment.Is64BitProcess) return SUNLinSolSolve_x64(linearSolver, A, x, b, tol);
515      else return SUNLinSolSolve_x86(linearSolver, A, x, b, tol);
516    }
517
518    public static int SUNLinSolFree(IntPtr linearSolver) {
519      if (Environment.Is64BitProcess) return SUNLinSolFree_x64(linearSolver);
520      else return SUNLinSolFree_x86(linearSolver);
521    }
522
523    #endregion
524
525    #region N_Vector
526    public static IntPtr N_VNew_Serial(long vec_length) {
527      if (Environment.Is64BitProcess) return N_VNew_Serial_x64(vec_length);
528      else return N_VNew_Serial_x86(vec_length);
529    }
530
531    public static void N_VDestroy_Serial(IntPtr vec) {
532      if (Environment.Is64BitProcess) N_VDestroy_Serial_x64(vec);
533      else N_VDestroy_Serial_x86(vec);
534    }
535
536    public static void N_VDestroyVectorArray_Serial(IntPtr vecArr, int count) {
537      if (Environment.Is64BitProcess) N_VDestroyVectorArray_Serial_x64(vecArr, count);
538      else N_VDestroyVectorArray_Serial_x86(vecArr, count);
539    }
540
541    public static void N_VPrint_Serial(IntPtr vec) {
542      if (Environment.Is64BitProcess) N_VPrint_Serial_x64(vec);
543      else N_VPrint_Serial_x86(vec);
544    }
545
546    public static void N_VConst_Serial(double c, IntPtr vec) {
547      if (Environment.Is64BitProcess) N_VConst_Serial_x64(c, vec);
548      else N_VConst_Serial_x86(c, vec);
549    }
550
551
552    public static double N_VL1Norm_Serial(IntPtr vec) {
553      if (Environment.Is64BitProcess) return N_VL1Norm_Serial_x64(vec);
554      else return N_VL1Norm_Serial_x86(vec);
555    }
556
557    public static IntPtr N_VMake_Serial(long vec_length, double[] v_data) {
558      if (Environment.Is64BitProcess) return N_VMake_Serial_x64(vec_length, v_data);
559      else return N_VMake_Serial_x86(vec_length, v_data);
560    }
561
562    ///  Performs the operation z = c*x
563    public static void N_VScale(double s,
564      IntPtr x, // N_Vector
565      IntPtr z // N_Vector
566      ) {
567      if (Environment.Is64BitProcess) N_VScale_x64(s, x, z);
568      else N_VScale_x86(s, x, z);
569    }
570
571    public unsafe static double* N_VGetArrayPointer_Serial(IntPtr vec) {
572      if (Environment.Is64BitProcess) return N_VGetArrayPointer_Serial_x64(vec);
573      else return N_VGetArrayPointer_Serial_x86(vec);
574    }
575
576    public static IntPtr N_VCloneVectorArray_Serial(int count, IntPtr vec) {
577      if (Environment.Is64BitProcess) return N_VCloneVectorArray_Serial_x64(count, vec);
578      else return N_VCloneVectorArray_Serial_x86(count, vec);
579    }
580
581
582    /*
583#define NV_CONTENT_S(v)  ( (N_VectorContent_Serial)(v->content) )
584
585#define NV_LENGTH_S(v)   ( NV_CONTENT_S(v)->length )
586
587#define NV_OWN_DATA_S(v) ( NV_CONTENT_S(v)->own_data )
588
589#define NV_DATA_S(v)     ( NV_CONTENT_S(v)->data )
590
591#define NV_Ith_S(v,i)    ( NV_DATA_S(v)[i] )
592*/
593    // methods for macros
594    // the vector is a pointer to a struct
595    // the struct has following layout
596    // {
597    //    long length,
598    //    int  own_data,   // flag to store whether the memory for the vector is managed internally
599    //    double* data,    // pointer to the actual elements
600    // }
601
602    private unsafe static IntPtr NV_CONTENT_S(IntPtr v) {
603      if (Environment.Is64BitProcess) {
604        return new IntPtr(*(Int64*)v.ToPointer());
605      } else {
606        return new IntPtr(*(Int32*)v.ToPointer());
607      }
608    }
609
610    public unsafe static long NV_LENGTH_S(IntPtr v) {
611      long length = *(long*)NV_CONTENT_S(v).ToPointer();
612      return length;
613    }
614    public unsafe static bool NV_OWN_DATA_S(IntPtr v) {
615      var content = (int*)NV_CONTENT_S(v).ToPointer();
616      int own_data = *(content + 2);
617      return own_data > 0;
618    }
619    public unsafe static double* NV_DATA_S(IntPtr v) {
620      if (Environment.Is64BitProcess) {
621        var content = (long*)NV_CONTENT_S(v).ToPointer();
622        double* data = (double*)*(content + 2); // different alignment or sizeof(int)?
623        return data;
624      } else {
625        var content = (int*)NV_CONTENT_S(v).ToPointer();
626        double* data = (double*)*(content + 3);
627        return data;
628      }
629    }
630    public unsafe static double NV_Get_Ith_S(IntPtr v, long i) {
631      return NV_DATA_S(v)[i];
632    }
633    public unsafe static void NV_Set_Ith_S(IntPtr v, long i, double val) {
634      NV_DATA_S(v)[i] = val;
635    }
636    #endregion
637
638
639    #region x86
640
641    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeCreate", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
642    // returns a void* to the cvodes memory block if successful otherwise NULL
643    private static extern IntPtr CVodeCreate_x86(MultistepMethod lmm, NonlinearSolverIteration iter);
644
645    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeInit", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
646    private static extern int CVodeInit_x86(
647      IntPtr cvode_mem, // pointer returned by CVodeCreate
648      CVRhsFunc f,
649      double t0, // realtype, the inital value of t
650      IntPtr y0 // N_Vector the initial value of y
651    );
652
653
654    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeSStolerances", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
655    private static extern int CVodeSStolerances_x86(
656      IntPtr cvode_mem,
657      double reltol,
658      double abstol
659      );
660
661    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVDlsSetLinearSolver", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
662    private static extern int CVDlsSetLinearSolver_x86(
663      IntPtr cvode_mem,
664      IntPtr linearSolver, // SUNLinearSolver
665      IntPtr j // SUNMatrix
666      );
667
668    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVDlsSetJacFn", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
669    private static extern int CVDlsSetJacFn_x86(
670      IntPtr cvode_mem,
671      CVDlsJacFunc jacFunc
672      );
673
674    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeSensInit", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
675    private static extern int CVodeSensInit_x86(
676      IntPtr cvode_mem,
677      int Ns, // number of parameters
678      int ism, // sensitivity solution method, CV_SIMULTANEOUS or CV_STAGGERED
679      CVSensRhsFn fS, // right hand side function which computes all sensitivity RHS at the same time
680      IntPtr yS0 // N_Vector
681      );
682
683    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeSensSStolerances", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
684    private static extern int CVodeSensSStolerances_x86(
685      IntPtr cvode_mem,
686      double reltol,
687      double[] abstol
688    );
689
690    /* Call CVodeSensEEtolerances to estimate tolerances for sensitivity
691   variables based on the rolerances supplied for states variables and
692   the scaling factor pbar */
693    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeSensEEtolerances", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
694    private static extern int CVodeSensEEtolerances_x86(
695      IntPtr cvode_mem
696      );
697
698    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeGetSens", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
699    private static extern int CVodeGetSens_x86(
700      IntPtr cvode_mem,
701      ref double tret,
702      IntPtr yS //N_Vector*, one vector for each parameter
703      );
704
705    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeQuadInit", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
706    private static extern int CVodeQuadInit_x86(
707      IntPtr cvode_mem,
708      CVQuadRhsFn qF,
709      IntPtr yQ0 // N_Vector, initial value of yQ
710      );
711
712    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeQuadInitB", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
713    private static extern int CVodeQuadInitB_x86(
714      IntPtr cvode_mem,
715      int indexB,
716      CVQuadRhsFnB rhsQB,
717      IntPtr yQB0 // N_Vector, initial value of yQB
718    );
719
720
721    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeAdjInit", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
722    private static extern int CVodeAdjInit_x86(
723      IntPtr cvode_mem,
724      int nd,  // integration steps between checkpoints
725      int interpType // either CV_POLYNOMIAL or CV_HERMITE
726      );
727
728    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeF", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
729    private static extern int CVodeF_x86(
730        IntPtr cvode_mem,
731        double t_out, // the next time at which a computed solution is desired
732        IntPtr y, // N_Vector, the solution vector y
733        ref double t_ret, // the time reached by the solver (output)
734        int itask, // CV_NORMAL or CV_ONE_STEP
735        ref int ncheck  // the number of internal checkpoints stored so far.
736      );
737
738    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeGetNumSteps", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
739    private static extern int CVodeGetNumSteps_x86(
740      IntPtr cvode_mem,
741      ref long ncheck  // the number of steps taken
742    );
743
744    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeCreateB", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
745    private static extern int CVodeCreateB_x86(
746        IntPtr cvode_mem,
747        MultistepMethod lmmB,
748        NonlinearSolverIteration iterB,
749        ref int which
750      );
751
752    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeInitB", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
753    private static extern int CVodeInitB_x86(
754        IntPtr cvode_mem,
755        int which,
756        CVRhsFuncB rightHandSide,
757        double tB0, // endpoint T where final conditions are provided (equal to endpoint of forward integration)
758        IntPtr yB0 // N_Vector inital value at t=tb0 of backward solution
759      );
760
761    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVode", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
762    private static extern int CVode_x86(
763      IntPtr cvode_mem,
764      double tout, // next time at which a solution is desired
765      IntPtr yout, // N_Vector
766      ref double tret, // the time reached by the solver (output)
767      int itask // flag indicating the job of the solver for the next step.
768      );
769
770
771    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeSStolerancesB", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
772    private static extern int CVodeSStolerancesB_x86(
773      IntPtr cvode_mem,
774      int which,
775      double reltol,
776      double abstol
777    );
778
779    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVDlsSetLinearSolverB", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
780    private static extern int CVDlsSetLinearSolverB_x86(
781      IntPtr cvode_mem,
782      int which,
783      IntPtr linearSolver, // SUNLinearSolver
784      IntPtr j // SUNMatrix
785    );
786
787    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVDlsSetJacFnB", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
788    private static extern int CVDlsSetJacFnB_x86(
789      IntPtr cvode_mem,
790      int which,
791      CVDlsJacFuncB jacFunc
792    );
793
794    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeB", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
795    private static extern int CVodeB_x86(
796      IntPtr cvode_mem,
797      double tout, // next time at which a solution is desired
798      int itask // flag indicating the job of the solver for the next step.
799    );
800
801    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeGetB", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
802    private static extern int CVodeGetB_x86(
803          IntPtr cvode_mem,
804          int which,
805          ref double tret,
806          IntPtr yB // the backward solution at time tret
807        );
808    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeGetAdjY", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
809    private static extern int CVodeGetAdjY_x86(
810          IntPtr cvode_mem,
811          double t,
812          IntPtr y // the forward solution y(t)
813        );
814
815    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeGetAdjCVodeBmem", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
816    private static extern IntPtr CVodeGetAdjCVodeBmem_x86(
817              IntPtr cvode_mem,
818              int which
819            );
820
821    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeGetQuadB", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
822    private static extern int CVodeGetQuadB_x86(
823        IntPtr cvode_mem,
824        int which,
825        ref double tret,
826        IntPtr yQB // N_Vector
827      );
828
829    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeSetErrHandlerFn", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
830    private static extern int CVodeSetErrHandlerFn_x86(IntPtr cvode_mem, CVErrHandlerFn errorHandler, IntPtr data);
831
832    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeFree", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
833    private static extern int CVodeFree_x86(ref IntPtr cvode_mem);
834
835    #region matrix
836    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "SUNDenseMatrix", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
837    private static extern IntPtr SUNDenseMatrix_x86(long m, long n);
838
839    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "SUNMatDestroy", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
840    private static extern void SUNMatDestroy_x86(IntPtr A);
841
842    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "SUNDenseMatrix_Data", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
843    private unsafe static extern double* SUNDenseMatrix_Data_x86(IntPtr matrix);
844
845    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "SUNDenseMatrix_Cols", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
846    private static extern long SUNDenseMatrix_Cols_x86(IntPtr matrix);
847
848    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "SUNDenseMatrix_Rows", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
849    private static extern long SUNDenseMatrix_Rows_x86(IntPtr matrix);
850
851
852    #endregion
853
854    #region linear solver
855    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "SUNDenseLinearSolver", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
856    private static extern IntPtr SUNDenseLinearSolver_x86(
857      IntPtr y, // N_Vector
858      IntPtr A // SUNMatrix
859      );
860
861    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "SUNLinSolInitialize", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
862    private static extern int SUNLinSolInitialize_x86(IntPtr linearSolver);
863
864    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "SUNLinSolSetup", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
865    private static extern int SUNLinSolSetup_x86(
866      IntPtr linearSolver,
867      IntPtr A // SUNMatrix
868      );
869
870    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "SUNLinSolSolve", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
871    private static extern int SUNLinSolSolve_x86(
872      IntPtr linearSolver,
873      IntPtr A, // SUNMatrix
874      IntPtr x, // N_Vector
875      IntPtr b, // N_Vector
876      double tol
877      );
878
879    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "SUNLinSolFree", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
880    private static extern int SUNLinSolFree_x86(IntPtr linearSolver);
881
882    #endregion
883
884    #region N_Vector
885    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "N_VNew_Serial", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
886    private static extern IntPtr N_VNew_Serial_x86(long vec_length);
887
888    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "N_VDestroy_Serial", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
889    private static extern void N_VDestroy_Serial_x86(IntPtr vec);
890
891    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "N_VDestroyVectorArray_Serial", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
892    private static extern void N_VDestroyVectorArray_Serial_x86(IntPtr vecArr, int count);
893
894    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "N_VPrint_Serial", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
895    private static extern void N_VPrint_Serial_x86(IntPtr vec);
896
897    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "N_VConst_Serial", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
898    private static extern void N_VConst_Serial_x86(double c, IntPtr vec);
899
900
901    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "N_VL1Norm_Serial", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
902    private static extern double N_VL1Norm_Serial_x86(IntPtr vec);
903
904    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "N_VMake_Serial", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
905    private static extern IntPtr N_VMake_Serial_x86(long vec_length, double[] v_data);
906
907    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "N_VScale", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
908    ///  Performs the operation z = c*x
909    private static extern void N_VScale_x86(double s,
910      IntPtr x, // N_Vector
911      IntPtr z // N_Vector
912      );
913
914    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "N_VMake_Serial", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
915    private unsafe static extern double* N_VGetArrayPointer_Serial_x86(IntPtr vec);
916
917    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "N_VCloneVectorArray_Serial", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
918    private static extern IntPtr N_VCloneVectorArray_Serial_x86(int count, IntPtr vec); // returns N_Vector* !
919
920    #endregion
921
922    #endregion
923
924    #region x64
925
926    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "CVodeCreate", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
927    // returns a void* to the cvodes memory block if successful otherwise NULL
928    private static extern IntPtr CVodeCreate_x64(MultistepMethod lmm, NonlinearSolverIteration iter);
929
930    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "CVodeInit", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
931    private static extern int CVodeInit_x64(
932      IntPtr cvode_mem, // pointer returned by CVodeCreate
933      CVRhsFunc f,
934      double t0, // realtype, the inital value of t
935      IntPtr y0 // N_Vector the initial value of y
936    );
937
938
939    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "CVodeSStolerances", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
940    private static extern int CVodeSStolerances_x64(
941      IntPtr cvode_mem,
942      double reltol,
943      double abstol
944      );
945
946    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "CVDlsSetLinearSolver", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
947    private static extern int CVDlsSetLinearSolver_x64(
948      IntPtr cvode_mem,
949      IntPtr linearSolver, // SUNLinearSolver
950      IntPtr j // SUNMatrix
951      );
952
953    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "CVDlsSetJacFn", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
954    private static extern int CVDlsSetJacFn_x64(
955      IntPtr cvode_mem,
956      CVDlsJacFunc jacFunc
957      );
958
959    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "CVodeSensInit", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
960    private static extern int CVodeSensInit_x64(
961      IntPtr cvode_mem,
962      int Ns, // number of parameters
963      int ism, // sensitivity solution method, CV_SIMULTANEOUS or CV_STAGGERED
964      CVSensRhsFn fS, // right hand side function which computes all sensitivity RHS at the same time
965      IntPtr yS0 // N_Vector
966      );
967
968    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "CVodeSensSStolerances", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
969    private static extern int CVodeSensSStolerances_x64(
970      IntPtr cvode_mem,
971      double reltol,
972      double[] abstol
973    );
974
975    /* Call CVodeSensEEtolerances to estimate tolerances for sensitivity
976   variables based on the rolerances supplied for states variables and
977   the scaling factor pbar */
978    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "CVodeSensEEtolerances", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
979    private static extern int CVodeSensEEtolerances_x64(
980      IntPtr cvode_mem
981      );
982
983    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "CVodeGetSens", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
984    private static extern int CVodeGetSens_x64(
985      IntPtr cvode_mem,
986      ref double tret,
987      IntPtr yS //N_Vector*, one vector for each parameter
988      );
989
990    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "CVodeQuadInit", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
991    private static extern int CVodeQuadInit_x64(
992      IntPtr cvode_mem,
993      CVQuadRhsFn qF,
994      IntPtr yQ0 // N_Vector, initial value of yQ
995      );
996
997    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "CVodeQuadInitB", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
998    private static extern int CVodeQuadInitB_x64(
999      IntPtr cvode_mem,
1000      int indexB,
1001      CVQuadRhsFnB rhsQB,
1002      IntPtr yQB0 // N_Vector, initial value of yQB
1003    );
1004
1005
1006    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "CVodeAdjInit", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1007    private static extern int CVodeAdjInit_x64(
1008      IntPtr cvode_mem,
1009      int nd,  // integration steps between checkpoints
1010      int interpType // either CV_POLYNOMIAL or CV_HERMITE
1011      );
1012
1013    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "CVodeF", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1014    private static extern int CVodeF_x64(
1015        IntPtr cvode_mem,
1016        double t_out, // the next time at which a computed solution is desired
1017        IntPtr y, // N_Vector, the solution vector y
1018        ref double t_ret, // the time reached by the solver (output)
1019        int itask, // CV_NORMAL or CV_ONE_STEP
1020        ref int ncheck  // the number of internal checkpoints stored so far.
1021      );
1022
1023    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "CVodeGetNumSteps", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1024    private static extern int CVodeGetNumSteps_x64(
1025      IntPtr cvode_mem,
1026      ref long ncheck  // the number of steps taken
1027    );
1028
1029    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "CVodeCreateB", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1030    private static extern int CVodeCreateB_x64(
1031        IntPtr cvode_mem,
1032        MultistepMethod lmmB,
1033        NonlinearSolverIteration iterB,
1034        ref int which
1035      );
1036
1037    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "CVodeInitB", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1038    private static extern int CVodeInitB_x64(
1039        IntPtr cvode_mem,
1040        int which,
1041        CVRhsFuncB rightHandSide,
1042        double tB0, // endpoint T where final conditions are provided (equal to endpoint of forward integration)
1043        IntPtr yB0 // N_Vector inital value at t=tb0 of backward solution
1044      );
1045
1046    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "CVode", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1047    private static extern int CVode_x64(
1048      IntPtr cvode_mem,
1049      double tout, // next time at which a solution is desired
1050      IntPtr yout, // N_Vector
1051      ref double tret, // the time reached by the solver (output)
1052      int itask // flag indicating the job of the solver for the next step.
1053      );
1054
1055
1056    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "CVodeSStolerancesB", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1057    private static extern int CVodeSStolerancesB_x64(
1058      IntPtr cvode_mem,
1059      int which,
1060      double reltol,
1061      double abstol
1062    );
1063
1064    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "CVDlsSetLinearSolverB", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1065    private static extern int CVDlsSetLinearSolverB_x64(
1066      IntPtr cvode_mem,
1067      int which,
1068      IntPtr linearSolver, // SUNLinearSolver
1069      IntPtr j // SUNMatrix
1070    );
1071
1072    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "CVDlsSetJacFnB", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1073    private static extern int CVDlsSetJacFnB_x64(
1074      IntPtr cvode_mem,
1075      int which,
1076      CVDlsJacFuncB jacFunc
1077    );
1078
1079    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "CVodeB", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1080    private static extern int CVodeB_x64(
1081      IntPtr cvode_mem,
1082      double tout, // next time at which a solution is desired
1083      int itask // flag indicating the job of the solver for the next step.
1084    );
1085
1086    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "CVodeGetB", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1087    private static extern int CVodeGetB_x64(
1088          IntPtr cvode_mem,
1089          int which,
1090          ref double tret,
1091          IntPtr yB // the backward solution at time tret
1092        );
1093    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "CVodeGetAdjY", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1094    private static extern int CVodeGetAdjY_x64(
1095          IntPtr cvode_mem,
1096          double t,
1097          IntPtr y // the forward solution y(t)
1098        );
1099
1100    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "CVodeGetAdjCVodeBmem", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1101    private static extern IntPtr CVodeGetAdjCVodeBmem_x64(
1102              IntPtr cvode_mem,
1103              int which
1104            );
1105
1106    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "CVodeGetQuadB", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1107    private static extern int CVodeGetQuadB_x64(
1108        IntPtr cvode_mem,
1109        int which,
1110        ref double tret,
1111        IntPtr yQB // N_Vector
1112      );
1113
1114    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "CVodeFree", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1115    private static extern int CVodeFree_x64(ref IntPtr cvode_mem);
1116
1117
1118    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "CVodeSetErrHandlerFn", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1119    private static extern int CVodeSetErrHandlerFn_x64(IntPtr cvode_mem, CVErrHandlerFn errorHandler, IntPtr data);
1120
1121    #region matrix
1122    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "SUNDenseMatrix", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1123    private static extern IntPtr SUNDenseMatrix_x64(long m, long n);
1124
1125    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "SUNMatDestroy", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1126    private static extern void SUNMatDestroy_x64(IntPtr A);
1127
1128    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "SUNDenseMatrix_Data", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1129    private unsafe static extern double* SUNDenseMatrix_Data_x64(IntPtr matrix);
1130
1131    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "SUNDenseMatrix_Cols", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1132    private static extern long SUNDenseMatrix_Cols_x64(IntPtr matrix);
1133
1134    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "SUNDenseMatrix_Rows", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1135    private static extern long SUNDenseMatrix_Rows_x64(IntPtr matrix);
1136
1137    #endregion
1138
1139    #region linear solver
1140    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "SUNDenseLinearSolver", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1141    private static extern IntPtr SUNDenseLinearSolver_x64(
1142      IntPtr y, // N_Vector
1143      IntPtr A // SUNMatrix
1144      );
1145
1146    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "SUNLinSolInitialize", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1147    private static extern int SUNLinSolInitialize_x64(IntPtr linearSolver);
1148
1149    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "SUNLinSolSetup", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1150    private static extern int SUNLinSolSetup_x64(
1151      IntPtr linearSolver,
1152      IntPtr A // SUNMatrix
1153      );
1154
1155    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "SUNLinSolSolve", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1156    private static extern int SUNLinSolSolve_x64(
1157      IntPtr linearSolver,
1158      IntPtr A, // SUNMatrix
1159      IntPtr x, // N_Vector
1160      IntPtr b, // N_Vector
1161      double tol
1162      );
1163
1164    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "SUNLinSolFree", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1165    private static extern int SUNLinSolFree_x64(IntPtr linearSolver);
1166
1167    #endregion
1168
1169    #region N_Vector
1170    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "N_VNew_Serial", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1171    private static extern IntPtr N_VNew_Serial_x64(long vec_length);
1172
1173    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "N_VDestroy_Serial", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1174    private static extern void N_VDestroy_Serial_x64(IntPtr vec);
1175
1176    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "N_VDestroyVectorArray_Serial", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1177    private static extern void N_VDestroyVectorArray_Serial_x64(IntPtr vecArr, int count);
1178
1179    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "N_VPrint_Serial", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1180    private static extern void N_VPrint_Serial_x64(IntPtr vec);
1181
1182    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "N_VConst_Serial", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1183    private static extern void N_VConst_Serial_x64(double c, IntPtr vec);
1184
1185
1186    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "N_VL1Norm_Serial", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1187    private static extern double N_VL1Norm_Serial_x64(IntPtr vec);
1188
1189    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "N_VMake_Serial", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1190    private static extern IntPtr N_VMake_Serial_x64(long vec_length, double[] v_data);
1191
1192    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "N_VScale", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1193    ///  Performs the operation z = c*x
1194    private static extern void N_VScale_x64(double s,
1195      IntPtr x, // N_Vector
1196      IntPtr z // N_Vector
1197      );
1198
1199    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "N_VMake_Serial", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1200    private unsafe static extern double* N_VGetArrayPointer_Serial_x64(IntPtr vec);
1201
1202    [DllImport("sundials_cvodes-x64.dll", EntryPoint = "N_VCloneVectorArray_Serial", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
1203    private static extern IntPtr N_VCloneVectorArray_Serial_x64(int count, IntPtr vec); // returns N_Vector* !
1204
1205    #endregion
1206
1207    #endregion
1208  }
1209}
Note: See TracBrowser for help on using the repository browser.