1 2 #ifndef lint 3 static char vcid[] = "$Id: snesj2.c,v 1.1 1996/10/09 19:38:25 bsmith Exp bsmith $"; 4 #endif 5 6 #include "src/mat/matimpl.h" /*I "mat.h" I*/ 7 #include "src/snes/snesimpl.h" /*I "snes.h" I*/ 8 9 10 /*@C 11 SNESDefaultComputeJacobianWithColoring 12 13 Input Parameters: 14 . snes - nonlinear solver object 15 . x1 - location at which to evaluate Jacobian 16 . ctx - MatFDColoring contex 17 18 Output Parameters: 19 . J - Jacobian matrix 20 . B - Jacobian preconditioner 21 . flag - flag indicating if the matrix nonzero structure has changed 22 23 .keywords: SNES, finite differences, Jacobian 24 25 .seealso: SNESSetJacobian(), SNESTestJacobian() 26 @*/ 27 int SNESDefaultComputeJacobianWithColoring(SNES snes,Vec x1,Mat *JJ,Mat *B,MatStructure *flag,void *ctx) 28 { 29 MatFDColoring color = (MatFDColoring) ctx; 30 Mat J = *B; 31 Vec jj1,jj2,x2; 32 int k, ierr,N,start,end,l,row,col,srow; 33 Scalar dx, mone = -1.0,*y,*scale = color->scale,*xx,*wscale = color->wscale; 34 double epsilon = color->error_rel, umin = color->umin; 35 MPI_Comm comm = color->comm; 36 37 ierr = MatZeroEntries(J); CHKERRQ(ierr); 38 if (!snes->nvwork) { 39 ierr = VecDuplicateVecs(x1,3,&snes->vwork); CHKERRQ(ierr); 40 snes->nvwork = 3; 41 PLogObjectParents(snes,3,snes->vwork); 42 } 43 jj1 = snes->vwork[0]; jj2 = snes->vwork[1]; x2 = snes->vwork[2]; 44 45 ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr); 46 ierr = VecGetSize(x1,&N); CHKERRQ(ierr); 47 ierr = VecGetArray(x1,&xx); CHKERRQ(ierr); 48 ierr = SNESComputeFunction(snes,x1,jj1); CHKERRQ(ierr); 49 50 PetscMemzero(wscale,N*sizeof(Scalar)); 51 /* 52 Loop over each color 53 */ 54 55 for (k=0; k<color->ncolors; k++) { 56 ierr = VecCopy(x1,x2); CHKERRQ(ierr); 57 /* 58 Loop over each column associated with color adding the 59 perturbation to the vector x2. 60 */ 61 for (l=0; l<color->ncolumns[k]; l++) { 62 col = color->columns[k][l]; /* column of the matrix we are probing for */ 63 dx = xx[col-start]; 64 #if !defined(PETSC_COMPLEX) 65 if (dx < umin && dx >= 0.0) dx = .1; 66 else if (dx < 0.0 && dx > -umin) dx = -.1; 67 #else 68 if (abs(dx) < umin && real(dx) >= 0.0) dx = .1; 69 else if (real(dx) < 0.0 && abs(dx) < umin) dx = -.1; 70 #endif 71 dx *= epsilon; 72 wscale[col] = 1.0/dx; 73 VecSetValues(x2,1,&col,&dx,ADD_VALUES); 74 } 75 VecRestoreArray(x1,&xx); 76 /* 77 Evaluate function at x1 + dx (here dx is a vector, of perturbations) 78 */ 79 ierr = SNESComputeFunction(snes,x2,jj2); CHKERRQ(ierr); 80 ierr = VecAXPY(&mone,jj1,jj2); CHKERRQ(ierr); 81 /* Communicate scale to all processors */ 82 #if !defined(PETSC_COMPLEX) 83 MPI_Allreduce(wscale,scale,N,MPI_DOUBLE,MPI_SUM,comm); 84 #else 85 MPI_Allreduce(wscale,scale,2*N,MPI_DOUBLE,MPI_SUM,comm); 86 #endif 87 /* 88 Loop over rows of vector putting results into Jacobian matrix 89 */ 90 VecGetArray(jj2,&y); 91 for (l=0; l<color->nrows[k]; l++) { 92 row = color->rows[k][l]; 93 col = color->columnsforrow[k][l]; 94 y[row] *= scale[col]; 95 srow = row + start; 96 ierr = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); 97 } 98 VecRestoreArray(jj2,&y); 99 } 100 ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 101 ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 102 *flag = SAME_NONZERO_PATTERN; 103 return 0; 104 } 105 106 107 108