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