xref: /petsc/src/snes/tutorials/ex3k.kokkos.cxx (revision 11486bccf1aeea1ca5536228f99d437b39bdaca6)
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 */
9*11486bccSBarry Smith typedef struct
10*11486bccSBarry Smith {
11c0558f20SBarry Smith   DM        da; /* distributed array */
12c0558f20SBarry Smith   Vec       F;  /* right-hand-side of PDE */
13c0558f20SBarry Smith   PetscReal h;  /* mesh spacing */
14c0558f20SBarry Smith } ApplicationCtx;
15c0558f20SBarry Smith 
16c0558f20SBarry Smith /* ------------------------------------------------------------------- */
17c0558f20SBarry Smith /*
18c0558f20SBarry Smith    FormInitialGuess - Computes initial guess.
19c0558f20SBarry Smith 
20c0558f20SBarry Smith    Input/Output Parameter:
21c0558f20SBarry Smith .  x - the solution vector
22c0558f20SBarry Smith */
23c0558f20SBarry Smith PetscErrorCode FormInitialGuess(Vec x)
24c0558f20SBarry Smith {
25c0558f20SBarry Smith   PetscScalar pfive = .50;
26c0558f20SBarry Smith 
27c0558f20SBarry Smith   PetscFunctionBeginUser;
289566063dSJacob Faibussowitsch   PetscCall(VecSet(x, pfive));
29c0558f20SBarry Smith   PetscFunctionReturn(0);
30c0558f20SBarry Smith }
31c0558f20SBarry Smith 
32c0558f20SBarry Smith /* ------------------------------------------------------------------- */
33c0558f20SBarry Smith /*
3441e469ddSJunchao Zhang    CpuFunction - Evaluates nonlinear function, F(x) on CPU
35c0558f20SBarry Smith 
36c0558f20SBarry Smith    Input Parameters:
37c0558f20SBarry Smith .  snes - the SNES context
38c0558f20SBarry Smith .  x - input vector
39c0558f20SBarry Smith .  ctx - optional user-defined context, as set by SNESSetFunction()
40c0558f20SBarry Smith 
41c0558f20SBarry Smith    Output Parameter:
4241e469ddSJunchao Zhang .  r - function vector
43c0558f20SBarry Smith 
44c0558f20SBarry Smith    Note:
45c0558f20SBarry Smith    The user-defined context can contain any application-specific
46c0558f20SBarry Smith    data needed for the function evaluation.
47c0558f20SBarry Smith */
4841e469ddSJunchao Zhang PetscErrorCode CpuFunction(SNES snes, Vec x, Vec r, void *ctx)
49c0558f20SBarry Smith {
50c0558f20SBarry Smith   ApplicationCtx *user = (ApplicationCtx *)ctx;
51c0558f20SBarry Smith   DM              da   = user->da;
5241e469ddSJunchao Zhang   PetscScalar    *X, *R, *F, d;
53c0558f20SBarry Smith   PetscInt        i, M, xs, xm;
5441e469ddSJunchao Zhang   Vec             xl;
55c0558f20SBarry Smith 
56c0558f20SBarry Smith   PetscFunctionBeginUser;
579566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &xl));
589566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocal(da, x, INSERT_VALUES, xl));
599566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, xl, &X));
609566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, r, &R));
619566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, user->F, &F));
62c0558f20SBarry Smith 
639566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL));
649566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, NULL, &M, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
65c0558f20SBarry Smith 
66c0558f20SBarry Smith   if (xs == 0) { /* left boundary */
6741e469ddSJunchao Zhang     R[0] = X[0];
68*11486bccSBarry Smith     xs++;
69*11486bccSBarry Smith     xm--;
70c0558f20SBarry Smith   }
71c0558f20SBarry Smith   if (xs + xm == M) { /* right boundary */
7241e469ddSJunchao Zhang     R[xs + xm - 1] = X[xs + xm - 1] - 1.0;
73c0558f20SBarry Smith     xm--;
74c0558f20SBarry Smith   }
75c0558f20SBarry Smith   d = 1.0 / (user->h * user->h);
7641e469ddSJunchao 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];
77c0558f20SBarry Smith 
789566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, xl, &X));
799566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, r, &R));
809566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, user->F, &F));
819566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &xl));
82c0558f20SBarry Smith   PetscFunctionReturn(0);
83c0558f20SBarry Smith }
84c0558f20SBarry Smith 
8541e469ddSJunchao Zhang using DefaultExecutionSpace            = Kokkos::DefaultExecutionSpace;
8641e469ddSJunchao Zhang using DefaultMemorySpace               = Kokkos::DefaultExecutionSpace::memory_space;
8741e469ddSJunchao Zhang using PetscScalarKokkosOffsetView      = Kokkos::Experimental::OffsetView<PetscScalar *, DefaultMemorySpace>;
8841e469ddSJunchao Zhang using ConstPetscScalarKokkosOffsetView = Kokkos::Experimental::OffsetView<const PetscScalar *, DefaultMemorySpace>;
8941e469ddSJunchao Zhang 
9041e469ddSJunchao Zhang PetscErrorCode KokkosFunction(SNES snes, Vec x, Vec r, void *ctx)
9141e469ddSJunchao Zhang {
9241e469ddSJunchao Zhang   ApplicationCtx                  *user = (ApplicationCtx *)ctx;
9341e469ddSJunchao Zhang   DM                               da   = user->da;
9441e469ddSJunchao Zhang   PetscScalar                      d;
9541e469ddSJunchao Zhang   PetscInt                         M;
9641e469ddSJunchao Zhang   Vec                              xl;
9741e469ddSJunchao Zhang   PetscScalarKokkosOffsetView      R;
9841e469ddSJunchao Zhang   ConstPetscScalarKokkosOffsetView X, F;
9941e469ddSJunchao Zhang 
10041e469ddSJunchao Zhang   PetscFunctionBeginUser;
1019566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &xl));
1029566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocal(da, x, INSERT_VALUES, xl));
10341e469ddSJunchao Zhang   d = 1.0 / (user->h * user->h);
1049566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, NULL, &M, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
1059566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetKokkosOffsetView(da, xl, &X));      /* read only */
1069566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetKokkosOffsetViewWrite(da, r, &R));  /* write only */
1079566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetKokkosOffsetView(da, user->F, &F)); /* read only */
108*11486bccSBarry Smith   Kokkos::parallel_for(
109*11486bccSBarry Smith     Kokkos::RangePolicy<>(R.begin(0), R.end(0)), KOKKOS_LAMBDA(int i) {
11041e469ddSJunchao Zhang       if (i == 0) R(0) = X(0);                                                 /* left boundary */
11141e469ddSJunchao Zhang       else if (i == M - 1) R(i) = X(i) - 1.0;                                  /* right boundary */
11241e469ddSJunchao Zhang       else R(i) = d * (X(i - 1) - 2.0 * X(i) + X(i + 1)) + X(i) * X(i) - F(i); /* interior */
11341e469ddSJunchao Zhang     });
1149566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreKokkosOffsetView(da, xl, &X));
1159566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreKokkosOffsetViewWrite(da, r, &R));
1169566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreKokkosOffsetView(da, user->F, &F));
1179566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &xl));
11841e469ddSJunchao Zhang   PetscFunctionReturn(0);
11941e469ddSJunchao Zhang }
12041e469ddSJunchao Zhang 
12141e469ddSJunchao Zhang PetscErrorCode StubFunction(SNES snes, Vec x, Vec r, void *ctx)
12241e469ddSJunchao Zhang {
12341e469ddSJunchao Zhang   ApplicationCtx *user = (ApplicationCtx *)ctx;
12441e469ddSJunchao Zhang   DM              da   = user->da;
12541e469ddSJunchao Zhang   Vec             rk;
12641e469ddSJunchao Zhang   PetscReal       norm = 0;
12741e469ddSJunchao Zhang 
12841e469ddSJunchao Zhang   PetscFunctionBeginUser;
1299566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(da, &rk));
1309566063dSJacob Faibussowitsch   PetscCall(CpuFunction(snes, x, r, ctx));
1319566063dSJacob Faibussowitsch   PetscCall(KokkosFunction(snes, x, rk, ctx));
1329566063dSJacob Faibussowitsch   PetscCall(VecAXPY(rk, -1.0, r));
1339566063dSJacob Faibussowitsch   PetscCall(VecNorm(rk, NORM_2, &norm));
1349566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(da, &rk));
13563a3b9bcSJacob Faibussowitsch   PetscCheck(norm <= 1e-6, PETSC_COMM_SELF, PETSC_ERR_PLIB, "KokkosFunction() different from CpuFunction() with a diff norm = %g", (double)norm);
13641e469ddSJunchao Zhang   PetscFunctionReturn(0);
13741e469ddSJunchao Zhang }
138c0558f20SBarry Smith /* ------------------------------------------------------------------- */
139c0558f20SBarry Smith /*
140c0558f20SBarry Smith    FormJacobian - Evaluates Jacobian matrix.
141c0558f20SBarry Smith 
142c0558f20SBarry Smith    Input Parameters:
143c0558f20SBarry Smith .  snes - the SNES context
144c0558f20SBarry Smith .  x - input vector
145c0558f20SBarry Smith .  dummy - optional user-defined context (not used here)
146c0558f20SBarry Smith 
147c0558f20SBarry Smith    Output Parameters:
148c0558f20SBarry Smith .  jac - Jacobian matrix
149c0558f20SBarry Smith .  B - optionally different preconditioning matrix
150c0558f20SBarry Smith .  flag - flag indicating matrix structure
151c0558f20SBarry Smith */
152c0558f20SBarry Smith PetscErrorCode FormJacobian(SNES snes, Vec x, Mat jac, Mat B, void *ctx)
153c0558f20SBarry Smith {
154c0558f20SBarry Smith   ApplicationCtx *user = (ApplicationCtx *)ctx;
155c0558f20SBarry Smith   PetscScalar    *xx, d, A[3];
156c0558f20SBarry Smith   PetscInt        i, j[3], M, xs, xm;
157c0558f20SBarry Smith   DM              da = user->da;
158c0558f20SBarry Smith 
159c0558f20SBarry Smith   PetscFunctionBeginUser;
160c0558f20SBarry Smith   /*
161c0558f20SBarry Smith      Get pointer to vector data
162c0558f20SBarry Smith   */
1639566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, x, &xx));
1649566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL));
165c0558f20SBarry Smith 
166c0558f20SBarry Smith   /*
167c0558f20SBarry Smith     Get range of locally owned matrix
168c0558f20SBarry Smith   */
1699566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, NULL, &M, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
170c0558f20SBarry Smith 
171c0558f20SBarry Smith   /*
172c0558f20SBarry Smith      Determine starting and ending local indices for interior grid points.
173c0558f20SBarry Smith      Set Jacobian entries for boundary points.
174c0558f20SBarry Smith   */
175c0558f20SBarry Smith 
176c0558f20SBarry Smith   if (xs == 0) { /* left boundary */
177*11486bccSBarry Smith     i    = 0;
178*11486bccSBarry Smith     A[0] = 1.0;
179c0558f20SBarry Smith 
1809566063dSJacob Faibussowitsch     PetscCall(MatSetValues(jac, 1, &i, 1, &i, A, INSERT_VALUES));
181*11486bccSBarry Smith     xs++;
182*11486bccSBarry Smith     xm--;
183c0558f20SBarry Smith   }
184c0558f20SBarry Smith   if (xs + xm == M) { /* right boundary */
185c0558f20SBarry Smith     i    = M - 1;
186c0558f20SBarry Smith     A[0] = 1.0;
1879566063dSJacob Faibussowitsch     PetscCall(MatSetValues(jac, 1, &i, 1, &i, A, INSERT_VALUES));
188c0558f20SBarry Smith     xm--;
189c0558f20SBarry Smith   }
190c0558f20SBarry Smith 
191c0558f20SBarry Smith   /*
192c0558f20SBarry Smith      Interior grid points
193c0558f20SBarry Smith       - Note that in this case we set all elements for a particular
194c0558f20SBarry Smith         row at once.
195c0558f20SBarry Smith   */
196c0558f20SBarry Smith   d = 1.0 / (user->h * user->h);
197c0558f20SBarry Smith   for (i = xs; i < xs + xm; i++) {
198*11486bccSBarry Smith     j[0] = i - 1;
199*11486bccSBarry Smith     j[1] = i;
200*11486bccSBarry Smith     j[2] = i + 1;
201*11486bccSBarry Smith     A[0] = A[2] = d;
202*11486bccSBarry Smith     A[1]        = -2.0 * d + 2.0 * xx[i];
2039566063dSJacob Faibussowitsch     PetscCall(MatSetValues(jac, 1, &i, 3, j, A, INSERT_VALUES));
204c0558f20SBarry Smith   }
205c0558f20SBarry Smith 
206c0558f20SBarry Smith   /*
207c0558f20SBarry Smith      Assemble matrix, using the 2-step process:
208c0558f20SBarry Smith        MatAssemblyBegin(), MatAssemblyEnd().
209c0558f20SBarry Smith      By placing code between these two statements, computations can be
210c0558f20SBarry Smith      done while messages are in transition.
211c0558f20SBarry Smith 
212c0558f20SBarry Smith      Also, restore vector.
213c0558f20SBarry Smith   */
214c0558f20SBarry Smith 
2159566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
2169566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, x, &xx));
2179566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
218c0558f20SBarry Smith 
219c0558f20SBarry Smith   PetscFunctionReturn(0);
220c0558f20SBarry Smith }
221c0558f20SBarry Smith 
22241e469ddSJunchao Zhang int main(int argc, char **argv)
22341e469ddSJunchao Zhang {
22441e469ddSJunchao Zhang   SNES           snes;       /* SNES context */
22541e469ddSJunchao Zhang   Mat            J;          /* Jacobian matrix */
22641e469ddSJunchao Zhang   ApplicationCtx ctx;        /* user-defined context */
22741e469ddSJunchao Zhang   Vec            x, r, U, F; /* vectors */
22841e469ddSJunchao Zhang   PetscScalar    none = -1.0;
22941e469ddSJunchao Zhang   PetscInt       its, N = 5, maxit, maxf;
23041e469ddSJunchao Zhang   PetscReal      abstol, rtol, stol, norm;
23141e469ddSJunchao Zhang   PetscBool      viewinitial = PETSC_FALSE;
23241e469ddSJunchao Zhang 
233327415f7SBarry Smith   PetscFunctionBeginUser;
2349566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
2359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &N, NULL));
23641e469ddSJunchao Zhang   ctx.h = 1.0 / (N - 1);
2379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-view_initial", &viewinitial, NULL));
23841e469ddSJunchao Zhang 
23941e469ddSJunchao Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
24041e469ddSJunchao Zhang      Create nonlinear solver context
24141e469ddSJunchao Zhang      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2429566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
24341e469ddSJunchao Zhang 
244c0558f20SBarry Smith   /*
24541e469ddSJunchao Zhang      Create distributed array (DMDA) to manage parallel grid and vectors
246c0558f20SBarry Smith   */
2479566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, N, 1, 1, NULL, &ctx.da));
2489566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(ctx.da));
2499566063dSJacob Faibussowitsch   PetscCall(DMSetUp(ctx.da));
25041e469ddSJunchao Zhang 
25141e469ddSJunchao Zhang   /*
25241e469ddSJunchao Zhang      Extract global and local vectors from DMDA; then duplicate for remaining
25341e469ddSJunchao Zhang      vectors that are the same types
25441e469ddSJunchao Zhang   */
2559566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(ctx.da, &x));
2569566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)x, "Approximate Solution"));
2579566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &r));
258*11486bccSBarry Smith   PetscCall(VecDuplicate(x, &F));
259*11486bccSBarry Smith   ctx.F = F;
2609566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)F, "Forcing function"));
2619566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &U));
2629566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)U, "Exact Solution"));
26341e469ddSJunchao Zhang 
26441e469ddSJunchao Zhang   /*
26541e469ddSJunchao Zhang      Set function evaluation routine and vector.  Whenever the nonlinear
26641e469ddSJunchao Zhang      solver needs to compute the nonlinear function, it will call this
26741e469ddSJunchao Zhang      routine.
26841e469ddSJunchao Zhang       - Note that the final routine argument is the user-defined
26941e469ddSJunchao Zhang         context that provides application-specific data for the
27041e469ddSJunchao Zhang         function evaluation routine.
27141e469ddSJunchao Zhang 
272a5b23f4aSJose E. Roman      At the beginning, one can use a stub function that checks the Kokkos version
27341e469ddSJunchao Zhang      against the CPU version to quickly expose errors.
2749566063dSJacob Faibussowitsch      PetscCall(SNESSetFunction(snes,r,StubFunction,&ctx));
27541e469ddSJunchao Zhang   */
2769566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, r, KokkosFunction, &ctx));
27741e469ddSJunchao Zhang 
27841e469ddSJunchao Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
27941e469ddSJunchao Zhang      Create matrix data structure; set Jacobian evaluation routine
28041e469ddSJunchao Zhang      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2819566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(ctx.da, &J));
28241e469ddSJunchao Zhang 
28341e469ddSJunchao Zhang   /*
28441e469ddSJunchao Zhang      Set Jacobian matrix data structure and default Jacobian evaluation
28541e469ddSJunchao Zhang      routine.  Whenever the nonlinear solver needs to compute the
28641e469ddSJunchao Zhang      Jacobian matrix, it will call this routine.
28741e469ddSJunchao Zhang       - Note that the final routine argument is the user-defined
28841e469ddSJunchao Zhang         context that provides application-specific data for the
28941e469ddSJunchao Zhang         Jacobian evaluation routine.
29041e469ddSJunchao Zhang   */
2919566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, J, J, FormJacobian, &ctx));
2929566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
2939566063dSJacob Faibussowitsch   PetscCall(SNESGetTolerances(snes, &abstol, &rtol, &stol, &maxit, &maxf));
29463a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "atol=%g, rtol=%g, stol=%g, maxit=%" PetscInt_FMT ", maxf=%" PetscInt_FMT "\n", (double)abstol, (double)rtol, (double)stol, maxit, maxf));
29541e469ddSJunchao Zhang 
29641e469ddSJunchao Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
29741e469ddSJunchao Zhang      Initialize application:
29841e469ddSJunchao Zhang      Store forcing function of PDE and exact solution
29941e469ddSJunchao Zhang    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3009000cb80SJunchao Zhang   {
3019000cb80SJunchao Zhang     PetscScalarKokkosOffsetView FF, UU;
3029566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetKokkosOffsetViewWrite(ctx.da, F, &FF));
3039566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetKokkosOffsetViewWrite(ctx.da, U, &UU));
304*11486bccSBarry Smith     Kokkos::parallel_for(
305*11486bccSBarry Smith       Kokkos::RangePolicy<>(FF.begin(0), FF.end(0)), KOKKOS_LAMBDA(int i) {
30641e469ddSJunchao Zhang         PetscReal xp = i * ctx.h;
30741e469ddSJunchao Zhang         FF(i)        = 6.0 * xp + pow(xp + 1.e-12, 6.0); /* +1.e-12 is to prevent 0^6 */
30841e469ddSJunchao Zhang         UU(i)        = xp * xp * xp;
30941e469ddSJunchao Zhang       });
3109566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreKokkosOffsetViewWrite(ctx.da, F, &FF));
3119566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreKokkosOffsetViewWrite(ctx.da, U, &UU));
3129000cb80SJunchao Zhang   }
31341e469ddSJunchao Zhang 
31441e469ddSJunchao Zhang   if (viewinitial) {
3159566063dSJacob Faibussowitsch     PetscCall(VecView(U, NULL));
3169566063dSJacob Faibussowitsch     PetscCall(VecView(F, NULL));
31741e469ddSJunchao Zhang   }
31841e469ddSJunchao Zhang 
31941e469ddSJunchao Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
32041e469ddSJunchao Zhang      Evaluate initial guess; then solve nonlinear system
32141e469ddSJunchao Zhang    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
32241e469ddSJunchao Zhang 
32341e469ddSJunchao Zhang   /*
32441e469ddSJunchao Zhang      Note: The user should initialize the vector, x, with the initial guess
32541e469ddSJunchao Zhang      for the nonlinear solver prior to calling SNESSolve().  In particular,
32641e469ddSJunchao Zhang      to employ an initial guess of zero, the user should explicitly set
32741e469ddSJunchao Zhang      this vector to zero by calling VecSet().
32841e469ddSJunchao Zhang   */
3299566063dSJacob Faibussowitsch   PetscCall(FormInitialGuess(x));
3309566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, x));
3319566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &its));
33263a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %" PetscInt_FMT "\n", its));
33341e469ddSJunchao Zhang 
33441e469ddSJunchao Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
33541e469ddSJunchao Zhang      Check solution and clean up
33641e469ddSJunchao Zhang    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
33741e469ddSJunchao Zhang   /*
33841e469ddSJunchao Zhang      Check the error
33941e469ddSJunchao Zhang   */
3409566063dSJacob Faibussowitsch   PetscCall(VecAXPY(x, none, U));
3419566063dSJacob Faibussowitsch   PetscCall(VecNorm(x, NORM_2, &norm));
34263a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Norm of error %g Iterations %" PetscInt_FMT "\n", (double)norm, its));
34341e469ddSJunchao Zhang 
34441e469ddSJunchao Zhang   /*
34541e469ddSJunchao Zhang      Free work space.  All PETSc objects should be destroyed when they
34641e469ddSJunchao Zhang      are no longer needed.
34741e469ddSJunchao Zhang   */
3489566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
3499566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
3509566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
3519566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&F));
3529566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
3539566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
3549566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&ctx.da));
3559566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
356b122ec5aSJacob Faibussowitsch   return 0;
35741e469ddSJunchao Zhang }
358c0558f20SBarry Smith 
359c0558f20SBarry Smith /*TEST
360c0558f20SBarry Smith 
361c0558f20SBarry Smith    build:
36241e469ddSJunchao Zhang      requires: kokkos_kernels
363c0558f20SBarry Smith 
364c0558f20SBarry Smith    test:
36587699089SJunchao Zhang      requires: kokkos_kernels !complex !single
366c0558f20SBarry Smith      nsize: 2
36741e469ddSJunchao Zhang      args: -dm_vec_type kokkos -dm_mat_type aijkokkos -view_initial -snes_monitor
368c0558f20SBarry Smith      output_file: output/ex3k_1.out
369c0558f20SBarry Smith 
370c0558f20SBarry Smith TEST*/
371