xref: /petsc/src/snes/mf/snesmfj.c (revision b1f0a012544317ce27ef4696e4e83969356d0f0f)
181e6777dSBarry Smith 
281e6777dSBarry Smith #ifndef lint
3*b1f0a012SBarry Smith static char vcid[] = "$Id: snesmfj.c,v 1.11 1995/06/23 12:41:49 bsmith Exp bsmith $";
481e6777dSBarry Smith #endif
581e6777dSBarry Smith 
6*b1f0a012SBarry Smith #include "draw.h"   /*I  "draw.h"   I*/
7*b1f0a012SBarry Smith #include "snes.h"   /*I  "snes.h"   I*/
881e6777dSBarry Smith 
939e2f89bSBarry Smith typedef struct {
1039e2f89bSBarry Smith   SNES snes;
1139e2f89bSBarry Smith   Vec  w;
1239e2f89bSBarry Smith } MFCtx_Private;
1339e2f89bSBarry Smith 
1439e2f89bSBarry Smith /*
1539e2f89bSBarry Smith     SNESMatrixFreeMult_Private - Default matrix free form of A*u.
1639e2f89bSBarry Smith 
1739e2f89bSBarry Smith */
1839e2f89bSBarry Smith int SNESMatrixFreeMult_Private(void *ptr,Vec dx,Vec y)
1939e2f89bSBarry Smith {
2039e2f89bSBarry Smith   MFCtx_Private *ctx = (MFCtx_Private* ) ptr;
2139e2f89bSBarry Smith   SNES          snes = ctx->snes;
2239e2f89bSBarry Smith   double        norm,epsilon = 1.e-8; /* assumes double precision */
2319a167f6SBarry Smith   Scalar        h,dot;
2419a167f6SBarry Smith   double        sum;
2539e2f89bSBarry Smith   Scalar        mone = -1.0;
2639e2f89bSBarry Smith   Vec           w = ctx->w,U,F;
2739e2f89bSBarry Smith   int           ierr;
2839e2f89bSBarry Smith 
2978b31e54SBarry Smith   ierr = SNESGetSolution(snes,&U); CHKERRQ(ierr);
3078b31e54SBarry Smith   ierr = SNESGetFunction(snes,&F); CHKERRQ(ierr);
3139e2f89bSBarry Smith   /* determine a "good" step size */
3239e2f89bSBarry Smith   VecDot(U,dx,&dot); VecASum(dx,&sum); VecNorm(dx,&norm);
33edd2f0e1SBarry Smith   if (sum == 0.0) {dot = 1.0; norm = 1.0;}
3419a167f6SBarry Smith #if defined(PETSC_COMPLEX)
3519a167f6SBarry Smith   else if (abs(dot) < 1.e-16*sum && real(dot) >= 0.0) dot = 1.e-16*sum;
3619a167f6SBarry Smith   else if (abs(dot) < 0.0 && real(dot) > 1.e-16*sum) dot = -1.e-16*sum;
3719a167f6SBarry Smith #else
38edd2f0e1SBarry Smith   else if (dot < 1.e-16*sum && dot >= 0.0) dot = 1.e-16*sum;
3939e2f89bSBarry Smith   else if (dot < 0.0 && dot > 1.e-16*sum) dot = -1.e-16*sum;
4019a167f6SBarry Smith #endif
4139e2f89bSBarry Smith   h = epsilon*dot/(norm*norm);
4239e2f89bSBarry Smith 
4339e2f89bSBarry Smith   /* evaluate function at F(x + dx) */
4439e2f89bSBarry Smith   VecWAXPY(&h,dx,U,w);
4578b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,y); CHKERRQ(ierr);
4639e2f89bSBarry Smith   VecAXPY(&mone,F,y);
4739e2f89bSBarry Smith   h = -1.0/h;
4839e2f89bSBarry Smith   VecScale(&h,y);
4939e2f89bSBarry Smith   return 0;
5039e2f89bSBarry Smith }
5181e6777dSBarry Smith /*@
525392566eSBarry Smith      SNESDefaultMatrixFreeMatCreate - Creates a matrix-free matrix
535392566eSBarry Smith          for use with SNES solver. You may use this matrix as
545392566eSBarry Smith          Jacobian argument for the routine SNESSetJacobian. This is
555392566eSBarry Smith          most useful when you are using finite differences for a
565392566eSBarry Smith          matrix free Newton method but explictly are forming a
575392566eSBarry Smith          preconditioner matrix.
585392566eSBarry Smith 
595392566eSBarry Smith   Input Parameters:
605392566eSBarry Smith .  x - vector where SNES solution is to be stored.
615392566eSBarry Smith 
625392566eSBarry Smith   Output Parameters:
635392566eSBarry Smith .  J - the matrix-free matrix
645392566eSBarry Smith 
655392566eSBarry Smith @*/
665392566eSBarry Smith int SNESDefaultMatrixFreeMatCreate(SNES snes,Vec x, Mat *J)
675392566eSBarry Smith {
685392566eSBarry Smith   MPI_Comm      comm;
695392566eSBarry Smith   MFCtx_Private *mfctx;
705392566eSBarry Smith   int           n,ierr;
715392566eSBarry Smith 
725392566eSBarry Smith   mfctx = (MFCtx_Private *) PETSCMALLOC(sizeof(MFCtx_Private));CHKPTRQ(mfctx);
735392566eSBarry Smith   mfctx->snes = snes;
745392566eSBarry Smith   ierr = VecDuplicate(x,&mfctx->w); CHKERRQ(ierr);
755392566eSBarry Smith   PetscObjectGetComm((PetscObject)x,&comm);
765392566eSBarry Smith   VecGetSize(x,&n);
775392566eSBarry Smith   ierr = MatShellCreate(comm,n,n,(void*)mfctx,J); CHKERRQ(ierr);
785392566eSBarry Smith   MatShellSetMult(*J,SNESMatrixFreeMult_Private);
795392566eSBarry Smith   return 0;
805392566eSBarry Smith }
815392566eSBarry Smith /*@
8239e2f89bSBarry Smith    SNESDefaultMatrixFreeComputeJacobian - Computes Jacobian using finite
8371dce806SLois Curfman McInnes    differences, matrix-free style.
8481e6777dSBarry Smith 
8581e6777dSBarry Smith    Input Parameters:
8681e6777dSBarry Smith .  x - compute Jacobian at this point
87ad960d00SLois Curfman McInnes .  ctx - application's function context, as set with SNESSetFunction()
8881e6777dSBarry Smith 
8981e6777dSBarry Smith    Output Parameters:
9081e6777dSBarry Smith .  J - Jacobian
9181e6777dSBarry Smith .  B - preconditioner, same as Jacobian
92ad960d00SLois Curfman McInnes .  flag - matrix flag
93ad960d00SLois Curfman McInnes 
94ad960d00SLois Curfman McInnes    Options Database Key:
95ad960d00SLois Curfman McInnes $  -snes_mf
9681e6777dSBarry Smith 
9771dce806SLois Curfman McInnes .keywords: SNES, finite differences, Jacobian
9881e6777dSBarry Smith 
9971dce806SLois Curfman McInnes .seealso: SNESSetJacobian(), SNESTestJacobian()
10081e6777dSBarry Smith @*/
10139e2f89bSBarry Smith int SNESDefaultMatrixFreeComputeJacobian(SNES snes, Vec x1,Mat *J,Mat *B,
102edd2f0e1SBarry Smith                                          MatStructure *flag,void *ctx)
10381e6777dSBarry Smith {
10439e2f89bSBarry Smith   int         n,ierr;
1057f223b93SBarry Smith   MatType     type;
1067f223b93SBarry Smith 
10739e2f89bSBarry Smith   VecGetSize(x1,&n);
1087f223b93SBarry Smith   if (*J) MatGetType(*J,&type);
1097f223b93SBarry Smith   if (!*J || type != MATSHELL) {
1107f223b93SBarry Smith     MPI_Comm    comm;
11139e2f89bSBarry Smith     MFCtx_Private *mfctx;
11239e2f89bSBarry Smith     /* first time in, therefore build datastructures */
11378b31e54SBarry Smith     mfctx = (MFCtx_Private *) PETSCMALLOC(sizeof(MFCtx_Private)); CHKPTRQ(mfctx);
11439e2f89bSBarry Smith     mfctx->snes = snes;
11578b31e54SBarry Smith     ierr = VecDuplicate(x1,&mfctx->w); CHKERRQ(ierr);
1167f223b93SBarry Smith     PetscObjectGetComm((PetscObject)x1,&comm);
11778b31e54SBarry Smith     ierr = MatShellCreate(comm,n,n,(void*)mfctx,J); CHKERRQ(ierr);
11839e2f89bSBarry Smith     MatShellSetMult(*J,SNESMatrixFreeMult_Private);
11939e2f89bSBarry Smith     *B = *J;
12081e6777dSBarry Smith   }
12181e6777dSBarry Smith   return 0;
12281e6777dSBarry Smith }
12381e6777dSBarry Smith 
124