1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*3a40ed3dSBarry Smith static char vcid[] = "$Id: snesj.c,v 1.44 1997/09/25 22:42:32 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 1211320018SBarry Smith Input Parameters: 13da88328cSLois Curfman McInnes . x1 - compute Jacobian at this point 14da88328cSLois Curfman McInnes . ctx - application's function context, as set with SNESSetFunction() 1511320018SBarry Smith 1611320018SBarry Smith Output Parameters: 17b4fc646aSLois Curfman McInnes . J - Jacobian matrix (not altered in this routine) 18b4fc646aSLois Curfman McInnes . B - newly computed Jacobian matrix to use with preconditioner (generally the same as J) 19b4fc646aSLois Curfman McInnes . flag - flag indicating whether the matrix sparsity structure has changed 2011320018SBarry Smith 21ad960d00SLois Curfman McInnes Options Database Key: 22ad960d00SLois Curfman McInnes $ -snes_fd 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 32b4fc646aSLois Curfman McInnes SNESDefaultComputeJacobianWithColoring(). 33b4fc646aSLois Curfman McInnes 345f3c43d9SLois Curfman McInnes .keywords: SNES, finite differences, Jacobian 355f3c43d9SLois Curfman McInnes 36b4fc646aSLois Curfman McInnes .seealso: SNESSetJacobian(), SNESDefaultComputeJacobianWithColoring() 3711320018SBarry Smith @*/ 385334005bSBarry Smith int SNESDefaultComputeJacobian(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx) 3911320018SBarry Smith { 4023242f5aSBarry Smith Vec j1,j2,x2; 4123242f5aSBarry Smith int i,ierr,N,start,end,j; 42bbb6d6a8SBarry Smith Scalar dx, mone = -1.0,*y,scale,*xx,wscale; 435334005bSBarry Smith double amax, epsilon = 1.e-8; /* assumes double precision */ 4454d73c9fSLois Curfman McInnes double dx_min = 1.e-16, dx_par = 1.e-1; 45bbb6d6a8SBarry Smith MPI_Comm comm; 460521c3abSLois Curfman McInnes int (*eval_fct)(SNES,Vec,Vec); 470521c3abSLois Curfman McInnes 48*3a40ed3dSBarry Smith PetscFunctionBegin; 49*3a40ed3dSBarry Smith if (snes->method_class == SNES_NONLINEAR_EQUATIONS) eval_fct = SNESComputeFunction; 50*3a40ed3dSBarry Smith else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) eval_fct = SNESComputeGradient; 51e3372554SBarry Smith else SETERRQ(1,0,"Invalid method class"); 5223242f5aSBarry Smith 53bbb6d6a8SBarry Smith PetscObjectGetComm((PetscObject)x1,&comm); 54b4fc646aSLois Curfman McInnes MatZeroEntries(*B); 55aa79bc6dSLois Curfman McInnes if (!snes->nvwork) { 56aa79bc6dSLois Curfman McInnes ierr = VecDuplicateVecs(x1,3,&snes->vwork); CHKERRQ(ierr); 57aa79bc6dSLois Curfman McInnes snes->nvwork = 3; 58aa79bc6dSLois Curfman McInnes PLogObjectParents(snes,3,snes->vwork); 59aa79bc6dSLois Curfman McInnes } 60aa79bc6dSLois Curfman McInnes j1 = snes->vwork[0]; j2 = snes->vwork[1]; x2 = snes->vwork[2]; 6123242f5aSBarry Smith 6278b31e54SBarry Smith ierr = VecGetSize(x1,&N); CHKERRQ(ierr); 6378b31e54SBarry Smith ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr); 6439e2f89bSBarry Smith VecGetArray(x1,&xx); 650521c3abSLois Curfman McInnes ierr = eval_fct(snes,x1,j1); CHKERRQ(ierr); 66c005e166SLois Curfman McInnes 67c005e166SLois Curfman McInnes /* Compute Jacobian approximation, 1 column at a time. 68c005e166SLois Curfman McInnes x1 = current iterate, j1 = F(x1) 69c005e166SLois Curfman McInnes x2 = perturbed iterate, j2 = F(x2) 70c005e166SLois Curfman McInnes */ 7139e2f89bSBarry Smith for ( i=0; i<N; i++ ) { 7278b31e54SBarry Smith ierr = VecCopy(x1,x2); CHKERRQ(ierr); 7323242f5aSBarry Smith if ( i>= start && i<end) { 7439e2f89bSBarry Smith dx = xx[i-start]; 75*3a40ed3dSBarry Smith #if !defined(USE_PETSC_COMPLEX) 7654d73c9fSLois Curfman McInnes if (dx < dx_min && dx >= 0.0) dx = dx_par; 7754d73c9fSLois Curfman McInnes else if (dx < 0.0 && dx > -dx_min) dx = -dx_par; 7819a167f6SBarry Smith #else 7954d73c9fSLois Curfman McInnes if (abs(dx) < dx_min && real(dx) >= 0.0) dx = dx_par; 8054d73c9fSLois Curfman McInnes else if (real(dx) < 0.0 && abs(dx) < dx_min) dx = -dx_par; 8119a167f6SBarry Smith #endif 8239e2f89bSBarry Smith dx *= epsilon; 8374f6f00dSLois Curfman McInnes wscale = 1.0/dx; 84dbb450caSBarry Smith VecSetValues(x2,1,&i,&dx,ADD_VALUES); 8511320018SBarry Smith } 86bbb6d6a8SBarry Smith else { 87bbb6d6a8SBarry Smith wscale = 0.0; 88bbb6d6a8SBarry Smith } 890521c3abSLois Curfman McInnes ierr = eval_fct(snes,x2,j2); CHKERRQ(ierr); 9078b31e54SBarry Smith ierr = VecAXPY(&mone,j1,j2); CHKERRQ(ierr); 91c005e166SLois Curfman McInnes /* Communicate scale to all processors */ 92*3a40ed3dSBarry Smith #if !defined(USE_PETSC_COMPLEX) 93bbb6d6a8SBarry Smith MPI_Allreduce(&wscale,&scale,1,MPI_DOUBLE,MPI_SUM,comm); 94bbb6d6a8SBarry Smith #else 95bbb6d6a8SBarry Smith MPI_Allreduce(&wscale,&scale,2,MPI_DOUBLE,MPI_SUM,comm); 96bbb6d6a8SBarry Smith #endif 9723242f5aSBarry Smith VecScale(&scale,j2); 9823242f5aSBarry Smith VecGetArray(j2,&y); 99cddf8d76SBarry Smith VecNorm(j2,NORM_INFINITY,&amax); 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 } 10523242f5aSBarry Smith VecRestoreArray(j2,&y); 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; 110*3a40ed3dSBarry 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 1180521c3abSLois Curfman McInnes Input Parameters: 1190521c3abSLois Curfman McInnes . x1 - compute Hessian at this point 1200521c3abSLois Curfman McInnes . ctx - application's gradient context, as set with SNESSetGradient() 1210521c3abSLois Curfman McInnes 1220521c3abSLois Curfman McInnes Output Parameters: 123b4fc646aSLois Curfman McInnes . J - Hessian matrix (not altered in this routine) 124b4fc646aSLois Curfman McInnes . B - newly computed Hessian matrix to use with preconditioner (generally the same as J) 125b4fc646aSLois Curfman McInnes . flag - flag indicating whether the matrix sparsity structure has changed 1260521c3abSLois Curfman McInnes 1270521c3abSLois Curfman McInnes Options Database Key: 1280521c3abSLois Curfman McInnes $ -snes_fd 1290521c3abSLois Curfman McInnes 1300521c3abSLois Curfman McInnes Notes: 1310521c3abSLois Curfman McInnes This routine is slow and expensive, and is not currently optimized 1320521c3abSLois Curfman McInnes to take advantage of sparsity in the problem. Although 1330521c3abSLois Curfman McInnes SNESDefaultComputeHessian() is not recommended for general use 1340521c3abSLois Curfman McInnes in large-scale applications, It can be useful in checking the 1350521c3abSLois Curfman McInnes correctness of a user-provided Hessian. 1360521c3abSLois Curfman McInnes 1370521c3abSLois Curfman McInnes .keywords: SNES, finite differences, Hessian 1380521c3abSLois Curfman McInnes 1399984d6fbSBarry Smith .seealso: SNESSetHessian() 1400521c3abSLois Curfman McInnes @*/ 1410521c3abSLois Curfman McInnes int SNESDefaultComputeHessian(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx) 1420521c3abSLois Curfman McInnes { 143*3a40ed3dSBarry Smith int ierr; 144*3a40ed3dSBarry Smith 145*3a40ed3dSBarry Smith PetscFunctionBegin; 146*3a40ed3dSBarry Smith ierr = SNESDefaultComputeJacobian(snes,x1,J,B,flag,ctx);CHKERRQ(ierr); 147*3a40ed3dSBarry Smith PetscFunctionReturn(0); 1480521c3abSLois Curfman McInnes } 149