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 PetscScalar pfive = .50; 25c0558f20SBarry Smith 26c0558f20SBarry Smith PetscFunctionBeginUser; 279566063dSJacob Faibussowitsch PetscCall(VecSet(x,pfive)); 28c0558f20SBarry Smith PetscFunctionReturn(0); 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 */ 4741e469ddSJunchao Zhang PetscErrorCode CpuFunction(SNES snes,Vec x,Vec r,void *ctx) 48c0558f20SBarry Smith { 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]; 67c0558f20SBarry Smith xs++;xm--; 68c0558f20SBarry Smith } 69c0558f20SBarry Smith if (xs+xm == M) { /* right boundary */ 7041e469ddSJunchao Zhang R[xs+xm-1] = X[xs+xm-1] - 1.0; 71c0558f20SBarry Smith xm--; 72c0558f20SBarry Smith } 73c0558f20SBarry Smith d = 1.0/(user->h*user->h); 7441e469ddSJunchao 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]; 75c0558f20SBarry Smith 769566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,xl,&X)); 779566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,r,&R)); 789566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,user->F,&F)); 799566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&xl)); 80c0558f20SBarry Smith PetscFunctionReturn(0); 81c0558f20SBarry Smith } 82c0558f20SBarry Smith 8341e469ddSJunchao Zhang using DefaultExecutionSpace = Kokkos::DefaultExecutionSpace; 8441e469ddSJunchao Zhang using DefaultMemorySpace = Kokkos::DefaultExecutionSpace::memory_space; 8541e469ddSJunchao Zhang using PetscScalarKokkosOffsetView = Kokkos::Experimental::OffsetView<PetscScalar*,DefaultMemorySpace>; 8641e469ddSJunchao Zhang using ConstPetscScalarKokkosOffsetView = Kokkos::Experimental::OffsetView<const PetscScalar*,DefaultMemorySpace>; 8741e469ddSJunchao Zhang 8841e469ddSJunchao Zhang PetscErrorCode KokkosFunction(SNES snes,Vec x,Vec r,void *ctx) 8941e469ddSJunchao Zhang { 9041e469ddSJunchao Zhang ApplicationCtx *user = (ApplicationCtx*) ctx; 9141e469ddSJunchao Zhang DM da = user->da; 9241e469ddSJunchao Zhang PetscScalar d; 9341e469ddSJunchao Zhang PetscInt M; 9441e469ddSJunchao Zhang Vec xl; 9541e469ddSJunchao Zhang PetscScalarKokkosOffsetView R; 9641e469ddSJunchao Zhang ConstPetscScalarKokkosOffsetView X,F; 9741e469ddSJunchao Zhang 9841e469ddSJunchao Zhang PetscFunctionBeginUser; 999566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&xl)); 1009566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocal(da,x,INSERT_VALUES,xl)); 10141e469ddSJunchao Zhang d = 1.0/(user->h*user->h); 1029566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,NULL,&M,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL)); 1039566063dSJacob Faibussowitsch PetscCall(DMDAVecGetKokkosOffsetView(da,xl,&X)); /* read only */ 1049566063dSJacob Faibussowitsch PetscCall(DMDAVecGetKokkosOffsetViewWrite(da,r,&R)); /* write only */ 1059566063dSJacob Faibussowitsch PetscCall(DMDAVecGetKokkosOffsetView(da,user->F,&F)); /* read only */ 10641e469ddSJunchao Zhang Kokkos:: parallel_for (Kokkos::RangePolicy<>(R.begin(0),R.end(0)),KOKKOS_LAMBDA (int i) { 10741e469ddSJunchao Zhang if (i == 0) R(0) = X(0); /* left boundary */ 10841e469ddSJunchao Zhang else if (i == M-1) R(i) = X(i) - 1.0; /* right boundary */ 10941e469ddSJunchao Zhang else R(i) = d*(X(i-1) - 2.0*X(i) + X(i+1)) + X(i)*X(i) - F(i); /* interior */ 11041e469ddSJunchao Zhang }); 1119566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreKokkosOffsetView(da,xl,&X)); 1129566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreKokkosOffsetViewWrite(da,r,&R)); 1139566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreKokkosOffsetView(da,user->F,&F)); 1149566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&xl)); 11541e469ddSJunchao Zhang PetscFunctionReturn(0); 11641e469ddSJunchao Zhang } 11741e469ddSJunchao Zhang 11841e469ddSJunchao Zhang PetscErrorCode StubFunction(SNES snes ,Vec x,Vec r,void *ctx) 11941e469ddSJunchao Zhang { 12041e469ddSJunchao Zhang ApplicationCtx *user = (ApplicationCtx*) ctx; 12141e469ddSJunchao Zhang DM da = user->da; 12241e469ddSJunchao Zhang Vec rk; 12341e469ddSJunchao Zhang PetscReal norm=0; 12441e469ddSJunchao Zhang 12541e469ddSJunchao Zhang PetscFunctionBeginUser; 1269566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(da,&rk)); 1279566063dSJacob Faibussowitsch PetscCall(CpuFunction(snes,x,r,ctx)); 1289566063dSJacob Faibussowitsch PetscCall(KokkosFunction(snes,x,rk,ctx)); 1299566063dSJacob Faibussowitsch PetscCall(VecAXPY(rk,-1.0,r)); 1309566063dSJacob Faibussowitsch PetscCall(VecNorm(rk,NORM_2,&norm)); 1319566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(da,&rk)); 13263a3b9bcSJacob Faibussowitsch PetscCheck(norm <= 1e-6,PETSC_COMM_SELF,PETSC_ERR_PLIB,"KokkosFunction() different from CpuFunction() with a diff norm = %g",(double)norm); 13341e469ddSJunchao Zhang PetscFunctionReturn(0); 13441e469ddSJunchao Zhang } 135c0558f20SBarry Smith /* ------------------------------------------------------------------- */ 136c0558f20SBarry Smith /* 137c0558f20SBarry Smith FormJacobian - Evaluates Jacobian matrix. 138c0558f20SBarry Smith 139c0558f20SBarry Smith Input Parameters: 140c0558f20SBarry Smith . snes - the SNES context 141c0558f20SBarry Smith . x - input vector 142c0558f20SBarry Smith . dummy - optional user-defined context (not used here) 143c0558f20SBarry Smith 144c0558f20SBarry Smith Output Parameters: 145c0558f20SBarry Smith . jac - Jacobian matrix 146c0558f20SBarry Smith . B - optionally different preconditioning matrix 147c0558f20SBarry Smith . flag - flag indicating matrix structure 148c0558f20SBarry Smith */ 149c0558f20SBarry Smith PetscErrorCode FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void *ctx) 150c0558f20SBarry Smith { 151c0558f20SBarry Smith ApplicationCtx *user = (ApplicationCtx*) ctx; 152c0558f20SBarry Smith PetscScalar *xx,d,A[3]; 153c0558f20SBarry Smith PetscInt i,j[3],M,xs,xm; 154c0558f20SBarry Smith DM da = user->da; 155c0558f20SBarry Smith 156c0558f20SBarry Smith PetscFunctionBeginUser; 157c0558f20SBarry Smith /* 158c0558f20SBarry Smith Get pointer to vector data 159c0558f20SBarry Smith */ 1609566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da,x,&xx)); 1619566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL)); 162c0558f20SBarry Smith 163c0558f20SBarry Smith /* 164c0558f20SBarry Smith Get range of locally owned matrix 165c0558f20SBarry Smith */ 1669566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,NULL,&M,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL)); 167c0558f20SBarry Smith 168c0558f20SBarry Smith /* 169c0558f20SBarry Smith Determine starting and ending local indices for interior grid points. 170c0558f20SBarry Smith Set Jacobian entries for boundary points. 171c0558f20SBarry Smith */ 172c0558f20SBarry Smith 173c0558f20SBarry Smith if (xs == 0) { /* left boundary */ 174c0558f20SBarry Smith i = 0; A[0] = 1.0; 175c0558f20SBarry Smith 1769566063dSJacob Faibussowitsch PetscCall(MatSetValues(jac,1,&i,1,&i,A,INSERT_VALUES)); 177c0558f20SBarry Smith xs++;xm--; 178c0558f20SBarry Smith } 179c0558f20SBarry Smith if (xs+xm == M) { /* right boundary */ 180c0558f20SBarry Smith i = M-1; 181c0558f20SBarry Smith A[0] = 1.0; 1829566063dSJacob Faibussowitsch PetscCall(MatSetValues(jac,1,&i,1,&i,A,INSERT_VALUES)); 183c0558f20SBarry Smith xm--; 184c0558f20SBarry Smith } 185c0558f20SBarry Smith 186c0558f20SBarry Smith /* 187c0558f20SBarry Smith Interior grid points 188c0558f20SBarry Smith - Note that in this case we set all elements for a particular 189c0558f20SBarry Smith row at once. 190c0558f20SBarry Smith */ 191c0558f20SBarry Smith d = 1.0/(user->h*user->h); 192c0558f20SBarry Smith for (i=xs; i<xs+xm; i++) { 193c0558f20SBarry Smith j[0] = i - 1; j[1] = i; j[2] = i + 1; 194c0558f20SBarry Smith A[0] = A[2] = d; A[1] = -2.0*d + 2.0*xx[i]; 1959566063dSJacob Faibussowitsch PetscCall(MatSetValues(jac,1,&i,3,j,A,INSERT_VALUES)); 196c0558f20SBarry Smith } 197c0558f20SBarry Smith 198c0558f20SBarry Smith /* 199c0558f20SBarry Smith Assemble matrix, using the 2-step process: 200c0558f20SBarry Smith MatAssemblyBegin(), MatAssemblyEnd(). 201c0558f20SBarry Smith By placing code between these two statements, computations can be 202c0558f20SBarry Smith done while messages are in transition. 203c0558f20SBarry Smith 204c0558f20SBarry Smith Also, restore vector. 205c0558f20SBarry Smith */ 206c0558f20SBarry Smith 2079566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); 2089566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da,x,&xx)); 2099566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); 210c0558f20SBarry Smith 211c0558f20SBarry Smith PetscFunctionReturn(0); 212c0558f20SBarry Smith } 213c0558f20SBarry Smith 21441e469ddSJunchao Zhang int main(int argc,char **argv) 21541e469ddSJunchao Zhang { 21641e469ddSJunchao Zhang SNES snes; /* SNES context */ 21741e469ddSJunchao Zhang Mat J; /* Jacobian matrix */ 21841e469ddSJunchao Zhang ApplicationCtx ctx; /* user-defined context */ 21941e469ddSJunchao Zhang Vec x,r,U,F; /* vectors */ 22041e469ddSJunchao Zhang PetscScalar none = -1.0; 22141e469ddSJunchao Zhang PetscInt its,N = 5,maxit,maxf; 22241e469ddSJunchao Zhang PetscReal abstol,rtol,stol,norm; 22341e469ddSJunchao Zhang PetscBool viewinitial = PETSC_FALSE; 22441e469ddSJunchao Zhang 225*327415f7SBarry Smith PetscFunctionBeginUser; 2269566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 2279566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL)); 22841e469ddSJunchao Zhang ctx.h = 1.0/(N-1); 2299566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-view_initial",&viewinitial,NULL)); 23041e469ddSJunchao Zhang 23141e469ddSJunchao Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 23241e469ddSJunchao Zhang Create nonlinear solver context 23341e469ddSJunchao Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2349566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes)); 23541e469ddSJunchao Zhang 236c0558f20SBarry Smith /* 23741e469ddSJunchao Zhang Create distributed array (DMDA) to manage parallel grid and vectors 238c0558f20SBarry Smith */ 2399566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,N,1,1,NULL,&ctx.da)); 2409566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(ctx.da)); 2419566063dSJacob Faibussowitsch PetscCall(DMSetUp(ctx.da)); 24241e469ddSJunchao Zhang 24341e469ddSJunchao Zhang /* 24441e469ddSJunchao Zhang Extract global and local vectors from DMDA; then duplicate for remaining 24541e469ddSJunchao Zhang vectors that are the same types 24641e469ddSJunchao Zhang */ 2479566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(ctx.da,&x)); 2489566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)x,"Approximate Solution")); 2499566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&r)); 2509566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&F)); ctx.F = F; 2519566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)F,"Forcing function")); 2529566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&U)); 2539566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)U,"Exact Solution")); 25441e469ddSJunchao Zhang 25541e469ddSJunchao Zhang /* 25641e469ddSJunchao Zhang Set function evaluation routine and vector. Whenever the nonlinear 25741e469ddSJunchao Zhang solver needs to compute the nonlinear function, it will call this 25841e469ddSJunchao Zhang routine. 25941e469ddSJunchao Zhang - Note that the final routine argument is the user-defined 26041e469ddSJunchao Zhang context that provides application-specific data for the 26141e469ddSJunchao Zhang function evaluation routine. 26241e469ddSJunchao Zhang 263a5b23f4aSJose E. Roman At the beginning, one can use a stub function that checks the Kokkos version 26441e469ddSJunchao Zhang against the CPU version to quickly expose errors. 2659566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes,r,StubFunction,&ctx)); 26641e469ddSJunchao Zhang */ 2679566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes,r,KokkosFunction,&ctx)); 26841e469ddSJunchao Zhang 26941e469ddSJunchao Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 27041e469ddSJunchao Zhang Create matrix data structure; set Jacobian evaluation routine 27141e469ddSJunchao Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2729566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(ctx.da,&J)); 27341e469ddSJunchao Zhang 27441e469ddSJunchao Zhang /* 27541e469ddSJunchao Zhang Set Jacobian matrix data structure and default Jacobian evaluation 27641e469ddSJunchao Zhang routine. Whenever the nonlinear solver needs to compute the 27741e469ddSJunchao Zhang Jacobian matrix, it will call this routine. 27841e469ddSJunchao Zhang - Note that the final routine argument is the user-defined 27941e469ddSJunchao Zhang context that provides application-specific data for the 28041e469ddSJunchao Zhang Jacobian evaluation routine. 28141e469ddSJunchao Zhang */ 2829566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes,J,J,FormJacobian,&ctx)); 2839566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 2849566063dSJacob Faibussowitsch PetscCall(SNESGetTolerances(snes,&abstol,&rtol,&stol,&maxit,&maxf)); 28563a3b9bcSJacob 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)); 28641e469ddSJunchao Zhang 28741e469ddSJunchao Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 28841e469ddSJunchao Zhang Initialize application: 28941e469ddSJunchao Zhang Store forcing function of PDE and exact solution 29041e469ddSJunchao Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2919000cb80SJunchao Zhang { 2929000cb80SJunchao Zhang PetscScalarKokkosOffsetView FF,UU; 2939566063dSJacob Faibussowitsch PetscCall(DMDAVecGetKokkosOffsetViewWrite(ctx.da,F,&FF)); 2949566063dSJacob Faibussowitsch PetscCall(DMDAVecGetKokkosOffsetViewWrite(ctx.da,U,&UU)); 29541e469ddSJunchao Zhang Kokkos:: parallel_for (Kokkos::RangePolicy<>(FF.begin(0),FF.end(0)),KOKKOS_LAMBDA (int i) { 29641e469ddSJunchao Zhang PetscReal xp = i*ctx.h; 29741e469ddSJunchao Zhang FF(i) = 6.0*xp + pow(xp+1.e-12,6.0); /* +1.e-12 is to prevent 0^6 */ 29841e469ddSJunchao Zhang UU(i) = xp*xp*xp; 29941e469ddSJunchao Zhang }); 3009566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreKokkosOffsetViewWrite(ctx.da,F,&FF)); 3019566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreKokkosOffsetViewWrite(ctx.da,U,&UU)); 3029000cb80SJunchao Zhang } 30341e469ddSJunchao Zhang 30441e469ddSJunchao Zhang if (viewinitial) { 3059566063dSJacob Faibussowitsch PetscCall(VecView(U,NULL)); 3069566063dSJacob Faibussowitsch PetscCall(VecView(F,NULL)); 30741e469ddSJunchao Zhang } 30841e469ddSJunchao Zhang 30941e469ddSJunchao Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 31041e469ddSJunchao Zhang Evaluate initial guess; then solve nonlinear system 31141e469ddSJunchao Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 31241e469ddSJunchao Zhang 31341e469ddSJunchao Zhang /* 31441e469ddSJunchao Zhang Note: The user should initialize the vector, x, with the initial guess 31541e469ddSJunchao Zhang for the nonlinear solver prior to calling SNESSolve(). In particular, 31641e469ddSJunchao Zhang to employ an initial guess of zero, the user should explicitly set 31741e469ddSJunchao Zhang this vector to zero by calling VecSet(). 31841e469ddSJunchao Zhang */ 3199566063dSJacob Faibussowitsch PetscCall(FormInitialGuess(x)); 3209566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes,NULL,x)); 3219566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes,&its)); 32263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %" PetscInt_FMT "\n",its)); 32341e469ddSJunchao Zhang 32441e469ddSJunchao Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 32541e469ddSJunchao Zhang Check solution and clean up 32641e469ddSJunchao Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 32741e469ddSJunchao Zhang /* 32841e469ddSJunchao Zhang Check the error 32941e469ddSJunchao Zhang */ 3309566063dSJacob Faibussowitsch PetscCall(VecAXPY(x,none,U)); 3319566063dSJacob Faibussowitsch PetscCall(VecNorm(x,NORM_2,&norm)); 33263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g Iterations %" PetscInt_FMT "\n",(double)norm,its)); 33341e469ddSJunchao Zhang 33441e469ddSJunchao Zhang /* 33541e469ddSJunchao Zhang Free work space. All PETSc objects should be destroyed when they 33641e469ddSJunchao Zhang are no longer needed. 33741e469ddSJunchao Zhang */ 3389566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 3399566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 3409566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 3419566063dSJacob Faibussowitsch PetscCall(VecDestroy(&F)); 3429566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 3439566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 3449566063dSJacob Faibussowitsch PetscCall(DMDestroy(&ctx.da)); 3459566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 346b122ec5aSJacob Faibussowitsch return 0; 34741e469ddSJunchao Zhang } 348c0558f20SBarry Smith 349c0558f20SBarry Smith /*TEST 350c0558f20SBarry Smith 351c0558f20SBarry Smith build: 35241e469ddSJunchao Zhang requires: kokkos_kernels 353c0558f20SBarry Smith 354c0558f20SBarry Smith test: 35587699089SJunchao Zhang requires: kokkos_kernels !complex !single 356c0558f20SBarry Smith nsize: 2 35741e469ddSJunchao Zhang args: -dm_vec_type kokkos -dm_mat_type aijkokkos -view_initial -snes_monitor 358c0558f20SBarry Smith output_file: output/ex3k_1.out 359c0558f20SBarry Smith 360c0558f20SBarry Smith TEST*/ 361