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