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 PetscErrorCode ierr; 37c4762a1bSJed Brown PetscReal lambda, mu; 38c4762a1bSJed Brown AppCtx *user; 39c4762a1bSJed Brown Vec out,work,y,x; 40c4762a1bSJed Brown Tao admm_tao,misfit; 41c4762a1bSJed Brown 42c4762a1bSJed Brown PetscFunctionBegin; 43c4762a1bSJed Brown user = NULL; 44c4762a1bSJed Brown mu = 0; 45c4762a1bSJed Brown ierr = TaoGetADMMParentTao(tao,&admm_tao);CHKERRQ(ierr); 46c4762a1bSJed Brown ierr = TaoADMMGetMisfitSubsolver(admm_tao, &misfit);CHKERRQ(ierr); 47c4762a1bSJed Brown ierr = TaoADMMGetSpectralPenalty(admm_tao,&mu);CHKERRQ(ierr); 483ec1f749SStefano Zampini ierr = TaoShellGetContext(tao,&user);CHKERRQ(ierr); 49c4762a1bSJed Brown 50c4762a1bSJed Brown lambda = user->lambda; 51c4762a1bSJed Brown work = user->workN; 52c4762a1bSJed Brown ierr = TaoGetSolutionVector(tao, &out);CHKERRQ(ierr); 53c4762a1bSJed Brown ierr = TaoGetSolutionVector(misfit, &x);CHKERRQ(ierr); 54c4762a1bSJed Brown ierr = TaoADMMGetDualVector(admm_tao, &y);CHKERRQ(ierr); 55c4762a1bSJed Brown 56c4762a1bSJed Brown /* Dx + y/mu */ 57c4762a1bSJed Brown ierr = MatMult(user->D,x,work);CHKERRQ(ierr); 58c4762a1bSJed Brown ierr = VecAXPY(work,1/mu,y);CHKERRQ(ierr); 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* soft thresholding */ 61c4762a1bSJed Brown ierr = TaoSoftThreshold(work, -lambda/mu, lambda/mu, out);CHKERRQ(ierr); 62c4762a1bSJed Brown PetscFunctionReturn(0); 63c4762a1bSJed Brown } 64c4762a1bSJed Brown 65c4762a1bSJed Brown /*------------------------------------------------------------*/ 66c4762a1bSJed Brown 67c4762a1bSJed Brown PetscErrorCode MisfitObjectiveAndGradient(Tao tao,Vec X,PetscReal *f,Vec g,void *ptr) 68c4762a1bSJed Brown { 69c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 70c4762a1bSJed Brown PetscErrorCode ierr; 71c4762a1bSJed Brown 72c4762a1bSJed Brown PetscFunctionBegin; 73c4762a1bSJed Brown /* Objective 0.5*||Ax-b||_2^2 */ 74c4762a1bSJed Brown ierr = MatMult(user->A,X,user->workM);CHKERRQ(ierr); 75c4762a1bSJed Brown ierr = VecAXPY(user->workM,-1,user->b);CHKERRQ(ierr); 76c4762a1bSJed Brown ierr = VecDot(user->workM,user->workM,f);CHKERRQ(ierr); 77c4762a1bSJed Brown *f *= 0.5; 78c4762a1bSJed Brown /* Gradient. ATAx-ATb */ 79c4762a1bSJed Brown ierr = MatMult(user->ATA,X,user->workN);CHKERRQ(ierr); 80c4762a1bSJed Brown ierr = MatMultTranspose(user->A,user->b,user->workN2);CHKERRQ(ierr); 81c4762a1bSJed Brown ierr = VecWAXPY(g,-1.,user->workN2,user->workN);CHKERRQ(ierr); 82c4762a1bSJed Brown PetscFunctionReturn(0); 83c4762a1bSJed Brown } 84c4762a1bSJed Brown 85c4762a1bSJed Brown /*------------------------------------------------------------*/ 86c4762a1bSJed Brown 87c4762a1bSJed Brown PetscErrorCode RegularizerObjectiveAndGradient1(Tao tao,Vec X,PetscReal *f_reg,Vec G_reg,void *ptr) 88c4762a1bSJed Brown { 89c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 90c4762a1bSJed Brown PetscErrorCode ierr; 91c4762a1bSJed Brown 92c4762a1bSJed Brown PetscFunctionBegin; 93c4762a1bSJed Brown /* compute regularizer objective 94c4762a1bSJed Brown * f = f + lambda*sum(sqrt(y.^2+epsilon^2) - epsilon), where y = D*x */ 95c4762a1bSJed Brown ierr = VecCopy(X,user->workN2);CHKERRQ(ierr); 96c4762a1bSJed Brown ierr = VecPow(user->workN2,2.);CHKERRQ(ierr); 97c4762a1bSJed Brown ierr = VecShift(user->workN2,user->eps*user->eps);CHKERRQ(ierr); 98c4762a1bSJed Brown ierr = VecSqrtAbs(user->workN2);CHKERRQ(ierr); 99c4762a1bSJed Brown ierr = VecCopy(user->workN2, user->workN3);CHKERRQ(ierr); 100c4762a1bSJed Brown ierr = VecShift(user->workN2,-user->eps);CHKERRQ(ierr); 101c4762a1bSJed Brown ierr = VecSum(user->workN2,f_reg);CHKERRQ(ierr); 102c4762a1bSJed Brown *f_reg *= user->lambda; 103c4762a1bSJed Brown /* compute regularizer gradient = lambda*x */ 104c4762a1bSJed Brown ierr = VecPointwiseDivide(G_reg,X,user->workN3);CHKERRQ(ierr); 105c4762a1bSJed Brown ierr = VecScale(G_reg,user->lambda);CHKERRQ(ierr); 106c4762a1bSJed Brown PetscFunctionReturn(0); 107c4762a1bSJed Brown } 108c4762a1bSJed Brown 109c4762a1bSJed Brown /*------------------------------------------------------------*/ 110c4762a1bSJed Brown 111c4762a1bSJed Brown PetscErrorCode RegularizerObjectiveAndGradient2(Tao tao,Vec X,PetscReal *f_reg,Vec G_reg,void *ptr) 112c4762a1bSJed Brown { 113c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 114c4762a1bSJed Brown PetscErrorCode ierr; 115c4762a1bSJed Brown PetscReal temp; 116c4762a1bSJed Brown 117c4762a1bSJed Brown PetscFunctionBegin; 118c4762a1bSJed Brown /* compute regularizer objective = lambda*|z|_2^2 */ 119c4762a1bSJed Brown ierr = VecDot(X,X,&temp);CHKERRQ(ierr); 120c4762a1bSJed Brown *f_reg = 0.5*user->lambda*temp; 121c4762a1bSJed Brown /* compute regularizer gradient = lambda*z */ 1221e1ea65dSPierre Jolivet ierr = VecCopy(X,G_reg);CHKERRQ(ierr); 123c4762a1bSJed Brown ierr = VecScale(G_reg,user->lambda);CHKERRQ(ierr); 124c4762a1bSJed Brown PetscFunctionReturn(0); 125c4762a1bSJed Brown } 126c4762a1bSJed Brown 127c4762a1bSJed Brown /*------------------------------------------------------------*/ 128c4762a1bSJed Brown 129c4762a1bSJed Brown static PetscErrorCode HessianMisfit(Tao tao, Vec x, Mat H, Mat Hpre, void *ptr) 130c4762a1bSJed Brown { 131c4762a1bSJed Brown PetscFunctionBegin; 132c4762a1bSJed Brown PetscFunctionReturn(0); 133c4762a1bSJed Brown } 134c4762a1bSJed Brown 135c4762a1bSJed Brown /*------------------------------------------------------------*/ 136c4762a1bSJed Brown 137c4762a1bSJed Brown static PetscErrorCode HessianReg(Tao tao, Vec x, Mat H, Mat Hpre, void *ptr) 138c4762a1bSJed Brown { 139c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 140c4762a1bSJed Brown PetscErrorCode ierr; 141c4762a1bSJed Brown 142c4762a1bSJed Brown PetscFunctionBegin; 143c4762a1bSJed Brown ierr = MatMult(user->D,x,user->workN);CHKERRQ(ierr); 144c4762a1bSJed Brown ierr = VecPow(user->workN2,2.);CHKERRQ(ierr); 145c4762a1bSJed Brown ierr = VecShift(user->workN2,user->eps*user->eps);CHKERRQ(ierr); 146c4762a1bSJed Brown ierr = VecSqrtAbs(user->workN2);CHKERRQ(ierr); 147c4762a1bSJed Brown ierr = VecShift(user->workN2,-user->eps);CHKERRQ(ierr); 148c4762a1bSJed Brown ierr = VecReciprocal(user->workN2);CHKERRQ(ierr); 149c4762a1bSJed Brown ierr = VecScale(user->workN2,user->eps*user->eps);CHKERRQ(ierr); 150c4762a1bSJed Brown ierr = MatDiagonalSet(H,user->workN2,INSERT_VALUES);CHKERRQ(ierr); 151c4762a1bSJed Brown PetscFunctionReturn(0); 152c4762a1bSJed Brown } 153c4762a1bSJed Brown 154c4762a1bSJed Brown /*------------------------------------------------------------*/ 155c4762a1bSJed Brown 156c4762a1bSJed Brown PetscErrorCode FullObjGrad(Tao tao,Vec X,PetscReal *f,Vec g,void *ptr) 157c4762a1bSJed Brown { 158c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 159c4762a1bSJed Brown PetscErrorCode ierr; 160c4762a1bSJed Brown PetscReal f_reg; 161c4762a1bSJed Brown 162c4762a1bSJed Brown PetscFunctionBegin; 163c4762a1bSJed Brown /* Objective 0.5*||Ax-b||_2^2 + lambda*||x||_2^2*/ 164c4762a1bSJed Brown ierr = MatMult(user->A,X,user->workM);CHKERRQ(ierr); 165c4762a1bSJed Brown ierr = VecAXPY(user->workM,-1,user->b);CHKERRQ(ierr); 166c4762a1bSJed Brown ierr = VecDot(user->workM,user->workM,f);CHKERRQ(ierr); 167c4762a1bSJed Brown ierr = VecNorm(X,NORM_2,&f_reg);CHKERRQ(ierr); 168c4762a1bSJed Brown *f *= 0.5; 169c4762a1bSJed Brown *f += user->lambda*f_reg*f_reg; 170c4762a1bSJed Brown /* Gradient. ATAx-ATb + 2*lambda*x */ 171c4762a1bSJed Brown ierr = MatMult(user->ATA,X,user->workN);CHKERRQ(ierr); 172c4762a1bSJed Brown ierr = MatMultTranspose(user->A,user->b,user->workN2);CHKERRQ(ierr); 173c4762a1bSJed Brown ierr = VecWAXPY(g,-1.,user->workN2,user->workN);CHKERRQ(ierr); 174c4762a1bSJed Brown ierr = VecAXPY(g,2*user->lambda,X);CHKERRQ(ierr); 175c4762a1bSJed Brown PetscFunctionReturn(0); 176c4762a1bSJed Brown } 177c4762a1bSJed Brown /*------------------------------------------------------------*/ 178c4762a1bSJed Brown 179c4762a1bSJed Brown static PetscErrorCode HessianFull(Tao tao, Vec x, Mat H, Mat Hpre, void *ptr) 180c4762a1bSJed Brown { 181c4762a1bSJed Brown PetscFunctionBegin; 182c4762a1bSJed Brown PetscFunctionReturn(0); 183c4762a1bSJed Brown } 184c4762a1bSJed Brown /*------------------------------------------------------------*/ 185c4762a1bSJed Brown 186c4762a1bSJed Brown PetscErrorCode InitializeUserData(AppCtx *user) 187c4762a1bSJed Brown { 188c4762a1bSJed 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". */ 189c4762a1bSJed Brown PetscViewer fd; /* used to load data from file */ 190c4762a1bSJed Brown PetscErrorCode ierr; 191c4762a1bSJed Brown PetscInt k,n; 192c4762a1bSJed Brown PetscScalar v; 193c4762a1bSJed Brown PetscFunctionBegin; 194c4762a1bSJed Brown 195c4762a1bSJed Brown /* Load the A matrix, b vector, and xGT vector from a binary file. */ 196c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,dataFile,FILE_MODE_READ,&fd);CHKERRQ(ierr); 197c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&user->A);CHKERRQ(ierr); 198c4762a1bSJed Brown ierr = MatSetType(user->A,MATAIJ);CHKERRQ(ierr); 199c4762a1bSJed Brown ierr = MatLoad(user->A,fd);CHKERRQ(ierr); 200c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&user->b);CHKERRQ(ierr); 201c4762a1bSJed Brown ierr = VecLoad(user->b,fd);CHKERRQ(ierr); 202c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&user->xGT);CHKERRQ(ierr); 203c4762a1bSJed Brown ierr = VecLoad(user->xGT,fd);CHKERRQ(ierr); 204c4762a1bSJed Brown ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr); 205c4762a1bSJed Brown 206c4762a1bSJed Brown ierr = MatGetSize(user->A,&user->M,&user->N);CHKERRQ(ierr); 207c4762a1bSJed Brown 208c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&user->D);CHKERRQ(ierr); 209c4762a1bSJed Brown ierr = MatSetSizes(user->D,PETSC_DECIDE,PETSC_DECIDE,user->N,user->N);CHKERRQ(ierr); 210c4762a1bSJed Brown ierr = MatSetFromOptions(user->D);CHKERRQ(ierr); 211c4762a1bSJed Brown ierr = MatSetUp(user->D);CHKERRQ(ierr); 212c4762a1bSJed Brown for (k=0; k<user->N; k++) { 213c4762a1bSJed Brown v = 1.0; 214c4762a1bSJed Brown n = k+1; 215c4762a1bSJed Brown if (k< user->N -1) { 216c4762a1bSJed Brown ierr = MatSetValues(user->D,1,&k,1,&n,&v,INSERT_VALUES);CHKERRQ(ierr); 217c4762a1bSJed Brown } 218c4762a1bSJed Brown v = -1.0; 219c4762a1bSJed Brown ierr = MatSetValues(user->D,1,&k,1,&k,&v,INSERT_VALUES);CHKERRQ(ierr); 220c4762a1bSJed Brown } 221c4762a1bSJed Brown ierr = MatAssemblyBegin(user->D,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 222c4762a1bSJed Brown ierr = MatAssemblyEnd(user->D,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 223c4762a1bSJed Brown 224c4762a1bSJed Brown ierr = MatTransposeMatMult(user->D,user->D,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&user->DTD);CHKERRQ(ierr); 225c4762a1bSJed Brown 226c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&user->Hz);CHKERRQ(ierr); 227c4762a1bSJed Brown ierr = MatSetSizes(user->Hz,PETSC_DECIDE,PETSC_DECIDE,user->N,user->N);CHKERRQ(ierr); 228c4762a1bSJed Brown ierr = MatSetFromOptions(user->Hz);CHKERRQ(ierr); 229c4762a1bSJed Brown ierr = MatSetUp(user->Hz);CHKERRQ(ierr); 230c4762a1bSJed Brown ierr = MatAssemblyBegin(user->Hz,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 231c4762a1bSJed Brown ierr = MatAssemblyEnd(user->Hz,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 232c4762a1bSJed Brown 233c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&(user->x));CHKERRQ(ierr); 234c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&(user->workM));CHKERRQ(ierr); 235c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&(user->workN));CHKERRQ(ierr); 236c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&(user->workN2));CHKERRQ(ierr); 237c4762a1bSJed Brown ierr = VecSetSizes(user->x,PETSC_DECIDE,user->N);CHKERRQ(ierr); 238c4762a1bSJed Brown ierr = VecSetSizes(user->workM,PETSC_DECIDE,user->M);CHKERRQ(ierr); 239c4762a1bSJed Brown ierr = VecSetSizes(user->workN,PETSC_DECIDE,user->N);CHKERRQ(ierr); 240c4762a1bSJed Brown ierr = VecSetSizes(user->workN2,PETSC_DECIDE,user->N);CHKERRQ(ierr); 241c4762a1bSJed Brown ierr = VecSetFromOptions(user->x);CHKERRQ(ierr); 242c4762a1bSJed Brown ierr = VecSetFromOptions(user->workM);CHKERRQ(ierr); 243c4762a1bSJed Brown ierr = VecSetFromOptions(user->workN);CHKERRQ(ierr); 244c4762a1bSJed Brown ierr = VecSetFromOptions(user->workN2);CHKERRQ(ierr); 245c4762a1bSJed Brown 246c4762a1bSJed Brown ierr = VecDuplicate(user->workN,&(user->workN3));CHKERRQ(ierr); 247c4762a1bSJed Brown ierr = VecDuplicate(user->x,&(user->xlb));CHKERRQ(ierr); 248c4762a1bSJed Brown ierr = VecDuplicate(user->x,&(user->xub));CHKERRQ(ierr); 249c4762a1bSJed Brown ierr = VecDuplicate(user->x,&(user->c));CHKERRQ(ierr); 250c4762a1bSJed Brown ierr = VecSet(user->xlb,0.0);CHKERRQ(ierr); 251c4762a1bSJed Brown ierr = VecSet(user->c,0.0);CHKERRQ(ierr); 252c4762a1bSJed Brown ierr = VecSet(user->xub,PETSC_INFINITY);CHKERRQ(ierr); 253c4762a1bSJed Brown 254c4762a1bSJed Brown ierr = MatTransposeMatMult(user->A,user->A, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &(user->ATA));CHKERRQ(ierr); 255c4762a1bSJed Brown ierr = MatTransposeMatMult(user->A,user->A, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &(user->Hx));CHKERRQ(ierr); 256c4762a1bSJed Brown ierr = MatTransposeMatMult(user->A,user->A, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &(user->HF));CHKERRQ(ierr); 257c4762a1bSJed Brown 258c4762a1bSJed Brown ierr = MatAssemblyBegin(user->ATA,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 259c4762a1bSJed Brown ierr = MatAssemblyEnd(user->ATA,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 260c4762a1bSJed Brown ierr = MatAssemblyBegin(user->Hx,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 261c4762a1bSJed Brown ierr = MatAssemblyEnd(user->Hx,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 262c4762a1bSJed Brown ierr = MatAssemblyBegin(user->HF,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 263c4762a1bSJed Brown ierr = MatAssemblyEnd(user->HF,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 264c4762a1bSJed Brown 265c4762a1bSJed Brown user->lambda = 1.e-8; 266c4762a1bSJed Brown user->eps = 1.e-3; 267c4762a1bSJed Brown user->reg = 2; 268c4762a1bSJed Brown user->mumin = 5.e-6; 269c4762a1bSJed Brown 270c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Configure separable objection example", "tomographyADMM.c");CHKERRQ(ierr); 271c4762a1bSJed Brown ierr = PetscOptionsInt("-reg","Regularization scheme for z solver (1,2)", "tomographyADMM.c", user->reg, &(user->reg), NULL);CHKERRQ(ierr); 272c4762a1bSJed Brown ierr = PetscOptionsReal("-lambda", "The regularization multiplier. 1 default", "tomographyADMM.c", user->lambda, &(user->lambda), NULL);CHKERRQ(ierr); 273c4762a1bSJed Brown ierr = PetscOptionsReal("-eps", "L1 norm epsilon padding", "tomographyADMM.c", user->eps, &(user->eps), NULL);CHKERRQ(ierr); 274c4762a1bSJed Brown ierr = PetscOptionsReal("-mumin", "Minimum value for ADMM spectral penalty", "tomographyADMM.c", user->mumin, &(user->mumin), NULL);CHKERRQ(ierr); 275c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 276c4762a1bSJed Brown PetscFunctionReturn(0); 277c4762a1bSJed Brown } 278c4762a1bSJed Brown 279c4762a1bSJed Brown /*------------------------------------------------------------*/ 280c4762a1bSJed Brown 281c4762a1bSJed Brown PetscErrorCode DestroyContext(AppCtx *user) 282c4762a1bSJed Brown { 283c4762a1bSJed Brown PetscErrorCode ierr; 284c4762a1bSJed Brown 285c4762a1bSJed Brown PetscFunctionBegin; 286c4762a1bSJed Brown ierr = MatDestroy(&user->A);CHKERRQ(ierr); 287c4762a1bSJed Brown ierr = MatDestroy(&user->ATA);CHKERRQ(ierr); 288c4762a1bSJed Brown ierr = MatDestroy(&user->Hx);CHKERRQ(ierr); 289c4762a1bSJed Brown ierr = MatDestroy(&user->Hz);CHKERRQ(ierr); 290c4762a1bSJed Brown ierr = MatDestroy(&user->HF);CHKERRQ(ierr); 291c4762a1bSJed Brown ierr = MatDestroy(&user->D);CHKERRQ(ierr); 292c4762a1bSJed Brown ierr = MatDestroy(&user->DTD);CHKERRQ(ierr); 293c4762a1bSJed Brown ierr = VecDestroy(&user->xGT);CHKERRQ(ierr); 294c4762a1bSJed Brown ierr = VecDestroy(&user->xlb);CHKERRQ(ierr); 295c4762a1bSJed Brown ierr = VecDestroy(&user->xub);CHKERRQ(ierr); 296c4762a1bSJed Brown ierr = VecDestroy(&user->b);CHKERRQ(ierr); 297c4762a1bSJed Brown ierr = VecDestroy(&user->x);CHKERRQ(ierr); 298c4762a1bSJed Brown ierr = VecDestroy(&user->c);CHKERRQ(ierr); 299c4762a1bSJed Brown ierr = VecDestroy(&user->workN3);CHKERRQ(ierr); 300c4762a1bSJed Brown ierr = VecDestroy(&user->workN2);CHKERRQ(ierr); 301c4762a1bSJed Brown ierr = VecDestroy(&user->workN);CHKERRQ(ierr); 302c4762a1bSJed Brown ierr = VecDestroy(&user->workM);CHKERRQ(ierr); 303c4762a1bSJed Brown PetscFunctionReturn(0); 304c4762a1bSJed Brown } 305c4762a1bSJed Brown 306c4762a1bSJed Brown /*------------------------------------------------------------*/ 307c4762a1bSJed Brown 308c4762a1bSJed Brown int main(int argc,char **argv) 309c4762a1bSJed Brown { 310c4762a1bSJed Brown PetscErrorCode ierr; 311c4762a1bSJed Brown Tao tao,misfit,reg; 312c4762a1bSJed Brown PetscReal v1,v2; 313c4762a1bSJed Brown AppCtx* user; 314c4762a1bSJed Brown PetscViewer fd; 315c4762a1bSJed Brown char resultFile[] = "tomographyResult_x"; 316c4762a1bSJed Brown 317c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 318c4762a1bSJed Brown ierr = PetscNew(&user);CHKERRQ(ierr); 319c4762a1bSJed Brown ierr = InitializeUserData(user);CHKERRQ(ierr); 320c4762a1bSJed Brown 321c4762a1bSJed Brown ierr = TaoCreate(PETSC_COMM_WORLD, &tao);CHKERRQ(ierr); 322c4762a1bSJed Brown ierr = TaoSetType(tao, TAOADMM);CHKERRQ(ierr); 323c4762a1bSJed Brown ierr = TaoSetInitialVector(tao, user->x);CHKERRQ(ierr); 324c4762a1bSJed Brown /* f(x) + g(x) for parent tao */ 325c4762a1bSJed Brown ierr = TaoADMMSetSpectralPenalty(tao,1.);CHKERRQ(ierr); 326c4762a1bSJed Brown ierr = TaoSetObjectiveAndGradientRoutine(tao, FullObjGrad, (void*)user);CHKERRQ(ierr); 327c4762a1bSJed Brown ierr = MatShift(user->HF,user->lambda);CHKERRQ(ierr); 328c4762a1bSJed Brown ierr = TaoSetHessianRoutine(tao, user->HF, user->HF, HessianFull, (void*)user);CHKERRQ(ierr); 329c4762a1bSJed Brown 330c4762a1bSJed Brown /* f(x) for misfit tao */ 331c4762a1bSJed Brown ierr = TaoADMMSetMisfitObjectiveAndGradientRoutine(tao, MisfitObjectiveAndGradient, (void*)user);CHKERRQ(ierr); 332c4762a1bSJed Brown ierr = TaoADMMSetMisfitHessianRoutine(tao, user->Hx, user->Hx, HessianMisfit, (void*)user);CHKERRQ(ierr); 333c4762a1bSJed Brown ierr = TaoADMMSetMisfitHessianChangeStatus(tao,PETSC_FALSE);CHKERRQ(ierr); 334c4762a1bSJed Brown ierr = TaoADMMSetMisfitConstraintJacobian(tao,user->D,user->D,NullJacobian,(void*)user);CHKERRQ(ierr); 335c4762a1bSJed Brown 336c4762a1bSJed Brown /* g(x) for regularizer tao */ 337c4762a1bSJed Brown if (user->reg == 1) { 338c4762a1bSJed Brown ierr = TaoADMMSetRegularizerObjectiveAndGradientRoutine(tao, RegularizerObjectiveAndGradient1, (void*)user);CHKERRQ(ierr); 339c4762a1bSJed Brown ierr = TaoADMMSetRegularizerHessianRoutine(tao, user->Hz, user->Hz, HessianReg, (void*)user);CHKERRQ(ierr); 340c4762a1bSJed Brown ierr = TaoADMMSetRegHessianChangeStatus(tao,PETSC_TRUE);CHKERRQ(ierr); 341c4762a1bSJed Brown } else if (user->reg == 2) { 342c4762a1bSJed Brown ierr = TaoADMMSetRegularizerObjectiveAndGradientRoutine(tao, RegularizerObjectiveAndGradient2, (void*)user);CHKERRQ(ierr); 343c4762a1bSJed Brown ierr = MatShift(user->Hz,1);CHKERRQ(ierr); 344c4762a1bSJed Brown ierr = MatScale(user->Hz,user->lambda);CHKERRQ(ierr); 345c4762a1bSJed Brown ierr = TaoADMMSetRegularizerHessianRoutine(tao, user->Hz, user->Hz, HessianMisfit, (void*)user);CHKERRQ(ierr); 346c4762a1bSJed Brown ierr = TaoADMMSetRegHessianChangeStatus(tao,PETSC_TRUE);CHKERRQ(ierr); 347d8185827SBarry Smith } else if (user->reg != 3) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_UNKNOWN_TYPE, "Incorrect Reg type"); /* TaoShell case */ 348c4762a1bSJed Brown 349c4762a1bSJed Brown /* Set type for the misfit solver */ 350c4762a1bSJed Brown ierr = TaoADMMGetMisfitSubsolver(tao, &misfit);CHKERRQ(ierr); 351c4762a1bSJed Brown ierr = TaoADMMGetRegularizationSubsolver(tao, ®);CHKERRQ(ierr); 352c4762a1bSJed Brown ierr = TaoSetType(misfit,TAONLS);CHKERRQ(ierr); 353c4762a1bSJed Brown if (user->reg == 3) { 354c4762a1bSJed Brown ierr = TaoSetType(reg,TAOSHELL);CHKERRQ(ierr); 355c4762a1bSJed Brown ierr = TaoShellSetContext(reg, (void*) user);CHKERRQ(ierr); 356c4762a1bSJed Brown ierr = TaoShellSetSolve(reg, TaoShellSolve_SoftThreshold);CHKERRQ(ierr); 357c4762a1bSJed Brown } else { 358c4762a1bSJed Brown ierr = TaoSetType(reg,TAONLS);CHKERRQ(ierr); 359c4762a1bSJed Brown } 360c4762a1bSJed Brown ierr = TaoSetVariableBounds(misfit,user->xlb,user->xub);CHKERRQ(ierr); 361c4762a1bSJed Brown 362c4762a1bSJed Brown /* Soft Thresholding solves the ADMM problem with the L1 regularizer lambda*||z||_1 and the x-z=0 constraint */ 363c4762a1bSJed Brown ierr = TaoADMMSetRegularizerCoefficient(tao, user->lambda);CHKERRQ(ierr); 364c4762a1bSJed Brown ierr = TaoADMMSetRegularizerConstraintJacobian(tao,NULL,NULL,NullJacobian,(void*)user);CHKERRQ(ierr); 365c4762a1bSJed Brown ierr = TaoADMMSetMinimumSpectralPenalty(tao,user->mumin);CHKERRQ(ierr); 366c4762a1bSJed Brown 367c4762a1bSJed Brown ierr = TaoADMMSetConstraintVectorRHS(tao,user->c);CHKERRQ(ierr); 368c4762a1bSJed Brown ierr = TaoSetFromOptions(tao);CHKERRQ(ierr); 369c4762a1bSJed Brown ierr = TaoSolve(tao);CHKERRQ(ierr); 370c4762a1bSJed Brown 371c4762a1bSJed Brown /* Save x (reconstruction of object) vector to a binary file, which maybe read from Matlab and convert to a 2D image for comparison. */ 372c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,resultFile,FILE_MODE_WRITE,&fd);CHKERRQ(ierr); 373c4762a1bSJed Brown ierr = VecView(user->x,fd);CHKERRQ(ierr); 374c4762a1bSJed Brown ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr); 375c4762a1bSJed Brown 376c4762a1bSJed Brown /* compute the error */ 377c4762a1bSJed Brown ierr = VecAXPY(user->x,-1,user->xGT);CHKERRQ(ierr); 378c4762a1bSJed Brown ierr = VecNorm(user->x,NORM_2,&v1);CHKERRQ(ierr); 379c4762a1bSJed Brown ierr = VecNorm(user->xGT,NORM_2,&v2);CHKERRQ(ierr); 380c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "relative reconstruction error: ||x-xGT||/||xGT|| = %6.4e.\n", (double)(v1/v2));CHKERRQ(ierr); 381c4762a1bSJed Brown 382c4762a1bSJed Brown /* Free TAO data structures */ 383c4762a1bSJed Brown ierr = TaoDestroy(&tao);CHKERRQ(ierr); 384c4762a1bSJed Brown ierr = DestroyContext(user);CHKERRQ(ierr); 385c4762a1bSJed Brown ierr = PetscFree(user);CHKERRQ(ierr); 386c4762a1bSJed Brown ierr = PetscFinalize(); 387c4762a1bSJed Brown return ierr; 388c4762a1bSJed Brown } 389c4762a1bSJed Brown 390c4762a1bSJed Brown /*TEST 391c4762a1bSJed Brown 392c4762a1bSJed Brown build: 393*dfd57a17SPierre Jolivet requires: !complex !single !__float128 !defined(PETSC_USE_64BIT_INDICES) 394c4762a1bSJed Brown 395c4762a1bSJed Brown test: 396c4762a1bSJed Brown suffix: 1 397c4762a1bSJed Brown localrunfiles: tomographyData_A_b_xGT 398c4762a1bSJed Brown args: -lambda 1.e-8 -tao_monitor -tao_type nls -tao_nls_pc_type icc 399c4762a1bSJed Brown 400c4762a1bSJed Brown test: 401c4762a1bSJed Brown suffix: 2 402c4762a1bSJed Brown localrunfiles: tomographyData_A_b_xGT 403c4762a1bSJed 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 404c4762a1bSJed Brown 405c4762a1bSJed Brown test: 406c4762a1bSJed Brown suffix: 3 407c4762a1bSJed Brown localrunfiles: tomographyData_A_b_xGT 408c4762a1bSJed 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 409c4762a1bSJed Brown 410c4762a1bSJed Brown test: 411c4762a1bSJed Brown suffix: 4 412c4762a1bSJed Brown localrunfiles: tomographyData_A_b_xGT 413c4762a1bSJed 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 414c4762a1bSJed Brown 415c4762a1bSJed Brown test: 416c4762a1bSJed Brown suffix: 5 417c4762a1bSJed Brown localrunfiles: tomographyData_A_b_xGT 418c4762a1bSJed 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 419c4762a1bSJed Brown 420c4762a1bSJed Brown test: 421c4762a1bSJed Brown suffix: 6 422c4762a1bSJed Brown localrunfiles: tomographyData_A_b_xGT 423c4762a1bSJed 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 424c4762a1bSJed Brown 425c4762a1bSJed Brown TEST*/ 426