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