1c4762a1bSJed Brown #include <petsctao.h> 2c4762a1bSJed Brown 3c4762a1bSJed Brown typedef struct { 4c4762a1bSJed Brown PetscInt n; /* Number of variables */ 5c4762a1bSJed Brown PetscInt m; /* Number of constraints */ 6c4762a1bSJed Brown PetscInt mx; /* grid points in each direction */ 7c4762a1bSJed Brown PetscInt nt; /* Number of time steps */ 8c4762a1bSJed Brown PetscInt ndata; /* Number of data points per sample */ 9c4762a1bSJed Brown IS s_is; 10c4762a1bSJed Brown IS d_is; 11c4762a1bSJed Brown VecScatter state_scatter; 12c4762a1bSJed Brown VecScatter design_scatter; 13c4762a1bSJed Brown VecScatter *uxi_scatter, *uyi_scatter, *ux_scatter, *uy_scatter, *ui_scatter; 14c4762a1bSJed Brown VecScatter *yi_scatter; 15c4762a1bSJed Brown 16c4762a1bSJed Brown Mat Js, Jd, JsBlockPrec, JsInv, JsBlock; 17c4762a1bSJed Brown PetscBool jformed, c_formed; 18c4762a1bSJed Brown 19c4762a1bSJed Brown PetscReal alpha; /* Regularization parameter */ 20c4762a1bSJed Brown PetscReal gamma; 21c4762a1bSJed Brown PetscReal ht; /* Time step */ 22c4762a1bSJed Brown PetscReal T; /* Final time */ 23c4762a1bSJed Brown Mat Q, QT; 24c4762a1bSJed Brown Mat L, LT; 25c4762a1bSJed Brown Mat Div, Divwork, Divxy[2]; 26c4762a1bSJed Brown Mat Grad, Gradxy[2]; 27c4762a1bSJed Brown Mat M; 28c4762a1bSJed Brown Mat *C, *Cwork; 29c4762a1bSJed Brown /* Mat Hs,Hd,Hsd; */ 30c4762a1bSJed Brown Vec q; 31c4762a1bSJed Brown Vec ur; /* reference */ 32c4762a1bSJed Brown 33c4762a1bSJed Brown Vec d; 34c4762a1bSJed Brown Vec dwork; 35c4762a1bSJed Brown 36c4762a1bSJed Brown Vec y; /* state variables */ 37c4762a1bSJed Brown Vec ywork; 38c4762a1bSJed Brown Vec ytrue; 39c4762a1bSJed Brown Vec *yi, *yiwork, *ziwork; 40c4762a1bSJed Brown Vec *uxi, *uyi, *uxiwork, *uyiwork, *ui, *uiwork; 41c4762a1bSJed Brown 42c4762a1bSJed Brown Vec u; /* design variables */ 43c4762a1bSJed Brown Vec uwork, vwork; 44c4762a1bSJed Brown Vec utrue; 45c4762a1bSJed Brown 46c4762a1bSJed Brown Vec js_diag; 47c4762a1bSJed Brown 48c4762a1bSJed Brown Vec c; /* constraint vector */ 49c4762a1bSJed Brown Vec cwork; 50c4762a1bSJed Brown 51c4762a1bSJed Brown Vec lwork; 52c4762a1bSJed Brown 53c4762a1bSJed Brown KSP solver; 54c4762a1bSJed Brown PC prec; 55c4762a1bSJed Brown PetscInt block_index; 56c4762a1bSJed Brown 57c4762a1bSJed Brown PetscInt ksp_its; 58c4762a1bSJed Brown PetscInt ksp_its_initial; 59c4762a1bSJed Brown } AppCtx; 60c4762a1bSJed Brown 61c4762a1bSJed Brown PetscErrorCode FormFunction(Tao, Vec, PetscReal *, void *); 62c4762a1bSJed Brown PetscErrorCode FormGradient(Tao, Vec, Vec, void *); 63c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *); 64c4762a1bSJed Brown PetscErrorCode FormJacobianState(Tao, Vec, Mat, Mat, Mat, void *); 65c4762a1bSJed Brown PetscErrorCode FormJacobianDesign(Tao, Vec, Mat, void *); 66c4762a1bSJed Brown PetscErrorCode FormConstraints(Tao, Vec, Vec, void *); 67c4762a1bSJed Brown PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void *); 68c4762a1bSJed Brown PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat); 69c4762a1bSJed Brown PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat); 70c4762a1bSJed Brown PetscErrorCode HyperbolicInitialize(AppCtx *user); 71c4762a1bSJed Brown PetscErrorCode HyperbolicDestroy(AppCtx *user); 72c4762a1bSJed Brown PetscErrorCode HyperbolicMonitor(Tao, void *); 73c4762a1bSJed Brown 74c4762a1bSJed Brown PetscErrorCode StateMatMult(Mat, Vec, Vec); 75c4762a1bSJed Brown PetscErrorCode StateMatBlockMult(Mat, Vec, Vec); 76c4762a1bSJed Brown PetscErrorCode StateMatBlockMultTranspose(Mat, Vec, Vec); 77c4762a1bSJed Brown PetscErrorCode StateMatMultTranspose(Mat, Vec, Vec); 78c4762a1bSJed Brown PetscErrorCode StateMatGetDiagonal(Mat, Vec); 79c4762a1bSJed Brown PetscErrorCode StateMatDuplicate(Mat, MatDuplicateOption, Mat *); 80c4762a1bSJed Brown PetscErrorCode StateMatInvMult(Mat, Vec, Vec); 81c4762a1bSJed Brown PetscErrorCode StateMatInvTransposeMult(Mat, Vec, Vec); 82c4762a1bSJed Brown PetscErrorCode StateMatBlockPrecMult(PC, Vec, Vec); 83c4762a1bSJed Brown 84c4762a1bSJed Brown PetscErrorCode DesignMatMult(Mat, Vec, Vec); 85c4762a1bSJed Brown PetscErrorCode DesignMatMultTranspose(Mat, Vec, Vec); 86c4762a1bSJed Brown 87c4762a1bSJed Brown PetscErrorCode Scatter_yi(Vec, Vec *, VecScatter *, PetscInt); /* y to y1,y2,...,y_nt */ 88c4762a1bSJed Brown PetscErrorCode Gather_yi(Vec, Vec *, VecScatter *, PetscInt); 89c4762a1bSJed Brown PetscErrorCode Scatter_uxi_uyi(Vec, Vec *, VecScatter *, Vec *, VecScatter *, PetscInt); /* u to ux_1,uy_1,ux_2,uy_2,...,u */ 90c4762a1bSJed Brown PetscErrorCode Gather_uxi_uyi(Vec, Vec *, VecScatter *, Vec *, VecScatter *, PetscInt); 91c4762a1bSJed Brown 92c4762a1bSJed Brown static char help[] = ""; 93c4762a1bSJed Brown 949371c9d4SSatish Balay int main(int argc, char **argv) { 95c4762a1bSJed Brown Vec x, x0; 96c4762a1bSJed Brown Tao tao; 97c4762a1bSJed Brown AppCtx user; 98c4762a1bSJed Brown IS is_allstate, is_alldesign; 99c4762a1bSJed Brown PetscInt lo, hi, hi2, lo2, ksp_old; 100c4762a1bSJed Brown PetscInt ntests = 1; 101c4762a1bSJed Brown PetscInt i; 102c4762a1bSJed Brown #if defined(PETSC_USE_LOG) 103c4762a1bSJed Brown PetscLogStage stages[1]; 104c4762a1bSJed Brown #endif 105c4762a1bSJed Brown 106327415f7SBarry Smith PetscFunctionBeginUser; 1079566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 108c4762a1bSJed Brown user.mx = 32; 109d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "hyperbolic example", NULL); 1109566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mx", "Number of grid points in each direction", "", user.mx, &user.mx, NULL)); 111c4762a1bSJed Brown user.nt = 16; 1129566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-nt", "Number of time steps", "", user.nt, &user.nt, NULL)); 113c4762a1bSJed Brown user.ndata = 64; 1149566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ndata", "Numbers of data points per sample", "", user.ndata, &user.ndata, NULL)); 115c4762a1bSJed Brown user.alpha = 10.0; 1169566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-alpha", "Regularization parameter", "", user.alpha, &user.alpha, NULL)); 117c4762a1bSJed Brown user.T = 1.0 / 32.0; 1189566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-Tfinal", "Final time", "", user.T, &user.T, NULL)); 1199566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ntests", "Number of times to repeat TaoSolve", "", ntests, &ntests, NULL)); 120d0609cedSBarry Smith PetscOptionsEnd(); 121c4762a1bSJed Brown 122c4762a1bSJed Brown user.m = user.mx * user.mx * user.nt; /* number of constraints */ 123c4762a1bSJed Brown user.n = user.mx * user.mx * 3 * user.nt; /* number of variables */ 124c4762a1bSJed Brown user.ht = user.T / user.nt; /* Time step */ 125c4762a1bSJed Brown user.gamma = user.T * user.ht / (user.mx * user.mx); 126c4762a1bSJed Brown 1279566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user.u)); 1289566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user.y)); 1299566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user.c)); 1309566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user.u, PETSC_DECIDE, user.n - user.m)); 1319566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user.y, PETSC_DECIDE, user.m)); 1329566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user.c, PETSC_DECIDE, user.m)); 1339566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user.u)); 1349566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user.y)); 1359566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user.c)); 136c4762a1bSJed Brown 137c4762a1bSJed Brown /* Create scatters for reduced spaces. 138c4762a1bSJed Brown If the state vector y and design vector u are partitioned as 139c4762a1bSJed Brown [y_1; y_2; ...; y_np] and [u_1; u_2; ...; u_np] (with np = # of processors), 140c4762a1bSJed Brown then the solution vector x is organized as 141c4762a1bSJed Brown [y_1; u_1; y_2; u_2; ...; y_np; u_np]. 142c4762a1bSJed Brown The index sets user.s_is and user.d_is correspond to the indices of the 143c4762a1bSJed Brown state and design variables owned by the current processor. 144c4762a1bSJed Brown */ 1459566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &x)); 146c4762a1bSJed Brown 1479566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user.y, &lo, &hi)); 1489566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user.u, &lo2, &hi2)); 149c4762a1bSJed Brown 1509566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_allstate)); 1519566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo + lo2, 1, &user.s_is)); 152c4762a1bSJed Brown 1539566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi2 - lo2, lo2, 1, &is_alldesign)); 1549566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi2 - lo2, hi + lo2, 1, &user.d_is)); 155c4762a1bSJed Brown 1569566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x, hi - lo + hi2 - lo2, user.n)); 1579566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x)); 158c4762a1bSJed Brown 1599566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(x, user.s_is, user.y, is_allstate, &user.state_scatter)); 1609566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(x, user.d_is, user.u, is_alldesign, &user.design_scatter)); 1619566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_alldesign)); 1629566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_allstate)); 163c4762a1bSJed Brown 164c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 1659566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao)); 1669566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOLCL)); 167c4762a1bSJed Brown 168c4762a1bSJed Brown /* Set up initial vectors and matrices */ 1699566063dSJacob Faibussowitsch PetscCall(HyperbolicInitialize(&user)); 170c4762a1bSJed Brown 1719566063dSJacob Faibussowitsch PetscCall(Gather(x, user.y, user.state_scatter, user.u, user.design_scatter)); 1729566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &x0)); 1739566063dSJacob Faibussowitsch PetscCall(VecCopy(x, x0)); 174c4762a1bSJed Brown 175c4762a1bSJed Brown /* Set solution vector with an initial guess */ 1769566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, x)); 1779566063dSJacob Faibussowitsch PetscCall(TaoSetObjective(tao, FormFunction, &user)); 1789566063dSJacob Faibussowitsch PetscCall(TaoSetGradient(tao, NULL, FormGradient, &user)); 1799566063dSJacob Faibussowitsch PetscCall(TaoSetConstraintsRoutine(tao, user.c, FormConstraints, &user)); 1809566063dSJacob Faibussowitsch PetscCall(TaoSetJacobianStateRoutine(tao, user.Js, user.Js, user.JsInv, FormJacobianState, &user)); 1819566063dSJacob Faibussowitsch PetscCall(TaoSetJacobianDesignRoutine(tao, user.Jd, FormJacobianDesign, &user)); 1829566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 1839566063dSJacob Faibussowitsch PetscCall(TaoSetStateDesignIS(tao, user.s_is, user.d_is)); 184c4762a1bSJed Brown 185c4762a1bSJed Brown /* SOLVE THE APPLICATION */ 1869566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("Trials", &stages[0])); 1879566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stages[0])); 188c4762a1bSJed Brown user.ksp_its_initial = user.ksp_its; 189c4762a1bSJed Brown ksp_old = user.ksp_its; 190c4762a1bSJed Brown for (i = 0; i < ntests; i++) { 1919566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 19263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "KSP Iterations = %" PetscInt_FMT "\n", user.ksp_its - ksp_old)); 1939566063dSJacob Faibussowitsch PetscCall(VecCopy(x0, x)); 1949566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, x)); 195c4762a1bSJed Brown } 1969566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 1979566063dSJacob Faibussowitsch PetscCall(PetscBarrier((PetscObject)x)); 1989566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "KSP iterations within initialization: ")); 19963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "\n", user.ksp_its_initial)); 20063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Total KSP iterations over %" PetscInt_FMT " trial(s): ", ntests)); 20163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "\n", user.ksp_its)); 2029566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "KSP iterations per trial: ")); 20363a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "\n", (user.ksp_its - user.ksp_its_initial) / ntests)); 204c4762a1bSJed Brown 2059566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 2069566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 2079566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x0)); 2089566063dSJacob Faibussowitsch PetscCall(HyperbolicDestroy(&user)); 2099566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 210b122ec5aSJacob Faibussowitsch return 0; 211c4762a1bSJed Brown } 212c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 213c4762a1bSJed Brown /* 214c4762a1bSJed Brown dwork = Qy - d 215c4762a1bSJed Brown lwork = L*(u-ur).^2 216c4762a1bSJed Brown f = 1/2 * (dwork.dork + alpha*y.lwork) 217c4762a1bSJed Brown */ 2189371c9d4SSatish Balay PetscErrorCode FormFunction(Tao tao, Vec X, PetscReal *f, void *ptr) { 219c4762a1bSJed Brown PetscReal d1 = 0, d2 = 0; 220c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 221c4762a1bSJed Brown 222c4762a1bSJed Brown PetscFunctionBegin; 2239566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter)); 2249566063dSJacob Faibussowitsch PetscCall(MatMult(user->Q, user->y, user->dwork)); 2259566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->dwork, -1.0, user->d)); 2269566063dSJacob Faibussowitsch PetscCall(VecDot(user->dwork, user->dwork, &d1)); 227c4762a1bSJed Brown 2289566063dSJacob Faibussowitsch PetscCall(VecWAXPY(user->uwork, -1.0, user->ur, user->u)); 2299566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uwork, user->uwork, user->uwork)); 2309566063dSJacob Faibussowitsch PetscCall(MatMult(user->L, user->uwork, user->lwork)); 2319566063dSJacob Faibussowitsch PetscCall(VecDot(user->y, user->lwork, &d2)); 232c4762a1bSJed Brown *f = 0.5 * (d1 + user->alpha * d2); 233c4762a1bSJed Brown PetscFunctionReturn(0); 234c4762a1bSJed Brown } 235c4762a1bSJed Brown 236c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 237c4762a1bSJed Brown /* 238c4762a1bSJed Brown state: g_s = Q' *(Qy - d) + 0.5*alpha*L*(u-ur).^2 239c4762a1bSJed Brown design: g_d = alpha*(L'y).*(u-ur) 240c4762a1bSJed Brown */ 2419371c9d4SSatish Balay PetscErrorCode FormGradient(Tao tao, Vec X, Vec G, void *ptr) { 242c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 243c4762a1bSJed Brown 244c4762a1bSJed Brown PetscFunctionBegin; 2459566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter)); 2469566063dSJacob Faibussowitsch PetscCall(MatMult(user->Q, user->y, user->dwork)); 2479566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->dwork, -1.0, user->d)); 248c4762a1bSJed Brown 2499566063dSJacob Faibussowitsch PetscCall(MatMult(user->QT, user->dwork, user->ywork)); 250c4762a1bSJed Brown 2519566063dSJacob Faibussowitsch PetscCall(MatMult(user->LT, user->y, user->uwork)); 2529566063dSJacob Faibussowitsch PetscCall(VecWAXPY(user->vwork, -1.0, user->ur, user->u)); 2539566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uwork, user->vwork, user->uwork)); 2549566063dSJacob Faibussowitsch PetscCall(VecScale(user->uwork, user->alpha)); 255c4762a1bSJed Brown 2569566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->vwork, user->vwork, user->vwork)); 2579566063dSJacob Faibussowitsch PetscCall(MatMult(user->L, user->vwork, user->lwork)); 2589566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->ywork, 0.5 * user->alpha, user->lwork)); 259c4762a1bSJed Brown 2609566063dSJacob Faibussowitsch PetscCall(Gather(G, user->ywork, user->state_scatter, user->uwork, user->design_scatter)); 261c4762a1bSJed Brown PetscFunctionReturn(0); 262c4762a1bSJed Brown } 263c4762a1bSJed Brown 2649371c9d4SSatish Balay PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr) { 265c4762a1bSJed Brown PetscReal d1, d2; 266c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 267c4762a1bSJed Brown 268c4762a1bSJed Brown PetscFunctionBegin; 2699566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter)); 2709566063dSJacob Faibussowitsch PetscCall(MatMult(user->Q, user->y, user->dwork)); 2719566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->dwork, -1.0, user->d)); 272c4762a1bSJed Brown 2739566063dSJacob Faibussowitsch PetscCall(MatMult(user->QT, user->dwork, user->ywork)); 274c4762a1bSJed Brown 2759566063dSJacob Faibussowitsch PetscCall(VecDot(user->dwork, user->dwork, &d1)); 276c4762a1bSJed Brown 2779566063dSJacob Faibussowitsch PetscCall(MatMult(user->LT, user->y, user->uwork)); 2789566063dSJacob Faibussowitsch PetscCall(VecWAXPY(user->vwork, -1.0, user->ur, user->u)); 2799566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uwork, user->vwork, user->uwork)); 2809566063dSJacob Faibussowitsch PetscCall(VecScale(user->uwork, user->alpha)); 281c4762a1bSJed Brown 2829566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->vwork, user->vwork, user->vwork)); 2839566063dSJacob Faibussowitsch PetscCall(MatMult(user->L, user->vwork, user->lwork)); 2849566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->ywork, 0.5 * user->alpha, user->lwork)); 285c4762a1bSJed Brown 2869566063dSJacob Faibussowitsch PetscCall(VecDot(user->y, user->lwork, &d2)); 287c4762a1bSJed Brown 288c4762a1bSJed Brown *f = 0.5 * (d1 + user->alpha * d2); 2899566063dSJacob Faibussowitsch PetscCall(Gather(G, user->ywork, user->state_scatter, user->uwork, user->design_scatter)); 290c4762a1bSJed Brown PetscFunctionReturn(0); 291c4762a1bSJed Brown } 292c4762a1bSJed Brown 293c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 294c4762a1bSJed Brown /* A 295c4762a1bSJed Brown MatShell object 296c4762a1bSJed Brown */ 2979371c9d4SSatish Balay PetscErrorCode FormJacobianState(Tao tao, Vec X, Mat J, Mat JPre, Mat JInv, void *ptr) { 298c4762a1bSJed Brown PetscInt i; 299c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 300c4762a1bSJed Brown 301c4762a1bSJed Brown PetscFunctionBegin; 3029566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter)); 3039566063dSJacob Faibussowitsch PetscCall(Scatter_yi(user->u, user->ui, user->ui_scatter, user->nt)); 3049566063dSJacob Faibussowitsch PetscCall(Scatter_uxi_uyi(user->u, user->uxi, user->uxi_scatter, user->uyi, user->uyi_scatter, user->nt)); 305c4762a1bSJed Brown for (i = 0; i < user->nt; i++) { 3069566063dSJacob Faibussowitsch PetscCall(MatCopy(user->Divxy[0], user->C[i], SUBSET_NONZERO_PATTERN)); 3079566063dSJacob Faibussowitsch PetscCall(MatCopy(user->Divxy[1], user->Cwork[i], SAME_NONZERO_PATTERN)); 308c4762a1bSJed Brown 3099566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(user->C[i], NULL, user->uxi[i])); 3109566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(user->Cwork[i], NULL, user->uyi[i])); 3119566063dSJacob Faibussowitsch PetscCall(MatAXPY(user->C[i], 1.0, user->Cwork[i], SUBSET_NONZERO_PATTERN)); 3129566063dSJacob Faibussowitsch PetscCall(MatScale(user->C[i], user->ht)); 3139566063dSJacob Faibussowitsch PetscCall(MatShift(user->C[i], 1.0)); 314c4762a1bSJed Brown } 315c4762a1bSJed Brown PetscFunctionReturn(0); 316c4762a1bSJed Brown } 317c4762a1bSJed Brown 318c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 319c4762a1bSJed Brown /* B */ 3209371c9d4SSatish Balay PetscErrorCode FormJacobianDesign(Tao tao, Vec X, Mat J, void *ptr) { 321c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 322c4762a1bSJed Brown 323c4762a1bSJed Brown PetscFunctionBegin; 3249566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter)); 325c4762a1bSJed Brown PetscFunctionReturn(0); 326c4762a1bSJed Brown } 327c4762a1bSJed Brown 3289371c9d4SSatish Balay PetscErrorCode StateMatMult(Mat J_shell, Vec X, Vec Y) { 329c4762a1bSJed Brown PetscInt i; 330c4762a1bSJed Brown AppCtx *user; 331c4762a1bSJed Brown 332c4762a1bSJed Brown PetscFunctionBegin; 3339566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user)); 3349566063dSJacob Faibussowitsch PetscCall(Scatter_yi(X, user->yi, user->yi_scatter, user->nt)); 335c4762a1bSJed Brown user->block_index = 0; 3369566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock, user->yi[0], user->yiwork[0])); 337c4762a1bSJed Brown 338c4762a1bSJed Brown for (i = 1; i < user->nt; i++) { 339c4762a1bSJed Brown user->block_index = i; 3409566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock, user->yi[i], user->yiwork[i])); 3419566063dSJacob Faibussowitsch PetscCall(MatMult(user->M, user->yi[i - 1], user->ziwork[i - 1])); 3429566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->yiwork[i], -1.0, user->ziwork[i - 1])); 343c4762a1bSJed Brown } 3449566063dSJacob Faibussowitsch PetscCall(Gather_yi(Y, user->yiwork, user->yi_scatter, user->nt)); 345c4762a1bSJed Brown PetscFunctionReturn(0); 346c4762a1bSJed Brown } 347c4762a1bSJed Brown 3489371c9d4SSatish Balay PetscErrorCode StateMatMultTranspose(Mat J_shell, Vec X, Vec Y) { 349c4762a1bSJed Brown PetscInt i; 350c4762a1bSJed Brown AppCtx *user; 351c4762a1bSJed Brown 352c4762a1bSJed Brown PetscFunctionBegin; 3539566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user)); 3549566063dSJacob Faibussowitsch PetscCall(Scatter_yi(X, user->yi, user->yi_scatter, user->nt)); 355c4762a1bSJed Brown 356c4762a1bSJed Brown for (i = 0; i < user->nt - 1; i++) { 357c4762a1bSJed Brown user->block_index = i; 3589566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->JsBlock, user->yi[i], user->yiwork[i])); 3599566063dSJacob Faibussowitsch PetscCall(MatMult(user->M, user->yi[i + 1], user->ziwork[i + 1])); 3609566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->yiwork[i], -1.0, user->ziwork[i + 1])); 361c4762a1bSJed Brown } 362c4762a1bSJed Brown 363c4762a1bSJed Brown i = user->nt - 1; 364c4762a1bSJed Brown user->block_index = i; 3659566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->JsBlock, user->yi[i], user->yiwork[i])); 3669566063dSJacob Faibussowitsch PetscCall(Gather_yi(Y, user->yiwork, user->yi_scatter, user->nt)); 367c4762a1bSJed Brown PetscFunctionReturn(0); 368c4762a1bSJed Brown } 369c4762a1bSJed Brown 3709371c9d4SSatish Balay PetscErrorCode StateMatBlockMult(Mat J_shell, Vec X, Vec Y) { 371c4762a1bSJed Brown PetscInt i; 372c4762a1bSJed Brown AppCtx *user; 373c4762a1bSJed Brown 374c4762a1bSJed Brown PetscFunctionBegin; 3759566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user)); 376c4762a1bSJed Brown i = user->block_index; 3779566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uxiwork[i], X, user->uxi[i])); 3789566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uyiwork[i], X, user->uyi[i])); 3799566063dSJacob Faibussowitsch PetscCall(Gather(user->uiwork[i], user->uxiwork[i], user->ux_scatter[i], user->uyiwork[i], user->uy_scatter[i])); 3809566063dSJacob Faibussowitsch PetscCall(MatMult(user->Div, user->uiwork[i], Y)); 3819566063dSJacob Faibussowitsch PetscCall(VecAYPX(Y, user->ht, X)); 382c4762a1bSJed Brown PetscFunctionReturn(0); 383c4762a1bSJed Brown } 384c4762a1bSJed Brown 3859371c9d4SSatish Balay PetscErrorCode StateMatBlockMultTranspose(Mat J_shell, Vec X, Vec Y) { 386c4762a1bSJed Brown PetscInt i; 387c4762a1bSJed Brown AppCtx *user; 388c4762a1bSJed Brown 389c4762a1bSJed Brown PetscFunctionBegin; 3909566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user)); 391c4762a1bSJed Brown i = user->block_index; 3929566063dSJacob Faibussowitsch PetscCall(MatMult(user->Grad, X, user->uiwork[i])); 3939566063dSJacob Faibussowitsch PetscCall(Scatter(user->uiwork[i], user->uxiwork[i], user->ux_scatter[i], user->uyiwork[i], user->uy_scatter[i])); 3949566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uxiwork[i], user->uxi[i], user->uxiwork[i])); 3959566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uyiwork[i], user->uyi[i], user->uyiwork[i])); 3969566063dSJacob Faibussowitsch PetscCall(VecWAXPY(Y, 1.0, user->uxiwork[i], user->uyiwork[i])); 3979566063dSJacob Faibussowitsch PetscCall(VecAYPX(Y, user->ht, X)); 398c4762a1bSJed Brown PetscFunctionReturn(0); 399c4762a1bSJed Brown } 400c4762a1bSJed Brown 4019371c9d4SSatish Balay PetscErrorCode DesignMatMult(Mat J_shell, Vec X, Vec Y) { 402c4762a1bSJed Brown PetscInt i; 403c4762a1bSJed Brown AppCtx *user; 404c4762a1bSJed Brown 405c4762a1bSJed Brown PetscFunctionBegin; 4069566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user)); 4079566063dSJacob Faibussowitsch PetscCall(Scatter_yi(user->y, user->yi, user->yi_scatter, user->nt)); 4089566063dSJacob Faibussowitsch PetscCall(Scatter_uxi_uyi(X, user->uxiwork, user->uxi_scatter, user->uyiwork, user->uyi_scatter, user->nt)); 409c4762a1bSJed Brown for (i = 0; i < user->nt; i++) { 4109566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uxiwork[i], user->yi[i], user->uxiwork[i])); 4119566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uyiwork[i], user->yi[i], user->uyiwork[i])); 4129566063dSJacob Faibussowitsch PetscCall(Gather(user->uiwork[i], user->uxiwork[i], user->ux_scatter[i], user->uyiwork[i], user->uy_scatter[i])); 4139566063dSJacob Faibussowitsch PetscCall(MatMult(user->Div, user->uiwork[i], user->ziwork[i])); 4149566063dSJacob Faibussowitsch PetscCall(VecScale(user->ziwork[i], user->ht)); 415c4762a1bSJed Brown } 4169566063dSJacob Faibussowitsch PetscCall(Gather_yi(Y, user->ziwork, user->yi_scatter, user->nt)); 417c4762a1bSJed Brown PetscFunctionReturn(0); 418c4762a1bSJed Brown } 419c4762a1bSJed Brown 4209371c9d4SSatish Balay PetscErrorCode DesignMatMultTranspose(Mat J_shell, Vec X, Vec Y) { 421c4762a1bSJed Brown PetscInt i; 422c4762a1bSJed Brown AppCtx *user; 423c4762a1bSJed Brown 424c4762a1bSJed Brown PetscFunctionBegin; 4259566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user)); 4269566063dSJacob Faibussowitsch PetscCall(Scatter_yi(user->y, user->yi, user->yi_scatter, user->nt)); 4279566063dSJacob Faibussowitsch PetscCall(Scatter_yi(X, user->yiwork, user->yi_scatter, user->nt)); 428c4762a1bSJed Brown for (i = 0; i < user->nt; i++) { 4299566063dSJacob Faibussowitsch PetscCall(MatMult(user->Grad, user->yiwork[i], user->uiwork[i])); 4309566063dSJacob Faibussowitsch PetscCall(Scatter(user->uiwork[i], user->uxiwork[i], user->ux_scatter[i], user->uyiwork[i], user->uy_scatter[i])); 4319566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uxiwork[i], user->yi[i], user->uxiwork[i])); 4329566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uyiwork[i], user->yi[i], user->uyiwork[i])); 4339566063dSJacob Faibussowitsch PetscCall(Gather(user->uiwork[i], user->uxiwork[i], user->ux_scatter[i], user->uyiwork[i], user->uy_scatter[i])); 4349566063dSJacob Faibussowitsch PetscCall(VecScale(user->uiwork[i], user->ht)); 435c4762a1bSJed Brown } 4369566063dSJacob Faibussowitsch PetscCall(Gather_yi(Y, user->uiwork, user->ui_scatter, user->nt)); 437c4762a1bSJed Brown PetscFunctionReturn(0); 438c4762a1bSJed Brown } 439c4762a1bSJed Brown 4409371c9d4SSatish Balay PetscErrorCode StateMatBlockPrecMult(PC PC_shell, Vec X, Vec Y) { 441c4762a1bSJed Brown PetscInt i; 442c4762a1bSJed Brown AppCtx *user; 443c4762a1bSJed Brown 444c4762a1bSJed Brown PetscFunctionBegin; 4459566063dSJacob Faibussowitsch PetscCall(PCShellGetContext(PC_shell, &user)); 446c4762a1bSJed Brown i = user->block_index; 447c4762a1bSJed Brown if (user->c_formed) { 4489566063dSJacob Faibussowitsch PetscCall(MatSOR(user->C[i], X, 1.0, (MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_SYMMETRIC_SWEEP), 0.0, 1, 1, Y)); 449c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not formed"); 450c4762a1bSJed Brown PetscFunctionReturn(0); 451c4762a1bSJed Brown } 452c4762a1bSJed Brown 4539371c9d4SSatish Balay PetscErrorCode StateMatBlockPrecMultTranspose(PC PC_shell, Vec X, Vec Y) { 454c4762a1bSJed Brown PetscInt i; 455c4762a1bSJed Brown AppCtx *user; 456c4762a1bSJed Brown 457c4762a1bSJed Brown PetscFunctionBegin; 4589566063dSJacob Faibussowitsch PetscCall(PCShellGetContext(PC_shell, &user)); 459c4762a1bSJed Brown 460c4762a1bSJed Brown i = user->block_index; 461c4762a1bSJed Brown if (user->c_formed) { 4629566063dSJacob Faibussowitsch PetscCall(MatSOR(user->C[i], X, 1.0, (MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_SYMMETRIC_SWEEP), 0.0, 1, 1, Y)); 463c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not formed"); 464c4762a1bSJed Brown PetscFunctionReturn(0); 465c4762a1bSJed Brown } 466c4762a1bSJed Brown 4679371c9d4SSatish Balay PetscErrorCode StateMatInvMult(Mat J_shell, Vec X, Vec Y) { 468c4762a1bSJed Brown AppCtx *user; 469c4762a1bSJed Brown PetscInt its, i; 470c4762a1bSJed Brown 471c4762a1bSJed Brown PetscFunctionBegin; 4729566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user)); 473c4762a1bSJed Brown 474c4762a1bSJed Brown if (Y == user->ytrue) { 475c4762a1bSJed Brown /* First solve is done using true solution to set up problem */ 4769566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(user->solver, 1e-4, 1e-20, PETSC_DEFAULT, PETSC_DEFAULT)); 477c4762a1bSJed Brown } else { 4789566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(user->solver, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 479c4762a1bSJed Brown } 4809566063dSJacob Faibussowitsch PetscCall(Scatter_yi(X, user->yi, user->yi_scatter, user->nt)); 4819566063dSJacob Faibussowitsch PetscCall(Scatter_yi(Y, user->yiwork, user->yi_scatter, user->nt)); 4829566063dSJacob Faibussowitsch PetscCall(Scatter_uxi_uyi(user->u, user->uxi, user->uxi_scatter, user->uyi, user->uyi_scatter, user->nt)); 483c4762a1bSJed Brown 484c4762a1bSJed Brown user->block_index = 0; 4859566063dSJacob Faibussowitsch PetscCall(KSPSolve(user->solver, user->yi[0], user->yiwork[0])); 486c4762a1bSJed Brown 4879566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(user->solver, &its)); 488c4762a1bSJed Brown user->ksp_its = user->ksp_its + its; 489c4762a1bSJed Brown for (i = 1; i < user->nt; i++) { 4909566063dSJacob Faibussowitsch PetscCall(MatMult(user->M, user->yiwork[i - 1], user->ziwork[i - 1])); 4919566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->yi[i], 1.0, user->ziwork[i - 1])); 492c4762a1bSJed Brown user->block_index = i; 4939566063dSJacob Faibussowitsch PetscCall(KSPSolve(user->solver, user->yi[i], user->yiwork[i])); 494c4762a1bSJed Brown 4959566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(user->solver, &its)); 496c4762a1bSJed Brown user->ksp_its = user->ksp_its + its; 497c4762a1bSJed Brown } 4989566063dSJacob Faibussowitsch PetscCall(Gather_yi(Y, user->yiwork, user->yi_scatter, user->nt)); 499c4762a1bSJed Brown PetscFunctionReturn(0); 500c4762a1bSJed Brown } 501c4762a1bSJed Brown 5029371c9d4SSatish Balay PetscErrorCode StateMatInvTransposeMult(Mat J_shell, Vec X, Vec Y) { 503c4762a1bSJed Brown AppCtx *user; 504c4762a1bSJed Brown PetscInt its, i; 505c4762a1bSJed Brown 506c4762a1bSJed Brown PetscFunctionBegin; 5079566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user)); 508c4762a1bSJed Brown 5099566063dSJacob Faibussowitsch PetscCall(Scatter_yi(X, user->yi, user->yi_scatter, user->nt)); 5109566063dSJacob Faibussowitsch PetscCall(Scatter_yi(Y, user->yiwork, user->yi_scatter, user->nt)); 5119566063dSJacob Faibussowitsch PetscCall(Scatter_uxi_uyi(user->u, user->uxi, user->uxi_scatter, user->uyi, user->uyi_scatter, user->nt)); 512c4762a1bSJed Brown 513c4762a1bSJed Brown i = user->nt - 1; 514c4762a1bSJed Brown user->block_index = i; 5159566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(user->solver, user->yi[i], user->yiwork[i])); 516c4762a1bSJed Brown 5179566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(user->solver, &its)); 518c4762a1bSJed Brown user->ksp_its = user->ksp_its + its; 519c4762a1bSJed Brown 520c4762a1bSJed Brown for (i = user->nt - 2; i >= 0; i--) { 5219566063dSJacob Faibussowitsch PetscCall(MatMult(user->M, user->yiwork[i + 1], user->ziwork[i + 1])); 5229566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->yi[i], 1.0, user->ziwork[i + 1])); 523c4762a1bSJed Brown user->block_index = i; 5249566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(user->solver, user->yi[i], user->yiwork[i])); 525c4762a1bSJed Brown 5269566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(user->solver, &its)); 527c4762a1bSJed Brown user->ksp_its = user->ksp_its + its; 528c4762a1bSJed Brown } 5299566063dSJacob Faibussowitsch PetscCall(Gather_yi(Y, user->yiwork, user->yi_scatter, user->nt)); 530c4762a1bSJed Brown PetscFunctionReturn(0); 531c4762a1bSJed Brown } 532c4762a1bSJed Brown 5339371c9d4SSatish Balay PetscErrorCode StateMatDuplicate(Mat J_shell, MatDuplicateOption opt, Mat *new_shell) { 534c4762a1bSJed Brown AppCtx *user; 535c4762a1bSJed Brown 536c4762a1bSJed Brown PetscFunctionBegin; 5379566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user)); 538c4762a1bSJed Brown 5399566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m, user->m, user, new_shell)); 5409566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*new_shell, MATOP_MULT, (void (*)(void))StateMatMult)); 5419566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*new_shell, MATOP_DUPLICATE, (void (*)(void))StateMatDuplicate)); 5429566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*new_shell, MATOP_MULT_TRANSPOSE, (void (*)(void))StateMatMultTranspose)); 5439566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*new_shell, MATOP_GET_DIAGONAL, (void (*)(void))StateMatGetDiagonal)); 544c4762a1bSJed Brown PetscFunctionReturn(0); 545c4762a1bSJed Brown } 546c4762a1bSJed Brown 5479371c9d4SSatish Balay PetscErrorCode StateMatGetDiagonal(Mat J_shell, Vec X) { 548c4762a1bSJed Brown AppCtx *user; 549c4762a1bSJed Brown 550c4762a1bSJed Brown PetscFunctionBegin; 5519566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user)); 5529566063dSJacob Faibussowitsch PetscCall(VecCopy(user->js_diag, X)); 553c4762a1bSJed Brown PetscFunctionReturn(0); 554c4762a1bSJed Brown } 555c4762a1bSJed Brown 5569371c9d4SSatish Balay PetscErrorCode FormConstraints(Tao tao, Vec X, Vec C, void *ptr) { 557c4762a1bSJed Brown /* con = Ay - q, A = [C(u1) 0 0 ... 0; 558c4762a1bSJed Brown -M C(u2) 0 ... 0; 559c4762a1bSJed Brown 0 -M C(u3) ... 0; 560c4762a1bSJed Brown ... ; 561c4762a1bSJed Brown 0 ... -M C(u_nt)] 562c4762a1bSJed Brown C(u) = eye + ht*Div*[diag(u1); diag(u2)] */ 563c4762a1bSJed Brown PetscInt i; 564c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 565c4762a1bSJed Brown 566c4762a1bSJed Brown PetscFunctionBegin; 5679566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter)); 5689566063dSJacob Faibussowitsch PetscCall(Scatter_yi(user->y, user->yi, user->yi_scatter, user->nt)); 5699566063dSJacob Faibussowitsch PetscCall(Scatter_uxi_uyi(user->u, user->uxi, user->uxi_scatter, user->uyi, user->uyi_scatter, user->nt)); 570c4762a1bSJed Brown 571c4762a1bSJed Brown user->block_index = 0; 5729566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock, user->yi[0], user->yiwork[0])); 573c4762a1bSJed Brown 574c4762a1bSJed Brown for (i = 1; i < user->nt; i++) { 575c4762a1bSJed Brown user->block_index = i; 5769566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock, user->yi[i], user->yiwork[i])); 5779566063dSJacob Faibussowitsch PetscCall(MatMult(user->M, user->yi[i - 1], user->ziwork[i - 1])); 5789566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->yiwork[i], -1.0, user->ziwork[i - 1])); 579c4762a1bSJed Brown } 580c4762a1bSJed Brown 5819566063dSJacob Faibussowitsch PetscCall(Gather_yi(C, user->yiwork, user->yi_scatter, user->nt)); 5829566063dSJacob Faibussowitsch PetscCall(VecAXPY(C, -1.0, user->q)); 583c4762a1bSJed Brown 584c4762a1bSJed Brown PetscFunctionReturn(0); 585c4762a1bSJed Brown } 586c4762a1bSJed Brown 5879371c9d4SSatish Balay PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat) { 588c4762a1bSJed Brown PetscFunctionBegin; 5899566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(s_scat, x, state, INSERT_VALUES, SCATTER_FORWARD)); 5909566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(s_scat, x, state, INSERT_VALUES, SCATTER_FORWARD)); 5919566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(d_scat, x, design, INSERT_VALUES, SCATTER_FORWARD)); 5929566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(d_scat, x, design, INSERT_VALUES, SCATTER_FORWARD)); 593c4762a1bSJed Brown PetscFunctionReturn(0); 594c4762a1bSJed Brown } 595c4762a1bSJed Brown 5969371c9d4SSatish Balay PetscErrorCode Scatter_uxi_uyi(Vec u, Vec *uxi, VecScatter *scatx, Vec *uyi, VecScatter *scaty, PetscInt nt) { 597c4762a1bSJed Brown PetscInt i; 598c4762a1bSJed Brown 599c4762a1bSJed Brown PetscFunctionBegin; 600c4762a1bSJed Brown for (i = 0; i < nt; i++) { 6019566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scatx[i], u, uxi[i], INSERT_VALUES, SCATTER_FORWARD)); 6029566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scatx[i], u, uxi[i], INSERT_VALUES, SCATTER_FORWARD)); 6039566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scaty[i], u, uyi[i], INSERT_VALUES, SCATTER_FORWARD)); 6049566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scaty[i], u, uyi[i], INSERT_VALUES, SCATTER_FORWARD)); 605c4762a1bSJed Brown } 606c4762a1bSJed Brown PetscFunctionReturn(0); 607c4762a1bSJed Brown } 608c4762a1bSJed Brown 6099371c9d4SSatish Balay PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat) { 610c4762a1bSJed Brown PetscFunctionBegin; 6119566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(s_scat, state, x, INSERT_VALUES, SCATTER_REVERSE)); 6129566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(s_scat, state, x, INSERT_VALUES, SCATTER_REVERSE)); 6139566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(d_scat, design, x, INSERT_VALUES, SCATTER_REVERSE)); 6149566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(d_scat, design, x, INSERT_VALUES, SCATTER_REVERSE)); 615c4762a1bSJed Brown PetscFunctionReturn(0); 616c4762a1bSJed Brown } 617c4762a1bSJed Brown 6189371c9d4SSatish Balay PetscErrorCode Gather_uxi_uyi(Vec u, Vec *uxi, VecScatter *scatx, Vec *uyi, VecScatter *scaty, PetscInt nt) { 619c4762a1bSJed Brown PetscInt i; 620c4762a1bSJed Brown 621c4762a1bSJed Brown PetscFunctionBegin; 622c4762a1bSJed Brown for (i = 0; i < nt; i++) { 6239566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scatx[i], uxi[i], u, INSERT_VALUES, SCATTER_REVERSE)); 6249566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scatx[i], uxi[i], u, INSERT_VALUES, SCATTER_REVERSE)); 6259566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scaty[i], uyi[i], u, INSERT_VALUES, SCATTER_REVERSE)); 6269566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scaty[i], uyi[i], u, INSERT_VALUES, SCATTER_REVERSE)); 627c4762a1bSJed Brown } 628c4762a1bSJed Brown PetscFunctionReturn(0); 629c4762a1bSJed Brown } 630c4762a1bSJed Brown 6319371c9d4SSatish Balay PetscErrorCode Scatter_yi(Vec y, Vec *yi, VecScatter *scat, PetscInt nt) { 632c4762a1bSJed Brown PetscInt i; 633c4762a1bSJed Brown 634c4762a1bSJed Brown PetscFunctionBegin; 635c4762a1bSJed Brown for (i = 0; i < nt; i++) { 6369566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scat[i], y, yi[i], INSERT_VALUES, SCATTER_FORWARD)); 6379566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scat[i], y, yi[i], INSERT_VALUES, SCATTER_FORWARD)); 638c4762a1bSJed Brown } 639c4762a1bSJed Brown PetscFunctionReturn(0); 640c4762a1bSJed Brown } 641c4762a1bSJed Brown 6429371c9d4SSatish Balay PetscErrorCode Gather_yi(Vec y, Vec *yi, VecScatter *scat, PetscInt nt) { 643c4762a1bSJed Brown PetscInt i; 644c4762a1bSJed Brown 645c4762a1bSJed Brown PetscFunctionBegin; 646c4762a1bSJed Brown for (i = 0; i < nt; i++) { 6479566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scat[i], yi[i], y, INSERT_VALUES, SCATTER_REVERSE)); 6489566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scat[i], yi[i], y, INSERT_VALUES, SCATTER_REVERSE)); 649c4762a1bSJed Brown } 650c4762a1bSJed Brown PetscFunctionReturn(0); 651c4762a1bSJed Brown } 652c4762a1bSJed Brown 6539371c9d4SSatish Balay PetscErrorCode HyperbolicInitialize(AppCtx *user) { 654c4762a1bSJed Brown PetscInt n, i, j, linear_index, istart, iend, iblock, lo, hi; 655c4762a1bSJed Brown Vec XX, YY, XXwork, YYwork, yi, uxi, ui, bc; 656c4762a1bSJed Brown PetscReal h, sum; 657c4762a1bSJed Brown PetscScalar hinv, neg_hinv, quarter = 0.25, one = 1.0, half_hinv, neg_half_hinv; 658c4762a1bSJed Brown PetscScalar vx, vy, zero = 0.0; 659c4762a1bSJed Brown IS is_from_y, is_to_yi, is_from_u, is_to_uxi, is_to_uyi; 660c4762a1bSJed Brown 661c4762a1bSJed Brown PetscFunctionBegin; 662c4762a1bSJed Brown user->jformed = PETSC_FALSE; 663c4762a1bSJed Brown user->c_formed = PETSC_FALSE; 664c4762a1bSJed Brown 665c4762a1bSJed Brown user->ksp_its = 0; 666c4762a1bSJed Brown user->ksp_its_initial = 0; 667c4762a1bSJed Brown 668c4762a1bSJed Brown n = user->mx * user->mx; 669c4762a1bSJed Brown 670c4762a1bSJed Brown h = 1.0 / user->mx; 671c4762a1bSJed Brown hinv = user->mx; 672c4762a1bSJed Brown neg_hinv = -hinv; 673c4762a1bSJed Brown half_hinv = hinv / 2.0; 674c4762a1bSJed Brown neg_half_hinv = neg_hinv / 2.0; 675c4762a1bSJed Brown 676c4762a1bSJed Brown /* Generate Grad matrix */ 6779566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Grad)); 6789566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Grad, PETSC_DECIDE, PETSC_DECIDE, 2 * n, n)); 6799566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Grad)); 6809566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->Grad, 3, NULL, 3, NULL)); 6819566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->Grad, 3, NULL)); 6829566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->Grad, &istart, &iend)); 683c4762a1bSJed Brown 684c4762a1bSJed Brown for (i = istart; i < iend; i++) { 685c4762a1bSJed Brown if (i < n) { 686c4762a1bSJed Brown iblock = i / user->mx; 687c4762a1bSJed Brown j = iblock * user->mx + ((i + user->mx - 1) % user->mx); 6889566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &half_hinv, INSERT_VALUES)); 689c4762a1bSJed Brown j = iblock * user->mx + ((i + 1) % user->mx); 6909566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &neg_half_hinv, INSERT_VALUES)); 691c4762a1bSJed Brown } 692c4762a1bSJed Brown if (i >= n) { 693c4762a1bSJed Brown j = (i - user->mx) % n; 6949566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &half_hinv, INSERT_VALUES)); 695c4762a1bSJed Brown j = (j + 2 * user->mx) % n; 6969566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &neg_half_hinv, INSERT_VALUES)); 697c4762a1bSJed Brown } 698c4762a1bSJed Brown } 699c4762a1bSJed Brown 7009566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Grad, MAT_FINAL_ASSEMBLY)); 7019566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Grad, MAT_FINAL_ASSEMBLY)); 702c4762a1bSJed Brown 7039566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Gradxy[0])); 7049566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Gradxy[0], PETSC_DECIDE, PETSC_DECIDE, n, n)); 7059566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Gradxy[0])); 7069566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->Gradxy[0], 3, NULL, 3, NULL)); 7079566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->Gradxy[0], 3, NULL)); 7089566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->Gradxy[0], &istart, &iend)); 709c4762a1bSJed Brown 710c4762a1bSJed Brown for (i = istart; i < iend; i++) { 711c4762a1bSJed Brown iblock = i / user->mx; 712c4762a1bSJed Brown j = iblock * user->mx + ((i + user->mx - 1) % user->mx); 7139566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Gradxy[0], 1, &i, 1, &j, &half_hinv, INSERT_VALUES)); 714c4762a1bSJed Brown j = iblock * user->mx + ((i + 1) % user->mx); 7159566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Gradxy[0], 1, &i, 1, &j, &neg_half_hinv, INSERT_VALUES)); 7169566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Gradxy[0], 1, &i, 1, &i, &zero, INSERT_VALUES)); 717c4762a1bSJed Brown } 7189566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Gradxy[0], MAT_FINAL_ASSEMBLY)); 7199566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Gradxy[0], MAT_FINAL_ASSEMBLY)); 720c4762a1bSJed Brown 7219566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Gradxy[1])); 7229566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Gradxy[1], PETSC_DECIDE, PETSC_DECIDE, n, n)); 7239566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Gradxy[1])); 7249566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->Gradxy[1], 3, NULL, 3, NULL)); 7259566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->Gradxy[1], 3, NULL)); 7269566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->Gradxy[1], &istart, &iend)); 727c4762a1bSJed Brown 728c4762a1bSJed Brown for (i = istart; i < iend; i++) { 729c4762a1bSJed Brown j = (i + n - user->mx) % n; 7309566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Gradxy[1], 1, &i, 1, &j, &half_hinv, INSERT_VALUES)); 731c4762a1bSJed Brown j = (j + 2 * user->mx) % n; 7329566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Gradxy[1], 1, &i, 1, &j, &neg_half_hinv, INSERT_VALUES)); 7339566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Gradxy[1], 1, &i, 1, &i, &zero, INSERT_VALUES)); 734c4762a1bSJed Brown } 7359566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Gradxy[1], MAT_FINAL_ASSEMBLY)); 7369566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Gradxy[1], MAT_FINAL_ASSEMBLY)); 737c4762a1bSJed Brown 738c4762a1bSJed Brown /* Generate Div matrix */ 7399566063dSJacob Faibussowitsch PetscCall(MatTranspose(user->Grad, MAT_INITIAL_MATRIX, &user->Div)); 7409566063dSJacob Faibussowitsch PetscCall(MatTranspose(user->Gradxy[0], MAT_INITIAL_MATRIX, &user->Divxy[0])); 7419566063dSJacob Faibussowitsch PetscCall(MatTranspose(user->Gradxy[1], MAT_INITIAL_MATRIX, &user->Divxy[1])); 742c4762a1bSJed Brown 743c4762a1bSJed Brown /* Off-diagonal averaging matrix */ 7449566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->M)); 7459566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->M, PETSC_DECIDE, PETSC_DECIDE, n, n)); 7469566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->M)); 7479566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->M, 4, NULL, 4, NULL)); 7489566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->M, 4, NULL)); 7499566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->M, &istart, &iend)); 750c4762a1bSJed Brown 751c4762a1bSJed Brown for (i = istart; i < iend; i++) { 752c4762a1bSJed Brown /* kron(Id,Av) */ 753c4762a1bSJed Brown iblock = i / user->mx; 754c4762a1bSJed Brown j = iblock * user->mx + ((i + user->mx - 1) % user->mx); 7559566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->M, 1, &i, 1, &j, &quarter, INSERT_VALUES)); 756c4762a1bSJed Brown j = iblock * user->mx + ((i + 1) % user->mx); 7579566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->M, 1, &i, 1, &j, &quarter, INSERT_VALUES)); 758c4762a1bSJed Brown 759c4762a1bSJed Brown /* kron(Av,Id) */ 760c4762a1bSJed Brown j = (i + user->mx) % n; 7619566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->M, 1, &i, 1, &j, &quarter, INSERT_VALUES)); 762c4762a1bSJed Brown j = (i + n - user->mx) % n; 7639566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->M, 1, &i, 1, &j, &quarter, INSERT_VALUES)); 764c4762a1bSJed Brown } 7659566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->M, MAT_FINAL_ASSEMBLY)); 7669566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->M, MAT_FINAL_ASSEMBLY)); 767c4762a1bSJed Brown 768c4762a1bSJed Brown /* Generate 2D grid */ 7699566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &XX)); 7709566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->q)); 7719566063dSJacob Faibussowitsch PetscCall(VecSetSizes(XX, PETSC_DECIDE, n)); 7729566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->q, PETSC_DECIDE, n * user->nt)); 7739566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(XX)); 7749566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->q)); 775c4762a1bSJed Brown 7769566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &YY)); 7779566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &XXwork)); 7789566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &YYwork)); 7799566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &user->d)); 7809566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &user->dwork)); 781c4762a1bSJed Brown 7829566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(XX, &istart, &iend)); 783c4762a1bSJed Brown for (linear_index = istart; linear_index < iend; linear_index++) { 784c4762a1bSJed Brown i = linear_index % user->mx; 785c4762a1bSJed Brown j = (linear_index - i) / user->mx; 786c4762a1bSJed Brown vx = h * (i + 0.5); 787c4762a1bSJed Brown vy = h * (j + 0.5); 7889566063dSJacob Faibussowitsch PetscCall(VecSetValues(XX, 1, &linear_index, &vx, INSERT_VALUES)); 7899566063dSJacob Faibussowitsch PetscCall(VecSetValues(YY, 1, &linear_index, &vy, INSERT_VALUES)); 790c4762a1bSJed Brown } 791c4762a1bSJed Brown 7929566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(XX)); 7939566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(XX)); 7949566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(YY)); 7959566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(YY)); 796c4762a1bSJed Brown 797c4762a1bSJed Brown /* Compute final density function yT 798c4762a1bSJed Brown yT = 1.0 + exp(-30*((x-0.25)^2+(y-0.25)^2)) + exp(-30*((x-0.75)^2+(y-0.75)^2)) 799c4762a1bSJed Brown yT = yT / (h^2*sum(yT)) */ 8009566063dSJacob Faibussowitsch PetscCall(VecCopy(XX, XXwork)); 8019566063dSJacob Faibussowitsch PetscCall(VecCopy(YY, YYwork)); 802c4762a1bSJed Brown 8039566063dSJacob Faibussowitsch PetscCall(VecShift(XXwork, -0.25)); 8049566063dSJacob Faibussowitsch PetscCall(VecShift(YYwork, -0.25)); 805c4762a1bSJed Brown 8069566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(XXwork, XXwork, XXwork)); 8079566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(YYwork, YYwork, YYwork)); 808c4762a1bSJed Brown 8099566063dSJacob Faibussowitsch PetscCall(VecCopy(XXwork, user->dwork)); 8109566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->dwork, 1.0, YYwork)); 8119566063dSJacob Faibussowitsch PetscCall(VecScale(user->dwork, -30.0)); 8129566063dSJacob Faibussowitsch PetscCall(VecExp(user->dwork)); 8139566063dSJacob Faibussowitsch PetscCall(VecCopy(user->dwork, user->d)); 814c4762a1bSJed Brown 8159566063dSJacob Faibussowitsch PetscCall(VecCopy(XX, XXwork)); 8169566063dSJacob Faibussowitsch PetscCall(VecCopy(YY, YYwork)); 817c4762a1bSJed Brown 8189566063dSJacob Faibussowitsch PetscCall(VecShift(XXwork, -0.75)); 8199566063dSJacob Faibussowitsch PetscCall(VecShift(YYwork, -0.75)); 820c4762a1bSJed Brown 8219566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(XXwork, XXwork, XXwork)); 8229566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(YYwork, YYwork, YYwork)); 823c4762a1bSJed Brown 8249566063dSJacob Faibussowitsch PetscCall(VecCopy(XXwork, user->dwork)); 8259566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->dwork, 1.0, YYwork)); 8269566063dSJacob Faibussowitsch PetscCall(VecScale(user->dwork, -30.0)); 8279566063dSJacob Faibussowitsch PetscCall(VecExp(user->dwork)); 828c4762a1bSJed Brown 8299566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->d, 1.0, user->dwork)); 8309566063dSJacob Faibussowitsch PetscCall(VecShift(user->d, 1.0)); 8319566063dSJacob Faibussowitsch PetscCall(VecSum(user->d, &sum)); 8329566063dSJacob Faibussowitsch PetscCall(VecScale(user->d, 1.0 / (h * h * sum))); 833c4762a1bSJed Brown 834c4762a1bSJed Brown /* Initial conditions of forward problem */ 8359566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &bc)); 8369566063dSJacob Faibussowitsch PetscCall(VecCopy(XX, XXwork)); 8379566063dSJacob Faibussowitsch PetscCall(VecCopy(YY, YYwork)); 838c4762a1bSJed Brown 8399566063dSJacob Faibussowitsch PetscCall(VecShift(XXwork, -0.5)); 8409566063dSJacob Faibussowitsch PetscCall(VecShift(YYwork, -0.5)); 841c4762a1bSJed Brown 8429566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(XXwork, XXwork, XXwork)); 8439566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(YYwork, YYwork, YYwork)); 844c4762a1bSJed Brown 8459566063dSJacob Faibussowitsch PetscCall(VecWAXPY(bc, 1.0, XXwork, YYwork)); 8469566063dSJacob Faibussowitsch PetscCall(VecScale(bc, -50.0)); 8479566063dSJacob Faibussowitsch PetscCall(VecExp(bc)); 8489566063dSJacob Faibussowitsch PetscCall(VecShift(bc, 1.0)); 8499566063dSJacob Faibussowitsch PetscCall(VecSum(bc, &sum)); 8509566063dSJacob Faibussowitsch PetscCall(VecScale(bc, 1.0 / (h * h * sum))); 851c4762a1bSJed Brown 852c4762a1bSJed Brown /* Create scatter from y to y_1,y_2,...,y_nt */ 853c4762a1bSJed Brown /* TODO: Reorder for better parallelism. (This will require reordering Q and L as well.) */ 8549566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->nt * user->mx * user->mx, &user->yi_scatter)); 8559566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &yi)); 8569566063dSJacob Faibussowitsch PetscCall(VecSetSizes(yi, PETSC_DECIDE, user->mx * user->mx)); 8579566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(yi)); 8589566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(yi, user->nt, &user->yi)); 8599566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(yi, user->nt, &user->yiwork)); 8609566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(yi, user->nt, &user->ziwork)); 861c4762a1bSJed Brown for (i = 0; i < user->nt; i++) { 8629566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->yi[i], &lo, &hi)); 8639566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_to_yi)); 8649566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo + i * user->mx * user->mx, 1, &is_from_y)); 8659566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(user->y, is_from_y, user->yi[i], is_to_yi, &user->yi_scatter[i])); 8669566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_to_yi)); 8679566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_from_y)); 868c4762a1bSJed Brown } 869c4762a1bSJed Brown 870c4762a1bSJed Brown /* Create scatter from u to ux_1,uy_1,ux_2,uy_2,...,ux_nt,uy_nt */ 871c4762a1bSJed Brown /* TODO: reorder for better parallelism */ 8729566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->nt * user->mx * user->mx, &user->uxi_scatter)); 8739566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->nt * user->mx * user->mx, &user->uyi_scatter)); 8749566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->nt * user->mx * user->mx, &user->ux_scatter)); 8759566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->nt * user->mx * user->mx, &user->uy_scatter)); 8769566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2 * user->nt * user->mx * user->mx, &user->ui_scatter)); 8779566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &uxi)); 8789566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &ui)); 8799566063dSJacob Faibussowitsch PetscCall(VecSetSizes(uxi, PETSC_DECIDE, user->mx * user->mx)); 8809566063dSJacob Faibussowitsch PetscCall(VecSetSizes(ui, PETSC_DECIDE, 2 * user->mx * user->mx)); 8819566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(uxi)); 8829566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(ui)); 8839566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(uxi, user->nt, &user->uxi)); 8849566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(uxi, user->nt, &user->uyi)); 8859566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(uxi, user->nt, &user->uxiwork)); 8869566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(uxi, user->nt, &user->uyiwork)); 8879566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ui, user->nt, &user->ui)); 8889566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ui, user->nt, &user->uiwork)); 889c4762a1bSJed Brown for (i = 0; i < user->nt; i++) { 8909566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->uxi[i], &lo, &hi)); 8919566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_to_uxi)); 8929566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo + 2 * i * user->mx * user->mx, 1, &is_from_u)); 8939566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(user->u, is_from_u, user->uxi[i], is_to_uxi, &user->uxi_scatter[i])); 894c4762a1bSJed Brown 8959566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_to_uxi)); 8969566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_from_u)); 897c4762a1bSJed Brown 8989566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->uyi[i], &lo, &hi)); 8999566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_to_uyi)); 9009566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo + (2 * i + 1) * user->mx * user->mx, 1, &is_from_u)); 9019566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(user->u, is_from_u, user->uyi[i], is_to_uyi, &user->uyi_scatter[i])); 902c4762a1bSJed Brown 9039566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_to_uyi)); 9049566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_from_u)); 905c4762a1bSJed Brown 9069566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->uxi[i], &lo, &hi)); 9079566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_to_uxi)); 9089566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_from_u)); 9099566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(user->ui[i], is_from_u, user->uxi[i], is_to_uxi, &user->ux_scatter[i])); 910c4762a1bSJed Brown 9119566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_to_uxi)); 9129566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_from_u)); 913c4762a1bSJed Brown 9149566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->uyi[i], &lo, &hi)); 9159566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_to_uyi)); 9169566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo + user->mx * user->mx, 1, &is_from_u)); 9179566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(user->ui[i], is_from_u, user->uyi[i], is_to_uyi, &user->uy_scatter[i])); 918c4762a1bSJed Brown 9199566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_to_uyi)); 9209566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_from_u)); 921c4762a1bSJed Brown 9229566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->ui[i], &lo, &hi)); 9239566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_to_uxi)); 9249566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo + 2 * i * user->mx * user->mx, 1, &is_from_u)); 9259566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(user->u, is_from_u, user->ui[i], is_to_uxi, &user->ui_scatter[i])); 926c4762a1bSJed Brown 9279566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_to_uxi)); 9289566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_from_u)); 929c4762a1bSJed Brown } 930c4762a1bSJed Brown 931c4762a1bSJed Brown /* RHS of forward problem */ 9329566063dSJacob Faibussowitsch PetscCall(MatMult(user->M, bc, user->yiwork[0])); 933*48a46eb9SPierre Jolivet for (i = 1; i < user->nt; i++) PetscCall(VecSet(user->yiwork[i], 0.0)); 9349566063dSJacob Faibussowitsch PetscCall(Gather_yi(user->q, user->yiwork, user->yi_scatter, user->nt)); 935c4762a1bSJed Brown 936c4762a1bSJed Brown /* Compute true velocity field utrue */ 9379566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->u, &user->utrue)); 938c4762a1bSJed Brown for (i = 0; i < user->nt; i++) { 9399566063dSJacob Faibussowitsch PetscCall(VecCopy(YY, user->uxi[i])); 9409566063dSJacob Faibussowitsch PetscCall(VecScale(user->uxi[i], 150.0 * i * user->ht)); 9419566063dSJacob Faibussowitsch PetscCall(VecCopy(XX, user->uyi[i])); 9429566063dSJacob Faibussowitsch PetscCall(VecShift(user->uyi[i], -10.0)); 9439566063dSJacob Faibussowitsch PetscCall(VecScale(user->uyi[i], 15.0 * i * user->ht)); 944c4762a1bSJed Brown } 9459566063dSJacob Faibussowitsch PetscCall(Gather_uxi_uyi(user->utrue, user->uxi, user->uxi_scatter, user->uyi, user->uyi_scatter, user->nt)); 946c4762a1bSJed Brown 947c4762a1bSJed Brown /* Initial guess and reference model */ 9489566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->utrue, &user->ur)); 949c4762a1bSJed Brown for (i = 0; i < user->nt; i++) { 9509566063dSJacob Faibussowitsch PetscCall(VecCopy(XX, user->uxi[i])); 9519566063dSJacob Faibussowitsch PetscCall(VecShift(user->uxi[i], i * user->ht)); 9529566063dSJacob Faibussowitsch PetscCall(VecCopy(YY, user->uyi[i])); 9539566063dSJacob Faibussowitsch PetscCall(VecShift(user->uyi[i], -i * user->ht)); 954c4762a1bSJed Brown } 9559566063dSJacob Faibussowitsch PetscCall(Gather_uxi_uyi(user->ur, user->uxi, user->uxi_scatter, user->uyi, user->uyi_scatter, user->nt)); 956c4762a1bSJed Brown 957c4762a1bSJed Brown /* Generate regularization matrix L */ 9589566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->LT)); 9599566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->LT, PETSC_DECIDE, PETSC_DECIDE, 2 * n * user->nt, n * user->nt)); 9609566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->LT)); 9619566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->LT, 1, NULL, 1, NULL)); 9629566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->LT, 1, NULL)); 9639566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->LT, &istart, &iend)); 964c4762a1bSJed Brown 965c4762a1bSJed Brown for (i = istart; i < iend; i++) { 966c4762a1bSJed Brown iblock = (i + n) / (2 * n); 967c4762a1bSJed Brown j = i - iblock * n; 9689566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->LT, 1, &i, 1, &j, &user->gamma, INSERT_VALUES)); 969c4762a1bSJed Brown } 970c4762a1bSJed Brown 9719566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->LT, MAT_FINAL_ASSEMBLY)); 9729566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->LT, MAT_FINAL_ASSEMBLY)); 973c4762a1bSJed Brown 9749566063dSJacob Faibussowitsch PetscCall(MatTranspose(user->LT, MAT_INITIAL_MATRIX, &user->L)); 975c4762a1bSJed Brown 976c4762a1bSJed Brown /* Build work vectors and matrices */ 9779566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->lwork)); 9789566063dSJacob Faibussowitsch PetscCall(VecSetType(user->lwork, VECMPI)); 9799566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->lwork, PETSC_DECIDE, user->m)); 9809566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->lwork)); 981c4762a1bSJed Brown 9829566063dSJacob Faibussowitsch PetscCall(MatDuplicate(user->Div, MAT_SHARE_NONZERO_PATTERN, &user->Divwork)); 983c4762a1bSJed Brown 9849566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->y, &user->ywork)); 9859566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->u, &user->uwork)); 9869566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->u, &user->vwork)); 9879566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->u, &user->js_diag)); 9889566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->c, &user->cwork)); 989c4762a1bSJed Brown 990c4762a1bSJed Brown /* Create matrix-free shell user->Js for computing A*x */ 9919566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m, user->m, user, &user->Js)); 9929566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->Js, MATOP_MULT, (void (*)(void))StateMatMult)); 9939566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->Js, MATOP_DUPLICATE, (void (*)(void))StateMatDuplicate)); 9949566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->Js, MATOP_MULT_TRANSPOSE, (void (*)(void))StateMatMultTranspose)); 9959566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->Js, MATOP_GET_DIAGONAL, (void (*)(void))StateMatGetDiagonal)); 996c4762a1bSJed Brown 997c4762a1bSJed Brown /* Diagonal blocks of user->Js */ 9989566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, n, n, user, &user->JsBlock)); 9999566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->JsBlock, MATOP_MULT, (void (*)(void))StateMatBlockMult)); 10009566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->JsBlock, MATOP_MULT_TRANSPOSE, (void (*)(void))StateMatBlockMultTranspose)); 1001c4762a1bSJed Brown 1002c4762a1bSJed Brown /* Create a matrix-free shell user->JsBlockPrec for computing (U+D)\D*(L+D)\x, where JsBlock = L+D+U, 1003c4762a1bSJed Brown D is diagonal, L is strictly lower triangular, and U is strictly upper triangular. 1004c4762a1bSJed Brown This is an SOR preconditioner for user->JsBlock. */ 10059566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, n, n, user, &user->JsBlockPrec)); 10069566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->JsBlockPrec, MATOP_MULT, (void (*)(void))StateMatBlockPrecMult)); 10079566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->JsBlockPrec, MATOP_MULT_TRANSPOSE, (void (*)(void))StateMatBlockPrecMultTranspose)); 1008c4762a1bSJed Brown 1009c4762a1bSJed Brown /* Create a matrix-free shell user->Jd for computing B*x */ 10109566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m, user->n - user->m, user, &user->Jd)); 10119566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->Jd, MATOP_MULT, (void (*)(void))DesignMatMult)); 10129566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->Jd, MATOP_MULT_TRANSPOSE, (void (*)(void))DesignMatMultTranspose)); 1013c4762a1bSJed Brown 1014c4762a1bSJed Brown /* User-defined routines for computing user->Js\x and user->Js^T\x*/ 10159566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m, user->m, user, &user->JsInv)); 10169566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->JsInv, MATOP_MULT, (void (*)(void))StateMatInvMult)); 10179566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->JsInv, MATOP_MULT_TRANSPOSE, (void (*)(void))StateMatInvTransposeMult)); 1018c4762a1bSJed Brown 1019c4762a1bSJed Brown /* Build matrices for SOR preconditioner */ 10209566063dSJacob Faibussowitsch PetscCall(Scatter_uxi_uyi(user->u, user->uxi, user->uxi_scatter, user->uyi, user->uyi_scatter, user->nt)); 10219566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(5 * n, &user->C)); 10229566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2 * n, &user->Cwork)); 1023c4762a1bSJed Brown for (i = 0; i < user->nt; i++) { 10249566063dSJacob Faibussowitsch PetscCall(MatDuplicate(user->Divxy[0], MAT_COPY_VALUES, &user->C[i])); 10259566063dSJacob Faibussowitsch PetscCall(MatDuplicate(user->Divxy[1], MAT_COPY_VALUES, &user->Cwork[i])); 1026c4762a1bSJed Brown 10279566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(user->C[i], NULL, user->uxi[i])); 10289566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(user->Cwork[i], NULL, user->uyi[i])); 10299566063dSJacob Faibussowitsch PetscCall(MatAXPY(user->C[i], 1.0, user->Cwork[i], DIFFERENT_NONZERO_PATTERN)); 10309566063dSJacob Faibussowitsch PetscCall(MatScale(user->C[i], user->ht)); 10319566063dSJacob Faibussowitsch PetscCall(MatShift(user->C[i], 1.0)); 1032c4762a1bSJed Brown } 1033c4762a1bSJed Brown 1034c4762a1bSJed Brown /* Solver options and tolerances */ 10359566063dSJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_WORLD, &user->solver)); 10369566063dSJacob Faibussowitsch PetscCall(KSPSetType(user->solver, KSPGMRES)); 10379566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(user->solver, user->JsBlock, user->JsBlockPrec)); 10389566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(user->solver, 1e-4, 1e-20, 1e3, 500)); 10399566063dSJacob Faibussowitsch /* PetscCall(KSPSetTolerances(user->solver,1e-8,1e-16,1e3,500)); */ 10409566063dSJacob Faibussowitsch PetscCall(KSPGetPC(user->solver, &user->prec)); 10419566063dSJacob Faibussowitsch PetscCall(PCSetType(user->prec, PCSHELL)); 1042c4762a1bSJed Brown 10439566063dSJacob Faibussowitsch PetscCall(PCShellSetApply(user->prec, StateMatBlockPrecMult)); 10449566063dSJacob Faibussowitsch PetscCall(PCShellSetApplyTranspose(user->prec, StateMatBlockPrecMultTranspose)); 10459566063dSJacob Faibussowitsch PetscCall(PCShellSetContext(user->prec, user)); 1046c4762a1bSJed Brown 1047c4762a1bSJed Brown /* Compute true state function yt given ut */ 10489566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->ytrue)); 10499566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->ytrue, PETSC_DECIDE, n * user->nt)); 10509566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->ytrue)); 1051c4762a1bSJed Brown user->c_formed = PETSC_TRUE; 10529566063dSJacob Faibussowitsch PetscCall(VecCopy(user->utrue, user->u)); /* Set u=utrue temporarily for StateMatInv */ 10539566063dSJacob Faibussowitsch PetscCall(VecSet(user->ytrue, 0.0)); /* Initial guess */ 10549566063dSJacob Faibussowitsch PetscCall(StateMatInvMult(user->Js, user->q, user->ytrue)); 10559566063dSJacob Faibussowitsch PetscCall(VecCopy(user->ur, user->u)); /* Reset u=ur */ 1056c4762a1bSJed Brown 1057c4762a1bSJed Brown /* Initial guess y0 for state given u0 */ 10589566063dSJacob Faibussowitsch PetscCall(StateMatInvMult(user->Js, user->q, user->y)); 1059c4762a1bSJed Brown 1060c4762a1bSJed Brown /* Data discretization */ 10619566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Q)); 10629566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Q, PETSC_DECIDE, PETSC_DECIDE, user->mx * user->mx, user->m)); 10639566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Q)); 10649566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->Q, 0, NULL, 1, NULL)); 10659566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->Q, 1, NULL)); 1066c4762a1bSJed Brown 10679566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->Q, &istart, &iend)); 1068c4762a1bSJed Brown 1069c4762a1bSJed Brown for (i = istart; i < iend; i++) { 1070c4762a1bSJed Brown j = i + user->m - user->mx * user->mx; 10719566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Q, 1, &i, 1, &j, &one, INSERT_VALUES)); 1072c4762a1bSJed Brown } 1073c4762a1bSJed Brown 10749566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Q, MAT_FINAL_ASSEMBLY)); 10759566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Q, MAT_FINAL_ASSEMBLY)); 1076c4762a1bSJed Brown 10779566063dSJacob Faibussowitsch PetscCall(MatTranspose(user->Q, MAT_INITIAL_MATRIX, &user->QT)); 1078c4762a1bSJed Brown 10799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&XX)); 10809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&YY)); 10819566063dSJacob Faibussowitsch PetscCall(VecDestroy(&XXwork)); 10829566063dSJacob Faibussowitsch PetscCall(VecDestroy(&YYwork)); 10839566063dSJacob Faibussowitsch PetscCall(VecDestroy(&yi)); 10849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&uxi)); 10859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ui)); 10869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bc)); 1087c4762a1bSJed Brown 1088c4762a1bSJed Brown /* Now that initial conditions have been set, let the user pass tolerance options to the KSP solver */ 10899566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(user->solver)); 1090c4762a1bSJed Brown PetscFunctionReturn(0); 1091c4762a1bSJed Brown } 1092c4762a1bSJed Brown 10939371c9d4SSatish Balay PetscErrorCode HyperbolicDestroy(AppCtx *user) { 1094c4762a1bSJed Brown PetscInt i; 1095c4762a1bSJed Brown 1096c4762a1bSJed Brown PetscFunctionBegin; 10979566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Q)); 10989566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->QT)); 10999566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Div)); 11009566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Divwork)); 11019566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Grad)); 11029566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->L)); 11039566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->LT)); 11049566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&user->solver)); 11059566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Js)); 11069566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Jd)); 11079566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->JsBlockPrec)); 11089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->JsInv)); 11099566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->JsBlock)); 11109566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Divxy[0])); 11119566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Divxy[1])); 11129566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Gradxy[0])); 11139566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Gradxy[1])); 11149566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->M)); 1115c4762a1bSJed Brown for (i = 0; i < user->nt; i++) { 11169566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->C[i])); 11179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Cwork[i])); 1118c4762a1bSJed Brown } 11199566063dSJacob Faibussowitsch PetscCall(PetscFree(user->C)); 11209566063dSJacob Faibussowitsch PetscCall(PetscFree(user->Cwork)); 11219566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->u)); 11229566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->uwork)); 11239566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->vwork)); 11249566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->utrue)); 11259566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->y)); 11269566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->ywork)); 11279566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->ytrue)); 11289566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->nt, &user->yi)); 11299566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->nt, &user->yiwork)); 11309566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->nt, &user->ziwork)); 11319566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->nt, &user->uxi)); 11329566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->nt, &user->uyi)); 11339566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->nt, &user->uxiwork)); 11349566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->nt, &user->uyiwork)); 11359566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->nt, &user->ui)); 11369566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->nt, &user->uiwork)); 11379566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->c)); 11389566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->cwork)); 11399566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->ur)); 11409566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->q)); 11419566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->d)); 11429566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->dwork)); 11439566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->lwork)); 11449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->js_diag)); 11459566063dSJacob Faibussowitsch PetscCall(ISDestroy(&user->s_is)); 11469566063dSJacob Faibussowitsch PetscCall(ISDestroy(&user->d_is)); 11479566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->state_scatter)); 11489566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->design_scatter)); 1149c4762a1bSJed Brown for (i = 0; i < user->nt; i++) { 11509566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->uxi_scatter[i])); 11519566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->uyi_scatter[i])); 11529566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->ux_scatter[i])); 11539566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->uy_scatter[i])); 11549566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->ui_scatter[i])); 11559566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->yi_scatter[i])); 1156c4762a1bSJed Brown } 11579566063dSJacob Faibussowitsch PetscCall(PetscFree(user->uxi_scatter)); 11589566063dSJacob Faibussowitsch PetscCall(PetscFree(user->uyi_scatter)); 11599566063dSJacob Faibussowitsch PetscCall(PetscFree(user->ux_scatter)); 11609566063dSJacob Faibussowitsch PetscCall(PetscFree(user->uy_scatter)); 11619566063dSJacob Faibussowitsch PetscCall(PetscFree(user->ui_scatter)); 11629566063dSJacob Faibussowitsch PetscCall(PetscFree(user->yi_scatter)); 1163c4762a1bSJed Brown PetscFunctionReturn(0); 1164c4762a1bSJed Brown } 1165c4762a1bSJed Brown 11669371c9d4SSatish Balay PetscErrorCode HyperbolicMonitor(Tao tao, void *ptr) { 1167c4762a1bSJed Brown Vec X; 1168c4762a1bSJed Brown PetscReal unorm, ynorm; 1169c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 1170c4762a1bSJed Brown 1171c4762a1bSJed Brown PetscFunctionBegin; 11729566063dSJacob Faibussowitsch PetscCall(TaoGetSolution(tao, &X)); 11739566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->ywork, user->state_scatter, user->uwork, user->design_scatter)); 11749566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->ywork, -1.0, user->ytrue)); 11759566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->uwork, -1.0, user->utrue)); 11769566063dSJacob Faibussowitsch PetscCall(VecNorm(user->uwork, NORM_2, &unorm)); 11779566063dSJacob Faibussowitsch PetscCall(VecNorm(user->ywork, NORM_2, &ynorm)); 11789566063dSJacob Faibussowitsch PetscCall(PetscPrintf(MPI_COMM_WORLD, "||u-ut||=%g ||y-yt||=%g\n", (double)unorm, (double)ynorm)); 1179c4762a1bSJed Brown PetscFunctionReturn(0); 1180c4762a1bSJed Brown } 1181c4762a1bSJed Brown 1182c4762a1bSJed Brown /*TEST 1183c4762a1bSJed Brown 1184c4762a1bSJed Brown build: 1185c4762a1bSJed Brown requires: !complex 1186c4762a1bSJed Brown 1187c4762a1bSJed Brown test: 1188c4762a1bSJed Brown requires: !single 1189c4762a1bSJed Brown args: -tao_cmonitor -tao_max_funcs 10 -tao_type lcl -tao_gatol 1.e-5 1190c4762a1bSJed Brown 1191c4762a1bSJed Brown test: 1192c4762a1bSJed Brown suffix: guess_pod 1193c4762a1bSJed Brown requires: !single 1194c4762a1bSJed Brown args: -tao_cmonitor -tao_max_funcs 10 -tao_type lcl -ksp_guess_type pod -tao_gatol 1.e-5 1195c4762a1bSJed Brown 1196c4762a1bSJed Brown TEST*/ 1197