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 PetscErrorCode ierr; 71 72 PetscFunctionBegin; 73 ierr = PetscOptionsHead("Linearly-Constrained Augmented Lagrangian Method for PDE-constrained optimization");CHKERRQ(ierr); 74 ierr = PetscOptionsReal("-tao_lcl_eps1","epsilon 1 tolerance","",lclP->eps1,&lclP->eps1,NULL);CHKERRQ(ierr); 75 ierr = PetscOptionsReal("-tao_lcl_eps2","epsilon 2 tolerance","",lclP->eps2,&lclP->eps2,NULL);CHKERRQ(ierr); 76 ierr = PetscOptionsReal("-tao_lcl_rho0","init value for rho","",lclP->rho0,&lclP->rho0,NULL);CHKERRQ(ierr); 77 ierr = PetscOptionsReal("-tao_lcl_rhomax","max value for rho","",lclP->rhomax,&lclP->rhomax,NULL);CHKERRQ(ierr); 78 lclP->phase2_niter = 1; 79 ierr = PetscOptionsInt("-tao_lcl_phase2_niter","Number of phase 2 iterations in LCL algorithm","",lclP->phase2_niter,&lclP->phase2_niter,NULL);CHKERRQ(ierr); 80 lclP->verbose = PETSC_FALSE; 81 ierr = PetscOptionsBool("-tao_lcl_verbose","Print verbose output","",lclP->verbose,&lclP->verbose,NULL);CHKERRQ(ierr); 82 lclP->tau[0] = lclP->tau[1] = lclP->tau[2] = lclP->tau[3] = 1.0e-4; 83 ierr = PetscOptionsReal("-tao_lcl_tola","Tolerance for first forward solve","",lclP->tau[0],&lclP->tau[0],NULL);CHKERRQ(ierr); 84 ierr = PetscOptionsReal("-tao_lcl_tolb","Tolerance for first adjoint solve","",lclP->tau[1],&lclP->tau[1],NULL);CHKERRQ(ierr); 85 ierr = PetscOptionsReal("-tao_lcl_tolc","Tolerance for second forward solve","",lclP->tau[2],&lclP->tau[2],NULL);CHKERRQ(ierr); 86 ierr = PetscOptionsReal("-tao_lcl_told","Tolerance for second adjoint solve","",lclP->tau[3],&lclP->tau[3],NULL);CHKERRQ(ierr); 87 ierr = PetscOptionsTail();CHKERRQ(ierr); 88 ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 89 PetscFunctionReturn(0); 90 } 91 92 #undef __FUNCT__ 93 #define __FUNCT__ "TaoView_LCL" 94 static PetscErrorCode TaoView_LCL(Tao tao, PetscViewer viewer) 95 { 96 return 0; 97 } 98 99 #undef __FUNCT__ 100 #define __FUNCT__ "TaoSetup_LCL" 101 static PetscErrorCode TaoSetup_LCL(Tao tao) 102 { 103 TAO_LCL *lclP = (TAO_LCL*)tao->data; 104 PetscInt lo, hi, nlocalstate, nlocaldesign; 105 PetscErrorCode ierr; 106 IS is_state, is_design; 107 108 PetscFunctionBegin; 109 if (!tao->state_is) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONGSTATE,"LCL Solver requires an initial state index set -- use TaoSetStateIS()"); 110 ierr = VecDuplicate(tao->solution, &tao->gradient);CHKERRQ(ierr); 111 ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr); 112 ierr = VecDuplicate(tao->solution, &lclP->W);CHKERRQ(ierr); 113 ierr = VecDuplicate(tao->solution, &lclP->X0);CHKERRQ(ierr); 114 ierr = VecDuplicate(tao->solution, &lclP->G0);CHKERRQ(ierr); 115 ierr = VecDuplicate(tao->solution, &lclP->GL);CHKERRQ(ierr); 116 ierr = VecDuplicate(tao->solution, &lclP->GAugL);CHKERRQ(ierr); 117 118 ierr = VecDuplicate(tao->constraints, &lclP->lamda);CHKERRQ(ierr); 119 ierr = VecDuplicate(tao->constraints, &lclP->WL);CHKERRQ(ierr); 120 ierr = VecDuplicate(tao->constraints, &lclP->lamda0);CHKERRQ(ierr); 121 ierr = VecDuplicate(tao->constraints, &lclP->con1);CHKERRQ(ierr); 122 123 ierr = VecSet(lclP->lamda,0.0);CHKERRQ(ierr); 124 125 ierr = VecGetSize(tao->solution, &lclP->n);CHKERRQ(ierr); 126 ierr = VecGetSize(tao->constraints, &lclP->m);CHKERRQ(ierr); 127 128 ierr = VecCreate(((PetscObject)tao)->comm,&lclP->U);CHKERRQ(ierr); 129 ierr = VecCreate(((PetscObject)tao)->comm,&lclP->V);CHKERRQ(ierr); 130 ierr = ISGetLocalSize(tao->state_is,&nlocalstate);CHKERRQ(ierr); 131 ierr = ISGetLocalSize(tao->design_is,&nlocaldesign);CHKERRQ(ierr); 132 ierr = VecSetSizes(lclP->U,nlocalstate,lclP->m);CHKERRQ(ierr); 133 ierr = VecSetSizes(lclP->V,nlocaldesign,lclP->n-lclP->m);CHKERRQ(ierr); 134 ierr = VecSetType(lclP->U,((PetscObject)(tao->solution))->type_name);CHKERRQ(ierr); 135 ierr = VecSetType(lclP->V,((PetscObject)(tao->solution))->type_name);CHKERRQ(ierr); 136 ierr = VecSetFromOptions(lclP->U);CHKERRQ(ierr); 137 ierr = VecSetFromOptions(lclP->V);CHKERRQ(ierr); 138 ierr = VecDuplicate(lclP->U,&lclP->DU);CHKERRQ(ierr); 139 ierr = VecDuplicate(lclP->U,&lclP->U0);CHKERRQ(ierr); 140 ierr = VecDuplicate(lclP->U,&lclP->GU);CHKERRQ(ierr); 141 ierr = VecDuplicate(lclP->U,&lclP->GU0);CHKERRQ(ierr); 142 ierr = VecDuplicate(lclP->U,&lclP->GAugL_U);CHKERRQ(ierr); 143 ierr = VecDuplicate(lclP->U,&lclP->GL_U);CHKERRQ(ierr); 144 ierr = VecDuplicate(lclP->U,&lclP->GAugL_U0);CHKERRQ(ierr); 145 ierr = VecDuplicate(lclP->U,&lclP->GL_U0);CHKERRQ(ierr); 146 ierr = VecDuplicate(lclP->U,&lclP->WU);CHKERRQ(ierr); 147 ierr = VecDuplicate(lclP->U,&lclP->r);CHKERRQ(ierr); 148 ierr = VecDuplicate(lclP->V,&lclP->V0);CHKERRQ(ierr); 149 ierr = VecDuplicate(lclP->V,&lclP->V1);CHKERRQ(ierr); 150 ierr = VecDuplicate(lclP->V,&lclP->DV);CHKERRQ(ierr); 151 ierr = VecDuplicate(lclP->V,&lclP->s);CHKERRQ(ierr); 152 ierr = VecDuplicate(lclP->V,&lclP->GV);CHKERRQ(ierr); 153 ierr = VecDuplicate(lclP->V,&lclP->GV0);CHKERRQ(ierr); 154 ierr = VecDuplicate(lclP->V,&lclP->dbar);CHKERRQ(ierr); 155 ierr = VecDuplicate(lclP->V,&lclP->GAugL_V);CHKERRQ(ierr); 156 ierr = VecDuplicate(lclP->V,&lclP->GL_V);CHKERRQ(ierr); 157 ierr = VecDuplicate(lclP->V,&lclP->GAugL_V0);CHKERRQ(ierr); 158 ierr = VecDuplicate(lclP->V,&lclP->GL_V0);CHKERRQ(ierr); 159 ierr = VecDuplicate(lclP->V,&lclP->WV);CHKERRQ(ierr); 160 ierr = VecDuplicate(lclP->V,&lclP->g1);CHKERRQ(ierr); 161 ierr = VecDuplicate(lclP->V,&lclP->g2);CHKERRQ(ierr); 162 163 /* create scatters for state, design subvecs */ 164 ierr = VecGetOwnershipRange(lclP->U,&lo,&hi);CHKERRQ(ierr); 165 ierr = ISCreateStride(((PetscObject)lclP->U)->comm,hi-lo,lo,1,&is_state);CHKERRQ(ierr); 166 ierr = VecGetOwnershipRange(lclP->V,&lo,&hi);CHKERRQ(ierr); 167 if (0) { 168 PetscInt sizeU,sizeV; 169 ierr = VecGetSize(lclP->U,&sizeU); 170 ierr = VecGetSize(lclP->V,&sizeV); 171 ierr = PetscPrintf(PETSC_COMM_WORLD,"size(U)=%D, size(V)=%D\n",sizeU,sizeV); 172 } 173 ierr = ISCreateStride(((PetscObject)lclP->V)->comm,hi-lo,lo,1,&is_design);CHKERRQ(ierr); 174 ierr = VecScatterCreate(tao->solution,tao->state_is,lclP->U,is_state,&lclP->state_scatter);CHKERRQ(ierr); 175 ierr = VecScatterCreate(tao->solution,tao->design_is,lclP->V,is_design,&lclP->design_scatter);CHKERRQ(ierr); 176 ierr = ISDestroy(&is_state);CHKERRQ(ierr); 177 ierr = ISDestroy(&is_design);CHKERRQ(ierr); 178 PetscFunctionReturn(0); 179 } 180 181 #undef __FUNCT__ 182 #define __FUNCT__ "TaoSolve_LCL" 183 static PetscErrorCode TaoSolve_LCL(Tao tao) 184 { 185 TAO_LCL *lclP = (TAO_LCL*)tao->data; 186 PetscInt iter=0,phase2_iter,nlocal,its; 187 TaoConvergedReason reason = TAO_CONTINUE_ITERATING; 188 TaoLineSearchConvergedReason ls_reason = TAOLINESEARCH_CONTINUE_ITERATING; 189 PetscReal step=1.0,f, descent, aldescent; 190 PetscReal cnorm, mnorm; 191 PetscReal adec,r2,rGL_U,rWU; 192 PetscBool set,pset,flag,pflag,symmetric; 193 PetscErrorCode ierr; 194 195 PetscFunctionBegin; 196 lclP->rho = lclP->rho0; 197 ierr = VecGetLocalSize(lclP->U,&nlocal);CHKERRQ(ierr); 198 ierr = VecGetLocalSize(lclP->V,&nlocal);CHKERRQ(ierr); 199 ierr = MatCreateLMVM(((PetscObject)tao)->comm,nlocal,lclP->n-lclP->m,&lclP->R);CHKERRQ(ierr); 200 ierr = MatLMVMAllocateVectors(lclP->R,lclP->V);CHKERRQ(ierr); 201 lclP->recompute_jacobian_flag = PETSC_TRUE; 202 203 /* Scatter to U,V */ 204 ierr = LCLScatter(lclP,tao->solution,lclP->U,lclP->V);CHKERRQ(ierr); 205 206 /* Evaluate Function, Gradient, Constraints, and Jacobian */ 207 ierr = TaoComputeObjectiveAndGradient(tao,tao->solution,&f,tao->gradient);CHKERRQ(ierr); 208 ierr = TaoComputeJacobianState(tao,tao->solution,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv);CHKERRQ(ierr); 209 ierr = TaoComputeJacobianDesign(tao,tao->solution,tao->jacobian_design);CHKERRQ(ierr); 210 ierr = TaoComputeConstraints(tao,tao->solution, tao->constraints);CHKERRQ(ierr); 211 212 /* Scatter gradient to GU,GV */ 213 ierr = LCLScatter(lclP,tao->gradient,lclP->GU,lclP->GV);CHKERRQ(ierr); 214 215 /* Evaluate Lagrangian function and gradient */ 216 /* p0 */ 217 ierr = VecSet(lclP->lamda,0.0);CHKERRQ(ierr); /* Initial guess in CG */ 218 ierr = MatIsSymmetricKnown(tao->jacobian_state,&set,&flag);CHKERRQ(ierr); 219 if (tao->jacobian_state_pre) { 220 ierr = MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag);CHKERRQ(ierr); 221 } else { 222 pset = pflag = PETSC_TRUE; 223 } 224 if (set && pset && flag && pflag) symmetric = PETSC_TRUE; 225 else symmetric = PETSC_FALSE; 226 227 lclP->solve_type = LCL_ADJOINT2; 228 if (tao->jacobian_state_inv) { 229 if (symmetric) { 230 ierr = MatMult(tao->jacobian_state_inv, lclP->GU, lclP->lamda);CHKERRQ(ierr); } else { 231 ierr = MatMultTranspose(tao->jacobian_state_inv, lclP->GU, lclP->lamda);CHKERRQ(ierr); 232 } 233 } else { 234 ierr = KSPSetOperators(tao->ksp, tao->jacobian_state, tao->jacobian_state_pre);CHKERRQ(ierr); 235 if (symmetric) { 236 ierr = KSPSolve(tao->ksp, lclP->GU, lclP->lamda);CHKERRQ(ierr); 237 } else { 238 ierr = KSPSolveTranspose(tao->ksp, lclP->GU, lclP->lamda);CHKERRQ(ierr); 239 } 240 ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 241 tao->ksp_its += its; 242 } 243 ierr = VecCopy(lclP->lamda,lclP->lamda0);CHKERRQ(ierr); 244 ierr = LCLComputeAugmentedLagrangianAndGradient(tao->linesearch,tao->solution,&lclP->aug,lclP->GAugL,tao);CHKERRQ(ierr); 245 246 ierr = LCLScatter(lclP,lclP->GL,lclP->GL_U,lclP->GL_V);CHKERRQ(ierr); 247 ierr = LCLScatter(lclP,lclP->GAugL,lclP->GAugL_U,lclP->GAugL_V);CHKERRQ(ierr); 248 249 /* Evaluate constraint norm */ 250 ierr = VecNorm(tao->constraints, NORM_2, &cnorm);CHKERRQ(ierr); 251 ierr = VecNorm(lclP->GAugL, NORM_2, &mnorm);CHKERRQ(ierr); 252 253 /* Monitor convergence */ 254 ierr = TaoMonitor(tao, iter,f,mnorm,cnorm,step,&reason);CHKERRQ(ierr); 255 256 while (reason == TAO_CONTINUE_ITERATING) { 257 /* Compute a descent direction for the linearly constrained subproblem 258 minimize f(u+du, v+dv) 259 s.t. A(u0,v0)du + B(u0,v0)dv = -g(u0,v0) */ 260 261 /* Store the points around the linearization */ 262 ierr = VecCopy(lclP->U, lclP->U0);CHKERRQ(ierr); 263 ierr = VecCopy(lclP->V, lclP->V0);CHKERRQ(ierr); 264 ierr = VecCopy(lclP->GU,lclP->GU0);CHKERRQ(ierr); 265 ierr = VecCopy(lclP->GV,lclP->GV0);CHKERRQ(ierr); 266 ierr = VecCopy(lclP->GAugL_U,lclP->GAugL_U0);CHKERRQ(ierr); 267 ierr = VecCopy(lclP->GAugL_V,lclP->GAugL_V0);CHKERRQ(ierr); 268 ierr = VecCopy(lclP->GL_U,lclP->GL_U0);CHKERRQ(ierr); 269 ierr = VecCopy(lclP->GL_V,lclP->GL_V0);CHKERRQ(ierr); 270 271 lclP->aug0 = lclP->aug; 272 lclP->lgn0 = lclP->lgn; 273 274 /* Given the design variables, we need to project the current iterate 275 onto the linearized constraint. We choose to fix the design variables 276 and solve the linear system for the state variables. The resulting 277 point is the Newton direction */ 278 279 /* Solve r = A\con */ 280 lclP->solve_type = LCL_FORWARD1; 281 ierr = VecSet(lclP->r,0.0);CHKERRQ(ierr); /* Initial guess in CG */ 282 283 if (tao->jacobian_state_inv) { 284 ierr = MatMult(tao->jacobian_state_inv, tao->constraints, lclP->r);CHKERRQ(ierr); 285 } else { 286 ierr = KSPSetOperators(tao->ksp, tao->jacobian_state, tao->jacobian_state_pre);CHKERRQ(ierr); 287 ierr = KSPSolve(tao->ksp, tao->constraints, lclP->r);CHKERRQ(ierr); 288 ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 289 tao->ksp_its += its; 290 } 291 292 /* Set design step direction dv to zero */ 293 ierr = VecSet(lclP->s, 0.0);CHKERRQ(ierr); 294 295 /* 296 Check sufficient descent for constraint merit function .5*||con||^2 297 con' Ak r >= eps1 ||r||^(2+eps2) 298 */ 299 300 /* Compute WU= Ak' * con */ 301 if (symmetric) { 302 ierr = MatMult(tao->jacobian_state,tao->constraints,lclP->WU);CHKERRQ(ierr); 303 } else { 304 ierr = MatMultTranspose(tao->jacobian_state,tao->constraints,lclP->WU);CHKERRQ(ierr); 305 } 306 /* Compute r * Ak' * con */ 307 ierr = VecDot(lclP->r,lclP->WU,&rWU);CHKERRQ(ierr); 308 309 /* compute ||r||^(2+eps2) */ 310 ierr = VecNorm(lclP->r,NORM_2,&r2);CHKERRQ(ierr); 311 r2 = PetscPowScalar(r2,2.0+lclP->eps2); 312 adec = lclP->eps1 * r2; 313 314 if (rWU < adec) { 315 ierr = PetscInfo(tao,"Newton direction not descent for constraint, feasibility phase required");CHKERRQ(ierr); 316 if (lclP->verbose) { 317 ierr = PetscPrintf(PETSC_COMM_WORLD,"Newton direction not descent for constraint: %g -- using steepest descent\n",(double)descent);CHKERRQ(ierr); 318 } 319 320 ierr = PetscInfo(tao,"Using steepest descent direction instead.\n");CHKERRQ(ierr); 321 ierr = VecSet(lclP->r,0.0);CHKERRQ(ierr); 322 ierr = VecAXPY(lclP->r,-1.0,lclP->WU);CHKERRQ(ierr); 323 ierr = VecDot(lclP->r,lclP->r,&rWU);CHKERRQ(ierr); 324 ierr = VecNorm(lclP->r,NORM_2,&r2);CHKERRQ(ierr); 325 r2 = PetscPowScalar(r2,2.0+lclP->eps2); 326 ierr = VecDot(lclP->r,lclP->GAugL_U,&descent);CHKERRQ(ierr); 327 adec = lclP->eps1 * r2; 328 } 329 330 331 /* 332 Check descent for aug. lagrangian 333 r' (GUk - Ak'*yk - rho*Ak'*con) <= -eps1 ||r||^(2+eps2) 334 GL_U = GUk - Ak'*yk 335 WU = Ak'*con 336 adec=eps1||r||^(2+eps2) 337 338 ==> 339 Check r'GL_U - rho*r'WU <= adec 340 */ 341 342 ierr = VecDot(lclP->r,lclP->GL_U,&rGL_U); 343 aldescent = rGL_U - lclP->rho*rWU; 344 if (aldescent > -adec) { 345 if (lclP->verbose) { 346 ierr = PetscPrintf(PETSC_COMM_WORLD," Newton direction not descent for augmented Lagrangian: %g",(double)aldescent);CHKERRQ(ierr); 347 } 348 ierr = PetscInfo1(tao,"Newton direction not descent for augmented Lagrangian: %g",(double)aldescent);CHKERRQ(ierr); 349 lclP->rho = (rGL_U - adec)/rWU; 350 if (lclP->rho > lclP->rhomax) { 351 lclP->rho = lclP->rhomax; 352 SETERRQ1(PETSC_COMM_WORLD,0,"rho=%g > rhomax, case not implemented. Increase rhomax (-tao_lcl_rhomax)",(double)lclP->rho); 353 } 354 if (lclP->verbose) { 355 ierr = PetscPrintf(PETSC_COMM_WORLD," Increasing penalty parameter to %g\n",(double)lclP->rho);CHKERRQ(ierr); 356 } 357 ierr = PetscInfo1(tao," Increasing penalty parameter to %g",(double)lclP->rho); 358 } 359 360 ierr = LCLComputeAugmentedLagrangianAndGradient(tao->linesearch,tao->solution,&lclP->aug,lclP->GAugL,tao);CHKERRQ(ierr); 361 ierr = LCLScatter(lclP,lclP->GL,lclP->GL_U,lclP->GL_V);CHKERRQ(ierr); 362 ierr = LCLScatter(lclP,lclP->GAugL,lclP->GAugL_U,lclP->GAugL_V);CHKERRQ(ierr); 363 364 /* We now minimize the augmented Lagrangian along the Newton direction */ 365 ierr = VecScale(lclP->r,-1.0);CHKERRQ(ierr); 366 ierr = LCLGather(lclP, lclP->r,lclP->s,tao->stepdirection);CHKERRQ(ierr); 367 ierr = VecScale(lclP->r,-1.0);CHKERRQ(ierr); 368 ierr = LCLGather(lclP, lclP->GAugL_U0, lclP->GAugL_V0, lclP->GAugL);CHKERRQ(ierr); 369 ierr = LCLGather(lclP, lclP->U0,lclP->V0,lclP->X0);CHKERRQ(ierr); 370 371 lclP->recompute_jacobian_flag = PETSC_TRUE; 372 373 ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr); 374 ierr = TaoLineSearchSetType(tao->linesearch, TAOLINESEARCHMT);CHKERRQ(ierr); 375 ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 376 ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &lclP->aug, lclP->GAugL, tao->stepdirection, &step, &ls_reason);CHKERRQ(ierr); 377 if (lclP->verbose) { 378 ierr = PetscPrintf(PETSC_COMM_WORLD,"Steplength = %g\n",(double)step);CHKERRQ(ierr); 379 } 380 381 ierr = LCLScatter(lclP,tao->solution,lclP->U,lclP->V);CHKERRQ(ierr); 382 ierr = TaoComputeObjectiveAndGradient(tao,tao->solution,&f,tao->gradient);CHKERRQ(ierr); 383 ierr = LCLScatter(lclP,tao->gradient,lclP->GU,lclP->GV);CHKERRQ(ierr); 384 385 ierr = LCLScatter(lclP,lclP->GAugL,lclP->GAugL_U,lclP->GAugL_V);CHKERRQ(ierr); 386 387 /* Check convergence */ 388 ierr = VecNorm(lclP->GAugL, NORM_2, &mnorm);CHKERRQ(ierr); 389 ierr = VecNorm(tao->constraints, NORM_2, &cnorm);CHKERRQ(ierr); 390 ierr = TaoMonitor(tao, iter,f,mnorm,cnorm,step,&reason);CHKERRQ(ierr); 391 if (reason != TAO_CONTINUE_ITERATING){ 392 break; 393 } 394 395 /* TODO: use a heuristic to choose how many iterations should be performed within phase 2 */ 396 for (phase2_iter=0; phase2_iter<lclP->phase2_niter; phase2_iter++){ 397 /* We now minimize the objective function starting from the fraction of 398 the Newton point accepted by applying one step of a reduced-space 399 method. The optimization problem is: 400 401 minimize f(u+du, v+dv) 402 s. t. A(u0,v0)du + B(u0,v0)du = -alpha g(u0,v0) 403 404 In particular, we have that 405 du = -inv(A)*(Bdv + alpha g) */ 406 407 ierr = TaoComputeJacobianState(tao,lclP->X0,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv);CHKERRQ(ierr); 408 ierr = TaoComputeJacobianDesign(tao,lclP->X0,tao->jacobian_design);CHKERRQ(ierr); 409 410 /* Store V and constraints */ 411 ierr = VecCopy(lclP->V, lclP->V1);CHKERRQ(ierr); 412 ierr = VecCopy(tao->constraints,lclP->con1);CHKERRQ(ierr); 413 414 /* Compute multipliers */ 415 /* p1 */ 416 ierr = VecSet(lclP->lamda,0.0);CHKERRQ(ierr); /* Initial guess in CG */ 417 lclP->solve_type = LCL_ADJOINT1; 418 ierr = MatIsSymmetricKnown(tao->jacobian_state,&set,&flag);CHKERRQ(ierr); 419 if (tao->jacobian_state_pre) { 420 ierr = MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag);CHKERRQ(ierr); 421 } else { 422 pset = pflag = PETSC_TRUE; 423 } 424 if (set && pset && flag && pflag) symmetric = PETSC_TRUE; 425 else symmetric = PETSC_FALSE; 426 427 if (tao->jacobian_state_inv) { 428 if (symmetric) { 429 ierr = MatMult(tao->jacobian_state_inv, lclP->GAugL_U, lclP->lamda);CHKERRQ(ierr); 430 } else { 431 ierr = MatMultTranspose(tao->jacobian_state_inv, lclP->GAugL_U, lclP->lamda);CHKERRQ(ierr); 432 } 433 } else { 434 if (symmetric) { 435 ierr = KSPSolve(tao->ksp, lclP->GAugL_U, lclP->lamda);CHKERRQ(ierr); 436 } else { 437 ierr = KSPSolveTranspose(tao->ksp, lclP->GAugL_U, lclP->lamda);CHKERRQ(ierr); 438 } 439 ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 440 tao->ksp_its += its; 441 } 442 ierr = MatMultTranspose(tao->jacobian_design,lclP->lamda,lclP->g1);CHKERRQ(ierr); 443 ierr = VecAXPY(lclP->g1,-1.0,lclP->GAugL_V);CHKERRQ(ierr); 444 445 /* Compute the limited-memory quasi-newton direction */ 446 if (iter > 0) { 447 ierr = MatLMVMSolve(lclP->R,lclP->g1,lclP->s);CHKERRQ(ierr); 448 ierr = VecDot(lclP->s,lclP->g1,&descent);CHKERRQ(ierr); 449 if (descent <= 0) { 450 if (lclP->verbose) { 451 ierr = PetscPrintf(PETSC_COMM_WORLD,"Reduced-space direction not descent: %g\n",(double)descent);CHKERRQ(ierr); 452 } 453 ierr = VecCopy(lclP->g1,lclP->s);CHKERRQ(ierr); 454 } 455 } else { 456 ierr = VecCopy(lclP->g1,lclP->s);CHKERRQ(ierr); 457 } 458 ierr = VecScale(lclP->g1,-1.0);CHKERRQ(ierr); 459 460 /* Recover the full space direction */ 461 ierr = MatMult(tao->jacobian_design,lclP->s,lclP->WU);CHKERRQ(ierr); 462 /* ierr = VecSet(lclP->r,0.0);CHKERRQ(ierr); */ /* Initial guess in CG */ 463 lclP->solve_type = LCL_FORWARD2; 464 if (tao->jacobian_state_inv) { 465 ierr = MatMult(tao->jacobian_state_inv,lclP->WU,lclP->r);CHKERRQ(ierr); 466 } else { 467 ierr = KSPSolve(tao->ksp, lclP->WU, lclP->r);CHKERRQ(ierr); 468 ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 469 tao->ksp_its += its; 470 } 471 472 /* We now minimize the augmented Lagrangian along the direction -r,s */ 473 ierr = VecScale(lclP->r, -1.0);CHKERRQ(ierr); 474 ierr = LCLGather(lclP,lclP->r,lclP->s,tao->stepdirection);CHKERRQ(ierr); 475 ierr = VecScale(lclP->r, -1.0);CHKERRQ(ierr); 476 lclP->recompute_jacobian_flag = PETSC_TRUE; 477 478 ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr); 479 ierr = TaoLineSearchSetType(tao->linesearch,TAOLINESEARCHMT);CHKERRQ(ierr); 480 ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 481 ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &lclP->aug, lclP->GAugL, tao->stepdirection,&step,&ls_reason);CHKERRQ(ierr); 482 if (lclP->verbose){ 483 ierr = PetscPrintf(PETSC_COMM_WORLD,"Reduced-space steplength = %g\n",(double)step);CHKERRQ(ierr); 484 } 485 486 ierr = LCLScatter(lclP,tao->solution,lclP->U,lclP->V);CHKERRQ(ierr); 487 ierr = LCLScatter(lclP,lclP->GL,lclP->GL_U,lclP->GL_V);CHKERRQ(ierr); 488 ierr = LCLScatter(lclP,lclP->GAugL,lclP->GAugL_U,lclP->GAugL_V);CHKERRQ(ierr); 489 ierr = TaoComputeObjectiveAndGradient(tao,tao->solution,&f,tao->gradient);CHKERRQ(ierr); 490 ierr = LCLScatter(lclP,tao->gradient,lclP->GU,lclP->GV);CHKERRQ(ierr); 491 492 /* Compute the reduced gradient at the new point */ 493 494 ierr = TaoComputeJacobianState(tao,lclP->X0,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv);CHKERRQ(ierr); 495 ierr = TaoComputeJacobianDesign(tao,lclP->X0,tao->jacobian_design);CHKERRQ(ierr); 496 497 /* p2 */ 498 /* Compute multipliers, using lamda-rho*con as an initial guess in PCG */ 499 if (phase2_iter==0){ 500 ierr = VecWAXPY(lclP->lamda,-lclP->rho,lclP->con1,lclP->lamda0);CHKERRQ(ierr); 501 } else { 502 ierr = VecAXPY(lclP->lamda,-lclP->rho,tao->constraints);CHKERRQ(ierr); 503 } 504 505 ierr = MatIsSymmetricKnown(tao->jacobian_state,&set,&flag);CHKERRQ(ierr); 506 if (tao->jacobian_state_pre) { 507 ierr = MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag);CHKERRQ(ierr); 508 } else { 509 pset = pflag = PETSC_TRUE; 510 } 511 if (set && pset && flag && pflag) symmetric = PETSC_TRUE; 512 else symmetric = PETSC_FALSE; 513 514 lclP->solve_type = LCL_ADJOINT2; 515 if (tao->jacobian_state_inv) { 516 if (symmetric) { 517 ierr = MatMult(tao->jacobian_state_inv, lclP->GU, lclP->lamda);CHKERRQ(ierr); 518 } else { 519 ierr = MatMultTranspose(tao->jacobian_state_inv, lclP->GU, lclP->lamda);CHKERRQ(ierr); 520 } 521 } else { 522 if (symmetric) { 523 ierr = KSPSolve(tao->ksp, lclP->GU, lclP->lamda);CHKERRQ(ierr); 524 } else { 525 ierr = KSPSolveTranspose(tao->ksp, lclP->GU, lclP->lamda);CHKERRQ(ierr); 526 } 527 ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 528 tao->ksp_its += its; 529 } 530 531 ierr = MatMultTranspose(tao->jacobian_design,lclP->lamda,lclP->g2);CHKERRQ(ierr); 532 ierr = VecAXPY(lclP->g2,-1.0,lclP->GV);CHKERRQ(ierr); 533 534 ierr = VecScale(lclP->g2,-1.0);CHKERRQ(ierr); 535 536 /* Update the quasi-newton approximation */ 537 if (phase2_iter >= 0){ 538 ierr = MatLMVMSetPrev(lclP->R,lclP->V1,lclP->g1); 539 } 540 ierr = MatLMVMUpdate(lclP->R,lclP->V,lclP->g2);CHKERRQ(ierr); 541 /* 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 */ 542 543 } 544 545 ierr = VecCopy(lclP->lamda,lclP->lamda0);CHKERRQ(ierr); 546 547 /* Evaluate Function, Gradient, Constraints, and Jacobian */ 548 ierr = TaoComputeObjectiveAndGradient(tao,tao->solution,&f,tao->gradient);CHKERRQ(ierr); 549 ierr = LCLScatter(lclP,tao->solution,lclP->U,lclP->V);CHKERRQ(ierr); 550 ierr = LCLScatter(lclP,tao->gradient,lclP->GU,lclP->GV);CHKERRQ(ierr); 551 552 ierr = TaoComputeJacobianState(tao,tao->solution,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv);CHKERRQ(ierr); 553 ierr = TaoComputeJacobianDesign(tao,tao->solution,tao->jacobian_design);CHKERRQ(ierr); 554 ierr = TaoComputeConstraints(tao,tao->solution, tao->constraints);CHKERRQ(ierr); 555 556 ierr = LCLComputeAugmentedLagrangianAndGradient(tao->linesearch,tao->solution,&lclP->aug,lclP->GAugL,tao);CHKERRQ(ierr); 557 558 ierr = VecNorm(lclP->GAugL, NORM_2, &mnorm);CHKERRQ(ierr); 559 560 /* Evaluate constraint norm */ 561 ierr = VecNorm(tao->constraints, NORM_2, &cnorm);CHKERRQ(ierr); 562 563 /* Monitor convergence */ 564 iter++; 565 ierr = TaoMonitor(tao, iter,f,mnorm,cnorm,step,&reason);CHKERRQ(ierr); 566 } 567 ierr = MatDestroy(&lclP->R);CHKERRQ(ierr); 568 PetscFunctionReturn(0); 569 } 570 571 /*MC 572 TAOLCL - linearly constrained lagrangian method for pde-constrained optimization 573 574 + -tao_lcl_eps1 - epsilon 1 tolerance 575 . -tao_lcl_eps2","epsilon 2 tolerance","",lclP->eps2,&lclP->eps2,NULL);CHKERRQ(ierr); 576 . -tao_lcl_rho0","init value for rho","",lclP->rho0,&lclP->rho0,NULL);CHKERRQ(ierr); 577 . -tao_lcl_rhomax","max value for rho","",lclP->rhomax,&lclP->rhomax,NULL);CHKERRQ(ierr); 578 . -tao_lcl_phase2_niter - Number of phase 2 iterations in LCL algorithm 579 . -tao_lcl_verbose - Print verbose output if True 580 . -tao_lcl_tola - Tolerance for first forward solve 581 . -tao_lcl_tolb - Tolerance for first adjoint solve 582 . -tao_lcl_tolc - Tolerance for second forward solve 583 - -tao_lcl_told - Tolerance for second adjoint solve 584 585 Level: beginner 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