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 26c4762a1bSJed Brown PetscErrorCode NullJacobian(Tao tao,Vec X,Mat J,Mat Jpre,void *ptr) 27c4762a1bSJed Brown { 28c4762a1bSJed Brown PetscFunctionBegin; 29c4762a1bSJed Brown PetscFunctionReturn(0); 30c4762a1bSJed Brown } 31c4762a1bSJed Brown 32c4762a1bSJed Brown /*------------------------------------------------------------*/ 33c4762a1bSJed Brown 34c4762a1bSJed Brown static PetscErrorCode TaoShellSolve_SoftThreshold(Tao tao) 35c4762a1bSJed Brown { 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; 445f80ce2aSJacob Faibussowitsch CHKERRQ(TaoGetADMMParentTao(tao,&admm_tao)); 455f80ce2aSJacob Faibussowitsch CHKERRQ(TaoADMMGetMisfitSubsolver(admm_tao, &misfit)); 465f80ce2aSJacob Faibussowitsch CHKERRQ(TaoADMMGetSpectralPenalty(admm_tao,&mu)); 475f80ce2aSJacob Faibussowitsch CHKERRQ(TaoShellGetContext(tao,&user)); 48c4762a1bSJed Brown 49c4762a1bSJed Brown lambda = user->lambda; 50c4762a1bSJed Brown work = user->workN; 515f80ce2aSJacob Faibussowitsch CHKERRQ(TaoGetSolution(tao, &out)); 525f80ce2aSJacob Faibussowitsch CHKERRQ(TaoGetSolution(misfit, &x)); 535f80ce2aSJacob Faibussowitsch CHKERRQ(TaoADMMGetDualVector(admm_tao, &y)); 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* Dx + y/mu */ 565f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(user->D,x,work)); 575f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(work,1/mu,y)); 58c4762a1bSJed Brown 59c4762a1bSJed Brown /* soft thresholding */ 605f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSoftThreshold(work, -lambda/mu, lambda/mu, out)); 61c4762a1bSJed Brown PetscFunctionReturn(0); 62c4762a1bSJed Brown } 63c4762a1bSJed Brown 64c4762a1bSJed Brown /*------------------------------------------------------------*/ 65c4762a1bSJed Brown 66c4762a1bSJed Brown PetscErrorCode MisfitObjectiveAndGradient(Tao tao,Vec X,PetscReal *f,Vec g,void *ptr) 67c4762a1bSJed Brown { 68c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 69c4762a1bSJed Brown 70c4762a1bSJed Brown PetscFunctionBegin; 71c4762a1bSJed Brown /* Objective 0.5*||Ax-b||_2^2 */ 725f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(user->A,X,user->workM)); 735f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(user->workM,-1,user->b)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(VecDot(user->workM,user->workM,f)); 75c4762a1bSJed Brown *f *= 0.5; 76c4762a1bSJed Brown /* Gradient. ATAx-ATb */ 775f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(user->ATA,X,user->workN)); 785f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(user->A,user->b,user->workN2)); 795f80ce2aSJacob Faibussowitsch CHKERRQ(VecWAXPY(g,-1.,user->workN2,user->workN)); 80c4762a1bSJed Brown PetscFunctionReturn(0); 81c4762a1bSJed Brown } 82c4762a1bSJed Brown 83c4762a1bSJed Brown /*------------------------------------------------------------*/ 84c4762a1bSJed Brown 85c4762a1bSJed Brown PetscErrorCode RegularizerObjectiveAndGradient1(Tao tao,Vec X,PetscReal *f_reg,Vec G_reg,void *ptr) 86c4762a1bSJed Brown { 87c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 88c4762a1bSJed Brown 89c4762a1bSJed Brown PetscFunctionBegin; 90c4762a1bSJed Brown /* compute regularizer objective 91c4762a1bSJed Brown * f = f + lambda*sum(sqrt(y.^2+epsilon^2) - epsilon), where y = D*x */ 925f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(X,user->workN2)); 935f80ce2aSJacob Faibussowitsch CHKERRQ(VecPow(user->workN2,2.)); 945f80ce2aSJacob Faibussowitsch CHKERRQ(VecShift(user->workN2,user->eps*user->eps)); 955f80ce2aSJacob Faibussowitsch CHKERRQ(VecSqrtAbs(user->workN2)); 965f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(user->workN2, user->workN3)); 975f80ce2aSJacob Faibussowitsch CHKERRQ(VecShift(user->workN2,-user->eps)); 985f80ce2aSJacob Faibussowitsch CHKERRQ(VecSum(user->workN2,f_reg)); 99c4762a1bSJed Brown *f_reg *= user->lambda; 100c4762a1bSJed Brown /* compute regularizer gradient = lambda*x */ 1015f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseDivide(G_reg,X,user->workN3)); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(G_reg,user->lambda)); 103c4762a1bSJed Brown PetscFunctionReturn(0); 104c4762a1bSJed Brown } 105c4762a1bSJed Brown 106c4762a1bSJed Brown /*------------------------------------------------------------*/ 107c4762a1bSJed Brown 108c4762a1bSJed Brown PetscErrorCode RegularizerObjectiveAndGradient2(Tao tao,Vec X,PetscReal *f_reg,Vec G_reg,void *ptr) 109c4762a1bSJed Brown { 110c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 111c4762a1bSJed Brown PetscReal temp; 112c4762a1bSJed Brown 113c4762a1bSJed Brown PetscFunctionBegin; 114c4762a1bSJed Brown /* compute regularizer objective = lambda*|z|_2^2 */ 1155f80ce2aSJacob Faibussowitsch CHKERRQ(VecDot(X,X,&temp)); 116c4762a1bSJed Brown *f_reg = 0.5*user->lambda*temp; 117c4762a1bSJed Brown /* compute regularizer gradient = lambda*z */ 1185f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(X,G_reg)); 1195f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(G_reg,user->lambda)); 120c4762a1bSJed Brown PetscFunctionReturn(0); 121c4762a1bSJed Brown } 122c4762a1bSJed Brown 123c4762a1bSJed Brown /*------------------------------------------------------------*/ 124c4762a1bSJed Brown 125c4762a1bSJed Brown static PetscErrorCode HessianMisfit(Tao tao, Vec x, Mat H, Mat Hpre, void *ptr) 126c4762a1bSJed Brown { 127c4762a1bSJed Brown PetscFunctionBegin; 128c4762a1bSJed Brown PetscFunctionReturn(0); 129c4762a1bSJed Brown } 130c4762a1bSJed Brown 131c4762a1bSJed Brown /*------------------------------------------------------------*/ 132c4762a1bSJed Brown 133c4762a1bSJed Brown static PetscErrorCode HessianReg(Tao tao, Vec x, Mat H, Mat Hpre, void *ptr) 134c4762a1bSJed Brown { 135c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 136c4762a1bSJed Brown 137c4762a1bSJed Brown PetscFunctionBegin; 1385f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(user->D,x,user->workN)); 1395f80ce2aSJacob Faibussowitsch CHKERRQ(VecPow(user->workN2,2.)); 1405f80ce2aSJacob Faibussowitsch CHKERRQ(VecShift(user->workN2,user->eps*user->eps)); 1415f80ce2aSJacob Faibussowitsch CHKERRQ(VecSqrtAbs(user->workN2)); 1425f80ce2aSJacob Faibussowitsch CHKERRQ(VecShift(user->workN2,-user->eps)); 1435f80ce2aSJacob Faibussowitsch CHKERRQ(VecReciprocal(user->workN2)); 1445f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(user->workN2,user->eps*user->eps)); 1455f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalSet(H,user->workN2,INSERT_VALUES)); 146c4762a1bSJed Brown PetscFunctionReturn(0); 147c4762a1bSJed Brown } 148c4762a1bSJed Brown 149c4762a1bSJed Brown /*------------------------------------------------------------*/ 150c4762a1bSJed Brown 151c4762a1bSJed Brown PetscErrorCode FullObjGrad(Tao tao,Vec X,PetscReal *f,Vec g,void *ptr) 152c4762a1bSJed Brown { 153c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 154c4762a1bSJed Brown PetscReal f_reg; 155c4762a1bSJed Brown 156c4762a1bSJed Brown PetscFunctionBegin; 157c4762a1bSJed Brown /* Objective 0.5*||Ax-b||_2^2 + lambda*||x||_2^2*/ 1585f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(user->A,X,user->workM)); 1595f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(user->workM,-1,user->b)); 1605f80ce2aSJacob Faibussowitsch CHKERRQ(VecDot(user->workM,user->workM,f)); 1615f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(X,NORM_2,&f_reg)); 162c4762a1bSJed Brown *f *= 0.5; 163c4762a1bSJed Brown *f += user->lambda*f_reg*f_reg; 164c4762a1bSJed Brown /* Gradient. ATAx-ATb + 2*lambda*x */ 1655f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(user->ATA,X,user->workN)); 1665f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(user->A,user->b,user->workN2)); 1675f80ce2aSJacob Faibussowitsch CHKERRQ(VecWAXPY(g,-1.,user->workN2,user->workN)); 1685f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(g,2*user->lambda,X)); 169c4762a1bSJed Brown PetscFunctionReturn(0); 170c4762a1bSJed Brown } 171c4762a1bSJed Brown /*------------------------------------------------------------*/ 172c4762a1bSJed Brown 173c4762a1bSJed Brown static PetscErrorCode HessianFull(Tao tao, Vec x, Mat H, Mat Hpre, void *ptr) 174c4762a1bSJed Brown { 175c4762a1bSJed Brown PetscFunctionBegin; 176c4762a1bSJed Brown PetscFunctionReturn(0); 177c4762a1bSJed Brown } 178c4762a1bSJed Brown /*------------------------------------------------------------*/ 179c4762a1bSJed Brown 180c4762a1bSJed Brown PetscErrorCode InitializeUserData(AppCtx *user) 181c4762a1bSJed Brown { 182c4762a1bSJed 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". */ 183c4762a1bSJed Brown PetscViewer fd; /* used to load data from file */ 184c4762a1bSJed Brown PetscErrorCode ierr; 185c4762a1bSJed Brown PetscInt k,n; 186c4762a1bSJed Brown PetscScalar v; 187c4762a1bSJed Brown PetscFunctionBegin; 188c4762a1bSJed Brown 189c4762a1bSJed Brown /* Load the A matrix, b vector, and xGT vector from a binary file. */ 1905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,dataFile,FILE_MODE_READ,&fd)); 1915f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&user->A)); 1925f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(user->A,MATAIJ)); 1935f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(user->A,fd)); 1945f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&user->b)); 1955f80ce2aSJacob Faibussowitsch CHKERRQ(VecLoad(user->b,fd)); 1965f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&user->xGT)); 1975f80ce2aSJacob Faibussowitsch CHKERRQ(VecLoad(user->xGT,fd)); 1985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&fd)); 199c4762a1bSJed Brown 2005f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(user->A,&user->M,&user->N)); 201c4762a1bSJed Brown 2025f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&user->D)); 2035f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(user->D,PETSC_DECIDE,PETSC_DECIDE,user->N,user->N)); 2045f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(user->D)); 2055f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(user->D)); 206c4762a1bSJed Brown for (k=0; k<user->N; k++) { 207c4762a1bSJed Brown v = 1.0; 208c4762a1bSJed Brown n = k+1; 209c4762a1bSJed Brown if (k< user->N -1) { 2105f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(user->D,1,&k,1,&n,&v,INSERT_VALUES)); 211c4762a1bSJed Brown } 212c4762a1bSJed Brown v = -1.0; 2135f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(user->D,1,&k,1,&k,&v,INSERT_VALUES)); 214c4762a1bSJed Brown } 2155f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(user->D,MAT_FINAL_ASSEMBLY)); 2165f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(user->D,MAT_FINAL_ASSEMBLY)); 217c4762a1bSJed Brown 2185f80ce2aSJacob Faibussowitsch CHKERRQ(MatTransposeMatMult(user->D,user->D,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&user->DTD)); 219c4762a1bSJed Brown 2205f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&user->Hz)); 2215f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(user->Hz,PETSC_DECIDE,PETSC_DECIDE,user->N,user->N)); 2225f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(user->Hz)); 2235f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(user->Hz)); 2245f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(user->Hz,MAT_FINAL_ASSEMBLY)); 2255f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(user->Hz,MAT_FINAL_ASSEMBLY)); 226c4762a1bSJed Brown 2275f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&(user->x))); 2285f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&(user->workM))); 2295f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&(user->workN))); 2305f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&(user->workN2))); 2315f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(user->x,PETSC_DECIDE,user->N)); 2325f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(user->workM,PETSC_DECIDE,user->M)); 2335f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(user->workN,PETSC_DECIDE,user->N)); 2345f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(user->workN2,PETSC_DECIDE,user->N)); 2355f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(user->x)); 2365f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(user->workM)); 2375f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(user->workN)); 2385f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(user->workN2)); 239c4762a1bSJed Brown 2405f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(user->workN,&(user->workN3))); 2415f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(user->x,&(user->xlb))); 2425f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(user->x,&(user->xub))); 2435f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(user->x,&(user->c))); 2445f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(user->xlb,0.0)); 2455f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(user->c,0.0)); 2465f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(user->xub,PETSC_INFINITY)); 247c4762a1bSJed Brown 2485f80ce2aSJacob Faibussowitsch CHKERRQ(MatTransposeMatMult(user->A,user->A, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &(user->ATA))); 2495f80ce2aSJacob Faibussowitsch CHKERRQ(MatTransposeMatMult(user->A,user->A, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &(user->Hx))); 2505f80ce2aSJacob Faibussowitsch CHKERRQ(MatTransposeMatMult(user->A,user->A, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &(user->HF))); 251c4762a1bSJed Brown 2525f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(user->ATA,MAT_FINAL_ASSEMBLY)); 2535f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(user->ATA,MAT_FINAL_ASSEMBLY)); 2545f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(user->Hx,MAT_FINAL_ASSEMBLY)); 2555f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(user->Hx,MAT_FINAL_ASSEMBLY)); 2565f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(user->HF,MAT_FINAL_ASSEMBLY)); 2575f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(user->HF,MAT_FINAL_ASSEMBLY)); 258c4762a1bSJed Brown 259c4762a1bSJed Brown user->lambda = 1.e-8; 260c4762a1bSJed Brown user->eps = 1.e-3; 261c4762a1bSJed Brown user->reg = 2; 262c4762a1bSJed Brown user->mumin = 5.e-6; 263c4762a1bSJed Brown 264c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Configure separable objection example", "tomographyADMM.c");CHKERRQ(ierr); 2655f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-reg","Regularization scheme for z solver (1,2)", "tomographyADMM.c", user->reg, &(user->reg), NULL)); 2665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-lambda", "The regularization multiplier. 1 default", "tomographyADMM.c", user->lambda, &(user->lambda), NULL)); 2675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-eps", "L1 norm epsilon padding", "tomographyADMM.c", user->eps, &(user->eps), NULL)); 2685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-mumin", "Minimum value for ADMM spectral penalty", "tomographyADMM.c", user->mumin, &(user->mumin), NULL)); 269c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 270c4762a1bSJed Brown PetscFunctionReturn(0); 271c4762a1bSJed Brown } 272c4762a1bSJed Brown 273c4762a1bSJed Brown /*------------------------------------------------------------*/ 274c4762a1bSJed Brown 275c4762a1bSJed Brown PetscErrorCode DestroyContext(AppCtx *user) 276c4762a1bSJed Brown { 277c4762a1bSJed Brown PetscFunctionBegin; 2785f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&user->A)); 2795f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&user->ATA)); 2805f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&user->Hx)); 2815f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&user->Hz)); 2825f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&user->HF)); 2835f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&user->D)); 2845f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&user->DTD)); 2855f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&user->xGT)); 2865f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&user->xlb)); 2875f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&user->xub)); 2885f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&user->b)); 2895f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&user->x)); 2905f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&user->c)); 2915f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&user->workN3)); 2925f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&user->workN2)); 2935f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&user->workN)); 2945f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&user->workM)); 295c4762a1bSJed Brown PetscFunctionReturn(0); 296c4762a1bSJed Brown } 297c4762a1bSJed Brown 298c4762a1bSJed Brown /*------------------------------------------------------------*/ 299c4762a1bSJed Brown 300c4762a1bSJed Brown int main(int argc,char **argv) 301c4762a1bSJed Brown { 302c4762a1bSJed Brown Tao tao,misfit,reg; 303c4762a1bSJed Brown PetscReal v1,v2; 304c4762a1bSJed Brown AppCtx* user; 305c4762a1bSJed Brown PetscViewer fd; 306c4762a1bSJed Brown char resultFile[] = "tomographyResult_x"; 307c4762a1bSJed Brown 308*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 3095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNew(&user)); 3105f80ce2aSJacob Faibussowitsch CHKERRQ(InitializeUserData(user)); 311c4762a1bSJed Brown 3125f80ce2aSJacob Faibussowitsch CHKERRQ(TaoCreate(PETSC_COMM_WORLD, &tao)); 3135f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetType(tao, TAOADMM)); 3145f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetSolution(tao, user->x)); 315c4762a1bSJed Brown /* f(x) + g(x) for parent tao */ 3165f80ce2aSJacob Faibussowitsch CHKERRQ(TaoADMMSetSpectralPenalty(tao,1.)); 3175f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetObjectiveAndGradient(tao,NULL, FullObjGrad, (void*)user)); 3185f80ce2aSJacob Faibussowitsch CHKERRQ(MatShift(user->HF,user->lambda)); 3195f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetHessian(tao, user->HF, user->HF, HessianFull, (void*)user)); 320c4762a1bSJed Brown 321c4762a1bSJed Brown /* f(x) for misfit tao */ 3225f80ce2aSJacob Faibussowitsch CHKERRQ(TaoADMMSetMisfitObjectiveAndGradientRoutine(tao, MisfitObjectiveAndGradient, (void*)user)); 3235f80ce2aSJacob Faibussowitsch CHKERRQ(TaoADMMSetMisfitHessianRoutine(tao, user->Hx, user->Hx, HessianMisfit, (void*)user)); 3245f80ce2aSJacob Faibussowitsch CHKERRQ(TaoADMMSetMisfitHessianChangeStatus(tao,PETSC_FALSE)); 3255f80ce2aSJacob Faibussowitsch CHKERRQ(TaoADMMSetMisfitConstraintJacobian(tao,user->D,user->D,NullJacobian,(void*)user)); 326c4762a1bSJed Brown 327c4762a1bSJed Brown /* g(x) for regularizer tao */ 328c4762a1bSJed Brown if (user->reg == 1) { 3295f80ce2aSJacob Faibussowitsch CHKERRQ(TaoADMMSetRegularizerObjectiveAndGradientRoutine(tao, RegularizerObjectiveAndGradient1, (void*)user)); 3305f80ce2aSJacob Faibussowitsch CHKERRQ(TaoADMMSetRegularizerHessianRoutine(tao, user->Hz, user->Hz, HessianReg, (void*)user)); 3315f80ce2aSJacob Faibussowitsch CHKERRQ(TaoADMMSetRegHessianChangeStatus(tao,PETSC_TRUE)); 332c4762a1bSJed Brown } else if (user->reg == 2) { 3335f80ce2aSJacob Faibussowitsch CHKERRQ(TaoADMMSetRegularizerObjectiveAndGradientRoutine(tao, RegularizerObjectiveAndGradient2, (void*)user)); 3345f80ce2aSJacob Faibussowitsch CHKERRQ(MatShift(user->Hz,1)); 3355f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(user->Hz,user->lambda)); 3365f80ce2aSJacob Faibussowitsch CHKERRQ(TaoADMMSetRegularizerHessianRoutine(tao, user->Hz, user->Hz, HessianMisfit, (void*)user)); 3375f80ce2aSJacob Faibussowitsch CHKERRQ(TaoADMMSetRegHessianChangeStatus(tao,PETSC_TRUE)); 3383c859ba3SBarry Smith } else PetscCheck(user->reg == 3,PETSC_COMM_WORLD, PETSC_ERR_ARG_UNKNOWN_TYPE, "Incorrect Reg type"); /* TaoShell case */ 339c4762a1bSJed Brown 340c4762a1bSJed Brown /* Set type for the misfit solver */ 3415f80ce2aSJacob Faibussowitsch CHKERRQ(TaoADMMGetMisfitSubsolver(tao, &misfit)); 3425f80ce2aSJacob Faibussowitsch CHKERRQ(TaoADMMGetRegularizationSubsolver(tao, ®)); 3435f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetType(misfit,TAONLS)); 344c4762a1bSJed Brown if (user->reg == 3) { 3455f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetType(reg,TAOSHELL)); 3465f80ce2aSJacob Faibussowitsch CHKERRQ(TaoShellSetContext(reg, (void*) user)); 3475f80ce2aSJacob Faibussowitsch CHKERRQ(TaoShellSetSolve(reg, TaoShellSolve_SoftThreshold)); 348c4762a1bSJed Brown } else { 3495f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetType(reg,TAONLS)); 350c4762a1bSJed Brown } 3515f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetVariableBounds(misfit,user->xlb,user->xub)); 352c4762a1bSJed Brown 353c4762a1bSJed Brown /* Soft Thresholding solves the ADMM problem with the L1 regularizer lambda*||z||_1 and the x-z=0 constraint */ 3545f80ce2aSJacob Faibussowitsch CHKERRQ(TaoADMMSetRegularizerCoefficient(tao, user->lambda)); 3555f80ce2aSJacob Faibussowitsch CHKERRQ(TaoADMMSetRegularizerConstraintJacobian(tao,NULL,NULL,NullJacobian,(void*)user)); 3565f80ce2aSJacob Faibussowitsch CHKERRQ(TaoADMMSetMinimumSpectralPenalty(tao,user->mumin)); 357c4762a1bSJed Brown 3585f80ce2aSJacob Faibussowitsch CHKERRQ(TaoADMMSetConstraintVectorRHS(tao,user->c)); 3595f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetFromOptions(tao)); 3605f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSolve(tao)); 361c4762a1bSJed Brown 362c4762a1bSJed Brown /* Save x (reconstruction of object) vector to a binary file, which maybe read from Matlab and convert to a 2D image for comparison. */ 3635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,resultFile,FILE_MODE_WRITE,&fd)); 3645f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(user->x,fd)); 3655f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&fd)); 366c4762a1bSJed Brown 367c4762a1bSJed Brown /* compute the error */ 3685f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(user->x,-1,user->xGT)); 3695f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(user->x,NORM_2,&v1)); 3705f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(user->xGT,NORM_2,&v2)); 3715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "relative reconstruction error: ||x-xGT||/||xGT|| = %6.4e.\n", (double)(v1/v2))); 372c4762a1bSJed Brown 373c4762a1bSJed Brown /* Free TAO data structures */ 3745f80ce2aSJacob Faibussowitsch CHKERRQ(TaoDestroy(&tao)); 3755f80ce2aSJacob Faibussowitsch CHKERRQ(DestroyContext(user)); 3765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(user)); 377*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 378*b122ec5aSJacob Faibussowitsch return 0; 379c4762a1bSJed Brown } 380c4762a1bSJed Brown 381c4762a1bSJed Brown /*TEST 382c4762a1bSJed Brown 383c4762a1bSJed Brown build: 384dfd57a17SPierre Jolivet requires: !complex !single !__float128 !defined(PETSC_USE_64BIT_INDICES) 385c4762a1bSJed Brown 386c4762a1bSJed Brown test: 387c4762a1bSJed Brown suffix: 1 388c4762a1bSJed Brown localrunfiles: tomographyData_A_b_xGT 389c4762a1bSJed Brown args: -lambda 1.e-8 -tao_monitor -tao_type nls -tao_nls_pc_type icc 390c4762a1bSJed Brown 391c4762a1bSJed Brown test: 392c4762a1bSJed Brown suffix: 2 393c4762a1bSJed Brown localrunfiles: tomographyData_A_b_xGT 394c4762a1bSJed 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 395c4762a1bSJed Brown 396c4762a1bSJed Brown test: 397c4762a1bSJed Brown suffix: 3 398c4762a1bSJed Brown localrunfiles: tomographyData_A_b_xGT 399c4762a1bSJed 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 400c4762a1bSJed Brown 401c4762a1bSJed Brown test: 402c4762a1bSJed Brown suffix: 4 403c4762a1bSJed Brown localrunfiles: tomographyData_A_b_xGT 404c4762a1bSJed 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 405c4762a1bSJed Brown 406c4762a1bSJed Brown test: 407c4762a1bSJed Brown suffix: 5 408c4762a1bSJed Brown localrunfiles: tomographyData_A_b_xGT 409c4762a1bSJed 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 410c4762a1bSJed Brown 411c4762a1bSJed Brown test: 412c4762a1bSJed Brown suffix: 6 413c4762a1bSJed Brown localrunfiles: tomographyData_A_b_xGT 414c4762a1bSJed 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 415c4762a1bSJed Brown 416c4762a1bSJed Brown TEST*/ 417