1*b0a32e0cSBarry Smith /*$Id: maij.c,v 1.12 2000/12/01 03:12:25 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__ 35*b0a32e0cSBarry Smith #define __FUNC__ "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__ 59*b0a32e0cSBarry Smith #define __FUNC__ "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__ 72*b0a32e0cSBarry Smith #define __FUNC__ "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 PetscFunctionReturn(0); 844cb9d9b8SBarry Smith } 854cb9d9b8SBarry Smith 864cb9d9b8SBarry Smith #undef __FUNC__ 87*b0a32e0cSBarry Smith #define __FUNC__ "MatDestroy_MPIMAIJ" 884cb9d9b8SBarry Smith int MatDestroy_MPIMAIJ(Mat A) 894cb9d9b8SBarry Smith { 904cb9d9b8SBarry Smith int ierr; 914cb9d9b8SBarry Smith Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data; 924cb9d9b8SBarry Smith 934cb9d9b8SBarry Smith PetscFunctionBegin; 944cb9d9b8SBarry Smith if (b->AIJ) { 954cb9d9b8SBarry Smith ierr = MatDestroy(b->AIJ);CHKERRQ(ierr); 964cb9d9b8SBarry Smith } 97f4a53059SBarry Smith if (b->OAIJ) { 98f4a53059SBarry Smith ierr = MatDestroy(b->OAIJ);CHKERRQ(ierr); 99f4a53059SBarry Smith } 100b9b97703SBarry Smith if (b->A) { 101b9b97703SBarry Smith ierr = MatDestroy(b->A);CHKERRQ(ierr); 102b9b97703SBarry Smith } 103f4a53059SBarry Smith if (b->ctx) { 104f4a53059SBarry Smith ierr = VecScatterDestroy(b->ctx);CHKERRQ(ierr); 105f4a53059SBarry Smith } 106f4a53059SBarry Smith if (b->w) { 107f4a53059SBarry Smith ierr = VecDestroy(b->w);CHKERRQ(ierr); 108f4a53059SBarry Smith } 109cd3bbe55SBarry Smith ierr = PetscFree(b);CHKERRQ(ierr); 110b9b97703SBarry Smith PetscFunctionReturn(0); 111b9b97703SBarry Smith } 112b9b97703SBarry Smith 11382b1193eSBarry Smith EXTERN_C_BEGIN 11482b1193eSBarry Smith #undef __FUNC__ 115*b0a32e0cSBarry Smith #define __FUNC__ "MatCreate_MAIJ" 116f4a53059SBarry Smith int MatCreate_MAIJ(Mat A) 11782b1193eSBarry Smith { 118cd3bbe55SBarry Smith int ierr; 1194cb9d9b8SBarry Smith Mat_MPIMAIJ *b; 12082b1193eSBarry Smith 12182b1193eSBarry Smith PetscFunctionBegin; 122*b0a32e0cSBarry Smith ierr = PetscNew(Mat_MPIMAIJ,&b);CHKERRQ(ierr); 123*b0a32e0cSBarry Smith A->data = (void*)b; 1244cb9d9b8SBarry Smith ierr = PetscMemzero(b,sizeof(Mat_MPIMAIJ));CHKERRQ(ierr); 125cd3bbe55SBarry Smith ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 126cd3bbe55SBarry Smith A->factor = 0; 127cd3bbe55SBarry Smith A->mapping = 0; 128f4a53059SBarry Smith 129cd3bbe55SBarry Smith b->AIJ = 0; 130cd3bbe55SBarry Smith b->dof = 0; 131f4a53059SBarry Smith b->OAIJ = 0; 132f4a53059SBarry Smith b->ctx = 0; 133f4a53059SBarry Smith b->w = 0; 13482b1193eSBarry Smith PetscFunctionReturn(0); 13582b1193eSBarry Smith } 13682b1193eSBarry Smith EXTERN_C_END 13782b1193eSBarry Smith 138cd3bbe55SBarry Smith /* --------------------------------------------------------------------------------------*/ 13982b1193eSBarry Smith #undef __FUNC__ 140cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMult_SeqMAIJ_2"></a>*/"MatMult_SeqMAIJ_2" 141cd3bbe55SBarry Smith int MatMult_SeqMAIJ_2(Mat A,Vec xx,Vec yy) 14282b1193eSBarry Smith { 1434cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 144bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 145bcc973b7SBarry Smith Scalar *x,*y,*v,sum1, sum2; 146273d9f13SBarry Smith int ierr,m = b->AIJ->m,*idx,shift = a->indexshift,*ii; 147bcc973b7SBarry Smith int n,i,jrow,j; 14882b1193eSBarry Smith 149bcc973b7SBarry Smith PetscFunctionBegin; 150bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 151bcc973b7SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 152bcc973b7SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 153bcc973b7SBarry Smith idx = a->j; 154bcc973b7SBarry Smith v = a->a; 155bcc973b7SBarry Smith ii = a->i; 156bcc973b7SBarry Smith 157bcc973b7SBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 158bcc973b7SBarry Smith idx += shift; 159bcc973b7SBarry Smith for (i=0; i<m; i++) { 160bcc973b7SBarry Smith jrow = ii[i]; 161bcc973b7SBarry Smith n = ii[i+1] - jrow; 162bcc973b7SBarry Smith sum1 = 0.0; 163bcc973b7SBarry Smith sum2 = 0.0; 164bcc973b7SBarry Smith for (j=0; j<n; j++) { 165bcc973b7SBarry Smith sum1 += v[jrow]*x[2*idx[jrow]]; 166bcc973b7SBarry Smith sum2 += v[jrow]*x[2*idx[jrow]+1]; 167bcc973b7SBarry Smith jrow++; 168bcc973b7SBarry Smith } 169bcc973b7SBarry Smith y[2*i] = sum1; 170bcc973b7SBarry Smith y[2*i+1] = sum2; 171bcc973b7SBarry Smith } 172bcc973b7SBarry Smith 173*b0a32e0cSBarry Smith PetscLogFlops(4*a->nz - 2*m); 174bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 175bcc973b7SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 17682b1193eSBarry Smith PetscFunctionReturn(0); 17782b1193eSBarry Smith } 178bcc973b7SBarry Smith 17982b1193eSBarry Smith #undef __FUNC__ 180cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_SeqMAIJ_2"></a>*/"MatMultTranspose_SeqMAIJ_2" 181cd3bbe55SBarry Smith int MatMultTranspose_SeqMAIJ_2(Mat A,Vec xx,Vec yy) 18282b1193eSBarry Smith { 1834cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 184bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 185bcc973b7SBarry Smith Scalar *x,*y,*v,alpha1,alpha2,zero = 0.0; 186273d9f13SBarry Smith int ierr,m = b->AIJ->m,n,i,*idx,shift = a->indexshift; 18782b1193eSBarry Smith 188bcc973b7SBarry Smith PetscFunctionBegin; 189bcc973b7SBarry Smith ierr = VecSet(&zero,yy);CHKERRQ(ierr); 190bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 191bcc973b7SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 192bcc973b7SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 193bcc973b7SBarry Smith for (i=0; i<m; i++) { 194bcc973b7SBarry Smith idx = a->j + a->i[i] + shift; 195bcc973b7SBarry Smith v = a->a + a->i[i] + shift; 196bcc973b7SBarry Smith n = a->i[i+1] - a->i[i]; 197bcc973b7SBarry Smith alpha1 = x[2*i]; 198bcc973b7SBarry Smith alpha2 = x[2*i+1]; 199bcc973b7SBarry Smith while (n-->0) {y[2*(*idx)] += alpha1*(*v); y[2*(*idx)+1] += alpha2*(*v); idx++; v++;} 200bcc973b7SBarry Smith } 201*b0a32e0cSBarry Smith PetscLogFlops(4*a->nz - 2*b->AIJ->n); 202bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 203bcc973b7SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 20482b1193eSBarry Smith PetscFunctionReturn(0); 20582b1193eSBarry Smith } 206bcc973b7SBarry Smith 20782b1193eSBarry Smith #undef __FUNC__ 208cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMultAdd_SeqMAIJ_2"></a>*/"MatMultAdd_SeqMAIJ_2" 209cd3bbe55SBarry Smith int MatMultAdd_SeqMAIJ_2(Mat A,Vec xx,Vec yy,Vec zz) 21082b1193eSBarry Smith { 2114cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 212bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 213bcc973b7SBarry Smith Scalar *x,*y,*v,sum1, sum2; 214273d9f13SBarry Smith int ierr,m = b->AIJ->m,*idx,shift = a->indexshift,*ii; 215bcc973b7SBarry Smith int n,i,jrow,j; 21682b1193eSBarry Smith 217bcc973b7SBarry Smith PetscFunctionBegin; 218f9fae5adSBarry Smith if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);} 219bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 220f9fae5adSBarry Smith ierr = VecGetArray(zz,&y);CHKERRQ(ierr); 221bcc973b7SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 222bcc973b7SBarry Smith idx = a->j; 223bcc973b7SBarry Smith v = a->a; 224bcc973b7SBarry Smith ii = a->i; 225bcc973b7SBarry Smith 226bcc973b7SBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 227bcc973b7SBarry Smith idx += shift; 228bcc973b7SBarry Smith for (i=0; i<m; i++) { 229bcc973b7SBarry Smith jrow = ii[i]; 230bcc973b7SBarry Smith n = ii[i+1] - jrow; 231bcc973b7SBarry Smith sum1 = 0.0; 232bcc973b7SBarry Smith sum2 = 0.0; 233bcc973b7SBarry Smith for (j=0; j<n; j++) { 234bcc973b7SBarry Smith sum1 += v[jrow]*x[2*idx[jrow]]; 235bcc973b7SBarry Smith sum2 += v[jrow]*x[2*idx[jrow]+1]; 236bcc973b7SBarry Smith jrow++; 237bcc973b7SBarry Smith } 238bcc973b7SBarry Smith y[2*i] += sum1; 239bcc973b7SBarry Smith y[2*i+1] += sum2; 240bcc973b7SBarry Smith } 241bcc973b7SBarry Smith 242*b0a32e0cSBarry Smith PetscLogFlops(4*a->nz - 2*m); 243bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 244f9fae5adSBarry Smith ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr); 245bcc973b7SBarry Smith PetscFunctionReturn(0); 24682b1193eSBarry Smith } 24782b1193eSBarry Smith #undef __FUNC__ 248cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_SeqMAIJ_2"></a>*/"MatMultTransposeAdd_SeqMAIJ_2" 249cd3bbe55SBarry Smith int MatMultTransposeAdd_SeqMAIJ_2(Mat A,Vec xx,Vec yy,Vec zz) 25082b1193eSBarry Smith { 2514cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 252bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 253bcc973b7SBarry Smith Scalar *x,*y,*v,alpha1,alpha2; 254273d9f13SBarry Smith int ierr,m = b->AIJ->m,n,i,*idx,shift = a->indexshift; 25582b1193eSBarry Smith 256bcc973b7SBarry Smith PetscFunctionBegin; 257f9fae5adSBarry Smith if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);} 258bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 259f9fae5adSBarry Smith ierr = VecGetArray(zz,&y);CHKERRQ(ierr); 260bcc973b7SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 261bcc973b7SBarry Smith for (i=0; i<m; i++) { 262bcc973b7SBarry Smith idx = a->j + a->i[i] + shift; 263bcc973b7SBarry Smith v = a->a + a->i[i] + shift; 264bcc973b7SBarry Smith n = a->i[i+1] - a->i[i]; 265bcc973b7SBarry Smith alpha1 = x[2*i]; 266bcc973b7SBarry Smith alpha2 = x[2*i+1]; 267bcc973b7SBarry Smith while (n-->0) {y[2*(*idx)] += alpha1*(*v); y[2*(*idx)+1] += alpha2*(*v); idx++; v++;} 268bcc973b7SBarry Smith } 269*b0a32e0cSBarry Smith PetscLogFlops(4*a->nz - 2*b->AIJ->n); 270bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 271f9fae5adSBarry Smith ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr); 272bcc973b7SBarry Smith PetscFunctionReturn(0); 27382b1193eSBarry Smith } 274cd3bbe55SBarry Smith /* --------------------------------------------------------------------------------------*/ 275bcc973b7SBarry Smith #undef __FUNC__ 276bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMult_SeqMAIJ_3"></a>*/"MatMult_SeqMAIJ_3" 277bcc973b7SBarry Smith int MatMult_SeqMAIJ_3(Mat A,Vec xx,Vec yy) 278bcc973b7SBarry Smith { 2794cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 280bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 281bcc973b7SBarry Smith Scalar *x,*y,*v,sum1, sum2, sum3; 282273d9f13SBarry Smith int ierr,m = b->AIJ->m,*idx,shift = a->indexshift,*ii; 283bcc973b7SBarry Smith int n,i,jrow,j; 28482b1193eSBarry Smith 285bcc973b7SBarry Smith PetscFunctionBegin; 286bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 287bcc973b7SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 288bcc973b7SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 289bcc973b7SBarry Smith idx = a->j; 290bcc973b7SBarry Smith v = a->a; 291bcc973b7SBarry Smith ii = a->i; 292bcc973b7SBarry Smith 293bcc973b7SBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 294bcc973b7SBarry Smith idx += shift; 295bcc973b7SBarry Smith for (i=0; i<m; i++) { 296bcc973b7SBarry Smith jrow = ii[i]; 297bcc973b7SBarry Smith n = ii[i+1] - jrow; 298bcc973b7SBarry Smith sum1 = 0.0; 299bcc973b7SBarry Smith sum2 = 0.0; 300bcc973b7SBarry Smith sum3 = 0.0; 301bcc973b7SBarry Smith for (j=0; j<n; j++) { 302bcc973b7SBarry Smith sum1 += v[jrow]*x[3*idx[jrow]]; 303bcc973b7SBarry Smith sum2 += v[jrow]*x[3*idx[jrow]+1]; 304bcc973b7SBarry Smith sum3 += v[jrow]*x[3*idx[jrow]+2]; 305bcc973b7SBarry Smith jrow++; 306bcc973b7SBarry Smith } 307bcc973b7SBarry Smith y[3*i] = sum1; 308bcc973b7SBarry Smith y[3*i+1] = sum2; 309bcc973b7SBarry Smith y[3*i+2] = sum3; 310bcc973b7SBarry Smith } 311bcc973b7SBarry Smith 312*b0a32e0cSBarry Smith PetscLogFlops(6*a->nz - 3*m); 313bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 314bcc973b7SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 315bcc973b7SBarry Smith PetscFunctionReturn(0); 316bcc973b7SBarry Smith } 317bcc973b7SBarry Smith 318bcc973b7SBarry Smith #undef __FUNC__ 319bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_SeqMAIJ_3"></a>*/"MatMultTranspose_SeqMAIJ_3" 320bcc973b7SBarry Smith int MatMultTranspose_SeqMAIJ_3(Mat A,Vec xx,Vec yy) 321bcc973b7SBarry Smith { 3224cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 323bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 324bcc973b7SBarry Smith Scalar *x,*y,*v,alpha1,alpha2,alpha3,zero = 0.0; 325273d9f13SBarry Smith int ierr,m = b->AIJ->m,n,i,*idx,shift = a->indexshift; 326bcc973b7SBarry Smith 327bcc973b7SBarry Smith PetscFunctionBegin; 328bcc973b7SBarry Smith ierr = VecSet(&zero,yy);CHKERRQ(ierr); 329bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 330bcc973b7SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 331bcc973b7SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 332bcc973b7SBarry Smith for (i=0; i<m; i++) { 333bcc973b7SBarry Smith idx = a->j + a->i[i] + shift; 334bcc973b7SBarry Smith v = a->a + a->i[i] + shift; 335bcc973b7SBarry Smith n = a->i[i+1] - a->i[i]; 336bcc973b7SBarry Smith alpha1 = x[3*i]; 337bcc973b7SBarry Smith alpha2 = x[3*i+1]; 338bcc973b7SBarry Smith alpha3 = x[3*i+2]; 339bcc973b7SBarry Smith while (n-->0) { 340bcc973b7SBarry Smith y[3*(*idx)] += alpha1*(*v); 341bcc973b7SBarry Smith y[3*(*idx)+1] += alpha2*(*v); 342bcc973b7SBarry Smith y[3*(*idx)+2] += alpha3*(*v); 343bcc973b7SBarry Smith idx++; v++; 344bcc973b7SBarry Smith } 345bcc973b7SBarry Smith } 346*b0a32e0cSBarry Smith PetscLogFlops(6*a->nz - 3*b->AIJ->n); 347bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 348bcc973b7SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 349bcc973b7SBarry Smith PetscFunctionReturn(0); 350bcc973b7SBarry Smith } 351bcc973b7SBarry Smith 352bcc973b7SBarry Smith #undef __FUNC__ 353bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultAdd_SeqMAIJ_3"></a>*/"MatMultAdd_SeqMAIJ_3" 354bcc973b7SBarry Smith int MatMultAdd_SeqMAIJ_3(Mat A,Vec xx,Vec yy,Vec zz) 355bcc973b7SBarry Smith { 3564cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 357bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 358bcc973b7SBarry Smith Scalar *x,*y,*v,sum1, sum2, sum3; 359273d9f13SBarry Smith int ierr,m = b->AIJ->m,*idx,shift = a->indexshift,*ii; 360bcc973b7SBarry Smith int n,i,jrow,j; 361bcc973b7SBarry Smith 362bcc973b7SBarry Smith PetscFunctionBegin; 363f9fae5adSBarry Smith if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);} 364bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 365f9fae5adSBarry Smith ierr = VecGetArray(zz,&y);CHKERRQ(ierr); 366bcc973b7SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 367bcc973b7SBarry Smith idx = a->j; 368bcc973b7SBarry Smith v = a->a; 369bcc973b7SBarry Smith ii = a->i; 370bcc973b7SBarry Smith 371bcc973b7SBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 372bcc973b7SBarry Smith idx += shift; 373bcc973b7SBarry Smith for (i=0; i<m; i++) { 374bcc973b7SBarry Smith jrow = ii[i]; 375bcc973b7SBarry Smith n = ii[i+1] - jrow; 376bcc973b7SBarry Smith sum1 = 0.0; 377bcc973b7SBarry Smith sum2 = 0.0; 378bcc973b7SBarry Smith sum3 = 0.0; 379bcc973b7SBarry Smith for (j=0; j<n; j++) { 380bcc973b7SBarry Smith sum1 += v[jrow]*x[3*idx[jrow]]; 381bcc973b7SBarry Smith sum2 += v[jrow]*x[3*idx[jrow]+1]; 382bcc973b7SBarry Smith sum3 += v[jrow]*x[3*idx[jrow]+2]; 383bcc973b7SBarry Smith jrow++; 384bcc973b7SBarry Smith } 385bcc973b7SBarry Smith y[3*i] += sum1; 386bcc973b7SBarry Smith y[3*i+1] += sum2; 387bcc973b7SBarry Smith y[3*i+2] += sum3; 388bcc973b7SBarry Smith } 389bcc973b7SBarry Smith 390*b0a32e0cSBarry Smith PetscLogFlops(6*a->nz); 391bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 392f9fae5adSBarry Smith ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr); 393bcc973b7SBarry Smith PetscFunctionReturn(0); 394bcc973b7SBarry Smith } 395bcc973b7SBarry Smith #undef __FUNC__ 396bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_SeqMAIJ_3"></a>*/"MatMultTransposeAdd_SeqMAIJ_3" 397bcc973b7SBarry Smith int MatMultTransposeAdd_SeqMAIJ_3(Mat A,Vec xx,Vec yy,Vec zz) 398bcc973b7SBarry Smith { 3994cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 400bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 401bcc973b7SBarry Smith Scalar *x,*y,*v,alpha1,alpha2,alpha3; 402273d9f13SBarry Smith int ierr,m = b->AIJ->m,n,i,*idx,shift = a->indexshift; 403bcc973b7SBarry Smith 404bcc973b7SBarry Smith PetscFunctionBegin; 405f9fae5adSBarry Smith if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);} 406bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 407f9fae5adSBarry Smith ierr = VecGetArray(zz,&y);CHKERRQ(ierr); 408bcc973b7SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 409bcc973b7SBarry Smith for (i=0; i<m; i++) { 410bcc973b7SBarry Smith idx = a->j + a->i[i] + shift; 411bcc973b7SBarry Smith v = a->a + a->i[i] + shift; 412bcc973b7SBarry Smith n = a->i[i+1] - a->i[i]; 413bcc973b7SBarry Smith alpha1 = x[3*i]; 414bcc973b7SBarry Smith alpha2 = x[3*i+1]; 415bcc973b7SBarry Smith alpha3 = x[3*i+2]; 416bcc973b7SBarry Smith while (n-->0) { 417bcc973b7SBarry Smith y[3*(*idx)] += alpha1*(*v); 418bcc973b7SBarry Smith y[3*(*idx)+1] += alpha2*(*v); 419bcc973b7SBarry Smith y[3*(*idx)+2] += alpha3*(*v); 420bcc973b7SBarry Smith idx++; v++; 421bcc973b7SBarry Smith } 422bcc973b7SBarry Smith } 423*b0a32e0cSBarry Smith PetscLogFlops(6*a->nz); 424bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 425f9fae5adSBarry Smith ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr); 426bcc973b7SBarry Smith PetscFunctionReturn(0); 427bcc973b7SBarry Smith } 428bcc973b7SBarry Smith 429bcc973b7SBarry Smith /* ------------------------------------------------------------------------------*/ 430bcc973b7SBarry Smith #undef __FUNC__ 431bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMult_SeqMAIJ_4"></a>*/"MatMult_SeqMAIJ_4" 432bcc973b7SBarry Smith int MatMult_SeqMAIJ_4(Mat A,Vec xx,Vec yy) 433bcc973b7SBarry Smith { 4344cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 435bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 436bcc973b7SBarry Smith Scalar *x,*y,*v,sum1, sum2, sum3, sum4; 437273d9f13SBarry Smith int ierr,m = b->AIJ->m,*idx,shift = a->indexshift,*ii; 438bcc973b7SBarry Smith int n,i,jrow,j; 439bcc973b7SBarry Smith 440bcc973b7SBarry Smith PetscFunctionBegin; 441bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 442bcc973b7SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 443bcc973b7SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 444bcc973b7SBarry Smith idx = a->j; 445bcc973b7SBarry Smith v = a->a; 446bcc973b7SBarry Smith ii = a->i; 447bcc973b7SBarry Smith 448bcc973b7SBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 449bcc973b7SBarry Smith idx += shift; 450bcc973b7SBarry Smith for (i=0; i<m; i++) { 451bcc973b7SBarry Smith jrow = ii[i]; 452bcc973b7SBarry Smith n = ii[i+1] - jrow; 453bcc973b7SBarry Smith sum1 = 0.0; 454bcc973b7SBarry Smith sum2 = 0.0; 455bcc973b7SBarry Smith sum3 = 0.0; 456bcc973b7SBarry Smith sum4 = 0.0; 457bcc973b7SBarry Smith for (j=0; j<n; j++) { 458bcc973b7SBarry Smith sum1 += v[jrow]*x[4*idx[jrow]]; 459bcc973b7SBarry Smith sum2 += v[jrow]*x[4*idx[jrow]+1]; 460bcc973b7SBarry Smith sum3 += v[jrow]*x[4*idx[jrow]+2]; 461bcc973b7SBarry Smith sum4 += v[jrow]*x[4*idx[jrow]+3]; 462bcc973b7SBarry Smith jrow++; 463bcc973b7SBarry Smith } 464bcc973b7SBarry Smith y[4*i] = sum1; 465bcc973b7SBarry Smith y[4*i+1] = sum2; 466bcc973b7SBarry Smith y[4*i+2] = sum3; 467bcc973b7SBarry Smith y[4*i+3] = sum4; 468bcc973b7SBarry Smith } 469bcc973b7SBarry Smith 470*b0a32e0cSBarry Smith PetscLogFlops(8*a->nz - 4*m); 471bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 472bcc973b7SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 473bcc973b7SBarry Smith PetscFunctionReturn(0); 474bcc973b7SBarry Smith } 475bcc973b7SBarry Smith 476bcc973b7SBarry Smith #undef __FUNC__ 477bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_SeqMAIJ_4"></a>*/"MatMultTranspose_SeqMAIJ_4" 478bcc973b7SBarry Smith int MatMultTranspose_SeqMAIJ_4(Mat A,Vec xx,Vec yy) 479bcc973b7SBarry Smith { 4804cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 481bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 482bcc973b7SBarry Smith Scalar *x,*y,*v,alpha1,alpha2,alpha3,alpha4,zero = 0.0; 483273d9f13SBarry Smith int ierr,m = b->AIJ->m,n,i,*idx,shift = a->indexshift; 484bcc973b7SBarry Smith 485bcc973b7SBarry Smith PetscFunctionBegin; 486bcc973b7SBarry Smith ierr = VecSet(&zero,yy);CHKERRQ(ierr); 487bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 488bcc973b7SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 489bcc973b7SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 490bcc973b7SBarry Smith for (i=0; i<m; i++) { 491bcc973b7SBarry Smith idx = a->j + a->i[i] + shift; 492bcc973b7SBarry Smith v = a->a + a->i[i] + shift; 493bcc973b7SBarry Smith n = a->i[i+1] - a->i[i]; 494bcc973b7SBarry Smith alpha1 = x[4*i]; 495bcc973b7SBarry Smith alpha2 = x[4*i+1]; 496bcc973b7SBarry Smith alpha3 = x[4*i+2]; 497bcc973b7SBarry Smith alpha4 = x[4*i+3]; 498bcc973b7SBarry Smith while (n-->0) { 499bcc973b7SBarry Smith y[4*(*idx)] += alpha1*(*v); 500bcc973b7SBarry Smith y[4*(*idx)+1] += alpha2*(*v); 501bcc973b7SBarry Smith y[4*(*idx)+2] += alpha3*(*v); 502bcc973b7SBarry Smith y[4*(*idx)+3] += alpha4*(*v); 503bcc973b7SBarry Smith idx++; v++; 504bcc973b7SBarry Smith } 505bcc973b7SBarry Smith } 506*b0a32e0cSBarry Smith PetscLogFlops(8*a->nz - 4*b->AIJ->n); 507bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 508bcc973b7SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 509bcc973b7SBarry Smith PetscFunctionReturn(0); 510bcc973b7SBarry Smith } 511bcc973b7SBarry Smith 512bcc973b7SBarry Smith #undef __FUNC__ 513bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultAdd_SeqMAIJ_4"></a>*/"MatMultAdd_SeqMAIJ_4" 514bcc973b7SBarry Smith int MatMultAdd_SeqMAIJ_4(Mat A,Vec xx,Vec yy,Vec zz) 515bcc973b7SBarry Smith { 5164cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 517f9fae5adSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 518f9fae5adSBarry Smith Scalar *x,*y,*v,sum1, sum2, sum3, sum4; 519273d9f13SBarry Smith int ierr,m = b->AIJ->m,*idx,shift = a->indexshift,*ii; 520f9fae5adSBarry Smith int n,i,jrow,j; 521f9fae5adSBarry Smith 522f9fae5adSBarry Smith PetscFunctionBegin; 523f9fae5adSBarry Smith if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);} 524f9fae5adSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 525f9fae5adSBarry Smith ierr = VecGetArray(zz,&y);CHKERRQ(ierr); 526f9fae5adSBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 527f9fae5adSBarry Smith idx = a->j; 528f9fae5adSBarry Smith v = a->a; 529f9fae5adSBarry Smith ii = a->i; 530f9fae5adSBarry Smith 531f9fae5adSBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 532f9fae5adSBarry Smith idx += shift; 533f9fae5adSBarry Smith for (i=0; i<m; i++) { 534f9fae5adSBarry Smith jrow = ii[i]; 535f9fae5adSBarry Smith n = ii[i+1] - jrow; 536f9fae5adSBarry Smith sum1 = 0.0; 537f9fae5adSBarry Smith sum2 = 0.0; 538f9fae5adSBarry Smith sum3 = 0.0; 539f9fae5adSBarry Smith sum4 = 0.0; 540f9fae5adSBarry Smith for (j=0; j<n; j++) { 541f9fae5adSBarry Smith sum1 += v[jrow]*x[4*idx[jrow]]; 542f9fae5adSBarry Smith sum2 += v[jrow]*x[4*idx[jrow]+1]; 543f9fae5adSBarry Smith sum3 += v[jrow]*x[4*idx[jrow]+2]; 544f9fae5adSBarry Smith sum4 += v[jrow]*x[4*idx[jrow]+3]; 545f9fae5adSBarry Smith jrow++; 546f9fae5adSBarry Smith } 547f9fae5adSBarry Smith y[4*i] += sum1; 548f9fae5adSBarry Smith y[4*i+1] += sum2; 549f9fae5adSBarry Smith y[4*i+2] += sum3; 550f9fae5adSBarry Smith y[4*i+3] += sum4; 551f9fae5adSBarry Smith } 552f9fae5adSBarry Smith 553*b0a32e0cSBarry Smith PetscLogFlops(8*a->nz - 4*m); 554f9fae5adSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 555f9fae5adSBarry Smith ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr); 556f9fae5adSBarry Smith PetscFunctionReturn(0); 557bcc973b7SBarry Smith } 558bcc973b7SBarry Smith #undef __FUNC__ 559bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_SeqMAIJ_4"></a>*/"MatMultTransposeAdd_SeqMAIJ_4" 560bcc973b7SBarry Smith int MatMultTransposeAdd_SeqMAIJ_4(Mat A,Vec xx,Vec yy,Vec zz) 561bcc973b7SBarry Smith { 5624cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 563f9fae5adSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 564f9fae5adSBarry Smith Scalar *x,*y,*v,alpha1,alpha2,alpha3,alpha4; 565273d9f13SBarry Smith int ierr,m = b->AIJ->m,n,i,*idx,shift = a->indexshift; 566f9fae5adSBarry Smith 567f9fae5adSBarry Smith PetscFunctionBegin; 568f9fae5adSBarry Smith if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);} 569f9fae5adSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 570f9fae5adSBarry Smith ierr = VecGetArray(zz,&y);CHKERRQ(ierr); 571f9fae5adSBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 572f9fae5adSBarry Smith for (i=0; i<m; i++) { 573f9fae5adSBarry Smith idx = a->j + a->i[i] + shift; 574f9fae5adSBarry Smith v = a->a + a->i[i] + shift; 575f9fae5adSBarry Smith n = a->i[i+1] - a->i[i]; 576f9fae5adSBarry Smith alpha1 = x[4*i]; 577f9fae5adSBarry Smith alpha2 = x[4*i+1]; 578f9fae5adSBarry Smith alpha3 = x[4*i+2]; 579f9fae5adSBarry Smith alpha4 = x[4*i+3]; 580f9fae5adSBarry Smith while (n-->0) { 581f9fae5adSBarry Smith y[4*(*idx)] += alpha1*(*v); 582f9fae5adSBarry Smith y[4*(*idx)+1] += alpha2*(*v); 583f9fae5adSBarry Smith y[4*(*idx)+2] += alpha3*(*v); 584f9fae5adSBarry Smith y[4*(*idx)+3] += alpha4*(*v); 585f9fae5adSBarry Smith idx++; v++; 586f9fae5adSBarry Smith } 587f9fae5adSBarry Smith } 588*b0a32e0cSBarry Smith PetscLogFlops(8*a->nz - 4*b->AIJ->n); 589f9fae5adSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 590f9fae5adSBarry Smith ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr); 591f9fae5adSBarry Smith PetscFunctionReturn(0); 592f9fae5adSBarry Smith } 593f9fae5adSBarry Smith /* ------------------------------------------------------------------------------*/ 5946cd98798SBarry Smith 595f9fae5adSBarry Smith #undef __FUNC__ 596f9fae5adSBarry Smith #define __FUNC__ /*<a name="MatMult_SeqMAIJ_5"></a>*/"MatMult_SeqMAIJ_5" 597f9fae5adSBarry Smith int MatMult_SeqMAIJ_5(Mat A,Vec xx,Vec yy) 598f9fae5adSBarry Smith { 5994cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 600f9fae5adSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 601f9fae5adSBarry Smith Scalar *x,*y,*v,sum1, sum2, sum3, sum4, sum5; 602273d9f13SBarry Smith int ierr,m = b->AIJ->m,*idx,shift = a->indexshift,*ii; 603f9fae5adSBarry Smith int n,i,jrow,j; 604f9fae5adSBarry Smith 605f9fae5adSBarry Smith PetscFunctionBegin; 606f9fae5adSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 607f9fae5adSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 608f9fae5adSBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 609f9fae5adSBarry Smith idx = a->j; 610f9fae5adSBarry Smith v = a->a; 611f9fae5adSBarry Smith ii = a->i; 612f9fae5adSBarry Smith 613f9fae5adSBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 614f9fae5adSBarry Smith idx += shift; 615f9fae5adSBarry Smith for (i=0; i<m; i++) { 616f9fae5adSBarry Smith jrow = ii[i]; 617f9fae5adSBarry Smith n = ii[i+1] - jrow; 618f9fae5adSBarry Smith sum1 = 0.0; 619f9fae5adSBarry Smith sum2 = 0.0; 620f9fae5adSBarry Smith sum3 = 0.0; 621f9fae5adSBarry Smith sum4 = 0.0; 622f9fae5adSBarry Smith sum5 = 0.0; 623f9fae5adSBarry Smith for (j=0; j<n; j++) { 624f9fae5adSBarry Smith sum1 += v[jrow]*x[5*idx[jrow]]; 625f9fae5adSBarry Smith sum2 += v[jrow]*x[5*idx[jrow]+1]; 626f9fae5adSBarry Smith sum3 += v[jrow]*x[5*idx[jrow]+2]; 627f9fae5adSBarry Smith sum4 += v[jrow]*x[5*idx[jrow]+3]; 628f9fae5adSBarry Smith sum5 += v[jrow]*x[5*idx[jrow]+4]; 629f9fae5adSBarry Smith jrow++; 630f9fae5adSBarry Smith } 631f9fae5adSBarry Smith y[5*i] = sum1; 632f9fae5adSBarry Smith y[5*i+1] = sum2; 633f9fae5adSBarry Smith y[5*i+2] = sum3; 634f9fae5adSBarry Smith y[5*i+3] = sum4; 635f9fae5adSBarry Smith y[5*i+4] = sum5; 636f9fae5adSBarry Smith } 637f9fae5adSBarry Smith 638*b0a32e0cSBarry Smith PetscLogFlops(10*a->nz - 5*m); 639f9fae5adSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 640f9fae5adSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 641f9fae5adSBarry Smith PetscFunctionReturn(0); 642f9fae5adSBarry Smith } 643f9fae5adSBarry Smith 644f9fae5adSBarry Smith #undef __FUNC__ 645f9fae5adSBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_SeqMAIJ_5"></a>*/"MatMultTranspose_SeqMAIJ_5" 646f9fae5adSBarry Smith int MatMultTranspose_SeqMAIJ_5(Mat A,Vec xx,Vec yy) 647f9fae5adSBarry Smith { 6484cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 649f9fae5adSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 650f9fae5adSBarry Smith Scalar *x,*y,*v,alpha1,alpha2,alpha3,alpha4,alpha5,zero = 0.0; 651273d9f13SBarry Smith int ierr,m = b->AIJ->m,n,i,*idx,shift = a->indexshift; 652f9fae5adSBarry Smith 653f9fae5adSBarry Smith PetscFunctionBegin; 654f9fae5adSBarry Smith ierr = VecSet(&zero,yy);CHKERRQ(ierr); 655f9fae5adSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 656f9fae5adSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 657f9fae5adSBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 658f9fae5adSBarry Smith for (i=0; i<m; i++) { 659f9fae5adSBarry Smith idx = a->j + a->i[i] + shift; 660f9fae5adSBarry Smith v = a->a + a->i[i] + shift; 661f9fae5adSBarry Smith n = a->i[i+1] - a->i[i]; 662f9fae5adSBarry Smith alpha1 = x[5*i]; 663f9fae5adSBarry Smith alpha2 = x[5*i+1]; 664f9fae5adSBarry Smith alpha3 = x[5*i+2]; 665f9fae5adSBarry Smith alpha4 = x[5*i+3]; 666f9fae5adSBarry Smith alpha5 = x[5*i+4]; 667f9fae5adSBarry Smith while (n-->0) { 668f9fae5adSBarry Smith y[5*(*idx)] += alpha1*(*v); 669f9fae5adSBarry Smith y[5*(*idx)+1] += alpha2*(*v); 670f9fae5adSBarry Smith y[5*(*idx)+2] += alpha3*(*v); 671f9fae5adSBarry Smith y[5*(*idx)+3] += alpha4*(*v); 672f9fae5adSBarry Smith y[5*(*idx)+4] += alpha5*(*v); 673f9fae5adSBarry Smith idx++; v++; 674f9fae5adSBarry Smith } 675f9fae5adSBarry Smith } 676*b0a32e0cSBarry Smith PetscLogFlops(10*a->nz - 5*b->AIJ->n); 677f9fae5adSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 678f9fae5adSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 679f9fae5adSBarry Smith PetscFunctionReturn(0); 680f9fae5adSBarry Smith } 681f9fae5adSBarry Smith 682f9fae5adSBarry Smith #undef __FUNC__ 683f9fae5adSBarry Smith #define __FUNC__ /*<a name="MatMultAdd_SeqMAIJ_5"></a>*/"MatMultAdd_SeqMAIJ_5" 684f9fae5adSBarry Smith int MatMultAdd_SeqMAIJ_5(Mat A,Vec xx,Vec yy,Vec zz) 685f9fae5adSBarry Smith { 6864cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 687f9fae5adSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 688f9fae5adSBarry Smith Scalar *x,*y,*v,sum1, sum2, sum3, sum4, sum5; 689273d9f13SBarry Smith int ierr,m = b->AIJ->m,*idx,shift = a->indexshift,*ii; 690f9fae5adSBarry Smith int n,i,jrow,j; 691f9fae5adSBarry Smith 692f9fae5adSBarry Smith PetscFunctionBegin; 693f9fae5adSBarry Smith if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);} 694f9fae5adSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 695f9fae5adSBarry Smith ierr = VecGetArray(zz,&y);CHKERRQ(ierr); 696f9fae5adSBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 697f9fae5adSBarry Smith idx = a->j; 698f9fae5adSBarry Smith v = a->a; 699f9fae5adSBarry Smith ii = a->i; 700f9fae5adSBarry Smith 701f9fae5adSBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 702f9fae5adSBarry Smith idx += shift; 703f9fae5adSBarry Smith for (i=0; i<m; i++) { 704f9fae5adSBarry Smith jrow = ii[i]; 705f9fae5adSBarry Smith n = ii[i+1] - jrow; 706f9fae5adSBarry Smith sum1 = 0.0; 707f9fae5adSBarry Smith sum2 = 0.0; 708f9fae5adSBarry Smith sum3 = 0.0; 709f9fae5adSBarry Smith sum4 = 0.0; 710f9fae5adSBarry Smith sum5 = 0.0; 711f9fae5adSBarry Smith for (j=0; j<n; j++) { 712f9fae5adSBarry Smith sum1 += v[jrow]*x[5*idx[jrow]]; 713f9fae5adSBarry Smith sum2 += v[jrow]*x[5*idx[jrow]+1]; 714f9fae5adSBarry Smith sum3 += v[jrow]*x[5*idx[jrow]+2]; 715f9fae5adSBarry Smith sum4 += v[jrow]*x[5*idx[jrow]+3]; 716f9fae5adSBarry Smith sum5 += v[jrow]*x[5*idx[jrow]+4]; 717f9fae5adSBarry Smith jrow++; 718f9fae5adSBarry Smith } 719f9fae5adSBarry Smith y[5*i] += sum1; 720f9fae5adSBarry Smith y[5*i+1] += sum2; 721f9fae5adSBarry Smith y[5*i+2] += sum3; 722f9fae5adSBarry Smith y[5*i+3] += sum4; 723f9fae5adSBarry Smith y[5*i+4] += sum5; 724f9fae5adSBarry Smith } 725f9fae5adSBarry Smith 726*b0a32e0cSBarry Smith PetscLogFlops(10*a->nz); 727f9fae5adSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 728f9fae5adSBarry Smith ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr); 729f9fae5adSBarry Smith PetscFunctionReturn(0); 730f9fae5adSBarry Smith } 731f9fae5adSBarry Smith 732f9fae5adSBarry Smith #undef __FUNC__ 733f9fae5adSBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_SeqMAIJ_5"></a>*/"MatMultTransposeAdd_SeqMAIJ_5" 734f9fae5adSBarry Smith int MatMultTransposeAdd_SeqMAIJ_5(Mat A,Vec xx,Vec yy,Vec zz) 735f9fae5adSBarry Smith { 7364cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 737f9fae5adSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 738f9fae5adSBarry Smith Scalar *x,*y,*v,alpha1,alpha2,alpha3,alpha4,alpha5; 739273d9f13SBarry Smith int ierr,m = b->AIJ->m,n,i,*idx,shift = a->indexshift; 740f9fae5adSBarry Smith 741f9fae5adSBarry Smith PetscFunctionBegin; 742f9fae5adSBarry Smith if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);} 743f9fae5adSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 744f9fae5adSBarry Smith ierr = VecGetArray(zz,&y);CHKERRQ(ierr); 745f9fae5adSBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 746f9fae5adSBarry Smith for (i=0; i<m; i++) { 747f9fae5adSBarry Smith idx = a->j + a->i[i] + shift; 748f9fae5adSBarry Smith v = a->a + a->i[i] + shift; 749f9fae5adSBarry Smith n = a->i[i+1] - a->i[i]; 750f9fae5adSBarry Smith alpha1 = x[5*i]; 751f9fae5adSBarry Smith alpha2 = x[5*i+1]; 752f9fae5adSBarry Smith alpha3 = x[5*i+2]; 753f9fae5adSBarry Smith alpha4 = x[5*i+3]; 754f9fae5adSBarry Smith alpha5 = x[5*i+4]; 755f9fae5adSBarry Smith while (n-->0) { 756f9fae5adSBarry Smith y[5*(*idx)] += alpha1*(*v); 757f9fae5adSBarry Smith y[5*(*idx)+1] += alpha2*(*v); 758f9fae5adSBarry Smith y[5*(*idx)+2] += alpha3*(*v); 759f9fae5adSBarry Smith y[5*(*idx)+3] += alpha4*(*v); 760f9fae5adSBarry Smith y[5*(*idx)+4] += alpha5*(*v); 761f9fae5adSBarry Smith idx++; v++; 762f9fae5adSBarry Smith } 763f9fae5adSBarry Smith } 764*b0a32e0cSBarry Smith PetscLogFlops(10*a->nz); 765f9fae5adSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 766f9fae5adSBarry Smith ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr); 767f9fae5adSBarry Smith PetscFunctionReturn(0); 768bcc973b7SBarry Smith } 769bcc973b7SBarry Smith 7706cd98798SBarry Smith /* ------------------------------------------------------------------------------*/ 7716cd98798SBarry Smith #undef __FUNC__ 7726cd98798SBarry Smith #define __FUNC__ /*<a name="MatMult_SeqMAIJ_6"></a>*/"MatMult_SeqMAIJ_6" 7736cd98798SBarry Smith int MatMult_SeqMAIJ_6(Mat A,Vec xx,Vec yy) 7746cd98798SBarry Smith { 7756cd98798SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 7766cd98798SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 7776cd98798SBarry Smith Scalar *x,*y,*v,sum1, sum2, sum3, sum4, sum5, sum6; 7786cd98798SBarry Smith int ierr,m = b->AIJ->m,*idx,shift = a->indexshift,*ii; 7796cd98798SBarry Smith int n,i,jrow,j; 7806cd98798SBarry Smith 7816cd98798SBarry Smith PetscFunctionBegin; 7826cd98798SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 7836cd98798SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 7846cd98798SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 7856cd98798SBarry Smith idx = a->j; 7866cd98798SBarry Smith v = a->a; 7876cd98798SBarry Smith ii = a->i; 7886cd98798SBarry Smith 7896cd98798SBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 7906cd98798SBarry Smith idx += shift; 7916cd98798SBarry Smith for (i=0; i<m; i++) { 7926cd98798SBarry Smith jrow = ii[i]; 7936cd98798SBarry Smith n = ii[i+1] - jrow; 7946cd98798SBarry Smith sum1 = 0.0; 7956cd98798SBarry Smith sum2 = 0.0; 7966cd98798SBarry Smith sum3 = 0.0; 7976cd98798SBarry Smith sum4 = 0.0; 7986cd98798SBarry Smith sum5 = 0.0; 7996cd98798SBarry Smith sum6 = 0.0; 8006cd98798SBarry Smith for (j=0; j<n; j++) { 8016cd98798SBarry Smith sum1 += v[jrow]*x[6*idx[jrow]]; 8026cd98798SBarry Smith sum2 += v[jrow]*x[6*idx[jrow]+1]; 8036cd98798SBarry Smith sum3 += v[jrow]*x[6*idx[jrow]+2]; 8046cd98798SBarry Smith sum4 += v[jrow]*x[6*idx[jrow]+3]; 8056cd98798SBarry Smith sum5 += v[jrow]*x[6*idx[jrow]+4]; 8066cd98798SBarry Smith sum6 += v[jrow]*x[6*idx[jrow]+5]; 8076cd98798SBarry Smith jrow++; 8086cd98798SBarry Smith } 8096cd98798SBarry Smith y[6*i] = sum1; 8106cd98798SBarry Smith y[6*i+1] = sum2; 8116cd98798SBarry Smith y[6*i+2] = sum3; 8126cd98798SBarry Smith y[6*i+3] = sum4; 8136cd98798SBarry Smith y[6*i+4] = sum5; 8146cd98798SBarry Smith y[6*i+5] = sum6; 8156cd98798SBarry Smith } 8166cd98798SBarry Smith 817*b0a32e0cSBarry Smith PetscLogFlops(12*a->nz - 6*m); 8186cd98798SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 8196cd98798SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 8206cd98798SBarry Smith PetscFunctionReturn(0); 8216cd98798SBarry Smith } 8226cd98798SBarry Smith 8236cd98798SBarry Smith #undef __FUNC__ 8246cd98798SBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_SeqMAIJ_6"></a>*/"MatMultTranspose_SeqMAIJ_6" 8256cd98798SBarry Smith int MatMultTranspose_SeqMAIJ_6(Mat A,Vec xx,Vec yy) 8266cd98798SBarry Smith { 8276cd98798SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 8286cd98798SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 8296cd98798SBarry Smith Scalar *x,*y,*v,alpha1,alpha2,alpha3,alpha4,alpha5,alpha6,zero = 0.0; 8306cd98798SBarry Smith int ierr,m = b->AIJ->m,n,i,*idx,shift = a->indexshift; 8316cd98798SBarry Smith 8326cd98798SBarry Smith PetscFunctionBegin; 8336cd98798SBarry Smith ierr = VecSet(&zero,yy);CHKERRQ(ierr); 8346cd98798SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 8356cd98798SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 8366cd98798SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 8376cd98798SBarry Smith for (i=0; i<m; i++) { 8386cd98798SBarry Smith idx = a->j + a->i[i] + shift; 8396cd98798SBarry Smith v = a->a + a->i[i] + shift; 8406cd98798SBarry Smith n = a->i[i+1] - a->i[i]; 8416cd98798SBarry Smith alpha1 = x[6*i]; 8426cd98798SBarry Smith alpha2 = x[6*i+1]; 8436cd98798SBarry Smith alpha3 = x[6*i+2]; 8446cd98798SBarry Smith alpha4 = x[6*i+3]; 8456cd98798SBarry Smith alpha5 = x[6*i+4]; 8466cd98798SBarry Smith alpha6 = x[6*i+5]; 8476cd98798SBarry Smith while (n-->0) { 8486cd98798SBarry Smith y[6*(*idx)] += alpha1*(*v); 8496cd98798SBarry Smith y[6*(*idx)+1] += alpha2*(*v); 8506cd98798SBarry Smith y[6*(*idx)+2] += alpha3*(*v); 8516cd98798SBarry Smith y[6*(*idx)+3] += alpha4*(*v); 8526cd98798SBarry Smith y[6*(*idx)+4] += alpha5*(*v); 8536cd98798SBarry Smith y[6*(*idx)+5] += alpha6*(*v); 8546cd98798SBarry Smith idx++; v++; 8556cd98798SBarry Smith } 8566cd98798SBarry Smith } 857*b0a32e0cSBarry Smith PetscLogFlops(12*a->nz - 6*b->AIJ->n); 8586cd98798SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 8596cd98798SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 8606cd98798SBarry Smith PetscFunctionReturn(0); 8616cd98798SBarry Smith } 8626cd98798SBarry Smith 8636cd98798SBarry Smith #undef __FUNC__ 8646cd98798SBarry Smith #define __FUNC__ /*<a name="MatMultAdd_SeqMAIJ_6"></a>*/"MatMultAdd_SeqMAIJ_6" 8656cd98798SBarry Smith int MatMultAdd_SeqMAIJ_6(Mat A,Vec xx,Vec yy,Vec zz) 8666cd98798SBarry Smith { 8676cd98798SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 8686cd98798SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 8696cd98798SBarry Smith Scalar *x,*y,*v,sum1, sum2, sum3, sum4, sum5, sum6; 8706cd98798SBarry Smith int ierr,m = b->AIJ->m,*idx,shift = a->indexshift,*ii; 8716cd98798SBarry Smith int n,i,jrow,j; 8726cd98798SBarry Smith 8736cd98798SBarry Smith PetscFunctionBegin; 8746cd98798SBarry Smith if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);} 8756cd98798SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 8766cd98798SBarry Smith ierr = VecGetArray(zz,&y);CHKERRQ(ierr); 8776cd98798SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 8786cd98798SBarry Smith idx = a->j; 8796cd98798SBarry Smith v = a->a; 8806cd98798SBarry Smith ii = a->i; 8816cd98798SBarry Smith 8826cd98798SBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 8836cd98798SBarry Smith idx += shift; 8846cd98798SBarry Smith for (i=0; i<m; i++) { 8856cd98798SBarry Smith jrow = ii[i]; 8866cd98798SBarry Smith n = ii[i+1] - jrow; 8876cd98798SBarry Smith sum1 = 0.0; 8886cd98798SBarry Smith sum2 = 0.0; 8896cd98798SBarry Smith sum3 = 0.0; 8906cd98798SBarry Smith sum4 = 0.0; 8916cd98798SBarry Smith sum5 = 0.0; 8926cd98798SBarry Smith sum6 = 0.0; 8936cd98798SBarry Smith for (j=0; j<n; j++) { 8946cd98798SBarry Smith sum1 += v[jrow]*x[6*idx[jrow]]; 8956cd98798SBarry Smith sum2 += v[jrow]*x[6*idx[jrow]+1]; 8966cd98798SBarry Smith sum3 += v[jrow]*x[6*idx[jrow]+2]; 8976cd98798SBarry Smith sum4 += v[jrow]*x[6*idx[jrow]+3]; 8986cd98798SBarry Smith sum5 += v[jrow]*x[6*idx[jrow]+4]; 8996cd98798SBarry Smith sum6 += v[jrow]*x[6*idx[jrow]+5]; 9006cd98798SBarry Smith jrow++; 9016cd98798SBarry Smith } 9026cd98798SBarry Smith y[6*i] += sum1; 9036cd98798SBarry Smith y[6*i+1] += sum2; 9046cd98798SBarry Smith y[6*i+2] += sum3; 9056cd98798SBarry Smith y[6*i+3] += sum4; 9066cd98798SBarry Smith y[6*i+4] += sum5; 9076cd98798SBarry Smith y[6*i+5] += sum6; 9086cd98798SBarry Smith } 9096cd98798SBarry Smith 910*b0a32e0cSBarry Smith PetscLogFlops(12*a->nz); 9116cd98798SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 9126cd98798SBarry Smith ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr); 9136cd98798SBarry Smith PetscFunctionReturn(0); 9146cd98798SBarry Smith } 9156cd98798SBarry Smith 9166cd98798SBarry Smith #undef __FUNC__ 9176cd98798SBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_SeqMAIJ_6"></a>*/"MatMultTransposeAdd_SeqMAIJ_6" 9186cd98798SBarry Smith int MatMultTransposeAdd_SeqMAIJ_6(Mat A,Vec xx,Vec yy,Vec zz) 9196cd98798SBarry Smith { 9206cd98798SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 9216cd98798SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 9226cd98798SBarry Smith Scalar *x,*y,*v,alpha1,alpha2,alpha3,alpha4,alpha5,alpha6; 9236cd98798SBarry Smith int ierr,m = b->AIJ->m,n,i,*idx,shift = a->indexshift; 9246cd98798SBarry Smith 9256cd98798SBarry Smith PetscFunctionBegin; 9266cd98798SBarry Smith if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);} 9276cd98798SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 9286cd98798SBarry Smith ierr = VecGetArray(zz,&y);CHKERRQ(ierr); 9296cd98798SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 9306cd98798SBarry Smith for (i=0; i<m; i++) { 9316cd98798SBarry Smith idx = a->j + a->i[i] + shift; 9326cd98798SBarry Smith v = a->a + a->i[i] + shift; 9336cd98798SBarry Smith n = a->i[i+1] - a->i[i]; 9346cd98798SBarry Smith alpha1 = x[6*i]; 9356cd98798SBarry Smith alpha2 = x[6*i+1]; 9366cd98798SBarry Smith alpha3 = x[6*i+2]; 9376cd98798SBarry Smith alpha4 = x[6*i+3]; 9386cd98798SBarry Smith alpha5 = x[6*i+4]; 9396cd98798SBarry Smith alpha6 = x[6*i+5]; 9406cd98798SBarry Smith while (n-->0) { 9416cd98798SBarry Smith y[6*(*idx)] += alpha1*(*v); 9426cd98798SBarry Smith y[6*(*idx)+1] += alpha2*(*v); 9436cd98798SBarry Smith y[6*(*idx)+2] += alpha3*(*v); 9446cd98798SBarry Smith y[6*(*idx)+3] += alpha4*(*v); 9456cd98798SBarry Smith y[6*(*idx)+4] += alpha5*(*v); 9466cd98798SBarry Smith y[6*(*idx)+5] += alpha6*(*v); 9476cd98798SBarry Smith idx++; v++; 9486cd98798SBarry Smith } 9496cd98798SBarry Smith } 950*b0a32e0cSBarry Smith PetscLogFlops(12*a->nz); 9516cd98798SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 9526cd98798SBarry Smith ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr); 9536cd98798SBarry Smith PetscFunctionReturn(0); 9546cd98798SBarry Smith } 9556cd98798SBarry Smith 956f4a53059SBarry Smith /*===================================================================================*/ 957f4a53059SBarry Smith #undef __FUNC__ 958*b0a32e0cSBarry Smith #define __FUNC__ "MatMult_MPIMAIJ_dof" 959f4a53059SBarry Smith int MatMult_MPIMAIJ_dof(Mat A,Vec xx,Vec yy) 960f4a53059SBarry Smith { 9614cb9d9b8SBarry Smith Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data; 962f4a53059SBarry Smith int ierr; 963f4a53059SBarry Smith PetscFunctionBegin; 964f4a53059SBarry Smith 965f4a53059SBarry Smith /* start the scatter */ 966f4a53059SBarry Smith ierr = VecScatterBegin(xx,b->w,INSERT_VALUES,SCATTER_FORWARD,b->ctx);CHKERRQ(ierr); 9674cb9d9b8SBarry Smith ierr = (*b->AIJ->ops->mult)(b->AIJ,xx,yy);CHKERRQ(ierr); 968f4a53059SBarry Smith ierr = VecScatterEnd(xx,b->w,INSERT_VALUES,SCATTER_FORWARD,b->ctx);CHKERRQ(ierr); 9694cb9d9b8SBarry Smith ierr = (*b->OAIJ->ops->multadd)(b->OAIJ,b->w,yy,yy);CHKERRQ(ierr); 970f4a53059SBarry Smith PetscFunctionReturn(0); 971f4a53059SBarry Smith } 972f4a53059SBarry Smith 9734cb9d9b8SBarry Smith #undef __FUNC__ 974*b0a32e0cSBarry Smith #define __FUNC__ "MatMultTranspose_MPIMAIJ_dof" 9754cb9d9b8SBarry Smith int MatMultTranspose_MPIMAIJ_dof(Mat A,Vec xx,Vec yy) 9764cb9d9b8SBarry Smith { 9774cb9d9b8SBarry Smith Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data; 9784cb9d9b8SBarry Smith int ierr; 9794cb9d9b8SBarry Smith PetscFunctionBegin; 9804cb9d9b8SBarry Smith ierr = (*b->OAIJ->ops->multtranspose)(b->OAIJ,xx,b->w);CHKERRQ(ierr); 9814cb9d9b8SBarry Smith ierr = VecScatterBegin(b->w,yy,ADD_VALUES,SCATTER_REVERSE,b->ctx);CHKERRQ(ierr); 9824cb9d9b8SBarry Smith ierr = (*b->AIJ->ops->multtranspose)(b->AIJ,xx,yy);CHKERRQ(ierr); 9834cb9d9b8SBarry Smith ierr = VecScatterEnd(b->w,yy,ADD_VALUES,SCATTER_REVERSE,b->ctx);CHKERRQ(ierr); 9844cb9d9b8SBarry Smith PetscFunctionReturn(0); 9854cb9d9b8SBarry Smith } 9864cb9d9b8SBarry Smith 9874cb9d9b8SBarry Smith #undef __FUNC__ 988*b0a32e0cSBarry Smith #define __FUNC__ "MatMultAdd_MPIMAIJ_dof" 989d72c5c08SBarry Smith int MatMultAdd_MPIMAIJ_dof(Mat A,Vec xx,Vec yy,Vec zz) 9904cb9d9b8SBarry Smith { 9914cb9d9b8SBarry Smith Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data; 9924cb9d9b8SBarry Smith int ierr; 9934cb9d9b8SBarry Smith PetscFunctionBegin; 9944cb9d9b8SBarry Smith 9954cb9d9b8SBarry Smith /* start the scatter */ 9964cb9d9b8SBarry Smith ierr = VecScatterBegin(xx,b->w,INSERT_VALUES,SCATTER_FORWARD,b->ctx);CHKERRQ(ierr); 997d72c5c08SBarry Smith ierr = (*b->AIJ->ops->multadd)(b->AIJ,xx,yy,zz);CHKERRQ(ierr); 9984cb9d9b8SBarry Smith ierr = VecScatterEnd(xx,b->w,INSERT_VALUES,SCATTER_FORWARD,b->ctx);CHKERRQ(ierr); 999d72c5c08SBarry Smith ierr = (*b->OAIJ->ops->multadd)(b->OAIJ,b->w,yy,zz);CHKERRQ(ierr); 10004cb9d9b8SBarry Smith PetscFunctionReturn(0); 10014cb9d9b8SBarry Smith } 10024cb9d9b8SBarry Smith 10034cb9d9b8SBarry Smith #undef __FUNC__ 1004*b0a32e0cSBarry Smith #define __FUNC__ "MatMultTransposeAdd_MPIMAIJ_dof" 1005d72c5c08SBarry Smith int MatMultTransposeAdd_MPIMAIJ_dof(Mat A,Vec xx,Vec yy,Vec zz) 10064cb9d9b8SBarry Smith { 10074cb9d9b8SBarry Smith Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data; 10084cb9d9b8SBarry Smith int ierr; 10094cb9d9b8SBarry Smith PetscFunctionBegin; 10104cb9d9b8SBarry Smith ierr = (*b->OAIJ->ops->multtranspose)(b->OAIJ,xx,b->w);CHKERRQ(ierr); 1011d72c5c08SBarry Smith ierr = VecScatterBegin(b->w,zz,ADD_VALUES,SCATTER_REVERSE,b->ctx);CHKERRQ(ierr); 1012d72c5c08SBarry Smith ierr = (*b->AIJ->ops->multtransposeadd)(b->AIJ,xx,yy,zz);CHKERRQ(ierr); 1013d72c5c08SBarry Smith ierr = VecScatterEnd(b->w,zz,ADD_VALUES,SCATTER_REVERSE,b->ctx);CHKERRQ(ierr); 10144cb9d9b8SBarry Smith PetscFunctionReturn(0); 10154cb9d9b8SBarry Smith } 10164cb9d9b8SBarry Smith 1017bcc973b7SBarry Smith /* ---------------------------------------------------------------------------------- */ 101882b1193eSBarry Smith #undef __FUNC__ 1019*b0a32e0cSBarry Smith #define __FUNC__ "MatCreateMAIJ" 1020cd3bbe55SBarry Smith int MatCreateMAIJ(Mat A,int dof,Mat *maij) 102182b1193eSBarry Smith { 1022f4a53059SBarry Smith int ierr,size,n; 10234cb9d9b8SBarry Smith Mat_MPIMAIJ *b; 102482b1193eSBarry Smith Mat B; 102582b1193eSBarry Smith 102682b1193eSBarry Smith PetscFunctionBegin; 1027d72c5c08SBarry Smith ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 1028d72c5c08SBarry Smith 1029d72c5c08SBarry Smith if (dof == 1) { 1030d72c5c08SBarry Smith *maij = A; 1031d72c5c08SBarry Smith } else { 103283903688SBarry Smith ierr = MatCreate(A->comm,dof*A->m,dof*A->n,dof*A->M,dof*A->N,&B);CHKERRQ(ierr); 1033cd3bbe55SBarry Smith B->assembled = PETSC_TRUE; 1034d72c5c08SBarry Smith 1035f4a53059SBarry Smith ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 1036f4a53059SBarry Smith if (size == 1) { 1037b9b97703SBarry Smith ierr = MatSetType(B,MATSEQMAIJ);CHKERRQ(ierr); 10384cb9d9b8SBarry Smith B->ops->destroy = MatDestroy_SeqMAIJ; 1039b9b97703SBarry Smith b = (Mat_MPIMAIJ*)B->data; 1040b9b97703SBarry Smith b->dof = dof; 10414cb9d9b8SBarry Smith b->AIJ = A; 1042d72c5c08SBarry Smith if (dof == 2) { 1043cd3bbe55SBarry Smith B->ops->mult = MatMult_SeqMAIJ_2; 1044cd3bbe55SBarry Smith B->ops->multadd = MatMultAdd_SeqMAIJ_2; 1045cd3bbe55SBarry Smith B->ops->multtranspose = MatMultTranspose_SeqMAIJ_2; 1046cd3bbe55SBarry Smith B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_2; 1047bcc973b7SBarry Smith } else if (dof == 3) { 1048bcc973b7SBarry Smith B->ops->mult = MatMult_SeqMAIJ_3; 1049bcc973b7SBarry Smith B->ops->multadd = MatMultAdd_SeqMAIJ_3; 1050bcc973b7SBarry Smith B->ops->multtranspose = MatMultTranspose_SeqMAIJ_3; 1051bcc973b7SBarry Smith B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_3; 1052bcc973b7SBarry Smith } else if (dof == 4) { 1053bcc973b7SBarry Smith B->ops->mult = MatMult_SeqMAIJ_4; 1054bcc973b7SBarry Smith B->ops->multadd = MatMultAdd_SeqMAIJ_4; 1055bcc973b7SBarry Smith B->ops->multtranspose = MatMultTranspose_SeqMAIJ_4; 1056bcc973b7SBarry Smith B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_4; 1057f9fae5adSBarry Smith } else if (dof == 5) { 1058f9fae5adSBarry Smith B->ops->mult = MatMult_SeqMAIJ_5; 1059f9fae5adSBarry Smith B->ops->multadd = MatMultAdd_SeqMAIJ_5; 1060f9fae5adSBarry Smith B->ops->multtranspose = MatMultTranspose_SeqMAIJ_5; 1061f9fae5adSBarry Smith B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_5; 10626cd98798SBarry Smith } else if (dof == 6) { 10636cd98798SBarry Smith B->ops->mult = MatMult_SeqMAIJ_6; 10646cd98798SBarry Smith B->ops->multadd = MatMultAdd_SeqMAIJ_6; 10656cd98798SBarry Smith B->ops->multtranspose = MatMultTranspose_SeqMAIJ_6; 10666cd98798SBarry Smith B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_6; 106782b1193eSBarry Smith } else { 106829bbc08cSBarry Smith SETERRQ1(1,"Cannot handle a dof of %d\n",dof); 106982b1193eSBarry Smith } 1070f4a53059SBarry Smith } else { 1071f4a53059SBarry Smith Mat_MPIAIJ *mpiaij = (Mat_MPIAIJ *)A->data; 1072f4a53059SBarry Smith IS from,to; 1073f4a53059SBarry Smith Vec gvec; 1074f4a53059SBarry Smith int *garray,i; 1075f4a53059SBarry Smith 1076b9b97703SBarry Smith ierr = MatSetType(B,MATMPIMAIJ);CHKERRQ(ierr); 1077d72c5c08SBarry Smith B->ops->destroy = MatDestroy_MPIMAIJ; 1078b9b97703SBarry Smith b = (Mat_MPIMAIJ*)B->data; 1079b9b97703SBarry Smith b->dof = dof; 1080b9b97703SBarry Smith b->A = A; 10814cb9d9b8SBarry Smith ierr = MatCreateMAIJ(mpiaij->A,dof,&b->AIJ);CHKERRQ(ierr); 10824cb9d9b8SBarry Smith ierr = MatCreateMAIJ(mpiaij->B,dof,&b->OAIJ);CHKERRQ(ierr); 10834cb9d9b8SBarry Smith 1084f4a53059SBarry Smith ierr = VecGetSize(mpiaij->lvec,&n);CHKERRQ(ierr); 1085f4a53059SBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,n*dof,&b->w);CHKERRQ(ierr); 1086f4a53059SBarry Smith 1087f4a53059SBarry Smith /* create two temporary Index sets for build scatter gather */ 1088*b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),& garray );CHKERRQ(ierr); 1089f4a53059SBarry Smith for (i=0; i<n; i++) garray[i] = dof*mpiaij->garray[i]; 1090f4a53059SBarry Smith ierr = ISCreateBlock(A->comm,dof,n,garray,&from);CHKERRQ(ierr); 1091f4a53059SBarry Smith ierr = PetscFree(garray);CHKERRQ(ierr); 1092f4a53059SBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,n*dof,0,1,&to);CHKERRQ(ierr); 1093f4a53059SBarry Smith 1094f4a53059SBarry Smith /* create temporary global vector to generate scatter context */ 1095f4a53059SBarry Smith ierr = VecCreateMPI(A->comm,dof*A->n,dof*A->N,&gvec);CHKERRQ(ierr); 1096f4a53059SBarry Smith 1097f4a53059SBarry Smith /* generate the scatter context */ 1098f4a53059SBarry Smith ierr = VecScatterCreate(gvec,from,b->w,to,&b->ctx);CHKERRQ(ierr); 1099f4a53059SBarry Smith 1100f4a53059SBarry Smith ierr = ISDestroy(from);CHKERRQ(ierr); 1101f4a53059SBarry Smith ierr = ISDestroy(to);CHKERRQ(ierr); 1102f4a53059SBarry Smith ierr = VecDestroy(gvec);CHKERRQ(ierr); 1103f4a53059SBarry Smith 1104f4a53059SBarry Smith B->ops->mult = MatMult_MPIMAIJ_dof; 11054cb9d9b8SBarry Smith B->ops->multtranspose = MatMultTranspose_MPIMAIJ_dof; 11064cb9d9b8SBarry Smith B->ops->multadd = MatMultAdd_MPIMAIJ_dof; 11074cb9d9b8SBarry Smith B->ops->multtransposeadd = MatMultTransposeAdd_MPIMAIJ_dof; 1108f4a53059SBarry Smith } 1109cd3bbe55SBarry Smith *maij = B; 1110d72c5c08SBarry Smith } 111182b1193eSBarry Smith PetscFunctionReturn(0); 111282b1193eSBarry Smith } 111382b1193eSBarry Smith 111482b1193eSBarry Smith 111582b1193eSBarry Smith 111682b1193eSBarry Smith 111782b1193eSBarry Smith 111882b1193eSBarry Smith 111982b1193eSBarry Smith 112082b1193eSBarry Smith 112182b1193eSBarry Smith 112282b1193eSBarry Smith 112382b1193eSBarry Smith 112482b1193eSBarry Smith 1125