xref: /petsc/src/snes/tutorials/ex3k.kokkos.cxx (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
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   PetscErrorCode ierr;
25c0558f20SBarry Smith   PetscScalar    pfive = .50;
26c0558f20SBarry Smith 
27c0558f20SBarry Smith   PetscFunctionBeginUser;
28c0558f20SBarry Smith   ierr = VecSet(x,pfive);CHKERRQ(ierr);
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   PetscErrorCode ierr;
54c0558f20SBarry Smith   PetscInt       i,M,xs,xm;
5541e469ddSJunchao Zhang   Vec            xl;
56c0558f20SBarry Smith 
57c0558f20SBarry Smith   PetscFunctionBeginUser;
5841e469ddSJunchao Zhang   ierr = DMGetLocalVector(da,&xl);CHKERRQ(ierr);
5941e469ddSJunchao Zhang   ierr = DMGlobalToLocal(da,x,INSERT_VALUES,xl);CHKERRQ(ierr);
6041e469ddSJunchao Zhang   ierr = DMDAVecGetArray(da,xl,&X);CHKERRQ(ierr);
6141e469ddSJunchao Zhang   ierr = DMDAVecGetArray(da,r,&R);CHKERRQ(ierr);
6241e469ddSJunchao Zhang   ierr = DMDAVecGetArray(da,user->F,&F);CHKERRQ(ierr);
63c0558f20SBarry Smith 
64c0558f20SBarry Smith   ierr = DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr);
65c0558f20SBarry Smith   ierr = DMDAGetInfo(da,NULL,&M,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
66c0558f20SBarry Smith 
67c0558f20SBarry Smith   if (xs == 0) { /* left boundary */
6841e469ddSJunchao Zhang     R[0] = X[0];
69c0558f20SBarry Smith     xs++;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 
7841e469ddSJunchao Zhang   ierr = DMDAVecRestoreArray(da,xl,&X);CHKERRQ(ierr);
7941e469ddSJunchao Zhang   ierr = DMDAVecRestoreArray(da,r,&R);CHKERRQ(ierr);
8041e469ddSJunchao Zhang   ierr = DMDAVecRestoreArray(da,user->F,&F);CHKERRQ(ierr);
8141e469ddSJunchao Zhang   ierr = DMRestoreLocalVector(da,&xl);CHKERRQ(ierr);
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   PetscErrorCode                       ierr;
9341e469ddSJunchao Zhang   ApplicationCtx                       *user = (ApplicationCtx*) ctx;
9441e469ddSJunchao Zhang   DM                                   da = user->da;
9541e469ddSJunchao Zhang   PetscScalar                          d;
9641e469ddSJunchao Zhang   PetscInt                             M;
9741e469ddSJunchao Zhang   Vec                                  xl;
9841e469ddSJunchao Zhang   PetscScalarKokkosOffsetView          R;
9941e469ddSJunchao Zhang   ConstPetscScalarKokkosOffsetView     X,F;
10041e469ddSJunchao Zhang 
10141e469ddSJunchao Zhang   PetscFunctionBeginUser;
10241e469ddSJunchao Zhang   ierr = DMGetLocalVector(da,&xl);CHKERRQ(ierr);
10341e469ddSJunchao Zhang   ierr = DMGlobalToLocal(da,x,INSERT_VALUES,xl);CHKERRQ(ierr);
10441e469ddSJunchao Zhang   d    = 1.0/(user->h*user->h);
10541e469ddSJunchao Zhang   ierr = DMDAGetInfo(da,NULL,&M,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
10641e469ddSJunchao Zhang   ierr = DMDAVecGetKokkosOffsetView(da,xl,&X);CHKERRQ(ierr); /* read only */
10741e469ddSJunchao Zhang   ierr = DMDAVecGetKokkosOffsetViewWrite(da,r,&R);CHKERRQ(ierr); /* write only */
10841e469ddSJunchao Zhang   ierr = DMDAVecGetKokkosOffsetView(da,user->F,&F);CHKERRQ(ierr); /* read only */
10941e469ddSJunchao Zhang   Kokkos:: parallel_for (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   });
11441e469ddSJunchao Zhang   ierr = DMDAVecRestoreKokkosOffsetView(da,xl,&X);CHKERRQ(ierr);
11541e469ddSJunchao Zhang   ierr = DMDAVecRestoreKokkosOffsetViewWrite(da,r,&R);CHKERRQ(ierr);
11641e469ddSJunchao Zhang   ierr = DMDAVecRestoreKokkosOffsetView(da,user->F,&F);CHKERRQ(ierr);
11741e469ddSJunchao Zhang   ierr = DMRestoreLocalVector(da,&xl);CHKERRQ(ierr);
11841e469ddSJunchao Zhang   PetscFunctionReturn(0);
11941e469ddSJunchao Zhang }
12041e469ddSJunchao Zhang 
12141e469ddSJunchao Zhang PetscErrorCode StubFunction(SNES snes ,Vec x,Vec r,void *ctx)
12241e469ddSJunchao Zhang {
12341e469ddSJunchao Zhang   PetscErrorCode                       ierr;
12441e469ddSJunchao Zhang   ApplicationCtx                       *user = (ApplicationCtx*) ctx;
12541e469ddSJunchao Zhang   DM                                   da = user->da;
12641e469ddSJunchao Zhang   Vec                                  rk;
12741e469ddSJunchao Zhang   PetscReal                            norm=0;
12841e469ddSJunchao Zhang 
12941e469ddSJunchao Zhang   PetscFunctionBeginUser;
13041e469ddSJunchao Zhang   ierr = DMGetGlobalVector(da,&rk);CHKERRQ(ierr);
13141e469ddSJunchao Zhang   ierr = CpuFunction(snes,x,r,ctx);CHKERRQ(ierr);
13241e469ddSJunchao Zhang   ierr = KokkosFunction(snes,x,rk,ctx);CHKERRQ(ierr);
13341e469ddSJunchao Zhang   ierr = VecAXPY(rk,-1.0,r);CHKERRQ(ierr);
13441e469ddSJunchao Zhang   ierr = VecNorm(rk,NORM_2,&norm);CHKERRQ(ierr);
13541e469ddSJunchao Zhang   ierr = DMRestoreGlobalVector(da,&rk);CHKERRQ(ierr);
136*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(norm > 1e-6,PETSC_COMM_SELF,PETSC_ERR_PLIB,"KokkosFunction() different from CpuFunction() with a diff norm = %g",norm);
13741e469ddSJunchao Zhang   PetscFunctionReturn(0);
13841e469ddSJunchao Zhang }
139c0558f20SBarry Smith /* ------------------------------------------------------------------- */
140c0558f20SBarry Smith /*
141c0558f20SBarry Smith    FormJacobian - Evaluates Jacobian matrix.
142c0558f20SBarry Smith 
143c0558f20SBarry Smith    Input Parameters:
144c0558f20SBarry Smith .  snes - the SNES context
145c0558f20SBarry Smith .  x - input vector
146c0558f20SBarry Smith .  dummy - optional user-defined context (not used here)
147c0558f20SBarry Smith 
148c0558f20SBarry Smith    Output Parameters:
149c0558f20SBarry Smith .  jac - Jacobian matrix
150c0558f20SBarry Smith .  B - optionally different preconditioning matrix
151c0558f20SBarry Smith .  flag - flag indicating matrix structure
152c0558f20SBarry Smith */
153c0558f20SBarry Smith PetscErrorCode FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void *ctx)
154c0558f20SBarry Smith {
155c0558f20SBarry Smith   ApplicationCtx *user = (ApplicationCtx*) ctx;
156c0558f20SBarry Smith   PetscScalar    *xx,d,A[3];
157c0558f20SBarry Smith   PetscErrorCode ierr;
158c0558f20SBarry Smith   PetscInt       i,j[3],M,xs,xm;
159c0558f20SBarry Smith   DM             da = user->da;
160c0558f20SBarry Smith 
161c0558f20SBarry Smith   PetscFunctionBeginUser;
162c0558f20SBarry Smith   /*
163c0558f20SBarry Smith      Get pointer to vector data
164c0558f20SBarry Smith   */
165c0558f20SBarry Smith   ierr = DMDAVecGetArrayRead(da,x,&xx);CHKERRQ(ierr);
166c0558f20SBarry Smith   ierr = DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr);
167c0558f20SBarry Smith 
168c0558f20SBarry Smith   /*
169c0558f20SBarry Smith     Get range of locally owned matrix
170c0558f20SBarry Smith   */
171c0558f20SBarry Smith   ierr = DMDAGetInfo(da,NULL,&M,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
172c0558f20SBarry Smith 
173c0558f20SBarry Smith   /*
174c0558f20SBarry Smith      Determine starting and ending local indices for interior grid points.
175c0558f20SBarry Smith      Set Jacobian entries for boundary points.
176c0558f20SBarry Smith   */
177c0558f20SBarry Smith 
178c0558f20SBarry Smith   if (xs == 0) {  /* left boundary */
179c0558f20SBarry Smith     i = 0; A[0] = 1.0;
180c0558f20SBarry Smith 
181c0558f20SBarry Smith     ierr = MatSetValues(jac,1,&i,1,&i,A,INSERT_VALUES);CHKERRQ(ierr);
182c0558f20SBarry Smith     xs++;xm--;
183c0558f20SBarry Smith   }
184c0558f20SBarry Smith   if (xs+xm == M) { /* right boundary */
185c0558f20SBarry Smith     i    = M-1;
186c0558f20SBarry Smith     A[0] = 1.0;
187c0558f20SBarry Smith     ierr = MatSetValues(jac,1,&i,1,&i,A,INSERT_VALUES);CHKERRQ(ierr);
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++) {
198c0558f20SBarry Smith     j[0] = i - 1; j[1] = i; j[2] = i + 1;
199c0558f20SBarry Smith     A[0] = A[2] = d; A[1] = -2.0*d + 2.0*xx[i];
200c0558f20SBarry Smith     ierr = MatSetValues(jac,1,&i,3,j,A,INSERT_VALUES);CHKERRQ(ierr);
201c0558f20SBarry Smith   }
202c0558f20SBarry Smith 
203c0558f20SBarry Smith   /*
204c0558f20SBarry Smith      Assemble matrix, using the 2-step process:
205c0558f20SBarry Smith        MatAssemblyBegin(), MatAssemblyEnd().
206c0558f20SBarry Smith      By placing code between these two statements, computations can be
207c0558f20SBarry Smith      done while messages are in transition.
208c0558f20SBarry Smith 
209c0558f20SBarry Smith      Also, restore vector.
210c0558f20SBarry Smith   */
211c0558f20SBarry Smith 
212c0558f20SBarry Smith   ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
213c0558f20SBarry Smith   ierr = DMDAVecRestoreArrayRead(da,x,&xx);CHKERRQ(ierr);
214c0558f20SBarry Smith   ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
215c0558f20SBarry Smith 
216c0558f20SBarry Smith   PetscFunctionReturn(0);
217c0558f20SBarry Smith }
218c0558f20SBarry Smith 
21941e469ddSJunchao Zhang int main(int argc,char **argv)
22041e469ddSJunchao Zhang {
22141e469ddSJunchao Zhang   SNES                        snes;                 /* SNES context */
22241e469ddSJunchao Zhang   Mat                         J;                    /* Jacobian matrix */
22341e469ddSJunchao Zhang   ApplicationCtx              ctx;                  /* user-defined context */
22441e469ddSJunchao Zhang   Vec                         x,r,U,F;              /* vectors */
22541e469ddSJunchao Zhang   PetscScalar                 none = -1.0;
22641e469ddSJunchao Zhang   PetscErrorCode              ierr;
22741e469ddSJunchao Zhang   PetscInt                    its,N = 5,maxit,maxf;
22841e469ddSJunchao Zhang   PetscReal                   abstol,rtol,stol,norm;
22941e469ddSJunchao Zhang   PetscBool                   viewinitial = PETSC_FALSE;
23041e469ddSJunchao Zhang 
23141e469ddSJunchao Zhang   ierr  = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
23241e469ddSJunchao Zhang   ierr  = PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL);CHKERRQ(ierr);
23341e469ddSJunchao Zhang   ctx.h = 1.0/(N-1);
23441e469ddSJunchao Zhang   ierr  = PetscOptionsGetBool(NULL,NULL,"-view_initial",&viewinitial,NULL);CHKERRQ(ierr);
23541e469ddSJunchao Zhang 
23641e469ddSJunchao Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
23741e469ddSJunchao Zhang      Create nonlinear solver context
23841e469ddSJunchao Zhang      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
23941e469ddSJunchao Zhang   ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr);
24041e469ddSJunchao Zhang 
241c0558f20SBarry Smith   /*
24241e469ddSJunchao Zhang      Create distributed array (DMDA) to manage parallel grid and vectors
243c0558f20SBarry Smith   */
24441e469ddSJunchao Zhang   ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,N,1,1,NULL,&ctx.da);CHKERRQ(ierr);
24541e469ddSJunchao Zhang   ierr = DMSetFromOptions(ctx.da);CHKERRQ(ierr);
24641e469ddSJunchao Zhang   ierr = DMSetUp(ctx.da);CHKERRQ(ierr);
24741e469ddSJunchao Zhang 
24841e469ddSJunchao Zhang   /*
24941e469ddSJunchao Zhang      Extract global and local vectors from DMDA; then duplicate for remaining
25041e469ddSJunchao Zhang      vectors that are the same types
25141e469ddSJunchao Zhang   */
25241e469ddSJunchao Zhang   ierr = DMCreateGlobalVector(ctx.da,&x);CHKERRQ(ierr);
25341e469ddSJunchao Zhang   ierr = PetscObjectSetName((PetscObject)x,"Approximate Solution");CHKERRQ(ierr);
25441e469ddSJunchao Zhang   ierr = VecDuplicate(x,&r);CHKERRQ(ierr);
25541e469ddSJunchao Zhang   ierr = VecDuplicate(x,&F);CHKERRQ(ierr); ctx.F = F;
25641e469ddSJunchao Zhang   ierr = PetscObjectSetName((PetscObject)F,"Forcing function");CHKERRQ(ierr);
25741e469ddSJunchao Zhang   ierr = VecDuplicate(x,&U);CHKERRQ(ierr);
25841e469ddSJunchao Zhang   ierr = PetscObjectSetName((PetscObject)U,"Exact Solution");CHKERRQ(ierr);
25941e469ddSJunchao Zhang 
26041e469ddSJunchao Zhang   /*
26141e469ddSJunchao Zhang      Set function evaluation routine and vector.  Whenever the nonlinear
26241e469ddSJunchao Zhang      solver needs to compute the nonlinear function, it will call this
26341e469ddSJunchao Zhang      routine.
26441e469ddSJunchao Zhang       - Note that the final routine argument is the user-defined
26541e469ddSJunchao Zhang         context that provides application-specific data for the
26641e469ddSJunchao Zhang         function evaluation routine.
26741e469ddSJunchao Zhang 
268a5b23f4aSJose E. Roman      At the beginning, one can use a stub function that checks the Kokkos version
26941e469ddSJunchao Zhang      against the CPU version to quickly expose errors.
27041e469ddSJunchao Zhang      ierr = SNESSetFunction(snes,r,StubFunction,&ctx);CHKERRQ(ierr);
27141e469ddSJunchao Zhang   */
27241e469ddSJunchao Zhang   ierr = SNESSetFunction(snes,r,KokkosFunction,&ctx);CHKERRQ(ierr);
27341e469ddSJunchao Zhang 
27441e469ddSJunchao Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
27541e469ddSJunchao Zhang      Create matrix data structure; set Jacobian evaluation routine
27641e469ddSJunchao Zhang      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
27741e469ddSJunchao Zhang   ierr = DMCreateMatrix(ctx.da,&J);CHKERRQ(ierr);
27841e469ddSJunchao Zhang 
27941e469ddSJunchao Zhang   /*
28041e469ddSJunchao Zhang      Set Jacobian matrix data structure and default Jacobian evaluation
28141e469ddSJunchao Zhang      routine.  Whenever the nonlinear solver needs to compute the
28241e469ddSJunchao Zhang      Jacobian matrix, it will call this routine.
28341e469ddSJunchao Zhang       - Note that the final routine argument is the user-defined
28441e469ddSJunchao Zhang         context that provides application-specific data for the
28541e469ddSJunchao Zhang         Jacobian evaluation routine.
28641e469ddSJunchao Zhang   */
28741e469ddSJunchao Zhang   ierr = SNESSetJacobian(snes,J,J,FormJacobian,&ctx);CHKERRQ(ierr);
28841e469ddSJunchao Zhang   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
28941e469ddSJunchao Zhang   ierr = SNESGetTolerances(snes,&abstol,&rtol,&stol,&maxit,&maxf);CHKERRQ(ierr);
29041e469ddSJunchao Zhang   ierr = PetscPrintf(PETSC_COMM_WORLD,"atol=%g, rtol=%g, stol=%g, maxit=%D, maxf=%D\n",(double)abstol,(double)rtol,(double)stol,maxit,maxf);CHKERRQ(ierr);
29141e469ddSJunchao Zhang 
29241e469ddSJunchao Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
29341e469ddSJunchao Zhang      Initialize application:
29441e469ddSJunchao Zhang      Store forcing function of PDE and exact solution
29541e469ddSJunchao Zhang    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2969000cb80SJunchao Zhang   {
2979000cb80SJunchao Zhang     PetscScalarKokkosOffsetView FF,UU;
29841e469ddSJunchao Zhang     ierr = DMDAVecGetKokkosOffsetViewWrite(ctx.da,F,&FF);CHKERRQ(ierr);
29941e469ddSJunchao Zhang     ierr = DMDAVecGetKokkosOffsetViewWrite(ctx.da,U,&UU);CHKERRQ(ierr);
30041e469ddSJunchao Zhang     Kokkos:: parallel_for (Kokkos::RangePolicy<>(FF.begin(0),FF.end(0)),KOKKOS_LAMBDA (int i) {
30141e469ddSJunchao Zhang       PetscReal xp = i*ctx.h;
30241e469ddSJunchao Zhang       FF(i) = 6.0*xp + pow(xp+1.e-12,6.0); /* +1.e-12 is to prevent 0^6 */
30341e469ddSJunchao Zhang       UU(i) = xp*xp*xp;
30441e469ddSJunchao Zhang     });
30541e469ddSJunchao Zhang     ierr = DMDAVecRestoreKokkosOffsetViewWrite(ctx.da,F,&FF);CHKERRQ(ierr);
30641e469ddSJunchao Zhang     ierr = DMDAVecRestoreKokkosOffsetViewWrite(ctx.da,U,&UU);CHKERRQ(ierr);
3079000cb80SJunchao Zhang   }
30841e469ddSJunchao Zhang 
30941e469ddSJunchao Zhang   if (viewinitial) {
31041e469ddSJunchao Zhang     ierr = VecView(U,NULL);CHKERRQ(ierr);
31141e469ddSJunchao Zhang     ierr = VecView(F,NULL);CHKERRQ(ierr);
31241e469ddSJunchao Zhang   }
31341e469ddSJunchao Zhang 
31441e469ddSJunchao Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
31541e469ddSJunchao Zhang      Evaluate initial guess; then solve nonlinear system
31641e469ddSJunchao Zhang    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
31741e469ddSJunchao Zhang 
31841e469ddSJunchao Zhang   /*
31941e469ddSJunchao Zhang      Note: The user should initialize the vector, x, with the initial guess
32041e469ddSJunchao Zhang      for the nonlinear solver prior to calling SNESSolve().  In particular,
32141e469ddSJunchao Zhang      to employ an initial guess of zero, the user should explicitly set
32241e469ddSJunchao Zhang      this vector to zero by calling VecSet().
32341e469ddSJunchao Zhang   */
32441e469ddSJunchao Zhang   ierr = FormInitialGuess(x);CHKERRQ(ierr);
32541e469ddSJunchao Zhang   ierr = SNESSolve(snes,NULL,x);CHKERRQ(ierr);
32641e469ddSJunchao Zhang   ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
32741e469ddSJunchao Zhang   ierr = PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %D\n",its);CHKERRQ(ierr);
32841e469ddSJunchao Zhang 
32941e469ddSJunchao Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
33041e469ddSJunchao Zhang      Check solution and clean up
33141e469ddSJunchao Zhang    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
33241e469ddSJunchao Zhang   /*
33341e469ddSJunchao Zhang      Check the error
33441e469ddSJunchao Zhang   */
33541e469ddSJunchao Zhang   ierr = VecAXPY(x,none,U);CHKERRQ(ierr);
33641e469ddSJunchao Zhang   ierr = VecNorm(x,NORM_2,&norm);CHKERRQ(ierr);
33741e469ddSJunchao Zhang   ierr = PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g Iterations %D\n",(double)norm,its);CHKERRQ(ierr);
33841e469ddSJunchao Zhang 
33941e469ddSJunchao Zhang   /*
34041e469ddSJunchao Zhang      Free work space.  All PETSc objects should be destroyed when they
34141e469ddSJunchao Zhang      are no longer needed.
34241e469ddSJunchao Zhang   */
34341e469ddSJunchao Zhang   ierr = VecDestroy(&x);CHKERRQ(ierr);
34441e469ddSJunchao Zhang   ierr = VecDestroy(&r);CHKERRQ(ierr);
34541e469ddSJunchao Zhang   ierr = VecDestroy(&U);CHKERRQ(ierr);
34641e469ddSJunchao Zhang   ierr = VecDestroy(&F);CHKERRQ(ierr);
34741e469ddSJunchao Zhang   ierr = MatDestroy(&J);CHKERRQ(ierr);
34841e469ddSJunchao Zhang   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
34941e469ddSJunchao Zhang   ierr = DMDestroy(&ctx.da);CHKERRQ(ierr);
35041e469ddSJunchao Zhang   ierr = PetscFinalize();
35141e469ddSJunchao Zhang   return ierr;
35241e469ddSJunchao Zhang }
353c0558f20SBarry Smith 
354c0558f20SBarry Smith /*TEST
355c0558f20SBarry Smith 
356c0558f20SBarry Smith    build:
35741e469ddSJunchao Zhang      requires: kokkos_kernels
358c0558f20SBarry Smith 
359c0558f20SBarry Smith    test:
36087699089SJunchao Zhang      requires: kokkos_kernels !complex !single
361c0558f20SBarry Smith      nsize: 2
36241e469ddSJunchao Zhang      args: -dm_vec_type kokkos -dm_mat_type aijkokkos -view_initial -snes_monitor
363c0558f20SBarry Smith      output_file: output/ex3k_1.out
364c0558f20SBarry Smith 
365c0558f20SBarry Smith TEST*/
366