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