1c4762a1bSJed Brown /* Program usage: mpiexec -n 1 maros1 [-help] [all TAO options] */ 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* ---------------------------------------------------------------------- 4c4762a1bSJed Brown TODO Explain maros example 5c4762a1bSJed Brown ---------------------------------------------------------------------- */ 6c4762a1bSJed Brown 7c4762a1bSJed Brown #include <petsctao.h> 8c4762a1bSJed Brown 9c4762a1bSJed Brown static char help[]=""; 10c4762a1bSJed Brown 11c4762a1bSJed Brown /*T 12c4762a1bSJed Brown Concepts: TAO^Solving an unconstrained minimization problem 13c4762a1bSJed Brown Routines: TaoCreate(); TaoSetType(); 14a82e8c82SStefano Zampini Routines: TaoSetSolution(); 15a82e8c82SStefano Zampini Routines: TaoSetObjectiveAndGradient(); 16c4762a1bSJed Brown Routines: TaoSetEqualityConstraintsRoutine(); 17c4762a1bSJed Brown Routines: TaoSetInequalityConstraintsRoutine(); 18c4762a1bSJed Brown Routines: TaoSetEqualityJacobianRoutine(); 19c4762a1bSJed Brown Routines: TaoSetInequalityJacobianRoutine(); 20a82e8c82SStefano Zampini Routines: TaoSetHessian(); TaoSetFromOptions(); 21c4762a1bSJed Brown Routines: TaoGetKSP(); TaoSolve(); 22c4762a1bSJed Brown Routines: TaoGetConvergedReason();TaoDestroy(); 23c4762a1bSJed Brown Processors: 1 24c4762a1bSJed Brown T*/ 25c4762a1bSJed Brown 26c4762a1bSJed Brown /* 27c4762a1bSJed Brown User-defined application context - contains data needed by the 28c4762a1bSJed Brown application-provided call-back routines, FormFunction(), 29c4762a1bSJed Brown FormGradient(), and FormHessian(). 30c4762a1bSJed Brown */ 31c4762a1bSJed Brown 32c4762a1bSJed Brown /* 33c4762a1bSJed Brown x,d in R^n 34c4762a1bSJed Brown f in R 35c4762a1bSJed Brown bin in R^mi 36c4762a1bSJed Brown beq in R^me 37c4762a1bSJed Brown Aeq in R^(me x n) 38c4762a1bSJed Brown Ain in R^(mi x n) 39c4762a1bSJed Brown H in R^(n x n) 40c4762a1bSJed Brown min f=(1/2)*x'*H*x + d'*x 41c4762a1bSJed Brown s.t. Aeq*x == beq 42c4762a1bSJed Brown Ain*x >= bin 43c4762a1bSJed Brown */ 44c4762a1bSJed Brown typedef struct { 45c4762a1bSJed Brown char name[32]; 46c4762a1bSJed Brown PetscInt n; /* Length x */ 47c4762a1bSJed Brown PetscInt me; /* number of equality constraints */ 48c4762a1bSJed Brown PetscInt mi; /* number of inequality constraints */ 49c4762a1bSJed Brown PetscInt m; /* me+mi */ 50c4762a1bSJed Brown Mat Aeq,Ain,H; 51c4762a1bSJed Brown Vec beq,bin,d; 52c4762a1bSJed Brown } AppCtx; 53c4762a1bSJed Brown 54c4762a1bSJed Brown /* -------- User-defined Routines --------- */ 55c4762a1bSJed Brown 56c4762a1bSJed Brown PetscErrorCode InitializeProblem(AppCtx*); 57c4762a1bSJed Brown PetscErrorCode DestroyProblem(AppCtx *); 58c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal *,Vec,void *); 59c4762a1bSJed Brown PetscErrorCode FormHessian(Tao,Vec,Mat,Mat, void*); 60c4762a1bSJed Brown PetscErrorCode FormInequalityConstraints(Tao,Vec,Vec,void*); 61c4762a1bSJed Brown PetscErrorCode FormEqualityConstraints(Tao,Vec,Vec,void*); 62c4762a1bSJed Brown PetscErrorCode FormInequalityJacobian(Tao,Vec,Mat,Mat, void*); 63c4762a1bSJed Brown PetscErrorCode FormEqualityJacobian(Tao,Vec,Mat,Mat, void*); 64c4762a1bSJed Brown 65c4762a1bSJed Brown PetscErrorCode main(int argc,char **argv) 66c4762a1bSJed Brown { 67c4762a1bSJed Brown PetscErrorCode ierr; /* used to check for functions returning nonzeros */ 68c4762a1bSJed Brown PetscMPIInt size; 69c4762a1bSJed Brown Vec x; /* solution */ 70c4762a1bSJed Brown KSP ksp; 71c4762a1bSJed Brown PC pc; 72c4762a1bSJed Brown Vec ceq,cin; 73c4762a1bSJed Brown PetscBool flg; /* A return value when checking for use options */ 74c4762a1bSJed Brown Tao tao; /* Tao solver context */ 75c4762a1bSJed Brown TaoConvergedReason reason; 76c4762a1bSJed Brown AppCtx user; /* application context */ 77c4762a1bSJed Brown 78c4762a1bSJed Brown /* Initialize TAO,PETSc */ 79c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char *)0,help);if (ierr) return ierr; 80ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 81c4762a1bSJed Brown /* Specify default parameters for the problem, check for command-line overrides */ 82589a23caSBarry Smith ierr = PetscStrncpy(user.name,"HS21",sizeof(user.name));CHKERRQ(ierr); 83589a23caSBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-cutername",user.name,sizeof(user.name),&flg);CHKERRQ(ierr); 84c4762a1bSJed Brown 85c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\n---- MAROS Problem %s -----\n",user.name);CHKERRQ(ierr); 86c4762a1bSJed Brown ierr = InitializeProblem(&user);CHKERRQ(ierr); 87c4762a1bSJed Brown ierr = VecDuplicate(user.d,&x);CHKERRQ(ierr); 88c4762a1bSJed Brown ierr = VecDuplicate(user.beq,&ceq);CHKERRQ(ierr); 89c4762a1bSJed Brown ierr = VecDuplicate(user.bin,&cin);CHKERRQ(ierr); 90c4762a1bSJed Brown ierr = VecSet(x,1.0);CHKERRQ(ierr); 91c4762a1bSJed Brown 92c4762a1bSJed Brown ierr = TaoCreate(PETSC_COMM_WORLD,&tao);CHKERRQ(ierr); 93c4762a1bSJed Brown ierr = TaoSetType(tao,TAOIPM);CHKERRQ(ierr); 94a82e8c82SStefano Zampini ierr = TaoSetSolution(tao,x);CHKERRQ(ierr); 95a82e8c82SStefano Zampini ierr = TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void*)&user);CHKERRQ(ierr); 96c4762a1bSJed Brown ierr = TaoSetEqualityConstraintsRoutine(tao,ceq,FormEqualityConstraints,(void*)&user);CHKERRQ(ierr); 97c4762a1bSJed Brown ierr = TaoSetInequalityConstraintsRoutine(tao,cin,FormInequalityConstraints,(void*)&user);CHKERRQ(ierr); 98c4762a1bSJed Brown ierr = TaoSetInequalityBounds(tao,user.bin,NULL);CHKERRQ(ierr); 99c4762a1bSJed Brown ierr = TaoSetJacobianEqualityRoutine(tao,user.Aeq,user.Aeq,FormEqualityJacobian,(void*)&user);CHKERRQ(ierr); 100c4762a1bSJed Brown ierr = TaoSetJacobianInequalityRoutine(tao,user.Ain,user.Ain,FormInequalityJacobian,(void*)&user);CHKERRQ(ierr); 101a82e8c82SStefano Zampini ierr = TaoSetHessian(tao,user.H,user.H,FormHessian,(void*)&user);CHKERRQ(ierr); 102c4762a1bSJed Brown ierr = TaoGetKSP(tao,&ksp);CHKERRQ(ierr); 103c4762a1bSJed Brown ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 104c4762a1bSJed Brown ierr = PCSetType(pc,PCLU);CHKERRQ(ierr); 105c4762a1bSJed Brown /* 106c4762a1bSJed Brown This algorithm produces matrices with zeros along the diagonal therefore we need to use 107c4762a1bSJed Brown SuperLU which does partial pivoting 108c4762a1bSJed Brown */ 109c4762a1bSJed Brown ierr = PCFactorSetMatSolverType(pc,MATSOLVERSUPERLU);CHKERRQ(ierr); 110c4762a1bSJed Brown ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 111c4762a1bSJed Brown ierr = TaoSetTolerances(tao,0,0,0);CHKERRQ(ierr); 112c4762a1bSJed Brown 113c4762a1bSJed Brown ierr = TaoSetFromOptions(tao);CHKERRQ(ierr); 114c4762a1bSJed Brown ierr = TaoSolve(tao);CHKERRQ(ierr); 115c4762a1bSJed Brown ierr = TaoGetConvergedReason(tao,&reason);CHKERRQ(ierr); 116c4762a1bSJed Brown if (reason < 0) { 117c4762a1bSJed Brown ierr = PetscPrintf(MPI_COMM_WORLD, "TAO failed to converge due to %s.\n",TaoConvergedReasons[reason]);CHKERRQ(ierr); 118c4762a1bSJed Brown } else { 119c4762a1bSJed Brown ierr = PetscPrintf(MPI_COMM_WORLD, "Optimization completed with status %s.\n",TaoConvergedReasons[reason]);CHKERRQ(ierr); 120c4762a1bSJed Brown } 121c4762a1bSJed Brown 122c4762a1bSJed Brown ierr = DestroyProblem(&user);CHKERRQ(ierr); 123c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 124c4762a1bSJed Brown ierr = VecDestroy(&ceq);CHKERRQ(ierr); 125c4762a1bSJed Brown ierr = VecDestroy(&cin);CHKERRQ(ierr); 126c4762a1bSJed Brown ierr = TaoDestroy(&tao);CHKERRQ(ierr); 127c4762a1bSJed Brown 128c4762a1bSJed Brown ierr = PetscFinalize(); 129c4762a1bSJed Brown return ierr; 130c4762a1bSJed Brown } 131c4762a1bSJed Brown 132c4762a1bSJed Brown PetscErrorCode InitializeProblem(AppCtx *user) 133c4762a1bSJed Brown { 134c4762a1bSJed Brown PetscErrorCode ierr; 135c4762a1bSJed Brown PetscViewer loader; 136c4762a1bSJed Brown MPI_Comm comm; 137c4762a1bSJed Brown PetscInt nrows,ncols,i; 138c4762a1bSJed Brown PetscScalar one=1.0; 139c4762a1bSJed Brown char filebase[128]; 140c4762a1bSJed Brown char filename[128]; 141c4762a1bSJed Brown 142c4762a1bSJed Brown PetscFunctionBegin; 143c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 144c4762a1bSJed Brown ierr = PetscStrncpy(filebase,user->name,sizeof(filebase));CHKERRQ(ierr); 145c4762a1bSJed Brown ierr = PetscStrlcat(filebase,"/",sizeof(filebase));CHKERRQ(ierr); 146c4762a1bSJed Brown ierr = PetscStrncpy(filename,filebase,sizeof(filename));CHKERRQ(ierr); 147c4762a1bSJed Brown ierr = PetscStrlcat(filename,"f",sizeof(filename));CHKERRQ(ierr); 148c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(comm,filename,FILE_MODE_READ,&loader);CHKERRQ(ierr); 149c4762a1bSJed Brown 150c4762a1bSJed Brown ierr = VecCreate(comm,&user->d);CHKERRQ(ierr); 151c4762a1bSJed Brown ierr = VecLoad(user->d,loader);CHKERRQ(ierr); 152c4762a1bSJed Brown ierr = PetscViewerDestroy(&loader);CHKERRQ(ierr); 153c4762a1bSJed Brown ierr = VecGetSize(user->d,&nrows);CHKERRQ(ierr); 154c4762a1bSJed Brown ierr = VecSetFromOptions(user->d);CHKERRQ(ierr); 155c4762a1bSJed Brown user->n = nrows; 156c4762a1bSJed Brown 157c4762a1bSJed Brown ierr = PetscStrncpy(filename,filebase,sizeof(filename));CHKERRQ(ierr); 158c4762a1bSJed Brown ierr = PetscStrlcat(filename,"H",sizeof(filename));CHKERRQ(ierr); 159c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(comm,filename,FILE_MODE_READ,&loader);CHKERRQ(ierr); 160c4762a1bSJed Brown 161c4762a1bSJed Brown ierr = MatCreate(comm,&user->H);CHKERRQ(ierr); 162c4762a1bSJed Brown ierr = MatSetSizes(user->H,PETSC_DECIDE,PETSC_DECIDE,nrows,nrows);CHKERRQ(ierr); 163c4762a1bSJed Brown ierr = MatLoad(user->H,loader);CHKERRQ(ierr); 164c4762a1bSJed Brown ierr = PetscViewerDestroy(&loader);CHKERRQ(ierr); 165c4762a1bSJed Brown ierr = MatGetSize(user->H,&nrows,&ncols);CHKERRQ(ierr); 166*3c859ba3SBarry Smith PetscCheck(nrows == user->n,comm,PETSC_ERR_ARG_SIZ,"H: nrows != n"); 167*3c859ba3SBarry Smith PetscCheck(ncols == user->n,comm,PETSC_ERR_ARG_SIZ,"H: ncols != n"); 168c4762a1bSJed Brown ierr = MatSetFromOptions(user->H);CHKERRQ(ierr); 169c4762a1bSJed Brown 170c4762a1bSJed Brown ierr = PetscStrncpy(filename,filebase,sizeof(filename));CHKERRQ(ierr); 171c4762a1bSJed Brown ierr = PetscStrlcat(filename,"Aeq",sizeof(filename));CHKERRQ(ierr); 1721e1ea65dSPierre Jolivet ierr = PetscViewerBinaryOpen(comm,filename,FILE_MODE_READ,&loader);CHKERRQ(ierr); 173c4762a1bSJed Brown if (ierr) { 174c4762a1bSJed Brown user->Aeq = NULL; 175c4762a1bSJed Brown user->me = 0; 176c4762a1bSJed Brown } else { 177c4762a1bSJed Brown ierr = MatCreate(comm,&user->Aeq);CHKERRQ(ierr); 178c4762a1bSJed Brown ierr = MatLoad(user->Aeq,loader);CHKERRQ(ierr); 179c4762a1bSJed Brown ierr = PetscViewerDestroy(&loader);CHKERRQ(ierr); 180c4762a1bSJed Brown ierr = MatGetSize(user->Aeq,&nrows,&ncols);CHKERRQ(ierr); 181*3c859ba3SBarry Smith PetscCheck(ncols == user->n,comm,PETSC_ERR_ARG_SIZ,"Aeq ncols != H nrows"); 182c4762a1bSJed Brown ierr = MatSetFromOptions(user->Aeq);CHKERRQ(ierr); 183c4762a1bSJed Brown user->me = nrows; 184c4762a1bSJed Brown } 185c4762a1bSJed Brown 186c4762a1bSJed Brown ierr = PetscStrncpy(filename,filebase,sizeof(filename));CHKERRQ(ierr); 187c4762a1bSJed Brown ierr = PetscStrlcat(filename,"Beq",sizeof(filename));CHKERRQ(ierr); 188c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(comm,filename,FILE_MODE_READ,&loader);CHKERRQ(ierr); 189c4762a1bSJed Brown if (ierr) { 190c4762a1bSJed Brown user->beq = 0; 191c4762a1bSJed Brown } else { 192c4762a1bSJed Brown ierr = VecCreate(comm,&user->beq);CHKERRQ(ierr); 193c4762a1bSJed Brown ierr = VecLoad(user->beq,loader);CHKERRQ(ierr); 194c4762a1bSJed Brown ierr = PetscViewerDestroy(&loader);CHKERRQ(ierr); 195c4762a1bSJed Brown ierr = VecGetSize(user->beq,&nrows);CHKERRQ(ierr); 196*3c859ba3SBarry Smith PetscCheck(nrows == user->me,comm,PETSC_ERR_ARG_SIZ,"Aeq nrows != Beq n"); 197c4762a1bSJed Brown ierr = VecSetFromOptions(user->beq);CHKERRQ(ierr); 198c4762a1bSJed Brown } 199c4762a1bSJed Brown 200c4762a1bSJed Brown user->mi = user->n; 201c4762a1bSJed Brown /* Ain = eye(n,n) */ 202c4762a1bSJed Brown ierr = MatCreate(comm,&user->Ain);CHKERRQ(ierr); 203c4762a1bSJed Brown ierr = MatSetType(user->Ain,MATAIJ);CHKERRQ(ierr); 204c4762a1bSJed Brown ierr = MatSetSizes(user->Ain,PETSC_DECIDE,PETSC_DECIDE,user->mi,user->mi);CHKERRQ(ierr); 205c4762a1bSJed Brown 206c4762a1bSJed Brown ierr = MatMPIAIJSetPreallocation(user->Ain,1,NULL,0,NULL);CHKERRQ(ierr); 207c4762a1bSJed Brown ierr = MatSeqAIJSetPreallocation(user->Ain,1,NULL);CHKERRQ(ierr); 208c4762a1bSJed Brown 209c4762a1bSJed Brown for (i=0;i<user->mi;i++) { 210c4762a1bSJed Brown ierr = MatSetValues(user->Ain,1,&i,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr); 211c4762a1bSJed Brown } 212c4762a1bSJed Brown ierr = MatAssemblyBegin(user->Ain,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 213c4762a1bSJed Brown ierr = MatAssemblyEnd(user->Ain,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 214c4762a1bSJed Brown ierr = MatSetFromOptions(user->Ain);CHKERRQ(ierr); 215c4762a1bSJed Brown 216c4762a1bSJed Brown /* bin = [0,0 ... 0]' */ 217c4762a1bSJed Brown ierr = VecCreate(comm,&user->bin);CHKERRQ(ierr); 218c4762a1bSJed Brown ierr = VecSetType(user->bin,VECMPI);CHKERRQ(ierr); 219c4762a1bSJed Brown ierr = VecSetSizes(user->bin,PETSC_DECIDE,user->mi);CHKERRQ(ierr); 220c4762a1bSJed Brown ierr = VecSet(user->bin,0.0);CHKERRQ(ierr); 221c4762a1bSJed Brown ierr = VecSetFromOptions(user->bin);CHKERRQ(ierr); 222c4762a1bSJed Brown user->m = user->me + user->mi; 223c4762a1bSJed Brown PetscFunctionReturn(0); 224c4762a1bSJed Brown } 225c4762a1bSJed Brown 226c4762a1bSJed Brown PetscErrorCode DestroyProblem(AppCtx *user) 227c4762a1bSJed Brown { 228c4762a1bSJed Brown PetscErrorCode ierr; 229c4762a1bSJed Brown 230c4762a1bSJed Brown PetscFunctionBegin; 231c4762a1bSJed Brown ierr = MatDestroy(&user->H);CHKERRQ(ierr); 232c4762a1bSJed Brown ierr = MatDestroy(&user->Aeq);CHKERRQ(ierr); 233c4762a1bSJed Brown ierr = MatDestroy(&user->Ain);CHKERRQ(ierr); 234c4762a1bSJed Brown ierr = VecDestroy(&user->beq);CHKERRQ(ierr); 235c4762a1bSJed Brown ierr = VecDestroy(&user->bin);CHKERRQ(ierr); 236c4762a1bSJed Brown ierr = VecDestroy(&user->d);CHKERRQ(ierr); 237c4762a1bSJed Brown PetscFunctionReturn(0); 238c4762a1bSJed Brown } 239c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao, Vec x, PetscReal *f, Vec g, void *ctx) 240c4762a1bSJed Brown { 241c4762a1bSJed Brown AppCtx *user = (AppCtx*)ctx; 242c4762a1bSJed Brown PetscScalar xtHx; 243c4762a1bSJed Brown PetscErrorCode ierr; 244c4762a1bSJed Brown 245c4762a1bSJed Brown PetscFunctionBegin; 246c4762a1bSJed Brown ierr = MatMult(user->H,x,g);CHKERRQ(ierr); 247c4762a1bSJed Brown ierr = VecDot(x,g,&xtHx);CHKERRQ(ierr); 248c4762a1bSJed Brown ierr = VecDot(x,user->d,f);CHKERRQ(ierr); 249c4762a1bSJed Brown *f += 0.5*xtHx; 250c4762a1bSJed Brown ierr = VecAXPY(g,1.0,user->d);CHKERRQ(ierr); 251c4762a1bSJed Brown PetscFunctionReturn(0); 252c4762a1bSJed Brown } 253c4762a1bSJed Brown 254c4762a1bSJed Brown PetscErrorCode FormHessian(Tao tao, Vec x, Mat H, Mat Hpre, void *ctx) 255c4762a1bSJed Brown { 256c4762a1bSJed Brown PetscFunctionBegin; 257c4762a1bSJed Brown PetscFunctionReturn(0); 258c4762a1bSJed Brown } 259c4762a1bSJed Brown 260c4762a1bSJed Brown PetscErrorCode FormInequalityConstraints(Tao tao, Vec x, Vec ci, void *ctx) 261c4762a1bSJed Brown { 262c4762a1bSJed Brown AppCtx *user = (AppCtx*)ctx; 263c4762a1bSJed Brown PetscErrorCode ierr; 264c4762a1bSJed Brown 265c4762a1bSJed Brown PetscFunctionBegin; 266c4762a1bSJed Brown ierr = MatMult(user->Ain,x,ci);CHKERRQ(ierr); 267c4762a1bSJed Brown PetscFunctionReturn(0); 268c4762a1bSJed Brown } 269c4762a1bSJed Brown 270c4762a1bSJed Brown PetscErrorCode FormEqualityConstraints(Tao tao, Vec x, Vec ce,void *ctx) 271c4762a1bSJed Brown { 272c4762a1bSJed Brown AppCtx *user = (AppCtx*)ctx; 273c4762a1bSJed Brown PetscErrorCode ierr; 274c4762a1bSJed Brown 275c4762a1bSJed Brown PetscFunctionBegin; 276c4762a1bSJed Brown ierr = MatMult(user->Aeq,x,ce);CHKERRQ(ierr); 277c4762a1bSJed Brown ierr = VecAXPY(ce,-1.0,user->beq);CHKERRQ(ierr); 278c4762a1bSJed Brown PetscFunctionReturn(0); 279c4762a1bSJed Brown } 280c4762a1bSJed Brown 281c4762a1bSJed Brown PetscErrorCode FormInequalityJacobian(Tao tao, Vec x, Mat JI, Mat JIpre, void *ctx) 282c4762a1bSJed Brown { 283c4762a1bSJed Brown PetscFunctionBegin; 284c4762a1bSJed Brown PetscFunctionReturn(0); 285c4762a1bSJed Brown } 286c4762a1bSJed Brown 287c4762a1bSJed Brown PetscErrorCode FormEqualityJacobian(Tao tao, Vec x, Mat JE, Mat JEpre, void *ctx) 288c4762a1bSJed Brown { 289c4762a1bSJed Brown PetscFunctionBegin; 290c4762a1bSJed Brown PetscFunctionReturn(0); 291c4762a1bSJed Brown } 292c4762a1bSJed Brown 293c4762a1bSJed Brown /*TEST 294c4762a1bSJed Brown 295c4762a1bSJed Brown build: 296c4762a1bSJed Brown requires: !complex 297c4762a1bSJed Brown 298c4762a1bSJed Brown test: 299c4762a1bSJed Brown requires: superlu 300c4762a1bSJed Brown localrunfiles: HS21 301c4762a1bSJed Brown 302c4762a1bSJed Brown TEST*/ 303