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