181e6777dSBarry Smith 281e6777dSBarry Smith #ifndef lint 3*65afa06eSLois Curfman McInnes static char vcid[] = "$Id: snesmfj.c,v 1.13 1995/08/02 04:18:56 bsmith Exp curfman $"; 481e6777dSBarry Smith #endif 581e6777dSBarry Smith 6b1f0a012SBarry Smith #include "draw.h" /*I "draw.h" I*/ 7b1f0a012SBarry Smith #include "snes.h" /*I "snes.h" I*/ 8b9fa9cd0SBarry Smith #include "ptscimpl.h" 9b9fa9cd0SBarry Smith #include "plog.h" 1081e6777dSBarry Smith 1139e2f89bSBarry Smith typedef struct { 1239e2f89bSBarry Smith SNES snes; 1339e2f89bSBarry Smith Vec w; 1439e2f89bSBarry Smith } MFCtx_Private; 1539e2f89bSBarry Smith 16b9fa9cd0SBarry Smith int SNESMatrixFreeDestroy_Private(void *ptr) 17b9fa9cd0SBarry Smith { 18b9fa9cd0SBarry Smith int ierr; 19b9fa9cd0SBarry Smith MFCtx_Private *ctx = (MFCtx_Private* ) ptr; 20b9fa9cd0SBarry Smith ierr = VecDestroy(ctx->w); CHKERRQ(ierr); 21b9fa9cd0SBarry Smith PETSCFREE(ptr); 22b9fa9cd0SBarry Smith return 0; 23b9fa9cd0SBarry Smith } 2439e2f89bSBarry Smith /* 2539e2f89bSBarry Smith SNESMatrixFreeMult_Private - Default matrix free form of A*u. 2639e2f89bSBarry Smith 2739e2f89bSBarry Smith */ 2839e2f89bSBarry Smith int SNESMatrixFreeMult_Private(void *ptr,Vec dx,Vec y) 2939e2f89bSBarry Smith { 3039e2f89bSBarry Smith MFCtx_Private *ctx = (MFCtx_Private* ) ptr; 3139e2f89bSBarry Smith SNES snes = ctx->snes; 3239e2f89bSBarry Smith double norm,epsilon = 1.e-8; /* assumes double precision */ 3319a167f6SBarry Smith Scalar h,dot; 3419a167f6SBarry Smith double sum; 3539e2f89bSBarry Smith Scalar mone = -1.0; 3639e2f89bSBarry Smith Vec w = ctx->w,U,F; 3739e2f89bSBarry Smith int ierr; 3839e2f89bSBarry Smith 3978b31e54SBarry Smith ierr = SNESGetSolution(snes,&U); CHKERRQ(ierr); 4078b31e54SBarry Smith ierr = SNESGetFunction(snes,&F); CHKERRQ(ierr); 4139e2f89bSBarry Smith /* determine a "good" step size */ 4239e2f89bSBarry Smith VecDot(U,dx,&dot); VecASum(dx,&sum); VecNorm(dx,&norm); 43edd2f0e1SBarry Smith if (sum == 0.0) {dot = 1.0; norm = 1.0;} 4419a167f6SBarry Smith #if defined(PETSC_COMPLEX) 4519a167f6SBarry Smith else if (abs(dot) < 1.e-16*sum && real(dot) >= 0.0) dot = 1.e-16*sum; 4619a167f6SBarry Smith else if (abs(dot) < 0.0 && real(dot) > 1.e-16*sum) dot = -1.e-16*sum; 4719a167f6SBarry Smith #else 48edd2f0e1SBarry Smith else if (dot < 1.e-16*sum && dot >= 0.0) dot = 1.e-16*sum; 4939e2f89bSBarry Smith else if (dot < 0.0 && dot > 1.e-16*sum) dot = -1.e-16*sum; 5019a167f6SBarry Smith #endif 5139e2f89bSBarry Smith h = epsilon*dot/(norm*norm); 5239e2f89bSBarry Smith 5339e2f89bSBarry Smith /* evaluate function at F(x + dx) */ 5439e2f89bSBarry Smith VecWAXPY(&h,dx,U,w); 5578b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,y); CHKERRQ(ierr); 5639e2f89bSBarry Smith VecAXPY(&mone,F,y); 5739e2f89bSBarry Smith h = -1.0/h; 5839e2f89bSBarry Smith VecScale(&h,y); 5939e2f89bSBarry Smith return 0; 6039e2f89bSBarry Smith } 6181e6777dSBarry Smith /*@ 625392566eSBarry Smith SNESDefaultMatrixFreeMatCreate - Creates a matrix-free matrix 63*65afa06eSLois Curfman McInnes context for use with a SNES solver. You can use this matrix as the 64*65afa06eSLois Curfman McInnes Jacobian argument for the routine SNESSetJacobian(). 655392566eSBarry Smith 665392566eSBarry Smith Input Parameters: 675392566eSBarry Smith . x - vector where SNES solution is to be stored. 685392566eSBarry Smith 695392566eSBarry Smith Output Parameters: 705392566eSBarry Smith . J - the matrix-free matrix 715392566eSBarry Smith 72*65afa06eSLois Curfman McInnes Notes: 73*65afa06eSLois Curfman McInnes The matrix-free matrix context merely contains the function pointers 74*65afa06eSLois Curfman McInnes and work space for performing finite difference approximations of 75*65afa06eSLois Curfman McInnes matrix operations such as matrix-vector products. 76*65afa06eSLois Curfman McInnes 77*65afa06eSLois Curfman McInnes The user should call MatDestroy() when finished with the matrix-free 78*65afa06eSLois Curfman McInnes matrix context. 79*65afa06eSLois Curfman McInnes 80*65afa06eSLois Curfman McInnes .keywords: SNES, default, matrix-free, create, matrix 81*65afa06eSLois Curfman McInnes 82*65afa06eSLois Curfman McInnes .seealso: MatDestroy() 835392566eSBarry Smith @*/ 845392566eSBarry Smith int SNESDefaultMatrixFreeMatCreate(SNES snes,Vec x, Mat *J) 855392566eSBarry Smith { 865392566eSBarry Smith MPI_Comm comm; 875392566eSBarry Smith MFCtx_Private *mfctx; 885392566eSBarry Smith int n,ierr; 895392566eSBarry Smith 905392566eSBarry Smith mfctx = (MFCtx_Private *) PETSCMALLOC(sizeof(MFCtx_Private));CHKPTRQ(mfctx); 915392566eSBarry Smith mfctx->snes = snes; 925392566eSBarry Smith ierr = VecDuplicate(x,&mfctx->w); CHKERRQ(ierr); 935392566eSBarry Smith PetscObjectGetComm((PetscObject)x,&comm); 945392566eSBarry Smith VecGetSize(x,&n); 955392566eSBarry Smith ierr = MatShellCreate(comm,n,n,(void*)mfctx,J); CHKERRQ(ierr); 965392566eSBarry Smith MatShellSetMult(*J,SNESMatrixFreeMult_Private); 97b9fa9cd0SBarry Smith MatShellSetDestroy(*J,SNESMatrixFreeDestroy_Private); 98b9fa9cd0SBarry Smith PLogObjectParent(*J,mfctx->w); 9981e6777dSBarry Smith return 0; 10081e6777dSBarry Smith } 10181e6777dSBarry Smith 102