xref: /petsc/src/snes/tutorials/ex3k.kokkos.cxx (revision 41e469dd2095110404231de8e03e0545f30cf604)
1c0558f20SBarry Smith static char help[] = "Newton methods to solve u'' + u^{2} = f in parallel. Uses Kokkos\n\\n";
2c0558f20SBarry Smith 
3*41e469ddSJunchao 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 /*
34*41e469ddSJunchao 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:
42*41e469ddSJunchao 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 */
48*41e469ddSJunchao 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;
52*41e469ddSJunchao Zhang   PetscScalar    *X,*R,*F,d;
53c0558f20SBarry Smith   PetscErrorCode ierr;
54c0558f20SBarry Smith   PetscInt       i,M,xs,xm;
55*41e469ddSJunchao Zhang   Vec            xl;
56c0558f20SBarry Smith 
57c0558f20SBarry Smith   PetscFunctionBeginUser;
58*41e469ddSJunchao Zhang   ierr = DMGetLocalVector(da,&xl);CHKERRQ(ierr);
59*41e469ddSJunchao Zhang   ierr = DMGlobalToLocal(da,x,INSERT_VALUES,xl);CHKERRQ(ierr);
60*41e469ddSJunchao Zhang   ierr = DMDAVecGetArray(da,xl,&X);CHKERRQ(ierr);
61*41e469ddSJunchao Zhang   ierr = DMDAVecGetArray(da,r,&R);CHKERRQ(ierr);
62*41e469ddSJunchao 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 */
68*41e469ddSJunchao Zhang     R[0] = X[0];
69c0558f20SBarry Smith     xs++;xm--;
70c0558f20SBarry Smith   }
71c0558f20SBarry Smith   if (xs+xm == M) {  /* right boundary */
72*41e469ddSJunchao 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);
76*41e469ddSJunchao 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 
78*41e469ddSJunchao Zhang   ierr = DMDAVecRestoreArray(da,xl,&X);CHKERRQ(ierr);
79*41e469ddSJunchao Zhang   ierr = DMDAVecRestoreArray(da,r,&R);CHKERRQ(ierr);
80*41e469ddSJunchao Zhang   ierr = DMDAVecRestoreArray(da,user->F,&F);CHKERRQ(ierr);
81*41e469ddSJunchao Zhang   ierr = DMRestoreLocalVector(da,&xl);CHKERRQ(ierr);
82c0558f20SBarry Smith   PetscFunctionReturn(0);
83c0558f20SBarry Smith }
84c0558f20SBarry Smith 
85*41e469ddSJunchao Zhang using DefaultExecutionSpace             = Kokkos::DefaultExecutionSpace;
86*41e469ddSJunchao Zhang using DefaultMemorySpace                = Kokkos::DefaultExecutionSpace::memory_space;
87*41e469ddSJunchao Zhang using PetscScalarKokkosOffsetView       = Kokkos::Experimental::OffsetView<PetscScalar*,DefaultMemorySpace>;
88*41e469ddSJunchao Zhang using ConstPetscScalarKokkosOffsetView  = Kokkos::Experimental::OffsetView<const PetscScalar*,DefaultMemorySpace>;
89*41e469ddSJunchao Zhang 
90*41e469ddSJunchao Zhang PetscErrorCode KokkosFunction(SNES snes,Vec x,Vec r,void *ctx)
91*41e469ddSJunchao Zhang {
92*41e469ddSJunchao Zhang   PetscErrorCode                       ierr;
93*41e469ddSJunchao Zhang   ApplicationCtx                       *user = (ApplicationCtx*) ctx;
94*41e469ddSJunchao Zhang   DM                                   da = user->da;
95*41e469ddSJunchao Zhang   PetscScalar                          d;
96*41e469ddSJunchao Zhang   PetscInt                             M;
97*41e469ddSJunchao Zhang   Vec                                  xl;
98*41e469ddSJunchao Zhang   PetscScalarKokkosOffsetView          R;
99*41e469ddSJunchao Zhang   ConstPetscScalarKokkosOffsetView     X,F;
100*41e469ddSJunchao Zhang 
101*41e469ddSJunchao Zhang   PetscFunctionBeginUser;
102*41e469ddSJunchao Zhang   ierr = DMGetLocalVector(da,&xl);CHKERRQ(ierr);
103*41e469ddSJunchao Zhang   ierr = DMGlobalToLocal(da,x,INSERT_VALUES,xl);CHKERRQ(ierr);
104*41e469ddSJunchao Zhang   d    = 1.0/(user->h*user->h);
105*41e469ddSJunchao Zhang   ierr = DMDAGetInfo(da,NULL,&M,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
106*41e469ddSJunchao Zhang   ierr = DMDAVecGetKokkosOffsetView(da,xl,&X);CHKERRQ(ierr); /* read only */
107*41e469ddSJunchao Zhang   ierr = DMDAVecGetKokkosOffsetViewWrite(da,r,&R);CHKERRQ(ierr); /* write only */
108*41e469ddSJunchao Zhang   ierr = DMDAVecGetKokkosOffsetView(da,user->F,&F);CHKERRQ(ierr); /* read only */
109*41e469ddSJunchao Zhang   Kokkos:: parallel_for (Kokkos::RangePolicy<>(R.begin(0),R.end(0)),KOKKOS_LAMBDA (int i) {
110*41e469ddSJunchao Zhang     if (i == 0)        R(0) = X(0);        /* left boundary */
111*41e469ddSJunchao Zhang     else if (i == M-1) R(i) = X(i) - 1.0;  /* right boundary */
112*41e469ddSJunchao Zhang     else               R(i) = d*(X(i-1) - 2.0*X(i) + X(i+1)) + X(i)*X(i) - F(i); /* interior */
113*41e469ddSJunchao Zhang   });
114*41e469ddSJunchao Zhang   ierr = DMDAVecRestoreKokkosOffsetView(da,xl,&X);CHKERRQ(ierr);
115*41e469ddSJunchao Zhang   ierr = DMDAVecRestoreKokkosOffsetViewWrite(da,r,&R);CHKERRQ(ierr);
116*41e469ddSJunchao Zhang   ierr = DMDAVecRestoreKokkosOffsetView(da,user->F,&F);CHKERRQ(ierr);
117*41e469ddSJunchao Zhang   ierr = DMRestoreLocalVector(da,&xl);CHKERRQ(ierr);
118*41e469ddSJunchao Zhang   PetscFunctionReturn(0);
119*41e469ddSJunchao Zhang }
120*41e469ddSJunchao Zhang 
121*41e469ddSJunchao Zhang PetscErrorCode StubFunction(SNES snes ,Vec x,Vec r,void *ctx)
122*41e469ddSJunchao Zhang {
123*41e469ddSJunchao Zhang   PetscErrorCode                       ierr;
124*41e469ddSJunchao Zhang   ApplicationCtx                       *user = (ApplicationCtx*) ctx;
125*41e469ddSJunchao Zhang   DM                                   da = user->da;
126*41e469ddSJunchao Zhang   Vec                                  rk;
127*41e469ddSJunchao Zhang   PetscReal                            norm=0;
128*41e469ddSJunchao Zhang 
129*41e469ddSJunchao Zhang   PetscFunctionBeginUser;
130*41e469ddSJunchao Zhang   ierr = DMGetGlobalVector(da,&rk);CHKERRQ(ierr);
131*41e469ddSJunchao Zhang   ierr = CpuFunction(snes,x,r,ctx);CHKERRQ(ierr);
132*41e469ddSJunchao Zhang   ierr = KokkosFunction(snes,x,rk,ctx);CHKERRQ(ierr);
133*41e469ddSJunchao Zhang   ierr = VecAXPY(rk,-1.0,r);CHKERRQ(ierr);
134*41e469ddSJunchao Zhang   ierr = VecNorm(rk,NORM_2,&norm);CHKERRQ(ierr);
135*41e469ddSJunchao Zhang   ierr = DMRestoreGlobalVector(da,&rk);CHKERRQ(ierr);
136*41e469ddSJunchao Zhang   if (norm > 1e-6) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"KokkosFunction() different from CpuFunction() with a diff norm = %g\n",norm);
137*41e469ddSJunchao Zhang   PetscFunctionReturn(0);
138*41e469ddSJunchao 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 
219*41e469ddSJunchao Zhang int main(int argc,char **argv)
220*41e469ddSJunchao Zhang {
221*41e469ddSJunchao Zhang   SNES                        snes;                 /* SNES context */
222*41e469ddSJunchao Zhang   Mat                         J;                    /* Jacobian matrix */
223*41e469ddSJunchao Zhang   ApplicationCtx              ctx;                  /* user-defined context */
224*41e469ddSJunchao Zhang   Vec                         x,r,U,F;              /* vectors */
225*41e469ddSJunchao Zhang   PetscScalar                 none = -1.0;
226*41e469ddSJunchao Zhang   PetscErrorCode              ierr;
227*41e469ddSJunchao Zhang   PetscInt                    its,N = 5,maxit,maxf;
228*41e469ddSJunchao Zhang   PetscReal                   abstol,rtol,stol,norm;
229*41e469ddSJunchao Zhang   PetscBool                   viewinitial = PETSC_FALSE;
230*41e469ddSJunchao Zhang   PetscScalarKokkosOffsetView FF,UU;
231*41e469ddSJunchao Zhang 
232*41e469ddSJunchao Zhang   ierr  = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
233*41e469ddSJunchao Zhang   ierr  = PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL);CHKERRQ(ierr);
234*41e469ddSJunchao Zhang   ctx.h = 1.0/(N-1);
235*41e469ddSJunchao Zhang   ierr  = PetscOptionsGetBool(NULL,NULL,"-view_initial",&viewinitial,NULL);CHKERRQ(ierr);
236*41e469ddSJunchao Zhang 
237*41e469ddSJunchao Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
238*41e469ddSJunchao Zhang      Create nonlinear solver context
239*41e469ddSJunchao Zhang      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
240*41e469ddSJunchao Zhang   ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr);
241*41e469ddSJunchao Zhang 
242c0558f20SBarry Smith   /*
243*41e469ddSJunchao Zhang      Create distributed array (DMDA) to manage parallel grid and vectors
244c0558f20SBarry Smith   */
245*41e469ddSJunchao Zhang   ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,N,1,1,NULL,&ctx.da);CHKERRQ(ierr);
246*41e469ddSJunchao Zhang   ierr = DMSetFromOptions(ctx.da);CHKERRQ(ierr);
247*41e469ddSJunchao Zhang   ierr = DMSetUp(ctx.da);CHKERRQ(ierr);
248*41e469ddSJunchao Zhang 
249*41e469ddSJunchao Zhang   /*
250*41e469ddSJunchao Zhang      Extract global and local vectors from DMDA; then duplicate for remaining
251*41e469ddSJunchao Zhang      vectors that are the same types
252*41e469ddSJunchao Zhang   */
253*41e469ddSJunchao Zhang   ierr = DMCreateGlobalVector(ctx.da,&x);CHKERRQ(ierr);
254*41e469ddSJunchao Zhang   ierr = PetscObjectSetName((PetscObject)x,"Approximate Solution");CHKERRQ(ierr);
255*41e469ddSJunchao Zhang   ierr = VecDuplicate(x,&r);CHKERRQ(ierr);
256*41e469ddSJunchao Zhang   ierr = VecDuplicate(x,&F);CHKERRQ(ierr); ctx.F = F;
257*41e469ddSJunchao Zhang   ierr = PetscObjectSetName((PetscObject)F,"Forcing function");CHKERRQ(ierr);
258*41e469ddSJunchao Zhang   ierr = VecDuplicate(x,&U);CHKERRQ(ierr);
259*41e469ddSJunchao Zhang   ierr = PetscObjectSetName((PetscObject)U,"Exact Solution");CHKERRQ(ierr);
260*41e469ddSJunchao Zhang 
261*41e469ddSJunchao Zhang   /*
262*41e469ddSJunchao Zhang      Set function evaluation routine and vector.  Whenever the nonlinear
263*41e469ddSJunchao Zhang      solver needs to compute the nonlinear function, it will call this
264*41e469ddSJunchao Zhang      routine.
265*41e469ddSJunchao Zhang       - Note that the final routine argument is the user-defined
266*41e469ddSJunchao Zhang         context that provides application-specific data for the
267*41e469ddSJunchao Zhang         function evaluation routine.
268*41e469ddSJunchao Zhang 
269*41e469ddSJunchao Zhang      At the begining, one can use a stub function that checks the Kokkos version
270*41e469ddSJunchao Zhang      against the CPU version to quickly expose errors.
271*41e469ddSJunchao Zhang      ierr = SNESSetFunction(snes,r,StubFunction,&ctx);CHKERRQ(ierr);
272*41e469ddSJunchao Zhang   */
273*41e469ddSJunchao Zhang   ierr = SNESSetFunction(snes,r,KokkosFunction,&ctx);CHKERRQ(ierr);
274*41e469ddSJunchao Zhang 
275*41e469ddSJunchao Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
276*41e469ddSJunchao Zhang      Create matrix data structure; set Jacobian evaluation routine
277*41e469ddSJunchao Zhang      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
278*41e469ddSJunchao Zhang   ierr = DMCreateMatrix(ctx.da,&J);CHKERRQ(ierr);
279*41e469ddSJunchao Zhang 
280*41e469ddSJunchao Zhang   /*
281*41e469ddSJunchao Zhang      Set Jacobian matrix data structure and default Jacobian evaluation
282*41e469ddSJunchao Zhang      routine.  Whenever the nonlinear solver needs to compute the
283*41e469ddSJunchao Zhang      Jacobian matrix, it will call this routine.
284*41e469ddSJunchao Zhang       - Note that the final routine argument is the user-defined
285*41e469ddSJunchao Zhang         context that provides application-specific data for the
286*41e469ddSJunchao Zhang         Jacobian evaluation routine.
287*41e469ddSJunchao Zhang   */
288*41e469ddSJunchao Zhang   ierr = SNESSetJacobian(snes,J,J,FormJacobian,&ctx);CHKERRQ(ierr);
289*41e469ddSJunchao Zhang   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
290*41e469ddSJunchao Zhang   ierr = SNESGetTolerances(snes,&abstol,&rtol,&stol,&maxit,&maxf);CHKERRQ(ierr);
291*41e469ddSJunchao 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);
292*41e469ddSJunchao Zhang 
293*41e469ddSJunchao Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
294*41e469ddSJunchao Zhang      Initialize application:
295*41e469ddSJunchao Zhang      Store forcing function of PDE and exact solution
296*41e469ddSJunchao Zhang    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
297*41e469ddSJunchao Zhang   ierr = DMDAVecGetKokkosOffsetViewWrite(ctx.da,F,&FF);CHKERRQ(ierr);
298*41e469ddSJunchao Zhang   ierr = DMDAVecGetKokkosOffsetViewWrite(ctx.da,U,&UU);CHKERRQ(ierr);
299*41e469ddSJunchao Zhang   Kokkos:: parallel_for (Kokkos::RangePolicy<>(FF.begin(0),FF.end(0)),KOKKOS_LAMBDA (int i) {
300*41e469ddSJunchao Zhang     PetscReal xp = i*ctx.h;
301*41e469ddSJunchao Zhang     FF(i) = 6.0*xp + pow(xp+1.e-12,6.0); /* +1.e-12 is to prevent 0^6 */
302*41e469ddSJunchao Zhang     UU(i) = xp*xp*xp;
303*41e469ddSJunchao Zhang   });
304*41e469ddSJunchao Zhang   ierr = DMDAVecRestoreKokkosOffsetViewWrite(ctx.da,F,&FF);CHKERRQ(ierr);
305*41e469ddSJunchao Zhang   ierr = DMDAVecRestoreKokkosOffsetViewWrite(ctx.da,U,&UU);CHKERRQ(ierr);
306*41e469ddSJunchao Zhang 
307*41e469ddSJunchao Zhang   if (viewinitial) {
308*41e469ddSJunchao Zhang     ierr = VecView(U,NULL);CHKERRQ(ierr);
309*41e469ddSJunchao Zhang     ierr = VecView(F,NULL);CHKERRQ(ierr);
310*41e469ddSJunchao Zhang   }
311*41e469ddSJunchao Zhang 
312*41e469ddSJunchao Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
313*41e469ddSJunchao Zhang      Evaluate initial guess; then solve nonlinear system
314*41e469ddSJunchao Zhang    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
315*41e469ddSJunchao Zhang 
316*41e469ddSJunchao Zhang   /*
317*41e469ddSJunchao Zhang      Note: The user should initialize the vector, x, with the initial guess
318*41e469ddSJunchao Zhang      for the nonlinear solver prior to calling SNESSolve().  In particular,
319*41e469ddSJunchao Zhang      to employ an initial guess of zero, the user should explicitly set
320*41e469ddSJunchao Zhang      this vector to zero by calling VecSet().
321*41e469ddSJunchao Zhang   */
322*41e469ddSJunchao Zhang   ierr = FormInitialGuess(x);CHKERRQ(ierr);
323*41e469ddSJunchao Zhang   ierr = SNESSolve(snes,NULL,x);CHKERRQ(ierr);
324*41e469ddSJunchao Zhang   ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
325*41e469ddSJunchao Zhang   ierr = PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %D\n",its);CHKERRQ(ierr);
326*41e469ddSJunchao Zhang 
327*41e469ddSJunchao Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
328*41e469ddSJunchao Zhang      Check solution and clean up
329*41e469ddSJunchao Zhang    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
330*41e469ddSJunchao Zhang   /*
331*41e469ddSJunchao Zhang      Check the error
332*41e469ddSJunchao Zhang   */
333*41e469ddSJunchao Zhang   ierr = VecAXPY(x,none,U);CHKERRQ(ierr);
334*41e469ddSJunchao Zhang   ierr = VecNorm(x,NORM_2,&norm);CHKERRQ(ierr);
335*41e469ddSJunchao Zhang   ierr = PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g Iterations %D\n",(double)norm,its);CHKERRQ(ierr);
336*41e469ddSJunchao Zhang 
337*41e469ddSJunchao Zhang   /*
338*41e469ddSJunchao Zhang      Free work space.  All PETSc objects should be destroyed when they
339*41e469ddSJunchao Zhang      are no longer needed.
340*41e469ddSJunchao Zhang   */
341*41e469ddSJunchao Zhang   ierr = VecDestroy(&x);CHKERRQ(ierr);
342*41e469ddSJunchao Zhang   ierr = VecDestroy(&r);CHKERRQ(ierr);
343*41e469ddSJunchao Zhang   ierr = VecDestroy(&U);CHKERRQ(ierr);
344*41e469ddSJunchao Zhang   ierr = VecDestroy(&F);CHKERRQ(ierr);
345*41e469ddSJunchao Zhang   ierr = MatDestroy(&J);CHKERRQ(ierr);
346*41e469ddSJunchao Zhang   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
347*41e469ddSJunchao Zhang   ierr = DMDestroy(&ctx.da);CHKERRQ(ierr);
348*41e469ddSJunchao Zhang   ierr = PetscFinalize();
349*41e469ddSJunchao Zhang   return ierr;
350*41e469ddSJunchao Zhang }
351c0558f20SBarry Smith 
352c0558f20SBarry Smith /*TEST
353c0558f20SBarry Smith 
354c0558f20SBarry Smith    build:
355*41e469ddSJunchao Zhang      requires: kokkos_kernels
356c0558f20SBarry Smith 
357c0558f20SBarry Smith    test:
358*41e469ddSJunchao Zhang      requires: kokkos_kernels !complex !single cuda
359c0558f20SBarry Smith      nsize: 2
360*41e469ddSJunchao Zhang      args: -dm_vec_type kokkos -dm_mat_type aijkokkos -view_initial -snes_monitor
361c0558f20SBarry Smith      output_file: output/ex3k_1.out
362c0558f20SBarry Smith 
363c0558f20SBarry Smith TEST*/