1*bcc973b7SBarry Smith /*$Id: maij.c,v 1.2 2000/05/27 03:52:46 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*) 1582b1193eSBarry Smith */ 1682b1193eSBarry Smith 1782b1193eSBarry Smith #include "src/mat/impls/aij/seq/aij.h" 1882b1193eSBarry Smith 19cd3bbe55SBarry Smith typedef struct { 20cd3bbe55SBarry Smith int dof; /* number of components */ 21cd3bbe55SBarry Smith Mat AIJ; /* representation of interpolation for one component */ 22cd3bbe55SBarry Smith } Mat_SeqMAIJ; 2382b1193eSBarry Smith 2482b1193eSBarry Smith #undef __FUNC__ 25cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatDestroy_SeqMAIJ"></a>*/"MatDestroy_SeqMAIJ" 26cd3bbe55SBarry Smith int MatDestroy_SeqMAIJ(Mat A) 2782b1193eSBarry Smith { 2882b1193eSBarry Smith int ierr; 29cd3bbe55SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 3082b1193eSBarry Smith 3182b1193eSBarry Smith PetscFunctionBegin; 32cd3bbe55SBarry Smith if (b->AIJ) { 33cd3bbe55SBarry Smith ierr = MatDestroy(b->AIJ);CHKERRQ(ierr); 3482b1193eSBarry Smith } 35cd3bbe55SBarry Smith ierr = PetscFree(b);CHKERRQ(ierr); 3682b1193eSBarry Smith PetscHeaderDestroy(A); 3782b1193eSBarry Smith PetscFunctionReturn(0); 3882b1193eSBarry Smith } 3982b1193eSBarry Smith 4082b1193eSBarry Smith EXTERN_C_BEGIN 4182b1193eSBarry Smith #undef __FUNC__ 42cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatCreate_SeqMAIJ"></a>*/"MatCreate_SeqMAIJ" 43cd3bbe55SBarry Smith int MatCreate_SeqMAIJ(Mat A) 4482b1193eSBarry Smith { 45cd3bbe55SBarry Smith int ierr; 46cd3bbe55SBarry Smith Mat_SeqMAIJ *b; 4782b1193eSBarry Smith 4882b1193eSBarry Smith PetscFunctionBegin; 49cd3bbe55SBarry Smith A->data = (void*)(b = PetscNew(Mat_SeqMAIJ));CHKPTRQ(b); 50cd3bbe55SBarry Smith ierr = PetscMemzero(b,sizeof(Mat_SeqMAIJ));CHKERRQ(ierr); 51cd3bbe55SBarry Smith ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 52cd3bbe55SBarry Smith A->factor = 0; 53cd3bbe55SBarry Smith A->mapping = 0; 54cd3bbe55SBarry Smith b->AIJ = 0; 55cd3bbe55SBarry Smith b->dof = 0; 5682b1193eSBarry Smith PetscFunctionReturn(0); 5782b1193eSBarry Smith } 5882b1193eSBarry Smith EXTERN_C_END 5982b1193eSBarry Smith 60*bcc973b7SBarry Smith /* ----------------------------------------------------------------------------------*/ 61cd3bbe55SBarry Smith EXTERN int MatMult_SeqAIJ(Mat,Vec,Vec); 62cd3bbe55SBarry Smith EXTERN int MatMultTranspose_SeqAIJ(Mat,Vec,Vec); 63cd3bbe55SBarry Smith EXTERN int MatMultTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec); 64cd3bbe55SBarry Smith EXTERN int matmulttransposeadd_seqaijMatMultAdd_SeqAIJ(Mat,Vec,Vec,Vec); 65cd3bbe55SBarry Smith 6682b1193eSBarry Smith #undef __FUNC__ 67cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMult_SeqMAIJ_1"></a>*/"MatMult_SeqMAIJ_1" 68cd3bbe55SBarry Smith int MatMult_SeqMAIJ_1(Mat A,Vec xx,Vec yy) 6982b1193eSBarry Smith { 70cd3bbe55SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 71cd3bbe55SBarry Smith int ierr; 7282b1193eSBarry Smith PetscFunctionBegin; 73cd3bbe55SBarry Smith ierr = MatMult_SeqAIJ(b->AIJ,xx,yy); 74cd3bbe55SBarry Smith PetscFunctionReturn(0); 7582b1193eSBarry Smith } 76cd3bbe55SBarry Smith #undef __FUNC__ 77cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_SeqMAIJ_1"></a>*/"MatMultTranspose_SeqMAIJ_1" 78cd3bbe55SBarry Smith int MatMultTranspose_SeqMAIJ_1(Mat A,Vec xx,Vec yy) 79cd3bbe55SBarry Smith { 80cd3bbe55SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 81cd3bbe55SBarry Smith int ierr; 82cd3bbe55SBarry Smith PetscFunctionBegin; 83cd3bbe55SBarry Smith ierr = MatMultTranspose_SeqAIJ(b->AIJ,xx,yy); 84cd3bbe55SBarry Smith PetscFunctionReturn(0); 85cd3bbe55SBarry Smith } 86cd3bbe55SBarry Smith #undef __FUNC__ 87cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMultAdd_SeqMAIJ_1"></a>*/"MatMultAdd_SeqMAIJ_1" 88cd3bbe55SBarry Smith int MatMultAdd_SeqMAIJ_1(Mat A,Vec xx,Vec yy,Vec zz) 89cd3bbe55SBarry Smith { 90cd3bbe55SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 91cd3bbe55SBarry Smith int ierr; 92cd3bbe55SBarry Smith PetscFunctionBegin; 93cd3bbe55SBarry Smith ierr = MatMultAdd_SeqAIJ(b->AIJ,xx,yy,zz); 94cd3bbe55SBarry Smith PetscFunctionReturn(0); 95cd3bbe55SBarry Smith } 96cd3bbe55SBarry Smith #undef __FUNC__ 97cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_SeqMAIJ_1"></a>*/"MatMultTransposeAdd_SeqMAIJ_1" 98cd3bbe55SBarry Smith int MatMultTransposeAdd_SeqMAIJ_1(Mat A,Vec xx,Vec yy,Vec zz) 99cd3bbe55SBarry Smith { 100cd3bbe55SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 101cd3bbe55SBarry Smith int ierr; 102cd3bbe55SBarry Smith PetscFunctionBegin; 103cd3bbe55SBarry Smith ierr = MatMultTransposeAdd_SeqAIJ(b->AIJ,xx,yy,zz); 10482b1193eSBarry Smith PetscFunctionReturn(0); 10582b1193eSBarry Smith } 10682b1193eSBarry Smith 107cd3bbe55SBarry Smith /* --------------------------------------------------------------------------------------*/ 10882b1193eSBarry Smith #undef __FUNC__ 109cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMult_SeqMAIJ_2"></a>*/"MatMult_SeqMAIJ_2" 110cd3bbe55SBarry Smith int MatMult_SeqMAIJ_2(Mat A,Vec xx,Vec yy) 11182b1193eSBarry Smith { 112cd3bbe55SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 113*bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 114*bcc973b7SBarry Smith Scalar *x,*y,*v,sum1, sum2; 115*bcc973b7SBarry Smith int ierr,m = a->m,*idx,shift = a->indexshift,*ii; 116*bcc973b7SBarry Smith int n,i,jrow,j; 11782b1193eSBarry Smith 118*bcc973b7SBarry Smith PetscFunctionBegin; 119*bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 120*bcc973b7SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 121*bcc973b7SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 122*bcc973b7SBarry Smith idx = a->j; 123*bcc973b7SBarry Smith v = a->a; 124*bcc973b7SBarry Smith ii = a->i; 125*bcc973b7SBarry Smith 126*bcc973b7SBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 127*bcc973b7SBarry Smith idx += shift; 128*bcc973b7SBarry Smith for (i=0; i<m; i++) { 129*bcc973b7SBarry Smith jrow = ii[i]; 130*bcc973b7SBarry Smith n = ii[i+1] - jrow; 131*bcc973b7SBarry Smith sum1 = 0.0; 132*bcc973b7SBarry Smith sum2 = 0.0; 133*bcc973b7SBarry Smith for (j=0; j<n; j++) { 134*bcc973b7SBarry Smith sum1 += v[jrow]*x[2*idx[jrow]]; 135*bcc973b7SBarry Smith sum2 += v[jrow]*x[2*idx[jrow]+1]; 136*bcc973b7SBarry Smith jrow++; 137*bcc973b7SBarry Smith } 138*bcc973b7SBarry Smith y[2*i] = sum1; 139*bcc973b7SBarry Smith y[2*i+1] = sum2; 140*bcc973b7SBarry Smith } 141*bcc973b7SBarry Smith 142*bcc973b7SBarry Smith PLogFlops(4*a->nz - 2*m); 143*bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 144*bcc973b7SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 14582b1193eSBarry Smith PetscFunctionReturn(0); 14682b1193eSBarry Smith } 147*bcc973b7SBarry Smith 14882b1193eSBarry Smith #undef __FUNC__ 149cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_SeqMAIJ_2"></a>*/"MatMultTranspose_SeqMAIJ_2" 150cd3bbe55SBarry Smith int MatMultTranspose_SeqMAIJ_2(Mat A,Vec xx,Vec yy) 15182b1193eSBarry Smith { 152cd3bbe55SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 153*bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 154*bcc973b7SBarry Smith Scalar *x,*y,*v,alpha1,alpha2,zero = 0.0; 155*bcc973b7SBarry Smith int ierr,m = a->m,n,i,*idx,shift = a->indexshift; 15682b1193eSBarry Smith 157*bcc973b7SBarry Smith PetscFunctionBegin; 158*bcc973b7SBarry Smith ierr = VecSet(&zero,yy);CHKERRQ(ierr); 159*bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 160*bcc973b7SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 161*bcc973b7SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 162*bcc973b7SBarry Smith for (i=0; i<m; i++) { 163*bcc973b7SBarry Smith idx = a->j + a->i[i] + shift; 164*bcc973b7SBarry Smith v = a->a + a->i[i] + shift; 165*bcc973b7SBarry Smith n = a->i[i+1] - a->i[i]; 166*bcc973b7SBarry Smith alpha1 = x[2*i]; 167*bcc973b7SBarry Smith alpha2 = x[2*i+1]; 168*bcc973b7SBarry Smith while (n-->0) {y[2*(*idx)] += alpha1*(*v); y[2*(*idx)+1] += alpha2*(*v); idx++; v++;} 169*bcc973b7SBarry Smith } 170*bcc973b7SBarry Smith PLogFlops(4*a->nz - 2*a->n); 171*bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 172*bcc973b7SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 17382b1193eSBarry Smith PetscFunctionReturn(0); 17482b1193eSBarry Smith } 175*bcc973b7SBarry Smith 17682b1193eSBarry Smith #undef __FUNC__ 177cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMultAdd_SeqMAIJ_2"></a>*/"MatMultAdd_SeqMAIJ_2" 178cd3bbe55SBarry Smith int MatMultAdd_SeqMAIJ_2(Mat A,Vec xx,Vec yy,Vec zz) 17982b1193eSBarry Smith { 180cd3bbe55SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 181*bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 182*bcc973b7SBarry Smith Scalar *x,*y,*v,sum1, sum2; 183*bcc973b7SBarry Smith int ierr,m = a->m,*idx,shift = a->indexshift,*ii; 184*bcc973b7SBarry Smith int n,i,jrow,j; 18582b1193eSBarry Smith 186*bcc973b7SBarry Smith PetscFunctionBegin; 187*bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 188*bcc973b7SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 189*bcc973b7SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 190*bcc973b7SBarry Smith idx = a->j; 191*bcc973b7SBarry Smith v = a->a; 192*bcc973b7SBarry Smith ii = a->i; 193*bcc973b7SBarry Smith 194*bcc973b7SBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 195*bcc973b7SBarry Smith idx += shift; 196*bcc973b7SBarry Smith for (i=0; i<m; i++) { 197*bcc973b7SBarry Smith jrow = ii[i]; 198*bcc973b7SBarry Smith n = ii[i+1] - jrow; 199*bcc973b7SBarry Smith sum1 = 0.0; 200*bcc973b7SBarry Smith sum2 = 0.0; 201*bcc973b7SBarry Smith for (j=0; j<n; j++) { 202*bcc973b7SBarry Smith sum1 += v[jrow]*x[2*idx[jrow]]; 203*bcc973b7SBarry Smith sum2 += v[jrow]*x[2*idx[jrow]+1]; 204*bcc973b7SBarry Smith jrow++; 205*bcc973b7SBarry Smith } 206*bcc973b7SBarry Smith y[2*i] += sum1; 207*bcc973b7SBarry Smith y[2*i+1] += sum2; 208*bcc973b7SBarry Smith } 209*bcc973b7SBarry Smith 210*bcc973b7SBarry Smith PLogFlops(4*a->nz - 2*m); 211*bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 212*bcc973b7SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 213*bcc973b7SBarry Smith PetscFunctionReturn(0); 21482b1193eSBarry Smith PetscFunctionReturn(0); 21582b1193eSBarry Smith } 21682b1193eSBarry Smith #undef __FUNC__ 217cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_SeqMAIJ_2"></a>*/"MatMultTransposeAdd_SeqMAIJ_2" 218cd3bbe55SBarry Smith int MatMultTransposeAdd_SeqMAIJ_2(Mat A,Vec xx,Vec yy,Vec zz) 21982b1193eSBarry Smith { 220cd3bbe55SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 221*bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 222*bcc973b7SBarry Smith Scalar *x,*y,*v,alpha1,alpha2; 223*bcc973b7SBarry Smith int ierr,m = a->m,n,i,*idx,shift = a->indexshift; 22482b1193eSBarry Smith 225*bcc973b7SBarry Smith PetscFunctionBegin; 226*bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 227*bcc973b7SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 228*bcc973b7SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 229*bcc973b7SBarry Smith for (i=0; i<m; i++) { 230*bcc973b7SBarry Smith idx = a->j + a->i[i] + shift; 231*bcc973b7SBarry Smith v = a->a + a->i[i] + shift; 232*bcc973b7SBarry Smith n = a->i[i+1] - a->i[i]; 233*bcc973b7SBarry Smith alpha1 = x[2*i]; 234*bcc973b7SBarry Smith alpha2 = x[2*i+1]; 235*bcc973b7SBarry Smith while (n-->0) {y[2*(*idx)] += alpha1*(*v); y[2*(*idx)+1] += alpha2*(*v); idx++; v++;} 236*bcc973b7SBarry Smith } 237*bcc973b7SBarry Smith PLogFlops(4*a->nz - 2*a->n); 238*bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 239*bcc973b7SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 240*bcc973b7SBarry Smith PetscFunctionReturn(0); 24182b1193eSBarry Smith PetscFunctionReturn(0); 24282b1193eSBarry Smith } 243cd3bbe55SBarry Smith /* --------------------------------------------------------------------------------------*/ 244*bcc973b7SBarry Smith #undef __FUNC__ 245*bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMult_SeqMAIJ_3"></a>*/"MatMult_SeqMAIJ_3" 246*bcc973b7SBarry Smith int MatMult_SeqMAIJ_3(Mat A,Vec xx,Vec yy) 247*bcc973b7SBarry Smith { 248*bcc973b7SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 249*bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 250*bcc973b7SBarry Smith Scalar *x,*y,*v,sum1, sum2, sum3; 251*bcc973b7SBarry Smith int ierr,m = a->m,*idx,shift = a->indexshift,*ii; 252*bcc973b7SBarry Smith int n,i,jrow,j; 25382b1193eSBarry Smith 254*bcc973b7SBarry Smith PetscFunctionBegin; 255*bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 256*bcc973b7SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 257*bcc973b7SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 258*bcc973b7SBarry Smith idx = a->j; 259*bcc973b7SBarry Smith v = a->a; 260*bcc973b7SBarry Smith ii = a->i; 261*bcc973b7SBarry Smith 262*bcc973b7SBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 263*bcc973b7SBarry Smith idx += shift; 264*bcc973b7SBarry Smith for (i=0; i<m; i++) { 265*bcc973b7SBarry Smith jrow = ii[i]; 266*bcc973b7SBarry Smith n = ii[i+1] - jrow; 267*bcc973b7SBarry Smith sum1 = 0.0; 268*bcc973b7SBarry Smith sum2 = 0.0; 269*bcc973b7SBarry Smith sum3 = 0.0; 270*bcc973b7SBarry Smith for (j=0; j<n; j++) { 271*bcc973b7SBarry Smith sum1 += v[jrow]*x[3*idx[jrow]]; 272*bcc973b7SBarry Smith sum2 += v[jrow]*x[3*idx[jrow]+1]; 273*bcc973b7SBarry Smith sum3 += v[jrow]*x[3*idx[jrow]+2]; 274*bcc973b7SBarry Smith jrow++; 275*bcc973b7SBarry Smith } 276*bcc973b7SBarry Smith y[3*i] = sum1; 277*bcc973b7SBarry Smith y[3*i+1] = sum2; 278*bcc973b7SBarry Smith y[3*i+2] = sum3; 279*bcc973b7SBarry Smith } 280*bcc973b7SBarry Smith 281*bcc973b7SBarry Smith PLogFlops(6*a->nz - 3*m); 282*bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 283*bcc973b7SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 284*bcc973b7SBarry Smith PetscFunctionReturn(0); 285*bcc973b7SBarry Smith } 286*bcc973b7SBarry Smith 287*bcc973b7SBarry Smith #undef __FUNC__ 288*bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_SeqMAIJ_3"></a>*/"MatMultTranspose_SeqMAIJ_3" 289*bcc973b7SBarry Smith int MatMultTranspose_SeqMAIJ_3(Mat A,Vec xx,Vec yy) 290*bcc973b7SBarry Smith { 291*bcc973b7SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 292*bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 293*bcc973b7SBarry Smith Scalar *x,*y,*v,alpha1,alpha2,alpha3,zero = 0.0; 294*bcc973b7SBarry Smith int ierr,m = a->m,n,i,*idx,shift = a->indexshift; 295*bcc973b7SBarry Smith 296*bcc973b7SBarry Smith PetscFunctionBegin; 297*bcc973b7SBarry Smith ierr = VecSet(&zero,yy);CHKERRQ(ierr); 298*bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 299*bcc973b7SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 300*bcc973b7SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 301*bcc973b7SBarry Smith for (i=0; i<m; i++) { 302*bcc973b7SBarry Smith idx = a->j + a->i[i] + shift; 303*bcc973b7SBarry Smith v = a->a + a->i[i] + shift; 304*bcc973b7SBarry Smith n = a->i[i+1] - a->i[i]; 305*bcc973b7SBarry Smith alpha1 = x[3*i]; 306*bcc973b7SBarry Smith alpha2 = x[3*i+1]; 307*bcc973b7SBarry Smith alpha3 = x[3*i+2]; 308*bcc973b7SBarry Smith while (n-->0) { 309*bcc973b7SBarry Smith y[3*(*idx)] += alpha1*(*v); 310*bcc973b7SBarry Smith y[3*(*idx)+1] += alpha2*(*v); 311*bcc973b7SBarry Smith y[3*(*idx)+2] += alpha3*(*v); 312*bcc973b7SBarry Smith idx++; v++; 313*bcc973b7SBarry Smith } 314*bcc973b7SBarry Smith } 315*bcc973b7SBarry Smith PLogFlops(6*a->nz - 3*a->n); 316*bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 317*bcc973b7SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 318*bcc973b7SBarry Smith PetscFunctionReturn(0); 319*bcc973b7SBarry Smith } 320*bcc973b7SBarry Smith 321*bcc973b7SBarry Smith #undef __FUNC__ 322*bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultAdd_SeqMAIJ_3"></a>*/"MatMultAdd_SeqMAIJ_3" 323*bcc973b7SBarry Smith int MatMultAdd_SeqMAIJ_3(Mat A,Vec xx,Vec yy,Vec zz) 324*bcc973b7SBarry Smith { 325*bcc973b7SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 326*bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 327*bcc973b7SBarry Smith Scalar *x,*y,*v,sum1, sum2, sum3; 328*bcc973b7SBarry Smith int ierr,m = a->m,*idx,shift = a->indexshift,*ii; 329*bcc973b7SBarry Smith int n,i,jrow,j; 330*bcc973b7SBarry Smith 331*bcc973b7SBarry Smith PetscFunctionBegin; 332*bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 333*bcc973b7SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 334*bcc973b7SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 335*bcc973b7SBarry Smith idx = a->j; 336*bcc973b7SBarry Smith v = a->a; 337*bcc973b7SBarry Smith ii = a->i; 338*bcc973b7SBarry Smith 339*bcc973b7SBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 340*bcc973b7SBarry Smith idx += shift; 341*bcc973b7SBarry Smith for (i=0; i<m; i++) { 342*bcc973b7SBarry Smith jrow = ii[i]; 343*bcc973b7SBarry Smith n = ii[i+1] - jrow; 344*bcc973b7SBarry Smith sum1 = 0.0; 345*bcc973b7SBarry Smith sum2 = 0.0; 346*bcc973b7SBarry Smith sum3 = 0.0; 347*bcc973b7SBarry Smith for (j=0; j<n; j++) { 348*bcc973b7SBarry Smith sum1 += v[jrow]*x[3*idx[jrow]]; 349*bcc973b7SBarry Smith sum2 += v[jrow]*x[3*idx[jrow]+1]; 350*bcc973b7SBarry Smith sum3 += v[jrow]*x[3*idx[jrow]+2]; 351*bcc973b7SBarry Smith jrow++; 352*bcc973b7SBarry Smith } 353*bcc973b7SBarry Smith y[3*i] += sum1; 354*bcc973b7SBarry Smith y[3*i+1] += sum2; 355*bcc973b7SBarry Smith y[3*i+2] += sum3; 356*bcc973b7SBarry Smith } 357*bcc973b7SBarry Smith 358*bcc973b7SBarry Smith PLogFlops(6*a->nz); 359*bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 360*bcc973b7SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 361*bcc973b7SBarry Smith PetscFunctionReturn(0); 362*bcc973b7SBarry Smith } 363*bcc973b7SBarry Smith #undef __FUNC__ 364*bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_SeqMAIJ_3"></a>*/"MatMultTransposeAdd_SeqMAIJ_3" 365*bcc973b7SBarry Smith int MatMultTransposeAdd_SeqMAIJ_3(Mat A,Vec xx,Vec yy,Vec zz) 366*bcc973b7SBarry Smith { 367*bcc973b7SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 368*bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 369*bcc973b7SBarry Smith Scalar *x,*y,*v,alpha1,alpha2,alpha3; 370*bcc973b7SBarry Smith int ierr,m = a->m,n,i,*idx,shift = a->indexshift; 371*bcc973b7SBarry Smith 372*bcc973b7SBarry Smith PetscFunctionBegin; 373*bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 374*bcc973b7SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 375*bcc973b7SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 376*bcc973b7SBarry Smith for (i=0; i<m; i++) { 377*bcc973b7SBarry Smith idx = a->j + a->i[i] + shift; 378*bcc973b7SBarry Smith v = a->a + a->i[i] + shift; 379*bcc973b7SBarry Smith n = a->i[i+1] - a->i[i]; 380*bcc973b7SBarry Smith alpha1 = x[3*i]; 381*bcc973b7SBarry Smith alpha2 = x[3*i+1]; 382*bcc973b7SBarry Smith alpha3 = x[3*i+2]; 383*bcc973b7SBarry Smith while (n-->0) { 384*bcc973b7SBarry Smith y[3*(*idx)] += alpha1*(*v); 385*bcc973b7SBarry Smith y[3*(*idx)+1] += alpha2*(*v); 386*bcc973b7SBarry Smith y[3*(*idx)+2] += alpha3*(*v); 387*bcc973b7SBarry Smith idx++; v++; 388*bcc973b7SBarry Smith } 389*bcc973b7SBarry Smith } 390*bcc973b7SBarry Smith PLogFlops(6*a->nz); 391*bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 392*bcc973b7SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 393*bcc973b7SBarry Smith PetscFunctionReturn(0); 394*bcc973b7SBarry Smith } 395*bcc973b7SBarry Smith 396*bcc973b7SBarry Smith /* ------------------------------------------------------------------------------*/ 397*bcc973b7SBarry Smith #undef __FUNC__ 398*bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMult_SeqMAIJ_4"></a>*/"MatMult_SeqMAIJ_4" 399*bcc973b7SBarry Smith int MatMult_SeqMAIJ_4(Mat A,Vec xx,Vec yy) 400*bcc973b7SBarry Smith { 401*bcc973b7SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 402*bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 403*bcc973b7SBarry Smith Scalar *x,*y,*v,sum1, sum2, sum3, sum4; 404*bcc973b7SBarry Smith int ierr,m = a->m,*idx,shift = a->indexshift,*ii; 405*bcc973b7SBarry Smith int n,i,jrow,j; 406*bcc973b7SBarry Smith 407*bcc973b7SBarry Smith PetscFunctionBegin; 408*bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 409*bcc973b7SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 410*bcc973b7SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 411*bcc973b7SBarry Smith idx = a->j; 412*bcc973b7SBarry Smith v = a->a; 413*bcc973b7SBarry Smith ii = a->i; 414*bcc973b7SBarry Smith 415*bcc973b7SBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 416*bcc973b7SBarry Smith idx += shift; 417*bcc973b7SBarry Smith for (i=0; i<m; i++) { 418*bcc973b7SBarry Smith jrow = ii[i]; 419*bcc973b7SBarry Smith n = ii[i+1] - jrow; 420*bcc973b7SBarry Smith sum1 = 0.0; 421*bcc973b7SBarry Smith sum2 = 0.0; 422*bcc973b7SBarry Smith sum3 = 0.0; 423*bcc973b7SBarry Smith sum4 = 0.0; 424*bcc973b7SBarry Smith for (j=0; j<n; j++) { 425*bcc973b7SBarry Smith sum1 += v[jrow]*x[4*idx[jrow]]; 426*bcc973b7SBarry Smith sum2 += v[jrow]*x[4*idx[jrow]+1]; 427*bcc973b7SBarry Smith sum3 += v[jrow]*x[4*idx[jrow]+2]; 428*bcc973b7SBarry Smith sum4 += v[jrow]*x[4*idx[jrow]+3]; 429*bcc973b7SBarry Smith jrow++; 430*bcc973b7SBarry Smith } 431*bcc973b7SBarry Smith y[4*i] = sum1; 432*bcc973b7SBarry Smith y[4*i+1] = sum2; 433*bcc973b7SBarry Smith y[4*i+2] = sum3; 434*bcc973b7SBarry Smith y[4*i+3] = sum4; 435*bcc973b7SBarry Smith } 436*bcc973b7SBarry Smith 437*bcc973b7SBarry Smith PLogFlops(8*a->nz - 4*m); 438*bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 439*bcc973b7SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 440*bcc973b7SBarry Smith PetscFunctionReturn(0); 441*bcc973b7SBarry Smith } 442*bcc973b7SBarry Smith 443*bcc973b7SBarry Smith #undef __FUNC__ 444*bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_SeqMAIJ_4"></a>*/"MatMultTranspose_SeqMAIJ_4" 445*bcc973b7SBarry Smith int MatMultTranspose_SeqMAIJ_4(Mat A,Vec xx,Vec yy) 446*bcc973b7SBarry Smith { 447*bcc973b7SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 448*bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 449*bcc973b7SBarry Smith Scalar *x,*y,*v,alpha1,alpha2,alpha3,alpha4,zero = 0.0; 450*bcc973b7SBarry Smith int ierr,m = a->m,n,i,*idx,shift = a->indexshift; 451*bcc973b7SBarry Smith 452*bcc973b7SBarry Smith PetscFunctionBegin; 453*bcc973b7SBarry Smith ierr = VecSet(&zero,yy);CHKERRQ(ierr); 454*bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 455*bcc973b7SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 456*bcc973b7SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 457*bcc973b7SBarry Smith for (i=0; i<m; i++) { 458*bcc973b7SBarry Smith idx = a->j + a->i[i] + shift; 459*bcc973b7SBarry Smith v = a->a + a->i[i] + shift; 460*bcc973b7SBarry Smith n = a->i[i+1] - a->i[i]; 461*bcc973b7SBarry Smith alpha1 = x[4*i]; 462*bcc973b7SBarry Smith alpha2 = x[4*i+1]; 463*bcc973b7SBarry Smith alpha3 = x[4*i+2]; 464*bcc973b7SBarry Smith alpha4 = x[4*i+3]; 465*bcc973b7SBarry Smith while (n-->0) { 466*bcc973b7SBarry Smith y[4*(*idx)] += alpha1*(*v); 467*bcc973b7SBarry Smith y[4*(*idx)+1] += alpha2*(*v); 468*bcc973b7SBarry Smith y[4*(*idx)+2] += alpha3*(*v); 469*bcc973b7SBarry Smith y[4*(*idx)+3] += alpha4*(*v); 470*bcc973b7SBarry Smith idx++; v++; 471*bcc973b7SBarry Smith } 472*bcc973b7SBarry Smith } 473*bcc973b7SBarry Smith PLogFlops(8*a->nz - 4*a->n); 474*bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 475*bcc973b7SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 476*bcc973b7SBarry Smith PetscFunctionReturn(0); 477*bcc973b7SBarry Smith } 478*bcc973b7SBarry Smith 479*bcc973b7SBarry Smith #undef __FUNC__ 480*bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultAdd_SeqMAIJ_4"></a>*/"MatMultAdd_SeqMAIJ_4" 481*bcc973b7SBarry Smith int MatMultAdd_SeqMAIJ_4(Mat A,Vec xx,Vec yy,Vec zz) 482*bcc973b7SBarry Smith { 483*bcc973b7SBarry Smith return 0; 484*bcc973b7SBarry Smith } 485*bcc973b7SBarry Smith #undef __FUNC__ 486*bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_SeqMAIJ_4"></a>*/"MatMultTransposeAdd_SeqMAIJ_4" 487*bcc973b7SBarry Smith int MatMultTransposeAdd_SeqMAIJ_4(Mat A,Vec xx,Vec yy,Vec zz) 488*bcc973b7SBarry Smith { 489*bcc973b7SBarry Smith return 0; 490*bcc973b7SBarry Smith } 491*bcc973b7SBarry Smith 492*bcc973b7SBarry Smith /* ---------------------------------------------------------------------------------- */ 49382b1193eSBarry Smith #undef __FUNC__ 494cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatCreateMAIJ"></a>*/"MatCreateMAIJ" 495cd3bbe55SBarry Smith int MatCreateMAIJ(Mat A,int dof,Mat *maij) 49682b1193eSBarry Smith { 497cd3bbe55SBarry Smith int ierr; 498cd3bbe55SBarry Smith Mat_SeqMAIJ *b; 49982b1193eSBarry Smith Mat B; 50082b1193eSBarry Smith 50182b1193eSBarry Smith PetscFunctionBegin; 502cd3bbe55SBarry Smith ierr = MATCreate(A->comm,dof*A->m,dof*A->n,dof*A->M,dof*A->N,&B);CHKERRQ(ierr); 503cd3bbe55SBarry Smith ierr = MatSetType(B,MATSEQMAIJ);CHKERRQ(ierr); 50482b1193eSBarry Smith 505cd3bbe55SBarry Smith B->assembled = PETSC_TRUE; 506cd3bbe55SBarry Smith B->ops->destroy = MatDestroy_SeqMAIJ; 507cd3bbe55SBarry Smith b = (Mat_SeqMAIJ*)B->data; 50882b1193eSBarry Smith 509cd3bbe55SBarry Smith b->AIJ = A; 510cd3bbe55SBarry Smith b->dof = dof; 511cd3bbe55SBarry Smith ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 512cd3bbe55SBarry Smith if (dof == 1) { 513cd3bbe55SBarry Smith B->ops->mult = MatMult_SeqMAIJ_1; 514cd3bbe55SBarry Smith B->ops->multadd = MatMultAdd_SeqMAIJ_1; 515cd3bbe55SBarry Smith B->ops->multtranspose = MatMultTranspose_SeqMAIJ_1; 516cd3bbe55SBarry Smith B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_1; 517cd3bbe55SBarry Smith } else if (dof == 2) { 518cd3bbe55SBarry Smith B->ops->mult = MatMult_SeqMAIJ_2; 519cd3bbe55SBarry Smith B->ops->multadd = MatMultAdd_SeqMAIJ_2; 520cd3bbe55SBarry Smith B->ops->multtranspose = MatMultTranspose_SeqMAIJ_2; 521cd3bbe55SBarry Smith B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_2; 522*bcc973b7SBarry Smith } else if (dof == 3) { 523*bcc973b7SBarry Smith B->ops->mult = MatMult_SeqMAIJ_3; 524*bcc973b7SBarry Smith B->ops->multadd = MatMultAdd_SeqMAIJ_3; 525*bcc973b7SBarry Smith B->ops->multtranspose = MatMultTranspose_SeqMAIJ_3; 526*bcc973b7SBarry Smith B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_3; 527*bcc973b7SBarry Smith } else if (dof == 4) { 528*bcc973b7SBarry Smith B->ops->mult = MatMult_SeqMAIJ_4; 529*bcc973b7SBarry Smith B->ops->multadd = MatMultAdd_SeqMAIJ_4; 530*bcc973b7SBarry Smith B->ops->multtranspose = MatMultTranspose_SeqMAIJ_4; 531*bcc973b7SBarry Smith B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_4; 53282b1193eSBarry Smith } else { 533cd3bbe55SBarry Smith SETERRQ1(1,1,"Cannot handle a dof of %d\n",dof); 53482b1193eSBarry Smith } 535cd3bbe55SBarry Smith *maij = B; 53682b1193eSBarry Smith PetscFunctionReturn(0); 53782b1193eSBarry Smith } 53882b1193eSBarry Smith 53982b1193eSBarry Smith 54082b1193eSBarry Smith 54182b1193eSBarry Smith 54282b1193eSBarry Smith 54382b1193eSBarry Smith 54482b1193eSBarry Smith 54582b1193eSBarry Smith 54682b1193eSBarry Smith 54782b1193eSBarry Smith 54882b1193eSBarry Smith 54982b1193eSBarry Smith 550