xref: /petsc/src/snes/tutorials/ex55.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c5566c22SJunchao Zhang static char help[] = "Copy of ex5.c\n";
2c5566c22SJunchao Zhang 
3c5566c22SJunchao Zhang /* ------------------------------------------------------------------------
4c5566c22SJunchao Zhang 
5c5566c22SJunchao Zhang   Copy of ex5.c.
6c5566c22SJunchao Zhang   Once petsc test harness supports conditional linking, we can remove this duplicate.
7c5566c22SJunchao Zhang   See https://gitlab.com/petsc/petsc/-/issues/1173
8c5566c22SJunchao Zhang   ------------------------------------------------------------------------- */
9c5566c22SJunchao Zhang 
10c5566c22SJunchao Zhang /*
11c5566c22SJunchao Zhang    Include "petscdmda.h" so that we can use distributed arrays (DMDAs).
12c5566c22SJunchao Zhang    Include "petscsnes.h" so that we can use SNES solvers.  Note that this
13c5566c22SJunchao Zhang */
14c5566c22SJunchao Zhang #include <petscdm.h>
15c5566c22SJunchao Zhang #include <petscdmda.h>
16c5566c22SJunchao Zhang #include <petscsnes.h>
17c5566c22SJunchao Zhang #include <petscmatlab.h>
18c5566c22SJunchao Zhang #include <petsc/private/snesimpl.h> /* For SNES_Solve event */
19c5566c22SJunchao Zhang #include "ex55.h"
20c5566c22SJunchao Zhang 
21c5566c22SJunchao Zhang /* ------------------------------------------------------------------- */
22c5566c22SJunchao Zhang /*
23c5566c22SJunchao Zhang    FormInitialGuess - Forms initial approximation.
24c5566c22SJunchao Zhang 
25c5566c22SJunchao Zhang    Input Parameters:
26c5566c22SJunchao Zhang    da - The DM
27c5566c22SJunchao Zhang    user - user-defined application context
28c5566c22SJunchao Zhang 
29c5566c22SJunchao Zhang    Output Parameter:
30c5566c22SJunchao Zhang    X - vector
31c5566c22SJunchao Zhang  */
32c5566c22SJunchao Zhang static PetscErrorCode FormInitialGuess(DM da,AppCtx *user,Vec X)
33c5566c22SJunchao Zhang {
34c5566c22SJunchao Zhang   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
35c5566c22SJunchao Zhang   PetscReal      lambda,temp1,temp,hx,hy;
36c5566c22SJunchao Zhang   PetscScalar    **x;
37c5566c22SJunchao Zhang 
38c5566c22SJunchao Zhang   PetscFunctionBeginUser;
39c5566c22SJunchao Zhang   PetscCall(DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE));
40c5566c22SJunchao Zhang 
41c5566c22SJunchao Zhang   lambda = user->param;
42c5566c22SJunchao Zhang   hx     = 1.0/(PetscReal)(Mx-1);
43c5566c22SJunchao Zhang   hy     = 1.0/(PetscReal)(My-1);
44c5566c22SJunchao Zhang   temp1  = lambda/(lambda + 1.0);
45c5566c22SJunchao Zhang 
46c5566c22SJunchao Zhang   /*
47c5566c22SJunchao Zhang      Get a pointer to vector data.
48c5566c22SJunchao Zhang        - For default PETSc vectors, VecGetArray() returns a pointer to
49c5566c22SJunchao Zhang          the data array.  Otherwise, the routine is implementation dependent.
50c5566c22SJunchao Zhang        - You MUST call VecRestoreArray() when you no longer need access to
51c5566c22SJunchao Zhang          the array.
52c5566c22SJunchao Zhang   */
53c5566c22SJunchao Zhang   PetscCall(DMDAVecGetArray(da,X,&x));
54c5566c22SJunchao Zhang 
55c5566c22SJunchao Zhang   /*
56c5566c22SJunchao Zhang      Get local grid boundaries (for 2-dimensional DMDA):
57c5566c22SJunchao Zhang        xs, ys   - starting grid indices (no ghost points)
58c5566c22SJunchao Zhang        xm, ym   - widths of local grid (no ghost points)
59c5566c22SJunchao Zhang 
60c5566c22SJunchao Zhang   */
61c5566c22SJunchao Zhang   PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
62c5566c22SJunchao Zhang 
63c5566c22SJunchao Zhang   /*
64c5566c22SJunchao Zhang      Compute initial guess over the locally owned part of the grid
65c5566c22SJunchao Zhang   */
66c5566c22SJunchao Zhang   for (j=ys; j<ys+ym; j++) {
67c5566c22SJunchao Zhang     temp = (PetscReal)(PetscMin(j,My-j-1))*hy;
68c5566c22SJunchao Zhang     for (i=xs; i<xs+xm; i++) {
69c5566c22SJunchao Zhang       if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
70c5566c22SJunchao Zhang         /* boundary conditions are all zero Dirichlet */
71c5566c22SJunchao Zhang         x[j][i] = 0.0;
72c5566c22SJunchao Zhang       } else {
73c5566c22SJunchao Zhang         x[j][i] = temp1*PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i,Mx-i-1))*hx,temp));
74c5566c22SJunchao Zhang       }
75c5566c22SJunchao Zhang     }
76c5566c22SJunchao Zhang   }
77c5566c22SJunchao Zhang 
78c5566c22SJunchao Zhang   /*
79c5566c22SJunchao Zhang      Restore vector
80c5566c22SJunchao Zhang   */
81c5566c22SJunchao Zhang   PetscCall(DMDAVecRestoreArray(da,X,&x));
82c5566c22SJunchao Zhang   PetscFunctionReturn(0);
83c5566c22SJunchao Zhang }
84c5566c22SJunchao Zhang 
85c5566c22SJunchao Zhang /*
86c5566c22SJunchao Zhang   FormExactSolution - Forms MMS solution
87c5566c22SJunchao Zhang 
88c5566c22SJunchao Zhang   Input Parameters:
89c5566c22SJunchao Zhang   da - The DM
90c5566c22SJunchao Zhang   user - user-defined application context
91c5566c22SJunchao Zhang 
92c5566c22SJunchao Zhang   Output Parameter:
93c5566c22SJunchao Zhang   X - vector
94c5566c22SJunchao Zhang  */
95c5566c22SJunchao Zhang static PetscErrorCode FormExactSolution(DM da, AppCtx *user, Vec U)
96c5566c22SJunchao Zhang {
97c5566c22SJunchao Zhang   DM             coordDA;
98c5566c22SJunchao Zhang   Vec            coordinates;
99c5566c22SJunchao Zhang   DMDACoor2d   **coords;
100c5566c22SJunchao Zhang   PetscScalar  **u;
101c5566c22SJunchao Zhang   PetscInt       xs, ys, xm, ym, i, j;
102c5566c22SJunchao Zhang 
103c5566c22SJunchao Zhang   PetscFunctionBeginUser;
104c5566c22SJunchao Zhang   PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
105c5566c22SJunchao Zhang   PetscCall(DMGetCoordinateDM(da, &coordDA));
106c5566c22SJunchao Zhang   PetscCall(DMGetCoordinates(da, &coordinates));
107c5566c22SJunchao Zhang   PetscCall(DMDAVecGetArray(coordDA, coordinates, &coords));
108c5566c22SJunchao Zhang   PetscCall(DMDAVecGetArray(da, U, &u));
109c5566c22SJunchao Zhang   for (j = ys; j < ys+ym; ++j) {
110c5566c22SJunchao Zhang     for (i = xs; i < xs+xm; ++i) {
111c5566c22SJunchao Zhang       user->mms_solution(user,&coords[j][i],&u[j][i]);
112c5566c22SJunchao Zhang     }
113c5566c22SJunchao Zhang   }
114c5566c22SJunchao Zhang   PetscCall(DMDAVecRestoreArray(da, U, &u));
115c5566c22SJunchao Zhang   PetscCall(DMDAVecRestoreArray(coordDA, coordinates, &coords));
116c5566c22SJunchao Zhang   PetscFunctionReturn(0);
117c5566c22SJunchao Zhang }
118c5566c22SJunchao Zhang 
119c5566c22SJunchao Zhang static PetscErrorCode ZeroBCSolution(AppCtx *user,const DMDACoor2d *c,PetscScalar *u)
120c5566c22SJunchao Zhang {
121c5566c22SJunchao Zhang   u[0] = 0.;
122c5566c22SJunchao Zhang   return 0;
123c5566c22SJunchao Zhang }
124c5566c22SJunchao Zhang 
125c5566c22SJunchao Zhang /* The functions below evaluate the MMS solution u(x,y) and associated forcing
126c5566c22SJunchao Zhang 
127c5566c22SJunchao Zhang      f(x,y) = -u_xx - u_yy - lambda exp(u)
128c5566c22SJunchao Zhang 
129c5566c22SJunchao Zhang   such that u(x,y) is an exact solution with f(x,y) as the right hand side forcing term.
130c5566c22SJunchao Zhang  */
131c5566c22SJunchao Zhang static PetscErrorCode MMSSolution1(AppCtx *user,const DMDACoor2d *c,PetscScalar *u)
132c5566c22SJunchao Zhang {
133c5566c22SJunchao Zhang   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
134c5566c22SJunchao Zhang   u[0] = x*(1 - x)*y*(1 - y);
135c5566c22SJunchao Zhang   PetscLogFlops(5);
136c5566c22SJunchao Zhang   return 0;
137c5566c22SJunchao Zhang }
138c5566c22SJunchao Zhang static PetscErrorCode MMSForcing1(AppCtx *user,const DMDACoor2d *c,PetscScalar *f)
139c5566c22SJunchao Zhang {
140c5566c22SJunchao Zhang   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
141c5566c22SJunchao Zhang   f[0] = 2*x*(1 - x) + 2*y*(1 - y) - user->param*PetscExpReal(x*(1 - x)*y*(1 - y));
142c5566c22SJunchao Zhang   return 0;
143c5566c22SJunchao Zhang }
144c5566c22SJunchao Zhang 
145c5566c22SJunchao Zhang static PetscErrorCode MMSSolution2(AppCtx *user,const DMDACoor2d *c,PetscScalar *u)
146c5566c22SJunchao Zhang {
147c5566c22SJunchao Zhang   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
148c5566c22SJunchao Zhang   u[0] = PetscSinReal(PETSC_PI*x)*PetscSinReal(PETSC_PI*y);
149c5566c22SJunchao Zhang   PetscLogFlops(5);
150c5566c22SJunchao Zhang   return 0;
151c5566c22SJunchao Zhang }
152c5566c22SJunchao Zhang static PetscErrorCode MMSForcing2(AppCtx *user,const DMDACoor2d *c,PetscScalar *f)
153c5566c22SJunchao Zhang {
154c5566c22SJunchao Zhang   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
155c5566c22SJunchao Zhang   f[0] = 2*PetscSqr(PETSC_PI)*PetscSinReal(PETSC_PI*x)*PetscSinReal(PETSC_PI*y) - user->param*PetscExpReal(PetscSinReal(PETSC_PI*x)*PetscSinReal(PETSC_PI*y));
156c5566c22SJunchao Zhang   return 0;
157c5566c22SJunchao Zhang }
158c5566c22SJunchao Zhang 
159c5566c22SJunchao Zhang static PetscErrorCode MMSSolution3(AppCtx *user,const DMDACoor2d *c,PetscScalar *u)
160c5566c22SJunchao Zhang {
161c5566c22SJunchao Zhang   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
162c5566c22SJunchao Zhang   u[0] = PetscSinReal(user->m*PETSC_PI*x*(1-y))*PetscSinReal(user->n*PETSC_PI*y*(1-x));
163c5566c22SJunchao Zhang   PetscLogFlops(5);
164c5566c22SJunchao Zhang   return 0;
165c5566c22SJunchao Zhang }
166c5566c22SJunchao Zhang static PetscErrorCode MMSForcing3(AppCtx *user,const DMDACoor2d *c,PetscScalar *f)
167c5566c22SJunchao Zhang {
168c5566c22SJunchao Zhang   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
169c5566c22SJunchao Zhang   PetscReal m = user->m, n = user->n, lambda = user->param;
170c5566c22SJunchao Zhang   f[0] = (-(PetscExpReal(PetscSinReal(m*PETSC_PI*x*(1 - y))*PetscSinReal(n*PETSC_PI*(1 - x)*y))*lambda)
171c5566c22SJunchao Zhang           + PetscSqr(PETSC_PI)*(-2*m*n*((-1 + x)*x + (-1 + y)*y)*PetscCosReal(m*PETSC_PI*x*(-1 + y))*PetscCosReal(n*PETSC_PI*(-1 + x)*y)
172c5566c22SJunchao Zhang                                 + (PetscSqr(m)*(PetscSqr(x) + PetscSqr(-1 + y)) + PetscSqr(n)*(PetscSqr(-1 + x) + PetscSqr(y)))
173c5566c22SJunchao Zhang                                 *PetscSinReal(m*PETSC_PI*x*(-1 + y))*PetscSinReal(n*PETSC_PI*(-1 + x)*y)));
174c5566c22SJunchao Zhang   return 0;
175c5566c22SJunchao Zhang }
176c5566c22SJunchao Zhang 
177c5566c22SJunchao Zhang static PetscErrorCode MMSSolution4(AppCtx *user,const DMDACoor2d *c,PetscScalar *u)
178c5566c22SJunchao Zhang {
179c5566c22SJunchao Zhang   const PetscReal Lx = 1.,Ly = 1.;
180c5566c22SJunchao Zhang   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
181c5566c22SJunchao Zhang   u[0] = (PetscPowReal(x,4)-PetscSqr(Lx)*PetscSqr(x))*(PetscPowReal(y,4)-PetscSqr(Ly)*PetscSqr(y));
182c5566c22SJunchao Zhang   PetscLogFlops(9);
183c5566c22SJunchao Zhang   return 0;
184c5566c22SJunchao Zhang }
185c5566c22SJunchao Zhang static PetscErrorCode MMSForcing4(AppCtx *user,const DMDACoor2d *c,PetscScalar *f)
186c5566c22SJunchao Zhang {
187c5566c22SJunchao Zhang   const PetscReal Lx = 1.,Ly = 1.;
188c5566c22SJunchao Zhang   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
189c5566c22SJunchao Zhang   f[0] = (2*PetscSqr(x)*(PetscSqr(x)-PetscSqr(Lx))*(PetscSqr(Ly)-6*PetscSqr(y))
190c5566c22SJunchao Zhang           + 2*PetscSqr(y)*(PetscSqr(Lx)-6*PetscSqr(x))*(PetscSqr(y)-PetscSqr(Ly))
191c5566c22SJunchao Zhang           - user->param*PetscExpReal((PetscPowReal(x,4)-PetscSqr(Lx)*PetscSqr(x))*(PetscPowReal(y,4)-PetscSqr(Ly)*PetscSqr(y))));
192c5566c22SJunchao Zhang   return 0;
193c5566c22SJunchao Zhang }
194c5566c22SJunchao Zhang 
195c5566c22SJunchao Zhang /* ------------------------------------------------------------------- */
196c5566c22SJunchao Zhang /*
197c5566c22SJunchao Zhang    FormFunctionLocal - Evaluates nonlinear function, F(x) on local process patch
198c5566c22SJunchao Zhang 
199c5566c22SJunchao Zhang  */
200c5566c22SJunchao Zhang static PetscErrorCode FormFunctionLocal(DMDALocalInfo *info,PetscScalar **x,PetscScalar **f,AppCtx *user)
201c5566c22SJunchao Zhang {
202c5566c22SJunchao Zhang   PetscInt       i,j;
203c5566c22SJunchao Zhang   PetscReal      lambda,hx,hy,hxdhy,hydhx;
204c5566c22SJunchao Zhang   PetscScalar    u,ue,uw,un,us,uxx,uyy,mms_solution,mms_forcing;
205c5566c22SJunchao Zhang   DMDACoor2d     c;
206c5566c22SJunchao Zhang 
207c5566c22SJunchao Zhang   PetscFunctionBeginUser;
208c5566c22SJunchao Zhang   lambda = user->param;
209c5566c22SJunchao Zhang   hx     = 1.0/(PetscReal)(info->mx-1);
210c5566c22SJunchao Zhang   hy     = 1.0/(PetscReal)(info->my-1);
211c5566c22SJunchao Zhang   hxdhy  = hx/hy;
212c5566c22SJunchao Zhang   hydhx  = hy/hx;
213c5566c22SJunchao Zhang   /*
214c5566c22SJunchao Zhang      Compute function over the locally owned part of the grid
215c5566c22SJunchao Zhang   */
216c5566c22SJunchao Zhang   for (j=info->ys; j<info->ys+info->ym; j++) {
217c5566c22SJunchao Zhang     for (i=info->xs; i<info->xs+info->xm; i++) {
218c5566c22SJunchao Zhang       if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) {
219c5566c22SJunchao Zhang         c.x = i*hx; c.y = j*hy;
220c5566c22SJunchao Zhang         PetscCall(user->mms_solution(user,&c,&mms_solution));
221c5566c22SJunchao Zhang         f[j][i] = 2.0*(hydhx+hxdhy)*(x[j][i] - mms_solution);
222c5566c22SJunchao Zhang       } else {
223c5566c22SJunchao Zhang         u  = x[j][i];
224c5566c22SJunchao Zhang         uw = x[j][i-1];
225c5566c22SJunchao Zhang         ue = x[j][i+1];
226c5566c22SJunchao Zhang         un = x[j-1][i];
227c5566c22SJunchao Zhang         us = x[j+1][i];
228c5566c22SJunchao Zhang 
229c5566c22SJunchao Zhang         /* Enforce boundary conditions at neighboring points -- setting these values causes the Jacobian to be symmetric. */
230c5566c22SJunchao Zhang         if (i-1 == 0) {c.x = (i-1)*hx; c.y = j*hy; PetscCall(user->mms_solution(user,&c,&uw));}
231c5566c22SJunchao Zhang         if (i+1 == info->mx-1) {c.x = (i+1)*hx; c.y = j*hy; PetscCall(user->mms_solution(user,&c,&ue));}
232c5566c22SJunchao Zhang         if (j-1 == 0) {c.x = i*hx; c.y = (j-1)*hy; PetscCall(user->mms_solution(user,&c,&un));}
233c5566c22SJunchao Zhang         if (j+1 == info->my-1) {c.x = i*hx; c.y = (j+1)*hy; PetscCall(user->mms_solution(user,&c,&us));}
234c5566c22SJunchao Zhang 
235c5566c22SJunchao Zhang         uxx     = (2.0*u - uw - ue)*hydhx;
236c5566c22SJunchao Zhang         uyy     = (2.0*u - un - us)*hxdhy;
237c5566c22SJunchao Zhang         mms_forcing = 0;
238c5566c22SJunchao Zhang         c.x = i*hx; c.y = j*hy;
239c5566c22SJunchao Zhang         if (user->mms_forcing) PetscCall(user->mms_forcing(user,&c,&mms_forcing));
240c5566c22SJunchao Zhang         f[j][i] = uxx + uyy - hx*hy*(lambda*PetscExpScalar(u) + mms_forcing);
241c5566c22SJunchao Zhang       }
242c5566c22SJunchao Zhang     }
243c5566c22SJunchao Zhang   }
244c5566c22SJunchao Zhang   PetscCall(PetscLogFlops(11.0*info->ym*info->xm));
245c5566c22SJunchao Zhang   PetscFunctionReturn(0);
246c5566c22SJunchao Zhang }
247c5566c22SJunchao Zhang 
248c5566c22SJunchao Zhang /* FormObjectiveLocal - Evaluates nonlinear function, F(x) on local process patch */
249c5566c22SJunchao Zhang static PetscErrorCode FormObjectiveLocal(DMDALocalInfo *info,PetscScalar **x,PetscReal *obj,AppCtx *user)
250c5566c22SJunchao Zhang {
251c5566c22SJunchao Zhang   PetscInt       i,j;
252c5566c22SJunchao Zhang   PetscReal      lambda,hx,hy,hxdhy,hydhx,sc,lobj=0;
253c5566c22SJunchao Zhang   PetscScalar    u,ue,uw,un,us,uxux,uyuy;
254c5566c22SJunchao Zhang   MPI_Comm       comm;
255c5566c22SJunchao Zhang 
256c5566c22SJunchao Zhang   PetscFunctionBeginUser;
257c5566c22SJunchao Zhang   *obj   = 0;
258c5566c22SJunchao Zhang   PetscCall(PetscObjectGetComm((PetscObject)info->da,&comm));
259c5566c22SJunchao Zhang   lambda = user->param;
260c5566c22SJunchao Zhang   hx     = 1.0/(PetscReal)(info->mx-1);
261c5566c22SJunchao Zhang   hy     = 1.0/(PetscReal)(info->my-1);
262c5566c22SJunchao Zhang   sc     = hx*hy*lambda;
263c5566c22SJunchao Zhang   hxdhy  = hx/hy;
264c5566c22SJunchao Zhang   hydhx  = hy/hx;
265c5566c22SJunchao Zhang   /*
266c5566c22SJunchao Zhang      Compute function over the locally owned part of the grid
267c5566c22SJunchao Zhang   */
268c5566c22SJunchao Zhang   for (j=info->ys; j<info->ys+info->ym; j++) {
269c5566c22SJunchao Zhang     for (i=info->xs; i<info->xs+info->xm; i++) {
270c5566c22SJunchao Zhang       if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) {
271c5566c22SJunchao Zhang         lobj += PetscRealPart((hydhx + hxdhy)*x[j][i]*x[j][i]);
272c5566c22SJunchao Zhang       } else {
273c5566c22SJunchao Zhang         u  = x[j][i];
274c5566c22SJunchao Zhang         uw = x[j][i-1];
275c5566c22SJunchao Zhang         ue = x[j][i+1];
276c5566c22SJunchao Zhang         un = x[j-1][i];
277c5566c22SJunchao Zhang         us = x[j+1][i];
278c5566c22SJunchao Zhang 
279c5566c22SJunchao Zhang         if (i-1 == 0) uw = 0.;
280c5566c22SJunchao Zhang         if (i+1 == info->mx-1) ue = 0.;
281c5566c22SJunchao Zhang         if (j-1 == 0) un = 0.;
282c5566c22SJunchao Zhang         if (j+1 == info->my-1) us = 0.;
283c5566c22SJunchao Zhang 
284c5566c22SJunchao Zhang         /* F[u] = 1/2\int_{\omega}\nabla^2u(x)*u(x)*dx */
285c5566c22SJunchao Zhang 
286c5566c22SJunchao Zhang         uxux = u*(2.*u - ue - uw)*hydhx;
287c5566c22SJunchao Zhang         uyuy = u*(2.*u - un - us)*hxdhy;
288c5566c22SJunchao Zhang 
289c5566c22SJunchao Zhang         lobj += PetscRealPart(0.5*(uxux + uyuy) - sc*PetscExpScalar(u));
290c5566c22SJunchao Zhang       }
291c5566c22SJunchao Zhang     }
292c5566c22SJunchao Zhang   }
293c5566c22SJunchao Zhang   PetscCall(PetscLogFlops(12.0*info->ym*info->xm));
294c5566c22SJunchao Zhang   PetscCallMPI(MPI_Allreduce(&lobj,obj,1,MPIU_REAL,MPIU_SUM,comm));
295c5566c22SJunchao Zhang   PetscFunctionReturn(0);
296c5566c22SJunchao Zhang }
297c5566c22SJunchao Zhang 
298c5566c22SJunchao Zhang /*
299c5566c22SJunchao Zhang    FormJacobianLocal - Evaluates Jacobian matrix on local process patch
300c5566c22SJunchao Zhang */
301c5566c22SJunchao Zhang static PetscErrorCode FormJacobianLocal(DMDALocalInfo *info,PetscScalar **x,Mat jac,Mat jacpre,AppCtx *user)
302c5566c22SJunchao Zhang {
303c5566c22SJunchao Zhang   PetscInt       i,j,k;
304c5566c22SJunchao Zhang   MatStencil     col[5],row;
305c5566c22SJunchao Zhang   PetscScalar    lambda,v[5],hx,hy,hxdhy,hydhx,sc;
306c5566c22SJunchao Zhang   DM             coordDA;
307c5566c22SJunchao Zhang   Vec            coordinates;
308c5566c22SJunchao Zhang   DMDACoor2d   **coords;
309c5566c22SJunchao Zhang 
310c5566c22SJunchao Zhang   PetscFunctionBeginUser;
311c5566c22SJunchao Zhang   lambda = user->param;
312c5566c22SJunchao Zhang   /* Extract coordinates */
313c5566c22SJunchao Zhang   PetscCall(DMGetCoordinateDM(info->da, &coordDA));
314c5566c22SJunchao Zhang   PetscCall(DMGetCoordinates(info->da, &coordinates));
315c5566c22SJunchao Zhang   PetscCall(DMDAVecGetArray(coordDA, coordinates, &coords));
316c5566c22SJunchao Zhang   hx     = info->xm > 1 ? PetscRealPart(coords[info->ys][info->xs+1].x) - PetscRealPart(coords[info->ys][info->xs].x) : 1.0;
317c5566c22SJunchao Zhang   hy     = info->ym > 1 ? PetscRealPart(coords[info->ys+1][info->xs].y) - PetscRealPart(coords[info->ys][info->xs].y) : 1.0;
318c5566c22SJunchao Zhang   PetscCall(DMDAVecRestoreArray(coordDA, coordinates, &coords));
319c5566c22SJunchao Zhang   hxdhy  = hx/hy;
320c5566c22SJunchao Zhang   hydhx  = hy/hx;
321c5566c22SJunchao Zhang   sc     = hx*hy*lambda;
322c5566c22SJunchao Zhang 
323c5566c22SJunchao Zhang   /*
324c5566c22SJunchao Zhang      Compute entries for the locally owned part of the Jacobian.
325c5566c22SJunchao Zhang       - Currently, all PETSc parallel matrix formats are partitioned by
326c5566c22SJunchao Zhang         contiguous chunks of rows across the processors.
327c5566c22SJunchao Zhang       - Each processor needs to insert only elements that it owns
328c5566c22SJunchao Zhang         locally (but any non-local elements will be sent to the
329c5566c22SJunchao Zhang         appropriate processor during matrix assembly).
330c5566c22SJunchao Zhang       - Here, we set all entries for a particular row at once.
331c5566c22SJunchao Zhang       - We can set matrix entries either using either
332c5566c22SJunchao Zhang         MatSetValuesLocal() or MatSetValues(), as discussed above.
333c5566c22SJunchao Zhang   */
334c5566c22SJunchao Zhang   for (j=info->ys; j<info->ys+info->ym; j++) {
335c5566c22SJunchao Zhang     for (i=info->xs; i<info->xs+info->xm; i++) {
336c5566c22SJunchao Zhang       row.j = j; row.i = i;
337c5566c22SJunchao Zhang       /* boundary points */
338c5566c22SJunchao Zhang       if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) {
339c5566c22SJunchao Zhang         v[0] =  2.0*(hydhx + hxdhy);
340c5566c22SJunchao Zhang         PetscCall(MatSetValuesStencil(jacpre,1,&row,1,&row,v,INSERT_VALUES));
341c5566c22SJunchao Zhang       } else {
342c5566c22SJunchao Zhang         k = 0;
343c5566c22SJunchao Zhang         /* interior grid points */
344c5566c22SJunchao Zhang         if (j-1 != 0) {
345c5566c22SJunchao Zhang           v[k]     = -hxdhy;
346c5566c22SJunchao Zhang           col[k].j = j - 1; col[k].i = i;
347c5566c22SJunchao Zhang           k++;
348c5566c22SJunchao Zhang         }
349c5566c22SJunchao Zhang         if (i-1 != 0) {
350c5566c22SJunchao Zhang           v[k]     = -hydhx;
351c5566c22SJunchao Zhang           col[k].j = j;     col[k].i = i-1;
352c5566c22SJunchao Zhang           k++;
353c5566c22SJunchao Zhang         }
354c5566c22SJunchao Zhang 
355c5566c22SJunchao Zhang         v[k] = 2.0*(hydhx + hxdhy) - sc*PetscExpScalar(x[j][i]); col[k].j = row.j; col[k].i = row.i; k++;
356c5566c22SJunchao Zhang 
357c5566c22SJunchao Zhang         if (i+1 != info->mx-1) {
358c5566c22SJunchao Zhang           v[k]     = -hydhx;
359c5566c22SJunchao Zhang           col[k].j = j;     col[k].i = i+1;
360c5566c22SJunchao Zhang           k++;
361c5566c22SJunchao Zhang         }
362c5566c22SJunchao Zhang         if (j+1 != info->mx-1) {
363c5566c22SJunchao Zhang           v[k]     = -hxdhy;
364c5566c22SJunchao Zhang           col[k].j = j + 1; col[k].i = i;
365c5566c22SJunchao Zhang           k++;
366c5566c22SJunchao Zhang         }
367c5566c22SJunchao Zhang         PetscCall(MatSetValuesStencil(jacpre,1,&row,k,col,v,INSERT_VALUES));
368c5566c22SJunchao Zhang       }
369c5566c22SJunchao Zhang     }
370c5566c22SJunchao Zhang   }
371c5566c22SJunchao Zhang 
372c5566c22SJunchao Zhang   /*
373c5566c22SJunchao Zhang      Assemble matrix, using the 2-step process:
374c5566c22SJunchao Zhang        MatAssemblyBegin(), MatAssemblyEnd().
375c5566c22SJunchao Zhang   */
376c5566c22SJunchao Zhang   PetscCall(MatAssemblyBegin(jacpre,MAT_FINAL_ASSEMBLY));
377c5566c22SJunchao Zhang   PetscCall(MatAssemblyEnd(jacpre,MAT_FINAL_ASSEMBLY));
378c5566c22SJunchao Zhang 
379c5566c22SJunchao Zhang   /*
380c5566c22SJunchao Zhang      Tell the matrix we will never add a new nonzero location to the
381c5566c22SJunchao Zhang      matrix. If we do, it will generate an error.
382c5566c22SJunchao Zhang   */
383c5566c22SJunchao Zhang   PetscCall(MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE));
384c5566c22SJunchao Zhang   PetscFunctionReturn(0);
385c5566c22SJunchao Zhang }
386c5566c22SJunchao Zhang 
387c5566c22SJunchao Zhang static PetscErrorCode FormFunctionMatlab(SNES snes,Vec X,Vec F,void *ptr)
388c5566c22SJunchao Zhang {
389c5566c22SJunchao Zhang #if PetscDefined(HAVE_MATLAB_ENGINE)
390c5566c22SJunchao Zhang   AppCtx         *user = (AppCtx*)ptr;
391c5566c22SJunchao Zhang   PetscInt       Mx,My;
392c5566c22SJunchao Zhang   PetscReal      lambda,hx,hy;
393c5566c22SJunchao Zhang   Vec            localX,localF;
394c5566c22SJunchao Zhang   MPI_Comm       comm;
395c5566c22SJunchao Zhang   DM             da;
396c5566c22SJunchao Zhang 
397c5566c22SJunchao Zhang   PetscFunctionBeginUser;
398c5566c22SJunchao Zhang   PetscCall(SNESGetDM(snes,&da));
399c5566c22SJunchao Zhang   PetscCall(DMGetLocalVector(da,&localX));
400c5566c22SJunchao Zhang   PetscCall(DMGetLocalVector(da,&localF));
401c5566c22SJunchao Zhang   PetscCall(PetscObjectSetName((PetscObject)localX,"localX"));
402c5566c22SJunchao Zhang   PetscCall(PetscObjectSetName((PetscObject)localF,"localF"));
403c5566c22SJunchao Zhang   PetscCall(DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE));
404c5566c22SJunchao Zhang 
405c5566c22SJunchao Zhang   lambda = user->param;
406c5566c22SJunchao Zhang   hx     = 1.0/(PetscReal)(Mx-1);
407c5566c22SJunchao Zhang   hy     = 1.0/(PetscReal)(My-1);
408c5566c22SJunchao Zhang 
409c5566c22SJunchao Zhang   PetscCall(PetscObjectGetComm((PetscObject)snes,&comm));
410c5566c22SJunchao Zhang   /*
411c5566c22SJunchao Zhang      Scatter ghost points to local vector,using the 2-step process
412c5566c22SJunchao Zhang         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
413c5566c22SJunchao Zhang      By placing code between these two statements, computations can be
414c5566c22SJunchao Zhang      done while messages are in transition.
415c5566c22SJunchao Zhang   */
416c5566c22SJunchao Zhang   PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX));
417c5566c22SJunchao Zhang   PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX));
418c5566c22SJunchao Zhang   PetscCall(PetscMatlabEnginePut(PETSC_MATLAB_ENGINE_(comm),(PetscObject)localX));
419c5566c22SJunchao Zhang   PetscCall(PetscMatlabEngineEvaluate(PETSC_MATLAB_ENGINE_(comm),"localF=ex5m(localX,%18.16e,%18.16e,%18.16e)",(double)hx,(double)hy,(double)lambda));
420c5566c22SJunchao Zhang   PetscCall(PetscMatlabEngineGet(PETSC_MATLAB_ENGINE_(comm),(PetscObject)localF));
421c5566c22SJunchao Zhang 
422c5566c22SJunchao Zhang   /*
423c5566c22SJunchao Zhang      Insert values into global vector
424c5566c22SJunchao Zhang   */
425c5566c22SJunchao Zhang   PetscCall(DMLocalToGlobalBegin(da,localF,INSERT_VALUES,F));
426c5566c22SJunchao Zhang   PetscCall(DMLocalToGlobalEnd(da,localF,INSERT_VALUES,F));
427c5566c22SJunchao Zhang   PetscCall(DMRestoreLocalVector(da,&localX));
428c5566c22SJunchao Zhang   PetscCall(DMRestoreLocalVector(da,&localF));
429c5566c22SJunchao Zhang   PetscFunctionReturn(0);
430c5566c22SJunchao Zhang #else
431c5566c22SJunchao Zhang     return 0;                     /* Never called */
432c5566c22SJunchao Zhang #endif
433c5566c22SJunchao Zhang }
434c5566c22SJunchao Zhang 
435c5566c22SJunchao Zhang /* ------------------------------------------------------------------- */
436c5566c22SJunchao Zhang /*
437c5566c22SJunchao Zhang       Applies some sweeps on nonlinear Gauss-Seidel on each process
438c5566c22SJunchao Zhang 
439c5566c22SJunchao Zhang  */
440c5566c22SJunchao Zhang static PetscErrorCode NonlinearGS(SNES snes,Vec X, Vec B, void *ctx)
441c5566c22SJunchao Zhang {
442c5566c22SJunchao Zhang   PetscInt       i,j,k,Mx,My,xs,ys,xm,ym,its,tot_its,sweeps,l;
443c5566c22SJunchao Zhang   PetscReal      lambda,hx,hy,hxdhy,hydhx,sc;
444c5566c22SJunchao Zhang   PetscScalar    **x,**b,bij,F,F0=0,J,u,un,us,ue,eu,uw,uxx,uyy,y;
445c5566c22SJunchao Zhang   PetscReal      atol,rtol,stol;
446c5566c22SJunchao Zhang   DM             da;
447c5566c22SJunchao Zhang   AppCtx         *user;
448c5566c22SJunchao Zhang   Vec            localX,localB;
449c5566c22SJunchao Zhang 
450c5566c22SJunchao Zhang   PetscFunctionBeginUser;
451c5566c22SJunchao Zhang   tot_its = 0;
452c5566c22SJunchao Zhang   PetscCall(SNESNGSGetSweeps(snes,&sweeps));
453c5566c22SJunchao Zhang   PetscCall(SNESNGSGetTolerances(snes,&atol,&rtol,&stol,&its));
454c5566c22SJunchao Zhang   PetscCall(SNESGetDM(snes,&da));
455c5566c22SJunchao Zhang   PetscCall(DMGetApplicationContext(da,&user));
456c5566c22SJunchao Zhang 
457c5566c22SJunchao Zhang   PetscCall(DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE));
458c5566c22SJunchao Zhang 
459c5566c22SJunchao Zhang   lambda = user->param;
460c5566c22SJunchao Zhang   hx     = 1.0/(PetscReal)(Mx-1);
461c5566c22SJunchao Zhang   hy     = 1.0/(PetscReal)(My-1);
462c5566c22SJunchao Zhang   sc     = hx*hy*lambda;
463c5566c22SJunchao Zhang   hxdhy  = hx/hy;
464c5566c22SJunchao Zhang   hydhx  = hy/hx;
465c5566c22SJunchao Zhang 
466c5566c22SJunchao Zhang   PetscCall(DMGetLocalVector(da,&localX));
467c5566c22SJunchao Zhang   if (B) {
468c5566c22SJunchao Zhang     PetscCall(DMGetLocalVector(da,&localB));
469c5566c22SJunchao Zhang   }
470c5566c22SJunchao Zhang   for (l=0; l<sweeps; l++) {
471c5566c22SJunchao Zhang 
472c5566c22SJunchao Zhang     PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX));
473c5566c22SJunchao Zhang     PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX));
474c5566c22SJunchao Zhang     if (B) {
475c5566c22SJunchao Zhang       PetscCall(DMGlobalToLocalBegin(da,B,INSERT_VALUES,localB));
476c5566c22SJunchao Zhang       PetscCall(DMGlobalToLocalEnd(da,B,INSERT_VALUES,localB));
477c5566c22SJunchao Zhang     }
478c5566c22SJunchao Zhang     /*
479c5566c22SJunchao Zhang      Get a pointer to vector data.
480c5566c22SJunchao Zhang      - For default PETSc vectors, VecGetArray() returns a pointer to
481c5566c22SJunchao Zhang      the data array.  Otherwise, the routine is implementation dependent.
482c5566c22SJunchao Zhang      - You MUST call VecRestoreArray() when you no longer need access to
483c5566c22SJunchao Zhang      the array.
484c5566c22SJunchao Zhang      */
485c5566c22SJunchao Zhang     PetscCall(DMDAVecGetArray(da,localX,&x));
486c5566c22SJunchao Zhang     if (B) PetscCall(DMDAVecGetArray(da,localB,&b));
487c5566c22SJunchao Zhang     /*
488c5566c22SJunchao Zhang      Get local grid boundaries (for 2-dimensional DMDA):
489c5566c22SJunchao Zhang      xs, ys   - starting grid indices (no ghost points)
490c5566c22SJunchao Zhang      xm, ym   - widths of local grid (no ghost points)
491c5566c22SJunchao Zhang      */
492c5566c22SJunchao Zhang     PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
493c5566c22SJunchao Zhang 
494c5566c22SJunchao Zhang     for (j=ys; j<ys+ym; j++) {
495c5566c22SJunchao Zhang       for (i=xs; i<xs+xm; i++) {
496c5566c22SJunchao Zhang         if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
497c5566c22SJunchao Zhang           /* boundary conditions are all zero Dirichlet */
498c5566c22SJunchao Zhang           x[j][i] = 0.0;
499c5566c22SJunchao Zhang         } else {
500c5566c22SJunchao Zhang           if (B) bij = b[j][i];
501c5566c22SJunchao Zhang           else   bij = 0.;
502c5566c22SJunchao Zhang 
503c5566c22SJunchao Zhang           u  = x[j][i];
504c5566c22SJunchao Zhang           un = x[j-1][i];
505c5566c22SJunchao Zhang           us = x[j+1][i];
506c5566c22SJunchao Zhang           ue = x[j][i-1];
507c5566c22SJunchao Zhang           uw = x[j][i+1];
508c5566c22SJunchao Zhang 
509c5566c22SJunchao Zhang           for (k=0; k<its; k++) {
510c5566c22SJunchao Zhang             eu  = PetscExpScalar(u);
511c5566c22SJunchao Zhang             uxx = (2.0*u - ue - uw)*hydhx;
512c5566c22SJunchao Zhang             uyy = (2.0*u - un - us)*hxdhy;
513c5566c22SJunchao Zhang             F   = uxx + uyy - sc*eu - bij;
514c5566c22SJunchao Zhang             if (k == 0) F0 = F;
515c5566c22SJunchao Zhang             J  = 2.0*(hydhx + hxdhy) - sc*eu;
516c5566c22SJunchao Zhang             y  = F/J;
517c5566c22SJunchao Zhang             u -= y;
518c5566c22SJunchao Zhang             tot_its++;
519c5566c22SJunchao Zhang 
520c5566c22SJunchao Zhang             if (atol > PetscAbsReal(PetscRealPart(F)) ||
521c5566c22SJunchao Zhang                 rtol*PetscAbsReal(PetscRealPart(F0)) > PetscAbsReal(PetscRealPart(F)) ||
522c5566c22SJunchao Zhang                 stol*PetscAbsReal(PetscRealPart(u)) > PetscAbsReal(PetscRealPart(y))) {
523c5566c22SJunchao Zhang               break;
524c5566c22SJunchao Zhang             }
525c5566c22SJunchao Zhang           }
526c5566c22SJunchao Zhang           x[j][i] = u;
527c5566c22SJunchao Zhang         }
528c5566c22SJunchao Zhang       }
529c5566c22SJunchao Zhang     }
530c5566c22SJunchao Zhang     /*
531c5566c22SJunchao Zhang      Restore vector
532c5566c22SJunchao Zhang      */
533c5566c22SJunchao Zhang     PetscCall(DMDAVecRestoreArray(da,localX,&x));
534c5566c22SJunchao Zhang     PetscCall(DMLocalToGlobalBegin(da,localX,INSERT_VALUES,X));
535c5566c22SJunchao Zhang     PetscCall(DMLocalToGlobalEnd(da,localX,INSERT_VALUES,X));
536c5566c22SJunchao Zhang   }
537c5566c22SJunchao Zhang   PetscCall(PetscLogFlops(tot_its*(21.0)));
538c5566c22SJunchao Zhang   PetscCall(DMRestoreLocalVector(da,&localX));
539c5566c22SJunchao Zhang   if (B) {
540c5566c22SJunchao Zhang     PetscCall(DMDAVecRestoreArray(da,localB,&b));
541c5566c22SJunchao Zhang     PetscCall(DMRestoreLocalVector(da,&localB));
542c5566c22SJunchao Zhang   }
543c5566c22SJunchao Zhang   PetscFunctionReturn(0);
544c5566c22SJunchao Zhang }
545c5566c22SJunchao Zhang 
546c5566c22SJunchao Zhang int main(int argc,char **argv)
547c5566c22SJunchao Zhang {
548c5566c22SJunchao Zhang   SNES           snes;                         /* nonlinear solver */
549c5566c22SJunchao Zhang   Vec            x;                            /* solution vector */
550c5566c22SJunchao Zhang   AppCtx         user;                         /* user-defined work context */
551c5566c22SJunchao Zhang   PetscInt       its;                          /* iterations for convergence */
552c5566c22SJunchao Zhang   PetscReal      bratu_lambda_max = 6.81;
553c5566c22SJunchao Zhang   PetscReal      bratu_lambda_min = 0.;
554c5566c22SJunchao Zhang   PetscInt       MMS              = 1;
555c5566c22SJunchao Zhang   PetscBool      flg              = PETSC_FALSE,setMMS;
556c5566c22SJunchao Zhang   DM             da;
557c5566c22SJunchao Zhang   Vec            r                = NULL;
558c5566c22SJunchao Zhang   KSP            ksp;
559c5566c22SJunchao Zhang   PetscInt       lits,slits;
560c5566c22SJunchao Zhang   PetscBool      useKokkos        = PETSC_FALSE;
561c5566c22SJunchao Zhang 
562c5566c22SJunchao Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
563c5566c22SJunchao Zhang      Initialize program
564c5566c22SJunchao Zhang      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
565c5566c22SJunchao Zhang 
566*327415f7SBarry Smith   PetscFunctionBeginUser;
567c5566c22SJunchao Zhang   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
568c5566c22SJunchao Zhang 
569c5566c22SJunchao Zhang   PetscCall(PetscOptionsGetBool(NULL,NULL,"-use_kokkos",&useKokkos,NULL));
570c5566c22SJunchao Zhang 
571c5566c22SJunchao Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
572c5566c22SJunchao Zhang      Initialize problem parameters
573c5566c22SJunchao Zhang   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
574c5566c22SJunchao Zhang   user.param = 6.0;
575c5566c22SJunchao Zhang   PetscCall(PetscOptionsGetReal(NULL,NULL,"-par",&user.param,NULL));
576c5566c22SJunchao Zhang   PetscCheck(user.param <= bratu_lambda_max && user.param >= bratu_lambda_min,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lambda, %g, is out of range, [%g, %g]", (double)user.param, (double)bratu_lambda_min, (double)bratu_lambda_max);
577c5566c22SJunchao Zhang   PetscCall(PetscOptionsGetInt(NULL,NULL,"-mms",&MMS,&setMMS));
578c5566c22SJunchao Zhang   if (MMS == 3) {
579c5566c22SJunchao Zhang     PetscInt mPar = 2, nPar = 1;
580c5566c22SJunchao Zhang     PetscCall(PetscOptionsGetInt(NULL,NULL,"-m_par",&mPar,NULL));
581c5566c22SJunchao Zhang     PetscCall(PetscOptionsGetInt(NULL,NULL,"-n_par",&nPar,NULL));
582c5566c22SJunchao Zhang     user.m = PetscPowInt(2,mPar);
583c5566c22SJunchao Zhang     user.n = PetscPowInt(2,nPar);
584c5566c22SJunchao Zhang   }
585c5566c22SJunchao Zhang 
586c5566c22SJunchao Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
587c5566c22SJunchao Zhang      Create nonlinear solver context
588c5566c22SJunchao Zhang      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
589c5566c22SJunchao Zhang   PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes));
590c5566c22SJunchao Zhang   PetscCall(SNESSetCountersReset(snes,PETSC_FALSE));
591c5566c22SJunchao Zhang   PetscCall(SNESSetNGS(snes, NonlinearGS, NULL));
592c5566c22SJunchao Zhang 
593c5566c22SJunchao Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
594c5566c22SJunchao Zhang      Create distributed array (DMDA) to manage parallel grid and vectors
595c5566c22SJunchao Zhang   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
596c5566c22SJunchao Zhang   PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,4,4,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da));
597c5566c22SJunchao Zhang   PetscCall(DMSetFromOptions(da));
598c5566c22SJunchao Zhang   PetscCall(DMSetUp(da));
599c5566c22SJunchao Zhang   PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
600c5566c22SJunchao Zhang   PetscCall(DMSetApplicationContext(da,&user));
601c5566c22SJunchao Zhang   PetscCall(SNESSetDM(snes,da));
602c5566c22SJunchao Zhang   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
603c5566c22SJunchao Zhang      Extract global vectors from DMDA; then duplicate for remaining
604c5566c22SJunchao Zhang      vectors that are the same types
605c5566c22SJunchao Zhang    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
606c5566c22SJunchao Zhang   PetscCall(DMCreateGlobalVector(da,&x));
607c5566c22SJunchao Zhang 
608c5566c22SJunchao Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
609c5566c22SJunchao Zhang      Set local function evaluation routine
610c5566c22SJunchao Zhang   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
611c5566c22SJunchao Zhang   switch (MMS) {
612c5566c22SJunchao Zhang   case 0: user.mms_solution = ZeroBCSolution; user.mms_forcing = NULL; break;
613c5566c22SJunchao Zhang   case 1: user.mms_solution = MMSSolution1; user.mms_forcing = MMSForcing1; break;
614c5566c22SJunchao Zhang   case 2: user.mms_solution = MMSSolution2; user.mms_forcing = MMSForcing2; break;
615c5566c22SJunchao Zhang   case 3: user.mms_solution = MMSSolution3; user.mms_forcing = MMSForcing3; break;
616c5566c22SJunchao Zhang   case 4: user.mms_solution = MMSSolution4; user.mms_forcing = MMSForcing4; break;
617c5566c22SJunchao Zhang   default: SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Unknown MMS type %" PetscInt_FMT,MMS);
618c5566c22SJunchao Zhang   }
619c5566c22SJunchao Zhang 
620c5566c22SJunchao Zhang   if (useKokkos) {
621c5566c22SJunchao Zhang     PetscCheck(MMS == 1,PETSC_COMM_WORLD,PETSC_ERR_USER,"FormFunctionLocalVec_Kokkos only works with MMS 1");
622c5566c22SJunchao Zhang     PetscCall(DMDASNESSetFunctionLocalVec(da,INSERT_VALUES,(DMDASNESFunctionVec)FormFunctionLocalVec,&user));
623c5566c22SJunchao Zhang   } else {
624c5566c22SJunchao Zhang     PetscCall(DMDASNESSetFunctionLocal(da,INSERT_VALUES,(DMDASNESFunction)FormFunctionLocal,&user));
625c5566c22SJunchao Zhang   }
626c5566c22SJunchao Zhang 
627c5566c22SJunchao Zhang   PetscCall(PetscOptionsGetBool(NULL,NULL,"-fd",&flg,NULL));
628c5566c22SJunchao Zhang   if (!flg) {
629c5566c22SJunchao Zhang     if (useKokkos) PetscCall(DMDASNESSetJacobianLocalVec(da,(DMDASNESJacobianVec)FormJacobianLocalVec,&user));
630c5566c22SJunchao Zhang     else PetscCall(DMDASNESSetJacobianLocal(da,(DMDASNESJacobian)FormJacobianLocal,&user));
631c5566c22SJunchao Zhang   }
632c5566c22SJunchao Zhang 
633c5566c22SJunchao Zhang   PetscCall(PetscOptionsGetBool(NULL,NULL,"-obj",&flg,NULL));
634c5566c22SJunchao Zhang   if (flg) {
635c5566c22SJunchao Zhang     if (useKokkos) PetscCall(DMDASNESSetObjectiveLocalVec(da,(DMDASNESObjectiveVec)FormObjectiveLocalVec,&user));
636c5566c22SJunchao Zhang     else PetscCall(DMDASNESSetObjectiveLocal(da,(DMDASNESObjective)FormObjectiveLocal,&user));
637c5566c22SJunchao Zhang   }
638c5566c22SJunchao Zhang 
639c5566c22SJunchao Zhang   if (PetscDefined(HAVE_MATLAB_ENGINE)) {
640c5566c22SJunchao Zhang     PetscBool matlab_function = PETSC_FALSE;
641c5566c22SJunchao Zhang     PetscCall(PetscOptionsGetBool(NULL,NULL,"-matlab_function",&matlab_function,0));
642c5566c22SJunchao Zhang     if (matlab_function) {
643c5566c22SJunchao Zhang       PetscCall(VecDuplicate(x,&r));
644c5566c22SJunchao Zhang       PetscCall(SNESSetFunction(snes,r,FormFunctionMatlab,&user));
645c5566c22SJunchao Zhang     }
646c5566c22SJunchao Zhang   }
647c5566c22SJunchao Zhang 
648c5566c22SJunchao Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
649c5566c22SJunchao Zhang      Customize nonlinear solver; set runtime options
650c5566c22SJunchao Zhang    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
651c5566c22SJunchao Zhang   PetscCall(SNESSetFromOptions(snes));
652c5566c22SJunchao Zhang 
653c5566c22SJunchao Zhang   PetscCall(FormInitialGuess(da,&user,x));
654c5566c22SJunchao Zhang 
655c5566c22SJunchao Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
656c5566c22SJunchao Zhang      Solve nonlinear system
657c5566c22SJunchao Zhang      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
658c5566c22SJunchao Zhang   PetscCall(SNESSolve(snes,NULL,x));
659c5566c22SJunchao Zhang   PetscCall(SNESGetIterationNumber(snes,&its));
660c5566c22SJunchao Zhang 
661c5566c22SJunchao Zhang   PetscCall(SNESGetLinearSolveIterations(snes,&slits));
662c5566c22SJunchao Zhang   PetscCall(SNESGetKSP(snes,&ksp));
663c5566c22SJunchao Zhang   PetscCall(KSPGetTotalIterations(ksp,&lits));
664c5566c22SJunchao Zhang   PetscCheck(lits == slits,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Number of total linear iterations reported by SNES %" PetscInt_FMT " does not match reported by KSP %" PetscInt_FMT,slits,lits);
665c5566c22SJunchao Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
666c5566c22SJunchao Zhang    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
667c5566c22SJunchao Zhang 
668c5566c22SJunchao Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
669c5566c22SJunchao Zhang      If using MMS, check the l_2 error
670c5566c22SJunchao Zhang    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
671c5566c22SJunchao Zhang   if (setMMS) {
672c5566c22SJunchao Zhang     Vec       e;
673c5566c22SJunchao Zhang     PetscReal errorl2, errorinf;
674c5566c22SJunchao Zhang     PetscInt  N;
675c5566c22SJunchao Zhang 
676c5566c22SJunchao Zhang     PetscCall(VecDuplicate(x, &e));
677c5566c22SJunchao Zhang     PetscCall(PetscObjectViewFromOptions((PetscObject) x, NULL, "-sol_view"));
678c5566c22SJunchao Zhang     PetscCall(FormExactSolution(da, &user, e));
679c5566c22SJunchao Zhang     PetscCall(PetscObjectViewFromOptions((PetscObject) e, NULL, "-exact_view"));
680c5566c22SJunchao Zhang     PetscCall(VecAXPY(e, -1.0, x));
681c5566c22SJunchao Zhang     PetscCall(PetscObjectViewFromOptions((PetscObject) e, NULL, "-error_view"));
682c5566c22SJunchao Zhang     PetscCall(VecNorm(e, NORM_2, &errorl2));
683c5566c22SJunchao Zhang     PetscCall(VecNorm(e, NORM_INFINITY, &errorinf));
684c5566c22SJunchao Zhang     PetscCall(VecGetSize(e, &N));
685c5566c22SJunchao Zhang     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "N: %" PetscInt_FMT " error L2 %g inf %g\n", N, (double)(errorl2/PetscSqrtReal((PetscReal)N)), (double) errorinf));
686c5566c22SJunchao Zhang     PetscCall(VecDestroy(&e));
687c5566c22SJunchao Zhang     PetscCall(PetscLogEventSetDof(SNES_Solve, 0, N));
688c5566c22SJunchao Zhang     PetscCall(PetscLogEventSetError(SNES_Solve, 0, errorl2/PetscSqrtReal(N)));
689c5566c22SJunchao Zhang   }
690c5566c22SJunchao Zhang 
691c5566c22SJunchao Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
692c5566c22SJunchao Zhang      Free work space.  All PETSc objects should be destroyed when they
693c5566c22SJunchao Zhang      are no longer needed.
694c5566c22SJunchao Zhang    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
695c5566c22SJunchao Zhang   PetscCall(VecDestroy(&r));
696c5566c22SJunchao Zhang   PetscCall(VecDestroy(&x));
697c5566c22SJunchao Zhang   PetscCall(SNESDestroy(&snes));
698c5566c22SJunchao Zhang   PetscCall(DMDestroy(&da));
699c5566c22SJunchao Zhang   PetscCall(PetscFinalize());
700c5566c22SJunchao Zhang   return 0;
701c5566c22SJunchao Zhang }
702c5566c22SJunchao Zhang 
703c5566c22SJunchao Zhang /*TEST
704c5566c22SJunchao Zhang   build:
705c5566c22SJunchao Zhang     requires: kokkos_kernels
706c5566c22SJunchao Zhang     depends: ex55k.kokkos.cxx
707c5566c22SJunchao Zhang 
708c5566c22SJunchao Zhang   testset:
709c5566c22SJunchao Zhang     output_file: output/ex55_asm_0.out
710c5566c22SJunchao Zhang     requires: !single
711c5566c22SJunchao Zhang     args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 2 -pc_asm_overlap 0 -pc_asm_local_type additive -sub_pc_type lu
712c5566c22SJunchao Zhang     filter: grep -v "type"
713c5566c22SJunchao Zhang 
714c5566c22SJunchao Zhang     test:
715c5566c22SJunchao Zhang       suffix: asm_0
716c5566c22SJunchao Zhang     test:
717c5566c22SJunchao Zhang       suffix: asm_0_kok
718c5566c22SJunchao Zhang       args: -use_kokkos 1 -dm_mat_type aijkokkos -dm_vec_type kokkos
719c5566c22SJunchao Zhang 
720c5566c22SJunchao Zhang TEST*/
721