181e6777dSBarry Smith 281e6777dSBarry Smith #ifndef lint 3*0452661fSBarry Smith static char vcid[] = "$Id: snesmfj.c,v 1.21 1995/11/01 19:12:03 bsmith Exp bsmith $"; 481e6777dSBarry Smith #endif 581e6777dSBarry Smith 6b1f0a012SBarry Smith #include "draw.h" /*I "draw.h" I*/ 7b1f0a012SBarry Smith #include "snes.h" /*I "snes.h" I*/ 881e6777dSBarry Smith 950361f65SLois Curfman McInnes typedef struct { /* default context for matrix-free SNES */ 1039e2f89bSBarry Smith SNES snes; 1139e2f89bSBarry Smith Vec w; 1239e2f89bSBarry Smith } MFCtx_Private; 1339e2f89bSBarry Smith 14b9fa9cd0SBarry Smith int SNESMatrixFreeDestroy_Private(void *ptr) 15b9fa9cd0SBarry Smith { 16b9fa9cd0SBarry Smith int ierr; 17b9fa9cd0SBarry Smith MFCtx_Private *ctx = (MFCtx_Private* ) ptr; 18b9fa9cd0SBarry Smith ierr = VecDestroy(ctx->w); CHKERRQ(ierr); 19*0452661fSBarry Smith PetscFree(ptr); 20b9fa9cd0SBarry Smith return 0; 21b9fa9cd0SBarry Smith } 2250361f65SLois Curfman McInnes 2339e2f89bSBarry Smith /* 2439e2f89bSBarry Smith SNESMatrixFreeMult_Private - Default matrix free form of A*u. 2539e2f89bSBarry Smith */ 2639e2f89bSBarry Smith int SNESMatrixFreeMult_Private(void *ptr,Vec dx,Vec y) 2739e2f89bSBarry Smith { 2839e2f89bSBarry Smith MFCtx_Private *ctx = (MFCtx_Private* ) ptr; 2939e2f89bSBarry Smith SNES snes = ctx->snes; 3039e2f89bSBarry Smith double norm,epsilon = 1.e-8; /* assumes double precision */ 3119a167f6SBarry Smith Scalar h,dot; 3219a167f6SBarry Smith double sum; 3339e2f89bSBarry Smith Scalar mone = -1.0; 3439e2f89bSBarry Smith Vec w = ctx->w,U,F; 3539e2f89bSBarry Smith int ierr; 3639e2f89bSBarry Smith 3778b31e54SBarry Smith ierr = SNESGetSolution(snes,&U); CHKERRQ(ierr); 3878b31e54SBarry Smith ierr = SNESGetFunction(snes,&F); CHKERRQ(ierr); 3950361f65SLois Curfman McInnes 4050361f65SLois Curfman McInnes /* Determine a "good" step size */ 4150361f65SLois Curfman McInnes ierr = VecDot(U,dx,&dot); CHKERRQ(ierr); 42cddf8d76SBarry Smith ierr = VecNorm(dx,NORM_1,&sum); CHKERRQ(ierr); 43cddf8d76SBarry Smith ierr = VecNorm(dx,NORM_2,&norm); CHKERRQ(ierr); 44edd2f0e1SBarry Smith if (sum == 0.0) {dot = 1.0; norm = 1.0;} 4519a167f6SBarry Smith #if defined(PETSC_COMPLEX) 4619a167f6SBarry Smith else if (abs(dot) < 1.e-16*sum && real(dot) >= 0.0) dot = 1.e-16*sum; 4719a167f6SBarry Smith else if (abs(dot) < 0.0 && real(dot) > 1.e-16*sum) dot = -1.e-16*sum; 4819a167f6SBarry Smith #else 49edd2f0e1SBarry Smith else if (dot < 1.e-16*sum && dot >= 0.0) dot = 1.e-16*sum; 5039e2f89bSBarry Smith else if (dot < 0.0 && dot > 1.e-16*sum) dot = -1.e-16*sum; 5119a167f6SBarry Smith #endif 5239e2f89bSBarry Smith h = epsilon*dot/(norm*norm); 5339e2f89bSBarry Smith 5450361f65SLois Curfman McInnes /* Evaluate function at F(x + dx) */ 55195ee392SLois Curfman McInnes ierr = VecWAXPY(&h,dx,U,w); CHKERRQ(ierr); 5678b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,y); CHKERRQ(ierr); 57195ee392SLois Curfman McInnes ierr = VecAXPY(&mone,F,y); CHKERRQ(ierr); 5839e2f89bSBarry Smith h = -1.0/h; 59195ee392SLois Curfman McInnes ierr = VecScale(&h,y); CHKERRQ(ierr); 6039e2f89bSBarry Smith return 0; 6139e2f89bSBarry Smith } 624b828684SBarry Smith /*@C 635392566eSBarry Smith SNESDefaultMatrixFreeMatCreate - Creates a matrix-free matrix 6450361f65SLois Curfman McInnes context for use with a SNES solver. This matrix can be used as 6550361f65SLois Curfman McInnes the Jacobian argument for the routine SNESSetJacobian(). 665392566eSBarry Smith 675392566eSBarry Smith Input Parameters: 685392566eSBarry Smith . x - vector where SNES solution is to be stored. 695392566eSBarry Smith 705392566eSBarry Smith Output Parameters: 715392566eSBarry Smith . J - the matrix-free matrix 725392566eSBarry Smith 7365afa06eSLois Curfman McInnes Notes: 7465afa06eSLois Curfman McInnes The matrix-free matrix context merely contains the function pointers 7565afa06eSLois Curfman McInnes and work space for performing finite difference approximations of 7665afa06eSLois Curfman McInnes matrix operations such as matrix-vector products. 7765afa06eSLois Curfman McInnes 7865afa06eSLois Curfman McInnes The user should call MatDestroy() when finished with the matrix-free 7965afa06eSLois Curfman McInnes matrix context. 8065afa06eSLois Curfman McInnes 8165afa06eSLois Curfman McInnes .keywords: SNES, default, matrix-free, create, matrix 8265afa06eSLois Curfman McInnes 8365afa06eSLois Curfman McInnes .seealso: MatDestroy() 845392566eSBarry Smith @*/ 855392566eSBarry Smith int SNESDefaultMatrixFreeMatCreate(SNES snes,Vec x, Mat *J) 865392566eSBarry Smith { 875392566eSBarry Smith MPI_Comm comm; 885392566eSBarry Smith MFCtx_Private *mfctx; 895392566eSBarry Smith int n,ierr; 905392566eSBarry Smith 91*0452661fSBarry Smith mfctx = (MFCtx_Private *) PetscMalloc(sizeof(MFCtx_Private)); CHKPTRQ(mfctx); 92464493b3SBarry Smith PLogObjectMemory(snes,sizeof(MFCtx_Private)); 935392566eSBarry Smith mfctx->snes = snes; 945392566eSBarry Smith ierr = VecDuplicate(x,&mfctx->w); CHKERRQ(ierr); 95195ee392SLois Curfman McInnes ierr = PetscObjectGetComm((PetscObject)x,&comm); CHKERRQ(ierr); 96195ee392SLois Curfman McInnes ierr = VecGetSize(x,&n); CHKERRQ(ierr); 975392566eSBarry Smith ierr = MatShellCreate(comm,n,n,(void*)mfctx,J); CHKERRQ(ierr); 98195ee392SLois Curfman McInnes ierr = MatShellSetMult(*J,SNESMatrixFreeMult_Private); CHKERRQ(ierr); 99195ee392SLois Curfman McInnes ierr = MatShellSetDestroy(*J,SNESMatrixFreeDestroy_Private); CHKERRQ(ierr); 100b9fa9cd0SBarry Smith PLogObjectParent(*J,mfctx->w); 101d370d78aSBarry Smith PLogObjectParent(snes,*J); 10281e6777dSBarry Smith return 0; 10381e6777dSBarry Smith } 10481e6777dSBarry Smith 105