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