1c4762a1bSJed Brown /* XH: todo add cs1f.F90 and asjust makefile */
2c4762a1bSJed Brown /*
3c4762a1bSJed Brown Include "petsctao.h" so that we can use TAO solvers. Note that this
4c4762a1bSJed Brown file automatically includes libraries such as:
5c4762a1bSJed Brown petsc.h - base PETSc routines petscvec.h - vectors
6a5b23f4aSJose E. Roman petscsys.h - system routines petscmat.h - matrices
7c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods
8c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners
9c4762a1bSJed Brown
10c4762a1bSJed Brown */
11c4762a1bSJed Brown
12c4762a1bSJed Brown #include <petsctao.h>
13c4762a1bSJed Brown
14c4762a1bSJed Brown /*
15c4762a1bSJed Brown Description: Compressive sensing test example 1.
16c4762a1bSJed Brown 0.5*||Ax-b||^2 + lambda*||D*x||_1
17c4762a1bSJed Brown Xiang Huang: Nov 19, 2018
18c4762a1bSJed Brown
19c4762a1bSJed Brown Reference: None
20c4762a1bSJed Brown */
21c4762a1bSJed Brown
22c4762a1bSJed Brown static char help[] = "Finds the least-squares solution to the under constraint linear model Ax = b, with L1-norm regularizer. \n\
23c4762a1bSJed Brown A is a M*N real matrix (M<N), x is sparse. \n\
24c4762a1bSJed Brown We find the sparse solution by solving 0.5*||Ax-b||^2 + lambda*||D*x||_1, where lambda (by default 1e-4) is a user specified weight.\n\
25c4762a1bSJed Brown D is the K*N transform matrix so that D*x is sparse. By default D is identity matrix, so that D*x = x.\n";
26c4762a1bSJed Brown
27c4762a1bSJed Brown #define M 3
28c4762a1bSJed Brown #define N 5
29c4762a1bSJed Brown #define K 4
30c4762a1bSJed Brown
31c0b7dd19SHansol Suh typedef enum {
32c0b7dd19SHansol Suh TEST_L1DICT,
33c0b7dd19SHansol Suh TEST_LM,
34c0b7dd19SHansol Suh TEST_NONE
35c0b7dd19SHansol Suh } TestType;
36c0b7dd19SHansol Suh
37c4762a1bSJed Brown /* User-defined application context */
38c4762a1bSJed Brown typedef struct {
39c4762a1bSJed Brown /* Working space. linear least square: f(x) = A*x - b */
40c4762a1bSJed Brown PetscReal A[M][N]; /* array of coefficients */
41c4762a1bSJed Brown PetscReal b[M]; /* array of observations */
42c4762a1bSJed Brown PetscReal xGT[M]; /* array of ground truth object, which can be used to compare the reconstruction result */
43c4762a1bSJed Brown PetscReal D[K][N]; /* array of coefficients for 0.5*||Ax-b||^2 + lambda*||D*x||_1 */
44c4762a1bSJed Brown PetscReal J[M][N]; /* dense jacobian matrix array. For linear least square, J = A. For nonlinear least square, it is different from A */
45c4762a1bSJed Brown PetscInt idm[M]; /* Matrix row, column indices for jacobian and dictionary */
46c4762a1bSJed Brown PetscInt idn[N];
47c4762a1bSJed Brown PetscInt idk[K];
48c0b7dd19SHansol Suh TestType tType;
49c0b7dd19SHansol Suh PetscBool view_sol;
50c4762a1bSJed Brown } AppCtx;
51c4762a1bSJed Brown
52c4762a1bSJed Brown /* User provided Routines */
53c4762a1bSJed Brown PetscErrorCode InitializeUserData(AppCtx *);
54c4762a1bSJed Brown PetscErrorCode FormStartingPoint(Vec);
55c4762a1bSJed Brown PetscErrorCode FormDictionaryMatrix(Mat, AppCtx *);
56c4762a1bSJed Brown PetscErrorCode EvaluateFunction(Tao, Vec, Vec, void *);
57c4762a1bSJed Brown PetscErrorCode EvaluateJacobian(Tao, Vec, Mat, Mat, void *);
58c4762a1bSJed Brown
SetTaoOptionsFromUserOptions(Tao tao,AppCtx * ctx)59c0b7dd19SHansol Suh static PetscErrorCode SetTaoOptionsFromUserOptions(Tao tao, AppCtx *ctx)
60c0b7dd19SHansol Suh {
61c0b7dd19SHansol Suh PetscBool isbrgn;
62c0b7dd19SHansol Suh
63c0b7dd19SHansol Suh PetscFunctionBeginUser;
64c0b7dd19SHansol Suh PetscCall(PetscObjectTypeCompare((PetscObject)tao, TAOBRGN, &isbrgn));
65c0b7dd19SHansol Suh if (isbrgn) {
66c0b7dd19SHansol Suh switch (ctx->tType) {
67c0b7dd19SHansol Suh case TEST_LM:
68c0b7dd19SHansol Suh PetscCall(TaoBRGNSetRegularizationType(tao, TAOBRGN_REGULARIZATION_LM));
69c0b7dd19SHansol Suh break;
70c0b7dd19SHansol Suh case TEST_L1DICT:
71c0b7dd19SHansol Suh PetscCall(TaoBRGNSetRegularizationType(tao, TAOBRGN_REGULARIZATION_L1DICT));
72c0b7dd19SHansol Suh PetscCall(TaoBRGNSetRegularizerWeight(tao, 0.0001));
73c0b7dd19SHansol Suh PetscCall(TaoBRGNSetL1SmoothEpsilon(tao, 1.e-6));
74c0b7dd19SHansol Suh break;
75c0b7dd19SHansol Suh case TEST_NONE:
76c0b7dd19SHansol Suh default:
77c0b7dd19SHansol Suh break;
78c0b7dd19SHansol Suh }
79c0b7dd19SHansol Suh }
80c0b7dd19SHansol Suh PetscFunctionReturn(PETSC_SUCCESS);
81c0b7dd19SHansol Suh }
82c0b7dd19SHansol Suh
TestOutType(Tao tao,AppCtx * ctx)83c0b7dd19SHansol Suh static PetscErrorCode TestOutType(Tao tao, AppCtx *ctx)
84c0b7dd19SHansol Suh {
85c0b7dd19SHansol Suh PetscBool isbrgn;
86c0b7dd19SHansol Suh
87c0b7dd19SHansol Suh PetscFunctionBeginUser;
88c0b7dd19SHansol Suh PetscCall(PetscObjectTypeCompare((PetscObject)tao, TAOBRGN, &isbrgn));
89c0b7dd19SHansol Suh if (isbrgn) {
90c0b7dd19SHansol Suh TaoBRGNRegularizationType type;
91c0b7dd19SHansol Suh
92c0b7dd19SHansol Suh PetscCall(TaoBRGNGetRegularizationType(tao, &type));
93c0b7dd19SHansol Suh switch (ctx->tType) {
94c0b7dd19SHansol Suh case TEST_LM:
95c0b7dd19SHansol Suh PetscCheck(type == TAOBRGN_REGULARIZATION_LM, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_NOTSAMETYPE, "BRGN Regularization type is not LM!");
96c0b7dd19SHansol Suh break;
97c0b7dd19SHansol Suh case TEST_L1DICT:
98c0b7dd19SHansol Suh PetscCheck(type == TAOBRGN_REGULARIZATION_L1DICT, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_NOTSAMETYPE, "BRGN Regularization type is not L1DICT!");
99c0b7dd19SHansol Suh break;
100c0b7dd19SHansol Suh case TEST_NONE:
101c0b7dd19SHansol Suh default:
102c0b7dd19SHansol Suh break;
103c0b7dd19SHansol Suh }
104c0b7dd19SHansol Suh }
105c0b7dd19SHansol Suh PetscFunctionReturn(PETSC_SUCCESS);
106c0b7dd19SHansol Suh }
107c0b7dd19SHansol Suh
ProcessOptions(MPI_Comm comm,AppCtx * ctx)108c0b7dd19SHansol Suh static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *ctx)
109c0b7dd19SHansol Suh {
110c0b7dd19SHansol Suh const char *testTypes[3] = {"l1dict", "lm", "none"};
111c0b7dd19SHansol Suh PetscInt run;
112c0b7dd19SHansol Suh
113c0b7dd19SHansol Suh PetscFunctionBeginUser;
114c0b7dd19SHansol Suh ctx->tType = TEST_NONE;
115c0b7dd19SHansol Suh ctx->view_sol = PETSC_TRUE;
116c0b7dd19SHansol Suh PetscOptionsBegin(comm, "", "Least squares coverage", "");
117c0b7dd19SHansol Suh PetscCall(PetscOptionsBool("-view_sol", "Penalize deviation from both goals", "cs1.c", ctx->view_sol, &ctx->view_sol, NULL));
118c0b7dd19SHansol Suh run = ctx->tType;
119c0b7dd19SHansol Suh PetscCall(PetscOptionsEList("-test_type", "The coverage test type", "cs1.c", testTypes, 3, testTypes[ctx->tType], &run, NULL));
120c0b7dd19SHansol Suh ctx->tType = (TestType)run;
121c0b7dd19SHansol Suh PetscOptionsEnd();
122c0b7dd19SHansol Suh PetscFunctionReturn(PETSC_SUCCESS);
123c0b7dd19SHansol Suh }
124c0b7dd19SHansol Suh
125c4762a1bSJed Brown /*--------------------------------------------------------------------*/
main(int argc,char ** argv)126d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
127d71ae5a4SJacob Faibussowitsch {
128c4762a1bSJed Brown Vec x, f; /* solution, function f(x) = A*x-b */
129c4762a1bSJed Brown Mat J, D; /* Jacobian matrix, Transform matrix */
130c4762a1bSJed Brown Tao tao; /* Tao solver context */
131c4762a1bSJed Brown PetscInt i; /* iteration information */
132c4762a1bSJed Brown PetscReal hist[100], resid[100];
133c4762a1bSJed Brown PetscInt lits[100];
134c4762a1bSJed Brown AppCtx user; /* user-defined work context */
135c4762a1bSJed Brown
136327415f7SBarry Smith PetscFunctionBeginUser;
137c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
138c0b7dd19SHansol Suh PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
139c4762a1bSJed Brown /* Allocate solution and vector function vectors */
1409566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, N, &x));
1419566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, M, &f));
142c4762a1bSJed Brown
143c4762a1bSJed Brown /* Allocate Jacobian and Dictionary matrix. */
1449566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, M, N, NULL, &J));
1459566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, K, N, NULL, &D)); /* XH: TODO: dense -> sparse/dense/shell etc, do it on fly */
146c4762a1bSJed Brown
147c4762a1bSJed Brown for (i = 0; i < M; i++) user.idm[i] = i;
148c4762a1bSJed Brown for (i = 0; i < N; i++) user.idn[i] = i;
149c4762a1bSJed Brown for (i = 0; i < K; i++) user.idk[i] = i;
150c4762a1bSJed Brown
151c4762a1bSJed Brown /* Create TAO solver and set desired solution method */
1529566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_SELF, &tao));
1539566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOBRGN));
154c4762a1bSJed Brown
155c4762a1bSJed Brown /* User set application context: A, D matrice, and b vector. */
1569566063dSJacob Faibussowitsch PetscCall(InitializeUserData(&user));
157c4762a1bSJed Brown
158c4762a1bSJed Brown /* Set initial guess */
1599566063dSJacob Faibussowitsch PetscCall(FormStartingPoint(x));
160c4762a1bSJed Brown
161c4762a1bSJed Brown /* Fill the content of matrix D from user application Context */
1629566063dSJacob Faibussowitsch PetscCall(FormDictionaryMatrix(D, &user));
163c4762a1bSJed Brown
164c0b7dd19SHansol Suh /* If needed, set options via function for testing purpose */
165c0b7dd19SHansol Suh PetscCall(SetTaoOptionsFromUserOptions(tao, &user));
166c4762a1bSJed Brown /* Bind x to tao->solution. */
1679566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, x));
168c4762a1bSJed Brown /* Bind D to tao->data->D */
1699566063dSJacob Faibussowitsch PetscCall(TaoBRGNSetDictionaryMatrix(tao, D));
170c4762a1bSJed Brown
171c4762a1bSJed Brown /* Set the function and Jacobian routines. */
1729566063dSJacob Faibussowitsch PetscCall(TaoSetResidualRoutine(tao, f, EvaluateFunction, (void *)&user));
1739566063dSJacob Faibussowitsch PetscCall(TaoSetJacobianResidualRoutine(tao, J, J, EvaluateJacobian, (void *)&user));
174c4762a1bSJed Brown
175c4762a1bSJed Brown /* Check for any TAO command line arguments */
1769566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao));
177c4762a1bSJed Brown
1789566063dSJacob Faibussowitsch PetscCall(TaoSetConvergenceHistory(tao, hist, resid, 0, lits, 100, PETSC_TRUE));
179c4762a1bSJed Brown
180c4762a1bSJed Brown /* Perform the Solve */
1819566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao));
182c4762a1bSJed Brown
183c4762a1bSJed Brown /* XH: Debug: View the result, function and Jacobian. */
184c0b7dd19SHansol Suh if (user.view_sol) {
1859566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "-------- result x, residual f=A*x-b, and Jacobian=A. -------- \n"));
1869566063dSJacob Faibussowitsch PetscCall(VecView(x, PETSC_VIEWER_STDOUT_SELF));
1879566063dSJacob Faibussowitsch PetscCall(VecView(f, PETSC_VIEWER_STDOUT_SELF));
1889566063dSJacob Faibussowitsch PetscCall(MatView(J, PETSC_VIEWER_STDOUT_SELF));
1899566063dSJacob Faibussowitsch PetscCall(MatView(D, PETSC_VIEWER_STDOUT_SELF));
190c0b7dd19SHansol Suh }
191c0b7dd19SHansol Suh PetscCall(TestOutType(tao, &user));
192c4762a1bSJed Brown
193c4762a1bSJed Brown /* Free TAO data structures */
1949566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao));
195c4762a1bSJed Brown
196c4762a1bSJed Brown /* Free PETSc data structures */
1979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
1989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&f));
1999566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J));
2009566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D));
201c4762a1bSJed Brown
2029566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
203b122ec5aSJacob Faibussowitsch return 0;
204c4762a1bSJed Brown }
205c4762a1bSJed Brown
206c4762a1bSJed Brown /*--------------------------------------------------------------------*/
EvaluateFunction(Tao tao,Vec X,Vec F,void * ptr)207d71ae5a4SJacob Faibussowitsch PetscErrorCode EvaluateFunction(Tao tao, Vec X, Vec F, void *ptr)
208d71ae5a4SJacob Faibussowitsch {
209c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
210c4762a1bSJed Brown PetscInt m, n;
211c4762a1bSJed Brown const PetscReal *x;
212c4762a1bSJed Brown PetscReal *b = user->b, *f;
213c4762a1bSJed Brown
214c4762a1bSJed Brown PetscFunctionBegin;
2159566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x));
2169566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f));
217c4762a1bSJed Brown
218a5b23f4aSJose E. Roman /* Even for linear least square, we do not direct use matrix operation f = A*x - b now, just for future modification and compatibility for nonlinear least square */
219c4762a1bSJed Brown for (m = 0; m < M; m++) {
220c4762a1bSJed Brown f[m] = -b[m];
221ad540459SPierre Jolivet for (n = 0; n < N; n++) f[m] += user->A[m][n] * x[n];
222c4762a1bSJed Brown }
2239566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x));
2249566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f));
2253ba16761SJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * M * N));
2263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
227c4762a1bSJed Brown }
228c4762a1bSJed Brown
229c4762a1bSJed Brown /*------------------------------------------------------------*/
230c4762a1bSJed Brown /* J[m][n] = df[m]/dx[n] */
EvaluateJacobian(Tao tao,Vec X,Mat J,Mat Jpre,void * ptr)231d71ae5a4SJacob Faibussowitsch PetscErrorCode EvaluateJacobian(Tao tao, Vec X, Mat J, Mat Jpre, void *ptr)
232d71ae5a4SJacob Faibussowitsch {
233c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
234c4762a1bSJed Brown PetscInt m, n;
235c4762a1bSJed Brown const PetscReal *x;
236c4762a1bSJed Brown
237c4762a1bSJed Brown PetscFunctionBegin;
2389566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); /* not used for linear least square, but keep for future nonlinear least square) */
239c4762a1bSJed Brown /* XH: TODO: For linear least square, we can just set J=A fixed once, instead of keep update it! Maybe just create a function getFixedJacobian?
240c4762a1bSJed Brown For nonlinear least square, we require x to compute J, keep codes here for future nonlinear least square*/
241c4762a1bSJed Brown for (m = 0; m < M; ++m) {
242ad540459SPierre Jolivet for (n = 0; n < N; ++n) user->J[m][n] = user->A[m][n];
243c4762a1bSJed Brown }
244c4762a1bSJed Brown
2459566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, M, user->idm, N, user->idn, (PetscReal *)user->J, INSERT_VALUES));
2469566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
2479566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
248c4762a1bSJed Brown
2499566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); /* not used for linear least square, but keep for future nonlinear least square) */
2503ba16761SJacob Faibussowitsch PetscCall(PetscLogFlops(0)); /* 0 for linear least square, >0 for nonlinear least square */
2513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
252c4762a1bSJed Brown }
253c4762a1bSJed Brown
254c4762a1bSJed Brown /* ------------------------------------------------------------ */
255c4762a1bSJed Brown /* Currently fixed matrix, in future may be dynamic for D(x)? */
FormDictionaryMatrix(Mat D,AppCtx * user)256d71ae5a4SJacob Faibussowitsch PetscErrorCode FormDictionaryMatrix(Mat D, AppCtx *user)
257d71ae5a4SJacob Faibussowitsch {
258c4762a1bSJed Brown PetscFunctionBegin;
2599566063dSJacob Faibussowitsch PetscCall(MatSetValues(D, K, user->idk, N, user->idn, (PetscReal *)user->D, INSERT_VALUES));
2609566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(D, MAT_FINAL_ASSEMBLY));
2619566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(D, MAT_FINAL_ASSEMBLY));
262c4762a1bSJed Brown
2633ba16761SJacob Faibussowitsch PetscCall(PetscLogFlops(0)); /* 0 for fixed dictionary matrix, >0 for varying dictionary matrix */
2643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
265c4762a1bSJed Brown }
266c4762a1bSJed Brown
267c4762a1bSJed Brown /* ------------------------------------------------------------ */
FormStartingPoint(Vec X)268d71ae5a4SJacob Faibussowitsch PetscErrorCode FormStartingPoint(Vec X)
269d71ae5a4SJacob Faibussowitsch {
270c4762a1bSJed Brown PetscFunctionBegin;
2719566063dSJacob Faibussowitsch PetscCall(VecSet(X, 0.0));
2723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
273c4762a1bSJed Brown }
274c4762a1bSJed Brown
275c4762a1bSJed Brown /* ---------------------------------------------------------------------- */
InitializeUserData(AppCtx * user)276d71ae5a4SJacob Faibussowitsch PetscErrorCode InitializeUserData(AppCtx *user)
277d71ae5a4SJacob Faibussowitsch {
278da81f932SPierre Jolivet PetscReal *b = user->b; /* **A=user->A, but we don't know the dimension of A in this way, how to fix? */
279c4762a1bSJed Brown PetscInt m, n, k; /* loop index for M,N,K dimension. */
280c4762a1bSJed Brown
281c4762a1bSJed Brown PetscFunctionBegin;
282c4762a1bSJed Brown /* b = A*x while x = [0;0;1;0;0] here*/
283c4762a1bSJed Brown m = 0;
284c4762a1bSJed Brown b[m++] = 0.28;
285c4762a1bSJed Brown b[m++] = 0.55;
286c4762a1bSJed Brown b[m++] = 0.96;
287c4762a1bSJed Brown
28821afe8ebSBarry Smith /* MATLAB generated random matrix, uniformly distributed in [0,1] with 2 digits accuracy. rng(0); A = rand(M, N); A = round(A*100)/100;
289c4762a1bSJed Brown A = [0.81 0.91 0.28 0.96 0.96
290c4762a1bSJed Brown 0.91 0.63 0.55 0.16 0.49
291c4762a1bSJed Brown 0.13 0.10 0.96 0.97 0.80]
292c4762a1bSJed Brown */
2939371c9d4SSatish Balay m = 0;
2949371c9d4SSatish Balay n = 0;
2959371c9d4SSatish Balay user->A[m][n++] = 0.81;
2969371c9d4SSatish Balay user->A[m][n++] = 0.91;
2979371c9d4SSatish Balay user->A[m][n++] = 0.28;
2989371c9d4SSatish Balay user->A[m][n++] = 0.96;
2999371c9d4SSatish Balay user->A[m][n++] = 0.96;
3009371c9d4SSatish Balay ++m;
3019371c9d4SSatish Balay n = 0;
3029371c9d4SSatish Balay user->A[m][n++] = 0.91;
3039371c9d4SSatish Balay user->A[m][n++] = 0.63;
3049371c9d4SSatish Balay user->A[m][n++] = 0.55;
3059371c9d4SSatish Balay user->A[m][n++] = 0.16;
3069371c9d4SSatish Balay user->A[m][n++] = 0.49;
3079371c9d4SSatish Balay ++m;
3089371c9d4SSatish Balay n = 0;
3099371c9d4SSatish Balay user->A[m][n++] = 0.13;
3109371c9d4SSatish Balay user->A[m][n++] = 0.10;
3119371c9d4SSatish Balay user->A[m][n++] = 0.96;
3129371c9d4SSatish Balay user->A[m][n++] = 0.97;
3139371c9d4SSatish Balay user->A[m][n++] = 0.80;
314c4762a1bSJed Brown
315c4762a1bSJed Brown /* initialize to 0 */
316c4762a1bSJed Brown for (k = 0; k < K; k++) {
317ad540459SPierre Jolivet for (n = 0; n < N; n++) user->D[k][n] = 0.0;
318c4762a1bSJed Brown }
319c4762a1bSJed Brown /* Choice I: set D to identity matrix of size N*N for testing */
320c4762a1bSJed Brown /* for (k=0; k<K; k++) user->D[k][k] = 1.0; */
321c4762a1bSJed Brown /* Choice II: set D to Backward difference matrix of size (N-1)*N, with zero extended boundary assumption */
322c4762a1bSJed Brown for (k = 0; k < K; k++) {
323c4762a1bSJed Brown user->D[k][k] = -1.0;
324c4762a1bSJed Brown user->D[k][k + 1] = 1.0;
325c4762a1bSJed Brown }
3263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
327c4762a1bSJed Brown }
328c4762a1bSJed Brown
329c4762a1bSJed Brown /*TEST
330c4762a1bSJed Brown
331c4762a1bSJed Brown build:
332c0b7dd19SHansol Suh requires: !complex !single !quad !defined(PETSC_USE_64BIT_INDICES) !__float128
333c4762a1bSJed Brown
334c4762a1bSJed Brown test:
335c4762a1bSJed Brown localrunfiles: cs1Data_A_b_xGT
33610978b7dSBarry Smith args: -tao_monitor_short -tao_max_it 100 -tao_type pounders -tao_gatol 1.e-6
337c4762a1bSJed Brown
338c4762a1bSJed Brown test:
339c4762a1bSJed Brown suffix: 2
340c4762a1bSJed Brown localrunfiles: cs1Data_A_b_xGT
341*a336c150SZach Atkins args: -tao_monitor -tao_max_it 100 -tao_type brgn -tao_brgn_regularization_type l2prox -tao_brgn_regularizer_weight 1e-8 -tao_gatol 1.e-6 -tao_brgn_subsolver_tao_bnk_ksp_converged_reason -tao_brgn_subsolver_tao_monitor
342c4762a1bSJed Brown
343c4762a1bSJed Brown test:
344c4762a1bSJed Brown suffix: 3
345c4762a1bSJed Brown localrunfiles: cs1Data_A_b_xGT
346*a336c150SZach Atkins args: -tao_monitor -tao_max_it 100 -tao_type brgn -tao_brgn_regularization_type l1dict -tao_brgn_regularizer_weight 1e-8 -tao_brgn_l1_smooth_epsilon 1e-6 -tao_gatol 1.e-6 -tao_brgn_subsolver_tao_monitor
347c4762a1bSJed Brown
348c4762a1bSJed Brown test:
349c4762a1bSJed Brown suffix: 4
350c4762a1bSJed Brown localrunfiles: cs1Data_A_b_xGT
351*a336c150SZach Atkins args: -tao_monitor -tao_max_it 100 -tao_type brgn -tao_brgn_regularization_type l2pure -tao_brgn_regularizer_weight 1e-8 -tao_gatol 1.e-6 -tao_brgn_subsolver_tao_monitor
352c4762a1bSJed Brown
353cd1c4666STristan Konolige test:
354cd1c4666STristan Konolige suffix: 5
355cd1c4666STristan Konolige localrunfiles: cs1Data_A_b_xGT
356*a336c150SZach Atkins args: -tao_monitor -tao_max_it 100 -tao_type brgn -tao_brgn_regularization_type lm -tao_gatol 1.e-6 -tao_brgn_subsolver_tao_type bnls -tao_brgn_subsolver_tao_monitor
357cd1c4666STristan Konolige
358c0b7dd19SHansol Suh test:
359c0b7dd19SHansol Suh suffix: view_lm
360c0b7dd19SHansol Suh localrunfiles: cs1Data_A_b_xGT
361c0b7dd19SHansol Suh args: -tao_type brgn -test_type lm -tao_gatol 1.e-6 -view_sol 0 -tao_view
362c0b7dd19SHansol Suh
363c0b7dd19SHansol Suh test:
364c0b7dd19SHansol Suh suffix: view_l1dict
365c0b7dd19SHansol Suh localrunfiles: cs1Data_A_b_xGT
366c0b7dd19SHansol Suh args: -tao_type brgn -test_type l1dict -tao_gatol 1.e-6 -view_sol 0 -tao_view
367c0b7dd19SHansol Suh
368c4762a1bSJed Brown TEST*/
369