xref: /petsc/src/tao/constrained/tutorials/maros.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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;
80*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
81c4762a1bSJed Brown   /* Specify default parameters for the problem, check for command-line overrides */
82*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrncpy(user.name,"HS21",sizeof(user.name)));
83*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-cutername",user.name,sizeof(user.name),&flg));
84c4762a1bSJed Brown 
85*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\n---- MAROS Problem %s -----\n",user.name));
86*5f80ce2aSJacob Faibussowitsch   CHKERRQ(InitializeProblem(&user));
87*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(user.d,&x));
88*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(user.beq,&ceq));
89*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(user.bin,&cin));
90*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(x,1.0));
91c4762a1bSJed Brown 
92*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoCreate(PETSC_COMM_WORLD,&tao));
93*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetType(tao,TAOIPM));
94*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetSolution(tao,x));
95*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void*)&user));
96*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetEqualityConstraintsRoutine(tao,ceq,FormEqualityConstraints,(void*)&user));
97*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetInequalityConstraintsRoutine(tao,cin,FormInequalityConstraints,(void*)&user));
98*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetInequalityBounds(tao,user.bin,NULL));
99*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetJacobianEqualityRoutine(tao,user.Aeq,user.Aeq,FormEqualityJacobian,(void*)&user));
100*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetJacobianInequalityRoutine(tao,user.Ain,user.Ain,FormInequalityJacobian,(void*)&user));
101*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetHessian(tao,user.H,user.H,FormHessian,(void*)&user));
102*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoGetKSP(tao,&ksp));
103*5f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPGetPC(ksp,&pc));
104*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PCSetType(pc,PCLU));
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   */
109*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PCFactorSetMatSolverType(pc,MATSOLVERSUPERLU));
110*5f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetType(ksp,KSPPREONLY));
111*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetTolerances(tao,0,0,0));
112c4762a1bSJed Brown 
113*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetFromOptions(tao));
114*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSolve(tao));
115*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoGetConvergedReason(tao,&reason));
116c4762a1bSJed Brown   if (reason < 0) {
117*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(MPI_COMM_WORLD, "TAO failed to converge due to %s.\n",TaoConvergedReasons[reason]));
118c4762a1bSJed Brown   } else {
119*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(MPI_COMM_WORLD, "Optimization completed with status %s.\n",TaoConvergedReasons[reason]));
120c4762a1bSJed Brown   }
121c4762a1bSJed Brown 
122*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DestroyProblem(&user));
123*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
124*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&ceq));
125*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&cin));
126*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoDestroy(&tao));
127c4762a1bSJed Brown 
128c4762a1bSJed Brown   ierr = PetscFinalize();
129c4762a1bSJed Brown   return ierr;
130c4762a1bSJed Brown }
131c4762a1bSJed Brown 
132c4762a1bSJed Brown PetscErrorCode InitializeProblem(AppCtx *user)
133c4762a1bSJed Brown {
134c4762a1bSJed Brown   PetscViewer loader;
135c4762a1bSJed Brown   MPI_Comm    comm;
136c4762a1bSJed Brown   PetscInt    nrows,ncols,i;
137c4762a1bSJed Brown   PetscScalar one = 1.0;
138c4762a1bSJed Brown   char        filebase[128];
139c4762a1bSJed Brown   char        filename[128];
140c4762a1bSJed Brown 
141c4762a1bSJed Brown   PetscFunctionBegin;
142c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
143*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrncpy(filebase,user->name,sizeof(filebase)));
144*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrlcat(filebase,"/",sizeof(filebase)));
145*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrncpy(filename,filebase,sizeof(filename)));
146*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrlcat(filename,"f",sizeof(filename)));
147*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(comm,filename,FILE_MODE_READ,&loader));
148c4762a1bSJed Brown 
149*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(comm,&user->d));
150*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecLoad(user->d,loader));
151*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&loader));
152*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(user->d,&nrows));
153*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(user->d));
154c4762a1bSJed Brown   user->n = nrows;
155c4762a1bSJed Brown 
156*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrncpy(filename,filebase,sizeof(filename)));
157*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrlcat(filename,"H",sizeof(filename)));
158*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(comm,filename,FILE_MODE_READ,&loader));
159c4762a1bSJed Brown 
160*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(comm,&user->H));
161*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(user->H,PETSC_DECIDE,PETSC_DECIDE,nrows,nrows));
162*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(user->H,loader));
163*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&loader));
164*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetSize(user->H,&nrows,&ncols));
1653c859ba3SBarry Smith   PetscCheck(nrows == user->n,comm,PETSC_ERR_ARG_SIZ,"H: nrows != n");
1663c859ba3SBarry Smith   PetscCheck(ncols == user->n,comm,PETSC_ERR_ARG_SIZ,"H: ncols != n");
167*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(user->H));
168c4762a1bSJed Brown 
169*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrncpy(filename,filebase,sizeof(filename)));
170*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrlcat(filename,"Aeq",sizeof(filename)));
171*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(comm,filename,FILE_MODE_READ,&loader));
172*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(comm,&user->Aeq));
173*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(user->Aeq,loader));
174*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&loader));
175*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetSize(user->Aeq,&nrows,&ncols));
1763c859ba3SBarry Smith   PetscCheck(ncols == user->n,comm,PETSC_ERR_ARG_SIZ,"Aeq ncols != H nrows");
177*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(user->Aeq));
178c4762a1bSJed Brown   user->me = nrows;
179c4762a1bSJed Brown 
180*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrncpy(filename,filebase,sizeof(filename)));
181*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrlcat(filename,"Beq",sizeof(filename)));
182*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(comm,filename,FILE_MODE_READ,&loader));
183*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(comm,&user->beq));
184*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecLoad(user->beq,loader));
185*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&loader));
186*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(user->beq,&nrows));
1873c859ba3SBarry Smith   PetscCheck(nrows == user->me,comm,PETSC_ERR_ARG_SIZ,"Aeq nrows != Beq n");
188*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(user->beq));
189c4762a1bSJed Brown 
190c4762a1bSJed Brown   user->mi = user->n;
191c4762a1bSJed Brown   /* Ain = eye(n,n) */
192*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(comm,&user->Ain));
193*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(user->Ain,MATAIJ));
194*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(user->Ain,PETSC_DECIDE,PETSC_DECIDE,user->mi,user->mi));
195c4762a1bSJed Brown 
196*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMPIAIJSetPreallocation(user->Ain,1,NULL,0,NULL));
197*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqAIJSetPreallocation(user->Ain,1,NULL));
198c4762a1bSJed Brown 
199*5f80ce2aSJacob Faibussowitsch   for (i=0;i<user->mi;i++) CHKERRQ(MatSetValues(user->Ain,1,&i,1,&i,&one,INSERT_VALUES));
200*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(user->Ain,MAT_FINAL_ASSEMBLY));
201*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(user->Ain,MAT_FINAL_ASSEMBLY));
202*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(user->Ain));
203c4762a1bSJed Brown 
204c4762a1bSJed Brown   /* bin = [0,0 ... 0]' */
205*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(comm,&user->bin));
206*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetType(user->bin,VECMPI));
207*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(user->bin,PETSC_DECIDE,user->mi));
208*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(user->bin,0.0));
209*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(user->bin));
210c4762a1bSJed Brown   user->m = user->me + user->mi;
211c4762a1bSJed Brown   PetscFunctionReturn(0);
212c4762a1bSJed Brown }
213c4762a1bSJed Brown 
214c4762a1bSJed Brown PetscErrorCode DestroyProblem(AppCtx *user)
215c4762a1bSJed Brown {
216c4762a1bSJed Brown   PetscFunctionBegin;
217*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user->H));
218*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user->Aeq));
219*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user->Ain));
220*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user->beq));
221*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user->bin));
222*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user->d));
223c4762a1bSJed Brown   PetscFunctionReturn(0);
224c4762a1bSJed Brown }
225*5f80ce2aSJacob Faibussowitsch 
226c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao, Vec x, PetscReal *f, Vec g, void *ctx)
227c4762a1bSJed Brown {
228c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ctx;
229c4762a1bSJed Brown   PetscScalar    xtHx;
230c4762a1bSJed Brown 
231c4762a1bSJed Brown   PetscFunctionBegin;
232*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(user->H,x,g));
233*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDot(x,g,&xtHx));
234*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDot(x,user->d,f));
235c4762a1bSJed Brown   *f += 0.5*xtHx;
236*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(g,1.0,user->d));
237c4762a1bSJed Brown   PetscFunctionReturn(0);
238c4762a1bSJed Brown }
239c4762a1bSJed Brown 
240c4762a1bSJed Brown PetscErrorCode FormHessian(Tao tao, Vec x, Mat H, Mat Hpre, void *ctx)
241c4762a1bSJed Brown {
242c4762a1bSJed Brown   PetscFunctionBegin;
243c4762a1bSJed Brown   PetscFunctionReturn(0);
244c4762a1bSJed Brown }
245c4762a1bSJed Brown 
246c4762a1bSJed Brown PetscErrorCode FormInequalityConstraints(Tao tao, Vec x, Vec ci, void *ctx)
247c4762a1bSJed Brown {
248c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ctx;
249c4762a1bSJed Brown 
250c4762a1bSJed Brown   PetscFunctionBegin;
251*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(user->Ain,x,ci));
252c4762a1bSJed Brown   PetscFunctionReturn(0);
253c4762a1bSJed Brown }
254c4762a1bSJed Brown 
255c4762a1bSJed Brown PetscErrorCode FormEqualityConstraints(Tao tao, Vec x, Vec ce,void *ctx)
256c4762a1bSJed Brown {
257c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ctx;
258c4762a1bSJed Brown 
259c4762a1bSJed Brown   PetscFunctionBegin;
260*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(user->Aeq,x,ce));
261*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(ce,-1.0,user->beq));
262c4762a1bSJed Brown   PetscFunctionReturn(0);
263c4762a1bSJed Brown }
264c4762a1bSJed Brown 
265c4762a1bSJed Brown PetscErrorCode FormInequalityJacobian(Tao tao, Vec x, Mat JI, Mat JIpre,  void *ctx)
266c4762a1bSJed Brown {
267c4762a1bSJed Brown   PetscFunctionBegin;
268c4762a1bSJed Brown   PetscFunctionReturn(0);
269c4762a1bSJed Brown }
270c4762a1bSJed Brown 
271c4762a1bSJed Brown PetscErrorCode FormEqualityJacobian(Tao tao, Vec x, Mat JE, Mat JEpre, void *ctx)
272c4762a1bSJed Brown {
273c4762a1bSJed Brown   PetscFunctionBegin;
274c4762a1bSJed Brown   PetscFunctionReturn(0);
275c4762a1bSJed Brown }
276c4762a1bSJed Brown 
277c4762a1bSJed Brown /*TEST
278c4762a1bSJed Brown 
279c4762a1bSJed Brown    build:
280c4762a1bSJed Brown       requires: !complex
281c4762a1bSJed Brown 
282c4762a1bSJed Brown    test:
283c4762a1bSJed Brown       requires: superlu
284c4762a1bSJed Brown       localrunfiles: HS21
285c4762a1bSJed Brown 
286c4762a1bSJed Brown TEST*/
287