1c4762a1bSJed Brown #include <petsctao.h> 2c4762a1bSJed Brown /* 3c4762a1bSJed Brown Description: ADMM tomography reconstruction example . 4c4762a1bSJed Brown 0.5*||Ax-b||^2 + lambda*g(x) 5c4762a1bSJed Brown Reference: BRGN Tomography Example 6c4762a1bSJed Brown */ 7c4762a1bSJed Brown 8c4762a1bSJed Brown static char help[] = "Finds the ADMM solution to the under constraint linear model Ax = b, with regularizer. \n\ 9c4762a1bSJed Brown A is a M*N real matrix (M<N), x is sparse. A good regularizer is an L1 regularizer. \n\ 10c4762a1bSJed Brown We first split the operator into 0.5*||Ax-b||^2, f(x), and lambda*||x||_1, g(z), where lambda is user specified weight. \n\ 11c4762a1bSJed Brown g(z) could be either ||z||_1, or ||z||_2^2. Default closed form solution for NORM1 would be soft-threshold, which is \n\ 12c4762a1bSJed Brown natively supported in admm.c with -tao_admm_regularizer_type soft-threshold. Or user can use regular TAO solver for \n\ 13c4762a1bSJed Brown either NORM1 or NORM2 or TAOSHELL, with -reg {1,2,3} \n\ 14c4762a1bSJed Brown Then, we augment both f and g, and solve it via ADMM. \n\ 15c4762a1bSJed Brown D is the M*N transform matrix so that D*x is sparse. \n"; 16c4762a1bSJed Brown 17c4762a1bSJed Brown typedef struct { 18c4762a1bSJed Brown PetscInt M, N, K, reg; 19c4762a1bSJed Brown PetscReal lambda, eps, mumin; 20c4762a1bSJed Brown Mat A, ATA, H, Hx, D, Hz, DTD, HF; 21c4762a1bSJed Brown Vec c, xlb, xub, x, b, workM, workN, workN2, workN3, xGT; /* observation b, ground truth xGT, the lower bound and upper bound of x*/ 22c4762a1bSJed Brown } AppCtx; 23c4762a1bSJed Brown 24c4762a1bSJed Brown /*------------------------------------------------------------*/ 25c4762a1bSJed Brown 269371c9d4SSatish Balay PetscErrorCode NullJacobian(Tao tao, Vec X, Mat J, Mat Jpre, void *ptr) { 27c4762a1bSJed Brown PetscFunctionBegin; 28c4762a1bSJed Brown PetscFunctionReturn(0); 29c4762a1bSJed Brown } 30c4762a1bSJed Brown 31c4762a1bSJed Brown /*------------------------------------------------------------*/ 32c4762a1bSJed Brown 339371c9d4SSatish Balay static PetscErrorCode TaoShellSolve_SoftThreshold(Tao tao) { 34c4762a1bSJed Brown PetscReal lambda, mu; 35c4762a1bSJed Brown AppCtx *user; 36c4762a1bSJed Brown Vec out, work, y, x; 37c4762a1bSJed Brown Tao admm_tao, misfit; 38c4762a1bSJed Brown 39c4762a1bSJed Brown PetscFunctionBegin; 40c4762a1bSJed Brown user = NULL; 41c4762a1bSJed Brown mu = 0; 429566063dSJacob Faibussowitsch PetscCall(TaoGetADMMParentTao(tao, &admm_tao)); 439566063dSJacob Faibussowitsch PetscCall(TaoADMMGetMisfitSubsolver(admm_tao, &misfit)); 449566063dSJacob Faibussowitsch PetscCall(TaoADMMGetSpectralPenalty(admm_tao, &mu)); 459566063dSJacob Faibussowitsch PetscCall(TaoShellGetContext(tao, &user)); 46c4762a1bSJed Brown 47c4762a1bSJed Brown lambda = user->lambda; 48c4762a1bSJed Brown work = user->workN; 499566063dSJacob Faibussowitsch PetscCall(TaoGetSolution(tao, &out)); 509566063dSJacob Faibussowitsch PetscCall(TaoGetSolution(misfit, &x)); 519566063dSJacob Faibussowitsch PetscCall(TaoADMMGetDualVector(admm_tao, &y)); 52c4762a1bSJed Brown 53c4762a1bSJed Brown /* Dx + y/mu */ 549566063dSJacob Faibussowitsch PetscCall(MatMult(user->D, x, work)); 559566063dSJacob Faibussowitsch PetscCall(VecAXPY(work, 1 / mu, y)); 56c4762a1bSJed Brown 57c4762a1bSJed Brown /* soft thresholding */ 589566063dSJacob Faibussowitsch PetscCall(TaoSoftThreshold(work, -lambda / mu, lambda / mu, out)); 59c4762a1bSJed Brown PetscFunctionReturn(0); 60c4762a1bSJed Brown } 61c4762a1bSJed Brown 62c4762a1bSJed Brown /*------------------------------------------------------------*/ 63c4762a1bSJed Brown 649371c9d4SSatish Balay PetscErrorCode MisfitObjectiveAndGradient(Tao tao, Vec X, PetscReal *f, Vec g, void *ptr) { 65c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 66c4762a1bSJed Brown 67c4762a1bSJed Brown PetscFunctionBegin; 68c4762a1bSJed Brown /* Objective 0.5*||Ax-b||_2^2 */ 699566063dSJacob Faibussowitsch PetscCall(MatMult(user->A, X, user->workM)); 709566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->workM, -1, user->b)); 719566063dSJacob Faibussowitsch PetscCall(VecDot(user->workM, user->workM, f)); 72c4762a1bSJed Brown *f *= 0.5; 73c4762a1bSJed Brown /* Gradient. ATAx-ATb */ 749566063dSJacob Faibussowitsch PetscCall(MatMult(user->ATA, X, user->workN)); 759566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->A, user->b, user->workN2)); 769566063dSJacob Faibussowitsch PetscCall(VecWAXPY(g, -1., user->workN2, user->workN)); 77c4762a1bSJed Brown PetscFunctionReturn(0); 78c4762a1bSJed Brown } 79c4762a1bSJed Brown 80c4762a1bSJed Brown /*------------------------------------------------------------*/ 81c4762a1bSJed Brown 829371c9d4SSatish Balay PetscErrorCode RegularizerObjectiveAndGradient1(Tao tao, Vec X, PetscReal *f_reg, Vec G_reg, void *ptr) { 83c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 84c4762a1bSJed Brown 85c4762a1bSJed Brown PetscFunctionBegin; 86c4762a1bSJed Brown /* compute regularizer objective 87c4762a1bSJed Brown * f = f + lambda*sum(sqrt(y.^2+epsilon^2) - epsilon), where y = D*x */ 889566063dSJacob Faibussowitsch PetscCall(VecCopy(X, user->workN2)); 899566063dSJacob Faibussowitsch PetscCall(VecPow(user->workN2, 2.)); 909566063dSJacob Faibussowitsch PetscCall(VecShift(user->workN2, user->eps * user->eps)); 919566063dSJacob Faibussowitsch PetscCall(VecSqrtAbs(user->workN2)); 929566063dSJacob Faibussowitsch PetscCall(VecCopy(user->workN2, user->workN3)); 939566063dSJacob Faibussowitsch PetscCall(VecShift(user->workN2, -user->eps)); 949566063dSJacob Faibussowitsch PetscCall(VecSum(user->workN2, f_reg)); 95c4762a1bSJed Brown *f_reg *= user->lambda; 96c4762a1bSJed Brown /* compute regularizer gradient = lambda*x */ 979566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(G_reg, X, user->workN3)); 989566063dSJacob Faibussowitsch PetscCall(VecScale(G_reg, user->lambda)); 99c4762a1bSJed Brown PetscFunctionReturn(0); 100c4762a1bSJed Brown } 101c4762a1bSJed Brown 102c4762a1bSJed Brown /*------------------------------------------------------------*/ 103c4762a1bSJed Brown 1049371c9d4SSatish Balay PetscErrorCode RegularizerObjectiveAndGradient2(Tao tao, Vec X, PetscReal *f_reg, Vec G_reg, void *ptr) { 105c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 106c4762a1bSJed Brown PetscReal temp; 107c4762a1bSJed Brown 108c4762a1bSJed Brown PetscFunctionBegin; 109c4762a1bSJed Brown /* compute regularizer objective = lambda*|z|_2^2 */ 1109566063dSJacob Faibussowitsch PetscCall(VecDot(X, X, &temp)); 111c4762a1bSJed Brown *f_reg = 0.5 * user->lambda * temp; 112c4762a1bSJed Brown /* compute regularizer gradient = lambda*z */ 1139566063dSJacob Faibussowitsch PetscCall(VecCopy(X, G_reg)); 1149566063dSJacob Faibussowitsch PetscCall(VecScale(G_reg, user->lambda)); 115c4762a1bSJed Brown PetscFunctionReturn(0); 116c4762a1bSJed Brown } 117c4762a1bSJed Brown 118c4762a1bSJed Brown /*------------------------------------------------------------*/ 119c4762a1bSJed Brown 1209371c9d4SSatish Balay static PetscErrorCode HessianMisfit(Tao tao, Vec x, Mat H, Mat Hpre, void *ptr) { 121c4762a1bSJed Brown PetscFunctionBegin; 122c4762a1bSJed Brown PetscFunctionReturn(0); 123c4762a1bSJed Brown } 124c4762a1bSJed Brown 125c4762a1bSJed Brown /*------------------------------------------------------------*/ 126c4762a1bSJed Brown 1279371c9d4SSatish Balay static PetscErrorCode HessianReg(Tao tao, Vec x, Mat H, Mat Hpre, void *ptr) { 128c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 129c4762a1bSJed Brown 130c4762a1bSJed Brown PetscFunctionBegin; 1319566063dSJacob Faibussowitsch PetscCall(MatMult(user->D, x, user->workN)); 1329566063dSJacob Faibussowitsch PetscCall(VecPow(user->workN2, 2.)); 1339566063dSJacob Faibussowitsch PetscCall(VecShift(user->workN2, user->eps * user->eps)); 1349566063dSJacob Faibussowitsch PetscCall(VecSqrtAbs(user->workN2)); 1359566063dSJacob Faibussowitsch PetscCall(VecShift(user->workN2, -user->eps)); 1369566063dSJacob Faibussowitsch PetscCall(VecReciprocal(user->workN2)); 1379566063dSJacob Faibussowitsch PetscCall(VecScale(user->workN2, user->eps * user->eps)); 1389566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(H, user->workN2, INSERT_VALUES)); 139c4762a1bSJed Brown PetscFunctionReturn(0); 140c4762a1bSJed Brown } 141c4762a1bSJed Brown 142c4762a1bSJed Brown /*------------------------------------------------------------*/ 143c4762a1bSJed Brown 1449371c9d4SSatish Balay PetscErrorCode FullObjGrad(Tao tao, Vec X, PetscReal *f, Vec g, void *ptr) { 145c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 146c4762a1bSJed Brown PetscReal f_reg; 147c4762a1bSJed Brown 148c4762a1bSJed Brown PetscFunctionBegin; 149c4762a1bSJed Brown /* Objective 0.5*||Ax-b||_2^2 + lambda*||x||_2^2*/ 1509566063dSJacob Faibussowitsch PetscCall(MatMult(user->A, X, user->workM)); 1519566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->workM, -1, user->b)); 1529566063dSJacob Faibussowitsch PetscCall(VecDot(user->workM, user->workM, f)); 1539566063dSJacob Faibussowitsch PetscCall(VecNorm(X, NORM_2, &f_reg)); 154c4762a1bSJed Brown *f *= 0.5; 155c4762a1bSJed Brown *f += user->lambda * f_reg * f_reg; 156c4762a1bSJed Brown /* Gradient. ATAx-ATb + 2*lambda*x */ 1579566063dSJacob Faibussowitsch PetscCall(MatMult(user->ATA, X, user->workN)); 1589566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->A, user->b, user->workN2)); 1599566063dSJacob Faibussowitsch PetscCall(VecWAXPY(g, -1., user->workN2, user->workN)); 1609566063dSJacob Faibussowitsch PetscCall(VecAXPY(g, 2 * user->lambda, X)); 161c4762a1bSJed Brown PetscFunctionReturn(0); 162c4762a1bSJed Brown } 163c4762a1bSJed Brown /*------------------------------------------------------------*/ 164c4762a1bSJed Brown 1659371c9d4SSatish Balay static PetscErrorCode HessianFull(Tao tao, Vec x, Mat H, Mat Hpre, void *ptr) { 166c4762a1bSJed Brown PetscFunctionBegin; 167c4762a1bSJed Brown PetscFunctionReturn(0); 168c4762a1bSJed Brown } 169c4762a1bSJed Brown /*------------------------------------------------------------*/ 170c4762a1bSJed Brown 1719371c9d4SSatish Balay PetscErrorCode InitializeUserData(AppCtx *user) { 172c4762a1bSJed Brown char dataFile[] = "tomographyData_A_b_xGT"; /* Matrix A and vectors b, xGT(ground truth) binary files generated by Matlab. Debug: change from "tomographyData_A_b_xGT" to "cs1Data_A_b_xGT". */ 173c4762a1bSJed Brown PetscViewer fd; /* used to load data from file */ 174c4762a1bSJed Brown PetscInt k, n; 175c4762a1bSJed Brown PetscScalar v; 176c4762a1bSJed Brown 177d0609cedSBarry Smith PetscFunctionBegin; 178c4762a1bSJed Brown /* Load the A matrix, b vector, and xGT vector from a binary file. */ 1799566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, dataFile, FILE_MODE_READ, &fd)); 1809566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->A)); 1819566063dSJacob Faibussowitsch PetscCall(MatSetType(user->A, MATAIJ)); 1829566063dSJacob Faibussowitsch PetscCall(MatLoad(user->A, fd)); 1839566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->b)); 1849566063dSJacob Faibussowitsch PetscCall(VecLoad(user->b, fd)); 1859566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->xGT)); 1869566063dSJacob Faibussowitsch PetscCall(VecLoad(user->xGT, fd)); 1879566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd)); 188c4762a1bSJed Brown 1899566063dSJacob Faibussowitsch PetscCall(MatGetSize(user->A, &user->M, &user->N)); 190c4762a1bSJed Brown 1919566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->D)); 1929566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->D, PETSC_DECIDE, PETSC_DECIDE, user->N, user->N)); 1939566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->D)); 1949566063dSJacob Faibussowitsch PetscCall(MatSetUp(user->D)); 195c4762a1bSJed Brown for (k = 0; k < user->N; k++) { 196c4762a1bSJed Brown v = 1.0; 197c4762a1bSJed Brown n = k + 1; 198*48a46eb9SPierre Jolivet if (k < user->N - 1) PetscCall(MatSetValues(user->D, 1, &k, 1, &n, &v, INSERT_VALUES)); 199c4762a1bSJed Brown v = -1.0; 2009566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->D, 1, &k, 1, &k, &v, INSERT_VALUES)); 201c4762a1bSJed Brown } 2029566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->D, MAT_FINAL_ASSEMBLY)); 2039566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->D, MAT_FINAL_ASSEMBLY)); 204c4762a1bSJed Brown 2059566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMult(user->D, user->D, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &user->DTD)); 206c4762a1bSJed Brown 2079566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Hz)); 2089566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Hz, PETSC_DECIDE, PETSC_DECIDE, user->N, user->N)); 2099566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Hz)); 2109566063dSJacob Faibussowitsch PetscCall(MatSetUp(user->Hz)); 2119566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Hz, MAT_FINAL_ASSEMBLY)); 2129566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Hz, MAT_FINAL_ASSEMBLY)); 213c4762a1bSJed Brown 2149566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &(user->x))); 2159566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &(user->workM))); 2169566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &(user->workN))); 2179566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &(user->workN2))); 2189566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->x, PETSC_DECIDE, user->N)); 2199566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->workM, PETSC_DECIDE, user->M)); 2209566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->workN, PETSC_DECIDE, user->N)); 2219566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->workN2, PETSC_DECIDE, user->N)); 2229566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->x)); 2239566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->workM)); 2249566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->workN)); 2259566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->workN2)); 226c4762a1bSJed Brown 2279566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->workN, &(user->workN3))); 2289566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->x, &(user->xlb))); 2299566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->x, &(user->xub))); 2309566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->x, &(user->c))); 2319566063dSJacob Faibussowitsch PetscCall(VecSet(user->xlb, 0.0)); 2329566063dSJacob Faibussowitsch PetscCall(VecSet(user->c, 0.0)); 2339566063dSJacob Faibussowitsch PetscCall(VecSet(user->xub, PETSC_INFINITY)); 234c4762a1bSJed Brown 2359566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMult(user->A, user->A, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &(user->ATA))); 2369566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMult(user->A, user->A, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &(user->Hx))); 2379566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMult(user->A, user->A, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &(user->HF))); 238c4762a1bSJed Brown 2399566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->ATA, MAT_FINAL_ASSEMBLY)); 2409566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->ATA, MAT_FINAL_ASSEMBLY)); 2419566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Hx, MAT_FINAL_ASSEMBLY)); 2429566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Hx, MAT_FINAL_ASSEMBLY)); 2439566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->HF, MAT_FINAL_ASSEMBLY)); 2449566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->HF, MAT_FINAL_ASSEMBLY)); 245c4762a1bSJed Brown 246c4762a1bSJed Brown user->lambda = 1.e-8; 247c4762a1bSJed Brown user->eps = 1.e-3; 248c4762a1bSJed Brown user->reg = 2; 249c4762a1bSJed Brown user->mumin = 5.e-6; 250c4762a1bSJed Brown 251d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Configure separable objection example", "tomographyADMM.c"); 2529566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-reg", "Regularization scheme for z solver (1,2)", "tomographyADMM.c", user->reg, &(user->reg), NULL)); 2539566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-lambda", "The regularization multiplier. 1 default", "tomographyADMM.c", user->lambda, &(user->lambda), NULL)); 2549566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-eps", "L1 norm epsilon padding", "tomographyADMM.c", user->eps, &(user->eps), NULL)); 2559566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mumin", "Minimum value for ADMM spectral penalty", "tomographyADMM.c", user->mumin, &(user->mumin), NULL)); 256d0609cedSBarry Smith PetscOptionsEnd(); 257c4762a1bSJed Brown PetscFunctionReturn(0); 258c4762a1bSJed Brown } 259c4762a1bSJed Brown 260c4762a1bSJed Brown /*------------------------------------------------------------*/ 261c4762a1bSJed Brown 2629371c9d4SSatish Balay PetscErrorCode DestroyContext(AppCtx *user) { 263c4762a1bSJed Brown PetscFunctionBegin; 2649566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->A)); 2659566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->ATA)); 2669566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Hx)); 2679566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Hz)); 2689566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->HF)); 2699566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->D)); 2709566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->DTD)); 2719566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->xGT)); 2729566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->xlb)); 2739566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->xub)); 2749566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->b)); 2759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->x)); 2769566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->c)); 2779566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->workN3)); 2789566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->workN2)); 2799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->workN)); 2809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->workM)); 281c4762a1bSJed Brown PetscFunctionReturn(0); 282c4762a1bSJed Brown } 283c4762a1bSJed Brown 284c4762a1bSJed Brown /*------------------------------------------------------------*/ 285c4762a1bSJed Brown 2869371c9d4SSatish Balay int main(int argc, char **argv) { 287c4762a1bSJed Brown Tao tao, misfit, reg; 288c4762a1bSJed Brown PetscReal v1, v2; 289c4762a1bSJed Brown AppCtx *user; 290c4762a1bSJed Brown PetscViewer fd; 291c4762a1bSJed Brown char resultFile[] = "tomographyResult_x"; 292c4762a1bSJed Brown 293327415f7SBarry Smith PetscFunctionBeginUser; 2949566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 2959566063dSJacob Faibussowitsch PetscCall(PetscNew(&user)); 2969566063dSJacob Faibussowitsch PetscCall(InitializeUserData(user)); 297c4762a1bSJed Brown 2989566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao)); 2999566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOADMM)); 3009566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, user->x)); 301c4762a1bSJed Brown /* f(x) + g(x) for parent tao */ 3029566063dSJacob Faibussowitsch PetscCall(TaoADMMSetSpectralPenalty(tao, 1.)); 3039566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FullObjGrad, (void *)user)); 3049566063dSJacob Faibussowitsch PetscCall(MatShift(user->HF, user->lambda)); 3059566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao, user->HF, user->HF, HessianFull, (void *)user)); 306c4762a1bSJed Brown 307c4762a1bSJed Brown /* f(x) for misfit tao */ 3089566063dSJacob Faibussowitsch PetscCall(TaoADMMSetMisfitObjectiveAndGradientRoutine(tao, MisfitObjectiveAndGradient, (void *)user)); 3099566063dSJacob Faibussowitsch PetscCall(TaoADMMSetMisfitHessianRoutine(tao, user->Hx, user->Hx, HessianMisfit, (void *)user)); 3109566063dSJacob Faibussowitsch PetscCall(TaoADMMSetMisfitHessianChangeStatus(tao, PETSC_FALSE)); 3119566063dSJacob Faibussowitsch PetscCall(TaoADMMSetMisfitConstraintJacobian(tao, user->D, user->D, NullJacobian, (void *)user)); 312c4762a1bSJed Brown 313c4762a1bSJed Brown /* g(x) for regularizer tao */ 314c4762a1bSJed Brown if (user->reg == 1) { 3159566063dSJacob Faibussowitsch PetscCall(TaoADMMSetRegularizerObjectiveAndGradientRoutine(tao, RegularizerObjectiveAndGradient1, (void *)user)); 3169566063dSJacob Faibussowitsch PetscCall(TaoADMMSetRegularizerHessianRoutine(tao, user->Hz, user->Hz, HessianReg, (void *)user)); 3179566063dSJacob Faibussowitsch PetscCall(TaoADMMSetRegHessianChangeStatus(tao, PETSC_TRUE)); 318c4762a1bSJed Brown } else if (user->reg == 2) { 3199566063dSJacob Faibussowitsch PetscCall(TaoADMMSetRegularizerObjectiveAndGradientRoutine(tao, RegularizerObjectiveAndGradient2, (void *)user)); 3209566063dSJacob Faibussowitsch PetscCall(MatShift(user->Hz, 1)); 3219566063dSJacob Faibussowitsch PetscCall(MatScale(user->Hz, user->lambda)); 3229566063dSJacob Faibussowitsch PetscCall(TaoADMMSetRegularizerHessianRoutine(tao, user->Hz, user->Hz, HessianMisfit, (void *)user)); 3239566063dSJacob Faibussowitsch PetscCall(TaoADMMSetRegHessianChangeStatus(tao, PETSC_TRUE)); 3243c859ba3SBarry Smith } else PetscCheck(user->reg == 3, PETSC_COMM_WORLD, PETSC_ERR_ARG_UNKNOWN_TYPE, "Incorrect Reg type"); /* TaoShell case */ 325c4762a1bSJed Brown 326c4762a1bSJed Brown /* Set type for the misfit solver */ 3279566063dSJacob Faibussowitsch PetscCall(TaoADMMGetMisfitSubsolver(tao, &misfit)); 3289566063dSJacob Faibussowitsch PetscCall(TaoADMMGetRegularizationSubsolver(tao, ®)); 3299566063dSJacob Faibussowitsch PetscCall(TaoSetType(misfit, TAONLS)); 330c4762a1bSJed Brown if (user->reg == 3) { 3319566063dSJacob Faibussowitsch PetscCall(TaoSetType(reg, TAOSHELL)); 3329566063dSJacob Faibussowitsch PetscCall(TaoShellSetContext(reg, (void *)user)); 3339566063dSJacob Faibussowitsch PetscCall(TaoShellSetSolve(reg, TaoShellSolve_SoftThreshold)); 334c4762a1bSJed Brown } else { 3359566063dSJacob Faibussowitsch PetscCall(TaoSetType(reg, TAONLS)); 336c4762a1bSJed Brown } 3379566063dSJacob Faibussowitsch PetscCall(TaoSetVariableBounds(misfit, user->xlb, user->xub)); 338c4762a1bSJed Brown 339c4762a1bSJed Brown /* Soft Thresholding solves the ADMM problem with the L1 regularizer lambda*||z||_1 and the x-z=0 constraint */ 3409566063dSJacob Faibussowitsch PetscCall(TaoADMMSetRegularizerCoefficient(tao, user->lambda)); 3419566063dSJacob Faibussowitsch PetscCall(TaoADMMSetRegularizerConstraintJacobian(tao, NULL, NULL, NullJacobian, (void *)user)); 3429566063dSJacob Faibussowitsch PetscCall(TaoADMMSetMinimumSpectralPenalty(tao, user->mumin)); 343c4762a1bSJed Brown 3449566063dSJacob Faibussowitsch PetscCall(TaoADMMSetConstraintVectorRHS(tao, user->c)); 3459566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 3469566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 347c4762a1bSJed Brown 348c4762a1bSJed Brown /* Save x (reconstruction of object) vector to a binary file, which maybe read from Matlab and convert to a 2D image for comparison. */ 3499566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, resultFile, FILE_MODE_WRITE, &fd)); 3509566063dSJacob Faibussowitsch PetscCall(VecView(user->x, fd)); 3519566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd)); 352c4762a1bSJed Brown 353c4762a1bSJed Brown /* compute the error */ 3549566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->x, -1, user->xGT)); 3559566063dSJacob Faibussowitsch PetscCall(VecNorm(user->x, NORM_2, &v1)); 3569566063dSJacob Faibussowitsch PetscCall(VecNorm(user->xGT, NORM_2, &v2)); 3579566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "relative reconstruction error: ||x-xGT||/||xGT|| = %6.4e.\n", (double)(v1 / v2))); 358c4762a1bSJed Brown 359c4762a1bSJed Brown /* Free TAO data structures */ 3609566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 3619566063dSJacob Faibussowitsch PetscCall(DestroyContext(user)); 3629566063dSJacob Faibussowitsch PetscCall(PetscFree(user)); 3639566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 364b122ec5aSJacob Faibussowitsch return 0; 365c4762a1bSJed Brown } 366c4762a1bSJed Brown 367c4762a1bSJed Brown /*TEST 368c4762a1bSJed Brown 369c4762a1bSJed Brown build: 370dfd57a17SPierre Jolivet requires: !complex !single !__float128 !defined(PETSC_USE_64BIT_INDICES) 371c4762a1bSJed Brown 372c4762a1bSJed Brown test: 373c4762a1bSJed Brown suffix: 1 374c4762a1bSJed Brown localrunfiles: tomographyData_A_b_xGT 375c4762a1bSJed Brown args: -lambda 1.e-8 -tao_monitor -tao_type nls -tao_nls_pc_type icc 376c4762a1bSJed Brown 377c4762a1bSJed Brown test: 378c4762a1bSJed Brown suffix: 2 379c4762a1bSJed Brown localrunfiles: tomographyData_A_b_xGT 380c4762a1bSJed Brown args: -reg 2 -lambda 1.e-8 -tao_admm_dual_update update_basic -tao_admm_regularizer_type regularizer_user -tao_max_it 20 -tao_monitor -tao_admm_tolerance_update_factor 1.e-8 -misfit_tao_nls_pc_type icc -misfit_tao_monitor -reg_tao_monitor 381c4762a1bSJed Brown 382c4762a1bSJed Brown test: 383c4762a1bSJed Brown suffix: 3 384c4762a1bSJed Brown localrunfiles: tomographyData_A_b_xGT 385c4762a1bSJed Brown args: -lambda 1.e-8 -tao_admm_dual_update update_basic -tao_admm_regularizer_type regularizer_soft_thresh -tao_max_it 20 -tao_monitor -tao_admm_tolerance_update_factor 1.e-8 -misfit_tao_nls_pc_type icc -misfit_tao_monitor 386c4762a1bSJed Brown 387c4762a1bSJed Brown test: 388c4762a1bSJed Brown suffix: 4 389c4762a1bSJed Brown localrunfiles: tomographyData_A_b_xGT 390c4762a1bSJed Brown args: -lambda 1.e-8 -tao_admm_dual_update update_adaptive -tao_admm_regularizer_type regularizer_soft_thresh -tao_max_it 20 -tao_monitor -misfit_tao_monitor -misfit_tao_nls_pc_type icc 391c4762a1bSJed Brown 392c4762a1bSJed Brown test: 393c4762a1bSJed Brown suffix: 5 394c4762a1bSJed Brown localrunfiles: tomographyData_A_b_xGT 395c4762a1bSJed Brown args: -reg 2 -lambda 1.e-8 -tao_admm_dual_update update_adaptive -tao_admm_regularizer_type regularizer_user -tao_max_it 20 -tao_monitor -tao_admm_tolerance_update_factor 1.e-8 -misfit_tao_monitor -reg_tao_monitor -misfit_tao_nls_pc_type icc 396c4762a1bSJed Brown 397c4762a1bSJed Brown test: 398c4762a1bSJed Brown suffix: 6 399c4762a1bSJed Brown localrunfiles: tomographyData_A_b_xGT 400c4762a1bSJed Brown args: -reg 3 -lambda 1.e-8 -tao_admm_dual_update update_adaptive -tao_admm_regularizer_type regularizer_user -tao_max_it 20 -tao_monitor -tao_admm_tolerance_update_factor 1.e-8 -misfit_tao_monitor -reg_tao_monitor -misfit_tao_nls_pc_type icc 401c4762a1bSJed Brown 402c4762a1bSJed Brown TEST*/ 403