Actual source code: ex10.c
1: /*$Id: ex10.c,v 1.29 2001/09/11 16:34:10 bsmith Exp $*/
3: /*
4: Program usage: mpirun -np <procs> usg [-help] [all PETSc options]
5: */
7: #if !defined(PETSC_USE_COMPLEX)
9: static char help[] = "An Unstructured Grid Example.n
10: This example demonstrates how to solve a nonlinear system in paralleln
11: with SNES for an unstructured mesh. The mesh and partitioning informationn
12: is read in an application defined ordering,which is later transformedn
13: into another convenient ordering (called the local ordering). The localn
14: ordering, apart from being efficient on cpu cycles and memory, allowsn
15: the use of the SPMD model of parallel programming. After partitioningn
16: is done, scatters are created between local (sequential)and globaln
17: (distributed) vectors. Finally, we set up the nonlinear solver contextn
18: in the usual way as a structured grid (seen
19: petsc/src/snes/examples/tutorials/ex5.c).n
20: The command line options include:n
21: -vert <Nv>, where Nv is the global number of nodesn
22: -elem <Ne>, where Ne is the global number of elementsn
23: -nl_par <lambda>, where lambda is the multiplier for the non linear term (u*u) termn
24: -lin_par <alpha>, where alpha is the multiplier for the linear term (u) n";
26: /*T
27: Concepts: SNES^unstructured grid
28: Concepts: AO^application to PETSc ordering or vice versa;
29: Concepts: VecScatter^using vector scatter operations;
30: Processors: n
31: T*/
33: /* ------------------------------------------------------------------------
35: PDE Solved : L(u) + lambda*u*u + alpha*u = 0 where L(u) is the Laplacian.
37: The Laplacian is approximated in the following way: each edge is given a weight
38: of one meaning that the diagonal term will have the weight equal to the degree
39: of a node. The off diagonal terms will get a weight of -1.
41: -----------------------------------------------------------------------*/
43: /*
44: Include petscao.h so that we can use AO (Application Ordering) object's services.
45: Include "petscsnes.h" so that we can use SNES solvers. Note that this
46: file automatically includes:
47: petsc.h - base PETSc routines petscvec.h - vectors
48: petscsys.h - system routines petscmat.h - matrices
49: petscis.h - index sets petscksp.h - Krylov subspace methods
50: petscviewer.h - viewers petscpc.h - preconditioners
51: petscsles.h - linear solvers
52: */
53: #include "petscao.h"
54: #include "petscsnes.h"
57: #define MAX_ELEM 500 /* Maximum number of elements */
58: #define MAX_VERT 100 /* Maximum number of vertices */
59: #define MAX_VERT_ELEM 3 /* Vertices per element */
61: /*
62: Application-defined context for problem specific data
63: */
64: typedef struct {
65: int Nvglobal,Nvlocal; /* global and local number of vertices */
66: int Neglobal,Nelocal; /* global and local number of vertices */
67: int AdjM[MAX_VERT][50]; /* adjacency list of a vertex */
68: int itot[MAX_VERT]; /* total number of neighbors for a vertex */
69: int icv[MAX_ELEM][MAX_VERT_ELEM]; /* vertices belonging to an element */
70: int v2p[MAX_VERT]; /* processor number for a vertex */
71: int *locInd,*gloInd; /* local and global orderings for a node */
72: Vec localX,localF; /* local solution (u) and f(u) vectors */
73: PetscReal non_lin_param; /* nonlinear parameter for the PDE */
74: PetscReal lin_param; /* linear parameter for the PDE */
75: VecScatter scatter; /* scatter context for the local and
76: distributed vectors */
77: } AppCtx;
79: /*
80: User-defined routines
81: */
82: int FormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*),
83: FormFunction(SNES,Vec,Vec,void*),
84: FormInitialGuess(AppCtx*,Vec);
86: int main(int argc,char **argv)
87: {
88: SNES snes; /* SNES context */
89: SNESType type = SNESLS; /* default nonlinear solution method */
90: Vec x,r; /* solution, residual vectors */
91: Mat Jac; /* Jacobian matrix */
92: AppCtx user; /* user-defined application context */
93: AO ao; /* Application Ordering object */
94: IS isglobal,islocal; /* global and local index sets */
95: int rank,size; /* rank of a process, number of processors */
96: int rstart; /* starting index of PETSc ordering for a processor */
97: int nfails; /* number of unsuccessful Newton steps */
98: int bs = 1; /* block size for multicomponent systems */
99: int nvertices; /* number of local plus ghost nodes of a processor */
100: int *pordering; /* PETSc ordering */
101: int *vertices; /* list of all vertices (incl. ghost ones)
102: on a processor */
103: int *verticesmask,*svertices;
104: int *tmp;
105: int i,j,jstart,inode,nb,nbrs,Nvneighborstotal = 0;
106: int ierr,its,N;
107: PetscScalar *xx;
108: char str[256],form[256],part_name[256];
109: FILE *fptr,*fptr1;
110: ISLocalToGlobalMapping isl2g;
111: #if defined (UNUSED_VARIABLES)
112: PetscDraw draw; /* drawing context */
113: PetscScalar *ff,*gg;
114: PetscReal tiny = 1.0e-10,zero = 0.0,one = 1.0,big = 1.0e+10;
115: int *tmp1,*tmp2;
116: #endif
117: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
118: Initialize program
119: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
121: PetscInitialize(&argc,&argv,"options.inf",help);
122: MPI_Comm_rank(MPI_COMM_WORLD,&rank);
123: MPI_Comm_size(MPI_COMM_WORLD,&size);
125: /* The current input file options.inf is for 2 proc run only */
126: if (size != 2) SETERRQ(1,"This Example currently runs on 2 procs only.");
128: /*
129: Initialize problem parameters
130: */
131: user.Nvglobal = 16; /*Global # of vertices */
132: user.Neglobal = 18; /*Global # of elements */
133: PetscOptionsGetInt(PETSC_NULL,"-vert",&user.Nvglobal,PETSC_NULL);
134: PetscOptionsGetInt(PETSC_NULL,"-elem",&user.Neglobal,PETSC_NULL);
135: user.non_lin_param = 0.06;
136: PetscOptionsGetReal(PETSC_NULL,"-nl_par",&user.non_lin_param,PETSC_NULL);
137: user.lin_param = -1.0;
138: PetscOptionsGetReal(PETSC_NULL,"-lin_par",&user.lin_param,PETSC_NULL);
139: user.Nvlocal = 0;
140: user.Nelocal = 0;
142: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
143: Read the mesh and partitioning information
144: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
145:
146: /*
147: Read the mesh and partitioning information from 'adj.in'.
148: The file format is as follows.
149: For each line the first entry is the processor rank where the
150: current node belongs. The second entry is the number of
151: neighbors of a node. The rest of the line is the adjacency
152: list of a node. Currently this file is set up to work on two
153: processors.
155: This is not a very good example of reading input. In the future,
156: we will put an example that shows the style that should be
157: used in a real application, where partitioning will be done
158: dynamically by calling partitioning routines (at present, we have
159: a ready interface to ParMeTiS).
160: */
161: fptr = fopen("adj.in","r");
162: if (!fptr) {
163: SETERRQ(0,"Could not open adj.in")
164: }
165:
166: /*
167: Each processor writes to the file output.<rank> where rank is the
168: processor's rank.
169: */
170: sprintf(part_name,"output.%d",rank);
171: fptr1 = fopen(part_name,"w");
172: if (!fptr1) {
173: SETERRQ(0,"Could no open output file");
174: }
175: PetscMalloc(user.Nvglobal*sizeof(int),&user.gloInd);
176: fprintf(fptr1,"Rank is %dn",rank);
177: for (inode = 0; inode < user.Nvglobal; inode++) {
178: fgets(str,256,fptr);
179: sscanf(str,"%d",&user.v2p[inode]);
180: if (user.v2p[inode] == rank) {
181: fprintf(fptr1,"Node %d belongs to processor %dn",inode,user.v2p[inode]);
182: user.gloInd[user.Nvlocal] = inode;
183: sscanf(str,"%*d %d",&nbrs);
184: fprintf(fptr1,"Number of neighbors for the vertex %d is %dn",inode,nbrs);
185: user.itot[user.Nvlocal] = nbrs;
186: Nvneighborstotal += nbrs;
187: for (i = 0; i < user.itot[user.Nvlocal]; i++){
188: form[0]='0';
189: for (j=0; j < i+2; j++){
190: PetscStrcat(form,"%*d ");
191: }
192: PetscStrcat(form,"%d");
193: sscanf(str,form,&user.AdjM[user.Nvlocal][i]);
194: fprintf(fptr1,"%d ",user.AdjM[user.Nvlocal][i]);
195: }
196: fprintf(fptr1,"n");
197: user.Nvlocal++;
198: }
199: }
200: fprintf(fptr1,"Total # of Local Vertices is %d n",user.Nvlocal);
202: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
203: Create different orderings
204: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
206: /*
207: Create the local ordering list for vertices. First a list using the PETSc global
208: ordering is created. Then we use the AO object to get the PETSc-to-application and
209: application-to-PETSc mappings. Each vertex also gets a local index (stored in the
210: locInd array).
211: */
212: MPI_Scan(&user.Nvlocal,&rstart,1,MPI_INT,MPI_SUM,MPI_COMM_WORLD);
213: rstart -= user.Nvlocal;
214: PetscMalloc(user.Nvlocal*sizeof(int),&pordering);
216: for (i=0; i < user.Nvlocal; i++) {
217: pordering[i] = rstart + i;
218: }
220: /*
221: Create the AO object
222: */
223: AOCreateBasic(MPI_COMM_WORLD,user.Nvlocal,user.gloInd,pordering,&ao);
224: PetscFree(pordering);
225:
226: /*
227: Keep the global indices for later use
228: */
229: PetscMalloc(user.Nvlocal*sizeof(int),&user.locInd);
230: PetscMalloc(Nvneighborstotal*sizeof(int),&tmp);
231:
232: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
233: Demonstrate the use of AO functionality
234: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
236: fprintf(fptr1,"Before AOApplicationToPetsc, local indices are : n");
237: for (i=0; i < user.Nvlocal; i++) {
238: fprintf(fptr1," %d ",user.gloInd[i]);
239: user.locInd[i] = user.gloInd[i];
240: }
241: fprintf(fptr1,"n");
242: jstart = 0;
243: for (i=0; i < user.Nvlocal; i++) {
244: fprintf(fptr1,"Neghbors of local vertex %d are : ",user.gloInd[i]);
245: for (j=0; j < user.itot[i]; j++) {
246: fprintf(fptr1,"%d ",user.AdjM[i][j]);
247: tmp[j + jstart] = user.AdjM[i][j];
248: }
249: jstart += user.itot[i];
250: fprintf(fptr1,"n");
251: }
253: /*
254: Now map the vlocal and neighbor lists to the PETSc ordering
255: */
256: AOApplicationToPetsc(ao,user.Nvlocal,user.locInd);
257: AOApplicationToPetsc(ao,Nvneighborstotal,tmp);
258:
259: fprintf(fptr1,"After AOApplicationToPetsc, local indices are : n");
260: for (i=0; i < user.Nvlocal; i++) {
261: fprintf(fptr1," %d ",user.locInd[i]);
262: }
263: fprintf(fptr1,"n");
265: jstart = 0;
266: for (i=0; i < user.Nvlocal; i++) {
267: fprintf(fptr1,"Neghbors of local vertex %d are : ",user.locInd[i]);
268: for (j=0; j < user.itot[i]; j++) {
269: user.AdjM[i][j] = tmp[j+jstart];
270: fprintf(fptr1,"%d ",user.AdjM[i][j]);
271: }
272: jstart += user.itot[i];
273: fprintf(fptr1,"n");
274: }
276: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
277: Extract the ghost vertex information for each processor
278: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
279: /*
280: Next, we need to generate a list of vertices required for this processor
281: and a local numbering scheme for all vertices required on this processor.
282: vertices - integer array of all vertices needed on this processor in PETSc
283: global numbering; this list consists of first the "locally owned"
284: vertices followed by the ghost vertices.
285: verticesmask - integer array that for each global vertex lists its local
286: vertex number (in vertices) + 1. If the global vertex is not
287: represented on this processor, then the corresponding
288: entry in verticesmask is zero
289:
290: Note: vertices and verticesmask are both Nvglobal in length; this may
291: sound terribly non-scalable, but in fact is not so bad for a reasonable
292: number of processors. Importantly, it allows us to use NO SEARCHING
293: in setting up the data structures.
294: */
295: ierr = PetscMalloc(user.Nvglobal*sizeof(int),&vertices);
296: ierr = PetscMalloc(user.Nvglobal*sizeof(int),&verticesmask);
297: ierr = PetscMemzero(verticesmask,user.Nvglobal*sizeof(int));
298: nvertices = 0;
299:
300: /*
301: First load "owned vertices" into list
302: */
303: for (i=0; i < user.Nvlocal; i++) {
304: vertices[nvertices++] = user.locInd[i];
305: verticesmask[user.locInd[i]] = nvertices;
306: }
307:
308: /*
309: Now load ghost vertices into list
310: */
311: for (i=0; i < user.Nvlocal; i++) {
312: for (j=0; j < user.itot[i]; j++) {
313: nb = user.AdjM[i][j];
314: if (!verticesmask[nb]) {
315: vertices[nvertices++] = nb;
316: verticesmask[nb] = nvertices;
317: }
318: }
319: }
321: fprintf(fptr1,"n");
322: fprintf(fptr1,"The array vertices is :n");
323: for (i=0; i < nvertices; i++) {
324: fprintf(fptr1,"%d ",vertices[i]);
325: }
326: fprintf(fptr1,"n");
327:
328: /*
329: Map the vertices listed in the neighbors to the local numbering from
330: the global ordering that they contained initially.
331: */
332: fprintf(fptr1,"n");
333: fprintf(fptr1,"After mapping neighbors in the local contiguous orderingn");
334: for (i=0; i<user.Nvlocal; i++) {
335: fprintf(fptr1,"Neghbors of local vertex %d are :n",i);
336: for (j = 0; j < user.itot[i]; j++) {
337: nb = user.AdjM[i][j];
338: user.AdjM[i][j] = verticesmask[nb] - 1;
339: fprintf(fptr1,"%d ",user.AdjM[i][j]);
340: }
341: fprintf(fptr1,"n");
342: }
344: N = user.Nvglobal;
345:
346: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
347: Create vector and matrix data structures
348: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
350: /*
351: Create vector data structures
352: */
353: VecCreate(MPI_COMM_WORLD,&x);
354: VecSetSizes(x,user.Nvlocal,N);
355: VecSetFromOptions(x);
356: VecDuplicate(x,&r);
357: VecCreateSeq(MPI_COMM_SELF,bs*nvertices,&user.localX);
358: VecDuplicate(user.localX,&user.localF);
360: /*
361: Create the scatter between the global representation and the
362: local representation
363: */
364: ISCreateStride(MPI_COMM_SELF,bs*nvertices,0,1,&islocal);
365: PetscMalloc(nvertices*sizeof(int),&svertices);
366: for (i=0; i<nvertices; i++) svertices[i] = bs*vertices[i];
367: ISCreateBlock(MPI_COMM_SELF,bs,nvertices,svertices,&isglobal);
368: PetscFree(svertices);
369: VecScatterCreate(x,isglobal,user.localX,islocal,&user.scatter);
371: /*
372: Create matrix data structure; Just to keep the example simple, we have not done any
373: preallocation of memory for the matrix. In real application code with big matrices,
374: preallocation should always be done to expedite the matrix creation.
375: */
376: MatCreate(MPI_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,N,N,&Jac);
377: MatSetFromOptions(Jac);
379: /*
380: The following routine allows us to set the matrix values in local ordering
381: */
382: ISLocalToGlobalMappingCreate(MPI_COMM_SELF,bs*nvertices,vertices,&isl2g);
383: MatSetLocalToGlobalMapping(Jac,isl2g);
385: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
386: Create nonlinear solver context
387: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
389: SNESCreate(MPI_COMM_WORLD,&snes);
390: SNESSetType(snes,type);
392: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
393: Set routines for function and Jacobian evaluation
394: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
396: FormInitialGuess(&user,x);
397: SNESSetFunction(snes,r,FormFunction,(void *)&user);
398: SNESSetJacobian(snes,Jac,Jac,FormJacobian,(void *)&user);
400: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
401: Customize nonlinear solver; set runtime options
402: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
404: SNESSetFromOptions(snes);
406: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
407: Evaluate initial guess; then solve nonlinear system
408: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
410: /*
411: Note: The user should initialize the vector, x, with the initial guess
412: for the nonlinear solver prior to calling SNESSolve(). In particular,
413: to employ an initial guess of zero, the user should explicitly set
414: this vector to zero by calling VecSet().
415: */
416: FormInitialGuess(&user,x);
418: /*
419: Print the initial guess
420: */
421: VecGetArray(x,&xx);
422: for (inode = 0; inode < user.Nvlocal; inode++)
423: fprintf(fptr1,"Initial Solution at node %d is %f n",inode,xx[inode]);
424: VecRestoreArray(x,&xx);
426: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
427: Now solve the nonlinear system
428: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
430: SNESSolve(snes,x,&its);
431: SNESGetNumberUnsuccessfulSteps(snes,&nfails);
432:
433: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
434: Print the output : solution vector and other information
435: Each processor writes to the file output.<rank> where rank is the
436: processor's rank.
437: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
439: VecGetArray(x,&xx);
440: for (inode = 0; inode < user.Nvlocal; inode++)
441: fprintf(fptr1,"Solution at node %d is %f n",inode,xx[inode]);
442: VecRestoreArray(x,&xx);
443: fclose(fptr1);
444: PetscPrintf(MPI_COMM_WORLD,"number of Newton iterations = %d, ",its);
445: PetscPrintf(MPI_COMM_WORLD,"number of unsuccessful steps = %dn",nfails);
447: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
448: Free work space. All PETSc objects should be destroyed when they
449: are no longer needed.
450: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
452: VecDestroy(x);
453: VecDestroy(r);
454: VecDestroy(user.localX);
455: VecDestroy(user.localF);
456: MatDestroy(Jac); SNESDestroy(snes);
457: /*PetscDrawDestroy(draw);*/
458: PetscFinalize();
460: return 0;
461: }
462: /* -------------------- Form initial approximation ----------------- */
464: /*
465: FormInitialGuess - Forms initial approximation.
467: Input Parameters:
468: user - user-defined application context
469: X - vector
471: Output Parameter:
472: X - vector
473: */
474: int FormInitialGuess(AppCtx *user,Vec X)
475: {
476: int i,Nvlocal,ierr;
477: int *gloInd;
478: PetscScalar *x;
479: #if defined (UNUSED_VARIABLES)
480: PetscReal temp1,temp,hx,hy,hxdhy,hydhx,sc;
481: int Neglobal,Nvglobal,j,row;
482: PetscReal alpha,lambda;
484: Nvglobal = user->Nvglobal;
485: Neglobal = user->Neglobal;
486: lambda = user->non_lin_param;
487: alpha = user->lin_param;
488: #endif
490: Nvlocal = user->Nvlocal;
491: gloInd = user->gloInd;
493: /*
494: Get a pointer to vector data.
495: - For default PETSc vectors, VecGetArray() returns a pointer to
496: the data array. Otherwise, the routine is implementation dependent.
497: - You MUST call VecRestoreArray() when you no longer need access to
498: the array.
499: */
500: VecGetArray(X,&x);
502: /*
503: Compute initial guess over the locally owned part of the grid
504: */
505: for (i=0; i < Nvlocal; i++) {
506: x[i] = (PetscReal)gloInd[i];
507: }
509: /*
510: Restore vector
511: */
512: VecRestoreArray(X,&x);
513: return 0;
514: }
515: /* -------------------- Evaluate Function F(x) --------------------- */
516: /*
517: FormFunction - Evaluates nonlinear function, F(x).
519: Input Parameters:
520: . snes - the SNES context
521: . X - input vector
522: . ptr - optional user-defined context, as set by SNESSetFunction()
524: Output Parameter:
525: . F - function vector
526: */
527: int FormFunction(SNES snes,Vec X,Vec F,void *ptr)
528: {
529: AppCtx *user = (AppCtx*)ptr;
530: int ierr,i,j,Nvlocal;
531: PetscReal alpha,lambda;
532: PetscScalar *x,*f;
533: VecScatter scatter;
534: Vec localX = user->localX;
535: #if defined (UNUSED_VARIABLES)
536: PetscScalar ut,ub,ul,ur,u,*g,sc,uyy,uxx;
537: PetscReal hx,hy,hxdhy,hydhx;
538: PetscReal two = 2.0,one = 1.0;
539: int Nvglobal,Neglobal,row;
540: int *gloInd;
542: Nvglobal = user->Nvglobal;
543: Neglobal = user->Neglobal;
544: gloInd = user->gloInd;
545: #endif
547: Nvlocal = user->Nvlocal;
548: lambda = user->non_lin_param;
549: alpha = user->lin_param;
550: scatter = user->scatter;
552: /*
553: PDE : L(u) + lambda*u*u +alpha*u = 0 where L(u) is the approximate Laplacian as
554: described in the beginning of this code
555:
556: First scatter the distributed vector X into local vector localX (that includes
557: values for ghost nodes. If we wish,we can put some other work between
558: VecScatterBegin() and VecScatterEnd() to overlap the communication with
559: computation.
560: */
561: VecScatterBegin(X,localX,INSERT_VALUES,SCATTER_FORWARD,scatter);
562: VecScatterEnd(X,localX,INSERT_VALUES,SCATTER_FORWARD,scatter);
564: /*
565: Get pointers to vector data
566: */
567: VecGetArray(localX,&x);
568: VecGetArray(F,&f);
570: /*
571: Now compute the f(x). As mentioned earlier, the computed Laplacian is just an
572: approximate one chosen for illustrative purpose only. Another point to notice
573: is that this is a local (completly parallel) calculation. In practical application
574: codes, function calculation time is a dominat portion of the overall execution time.
575: */
576: for (i=0; i < Nvlocal; i++) {
577: f[i] = (user->itot[i] - alpha)*x[i] - lambda*x[i]*x[i];
578: for (j = 0; j < user->itot[i]; j++) {
579: f[i] -= x[user->AdjM[i][j]];
580: }
581: }
583: /*
584: Restore vectors
585: */
586: VecRestoreArray(localX,&x);
587: VecRestoreArray(F,&f);
588: /*VecView(F,PETSC_VIEWER_STDOUT_WORLD);*/
590: return 0;
591: }
593: /* -------------------- Evaluate Jacobian F'(x) -------------------- */
594: /*
595: FormJacobian - Evaluates Jacobian matrix.
597: Input Parameters:
598: . snes - the SNES context
599: . X - input vector
600: . ptr - optional user-defined context, as set by SNESSetJacobian()
602: Output Parameters:
603: . A - Jacobian matrix
604: . B - optionally different preconditioning matrix
605: . flag - flag indicating matrix structure
607: */
608: int FormJacobian(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
609: {
610: AppCtx *user = (AppCtx*)ptr;
611: Mat jac = *B;
612: int i,j,Nvlocal,col[50],ierr;
613: PetscScalar alpha,lambda,value[50];
614: Vec localX = user->localX;
615: VecScatter scatter;
616: PetscScalar *x;
617: #if defined (UNUSED_VARIABLES)
618: PetscScalar two = 2.0,one = 1.0;
619: int row,Nvglobal,Neglobal;
620: int *gloInd;
622: Nvglobal = user->Nvglobal;
623: Neglobal = user->Neglobal;
624: gloInd = user->gloInd;
625: #endif
626:
627: /*printf("Entering into FormJacobian n");*/
628: Nvlocal = user->Nvlocal;
629: lambda = user->non_lin_param;
630: alpha = user->lin_param;
631: scatter = user->scatter;
633: /*
634: PDE : L(u) + lambda*u*u +alpha*u = 0 where L(u) is the approximate Laplacian as
635: described in the beginning of this code
636:
637: First scatter the distributed vector X into local vector localX (that includes
638: values for ghost nodes. If we wish, we can put some other work between
639: VecScatterBegin() and VecScatterEnd() to overlap the communication with
640: computation.
641: */
642: VecScatterBegin(X,localX,INSERT_VALUES,SCATTER_FORWARD,scatter);
643: VecScatterEnd(X,localX,INSERT_VALUES,SCATTER_FORWARD,scatter);
644:
645: /*
646: Get pointer to vector data
647: */
648: VecGetArray(localX,&x);
650: for (i=0; i < Nvlocal; i++) {
651: col[0] = i;
652: value[0] = user->itot[i] - 2.0*lambda*x[i] - alpha;
653: for (j = 0; j < user->itot[i]; j++) {
654: col[j+1] = user->AdjM[i][j];
655: value[j+1] = -1.0;
656: }
658: /*
659: Set the matrix values in the local ordering. Note that in order to use this
660: feature we must call the routine MatSetLocalToGlobalMapping() after the
661: matrix has been created.
662: */
663: MatSetValuesLocal(jac,1,&i,1+user->itot[i],col,value,INSERT_VALUES);
664: }
666: /*
667: Assemble matrix, using the 2-step process:
668: MatAssemblyBegin(), MatAssemblyEnd().
669: Between these two calls, the pointer to vector data has been restored to
670: demonstrate the use of overlapping communicationn with computation.
671: */
672: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
673: VecRestoreArray(localX,&x);
674: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
676: /*
677: Set flag to indicate that the Jacobian matrix retains an identical
678: nonzero structure throughout all nonlinear iterations (although the
679: values of the entries change). Thus, we can save some work in setting
680: up the preconditioner (e.g., no need to redo symbolic factorization for
681: ILU/ICC preconditioners).
682: - If the nonzero structure of the matrix is different during
683: successive linear solves, then the flag DIFFERENT_NONZERO_PATTERN
684: must be used instead. If you are unsure whether the matrix
685: structure has changed or not, use the flag DIFFERENT_NONZERO_PATTERN.
686: - Caution: If you specify SAME_NONZERO_PATTERN, PETSc
687: believes your assertion and does not check the structure
688: of the matrix. If you erroneously claim that the structure
689: is the same when it actually is not, the new preconditioner
690: will not function correctly. Thus, use this optimization
691: feature with caution!
692: */
693: *flag = SAME_NONZERO_PATTERN;
695: /*
696: Tell the matrix we will never add a new nonzero location to the
697: matrix. If we do, it will generate an error.
698: */
699: MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR);
700: /* MatView(jac,PETSC_VIEWER_STDOUT_SELF); */
701: return 0;
702: }
703: #else
704: #include <stdio.h>
705: int main(int argc,char **args)
706: {
707: fprintf(stdout,"This example does not work for complex numbers.n");
708: return 0;
709: }
710: #endif