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