1*367daffbSBarry Smith #include <petscmat.h> /*I "mat.h" I*/ 2e0877f53SBarry Smith 3e0877f53SBarry Smith typedef struct{ 4e0877f53SBarry Smith Mat A; 5e0877f53SBarry Smith Vec D1; 6e0877f53SBarry Smith Vec D2; 7e0877f53SBarry Smith Vec W; 8e0877f53SBarry Smith Vec W2; 9e0877f53SBarry Smith Vec ADADiag; 10e0877f53SBarry Smith PetscInt GotDiag; 11e0877f53SBarry Smith } _p_TaoMatADACtx; 12e0877f53SBarry Smith typedef _p_TaoMatADACtx* TaoMatADACtx; 13e0877f53SBarry Smith 14e0877f53SBarry Smith extern PetscErrorCode MatCreateADA(Mat,Vec, Vec, Mat*); 15e0877f53SBarry Smith 16e0877f53SBarry Smith #undef __FUNCT__ 17e0877f53SBarry Smith #define __FUNCT__ "MatMult_ADA" 18e0877f53SBarry Smith static PetscErrorCode MatMult_ADA(Mat mat,Vec a,Vec y) 19e0877f53SBarry Smith { 20e0877f53SBarry Smith TaoMatADACtx ctx; 21e0877f53SBarry Smith PetscReal one = 1.0; 22e0877f53SBarry Smith PetscErrorCode ierr; 23e0877f53SBarry Smith 24e0877f53SBarry Smith PetscFunctionBegin; 25e0877f53SBarry Smith ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 26e0877f53SBarry Smith ierr = MatMult(ctx->A,a,ctx->W);CHKERRQ(ierr); 27e0877f53SBarry Smith if (ctx->D1){ 28e0877f53SBarry Smith ierr = VecPointwiseMult(ctx->W,ctx->D1,ctx->W);CHKERRQ(ierr); 29e0877f53SBarry Smith } 30e0877f53SBarry Smith ierr = MatMultTranspose(ctx->A,ctx->W,y);CHKERRQ(ierr); 31e0877f53SBarry Smith if (ctx->D2){ 32e0877f53SBarry Smith ierr = VecPointwiseMult(ctx->W2, ctx->D2, a);CHKERRQ(ierr); 33e0877f53SBarry Smith ierr = VecAXPY(y, one, ctx->W2);CHKERRQ(ierr); 34e0877f53SBarry Smith } 35e0877f53SBarry Smith PetscFunctionReturn(0); 36e0877f53SBarry Smith } 37e0877f53SBarry Smith 38e0877f53SBarry Smith #undef __FUNCT__ 39e0877f53SBarry Smith #define __FUNCT__ "MatMultTranspose_ADA" 40e0877f53SBarry Smith static PetscErrorCode MatMultTranspose_ADA(Mat mat,Vec a,Vec y) 41e0877f53SBarry Smith { 42e0877f53SBarry Smith PetscErrorCode ierr; 43e0877f53SBarry Smith 44e0877f53SBarry Smith PetscFunctionBegin; 45e0877f53SBarry Smith ierr = MatMult_ADA(mat,a,y);CHKERRQ(ierr); 46e0877f53SBarry Smith PetscFunctionReturn(0); 47e0877f53SBarry Smith } 48e0877f53SBarry Smith 49e0877f53SBarry Smith #undef __FUNCT__ 50e0877f53SBarry Smith #define __FUNCT__ "MatDiagonalSet_ADA" 51e0877f53SBarry Smith static PetscErrorCode MatDiagonalSet_ADA(Vec D, Mat M) 52e0877f53SBarry Smith { 53e0877f53SBarry Smith TaoMatADACtx ctx; 54e0877f53SBarry Smith PetscReal zero=0.0,one = 1.0; 55e0877f53SBarry Smith PetscErrorCode ierr; 56e0877f53SBarry Smith 57e0877f53SBarry Smith PetscFunctionBegin; 58e0877f53SBarry Smith ierr = MatShellGetContext(M,(void **)&ctx);CHKERRQ(ierr); 59e0877f53SBarry Smith if (!ctx->D2){ 60e0877f53SBarry Smith ierr = VecDuplicate(D,&ctx->D2);CHKERRQ(ierr); 61e0877f53SBarry Smith ierr = VecSet(ctx->D2, zero);CHKERRQ(ierr); 62e0877f53SBarry Smith } 63e0877f53SBarry Smith ierr = VecAXPY(ctx->D2, one, D);CHKERRQ(ierr); 64e0877f53SBarry Smith PetscFunctionReturn(0); 65e0877f53SBarry Smith } 66e0877f53SBarry Smith 67e0877f53SBarry Smith #undef __FUNCT__ 68e0877f53SBarry Smith #define __FUNCT__ "MatDestroy_ADA" 69e0877f53SBarry Smith static PetscErrorCode MatDestroy_ADA(Mat mat) 70e0877f53SBarry Smith { 71e0877f53SBarry Smith PetscErrorCode ierr; 72e0877f53SBarry Smith TaoMatADACtx ctx; 73e0877f53SBarry Smith 74e0877f53SBarry Smith PetscFunctionBegin; 75e0877f53SBarry Smith ierr=MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 76e0877f53SBarry Smith ierr=VecDestroy(&ctx->W);CHKERRQ(ierr); 77e0877f53SBarry Smith ierr=VecDestroy(&ctx->W2);CHKERRQ(ierr); 78e0877f53SBarry Smith ierr=VecDestroy(&ctx->ADADiag);CHKERRQ(ierr); 79e0877f53SBarry Smith ierr=MatDestroy(&ctx->A);CHKERRQ(ierr); 80e0877f53SBarry Smith ierr=VecDestroy(&ctx->D1);CHKERRQ(ierr); 81e0877f53SBarry Smith ierr=VecDestroy(&ctx->D2);CHKERRQ(ierr); 82e0877f53SBarry Smith ierr = PetscFree(ctx);CHKERRQ(ierr); 83e0877f53SBarry Smith PetscFunctionReturn(0); 84e0877f53SBarry Smith } 85e0877f53SBarry Smith 86e0877f53SBarry Smith #undef __FUNCT__ 87e0877f53SBarry Smith #define __FUNCT__ "MatView_ADA" 88e0877f53SBarry Smith static PetscErrorCode MatView_ADA(Mat mat,PetscViewer viewer) 89e0877f53SBarry Smith { 90e0877f53SBarry Smith PetscFunctionBegin; 91e0877f53SBarry Smith PetscFunctionReturn(0); 92e0877f53SBarry Smith } 93e0877f53SBarry Smith 94e0877f53SBarry Smith #undef __FUNCT__ 95e0877f53SBarry Smith #define __FUNCT__ "MatShift_ADA" 96e0877f53SBarry Smith static PetscErrorCode MatShift_ADA(Mat Y, PetscReal a) 97e0877f53SBarry Smith { 98e0877f53SBarry Smith PetscErrorCode ierr; 99e0877f53SBarry Smith TaoMatADACtx ctx; 100e0877f53SBarry Smith 101e0877f53SBarry Smith PetscFunctionBegin; 102e0877f53SBarry Smith ierr = MatShellGetContext(Y,(void **)&ctx);CHKERRQ(ierr); 103e0877f53SBarry Smith ierr = VecShift(ctx->D2,a);CHKERRQ(ierr); 104e0877f53SBarry Smith PetscFunctionReturn(0); 105e0877f53SBarry Smith } 106e0877f53SBarry Smith 107e0877f53SBarry Smith #undef __FUNCT__ 108e0877f53SBarry Smith #define __FUNCT__ "MatDuplicate_ADA" 109e0877f53SBarry Smith static PetscErrorCode MatDuplicate_ADA(Mat mat,MatDuplicateOption op,Mat *M) 110e0877f53SBarry Smith { 111e0877f53SBarry Smith PetscErrorCode ierr; 112e0877f53SBarry Smith TaoMatADACtx ctx; 113e0877f53SBarry Smith Mat A2; 114e0877f53SBarry Smith Vec D1b=NULL,D2b; 115e0877f53SBarry Smith 116e0877f53SBarry Smith PetscFunctionBegin; 117e0877f53SBarry Smith ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 118e0877f53SBarry Smith ierr = MatDuplicate(ctx->A,op,&A2);CHKERRQ(ierr); 119e0877f53SBarry Smith if (ctx->D1){ 120e0877f53SBarry Smith ierr = VecDuplicate(ctx->D1,&D1b);CHKERRQ(ierr); 121e0877f53SBarry Smith ierr = VecCopy(ctx->D1,D1b);CHKERRQ(ierr); 122e0877f53SBarry Smith } 123e0877f53SBarry Smith ierr = VecDuplicate(ctx->D2,&D2b);CHKERRQ(ierr); 124e0877f53SBarry Smith ierr = VecCopy(ctx->D2,D2b);CHKERRQ(ierr); 125e0877f53SBarry Smith ierr = MatCreateADA(A2,D1b,D2b,M);CHKERRQ(ierr); 126e0877f53SBarry Smith if (ctx->D1){ 127e0877f53SBarry Smith ierr = PetscObjectDereference((PetscObject)D1b);CHKERRQ(ierr); 128e0877f53SBarry Smith } 129e0877f53SBarry Smith ierr = PetscObjectDereference((PetscObject)D2b);CHKERRQ(ierr); 130e0877f53SBarry Smith ierr = PetscObjectDereference((PetscObject)A2);CHKERRQ(ierr); 131e0877f53SBarry Smith PetscFunctionReturn(0); 132e0877f53SBarry Smith } 133e0877f53SBarry Smith 134e0877f53SBarry Smith #undef __FUNCT__ 135e0877f53SBarry Smith #define __FUNCT__ "MatEqual_ADA" 136e0877f53SBarry Smith static PetscErrorCode MatEqual_ADA(Mat A,Mat B,PetscBool *flg) 137e0877f53SBarry Smith { 138e0877f53SBarry Smith PetscErrorCode ierr; 139e0877f53SBarry Smith TaoMatADACtx ctx1,ctx2; 140e0877f53SBarry Smith 141e0877f53SBarry Smith PetscFunctionBegin; 142e0877f53SBarry Smith ierr = MatShellGetContext(A,(void **)&ctx1);CHKERRQ(ierr); 143e0877f53SBarry Smith ierr = MatShellGetContext(B,(void **)&ctx2);CHKERRQ(ierr); 144e0877f53SBarry Smith ierr = VecEqual(ctx1->D2,ctx2->D2,flg);CHKERRQ(ierr); 145e0877f53SBarry Smith if (*flg==PETSC_TRUE){ 146e0877f53SBarry Smith ierr = VecEqual(ctx1->D1,ctx2->D1,flg);CHKERRQ(ierr); 147e0877f53SBarry Smith } 148e0877f53SBarry Smith if (*flg==PETSC_TRUE){ 149e0877f53SBarry Smith ierr = MatEqual(ctx1->A,ctx2->A,flg);CHKERRQ(ierr); 150e0877f53SBarry Smith } 151e0877f53SBarry Smith PetscFunctionReturn(0); 152e0877f53SBarry Smith } 153e0877f53SBarry Smith 154e0877f53SBarry Smith #undef __FUNCT__ 155e0877f53SBarry Smith #define __FUNCT__ "MatScale_ADA" 156e0877f53SBarry Smith static PetscErrorCode MatScale_ADA(Mat mat, PetscReal a) 157e0877f53SBarry Smith { 158e0877f53SBarry Smith PetscErrorCode ierr; 159e0877f53SBarry Smith TaoMatADACtx ctx; 160e0877f53SBarry Smith 161e0877f53SBarry Smith PetscFunctionBegin; 162e0877f53SBarry Smith ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 163e0877f53SBarry Smith ierr = VecScale(ctx->D1,a);CHKERRQ(ierr); 164e0877f53SBarry Smith if (ctx->D2){ 165e0877f53SBarry Smith ierr = VecScale(ctx->D2,a);CHKERRQ(ierr); 166e0877f53SBarry Smith } 167e0877f53SBarry Smith PetscFunctionReturn(0); 168e0877f53SBarry Smith } 169e0877f53SBarry Smith 170e0877f53SBarry Smith #undef __FUNCT__ 171e0877f53SBarry Smith #define __FUNCT__ "MatTranspose_ADA" 172e0877f53SBarry Smith static PetscErrorCode MatTranspose_ADA(Mat mat,Mat *B) 173e0877f53SBarry Smith { 174e0877f53SBarry Smith PetscErrorCode ierr; 175e0877f53SBarry Smith TaoMatADACtx ctx; 176e0877f53SBarry Smith 177e0877f53SBarry Smith PetscFunctionBegin; 178e0877f53SBarry Smith if (*B){ 179e0877f53SBarry Smith ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 180e0877f53SBarry Smith ierr = MatDuplicate(mat,MAT_COPY_VALUES,B);CHKERRQ(ierr); 181e0877f53SBarry Smith } 182e0877f53SBarry Smith PetscFunctionReturn(0); 183e0877f53SBarry Smith } 184e0877f53SBarry Smith 185e0877f53SBarry Smith #undef __FUNCT__ 186e0877f53SBarry Smith #define __FUNCT__ "MatADAComputeDiagonal" 187e0877f53SBarry Smith PetscErrorCode MatADAComputeDiagonal(Mat mat) 188e0877f53SBarry Smith { 189e0877f53SBarry Smith PetscErrorCode ierr; 190e0877f53SBarry Smith PetscInt i,m,n,low,high; 191e0877f53SBarry Smith PetscScalar *dtemp,*dptr; 192e0877f53SBarry Smith TaoMatADACtx ctx; 193e0877f53SBarry Smith 194e0877f53SBarry Smith PetscFunctionBegin; 195e0877f53SBarry Smith ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 196e0877f53SBarry Smith ierr = MatGetOwnershipRange(mat, &low, &high);CHKERRQ(ierr); 197e0877f53SBarry Smith ierr = MatGetSize(mat,&m,&n);CHKERRQ(ierr); 198e0877f53SBarry Smith 199e0877f53SBarry Smith ierr = PetscMalloc1(n,&dtemp );CHKERRQ(ierr); 200e0877f53SBarry Smith 201e0877f53SBarry Smith for (i=0; i<n; i++){ 202e0877f53SBarry Smith ierr = MatGetColumnVector(ctx->A, ctx->W, i);CHKERRQ(ierr); 203e0877f53SBarry Smith ierr = VecPointwiseMult(ctx->W,ctx->W,ctx->W);CHKERRQ(ierr); 204e0877f53SBarry Smith ierr = VecDotBegin(ctx->D1, ctx->W,dtemp+i);CHKERRQ(ierr); 205e0877f53SBarry Smith } 206e0877f53SBarry Smith for (i=0; i<n; i++){ 207e0877f53SBarry Smith ierr = VecDotEnd(ctx->D1, ctx->W,dtemp+i);CHKERRQ(ierr); 208e0877f53SBarry Smith } 209e0877f53SBarry Smith 210e0877f53SBarry Smith ierr = VecGetArray(ctx->ADADiag,&dptr);CHKERRQ(ierr); 211e0877f53SBarry Smith for (i=low; i<high; i++){ 212e0877f53SBarry Smith dptr[i-low]= dtemp[i]; 213e0877f53SBarry Smith } 214e0877f53SBarry Smith ierr = VecRestoreArray(ctx->ADADiag,&dptr);CHKERRQ(ierr); 215e0877f53SBarry Smith ierr = PetscFree(dtemp);CHKERRQ(ierr); 216e0877f53SBarry Smith PetscFunctionReturn(0); 217e0877f53SBarry Smith } 218e0877f53SBarry Smith 219e0877f53SBarry Smith #undef __FUNCT__ 220e0877f53SBarry Smith #define __FUNCT__ "MatGetDiagonal_ADA" 221e0877f53SBarry Smith static PetscErrorCode MatGetDiagonal_ADA(Mat mat,Vec v) 222e0877f53SBarry Smith { 223e0877f53SBarry Smith PetscErrorCode ierr; 224e0877f53SBarry Smith PetscReal one=1.0; 225e0877f53SBarry Smith TaoMatADACtx ctx; 226e0877f53SBarry Smith 227e0877f53SBarry Smith PetscFunctionBegin; 228e0877f53SBarry Smith ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 229e0877f53SBarry Smith ierr = MatADAComputeDiagonal(mat);CHKERRQ(ierr); 230e0877f53SBarry Smith ierr=VecCopy(ctx->ADADiag,v);CHKERRQ(ierr); 231e0877f53SBarry Smith if (ctx->D2){ 232e0877f53SBarry Smith ierr=VecAXPY(v, one, ctx->D2);CHKERRQ(ierr); 233e0877f53SBarry Smith } 234e0877f53SBarry Smith PetscFunctionReturn(0); 235e0877f53SBarry Smith } 236e0877f53SBarry Smith 237e0877f53SBarry Smith #undef __FUNCT__ 238e0877f53SBarry Smith #define __FUNCT__ "MatGetSubMatrix_ADA" 239e0877f53SBarry Smith static PetscErrorCode MatGetSubMatrix_ADA(Mat mat,IS isrow,IS iscol,MatReuse cll, Mat *newmat) 240e0877f53SBarry Smith { 241e0877f53SBarry Smith PetscErrorCode ierr; 242e0877f53SBarry Smith PetscInt low,high; 243e0877f53SBarry Smith IS ISrow; 244e0877f53SBarry Smith Vec D1,D2; 245e0877f53SBarry Smith Mat Atemp; 246e0877f53SBarry Smith TaoMatADACtx ctx; 247e0877f53SBarry Smith PetscBool isequal; 248e0877f53SBarry Smith 249e0877f53SBarry Smith PetscFunctionBegin; 250e0877f53SBarry Smith ierr = ISEqual(isrow,iscol,&isequal);CHKERRQ(ierr); 251e0877f53SBarry Smith if (!isequal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only for idential column and row indices"); 252e0877f53SBarry Smith ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 253e0877f53SBarry Smith 254e0877f53SBarry Smith ierr = MatGetOwnershipRange(ctx->A,&low,&high);CHKERRQ(ierr); 255*367daffbSBarry Smith ierr = ISCreateStride(PetscObjectComm((PetscObject)mat),high-low,low,1,&ISrow);CHKERRQ(ierr); 256e0877f53SBarry Smith ierr = MatGetSubMatrix(ctx->A,ISrow,iscol,cll,&Atemp);CHKERRQ(ierr); 257e0877f53SBarry Smith ierr = ISDestroy(&ISrow);CHKERRQ(ierr); 258e0877f53SBarry Smith 259e0877f53SBarry Smith if (ctx->D1){ 260e0877f53SBarry Smith ierr=VecDuplicate(ctx->D1,&D1);CHKERRQ(ierr); 261e0877f53SBarry Smith ierr=VecCopy(ctx->D1,D1);CHKERRQ(ierr); 262e0877f53SBarry Smith } else { 263e0877f53SBarry Smith D1 = NULL; 264e0877f53SBarry Smith } 265e0877f53SBarry Smith 266e0877f53SBarry Smith if (ctx->D2){ 267e0877f53SBarry Smith Vec D2sub; 268e0877f53SBarry Smith 269e0877f53SBarry Smith ierr=VecGetSubVector(ctx->D2,isrow,&D2sub);CHKERRQ(ierr); 270e0877f53SBarry Smith ierr=VecDuplicate(D2sub,&D2);CHKERRQ(ierr); 271e0877f53SBarry Smith ierr=VecCopy(D2sub,D2);CHKERRQ(ierr); 272e0877f53SBarry Smith ierr=VecRestoreSubVector(ctx->D2,isrow,&D2sub);CHKERRQ(ierr); 273e0877f53SBarry Smith } else { 274e0877f53SBarry Smith D2 = NULL; 275e0877f53SBarry Smith } 276e0877f53SBarry Smith 277e0877f53SBarry Smith ierr = MatCreateADA(Atemp,D1,D2,newmat);CHKERRQ(ierr); 278e0877f53SBarry Smith ierr = MatShellGetContext(*newmat,(void **)&ctx);CHKERRQ(ierr); 279e0877f53SBarry Smith ierr = PetscObjectDereference((PetscObject)Atemp);CHKERRQ(ierr); 280e0877f53SBarry Smith if (ctx->D1){ 281e0877f53SBarry Smith ierr = PetscObjectDereference((PetscObject)D1);CHKERRQ(ierr); 282e0877f53SBarry Smith } 283e0877f53SBarry Smith if (ctx->D2){ 284e0877f53SBarry Smith ierr = PetscObjectDereference((PetscObject)D2);CHKERRQ(ierr); 285e0877f53SBarry Smith } 286e0877f53SBarry Smith PetscFunctionReturn(0); 287e0877f53SBarry Smith } 288e0877f53SBarry Smith 289e0877f53SBarry Smith #undef __FUNCT__ 290e0877f53SBarry Smith #define __FUNCT__ "MatGetSubMatrices_ADA" 291e0877f53SBarry Smith static PetscErrorCode MatGetSubMatrices_ADA(Mat A,PetscInt n, IS *irow,IS *icol,MatReuse scall,Mat **B) 292e0877f53SBarry Smith { 293e0877f53SBarry Smith PetscErrorCode ierr; 294e0877f53SBarry Smith PetscInt i; 295e0877f53SBarry Smith 296e0877f53SBarry Smith PetscFunctionBegin; 297e0877f53SBarry Smith if (scall == MAT_INITIAL_MATRIX) { 298e0877f53SBarry Smith ierr = PetscMalloc1(n+1,B );CHKERRQ(ierr); 299e0877f53SBarry Smith } 300e0877f53SBarry Smith for ( i=0; i<n; i++ ) { 301e0877f53SBarry Smith ierr = MatGetSubMatrix_ADA(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr); 302e0877f53SBarry Smith } 303e0877f53SBarry Smith PetscFunctionReturn(0); 304e0877f53SBarry Smith } 305e0877f53SBarry Smith 306e0877f53SBarry Smith #undef __FUNCT__ 307e0877f53SBarry Smith #define __FUNCT__ "MatGetColumnVector_ADA" 308e0877f53SBarry Smith static PetscErrorCode MatGetColumnVector_ADA(Mat mat,Vec Y, PetscInt col) 309e0877f53SBarry Smith { 310e0877f53SBarry Smith PetscErrorCode ierr; 311e0877f53SBarry Smith PetscInt low,high; 312e0877f53SBarry Smith PetscScalar zero=0.0,one=1.0; 313e0877f53SBarry Smith 314e0877f53SBarry Smith PetscFunctionBegin; 315e0877f53SBarry Smith ierr=VecSet(Y, zero);CHKERRQ(ierr); 316e0877f53SBarry Smith ierr=VecGetOwnershipRange(Y,&low,&high);CHKERRQ(ierr); 317e0877f53SBarry Smith if (col>=low && col<high){ 318e0877f53SBarry Smith ierr=VecSetValue(Y,col,one,INSERT_VALUES);CHKERRQ(ierr); 319e0877f53SBarry Smith } 320e0877f53SBarry Smith ierr=VecAssemblyBegin(Y);CHKERRQ(ierr); 321e0877f53SBarry Smith ierr=VecAssemblyEnd(Y);CHKERRQ(ierr); 322e0877f53SBarry Smith ierr=MatMult_ADA(mat,Y,Y);CHKERRQ(ierr); 323e0877f53SBarry Smith PetscFunctionReturn(0); 324e0877f53SBarry Smith } 325e0877f53SBarry Smith 326e0877f53SBarry Smith #undef __FUNCT__ 327e0877f53SBarry Smith #define __FUNCT__ "MatConvert_ADA" 328cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_ADA(Mat mat,MatType newtype,Mat *NewMat) 329e0877f53SBarry Smith { 330e0877f53SBarry Smith PetscErrorCode ierr; 331e0877f53SBarry Smith PetscMPIInt size; 332e0877f53SBarry Smith PetscBool sametype, issame, isdense, isseqdense; 333e0877f53SBarry Smith TaoMatADACtx ctx; 334e0877f53SBarry Smith 335e0877f53SBarry Smith PetscFunctionBegin; 336e0877f53SBarry Smith ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 337*367daffbSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr); 338e0877f53SBarry Smith 339e0877f53SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)mat,newtype,&sametype);CHKERRQ(ierr); 340e0877f53SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)mat,MATSAME,&issame);CHKERRQ(ierr); 341e0877f53SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPIDENSE,&isdense);CHKERRQ(ierr); 342e0877f53SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQDENSE,&isseqdense);CHKERRQ(ierr); 343e0877f53SBarry Smith 344e0877f53SBarry Smith if (sametype || issame) { 345e0877f53SBarry Smith ierr=MatDuplicate(mat,MAT_COPY_VALUES,NewMat);CHKERRQ(ierr); 346e0877f53SBarry Smith } else if (isdense) { 347e0877f53SBarry Smith PetscInt i,j,low,high,m,n,M,N; 348e0877f53SBarry Smith const PetscScalar *dptr; 349e0877f53SBarry Smith Vec X; 350e0877f53SBarry Smith 351e0877f53SBarry Smith ierr = VecDuplicate(ctx->D2,&X);CHKERRQ(ierr); 352e0877f53SBarry Smith ierr=MatGetSize(mat,&M,&N);CHKERRQ(ierr); 353e0877f53SBarry Smith ierr=MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 354*367daffbSBarry Smith ierr = MatCreateDense(PetscObjectComm((PetscObject)mat),m,m,N,N,NULL,NewMat);CHKERRQ(ierr); 355e0877f53SBarry Smith ierr = MatGetOwnershipRange(*NewMat,&low,&high);CHKERRQ(ierr); 356e0877f53SBarry Smith for (i=0;i<M;i++){ 357e0877f53SBarry Smith ierr = MatGetColumnVector_ADA(mat,X,i);CHKERRQ(ierr); 358e0877f53SBarry Smith ierr = VecGetArrayRead(X,&dptr);CHKERRQ(ierr); 359e0877f53SBarry Smith for (j=0; j<high-low; j++){ 360e0877f53SBarry Smith ierr = MatSetValue(*NewMat,low+j,i,dptr[j],INSERT_VALUES);CHKERRQ(ierr); 361e0877f53SBarry Smith } 362e0877f53SBarry Smith ierr = VecRestoreArrayRead(X,&dptr);CHKERRQ(ierr); 363e0877f53SBarry Smith } 364e0877f53SBarry Smith ierr = MatAssemblyBegin(*NewMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 365e0877f53SBarry Smith ierr = MatAssemblyEnd(*NewMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 366e0877f53SBarry Smith ierr = VecDestroy(&X);CHKERRQ(ierr); 367e0877f53SBarry Smith } else if (isseqdense && size==1){ 368e0877f53SBarry Smith PetscInt i,j,low,high,m,n,M,N; 369e0877f53SBarry Smith const PetscScalar *dptr; 370e0877f53SBarry Smith Vec X; 371e0877f53SBarry Smith 372e0877f53SBarry Smith ierr = VecDuplicate(ctx->D2,&X);CHKERRQ(ierr); 373e0877f53SBarry Smith ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 374e0877f53SBarry Smith ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 375*367daffbSBarry Smith ierr = MatCreateSeqDense(PetscObjectComm((PetscObject)mat),N,N,NULL,NewMat);CHKERRQ(ierr); 376e0877f53SBarry Smith ierr = MatGetOwnershipRange(*NewMat,&low,&high);CHKERRQ(ierr); 377e0877f53SBarry Smith for (i=0;i<M;i++){ 378e0877f53SBarry Smith ierr = MatGetColumnVector_ADA(mat,X,i);CHKERRQ(ierr); 379e0877f53SBarry Smith ierr = VecGetArrayRead(X,&dptr);CHKERRQ(ierr); 380e0877f53SBarry Smith for (j=0; j<high-low; j++){ 381e0877f53SBarry Smith ierr = MatSetValue(*NewMat,low+j,i,dptr[j],INSERT_VALUES);CHKERRQ(ierr); 382e0877f53SBarry Smith } 383e0877f53SBarry Smith ierr = VecRestoreArrayRead(X,&dptr);CHKERRQ(ierr); 384e0877f53SBarry Smith } 385e0877f53SBarry Smith ierr = MatAssemblyBegin(*NewMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 386e0877f53SBarry Smith ierr = MatAssemblyEnd(*NewMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 387e0877f53SBarry Smith ierr = VecDestroy(&X);CHKERRQ(ierr); 388e0877f53SBarry Smith } else SETERRQ(PETSC_COMM_SELF,1,"No support to convert objects to that type"); 389e0877f53SBarry Smith PetscFunctionReturn(0); 390e0877f53SBarry Smith } 391e0877f53SBarry Smith 392e0877f53SBarry Smith #undef __FUNCT__ 393e0877f53SBarry Smith #define __FUNCT__ "MatNorm_ADA" 394e0877f53SBarry Smith static PetscErrorCode MatNorm_ADA(Mat mat,NormType type,PetscReal *norm) 395e0877f53SBarry Smith { 396e0877f53SBarry Smith PetscErrorCode ierr; 397e0877f53SBarry Smith TaoMatADACtx ctx; 398e0877f53SBarry Smith 399e0877f53SBarry Smith PetscFunctionBegin; 400e0877f53SBarry Smith ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 401e0877f53SBarry Smith if (type == NORM_FROBENIUS) { 402e0877f53SBarry Smith *norm = 1.0; 403e0877f53SBarry Smith } else if (type == NORM_1 || type == NORM_INFINITY) { 404e0877f53SBarry Smith *norm = 1.0; 405e0877f53SBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No two norm"); 406e0877f53SBarry Smith PetscFunctionReturn(0); 407e0877f53SBarry Smith } 408a7e14dcfSSatish Balay 409a7e14dcfSSatish Balay #undef __FUNCT__ 410a7e14dcfSSatish Balay #define __FUNCT__ "MatCreateADA" 411a7e14dcfSSatish Balay /*@C 412a7e14dcfSSatish Balay MatCreateADA - Creates a matrix M=A^T D1 A + D2 where D1, D2 are diagonal 413a7e14dcfSSatish Balay 414a7e14dcfSSatish Balay Collective on matrix 415a7e14dcfSSatish Balay 416a7e14dcfSSatish Balay Input Parameters: 417a7e14dcfSSatish Balay + mat - matrix of arbitrary type 418a7e14dcfSSatish Balay . d1 - A vector with diagonal elements of D1 419a7e14dcfSSatish Balay - d2 - A vector 420a7e14dcfSSatish Balay 421a7e14dcfSSatish Balay Output Parameters: 422a7e14dcfSSatish Balay . J - New matrix whose operations are defined in terms of mat, D1, and D2. 423a7e14dcfSSatish Balay 424a7e14dcfSSatish Balay Notes: 425a7e14dcfSSatish Balay The user provides the input data and is responsible for destroying 426a7e14dcfSSatish Balay this data after matrix J has been destroyed. 427a7e14dcfSSatish Balay The operation MatMult(A,D2,D1) must be well defined. 428a7e14dcfSSatish Balay Before calling the operation MatGetDiagonal(), the function 429a7e14dcfSSatish Balay MatADAComputeDiagonal() must be called. The matrices A and D1 must 430a7e14dcfSSatish Balay be the same during calls to MatADAComputeDiagonal() and 431a7e14dcfSSatish Balay MatGetDiagonal(). 432a7e14dcfSSatish Balay 433a7e14dcfSSatish Balay Level: developer 434a7e14dcfSSatish Balay 435a7e14dcfSSatish Balay .seealso: MatCreate() 436a7e14dcfSSatish Balay @*/ 437a7e14dcfSSatish Balay PetscErrorCode MatCreateADA(Mat mat,Vec d1, Vec d2, Mat *J) 438a7e14dcfSSatish Balay { 439*367daffbSBarry Smith MPI_Comm comm=PetscObjectComm((PetscObject)mat); 440a7e14dcfSSatish Balay TaoMatADACtx ctx; 441a7e14dcfSSatish Balay PetscErrorCode ierr; 442a7e14dcfSSatish Balay PetscInt nloc,n; 443a7e14dcfSSatish Balay 444a7e14dcfSSatish Balay PetscFunctionBegin; 4453c9e27cfSGeoffrey Irving ierr = PetscNew(&ctx);CHKERRQ(ierr); 446a7e14dcfSSatish Balay ctx->A=mat; 447a7e14dcfSSatish Balay ctx->D1=d1; 448a7e14dcfSSatish Balay ctx->D2=d2; 449a7e14dcfSSatish Balay if (d1){ 450a7e14dcfSSatish Balay ierr = VecDuplicate(d1,&ctx->W);CHKERRQ(ierr); 451a7e14dcfSSatish Balay ierr = PetscObjectReference((PetscObject)d1);CHKERRQ(ierr); 452a7e14dcfSSatish Balay } else { 4530e660641SBarry Smith ctx->W = NULL; 454a7e14dcfSSatish Balay } 455a7e14dcfSSatish Balay if (d2){ 456a7e14dcfSSatish Balay ierr = VecDuplicate(d2,&ctx->W2);CHKERRQ(ierr); 457a7e14dcfSSatish Balay ierr = VecDuplicate(d2,&ctx->ADADiag);CHKERRQ(ierr); 458a7e14dcfSSatish Balay ierr = PetscObjectReference((PetscObject)d2);CHKERRQ(ierr); 459a7e14dcfSSatish Balay } else { 4600e660641SBarry Smith ctx->W2 = NULL; 4610e660641SBarry Smith ctx->ADADiag = NULL; 462a7e14dcfSSatish Balay } 463a7e14dcfSSatish Balay 464a7e14dcfSSatish Balay ctx->GotDiag=0; 465a7e14dcfSSatish Balay ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 466a7e14dcfSSatish Balay 467a7e14dcfSSatish Balay ierr=VecGetLocalSize(d2,&nloc);CHKERRQ(ierr); 468a7e14dcfSSatish Balay ierr=VecGetSize(d2,&n);CHKERRQ(ierr); 469a7e14dcfSSatish Balay 470a7e14dcfSSatish Balay ierr = MatCreateShell(comm,nloc,nloc,n,n,ctx,J);CHKERRQ(ierr); 471a7e14dcfSSatish Balay 472a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_MULT,(void(*)(void))MatMult_ADA);CHKERRQ(ierr); 473a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void(*)(void))MatDestroy_ADA);CHKERRQ(ierr); 474a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_VIEW,(void(*)(void))MatView_ADA);CHKERRQ(ierr); 475a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_ADA);CHKERRQ(ierr); 476a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_DIAGONAL_SET,(void(*)(void))MatDiagonalSet_ADA);CHKERRQ(ierr); 477a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_SHIFT,(void(*)(void))MatShift_ADA);CHKERRQ(ierr); 478a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_EQUAL,(void(*)(void))MatEqual_ADA);CHKERRQ(ierr); 479a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_SCALE,(void(*)(void))MatScale_ADA);CHKERRQ(ierr); 480a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_TRANSPOSE,(void(*)(void))MatTranspose_ADA);CHKERRQ(ierr); 481a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_ADA);CHKERRQ(ierr); 482a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_GET_SUBMATRICES,(void(*)(void))MatGetSubMatrices_ADA);CHKERRQ(ierr); 483a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_NORM,(void(*)(void))MatNorm_ADA);CHKERRQ(ierr); 484a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_DUPLICATE,(void(*)(void))MatDuplicate_ADA);CHKERRQ(ierr); 485a7e14dcfSSatish Balay ierr = MatShellSetOperation(*J,MATOP_GET_SUBMATRIX,(void(*)(void))MatGetSubMatrix_ADA);CHKERRQ(ierr); 486a7e14dcfSSatish Balay 487a7e14dcfSSatish Balay ierr = PetscLogObjectParent((PetscObject)(*J),(PetscObject)ctx->W);CHKERRQ(ierr); 488a7e14dcfSSatish Balay ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)(*J));CHKERRQ(ierr); 489a7e14dcfSSatish Balay 490a7e14dcfSSatish Balay ierr = MatSetOption(*J,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr); 491a7e14dcfSSatish Balay PetscFunctionReturn(0); 492a7e14dcfSSatish Balay } 493