111320018SBarry Smith 211320018SBarry Smith #ifndef lint 3*aa79bc6dSLois Curfman McInnes static char vcid[] = "$Id: snesj.c,v 1.27 1996/02/15 02:06:09 curfman Exp curfman $"; 411320018SBarry Smith #endif 511320018SBarry Smith 6b1f0a012SBarry Smith #include "draw.h" /*I "draw.h" I*/ 70521c3abSLois Curfman McInnes #include "snesimpl.h" /*I "snes.h" I*/ 811320018SBarry Smith 94b828684SBarry Smith /*@C 105f3c43d9SLois Curfman McInnes SNESDefaultComputeJacobian - Computes the Jacobian using finite 115f3c43d9SLois Curfman McInnes differences. 1211320018SBarry Smith 1311320018SBarry Smith Input Parameters: 14da88328cSLois Curfman McInnes . x1 - compute Jacobian at this point 15da88328cSLois Curfman McInnes . ctx - application's function context, as set with SNESSetFunction() 1611320018SBarry Smith 1711320018SBarry Smith Output Parameters: 1823242f5aSBarry Smith . J - Jacobian 1923242f5aSBarry Smith . B - preconditioner, same as Jacobian 20da88328cSLois Curfman McInnes . flag - matrix flag 2111320018SBarry Smith 22ad960d00SLois Curfman McInnes Options Database Key: 23ad960d00SLois Curfman McInnes $ -snes_fd 24ad960d00SLois Curfman McInnes 255f3c43d9SLois Curfman McInnes Notes: 265f3c43d9SLois Curfman McInnes This routine is slow and expensive, and is not currently optimized 275f3c43d9SLois Curfman McInnes to take advantage of sparsity in the problem. Although 285f3c43d9SLois Curfman McInnes SNESDefaultComputeJacobian() is not recommended for general use 295f3c43d9SLois Curfman McInnes in large-scale applications, It can be useful in checking the 305f3c43d9SLois Curfman McInnes correctness of a user-provided Jacobian. 3111320018SBarry Smith 325f3c43d9SLois Curfman McInnes .keywords: SNES, finite differences, Jacobian 335f3c43d9SLois Curfman McInnes 345f3c43d9SLois Curfman McInnes .seealso: SNESSetJacobian(), SNESTestJacobian() 3511320018SBarry Smith @*/ 365334005bSBarry Smith int SNESDefaultComputeJacobian(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx) 3711320018SBarry Smith { 3823242f5aSBarry Smith Vec j1,j2,x2; 3923242f5aSBarry Smith int i,ierr,N,start,end,j; 40bbb6d6a8SBarry Smith Scalar dx, mone = -1.0,*y,scale,*xx,wscale; 415334005bSBarry Smith double amax, epsilon = 1.e-8; /* assumes double precision */ 42bbb6d6a8SBarry Smith MPI_Comm comm; 430521c3abSLois Curfman McInnes int (*eval_fct)(SNES,Vec,Vec); 440521c3abSLois Curfman McInnes 45a86d99e1SLois Curfman McInnes if (snes->method_class == SNES_NONLINEAR_EQUATIONS) 460521c3abSLois Curfman McInnes eval_fct = SNESComputeFunction; 47a86d99e1SLois Curfman McInnes else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) 480521c3abSLois Curfman McInnes eval_fct = SNESComputeGradient; 490521c3abSLois Curfman McInnes else SETERRQ(1,"SNESDefaultComputeJacobian: Invalid method class"); 5023242f5aSBarry Smith 51bbb6d6a8SBarry Smith PetscObjectGetComm((PetscObject)x1,&comm); 52336a5e98SBarry Smith MatZeroEntries(*J); 53*aa79bc6dSLois Curfman McInnes if (!snes->nvwork) { 54*aa79bc6dSLois Curfman McInnes ierr = VecDuplicateVecs(x1,3,&snes->vwork); CHKERRQ(ierr); 55*aa79bc6dSLois Curfman McInnes snes->nvwork = 3; 56*aa79bc6dSLois Curfman McInnes PLogObjectParents(snes,3,snes->vwork); 57*aa79bc6dSLois Curfman McInnes } 58*aa79bc6dSLois Curfman McInnes j1 = snes->vwork[0]; j2 = snes->vwork[1]; x2 = snes->vwork[2]; 5923242f5aSBarry Smith 6078b31e54SBarry Smith ierr = VecGetSize(x1,&N); CHKERRQ(ierr); 6178b31e54SBarry Smith ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr); 6239e2f89bSBarry Smith VecGetArray(x1,&xx); 630521c3abSLois Curfman McInnes ierr = eval_fct(snes,x1,j1); CHKERRQ(ierr); 6439e2f89bSBarry Smith for ( i=0; i<N; i++ ) { 6578b31e54SBarry Smith ierr = VecCopy(x1,x2); CHKERRQ(ierr); 6623242f5aSBarry Smith if ( i>= start && i<end) { 6739e2f89bSBarry Smith dx = xx[i-start]; 6819a167f6SBarry Smith #if !defined(PETSC_COMPLEX) 69336a5e98SBarry Smith if (dx < 1.e-16 && dx >= 0.0) dx = 1.e-1; 70336a5e98SBarry Smith else if (dx < 0.0 && dx > -1.e-16) dx = -1.e-1; 7119a167f6SBarry Smith #else 7219a167f6SBarry Smith if (abs(dx) < 1.e-16 && real(dx) >= 0.0) dx = 1.e-1; 73ec12a082SBarry Smith else if (real(dx) < 0.0 && abs(dx) < 1.e-16) dx = -1.e-1; 7419a167f6SBarry Smith #endif 7539e2f89bSBarry Smith dx *= epsilon; 76bbb6d6a8SBarry Smith wscale = -1.0/dx; 77dbb450caSBarry Smith VecSetValues(x2,1,&i,&dx,ADD_VALUES); 7811320018SBarry Smith } 79bbb6d6a8SBarry Smith else { 80bbb6d6a8SBarry Smith wscale = 0.0; 81bbb6d6a8SBarry Smith } 820521c3abSLois Curfman McInnes ierr = eval_fct(snes,x2,j2); CHKERRQ(ierr); 8378b31e54SBarry Smith ierr = VecAXPY(&mone,j1,j2); CHKERRQ(ierr); 84bbb6d6a8SBarry Smith /* communicate scale to all processors */ 85bbb6d6a8SBarry Smith #if !defined(PETSC_COMPLEX) 86bbb6d6a8SBarry Smith MPI_Allreduce(&wscale,&scale,1,MPI_DOUBLE,MPI_SUM,comm); 87bbb6d6a8SBarry Smith #else 88bbb6d6a8SBarry Smith MPI_Allreduce(&wscale,&scale,2,MPI_DOUBLE,MPI_SUM,comm); 89bbb6d6a8SBarry Smith #endif 905334005bSBarry Smith scale = -scale; 9123242f5aSBarry Smith VecScale(&scale,j2); 9223242f5aSBarry Smith VecGetArray(j2,&y); 93cddf8d76SBarry Smith VecNorm(j2,NORM_INFINITY,&amax); amax *= 1.e-14; 9423242f5aSBarry Smith for ( j=start; j<end; j++ ) { 95cddf8d76SBarry Smith if (PetscAbsScalar(y[j-start]) > amax) { 96dbb450caSBarry Smith ierr = MatSetValues(*J,1,&j,1,&i,y+j-start,INSERT_VALUES); CHKERRQ(ierr); 9723242f5aSBarry Smith } 9823242f5aSBarry Smith } 9923242f5aSBarry Smith VecRestoreArray(j2,&y); 10023242f5aSBarry Smith } 101ee50ffe9SBarry Smith MatAssemblyBegin(*J,FINAL_ASSEMBLY); 102ee50ffe9SBarry Smith MatAssemblyEnd(*J,FINAL_ASSEMBLY); 10311320018SBarry Smith return 0; 10411320018SBarry Smith } 10511320018SBarry Smith 1060521c3abSLois Curfman McInnes /*@C 1070521c3abSLois Curfman McInnes SNESDefaultComputeHessian - Computes the Hessian using finite 1080521c3abSLois Curfman McInnes differences. 1090521c3abSLois Curfman McInnes 1100521c3abSLois Curfman McInnes Input Parameters: 1110521c3abSLois Curfman McInnes . x1 - compute Hessian at this point 1120521c3abSLois Curfman McInnes . ctx - application's gradient context, as set with SNESSetGradient() 1130521c3abSLois Curfman McInnes 1140521c3abSLois Curfman McInnes Output Parameters: 1150521c3abSLois Curfman McInnes . J - Hessian 1160521c3abSLois Curfman McInnes . B - preconditioner, same as Hessian 1170521c3abSLois Curfman McInnes . flag - matrix flag 1180521c3abSLois Curfman McInnes 1190521c3abSLois Curfman McInnes Options Database Key: 1200521c3abSLois Curfman McInnes $ -snes_fd 1210521c3abSLois Curfman McInnes 1220521c3abSLois Curfman McInnes Notes: 1230521c3abSLois Curfman McInnes This routine is slow and expensive, and is not currently optimized 1240521c3abSLois Curfman McInnes to take advantage of sparsity in the problem. Although 1250521c3abSLois Curfman McInnes SNESDefaultComputeHessian() is not recommended for general use 1260521c3abSLois Curfman McInnes in large-scale applications, It can be useful in checking the 1270521c3abSLois Curfman McInnes correctness of a user-provided Hessian. 1280521c3abSLois Curfman McInnes 1290521c3abSLois Curfman McInnes .keywords: SNES, finite differences, Hessian 1300521c3abSLois Curfman McInnes 1310521c3abSLois Curfman McInnes .seealso: SNESSetHessian(), SNESTestHessian() 1320521c3abSLois Curfman McInnes @*/ 1330521c3abSLois Curfman McInnes int SNESDefaultComputeHessian(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx) 1340521c3abSLois Curfman McInnes { 1350521c3abSLois Curfman McInnes return SNESDefaultComputeJacobian(snes,x1,J,B,flag,ctx); 1360521c3abSLois Curfman McInnes } 137