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