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 72c4762a1bSJed Brown int main(int argc,char **argv) 73c4762a1bSJed Brown { 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 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 */ 163c4762a1bSJed Brown sprintf(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"); 170c4762a1bSJed Brown sscanf(str,"%d",&dtmp);user.v2p[inode] = dtmp; 171c4762a1bSJed Brown if (user.v2p[inode] == rank) { 172*63a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,fptr1,"Node %" PetscInt_FMT " belongs to processor %" PetscInt_FMT "\n",inode,user.v2p[inode])); 173c4762a1bSJed Brown 174c4762a1bSJed Brown user.gloInd[user.Nvlocal] = inode; 175c4762a1bSJed Brown sscanf(str,"%*d %d",&dtmp); 176c4762a1bSJed Brown nbrs = dtmp; 177*63a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,fptr1,"Number of neighbors for the vertex %" PetscInt_FMT " is %" PetscInt_FMT "\n",inode,nbrs)); 178c4762a1bSJed Brown 179c4762a1bSJed Brown user.itot[user.Nvlocal] = nbrs; 180c4762a1bSJed Brown Nvneighborstotal += nbrs; 181c4762a1bSJed Brown for (i = 0; i < user.itot[user.Nvlocal]; i++) { 182c4762a1bSJed Brown form[0]='\0'; 183c4762a1bSJed Brown for (j=0; j < i+2; j++) { 184b8abcfdeSJacob Faibussowitsch PetscCall(PetscStrlcat(form,"%*d ",sizeof(form))); 185c4762a1bSJed Brown } 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 191*63a3b9bcSJacob 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 } 197*63a3b9bcSJacob 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++) { 233*63a3b9bcSJacob 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++) { 240*63a3b9bcSJacob 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++) { 242*63a3b9bcSJacob 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")); 258c4762a1bSJed Brown for (i=0; i < user.Nvlocal; i++) { 259*63a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,fptr1," %" PetscInt_FMT " ",user.locInd[i])); 260c4762a1bSJed Brown } 261b8abcfdeSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,fptr1,"\n")); 262c4762a1bSJed Brown 263c4762a1bSJed Brown jstart = 0; 264c4762a1bSJed Brown for (i=0; i < user.Nvlocal; i++) { 265*63a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,fptr1,"Neghbors of local vertex %" PetscInt_FMT " are : ",user.locInd[i])); 266c4762a1bSJed Brown for (j=0; j < user.itot[i]; j++) { 267c4762a1bSJed Brown user.AdjM[i][j] = tmp[j+jstart]; 268c4762a1bSJed Brown 269*63a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,fptr1,"%" PetscInt_FMT " ",user.AdjM[i][j])); 270c4762a1bSJed Brown } 271c4762a1bSJed Brown jstart += user.itot[i]; 272b8abcfdeSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,fptr1,"\n")); 273c4762a1bSJed Brown } 274c4762a1bSJed Brown 275c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 276c4762a1bSJed Brown Extract the ghost vertex information for each processor 277c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 278c4762a1bSJed Brown /* 279c4762a1bSJed Brown Next, we need to generate a list of vertices required for this processor 280c4762a1bSJed Brown and a local numbering scheme for all vertices required on this processor. 281c4762a1bSJed Brown vertices - integer array of all vertices needed on this processor in PETSc 282c4762a1bSJed Brown global numbering; this list consists of first the "locally owned" 283c4762a1bSJed Brown vertices followed by the ghost vertices. 284c4762a1bSJed Brown verticesmask - integer array that for each global vertex lists its local 285c4762a1bSJed Brown vertex number (in vertices) + 1. If the global vertex is not 286c4762a1bSJed Brown represented on this processor, then the corresponding 287c4762a1bSJed Brown entry in verticesmask is zero 288c4762a1bSJed Brown 289c4762a1bSJed Brown Note: vertices and verticesmask are both Nvglobal in length; this may 290c4762a1bSJed Brown sound terribly non-scalable, but in fact is not so bad for a reasonable 291c4762a1bSJed Brown number of processors. Importantly, it allows us to use NO SEARCHING 292c4762a1bSJed Brown in setting up the data structures. 293c4762a1bSJed Brown */ 294b8abcfdeSJacob Faibussowitsch PetscCall(PetscMalloc1(user.Nvglobal,&vertices)); 295b8abcfdeSJacob Faibussowitsch PetscCall(PetscCalloc1(user.Nvglobal,&verticesmask)); 296c4762a1bSJed Brown nvertices = 0; 297c4762a1bSJed Brown 298c4762a1bSJed Brown /* 299c4762a1bSJed Brown First load "owned vertices" into list 300c4762a1bSJed Brown */ 301c4762a1bSJed Brown for (i=0; i < user.Nvlocal; i++) { 302c4762a1bSJed Brown vertices[nvertices++] = user.locInd[i]; 303c4762a1bSJed Brown verticesmask[user.locInd[i]] = nvertices; 304c4762a1bSJed Brown } 305c4762a1bSJed Brown 306c4762a1bSJed Brown /* 307c4762a1bSJed Brown Now load ghost vertices into list 308c4762a1bSJed Brown */ 309c4762a1bSJed Brown for (i=0; i < user.Nvlocal; i++) { 310c4762a1bSJed Brown for (j=0; j < user.itot[i]; j++) { 311c4762a1bSJed Brown nb = user.AdjM[i][j]; 312c4762a1bSJed Brown if (!verticesmask[nb]) { 313c4762a1bSJed Brown vertices[nvertices++] = nb; 314c4762a1bSJed Brown verticesmask[nb] = nvertices; 315c4762a1bSJed Brown } 316c4762a1bSJed Brown } 317c4762a1bSJed Brown } 318c4762a1bSJed Brown 319b8abcfdeSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,fptr1,"\n")); 320b8abcfdeSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,fptr1,"The array vertices is :\n")); 321c4762a1bSJed Brown for (i=0; i < nvertices; i++) { 322*63a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,fptr1,"%" PetscInt_FMT " ",vertices[i])); 323c4762a1bSJed Brown } 324b8abcfdeSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,fptr1,"\n")); 325c4762a1bSJed Brown 326c4762a1bSJed Brown /* 327c4762a1bSJed Brown Map the vertices listed in the neighbors to the local numbering from 328c4762a1bSJed Brown the global ordering that they contained initially. 329c4762a1bSJed Brown */ 330b8abcfdeSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,fptr1,"\n")); 331b8abcfdeSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,fptr1,"After mapping neighbors in the local contiguous ordering\n")); 332c4762a1bSJed Brown for (i=0; i<user.Nvlocal; i++) { 333*63a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,fptr1,"Neghbors of local vertex %" PetscInt_FMT " are :\n",i)); 334c4762a1bSJed Brown for (j = 0; j < user.itot[i]; j++) { 335c4762a1bSJed Brown nb = user.AdjM[i][j]; 336c4762a1bSJed Brown user.AdjM[i][j] = verticesmask[nb] - 1; 337c4762a1bSJed Brown 338*63a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,fptr1,"%" PetscInt_FMT " ",user.AdjM[i][j])); 339c4762a1bSJed Brown } 340b8abcfdeSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,fptr1,"\n")); 341c4762a1bSJed Brown } 342c4762a1bSJed Brown 343c4762a1bSJed Brown N = user.Nvglobal; 344c4762a1bSJed Brown 345c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 346c4762a1bSJed Brown Create vector and matrix data structures 347c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 348c4762a1bSJed Brown 349c4762a1bSJed Brown /* 350c4762a1bSJed Brown Create vector data structures 351c4762a1bSJed Brown */ 352b8abcfdeSJacob Faibussowitsch PetscCall(VecCreate(MPI_COMM_WORLD,&x)); 353b8abcfdeSJacob Faibussowitsch PetscCall(VecSetSizes(x,user.Nvlocal,N)); 354b8abcfdeSJacob Faibussowitsch PetscCall(VecSetFromOptions(x)); 355b8abcfdeSJacob Faibussowitsch PetscCall(VecDuplicate(x,&r)); 356b8abcfdeSJacob Faibussowitsch PetscCall(VecCreateSeq(MPI_COMM_SELF,bs*nvertices,&user.localX)); 357b8abcfdeSJacob Faibussowitsch PetscCall(VecDuplicate(user.localX,&user.localF)); 358c4762a1bSJed Brown 359c4762a1bSJed Brown /* 360c4762a1bSJed Brown Create the scatter between the global representation and the 361c4762a1bSJed Brown local representation 362c4762a1bSJed Brown */ 363b8abcfdeSJacob Faibussowitsch PetscCall(ISCreateStride(MPI_COMM_SELF,bs*nvertices,0,1,&islocal)); 364b8abcfdeSJacob Faibussowitsch PetscCall(ISCreateBlock(MPI_COMM_SELF,bs,nvertices,vertices,PETSC_COPY_VALUES,&isglobal)); 365b8abcfdeSJacob Faibussowitsch PetscCall(VecScatterCreate(x,isglobal,user.localX,islocal,&user.scatter)); 366b8abcfdeSJacob Faibussowitsch PetscCall(ISDestroy(&isglobal)); 367b8abcfdeSJacob Faibussowitsch PetscCall(ISDestroy(&islocal)); 368c4762a1bSJed Brown 369c4762a1bSJed Brown /* 370c4762a1bSJed Brown Create matrix data structure; Just to keep the example simple, we have not done any 371c4762a1bSJed Brown preallocation of memory for the matrix. In real application code with big matrices, 372c4762a1bSJed Brown preallocation should always be done to expedite the matrix creation. 373c4762a1bSJed Brown */ 374b8abcfdeSJacob Faibussowitsch PetscCall(MatCreate(MPI_COMM_WORLD,&Jac)); 375b8abcfdeSJacob Faibussowitsch PetscCall(MatSetSizes(Jac,PETSC_DECIDE,PETSC_DECIDE,N,N)); 376b8abcfdeSJacob Faibussowitsch PetscCall(MatSetFromOptions(Jac)); 377b8abcfdeSJacob Faibussowitsch PetscCall(MatSetUp(Jac)); 378c4762a1bSJed Brown 379c4762a1bSJed Brown /* 380c4762a1bSJed Brown The following routine allows us to set the matrix values in local ordering 381c4762a1bSJed Brown */ 382b8abcfdeSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(MPI_COMM_SELF,bs,nvertices,vertices,PETSC_COPY_VALUES,&isl2g)); 383b8abcfdeSJacob Faibussowitsch PetscCall(MatSetLocalToGlobalMapping(Jac,isl2g,isl2g)); 384c4762a1bSJed Brown 385c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 386c4762a1bSJed Brown Create nonlinear solver context 387c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 388c4762a1bSJed Brown 389b8abcfdeSJacob Faibussowitsch PetscCall(SNESCreate(MPI_COMM_WORLD,&snes)); 390b8abcfdeSJacob Faibussowitsch PetscCall(SNESSetType(snes,type)); 391c4762a1bSJed Brown 392c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 393c4762a1bSJed Brown Set routines for function and Jacobian evaluation 394c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 395b8abcfdeSJacob Faibussowitsch PetscCall(SNESSetFunction(snes,r,FormFunction,(void*)&user)); 396c4762a1bSJed Brown 397b8abcfdeSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-fd_jacobian_coloring",&fd_jacobian_coloring,0)); 398c4762a1bSJed Brown if (!fd_jacobian_coloring) { 399b8abcfdeSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes,Jac,Jac,FormJacobian,(void*)&user)); 400c4762a1bSJed Brown } else { /* Use matfdcoloring */ 401c4762a1bSJed Brown ISColoring iscoloring; 402c4762a1bSJed Brown MatColoring mc; 403c4762a1bSJed Brown 404c4762a1bSJed Brown /* Get the data structure of Jac */ 405b8abcfdeSJacob Faibussowitsch PetscCall(FormJacobian(snes,x,Jac,Jac,&user)); 406c4762a1bSJed Brown /* Create coloring context */ 407b8abcfdeSJacob Faibussowitsch PetscCall(MatColoringCreate(Jac,&mc)); 408b8abcfdeSJacob Faibussowitsch PetscCall(MatColoringSetType(mc,MATCOLORINGSL)); 409b8abcfdeSJacob Faibussowitsch PetscCall(MatColoringSetFromOptions(mc)); 410b8abcfdeSJacob Faibussowitsch PetscCall(MatColoringApply(mc,&iscoloring)); 411b8abcfdeSJacob Faibussowitsch PetscCall(MatColoringDestroy(&mc)); 412b8abcfdeSJacob Faibussowitsch PetscCall(MatFDColoringCreate(Jac,iscoloring,&matfdcoloring)); 413b8abcfdeSJacob Faibussowitsch PetscCall(MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))FormFunction,&user)); 414b8abcfdeSJacob Faibussowitsch PetscCall(MatFDColoringSetFromOptions(matfdcoloring)); 415b8abcfdeSJacob Faibussowitsch PetscCall(MatFDColoringSetUp(Jac,iscoloring,matfdcoloring)); 416b8abcfdeSJacob Faibussowitsch /* PetscCall(MatFDColoringView(matfdcoloring,PETSC_VIEWER_STDOUT_WORLD)); */ 417b8abcfdeSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes,Jac,Jac,SNESComputeJacobianDefaultColor,matfdcoloring)); 418b8abcfdeSJacob Faibussowitsch PetscCall(ISColoringDestroy(&iscoloring)); 419c4762a1bSJed Brown } 420c4762a1bSJed Brown 421c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 422c4762a1bSJed Brown Customize nonlinear solver; set runtime options 423c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 424c4762a1bSJed Brown 425b8abcfdeSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 426c4762a1bSJed Brown 427c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 428c4762a1bSJed Brown Evaluate initial guess; then solve nonlinear system 429c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 430c4762a1bSJed Brown 431c4762a1bSJed Brown /* 432c4762a1bSJed Brown Note: The user should initialize the vector, x, with the initial guess 433c4762a1bSJed Brown for the nonlinear solver prior to calling SNESSolve(). In particular, 434c4762a1bSJed Brown to employ an initial guess of zero, the user should explicitly set 435c4762a1bSJed Brown this vector to zero by calling VecSet(). 436c4762a1bSJed Brown */ 437b8abcfdeSJacob Faibussowitsch PetscCall(FormInitialGuess(&user,x)); 438c4762a1bSJed Brown 439c4762a1bSJed Brown /* 440c4762a1bSJed Brown Print the initial guess 441c4762a1bSJed Brown */ 442b8abcfdeSJacob Faibussowitsch PetscCall(VecGetArray(x,&xx)); 443c4762a1bSJed Brown for (inode = 0; inode < user.Nvlocal; inode++) { 444*63a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,fptr1,"Initial Solution at node %" PetscInt_FMT " is %f \n",inode,(double)PetscRealPart(xx[inode]))); 445c4762a1bSJed Brown } 446b8abcfdeSJacob Faibussowitsch PetscCall(VecRestoreArray(x,&xx)); 447c4762a1bSJed Brown 448c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 449c4762a1bSJed Brown Now solve the nonlinear system 450c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 451c4762a1bSJed Brown 452b8abcfdeSJacob Faibussowitsch PetscCall(SNESSolve(snes,NULL,x)); 453b8abcfdeSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes,&its)); 454b8abcfdeSJacob Faibussowitsch PetscCall(SNESGetNonlinearStepFailures(snes,&nfails)); 455c4762a1bSJed Brown 456c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 457c4762a1bSJed Brown Print the output : solution vector and other information 458c4762a1bSJed Brown Each processor writes to the file output.<rank> where rank is the 459c4762a1bSJed Brown processor's rank. 460c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 461c4762a1bSJed Brown 462b8abcfdeSJacob Faibussowitsch PetscCall(VecGetArray(x,&xx)); 463c4762a1bSJed Brown for (inode = 0; inode < user.Nvlocal; inode++) { 464*63a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_SELF,fptr1,"Solution at node %" PetscInt_FMT " is %f \n",inode,(double)PetscRealPart(xx[inode]))); 465c4762a1bSJed Brown } 466b8abcfdeSJacob Faibussowitsch PetscCall(VecRestoreArray(x,&xx)); 467c4762a1bSJed Brown fclose(fptr1); 468*63a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(MPI_COMM_WORLD,"number of SNES iterations = %" PetscInt_FMT ", ",its)); 469*63a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(MPI_COMM_WORLD,"number of unsuccessful steps = %" PetscInt_FMT "\n",nfails)); 470c4762a1bSJed Brown 471c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 472c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 473c4762a1bSJed Brown are no longer needed. 474c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 475b8abcfdeSJacob Faibussowitsch PetscCall(PetscFree(user.gloInd)); 476b8abcfdeSJacob Faibussowitsch PetscCall(PetscFree(user.locInd)); 477b8abcfdeSJacob Faibussowitsch PetscCall(PetscFree(vertices)); 478b8abcfdeSJacob Faibussowitsch PetscCall(PetscFree(verticesmask)); 479b8abcfdeSJacob Faibussowitsch PetscCall(PetscFree(tmp)); 480b8abcfdeSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user.scatter)); 481b8abcfdeSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&isl2g)); 482b8abcfdeSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 483b8abcfdeSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 484b8abcfdeSJacob Faibussowitsch PetscCall(VecDestroy(&user.localX)); 485b8abcfdeSJacob Faibussowitsch PetscCall(VecDestroy(&user.localF)); 486b8abcfdeSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 487b8abcfdeSJacob Faibussowitsch PetscCall(MatDestroy(&Jac)); 488b8abcfdeSJacob Faibussowitsch /* PetscCall(PetscDrawDestroy(draw));*/ 489b8abcfdeSJacob Faibussowitsch if (fd_jacobian_coloring) PetscCall(MatFDColoringDestroy(&matfdcoloring)); 490b8abcfdeSJacob Faibussowitsch PetscCall(PetscFinalize()); 491b8abcfdeSJacob Faibussowitsch return 0; 492c4762a1bSJed Brown } 493c4762a1bSJed Brown /* -------------------- Form initial approximation ----------------- */ 494c4762a1bSJed Brown 495c4762a1bSJed Brown /* 496c4762a1bSJed Brown FormInitialGuess - Forms initial approximation. 497c4762a1bSJed Brown 498c4762a1bSJed Brown Input Parameters: 499c4762a1bSJed Brown user - user-defined application context 500c4762a1bSJed Brown X - vector 501c4762a1bSJed Brown 502c4762a1bSJed Brown Output Parameter: 503c4762a1bSJed Brown X - vector 504c4762a1bSJed Brown */ 505c4762a1bSJed Brown PetscErrorCode FormInitialGuess(AppCtx *user,Vec X) 506c4762a1bSJed Brown { 507b8abcfdeSJacob Faibussowitsch PetscInt i,Nvlocal; 508c4762a1bSJed Brown PetscInt *gloInd; 509c4762a1bSJed Brown PetscScalar *x; 510c4762a1bSJed Brown #if defined(UNUSED_VARIABLES) 511c4762a1bSJed Brown PetscReal temp1,temp,hx,hy,hxdhy,hydhx,sc; 512c4762a1bSJed Brown PetscInt Neglobal,Nvglobal,j,row; 513c4762a1bSJed Brown PetscReal alpha,lambda; 514c4762a1bSJed Brown 515c4762a1bSJed Brown Nvglobal = user->Nvglobal; 516c4762a1bSJed Brown Neglobal = user->Neglobal; 517c4762a1bSJed Brown lambda = user->non_lin_param; 518c4762a1bSJed Brown alpha = user->lin_param; 519c4762a1bSJed Brown #endif 520c4762a1bSJed Brown 521c4762a1bSJed Brown Nvlocal = user->Nvlocal; 522c4762a1bSJed Brown gloInd = user->gloInd; 523c4762a1bSJed Brown 524c4762a1bSJed Brown /* 525c4762a1bSJed Brown Get a pointer to vector data. 526c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 527c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 528c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 529c4762a1bSJed Brown the array. 530c4762a1bSJed Brown */ 531b8abcfdeSJacob Faibussowitsch PetscCall(VecGetArray(X,&x)); 532c4762a1bSJed Brown 533c4762a1bSJed Brown /* 534c4762a1bSJed Brown Compute initial guess over the locally owned part of the grid 535c4762a1bSJed Brown */ 536c4762a1bSJed Brown for (i=0; i < Nvlocal; i++) x[i] = (PetscReal)gloInd[i]; 537c4762a1bSJed Brown 538c4762a1bSJed Brown /* 539c4762a1bSJed Brown Restore vector 540c4762a1bSJed Brown */ 541b8abcfdeSJacob Faibussowitsch PetscCall(VecRestoreArray(X,&x)); 542c4762a1bSJed Brown return 0; 543c4762a1bSJed Brown } 544c4762a1bSJed Brown /* -------------------- Evaluate Function F(x) --------------------- */ 545c4762a1bSJed Brown /* 546c4762a1bSJed Brown FormFunction - Evaluates nonlinear function, F(x). 547c4762a1bSJed Brown 548c4762a1bSJed Brown Input Parameters: 549c4762a1bSJed Brown . snes - the SNES context 550c4762a1bSJed Brown . X - input vector 551c4762a1bSJed Brown . ptr - optional user-defined context, as set by SNESSetFunction() 552c4762a1bSJed Brown 553c4762a1bSJed Brown Output Parameter: 554c4762a1bSJed Brown . F - function vector 555c4762a1bSJed Brown */ 556c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec X,Vec F,void *ptr) 557c4762a1bSJed Brown { 558c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 559c4762a1bSJed Brown PetscInt i,j,Nvlocal; 560c4762a1bSJed Brown PetscReal alpha,lambda; 561c4762a1bSJed Brown PetscScalar *x,*f; 562c4762a1bSJed Brown VecScatter scatter; 563c4762a1bSJed Brown Vec localX = user->localX; 564c4762a1bSJed Brown #if defined(UNUSED_VARIABLES) 565c4762a1bSJed Brown PetscScalar ut,ub,ul,ur,u,*g,sc,uyy,uxx; 566c4762a1bSJed Brown PetscReal hx,hy,hxdhy,hydhx; 567c4762a1bSJed Brown PetscReal two = 2.0,one = 1.0; 568c4762a1bSJed Brown PetscInt Nvglobal,Neglobal,row; 569c4762a1bSJed Brown PetscInt *gloInd; 570c4762a1bSJed Brown 571c4762a1bSJed Brown Nvglobal = user->Nvglobal; 572c4762a1bSJed Brown Neglobal = user->Neglobal; 573c4762a1bSJed Brown gloInd = user->gloInd; 574c4762a1bSJed Brown #endif 575c4762a1bSJed Brown 576c4762a1bSJed Brown Nvlocal = user->Nvlocal; 577c4762a1bSJed Brown lambda = user->non_lin_param; 578c4762a1bSJed Brown alpha = user->lin_param; 579c4762a1bSJed Brown scatter = user->scatter; 580c4762a1bSJed Brown 581c4762a1bSJed Brown /* 582c4762a1bSJed Brown PDE : L(u) + lambda*u*u +alpha*u = 0 where L(u) is the approximate Laplacian as 583c4762a1bSJed Brown described in the beginning of this code 584c4762a1bSJed Brown 585c4762a1bSJed Brown First scatter the distributed vector X into local vector localX (that includes 586c4762a1bSJed Brown values for ghost nodes. If we wish,we can put some other work between 587c4762a1bSJed Brown VecScatterBegin() and VecScatterEnd() to overlap the communication with 588c4762a1bSJed Brown computation. 589c4762a1bSJed Brown */ 590b8abcfdeSJacob Faibussowitsch PetscCall(VecScatterBegin(scatter,X,localX,INSERT_VALUES,SCATTER_FORWARD)); 591b8abcfdeSJacob Faibussowitsch PetscCall(VecScatterEnd(scatter,X,localX,INSERT_VALUES,SCATTER_FORWARD)); 592c4762a1bSJed Brown 593c4762a1bSJed Brown /* 594c4762a1bSJed Brown Get pointers to vector data 595c4762a1bSJed Brown */ 596b8abcfdeSJacob Faibussowitsch PetscCall(VecGetArray(localX,&x)); 597b8abcfdeSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 598c4762a1bSJed Brown 599c4762a1bSJed Brown /* 600c4762a1bSJed Brown Now compute the f(x). As mentioned earlier, the computed Laplacian is just an 601c4762a1bSJed Brown approximate one chosen for illustrative purpose only. Another point to notice 602c4762a1bSJed Brown is that this is a local (completly parallel) calculation. In practical application 603c4762a1bSJed Brown codes, function calculation time is a dominat portion of the overall execution time. 604c4762a1bSJed Brown */ 605c4762a1bSJed Brown for (i=0; i < Nvlocal; i++) { 606c4762a1bSJed Brown f[i] = (user->itot[i] - alpha)*x[i] - lambda*x[i]*x[i]; 607c4762a1bSJed Brown for (j = 0; j < user->itot[i]; j++) f[i] -= x[user->AdjM[i][j]]; 608c4762a1bSJed Brown } 609c4762a1bSJed Brown 610c4762a1bSJed Brown /* 611c4762a1bSJed Brown Restore vectors 612c4762a1bSJed Brown */ 613b8abcfdeSJacob Faibussowitsch PetscCall(VecRestoreArray(localX,&x)); 614b8abcfdeSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 615b8abcfdeSJacob Faibussowitsch /*PetscCall(VecView(F,PETSC_VIEWER_STDOUT_WORLD));*/ 616c4762a1bSJed Brown 617c4762a1bSJed Brown return 0; 618c4762a1bSJed Brown } 619c4762a1bSJed Brown 620c4762a1bSJed Brown /* -------------------- Evaluate Jacobian F'(x) -------------------- */ 621c4762a1bSJed Brown /* 622c4762a1bSJed Brown FormJacobian - Evaluates Jacobian matrix. 623c4762a1bSJed Brown 624c4762a1bSJed Brown Input Parameters: 625c4762a1bSJed Brown . snes - the SNES context 626c4762a1bSJed Brown . X - input vector 627c4762a1bSJed Brown . ptr - optional user-defined context, as set by SNESSetJacobian() 628c4762a1bSJed Brown 629c4762a1bSJed Brown Output Parameters: 630c4762a1bSJed Brown . A - Jacobian matrix 631c4762a1bSJed Brown . B - optionally different preconditioning matrix 632c4762a1bSJed Brown . flag - flag indicating matrix structure 633c4762a1bSJed Brown 634c4762a1bSJed Brown */ 635c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes,Vec X,Mat J,Mat jac,void *ptr) 636c4762a1bSJed Brown { 637c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 638b8abcfdeSJacob Faibussowitsch PetscInt i,j,Nvlocal,col[50]; 639c4762a1bSJed Brown PetscScalar alpha,lambda,value[50]; 640c4762a1bSJed Brown Vec localX = user->localX; 641c4762a1bSJed Brown VecScatter scatter; 642c4762a1bSJed Brown PetscScalar *x; 643c4762a1bSJed Brown #if defined(UNUSED_VARIABLES) 644c4762a1bSJed Brown PetscScalar two = 2.0,one = 1.0; 645c4762a1bSJed Brown PetscInt row,Nvglobal,Neglobal; 646c4762a1bSJed Brown PetscInt *gloInd; 647c4762a1bSJed Brown 648c4762a1bSJed Brown Nvglobal = user->Nvglobal; 649c4762a1bSJed Brown Neglobal = user->Neglobal; 650c4762a1bSJed Brown gloInd = user->gloInd; 651c4762a1bSJed Brown #endif 652c4762a1bSJed Brown 653c4762a1bSJed Brown /*printf("Entering into FormJacobian \n");*/ 654c4762a1bSJed Brown Nvlocal = user->Nvlocal; 655c4762a1bSJed Brown lambda = user->non_lin_param; 656c4762a1bSJed Brown alpha = user->lin_param; 657c4762a1bSJed Brown scatter = user->scatter; 658c4762a1bSJed Brown 659c4762a1bSJed Brown /* 660c4762a1bSJed Brown PDE : L(u) + lambda*u*u +alpha*u = 0 where L(u) is the approximate Laplacian as 661c4762a1bSJed Brown described in the beginning of this code 662c4762a1bSJed Brown 663c4762a1bSJed Brown First scatter the distributed vector X into local vector localX (that includes 664c4762a1bSJed Brown values for ghost nodes. If we wish, we can put some other work between 665c4762a1bSJed Brown VecScatterBegin() and VecScatterEnd() to overlap the communication with 666c4762a1bSJed Brown computation. 667c4762a1bSJed Brown */ 668b8abcfdeSJacob Faibussowitsch PetscCall(VecScatterBegin(scatter,X,localX,INSERT_VALUES,SCATTER_FORWARD)); 669b8abcfdeSJacob Faibussowitsch PetscCall(VecScatterEnd(scatter,X,localX,INSERT_VALUES,SCATTER_FORWARD)); 670c4762a1bSJed Brown 671c4762a1bSJed Brown /* 672c4762a1bSJed Brown Get pointer to vector data 673c4762a1bSJed Brown */ 674b8abcfdeSJacob Faibussowitsch PetscCall(VecGetArray(localX,&x)); 675c4762a1bSJed Brown 676c4762a1bSJed Brown for (i=0; i < Nvlocal; i++) { 677c4762a1bSJed Brown col[0] = i; 678c4762a1bSJed Brown value[0] = user->itot[i] - 2.0*lambda*x[i] - alpha; 679c4762a1bSJed Brown for (j = 0; j < user->itot[i]; j++) { 680c4762a1bSJed Brown col[j+1] = user->AdjM[i][j]; 681c4762a1bSJed Brown value[j+1] = -1.0; 682c4762a1bSJed Brown } 683c4762a1bSJed Brown 684c4762a1bSJed Brown /* 685c4762a1bSJed Brown Set the matrix values in the local ordering. Note that in order to use this 686c4762a1bSJed Brown feature we must call the routine MatSetLocalToGlobalMapping() after the 687c4762a1bSJed Brown matrix has been created. 688c4762a1bSJed Brown */ 689b8abcfdeSJacob Faibussowitsch PetscCall(MatSetValuesLocal(jac,1,&i,1+user->itot[i],col,value,INSERT_VALUES)); 690c4762a1bSJed Brown } 691c4762a1bSJed Brown 692c4762a1bSJed Brown /* 693c4762a1bSJed Brown Assemble matrix, using the 2-step process: 694c4762a1bSJed Brown MatAssemblyBegin(), MatAssemblyEnd(). 695c4762a1bSJed Brown Between these two calls, the pointer to vector data has been restored to 696c4762a1bSJed Brown demonstrate the use of overlapping communicationn with computation. 697c4762a1bSJed Brown */ 698b8abcfdeSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); 699b8abcfdeSJacob Faibussowitsch PetscCall(VecRestoreArray(localX,&x)); 700b8abcfdeSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); 701c4762a1bSJed Brown 702c4762a1bSJed Brown /* 703c4762a1bSJed Brown Tell the matrix we will never add a new nonzero location to the 704c4762a1bSJed Brown matrix. If we do, it will generate an error. 705c4762a1bSJed Brown */ 706b8abcfdeSJacob Faibussowitsch PetscCall(MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE)); 707c4762a1bSJed Brown /* MatView(jac,PETSC_VIEWER_STDOUT_SELF); */ 708c4762a1bSJed Brown return 0; 709c4762a1bSJed Brown } 710c4762a1bSJed Brown 711c4762a1bSJed Brown /*TEST 712c4762a1bSJed Brown 713c4762a1bSJed Brown build: 714c4762a1bSJed Brown requires: !complex 715c4762a1bSJed Brown 716c4762a1bSJed Brown test: 717c4762a1bSJed Brown nsize: 2 718c4762a1bSJed Brown args: -snes_monitor_short 719c4762a1bSJed Brown localrunfiles: options.inf adj.in 720c4762a1bSJed Brown 721c4762a1bSJed Brown test: 722c4762a1bSJed Brown suffix: 2 723c4762a1bSJed Brown nsize: 2 724c4762a1bSJed Brown args: -snes_monitor_short -fd_jacobian_coloring 725c4762a1bSJed Brown localrunfiles: options.inf adj.in 726c4762a1bSJed Brown 727c4762a1bSJed Brown TEST*/ 728