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 */ 99371c9d4SSatish Balay typedef struct { 10c0558f20SBarry Smith DM da; /* distributed array */ 11dd8e379bSPierre Jolivet 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 */ 22d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialGuess(Vec x) 23d71ae5a4SJacob Faibussowitsch { 24c0558f20SBarry Smith PetscScalar pfive = .50; 25c0558f20SBarry Smith 26c0558f20SBarry Smith PetscFunctionBeginUser; 279566063dSJacob Faibussowitsch PetscCall(VecSet(x, pfive)); 283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 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 */ 47d71ae5a4SJacob Faibussowitsch PetscErrorCode CpuFunction(SNES snes, Vec x, Vec r, void *ctx) 48d71ae5a4SJacob Faibussowitsch { 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; 569566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &xl)); 579566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocal(da, x, INSERT_VALUES, xl)); 589566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, xl, &X)); 599566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, r, &R)); 609566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, user->F, &F)); 61c0558f20SBarry Smith 629566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL)); 639566063dSJacob Faibussowitsch PetscCall(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]; 6711486bccSBarry Smith xs++; 6811486bccSBarry Smith xm--; 69c0558f20SBarry Smith } 70c0558f20SBarry Smith if (xs + xm == M) { /* right boundary */ 7141e469ddSJunchao Zhang R[xs + xm - 1] = X[xs + xm - 1] - 1.0; 72c0558f20SBarry Smith xm--; 73c0558f20SBarry Smith } 74c0558f20SBarry Smith d = 1.0 / (user->h * user->h); 7541e469ddSJunchao 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]; 76c0558f20SBarry Smith 779566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, xl, &X)); 789566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, r, &R)); 799566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, user->F, &F)); 809566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &xl)); 813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 82c0558f20SBarry Smith } 83c0558f20SBarry Smith 8441e469ddSJunchao Zhang using DefaultExecutionSpace = Kokkos::DefaultExecutionSpace; 8541e469ddSJunchao Zhang using DefaultMemorySpace = Kokkos::DefaultExecutionSpace::memory_space; 8641e469ddSJunchao Zhang using PetscScalarKokkosOffsetView = Kokkos::Experimental::OffsetView<PetscScalar *, DefaultMemorySpace>; 8741e469ddSJunchao Zhang using ConstPetscScalarKokkosOffsetView = Kokkos::Experimental::OffsetView<const PetscScalar *, DefaultMemorySpace>; 8841e469ddSJunchao Zhang 89d71ae5a4SJacob Faibussowitsch PetscErrorCode KokkosFunction(SNES snes, Vec x, Vec r, void *ctx) 90d71ae5a4SJacob Faibussowitsch { 9141e469ddSJunchao Zhang ApplicationCtx *user = (ApplicationCtx *)ctx; 9241e469ddSJunchao Zhang DM da = user->da; 9341e469ddSJunchao Zhang PetscScalar d; 9441e469ddSJunchao Zhang PetscInt M; 9541e469ddSJunchao Zhang Vec xl; 9641e469ddSJunchao Zhang PetscScalarKokkosOffsetView R; 9741e469ddSJunchao Zhang ConstPetscScalarKokkosOffsetView X, F; 9841e469ddSJunchao Zhang 9941e469ddSJunchao Zhang PetscFunctionBeginUser; 1009566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &xl)); 1019566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocal(da, x, INSERT_VALUES, xl)); 10241e469ddSJunchao Zhang d = 1.0 / (user->h * user->h); 1039566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, NULL, &M, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 1049566063dSJacob Faibussowitsch PetscCall(DMDAVecGetKokkosOffsetView(da, xl, &X)); /* read only */ 1059566063dSJacob Faibussowitsch PetscCall(DMDAVecGetKokkosOffsetViewWrite(da, r, &R)); /* write only */ 1069566063dSJacob Faibussowitsch PetscCall(DMDAVecGetKokkosOffsetView(da, user->F, &F)); /* read only */ 10711486bccSBarry Smith Kokkos::parallel_for( 10811486bccSBarry Smith Kokkos::RangePolicy<>(R.begin(0), R.end(0)), KOKKOS_LAMBDA(int i) { 10941e469ddSJunchao Zhang if (i == 0) R(0) = X(0); /* left boundary */ 11041e469ddSJunchao Zhang else if (i == M - 1) R(i) = X(i) - 1.0; /* right boundary */ 11141e469ddSJunchao Zhang else R(i) = d * (X(i - 1) - 2.0 * X(i) + X(i + 1)) + X(i) * X(i) - F(i); /* interior */ 11241e469ddSJunchao Zhang }); 1139566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreKokkosOffsetView(da, xl, &X)); 1149566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreKokkosOffsetViewWrite(da, r, &R)); 1159566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreKokkosOffsetView(da, user->F, &F)); 1169566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &xl)); 1173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11841e469ddSJunchao Zhang } 11941e469ddSJunchao Zhang 120d71ae5a4SJacob Faibussowitsch PetscErrorCode StubFunction(SNES snes, Vec x, Vec r, void *ctx) 121d71ae5a4SJacob Faibussowitsch { 12241e469ddSJunchao Zhang ApplicationCtx *user = (ApplicationCtx *)ctx; 12341e469ddSJunchao Zhang DM da = user->da; 12441e469ddSJunchao Zhang Vec rk; 12541e469ddSJunchao Zhang PetscReal norm = 0; 12641e469ddSJunchao Zhang 12741e469ddSJunchao Zhang PetscFunctionBeginUser; 1289566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(da, &rk)); 1299566063dSJacob Faibussowitsch PetscCall(CpuFunction(snes, x, r, ctx)); 1309566063dSJacob Faibussowitsch PetscCall(KokkosFunction(snes, x, rk, ctx)); 1319566063dSJacob Faibussowitsch PetscCall(VecAXPY(rk, -1.0, r)); 1329566063dSJacob Faibussowitsch PetscCall(VecNorm(rk, NORM_2, &norm)); 1339566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(da, &rk)); 13463a3b9bcSJacob Faibussowitsch PetscCheck(norm <= 1e-6, PETSC_COMM_SELF, PETSC_ERR_PLIB, "KokkosFunction() different from CpuFunction() with a diff norm = %g", (double)norm); 1353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13641e469ddSJunchao Zhang } 137c0558f20SBarry Smith /* ------------------------------------------------------------------- */ 138c0558f20SBarry Smith /* 139c0558f20SBarry Smith FormJacobian - Evaluates Jacobian matrix. 140c0558f20SBarry Smith 141c0558f20SBarry Smith Input Parameters: 142c0558f20SBarry Smith . snes - the SNES context 143c0558f20SBarry Smith . x - input vector 144c0558f20SBarry Smith . dummy - optional user-defined context (not used here) 145c0558f20SBarry Smith 146c0558f20SBarry Smith Output Parameters: 147c0558f20SBarry Smith . jac - Jacobian matrix 148c0558f20SBarry Smith . B - optionally different preconditioning matrix 149*0b4b7b1cSBarry Smith 150c0558f20SBarry Smith */ 151d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobian(SNES snes, Vec x, Mat jac, Mat B, void *ctx) 152d71ae5a4SJacob Faibussowitsch { 153c0558f20SBarry Smith ApplicationCtx *user = (ApplicationCtx *)ctx; 154c0558f20SBarry Smith PetscScalar *xx, d, A[3]; 155c0558f20SBarry Smith PetscInt i, j[3], M, xs, xm; 156c0558f20SBarry Smith DM da = user->da; 157c0558f20SBarry Smith 158c0558f20SBarry Smith PetscFunctionBeginUser; 159c0558f20SBarry Smith /* 160c0558f20SBarry Smith Get pointer to vector data 161c0558f20SBarry Smith */ 1629566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, x, &xx)); 1639566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL)); 164c0558f20SBarry Smith 165c0558f20SBarry Smith /* 166c0558f20SBarry Smith Get range of locally owned matrix 167c0558f20SBarry Smith */ 1689566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, NULL, &M, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 169c0558f20SBarry Smith 170c0558f20SBarry Smith /* 171c0558f20SBarry Smith Determine starting and ending local indices for interior grid points. 172c0558f20SBarry Smith Set Jacobian entries for boundary points. 173c0558f20SBarry Smith */ 174c0558f20SBarry Smith 175c0558f20SBarry Smith if (xs == 0) { /* left boundary */ 17611486bccSBarry Smith i = 0; 17711486bccSBarry Smith A[0] = 1.0; 178c0558f20SBarry Smith 1799566063dSJacob Faibussowitsch PetscCall(MatSetValues(jac, 1, &i, 1, &i, A, INSERT_VALUES)); 18011486bccSBarry Smith xs++; 18111486bccSBarry Smith xm--; 182c0558f20SBarry Smith } 183c0558f20SBarry Smith if (xs + xm == M) { /* right boundary */ 184c0558f20SBarry Smith i = M - 1; 185c0558f20SBarry Smith A[0] = 1.0; 1869566063dSJacob Faibussowitsch PetscCall(MatSetValues(jac, 1, &i, 1, &i, A, INSERT_VALUES)); 187c0558f20SBarry Smith xm--; 188c0558f20SBarry Smith } 189c0558f20SBarry Smith 190c0558f20SBarry Smith /* 191c0558f20SBarry Smith Interior grid points 192c0558f20SBarry Smith - Note that in this case we set all elements for a particular 193c0558f20SBarry Smith row at once. 194c0558f20SBarry Smith */ 195c0558f20SBarry Smith d = 1.0 / (user->h * user->h); 196c0558f20SBarry Smith for (i = xs; i < xs + xm; i++) { 19711486bccSBarry Smith j[0] = i - 1; 19811486bccSBarry Smith j[1] = i; 19911486bccSBarry Smith j[2] = i + 1; 20011486bccSBarry Smith A[0] = A[2] = d; 20111486bccSBarry Smith A[1] = -2.0 * d + 2.0 * xx[i]; 2029566063dSJacob Faibussowitsch PetscCall(MatSetValues(jac, 1, &i, 3, j, A, INSERT_VALUES)); 203c0558f20SBarry Smith } 204c0558f20SBarry Smith 205c0558f20SBarry Smith /* 206c0558f20SBarry Smith Assemble matrix, using the 2-step process: 207c0558f20SBarry Smith MatAssemblyBegin(), MatAssemblyEnd(). 208c0558f20SBarry Smith By placing code between these two statements, computations can be 209c0558f20SBarry Smith done while messages are in transition. 210c0558f20SBarry Smith 211c0558f20SBarry Smith Also, restore vector. 212c0558f20SBarry Smith */ 213c0558f20SBarry Smith 2149566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY)); 2159566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, x, &xx)); 2169566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY)); 2173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 218c0558f20SBarry Smith } 219c0558f20SBarry Smith 220d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 221d71ae5a4SJacob Faibussowitsch { 22241e469ddSJunchao Zhang SNES snes; /* SNES context */ 22341e469ddSJunchao Zhang Mat J; /* Jacobian matrix */ 22441e469ddSJunchao Zhang ApplicationCtx ctx; /* user-defined context */ 22541e469ddSJunchao Zhang Vec x, r, U, F; /* vectors */ 22641e469ddSJunchao Zhang PetscScalar none = -1.0; 22741e469ddSJunchao Zhang PetscInt its, N = 5, maxit, maxf; 22841e469ddSJunchao Zhang PetscReal abstol, rtol, stol, norm; 22941e469ddSJunchao Zhang PetscBool viewinitial = PETSC_FALSE; 23041e469ddSJunchao Zhang 231327415f7SBarry Smith PetscFunctionBeginUser; 232c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, nullptr, help)); 2339566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &N, NULL)); 23441e469ddSJunchao Zhang ctx.h = 1.0 / (N - 1); 2359566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-view_initial", &viewinitial, NULL)); 23641e469ddSJunchao Zhang 23741e469ddSJunchao Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 23841e469ddSJunchao Zhang Create nonlinear solver context 23941e469ddSJunchao Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2409566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 24141e469ddSJunchao Zhang 242c0558f20SBarry Smith /* 24341e469ddSJunchao Zhang Create distributed array (DMDA) to manage parallel grid and vectors 244c0558f20SBarry Smith */ 2459566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, N, 1, 1, NULL, &ctx.da)); 2469566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(ctx.da)); 2479566063dSJacob Faibussowitsch PetscCall(DMSetUp(ctx.da)); 24841e469ddSJunchao Zhang 24941e469ddSJunchao Zhang /* 25041e469ddSJunchao Zhang Extract global and local vectors from DMDA; then duplicate for remaining 25141e469ddSJunchao Zhang vectors that are the same types 25241e469ddSJunchao Zhang */ 2539566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(ctx.da, &x)); 2549566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)x, "Approximate Solution")); 2559566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &r)); 25611486bccSBarry Smith PetscCall(VecDuplicate(x, &F)); 25711486bccSBarry Smith ctx.F = F; 2589566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)F, "Forcing function")); 2599566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &U)); 2609566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)U, "Exact Solution")); 26141e469ddSJunchao Zhang 26241e469ddSJunchao Zhang /* 26341e469ddSJunchao Zhang Set function evaluation routine and vector. Whenever the nonlinear 26441e469ddSJunchao Zhang solver needs to compute the nonlinear function, it will call this 26541e469ddSJunchao Zhang routine. 26641e469ddSJunchao Zhang - Note that the final routine argument is the user-defined 26741e469ddSJunchao Zhang context that provides application-specific data for the 26841e469ddSJunchao Zhang function evaluation routine. 26941e469ddSJunchao Zhang 270a5b23f4aSJose E. Roman At the beginning, one can use a stub function that checks the Kokkos version 27141e469ddSJunchao Zhang against the CPU version to quickly expose errors. 2729566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes,r,StubFunction,&ctx)); 27341e469ddSJunchao Zhang */ 2749566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, r, KokkosFunction, &ctx)); 27541e469ddSJunchao Zhang 27641e469ddSJunchao Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 27741e469ddSJunchao Zhang Create matrix data structure; set Jacobian evaluation routine 27841e469ddSJunchao Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2799566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(ctx.da, &J)); 28041e469ddSJunchao Zhang 28141e469ddSJunchao Zhang /* 28241e469ddSJunchao Zhang Set Jacobian matrix data structure and default Jacobian evaluation 28341e469ddSJunchao Zhang routine. Whenever the nonlinear solver needs to compute the 28441e469ddSJunchao Zhang Jacobian matrix, it will call this routine. 28541e469ddSJunchao Zhang - Note that the final routine argument is the user-defined 28641e469ddSJunchao Zhang context that provides application-specific data for the 28741e469ddSJunchao Zhang Jacobian evaluation routine. 28841e469ddSJunchao Zhang */ 2899566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, J, J, FormJacobian, &ctx)); 2909566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 2919566063dSJacob Faibussowitsch PetscCall(SNESGetTolerances(snes, &abstol, &rtol, &stol, &maxit, &maxf)); 29263a3b9bcSJacob 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)); 29341e469ddSJunchao Zhang 29441e469ddSJunchao Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 29541e469ddSJunchao Zhang Initialize application: 29641e469ddSJunchao Zhang Store forcing function of PDE and exact solution 29741e469ddSJunchao Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2989000cb80SJunchao Zhang { 2999000cb80SJunchao Zhang PetscScalarKokkosOffsetView FF, UU; 3009566063dSJacob Faibussowitsch PetscCall(DMDAVecGetKokkosOffsetViewWrite(ctx.da, F, &FF)); 3019566063dSJacob Faibussowitsch PetscCall(DMDAVecGetKokkosOffsetViewWrite(ctx.da, U, &UU)); 30211486bccSBarry Smith Kokkos::parallel_for( 30311486bccSBarry Smith Kokkos::RangePolicy<>(FF.begin(0), FF.end(0)), KOKKOS_LAMBDA(int i) { 30441e469ddSJunchao Zhang PetscReal xp = i * ctx.h; 30541e469ddSJunchao Zhang FF(i) = 6.0 * xp + pow(xp + 1.e-12, 6.0); /* +1.e-12 is to prevent 0^6 */ 30641e469ddSJunchao Zhang UU(i) = xp * xp * xp; 30741e469ddSJunchao Zhang }); 3089566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreKokkosOffsetViewWrite(ctx.da, F, &FF)); 3099566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreKokkosOffsetViewWrite(ctx.da, U, &UU)); 3109000cb80SJunchao Zhang } 31141e469ddSJunchao Zhang 31241e469ddSJunchao Zhang if (viewinitial) { 3139566063dSJacob Faibussowitsch PetscCall(VecView(U, NULL)); 3149566063dSJacob Faibussowitsch PetscCall(VecView(F, NULL)); 31541e469ddSJunchao Zhang } 31641e469ddSJunchao Zhang 31741e469ddSJunchao Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 31841e469ddSJunchao Zhang Evaluate initial guess; then solve nonlinear system 31941e469ddSJunchao Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 32041e469ddSJunchao Zhang 32141e469ddSJunchao Zhang /* 32241e469ddSJunchao Zhang Note: The user should initialize the vector, x, with the initial guess 32341e469ddSJunchao Zhang for the nonlinear solver prior to calling SNESSolve(). In particular, 32441e469ddSJunchao Zhang to employ an initial guess of zero, the user should explicitly set 32541e469ddSJunchao Zhang this vector to zero by calling VecSet(). 32641e469ddSJunchao Zhang */ 3279566063dSJacob Faibussowitsch PetscCall(FormInitialGuess(x)); 3289566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, x)); 3299566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &its)); 33063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %" PetscInt_FMT "\n", its)); 33141e469ddSJunchao Zhang 33241e469ddSJunchao Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 33341e469ddSJunchao Zhang Check solution and clean up 33441e469ddSJunchao Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 33541e469ddSJunchao Zhang /* 33641e469ddSJunchao Zhang Check the error 33741e469ddSJunchao Zhang */ 3389566063dSJacob Faibussowitsch PetscCall(VecAXPY(x, none, U)); 3399566063dSJacob Faibussowitsch PetscCall(VecNorm(x, NORM_2, &norm)); 34063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Norm of error %g Iterations %" PetscInt_FMT "\n", (double)norm, its)); 34141e469ddSJunchao Zhang 34241e469ddSJunchao Zhang /* 34341e469ddSJunchao Zhang Free work space. All PETSc objects should be destroyed when they 34441e469ddSJunchao Zhang are no longer needed. 34541e469ddSJunchao Zhang */ 3469566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 3479566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 3489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 3499566063dSJacob Faibussowitsch PetscCall(VecDestroy(&F)); 3509566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 3519566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 3529566063dSJacob Faibussowitsch PetscCall(DMDestroy(&ctx.da)); 3539566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 354b122ec5aSJacob Faibussowitsch return 0; 35541e469ddSJunchao Zhang } 356c0558f20SBarry Smith 357c0558f20SBarry Smith /*TEST 358c0558f20SBarry Smith 359c0558f20SBarry Smith build: 36041e469ddSJunchao Zhang requires: kokkos_kernels 361c0558f20SBarry Smith 362c0558f20SBarry Smith test: 36387699089SJunchao Zhang requires: kokkos_kernels !complex !single 364c0558f20SBarry Smith nsize: 2 36541e469ddSJunchao Zhang args: -dm_vec_type kokkos -dm_mat_type aijkokkos -view_initial -snes_monitor 366c0558f20SBarry Smith output_file: output/ex3k_1.out 367c0558f20SBarry Smith 368c0558f20SBarry Smith TEST*/ 369