/*$Id: maij.c,v 1.6 2000/06/27 17:46:45 bsmith Exp bsmith $*/ /* Defines the basic matrix operations for the MAIJ matrix storage format. This format is used for restriction and interpolation operations for multicomponent problems. It interpolates each component the same way independently. We provide: MatMult() MatMultTranspose() MatMultTransposeAdd() MatMultAdd() and MatCreateMAIJ(Mat,dof,Mat*) This single directory handles both the sequential and parallel codes */ #include "src/mat/impls/aij/mpi/mpiaij.h" typedef struct { int dof; /* number of components */ Mat AIJ; /* representation of interpolation for one component */ } Mat_SeqMAIJ; typedef struct { int dof; /* number of components */ Mat AIJ,OAIJ; /* representation of interpolation for one component */ Mat A; VecScatter ctx; /* update ghost points for parallel case */ Vec w; /* work space for ghost values for parallel case */ } Mat_MPIMAIJ; #undef __FUNC__ #define __FUNC__ /**/"MatDestroy_SeqMAIJ" int MatDestroy_SeqMAIJ(Mat A) { int ierr; Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; PetscFunctionBegin; if (b->AIJ) { ierr = MatDestroy(b->AIJ);CHKERRQ(ierr); } ierr = PetscFree(b);CHKERRQ(ierr); PetscHeaderDestroy(A); PetscFunctionReturn(0); } #undef __FUNC__ #define __FUNC__ /**/"MatDestroy_MPIMAIJ" int MatDestroy_MPIMAIJ(Mat A) { int ierr; Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data; PetscFunctionBegin; if (b->A) { ierr = MatDestroy(b->A);CHKERRQ(ierr); } if (b->AIJ) { ierr = MatDestroy(b->AIJ);CHKERRQ(ierr); } if (b->OAIJ) { ierr = MatDestroy(b->OAIJ);CHKERRQ(ierr); } if (b->ctx) { ierr = VecScatterDestroy(b->ctx);CHKERRQ(ierr); } if (b->w) { ierr = VecDestroy(b->w);CHKERRQ(ierr); } ierr = PetscFree(b);CHKERRQ(ierr); PetscHeaderDestroy(A); PetscFunctionReturn(0); } EXTERN_C_BEGIN #undef __FUNC__ #define __FUNC__ /**/"MatCreate_MAIJ" int MatCreate_MAIJ(Mat A) { int ierr; Mat_MPIMAIJ *b; PetscFunctionBegin; A->data = (void*)(b = PetscNew(Mat_MPIMAIJ));CHKPTRQ(b); ierr = PetscMemzero(b,sizeof(Mat_MPIMAIJ));CHKERRQ(ierr); ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); A->factor = 0; A->mapping = 0; b->AIJ = 0; b->dof = 0; b->OAIJ = 0; b->ctx = 0; b->w = 0; PetscFunctionReturn(0); } EXTERN_C_END /* --------------------------------------------------------------------------------------*/ #undef __FUNC__ #define __FUNC__ /**/"MatMult_SeqMAIJ_2" int MatMult_SeqMAIJ_2(Mat A,Vec xx,Vec yy) { Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; Scalar *x,*y,*v,sum1, sum2; int ierr,m = a->m,*idx,shift = a->indexshift,*ii; int n,i,jrow,j; PetscFunctionBegin; ierr = VecGetArray(xx,&x);CHKERRQ(ierr); ierr = VecGetArray(yy,&y);CHKERRQ(ierr); x = x + shift; /* shift for Fortran start by 1 indexing */ idx = a->j; v = a->a; ii = a->i; v += shift; /* shift for Fortran start by 1 indexing */ idx += shift; for (i=0; inz - 2*m); ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNC__ #define __FUNC__ /**/"MatMultTranspose_SeqMAIJ_2" int MatMultTranspose_SeqMAIJ_2(Mat A,Vec xx,Vec yy) { Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; Scalar *x,*y,*v,alpha1,alpha2,zero = 0.0; int ierr,m = a->m,n,i,*idx,shift = a->indexshift; PetscFunctionBegin; ierr = VecSet(&zero,yy);CHKERRQ(ierr); ierr = VecGetArray(xx,&x);CHKERRQ(ierr); ierr = VecGetArray(yy,&y);CHKERRQ(ierr); y = y + shift; /* shift for Fortran start by 1 indexing */ for (i=0; ij + a->i[i] + shift; v = a->a + a->i[i] + shift; n = a->i[i+1] - a->i[i]; alpha1 = x[2*i]; alpha2 = x[2*i+1]; while (n-->0) {y[2*(*idx)] += alpha1*(*v); y[2*(*idx)+1] += alpha2*(*v); idx++; v++;} } PLogFlops(4*a->nz - 2*a->n); ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNC__ #define __FUNC__ /**/"MatMultAdd_SeqMAIJ_2" int MatMultAdd_SeqMAIJ_2(Mat A,Vec xx,Vec yy,Vec zz) { Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; Scalar *x,*y,*v,sum1, sum2; int ierr,m = a->m,*idx,shift = a->indexshift,*ii; int n,i,jrow,j; PetscFunctionBegin; if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);} ierr = VecGetArray(xx,&x);CHKERRQ(ierr); ierr = VecGetArray(zz,&y);CHKERRQ(ierr); x = x + shift; /* shift for Fortran start by 1 indexing */ idx = a->j; v = a->a; ii = a->i; v += shift; /* shift for Fortran start by 1 indexing */ idx += shift; for (i=0; inz - 2*m); ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNC__ #define __FUNC__ /**/"MatMultTransposeAdd_SeqMAIJ_2" int MatMultTransposeAdd_SeqMAIJ_2(Mat A,Vec xx,Vec yy,Vec zz) { Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; Scalar *x,*y,*v,alpha1,alpha2; int ierr,m = a->m,n,i,*idx,shift = a->indexshift; PetscFunctionBegin; if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);} ierr = VecGetArray(xx,&x);CHKERRQ(ierr); ierr = VecGetArray(zz,&y);CHKERRQ(ierr); y = y + shift; /* shift for Fortran start by 1 indexing */ for (i=0; ij + a->i[i] + shift; v = a->a + a->i[i] + shift; n = a->i[i+1] - a->i[i]; alpha1 = x[2*i]; alpha2 = x[2*i+1]; while (n-->0) {y[2*(*idx)] += alpha1*(*v); y[2*(*idx)+1] += alpha2*(*v); idx++; v++;} } PLogFlops(4*a->nz - 2*a->n); ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr); PetscFunctionReturn(0); } /* --------------------------------------------------------------------------------------*/ #undef __FUNC__ #define __FUNC__ /**/"MatMult_SeqMAIJ_3" int MatMult_SeqMAIJ_3(Mat A,Vec xx,Vec yy) { Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; Scalar *x,*y,*v,sum1, sum2, sum3; int ierr,m = a->m,*idx,shift = a->indexshift,*ii; int n,i,jrow,j; PetscFunctionBegin; ierr = VecGetArray(xx,&x);CHKERRQ(ierr); ierr = VecGetArray(yy,&y);CHKERRQ(ierr); x = x + shift; /* shift for Fortran start by 1 indexing */ idx = a->j; v = a->a; ii = a->i; v += shift; /* shift for Fortran start by 1 indexing */ idx += shift; for (i=0; inz - 3*m); ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNC__ #define __FUNC__ /**/"MatMultTranspose_SeqMAIJ_3" int MatMultTranspose_SeqMAIJ_3(Mat A,Vec xx,Vec yy) { Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; Scalar *x,*y,*v,alpha1,alpha2,alpha3,zero = 0.0; int ierr,m = a->m,n,i,*idx,shift = a->indexshift; PetscFunctionBegin; ierr = VecSet(&zero,yy);CHKERRQ(ierr); ierr = VecGetArray(xx,&x);CHKERRQ(ierr); ierr = VecGetArray(yy,&y);CHKERRQ(ierr); y = y + shift; /* shift for Fortran start by 1 indexing */ for (i=0; ij + a->i[i] + shift; v = a->a + a->i[i] + shift; n = a->i[i+1] - a->i[i]; alpha1 = x[3*i]; alpha2 = x[3*i+1]; alpha3 = x[3*i+2]; while (n-->0) { y[3*(*idx)] += alpha1*(*v); y[3*(*idx)+1] += alpha2*(*v); y[3*(*idx)+2] += alpha3*(*v); idx++; v++; } } PLogFlops(6*a->nz - 3*a->n); ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNC__ #define __FUNC__ /**/"MatMultAdd_SeqMAIJ_3" int MatMultAdd_SeqMAIJ_3(Mat A,Vec xx,Vec yy,Vec zz) { Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; Scalar *x,*y,*v,sum1, sum2, sum3; int ierr,m = a->m,*idx,shift = a->indexshift,*ii; int n,i,jrow,j; PetscFunctionBegin; if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);} ierr = VecGetArray(xx,&x);CHKERRQ(ierr); ierr = VecGetArray(zz,&y);CHKERRQ(ierr); x = x + shift; /* shift for Fortran start by 1 indexing */ idx = a->j; v = a->a; ii = a->i; v += shift; /* shift for Fortran start by 1 indexing */ idx += shift; for (i=0; inz); ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNC__ #define __FUNC__ /**/"MatMultTransposeAdd_SeqMAIJ_3" int MatMultTransposeAdd_SeqMAIJ_3(Mat A,Vec xx,Vec yy,Vec zz) { Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; Scalar *x,*y,*v,alpha1,alpha2,alpha3; int ierr,m = a->m,n,i,*idx,shift = a->indexshift; PetscFunctionBegin; if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);} ierr = VecGetArray(xx,&x);CHKERRQ(ierr); ierr = VecGetArray(zz,&y);CHKERRQ(ierr); y = y + shift; /* shift for Fortran start by 1 indexing */ for (i=0; ij + a->i[i] + shift; v = a->a + a->i[i] + shift; n = a->i[i+1] - a->i[i]; alpha1 = x[3*i]; alpha2 = x[3*i+1]; alpha3 = x[3*i+2]; while (n-->0) { y[3*(*idx)] += alpha1*(*v); y[3*(*idx)+1] += alpha2*(*v); y[3*(*idx)+2] += alpha3*(*v); idx++; v++; } } PLogFlops(6*a->nz); ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr); PetscFunctionReturn(0); } /* ------------------------------------------------------------------------------*/ #undef __FUNC__ #define __FUNC__ /**/"MatMult_SeqMAIJ_4" int MatMult_SeqMAIJ_4(Mat A,Vec xx,Vec yy) { Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; Scalar *x,*y,*v,sum1, sum2, sum3, sum4; int ierr,m = a->m,*idx,shift = a->indexshift,*ii; int n,i,jrow,j; PetscFunctionBegin; ierr = VecGetArray(xx,&x);CHKERRQ(ierr); ierr = VecGetArray(yy,&y);CHKERRQ(ierr); x = x + shift; /* shift for Fortran start by 1 indexing */ idx = a->j; v = a->a; ii = a->i; v += shift; /* shift for Fortran start by 1 indexing */ idx += shift; for (i=0; inz - 4*m); ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNC__ #define __FUNC__ /**/"MatMultTranspose_SeqMAIJ_4" int MatMultTranspose_SeqMAIJ_4(Mat A,Vec xx,Vec yy) { Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; Scalar *x,*y,*v,alpha1,alpha2,alpha3,alpha4,zero = 0.0; int ierr,m = a->m,n,i,*idx,shift = a->indexshift; PetscFunctionBegin; ierr = VecSet(&zero,yy);CHKERRQ(ierr); ierr = VecGetArray(xx,&x);CHKERRQ(ierr); ierr = VecGetArray(yy,&y);CHKERRQ(ierr); y = y + shift; /* shift for Fortran start by 1 indexing */ for (i=0; ij + a->i[i] + shift; v = a->a + a->i[i] + shift; n = a->i[i+1] - a->i[i]; alpha1 = x[4*i]; alpha2 = x[4*i+1]; alpha3 = x[4*i+2]; alpha4 = x[4*i+3]; while (n-->0) { y[4*(*idx)] += alpha1*(*v); y[4*(*idx)+1] += alpha2*(*v); y[4*(*idx)+2] += alpha3*(*v); y[4*(*idx)+3] += alpha4*(*v); idx++; v++; } } PLogFlops(8*a->nz - 4*a->n); ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNC__ #define __FUNC__ /**/"MatMultAdd_SeqMAIJ_4" int MatMultAdd_SeqMAIJ_4(Mat A,Vec xx,Vec yy,Vec zz) { Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; Scalar *x,*y,*v,sum1, sum2, sum3, sum4; int ierr,m = a->m,*idx,shift = a->indexshift,*ii; int n,i,jrow,j; PetscFunctionBegin; if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);} ierr = VecGetArray(xx,&x);CHKERRQ(ierr); ierr = VecGetArray(zz,&y);CHKERRQ(ierr); x = x + shift; /* shift for Fortran start by 1 indexing */ idx = a->j; v = a->a; ii = a->i; v += shift; /* shift for Fortran start by 1 indexing */ idx += shift; for (i=0; inz - 4*m); ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNC__ #define __FUNC__ /**/"MatMultTransposeAdd_SeqMAIJ_4" int MatMultTransposeAdd_SeqMAIJ_4(Mat A,Vec xx,Vec yy,Vec zz) { Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; Scalar *x,*y,*v,alpha1,alpha2,alpha3,alpha4; int ierr,m = a->m,n,i,*idx,shift = a->indexshift; PetscFunctionBegin; if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);} ierr = VecGetArray(xx,&x);CHKERRQ(ierr); ierr = VecGetArray(zz,&y);CHKERRQ(ierr); y = y + shift; /* shift for Fortran start by 1 indexing */ for (i=0; ij + a->i[i] + shift; v = a->a + a->i[i] + shift; n = a->i[i+1] - a->i[i]; alpha1 = x[4*i]; alpha2 = x[4*i+1]; alpha3 = x[4*i+2]; alpha4 = x[4*i+3]; while (n-->0) { y[4*(*idx)] += alpha1*(*v); y[4*(*idx)+1] += alpha2*(*v); y[4*(*idx)+2] += alpha3*(*v); y[4*(*idx)+3] += alpha4*(*v); idx++; v++; } } PLogFlops(8*a->nz - 4*a->n); ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr); PetscFunctionReturn(0); } /* ------------------------------------------------------------------------------*/ #undef __FUNC__ #define __FUNC__ /**/"MatMult_SeqMAIJ_5" int MatMult_SeqMAIJ_5(Mat A,Vec xx,Vec yy) { Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; Scalar *x,*y,*v,sum1, sum2, sum3, sum4, sum5; int ierr,m = a->m,*idx,shift = a->indexshift,*ii; int n,i,jrow,j; PetscFunctionBegin; ierr = VecGetArray(xx,&x);CHKERRQ(ierr); ierr = VecGetArray(yy,&y);CHKERRQ(ierr); x = x + shift; /* shift for Fortran start by 1 indexing */ idx = a->j; v = a->a; ii = a->i; v += shift; /* shift for Fortran start by 1 indexing */ idx += shift; for (i=0; inz - 5*m); ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNC__ #define __FUNC__ /**/"MatMultTranspose_SeqMAIJ_5" int MatMultTranspose_SeqMAIJ_5(Mat A,Vec xx,Vec yy) { Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; Scalar *x,*y,*v,alpha1,alpha2,alpha3,alpha4,alpha5,zero = 0.0; int ierr,m = a->m,n,i,*idx,shift = a->indexshift; PetscFunctionBegin; ierr = VecSet(&zero,yy);CHKERRQ(ierr); ierr = VecGetArray(xx,&x);CHKERRQ(ierr); ierr = VecGetArray(yy,&y);CHKERRQ(ierr); y = y + shift; /* shift for Fortran start by 1 indexing */ for (i=0; ij + a->i[i] + shift; v = a->a + a->i[i] + shift; n = a->i[i+1] - a->i[i]; alpha1 = x[5*i]; alpha2 = x[5*i+1]; alpha3 = x[5*i+2]; alpha4 = x[5*i+3]; alpha5 = x[5*i+4]; while (n-->0) { y[5*(*idx)] += alpha1*(*v); y[5*(*idx)+1] += alpha2*(*v); y[5*(*idx)+2] += alpha3*(*v); y[5*(*idx)+3] += alpha4*(*v); y[5*(*idx)+4] += alpha5*(*v); idx++; v++; } } PLogFlops(10*a->nz - 5*a->n); ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNC__ #define __FUNC__ /**/"MatMultAdd_SeqMAIJ_5" int MatMultAdd_SeqMAIJ_5(Mat A,Vec xx,Vec yy,Vec zz) { Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; Scalar *x,*y,*v,sum1, sum2, sum3, sum4, sum5; int ierr,m = a->m,*idx,shift = a->indexshift,*ii; int n,i,jrow,j; PetscFunctionBegin; if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);} ierr = VecGetArray(xx,&x);CHKERRQ(ierr); ierr = VecGetArray(zz,&y);CHKERRQ(ierr); x = x + shift; /* shift for Fortran start by 1 indexing */ idx = a->j; v = a->a; ii = a->i; v += shift; /* shift for Fortran start by 1 indexing */ idx += shift; for (i=0; inz); ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNC__ #define __FUNC__ /**/"MatMultTransposeAdd_SeqMAIJ_5" int MatMultTransposeAdd_SeqMAIJ_5(Mat A,Vec xx,Vec yy,Vec zz) { Mat_SeqMAIJ *b = (Mat_SeqMAIJ*)A->data; Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; Scalar *x,*y,*v,alpha1,alpha2,alpha3,alpha4,alpha5; int ierr,m = a->m,n,i,*idx,shift = a->indexshift; PetscFunctionBegin; if (yy != zz) {ierr = VecCopy(yy,zz);CHKERRQ(ierr);} ierr = VecGetArray(xx,&x);CHKERRQ(ierr); ierr = VecGetArray(zz,&y);CHKERRQ(ierr); y = y + shift; /* shift for Fortran start by 1 indexing */ for (i=0; ij + a->i[i] + shift; v = a->a + a->i[i] + shift; n = a->i[i+1] - a->i[i]; alpha1 = x[5*i]; alpha2 = x[5*i+1]; alpha3 = x[5*i+2]; alpha4 = x[5*i+3]; alpha5 = x[5*i+4]; while (n-->0) { y[5*(*idx)] += alpha1*(*v); y[5*(*idx)+1] += alpha2*(*v); y[5*(*idx)+2] += alpha3*(*v); y[5*(*idx)+3] += alpha4*(*v); y[5*(*idx)+4] += alpha5*(*v); idx++; v++; } } PLogFlops(10*a->nz); ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr); PetscFunctionReturn(0); } /*===================================================================================*/ #undef __FUNC__ #define __FUNC__ /**/"MatMult_MPIMAIJ_dof" int MatMult_MPIMAIJ_dof(Mat A,Vec xx,Vec yy) { Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data; int ierr; PetscFunctionBegin; /* start the scatter */ ierr = VecScatterBegin(xx,b->w,INSERT_VALUES,SCATTER_FORWARD,b->ctx);CHKERRQ(ierr); ierr = (*b->AIJ->ops->mult)(b->AIJ,xx,yy);CHKERRQ(ierr); ierr = VecScatterEnd(xx,b->w,INSERT_VALUES,SCATTER_FORWARD,b->ctx);CHKERRQ(ierr); ierr = (*b->OAIJ->ops->multadd)(b->OAIJ,b->w,yy,yy);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNC__ #define __FUNC__ /**/"MatMultTranspose_MPIMAIJ_dof" int MatMultTranspose_MPIMAIJ_dof(Mat A,Vec xx,Vec yy) { Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data; int ierr; PetscFunctionBegin; ierr = (*b->OAIJ->ops->multtranspose)(b->OAIJ,xx,b->w);CHKERRQ(ierr); ierr = VecScatterBegin(b->w,yy,ADD_VALUES,SCATTER_REVERSE,b->ctx);CHKERRQ(ierr); ierr = (*b->AIJ->ops->multtranspose)(b->AIJ,xx,yy);CHKERRQ(ierr); ierr = VecScatterEnd(b->w,yy,ADD_VALUES,SCATTER_REVERSE,b->ctx);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNC__ #define __FUNC__ /**/"MatMultAdd_MPIMAIJ_dof" int MatMultAdd_MPIMAIJ_dof(Mat A,Vec xx,Vec yy,Vec zz) { Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data; int ierr; PetscFunctionBegin; /* start the scatter */ ierr = VecScatterBegin(xx,b->w,INSERT_VALUES,SCATTER_FORWARD,b->ctx);CHKERRQ(ierr); ierr = (*b->AIJ->ops->multadd)(b->AIJ,xx,yy,zz);CHKERRQ(ierr); ierr = VecScatterEnd(xx,b->w,INSERT_VALUES,SCATTER_FORWARD,b->ctx);CHKERRQ(ierr); ierr = (*b->OAIJ->ops->multadd)(b->OAIJ,b->w,yy,zz);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNC__ #define __FUNC__ /**/"MatMultTransposeAdd_MPIMAIJ_dof" int MatMultTransposeAdd_MPIMAIJ_dof(Mat A,Vec xx,Vec yy,Vec zz) { Mat_MPIMAIJ *b = (Mat_MPIMAIJ*)A->data; int ierr; PetscFunctionBegin; ierr = (*b->OAIJ->ops->multtranspose)(b->OAIJ,xx,b->w);CHKERRQ(ierr); ierr = VecScatterBegin(b->w,zz,ADD_VALUES,SCATTER_REVERSE,b->ctx);CHKERRQ(ierr); ierr = (*b->AIJ->ops->multtransposeadd)(b->AIJ,xx,yy,zz);CHKERRQ(ierr); ierr = VecScatterEnd(b->w,zz,ADD_VALUES,SCATTER_REVERSE,b->ctx);CHKERRQ(ierr); PetscFunctionReturn(0); } /* ---------------------------------------------------------------------------------- */ #undef __FUNC__ #define __FUNC__ /**/"MatCreateMAIJ" int MatCreateMAIJ(Mat A,int dof,Mat *maij) { int ierr,size,n; Mat_MPIMAIJ *b; Mat B; PetscFunctionBegin; ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); if (dof == 1) { *maij = A; } else { ierr = MATCreate(A->comm,dof*A->m,dof*A->n,dof*A->M,dof*A->N,&B);CHKERRQ(ierr); ierr = MatSetType(B,MATMAIJ);CHKERRQ(ierr); B->assembled = PETSC_TRUE; b = (Mat_MPIMAIJ*)B->data; b->dof = dof; ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); if (size == 1) { B->ops->destroy = MatDestroy_SeqMAIJ; b->AIJ = A; if (dof == 2) { B->ops->mult = MatMult_SeqMAIJ_2; B->ops->multadd = MatMultAdd_SeqMAIJ_2; B->ops->multtranspose = MatMultTranspose_SeqMAIJ_2; B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_2; } else if (dof == 3) { B->ops->mult = MatMult_SeqMAIJ_3; B->ops->multadd = MatMultAdd_SeqMAIJ_3; B->ops->multtranspose = MatMultTranspose_SeqMAIJ_3; B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_3; } else if (dof == 4) { B->ops->mult = MatMult_SeqMAIJ_4; B->ops->multadd = MatMultAdd_SeqMAIJ_4; B->ops->multtranspose = MatMultTranspose_SeqMAIJ_4; B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_4; } else if (dof == 5) { B->ops->mult = MatMult_SeqMAIJ_5; B->ops->multadd = MatMultAdd_SeqMAIJ_5; B->ops->multtranspose = MatMultTranspose_SeqMAIJ_5; B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_5; } else { SETERRQ1(1,1,"Cannot handle a dof of %d\n",dof); } } else { Mat_MPIAIJ *mpiaij = (Mat_MPIAIJ *)A->data; IS from,to; Vec gvec; int *garray,i; b->A = A; B->ops->destroy = MatDestroy_MPIMAIJ; ierr = MatCreateMAIJ(mpiaij->A,dof,&b->AIJ);CHKERRQ(ierr); ierr = MatCreateMAIJ(mpiaij->B,dof,&b->OAIJ);CHKERRQ(ierr); ierr = VecGetSize(mpiaij->lvec,&n);CHKERRQ(ierr); ierr = VecCreateSeq(PETSC_COMM_SELF,n*dof,&b->w);CHKERRQ(ierr); /* create two temporary Index sets for build scatter gather */ garray = (int*)PetscMalloc((n+1)*sizeof(int));CHKPTRQ(garray); for (i=0; igarray[i]; ierr = ISCreateBlock(A->comm,dof,n,garray,&from);CHKERRQ(ierr); ierr = PetscFree(garray);CHKERRQ(ierr); ierr = ISCreateStride(PETSC_COMM_SELF,n*dof,0,1,&to);CHKERRQ(ierr); /* create temporary global vector to generate scatter context */ ierr = VecCreateMPI(A->comm,dof*A->n,dof*A->N,&gvec);CHKERRQ(ierr); /* generate the scatter context */ ierr = VecScatterCreate(gvec,from,b->w,to,&b->ctx);CHKERRQ(ierr); ierr = ISDestroy(from);CHKERRQ(ierr); ierr = ISDestroy(to);CHKERRQ(ierr); ierr = VecDestroy(gvec);CHKERRQ(ierr); B->ops->mult = MatMult_MPIMAIJ_dof; B->ops->multtranspose = MatMultTranspose_MPIMAIJ_dof; B->ops->multadd = MatMultAdd_MPIMAIJ_dof; B->ops->multtransposeadd = MatMultTransposeAdd_MPIMAIJ_dof; } *maij = B; } PetscFunctionReturn(0); }