1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*6831982aSBarry Smith static char vcid[] = "$Id: snesj.c,v 1.60 1999/10/22 21:16:38 curfman Exp bsmith $"; 311320018SBarry Smith #endif 411320018SBarry Smith 570f55243SBarry Smith #include "src/snes/snesimpl.h" /*I "snes.h" I*/ 611320018SBarry Smith 75615d1e5SSatish Balay #undef __FUNC__ 85615d1e5SSatish Balay #define __FUNC__ "SNESDefaultComputeJacobian" 94b828684SBarry Smith /*@C 1084a7bf8dSLois Curfman McInnes SNESDefaultComputeJacobian - Computes the Jacobian using finite differences. 1111320018SBarry Smith 12fee21e36SBarry Smith Collective on SNES 13fee21e36SBarry Smith 14c7afd0dbSLois Curfman McInnes Input Parameters: 15c7afd0dbSLois Curfman McInnes + x1 - compute Jacobian at this point 16c7afd0dbSLois Curfman McInnes - ctx - application's function context, as set with SNESSetFunction() 17c7afd0dbSLois Curfman McInnes 18c7afd0dbSLois Curfman McInnes Output Parameters: 19c7afd0dbSLois Curfman McInnes + J - Jacobian matrix (not altered in this routine) 20c7afd0dbSLois Curfman McInnes . B - newly computed Jacobian matrix to use with preconditioner (generally the same as J) 21c7afd0dbSLois Curfman McInnes - flag - flag indicating whether the matrix sparsity structure has changed 22c7afd0dbSLois Curfman McInnes 23ad960d00SLois Curfman McInnes Options Database Key: 24c7afd0dbSLois Curfman McInnes . -snes_fd - Activates SNESDefaultComputeJacobian() 25ad960d00SLois Curfman McInnes 265f3c43d9SLois Curfman McInnes Notes: 275f3c43d9SLois Curfman McInnes This routine is slow and expensive, and is not currently optimized 285f3c43d9SLois Curfman McInnes to take advantage of sparsity in the problem. Although 295f3c43d9SLois Curfman McInnes SNESDefaultComputeJacobian() is not recommended for general use 305f3c43d9SLois Curfman McInnes in large-scale applications, It can be useful in checking the 315f3c43d9SLois Curfman McInnes correctness of a user-provided Jacobian. 3211320018SBarry Smith 33b4fc646aSLois Curfman McInnes An alternative routine that uses coloring to explot matrix sparsity is 342d0c0e3bSBarry Smith SNESDefaultComputeJacobianColor(). 35b4fc646aSLois Curfman McInnes 3636851e7fSLois Curfman McInnes Level: intermediate 3736851e7fSLois Curfman McInnes 385f3c43d9SLois Curfman McInnes .keywords: SNES, finite differences, Jacobian 395f3c43d9SLois Curfman McInnes 402d0c0e3bSBarry Smith .seealso: SNESSetJacobian(), SNESDefaultComputeJacobianColor() 4111320018SBarry Smith @*/ 422d0c0e3bSBarry Smith int SNESDefaultComputeJacobian(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx) 4311320018SBarry Smith { 4488c956adSLois Curfman McInnes Vec j1a,j2a,x2; 4523242f5aSBarry Smith int i,ierr,N,start,end,j; 46bbb6d6a8SBarry Smith Scalar dx, mone = -1.0,*y,scale,*xx,wscale; 475334005bSBarry Smith double amax, epsilon = 1.e-8; /* assumes double precision */ 4854d73c9fSLois Curfman McInnes double dx_min = 1.e-16, dx_par = 1.e-1; 49bbb6d6a8SBarry Smith MPI_Comm comm; 50a305c92eSSatish Balay int (*eval_fct)(SNES,Vec,Vec)=0; 510521c3abSLois Curfman McInnes 523a40ed3dSBarry Smith PetscFunctionBegin; 533a40ed3dSBarry Smith if (snes->method_class == SNES_NONLINEAR_EQUATIONS) eval_fct = SNESComputeFunction; 543a40ed3dSBarry Smith else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) eval_fct = SNESComputeGradient; 55a8c6a408SBarry Smith else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Invalid method class"); 5623242f5aSBarry Smith 572d0c0e3bSBarry Smith ierr = PetscObjectGetComm((PetscObject)x1,&comm);CHKERRQ(ierr); 582d0c0e3bSBarry Smith ierr = MatZeroEntries(*B);CHKERRQ(ierr); 59aa79bc6dSLois Curfman McInnes if (!snes->nvwork) { 60aa79bc6dSLois Curfman McInnes ierr = VecDuplicateVecs(x1,3,&snes->vwork);CHKERRQ(ierr); 61aa79bc6dSLois Curfman McInnes snes->nvwork = 3; 62aa79bc6dSLois Curfman McInnes PLogObjectParents(snes,3,snes->vwork); 63aa79bc6dSLois Curfman McInnes } 6488c956adSLois Curfman McInnes j1a = snes->vwork[0]; j2a = snes->vwork[1]; x2 = snes->vwork[2]; 6523242f5aSBarry Smith 6678b31e54SBarry Smith ierr = VecGetSize(x1,&N);CHKERRQ(ierr); 6778b31e54SBarry Smith ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr); 6888c956adSLois Curfman McInnes ierr = eval_fct(snes,x1,j1a);CHKERRQ(ierr); 69c005e166SLois Curfman McInnes 70c005e166SLois Curfman McInnes /* Compute Jacobian approximation, 1 column at a time. 7188c956adSLois Curfman McInnes x1 = current iterate, j1a = F(x1) 7288c956adSLois Curfman McInnes x2 = perturbed iterate, j2a = F(x2) 73c005e166SLois Curfman McInnes */ 7439e2f89bSBarry Smith for ( i=0; i<N; i++ ) { 7578b31e54SBarry Smith ierr = VecCopy(x1,x2);CHKERRQ(ierr); 7623242f5aSBarry Smith if ( i>= start && i<end) { 772e8a6d31SBarry Smith ierr = VecGetArray(x1,&xx);CHKERRQ(ierr); 7839e2f89bSBarry Smith dx = xx[i-start]; 792e8a6d31SBarry Smith ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 80aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX) 8154d73c9fSLois Curfman McInnes if (dx < dx_min && dx >= 0.0) dx = dx_par; 8254d73c9fSLois Curfman McInnes else if (dx < 0.0 && dx > -dx_min) dx = -dx_par; 8319a167f6SBarry Smith #else 84e20fef11SSatish Balay if (PetscAbsScalar(dx) < dx_min && PetscReal(dx) >= 0.0) dx = dx_par; 85e20fef11SSatish Balay else if (PetscReal(dx) < 0.0 && PetscAbsScalar(dx) < dx_min) dx = -dx_par; 8619a167f6SBarry Smith #endif 8739e2f89bSBarry Smith dx *= epsilon; 8874f6f00dSLois Curfman McInnes wscale = 1.0/dx; 892e8a6d31SBarry Smith ierr = VecSetValues(x2,1,&i,&dx,ADD_VALUES);CHKERRQ(ierr); 906c783aadSBarry Smith } else { 91bbb6d6a8SBarry Smith wscale = 0.0; 92bbb6d6a8SBarry Smith } 9388c956adSLois Curfman McInnes ierr = eval_fct(snes,x2,j2a);CHKERRQ(ierr); 9488c956adSLois Curfman McInnes ierr = VecAXPY(&mone,j1a,j2a);CHKERRQ(ierr); 95c005e166SLois Curfman McInnes /* Communicate scale to all processors */ 96*6831982aSBarry Smith ierr = MPI_Allreduce(&wscale,&scale,1,MPIU_SCALAR,PetscSum_Op,comm);CHKERRQ(ierr); 972e8a6d31SBarry Smith ierr = VecScale(&scale,j2a);CHKERRQ(ierr); 982e8a6d31SBarry Smith ierr = VecGetArray(j2a,&y);CHKERRQ(ierr); 992e8a6d31SBarry Smith ierr = VecNorm(j2a,NORM_INFINITY,&amax);CHKERRQ(ierr); amax *= 1.e-14; 10023242f5aSBarry Smith for ( j=start; j<end; j++ ) { 101cddf8d76SBarry Smith if (PetscAbsScalar(y[j-start]) > amax) { 102b4fc646aSLois Curfman McInnes ierr = MatSetValues(*B,1,&j,1,&i,y+j-start,INSERT_VALUES);CHKERRQ(ierr); 10323242f5aSBarry Smith } 10423242f5aSBarry Smith } 1052e8a6d31SBarry Smith ierr = VecRestoreArray(j2a,&y);CHKERRQ(ierr); 10623242f5aSBarry Smith } 107b4fc646aSLois Curfman McInnes ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 108b4fc646aSLois Curfman McInnes ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 109595b82d2SBarry Smith *flag = DIFFERENT_NONZERO_PATTERN; 1103a40ed3dSBarry Smith PetscFunctionReturn(0); 11111320018SBarry Smith } 11211320018SBarry Smith 1135615d1e5SSatish Balay #undef __FUNC__ 1145615d1e5SSatish Balay #define __FUNC__ "SNESDefaultComputeHessian" 1150521c3abSLois Curfman McInnes /*@C 11684a7bf8dSLois Curfman McInnes SNESDefaultComputeHessian - Computes the Hessian using finite differences. 1170521c3abSLois Curfman McInnes 118fee21e36SBarry Smith Collective on SNES 119fee21e36SBarry Smith 120c7afd0dbSLois Curfman McInnes Input Parameters: 121c7afd0dbSLois Curfman McInnes + x1 - compute Hessian at this point 122c7afd0dbSLois Curfman McInnes - ctx - application's gradient context, as set with SNESSetGradient() 123c7afd0dbSLois Curfman McInnes 124c7afd0dbSLois Curfman McInnes Output Parameters: 125c7afd0dbSLois Curfman McInnes + J - Hessian matrix (not altered in this routine) 126c7afd0dbSLois Curfman McInnes . B - newly computed Hessian matrix to use with preconditioner (generally the same as J) 127c7afd0dbSLois Curfman McInnes - flag - flag indicating whether the matrix sparsity structure has changed 128c7afd0dbSLois Curfman McInnes 1290521c3abSLois Curfman McInnes Options Database Key: 130c7afd0dbSLois Curfman McInnes $ -snes_fd - Activates SNESDefaultComputeHessian() 1310521c3abSLois Curfman McInnes 13215091d37SBarry Smith 13315091d37SBarry Smith Level: intermediate 13415091d37SBarry Smith 1350521c3abSLois Curfman McInnes Notes: 1360521c3abSLois Curfman McInnes This routine is slow and expensive, and is not currently optimized 1370521c3abSLois Curfman McInnes to take advantage of sparsity in the problem. Although 1380521c3abSLois Curfman McInnes SNESDefaultComputeHessian() is not recommended for general use 1390521c3abSLois Curfman McInnes in large-scale applications, It can be useful in checking the 1400521c3abSLois Curfman McInnes correctness of a user-provided Hessian. 1410521c3abSLois Curfman McInnes 1420521c3abSLois Curfman McInnes .keywords: SNES, finite differences, Hessian 1430521c3abSLois Curfman McInnes 1449984d6fbSBarry Smith .seealso: SNESSetHessian() 1450521c3abSLois Curfman McInnes @*/ 146c7afd0dbSLois Curfman McInnes int SNESDefaultComputeHessian(SNES snes,Vec x1,Mat *J,Mat *B, 147c7afd0dbSLois Curfman McInnes MatStructure *flag,void *ctx) 1480521c3abSLois Curfman McInnes { 1493a40ed3dSBarry Smith int ierr; 1503a40ed3dSBarry Smith 1513a40ed3dSBarry Smith PetscFunctionBegin; 1523a40ed3dSBarry Smith ierr = SNESDefaultComputeJacobian(snes,x1,J,B,flag,ctx);CHKERRQ(ierr); 1533a40ed3dSBarry Smith PetscFunctionReturn(0); 1540521c3abSLois Curfman McInnes } 155bd45824fSLois Curfman McInnes 156bd45824fSLois Curfman McInnes #undef __FUNC__ 157bd45824fSLois Curfman McInnes #define __FUNC__ "SNESDefaultComputeHessianColor" 158bd45824fSLois Curfman McInnes /*@C 159bd45824fSLois Curfman McInnes SNESDefaultComputeHessianColor - Computes the Hessian using colored finite differences. 160bd45824fSLois Curfman McInnes 161bd45824fSLois Curfman McInnes Collective on SNES 162bd45824fSLois Curfman McInnes 163bd45824fSLois Curfman McInnes Input Parameters: 164bd45824fSLois Curfman McInnes + x1 - compute Hessian at this point 165bd45824fSLois Curfman McInnes - ctx - application's gradient context, as set with SNESSetGradient() 166bd45824fSLois Curfman McInnes 167bd45824fSLois Curfman McInnes Output Parameters: 168bd45824fSLois Curfman McInnes + J - Hessian matrix (not altered in this routine) 169bd45824fSLois Curfman McInnes . B - newly computed Hessian matrix to use with preconditioner (generally the same as J) 170bd45824fSLois Curfman McInnes - flag - flag indicating whether the matrix sparsity structure has changed 171bd45824fSLois Curfman McInnes 172bd45824fSLois Curfman McInnes Options Database Keys: 173bd45824fSLois Curfman McInnes . -mat_fd_coloring_freq <freq> - Activates SNESDefaultComputeJacobianColor() 174bd45824fSLois Curfman McInnes 175bd45824fSLois Curfman McInnes Level: intermediate 176bd45824fSLois Curfman McInnes 177bd45824fSLois Curfman McInnes .keywords: SNES, finite differences, Hessian, coloring, sparse 178bd45824fSLois Curfman McInnes 179bd45824fSLois Curfman McInnes .seealso: SNESSetHessian() 180bd45824fSLois Curfman McInnes @*/ 181bd45824fSLois Curfman McInnes int SNESDefaultComputeHessianColor(SNES snes,Vec x1,Mat *J,Mat *B, 182bd45824fSLois Curfman McInnes MatStructure *flag,void *ctx) 183bd45824fSLois Curfman McInnes { 184bd45824fSLois Curfman McInnes int ierr; 185bd45824fSLois Curfman McInnes 186bd45824fSLois Curfman McInnes PetscFunctionBegin; 187bd45824fSLois Curfman McInnes ierr = SNESDefaultComputeJacobianColor(snes,x1,J,B,flag,ctx);CHKERRQ(ierr); 188bd45824fSLois Curfman McInnes PetscFunctionReturn(0); 189bd45824fSLois Curfman McInnes } 190bd45824fSLois Curfman McInnes 191