1c4762a1bSJed Brown #include <petsc/private/taoimpl.h> 2c4762a1bSJed Brown 3c4762a1bSJed Brown typedef struct { 4c4762a1bSJed Brown PetscInt n; /* Number of variables */ 5c4762a1bSJed Brown PetscInt m; /* Number of constraints per time step */ 6c4762a1bSJed Brown PetscInt mx; /* grid points in each direction */ 7c4762a1bSJed Brown PetscInt nt; /* Number of time steps; as of now, must be divisible by 8 */ 8c4762a1bSJed Brown PetscInt ndata; /* Number of data points per sample */ 9c4762a1bSJed Brown PetscInt ns; /* Number of samples */ 10c4762a1bSJed Brown PetscInt *sample_times; /* Times of samples */ 11c4762a1bSJed Brown IS s_is; 12c4762a1bSJed Brown IS d_is; 13c4762a1bSJed Brown 14c4762a1bSJed Brown VecScatter state_scatter; 15c4762a1bSJed Brown VecScatter design_scatter; 16c4762a1bSJed Brown VecScatter *yi_scatter; 17c4762a1bSJed Brown VecScatter *di_scatter; 18c4762a1bSJed Brown 19c4762a1bSJed Brown Mat Js,Jd,JsBlockPrec,JsInv,JsBlock; 20c4762a1bSJed Brown PetscBool jformed,dsg_formed; 21c4762a1bSJed Brown 22c4762a1bSJed Brown PetscReal alpha; /* Regularization parameter */ 23c4762a1bSJed Brown PetscReal beta; /* Weight attributed to ||u||^2 in regularization functional */ 24c4762a1bSJed Brown PetscReal noise; /* Amount of noise to add to data */ 25c4762a1bSJed Brown PetscReal ht; /* Time step */ 26c4762a1bSJed Brown 27c4762a1bSJed Brown Mat Qblock,QblockT; 28c4762a1bSJed Brown Mat L,LT; 29c4762a1bSJed Brown Mat Div,Divwork; 30c4762a1bSJed Brown Mat Grad; 31c4762a1bSJed Brown Mat Av,Avwork,AvT; 32c4762a1bSJed Brown Mat DSG; 33c4762a1bSJed Brown Vec q; 34c4762a1bSJed Brown Vec ur; /* reference */ 35c4762a1bSJed Brown 36c4762a1bSJed Brown Vec d; 37c4762a1bSJed Brown Vec dwork; 38c4762a1bSJed Brown Vec *di; 39c4762a1bSJed Brown 40c4762a1bSJed Brown Vec y; /* state variables */ 41c4762a1bSJed Brown Vec ywork; 42c4762a1bSJed Brown 43c4762a1bSJed Brown Vec ytrue; 44c4762a1bSJed Brown Vec *yi,*yiwork; 45c4762a1bSJed Brown 46c4762a1bSJed Brown Vec u; /* design variables */ 47c4762a1bSJed Brown Vec uwork; 48c4762a1bSJed Brown 49c4762a1bSJed Brown Vec utrue; 50c4762a1bSJed Brown Vec js_diag; 51c4762a1bSJed Brown Vec c; /* constraint vector */ 52c4762a1bSJed Brown Vec cwork; 53c4762a1bSJed Brown 54c4762a1bSJed Brown Vec lwork; 55c4762a1bSJed Brown Vec S; 56c4762a1bSJed Brown Vec Rwork,Swork,Twork; 57c4762a1bSJed Brown Vec Av_u; 58c4762a1bSJed Brown 59c4762a1bSJed Brown KSP solver; 60c4762a1bSJed Brown PC prec; 61c4762a1bSJed Brown 62c4762a1bSJed Brown PetscInt ksp_its; 63c4762a1bSJed Brown PetscInt ksp_its_initial; 64c4762a1bSJed Brown } AppCtx; 65c4762a1bSJed Brown 66c4762a1bSJed Brown PetscErrorCode FormFunction(Tao, Vec, PetscReal*, void*); 67c4762a1bSJed Brown PetscErrorCode FormGradient(Tao, Vec, Vec, void*); 68c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal*, Vec, void*); 69c4762a1bSJed Brown PetscErrorCode FormJacobianState(Tao, Vec, Mat, Mat, Mat, void*); 70c4762a1bSJed Brown PetscErrorCode FormJacobianDesign(Tao, Vec, Mat, void*); 71c4762a1bSJed Brown PetscErrorCode FormConstraints(Tao, Vec, Vec, void*); 72c4762a1bSJed Brown PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void*); 73c4762a1bSJed Brown PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat); 74c4762a1bSJed Brown PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat); 75c4762a1bSJed Brown PetscErrorCode ParabolicInitialize(AppCtx *user); 76c4762a1bSJed Brown PetscErrorCode ParabolicDestroy(AppCtx *user); 77c4762a1bSJed Brown PetscErrorCode ParabolicMonitor(Tao, void*); 78c4762a1bSJed Brown 79c4762a1bSJed Brown PetscErrorCode StateMatMult(Mat,Vec,Vec); 80c4762a1bSJed Brown PetscErrorCode StateMatBlockMult(Mat,Vec,Vec); 81c4762a1bSJed Brown PetscErrorCode StateMatMultTranspose(Mat,Vec,Vec); 82c4762a1bSJed Brown PetscErrorCode StateMatGetDiagonal(Mat,Vec); 83c4762a1bSJed Brown PetscErrorCode StateMatDuplicate(Mat,MatDuplicateOption,Mat*); 84c4762a1bSJed Brown PetscErrorCode StateMatInvMult(Mat,Vec,Vec); 85c4762a1bSJed Brown PetscErrorCode StateMatInvTransposeMult(Mat,Vec,Vec); 86c4762a1bSJed Brown PetscErrorCode StateMatBlockPrecMult(PC,Vec,Vec); 87c4762a1bSJed Brown 88c4762a1bSJed Brown PetscErrorCode DesignMatMult(Mat,Vec,Vec); 89c4762a1bSJed Brown PetscErrorCode DesignMatMultTranspose(Mat,Vec,Vec); 90c4762a1bSJed Brown 91c4762a1bSJed Brown PetscErrorCode Gather_i(Vec,Vec*,VecScatter*,PetscInt); 92c4762a1bSJed Brown PetscErrorCode Scatter_i(Vec,Vec*,VecScatter*,PetscInt); 93c4762a1bSJed Brown 94c4762a1bSJed Brown static char help[]=""; 95c4762a1bSJed Brown 96c4762a1bSJed Brown int main(int argc, char **argv) 97c4762a1bSJed Brown { 98c4762a1bSJed Brown Vec x,x0; 99c4762a1bSJed Brown Tao tao; 100c4762a1bSJed Brown AppCtx user; 101c4762a1bSJed Brown IS is_allstate,is_alldesign; 102c4762a1bSJed Brown PetscInt lo,hi,hi2,lo2,ksp_old; 103c4762a1bSJed Brown PetscInt ntests = 1; 104c4762a1bSJed Brown PetscInt i; 105c4762a1bSJed Brown #if defined(PETSC_USE_LOG) 106c4762a1bSJed Brown PetscLogStage stages[1]; 107c4762a1bSJed Brown #endif 108c4762a1bSJed Brown 1099566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char*)0,help)); 110c4762a1bSJed Brown user.mx = 8; 111*d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"parabolic example",NULL); 1129566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mx","Number of grid points in each direction","",user.mx,&user.mx,NULL)); 113c4762a1bSJed Brown user.nt = 8; 1149566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-nt","Number of time steps","",user.nt,&user.nt,NULL)); 115c4762a1bSJed Brown user.ndata = 64; 1169566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ndata","Numbers of data points per sample","",user.ndata,&user.ndata,NULL)); 117c4762a1bSJed Brown user.ns = 8; 1189566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ns","Number of samples","",user.ns,&user.ns,NULL)); 119c4762a1bSJed Brown user.alpha = 1.0; 1209566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-alpha","Regularization parameter","",user.alpha,&user.alpha,NULL)); 121c4762a1bSJed Brown user.beta = 0.01; 1229566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-beta","Weight attributed to ||u||^2 in regularization functional","",user.beta,&user.beta,NULL)); 123c4762a1bSJed Brown user.noise = 0.01; 1249566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-noise","Amount of noise to add to data","",user.noise,&user.noise,NULL)); 1259566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ntests","Number of times to repeat TaoSolve","",ntests,&ntests,NULL)); 126*d0609cedSBarry Smith PetscOptionsEnd(); 127c4762a1bSJed Brown 128c4762a1bSJed Brown user.m = user.mx*user.mx*user.mx; /* number of constraints per time step */ 129c4762a1bSJed Brown user.n = user.m*(user.nt+1); /* number of variables */ 130c4762a1bSJed Brown user.ht = (PetscReal)1/user.nt; 131c4762a1bSJed Brown 1329566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&user.u)); 1339566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&user.y)); 1349566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&user.c)); 1359566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user.u,PETSC_DECIDE,user.n-user.m*user.nt)); 1369566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user.y,PETSC_DECIDE,user.m*user.nt)); 1379566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user.c,PETSC_DECIDE,user.m*user.nt)); 1389566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user.u)); 1399566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user.y)); 1409566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user.c)); 141c4762a1bSJed Brown 142c4762a1bSJed Brown /* Create scatters for reduced spaces. 143c4762a1bSJed Brown If the state vector y and design vector u are partitioned as 144c4762a1bSJed Brown [y_1; y_2; ...; y_np] and [u_1; u_2; ...; u_np] (with np = # of processors), 145c4762a1bSJed Brown then the solution vector x is organized as 146c4762a1bSJed Brown [y_1; u_1; y_2; u_2; ...; y_np; u_np]. 147c4762a1bSJed Brown The index sets user.s_is and user.d_is correspond to the indices of the 148c4762a1bSJed Brown state and design variables owned by the current processor. 149c4762a1bSJed Brown */ 1509566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&x)); 151c4762a1bSJed Brown 1529566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user.y,&lo,&hi)); 1539566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user.u,&lo2,&hi2)); 154c4762a1bSJed Brown 1559566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_allstate)); 1569566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+lo2,1,&user.s_is)); 157c4762a1bSJed Brown 1589566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,hi2-lo2,lo2,1,&is_alldesign)); 1599566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,hi2-lo2,hi+lo2,1,&user.d_is)); 160c4762a1bSJed Brown 1619566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x,hi-lo+hi2-lo2,user.n)); 1629566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x)); 163c4762a1bSJed Brown 1649566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(x,user.s_is,user.y,is_allstate,&user.state_scatter)); 1659566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(x,user.d_is,user.u,is_alldesign,&user.design_scatter)); 1669566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_alldesign)); 1679566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_allstate)); 168c4762a1bSJed Brown 169c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 1709566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD,&tao)); 1719566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao,TAOLCL)); 172c4762a1bSJed Brown 173c4762a1bSJed Brown /* Set up initial vectors and matrices */ 1749566063dSJacob Faibussowitsch PetscCall(ParabolicInitialize(&user)); 175c4762a1bSJed Brown 1769566063dSJacob Faibussowitsch PetscCall(Gather(x,user.y,user.state_scatter,user.u,user.design_scatter)); 1779566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&x0)); 1789566063dSJacob Faibussowitsch PetscCall(VecCopy(x,x0)); 179c4762a1bSJed Brown 180c4762a1bSJed Brown /* Set solution vector with an initial guess */ 1819566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao,x)); 1829566063dSJacob Faibussowitsch PetscCall(TaoSetObjective(tao, FormFunction, &user)); 1839566063dSJacob Faibussowitsch PetscCall(TaoSetGradient(tao, NULL, FormGradient, &user)); 1849566063dSJacob Faibussowitsch PetscCall(TaoSetConstraintsRoutine(tao, user.c, FormConstraints, &user)); 185c4762a1bSJed Brown 1869566063dSJacob Faibussowitsch PetscCall(TaoSetJacobianStateRoutine(tao, user.Js, user.JsBlockPrec, user.JsInv, FormJacobianState, &user)); 1879566063dSJacob Faibussowitsch PetscCall(TaoSetJacobianDesignRoutine(tao, user.Jd, FormJacobianDesign, &user)); 188c4762a1bSJed Brown 1899566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 1909566063dSJacob Faibussowitsch PetscCall(TaoSetStateDesignIS(tao,user.s_is,user.d_is)); 191c4762a1bSJed Brown 192c4762a1bSJed Brown /* SOLVE THE APPLICATION */ 1939566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("Trials",&stages[0])); 1949566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stages[0])); 195c4762a1bSJed Brown user.ksp_its_initial = user.ksp_its; 196c4762a1bSJed Brown ksp_old = user.ksp_its; 197c4762a1bSJed Brown for (i=0; i<ntests; i++) { 1989566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 1999566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"KSP Iterations = %D\n",user.ksp_its-ksp_old)); 2009566063dSJacob Faibussowitsch PetscCall(VecCopy(x0,x)); 2019566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao,x)); 202c4762a1bSJed Brown } 2039566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 2049566063dSJacob Faibussowitsch PetscCall(PetscBarrier((PetscObject)x)); 2059566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"KSP iterations within initialization: ")); 2069566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%D\n",user.ksp_its_initial)); 2079566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Total KSP iterations over %D trial(s): ",ntests)); 2089566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%D\n",user.ksp_its)); 2099566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"KSP iterations per trial: ")); 2109566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%D\n",(user.ksp_its-user.ksp_its_initial)/ntests)); 211c4762a1bSJed Brown 2129566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 2139566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 2149566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x0)); 2159566063dSJacob Faibussowitsch PetscCall(ParabolicDestroy(&user)); 216c4762a1bSJed Brown 2179566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 218b122ec5aSJacob Faibussowitsch return 0; 219c4762a1bSJed Brown } 220c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 221c4762a1bSJed Brown /* 222c4762a1bSJed Brown dwork = Qy - d 223c4762a1bSJed Brown lwork = L*(u-ur) 224c4762a1bSJed Brown f = 1/2 * (dwork.dork + alpha*lwork.lwork) 225c4762a1bSJed Brown */ 226c4762a1bSJed Brown PetscErrorCode FormFunction(Tao tao,Vec X,PetscReal *f,void *ptr) 227c4762a1bSJed Brown { 228c4762a1bSJed Brown PetscReal d1=0,d2=0; 229c4762a1bSJed Brown PetscInt i,j; 230c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 231c4762a1bSJed Brown 232c4762a1bSJed Brown PetscFunctionBegin; 2339566063dSJacob Faibussowitsch PetscCall(Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter)); 2349566063dSJacob Faibussowitsch PetscCall(Scatter_i(user->y,user->yi,user->yi_scatter,user->nt)); 235c4762a1bSJed Brown for (j=0; j<user->ns; j++) { 236c4762a1bSJed Brown i = user->sample_times[j]; 2379566063dSJacob Faibussowitsch PetscCall(MatMult(user->Qblock,user->yi[i],user->di[j])); 238c4762a1bSJed Brown } 2399566063dSJacob Faibussowitsch PetscCall(Gather_i(user->dwork,user->di,user->di_scatter,user->ns)); 2409566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->dwork,-1.0,user->d)); 2419566063dSJacob Faibussowitsch PetscCall(VecDot(user->dwork,user->dwork,&d1)); 242c4762a1bSJed Brown 2439566063dSJacob Faibussowitsch PetscCall(VecWAXPY(user->uwork,-1.0,user->ur,user->u)); 2449566063dSJacob Faibussowitsch PetscCall(MatMult(user->L,user->uwork,user->lwork)); 2459566063dSJacob Faibussowitsch PetscCall(VecDot(user->lwork,user->lwork,&d2)); 246c4762a1bSJed Brown 247c4762a1bSJed Brown *f = 0.5 * (d1 + user->alpha*d2); 248c4762a1bSJed Brown PetscFunctionReturn(0); 249c4762a1bSJed Brown } 250c4762a1bSJed Brown 251c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 252c4762a1bSJed Brown /* 253c4762a1bSJed Brown state: g_s = Q' *(Qy - d) 254c4762a1bSJed Brown design: g_d = alpha*L'*L*(u-ur) 255c4762a1bSJed Brown */ 256c4762a1bSJed Brown PetscErrorCode FormGradient(Tao tao,Vec X,Vec G,void *ptr) 257c4762a1bSJed Brown { 258c4762a1bSJed Brown PetscInt i,j; 259c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 260c4762a1bSJed Brown 261c4762a1bSJed Brown PetscFunctionBegin; 2629566063dSJacob Faibussowitsch PetscCall(Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter)); 2639566063dSJacob Faibussowitsch PetscCall(Scatter_i(user->y,user->yi,user->yi_scatter,user->nt)); 264c4762a1bSJed Brown for (j=0; j<user->ns; j++) { 265c4762a1bSJed Brown i = user->sample_times[j]; 2669566063dSJacob Faibussowitsch PetscCall(MatMult(user->Qblock,user->yi[i],user->di[j])); 267c4762a1bSJed Brown } 2689566063dSJacob Faibussowitsch PetscCall(Gather_i(user->dwork,user->di,user->di_scatter,user->ns)); 2699566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->dwork,-1.0,user->d)); 2709566063dSJacob Faibussowitsch PetscCall(Scatter_i(user->dwork,user->di,user->di_scatter,user->ns)); 2719566063dSJacob Faibussowitsch PetscCall(VecSet(user->ywork,0.0)); 2729566063dSJacob Faibussowitsch PetscCall(Scatter_i(user->ywork,user->yiwork,user->yi_scatter,user->nt)); 273c4762a1bSJed Brown for (j=0; j<user->ns; j++) { 274c4762a1bSJed Brown i = user->sample_times[j]; 2759566063dSJacob Faibussowitsch PetscCall(MatMult(user->QblockT,user->di[j],user->yiwork[i])); 276c4762a1bSJed Brown } 2779566063dSJacob Faibussowitsch PetscCall(Gather_i(user->ywork,user->yiwork,user->yi_scatter,user->nt)); 278c4762a1bSJed Brown 2799566063dSJacob Faibussowitsch PetscCall(VecWAXPY(user->uwork,-1.0,user->ur,user->u)); 2809566063dSJacob Faibussowitsch PetscCall(MatMult(user->L,user->uwork,user->lwork)); 2819566063dSJacob Faibussowitsch PetscCall(MatMult(user->LT,user->lwork,user->uwork)); 2829566063dSJacob Faibussowitsch PetscCall(VecScale(user->uwork, user->alpha)); 2839566063dSJacob Faibussowitsch PetscCall(Gather(G,user->ywork,user->state_scatter,user->uwork,user->design_scatter)); 284c4762a1bSJed Brown PetscFunctionReturn(0); 285c4762a1bSJed Brown } 286c4762a1bSJed Brown 287c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr) 288c4762a1bSJed Brown { 289c4762a1bSJed Brown PetscReal d1,d2; 290c4762a1bSJed Brown PetscInt i,j; 291c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 292c4762a1bSJed Brown 293c4762a1bSJed Brown PetscFunctionBegin; 2949566063dSJacob Faibussowitsch PetscCall(Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter)); 2959566063dSJacob Faibussowitsch PetscCall(Scatter_i(user->y,user->yi,user->yi_scatter,user->nt)); 296c4762a1bSJed Brown for (j=0; j<user->ns; j++) { 297c4762a1bSJed Brown i = user->sample_times[j]; 2989566063dSJacob Faibussowitsch PetscCall(MatMult(user->Qblock,user->yi[i],user->di[j])); 299c4762a1bSJed Brown } 3009566063dSJacob Faibussowitsch PetscCall(Gather_i(user->dwork,user->di,user->di_scatter,user->ns)); 3019566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->dwork,-1.0,user->d)); 3029566063dSJacob Faibussowitsch PetscCall(VecDot(user->dwork,user->dwork,&d1)); 3039566063dSJacob Faibussowitsch PetscCall(Scatter_i(user->dwork,user->di,user->di_scatter,user->ns)); 3049566063dSJacob Faibussowitsch PetscCall(VecSet(user->ywork,0.0)); 3059566063dSJacob Faibussowitsch PetscCall(Scatter_i(user->ywork,user->yiwork,user->yi_scatter,user->nt)); 306c4762a1bSJed Brown for (j=0; j<user->ns; j++) { 307c4762a1bSJed Brown i = user->sample_times[j]; 3089566063dSJacob Faibussowitsch PetscCall(MatMult(user->QblockT,user->di[j],user->yiwork[i])); 309c4762a1bSJed Brown } 3109566063dSJacob Faibussowitsch PetscCall(Gather_i(user->ywork,user->yiwork,user->yi_scatter,user->nt)); 311c4762a1bSJed Brown 3129566063dSJacob Faibussowitsch PetscCall(VecWAXPY(user->uwork,-1.0,user->ur,user->u)); 3139566063dSJacob Faibussowitsch PetscCall(MatMult(user->L,user->uwork,user->lwork)); 3149566063dSJacob Faibussowitsch PetscCall(VecDot(user->lwork,user->lwork,&d2)); 3159566063dSJacob Faibussowitsch PetscCall(MatMult(user->LT,user->lwork,user->uwork)); 3169566063dSJacob Faibussowitsch PetscCall(VecScale(user->uwork, user->alpha)); 317c4762a1bSJed Brown *f = 0.5 * (d1 + user->alpha*d2); 318c4762a1bSJed Brown 3199566063dSJacob Faibussowitsch PetscCall(Gather(G,user->ywork,user->state_scatter,user->uwork,user->design_scatter)); 320c4762a1bSJed Brown PetscFunctionReturn(0); 321c4762a1bSJed Brown } 322c4762a1bSJed Brown 323c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 324c4762a1bSJed Brown /* A 325c4762a1bSJed Brown MatShell object 326c4762a1bSJed Brown */ 327c4762a1bSJed Brown PetscErrorCode FormJacobianState(Tao tao, Vec X, Mat J, Mat JPre, Mat JInv, void *ptr) 328c4762a1bSJed Brown { 329c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 330c4762a1bSJed Brown 331c4762a1bSJed Brown PetscFunctionBegin; 3329566063dSJacob Faibussowitsch PetscCall(Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter)); 3339566063dSJacob Faibussowitsch PetscCall(VecSet(user->uwork,0)); 3349566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->uwork,-1.0,user->u)); 3359566063dSJacob Faibussowitsch PetscCall(VecExp(user->uwork)); 3369566063dSJacob Faibussowitsch PetscCall(MatMult(user->Av,user->uwork,user->Av_u)); 3379566063dSJacob Faibussowitsch PetscCall(VecCopy(user->Av_u,user->Swork)); 3389566063dSJacob Faibussowitsch PetscCall(VecReciprocal(user->Swork)); 3399566063dSJacob Faibussowitsch PetscCall(MatCopy(user->Div,user->Divwork,SAME_NONZERO_PATTERN)); 3409566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(user->Divwork,NULL,user->Swork)); 341c4762a1bSJed Brown if (user->dsg_formed) { 3429566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(user->DSG)); 343c4762a1bSJed Brown } else { 3449566063dSJacob Faibussowitsch PetscCall(MatMatMult(user->Divwork,user->Grad,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&user->DSG)); 345c4762a1bSJed Brown user->dsg_formed = PETSC_TRUE; 346c4762a1bSJed Brown } 347c4762a1bSJed Brown 348c4762a1bSJed Brown /* B = speye(nx^3) + ht*DSG; */ 3499566063dSJacob Faibussowitsch PetscCall(MatScale(user->DSG,user->ht)); 3509566063dSJacob Faibussowitsch PetscCall(MatShift(user->DSG,1.0)); 351c4762a1bSJed Brown PetscFunctionReturn(0); 352c4762a1bSJed Brown } 353c4762a1bSJed Brown 354c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 355c4762a1bSJed Brown /* B */ 356c4762a1bSJed Brown PetscErrorCode FormJacobianDesign(Tao tao, Vec X, Mat J, void *ptr) 357c4762a1bSJed Brown { 358c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 359c4762a1bSJed Brown 360c4762a1bSJed Brown PetscFunctionBegin; 3619566063dSJacob Faibussowitsch PetscCall(Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter)); 362c4762a1bSJed Brown PetscFunctionReturn(0); 363c4762a1bSJed Brown } 364c4762a1bSJed Brown 365c4762a1bSJed Brown PetscErrorCode StateMatMult(Mat J_shell, Vec X, Vec Y) 366c4762a1bSJed Brown { 367c4762a1bSJed Brown PetscInt i; 368c4762a1bSJed Brown AppCtx *user; 369c4762a1bSJed Brown 370c4762a1bSJed Brown PetscFunctionBegin; 3719566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell,&user)); 3729566063dSJacob Faibussowitsch PetscCall(Scatter_i(X,user->yi,user->yi_scatter,user->nt)); 3739566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock,user->yi[0],user->yiwork[0])); 374c4762a1bSJed Brown for (i=1; i<user->nt; i++) { 3759566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock,user->yi[i],user->yiwork[i])); 3769566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->yiwork[i],-1.0,user->yi[i-1])); 377c4762a1bSJed Brown } 3789566063dSJacob Faibussowitsch PetscCall(Gather_i(Y,user->yiwork,user->yi_scatter,user->nt)); 379c4762a1bSJed Brown PetscFunctionReturn(0); 380c4762a1bSJed Brown } 381c4762a1bSJed Brown 382c4762a1bSJed Brown PetscErrorCode StateMatMultTranspose(Mat J_shell, Vec X, Vec Y) 383c4762a1bSJed Brown { 384c4762a1bSJed Brown PetscInt i; 385c4762a1bSJed Brown AppCtx *user; 386c4762a1bSJed Brown 387c4762a1bSJed Brown PetscFunctionBegin; 3889566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell,&user)); 3899566063dSJacob Faibussowitsch PetscCall(Scatter_i(X,user->yi,user->yi_scatter,user->nt)); 390c4762a1bSJed Brown for (i=0; i<user->nt-1; i++) { 3919566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock,user->yi[i],user->yiwork[i])); 3929566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->yiwork[i],-1.0,user->yi[i+1])); 393c4762a1bSJed Brown } 394c4762a1bSJed Brown i = user->nt-1; 3959566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock,user->yi[i],user->yiwork[i])); 3969566063dSJacob Faibussowitsch PetscCall(Gather_i(Y,user->yiwork,user->yi_scatter,user->nt)); 397c4762a1bSJed Brown PetscFunctionReturn(0); 398c4762a1bSJed Brown } 399c4762a1bSJed Brown 400c4762a1bSJed Brown PetscErrorCode StateMatBlockMult(Mat J_shell, Vec X, Vec Y) 401c4762a1bSJed Brown { 402c4762a1bSJed Brown AppCtx *user; 403c4762a1bSJed Brown 404c4762a1bSJed Brown PetscFunctionBegin; 4059566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell,&user)); 4069566063dSJacob Faibussowitsch PetscCall(MatMult(user->Grad,X,user->Swork)); 4079566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(user->Swork,user->Swork,user->Av_u)); 4089566063dSJacob Faibussowitsch PetscCall(MatMult(user->Div,user->Swork,Y)); 4099566063dSJacob Faibussowitsch PetscCall(VecAYPX(Y,user->ht,X)); 410c4762a1bSJed Brown PetscFunctionReturn(0); 411c4762a1bSJed Brown } 412c4762a1bSJed Brown 413c4762a1bSJed Brown PetscErrorCode DesignMatMult(Mat J_shell, Vec X, Vec Y) 414c4762a1bSJed Brown { 415c4762a1bSJed Brown PetscInt i; 416c4762a1bSJed Brown AppCtx *user; 417c4762a1bSJed Brown 418c4762a1bSJed Brown PetscFunctionBegin; 4199566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell,&user)); 420c4762a1bSJed Brown 421c4762a1bSJed Brown /* sdiag(1./v) */ 4229566063dSJacob Faibussowitsch PetscCall(VecSet(user->uwork,0)); 4239566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->uwork,-1.0,user->u)); 4249566063dSJacob Faibussowitsch PetscCall(VecExp(user->uwork)); 425c4762a1bSJed Brown 426c4762a1bSJed Brown /* sdiag(1./((Av*(1./v)).^2)) */ 4279566063dSJacob Faibussowitsch PetscCall(MatMult(user->Av,user->uwork,user->Swork)); 4289566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Swork,user->Swork,user->Swork)); 4299566063dSJacob Faibussowitsch PetscCall(VecReciprocal(user->Swork)); 430c4762a1bSJed Brown 431c4762a1bSJed Brown /* (Av * (sdiag(1./v) * b)) */ 4329566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uwork,user->uwork,X)); 4339566063dSJacob Faibussowitsch PetscCall(MatMult(user->Av,user->uwork,user->Twork)); 434c4762a1bSJed Brown 435c4762a1bSJed Brown /* (sdiag(1./((Av*(1./v)).^2)) * (Av * (sdiag(1./v) * b))) */ 4369566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Swork,user->Twork,user->Swork)); 437c4762a1bSJed Brown 4389566063dSJacob Faibussowitsch PetscCall(Scatter_i(user->y,user->yi,user->yi_scatter,user->nt)); 439c4762a1bSJed Brown for (i=0; i<user->nt; i++) { 440c4762a1bSJed Brown /* (sdiag(Grad*y(:,i)) */ 4419566063dSJacob Faibussowitsch PetscCall(MatMult(user->Grad,user->yi[i],user->Twork)); 442c4762a1bSJed Brown 443c4762a1bSJed Brown /* ht * Div * (sdiag(Grad*y(:,i)) * (sdiag(1./((Av*(1./v)).^2)) * (Av * (sdiag(1./v) * b)))) */ 4449566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Twork,user->Twork,user->Swork)); 4459566063dSJacob Faibussowitsch PetscCall(MatMult(user->Div,user->Twork,user->yiwork[i])); 4469566063dSJacob Faibussowitsch PetscCall(VecScale(user->yiwork[i],user->ht)); 447c4762a1bSJed Brown } 4489566063dSJacob Faibussowitsch PetscCall(Gather_i(Y,user->yiwork,user->yi_scatter,user->nt)); 449c4762a1bSJed Brown 450c4762a1bSJed Brown PetscFunctionReturn(0); 451c4762a1bSJed Brown } 452c4762a1bSJed Brown 453c4762a1bSJed Brown PetscErrorCode DesignMatMultTranspose(Mat J_shell, Vec X, Vec Y) 454c4762a1bSJed Brown { 455c4762a1bSJed Brown PetscInt i; 456c4762a1bSJed Brown AppCtx *user; 457c4762a1bSJed Brown 458c4762a1bSJed Brown PetscFunctionBegin; 4599566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell,&user)); 460c4762a1bSJed Brown 461c4762a1bSJed Brown /* sdiag(1./((Av*(1./v)).^2)) */ 4629566063dSJacob Faibussowitsch PetscCall(VecSet(user->uwork,0)); 4639566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->uwork,-1.0,user->u)); 4649566063dSJacob Faibussowitsch PetscCall(VecExp(user->uwork)); 4659566063dSJacob Faibussowitsch PetscCall(MatMult(user->Av,user->uwork,user->Rwork)); 4669566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Rwork,user->Rwork,user->Rwork)); 4679566063dSJacob Faibussowitsch PetscCall(VecReciprocal(user->Rwork)); 468c4762a1bSJed Brown 4699566063dSJacob Faibussowitsch PetscCall(VecSet(Y,0.0)); 4709566063dSJacob Faibussowitsch PetscCall(Scatter_i(user->y,user->yi,user->yi_scatter,user->nt)); 4719566063dSJacob Faibussowitsch PetscCall(Scatter_i(X,user->yiwork,user->yi_scatter,user->nt)); 472c4762a1bSJed Brown for (i=0; i<user->nt; i++) { 473c4762a1bSJed Brown /* (Div' * b(:,i)) */ 4749566063dSJacob Faibussowitsch PetscCall(MatMult(user->Grad,user->yiwork[i],user->Swork)); 475c4762a1bSJed Brown 476c4762a1bSJed Brown /* sdiag(Grad*y(:,i)) */ 4779566063dSJacob Faibussowitsch PetscCall(MatMult(user->Grad,user->yi[i],user->Twork)); 478c4762a1bSJed Brown 479c4762a1bSJed Brown /* (sdiag(Grad*y(:,i)) * (Div' * b(:,i))) */ 4809566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Twork,user->Swork,user->Twork)); 481c4762a1bSJed Brown 482c4762a1bSJed Brown /* (sdiag(1./((Av*(1./v)).^2)) * (sdiag(Grad*y(:,i)) * (Div' * b(:,i)))) */ 4839566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Twork,user->Rwork,user->Twork)); 484c4762a1bSJed Brown 485c4762a1bSJed Brown /* (Av' * (sdiag(1./((Av*(1./v)).^2)) * (sdiag(Grad*y(:,i)) * (Div' * b(:,i))))) */ 4869566063dSJacob Faibussowitsch PetscCall(MatMult(user->AvT,user->Twork,user->yiwork[i])); 487c4762a1bSJed Brown 488c4762a1bSJed Brown /* sdiag(1./v) * (Av' * (sdiag(1./((Av*(1./v)).^2)) * (sdiag(Grad*y(:,i)) * (Div' * b(:,i))))) */ 4899566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->yiwork[i],user->uwork,user->yiwork[i])); 4909566063dSJacob Faibussowitsch PetscCall(VecAXPY(Y,user->ht,user->yiwork[i])); 491c4762a1bSJed Brown } 492c4762a1bSJed Brown PetscFunctionReturn(0); 493c4762a1bSJed Brown } 494c4762a1bSJed Brown 495c4762a1bSJed Brown PetscErrorCode StateMatBlockPrecMult(PC PC_shell, Vec X, Vec Y) 496c4762a1bSJed Brown { 497c4762a1bSJed Brown AppCtx *user; 498c4762a1bSJed Brown 499c4762a1bSJed Brown PetscFunctionBegin; 5009566063dSJacob Faibussowitsch PetscCall(PCShellGetContext(PC_shell,&user)); 501c4762a1bSJed Brown 502c4762a1bSJed Brown if (user->dsg_formed) { 5039566063dSJacob Faibussowitsch PetscCall(MatSOR(user->DSG,X,1.0,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_SYMMETRIC_SWEEP),0.0,1,1,Y)); 504c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"DSG not formed"); 505c4762a1bSJed Brown PetscFunctionReturn(0); 506c4762a1bSJed Brown } 507c4762a1bSJed Brown 508c4762a1bSJed Brown PetscErrorCode StateMatInvMult(Mat J_shell, Vec X, Vec Y) 509c4762a1bSJed Brown { 510c4762a1bSJed Brown AppCtx *user; 511c4762a1bSJed Brown PetscInt its,i; 512c4762a1bSJed Brown 513c4762a1bSJed Brown PetscFunctionBegin; 5149566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell,&user)); 515c4762a1bSJed Brown 516c4762a1bSJed Brown if (Y == user->ytrue) { 517c4762a1bSJed Brown /* First solve is done with true solution to set up problem */ 5189566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(user->solver,1e-8,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT)); 519c4762a1bSJed Brown } else { 5209566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(user->solver,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT)); 521c4762a1bSJed Brown } 522c4762a1bSJed Brown 5239566063dSJacob Faibussowitsch PetscCall(Scatter_i(X,user->yi,user->yi_scatter,user->nt)); 5249566063dSJacob Faibussowitsch PetscCall(KSPSolve(user->solver,user->yi[0],user->yiwork[0])); 5259566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(user->solver,&its)); 526c4762a1bSJed Brown user->ksp_its = user->ksp_its + its; 527c4762a1bSJed Brown 528c4762a1bSJed Brown for (i=1; i<user->nt; i++) { 5299566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->yi[i],1.0,user->yiwork[i-1])); 5309566063dSJacob Faibussowitsch PetscCall(KSPSolve(user->solver,user->yi[i],user->yiwork[i])); 5319566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(user->solver,&its)); 532c4762a1bSJed Brown user->ksp_its = user->ksp_its + its; 533c4762a1bSJed Brown } 5349566063dSJacob Faibussowitsch PetscCall(Gather_i(Y,user->yiwork,user->yi_scatter,user->nt)); 535c4762a1bSJed Brown PetscFunctionReturn(0); 536c4762a1bSJed Brown } 537c4762a1bSJed Brown 538c4762a1bSJed Brown PetscErrorCode StateMatInvTransposeMult(Mat J_shell, Vec X, Vec Y) 539c4762a1bSJed Brown { 540c4762a1bSJed Brown AppCtx *user; 541c4762a1bSJed Brown PetscInt its,i; 542c4762a1bSJed Brown 543c4762a1bSJed Brown PetscFunctionBegin; 5449566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell,&user)); 545c4762a1bSJed Brown 5469566063dSJacob Faibussowitsch PetscCall(Scatter_i(X,user->yi,user->yi_scatter,user->nt)); 547c4762a1bSJed Brown 548c4762a1bSJed Brown i = user->nt - 1; 5499566063dSJacob Faibussowitsch PetscCall(KSPSolve(user->solver,user->yi[i],user->yiwork[i])); 550c4762a1bSJed Brown 5519566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(user->solver,&its)); 552c4762a1bSJed Brown user->ksp_its = user->ksp_its + its; 553c4762a1bSJed Brown 554c4762a1bSJed Brown for (i=user->nt-2; i>=0; i--) { 5559566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->yi[i],1.0,user->yiwork[i+1])); 5569566063dSJacob Faibussowitsch PetscCall(KSPSolve(user->solver,user->yi[i],user->yiwork[i])); 557c4762a1bSJed Brown 5589566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(user->solver,&its)); 559c4762a1bSJed Brown user->ksp_its = user->ksp_its + its; 560c4762a1bSJed Brown 561c4762a1bSJed Brown } 562c4762a1bSJed Brown 5639566063dSJacob Faibussowitsch PetscCall(Gather_i(Y,user->yiwork,user->yi_scatter,user->nt)); 564c4762a1bSJed Brown PetscFunctionReturn(0); 565c4762a1bSJed Brown } 566c4762a1bSJed Brown 567c4762a1bSJed Brown PetscErrorCode StateMatDuplicate(Mat J_shell, MatDuplicateOption opt, Mat *new_shell) 568c4762a1bSJed Brown { 569c4762a1bSJed Brown AppCtx *user; 570c4762a1bSJed Brown 571c4762a1bSJed Brown PetscFunctionBegin; 5729566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell,&user)); 573c4762a1bSJed Brown 5749566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,new_shell)); 5759566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*new_shell,MATOP_MULT,(void(*)(void))StateMatMult)); 5769566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*new_shell,MATOP_DUPLICATE,(void(*)(void))StateMatDuplicate)); 5779566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*new_shell,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatMultTranspose)); 5789566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*new_shell,MATOP_GET_DIAGONAL,(void(*)(void))StateMatGetDiagonal)); 579c4762a1bSJed Brown PetscFunctionReturn(0); 580c4762a1bSJed Brown } 581c4762a1bSJed Brown 582c4762a1bSJed Brown PetscErrorCode StateMatGetDiagonal(Mat J_shell, Vec X) 583c4762a1bSJed Brown { 584c4762a1bSJed Brown AppCtx *user; 585c4762a1bSJed Brown 586c4762a1bSJed Brown PetscFunctionBegin; 5879566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell,&user)); 5889566063dSJacob Faibussowitsch PetscCall(VecCopy(user->js_diag,X)); 589c4762a1bSJed Brown PetscFunctionReturn(0); 590c4762a1bSJed Brown 591c4762a1bSJed Brown } 592c4762a1bSJed Brown 593c4762a1bSJed Brown PetscErrorCode FormConstraints(Tao tao, Vec X, Vec C, void *ptr) 594c4762a1bSJed Brown { 595c4762a1bSJed Brown /* con = Ay - q, A = [B 0 0 ... 0; 596c4762a1bSJed Brown -I B 0 ... 0; 597c4762a1bSJed Brown 0 -I B ... 0; 598c4762a1bSJed Brown ... ; 599c4762a1bSJed Brown 0 ... -I B] 600c4762a1bSJed Brown B = ht * Div * Sigma * Grad + eye */ 601c4762a1bSJed Brown PetscInt i; 602c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 603c4762a1bSJed Brown 604c4762a1bSJed Brown PetscFunctionBegin; 6059566063dSJacob Faibussowitsch PetscCall(Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter)); 6069566063dSJacob Faibussowitsch PetscCall(Scatter_i(user->y,user->yi,user->yi_scatter,user->nt)); 6079566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock,user->yi[0],user->yiwork[0])); 608c4762a1bSJed Brown for (i=1; i<user->nt; i++) { 6099566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock,user->yi[i],user->yiwork[i])); 6109566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->yiwork[i],-1.0,user->yi[i-1])); 611c4762a1bSJed Brown } 6129566063dSJacob Faibussowitsch PetscCall(Gather_i(C,user->yiwork,user->yi_scatter,user->nt)); 6139566063dSJacob Faibussowitsch PetscCall(VecAXPY(C,-1.0,user->q)); 614c4762a1bSJed Brown PetscFunctionReturn(0); 615c4762a1bSJed Brown } 616c4762a1bSJed Brown 617c4762a1bSJed Brown PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat) 618c4762a1bSJed Brown { 619c4762a1bSJed Brown PetscFunctionBegin; 6209566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(s_scat,x,state,INSERT_VALUES,SCATTER_FORWARD)); 6219566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(s_scat,x,state,INSERT_VALUES,SCATTER_FORWARD)); 6229566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(d_scat,x,design,INSERT_VALUES,SCATTER_FORWARD)); 6239566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(d_scat,x,design,INSERT_VALUES,SCATTER_FORWARD)); 624c4762a1bSJed Brown PetscFunctionReturn(0); 625c4762a1bSJed Brown } 626c4762a1bSJed Brown 627c4762a1bSJed Brown PetscErrorCode Scatter_i(Vec y, Vec *yi, VecScatter *scat, PetscInt nt) 628c4762a1bSJed Brown { 629c4762a1bSJed Brown PetscInt i; 630c4762a1bSJed Brown 631c4762a1bSJed Brown PetscFunctionBegin; 632c4762a1bSJed Brown for (i=0; i<nt; i++) { 6339566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scat[i],y,yi[i],INSERT_VALUES,SCATTER_FORWARD)); 6349566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scat[i],y,yi[i],INSERT_VALUES,SCATTER_FORWARD)); 635c4762a1bSJed Brown } 636c4762a1bSJed Brown PetscFunctionReturn(0); 637c4762a1bSJed Brown } 638c4762a1bSJed Brown 639c4762a1bSJed Brown PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat) 640c4762a1bSJed Brown { 641c4762a1bSJed Brown PetscFunctionBegin; 6429566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(s_scat,state,x,INSERT_VALUES,SCATTER_REVERSE)); 6439566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(s_scat,state,x,INSERT_VALUES,SCATTER_REVERSE)); 6449566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(d_scat,design,x,INSERT_VALUES,SCATTER_REVERSE)); 6459566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(d_scat,design,x,INSERT_VALUES,SCATTER_REVERSE)); 646c4762a1bSJed Brown PetscFunctionReturn(0); 647c4762a1bSJed Brown } 648c4762a1bSJed Brown 649c4762a1bSJed Brown PetscErrorCode Gather_i(Vec y, Vec *yi, VecScatter *scat, PetscInt nt) 650c4762a1bSJed Brown { 651c4762a1bSJed Brown PetscInt i; 652c4762a1bSJed Brown 653c4762a1bSJed Brown PetscFunctionBegin; 654c4762a1bSJed Brown for (i=0; i<nt; i++) { 6559566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scat[i],yi[i],y,INSERT_VALUES,SCATTER_REVERSE)); 6569566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scat[i],yi[i],y,INSERT_VALUES,SCATTER_REVERSE)); 657c4762a1bSJed Brown } 658c4762a1bSJed Brown PetscFunctionReturn(0); 659c4762a1bSJed Brown } 660c4762a1bSJed Brown 661c4762a1bSJed Brown PetscErrorCode ParabolicInitialize(AppCtx *user) 662c4762a1bSJed Brown { 663c4762a1bSJed Brown PetscInt m,n,i,j,k,linear_index,istart,iend,iblock,lo,hi,lo2,hi2; 664c4762a1bSJed Brown Vec XX,YY,ZZ,XXwork,YYwork,ZZwork,UTwork,yi,di,bc; 665c4762a1bSJed Brown PetscReal *x, *y, *z; 666c4762a1bSJed Brown PetscReal h,stime; 667c4762a1bSJed Brown PetscScalar hinv,neg_hinv,half = 0.5,sqrt_beta; 668c4762a1bSJed Brown PetscInt im,indx1,indx2,indy1,indy2,indz1,indz2,nx,ny,nz; 669c4762a1bSJed Brown PetscReal xri,yri,zri,xim,yim,zim,dx1,dx2,dy1,dy2,dz1,dz2,Dx,Dy,Dz; 670c4762a1bSJed Brown PetscScalar v,vx,vy,vz; 671c4762a1bSJed Brown IS is_from_y,is_to_yi,is_from_d,is_to_di; 672c4762a1bSJed Brown /* Data locations */ 673c4762a1bSJed Brown PetscScalar xr[64] = {0.4970, 0.8498, 0.7814, 0.6268, 0.7782, 0.6402, 0.3617, 0.3160, 674c4762a1bSJed Brown 0.3610, 0.5298, 0.6987, 0.3331, 0.7962, 0.5596, 0.3866, 0.6774, 675c4762a1bSJed Brown 0.5407, 0.4518, 0.6702, 0.6061, 0.7580, 0.8997, 0.5198, 0.8326, 676c4762a1bSJed Brown 0.2138, 0.9198, 0.3000, 0.2833, 0.8288, 0.7076, 0.1820, 0.0728, 677c4762a1bSJed Brown 0.8447, 0.2367, 0.3239, 0.6413, 0.3114, 0.4731, 0.1192, 0.9273, 678c4762a1bSJed Brown 0.5724, 0.4331, 0.5136, 0.3547, 0.4413, 0.2602, 0.5698, 0.7278, 679c4762a1bSJed Brown 0.5261, 0.6230, 0.2454, 0.3948, 0.7479, 0.6582, 0.4660, 0.5594, 680c4762a1bSJed Brown 0.7574, 0.1143, 0.5900, 0.1065, 0.4260, 0.3294, 0.8276, 0.0756}; 681c4762a1bSJed Brown 682c4762a1bSJed Brown PetscScalar yr[64] = {0.7345, 0.9120, 0.9288, 0.7528, 0.4463, 0.4985, 0.2497, 0.6256, 683c4762a1bSJed Brown 0.3425, 0.9026, 0.6983, 0.4230, 0.7140, 0.2970, 0.4474, 0.8792, 684c4762a1bSJed Brown 0.6604, 0.2485, 0.7968, 0.6127, 0.1796, 0.2437, 0.5938, 0.6137, 685c4762a1bSJed Brown 0.3867, 0.5658, 0.4575, 0.1009, 0.0863, 0.3361, 0.0738, 0.3985, 686c4762a1bSJed Brown 0.6602, 0.1437, 0.0934, 0.5983, 0.5950, 0.0763, 0.0768, 0.2288, 687c4762a1bSJed Brown 0.5761, 0.1129, 0.3841, 0.6150, 0.6904, 0.6686, 0.1361, 0.4601, 688c4762a1bSJed Brown 0.4491, 0.3716, 0.1969, 0.6537, 0.6743, 0.6991, 0.4811, 0.5480, 689c4762a1bSJed Brown 0.1684, 0.4569, 0.6889, 0.8437, 0.3015, 0.2854, 0.8199, 0.2658}; 690c4762a1bSJed Brown 691c4762a1bSJed Brown PetscScalar zr[64] = {0.7668, 0.8573, 0.2654, 0.2719, 0.1060, 0.1311, 0.6232, 0.2295, 692c4762a1bSJed Brown 0.8009, 0.2147, 0.2119, 0.9325, 0.4473, 0.3600, 0.3374, 0.3819, 693c4762a1bSJed Brown 0.4066, 0.5801, 0.1673, 0.0959, 0.4638, 0.8236, 0.8800, 0.2939, 694c4762a1bSJed Brown 0.2028, 0.8262, 0.2706, 0.6276, 0.9085, 0.6443, 0.8241, 0.0712, 695c4762a1bSJed Brown 0.1824, 0.7789, 0.4389, 0.8415, 0.7055, 0.6639, 0.3653, 0.2078, 696c4762a1bSJed Brown 0.1987, 0.2297, 0.4321, 0.8115, 0.4915, 0.7764, 0.4657, 0.4627, 697c4762a1bSJed Brown 0.4569, 0.4232, 0.8514, 0.0674, 0.3227, 0.1055, 0.6690, 0.6313, 698c4762a1bSJed Brown 0.9226, 0.5461, 0.4126, 0.2364, 0.6096, 0.7042, 0.3914, 0.0711}; 699c4762a1bSJed Brown 700c4762a1bSJed Brown PetscFunctionBegin; 7019566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->mx,&x)); 7029566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->mx,&y)); 7039566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->mx,&z)); 704c4762a1bSJed Brown user->jformed = PETSC_FALSE; 705c4762a1bSJed Brown user->dsg_formed = PETSC_FALSE; 706c4762a1bSJed Brown 707c4762a1bSJed Brown n = user->mx * user->mx * user->mx; 708c4762a1bSJed Brown m = 3 * user->mx * user->mx * (user->mx-1); 709c4762a1bSJed Brown sqrt_beta = PetscSqrtScalar(user->beta); 710c4762a1bSJed Brown 711c4762a1bSJed Brown user->ksp_its = 0; 712c4762a1bSJed Brown user->ksp_its_initial = 0; 713c4762a1bSJed Brown 714c4762a1bSJed Brown stime = (PetscReal)user->nt/user->ns; 7159566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->ns,&user->sample_times)); 716c4762a1bSJed Brown for (i=0; i<user->ns; i++) { 717c4762a1bSJed Brown user->sample_times[i] = (PetscInt)(stime*((PetscReal)i+1.0)-0.5); 718c4762a1bSJed Brown } 719c4762a1bSJed Brown 7209566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&XX)); 7219566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&user->q)); 7229566063dSJacob Faibussowitsch PetscCall(VecSetSizes(XX,PETSC_DECIDE,n)); 7239566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->q,PETSC_DECIDE,n*user->nt)); 7249566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(XX)); 7259566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->q)); 726c4762a1bSJed Brown 7279566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX,&YY)); 7289566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX,&ZZ)); 7299566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX,&XXwork)); 7309566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX,&YYwork)); 7319566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX,&ZZwork)); 7329566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX,&UTwork)); 7339566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX,&user->utrue)); 7349566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX,&bc)); 735c4762a1bSJed Brown 736c4762a1bSJed Brown /* Generate 3D grid, and collect ns (1<=ns<=8) right-hand-side vectors into user->q */ 737c4762a1bSJed Brown h = 1.0/user->mx; 738c4762a1bSJed Brown hinv = user->mx; 739c4762a1bSJed Brown neg_hinv = -hinv; 740c4762a1bSJed Brown 7419566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(XX,&istart,&iend)); 742c4762a1bSJed Brown for (linear_index=istart; linear_index<iend; linear_index++) { 743c4762a1bSJed Brown i = linear_index % user->mx; 744c4762a1bSJed Brown j = ((linear_index-i)/user->mx) % user->mx; 745c4762a1bSJed Brown k = ((linear_index-i)/user->mx-j) / user->mx; 746c4762a1bSJed Brown vx = h*(i+0.5); 747c4762a1bSJed Brown vy = h*(j+0.5); 748c4762a1bSJed Brown vz = h*(k+0.5); 7499566063dSJacob Faibussowitsch PetscCall(VecSetValues(XX,1,&linear_index,&vx,INSERT_VALUES)); 7509566063dSJacob Faibussowitsch PetscCall(VecSetValues(YY,1,&linear_index,&vy,INSERT_VALUES)); 7519566063dSJacob Faibussowitsch PetscCall(VecSetValues(ZZ,1,&linear_index,&vz,INSERT_VALUES)); 752c4762a1bSJed Brown if ((vx<0.6) && (vx>0.4) && (vy<0.6) && (vy>0.4) && (vy<0.6) && (vz<0.6) && (vz>0.4)) { 753c4762a1bSJed Brown v = 1000.0; 7549566063dSJacob Faibussowitsch PetscCall(VecSetValues(bc,1,&linear_index,&v,INSERT_VALUES)); 755c4762a1bSJed Brown } 756c4762a1bSJed Brown } 757c4762a1bSJed Brown 7589566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(XX)); 7599566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(XX)); 7609566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(YY)); 7619566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(YY)); 7629566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(ZZ)); 7639566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(ZZ)); 7649566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(bc)); 7659566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(bc)); 766c4762a1bSJed Brown 767c4762a1bSJed Brown /* Compute true parameter function 768c4762a1bSJed Brown ut = 0.5 + exp(-10*((x-0.5)^2+(y-0.5)^2+(z-0.5)^2)) */ 7699566063dSJacob Faibussowitsch PetscCall(VecCopy(XX,XXwork)); 7709566063dSJacob Faibussowitsch PetscCall(VecCopy(YY,YYwork)); 7719566063dSJacob Faibussowitsch PetscCall(VecCopy(ZZ,ZZwork)); 772c4762a1bSJed Brown 7739566063dSJacob Faibussowitsch PetscCall(VecShift(XXwork,-0.5)); 7749566063dSJacob Faibussowitsch PetscCall(VecShift(YYwork,-0.5)); 7759566063dSJacob Faibussowitsch PetscCall(VecShift(ZZwork,-0.5)); 776c4762a1bSJed Brown 7779566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(XXwork,XXwork,XXwork)); 7789566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(YYwork,YYwork,YYwork)); 7799566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(ZZwork,ZZwork,ZZwork)); 780c4762a1bSJed Brown 7819566063dSJacob Faibussowitsch PetscCall(VecCopy(XXwork,user->utrue)); 7829566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->utrue,1.0,YYwork)); 7839566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->utrue,1.0,ZZwork)); 7849566063dSJacob Faibussowitsch PetscCall(VecScale(user->utrue,-10.0)); 7859566063dSJacob Faibussowitsch PetscCall(VecExp(user->utrue)); 7869566063dSJacob Faibussowitsch PetscCall(VecShift(user->utrue,0.5)); 787c4762a1bSJed Brown 7889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&XX)); 7899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&YY)); 7909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ZZ)); 7919566063dSJacob Faibussowitsch PetscCall(VecDestroy(&XXwork)); 7929566063dSJacob Faibussowitsch PetscCall(VecDestroy(&YYwork)); 7939566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ZZwork)); 7949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&UTwork)); 795c4762a1bSJed Brown 796c4762a1bSJed Brown /* Initial guess and reference model */ 7979566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->utrue,&user->ur)); 7989566063dSJacob Faibussowitsch PetscCall(VecSet(user->ur,0.0)); 799c4762a1bSJed Brown 800c4762a1bSJed Brown /* Generate Grad matrix */ 8019566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&user->Grad)); 8029566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Grad,PETSC_DECIDE,PETSC_DECIDE,m,n)); 8039566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Grad)); 8049566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->Grad,2,NULL,2,NULL)); 8059566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->Grad,2,NULL)); 8069566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->Grad,&istart,&iend)); 807c4762a1bSJed Brown 808c4762a1bSJed Brown for (i=istart; i<iend; i++) { 809c4762a1bSJed Brown if (i<m/3) { 810c4762a1bSJed Brown iblock = i / (user->mx-1); 811c4762a1bSJed Brown j = iblock*user->mx + (i % (user->mx-1)); 8129566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad,1,&i,1,&j,&neg_hinv,INSERT_VALUES)); 813c4762a1bSJed Brown j = j+1; 8149566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad,1,&i,1,&j,&hinv,INSERT_VALUES)); 815c4762a1bSJed Brown } 816c4762a1bSJed Brown if (i>=m/3 && i<2*m/3) { 817c4762a1bSJed Brown iblock = (i-m/3) / (user->mx*(user->mx-1)); 818c4762a1bSJed Brown j = iblock*user->mx*user->mx + ((i-m/3) % (user->mx*(user->mx-1))); 8199566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad,1,&i,1,&j,&neg_hinv,INSERT_VALUES)); 820c4762a1bSJed Brown j = j + user->mx; 8219566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad,1,&i,1,&j,&hinv,INSERT_VALUES)); 822c4762a1bSJed Brown } 823c4762a1bSJed Brown if (i>=2*m/3) { 824c4762a1bSJed Brown j = i-2*m/3; 8259566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad,1,&i,1,&j,&neg_hinv,INSERT_VALUES)); 826c4762a1bSJed Brown j = j + user->mx*user->mx; 8279566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad,1,&i,1,&j,&hinv,INSERT_VALUES)); 828c4762a1bSJed Brown } 829c4762a1bSJed Brown } 830c4762a1bSJed Brown 8319566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Grad,MAT_FINAL_ASSEMBLY)); 8329566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Grad,MAT_FINAL_ASSEMBLY)); 833c4762a1bSJed Brown 834c4762a1bSJed Brown /* Generate arithmetic averaging matrix Av */ 8359566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&user->Av)); 8369566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Av,PETSC_DECIDE,PETSC_DECIDE,m,n)); 8379566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Av)); 8389566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->Av,2,NULL,2,NULL)); 8399566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->Av,2,NULL)); 8409566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->Av,&istart,&iend)); 841c4762a1bSJed Brown 842c4762a1bSJed Brown for (i=istart; i<iend; i++) { 843c4762a1bSJed Brown if (i<m/3) { 844c4762a1bSJed Brown iblock = i / (user->mx-1); 845c4762a1bSJed Brown j = iblock*user->mx + (i % (user->mx-1)); 8469566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES)); 847c4762a1bSJed Brown j = j+1; 8489566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES)); 849c4762a1bSJed Brown } 850c4762a1bSJed Brown if (i>=m/3 && i<2*m/3) { 851c4762a1bSJed Brown iblock = (i-m/3) / (user->mx*(user->mx-1)); 852c4762a1bSJed Brown j = iblock*user->mx*user->mx + ((i-m/3) % (user->mx*(user->mx-1))); 8539566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES)); 854c4762a1bSJed Brown j = j + user->mx; 8559566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES)); 856c4762a1bSJed Brown } 857c4762a1bSJed Brown if (i>=2*m/3) { 858c4762a1bSJed Brown j = i-2*m/3; 8599566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES)); 860c4762a1bSJed Brown j = j + user->mx*user->mx; 8619566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES)); 862c4762a1bSJed Brown } 863c4762a1bSJed Brown } 864c4762a1bSJed Brown 8659566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Av,MAT_FINAL_ASSEMBLY)); 8669566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Av,MAT_FINAL_ASSEMBLY)); 867c4762a1bSJed Brown 868c4762a1bSJed Brown /* Generate transpose of averaging matrix Av */ 8699566063dSJacob Faibussowitsch PetscCall(MatTranspose(user->Av,MAT_INITIAL_MATRIX,&user->AvT)); 870c4762a1bSJed Brown 8719566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&user->L)); 8729566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->L,PETSC_DECIDE,PETSC_DECIDE,m+n,n)); 8739566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->L)); 8749566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->L,2,NULL,2,NULL)); 8759566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->L,2,NULL)); 8769566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->L,&istart,&iend)); 877c4762a1bSJed Brown 878c4762a1bSJed Brown for (i=istart; i<iend; i++) { 879c4762a1bSJed Brown if (i<m/3) { 880c4762a1bSJed Brown iblock = i / (user->mx-1); 881c4762a1bSJed Brown j = iblock*user->mx + (i % (user->mx-1)); 8829566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L,1,&i,1,&j,&neg_hinv,INSERT_VALUES)); 883c4762a1bSJed Brown j = j+1; 8849566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L,1,&i,1,&j,&hinv,INSERT_VALUES)); 885c4762a1bSJed Brown } 886c4762a1bSJed Brown if (i>=m/3 && i<2*m/3) { 887c4762a1bSJed Brown iblock = (i-m/3) / (user->mx*(user->mx-1)); 888c4762a1bSJed Brown j = iblock*user->mx*user->mx + ((i-m/3) % (user->mx*(user->mx-1))); 8899566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L,1,&i,1,&j,&neg_hinv,INSERT_VALUES)); 890c4762a1bSJed Brown j = j + user->mx; 8919566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L,1,&i,1,&j,&hinv,INSERT_VALUES)); 892c4762a1bSJed Brown } 893c4762a1bSJed Brown if (i>=2*m/3 && i<m) { 894c4762a1bSJed Brown j = i-2*m/3; 8959566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L,1,&i,1,&j,&neg_hinv,INSERT_VALUES)); 896c4762a1bSJed Brown j = j + user->mx*user->mx; 8979566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L,1,&i,1,&j,&hinv,INSERT_VALUES)); 898c4762a1bSJed Brown } 899c4762a1bSJed Brown if (i>=m) { 900c4762a1bSJed Brown j = i - m; 9019566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L,1,&i,1,&j,&sqrt_beta,INSERT_VALUES)); 902c4762a1bSJed Brown } 903c4762a1bSJed Brown } 904c4762a1bSJed Brown 9059566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->L,MAT_FINAL_ASSEMBLY)); 9069566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->L,MAT_FINAL_ASSEMBLY)); 907c4762a1bSJed Brown 9089566063dSJacob Faibussowitsch PetscCall(MatScale(user->L,PetscPowScalar(h,1.5))); 909c4762a1bSJed Brown 910c4762a1bSJed Brown /* Generate Div matrix */ 9119566063dSJacob Faibussowitsch PetscCall(MatTranspose(user->Grad,MAT_INITIAL_MATRIX,&user->Div)); 912c4762a1bSJed Brown 913c4762a1bSJed Brown /* Build work vectors and matrices */ 9149566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&user->S)); 9159566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->S, PETSC_DECIDE, user->mx*user->mx*(user->mx-1)*3)); 9169566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->S)); 917c4762a1bSJed Brown 9189566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&user->lwork)); 9199566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->lwork,PETSC_DECIDE,m+user->mx*user->mx*user->mx)); 9209566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->lwork)); 921c4762a1bSJed Brown 9229566063dSJacob Faibussowitsch PetscCall(MatDuplicate(user->Div,MAT_SHARE_NONZERO_PATTERN,&user->Divwork)); 9239566063dSJacob Faibussowitsch PetscCall(MatDuplicate(user->Av,MAT_SHARE_NONZERO_PATTERN,&user->Avwork)); 924c4762a1bSJed Brown 9259566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&user->d)); 9269566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->d,PETSC_DECIDE,user->ndata*user->nt)); 9279566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->d)); 928c4762a1bSJed Brown 9299566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->S,&user->Swork)); 9309566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->S,&user->Av_u)); 9319566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->S,&user->Twork)); 9329566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->S,&user->Rwork)); 9339566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->y,&user->ywork)); 9349566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->u,&user->uwork)); 9359566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->u,&user->js_diag)); 9369566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->c,&user->cwork)); 9379566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->d,&user->dwork)); 938c4762a1bSJed Brown 939c4762a1bSJed Brown /* Create matrix-free shell user->Js for computing A*x */ 9409566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m*user->nt,user->m*user->nt,user,&user->Js)); 9419566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->Js,MATOP_MULT,(void(*)(void))StateMatMult)); 9429566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->Js,MATOP_DUPLICATE,(void(*)(void))StateMatDuplicate)); 9439566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->Js,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatMultTranspose)); 9449566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->Js,MATOP_GET_DIAGONAL,(void(*)(void))StateMatGetDiagonal)); 945c4762a1bSJed Brown 946c4762a1bSJed Brown /* Diagonal blocks of user->Js */ 9479566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,&user->JsBlock)); 9489566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->JsBlock,MATOP_MULT,(void(*)(void))StateMatBlockMult)); 949c4762a1bSJed Brown /* Blocks are symmetric */ 9509566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->JsBlock,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatBlockMult)); 951c4762a1bSJed Brown 952c4762a1bSJed Brown /* Create a matrix-free shell user->JsBlockPrec for computing (U+D)\D*(L+D)\x, where JsBlock = L+D+U, 953c4762a1bSJed Brown D is diagonal, L is strictly lower triangular, and U is strictly upper triangular. 954c4762a1bSJed Brown This is an SSOR preconditioner for user->JsBlock. */ 9559566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,&user->JsBlockPrec)); 9569566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->JsBlockPrec,MATOP_MULT,(void(*)(void))StateMatBlockPrecMult)); 957c4762a1bSJed Brown /* JsBlockPrec is symmetric */ 9589566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->JsBlockPrec,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatBlockPrecMult)); 9599566063dSJacob Faibussowitsch PetscCall(MatSetOption(user->JsBlockPrec,MAT_SYMMETRY_ETERNAL,PETSC_TRUE)); 960c4762a1bSJed Brown 961c4762a1bSJed Brown /* Create a matrix-free shell user->Jd for computing B*x */ 9629566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m*user->nt,user->m,user,&user->Jd)); 9639566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->Jd,MATOP_MULT,(void(*)(void))DesignMatMult)); 9649566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->Jd,MATOP_MULT_TRANSPOSE,(void(*)(void))DesignMatMultTranspose)); 965c4762a1bSJed Brown 966c4762a1bSJed Brown /* User-defined routines for computing user->Js\x and user->Js^T\x*/ 9679566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m*user->nt,user->m*user->nt,user,&user->JsInv)); 9689566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->JsInv,MATOP_MULT,(void(*)(void))StateMatInvMult)); 9699566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->JsInv,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatInvTransposeMult)); 970c4762a1bSJed Brown 971c4762a1bSJed Brown /* Solver options and tolerances */ 9729566063dSJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_WORLD,&user->solver)); 9739566063dSJacob Faibussowitsch PetscCall(KSPSetType(user->solver,KSPCG)); 9749566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(user->solver,user->JsBlock,user->JsBlockPrec)); 9759566063dSJacob Faibussowitsch PetscCall(KSPSetInitialGuessNonzero(user->solver,PETSC_FALSE)); 9769566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(user->solver,1e-4,1e-20,1e3,500)); 9779566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(user->solver)); 9789566063dSJacob Faibussowitsch PetscCall(KSPGetPC(user->solver,&user->prec)); 9799566063dSJacob Faibussowitsch PetscCall(PCSetType(user->prec,PCSHELL)); 980c4762a1bSJed Brown 9819566063dSJacob Faibussowitsch PetscCall(PCShellSetApply(user->prec,StateMatBlockPrecMult)); 9829566063dSJacob Faibussowitsch PetscCall(PCShellSetApplyTranspose(user->prec,StateMatBlockPrecMult)); 9839566063dSJacob Faibussowitsch PetscCall(PCShellSetContext(user->prec,user)); 984c4762a1bSJed Brown 985c4762a1bSJed Brown /* Create scatter from y to y_1,y_2,...,y_nt */ 9869566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->nt*user->m,&user->yi_scatter)); 9879566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&yi)); 9889566063dSJacob Faibussowitsch PetscCall(VecSetSizes(yi,PETSC_DECIDE,user->mx*user->mx*user->mx)); 9899566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(yi)); 9909566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(yi,user->nt,&user->yi)); 9919566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(yi,user->nt,&user->yiwork)); 992c4762a1bSJed Brown 9939566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->y,&lo2,&hi2)); 994c4762a1bSJed Brown istart = 0; 995c4762a1bSJed Brown for (i=0; i<user->nt; i++) { 9969566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->yi[i],&lo,&hi)); 9979566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_yi)); 9989566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo2+istart,1,&is_from_y)); 9999566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(user->y,is_from_y,user->yi[i],is_to_yi,&user->yi_scatter[i])); 1000c4762a1bSJed Brown istart = istart + hi-lo; 10019566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_to_yi)); 10029566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_from_y)); 1003c4762a1bSJed Brown } 10049566063dSJacob Faibussowitsch PetscCall(VecDestroy(&yi)); 1005c4762a1bSJed Brown 1006c4762a1bSJed Brown /* Create scatter from d to d_1,d_2,...,d_ns */ 10079566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->ns*user->ndata,&user->di_scatter)); 10089566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&di)); 10099566063dSJacob Faibussowitsch PetscCall(VecSetSizes(di,PETSC_DECIDE,user->ndata)); 10109566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(di)); 10119566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(di,user->ns,&user->di)); 10129566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->d,&lo2,&hi2)); 1013c4762a1bSJed Brown istart = 0; 1014c4762a1bSJed Brown for (i=0; i<user->ns; i++) { 10159566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->di[i],&lo,&hi)); 10169566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_di)); 10179566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo2+istart,1,&is_from_d)); 10189566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(user->d,is_from_d,user->di[i],is_to_di,&user->di_scatter[i])); 1019c4762a1bSJed Brown istart = istart + hi-lo; 10209566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_to_di)); 10219566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_from_d)); 1022c4762a1bSJed Brown } 10239566063dSJacob Faibussowitsch PetscCall(VecDestroy(&di)); 1024c4762a1bSJed Brown 1025c4762a1bSJed Brown /* Assemble RHS of forward problem */ 10269566063dSJacob Faibussowitsch PetscCall(VecCopy(bc,user->yiwork[0])); 1027c4762a1bSJed Brown for (i=1; i<user->nt; i++) { 10289566063dSJacob Faibussowitsch PetscCall(VecSet(user->yiwork[i],0.0)); 1029c4762a1bSJed Brown } 10309566063dSJacob Faibussowitsch PetscCall(Gather_i(user->q,user->yiwork,user->yi_scatter,user->nt)); 10319566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bc)); 1032c4762a1bSJed Brown 1033c4762a1bSJed Brown /* Compute true state function yt given ut */ 10349566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&user->ytrue)); 10359566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->ytrue,PETSC_DECIDE,n*user->nt)); 10369566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->ytrue)); 1037c4762a1bSJed Brown 1038c4762a1bSJed Brown /* First compute Av_u = Av*exp(-u) */ 10399566063dSJacob Faibussowitsch PetscCall(VecSet(user->uwork,0)); 10409566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->uwork,-1.0,user->utrue)); /* Note: user->utrue */ 10419566063dSJacob Faibussowitsch PetscCall(VecExp(user->uwork)); 10429566063dSJacob Faibussowitsch PetscCall(MatMult(user->Av,user->uwork,user->Av_u)); 1043c4762a1bSJed Brown 1044c20d7725SJed Brown /* Symbolic DSG = Div * Grad */ 10459566063dSJacob Faibussowitsch PetscCall(MatProductCreate(user->Div,user->Grad,NULL,&user->DSG)); 10469566063dSJacob Faibussowitsch PetscCall(MatProductSetType(user->DSG,MATPRODUCT_AB)); 10479566063dSJacob Faibussowitsch PetscCall(MatProductSetAlgorithm(user->DSG,"default")); 10489566063dSJacob Faibussowitsch PetscCall(MatProductSetFill(user->DSG,PETSC_DEFAULT)); 10499566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(user->DSG)); 10509566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(user->DSG)); 1051c20d7725SJed Brown 1052c4762a1bSJed Brown user->dsg_formed = PETSC_TRUE; 1053c4762a1bSJed Brown 1054c20d7725SJed Brown /* Next form DSG = Div*Grad */ 10559566063dSJacob Faibussowitsch PetscCall(MatCopy(user->Div,user->Divwork,SAME_NONZERO_PATTERN)); 10569566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(user->Divwork,NULL,user->Av_u)); 1057c4762a1bSJed Brown if (user->dsg_formed) { 10589566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(user->DSG)); 1059c4762a1bSJed Brown } else { 10609566063dSJacob Faibussowitsch PetscCall(MatMatMult(user->Div,user->Grad,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&user->DSG)); 1061c4762a1bSJed Brown user->dsg_formed = PETSC_TRUE; 1062c4762a1bSJed Brown } 1063c4762a1bSJed Brown /* B = speye(nx^3) + ht*DSG; */ 10649566063dSJacob Faibussowitsch PetscCall(MatScale(user->DSG,user->ht)); 10659566063dSJacob Faibussowitsch PetscCall(MatShift(user->DSG,1.0)); 1066c4762a1bSJed Brown 1067c4762a1bSJed Brown /* Now solve for ytrue */ 10689566063dSJacob Faibussowitsch PetscCall(StateMatInvMult(user->Js,user->q,user->ytrue)); 1069c4762a1bSJed Brown 1070c4762a1bSJed Brown /* Initial guess y0 for state given u0 */ 1071c4762a1bSJed Brown 1072c4762a1bSJed Brown /* First compute Av_u = Av*exp(-u) */ 10739566063dSJacob Faibussowitsch PetscCall(VecSet(user->uwork,0)); 10749566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->uwork,-1.0,user->u)); /* Note: user->u */ 10759566063dSJacob Faibussowitsch PetscCall(VecExp(user->uwork)); 10769566063dSJacob Faibussowitsch PetscCall(MatMult(user->Av,user->uwork,user->Av_u)); 1077c4762a1bSJed Brown 1078c4762a1bSJed Brown /* Next form DSG = Div*S*Grad */ 10799566063dSJacob Faibussowitsch PetscCall(MatCopy(user->Div,user->Divwork,SAME_NONZERO_PATTERN)); 10809566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(user->Divwork,NULL,user->Av_u)); 1081c4762a1bSJed Brown if (user->dsg_formed) { 10829566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(user->DSG)); 1083c4762a1bSJed Brown } else { 10849566063dSJacob Faibussowitsch PetscCall(MatMatMult(user->Div,user->Grad,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&user->DSG)); 1085c20d7725SJed Brown 1086c4762a1bSJed Brown user->dsg_formed = PETSC_TRUE; 1087c4762a1bSJed Brown } 1088c4762a1bSJed Brown /* B = speye(nx^3) + ht*DSG; */ 10899566063dSJacob Faibussowitsch PetscCall(MatScale(user->DSG,user->ht)); 10909566063dSJacob Faibussowitsch PetscCall(MatShift(user->DSG,1.0)); 1091c4762a1bSJed Brown 1092c4762a1bSJed Brown /* Now solve for y */ 10939566063dSJacob Faibussowitsch PetscCall(StateMatInvMult(user->Js,user->q,user->y)); 1094c4762a1bSJed Brown 1095c4762a1bSJed Brown /* Construct projection matrix Q, a block diagonal matrix consisting of nt copies of Qblock along the diagonal */ 10969566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&user->Qblock)); 10979566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Qblock,PETSC_DECIDE,PETSC_DECIDE,user->ndata,n)); 10989566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Qblock)); 10999566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->Qblock,8,NULL,8,NULL)); 11009566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->Qblock,8,NULL)); 1101c4762a1bSJed Brown 1102c4762a1bSJed Brown for (i=0; i<user->mx; i++) { 1103c4762a1bSJed Brown x[i] = h*(i+0.5); 1104c4762a1bSJed Brown y[i] = h*(i+0.5); 1105c4762a1bSJed Brown z[i] = h*(i+0.5); 1106c4762a1bSJed Brown } 1107c4762a1bSJed Brown 11089566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->Qblock,&istart,&iend)); 1109c4762a1bSJed Brown nx = user->mx; ny = user->mx; nz = user->mx; 1110c4762a1bSJed Brown for (i=istart; i<iend; i++) { 1111c4762a1bSJed Brown xri = xr[i]; 1112c4762a1bSJed Brown im = 0; 1113c4762a1bSJed Brown xim = x[im]; 1114c4762a1bSJed Brown while (xri>xim && im<nx) { 1115c4762a1bSJed Brown im = im+1; 1116c4762a1bSJed Brown xim = x[im]; 1117c4762a1bSJed Brown } 1118c4762a1bSJed Brown indx1 = im-1; 1119c4762a1bSJed Brown indx2 = im; 1120c4762a1bSJed Brown dx1 = xri - x[indx1]; 1121c4762a1bSJed Brown dx2 = x[indx2] - xri; 1122c4762a1bSJed Brown 1123c4762a1bSJed Brown yri = yr[i]; 1124c4762a1bSJed Brown im = 0; 1125c4762a1bSJed Brown yim = y[im]; 1126c4762a1bSJed Brown while (yri>yim && im<ny) { 1127c4762a1bSJed Brown im = im+1; 1128c4762a1bSJed Brown yim = y[im]; 1129c4762a1bSJed Brown } 1130c4762a1bSJed Brown indy1 = im-1; 1131c4762a1bSJed Brown indy2 = im; 1132c4762a1bSJed Brown dy1 = yri - y[indy1]; 1133c4762a1bSJed Brown dy2 = y[indy2] - yri; 1134c4762a1bSJed Brown 1135c4762a1bSJed Brown zri = zr[i]; 1136c4762a1bSJed Brown im = 0; 1137c4762a1bSJed Brown zim = z[im]; 1138c4762a1bSJed Brown while (zri>zim && im<nz) { 1139c4762a1bSJed Brown im = im+1; 1140c4762a1bSJed Brown zim = z[im]; 1141c4762a1bSJed Brown } 1142c4762a1bSJed Brown indz1 = im-1; 1143c4762a1bSJed Brown indz2 = im; 1144c4762a1bSJed Brown dz1 = zri - z[indz1]; 1145c4762a1bSJed Brown dz2 = z[indz2] - zri; 1146c4762a1bSJed Brown 1147c4762a1bSJed Brown Dx = x[indx2] - x[indx1]; 1148c4762a1bSJed Brown Dy = y[indy2] - y[indy1]; 1149c4762a1bSJed Brown Dz = z[indz2] - z[indz1]; 1150c4762a1bSJed Brown 1151c4762a1bSJed Brown j = indx1 + indy1*nx + indz1*nx*ny; 1152c4762a1bSJed Brown v = (1-dx1/Dx)*(1-dy1/Dy)*(1-dz1/Dz); 11539566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES)); 1154c4762a1bSJed Brown 1155c4762a1bSJed Brown j = indx1 + indy1*nx + indz2*nx*ny; 1156c4762a1bSJed Brown v = (1-dx1/Dx)*(1-dy1/Dy)*(1-dz2/Dz); 11579566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES)); 1158c4762a1bSJed Brown 1159c4762a1bSJed Brown j = indx1 + indy2*nx + indz1*nx*ny; 1160c4762a1bSJed Brown v = (1-dx1/Dx)*(1-dy2/Dy)*(1-dz1/Dz); 11619566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES)); 1162c4762a1bSJed Brown 1163c4762a1bSJed Brown j = indx1 + indy2*nx + indz2*nx*ny; 1164c4762a1bSJed Brown v = (1-dx1/Dx)*(1-dy2/Dy)*(1-dz2/Dz); 11659566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES)); 1166c4762a1bSJed Brown 1167c4762a1bSJed Brown j = indx2 + indy1*nx + indz1*nx*ny; 1168c4762a1bSJed Brown v = (1-dx2/Dx)*(1-dy1/Dy)*(1-dz1/Dz); 11699566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES)); 1170c4762a1bSJed Brown 1171c4762a1bSJed Brown j = indx2 + indy1*nx + indz2*nx*ny; 1172c4762a1bSJed Brown v = (1-dx2/Dx)*(1-dy1/Dy)*(1-dz2/Dz); 11739566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES)); 1174c4762a1bSJed Brown 1175c4762a1bSJed Brown j = indx2 + indy2*nx + indz1*nx*ny; 1176c4762a1bSJed Brown v = (1-dx2/Dx)*(1-dy2/Dy)*(1-dz1/Dz); 11779566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES)); 1178c4762a1bSJed Brown 1179c4762a1bSJed Brown j = indx2 + indy2*nx + indz2*nx*ny; 1180c4762a1bSJed Brown v = (1-dx2/Dx)*(1-dy2/Dy)*(1-dz2/Dz); 11819566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES)); 1182c4762a1bSJed Brown } 11839566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Qblock,MAT_FINAL_ASSEMBLY)); 11849566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Qblock,MAT_FINAL_ASSEMBLY)); 1185c4762a1bSJed Brown 11869566063dSJacob Faibussowitsch PetscCall(MatTranspose(user->Qblock,MAT_INITIAL_MATRIX,&user->QblockT)); 11879566063dSJacob Faibussowitsch PetscCall(MatTranspose(user->L,MAT_INITIAL_MATRIX,&user->LT)); 1188c4762a1bSJed Brown 1189c4762a1bSJed Brown /* Add noise to the measurement data */ 11909566063dSJacob Faibussowitsch PetscCall(VecSet(user->ywork,1.0)); 11919566063dSJacob Faibussowitsch PetscCall(VecAYPX(user->ywork,user->noise,user->ytrue)); 11929566063dSJacob Faibussowitsch PetscCall(Scatter_i(user->ywork,user->yiwork,user->yi_scatter,user->nt)); 1193c4762a1bSJed Brown for (j=0; j<user->ns; j++) { 1194c4762a1bSJed Brown i = user->sample_times[j]; 11959566063dSJacob Faibussowitsch PetscCall(MatMult(user->Qblock,user->yiwork[i],user->di[j])); 1196c4762a1bSJed Brown } 11979566063dSJacob Faibussowitsch PetscCall(Gather_i(user->d,user->di,user->di_scatter,user->ns)); 1198c4762a1bSJed Brown 1199c4762a1bSJed Brown /* Now that initial conditions have been set, let the user pass tolerance options to the KSP solver */ 12009566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(user->solver)); 12019566063dSJacob Faibussowitsch PetscCall(PetscFree(x)); 12029566063dSJacob Faibussowitsch PetscCall(PetscFree(y)); 12039566063dSJacob Faibussowitsch PetscCall(PetscFree(z)); 1204c4762a1bSJed Brown PetscFunctionReturn(0); 1205c4762a1bSJed Brown } 1206c4762a1bSJed Brown 1207c4762a1bSJed Brown PetscErrorCode ParabolicDestroy(AppCtx *user) 1208c4762a1bSJed Brown { 1209c4762a1bSJed Brown PetscInt i; 1210c4762a1bSJed Brown 1211c4762a1bSJed Brown PetscFunctionBegin; 12129566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Qblock)); 12139566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->QblockT)); 12149566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Div)); 12159566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Divwork)); 12169566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Grad)); 12179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Av)); 12189566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Avwork)); 12199566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->AvT)); 12209566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->DSG)); 12219566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->L)); 12229566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->LT)); 12239566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&user->solver)); 12249566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Js)); 12259566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Jd)); 12269566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->JsInv)); 12279566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->JsBlock)); 12289566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->JsBlockPrec)); 12299566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->u)); 12309566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->uwork)); 12319566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->utrue)); 12329566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->y)); 12339566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->ywork)); 12349566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->ytrue)); 12359566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->nt,&user->yi)); 12369566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->nt,&user->yiwork)); 12379566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->ns,&user->di)); 12389566063dSJacob Faibussowitsch PetscCall(PetscFree(user->yi)); 12399566063dSJacob Faibussowitsch PetscCall(PetscFree(user->yiwork)); 12409566063dSJacob Faibussowitsch PetscCall(PetscFree(user->di)); 12419566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->c)); 12429566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->cwork)); 12439566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->ur)); 12449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->q)); 12459566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->d)); 12469566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->dwork)); 12479566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->lwork)); 12489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->S)); 12499566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->Swork)); 12509566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->Av_u)); 12519566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->Twork)); 12529566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->Rwork)); 12539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->js_diag)); 12549566063dSJacob Faibussowitsch PetscCall(ISDestroy(&user->s_is)); 12559566063dSJacob Faibussowitsch PetscCall(ISDestroy(&user->d_is)); 12569566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->state_scatter)); 12579566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->design_scatter)); 1258c4762a1bSJed Brown for (i=0; i<user->nt; i++) { 12599566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->yi_scatter[i])); 1260c4762a1bSJed Brown } 1261c4762a1bSJed Brown for (i=0; i<user->ns; i++) { 12629566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->di_scatter[i])); 1263c4762a1bSJed Brown } 12649566063dSJacob Faibussowitsch PetscCall(PetscFree(user->yi_scatter)); 12659566063dSJacob Faibussowitsch PetscCall(PetscFree(user->di_scatter)); 12669566063dSJacob Faibussowitsch PetscCall(PetscFree(user->sample_times)); 1267c4762a1bSJed Brown PetscFunctionReturn(0); 1268c4762a1bSJed Brown } 1269c4762a1bSJed Brown 1270c4762a1bSJed Brown PetscErrorCode ParabolicMonitor(Tao tao, void *ptr) 1271c4762a1bSJed Brown { 1272c4762a1bSJed Brown Vec X; 1273c4762a1bSJed Brown PetscReal unorm,ynorm; 1274c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 1275c4762a1bSJed Brown 1276c4762a1bSJed Brown PetscFunctionBegin; 12779566063dSJacob Faibussowitsch PetscCall(TaoGetSolution(tao,&X)); 12789566063dSJacob Faibussowitsch PetscCall(Scatter(X,user->ywork,user->state_scatter,user->uwork,user->design_scatter)); 12799566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->ywork,-1.0,user->ytrue)); 12809566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->uwork,-1.0,user->utrue)); 12819566063dSJacob Faibussowitsch PetscCall(VecNorm(user->uwork,NORM_2,&unorm)); 12829566063dSJacob Faibussowitsch PetscCall(VecNorm(user->ywork,NORM_2,&ynorm)); 12839566063dSJacob Faibussowitsch PetscCall(PetscPrintf(MPI_COMM_WORLD, "||u-ut||=%g ||y-yt||=%g\n",(double)unorm,(double)ynorm)); 1284c4762a1bSJed Brown PetscFunctionReturn(0); 1285c4762a1bSJed Brown } 1286c4762a1bSJed Brown 1287c4762a1bSJed Brown /*TEST 1288c4762a1bSJed Brown 1289c4762a1bSJed Brown build: 1290c4762a1bSJed Brown requires: !complex 1291c4762a1bSJed Brown 1292c4762a1bSJed Brown test: 1293c4762a1bSJed Brown args: -tao_cmonitor -tao_type lcl -ns 1 -tao_gatol 1.e-4 -ksp_max_it 30 1294c4762a1bSJed Brown requires: !single 1295c4762a1bSJed Brown 1296c4762a1bSJed Brown TEST*/ 1297