1 #include <petsc/private/taoimpl.h> 2 3 typedef struct { 4 PetscInt n; /* Number of variables */ 5 PetscInt m; /* Number of constraints per time step */ 6 PetscInt mx; /* grid points in each direction */ 7 PetscInt nt; /* Number of time steps; as of now, must be divisible by 8 */ 8 PetscInt ndata; /* Number of data points per sample */ 9 PetscInt ns; /* Number of samples */ 10 PetscInt *sample_times; /* Times of samples */ 11 IS s_is; 12 IS d_is; 13 14 VecScatter state_scatter; 15 VecScatter design_scatter; 16 VecScatter *yi_scatter; 17 VecScatter *di_scatter; 18 19 Mat Js, Jd, JsBlockPrec, JsInv, JsBlock; 20 PetscBool jformed, dsg_formed; 21 22 PetscReal alpha; /* Regularization parameter */ 23 PetscReal beta; /* Weight attributed to ||u||^2 in regularization functional */ 24 PetscReal noise; /* Amount of noise to add to data */ 25 PetscReal ht; /* Time step */ 26 27 Mat Qblock, QblockT; 28 Mat L, LT; 29 Mat Div, Divwork; 30 Mat Grad; 31 Mat Av, Avwork, AvT; 32 Mat DSG; 33 Vec q; 34 Vec ur; /* reference */ 35 36 Vec d; 37 Vec dwork; 38 Vec *di; 39 40 Vec y; /* state variables */ 41 Vec ywork; 42 43 Vec ytrue; 44 Vec *yi, *yiwork; 45 46 Vec u; /* design variables */ 47 Vec uwork; 48 49 Vec utrue; 50 Vec js_diag; 51 Vec c; /* constraint vector */ 52 Vec cwork; 53 54 Vec lwork; 55 Vec S; 56 Vec Rwork, Swork, Twork; 57 Vec Av_u; 58 59 KSP solver; 60 PC prec; 61 62 PetscInt ksp_its; 63 PetscInt ksp_its_initial; 64 } AppCtx; 65 66 PetscErrorCode FormFunction(Tao, Vec, PetscReal *, void *); 67 PetscErrorCode FormGradient(Tao, Vec, Vec, void *); 68 PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *); 69 PetscErrorCode FormJacobianState(Tao, Vec, Mat, Mat, Mat, void *); 70 PetscErrorCode FormJacobianDesign(Tao, Vec, Mat, void *); 71 PetscErrorCode FormConstraints(Tao, Vec, Vec, void *); 72 PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void *); 73 PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat); 74 PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat); 75 PetscErrorCode ParabolicInitialize(AppCtx *user); 76 PetscErrorCode ParabolicDestroy(AppCtx *user); 77 PetscErrorCode ParabolicMonitor(Tao, void *); 78 79 PetscErrorCode StateMatMult(Mat, Vec, Vec); 80 PetscErrorCode StateMatBlockMult(Mat, Vec, Vec); 81 PetscErrorCode StateMatMultTranspose(Mat, Vec, Vec); 82 PetscErrorCode StateMatGetDiagonal(Mat, Vec); 83 PetscErrorCode StateMatDuplicate(Mat, MatDuplicateOption, Mat *); 84 PetscErrorCode StateMatInvMult(Mat, Vec, Vec); 85 PetscErrorCode StateMatInvTransposeMult(Mat, Vec, Vec); 86 PetscErrorCode StateMatBlockPrecMult(PC, Vec, Vec); 87 88 PetscErrorCode DesignMatMult(Mat, Vec, Vec); 89 PetscErrorCode DesignMatMultTranspose(Mat, Vec, Vec); 90 91 PetscErrorCode Gather_i(Vec, Vec *, VecScatter *, PetscInt); 92 PetscErrorCode Scatter_i(Vec, Vec *, VecScatter *, PetscInt); 93 94 static char help[] = ""; 95 96 int main(int argc, char **argv) { 97 Vec x, x0; 98 Tao tao; 99 AppCtx user; 100 IS is_allstate, is_alldesign; 101 PetscInt lo, hi, hi2, lo2, ksp_old; 102 PetscInt ntests = 1; 103 PetscInt i; 104 #if defined(PETSC_USE_LOG) 105 PetscLogStage stages[1]; 106 #endif 107 108 PetscFunctionBeginUser; 109 PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 110 user.mx = 8; 111 PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "parabolic example", NULL); 112 PetscCall(PetscOptionsInt("-mx", "Number of grid points in each direction", "", user.mx, &user.mx, NULL)); 113 user.nt = 8; 114 PetscCall(PetscOptionsInt("-nt", "Number of time steps", "", user.nt, &user.nt, NULL)); 115 user.ndata = 64; 116 PetscCall(PetscOptionsInt("-ndata", "Numbers of data points per sample", "", user.ndata, &user.ndata, NULL)); 117 user.ns = 8; 118 PetscCall(PetscOptionsInt("-ns", "Number of samples", "", user.ns, &user.ns, NULL)); 119 user.alpha = 1.0; 120 PetscCall(PetscOptionsReal("-alpha", "Regularization parameter", "", user.alpha, &user.alpha, NULL)); 121 user.beta = 0.01; 122 PetscCall(PetscOptionsReal("-beta", "Weight attributed to ||u||^2 in regularization functional", "", user.beta, &user.beta, NULL)); 123 user.noise = 0.01; 124 PetscCall(PetscOptionsReal("-noise", "Amount of noise to add to data", "", user.noise, &user.noise, NULL)); 125 PetscCall(PetscOptionsInt("-ntests", "Number of times to repeat TaoSolve", "", ntests, &ntests, NULL)); 126 PetscOptionsEnd(); 127 128 user.m = user.mx * user.mx * user.mx; /* number of constraints per time step */ 129 user.n = user.m * (user.nt + 1); /* number of variables */ 130 user.ht = (PetscReal)1 / user.nt; 131 132 PetscCall(VecCreate(PETSC_COMM_WORLD, &user.u)); 133 PetscCall(VecCreate(PETSC_COMM_WORLD, &user.y)); 134 PetscCall(VecCreate(PETSC_COMM_WORLD, &user.c)); 135 PetscCall(VecSetSizes(user.u, PETSC_DECIDE, user.n - user.m * user.nt)); 136 PetscCall(VecSetSizes(user.y, PETSC_DECIDE, user.m * user.nt)); 137 PetscCall(VecSetSizes(user.c, PETSC_DECIDE, user.m * user.nt)); 138 PetscCall(VecSetFromOptions(user.u)); 139 PetscCall(VecSetFromOptions(user.y)); 140 PetscCall(VecSetFromOptions(user.c)); 141 142 /* Create scatters for reduced spaces. 143 If the state vector y and design vector u are partitioned as 144 [y_1; y_2; ...; y_np] and [u_1; u_2; ...; u_np] (with np = # of processors), 145 then the solution vector x is organized as 146 [y_1; u_1; y_2; u_2; ...; y_np; u_np]. 147 The index sets user.s_is and user.d_is correspond to the indices of the 148 state and design variables owned by the current processor. 149 */ 150 PetscCall(VecCreate(PETSC_COMM_WORLD, &x)); 151 152 PetscCall(VecGetOwnershipRange(user.y, &lo, &hi)); 153 PetscCall(VecGetOwnershipRange(user.u, &lo2, &hi2)); 154 155 PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_allstate)); 156 PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo + lo2, 1, &user.s_is)); 157 158 PetscCall(ISCreateStride(PETSC_COMM_SELF, hi2 - lo2, lo2, 1, &is_alldesign)); 159 PetscCall(ISCreateStride(PETSC_COMM_SELF, hi2 - lo2, hi + lo2, 1, &user.d_is)); 160 161 PetscCall(VecSetSizes(x, hi - lo + hi2 - lo2, user.n)); 162 PetscCall(VecSetFromOptions(x)); 163 164 PetscCall(VecScatterCreate(x, user.s_is, user.y, is_allstate, &user.state_scatter)); 165 PetscCall(VecScatterCreate(x, user.d_is, user.u, is_alldesign, &user.design_scatter)); 166 PetscCall(ISDestroy(&is_alldesign)); 167 PetscCall(ISDestroy(&is_allstate)); 168 169 /* Create TAO solver and set desired solution method */ 170 PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao)); 171 PetscCall(TaoSetType(tao, TAOLCL)); 172 173 /* Set up initial vectors and matrices */ 174 PetscCall(ParabolicInitialize(&user)); 175 176 PetscCall(Gather(x, user.y, user.state_scatter, user.u, user.design_scatter)); 177 PetscCall(VecDuplicate(x, &x0)); 178 PetscCall(VecCopy(x, x0)); 179 180 /* Set solution vector with an initial guess */ 181 PetscCall(TaoSetSolution(tao, x)); 182 PetscCall(TaoSetObjective(tao, FormFunction, &user)); 183 PetscCall(TaoSetGradient(tao, NULL, FormGradient, &user)); 184 PetscCall(TaoSetConstraintsRoutine(tao, user.c, FormConstraints, &user)); 185 186 PetscCall(TaoSetJacobianStateRoutine(tao, user.Js, user.JsBlockPrec, user.JsInv, FormJacobianState, &user)); 187 PetscCall(TaoSetJacobianDesignRoutine(tao, user.Jd, FormJacobianDesign, &user)); 188 189 PetscCall(TaoSetFromOptions(tao)); 190 PetscCall(TaoSetStateDesignIS(tao, user.s_is, user.d_is)); 191 192 /* SOLVE THE APPLICATION */ 193 PetscCall(PetscLogStageRegister("Trials", &stages[0])); 194 PetscCall(PetscLogStagePush(stages[0])); 195 user.ksp_its_initial = user.ksp_its; 196 ksp_old = user.ksp_its; 197 for (i = 0; i < ntests; i++) { 198 PetscCall(TaoSolve(tao)); 199 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "KSP Iterations = %" PetscInt_FMT "\n", user.ksp_its - ksp_old)); 200 PetscCall(VecCopy(x0, x)); 201 PetscCall(TaoSetSolution(tao, x)); 202 } 203 PetscCall(PetscLogStagePop()); 204 PetscCall(PetscBarrier((PetscObject)x)); 205 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "KSP iterations within initialization: ")); 206 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "\n", user.ksp_its_initial)); 207 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Total KSP iterations over %" PetscInt_FMT " trial(s): ", ntests)); 208 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "\n", user.ksp_its)); 209 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "KSP iterations per trial: ")); 210 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "\n", (user.ksp_its - user.ksp_its_initial) / ntests)); 211 212 PetscCall(TaoDestroy(&tao)); 213 PetscCall(VecDestroy(&x)); 214 PetscCall(VecDestroy(&x0)); 215 PetscCall(ParabolicDestroy(&user)); 216 217 PetscCall(PetscFinalize()); 218 return 0; 219 } 220 /* ------------------------------------------------------------------- */ 221 /* 222 dwork = Qy - d 223 lwork = L*(u-ur) 224 f = 1/2 * (dwork.dork + alpha*lwork.lwork) 225 */ 226 PetscErrorCode FormFunction(Tao tao, Vec X, PetscReal *f, void *ptr) { 227 PetscReal d1 = 0, d2 = 0; 228 PetscInt i, j; 229 AppCtx *user = (AppCtx *)ptr; 230 231 PetscFunctionBegin; 232 PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter)); 233 PetscCall(Scatter_i(user->y, user->yi, user->yi_scatter, user->nt)); 234 for (j = 0; j < user->ns; j++) { 235 i = user->sample_times[j]; 236 PetscCall(MatMult(user->Qblock, user->yi[i], user->di[j])); 237 } 238 PetscCall(Gather_i(user->dwork, user->di, user->di_scatter, user->ns)); 239 PetscCall(VecAXPY(user->dwork, -1.0, user->d)); 240 PetscCall(VecDot(user->dwork, user->dwork, &d1)); 241 242 PetscCall(VecWAXPY(user->uwork, -1.0, user->ur, user->u)); 243 PetscCall(MatMult(user->L, user->uwork, user->lwork)); 244 PetscCall(VecDot(user->lwork, user->lwork, &d2)); 245 246 *f = 0.5 * (d1 + user->alpha * d2); 247 PetscFunctionReturn(0); 248 } 249 250 /* ------------------------------------------------------------------- */ 251 /* 252 state: g_s = Q' *(Qy - d) 253 design: g_d = alpha*L'*L*(u-ur) 254 */ 255 PetscErrorCode FormGradient(Tao tao, Vec X, Vec G, void *ptr) { 256 PetscInt i, j; 257 AppCtx *user = (AppCtx *)ptr; 258 259 PetscFunctionBegin; 260 PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter)); 261 PetscCall(Scatter_i(user->y, user->yi, user->yi_scatter, user->nt)); 262 for (j = 0; j < user->ns; j++) { 263 i = user->sample_times[j]; 264 PetscCall(MatMult(user->Qblock, user->yi[i], user->di[j])); 265 } 266 PetscCall(Gather_i(user->dwork, user->di, user->di_scatter, user->ns)); 267 PetscCall(VecAXPY(user->dwork, -1.0, user->d)); 268 PetscCall(Scatter_i(user->dwork, user->di, user->di_scatter, user->ns)); 269 PetscCall(VecSet(user->ywork, 0.0)); 270 PetscCall(Scatter_i(user->ywork, user->yiwork, user->yi_scatter, user->nt)); 271 for (j = 0; j < user->ns; j++) { 272 i = user->sample_times[j]; 273 PetscCall(MatMult(user->QblockT, user->di[j], user->yiwork[i])); 274 } 275 PetscCall(Gather_i(user->ywork, user->yiwork, user->yi_scatter, user->nt)); 276 277 PetscCall(VecWAXPY(user->uwork, -1.0, user->ur, user->u)); 278 PetscCall(MatMult(user->L, user->uwork, user->lwork)); 279 PetscCall(MatMult(user->LT, user->lwork, user->uwork)); 280 PetscCall(VecScale(user->uwork, user->alpha)); 281 PetscCall(Gather(G, user->ywork, user->state_scatter, user->uwork, user->design_scatter)); 282 PetscFunctionReturn(0); 283 } 284 285 PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr) { 286 PetscReal d1, d2; 287 PetscInt i, j; 288 AppCtx *user = (AppCtx *)ptr; 289 290 PetscFunctionBegin; 291 PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter)); 292 PetscCall(Scatter_i(user->y, user->yi, user->yi_scatter, user->nt)); 293 for (j = 0; j < user->ns; j++) { 294 i = user->sample_times[j]; 295 PetscCall(MatMult(user->Qblock, user->yi[i], user->di[j])); 296 } 297 PetscCall(Gather_i(user->dwork, user->di, user->di_scatter, user->ns)); 298 PetscCall(VecAXPY(user->dwork, -1.0, user->d)); 299 PetscCall(VecDot(user->dwork, user->dwork, &d1)); 300 PetscCall(Scatter_i(user->dwork, user->di, user->di_scatter, user->ns)); 301 PetscCall(VecSet(user->ywork, 0.0)); 302 PetscCall(Scatter_i(user->ywork, user->yiwork, user->yi_scatter, user->nt)); 303 for (j = 0; j < user->ns; j++) { 304 i = user->sample_times[j]; 305 PetscCall(MatMult(user->QblockT, user->di[j], user->yiwork[i])); 306 } 307 PetscCall(Gather_i(user->ywork, user->yiwork, user->yi_scatter, user->nt)); 308 309 PetscCall(VecWAXPY(user->uwork, -1.0, user->ur, user->u)); 310 PetscCall(MatMult(user->L, user->uwork, user->lwork)); 311 PetscCall(VecDot(user->lwork, user->lwork, &d2)); 312 PetscCall(MatMult(user->LT, user->lwork, user->uwork)); 313 PetscCall(VecScale(user->uwork, user->alpha)); 314 *f = 0.5 * (d1 + user->alpha * d2); 315 316 PetscCall(Gather(G, user->ywork, user->state_scatter, user->uwork, user->design_scatter)); 317 PetscFunctionReturn(0); 318 } 319 320 /* ------------------------------------------------------------------- */ 321 /* A 322 MatShell object 323 */ 324 PetscErrorCode FormJacobianState(Tao tao, Vec X, Mat J, Mat JPre, Mat JInv, void *ptr) { 325 AppCtx *user = (AppCtx *)ptr; 326 327 PetscFunctionBegin; 328 PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter)); 329 PetscCall(VecSet(user->uwork, 0)); 330 PetscCall(VecAXPY(user->uwork, -1.0, user->u)); 331 PetscCall(VecExp(user->uwork)); 332 PetscCall(MatMult(user->Av, user->uwork, user->Av_u)); 333 PetscCall(VecCopy(user->Av_u, user->Swork)); 334 PetscCall(VecReciprocal(user->Swork)); 335 PetscCall(MatCopy(user->Div, user->Divwork, SAME_NONZERO_PATTERN)); 336 PetscCall(MatDiagonalScale(user->Divwork, NULL, user->Swork)); 337 if (user->dsg_formed) { 338 PetscCall(MatProductNumeric(user->DSG)); 339 } else { 340 PetscCall(MatMatMult(user->Divwork, user->Grad, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &user->DSG)); 341 user->dsg_formed = PETSC_TRUE; 342 } 343 344 /* B = speye(nx^3) + ht*DSG; */ 345 PetscCall(MatScale(user->DSG, user->ht)); 346 PetscCall(MatShift(user->DSG, 1.0)); 347 PetscFunctionReturn(0); 348 } 349 350 /* ------------------------------------------------------------------- */ 351 /* B */ 352 PetscErrorCode FormJacobianDesign(Tao tao, Vec X, Mat J, void *ptr) { 353 AppCtx *user = (AppCtx *)ptr; 354 355 PetscFunctionBegin; 356 PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter)); 357 PetscFunctionReturn(0); 358 } 359 360 PetscErrorCode StateMatMult(Mat J_shell, Vec X, Vec Y) { 361 PetscInt i; 362 AppCtx *user; 363 364 PetscFunctionBegin; 365 PetscCall(MatShellGetContext(J_shell, &user)); 366 PetscCall(Scatter_i(X, user->yi, user->yi_scatter, user->nt)); 367 PetscCall(MatMult(user->JsBlock, user->yi[0], user->yiwork[0])); 368 for (i = 1; i < user->nt; i++) { 369 PetscCall(MatMult(user->JsBlock, user->yi[i], user->yiwork[i])); 370 PetscCall(VecAXPY(user->yiwork[i], -1.0, user->yi[i - 1])); 371 } 372 PetscCall(Gather_i(Y, user->yiwork, user->yi_scatter, user->nt)); 373 PetscFunctionReturn(0); 374 } 375 376 PetscErrorCode StateMatMultTranspose(Mat J_shell, Vec X, Vec Y) { 377 PetscInt i; 378 AppCtx *user; 379 380 PetscFunctionBegin; 381 PetscCall(MatShellGetContext(J_shell, &user)); 382 PetscCall(Scatter_i(X, user->yi, user->yi_scatter, user->nt)); 383 for (i = 0; i < user->nt - 1; i++) { 384 PetscCall(MatMult(user->JsBlock, user->yi[i], user->yiwork[i])); 385 PetscCall(VecAXPY(user->yiwork[i], -1.0, user->yi[i + 1])); 386 } 387 i = user->nt - 1; 388 PetscCall(MatMult(user->JsBlock, user->yi[i], user->yiwork[i])); 389 PetscCall(Gather_i(Y, user->yiwork, user->yi_scatter, user->nt)); 390 PetscFunctionReturn(0); 391 } 392 393 PetscErrorCode StateMatBlockMult(Mat J_shell, Vec X, Vec Y) { 394 AppCtx *user; 395 396 PetscFunctionBegin; 397 PetscCall(MatShellGetContext(J_shell, &user)); 398 PetscCall(MatMult(user->Grad, X, user->Swork)); 399 PetscCall(VecPointwiseDivide(user->Swork, user->Swork, user->Av_u)); 400 PetscCall(MatMult(user->Div, user->Swork, Y)); 401 PetscCall(VecAYPX(Y, user->ht, X)); 402 PetscFunctionReturn(0); 403 } 404 405 PetscErrorCode DesignMatMult(Mat J_shell, Vec X, Vec Y) { 406 PetscInt i; 407 AppCtx *user; 408 409 PetscFunctionBegin; 410 PetscCall(MatShellGetContext(J_shell, &user)); 411 412 /* sdiag(1./v) */ 413 PetscCall(VecSet(user->uwork, 0)); 414 PetscCall(VecAXPY(user->uwork, -1.0, user->u)); 415 PetscCall(VecExp(user->uwork)); 416 417 /* sdiag(1./((Av*(1./v)).^2)) */ 418 PetscCall(MatMult(user->Av, user->uwork, user->Swork)); 419 PetscCall(VecPointwiseMult(user->Swork, user->Swork, user->Swork)); 420 PetscCall(VecReciprocal(user->Swork)); 421 422 /* (Av * (sdiag(1./v) * b)) */ 423 PetscCall(VecPointwiseMult(user->uwork, user->uwork, X)); 424 PetscCall(MatMult(user->Av, user->uwork, user->Twork)); 425 426 /* (sdiag(1./((Av*(1./v)).^2)) * (Av * (sdiag(1./v) * b))) */ 427 PetscCall(VecPointwiseMult(user->Swork, user->Twork, user->Swork)); 428 429 PetscCall(Scatter_i(user->y, user->yi, user->yi_scatter, user->nt)); 430 for (i = 0; i < user->nt; i++) { 431 /* (sdiag(Grad*y(:,i)) */ 432 PetscCall(MatMult(user->Grad, user->yi[i], user->Twork)); 433 434 /* ht * Div * (sdiag(Grad*y(:,i)) * (sdiag(1./((Av*(1./v)).^2)) * (Av * (sdiag(1./v) * b)))) */ 435 PetscCall(VecPointwiseMult(user->Twork, user->Twork, user->Swork)); 436 PetscCall(MatMult(user->Div, user->Twork, user->yiwork[i])); 437 PetscCall(VecScale(user->yiwork[i], user->ht)); 438 } 439 PetscCall(Gather_i(Y, user->yiwork, user->yi_scatter, user->nt)); 440 441 PetscFunctionReturn(0); 442 } 443 444 PetscErrorCode DesignMatMultTranspose(Mat J_shell, Vec X, Vec Y) { 445 PetscInt i; 446 AppCtx *user; 447 448 PetscFunctionBegin; 449 PetscCall(MatShellGetContext(J_shell, &user)); 450 451 /* sdiag(1./((Av*(1./v)).^2)) */ 452 PetscCall(VecSet(user->uwork, 0)); 453 PetscCall(VecAXPY(user->uwork, -1.0, user->u)); 454 PetscCall(VecExp(user->uwork)); 455 PetscCall(MatMult(user->Av, user->uwork, user->Rwork)); 456 PetscCall(VecPointwiseMult(user->Rwork, user->Rwork, user->Rwork)); 457 PetscCall(VecReciprocal(user->Rwork)); 458 459 PetscCall(VecSet(Y, 0.0)); 460 PetscCall(Scatter_i(user->y, user->yi, user->yi_scatter, user->nt)); 461 PetscCall(Scatter_i(X, user->yiwork, user->yi_scatter, user->nt)); 462 for (i = 0; i < user->nt; i++) { 463 /* (Div' * b(:,i)) */ 464 PetscCall(MatMult(user->Grad, user->yiwork[i], user->Swork)); 465 466 /* sdiag(Grad*y(:,i)) */ 467 PetscCall(MatMult(user->Grad, user->yi[i], user->Twork)); 468 469 /* (sdiag(Grad*y(:,i)) * (Div' * b(:,i))) */ 470 PetscCall(VecPointwiseMult(user->Twork, user->Swork, user->Twork)); 471 472 /* (sdiag(1./((Av*(1./v)).^2)) * (sdiag(Grad*y(:,i)) * (Div' * b(:,i)))) */ 473 PetscCall(VecPointwiseMult(user->Twork, user->Rwork, user->Twork)); 474 475 /* (Av' * (sdiag(1./((Av*(1./v)).^2)) * (sdiag(Grad*y(:,i)) * (Div' * b(:,i))))) */ 476 PetscCall(MatMult(user->AvT, user->Twork, user->yiwork[i])); 477 478 /* sdiag(1./v) * (Av' * (sdiag(1./((Av*(1./v)).^2)) * (sdiag(Grad*y(:,i)) * (Div' * b(:,i))))) */ 479 PetscCall(VecPointwiseMult(user->yiwork[i], user->uwork, user->yiwork[i])); 480 PetscCall(VecAXPY(Y, user->ht, user->yiwork[i])); 481 } 482 PetscFunctionReturn(0); 483 } 484 485 PetscErrorCode StateMatBlockPrecMult(PC PC_shell, Vec X, Vec Y) { 486 AppCtx *user; 487 488 PetscFunctionBegin; 489 PetscCall(PCShellGetContext(PC_shell, &user)); 490 491 if (user->dsg_formed) { 492 PetscCall(MatSOR(user->DSG, X, 1.0, (MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_SYMMETRIC_SWEEP), 0.0, 1, 1, Y)); 493 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "DSG not formed"); 494 PetscFunctionReturn(0); 495 } 496 497 PetscErrorCode StateMatInvMult(Mat J_shell, Vec X, Vec Y) { 498 AppCtx *user; 499 PetscInt its, i; 500 501 PetscFunctionBegin; 502 PetscCall(MatShellGetContext(J_shell, &user)); 503 504 if (Y == user->ytrue) { 505 /* First solve is done with true solution to set up problem */ 506 PetscCall(KSPSetTolerances(user->solver, 1e-8, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 507 } else { 508 PetscCall(KSPSetTolerances(user->solver, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 509 } 510 511 PetscCall(Scatter_i(X, user->yi, user->yi_scatter, user->nt)); 512 PetscCall(KSPSolve(user->solver, user->yi[0], user->yiwork[0])); 513 PetscCall(KSPGetIterationNumber(user->solver, &its)); 514 user->ksp_its = user->ksp_its + its; 515 516 for (i = 1; i < user->nt; i++) { 517 PetscCall(VecAXPY(user->yi[i], 1.0, user->yiwork[i - 1])); 518 PetscCall(KSPSolve(user->solver, user->yi[i], user->yiwork[i])); 519 PetscCall(KSPGetIterationNumber(user->solver, &its)); 520 user->ksp_its = user->ksp_its + its; 521 } 522 PetscCall(Gather_i(Y, user->yiwork, user->yi_scatter, user->nt)); 523 PetscFunctionReturn(0); 524 } 525 526 PetscErrorCode StateMatInvTransposeMult(Mat J_shell, Vec X, Vec Y) { 527 AppCtx *user; 528 PetscInt its, i; 529 530 PetscFunctionBegin; 531 PetscCall(MatShellGetContext(J_shell, &user)); 532 533 PetscCall(Scatter_i(X, user->yi, user->yi_scatter, user->nt)); 534 535 i = user->nt - 1; 536 PetscCall(KSPSolve(user->solver, user->yi[i], user->yiwork[i])); 537 538 PetscCall(KSPGetIterationNumber(user->solver, &its)); 539 user->ksp_its = user->ksp_its + its; 540 541 for (i = user->nt - 2; i >= 0; i--) { 542 PetscCall(VecAXPY(user->yi[i], 1.0, user->yiwork[i + 1])); 543 PetscCall(KSPSolve(user->solver, user->yi[i], user->yiwork[i])); 544 545 PetscCall(KSPGetIterationNumber(user->solver, &its)); 546 user->ksp_its = user->ksp_its + its; 547 } 548 549 PetscCall(Gather_i(Y, user->yiwork, user->yi_scatter, user->nt)); 550 PetscFunctionReturn(0); 551 } 552 553 PetscErrorCode StateMatDuplicate(Mat J_shell, MatDuplicateOption opt, Mat *new_shell) { 554 AppCtx *user; 555 556 PetscFunctionBegin; 557 PetscCall(MatShellGetContext(J_shell, &user)); 558 559 PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m, user->m, user, new_shell)); 560 PetscCall(MatShellSetOperation(*new_shell, MATOP_MULT, (void (*)(void))StateMatMult)); 561 PetscCall(MatShellSetOperation(*new_shell, MATOP_DUPLICATE, (void (*)(void))StateMatDuplicate)); 562 PetscCall(MatShellSetOperation(*new_shell, MATOP_MULT_TRANSPOSE, (void (*)(void))StateMatMultTranspose)); 563 PetscCall(MatShellSetOperation(*new_shell, MATOP_GET_DIAGONAL, (void (*)(void))StateMatGetDiagonal)); 564 PetscFunctionReturn(0); 565 } 566 567 PetscErrorCode StateMatGetDiagonal(Mat J_shell, Vec X) { 568 AppCtx *user; 569 570 PetscFunctionBegin; 571 PetscCall(MatShellGetContext(J_shell, &user)); 572 PetscCall(VecCopy(user->js_diag, X)); 573 PetscFunctionReturn(0); 574 } 575 576 PetscErrorCode FormConstraints(Tao tao, Vec X, Vec C, void *ptr) { 577 /* con = Ay - q, A = [B 0 0 ... 0; 578 -I B 0 ... 0; 579 0 -I B ... 0; 580 ... ; 581 0 ... -I B] 582 B = ht * Div * Sigma * Grad + eye */ 583 PetscInt i; 584 AppCtx *user = (AppCtx *)ptr; 585 586 PetscFunctionBegin; 587 PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter)); 588 PetscCall(Scatter_i(user->y, user->yi, user->yi_scatter, user->nt)); 589 PetscCall(MatMult(user->JsBlock, user->yi[0], user->yiwork[0])); 590 for (i = 1; i < user->nt; i++) { 591 PetscCall(MatMult(user->JsBlock, user->yi[i], user->yiwork[i])); 592 PetscCall(VecAXPY(user->yiwork[i], -1.0, user->yi[i - 1])); 593 } 594 PetscCall(Gather_i(C, user->yiwork, user->yi_scatter, user->nt)); 595 PetscCall(VecAXPY(C, -1.0, user->q)); 596 PetscFunctionReturn(0); 597 } 598 599 PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat) { 600 PetscFunctionBegin; 601 PetscCall(VecScatterBegin(s_scat, x, state, INSERT_VALUES, SCATTER_FORWARD)); 602 PetscCall(VecScatterEnd(s_scat, x, state, INSERT_VALUES, SCATTER_FORWARD)); 603 PetscCall(VecScatterBegin(d_scat, x, design, INSERT_VALUES, SCATTER_FORWARD)); 604 PetscCall(VecScatterEnd(d_scat, x, design, INSERT_VALUES, SCATTER_FORWARD)); 605 PetscFunctionReturn(0); 606 } 607 608 PetscErrorCode Scatter_i(Vec y, Vec *yi, VecScatter *scat, PetscInt nt) { 609 PetscInt i; 610 611 PetscFunctionBegin; 612 for (i = 0; i < nt; i++) { 613 PetscCall(VecScatterBegin(scat[i], y, yi[i], INSERT_VALUES, SCATTER_FORWARD)); 614 PetscCall(VecScatterEnd(scat[i], y, yi[i], INSERT_VALUES, SCATTER_FORWARD)); 615 } 616 PetscFunctionReturn(0); 617 } 618 619 PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat) { 620 PetscFunctionBegin; 621 PetscCall(VecScatterBegin(s_scat, state, x, INSERT_VALUES, SCATTER_REVERSE)); 622 PetscCall(VecScatterEnd(s_scat, state, x, INSERT_VALUES, SCATTER_REVERSE)); 623 PetscCall(VecScatterBegin(d_scat, design, x, INSERT_VALUES, SCATTER_REVERSE)); 624 PetscCall(VecScatterEnd(d_scat, design, x, INSERT_VALUES, SCATTER_REVERSE)); 625 PetscFunctionReturn(0); 626 } 627 628 PetscErrorCode Gather_i(Vec y, Vec *yi, VecScatter *scat, PetscInt nt) { 629 PetscInt i; 630 631 PetscFunctionBegin; 632 for (i = 0; i < nt; i++) { 633 PetscCall(VecScatterBegin(scat[i], yi[i], y, INSERT_VALUES, SCATTER_REVERSE)); 634 PetscCall(VecScatterEnd(scat[i], yi[i], y, INSERT_VALUES, SCATTER_REVERSE)); 635 } 636 PetscFunctionReturn(0); 637 } 638 639 PetscErrorCode ParabolicInitialize(AppCtx *user) { 640 PetscInt m, n, i, j, k, linear_index, istart, iend, iblock, lo, hi, lo2, hi2; 641 Vec XX, YY, ZZ, XXwork, YYwork, ZZwork, UTwork, yi, di, bc; 642 PetscReal *x, *y, *z; 643 PetscReal h, stime; 644 PetscScalar hinv, neg_hinv, half = 0.5, sqrt_beta; 645 PetscInt im, indx1, indx2, indy1, indy2, indz1, indz2, nx, ny, nz; 646 PetscReal xri, yri, zri, xim, yim, zim, dx1, dx2, dy1, dy2, dz1, dz2, Dx, Dy, Dz; 647 PetscScalar v, vx, vy, vz; 648 IS is_from_y, is_to_yi, is_from_d, is_to_di; 649 /* Data locations */ 650 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, 651 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, 652 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}; 653 654 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, 655 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, 656 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}; 657 658 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, 659 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, 660 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}; 661 662 PetscFunctionBegin; 663 PetscCall(PetscMalloc1(user->mx, &x)); 664 PetscCall(PetscMalloc1(user->mx, &y)); 665 PetscCall(PetscMalloc1(user->mx, &z)); 666 user->jformed = PETSC_FALSE; 667 user->dsg_formed = PETSC_FALSE; 668 669 n = user->mx * user->mx * user->mx; 670 m = 3 * user->mx * user->mx * (user->mx - 1); 671 sqrt_beta = PetscSqrtScalar(user->beta); 672 673 user->ksp_its = 0; 674 user->ksp_its_initial = 0; 675 676 stime = (PetscReal)user->nt / user->ns; 677 PetscCall(PetscMalloc1(user->ns, &user->sample_times)); 678 for (i = 0; i < user->ns; i++) user->sample_times[i] = (PetscInt)(stime * ((PetscReal)i + 1.0) - 0.5); 679 680 PetscCall(VecCreate(PETSC_COMM_WORLD, &XX)); 681 PetscCall(VecCreate(PETSC_COMM_WORLD, &user->q)); 682 PetscCall(VecSetSizes(XX, PETSC_DECIDE, n)); 683 PetscCall(VecSetSizes(user->q, PETSC_DECIDE, n * user->nt)); 684 PetscCall(VecSetFromOptions(XX)); 685 PetscCall(VecSetFromOptions(user->q)); 686 687 PetscCall(VecDuplicate(XX, &YY)); 688 PetscCall(VecDuplicate(XX, &ZZ)); 689 PetscCall(VecDuplicate(XX, &XXwork)); 690 PetscCall(VecDuplicate(XX, &YYwork)); 691 PetscCall(VecDuplicate(XX, &ZZwork)); 692 PetscCall(VecDuplicate(XX, &UTwork)); 693 PetscCall(VecDuplicate(XX, &user->utrue)); 694 PetscCall(VecDuplicate(XX, &bc)); 695 696 /* Generate 3D grid, and collect ns (1<=ns<=8) right-hand-side vectors into user->q */ 697 h = 1.0 / user->mx; 698 hinv = user->mx; 699 neg_hinv = -hinv; 700 701 PetscCall(VecGetOwnershipRange(XX, &istart, &iend)); 702 for (linear_index = istart; linear_index < iend; linear_index++) { 703 i = linear_index % user->mx; 704 j = ((linear_index - i) / user->mx) % user->mx; 705 k = ((linear_index - i) / user->mx - j) / user->mx; 706 vx = h * (i + 0.5); 707 vy = h * (j + 0.5); 708 vz = h * (k + 0.5); 709 PetscCall(VecSetValues(XX, 1, &linear_index, &vx, INSERT_VALUES)); 710 PetscCall(VecSetValues(YY, 1, &linear_index, &vy, INSERT_VALUES)); 711 PetscCall(VecSetValues(ZZ, 1, &linear_index, &vz, INSERT_VALUES)); 712 if ((vx < 0.6) && (vx > 0.4) && (vy < 0.6) && (vy > 0.4) && (vy < 0.6) && (vz < 0.6) && (vz > 0.4)) { 713 v = 1000.0; 714 PetscCall(VecSetValues(bc, 1, &linear_index, &v, INSERT_VALUES)); 715 } 716 } 717 718 PetscCall(VecAssemblyBegin(XX)); 719 PetscCall(VecAssemblyEnd(XX)); 720 PetscCall(VecAssemblyBegin(YY)); 721 PetscCall(VecAssemblyEnd(YY)); 722 PetscCall(VecAssemblyBegin(ZZ)); 723 PetscCall(VecAssemblyEnd(ZZ)); 724 PetscCall(VecAssemblyBegin(bc)); 725 PetscCall(VecAssemblyEnd(bc)); 726 727 /* Compute true parameter function 728 ut = 0.5 + exp(-10*((x-0.5)^2+(y-0.5)^2+(z-0.5)^2)) */ 729 PetscCall(VecCopy(XX, XXwork)); 730 PetscCall(VecCopy(YY, YYwork)); 731 PetscCall(VecCopy(ZZ, ZZwork)); 732 733 PetscCall(VecShift(XXwork, -0.5)); 734 PetscCall(VecShift(YYwork, -0.5)); 735 PetscCall(VecShift(ZZwork, -0.5)); 736 737 PetscCall(VecPointwiseMult(XXwork, XXwork, XXwork)); 738 PetscCall(VecPointwiseMult(YYwork, YYwork, YYwork)); 739 PetscCall(VecPointwiseMult(ZZwork, ZZwork, ZZwork)); 740 741 PetscCall(VecCopy(XXwork, user->utrue)); 742 PetscCall(VecAXPY(user->utrue, 1.0, YYwork)); 743 PetscCall(VecAXPY(user->utrue, 1.0, ZZwork)); 744 PetscCall(VecScale(user->utrue, -10.0)); 745 PetscCall(VecExp(user->utrue)); 746 PetscCall(VecShift(user->utrue, 0.5)); 747 748 PetscCall(VecDestroy(&XX)); 749 PetscCall(VecDestroy(&YY)); 750 PetscCall(VecDestroy(&ZZ)); 751 PetscCall(VecDestroy(&XXwork)); 752 PetscCall(VecDestroy(&YYwork)); 753 PetscCall(VecDestroy(&ZZwork)); 754 PetscCall(VecDestroy(&UTwork)); 755 756 /* Initial guess and reference model */ 757 PetscCall(VecDuplicate(user->utrue, &user->ur)); 758 PetscCall(VecSet(user->ur, 0.0)); 759 760 /* Generate Grad matrix */ 761 PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Grad)); 762 PetscCall(MatSetSizes(user->Grad, PETSC_DECIDE, PETSC_DECIDE, m, n)); 763 PetscCall(MatSetFromOptions(user->Grad)); 764 PetscCall(MatMPIAIJSetPreallocation(user->Grad, 2, NULL, 2, NULL)); 765 PetscCall(MatSeqAIJSetPreallocation(user->Grad, 2, NULL)); 766 PetscCall(MatGetOwnershipRange(user->Grad, &istart, &iend)); 767 768 for (i = istart; i < iend; i++) { 769 if (i < m / 3) { 770 iblock = i / (user->mx - 1); 771 j = iblock * user->mx + (i % (user->mx - 1)); 772 PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES)); 773 j = j + 1; 774 PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &hinv, INSERT_VALUES)); 775 } 776 if (i >= m / 3 && i < 2 * m / 3) { 777 iblock = (i - m / 3) / (user->mx * (user->mx - 1)); 778 j = iblock * user->mx * user->mx + ((i - m / 3) % (user->mx * (user->mx - 1))); 779 PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES)); 780 j = j + user->mx; 781 PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &hinv, INSERT_VALUES)); 782 } 783 if (i >= 2 * m / 3) { 784 j = i - 2 * m / 3; 785 PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES)); 786 j = j + user->mx * user->mx; 787 PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &hinv, INSERT_VALUES)); 788 } 789 } 790 791 PetscCall(MatAssemblyBegin(user->Grad, MAT_FINAL_ASSEMBLY)); 792 PetscCall(MatAssemblyEnd(user->Grad, MAT_FINAL_ASSEMBLY)); 793 794 /* Generate arithmetic averaging matrix Av */ 795 PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Av)); 796 PetscCall(MatSetSizes(user->Av, PETSC_DECIDE, PETSC_DECIDE, m, n)); 797 PetscCall(MatSetFromOptions(user->Av)); 798 PetscCall(MatMPIAIJSetPreallocation(user->Av, 2, NULL, 2, NULL)); 799 PetscCall(MatSeqAIJSetPreallocation(user->Av, 2, NULL)); 800 PetscCall(MatGetOwnershipRange(user->Av, &istart, &iend)); 801 802 for (i = istart; i < iend; i++) { 803 if (i < m / 3) { 804 iblock = i / (user->mx - 1); 805 j = iblock * user->mx + (i % (user->mx - 1)); 806 PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES)); 807 j = j + 1; 808 PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES)); 809 } 810 if (i >= m / 3 && i < 2 * m / 3) { 811 iblock = (i - m / 3) / (user->mx * (user->mx - 1)); 812 j = iblock * user->mx * user->mx + ((i - m / 3) % (user->mx * (user->mx - 1))); 813 PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES)); 814 j = j + user->mx; 815 PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES)); 816 } 817 if (i >= 2 * m / 3) { 818 j = i - 2 * m / 3; 819 PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES)); 820 j = j + user->mx * user->mx; 821 PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES)); 822 } 823 } 824 825 PetscCall(MatAssemblyBegin(user->Av, MAT_FINAL_ASSEMBLY)); 826 PetscCall(MatAssemblyEnd(user->Av, MAT_FINAL_ASSEMBLY)); 827 828 /* Generate transpose of averaging matrix Av */ 829 PetscCall(MatTranspose(user->Av, MAT_INITIAL_MATRIX, &user->AvT)); 830 831 PetscCall(MatCreate(PETSC_COMM_WORLD, &user->L)); 832 PetscCall(MatSetSizes(user->L, PETSC_DECIDE, PETSC_DECIDE, m + n, n)); 833 PetscCall(MatSetFromOptions(user->L)); 834 PetscCall(MatMPIAIJSetPreallocation(user->L, 2, NULL, 2, NULL)); 835 PetscCall(MatSeqAIJSetPreallocation(user->L, 2, NULL)); 836 PetscCall(MatGetOwnershipRange(user->L, &istart, &iend)); 837 838 for (i = istart; i < iend; i++) { 839 if (i < m / 3) { 840 iblock = i / (user->mx - 1); 841 j = iblock * user->mx + (i % (user->mx - 1)); 842 PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES)); 843 j = j + 1; 844 PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &hinv, INSERT_VALUES)); 845 } 846 if (i >= m / 3 && i < 2 * m / 3) { 847 iblock = (i - m / 3) / (user->mx * (user->mx - 1)); 848 j = iblock * user->mx * user->mx + ((i - m / 3) % (user->mx * (user->mx - 1))); 849 PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES)); 850 j = j + user->mx; 851 PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &hinv, INSERT_VALUES)); 852 } 853 if (i >= 2 * m / 3 && i < m) { 854 j = i - 2 * m / 3; 855 PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES)); 856 j = j + user->mx * user->mx; 857 PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &hinv, INSERT_VALUES)); 858 } 859 if (i >= m) { 860 j = i - m; 861 PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &sqrt_beta, INSERT_VALUES)); 862 } 863 } 864 865 PetscCall(MatAssemblyBegin(user->L, MAT_FINAL_ASSEMBLY)); 866 PetscCall(MatAssemblyEnd(user->L, MAT_FINAL_ASSEMBLY)); 867 868 PetscCall(MatScale(user->L, PetscPowScalar(h, 1.5))); 869 870 /* Generate Div matrix */ 871 PetscCall(MatTranspose(user->Grad, MAT_INITIAL_MATRIX, &user->Div)); 872 873 /* Build work vectors and matrices */ 874 PetscCall(VecCreate(PETSC_COMM_WORLD, &user->S)); 875 PetscCall(VecSetSizes(user->S, PETSC_DECIDE, user->mx * user->mx * (user->mx - 1) * 3)); 876 PetscCall(VecSetFromOptions(user->S)); 877 878 PetscCall(VecCreate(PETSC_COMM_WORLD, &user->lwork)); 879 PetscCall(VecSetSizes(user->lwork, PETSC_DECIDE, m + user->mx * user->mx * user->mx)); 880 PetscCall(VecSetFromOptions(user->lwork)); 881 882 PetscCall(MatDuplicate(user->Div, MAT_SHARE_NONZERO_PATTERN, &user->Divwork)); 883 PetscCall(MatDuplicate(user->Av, MAT_SHARE_NONZERO_PATTERN, &user->Avwork)); 884 885 PetscCall(VecCreate(PETSC_COMM_WORLD, &user->d)); 886 PetscCall(VecSetSizes(user->d, PETSC_DECIDE, user->ndata * user->nt)); 887 PetscCall(VecSetFromOptions(user->d)); 888 889 PetscCall(VecDuplicate(user->S, &user->Swork)); 890 PetscCall(VecDuplicate(user->S, &user->Av_u)); 891 PetscCall(VecDuplicate(user->S, &user->Twork)); 892 PetscCall(VecDuplicate(user->S, &user->Rwork)); 893 PetscCall(VecDuplicate(user->y, &user->ywork)); 894 PetscCall(VecDuplicate(user->u, &user->uwork)); 895 PetscCall(VecDuplicate(user->u, &user->js_diag)); 896 PetscCall(VecDuplicate(user->c, &user->cwork)); 897 PetscCall(VecDuplicate(user->d, &user->dwork)); 898 899 /* Create matrix-free shell user->Js for computing A*x */ 900 PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m * user->nt, user->m * user->nt, user, &user->Js)); 901 PetscCall(MatShellSetOperation(user->Js, MATOP_MULT, (void (*)(void))StateMatMult)); 902 PetscCall(MatShellSetOperation(user->Js, MATOP_DUPLICATE, (void (*)(void))StateMatDuplicate)); 903 PetscCall(MatShellSetOperation(user->Js, MATOP_MULT_TRANSPOSE, (void (*)(void))StateMatMultTranspose)); 904 PetscCall(MatShellSetOperation(user->Js, MATOP_GET_DIAGONAL, (void (*)(void))StateMatGetDiagonal)); 905 906 /* Diagonal blocks of user->Js */ 907 PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m, user->m, user, &user->JsBlock)); 908 PetscCall(MatShellSetOperation(user->JsBlock, MATOP_MULT, (void (*)(void))StateMatBlockMult)); 909 /* Blocks are symmetric */ 910 PetscCall(MatShellSetOperation(user->JsBlock, MATOP_MULT_TRANSPOSE, (void (*)(void))StateMatBlockMult)); 911 912 /* Create a matrix-free shell user->JsBlockPrec for computing (U+D)\D*(L+D)\x, where JsBlock = L+D+U, 913 D is diagonal, L is strictly lower triangular, and U is strictly upper triangular. 914 This is an SSOR preconditioner for user->JsBlock. */ 915 PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m, user->m, user, &user->JsBlockPrec)); 916 PetscCall(MatShellSetOperation(user->JsBlockPrec, MATOP_MULT, (void (*)(void))StateMatBlockPrecMult)); 917 /* JsBlockPrec is symmetric */ 918 PetscCall(MatShellSetOperation(user->JsBlockPrec, MATOP_MULT_TRANSPOSE, (void (*)(void))StateMatBlockPrecMult)); 919 PetscCall(MatSetOption(user->JsBlockPrec, MAT_SYMMETRIC, PETSC_TRUE)); 920 PetscCall(MatSetOption(user->JsBlockPrec, MAT_SYMMETRY_ETERNAL, PETSC_TRUE)); 921 922 /* Create a matrix-free shell user->Jd for computing B*x */ 923 PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m * user->nt, user->m, user, &user->Jd)); 924 PetscCall(MatShellSetOperation(user->Jd, MATOP_MULT, (void (*)(void))DesignMatMult)); 925 PetscCall(MatShellSetOperation(user->Jd, MATOP_MULT_TRANSPOSE, (void (*)(void))DesignMatMultTranspose)); 926 927 /* User-defined routines for computing user->Js\x and user->Js^T\x*/ 928 PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, user->m * user->nt, user->m * user->nt, user, &user->JsInv)); 929 PetscCall(MatShellSetOperation(user->JsInv, MATOP_MULT, (void (*)(void))StateMatInvMult)); 930 PetscCall(MatShellSetOperation(user->JsInv, MATOP_MULT_TRANSPOSE, (void (*)(void))StateMatInvTransposeMult)); 931 932 /* Solver options and tolerances */ 933 PetscCall(KSPCreate(PETSC_COMM_WORLD, &user->solver)); 934 PetscCall(KSPSetType(user->solver, KSPCG)); 935 PetscCall(KSPSetOperators(user->solver, user->JsBlock, user->JsBlockPrec)); 936 PetscCall(KSPSetInitialGuessNonzero(user->solver, PETSC_FALSE)); 937 PetscCall(KSPSetTolerances(user->solver, 1e-4, 1e-20, 1e3, 500)); 938 PetscCall(KSPSetFromOptions(user->solver)); 939 PetscCall(KSPGetPC(user->solver, &user->prec)); 940 PetscCall(PCSetType(user->prec, PCSHELL)); 941 942 PetscCall(PCShellSetApply(user->prec, StateMatBlockPrecMult)); 943 PetscCall(PCShellSetApplyTranspose(user->prec, StateMatBlockPrecMult)); 944 PetscCall(PCShellSetContext(user->prec, user)); 945 946 /* Create scatter from y to y_1,y_2,...,y_nt */ 947 PetscCall(PetscMalloc1(user->nt * user->m, &user->yi_scatter)); 948 PetscCall(VecCreate(PETSC_COMM_WORLD, &yi)); 949 PetscCall(VecSetSizes(yi, PETSC_DECIDE, user->mx * user->mx * user->mx)); 950 PetscCall(VecSetFromOptions(yi)); 951 PetscCall(VecDuplicateVecs(yi, user->nt, &user->yi)); 952 PetscCall(VecDuplicateVecs(yi, user->nt, &user->yiwork)); 953 954 PetscCall(VecGetOwnershipRange(user->y, &lo2, &hi2)); 955 istart = 0; 956 for (i = 0; i < user->nt; i++) { 957 PetscCall(VecGetOwnershipRange(user->yi[i], &lo, &hi)); 958 PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_to_yi)); 959 PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo2 + istart, 1, &is_from_y)); 960 PetscCall(VecScatterCreate(user->y, is_from_y, user->yi[i], is_to_yi, &user->yi_scatter[i])); 961 istart = istart + hi - lo; 962 PetscCall(ISDestroy(&is_to_yi)); 963 PetscCall(ISDestroy(&is_from_y)); 964 } 965 PetscCall(VecDestroy(&yi)); 966 967 /* Create scatter from d to d_1,d_2,...,d_ns */ 968 PetscCall(PetscMalloc1(user->ns * user->ndata, &user->di_scatter)); 969 PetscCall(VecCreate(PETSC_COMM_WORLD, &di)); 970 PetscCall(VecSetSizes(di, PETSC_DECIDE, user->ndata)); 971 PetscCall(VecSetFromOptions(di)); 972 PetscCall(VecDuplicateVecs(di, user->ns, &user->di)); 973 PetscCall(VecGetOwnershipRange(user->d, &lo2, &hi2)); 974 istart = 0; 975 for (i = 0; i < user->ns; i++) { 976 PetscCall(VecGetOwnershipRange(user->di[i], &lo, &hi)); 977 PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_to_di)); 978 PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo2 + istart, 1, &is_from_d)); 979 PetscCall(VecScatterCreate(user->d, is_from_d, user->di[i], is_to_di, &user->di_scatter[i])); 980 istart = istart + hi - lo; 981 PetscCall(ISDestroy(&is_to_di)); 982 PetscCall(ISDestroy(&is_from_d)); 983 } 984 PetscCall(VecDestroy(&di)); 985 986 /* Assemble RHS of forward problem */ 987 PetscCall(VecCopy(bc, user->yiwork[0])); 988 for (i = 1; i < user->nt; i++) PetscCall(VecSet(user->yiwork[i], 0.0)); 989 PetscCall(Gather_i(user->q, user->yiwork, user->yi_scatter, user->nt)); 990 PetscCall(VecDestroy(&bc)); 991 992 /* Compute true state function yt given ut */ 993 PetscCall(VecCreate(PETSC_COMM_WORLD, &user->ytrue)); 994 PetscCall(VecSetSizes(user->ytrue, PETSC_DECIDE, n * user->nt)); 995 PetscCall(VecSetFromOptions(user->ytrue)); 996 997 /* First compute Av_u = Av*exp(-u) */ 998 PetscCall(VecSet(user->uwork, 0)); 999 PetscCall(VecAXPY(user->uwork, -1.0, user->utrue)); /* Note: user->utrue */ 1000 PetscCall(VecExp(user->uwork)); 1001 PetscCall(MatMult(user->Av, user->uwork, user->Av_u)); 1002 1003 /* Symbolic DSG = Div * Grad */ 1004 PetscCall(MatProductCreate(user->Div, user->Grad, NULL, &user->DSG)); 1005 PetscCall(MatProductSetType(user->DSG, MATPRODUCT_AB)); 1006 PetscCall(MatProductSetAlgorithm(user->DSG, "default")); 1007 PetscCall(MatProductSetFill(user->DSG, PETSC_DEFAULT)); 1008 PetscCall(MatProductSetFromOptions(user->DSG)); 1009 PetscCall(MatProductSymbolic(user->DSG)); 1010 1011 user->dsg_formed = PETSC_TRUE; 1012 1013 /* Next form DSG = Div*Grad */ 1014 PetscCall(MatCopy(user->Div, user->Divwork, SAME_NONZERO_PATTERN)); 1015 PetscCall(MatDiagonalScale(user->Divwork, NULL, user->Av_u)); 1016 if (user->dsg_formed) { 1017 PetscCall(MatProductNumeric(user->DSG)); 1018 } else { 1019 PetscCall(MatMatMult(user->Div, user->Grad, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &user->DSG)); 1020 user->dsg_formed = PETSC_TRUE; 1021 } 1022 /* B = speye(nx^3) + ht*DSG; */ 1023 PetscCall(MatScale(user->DSG, user->ht)); 1024 PetscCall(MatShift(user->DSG, 1.0)); 1025 1026 /* Now solve for ytrue */ 1027 PetscCall(StateMatInvMult(user->Js, user->q, user->ytrue)); 1028 1029 /* Initial guess y0 for state given u0 */ 1030 1031 /* First compute Av_u = Av*exp(-u) */ 1032 PetscCall(VecSet(user->uwork, 0)); 1033 PetscCall(VecAXPY(user->uwork, -1.0, user->u)); /* Note: user->u */ 1034 PetscCall(VecExp(user->uwork)); 1035 PetscCall(MatMult(user->Av, user->uwork, user->Av_u)); 1036 1037 /* Next form DSG = Div*S*Grad */ 1038 PetscCall(MatCopy(user->Div, user->Divwork, SAME_NONZERO_PATTERN)); 1039 PetscCall(MatDiagonalScale(user->Divwork, NULL, user->Av_u)); 1040 if (user->dsg_formed) { 1041 PetscCall(MatProductNumeric(user->DSG)); 1042 } else { 1043 PetscCall(MatMatMult(user->Div, user->Grad, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &user->DSG)); 1044 1045 user->dsg_formed = PETSC_TRUE; 1046 } 1047 /* B = speye(nx^3) + ht*DSG; */ 1048 PetscCall(MatScale(user->DSG, user->ht)); 1049 PetscCall(MatShift(user->DSG, 1.0)); 1050 1051 /* Now solve for y */ 1052 PetscCall(StateMatInvMult(user->Js, user->q, user->y)); 1053 1054 /* Construct projection matrix Q, a block diagonal matrix consisting of nt copies of Qblock along the diagonal */ 1055 PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Qblock)); 1056 PetscCall(MatSetSizes(user->Qblock, PETSC_DECIDE, PETSC_DECIDE, user->ndata, n)); 1057 PetscCall(MatSetFromOptions(user->Qblock)); 1058 PetscCall(MatMPIAIJSetPreallocation(user->Qblock, 8, NULL, 8, NULL)); 1059 PetscCall(MatSeqAIJSetPreallocation(user->Qblock, 8, NULL)); 1060 1061 for (i = 0; i < user->mx; i++) { 1062 x[i] = h * (i + 0.5); 1063 y[i] = h * (i + 0.5); 1064 z[i] = h * (i + 0.5); 1065 } 1066 1067 PetscCall(MatGetOwnershipRange(user->Qblock, &istart, &iend)); 1068 nx = user->mx; 1069 ny = user->mx; 1070 nz = user->mx; 1071 for (i = istart; i < iend; i++) { 1072 xri = xr[i]; 1073 im = 0; 1074 xim = x[im]; 1075 while (xri > xim && im < nx) { 1076 im = im + 1; 1077 xim = x[im]; 1078 } 1079 indx1 = im - 1; 1080 indx2 = im; 1081 dx1 = xri - x[indx1]; 1082 dx2 = x[indx2] - xri; 1083 1084 yri = yr[i]; 1085 im = 0; 1086 yim = y[im]; 1087 while (yri > yim && im < ny) { 1088 im = im + 1; 1089 yim = y[im]; 1090 } 1091 indy1 = im - 1; 1092 indy2 = im; 1093 dy1 = yri - y[indy1]; 1094 dy2 = y[indy2] - yri; 1095 1096 zri = zr[i]; 1097 im = 0; 1098 zim = z[im]; 1099 while (zri > zim && im < nz) { 1100 im = im + 1; 1101 zim = z[im]; 1102 } 1103 indz1 = im - 1; 1104 indz2 = im; 1105 dz1 = zri - z[indz1]; 1106 dz2 = z[indz2] - zri; 1107 1108 Dx = x[indx2] - x[indx1]; 1109 Dy = y[indy2] - y[indy1]; 1110 Dz = z[indz2] - z[indz1]; 1111 1112 j = indx1 + indy1 * nx + indz1 * nx * ny; 1113 v = (1 - dx1 / Dx) * (1 - dy1 / Dy) * (1 - dz1 / Dz); 1114 PetscCall(MatSetValues(user->Qblock, 1, &i, 1, &j, &v, INSERT_VALUES)); 1115 1116 j = indx1 + indy1 * nx + indz2 * nx * ny; 1117 v = (1 - dx1 / Dx) * (1 - dy1 / Dy) * (1 - dz2 / Dz); 1118 PetscCall(MatSetValues(user->Qblock, 1, &i, 1, &j, &v, INSERT_VALUES)); 1119 1120 j = indx1 + indy2 * nx + indz1 * nx * ny; 1121 v = (1 - dx1 / Dx) * (1 - dy2 / Dy) * (1 - dz1 / Dz); 1122 PetscCall(MatSetValues(user->Qblock, 1, &i, 1, &j, &v, INSERT_VALUES)); 1123 1124 j = indx1 + indy2 * nx + indz2 * nx * ny; 1125 v = (1 - dx1 / Dx) * (1 - dy2 / Dy) * (1 - dz2 / Dz); 1126 PetscCall(MatSetValues(user->Qblock, 1, &i, 1, &j, &v, INSERT_VALUES)); 1127 1128 j = indx2 + indy1 * nx + indz1 * nx * ny; 1129 v = (1 - dx2 / Dx) * (1 - dy1 / Dy) * (1 - dz1 / Dz); 1130 PetscCall(MatSetValues(user->Qblock, 1, &i, 1, &j, &v, INSERT_VALUES)); 1131 1132 j = indx2 + indy1 * nx + indz2 * nx * ny; 1133 v = (1 - dx2 / Dx) * (1 - dy1 / Dy) * (1 - dz2 / Dz); 1134 PetscCall(MatSetValues(user->Qblock, 1, &i, 1, &j, &v, INSERT_VALUES)); 1135 1136 j = indx2 + indy2 * nx + indz1 * nx * ny; 1137 v = (1 - dx2 / Dx) * (1 - dy2 / Dy) * (1 - dz1 / Dz); 1138 PetscCall(MatSetValues(user->Qblock, 1, &i, 1, &j, &v, INSERT_VALUES)); 1139 1140 j = indx2 + indy2 * nx + indz2 * nx * ny; 1141 v = (1 - dx2 / Dx) * (1 - dy2 / Dy) * (1 - dz2 / Dz); 1142 PetscCall(MatSetValues(user->Qblock, 1, &i, 1, &j, &v, INSERT_VALUES)); 1143 } 1144 PetscCall(MatAssemblyBegin(user->Qblock, MAT_FINAL_ASSEMBLY)); 1145 PetscCall(MatAssemblyEnd(user->Qblock, MAT_FINAL_ASSEMBLY)); 1146 1147 PetscCall(MatTranspose(user->Qblock, MAT_INITIAL_MATRIX, &user->QblockT)); 1148 PetscCall(MatTranspose(user->L, MAT_INITIAL_MATRIX, &user->LT)); 1149 1150 /* Add noise to the measurement data */ 1151 PetscCall(VecSet(user->ywork, 1.0)); 1152 PetscCall(VecAYPX(user->ywork, user->noise, user->ytrue)); 1153 PetscCall(Scatter_i(user->ywork, user->yiwork, user->yi_scatter, user->nt)); 1154 for (j = 0; j < user->ns; j++) { 1155 i = user->sample_times[j]; 1156 PetscCall(MatMult(user->Qblock, user->yiwork[i], user->di[j])); 1157 } 1158 PetscCall(Gather_i(user->d, user->di, user->di_scatter, user->ns)); 1159 1160 /* Now that initial conditions have been set, let the user pass tolerance options to the KSP solver */ 1161 PetscCall(KSPSetFromOptions(user->solver)); 1162 PetscCall(PetscFree(x)); 1163 PetscCall(PetscFree(y)); 1164 PetscCall(PetscFree(z)); 1165 PetscFunctionReturn(0); 1166 } 1167 1168 PetscErrorCode ParabolicDestroy(AppCtx *user) { 1169 PetscInt i; 1170 1171 PetscFunctionBegin; 1172 PetscCall(MatDestroy(&user->Qblock)); 1173 PetscCall(MatDestroy(&user->QblockT)); 1174 PetscCall(MatDestroy(&user->Div)); 1175 PetscCall(MatDestroy(&user->Divwork)); 1176 PetscCall(MatDestroy(&user->Grad)); 1177 PetscCall(MatDestroy(&user->Av)); 1178 PetscCall(MatDestroy(&user->Avwork)); 1179 PetscCall(MatDestroy(&user->AvT)); 1180 PetscCall(MatDestroy(&user->DSG)); 1181 PetscCall(MatDestroy(&user->L)); 1182 PetscCall(MatDestroy(&user->LT)); 1183 PetscCall(KSPDestroy(&user->solver)); 1184 PetscCall(MatDestroy(&user->Js)); 1185 PetscCall(MatDestroy(&user->Jd)); 1186 PetscCall(MatDestroy(&user->JsInv)); 1187 PetscCall(MatDestroy(&user->JsBlock)); 1188 PetscCall(MatDestroy(&user->JsBlockPrec)); 1189 PetscCall(VecDestroy(&user->u)); 1190 PetscCall(VecDestroy(&user->uwork)); 1191 PetscCall(VecDestroy(&user->utrue)); 1192 PetscCall(VecDestroy(&user->y)); 1193 PetscCall(VecDestroy(&user->ywork)); 1194 PetscCall(VecDestroy(&user->ytrue)); 1195 PetscCall(VecDestroyVecs(user->nt, &user->yi)); 1196 PetscCall(VecDestroyVecs(user->nt, &user->yiwork)); 1197 PetscCall(VecDestroyVecs(user->ns, &user->di)); 1198 PetscCall(PetscFree(user->yi)); 1199 PetscCall(PetscFree(user->yiwork)); 1200 PetscCall(PetscFree(user->di)); 1201 PetscCall(VecDestroy(&user->c)); 1202 PetscCall(VecDestroy(&user->cwork)); 1203 PetscCall(VecDestroy(&user->ur)); 1204 PetscCall(VecDestroy(&user->q)); 1205 PetscCall(VecDestroy(&user->d)); 1206 PetscCall(VecDestroy(&user->dwork)); 1207 PetscCall(VecDestroy(&user->lwork)); 1208 PetscCall(VecDestroy(&user->S)); 1209 PetscCall(VecDestroy(&user->Swork)); 1210 PetscCall(VecDestroy(&user->Av_u)); 1211 PetscCall(VecDestroy(&user->Twork)); 1212 PetscCall(VecDestroy(&user->Rwork)); 1213 PetscCall(VecDestroy(&user->js_diag)); 1214 PetscCall(ISDestroy(&user->s_is)); 1215 PetscCall(ISDestroy(&user->d_is)); 1216 PetscCall(VecScatterDestroy(&user->state_scatter)); 1217 PetscCall(VecScatterDestroy(&user->design_scatter)); 1218 for (i = 0; i < user->nt; i++) PetscCall(VecScatterDestroy(&user->yi_scatter[i])); 1219 for (i = 0; i < user->ns; i++) PetscCall(VecScatterDestroy(&user->di_scatter[i])); 1220 PetscCall(PetscFree(user->yi_scatter)); 1221 PetscCall(PetscFree(user->di_scatter)); 1222 PetscCall(PetscFree(user->sample_times)); 1223 PetscFunctionReturn(0); 1224 } 1225 1226 PetscErrorCode ParabolicMonitor(Tao tao, void *ptr) { 1227 Vec X; 1228 PetscReal unorm, ynorm; 1229 AppCtx *user = (AppCtx *)ptr; 1230 1231 PetscFunctionBegin; 1232 PetscCall(TaoGetSolution(tao, &X)); 1233 PetscCall(Scatter(X, user->ywork, user->state_scatter, user->uwork, user->design_scatter)); 1234 PetscCall(VecAXPY(user->ywork, -1.0, user->ytrue)); 1235 PetscCall(VecAXPY(user->uwork, -1.0, user->utrue)); 1236 PetscCall(VecNorm(user->uwork, NORM_2, &unorm)); 1237 PetscCall(VecNorm(user->ywork, NORM_2, &ynorm)); 1238 PetscCall(PetscPrintf(MPI_COMM_WORLD, "||u-ut||=%g ||y-yt||=%g\n", (double)unorm, (double)ynorm)); 1239 PetscFunctionReturn(0); 1240 } 1241 1242 /*TEST 1243 1244 build: 1245 requires: !complex 1246 1247 test: 1248 args: -tao_cmonitor -tao_type lcl -ns 1 -tao_gatol 1.e-4 -ksp_max_it 30 1249 requires: !single 1250 1251 TEST*/ 1252