xref: /petsc/src/snes/tutorials/ex3k.kokkos.cxx (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c0558f20SBarry Smith static char help[] = "Newton methods to solve u'' + u^{2} = f in parallel. Uses Kokkos\n\\n";
2c0558f20SBarry Smith 
341e469ddSJunchao Zhang #include <petscdmda_kokkos.hpp>
4c0558f20SBarry Smith #include <petscsnes.h>
5c0558f20SBarry Smith 
6c0558f20SBarry Smith /*
7c0558f20SBarry Smith    User-defined application context
8c0558f20SBarry Smith */
9c0558f20SBarry Smith typedef struct {
10c0558f20SBarry Smith   DM          da;      /* distributed array */
11c0558f20SBarry Smith   Vec         F;       /* right-hand-side of PDE */
12c0558f20SBarry Smith   PetscReal   h;       /* mesh spacing */
13c0558f20SBarry Smith } ApplicationCtx;
14c0558f20SBarry Smith 
15c0558f20SBarry Smith /* ------------------------------------------------------------------- */
16c0558f20SBarry Smith /*
17c0558f20SBarry Smith    FormInitialGuess - Computes initial guess.
18c0558f20SBarry Smith 
19c0558f20SBarry Smith    Input/Output Parameter:
20c0558f20SBarry Smith .  x - the solution vector
21c0558f20SBarry Smith */
22c0558f20SBarry Smith PetscErrorCode FormInitialGuess(Vec x)
23c0558f20SBarry Smith {
24c0558f20SBarry Smith   PetscScalar    pfive = .50;
25c0558f20SBarry Smith 
26c0558f20SBarry Smith   PetscFunctionBeginUser;
275f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(x,pfive));
28c0558f20SBarry Smith   PetscFunctionReturn(0);
29c0558f20SBarry Smith }
30c0558f20SBarry Smith 
31c0558f20SBarry Smith /* ------------------------------------------------------------------- */
32c0558f20SBarry Smith /*
3341e469ddSJunchao Zhang    CpuFunction - Evaluates nonlinear function, F(x) on CPU
34c0558f20SBarry Smith 
35c0558f20SBarry Smith    Input Parameters:
36c0558f20SBarry Smith .  snes - the SNES context
37c0558f20SBarry Smith .  x - input vector
38c0558f20SBarry Smith .  ctx - optional user-defined context, as set by SNESSetFunction()
39c0558f20SBarry Smith 
40c0558f20SBarry Smith    Output Parameter:
4141e469ddSJunchao Zhang .  r - function vector
42c0558f20SBarry Smith 
43c0558f20SBarry Smith    Note:
44c0558f20SBarry Smith    The user-defined context can contain any application-specific
45c0558f20SBarry Smith    data needed for the function evaluation.
46c0558f20SBarry Smith */
4741e469ddSJunchao Zhang PetscErrorCode CpuFunction(SNES snes,Vec x,Vec r,void *ctx)
48c0558f20SBarry Smith {
49c0558f20SBarry Smith   ApplicationCtx *user = (ApplicationCtx*) ctx;
50c0558f20SBarry Smith   DM             da    = user->da;
5141e469ddSJunchao Zhang   PetscScalar    *X,*R,*F,d;
52c0558f20SBarry Smith   PetscInt       i,M,xs,xm;
5341e469ddSJunchao Zhang   Vec            xl;
54c0558f20SBarry Smith 
55c0558f20SBarry Smith   PetscFunctionBeginUser;
565f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&xl));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocal(da,x,INSERT_VALUES,xl));
585f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,xl,&X));
595f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,r,&R));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,user->F,&F));
61c0558f20SBarry Smith 
625f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL));
635f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,NULL,&M,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL));
64c0558f20SBarry Smith 
65c0558f20SBarry Smith   if (xs == 0) { /* left boundary */
6641e469ddSJunchao Zhang     R[0] = X[0];
67c0558f20SBarry Smith     xs++;xm--;
68c0558f20SBarry Smith   }
69c0558f20SBarry Smith   if (xs+xm == M) {  /* right boundary */
7041e469ddSJunchao Zhang     R[xs+xm-1] = X[xs+xm-1] - 1.0;
71c0558f20SBarry Smith     xm--;
72c0558f20SBarry Smith   }
73c0558f20SBarry Smith   d = 1.0/(user->h*user->h);
7441e469ddSJunchao Zhang   for (i=xs; i<xs+xm; i++) R[i] = d*(X[i-1] - 2.0*X[i] + X[i+1]) + X[i]*X[i] - F[i];
75c0558f20SBarry Smith 
765f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,xl,&X));
775f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,r,&R));
785f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,user->F,&F));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&xl));
80c0558f20SBarry Smith   PetscFunctionReturn(0);
81c0558f20SBarry Smith }
82c0558f20SBarry Smith 
8341e469ddSJunchao Zhang using DefaultExecutionSpace             = Kokkos::DefaultExecutionSpace;
8441e469ddSJunchao Zhang using DefaultMemorySpace                = Kokkos::DefaultExecutionSpace::memory_space;
8541e469ddSJunchao Zhang using PetscScalarKokkosOffsetView       = Kokkos::Experimental::OffsetView<PetscScalar*,DefaultMemorySpace>;
8641e469ddSJunchao Zhang using ConstPetscScalarKokkosOffsetView  = Kokkos::Experimental::OffsetView<const PetscScalar*,DefaultMemorySpace>;
8741e469ddSJunchao Zhang 
8841e469ddSJunchao Zhang PetscErrorCode KokkosFunction(SNES snes,Vec x,Vec r,void *ctx)
8941e469ddSJunchao Zhang {
9041e469ddSJunchao Zhang   ApplicationCtx                       *user = (ApplicationCtx*) ctx;
9141e469ddSJunchao Zhang   DM                                   da = user->da;
9241e469ddSJunchao Zhang   PetscScalar                          d;
9341e469ddSJunchao Zhang   PetscInt                             M;
9441e469ddSJunchao Zhang   Vec                                  xl;
9541e469ddSJunchao Zhang   PetscScalarKokkosOffsetView          R;
9641e469ddSJunchao Zhang   ConstPetscScalarKokkosOffsetView     X,F;
9741e469ddSJunchao Zhang 
9841e469ddSJunchao Zhang   PetscFunctionBeginUser;
995f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&xl));
1005f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocal(da,x,INSERT_VALUES,xl));
10141e469ddSJunchao Zhang   d    = 1.0/(user->h*user->h);
1025f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,NULL,&M,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL));
1035f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetKokkosOffsetView(da,xl,&X)); /* read only */
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetKokkosOffsetViewWrite(da,r,&R)); /* write only */
1055f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetKokkosOffsetView(da,user->F,&F)); /* read only */
10641e469ddSJunchao Zhang   Kokkos:: parallel_for (Kokkos::RangePolicy<>(R.begin(0),R.end(0)),KOKKOS_LAMBDA (int i) {
10741e469ddSJunchao Zhang     if (i == 0)        R(0) = X(0);        /* left boundary */
10841e469ddSJunchao Zhang     else if (i == M-1) R(i) = X(i) - 1.0;  /* right boundary */
10941e469ddSJunchao Zhang     else               R(i) = d*(X(i-1) - 2.0*X(i) + X(i+1)) + X(i)*X(i) - F(i); /* interior */
11041e469ddSJunchao Zhang   });
1115f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreKokkosOffsetView(da,xl,&X));
1125f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreKokkosOffsetViewWrite(da,r,&R));
1135f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreKokkosOffsetView(da,user->F,&F));
1145f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&xl));
11541e469ddSJunchao Zhang   PetscFunctionReturn(0);
11641e469ddSJunchao Zhang }
11741e469ddSJunchao Zhang 
11841e469ddSJunchao Zhang PetscErrorCode StubFunction(SNES snes ,Vec x,Vec r,void *ctx)
11941e469ddSJunchao Zhang {
12041e469ddSJunchao Zhang   ApplicationCtx                       *user = (ApplicationCtx*) ctx;
12141e469ddSJunchao Zhang   DM                                   da = user->da;
12241e469ddSJunchao Zhang   Vec                                  rk;
12341e469ddSJunchao Zhang   PetscReal                            norm=0;
12441e469ddSJunchao Zhang 
12541e469ddSJunchao Zhang   PetscFunctionBeginUser;
1265f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetGlobalVector(da,&rk));
1275f80ce2aSJacob Faibussowitsch   CHKERRQ(CpuFunction(snes,x,r,ctx));
1285f80ce2aSJacob Faibussowitsch   CHKERRQ(KokkosFunction(snes,x,rk,ctx));
1295f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(rk,-1.0,r));
1305f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(rk,NORM_2,&norm));
1315f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreGlobalVector(da,&rk));
1322c71b3e2SJacob Faibussowitsch   PetscCheckFalse(norm > 1e-6,PETSC_COMM_SELF,PETSC_ERR_PLIB,"KokkosFunction() different from CpuFunction() with a diff norm = %g",norm);
13341e469ddSJunchao Zhang   PetscFunctionReturn(0);
13441e469ddSJunchao Zhang }
135c0558f20SBarry Smith /* ------------------------------------------------------------------- */
136c0558f20SBarry Smith /*
137c0558f20SBarry Smith    FormJacobian - Evaluates Jacobian matrix.
138c0558f20SBarry Smith 
139c0558f20SBarry Smith    Input Parameters:
140c0558f20SBarry Smith .  snes - the SNES context
141c0558f20SBarry Smith .  x - input vector
142c0558f20SBarry Smith .  dummy - optional user-defined context (not used here)
143c0558f20SBarry Smith 
144c0558f20SBarry Smith    Output Parameters:
145c0558f20SBarry Smith .  jac - Jacobian matrix
146c0558f20SBarry Smith .  B - optionally different preconditioning matrix
147c0558f20SBarry Smith .  flag - flag indicating matrix structure
148c0558f20SBarry Smith */
149c0558f20SBarry Smith PetscErrorCode FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void *ctx)
150c0558f20SBarry Smith {
151c0558f20SBarry Smith   ApplicationCtx *user = (ApplicationCtx*) ctx;
152c0558f20SBarry Smith   PetscScalar    *xx,d,A[3];
153c0558f20SBarry Smith   PetscInt       i,j[3],M,xs,xm;
154c0558f20SBarry Smith   DM             da = user->da;
155c0558f20SBarry Smith 
156c0558f20SBarry Smith   PetscFunctionBeginUser;
157c0558f20SBarry Smith   /*
158c0558f20SBarry Smith      Get pointer to vector data
159c0558f20SBarry Smith   */
1605f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArrayRead(da,x,&xx));
1615f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL));
162c0558f20SBarry Smith 
163c0558f20SBarry Smith   /*
164c0558f20SBarry Smith     Get range of locally owned matrix
165c0558f20SBarry Smith   */
1665f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,NULL,&M,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL));
167c0558f20SBarry Smith 
168c0558f20SBarry Smith   /*
169c0558f20SBarry Smith      Determine starting and ending local indices for interior grid points.
170c0558f20SBarry Smith      Set Jacobian entries for boundary points.
171c0558f20SBarry Smith   */
172c0558f20SBarry Smith 
173c0558f20SBarry Smith   if (xs == 0) {  /* left boundary */
174c0558f20SBarry Smith     i = 0; A[0] = 1.0;
175c0558f20SBarry Smith 
1765f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(jac,1,&i,1,&i,A,INSERT_VALUES));
177c0558f20SBarry Smith     xs++;xm--;
178c0558f20SBarry Smith   }
179c0558f20SBarry Smith   if (xs+xm == M) { /* right boundary */
180c0558f20SBarry Smith     i    = M-1;
181c0558f20SBarry Smith     A[0] = 1.0;
1825f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(jac,1,&i,1,&i,A,INSERT_VALUES));
183c0558f20SBarry Smith     xm--;
184c0558f20SBarry Smith   }
185c0558f20SBarry Smith 
186c0558f20SBarry Smith   /*
187c0558f20SBarry Smith      Interior grid points
188c0558f20SBarry Smith       - Note that in this case we set all elements for a particular
189c0558f20SBarry Smith         row at once.
190c0558f20SBarry Smith   */
191c0558f20SBarry Smith   d = 1.0/(user->h*user->h);
192c0558f20SBarry Smith   for (i=xs; i<xs+xm; i++) {
193c0558f20SBarry Smith     j[0] = i - 1; j[1] = i; j[2] = i + 1;
194c0558f20SBarry Smith     A[0] = A[2] = d; A[1] = -2.0*d + 2.0*xx[i];
1955f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(jac,1,&i,3,j,A,INSERT_VALUES));
196c0558f20SBarry Smith   }
197c0558f20SBarry Smith 
198c0558f20SBarry Smith   /*
199c0558f20SBarry Smith      Assemble matrix, using the 2-step process:
200c0558f20SBarry Smith        MatAssemblyBegin(), MatAssemblyEnd().
201c0558f20SBarry Smith      By placing code between these two statements, computations can be
202c0558f20SBarry Smith      done while messages are in transition.
203c0558f20SBarry Smith 
204c0558f20SBarry Smith      Also, restore vector.
205c0558f20SBarry Smith   */
206c0558f20SBarry Smith 
2075f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
2085f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayRead(da,x,&xx));
2095f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
210c0558f20SBarry Smith 
211c0558f20SBarry Smith   PetscFunctionReturn(0);
212c0558f20SBarry Smith }
213c0558f20SBarry Smith 
21441e469ddSJunchao Zhang int main(int argc,char **argv)
21541e469ddSJunchao Zhang {
21641e469ddSJunchao Zhang   SNES           snes;          /* SNES context */
21741e469ddSJunchao Zhang   Mat            J;             /* Jacobian matrix */
21841e469ddSJunchao Zhang   ApplicationCtx ctx;           /* user-defined context */
21941e469ddSJunchao Zhang   Vec            x,r,U,F;       /* vectors */
22041e469ddSJunchao Zhang   PetscScalar    none        = -1.0;
22141e469ddSJunchao Zhang   PetscInt       its,N       = 5,maxit,maxf;
22241e469ddSJunchao Zhang   PetscReal      abstol,rtol,stol,norm;
22341e469ddSJunchao Zhang   PetscBool      viewinitial = PETSC_FALSE;
22441e469ddSJunchao Zhang 
225*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
2265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL));
22741e469ddSJunchao Zhang   ctx.h = 1.0/(N-1);
2285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-view_initial",&viewinitial,NULL));
22941e469ddSJunchao Zhang 
23041e469ddSJunchao Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
23141e469ddSJunchao Zhang      Create nonlinear solver context
23241e469ddSJunchao Zhang      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2335f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESCreate(PETSC_COMM_WORLD,&snes));
23441e469ddSJunchao Zhang 
235c0558f20SBarry Smith   /*
23641e469ddSJunchao Zhang      Create distributed array (DMDA) to manage parallel grid and vectors
237c0558f20SBarry Smith   */
2385f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,N,1,1,NULL,&ctx.da));
2395f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(ctx.da));
2405f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(ctx.da));
24141e469ddSJunchao Zhang 
24241e469ddSJunchao Zhang   /*
24341e469ddSJunchao Zhang      Extract global and local vectors from DMDA; then duplicate for remaining
24441e469ddSJunchao Zhang      vectors that are the same types
24541e469ddSJunchao Zhang   */
2465f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(ctx.da,&x));
2475f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)x,"Approximate Solution"));
2485f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&r));
2495f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&F)); ctx.F = F;
2505f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)F,"Forcing function"));
2515f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&U));
2525f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)U,"Exact Solution"));
25341e469ddSJunchao Zhang 
25441e469ddSJunchao Zhang   /*
25541e469ddSJunchao Zhang      Set function evaluation routine and vector.  Whenever the nonlinear
25641e469ddSJunchao Zhang      solver needs to compute the nonlinear function, it will call this
25741e469ddSJunchao Zhang      routine.
25841e469ddSJunchao Zhang       - Note that the final routine argument is the user-defined
25941e469ddSJunchao Zhang         context that provides application-specific data for the
26041e469ddSJunchao Zhang         function evaluation routine.
26141e469ddSJunchao Zhang 
262a5b23f4aSJose E. Roman      At the beginning, one can use a stub function that checks the Kokkos version
26341e469ddSJunchao Zhang      against the CPU version to quickly expose errors.
2645f80ce2aSJacob Faibussowitsch      CHKERRQ(SNESSetFunction(snes,r,StubFunction,&ctx));
26541e469ddSJunchao Zhang   */
2665f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFunction(snes,r,KokkosFunction,&ctx));
26741e469ddSJunchao Zhang 
26841e469ddSJunchao Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
26941e469ddSJunchao Zhang      Create matrix data structure; set Jacobian evaluation routine
27041e469ddSJunchao Zhang      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2715f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateMatrix(ctx.da,&J));
27241e469ddSJunchao Zhang 
27341e469ddSJunchao Zhang   /*
27441e469ddSJunchao Zhang      Set Jacobian matrix data structure and default Jacobian evaluation
27541e469ddSJunchao Zhang      routine.  Whenever the nonlinear solver needs to compute the
27641e469ddSJunchao Zhang      Jacobian matrix, it will call this routine.
27741e469ddSJunchao Zhang       - Note that the final routine argument is the user-defined
27841e469ddSJunchao Zhang         context that provides application-specific data for the
27941e469ddSJunchao Zhang         Jacobian evaluation routine.
28041e469ddSJunchao Zhang   */
2815f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetJacobian(snes,J,J,FormJacobian,&ctx));
2825f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFromOptions(snes));
2835f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetTolerances(snes,&abstol,&rtol,&stol,&maxit,&maxf));
2845f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"atol=%g, rtol=%g, stol=%g, maxit=%D, maxf=%D\n",(double)abstol,(double)rtol,(double)stol,maxit,maxf));
28541e469ddSJunchao Zhang 
28641e469ddSJunchao Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
28741e469ddSJunchao Zhang      Initialize application:
28841e469ddSJunchao Zhang      Store forcing function of PDE and exact solution
28941e469ddSJunchao Zhang    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2909000cb80SJunchao Zhang   {
2919000cb80SJunchao Zhang     PetscScalarKokkosOffsetView FF,UU;
2925f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecGetKokkosOffsetViewWrite(ctx.da,F,&FF));
2935f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecGetKokkosOffsetViewWrite(ctx.da,U,&UU));
29441e469ddSJunchao Zhang     Kokkos:: parallel_for (Kokkos::RangePolicy<>(FF.begin(0),FF.end(0)),KOKKOS_LAMBDA (int i) {
29541e469ddSJunchao Zhang       PetscReal xp = i*ctx.h;
29641e469ddSJunchao Zhang       FF(i) = 6.0*xp + pow(xp+1.e-12,6.0); /* +1.e-12 is to prevent 0^6 */
29741e469ddSJunchao Zhang       UU(i) = xp*xp*xp;
29841e469ddSJunchao Zhang     });
2995f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecRestoreKokkosOffsetViewWrite(ctx.da,F,&FF));
3005f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecRestoreKokkosOffsetViewWrite(ctx.da,U,&UU));
3019000cb80SJunchao Zhang   }
30241e469ddSJunchao Zhang 
30341e469ddSJunchao Zhang   if (viewinitial) {
3045f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(U,NULL));
3055f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(F,NULL));
30641e469ddSJunchao Zhang   }
30741e469ddSJunchao Zhang 
30841e469ddSJunchao Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
30941e469ddSJunchao Zhang      Evaluate initial guess; then solve nonlinear system
31041e469ddSJunchao Zhang    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
31141e469ddSJunchao Zhang 
31241e469ddSJunchao Zhang   /*
31341e469ddSJunchao Zhang      Note: The user should initialize the vector, x, with the initial guess
31441e469ddSJunchao Zhang      for the nonlinear solver prior to calling SNESSolve().  In particular,
31541e469ddSJunchao Zhang      to employ an initial guess of zero, the user should explicitly set
31641e469ddSJunchao Zhang      this vector to zero by calling VecSet().
31741e469ddSJunchao Zhang   */
3185f80ce2aSJacob Faibussowitsch   CHKERRQ(FormInitialGuess(x));
3195f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSolve(snes,NULL,x));
3205f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetIterationNumber(snes,&its));
3215f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %D\n",its));
32241e469ddSJunchao Zhang 
32341e469ddSJunchao Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
32441e469ddSJunchao Zhang      Check solution and clean up
32541e469ddSJunchao Zhang    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
32641e469ddSJunchao Zhang   /*
32741e469ddSJunchao Zhang      Check the error
32841e469ddSJunchao Zhang   */
3295f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(x,none,U));
3305f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(x,NORM_2,&norm));
3315f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g Iterations %D\n",(double)norm,its));
33241e469ddSJunchao Zhang 
33341e469ddSJunchao Zhang   /*
33441e469ddSJunchao Zhang      Free work space.  All PETSc objects should be destroyed when they
33541e469ddSJunchao Zhang      are no longer needed.
33641e469ddSJunchao Zhang   */
3375f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
3385f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&r));
3395f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&U));
3405f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&F));
3415f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&J));
3425f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESDestroy(&snes));
3435f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&ctx.da));
344*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
345*b122ec5aSJacob Faibussowitsch   return 0;
34641e469ddSJunchao Zhang }
347c0558f20SBarry Smith 
348c0558f20SBarry Smith /*TEST
349c0558f20SBarry Smith 
350c0558f20SBarry Smith    build:
35141e469ddSJunchao Zhang      requires: kokkos_kernels
352c0558f20SBarry Smith 
353c0558f20SBarry Smith    test:
35487699089SJunchao Zhang      requires: kokkos_kernels !complex !single
355c0558f20SBarry Smith      nsize: 2
35641e469ddSJunchao Zhang      args: -dm_vec_type kokkos -dm_mat_type aijkokkos -view_initial -snes_monitor
357c0558f20SBarry Smith      output_file: output/ex3k_1.out
358c0558f20SBarry Smith 
359c0558f20SBarry Smith TEST*/
360