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