1*e090d566SSatish Balay /*$Id: snesj.c,v 1.65 2000/04/12 04:25:27 bsmith Exp balay $*/ 211320018SBarry Smith 3*e090d566SSatish Balay #include "src/snes/snesimpl.h" /*I "petscsnes.h" I*/ 411320018SBarry Smith 55615d1e5SSatish Balay #undef __FUNC__ 6b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESDefaultComputeJacobian" 74b828684SBarry Smith /*@C 884a7bf8dSLois Curfman McInnes SNESDefaultComputeJacobian - Computes the Jacobian using finite differences. 911320018SBarry Smith 10fee21e36SBarry Smith Collective on SNES 11fee21e36SBarry Smith 12c7afd0dbSLois Curfman McInnes Input Parameters: 13c7afd0dbSLois Curfman McInnes + x1 - compute Jacobian at this point 14c7afd0dbSLois Curfman McInnes - ctx - application's function context, as set with SNESSetFunction() 15c7afd0dbSLois Curfman McInnes 16c7afd0dbSLois Curfman McInnes Output Parameters: 17c7afd0dbSLois Curfman McInnes + J - Jacobian matrix (not altered in this routine) 18c7afd0dbSLois Curfman McInnes . B - newly computed Jacobian matrix to use with preconditioner (generally the same as J) 19c7afd0dbSLois Curfman McInnes - flag - flag indicating whether the matrix sparsity structure has changed 20c7afd0dbSLois Curfman McInnes 21ad960d00SLois Curfman McInnes Options Database Key: 22c7afd0dbSLois Curfman McInnes . -snes_fd - Activates SNESDefaultComputeJacobian() 23ad960d00SLois Curfman McInnes 245f3c43d9SLois Curfman McInnes Notes: 255f3c43d9SLois Curfman McInnes This routine is slow and expensive, and is not currently optimized 265f3c43d9SLois Curfman McInnes to take advantage of sparsity in the problem. Although 275f3c43d9SLois Curfman McInnes SNESDefaultComputeJacobian() is not recommended for general use 285f3c43d9SLois Curfman McInnes in large-scale applications, It can be useful in checking the 295f3c43d9SLois Curfman McInnes correctness of a user-provided Jacobian. 3011320018SBarry Smith 31b4fc646aSLois Curfman McInnes An alternative routine that uses coloring to explot matrix sparsity is 322d0c0e3bSBarry Smith SNESDefaultComputeJacobianColor(). 33b4fc646aSLois Curfman McInnes 3436851e7fSLois Curfman McInnes Level: intermediate 3536851e7fSLois Curfman McInnes 365f3c43d9SLois Curfman McInnes .keywords: SNES, finite differences, Jacobian 375f3c43d9SLois Curfman McInnes 382d0c0e3bSBarry Smith .seealso: SNESSetJacobian(), SNESDefaultComputeJacobianColor() 3911320018SBarry Smith @*/ 402d0c0e3bSBarry Smith int SNESDefaultComputeJacobian(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx) 4111320018SBarry Smith { 4288c956adSLois Curfman McInnes Vec j1a,j2a,x2; 4323242f5aSBarry Smith int i,ierr,N,start,end,j; 44bbb6d6a8SBarry Smith Scalar dx,mone = -1.0,*y,scale,*xx,wscale; 45329f5518SBarry Smith PetscReal amax,epsilon = 1.e-8; /* assumes PetscReal precision */ 46329f5518SBarry Smith PetscReal dx_min = 1.e-16,dx_par = 1.e-1; 47bbb6d6a8SBarry Smith MPI_Comm comm; 48a305c92eSSatish Balay int (*eval_fct)(SNES,Vec,Vec)=0; 490521c3abSLois Curfman McInnes 503a40ed3dSBarry Smith PetscFunctionBegin; 513a40ed3dSBarry Smith if (snes->method_class == SNES_NONLINEAR_EQUATIONS) eval_fct = SNESComputeFunction; 523a40ed3dSBarry Smith else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) eval_fct = SNESComputeGradient; 53a8c6a408SBarry Smith else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Invalid method class"); 5423242f5aSBarry Smith 552d0c0e3bSBarry Smith ierr = PetscObjectGetComm((PetscObject)x1,&comm);CHKERRQ(ierr); 562d0c0e3bSBarry Smith ierr = MatZeroEntries(*B);CHKERRQ(ierr); 57aa79bc6dSLois Curfman McInnes if (!snes->nvwork) { 58aa79bc6dSLois Curfman McInnes ierr = VecDuplicateVecs(x1,3,&snes->vwork);CHKERRQ(ierr); 59aa79bc6dSLois Curfman McInnes snes->nvwork = 3; 60aa79bc6dSLois Curfman McInnes PLogObjectParents(snes,3,snes->vwork); 61aa79bc6dSLois Curfman McInnes } 6288c956adSLois Curfman McInnes j1a = snes->vwork[0]; j2a = snes->vwork[1]; x2 = snes->vwork[2]; 6323242f5aSBarry Smith 6478b31e54SBarry Smith ierr = VecGetSize(x1,&N);CHKERRQ(ierr); 6578b31e54SBarry Smith ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr); 6688c956adSLois Curfman McInnes ierr = eval_fct(snes,x1,j1a);CHKERRQ(ierr); 67c005e166SLois Curfman McInnes 68c005e166SLois Curfman McInnes /* Compute Jacobian approximation, 1 column at a time. 6988c956adSLois Curfman McInnes x1 = current iterate, j1a = F(x1) 7088c956adSLois Curfman McInnes x2 = perturbed iterate, j2a = F(x2) 71c005e166SLois Curfman McInnes */ 7239e2f89bSBarry Smith for (i=0; i<N; i++) { 7378b31e54SBarry Smith ierr = VecCopy(x1,x2);CHKERRQ(ierr); 7423242f5aSBarry Smith if (i>= start && i<end) { 752e8a6d31SBarry Smith ierr = VecGetArray(x1,&xx);CHKERRQ(ierr); 7639e2f89bSBarry Smith dx = xx[i-start]; 772e8a6d31SBarry Smith ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 78aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX) 7954d73c9fSLois Curfman McInnes if (dx < dx_min && dx >= 0.0) dx = dx_par; 8054d73c9fSLois Curfman McInnes else if (dx < 0.0 && dx > -dx_min) dx = -dx_par; 8119a167f6SBarry Smith #else 82329f5518SBarry Smith if (PetscAbsScalar(dx) < dx_min && PetscRealPart(dx) >= 0.0) dx = dx_par; 83329f5518SBarry Smith else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < dx_min) dx = -dx_par; 8419a167f6SBarry Smith #endif 8539e2f89bSBarry Smith dx *= epsilon; 8674f6f00dSLois Curfman McInnes wscale = 1.0/dx; 872e8a6d31SBarry Smith ierr = VecSetValues(x2,1,&i,&dx,ADD_VALUES);CHKERRQ(ierr); 886c783aadSBarry Smith } else { 89bbb6d6a8SBarry Smith wscale = 0.0; 90bbb6d6a8SBarry Smith } 9188c956adSLois Curfman McInnes ierr = eval_fct(snes,x2,j2a);CHKERRQ(ierr); 9288c956adSLois Curfman McInnes ierr = VecAXPY(&mone,j1a,j2a);CHKERRQ(ierr); 93c005e166SLois Curfman McInnes /* Communicate scale to all processors */ 946831982aSBarry Smith ierr = MPI_Allreduce(&wscale,&scale,1,MPIU_SCALAR,PetscSum_Op,comm);CHKERRQ(ierr); 952e8a6d31SBarry Smith ierr = VecScale(&scale,j2a);CHKERRQ(ierr); 962e8a6d31SBarry Smith ierr = VecGetArray(j2a,&y);CHKERRQ(ierr); 972e8a6d31SBarry Smith ierr = VecNorm(j2a,NORM_INFINITY,&amax);CHKERRQ(ierr); amax *= 1.e-14; 9823242f5aSBarry Smith for (j=start; j<end; j++) { 99cddf8d76SBarry Smith if (PetscAbsScalar(y[j-start]) > amax) { 100b4fc646aSLois Curfman McInnes ierr = MatSetValues(*B,1,&j,1,&i,y+j-start,INSERT_VALUES);CHKERRQ(ierr); 10123242f5aSBarry Smith } 10223242f5aSBarry Smith } 1032e8a6d31SBarry Smith ierr = VecRestoreArray(j2a,&y);CHKERRQ(ierr); 10423242f5aSBarry Smith } 105b4fc646aSLois Curfman McInnes ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 106b4fc646aSLois Curfman McInnes ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 107595b82d2SBarry Smith *flag = DIFFERENT_NONZERO_PATTERN; 1083a40ed3dSBarry Smith PetscFunctionReturn(0); 10911320018SBarry Smith } 11011320018SBarry Smith 1115615d1e5SSatish Balay #undef __FUNC__ 112b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESDefaultComputeHessian" 1130521c3abSLois Curfman McInnes /*@C 11484a7bf8dSLois Curfman McInnes SNESDefaultComputeHessian - Computes the Hessian using finite differences. 1150521c3abSLois Curfman McInnes 116fee21e36SBarry Smith Collective on SNES 117fee21e36SBarry Smith 118c7afd0dbSLois Curfman McInnes Input Parameters: 119c7afd0dbSLois Curfman McInnes + x1 - compute Hessian at this point 120c7afd0dbSLois Curfman McInnes - ctx - application's gradient context, as set with SNESSetGradient() 121c7afd0dbSLois Curfman McInnes 122c7afd0dbSLois Curfman McInnes Output Parameters: 123c7afd0dbSLois Curfman McInnes + J - Hessian matrix (not altered in this routine) 124c7afd0dbSLois Curfman McInnes . B - newly computed Hessian matrix to use with preconditioner (generally the same as J) 125c7afd0dbSLois Curfman McInnes - flag - flag indicating whether the matrix sparsity structure has changed 126c7afd0dbSLois Curfman McInnes 1270521c3abSLois Curfman McInnes Options Database Key: 128c7afd0dbSLois Curfman McInnes $ -snes_fd - Activates SNESDefaultComputeHessian() 1290521c3abSLois Curfman McInnes 13015091d37SBarry Smith 13115091d37SBarry Smith Level: intermediate 13215091d37SBarry Smith 1330521c3abSLois Curfman McInnes Notes: 1340521c3abSLois Curfman McInnes This routine is slow and expensive, and is not currently optimized 1350521c3abSLois Curfman McInnes to take advantage of sparsity in the problem. Although 1360521c3abSLois Curfman McInnes SNESDefaultComputeHessian() is not recommended for general use 1370521c3abSLois Curfman McInnes in large-scale applications, It can be useful in checking the 1380521c3abSLois Curfman McInnes correctness of a user-provided Hessian. 1390521c3abSLois Curfman McInnes 1400521c3abSLois Curfman McInnes .keywords: SNES, finite differences, Hessian 1410521c3abSLois Curfman McInnes 1429984d6fbSBarry Smith .seealso: SNESSetHessian() 1430521c3abSLois Curfman McInnes @*/ 144329f5518SBarry Smith int SNESDefaultComputeHessian(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx) 1450521c3abSLois Curfman McInnes { 1463a40ed3dSBarry Smith int ierr; 1473a40ed3dSBarry Smith 1483a40ed3dSBarry Smith PetscFunctionBegin; 1493a40ed3dSBarry Smith ierr = SNESDefaultComputeJacobian(snes,x1,J,B,flag,ctx);CHKERRQ(ierr); 1503a40ed3dSBarry Smith PetscFunctionReturn(0); 1510521c3abSLois Curfman McInnes } 152bd45824fSLois Curfman McInnes 153bd45824fSLois Curfman McInnes #undef __FUNC__ 154b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESDefaultComputeHessianColor" 155bd45824fSLois Curfman McInnes /*@C 156bd45824fSLois Curfman McInnes SNESDefaultComputeHessianColor - Computes the Hessian using colored finite differences. 157bd45824fSLois Curfman McInnes 158bd45824fSLois Curfman McInnes Collective on SNES 159bd45824fSLois Curfman McInnes 160bd45824fSLois Curfman McInnes Input Parameters: 161bd45824fSLois Curfman McInnes + x1 - compute Hessian at this point 162bd45824fSLois Curfman McInnes - ctx - application's gradient context, as set with SNESSetGradient() 163bd45824fSLois Curfman McInnes 164bd45824fSLois Curfman McInnes Output Parameters: 165bd45824fSLois Curfman McInnes + J - Hessian matrix (not altered in this routine) 166bd45824fSLois Curfman McInnes . B - newly computed Hessian matrix to use with preconditioner (generally the same as J) 167bd45824fSLois Curfman McInnes - flag - flag indicating whether the matrix sparsity structure has changed 168bd45824fSLois Curfman McInnes 169bd45824fSLois Curfman McInnes Options Database Keys: 170bd45824fSLois Curfman McInnes . -mat_fd_coloring_freq <freq> - Activates SNESDefaultComputeJacobianColor() 171bd45824fSLois Curfman McInnes 172bd45824fSLois Curfman McInnes Level: intermediate 173bd45824fSLois Curfman McInnes 174bd45824fSLois Curfman McInnes .keywords: SNES, finite differences, Hessian, coloring, sparse 175bd45824fSLois Curfman McInnes 176bd45824fSLois Curfman McInnes .seealso: SNESSetHessian() 177bd45824fSLois Curfman McInnes @*/ 178329f5518SBarry Smith int SNESDefaultComputeHessianColor(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx) 179bd45824fSLois Curfman McInnes { 180bd45824fSLois Curfman McInnes int ierr; 181bd45824fSLois Curfman McInnes 182bd45824fSLois Curfman McInnes PetscFunctionBegin; 183bd45824fSLois Curfman McInnes ierr = SNESDefaultComputeJacobianColor(snes,x1,J,B,flag,ctx);CHKERRQ(ierr); 184bd45824fSLois Curfman McInnes PetscFunctionReturn(0); 185bd45824fSLois Curfman McInnes } 186bd45824fSLois Curfman McInnes 187