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 94*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 95*d71ae5a4SJacob Faibussowitsch { 96c4762a1bSJed Brown Vec x, x0; 97c4762a1bSJed Brown Tao tao; 98c4762a1bSJed Brown AppCtx user; 99c4762a1bSJed Brown IS is_allstate, is_alldesign; 100c4762a1bSJed Brown PetscInt lo, hi, hi2, lo2, ksp_old; 101c4762a1bSJed Brown PetscInt ntests = 1; 102c4762a1bSJed Brown PetscInt i; 103c4762a1bSJed Brown #if defined(PETSC_USE_LOG) 104c4762a1bSJed Brown PetscLogStage stages[1]; 105c4762a1bSJed Brown #endif 106c4762a1bSJed Brown 107327415f7SBarry Smith PetscFunctionBeginUser; 1089566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 109c4762a1bSJed Brown user.mx = 32; 110d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "hyperbolic example", NULL); 1119566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mx", "Number of grid points in each direction", "", user.mx, &user.mx, NULL)); 112c4762a1bSJed Brown user.nt = 16; 1139566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-nt", "Number of time steps", "", user.nt, &user.nt, NULL)); 114c4762a1bSJed Brown user.ndata = 64; 1159566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ndata", "Numbers of data points per sample", "", user.ndata, &user.ndata, NULL)); 116c4762a1bSJed Brown user.alpha = 10.0; 1179566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-alpha", "Regularization parameter", "", user.alpha, &user.alpha, NULL)); 118c4762a1bSJed Brown user.T = 1.0 / 32.0; 1199566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-Tfinal", "Final time", "", user.T, &user.T, NULL)); 1209566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ntests", "Number of times to repeat TaoSolve", "", ntests, &ntests, NULL)); 121d0609cedSBarry Smith PetscOptionsEnd(); 122c4762a1bSJed Brown 123c4762a1bSJed Brown user.m = user.mx * user.mx * user.nt; /* number of constraints */ 124c4762a1bSJed Brown user.n = user.mx * user.mx * 3 * user.nt; /* number of variables */ 125c4762a1bSJed Brown user.ht = user.T / user.nt; /* Time step */ 126c4762a1bSJed Brown user.gamma = user.T * user.ht / (user.mx * user.mx); 127c4762a1bSJed Brown 1289566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user.u)); 1299566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user.y)); 1309566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user.c)); 1319566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user.u, PETSC_DECIDE, user.n - user.m)); 1329566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user.y, PETSC_DECIDE, user.m)); 1339566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user.c, PETSC_DECIDE, user.m)); 1349566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user.u)); 1359566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user.y)); 1369566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user.c)); 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* Create scatters for reduced spaces. 139c4762a1bSJed Brown If the state vector y and design vector u are partitioned as 140c4762a1bSJed Brown [y_1; y_2; ...; y_np] and [u_1; u_2; ...; u_np] (with np = # of processors), 141c4762a1bSJed Brown then the solution vector x is organized as 142c4762a1bSJed Brown [y_1; u_1; y_2; u_2; ...; y_np; u_np]. 143c4762a1bSJed Brown The index sets user.s_is and user.d_is correspond to the indices of the 144c4762a1bSJed Brown state and design variables owned by the current processor. 145c4762a1bSJed Brown */ 1469566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &x)); 147c4762a1bSJed Brown 1489566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user.y, &lo, &hi)); 1499566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user.u, &lo2, &hi2)); 150c4762a1bSJed Brown 1519566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_allstate)); 1529566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo + lo2, 1, &user.s_is)); 153c4762a1bSJed Brown 1549566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi2 - lo2, lo2, 1, &is_alldesign)); 1559566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi2 - lo2, hi + lo2, 1, &user.d_is)); 156c4762a1bSJed Brown 1579566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x, hi - lo + hi2 - lo2, user.n)); 1589566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x)); 159c4762a1bSJed Brown 1609566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(x, user.s_is, user.y, is_allstate, &user.state_scatter)); 1619566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(x, user.d_is, user.u, is_alldesign, &user.design_scatter)); 1629566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_alldesign)); 1639566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_allstate)); 164c4762a1bSJed Brown 165c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 1669566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao)); 1679566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOLCL)); 168c4762a1bSJed Brown 169c4762a1bSJed Brown /* Set up initial vectors and matrices */ 1709566063dSJacob Faibussowitsch PetscCall(HyperbolicInitialize(&user)); 171c4762a1bSJed Brown 1729566063dSJacob Faibussowitsch PetscCall(Gather(x, user.y, user.state_scatter, user.u, user.design_scatter)); 1739566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &x0)); 1749566063dSJacob Faibussowitsch PetscCall(VecCopy(x, x0)); 175c4762a1bSJed Brown 176c4762a1bSJed Brown /* Set solution vector with an initial guess */ 1779566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, x)); 1789566063dSJacob Faibussowitsch PetscCall(TaoSetObjective(tao, FormFunction, &user)); 1799566063dSJacob Faibussowitsch PetscCall(TaoSetGradient(tao, NULL, FormGradient, &user)); 1809566063dSJacob Faibussowitsch PetscCall(TaoSetConstraintsRoutine(tao, user.c, FormConstraints, &user)); 1819566063dSJacob Faibussowitsch PetscCall(TaoSetJacobianStateRoutine(tao, user.Js, user.Js, user.JsInv, FormJacobianState, &user)); 1829566063dSJacob Faibussowitsch PetscCall(TaoSetJacobianDesignRoutine(tao, user.Jd, FormJacobianDesign, &user)); 1839566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 1849566063dSJacob Faibussowitsch PetscCall(TaoSetStateDesignIS(tao, user.s_is, user.d_is)); 185c4762a1bSJed Brown 186c4762a1bSJed Brown /* SOLVE THE APPLICATION */ 1879566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("Trials", &stages[0])); 1889566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stages[0])); 189c4762a1bSJed Brown user.ksp_its_initial = user.ksp_its; 190c4762a1bSJed Brown ksp_old = user.ksp_its; 191c4762a1bSJed Brown for (i = 0; i < ntests; i++) { 1929566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 19363a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "KSP Iterations = %" PetscInt_FMT "\n", user.ksp_its - ksp_old)); 1949566063dSJacob Faibussowitsch PetscCall(VecCopy(x0, x)); 1959566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, x)); 196c4762a1bSJed Brown } 1979566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 1989566063dSJacob Faibussowitsch PetscCall(PetscBarrier((PetscObject)x)); 1999566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "KSP iterations within initialization: ")); 20063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "\n", user.ksp_its_initial)); 20163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Total KSP iterations over %" PetscInt_FMT " trial(s): ", ntests)); 20263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "\n", user.ksp_its)); 2039566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "KSP iterations per trial: ")); 20463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "\n", (user.ksp_its - user.ksp_its_initial) / ntests)); 205c4762a1bSJed Brown 2069566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 2079566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 2089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x0)); 2099566063dSJacob Faibussowitsch PetscCall(HyperbolicDestroy(&user)); 2109566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 211b122ec5aSJacob Faibussowitsch return 0; 212c4762a1bSJed Brown } 213c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 214c4762a1bSJed Brown /* 215c4762a1bSJed Brown dwork = Qy - d 216c4762a1bSJed Brown lwork = L*(u-ur).^2 217c4762a1bSJed Brown f = 1/2 * (dwork.dork + alpha*y.lwork) 218c4762a1bSJed Brown */ 219*d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(Tao tao, Vec X, PetscReal *f, void *ptr) 220*d71ae5a4SJacob Faibussowitsch { 221c4762a1bSJed Brown PetscReal d1 = 0, d2 = 0; 222c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 223c4762a1bSJed Brown 224c4762a1bSJed Brown PetscFunctionBegin; 2259566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter)); 2269566063dSJacob Faibussowitsch PetscCall(MatMult(user->Q, user->y, user->dwork)); 2279566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->dwork, -1.0, user->d)); 2289566063dSJacob Faibussowitsch PetscCall(VecDot(user->dwork, user->dwork, &d1)); 229c4762a1bSJed Brown 2309566063dSJacob Faibussowitsch PetscCall(VecWAXPY(user->uwork, -1.0, user->ur, user->u)); 2319566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uwork, user->uwork, user->uwork)); 2329566063dSJacob Faibussowitsch PetscCall(MatMult(user->L, user->uwork, user->lwork)); 2339566063dSJacob Faibussowitsch PetscCall(VecDot(user->y, user->lwork, &d2)); 234c4762a1bSJed Brown *f = 0.5 * (d1 + user->alpha * d2); 235c4762a1bSJed Brown PetscFunctionReturn(0); 236c4762a1bSJed Brown } 237c4762a1bSJed Brown 238c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 239c4762a1bSJed Brown /* 240c4762a1bSJed Brown state: g_s = Q' *(Qy - d) + 0.5*alpha*L*(u-ur).^2 241c4762a1bSJed Brown design: g_d = alpha*(L'y).*(u-ur) 242c4762a1bSJed Brown */ 243*d71ae5a4SJacob Faibussowitsch PetscErrorCode FormGradient(Tao tao, Vec X, Vec G, void *ptr) 244*d71ae5a4SJacob Faibussowitsch { 245c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 246c4762a1bSJed Brown 247c4762a1bSJed Brown PetscFunctionBegin; 2489566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter)); 2499566063dSJacob Faibussowitsch PetscCall(MatMult(user->Q, user->y, user->dwork)); 2509566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->dwork, -1.0, user->d)); 251c4762a1bSJed Brown 2529566063dSJacob Faibussowitsch PetscCall(MatMult(user->QT, user->dwork, user->ywork)); 253c4762a1bSJed Brown 2549566063dSJacob Faibussowitsch PetscCall(MatMult(user->LT, user->y, user->uwork)); 2559566063dSJacob Faibussowitsch PetscCall(VecWAXPY(user->vwork, -1.0, user->ur, user->u)); 2569566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uwork, user->vwork, user->uwork)); 2579566063dSJacob Faibussowitsch PetscCall(VecScale(user->uwork, user->alpha)); 258c4762a1bSJed Brown 2599566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->vwork, user->vwork, user->vwork)); 2609566063dSJacob Faibussowitsch PetscCall(MatMult(user->L, user->vwork, user->lwork)); 2619566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->ywork, 0.5 * user->alpha, user->lwork)); 262c4762a1bSJed Brown 2639566063dSJacob Faibussowitsch PetscCall(Gather(G, user->ywork, user->state_scatter, user->uwork, user->design_scatter)); 264c4762a1bSJed Brown PetscFunctionReturn(0); 265c4762a1bSJed Brown } 266c4762a1bSJed Brown 267*d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr) 268*d71ae5a4SJacob Faibussowitsch { 269c4762a1bSJed Brown PetscReal d1, d2; 270c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 271c4762a1bSJed Brown 272c4762a1bSJed Brown PetscFunctionBegin; 2739566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter)); 2749566063dSJacob Faibussowitsch PetscCall(MatMult(user->Q, user->y, user->dwork)); 2759566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->dwork, -1.0, user->d)); 276c4762a1bSJed Brown 2779566063dSJacob Faibussowitsch PetscCall(MatMult(user->QT, user->dwork, user->ywork)); 278c4762a1bSJed Brown 2799566063dSJacob Faibussowitsch PetscCall(VecDot(user->dwork, user->dwork, &d1)); 280c4762a1bSJed Brown 2819566063dSJacob Faibussowitsch PetscCall(MatMult(user->LT, user->y, user->uwork)); 2829566063dSJacob Faibussowitsch PetscCall(VecWAXPY(user->vwork, -1.0, user->ur, user->u)); 2839566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uwork, user->vwork, user->uwork)); 2849566063dSJacob Faibussowitsch PetscCall(VecScale(user->uwork, user->alpha)); 285c4762a1bSJed Brown 2869566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->vwork, user->vwork, user->vwork)); 2879566063dSJacob Faibussowitsch PetscCall(MatMult(user->L, user->vwork, user->lwork)); 2889566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->ywork, 0.5 * user->alpha, user->lwork)); 289c4762a1bSJed Brown 2909566063dSJacob Faibussowitsch PetscCall(VecDot(user->y, user->lwork, &d2)); 291c4762a1bSJed Brown 292c4762a1bSJed Brown *f = 0.5 * (d1 + user->alpha * d2); 2939566063dSJacob Faibussowitsch PetscCall(Gather(G, user->ywork, user->state_scatter, user->uwork, user->design_scatter)); 294c4762a1bSJed Brown PetscFunctionReturn(0); 295c4762a1bSJed Brown } 296c4762a1bSJed Brown 297c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 298c4762a1bSJed Brown /* A 299c4762a1bSJed Brown MatShell object 300c4762a1bSJed Brown */ 301*d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobianState(Tao tao, Vec X, Mat J, Mat JPre, Mat JInv, void *ptr) 302*d71ae5a4SJacob Faibussowitsch { 303c4762a1bSJed Brown PetscInt i; 304c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 305c4762a1bSJed Brown 306c4762a1bSJed Brown PetscFunctionBegin; 3079566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter)); 3089566063dSJacob Faibussowitsch PetscCall(Scatter_yi(user->u, user->ui, user->ui_scatter, user->nt)); 3099566063dSJacob Faibussowitsch PetscCall(Scatter_uxi_uyi(user->u, user->uxi, user->uxi_scatter, user->uyi, user->uyi_scatter, user->nt)); 310c4762a1bSJed Brown for (i = 0; i < user->nt; i++) { 3119566063dSJacob Faibussowitsch PetscCall(MatCopy(user->Divxy[0], user->C[i], SUBSET_NONZERO_PATTERN)); 3129566063dSJacob Faibussowitsch PetscCall(MatCopy(user->Divxy[1], user->Cwork[i], SAME_NONZERO_PATTERN)); 313c4762a1bSJed Brown 3149566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(user->C[i], NULL, user->uxi[i])); 3159566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(user->Cwork[i], NULL, user->uyi[i])); 3169566063dSJacob Faibussowitsch PetscCall(MatAXPY(user->C[i], 1.0, user->Cwork[i], SUBSET_NONZERO_PATTERN)); 3179566063dSJacob Faibussowitsch PetscCall(MatScale(user->C[i], user->ht)); 3189566063dSJacob Faibussowitsch PetscCall(MatShift(user->C[i], 1.0)); 319c4762a1bSJed Brown } 320c4762a1bSJed Brown PetscFunctionReturn(0); 321c4762a1bSJed Brown } 322c4762a1bSJed Brown 323c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 324c4762a1bSJed Brown /* B */ 325*d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobianDesign(Tao tao, Vec X, Mat J, void *ptr) 326*d71ae5a4SJacob Faibussowitsch { 327c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 328c4762a1bSJed Brown 329c4762a1bSJed Brown PetscFunctionBegin; 3309566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter)); 331c4762a1bSJed Brown PetscFunctionReturn(0); 332c4762a1bSJed Brown } 333c4762a1bSJed Brown 334*d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatMult(Mat J_shell, Vec X, Vec Y) 335*d71ae5a4SJacob Faibussowitsch { 336c4762a1bSJed Brown PetscInt i; 337c4762a1bSJed Brown AppCtx *user; 338c4762a1bSJed Brown 339c4762a1bSJed Brown PetscFunctionBegin; 3409566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user)); 3419566063dSJacob Faibussowitsch PetscCall(Scatter_yi(X, user->yi, user->yi_scatter, user->nt)); 342c4762a1bSJed Brown user->block_index = 0; 3439566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock, user->yi[0], user->yiwork[0])); 344c4762a1bSJed Brown 345c4762a1bSJed Brown for (i = 1; i < user->nt; i++) { 346c4762a1bSJed Brown user->block_index = i; 3479566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock, user->yi[i], user->yiwork[i])); 3489566063dSJacob Faibussowitsch PetscCall(MatMult(user->M, user->yi[i - 1], user->ziwork[i - 1])); 3499566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->yiwork[i], -1.0, user->ziwork[i - 1])); 350c4762a1bSJed Brown } 3519566063dSJacob Faibussowitsch PetscCall(Gather_yi(Y, user->yiwork, user->yi_scatter, user->nt)); 352c4762a1bSJed Brown PetscFunctionReturn(0); 353c4762a1bSJed Brown } 354c4762a1bSJed Brown 355*d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatMultTranspose(Mat J_shell, Vec X, Vec Y) 356*d71ae5a4SJacob Faibussowitsch { 357c4762a1bSJed Brown PetscInt i; 358c4762a1bSJed Brown AppCtx *user; 359c4762a1bSJed Brown 360c4762a1bSJed Brown PetscFunctionBegin; 3619566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user)); 3629566063dSJacob Faibussowitsch PetscCall(Scatter_yi(X, user->yi, user->yi_scatter, user->nt)); 363c4762a1bSJed Brown 364c4762a1bSJed Brown for (i = 0; i < user->nt - 1; i++) { 365c4762a1bSJed Brown user->block_index = i; 3669566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->JsBlock, user->yi[i], user->yiwork[i])); 3679566063dSJacob Faibussowitsch PetscCall(MatMult(user->M, user->yi[i + 1], user->ziwork[i + 1])); 3689566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->yiwork[i], -1.0, user->ziwork[i + 1])); 369c4762a1bSJed Brown } 370c4762a1bSJed Brown 371c4762a1bSJed Brown i = user->nt - 1; 372c4762a1bSJed Brown user->block_index = i; 3739566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->JsBlock, user->yi[i], user->yiwork[i])); 3749566063dSJacob Faibussowitsch PetscCall(Gather_yi(Y, user->yiwork, user->yi_scatter, user->nt)); 375c4762a1bSJed Brown PetscFunctionReturn(0); 376c4762a1bSJed Brown } 377c4762a1bSJed Brown 378*d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatBlockMult(Mat J_shell, Vec X, Vec Y) 379*d71ae5a4SJacob Faibussowitsch { 380c4762a1bSJed Brown PetscInt i; 381c4762a1bSJed Brown AppCtx *user; 382c4762a1bSJed Brown 383c4762a1bSJed Brown PetscFunctionBegin; 3849566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user)); 385c4762a1bSJed Brown i = user->block_index; 3869566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uxiwork[i], X, user->uxi[i])); 3879566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uyiwork[i], X, user->uyi[i])); 3889566063dSJacob Faibussowitsch PetscCall(Gather(user->uiwork[i], user->uxiwork[i], user->ux_scatter[i], user->uyiwork[i], user->uy_scatter[i])); 3899566063dSJacob Faibussowitsch PetscCall(MatMult(user->Div, user->uiwork[i], Y)); 3909566063dSJacob Faibussowitsch PetscCall(VecAYPX(Y, user->ht, X)); 391c4762a1bSJed Brown PetscFunctionReturn(0); 392c4762a1bSJed Brown } 393c4762a1bSJed Brown 394*d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatBlockMultTranspose(Mat J_shell, Vec X, Vec Y) 395*d71ae5a4SJacob Faibussowitsch { 396c4762a1bSJed Brown PetscInt i; 397c4762a1bSJed Brown AppCtx *user; 398c4762a1bSJed Brown 399c4762a1bSJed Brown PetscFunctionBegin; 4009566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user)); 401c4762a1bSJed Brown i = user->block_index; 4029566063dSJacob Faibussowitsch PetscCall(MatMult(user->Grad, X, user->uiwork[i])); 4039566063dSJacob Faibussowitsch PetscCall(Scatter(user->uiwork[i], user->uxiwork[i], user->ux_scatter[i], user->uyiwork[i], user->uy_scatter[i])); 4049566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uxiwork[i], user->uxi[i], user->uxiwork[i])); 4059566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uyiwork[i], user->uyi[i], user->uyiwork[i])); 4069566063dSJacob Faibussowitsch PetscCall(VecWAXPY(Y, 1.0, user->uxiwork[i], user->uyiwork[i])); 4079566063dSJacob Faibussowitsch PetscCall(VecAYPX(Y, user->ht, X)); 408c4762a1bSJed Brown PetscFunctionReturn(0); 409c4762a1bSJed Brown } 410c4762a1bSJed Brown 411*d71ae5a4SJacob Faibussowitsch PetscErrorCode DesignMatMult(Mat J_shell, Vec X, Vec Y) 412*d71ae5a4SJacob Faibussowitsch { 413c4762a1bSJed Brown PetscInt i; 414c4762a1bSJed Brown AppCtx *user; 415c4762a1bSJed Brown 416c4762a1bSJed Brown PetscFunctionBegin; 4179566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user)); 4189566063dSJacob Faibussowitsch PetscCall(Scatter_yi(user->y, user->yi, user->yi_scatter, user->nt)); 4199566063dSJacob Faibussowitsch PetscCall(Scatter_uxi_uyi(X, user->uxiwork, user->uxi_scatter, user->uyiwork, user->uyi_scatter, user->nt)); 420c4762a1bSJed Brown for (i = 0; i < user->nt; i++) { 4219566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uxiwork[i], user->yi[i], user->uxiwork[i])); 4229566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uyiwork[i], user->yi[i], user->uyiwork[i])); 4239566063dSJacob Faibussowitsch PetscCall(Gather(user->uiwork[i], user->uxiwork[i], user->ux_scatter[i], user->uyiwork[i], user->uy_scatter[i])); 4249566063dSJacob Faibussowitsch PetscCall(MatMult(user->Div, user->uiwork[i], user->ziwork[i])); 4259566063dSJacob Faibussowitsch PetscCall(VecScale(user->ziwork[i], user->ht)); 426c4762a1bSJed Brown } 4279566063dSJacob Faibussowitsch PetscCall(Gather_yi(Y, user->ziwork, user->yi_scatter, user->nt)); 428c4762a1bSJed Brown PetscFunctionReturn(0); 429c4762a1bSJed Brown } 430c4762a1bSJed Brown 431*d71ae5a4SJacob Faibussowitsch PetscErrorCode DesignMatMultTranspose(Mat J_shell, Vec X, Vec Y) 432*d71ae5a4SJacob Faibussowitsch { 433c4762a1bSJed Brown PetscInt i; 434c4762a1bSJed Brown AppCtx *user; 435c4762a1bSJed Brown 436c4762a1bSJed Brown PetscFunctionBegin; 4379566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user)); 4389566063dSJacob Faibussowitsch PetscCall(Scatter_yi(user->y, user->yi, user->yi_scatter, user->nt)); 4399566063dSJacob Faibussowitsch PetscCall(Scatter_yi(X, user->yiwork, user->yi_scatter, user->nt)); 440c4762a1bSJed Brown for (i = 0; i < user->nt; i++) { 4419566063dSJacob Faibussowitsch PetscCall(MatMult(user->Grad, user->yiwork[i], user->uiwork[i])); 4429566063dSJacob Faibussowitsch PetscCall(Scatter(user->uiwork[i], user->uxiwork[i], user->ux_scatter[i], user->uyiwork[i], user->uy_scatter[i])); 4439566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uxiwork[i], user->yi[i], user->uxiwork[i])); 4449566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uyiwork[i], user->yi[i], user->uyiwork[i])); 4459566063dSJacob Faibussowitsch PetscCall(Gather(user->uiwork[i], user->uxiwork[i], user->ux_scatter[i], user->uyiwork[i], user->uy_scatter[i])); 4469566063dSJacob Faibussowitsch PetscCall(VecScale(user->uiwork[i], user->ht)); 447c4762a1bSJed Brown } 4489566063dSJacob Faibussowitsch PetscCall(Gather_yi(Y, user->uiwork, user->ui_scatter, user->nt)); 449c4762a1bSJed Brown PetscFunctionReturn(0); 450c4762a1bSJed Brown } 451c4762a1bSJed Brown 452*d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatBlockPrecMult(PC PC_shell, Vec X, Vec Y) 453*d71ae5a4SJacob Faibussowitsch { 454c4762a1bSJed Brown PetscInt i; 455c4762a1bSJed Brown AppCtx *user; 456c4762a1bSJed Brown 457c4762a1bSJed Brown PetscFunctionBegin; 4589566063dSJacob Faibussowitsch PetscCall(PCShellGetContext(PC_shell, &user)); 459c4762a1bSJed Brown i = user->block_index; 460c4762a1bSJed Brown if (user->c_formed) { 4619566063dSJacob Faibussowitsch PetscCall(MatSOR(user->C[i], X, 1.0, (MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_SYMMETRIC_SWEEP), 0.0, 1, 1, Y)); 462c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not formed"); 463c4762a1bSJed Brown PetscFunctionReturn(0); 464c4762a1bSJed Brown } 465c4762a1bSJed Brown 466*d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatBlockPrecMultTranspose(PC PC_shell, Vec X, Vec Y) 467*d71ae5a4SJacob Faibussowitsch { 468c4762a1bSJed Brown PetscInt i; 469c4762a1bSJed Brown AppCtx *user; 470c4762a1bSJed Brown 471c4762a1bSJed Brown PetscFunctionBegin; 4729566063dSJacob Faibussowitsch PetscCall(PCShellGetContext(PC_shell, &user)); 473c4762a1bSJed Brown 474c4762a1bSJed Brown i = user->block_index; 475c4762a1bSJed Brown if (user->c_formed) { 4769566063dSJacob Faibussowitsch PetscCall(MatSOR(user->C[i], X, 1.0, (MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_SYMMETRIC_SWEEP), 0.0, 1, 1, Y)); 477c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not formed"); 478c4762a1bSJed Brown PetscFunctionReturn(0); 479c4762a1bSJed Brown } 480c4762a1bSJed Brown 481*d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatInvMult(Mat J_shell, Vec X, Vec Y) 482*d71ae5a4SJacob Faibussowitsch { 483c4762a1bSJed Brown AppCtx *user; 484c4762a1bSJed Brown PetscInt its, i; 485c4762a1bSJed Brown 486c4762a1bSJed Brown PetscFunctionBegin; 4879566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user)); 488c4762a1bSJed Brown 489c4762a1bSJed Brown if (Y == user->ytrue) { 490c4762a1bSJed Brown /* First solve is done using true solution to set up problem */ 4919566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(user->solver, 1e-4, 1e-20, PETSC_DEFAULT, PETSC_DEFAULT)); 492c4762a1bSJed Brown } else { 4939566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(user->solver, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 494c4762a1bSJed Brown } 4959566063dSJacob Faibussowitsch PetscCall(Scatter_yi(X, user->yi, user->yi_scatter, user->nt)); 4969566063dSJacob Faibussowitsch PetscCall(Scatter_yi(Y, user->yiwork, user->yi_scatter, user->nt)); 4979566063dSJacob Faibussowitsch PetscCall(Scatter_uxi_uyi(user->u, user->uxi, user->uxi_scatter, user->uyi, user->uyi_scatter, user->nt)); 498c4762a1bSJed Brown 499c4762a1bSJed Brown user->block_index = 0; 5009566063dSJacob Faibussowitsch PetscCall(KSPSolve(user->solver, user->yi[0], user->yiwork[0])); 501c4762a1bSJed Brown 5029566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(user->solver, &its)); 503c4762a1bSJed Brown user->ksp_its = user->ksp_its + its; 504c4762a1bSJed Brown for (i = 1; i < user->nt; i++) { 5059566063dSJacob Faibussowitsch PetscCall(MatMult(user->M, user->yiwork[i - 1], user->ziwork[i - 1])); 5069566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->yi[i], 1.0, user->ziwork[i - 1])); 507c4762a1bSJed Brown user->block_index = i; 5089566063dSJacob Faibussowitsch PetscCall(KSPSolve(user->solver, user->yi[i], user->yiwork[i])); 509c4762a1bSJed Brown 5109566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(user->solver, &its)); 511c4762a1bSJed Brown user->ksp_its = user->ksp_its + its; 512c4762a1bSJed Brown } 5139566063dSJacob Faibussowitsch PetscCall(Gather_yi(Y, user->yiwork, user->yi_scatter, user->nt)); 514c4762a1bSJed Brown PetscFunctionReturn(0); 515c4762a1bSJed Brown } 516c4762a1bSJed Brown 517*d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatInvTransposeMult(Mat J_shell, Vec X, Vec Y) 518*d71ae5a4SJacob Faibussowitsch { 519c4762a1bSJed Brown AppCtx *user; 520c4762a1bSJed Brown PetscInt its, i; 521c4762a1bSJed Brown 522c4762a1bSJed Brown PetscFunctionBegin; 5239566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user)); 524c4762a1bSJed Brown 5259566063dSJacob Faibussowitsch PetscCall(Scatter_yi(X, user->yi, user->yi_scatter, user->nt)); 5269566063dSJacob Faibussowitsch PetscCall(Scatter_yi(Y, user->yiwork, user->yi_scatter, user->nt)); 5279566063dSJacob Faibussowitsch PetscCall(Scatter_uxi_uyi(user->u, user->uxi, user->uxi_scatter, user->uyi, user->uyi_scatter, user->nt)); 528c4762a1bSJed Brown 529c4762a1bSJed Brown i = user->nt - 1; 530c4762a1bSJed Brown user->block_index = i; 5319566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(user->solver, user->yi[i], user->yiwork[i])); 532c4762a1bSJed Brown 5339566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(user->solver, &its)); 534c4762a1bSJed Brown user->ksp_its = user->ksp_its + its; 535c4762a1bSJed Brown 536c4762a1bSJed Brown for (i = user->nt - 2; i >= 0; i--) { 5379566063dSJacob Faibussowitsch PetscCall(MatMult(user->M, user->yiwork[i + 1], user->ziwork[i + 1])); 5389566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->yi[i], 1.0, user->ziwork[i + 1])); 539c4762a1bSJed Brown user->block_index = i; 5409566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(user->solver, user->yi[i], user->yiwork[i])); 541c4762a1bSJed Brown 5429566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(user->solver, &its)); 543c4762a1bSJed Brown user->ksp_its = user->ksp_its + its; 544c4762a1bSJed Brown } 5459566063dSJacob Faibussowitsch PetscCall(Gather_yi(Y, user->yiwork, user->yi_scatter, user->nt)); 546c4762a1bSJed Brown PetscFunctionReturn(0); 547c4762a1bSJed Brown } 548c4762a1bSJed Brown 549*d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatDuplicate(Mat J_shell, MatDuplicateOption opt, Mat *new_shell) 550*d71ae5a4SJacob Faibussowitsch { 551c4762a1bSJed Brown AppCtx *user; 552c4762a1bSJed Brown 553c4762a1bSJed Brown PetscFunctionBegin; 5549566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user)); 555c4762a1bSJed Brown 5569566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m, user->m, user, new_shell)); 5579566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*new_shell, MATOP_MULT, (void (*)(void))StateMatMult)); 5589566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*new_shell, MATOP_DUPLICATE, (void (*)(void))StateMatDuplicate)); 5599566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*new_shell, MATOP_MULT_TRANSPOSE, (void (*)(void))StateMatMultTranspose)); 5609566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*new_shell, MATOP_GET_DIAGONAL, (void (*)(void))StateMatGetDiagonal)); 561c4762a1bSJed Brown PetscFunctionReturn(0); 562c4762a1bSJed Brown } 563c4762a1bSJed Brown 564*d71ae5a4SJacob Faibussowitsch PetscErrorCode StateMatGetDiagonal(Mat J_shell, Vec X) 565*d71ae5a4SJacob Faibussowitsch { 566c4762a1bSJed Brown AppCtx *user; 567c4762a1bSJed Brown 568c4762a1bSJed Brown PetscFunctionBegin; 5699566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user)); 5709566063dSJacob Faibussowitsch PetscCall(VecCopy(user->js_diag, X)); 571c4762a1bSJed Brown PetscFunctionReturn(0); 572c4762a1bSJed Brown } 573c4762a1bSJed Brown 574*d71ae5a4SJacob Faibussowitsch PetscErrorCode FormConstraints(Tao tao, Vec X, Vec C, void *ptr) 575*d71ae5a4SJacob Faibussowitsch { 576c4762a1bSJed Brown /* con = Ay - q, A = [C(u1) 0 0 ... 0; 577c4762a1bSJed Brown -M C(u2) 0 ... 0; 578c4762a1bSJed Brown 0 -M C(u3) ... 0; 579c4762a1bSJed Brown ... ; 580c4762a1bSJed Brown 0 ... -M C(u_nt)] 581c4762a1bSJed Brown C(u) = eye + ht*Div*[diag(u1); diag(u2)] */ 582c4762a1bSJed Brown PetscInt i; 583c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 584c4762a1bSJed Brown 585c4762a1bSJed Brown PetscFunctionBegin; 5869566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter)); 5879566063dSJacob Faibussowitsch PetscCall(Scatter_yi(user->y, user->yi, user->yi_scatter, user->nt)); 5889566063dSJacob Faibussowitsch PetscCall(Scatter_uxi_uyi(user->u, user->uxi, user->uxi_scatter, user->uyi, user->uyi_scatter, user->nt)); 589c4762a1bSJed Brown 590c4762a1bSJed Brown user->block_index = 0; 5919566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock, user->yi[0], user->yiwork[0])); 592c4762a1bSJed Brown 593c4762a1bSJed Brown for (i = 1; i < user->nt; i++) { 594c4762a1bSJed Brown user->block_index = i; 5959566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock, user->yi[i], user->yiwork[i])); 5969566063dSJacob Faibussowitsch PetscCall(MatMult(user->M, user->yi[i - 1], user->ziwork[i - 1])); 5979566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->yiwork[i], -1.0, user->ziwork[i - 1])); 598c4762a1bSJed Brown } 599c4762a1bSJed Brown 6009566063dSJacob Faibussowitsch PetscCall(Gather_yi(C, user->yiwork, user->yi_scatter, user->nt)); 6019566063dSJacob Faibussowitsch PetscCall(VecAXPY(C, -1.0, user->q)); 602c4762a1bSJed Brown 603c4762a1bSJed Brown PetscFunctionReturn(0); 604c4762a1bSJed Brown } 605c4762a1bSJed Brown 606*d71ae5a4SJacob Faibussowitsch PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat) 607*d71ae5a4SJacob Faibussowitsch { 608c4762a1bSJed Brown PetscFunctionBegin; 6099566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(s_scat, x, state, INSERT_VALUES, SCATTER_FORWARD)); 6109566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(s_scat, x, state, INSERT_VALUES, SCATTER_FORWARD)); 6119566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(d_scat, x, design, INSERT_VALUES, SCATTER_FORWARD)); 6129566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(d_scat, x, design, INSERT_VALUES, SCATTER_FORWARD)); 613c4762a1bSJed Brown PetscFunctionReturn(0); 614c4762a1bSJed Brown } 615c4762a1bSJed Brown 616*d71ae5a4SJacob Faibussowitsch PetscErrorCode Scatter_uxi_uyi(Vec u, Vec *uxi, VecScatter *scatx, Vec *uyi, VecScatter *scaty, PetscInt nt) 617*d71ae5a4SJacob Faibussowitsch { 618c4762a1bSJed Brown PetscInt i; 619c4762a1bSJed Brown 620c4762a1bSJed Brown PetscFunctionBegin; 621c4762a1bSJed Brown for (i = 0; i < nt; i++) { 6229566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scatx[i], u, uxi[i], INSERT_VALUES, SCATTER_FORWARD)); 6239566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scatx[i], u, uxi[i], INSERT_VALUES, SCATTER_FORWARD)); 6249566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scaty[i], u, uyi[i], INSERT_VALUES, SCATTER_FORWARD)); 6259566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scaty[i], u, uyi[i], INSERT_VALUES, SCATTER_FORWARD)); 626c4762a1bSJed Brown } 627c4762a1bSJed Brown PetscFunctionReturn(0); 628c4762a1bSJed Brown } 629c4762a1bSJed Brown 630*d71ae5a4SJacob Faibussowitsch PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat) 631*d71ae5a4SJacob Faibussowitsch { 632c4762a1bSJed Brown PetscFunctionBegin; 6339566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(s_scat, state, x, INSERT_VALUES, SCATTER_REVERSE)); 6349566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(s_scat, state, x, INSERT_VALUES, SCATTER_REVERSE)); 6359566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(d_scat, design, x, INSERT_VALUES, SCATTER_REVERSE)); 6369566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(d_scat, design, x, INSERT_VALUES, SCATTER_REVERSE)); 637c4762a1bSJed Brown PetscFunctionReturn(0); 638c4762a1bSJed Brown } 639c4762a1bSJed Brown 640*d71ae5a4SJacob Faibussowitsch PetscErrorCode Gather_uxi_uyi(Vec u, Vec *uxi, VecScatter *scatx, Vec *uyi, VecScatter *scaty, PetscInt nt) 641*d71ae5a4SJacob Faibussowitsch { 642c4762a1bSJed Brown PetscInt i; 643c4762a1bSJed Brown 644c4762a1bSJed Brown PetscFunctionBegin; 645c4762a1bSJed Brown for (i = 0; i < nt; i++) { 6469566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scatx[i], uxi[i], u, INSERT_VALUES, SCATTER_REVERSE)); 6479566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scatx[i], uxi[i], u, INSERT_VALUES, SCATTER_REVERSE)); 6489566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scaty[i], uyi[i], u, INSERT_VALUES, SCATTER_REVERSE)); 6499566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scaty[i], uyi[i], u, INSERT_VALUES, SCATTER_REVERSE)); 650c4762a1bSJed Brown } 651c4762a1bSJed Brown PetscFunctionReturn(0); 652c4762a1bSJed Brown } 653c4762a1bSJed Brown 654*d71ae5a4SJacob Faibussowitsch PetscErrorCode Scatter_yi(Vec y, Vec *yi, VecScatter *scat, PetscInt nt) 655*d71ae5a4SJacob Faibussowitsch { 656c4762a1bSJed Brown PetscInt i; 657c4762a1bSJed Brown 658c4762a1bSJed Brown PetscFunctionBegin; 659c4762a1bSJed Brown for (i = 0; i < nt; i++) { 6609566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scat[i], y, yi[i], INSERT_VALUES, SCATTER_FORWARD)); 6619566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scat[i], y, yi[i], INSERT_VALUES, SCATTER_FORWARD)); 662c4762a1bSJed Brown } 663c4762a1bSJed Brown PetscFunctionReturn(0); 664c4762a1bSJed Brown } 665c4762a1bSJed Brown 666*d71ae5a4SJacob Faibussowitsch PetscErrorCode Gather_yi(Vec y, Vec *yi, VecScatter *scat, PetscInt nt) 667*d71ae5a4SJacob Faibussowitsch { 668c4762a1bSJed Brown PetscInt i; 669c4762a1bSJed Brown 670c4762a1bSJed Brown PetscFunctionBegin; 671c4762a1bSJed Brown for (i = 0; i < nt; i++) { 6729566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scat[i], yi[i], y, INSERT_VALUES, SCATTER_REVERSE)); 6739566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scat[i], yi[i], y, INSERT_VALUES, SCATTER_REVERSE)); 674c4762a1bSJed Brown } 675c4762a1bSJed Brown PetscFunctionReturn(0); 676c4762a1bSJed Brown } 677c4762a1bSJed Brown 678*d71ae5a4SJacob Faibussowitsch PetscErrorCode HyperbolicInitialize(AppCtx *user) 679*d71ae5a4SJacob Faibussowitsch { 680c4762a1bSJed Brown PetscInt n, i, j, linear_index, istart, iend, iblock, lo, hi; 681c4762a1bSJed Brown Vec XX, YY, XXwork, YYwork, yi, uxi, ui, bc; 682c4762a1bSJed Brown PetscReal h, sum; 683c4762a1bSJed Brown PetscScalar hinv, neg_hinv, quarter = 0.25, one = 1.0, half_hinv, neg_half_hinv; 684c4762a1bSJed Brown PetscScalar vx, vy, zero = 0.0; 685c4762a1bSJed Brown IS is_from_y, is_to_yi, is_from_u, is_to_uxi, is_to_uyi; 686c4762a1bSJed Brown 687c4762a1bSJed Brown PetscFunctionBegin; 688c4762a1bSJed Brown user->jformed = PETSC_FALSE; 689c4762a1bSJed Brown user->c_formed = PETSC_FALSE; 690c4762a1bSJed Brown 691c4762a1bSJed Brown user->ksp_its = 0; 692c4762a1bSJed Brown user->ksp_its_initial = 0; 693c4762a1bSJed Brown 694c4762a1bSJed Brown n = user->mx * user->mx; 695c4762a1bSJed Brown 696c4762a1bSJed Brown h = 1.0 / user->mx; 697c4762a1bSJed Brown hinv = user->mx; 698c4762a1bSJed Brown neg_hinv = -hinv; 699c4762a1bSJed Brown half_hinv = hinv / 2.0; 700c4762a1bSJed Brown neg_half_hinv = neg_hinv / 2.0; 701c4762a1bSJed Brown 702c4762a1bSJed Brown /* Generate Grad matrix */ 7039566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Grad)); 7049566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Grad, PETSC_DECIDE, PETSC_DECIDE, 2 * n, n)); 7059566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Grad)); 7069566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->Grad, 3, NULL, 3, NULL)); 7079566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->Grad, 3, NULL)); 7089566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->Grad, &istart, &iend)); 709c4762a1bSJed Brown 710c4762a1bSJed Brown for (i = istart; i < iend; i++) { 711c4762a1bSJed Brown if (i < n) { 712c4762a1bSJed Brown iblock = i / user->mx; 713c4762a1bSJed Brown j = iblock * user->mx + ((i + user->mx - 1) % user->mx); 7149566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &half_hinv, INSERT_VALUES)); 715c4762a1bSJed Brown j = iblock * user->mx + ((i + 1) % user->mx); 7169566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &neg_half_hinv, INSERT_VALUES)); 717c4762a1bSJed Brown } 718c4762a1bSJed Brown if (i >= n) { 719c4762a1bSJed Brown j = (i - user->mx) % n; 7209566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &half_hinv, INSERT_VALUES)); 721c4762a1bSJed Brown j = (j + 2 * user->mx) % n; 7229566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &neg_half_hinv, INSERT_VALUES)); 723c4762a1bSJed Brown } 724c4762a1bSJed Brown } 725c4762a1bSJed Brown 7269566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Grad, MAT_FINAL_ASSEMBLY)); 7279566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Grad, MAT_FINAL_ASSEMBLY)); 728c4762a1bSJed Brown 7299566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Gradxy[0])); 7309566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Gradxy[0], PETSC_DECIDE, PETSC_DECIDE, n, n)); 7319566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Gradxy[0])); 7329566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->Gradxy[0], 3, NULL, 3, NULL)); 7339566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->Gradxy[0], 3, NULL)); 7349566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->Gradxy[0], &istart, &iend)); 735c4762a1bSJed Brown 736c4762a1bSJed Brown for (i = istart; i < iend; i++) { 737c4762a1bSJed Brown iblock = i / user->mx; 738c4762a1bSJed Brown j = iblock * user->mx + ((i + user->mx - 1) % user->mx); 7399566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Gradxy[0], 1, &i, 1, &j, &half_hinv, INSERT_VALUES)); 740c4762a1bSJed Brown j = iblock * user->mx + ((i + 1) % user->mx); 7419566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Gradxy[0], 1, &i, 1, &j, &neg_half_hinv, INSERT_VALUES)); 7429566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Gradxy[0], 1, &i, 1, &i, &zero, INSERT_VALUES)); 743c4762a1bSJed Brown } 7449566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Gradxy[0], MAT_FINAL_ASSEMBLY)); 7459566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Gradxy[0], MAT_FINAL_ASSEMBLY)); 746c4762a1bSJed Brown 7479566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Gradxy[1])); 7489566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Gradxy[1], PETSC_DECIDE, PETSC_DECIDE, n, n)); 7499566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Gradxy[1])); 7509566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->Gradxy[1], 3, NULL, 3, NULL)); 7519566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->Gradxy[1], 3, NULL)); 7529566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->Gradxy[1], &istart, &iend)); 753c4762a1bSJed Brown 754c4762a1bSJed Brown for (i = istart; i < iend; i++) { 755c4762a1bSJed Brown j = (i + n - user->mx) % n; 7569566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Gradxy[1], 1, &i, 1, &j, &half_hinv, INSERT_VALUES)); 757c4762a1bSJed Brown j = (j + 2 * user->mx) % n; 7589566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Gradxy[1], 1, &i, 1, &j, &neg_half_hinv, INSERT_VALUES)); 7599566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Gradxy[1], 1, &i, 1, &i, &zero, INSERT_VALUES)); 760c4762a1bSJed Brown } 7619566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Gradxy[1], MAT_FINAL_ASSEMBLY)); 7629566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Gradxy[1], MAT_FINAL_ASSEMBLY)); 763c4762a1bSJed Brown 764c4762a1bSJed Brown /* Generate Div matrix */ 7659566063dSJacob Faibussowitsch PetscCall(MatTranspose(user->Grad, MAT_INITIAL_MATRIX, &user->Div)); 7669566063dSJacob Faibussowitsch PetscCall(MatTranspose(user->Gradxy[0], MAT_INITIAL_MATRIX, &user->Divxy[0])); 7679566063dSJacob Faibussowitsch PetscCall(MatTranspose(user->Gradxy[1], MAT_INITIAL_MATRIX, &user->Divxy[1])); 768c4762a1bSJed Brown 769c4762a1bSJed Brown /* Off-diagonal averaging matrix */ 7709566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->M)); 7719566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->M, PETSC_DECIDE, PETSC_DECIDE, n, n)); 7729566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->M)); 7739566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->M, 4, NULL, 4, NULL)); 7749566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->M, 4, NULL)); 7759566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->M, &istart, &iend)); 776c4762a1bSJed Brown 777c4762a1bSJed Brown for (i = istart; i < iend; i++) { 778c4762a1bSJed Brown /* kron(Id,Av) */ 779c4762a1bSJed Brown iblock = i / user->mx; 780c4762a1bSJed Brown j = iblock * user->mx + ((i + user->mx - 1) % user->mx); 7819566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->M, 1, &i, 1, &j, &quarter, INSERT_VALUES)); 782c4762a1bSJed Brown j = iblock * user->mx + ((i + 1) % user->mx); 7839566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->M, 1, &i, 1, &j, &quarter, INSERT_VALUES)); 784c4762a1bSJed Brown 785c4762a1bSJed Brown /* kron(Av,Id) */ 786c4762a1bSJed Brown j = (i + user->mx) % n; 7879566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->M, 1, &i, 1, &j, &quarter, INSERT_VALUES)); 788c4762a1bSJed Brown j = (i + n - user->mx) % n; 7899566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->M, 1, &i, 1, &j, &quarter, INSERT_VALUES)); 790c4762a1bSJed Brown } 7919566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->M, MAT_FINAL_ASSEMBLY)); 7929566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->M, MAT_FINAL_ASSEMBLY)); 793c4762a1bSJed Brown 794c4762a1bSJed Brown /* Generate 2D grid */ 7959566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &XX)); 7969566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->q)); 7979566063dSJacob Faibussowitsch PetscCall(VecSetSizes(XX, PETSC_DECIDE, n)); 7989566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->q, PETSC_DECIDE, n * user->nt)); 7999566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(XX)); 8009566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->q)); 801c4762a1bSJed Brown 8029566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &YY)); 8039566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &XXwork)); 8049566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &YYwork)); 8059566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &user->d)); 8069566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &user->dwork)); 807c4762a1bSJed Brown 8089566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(XX, &istart, &iend)); 809c4762a1bSJed Brown for (linear_index = istart; linear_index < iend; linear_index++) { 810c4762a1bSJed Brown i = linear_index % user->mx; 811c4762a1bSJed Brown j = (linear_index - i) / user->mx; 812c4762a1bSJed Brown vx = h * (i + 0.5); 813c4762a1bSJed Brown vy = h * (j + 0.5); 8149566063dSJacob Faibussowitsch PetscCall(VecSetValues(XX, 1, &linear_index, &vx, INSERT_VALUES)); 8159566063dSJacob Faibussowitsch PetscCall(VecSetValues(YY, 1, &linear_index, &vy, INSERT_VALUES)); 816c4762a1bSJed Brown } 817c4762a1bSJed Brown 8189566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(XX)); 8199566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(XX)); 8209566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(YY)); 8219566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(YY)); 822c4762a1bSJed Brown 823c4762a1bSJed Brown /* Compute final density function yT 824c4762a1bSJed Brown yT = 1.0 + exp(-30*((x-0.25)^2+(y-0.25)^2)) + exp(-30*((x-0.75)^2+(y-0.75)^2)) 825c4762a1bSJed Brown yT = yT / (h^2*sum(yT)) */ 8269566063dSJacob Faibussowitsch PetscCall(VecCopy(XX, XXwork)); 8279566063dSJacob Faibussowitsch PetscCall(VecCopy(YY, YYwork)); 828c4762a1bSJed Brown 8299566063dSJacob Faibussowitsch PetscCall(VecShift(XXwork, -0.25)); 8309566063dSJacob Faibussowitsch PetscCall(VecShift(YYwork, -0.25)); 831c4762a1bSJed Brown 8329566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(XXwork, XXwork, XXwork)); 8339566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(YYwork, YYwork, YYwork)); 834c4762a1bSJed Brown 8359566063dSJacob Faibussowitsch PetscCall(VecCopy(XXwork, user->dwork)); 8369566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->dwork, 1.0, YYwork)); 8379566063dSJacob Faibussowitsch PetscCall(VecScale(user->dwork, -30.0)); 8389566063dSJacob Faibussowitsch PetscCall(VecExp(user->dwork)); 8399566063dSJacob Faibussowitsch PetscCall(VecCopy(user->dwork, user->d)); 840c4762a1bSJed Brown 8419566063dSJacob Faibussowitsch PetscCall(VecCopy(XX, XXwork)); 8429566063dSJacob Faibussowitsch PetscCall(VecCopy(YY, YYwork)); 843c4762a1bSJed Brown 8449566063dSJacob Faibussowitsch PetscCall(VecShift(XXwork, -0.75)); 8459566063dSJacob Faibussowitsch PetscCall(VecShift(YYwork, -0.75)); 846c4762a1bSJed Brown 8479566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(XXwork, XXwork, XXwork)); 8489566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(YYwork, YYwork, YYwork)); 849c4762a1bSJed Brown 8509566063dSJacob Faibussowitsch PetscCall(VecCopy(XXwork, user->dwork)); 8519566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->dwork, 1.0, YYwork)); 8529566063dSJacob Faibussowitsch PetscCall(VecScale(user->dwork, -30.0)); 8539566063dSJacob Faibussowitsch PetscCall(VecExp(user->dwork)); 854c4762a1bSJed Brown 8559566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->d, 1.0, user->dwork)); 8569566063dSJacob Faibussowitsch PetscCall(VecShift(user->d, 1.0)); 8579566063dSJacob Faibussowitsch PetscCall(VecSum(user->d, &sum)); 8589566063dSJacob Faibussowitsch PetscCall(VecScale(user->d, 1.0 / (h * h * sum))); 859c4762a1bSJed Brown 860c4762a1bSJed Brown /* Initial conditions of forward problem */ 8619566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &bc)); 8629566063dSJacob Faibussowitsch PetscCall(VecCopy(XX, XXwork)); 8639566063dSJacob Faibussowitsch PetscCall(VecCopy(YY, YYwork)); 864c4762a1bSJed Brown 8659566063dSJacob Faibussowitsch PetscCall(VecShift(XXwork, -0.5)); 8669566063dSJacob Faibussowitsch PetscCall(VecShift(YYwork, -0.5)); 867c4762a1bSJed Brown 8689566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(XXwork, XXwork, XXwork)); 8699566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(YYwork, YYwork, YYwork)); 870c4762a1bSJed Brown 8719566063dSJacob Faibussowitsch PetscCall(VecWAXPY(bc, 1.0, XXwork, YYwork)); 8729566063dSJacob Faibussowitsch PetscCall(VecScale(bc, -50.0)); 8739566063dSJacob Faibussowitsch PetscCall(VecExp(bc)); 8749566063dSJacob Faibussowitsch PetscCall(VecShift(bc, 1.0)); 8759566063dSJacob Faibussowitsch PetscCall(VecSum(bc, &sum)); 8769566063dSJacob Faibussowitsch PetscCall(VecScale(bc, 1.0 / (h * h * sum))); 877c4762a1bSJed Brown 878c4762a1bSJed Brown /* Create scatter from y to y_1,y_2,...,y_nt */ 879c4762a1bSJed Brown /* TODO: Reorder for better parallelism. (This will require reordering Q and L as well.) */ 8809566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->nt * user->mx * user->mx, &user->yi_scatter)); 8819566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &yi)); 8829566063dSJacob Faibussowitsch PetscCall(VecSetSizes(yi, PETSC_DECIDE, user->mx * user->mx)); 8839566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(yi)); 8849566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(yi, user->nt, &user->yi)); 8859566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(yi, user->nt, &user->yiwork)); 8869566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(yi, user->nt, &user->ziwork)); 887c4762a1bSJed Brown for (i = 0; i < user->nt; i++) { 8889566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->yi[i], &lo, &hi)); 8899566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_to_yi)); 8909566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo + i * user->mx * user->mx, 1, &is_from_y)); 8919566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(user->y, is_from_y, user->yi[i], is_to_yi, &user->yi_scatter[i])); 8929566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_to_yi)); 8939566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_from_y)); 894c4762a1bSJed Brown } 895c4762a1bSJed Brown 896c4762a1bSJed Brown /* Create scatter from u to ux_1,uy_1,ux_2,uy_2,...,ux_nt,uy_nt */ 897c4762a1bSJed Brown /* TODO: reorder for better parallelism */ 8989566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->nt * user->mx * user->mx, &user->uxi_scatter)); 8999566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->nt * user->mx * user->mx, &user->uyi_scatter)); 9009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->nt * user->mx * user->mx, &user->ux_scatter)); 9019566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->nt * user->mx * user->mx, &user->uy_scatter)); 9029566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2 * user->nt * user->mx * user->mx, &user->ui_scatter)); 9039566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &uxi)); 9049566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &ui)); 9059566063dSJacob Faibussowitsch PetscCall(VecSetSizes(uxi, PETSC_DECIDE, user->mx * user->mx)); 9069566063dSJacob Faibussowitsch PetscCall(VecSetSizes(ui, PETSC_DECIDE, 2 * user->mx * user->mx)); 9079566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(uxi)); 9089566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(ui)); 9099566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(uxi, user->nt, &user->uxi)); 9109566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(uxi, user->nt, &user->uyi)); 9119566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(uxi, user->nt, &user->uxiwork)); 9129566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(uxi, user->nt, &user->uyiwork)); 9139566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ui, user->nt, &user->ui)); 9149566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ui, user->nt, &user->uiwork)); 915c4762a1bSJed Brown for (i = 0; i < user->nt; i++) { 9169566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->uxi[i], &lo, &hi)); 9179566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_to_uxi)); 9189566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo + 2 * i * user->mx * user->mx, 1, &is_from_u)); 9199566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(user->u, is_from_u, user->uxi[i], is_to_uxi, &user->uxi_scatter[i])); 920c4762a1bSJed Brown 9219566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_to_uxi)); 9229566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_from_u)); 923c4762a1bSJed Brown 9249566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->uyi[i], &lo, &hi)); 9259566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_to_uyi)); 9269566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo + (2 * i + 1) * user->mx * user->mx, 1, &is_from_u)); 9279566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(user->u, is_from_u, user->uyi[i], is_to_uyi, &user->uyi_scatter[i])); 928c4762a1bSJed Brown 9299566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_to_uyi)); 9309566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_from_u)); 931c4762a1bSJed Brown 9329566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->uxi[i], &lo, &hi)); 9339566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_to_uxi)); 9349566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_from_u)); 9359566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(user->ui[i], is_from_u, user->uxi[i], is_to_uxi, &user->ux_scatter[i])); 936c4762a1bSJed Brown 9379566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_to_uxi)); 9389566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_from_u)); 939c4762a1bSJed Brown 9409566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->uyi[i], &lo, &hi)); 9419566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_to_uyi)); 9429566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo + user->mx * user->mx, 1, &is_from_u)); 9439566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(user->ui[i], is_from_u, user->uyi[i], is_to_uyi, &user->uy_scatter[i])); 944c4762a1bSJed Brown 9459566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_to_uyi)); 9469566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_from_u)); 947c4762a1bSJed Brown 9489566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->ui[i], &lo, &hi)); 9499566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_to_uxi)); 9509566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo + 2 * i * user->mx * user->mx, 1, &is_from_u)); 9519566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(user->u, is_from_u, user->ui[i], is_to_uxi, &user->ui_scatter[i])); 952c4762a1bSJed Brown 9539566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_to_uxi)); 9549566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_from_u)); 955c4762a1bSJed Brown } 956c4762a1bSJed Brown 957c4762a1bSJed Brown /* RHS of forward problem */ 9589566063dSJacob Faibussowitsch PetscCall(MatMult(user->M, bc, user->yiwork[0])); 95948a46eb9SPierre Jolivet for (i = 1; i < user->nt; i++) PetscCall(VecSet(user->yiwork[i], 0.0)); 9609566063dSJacob Faibussowitsch PetscCall(Gather_yi(user->q, user->yiwork, user->yi_scatter, user->nt)); 961c4762a1bSJed Brown 962c4762a1bSJed Brown /* Compute true velocity field utrue */ 9639566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->u, &user->utrue)); 964c4762a1bSJed Brown for (i = 0; i < user->nt; i++) { 9659566063dSJacob Faibussowitsch PetscCall(VecCopy(YY, user->uxi[i])); 9669566063dSJacob Faibussowitsch PetscCall(VecScale(user->uxi[i], 150.0 * i * user->ht)); 9679566063dSJacob Faibussowitsch PetscCall(VecCopy(XX, user->uyi[i])); 9689566063dSJacob Faibussowitsch PetscCall(VecShift(user->uyi[i], -10.0)); 9699566063dSJacob Faibussowitsch PetscCall(VecScale(user->uyi[i], 15.0 * i * user->ht)); 970c4762a1bSJed Brown } 9719566063dSJacob Faibussowitsch PetscCall(Gather_uxi_uyi(user->utrue, user->uxi, user->uxi_scatter, user->uyi, user->uyi_scatter, user->nt)); 972c4762a1bSJed Brown 973c4762a1bSJed Brown /* Initial guess and reference model */ 9749566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->utrue, &user->ur)); 975c4762a1bSJed Brown for (i = 0; i < user->nt; i++) { 9769566063dSJacob Faibussowitsch PetscCall(VecCopy(XX, user->uxi[i])); 9779566063dSJacob Faibussowitsch PetscCall(VecShift(user->uxi[i], i * user->ht)); 9789566063dSJacob Faibussowitsch PetscCall(VecCopy(YY, user->uyi[i])); 9799566063dSJacob Faibussowitsch PetscCall(VecShift(user->uyi[i], -i * user->ht)); 980c4762a1bSJed Brown } 9819566063dSJacob Faibussowitsch PetscCall(Gather_uxi_uyi(user->ur, user->uxi, user->uxi_scatter, user->uyi, user->uyi_scatter, user->nt)); 982c4762a1bSJed Brown 983c4762a1bSJed Brown /* Generate regularization matrix L */ 9849566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->LT)); 9859566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->LT, PETSC_DECIDE, PETSC_DECIDE, 2 * n * user->nt, n * user->nt)); 9869566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->LT)); 9879566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->LT, 1, NULL, 1, NULL)); 9889566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->LT, 1, NULL)); 9899566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->LT, &istart, &iend)); 990c4762a1bSJed Brown 991c4762a1bSJed Brown for (i = istart; i < iend; i++) { 992c4762a1bSJed Brown iblock = (i + n) / (2 * n); 993c4762a1bSJed Brown j = i - iblock * n; 9949566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->LT, 1, &i, 1, &j, &user->gamma, INSERT_VALUES)); 995c4762a1bSJed Brown } 996c4762a1bSJed Brown 9979566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->LT, MAT_FINAL_ASSEMBLY)); 9989566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->LT, MAT_FINAL_ASSEMBLY)); 999c4762a1bSJed Brown 10009566063dSJacob Faibussowitsch PetscCall(MatTranspose(user->LT, MAT_INITIAL_MATRIX, &user->L)); 1001c4762a1bSJed Brown 1002c4762a1bSJed Brown /* Build work vectors and matrices */ 10039566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->lwork)); 10049566063dSJacob Faibussowitsch PetscCall(VecSetType(user->lwork, VECMPI)); 10059566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->lwork, PETSC_DECIDE, user->m)); 10069566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->lwork)); 1007c4762a1bSJed Brown 10089566063dSJacob Faibussowitsch PetscCall(MatDuplicate(user->Div, MAT_SHARE_NONZERO_PATTERN, &user->Divwork)); 1009c4762a1bSJed Brown 10109566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->y, &user->ywork)); 10119566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->u, &user->uwork)); 10129566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->u, &user->vwork)); 10139566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->u, &user->js_diag)); 10149566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->c, &user->cwork)); 1015c4762a1bSJed Brown 1016c4762a1bSJed Brown /* Create matrix-free shell user->Js for computing A*x */ 10179566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m, user->m, user, &user->Js)); 10189566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->Js, MATOP_MULT, (void (*)(void))StateMatMult)); 10199566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->Js, MATOP_DUPLICATE, (void (*)(void))StateMatDuplicate)); 10209566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->Js, MATOP_MULT_TRANSPOSE, (void (*)(void))StateMatMultTranspose)); 10219566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->Js, MATOP_GET_DIAGONAL, (void (*)(void))StateMatGetDiagonal)); 1022c4762a1bSJed Brown 1023c4762a1bSJed Brown /* Diagonal blocks of user->Js */ 10249566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, n, n, user, &user->JsBlock)); 10259566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->JsBlock, MATOP_MULT, (void (*)(void))StateMatBlockMult)); 10269566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->JsBlock, MATOP_MULT_TRANSPOSE, (void (*)(void))StateMatBlockMultTranspose)); 1027c4762a1bSJed Brown 1028c4762a1bSJed Brown /* Create a matrix-free shell user->JsBlockPrec for computing (U+D)\D*(L+D)\x, where JsBlock = L+D+U, 1029c4762a1bSJed Brown D is diagonal, L is strictly lower triangular, and U is strictly upper triangular. 1030c4762a1bSJed Brown This is an SOR preconditioner for user->JsBlock. */ 10319566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, n, n, user, &user->JsBlockPrec)); 10329566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->JsBlockPrec, MATOP_MULT, (void (*)(void))StateMatBlockPrecMult)); 10339566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->JsBlockPrec, MATOP_MULT_TRANSPOSE, (void (*)(void))StateMatBlockPrecMultTranspose)); 1034c4762a1bSJed Brown 1035c4762a1bSJed Brown /* Create a matrix-free shell user->Jd for computing B*x */ 10369566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m, user->n - user->m, user, &user->Jd)); 10379566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->Jd, MATOP_MULT, (void (*)(void))DesignMatMult)); 10389566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->Jd, MATOP_MULT_TRANSPOSE, (void (*)(void))DesignMatMultTranspose)); 1039c4762a1bSJed Brown 1040c4762a1bSJed Brown /* User-defined routines for computing user->Js\x and user->Js^T\x*/ 10419566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m, user->m, user, &user->JsInv)); 10429566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->JsInv, MATOP_MULT, (void (*)(void))StateMatInvMult)); 10439566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->JsInv, MATOP_MULT_TRANSPOSE, (void (*)(void))StateMatInvTransposeMult)); 1044c4762a1bSJed Brown 1045c4762a1bSJed Brown /* Build matrices for SOR preconditioner */ 10469566063dSJacob Faibussowitsch PetscCall(Scatter_uxi_uyi(user->u, user->uxi, user->uxi_scatter, user->uyi, user->uyi_scatter, user->nt)); 10479566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(5 * n, &user->C)); 10489566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2 * n, &user->Cwork)); 1049c4762a1bSJed Brown for (i = 0; i < user->nt; i++) { 10509566063dSJacob Faibussowitsch PetscCall(MatDuplicate(user->Divxy[0], MAT_COPY_VALUES, &user->C[i])); 10519566063dSJacob Faibussowitsch PetscCall(MatDuplicate(user->Divxy[1], MAT_COPY_VALUES, &user->Cwork[i])); 1052c4762a1bSJed Brown 10539566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(user->C[i], NULL, user->uxi[i])); 10549566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(user->Cwork[i], NULL, user->uyi[i])); 10559566063dSJacob Faibussowitsch PetscCall(MatAXPY(user->C[i], 1.0, user->Cwork[i], DIFFERENT_NONZERO_PATTERN)); 10569566063dSJacob Faibussowitsch PetscCall(MatScale(user->C[i], user->ht)); 10579566063dSJacob Faibussowitsch PetscCall(MatShift(user->C[i], 1.0)); 1058c4762a1bSJed Brown } 1059c4762a1bSJed Brown 1060c4762a1bSJed Brown /* Solver options and tolerances */ 10619566063dSJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_WORLD, &user->solver)); 10629566063dSJacob Faibussowitsch PetscCall(KSPSetType(user->solver, KSPGMRES)); 10639566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(user->solver, user->JsBlock, user->JsBlockPrec)); 10649566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(user->solver, 1e-4, 1e-20, 1e3, 500)); 10659566063dSJacob Faibussowitsch /* PetscCall(KSPSetTolerances(user->solver,1e-8,1e-16,1e3,500)); */ 10669566063dSJacob Faibussowitsch PetscCall(KSPGetPC(user->solver, &user->prec)); 10679566063dSJacob Faibussowitsch PetscCall(PCSetType(user->prec, PCSHELL)); 1068c4762a1bSJed Brown 10699566063dSJacob Faibussowitsch PetscCall(PCShellSetApply(user->prec, StateMatBlockPrecMult)); 10709566063dSJacob Faibussowitsch PetscCall(PCShellSetApplyTranspose(user->prec, StateMatBlockPrecMultTranspose)); 10719566063dSJacob Faibussowitsch PetscCall(PCShellSetContext(user->prec, user)); 1072c4762a1bSJed Brown 1073c4762a1bSJed Brown /* Compute true state function yt given ut */ 10749566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->ytrue)); 10759566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->ytrue, PETSC_DECIDE, n * user->nt)); 10769566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->ytrue)); 1077c4762a1bSJed Brown user->c_formed = PETSC_TRUE; 10789566063dSJacob Faibussowitsch PetscCall(VecCopy(user->utrue, user->u)); /* Set u=utrue temporarily for StateMatInv */ 10799566063dSJacob Faibussowitsch PetscCall(VecSet(user->ytrue, 0.0)); /* Initial guess */ 10809566063dSJacob Faibussowitsch PetscCall(StateMatInvMult(user->Js, user->q, user->ytrue)); 10819566063dSJacob Faibussowitsch PetscCall(VecCopy(user->ur, user->u)); /* Reset u=ur */ 1082c4762a1bSJed Brown 1083c4762a1bSJed Brown /* Initial guess y0 for state given u0 */ 10849566063dSJacob Faibussowitsch PetscCall(StateMatInvMult(user->Js, user->q, user->y)); 1085c4762a1bSJed Brown 1086c4762a1bSJed Brown /* Data discretization */ 10879566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Q)); 10889566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Q, PETSC_DECIDE, PETSC_DECIDE, user->mx * user->mx, user->m)); 10899566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Q)); 10909566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->Q, 0, NULL, 1, NULL)); 10919566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->Q, 1, NULL)); 1092c4762a1bSJed Brown 10939566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->Q, &istart, &iend)); 1094c4762a1bSJed Brown 1095c4762a1bSJed Brown for (i = istart; i < iend; i++) { 1096c4762a1bSJed Brown j = i + user->m - user->mx * user->mx; 10979566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Q, 1, &i, 1, &j, &one, INSERT_VALUES)); 1098c4762a1bSJed Brown } 1099c4762a1bSJed Brown 11009566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Q, MAT_FINAL_ASSEMBLY)); 11019566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Q, MAT_FINAL_ASSEMBLY)); 1102c4762a1bSJed Brown 11039566063dSJacob Faibussowitsch PetscCall(MatTranspose(user->Q, MAT_INITIAL_MATRIX, &user->QT)); 1104c4762a1bSJed Brown 11059566063dSJacob Faibussowitsch PetscCall(VecDestroy(&XX)); 11069566063dSJacob Faibussowitsch PetscCall(VecDestroy(&YY)); 11079566063dSJacob Faibussowitsch PetscCall(VecDestroy(&XXwork)); 11089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&YYwork)); 11099566063dSJacob Faibussowitsch PetscCall(VecDestroy(&yi)); 11109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&uxi)); 11119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ui)); 11129566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bc)); 1113c4762a1bSJed Brown 1114c4762a1bSJed Brown /* Now that initial conditions have been set, let the user pass tolerance options to the KSP solver */ 11159566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(user->solver)); 1116c4762a1bSJed Brown PetscFunctionReturn(0); 1117c4762a1bSJed Brown } 1118c4762a1bSJed Brown 1119*d71ae5a4SJacob Faibussowitsch PetscErrorCode HyperbolicDestroy(AppCtx *user) 1120*d71ae5a4SJacob Faibussowitsch { 1121c4762a1bSJed Brown PetscInt i; 1122c4762a1bSJed Brown 1123c4762a1bSJed Brown PetscFunctionBegin; 11249566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Q)); 11259566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->QT)); 11269566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Div)); 11279566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Divwork)); 11289566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Grad)); 11299566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->L)); 11309566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->LT)); 11319566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&user->solver)); 11329566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Js)); 11339566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Jd)); 11349566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->JsBlockPrec)); 11359566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->JsInv)); 11369566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->JsBlock)); 11379566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Divxy[0])); 11389566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Divxy[1])); 11399566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Gradxy[0])); 11409566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Gradxy[1])); 11419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->M)); 1142c4762a1bSJed Brown for (i = 0; i < user->nt; i++) { 11439566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->C[i])); 11449566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Cwork[i])); 1145c4762a1bSJed Brown } 11469566063dSJacob Faibussowitsch PetscCall(PetscFree(user->C)); 11479566063dSJacob Faibussowitsch PetscCall(PetscFree(user->Cwork)); 11489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->u)); 11499566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->uwork)); 11509566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->vwork)); 11519566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->utrue)); 11529566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->y)); 11539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->ywork)); 11549566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->ytrue)); 11559566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->nt, &user->yi)); 11569566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->nt, &user->yiwork)); 11579566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->nt, &user->ziwork)); 11589566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->nt, &user->uxi)); 11599566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->nt, &user->uyi)); 11609566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->nt, &user->uxiwork)); 11619566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->nt, &user->uyiwork)); 11629566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->nt, &user->ui)); 11639566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->nt, &user->uiwork)); 11649566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->c)); 11659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->cwork)); 11669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->ur)); 11679566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->q)); 11689566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->d)); 11699566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->dwork)); 11709566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->lwork)); 11719566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->js_diag)); 11729566063dSJacob Faibussowitsch PetscCall(ISDestroy(&user->s_is)); 11739566063dSJacob Faibussowitsch PetscCall(ISDestroy(&user->d_is)); 11749566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->state_scatter)); 11759566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->design_scatter)); 1176c4762a1bSJed Brown for (i = 0; i < user->nt; i++) { 11779566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->uxi_scatter[i])); 11789566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->uyi_scatter[i])); 11799566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->ux_scatter[i])); 11809566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->uy_scatter[i])); 11819566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->ui_scatter[i])); 11829566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->yi_scatter[i])); 1183c4762a1bSJed Brown } 11849566063dSJacob Faibussowitsch PetscCall(PetscFree(user->uxi_scatter)); 11859566063dSJacob Faibussowitsch PetscCall(PetscFree(user->uyi_scatter)); 11869566063dSJacob Faibussowitsch PetscCall(PetscFree(user->ux_scatter)); 11879566063dSJacob Faibussowitsch PetscCall(PetscFree(user->uy_scatter)); 11889566063dSJacob Faibussowitsch PetscCall(PetscFree(user->ui_scatter)); 11899566063dSJacob Faibussowitsch PetscCall(PetscFree(user->yi_scatter)); 1190c4762a1bSJed Brown PetscFunctionReturn(0); 1191c4762a1bSJed Brown } 1192c4762a1bSJed Brown 1193*d71ae5a4SJacob Faibussowitsch PetscErrorCode HyperbolicMonitor(Tao tao, void *ptr) 1194*d71ae5a4SJacob Faibussowitsch { 1195c4762a1bSJed Brown Vec X; 1196c4762a1bSJed Brown PetscReal unorm, ynorm; 1197c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 1198c4762a1bSJed Brown 1199c4762a1bSJed Brown PetscFunctionBegin; 12009566063dSJacob Faibussowitsch PetscCall(TaoGetSolution(tao, &X)); 12019566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->ywork, user->state_scatter, user->uwork, user->design_scatter)); 12029566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->ywork, -1.0, user->ytrue)); 12039566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->uwork, -1.0, user->utrue)); 12049566063dSJacob Faibussowitsch PetscCall(VecNorm(user->uwork, NORM_2, &unorm)); 12059566063dSJacob Faibussowitsch PetscCall(VecNorm(user->ywork, NORM_2, &ynorm)); 12069566063dSJacob Faibussowitsch PetscCall(PetscPrintf(MPI_COMM_WORLD, "||u-ut||=%g ||y-yt||=%g\n", (double)unorm, (double)ynorm)); 1207c4762a1bSJed Brown PetscFunctionReturn(0); 1208c4762a1bSJed Brown } 1209c4762a1bSJed Brown 1210c4762a1bSJed Brown /*TEST 1211c4762a1bSJed Brown 1212c4762a1bSJed Brown build: 1213c4762a1bSJed Brown requires: !complex 1214c4762a1bSJed Brown 1215c4762a1bSJed Brown test: 1216c4762a1bSJed Brown requires: !single 1217c4762a1bSJed Brown args: -tao_cmonitor -tao_max_funcs 10 -tao_type lcl -tao_gatol 1.e-5 1218c4762a1bSJed Brown 1219c4762a1bSJed Brown test: 1220c4762a1bSJed Brown suffix: guess_pod 1221c4762a1bSJed Brown requires: !single 1222c4762a1bSJed Brown args: -tao_cmonitor -tao_max_funcs 10 -tao_type lcl -ksp_guess_type pod -tao_gatol 1.e-5 1223c4762a1bSJed Brown 1224c4762a1bSJed Brown TEST*/ 1225