1*a7e14dcfSSatish Balay #include "ssls.h" 2*a7e14dcfSSatish Balay 3*a7e14dcfSSatish Balay #undef __FUNCT__ 4*a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetUp_SSILS" 5*a7e14dcfSSatish Balay PetscErrorCode TaoSetUp_SSILS(TaoSolver tao) 6*a7e14dcfSSatish Balay { 7*a7e14dcfSSatish Balay TAO_SSLS *ssls = (TAO_SSLS *)tao->data; 8*a7e14dcfSSatish Balay PetscErrorCode ierr; 9*a7e14dcfSSatish Balay 10*a7e14dcfSSatish Balay PetscFunctionBegin; 11*a7e14dcfSSatish Balay 12*a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&tao->gradient); CHKERRQ(ierr); 13*a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&tao->stepdirection); CHKERRQ(ierr); 14*a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&ssls->ff); CHKERRQ(ierr); 15*a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&ssls->dpsi); CHKERRQ(ierr); 16*a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&ssls->da); CHKERRQ(ierr); 17*a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&ssls->db); CHKERRQ(ierr); 18*a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&ssls->t1); CHKERRQ(ierr); 19*a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&ssls->t2); CHKERRQ(ierr); 20*a7e14dcfSSatish Balay PetscFunctionReturn(0); 21*a7e14dcfSSatish Balay } 22*a7e14dcfSSatish Balay 23*a7e14dcfSSatish Balay 24*a7e14dcfSSatish Balay #undef __FUNCT__ 25*a7e14dcfSSatish Balay #define __FUNCT__ "TaoDestroy_SSILS" 26*a7e14dcfSSatish Balay PetscErrorCode TaoDestroy_SSILS(TaoSolver tao) 27*a7e14dcfSSatish Balay { 28*a7e14dcfSSatish Balay TAO_SSLS *ssls = (TAO_SSLS *)tao->data; 29*a7e14dcfSSatish Balay PetscErrorCode ierr; 30*a7e14dcfSSatish Balay 31*a7e14dcfSSatish Balay PetscFunctionBegin; 32*a7e14dcfSSatish Balay 33*a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->ff); CHKERRQ(ierr); 34*a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->dpsi); CHKERRQ(ierr); 35*a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->da); CHKERRQ(ierr); 36*a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->db); CHKERRQ(ierr); 37*a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->t1); CHKERRQ(ierr); 38*a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->t2); CHKERRQ(ierr); 39*a7e14dcfSSatish Balay ierr = PetscFree(tao->data); CHKERRQ(ierr); 40*a7e14dcfSSatish Balay tao->data = PETSC_NULL; 41*a7e14dcfSSatish Balay PetscFunctionReturn(0); 42*a7e14dcfSSatish Balay } 43*a7e14dcfSSatish Balay 44*a7e14dcfSSatish Balay #undef __FUNCT__ 45*a7e14dcfSSatish Balay #define __FUNCT__ "TaoSolve_SSILS" 46*a7e14dcfSSatish Balay static PetscErrorCode TaoSolve_SSILS(TaoSolver tao) 47*a7e14dcfSSatish Balay { 48*a7e14dcfSSatish Balay TAO_SSLS *ssls = (TAO_SSLS *)tao->data; 49*a7e14dcfSSatish Balay PetscReal psi, ndpsi, normd, innerd, t=0; 50*a7e14dcfSSatish Balay PetscReal delta, rho; 51*a7e14dcfSSatish Balay PetscInt iter=0,kspits; 52*a7e14dcfSSatish Balay TaoSolverTerminationReason reason; 53*a7e14dcfSSatish Balay TaoLineSearchTerminationReason ls_reason; 54*a7e14dcfSSatish Balay PetscErrorCode ierr; 55*a7e14dcfSSatish Balay 56*a7e14dcfSSatish Balay PetscFunctionBegin; 57*a7e14dcfSSatish Balay 58*a7e14dcfSSatish Balay /* Assume that Setup has been called! 59*a7e14dcfSSatish Balay Set the structure for the Jacobian and create a linear solver. */ 60*a7e14dcfSSatish Balay 61*a7e14dcfSSatish Balay delta = ssls->delta; 62*a7e14dcfSSatish Balay rho = ssls->rho; 63*a7e14dcfSSatish Balay 64*a7e14dcfSSatish Balay ierr = TaoComputeVariableBounds(tao); CHKERRQ(ierr); 65*a7e14dcfSSatish Balay ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution); CHKERRQ(ierr); 66*a7e14dcfSSatish Balay ierr = TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,Tao_SSLS_FunctionGradient,tao); 67*a7e14dcfSSatish Balay CHKERRQ(ierr); 68*a7e14dcfSSatish Balay ierr = TaoLineSearchSetObjectiveRoutine(tao->linesearch,Tao_SSLS_Function,tao); CHKERRQ(ierr); 69*a7e14dcfSSatish Balay 70*a7e14dcfSSatish Balay /* Calculate the function value and fischer function value at the 71*a7e14dcfSSatish Balay current iterate */ 72*a7e14dcfSSatish Balay ierr = TaoLineSearchComputeObjectiveAndGradient(tao->linesearch,tao->solution,&psi,ssls->dpsi); CHKERRQ(ierr); 73*a7e14dcfSSatish Balay ierr = VecNorm(ssls->dpsi,NORM_2,&ndpsi); CHKERRQ(ierr); 74*a7e14dcfSSatish Balay 75*a7e14dcfSSatish Balay 76*a7e14dcfSSatish Balay while (1) { 77*a7e14dcfSSatish Balay ierr=PetscInfo3(tao, "iter: %D, merit: %G, ndpsi: %G\n", 78*a7e14dcfSSatish Balay iter, ssls->merit, ndpsi); CHKERRQ(ierr); 79*a7e14dcfSSatish Balay /* Check the termination criteria */ 80*a7e14dcfSSatish Balay ierr = TaoMonitor(tao,iter++,ssls->merit,ndpsi,0.0,t,&reason); 81*a7e14dcfSSatish Balay CHKERRQ(ierr); 82*a7e14dcfSSatish Balay if (reason!=TAO_CONTINUE_ITERATING) break; 83*a7e14dcfSSatish Balay 84*a7e14dcfSSatish Balay /* Calculate direction. (Really negative of newton direction. Therefore, 85*a7e14dcfSSatish Balay rest of the code uses -d.) */ 86*a7e14dcfSSatish Balay ierr = KSPSetOperators(tao->ksp,tao->jacobian,tao->jacobian_pre,ssls->matflag); CHKERRQ(ierr); 87*a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp,ssls->ff,tao->stepdirection); CHKERRQ(ierr); 88*a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&kspits); CHKERRQ(ierr); 89*a7e14dcfSSatish Balay tao->ksp_its+=kspits; 90*a7e14dcfSSatish Balay ierr = VecNorm(tao->stepdirection,NORM_2,&normd); CHKERRQ(ierr); 91*a7e14dcfSSatish Balay ierr = VecDot(tao->stepdirection,ssls->dpsi,&innerd); CHKERRQ(ierr); 92*a7e14dcfSSatish Balay 93*a7e14dcfSSatish Balay /* Make sure that we have a descent direction */ 94*a7e14dcfSSatish Balay if (innerd <= delta*pow(normd, rho)) { 95*a7e14dcfSSatish Balay ierr = PetscInfo(tao, "newton direction not descent\n"); CHKERRQ(ierr); 96*a7e14dcfSSatish Balay ierr = VecCopy(ssls->dpsi,tao->stepdirection); CHKERRQ(ierr); 97*a7e14dcfSSatish Balay ierr = VecDot(tao->stepdirection,ssls->dpsi,&innerd); CHKERRQ(ierr); 98*a7e14dcfSSatish Balay } 99*a7e14dcfSSatish Balay 100*a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection, -1.0); CHKERRQ(ierr); 101*a7e14dcfSSatish Balay innerd = -innerd; 102*a7e14dcfSSatish Balay 103*a7e14dcfSSatish Balay ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0); 104*a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch,tao->solution,&psi,ssls->dpsi,tao->stepdirection,&t,&ls_reason); CHKERRQ(ierr); 105*a7e14dcfSSatish Balay ierr = VecNorm(ssls->dpsi,NORM_2,&ndpsi); CHKERRQ(ierr); 106*a7e14dcfSSatish Balay } 107*a7e14dcfSSatish Balay 108*a7e14dcfSSatish Balay PetscFunctionReturn(0); 109*a7e14dcfSSatish Balay } 110*a7e14dcfSSatish Balay 111*a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 112*a7e14dcfSSatish Balay EXTERN_C_BEGIN 113*a7e14dcfSSatish Balay #undef __FUNCT__ 114*a7e14dcfSSatish Balay #define __FUNCT__ "TaoCreate_SSILS" 115*a7e14dcfSSatish Balay PetscErrorCode TaoCreate_SSILS(TaoSolver tao) 116*a7e14dcfSSatish Balay { 117*a7e14dcfSSatish Balay TAO_SSLS *ssls; 118*a7e14dcfSSatish Balay PetscErrorCode ierr; 119*a7e14dcfSSatish Balay const char *armijo_type = TAOLINESEARCH_ARMIJO; 120*a7e14dcfSSatish Balay 121*a7e14dcfSSatish Balay PetscFunctionBegin; 122*a7e14dcfSSatish Balay ierr = PetscNewLog(tao,TAO_SSLS,&ssls); CHKERRQ(ierr); 123*a7e14dcfSSatish Balay tao->data = (void*)ssls; 124*a7e14dcfSSatish Balay tao->ops->solve=TaoSolve_SSILS; 125*a7e14dcfSSatish Balay tao->ops->setup=TaoSetUp_SSILS; 126*a7e14dcfSSatish Balay tao->ops->view=TaoView_SSLS; 127*a7e14dcfSSatish Balay tao->ops->setfromoptions=TaoSetFromOptions_SSLS; 128*a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_SSILS; 129*a7e14dcfSSatish Balay 130*a7e14dcfSSatish Balay ssls->delta = 1e-10; 131*a7e14dcfSSatish Balay ssls->rho = 2.1; 132*a7e14dcfSSatish Balay 133*a7e14dcfSSatish Balay ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch); CHKERRQ(ierr); 134*a7e14dcfSSatish Balay ierr = TaoLineSearchSetType(tao->linesearch,armijo_type); CHKERRQ(ierr); 135*a7e14dcfSSatish Balay ierr = TaoLineSearchSetFromOptions(tao->linesearch); 136*a7e14dcfSSatish Balay /* Note: linesearch objective and objectivegradient routines are set in solve routine */ 137*a7e14dcfSSatish Balay ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp); CHKERRQ(ierr); 138*a7e14dcfSSatish Balay 139*a7e14dcfSSatish Balay tao->max_it = 2000; 140*a7e14dcfSSatish Balay tao->max_funcs = 4000; 141*a7e14dcfSSatish Balay tao->fatol = 0; tao->frtol = 0; tao->gttol=0; tao->grtol=0; 142*a7e14dcfSSatish Balay tao->gatol = 1.0e-16; 143*a7e14dcfSSatish Balay tao->fmin = 1.0e-8; 144*a7e14dcfSSatish Balay 145*a7e14dcfSSatish Balay PetscFunctionReturn(0); 146*a7e14dcfSSatish Balay } 147*a7e14dcfSSatish Balay EXTERN_C_END 148*a7e14dcfSSatish Balay 149