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