xref: /petsc/src/ts/tests/ex3.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
1c4762a1bSJed Brown static char help[] = "Solves 1D heat equation with FEM formulation.\n\
2c4762a1bSJed Brown Input arguments are\n\
3c4762a1bSJed Brown   -useAlhs: solve Alhs*U' =  (Arhs*U + g) \n\
4c4762a1bSJed Brown             otherwise, solve U' = inv(Alhs)*(Arhs*U + g) \n\n";
5c4762a1bSJed Brown 
6c4762a1bSJed Brown /*--------------------------------------------------------------------------
7c4762a1bSJed Brown   Solves 1D heat equation U_t = U_xx with FEM formulation:
8c4762a1bSJed Brown                           Alhs*U' = rhs (= Arhs*U + g)
9c4762a1bSJed Brown   We thank Chris Cox <clcox@clemson.edu> for contributing the original code
10c4762a1bSJed Brown ----------------------------------------------------------------------------*/
11c4762a1bSJed Brown 
12c4762a1bSJed Brown #include <petscksp.h>
13c4762a1bSJed Brown #include <petscts.h>
14c4762a1bSJed Brown 
15c4762a1bSJed Brown /* special variable - max size of all arrays  */
16c4762a1bSJed Brown #define num_z 10
17c4762a1bSJed Brown 
18c4762a1bSJed Brown /*
19c4762a1bSJed Brown    User-defined application context - contains data needed by the
20c4762a1bSJed Brown    application-provided call-back routines.
21c4762a1bSJed Brown */
22c4762a1bSJed Brown typedef struct {
23c4762a1bSJed Brown   Mat          Amat;             /* left hand side matrix */
24c4762a1bSJed Brown   Vec          ksp_rhs, ksp_sol; /* working vectors for formulating inv(Alhs)*(Arhs*U+g) */
256497c311SBarry Smith   PetscInt     max_probsz;       /* max size of the problem */
26c4762a1bSJed Brown   PetscBool    useAlhs;          /* flag (1 indicates solving Alhs*U' = Arhs*U+g */
276497c311SBarry Smith   PetscInt     nz;               /* total number of grid points */
28c4762a1bSJed Brown   PetscInt     m;                /* total number of interio grid points */
29c4762a1bSJed Brown   Vec          solution;         /* global exact ts solution vector */
30c4762a1bSJed Brown   PetscScalar *z;                /* array of grid points */
31c4762a1bSJed Brown   PetscBool    debug;            /* flag (1 indicates activation of debugging printouts) */
32c4762a1bSJed Brown } AppCtx;
33c4762a1bSJed Brown 
34c4762a1bSJed Brown extern PetscScalar    exact(PetscScalar, PetscReal);
35c4762a1bSJed Brown extern PetscErrorCode Monitor(TS, PetscInt, PetscReal, Vec, void *);
36c4762a1bSJed Brown extern PetscErrorCode Petsc_KSPSolve(AppCtx *);
37c4762a1bSJed Brown extern PetscScalar    bspl(PetscScalar *, PetscScalar, PetscInt, PetscInt, PetscInt[][2], PetscInt);
38c4762a1bSJed Brown extern PetscErrorCode femBg(PetscScalar[][3], PetscScalar *, PetscInt, PetscScalar *, PetscReal);
39c4762a1bSJed Brown extern PetscErrorCode femA(AppCtx *, PetscInt, PetscScalar *);
40c4762a1bSJed Brown extern PetscErrorCode rhs(AppCtx *, PetscScalar *, PetscInt, PetscScalar *, PetscReal);
41c4762a1bSJed Brown extern PetscErrorCode RHSfunction(TS, PetscReal, Vec, Vec, void *);
42c4762a1bSJed Brown 
main(int argc,char ** argv)43d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
44d71ae5a4SJacob Faibussowitsch {
45c4762a1bSJed Brown   PetscInt    i, m, nz, steps, max_steps, k, nphase = 1;
46c4762a1bSJed Brown   PetscScalar zInitial, zFinal, val, *z;
47c4762a1bSJed Brown   PetscReal   stepsz[4], T, ftime;
48c4762a1bSJed Brown   TS          ts;
49c4762a1bSJed Brown   SNES        snes;
50c4762a1bSJed Brown   Mat         Jmat;
51c4762a1bSJed Brown   AppCtx      appctx;   /* user-defined application context */
52c4762a1bSJed Brown   Vec         init_sol; /* ts solution vector */
53c4762a1bSJed Brown   PetscMPIInt size;
54c4762a1bSJed Brown 
55327415f7SBarry Smith   PetscFunctionBeginUser;
56c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
579566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
583c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_SELF, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only");
59c4762a1bSJed Brown 
60c4762a1bSJed Brown   /* initializations */
61c4762a1bSJed Brown   zInitial  = 0.0;
62c4762a1bSJed Brown   zFinal    = 1.0;
63c4762a1bSJed Brown   nz        = num_z;
64c4762a1bSJed Brown   m         = nz - 2;
65c4762a1bSJed Brown   appctx.nz = nz;
66300f1712SStefano Zampini   max_steps = 10000;
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   appctx.m          = m;
69c4762a1bSJed Brown   appctx.max_probsz = nz;
70c4762a1bSJed Brown   appctx.debug      = PETSC_FALSE;
71c4762a1bSJed Brown   appctx.useAlhs    = PETSC_FALSE;
72c4762a1bSJed Brown 
73d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "", "");
749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsName("-debug", NULL, NULL, &appctx.debug));
759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsName("-useAlhs", NULL, NULL, &appctx.useAlhs));
769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRangeInt("-nphase", NULL, NULL, nphase, &nphase, NULL, 1, 3));
77d0609cedSBarry Smith   PetscOptionsEnd();
78303a5415SBarry Smith   T = 0.014 / nphase;
79303a5415SBarry Smith 
80c4762a1bSJed Brown   /* create vector to hold ts solution */
81c4762a1bSJed Brown   /*-----------------------------------*/
829566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &init_sol));
839566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(init_sol, PETSC_DECIDE, m));
849566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(init_sol));
85c4762a1bSJed Brown 
86c4762a1bSJed Brown   /* create vector to hold true ts soln for comparison */
879566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(init_sol, &appctx.solution));
88c4762a1bSJed Brown 
89c4762a1bSJed Brown   /* create LHS matrix Amat */
90c4762a1bSJed Brown   /*------------------------*/
919566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD, m, m, 3, NULL, &appctx.Amat));
929566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(appctx.Amat));
939566063dSJacob Faibussowitsch   PetscCall(MatSetUp(appctx.Amat));
94c4762a1bSJed Brown   /* set space grid points - interio points only! */
959566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nz + 1, &z));
96c4762a1bSJed Brown   for (i = 0; i < nz; i++) z[i] = (i) * ((zFinal - zInitial) / (nz - 1));
97c4762a1bSJed Brown   appctx.z = z;
983ba16761SJacob Faibussowitsch   PetscCall(femA(&appctx, nz, z));
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   /* create the jacobian matrix */
101c4762a1bSJed Brown   /*----------------------------*/
1029566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &Jmat));
1039566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Jmat, PETSC_DECIDE, PETSC_DECIDE, m, m));
1049566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(Jmat));
1059566063dSJacob Faibussowitsch   PetscCall(MatSetUp(Jmat));
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   /* create working vectors for formulating rhs=inv(Alhs)*(Arhs*U + g) */
1089566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(init_sol, &appctx.ksp_rhs));
1099566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(init_sol, &appctx.ksp_sol));
110c4762a1bSJed Brown 
1112d4ee042Sprj-   /* set initial guess */
1122d4ee042Sprj-   /*-------------------*/
113c4762a1bSJed Brown   for (i = 0; i < nz - 2; i++) {
114c4762a1bSJed Brown     val = exact(z[i + 1], 0.0);
115300f1712SStefano Zampini     PetscCall(VecSetValue(init_sol, i, val, INSERT_VALUES));
116c4762a1bSJed Brown   }
1179566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(init_sol));
1189566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(init_sol));
119c4762a1bSJed Brown 
120c4762a1bSJed Brown   /*create a time-stepping context and set the problem type */
121c4762a1bSJed Brown   /*--------------------------------------------------------*/
1229566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1239566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   /* set time-step method */
1269566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSCN));
127c4762a1bSJed Brown 
128c4762a1bSJed Brown   /* Set optional user-defined monitoring routine */
1299566063dSJacob Faibussowitsch   PetscCall(TSMonitorSet(ts, Monitor, &appctx, NULL));
130dd8e379bSPierre Jolivet   /* set the right-hand side of U_t = RHSfunction(U,t) */
1319566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, NULL, (PetscErrorCode (*)(TS, PetscScalar, Vec, Vec, void *))RHSfunction, &appctx));
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   if (appctx.useAlhs) {
134c4762a1bSJed Brown     /* set the left hand side matrix of Amat*U_t = rhs(U,t) */
135c4762a1bSJed Brown 
136c4762a1bSJed Brown     /* Note: this approach is incompatible with the finite differenced Jacobian set below because we can't restore the
137c4762a1bSJed Brown      * Alhs matrix without making a copy.  Either finite difference the entire thing or use analytic Jacobians in both
138c4762a1bSJed Brown      * places.
139c4762a1bSJed Brown      */
1409566063dSJacob Faibussowitsch     PetscCall(TSSetIFunction(ts, NULL, TSComputeIFunctionLinear, &appctx));
1419566063dSJacob Faibussowitsch     PetscCall(TSSetIJacobian(ts, appctx.Amat, appctx.Amat, TSComputeIJacobianConstant, &appctx));
142c4762a1bSJed Brown   }
143c4762a1bSJed Brown 
144f0b74427SPierre Jolivet   /* use PETSc to compute the jacobian by finite differences */
1459566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
1469566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, Jmat, Jmat, SNESComputeJacobianDefault, NULL));
147c4762a1bSJed Brown 
148c4762a1bSJed Brown   /* get the command line options if there are any and set them */
1499566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
150c4762a1bSJed Brown 
151e808b789SPatrick Sanan #if defined(PETSC_HAVE_SUNDIALS2)
152c4762a1bSJed Brown   {
153c4762a1bSJed Brown     TSType    type;
154c4762a1bSJed Brown     PetscBool sundialstype = PETSC_FALSE;
1559566063dSJacob Faibussowitsch     PetscCall(TSGetType(ts, &type));
1569566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)ts, TSSUNDIALS, &sundialstype));
1573c633725SBarry Smith     PetscCheck(!sundialstype || !appctx.useAlhs, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot use Alhs formulation for TSSUNDIALS type");
158c4762a1bSJed Brown   }
159c4762a1bSJed Brown #endif
160c4762a1bSJed Brown   /* Sets the initial solution */
1619566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, init_sol));
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   stepsz[0] = 1.0 / (2.0 * (nz - 1) * (nz - 1)); /* (mesh_size)^2/2.0 */
164c4762a1bSJed Brown   ftime     = 0.0;
165c4762a1bSJed Brown   for (k = 0; k < nphase; k++) {
16663a3b9bcSJacob Faibussowitsch     if (nphase > 1) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Phase %" PetscInt_FMT " initial time %g, stepsz %g, duration: %g\n", k, (double)ftime, (double)stepsz[k], (double)((k + 1) * T)));
1679566063dSJacob Faibussowitsch     PetscCall(TSSetTime(ts, ftime));
1689566063dSJacob Faibussowitsch     PetscCall(TSSetTimeStep(ts, stepsz[k]));
1699566063dSJacob Faibussowitsch     PetscCall(TSSetMaxSteps(ts, max_steps));
1709566063dSJacob Faibussowitsch     PetscCall(TSSetMaxTime(ts, (k + 1) * T));
1719566063dSJacob Faibussowitsch     PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
172c4762a1bSJed Brown 
173c4762a1bSJed Brown     /* loop over time steps */
174c4762a1bSJed Brown     /*----------------------*/
1759566063dSJacob Faibussowitsch     PetscCall(TSSolve(ts, init_sol));
1769566063dSJacob Faibussowitsch     PetscCall(TSGetSolveTime(ts, &ftime));
1779566063dSJacob Faibussowitsch     PetscCall(TSGetStepNumber(ts, &steps));
178c4762a1bSJed Brown     stepsz[k + 1] = stepsz[k] * 1.5; /* change step size for the next phase */
179c4762a1bSJed Brown   }
180c4762a1bSJed Brown 
181c4762a1bSJed Brown   /* free space */
1829566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1839566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&appctx.Amat));
1849566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Jmat));
1859566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.ksp_rhs));
1869566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.ksp_sol));
1879566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&init_sol));
1889566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&appctx.solution));
1899566063dSJacob Faibussowitsch   PetscCall(PetscFree(z));
190c4762a1bSJed Brown 
1919566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
192b122ec5aSJacob Faibussowitsch   return 0;
193c4762a1bSJed Brown }
194c4762a1bSJed Brown 
195c4762a1bSJed Brown /*------------------------------------------------------------------------
196c4762a1bSJed Brown   Set exact solution
197c4762a1bSJed Brown   u(z,t) = sin(6*PI*z)*exp(-36.*PI*PI*t) + 3.*sin(2*PI*z)*exp(-4.*PI*PI*t)
198c4762a1bSJed Brown --------------------------------------------------------------------------*/
exact(PetscScalar z,PetscReal t)199d71ae5a4SJacob Faibussowitsch PetscScalar exact(PetscScalar z, PetscReal t)
200d71ae5a4SJacob Faibussowitsch {
201c4762a1bSJed Brown   PetscScalar val, ex1, ex2;
202c4762a1bSJed Brown 
203c4762a1bSJed Brown   ex1 = PetscExpReal(-36. * PETSC_PI * PETSC_PI * t);
204c4762a1bSJed Brown   ex2 = PetscExpReal(-4. * PETSC_PI * PETSC_PI * t);
205c4762a1bSJed Brown   val = PetscSinScalar(6 * PETSC_PI * z) * ex1 + 3. * PetscSinScalar(2 * PETSC_PI * z) * ex2;
206c4762a1bSJed Brown   return val;
207c4762a1bSJed Brown }
208c4762a1bSJed Brown 
209c4762a1bSJed Brown /*
210c4762a1bSJed Brown    Monitor - User-provided routine to monitor the solution computed at
211c4762a1bSJed Brown    each timestep.  This example plots the solution and computes the
212c4762a1bSJed Brown    error in two different norms.
213c4762a1bSJed Brown 
214c4762a1bSJed Brown    Input Parameters:
215c4762a1bSJed Brown    ts     - the timestep context
216c4762a1bSJed Brown    step   - the count of the current step (with 0 meaning the
217c4762a1bSJed Brown              initial condition)
218c4762a1bSJed Brown    time   - the current time
219c4762a1bSJed Brown    u      - the solution at this timestep
220c4762a1bSJed Brown    ctx    - the user-provided context for this monitoring routine.
221c4762a1bSJed Brown             In this case we use the application context which contains
222c4762a1bSJed Brown             information about the problem size, workspace and the exact
223c4762a1bSJed Brown             solution.
224c4762a1bSJed Brown */
Monitor(TS ts,PetscInt step,PetscReal time,Vec u,PetscCtx ctx)225*2a8381b2SBarry Smith PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal time, Vec u, PetscCtx ctx)
226d71ae5a4SJacob Faibussowitsch {
227c4762a1bSJed Brown   AppCtx      *appctx = (AppCtx *)ctx;
228c4762a1bSJed Brown   PetscInt     i, m = appctx->m;
229c4762a1bSJed Brown   PetscReal    norm_2, norm_max, h = 1.0 / (m + 1);
230c4762a1bSJed Brown   PetscScalar *u_exact;
231c4762a1bSJed Brown 
2323ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
233c4762a1bSJed Brown   /* Compute the exact solution */
2349566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(appctx->solution, &u_exact));
235c4762a1bSJed Brown   for (i = 0; i < m; i++) u_exact[i] = exact(appctx->z[i + 1], time);
2369566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(appctx->solution, &u_exact));
237c4762a1bSJed Brown 
238c4762a1bSJed Brown   /* Print debugging information if desired */
239c4762a1bSJed Brown   if (appctx->debug) {
2409566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Computed solution vector at time %g\n", (double)time));
2419566063dSJacob Faibussowitsch     PetscCall(VecView(u, PETSC_VIEWER_STDOUT_SELF));
2429566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Exact solution vector\n"));
2439566063dSJacob Faibussowitsch     PetscCall(VecView(appctx->solution, PETSC_VIEWER_STDOUT_SELF));
244c4762a1bSJed Brown   }
245c4762a1bSJed Brown 
246c4762a1bSJed Brown   /* Compute the 2-norm and max-norm of the error */
2479566063dSJacob Faibussowitsch   PetscCall(VecAXPY(appctx->solution, -1.0, u));
2489566063dSJacob Faibussowitsch   PetscCall(VecNorm(appctx->solution, NORM_2, &norm_2));
249c4762a1bSJed Brown 
250c4762a1bSJed Brown   norm_2 = PetscSqrtReal(h) * norm_2;
2519566063dSJacob Faibussowitsch   PetscCall(VecNorm(appctx->solution, NORM_MAX, &norm_max));
25263a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Timestep %" PetscInt_FMT ": time = %g, 2-norm error = %6.4f, max norm error = %6.4f\n", step, (double)time, (double)norm_2, (double)norm_max));
253c4762a1bSJed Brown 
254c4762a1bSJed Brown   /*
255c4762a1bSJed Brown      Print debugging information if desired
256c4762a1bSJed Brown   */
257c4762a1bSJed Brown   if (appctx->debug) {
2589566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error vector\n"));
2599566063dSJacob Faibussowitsch     PetscCall(VecView(appctx->solution, PETSC_VIEWER_STDOUT_SELF));
260c4762a1bSJed Brown   }
2613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
262c4762a1bSJed Brown }
263c4762a1bSJed Brown 
264c4762a1bSJed Brown /*%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
2650e3d61c9SBarry Smith       Function to solve a linear system using KSP
266c4762a1bSJed Brown %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%*/
267c4762a1bSJed Brown 
Petsc_KSPSolve(AppCtx * obj)268d71ae5a4SJacob Faibussowitsch PetscErrorCode Petsc_KSPSolve(AppCtx *obj)
269d71ae5a4SJacob Faibussowitsch {
270c4762a1bSJed Brown   KSP ksp;
271c4762a1bSJed Brown   PC  pc;
272c4762a1bSJed Brown 
2733ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
274c4762a1bSJed Brown   /*create the ksp context and set the operators,that is, associate the system matrix with it*/
2759566063dSJacob Faibussowitsch   PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
2769566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(ksp, obj->Amat, obj->Amat));
277c4762a1bSJed Brown 
278c4762a1bSJed Brown   /*get the preconditioner context, set its type and the tolerances*/
2799566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(ksp, &pc));
2809566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCLU));
28109cb0f53SBarry Smith   PetscCall(KSPSetTolerances(ksp, 1.e-7, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT));
282c4762a1bSJed Brown 
283c4762a1bSJed Brown   /*get the command line options if there are any and set them*/
2849566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(ksp));
285c4762a1bSJed Brown 
286c4762a1bSJed Brown   /*get the linear system (ksp) solve*/
2879566063dSJacob Faibussowitsch   PetscCall(KSPSolve(ksp, obj->ksp_rhs, obj->ksp_sol));
288c4762a1bSJed Brown 
2899566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&ksp));
2903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
291c4762a1bSJed Brown }
292c4762a1bSJed Brown 
293c4762a1bSJed Brown /***********************************************************************
2940e3d61c9SBarry Smith   Function to return value of basis function or derivative of basis function.
295c4762a1bSJed Brown  ***********************************************************************
2960e3d61c9SBarry Smith 
2970e3d61c9SBarry Smith         Arguments:
2980e3d61c9SBarry Smith           x       = array of xpoints or nodal values
2990e3d61c9SBarry Smith           xx      = point at which the basis function is to be
3000e3d61c9SBarry Smith                       evaluated.
3010e3d61c9SBarry Smith           il      = interval containing xx.
3020e3d61c9SBarry Smith           iq      = indicates which of the two basis functions in
3030e3d61c9SBarry Smith                       interval intrvl should be used
3040e3d61c9SBarry Smith           nll     = array containing the endpoints of each interval.
3050e3d61c9SBarry Smith           id      = If id ~= 2, the value of the basis function
3060e3d61c9SBarry Smith                       is calculated; if id = 2, the value of the
3070e3d61c9SBarry Smith                       derivative of the basis function is returned.
308c4762a1bSJed Brown  ***********************************************************************/
309c4762a1bSJed Brown 
bspl(PetscScalar * x,PetscScalar xx,PetscInt il,PetscInt iq,PetscInt nll[][2],PetscInt id)310d71ae5a4SJacob Faibussowitsch PetscScalar bspl(PetscScalar *x, PetscScalar xx, PetscInt il, PetscInt iq, PetscInt nll[][2], PetscInt id)
311d71ae5a4SJacob Faibussowitsch {
312c4762a1bSJed Brown   PetscScalar x1, x2, bfcn;
313c4762a1bSJed Brown   PetscInt    i1, i2, iq1, iq2;
314c4762a1bSJed Brown 
3150e3d61c9SBarry Smith   /* Determine which basis function in interval intrvl is to be used in */
316c4762a1bSJed Brown   iq1 = iq;
317c4762a1bSJed Brown   if (iq1 == 0) iq2 = 1;
318c4762a1bSJed Brown   else iq2 = 0;
319c4762a1bSJed Brown 
3200e3d61c9SBarry Smith   /*    Determine endpoint of the interval intrvl   */
321c4762a1bSJed Brown   i1 = nll[il][iq1];
322c4762a1bSJed Brown   i2 = nll[il][iq2];
323c4762a1bSJed Brown 
3240e3d61c9SBarry Smith   /*   Determine nodal values at the endpoints of the interval intrvl   */
325c4762a1bSJed Brown   x1 = x[i1];
326c4762a1bSJed Brown   x2 = x[i2];
327303a5415SBarry Smith 
3280e3d61c9SBarry Smith   /*   Evaluate basis function   */
329c4762a1bSJed Brown   if (id == 2) bfcn = (1.0) / (x1 - x2);
330c4762a1bSJed Brown   else bfcn = (xx - x2) / (x1 - x2);
331c4762a1bSJed Brown   return bfcn;
332c4762a1bSJed Brown }
333c4762a1bSJed Brown 
334c4762a1bSJed Brown /*---------------------------------------------------------
335c4762a1bSJed Brown   Function called by rhs function to get B and g
336c4762a1bSJed Brown ---------------------------------------------------------*/
femBg(PetscScalar btri[][3],PetscScalar * f,PetscInt nz,PetscScalar * z,PetscReal t)337d71ae5a4SJacob Faibussowitsch PetscErrorCode femBg(PetscScalar btri[][3], PetscScalar *f, PetscInt nz, PetscScalar *z, PetscReal t)
338d71ae5a4SJacob Faibussowitsch {
339c4762a1bSJed Brown   PetscInt    i, j, jj, il, ip, ipp, ipq, iq, iquad, iqq;
340c4762a1bSJed Brown   PetscInt    nli[num_z][2], indx[num_z];
341c4762a1bSJed Brown   PetscScalar dd, dl, zip, zipq, zz, b_z, bb_z, bij;
342c4762a1bSJed Brown   PetscScalar zquad[num_z][3], dlen[num_z], qdwt[3];
343c4762a1bSJed Brown 
3443ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
345c4762a1bSJed Brown   /*  initializing everything - btri and f are initialized in rhs.c  */
346c4762a1bSJed Brown   for (i = 0; i < nz; i++) {
347c4762a1bSJed Brown     nli[i][0]   = 0;
348c4762a1bSJed Brown     nli[i][1]   = 0;
349c4762a1bSJed Brown     indx[i]     = 0;
350c4762a1bSJed Brown     zquad[i][0] = 0.0;
351c4762a1bSJed Brown     zquad[i][1] = 0.0;
352c4762a1bSJed Brown     zquad[i][2] = 0.0;
353c4762a1bSJed Brown     dlen[i]     = 0.0;
354c4762a1bSJed Brown   } /*end for (i)*/
355c4762a1bSJed Brown 
356c4762a1bSJed Brown   /*  quadrature weights  */
357c4762a1bSJed Brown   qdwt[0] = 1.0 / 6.0;
358c4762a1bSJed Brown   qdwt[1] = 4.0 / 6.0;
359c4762a1bSJed Brown   qdwt[2] = 1.0 / 6.0;
360c4762a1bSJed Brown 
361c4762a1bSJed Brown   /* 1st and last nodes have Dirichlet boundary condition -
362c4762a1bSJed Brown      set indices there to -1 */
363c4762a1bSJed Brown 
364c4762a1bSJed Brown   for (i = 0; i < nz - 1; i++) indx[i] = i - 1;
365c4762a1bSJed Brown   indx[nz - 1] = -1;
366c4762a1bSJed Brown 
367c4762a1bSJed Brown   ipq = 0;
368c4762a1bSJed Brown   for (il = 0; il < nz - 1; il++) {
369c4762a1bSJed Brown     ip           = ipq;
370c4762a1bSJed Brown     ipq          = ip + 1;
371c4762a1bSJed Brown     zip          = z[ip];
372c4762a1bSJed Brown     zipq         = z[ipq];
373c4762a1bSJed Brown     dl           = zipq - zip;
374c4762a1bSJed Brown     zquad[il][0] = zip;
375c4762a1bSJed Brown     zquad[il][1] = (0.5) * (zip + zipq);
376c4762a1bSJed Brown     zquad[il][2] = zipq;
377c4762a1bSJed Brown     dlen[il]     = PetscAbsScalar(dl);
378c4762a1bSJed Brown     nli[il][0]   = ip;
379c4762a1bSJed Brown     nli[il][1]   = ipq;
380c4762a1bSJed Brown   }
381c4762a1bSJed Brown 
382c4762a1bSJed Brown   for (il = 0; il < nz - 1; il++) {
383c4762a1bSJed Brown     for (iquad = 0; iquad < 3; iquad++) {
384c4762a1bSJed Brown       dd = (dlen[il]) * (qdwt[iquad]);
385c4762a1bSJed Brown       zz = zquad[il][iquad];
386c4762a1bSJed Brown 
387c4762a1bSJed Brown       for (iq = 0; iq < 2; iq++) {
388c4762a1bSJed Brown         ip  = nli[il][iq];
389c4762a1bSJed Brown         b_z = bspl(z, zz, il, iq, nli, 2);
390c4762a1bSJed Brown         i   = indx[ip];
391c4762a1bSJed Brown 
392c4762a1bSJed Brown         if (i > -1) {
393c4762a1bSJed Brown           for (iqq = 0; iqq < 2; iqq++) {
394c4762a1bSJed Brown             ipp  = nli[il][iqq];
395c4762a1bSJed Brown             bb_z = bspl(z, zz, il, iqq, nli, 2);
396c4762a1bSJed Brown             j    = indx[ipp];
397c4762a1bSJed Brown             bij  = -b_z * bb_z;
398c4762a1bSJed Brown 
399c4762a1bSJed Brown             if (j > -1) {
400c4762a1bSJed Brown               jj = 1 + j - i;
401c4762a1bSJed Brown               btri[i][jj] += bij * dd;
402c4762a1bSJed Brown             } else {
403c4762a1bSJed Brown               f[i] += bij * dd * exact(z[ipp], t);
404c4762a1bSJed Brown               /* f[i] += 0.0; */
405c4762a1bSJed Brown               /* if (il==0 && j==-1) { */
406c4762a1bSJed Brown               /* f[i] += bij*dd*exact(zz,t); */
407c4762a1bSJed Brown               /* }*/ /*end if*/
408c4762a1bSJed Brown             } /*end else*/
409c4762a1bSJed Brown           } /*end for (iqq)*/
410c4762a1bSJed Brown         } /*end if (i>0)*/
411c4762a1bSJed Brown       } /*end for (iq)*/
412c4762a1bSJed Brown     } /*end for (iquad)*/
413c4762a1bSJed Brown   } /*end for (il)*/
4143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
415c4762a1bSJed Brown }
416c4762a1bSJed Brown 
femA(AppCtx * obj,PetscInt nz,PetscScalar * z)417d71ae5a4SJacob Faibussowitsch PetscErrorCode femA(AppCtx *obj, PetscInt nz, PetscScalar *z)
418d71ae5a4SJacob Faibussowitsch {
419c4762a1bSJed Brown   PetscInt    i, j, il, ip, ipp, ipq, iq, iquad, iqq;
420c4762a1bSJed Brown   PetscInt    nli[num_z][2], indx[num_z];
421c4762a1bSJed Brown   PetscScalar dd, dl, zip, zipq, zz, bb, bbb, aij;
422c4762a1bSJed Brown   PetscScalar rquad[num_z][3], dlen[num_z], qdwt[3], add_term;
423c4762a1bSJed Brown 
4243ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
425c4762a1bSJed Brown   /*  initializing everything  */
426c4762a1bSJed Brown   for (i = 0; i < nz; i++) {
427c4762a1bSJed Brown     nli[i][0]   = 0;
428c4762a1bSJed Brown     nli[i][1]   = 0;
429c4762a1bSJed Brown     indx[i]     = 0;
430c4762a1bSJed Brown     rquad[i][0] = 0.0;
431c4762a1bSJed Brown     rquad[i][1] = 0.0;
432c4762a1bSJed Brown     rquad[i][2] = 0.0;
433c4762a1bSJed Brown     dlen[i]     = 0.0;
434c4762a1bSJed Brown   } /*end for (i)*/
435c4762a1bSJed Brown 
436c4762a1bSJed Brown   /*  quadrature weights  */
437c4762a1bSJed Brown   qdwt[0] = 1.0 / 6.0;
438c4762a1bSJed Brown   qdwt[1] = 4.0 / 6.0;
439c4762a1bSJed Brown   qdwt[2] = 1.0 / 6.0;
440c4762a1bSJed Brown 
441c4762a1bSJed Brown   /* 1st and last nodes have Dirichlet boundary condition -
442c4762a1bSJed Brown      set indices there to -1 */
443c4762a1bSJed Brown 
444c4762a1bSJed Brown   for (i = 0; i < nz - 1; i++) indx[i] = i - 1;
445c4762a1bSJed Brown   indx[nz - 1] = -1;
446c4762a1bSJed Brown 
447c4762a1bSJed Brown   ipq = 0;
448c4762a1bSJed Brown 
449c4762a1bSJed Brown   for (il = 0; il < nz - 1; il++) {
450c4762a1bSJed Brown     ip           = ipq;
451c4762a1bSJed Brown     ipq          = ip + 1;
452c4762a1bSJed Brown     zip          = z[ip];
453c4762a1bSJed Brown     zipq         = z[ipq];
454c4762a1bSJed Brown     dl           = zipq - zip;
455c4762a1bSJed Brown     rquad[il][0] = zip;
456c4762a1bSJed Brown     rquad[il][1] = (0.5) * (zip + zipq);
457c4762a1bSJed Brown     rquad[il][2] = zipq;
458c4762a1bSJed Brown     dlen[il]     = PetscAbsScalar(dl);
459c4762a1bSJed Brown     nli[il][0]   = ip;
460c4762a1bSJed Brown     nli[il][1]   = ipq;
461c4762a1bSJed Brown   } /*end for (il)*/
462c4762a1bSJed Brown 
463c4762a1bSJed Brown   for (il = 0; il < nz - 1; il++) {
464c4762a1bSJed Brown     for (iquad = 0; iquad < 3; iquad++) {
465c4762a1bSJed Brown       dd = (dlen[il]) * (qdwt[iquad]);
466c4762a1bSJed Brown       zz = rquad[il][iquad];
467c4762a1bSJed Brown 
468c4762a1bSJed Brown       for (iq = 0; iq < 2; iq++) {
469c4762a1bSJed Brown         ip = nli[il][iq];
470c4762a1bSJed Brown         bb = bspl(z, zz, il, iq, nli, 1);
471c4762a1bSJed Brown         i  = indx[ip];
472c4762a1bSJed Brown         if (i > -1) {
473c4762a1bSJed Brown           for (iqq = 0; iqq < 2; iqq++) {
474c4762a1bSJed Brown             ipp = nli[il][iqq];
475c4762a1bSJed Brown             bbb = bspl(z, zz, il, iqq, nli, 1);
476c4762a1bSJed Brown             j   = indx[ipp];
477c4762a1bSJed Brown             aij = bb * bbb;
478c4762a1bSJed Brown             if (j > -1) {
479c4762a1bSJed Brown               add_term = aij * dd;
4809566063dSJacob Faibussowitsch               PetscCall(MatSetValue(obj->Amat, i, j, add_term, ADD_VALUES));
481c4762a1bSJed Brown             } /*endif*/
482c4762a1bSJed Brown           } /*end for (iqq)*/
483c4762a1bSJed Brown         } /*end if (i>0)*/
484c4762a1bSJed Brown       } /*end for (iq)*/
485c4762a1bSJed Brown     } /*end for (iquad)*/
486c4762a1bSJed Brown   } /*end for (il)*/
4879566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(obj->Amat, MAT_FINAL_ASSEMBLY));
4889566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(obj->Amat, MAT_FINAL_ASSEMBLY));
4893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
490c4762a1bSJed Brown }
491c4762a1bSJed Brown 
492c4762a1bSJed Brown /*---------------------------------------------------------
493c4762a1bSJed Brown         Function to fill the rhs vector with
494c4762a1bSJed Brown         By + g values ****
495c4762a1bSJed Brown ---------------------------------------------------------*/
rhs(AppCtx * obj,PetscScalar * y,PetscInt nz,PetscScalar * z,PetscReal t)496d71ae5a4SJacob Faibussowitsch PetscErrorCode rhs(AppCtx *obj, PetscScalar *y, PetscInt nz, PetscScalar *z, PetscReal t)
497d71ae5a4SJacob Faibussowitsch {
498c4762a1bSJed Brown   PetscInt    i, j, js, je, jj;
499c4762a1bSJed Brown   PetscScalar val, g[num_z], btri[num_z][3], add_term;
500c4762a1bSJed Brown 
5013ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
502c4762a1bSJed Brown   for (i = 0; i < nz - 2; i++) {
503c4762a1bSJed Brown     for (j = 0; j <= 2; j++) btri[i][j] = 0.0;
504c4762a1bSJed Brown     g[i] = 0.0;
505c4762a1bSJed Brown   }
506c4762a1bSJed Brown 
507c4762a1bSJed Brown   /*  call femBg to set the tri-diagonal b matrix and vector g  */
5083ba16761SJacob Faibussowitsch   PetscCall(femBg(btri, g, nz, z, t));
509c4762a1bSJed Brown 
510dd8e379bSPierre Jolivet   /*  setting the entries of the right-hand side vector  */
511c4762a1bSJed Brown   for (i = 0; i < nz - 2; i++) {
512c4762a1bSJed Brown     val = 0.0;
513c4762a1bSJed Brown     js  = 0;
514c4762a1bSJed Brown     if (i == 0) js = 1;
515c4762a1bSJed Brown     je = 2;
516c4762a1bSJed Brown     if (i == nz - 2) je = 1;
517c4762a1bSJed Brown 
518c4762a1bSJed Brown     for (jj = js; jj <= je; jj++) {
519c4762a1bSJed Brown       j = i + jj - 1;
520c4762a1bSJed Brown       val += (btri[i][jj]) * (y[j]);
521c4762a1bSJed Brown     }
522c4762a1bSJed Brown     add_term = val + g[i];
523300f1712SStefano Zampini     PetscCall(VecSetValue(obj->ksp_rhs, i, add_term, INSERT_VALUES));
524c4762a1bSJed Brown   }
5259566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(obj->ksp_rhs));
5269566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(obj->ksp_rhs));
5273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
528c4762a1bSJed Brown }
529c4762a1bSJed Brown 
530c4762a1bSJed Brown /*%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
531dd8e379bSPierre Jolivet %%   Function to form the right-hand side of the time-stepping problem.                       %%
532c4762a1bSJed Brown %% -------------------------------------------------------------------------------------------%%
533c4762a1bSJed Brown   if (useAlhs):
534c4762a1bSJed Brown     globalout = By+g
535c4762a1bSJed Brown   else if (!useAlhs):
536c4762a1bSJed Brown     globalout = f(y,t)=Ainv(By+g),
537c4762a1bSJed Brown       in which the ksp solver to transform the problem A*ydot=By+g
538c4762a1bSJed Brown       to the problem ydot=f(y,t)=inv(A)*(By+g)
539c4762a1bSJed Brown %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%*/
540c4762a1bSJed Brown 
RHSfunction(TS ts,PetscReal t,Vec globalin,Vec globalout,PetscCtx ctx)541*2a8381b2SBarry Smith PetscErrorCode RHSfunction(TS ts, PetscReal t, Vec globalin, Vec globalout, PetscCtx ctx)
542d71ae5a4SJacob Faibussowitsch {
543c4762a1bSJed Brown   AppCtx            *obj = (AppCtx *)ctx;
544c4762a1bSJed Brown   PetscScalar        soln[num_z];
545c4762a1bSJed Brown   const PetscScalar *soln_ptr;
546c4762a1bSJed Brown   PetscInt           i, nz = obj->nz;
547c4762a1bSJed Brown   PetscReal          time;
548c4762a1bSJed Brown 
5493ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
550c4762a1bSJed Brown   /* get the previous solution to compute updated system */
5519566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(globalin, &soln_ptr));
552c4762a1bSJed Brown   for (i = 0; i < num_z - 2; i++) soln[i] = soln_ptr[i];
5539566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(globalin, &soln_ptr));
554c4762a1bSJed Brown   soln[num_z - 1] = 0.0;
555c4762a1bSJed Brown   soln[num_z - 2] = 0.0;
556c4762a1bSJed Brown 
557c4762a1bSJed Brown   /* clear out the matrix and rhs for ksp to keep things straight */
5589566063dSJacob Faibussowitsch   PetscCall(VecSet(obj->ksp_rhs, (PetscScalar)0.0));
559c4762a1bSJed Brown 
560c4762a1bSJed Brown   time = t;
561c4762a1bSJed Brown   /* get the updated system */
5623ba16761SJacob Faibussowitsch   PetscCall(rhs(obj, soln, nz, obj->z, time)); /* setup of the By+g rhs */
563c4762a1bSJed Brown 
564c4762a1bSJed Brown   /* do a ksp solve to get the rhs for the ts problem */
565c4762a1bSJed Brown   if (obj->useAlhs) {
566c4762a1bSJed Brown     /* ksp_sol = ksp_rhs */
5679566063dSJacob Faibussowitsch     PetscCall(VecCopy(obj->ksp_rhs, globalout));
568c4762a1bSJed Brown   } else {
569c4762a1bSJed Brown     /* ksp_sol = inv(Amat)*ksp_rhs */
5709566063dSJacob Faibussowitsch     PetscCall(Petsc_KSPSolve(obj));
5719566063dSJacob Faibussowitsch     PetscCall(VecCopy(obj->ksp_sol, globalout));
572c4762a1bSJed Brown   }
5733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
574c4762a1bSJed Brown }
575c4762a1bSJed Brown 
576c4762a1bSJed Brown /*TEST
577c4762a1bSJed Brown 
578c4762a1bSJed Brown     build:
579c4762a1bSJed Brown       requires: !complex
580c4762a1bSJed Brown 
581c4762a1bSJed Brown     test:
582c4762a1bSJed Brown       suffix: euler
583c4762a1bSJed Brown       output_file: output/ex3.out
584c4762a1bSJed Brown 
585c4762a1bSJed Brown     test:
586c4762a1bSJed Brown       suffix: 2
587c4762a1bSJed Brown       args: -useAlhs
588c4762a1bSJed Brown       output_file: output/ex3.out
589c4762a1bSJed Brown       TODO: Broken because SNESComputeJacobianDefault is incompatible with TSComputeIJacobianConstant
590c4762a1bSJed Brown 
591c4762a1bSJed Brown TEST*/
592