1c4762a1bSJed Brown static char help[] = "Single-DOF oscillator formulated as a second-order system.\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscts.h> 4c4762a1bSJed Brown 5c4762a1bSJed Brown typedef struct { 6c4762a1bSJed Brown PetscReal Omega; /* natural frequency */ 7c4762a1bSJed Brown PetscReal Xi; /* damping coefficient */ 8c4762a1bSJed Brown PetscReal u0, v0; /* initial conditions */ 9220f924aSDavid Kamensky PetscBool use_pred; /* whether to use a predictor callback */ 10c4762a1bSJed Brown } UserParams; 11c4762a1bSJed Brown 12d71ae5a4SJacob Faibussowitsch static void Exact(PetscReal t, PetscReal omega, PetscReal xi, PetscReal u0, PetscReal v0, PetscReal *ut, PetscReal *vt) 13d71ae5a4SJacob Faibussowitsch { 14c4762a1bSJed Brown PetscReal u, v; 15c4762a1bSJed Brown if (xi < 1) { 16c4762a1bSJed Brown PetscReal a = xi * omega; 17303a5415SBarry Smith PetscReal w = PetscSqrtReal(1 - xi * xi) * omega; 18c4762a1bSJed Brown PetscReal C1 = (v0 + a * u0) / w; 19c4762a1bSJed Brown PetscReal C2 = u0; 20303a5415SBarry Smith u = PetscExpReal(-a * t) * (C1 * PetscSinReal(w * t) + C2 * PetscCosReal(w * t)); 21303a5415SBarry Smith v = (-a * PetscExpReal(-a * t) * (C1 * PetscSinReal(w * t) + C2 * PetscCosReal(w * t)) + w * PetscExpReal(-a * t) * (C1 * PetscCosReal(w * t) - C2 * PetscSinReal(w * t))); 22c4762a1bSJed Brown } else if (xi > 1) { 23303a5415SBarry Smith PetscReal w = PetscSqrtReal(xi * xi - 1) * omega; 24c4762a1bSJed Brown PetscReal C1 = (w * u0 + xi * u0 + v0) / (2 * w); 25c4762a1bSJed Brown PetscReal C2 = (w * u0 - xi * u0 - v0) / (2 * w); 26303a5415SBarry Smith u = C1 * PetscExpReal((-xi + w) * t) + C2 * PetscExpReal((-xi - w) * t); 27303a5415SBarry Smith v = C1 * (-xi + w) * PetscExpReal((-xi + w) * t) + C2 * (-xi - w) * PetscExpReal((-xi - w) * t); 28c4762a1bSJed Brown } else { 29c4762a1bSJed Brown PetscReal a = xi * omega; 30c4762a1bSJed Brown PetscReal C1 = v0 + a * u0; 31c4762a1bSJed Brown PetscReal C2 = u0; 32303a5415SBarry Smith u = (C1 * t + C2) * PetscExpReal(-a * t); 33303a5415SBarry Smith v = (C1 - a * (C1 * t + C2)) * PetscExpReal(-a * t); 34c4762a1bSJed Brown } 35c4762a1bSJed Brown if (ut) *ut = u; 36c4762a1bSJed Brown if (vt) *vt = v; 37c4762a1bSJed Brown } 38c4762a1bSJed Brown 39*2a8381b2SBarry Smith PetscErrorCode Solution(TS ts, PetscReal t, Vec X, PetscCtx ctx) 40d71ae5a4SJacob Faibussowitsch { 41c4762a1bSJed Brown UserParams *user = (UserParams *)ctx; 42c4762a1bSJed Brown PetscReal u, v; 43c4762a1bSJed Brown PetscScalar *x; 44c4762a1bSJed Brown 457510d9b0SBarry Smith PetscFunctionBeginUser; 46c4762a1bSJed Brown Exact(t, user->Omega, user->Xi, user->u0, user->v0, &u, &v); 479566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x)); 48c4762a1bSJed Brown x[0] = (PetscScalar)u; 499566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x)); 503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 51c4762a1bSJed Brown } 52c4762a1bSJed Brown 53*2a8381b2SBarry Smith PetscErrorCode Predictor(TS ts, Vec X0, Vec V0, Vec A0, Vec X1, PetscCtx ctx) 54220f924aSDavid Kamensky { 55220f924aSDavid Kamensky PetscReal dt, accel_fac; 56220f924aSDavid Kamensky 57220f924aSDavid Kamensky PetscFunctionBeginUser; 58220f924aSDavid Kamensky PetscCall(TSGetTimeStep(ts, &dt)); 59220f924aSDavid Kamensky accel_fac = 0.5 * dt * dt; 60220f924aSDavid Kamensky /* X1 = X0 + dt*V0 + accel_fac*A0 */ 61220f924aSDavid Kamensky PetscCall(VecCopy(X0, X1)); 62220f924aSDavid Kamensky PetscCall(VecAXPY(X1, dt, V0)); 63220f924aSDavid Kamensky PetscCall(VecAXPY(X1, accel_fac, A0)); 64220f924aSDavid Kamensky PetscFunctionReturn(PETSC_SUCCESS); 65220f924aSDavid Kamensky } 66220f924aSDavid Kamensky 67*2a8381b2SBarry Smith PetscErrorCode Residual1(TS ts, PetscReal t, Vec U, Vec A, Vec R, PetscCtx ctx) 68d71ae5a4SJacob Faibussowitsch { 69c4762a1bSJed Brown UserParams *user = (UserParams *)ctx; 70c4762a1bSJed Brown PetscReal Omega = user->Omega; 71c4762a1bSJed Brown const PetscScalar *u, *a; 72c4762a1bSJed Brown PetscScalar *r; 73c4762a1bSJed Brown 747510d9b0SBarry Smith PetscFunctionBeginUser; 759566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 769566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(A, &a)); 779566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(R, &r)); 78c4762a1bSJed Brown 79c4762a1bSJed Brown r[0] = a[0] + (Omega * Omega) * u[0]; 80c4762a1bSJed Brown 819566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 829566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(A, &a)); 839566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(R, &r)); 849566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(R)); 859566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(R)); 863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 87c4762a1bSJed Brown } 88c4762a1bSJed Brown 89*2a8381b2SBarry Smith PetscErrorCode Tangent1(TS ts, PetscReal t, Vec U, Vec A, PetscReal shiftA, Mat J, Mat P, PetscCtx ctx) 90d71ae5a4SJacob Faibussowitsch { 91c4762a1bSJed Brown UserParams *user = (UserParams *)ctx; 92c4762a1bSJed Brown PetscReal Omega = user->Omega; 93c4762a1bSJed Brown PetscReal T = 0; 94c4762a1bSJed Brown 957510d9b0SBarry Smith PetscFunctionBeginUser; 96c4762a1bSJed Brown T = shiftA + (Omega * Omega); 97c4762a1bSJed Brown 989566063dSJacob Faibussowitsch PetscCall(MatSetValue(P, 0, 0, T, INSERT_VALUES)); 999566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(P, MAT_FINAL_ASSEMBLY)); 1009566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(P, MAT_FINAL_ASSEMBLY)); 101c4762a1bSJed Brown if (J != P) { 1029566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 1039566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 104c4762a1bSJed Brown } 1053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 106c4762a1bSJed Brown } 107c4762a1bSJed Brown 108*2a8381b2SBarry Smith PetscErrorCode Residual2(TS ts, PetscReal t, Vec U, Vec V, Vec A, Vec R, PetscCtx ctx) 109d71ae5a4SJacob Faibussowitsch { 110c4762a1bSJed Brown UserParams *user = (UserParams *)ctx; 111c4762a1bSJed Brown PetscReal Omega = user->Omega, Xi = user->Xi; 112c4762a1bSJed Brown const PetscScalar *u, *v, *a; 113c4762a1bSJed Brown PetscScalar *r; 114c4762a1bSJed Brown 1157510d9b0SBarry Smith PetscFunctionBeginUser; 1169566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1179566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(V, &v)); 1189566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(A, &a)); 1199566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(R, &r)); 120c4762a1bSJed Brown 121c4762a1bSJed Brown r[0] = a[0] + (2 * Xi * Omega) * v[0] + (Omega * Omega) * u[0]; 122c4762a1bSJed Brown 1239566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1249566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(V, &v)); 1259566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(A, &a)); 1269566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(R, &r)); 1279566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(R)); 1289566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(R)); 1293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 130c4762a1bSJed Brown } 131c4762a1bSJed Brown 132*2a8381b2SBarry Smith PetscErrorCode Tangent2(TS ts, PetscReal t, Vec U, Vec V, Vec A, PetscReal shiftV, PetscReal shiftA, Mat J, Mat P, PetscCtx ctx) 133d71ae5a4SJacob Faibussowitsch { 134c4762a1bSJed Brown UserParams *user = (UserParams *)ctx; 135c4762a1bSJed Brown PetscReal Omega = user->Omega, Xi = user->Xi; 136c4762a1bSJed Brown PetscReal T = 0; 137c4762a1bSJed Brown 1387510d9b0SBarry Smith PetscFunctionBeginUser; 139c4762a1bSJed Brown T = shiftA + shiftV * (2 * Xi * Omega) + (Omega * Omega); 140c4762a1bSJed Brown 1419566063dSJacob Faibussowitsch PetscCall(MatSetValue(P, 0, 0, T, INSERT_VALUES)); 1429566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(P, MAT_FINAL_ASSEMBLY)); 1439566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(P, MAT_FINAL_ASSEMBLY)); 144c4762a1bSJed Brown if (J != P) { 1459566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 1469566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 147c4762a1bSJed Brown } 1483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 149c4762a1bSJed Brown } 150c4762a1bSJed Brown 151d71ae5a4SJacob Faibussowitsch int main(int argc, char *argv[]) 152d71ae5a4SJacob Faibussowitsch { 153c4762a1bSJed Brown PetscMPIInt size; 154c4762a1bSJed Brown TS ts; 155c4762a1bSJed Brown Vec R; 156c4762a1bSJed Brown Mat J; 157c4762a1bSJed Brown Vec U, V; 158c4762a1bSJed Brown PetscScalar *u, *v; 159220f924aSDavid Kamensky UserParams user = {/*Omega=*/1, /*Xi=*/0, /*u0=*/1, /*,v0=*/0, /*,use_pred=*/PETSC_FALSE}; 160c4762a1bSJed Brown 161327415f7SBarry Smith PetscFunctionBeginUser; 1629566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 1639566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 1643c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs"); 165c4762a1bSJed Brown 166d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_SELF, "", "ex43 options", ""); 1679566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-frequency", "Natual frequency", __FILE__, user.Omega, &user.Omega, NULL)); 1689566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-damping", "Damping coefficient", __FILE__, user.Xi, &user.Xi, NULL)); 1699566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-initial_u", "Initial displacement", __FILE__, user.u0, &user.u0, NULL)); 1709566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-initial_v", "Initial velocity", __FILE__, user.v0, &user.v0, NULL)); 171220f924aSDavid Kamensky PetscCall(PetscOptionsBool("-use_pred", "Use a custom predictor", __FILE__, user.use_pred, &user.use_pred, NULL)); 172d0609cedSBarry Smith PetscOptionsEnd(); 173c4762a1bSJed Brown 1749566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_SELF, &ts)); 1759566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSALPHA2)); 1769566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 5 * (2 * PETSC_PI))); 1779566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 1789566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 0.01)); 179c4762a1bSJed Brown 1809566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, 1, &R)); 1819566063dSJacob Faibussowitsch PetscCall(VecSetUp(R)); 1829566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, 1, 1, NULL, &J)); 1839566063dSJacob Faibussowitsch PetscCall(MatSetUp(J)); 184c4762a1bSJed Brown if (user.Xi) { 1859566063dSJacob Faibussowitsch PetscCall(TSSetI2Function(ts, R, Residual2, &user)); 1869566063dSJacob Faibussowitsch PetscCall(TSSetI2Jacobian(ts, J, J, Tangent2, &user)); 187c4762a1bSJed Brown } else { 1889566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, R, Residual1, &user)); 1899566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, J, J, Tangent1, &user)); 190c4762a1bSJed Brown } 1919566063dSJacob Faibussowitsch PetscCall(VecDestroy(&R)); 1929566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 1939566063dSJacob Faibussowitsch PetscCall(TSSetSolutionFunction(ts, Solution, &user)); 194c4762a1bSJed Brown 1959566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, 1, &U)); 1969566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, 1, &V)); 1979566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(U, &u)); 1989566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(V, &v)); 199c4762a1bSJed Brown u[0] = user.u0; 200c4762a1bSJed Brown v[0] = user.v0; 2019566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(U, &u)); 2029566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(V, &v)); 203c4762a1bSJed Brown 204220f924aSDavid Kamensky if (user.use_pred) PetscCall(TSAlpha2SetPredictor(ts, Predictor, NULL)); 205220f924aSDavid Kamensky 2069566063dSJacob Faibussowitsch PetscCall(TS2SetSolution(ts, U, V)); 2079566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 2089566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, NULL)); 209c4762a1bSJed Brown 2109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 2119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&V)); 2129566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 2139566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 214b122ec5aSJacob Faibussowitsch return 0; 215c4762a1bSJed Brown } 216c4762a1bSJed Brown 217c4762a1bSJed Brown /*TEST 218c4762a1bSJed Brown 219c4762a1bSJed Brown test: 220c4762a1bSJed Brown suffix: a 221220f924aSDavid Kamensky args: -ts_max_steps 10 -ts_view -use_pred 222c4762a1bSJed Brown requires: !single 223c4762a1bSJed Brown 224c4762a1bSJed Brown test: 225c4762a1bSJed Brown suffix: b 226c4762a1bSJed Brown args: -ts_max_steps 10 -ts_rtol 0 -ts_atol 1e-5 -ts_adapt_type basic -ts_adapt_monitor 227c4762a1bSJed Brown requires: !single 228c4762a1bSJed Brown 229c4762a1bSJed Brown TEST*/ 230