PETSc version 3.16.6
Fix/Edit manual page

MatH2OpusMapVec

Maps a vector between PETSc and H2Opus ordering

Synopsis


Input Parameters

A - the matrix
nativetopetsc - if true, maps from H2Opus ordering to PETSc ordering. If false, applies the reverse map
in - the vector to be mapped

Output Parameter

out - the newly created mapped vector

See Also

MatCreate(), MATH2OPUS, MatCreateH2OpusFromMat(), MatCreateH2OpusFromKernel()
*/ PetscErrorCode MatH2OpusMapVec(Mat A, PetscBool nativetopetsc, Vec in, Vec* out) { PetscBool ish2opus; Mat_H2OPUS *a = (Mat_H2OPUS*)A->data; PetscScalar *xin,*xout; PetscBool nm; PetscErrorCode ierr;

PetscFunctionBegin; PetscValidHeaderSpecific(A,MAT_CLASSID,1); PetscValidType(A,1); PetscValidLogicalCollectiveBool(A,nativetopetsc,2); PetscValidHeaderSpecific(in,VEC_CLASSID,3); PetscValidPointer(out,4); if (!A->assembled) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); ierr = PetscObjectTypeCompare((PetscObject)A,MATH2OPUS,&ish2opus);CHKERRQ(ierr); if (!ish2opus) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not for type %s",((PetscObject)A)->type_name); nm = a->nativemult; ierr = MatH2OpusSetNativeMult(A,(PetscBool)!nativetopetsc);CHKERRQ(ierr); ierr = MatCreateVecs(A,out,NULL);CHKERRQ(ierr); ierr = MatH2OpusSetNativeMult(A,nm);CHKERRQ(ierr); if (!a->sf) { /* same ordering */ ierr = VecCopy(in,*out);CHKERRQ(ierr); PetscFunctionReturn(0); } ierr = VecGetArrayRead(in,(const PetscScalar**)&xin);CHKERRQ(ierr); ierr = VecGetArrayWrite(*out,&xout);CHKERRQ(ierr); if (nativetopetsc) { ierr = PetscSFReduceBegin(a->sf,MPIU_SCALAR,xin,xout,MPI_REPLACE);CHKERRQ(ierr); ierr = PetscSFReduceEnd(a->sf,MPIU_SCALAR,xin,xout,MPI_REPLACE);CHKERRQ(ierr); } else { ierr = PetscSFBcastBegin(a->sf,MPIU_SCALAR,xin,xout,MPI_REPLACE);CHKERRQ(ierr); ierr = PetscSFBcastEnd(a->sf,MPIU_SCALAR,xin,xout,MPI_REPLACE);CHKERRQ(ierr); } ierr = VecRestoreArrayRead(in,(const PetscScalar**)&xin);CHKERRQ(ierr); ierr = VecRestoreArrayWrite(*out,&xout);CHKERRQ(ierr); PetscFunctionReturn(0); } #endif

Level

intermediate

Location

src/mat/impls/h2opus/math2opus.cu
Index of all Mat routines
Table of Contents for all manual pages
Index of all manual pages