1 #include <../src/tao/pde_constrained/impls/lcl/lcl.h> 2 #include <../src/tao/matrix/lmvmmat.h> 3 static PetscErrorCode LCLComputeLagrangianAndGradient(TaoLineSearch,Vec,PetscReal*,Vec,void*); 4 static PetscErrorCode LCLComputeAugmentedLagrangianAndGradient(TaoLineSearch,Vec,PetscReal*,Vec,void*); 5 static PetscErrorCode LCLScatter(TAO_LCL*,Vec,Vec,Vec); 6 static PetscErrorCode LCLGather(TAO_LCL*,Vec,Vec,Vec); 7 8 #undef __FUNCT__ 9 #define __FUNCT__ "TaoDestroy_LCL" 10 static PetscErrorCode TaoDestroy_LCL(Tao tao) 11 { 12 TAO_LCL *lclP = (TAO_LCL*)tao->data; 13 PetscErrorCode ierr; 14 15 PetscFunctionBegin; 16 if (tao->setupcalled) { 17 ierr = MatDestroy(&lclP->R);CHKERRQ(ierr); 18 ierr = VecDestroy(&lclP->lamda);CHKERRQ(ierr); 19 ierr = VecDestroy(&lclP->lamda0);CHKERRQ(ierr); 20 ierr = VecDestroy(&lclP->WL);CHKERRQ(ierr); 21 ierr = VecDestroy(&lclP->W);CHKERRQ(ierr); 22 ierr = VecDestroy(&lclP->X0);CHKERRQ(ierr); 23 ierr = VecDestroy(&lclP->G0);CHKERRQ(ierr); 24 ierr = VecDestroy(&lclP->GL);CHKERRQ(ierr); 25 ierr = VecDestroy(&lclP->GAugL);CHKERRQ(ierr); 26 ierr = VecDestroy(&lclP->dbar);CHKERRQ(ierr); 27 ierr = VecDestroy(&lclP->U);CHKERRQ(ierr); 28 ierr = VecDestroy(&lclP->U0);CHKERRQ(ierr); 29 ierr = VecDestroy(&lclP->V);CHKERRQ(ierr); 30 ierr = VecDestroy(&lclP->V0);CHKERRQ(ierr); 31 ierr = VecDestroy(&lclP->V1);CHKERRQ(ierr); 32 ierr = VecDestroy(&lclP->GU);CHKERRQ(ierr); 33 ierr = VecDestroy(&lclP->GV);CHKERRQ(ierr); 34 ierr = VecDestroy(&lclP->GU0);CHKERRQ(ierr); 35 ierr = VecDestroy(&lclP->GV0);CHKERRQ(ierr); 36 ierr = VecDestroy(&lclP->GL_U);CHKERRQ(ierr); 37 ierr = VecDestroy(&lclP->GL_V);CHKERRQ(ierr); 38 ierr = VecDestroy(&lclP->GAugL_U);CHKERRQ(ierr); 39 ierr = VecDestroy(&lclP->GAugL_V);CHKERRQ(ierr); 40 ierr = VecDestroy(&lclP->GL_U0);CHKERRQ(ierr); 41 ierr = VecDestroy(&lclP->GL_V0);CHKERRQ(ierr); 42 ierr = VecDestroy(&lclP->GAugL_U0);CHKERRQ(ierr); 43 ierr = VecDestroy(&lclP->GAugL_V0);CHKERRQ(ierr); 44 ierr = VecDestroy(&lclP->DU);CHKERRQ(ierr); 45 ierr = VecDestroy(&lclP->DV);CHKERRQ(ierr); 46 ierr = VecDestroy(&lclP->WU);CHKERRQ(ierr); 47 ierr = VecDestroy(&lclP->WV);CHKERRQ(ierr); 48 ierr = VecDestroy(&lclP->g1);CHKERRQ(ierr); 49 ierr = VecDestroy(&lclP->g2);CHKERRQ(ierr); 50 ierr = VecDestroy(&lclP->con1);CHKERRQ(ierr); 51 52 ierr = VecDestroy(&lclP->r);CHKERRQ(ierr); 53 ierr = VecDestroy(&lclP->s);CHKERRQ(ierr); 54 55 ierr = ISDestroy(&tao->state_is);CHKERRQ(ierr); 56 ierr = ISDestroy(&tao->design_is);CHKERRQ(ierr); 57 58 ierr = VecScatterDestroy(&lclP->state_scatter);CHKERRQ(ierr); 59 ierr = VecScatterDestroy(&lclP->design_scatter);CHKERRQ(ierr); 60 } 61 ierr = PetscFree(tao->data); 62 PetscFunctionReturn(0); 63 } 64 65 #undef __FUNCT__ 66 #define __FUNCT__ "TaoSetFromOptions_LCL" 67 static PetscErrorCode TaoSetFromOptions_LCL(Tao tao) 68 { 69 TAO_LCL *lclP = (TAO_LCL*)tao->data; 70 PetscBool flg; 71 PetscErrorCode ierr; 72 73 PetscFunctionBegin; 74 ierr = PetscOptionsHead("Linearly-Constrained Augmented Lagrangian Method for PDE-constrained optimization");CHKERRQ(ierr); 75 ierr = PetscOptionsReal("-tao_lcl_eps1","epsilon 1 tolerance","",lclP->eps1,&lclP->eps1,&flg);CHKERRQ(ierr); 76 ierr = PetscOptionsReal("-tao_lcl_eps2","epsilon 2 tolerance","",lclP->eps2,&lclP->eps2,&flg);CHKERRQ(ierr); 77 ierr = PetscOptionsReal("-tao_lcl_rho0","init value for rho","",lclP->rho0,&lclP->rho0,&flg);CHKERRQ(ierr); 78 ierr = PetscOptionsReal("-tao_lcl_rhomax","max value for rho","",lclP->rhomax,&lclP->rhomax,&flg);CHKERRQ(ierr); 79 lclP->phase2_niter = 1; 80 ierr = PetscOptionsInt("-tao_lcl_phase2_niter","Number of phase 2 iterations in LCL algorithm","",lclP->phase2_niter,&lclP->phase2_niter,&flg);CHKERRQ(ierr); 81 lclP->verbose = PETSC_FALSE; 82 ierr = PetscOptionsBool("-tao_lcl_verbose","Print verbose output","",lclP->verbose,&lclP->verbose,&flg);CHKERRQ(ierr); 83 lclP->tau[0] = lclP->tau[1] = lclP->tau[2] = lclP->tau[3] = 1.0e-4; 84 ierr = PetscOptionsReal("-tao_lcl_tola","Tolerance for first forward solve","",lclP->tau[0],&lclP->tau[0],&flg);CHKERRQ(ierr); 85 ierr = PetscOptionsReal("-tao_lcl_tolb","Tolerance for first adjoint solve","",lclP->tau[1],&lclP->tau[1],&flg);CHKERRQ(ierr); 86 ierr = PetscOptionsReal("-tao_lcl_tolc","Tolerance for second forward solve","",lclP->tau[2],&lclP->tau[2],&flg);CHKERRQ(ierr); 87 ierr = PetscOptionsReal("-tao_lcl_told","Tolerance for second adjoint solve","",lclP->tau[3],&lclP->tau[3],&flg);CHKERRQ(ierr); 88 ierr = PetscOptionsTail();CHKERRQ(ierr); 89 ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 90 PetscFunctionReturn(0); 91 } 92 93 #undef __FUNCT__ 94 #define __FUNCT__ "TaoView_LCL" 95 static PetscErrorCode TaoView_LCL(Tao tao, PetscViewer viewer) 96 { 97 return 0; 98 } 99 100 #undef __FUNCT__ 101 #define __FUNCT__ "TaoSetup_LCL" 102 static PetscErrorCode TaoSetup_LCL(Tao tao) 103 { 104 TAO_LCL *lclP = (TAO_LCL*)tao->data; 105 PetscInt lo, hi, nlocalstate, nlocaldesign; 106 PetscErrorCode ierr; 107 IS is_state, is_design; 108 109 PetscFunctionBegin; 110 if (!tao->state_is) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONGSTATE,"LCL Solver requires an initial state index set -- use TaoSetStateIS()"); 111 ierr = VecDuplicate(tao->solution, &tao->gradient);CHKERRQ(ierr); 112 ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr); 113 ierr = VecDuplicate(tao->solution, &lclP->W);CHKERRQ(ierr); 114 ierr = VecDuplicate(tao->solution, &lclP->X0);CHKERRQ(ierr); 115 ierr = VecDuplicate(tao->solution, &lclP->G0);CHKERRQ(ierr); 116 ierr = VecDuplicate(tao->solution, &lclP->GL);CHKERRQ(ierr); 117 ierr = VecDuplicate(tao->solution, &lclP->GAugL);CHKERRQ(ierr); 118 119 ierr = VecDuplicate(tao->constraints, &lclP->lamda);CHKERRQ(ierr); 120 ierr = VecDuplicate(tao->constraints, &lclP->WL);CHKERRQ(ierr); 121 ierr = VecDuplicate(tao->constraints, &lclP->lamda0);CHKERRQ(ierr); 122 ierr = VecDuplicate(tao->constraints, &lclP->con1);CHKERRQ(ierr); 123 124 ierr = VecSet(lclP->lamda,0.0);CHKERRQ(ierr); 125 126 ierr = VecGetSize(tao->solution, &lclP->n);CHKERRQ(ierr); 127 ierr = VecGetSize(tao->constraints, &lclP->m);CHKERRQ(ierr); 128 129 ierr = VecCreate(((PetscObject)tao)->comm,&lclP->U);CHKERRQ(ierr); 130 ierr = VecCreate(((PetscObject)tao)->comm,&lclP->V);CHKERRQ(ierr); 131 ierr = ISGetLocalSize(tao->state_is,&nlocalstate);CHKERRQ(ierr); 132 ierr = ISGetLocalSize(tao->design_is,&nlocaldesign);CHKERRQ(ierr); 133 ierr = VecSetSizes(lclP->U,nlocalstate,lclP->m);CHKERRQ(ierr); 134 ierr = VecSetSizes(lclP->V,nlocaldesign,lclP->n-lclP->m);CHKERRQ(ierr); 135 ierr = VecSetType(lclP->U,((PetscObject)(tao->solution))->type_name);CHKERRQ(ierr); 136 ierr = VecSetType(lclP->V,((PetscObject)(tao->solution))->type_name);CHKERRQ(ierr); 137 ierr = VecSetFromOptions(lclP->U);CHKERRQ(ierr); 138 ierr = VecSetFromOptions(lclP->V);CHKERRQ(ierr); 139 ierr = VecDuplicate(lclP->U,&lclP->DU);CHKERRQ(ierr); 140 ierr = VecDuplicate(lclP->U,&lclP->U0);CHKERRQ(ierr); 141 ierr = VecDuplicate(lclP->U,&lclP->GU);CHKERRQ(ierr); 142 ierr = VecDuplicate(lclP->U,&lclP->GU0);CHKERRQ(ierr); 143 ierr = VecDuplicate(lclP->U,&lclP->GAugL_U);CHKERRQ(ierr); 144 ierr = VecDuplicate(lclP->U,&lclP->GL_U);CHKERRQ(ierr); 145 ierr = VecDuplicate(lclP->U,&lclP->GAugL_U0);CHKERRQ(ierr); 146 ierr = VecDuplicate(lclP->U,&lclP->GL_U0);CHKERRQ(ierr); 147 ierr = VecDuplicate(lclP->U,&lclP->WU);CHKERRQ(ierr); 148 ierr = VecDuplicate(lclP->U,&lclP->r);CHKERRQ(ierr); 149 ierr = VecDuplicate(lclP->V,&lclP->V0);CHKERRQ(ierr); 150 ierr = VecDuplicate(lclP->V,&lclP->V1);CHKERRQ(ierr); 151 ierr = VecDuplicate(lclP->V,&lclP->DV);CHKERRQ(ierr); 152 ierr = VecDuplicate(lclP->V,&lclP->s);CHKERRQ(ierr); 153 ierr = VecDuplicate(lclP->V,&lclP->GV);CHKERRQ(ierr); 154 ierr = VecDuplicate(lclP->V,&lclP->GV0);CHKERRQ(ierr); 155 ierr = VecDuplicate(lclP->V,&lclP->dbar);CHKERRQ(ierr); 156 ierr = VecDuplicate(lclP->V,&lclP->GAugL_V);CHKERRQ(ierr); 157 ierr = VecDuplicate(lclP->V,&lclP->GL_V);CHKERRQ(ierr); 158 ierr = VecDuplicate(lclP->V,&lclP->GAugL_V0);CHKERRQ(ierr); 159 ierr = VecDuplicate(lclP->V,&lclP->GL_V0);CHKERRQ(ierr); 160 ierr = VecDuplicate(lclP->V,&lclP->WV);CHKERRQ(ierr); 161 ierr = VecDuplicate(lclP->V,&lclP->g1);CHKERRQ(ierr); 162 ierr = VecDuplicate(lclP->V,&lclP->g2);CHKERRQ(ierr); 163 164 /* create scatters for state, design subvecs */ 165 ierr = VecGetOwnershipRange(lclP->U,&lo,&hi);CHKERRQ(ierr); 166 ierr = ISCreateStride(((PetscObject)lclP->U)->comm,hi-lo,lo,1,&is_state);CHKERRQ(ierr); 167 ierr = VecGetOwnershipRange(lclP->V,&lo,&hi);CHKERRQ(ierr); 168 if (0) { 169 PetscInt sizeU,sizeV; 170 ierr = VecGetSize(lclP->U,&sizeU); 171 ierr = VecGetSize(lclP->V,&sizeV); 172 ierr = PetscPrintf(PETSC_COMM_WORLD,"size(U)=%D, size(V)=%D\n",sizeU,sizeV); 173 } 174 ierr = ISCreateStride(((PetscObject)lclP->V)->comm,hi-lo,lo,1,&is_design);CHKERRQ(ierr); 175 ierr = VecScatterCreate(tao->solution,tao->state_is,lclP->U,is_state,&lclP->state_scatter);CHKERRQ(ierr); 176 ierr = VecScatterCreate(tao->solution,tao->design_is,lclP->V,is_design,&lclP->design_scatter);CHKERRQ(ierr); 177 ierr = ISDestroy(&is_state);CHKERRQ(ierr); 178 ierr = ISDestroy(&is_design);CHKERRQ(ierr); 179 PetscFunctionReturn(0); 180 } 181 182 #undef __FUNCT__ 183 #define __FUNCT__ "TaoSolve_LCL" 184 static PetscErrorCode TaoSolve_LCL(Tao tao) 185 { 186 TAO_LCL *lclP = (TAO_LCL*)tao->data; 187 PetscInt iter=0,phase2_iter,nlocal,its; 188 TaoConvergedReason reason = TAO_CONTINUE_ITERATING; 189 TaoLineSearchConvergedReason ls_reason = TAOLINESEARCH_CONTINUE_ITERATING; 190 PetscReal step=1.0,f, descent, aldescent; 191 PetscReal cnorm, mnorm; 192 PetscReal adec,r2,rGL_U,rWU; 193 PetscBool set,pset,flag,pflag,symmetric; 194 PetscErrorCode ierr; 195 196 PetscFunctionBegin; 197 lclP->rho = lclP->rho0; 198 ierr = VecGetLocalSize(lclP->U,&nlocal);CHKERRQ(ierr); 199 ierr = VecGetLocalSize(lclP->V,&nlocal);CHKERRQ(ierr); 200 ierr = MatCreateLMVM(((PetscObject)tao)->comm,nlocal,lclP->n-lclP->m,&lclP->R);CHKERRQ(ierr); 201 ierr = MatLMVMAllocateVectors(lclP->R,lclP->V);CHKERRQ(ierr); 202 lclP->recompute_jacobian_flag = PETSC_TRUE; 203 204 /* Scatter to U,V */ 205 ierr = LCLScatter(lclP,tao->solution,lclP->U,lclP->V);CHKERRQ(ierr); 206 207 /* Evaluate Function, Gradient, Constraints, and Jacobian */ 208 ierr = TaoComputeObjectiveAndGradient(tao,tao->solution,&f,tao->gradient);CHKERRQ(ierr); 209 ierr = TaoComputeJacobianState(tao,tao->solution,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv);CHKERRQ(ierr); 210 ierr = TaoComputeJacobianDesign(tao,tao->solution,tao->jacobian_design);CHKERRQ(ierr); 211 ierr = TaoComputeConstraints(tao,tao->solution, tao->constraints);CHKERRQ(ierr); 212 213 /* Scatter gradient to GU,GV */ 214 ierr = LCLScatter(lclP,tao->gradient,lclP->GU,lclP->GV);CHKERRQ(ierr); 215 216 /* Evaluate Lagrangian function and gradient */ 217 /* p0 */ 218 ierr = VecSet(lclP->lamda,0.0);CHKERRQ(ierr); /* Initial guess in CG */ 219 ierr = MatIsSymmetricKnown(tao->jacobian_state,&set,&flag);CHKERRQ(ierr); 220 if (tao->jacobian_state_pre) { 221 ierr = MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag);CHKERRQ(ierr); 222 } else { 223 pset = pflag = PETSC_TRUE; 224 } 225 if (set && pset && flag && pflag) symmetric = PETSC_TRUE; 226 else symmetric = PETSC_FALSE; 227 228 lclP->solve_type = LCL_ADJOINT2; 229 if (tao->jacobian_state_inv) { 230 if (symmetric) { 231 ierr = MatMult(tao->jacobian_state_inv, lclP->GU, lclP->lamda);CHKERRQ(ierr); } else { 232 ierr = MatMultTranspose(tao->jacobian_state_inv, lclP->GU, lclP->lamda);CHKERRQ(ierr); 233 } 234 } else { 235 ierr = KSPSetOperators(tao->ksp, tao->jacobian_state, tao->jacobian_state_pre);CHKERRQ(ierr); 236 if (symmetric) { 237 ierr = KSPSolve(tao->ksp, lclP->GU, lclP->lamda);CHKERRQ(ierr); 238 } else { 239 ierr = KSPSolveTranspose(tao->ksp, lclP->GU, lclP->lamda);CHKERRQ(ierr); 240 } 241 ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 242 tao->ksp_its += its; 243 } 244 ierr = VecCopy(lclP->lamda,lclP->lamda0);CHKERRQ(ierr); 245 ierr = LCLComputeAugmentedLagrangianAndGradient(tao->linesearch,tao->solution,&lclP->aug,lclP->GAugL,tao);CHKERRQ(ierr); 246 247 ierr = LCLScatter(lclP,lclP->GL,lclP->GL_U,lclP->GL_V);CHKERRQ(ierr); 248 ierr = LCLScatter(lclP,lclP->GAugL,lclP->GAugL_U,lclP->GAugL_V);CHKERRQ(ierr); 249 250 /* Evaluate constraint norm */ 251 ierr = VecNorm(tao->constraints, NORM_2, &cnorm);CHKERRQ(ierr); 252 ierr = VecNorm(lclP->GAugL, NORM_2, &mnorm);CHKERRQ(ierr); 253 254 /* Monitor convergence */ 255 ierr = TaoMonitor(tao, iter,f,mnorm,cnorm,step,&reason);CHKERRQ(ierr); 256 257 while (reason == TAO_CONTINUE_ITERATING) { 258 /* Compute a descent direction for the linearly constrained subproblem 259 minimize f(u+du, v+dv) 260 s.t. A(u0,v0)du + B(u0,v0)dv = -g(u0,v0) */ 261 262 /* Store the points around the linearization */ 263 ierr = VecCopy(lclP->U, lclP->U0);CHKERRQ(ierr); 264 ierr = VecCopy(lclP->V, lclP->V0);CHKERRQ(ierr); 265 ierr = VecCopy(lclP->GU,lclP->GU0);CHKERRQ(ierr); 266 ierr = VecCopy(lclP->GV,lclP->GV0);CHKERRQ(ierr); 267 ierr = VecCopy(lclP->GAugL_U,lclP->GAugL_U0);CHKERRQ(ierr); 268 ierr = VecCopy(lclP->GAugL_V,lclP->GAugL_V0);CHKERRQ(ierr); 269 ierr = VecCopy(lclP->GL_U,lclP->GL_U0);CHKERRQ(ierr); 270 ierr = VecCopy(lclP->GL_V,lclP->GL_V0);CHKERRQ(ierr); 271 272 lclP->aug0 = lclP->aug; 273 lclP->lgn0 = lclP->lgn; 274 275 /* Given the design variables, we need to project the current iterate 276 onto the linearized constraint. We choose to fix the design variables 277 and solve the linear system for the state variables. The resulting 278 point is the Newton direction */ 279 280 /* Solve r = A\con */ 281 lclP->solve_type = LCL_FORWARD1; 282 ierr = VecSet(lclP->r,0.0);CHKERRQ(ierr); /* Initial guess in CG */ 283 284 if (tao->jacobian_state_inv) { 285 ierr = MatMult(tao->jacobian_state_inv, tao->constraints, lclP->r);CHKERRQ(ierr); 286 } else { 287 ierr = KSPSetOperators(tao->ksp, tao->jacobian_state, tao->jacobian_state_pre);CHKERRQ(ierr); 288 ierr = KSPSolve(tao->ksp, tao->constraints, lclP->r);CHKERRQ(ierr); 289 ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 290 tao->ksp_its += its; 291 } 292 293 /* Set design step direction dv to zero */ 294 ierr = VecSet(lclP->s, 0.0);CHKERRQ(ierr); 295 296 /* 297 Check sufficient descent for constraint merit function .5*||con||^2 298 con' Ak r >= eps1 ||r||^(2+eps2) 299 */ 300 301 /* Compute WU= Ak' * con */ 302 if (symmetric) { 303 ierr = MatMult(tao->jacobian_state,tao->constraints,lclP->WU);CHKERRQ(ierr); 304 } else { 305 ierr = MatMultTranspose(tao->jacobian_state,tao->constraints,lclP->WU);CHKERRQ(ierr); 306 } 307 /* Compute r * Ak' * con */ 308 ierr = VecDot(lclP->r,lclP->WU,&rWU);CHKERRQ(ierr); 309 310 /* compute ||r||^(2+eps2) */ 311 ierr = VecNorm(lclP->r,NORM_2,&r2);CHKERRQ(ierr); 312 r2 = PetscPowScalar(r2,2.0+lclP->eps2); 313 adec = lclP->eps1 * r2; 314 315 if (rWU < adec) { 316 ierr = PetscInfo(tao,"Newton direction not descent for constraint, feasibility phase required");CHKERRQ(ierr); 317 if (lclP->verbose) { 318 ierr = PetscPrintf(PETSC_COMM_WORLD,"Newton direction not descent for constraint: %g -- using steepest descent\n",(double)descent);CHKERRQ(ierr); 319 } 320 321 ierr = PetscInfo(tao,"Using steepest descent direction instead.\n");CHKERRQ(ierr); 322 ierr = VecSet(lclP->r,0.0);CHKERRQ(ierr); 323 ierr = VecAXPY(lclP->r,-1.0,lclP->WU);CHKERRQ(ierr); 324 ierr = VecDot(lclP->r,lclP->r,&rWU);CHKERRQ(ierr); 325 ierr = VecNorm(lclP->r,NORM_2,&r2);CHKERRQ(ierr); 326 r2 = PetscPowScalar(r2,2.0+lclP->eps2); 327 ierr = VecDot(lclP->r,lclP->GAugL_U,&descent);CHKERRQ(ierr); 328 adec = lclP->eps1 * r2; 329 } 330 331 332 /* 333 Check descent for aug. lagrangian 334 r' (GUk - Ak'*yk - rho*Ak'*con) <= -eps1 ||r||^(2+eps2) 335 GL_U = GUk - Ak'*yk 336 WU = Ak'*con 337 adec=eps1||r||^(2+eps2) 338 339 ==> 340 Check r'GL_U - rho*r'WU <= adec 341 */ 342 343 ierr = VecDot(lclP->r,lclP->GL_U,&rGL_U); 344 aldescent = rGL_U - lclP->rho*rWU; 345 if (aldescent > -adec) { 346 if (lclP->verbose) { 347 ierr = PetscPrintf(PETSC_COMM_WORLD," Newton direction not descent for augmented Lagrangian: %g",(double)aldescent);CHKERRQ(ierr); 348 } 349 ierr = PetscInfo1(tao,"Newton direction not descent for augmented Lagrangian: %g",(double)aldescent);CHKERRQ(ierr); 350 lclP->rho = (rGL_U - adec)/rWU; 351 if (lclP->rho > lclP->rhomax) { 352 lclP->rho = lclP->rhomax; 353 SETERRQ1(PETSC_COMM_WORLD,0,"rho=%g > rhomax, case not implemented. Increase rhomax (-tao_lcl_rhomax)",(double)lclP->rho); 354 } 355 if (lclP->verbose) { 356 ierr = PetscPrintf(PETSC_COMM_WORLD," Increasing penalty parameter to %g\n",(double)lclP->rho);CHKERRQ(ierr); 357 } 358 ierr = PetscInfo1(tao," Increasing penalty parameter to %g",(double)lclP->rho); 359 } 360 361 ierr = LCLComputeAugmentedLagrangianAndGradient(tao->linesearch,tao->solution,&lclP->aug,lclP->GAugL,tao);CHKERRQ(ierr); 362 ierr = LCLScatter(lclP,lclP->GL,lclP->GL_U,lclP->GL_V);CHKERRQ(ierr); 363 ierr = LCLScatter(lclP,lclP->GAugL,lclP->GAugL_U,lclP->GAugL_V);CHKERRQ(ierr); 364 365 /* We now minimize the augmented Lagrangian along the Newton direction */ 366 ierr = VecScale(lclP->r,-1.0);CHKERRQ(ierr); 367 ierr = LCLGather(lclP, lclP->r,lclP->s,tao->stepdirection);CHKERRQ(ierr); 368 ierr = VecScale(lclP->r,-1.0);CHKERRQ(ierr); 369 ierr = LCLGather(lclP, lclP->GAugL_U0, lclP->GAugL_V0, lclP->GAugL);CHKERRQ(ierr); 370 ierr = LCLGather(lclP, lclP->U0,lclP->V0,lclP->X0);CHKERRQ(ierr); 371 372 lclP->recompute_jacobian_flag = PETSC_TRUE; 373 374 ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr); 375 ierr = TaoLineSearchSetType(tao->linesearch, TAOLINESEARCHMT);CHKERRQ(ierr); 376 ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 377 ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &lclP->aug, lclP->GAugL, tao->stepdirection, &step, &ls_reason);CHKERRQ(ierr); 378 if (lclP->verbose) { 379 ierr = PetscPrintf(PETSC_COMM_WORLD,"Steplength = %g\n",(double)step);CHKERRQ(ierr); 380 } 381 382 ierr = LCLScatter(lclP,tao->solution,lclP->U,lclP->V);CHKERRQ(ierr); 383 ierr = TaoComputeObjectiveAndGradient(tao,tao->solution,&f,tao->gradient);CHKERRQ(ierr); 384 ierr = LCLScatter(lclP,tao->gradient,lclP->GU,lclP->GV);CHKERRQ(ierr); 385 386 ierr = LCLScatter(lclP,lclP->GAugL,lclP->GAugL_U,lclP->GAugL_V);CHKERRQ(ierr); 387 388 /* Check convergence */ 389 ierr = VecNorm(lclP->GAugL, NORM_2, &mnorm);CHKERRQ(ierr); 390 ierr = VecNorm(tao->constraints, NORM_2, &cnorm);CHKERRQ(ierr); 391 ierr = TaoMonitor(tao, iter,f,mnorm,cnorm,step,&reason);CHKERRQ(ierr); 392 if (reason != TAO_CONTINUE_ITERATING){ 393 break; 394 } 395 396 /* TODO: use a heuristic to choose how many iterations should be performed within phase 2 */ 397 for (phase2_iter=0; phase2_iter<lclP->phase2_niter; phase2_iter++){ 398 /* We now minimize the objective function starting from the fraction of 399 the Newton point accepted by applying one step of a reduced-space 400 method. The optimization problem is: 401 402 minimize f(u+du, v+dv) 403 s. t. A(u0,v0)du + B(u0,v0)du = -alpha g(u0,v0) 404 405 In particular, we have that 406 du = -inv(A)*(Bdv + alpha g) */ 407 408 ierr = TaoComputeJacobianState(tao,lclP->X0,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv);CHKERRQ(ierr); 409 ierr = TaoComputeJacobianDesign(tao,lclP->X0,tao->jacobian_design);CHKERRQ(ierr); 410 411 /* Store V and constraints */ 412 ierr = VecCopy(lclP->V, lclP->V1);CHKERRQ(ierr); 413 ierr = VecCopy(tao->constraints,lclP->con1);CHKERRQ(ierr); 414 415 /* Compute multipliers */ 416 /* p1 */ 417 ierr = VecSet(lclP->lamda,0.0);CHKERRQ(ierr); /* Initial guess in CG */ 418 lclP->solve_type = LCL_ADJOINT1; 419 ierr = MatIsSymmetricKnown(tao->jacobian_state,&set,&flag);CHKERRQ(ierr); 420 if (tao->jacobian_state_pre) { 421 ierr = MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag);CHKERRQ(ierr); 422 } else { 423 pset = pflag = PETSC_TRUE; 424 } 425 if (set && pset && flag && pflag) symmetric = PETSC_TRUE; 426 else symmetric = PETSC_FALSE; 427 428 if (tao->jacobian_state_inv) { 429 if (symmetric) { 430 ierr = MatMult(tao->jacobian_state_inv, lclP->GAugL_U, lclP->lamda);CHKERRQ(ierr); 431 } else { 432 ierr = MatMultTranspose(tao->jacobian_state_inv, lclP->GAugL_U, lclP->lamda);CHKERRQ(ierr); 433 } 434 } else { 435 if (symmetric) { 436 ierr = KSPSolve(tao->ksp, lclP->GAugL_U, lclP->lamda);CHKERRQ(ierr); 437 } else { 438 ierr = KSPSolveTranspose(tao->ksp, lclP->GAugL_U, lclP->lamda);CHKERRQ(ierr); 439 } 440 ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 441 tao->ksp_its += its; 442 } 443 ierr = MatMultTranspose(tao->jacobian_design,lclP->lamda,lclP->g1);CHKERRQ(ierr); 444 ierr = VecAXPY(lclP->g1,-1.0,lclP->GAugL_V);CHKERRQ(ierr); 445 446 /* Compute the limited-memory quasi-newton direction */ 447 if (iter > 0) { 448 ierr = MatLMVMSolve(lclP->R,lclP->g1,lclP->s);CHKERRQ(ierr); 449 ierr = VecDot(lclP->s,lclP->g1,&descent);CHKERRQ(ierr); 450 if (descent <= 0) { 451 if (lclP->verbose) { 452 ierr = PetscPrintf(PETSC_COMM_WORLD,"Reduced-space direction not descent: %g\n",(double)descent);CHKERRQ(ierr); 453 } 454 ierr = VecCopy(lclP->g1,lclP->s);CHKERRQ(ierr); 455 } 456 } else { 457 ierr = VecCopy(lclP->g1,lclP->s);CHKERRQ(ierr); 458 } 459 ierr = VecScale(lclP->g1,-1.0);CHKERRQ(ierr); 460 461 /* Recover the full space direction */ 462 ierr = MatMult(tao->jacobian_design,lclP->s,lclP->WU);CHKERRQ(ierr); 463 /* ierr = VecSet(lclP->r,0.0);CHKERRQ(ierr); */ /* Initial guess in CG */ 464 lclP->solve_type = LCL_FORWARD2; 465 if (tao->jacobian_state_inv) { 466 ierr = MatMult(tao->jacobian_state_inv,lclP->WU,lclP->r);CHKERRQ(ierr); 467 } else { 468 ierr = KSPSolve(tao->ksp, lclP->WU, lclP->r);CHKERRQ(ierr); 469 ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 470 tao->ksp_its += its; 471 } 472 473 /* We now minimize the augmented Lagrangian along the direction -r,s */ 474 ierr = VecScale(lclP->r, -1.0);CHKERRQ(ierr); 475 ierr = LCLGather(lclP,lclP->r,lclP->s,tao->stepdirection);CHKERRQ(ierr); 476 ierr = VecScale(lclP->r, -1.0);CHKERRQ(ierr); 477 lclP->recompute_jacobian_flag = PETSC_TRUE; 478 479 ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr); 480 ierr = TaoLineSearchSetType(tao->linesearch,TAOLINESEARCHMT);CHKERRQ(ierr); 481 ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 482 ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &lclP->aug, lclP->GAugL, tao->stepdirection,&step,&ls_reason);CHKERRQ(ierr); 483 if (lclP->verbose){ 484 ierr = PetscPrintf(PETSC_COMM_WORLD,"Reduced-space steplength = %g\n",(double)step);CHKERRQ(ierr); 485 } 486 487 ierr = LCLScatter(lclP,tao->solution,lclP->U,lclP->V);CHKERRQ(ierr); 488 ierr = LCLScatter(lclP,lclP->GL,lclP->GL_U,lclP->GL_V);CHKERRQ(ierr); 489 ierr = LCLScatter(lclP,lclP->GAugL,lclP->GAugL_U,lclP->GAugL_V);CHKERRQ(ierr); 490 ierr = TaoComputeObjectiveAndGradient(tao,tao->solution,&f,tao->gradient);CHKERRQ(ierr); 491 ierr = LCLScatter(lclP,tao->gradient,lclP->GU,lclP->GV);CHKERRQ(ierr); 492 493 /* Compute the reduced gradient at the new point */ 494 495 ierr = TaoComputeJacobianState(tao,lclP->X0,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv);CHKERRQ(ierr); 496 ierr = TaoComputeJacobianDesign(tao,lclP->X0,tao->jacobian_design);CHKERRQ(ierr); 497 498 /* p2 */ 499 /* Compute multipliers, using lamda-rho*con as an initial guess in PCG */ 500 if (phase2_iter==0){ 501 ierr = VecWAXPY(lclP->lamda,-lclP->rho,lclP->con1,lclP->lamda0);CHKERRQ(ierr); 502 } else { 503 ierr = VecAXPY(lclP->lamda,-lclP->rho,tao->constraints);CHKERRQ(ierr); 504 } 505 506 ierr = MatIsSymmetricKnown(tao->jacobian_state,&set,&flag);CHKERRQ(ierr); 507 if (tao->jacobian_state_pre) { 508 ierr = MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag);CHKERRQ(ierr); 509 } else { 510 pset = pflag = PETSC_TRUE; 511 } 512 if (set && pset && flag && pflag) symmetric = PETSC_TRUE; 513 else symmetric = PETSC_FALSE; 514 515 lclP->solve_type = LCL_ADJOINT2; 516 if (tao->jacobian_state_inv) { 517 if (symmetric) { 518 ierr = MatMult(tao->jacobian_state_inv, lclP->GU, lclP->lamda);CHKERRQ(ierr); 519 } else { 520 ierr = MatMultTranspose(tao->jacobian_state_inv, lclP->GU, lclP->lamda);CHKERRQ(ierr); 521 } 522 } else { 523 if (symmetric) { 524 ierr = KSPSolve(tao->ksp, lclP->GU, lclP->lamda);CHKERRQ(ierr); 525 } else { 526 ierr = KSPSolveTranspose(tao->ksp, lclP->GU, lclP->lamda);CHKERRQ(ierr); 527 } 528 ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 529 tao->ksp_its += its; 530 } 531 532 ierr = MatMultTranspose(tao->jacobian_design,lclP->lamda,lclP->g2);CHKERRQ(ierr); 533 ierr = VecAXPY(lclP->g2,-1.0,lclP->GV);CHKERRQ(ierr); 534 535 ierr = VecScale(lclP->g2,-1.0);CHKERRQ(ierr); 536 537 /* Update the quasi-newton approximation */ 538 if (phase2_iter >= 0){ 539 ierr = MatLMVMSetPrev(lclP->R,lclP->V1,lclP->g1); 540 } 541 ierr = MatLMVMUpdate(lclP->R,lclP->V,lclP->g2);CHKERRQ(ierr); 542 /* Use "-tao_ls_type gpcg -tao_ls_ftol 0 -tao_lmm_broyden_phi 0.0 -tao_lmm_scale_type scalar" to obtain agreement with Matlab code */ 543 544 } 545 546 ierr = VecCopy(lclP->lamda,lclP->lamda0);CHKERRQ(ierr); 547 548 /* Evaluate Function, Gradient, Constraints, and Jacobian */ 549 ierr = TaoComputeObjectiveAndGradient(tao,tao->solution,&f,tao->gradient);CHKERRQ(ierr); 550 ierr = LCLScatter(lclP,tao->solution,lclP->U,lclP->V);CHKERRQ(ierr); 551 ierr = LCLScatter(lclP,tao->gradient,lclP->GU,lclP->GV);CHKERRQ(ierr); 552 553 ierr = TaoComputeJacobianState(tao,tao->solution,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv);CHKERRQ(ierr); 554 ierr = TaoComputeJacobianDesign(tao,tao->solution,tao->jacobian_design);CHKERRQ(ierr); 555 ierr = TaoComputeConstraints(tao,tao->solution, tao->constraints);CHKERRQ(ierr); 556 557 ierr = LCLComputeAugmentedLagrangianAndGradient(tao->linesearch,tao->solution,&lclP->aug,lclP->GAugL,tao);CHKERRQ(ierr); 558 559 ierr = VecNorm(lclP->GAugL, NORM_2, &mnorm);CHKERRQ(ierr); 560 561 /* Evaluate constraint norm */ 562 ierr = VecNorm(tao->constraints, NORM_2, &cnorm);CHKERRQ(ierr); 563 564 /* Monitor convergence */ 565 iter++; 566 ierr = TaoMonitor(tao, iter,f,mnorm,cnorm,step,&reason);CHKERRQ(ierr); 567 } 568 ierr = MatDestroy(&lclP->R);CHKERRQ(ierr); 569 PetscFunctionReturn(0); 570 } 571 572 /*MC 573 TAOLCL - linearly constrained lagrangian method for pde-constrained optimization 574 575 + -tao_lcl_eps1 - epsilon 1 tolerance 576 . -tao_lcl_eps2","epsilon 2 tolerance","",lclP->eps2,&lclP->eps2,&flg);CHKERRQ(ierr); 577 . -tao_lcl_rho0","init value for rho","",lclP->rho0,&lclP->rho0,&flg);CHKERRQ(ierr); 578 . -tao_lcl_rhomax","max value for rho","",lclP->rhomax,&lclP->rhomax,&flg);CHKERRQ(ierr); 579 . -tao_lcl_phase2_niter - Number of phase 2 iterations in LCL algorithm 580 . -tao_lcl_verbose - Print verbose output if True 581 . -tao_lcl_tola - Tolerance for first forward solve 582 . -tao_lcl_tolb - Tolerance for first adjoint solve 583 . -tao_lcl_tolc - Tolerance for second forward solve 584 - -tao_lcl_told - Tolerance for second adjoint solve 585 586 M*/ 587 EXTERN_C_BEGIN 588 #undef __FUNCT__ 589 #define __FUNCT__ "TaoCreate_LCL" 590 PetscErrorCode TaoCreate_LCL(Tao tao) 591 { 592 TAO_LCL *lclP; 593 PetscErrorCode ierr; 594 const char *morethuente_type = TAOLINESEARCHMT; 595 596 PetscFunctionBegin; 597 tao->ops->setup = TaoSetup_LCL; 598 tao->ops->solve = TaoSolve_LCL; 599 tao->ops->view = TaoView_LCL; 600 tao->ops->setfromoptions = TaoSetFromOptions_LCL; 601 tao->ops->destroy = TaoDestroy_LCL; 602 ierr = PetscNewLog(tao,&lclP);CHKERRQ(ierr); 603 tao->data = (void*)lclP; 604 605 tao->max_it=200; 606 #if defined(PETSC_USE_REAL_SINGLE) 607 tao->fatol=1e-5; 608 tao->frtol=1e-5; 609 #else 610 tao->fatol=1e-8; 611 tao->frtol=1e-8; 612 #endif 613 tao->catol=1e-4; 614 tao->crtol=1e-4; 615 tao->gttol=1e-4; 616 tao->gatol=1e-4; 617 tao->grtol=1e-4; 618 lclP->rho0 = 1.0e-4; 619 lclP->rhomax=1e5; 620 lclP->eps1 = 1.0e-8; 621 lclP->eps2 = 0.0; 622 lclP->solve_type=2; 623 lclP->tau[0] = lclP->tau[1] = lclP->tau[2] = lclP->tau[3] = 1.0e-4; 624 ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr); 625 ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr); 626 627 ierr = TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,LCLComputeAugmentedLagrangianAndGradient, tao);CHKERRQ(ierr); 628 ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr); 629 ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 630 PetscFunctionReturn(0); 631 } 632 EXTERN_C_END 633 634 #undef __FUNCT__ 635 #define __FUNCT__ "LCLComputeLagrangianAndGradient" 636 static PetscErrorCode LCLComputeLagrangianAndGradient(TaoLineSearch ls, Vec X, PetscReal *f, Vec G, void *ptr) 637 { 638 Tao tao = (Tao)ptr; 639 TAO_LCL *lclP = (TAO_LCL*)tao->data; 640 PetscBool set,pset,flag,pflag,symmetric; 641 PetscReal cdotl; 642 PetscErrorCode ierr; 643 644 PetscFunctionBegin; 645 ierr = TaoComputeObjectiveAndGradient(tao,X,f,G);CHKERRQ(ierr); 646 ierr = LCLScatter(lclP,G,lclP->GU,lclP->GV);CHKERRQ(ierr); 647 if (lclP->recompute_jacobian_flag) { 648 ierr = TaoComputeJacobianState(tao,X,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv);CHKERRQ(ierr); 649 ierr = TaoComputeJacobianDesign(tao,X,tao->jacobian_design);CHKERRQ(ierr); 650 } 651 ierr = TaoComputeConstraints(tao,X, tao->constraints);CHKERRQ(ierr); 652 ierr = MatIsSymmetricKnown(tao->jacobian_state,&set,&flag);CHKERRQ(ierr); 653 if (tao->jacobian_state_pre) { 654 ierr = MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag);CHKERRQ(ierr); 655 } else { 656 pset = pflag = PETSC_TRUE; 657 } 658 if (set && pset && flag && pflag) symmetric = PETSC_TRUE; 659 else symmetric = PETSC_FALSE; 660 661 ierr = VecDot(lclP->lamda0, tao->constraints, &cdotl);CHKERRQ(ierr); 662 lclP->lgn = *f - cdotl; 663 664 /* Gradient of Lagrangian GL = G - J' * lamda */ 665 /* WU = A' * WL 666 WV = B' * WL */ 667 if (symmetric) { 668 ierr = MatMult(tao->jacobian_state,lclP->lamda0,lclP->GL_U);CHKERRQ(ierr); 669 } else { 670 ierr = MatMultTranspose(tao->jacobian_state,lclP->lamda0,lclP->GL_U);CHKERRQ(ierr); 671 } 672 ierr = MatMultTranspose(tao->jacobian_design,lclP->lamda0,lclP->GL_V);CHKERRQ(ierr); 673 ierr = VecScale(lclP->GL_U,-1.0);CHKERRQ(ierr); 674 ierr = VecScale(lclP->GL_V,-1.0);CHKERRQ(ierr); 675 ierr = VecAXPY(lclP->GL_U,1.0,lclP->GU);CHKERRQ(ierr); 676 ierr = VecAXPY(lclP->GL_V,1.0,lclP->GV);CHKERRQ(ierr); 677 ierr = LCLGather(lclP,lclP->GL_U,lclP->GL_V,lclP->GL);CHKERRQ(ierr); 678 679 f[0] = lclP->lgn; 680 ierr = VecCopy(lclP->GL,G); 681 PetscFunctionReturn(0); 682 } 683 684 #undef __FUNCT__ 685 #define __FUNCT__ "LCLComputeAugmentedLagrangianAndGradient" 686 static PetscErrorCode LCLComputeAugmentedLagrangianAndGradient(TaoLineSearch ls, Vec X, PetscReal *f, Vec G, void *ptr) 687 { 688 Tao tao = (Tao)ptr; 689 TAO_LCL *lclP = (TAO_LCL*)tao->data; 690 PetscReal con2; 691 PetscBool flag,pflag,set,pset,symmetric; 692 PetscErrorCode ierr; 693 694 PetscFunctionBegin; 695 ierr = LCLComputeLagrangianAndGradient(tao->linesearch,X,f,G,tao);CHKERRQ(ierr); 696 ierr = LCLScatter(lclP,G,lclP->GL_U,lclP->GL_V);CHKERRQ(ierr); 697 ierr = VecDot(tao->constraints,tao->constraints,&con2);CHKERRQ(ierr); 698 lclP->aug = lclP->lgn + 0.5*lclP->rho*con2; 699 700 /* Gradient of Aug. Lagrangian GAugL = GL + rho * J' c */ 701 /* WU = A' * c 702 WV = B' * c */ 703 ierr = MatIsSymmetricKnown(tao->jacobian_state,&set,&flag);CHKERRQ(ierr); 704 if (tao->jacobian_state_pre) { 705 ierr = MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag);CHKERRQ(ierr); 706 } else { 707 pset = pflag = PETSC_TRUE; 708 } 709 if (set && pset && flag && pflag) symmetric = PETSC_TRUE; 710 else symmetric = PETSC_FALSE; 711 712 if (symmetric) { 713 ierr = MatMult(tao->jacobian_state,tao->constraints,lclP->GAugL_U);CHKERRQ(ierr); 714 } else { 715 ierr = MatMultTranspose(tao->jacobian_state,tao->constraints,lclP->GAugL_U);CHKERRQ(ierr); 716 } 717 718 ierr = MatMultTranspose(tao->jacobian_design,tao->constraints,lclP->GAugL_V);CHKERRQ(ierr); 719 ierr = VecAYPX(lclP->GAugL_U,lclP->rho,lclP->GL_U);CHKERRQ(ierr); 720 ierr = VecAYPX(lclP->GAugL_V,lclP->rho,lclP->GL_V);CHKERRQ(ierr); 721 ierr = LCLGather(lclP,lclP->GAugL_U,lclP->GAugL_V,lclP->GAugL);CHKERRQ(ierr); 722 723 f[0] = lclP->aug; 724 ierr = VecCopy(lclP->GAugL,G);CHKERRQ(ierr); 725 PetscFunctionReturn(0); 726 } 727 728 #undef __FUNCT__ 729 #define __FUNCT__ "LCLGather" 730 PetscErrorCode LCLGather(TAO_LCL *lclP, Vec u, Vec v, Vec x) 731 { 732 PetscErrorCode ierr; 733 PetscFunctionBegin; 734 ierr = VecScatterBegin(lclP->state_scatter, u, x, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 735 ierr = VecScatterEnd(lclP->state_scatter, u, x, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 736 ierr = VecScatterBegin(lclP->design_scatter, v, x, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 737 ierr = VecScatterEnd(lclP->design_scatter, v, x, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 738 PetscFunctionReturn(0); 739 740 } 741 #undef __FUNCT__ 742 #define __FUNCT__ "LCLScatter" 743 PetscErrorCode LCLScatter(TAO_LCL *lclP, Vec x, Vec u, Vec v) 744 { 745 PetscErrorCode ierr; 746 PetscFunctionBegin; 747 ierr = VecScatterBegin(lclP->state_scatter, x, u, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 748 ierr = VecScatterEnd(lclP->state_scatter, x, u, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 749 ierr = VecScatterBegin(lclP->design_scatter, x, v, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 750 ierr = VecScatterEnd(lclP->design_scatter, x, v, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 751 PetscFunctionReturn(0); 752 753 } 754