1*414ab8ebSDavid Kamensky static char help[] = "Constant velocity check with 1st-order generalized-alpha.\n"; 2*414ab8ebSDavid Kamensky 3*414ab8ebSDavid Kamensky #include <petscts.h> 4*414ab8ebSDavid Kamensky 5*414ab8ebSDavid Kamensky typedef struct { 6*414ab8ebSDavid Kamensky PetscReal v0; /* contant velocity */ 7*414ab8ebSDavid Kamensky PetscReal u0; /* initial condition */ 8*414ab8ebSDavid Kamensky PetscReal radius; /* spectral radius of integrator */ 9*414ab8ebSDavid Kamensky } UserParams; 10*414ab8ebSDavid Kamensky 11*414ab8ebSDavid Kamensky static void Exact(PetscReal t, PetscReal v0, PetscReal u0, PetscScalar *ut) 12*414ab8ebSDavid Kamensky { 13*414ab8ebSDavid Kamensky if (ut) *ut = u0 + v0 * t; 14*414ab8ebSDavid Kamensky } 15*414ab8ebSDavid Kamensky 16*414ab8ebSDavid Kamensky PetscErrorCode Residual(TS ts, PetscReal t, Vec U, Vec V, Vec R, void *ctx) 17*414ab8ebSDavid Kamensky { 18*414ab8ebSDavid Kamensky UserParams *user = (UserParams *)ctx; 19*414ab8ebSDavid Kamensky const PetscScalar *v; 20*414ab8ebSDavid Kamensky PetscScalar *r; 21*414ab8ebSDavid Kamensky 22*414ab8ebSDavid Kamensky PetscFunctionBeginUser; 23*414ab8ebSDavid Kamensky PetscCall(VecGetArrayRead(V, &v)); 24*414ab8ebSDavid Kamensky PetscCall(VecGetArrayWrite(R, &r)); 25*414ab8ebSDavid Kamensky 26*414ab8ebSDavid Kamensky r[0] = v[0] - user->v0; 27*414ab8ebSDavid Kamensky 28*414ab8ebSDavid Kamensky PetscCall(VecRestoreArrayRead(V, &v)); 29*414ab8ebSDavid Kamensky PetscCall(VecRestoreArrayWrite(R, &r)); 30*414ab8ebSDavid Kamensky PetscFunctionReturn(PETSC_SUCCESS); 31*414ab8ebSDavid Kamensky } 32*414ab8ebSDavid Kamensky 33*414ab8ebSDavid Kamensky PetscErrorCode Tangent(TS ts, PetscReal t, Vec U, Vec V, PetscReal shiftV, Mat J, Mat P, void *ctx) 34*414ab8ebSDavid Kamensky { 35*414ab8ebSDavid Kamensky PetscReal T = 0; 36*414ab8ebSDavid Kamensky 37*414ab8ebSDavid Kamensky PetscFunctionBeginUser; 38*414ab8ebSDavid Kamensky T = shiftV; 39*414ab8ebSDavid Kamensky 40*414ab8ebSDavid Kamensky PetscCall(MatSetValue(P, 0, 0, T, INSERT_VALUES)); 41*414ab8ebSDavid Kamensky PetscCall(MatAssemblyBegin(P, MAT_FINAL_ASSEMBLY)); 42*414ab8ebSDavid Kamensky PetscCall(MatAssemblyEnd(P, MAT_FINAL_ASSEMBLY)); 43*414ab8ebSDavid Kamensky if (J != P) { 44*414ab8ebSDavid Kamensky PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 45*414ab8ebSDavid Kamensky PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 46*414ab8ebSDavid Kamensky } 47*414ab8ebSDavid Kamensky PetscFunctionReturn(PETSC_SUCCESS); 48*414ab8ebSDavid Kamensky } 49*414ab8ebSDavid Kamensky 50*414ab8ebSDavid Kamensky int main(int argc, char *argv[]) 51*414ab8ebSDavid Kamensky { 52*414ab8ebSDavid Kamensky PetscMPIInt size; 53*414ab8ebSDavid Kamensky TS ts; 54*414ab8ebSDavid Kamensky Vec R; 55*414ab8ebSDavid Kamensky Mat J; 56*414ab8ebSDavid Kamensky Vec U; 57*414ab8ebSDavid Kamensky PetscScalar *u, u_exact; 58*414ab8ebSDavid Kamensky PetscReal u_err; 59*414ab8ebSDavid Kamensky PetscReal atol = 1e-15; 60*414ab8ebSDavid Kamensky PetscReal t_final = 3.0; 61*414ab8ebSDavid Kamensky PetscInt n_step = 8; 62*414ab8ebSDavid Kamensky UserParams user = {/*v0=*/1, /*u0=*/1, /*radius=*/0.0}; 63*414ab8ebSDavid Kamensky 64*414ab8ebSDavid Kamensky PetscFunctionBeginUser; 65*414ab8ebSDavid Kamensky PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 66*414ab8ebSDavid Kamensky PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 67*414ab8ebSDavid Kamensky PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs"); 68*414ab8ebSDavid Kamensky 69*414ab8ebSDavid Kamensky PetscOptionsBegin(PETSC_COMM_SELF, "", "ex81 options", ""); 70*414ab8ebSDavid Kamensky PetscCall(PetscOptionsReal("-velocity", "Constant velocity", __FILE__, user.v0, &user.v0, NULL)); 71*414ab8ebSDavid Kamensky PetscCall(PetscOptionsReal("-initial_u", "Initial displacement", __FILE__, user.u0, &user.u0, NULL)); 72*414ab8ebSDavid Kamensky PetscCall(PetscOptionsReal("-radius", "Spectral radius", __FILE__, user.radius, &user.radius, NULL)); 73*414ab8ebSDavid Kamensky PetscOptionsEnd(); 74*414ab8ebSDavid Kamensky 75*414ab8ebSDavid Kamensky PetscCall(TSCreate(PETSC_COMM_SELF, &ts)); 76*414ab8ebSDavid Kamensky PetscCall(TSSetType(ts, TSALPHA)); 77*414ab8ebSDavid Kamensky PetscCall(TSSetMaxTime(ts, t_final)); 78*414ab8ebSDavid Kamensky PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 79*414ab8ebSDavid Kamensky PetscCall(TSSetTimeStep(ts, t_final / n_step)); 80*414ab8ebSDavid Kamensky PetscCall(TSAlphaSetRadius(ts, user.radius)); 81*414ab8ebSDavid Kamensky 82*414ab8ebSDavid Kamensky PetscCall(VecCreateSeq(PETSC_COMM_SELF, 1, &R)); 83*414ab8ebSDavid Kamensky PetscCall(VecSetUp(R)); 84*414ab8ebSDavid Kamensky PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, 1, 1, NULL, &J)); 85*414ab8ebSDavid Kamensky PetscCall(MatSetUp(J)); 86*414ab8ebSDavid Kamensky PetscCall(TSSetIFunction(ts, R, Residual, &user)); 87*414ab8ebSDavid Kamensky PetscCall(TSSetIJacobian(ts, J, J, Tangent, &user)); 88*414ab8ebSDavid Kamensky PetscCall(VecDestroy(&R)); 89*414ab8ebSDavid Kamensky PetscCall(MatDestroy(&J)); 90*414ab8ebSDavid Kamensky 91*414ab8ebSDavid Kamensky PetscCall(VecCreateSeq(PETSC_COMM_SELF, 1, &U)); 92*414ab8ebSDavid Kamensky PetscCall(VecGetArrayWrite(U, &u)); 93*414ab8ebSDavid Kamensky u[0] = user.u0; 94*414ab8ebSDavid Kamensky PetscCall(VecRestoreArrayWrite(U, &u)); 95*414ab8ebSDavid Kamensky 96*414ab8ebSDavid Kamensky PetscCall(TSSetSolution(ts, U)); 97*414ab8ebSDavid Kamensky PetscCall(TSSetFromOptions(ts)); 98*414ab8ebSDavid Kamensky PetscCall(TSSolve(ts, NULL)); 99*414ab8ebSDavid Kamensky 100*414ab8ebSDavid Kamensky PetscCall(VecGetArray(U, &u)); 101*414ab8ebSDavid Kamensky Exact(t_final, user.v0, user.u0, &u_exact); 102*414ab8ebSDavid Kamensky u_err = PetscAbsScalar(u[0] - u_exact); 103*414ab8ebSDavid Kamensky PetscCall(PetscPrintf(PETSC_COMM_WORLD, "u(t=%g) = %g (error = %g)\n", (double)t_final, (double)PetscRealPart(u[0]), (double)u_err)); 104*414ab8ebSDavid Kamensky PetscCheck(u_err < atol, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Inexact displacement."); 105*414ab8ebSDavid Kamensky PetscCall(VecRestoreArray(U, &u)); 106*414ab8ebSDavid Kamensky 107*414ab8ebSDavid Kamensky PetscCall(VecDestroy(&U)); 108*414ab8ebSDavid Kamensky PetscCall(TSDestroy(&ts)); 109*414ab8ebSDavid Kamensky PetscCall(PetscFinalize()); 110*414ab8ebSDavid Kamensky return 0; 111*414ab8ebSDavid Kamensky } 112*414ab8ebSDavid Kamensky 113*414ab8ebSDavid Kamensky /*TEST 114*414ab8ebSDavid Kamensky 115*414ab8ebSDavid Kamensky test: 116*414ab8ebSDavid Kamensky suffix: a 117*414ab8ebSDavid Kamensky args: -radius 0.0 118*414ab8ebSDavid Kamensky requires: !single 119*414ab8ebSDavid Kamensky 120*414ab8ebSDavid Kamensky test: 121*414ab8ebSDavid Kamensky suffix: b 122*414ab8ebSDavid Kamensky args: -radius 0.5 123*414ab8ebSDavid Kamensky requires: !single 124*414ab8ebSDavid Kamensky 125*414ab8ebSDavid Kamensky test: 126*414ab8ebSDavid Kamensky suffix: c 127*414ab8ebSDavid Kamensky args: -radius 1.0 128*414ab8ebSDavid Kamensky requires: !single 129*414ab8ebSDavid Kamensky 130*414ab8ebSDavid Kamensky TEST*/ 131