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 109*327415f7SBarry Smith PetscFunctionBeginUser; 1109566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char*)0,help)); 111c4762a1bSJed Brown user.mx = 8; 112d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"parabolic example",NULL); 1139566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mx","Number of grid points in each direction","",user.mx,&user.mx,NULL)); 114c4762a1bSJed Brown user.nt = 8; 1159566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-nt","Number of time steps","",user.nt,&user.nt,NULL)); 116c4762a1bSJed Brown user.ndata = 64; 1179566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ndata","Numbers of data points per sample","",user.ndata,&user.ndata,NULL)); 118c4762a1bSJed Brown user.ns = 8; 1199566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ns","Number of samples","",user.ns,&user.ns,NULL)); 120c4762a1bSJed Brown user.alpha = 1.0; 1219566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-alpha","Regularization parameter","",user.alpha,&user.alpha,NULL)); 122c4762a1bSJed Brown user.beta = 0.01; 1239566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-beta","Weight attributed to ||u||^2 in regularization functional","",user.beta,&user.beta,NULL)); 124c4762a1bSJed Brown user.noise = 0.01; 1259566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-noise","Amount of noise to add to data","",user.noise,&user.noise,NULL)); 1269566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ntests","Number of times to repeat TaoSolve","",ntests,&ntests,NULL)); 127d0609cedSBarry Smith PetscOptionsEnd(); 128c4762a1bSJed Brown 129c4762a1bSJed Brown user.m = user.mx*user.mx*user.mx; /* number of constraints per time step */ 130c4762a1bSJed Brown user.n = user.m*(user.nt+1); /* number of variables */ 131c4762a1bSJed Brown user.ht = (PetscReal)1/user.nt; 132c4762a1bSJed Brown 1339566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&user.u)); 1349566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&user.y)); 1359566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&user.c)); 1369566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user.u,PETSC_DECIDE,user.n-user.m*user.nt)); 1379566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user.y,PETSC_DECIDE,user.m*user.nt)); 1389566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user.c,PETSC_DECIDE,user.m*user.nt)); 1399566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user.u)); 1409566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user.y)); 1419566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user.c)); 142c4762a1bSJed Brown 143c4762a1bSJed Brown /* Create scatters for reduced spaces. 144c4762a1bSJed Brown If the state vector y and design vector u are partitioned as 145c4762a1bSJed Brown [y_1; y_2; ...; y_np] and [u_1; u_2; ...; u_np] (with np = # of processors), 146c4762a1bSJed Brown then the solution vector x is organized as 147c4762a1bSJed Brown [y_1; u_1; y_2; u_2; ...; y_np; u_np]. 148c4762a1bSJed Brown The index sets user.s_is and user.d_is correspond to the indices of the 149c4762a1bSJed Brown state and design variables owned by the current processor. 150c4762a1bSJed Brown */ 1519566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&x)); 152c4762a1bSJed Brown 1539566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user.y,&lo,&hi)); 1549566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user.u,&lo2,&hi2)); 155c4762a1bSJed Brown 1569566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_allstate)); 1579566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+lo2,1,&user.s_is)); 158c4762a1bSJed Brown 1599566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,hi2-lo2,lo2,1,&is_alldesign)); 1609566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,hi2-lo2,hi+lo2,1,&user.d_is)); 161c4762a1bSJed Brown 1629566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x,hi-lo+hi2-lo2,user.n)); 1639566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x)); 164c4762a1bSJed Brown 1659566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(x,user.s_is,user.y,is_allstate,&user.state_scatter)); 1669566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(x,user.d_is,user.u,is_alldesign,&user.design_scatter)); 1679566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_alldesign)); 1689566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_allstate)); 169c4762a1bSJed Brown 170c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 1719566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD,&tao)); 1729566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao,TAOLCL)); 173c4762a1bSJed Brown 174c4762a1bSJed Brown /* Set up initial vectors and matrices */ 1759566063dSJacob Faibussowitsch PetscCall(ParabolicInitialize(&user)); 176c4762a1bSJed Brown 1779566063dSJacob Faibussowitsch PetscCall(Gather(x,user.y,user.state_scatter,user.u,user.design_scatter)); 1789566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&x0)); 1799566063dSJacob Faibussowitsch PetscCall(VecCopy(x,x0)); 180c4762a1bSJed Brown 181c4762a1bSJed Brown /* Set solution vector with an initial guess */ 1829566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao,x)); 1839566063dSJacob Faibussowitsch PetscCall(TaoSetObjective(tao, FormFunction, &user)); 1849566063dSJacob Faibussowitsch PetscCall(TaoSetGradient(tao, NULL, FormGradient, &user)); 1859566063dSJacob Faibussowitsch PetscCall(TaoSetConstraintsRoutine(tao, user.c, FormConstraints, &user)); 186c4762a1bSJed Brown 1879566063dSJacob Faibussowitsch PetscCall(TaoSetJacobianStateRoutine(tao, user.Js, user.JsBlockPrec, user.JsInv, FormJacobianState, &user)); 1889566063dSJacob Faibussowitsch PetscCall(TaoSetJacobianDesignRoutine(tao, user.Jd, FormJacobianDesign, &user)); 189c4762a1bSJed Brown 1909566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 1919566063dSJacob Faibussowitsch PetscCall(TaoSetStateDesignIS(tao,user.s_is,user.d_is)); 192c4762a1bSJed Brown 193c4762a1bSJed Brown /* SOLVE THE APPLICATION */ 1949566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("Trials",&stages[0])); 1959566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stages[0])); 196c4762a1bSJed Brown user.ksp_its_initial = user.ksp_its; 197c4762a1bSJed Brown ksp_old = user.ksp_its; 198c4762a1bSJed Brown for (i=0; i<ntests; i++) { 1999566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 20063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"KSP Iterations = %" PetscInt_FMT "\n",user.ksp_its-ksp_old)); 2019566063dSJacob Faibussowitsch PetscCall(VecCopy(x0,x)); 2029566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao,x)); 203c4762a1bSJed Brown } 2049566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 2059566063dSJacob Faibussowitsch PetscCall(PetscBarrier((PetscObject)x)); 2069566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"KSP iterations within initialization: ")); 20763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%" PetscInt_FMT "\n",user.ksp_its_initial)); 20863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Total KSP iterations over %" PetscInt_FMT " trial(s): ",ntests)); 20963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%" PetscInt_FMT "\n",user.ksp_its)); 2109566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"KSP iterations per trial: ")); 21163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%" PetscInt_FMT "\n",(user.ksp_its-user.ksp_its_initial)/ntests)); 212c4762a1bSJed Brown 2139566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 2149566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 2159566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x0)); 2169566063dSJacob Faibussowitsch PetscCall(ParabolicDestroy(&user)); 217c4762a1bSJed Brown 2189566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 219b122ec5aSJacob Faibussowitsch return 0; 220c4762a1bSJed Brown } 221c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 222c4762a1bSJed Brown /* 223c4762a1bSJed Brown dwork = Qy - d 224c4762a1bSJed Brown lwork = L*(u-ur) 225c4762a1bSJed Brown f = 1/2 * (dwork.dork + alpha*lwork.lwork) 226c4762a1bSJed Brown */ 227c4762a1bSJed Brown PetscErrorCode FormFunction(Tao tao,Vec X,PetscReal *f,void *ptr) 228c4762a1bSJed Brown { 229c4762a1bSJed Brown PetscReal d1=0,d2=0; 230c4762a1bSJed Brown PetscInt i,j; 231c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 232c4762a1bSJed Brown 233c4762a1bSJed Brown PetscFunctionBegin; 2349566063dSJacob Faibussowitsch PetscCall(Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter)); 2359566063dSJacob Faibussowitsch PetscCall(Scatter_i(user->y,user->yi,user->yi_scatter,user->nt)); 236c4762a1bSJed Brown for (j=0; j<user->ns; j++) { 237c4762a1bSJed Brown i = user->sample_times[j]; 2389566063dSJacob Faibussowitsch PetscCall(MatMult(user->Qblock,user->yi[i],user->di[j])); 239c4762a1bSJed Brown } 2409566063dSJacob Faibussowitsch PetscCall(Gather_i(user->dwork,user->di,user->di_scatter,user->ns)); 2419566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->dwork,-1.0,user->d)); 2429566063dSJacob Faibussowitsch PetscCall(VecDot(user->dwork,user->dwork,&d1)); 243c4762a1bSJed Brown 2449566063dSJacob Faibussowitsch PetscCall(VecWAXPY(user->uwork,-1.0,user->ur,user->u)); 2459566063dSJacob Faibussowitsch PetscCall(MatMult(user->L,user->uwork,user->lwork)); 2469566063dSJacob Faibussowitsch PetscCall(VecDot(user->lwork,user->lwork,&d2)); 247c4762a1bSJed Brown 248c4762a1bSJed Brown *f = 0.5 * (d1 + user->alpha*d2); 249c4762a1bSJed Brown PetscFunctionReturn(0); 250c4762a1bSJed Brown } 251c4762a1bSJed Brown 252c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 253c4762a1bSJed Brown /* 254c4762a1bSJed Brown state: g_s = Q' *(Qy - d) 255c4762a1bSJed Brown design: g_d = alpha*L'*L*(u-ur) 256c4762a1bSJed Brown */ 257c4762a1bSJed Brown PetscErrorCode FormGradient(Tao tao,Vec X,Vec G,void *ptr) 258c4762a1bSJed Brown { 259c4762a1bSJed Brown PetscInt i,j; 260c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 261c4762a1bSJed Brown 262c4762a1bSJed Brown PetscFunctionBegin; 2639566063dSJacob Faibussowitsch PetscCall(Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter)); 2649566063dSJacob Faibussowitsch PetscCall(Scatter_i(user->y,user->yi,user->yi_scatter,user->nt)); 265c4762a1bSJed Brown for (j=0; j<user->ns; j++) { 266c4762a1bSJed Brown i = user->sample_times[j]; 2679566063dSJacob Faibussowitsch PetscCall(MatMult(user->Qblock,user->yi[i],user->di[j])); 268c4762a1bSJed Brown } 2699566063dSJacob Faibussowitsch PetscCall(Gather_i(user->dwork,user->di,user->di_scatter,user->ns)); 2709566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->dwork,-1.0,user->d)); 2719566063dSJacob Faibussowitsch PetscCall(Scatter_i(user->dwork,user->di,user->di_scatter,user->ns)); 2729566063dSJacob Faibussowitsch PetscCall(VecSet(user->ywork,0.0)); 2739566063dSJacob Faibussowitsch PetscCall(Scatter_i(user->ywork,user->yiwork,user->yi_scatter,user->nt)); 274c4762a1bSJed Brown for (j=0; j<user->ns; j++) { 275c4762a1bSJed Brown i = user->sample_times[j]; 2769566063dSJacob Faibussowitsch PetscCall(MatMult(user->QblockT,user->di[j],user->yiwork[i])); 277c4762a1bSJed Brown } 2789566063dSJacob Faibussowitsch PetscCall(Gather_i(user->ywork,user->yiwork,user->yi_scatter,user->nt)); 279c4762a1bSJed Brown 2809566063dSJacob Faibussowitsch PetscCall(VecWAXPY(user->uwork,-1.0,user->ur,user->u)); 2819566063dSJacob Faibussowitsch PetscCall(MatMult(user->L,user->uwork,user->lwork)); 2829566063dSJacob Faibussowitsch PetscCall(MatMult(user->LT,user->lwork,user->uwork)); 2839566063dSJacob Faibussowitsch PetscCall(VecScale(user->uwork, user->alpha)); 2849566063dSJacob Faibussowitsch PetscCall(Gather(G,user->ywork,user->state_scatter,user->uwork,user->design_scatter)); 285c4762a1bSJed Brown PetscFunctionReturn(0); 286c4762a1bSJed Brown } 287c4762a1bSJed Brown 288c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr) 289c4762a1bSJed Brown { 290c4762a1bSJed Brown PetscReal d1,d2; 291c4762a1bSJed Brown PetscInt i,j; 292c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 293c4762a1bSJed Brown 294c4762a1bSJed Brown PetscFunctionBegin; 2959566063dSJacob Faibussowitsch PetscCall(Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter)); 2969566063dSJacob Faibussowitsch PetscCall(Scatter_i(user->y,user->yi,user->yi_scatter,user->nt)); 297c4762a1bSJed Brown for (j=0; j<user->ns; j++) { 298c4762a1bSJed Brown i = user->sample_times[j]; 2999566063dSJacob Faibussowitsch PetscCall(MatMult(user->Qblock,user->yi[i],user->di[j])); 300c4762a1bSJed Brown } 3019566063dSJacob Faibussowitsch PetscCall(Gather_i(user->dwork,user->di,user->di_scatter,user->ns)); 3029566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->dwork,-1.0,user->d)); 3039566063dSJacob Faibussowitsch PetscCall(VecDot(user->dwork,user->dwork,&d1)); 3049566063dSJacob Faibussowitsch PetscCall(Scatter_i(user->dwork,user->di,user->di_scatter,user->ns)); 3059566063dSJacob Faibussowitsch PetscCall(VecSet(user->ywork,0.0)); 3069566063dSJacob Faibussowitsch PetscCall(Scatter_i(user->ywork,user->yiwork,user->yi_scatter,user->nt)); 307c4762a1bSJed Brown for (j=0; j<user->ns; j++) { 308c4762a1bSJed Brown i = user->sample_times[j]; 3099566063dSJacob Faibussowitsch PetscCall(MatMult(user->QblockT,user->di[j],user->yiwork[i])); 310c4762a1bSJed Brown } 3119566063dSJacob Faibussowitsch PetscCall(Gather_i(user->ywork,user->yiwork,user->yi_scatter,user->nt)); 312c4762a1bSJed Brown 3139566063dSJacob Faibussowitsch PetscCall(VecWAXPY(user->uwork,-1.0,user->ur,user->u)); 3149566063dSJacob Faibussowitsch PetscCall(MatMult(user->L,user->uwork,user->lwork)); 3159566063dSJacob Faibussowitsch PetscCall(VecDot(user->lwork,user->lwork,&d2)); 3169566063dSJacob Faibussowitsch PetscCall(MatMult(user->LT,user->lwork,user->uwork)); 3179566063dSJacob Faibussowitsch PetscCall(VecScale(user->uwork, user->alpha)); 318c4762a1bSJed Brown *f = 0.5 * (d1 + user->alpha*d2); 319c4762a1bSJed Brown 3209566063dSJacob Faibussowitsch PetscCall(Gather(G,user->ywork,user->state_scatter,user->uwork,user->design_scatter)); 321c4762a1bSJed Brown PetscFunctionReturn(0); 322c4762a1bSJed Brown } 323c4762a1bSJed Brown 324c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 325c4762a1bSJed Brown /* A 326c4762a1bSJed Brown MatShell object 327c4762a1bSJed Brown */ 328c4762a1bSJed Brown PetscErrorCode FormJacobianState(Tao tao, Vec X, Mat J, Mat JPre, Mat JInv, void *ptr) 329c4762a1bSJed Brown { 330c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 331c4762a1bSJed Brown 332c4762a1bSJed Brown PetscFunctionBegin; 3339566063dSJacob Faibussowitsch PetscCall(Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter)); 3349566063dSJacob Faibussowitsch PetscCall(VecSet(user->uwork,0)); 3359566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->uwork,-1.0,user->u)); 3369566063dSJacob Faibussowitsch PetscCall(VecExp(user->uwork)); 3379566063dSJacob Faibussowitsch PetscCall(MatMult(user->Av,user->uwork,user->Av_u)); 3389566063dSJacob Faibussowitsch PetscCall(VecCopy(user->Av_u,user->Swork)); 3399566063dSJacob Faibussowitsch PetscCall(VecReciprocal(user->Swork)); 3409566063dSJacob Faibussowitsch PetscCall(MatCopy(user->Div,user->Divwork,SAME_NONZERO_PATTERN)); 3419566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(user->Divwork,NULL,user->Swork)); 342c4762a1bSJed Brown if (user->dsg_formed) { 3439566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(user->DSG)); 344c4762a1bSJed Brown } else { 3459566063dSJacob Faibussowitsch PetscCall(MatMatMult(user->Divwork,user->Grad,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&user->DSG)); 346c4762a1bSJed Brown user->dsg_formed = PETSC_TRUE; 347c4762a1bSJed Brown } 348c4762a1bSJed Brown 349c4762a1bSJed Brown /* B = speye(nx^3) + ht*DSG; */ 3509566063dSJacob Faibussowitsch PetscCall(MatScale(user->DSG,user->ht)); 3519566063dSJacob Faibussowitsch PetscCall(MatShift(user->DSG,1.0)); 352c4762a1bSJed Brown PetscFunctionReturn(0); 353c4762a1bSJed Brown } 354c4762a1bSJed Brown 355c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 356c4762a1bSJed Brown /* B */ 357c4762a1bSJed Brown PetscErrorCode FormJacobianDesign(Tao tao, Vec X, Mat J, void *ptr) 358c4762a1bSJed Brown { 359c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 360c4762a1bSJed Brown 361c4762a1bSJed Brown PetscFunctionBegin; 3629566063dSJacob Faibussowitsch PetscCall(Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter)); 363c4762a1bSJed Brown PetscFunctionReturn(0); 364c4762a1bSJed Brown } 365c4762a1bSJed Brown 366c4762a1bSJed Brown PetscErrorCode StateMatMult(Mat J_shell, Vec X, Vec Y) 367c4762a1bSJed Brown { 368c4762a1bSJed Brown PetscInt i; 369c4762a1bSJed Brown AppCtx *user; 370c4762a1bSJed Brown 371c4762a1bSJed Brown PetscFunctionBegin; 3729566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell,&user)); 3739566063dSJacob Faibussowitsch PetscCall(Scatter_i(X,user->yi,user->yi_scatter,user->nt)); 3749566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock,user->yi[0],user->yiwork[0])); 375c4762a1bSJed Brown for (i=1; i<user->nt; i++) { 3769566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock,user->yi[i],user->yiwork[i])); 3779566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->yiwork[i],-1.0,user->yi[i-1])); 378c4762a1bSJed Brown } 3799566063dSJacob Faibussowitsch PetscCall(Gather_i(Y,user->yiwork,user->yi_scatter,user->nt)); 380c4762a1bSJed Brown PetscFunctionReturn(0); 381c4762a1bSJed Brown } 382c4762a1bSJed Brown 383c4762a1bSJed Brown PetscErrorCode StateMatMultTranspose(Mat J_shell, Vec X, Vec Y) 384c4762a1bSJed Brown { 385c4762a1bSJed Brown PetscInt i; 386c4762a1bSJed Brown AppCtx *user; 387c4762a1bSJed Brown 388c4762a1bSJed Brown PetscFunctionBegin; 3899566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell,&user)); 3909566063dSJacob Faibussowitsch PetscCall(Scatter_i(X,user->yi,user->yi_scatter,user->nt)); 391c4762a1bSJed Brown for (i=0; i<user->nt-1; i++) { 3929566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock,user->yi[i],user->yiwork[i])); 3939566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->yiwork[i],-1.0,user->yi[i+1])); 394c4762a1bSJed Brown } 395c4762a1bSJed Brown i = user->nt-1; 3969566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock,user->yi[i],user->yiwork[i])); 3979566063dSJacob Faibussowitsch PetscCall(Gather_i(Y,user->yiwork,user->yi_scatter,user->nt)); 398c4762a1bSJed Brown PetscFunctionReturn(0); 399c4762a1bSJed Brown } 400c4762a1bSJed Brown 401c4762a1bSJed Brown PetscErrorCode StateMatBlockMult(Mat J_shell, Vec X, Vec Y) 402c4762a1bSJed Brown { 403c4762a1bSJed Brown AppCtx *user; 404c4762a1bSJed Brown 405c4762a1bSJed Brown PetscFunctionBegin; 4069566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell,&user)); 4079566063dSJacob Faibussowitsch PetscCall(MatMult(user->Grad,X,user->Swork)); 4089566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(user->Swork,user->Swork,user->Av_u)); 4099566063dSJacob Faibussowitsch PetscCall(MatMult(user->Div,user->Swork,Y)); 4109566063dSJacob Faibussowitsch PetscCall(VecAYPX(Y,user->ht,X)); 411c4762a1bSJed Brown PetscFunctionReturn(0); 412c4762a1bSJed Brown } 413c4762a1bSJed Brown 414c4762a1bSJed Brown PetscErrorCode DesignMatMult(Mat J_shell, Vec X, Vec Y) 415c4762a1bSJed Brown { 416c4762a1bSJed Brown PetscInt i; 417c4762a1bSJed Brown AppCtx *user; 418c4762a1bSJed Brown 419c4762a1bSJed Brown PetscFunctionBegin; 4209566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell,&user)); 421c4762a1bSJed Brown 422c4762a1bSJed Brown /* sdiag(1./v) */ 4239566063dSJacob Faibussowitsch PetscCall(VecSet(user->uwork,0)); 4249566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->uwork,-1.0,user->u)); 4259566063dSJacob Faibussowitsch PetscCall(VecExp(user->uwork)); 426c4762a1bSJed Brown 427c4762a1bSJed Brown /* sdiag(1./((Av*(1./v)).^2)) */ 4289566063dSJacob Faibussowitsch PetscCall(MatMult(user->Av,user->uwork,user->Swork)); 4299566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Swork,user->Swork,user->Swork)); 4309566063dSJacob Faibussowitsch PetscCall(VecReciprocal(user->Swork)); 431c4762a1bSJed Brown 432c4762a1bSJed Brown /* (Av * (sdiag(1./v) * b)) */ 4339566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->uwork,user->uwork,X)); 4349566063dSJacob Faibussowitsch PetscCall(MatMult(user->Av,user->uwork,user->Twork)); 435c4762a1bSJed Brown 436c4762a1bSJed Brown /* (sdiag(1./((Av*(1./v)).^2)) * (Av * (sdiag(1./v) * b))) */ 4379566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Swork,user->Twork,user->Swork)); 438c4762a1bSJed Brown 4399566063dSJacob Faibussowitsch PetscCall(Scatter_i(user->y,user->yi,user->yi_scatter,user->nt)); 440c4762a1bSJed Brown for (i=0; i<user->nt; i++) { 441c4762a1bSJed Brown /* (sdiag(Grad*y(:,i)) */ 4429566063dSJacob Faibussowitsch PetscCall(MatMult(user->Grad,user->yi[i],user->Twork)); 443c4762a1bSJed Brown 444c4762a1bSJed Brown /* ht * Div * (sdiag(Grad*y(:,i)) * (sdiag(1./((Av*(1./v)).^2)) * (Av * (sdiag(1./v) * b)))) */ 4459566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Twork,user->Twork,user->Swork)); 4469566063dSJacob Faibussowitsch PetscCall(MatMult(user->Div,user->Twork,user->yiwork[i])); 4479566063dSJacob Faibussowitsch PetscCall(VecScale(user->yiwork[i],user->ht)); 448c4762a1bSJed Brown } 4499566063dSJacob Faibussowitsch PetscCall(Gather_i(Y,user->yiwork,user->yi_scatter,user->nt)); 450c4762a1bSJed Brown 451c4762a1bSJed Brown PetscFunctionReturn(0); 452c4762a1bSJed Brown } 453c4762a1bSJed Brown 454c4762a1bSJed Brown PetscErrorCode DesignMatMultTranspose(Mat J_shell, Vec X, Vec Y) 455c4762a1bSJed Brown { 456c4762a1bSJed Brown PetscInt i; 457c4762a1bSJed Brown AppCtx *user; 458c4762a1bSJed Brown 459c4762a1bSJed Brown PetscFunctionBegin; 4609566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell,&user)); 461c4762a1bSJed Brown 462c4762a1bSJed Brown /* sdiag(1./((Av*(1./v)).^2)) */ 4639566063dSJacob Faibussowitsch PetscCall(VecSet(user->uwork,0)); 4649566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->uwork,-1.0,user->u)); 4659566063dSJacob Faibussowitsch PetscCall(VecExp(user->uwork)); 4669566063dSJacob Faibussowitsch PetscCall(MatMult(user->Av,user->uwork,user->Rwork)); 4679566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Rwork,user->Rwork,user->Rwork)); 4689566063dSJacob Faibussowitsch PetscCall(VecReciprocal(user->Rwork)); 469c4762a1bSJed Brown 4709566063dSJacob Faibussowitsch PetscCall(VecSet(Y,0.0)); 4719566063dSJacob Faibussowitsch PetscCall(Scatter_i(user->y,user->yi,user->yi_scatter,user->nt)); 4729566063dSJacob Faibussowitsch PetscCall(Scatter_i(X,user->yiwork,user->yi_scatter,user->nt)); 473c4762a1bSJed Brown for (i=0; i<user->nt; i++) { 474c4762a1bSJed Brown /* (Div' * b(:,i)) */ 4759566063dSJacob Faibussowitsch PetscCall(MatMult(user->Grad,user->yiwork[i],user->Swork)); 476c4762a1bSJed Brown 477c4762a1bSJed Brown /* sdiag(Grad*y(:,i)) */ 4789566063dSJacob Faibussowitsch PetscCall(MatMult(user->Grad,user->yi[i],user->Twork)); 479c4762a1bSJed Brown 480c4762a1bSJed Brown /* (sdiag(Grad*y(:,i)) * (Div' * b(:,i))) */ 4819566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Twork,user->Swork,user->Twork)); 482c4762a1bSJed Brown 483c4762a1bSJed Brown /* (sdiag(1./((Av*(1./v)).^2)) * (sdiag(Grad*y(:,i)) * (Div' * b(:,i)))) */ 4849566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->Twork,user->Rwork,user->Twork)); 485c4762a1bSJed Brown 486c4762a1bSJed Brown /* (Av' * (sdiag(1./((Av*(1./v)).^2)) * (sdiag(Grad*y(:,i)) * (Div' * b(:,i))))) */ 4879566063dSJacob Faibussowitsch PetscCall(MatMult(user->AvT,user->Twork,user->yiwork[i])); 488c4762a1bSJed Brown 489c4762a1bSJed Brown /* sdiag(1./v) * (Av' * (sdiag(1./((Av*(1./v)).^2)) * (sdiag(Grad*y(:,i)) * (Div' * b(:,i))))) */ 4909566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(user->yiwork[i],user->uwork,user->yiwork[i])); 4919566063dSJacob Faibussowitsch PetscCall(VecAXPY(Y,user->ht,user->yiwork[i])); 492c4762a1bSJed Brown } 493c4762a1bSJed Brown PetscFunctionReturn(0); 494c4762a1bSJed Brown } 495c4762a1bSJed Brown 496c4762a1bSJed Brown PetscErrorCode StateMatBlockPrecMult(PC PC_shell, Vec X, Vec Y) 497c4762a1bSJed Brown { 498c4762a1bSJed Brown AppCtx *user; 499c4762a1bSJed Brown 500c4762a1bSJed Brown PetscFunctionBegin; 5019566063dSJacob Faibussowitsch PetscCall(PCShellGetContext(PC_shell,&user)); 502c4762a1bSJed Brown 503c4762a1bSJed Brown if (user->dsg_formed) { 5049566063dSJacob Faibussowitsch PetscCall(MatSOR(user->DSG,X,1.0,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_SYMMETRIC_SWEEP),0.0,1,1,Y)); 505c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"DSG not formed"); 506c4762a1bSJed Brown PetscFunctionReturn(0); 507c4762a1bSJed Brown } 508c4762a1bSJed Brown 509c4762a1bSJed Brown PetscErrorCode StateMatInvMult(Mat J_shell, Vec X, Vec Y) 510c4762a1bSJed Brown { 511c4762a1bSJed Brown AppCtx *user; 512c4762a1bSJed Brown PetscInt its,i; 513c4762a1bSJed Brown 514c4762a1bSJed Brown PetscFunctionBegin; 5159566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell,&user)); 516c4762a1bSJed Brown 517c4762a1bSJed Brown if (Y == user->ytrue) { 518c4762a1bSJed Brown /* First solve is done with true solution to set up problem */ 5199566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(user->solver,1e-8,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT)); 520c4762a1bSJed Brown } else { 5219566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(user->solver,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT)); 522c4762a1bSJed Brown } 523c4762a1bSJed Brown 5249566063dSJacob Faibussowitsch PetscCall(Scatter_i(X,user->yi,user->yi_scatter,user->nt)); 5259566063dSJacob Faibussowitsch PetscCall(KSPSolve(user->solver,user->yi[0],user->yiwork[0])); 5269566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(user->solver,&its)); 527c4762a1bSJed Brown user->ksp_its = user->ksp_its + its; 528c4762a1bSJed Brown 529c4762a1bSJed Brown for (i=1; i<user->nt; i++) { 5309566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->yi[i],1.0,user->yiwork[i-1])); 5319566063dSJacob Faibussowitsch PetscCall(KSPSolve(user->solver,user->yi[i],user->yiwork[i])); 5329566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(user->solver,&its)); 533c4762a1bSJed Brown user->ksp_its = user->ksp_its + its; 534c4762a1bSJed Brown } 5359566063dSJacob Faibussowitsch PetscCall(Gather_i(Y,user->yiwork,user->yi_scatter,user->nt)); 536c4762a1bSJed Brown PetscFunctionReturn(0); 537c4762a1bSJed Brown } 538c4762a1bSJed Brown 539c4762a1bSJed Brown PetscErrorCode StateMatInvTransposeMult(Mat J_shell, Vec X, Vec Y) 540c4762a1bSJed Brown { 541c4762a1bSJed Brown AppCtx *user; 542c4762a1bSJed Brown PetscInt its,i; 543c4762a1bSJed Brown 544c4762a1bSJed Brown PetscFunctionBegin; 5459566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell,&user)); 546c4762a1bSJed Brown 5479566063dSJacob Faibussowitsch PetscCall(Scatter_i(X,user->yi,user->yi_scatter,user->nt)); 548c4762a1bSJed Brown 549c4762a1bSJed Brown i = user->nt - 1; 5509566063dSJacob Faibussowitsch PetscCall(KSPSolve(user->solver,user->yi[i],user->yiwork[i])); 551c4762a1bSJed Brown 5529566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(user->solver,&its)); 553c4762a1bSJed Brown user->ksp_its = user->ksp_its + its; 554c4762a1bSJed Brown 555c4762a1bSJed Brown for (i=user->nt-2; i>=0; i--) { 5569566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->yi[i],1.0,user->yiwork[i+1])); 5579566063dSJacob Faibussowitsch PetscCall(KSPSolve(user->solver,user->yi[i],user->yiwork[i])); 558c4762a1bSJed Brown 5599566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(user->solver,&its)); 560c4762a1bSJed Brown user->ksp_its = user->ksp_its + its; 561c4762a1bSJed Brown 562c4762a1bSJed Brown } 563c4762a1bSJed Brown 5649566063dSJacob Faibussowitsch PetscCall(Gather_i(Y,user->yiwork,user->yi_scatter,user->nt)); 565c4762a1bSJed Brown PetscFunctionReturn(0); 566c4762a1bSJed Brown } 567c4762a1bSJed Brown 568c4762a1bSJed Brown PetscErrorCode StateMatDuplicate(Mat J_shell, MatDuplicateOption opt, Mat *new_shell) 569c4762a1bSJed Brown { 570c4762a1bSJed Brown AppCtx *user; 571c4762a1bSJed Brown 572c4762a1bSJed Brown PetscFunctionBegin; 5739566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell,&user)); 574c4762a1bSJed Brown 5759566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,new_shell)); 5769566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*new_shell,MATOP_MULT,(void(*)(void))StateMatMult)); 5779566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*new_shell,MATOP_DUPLICATE,(void(*)(void))StateMatDuplicate)); 5789566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*new_shell,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatMultTranspose)); 5799566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*new_shell,MATOP_GET_DIAGONAL,(void(*)(void))StateMatGetDiagonal)); 580c4762a1bSJed Brown PetscFunctionReturn(0); 581c4762a1bSJed Brown } 582c4762a1bSJed Brown 583c4762a1bSJed Brown PetscErrorCode StateMatGetDiagonal(Mat J_shell, Vec X) 584c4762a1bSJed Brown { 585c4762a1bSJed Brown AppCtx *user; 586c4762a1bSJed Brown 587c4762a1bSJed Brown PetscFunctionBegin; 5889566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J_shell,&user)); 5899566063dSJacob Faibussowitsch PetscCall(VecCopy(user->js_diag,X)); 590c4762a1bSJed Brown PetscFunctionReturn(0); 591c4762a1bSJed Brown 592c4762a1bSJed Brown } 593c4762a1bSJed Brown 594c4762a1bSJed Brown PetscErrorCode FormConstraints(Tao tao, Vec X, Vec C, void *ptr) 595c4762a1bSJed Brown { 596c4762a1bSJed Brown /* con = Ay - q, A = [B 0 0 ... 0; 597c4762a1bSJed Brown -I B 0 ... 0; 598c4762a1bSJed Brown 0 -I B ... 0; 599c4762a1bSJed Brown ... ; 600c4762a1bSJed Brown 0 ... -I B] 601c4762a1bSJed Brown B = ht * Div * Sigma * Grad + eye */ 602c4762a1bSJed Brown PetscInt i; 603c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 604c4762a1bSJed Brown 605c4762a1bSJed Brown PetscFunctionBegin; 6069566063dSJacob Faibussowitsch PetscCall(Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter)); 6079566063dSJacob Faibussowitsch PetscCall(Scatter_i(user->y,user->yi,user->yi_scatter,user->nt)); 6089566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock,user->yi[0],user->yiwork[0])); 609c4762a1bSJed Brown for (i=1; i<user->nt; i++) { 6109566063dSJacob Faibussowitsch PetscCall(MatMult(user->JsBlock,user->yi[i],user->yiwork[i])); 6119566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->yiwork[i],-1.0,user->yi[i-1])); 612c4762a1bSJed Brown } 6139566063dSJacob Faibussowitsch PetscCall(Gather_i(C,user->yiwork,user->yi_scatter,user->nt)); 6149566063dSJacob Faibussowitsch PetscCall(VecAXPY(C,-1.0,user->q)); 615c4762a1bSJed Brown PetscFunctionReturn(0); 616c4762a1bSJed Brown } 617c4762a1bSJed Brown 618c4762a1bSJed Brown PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat) 619c4762a1bSJed Brown { 620c4762a1bSJed Brown PetscFunctionBegin; 6219566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(s_scat,x,state,INSERT_VALUES,SCATTER_FORWARD)); 6229566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(s_scat,x,state,INSERT_VALUES,SCATTER_FORWARD)); 6239566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(d_scat,x,design,INSERT_VALUES,SCATTER_FORWARD)); 6249566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(d_scat,x,design,INSERT_VALUES,SCATTER_FORWARD)); 625c4762a1bSJed Brown PetscFunctionReturn(0); 626c4762a1bSJed Brown } 627c4762a1bSJed Brown 628c4762a1bSJed Brown PetscErrorCode Scatter_i(Vec y, Vec *yi, VecScatter *scat, PetscInt nt) 629c4762a1bSJed Brown { 630c4762a1bSJed Brown PetscInt i; 631c4762a1bSJed Brown 632c4762a1bSJed Brown PetscFunctionBegin; 633c4762a1bSJed Brown for (i=0; i<nt; i++) { 6349566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scat[i],y,yi[i],INSERT_VALUES,SCATTER_FORWARD)); 6359566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scat[i],y,yi[i],INSERT_VALUES,SCATTER_FORWARD)); 636c4762a1bSJed Brown } 637c4762a1bSJed Brown PetscFunctionReturn(0); 638c4762a1bSJed Brown } 639c4762a1bSJed Brown 640c4762a1bSJed Brown PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat) 641c4762a1bSJed Brown { 642c4762a1bSJed Brown PetscFunctionBegin; 6439566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(s_scat,state,x,INSERT_VALUES,SCATTER_REVERSE)); 6449566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(s_scat,state,x,INSERT_VALUES,SCATTER_REVERSE)); 6459566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(d_scat,design,x,INSERT_VALUES,SCATTER_REVERSE)); 6469566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(d_scat,design,x,INSERT_VALUES,SCATTER_REVERSE)); 647c4762a1bSJed Brown PetscFunctionReturn(0); 648c4762a1bSJed Brown } 649c4762a1bSJed Brown 650c4762a1bSJed Brown PetscErrorCode Gather_i(Vec y, Vec *yi, VecScatter *scat, PetscInt nt) 651c4762a1bSJed Brown { 652c4762a1bSJed Brown PetscInt i; 653c4762a1bSJed Brown 654c4762a1bSJed Brown PetscFunctionBegin; 655c4762a1bSJed Brown for (i=0; i<nt; i++) { 6569566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scat[i],yi[i],y,INSERT_VALUES,SCATTER_REVERSE)); 6579566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scat[i],yi[i],y,INSERT_VALUES,SCATTER_REVERSE)); 658c4762a1bSJed Brown } 659c4762a1bSJed Brown PetscFunctionReturn(0); 660c4762a1bSJed Brown } 661c4762a1bSJed Brown 662c4762a1bSJed Brown PetscErrorCode ParabolicInitialize(AppCtx *user) 663c4762a1bSJed Brown { 664c4762a1bSJed Brown PetscInt m,n,i,j,k,linear_index,istart,iend,iblock,lo,hi,lo2,hi2; 665c4762a1bSJed Brown Vec XX,YY,ZZ,XXwork,YYwork,ZZwork,UTwork,yi,di,bc; 666c4762a1bSJed Brown PetscReal *x, *y, *z; 667c4762a1bSJed Brown PetscReal h,stime; 668c4762a1bSJed Brown PetscScalar hinv,neg_hinv,half = 0.5,sqrt_beta; 669c4762a1bSJed Brown PetscInt im,indx1,indx2,indy1,indy2,indz1,indz2,nx,ny,nz; 670c4762a1bSJed Brown PetscReal xri,yri,zri,xim,yim,zim,dx1,dx2,dy1,dy2,dz1,dz2,Dx,Dy,Dz; 671c4762a1bSJed Brown PetscScalar v,vx,vy,vz; 672c4762a1bSJed Brown IS is_from_y,is_to_yi,is_from_d,is_to_di; 673c4762a1bSJed Brown /* Data locations */ 674c4762a1bSJed Brown PetscScalar xr[64] = {0.4970, 0.8498, 0.7814, 0.6268, 0.7782, 0.6402, 0.3617, 0.3160, 675c4762a1bSJed Brown 0.3610, 0.5298, 0.6987, 0.3331, 0.7962, 0.5596, 0.3866, 0.6774, 676c4762a1bSJed Brown 0.5407, 0.4518, 0.6702, 0.6061, 0.7580, 0.8997, 0.5198, 0.8326, 677c4762a1bSJed Brown 0.2138, 0.9198, 0.3000, 0.2833, 0.8288, 0.7076, 0.1820, 0.0728, 678c4762a1bSJed Brown 0.8447, 0.2367, 0.3239, 0.6413, 0.3114, 0.4731, 0.1192, 0.9273, 679c4762a1bSJed Brown 0.5724, 0.4331, 0.5136, 0.3547, 0.4413, 0.2602, 0.5698, 0.7278, 680c4762a1bSJed Brown 0.5261, 0.6230, 0.2454, 0.3948, 0.7479, 0.6582, 0.4660, 0.5594, 681c4762a1bSJed Brown 0.7574, 0.1143, 0.5900, 0.1065, 0.4260, 0.3294, 0.8276, 0.0756}; 682c4762a1bSJed Brown 683c4762a1bSJed Brown PetscScalar yr[64] = {0.7345, 0.9120, 0.9288, 0.7528, 0.4463, 0.4985, 0.2497, 0.6256, 684c4762a1bSJed Brown 0.3425, 0.9026, 0.6983, 0.4230, 0.7140, 0.2970, 0.4474, 0.8792, 685c4762a1bSJed Brown 0.6604, 0.2485, 0.7968, 0.6127, 0.1796, 0.2437, 0.5938, 0.6137, 686c4762a1bSJed Brown 0.3867, 0.5658, 0.4575, 0.1009, 0.0863, 0.3361, 0.0738, 0.3985, 687c4762a1bSJed Brown 0.6602, 0.1437, 0.0934, 0.5983, 0.5950, 0.0763, 0.0768, 0.2288, 688c4762a1bSJed Brown 0.5761, 0.1129, 0.3841, 0.6150, 0.6904, 0.6686, 0.1361, 0.4601, 689c4762a1bSJed Brown 0.4491, 0.3716, 0.1969, 0.6537, 0.6743, 0.6991, 0.4811, 0.5480, 690c4762a1bSJed Brown 0.1684, 0.4569, 0.6889, 0.8437, 0.3015, 0.2854, 0.8199, 0.2658}; 691c4762a1bSJed Brown 692c4762a1bSJed Brown PetscScalar zr[64] = {0.7668, 0.8573, 0.2654, 0.2719, 0.1060, 0.1311, 0.6232, 0.2295, 693c4762a1bSJed Brown 0.8009, 0.2147, 0.2119, 0.9325, 0.4473, 0.3600, 0.3374, 0.3819, 694c4762a1bSJed Brown 0.4066, 0.5801, 0.1673, 0.0959, 0.4638, 0.8236, 0.8800, 0.2939, 695c4762a1bSJed Brown 0.2028, 0.8262, 0.2706, 0.6276, 0.9085, 0.6443, 0.8241, 0.0712, 696c4762a1bSJed Brown 0.1824, 0.7789, 0.4389, 0.8415, 0.7055, 0.6639, 0.3653, 0.2078, 697c4762a1bSJed Brown 0.1987, 0.2297, 0.4321, 0.8115, 0.4915, 0.7764, 0.4657, 0.4627, 698c4762a1bSJed Brown 0.4569, 0.4232, 0.8514, 0.0674, 0.3227, 0.1055, 0.6690, 0.6313, 699c4762a1bSJed Brown 0.9226, 0.5461, 0.4126, 0.2364, 0.6096, 0.7042, 0.3914, 0.0711}; 700c4762a1bSJed Brown 701c4762a1bSJed Brown PetscFunctionBegin; 7029566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->mx,&x)); 7039566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->mx,&y)); 7049566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->mx,&z)); 705c4762a1bSJed Brown user->jformed = PETSC_FALSE; 706c4762a1bSJed Brown user->dsg_formed = PETSC_FALSE; 707c4762a1bSJed Brown 708c4762a1bSJed Brown n = user->mx * user->mx * user->mx; 709c4762a1bSJed Brown m = 3 * user->mx * user->mx * (user->mx-1); 710c4762a1bSJed Brown sqrt_beta = PetscSqrtScalar(user->beta); 711c4762a1bSJed Brown 712c4762a1bSJed Brown user->ksp_its = 0; 713c4762a1bSJed Brown user->ksp_its_initial = 0; 714c4762a1bSJed Brown 715c4762a1bSJed Brown stime = (PetscReal)user->nt/user->ns; 7169566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->ns,&user->sample_times)); 717c4762a1bSJed Brown for (i=0; i<user->ns; i++) { 718c4762a1bSJed Brown user->sample_times[i] = (PetscInt)(stime*((PetscReal)i+1.0)-0.5); 719c4762a1bSJed Brown } 720c4762a1bSJed Brown 7219566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&XX)); 7229566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&user->q)); 7239566063dSJacob Faibussowitsch PetscCall(VecSetSizes(XX,PETSC_DECIDE,n)); 7249566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->q,PETSC_DECIDE,n*user->nt)); 7259566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(XX)); 7269566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->q)); 727c4762a1bSJed Brown 7289566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX,&YY)); 7299566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX,&ZZ)); 7309566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX,&XXwork)); 7319566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX,&YYwork)); 7329566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX,&ZZwork)); 7339566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX,&UTwork)); 7349566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX,&user->utrue)); 7359566063dSJacob Faibussowitsch PetscCall(VecDuplicate(XX,&bc)); 736c4762a1bSJed Brown 737c4762a1bSJed Brown /* Generate 3D grid, and collect ns (1<=ns<=8) right-hand-side vectors into user->q */ 738c4762a1bSJed Brown h = 1.0/user->mx; 739c4762a1bSJed Brown hinv = user->mx; 740c4762a1bSJed Brown neg_hinv = -hinv; 741c4762a1bSJed Brown 7429566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(XX,&istart,&iend)); 743c4762a1bSJed Brown for (linear_index=istart; linear_index<iend; linear_index++) { 744c4762a1bSJed Brown i = linear_index % user->mx; 745c4762a1bSJed Brown j = ((linear_index-i)/user->mx) % user->mx; 746c4762a1bSJed Brown k = ((linear_index-i)/user->mx-j) / user->mx; 747c4762a1bSJed Brown vx = h*(i+0.5); 748c4762a1bSJed Brown vy = h*(j+0.5); 749c4762a1bSJed Brown vz = h*(k+0.5); 7509566063dSJacob Faibussowitsch PetscCall(VecSetValues(XX,1,&linear_index,&vx,INSERT_VALUES)); 7519566063dSJacob Faibussowitsch PetscCall(VecSetValues(YY,1,&linear_index,&vy,INSERT_VALUES)); 7529566063dSJacob Faibussowitsch PetscCall(VecSetValues(ZZ,1,&linear_index,&vz,INSERT_VALUES)); 753c4762a1bSJed Brown if ((vx<0.6) && (vx>0.4) && (vy<0.6) && (vy>0.4) && (vy<0.6) && (vz<0.6) && (vz>0.4)) { 754c4762a1bSJed Brown v = 1000.0; 7559566063dSJacob Faibussowitsch PetscCall(VecSetValues(bc,1,&linear_index,&v,INSERT_VALUES)); 756c4762a1bSJed Brown } 757c4762a1bSJed Brown } 758c4762a1bSJed Brown 7599566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(XX)); 7609566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(XX)); 7619566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(YY)); 7629566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(YY)); 7639566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(ZZ)); 7649566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(ZZ)); 7659566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(bc)); 7669566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(bc)); 767c4762a1bSJed Brown 768c4762a1bSJed Brown /* Compute true parameter function 769c4762a1bSJed Brown ut = 0.5 + exp(-10*((x-0.5)^2+(y-0.5)^2+(z-0.5)^2)) */ 7709566063dSJacob Faibussowitsch PetscCall(VecCopy(XX,XXwork)); 7719566063dSJacob Faibussowitsch PetscCall(VecCopy(YY,YYwork)); 7729566063dSJacob Faibussowitsch PetscCall(VecCopy(ZZ,ZZwork)); 773c4762a1bSJed Brown 7749566063dSJacob Faibussowitsch PetscCall(VecShift(XXwork,-0.5)); 7759566063dSJacob Faibussowitsch PetscCall(VecShift(YYwork,-0.5)); 7769566063dSJacob Faibussowitsch PetscCall(VecShift(ZZwork,-0.5)); 777c4762a1bSJed Brown 7789566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(XXwork,XXwork,XXwork)); 7799566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(YYwork,YYwork,YYwork)); 7809566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(ZZwork,ZZwork,ZZwork)); 781c4762a1bSJed Brown 7829566063dSJacob Faibussowitsch PetscCall(VecCopy(XXwork,user->utrue)); 7839566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->utrue,1.0,YYwork)); 7849566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->utrue,1.0,ZZwork)); 7859566063dSJacob Faibussowitsch PetscCall(VecScale(user->utrue,-10.0)); 7869566063dSJacob Faibussowitsch PetscCall(VecExp(user->utrue)); 7879566063dSJacob Faibussowitsch PetscCall(VecShift(user->utrue,0.5)); 788c4762a1bSJed Brown 7899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&XX)); 7909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&YY)); 7919566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ZZ)); 7929566063dSJacob Faibussowitsch PetscCall(VecDestroy(&XXwork)); 7939566063dSJacob Faibussowitsch PetscCall(VecDestroy(&YYwork)); 7949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ZZwork)); 7959566063dSJacob Faibussowitsch PetscCall(VecDestroy(&UTwork)); 796c4762a1bSJed Brown 797c4762a1bSJed Brown /* Initial guess and reference model */ 7989566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->utrue,&user->ur)); 7999566063dSJacob Faibussowitsch PetscCall(VecSet(user->ur,0.0)); 800c4762a1bSJed Brown 801c4762a1bSJed Brown /* Generate Grad matrix */ 8029566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&user->Grad)); 8039566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Grad,PETSC_DECIDE,PETSC_DECIDE,m,n)); 8049566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Grad)); 8059566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->Grad,2,NULL,2,NULL)); 8069566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->Grad,2,NULL)); 8079566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->Grad,&istart,&iend)); 808c4762a1bSJed Brown 809c4762a1bSJed Brown for (i=istart; i<iend; i++) { 810c4762a1bSJed Brown if (i<m/3) { 811c4762a1bSJed Brown iblock = i / (user->mx-1); 812c4762a1bSJed Brown j = iblock*user->mx + (i % (user->mx-1)); 8139566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad,1,&i,1,&j,&neg_hinv,INSERT_VALUES)); 814c4762a1bSJed Brown j = j+1; 8159566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad,1,&i,1,&j,&hinv,INSERT_VALUES)); 816c4762a1bSJed Brown } 817c4762a1bSJed Brown if (i>=m/3 && i<2*m/3) { 818c4762a1bSJed Brown iblock = (i-m/3) / (user->mx*(user->mx-1)); 819c4762a1bSJed Brown j = iblock*user->mx*user->mx + ((i-m/3) % (user->mx*(user->mx-1))); 8209566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad,1,&i,1,&j,&neg_hinv,INSERT_VALUES)); 821c4762a1bSJed Brown j = j + user->mx; 8229566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad,1,&i,1,&j,&hinv,INSERT_VALUES)); 823c4762a1bSJed Brown } 824c4762a1bSJed Brown if (i>=2*m/3) { 825c4762a1bSJed Brown j = i-2*m/3; 8269566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad,1,&i,1,&j,&neg_hinv,INSERT_VALUES)); 827c4762a1bSJed Brown j = j + user->mx*user->mx; 8289566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Grad,1,&i,1,&j,&hinv,INSERT_VALUES)); 829c4762a1bSJed Brown } 830c4762a1bSJed Brown } 831c4762a1bSJed Brown 8329566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Grad,MAT_FINAL_ASSEMBLY)); 8339566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Grad,MAT_FINAL_ASSEMBLY)); 834c4762a1bSJed Brown 835c4762a1bSJed Brown /* Generate arithmetic averaging matrix Av */ 8369566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&user->Av)); 8379566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Av,PETSC_DECIDE,PETSC_DECIDE,m,n)); 8389566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Av)); 8399566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->Av,2,NULL,2,NULL)); 8409566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->Av,2,NULL)); 8419566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->Av,&istart,&iend)); 842c4762a1bSJed Brown 843c4762a1bSJed Brown for (i=istart; i<iend; i++) { 844c4762a1bSJed Brown if (i<m/3) { 845c4762a1bSJed Brown iblock = i / (user->mx-1); 846c4762a1bSJed Brown j = iblock*user->mx + (i % (user->mx-1)); 8479566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES)); 848c4762a1bSJed Brown j = j+1; 8499566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES)); 850c4762a1bSJed Brown } 851c4762a1bSJed Brown if (i>=m/3 && i<2*m/3) { 852c4762a1bSJed Brown iblock = (i-m/3) / (user->mx*(user->mx-1)); 853c4762a1bSJed Brown j = iblock*user->mx*user->mx + ((i-m/3) % (user->mx*(user->mx-1))); 8549566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES)); 855c4762a1bSJed Brown j = j + user->mx; 8569566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES)); 857c4762a1bSJed Brown } 858c4762a1bSJed Brown if (i>=2*m/3) { 859c4762a1bSJed Brown j = i-2*m/3; 8609566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES)); 861c4762a1bSJed Brown j = j + user->mx*user->mx; 8629566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES)); 863c4762a1bSJed Brown } 864c4762a1bSJed Brown } 865c4762a1bSJed Brown 8669566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Av,MAT_FINAL_ASSEMBLY)); 8679566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Av,MAT_FINAL_ASSEMBLY)); 868c4762a1bSJed Brown 869c4762a1bSJed Brown /* Generate transpose of averaging matrix Av */ 8709566063dSJacob Faibussowitsch PetscCall(MatTranspose(user->Av,MAT_INITIAL_MATRIX,&user->AvT)); 871c4762a1bSJed Brown 8729566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&user->L)); 8739566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->L,PETSC_DECIDE,PETSC_DECIDE,m+n,n)); 8749566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->L)); 8759566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->L,2,NULL,2,NULL)); 8769566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->L,2,NULL)); 8779566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->L,&istart,&iend)); 878c4762a1bSJed Brown 879c4762a1bSJed Brown for (i=istart; i<iend; i++) { 880c4762a1bSJed Brown if (i<m/3) { 881c4762a1bSJed Brown iblock = i / (user->mx-1); 882c4762a1bSJed Brown j = iblock*user->mx + (i % (user->mx-1)); 8839566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L,1,&i,1,&j,&neg_hinv,INSERT_VALUES)); 884c4762a1bSJed Brown j = j+1; 8859566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L,1,&i,1,&j,&hinv,INSERT_VALUES)); 886c4762a1bSJed Brown } 887c4762a1bSJed Brown if (i>=m/3 && i<2*m/3) { 888c4762a1bSJed Brown iblock = (i-m/3) / (user->mx*(user->mx-1)); 889c4762a1bSJed Brown j = iblock*user->mx*user->mx + ((i-m/3) % (user->mx*(user->mx-1))); 8909566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L,1,&i,1,&j,&neg_hinv,INSERT_VALUES)); 891c4762a1bSJed Brown j = j + user->mx; 8929566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L,1,&i,1,&j,&hinv,INSERT_VALUES)); 893c4762a1bSJed Brown } 894c4762a1bSJed Brown if (i>=2*m/3 && i<m) { 895c4762a1bSJed Brown j = i-2*m/3; 8969566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L,1,&i,1,&j,&neg_hinv,INSERT_VALUES)); 897c4762a1bSJed Brown j = j + user->mx*user->mx; 8989566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L,1,&i,1,&j,&hinv,INSERT_VALUES)); 899c4762a1bSJed Brown } 900c4762a1bSJed Brown if (i>=m) { 901c4762a1bSJed Brown j = i - m; 9029566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->L,1,&i,1,&j,&sqrt_beta,INSERT_VALUES)); 903c4762a1bSJed Brown } 904c4762a1bSJed Brown } 905c4762a1bSJed Brown 9069566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->L,MAT_FINAL_ASSEMBLY)); 9079566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->L,MAT_FINAL_ASSEMBLY)); 908c4762a1bSJed Brown 9099566063dSJacob Faibussowitsch PetscCall(MatScale(user->L,PetscPowScalar(h,1.5))); 910c4762a1bSJed Brown 911c4762a1bSJed Brown /* Generate Div matrix */ 9129566063dSJacob Faibussowitsch PetscCall(MatTranspose(user->Grad,MAT_INITIAL_MATRIX,&user->Div)); 913c4762a1bSJed Brown 914c4762a1bSJed Brown /* Build work vectors and matrices */ 9159566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&user->S)); 9169566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->S, PETSC_DECIDE, user->mx*user->mx*(user->mx-1)*3)); 9179566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->S)); 918c4762a1bSJed Brown 9199566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&user->lwork)); 9209566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->lwork,PETSC_DECIDE,m+user->mx*user->mx*user->mx)); 9219566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->lwork)); 922c4762a1bSJed Brown 9239566063dSJacob Faibussowitsch PetscCall(MatDuplicate(user->Div,MAT_SHARE_NONZERO_PATTERN,&user->Divwork)); 9249566063dSJacob Faibussowitsch PetscCall(MatDuplicate(user->Av,MAT_SHARE_NONZERO_PATTERN,&user->Avwork)); 925c4762a1bSJed Brown 9269566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&user->d)); 9279566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->d,PETSC_DECIDE,user->ndata*user->nt)); 9289566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->d)); 929c4762a1bSJed Brown 9309566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->S,&user->Swork)); 9319566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->S,&user->Av_u)); 9329566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->S,&user->Twork)); 9339566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->S,&user->Rwork)); 9349566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->y,&user->ywork)); 9359566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->u,&user->uwork)); 9369566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->u,&user->js_diag)); 9379566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->c,&user->cwork)); 9389566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->d,&user->dwork)); 939c4762a1bSJed Brown 940c4762a1bSJed Brown /* Create matrix-free shell user->Js for computing A*x */ 9419566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m*user->nt,user->m*user->nt,user,&user->Js)); 9429566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->Js,MATOP_MULT,(void(*)(void))StateMatMult)); 9439566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->Js,MATOP_DUPLICATE,(void(*)(void))StateMatDuplicate)); 9449566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->Js,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatMultTranspose)); 9459566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->Js,MATOP_GET_DIAGONAL,(void(*)(void))StateMatGetDiagonal)); 946c4762a1bSJed Brown 947c4762a1bSJed Brown /* Diagonal blocks of user->Js */ 9489566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,&user->JsBlock)); 9499566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->JsBlock,MATOP_MULT,(void(*)(void))StateMatBlockMult)); 950c4762a1bSJed Brown /* Blocks are symmetric */ 9519566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->JsBlock,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatBlockMult)); 952c4762a1bSJed Brown 953c4762a1bSJed Brown /* Create a matrix-free shell user->JsBlockPrec for computing (U+D)\D*(L+D)\x, where JsBlock = L+D+U, 954c4762a1bSJed Brown D is diagonal, L is strictly lower triangular, and U is strictly upper triangular. 955c4762a1bSJed Brown This is an SSOR preconditioner for user->JsBlock. */ 9569566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,&user->JsBlockPrec)); 9579566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->JsBlockPrec,MATOP_MULT,(void(*)(void))StateMatBlockPrecMult)); 958c4762a1bSJed Brown /* JsBlockPrec is symmetric */ 9599566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->JsBlockPrec,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatBlockPrecMult)); 9609566063dSJacob Faibussowitsch PetscCall(MatSetOption(user->JsBlockPrec,MAT_SYMMETRY_ETERNAL,PETSC_TRUE)); 961c4762a1bSJed Brown 962c4762a1bSJed Brown /* Create a matrix-free shell user->Jd for computing B*x */ 9639566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m*user->nt,user->m,user,&user->Jd)); 9649566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->Jd,MATOP_MULT,(void(*)(void))DesignMatMult)); 9659566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->Jd,MATOP_MULT_TRANSPOSE,(void(*)(void))DesignMatMultTranspose)); 966c4762a1bSJed Brown 967c4762a1bSJed Brown /* User-defined routines for computing user->Js\x and user->Js^T\x*/ 9689566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m*user->nt,user->m*user->nt,user,&user->JsInv)); 9699566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->JsInv,MATOP_MULT,(void(*)(void))StateMatInvMult)); 9709566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(user->JsInv,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatInvTransposeMult)); 971c4762a1bSJed Brown 972c4762a1bSJed Brown /* Solver options and tolerances */ 9739566063dSJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_WORLD,&user->solver)); 9749566063dSJacob Faibussowitsch PetscCall(KSPSetType(user->solver,KSPCG)); 9759566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(user->solver,user->JsBlock,user->JsBlockPrec)); 9769566063dSJacob Faibussowitsch PetscCall(KSPSetInitialGuessNonzero(user->solver,PETSC_FALSE)); 9779566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(user->solver,1e-4,1e-20,1e3,500)); 9789566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(user->solver)); 9799566063dSJacob Faibussowitsch PetscCall(KSPGetPC(user->solver,&user->prec)); 9809566063dSJacob Faibussowitsch PetscCall(PCSetType(user->prec,PCSHELL)); 981c4762a1bSJed Brown 9829566063dSJacob Faibussowitsch PetscCall(PCShellSetApply(user->prec,StateMatBlockPrecMult)); 9839566063dSJacob Faibussowitsch PetscCall(PCShellSetApplyTranspose(user->prec,StateMatBlockPrecMult)); 9849566063dSJacob Faibussowitsch PetscCall(PCShellSetContext(user->prec,user)); 985c4762a1bSJed Brown 986c4762a1bSJed Brown /* Create scatter from y to y_1,y_2,...,y_nt */ 9879566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->nt*user->m,&user->yi_scatter)); 9889566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&yi)); 9899566063dSJacob Faibussowitsch PetscCall(VecSetSizes(yi,PETSC_DECIDE,user->mx*user->mx*user->mx)); 9909566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(yi)); 9919566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(yi,user->nt,&user->yi)); 9929566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(yi,user->nt,&user->yiwork)); 993c4762a1bSJed Brown 9949566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->y,&lo2,&hi2)); 995c4762a1bSJed Brown istart = 0; 996c4762a1bSJed Brown for (i=0; i<user->nt; i++) { 9979566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->yi[i],&lo,&hi)); 9989566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_yi)); 9999566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo2+istart,1,&is_from_y)); 10009566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(user->y,is_from_y,user->yi[i],is_to_yi,&user->yi_scatter[i])); 1001c4762a1bSJed Brown istart = istart + hi-lo; 10029566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_to_yi)); 10039566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_from_y)); 1004c4762a1bSJed Brown } 10059566063dSJacob Faibussowitsch PetscCall(VecDestroy(&yi)); 1006c4762a1bSJed Brown 1007c4762a1bSJed Brown /* Create scatter from d to d_1,d_2,...,d_ns */ 10089566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->ns*user->ndata,&user->di_scatter)); 10099566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&di)); 10109566063dSJacob Faibussowitsch PetscCall(VecSetSizes(di,PETSC_DECIDE,user->ndata)); 10119566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(di)); 10129566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(di,user->ns,&user->di)); 10139566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->d,&lo2,&hi2)); 1014c4762a1bSJed Brown istart = 0; 1015c4762a1bSJed Brown for (i=0; i<user->ns; i++) { 10169566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(user->di[i],&lo,&hi)); 10179566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_di)); 10189566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,hi-lo,lo2+istart,1,&is_from_d)); 10199566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(user->d,is_from_d,user->di[i],is_to_di,&user->di_scatter[i])); 1020c4762a1bSJed Brown istart = istart + hi-lo; 10219566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_to_di)); 10229566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_from_d)); 1023c4762a1bSJed Brown } 10249566063dSJacob Faibussowitsch PetscCall(VecDestroy(&di)); 1025c4762a1bSJed Brown 1026c4762a1bSJed Brown /* Assemble RHS of forward problem */ 10279566063dSJacob Faibussowitsch PetscCall(VecCopy(bc,user->yiwork[0])); 1028c4762a1bSJed Brown for (i=1; i<user->nt; i++) { 10299566063dSJacob Faibussowitsch PetscCall(VecSet(user->yiwork[i],0.0)); 1030c4762a1bSJed Brown } 10319566063dSJacob Faibussowitsch PetscCall(Gather_i(user->q,user->yiwork,user->yi_scatter,user->nt)); 10329566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bc)); 1033c4762a1bSJed Brown 1034c4762a1bSJed Brown /* Compute true state function yt given ut */ 10359566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&user->ytrue)); 10369566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->ytrue,PETSC_DECIDE,n*user->nt)); 10379566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->ytrue)); 1038c4762a1bSJed Brown 1039c4762a1bSJed Brown /* First compute Av_u = Av*exp(-u) */ 10409566063dSJacob Faibussowitsch PetscCall(VecSet(user->uwork,0)); 10419566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->uwork,-1.0,user->utrue)); /* Note: user->utrue */ 10429566063dSJacob Faibussowitsch PetscCall(VecExp(user->uwork)); 10439566063dSJacob Faibussowitsch PetscCall(MatMult(user->Av,user->uwork,user->Av_u)); 1044c4762a1bSJed Brown 1045c20d7725SJed Brown /* Symbolic DSG = Div * Grad */ 10469566063dSJacob Faibussowitsch PetscCall(MatProductCreate(user->Div,user->Grad,NULL,&user->DSG)); 10479566063dSJacob Faibussowitsch PetscCall(MatProductSetType(user->DSG,MATPRODUCT_AB)); 10489566063dSJacob Faibussowitsch PetscCall(MatProductSetAlgorithm(user->DSG,"default")); 10499566063dSJacob Faibussowitsch PetscCall(MatProductSetFill(user->DSG,PETSC_DEFAULT)); 10509566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(user->DSG)); 10519566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(user->DSG)); 1052c20d7725SJed Brown 1053c4762a1bSJed Brown user->dsg_formed = PETSC_TRUE; 1054c4762a1bSJed Brown 1055c20d7725SJed Brown /* Next form DSG = Div*Grad */ 10569566063dSJacob Faibussowitsch PetscCall(MatCopy(user->Div,user->Divwork,SAME_NONZERO_PATTERN)); 10579566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(user->Divwork,NULL,user->Av_u)); 1058c4762a1bSJed Brown if (user->dsg_formed) { 10599566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(user->DSG)); 1060c4762a1bSJed Brown } else { 10619566063dSJacob Faibussowitsch PetscCall(MatMatMult(user->Div,user->Grad,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&user->DSG)); 1062c4762a1bSJed Brown user->dsg_formed = PETSC_TRUE; 1063c4762a1bSJed Brown } 1064c4762a1bSJed Brown /* B = speye(nx^3) + ht*DSG; */ 10659566063dSJacob Faibussowitsch PetscCall(MatScale(user->DSG,user->ht)); 10669566063dSJacob Faibussowitsch PetscCall(MatShift(user->DSG,1.0)); 1067c4762a1bSJed Brown 1068c4762a1bSJed Brown /* Now solve for ytrue */ 10699566063dSJacob Faibussowitsch PetscCall(StateMatInvMult(user->Js,user->q,user->ytrue)); 1070c4762a1bSJed Brown 1071c4762a1bSJed Brown /* Initial guess y0 for state given u0 */ 1072c4762a1bSJed Brown 1073c4762a1bSJed Brown /* First compute Av_u = Av*exp(-u) */ 10749566063dSJacob Faibussowitsch PetscCall(VecSet(user->uwork,0)); 10759566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->uwork,-1.0,user->u)); /* Note: user->u */ 10769566063dSJacob Faibussowitsch PetscCall(VecExp(user->uwork)); 10779566063dSJacob Faibussowitsch PetscCall(MatMult(user->Av,user->uwork,user->Av_u)); 1078c4762a1bSJed Brown 1079c4762a1bSJed Brown /* Next form DSG = Div*S*Grad */ 10809566063dSJacob Faibussowitsch PetscCall(MatCopy(user->Div,user->Divwork,SAME_NONZERO_PATTERN)); 10819566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(user->Divwork,NULL,user->Av_u)); 1082c4762a1bSJed Brown if (user->dsg_formed) { 10839566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(user->DSG)); 1084c4762a1bSJed Brown } else { 10859566063dSJacob Faibussowitsch PetscCall(MatMatMult(user->Div,user->Grad,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&user->DSG)); 1086c20d7725SJed Brown 1087c4762a1bSJed Brown user->dsg_formed = PETSC_TRUE; 1088c4762a1bSJed Brown } 1089c4762a1bSJed Brown /* B = speye(nx^3) + ht*DSG; */ 10909566063dSJacob Faibussowitsch PetscCall(MatScale(user->DSG,user->ht)); 10919566063dSJacob Faibussowitsch PetscCall(MatShift(user->DSG,1.0)); 1092c4762a1bSJed Brown 1093c4762a1bSJed Brown /* Now solve for y */ 10949566063dSJacob Faibussowitsch PetscCall(StateMatInvMult(user->Js,user->q,user->y)); 1095c4762a1bSJed Brown 1096c4762a1bSJed Brown /* Construct projection matrix Q, a block diagonal matrix consisting of nt copies of Qblock along the diagonal */ 10979566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&user->Qblock)); 10989566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Qblock,PETSC_DECIDE,PETSC_DECIDE,user->ndata,n)); 10999566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Qblock)); 11009566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(user->Qblock,8,NULL,8,NULL)); 11019566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(user->Qblock,8,NULL)); 1102c4762a1bSJed Brown 1103c4762a1bSJed Brown for (i=0; i<user->mx; i++) { 1104c4762a1bSJed Brown x[i] = h*(i+0.5); 1105c4762a1bSJed Brown y[i] = h*(i+0.5); 1106c4762a1bSJed Brown z[i] = h*(i+0.5); 1107c4762a1bSJed Brown } 1108c4762a1bSJed Brown 11099566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(user->Qblock,&istart,&iend)); 1110c4762a1bSJed Brown nx = user->mx; ny = user->mx; nz = user->mx; 1111c4762a1bSJed Brown for (i=istart; i<iend; i++) { 1112c4762a1bSJed Brown xri = xr[i]; 1113c4762a1bSJed Brown im = 0; 1114c4762a1bSJed Brown xim = x[im]; 1115c4762a1bSJed Brown while (xri>xim && im<nx) { 1116c4762a1bSJed Brown im = im+1; 1117c4762a1bSJed Brown xim = x[im]; 1118c4762a1bSJed Brown } 1119c4762a1bSJed Brown indx1 = im-1; 1120c4762a1bSJed Brown indx2 = im; 1121c4762a1bSJed Brown dx1 = xri - x[indx1]; 1122c4762a1bSJed Brown dx2 = x[indx2] - xri; 1123c4762a1bSJed Brown 1124c4762a1bSJed Brown yri = yr[i]; 1125c4762a1bSJed Brown im = 0; 1126c4762a1bSJed Brown yim = y[im]; 1127c4762a1bSJed Brown while (yri>yim && im<ny) { 1128c4762a1bSJed Brown im = im+1; 1129c4762a1bSJed Brown yim = y[im]; 1130c4762a1bSJed Brown } 1131c4762a1bSJed Brown indy1 = im-1; 1132c4762a1bSJed Brown indy2 = im; 1133c4762a1bSJed Brown dy1 = yri - y[indy1]; 1134c4762a1bSJed Brown dy2 = y[indy2] - yri; 1135c4762a1bSJed Brown 1136c4762a1bSJed Brown zri = zr[i]; 1137c4762a1bSJed Brown im = 0; 1138c4762a1bSJed Brown zim = z[im]; 1139c4762a1bSJed Brown while (zri>zim && im<nz) { 1140c4762a1bSJed Brown im = im+1; 1141c4762a1bSJed Brown zim = z[im]; 1142c4762a1bSJed Brown } 1143c4762a1bSJed Brown indz1 = im-1; 1144c4762a1bSJed Brown indz2 = im; 1145c4762a1bSJed Brown dz1 = zri - z[indz1]; 1146c4762a1bSJed Brown dz2 = z[indz2] - zri; 1147c4762a1bSJed Brown 1148c4762a1bSJed Brown Dx = x[indx2] - x[indx1]; 1149c4762a1bSJed Brown Dy = y[indy2] - y[indy1]; 1150c4762a1bSJed Brown Dz = z[indz2] - z[indz1]; 1151c4762a1bSJed Brown 1152c4762a1bSJed Brown j = indx1 + indy1*nx + indz1*nx*ny; 1153c4762a1bSJed Brown v = (1-dx1/Dx)*(1-dy1/Dy)*(1-dz1/Dz); 11549566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES)); 1155c4762a1bSJed Brown 1156c4762a1bSJed Brown j = indx1 + indy1*nx + indz2*nx*ny; 1157c4762a1bSJed Brown v = (1-dx1/Dx)*(1-dy1/Dy)*(1-dz2/Dz); 11589566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES)); 1159c4762a1bSJed Brown 1160c4762a1bSJed Brown j = indx1 + indy2*nx + indz1*nx*ny; 1161c4762a1bSJed Brown v = (1-dx1/Dx)*(1-dy2/Dy)*(1-dz1/Dz); 11629566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES)); 1163c4762a1bSJed Brown 1164c4762a1bSJed Brown j = indx1 + indy2*nx + indz2*nx*ny; 1165c4762a1bSJed Brown v = (1-dx1/Dx)*(1-dy2/Dy)*(1-dz2/Dz); 11669566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES)); 1167c4762a1bSJed Brown 1168c4762a1bSJed Brown j = indx2 + indy1*nx + indz1*nx*ny; 1169c4762a1bSJed Brown v = (1-dx2/Dx)*(1-dy1/Dy)*(1-dz1/Dz); 11709566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES)); 1171c4762a1bSJed Brown 1172c4762a1bSJed Brown j = indx2 + indy1*nx + indz2*nx*ny; 1173c4762a1bSJed Brown v = (1-dx2/Dx)*(1-dy1/Dy)*(1-dz2/Dz); 11749566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES)); 1175c4762a1bSJed Brown 1176c4762a1bSJed Brown j = indx2 + indy2*nx + indz1*nx*ny; 1177c4762a1bSJed Brown v = (1-dx2/Dx)*(1-dy2/Dy)*(1-dz1/Dz); 11789566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES)); 1179c4762a1bSJed Brown 1180c4762a1bSJed Brown j = indx2 + indy2*nx + indz2*nx*ny; 1181c4762a1bSJed Brown v = (1-dx2/Dx)*(1-dy2/Dy)*(1-dz2/Dz); 11829566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES)); 1183c4762a1bSJed Brown } 11849566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Qblock,MAT_FINAL_ASSEMBLY)); 11859566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Qblock,MAT_FINAL_ASSEMBLY)); 1186c4762a1bSJed Brown 11879566063dSJacob Faibussowitsch PetscCall(MatTranspose(user->Qblock,MAT_INITIAL_MATRIX,&user->QblockT)); 11889566063dSJacob Faibussowitsch PetscCall(MatTranspose(user->L,MAT_INITIAL_MATRIX,&user->LT)); 1189c4762a1bSJed Brown 1190c4762a1bSJed Brown /* Add noise to the measurement data */ 11919566063dSJacob Faibussowitsch PetscCall(VecSet(user->ywork,1.0)); 11929566063dSJacob Faibussowitsch PetscCall(VecAYPX(user->ywork,user->noise,user->ytrue)); 11939566063dSJacob Faibussowitsch PetscCall(Scatter_i(user->ywork,user->yiwork,user->yi_scatter,user->nt)); 1194c4762a1bSJed Brown for (j=0; j<user->ns; j++) { 1195c4762a1bSJed Brown i = user->sample_times[j]; 11969566063dSJacob Faibussowitsch PetscCall(MatMult(user->Qblock,user->yiwork[i],user->di[j])); 1197c4762a1bSJed Brown } 11989566063dSJacob Faibussowitsch PetscCall(Gather_i(user->d,user->di,user->di_scatter,user->ns)); 1199c4762a1bSJed Brown 1200c4762a1bSJed Brown /* Now that initial conditions have been set, let the user pass tolerance options to the KSP solver */ 12019566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(user->solver)); 12029566063dSJacob Faibussowitsch PetscCall(PetscFree(x)); 12039566063dSJacob Faibussowitsch PetscCall(PetscFree(y)); 12049566063dSJacob Faibussowitsch PetscCall(PetscFree(z)); 1205c4762a1bSJed Brown PetscFunctionReturn(0); 1206c4762a1bSJed Brown } 1207c4762a1bSJed Brown 1208c4762a1bSJed Brown PetscErrorCode ParabolicDestroy(AppCtx *user) 1209c4762a1bSJed Brown { 1210c4762a1bSJed Brown PetscInt i; 1211c4762a1bSJed Brown 1212c4762a1bSJed Brown PetscFunctionBegin; 12139566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Qblock)); 12149566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->QblockT)); 12159566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Div)); 12169566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Divwork)); 12179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Grad)); 12189566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Av)); 12199566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Avwork)); 12209566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->AvT)); 12219566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->DSG)); 12229566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->L)); 12239566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->LT)); 12249566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&user->solver)); 12259566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Js)); 12269566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Jd)); 12279566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->JsInv)); 12289566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->JsBlock)); 12299566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->JsBlockPrec)); 12309566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->u)); 12319566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->uwork)); 12329566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->utrue)); 12339566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->y)); 12349566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->ywork)); 12359566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->ytrue)); 12369566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->nt,&user->yi)); 12379566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->nt,&user->yiwork)); 12389566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(user->ns,&user->di)); 12399566063dSJacob Faibussowitsch PetscCall(PetscFree(user->yi)); 12409566063dSJacob Faibussowitsch PetscCall(PetscFree(user->yiwork)); 12419566063dSJacob Faibussowitsch PetscCall(PetscFree(user->di)); 12429566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->c)); 12439566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->cwork)); 12449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->ur)); 12459566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->q)); 12469566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->d)); 12479566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->dwork)); 12489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->lwork)); 12499566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->S)); 12509566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->Swork)); 12519566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->Av_u)); 12529566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->Twork)); 12539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->Rwork)); 12549566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->js_diag)); 12559566063dSJacob Faibussowitsch PetscCall(ISDestroy(&user->s_is)); 12569566063dSJacob Faibussowitsch PetscCall(ISDestroy(&user->d_is)); 12579566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->state_scatter)); 12589566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->design_scatter)); 1259c4762a1bSJed Brown for (i=0; i<user->nt; i++) { 12609566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->yi_scatter[i])); 1261c4762a1bSJed Brown } 1262c4762a1bSJed Brown for (i=0; i<user->ns; i++) { 12639566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&user->di_scatter[i])); 1264c4762a1bSJed Brown } 12659566063dSJacob Faibussowitsch PetscCall(PetscFree(user->yi_scatter)); 12669566063dSJacob Faibussowitsch PetscCall(PetscFree(user->di_scatter)); 12679566063dSJacob Faibussowitsch PetscCall(PetscFree(user->sample_times)); 1268c4762a1bSJed Brown PetscFunctionReturn(0); 1269c4762a1bSJed Brown } 1270c4762a1bSJed Brown 1271c4762a1bSJed Brown PetscErrorCode ParabolicMonitor(Tao tao, void *ptr) 1272c4762a1bSJed Brown { 1273c4762a1bSJed Brown Vec X; 1274c4762a1bSJed Brown PetscReal unorm,ynorm; 1275c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 1276c4762a1bSJed Brown 1277c4762a1bSJed Brown PetscFunctionBegin; 12789566063dSJacob Faibussowitsch PetscCall(TaoGetSolution(tao,&X)); 12799566063dSJacob Faibussowitsch PetscCall(Scatter(X,user->ywork,user->state_scatter,user->uwork,user->design_scatter)); 12809566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->ywork,-1.0,user->ytrue)); 12819566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->uwork,-1.0,user->utrue)); 12829566063dSJacob Faibussowitsch PetscCall(VecNorm(user->uwork,NORM_2,&unorm)); 12839566063dSJacob Faibussowitsch PetscCall(VecNorm(user->ywork,NORM_2,&ynorm)); 12849566063dSJacob Faibussowitsch PetscCall(PetscPrintf(MPI_COMM_WORLD, "||u-ut||=%g ||y-yt||=%g\n",(double)unorm,(double)ynorm)); 1285c4762a1bSJed Brown PetscFunctionReturn(0); 1286c4762a1bSJed Brown } 1287c4762a1bSJed Brown 1288c4762a1bSJed Brown /*TEST 1289c4762a1bSJed Brown 1290c4762a1bSJed Brown build: 1291c4762a1bSJed Brown requires: !complex 1292c4762a1bSJed Brown 1293c4762a1bSJed Brown test: 1294c4762a1bSJed Brown args: -tao_cmonitor -tao_type lcl -ns 1 -tao_gatol 1.e-4 -ksp_max_it 30 1295c4762a1bSJed Brown requires: !single 1296c4762a1bSJed Brown 1297c4762a1bSJed Brown TEST*/ 1298