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