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