1*273d9f13SBarry Smith /*$Id: maij.c,v 1.10 2000/10/10 19:24:12 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 PetscFunctionReturn(0); 844cb9d9b8SBarry Smith } 854cb9d9b8SBarry Smith 864cb9d9b8SBarry Smith #undef __FUNC__ 874cb9d9b8SBarry Smith #define __FUNC__ /*<a name="MatDestroy_MPIMAIJ"></a>*/"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__ 115f4a53059SBarry Smith #define __FUNC__ /*<a name="MatCreate_MAIJ"></a>*/"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; 1224cb9d9b8SBarry Smith A->data = (void*)(b = PetscNew(Mat_MPIMAIJ));CHKPTRQ(b); 1234cb9d9b8SBarry Smith ierr = PetscMemzero(b,sizeof(Mat_MPIMAIJ));CHKERRQ(ierr); 124cd3bbe55SBarry Smith ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 125cd3bbe55SBarry Smith A->factor = 0; 126cd3bbe55SBarry Smith A->mapping = 0; 127f4a53059SBarry Smith 128cd3bbe55SBarry Smith b->AIJ = 0; 129cd3bbe55SBarry Smith b->dof = 0; 130f4a53059SBarry Smith b->OAIJ = 0; 131f4a53059SBarry Smith b->ctx = 0; 132f4a53059SBarry Smith b->w = 0; 13382b1193eSBarry Smith PetscFunctionReturn(0); 13482b1193eSBarry Smith } 13582b1193eSBarry Smith EXTERN_C_END 13682b1193eSBarry Smith 137cd3bbe55SBarry Smith /* --------------------------------------------------------------------------------------*/ 13882b1193eSBarry Smith #undef __FUNC__ 139cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMult_SeqMAIJ_2"></a>*/"MatMult_SeqMAIJ_2" 140cd3bbe55SBarry Smith int MatMult_SeqMAIJ_2(Mat A,Vec xx,Vec yy) 14182b1193eSBarry Smith { 1424cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 143bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 144bcc973b7SBarry Smith Scalar *x,*y,*v,sum1, sum2; 145*273d9f13SBarry Smith int ierr,m = b->AIJ->m,*idx,shift = a->indexshift,*ii; 146bcc973b7SBarry Smith int n,i,jrow,j; 14782b1193eSBarry Smith 148bcc973b7SBarry Smith PetscFunctionBegin; 149bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 150bcc973b7SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 151bcc973b7SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 152bcc973b7SBarry Smith idx = a->j; 153bcc973b7SBarry Smith v = a->a; 154bcc973b7SBarry Smith ii = a->i; 155bcc973b7SBarry Smith 156bcc973b7SBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 157bcc973b7SBarry Smith idx += shift; 158bcc973b7SBarry Smith for (i=0; i<m; i++) { 159bcc973b7SBarry Smith jrow = ii[i]; 160bcc973b7SBarry Smith n = ii[i+1] - jrow; 161bcc973b7SBarry Smith sum1 = 0.0; 162bcc973b7SBarry Smith sum2 = 0.0; 163bcc973b7SBarry Smith for (j=0; j<n; j++) { 164bcc973b7SBarry Smith sum1 += v[jrow]*x[2*idx[jrow]]; 165bcc973b7SBarry Smith sum2 += v[jrow]*x[2*idx[jrow]+1]; 166bcc973b7SBarry Smith jrow++; 167bcc973b7SBarry Smith } 168bcc973b7SBarry Smith y[2*i] = sum1; 169bcc973b7SBarry Smith y[2*i+1] = sum2; 170bcc973b7SBarry Smith } 171bcc973b7SBarry Smith 172bcc973b7SBarry Smith PLogFlops(4*a->nz - 2*m); 173bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 174bcc973b7SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 17582b1193eSBarry Smith PetscFunctionReturn(0); 17682b1193eSBarry Smith } 177bcc973b7SBarry Smith 17882b1193eSBarry Smith #undef __FUNC__ 179cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_SeqMAIJ_2"></a>*/"MatMultTranspose_SeqMAIJ_2" 180cd3bbe55SBarry Smith int MatMultTranspose_SeqMAIJ_2(Mat A,Vec xx,Vec yy) 18182b1193eSBarry Smith { 1824cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 183bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 184bcc973b7SBarry Smith Scalar *x,*y,*v,alpha1,alpha2,zero = 0.0; 185*273d9f13SBarry Smith int ierr,m = b->AIJ->m,n,i,*idx,shift = a->indexshift; 18682b1193eSBarry Smith 187bcc973b7SBarry Smith PetscFunctionBegin; 188bcc973b7SBarry Smith ierr = VecSet(&zero,yy);CHKERRQ(ierr); 189bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 190bcc973b7SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 191bcc973b7SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 192bcc973b7SBarry Smith for (i=0; i<m; i++) { 193bcc973b7SBarry Smith idx = a->j + a->i[i] + shift; 194bcc973b7SBarry Smith v = a->a + a->i[i] + shift; 195bcc973b7SBarry Smith n = a->i[i+1] - a->i[i]; 196bcc973b7SBarry Smith alpha1 = x[2*i]; 197bcc973b7SBarry Smith alpha2 = x[2*i+1]; 198bcc973b7SBarry Smith while (n-->0) {y[2*(*idx)] += alpha1*(*v); y[2*(*idx)+1] += alpha2*(*v); idx++; v++;} 199bcc973b7SBarry Smith } 200*273d9f13SBarry Smith PLogFlops(4*a->nz - 2*b->AIJ->n); 201bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 202bcc973b7SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 20382b1193eSBarry Smith PetscFunctionReturn(0); 20482b1193eSBarry Smith } 205bcc973b7SBarry Smith 20682b1193eSBarry Smith #undef __FUNC__ 207cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMultAdd_SeqMAIJ_2"></a>*/"MatMultAdd_SeqMAIJ_2" 208cd3bbe55SBarry Smith int MatMultAdd_SeqMAIJ_2(Mat A,Vec xx,Vec yy,Vec zz) 20982b1193eSBarry Smith { 2104cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 211bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 212bcc973b7SBarry Smith Scalar *x,*y,*v,sum1, sum2; 213*273d9f13SBarry Smith int ierr,m = b->AIJ->m,*idx,shift = a->indexshift,*ii; 214bcc973b7SBarry Smith int n,i,jrow,j; 21582b1193eSBarry Smith 216bcc973b7SBarry Smith PetscFunctionBegin; 217f9fae5adSBarry Smith if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);} 218bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 219f9fae5adSBarry Smith ierr = VecGetArray(zz,&y);CHKERRQ(ierr); 220bcc973b7SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 221bcc973b7SBarry Smith idx = a->j; 222bcc973b7SBarry Smith v = a->a; 223bcc973b7SBarry Smith ii = a->i; 224bcc973b7SBarry Smith 225bcc973b7SBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 226bcc973b7SBarry Smith idx += shift; 227bcc973b7SBarry Smith for (i=0; i<m; i++) { 228bcc973b7SBarry Smith jrow = ii[i]; 229bcc973b7SBarry Smith n = ii[i+1] - jrow; 230bcc973b7SBarry Smith sum1 = 0.0; 231bcc973b7SBarry Smith sum2 = 0.0; 232bcc973b7SBarry Smith for (j=0; j<n; j++) { 233bcc973b7SBarry Smith sum1 += v[jrow]*x[2*idx[jrow]]; 234bcc973b7SBarry Smith sum2 += v[jrow]*x[2*idx[jrow]+1]; 235bcc973b7SBarry Smith jrow++; 236bcc973b7SBarry Smith } 237bcc973b7SBarry Smith y[2*i] += sum1; 238bcc973b7SBarry Smith y[2*i+1] += sum2; 239bcc973b7SBarry Smith } 240bcc973b7SBarry Smith 241bcc973b7SBarry Smith PLogFlops(4*a->nz - 2*m); 242bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 243f9fae5adSBarry Smith ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr); 244bcc973b7SBarry Smith PetscFunctionReturn(0); 24582b1193eSBarry Smith } 24682b1193eSBarry Smith #undef __FUNC__ 247cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_SeqMAIJ_2"></a>*/"MatMultTransposeAdd_SeqMAIJ_2" 248cd3bbe55SBarry Smith int MatMultTransposeAdd_SeqMAIJ_2(Mat A,Vec xx,Vec yy,Vec zz) 24982b1193eSBarry Smith { 2504cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 251bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 252bcc973b7SBarry Smith Scalar *x,*y,*v,alpha1,alpha2; 253*273d9f13SBarry Smith int ierr,m = b->AIJ->m,n,i,*idx,shift = a->indexshift; 25482b1193eSBarry Smith 255bcc973b7SBarry Smith PetscFunctionBegin; 256f9fae5adSBarry Smith if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);} 257bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 258f9fae5adSBarry Smith ierr = VecGetArray(zz,&y);CHKERRQ(ierr); 259bcc973b7SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 260bcc973b7SBarry Smith for (i=0; i<m; i++) { 261bcc973b7SBarry Smith idx = a->j + a->i[i] + shift; 262bcc973b7SBarry Smith v = a->a + a->i[i] + shift; 263bcc973b7SBarry Smith n = a->i[i+1] - a->i[i]; 264bcc973b7SBarry Smith alpha1 = x[2*i]; 265bcc973b7SBarry Smith alpha2 = x[2*i+1]; 266bcc973b7SBarry Smith while (n-->0) {y[2*(*idx)] += alpha1*(*v); y[2*(*idx)+1] += alpha2*(*v); idx++; v++;} 267bcc973b7SBarry Smith } 268*273d9f13SBarry Smith PLogFlops(4*a->nz - 2*b->AIJ->n); 269bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 270f9fae5adSBarry Smith ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr); 271bcc973b7SBarry Smith PetscFunctionReturn(0); 27282b1193eSBarry Smith } 273cd3bbe55SBarry Smith /* --------------------------------------------------------------------------------------*/ 274bcc973b7SBarry Smith #undef __FUNC__ 275bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMult_SeqMAIJ_3"></a>*/"MatMult_SeqMAIJ_3" 276bcc973b7SBarry Smith int MatMult_SeqMAIJ_3(Mat A,Vec xx,Vec yy) 277bcc973b7SBarry Smith { 2784cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 279bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 280bcc973b7SBarry Smith Scalar *x,*y,*v,sum1, sum2, sum3; 281*273d9f13SBarry Smith int ierr,m = b->AIJ->m,*idx,shift = a->indexshift,*ii; 282bcc973b7SBarry Smith int n,i,jrow,j; 28382b1193eSBarry Smith 284bcc973b7SBarry Smith PetscFunctionBegin; 285bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 286bcc973b7SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 287bcc973b7SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 288bcc973b7SBarry Smith idx = a->j; 289bcc973b7SBarry Smith v = a->a; 290bcc973b7SBarry Smith ii = a->i; 291bcc973b7SBarry Smith 292bcc973b7SBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 293bcc973b7SBarry Smith idx += shift; 294bcc973b7SBarry Smith for (i=0; i<m; i++) { 295bcc973b7SBarry Smith jrow = ii[i]; 296bcc973b7SBarry Smith n = ii[i+1] - jrow; 297bcc973b7SBarry Smith sum1 = 0.0; 298bcc973b7SBarry Smith sum2 = 0.0; 299bcc973b7SBarry Smith sum3 = 0.0; 300bcc973b7SBarry Smith for (j=0; j<n; j++) { 301bcc973b7SBarry Smith sum1 += v[jrow]*x[3*idx[jrow]]; 302bcc973b7SBarry Smith sum2 += v[jrow]*x[3*idx[jrow]+1]; 303bcc973b7SBarry Smith sum3 += v[jrow]*x[3*idx[jrow]+2]; 304bcc973b7SBarry Smith jrow++; 305bcc973b7SBarry Smith } 306bcc973b7SBarry Smith y[3*i] = sum1; 307bcc973b7SBarry Smith y[3*i+1] = sum2; 308bcc973b7SBarry Smith y[3*i+2] = sum3; 309bcc973b7SBarry Smith } 310bcc973b7SBarry Smith 311bcc973b7SBarry Smith PLogFlops(6*a->nz - 3*m); 312bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 313bcc973b7SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 314bcc973b7SBarry Smith PetscFunctionReturn(0); 315bcc973b7SBarry Smith } 316bcc973b7SBarry Smith 317bcc973b7SBarry Smith #undef __FUNC__ 318bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_SeqMAIJ_3"></a>*/"MatMultTranspose_SeqMAIJ_3" 319bcc973b7SBarry Smith int MatMultTranspose_SeqMAIJ_3(Mat A,Vec xx,Vec yy) 320bcc973b7SBarry Smith { 3214cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 322bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 323bcc973b7SBarry Smith Scalar *x,*y,*v,alpha1,alpha2,alpha3,zero = 0.0; 324*273d9f13SBarry Smith int ierr,m = b->AIJ->m,n,i,*idx,shift = a->indexshift; 325bcc973b7SBarry Smith 326bcc973b7SBarry Smith PetscFunctionBegin; 327bcc973b7SBarry Smith ierr = VecSet(&zero,yy);CHKERRQ(ierr); 328bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 329bcc973b7SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 330bcc973b7SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 331bcc973b7SBarry Smith for (i=0; i<m; i++) { 332bcc973b7SBarry Smith idx = a->j + a->i[i] + shift; 333bcc973b7SBarry Smith v = a->a + a->i[i] + shift; 334bcc973b7SBarry Smith n = a->i[i+1] - a->i[i]; 335bcc973b7SBarry Smith alpha1 = x[3*i]; 336bcc973b7SBarry Smith alpha2 = x[3*i+1]; 337bcc973b7SBarry Smith alpha3 = x[3*i+2]; 338bcc973b7SBarry Smith while (n-->0) { 339bcc973b7SBarry Smith y[3*(*idx)] += alpha1*(*v); 340bcc973b7SBarry Smith y[3*(*idx)+1] += alpha2*(*v); 341bcc973b7SBarry Smith y[3*(*idx)+2] += alpha3*(*v); 342bcc973b7SBarry Smith idx++; v++; 343bcc973b7SBarry Smith } 344bcc973b7SBarry Smith } 345*273d9f13SBarry Smith PLogFlops(6*a->nz - 3*b->AIJ->n); 346bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 347bcc973b7SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 348bcc973b7SBarry Smith PetscFunctionReturn(0); 349bcc973b7SBarry Smith } 350bcc973b7SBarry Smith 351bcc973b7SBarry Smith #undef __FUNC__ 352bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultAdd_SeqMAIJ_3"></a>*/"MatMultAdd_SeqMAIJ_3" 353bcc973b7SBarry Smith int MatMultAdd_SeqMAIJ_3(Mat A,Vec xx,Vec yy,Vec zz) 354bcc973b7SBarry Smith { 3554cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 356bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 357bcc973b7SBarry Smith Scalar *x,*y,*v,sum1, sum2, sum3; 358*273d9f13SBarry Smith int ierr,m = b->AIJ->m,*idx,shift = a->indexshift,*ii; 359bcc973b7SBarry Smith int n,i,jrow,j; 360bcc973b7SBarry Smith 361bcc973b7SBarry Smith PetscFunctionBegin; 362f9fae5adSBarry Smith if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);} 363bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 364f9fae5adSBarry Smith ierr = VecGetArray(zz,&y);CHKERRQ(ierr); 365bcc973b7SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 366bcc973b7SBarry Smith idx = a->j; 367bcc973b7SBarry Smith v = a->a; 368bcc973b7SBarry Smith ii = a->i; 369bcc973b7SBarry Smith 370bcc973b7SBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 371bcc973b7SBarry Smith idx += shift; 372bcc973b7SBarry Smith for (i=0; i<m; i++) { 373bcc973b7SBarry Smith jrow = ii[i]; 374bcc973b7SBarry Smith n = ii[i+1] - jrow; 375bcc973b7SBarry Smith sum1 = 0.0; 376bcc973b7SBarry Smith sum2 = 0.0; 377bcc973b7SBarry Smith sum3 = 0.0; 378bcc973b7SBarry Smith for (j=0; j<n; j++) { 379bcc973b7SBarry Smith sum1 += v[jrow]*x[3*idx[jrow]]; 380bcc973b7SBarry Smith sum2 += v[jrow]*x[3*idx[jrow]+1]; 381bcc973b7SBarry Smith sum3 += v[jrow]*x[3*idx[jrow]+2]; 382bcc973b7SBarry Smith jrow++; 383bcc973b7SBarry Smith } 384bcc973b7SBarry Smith y[3*i] += sum1; 385bcc973b7SBarry Smith y[3*i+1] += sum2; 386bcc973b7SBarry Smith y[3*i+2] += sum3; 387bcc973b7SBarry Smith } 388bcc973b7SBarry Smith 389bcc973b7SBarry Smith PLogFlops(6*a->nz); 390bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 391f9fae5adSBarry Smith ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr); 392bcc973b7SBarry Smith PetscFunctionReturn(0); 393bcc973b7SBarry Smith } 394bcc973b7SBarry Smith #undef __FUNC__ 395bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_SeqMAIJ_3"></a>*/"MatMultTransposeAdd_SeqMAIJ_3" 396bcc973b7SBarry Smith int MatMultTransposeAdd_SeqMAIJ_3(Mat A,Vec xx,Vec yy,Vec zz) 397bcc973b7SBarry Smith { 3984cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 399bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 400bcc973b7SBarry Smith Scalar *x,*y,*v,alpha1,alpha2,alpha3; 401*273d9f13SBarry Smith int ierr,m = b->AIJ->m,n,i,*idx,shift = a->indexshift; 402bcc973b7SBarry Smith 403bcc973b7SBarry Smith PetscFunctionBegin; 404f9fae5adSBarry Smith if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);} 405bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 406f9fae5adSBarry Smith ierr = VecGetArray(zz,&y);CHKERRQ(ierr); 407bcc973b7SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 408bcc973b7SBarry Smith for (i=0; i<m; i++) { 409bcc973b7SBarry Smith idx = a->j + a->i[i] + shift; 410bcc973b7SBarry Smith v = a->a + a->i[i] + shift; 411bcc973b7SBarry Smith n = a->i[i+1] - a->i[i]; 412bcc973b7SBarry Smith alpha1 = x[3*i]; 413bcc973b7SBarry Smith alpha2 = x[3*i+1]; 414bcc973b7SBarry Smith alpha3 = x[3*i+2]; 415bcc973b7SBarry Smith while (n-->0) { 416bcc973b7SBarry Smith y[3*(*idx)] += alpha1*(*v); 417bcc973b7SBarry Smith y[3*(*idx)+1] += alpha2*(*v); 418bcc973b7SBarry Smith y[3*(*idx)+2] += alpha3*(*v); 419bcc973b7SBarry Smith idx++; v++; 420bcc973b7SBarry Smith } 421bcc973b7SBarry Smith } 422bcc973b7SBarry Smith PLogFlops(6*a->nz); 423bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 424f9fae5adSBarry Smith ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr); 425bcc973b7SBarry Smith PetscFunctionReturn(0); 426bcc973b7SBarry Smith } 427bcc973b7SBarry Smith 428bcc973b7SBarry Smith /* ------------------------------------------------------------------------------*/ 429bcc973b7SBarry Smith #undef __FUNC__ 430bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMult_SeqMAIJ_4"></a>*/"MatMult_SeqMAIJ_4" 431bcc973b7SBarry Smith int MatMult_SeqMAIJ_4(Mat A,Vec xx,Vec yy) 432bcc973b7SBarry Smith { 4334cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 434bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 435bcc973b7SBarry Smith Scalar *x,*y,*v,sum1, sum2, sum3, sum4; 436*273d9f13SBarry Smith int ierr,m = b->AIJ->m,*idx,shift = a->indexshift,*ii; 437bcc973b7SBarry Smith int n,i,jrow,j; 438bcc973b7SBarry Smith 439bcc973b7SBarry Smith PetscFunctionBegin; 440bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 441bcc973b7SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 442bcc973b7SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 443bcc973b7SBarry Smith idx = a->j; 444bcc973b7SBarry Smith v = a->a; 445bcc973b7SBarry Smith ii = a->i; 446bcc973b7SBarry Smith 447bcc973b7SBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 448bcc973b7SBarry Smith idx += shift; 449bcc973b7SBarry Smith for (i=0; i<m; i++) { 450bcc973b7SBarry Smith jrow = ii[i]; 451bcc973b7SBarry Smith n = ii[i+1] - jrow; 452bcc973b7SBarry Smith sum1 = 0.0; 453bcc973b7SBarry Smith sum2 = 0.0; 454bcc973b7SBarry Smith sum3 = 0.0; 455bcc973b7SBarry Smith sum4 = 0.0; 456bcc973b7SBarry Smith for (j=0; j<n; j++) { 457bcc973b7SBarry Smith sum1 += v[jrow]*x[4*idx[jrow]]; 458bcc973b7SBarry Smith sum2 += v[jrow]*x[4*idx[jrow]+1]; 459bcc973b7SBarry Smith sum3 += v[jrow]*x[4*idx[jrow]+2]; 460bcc973b7SBarry Smith sum4 += v[jrow]*x[4*idx[jrow]+3]; 461bcc973b7SBarry Smith jrow++; 462bcc973b7SBarry Smith } 463bcc973b7SBarry Smith y[4*i] = sum1; 464bcc973b7SBarry Smith y[4*i+1] = sum2; 465bcc973b7SBarry Smith y[4*i+2] = sum3; 466bcc973b7SBarry Smith y[4*i+3] = sum4; 467bcc973b7SBarry Smith } 468bcc973b7SBarry Smith 469bcc973b7SBarry Smith PLogFlops(8*a->nz - 4*m); 470bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 471bcc973b7SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 472bcc973b7SBarry Smith PetscFunctionReturn(0); 473bcc973b7SBarry Smith } 474bcc973b7SBarry Smith 475bcc973b7SBarry Smith #undef __FUNC__ 476bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_SeqMAIJ_4"></a>*/"MatMultTranspose_SeqMAIJ_4" 477bcc973b7SBarry Smith int MatMultTranspose_SeqMAIJ_4(Mat A,Vec xx,Vec yy) 478bcc973b7SBarry Smith { 4794cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 480bcc973b7SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 481bcc973b7SBarry Smith Scalar *x,*y,*v,alpha1,alpha2,alpha3,alpha4,zero = 0.0; 482*273d9f13SBarry Smith int ierr,m = b->AIJ->m,n,i,*idx,shift = a->indexshift; 483bcc973b7SBarry Smith 484bcc973b7SBarry Smith PetscFunctionBegin; 485bcc973b7SBarry Smith ierr = VecSet(&zero,yy);CHKERRQ(ierr); 486bcc973b7SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 487bcc973b7SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 488bcc973b7SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 489bcc973b7SBarry Smith for (i=0; i<m; i++) { 490bcc973b7SBarry Smith idx = a->j + a->i[i] + shift; 491bcc973b7SBarry Smith v = a->a + a->i[i] + shift; 492bcc973b7SBarry Smith n = a->i[i+1] - a->i[i]; 493bcc973b7SBarry Smith alpha1 = x[4*i]; 494bcc973b7SBarry Smith alpha2 = x[4*i+1]; 495bcc973b7SBarry Smith alpha3 = x[4*i+2]; 496bcc973b7SBarry Smith alpha4 = x[4*i+3]; 497bcc973b7SBarry Smith while (n-->0) { 498bcc973b7SBarry Smith y[4*(*idx)] += alpha1*(*v); 499bcc973b7SBarry Smith y[4*(*idx)+1] += alpha2*(*v); 500bcc973b7SBarry Smith y[4*(*idx)+2] += alpha3*(*v); 501bcc973b7SBarry Smith y[4*(*idx)+3] += alpha4*(*v); 502bcc973b7SBarry Smith idx++; v++; 503bcc973b7SBarry Smith } 504bcc973b7SBarry Smith } 505*273d9f13SBarry Smith PLogFlops(8*a->nz - 4*b->AIJ->n); 506bcc973b7SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 507bcc973b7SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 508bcc973b7SBarry Smith PetscFunctionReturn(0); 509bcc973b7SBarry Smith } 510bcc973b7SBarry Smith 511bcc973b7SBarry Smith #undef __FUNC__ 512bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultAdd_SeqMAIJ_4"></a>*/"MatMultAdd_SeqMAIJ_4" 513bcc973b7SBarry Smith int MatMultAdd_SeqMAIJ_4(Mat A,Vec xx,Vec yy,Vec zz) 514bcc973b7SBarry Smith { 5154cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 516f9fae5adSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 517f9fae5adSBarry Smith Scalar *x,*y,*v,sum1, sum2, sum3, sum4; 518*273d9f13SBarry Smith int ierr,m = b->AIJ->m,*idx,shift = a->indexshift,*ii; 519f9fae5adSBarry Smith int n,i,jrow,j; 520f9fae5adSBarry Smith 521f9fae5adSBarry Smith PetscFunctionBegin; 522f9fae5adSBarry Smith if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);} 523f9fae5adSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 524f9fae5adSBarry Smith ierr = VecGetArray(zz,&y);CHKERRQ(ierr); 525f9fae5adSBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 526f9fae5adSBarry Smith idx = a->j; 527f9fae5adSBarry Smith v = a->a; 528f9fae5adSBarry Smith ii = a->i; 529f9fae5adSBarry Smith 530f9fae5adSBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 531f9fae5adSBarry Smith idx += shift; 532f9fae5adSBarry Smith for (i=0; i<m; i++) { 533f9fae5adSBarry Smith jrow = ii[i]; 534f9fae5adSBarry Smith n = ii[i+1] - jrow; 535f9fae5adSBarry Smith sum1 = 0.0; 536f9fae5adSBarry Smith sum2 = 0.0; 537f9fae5adSBarry Smith sum3 = 0.0; 538f9fae5adSBarry Smith sum4 = 0.0; 539f9fae5adSBarry Smith for (j=0; j<n; j++) { 540f9fae5adSBarry Smith sum1 += v[jrow]*x[4*idx[jrow]]; 541f9fae5adSBarry Smith sum2 += v[jrow]*x[4*idx[jrow]+1]; 542f9fae5adSBarry Smith sum3 += v[jrow]*x[4*idx[jrow]+2]; 543f9fae5adSBarry Smith sum4 += v[jrow]*x[4*idx[jrow]+3]; 544f9fae5adSBarry Smith jrow++; 545f9fae5adSBarry Smith } 546f9fae5adSBarry Smith y[4*i] += sum1; 547f9fae5adSBarry Smith y[4*i+1] += sum2; 548f9fae5adSBarry Smith y[4*i+2] += sum3; 549f9fae5adSBarry Smith y[4*i+3] += sum4; 550f9fae5adSBarry Smith } 551f9fae5adSBarry Smith 552f9fae5adSBarry Smith PLogFlops(8*a->nz - 4*m); 553f9fae5adSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 554f9fae5adSBarry Smith ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr); 555f9fae5adSBarry Smith PetscFunctionReturn(0); 556bcc973b7SBarry Smith } 557bcc973b7SBarry Smith #undef __FUNC__ 558bcc973b7SBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_SeqMAIJ_4"></a>*/"MatMultTransposeAdd_SeqMAIJ_4" 559bcc973b7SBarry Smith int MatMultTransposeAdd_SeqMAIJ_4(Mat A,Vec xx,Vec yy,Vec zz) 560bcc973b7SBarry Smith { 5614cb9d9b8SBarry Smith Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; 562f9fae5adSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 563f9fae5adSBarry Smith Scalar *x,*y,*v,alpha1,alpha2,alpha3,alpha4; 564*273d9f13SBarry Smith int ierr,m = b->AIJ->m,n,i,*idx,shift = a->indexshift; 565f9fae5adSBarry Smith 566f9fae5adSBarry Smith PetscFunctionBegin; 567f9fae5adSBarry Smith if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);} 568f9fae5adSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 569f9fae5adSBarry Smith ierr = VecGetArray(zz,&y);CHKERRQ(ierr); 570f9fae5adSBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 571f9fae5adSBarry Smith for (i=0; i<m; i++) { 572f9fae5adSBarry Smith idx = a->j + a->i[i] + shift; 573f9fae5adSBarry Smith v = a->a + a->i[i] + shift; 574f9fae5adSBarry Smith n = a->i[i+1] - a->i[i]; 575f9fae5adSBarry Smith alpha1 = x[4*i]; 576f9fae5adSBarry Smith alpha2 = x[4*i+1]; 577f9fae5adSBarry Smith alpha3 = x[4*i+2]; 578f9fae5adSBarry Smith alpha4 = x[4*i+3]; 579f9fae5adSBarry Smith while (n-->0) { 580f9fae5adSBarry Smith y[4*(*idx)] += alpha1*(*v); 581f9fae5adSBarry Smith y[4*(*idx)+1] += alpha2*(*v); 582f9fae5adSBarry Smith y[4*(*idx)+2] += alpha3*(*v); 583f9fae5adSBarry Smith y[4*(*idx)+3] += alpha4*(*v); 584f9fae5adSBarry Smith idx++; v++; 585f9fae5adSBarry Smith } 586f9fae5adSBarry Smith } 587*273d9f13SBarry Smith PLogFlops(8*a->nz - 4*b->AIJ->n); 588f9fae5adSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 589f9fae5adSBarry Smith ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr); 590f9fae5adSBarry Smith PetscFunctionReturn(0); 591f9fae5adSBarry Smith 592f9fae5adSBarry Smith } 593f9fae5adSBarry Smith 594f9fae5adSBarry 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; 602*273d9f13SBarry 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 638f9fae5adSBarry Smith PLogFlops(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; 651*273d9f13SBarry 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*273d9f13SBarry Smith PLogFlops(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; 689*273d9f13SBarry 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 726f9fae5adSBarry Smith PLogFlops(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; 739*273d9f13SBarry 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 } 764f9fae5adSBarry Smith PLogFlops(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 770f4a53059SBarry Smith /*===================================================================================*/ 771f4a53059SBarry Smith #undef __FUNC__ 772f4a53059SBarry Smith #define __FUNC__ /*<a name="MatMult_MPIMAIJ_dof"></a>*/"MatMult_MPIMAIJ_dof" 773f4a53059SBarry Smith int MatMult_MPIMAIJ_dof(Mat A,Vec xx,Vec yy) 774f4a53059SBarry Smith { 7754cb9d9b8SBarry Smith Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data; 776f4a53059SBarry Smith int ierr; 777f4a53059SBarry Smith PetscFunctionBegin; 778f4a53059SBarry Smith 779f4a53059SBarry Smith /* start the scatter */ 780f4a53059SBarry Smith ierr = VecScatterBegin(xx,b->w,INSERT_VALUES,SCATTER_FORWARD,b->ctx);CHKERRQ(ierr); 7814cb9d9b8SBarry Smith ierr = (*b->AIJ->ops->mult)(b->AIJ,xx,yy);CHKERRQ(ierr); 782f4a53059SBarry Smith ierr = VecScatterEnd(xx,b->w,INSERT_VALUES,SCATTER_FORWARD,b->ctx);CHKERRQ(ierr); 7834cb9d9b8SBarry Smith ierr = (*b->OAIJ->ops->multadd)(b->OAIJ,b->w,yy,yy);CHKERRQ(ierr); 784f4a53059SBarry Smith PetscFunctionReturn(0); 785f4a53059SBarry Smith } 786f4a53059SBarry Smith 7874cb9d9b8SBarry Smith #undef __FUNC__ 7884cb9d9b8SBarry Smith #define __FUNC__ /*<a name="MatMultTranspose_MPIMAIJ_dof"></a>*/"MatMultTranspose_MPIMAIJ_dof" 7894cb9d9b8SBarry Smith int MatMultTranspose_MPIMAIJ_dof(Mat A,Vec xx,Vec yy) 7904cb9d9b8SBarry Smith { 7914cb9d9b8SBarry Smith Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data; 7924cb9d9b8SBarry Smith int ierr; 7934cb9d9b8SBarry Smith PetscFunctionBegin; 7944cb9d9b8SBarry Smith ierr = (*b->OAIJ->ops->multtranspose)(b->OAIJ,xx,b->w);CHKERRQ(ierr); 7954cb9d9b8SBarry Smith ierr = VecScatterBegin(b->w,yy,ADD_VALUES,SCATTER_REVERSE,b->ctx);CHKERRQ(ierr); 7964cb9d9b8SBarry Smith ierr = (*b->AIJ->ops->multtranspose)(b->AIJ,xx,yy);CHKERRQ(ierr); 7974cb9d9b8SBarry Smith ierr = VecScatterEnd(b->w,yy,ADD_VALUES,SCATTER_REVERSE,b->ctx);CHKERRQ(ierr); 7984cb9d9b8SBarry Smith PetscFunctionReturn(0); 7994cb9d9b8SBarry Smith } 8004cb9d9b8SBarry Smith 8014cb9d9b8SBarry Smith #undef __FUNC__ 8024cb9d9b8SBarry Smith #define __FUNC__ /*<a name="MatMultAdd_MPIMAIJ_dof"></a>*/"MatMultAdd_MPIMAIJ_dof" 803d72c5c08SBarry Smith int MatMultAdd_MPIMAIJ_dof(Mat A,Vec xx,Vec yy,Vec zz) 8044cb9d9b8SBarry Smith { 8054cb9d9b8SBarry Smith Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data; 8064cb9d9b8SBarry Smith int ierr; 8074cb9d9b8SBarry Smith PetscFunctionBegin; 8084cb9d9b8SBarry Smith 8094cb9d9b8SBarry Smith /* start the scatter */ 8104cb9d9b8SBarry Smith ierr = VecScatterBegin(xx,b->w,INSERT_VALUES,SCATTER_FORWARD,b->ctx);CHKERRQ(ierr); 811d72c5c08SBarry Smith ierr = (*b->AIJ->ops->multadd)(b->AIJ,xx,yy,zz);CHKERRQ(ierr); 8124cb9d9b8SBarry Smith ierr = VecScatterEnd(xx,b->w,INSERT_VALUES,SCATTER_FORWARD,b->ctx);CHKERRQ(ierr); 813d72c5c08SBarry Smith ierr = (*b->OAIJ->ops->multadd)(b->OAIJ,b->w,yy,zz);CHKERRQ(ierr); 8144cb9d9b8SBarry Smith PetscFunctionReturn(0); 8154cb9d9b8SBarry Smith } 8164cb9d9b8SBarry Smith 8174cb9d9b8SBarry Smith #undef __FUNC__ 8184cb9d9b8SBarry Smith #define __FUNC__ /*<a name="MatMultTransposeAdd_MPIMAIJ_dof"></a>*/"MatMultTransposeAdd_MPIMAIJ_dof" 819d72c5c08SBarry Smith int MatMultTransposeAdd_MPIMAIJ_dof(Mat A,Vec xx,Vec yy,Vec zz) 8204cb9d9b8SBarry Smith { 8214cb9d9b8SBarry Smith Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data; 8224cb9d9b8SBarry Smith int ierr; 8234cb9d9b8SBarry Smith PetscFunctionBegin; 8244cb9d9b8SBarry Smith ierr = (*b->OAIJ->ops->multtranspose)(b->OAIJ,xx,b->w);CHKERRQ(ierr); 825d72c5c08SBarry Smith ierr = VecScatterBegin(b->w,zz,ADD_VALUES,SCATTER_REVERSE,b->ctx);CHKERRQ(ierr); 826d72c5c08SBarry Smith ierr = (*b->AIJ->ops->multtransposeadd)(b->AIJ,xx,yy,zz);CHKERRQ(ierr); 827d72c5c08SBarry Smith ierr = VecScatterEnd(b->w,zz,ADD_VALUES,SCATTER_REVERSE,b->ctx);CHKERRQ(ierr); 8284cb9d9b8SBarry Smith PetscFunctionReturn(0); 8294cb9d9b8SBarry Smith } 8304cb9d9b8SBarry Smith 831bcc973b7SBarry Smith /* ---------------------------------------------------------------------------------- */ 83282b1193eSBarry Smith #undef __FUNC__ 833cd3bbe55SBarry Smith #define __FUNC__ /*<a name="MatCreateMAIJ"></a>*/"MatCreateMAIJ" 834cd3bbe55SBarry Smith int MatCreateMAIJ(Mat A,int dof,Mat *maij) 83582b1193eSBarry Smith { 836f4a53059SBarry Smith int ierr,size,n; 8374cb9d9b8SBarry Smith Mat_MPIMAIJ *b; 83882b1193eSBarry Smith Mat B; 83982b1193eSBarry Smith 84082b1193eSBarry Smith PetscFunctionBegin; 841d72c5c08SBarry Smith ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 842d72c5c08SBarry Smith 843d72c5c08SBarry Smith if (dof == 1) { 844d72c5c08SBarry Smith *maij = A; 845d72c5c08SBarry Smith } else { 84683903688SBarry Smith ierr = MatCreate(A->comm,dof*A->m,dof*A->n,dof*A->M,dof*A->N,&B);CHKERRQ(ierr); 847cd3bbe55SBarry Smith B->assembled = PETSC_TRUE; 848d72c5c08SBarry Smith 849f4a53059SBarry Smith ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 850f4a53059SBarry Smith if (size == 1) { 851b9b97703SBarry Smith ierr = MatSetType(B,MATSEQMAIJ);CHKERRQ(ierr); 8524cb9d9b8SBarry Smith B->ops->destroy = MatDestroy_SeqMAIJ; 853b9b97703SBarry Smith b = (Mat_MPIMAIJ*)B->data; 854b9b97703SBarry Smith b->dof = dof; 8554cb9d9b8SBarry Smith b->AIJ = A; 856d72c5c08SBarry Smith if (dof == 2) { 857cd3bbe55SBarry Smith B->ops->mult = MatMult_SeqMAIJ_2; 858cd3bbe55SBarry Smith B->ops->multadd = MatMultAdd_SeqMAIJ_2; 859cd3bbe55SBarry Smith B->ops->multtranspose = MatMultTranspose_SeqMAIJ_2; 860cd3bbe55SBarry Smith B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_2; 861bcc973b7SBarry Smith } else if (dof == 3) { 862bcc973b7SBarry Smith B->ops->mult = MatMult_SeqMAIJ_3; 863bcc973b7SBarry Smith B->ops->multadd = MatMultAdd_SeqMAIJ_3; 864bcc973b7SBarry Smith B->ops->multtranspose = MatMultTranspose_SeqMAIJ_3; 865bcc973b7SBarry Smith B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_3; 866bcc973b7SBarry Smith } else if (dof == 4) { 867bcc973b7SBarry Smith B->ops->mult = MatMult_SeqMAIJ_4; 868bcc973b7SBarry Smith B->ops->multadd = MatMultAdd_SeqMAIJ_4; 869bcc973b7SBarry Smith B->ops->multtranspose = MatMultTranspose_SeqMAIJ_4; 870bcc973b7SBarry Smith B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_4; 871f9fae5adSBarry Smith } else if (dof == 5) { 872f9fae5adSBarry Smith B->ops->mult = MatMult_SeqMAIJ_5; 873f9fae5adSBarry Smith B->ops->multadd = MatMultAdd_SeqMAIJ_5; 874f9fae5adSBarry Smith B->ops->multtranspose = MatMultTranspose_SeqMAIJ_5; 875f9fae5adSBarry Smith B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_5; 87682b1193eSBarry Smith } else { 87729bbc08cSBarry Smith SETERRQ1(1,"Cannot handle a dof of %d\n",dof); 87882b1193eSBarry Smith } 879f4a53059SBarry Smith } else { 880f4a53059SBarry Smith Mat_MPIAIJ *mpiaij = (Mat_MPIAIJ *)A->data; 881f4a53059SBarry Smith IS from,to; 882f4a53059SBarry Smith Vec gvec; 883f4a53059SBarry Smith int *garray,i; 884f4a53059SBarry Smith 885b9b97703SBarry Smith ierr = MatSetType(B,MATMPIMAIJ);CHKERRQ(ierr); 886d72c5c08SBarry Smith B->ops->destroy = MatDestroy_MPIMAIJ; 887b9b97703SBarry Smith b = (Mat_MPIMAIJ*)B->data; 888b9b97703SBarry Smith b->dof = dof; 889b9b97703SBarry Smith b->A = A; 8904cb9d9b8SBarry Smith ierr = MatCreateMAIJ(mpiaij->A,dof,&b->AIJ);CHKERRQ(ierr); 8914cb9d9b8SBarry Smith ierr = MatCreateMAIJ(mpiaij->B,dof,&b->OAIJ);CHKERRQ(ierr); 8924cb9d9b8SBarry Smith 893f4a53059SBarry Smith ierr = VecGetSize(mpiaij->lvec,&n);CHKERRQ(ierr); 894f4a53059SBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,n*dof,&b->w);CHKERRQ(ierr); 895f4a53059SBarry Smith 896f4a53059SBarry Smith /* create two temporary Index sets for build scatter gather */ 897f4a53059SBarry Smith garray = (int*)PetscMalloc((n+1)*sizeof(int));CHKPTRQ(garray); 898f4a53059SBarry Smith for (i=0; i<n; i++) garray[i] = dof*mpiaij->garray[i]; 899f4a53059SBarry Smith ierr = ISCreateBlock(A->comm,dof,n,garray,&from);CHKERRQ(ierr); 900f4a53059SBarry Smith ierr = PetscFree(garray);CHKERRQ(ierr); 901f4a53059SBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,n*dof,0,1,&to);CHKERRQ(ierr); 902f4a53059SBarry Smith 903f4a53059SBarry Smith /* create temporary global vector to generate scatter context */ 904f4a53059SBarry Smith ierr = VecCreateMPI(A->comm,dof*A->n,dof*A->N,&gvec);CHKERRQ(ierr); 905f4a53059SBarry Smith 906f4a53059SBarry Smith /* generate the scatter context */ 907f4a53059SBarry Smith ierr = VecScatterCreate(gvec,from,b->w,to,&b->ctx);CHKERRQ(ierr); 908f4a53059SBarry Smith 909f4a53059SBarry Smith ierr = ISDestroy(from);CHKERRQ(ierr); 910f4a53059SBarry Smith ierr = ISDestroy(to);CHKERRQ(ierr); 911f4a53059SBarry Smith ierr = VecDestroy(gvec);CHKERRQ(ierr); 912f4a53059SBarry Smith 913f4a53059SBarry Smith B->ops->mult = MatMult_MPIMAIJ_dof; 9144cb9d9b8SBarry Smith B->ops->multtranspose = MatMultTranspose_MPIMAIJ_dof; 9154cb9d9b8SBarry Smith B->ops->multadd = MatMultAdd_MPIMAIJ_dof; 9164cb9d9b8SBarry Smith B->ops->multtransposeadd = MatMultTransposeAdd_MPIMAIJ_dof; 917f4a53059SBarry Smith } 918cd3bbe55SBarry Smith *maij = B; 919d72c5c08SBarry Smith } 92082b1193eSBarry Smith PetscFunctionReturn(0); 92182b1193eSBarry Smith } 92282b1193eSBarry Smith 92382b1193eSBarry Smith 92482b1193eSBarry Smith 92582b1193eSBarry Smith 92682b1193eSBarry Smith 92782b1193eSBarry Smith 92882b1193eSBarry Smith 92982b1193eSBarry Smith 93082b1193eSBarry Smith 93182b1193eSBarry Smith 93282b1193eSBarry Smith 93382b1193eSBarry Smith 934