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