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