Free cookie consent management tool by TermsFeed Policy Generator

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

Last change on this file since 16254 was 16254, checked in by gkronber, 6 years ago

#2925: made usage of x86 native library explicit.

File size: 20.1 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    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeCreate", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
177    // returns a void* to the cvodes memory block if successful otherwise NULL
178    public static extern IntPtr CVodeCreate(MultistepMethod lmm, NonlinearSolverIteration iter);
179
180    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeInit", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
181    public static extern int CVodeInit(
182      IntPtr cvode_mem, // pointer returned by CVodeCreate
183      CVRhsFunc f,
184      double t0, // realtype, the inital value of t
185      IntPtr y0 // N_Vector the initial value of y
186    );
187
188
189    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeSStolerances", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
190    public static extern int CVodeSStolerances(
191      IntPtr cvode_mem,
192      double reltol,
193      double abstol
194      );
195
196    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVDlsSetLinearSolver", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
197    public static extern int CVDlsSetLinearSolver(
198      IntPtr cvode_mem,
199      IntPtr linearSolver, // SUNLinearSolver
200      IntPtr j // SUNMatrix
201      );
202
203    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVDlsSetJacFn", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
204    public static extern int CVDlsSetJacFn(
205      IntPtr cvode_mem,
206      CVDlsJacFunc jacFunc
207      );
208
209    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeSensInit", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
210    public static extern int CVodeSensInit(
211      IntPtr cvode_mem,
212      int Ns, // number of parameters
213      int ism, // sensitivity solution method, CV_SIMULTANEOUS or CV_STAGGERED
214      CVSensRhsFn fS, // right hand side function which computes all sensitivity RHS at the same time
215      IntPtr yS0 // N_Vector
216      );
217
218    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeSensSStolerances", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
219    public static extern int CVodeSensSStolerances(
220      IntPtr cvode_mem,
221      double reltol,
222      double[] abstol
223    );
224
225    /* Call CVodeSensEEtolerances to estimate tolerances for sensitivity
226   variables based on the rolerances supplied for states variables and
227   the scaling factor pbar */
228    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeSensEEtolerances", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
229    public static extern int CVodeSensEEtolerances(
230      IntPtr cvode_mem
231      );
232
233    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeGetSens", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
234    public static extern int CVodeGetSens(
235      IntPtr cvode_mem,
236      ref double tret,
237      IntPtr yS //N_Vector*, one vector for each parameter
238      );
239
240    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeQuadInit", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
241    public static extern int CVodeQuadInit(
242      IntPtr cvode_mem,
243      CVQuadRhsFn qF,
244      IntPtr yQ0 // N_Vector, initial value of yQ
245      );
246
247    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeQuadInitB", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
248    public static extern int CVodeQuadInitB(
249      IntPtr cvode_mem,
250      int indexB,
251      CVQuadRhsFnB rhsQB,
252      IntPtr yQB0 // N_Vector, initial value of yQB
253    );
254
255
256    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeAdjInit", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
257    public static extern int CVodeAdjInit(
258      IntPtr cvode_mem,
259      int nd,  // integration steps between checkpoints
260      int interpType // either CV_POLYNOMIAL or CV_HERMITE
261      );
262
263    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeF", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
264    public static extern int CVodeF(
265        IntPtr cvode_mem,
266        double t_out, // the next time at which a computed solution is desired
267        IntPtr y, // N_Vector, the solution vector y
268        ref double t_ret, // the time reached by the solver (output)
269        int itask, // CV_NORMAL or CV_ONE_STEP
270        ref int ncheck  // the number of internal checkpoints stored so far.
271      );
272
273    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeGetNumSteps", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
274    public static extern int CVodeGetNumSteps(
275      IntPtr cvode_mem,
276      ref long ncheck  // the number of steps taken
277    );
278
279    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeCreateB", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
280    public static extern int CVodeCreateB(
281        IntPtr cvode_mem,
282        MultistepMethod lmmB,
283        NonlinearSolverIteration iterB,
284        ref int which
285      );
286
287    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeInitB", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
288    public static extern int CVodeInitB(
289        IntPtr cvode_mem,
290        int which,
291        CVRhsFuncB rightHandSide,
292        double tB0, // endpoint T where final conditions are provided (equal to endpoint of forward integration)
293        IntPtr yB0 // N_Vector inital value at t=tb0 of backward solution
294      );
295
296    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVode", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
297    public static extern int CVode(
298      IntPtr cvode_mem,
299      double tout, // next time at which a solution is desired
300      IntPtr yout, // N_Vector
301      ref double tret, // the time reached by the solver (output)
302      int itask // flag indicating the job of the solver for the next step.
303      );
304
305
306    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeSStolerancesB", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
307    public static extern int CVodeSStolerancesB(
308      IntPtr cvode_mem,
309      int which,
310      double reltol,
311      double abstol
312    );
313
314    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVDlsSetLinearSolverB", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
315    public static extern int CVDlsSetLinearSolverB(
316      IntPtr cvode_mem,
317      int which,
318      IntPtr linearSolver, // SUNLinearSolver
319      IntPtr j // SUNMatrix
320    );
321
322    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVDlsSetJacFnB", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
323    public static extern int CVDlsSetJacFnB(
324      IntPtr cvode_mem,
325      int which,
326      CVDlsJacFuncB jacFunc
327    );
328
329    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeB", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
330    public static extern int CVodeB(
331      IntPtr cvode_mem,
332      double tout, // next time at which a solution is desired
333      int itask // flag indicating the job of the solver for the next step.
334    );
335
336    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeGetB", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
337    public static extern int CVodeGetB(
338          IntPtr cvode_mem,
339          int which,
340          ref double tret,
341          IntPtr yB // the backward solution at time tret
342        );
343    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeGetAdjY", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
344    public static extern int CVodeGetAdjY(
345          IntPtr cvode_mem,
346          double t,
347          IntPtr y // the forward solution y(t)
348        );
349
350    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeGetAdjCVodeBmem", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
351    public static extern IntPtr CVodeGetAdjCVodeBmem(
352              IntPtr cvode_mem,
353              int which
354            );
355
356    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeGetQuadB", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
357    public static extern int CVodeGetQuadB(
358        IntPtr cvode_mem,
359        int which,
360        ref double tret,
361        IntPtr yQB // N_Vector
362      );
363
364    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "CVodeFree", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
365
366    public static extern int CVodeFree(ref IntPtr cvode_mem);
367
368    #region matrix
369    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "SUNDenseMatrix", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
370    public static extern IntPtr SUNDenseMatrix(long m, long n);
371
372    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "SUNMatDestroy", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
373    public static extern void SUNMatDestroy(IntPtr A);
374
375    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "SUNDenseMatrix_Data", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
376    public unsafe static extern double* SUNDenseMatrix_Data(IntPtr matrix);
377
378    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "SUNDenseMatrix_Cols", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
379    public static extern long SUNDenseMatrix_Cols(IntPtr matrix);
380
381    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "SUNDenseMatrix_Rows", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
382    public static extern long SUNDenseMatrix_Rows(IntPtr matrix);
383
384
385    public unsafe static double SUNDenseMatrix_Get(IntPtr mat, long i, long j) {
386      long M = SUNDenseMatrix_Rows(mat);
387      // the(i, j)th element of A(with 0 <= i<M and 0 <= j<N) is given by(A->data)[j*M+i]
388      return SUNDenseMatrix_Data(mat)[j * M + i];
389    }
390
391    public unsafe static void SUNDenseMatrix_Set(IntPtr mat, long i, long j, double val) {
392      long M = SUNDenseMatrix_Rows(mat);
393      // the(i, j)th element of A(with 0 <= i<M and 0 <= j<N) is given by(A->data)[j*M+i]
394      SUNDenseMatrix_Data(mat)[j * M + i] = val;
395    }
396
397    #endregion
398
399    #region linear solver
400    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "SUNDenseLinearSolver", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
401    public static extern IntPtr SUNDenseLinearSolver(
402      IntPtr y, // N_Vector
403      IntPtr A // SUNMatrix
404      );
405
406    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "SUNLinSolInitialize", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
407    public static extern int SUNLinSolInitialize(IntPtr linearSolver);
408
409    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "SUNLinSolSetup", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
410    public static extern int SUNLinSolSetup(
411      IntPtr linearSolver,
412      IntPtr A // SUNMatrix
413      );
414
415    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "SUNLinSolSolve", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
416    public static extern int SUNLinSolSolve(
417      IntPtr linearSolver,
418      IntPtr A, // SUNMatrix
419      IntPtr x, // N_Vector
420      IntPtr b, // N_Vector
421      double tol
422      );
423
424    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "SUNLinSolFree", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
425    public static extern int SUNLinSolFree(IntPtr linearSolver);
426
427    #endregion
428
429    #region N_Vector
430    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "N_VNew_Serial", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
431    public static extern IntPtr N_VNew_Serial(long vec_length);
432
433    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "N_VDestroy_Serial", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
434    public static extern void N_VDestroy_Serial(IntPtr vec);
435
436    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "N_VDestroyVectorArray_Serial", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
437    public static extern void N_VDestroyVectorArray_Serial(IntPtr vecArr, int count);
438
439    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "N_VPrint_Serial", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
440    public static extern void N_VPrint_Serial(IntPtr vec);
441
442    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "N_VConst_Serial", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
443    public static extern void N_VConst_Serial(double c, IntPtr vec);
444
445
446    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "N_VL1Norm_Serial", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
447    public static extern double N_VL1Norm_Serial(IntPtr vec);
448
449    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "N_VMake_Serial", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
450    public static extern IntPtr N_VMake_Serial(long vec_length, double[] v_data);
451
452    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "N_VScale", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
453    ///  Performs the operation z = c*x
454    public static extern void N_VScale(double s,
455      IntPtr x, // N_Vector
456      IntPtr z // N_Vector
457      );
458
459    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "N_VMake_Serial", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
460    public unsafe static extern double* N_VGetArrayPointer_Serial(IntPtr vec);
461
462    [DllImport("sundials_cvodes-x86.dll", EntryPoint = "N_VCloneVectorArray_Serial", ExactSpelling = true, CallingConvention = CallingConvention.Cdecl)]
463    public static extern IntPtr N_VCloneVectorArray_Serial(int count, IntPtr vec); // returns N_Vector* !
464
465
466    /*
467#define NV_CONTENT_S(v)  ( (N_VectorContent_Serial)(v->content) )
468
469#define NV_LENGTH_S(v)   ( NV_CONTENT_S(v)->length )
470
471#define NV_OWN_DATA_S(v) ( NV_CONTENT_S(v)->own_data )
472
473#define NV_DATA_S(v)     ( NV_CONTENT_S(v)->data )
474
475#define NV_Ith_S(v,i)    ( NV_DATA_S(v)[i] )
476*/
477    // methods for macros
478    public unsafe static int* NV_CONTENT_S(IntPtr v) {
479      int* content = (int*)*(int*)v.ToPointer();
480      return content;
481    }
482
483    public unsafe static long NV_LENGTH_S(IntPtr v) {
484      long length = *NV_CONTENT_S(v);
485      return length;
486    }
487    public unsafe static bool NV_OWN_DATA_S(IntPtr v) {
488      var content = NV_CONTENT_S(v);
489      int own_data = *(content + 2);
490      return own_data > 0;
491    }
492    public unsafe static double* NV_DATA_S(IntPtr v) {
493      var content = NV_CONTENT_S(v);
494      double* data = (double*)*(content + 3);
495      return data;
496    }
497    public unsafe static double NV_Get_Ith_S(IntPtr v, long i) {
498      return NV_DATA_S(v)[i];
499    }
500    public unsafe static void NV_Set_Ith_S(IntPtr v, long i, double val) {
501      NV_DATA_S(v)[i] = val;
502    }
503    #endregion
504
505
506  }
507}
Note: See TracBrowser for help on using the repository browser.