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