1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*c7afd0dbSLois Curfman McInnes static char vcid[] = "$Id: snesj.c,v 1.49 1998/04/20 18:09:47 curfman Exp curfman $"; 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 12fee21e36SBarry Smith Collective on SNES 13fee21e36SBarry Smith 14*c7afd0dbSLois Curfman McInnes Input Parameters: 15*c7afd0dbSLois Curfman McInnes + x1 - compute Jacobian at this point 16*c7afd0dbSLois Curfman McInnes - ctx - application's function context, as set with SNESSetFunction() 17*c7afd0dbSLois Curfman McInnes 18*c7afd0dbSLois Curfman McInnes Output Parameters: 19*c7afd0dbSLois Curfman McInnes + J - Jacobian matrix (not altered in this routine) 20*c7afd0dbSLois Curfman McInnes . B - newly computed Jacobian matrix to use with preconditioner (generally the same as J) 21*c7afd0dbSLois Curfman McInnes - flag - flag indicating whether the matrix sparsity structure has changed 22*c7afd0dbSLois Curfman McInnes 23ad960d00SLois Curfman McInnes Options Database Key: 24*c7afd0dbSLois Curfman McInnes . -snes_fd - Activates SNESDefaultComputeJacobian() 25ad960d00SLois Curfman McInnes 265f3c43d9SLois Curfman McInnes Notes: 275f3c43d9SLois Curfman McInnes This routine is slow and expensive, and is not currently optimized 285f3c43d9SLois Curfman McInnes to take advantage of sparsity in the problem. Although 295f3c43d9SLois Curfman McInnes SNESDefaultComputeJacobian() is not recommended for general use 305f3c43d9SLois Curfman McInnes in large-scale applications, It can be useful in checking the 315f3c43d9SLois Curfman McInnes correctness of a user-provided Jacobian. 3211320018SBarry Smith 33b4fc646aSLois Curfman McInnes An alternative routine that uses coloring to explot matrix sparsity is 34b4fc646aSLois Curfman McInnes SNESDefaultComputeJacobianWithColoring(). 35b4fc646aSLois Curfman McInnes 365f3c43d9SLois Curfman McInnes .keywords: SNES, finite differences, Jacobian 375f3c43d9SLois Curfman McInnes 38b4fc646aSLois Curfman McInnes .seealso: SNESSetJacobian(), SNESDefaultComputeJacobianWithColoring() 3911320018SBarry Smith @*/ 40*c7afd0dbSLois Curfman McInnes int SNESDefaultComputeJacobian(SNES snes,Vec x1,Mat *J,Mat *B, 41*c7afd0dbSLois Curfman McInnes MatStructure *flag,void *ctx) 4211320018SBarry Smith { 4388c956adSLois Curfman McInnes Vec j1a,j2a,x2; 4423242f5aSBarry Smith int i,ierr,N,start,end,j; 45bbb6d6a8SBarry Smith Scalar dx, mone = -1.0,*y,scale,*xx,wscale; 465334005bSBarry Smith double amax, epsilon = 1.e-8; /* assumes double precision */ 4754d73c9fSLois Curfman McInnes double dx_min = 1.e-16, dx_par = 1.e-1; 48bbb6d6a8SBarry Smith MPI_Comm comm; 490521c3abSLois Curfman McInnes int (*eval_fct)(SNES,Vec,Vec); 500521c3abSLois Curfman McInnes 513a40ed3dSBarry Smith PetscFunctionBegin; 523a40ed3dSBarry Smith if (snes->method_class == SNES_NONLINEAR_EQUATIONS) eval_fct = SNESComputeFunction; 533a40ed3dSBarry Smith else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) eval_fct = SNESComputeGradient; 54a8c6a408SBarry Smith else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Invalid method class"); 5523242f5aSBarry Smith 56bbb6d6a8SBarry Smith PetscObjectGetComm((PetscObject)x1,&comm); 57b4fc646aSLois Curfman McInnes MatZeroEntries(*B); 58aa79bc6dSLois Curfman McInnes if (!snes->nvwork) { 59aa79bc6dSLois Curfman McInnes ierr = VecDuplicateVecs(x1,3,&snes->vwork); CHKERRQ(ierr); 60aa79bc6dSLois Curfman McInnes snes->nvwork = 3; 61aa79bc6dSLois Curfman McInnes PLogObjectParents(snes,3,snes->vwork); 62aa79bc6dSLois Curfman McInnes } 6388c956adSLois Curfman McInnes j1a = snes->vwork[0]; j2a = snes->vwork[1]; x2 = snes->vwork[2]; 6423242f5aSBarry Smith 6578b31e54SBarry Smith ierr = VecGetSize(x1,&N); CHKERRQ(ierr); 6678b31e54SBarry Smith ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr); 6739e2f89bSBarry Smith VecGetArray(x1,&xx); 6888c956adSLois Curfman McInnes ierr = eval_fct(snes,x1,j1a); CHKERRQ(ierr); 69c005e166SLois Curfman McInnes 70c005e166SLois Curfman McInnes /* Compute Jacobian approximation, 1 column at a time. 7188c956adSLois Curfman McInnes x1 = current iterate, j1a = F(x1) 7288c956adSLois Curfman McInnes x2 = perturbed iterate, j2a = F(x2) 73c005e166SLois Curfman McInnes */ 7439e2f89bSBarry Smith for ( i=0; i<N; i++ ) { 7578b31e54SBarry Smith ierr = VecCopy(x1,x2); CHKERRQ(ierr); 7623242f5aSBarry Smith if ( i>= start && i<end) { 7739e2f89bSBarry Smith dx = xx[i-start]; 783a40ed3dSBarry Smith #if !defined(USE_PETSC_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 8254d73c9fSLois Curfman McInnes if (abs(dx) < dx_min && real(dx) >= 0.0) dx = dx_par; 8354d73c9fSLois Curfman McInnes else if (real(dx) < 0.0 && abs(dx) < dx_min) dx = -dx_par; 8419a167f6SBarry Smith #endif 8539e2f89bSBarry Smith dx *= epsilon; 8674f6f00dSLois Curfman McInnes wscale = 1.0/dx; 87dbb450caSBarry Smith VecSetValues(x2,1,&i,&dx,ADD_VALUES); 8811320018SBarry Smith } 89bbb6d6a8SBarry Smith else { 90bbb6d6a8SBarry Smith wscale = 0.0; 91bbb6d6a8SBarry Smith } 9288c956adSLois Curfman McInnes ierr = eval_fct(snes,x2,j2a); CHKERRQ(ierr); 9388c956adSLois Curfman McInnes ierr = VecAXPY(&mone,j1a,j2a); CHKERRQ(ierr); 94c005e166SLois Curfman McInnes /* Communicate scale to all processors */ 953a40ed3dSBarry Smith #if !defined(USE_PETSC_COMPLEX) 96ca161407SBarry Smith ierr = MPI_Allreduce(&wscale,&scale,1,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr); 97bbb6d6a8SBarry Smith #else 98ca161407SBarry Smith ierr = MPI_Allreduce(&wscale,&scale,2,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr); 99bbb6d6a8SBarry Smith #endif 10088c956adSLois Curfman McInnes VecScale(&scale,j2a); 10188c956adSLois Curfman McInnes VecGetArray(j2a,&y); 10288c956adSLois Curfman McInnes VecNorm(j2a,NORM_INFINITY,&amax); amax *= 1.e-14; 10323242f5aSBarry Smith for ( j=start; j<end; j++ ) { 104cddf8d76SBarry Smith if (PetscAbsScalar(y[j-start]) > amax) { 105b4fc646aSLois Curfman McInnes ierr = MatSetValues(*B,1,&j,1,&i,y+j-start,INSERT_VALUES); CHKERRQ(ierr); 10623242f5aSBarry Smith } 10723242f5aSBarry Smith } 10888c956adSLois Curfman McInnes VecRestoreArray(j2a,&y); 10923242f5aSBarry Smith } 110b4fc646aSLois Curfman McInnes ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 111b4fc646aSLois Curfman McInnes ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 112595b82d2SBarry Smith *flag = DIFFERENT_NONZERO_PATTERN; 1133a40ed3dSBarry Smith PetscFunctionReturn(0); 11411320018SBarry Smith } 11511320018SBarry Smith 1165615d1e5SSatish Balay #undef __FUNC__ 1175615d1e5SSatish Balay #define __FUNC__ "SNESDefaultComputeHessian" 1180521c3abSLois Curfman McInnes /*@C 11984a7bf8dSLois Curfman McInnes SNESDefaultComputeHessian - Computes the Hessian using finite differences. 1200521c3abSLois Curfman McInnes 121fee21e36SBarry Smith Collective on SNES 122fee21e36SBarry Smith 123*c7afd0dbSLois Curfman McInnes Input Parameters: 124*c7afd0dbSLois Curfman McInnes + x1 - compute Hessian at this point 125*c7afd0dbSLois Curfman McInnes - ctx - application's gradient context, as set with SNESSetGradient() 126*c7afd0dbSLois Curfman McInnes 127*c7afd0dbSLois Curfman McInnes Output Parameters: 128*c7afd0dbSLois Curfman McInnes + J - Hessian matrix (not altered in this routine) 129*c7afd0dbSLois Curfman McInnes . B - newly computed Hessian matrix to use with preconditioner (generally the same as J) 130*c7afd0dbSLois Curfman McInnes - flag - flag indicating whether the matrix sparsity structure has changed 131*c7afd0dbSLois Curfman McInnes 1320521c3abSLois Curfman McInnes Options Database Key: 133*c7afd0dbSLois Curfman McInnes $ -snes_fd - Activates SNESDefaultComputeHessian() 1340521c3abSLois Curfman McInnes 1350521c3abSLois Curfman McInnes Notes: 1360521c3abSLois Curfman McInnes This routine is slow and expensive, and is not currently optimized 1370521c3abSLois Curfman McInnes to take advantage of sparsity in the problem. Although 1380521c3abSLois Curfman McInnes SNESDefaultComputeHessian() is not recommended for general use 1390521c3abSLois Curfman McInnes in large-scale applications, It can be useful in checking the 1400521c3abSLois Curfman McInnes correctness of a user-provided Hessian. 1410521c3abSLois Curfman McInnes 1420521c3abSLois Curfman McInnes .keywords: SNES, finite differences, Hessian 1430521c3abSLois Curfman McInnes 1449984d6fbSBarry Smith .seealso: SNESSetHessian() 1450521c3abSLois Curfman McInnes @*/ 146*c7afd0dbSLois Curfman McInnes int SNESDefaultComputeHessian(SNES snes,Vec x1,Mat *J,Mat *B, 147*c7afd0dbSLois Curfman McInnes MatStructure *flag,void *ctx) 1480521c3abSLois Curfman McInnes { 1493a40ed3dSBarry Smith int ierr; 1503a40ed3dSBarry Smith 1513a40ed3dSBarry Smith PetscFunctionBegin; 1523a40ed3dSBarry Smith ierr = SNESDefaultComputeJacobian(snes,x1,J,B,flag,ctx);CHKERRQ(ierr); 1533a40ed3dSBarry Smith PetscFunctionReturn(0); 1540521c3abSLois Curfman McInnes } 155