[16225] | 1 | using System;
|
---|
| 2 | using HeuristicLab.Problems.DynamicalSystemsModelling;
|
---|
| 3 | using Microsoft.VisualStudio.TestTools.UnitTesting;
|
---|
| 4 |
|
---|
| 5 | namespace AutoDiffForDynamicalModelsTest {
|
---|
| 6 | [TestClass]
|
---|
| 7 | public class TestCvodes {
|
---|
| 8 | [TestMethod]
|
---|
[16246] | 9 | public unsafe void TestCvodesASAMethod() {
|
---|
[16225] | 10 | // test vectors
|
---|
| 11 | var arr = new double[] { 3.14, 2.71 };
|
---|
[16248] | 12 | var vec = CVODES.N_VMake_Serial(2, arr);
|
---|
| 13 | CVODES.N_VDestroy_Serial(vec);
|
---|
[16225] | 14 |
|
---|
[16248] | 15 | vec = CVODES.N_VNew_Serial(10);
|
---|
[16225] | 16 |
|
---|
| 17 | unsafe {
|
---|
| 18 | int* content = (int*)*(int*)vec.ToPointer();
|
---|
| 19 | long length = *content;
|
---|
| 20 | Console.WriteLine(*content);
|
---|
| 21 | int own_data = *(content + 2);
|
---|
| 22 | Console.WriteLine(own_data);
|
---|
[16245] | 23 | double* data = (double*)*(content + 3);
|
---|
[16225] | 24 | double data0 = *data;
|
---|
| 25 | }
|
---|
[16248] | 26 | CVODES.N_VConst_Serial(2.0, vec);
|
---|
| 27 | CVODES.N_VPrint_Serial(vec);
|
---|
| 28 | Assert.AreEqual(20, CVODES.N_VL1Norm_Serial(vec));
|
---|
| 29 | CVODES.N_VDestroy_Serial(vec);
|
---|
[16225] | 30 |
|
---|
| 31 |
|
---|
[16227] | 32 | // linear oscillator
|
---|
| 33 | int numberOfEquations = 2;
|
---|
[16248] | 34 | var y = CVODES.N_VNew_Serial(numberOfEquations);
|
---|
[16226] | 35 | // y must be initialized before calling CVodeInit
|
---|
[16248] | 36 | // CVODES.N_VConst_Serial(100.0, y);
|
---|
| 37 | CVODES.NV_Set_Ith_S(y, 0, 0.5); // x
|
---|
| 38 | CVODES.NV_Set_Ith_S(y, 1, 1); // v
|
---|
[16225] | 39 |
|
---|
[16248] | 40 | var cvode_mem = CVODES.CVodeCreate(CVODES.MultistepMethod.CV_ADAMS, CVODES.NonlinearSolverIteration.CV_FUNCTIONAL);
|
---|
[16225] | 41 |
|
---|
[16248] | 42 | var flag = CVODES.CVodeInit(cvode_mem, F, 0.0, y);
|
---|
| 43 | Assert.AreEqual(CVODES.CV_SUCCESS, flag);
|
---|
[16225] | 44 |
|
---|
[16248] | 45 | flag = CVODES.CVodeSStolerances(cvode_mem, 1E-4, 1.0);
|
---|
| 46 | Assert.AreEqual(CVODES.CV_SUCCESS, flag);
|
---|
[16225] | 47 |
|
---|
[16248] | 48 | var A = CVODES.SUNDenseMatrix(numberOfEquations, numberOfEquations);
|
---|
[16225] | 49 | Assert.AreNotSame(A, IntPtr.Zero);
|
---|
| 50 |
|
---|
| 51 |
|
---|
[16248] | 52 | var linearSolver = CVODES.SUNDenseLinearSolver(y, A);
|
---|
[16225] | 53 | Assert.AreNotSame(linearSolver, IntPtr.Zero);
|
---|
| 54 |
|
---|
[16248] | 55 | flag = CVODES.CVDlsSetLinearSolver(cvode_mem, linearSolver, A);
|
---|
| 56 | Assert.AreEqual(CVODES.CV_SUCCESS, flag);
|
---|
[16225] | 57 |
|
---|
[16248] | 58 | // flag = CVODES.CVDlsSetJacFn(cvode_mem, JacF);
|
---|
| 59 | // Assert.AreEqual(CVODES.CV_SUCCESS, flag);
|
---|
[16245] | 60 |
|
---|
| 61 | // var ns = 1; // number of parameters
|
---|
| 62 | // var p = new double[1]; // set as user-data
|
---|
[16248] | 63 | // var yS0 = CVODES.N_VCloneVectorArray_Serial(ns, y); // clone the output vector for each parameter, TODO: free
|
---|
[16245] | 64 | // for(int i=0;i<ns;i++) {
|
---|
| 65 | //
|
---|
| 66 | // }
|
---|
| 67 |
|
---|
| 68 |
|
---|
[16248] | 69 | var q = CVODES.N_VNew_Serial(1);
|
---|
| 70 | CVODES.NV_Set_Ith_S(q, 0, 0.0);
|
---|
| 71 | flag = CVODES.CVodeQuadInit(cvode_mem, FQ, q);
|
---|
| 72 | Assert.AreEqual(CVODES.CV_SUCCESS, flag);
|
---|
[16225] | 73 |
|
---|
| 74 |
|
---|
[16245] | 75 | int steps = 10;
|
---|
[16248] | 76 | flag = CVODES.CVodeAdjInit(cvode_mem, steps, CVODES.CV_HERMITE);
|
---|
| 77 | Assert.AreEqual(CVODES.CV_SUCCESS, flag);
|
---|
[16245] | 78 |
|
---|
| 79 | /* step by step forward integration
|
---|
[16225] | 80 | double t = 0.0;
|
---|
| 81 | double tout = 0.1; // first output time
|
---|
| 82 | int nout = 100; // number of output times
|
---|
| 83 | for (int iout = 0; iout < nout; iout++) {
|
---|
[16248] | 84 | flag = CVODES.CVode(cvode_mem, tout, y, ref t, CVODES.CV_NORMAL);
|
---|
| 85 | Assert.AreEqual(CVODES.CV_SUCCESS, flag);
|
---|
| 86 | Console.WriteLine("{0} {1} {2}", t, CVODES.NV_Get_Ith_S(y, 0), CVODES.NV_Get_Ith_S(y, 1));
|
---|
| 87 | // CVODES.N_VPrint_Serial(y);
|
---|
[16225] | 88 | tout += 0.1;
|
---|
| 89 | }
|
---|
[16245] | 90 | */
|
---|
[16225] | 91 |
|
---|
[16245] | 92 | // complete forward integration
|
---|
| 93 | double tout = 100.0; // last time
|
---|
| 94 | double time = 0.0;
|
---|
| 95 | int ncheck = 0; // number of checkpoints
|
---|
[16248] | 96 | flag = CVODES.CVodeF(cvode_mem, tout, y, ref time, CVODES.CV_NORMAL, ref ncheck);
|
---|
| 97 | Assert.AreEqual(CVODES.CV_SUCCESS, flag);
|
---|
[16245] | 98 |
|
---|
| 99 | long numSteps = 0;
|
---|
[16248] | 100 | flag = CVODES.CVodeGetNumSteps(cvode_mem, ref numSteps);
|
---|
| 101 | Assert.AreEqual(CVODES.CV_SUCCESS, flag);
|
---|
[16245] | 102 |
|
---|
[16248] | 103 | var yB = CVODES.N_VNew_Serial(numberOfEquations);
|
---|
| 104 | CVODES.N_VConst_Serial(0.0, yB);
|
---|
[16245] | 105 |
|
---|
| 106 | int numberOfParameters = 2;
|
---|
[16248] | 107 | var qB = CVODES.N_VNew_Serial(numberOfParameters);
|
---|
| 108 | CVODES.N_VConst_Serial(0.0, qB);
|
---|
[16245] | 109 |
|
---|
| 110 | int indexB = 0;
|
---|
[16248] | 111 | flag = CVODES.CVodeCreateB(cvode_mem, CVODES.MultistepMethod.CV_BDF, CVODES.NonlinearSolverIteration.CV_NEWTON, ref indexB);
|
---|
| 112 | Assert.AreEqual(CVODES.CV_SUCCESS, flag);
|
---|
[16245] | 113 |
|
---|
| 114 | var TB1 = tout;
|
---|
[16248] | 115 | flag = CVODES.CVodeInitB(cvode_mem, indexB, FB, TB1, yB);
|
---|
| 116 | Assert.AreEqual(CVODES.CV_SUCCESS, flag);
|
---|
[16245] | 117 |
|
---|
| 118 | var relTolB = 1E-6;
|
---|
| 119 | var absTolB = 1E-8;
|
---|
[16248] | 120 | flag = CVODES.CVodeSStolerancesB(cvode_mem, indexB, relTolB, absTolB);
|
---|
| 121 | Assert.AreEqual(CVODES.CV_SUCCESS, flag);
|
---|
[16245] | 122 |
|
---|
[16248] | 123 | var AB = CVODES.SUNDenseMatrix(numberOfEquations, numberOfEquations);
|
---|
[16245] | 124 | Assert.AreNotSame(linearSolver, IntPtr.Zero);
|
---|
| 125 |
|
---|
[16248] | 126 | var lsB = CVODES.SUNDenseLinearSolver(yB, AB);
|
---|
[16245] | 127 | Assert.AreNotSame(linearSolver, IntPtr.Zero);
|
---|
| 128 |
|
---|
[16248] | 129 | flag = CVODES.CVDlsSetLinearSolverB(cvode_mem, indexB, lsB, AB);
|
---|
| 130 | Assert.AreEqual(CVODES.CV_SUCCESS, flag);
|
---|
[16245] | 131 |
|
---|
[16248] | 132 | flag = CVODES.CVDlsSetJacFnB(cvode_mem, indexB, JacFB);
|
---|
| 133 | Assert.AreEqual(CVODES.CV_SUCCESS, flag);
|
---|
[16245] | 134 |
|
---|
[16248] | 135 | flag = CVODES.CVodeQuadInitB(cvode_mem, indexB, FQB, qB);
|
---|
| 136 | Assert.AreEqual(CVODES.CV_SUCCESS, flag);
|
---|
[16245] | 137 |
|
---|
| 138 | /* First get results at t = TBout1 */
|
---|
| 139 |
|
---|
[16248] | 140 | /* Call CVodeB to integrate the backward ODE CVODES. */
|
---|
[16245] | 141 | var tBackOut = 50.0;
|
---|
[16248] | 142 | flag = CVODES.CVodeB(cvode_mem, tBackOut, CVODES.CV_NORMAL);
|
---|
| 143 | Assert.AreEqual(CVODES.CV_SUCCESS, flag);
|
---|
[16245] | 144 |
|
---|
[16248] | 145 | /* Call CVodeGetB to get yB of the backward ODE CVODES. */
|
---|
| 146 | flag = CVODES.CVodeGetB(cvode_mem, indexB, ref time, yB);
|
---|
| 147 | Assert.AreEqual(CVODES.CV_SUCCESS, flag);
|
---|
[16245] | 148 |
|
---|
| 149 | /* Call CVodeGetAdjY to get the interpolated value of the forward solution
|
---|
| 150 | y during a backward integration. */
|
---|
[16248] | 151 | flag = CVODES.CVodeGetAdjY(cvode_mem, tBackOut, y);
|
---|
| 152 | Assert.AreEqual(CVODES.CV_SUCCESS, flag);
|
---|
[16245] | 153 |
|
---|
| 154 | /* Then at t = T0 */
|
---|
| 155 |
|
---|
| 156 | double t0 = 0.0;
|
---|
[16248] | 157 | flag = CVODES.CVodeB(cvode_mem, t0, CVODES.CV_NORMAL);
|
---|
| 158 | Assert.AreEqual(CVODES.CV_SUCCESS, flag);
|
---|
[16245] | 159 |
|
---|
| 160 | long nstB = 0;
|
---|
[16248] | 161 | CVODES.CVodeGetNumSteps(CVODES.CVodeGetAdjCVodeBmem(cvode_mem, indexB), ref nstB);
|
---|
[16245] | 162 |
|
---|
[16248] | 163 | flag = CVODES.CVodeGetB(cvode_mem, indexB, ref time, yB);
|
---|
| 164 | Assert.AreEqual(CVODES.CV_SUCCESS, flag);
|
---|
[16245] | 165 |
|
---|
| 166 |
|
---|
[16248] | 167 | flag = CVODES.CVodeGetQuadB(cvode_mem, indexB, ref time, qB);
|
---|
| 168 | Assert.AreEqual(CVODES.CV_SUCCESS, flag);
|
---|
[16245] | 169 |
|
---|
[16248] | 170 | flag = CVODES.CVodeGetAdjY(cvode_mem, t0, y);
|
---|
| 171 | Assert.AreEqual(CVODES.CV_SUCCESS, flag);
|
---|
[16245] | 172 |
|
---|
[16248] | 173 | CVODES.N_VDestroy_Serial(y);
|
---|
[16253] | 174 | CVODES.CVodeFree(ref cvode_mem);
|
---|
[16248] | 175 | CVODES.SUNLinSolFree(linearSolver);
|
---|
| 176 | CVODES.SUNMatDestroy(A);
|
---|
[16225] | 177 | }
|
---|
| 178 |
|
---|
[16246] | 179 | [TestMethod]
|
---|
| 180 | public void TestCvodesFSAMethod() {
|
---|
[16245] | 181 |
|
---|
[16246] | 182 | // linear oscillator
|
---|
| 183 | int numberOfEquations = 2;
|
---|
[16248] | 184 | var y = CVODES.N_VNew_Serial(numberOfEquations);
|
---|
[16246] | 185 | // y must be initialized before calling CVodeInit
|
---|
[16248] | 186 | // CVODES.N_VConst_Serial(100.0, y);
|
---|
| 187 | CVODES.NV_Set_Ith_S(y, 0, 0.5); // x
|
---|
| 188 | CVODES.NV_Set_Ith_S(y, 1, 1); // v
|
---|
[16225] | 189 |
|
---|
[16253] | 190 | var cvode_mem = CVODES.CVodeCreate(CVODES.MultistepMethod.CV_BDF, CVODES.NonlinearSolverIteration.CV_NEWTON);
|
---|
[16246] | 191 |
|
---|
[16248] | 192 | var flag = CVODES.CVodeInit(cvode_mem, F, 0.0, y);
|
---|
| 193 | Assert.AreEqual(CVODES.CV_SUCCESS, flag);
|
---|
[16246] | 194 |
|
---|
[16248] | 195 | flag = CVODES.CVodeSStolerances(cvode_mem, 1E-4, 1.0);
|
---|
| 196 | Assert.AreEqual(CVODES.CV_SUCCESS, flag);
|
---|
[16246] | 197 |
|
---|
[16248] | 198 | var A = CVODES.SUNDenseMatrix(numberOfEquations, numberOfEquations);
|
---|
[16246] | 199 | Assert.AreNotSame(A, IntPtr.Zero);
|
---|
[16248] | 200 | Console.WriteLine(CVODES.SUNDenseMatrix_Get(A, 0, 0));
|
---|
[16246] | 201 |
|
---|
| 202 |
|
---|
[16248] | 203 | var linearSolver = CVODES.SUNDenseLinearSolver(y, A);
|
---|
[16246] | 204 | Assert.AreNotSame(linearSolver, IntPtr.Zero);
|
---|
| 205 |
|
---|
[16248] | 206 | flag = CVODES.CVDlsSetLinearSolver(cvode_mem, linearSolver, A);
|
---|
| 207 | Assert.AreEqual(CVODES.CV_SUCCESS, flag);
|
---|
[16246] | 208 |
|
---|
[16248] | 209 | flag = CVODES.CVDlsSetJacFn(cvode_mem, JacF);
|
---|
| 210 | Assert.AreEqual(CVODES.CV_SUCCESS, flag);
|
---|
[16246] | 211 |
|
---|
| 212 | var ns = 1; // number of parameters
|
---|
| 213 | var p = new double[1]; // set as user-data
|
---|
[16248] | 214 | var yS0 = CVODES.N_VCloneVectorArray_Serial(ns, y); // clone the output vector for each parameter, TODO: free
|
---|
[16246] | 215 | unsafe {
|
---|
| 216 | for (int i = 0; i < ns; i++) {
|
---|
| 217 | var yS0_i = *((IntPtr*)yS0.ToPointer() + i);
|
---|
[16248] | 218 | CVODES.N_VConst_Serial(0.0, yS0_i);
|
---|
[16246] | 219 | }
|
---|
| 220 | }
|
---|
| 221 |
|
---|
| 222 |
|
---|
[16248] | 223 | flag = CVODES.CVodeSensInit(cvode_mem, ns, CVODES.CV_SIMULTANEOUS, FS, yS0);
|
---|
| 224 | Assert.AreEqual(CVODES.CV_SUCCESS, flag);
|
---|
[16246] | 225 |
|
---|
[16248] | 226 | // flag = CVODES.CVodeSensSStolerances(cvode_mem, reltolS, abstolS);
|
---|
| 227 | // Assert.AreEqual(CVODES.CV_SUCCESS, flag);
|
---|
| 228 | flag = CVODES.CVodeSensEEtolerances(cvode_mem);
|
---|
| 229 | Assert.AreEqual(CVODES.CV_SUCCESS, flag);
|
---|
[16246] | 230 |
|
---|
| 231 |
|
---|
[16248] | 232 | // var q = CVODES.N_VNew_Serial(1);
|
---|
| 233 | // CVODES.NV_Set_Ith_S(q, 0, 0.0);
|
---|
| 234 | // flag = CVODES.CVodeQuadInit(cvode_mem, FQ, q);
|
---|
| 235 | // Assert.AreEqual(CVODES.CV_SUCCESS, flag);
|
---|
[16246] | 236 |
|
---|
| 237 |
|
---|
| 238 | // int steps = 10;
|
---|
[16248] | 239 | // flag = CVODES.CVodeAdjInit(cvode_mem, steps, CVODES.CV_HERMITE);
|
---|
| 240 | // Assert.AreEqual(CVODES.CV_SUCCESS, flag);
|
---|
[16246] | 241 |
|
---|
| 242 | // step by step forward integration
|
---|
| 243 | double t = 0.0;
|
---|
| 244 | double tout = 0.1; // first output time
|
---|
| 245 | int nout = 1000; // number of output times
|
---|
| 246 | for (int iout = 0; iout < nout; iout++) {
|
---|
[16248] | 247 | flag = CVODES.CVode(cvode_mem, tout, y, ref t, CVODES.CV_NORMAL);
|
---|
| 248 | Assert.AreEqual(CVODES.CV_SUCCESS, flag);
|
---|
[16246] | 249 |
|
---|
| 250 | // get sensitivities
|
---|
[16248] | 251 | flag = CVODES.CVodeGetSens(cvode_mem, ref t, yS0);
|
---|
| 252 | Assert.AreEqual(CVODES.CV_SUCCESS, flag);
|
---|
[16246] | 253 | unsafe {
|
---|
| 254 | var ySP0 = *(IntPtr*)yS0.ToPointer(); // parameter of sensitivities for first component
|
---|
[16248] | 255 | Console.WriteLine("{0} {1} {2} {3} {4}", t, CVODES.NV_Get_Ith_S(y, 0), CVODES.NV_Get_Ith_S(y, 1), CVODES.NV_Get_Ith_S(ySP0, 0), CVODES.NV_Get_Ith_S(ySP0, 1));
|
---|
[16246] | 256 | }
|
---|
| 257 |
|
---|
| 258 | tout += 0.1;
|
---|
| 259 | }
|
---|
| 260 |
|
---|
| 261 |
|
---|
| 262 |
|
---|
[16248] | 263 | CVODES.N_VDestroy_Serial(y);
|
---|
[16253] | 264 | CVODES.CVodeFree(ref cvode_mem);
|
---|
[16248] | 265 | CVODES.SUNLinSolFree(linearSolver);
|
---|
| 266 | CVODES.SUNMatDestroy(A);
|
---|
[16245] | 267 | }
|
---|
| 268 |
|
---|
[16246] | 269 |
|
---|
| 270 |
|
---|
[16225] | 271 | public static int F(
|
---|
| 272 | double t, // realtype
|
---|
| 273 | IntPtr y, // N_Vector
|
---|
| 274 | IntPtr ydot, // N_Vector
|
---|
| 275 | IntPtr user_data) {
|
---|
| 276 |
|
---|
[16246] | 277 | // the system
|
---|
| 278 | // ∂y1/∂t = y2
|
---|
| 279 | // ∂y2/∂t = -0.3 y1
|
---|
[16248] | 280 | CVODES.NV_Set_Ith_S(ydot, 0, CVODES.NV_Get_Ith_S(y, 1));
|
---|
| 281 | CVODES.NV_Set_Ith_S(ydot, 1, -0.3 * CVODES.NV_Get_Ith_S(y, 0));
|
---|
[16225] | 282 | return 0;
|
---|
| 283 | ;
|
---|
| 284 | }
|
---|
[16246] | 285 |
|
---|
| 286 | // computes ∂f/∂y or an approximation
|
---|
| 287 | private int JacF(double t, IntPtr y, IntPtr fy, IntPtr Jac, IntPtr user_data, IntPtr tmp1, IntPtr tmp2, IntPtr tmp3) {
|
---|
| 288 | // t is current time
|
---|
| 289 | // y is current y
|
---|
| 290 | // fy is current f(t,y)
|
---|
| 291 | // Jac is output
|
---|
| 292 |
|
---|
| 293 | // Matrix:
|
---|
| 294 | // ∂f1/∂y1 ∂f1/∂y2
|
---|
| 295 | // ∂f2/∂y1 ∂f2/∂y2
|
---|
| 296 |
|
---|
[16248] | 297 | CVODES.SUNDenseMatrix_Set(Jac, 0, 0, 0.0);
|
---|
| 298 | CVODES.SUNDenseMatrix_Set(Jac, 0, 1, 1.0);
|
---|
| 299 | CVODES.SUNDenseMatrix_Set(Jac, 1, 0, -0.3);
|
---|
| 300 | CVODES.SUNDenseMatrix_Set(Jac, 1, 1, 0.0);
|
---|
[16246] | 301 | return 0;
|
---|
| 302 | }
|
---|
| 303 |
|
---|
| 304 | private int JacFB(double t, IntPtr y, IntPtr yB, IntPtr fyB, IntPtr Jac, IntPtr user_data, IntPtr tmp1, IntPtr tmp2, IntPtr tmp3) {
|
---|
| 305 | throw new NotImplementedException();
|
---|
| 306 | return 0;
|
---|
| 307 | }
|
---|
| 308 |
|
---|
| 309 | // sensitivity function must compute the vector (∂f/∂y)s_i(t) + ∂f/∂p_i and store in ySdot[i]
|
---|
| 310 | public static int FS(
|
---|
| 311 | int Ns,
|
---|
| 312 | double t,
|
---|
| 313 | IntPtr y, // N_Vector
|
---|
| 314 | IntPtr ydot, // N_Vector
|
---|
| 315 | IntPtr yS, // N_Vector* one vector for each parameter
|
---|
| 316 | IntPtr ySdot, // N_Vector* one vector for each parameter
|
---|
| 317 | IntPtr user_data,
|
---|
| 318 | IntPtr tmp1, // N_Vector
|
---|
| 319 | IntPtr tmp2 // N_Vector
|
---|
| 320 | ) {
|
---|
| 321 |
|
---|
| 322 | // (∂f /∂y)s_i(t) + ∂f /∂p_i
|
---|
[16248] | 323 | var y1 = CVODES.NV_Get_Ith_S(y, 0); var y2 = CVODES.NV_Get_Ith_S(y, 1);
|
---|
[16246] | 324 |
|
---|
| 325 | unsafe {
|
---|
| 326 | var yS_p0 = *(IntPtr*)yS.ToPointer();
|
---|
[16248] | 327 | var s1 = CVODES.NV_Get_Ith_S(yS_p0, 0); var s2 = CVODES.NV_Get_Ith_S(yS_p0, 1);
|
---|
[16246] | 328 |
|
---|
| 329 |
|
---|
[16248] | 330 | // (∂f /∂y)s_i(t) == Jac s_i(t)
|
---|
[16246] | 331 | var sd1 = s2;
|
---|
| 332 | var sd2 = -0.3 * s1;
|
---|
| 333 |
|
---|
| 334 | var s_p0 = *(IntPtr*)ySdot.ToPointer();
|
---|
[16248] | 335 | CVODES.NV_Set_Ith_S(s_p0, 0, sd1 + 0.0);
|
---|
| 336 | CVODES.NV_Set_Ith_S(s_p0, 1, sd2 - 0.3);
|
---|
[16246] | 337 | }
|
---|
| 338 |
|
---|
| 339 | return 0;
|
---|
| 340 | }
|
---|
| 341 |
|
---|
[16245] | 342 | public static int FB(
|
---|
| 343 | double t, // realtype
|
---|
| 344 | IntPtr y, // N_Vector
|
---|
| 345 | IntPtr yB, // N_Vector
|
---|
| 346 | IntPtr yBdot, // N_Vector
|
---|
| 347 | IntPtr user_data) {
|
---|
| 348 |
|
---|
| 349 | // y' = f(y,p)
|
---|
| 350 |
|
---|
| 351 | // yB = λ
|
---|
| 352 | // λ' = -(∂f / ∂y)^T λ - (∂g / ∂y)
|
---|
| 353 | // the goal is to find dG/dp where G = integrate( g(y,t,p), t=0, t=T) )
|
---|
| 354 |
|
---|
| 355 | // for F above:
|
---|
| 356 | // ∂f / ∂ =
|
---|
| 357 | // 0.0 1.0
|
---|
| 358 | // -0.3 0.0
|
---|
| 359 |
|
---|
| 360 | // we have no g!?
|
---|
| 361 | // therefore λ' =
|
---|
| 362 | // λ1' = 0.3 λ2 - 0.0
|
---|
| 363 | // λ2' = 1.0 λ1 - 0.0
|
---|
| 364 |
|
---|
[16248] | 365 | CVODES.NV_Set_Ith_S(yBdot, 0, 0.3 * CVODES.NV_Get_Ith_S(yB, 1));
|
---|
| 366 | CVODES.NV_Set_Ith_S(yBdot, 1, CVODES.NV_Get_Ith_S(yB, 0));
|
---|
[16245] | 367 |
|
---|
| 368 | return 0;
|
---|
| 369 | ;
|
---|
| 370 | }
|
---|
| 371 |
|
---|
| 372 | private int FQ(double t, IntPtr y, IntPtr yQdot, IntPtr user_data) {
|
---|
| 373 | // TODO: squared error
|
---|
[16248] | 374 | CVODES.NV_Set_Ith_S(yQdot, 0, CVODES.NV_Get_Ith_S(y, 2));
|
---|
[16245] | 375 | return 0;
|
---|
| 376 | }
|
---|
| 377 |
|
---|
| 378 |
|
---|
| 379 | private int FQB(double t, IntPtr y, IntPtr yB, IntPtr qBdot, IntPtr user_data) {
|
---|
| 380 | throw new NotImplementedException();
|
---|
| 381 | }
|
---|
[16225] | 382 | }
|
---|
| 383 | }
|
---|