/* Evaluate analytical Jacobian matrix. When debugging Jacobian: ./marine -noscale -dx 60000 -mat_view -exactinit |tail -n 18 versus ./marine -noscale -dx 60000 -mat_view -exactinit -snes_fd |tail -n 18 and do ./marine -noscale -dx 60000 -snes_type test ./marine -noscale -dx 60000 -snes_type test -snes_test_display */ PetscErrorCode JacobianMatrixLocal(DMDALocalInfo *info, Node *Hu, Mat A, Mat jac, MatStructure *str, AppCtx *user) { PetscErrorCode ierr; PetscReal rg = user->rho * user->g, dx = user->dx, Hg = user->zocean * (user->rhow / user->rho); PetscReal *Bstag, Hl, ul, Fl, Fr, Gl, Gr, hr, hl; PetscInt i, Mx = info->mx; Vec locBstag; PetscReal v[6]; MatStencil row,col[6]; PetscFunctionBegin; ierr = DMGetLocalVector(user->stagda,&locBstag);CHKERRQ(ierr); /* do NOT destroy it */ ierr = DMGlobalToLocalBegin(user->stagda,user->Bstag,INSERT_VALUES,locBstag); CHKERRQ(ierr); ierr = DMGlobalToLocalEnd(user->stagda,user->Bstag,INSERT_VALUES,locBstag); CHKERRQ(ierr); ierr = DMDAVecGetArray(user->stagda,locBstag,&Bstag);CHKERRQ(ierr); for (i=info->xs; i<info->xs+info->xm; i++) { /* MASS CONT */ row.i = i; row.c = 0; /* "row.c=0" means H component of f */ if (i == 0) { col[0].i = i; col[0].c = 0; v[0] = 1.0; ierr = MatSetValuesStencil(jac,1,&row,1,col,v,INSERT_VALUES);CHKERRQ(ierr); } else { col[0].i = i; col[0].c = 0; v[0] = Hu[i].u / dx; col[1].i = i; col[1].c = 1; v[1] = Hu[i].H / dx; if (i > 1) { col[2].i = i-1; col[2].c = 0; v[2] = - Hu[i-1].u / dx; col[3].i = i-1; col[3].c = 1; v[3] = - Hu[i-1].H / dx; } ierr = MatSetValuesStencil(jac,1,&row,(i > 1) ? 4 : 2,col,v,INSERT_VALUES);CHKERRQ(ierr); } /* done with MASS CONT */ /* SSA */ row.i = i; row.c = 1; /* "row.c=1" means u component of f */ if (i == 0) { col[0].i = i; col[0].c = 1; v[0] = 1.0; ierr = MatSetValuesStencil(jac,1,&row,1,col,v,INSERT_VALUES);CHKERRQ(ierr); } else if (i == Mx - 1) { col[0].i = i-1; col[0].c = 0; v[0] = 0.25 * user->omega * rg; col[1].i = i; col[1].c = 0; v[1] = 0.25 * user->omega * rg; Gl = GSR(dx,user->epsilon,user->n,Hu[i-1].u,Hu[i].u); col[2].i = i-1; col[2].c = 1; v[2] = 2.0 * Bstag[i-1] * Gl / dx; col[3].i = i; col[3].c = 1; v[3] = - 2.0 * Bstag[i-1] * Gl / dx; ierr = MatSetValuesStencil(jac,1,&row,4,col,v,INSERT_VALUES);CHKERRQ(ierr); } else { ul = (i == 1) ? user->ua : Hu[i-1].u; Hl = (i == 1) ? user->Ha : Hu[i-1].H; Fl = GetFSR(dx,user->epsilon,user->n,ul,Hu[i].u); Fr = GetFSR(dx,user->epsilon,user->n,Hu[i].u,Hu[i+1].u); Gl = GSR(dx,user->epsilon,user->n,ul,Hu[i].u); Gr = GSR(dx,user->epsilon,user->n,Hu[i].u,Hu[i+1].u); // df^u / dH col[0].i = i-1; col[0].c = 0; if (i == 1) v[0] = 0.0; else { v[0] = ( - Bstag[i-1] * Fl ) / dx - rg * Hu[i].H * ( - dsurfdH(Hl,Hg,user->omega,0.0) ) / (2.0 * dx); } col[1].i = i; col[1].c = 0; hl = getsurf(Hl,Hg,user->omega,user->zocean,0.0); hr = getsurf(Hu[i+1].H,Hg,user->omega,user->zocean,0.0); v[1] = ( Bstag[i] * Fr - Bstag[i-1] * Fl ) / dx - user->k * rg * GLREG(Hu[i].H,Hg,0.0) * Hu[i].u - rg * (hr - hl) / (2.0 * dx); col[2].i = i+1; col[2].c = 0; v[2] = ( Bstag[i] * Fr ) / dx - rg * Hu[i].H * ( dsurfdH(Hu[i+1].H,Hg,user->omega,0.0) ) / (2.0 * dx); // df^u / du col[3].i = i-1; col[3].c = 1; v[3] = ( - Bstag[i-1] * (Hl + Hu[i].H) * (-1.0/dx) * Gl ) / dx; col[4].i = i; col[4].c = 1; v[4] = ( Bstag[i] * (Hu[i].H + Hu[i+1].H) * (-1.0/dx) * Gr - Bstag[i-1] * (Hl + Hu[i].H) * (1.0/dx) * Gl ) / dx - user->k * rg * Hu[i].H * GLREG(Hu[i].H,Hg,0.0); col[5].i = i+1; col[5].c = 1; v[5] = ( Bstag[i] * (Hu[i].H + Hu[i+1].H) * (1.0/dx) * Gr ) / dx; ierr = MatSetValuesStencil(jac,1,&row,6,col,v,INSERT_VALUES);CHKERRQ(ierr); } } ierr = DMDAVecRestoreArray(user->stagda,locBstag,&Bstag);CHKERRQ(ierr); ierr = DMRestoreLocalVector(user->stagda,&locBstag);CHKERRQ(ierr); /* assemble matrix, using the 2-step process */ ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); if (A != jac) { ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); } *str = SAME_NONZERO_PATTERN; /* tell matrix we will never add a new nonzero location; if we do then gives error */ ierr = MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); PetscFunctionReturn(0); }
static PetscErrorCode BodFunctionLocal(DALocalInfo *info, Node *Hu, Node *f, AppCtx *user) { PetscErrorCode ierr; PetscReal rg = user->rho * user->g, dx = user->dx; PetscScalar *M, *Bstag, *beta, duH, ul, u, ur, dHdx, Fl, Fr, Tl, Tr; PetscInt i, Mx = info->mx; Vec locBstag; PetscFunctionBegin; /* we need stencil width on Bstag (but not for M, beta) */ ierr = DAGetLocalVector(user->scalarda,&locBstag);CHKERRQ(ierr); /* do NOT destroy it */ ierr = DAGlobalToLocalBegin(user->scalarda,user->Bstag,INSERT_VALUES,locBstag); CHKERRQ(ierr); ierr = DAGlobalToLocalEnd(user->scalarda,user->Bstag,INSERT_VALUES,locBstag); CHKERRQ(ierr); ierr = DAVecGetArray(user->scalarda,locBstag,&Bstag);CHKERRQ(ierr); ierr = DAVecGetArray(user->scalarda,user->M,&M);CHKERRQ(ierr); ierr = DAVecGetArray(user->scalarda,user->beta,&beta);CHKERRQ(ierr); for (i = info->xs; i < info->xs + info->xm; i++) { /* MASS CONT */ if (i == 0) { /* residual at left-most point is Dirichlet cond. */ f[0].H = Hu[0].H - user->H0; } else { /* centered difference IS UNSTABLE: duH = Hu[i+1].u * Hu[i+1].H - ( (i == 1) ? 0.0 : Hu[i-1].u * Hu[i-1].H ); duH *= 0.5; */ if (user->upwind1) { /* 1st-order upwind; leftward difference because u > 0 (because dH/dx < 0) */ duH = Hu[i].u * Hu[i].H - ( (i == 1) ? 0.0 : Hu[i-1].u * Hu[i-1].H ); } else { /* 2nd-order upwind; see Beam-Warming discussion in R. LeVeque, "Finite Volume ..." */ if (i == 1) { /* use PDE M - (uH)_x = 0 to get quadratic poly, then diff that */ duH = - dx * M[0] + 2.0 * Hu[1].u * Hu[1].H; } else { duH = 3.0 * Hu[i].u * Hu[i].H - 4.0 * Hu[i-1].u * Hu[i-1].H; /* if i == 2 then u=0 so uH = 0 */ if (i >= 3) duH += Hu[i-2].u * Hu[i-2].H; duH *= 0.5; } } f[i].H = dx * M[i] - duH; } /* SSA */ if (i == 0) { /* residual at left-most point is Dirichlet cond. */ f[0].u = Hu[0].u - 0.0; } else { /* residual: SSA eqn */ /* consecutive values of u */ ul = (i == 1) ? 0.0 : Hu[i-1].u; u = Hu[i].u; ur = (i == Mx-1) ? -1.1e30 : Hu[i+1].u; /* surface slope */ if (i == 1) { dHdx = (Hu[i+1].H - user->H0) / (2.0 * dx); } else if (i == Mx-1) { /* nearly 2nd-order global convergence seems to occur even with this: dHdx = (Hu[i].H - Hu[i-1].H) / dx; */ dHdx = (3.0*Hu[i].H - 4.0*Hu[i-1].H + Hu[i-2].H) / (2.0 * dx); } else { /* generic case */ dHdx = (Hu[i+1].H - Hu[i-1].H) / (2.0 * dx); } /* vertically-integrated longitudinal stress */ Fl = GetFSR(dx,user->epsilon,user->n, ul,u); if (i == Mx-1) { Tl = 2.0 * (Hu[i-1].H + Hu[i].H) * Bstag[i-1] * Fl; Tr = (1.0 - user->rho / user->rhow) * user->rho * user->g * Hu[i].H * Hu[i].H; /* exact value: Tr = 2.0 * user->Txc; */ } else { Fr = GetFSR(dx,user->epsilon,user->n, u,ur); Tl = (Hu[i-1].H + Hu[i].H) * Bstag[i-1] * Fl; Tr = (Hu[i].H + Hu[i+1].H) * Bstag[i] * Fr; } f[i].u = (Tr - Tl) - dx * beta[i] * u - dx * rg * Hu[i].H * dHdx; /* SSA */ } } ierr = DAVecRestoreArray(user->scalarda,locBstag,&Bstag);CHKERRQ(ierr); ierr = DAVecRestoreArray(user->scalarda,user->M,&M);CHKERRQ(ierr); ierr = DAVecRestoreArray(user->scalarda,user->beta,&beta);CHKERRQ(ierr); ierr = DARestoreLocalVector(user->scalarda,&locBstag);CHKERRQ(ierr); PetscFunctionReturn(0); }
/* Evaluate residual f at current iterate Hu. */ PetscErrorCode FunctionLocal(DMDALocalInfo *info, Node *Hu, Node *f, AppCtx *user) { PetscErrorCode ierr; PetscReal rg = user->rho * user->g, dx = user->dx, Hg = user->zocean * (user->rhow / user->rho); PetscReal *Mstag, *Bstag, duH, ul, Fl, Fr, Tl, Tr, Hl, hl, hr, beta; PetscInt i, Mx = info->mx; Vec locBstag, locMstag; PetscFunctionBegin; /* we need stencil width on Bstag and Mstag */ ierr = DMGetLocalVector(user->stagda,&locBstag);CHKERRQ(ierr); /* do NOT destroy it */ ierr = DMGlobalToLocalBegin(user->stagda,user->Bstag,INSERT_VALUES,locBstag); CHKERRQ(ierr); ierr = DMGlobalToLocalEnd(user->stagda,user->Bstag,INSERT_VALUES,locBstag); CHKERRQ(ierr); ierr = DMGetLocalVector(user->stagda,&locMstag);CHKERRQ(ierr); /* do NOT destroy it */ ierr = DMGlobalToLocalBegin(user->stagda,user->Mstag,INSERT_VALUES,locMstag); CHKERRQ(ierr); ierr = DMGlobalToLocalEnd(user->stagda,user->Mstag,INSERT_VALUES,locMstag); CHKERRQ(ierr); ierr = DMDAVecGetArray(user->stagda,locBstag,&Bstag);CHKERRQ(ierr); ierr = DMDAVecGetArray(user->stagda,locMstag,&Mstag);CHKERRQ(ierr); for (i = info->xs; i < info->xs + info->xm; i++) { /* MASS CONT */ if (i == 0) { /* residual at left-most point is Dirichlet cond. */ f[0].H = Hu[0].H - user->Ha; } else { /* centered */ if (i == 1) { duH = Hu[i].u * Hu[i].H - user->ua * user->Ha; } else { duH = Hu[i].u * Hu[i].H - Hu[i-1].u * Hu[i-1].H; } /**** residual for mass continuity ****/ f[i].H = duH/dx - Mstag[i-1]; } /* SSA */ if (i == 0) { /* residual at left-most point is Dirichlet condition */ f[0].u = Hu[0].u - user->ua; } else if (i == Mx-1) { /* residual at right-most point is calving condition */ Fr = GetFSR(dx,user->epsilon,user->n, Hu[i-1].u,Hu[i].u); f[i].u = 0.25 * user->omega * rg * (Hu[i-1].H + Hu[i].H) - 2.0 * Bstag[i-1] * Fr; } else { /* T = vertically-integrated longitudinal stress */ ul = (i == 1) ? user->ua : Hu[i-1].u; Hl = (i == 1) ? user->Ha : Hu[i-1].H; Fl = GetFSR(dx,user->epsilon,user->n, ul, Hu[i].u); Fr = GetFSR(dx,user->epsilon,user->n, Hu[i].u,Hu[i+1].u); Tl = Bstag[i-1] * (Hl + Hu[i].H ) * Fl; Tr = Bstag[i] * (Hu[i].H + Hu[i+1].H) * Fr; /* sliding coefficient */ beta = user->k * rg * Hu[i].H * GLREG(Hu[i].H,Hg,0.0); /* surface elevations */ hl = getsurf(Hl,Hg,user->omega,user->zocean,0.0); hr = getsurf(Hu[i+1].H,Hg,user->omega,user->zocean,0.0); /**** residual for SSA ****/ f[i].u = (Tr - Tl) / dx - beta * Hu[i].u - rg * Hu[i].H * (hr - hl) / (2.0 * dx); } } ierr = DMDAVecRestoreArray(user->stagda,locBstag,&Bstag);CHKERRQ(ierr); ierr = DMDAVecRestoreArray(user->stagda,locMstag,&Mstag);CHKERRQ(ierr); ierr = DMRestoreLocalVector(user->stagda,&locBstag);CHKERRQ(ierr); ierr = DMRestoreLocalVector(user->stagda,&locMstag);CHKERRQ(ierr); PetscFunctionReturn(0); }