111320018SBarry Smith 211320018SBarry Smith #ifndef lint 3*54d73c9fSLois Curfman McInnes static char vcid[] = "$Id: snesj.c,v 1.34 1996/08/27 21:07:44 bsmith Exp curfman $"; 411320018SBarry Smith #endif 511320018SBarry Smith 6b1f0a012SBarry Smith #include "draw.h" /*I "draw.h" I*/ 770f55243SBarry Smith #include "src/snes/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 */ 42*54d73c9fSLois Curfman McInnes double dx_min = 1.e-16, dx_par = 1.e-1; 43bbb6d6a8SBarry Smith MPI_Comm comm; 440521c3abSLois Curfman McInnes int (*eval_fct)(SNES,Vec,Vec); 450521c3abSLois Curfman McInnes 46a86d99e1SLois Curfman McInnes if (snes->method_class == SNES_NONLINEAR_EQUATIONS) 470521c3abSLois Curfman McInnes eval_fct = SNESComputeFunction; 48a86d99e1SLois Curfman McInnes else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) 490521c3abSLois Curfman McInnes eval_fct = SNESComputeGradient; 500521c3abSLois Curfman McInnes else SETERRQ(1,"SNESDefaultComputeJacobian: Invalid method class"); 5123242f5aSBarry Smith 52bbb6d6a8SBarry Smith PetscObjectGetComm((PetscObject)x1,&comm); 53336a5e98SBarry Smith MatZeroEntries(*J); 54aa79bc6dSLois Curfman McInnes if (!snes->nvwork) { 55aa79bc6dSLois Curfman McInnes ierr = VecDuplicateVecs(x1,3,&snes->vwork); CHKERRQ(ierr); 56aa79bc6dSLois Curfman McInnes snes->nvwork = 3; 57aa79bc6dSLois Curfman McInnes PLogObjectParents(snes,3,snes->vwork); 58aa79bc6dSLois Curfman McInnes } 59aa79bc6dSLois Curfman McInnes j1 = snes->vwork[0]; j2 = snes->vwork[1]; x2 = snes->vwork[2]; 6023242f5aSBarry Smith 6178b31e54SBarry Smith ierr = VecGetSize(x1,&N); CHKERRQ(ierr); 6278b31e54SBarry Smith ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr); 6339e2f89bSBarry Smith VecGetArray(x1,&xx); 640521c3abSLois Curfman McInnes ierr = eval_fct(snes,x1,j1); CHKERRQ(ierr); 65c005e166SLois Curfman McInnes 66c005e166SLois Curfman McInnes /* Compute Jacobian approximation, 1 column at a time. 67c005e166SLois Curfman McInnes x1 = current iterate, j1 = F(x1) 68c005e166SLois Curfman McInnes x2 = perturbed iterate, j2 = F(x2) 69c005e166SLois Curfman McInnes */ 7039e2f89bSBarry Smith for ( i=0; i<N; i++ ) { 7178b31e54SBarry Smith ierr = VecCopy(x1,x2); CHKERRQ(ierr); 7223242f5aSBarry Smith if ( i>= start && i<end) { 7339e2f89bSBarry Smith dx = xx[i-start]; 7419a167f6SBarry Smith #if !defined(PETSC_COMPLEX) 75*54d73c9fSLois Curfman McInnes if (dx < dx_min && dx >= 0.0) dx = dx_par; 76*54d73c9fSLois Curfman McInnes else if (dx < 0.0 && dx > -dx_min) dx = -dx_par; 7719a167f6SBarry Smith #else 78*54d73c9fSLois Curfman McInnes if (abs(dx) < dx_min && real(dx) >= 0.0) dx = dx_par; 79*54d73c9fSLois Curfman McInnes else if (real(dx) < 0.0 && abs(dx) < dx_min) dx = -dx_par; 8019a167f6SBarry Smith #endif 8139e2f89bSBarry Smith dx *= epsilon; 8274f6f00dSLois Curfman McInnes wscale = 1.0/dx; 83dbb450caSBarry Smith VecSetValues(x2,1,&i,&dx,ADD_VALUES); 8411320018SBarry Smith } 85bbb6d6a8SBarry Smith else { 86bbb6d6a8SBarry Smith wscale = 0.0; 87bbb6d6a8SBarry Smith } 880521c3abSLois Curfman McInnes ierr = eval_fct(snes,x2,j2); CHKERRQ(ierr); 8978b31e54SBarry Smith ierr = VecAXPY(&mone,j1,j2); CHKERRQ(ierr); 90c005e166SLois Curfman McInnes /* Communicate scale to all processors */ 91bbb6d6a8SBarry Smith #if !defined(PETSC_COMPLEX) 92bbb6d6a8SBarry Smith MPI_Allreduce(&wscale,&scale,1,MPI_DOUBLE,MPI_SUM,comm); 93bbb6d6a8SBarry Smith #else 94bbb6d6a8SBarry Smith MPI_Allreduce(&wscale,&scale,2,MPI_DOUBLE,MPI_SUM,comm); 95bbb6d6a8SBarry Smith #endif 9623242f5aSBarry Smith VecScale(&scale,j2); 9723242f5aSBarry Smith VecGetArray(j2,&y); 98cddf8d76SBarry Smith VecNorm(j2,NORM_INFINITY,&amax); amax *= 1.e-14; 9923242f5aSBarry Smith for ( j=start; j<end; j++ ) { 100cddf8d76SBarry Smith if (PetscAbsScalar(y[j-start]) > amax) { 101dbb450caSBarry Smith ierr = MatSetValues(*J,1,&j,1,&i,y+j-start,INSERT_VALUES); CHKERRQ(ierr); 10223242f5aSBarry Smith } 10323242f5aSBarry Smith } 10423242f5aSBarry Smith VecRestoreArray(j2,&y); 10523242f5aSBarry Smith } 1060989deb7SLois Curfman McInnes ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 1070989deb7SLois Curfman McInnes ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 108595b82d2SBarry Smith *flag = DIFFERENT_NONZERO_PATTERN; 10911320018SBarry Smith return 0; 11011320018SBarry Smith } 11111320018SBarry Smith 1120521c3abSLois Curfman McInnes /*@C 1130521c3abSLois Curfman McInnes SNESDefaultComputeHessian - Computes the Hessian using finite 1140521c3abSLois Curfman McInnes differences. 1150521c3abSLois Curfman McInnes 1160521c3abSLois Curfman McInnes Input Parameters: 1170521c3abSLois Curfman McInnes . x1 - compute Hessian at this point 1180521c3abSLois Curfman McInnes . ctx - application's gradient context, as set with SNESSetGradient() 1190521c3abSLois Curfman McInnes 1200521c3abSLois Curfman McInnes Output Parameters: 1210521c3abSLois Curfman McInnes . J - Hessian 1220521c3abSLois Curfman McInnes . B - preconditioner, same as Hessian 1230521c3abSLois Curfman McInnes . flag - matrix flag 1240521c3abSLois Curfman McInnes 1250521c3abSLois Curfman McInnes Options Database Key: 1260521c3abSLois Curfman McInnes $ -snes_fd 1270521c3abSLois Curfman McInnes 1280521c3abSLois Curfman McInnes Notes: 1290521c3abSLois Curfman McInnes This routine is slow and expensive, and is not currently optimized 1300521c3abSLois Curfman McInnes to take advantage of sparsity in the problem. Although 1310521c3abSLois Curfman McInnes SNESDefaultComputeHessian() is not recommended for general use 1320521c3abSLois Curfman McInnes in large-scale applications, It can be useful in checking the 1330521c3abSLois Curfman McInnes correctness of a user-provided Hessian. 1340521c3abSLois Curfman McInnes 1350521c3abSLois Curfman McInnes .keywords: SNES, finite differences, Hessian 1360521c3abSLois Curfman McInnes 1370521c3abSLois Curfman McInnes .seealso: SNESSetHessian(), SNESTestHessian() 1380521c3abSLois Curfman McInnes @*/ 1390521c3abSLois Curfman McInnes int SNESDefaultComputeHessian(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx) 1400521c3abSLois Curfman McInnes { 1410521c3abSLois Curfman McInnes return SNESDefaultComputeJacobian(snes,x1,J,B,flag,ctx); 1420521c3abSLois Curfman McInnes } 143