Lines Matching +full:test +full:- +full:experimental
7 User-defined application context
11 Vec F; /* right-hand side of PDE */
15 /* ------------------------------------------------------------------- */
17 FormInitialGuess - Computes initial guess.
20 . x - the solution vector
31 /* ------------------------------------------------------------------- */
33 CpuFunction - Evaluates nonlinear function, F(x) on CPU
36 . snes - the SNES context
37 . x - input vector
38 . ctx - optional user-defined context, as set by SNESSetFunction()
41 . r - function vector
44 The user-defined context can contain any application-specific
50 DM da = user->da; in CpuFunction()
60 PetscCall(DMDAVecGetArray(da, user->F, &F)); in CpuFunction()
68 xm--; in CpuFunction()
71 R[xs + xm - 1] = X[xs + xm - 1] - 1.0; in CpuFunction()
72 xm--; in CpuFunction()
74 d = 1.0 / (user->h * user->h); in CpuFunction()
75 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]; in CpuFunction()
79 PetscCall(DMDAVecRestoreArray(da, user->F, &F)); in CpuFunction()
86 using PetscScalarKokkosOffsetView = Kokkos::Experimental::OffsetView<PetscScalar *, DefaultMem…
87 using ConstPetscScalarKokkosOffsetView = Kokkos::Experimental::OffsetView<const PetscScalar *, Defa…
92 DM da = user->da; in KokkosFunction()
102 d = 1.0 / (user->h * user->h); in KokkosFunction()
106 PetscCall(DMDAVecGetKokkosOffsetView(da, user->F, &F)); /* read only */ in KokkosFunction()
110 else if (i == M - 1) R(i) = X(i) - 1.0; /* right boundary */ in KokkosFunction()
111 else R(i) = d * (X(i - 1) - 2.0 * X(i) + X(i + 1)) + X(i) * X(i) - F(i); /* interior */ in KokkosFunction()
115 PetscCall(DMDAVecRestoreKokkosOffsetView(da, user->F, &F)); in KokkosFunction()
123 DM da = user->da; in StubFunction()
131 PetscCall(VecAXPY(rk, -1.0, r)); in StubFunction()
134 …PetscCheck(norm <= 1e-6, PETSC_COMM_SELF, PETSC_ERR_PLIB, "KokkosFunction() different from CpuFunc… in StubFunction()
137 /* ------------------------------------------------------------------- */
139 FormJacobian - Evaluates Jacobian matrix.
142 . snes - the SNES context
143 . x - input vector
144 . dummy - optional user-defined context (not used here)
147 . jac - Jacobian matrix
148 . B - optionally different matrix used to construct the preconditioner
156 DM da = user->da; in FormJacobian()
181 xm--; in FormJacobian()
184 i = M - 1; in FormJacobian()
187 xm--; in FormJacobian()
192 - Note that in this case we set all elements for a particular in FormJacobian()
195 d = 1.0 / (user->h * user->h); in FormJacobian()
197 j[0] = i - 1; in FormJacobian()
201 A[1] = -2.0 * d + 2.0 * xx[i]; in FormJacobian()
206 Assemble matrix, using the 2-step process: in FormJacobian()
224 ApplicationCtx ctx; /* user-defined context */ in main()
226 PetscScalar none = -1.0; in main()
233 PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &N, NULL)); in main()
234 ctx.h = 1.0 / (N - 1); in main()
235 PetscCall(PetscOptionsGetBool(NULL, NULL, "-view_initial", &viewinitial, NULL)); in main()
237 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - in main()
239 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ in main()
266 - Note that the final routine argument is the user-defined in main()
267 context that provides application-specific data for the in main()
276 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - in main()
278 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ in main()
285 - Note that the final routine argument is the user-defined in main()
286 context that provides application-specific data for the in main()
294 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - in main()
297 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ in main()
305 FF(i) = 6.0 * xp + pow(xp + 1.e-12, 6.0); /* +1.e-12 is to prevent 0^6 */ in main()
317 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - in main()
319 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ in main()
332 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - in main()
334 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ in main()
357 /*TEST
362 test:
365 args: -dm_vec_type kokkos -dm_mat_type aijkokkos -view_initial -snes_monitor
368 TEST*/