1*d27334e2SStefano Zampini static char help[] = "Solves a DAE with a non-trivial mass matrix. \n\n"; 2*d27334e2SStefano Zampini /* 3*d27334e2SStefano Zampini Solves: 4*d27334e2SStefano Zampini U * dU/dt = U*V 5*d27334e2SStefano Zampini U - V = 0 6*d27334e2SStefano Zampini 7*d27334e2SStefano Zampini that can be rewritten in implicit form F = 0, with F as 8*d27334e2SStefano Zampini x[0] * xdot[0] - x[0] * x[1] 9*d27334e2SStefano Zampini F(t,x,xdot) = 10*d27334e2SStefano Zampini x[0] - x[1] 11*d27334e2SStefano Zampini It is equivalent to solve dU/dt = U, U = U0 with solution U = U0 * exp(tfinal) 12*d27334e2SStefano Zampini */ 13*d27334e2SStefano Zampini 14*d27334e2SStefano Zampini #include <petscts.h> 15*d27334e2SStefano Zampini 16*d27334e2SStefano Zampini PetscErrorCode IFunction(TS, PetscReal, Vec, Vec, Vec, void *); 17*d27334e2SStefano Zampini PetscErrorCode IJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *); 18*d27334e2SStefano Zampini 19*d27334e2SStefano Zampini int main(int argc, char **argv) 20*d27334e2SStefano Zampini { 21*d27334e2SStefano Zampini TS ts; 22*d27334e2SStefano Zampini Vec x; 23*d27334e2SStefano Zampini PetscBool dae = PETSC_TRUE; 24*d27334e2SStefano Zampini 25*d27334e2SStefano Zampini PetscFunctionBeginUser; 26*d27334e2SStefano Zampini PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 27*d27334e2SStefano Zampini PetscCall(PetscOptionsGetBool(NULL, NULL, "-dae", &dae, NULL)); 28*d27334e2SStefano Zampini 29*d27334e2SStefano Zampini PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 30*d27334e2SStefano Zampini PetscCall(TSSetIFunction(ts, NULL, IFunction, &dae)); 31*d27334e2SStefano Zampini PetscCall(TSSetIJacobian(ts, NULL, NULL, IJacobian, &dae)); 32*d27334e2SStefano Zampini 33*d27334e2SStefano Zampini PetscCall(VecCreate(PETSC_COMM_WORLD, &x)); 34*d27334e2SStefano Zampini PetscCall(VecSetSizes(x, 2, PETSC_DECIDE)); 35*d27334e2SStefano Zampini PetscCall(VecSetFromOptions(x)); 36*d27334e2SStefano Zampini PetscCall(VecSetUp(x)); 37*d27334e2SStefano Zampini PetscCall(VecSet(x, 0.5)); 38*d27334e2SStefano Zampini PetscCall(TSSetSolution(ts, x)); 39*d27334e2SStefano Zampini PetscCall(VecDestroy(&x)); 40*d27334e2SStefano Zampini 41*d27334e2SStefano Zampini PetscCall(TSSetTimeStep(ts, 1.0 / 16.0)); 42*d27334e2SStefano Zampini PetscCall(TSSetMaxTime(ts, 1.0)); 43*d27334e2SStefano Zampini PetscCall(TSSetFromOptions(ts)); 44*d27334e2SStefano Zampini PetscCall(TSSolve(ts, NULL)); 45*d27334e2SStefano Zampini 46*d27334e2SStefano Zampini PetscCall(TSDestroy(&ts)); 47*d27334e2SStefano Zampini PetscCall(PetscFinalize()); 48*d27334e2SStefano Zampini return 0; 49*d27334e2SStefano Zampini } 50*d27334e2SStefano Zampini 51*d27334e2SStefano Zampini PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx) 52*d27334e2SStefano Zampini { 53*d27334e2SStefano Zampini const PetscScalar *xdot, *x; 54*d27334e2SStefano Zampini PetscScalar *f; 55*d27334e2SStefano Zampini PetscBool dae = *(PetscBool *)(ctx); 56*d27334e2SStefano Zampini 57*d27334e2SStefano Zampini PetscFunctionBeginUser; 58*d27334e2SStefano Zampini PetscCall(VecGetArrayRead(Xdot, &xdot)); 59*d27334e2SStefano Zampini PetscCall(VecGetArrayRead(X, &x)); 60*d27334e2SStefano Zampini PetscCall(VecGetArrayWrite(F, &f)); 61*d27334e2SStefano Zampini if (dae) { 62*d27334e2SStefano Zampini f[0] = x[0] * xdot[0] - x[0] * x[1]; 63*d27334e2SStefano Zampini f[1] = x[0] - x[1]; 64*d27334e2SStefano Zampini } else { 65*d27334e2SStefano Zampini f[0] = xdot[0] - x[0]; 66*d27334e2SStefano Zampini f[1] = xdot[1] - x[1]; 67*d27334e2SStefano Zampini } 68*d27334e2SStefano Zampini PetscCall(VecRestoreArrayRead(Xdot, &xdot)); 69*d27334e2SStefano Zampini PetscCall(VecRestoreArrayRead(X, &x)); 70*d27334e2SStefano Zampini PetscCall(VecRestoreArrayWrite(F, &f)); 71*d27334e2SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 72*d27334e2SStefano Zampini } 73*d27334e2SStefano Zampini 74*d27334e2SStefano Zampini PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal shift, Mat A, Mat B, void *ctx) 75*d27334e2SStefano Zampini { 76*d27334e2SStefano Zampini const PetscScalar *xdot, *x; 77*d27334e2SStefano Zampini PetscBool dae = *(PetscBool *)(ctx); 78*d27334e2SStefano Zampini 79*d27334e2SStefano Zampini PetscFunctionBeginUser; 80*d27334e2SStefano Zampini PetscCall(VecGetArrayRead(Xdot, &xdot)); 81*d27334e2SStefano Zampini PetscCall(VecGetArrayRead(X, &x)); 82*d27334e2SStefano Zampini if (dae) { 83*d27334e2SStefano Zampini PetscCall(MatSetValue(B, 0, 0, shift * x[0] + xdot[0] - x[1], INSERT_VALUES)); 84*d27334e2SStefano Zampini PetscCall(MatSetValue(B, 0, 1, -x[0], INSERT_VALUES)); 85*d27334e2SStefano Zampini PetscCall(MatSetValue(B, 1, 0, 1.0, INSERT_VALUES)); 86*d27334e2SStefano Zampini PetscCall(MatSetValue(B, 1, 1, -1.0, INSERT_VALUES)); 87*d27334e2SStefano Zampini } else { 88*d27334e2SStefano Zampini PetscCall(MatZeroEntries(B)); 89*d27334e2SStefano Zampini PetscCall(MatShift(B, shift)); 90*d27334e2SStefano Zampini } 91*d27334e2SStefano Zampini PetscCall(VecRestoreArrayRead(X, &x)); 92*d27334e2SStefano Zampini PetscCall(VecRestoreArrayRead(Xdot, &xdot)); 93*d27334e2SStefano Zampini PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 94*d27334e2SStefano Zampini PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 95*d27334e2SStefano Zampini if (A != B) { 96*d27334e2SStefano Zampini PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 97*d27334e2SStefano Zampini PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 98*d27334e2SStefano Zampini } 99*d27334e2SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 100*d27334e2SStefano Zampini } 101*d27334e2SStefano Zampini 102*d27334e2SStefano Zampini /*TEST 103*d27334e2SStefano Zampini 104*d27334e2SStefano Zampini testset: 105*d27334e2SStefano Zampini args: -ts_view_solution -ts_max_steps 10 -ts_dt 0.1 -ts_view_solution -ts_adapt_type {{none basic}} -ts_exact_final_time matchstep 106*d27334e2SStefano Zampini 107*d27334e2SStefano Zampini test: 108*d27334e2SStefano Zampini output_file: output/ex18_1.out 109*d27334e2SStefano Zampini suffix: bdf 110*d27334e2SStefano Zampini args: -ts_type bdf 111*d27334e2SStefano Zampini 112*d27334e2SStefano Zampini test: 113*d27334e2SStefano Zampini output_file: output/ex18_1.out 114*d27334e2SStefano Zampini suffix: dirk 115*d27334e2SStefano Zampini args: -dae {{0 1}} -ts_type dirk -ts_dirk_type {{s212 es212 es213sal es324sal es325sal 657a es648sa 658a s659a 7510sal es7510sa 759a s7511sal 8614a 8616sal es8516sal}} 116*d27334e2SStefano Zampini 117*d27334e2SStefano Zampini TEST*/ 118