1*29bbc08cSBarry Smith /*$Id: maij.c,v 1.8 2000/07/10 03:39:48 bsmith Exp bsmith $*/ 282b1193eSBarry Smith /* 3cd3bbe55SBarry Smith Defines the basic matrix operations for the MAIJ matrix storage format. 4cd3bbe55SBarry Smith This format is used for restriction and interpolation operations for 5cd3bbe55SBarry Smith multicomponent problems. It interpolates each component the same way 6cd3bbe55SBarry Smith independently. 7cd3bbe55SBarry Smith 8cd3bbe55SBarry Smith We provide: 9cd3bbe55SBarry Smith MatMult() 10cd3bbe55SBarry Smith MatMultTranspose() 11cd3bbe55SBarry Smith MatMultTransposeAdd() 12cd3bbe55SBarry Smith MatMultAdd() 13cd3bbe55SBarry Smith and 14cd3bbe55SBarry Smith MatCreateMAIJ(Mat,dof,Mat*) 15f4a53059SBarry Smith 16f4a53059SBarry Smith This single directory handles both the sequential and parallel codes 1782b1193eSBarry Smith */ 1882b1193eSBarry Smith 19f4a53059SBarry Smith #include "src/mat/impls/aij/mpi/mpiaij.h" 2082b1193eSBarry Smith 21cd3bbe55SBarry Smith typedef struct { 22cd3bbe55SBarry Smith int dof; /* number of components */ 234cb9d9b8SBarry Smith Mat AIJ; /* representation of interpolation for one component */ 244cb9d9b8SBarry Smith } Mat_SeqMAIJ; 254cb9d9b8SBarry Smith 264cb9d9b8SBarry Smith typedef struct { 274cb9d9b8SBarry Smith int dof; /* number of components */ 28f4a53059SBarry Smith Mat AIJ,OAIJ; /* representation of interpolation for one component */ 294cb9d9b8SBarry Smith Mat A; 30f4a53059SBarry Smith VecScatter ctx; /* update ghost points for parallel case */ 31f4a53059SBarry Smith Vec w; /* work space for ghost values for parallel case */ 324cb9d9b8SBarry Smith } Mat_MPIMAIJ; 3382b1193eSBarry Smith 3482b1193eSBarry Smith #undef __FUNC__ 35b9b97703SBarry Smith #define __FUNC__ /*<a name="MatMAIJGetAIJ"></a>*/"MatMAIJGetAIJ" 36b9b97703SBarry Smith int MatMAIJGetAIJ(Mat A,Mat *B) 37b9b97703SBarry Smith { 38b9b97703SBarry Smith int ierr; 39b9b97703SBarry Smith PetscTruth ismpimaij,isseqmaij; 40b9b97703SBarry Smith 41b9b97703SBarry Smith PetscFunctionBegin; 42b9b97703SBarry Smith ierr = PetscTypeCompare((PetscObject)A,MATMPIMAIJ,&ismpimaij);CHKERRQ(ierr); 43b9b97703SBarry Smith ierr = PetscTypeCompare((PetscObject)A,MATSEQMAIJ,&isseqmaij);CHKERRQ(ierr); 44b9b97703SBarry Smith if (ismpimaij) { 45b9b97703SBarry Smith Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data; 46b9b97703SBarry Smith 47b9b97703SBarry Smith *B = b->A; 48b9b97703SBarry Smith } else if (isseqmaij) { 49b9b97703SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 50b9b97703SBarry Smith 51b9b97703SBarry Smith *B = b->AIJ; 52b9b97703SBarry Smith } else { 53b9b97703SBarry Smith *B = A; 54b9b97703SBarry Smith } 55b9b97703SBarry Smith PetscFunctionReturn(0); 56b9b97703SBarry Smith } 57b9b97703SBarry Smith 58b9b97703SBarry Smith #undef __FUNC__ 59b9b97703SBarry Smith #define __FUNC__ /*<a name="MatMAIJRedimension"></a>*/"MatMAIJRedimension" 60b9b97703SBarry Smith int MatMAIJRedimension(Mat A,int dof,Mat *B) 61b9b97703SBarry Smith { 62b9b97703SBarry Smith int ierr; 63b9b97703SBarry Smith Mat Aij; 64b9b97703SBarry Smith 65b9b97703SBarry Smith PetscFunctionBegin; 66b9b97703SBarry Smith ierr = MatMAIJGetAIJ(A,&Aij);CHKERRQ(ierr); 67b9b97703SBarry Smith ierr = MatCreateMAIJ(Aij,dof,B);CHKERRQ(ierr); 68b9b97703SBarry Smith PetscFunctionReturn(0); 69b9b97703SBarry Smith } 70b9b97703SBarry Smith 71b9b97703SBarry Smith #undef __FUNC__ 724cb9d9b8SBarry Smith #define __FUNC__ /*<a name="MatDestroy_SeqMAIJ"></a>*/"MatDestroy_SeqMAIJ" 734cb9d9b8SBarry Smith int MatDestroy_SeqMAIJ(Mat A) 7482b1193eSBarry Smith { 7582b1193eSBarry Smith int ierr; 764cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 7782b1193eSBarry Smith 7882b1193eSBarry Smith PetscFunctionBegin; 79cd3bbe55SBarry Smith if (b->AIJ) { 80cd3bbe55SBarry Smith ierr = MatDestroy(b->AIJ);CHKERRQ(ierr); 8182b1193eSBarry Smith } 824cb9d9b8SBarry Smith ierr = PetscFree(b);CHKERRQ(ierr); 834cb9d9b8SBarry Smith PetscHeaderDestroy(A); 844cb9d9b8SBarry Smith PetscFunctionReturn(0); 854cb9d9b8SBarry Smith } 864cb9d9b8SBarry Smith 874cb9d9b8SBarry Smith #undef __FUNC__ 884cb9d9b8SBarry Smith #define __FUNC__ /*<a name="MatDestroy_MPIMAIJ"></a>*/"MatDestroy_MPIMAIJ" 894cb9d9b8SBarry Smith int MatDestroy_MPIMAIJ(Mat A) 904cb9d9b8SBarry Smith { 914cb9d9b8SBarry Smith int ierr; 924cb9d9b8SBarry Smith Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data; 934cb9d9b8SBarry Smith 944cb9d9b8SBarry Smith PetscFunctionBegin; 954cb9d9b8SBarry Smith if (b->AIJ) { 964cb9d9b8SBarry Smith ierr = MatDestroy(b->AIJ);CHKERRQ(ierr); 974cb9d9b8SBarry Smith } 98f4a53059SBarry Smith if (b->OAIJ) { 99f4a53059SBarry Smith ierr = MatDestroy(b->OAIJ);CHKERRQ(ierr); 100f4a53059SBarry Smith } 101b9b97703SBarry Smith if (b->A) { 102b9b97703SBarry Smith ierr = MatDestroy(b->A);CHKERRQ(ierr); 103b9b97703SBarry Smith } 104f4a53059SBarry Smith if (b->ctx) { 105f4a53059SBarry Smith ierr = VecScatterDestroy(b->ctx);CHKERRQ(ierr); 106f4a53059SBarry Smith } 107f4a53059SBarry Smith if (b->w) { 108f4a53059SBarry Smith ierr = VecDestroy(b->w);CHKERRQ(ierr); 109f4a53059SBarry Smith } 110cd3bbe55SBarry Smith ierr = PetscFree(b);CHKERRQ(ierr); 11182b1193eSBarry Smith PetscHeaderDestroy(A); 11282b1193eSBarry Smith PetscFunctionReturn(0); 11382b1193eSBarry Smith } 11482b1193eSBarry Smith 115b9b97703SBarry Smith #undef __FUNC__ 116b9b97703SBarry Smith #define __FUNC__ /*<a name="MatGetSize_MAIJ"></a>*/"MatGetSize_MAIJ" 117b9b97703SBarry Smith static int MatGetSize_MAIJ(Mat A,int *m,int *n) 118b9b97703SBarry Smith { 119b9b97703SBarry Smith PetscFunctionBegin; 120b9b97703SBarry Smith if (m) *m = A->M; 121b9b97703SBarry Smith if (n) *n = A->N; 122b9b97703SBarry Smith PetscFunctionReturn(0); 123b9b97703SBarry Smith } 124b9b97703SBarry Smith 12582b1193eSBarry Smith EXTERN_C_BEGIN 12682b1193eSBarry Smith #undef __FUNC__ 127f4a53059SBarry Smith #define __FUNC__ /*<a name="MatCreate_MAIJ"></a>*/"MatCreate_MAIJ" 128f4a53059SBarry Smith int MatCreate_MAIJ(Mat A) 12982b1193eSBarry Smith { 130cd3bbe55SBarry Smith int ierr; 1314cb9d9b8SBarry Smith Mat_MPIMAIJ *b; 13282b1193eSBarry Smith 13382b1193eSBarry Smith PetscFunctionBegin; 1344cb9d9b8SBarry Smith A->data = (void*)(b = PetscNew(Mat_MPIMAIJ));CHKPTRQ(b); 1354cb9d9b8SBarry Smith ierr = PetscMemzero(b,sizeof(Mat_MPIMAIJ));CHKERRQ(ierr); 136cd3bbe55SBarry Smith ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 137b9b97703SBarry Smith A->ops->getsize = MatGetSize_MAIJ; 138cd3bbe55SBarry Smith A->factor = 0; 139cd3bbe55SBarry Smith A->mapping = 0; 140f4a53059SBarry Smith 141cd3bbe55SBarry Smith b->AIJ = 0; 142cd3bbe55SBarry Smith b->dof = 0; 143f4a53059SBarry Smith b->OAIJ = 0; 144f4a53059SBarry Smith b->ctx = 0; 145f4a53059SBarry Smith b->w = 0; 14682b1193eSBarry Smith PetscFunctionReturn(0); 14782b1193eSBarry Smith } 14882b1193eSBarry Smith EXTERN_C_END 14982b1193eSBarry Smith 150cd3bbe55SBarry Smith /* --------------------------------------------------------------------------------------*/ 15182b1193eSBarry Smith #undef __FUNC__ 152cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMult_SeqMAIJ_2"></a>*/"MatMult_SeqMAIJ_2" 153cd3bbe55SBarry Smith int MatMult_SeqMAIJ_2(Mat A,Vec xx,Vec yy) 15482b1193eSBarry Smith { 1554cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 156bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 157bcc973b7SBarry Smith Scalar *x,*y,*v,sum1, sum2; 158bcc973b7SBarry Smith int ierr,m = a->m,*idx,shift = a->indexshift,*ii; 159bcc973b7SBarry Smith int n,i,jrow,j; 16082b1193eSBarry Smith 161bcc973b7SBarry Smith PetscFunctionBegin; 162bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 163bcc973b7SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 164bcc973b7SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 165bcc973b7SBarry Smith idx = a->j; 166bcc973b7SBarry Smith v = a->a; 167bcc973b7SBarry Smith ii = a->i; 168bcc973b7SBarry Smith 169bcc973b7SBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 170bcc973b7SBarry Smith idx += shift; 171bcc973b7SBarry Smith for (i=0; i<m; i++) { 172bcc973b7SBarry Smith jrow = ii[i]; 173bcc973b7SBarry Smith n = ii[i+1] - jrow; 174bcc973b7SBarry Smith sum1 = 0.0; 175bcc973b7SBarry Smith sum2 = 0.0; 176bcc973b7SBarry Smith for (j=0; j<n; j++) { 177bcc973b7SBarry Smith sum1 += v[jrow]*x[2*idx[jrow]]; 178bcc973b7SBarry Smith sum2 += v[jrow]*x[2*idx[jrow]+1]; 179bcc973b7SBarry Smith jrow++; 180bcc973b7SBarry Smith } 181bcc973b7SBarry Smith y[2*i] = sum1; 182bcc973b7SBarry Smith y[2*i+1] = sum2; 183bcc973b7SBarry Smith } 184bcc973b7SBarry Smith 185bcc973b7SBarry Smith PLogFlops(4*a->nz - 2*m); 186bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 187bcc973b7SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 18882b1193eSBarry Smith PetscFunctionReturn(0); 18982b1193eSBarry Smith } 190bcc973b7SBarry Smith 19182b1193eSBarry Smith #undef __FUNC__ 192cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_SeqMAIJ_2"></a>*/"MatMultTranspose_SeqMAIJ_2" 193cd3bbe55SBarry Smith int MatMultTranspose_SeqMAIJ_2(Mat A,Vec xx,Vec yy) 19482b1193eSBarry Smith { 1954cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 196bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 197bcc973b7SBarry Smith Scalar *x,*y,*v,alpha1,alpha2,zero = 0.0; 198bcc973b7SBarry Smith int ierr,m = a->m,n,i,*idx,shift = a->indexshift; 19982b1193eSBarry Smith 200bcc973b7SBarry Smith PetscFunctionBegin; 201bcc973b7SBarry Smith ierr = VecSet(&zero,yy);CHKERRQ(ierr); 202bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 203bcc973b7SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 204bcc973b7SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 205bcc973b7SBarry Smith for (i=0; i<m; i++) { 206bcc973b7SBarry Smith idx = a->j + a->i[i] + shift; 207bcc973b7SBarry Smith v = a->a + a->i[i] + shift; 208bcc973b7SBarry Smith n = a->i[i+1] - a->i[i]; 209bcc973b7SBarry Smith alpha1 = x[2*i]; 210bcc973b7SBarry Smith alpha2 = x[2*i+1]; 211bcc973b7SBarry Smith while (n-->0) {y[2*(*idx)] += alpha1*(*v); y[2*(*idx)+1] += alpha2*(*v); idx++; v++;} 212bcc973b7SBarry Smith } 213bcc973b7SBarry Smith PLogFlops(4*a->nz - 2*a->n); 214bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 215bcc973b7SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 21682b1193eSBarry Smith PetscFunctionReturn(0); 21782b1193eSBarry Smith } 218bcc973b7SBarry Smith 21982b1193eSBarry Smith #undef __FUNC__ 220cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMultAdd_SeqMAIJ_2"></a>*/"MatMultAdd_SeqMAIJ_2" 221cd3bbe55SBarry Smith int MatMultAdd_SeqMAIJ_2(Mat A,Vec xx,Vec yy,Vec zz) 22282b1193eSBarry Smith { 2234cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 224bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 225bcc973b7SBarry Smith Scalar *x,*y,*v,sum1, sum2; 226bcc973b7SBarry Smith int ierr,m = a->m,*idx,shift = a->indexshift,*ii; 227bcc973b7SBarry Smith int n,i,jrow,j; 22882b1193eSBarry Smith 229bcc973b7SBarry Smith PetscFunctionBegin; 230f9fae5adSBarry Smith if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);} 231bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 232f9fae5adSBarry Smith ierr = VecGetArray(zz,&y);CHKERRQ(ierr); 233bcc973b7SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 234bcc973b7SBarry Smith idx = a->j; 235bcc973b7SBarry Smith v = a->a; 236bcc973b7SBarry Smith ii = a->i; 237bcc973b7SBarry Smith 238bcc973b7SBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 239bcc973b7SBarry Smith idx += shift; 240bcc973b7SBarry Smith for (i=0; i<m; i++) { 241bcc973b7SBarry Smith jrow = ii[i]; 242bcc973b7SBarry Smith n = ii[i+1] - jrow; 243bcc973b7SBarry Smith sum1 = 0.0; 244bcc973b7SBarry Smith sum2 = 0.0; 245bcc973b7SBarry Smith for (j=0; j<n; j++) { 246bcc973b7SBarry Smith sum1 += v[jrow]*x[2*idx[jrow]]; 247bcc973b7SBarry Smith sum2 += v[jrow]*x[2*idx[jrow]+1]; 248bcc973b7SBarry Smith jrow++; 249bcc973b7SBarry Smith } 250bcc973b7SBarry Smith y[2*i] += sum1; 251bcc973b7SBarry Smith y[2*i+1] += sum2; 252bcc973b7SBarry Smith } 253bcc973b7SBarry Smith 254bcc973b7SBarry Smith PLogFlops(4*a->nz - 2*m); 255bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 256f9fae5adSBarry Smith ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr); 257bcc973b7SBarry Smith PetscFunctionReturn(0); 25882b1193eSBarry Smith } 25982b1193eSBarry Smith #undef __FUNC__ 260cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_SeqMAIJ_2"></a>*/"MatMultTransposeAdd_SeqMAIJ_2" 261cd3bbe55SBarry Smith int MatMultTransposeAdd_SeqMAIJ_2(Mat A,Vec xx,Vec yy,Vec zz) 26282b1193eSBarry Smith { 2634cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 264bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 265bcc973b7SBarry Smith Scalar *x,*y,*v,alpha1,alpha2; 266bcc973b7SBarry Smith int ierr,m = a->m,n,i,*idx,shift = a->indexshift; 26782b1193eSBarry Smith 268bcc973b7SBarry Smith PetscFunctionBegin; 269f9fae5adSBarry Smith if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);} 270bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 271f9fae5adSBarry Smith ierr = VecGetArray(zz,&y);CHKERRQ(ierr); 272bcc973b7SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 273bcc973b7SBarry Smith for (i=0; i<m; i++) { 274bcc973b7SBarry Smith idx = a->j + a->i[i] + shift; 275bcc973b7SBarry Smith v = a->a + a->i[i] + shift; 276bcc973b7SBarry Smith n = a->i[i+1] - a->i[i]; 277bcc973b7SBarry Smith alpha1 = x[2*i]; 278bcc973b7SBarry Smith alpha2 = x[2*i+1]; 279bcc973b7SBarry Smith while (n-->0) {y[2*(*idx)] += alpha1*(*v); y[2*(*idx)+1] += alpha2*(*v); idx++; v++;} 280bcc973b7SBarry Smith } 281bcc973b7SBarry Smith PLogFlops(4*a->nz - 2*a->n); 282bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 283f9fae5adSBarry Smith ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr); 284bcc973b7SBarry Smith PetscFunctionReturn(0); 28582b1193eSBarry Smith } 286cd3bbe55SBarry Smith /* --------------------------------------------------------------------------------------*/ 287bcc973b7SBarry Smith #undef __FUNC__ 288bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMult_SeqMAIJ_3"></a>*/"MatMult_SeqMAIJ_3" 289bcc973b7SBarry Smith int MatMult_SeqMAIJ_3(Mat A,Vec xx,Vec yy) 290bcc973b7SBarry Smith { 2914cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 292bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 293bcc973b7SBarry Smith Scalar *x,*y,*v,sum1, sum2, sum3; 294bcc973b7SBarry Smith int ierr,m = a->m,*idx,shift = a->indexshift,*ii; 295bcc973b7SBarry Smith int n,i,jrow,j; 29682b1193eSBarry Smith 297bcc973b7SBarry Smith PetscFunctionBegin; 298bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 299bcc973b7SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 300bcc973b7SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 301bcc973b7SBarry Smith idx = a->j; 302bcc973b7SBarry Smith v = a->a; 303bcc973b7SBarry Smith ii = a->i; 304bcc973b7SBarry Smith 305bcc973b7SBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 306bcc973b7SBarry Smith idx += shift; 307bcc973b7SBarry Smith for (i=0; i<m; i++) { 308bcc973b7SBarry Smith jrow = ii[i]; 309bcc973b7SBarry Smith n = ii[i+1] - jrow; 310bcc973b7SBarry Smith sum1 = 0.0; 311bcc973b7SBarry Smith sum2 = 0.0; 312bcc973b7SBarry Smith sum3 = 0.0; 313bcc973b7SBarry Smith for (j=0; j<n; j++) { 314bcc973b7SBarry Smith sum1 += v[jrow]*x[3*idx[jrow]]; 315bcc973b7SBarry Smith sum2 += v[jrow]*x[3*idx[jrow]+1]; 316bcc973b7SBarry Smith sum3 += v[jrow]*x[3*idx[jrow]+2]; 317bcc973b7SBarry Smith jrow++; 318bcc973b7SBarry Smith } 319bcc973b7SBarry Smith y[3*i] = sum1; 320bcc973b7SBarry Smith y[3*i+1] = sum2; 321bcc973b7SBarry Smith y[3*i+2] = sum3; 322bcc973b7SBarry Smith } 323bcc973b7SBarry Smith 324bcc973b7SBarry Smith PLogFlops(6*a->nz - 3*m); 325bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 326bcc973b7SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 327bcc973b7SBarry Smith PetscFunctionReturn(0); 328bcc973b7SBarry Smith } 329bcc973b7SBarry Smith 330bcc973b7SBarry Smith #undef __FUNC__ 331bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_SeqMAIJ_3"></a>*/"MatMultTranspose_SeqMAIJ_3" 332bcc973b7SBarry Smith int MatMultTranspose_SeqMAIJ_3(Mat A,Vec xx,Vec yy) 333bcc973b7SBarry Smith { 3344cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 335bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 336bcc973b7SBarry Smith Scalar *x,*y,*v,alpha1,alpha2,alpha3,zero = 0.0; 337bcc973b7SBarry Smith int ierr,m = a->m,n,i,*idx,shift = a->indexshift; 338bcc973b7SBarry Smith 339bcc973b7SBarry Smith PetscFunctionBegin; 340bcc973b7SBarry Smith ierr = VecSet(&zero,yy);CHKERRQ(ierr); 341bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 342bcc973b7SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 343bcc973b7SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 344bcc973b7SBarry Smith for (i=0; i<m; i++) { 345bcc973b7SBarry Smith idx = a->j + a->i[i] + shift; 346bcc973b7SBarry Smith v = a->a + a->i[i] + shift; 347bcc973b7SBarry Smith n = a->i[i+1] - a->i[i]; 348bcc973b7SBarry Smith alpha1 = x[3*i]; 349bcc973b7SBarry Smith alpha2 = x[3*i+1]; 350bcc973b7SBarry Smith alpha3 = x[3*i+2]; 351bcc973b7SBarry Smith while (n-->0) { 352bcc973b7SBarry Smith y[3*(*idx)] += alpha1*(*v); 353bcc973b7SBarry Smith y[3*(*idx)+1] += alpha2*(*v); 354bcc973b7SBarry Smith y[3*(*idx)+2] += alpha3*(*v); 355bcc973b7SBarry Smith idx++; v++; 356bcc973b7SBarry Smith } 357bcc973b7SBarry Smith } 358bcc973b7SBarry Smith PLogFlops(6*a->nz - 3*a->n); 359bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 360bcc973b7SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 361bcc973b7SBarry Smith PetscFunctionReturn(0); 362bcc973b7SBarry Smith } 363bcc973b7SBarry Smith 364bcc973b7SBarry Smith #undef __FUNC__ 365bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultAdd_SeqMAIJ_3"></a>*/"MatMultAdd_SeqMAIJ_3" 366bcc973b7SBarry Smith int MatMultAdd_SeqMAIJ_3(Mat A,Vec xx,Vec yy,Vec zz) 367bcc973b7SBarry Smith { 3684cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 369bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 370bcc973b7SBarry Smith Scalar *x,*y,*v,sum1, sum2, sum3; 371bcc973b7SBarry Smith int ierr,m = a->m,*idx,shift = a->indexshift,*ii; 372bcc973b7SBarry Smith int n,i,jrow,j; 373bcc973b7SBarry Smith 374bcc973b7SBarry Smith PetscFunctionBegin; 375f9fae5adSBarry Smith if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);} 376bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 377f9fae5adSBarry Smith ierr = VecGetArray(zz,&y);CHKERRQ(ierr); 378bcc973b7SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 379bcc973b7SBarry Smith idx = a->j; 380bcc973b7SBarry Smith v = a->a; 381bcc973b7SBarry Smith ii = a->i; 382bcc973b7SBarry Smith 383bcc973b7SBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 384bcc973b7SBarry Smith idx += shift; 385bcc973b7SBarry Smith for (i=0; i<m; i++) { 386bcc973b7SBarry Smith jrow = ii[i]; 387bcc973b7SBarry Smith n = ii[i+1] - jrow; 388bcc973b7SBarry Smith sum1 = 0.0; 389bcc973b7SBarry Smith sum2 = 0.0; 390bcc973b7SBarry Smith sum3 = 0.0; 391bcc973b7SBarry Smith for (j=0; j<n; j++) { 392bcc973b7SBarry Smith sum1 += v[jrow]*x[3*idx[jrow]]; 393bcc973b7SBarry Smith sum2 += v[jrow]*x[3*idx[jrow]+1]; 394bcc973b7SBarry Smith sum3 += v[jrow]*x[3*idx[jrow]+2]; 395bcc973b7SBarry Smith jrow++; 396bcc973b7SBarry Smith } 397bcc973b7SBarry Smith y[3*i] += sum1; 398bcc973b7SBarry Smith y[3*i+1] += sum2; 399bcc973b7SBarry Smith y[3*i+2] += sum3; 400bcc973b7SBarry Smith } 401bcc973b7SBarry Smith 402bcc973b7SBarry Smith PLogFlops(6*a->nz); 403bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 404f9fae5adSBarry Smith ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr); 405bcc973b7SBarry Smith PetscFunctionReturn(0); 406bcc973b7SBarry Smith } 407bcc973b7SBarry Smith #undef __FUNC__ 408bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_SeqMAIJ_3"></a>*/"MatMultTransposeAdd_SeqMAIJ_3" 409bcc973b7SBarry Smith int MatMultTransposeAdd_SeqMAIJ_3(Mat A,Vec xx,Vec yy,Vec zz) 410bcc973b7SBarry Smith { 4114cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 412bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 413bcc973b7SBarry Smith Scalar *x,*y,*v,alpha1,alpha2,alpha3; 414bcc973b7SBarry Smith int ierr,m = a->m,n,i,*idx,shift = a->indexshift; 415bcc973b7SBarry Smith 416bcc973b7SBarry Smith PetscFunctionBegin; 417f9fae5adSBarry Smith if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);} 418bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 419f9fae5adSBarry Smith ierr = VecGetArray(zz,&y);CHKERRQ(ierr); 420bcc973b7SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 421bcc973b7SBarry Smith for (i=0; i<m; i++) { 422bcc973b7SBarry Smith idx = a->j + a->i[i] + shift; 423bcc973b7SBarry Smith v = a->a + a->i[i] + shift; 424bcc973b7SBarry Smith n = a->i[i+1] - a->i[i]; 425bcc973b7SBarry Smith alpha1 = x[3*i]; 426bcc973b7SBarry Smith alpha2 = x[3*i+1]; 427bcc973b7SBarry Smith alpha3 = x[3*i+2]; 428bcc973b7SBarry Smith while (n-->0) { 429bcc973b7SBarry Smith y[3*(*idx)] += alpha1*(*v); 430bcc973b7SBarry Smith y[3*(*idx)+1] += alpha2*(*v); 431bcc973b7SBarry Smith y[3*(*idx)+2] += alpha3*(*v); 432bcc973b7SBarry Smith idx++; v++; 433bcc973b7SBarry Smith } 434bcc973b7SBarry Smith } 435bcc973b7SBarry Smith PLogFlops(6*a->nz); 436bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 437f9fae5adSBarry Smith ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr); 438bcc973b7SBarry Smith PetscFunctionReturn(0); 439bcc973b7SBarry Smith } 440bcc973b7SBarry Smith 441bcc973b7SBarry Smith /* ------------------------------------------------------------------------------*/ 442bcc973b7SBarry Smith #undef __FUNC__ 443bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMult_SeqMAIJ_4"></a>*/"MatMult_SeqMAIJ_4" 444bcc973b7SBarry Smith int MatMult_SeqMAIJ_4(Mat A,Vec xx,Vec yy) 445bcc973b7SBarry Smith { 4464cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 447bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 448bcc973b7SBarry Smith Scalar *x,*y,*v,sum1, sum2, sum3, sum4; 449bcc973b7SBarry Smith int ierr,m = a->m,*idx,shift = a->indexshift,*ii; 450bcc973b7SBarry Smith int n,i,jrow,j; 451bcc973b7SBarry Smith 452bcc973b7SBarry Smith PetscFunctionBegin; 453bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 454bcc973b7SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 455bcc973b7SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 456bcc973b7SBarry Smith idx = a->j; 457bcc973b7SBarry Smith v = a->a; 458bcc973b7SBarry Smith ii = a->i; 459bcc973b7SBarry Smith 460bcc973b7SBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 461bcc973b7SBarry Smith idx += shift; 462bcc973b7SBarry Smith for (i=0; i<m; i++) { 463bcc973b7SBarry Smith jrow = ii[i]; 464bcc973b7SBarry Smith n = ii[i+1] - jrow; 465bcc973b7SBarry Smith sum1 = 0.0; 466bcc973b7SBarry Smith sum2 = 0.0; 467bcc973b7SBarry Smith sum3 = 0.0; 468bcc973b7SBarry Smith sum4 = 0.0; 469bcc973b7SBarry Smith for (j=0; j<n; j++) { 470bcc973b7SBarry Smith sum1 += v[jrow]*x[4*idx[jrow]]; 471bcc973b7SBarry Smith sum2 += v[jrow]*x[4*idx[jrow]+1]; 472bcc973b7SBarry Smith sum3 += v[jrow]*x[4*idx[jrow]+2]; 473bcc973b7SBarry Smith sum4 += v[jrow]*x[4*idx[jrow]+3]; 474bcc973b7SBarry Smith jrow++; 475bcc973b7SBarry Smith } 476bcc973b7SBarry Smith y[4*i] = sum1; 477bcc973b7SBarry Smith y[4*i+1] = sum2; 478bcc973b7SBarry Smith y[4*i+2] = sum3; 479bcc973b7SBarry Smith y[4*i+3] = sum4; 480bcc973b7SBarry Smith } 481bcc973b7SBarry Smith 482bcc973b7SBarry Smith PLogFlops(8*a->nz - 4*m); 483bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 484bcc973b7SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 485bcc973b7SBarry Smith PetscFunctionReturn(0); 486bcc973b7SBarry Smith } 487bcc973b7SBarry Smith 488bcc973b7SBarry Smith #undef __FUNC__ 489bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_SeqMAIJ_4"></a>*/"MatMultTranspose_SeqMAIJ_4" 490bcc973b7SBarry Smith int MatMultTranspose_SeqMAIJ_4(Mat A,Vec xx,Vec yy) 491bcc973b7SBarry Smith { 4924cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 493bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 494bcc973b7SBarry Smith Scalar *x,*y,*v,alpha1,alpha2,alpha3,alpha4,zero = 0.0; 495bcc973b7SBarry Smith int ierr,m = a->m,n,i,*idx,shift = a->indexshift; 496bcc973b7SBarry Smith 497bcc973b7SBarry Smith PetscFunctionBegin; 498bcc973b7SBarry Smith ierr = VecSet(&zero,yy);CHKERRQ(ierr); 499bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 500bcc973b7SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 501bcc973b7SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 502bcc973b7SBarry Smith for (i=0; i<m; i++) { 503bcc973b7SBarry Smith idx = a->j + a->i[i] + shift; 504bcc973b7SBarry Smith v = a->a + a->i[i] + shift; 505bcc973b7SBarry Smith n = a->i[i+1] - a->i[i]; 506bcc973b7SBarry Smith alpha1 = x[4*i]; 507bcc973b7SBarry Smith alpha2 = x[4*i+1]; 508bcc973b7SBarry Smith alpha3 = x[4*i+2]; 509bcc973b7SBarry Smith alpha4 = x[4*i+3]; 510bcc973b7SBarry Smith while (n-->0) { 511bcc973b7SBarry Smith y[4*(*idx)] += alpha1*(*v); 512bcc973b7SBarry Smith y[4*(*idx)+1] += alpha2*(*v); 513bcc973b7SBarry Smith y[4*(*idx)+2] += alpha3*(*v); 514bcc973b7SBarry Smith y[4*(*idx)+3] += alpha4*(*v); 515bcc973b7SBarry Smith idx++; v++; 516bcc973b7SBarry Smith } 517bcc973b7SBarry Smith } 518bcc973b7SBarry Smith PLogFlops(8*a->nz - 4*a->n); 519bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 520bcc973b7SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 521bcc973b7SBarry Smith PetscFunctionReturn(0); 522bcc973b7SBarry Smith } 523bcc973b7SBarry Smith 524bcc973b7SBarry Smith #undef __FUNC__ 525bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultAdd_SeqMAIJ_4"></a>*/"MatMultAdd_SeqMAIJ_4" 526bcc973b7SBarry Smith int MatMultAdd_SeqMAIJ_4(Mat A,Vec xx,Vec yy,Vec zz) 527bcc973b7SBarry Smith { 5284cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 529f9fae5adSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 530f9fae5adSBarry Smith Scalar *x,*y,*v,sum1, sum2, sum3, sum4; 531f9fae5adSBarry Smith int ierr,m = a->m,*idx,shift = a->indexshift,*ii; 532f9fae5adSBarry Smith int n,i,jrow,j; 533f9fae5adSBarry Smith 534f9fae5adSBarry Smith PetscFunctionBegin; 535f9fae5adSBarry Smith if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);} 536f9fae5adSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 537f9fae5adSBarry Smith ierr = VecGetArray(zz,&y);CHKERRQ(ierr); 538f9fae5adSBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 539f9fae5adSBarry Smith idx = a->j; 540f9fae5adSBarry Smith v = a->a; 541f9fae5adSBarry Smith ii = a->i; 542f9fae5adSBarry Smith 543f9fae5adSBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 544f9fae5adSBarry Smith idx += shift; 545f9fae5adSBarry Smith for (i=0; i<m; i++) { 546f9fae5adSBarry Smith jrow = ii[i]; 547f9fae5adSBarry Smith n = ii[i+1] - jrow; 548f9fae5adSBarry Smith sum1 = 0.0; 549f9fae5adSBarry Smith sum2 = 0.0; 550f9fae5adSBarry Smith sum3 = 0.0; 551f9fae5adSBarry Smith sum4 = 0.0; 552f9fae5adSBarry Smith for (j=0; j<n; j++) { 553f9fae5adSBarry Smith sum1 += v[jrow]*x[4*idx[jrow]]; 554f9fae5adSBarry Smith sum2 += v[jrow]*x[4*idx[jrow]+1]; 555f9fae5adSBarry Smith sum3 += v[jrow]*x[4*idx[jrow]+2]; 556f9fae5adSBarry Smith sum4 += v[jrow]*x[4*idx[jrow]+3]; 557f9fae5adSBarry Smith jrow++; 558f9fae5adSBarry Smith } 559f9fae5adSBarry Smith y[4*i] += sum1; 560f9fae5adSBarry Smith y[4*i+1] += sum2; 561f9fae5adSBarry Smith y[4*i+2] += sum3; 562f9fae5adSBarry Smith y[4*i+3] += sum4; 563f9fae5adSBarry Smith } 564f9fae5adSBarry Smith 565f9fae5adSBarry Smith PLogFlops(8*a->nz - 4*m); 566f9fae5adSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 567f9fae5adSBarry Smith ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr); 568f9fae5adSBarry Smith PetscFunctionReturn(0); 569bcc973b7SBarry Smith } 570bcc973b7SBarry Smith #undef __FUNC__ 571bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_SeqMAIJ_4"></a>*/"MatMultTransposeAdd_SeqMAIJ_4" 572bcc973b7SBarry Smith int MatMultTransposeAdd_SeqMAIJ_4(Mat A,Vec xx,Vec yy,Vec zz) 573bcc973b7SBarry Smith { 5744cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 575f9fae5adSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 576f9fae5adSBarry Smith Scalar *x,*y,*v,alpha1,alpha2,alpha3,alpha4; 577f9fae5adSBarry Smith int ierr,m = a->m,n,i,*idx,shift = a->indexshift; 578f9fae5adSBarry Smith 579f9fae5adSBarry Smith PetscFunctionBegin; 580f9fae5adSBarry Smith if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);} 581f9fae5adSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 582f9fae5adSBarry Smith ierr = VecGetArray(zz,&y);CHKERRQ(ierr); 583f9fae5adSBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 584f9fae5adSBarry Smith for (i=0; i<m; i++) { 585f9fae5adSBarry Smith idx = a->j + a->i[i] + shift; 586f9fae5adSBarry Smith v = a->a + a->i[i] + shift; 587f9fae5adSBarry Smith n = a->i[i+1] - a->i[i]; 588f9fae5adSBarry Smith alpha1 = x[4*i]; 589f9fae5adSBarry Smith alpha2 = x[4*i+1]; 590f9fae5adSBarry Smith alpha3 = x[4*i+2]; 591f9fae5adSBarry Smith alpha4 = x[4*i+3]; 592f9fae5adSBarry Smith while (n-->0) { 593f9fae5adSBarry Smith y[4*(*idx)] += alpha1*(*v); 594f9fae5adSBarry Smith y[4*(*idx)+1] += alpha2*(*v); 595f9fae5adSBarry Smith y[4*(*idx)+2] += alpha3*(*v); 596f9fae5adSBarry Smith y[4*(*idx)+3] += alpha4*(*v); 597f9fae5adSBarry Smith idx++; v++; 598f9fae5adSBarry Smith } 599f9fae5adSBarry Smith } 600f9fae5adSBarry Smith PLogFlops(8*a->nz - 4*a->n); 601f9fae5adSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 602f9fae5adSBarry Smith ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr); 603f9fae5adSBarry Smith PetscFunctionReturn(0); 604f9fae5adSBarry Smith 605f9fae5adSBarry Smith } 606f9fae5adSBarry Smith 607f9fae5adSBarry Smith /* ------------------------------------------------------------------------------*/ 608f9fae5adSBarry Smith #undef __FUNC__ 609f9fae5adSBarry Smith #define __FUNC__ /*<a name="MatMult_SeqMAIJ_5"></a>*/"MatMult_SeqMAIJ_5" 610f9fae5adSBarry Smith int MatMult_SeqMAIJ_5(Mat A,Vec xx,Vec yy) 611f9fae5adSBarry Smith { 6124cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 613f9fae5adSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 614f9fae5adSBarry Smith Scalar *x,*y,*v,sum1, sum2, sum3, sum4, sum5; 615f9fae5adSBarry Smith int ierr,m = a->m,*idx,shift = a->indexshift,*ii; 616f9fae5adSBarry Smith int n,i,jrow,j; 617f9fae5adSBarry Smith 618f9fae5adSBarry Smith PetscFunctionBegin; 619f9fae5adSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 620f9fae5adSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 621f9fae5adSBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 622f9fae5adSBarry Smith idx = a->j; 623f9fae5adSBarry Smith v = a->a; 624f9fae5adSBarry Smith ii = a->i; 625f9fae5adSBarry Smith 626f9fae5adSBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 627f9fae5adSBarry Smith idx += shift; 628f9fae5adSBarry Smith for (i=0; i<m; i++) { 629f9fae5adSBarry Smith jrow = ii[i]; 630f9fae5adSBarry Smith n = ii[i+1] - jrow; 631f9fae5adSBarry Smith sum1 = 0.0; 632f9fae5adSBarry Smith sum2 = 0.0; 633f9fae5adSBarry Smith sum3 = 0.0; 634f9fae5adSBarry Smith sum4 = 0.0; 635f9fae5adSBarry Smith sum5 = 0.0; 636f9fae5adSBarry Smith for (j=0; j<n; j++) { 637f9fae5adSBarry Smith sum1 += v[jrow]*x[5*idx[jrow]]; 638f9fae5adSBarry Smith sum2 += v[jrow]*x[5*idx[jrow]+1]; 639f9fae5adSBarry Smith sum3 += v[jrow]*x[5*idx[jrow]+2]; 640f9fae5adSBarry Smith sum4 += v[jrow]*x[5*idx[jrow]+3]; 641f9fae5adSBarry Smith sum5 += v[jrow]*x[5*idx[jrow]+4]; 642f9fae5adSBarry Smith jrow++; 643f9fae5adSBarry Smith } 644f9fae5adSBarry Smith y[5*i] = sum1; 645f9fae5adSBarry Smith y[5*i+1] = sum2; 646f9fae5adSBarry Smith y[5*i+2] = sum3; 647f9fae5adSBarry Smith y[5*i+3] = sum4; 648f9fae5adSBarry Smith y[5*i+4] = sum5; 649f9fae5adSBarry Smith } 650f9fae5adSBarry Smith 651f9fae5adSBarry Smith PLogFlops(10*a->nz - 5*m); 652f9fae5adSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 653f9fae5adSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 654f9fae5adSBarry Smith PetscFunctionReturn(0); 655f9fae5adSBarry Smith } 656f9fae5adSBarry Smith 657f9fae5adSBarry Smith #undef __FUNC__ 658f9fae5adSBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_SeqMAIJ_5"></a>*/"MatMultTranspose_SeqMAIJ_5" 659f9fae5adSBarry Smith int MatMultTranspose_SeqMAIJ_5(Mat A,Vec xx,Vec yy) 660f9fae5adSBarry Smith { 6614cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 662f9fae5adSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 663f9fae5adSBarry Smith Scalar *x,*y,*v,alpha1,alpha2,alpha3,alpha4,alpha5,zero = 0.0; 664f9fae5adSBarry Smith int ierr,m = a->m,n,i,*idx,shift = a->indexshift; 665f9fae5adSBarry Smith 666f9fae5adSBarry Smith PetscFunctionBegin; 667f9fae5adSBarry Smith ierr = VecSet(&zero,yy);CHKERRQ(ierr); 668f9fae5adSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 669f9fae5adSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 670f9fae5adSBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 671f9fae5adSBarry Smith for (i=0; i<m; i++) { 672f9fae5adSBarry Smith idx = a->j + a->i[i] + shift; 673f9fae5adSBarry Smith v = a->a + a->i[i] + shift; 674f9fae5adSBarry Smith n = a->i[i+1] - a->i[i]; 675f9fae5adSBarry Smith alpha1 = x[5*i]; 676f9fae5adSBarry Smith alpha2 = x[5*i+1]; 677f9fae5adSBarry Smith alpha3 = x[5*i+2]; 678f9fae5adSBarry Smith alpha4 = x[5*i+3]; 679f9fae5adSBarry Smith alpha5 = x[5*i+4]; 680f9fae5adSBarry Smith while (n-->0) { 681f9fae5adSBarry Smith y[5*(*idx)] += alpha1*(*v); 682f9fae5adSBarry Smith y[5*(*idx)+1] += alpha2*(*v); 683f9fae5adSBarry Smith y[5*(*idx)+2] += alpha3*(*v); 684f9fae5adSBarry Smith y[5*(*idx)+3] += alpha4*(*v); 685f9fae5adSBarry Smith y[5*(*idx)+4] += alpha5*(*v); 686f9fae5adSBarry Smith idx++; v++; 687f9fae5adSBarry Smith } 688f9fae5adSBarry Smith } 689f9fae5adSBarry Smith PLogFlops(10*a->nz - 5*a->n); 690f9fae5adSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 691f9fae5adSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 692f9fae5adSBarry Smith PetscFunctionReturn(0); 693f9fae5adSBarry Smith } 694f9fae5adSBarry Smith 695f9fae5adSBarry Smith #undef __FUNC__ 696f9fae5adSBarry Smith #define __FUNC__ /*<a name="MatMultAdd_SeqMAIJ_5"></a>*/"MatMultAdd_SeqMAIJ_5" 697f9fae5adSBarry Smith int MatMultAdd_SeqMAIJ_5(Mat A,Vec xx,Vec yy,Vec zz) 698f9fae5adSBarry Smith { 6994cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 700f9fae5adSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 701f9fae5adSBarry Smith Scalar *x,*y,*v,sum1, sum2, sum3, sum4, sum5; 702f9fae5adSBarry Smith int ierr,m = a->m,*idx,shift = a->indexshift,*ii; 703f9fae5adSBarry Smith int n,i,jrow,j; 704f9fae5adSBarry Smith 705f9fae5adSBarry Smith PetscFunctionBegin; 706f9fae5adSBarry Smith if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);} 707f9fae5adSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 708f9fae5adSBarry Smith ierr = VecGetArray(zz,&y);CHKERRQ(ierr); 709f9fae5adSBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 710f9fae5adSBarry Smith idx = a->j; 711f9fae5adSBarry Smith v = a->a; 712f9fae5adSBarry Smith ii = a->i; 713f9fae5adSBarry Smith 714f9fae5adSBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 715f9fae5adSBarry Smith idx += shift; 716f9fae5adSBarry Smith for (i=0; i<m; i++) { 717f9fae5adSBarry Smith jrow = ii[i]; 718f9fae5adSBarry Smith n = ii[i+1] - jrow; 719f9fae5adSBarry Smith sum1 = 0.0; 720f9fae5adSBarry Smith sum2 = 0.0; 721f9fae5adSBarry Smith sum3 = 0.0; 722f9fae5adSBarry Smith sum4 = 0.0; 723f9fae5adSBarry Smith sum5 = 0.0; 724f9fae5adSBarry Smith for (j=0; j<n; j++) { 725f9fae5adSBarry Smith sum1 += v[jrow]*x[5*idx[jrow]]; 726f9fae5adSBarry Smith sum2 += v[jrow]*x[5*idx[jrow]+1]; 727f9fae5adSBarry Smith sum3 += v[jrow]*x[5*idx[jrow]+2]; 728f9fae5adSBarry Smith sum4 += v[jrow]*x[5*idx[jrow]+3]; 729f9fae5adSBarry Smith sum5 += v[jrow]*x[5*idx[jrow]+4]; 730f9fae5adSBarry Smith jrow++; 731f9fae5adSBarry Smith } 732f9fae5adSBarry Smith y[5*i] += sum1; 733f9fae5adSBarry Smith y[5*i+1] += sum2; 734f9fae5adSBarry Smith y[5*i+2] += sum3; 735f9fae5adSBarry Smith y[5*i+3] += sum4; 736f9fae5adSBarry Smith y[5*i+4] += sum5; 737f9fae5adSBarry Smith } 738f9fae5adSBarry Smith 739f9fae5adSBarry Smith PLogFlops(10*a->nz); 740f9fae5adSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 741f9fae5adSBarry Smith ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr); 742f9fae5adSBarry Smith PetscFunctionReturn(0); 743f9fae5adSBarry Smith } 744f9fae5adSBarry Smith 745f9fae5adSBarry Smith #undef __FUNC__ 746f9fae5adSBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_SeqMAIJ_5"></a>*/"MatMultTransposeAdd_SeqMAIJ_5" 747f9fae5adSBarry Smith int MatMultTransposeAdd_SeqMAIJ_5(Mat A,Vec xx,Vec yy,Vec zz) 748f9fae5adSBarry Smith { 7494cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 750f9fae5adSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 751f9fae5adSBarry Smith Scalar *x,*y,*v,alpha1,alpha2,alpha3,alpha4,alpha5; 752f9fae5adSBarry Smith int ierr,m = a->m,n,i,*idx,shift = a->indexshift; 753f9fae5adSBarry Smith 754f9fae5adSBarry Smith PetscFunctionBegin; 755f9fae5adSBarry Smith if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);} 756f9fae5adSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 757f9fae5adSBarry Smith ierr = VecGetArray(zz,&y);CHKERRQ(ierr); 758f9fae5adSBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 759f9fae5adSBarry Smith for (i=0; i<m; i++) { 760f9fae5adSBarry Smith idx = a->j + a->i[i] + shift; 761f9fae5adSBarry Smith v = a->a + a->i[i] + shift; 762f9fae5adSBarry Smith n = a->i[i+1] - a->i[i]; 763f9fae5adSBarry Smith alpha1 = x[5*i]; 764f9fae5adSBarry Smith alpha2 = x[5*i+1]; 765f9fae5adSBarry Smith alpha3 = x[5*i+2]; 766f9fae5adSBarry Smith alpha4 = x[5*i+3]; 767f9fae5adSBarry Smith alpha5 = x[5*i+4]; 768f9fae5adSBarry Smith while (n-->0) { 769f9fae5adSBarry Smith y[5*(*idx)] += alpha1*(*v); 770f9fae5adSBarry Smith y[5*(*idx)+1] += alpha2*(*v); 771f9fae5adSBarry Smith y[5*(*idx)+2] += alpha3*(*v); 772f9fae5adSBarry Smith y[5*(*idx)+3] += alpha4*(*v); 773f9fae5adSBarry Smith y[5*(*idx)+4] += alpha5*(*v); 774f9fae5adSBarry Smith idx++; v++; 775f9fae5adSBarry Smith } 776f9fae5adSBarry Smith } 777f9fae5adSBarry Smith PLogFlops(10*a->nz); 778f9fae5adSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 779f9fae5adSBarry Smith ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr); 780f9fae5adSBarry Smith PetscFunctionReturn(0); 781bcc973b7SBarry Smith } 782bcc973b7SBarry Smith 783f4a53059SBarry Smith /*===================================================================================*/ 784f4a53059SBarry Smith #undef __FUNC__ 785f4a53059SBarry Smith #define __FUNC__ /*<a name="MatMult_MPIMAIJ_dof"></a>*/"MatMult_MPIMAIJ_dof" 786f4a53059SBarry Smith int MatMult_MPIMAIJ_dof(Mat A,Vec xx,Vec yy) 787f4a53059SBarry Smith { 7884cb9d9b8SBarry Smith Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data; 789f4a53059SBarry Smith int ierr; 790f4a53059SBarry Smith PetscFunctionBegin; 791f4a53059SBarry Smith 792f4a53059SBarry Smith /* start the scatter */ 793f4a53059SBarry Smith ierr = VecScatterBegin(xx,b->w,INSERT_VALUES,SCATTER_FORWARD,b->ctx);CHKERRQ(ierr); 7944cb9d9b8SBarry Smith ierr = (*b->AIJ->ops->mult)(b->AIJ,xx,yy);CHKERRQ(ierr); 795f4a53059SBarry Smith ierr = VecScatterEnd(xx,b->w,INSERT_VALUES,SCATTER_FORWARD,b->ctx);CHKERRQ(ierr); 7964cb9d9b8SBarry Smith ierr = (*b->OAIJ->ops->multadd)(b->OAIJ,b->w,yy,yy);CHKERRQ(ierr); 797f4a53059SBarry Smith PetscFunctionReturn(0); 798f4a53059SBarry Smith } 799f4a53059SBarry Smith 8004cb9d9b8SBarry Smith #undef __FUNC__ 8014cb9d9b8SBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_MPIMAIJ_dof"></a>*/"MatMultTranspose_MPIMAIJ_dof" 8024cb9d9b8SBarry Smith int MatMultTranspose_MPIMAIJ_dof(Mat A,Vec xx,Vec yy) 8034cb9d9b8SBarry Smith { 8044cb9d9b8SBarry Smith Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data; 8054cb9d9b8SBarry Smith int ierr; 8064cb9d9b8SBarry Smith PetscFunctionBegin; 8074cb9d9b8SBarry Smith ierr = (*b->OAIJ->ops->multtranspose)(b->OAIJ,xx,b->w);CHKERRQ(ierr); 8084cb9d9b8SBarry Smith ierr = VecScatterBegin(b->w,yy,ADD_VALUES,SCATTER_REVERSE,b->ctx);CHKERRQ(ierr); 8094cb9d9b8SBarry Smith ierr = (*b->AIJ->ops->multtranspose)(b->AIJ,xx,yy);CHKERRQ(ierr); 8104cb9d9b8SBarry Smith ierr = VecScatterEnd(b->w,yy,ADD_VALUES,SCATTER_REVERSE,b->ctx);CHKERRQ(ierr); 8114cb9d9b8SBarry Smith PetscFunctionReturn(0); 8124cb9d9b8SBarry Smith } 8134cb9d9b8SBarry Smith 8144cb9d9b8SBarry Smith #undef __FUNC__ 8154cb9d9b8SBarry Smith #define __FUNC__ /*<a name="MatMultAdd_MPIMAIJ_dof"></a>*/"MatMultAdd_MPIMAIJ_dof" 816d72c5c08SBarry Smith int MatMultAdd_MPIMAIJ_dof(Mat A,Vec xx,Vec yy,Vec zz) 8174cb9d9b8SBarry Smith { 8184cb9d9b8SBarry Smith Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data; 8194cb9d9b8SBarry Smith int ierr; 8204cb9d9b8SBarry Smith PetscFunctionBegin; 8214cb9d9b8SBarry Smith 8224cb9d9b8SBarry Smith /* start the scatter */ 8234cb9d9b8SBarry Smith ierr = VecScatterBegin(xx,b->w,INSERT_VALUES,SCATTER_FORWARD,b->ctx);CHKERRQ(ierr); 824d72c5c08SBarry Smith ierr = (*b->AIJ->ops->multadd)(b->AIJ,xx,yy,zz);CHKERRQ(ierr); 8254cb9d9b8SBarry Smith ierr = VecScatterEnd(xx,b->w,INSERT_VALUES,SCATTER_FORWARD,b->ctx);CHKERRQ(ierr); 826d72c5c08SBarry Smith ierr = (*b->OAIJ->ops->multadd)(b->OAIJ,b->w,yy,zz);CHKERRQ(ierr); 8274cb9d9b8SBarry Smith PetscFunctionReturn(0); 8284cb9d9b8SBarry Smith } 8294cb9d9b8SBarry Smith 8304cb9d9b8SBarry Smith #undef __FUNC__ 8314cb9d9b8SBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_MPIMAIJ_dof"></a>*/"MatMultTransposeAdd_MPIMAIJ_dof" 832d72c5c08SBarry Smith int MatMultTransposeAdd_MPIMAIJ_dof(Mat A,Vec xx,Vec yy,Vec zz) 8334cb9d9b8SBarry Smith { 8344cb9d9b8SBarry Smith Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data; 8354cb9d9b8SBarry Smith int ierr; 8364cb9d9b8SBarry Smith PetscFunctionBegin; 8374cb9d9b8SBarry Smith ierr = (*b->OAIJ->ops->multtranspose)(b->OAIJ,xx,b->w);CHKERRQ(ierr); 838d72c5c08SBarry Smith ierr = VecScatterBegin(b->w,zz,ADD_VALUES,SCATTER_REVERSE,b->ctx);CHKERRQ(ierr); 839d72c5c08SBarry Smith ierr = (*b->AIJ->ops->multtransposeadd)(b->AIJ,xx,yy,zz);CHKERRQ(ierr); 840d72c5c08SBarry Smith ierr = VecScatterEnd(b->w,zz,ADD_VALUES,SCATTER_REVERSE,b->ctx);CHKERRQ(ierr); 8414cb9d9b8SBarry Smith PetscFunctionReturn(0); 8424cb9d9b8SBarry Smith } 8434cb9d9b8SBarry Smith 844bcc973b7SBarry Smith /* ---------------------------------------------------------------------------------- */ 84582b1193eSBarry Smith #undef __FUNC__ 846cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatCreateMAIJ"></a>*/"MatCreateMAIJ" 847cd3bbe55SBarry Smith int MatCreateMAIJ(Mat A,int dof,Mat *maij) 84882b1193eSBarry Smith { 849f4a53059SBarry Smith int ierr,size,n; 8504cb9d9b8SBarry Smith Mat_MPIMAIJ *b; 85182b1193eSBarry Smith Mat B; 85282b1193eSBarry Smith 85382b1193eSBarry Smith PetscFunctionBegin; 854d72c5c08SBarry Smith ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 855d72c5c08SBarry Smith 856d72c5c08SBarry Smith if (dof == 1) { 857d72c5c08SBarry Smith *maij = A; 858d72c5c08SBarry Smith } else { 859cd3bbe55SBarry Smith ierr = MATCreate(A->comm,dof*A->m,dof*A->n,dof*A->M,dof*A->N,&B);CHKERRQ(ierr); 860cd3bbe55SBarry Smith B->assembled = PETSC_TRUE; 861d72c5c08SBarry Smith 862f4a53059SBarry Smith ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 863f4a53059SBarry Smith if (size == 1) { 864b9b97703SBarry Smith ierr = MatSetType(B,MATSEQMAIJ);CHKERRQ(ierr); 8654cb9d9b8SBarry Smith B->ops->destroy = MatDestroy_SeqMAIJ; 866b9b97703SBarry Smith b = (Mat_MPIMAIJ*)B->data; 867b9b97703SBarry Smith b->dof = dof; 8684cb9d9b8SBarry Smith b->AIJ = A; 869d72c5c08SBarry Smith if (dof == 2) { 870cd3bbe55SBarry Smith B->ops->mult = MatMult_SeqMAIJ_2; 871cd3bbe55SBarry Smith B->ops->multadd = MatMultAdd_SeqMAIJ_2; 872cd3bbe55SBarry Smith B->ops->multtranspose = MatMultTranspose_SeqMAIJ_2; 873cd3bbe55SBarry Smith B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_2; 874bcc973b7SBarry Smith } else if (dof == 3) { 875bcc973b7SBarry Smith B->ops->mult = MatMult_SeqMAIJ_3; 876bcc973b7SBarry Smith B->ops->multadd = MatMultAdd_SeqMAIJ_3; 877bcc973b7SBarry Smith B->ops->multtranspose = MatMultTranspose_SeqMAIJ_3; 878bcc973b7SBarry Smith B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_3; 879bcc973b7SBarry Smith } else if (dof == 4) { 880bcc973b7SBarry Smith B->ops->mult = MatMult_SeqMAIJ_4; 881bcc973b7SBarry Smith B->ops->multadd = MatMultAdd_SeqMAIJ_4; 882bcc973b7SBarry Smith B->ops->multtranspose = MatMultTranspose_SeqMAIJ_4; 883bcc973b7SBarry Smith B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_4; 884f9fae5adSBarry Smith } else if (dof == 5) { 885f9fae5adSBarry Smith B->ops->mult = MatMult_SeqMAIJ_5; 886f9fae5adSBarry Smith B->ops->multadd = MatMultAdd_SeqMAIJ_5; 887f9fae5adSBarry Smith B->ops->multtranspose = MatMultTranspose_SeqMAIJ_5; 888f9fae5adSBarry Smith B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_5; 88982b1193eSBarry Smith } else { 890*29bbc08cSBarry Smith SETERRQ1(1,"Cannot handle a dof of %d\n",dof); 89182b1193eSBarry Smith } 892f4a53059SBarry Smith } else { 893f4a53059SBarry Smith Mat_MPIAIJ *mpiaij = (Mat_MPIAIJ *)A->data; 894f4a53059SBarry Smith IS from,to; 895f4a53059SBarry Smith Vec gvec; 896f4a53059SBarry Smith int *garray,i; 897f4a53059SBarry Smith 898b9b97703SBarry Smith ierr = MatSetType(B,MATMPIMAIJ);CHKERRQ(ierr); 899d72c5c08SBarry Smith B->ops->destroy = MatDestroy_MPIMAIJ; 900b9b97703SBarry Smith b = (Mat_MPIMAIJ*)B->data; 901b9b97703SBarry Smith b->dof = dof; 902b9b97703SBarry Smith b->A = A; 9034cb9d9b8SBarry Smith ierr = MatCreateMAIJ(mpiaij->A,dof,&b->AIJ);CHKERRQ(ierr); 9044cb9d9b8SBarry Smith ierr = MatCreateMAIJ(mpiaij->B,dof,&b->OAIJ);CHKERRQ(ierr); 9054cb9d9b8SBarry Smith 906f4a53059SBarry Smith ierr = VecGetSize(mpiaij->lvec,&n);CHKERRQ(ierr); 907f4a53059SBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,n*dof,&b->w);CHKERRQ(ierr); 908f4a53059SBarry Smith 909f4a53059SBarry Smith /* create two temporary Index sets for build scatter gather */ 910f4a53059SBarry Smith garray = (int*)PetscMalloc((n+1)*sizeof(int));CHKPTRQ(garray); 911f4a53059SBarry Smith for (i=0; i<n; i++) garray[i] = dof*mpiaij->garray[i]; 912f4a53059SBarry Smith ierr = ISCreateBlock(A->comm,dof,n,garray,&from);CHKERRQ(ierr); 913f4a53059SBarry Smith ierr = PetscFree(garray);CHKERRQ(ierr); 914f4a53059SBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,n*dof,0,1,&to);CHKERRQ(ierr); 915f4a53059SBarry Smith 916f4a53059SBarry Smith /* create temporary global vector to generate scatter context */ 917f4a53059SBarry Smith ierr = VecCreateMPI(A->comm,dof*A->n,dof*A->N,&gvec);CHKERRQ(ierr); 918f4a53059SBarry Smith 919f4a53059SBarry Smith /* generate the scatter context */ 920f4a53059SBarry Smith ierr = VecScatterCreate(gvec,from,b->w,to,&b->ctx);CHKERRQ(ierr); 921f4a53059SBarry Smith 922f4a53059SBarry Smith ierr = ISDestroy(from);CHKERRQ(ierr); 923f4a53059SBarry Smith ierr = ISDestroy(to);CHKERRQ(ierr); 924f4a53059SBarry Smith ierr = VecDestroy(gvec);CHKERRQ(ierr); 925f4a53059SBarry Smith 926f4a53059SBarry Smith B->ops->mult = MatMult_MPIMAIJ_dof; 9274cb9d9b8SBarry Smith B->ops->multtranspose = MatMultTranspose_MPIMAIJ_dof; 9284cb9d9b8SBarry Smith B->ops->multadd = MatMultAdd_MPIMAIJ_dof; 9294cb9d9b8SBarry Smith B->ops->multtransposeadd = MatMultTransposeAdd_MPIMAIJ_dof; 930f4a53059SBarry Smith } 931cd3bbe55SBarry Smith *maij = B; 932d72c5c08SBarry Smith } 93382b1193eSBarry Smith PetscFunctionReturn(0); 93482b1193eSBarry Smith } 93582b1193eSBarry Smith 93682b1193eSBarry Smith 93782b1193eSBarry Smith 93882b1193eSBarry Smith 93982b1193eSBarry Smith 94082b1193eSBarry Smith 94182b1193eSBarry Smith 94282b1193eSBarry Smith 94382b1193eSBarry Smith 94482b1193eSBarry Smith 94582b1193eSBarry Smith 94682b1193eSBarry Smith 947