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 26d71ae5a4SJacob Faibussowitsch PetscErrorCode NullJacobian(Tao tao, Vec X, Mat J, Mat Jpre, void *ptr) 27d71ae5a4SJacob Faibussowitsch { 28c4762a1bSJed Brown PetscFunctionBegin; 293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 30c4762a1bSJed Brown } 31c4762a1bSJed Brown 32c4762a1bSJed Brown /*------------------------------------------------------------*/ 33c4762a1bSJed Brown 34d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoShellSolve_SoftThreshold(Tao tao) 35d71ae5a4SJacob Faibussowitsch { 36c4762a1bSJed Brown PetscReal lambda, mu; 37c4762a1bSJed Brown AppCtx *user; 38c4762a1bSJed Brown Vec out, work, y, x; 39c4762a1bSJed Brown Tao admm_tao, misfit; 40c4762a1bSJed Brown 41c4762a1bSJed Brown PetscFunctionBegin; 42c4762a1bSJed Brown user = NULL; 43c4762a1bSJed Brown mu = 0; 449566063dSJacob Faibussowitsch PetscCall(TaoGetADMMParentTao(tao, &admm_tao)); 459566063dSJacob Faibussowitsch PetscCall(TaoADMMGetMisfitSubsolver(admm_tao, &misfit)); 469566063dSJacob Faibussowitsch PetscCall(TaoADMMGetSpectralPenalty(admm_tao, &mu)); 479566063dSJacob Faibussowitsch PetscCall(TaoShellGetContext(tao, &user)); 48b4623fecSHansol Suh PetscCall(TaoADMMGetRegularizerCoefficient(admm_tao, &lambda)); 49c4762a1bSJed Brown 50c4762a1bSJed Brown work = user->workN; 519566063dSJacob Faibussowitsch PetscCall(TaoGetSolution(tao, &out)); 529566063dSJacob Faibussowitsch PetscCall(TaoGetSolution(misfit, &x)); 539566063dSJacob Faibussowitsch PetscCall(TaoADMMGetDualVector(admm_tao, &y)); 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* Dx + y/mu */ 569566063dSJacob Faibussowitsch PetscCall(MatMult(user->D, x, work)); 579566063dSJacob Faibussowitsch PetscCall(VecAXPY(work, 1 / mu, y)); 58c4762a1bSJed Brown 59c4762a1bSJed Brown /* soft thresholding */ 609566063dSJacob Faibussowitsch PetscCall(TaoSoftThreshold(work, -lambda / mu, lambda / mu, out)); 613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 62c4762a1bSJed Brown } 63c4762a1bSJed Brown 64c4762a1bSJed Brown /*------------------------------------------------------------*/ 65c4762a1bSJed Brown 66d71ae5a4SJacob Faibussowitsch PetscErrorCode MisfitObjectiveAndGradient(Tao tao, Vec X, PetscReal *f, Vec g, void *ptr) 67d71ae5a4SJacob Faibussowitsch { 68c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 69c4762a1bSJed Brown 70c4762a1bSJed Brown PetscFunctionBegin; 71c4762a1bSJed Brown /* Objective 0.5*||Ax-b||_2^2 */ 729566063dSJacob Faibussowitsch PetscCall(MatMult(user->A, X, user->workM)); 739566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->workM, -1, user->b)); 749566063dSJacob Faibussowitsch PetscCall(VecDot(user->workM, user->workM, f)); 75c4762a1bSJed Brown *f *= 0.5; 76c4762a1bSJed Brown /* Gradient. ATAx-ATb */ 779566063dSJacob Faibussowitsch PetscCall(MatMult(user->ATA, X, user->workN)); 789566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->A, user->b, user->workN2)); 799566063dSJacob Faibussowitsch PetscCall(VecWAXPY(g, -1., user->workN2, user->workN)); 803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 81c4762a1bSJed Brown } 82c4762a1bSJed Brown 83c4762a1bSJed Brown /*------------------------------------------------------------*/ 84c4762a1bSJed Brown 85d71ae5a4SJacob Faibussowitsch PetscErrorCode RegularizerObjectiveAndGradient1(Tao tao, Vec X, PetscReal *f_reg, Vec G_reg, void *ptr) 86d71ae5a4SJacob Faibussowitsch { 87c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 88b4623fecSHansol Suh PetscReal lambda; 89b4623fecSHansol Suh Tao admm_tao; 90c4762a1bSJed Brown 91c4762a1bSJed Brown PetscFunctionBegin; 92c4762a1bSJed Brown /* compute regularizer objective 93c4762a1bSJed Brown * f = f + lambda*sum(sqrt(y.^2+epsilon^2) - epsilon), where y = D*x */ 949566063dSJacob Faibussowitsch PetscCall(VecCopy(X, user->workN2)); 959566063dSJacob Faibussowitsch PetscCall(VecPow(user->workN2, 2.)); 969566063dSJacob Faibussowitsch PetscCall(VecShift(user->workN2, user->eps * user->eps)); 979566063dSJacob Faibussowitsch PetscCall(VecSqrtAbs(user->workN2)); 989566063dSJacob Faibussowitsch PetscCall(VecCopy(user->workN2, user->workN3)); 999566063dSJacob Faibussowitsch PetscCall(VecShift(user->workN2, -user->eps)); 1009566063dSJacob Faibussowitsch PetscCall(VecSum(user->workN2, f_reg)); 101b4623fecSHansol Suh PetscCall(TaoGetADMMParentTao(tao, &admm_tao)); 102b4623fecSHansol Suh PetscCall(TaoADMMGetRegularizerCoefficient(admm_tao, &lambda)); 103b4623fecSHansol Suh *f_reg *= lambda; 104c4762a1bSJed Brown /* compute regularizer gradient = lambda*x */ 1059566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(G_reg, X, user->workN3)); 106b4623fecSHansol Suh PetscCall(VecScale(G_reg, lambda)); 1073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 108c4762a1bSJed Brown } 109c4762a1bSJed Brown 110c4762a1bSJed Brown /*------------------------------------------------------------*/ 111c4762a1bSJed Brown 112d71ae5a4SJacob Faibussowitsch PetscErrorCode RegularizerObjectiveAndGradient2(Tao tao, Vec X, PetscReal *f_reg, Vec G_reg, void *ptr) 113d71ae5a4SJacob Faibussowitsch { 114b4623fecSHansol Suh PetscReal temp, lambda; 115b4623fecSHansol Suh Tao admm_tao; 116c4762a1bSJed Brown 117c4762a1bSJed Brown PetscFunctionBegin; 118c4762a1bSJed Brown /* compute regularizer objective = lambda*|z|_2^2 */ 1199566063dSJacob Faibussowitsch PetscCall(VecDot(X, X, &temp)); 120b4623fecSHansol Suh PetscCall(TaoGetADMMParentTao(tao, &admm_tao)); 121b4623fecSHansol Suh PetscCall(TaoADMMGetRegularizerCoefficient(admm_tao, &lambda)); 122b4623fecSHansol Suh *f_reg = 0.5 * lambda * temp; 123c4762a1bSJed Brown /* compute regularizer gradient = lambda*z */ 1249566063dSJacob Faibussowitsch PetscCall(VecCopy(X, G_reg)); 125b4623fecSHansol Suh PetscCall(VecScale(G_reg, lambda)); 1263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 127c4762a1bSJed Brown } 128c4762a1bSJed Brown 129c4762a1bSJed Brown /*------------------------------------------------------------*/ 130c4762a1bSJed Brown 131d71ae5a4SJacob Faibussowitsch static PetscErrorCode HessianMisfit(Tao tao, Vec x, Mat H, Mat Hpre, void *ptr) 132d71ae5a4SJacob Faibussowitsch { 133c4762a1bSJed Brown PetscFunctionBegin; 1343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 135c4762a1bSJed Brown } 136c4762a1bSJed Brown 137c4762a1bSJed Brown /*------------------------------------------------------------*/ 138c4762a1bSJed Brown 139d71ae5a4SJacob Faibussowitsch static PetscErrorCode HessianReg(Tao tao, Vec x, Mat H, Mat Hpre, void *ptr) 140d71ae5a4SJacob Faibussowitsch { 141c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 142c4762a1bSJed Brown 143c4762a1bSJed Brown PetscFunctionBegin; 1449566063dSJacob Faibussowitsch PetscCall(MatMult(user->D, x, user->workN)); 1459566063dSJacob Faibussowitsch PetscCall(VecPow(user->workN2, 2.)); 1469566063dSJacob Faibussowitsch PetscCall(VecShift(user->workN2, user->eps * user->eps)); 1479566063dSJacob Faibussowitsch PetscCall(VecSqrtAbs(user->workN2)); 1489566063dSJacob Faibussowitsch PetscCall(VecShift(user->workN2, -user->eps)); 1499566063dSJacob Faibussowitsch PetscCall(VecReciprocal(user->workN2)); 1509566063dSJacob Faibussowitsch PetscCall(VecScale(user->workN2, user->eps * user->eps)); 1519566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(H, user->workN2, INSERT_VALUES)); 1523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 153c4762a1bSJed Brown } 154c4762a1bSJed Brown 155c4762a1bSJed Brown /*------------------------------------------------------------*/ 156c4762a1bSJed Brown 157d71ae5a4SJacob Faibussowitsch PetscErrorCode FullObjGrad(Tao tao, Vec X, PetscReal *f, Vec g, void *ptr) 158d71ae5a4SJacob Faibussowitsch { 159c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 160b4623fecSHansol Suh PetscReal f_reg, lambda; 161b4623fecSHansol Suh PetscBool is_admm; 162c4762a1bSJed Brown 163c4762a1bSJed Brown PetscFunctionBegin; 164b4623fecSHansol Suh /* Objective 0.5*||Ax-b||_2^2 + lambda*||x||_{1,2}^2*/ 1659566063dSJacob Faibussowitsch PetscCall(MatMult(user->A, X, user->workM)); 1669566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->workM, -1, user->b)); 1679566063dSJacob Faibussowitsch PetscCall(VecDot(user->workM, user->workM, f)); 168b4623fecSHansol Suh if (user->reg == 1) { 169b4623fecSHansol Suh PetscCall(VecNorm(X, NORM_1, &f_reg)); 170b4623fecSHansol Suh } else { 1719566063dSJacob Faibussowitsch PetscCall(VecNorm(X, NORM_2, &f_reg)); 172b4623fecSHansol Suh } 173b4623fecSHansol Suh PetscCall(PetscObjectTypeCompare((PetscObject)tao, TAOADMM, &is_admm)); 174b4623fecSHansol Suh if (is_admm) { 175b4623fecSHansol Suh PetscCall(TaoADMMGetRegularizerCoefficient(tao, &lambda)); 176b4623fecSHansol Suh } else { 177b4623fecSHansol Suh lambda = user->lambda; 178b4623fecSHansol Suh } 179c4762a1bSJed Brown *f *= 0.5; 180b4623fecSHansol Suh *f += lambda * f_reg * f_reg; 181c4762a1bSJed Brown /* Gradient. ATAx-ATb + 2*lambda*x */ 1829566063dSJacob Faibussowitsch PetscCall(MatMult(user->ATA, X, user->workN)); 1839566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->A, user->b, user->workN2)); 1849566063dSJacob Faibussowitsch PetscCall(VecWAXPY(g, -1., user->workN2, user->workN)); 185b4623fecSHansol Suh PetscCall(VecAXPY(g, 2 * lambda, X)); 1863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 187c4762a1bSJed Brown } 188c4762a1bSJed Brown /*------------------------------------------------------------*/ 189c4762a1bSJed Brown 190d71ae5a4SJacob Faibussowitsch static PetscErrorCode HessianFull(Tao tao, Vec x, Mat H, Mat Hpre, void *ptr) 191d71ae5a4SJacob Faibussowitsch { 192c4762a1bSJed Brown PetscFunctionBegin; 1933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 194c4762a1bSJed Brown } 195c4762a1bSJed Brown /*------------------------------------------------------------*/ 196c4762a1bSJed Brown 197d71ae5a4SJacob Faibussowitsch PetscErrorCode InitializeUserData(AppCtx *user) 198d71ae5a4SJacob Faibussowitsch { 199c4762a1bSJed 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". */ 200c4762a1bSJed Brown PetscViewer fd; /* used to load data from file */ 201c4762a1bSJed Brown PetscInt k, n; 202c4762a1bSJed Brown PetscScalar v; 203c4762a1bSJed Brown 204d0609cedSBarry Smith PetscFunctionBegin; 205c4762a1bSJed Brown /* Load the A matrix, b vector, and xGT vector from a binary file. */ 2069566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, dataFile, FILE_MODE_READ, &fd)); 2079566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->A)); 2089566063dSJacob Faibussowitsch PetscCall(MatSetType(user->A, MATAIJ)); 2099566063dSJacob Faibussowitsch PetscCall(MatLoad(user->A, fd)); 2109566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->b)); 2119566063dSJacob Faibussowitsch PetscCall(VecLoad(user->b, fd)); 2129566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user->xGT)); 2139566063dSJacob Faibussowitsch PetscCall(VecLoad(user->xGT, fd)); 2149566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd)); 215c4762a1bSJed Brown 2169566063dSJacob Faibussowitsch PetscCall(MatGetSize(user->A, &user->M, &user->N)); 217c4762a1bSJed Brown 2189566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->D)); 2199566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->D, PETSC_DECIDE, PETSC_DECIDE, user->N, user->N)); 2209566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->D)); 2219566063dSJacob Faibussowitsch PetscCall(MatSetUp(user->D)); 222c4762a1bSJed Brown for (k = 0; k < user->N; k++) { 223c4762a1bSJed Brown v = 1.0; 224c4762a1bSJed Brown n = k + 1; 22548a46eb9SPierre Jolivet if (k < user->N - 1) PetscCall(MatSetValues(user->D, 1, &k, 1, &n, &v, INSERT_VALUES)); 226c4762a1bSJed Brown v = -1.0; 2279566063dSJacob Faibussowitsch PetscCall(MatSetValues(user->D, 1, &k, 1, &k, &v, INSERT_VALUES)); 228c4762a1bSJed Brown } 2299566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->D, MAT_FINAL_ASSEMBLY)); 2309566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->D, MAT_FINAL_ASSEMBLY)); 231c4762a1bSJed Brown 2329566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMult(user->D, user->D, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &user->DTD)); 233c4762a1bSJed Brown 2349566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Hz)); 2359566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user->Hz, PETSC_DECIDE, PETSC_DECIDE, user->N, user->N)); 2369566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user->Hz)); 2379566063dSJacob Faibussowitsch PetscCall(MatSetUp(user->Hz)); 2389566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Hz, MAT_FINAL_ASSEMBLY)); 2399566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Hz, MAT_FINAL_ASSEMBLY)); 240c4762a1bSJed Brown 241*f4f49eeaSPierre Jolivet PetscCall(VecCreate(PETSC_COMM_WORLD, &user->x)); 242*f4f49eeaSPierre Jolivet PetscCall(VecCreate(PETSC_COMM_WORLD, &user->workM)); 243*f4f49eeaSPierre Jolivet PetscCall(VecCreate(PETSC_COMM_WORLD, &user->workN)); 244*f4f49eeaSPierre Jolivet PetscCall(VecCreate(PETSC_COMM_WORLD, &user->workN2)); 2459566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->x, PETSC_DECIDE, user->N)); 2469566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->workM, PETSC_DECIDE, user->M)); 2479566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->workN, PETSC_DECIDE, user->N)); 2489566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user->workN2, PETSC_DECIDE, user->N)); 2499566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->x)); 2509566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->workM)); 2519566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->workN)); 2529566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(user->workN2)); 253c4762a1bSJed Brown 254*f4f49eeaSPierre Jolivet PetscCall(VecDuplicate(user->workN, &user->workN3)); 255*f4f49eeaSPierre Jolivet PetscCall(VecDuplicate(user->x, &user->xlb)); 256*f4f49eeaSPierre Jolivet PetscCall(VecDuplicate(user->x, &user->xub)); 257*f4f49eeaSPierre Jolivet PetscCall(VecDuplicate(user->x, &user->c)); 2589566063dSJacob Faibussowitsch PetscCall(VecSet(user->xlb, 0.0)); 2599566063dSJacob Faibussowitsch PetscCall(VecSet(user->c, 0.0)); 2609566063dSJacob Faibussowitsch PetscCall(VecSet(user->xub, PETSC_INFINITY)); 261c4762a1bSJed Brown 262*f4f49eeaSPierre Jolivet PetscCall(MatTransposeMatMult(user->A, user->A, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &user->ATA)); 263*f4f49eeaSPierre Jolivet PetscCall(MatTransposeMatMult(user->A, user->A, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &user->Hx)); 264*f4f49eeaSPierre Jolivet PetscCall(MatTransposeMatMult(user->A, user->A, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &user->HF)); 265c4762a1bSJed Brown 2669566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->ATA, MAT_FINAL_ASSEMBLY)); 2679566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->ATA, MAT_FINAL_ASSEMBLY)); 2689566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->Hx, MAT_FINAL_ASSEMBLY)); 2699566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->Hx, MAT_FINAL_ASSEMBLY)); 2709566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user->HF, MAT_FINAL_ASSEMBLY)); 2719566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user->HF, MAT_FINAL_ASSEMBLY)); 272c4762a1bSJed Brown 273c4762a1bSJed Brown user->lambda = 1.e-8; 274c4762a1bSJed Brown user->eps = 1.e-3; 275c4762a1bSJed Brown user->reg = 2; 276c4762a1bSJed Brown user->mumin = 5.e-6; 277c4762a1bSJed Brown 278d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Configure separable objection example", "tomographyADMM.c"); 279*f4f49eeaSPierre Jolivet PetscCall(PetscOptionsInt("-reg", "Regularization scheme for z solver (1,2)", "tomographyADMM.c", user->reg, &user->reg, NULL)); 280*f4f49eeaSPierre Jolivet PetscCall(PetscOptionsReal("-lambda", "The regularization multiplier. 1 default", "tomographyADMM.c", user->lambda, &user->lambda, NULL)); 281*f4f49eeaSPierre Jolivet PetscCall(PetscOptionsReal("-eps", "L1 norm epsilon padding", "tomographyADMM.c", user->eps, &user->eps, NULL)); 282*f4f49eeaSPierre Jolivet PetscCall(PetscOptionsReal("-mumin", "Minimum value for ADMM spectral penalty", "tomographyADMM.c", user->mumin, &user->mumin, NULL)); 283d0609cedSBarry Smith PetscOptionsEnd(); 2843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 285c4762a1bSJed Brown } 286c4762a1bSJed Brown 287c4762a1bSJed Brown /*------------------------------------------------------------*/ 288c4762a1bSJed Brown 289d71ae5a4SJacob Faibussowitsch PetscErrorCode DestroyContext(AppCtx *user) 290d71ae5a4SJacob Faibussowitsch { 291c4762a1bSJed Brown PetscFunctionBegin; 2929566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->A)); 2939566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->ATA)); 2949566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Hx)); 2959566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->Hz)); 2969566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->HF)); 2979566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->D)); 2989566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->DTD)); 2999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->xGT)); 3009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->xlb)); 3019566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->xub)); 3029566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->b)); 3039566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->x)); 3049566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->c)); 3059566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->workN3)); 3069566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->workN2)); 3079566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->workN)); 3089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->workM)); 3093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 310c4762a1bSJed Brown } 311c4762a1bSJed Brown 312c4762a1bSJed Brown /*------------------------------------------------------------*/ 313c4762a1bSJed Brown 314d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 315d71ae5a4SJacob Faibussowitsch { 316c4762a1bSJed Brown Tao tao, misfit, reg; 317c4762a1bSJed Brown PetscReal v1, v2; 318c4762a1bSJed Brown AppCtx *user; 319c4762a1bSJed Brown PetscViewer fd; 320c4762a1bSJed Brown char resultFile[] = "tomographyResult_x"; 321c4762a1bSJed Brown 322327415f7SBarry Smith PetscFunctionBeginUser; 3239566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 3249566063dSJacob Faibussowitsch PetscCall(PetscNew(&user)); 3259566063dSJacob Faibussowitsch PetscCall(InitializeUserData(user)); 326c4762a1bSJed Brown 3279566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao)); 3289566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOADMM)); 3299566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, user->x)); 330c4762a1bSJed Brown /* f(x) + g(x) for parent tao */ 3319566063dSJacob Faibussowitsch PetscCall(TaoADMMSetSpectralPenalty(tao, 1.)); 3329566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FullObjGrad, (void *)user)); 3339566063dSJacob Faibussowitsch PetscCall(MatShift(user->HF, user->lambda)); 3349566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao, user->HF, user->HF, HessianFull, (void *)user)); 335c4762a1bSJed Brown 336c4762a1bSJed Brown /* f(x) for misfit tao */ 3379566063dSJacob Faibussowitsch PetscCall(TaoADMMSetMisfitObjectiveAndGradientRoutine(tao, MisfitObjectiveAndGradient, (void *)user)); 3389566063dSJacob Faibussowitsch PetscCall(TaoADMMSetMisfitHessianRoutine(tao, user->Hx, user->Hx, HessianMisfit, (void *)user)); 3399566063dSJacob Faibussowitsch PetscCall(TaoADMMSetMisfitHessianChangeStatus(tao, PETSC_FALSE)); 3409566063dSJacob Faibussowitsch PetscCall(TaoADMMSetMisfitConstraintJacobian(tao, user->D, user->D, NullJacobian, (void *)user)); 341c4762a1bSJed Brown 342c4762a1bSJed Brown /* g(x) for regularizer tao */ 343c4762a1bSJed Brown if (user->reg == 1) { 3449566063dSJacob Faibussowitsch PetscCall(TaoADMMSetRegularizerObjectiveAndGradientRoutine(tao, RegularizerObjectiveAndGradient1, (void *)user)); 3459566063dSJacob Faibussowitsch PetscCall(TaoADMMSetRegularizerHessianRoutine(tao, user->Hz, user->Hz, HessianReg, (void *)user)); 3469566063dSJacob Faibussowitsch PetscCall(TaoADMMSetRegHessianChangeStatus(tao, PETSC_TRUE)); 347c4762a1bSJed Brown } else if (user->reg == 2) { 3489566063dSJacob Faibussowitsch PetscCall(TaoADMMSetRegularizerObjectiveAndGradientRoutine(tao, RegularizerObjectiveAndGradient2, (void *)user)); 3499566063dSJacob Faibussowitsch PetscCall(MatShift(user->Hz, 1)); 3509566063dSJacob Faibussowitsch PetscCall(MatScale(user->Hz, user->lambda)); 3519566063dSJacob Faibussowitsch PetscCall(TaoADMMSetRegularizerHessianRoutine(tao, user->Hz, user->Hz, HessianMisfit, (void *)user)); 3529566063dSJacob Faibussowitsch PetscCall(TaoADMMSetRegHessianChangeStatus(tao, PETSC_TRUE)); 3533c859ba3SBarry Smith } else PetscCheck(user->reg == 3, PETSC_COMM_WORLD, PETSC_ERR_ARG_UNKNOWN_TYPE, "Incorrect Reg type"); /* TaoShell case */ 354c4762a1bSJed Brown 355c4762a1bSJed Brown /* Set type for the misfit solver */ 3569566063dSJacob Faibussowitsch PetscCall(TaoADMMGetMisfitSubsolver(tao, &misfit)); 3579566063dSJacob Faibussowitsch PetscCall(TaoADMMGetRegularizationSubsolver(tao, ®)); 3589566063dSJacob Faibussowitsch PetscCall(TaoSetType(misfit, TAONLS)); 359c4762a1bSJed Brown if (user->reg == 3) { 3609566063dSJacob Faibussowitsch PetscCall(TaoSetType(reg, TAOSHELL)); 3619566063dSJacob Faibussowitsch PetscCall(TaoShellSetContext(reg, (void *)user)); 3629566063dSJacob Faibussowitsch PetscCall(TaoShellSetSolve(reg, TaoShellSolve_SoftThreshold)); 363c4762a1bSJed Brown } else { 3649566063dSJacob Faibussowitsch PetscCall(TaoSetType(reg, TAONLS)); 365c4762a1bSJed Brown } 3669566063dSJacob Faibussowitsch PetscCall(TaoSetVariableBounds(misfit, user->xlb, user->xub)); 367c4762a1bSJed Brown 368c4762a1bSJed Brown /* Soft Thresholding solves the ADMM problem with the L1 regularizer lambda*||z||_1 and the x-z=0 constraint */ 3699566063dSJacob Faibussowitsch PetscCall(TaoADMMSetRegularizerCoefficient(tao, user->lambda)); 3709566063dSJacob Faibussowitsch PetscCall(TaoADMMSetRegularizerConstraintJacobian(tao, NULL, NULL, NullJacobian, (void *)user)); 3719566063dSJacob Faibussowitsch PetscCall(TaoADMMSetMinimumSpectralPenalty(tao, user->mumin)); 372c4762a1bSJed Brown 3739566063dSJacob Faibussowitsch PetscCall(TaoADMMSetConstraintVectorRHS(tao, user->c)); 3749566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 3759566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 376c4762a1bSJed Brown 37721afe8ebSBarry Smith /* Save x (reconstruction of object) vector to a binary file, which maybe read from MATLAB and convert to a 2D image for comparison. */ 3789566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, resultFile, FILE_MODE_WRITE, &fd)); 3799566063dSJacob Faibussowitsch PetscCall(VecView(user->x, fd)); 3809566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd)); 381c4762a1bSJed Brown 382c4762a1bSJed Brown /* compute the error */ 3839566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->x, -1, user->xGT)); 3849566063dSJacob Faibussowitsch PetscCall(VecNorm(user->x, NORM_2, &v1)); 3859566063dSJacob Faibussowitsch PetscCall(VecNorm(user->xGT, NORM_2, &v2)); 3869566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "relative reconstruction error: ||x-xGT||/||xGT|| = %6.4e.\n", (double)(v1 / v2))); 387c4762a1bSJed Brown 388c4762a1bSJed Brown /* Free TAO data structures */ 3899566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 3909566063dSJacob Faibussowitsch PetscCall(DestroyContext(user)); 3919566063dSJacob Faibussowitsch PetscCall(PetscFree(user)); 3929566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 393b122ec5aSJacob Faibussowitsch return 0; 394c4762a1bSJed Brown } 395c4762a1bSJed Brown 396c4762a1bSJed Brown /*TEST 397c4762a1bSJed Brown 398c4762a1bSJed Brown build: 399dfd57a17SPierre Jolivet requires: !complex !single !__float128 !defined(PETSC_USE_64BIT_INDICES) 400c4762a1bSJed Brown 401c4762a1bSJed Brown test: 402c4762a1bSJed Brown suffix: 1 403c4762a1bSJed Brown localrunfiles: tomographyData_A_b_xGT 404c4762a1bSJed Brown args: -lambda 1.e-8 -tao_monitor -tao_type nls -tao_nls_pc_type icc 405c4762a1bSJed Brown 406c4762a1bSJed Brown test: 407c4762a1bSJed Brown suffix: 2 408c4762a1bSJed Brown localrunfiles: tomographyData_A_b_xGT 409c4762a1bSJed 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 410c4762a1bSJed Brown 411c4762a1bSJed Brown test: 412c4762a1bSJed Brown suffix: 3 413c4762a1bSJed Brown localrunfiles: tomographyData_A_b_xGT 414c4762a1bSJed 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 415c4762a1bSJed Brown 416c4762a1bSJed Brown test: 417c4762a1bSJed Brown suffix: 4 418c4762a1bSJed Brown localrunfiles: tomographyData_A_b_xGT 419c4762a1bSJed 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 420c4762a1bSJed Brown 421c4762a1bSJed Brown test: 422c4762a1bSJed Brown suffix: 5 423c4762a1bSJed Brown localrunfiles: tomographyData_A_b_xGT 424c4762a1bSJed 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 425c4762a1bSJed Brown 426c4762a1bSJed Brown test: 427c4762a1bSJed Brown suffix: 6 428c4762a1bSJed Brown localrunfiles: tomographyData_A_b_xGT 429c4762a1bSJed 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 430c4762a1bSJed Brown 431c4762a1bSJed Brown TEST*/ 432