LCOV - code coverage report
Current view: top level - sys/classes/bv/impls/mat - bvmat.c (source / functions) Hit Total Coverage
Test: SLEPc Lines: 326 332 98.2 %
Date: 2024-11-21 00:40:22 Functions: 24 24 100.0 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /*
       2             :    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
       3             :    SLEPc - Scalable Library for Eigenvalue Problem Computations
       4             :    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
       5             : 
       6             :    This file is part of SLEPc.
       7             :    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
       8             :    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
       9             : */
      10             : /*
      11             :    BV implemented with a dense Mat
      12             : */
      13             : 
      14             : #include <slepc/private/bvimpl.h>
      15             : #include "bvmat.h"
      16             : 
      17        6400 : static PetscErrorCode BVMult_Mat(BV Y,PetscScalar alpha,PetscScalar beta,BV X,Mat Q)
      18             : {
      19        6400 :   BV_MAT            *y = (BV_MAT*)Y->data,*x = (BV_MAT*)X->data;
      20        6400 :   PetscScalar       *py;
      21        6400 :   const PetscScalar *px,*q;
      22        6400 :   PetscInt          ldq;
      23             : 
      24        6400 :   PetscFunctionBegin;
      25        6400 :   PetscCall(MatDenseGetArrayRead(x->A,&px));
      26        6400 :   PetscCall(MatDenseGetArray(y->A,&py));
      27        6400 :   if (Q) {
      28        4948 :     PetscCall(MatDenseGetLDA(Q,&ldq));
      29        4948 :     PetscCall(MatDenseGetArrayRead(Q,&q));
      30        4948 :     PetscCall(BVMult_BLAS_Private(Y,Y->n,Y->k-Y->l,X->k-X->l,alpha,px+(X->nc+X->l)*X->ld,X->ld,q+Y->l*ldq+X->l,ldq,beta,py+(Y->nc+Y->l)*Y->ld,Y->ld));
      31        4948 :     PetscCall(MatDenseRestoreArrayRead(Q,&q));
      32        1452 :   } else PetscCall(BVAXPY_BLAS_Private(Y,Y->n,Y->k-Y->l,alpha,px+(X->nc+X->l)*X->ld,X->ld,beta,py+(Y->nc+Y->l)*Y->ld,Y->ld));
      33        6400 :   PetscCall(MatDenseRestoreArrayRead(x->A,&px));
      34        6400 :   PetscCall(MatDenseRestoreArray(y->A,&py));
      35        6400 :   PetscFunctionReturn(PETSC_SUCCESS);
      36             : }
      37             : 
      38      719194 : static PetscErrorCode BVMultVec_Mat(BV X,PetscScalar alpha,PetscScalar beta,Vec y,PetscScalar *q)
      39             : {
      40      719194 :   BV_MAT            *x = (BV_MAT*)X->data;
      41      719194 :   PetscScalar       *py,*qq=q;
      42      719194 :   const PetscScalar *px;
      43             : 
      44      719194 :   PetscFunctionBegin;
      45      719194 :   PetscCall(MatDenseGetArrayRead(x->A,&px));
      46      719194 :   PetscCall(VecGetArray(y,&py));
      47      719194 :   if (!q) PetscCall(VecGetArray(X->buffer,&qq));
      48      719194 :   PetscCall(BVMultVec_BLAS_Private(X,X->n,X->k-X->l,alpha,px+(X->nc+X->l)*X->ld,X->ld,qq,beta,py));
      49      719194 :   if (!q) PetscCall(VecRestoreArray(X->buffer,&qq));
      50      719194 :   PetscCall(MatDenseRestoreArrayRead(x->A,&px));
      51      719194 :   PetscCall(VecRestoreArray(y,&py));
      52      719194 :   PetscFunctionReturn(PETSC_SUCCESS);
      53             : }
      54             : 
      55       24216 : static PetscErrorCode BVMultInPlace_Mat(BV V,Mat Q,PetscInt s,PetscInt e)
      56             : {
      57       24216 :   BV_MAT            *ctx = (BV_MAT*)V->data;
      58       24216 :   PetscScalar       *pv;
      59       24216 :   const PetscScalar *q;
      60       24216 :   PetscInt          ldq;
      61             : 
      62       24216 :   PetscFunctionBegin;
      63       24216 :   if (s>=e || !V->n) PetscFunctionReturn(PETSC_SUCCESS);
      64       22376 :   PetscCall(MatDenseGetLDA(Q,&ldq));
      65       22376 :   PetscCall(MatDenseGetArray(ctx->A,&pv));
      66       22376 :   PetscCall(MatDenseGetArrayRead(Q,&q));
      67       22376 :   PetscCall(BVMultInPlace_BLAS_Private(V,V->n,V->k-V->l,s-V->l,e-V->l,pv+(V->nc+V->l)*V->ld,V->ld,q+V->l*ldq+V->l,ldq,PETSC_FALSE));
      68       22376 :   PetscCall(MatDenseRestoreArrayRead(Q,&q));
      69       22376 :   PetscCall(MatDenseRestoreArray(ctx->A,&pv));
      70       22376 :   PetscFunctionReturn(PETSC_SUCCESS);
      71             : }
      72             : 
      73           1 : static PetscErrorCode BVMultInPlaceHermitianTranspose_Mat(BV V,Mat Q,PetscInt s,PetscInt e)
      74             : {
      75           1 :   BV_MAT            *ctx = (BV_MAT*)V->data;
      76           1 :   PetscScalar       *pv;
      77           1 :   const PetscScalar *q;
      78           1 :   PetscInt          ldq;
      79             : 
      80           1 :   PetscFunctionBegin;
      81           1 :   if (s>=e || !V->n) PetscFunctionReturn(PETSC_SUCCESS);
      82           1 :   PetscCall(MatDenseGetLDA(Q,&ldq));
      83           1 :   PetscCall(MatDenseGetArray(ctx->A,&pv));
      84           1 :   PetscCall(MatDenseGetArrayRead(Q,&q));
      85           1 :   PetscCall(BVMultInPlace_BLAS_Private(V,V->n,V->k-V->l,s-V->l,e-V->l,pv+(V->nc+V->l)*V->ld,V->ld,q+V->l*ldq+V->l,ldq,PETSC_TRUE));
      86           1 :   PetscCall(MatDenseRestoreArrayRead(Q,&q));
      87           1 :   PetscCall(MatDenseRestoreArray(ctx->A,&pv));
      88           1 :   PetscFunctionReturn(PETSC_SUCCESS);
      89             : }
      90             : 
      91       29604 : static PetscErrorCode BVDot_Mat(BV X,BV Y,Mat M)
      92             : {
      93       29604 :   BV_MAT            *x = (BV_MAT*)X->data,*y = (BV_MAT*)Y->data;
      94       29604 :   PetscScalar       *m;
      95       29604 :   const PetscScalar *px,*py;
      96       29604 :   PetscInt          ldm;
      97             : 
      98       29604 :   PetscFunctionBegin;
      99       29604 :   PetscCall(MatDenseGetLDA(M,&ldm));
     100       29604 :   PetscCall(MatDenseGetArrayRead(x->A,&px));
     101       29604 :   PetscCall(MatDenseGetArrayRead(y->A,&py));
     102       29604 :   PetscCall(MatDenseGetArray(M,&m));
     103       29604 :   PetscCall(BVDot_BLAS_Private(X,Y->k-Y->l,X->k-X->l,X->n,py+(Y->nc+Y->l)*Y->ld,Y->ld,px+(X->nc+X->l)*X->ld,X->ld,m+X->l*ldm+Y->l,ldm,x->mpi));
     104       29604 :   PetscCall(MatDenseRestoreArray(M,&m));
     105       29604 :   PetscCall(MatDenseRestoreArrayRead(x->A,&px));
     106       29604 :   PetscCall(MatDenseRestoreArrayRead(y->A,&py));
     107       29604 :   PetscFunctionReturn(PETSC_SUCCESS);
     108             : }
     109             : 
     110      415698 : static PetscErrorCode BVDotVec_Mat(BV X,Vec y,PetscScalar *q)
     111             : {
     112      415698 :   BV_MAT            *x = (BV_MAT*)X->data;
     113      415698 :   PetscScalar       *qq=q;
     114      415698 :   const PetscScalar *px,*py;
     115      415698 :   Vec               z = y;
     116             : 
     117      415698 :   PetscFunctionBegin;
     118      415698 :   if (PetscUnlikely(X->matrix)) {
     119       47463 :     PetscCall(BV_IPMatMult(X,y));
     120       47463 :     z = X->Bx;
     121             :   }
     122      415698 :   PetscCall(MatDenseGetArrayRead(x->A,&px));
     123      415698 :   PetscCall(VecGetArrayRead(z,&py));
     124      415698 :   if (!q) PetscCall(VecGetArray(X->buffer,&qq));
     125      415698 :   PetscCall(BVDotVec_BLAS_Private(X,X->n,X->k-X->l,px+(X->nc+X->l)*X->ld,X->ld,py,qq,x->mpi));
     126      415698 :   if (!q) PetscCall(VecRestoreArray(X->buffer,&qq));
     127      415698 :   PetscCall(VecRestoreArrayRead(z,&py));
     128      415698 :   PetscCall(MatDenseRestoreArrayRead(x->A,&px));
     129      415698 :   PetscFunctionReturn(PETSC_SUCCESS);
     130             : }
     131             : 
     132        3642 : static PetscErrorCode BVDotVec_Local_Mat(BV X,Vec y,PetscScalar *m)
     133             : {
     134        3642 :   BV_MAT            *x = (BV_MAT*)X->data;
     135        3642 :   const PetscScalar *px,*py;
     136        3642 :   Vec               z = y;
     137             : 
     138        3642 :   PetscFunctionBegin;
     139        3642 :   if (PetscUnlikely(X->matrix)) {
     140           0 :     PetscCall(BV_IPMatMult(X,y));
     141           0 :     z = X->Bx;
     142             :   }
     143        3642 :   PetscCall(MatDenseGetArrayRead(x->A,&px));
     144        3642 :   PetscCall(VecGetArrayRead(z,&py));
     145        3642 :   PetscCall(BVDotVec_BLAS_Private(X,X->n,X->k-X->l,px+(X->nc+X->l)*X->ld,X->ld,py,m,PETSC_FALSE));
     146        3642 :   PetscCall(VecRestoreArrayRead(z,&py));
     147        3642 :   PetscCall(MatDenseRestoreArrayRead(x->A,&px));
     148        3642 :   PetscFunctionReturn(PETSC_SUCCESS);
     149             : }
     150             : 
     151      179187 : static PetscErrorCode BVScale_Mat(BV bv,PetscInt j,PetscScalar alpha)
     152             : {
     153      179187 :   BV_MAT         *ctx = (BV_MAT*)bv->data;
     154      179187 :   PetscScalar    *array;
     155             : 
     156      179187 :   PetscFunctionBegin;
     157      179187 :   if (!bv->n) PetscFunctionReturn(PETSC_SUCCESS);
     158      179187 :   PetscCall(MatDenseGetArray(ctx->A,&array));
     159      179187 :   if (PetscUnlikely(j<0)) PetscCall(BVScale_BLAS_Private(bv,(bv->k-bv->l)*bv->ld,array+(bv->nc+bv->l)*bv->ld,alpha));
     160      178789 :   else PetscCall(BVScale_BLAS_Private(bv,bv->n,array+(bv->nc+j)*bv->ld,alpha));
     161      179187 :   PetscCall(MatDenseRestoreArray(ctx->A,&array));
     162      179187 :   PetscFunctionReturn(PETSC_SUCCESS);
     163             : }
     164             : 
     165        8279 : static PetscErrorCode BVNorm_Mat(BV bv,PetscInt j,NormType type,PetscReal *val)
     166             : {
     167        8279 :   BV_MAT            *ctx = (BV_MAT*)bv->data;
     168        8279 :   const PetscScalar *array;
     169             : 
     170        8279 :   PetscFunctionBegin;
     171        8279 :   PetscCall(MatDenseGetArrayRead(ctx->A,&array));
     172        8279 :   if (PetscUnlikely(j<0)) PetscCall(BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->ld,bv->ld,type,val,ctx->mpi));
     173        8209 :   else PetscCall(BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->ld,bv->ld,type,val,ctx->mpi));
     174        8279 :   PetscCall(MatDenseRestoreArrayRead(ctx->A,&array));
     175        8279 :   PetscFunctionReturn(PETSC_SUCCESS);
     176             : }
     177             : 
     178        2543 : static PetscErrorCode BVNorm_Local_Mat(BV bv,PetscInt j,NormType type,PetscReal *val)
     179             : {
     180        2543 :   BV_MAT            *ctx = (BV_MAT*)bv->data;
     181        2543 :   const PetscScalar *array;
     182             : 
     183        2543 :   PetscFunctionBegin;
     184        2543 :   PetscCall(MatDenseGetArrayRead(ctx->A,&array));
     185        2543 :   if (PetscUnlikely(j<0)) PetscCall(BVNorm_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->ld,bv->ld,type,val,PETSC_FALSE));
     186        2543 :   else PetscCall(BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->ld,bv->ld,type,val,PETSC_FALSE));
     187        2543 :   PetscCall(MatDenseRestoreArrayRead(ctx->A,&array));
     188        2543 :   PetscFunctionReturn(PETSC_SUCCESS);
     189             : }
     190             : 
     191        3882 : static PetscErrorCode BVNormalize_Mat(BV bv,PetscScalar *eigi)
     192             : {
     193        3882 :   BV_MAT         *ctx = (BV_MAT*)bv->data;
     194        3882 :   PetscScalar    *array,*wi=NULL;
     195             : 
     196        3882 :   PetscFunctionBegin;
     197        3882 :   PetscCall(MatDenseGetArray(ctx->A,&array));
     198        3882 :   if (eigi) wi = eigi+bv->l;
     199        3882 :   PetscCall(BVNormalize_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->ld,bv->ld,wi,ctx->mpi));
     200        3882 :   PetscCall(MatDenseRestoreArray(ctx->A,&array));
     201        3882 :   PetscFunctionReturn(PETSC_SUCCESS);
     202             : }
     203             : 
     204       10783 : static PetscErrorCode BVMatMult_Mat(BV V,Mat A,BV W)
     205             : {
     206       10783 :   PetscInt       j;
     207       10783 :   Mat            Vmat,Wmat;
     208       10783 :   Vec            vv,ww;
     209             : 
     210       10783 :   PetscFunctionBegin;
     211       10783 :   if (V->vmm) {
     212       10780 :     PetscCall(BVGetMat(V,&Vmat));
     213       10780 :     PetscCall(BVGetMat(W,&Wmat));
     214       10780 :     PetscCall(MatProductCreateWithMat(A,Vmat,NULL,Wmat));
     215       10780 :     PetscCall(MatProductSetType(Wmat,MATPRODUCT_AB));
     216       10780 :     PetscCall(MatProductSetFromOptions(Wmat));
     217       10780 :     PetscCall(MatProductSymbolic(Wmat));
     218       10780 :     PetscCall(MatProductNumeric(Wmat));
     219       10780 :     PetscCall(MatProductClear(Wmat));
     220       10780 :     PetscCall(BVRestoreMat(V,&Vmat));
     221       10780 :     PetscCall(BVRestoreMat(W,&Wmat));
     222             :   } else {
     223          26 :     for (j=0;j<V->k-V->l;j++) {
     224          23 :       PetscCall(BVGetColumn(V,V->l+j,&vv));
     225          23 :       PetscCall(BVGetColumn(W,W->l+j,&ww));
     226          23 :       PetscCall(MatMult(A,vv,ww));
     227          23 :       PetscCall(BVRestoreColumn(V,V->l+j,&vv));
     228          23 :       PetscCall(BVRestoreColumn(W,W->l+j,&ww));
     229             :     }
     230             :   }
     231       10783 :   PetscFunctionReturn(PETSC_SUCCESS);
     232             : }
     233             : 
     234       12816 : static PetscErrorCode BVCopy_Mat(BV V,BV W)
     235             : {
     236       12816 :   BV_MAT            *v = (BV_MAT*)V->data,*w = (BV_MAT*)W->data;
     237       12816 :   const PetscScalar *pv;
     238       12816 :   PetscScalar       *pw;
     239       12816 :   PetscInt          j;
     240             : 
     241       12816 :   PetscFunctionBegin;
     242       12816 :   PetscCall(MatDenseGetArrayRead(v->A,&pv));
     243       12816 :   PetscCall(MatDenseGetArray(w->A,&pw));
     244       82493 :   for (j=0;j<V->k-V->l;j++) PetscCall(PetscArraycpy(pw+(W->nc+W->l+j)*W->ld,pv+(V->nc+V->l+j)*V->ld,V->n));
     245       12816 :   PetscCall(MatDenseRestoreArrayRead(v->A,&pv));
     246       12816 :   PetscCall(MatDenseRestoreArray(w->A,&pw));
     247       12816 :   PetscFunctionReturn(PETSC_SUCCESS);
     248             : }
     249             : 
     250        6065 : static PetscErrorCode BVCopyColumn_Mat(BV V,PetscInt j,PetscInt i)
     251             : {
     252        6065 :   BV_MAT         *v = (BV_MAT*)V->data;
     253        6065 :   PetscScalar    *pv;
     254             : 
     255        6065 :   PetscFunctionBegin;
     256        6065 :   PetscCall(MatDenseGetArray(v->A,&pv));
     257        6065 :   PetscCall(PetscArraycpy(pv+(V->nc+i)*V->ld,pv+(V->nc+j)*V->ld,V->n));
     258        6065 :   PetscCall(MatDenseRestoreArray(v->A,&pv));
     259        6065 :   PetscFunctionReturn(PETSC_SUCCESS);
     260             : }
     261             : 
     262         164 : static PetscErrorCode BVResize_Mat(BV bv,PetscInt m,PetscBool copy)
     263             : {
     264         164 :   BV_MAT            *ctx = (BV_MAT*)bv->data;
     265         164 :   Mat               A,Msrc,Mdst;
     266         164 :   char              str[50];
     267             : 
     268         164 :   PetscFunctionBegin;
     269         164 :   PetscCall(MatCreateDenseFromVecType(PetscObjectComm((PetscObject)bv),bv->vtype,bv->n,PETSC_DECIDE,bv->N,m,bv->ld,NULL,&A));
     270         164 :   if (((PetscObject)bv)->name) {
     271          18 :     PetscCall(PetscSNPrintf(str,sizeof(str),"%s_0",((PetscObject)bv)->name));
     272          18 :     PetscCall(PetscObjectSetName((PetscObject)A,str));
     273             :   }
     274         164 :   if (copy) {
     275          19 :     PetscCall(MatDenseGetSubMatrix(ctx->A,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,PetscMin(m,bv->m),&Msrc));
     276          19 :     PetscCall(MatDenseGetSubMatrix(A,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,PetscMin(m,bv->m),&Mdst));
     277          19 :     PetscCall(MatCopy(Msrc,Mdst,SAME_NONZERO_PATTERN));
     278          19 :     PetscCall(MatDenseRestoreSubMatrix(ctx->A,&Msrc));
     279          19 :     PetscCall(MatDenseRestoreSubMatrix(A,&Mdst));
     280             :   }
     281         164 :   PetscCall(MatDestroy(&ctx->A));
     282         164 :   ctx->A = A;
     283         164 :   PetscFunctionReturn(PETSC_SUCCESS);
     284             : }
     285             : 
     286     1200414 : static PetscErrorCode BVGetColumn_Mat(BV bv,PetscInt j,Vec *v)
     287             : {
     288     1200414 :   BV_MAT         *ctx = (BV_MAT*)bv->data;
     289     1200414 :   PetscScalar    *pA;
     290     1200414 :   PetscInt       l;
     291             : 
     292     1200414 :   PetscFunctionBegin;
     293     1200414 :   l = BVAvailableVec;
     294     1200414 :   PetscCall(MatDenseGetArray(ctx->A,&pA));
     295     1200414 :   PetscCall(VecPlaceArray(bv->cv[l],pA+(bv->nc+j)*bv->ld));
     296     1200414 :   PetscFunctionReturn(PETSC_SUCCESS);
     297             : }
     298             : 
     299     1200414 : static PetscErrorCode BVRestoreColumn_Mat(BV bv,PetscInt j,Vec *v)
     300             : {
     301     1200414 :   BV_MAT         *ctx = (BV_MAT*)bv->data;
     302     1200414 :   PetscScalar    *pA;
     303     1200414 :   PetscInt       l;
     304             : 
     305     1200414 :   PetscFunctionBegin;
     306     1200414 :   l = (j==bv->ci[0])? 0: 1;
     307     1200414 :   PetscCall(VecResetArray(bv->cv[l]));
     308     1200414 :   PetscCall(MatDenseRestoreArray(ctx->A,&pA));
     309     1200414 :   PetscFunctionReturn(PETSC_SUCCESS);
     310             : }
     311             : 
     312       25975 : static PetscErrorCode BVGetArray_Mat(BV bv,PetscScalar **a)
     313             : {
     314       25975 :   BV_MAT         *ctx = (BV_MAT*)bv->data;
     315             : 
     316       25975 :   PetscFunctionBegin;
     317       25975 :   PetscCall(MatDenseGetArray(ctx->A,a));
     318       25975 :   PetscFunctionReturn(PETSC_SUCCESS);
     319             : }
     320             : 
     321       25975 : static PetscErrorCode BVRestoreArray_Mat(BV bv,PetscScalar **a)
     322             : {
     323       25975 :   BV_MAT         *ctx = (BV_MAT*)bv->data;
     324             : 
     325       25975 :   PetscFunctionBegin;
     326       25975 :   if (a) PetscCall(MatDenseRestoreArray(ctx->A,a));
     327       25975 :   PetscFunctionReturn(PETSC_SUCCESS);
     328             : }
     329             : 
     330         160 : static PetscErrorCode BVGetArrayRead_Mat(BV bv,const PetscScalar **a)
     331             : {
     332         160 :   BV_MAT         *ctx = (BV_MAT*)bv->data;
     333             : 
     334         160 :   PetscFunctionBegin;
     335         160 :   PetscCall(MatDenseGetArrayRead(ctx->A,a));
     336         160 :   PetscFunctionReturn(PETSC_SUCCESS);
     337             : }
     338             : 
     339         160 : static PetscErrorCode BVRestoreArrayRead_Mat(BV bv,const PetscScalar **a)
     340             : {
     341         160 :   BV_MAT         *ctx = (BV_MAT*)bv->data;
     342             : 
     343         160 :   PetscFunctionBegin;
     344         160 :   if (a) PetscCall(MatDenseRestoreArrayRead(ctx->A,a));
     345         160 :   PetscFunctionReturn(PETSC_SUCCESS);
     346             : }
     347             : 
     348          28 : static PetscErrorCode BVView_Mat(BV bv,PetscViewer viewer)
     349             : {
     350          28 :   Mat               A;
     351          28 :   PetscViewerFormat format;
     352          28 :   PetscBool         isascii;
     353          28 :   const char        *bvname,*name;
     354             : 
     355          28 :   PetscFunctionBegin;
     356          28 :   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
     357          28 :   if (isascii) {
     358          28 :     PetscCall(PetscViewerGetFormat(viewer,&format));
     359          28 :     if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) PetscFunctionReturn(PETSC_SUCCESS);
     360             :   }
     361          10 :   PetscCall(BVGetMat(bv,&A));
     362          10 :   PetscCall(MatView(A,viewer));
     363          10 :   if (format == PETSC_VIEWER_ASCII_MATLAB) {
     364           0 :     PetscCall(PetscObjectGetName((PetscObject)A,&name));
     365           0 :     PetscCall(PetscObjectGetName((PetscObject)bv,&bvname));
     366           0 :     PetscCall(PetscViewerASCIIPrintf(viewer,"%s=%s;clear %s\n",bvname,name,name));
     367           0 :     if (bv->nc) PetscCall(PetscViewerASCIIPrintf(viewer,"%s=%s(:,%" PetscInt_FMT ":end);\n",bvname,bvname,bv->nc+1));
     368             :   }
     369          10 :   PetscCall(BVRestoreMat(bv,&A));
     370          10 :   PetscFunctionReturn(PETSC_SUCCESS);
     371             : }
     372             : 
     373       14049 : static PetscErrorCode BVDestroy_Mat(BV bv)
     374             : {
     375       14049 :   BV_MAT         *ctx = (BV_MAT*)bv->data;
     376             : 
     377       14049 :   PetscFunctionBegin;
     378       14049 :   PetscCall(MatDestroy(&ctx->A));
     379       14049 :   PetscCall(VecDestroy(&bv->cv[0]));
     380       14049 :   PetscCall(VecDestroy(&bv->cv[1]));
     381       14049 :   PetscCall(PetscFree(bv->data));
     382       14049 :   PetscFunctionReturn(PETSC_SUCCESS);
     383             : }
     384             : 
     385       14049 : SLEPC_EXTERN PetscErrorCode BVCreate_Mat(BV bv)
     386             : {
     387       14049 :   BV_MAT         *ctx;
     388       14049 :   PetscInt       nloc,lsplit;
     389       14049 :   PetscBool      seq;
     390       14049 :   char           str[50];
     391       14049 :   PetscScalar    *array,*ptr=NULL;
     392       14049 :   BV             parent;
     393       14049 :   Mat            Apar;
     394             : 
     395       14049 :   PetscFunctionBegin;
     396       14049 :   PetscCall(PetscNew(&ctx));
     397       14049 :   bv->data = (void*)ctx;
     398             : 
     399       14049 :   PetscCall(PetscStrcmpAny(bv->vtype,&bv->cuda,VECSEQCUDA,VECMPICUDA,""));
     400       14049 :   PetscCall(PetscStrcmpAny(bv->vtype,&bv->hip,VECSEQHIP,VECMPIHIP,""));
     401       14049 :   PetscCall(PetscStrcmpAny(bv->vtype,&ctx->mpi,VECMPI,VECMPICUDA,VECMPIHIP,""));
     402             : 
     403       14049 :   PetscCall(PetscStrcmp(bv->vtype,VECSEQ,&seq));
     404       14049 :   PetscCheck(seq || ctx->mpi || bv->cuda || bv->hip,PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"BVMAT does not support the requested vector type: %s",bv->vtype);
     405             : 
     406       14049 :   PetscCall(PetscLayoutGetLocalSize(bv->map,&nloc));
     407       14049 :   PetscCall(BV_SetDefaultLD(bv,nloc));
     408             : 
     409       14049 :   if (PetscUnlikely(bv->issplit)) {
     410             :     /* split BV: share the memory of the parent BV */
     411         120 :     parent = bv->splitparent;
     412         120 :     lsplit = parent->lsplit;
     413         120 :     Apar = ((BV_MAT*)parent->data)->A;
     414         120 :     if (bv->cuda) {
     415             : #if defined(PETSC_HAVE_CUDA)
     416             :       PetscCall(MatDenseCUDAGetArray(Apar,&array));
     417             :       if (bv->issplit>0) ptr = (bv->issplit==1)? array: array+lsplit*bv->ld;
     418             :       else ptr = (bv->issplit==1)? array: array-lsplit;
     419             :       PetscCall(MatDenseCUDARestoreArray(Apar,&array));
     420             : #endif
     421         120 :     } else if (bv->hip) {
     422             : #if defined(PETSC_HAVE_HIP)
     423             :       PetscCall(MatDenseHIPGetArray(Apar,&array));
     424             :       if (bv->issplit>0) ptr = (bv->issplit==1)? array: array+lsplit*bv->ld;
     425             :       else ptr = (bv->issplit==1)? array: array-lsplit;
     426             :       PetscCall(MatDenseHIPRestoreArray(Apar,&array));
     427             : #endif
     428             :     } else {
     429         120 :       PetscCall(MatDenseGetArray(Apar,&array));
     430         120 :       if (bv->issplit>0) ptr = (bv->issplit==1)? array: array+lsplit*bv->ld;
     431          40 :       else ptr = (bv->issplit==-1)? array: array-lsplit;
     432         120 :       PetscCall(MatDenseRestoreArray(Apar,&array));
     433             :     }
     434             :   }
     435             : 
     436       14049 :   PetscCall(MatCreateDenseFromVecType(PetscObjectComm((PetscObject)bv),bv->vtype,nloc,PETSC_DECIDE,bv->N,bv->m,bv->ld,ptr,&ctx->A));
     437       14049 :   if (((PetscObject)bv)->name) {
     438         118 :     PetscCall(PetscSNPrintf(str,sizeof(str),"%s_0",((PetscObject)bv)->name));
     439         118 :     PetscCall(PetscObjectSetName((PetscObject)ctx->A,str));
     440             :   }
     441             : 
     442       14049 :   if (PetscUnlikely(bv->Acreate)) {
     443         144 :     PetscCall(MatConvert(bv->Acreate,bv->cuda?MATDENSECUDA:bv->hip?MATDENSEHIP:MATDENSE,MAT_REUSE_MATRIX,&ctx->A));
     444          72 :     PetscCall(MatDestroy(&bv->Acreate));
     445             :   }
     446             : 
     447       14049 :   PetscCall(BVCreateVecEmpty(bv,&bv->cv[0]));
     448       14049 :   PetscCall(BVCreateVecEmpty(bv,&bv->cv[1]));
     449             : 
     450       14049 :   if (bv->cuda) {
     451             : #if defined(PETSC_HAVE_CUDA)
     452             :     bv->ops->mult             = BVMult_Mat_CUDA;
     453             :     bv->ops->multvec          = BVMultVec_Mat_CUDA;
     454             :     bv->ops->multinplace      = BVMultInPlace_Mat_CUDA;
     455             :     bv->ops->multinplacetrans = BVMultInPlaceHermitianTranspose_Mat_CUDA;
     456             :     bv->ops->dot              = BVDot_Mat_CUDA;
     457             :     bv->ops->dotvec           = BVDotVec_Mat_CUDA;
     458             :     bv->ops->dotvec_local     = BVDotVec_Local_Mat_CUDA;
     459             :     bv->ops->scale            = BVScale_Mat_CUDA;
     460             :     bv->ops->norm             = BVNorm_Mat_CUDA;
     461             :     bv->ops->norm_local       = BVNorm_Local_Mat_CUDA;
     462             :     bv->ops->normalize        = BVNormalize_Mat_CUDA;
     463             :     bv->ops->matmult          = BVMatMult_Mat_CUDA;
     464             :     bv->ops->copy             = BVCopy_Mat_CUDA;
     465             :     bv->ops->copycolumn       = BVCopyColumn_Mat_CUDA;
     466             :     bv->ops->getcolumn        = BVGetColumn_Mat_CUDA;
     467             :     bv->ops->restorecolumn    = BVRestoreColumn_Mat_CUDA;
     468             :     bv->ops->restoresplit     = BVRestoreSplit_Mat_CUDA;
     469             :     bv->ops->restoresplitrows = BVRestoreSplitRows_Mat_CUDA;
     470             :     bv->ops->getmat           = BVGetMat_Mat_CUDA;
     471             :     bv->ops->restoremat       = BVRestoreMat_Mat_CUDA;
     472             : #endif
     473       14049 :   } else if (bv->hip) {
     474             : #if defined(PETSC_HAVE_HIP)
     475             :     bv->ops->mult             = BVMult_Mat_HIP;
     476             :     bv->ops->multvec          = BVMultVec_Mat_HIP;
     477             :     bv->ops->multinplace      = BVMultInPlace_Mat_HIP;
     478             :     bv->ops->multinplacetrans = BVMultInPlaceHermitianTranspose_Mat_HIP;
     479             :     bv->ops->dot              = BVDot_Mat_HIP;
     480             :     bv->ops->dotvec           = BVDotVec_Mat_HIP;
     481             :     bv->ops->dotvec_local     = BVDotVec_Local_Mat_HIP;
     482             :     bv->ops->scale            = BVScale_Mat_HIP;
     483             :     bv->ops->norm             = BVNorm_Mat_HIP;
     484             :     bv->ops->norm_local       = BVNorm_Local_Mat_HIP;
     485             :     bv->ops->normalize        = BVNormalize_Mat_HIP;
     486             :     bv->ops->matmult          = BVMatMult_Mat_HIP;
     487             :     bv->ops->copy             = BVCopy_Mat_HIP;
     488             :     bv->ops->copycolumn       = BVCopyColumn_Mat_HIP;
     489             :     bv->ops->getcolumn        = BVGetColumn_Mat_HIP;
     490             :     bv->ops->restorecolumn    = BVRestoreColumn_Mat_HIP;
     491             :     bv->ops->restoresplit     = BVRestoreSplit_Mat_HIP;
     492             :     bv->ops->restoresplitrows = BVRestoreSplitRows_Mat_HIP;
     493             :     bv->ops->getmat           = BVGetMat_Mat_HIP;
     494             :     bv->ops->restoremat       = BVRestoreMat_Mat_HIP;
     495             : #endif
     496             :   } else {
     497       14049 :     bv->ops->mult             = BVMult_Mat;
     498       14049 :     bv->ops->multvec          = BVMultVec_Mat;
     499       14049 :     bv->ops->multinplace      = BVMultInPlace_Mat;
     500       14049 :     bv->ops->multinplacetrans = BVMultInPlaceHermitianTranspose_Mat;
     501       14049 :     bv->ops->dot              = BVDot_Mat;
     502       14049 :     bv->ops->dotvec           = BVDotVec_Mat;
     503       14049 :     bv->ops->dotvec_local     = BVDotVec_Local_Mat;
     504       14049 :     bv->ops->scale            = BVScale_Mat;
     505       14049 :     bv->ops->norm             = BVNorm_Mat;
     506       14049 :     bv->ops->norm_local       = BVNorm_Local_Mat;
     507       14049 :     bv->ops->normalize        = BVNormalize_Mat;
     508       14049 :     bv->ops->matmult          = BVMatMult_Mat;
     509       14049 :     bv->ops->copy             = BVCopy_Mat;
     510       14049 :     bv->ops->copycolumn       = BVCopyColumn_Mat;
     511       14049 :     bv->ops->getcolumn        = BVGetColumn_Mat;
     512       14049 :     bv->ops->restorecolumn    = BVRestoreColumn_Mat;
     513       14049 :     bv->ops->getmat           = BVGetMat_Default;
     514       14049 :     bv->ops->restoremat       = BVRestoreMat_Default;
     515             :   }
     516       14049 :   bv->ops->resize           = BVResize_Mat;
     517       14049 :   bv->ops->getarray         = BVGetArray_Mat;
     518       14049 :   bv->ops->restorearray     = BVRestoreArray_Mat;
     519       14049 :   bv->ops->getarrayread     = BVGetArrayRead_Mat;
     520       14049 :   bv->ops->restorearrayread = BVRestoreArrayRead_Mat;
     521       14049 :   bv->ops->destroy          = BVDestroy_Mat;
     522       14049 :   bv->ops->view             = BVView_Mat;
     523       14049 :   PetscFunctionReturn(PETSC_SUCCESS);
     524             : }

Generated by: LCOV version 1.14