1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*36851e7fSLois Curfman McInnes static char vcid[] = "$Id: snesj.c,v 1.54 1998/12/09 16:06:57 balay 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 14c7afd0dbSLois Curfman McInnes Input Parameters: 15c7afd0dbSLois Curfman McInnes + x1 - compute Jacobian at this point 16c7afd0dbSLois Curfman McInnes - ctx - application's function context, as set with SNESSetFunction() 17c7afd0dbSLois Curfman McInnes 18c7afd0dbSLois Curfman McInnes Output Parameters: 19c7afd0dbSLois Curfman McInnes + J - Jacobian matrix (not altered in this routine) 20c7afd0dbSLois Curfman McInnes . B - newly computed Jacobian matrix to use with preconditioner (generally the same as J) 21c7afd0dbSLois Curfman McInnes - flag - flag indicating whether the matrix sparsity structure has changed 22c7afd0dbSLois Curfman McInnes 23ad960d00SLois Curfman McInnes Options Database Key: 24c7afd0dbSLois 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 36*36851e7fSLois Curfman McInnes Level: intermediate 37*36851e7fSLois Curfman McInnes 385f3c43d9SLois Curfman McInnes .keywords: SNES, finite differences, Jacobian 395f3c43d9SLois Curfman McInnes 40b4fc646aSLois Curfman McInnes .seealso: SNESSetJacobian(), SNESDefaultComputeJacobianWithColoring() 4111320018SBarry Smith @*/ 42c7afd0dbSLois Curfman McInnes int SNESDefaultComputeJacobian(SNES snes,Vec x1,Mat *J,Mat *B, 43c7afd0dbSLois Curfman McInnes MatStructure *flag,void *ctx) 4411320018SBarry Smith { 4588c956adSLois Curfman McInnes Vec j1a,j2a,x2; 4623242f5aSBarry Smith int i,ierr,N,start,end,j; 47bbb6d6a8SBarry Smith Scalar dx, mone = -1.0,*y,scale,*xx,wscale; 485334005bSBarry Smith double amax, epsilon = 1.e-8; /* assumes double precision */ 4954d73c9fSLois Curfman McInnes double dx_min = 1.e-16, dx_par = 1.e-1; 50bbb6d6a8SBarry Smith MPI_Comm comm; 51a305c92eSSatish Balay int (*eval_fct)(SNES,Vec,Vec)=0; 520521c3abSLois Curfman McInnes 533a40ed3dSBarry Smith PetscFunctionBegin; 543a40ed3dSBarry Smith if (snes->method_class == SNES_NONLINEAR_EQUATIONS) eval_fct = SNESComputeFunction; 553a40ed3dSBarry Smith else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) eval_fct = SNESComputeGradient; 56a8c6a408SBarry Smith else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Invalid method class"); 5723242f5aSBarry Smith 58bbb6d6a8SBarry Smith PetscObjectGetComm((PetscObject)x1,&comm); 59b4fc646aSLois Curfman McInnes MatZeroEntries(*B); 60aa79bc6dSLois Curfman McInnes if (!snes->nvwork) { 61aa79bc6dSLois Curfman McInnes ierr = VecDuplicateVecs(x1,3,&snes->vwork); CHKERRQ(ierr); 62aa79bc6dSLois Curfman McInnes snes->nvwork = 3; 63aa79bc6dSLois Curfman McInnes PLogObjectParents(snes,3,snes->vwork); 64aa79bc6dSLois Curfman McInnes } 6588c956adSLois Curfman McInnes j1a = snes->vwork[0]; j2a = snes->vwork[1]; x2 = snes->vwork[2]; 6623242f5aSBarry Smith 6778b31e54SBarry Smith ierr = VecGetSize(x1,&N); CHKERRQ(ierr); 6878b31e54SBarry Smith ierr = VecGetOwnershipRange(x1,&start,&end); CHKERRQ(ierr); 6988c956adSLois Curfman McInnes ierr = eval_fct(snes,x1,j1a); CHKERRQ(ierr); 70c005e166SLois Curfman McInnes 71c005e166SLois Curfman McInnes /* Compute Jacobian approximation, 1 column at a time. 7288c956adSLois Curfman McInnes x1 = current iterate, j1a = F(x1) 7388c956adSLois Curfman McInnes x2 = perturbed iterate, j2a = F(x2) 74c005e166SLois Curfman McInnes */ 7539e2f89bSBarry Smith for ( i=0; i<N; i++ ) { 7678b31e54SBarry Smith ierr = VecCopy(x1,x2); CHKERRQ(ierr); 7723242f5aSBarry Smith if ( i>= start && i<end) { 782e8a6d31SBarry Smith ierr = VecGetArray(x1,&xx);CHKERRQ(ierr); 7939e2f89bSBarry Smith dx = xx[i-start]; 802e8a6d31SBarry Smith ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 813a40ed3dSBarry Smith #if !defined(USE_PETSC_COMPLEX) 8254d73c9fSLois Curfman McInnes if (dx < dx_min && dx >= 0.0) dx = dx_par; 8354d73c9fSLois Curfman McInnes else if (dx < 0.0 && dx > -dx_min) dx = -dx_par; 8419a167f6SBarry Smith #else 85e20fef11SSatish Balay if (PetscAbsScalar(dx) < dx_min && PetscReal(dx) >= 0.0) dx = dx_par; 86e20fef11SSatish Balay else if (PetscReal(dx) < 0.0 && PetscAbsScalar(dx) < dx_min) dx = -dx_par; 8719a167f6SBarry Smith #endif 8839e2f89bSBarry Smith dx *= epsilon; 8974f6f00dSLois Curfman McInnes wscale = 1.0/dx; 902e8a6d31SBarry Smith ierr = VecSetValues(x2,1,&i,&dx,ADD_VALUES); CHKERRQ(ierr); 916c783aadSBarry Smith } else { 92bbb6d6a8SBarry Smith wscale = 0.0; 93bbb6d6a8SBarry Smith } 9488c956adSLois Curfman McInnes ierr = eval_fct(snes,x2,j2a); CHKERRQ(ierr); 9588c956adSLois Curfman McInnes ierr = VecAXPY(&mone,j1a,j2a); CHKERRQ(ierr); 96c005e166SLois Curfman McInnes /* Communicate scale to all processors */ 973a40ed3dSBarry Smith #if !defined(USE_PETSC_COMPLEX) 98ca161407SBarry Smith ierr = MPI_Allreduce(&wscale,&scale,1,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr); 99bbb6d6a8SBarry Smith #else 100ca161407SBarry Smith ierr = MPI_Allreduce(&wscale,&scale,2,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr); 101bbb6d6a8SBarry Smith #endif 1022e8a6d31SBarry Smith ierr = VecScale(&scale,j2a);CHKERRQ(ierr); 1032e8a6d31SBarry Smith ierr = VecGetArray(j2a,&y);CHKERRQ(ierr); 1042e8a6d31SBarry Smith ierr = VecNorm(j2a,NORM_INFINITY,&amax);CHKERRQ(ierr); amax *= 1.e-14; 10523242f5aSBarry Smith for ( j=start; j<end; j++ ) { 106cddf8d76SBarry Smith if (PetscAbsScalar(y[j-start]) > amax) { 107b4fc646aSLois Curfman McInnes ierr = MatSetValues(*B,1,&j,1,&i,y+j-start,INSERT_VALUES); CHKERRQ(ierr); 10823242f5aSBarry Smith } 10923242f5aSBarry Smith } 1102e8a6d31SBarry Smith ierr = VecRestoreArray(j2a,&y);CHKERRQ(ierr); 11123242f5aSBarry Smith } 112b4fc646aSLois Curfman McInnes ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 113b4fc646aSLois Curfman McInnes ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 114595b82d2SBarry Smith *flag = DIFFERENT_NONZERO_PATTERN; 1153a40ed3dSBarry Smith PetscFunctionReturn(0); 11611320018SBarry Smith } 11711320018SBarry Smith 1185615d1e5SSatish Balay #undef __FUNC__ 1195615d1e5SSatish Balay #define __FUNC__ "SNESDefaultComputeHessian" 1200521c3abSLois Curfman McInnes /*@C 12184a7bf8dSLois Curfman McInnes SNESDefaultComputeHessian - Computes the Hessian using finite differences. 1220521c3abSLois Curfman McInnes 123fee21e36SBarry Smith Collective on SNES 124fee21e36SBarry Smith 125c7afd0dbSLois Curfman McInnes Input Parameters: 126c7afd0dbSLois Curfman McInnes + x1 - compute Hessian at this point 127c7afd0dbSLois Curfman McInnes - ctx - application's gradient context, as set with SNESSetGradient() 128c7afd0dbSLois Curfman McInnes 129c7afd0dbSLois Curfman McInnes Output Parameters: 130c7afd0dbSLois Curfman McInnes + J - Hessian matrix (not altered in this routine) 131c7afd0dbSLois Curfman McInnes . B - newly computed Hessian matrix to use with preconditioner (generally the same as J) 132c7afd0dbSLois Curfman McInnes - flag - flag indicating whether the matrix sparsity structure has changed 133c7afd0dbSLois Curfman McInnes 1340521c3abSLois Curfman McInnes Options Database Key: 135c7afd0dbSLois Curfman McInnes $ -snes_fd - Activates SNESDefaultComputeHessian() 1360521c3abSLois Curfman McInnes 1370521c3abSLois Curfman McInnes Notes: 1380521c3abSLois Curfman McInnes This routine is slow and expensive, and is not currently optimized 1390521c3abSLois Curfman McInnes to take advantage of sparsity in the problem. Although 1400521c3abSLois Curfman McInnes SNESDefaultComputeHessian() is not recommended for general use 1410521c3abSLois Curfman McInnes in large-scale applications, It can be useful in checking the 1420521c3abSLois Curfman McInnes correctness of a user-provided Hessian. 1430521c3abSLois Curfman McInnes 1440521c3abSLois Curfman McInnes .keywords: SNES, finite differences, Hessian 1450521c3abSLois Curfman McInnes 1469984d6fbSBarry Smith .seealso: SNESSetHessian() 1470521c3abSLois Curfman McInnes @*/ 148c7afd0dbSLois Curfman McInnes int SNESDefaultComputeHessian(SNES snes,Vec x1,Mat *J,Mat *B, 149c7afd0dbSLois Curfman McInnes MatStructure *flag,void *ctx) 1500521c3abSLois Curfman McInnes { 1513a40ed3dSBarry Smith int ierr; 1523a40ed3dSBarry Smith 1533a40ed3dSBarry Smith PetscFunctionBegin; 1543a40ed3dSBarry Smith ierr = SNESDefaultComputeJacobian(snes,x1,J,B,flag,ctx);CHKERRQ(ierr); 1553a40ed3dSBarry Smith PetscFunctionReturn(0); 1560521c3abSLois Curfman McInnes } 157