1367daffbSBarry Smith #include <petscmat.h> /*I "mat.h" I*/ 2e0877f53SBarry Smith 30ebee16dSLisandro Dalcin PETSC_INTERN PetscErrorCode MatCreateADA(Mat,Vec, Vec, Mat*); 40ebee16dSLisandro Dalcin 5e0877f53SBarry Smith typedef struct { 6e0877f53SBarry Smith Mat A; 7e0877f53SBarry Smith Vec D1; 8e0877f53SBarry Smith Vec D2; 9e0877f53SBarry Smith Vec W; 10e0877f53SBarry Smith Vec W2; 11e0877f53SBarry Smith Vec ADADiag; 12e0877f53SBarry Smith PetscInt GotDiag; 13e0877f53SBarry Smith } _p_TaoMatADACtx; 14e0877f53SBarry Smith typedef _p_TaoMatADACtx* TaoMatADACtx; 15e0877f53SBarry Smith 16e0877f53SBarry Smith static PetscErrorCode MatMult_ADA(Mat mat,Vec a,Vec y) 17e0877f53SBarry Smith { 18e0877f53SBarry Smith TaoMatADACtx ctx; 19e0877f53SBarry Smith PetscReal one = 1.0; 20e0877f53SBarry Smith PetscErrorCode ierr; 21e0877f53SBarry Smith 22e0877f53SBarry Smith PetscFunctionBegin; 233ec1f749SStefano Zampini ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 24e0877f53SBarry Smith ierr = MatMult(ctx->A,a,ctx->W);CHKERRQ(ierr); 25e0877f53SBarry Smith if (ctx->D1) { 26e0877f53SBarry Smith ierr = VecPointwiseMult(ctx->W,ctx->D1,ctx->W);CHKERRQ(ierr); 27e0877f53SBarry Smith } 28e0877f53SBarry Smith ierr = MatMultTranspose(ctx->A,ctx->W,y);CHKERRQ(ierr); 29e0877f53SBarry Smith if (ctx->D2) { 30e0877f53SBarry Smith ierr = VecPointwiseMult(ctx->W2, ctx->D2, a);CHKERRQ(ierr); 31e0877f53SBarry Smith ierr = VecAXPY(y, one, ctx->W2);CHKERRQ(ierr); 32e0877f53SBarry Smith } 33e0877f53SBarry Smith PetscFunctionReturn(0); 34e0877f53SBarry Smith } 35e0877f53SBarry Smith 36e0877f53SBarry Smith static PetscErrorCode MatMultTranspose_ADA(Mat mat,Vec a,Vec y) 37e0877f53SBarry Smith { 38e0877f53SBarry Smith PetscErrorCode ierr; 39e0877f53SBarry Smith 40e0877f53SBarry Smith PetscFunctionBegin; 41e0877f53SBarry Smith ierr = MatMult_ADA(mat,a,y);CHKERRQ(ierr); 42e0877f53SBarry Smith PetscFunctionReturn(0); 43e0877f53SBarry Smith } 44e0877f53SBarry Smith 450c0fd78eSBarry Smith static PetscErrorCode MatDiagonalSet_ADA(Mat M,Vec D, InsertMode mode) 46e0877f53SBarry Smith { 47e0877f53SBarry Smith TaoMatADACtx ctx; 48e0877f53SBarry Smith PetscReal zero=0.0,one = 1.0; 49e0877f53SBarry Smith PetscErrorCode ierr; 50e0877f53SBarry Smith 51e0877f53SBarry Smith PetscFunctionBegin; 52*3c859ba3SBarry Smith PetscCheck(mode != INSERT_VALUES,PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"Cannot insert diagonal entries of this matrix type, can only add"); 533ec1f749SStefano Zampini ierr = MatShellGetContext(M,&ctx);CHKERRQ(ierr); 54e0877f53SBarry Smith if (!ctx->D2) { 55e0877f53SBarry Smith ierr = VecDuplicate(D,&ctx->D2);CHKERRQ(ierr); 56e0877f53SBarry Smith ierr = VecSet(ctx->D2, zero);CHKERRQ(ierr); 57e0877f53SBarry Smith } 58e0877f53SBarry Smith ierr = VecAXPY(ctx->D2, one, D);CHKERRQ(ierr); 59e0877f53SBarry Smith PetscFunctionReturn(0); 60e0877f53SBarry Smith } 61e0877f53SBarry Smith 62e0877f53SBarry Smith static PetscErrorCode MatDestroy_ADA(Mat mat) 63e0877f53SBarry Smith { 64e0877f53SBarry Smith PetscErrorCode ierr; 65e0877f53SBarry Smith TaoMatADACtx ctx; 66e0877f53SBarry Smith 67e0877f53SBarry Smith PetscFunctionBegin; 683ec1f749SStefano Zampini ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 69e0877f53SBarry Smith ierr = VecDestroy(&ctx->W);CHKERRQ(ierr); 70e0877f53SBarry Smith ierr = VecDestroy(&ctx->W2);CHKERRQ(ierr); 71e0877f53SBarry Smith ierr = VecDestroy(&ctx->ADADiag);CHKERRQ(ierr); 72e0877f53SBarry Smith ierr = MatDestroy(&ctx->A);CHKERRQ(ierr); 73e0877f53SBarry Smith ierr = VecDestroy(&ctx->D1);CHKERRQ(ierr); 74e0877f53SBarry Smith ierr = VecDestroy(&ctx->D2);CHKERRQ(ierr); 75e0877f53SBarry Smith ierr = PetscFree(ctx);CHKERRQ(ierr); 76e0877f53SBarry Smith PetscFunctionReturn(0); 77e0877f53SBarry Smith } 78e0877f53SBarry Smith 79e0877f53SBarry Smith static PetscErrorCode MatView_ADA(Mat mat,PetscViewer viewer) 80e0877f53SBarry Smith { 81e0877f53SBarry Smith PetscFunctionBegin; 82e0877f53SBarry Smith PetscFunctionReturn(0); 83e0877f53SBarry Smith } 84e0877f53SBarry Smith 85e0877f53SBarry Smith static PetscErrorCode MatShift_ADA(Mat Y, PetscReal a) 86e0877f53SBarry Smith { 87e0877f53SBarry Smith PetscErrorCode ierr; 88e0877f53SBarry Smith TaoMatADACtx ctx; 89e0877f53SBarry Smith 90e0877f53SBarry Smith PetscFunctionBegin; 913ec1f749SStefano Zampini ierr = MatShellGetContext(Y,&ctx);CHKERRQ(ierr); 92e0877f53SBarry Smith ierr = VecShift(ctx->D2,a);CHKERRQ(ierr); 93e0877f53SBarry Smith PetscFunctionReturn(0); 94e0877f53SBarry Smith } 95e0877f53SBarry Smith 96e0877f53SBarry Smith static PetscErrorCode MatDuplicate_ADA(Mat mat,MatDuplicateOption op,Mat *M) 97e0877f53SBarry Smith { 98e0877f53SBarry Smith PetscErrorCode ierr; 99e0877f53SBarry Smith TaoMatADACtx ctx; 100e0877f53SBarry Smith Mat A2; 101e0877f53SBarry Smith Vec D1b=NULL,D2b; 102e0877f53SBarry Smith 103e0877f53SBarry Smith PetscFunctionBegin; 1043ec1f749SStefano Zampini ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 105e0877f53SBarry Smith ierr = MatDuplicate(ctx->A,op,&A2);CHKERRQ(ierr); 106e0877f53SBarry Smith if (ctx->D1) { 107e0877f53SBarry Smith ierr = VecDuplicate(ctx->D1,&D1b);CHKERRQ(ierr); 108e0877f53SBarry Smith ierr = VecCopy(ctx->D1,D1b);CHKERRQ(ierr); 109e0877f53SBarry Smith } 110e0877f53SBarry Smith ierr = VecDuplicate(ctx->D2,&D2b);CHKERRQ(ierr); 111e0877f53SBarry Smith ierr = VecCopy(ctx->D2,D2b);CHKERRQ(ierr); 112e0877f53SBarry Smith ierr = MatCreateADA(A2,D1b,D2b,M);CHKERRQ(ierr); 113e0877f53SBarry Smith if (ctx->D1) { 114e0877f53SBarry Smith ierr = PetscObjectDereference((PetscObject)D1b);CHKERRQ(ierr); 115e0877f53SBarry Smith } 116e0877f53SBarry Smith ierr = PetscObjectDereference((PetscObject)D2b);CHKERRQ(ierr); 117e0877f53SBarry Smith ierr = PetscObjectDereference((PetscObject)A2);CHKERRQ(ierr); 118e0877f53SBarry Smith PetscFunctionReturn(0); 119e0877f53SBarry Smith } 120e0877f53SBarry Smith 121e0877f53SBarry Smith static PetscErrorCode MatEqual_ADA(Mat A,Mat B,PetscBool *flg) 122e0877f53SBarry Smith { 123e0877f53SBarry Smith PetscErrorCode ierr; 124e0877f53SBarry Smith TaoMatADACtx ctx1,ctx2; 125e0877f53SBarry Smith 126e0877f53SBarry Smith PetscFunctionBegin; 1273ec1f749SStefano Zampini ierr = MatShellGetContext(A,&ctx1);CHKERRQ(ierr); 1283ec1f749SStefano Zampini ierr = MatShellGetContext(B,&ctx2);CHKERRQ(ierr); 129e0877f53SBarry Smith ierr = VecEqual(ctx1->D2,ctx2->D2,flg);CHKERRQ(ierr); 130e0877f53SBarry Smith if (*flg==PETSC_TRUE) { 131e0877f53SBarry Smith ierr = VecEqual(ctx1->D1,ctx2->D1,flg);CHKERRQ(ierr); 132e0877f53SBarry Smith } 133e0877f53SBarry Smith if (*flg==PETSC_TRUE) { 134e0877f53SBarry Smith ierr = MatEqual(ctx1->A,ctx2->A,flg);CHKERRQ(ierr); 135e0877f53SBarry Smith } 136e0877f53SBarry Smith PetscFunctionReturn(0); 137e0877f53SBarry Smith } 138e0877f53SBarry Smith 139e0877f53SBarry Smith static PetscErrorCode MatScale_ADA(Mat mat, PetscReal a) 140e0877f53SBarry Smith { 141e0877f53SBarry Smith PetscErrorCode ierr; 142e0877f53SBarry Smith TaoMatADACtx ctx; 143e0877f53SBarry Smith 144e0877f53SBarry Smith PetscFunctionBegin; 1453ec1f749SStefano Zampini ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 146e0877f53SBarry Smith ierr = VecScale(ctx->D1,a);CHKERRQ(ierr); 147e0877f53SBarry Smith if (ctx->D2) { 148e0877f53SBarry Smith ierr = VecScale(ctx->D2,a);CHKERRQ(ierr); 149e0877f53SBarry Smith } 150e0877f53SBarry Smith PetscFunctionReturn(0); 151e0877f53SBarry Smith } 152e0877f53SBarry Smith 153cf37664fSBarry Smith static PetscErrorCode MatTranspose_ADA(Mat mat,MatReuse reuse,Mat *B) 154e0877f53SBarry Smith { 155e0877f53SBarry Smith PetscErrorCode ierr; 156e0877f53SBarry Smith TaoMatADACtx ctx; 157e0877f53SBarry Smith 158e0877f53SBarry Smith PetscFunctionBegin; 1593ec1f749SStefano Zampini ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 160cf37664fSBarry Smith if (reuse == MAT_INITIAL_MATRIX) { 161e0877f53SBarry Smith ierr = MatDuplicate(mat,MAT_COPY_VALUES,B);CHKERRQ(ierr); 162cf37664fSBarry Smith } else if (reuse == MAT_REUSE_MATRIX) { 163cf37664fSBarry Smith ierr = MatCopy(mat,*B,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 164cf37664fSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Does not support inplace transpose"); 165e0877f53SBarry Smith PetscFunctionReturn(0); 166e0877f53SBarry Smith } 167e0877f53SBarry Smith 1680c0fd78eSBarry Smith static PetscErrorCode MatADAComputeDiagonal(Mat mat) 169e0877f53SBarry Smith { 170e0877f53SBarry Smith PetscErrorCode ierr; 171e0877f53SBarry Smith PetscInt i,m,n,low,high; 172e0877f53SBarry Smith PetscScalar *dtemp,*dptr; 173e0877f53SBarry Smith TaoMatADACtx ctx; 174e0877f53SBarry Smith 175e0877f53SBarry Smith PetscFunctionBegin; 1763ec1f749SStefano Zampini ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 177e0877f53SBarry Smith ierr = MatGetOwnershipRange(mat, &low, &high);CHKERRQ(ierr); 178e0877f53SBarry Smith ierr = MatGetSize(mat,&m,&n);CHKERRQ(ierr); 179e0877f53SBarry Smith 180e0877f53SBarry Smith ierr = PetscMalloc1(n,&dtemp);CHKERRQ(ierr); 181e0877f53SBarry Smith for (i=0; i<n; i++) { 182e0877f53SBarry Smith ierr = MatGetColumnVector(ctx->A, ctx->W, i);CHKERRQ(ierr); 183e0877f53SBarry Smith ierr = VecPointwiseMult(ctx->W,ctx->W,ctx->W);CHKERRQ(ierr); 184e0877f53SBarry Smith ierr = VecDotBegin(ctx->D1, ctx->W,dtemp+i);CHKERRQ(ierr); 185e0877f53SBarry Smith } 186e0877f53SBarry Smith for (i=0; i<n; i++) { 187e0877f53SBarry Smith ierr = VecDotEnd(ctx->D1, ctx->W,dtemp+i);CHKERRQ(ierr); 188e0877f53SBarry Smith } 189e0877f53SBarry Smith 190e0877f53SBarry Smith ierr = VecGetArray(ctx->ADADiag,&dptr);CHKERRQ(ierr); 191e0877f53SBarry Smith for (i=low; i<high; i++) { 192e0877f53SBarry Smith dptr[i-low]= dtemp[i]; 193e0877f53SBarry Smith } 194e0877f53SBarry Smith ierr = VecRestoreArray(ctx->ADADiag,&dptr);CHKERRQ(ierr); 195e0877f53SBarry Smith ierr = PetscFree(dtemp);CHKERRQ(ierr); 196e0877f53SBarry Smith PetscFunctionReturn(0); 197e0877f53SBarry Smith } 198e0877f53SBarry Smith 199e0877f53SBarry Smith static PetscErrorCode MatGetDiagonal_ADA(Mat mat,Vec v) 200e0877f53SBarry Smith { 201e0877f53SBarry Smith PetscErrorCode ierr; 202e0877f53SBarry Smith PetscReal one=1.0; 203e0877f53SBarry Smith TaoMatADACtx ctx; 204e0877f53SBarry Smith 205e0877f53SBarry Smith PetscFunctionBegin; 2063ec1f749SStefano Zampini ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 207e0877f53SBarry Smith ierr = MatADAComputeDiagonal(mat);CHKERRQ(ierr); 208e0877f53SBarry Smith ierr = VecCopy(ctx->ADADiag,v);CHKERRQ(ierr); 209e0877f53SBarry Smith if (ctx->D2) { 210e0877f53SBarry Smith ierr = VecAXPY(v, one, ctx->D2);CHKERRQ(ierr); 211e0877f53SBarry Smith } 212e0877f53SBarry Smith PetscFunctionReturn(0); 213e0877f53SBarry Smith } 214e0877f53SBarry Smith 2157dae84e0SHong Zhang static PetscErrorCode MatCreateSubMatrix_ADA(Mat mat,IS isrow,IS iscol,MatReuse cll, Mat *newmat) 216e0877f53SBarry Smith { 217e0877f53SBarry Smith PetscErrorCode ierr; 218e0877f53SBarry Smith PetscInt low,high; 219e0877f53SBarry Smith IS ISrow; 220e0877f53SBarry Smith Vec D1,D2; 221e0877f53SBarry Smith Mat Atemp; 222e0877f53SBarry Smith TaoMatADACtx ctx; 223e0877f53SBarry Smith PetscBool isequal; 224e0877f53SBarry Smith 225e0877f53SBarry Smith PetscFunctionBegin; 226e0877f53SBarry Smith ierr = ISEqual(isrow,iscol,&isequal);CHKERRQ(ierr); 227*3c859ba3SBarry Smith PetscCheck(isequal,PETSC_COMM_SELF,PETSC_ERR_SUP,"Only for identical column and row indices"); 2283ec1f749SStefano Zampini ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 229e0877f53SBarry Smith 230e0877f53SBarry Smith ierr = MatGetOwnershipRange(ctx->A,&low,&high);CHKERRQ(ierr); 231367daffbSBarry Smith ierr = ISCreateStride(PetscObjectComm((PetscObject)mat),high-low,low,1,&ISrow);CHKERRQ(ierr); 2327dae84e0SHong Zhang ierr = MatCreateSubMatrix(ctx->A,ISrow,iscol,cll,&Atemp);CHKERRQ(ierr); 233e0877f53SBarry Smith ierr = ISDestroy(&ISrow);CHKERRQ(ierr); 234e0877f53SBarry Smith 235e0877f53SBarry Smith if (ctx->D1) { 236e0877f53SBarry Smith ierr = VecDuplicate(ctx->D1,&D1);CHKERRQ(ierr); 237e0877f53SBarry Smith ierr = VecCopy(ctx->D1,D1);CHKERRQ(ierr); 238e0877f53SBarry Smith } else { 239e0877f53SBarry Smith D1 = NULL; 240e0877f53SBarry Smith } 241e0877f53SBarry Smith 242e0877f53SBarry Smith if (ctx->D2) { 243e0877f53SBarry Smith Vec D2sub; 244e0877f53SBarry Smith 245e0877f53SBarry Smith ierr = VecGetSubVector(ctx->D2,isrow,&D2sub);CHKERRQ(ierr); 246e0877f53SBarry Smith ierr = VecDuplicate(D2sub,&D2);CHKERRQ(ierr); 247e0877f53SBarry Smith ierr = VecCopy(D2sub,D2);CHKERRQ(ierr); 248e0877f53SBarry Smith ierr = VecRestoreSubVector(ctx->D2,isrow,&D2sub);CHKERRQ(ierr); 249e0877f53SBarry Smith } else { 250e0877f53SBarry Smith D2 = NULL; 251e0877f53SBarry Smith } 252e0877f53SBarry Smith 253e0877f53SBarry Smith ierr = MatCreateADA(Atemp,D1,D2,newmat);CHKERRQ(ierr); 2543ec1f749SStefano Zampini ierr = MatShellGetContext(*newmat,&ctx);CHKERRQ(ierr); 255e0877f53SBarry Smith ierr = PetscObjectDereference((PetscObject)Atemp);CHKERRQ(ierr); 256e0877f53SBarry Smith if (ctx->D1) { 257e0877f53SBarry Smith ierr = PetscObjectDereference((PetscObject)D1);CHKERRQ(ierr); 258e0877f53SBarry Smith } 259e0877f53SBarry Smith if (ctx->D2) { 260e0877f53SBarry Smith ierr = PetscObjectDereference((PetscObject)D2);CHKERRQ(ierr); 261e0877f53SBarry Smith } 262e0877f53SBarry Smith PetscFunctionReturn(0); 263e0877f53SBarry Smith } 264e0877f53SBarry Smith 2657dae84e0SHong Zhang static PetscErrorCode MatCreateSubMatrices_ADA(Mat A,PetscInt n, IS *irow,IS *icol,MatReuse scall,Mat **B) 266e0877f53SBarry Smith { 267e0877f53SBarry Smith PetscErrorCode ierr; 268e0877f53SBarry Smith PetscInt i; 269e0877f53SBarry Smith 270e0877f53SBarry Smith PetscFunctionBegin; 271e0877f53SBarry Smith if (scall == MAT_INITIAL_MATRIX) { 272df750dc8SHong Zhang ierr = PetscCalloc1(n+1,B);CHKERRQ(ierr); 273e0877f53SBarry Smith } 274e0877f53SBarry Smith for (i=0; i<n; i++) { 2757dae84e0SHong Zhang ierr = MatCreateSubMatrix_ADA(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr); 276e0877f53SBarry Smith } 277e0877f53SBarry Smith PetscFunctionReturn(0); 278e0877f53SBarry Smith } 279e0877f53SBarry Smith 280e0877f53SBarry Smith static PetscErrorCode MatGetColumnVector_ADA(Mat mat,Vec Y, PetscInt col) 281e0877f53SBarry Smith { 282e0877f53SBarry Smith PetscErrorCode ierr; 283e0877f53SBarry Smith PetscInt low,high; 284e0877f53SBarry Smith PetscScalar zero=0.0,one=1.0; 285e0877f53SBarry Smith 286e0877f53SBarry Smith PetscFunctionBegin; 287e0877f53SBarry Smith ierr = VecSet(Y, zero);CHKERRQ(ierr); 288e0877f53SBarry Smith ierr = VecGetOwnershipRange(Y,&low,&high);CHKERRQ(ierr); 289e0877f53SBarry Smith if (col>=low && col<high) { 290e0877f53SBarry Smith ierr = VecSetValue(Y,col,one,INSERT_VALUES);CHKERRQ(ierr); 291e0877f53SBarry Smith } 292e0877f53SBarry Smith ierr = VecAssemblyBegin(Y);CHKERRQ(ierr); 293e0877f53SBarry Smith ierr = VecAssemblyEnd(Y);CHKERRQ(ierr); 294e0877f53SBarry Smith ierr = MatMult_ADA(mat,Y,Y);CHKERRQ(ierr); 295e0877f53SBarry Smith PetscFunctionReturn(0); 296e0877f53SBarry Smith } 297e0877f53SBarry Smith 298cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_ADA(Mat mat,MatType newtype,Mat *NewMat) 299e0877f53SBarry Smith { 300e0877f53SBarry Smith PetscErrorCode ierr; 301e0877f53SBarry Smith PetscMPIInt size; 302e0877f53SBarry Smith PetscBool sametype, issame, isdense, isseqdense; 303e0877f53SBarry Smith TaoMatADACtx ctx; 304e0877f53SBarry Smith 305e0877f53SBarry Smith PetscFunctionBegin; 3063ec1f749SStefano Zampini ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 307ffc4695bSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRMPI(ierr); 308e0877f53SBarry Smith 309e0877f53SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)mat,newtype,&sametype);CHKERRQ(ierr); 310e0877f53SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)mat,MATSAME,&issame);CHKERRQ(ierr); 311e0877f53SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPIDENSE,&isdense);CHKERRQ(ierr); 312e0877f53SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQDENSE,&isseqdense);CHKERRQ(ierr); 313e0877f53SBarry Smith 314e0877f53SBarry Smith if (sametype || issame) { 315e0877f53SBarry Smith ierr = MatDuplicate(mat,MAT_COPY_VALUES,NewMat);CHKERRQ(ierr); 316e0877f53SBarry Smith } else if (isdense) { 317e0877f53SBarry Smith PetscInt i,j,low,high,m,n,M,N; 318e0877f53SBarry Smith const PetscScalar *dptr; 319e0877f53SBarry Smith Vec X; 320e0877f53SBarry Smith 321e0877f53SBarry Smith ierr = VecDuplicate(ctx->D2,&X);CHKERRQ(ierr); 322e0877f53SBarry Smith ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 323e0877f53SBarry Smith ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 324367daffbSBarry Smith ierr = MatCreateDense(PetscObjectComm((PetscObject)mat),m,m,N,N,NULL,NewMat);CHKERRQ(ierr); 325e0877f53SBarry Smith ierr = MatGetOwnershipRange(*NewMat,&low,&high);CHKERRQ(ierr); 326e0877f53SBarry Smith for (i=0;i<M;i++) { 327e0877f53SBarry Smith ierr = MatGetColumnVector_ADA(mat,X,i);CHKERRQ(ierr); 328e0877f53SBarry Smith ierr = VecGetArrayRead(X,&dptr);CHKERRQ(ierr); 329e0877f53SBarry Smith for (j=0; j<high-low; j++) { 330e0877f53SBarry Smith ierr = MatSetValue(*NewMat,low+j,i,dptr[j],INSERT_VALUES);CHKERRQ(ierr); 331e0877f53SBarry Smith } 332e0877f53SBarry Smith ierr = VecRestoreArrayRead(X,&dptr);CHKERRQ(ierr); 333e0877f53SBarry Smith } 334e0877f53SBarry Smith ierr = MatAssemblyBegin(*NewMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 335e0877f53SBarry Smith ierr = MatAssemblyEnd(*NewMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 336e0877f53SBarry Smith ierr = VecDestroy(&X);CHKERRQ(ierr); 337e0877f53SBarry Smith } else if (isseqdense && size==1) { 338e0877f53SBarry Smith PetscInt i,j,low,high,m,n,M,N; 339e0877f53SBarry Smith const PetscScalar *dptr; 340e0877f53SBarry Smith Vec X; 341e0877f53SBarry Smith 342e0877f53SBarry Smith ierr = VecDuplicate(ctx->D2,&X);CHKERRQ(ierr); 343e0877f53SBarry Smith ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 344e0877f53SBarry Smith ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 345367daffbSBarry Smith ierr = MatCreateSeqDense(PetscObjectComm((PetscObject)mat),N,N,NULL,NewMat);CHKERRQ(ierr); 346e0877f53SBarry Smith ierr = MatGetOwnershipRange(*NewMat,&low,&high);CHKERRQ(ierr); 347e0877f53SBarry Smith for (i=0;i<M;i++) { 348e0877f53SBarry Smith ierr = MatGetColumnVector_ADA(mat,X,i);CHKERRQ(ierr); 349e0877f53SBarry Smith ierr = VecGetArrayRead(X,&dptr);CHKERRQ(ierr); 350e0877f53SBarry Smith for (j=0; j<high-low; j++) { 351e0877f53SBarry Smith ierr = MatSetValue(*NewMat,low+j,i,dptr[j],INSERT_VALUES);CHKERRQ(ierr); 352e0877f53SBarry Smith } 353e0877f53SBarry Smith ierr = VecRestoreArrayRead(X,&dptr);CHKERRQ(ierr); 354e0877f53SBarry Smith } 355e0877f53SBarry Smith ierr = MatAssemblyBegin(*NewMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 356e0877f53SBarry Smith ierr = MatAssemblyEnd(*NewMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 357e0877f53SBarry Smith ierr = VecDestroy(&X);CHKERRQ(ierr); 358691b26d3SBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"No support to convert objects to that type"); 359e0877f53SBarry Smith PetscFunctionReturn(0); 360e0877f53SBarry Smith } 361e0877f53SBarry Smith 362e0877f53SBarry Smith static PetscErrorCode MatNorm_ADA(Mat mat,NormType type,PetscReal *norm) 363e0877f53SBarry Smith { 364e0877f53SBarry Smith PetscErrorCode ierr; 365e0877f53SBarry Smith TaoMatADACtx ctx; 366e0877f53SBarry Smith 367e0877f53SBarry Smith PetscFunctionBegin; 3683ec1f749SStefano Zampini ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 369e0877f53SBarry Smith if (type == NORM_FROBENIUS) { 370e0877f53SBarry Smith *norm = 1.0; 371e0877f53SBarry Smith } else if (type == NORM_1 || type == NORM_INFINITY) { 372e0877f53SBarry Smith *norm = 1.0; 373e0877f53SBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No two norm"); 374e0877f53SBarry Smith PetscFunctionReturn(0); 375e0877f53SBarry Smith } 376a7e14dcfSSatish Balay 377a7e14dcfSSatish Balay /*@C 378a7e14dcfSSatish Balay MatCreateADA - Creates a matrix M=A^T D1 A + D2 where D1, D2 are diagonal 379a7e14dcfSSatish Balay 380a7e14dcfSSatish Balay Collective on matrix 381a7e14dcfSSatish Balay 382a7e14dcfSSatish Balay Input Parameters: 383a7e14dcfSSatish Balay + mat - matrix of arbitrary type 3840c0fd78eSBarry Smith . d1 - A vector defining a diagonal matrix 3850c0fd78eSBarry Smith - d2 - A vector defining a diagonal matrix 386a7e14dcfSSatish Balay 387a7e14dcfSSatish Balay Output Parameters: 388a7e14dcfSSatish Balay . J - New matrix whose operations are defined in terms of mat, D1, and D2. 389a7e14dcfSSatish Balay 390a7e14dcfSSatish Balay Notes: 391a7e14dcfSSatish Balay The user provides the input data and is responsible for destroying 392a7e14dcfSSatish Balay this data after matrix J has been destroyed. 393a7e14dcfSSatish Balay 394a7e14dcfSSatish Balay Level: developer 395a7e14dcfSSatish Balay 396a7e14dcfSSatish Balay .seealso: MatCreate() 397a7e14dcfSSatish Balay @*/ 398a7e14dcfSSatish Balay PetscErrorCode MatCreateADA(Mat mat,Vec d1, Vec d2, Mat *J) 399a7e14dcfSSatish Balay { 400367daffbSBarry Smith MPI_Comm comm = PetscObjectComm((PetscObject)mat); 401a7e14dcfSSatish Balay TaoMatADACtx ctx; 402a7e14dcfSSatish Balay PetscErrorCode ierr; 403a7e14dcfSSatish Balay PetscInt nloc,n; 404a7e14dcfSSatish Balay 405a7e14dcfSSatish Balay PetscFunctionBegin; 4063c9e27cfSGeoffrey Irving ierr = PetscNew(&ctx);CHKERRQ(ierr); 407a7e14dcfSSatish Balay ctx->A=mat; 408a7e14dcfSSatish Balay ctx->D1=d1; 409a7e14dcfSSatish Balay ctx->D2=d2; 410a7e14dcfSSatish Balay if (d1) { 411a7e14dcfSSatish Balay ierr = VecDuplicate(d1,&ctx->W);CHKERRQ(ierr); 412a7e14dcfSSatish Balay ierr = PetscObjectReference((PetscObject)d1);CHKERRQ(ierr); 413a7e14dcfSSatish Balay } else { 4140e660641SBarry Smith ctx->W = NULL; 415a7e14dcfSSatish Balay } 416a7e14dcfSSatish Balay if (d2) { 417a7e14dcfSSatish Balay ierr = VecDuplicate(d2,&ctx->W2);CHKERRQ(ierr); 418a7e14dcfSSatish Balay ierr = VecDuplicate(d2,&ctx->ADADiag);CHKERRQ(ierr); 419a7e14dcfSSatish Balay ierr = PetscObjectReference((PetscObject)d2);CHKERRQ(ierr); 420a7e14dcfSSatish Balay } else { 4210e660641SBarry Smith ctx->W2 = NULL; 4220e660641SBarry Smith ctx->ADADiag = NULL; 423a7e14dcfSSatish Balay } 424a7e14dcfSSatish Balay 425a7e14dcfSSatish Balay ctx->GotDiag = 0; 426a7e14dcfSSatish Balay ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 427a7e14dcfSSatish Balay 428a7e14dcfSSatish Balay ierr = VecGetLocalSize(d2,&nloc);CHKERRQ(ierr); 429a7e14dcfSSatish Balay ierr = VecGetSize(d2,&n);CHKERRQ(ierr); 430a7e14dcfSSatish Balay 431a7e14dcfSSatish Balay ierr = MatCreateShell(comm,nloc,nloc,n,n,ctx,J);CHKERRQ(ierr); 4320c0fd78eSBarry Smith ierr = MatShellSetManageScalingShifts(*J);CHKERRQ(ierr); 433a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_MULT,(void(*)(void))MatMult_ADA);CHKERRQ(ierr); 434a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void(*)(void))MatDestroy_ADA);CHKERRQ(ierr); 435a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_VIEW,(void(*)(void))MatView_ADA);CHKERRQ(ierr); 436a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_ADA);CHKERRQ(ierr); 437a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_DIAGONAL_SET,(void(*)(void))MatDiagonalSet_ADA);CHKERRQ(ierr); 438a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_SHIFT,(void(*)(void))MatShift_ADA);CHKERRQ(ierr); 439a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_EQUAL,(void(*)(void))MatEqual_ADA);CHKERRQ(ierr); 440a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_SCALE,(void(*)(void))MatScale_ADA);CHKERRQ(ierr); 441a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_TRANSPOSE,(void(*)(void))MatTranspose_ADA);CHKERRQ(ierr); 442a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_ADA);CHKERRQ(ierr); 4437dae84e0SHong Zhang ierr = MatShellSetOperation(*J,MATOP_CREATE_SUBMATRICES,(void(*)(void))MatCreateSubMatrices_ADA);CHKERRQ(ierr); 444a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_NORM,(void(*)(void))MatNorm_ADA);CHKERRQ(ierr); 445a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_DUPLICATE,(void(*)(void))MatDuplicate_ADA);CHKERRQ(ierr); 4467dae84e0SHong Zhang ierr = MatShellSetOperation(*J,MATOP_CREATE_SUBMATRIX,(void(*)(void))MatCreateSubMatrix_ADA);CHKERRQ(ierr); 447a7e14dcfSSatish Balay 448a7e14dcfSSatish Balay ierr = PetscLogObjectParent((PetscObject)(*J),(PetscObject)ctx->W);CHKERRQ(ierr); 449a7e14dcfSSatish Balay ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)(*J));CHKERRQ(ierr); 450a7e14dcfSSatish Balay 451a7e14dcfSSatish Balay ierr = MatSetOption(*J,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr); 452a7e14dcfSSatish Balay PetscFunctionReturn(0); 453a7e14dcfSSatish Balay } 454