111320018SBarry Smith 211320018SBarry Smith #ifndef lint 3*a86d99e1SLois Curfman McInnes static char vcid[] = "$Id: snesj.c,v 1.26 1996/02/05 19:32:38 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 45*a86d99e1SLois Curfman McInnes if (snes->method_class == SNES_NONLINEAR_EQUATIONS) 460521c3abSLois Curfman McInnes eval_fct = SNESComputeFunction; 47*a86d99e1SLois 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); 5378b31e54SBarry Smith ierr = VecDuplicate(x1,&j1); CHKERRQ(ierr); 5478b31e54SBarry Smith ierr = VecDuplicate(x1,&j2); CHKERRQ(ierr); 5578b31e54SBarry Smith ierr = VecDuplicate(x1,&x2); CHKERRQ(ierr); 56d370d78aSBarry Smith PLogObjectParent(snes,j1); PLogObjectParent(snes,j2); 57393d2d9aSLois Curfman McInnes PLogObjectParent(snes,x2); 5823242f5aSBarry Smith 5978b31e54SBarry Smith ierr = VecGetSize(x1,&N); CHKERRQ(ierr); 6078b31e54SBarry Smith ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr); 6139e2f89bSBarry Smith VecGetArray(x1,&xx); 620521c3abSLois Curfman McInnes ierr = eval_fct(snes,x1,j1); CHKERRQ(ierr); 6339e2f89bSBarry Smith for ( i=0; i<N; i++ ) { 6478b31e54SBarry Smith ierr = VecCopy(x1,x2); CHKERRQ(ierr); 6523242f5aSBarry Smith if ( i>= start && i<end) { 6639e2f89bSBarry Smith dx = xx[i-start]; 6719a167f6SBarry Smith #if !defined(PETSC_COMPLEX) 68336a5e98SBarry Smith if (dx < 1.e-16 && dx >= 0.0) dx = 1.e-1; 69336a5e98SBarry Smith else if (dx < 0.0 && dx > -1.e-16) dx = -1.e-1; 7019a167f6SBarry Smith #else 7119a167f6SBarry Smith if (abs(dx) < 1.e-16 && real(dx) >= 0.0) dx = 1.e-1; 72ec12a082SBarry Smith else if (real(dx) < 0.0 && abs(dx) < 1.e-16) dx = -1.e-1; 7319a167f6SBarry Smith #endif 7439e2f89bSBarry Smith dx *= epsilon; 75bbb6d6a8SBarry Smith wscale = -1.0/dx; 76dbb450caSBarry Smith VecSetValues(x2,1,&i,&dx,ADD_VALUES); 7711320018SBarry Smith } 78bbb6d6a8SBarry Smith else { 79bbb6d6a8SBarry Smith wscale = 0.0; 80bbb6d6a8SBarry Smith } 810521c3abSLois Curfman McInnes ierr = eval_fct(snes,x2,j2); CHKERRQ(ierr); 8278b31e54SBarry Smith ierr = VecAXPY(&mone,j1,j2); CHKERRQ(ierr); 83bbb6d6a8SBarry Smith /* communicate scale to all processors */ 84bbb6d6a8SBarry Smith #if !defined(PETSC_COMPLEX) 85bbb6d6a8SBarry Smith MPI_Allreduce(&wscale,&scale,1,MPI_DOUBLE,MPI_SUM,comm); 86bbb6d6a8SBarry Smith #else 87bbb6d6a8SBarry Smith MPI_Allreduce(&wscale,&scale,2,MPI_DOUBLE,MPI_SUM,comm); 88bbb6d6a8SBarry Smith #endif 895334005bSBarry Smith scale = -scale; 9023242f5aSBarry Smith VecScale(&scale,j2); 9123242f5aSBarry Smith VecGetArray(j2,&y); 92cddf8d76SBarry Smith VecNorm(j2,NORM_INFINITY,&amax); amax *= 1.e-14; 9323242f5aSBarry Smith for ( j=start; j<end; j++ ) { 94cddf8d76SBarry Smith if (PetscAbsScalar(y[j-start]) > amax) { 95dbb450caSBarry Smith ierr = MatSetValues(*J,1,&j,1,&i,y+j-start,INSERT_VALUES); CHKERRQ(ierr); 9623242f5aSBarry Smith } 9723242f5aSBarry Smith } 9823242f5aSBarry Smith VecRestoreArray(j2,&y); 9923242f5aSBarry Smith } 100ee50ffe9SBarry Smith MatAssemblyBegin(*J,FINAL_ASSEMBLY); 10123242f5aSBarry Smith VecDestroy(x2); VecDestroy(j1); VecDestroy(j2); 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