Actual source code: ex10.c
1: /*
2: Program usage: mpiexec -np <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: 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: 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: 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,*svertices;
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: PetscTruth 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(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) {
170: SETERRQ(0,"Could not open adj.in")
171: }
172:
173: /*
174: Each processor writes to the file output.<rank> where rank is the
175: processor's rank.
176: */
177: sprintf(part_name,"output.%d",rank);
178: fptr1 = fopen(part_name,"w");
179: if (!fptr1) {
180: SETERRQ(0,"Could no open output file");
181: }
182: PetscMalloc(user.Nvglobal*sizeof(PetscInt),&user.gloInd);
183: PetscFPrintf(PETSC_COMM_SELF,fptr1,"Rank is %D\n",rank);
184: for (inode = 0; inode < user.Nvglobal; inode++) {
185: fgets(str,256,fptr);
186: sscanf(str,"%d",&dtmp);user.v2p[inode] = dtmp;
187: if (user.v2p[inode] == rank) {
188: PetscFPrintf(PETSC_COMM_SELF,fptr1,"Node %D belongs to processor %D\n",inode,user.v2p[inode]);
189: user.gloInd[user.Nvlocal] = inode;
190: sscanf(str,"%*d %d",&dtmp);nbrs = dtmp;
191: PetscFPrintf(PETSC_COMM_SELF,fptr1,"Number of neighbors for the vertex %D is %D\n",inode,nbrs);
192: user.itot[user.Nvlocal] = nbrs;
193: Nvneighborstotal += nbrs;
194: for (i = 0; i < user.itot[user.Nvlocal]; i++){
195: form[0]='\0';
196: for (j=0; j < i+2; j++){
197: PetscStrcat(form,"%*d ");
198: }
199: PetscStrcat(form,"%d");
200: sscanf(str,form,&dtmp); user.AdjM[user.Nvlocal][i] = dtmp;
201: PetscFPrintf(PETSC_COMM_SELF,fptr1,"%D ",user.AdjM[user.Nvlocal][i]);
202: }
203: PetscFPrintf(PETSC_COMM_SELF,fptr1,"\n");
204: user.Nvlocal++;
205: }
206: }
207: PetscFPrintf(PETSC_COMM_SELF,fptr1,"Total # of Local Vertices is %D \n",user.Nvlocal);
209: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
210: Create different orderings
211: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
213: /*
214: Create the local ordering list for vertices. First a list using the PETSc global
215: ordering is created. Then we use the AO object to get the PETSc-to-application and
216: application-to-PETSc mappings. Each vertex also gets a local index (stored in the
217: locInd array).
218: */
219: MPI_Scan(&user.Nvlocal,&rstart,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD);
220: rstart -= user.Nvlocal;
221: PetscMalloc(user.Nvlocal*sizeof(PetscInt),&pordering);
223: for (i=0; i < user.Nvlocal; i++) {
224: pordering[i] = rstart + i;
225: }
227: /*
228: Create the AO object
229: */
230: AOCreateBasic(MPI_COMM_WORLD,user.Nvlocal,user.gloInd,pordering,&ao);
231: PetscFree(pordering);
232:
233: /*
234: Keep the global indices for later use
235: */
236: PetscMalloc(user.Nvlocal*sizeof(PetscInt),&user.locInd);
237: PetscMalloc(Nvneighborstotal*sizeof(PetscInt),&tmp);
238:
239: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
240: Demonstrate the use of AO functionality
241: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
243: PetscFPrintf(PETSC_COMM_SELF,fptr1,"Before AOApplicationToPetsc, local indices are : \n");
244: for (i=0; i < user.Nvlocal; i++) {
245: PetscFPrintf(PETSC_COMM_SELF,fptr1," %D ",user.gloInd[i]);
246: user.locInd[i] = user.gloInd[i];
247: }
248: PetscFPrintf(PETSC_COMM_SELF,fptr1,"\n");
249: jstart = 0;
250: for (i=0; i < user.Nvlocal; i++) {
251: PetscFPrintf(PETSC_COMM_SELF,fptr1,"Neghbors of local vertex %D are : ",user.gloInd[i]);
252: for (j=0; j < user.itot[i]; j++) {
253: PetscFPrintf(PETSC_COMM_SELF,fptr1,"%D ",user.AdjM[i][j]);
254: tmp[j + jstart] = user.AdjM[i][j];
255: }
256: jstart += user.itot[i];
257: PetscFPrintf(PETSC_COMM_SELF,fptr1,"\n");
258: }
260: /*
261: Now map the vlocal and neighbor lists to the PETSc ordering
262: */
263: AOApplicationToPetsc(ao,user.Nvlocal,user.locInd);
264: AOApplicationToPetsc(ao,Nvneighborstotal,tmp);
265: AODestroy(ao);
267: PetscFPrintf(PETSC_COMM_SELF,fptr1,"After AOApplicationToPetsc, local indices are : \n");
268: for (i=0; i < user.Nvlocal; i++) {
269: PetscFPrintf(PETSC_COMM_SELF,fptr1," %D ",user.locInd[i]);
270: }
271: PetscFPrintf(PETSC_COMM_SELF,fptr1,"\n");
273: jstart = 0;
274: for (i=0; i < user.Nvlocal; i++) {
275: PetscFPrintf(PETSC_COMM_SELF,fptr1,"Neghbors of local vertex %D are : ",user.locInd[i]);
276: for (j=0; j < user.itot[i]; j++) {
277: user.AdjM[i][j] = tmp[j+jstart];
278: PetscFPrintf(PETSC_COMM_SELF,fptr1,"%D ",user.AdjM[i][j]);
279: }
280: jstart += user.itot[i];
281: PetscFPrintf(PETSC_COMM_SELF,fptr1,"\n");
282: }
284: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
285: Extract the ghost vertex information for each processor
286: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
287: /*
288: Next, we need to generate a list of vertices required for this processor
289: and a local numbering scheme for all vertices required on this processor.
290: vertices - integer array of all vertices needed on this processor in PETSc
291: global numbering; this list consists of first the "locally owned"
292: vertices followed by the ghost vertices.
293: verticesmask - integer array that for each global vertex lists its local
294: vertex number (in vertices) + 1. If the global vertex is not
295: represented on this processor, then the corresponding
296: entry in verticesmask is zero
297:
298: Note: vertices and verticesmask are both Nvglobal in length; this may
299: sound terribly non-scalable, but in fact is not so bad for a reasonable
300: number of processors. Importantly, it allows us to use NO SEARCHING
301: in setting up the data structures.
302: */
303: PetscMalloc(user.Nvglobal*sizeof(PetscInt),&vertices);
304: PetscMalloc(user.Nvglobal*sizeof(PetscInt),&verticesmask);
305: PetscMemzero(verticesmask,user.Nvglobal*sizeof(PetscInt));
306: nvertices = 0;
307:
308: /*
309: First load "owned vertices" into list
310: */
311: for (i=0; i < user.Nvlocal; i++) {
312: vertices[nvertices++] = user.locInd[i];
313: verticesmask[user.locInd[i]] = nvertices;
314: }
315:
316: /*
317: Now load ghost vertices into list
318: */
319: for (i=0; i < user.Nvlocal; i++) {
320: for (j=0; j < user.itot[i]; j++) {
321: nb = user.AdjM[i][j];
322: if (!verticesmask[nb]) {
323: vertices[nvertices++] = nb;
324: verticesmask[nb] = nvertices;
325: }
326: }
327: }
329: PetscFPrintf(PETSC_COMM_SELF,fptr1,"\n");
330: PetscFPrintf(PETSC_COMM_SELF,fptr1,"The array vertices is :\n");
331: for (i=0; i < nvertices; i++) {
332: PetscFPrintf(PETSC_COMM_SELF,fptr1,"%D ",vertices[i]);
333: }
334: PetscFPrintf(PETSC_COMM_SELF,fptr1,"\n");
335:
336: /*
337: Map the vertices listed in the neighbors to the local numbering from
338: the global ordering that they contained initially.
339: */
340: PetscFPrintf(PETSC_COMM_SELF,fptr1,"\n");
341: PetscFPrintf(PETSC_COMM_SELF,fptr1,"After mapping neighbors in the local contiguous ordering\n");
342: for (i=0; i<user.Nvlocal; i++) {
343: PetscFPrintf(PETSC_COMM_SELF,fptr1,"Neghbors of local vertex %D are :\n",i);
344: for (j = 0; j < user.itot[i]; j++) {
345: nb = user.AdjM[i][j];
346: user.AdjM[i][j] = verticesmask[nb] - 1;
347: PetscFPrintf(PETSC_COMM_SELF,fptr1,"%D ",user.AdjM[i][j]);
348: }
349: PetscFPrintf(PETSC_COMM_SELF,fptr1,"\n");
350: }
352: N = user.Nvglobal;
353:
354: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
355: Create vector and matrix data structures
356: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
358: /*
359: Create vector data structures
360: */
361: VecCreate(MPI_COMM_WORLD,&x);
362: VecSetSizes(x,user.Nvlocal,N);
363: VecSetFromOptions(x);
364: VecDuplicate(x,&r);
365: VecCreateSeq(MPI_COMM_SELF,bs*nvertices,&user.localX);
366: VecDuplicate(user.localX,&user.localF);
368: /*
369: Create the scatter between the global representation and the
370: local representation
371: */
372: ISCreateStride(MPI_COMM_SELF,bs*nvertices,0,1,&islocal);
373: PetscMalloc(nvertices*sizeof(PetscInt),&svertices);
374: for (i=0; i<nvertices; i++) svertices[i] = bs*vertices[i];
375: ISCreateBlock(MPI_COMM_SELF,bs,nvertices,svertices,&isglobal);
376: PetscFree(svertices);
377: VecScatterCreate(x,isglobal,user.localX,islocal,&user.scatter);
378: ISDestroy(isglobal);
379: ISDestroy(islocal);
381: /*
382: Create matrix data structure; Just to keep the example simple, we have not done any
383: preallocation of memory for the matrix. In real application code with big matrices,
384: preallocation should always be done to expedite the matrix creation.
385: */
386: MatCreate(MPI_COMM_WORLD,&Jac);
387: MatSetSizes(Jac,PETSC_DECIDE,PETSC_DECIDE,N,N);
388: MatSetFromOptions(Jac);
390: /*
391: The following routine allows us to set the matrix values in local ordering
392: */
393: ISLocalToGlobalMappingCreate(MPI_COMM_SELF,bs*nvertices,vertices,&isl2g);
394: MatSetLocalToGlobalMapping(Jac,isl2g);
396: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
397: Create nonlinear solver context
398: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
400: SNESCreate(MPI_COMM_WORLD,&snes);
401: SNESSetType(snes,type);
402: FormInitialGuess(&user,x);
404: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
405: Set routines for function and Jacobian evaluation
406: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
407: SNESSetFunction(snes,r,FormFunction,(void*)&user);
409: PetscOptionsGetTruth(PETSC_NULL,"-fd_jacobian_coloring",&fd_jacobian_coloring,0);
410: if (!fd_jacobian_coloring){
411: SNESSetJacobian(snes,Jac,Jac,FormJacobian,(void*)&user);
412: } else { /* Use matfdcoloring */
413: ISColoring iscoloring;
414: MatStructure flag;
415: FormJacobian(snes,x,&Jac,&Jac,&flag,&user);
416: MatGetColoring(Jac,MATCOLORING_SL,&iscoloring);
417: MatFDColoringCreate(Jac,iscoloring,&matfdcoloring);
418: MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))FormFunction,&user);
419: MatFDColoringSetFromOptions(matfdcoloring);
420: /* MatFDColoringView(matfdcoloring,PETSC_VIEWER_STDOUT_WORLD); */
421: SNESSetJacobian(snes,Jac,Jac,SNESDefaultComputeJacobianColor,matfdcoloring);
422: ISColoringDestroy(iscoloring);
423: }
425: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
426: Customize nonlinear solver; set runtime options
427: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
429: SNESSetFromOptions(snes);
431: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
432: Evaluate initial guess; then solve nonlinear system
433: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
435: /*
436: Note: The user should initialize the vector, x, with the initial guess
437: for the nonlinear solver prior to calling SNESSolve(). In particular,
438: to employ an initial guess of zero, the user should explicitly set
439: this vector to zero by calling VecSet().
440: */
441: FormInitialGuess(&user,x);
443: /*
444: Print the initial guess
445: */
446: VecGetArray(x,&xx);
447: for (inode = 0; inode < user.Nvlocal; inode++)
448: PetscFPrintf(PETSC_COMM_SELF,fptr1,"Initial Solution at node %D is %f \n",inode,xx[inode]);
449: VecRestoreArray(x,&xx);
451: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
452: Now solve the nonlinear system
453: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
455: SNESSolve(snes,PETSC_NULL,x);
456: SNESGetIterationNumber(snes,&its);
457: SNESGetNonlinearStepFailures(snes,&nfails);
458:
459: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
460: Print the output : solution vector and other information
461: Each processor writes to the file output.<rank> where rank is the
462: processor's rank.
463: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
465: VecGetArray(x,&xx);
466: for (inode = 0; inode < user.Nvlocal; inode++)
467: PetscFPrintf(PETSC_COMM_SELF,fptr1,"Solution at node %D is %f \n",inode,xx[inode]);
468: VecRestoreArray(x,&xx);
469: fclose(fptr1);
470: PetscPrintf(MPI_COMM_WORLD,"number of Newton iterations = %D, ",its);
471: PetscPrintf(MPI_COMM_WORLD,"number of unsuccessful steps = %D\n",nfails);
473: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
474: Free work space. All PETSc objects should be destroyed when they
475: are no longer needed.
476: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
477: PetscFree(user.gloInd);
478: PetscFree(user.locInd);
479: PetscFree(vertices);
480: PetscFree(verticesmask);
481: PetscFree(tmp);
482: VecScatterDestroy(user.scatter);
483: ISLocalToGlobalMappingDestroy(isl2g);
484: VecDestroy(x);
485: VecDestroy(r);
486: VecDestroy(user.localX);
487: VecDestroy(user.localF);
488: MatDestroy(Jac); SNESDestroy(snes);
489: /*PetscDrawDestroy(draw);*/
490: if (fd_jacobian_coloring){
491: MatFDColoringDestroy(matfdcoloring);
492: }
493: PetscFinalize();
495: return 0;
496: }
499: /* -------------------- Form initial approximation ----------------- */
501: /*
502: FormInitialGuess - Forms initial approximation.
504: Input Parameters:
505: user - user-defined application context
506: X - vector
508: Output Parameter:
509: X - vector
510: */
511: PetscErrorCode FormInitialGuess(AppCtx *user,Vec X)
512: {
513: PetscInt i,Nvlocal,ierr;
514: PetscInt *gloInd;
515: PetscScalar *x;
516: #if defined (UNUSED_VARIABLES)
517: PetscReal temp1,temp,hx,hy,hxdhy,hydhx,sc;
518: PetscInt Neglobal,Nvglobal,j,row;
519: PetscReal alpha,lambda;
521: Nvglobal = user->Nvglobal;
522: Neglobal = user->Neglobal;
523: lambda = user->non_lin_param;
524: alpha = user->lin_param;
525: #endif
527: Nvlocal = user->Nvlocal;
528: gloInd = user->gloInd;
530: /*
531: Get a pointer to vector data.
532: - For default PETSc vectors, VecGetArray() returns a pointer to
533: the data array. Otherwise, the routine is implementation dependent.
534: - You MUST call VecRestoreArray() when you no longer need access to
535: the array.
536: */
537: VecGetArray(X,&x);
539: /*
540: Compute initial guess over the locally owned part of the grid
541: */
542: for (i=0; i < Nvlocal; i++) {
543: x[i] = (PetscReal)gloInd[i];
544: }
546: /*
547: Restore vector
548: */
549: VecRestoreArray(X,&x);
550: return 0;
551: }
554: /* -------------------- Evaluate Function F(x) --------------------- */
555: /*
556: FormFunction - Evaluates nonlinear function, F(x).
558: Input Parameters:
559: . snes - the SNES context
560: . X - input vector
561: . ptr - optional user-defined context, as set by SNESSetFunction()
563: Output Parameter:
564: . F - function vector
565: */
566: PetscErrorCode FormFunction(SNES snes,Vec X,Vec F,void *ptr)
567: {
568: AppCtx *user = (AppCtx*)ptr;
570: PetscInt i,j,Nvlocal;
571: PetscReal alpha,lambda;
572: PetscScalar *x,*f;
573: VecScatter scatter;
574: Vec localX = user->localX;
575: #if defined (UNUSED_VARIABLES)
576: PetscScalar ut,ub,ul,ur,u,*g,sc,uyy,uxx;
577: PetscReal hx,hy,hxdhy,hydhx;
578: PetscReal two = 2.0,one = 1.0;
579: PetscInt Nvglobal,Neglobal,row;
580: PetscInt *gloInd;
582: Nvglobal = user->Nvglobal;
583: Neglobal = user->Neglobal;
584: gloInd = user->gloInd;
585: #endif
587: Nvlocal = user->Nvlocal;
588: lambda = user->non_lin_param;
589: alpha = user->lin_param;
590: scatter = user->scatter;
592: /*
593: PDE : L(u) + lambda*u*u +alpha*u = 0 where L(u) is the approximate Laplacian as
594: described in the beginning of this code
595:
596: First scatter the distributed vector X into local vector localX (that includes
597: values for ghost nodes. If we wish,we can put some other work between
598: VecScatterBegin() and VecScatterEnd() to overlap the communication with
599: computation.
600: */
601: VecScatterBegin(scatter,X,localX,INSERT_VALUES,SCATTER_FORWARD);
602: VecScatterEnd(scatter,X,localX,INSERT_VALUES,SCATTER_FORWARD);
604: /*
605: Get pointers to vector data
606: */
607: VecGetArray(localX,&x);
608: VecGetArray(F,&f);
610: /*
611: Now compute the f(x). As mentioned earlier, the computed Laplacian is just an
612: approximate one chosen for illustrative purpose only. Another point to notice
613: is that this is a local (completly parallel) calculation. In practical application
614: codes, function calculation time is a dominat portion of the overall execution time.
615: */
616: for (i=0; i < Nvlocal; i++) {
617: f[i] = (user->itot[i] - alpha)*x[i] - lambda*x[i]*x[i];
618: for (j = 0; j < user->itot[i]; j++) {
619: f[i] -= x[user->AdjM[i][j]];
620: }
621: }
623: /*
624: Restore vectors
625: */
626: VecRestoreArray(localX,&x);
627: VecRestoreArray(F,&f);
628: /*VecView(F,PETSC_VIEWER_STDOUT_WORLD);*/
630: return 0;
631: }
635: /* -------------------- Evaluate Jacobian F'(x) -------------------- */
636: /*
637: FormJacobian - Evaluates Jacobian matrix.
639: Input Parameters:
640: . snes - the SNES context
641: . X - input vector
642: . ptr - optional user-defined context, as set by SNESSetJacobian()
644: Output Parameters:
645: . A - Jacobian matrix
646: . B - optionally different preconditioning matrix
647: . flag - flag indicating matrix structure
649: */
650: PetscErrorCode FormJacobian(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
651: {
652: AppCtx *user = (AppCtx*)ptr;
653: Mat jac = *B;
654: PetscInt i,j,Nvlocal,col[50],ierr;
655: PetscScalar alpha,lambda,value[50];
656: Vec localX = user->localX;
657: VecScatter scatter;
658: PetscScalar *x;
659: #if defined (UNUSED_VARIABLES)
660: PetscScalar two = 2.0,one = 1.0;
661: PetscInt row,Nvglobal,Neglobal;
662: PetscInt *gloInd;
664: Nvglobal = user->Nvglobal;
665: Neglobal = user->Neglobal;
666: gloInd = user->gloInd;
667: #endif
668:
669: /*printf("Entering into FormJacobian \n");*/
670: Nvlocal = user->Nvlocal;
671: lambda = user->non_lin_param;
672: alpha = user->lin_param;
673: scatter = user->scatter;
675: /*
676: PDE : L(u) + lambda*u*u +alpha*u = 0 where L(u) is the approximate Laplacian as
677: described in the beginning of this code
678:
679: First scatter the distributed vector X into local vector localX (that includes
680: values for ghost nodes. If we wish, we can put some other work between
681: VecScatterBegin() and VecScatterEnd() to overlap the communication with
682: computation.
683: */
684: VecScatterBegin(scatter,X,localX,INSERT_VALUES,SCATTER_FORWARD);
685: VecScatterEnd(scatter,X,localX,INSERT_VALUES,SCATTER_FORWARD);
686:
687: /*
688: Get pointer to vector data
689: */
690: VecGetArray(localX,&x);
692: for (i=0; i < Nvlocal; i++) {
693: col[0] = i;
694: value[0] = user->itot[i] - 2.0*lambda*x[i] - alpha;
695: for (j = 0; j < user->itot[i]; j++) {
696: col[j+1] = user->AdjM[i][j];
697: value[j+1] = -1.0;
698: }
700: /*
701: Set the matrix values in the local ordering. Note that in order to use this
702: feature we must call the routine MatSetLocalToGlobalMapping() after the
703: matrix has been created.
704: */
705: MatSetValuesLocal(jac,1,&i,1+user->itot[i],col,value,INSERT_VALUES);
706: }
708: /*
709: Assemble matrix, using the 2-step process:
710: MatAssemblyBegin(), MatAssemblyEnd().
711: Between these two calls, the pointer to vector data has been restored to
712: demonstrate the use of overlapping communicationn with computation.
713: */
714: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
715: VecRestoreArray(localX,&x);
716: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
718: /*
719: Set flag to indicate that the Jacobian matrix retains an identical
720: nonzero structure throughout all nonlinear iterations (although the
721: values of the entries change). Thus, we can save some work in setting
722: up the preconditioner (e.g., no need to redo symbolic factorization for
723: ILU/ICC preconditioners).
724: - If the nonzero structure of the matrix is different during
725: successive linear solves, then the flag DIFFERENT_NONZERO_PATTERN
726: must be used instead. If you are unsure whether the matrix
727: structure has changed or not, use the flag DIFFERENT_NONZERO_PATTERN.
728: - Caution: If you specify SAME_NONZERO_PATTERN, PETSc
729: believes your assertion and does not check the structure
730: of the matrix. If you erroneously claim that the structure
731: is the same when it actually is not, the new preconditioner
732: will not function correctly. Thus, use this optimization
733: feature with caution!
734: */
735: *flag = SAME_NONZERO_PATTERN;
737: /*
738: Tell the matrix we will never add a new nonzero location to the
739: matrix. If we do, it will generate an error.
740: */
741: MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
742: /* MatView(jac,PETSC_VIEWER_STDOUT_SELF); */
743: return 0;
744: }
745: #else
746: #include <stdio.h>
747: int main(int argc,char **args)
748: {
749: fprintf(stdout,"This example does not work for complex numbers.\n");
750: return 0;
751: }
752: #endif