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);CHKERRQ(ierr); 62 PetscFunctionReturn(0); 63 } 64 65 #undef __FUNCT__ 66 #define __FUNCT__ "TaoSetFromOptions_LCL" 67 static PetscErrorCode TaoSetFromOptions_LCL(PetscOptionItems *PetscOptionsObject,Tao tao) 68 { 69 TAO_LCL *lclP = (TAO_LCL*)tao->data; 70 PetscErrorCode ierr; 71 72 PetscFunctionBegin; 73 ierr = PetscOptionsHead(PetscOptionsObject,"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 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 tao->ksp_tot_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, tao->niter,f,mnorm,cnorm,step,&reason);CHKERRQ(ierr); 256 257 while (reason == TAO_CONTINUE_ITERATING) { 258 tao->ksp_its=0; 259 /* Compute a descent direction for the linearly constrained subproblem 260 minimize f(u+du, v+dv) 261 s.t. A(u0,v0)du + B(u0,v0)dv = -g(u0,v0) */ 262 263 /* Store the points around the linearization */ 264 ierr = VecCopy(lclP->U, lclP->U0);CHKERRQ(ierr); 265 ierr = VecCopy(lclP->V, lclP->V0);CHKERRQ(ierr); 266 ierr = VecCopy(lclP->GU,lclP->GU0);CHKERRQ(ierr); 267 ierr = VecCopy(lclP->GV,lclP->GV0);CHKERRQ(ierr); 268 ierr = VecCopy(lclP->GAugL_U,lclP->GAugL_U0);CHKERRQ(ierr); 269 ierr = VecCopy(lclP->GAugL_V,lclP->GAugL_V0);CHKERRQ(ierr); 270 ierr = VecCopy(lclP->GL_U,lclP->GL_U0);CHKERRQ(ierr); 271 ierr = VecCopy(lclP->GL_V,lclP->GL_V0);CHKERRQ(ierr); 272 273 lclP->aug0 = lclP->aug; 274 lclP->lgn0 = lclP->lgn; 275 276 /* Given the design variables, we need to project the current iterate 277 onto the linearized constraint. We choose to fix the design variables 278 and solve the linear system for the state variables. The resulting 279 point is the Newton direction */ 280 281 /* Solve r = A\con */ 282 lclP->solve_type = LCL_FORWARD1; 283 ierr = VecSet(lclP->r,0.0);CHKERRQ(ierr); /* Initial guess in CG */ 284 285 if (tao->jacobian_state_inv) { 286 ierr = MatMult(tao->jacobian_state_inv, tao->constraints, lclP->r);CHKERRQ(ierr); 287 } else { 288 ierr = KSPSetOperators(tao->ksp, tao->jacobian_state, tao->jacobian_state_pre);CHKERRQ(ierr); 289 ierr = KSPSolve(tao->ksp, tao->constraints, lclP->r);CHKERRQ(ierr); 290 ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 291 tao->ksp_its+=its; 292 tao->ksp_tot_its+=tao->ksp_its; 293 } 294 295 /* Set design step direction dv to zero */ 296 ierr = VecSet(lclP->s, 0.0);CHKERRQ(ierr); 297 298 /* 299 Check sufficient descent for constraint merit function .5*||con||^2 300 con' Ak r >= eps1 ||r||^(2+eps2) 301 */ 302 303 /* Compute WU= Ak' * con */ 304 if (symmetric) { 305 ierr = MatMult(tao->jacobian_state,tao->constraints,lclP->WU);CHKERRQ(ierr); 306 } else { 307 ierr = MatMultTranspose(tao->jacobian_state,tao->constraints,lclP->WU);CHKERRQ(ierr); 308 } 309 /* Compute r * Ak' * con */ 310 ierr = VecDot(lclP->r,lclP->WU,&rWU);CHKERRQ(ierr); 311 312 /* compute ||r||^(2+eps2) */ 313 ierr = VecNorm(lclP->r,NORM_2,&r2);CHKERRQ(ierr); 314 r2 = PetscPowScalar(r2,2.0+lclP->eps2); 315 adec = lclP->eps1 * r2; 316 317 if (rWU < adec) { 318 ierr = PetscInfo(tao,"Newton direction not descent for constraint, feasibility phase required\n");CHKERRQ(ierr); 319 if (lclP->verbose) { 320 ierr = PetscPrintf(PETSC_COMM_WORLD,"Newton direction not descent for constraint: %g -- using steepest descent\n",(double)descent);CHKERRQ(ierr); 321 } 322 323 ierr = PetscInfo(tao,"Using steepest descent direction instead.\n");CHKERRQ(ierr); 324 ierr = VecSet(lclP->r,0.0);CHKERRQ(ierr); 325 ierr = VecAXPY(lclP->r,-1.0,lclP->WU);CHKERRQ(ierr); 326 ierr = VecDot(lclP->r,lclP->r,&rWU);CHKERRQ(ierr); 327 ierr = VecNorm(lclP->r,NORM_2,&r2);CHKERRQ(ierr); 328 r2 = PetscPowScalar(r2,2.0+lclP->eps2); 329 ierr = VecDot(lclP->r,lclP->GAugL_U,&descent);CHKERRQ(ierr); 330 adec = lclP->eps1 * r2; 331 } 332 333 334 /* 335 Check descent for aug. lagrangian 336 r' (GUk - Ak'*yk - rho*Ak'*con) <= -eps1 ||r||^(2+eps2) 337 GL_U = GUk - Ak'*yk 338 WU = Ak'*con 339 adec=eps1||r||^(2+eps2) 340 341 ==> 342 Check r'GL_U - rho*r'WU <= adec 343 */ 344 345 ierr = VecDot(lclP->r,lclP->GL_U,&rGL_U);CHKERRQ(ierr); 346 aldescent = rGL_U - lclP->rho*rWU; 347 if (aldescent > -adec) { 348 if (lclP->verbose) { 349 ierr = PetscPrintf(PETSC_COMM_WORLD," Newton direction not descent for augmented Lagrangian: %g",(double)aldescent);CHKERRQ(ierr); 350 } 351 ierr = PetscInfo1(tao,"Newton direction not descent for augmented Lagrangian: %g\n",(double)aldescent);CHKERRQ(ierr); 352 lclP->rho = (rGL_U - adec)/rWU; 353 if (lclP->rho > lclP->rhomax) { 354 lclP->rho = lclP->rhomax; 355 SETERRQ1(PETSC_COMM_WORLD,0,"rho=%g > rhomax, case not implemented. Increase rhomax (-tao_lcl_rhomax)",(double)lclP->rho); 356 } 357 if (lclP->verbose) { 358 ierr = PetscPrintf(PETSC_COMM_WORLD," Increasing penalty parameter to %g\n",(double)lclP->rho);CHKERRQ(ierr); 359 } 360 ierr = PetscInfo1(tao," Increasing penalty parameter to %g\n",(double)lclP->rho);CHKERRQ(ierr); 361 } 362 363 ierr = LCLComputeAugmentedLagrangianAndGradient(tao->linesearch,tao->solution,&lclP->aug,lclP->GAugL,tao);CHKERRQ(ierr); 364 ierr = LCLScatter(lclP,lclP->GL,lclP->GL_U,lclP->GL_V);CHKERRQ(ierr); 365 ierr = LCLScatter(lclP,lclP->GAugL,lclP->GAugL_U,lclP->GAugL_V);CHKERRQ(ierr); 366 367 /* We now minimize the augmented Lagrangian along the Newton direction */ 368 ierr = VecScale(lclP->r,-1.0);CHKERRQ(ierr); 369 ierr = LCLGather(lclP, lclP->r,lclP->s,tao->stepdirection);CHKERRQ(ierr); 370 ierr = VecScale(lclP->r,-1.0);CHKERRQ(ierr); 371 ierr = LCLGather(lclP, lclP->GAugL_U0, lclP->GAugL_V0, lclP->GAugL);CHKERRQ(ierr); 372 ierr = LCLGather(lclP, lclP->U0,lclP->V0,lclP->X0);CHKERRQ(ierr); 373 374 lclP->recompute_jacobian_flag = PETSC_TRUE; 375 376 ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr); 377 ierr = TaoLineSearchSetType(tao->linesearch, TAOLINESEARCHMT);CHKERRQ(ierr); 378 ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 379 ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &lclP->aug, lclP->GAugL, tao->stepdirection, &step, &ls_reason);CHKERRQ(ierr); 380 if (lclP->verbose) { 381 ierr = PetscPrintf(PETSC_COMM_WORLD,"Steplength = %g\n",(double)step);CHKERRQ(ierr); 382 } 383 384 ierr = LCLScatter(lclP,tao->solution,lclP->U,lclP->V);CHKERRQ(ierr); 385 ierr = TaoComputeObjectiveAndGradient(tao,tao->solution,&f,tao->gradient);CHKERRQ(ierr); 386 ierr = LCLScatter(lclP,tao->gradient,lclP->GU,lclP->GV);CHKERRQ(ierr); 387 388 ierr = LCLScatter(lclP,lclP->GAugL,lclP->GAugL_U,lclP->GAugL_V);CHKERRQ(ierr); 389 390 /* Check convergence */ 391 ierr = VecNorm(lclP->GAugL, NORM_2, &mnorm);CHKERRQ(ierr); 392 ierr = VecNorm(tao->constraints, NORM_2, &cnorm);CHKERRQ(ierr); 393 ierr = TaoMonitor(tao,tao->niter,f,mnorm,cnorm,step,&reason);CHKERRQ(ierr); 394 if (reason != TAO_CONTINUE_ITERATING){ 395 break; 396 } 397 398 /* TODO: use a heuristic to choose how many iterations should be performed within phase 2 */ 399 for (phase2_iter=0; phase2_iter<lclP->phase2_niter; phase2_iter++){ 400 /* We now minimize the objective function starting from the fraction of 401 the Newton point accepted by applying one step of a reduced-space 402 method. The optimization problem is: 403 404 minimize f(u+du, v+dv) 405 s. t. A(u0,v0)du + B(u0,v0)du = -alpha g(u0,v0) 406 407 In particular, we have that 408 du = -inv(A)*(Bdv + alpha g) */ 409 410 ierr = TaoComputeJacobianState(tao,lclP->X0,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv);CHKERRQ(ierr); 411 ierr = TaoComputeJacobianDesign(tao,lclP->X0,tao->jacobian_design);CHKERRQ(ierr); 412 413 /* Store V and constraints */ 414 ierr = VecCopy(lclP->V, lclP->V1);CHKERRQ(ierr); 415 ierr = VecCopy(tao->constraints,lclP->con1);CHKERRQ(ierr); 416 417 /* Compute multipliers */ 418 /* p1 */ 419 ierr = VecSet(lclP->lamda,0.0);CHKERRQ(ierr); /* Initial guess in CG */ 420 lclP->solve_type = LCL_ADJOINT1; 421 ierr = MatIsSymmetricKnown(tao->jacobian_state,&set,&flag);CHKERRQ(ierr); 422 if (tao->jacobian_state_pre) { 423 ierr = MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag);CHKERRQ(ierr); 424 } else { 425 pset = pflag = PETSC_TRUE; 426 } 427 if (set && pset && flag && pflag) symmetric = PETSC_TRUE; 428 else symmetric = PETSC_FALSE; 429 430 if (tao->jacobian_state_inv) { 431 if (symmetric) { 432 ierr = MatMult(tao->jacobian_state_inv, lclP->GAugL_U, lclP->lamda);CHKERRQ(ierr); 433 } else { 434 ierr = MatMultTranspose(tao->jacobian_state_inv, lclP->GAugL_U, lclP->lamda);CHKERRQ(ierr); 435 } 436 } else { 437 if (symmetric) { 438 ierr = KSPSolve(tao->ksp, lclP->GAugL_U, lclP->lamda);CHKERRQ(ierr); 439 } else { 440 ierr = KSPSolveTranspose(tao->ksp, lclP->GAugL_U, lclP->lamda);CHKERRQ(ierr); 441 } 442 ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 443 tao->ksp_its+=its; 444 tao->ksp_tot_its+=its; 445 } 446 ierr = MatMultTranspose(tao->jacobian_design,lclP->lamda,lclP->g1);CHKERRQ(ierr); 447 ierr = VecAXPY(lclP->g1,-1.0,lclP->GAugL_V);CHKERRQ(ierr); 448 449 /* Compute the limited-memory quasi-newton direction */ 450 if (tao->niter > 0) { 451 ierr = MatLMVMSolve(lclP->R,lclP->g1,lclP->s);CHKERRQ(ierr); 452 ierr = VecDot(lclP->s,lclP->g1,&descent);CHKERRQ(ierr); 453 if (descent <= 0) { 454 if (lclP->verbose) { 455 ierr = PetscPrintf(PETSC_COMM_WORLD,"Reduced-space direction not descent: %g\n",(double)descent);CHKERRQ(ierr); 456 } 457 ierr = VecCopy(lclP->g1,lclP->s);CHKERRQ(ierr); 458 } 459 } else { 460 ierr = VecCopy(lclP->g1,lclP->s);CHKERRQ(ierr); 461 } 462 ierr = VecScale(lclP->g1,-1.0);CHKERRQ(ierr); 463 464 /* Recover the full space direction */ 465 ierr = MatMult(tao->jacobian_design,lclP->s,lclP->WU);CHKERRQ(ierr); 466 /* ierr = VecSet(lclP->r,0.0);CHKERRQ(ierr); */ /* Initial guess in CG */ 467 lclP->solve_type = LCL_FORWARD2; 468 if (tao->jacobian_state_inv) { 469 ierr = MatMult(tao->jacobian_state_inv,lclP->WU,lclP->r);CHKERRQ(ierr); 470 } else { 471 ierr = KSPSolve(tao->ksp, lclP->WU, lclP->r);CHKERRQ(ierr); 472 ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 473 tao->ksp_its += its; 474 tao->ksp_tot_its+=its; 475 } 476 477 /* We now minimize the augmented Lagrangian along the direction -r,s */ 478 ierr = VecScale(lclP->r, -1.0);CHKERRQ(ierr); 479 ierr = LCLGather(lclP,lclP->r,lclP->s,tao->stepdirection);CHKERRQ(ierr); 480 ierr = VecScale(lclP->r, -1.0);CHKERRQ(ierr); 481 lclP->recompute_jacobian_flag = PETSC_TRUE; 482 483 ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr); 484 ierr = TaoLineSearchSetType(tao->linesearch,TAOLINESEARCHMT);CHKERRQ(ierr); 485 ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 486 ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &lclP->aug, lclP->GAugL, tao->stepdirection,&step,&ls_reason);CHKERRQ(ierr); 487 if (lclP->verbose){ 488 ierr = PetscPrintf(PETSC_COMM_WORLD,"Reduced-space steplength = %g\n",(double)step);CHKERRQ(ierr); 489 } 490 491 ierr = LCLScatter(lclP,tao->solution,lclP->U,lclP->V);CHKERRQ(ierr); 492 ierr = LCLScatter(lclP,lclP->GL,lclP->GL_U,lclP->GL_V);CHKERRQ(ierr); 493 ierr = LCLScatter(lclP,lclP->GAugL,lclP->GAugL_U,lclP->GAugL_V);CHKERRQ(ierr); 494 ierr = TaoComputeObjectiveAndGradient(tao,tao->solution,&f,tao->gradient);CHKERRQ(ierr); 495 ierr = LCLScatter(lclP,tao->gradient,lclP->GU,lclP->GV);CHKERRQ(ierr); 496 497 /* Compute the reduced gradient at the new point */ 498 499 ierr = TaoComputeJacobianState(tao,lclP->X0,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv);CHKERRQ(ierr); 500 ierr = TaoComputeJacobianDesign(tao,lclP->X0,tao->jacobian_design);CHKERRQ(ierr); 501 502 /* p2 */ 503 /* Compute multipliers, using lamda-rho*con as an initial guess in PCG */ 504 if (phase2_iter==0){ 505 ierr = VecWAXPY(lclP->lamda,-lclP->rho,lclP->con1,lclP->lamda0);CHKERRQ(ierr); 506 } else { 507 ierr = VecAXPY(lclP->lamda,-lclP->rho,tao->constraints);CHKERRQ(ierr); 508 } 509 510 ierr = MatIsSymmetricKnown(tao->jacobian_state,&set,&flag);CHKERRQ(ierr); 511 if (tao->jacobian_state_pre) { 512 ierr = MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag);CHKERRQ(ierr); 513 } else { 514 pset = pflag = PETSC_TRUE; 515 } 516 if (set && pset && flag && pflag) symmetric = PETSC_TRUE; 517 else symmetric = PETSC_FALSE; 518 519 lclP->solve_type = LCL_ADJOINT2; 520 if (tao->jacobian_state_inv) { 521 if (symmetric) { 522 ierr = MatMult(tao->jacobian_state_inv, lclP->GU, lclP->lamda);CHKERRQ(ierr); 523 } else { 524 ierr = MatMultTranspose(tao->jacobian_state_inv, lclP->GU, lclP->lamda);CHKERRQ(ierr); 525 } 526 } else { 527 if (symmetric) { 528 ierr = KSPSolve(tao->ksp, lclP->GU, lclP->lamda);CHKERRQ(ierr); 529 } else { 530 ierr = KSPSolveTranspose(tao->ksp, lclP->GU, lclP->lamda);CHKERRQ(ierr); 531 } 532 ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 533 tao->ksp_its += its; 534 tao->ksp_tot_its += its; 535 } 536 537 ierr = MatMultTranspose(tao->jacobian_design,lclP->lamda,lclP->g2);CHKERRQ(ierr); 538 ierr = VecAXPY(lclP->g2,-1.0,lclP->GV);CHKERRQ(ierr); 539 540 ierr = VecScale(lclP->g2,-1.0);CHKERRQ(ierr); 541 542 /* Update the quasi-newton approximation */ 543 if (phase2_iter >= 0){ 544 ierr = MatLMVMSetPrev(lclP->R,lclP->V1,lclP->g1);CHKERRQ(ierr); 545 } 546 ierr = MatLMVMUpdate(lclP->R,lclP->V,lclP->g2);CHKERRQ(ierr); 547 /* 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 */ 548 549 } 550 551 ierr = VecCopy(lclP->lamda,lclP->lamda0);CHKERRQ(ierr); 552 553 /* Evaluate Function, Gradient, Constraints, and Jacobian */ 554 ierr = TaoComputeObjectiveAndGradient(tao,tao->solution,&f,tao->gradient);CHKERRQ(ierr); 555 ierr = LCLScatter(lclP,tao->solution,lclP->U,lclP->V);CHKERRQ(ierr); 556 ierr = LCLScatter(lclP,tao->gradient,lclP->GU,lclP->GV);CHKERRQ(ierr); 557 558 ierr = TaoComputeJacobianState(tao,tao->solution,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv);CHKERRQ(ierr); 559 ierr = TaoComputeJacobianDesign(tao,tao->solution,tao->jacobian_design);CHKERRQ(ierr); 560 ierr = TaoComputeConstraints(tao,tao->solution, tao->constraints);CHKERRQ(ierr); 561 562 ierr = LCLComputeAugmentedLagrangianAndGradient(tao->linesearch,tao->solution,&lclP->aug,lclP->GAugL,tao);CHKERRQ(ierr); 563 564 ierr = VecNorm(lclP->GAugL, NORM_2, &mnorm);CHKERRQ(ierr); 565 566 /* Evaluate constraint norm */ 567 ierr = VecNorm(tao->constraints, NORM_2, &cnorm);CHKERRQ(ierr); 568 569 /* Monitor convergence */ 570 tao->niter++; 571 ierr = TaoMonitor(tao, tao->niter,f,mnorm,cnorm,step,&reason);CHKERRQ(ierr); 572 } 573 ierr = MatDestroy(&lclP->R);CHKERRQ(ierr); 574 PetscFunctionReturn(0); 575 } 576 577 /*MC 578 TAOLCL - linearly constrained lagrangian method for pde-constrained optimization 579 580 + -tao_lcl_eps1 - epsilon 1 tolerance 581 . -tao_lcl_eps2","epsilon 2 tolerance","",lclP->eps2,&lclP->eps2,NULL);CHKERRQ(ierr); 582 . -tao_lcl_rho0","init value for rho","",lclP->rho0,&lclP->rho0,NULL);CHKERRQ(ierr); 583 . -tao_lcl_rhomax","max value for rho","",lclP->rhomax,&lclP->rhomax,NULL);CHKERRQ(ierr); 584 . -tao_lcl_phase2_niter - Number of phase 2 iterations in LCL algorithm 585 . -tao_lcl_verbose - Print verbose output if True 586 . -tao_lcl_tola - Tolerance for first forward solve 587 . -tao_lcl_tolb - Tolerance for first adjoint solve 588 . -tao_lcl_tolc - Tolerance for second forward solve 589 - -tao_lcl_told - Tolerance for second adjoint solve 590 591 Level: beginner 592 M*/ 593 #undef __FUNCT__ 594 #define __FUNCT__ "TaoCreate_LCL" 595 PETSC_EXTERN PetscErrorCode TaoCreate_LCL(Tao tao) 596 { 597 TAO_LCL *lclP; 598 PetscErrorCode ierr; 599 const char *morethuente_type = TAOLINESEARCHMT; 600 601 PetscFunctionBegin; 602 tao->ops->setup = TaoSetup_LCL; 603 tao->ops->solve = TaoSolve_LCL; 604 tao->ops->view = TaoView_LCL; 605 tao->ops->setfromoptions = TaoSetFromOptions_LCL; 606 tao->ops->destroy = TaoDestroy_LCL; 607 ierr = PetscNewLog(tao,&lclP);CHKERRQ(ierr); 608 tao->data = (void*)lclP; 609 610 /* Override default settings (unless already changed) */ 611 if (!tao->max_it_changed) tao->max_it = 200; 612 if (!tao->catol_changed) tao->catol = 1.0e-4; 613 if (!tao->gatol_changed) tao->gttol = 1.0e-4; 614 if (!tao->grtol_changed) tao->gttol = 1.0e-4; 615 if (!tao->gttol_changed) tao->gttol = 1.0e-4; 616 lclP->rho0 = 1.0e-4; 617 lclP->rhomax=1e5; 618 lclP->eps1 = 1.0e-8; 619 lclP->eps2 = 0.0; 620 lclP->solve_type=2; 621 lclP->tau[0] = lclP->tau[1] = lclP->tau[2] = lclP->tau[3] = 1.0e-4; 622 ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr); 623 ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr); 624 ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 625 626 ierr = TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,LCLComputeAugmentedLagrangianAndGradient, tao);CHKERRQ(ierr); 627 ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr); 628 ierr = KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix);CHKERRQ(ierr); 629 ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 630 PetscFunctionReturn(0); 631 } 632 633 #undef __FUNCT__ 634 #define __FUNCT__ "LCLComputeLagrangianAndGradient" 635 static PetscErrorCode LCLComputeLagrangianAndGradient(TaoLineSearch ls, Vec X, PetscReal *f, Vec G, void *ptr) 636 { 637 Tao tao = (Tao)ptr; 638 TAO_LCL *lclP = (TAO_LCL*)tao->data; 639 PetscBool set,pset,flag,pflag,symmetric; 640 PetscReal cdotl; 641 PetscErrorCode ierr; 642 643 PetscFunctionBegin; 644 ierr = TaoComputeObjectiveAndGradient(tao,X,f,G);CHKERRQ(ierr); 645 ierr = LCLScatter(lclP,G,lclP->GU,lclP->GV);CHKERRQ(ierr); 646 if (lclP->recompute_jacobian_flag) { 647 ierr = TaoComputeJacobianState(tao,X,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv);CHKERRQ(ierr); 648 ierr = TaoComputeJacobianDesign(tao,X,tao->jacobian_design);CHKERRQ(ierr); 649 } 650 ierr = TaoComputeConstraints(tao,X, tao->constraints);CHKERRQ(ierr); 651 ierr = MatIsSymmetricKnown(tao->jacobian_state,&set,&flag);CHKERRQ(ierr); 652 if (tao->jacobian_state_pre) { 653 ierr = MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag);CHKERRQ(ierr); 654 } else { 655 pset = pflag = PETSC_TRUE; 656 } 657 if (set && pset && flag && pflag) symmetric = PETSC_TRUE; 658 else symmetric = PETSC_FALSE; 659 660 ierr = VecDot(lclP->lamda0, tao->constraints, &cdotl);CHKERRQ(ierr); 661 lclP->lgn = *f - cdotl; 662 663 /* Gradient of Lagrangian GL = G - J' * lamda */ 664 /* WU = A' * WL 665 WV = B' * WL */ 666 if (symmetric) { 667 ierr = MatMult(tao->jacobian_state,lclP->lamda0,lclP->GL_U);CHKERRQ(ierr); 668 } else { 669 ierr = MatMultTranspose(tao->jacobian_state,lclP->lamda0,lclP->GL_U);CHKERRQ(ierr); 670 } 671 ierr = MatMultTranspose(tao->jacobian_design,lclP->lamda0,lclP->GL_V);CHKERRQ(ierr); 672 ierr = VecScale(lclP->GL_U,-1.0);CHKERRQ(ierr); 673 ierr = VecScale(lclP->GL_V,-1.0);CHKERRQ(ierr); 674 ierr = VecAXPY(lclP->GL_U,1.0,lclP->GU);CHKERRQ(ierr); 675 ierr = VecAXPY(lclP->GL_V,1.0,lclP->GV);CHKERRQ(ierr); 676 ierr = LCLGather(lclP,lclP->GL_U,lclP->GL_V,lclP->GL);CHKERRQ(ierr); 677 678 f[0] = lclP->lgn; 679 ierr = VecCopy(lclP->GL,G);CHKERRQ(ierr); 680 PetscFunctionReturn(0); 681 } 682 683 #undef __FUNCT__ 684 #define __FUNCT__ "LCLComputeAugmentedLagrangianAndGradient" 685 static PetscErrorCode LCLComputeAugmentedLagrangianAndGradient(TaoLineSearch ls, Vec X, PetscReal *f, Vec G, void *ptr) 686 { 687 Tao tao = (Tao)ptr; 688 TAO_LCL *lclP = (TAO_LCL*)tao->data; 689 PetscReal con2; 690 PetscBool flag,pflag,set,pset,symmetric; 691 PetscErrorCode ierr; 692 693 PetscFunctionBegin; 694 ierr = LCLComputeLagrangianAndGradient(tao->linesearch,X,f,G,tao);CHKERRQ(ierr); 695 ierr = LCLScatter(lclP,G,lclP->GL_U,lclP->GL_V);CHKERRQ(ierr); 696 ierr = VecDot(tao->constraints,tao->constraints,&con2);CHKERRQ(ierr); 697 lclP->aug = lclP->lgn + 0.5*lclP->rho*con2; 698 699 /* Gradient of Aug. Lagrangian GAugL = GL + rho * J' c */ 700 /* WU = A' * c 701 WV = B' * c */ 702 ierr = MatIsSymmetricKnown(tao->jacobian_state,&set,&flag);CHKERRQ(ierr); 703 if (tao->jacobian_state_pre) { 704 ierr = MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag);CHKERRQ(ierr); 705 } else { 706 pset = pflag = PETSC_TRUE; 707 } 708 if (set && pset && flag && pflag) symmetric = PETSC_TRUE; 709 else symmetric = PETSC_FALSE; 710 711 if (symmetric) { 712 ierr = MatMult(tao->jacobian_state,tao->constraints,lclP->GAugL_U);CHKERRQ(ierr); 713 } else { 714 ierr = MatMultTranspose(tao->jacobian_state,tao->constraints,lclP->GAugL_U);CHKERRQ(ierr); 715 } 716 717 ierr = MatMultTranspose(tao->jacobian_design,tao->constraints,lclP->GAugL_V);CHKERRQ(ierr); 718 ierr = VecAYPX(lclP->GAugL_U,lclP->rho,lclP->GL_U);CHKERRQ(ierr); 719 ierr = VecAYPX(lclP->GAugL_V,lclP->rho,lclP->GL_V);CHKERRQ(ierr); 720 ierr = LCLGather(lclP,lclP->GAugL_U,lclP->GAugL_V,lclP->GAugL);CHKERRQ(ierr); 721 722 f[0] = lclP->aug; 723 ierr = VecCopy(lclP->GAugL,G);CHKERRQ(ierr); 724 PetscFunctionReturn(0); 725 } 726 727 #undef __FUNCT__ 728 #define __FUNCT__ "LCLGather" 729 PetscErrorCode LCLGather(TAO_LCL *lclP, Vec u, Vec v, Vec x) 730 { 731 PetscErrorCode ierr; 732 PetscFunctionBegin; 733 ierr = VecScatterBegin(lclP->state_scatter, u, x, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 734 ierr = VecScatterEnd(lclP->state_scatter, u, x, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 735 ierr = VecScatterBegin(lclP->design_scatter, v, x, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 736 ierr = VecScatterEnd(lclP->design_scatter, v, x, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 737 PetscFunctionReturn(0); 738 739 } 740 #undef __FUNCT__ 741 #define __FUNCT__ "LCLScatter" 742 PetscErrorCode LCLScatter(TAO_LCL *lclP, Vec x, Vec u, Vec v) 743 { 744 PetscErrorCode ierr; 745 PetscFunctionBegin; 746 ierr = VecScatterBegin(lclP->state_scatter, x, u, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 747 ierr = VecScatterEnd(lclP->state_scatter, x, u, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 748 ierr = VecScatterBegin(lclP->design_scatter, x, v, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 749 ierr = VecScatterEnd(lclP->design_scatter, x, v, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 750 PetscFunctionReturn(0); 751 752 } 753