xref: /petsc/src/ts/tests/ex3.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Solves 1D heat equation with FEM formulation.\n\
3*c4762a1bSJed Brown Input arguments are\n\
4*c4762a1bSJed Brown   -useAlhs: solve Alhs*U' =  (Arhs*U + g) \n\
5*c4762a1bSJed Brown             otherwise, solve U' = inv(Alhs)*(Arhs*U + g) \n\n";
6*c4762a1bSJed Brown 
7*c4762a1bSJed Brown /*--------------------------------------------------------------------------
8*c4762a1bSJed Brown   Solves 1D heat equation U_t = U_xx with FEM formulation:
9*c4762a1bSJed Brown                           Alhs*U' = rhs (= Arhs*U + g)
10*c4762a1bSJed Brown   We thank Chris Cox <clcox@clemson.edu> for contributing the original code
11*c4762a1bSJed Brown ----------------------------------------------------------------------------*/
12*c4762a1bSJed Brown 
13*c4762a1bSJed Brown #include <petscksp.h>
14*c4762a1bSJed Brown #include <petscts.h>
15*c4762a1bSJed Brown 
16*c4762a1bSJed Brown /* special variable - max size of all arrays  */
17*c4762a1bSJed Brown #define num_z 10
18*c4762a1bSJed Brown 
19*c4762a1bSJed Brown /*
20*c4762a1bSJed Brown    User-defined application context - contains data needed by the
21*c4762a1bSJed Brown    application-provided call-back routines.
22*c4762a1bSJed Brown */
23*c4762a1bSJed Brown typedef struct {
24*c4762a1bSJed Brown   Mat         Amat;               /* left hand side matrix */
25*c4762a1bSJed Brown   Vec         ksp_rhs,ksp_sol;    /* working vectors for formulating inv(Alhs)*(Arhs*U+g) */
26*c4762a1bSJed Brown   int         max_probsz;         /* max size of the problem */
27*c4762a1bSJed Brown   PetscBool   useAlhs;            /* flag (1 indicates solving Alhs*U' = Arhs*U+g */
28*c4762a1bSJed Brown   int         nz;                 /* total number of grid points */
29*c4762a1bSJed Brown   PetscInt    m;                  /* total number of interio grid points */
30*c4762a1bSJed Brown   Vec         solution;           /* global exact ts solution vector */
31*c4762a1bSJed Brown   PetscScalar *z;                 /* array of grid points */
32*c4762a1bSJed Brown   PetscBool   debug;              /* flag (1 indicates activation of debugging printouts) */
33*c4762a1bSJed Brown } AppCtx;
34*c4762a1bSJed Brown 
35*c4762a1bSJed Brown extern PetscScalar exact(PetscScalar,PetscReal);
36*c4762a1bSJed Brown extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*);
37*c4762a1bSJed Brown extern PetscErrorCode Petsc_KSPSolve(AppCtx*);
38*c4762a1bSJed Brown extern PetscScalar bspl(PetscScalar*,PetscScalar,PetscInt,PetscInt,PetscInt[][2],PetscInt);
39*c4762a1bSJed Brown extern PetscErrorCode femBg(PetscScalar[][3],PetscScalar*,PetscInt,PetscScalar*,PetscReal);
40*c4762a1bSJed Brown extern PetscErrorCode femA(AppCtx*,PetscInt,PetscScalar*);
41*c4762a1bSJed Brown extern PetscErrorCode rhs(AppCtx*,PetscScalar*, PetscInt, PetscScalar*,PetscReal);
42*c4762a1bSJed Brown extern PetscErrorCode RHSfunction(TS,PetscReal,Vec,Vec,void*);
43*c4762a1bSJed Brown 
44*c4762a1bSJed Brown int main(int argc,char **argv)
45*c4762a1bSJed Brown {
46*c4762a1bSJed Brown   PetscInt       i,m,nz,steps,max_steps,k,nphase=1;
47*c4762a1bSJed Brown   PetscScalar    zInitial,zFinal,val,*z;
48*c4762a1bSJed Brown   PetscReal      stepsz[4],T,ftime;
49*c4762a1bSJed Brown   PetscErrorCode ierr;
50*c4762a1bSJed Brown   TS             ts;
51*c4762a1bSJed Brown   SNES           snes;
52*c4762a1bSJed Brown   Mat            Jmat;
53*c4762a1bSJed Brown   AppCtx         appctx;   /* user-defined application context */
54*c4762a1bSJed Brown   Vec            init_sol; /* ts solution vector */
55*c4762a1bSJed Brown   PetscMPIInt    size;
56*c4762a1bSJed Brown 
57*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
58*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
59*c4762a1bSJed Brown   if (size != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"This is a uniprocessor example only");
60*c4762a1bSJed Brown 
61*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-nphase",&nphase,NULL);CHKERRQ(ierr);
62*c4762a1bSJed Brown   if (nphase > 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"nphase must be an integer between 1 and 3");
63*c4762a1bSJed Brown 
64*c4762a1bSJed Brown   /* initializations */
65*c4762a1bSJed Brown   zInitial  = 0.0;
66*c4762a1bSJed Brown   zFinal    = 1.0;
67*c4762a1bSJed Brown   T         = 0.014/nphase;
68*c4762a1bSJed Brown   nz        = num_z;
69*c4762a1bSJed Brown   m         = nz-2;
70*c4762a1bSJed Brown   appctx.nz = nz;
71*c4762a1bSJed Brown   max_steps = (PetscInt)10000;
72*c4762a1bSJed Brown 
73*c4762a1bSJed Brown   appctx.m          = m;
74*c4762a1bSJed Brown   appctx.max_probsz = nz;
75*c4762a1bSJed Brown   appctx.debug      = PETSC_FALSE;
76*c4762a1bSJed Brown   appctx.useAlhs    = PETSC_FALSE;
77*c4762a1bSJed Brown 
78*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-debug",&appctx.debug);CHKERRQ(ierr);
79*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-useAlhs",&appctx.useAlhs);CHKERRQ(ierr);
80*c4762a1bSJed Brown 
81*c4762a1bSJed Brown   /* create vector to hold ts solution */
82*c4762a1bSJed Brown   /*-----------------------------------*/
83*c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD, &init_sol);CHKERRQ(ierr);
84*c4762a1bSJed Brown   ierr = VecSetSizes(init_sol, PETSC_DECIDE, m);CHKERRQ(ierr);
85*c4762a1bSJed Brown   ierr = VecSetFromOptions(init_sol);CHKERRQ(ierr);
86*c4762a1bSJed Brown 
87*c4762a1bSJed Brown   /* create vector to hold true ts soln for comparison */
88*c4762a1bSJed Brown   ierr = VecDuplicate(init_sol, &appctx.solution);CHKERRQ(ierr);
89*c4762a1bSJed Brown 
90*c4762a1bSJed Brown   /* create LHS matrix Amat */
91*c4762a1bSJed Brown   /*------------------------*/
92*c4762a1bSJed Brown   ierr = MatCreateSeqAIJ(PETSC_COMM_WORLD, m, m, 3, NULL, &appctx.Amat);CHKERRQ(ierr);
93*c4762a1bSJed Brown   ierr = MatSetFromOptions(appctx.Amat);CHKERRQ(ierr);
94*c4762a1bSJed Brown   ierr = MatSetUp(appctx.Amat);CHKERRQ(ierr);
95*c4762a1bSJed Brown   /* set space grid points - interio points only! */
96*c4762a1bSJed Brown   ierr = PetscMalloc1(nz+1,&z);CHKERRQ(ierr);
97*c4762a1bSJed Brown   for (i=0; i<nz; i++) z[i]=(i)*((zFinal-zInitial)/(nz-1));
98*c4762a1bSJed Brown   appctx.z = z;
99*c4762a1bSJed Brown   femA(&appctx,nz,z);
100*c4762a1bSJed Brown 
101*c4762a1bSJed Brown   /* create the jacobian matrix */
102*c4762a1bSJed Brown   /*----------------------------*/
103*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD, &Jmat);CHKERRQ(ierr);
104*c4762a1bSJed Brown   ierr = MatSetSizes(Jmat,PETSC_DECIDE,PETSC_DECIDE,m,m);CHKERRQ(ierr);
105*c4762a1bSJed Brown   ierr = MatSetFromOptions(Jmat);CHKERRQ(ierr);
106*c4762a1bSJed Brown   ierr = MatSetUp(Jmat);CHKERRQ(ierr);
107*c4762a1bSJed Brown 
108*c4762a1bSJed Brown   /* create working vectors for formulating rhs=inv(Alhs)*(Arhs*U + g) */
109*c4762a1bSJed Brown   ierr = VecDuplicate(init_sol,&appctx.ksp_rhs);CHKERRQ(ierr);
110*c4762a1bSJed Brown   ierr = VecDuplicate(init_sol,&appctx.ksp_sol);CHKERRQ(ierr);
111*c4762a1bSJed Brown 
112*c4762a1bSJed Brown   /* set intial guess */
113*c4762a1bSJed Brown   /*------------------*/
114*c4762a1bSJed Brown   for (i=0; i<nz-2; i++) {
115*c4762a1bSJed Brown     val  = exact(z[i+1], 0.0);
116*c4762a1bSJed Brown     ierr = VecSetValue(init_sol,i,(PetscScalar)val,INSERT_VALUES);CHKERRQ(ierr);
117*c4762a1bSJed Brown   }
118*c4762a1bSJed Brown   ierr = VecAssemblyBegin(init_sol);CHKERRQ(ierr);
119*c4762a1bSJed Brown   ierr = VecAssemblyEnd(init_sol);CHKERRQ(ierr);
120*c4762a1bSJed Brown 
121*c4762a1bSJed Brown   /*create a time-stepping context and set the problem type */
122*c4762a1bSJed Brown   /*--------------------------------------------------------*/
123*c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD, &ts);CHKERRQ(ierr);
124*c4762a1bSJed Brown   ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr);
125*c4762a1bSJed Brown 
126*c4762a1bSJed Brown   /* set time-step method */
127*c4762a1bSJed Brown   ierr = TSSetType(ts,TSCN);CHKERRQ(ierr);
128*c4762a1bSJed Brown 
129*c4762a1bSJed Brown   /* Set optional user-defined monitoring routine */
130*c4762a1bSJed Brown   ierr = TSMonitorSet(ts,Monitor,&appctx,NULL);CHKERRQ(ierr);
131*c4762a1bSJed Brown   /* set the right hand side of U_t = RHSfunction(U,t) */
132*c4762a1bSJed Brown   ierr = TSSetRHSFunction(ts,NULL,(PetscErrorCode (*)(TS,PetscScalar,Vec,Vec,void*))RHSfunction,&appctx);CHKERRQ(ierr);
133*c4762a1bSJed Brown 
134*c4762a1bSJed Brown   if (appctx.useAlhs) {
135*c4762a1bSJed Brown     /* set the left hand side matrix of Amat*U_t = rhs(U,t) */
136*c4762a1bSJed Brown 
137*c4762a1bSJed Brown     /* Note: this approach is incompatible with the finite differenced Jacobian set below because we can't restore the
138*c4762a1bSJed Brown      * Alhs matrix without making a copy.  Either finite difference the entire thing or use analytic Jacobians in both
139*c4762a1bSJed Brown      * places.
140*c4762a1bSJed Brown      */
141*c4762a1bSJed Brown     ierr = TSSetIFunction(ts,NULL,TSComputeIFunctionLinear,&appctx);CHKERRQ(ierr);
142*c4762a1bSJed Brown     ierr = TSSetIJacobian(ts,appctx.Amat,appctx.Amat,TSComputeIJacobianConstant,&appctx);CHKERRQ(ierr);
143*c4762a1bSJed Brown   }
144*c4762a1bSJed Brown 
145*c4762a1bSJed Brown   /* use petsc to compute the jacobian by finite differences */
146*c4762a1bSJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
147*c4762a1bSJed Brown   ierr = SNESSetJacobian(snes,Jmat,Jmat,SNESComputeJacobianDefault,NULL);CHKERRQ(ierr);
148*c4762a1bSJed Brown 
149*c4762a1bSJed Brown   /* get the command line options if there are any and set them */
150*c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
151*c4762a1bSJed Brown 
152*c4762a1bSJed Brown #if defined(PETSC_HAVE_SUNDIALS)
153*c4762a1bSJed Brown   {
154*c4762a1bSJed Brown     TSType    type;
155*c4762a1bSJed Brown     PetscBool sundialstype=PETSC_FALSE;
156*c4762a1bSJed Brown     ierr = TSGetType(ts,&type);CHKERRQ(ierr);
157*c4762a1bSJed Brown     ierr = PetscObjectTypeCompare((PetscObject)ts,TSSUNDIALS,&sundialstype);CHKERRQ(ierr);
158*c4762a1bSJed Brown     if (sundialstype && appctx.useAlhs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use Alhs formulation for TSSUNDIALS type");
159*c4762a1bSJed Brown   }
160*c4762a1bSJed Brown #endif
161*c4762a1bSJed Brown   /* Sets the initial solution */
162*c4762a1bSJed Brown   ierr = TSSetSolution(ts,init_sol);CHKERRQ(ierr);
163*c4762a1bSJed Brown 
164*c4762a1bSJed Brown   stepsz[0] = 1.0/(2.0*(nz-1)*(nz-1)); /* (mesh_size)^2/2.0 */
165*c4762a1bSJed Brown   ftime     = 0.0;
166*c4762a1bSJed Brown   for (k=0; k<nphase; k++) {
167*c4762a1bSJed Brown     if (nphase > 1) {ierr = PetscPrintf(PETSC_COMM_WORLD,"Phase %D initial time %g, stepsz %g, duration: %g\n",k,(double)ftime,(double)stepsz[k],(double)((k+1)*T));CHKERRQ(ierr);}
168*c4762a1bSJed Brown     ierr = TSSetTime(ts,ftime);CHKERRQ(ierr);
169*c4762a1bSJed Brown     ierr = TSSetTimeStep(ts,stepsz[k]);CHKERRQ(ierr);
170*c4762a1bSJed Brown     ierr = TSSetMaxSteps(ts,max_steps);CHKERRQ(ierr);
171*c4762a1bSJed Brown     ierr = TSSetMaxTime(ts,(k+1)*T);CHKERRQ(ierr);
172*c4762a1bSJed Brown     ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
173*c4762a1bSJed Brown 
174*c4762a1bSJed Brown     /* loop over time steps */
175*c4762a1bSJed Brown     /*----------------------*/
176*c4762a1bSJed Brown     ierr = TSSolve(ts,init_sol);CHKERRQ(ierr);
177*c4762a1bSJed Brown     ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr);
178*c4762a1bSJed Brown     ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr);
179*c4762a1bSJed Brown     stepsz[k+1] = stepsz[k]*1.5; /* change step size for the next phase */
180*c4762a1bSJed Brown   }
181*c4762a1bSJed Brown 
182*c4762a1bSJed Brown   /* free space */
183*c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
184*c4762a1bSJed Brown   ierr = MatDestroy(&appctx.Amat);CHKERRQ(ierr);
185*c4762a1bSJed Brown   ierr = MatDestroy(&Jmat);CHKERRQ(ierr);
186*c4762a1bSJed Brown   ierr = VecDestroy(&appctx.ksp_rhs);CHKERRQ(ierr);
187*c4762a1bSJed Brown   ierr = VecDestroy(&appctx.ksp_sol);CHKERRQ(ierr);
188*c4762a1bSJed Brown   ierr = VecDestroy(&init_sol);CHKERRQ(ierr);
189*c4762a1bSJed Brown   ierr = VecDestroy(&appctx.solution);CHKERRQ(ierr);
190*c4762a1bSJed Brown   ierr = PetscFree(z);CHKERRQ(ierr);
191*c4762a1bSJed Brown 
192*c4762a1bSJed Brown   ierr = PetscFinalize();
193*c4762a1bSJed Brown   return ierr;
194*c4762a1bSJed Brown }
195*c4762a1bSJed Brown 
196*c4762a1bSJed Brown /*------------------------------------------------------------------------
197*c4762a1bSJed Brown   Set exact solution
198*c4762a1bSJed Brown   u(z,t) = sin(6*PI*z)*exp(-36.*PI*PI*t) + 3.*sin(2*PI*z)*exp(-4.*PI*PI*t)
199*c4762a1bSJed Brown --------------------------------------------------------------------------*/
200*c4762a1bSJed Brown PetscScalar exact(PetscScalar z,PetscReal t)
201*c4762a1bSJed Brown {
202*c4762a1bSJed Brown   PetscScalar val, ex1, ex2;
203*c4762a1bSJed Brown 
204*c4762a1bSJed Brown   ex1 = PetscExpReal(-36.*PETSC_PI*PETSC_PI*t);
205*c4762a1bSJed Brown   ex2 = PetscExpReal(-4.*PETSC_PI*PETSC_PI*t);
206*c4762a1bSJed Brown   val = PetscSinScalar(6*PETSC_PI*z)*ex1 + 3.*PetscSinScalar(2*PETSC_PI*z)*ex2;
207*c4762a1bSJed Brown   return val;
208*c4762a1bSJed Brown }
209*c4762a1bSJed Brown 
210*c4762a1bSJed Brown /*
211*c4762a1bSJed Brown    Monitor - User-provided routine to monitor the solution computed at
212*c4762a1bSJed Brown    each timestep.  This example plots the solution and computes the
213*c4762a1bSJed Brown    error in two different norms.
214*c4762a1bSJed Brown 
215*c4762a1bSJed Brown    Input Parameters:
216*c4762a1bSJed Brown    ts     - the timestep context
217*c4762a1bSJed Brown    step   - the count of the current step (with 0 meaning the
218*c4762a1bSJed Brown              initial condition)
219*c4762a1bSJed Brown    time   - the current time
220*c4762a1bSJed Brown    u      - the solution at this timestep
221*c4762a1bSJed Brown    ctx    - the user-provided context for this monitoring routine.
222*c4762a1bSJed Brown             In this case we use the application context which contains
223*c4762a1bSJed Brown             information about the problem size, workspace and the exact
224*c4762a1bSJed Brown             solution.
225*c4762a1bSJed Brown */
226*c4762a1bSJed Brown PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec u,void *ctx)
227*c4762a1bSJed Brown {
228*c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ctx;
229*c4762a1bSJed Brown   PetscErrorCode ierr;
230*c4762a1bSJed Brown   PetscInt       i,m=appctx->m;
231*c4762a1bSJed Brown   PetscReal      norm_2,norm_max,h=1.0/(m+1);
232*c4762a1bSJed Brown   PetscScalar    *u_exact;
233*c4762a1bSJed Brown 
234*c4762a1bSJed Brown   /* Compute the exact solution */
235*c4762a1bSJed Brown   ierr = VecGetArray(appctx->solution,&u_exact);CHKERRQ(ierr);
236*c4762a1bSJed Brown   for (i=0; i<m; i++) u_exact[i] = exact(appctx->z[i+1],time);
237*c4762a1bSJed Brown   ierr = VecRestoreArray(appctx->solution,&u_exact);CHKERRQ(ierr);
238*c4762a1bSJed Brown 
239*c4762a1bSJed Brown   /* Print debugging information if desired */
240*c4762a1bSJed Brown   if (appctx->debug) {
241*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF,"Computed solution vector at time %g\n",(double)time);CHKERRQ(ierr);
242*c4762a1bSJed Brown     ierr = VecView(u,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
243*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF,"Exact solution vector\n");CHKERRQ(ierr);
244*c4762a1bSJed Brown     ierr = VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
245*c4762a1bSJed Brown   }
246*c4762a1bSJed Brown 
247*c4762a1bSJed Brown   /* Compute the 2-norm and max-norm of the error */
248*c4762a1bSJed Brown   ierr = VecAXPY(appctx->solution,-1.0,u);CHKERRQ(ierr);
249*c4762a1bSJed Brown   ierr = VecNorm(appctx->solution,NORM_2,&norm_2);CHKERRQ(ierr);
250*c4762a1bSJed Brown 
251*c4762a1bSJed Brown   norm_2 = PetscSqrtReal(h)*norm_2;
252*c4762a1bSJed Brown   ierr   = VecNorm(appctx->solution,NORM_MAX,&norm_max);CHKERRQ(ierr);
253*c4762a1bSJed Brown 
254*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF,"Timestep %D: time = %g, 2-norm error = %6.4f, max norm error = %6.4f\n",step,(double)time,(double)norm_2,(double)norm_max);CHKERRQ(ierr);
255*c4762a1bSJed Brown 
256*c4762a1bSJed Brown   /*
257*c4762a1bSJed Brown      Print debugging information if desired
258*c4762a1bSJed Brown   */
259*c4762a1bSJed Brown   if (appctx->debug) {
260*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF,"Error vector\n");CHKERRQ(ierr);
261*c4762a1bSJed Brown     ierr = VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
262*c4762a1bSJed Brown   }
263*c4762a1bSJed Brown   return 0;
264*c4762a1bSJed Brown }
265*c4762a1bSJed Brown 
266*c4762a1bSJed Brown /*%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
267*c4762a1bSJed Brown %%      Function to solve a linear system using KSP                                           %%
268*c4762a1bSJed Brown %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%*/
269*c4762a1bSJed Brown 
270*c4762a1bSJed Brown PetscErrorCode Petsc_KSPSolve(AppCtx *obj)
271*c4762a1bSJed Brown {
272*c4762a1bSJed Brown   PetscErrorCode ierr;
273*c4762a1bSJed Brown   KSP            ksp;
274*c4762a1bSJed Brown   PC             pc;
275*c4762a1bSJed Brown 
276*c4762a1bSJed Brown   /*create the ksp context and set the operators,that is, associate the system matrix with it*/
277*c4762a1bSJed Brown   ierr = KSPCreate(PETSC_COMM_WORLD,&ksp);CHKERRQ(ierr);
278*c4762a1bSJed Brown   ierr = KSPSetOperators(ksp,obj->Amat,obj->Amat);CHKERRQ(ierr);
279*c4762a1bSJed Brown 
280*c4762a1bSJed Brown   /*get the preconditioner context, set its type and the tolerances*/
281*c4762a1bSJed Brown   ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
282*c4762a1bSJed Brown   ierr = PCSetType(pc,PCLU);CHKERRQ(ierr);
283*c4762a1bSJed Brown   ierr = KSPSetTolerances(ksp,1.e-7,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
284*c4762a1bSJed Brown 
285*c4762a1bSJed Brown   /*get the command line options if there are any and set them*/
286*c4762a1bSJed Brown   ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
287*c4762a1bSJed Brown 
288*c4762a1bSJed Brown   /*get the linear system (ksp) solve*/
289*c4762a1bSJed Brown   ierr = KSPSolve(ksp,obj->ksp_rhs,obj->ksp_sol);CHKERRQ(ierr);
290*c4762a1bSJed Brown 
291*c4762a1bSJed Brown   KSPDestroy(&ksp);
292*c4762a1bSJed Brown   return 0;
293*c4762a1bSJed Brown }
294*c4762a1bSJed Brown 
295*c4762a1bSJed Brown /***********************************************************************
296*c4762a1bSJed Brown  * Function to return value of basis function or derivative of basis   *
297*c4762a1bSJed Brown  *              function.                                              *
298*c4762a1bSJed Brown  ***********************************************************************
299*c4762a1bSJed Brown  *                                                                     *
300*c4762a1bSJed Brown  *       Arguments:                                                    *
301*c4762a1bSJed Brown  *         x       = array of xpoints or nodal values                  *
302*c4762a1bSJed Brown  *         xx      = point at which the basis function is to be        *
303*c4762a1bSJed Brown  *                     evaluated.                                      *
304*c4762a1bSJed Brown  *         il      = interval containing xx.                           *
305*c4762a1bSJed Brown  *         iq      = indicates which of the two basis functions in     *
306*c4762a1bSJed Brown  *                     interval intrvl should be used                  *
307*c4762a1bSJed Brown  *         nll     = array containing the endpoints of each interval.  *
308*c4762a1bSJed Brown  *         id      = If id ~= 2, the value of the basis function       *
309*c4762a1bSJed Brown  *                     is calculated; if id = 2, the value of the      *
310*c4762a1bSJed Brown  *                     derivative of the basis function is returned.   *
311*c4762a1bSJed Brown  ***********************************************************************/
312*c4762a1bSJed Brown 
313*c4762a1bSJed Brown PetscScalar bspl(PetscScalar *x, PetscScalar xx,PetscInt il,PetscInt iq,PetscInt nll[][2],PetscInt id)
314*c4762a1bSJed Brown {
315*c4762a1bSJed Brown   PetscScalar x1,x2,bfcn;
316*c4762a1bSJed Brown   PetscInt    i1,i2,iq1,iq2;
317*c4762a1bSJed Brown 
318*c4762a1bSJed Brown   /*** Determine which basis function in interval intrvl is to be used in ***/
319*c4762a1bSJed Brown   iq1 = iq;
320*c4762a1bSJed Brown   if (iq1==0) iq2 = 1;
321*c4762a1bSJed Brown   else iq2 = 0;
322*c4762a1bSJed Brown 
323*c4762a1bSJed Brown   /***  Determine endpoint of the interval intrvl ***/
324*c4762a1bSJed Brown   i1=nll[il][iq1];
325*c4762a1bSJed Brown   i2=nll[il][iq2];
326*c4762a1bSJed Brown 
327*c4762a1bSJed Brown   /*** Determine nodal values at the endpoints of the interval intrvl ***/
328*c4762a1bSJed Brown   x1=x[i1];
329*c4762a1bSJed Brown   x2=x[i2];
330*c4762a1bSJed Brown   /* printf("x1=%g\tx2=%g\txx=%g\n",x1,x2,xx); */
331*c4762a1bSJed Brown   /*** Evaluate basis function ***/
332*c4762a1bSJed Brown   if (id == 2) bfcn=(1.0)/(x1-x2);
333*c4762a1bSJed Brown   else bfcn=(xx-x2)/(x1-x2);
334*c4762a1bSJed Brown   /* printf("bfcn=%g\n",bfcn); */
335*c4762a1bSJed Brown   return bfcn;
336*c4762a1bSJed Brown }
337*c4762a1bSJed Brown 
338*c4762a1bSJed Brown /*---------------------------------------------------------
339*c4762a1bSJed Brown   Function called by rhs function to get B and g
340*c4762a1bSJed Brown ---------------------------------------------------------*/
341*c4762a1bSJed Brown PetscErrorCode femBg(PetscScalar btri[][3],PetscScalar *f,PetscInt nz,PetscScalar *z, PetscReal t)
342*c4762a1bSJed Brown {
343*c4762a1bSJed Brown   PetscInt    i,j,jj,il,ip,ipp,ipq,iq,iquad,iqq;
344*c4762a1bSJed Brown   PetscInt    nli[num_z][2],indx[num_z];
345*c4762a1bSJed Brown   PetscScalar dd,dl,zip,zipq,zz,b_z,bb_z,bij;
346*c4762a1bSJed Brown   PetscScalar zquad[num_z][3],dlen[num_z],qdwt[3];
347*c4762a1bSJed Brown 
348*c4762a1bSJed Brown   /*  initializing everything - btri and f are initialized in rhs.c  */
349*c4762a1bSJed Brown   for (i=0; i < nz; i++) {
350*c4762a1bSJed Brown     nli[i][0]   = 0;
351*c4762a1bSJed Brown     nli[i][1]   = 0;
352*c4762a1bSJed Brown     indx[i]     = 0;
353*c4762a1bSJed Brown     zquad[i][0] = 0.0;
354*c4762a1bSJed Brown     zquad[i][1] = 0.0;
355*c4762a1bSJed Brown     zquad[i][2] = 0.0;
356*c4762a1bSJed Brown     dlen[i]     = 0.0;
357*c4762a1bSJed Brown   } /*end for (i)*/
358*c4762a1bSJed Brown 
359*c4762a1bSJed Brown   /*  quadrature weights  */
360*c4762a1bSJed Brown   qdwt[0] = 1.0/6.0;
361*c4762a1bSJed Brown   qdwt[1] = 4.0/6.0;
362*c4762a1bSJed Brown   qdwt[2] = 1.0/6.0;
363*c4762a1bSJed Brown 
364*c4762a1bSJed Brown   /* 1st and last nodes have Dirichlet boundary condition -
365*c4762a1bSJed Brown      set indices there to -1 */
366*c4762a1bSJed Brown 
367*c4762a1bSJed Brown   for (i=0; i < nz-1; i++) indx[i] = i-1;
368*c4762a1bSJed Brown   indx[nz-1] = -1;
369*c4762a1bSJed Brown 
370*c4762a1bSJed Brown   ipq = 0;
371*c4762a1bSJed Brown   for (il=0; il < nz-1; il++) {
372*c4762a1bSJed Brown     ip           = ipq;
373*c4762a1bSJed Brown     ipq          = ip+1;
374*c4762a1bSJed Brown     zip          = z[ip];
375*c4762a1bSJed Brown     zipq         = z[ipq];
376*c4762a1bSJed Brown     dl           = zipq-zip;
377*c4762a1bSJed Brown     zquad[il][0] = zip;
378*c4762a1bSJed Brown     zquad[il][1] = (0.5)*(zip+zipq);
379*c4762a1bSJed Brown     zquad[il][2] = zipq;
380*c4762a1bSJed Brown     dlen[il]     = PetscAbsScalar(dl);
381*c4762a1bSJed Brown     nli[il][0]   = ip;
382*c4762a1bSJed Brown     nli[il][1]   = ipq;
383*c4762a1bSJed Brown   }
384*c4762a1bSJed Brown 
385*c4762a1bSJed Brown   for (il=0; il < nz-1; il++) {
386*c4762a1bSJed Brown     for (iquad=0; iquad < 3; iquad++) {
387*c4762a1bSJed Brown       dd = (dlen[il])*(qdwt[iquad]);
388*c4762a1bSJed Brown       zz = zquad[il][iquad];
389*c4762a1bSJed Brown 
390*c4762a1bSJed Brown       for (iq=0; iq < 2; iq++) {
391*c4762a1bSJed Brown         ip  = nli[il][iq];
392*c4762a1bSJed Brown         b_z = bspl(z,zz,il,iq,nli,2);
393*c4762a1bSJed Brown         i   = indx[ip];
394*c4762a1bSJed Brown 
395*c4762a1bSJed Brown         if (i > -1) {
396*c4762a1bSJed Brown           for (iqq=0; iqq < 2; iqq++) {
397*c4762a1bSJed Brown             ipp  = nli[il][iqq];
398*c4762a1bSJed Brown             bb_z = bspl(z,zz,il,iqq,nli,2);
399*c4762a1bSJed Brown             j    = indx[ipp];
400*c4762a1bSJed Brown             bij  = -b_z*bb_z;
401*c4762a1bSJed Brown 
402*c4762a1bSJed Brown             if (j > -1) {
403*c4762a1bSJed Brown               jj = 1+j-i;
404*c4762a1bSJed Brown               btri[i][jj] += bij*dd;
405*c4762a1bSJed Brown             } else {
406*c4762a1bSJed Brown               f[i] += bij*dd*exact(z[ipp], t);
407*c4762a1bSJed Brown               /* f[i] += 0.0; */
408*c4762a1bSJed Brown               /* if (il==0 && j==-1) { */
409*c4762a1bSJed Brown               /* f[i] += bij*dd*exact(zz,t); */
410*c4762a1bSJed Brown               /* }*/ /*end if*/
411*c4762a1bSJed Brown             } /*end else*/
412*c4762a1bSJed Brown           } /*end for (iqq)*/
413*c4762a1bSJed Brown         } /*end if (i>0)*/
414*c4762a1bSJed Brown       } /*end for (iq)*/
415*c4762a1bSJed Brown     } /*end for (iquad)*/
416*c4762a1bSJed Brown   } /*end for (il)*/
417*c4762a1bSJed Brown   return 0;
418*c4762a1bSJed Brown }
419*c4762a1bSJed Brown 
420*c4762a1bSJed Brown PetscErrorCode femA(AppCtx *obj,PetscInt nz,PetscScalar *z)
421*c4762a1bSJed Brown {
422*c4762a1bSJed Brown   PetscInt       i,j,il,ip,ipp,ipq,iq,iquad,iqq;
423*c4762a1bSJed Brown   PetscInt       nli[num_z][2],indx[num_z];
424*c4762a1bSJed Brown   PetscScalar    dd,dl,zip,zipq,zz,bb,bbb,aij;
425*c4762a1bSJed Brown   PetscScalar    rquad[num_z][3],dlen[num_z],qdwt[3],add_term;
426*c4762a1bSJed Brown   PetscErrorCode ierr;
427*c4762a1bSJed Brown 
428*c4762a1bSJed Brown   /*  initializing everything  */
429*c4762a1bSJed Brown 
430*c4762a1bSJed Brown   for (i=0; i < nz; i++) {
431*c4762a1bSJed Brown     nli[i][0]   = 0;
432*c4762a1bSJed Brown     nli[i][1]   = 0;
433*c4762a1bSJed Brown     indx[i]     = 0;
434*c4762a1bSJed Brown     rquad[i][0] = 0.0;
435*c4762a1bSJed Brown     rquad[i][1] = 0.0;
436*c4762a1bSJed Brown     rquad[i][2] = 0.0;
437*c4762a1bSJed Brown     dlen[i]     = 0.0;
438*c4762a1bSJed Brown   } /*end for (i)*/
439*c4762a1bSJed Brown 
440*c4762a1bSJed Brown   /*  quadrature weights  */
441*c4762a1bSJed Brown   qdwt[0] = 1.0/6.0;
442*c4762a1bSJed Brown   qdwt[1] = 4.0/6.0;
443*c4762a1bSJed Brown   qdwt[2] = 1.0/6.0;
444*c4762a1bSJed Brown 
445*c4762a1bSJed Brown   /* 1st and last nodes have Dirichlet boundary condition -
446*c4762a1bSJed Brown      set indices there to -1 */
447*c4762a1bSJed Brown 
448*c4762a1bSJed Brown   for (i=0; i < nz-1; i++) indx[i]=i-1;
449*c4762a1bSJed Brown   indx[nz-1]=-1;
450*c4762a1bSJed Brown 
451*c4762a1bSJed Brown   ipq = 0;
452*c4762a1bSJed Brown 
453*c4762a1bSJed Brown   for (il=0; il < nz-1; il++) {
454*c4762a1bSJed Brown     ip           = ipq;
455*c4762a1bSJed Brown     ipq          = ip+1;
456*c4762a1bSJed Brown     zip          = z[ip];
457*c4762a1bSJed Brown     zipq         = z[ipq];
458*c4762a1bSJed Brown     dl           = zipq-zip;
459*c4762a1bSJed Brown     rquad[il][0] = zip;
460*c4762a1bSJed Brown     rquad[il][1] = (0.5)*(zip+zipq);
461*c4762a1bSJed Brown     rquad[il][2] = zipq;
462*c4762a1bSJed Brown     dlen[il]     = PetscAbsScalar(dl);
463*c4762a1bSJed Brown     nli[il][0]   = ip;
464*c4762a1bSJed Brown     nli[il][1]   = ipq;
465*c4762a1bSJed Brown   } /*end for (il)*/
466*c4762a1bSJed Brown 
467*c4762a1bSJed Brown   for (il=0; il < nz-1; il++) {
468*c4762a1bSJed Brown     for (iquad=0; iquad < 3; iquad++) {
469*c4762a1bSJed Brown       dd = (dlen[il])*(qdwt[iquad]);
470*c4762a1bSJed Brown       zz = rquad[il][iquad];
471*c4762a1bSJed Brown 
472*c4762a1bSJed Brown       for (iq=0; iq < 2; iq++) {
473*c4762a1bSJed Brown         ip = nli[il][iq];
474*c4762a1bSJed Brown         bb = bspl(z,zz,il,iq,nli,1);
475*c4762a1bSJed Brown         i = indx[ip];
476*c4762a1bSJed Brown         if (i > -1) {
477*c4762a1bSJed Brown           for (iqq=0; iqq < 2; iqq++) {
478*c4762a1bSJed Brown             ipp = nli[il][iqq];
479*c4762a1bSJed Brown             bbb = bspl(z,zz,il,iqq,nli,1);
480*c4762a1bSJed Brown             j = indx[ipp];
481*c4762a1bSJed Brown             aij = bb*bbb;
482*c4762a1bSJed Brown             if (j > -1) {
483*c4762a1bSJed Brown               add_term = aij*dd;
484*c4762a1bSJed Brown               ierr = MatSetValue(obj->Amat,i,j,add_term,ADD_VALUES);CHKERRQ(ierr);
485*c4762a1bSJed Brown             }/*endif*/
486*c4762a1bSJed Brown           } /*end for (iqq)*/
487*c4762a1bSJed Brown         } /*end if (i>0)*/
488*c4762a1bSJed Brown       } /*end for (iq)*/
489*c4762a1bSJed Brown     } /*end for (iquad)*/
490*c4762a1bSJed Brown   } /*end for (il)*/
491*c4762a1bSJed Brown   ierr = MatAssemblyBegin(obj->Amat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
492*c4762a1bSJed Brown   ierr = MatAssemblyEnd(obj->Amat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
493*c4762a1bSJed Brown   return 0;
494*c4762a1bSJed Brown }
495*c4762a1bSJed Brown 
496*c4762a1bSJed Brown /*---------------------------------------------------------
497*c4762a1bSJed Brown         Function to fill the rhs vector with
498*c4762a1bSJed Brown         By + g values ****
499*c4762a1bSJed Brown ---------------------------------------------------------*/
500*c4762a1bSJed Brown PetscErrorCode rhs(AppCtx *obj,PetscScalar *y, PetscInt nz, PetscScalar *z, PetscReal t)
501*c4762a1bSJed Brown {
502*c4762a1bSJed Brown   PetscInt       i,j,js,je,jj;
503*c4762a1bSJed Brown   PetscScalar    val,g[num_z],btri[num_z][3],add_term;
504*c4762a1bSJed Brown   PetscErrorCode ierr;
505*c4762a1bSJed Brown 
506*c4762a1bSJed Brown   for (i=0; i < nz-2; i++) {
507*c4762a1bSJed Brown     for (j=0; j <= 2; j++) btri[i][j]=0.0;
508*c4762a1bSJed Brown     g[i] = 0.0;
509*c4762a1bSJed Brown   }
510*c4762a1bSJed Brown 
511*c4762a1bSJed Brown   /*  call femBg to set the tri-diagonal b matrix and vector g  */
512*c4762a1bSJed Brown   femBg(btri,g,nz,z,t);
513*c4762a1bSJed Brown 
514*c4762a1bSJed Brown   /*  setting the entries of the right hand side vector  */
515*c4762a1bSJed Brown   for (i=0; i < nz-2; i++) {
516*c4762a1bSJed Brown     val = 0.0;
517*c4762a1bSJed Brown     js  = 0;
518*c4762a1bSJed Brown     if (i == 0) js = 1;
519*c4762a1bSJed Brown     je = 2;
520*c4762a1bSJed Brown     if (i == nz-2) je = 1;
521*c4762a1bSJed Brown 
522*c4762a1bSJed Brown     for (jj=js; jj <= je; jj++) {
523*c4762a1bSJed Brown       j    = i+jj-1;
524*c4762a1bSJed Brown       val += (btri[i][jj])*(y[j]);
525*c4762a1bSJed Brown     }
526*c4762a1bSJed Brown     add_term = val + g[i];
527*c4762a1bSJed Brown     ierr = VecSetValue(obj->ksp_rhs,(PetscInt)i,(PetscScalar)add_term,INSERT_VALUES);CHKERRQ(ierr);
528*c4762a1bSJed Brown   }
529*c4762a1bSJed Brown   ierr = VecAssemblyBegin(obj->ksp_rhs);CHKERRQ(ierr);
530*c4762a1bSJed Brown   ierr = VecAssemblyEnd(obj->ksp_rhs);CHKERRQ(ierr);
531*c4762a1bSJed Brown 
532*c4762a1bSJed Brown   /*  return to main driver function  */
533*c4762a1bSJed Brown   return 0;
534*c4762a1bSJed Brown }
535*c4762a1bSJed Brown 
536*c4762a1bSJed Brown /*%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
537*c4762a1bSJed Brown %%   Function to form the right hand side of the time-stepping problem.                       %%
538*c4762a1bSJed Brown %% -------------------------------------------------------------------------------------------%%
539*c4762a1bSJed Brown   if (useAlhs):
540*c4762a1bSJed Brown     globalout = By+g
541*c4762a1bSJed Brown   else if (!useAlhs):
542*c4762a1bSJed Brown     globalout = f(y,t)=Ainv(By+g),
543*c4762a1bSJed Brown       in which the ksp solver to transform the problem A*ydot=By+g
544*c4762a1bSJed Brown       to the problem ydot=f(y,t)=inv(A)*(By+g)
545*c4762a1bSJed Brown %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%*/
546*c4762a1bSJed Brown 
547*c4762a1bSJed Brown PetscErrorCode RHSfunction(TS ts,PetscReal t,Vec globalin,Vec globalout,void *ctx)
548*c4762a1bSJed Brown {
549*c4762a1bSJed Brown   PetscErrorCode    ierr;
550*c4762a1bSJed Brown   AppCtx            *obj = (AppCtx*)ctx;
551*c4762a1bSJed Brown   PetscScalar       soln[num_z];
552*c4762a1bSJed Brown   const PetscScalar *soln_ptr;
553*c4762a1bSJed Brown   PetscInt          i,nz=obj->nz;
554*c4762a1bSJed Brown   PetscReal         time;
555*c4762a1bSJed Brown 
556*c4762a1bSJed Brown   /* get the previous solution to compute updated system */
557*c4762a1bSJed Brown   ierr = VecGetArrayRead(globalin,&soln_ptr);CHKERRQ(ierr);
558*c4762a1bSJed Brown   for (i=0; i < num_z-2; i++) soln[i] = soln_ptr[i];
559*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(globalin,&soln_ptr);CHKERRQ(ierr);
560*c4762a1bSJed Brown   soln[num_z-1] = 0.0;
561*c4762a1bSJed Brown   soln[num_z-2] = 0.0;
562*c4762a1bSJed Brown 
563*c4762a1bSJed Brown   /* clear out the matrix and rhs for ksp to keep things straight */
564*c4762a1bSJed Brown   ierr = VecSet(obj->ksp_rhs,(PetscScalar)0.0);CHKERRQ(ierr);
565*c4762a1bSJed Brown 
566*c4762a1bSJed Brown   time = t;
567*c4762a1bSJed Brown   /* get the updated system */
568*c4762a1bSJed Brown   rhs(obj,soln,nz,obj->z,time); /* setup of the By+g rhs */
569*c4762a1bSJed Brown 
570*c4762a1bSJed Brown   /* do a ksp solve to get the rhs for the ts problem */
571*c4762a1bSJed Brown   if (obj->useAlhs) {
572*c4762a1bSJed Brown     /* ksp_sol = ksp_rhs */
573*c4762a1bSJed Brown     ierr = VecCopy(obj->ksp_rhs,globalout);CHKERRQ(ierr);
574*c4762a1bSJed Brown   } else {
575*c4762a1bSJed Brown     /* ksp_sol = inv(Amat)*ksp_rhs */
576*c4762a1bSJed Brown     ierr = Petsc_KSPSolve(obj);CHKERRQ(ierr);
577*c4762a1bSJed Brown     ierr = VecCopy(obj->ksp_sol,globalout);CHKERRQ(ierr);
578*c4762a1bSJed Brown   }
579*c4762a1bSJed Brown   return 0;
580*c4762a1bSJed Brown }
581*c4762a1bSJed Brown 
582*c4762a1bSJed Brown /*TEST
583*c4762a1bSJed Brown 
584*c4762a1bSJed Brown     build:
585*c4762a1bSJed Brown       requires: !complex
586*c4762a1bSJed Brown 
587*c4762a1bSJed Brown     test:
588*c4762a1bSJed Brown       suffix: euler
589*c4762a1bSJed Brown       output_file: output/ex3.out
590*c4762a1bSJed Brown 
591*c4762a1bSJed Brown     test:
592*c4762a1bSJed Brown       suffix: 2
593*c4762a1bSJed Brown       args:   -useAlhs
594*c4762a1bSJed Brown       output_file: output/ex3.out
595*c4762a1bSJed Brown       TODO: Broken because SNESComputeJacobianDefault is incompatible with TSComputeIJacobianConstant
596*c4762a1bSJed Brown 
597*c4762a1bSJed Brown TEST*/
598*c4762a1bSJed Brown 
599