LCOV - code coverage report
Current view: top level - sys/classes/bv/impls/mat - bvmat.c (source / functions) Hit Total Coverage
Test: SLEPc Lines: 320 327 97.9 %
Date: 2024-04-25 00:48:42 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        6522 : static PetscErrorCode BVMult_Mat(BV Y,PetscScalar alpha,PetscScalar beta,BV X,Mat Q)
      18             : {
      19        6522 :   BV_MAT            *y = (BV_MAT*)Y->data,*x = (BV_MAT*)X->data;
      20        6522 :   PetscScalar       *py;
      21        6522 :   const PetscScalar *px,*q;
      22        6522 :   PetscInt          ldq;
      23             : 
      24        6522 :   PetscFunctionBegin;
      25        6522 :   PetscCall(MatDenseGetArrayRead(x->A,&px));
      26        6522 :   PetscCall(MatDenseGetArray(y->A,&py));
      27        6522 :   if (Q) {
      28        5439 :     PetscCall(MatDenseGetLDA(Q,&ldq));
      29        5439 :     PetscCall(MatDenseGetArrayRead(Q,&q));
      30        5439 :     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        5439 :     PetscCall(MatDenseRestoreArrayRead(Q,&q));
      32        1083 :   } 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        6522 :   PetscCall(MatDenseRestoreArrayRead(x->A,&px));
      34        6522 :   PetscCall(MatDenseRestoreArray(y->A,&py));
      35        6522 :   PetscFunctionReturn(PETSC_SUCCESS);
      36             : }
      37             : 
      38      609383 : static PetscErrorCode BVMultVec_Mat(BV X,PetscScalar alpha,PetscScalar beta,Vec y,PetscScalar *q)
      39             : {
      40      609383 :   BV_MAT            *x = (BV_MAT*)X->data;
      41      609383 :   PetscScalar       *py,*qq=q;
      42      609383 :   const PetscScalar *px;
      43             : 
      44      609383 :   PetscFunctionBegin;
      45      609383 :   PetscCall(MatDenseGetArrayRead(x->A,&px));
      46      609383 :   PetscCall(VecGetArray(y,&py));
      47      609383 :   if (!q) PetscCall(VecGetArray(X->buffer,&qq));
      48      609383 :   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      609383 :   if (!q) PetscCall(VecRestoreArray(X->buffer,&qq));
      50      609383 :   PetscCall(MatDenseRestoreArrayRead(x->A,&px));
      51      609383 :   PetscCall(VecRestoreArray(y,&py));
      52      609383 :   PetscFunctionReturn(PETSC_SUCCESS);
      53             : }
      54             : 
      55       27386 : static PetscErrorCode BVMultInPlace_Mat(BV V,Mat Q,PetscInt s,PetscInt e)
      56             : {
      57       27386 :   BV_MAT            *ctx = (BV_MAT*)V->data;
      58       27386 :   PetscScalar       *pv;
      59       27386 :   const PetscScalar *q;
      60       27386 :   PetscInt          ldq;
      61             : 
      62       27386 :   PetscFunctionBegin;
      63       27386 :   if (s>=e || !V->n) PetscFunctionReturn(PETSC_SUCCESS);
      64       25573 :   PetscCall(MatDenseGetLDA(Q,&ldq));
      65       25573 :   PetscCall(MatDenseGetArray(ctx->A,&pv));
      66       25573 :   PetscCall(MatDenseGetArrayRead(Q,&q));
      67       25573 :   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       25573 :   PetscCall(MatDenseRestoreArrayRead(Q,&q));
      69       25573 :   PetscCall(MatDenseRestoreArray(ctx->A,&pv));
      70       25573 :   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       26133 : static PetscErrorCode BVDot_Mat(BV X,BV Y,Mat M)
      92             : {
      93       26133 :   BV_MAT            *x = (BV_MAT*)X->data,*y = (BV_MAT*)Y->data;
      94       26133 :   PetscScalar       *m;
      95       26133 :   const PetscScalar *px,*py;
      96       26133 :   PetscInt          ldm;
      97             : 
      98       26133 :   PetscFunctionBegin;
      99       26133 :   PetscCall(MatDenseGetLDA(M,&ldm));
     100       26133 :   PetscCall(MatDenseGetArrayRead(x->A,&px));
     101       26133 :   PetscCall(MatDenseGetArrayRead(y->A,&py));
     102       26133 :   PetscCall(MatDenseGetArray(M,&m));
     103       26133 :   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       26133 :   PetscCall(MatDenseRestoreArray(M,&m));
     105       26133 :   PetscCall(MatDenseRestoreArrayRead(x->A,&px));
     106       26133 :   PetscCall(MatDenseRestoreArrayRead(y->A,&py));
     107       26133 :   PetscFunctionReturn(PETSC_SUCCESS);
     108             : }
     109             : 
     110      466082 : static PetscErrorCode BVDotVec_Mat(BV X,Vec y,PetscScalar *q)
     111             : {
     112      466082 :   BV_MAT            *x = (BV_MAT*)X->data;
     113      466082 :   PetscScalar       *qq=q;
     114      466082 :   const PetscScalar *px,*py;
     115      466082 :   Vec               z = y;
     116             : 
     117      466082 :   PetscFunctionBegin;
     118      466082 :   if (PetscUnlikely(X->matrix)) {
     119      137166 :     PetscCall(BV_IPMatMult(X,y));
     120      137166 :     z = X->Bx;
     121             :   }
     122      466082 :   PetscCall(MatDenseGetArrayRead(x->A,&px));
     123      466082 :   PetscCall(VecGetArrayRead(z,&py));
     124      466082 :   if (!q) PetscCall(VecGetArray(X->buffer,&qq));
     125      466082 :   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      466082 :   if (!q) PetscCall(VecRestoreArray(X->buffer,&qq));
     127      466082 :   PetscCall(VecRestoreArrayRead(z,&py));
     128      466082 :   PetscCall(MatDenseRestoreArrayRead(x->A,&px));
     129      466082 :   PetscFunctionReturn(PETSC_SUCCESS);
     130             : }
     131             : 
     132        2485 : static PetscErrorCode BVDotVec_Local_Mat(BV X,Vec y,PetscScalar *m)
     133             : {
     134        2485 :   BV_MAT            *x = (BV_MAT*)X->data;
     135        2485 :   const PetscScalar *px,*py;
     136        2485 :   Vec               z = y;
     137             : 
     138        2485 :   PetscFunctionBegin;
     139        2485 :   if (PetscUnlikely(X->matrix)) {
     140           0 :     PetscCall(BV_IPMatMult(X,y));
     141           0 :     z = X->Bx;
     142             :   }
     143        2485 :   PetscCall(MatDenseGetArrayRead(x->A,&px));
     144        2485 :   PetscCall(VecGetArrayRead(z,&py));
     145        2485 :   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        2485 :   PetscCall(VecRestoreArrayRead(z,&py));
     147        2485 :   PetscCall(MatDenseRestoreArrayRead(x->A,&px));
     148        2485 :   PetscFunctionReturn(PETSC_SUCCESS);
     149             : }
     150             : 
     151      232144 : static PetscErrorCode BVScale_Mat(BV bv,PetscInt j,PetscScalar alpha)
     152             : {
     153      232144 :   BV_MAT         *ctx = (BV_MAT*)bv->data;
     154      232144 :   PetscScalar    *array;
     155             : 
     156      232144 :   PetscFunctionBegin;
     157      232144 :   if (!bv->n) PetscFunctionReturn(PETSC_SUCCESS);
     158      232144 :   PetscCall(MatDenseGetArray(ctx->A,&array));
     159      232144 :   if (PetscUnlikely(j<0)) PetscCall(BVScale_BLAS_Private(bv,(bv->k-bv->l)*bv->ld,array+(bv->nc+bv->l)*bv->ld,alpha));
     160      231853 :   else PetscCall(BVScale_BLAS_Private(bv,bv->n,array+(bv->nc+j)*bv->ld,alpha));
     161      232144 :   PetscCall(MatDenseRestoreArray(ctx->A,&array));
     162      232144 :   PetscFunctionReturn(PETSC_SUCCESS);
     163             : }
     164             : 
     165        5187 : static PetscErrorCode BVNorm_Mat(BV bv,PetscInt j,NormType type,PetscReal *val)
     166             : {
     167        5187 :   BV_MAT            *ctx = (BV_MAT*)bv->data;
     168        5187 :   const PetscScalar *array;
     169             : 
     170        5187 :   PetscFunctionBegin;
     171        5187 :   PetscCall(MatDenseGetArrayRead(ctx->A,&array));
     172        5187 :   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        5115 :   else PetscCall(BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->ld,bv->ld,type,val,ctx->mpi));
     174        5187 :   PetscCall(MatDenseRestoreArrayRead(ctx->A,&array));
     175        5187 :   PetscFunctionReturn(PETSC_SUCCESS);
     176             : }
     177             : 
     178        2308 : static PetscErrorCode BVNorm_Local_Mat(BV bv,PetscInt j,NormType type,PetscReal *val)
     179             : {
     180        2308 :   BV_MAT            *ctx = (BV_MAT*)bv->data;
     181        2308 :   const PetscScalar *array;
     182             : 
     183        2308 :   PetscFunctionBegin;
     184        2308 :   PetscCall(MatDenseGetArrayRead(ctx->A,&array));
     185        2308 :   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        2308 :   else PetscCall(BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->ld,bv->ld,type,val,PETSC_FALSE));
     187        2308 :   PetscCall(MatDenseRestoreArrayRead(ctx->A,&array));
     188        2308 :   PetscFunctionReturn(PETSC_SUCCESS);
     189             : }
     190             : 
     191        3484 : static PetscErrorCode BVNormalize_Mat(BV bv,PetscScalar *eigi)
     192             : {
     193        3484 :   BV_MAT         *ctx = (BV_MAT*)bv->data;
     194        3484 :   PetscScalar    *array,*wi=NULL;
     195             : 
     196        3484 :   PetscFunctionBegin;
     197        3484 :   PetscCall(MatDenseGetArray(ctx->A,&array));
     198        3484 :   if (eigi) wi = eigi+bv->l;
     199        3484 :   PetscCall(BVNormalize_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->ld,bv->ld,wi,ctx->mpi));
     200        3484 :   PetscCall(MatDenseRestoreArray(ctx->A,&array));
     201        3484 :   PetscFunctionReturn(PETSC_SUCCESS);
     202             : }
     203             : 
     204       11272 : static PetscErrorCode BVMatMult_Mat(BV V,Mat A,BV W)
     205             : {
     206       11272 :   PetscInt       j;
     207       11272 :   Mat            Vmat,Wmat;
     208       11272 :   Vec            vv,ww;
     209             : 
     210       11272 :   PetscFunctionBegin;
     211       11272 :   if (V->vmm) {
     212       11265 :     PetscCall(BVGetMat(V,&Vmat));
     213       11265 :     PetscCall(BVGetMat(W,&Wmat));
     214       11265 :     PetscCall(MatProductCreateWithMat(A,Vmat,NULL,Wmat));
     215       11265 :     PetscCall(MatProductSetType(Wmat,MATPRODUCT_AB));
     216       11265 :     PetscCall(MatProductSetFromOptions(Wmat));
     217       11265 :     PetscCall(MatProductSymbolic(Wmat));
     218       11265 :     PetscCall(MatProductNumeric(Wmat));
     219       11265 :     PetscCall(MatProductClear(Wmat));
     220       11265 :     PetscCall(BVRestoreMat(V,&Vmat));
     221       11265 :     PetscCall(BVRestoreMat(W,&Wmat));
     222             :   } else {
     223          50 :     for (j=0;j<V->k-V->l;j++) {
     224          43 :       PetscCall(BVGetColumn(V,V->l+j,&vv));
     225          43 :       PetscCall(BVGetColumn(W,W->l+j,&ww));
     226          43 :       PetscCall(MatMult(A,vv,ww));
     227          43 :       PetscCall(BVRestoreColumn(V,V->l+j,&vv));
     228          43 :       PetscCall(BVRestoreColumn(W,W->l+j,&ww));
     229             :     }
     230             :   }
     231       11272 :   PetscFunctionReturn(PETSC_SUCCESS);
     232             : }
     233             : 
     234       13482 : static PetscErrorCode BVCopy_Mat(BV V,BV W)
     235             : {
     236       13482 :   BV_MAT            *v = (BV_MAT*)V->data,*w = (BV_MAT*)W->data;
     237       13482 :   const PetscScalar *pv;
     238       13482 :   PetscScalar       *pw;
     239       13482 :   PetscInt          j;
     240             : 
     241       13482 :   PetscFunctionBegin;
     242       13482 :   PetscCall(MatDenseGetArrayRead(v->A,&pv));
     243       13482 :   PetscCall(MatDenseGetArray(w->A,&pw));
     244       80423 :   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       13482 :   PetscCall(MatDenseRestoreArrayRead(v->A,&pv));
     246       13482 :   PetscCall(MatDenseRestoreArray(w->A,&pw));
     247       13482 :   PetscFunctionReturn(PETSC_SUCCESS);
     248             : }
     249             : 
     250        8906 : static PetscErrorCode BVCopyColumn_Mat(BV V,PetscInt j,PetscInt i)
     251             : {
     252        8906 :   BV_MAT         *v = (BV_MAT*)V->data;
     253        8906 :   PetscScalar    *pv;
     254             : 
     255        8906 :   PetscFunctionBegin;
     256        8906 :   PetscCall(MatDenseGetArray(v->A,&pv));
     257        8906 :   PetscCall(PetscArraycpy(pv+(V->nc+i)*V->ld,pv+(V->nc+j)*V->ld,V->n));
     258        8906 :   PetscCall(MatDenseRestoreArray(v->A,&pv));
     259        8906 :   PetscFunctionReturn(PETSC_SUCCESS);
     260             : }
     261             : 
     262          95 : static PetscErrorCode BVResize_Mat(BV bv,PetscInt m,PetscBool copy)
     263             : {
     264          95 :   BV_MAT            *ctx = (BV_MAT*)bv->data;
     265          95 :   Mat               A,Msrc,Mdst;
     266          95 :   char              str[50];
     267             : 
     268          95 :   PetscFunctionBegin;
     269          95 :   PetscCall(MatCreateDenseFromVecType(PetscObjectComm((PetscObject)bv),bv->vtype,bv->n,PETSC_DECIDE,bv->N,m,bv->ld,NULL,&A));
     270          95 :   if (((PetscObject)bv)->name) {
     271          15 :     PetscCall(PetscSNPrintf(str,sizeof(str),"%s_0",((PetscObject)bv)->name));
     272          15 :     PetscCall(PetscObjectSetName((PetscObject)A,str));
     273             :   }
     274          95 :   if (copy) {
     275           5 :     PetscCall(MatDenseGetSubMatrix(ctx->A,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,PetscMin(m,bv->m),&Msrc));
     276           5 :     PetscCall(MatDenseGetSubMatrix(A,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,PetscMin(m,bv->m),&Mdst));
     277           5 :     PetscCall(MatCopy(Msrc,Mdst,SAME_NONZERO_PATTERN));
     278           5 :     PetscCall(MatDenseRestoreSubMatrix(ctx->A,&Msrc));
     279           5 :     PetscCall(MatDenseRestoreSubMatrix(A,&Mdst));
     280             :   }
     281          95 :   PetscCall(MatDestroy(&ctx->A));
     282          95 :   ctx->A = A;
     283          95 :   PetscFunctionReturn(PETSC_SUCCESS);
     284             : }
     285             : 
     286     1546048 : static PetscErrorCode BVGetColumn_Mat(BV bv,PetscInt j,Vec *v)
     287             : {
     288     1546048 :   BV_MAT         *ctx = (BV_MAT*)bv->data;
     289     1546048 :   PetscScalar    *pA;
     290     1546048 :   PetscInt       l;
     291             : 
     292     1546048 :   PetscFunctionBegin;
     293     1546048 :   l = BVAvailableVec;
     294     1546048 :   PetscCall(MatDenseGetArray(ctx->A,&pA));
     295     1546048 :   PetscCall(VecPlaceArray(bv->cv[l],pA+(bv->nc+j)*bv->ld));
     296     1546048 :   PetscFunctionReturn(PETSC_SUCCESS);
     297             : }
     298             : 
     299     1546048 : static PetscErrorCode BVRestoreColumn_Mat(BV bv,PetscInt j,Vec *v)
     300             : {
     301     1546048 :   BV_MAT         *ctx = (BV_MAT*)bv->data;
     302     1546048 :   PetscScalar    *pA;
     303     1546048 :   PetscInt       l;
     304             : 
     305     1546048 :   PetscFunctionBegin;
     306     1546048 :   l = (j==bv->ci[0])? 0: 1;
     307     1546048 :   PetscCall(VecResetArray(bv->cv[l]));
     308     1546048 :   PetscCall(MatDenseRestoreArray(ctx->A,&pA));
     309     1546048 :   PetscFunctionReturn(PETSC_SUCCESS);
     310             : }
     311             : 
     312       25293 : static PetscErrorCode BVGetArray_Mat(BV bv,PetscScalar **a)
     313             : {
     314       25293 :   BV_MAT         *ctx = (BV_MAT*)bv->data;
     315             : 
     316       25293 :   PetscFunctionBegin;
     317       25293 :   PetscCall(MatDenseGetArray(ctx->A,a));
     318       25293 :   PetscFunctionReturn(PETSC_SUCCESS);
     319             : }
     320             : 
     321       25293 : static PetscErrorCode BVRestoreArray_Mat(BV bv,PetscScalar **a)
     322             : {
     323       25293 :   BV_MAT         *ctx = (BV_MAT*)bv->data;
     324             : 
     325       25293 :   PetscFunctionBegin;
     326       25293 :   if (a) PetscCall(MatDenseRestoreArray(ctx->A,a));
     327       25293 :   PetscFunctionReturn(PETSC_SUCCESS);
     328             : }
     329             : 
     330         166 : static PetscErrorCode BVGetArrayRead_Mat(BV bv,const PetscScalar **a)
     331             : {
     332         166 :   BV_MAT         *ctx = (BV_MAT*)bv->data;
     333             : 
     334         166 :   PetscFunctionBegin;
     335         166 :   PetscCall(MatDenseGetArrayRead(ctx->A,a));
     336         166 :   PetscFunctionReturn(PETSC_SUCCESS);
     337             : }
     338             : 
     339         166 : static PetscErrorCode BVRestoreArrayRead_Mat(BV bv,const PetscScalar **a)
     340             : {
     341         166 :   BV_MAT         *ctx = (BV_MAT*)bv->data;
     342             : 
     343         166 :   PetscFunctionBegin;
     344         166 :   if (a) PetscCall(MatDenseRestoreArrayRead(ctx->A,a));
     345         166 :   PetscFunctionReturn(PETSC_SUCCESS);
     346             : }
     347             : 
     348          28 : static PetscErrorCode BVView_Mat(BV bv,PetscViewer viewer)
     349             : {
     350          28 :   BV_MAT            *ctx = (BV_MAT*)bv->data;
     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          10 :     PetscCall(MatView(ctx->A,viewer));
     361          10 :     if (format == PETSC_VIEWER_ASCII_MATLAB) {
     362           0 :       PetscCall(PetscObjectGetName((PetscObject)bv,&bvname));
     363           0 :       PetscCall(PetscObjectGetName((PetscObject)ctx->A,&name));
     364           0 :       PetscCall(PetscViewerASCIIPrintf(viewer,"%s=%s;clear %s\n",bvname,name,name));
     365           0 :       if (bv->nc) PetscCall(PetscViewerASCIIPrintf(viewer,"%s=%s(:,%" PetscInt_FMT ":end);\n",bvname,bvname,bv->nc+1));
     366             :     }
     367           0 :   } else PetscCall(MatView(ctx->A,viewer));
     368          10 :   PetscFunctionReturn(PETSC_SUCCESS);
     369             : }
     370             : 
     371       13046 : static PetscErrorCode BVDestroy_Mat(BV bv)
     372             : {
     373       13046 :   BV_MAT         *ctx = (BV_MAT*)bv->data;
     374             : 
     375       13046 :   PetscFunctionBegin;
     376       13046 :   PetscCall(MatDestroy(&ctx->A));
     377       13046 :   PetscCall(VecDestroy(&bv->cv[0]));
     378       13046 :   PetscCall(VecDestroy(&bv->cv[1]));
     379       13046 :   PetscCall(PetscFree(bv->data));
     380       13046 :   PetscFunctionReturn(PETSC_SUCCESS);
     381             : }
     382             : 
     383       13046 : SLEPC_EXTERN PetscErrorCode BVCreate_Mat(BV bv)
     384             : {
     385       13046 :   BV_MAT         *ctx;
     386       13046 :   PetscInt       nloc,lsplit;
     387       13046 :   PetscBool      seq;
     388       13046 :   char           str[50];
     389       13046 :   PetscScalar    *array,*ptr=NULL;
     390       13046 :   BV             parent;
     391       13046 :   Mat            Apar;
     392             : 
     393       13046 :   PetscFunctionBegin;
     394       13046 :   PetscCall(PetscNew(&ctx));
     395       13046 :   bv->data = (void*)ctx;
     396             : 
     397       13046 :   PetscCall(PetscStrcmpAny(bv->vtype,&bv->cuda,VECSEQCUDA,VECMPICUDA,""));
     398       13046 :   PetscCall(PetscStrcmpAny(bv->vtype,&ctx->mpi,VECMPI,VECMPICUDA,""));
     399             : 
     400       13046 :   PetscCall(PetscStrcmp(bv->vtype,VECSEQ,&seq));
     401       13046 :   PetscCheck(seq || ctx->mpi || bv->cuda,PetscObjectComm((PetscObject)bv),PETSC_ERR_SUP,"BVMAT does not support the requested vector type: %s",bv->vtype);
     402             : 
     403       13046 :   PetscCall(PetscLayoutGetLocalSize(bv->map,&nloc));
     404       13046 :   PetscCall(BV_SetDefaultLD(bv,nloc));
     405             : 
     406       13046 :   if (PetscUnlikely(bv->issplit)) {
     407             :     /* split BV: share the memory of the parent BV */
     408          80 :     parent = bv->splitparent;
     409          80 :     lsplit = parent->lsplit;
     410          80 :     Apar = ((BV_MAT*)parent->data)->A;
     411          80 :     if (bv->cuda) {
     412             : #if defined(PETSC_HAVE_CUDA)
     413             :       PetscCall(MatDenseCUDAGetArray(Apar,&array));
     414             :       ptr = (bv->issplit==1)? array: array+lsplit*bv->ld;
     415             :       PetscCall(MatDenseCUDARestoreArray(Apar,&array));
     416             : #endif
     417             :     } else {
     418          80 :       PetscCall(MatDenseGetArray(Apar,&array));
     419          80 :       ptr = (bv->issplit==1)? array: array+lsplit*bv->ld;
     420          80 :       PetscCall(MatDenseRestoreArray(Apar,&array));
     421             :     }
     422             :   }
     423             : 
     424       13046 :   PetscCall(MatCreateDenseFromVecType(PetscObjectComm((PetscObject)bv),bv->vtype,nloc,PETSC_DECIDE,bv->N,bv->m,bv->ld,ptr,&ctx->A));
     425       13046 :   if (((PetscObject)bv)->name) {
     426         122 :     PetscCall(PetscSNPrintf(str,sizeof(str),"%s_0",((PetscObject)bv)->name));
     427         122 :     PetscCall(PetscObjectSetName((PetscObject)ctx->A,str));
     428             :   }
     429             : 
     430       13046 :   if (PetscUnlikely(bv->Acreate)) {
     431         136 :     PetscCall(MatConvert(bv->Acreate,bv->cuda?MATDENSECUDA:MATDENSE,MAT_REUSE_MATRIX,&ctx->A));
     432          68 :     PetscCall(MatDestroy(&bv->Acreate));
     433             :   }
     434             : 
     435       13046 :   PetscCall(BVCreateVecEmpty(bv,&bv->cv[0]));
     436       13046 :   PetscCall(BVCreateVecEmpty(bv,&bv->cv[1]));
     437             : 
     438       13046 :   if (bv->cuda) {
     439             : #if defined(PETSC_HAVE_CUDA)
     440             :     bv->ops->mult             = BVMult_Mat_CUDA;
     441             :     bv->ops->multvec          = BVMultVec_Mat_CUDA;
     442             :     bv->ops->multinplace      = BVMultInPlace_Mat_CUDA;
     443             :     bv->ops->multinplacetrans = BVMultInPlaceHermitianTranspose_Mat_CUDA;
     444             :     bv->ops->dot              = BVDot_Mat_CUDA;
     445             :     bv->ops->dotvec           = BVDotVec_Mat_CUDA;
     446             :     bv->ops->dotvec_local     = BVDotVec_Local_Mat_CUDA;
     447             :     bv->ops->scale            = BVScale_Mat_CUDA;
     448             :     bv->ops->norm             = BVNorm_Mat_CUDA;
     449             :     bv->ops->norm_local       = BVNorm_Local_Mat_CUDA;
     450             :     bv->ops->normalize        = BVNormalize_Mat_CUDA;
     451             :     bv->ops->matmult          = BVMatMult_Mat_CUDA;
     452             :     bv->ops->copy             = BVCopy_Mat_CUDA;
     453             :     bv->ops->copycolumn       = BVCopyColumn_Mat_CUDA;
     454             :     bv->ops->getcolumn        = BVGetColumn_Mat_CUDA;
     455             :     bv->ops->restorecolumn    = BVRestoreColumn_Mat_CUDA;
     456             :     bv->ops->restoresplit     = BVRestoreSplit_Mat_CUDA;
     457             :     bv->ops->getmat           = BVGetMat_Mat_CUDA;
     458             :     bv->ops->restoremat       = BVRestoreMat_Mat_CUDA;
     459             : #endif
     460             :   } else {
     461       13046 :     bv->ops->mult             = BVMult_Mat;
     462       13046 :     bv->ops->multvec          = BVMultVec_Mat;
     463       13046 :     bv->ops->multinplace      = BVMultInPlace_Mat;
     464       13046 :     bv->ops->multinplacetrans = BVMultInPlaceHermitianTranspose_Mat;
     465       13046 :     bv->ops->dot              = BVDot_Mat;
     466       13046 :     bv->ops->dotvec           = BVDotVec_Mat;
     467       13046 :     bv->ops->dotvec_local     = BVDotVec_Local_Mat;
     468       13046 :     bv->ops->scale            = BVScale_Mat;
     469       13046 :     bv->ops->norm             = BVNorm_Mat;
     470       13046 :     bv->ops->norm_local       = BVNorm_Local_Mat;
     471       13046 :     bv->ops->normalize        = BVNormalize_Mat;
     472       13046 :     bv->ops->matmult          = BVMatMult_Mat;
     473       13046 :     bv->ops->copy             = BVCopy_Mat;
     474       13046 :     bv->ops->copycolumn       = BVCopyColumn_Mat;
     475       13046 :     bv->ops->getcolumn        = BVGetColumn_Mat;
     476       13046 :     bv->ops->restorecolumn    = BVRestoreColumn_Mat;
     477       13046 :     bv->ops->getmat           = BVGetMat_Default;
     478       13046 :     bv->ops->restoremat       = BVRestoreMat_Default;
     479             :   }
     480       13046 :   bv->ops->resize           = BVResize_Mat;
     481       13046 :   bv->ops->getarray         = BVGetArray_Mat;
     482       13046 :   bv->ops->restorearray     = BVRestoreArray_Mat;
     483       13046 :   bv->ops->getarrayread     = BVGetArrayRead_Mat;
     484       13046 :   bv->ops->restorearrayread = BVRestoreArrayRead_Mat;
     485       13046 :   bv->ops->destroy          = BVDestroy_Mat;
     486       13046 :   if (!ctx->mpi) bv->ops->view = BVView_Mat;
     487       13046 :   PetscFunctionReturn(PETSC_SUCCESS);
     488             : }

Generated by: LCOV version 1.14