1c4762a1bSJed Brown #include <petsc/private/taoimpl.h> 2c4762a1bSJed Brown 3c4762a1bSJed Brown typedef struct { 4c4762a1bSJed Brown PetscInt n; /* Number of total variables */ 5c4762a1bSJed Brown PetscInt m; /* Number of constraints */ 6c4762a1bSJed Brown PetscInt nstate; 7c4762a1bSJed Brown PetscInt ndesign; 8c4762a1bSJed Brown PetscInt mx; /* grid points in each direction */ 9c4762a1bSJed Brown PetscInt ns; /* Number of data samples (1<=ns<=8) 10c4762a1bSJed Brown Currently only ns=1 is supported */ 11c4762a1bSJed Brown PetscInt ndata; /* Number of data points per sample */ 12c4762a1bSJed Brown IS s_is; 13c4762a1bSJed Brown IS d_is; 14c4762a1bSJed Brown 15c4762a1bSJed Brown VecScatter state_scatter; 16c4762a1bSJed Brown VecScatter design_scatter; 17c4762a1bSJed Brown VecScatter *yi_scatter, *di_scatter; 18c4762a1bSJed Brown Vec suby, subq, subd; 19c4762a1bSJed Brown Mat Js, Jd, JsPrec, JsInv, JsBlock; 20c4762a1bSJed Brown 21c4762a1bSJed Brown PetscReal alpha; /* Regularization parameter */ 22c4762a1bSJed Brown PetscReal beta; /* Weight attributed to ||u||^2 in regularization functional */ 23c4762a1bSJed Brown PetscReal noise; /* Amount of noise to add to data */ 24c4762a1bSJed Brown PetscReal *ones; 25c4762a1bSJed Brown Mat Q; 26c4762a1bSJed Brown Mat MQ; 27c4762a1bSJed Brown Mat L; 28c4762a1bSJed Brown 29c4762a1bSJed Brown Mat Grad; 30c4762a1bSJed Brown Mat Av, Avwork; 31c4762a1bSJed Brown Mat Div, Divwork; 32c4762a1bSJed Brown Mat DSG; 33c4762a1bSJed Brown Mat Diag, Ones; 34c4762a1bSJed Brown 35c4762a1bSJed Brown Vec q; 36c4762a1bSJed Brown Vec ur; /* reference */ 37c4762a1bSJed Brown 38c4762a1bSJed Brown Vec d; 39c4762a1bSJed Brown Vec dwork; 40c4762a1bSJed Brown 41c4762a1bSJed Brown Vec x; /* super vec of y,u */ 42c4762a1bSJed Brown 43c4762a1bSJed Brown Vec y; /* state variables */ 44c4762a1bSJed Brown Vec ywork; 45c4762a1bSJed Brown 46c4762a1bSJed Brown Vec ytrue; 47c4762a1bSJed Brown 48c4762a1bSJed Brown Vec u; /* design variables */ 49c4762a1bSJed Brown Vec uwork; 50c4762a1bSJed Brown 51c4762a1bSJed Brown Vec utrue; 52c4762a1bSJed Brown 53c4762a1bSJed Brown Vec js_diag; 54c4762a1bSJed Brown 55c4762a1bSJed Brown Vec c; /* constraint vector */ 56c4762a1bSJed Brown Vec cwork; 57c4762a1bSJed Brown 58c4762a1bSJed Brown Vec lwork; 59c4762a1bSJed Brown Vec S; 60c4762a1bSJed Brown Vec Swork, Twork, Sdiag, Ywork; 61c4762a1bSJed Brown Vec Av_u; 62c4762a1bSJed Brown 63c4762a1bSJed Brown KSP solver; 64c4762a1bSJed Brown PC prec; 65c4762a1bSJed Brown 66c4762a1bSJed Brown PetscReal tola, tolb, tolc, told; 67c4762a1bSJed Brown PetscInt ksp_its; 68c4762a1bSJed Brown PetscInt ksp_its_initial; 69c4762a1bSJed Brown PetscLogStage stages[10]; 70c4762a1bSJed Brown PetscBool use_ptap; 71c4762a1bSJed Brown PetscBool use_lrc; 72c4762a1bSJed Brown } AppCtx; 73c4762a1bSJed Brown 74c4762a1bSJed Brown PetscErrorCode FormFunction(Tao, Vec, PetscReal *, void *); 75c4762a1bSJed Brown PetscErrorCode FormGradient(Tao, Vec, Vec, void *); 76c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *); 77c4762a1bSJed Brown PetscErrorCode FormJacobianState(Tao, Vec, Mat, Mat, Mat, void *); 78c4762a1bSJed Brown PetscErrorCode FormJacobianDesign(Tao, Vec, Mat, void *); 79c4762a1bSJed Brown PetscErrorCode FormConstraints(Tao, Vec, Vec, void *); 80c4762a1bSJed Brown PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void *); 81c4762a1bSJed Brown PetscErrorCode Gather(Vec, Vec, VecScatter, Vec, VecScatter); 82c4762a1bSJed Brown PetscErrorCode Scatter(Vec, Vec, VecScatter, Vec, VecScatter); 83c4762a1bSJed Brown PetscErrorCode EllipticInitialize(AppCtx *); 84c4762a1bSJed Brown PetscErrorCode EllipticDestroy(AppCtx *); 85c4762a1bSJed Brown PetscErrorCode EllipticMonitor(Tao, void *); 86c4762a1bSJed Brown 87c4762a1bSJed Brown PetscErrorCode StateBlockMatMult(Mat, Vec, Vec); 88c4762a1bSJed Brown PetscErrorCode StateMatMult(Mat, Vec, Vec); 89c4762a1bSJed Brown 90c4762a1bSJed Brown PetscErrorCode StateInvMatMult(Mat, Vec, Vec); 91c4762a1bSJed Brown PetscErrorCode DesignMatMult(Mat, Vec, Vec); 92c4762a1bSJed Brown PetscErrorCode DesignMatMultTranspose(Mat, Vec, Vec); 93c4762a1bSJed Brown 94c4762a1bSJed Brown PetscErrorCode QMatMult(Mat, Vec, Vec); 95c4762a1bSJed Brown PetscErrorCode QMatMultTranspose(Mat, Vec, Vec); 96c4762a1bSJed Brown 97c4762a1bSJed Brown static char help[] = ""; 98c4762a1bSJed Brown 99*9371c9d4SSatish Balay int main(int argc, char **argv) { 100c4762a1bSJed Brown Vec x0; 101c4762a1bSJed Brown Tao tao; 102c4762a1bSJed Brown AppCtx user; 103c4762a1bSJed Brown PetscInt ntests = 1; 104c4762a1bSJed Brown PetscInt i; 105c4762a1bSJed Brown 106327415f7SBarry Smith PetscFunctionBeginUser; 1079566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 108c4762a1bSJed Brown user.mx = 8; 109d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "elliptic example", NULL); 1109566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mx", "Number of grid points in each direction", "", user.mx, &user.mx, NULL)); 111c4762a1bSJed Brown user.ns = 6; 1129566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ns", "Number of data samples (1<=ns<=8)", "", user.ns, &user.ns, 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 = 0.1; 1169566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-alpha", "Regularization parameter", "", user.alpha, &user.alpha, NULL)); 117c4762a1bSJed Brown user.beta = 0.00001; 1189566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-beta", "Weight attributed to ||u||^2 in regularization functional", "", user.beta, &user.beta, NULL)); 119c4762a1bSJed Brown user.noise = 0.01; 1209566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-noise", "Amount of noise to add to data", "", user.noise, &user.noise, NULL)); 121c4762a1bSJed Brown 122c4762a1bSJed Brown user.use_ptap = PETSC_FALSE; 1239566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-use_ptap", "Use ptap matrix for DSG", "", user.use_ptap, &user.use_ptap, NULL)); 124c4762a1bSJed Brown user.use_lrc = PETSC_FALSE; 1259566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-use_lrc", "Use lrc matrix for Js", "", user.use_lrc, &user.use_lrc, NULL)); 1269566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ntests", "Number of times to repeat TaoSolve", "", ntests, &ntests, NULL)); 127d0609cedSBarry Smith PetscOptionsEnd(); 12876280437SVaclav Hapla 129c4762a1bSJed Brown user.m = user.ns * user.mx * user.mx * user.mx; /* number of constraints */ 130c4762a1bSJed Brown user.nstate = user.m; 131c4762a1bSJed Brown user.ndesign = user.mx * user.mx * user.mx; 132c4762a1bSJed Brown user.n = user.nstate + user.ndesign; /* number of variables */ 133c4762a1bSJed Brown 134c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 1359566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao)); 1369566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOLCL)); 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* Set up initial vectors and matrices */ 1399566063dSJacob Faibussowitsch PetscCall(EllipticInitialize(&user)); 140c4762a1bSJed Brown 1419566063dSJacob Faibussowitsch PetscCall(Gather(user.x, user.y, user.state_scatter, user.u, user.design_scatter)); 1429566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user.x, &x0)); 1439566063dSJacob Faibussowitsch PetscCall(VecCopy(user.x, x0)); 144c4762a1bSJed Brown 145c4762a1bSJed Brown /* Set solution vector with an initial guess */ 1469566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, user.x)); 1479566063dSJacob Faibussowitsch PetscCall(TaoSetObjective(tao, FormFunction, &user)); 1489566063dSJacob Faibussowitsch PetscCall(TaoSetGradient(tao, NULL, FormGradient, &user)); 1499566063dSJacob Faibussowitsch PetscCall(TaoSetConstraintsRoutine(tao, user.c, FormConstraints, &user)); 150c4762a1bSJed Brown 1519566063dSJacob Faibussowitsch PetscCall(TaoSetJacobianStateRoutine(tao, user.Js, NULL, user.JsInv, FormJacobianState, &user)); 1529566063dSJacob Faibussowitsch PetscCall(TaoSetJacobianDesignRoutine(tao, user.Jd, FormJacobianDesign, &user)); 153c4762a1bSJed Brown 1549566063dSJacob Faibussowitsch PetscCall(TaoSetStateDesignIS(tao, user.s_is, user.d_is)); 1559566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 156c4762a1bSJed Brown 157c4762a1bSJed Brown /* SOLVE THE APPLICATION */ 1589566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("Trials", &user.stages[1])); 1599566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(user.stages[1])); 160c4762a1bSJed Brown for (i = 0; i < ntests; i++) { 1619566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 16263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "KSP Iterations = %" PetscInt_FMT "\n", user.ksp_its)); 1639566063dSJacob Faibussowitsch PetscCall(VecCopy(x0, user.x)); 164c4762a1bSJed Brown } 1659566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 1669566063dSJacob Faibussowitsch PetscCall(PetscBarrier((PetscObject)user.x)); 1679566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "KSP iterations within initialization: ")); 16863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "\n", user.ksp_its_initial)); 169c4762a1bSJed Brown 1709566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 1719566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x0)); 1729566063dSJacob Faibussowitsch PetscCall(EllipticDestroy(&user)); 1739566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 174b122ec5aSJacob Faibussowitsch return 0; 175c4762a1bSJed Brown } 176c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 177c4762a1bSJed Brown /* 178c4762a1bSJed Brown dwork = Qy - d 179c4762a1bSJed Brown lwork = L*(u-ur) 180c4762a1bSJed Brown f = 1/2 * (dwork.dwork + alpha*lwork.lwork) 181c4762a1bSJed Brown */ 182*9371c9d4SSatish Balay PetscErrorCode FormFunction(Tao tao, Vec X, PetscReal *f, void *ptr) { 183c4762a1bSJed Brown PetscReal d1 = 0, d2 = 0; 184c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 185c4762a1bSJed Brown 186c4762a1bSJed Brown PetscFunctionBegin; 1879566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter)); 1889566063dSJacob Faibussowitsch PetscCall(MatMult(user->MQ, user->y, user->dwork)); 1899566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->dwork, -1.0, user->d)); 1909566063dSJacob Faibussowitsch PetscCall(VecDot(user->dwork, user->dwork, &d1)); 1919566063dSJacob Faibussowitsch PetscCall(VecWAXPY(user->uwork, -1.0, user->ur, user->u)); 1929566063dSJacob Faibussowitsch PetscCall(MatMult(user->L, user->uwork, user->lwork)); 1939566063dSJacob Faibussowitsch PetscCall(VecDot(user->lwork, user->lwork, &d2)); 194c4762a1bSJed Brown *f = 0.5 * (d1 + user->alpha * d2); 195c4762a1bSJed Brown PetscFunctionReturn(0); 196c4762a1bSJed Brown } 197c4762a1bSJed Brown 198c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 199c4762a1bSJed Brown /* 200c4762a1bSJed Brown state: g_s = Q' *(Qy - d) 201c4762a1bSJed Brown design: g_d = alpha*L'*L*(u-ur) 202c4762a1bSJed Brown */ 203*9371c9d4SSatish Balay PetscErrorCode FormGradient(Tao tao, Vec X, Vec G, void *ptr) { 204c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 205c4762a1bSJed Brown 206c4762a1bSJed Brown PetscFunctionBegin; 2079566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter)); 2089566063dSJacob Faibussowitsch PetscCall(MatMult(user->MQ, user->y, user->dwork)); 2099566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->dwork, -1.0, user->d)); 2109566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->MQ, user->dwork, user->ywork)); 2119566063dSJacob Faibussowitsch PetscCall(VecWAXPY(user->uwork, -1.0, user->ur, user->u)); 2129566063dSJacob Faibussowitsch PetscCall(MatMult(user->L, user->uwork, user->lwork)); 2139566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->L, user->lwork, user->uwork)); 2149566063dSJacob Faibussowitsch PetscCall(VecScale(user->uwork, user->alpha)); 2159566063dSJacob Faibussowitsch PetscCall(Gather(G, user->ywork, user->state_scatter, user->uwork, user->design_scatter)); 216c4762a1bSJed Brown PetscFunctionReturn(0); 217c4762a1bSJed Brown } 218c4762a1bSJed Brown 219*9371c9d4SSatish Balay PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr) { 220c4762a1bSJed Brown PetscReal d1, d2; 221c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 222c4762a1bSJed Brown 223c4762a1bSJed Brown PetscFunctionBegin; 2249566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter)); 2259566063dSJacob Faibussowitsch PetscCall(MatMult(user->MQ, user->y, user->dwork)); 2269566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->dwork, -1.0, user->d)); 2279566063dSJacob Faibussowitsch PetscCall(VecDot(user->dwork, user->dwork, &d1)); 2289566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->MQ, user->dwork, user->ywork)); 229c4762a1bSJed Brown 2309566063dSJacob Faibussowitsch PetscCall(VecWAXPY(user->uwork, -1.0, user->ur, user->u)); 2319566063dSJacob Faibussowitsch PetscCall(MatMult(user->L, user->uwork, user->lwork)); 2329566063dSJacob Faibussowitsch PetscCall(VecDot(user->lwork, user->lwork, &d2)); 2339566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->L, user->lwork, user->uwork)); 2349566063dSJacob Faibussowitsch PetscCall(VecScale(user->uwork, user->alpha)); 235c4762a1bSJed Brown *f = 0.5 * (d1 + user->alpha * d2); 2369566063dSJacob Faibussowitsch PetscCall(Gather(G, user->ywork, user->state_scatter, user->uwork, user->design_scatter)); 237c4762a1bSJed Brown PetscFunctionReturn(0); 238c4762a1bSJed Brown } 239c4762a1bSJed Brown 240c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 241c4762a1bSJed Brown /* A 242c4762a1bSJed Brown MatShell object 243c4762a1bSJed Brown */ 244*9371c9d4SSatish Balay PetscErrorCode FormJacobianState(Tao tao, Vec X, Mat J, Mat JPre, Mat JInv, void *ptr) { 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)); 249c4762a1bSJed Brown /* DSG = Div * (1/Av_u) * Grad */ 2509566063dSJacob Faibussowitsch PetscCall(VecSet(user->uwork, 0)); 2519566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->uwork, -1.0, user->u)); 2529566063dSJacob Faibussowitsch PetscCall(VecExp(user->uwork)); 2539566063dSJacob Faibussowitsch PetscCall(MatMult(user->Av, user->uwork, user->Av_u)); 2549566063dSJacob Faibussowitsch PetscCall(VecCopy(user->Av_u, user->Swork)); 2559566063dSJacob Faibussowitsch PetscCall(VecReciprocal(user->Swork)); 256c4762a1bSJed Brown if (user->use_ptap) { 2579566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(user->Diag, user->Swork, INSERT_VALUES)); 2589566063dSJacob Faibussowitsch PetscCall(MatPtAP(user->Diag, user->Grad, MAT_REUSE_MATRIX, 1.0, &user->DSG)); 259c4762a1bSJed Brown } else { 2609566063dSJacob Faibussowitsch PetscCall(MatCopy(user->Div, user->Divwork, SAME_NONZERO_PATTERN)); 2619566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(user->Divwork, NULL, user->Swork)); 2629566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(user->DSG)); 263c4762a1bSJed Brown } 264c4762a1bSJed Brown PetscFunctionReturn(0); 265c4762a1bSJed Brown } 266c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 267c4762a1bSJed Brown /* B */ 268*9371c9d4SSatish Balay PetscErrorCode FormJacobianDesign(Tao tao, Vec X, Mat J, void *ptr) { 269c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 270c4762a1bSJed Brown 271c4762a1bSJed Brown PetscFunctionBegin; 2729566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter)); 273c4762a1bSJed Brown PetscFunctionReturn(0); 274c4762a1bSJed Brown } 275c4762a1bSJed Brown 276*9371c9d4SSatish Balay PetscErrorCode StateBlockMatMult(Mat J_shell, Vec X, Vec Y) { 277c4762a1bSJed Brown PetscReal sum; 278c4762a1bSJed Brown AppCtx *user; 279c4762a1bSJed Brown 280c4762a1bSJed Brown PetscFunctionBegin; 2819566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user)); 2829566063dSJacob Faibussowitsch PetscCall(MatMult(user->DSG, X, Y)); 2839566063dSJacob Faibussowitsch PetscCall(VecSum(X, &sum)); 284c4762a1bSJed Brown sum /= user->ndesign; 2859566063dSJacob Faibussowitsch PetscCall(VecShift(Y, sum)); 286c4762a1bSJed Brown PetscFunctionReturn(0); 287c4762a1bSJed Brown } 288c4762a1bSJed Brown 289*9371c9d4SSatish Balay PetscErrorCode StateMatMult(Mat J_shell, Vec X, Vec Y) { 290c4762a1bSJed Brown PetscInt i; 291c4762a1bSJed Brown AppCtx *user; 292c4762a1bSJed Brown 293c4762a1bSJed Brown PetscFunctionBegin; 2949566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user)); 295c4762a1bSJed Brown if (user->ns == 1) { 2969566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock, X, Y)); 297c4762a1bSJed Brown } else { 298c4762a1bSJed Brown for (i = 0; i < user->ns; i++) { 2999566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->subq, user->yi_scatter[i], 0, 0)); 3009566063dSJacob Faibussowitsch PetscCall(Scatter(Y, user->suby, user->yi_scatter[i], 0, 0)); 3019566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock, user->subq, user->suby)); 3029566063dSJacob Faibussowitsch PetscCall(Gather(Y, user->suby, user->yi_scatter[i], 0, 0)); 303c4762a1bSJed Brown } 304c4762a1bSJed Brown } 305c4762a1bSJed Brown PetscFunctionReturn(0); 306c4762a1bSJed Brown } 307c4762a1bSJed Brown 308*9371c9d4SSatish Balay PetscErrorCode StateInvMatMult(Mat J_shell, Vec X, Vec Y) { 309c4762a1bSJed Brown PetscInt its, i; 310c4762a1bSJed Brown AppCtx *user; 311c4762a1bSJed Brown 312c4762a1bSJed Brown PetscFunctionBegin; 3139566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user)); 3149566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(user->solver, user->JsBlock, user->DSG)); 315c4762a1bSJed Brown if (Y == user->ytrue) { 316c4762a1bSJed Brown /* First solve is done using true solution to set up problem */ 3179566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(user->solver, 1e-8, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 318c4762a1bSJed Brown } else { 3199566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(user->solver, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 320c4762a1bSJed Brown } 321c4762a1bSJed Brown if (user->ns == 1) { 3229566063dSJacob Faibussowitsch PetscCall(KSPSolve(user->solver, X, Y)); 3239566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(user->solver, &its)); 324c4762a1bSJed Brown user->ksp_its += its; 325c4762a1bSJed Brown } else { 326c4762a1bSJed Brown for (i = 0; i < user->ns; i++) { 3279566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->subq, user->yi_scatter[i], 0, 0)); 3289566063dSJacob Faibussowitsch PetscCall(Scatter(Y, user->suby, user->yi_scatter[i], 0, 0)); 3299566063dSJacob Faibussowitsch PetscCall(KSPSolve(user->solver, user->subq, user->suby)); 3309566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(user->solver, &its)); 331c4762a1bSJed Brown user->ksp_its += its; 3329566063dSJacob Faibussowitsch PetscCall(Gather(Y, user->suby, user->yi_scatter[i], 0, 0)); 333c4762a1bSJed Brown } 334c4762a1bSJed Brown } 335c4762a1bSJed Brown PetscFunctionReturn(0); 336c4762a1bSJed Brown } 337*9371c9d4SSatish Balay PetscErrorCode QMatMult(Mat J_shell, Vec X, Vec Y) { 338c4762a1bSJed Brown AppCtx *user; 339c4762a1bSJed Brown PetscInt i; 340c4762a1bSJed Brown 341c4762a1bSJed Brown PetscFunctionBegin; 3429566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user)); 343c4762a1bSJed Brown if (user->ns == 1) { 3449566063dSJacob Faibussowitsch PetscCall(MatMult(user->Q, X, Y)); 345c4762a1bSJed Brown } else { 346c4762a1bSJed Brown for (i = 0; i < user->ns; i++) { 3479566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->subq, user->yi_scatter[i], 0, 0)); 3489566063dSJacob Faibussowitsch PetscCall(Scatter(Y, user->subd, user->di_scatter[i], 0, 0)); 3499566063dSJacob Faibussowitsch PetscCall(MatMult(user->Q, user->subq, user->subd)); 3509566063dSJacob Faibussowitsch PetscCall(Gather(Y, user->subd, user->di_scatter[i], 0, 0)); 351c4762a1bSJed Brown } 352c4762a1bSJed Brown } 353c4762a1bSJed Brown PetscFunctionReturn(0); 354c4762a1bSJed Brown } 355c4762a1bSJed Brown 356*9371c9d4SSatish Balay PetscErrorCode QMatMultTranspose(Mat J_shell, Vec X, Vec Y) { 357c4762a1bSJed Brown AppCtx *user; 358c4762a1bSJed Brown PetscInt i; 359c4762a1bSJed Brown 360c4762a1bSJed Brown PetscFunctionBegin; 3619566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user)); 362c4762a1bSJed Brown if (user->ns == 1) { 3639566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->Q, X, Y)); 364c4762a1bSJed Brown } else { 365c4762a1bSJed Brown for (i = 0; i < user->ns; i++) { 3669566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->subd, user->di_scatter[i], 0, 0)); 3679566063dSJacob Faibussowitsch PetscCall(Scatter(Y, user->suby, user->yi_scatter[i], 0, 0)); 3689566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->Q, user->subd, user->suby)); 3699566063dSJacob Faibussowitsch PetscCall(Gather(Y, user->suby, user->yi_scatter[i], 0, 0)); 370c4762a1bSJed Brown } 371c4762a1bSJed Brown } 372c4762a1bSJed Brown PetscFunctionReturn(0); 373c4762a1bSJed Brown } 374c4762a1bSJed Brown 375*9371c9d4SSatish Balay PetscErrorCode DesignMatMult(Mat J_shell, Vec X, Vec Y) { 376c4762a1bSJed Brown PetscInt i; 377c4762a1bSJed Brown AppCtx *user; 378c4762a1bSJed Brown 379c4762a1bSJed Brown PetscFunctionBegin; 3809566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user)); 381c4762a1bSJed Brown 382c4762a1bSJed Brown /* sdiag(1./v) */ 3839566063dSJacob Faibussowitsch PetscCall(VecSet(user->uwork, 0)); 3849566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->uwork, -1.0, user->u)); 3859566063dSJacob Faibussowitsch PetscCall(VecExp(user->uwork)); 386c4762a1bSJed Brown 387c4762a1bSJed Brown /* sdiag(1./((Av*(1./v)).^2)) */ 3889566063dSJacob Faibussowitsch PetscCall(MatMult(user->Av, user->uwork, user->Swork)); 3899566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Swork, user->Swork, user->Swork)); 3909566063dSJacob Faibussowitsch PetscCall(VecReciprocal(user->Swork)); 391c4762a1bSJed Brown 392c4762a1bSJed Brown /* (Av * (sdiag(1./v) * b)) */ 3939566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uwork, user->uwork, X)); 3949566063dSJacob Faibussowitsch PetscCall(MatMult(user->Av, user->uwork, user->Twork)); 395c4762a1bSJed Brown 396c4762a1bSJed Brown /* (sdiag(1./((Av*(1./v)).^2)) * (Av * (sdiag(1./v) * b))) */ 3979566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Swork, user->Twork, user->Swork)); 398c4762a1bSJed Brown 399c4762a1bSJed Brown if (user->ns == 1) { 400c4762a1bSJed Brown /* (sdiag(Grad*y(:,i)) */ 4019566063dSJacob Faibussowitsch PetscCall(MatMult(user->Grad, user->y, user->Twork)); 402c4762a1bSJed Brown 403c4762a1bSJed Brown /* Div * (sdiag(Grad*y(:,i)) * (sdiag(1./((Av*(1./v)).^2)) * (Av * (sdiag(1./v) * b)))) */ 4049566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Swork, user->Twork, user->Swork)); 4059566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->Grad, user->Swork, Y)); 406c4762a1bSJed Brown } else { 407c4762a1bSJed Brown for (i = 0; i < user->ns; i++) { 4089566063dSJacob Faibussowitsch PetscCall(Scatter(user->y, user->suby, user->yi_scatter[i], 0, 0)); 4099566063dSJacob Faibussowitsch PetscCall(Scatter(Y, user->subq, user->yi_scatter[i], 0, 0)); 410c4762a1bSJed Brown 4119566063dSJacob Faibussowitsch PetscCall(MatMult(user->Grad, user->suby, user->Twork)); 4129566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Twork, user->Twork, user->Swork)); 4139566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->Grad, user->Twork, user->subq)); 4149566063dSJacob Faibussowitsch PetscCall(Gather(user->y, user->suby, user->yi_scatter[i], 0, 0)); 4159566063dSJacob Faibussowitsch PetscCall(Gather(Y, user->subq, user->yi_scatter[i], 0, 0)); 416c4762a1bSJed Brown } 417c4762a1bSJed Brown } 418c4762a1bSJed Brown PetscFunctionReturn(0); 419c4762a1bSJed Brown } 420c4762a1bSJed Brown 421*9371c9d4SSatish Balay PetscErrorCode DesignMatMultTranspose(Mat J_shell, Vec X, Vec Y) { 422c4762a1bSJed Brown PetscInt i; 423c4762a1bSJed Brown AppCtx *user; 424c4762a1bSJed Brown 425c4762a1bSJed Brown PetscFunctionBegin; 4269566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell, &user)); 4279566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Y)); 428c4762a1bSJed Brown 429c4762a1bSJed Brown /* Sdiag = 1./((Av*(1./v)).^2) */ 4309566063dSJacob Faibussowitsch PetscCall(VecSet(user->uwork, 0)); 4319566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->uwork, -1.0, user->u)); 4329566063dSJacob Faibussowitsch PetscCall(VecExp(user->uwork)); 4339566063dSJacob Faibussowitsch PetscCall(MatMult(user->Av, user->uwork, user->Swork)); 4349566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Sdiag, user->Swork, user->Swork)); 4359566063dSJacob Faibussowitsch PetscCall(VecReciprocal(user->Sdiag)); 436c4762a1bSJed Brown 437c4762a1bSJed Brown for (i = 0; i < user->ns; i++) { 4389566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->subq, user->yi_scatter[i], 0, 0)); 4399566063dSJacob Faibussowitsch PetscCall(Scatter(user->y, user->suby, user->yi_scatter[i], 0, 0)); 440c4762a1bSJed Brown 441c4762a1bSJed Brown /* Swork = (Div' * b(:,i)) */ 4429566063dSJacob Faibussowitsch PetscCall(MatMult(user->Grad, user->subq, user->Swork)); 443c4762a1bSJed Brown 444c4762a1bSJed Brown /* Twork = Grad*y(:,i) */ 4459566063dSJacob Faibussowitsch PetscCall(MatMult(user->Grad, user->suby, user->Twork)); 446c4762a1bSJed Brown 447c4762a1bSJed Brown /* Twork = sdiag(Twork) * Swork */ 4489566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Twork, user->Swork, user->Twork)); 449c4762a1bSJed Brown 450c4762a1bSJed Brown /* Swork = pointwisemult(Sdiag,Twork) */ 4519566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Swork, user->Twork, user->Sdiag)); 452c4762a1bSJed Brown 453c4762a1bSJed Brown /* Ywork = Av' * Swork */ 4549566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->Av, user->Swork, user->Ywork)); 455c4762a1bSJed Brown 456c4762a1bSJed Brown /* Ywork = pointwisemult(uwork,Ywork) */ 4579566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Ywork, user->uwork, user->Ywork)); 4589566063dSJacob Faibussowitsch PetscCall(VecAXPY(Y, 1.0, user->Ywork)); 4599566063dSJacob Faibussowitsch PetscCall(Gather(user->y, user->suby, user->yi_scatter[i], 0, 0)); 460c4762a1bSJed Brown } 461c4762a1bSJed Brown PetscFunctionReturn(0); 462c4762a1bSJed Brown } 463c4762a1bSJed Brown 464*9371c9d4SSatish Balay PetscErrorCode FormConstraints(Tao tao, Vec X, Vec C, void *ptr) { 465c4762a1bSJed Brown /* C=Ay - q A = Div * Sigma * Grad + hx*hx*hx*ones(n,n) */ 466c4762a1bSJed Brown PetscReal sum; 467c4762a1bSJed Brown PetscInt i; 468c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 469c4762a1bSJed Brown 470c4762a1bSJed Brown PetscFunctionBegin; 4719566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter)); 472c4762a1bSJed Brown if (user->ns == 1) { 4739566063dSJacob Faibussowitsch PetscCall(MatMult(user->Grad, user->y, user->Swork)); 4749566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(user->Swork, user->Swork, user->Av_u)); 4759566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->Grad, user->Swork, C)); 4769566063dSJacob Faibussowitsch PetscCall(VecSum(user->y, &sum)); 477c4762a1bSJed Brown sum /= user->ndesign; 4789566063dSJacob Faibussowitsch PetscCall(VecShift(C, sum)); 479c4762a1bSJed Brown } else { 480c4762a1bSJed Brown for (i = 0; i < user->ns; i++) { 4819566063dSJacob Faibussowitsch PetscCall(Scatter(user->y, user->suby, user->yi_scatter[i], 0, 0)); 4829566063dSJacob Faibussowitsch PetscCall(Scatter(C, user->subq, user->yi_scatter[i], 0, 0)); 4839566063dSJacob Faibussowitsch PetscCall(MatMult(user->Grad, user->suby, user->Swork)); 4849566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(user->Swork, user->Swork, user->Av_u)); 4859566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->Grad, user->Swork, user->subq)); 486c4762a1bSJed Brown 4879566063dSJacob Faibussowitsch PetscCall(VecSum(user->suby, &sum)); 488c4762a1bSJed Brown sum /= user->ndesign; 4899566063dSJacob Faibussowitsch PetscCall(VecShift(user->subq, sum)); 490c4762a1bSJed Brown 4919566063dSJacob Faibussowitsch PetscCall(Gather(user->y, user->suby, user->yi_scatter[i], 0, 0)); 4929566063dSJacob Faibussowitsch PetscCall(Gather(C, user->subq, user->yi_scatter[i], 0, 0)); 493c4762a1bSJed Brown } 494c4762a1bSJed Brown } 4959566063dSJacob Faibussowitsch PetscCall(VecAXPY(C, -1.0, user->q)); 496c4762a1bSJed Brown PetscFunctionReturn(0); 497c4762a1bSJed Brown } 498c4762a1bSJed Brown 499*9371c9d4SSatish Balay PetscErrorCode Scatter(Vec x, Vec sub1, VecScatter scat1, Vec sub2, VecScatter scat2) { 500c4762a1bSJed Brown PetscFunctionBegin; 5019566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scat1, x, sub1, INSERT_VALUES, SCATTER_FORWARD)); 5029566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scat1, x, sub1, INSERT_VALUES, SCATTER_FORWARD)); 503c4762a1bSJed Brown if (sub2) { 5049566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scat2, x, sub2, INSERT_VALUES, SCATTER_FORWARD)); 5059566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scat2, x, sub2, INSERT_VALUES, SCATTER_FORWARD)); 506c4762a1bSJed Brown } 507c4762a1bSJed Brown PetscFunctionReturn(0); 508c4762a1bSJed Brown } 509c4762a1bSJed Brown 510*9371c9d4SSatish Balay PetscErrorCode Gather(Vec x, Vec sub1, VecScatter scat1, Vec sub2, VecScatter scat2) { 511c4762a1bSJed Brown PetscFunctionBegin; 5129566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scat1, sub1, x, INSERT_VALUES, SCATTER_REVERSE)); 5139566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scat1, sub1, x, INSERT_VALUES, SCATTER_REVERSE)); 514c4762a1bSJed Brown if (sub2) { 5159566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scat2, sub2, x, INSERT_VALUES, SCATTER_REVERSE)); 5169566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scat2, sub2, x, INSERT_VALUES, SCATTER_REVERSE)); 517c4762a1bSJed Brown } 518c4762a1bSJed Brown PetscFunctionReturn(0); 519c4762a1bSJed Brown } 520c4762a1bSJed Brown 521*9371c9d4SSatish Balay PetscErrorCode EllipticInitialize(AppCtx *user) { 522c4762a1bSJed Brown PetscInt m, n, i, j, k, l, linear_index, is, js, ks, ls, istart, iend, iblock; 523c4762a1bSJed Brown Vec XX, YY, ZZ, XXwork, YYwork, ZZwork, UTwork; 524c4762a1bSJed Brown PetscReal *x, *y, *z; 525c4762a1bSJed Brown PetscReal h, meanut; 526c4762a1bSJed Brown PetscScalar hinv, neg_hinv, half = 0.5, sqrt_beta; 527c4762a1bSJed Brown PetscInt im, indx1, indx2, indy1, indy2, indz1, indz2, nx, ny, nz; 528c4762a1bSJed Brown IS is_alldesign, is_allstate; 529c4762a1bSJed Brown IS is_from_d; 530c4762a1bSJed Brown IS is_from_y; 531c4762a1bSJed Brown PetscInt lo, hi, hi2, lo2, ysubnlocal, dsubnlocal; 532c4762a1bSJed Brown const PetscInt *ranges, *subranges; 533c4762a1bSJed Brown PetscMPIInt size; 534c4762a1bSJed Brown PetscReal xri, yri, zri, xim, yim, zim, dx1, dx2, dy1, dy2, dz1, dz2, Dx, Dy, Dz; 535c4762a1bSJed Brown PetscScalar v, vx, vy, vz; 536c4762a1bSJed Brown PetscInt offset, subindex, subvec, nrank, kk; 537c4762a1bSJed Brown 538*9371c9d4SSatish Balay PetscScalar xr[64] = {0.4970, 0.8498, 0.7814, 0.6268, 0.7782, 0.6402, 0.3617, 0.3160, 0.3610, 0.5298, 0.6987, 0.3331, 0.7962, 0.5596, 0.3866, 0.6774, 0.5407, 0.4518, 0.6702, 0.6061, 0.7580, 0.8997, 539*9371c9d4SSatish Balay 0.5198, 0.8326, 0.2138, 0.9198, 0.3000, 0.2833, 0.8288, 0.7076, 0.1820, 0.0728, 0.8447, 0.2367, 0.3239, 0.6413, 0.3114, 0.4731, 0.1192, 0.9273, 0.5724, 0.4331, 0.5136, 0.3547, 540*9371c9d4SSatish Balay 0.4413, 0.2602, 0.5698, 0.7278, 0.5261, 0.6230, 0.2454, 0.3948, 0.7479, 0.6582, 0.4660, 0.5594, 0.7574, 0.1143, 0.5900, 0.1065, 0.4260, 0.3294, 0.8276, 0.0756}; 541c4762a1bSJed Brown 542*9371c9d4SSatish Balay PetscScalar yr[64] = {0.7345, 0.9120, 0.9288, 0.7528, 0.4463, 0.4985, 0.2497, 0.6256, 0.3425, 0.9026, 0.6983, 0.4230, 0.7140, 0.2970, 0.4474, 0.8792, 0.6604, 0.2485, 0.7968, 0.6127, 0.1796, 0.2437, 543*9371c9d4SSatish Balay 0.5938, 0.6137, 0.3867, 0.5658, 0.4575, 0.1009, 0.0863, 0.3361, 0.0738, 0.3985, 0.6602, 0.1437, 0.0934, 0.5983, 0.5950, 0.0763, 0.0768, 0.2288, 0.5761, 0.1129, 0.3841, 0.6150, 544*9371c9d4SSatish Balay 0.6904, 0.6686, 0.1361, 0.4601, 0.4491, 0.3716, 0.1969, 0.6537, 0.6743, 0.6991, 0.4811, 0.5480, 0.1684, 0.4569, 0.6889, 0.8437, 0.3015, 0.2854, 0.8199, 0.2658}; 545c4762a1bSJed Brown 546*9371c9d4SSatish Balay PetscScalar zr[64] = {0.7668, 0.8573, 0.2654, 0.2719, 0.1060, 0.1311, 0.6232, 0.2295, 0.8009, 0.2147, 0.2119, 0.9325, 0.4473, 0.3600, 0.3374, 0.3819, 0.4066, 0.5801, 0.1673, 0.0959, 0.4638, 0.8236, 547*9371c9d4SSatish Balay 0.8800, 0.2939, 0.2028, 0.8262, 0.2706, 0.6276, 0.9085, 0.6443, 0.8241, 0.0712, 0.1824, 0.7789, 0.4389, 0.8415, 0.7055, 0.6639, 0.3653, 0.2078, 0.1987, 0.2297, 0.4321, 0.8115, 548*9371c9d4SSatish Balay 0.4915, 0.7764, 0.4657, 0.4627, 0.4569, 0.4232, 0.8514, 0.0674, 0.3227, 0.1055, 0.6690, 0.6313, 0.9226, 0.5461, 0.4126, 0.2364, 0.6096, 0.7042, 0.3914, 0.0711}; 549c4762a1bSJed Brown 550c4762a1bSJed Brown PetscFunctionBegin; 5519566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 5529566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("Elliptic Setup", &user->stages[0])); 5539566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(user->stages[0])); 554c4762a1bSJed Brown 555c4762a1bSJed Brown /* Create u,y,c,x */ 5569566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->u)); 5579566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->y)); 5589566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->c)); 5599566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->u, PETSC_DECIDE, user->ndesign)); 5609566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->u)); 5619566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(user->u, &ysubnlocal)); 5629566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->y, ysubnlocal * user->ns, user->nstate)); 5639566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->c, ysubnlocal * user->ns, user->m)); 5649566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->y)); 5659566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->c)); 566c4762a1bSJed Brown 567c4762a1bSJed Brown /* 568c4762a1bSJed Brown ******************************* 569c4762a1bSJed Brown Create scatters for x <-> y,u 570c4762a1bSJed Brown ******************************* 571c4762a1bSJed Brown 572c4762a1bSJed Brown If the state vector y and design vector u are partitioned as 573c4762a1bSJed Brown [y_1; y_2; ...; y_np] and [u_1; u_2; ...; u_np] (with np = # of processors), 574c4762a1bSJed Brown then the solution vector x is organized as 575c4762a1bSJed Brown [y_1; u_1; y_2; u_2; ...; y_np; u_np]. 576c4762a1bSJed Brown The index sets user->s_is and user->d_is correspond to the indices of the 577c4762a1bSJed Brown state and design variables owned by the current processor. 578c4762a1bSJed Brown */ 5799566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->x)); 580c4762a1bSJed Brown 5819566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->y, &lo, &hi)); 5829566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->u, &lo2, &hi2)); 583c4762a1bSJed Brown 5849566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_allstate)); 5859566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo + lo2, 1, &user->s_is)); 5869566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi2 - lo2, lo2, 1, &is_alldesign)); 5879566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi2 - lo2, hi + lo2, 1, &user->d_is)); 588c4762a1bSJed Brown 5899566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->x, hi - lo + hi2 - lo2, user->n)); 5909566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->x)); 591c4762a1bSJed Brown 5929566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(user->x, user->s_is, user->y, is_allstate, &user->state_scatter)); 5939566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(user->x, user->d_is, user->u, is_alldesign, &user->design_scatter)); 5949566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_alldesign)); 5959566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_allstate)); 596c4762a1bSJed Brown /* 597c4762a1bSJed Brown ******************************* 598c4762a1bSJed Brown Create scatter from y to y_1,y_2,...,y_ns 599c4762a1bSJed Brown ******************************* 600c4762a1bSJed Brown */ 6019566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->ns, &user->yi_scatter)); 6029566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->u, &user->suby)); 6039566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->u, &user->subq)); 604c4762a1bSJed Brown 6059566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->y, &lo2, &hi2)); 606c4762a1bSJed Brown istart = 0; 607c4762a1bSJed Brown for (i = 0; i < user->ns; i++) { 6089566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->suby, &lo, &hi)); 6099566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo2 + istart, 1, &is_from_y)); 6109566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(user->y, is_from_y, user->suby, NULL, &user->yi_scatter[i])); 611c4762a1bSJed Brown istart = istart + hi - lo; 6129566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_from_y)); 613c4762a1bSJed Brown } 614c4762a1bSJed Brown /* 615c4762a1bSJed Brown ******************************* 616c4762a1bSJed Brown Create scatter from d to d_1,d_2,...,d_ns 617c4762a1bSJed Brown ******************************* 618c4762a1bSJed Brown */ 6199566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->subd)); 6209566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->subd, PETSC_DECIDE, user->ndata)); 6219566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->subd)); 6229566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->d)); 6239566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(user->subd, &dsubnlocal)); 6249566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->d, dsubnlocal * user->ns, user->ndata * user->ns)); 6259566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->d)); 6269566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->ns, &user->di_scatter)); 627c4762a1bSJed Brown 6289566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->d, &lo2, &hi2)); 629c4762a1bSJed Brown istart = 0; 630c4762a1bSJed Brown for (i = 0; i < user->ns; i++) { 6319566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->subd, &lo, &hi)); 6329566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo2 + istart, 1, &is_from_d)); 6339566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(user->d, is_from_d, user->subd, NULL, &user->di_scatter[i])); 634c4762a1bSJed Brown istart = istart + hi - lo; 6359566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_from_d)); 636c4762a1bSJed Brown } 637c4762a1bSJed Brown 6389566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->mx, &x)); 6399566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->mx, &y)); 6409566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->mx, &z)); 641c4762a1bSJed Brown 642c4762a1bSJed Brown user->ksp_its = 0; 643c4762a1bSJed Brown user->ksp_its_initial = 0; 644c4762a1bSJed Brown 645c4762a1bSJed Brown n = user->mx * user->mx * user->mx; 646c4762a1bSJed Brown m = 3 * user->mx * user->mx * (user->mx - 1); 647c4762a1bSJed Brown sqrt_beta = PetscSqrtScalar(user->beta); 648c4762a1bSJed Brown 6499566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &XX)); 6509566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->q)); 6519566063dSJacob Faibussowitsch PetscCall(VecSetSizes(XX, ysubnlocal, n)); 6529566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->q, ysubnlocal * user->ns, user->m)); 6539566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(XX)); 6549566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->q)); 655c4762a1bSJed Brown 6569566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &YY)); 6579566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &ZZ)); 6589566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &XXwork)); 6599566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &YYwork)); 6609566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &ZZwork)); 6619566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &UTwork)); 6629566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX, &user->utrue)); 663c4762a1bSJed Brown 664c4762a1bSJed Brown /* map for striding q */ 6659566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRanges(user->q, &ranges)); 6669566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRanges(user->u, &subranges)); 667c4762a1bSJed Brown 6689566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->q, &lo2, &hi2)); 6699566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->u, &lo, &hi)); 670c4762a1bSJed Brown /* Generate 3D grid, and collect ns (1<=ns<=8) right-hand-side vectors into user->q */ 671c4762a1bSJed Brown h = 1.0 / user->mx; 672c4762a1bSJed Brown hinv = user->mx; 673c4762a1bSJed Brown neg_hinv = -hinv; 674c4762a1bSJed Brown 6759566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(XX, &istart, &iend)); 676c4762a1bSJed Brown for (linear_index = istart; linear_index < iend; linear_index++) { 677c4762a1bSJed Brown i = linear_index % user->mx; 678c4762a1bSJed Brown j = ((linear_index - i) / user->mx) % user->mx; 679c4762a1bSJed Brown k = ((linear_index - i) / user->mx - j) / user->mx; 680c4762a1bSJed Brown vx = h * (i + 0.5); 681c4762a1bSJed Brown vy = h * (j + 0.5); 682c4762a1bSJed Brown vz = h * (k + 0.5); 6839566063dSJacob Faibussowitsch PetscCall(VecSetValues(XX, 1, &linear_index, &vx, INSERT_VALUES)); 6849566063dSJacob Faibussowitsch PetscCall(VecSetValues(YY, 1, &linear_index, &vy, INSERT_VALUES)); 6859566063dSJacob Faibussowitsch PetscCall(VecSetValues(ZZ, 1, &linear_index, &vz, INSERT_VALUES)); 686c4762a1bSJed Brown for (is = 0; is < 2; is++) { 687c4762a1bSJed Brown for (js = 0; js < 2; js++) { 688c4762a1bSJed Brown for (ks = 0; ks < 2; ks++) { 689c4762a1bSJed Brown ls = is * 4 + js * 2 + ks; 690c4762a1bSJed Brown if (ls < user->ns) { 691c4762a1bSJed Brown l = ls * n + linear_index; 692c4762a1bSJed Brown /* remap */ 693c4762a1bSJed Brown subindex = l % n; 694c4762a1bSJed Brown subvec = l / n; 695c4762a1bSJed Brown nrank = 0; 696c4762a1bSJed Brown while (subindex >= subranges[nrank + 1]) nrank++; 697c4762a1bSJed Brown offset = subindex - subranges[nrank]; 698c4762a1bSJed Brown istart = 0; 699c4762a1bSJed Brown for (kk = 0; kk < nrank; kk++) istart += user->ns * (subranges[kk + 1] - subranges[kk]); 700c4762a1bSJed Brown istart += (subranges[nrank + 1] - subranges[nrank]) * subvec; 701c4762a1bSJed Brown l = istart + offset; 702c4762a1bSJed Brown v = 100 * PetscSinScalar(2 * PETSC_PI * (vx + 0.25 * is)) * PetscSinScalar(2 * PETSC_PI * (vy + 0.25 * js)) * PetscSinScalar(2 * PETSC_PI * (vz + 0.25 * ks)); 7039566063dSJacob Faibussowitsch PetscCall(VecSetValues(user->q, 1, &l, &v, INSERT_VALUES)); 704c4762a1bSJed Brown } 705c4762a1bSJed Brown } 706c4762a1bSJed Brown } 707c4762a1bSJed Brown } 708c4762a1bSJed Brown } 709c4762a1bSJed Brown 7109566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(XX)); 7119566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(XX)); 7129566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(YY)); 7139566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(YY)); 7149566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(ZZ)); 7159566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(ZZ)); 7169566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(user->q)); 7179566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(user->q)); 718c4762a1bSJed Brown 719c4762a1bSJed Brown /* Compute true parameter function 720c4762a1bSJed Brown ut = exp(-((x-0.25)^2+(y-0.25)^2+(z-0.25)^2)/0.05) - exp((x-0.75)^2-(y-0.75)^2-(z-0.75))^2/0.05) */ 7219566063dSJacob Faibussowitsch PetscCall(VecCopy(XX, XXwork)); 7229566063dSJacob Faibussowitsch PetscCall(VecCopy(YY, YYwork)); 7239566063dSJacob Faibussowitsch PetscCall(VecCopy(ZZ, ZZwork)); 724c4762a1bSJed Brown 7259566063dSJacob Faibussowitsch PetscCall(VecShift(XXwork, -0.25)); 7269566063dSJacob Faibussowitsch PetscCall(VecShift(YYwork, -0.25)); 7279566063dSJacob Faibussowitsch PetscCall(VecShift(ZZwork, -0.25)); 728c4762a1bSJed Brown 7299566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(XXwork, XXwork, XXwork)); 7309566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(YYwork, YYwork, YYwork)); 7319566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(ZZwork, ZZwork, ZZwork)); 732c4762a1bSJed Brown 7339566063dSJacob Faibussowitsch PetscCall(VecCopy(XXwork, UTwork)); 7349566063dSJacob Faibussowitsch PetscCall(VecAXPY(UTwork, 1.0, YYwork)); 7359566063dSJacob Faibussowitsch PetscCall(VecAXPY(UTwork, 1.0, ZZwork)); 7369566063dSJacob Faibussowitsch PetscCall(VecScale(UTwork, -20.0)); 7379566063dSJacob Faibussowitsch PetscCall(VecExp(UTwork)); 7389566063dSJacob Faibussowitsch PetscCall(VecCopy(UTwork, user->utrue)); 739c4762a1bSJed Brown 7409566063dSJacob Faibussowitsch PetscCall(VecCopy(XX, XXwork)); 7419566063dSJacob Faibussowitsch PetscCall(VecCopy(YY, YYwork)); 7429566063dSJacob Faibussowitsch PetscCall(VecCopy(ZZ, ZZwork)); 743c4762a1bSJed Brown 7449566063dSJacob Faibussowitsch PetscCall(VecShift(XXwork, -0.75)); 7459566063dSJacob Faibussowitsch PetscCall(VecShift(YYwork, -0.75)); 7469566063dSJacob Faibussowitsch PetscCall(VecShift(ZZwork, -0.75)); 747c4762a1bSJed Brown 7489566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(XXwork, XXwork, XXwork)); 7499566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(YYwork, YYwork, YYwork)); 7509566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(ZZwork, ZZwork, ZZwork)); 751c4762a1bSJed Brown 7529566063dSJacob Faibussowitsch PetscCall(VecCopy(XXwork, UTwork)); 7539566063dSJacob Faibussowitsch PetscCall(VecAXPY(UTwork, 1.0, YYwork)); 7549566063dSJacob Faibussowitsch PetscCall(VecAXPY(UTwork, 1.0, ZZwork)); 7559566063dSJacob Faibussowitsch PetscCall(VecScale(UTwork, -20.0)); 7569566063dSJacob Faibussowitsch PetscCall(VecExp(UTwork)); 757c4762a1bSJed Brown 7589566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->utrue, -1.0, UTwork)); 759c4762a1bSJed Brown 7609566063dSJacob Faibussowitsch PetscCall(VecDestroy(&XX)); 7619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&YY)); 7629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ZZ)); 7639566063dSJacob Faibussowitsch PetscCall(VecDestroy(&XXwork)); 7649566063dSJacob Faibussowitsch PetscCall(VecDestroy(&YYwork)); 7659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ZZwork)); 7669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&UTwork)); 767c4762a1bSJed Brown 768c4762a1bSJed Brown /* Initial guess and reference model */ 7699566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->utrue, &user->ur)); 7709566063dSJacob Faibussowitsch PetscCall(VecSum(user->utrue, &meanut)); 771c4762a1bSJed Brown meanut = meanut / n; 7729566063dSJacob Faibussowitsch PetscCall(VecSet(user->ur, meanut)); 7739566063dSJacob Faibussowitsch PetscCall(VecCopy(user->ur, user->u)); 774c4762a1bSJed Brown 775c4762a1bSJed Brown /* Generate Grad matrix */ 7769566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Grad)); 7779566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Grad, PETSC_DECIDE, ysubnlocal, m, n)); 7789566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Grad)); 7799566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->Grad, 2, NULL, 2, NULL)); 7809566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->Grad, 2, NULL)); 7819566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->Grad, &istart, &iend)); 782c4762a1bSJed Brown 783c4762a1bSJed Brown for (i = istart; i < iend; i++) { 784c4762a1bSJed Brown if (i < m / 3) { 785c4762a1bSJed Brown iblock = i / (user->mx - 1); 786c4762a1bSJed Brown j = iblock * user->mx + (i % (user->mx - 1)); 7879566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES)); 788c4762a1bSJed Brown j = j + 1; 7899566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &hinv, INSERT_VALUES)); 790c4762a1bSJed Brown } 791c4762a1bSJed Brown if (i >= m / 3 && i < 2 * m / 3) { 792c4762a1bSJed Brown iblock = (i - m / 3) / (user->mx * (user->mx - 1)); 793c4762a1bSJed Brown j = iblock * user->mx * user->mx + ((i - m / 3) % (user->mx * (user->mx - 1))); 7949566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES)); 795c4762a1bSJed Brown j = j + user->mx; 7969566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &hinv, INSERT_VALUES)); 797c4762a1bSJed Brown } 798c4762a1bSJed Brown if (i >= 2 * m / 3) { 799c4762a1bSJed Brown j = i - 2 * m / 3; 8009566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES)); 801c4762a1bSJed Brown j = j + user->mx * user->mx; 8029566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &hinv, INSERT_VALUES)); 803c4762a1bSJed Brown } 804c4762a1bSJed Brown } 805c4762a1bSJed Brown 8069566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Grad, MAT_FINAL_ASSEMBLY)); 8079566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Grad, MAT_FINAL_ASSEMBLY)); 808c4762a1bSJed Brown 809c4762a1bSJed Brown /* Generate arithmetic averaging matrix Av */ 8109566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Av)); 8119566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Av, PETSC_DECIDE, ysubnlocal, m, n)); 8129566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Av)); 8139566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->Av, 2, NULL, 2, NULL)); 8149566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->Av, 2, NULL)); 8159566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->Av, &istart, &iend)); 816c4762a1bSJed Brown 817c4762a1bSJed Brown for (i = istart; i < iend; i++) { 818c4762a1bSJed Brown if (i < m / 3) { 819c4762a1bSJed Brown iblock = i / (user->mx - 1); 820c4762a1bSJed Brown j = iblock * user->mx + (i % (user->mx - 1)); 8219566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES)); 822c4762a1bSJed Brown j = j + 1; 8239566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES)); 824c4762a1bSJed Brown } 825c4762a1bSJed Brown if (i >= m / 3 && i < 2 * m / 3) { 826c4762a1bSJed Brown iblock = (i - m / 3) / (user->mx * (user->mx - 1)); 827c4762a1bSJed Brown j = iblock * user->mx * user->mx + ((i - m / 3) % (user->mx * (user->mx - 1))); 8289566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES)); 829c4762a1bSJed Brown j = j + user->mx; 8309566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES)); 831c4762a1bSJed Brown } 832c4762a1bSJed Brown if (i >= 2 * m / 3) { 833c4762a1bSJed Brown j = i - 2 * m / 3; 8349566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES)); 835c4762a1bSJed Brown j = j + user->mx * user->mx; 8369566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES)); 837c4762a1bSJed Brown } 838c4762a1bSJed Brown } 839c4762a1bSJed Brown 8409566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Av, MAT_FINAL_ASSEMBLY)); 8419566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Av, MAT_FINAL_ASSEMBLY)); 842c4762a1bSJed Brown 8439566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->L)); 8449566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->L, PETSC_DECIDE, ysubnlocal, m + n, n)); 8459566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->L)); 8469566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->L, 2, NULL, 2, NULL)); 8479566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->L, 2, NULL)); 8489566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->L, &istart, &iend)); 849c4762a1bSJed Brown 850c4762a1bSJed Brown for (i = istart; i < iend; i++) { 851c4762a1bSJed Brown if (i < m / 3) { 852c4762a1bSJed Brown iblock = i / (user->mx - 1); 853c4762a1bSJed Brown j = iblock * user->mx + (i % (user->mx - 1)); 8549566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES)); 855c4762a1bSJed Brown j = j + 1; 8569566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &hinv, INSERT_VALUES)); 857c4762a1bSJed Brown } 858c4762a1bSJed Brown if (i >= m / 3 && i < 2 * m / 3) { 859c4762a1bSJed Brown iblock = (i - m / 3) / (user->mx * (user->mx - 1)); 860c4762a1bSJed Brown j = iblock * user->mx * user->mx + ((i - m / 3) % (user->mx * (user->mx - 1))); 8619566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES)); 862c4762a1bSJed Brown j = j + user->mx; 8639566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &hinv, INSERT_VALUES)); 864c4762a1bSJed Brown } 865c4762a1bSJed Brown if (i >= 2 * m / 3 && i < m) { 866c4762a1bSJed Brown j = i - 2 * m / 3; 8679566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES)); 868c4762a1bSJed Brown j = j + user->mx * user->mx; 8699566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &hinv, INSERT_VALUES)); 870c4762a1bSJed Brown } 871c4762a1bSJed Brown if (i >= m) { 872c4762a1bSJed Brown j = i - m; 8739566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &sqrt_beta, INSERT_VALUES)); 874c4762a1bSJed Brown } 875c4762a1bSJed Brown } 8769566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->L, MAT_FINAL_ASSEMBLY)); 8779566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->L, MAT_FINAL_ASSEMBLY)); 8789566063dSJacob Faibussowitsch PetscCall(MatScale(user->L, PetscPowScalar(h, 1.5))); 879c4762a1bSJed Brown 880c4762a1bSJed Brown /* Generate Div matrix */ 881c4762a1bSJed Brown if (!user->use_ptap) { 882c4762a1bSJed Brown /* Generate Div matrix */ 8839566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Div)); 8849566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Div, ysubnlocal, PETSC_DECIDE, n, m)); 8859566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Div)); 8869566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->Div, 4, NULL, 4, NULL)); 8879566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->Div, 6, NULL)); 8889566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->Grad, &istart, &iend)); 889c4762a1bSJed Brown 890c4762a1bSJed Brown for (i = istart; i < iend; i++) { 891c4762a1bSJed Brown if (i < m / 3) { 892c4762a1bSJed Brown iblock = i / (user->mx - 1); 893c4762a1bSJed Brown j = iblock * user->mx + (i % (user->mx - 1)); 8949566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Div, 1, &j, 1, &i, &neg_hinv, INSERT_VALUES)); 895c4762a1bSJed Brown j = j + 1; 8969566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Div, 1, &j, 1, &i, &hinv, INSERT_VALUES)); 897c4762a1bSJed Brown } 898c4762a1bSJed Brown if (i >= m / 3 && i < 2 * m / 3) { 899c4762a1bSJed Brown iblock = (i - m / 3) / (user->mx * (user->mx - 1)); 900c4762a1bSJed Brown j = iblock * user->mx * user->mx + ((i - m / 3) % (user->mx * (user->mx - 1))); 9019566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Div, 1, &j, 1, &i, &neg_hinv, INSERT_VALUES)); 902c4762a1bSJed Brown j = j + user->mx; 9039566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Div, 1, &j, 1, &i, &hinv, INSERT_VALUES)); 904c4762a1bSJed Brown } 905c4762a1bSJed Brown if (i >= 2 * m / 3) { 906c4762a1bSJed Brown j = i - 2 * m / 3; 9079566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Div, 1, &j, 1, &i, &neg_hinv, INSERT_VALUES)); 908c4762a1bSJed Brown j = j + user->mx * user->mx; 9099566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Div, 1, &j, 1, &i, &hinv, INSERT_VALUES)); 910c4762a1bSJed Brown } 911c4762a1bSJed Brown } 912c4762a1bSJed Brown 9139566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Div, MAT_FINAL_ASSEMBLY)); 9149566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Div, MAT_FINAL_ASSEMBLY)); 9159566063dSJacob Faibussowitsch PetscCall(MatDuplicate(user->Div, MAT_SHARE_NONZERO_PATTERN, &user->Divwork)); 916c4762a1bSJed Brown } else { 9179566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Diag)); 9189566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Diag, PETSC_DECIDE, PETSC_DECIDE, m, m)); 9199566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Diag)); 9209566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->Diag, 1, NULL, 0, NULL)); 9219566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->Diag, 1, NULL)); 922c4762a1bSJed Brown } 923c4762a1bSJed Brown 924c4762a1bSJed Brown /* Build work vectors and matrices */ 9259566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->S)); 9269566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->S, PETSC_DECIDE, m)); 9279566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->S)); 928c4762a1bSJed Brown 9299566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->lwork)); 9309566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->lwork, PETSC_DECIDE, m + user->mx * user->mx * user->mx)); 9319566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->lwork)); 932c4762a1bSJed Brown 9339566063dSJacob Faibussowitsch PetscCall(MatDuplicate(user->Av, MAT_SHARE_NONZERO_PATTERN, &user->Avwork)); 934c4762a1bSJed Brown 9359566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->S, &user->Swork)); 9369566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->S, &user->Sdiag)); 9379566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->S, &user->Av_u)); 9389566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->S, &user->Twork)); 9399566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->y, &user->ywork)); 9409566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->u, &user->Ywork)); 9419566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->u, &user->uwork)); 9429566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->u, &user->js_diag)); 9439566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->c, &user->cwork)); 9449566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->d, &user->dwork)); 945c4762a1bSJed Brown 946c4762a1bSJed Brown /* Create a matrix-free shell user->Jd for computing B*x */ 9479566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, ysubnlocal * user->ns, ysubnlocal, user->nstate, user->ndesign, user, &user->Jd)); 9489566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->Jd, MATOP_MULT, (void (*)(void))DesignMatMult)); 9499566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->Jd, MATOP_MULT_TRANSPOSE, (void (*)(void))DesignMatMultTranspose)); 950c4762a1bSJed Brown 951c4762a1bSJed Brown /* Compute true state function ytrue given utrue */ 9529566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->y, &user->ytrue)); 953c4762a1bSJed Brown 954c4762a1bSJed Brown /* First compute Av_u = Av*exp(-u) */ 9559566063dSJacob Faibussowitsch PetscCall(VecSet(user->uwork, 0)); 9569566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->uwork, -1.0, user->utrue)); /* Note: user->utrue */ 9579566063dSJacob Faibussowitsch PetscCall(VecExp(user->uwork)); 9589566063dSJacob Faibussowitsch PetscCall(MatMult(user->Av, user->uwork, user->Av_u)); 959c4762a1bSJed Brown 960c4762a1bSJed Brown /* Next form DSG = Div*S*Grad */ 9619566063dSJacob Faibussowitsch PetscCall(VecCopy(user->Av_u, user->Swork)); 9629566063dSJacob Faibussowitsch PetscCall(VecReciprocal(user->Swork)); 963c4762a1bSJed Brown if (user->use_ptap) { 9649566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(user->Diag, user->Swork, INSERT_VALUES)); 9659566063dSJacob Faibussowitsch PetscCall(MatPtAP(user->Diag, user->Grad, MAT_INITIAL_MATRIX, 1.0, &user->DSG)); 966c4762a1bSJed Brown } else { 9679566063dSJacob Faibussowitsch PetscCall(MatCopy(user->Div, user->Divwork, SAME_NONZERO_PATTERN)); 9689566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(user->Divwork, NULL, user->Swork)); 969c20d7725SJed Brown 9709566063dSJacob Faibussowitsch PetscCall(MatMatMult(user->Divwork, user->Grad, MAT_INITIAL_MATRIX, 1.0, &user->DSG)); 971c4762a1bSJed Brown } 972c4762a1bSJed Brown 9739566063dSJacob Faibussowitsch PetscCall(MatSetOption(user->DSG, MAT_SYMMETRIC, PETSC_TRUE)); 9749566063dSJacob Faibussowitsch PetscCall(MatSetOption(user->DSG, MAT_SYMMETRY_ETERNAL, PETSC_TRUE)); 975c4762a1bSJed Brown 976c4762a1bSJed Brown if (user->use_lrc == PETSC_TRUE) { 977c4762a1bSJed Brown v = PetscSqrtReal(1.0 / user->ndesign); 9789566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->ndesign, &user->ones)); 979c4762a1bSJed Brown 980*9371c9d4SSatish Balay for (i = 0; i < user->ndesign; i++) { user->ones[i] = v; } 9819566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD, ysubnlocal, PETSC_DECIDE, user->ndesign, 1, user->ones, &user->Ones)); 9829566063dSJacob Faibussowitsch PetscCall(MatCreateLRC(user->DSG, user->Ones, NULL, user->Ones, &user->JsBlock)); 9839566063dSJacob Faibussowitsch PetscCall(MatSetUp(user->JsBlock)); 984c4762a1bSJed Brown } else { 985c4762a1bSJed Brown /* Create matrix-free shell user->Js for computing (A + h^3*e*e^T)*x */ 9869566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, ysubnlocal, ysubnlocal, user->ndesign, user->ndesign, user, &user->JsBlock)); 9879566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->JsBlock, MATOP_MULT, (void (*)(void))StateBlockMatMult)); 9889566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->JsBlock, MATOP_MULT_TRANSPOSE, (void (*)(void))StateBlockMatMult)); 989c4762a1bSJed Brown } 9909566063dSJacob Faibussowitsch PetscCall(MatSetOption(user->JsBlock, MAT_SYMMETRIC, PETSC_TRUE)); 9919566063dSJacob Faibussowitsch PetscCall(MatSetOption(user->JsBlock, MAT_SYMMETRY_ETERNAL, PETSC_TRUE)); 9929566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, ysubnlocal * user->ns, ysubnlocal * user->ns, user->nstate, user->nstate, user, &user->Js)); 9939566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->Js, MATOP_MULT, (void (*)(void))StateMatMult)); 9949566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->Js, MATOP_MULT_TRANSPOSE, (void (*)(void))StateMatMult)); 9959566063dSJacob Faibussowitsch PetscCall(MatSetOption(user->Js, MAT_SYMMETRIC, PETSC_TRUE)); 9969566063dSJacob Faibussowitsch PetscCall(MatSetOption(user->Js, MAT_SYMMETRY_ETERNAL, PETSC_TRUE)); 997c4762a1bSJed Brown 9989566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, ysubnlocal * user->ns, ysubnlocal * user->ns, user->nstate, user->nstate, user, &user->JsInv)); 9999566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->JsInv, MATOP_MULT, (void (*)(void))StateInvMatMult)); 10009566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->JsInv, MATOP_MULT_TRANSPOSE, (void (*)(void))StateInvMatMult)); 10019566063dSJacob Faibussowitsch PetscCall(MatSetOption(user->JsInv, MAT_SYMMETRIC, PETSC_TRUE)); 10029566063dSJacob Faibussowitsch PetscCall(MatSetOption(user->JsInv, MAT_SYMMETRY_ETERNAL, PETSC_TRUE)); 1003c4762a1bSJed Brown 10049566063dSJacob Faibussowitsch PetscCall(MatSetOption(user->DSG, MAT_SYMMETRIC, PETSC_TRUE)); 10059566063dSJacob Faibussowitsch PetscCall(MatSetOption(user->DSG, MAT_SYMMETRY_ETERNAL, PETSC_TRUE)); 1006c4762a1bSJed Brown /* Now solve for ytrue */ 10079566063dSJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_WORLD, &user->solver)); 10089566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(user->solver)); 1009c4762a1bSJed Brown 10109566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(user->solver, user->JsBlock, user->DSG)); 1011c4762a1bSJed Brown 10129566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsInv, user->q, user->ytrue)); 1013c4762a1bSJed Brown /* First compute Av_u = Av*exp(-u) */ 10149566063dSJacob Faibussowitsch PetscCall(VecSet(user->uwork, 0)); 10159566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->uwork, -1.0, user->u)); /* Note: user->u */ 10169566063dSJacob Faibussowitsch PetscCall(VecExp(user->uwork)); 10179566063dSJacob Faibussowitsch PetscCall(MatMult(user->Av, user->uwork, user->Av_u)); 1018c4762a1bSJed Brown 1019c4762a1bSJed Brown /* Next update DSG = Div*S*Grad with user->u */ 10209566063dSJacob Faibussowitsch PetscCall(VecCopy(user->Av_u, user->Swork)); 10219566063dSJacob Faibussowitsch PetscCall(VecReciprocal(user->Swork)); 1022c4762a1bSJed Brown if (user->use_ptap) { 10239566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(user->Diag, user->Swork, INSERT_VALUES)); 10249566063dSJacob Faibussowitsch PetscCall(MatPtAP(user->Diag, user->Grad, MAT_REUSE_MATRIX, 1.0, &user->DSG)); 1025c4762a1bSJed Brown } else { 10269566063dSJacob Faibussowitsch PetscCall(MatCopy(user->Div, user->Divwork, SAME_NONZERO_PATTERN)); 10279566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(user->Divwork, NULL, user->Av_u)); 10289566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(user->DSG)); 1029c4762a1bSJed Brown } 1030c4762a1bSJed Brown 1031c4762a1bSJed Brown /* Now solve for y */ 1032c4762a1bSJed Brown 10339566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsInv, user->q, user->y)); 1034c4762a1bSJed Brown 1035c4762a1bSJed Brown user->ksp_its_initial = user->ksp_its; 1036c4762a1bSJed Brown user->ksp_its = 0; 1037c4762a1bSJed Brown /* Construct projection matrix Q (blocks) */ 10389566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Q)); 10399566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Q, dsubnlocal, ysubnlocal, user->ndata, user->ndesign)); 10409566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Q)); 10419566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->Q, 8, NULL, 8, NULL)); 10429566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->Q, 8, NULL)); 1043c4762a1bSJed Brown 1044c4762a1bSJed Brown for (i = 0; i < user->mx; i++) { 1045c4762a1bSJed Brown x[i] = h * (i + 0.5); 1046c4762a1bSJed Brown y[i] = h * (i + 0.5); 1047c4762a1bSJed Brown z[i] = h * (i + 0.5); 1048c4762a1bSJed Brown } 10499566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->Q, &istart, &iend)); 1050c4762a1bSJed Brown 1051*9371c9d4SSatish Balay nx = user->mx; 1052*9371c9d4SSatish Balay ny = user->mx; 1053*9371c9d4SSatish Balay nz = user->mx; 1054c4762a1bSJed Brown for (i = istart; i < iend; i++) { 1055c4762a1bSJed Brown xri = xr[i]; 1056c4762a1bSJed Brown im = 0; 1057c4762a1bSJed Brown xim = x[im]; 1058c4762a1bSJed Brown while (xri > xim && im < nx) { 1059c4762a1bSJed Brown im = im + 1; 1060c4762a1bSJed Brown xim = x[im]; 1061c4762a1bSJed Brown } 1062c4762a1bSJed Brown indx1 = im - 1; 1063c4762a1bSJed Brown indx2 = im; 1064c4762a1bSJed Brown dx1 = xri - x[indx1]; 1065c4762a1bSJed Brown dx2 = x[indx2] - xri; 1066c4762a1bSJed Brown 1067c4762a1bSJed Brown yri = yr[i]; 1068c4762a1bSJed Brown im = 0; 1069c4762a1bSJed Brown yim = y[im]; 1070c4762a1bSJed Brown while (yri > yim && im < ny) { 1071c4762a1bSJed Brown im = im + 1; 1072c4762a1bSJed Brown yim = y[im]; 1073c4762a1bSJed Brown } 1074c4762a1bSJed Brown indy1 = im - 1; 1075c4762a1bSJed Brown indy2 = im; 1076c4762a1bSJed Brown dy1 = yri - y[indy1]; 1077c4762a1bSJed Brown dy2 = y[indy2] - yri; 1078c4762a1bSJed Brown 1079c4762a1bSJed Brown zri = zr[i]; 1080c4762a1bSJed Brown im = 0; 1081c4762a1bSJed Brown zim = z[im]; 1082c4762a1bSJed Brown while (zri > zim && im < nz) { 1083c4762a1bSJed Brown im = im + 1; 1084c4762a1bSJed Brown zim = z[im]; 1085c4762a1bSJed Brown } 1086c4762a1bSJed Brown indz1 = im - 1; 1087c4762a1bSJed Brown indz2 = im; 1088c4762a1bSJed Brown dz1 = zri - z[indz1]; 1089c4762a1bSJed Brown dz2 = z[indz2] - zri; 1090c4762a1bSJed Brown 1091c4762a1bSJed Brown Dx = x[indx2] - x[indx1]; 1092c4762a1bSJed Brown Dy = y[indy2] - y[indy1]; 1093c4762a1bSJed Brown Dz = z[indz2] - z[indz1]; 1094c4762a1bSJed Brown 1095c4762a1bSJed Brown j = indx1 + indy1 * nx + indz1 * nx * ny; 1096c4762a1bSJed Brown v = (1 - dx1 / Dx) * (1 - dy1 / Dy) * (1 - dz1 / Dz); 10979566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Q, 1, &i, 1, &j, &v, INSERT_VALUES)); 1098c4762a1bSJed Brown 1099c4762a1bSJed Brown j = indx1 + indy1 * nx + indz2 * nx * ny; 1100c4762a1bSJed Brown v = (1 - dx1 / Dx) * (1 - dy1 / Dy) * (1 - dz2 / Dz); 11019566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Q, 1, &i, 1, &j, &v, INSERT_VALUES)); 1102c4762a1bSJed Brown 1103c4762a1bSJed Brown j = indx1 + indy2 * nx + indz1 * nx * ny; 1104c4762a1bSJed Brown v = (1 - dx1 / Dx) * (1 - dy2 / Dy) * (1 - dz1 / Dz); 11059566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Q, 1, &i, 1, &j, &v, INSERT_VALUES)); 1106c4762a1bSJed Brown 1107c4762a1bSJed Brown j = indx1 + indy2 * nx + indz2 * nx * ny; 1108c4762a1bSJed Brown v = (1 - dx1 / Dx) * (1 - dy2 / Dy) * (1 - dz2 / Dz); 11099566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Q, 1, &i, 1, &j, &v, INSERT_VALUES)); 1110c4762a1bSJed Brown 1111c4762a1bSJed Brown j = indx2 + indy1 * nx + indz1 * nx * ny; 1112c4762a1bSJed Brown v = (1 - dx2 / Dx) * (1 - dy1 / Dy) * (1 - dz1 / Dz); 11139566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Q, 1, &i, 1, &j, &v, INSERT_VALUES)); 1114c4762a1bSJed Brown 1115c4762a1bSJed Brown j = indx2 + indy1 * nx + indz2 * nx * ny; 1116c4762a1bSJed Brown v = (1 - dx2 / Dx) * (1 - dy1 / Dy) * (1 - dz2 / Dz); 11179566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Q, 1, &i, 1, &j, &v, INSERT_VALUES)); 1118c4762a1bSJed Brown 1119c4762a1bSJed Brown j = indx2 + indy2 * nx + indz1 * nx * ny; 1120c4762a1bSJed Brown v = (1 - dx2 / Dx) * (1 - dy2 / Dy) * (1 - dz1 / Dz); 11219566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Q, 1, &i, 1, &j, &v, INSERT_VALUES)); 1122c4762a1bSJed Brown 1123c4762a1bSJed Brown j = indx2 + indy2 * nx + indz2 * nx * ny; 1124c4762a1bSJed Brown v = (1 - dx2 / Dx) * (1 - dy2 / Dy) * (1 - dz2 / Dz); 11259566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Q, 1, &i, 1, &j, &v, INSERT_VALUES)); 1126c4762a1bSJed Brown } 1127c4762a1bSJed Brown 11289566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Q, MAT_FINAL_ASSEMBLY)); 11299566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Q, MAT_FINAL_ASSEMBLY)); 1130c4762a1bSJed Brown /* Create MQ (composed of blocks of Q */ 11319566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, dsubnlocal * user->ns, PETSC_DECIDE, user->ndata * user->ns, user->nstate, user, &user->MQ)); 11329566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->MQ, MATOP_MULT, (void (*)(void))QMatMult)); 11339566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->MQ, MATOP_MULT_TRANSPOSE, (void (*)(void))QMatMultTranspose)); 1134c4762a1bSJed Brown 1135c4762a1bSJed Brown /* Add noise to the measurement data */ 11369566063dSJacob Faibussowitsch PetscCall(VecSet(user->ywork, 1.0)); 11379566063dSJacob Faibussowitsch PetscCall(VecAYPX(user->ywork, user->noise, user->ytrue)); 11389566063dSJacob Faibussowitsch PetscCall(MatMult(user->MQ, user->ywork, user->d)); 1139c4762a1bSJed Brown 1140c4762a1bSJed Brown /* Now that initial conditions have been set, let the user pass tolerance options to the KSP solver */ 11419566063dSJacob Faibussowitsch PetscCall(PetscFree(x)); 11429566063dSJacob Faibussowitsch PetscCall(PetscFree(y)); 11439566063dSJacob Faibussowitsch PetscCall(PetscFree(z)); 11449566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 1145c4762a1bSJed Brown PetscFunctionReturn(0); 1146c4762a1bSJed Brown } 1147c4762a1bSJed Brown 1148*9371c9d4SSatish Balay PetscErrorCode EllipticDestroy(AppCtx *user) { 1149c4762a1bSJed Brown PetscInt i; 1150c4762a1bSJed Brown 1151c4762a1bSJed Brown PetscFunctionBegin; 11529566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->DSG)); 11539566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&user->solver)); 11549566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Q)); 11559566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->MQ)); 1156c4762a1bSJed Brown if (!user->use_ptap) { 11579566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Div)); 11589566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Divwork)); 1159c4762a1bSJed Brown } else { 11609566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Diag)); 1161c4762a1bSJed Brown } 1162*9371c9d4SSatish Balay if (user->use_lrc) { PetscCall(MatDestroy(&user->Ones)); } 1163c4762a1bSJed Brown 11649566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Grad)); 11659566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Av)); 11669566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Avwork)); 11679566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->L)); 11689566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Js)); 11699566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Jd)); 11709566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->JsBlock)); 11719566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->JsInv)); 1172c4762a1bSJed Brown 11739566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->x)); 11749566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->u)); 11759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->uwork)); 11769566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->utrue)); 11779566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->y)); 11789566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->ywork)); 11799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->ytrue)); 11809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->c)); 11819566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->cwork)); 11829566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->ur)); 11839566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->q)); 11849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->d)); 11859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->dwork)); 11869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->lwork)); 11879566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->S)); 11889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->Swork)); 11899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->Sdiag)); 11909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->Ywork)); 11919566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->Twork)); 11929566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->Av_u)); 11939566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->js_diag)); 11949566063dSJacob Faibussowitsch PetscCall(ISDestroy(&user->s_is)); 11959566063dSJacob Faibussowitsch PetscCall(ISDestroy(&user->d_is)); 11969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->suby)); 11979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->subd)); 11989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->subq)); 11999566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->state_scatter)); 12009566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->design_scatter)); 1201c4762a1bSJed Brown for (i = 0; i < user->ns; i++) { 12029566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->yi_scatter[i])); 12039566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->di_scatter[i])); 1204c4762a1bSJed Brown } 12059566063dSJacob Faibussowitsch PetscCall(PetscFree(user->yi_scatter)); 12069566063dSJacob Faibussowitsch PetscCall(PetscFree(user->di_scatter)); 1207c4762a1bSJed Brown if (user->use_lrc) { 12089566063dSJacob Faibussowitsch PetscCall(PetscFree(user->ones)); 12099566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Ones)); 1210c4762a1bSJed Brown } 1211c4762a1bSJed Brown PetscFunctionReturn(0); 1212c4762a1bSJed Brown } 1213c4762a1bSJed Brown 1214*9371c9d4SSatish Balay PetscErrorCode EllipticMonitor(Tao tao, void *ptr) { 1215c4762a1bSJed Brown Vec X; 1216c4762a1bSJed Brown PetscReal unorm, ynorm; 1217c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 1218c4762a1bSJed Brown 1219c4762a1bSJed Brown PetscFunctionBegin; 12209566063dSJacob Faibussowitsch PetscCall(TaoGetSolution(tao, &X)); 12219566063dSJacob Faibussowitsch PetscCall(Scatter(X, user->ywork, user->state_scatter, user->uwork, user->design_scatter)); 12229566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->ywork, -1.0, user->ytrue)); 12239566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->uwork, -1.0, user->utrue)); 12249566063dSJacob Faibussowitsch PetscCall(VecNorm(user->uwork, NORM_2, &unorm)); 12259566063dSJacob Faibussowitsch PetscCall(VecNorm(user->ywork, NORM_2, &ynorm)); 12269566063dSJacob Faibussowitsch PetscCall(PetscPrintf(MPI_COMM_WORLD, "||u-ut||=%g ||y-yt||=%g\n", (double)unorm, (double)ynorm)); 1227c4762a1bSJed Brown PetscFunctionReturn(0); 1228c4762a1bSJed Brown } 1229c4762a1bSJed Brown 1230c4762a1bSJed Brown /*TEST 1231c4762a1bSJed Brown 1232c4762a1bSJed Brown build: 1233c4762a1bSJed Brown requires: !complex 1234c4762a1bSJed Brown 1235c4762a1bSJed Brown test: 1236c4762a1bSJed Brown args: -tao_cmonitor -ns 1 -tao_type lcl -tao_gatol 1.e-3 -tao_max_it 11 1237c4762a1bSJed Brown requires: !single 1238c4762a1bSJed Brown 1239c4762a1bSJed Brown test: 1240c4762a1bSJed Brown suffix: 2 1241c4762a1bSJed Brown args: -tao_cmonitor -tao_type lcl -tao_max_it 11 -use_ptap -use_lrc -ns 1 -tao_gatol 1.e-3 1242c4762a1bSJed Brown requires: !single 1243c4762a1bSJed Brown 1244c4762a1bSJed Brown TEST*/ 1245