| Line | Branch | Exec | Source |
|---|---|---|---|
| 1 | 3553 | /* | |
| 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 (CUDA version) | ||
| 12 | */ | ||
| 13 | |||
| 14 | #include <slepc/private/bvimpl.h> | ||
| 15 | #include <slepccupmblas.h> | ||
| 16 | #include "../src/sys/classes/bv/impls/mat/bvmat.h" | ||
| 17 | |||
| 18 | 560 | PetscErrorCode BVMult_Mat_CUDA(BV Y,PetscScalar alpha,PetscScalar beta,BV X,Mat Q) | |
| 19 | { | ||
| 20 | 560 | BV_MAT *y = (BV_MAT*)Y->data,*x = (BV_MAT*)X->data; | |
| 21 | 560 | const PetscScalar *d_px,*d_A,*d_B,*d_q; | |
| 22 | 560 | PetscScalar *d_py,*d_C; | |
| 23 | 560 | PetscInt ldq; | |
| 24 | |||
| 25 | 560 | PetscFunctionBegin; | |
| 26 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
560 | if (!Y->n) PetscFunctionReturn(PETSC_SUCCESS); |
| 27 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
560 | PetscCall(MatDenseCUDAGetArrayRead(x->A,&d_px)); |
| 28 |
3/4✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
|
588 | if (beta==(PetscScalar)0.0) PetscCall(MatDenseCUDAGetArrayWrite(y->A,&d_py)); |
| 29 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
504 | else PetscCall(MatDenseCUDAGetArray(y->A,&d_py)); |
| 30 | 560 | d_A = d_px+(X->nc+X->l)*X->ld; | |
| 31 | 560 | d_C = d_py+(Y->nc+Y->l)*Y->ld; | |
| 32 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
|
560 | if (Q) { |
| 33 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
318 | PetscCall(MatDenseGetLDA(Q,&ldq)); |
| 34 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
318 | PetscCall(BV_MatDenseCUDAGetArrayRead(Y,Q,&d_q)); |
| 35 | 318 | d_B = d_q+Y->l*ldq+X->l; | |
| 36 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
318 | PetscCall(BVMult_BLAS_CUDA(Y,Y->n,Y->k-Y->l,X->k-X->l,alpha,d_A,X->ld,d_B,ldq,beta,d_C,Y->ld)); |
| 37 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
318 | PetscCall(BV_MatDenseCUDARestoreArrayRead(Y,Q,&d_q)); |
| 38 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
242 | } else PetscCall(BVAXPY_BLAS_CUDA(Y,Y->n,Y->k-Y->l,alpha,d_A,X->ld,beta,d_C,Y->ld)); |
| 39 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
560 | PetscCall(MatDenseCUDARestoreArrayRead(x->A,&d_px)); |
| 40 |
3/4✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
|
588 | if (beta==(PetscScalar)0.0) PetscCall(MatDenseCUDARestoreArrayWrite(y->A,&d_py)); |
| 41 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
504 | else PetscCall(MatDenseCUDARestoreArray(y->A,&d_py)); |
| 42 | PetscFunctionReturn(PETSC_SUCCESS); | ||
| 43 | } | ||
| 44 | |||
| 45 | 66220 | PetscErrorCode BVMultVec_Mat_CUDA(BV X,PetscScalar alpha,PetscScalar beta,Vec y,PetscScalar *q) | |
| 46 | { | ||
| 47 | 66220 | BV_MAT *x = (BV_MAT*)X->data; | |
| 48 | 66220 | PetscScalar *d_py,*d_q; | |
| 49 | 66220 | const PetscScalar *d_px; | |
| 50 | |||
| 51 | 66220 | PetscFunctionBegin; | |
| 52 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
66220 | PetscCall(MatDenseCUDAGetArrayRead(x->A,&d_px)); |
| 53 |
3/4✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
|
70098 | if (beta==(PetscScalar)0.0) PetscCall(VecCUDAGetArrayWrite(y,&d_py)); |
| 54 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
58563 | else PetscCall(VecCUDAGetArray(y,&d_py)); |
| 55 |
3/4✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
|
66220 | if (!q) PetscCall(VecCUDAGetArray(X->buffer,&d_q)); |
| 56 | else { | ||
| 57 | 34506 | PetscInt k=X->k-X->l; | |
| 58 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
34506 | PetscCallCUDA(cudaMalloc((void**)&d_q,k*sizeof(PetscScalar))); |
| 59 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
34506 | PetscCallCUDA(cudaMemcpy(d_q,q,k*sizeof(PetscScalar),cudaMemcpyHostToDevice)); |
| 60 | 34506 | PetscCall(PetscLogCpuToGpu(k*sizeof(PetscScalar))); | |
| 61 | } | ||
| 62 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
66220 | PetscCall(BVMultVec_BLAS_CUDA(X,X->n,X->k-X->l,alpha,d_px+(X->nc+X->l)*X->ld,X->ld,d_q,beta,d_py)); |
| 63 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
66220 | PetscCall(MatDenseCUDARestoreArrayRead(x->A,&d_px)); |
| 64 |
3/4✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
|
70098 | if (beta==(PetscScalar)0.0) PetscCall(VecCUDARestoreArrayWrite(y,&d_py)); |
| 65 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
58563 | else PetscCall(VecCUDARestoreArray(y,&d_py)); |
| 66 |
3/4✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
|
66220 | if (!q) PetscCall(VecCUDARestoreArray(X->buffer,&d_q)); |
| 67 |
1/4✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
34506 | else PetscCallCUDA(cudaFree(d_q)); |
| 68 | PetscFunctionReturn(PETSC_SUCCESS); | ||
| 69 | } | ||
| 70 | |||
| 71 | 3530 | PetscErrorCode BVMultInPlace_Mat_CUDA(BV V,Mat Q,PetscInt s,PetscInt e) | |
| 72 | { | ||
| 73 | 3530 | BV_MAT *ctx = (BV_MAT*)V->data; | |
| 74 | 3530 | PetscScalar *d_pv; | |
| 75 | 3530 | const PetscScalar *d_q; | |
| 76 | 3530 | PetscInt ldq; | |
| 77 | |||
| 78 | 3530 | PetscFunctionBegin; | |
| 79 |
3/4✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
|
3530 | if (s>=e || !V->n) PetscFunctionReturn(PETSC_SUCCESS); |
| 80 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
3263 | PetscCall(MatDenseGetLDA(Q,&ldq)); |
| 81 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
3263 | PetscCall(MatDenseCUDAGetArray(ctx->A,&d_pv)); |
| 82 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
3263 | PetscCall(BV_MatDenseCUDAGetArrayRead(V,Q,&d_q)); |
| 83 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
3263 | PetscCall(BVMultInPlace_BLAS_CUDA(V,V->n,V->k-V->l,s-V->l,e-V->l,d_pv+(V->nc+V->l)*V->ld,V->ld,d_q+V->l*ldq+V->l,ldq,PETSC_FALSE)); |
| 84 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
3263 | PetscCall(BV_MatDenseCUDARestoreArrayRead(V,Q,&d_q)); |
| 85 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
3263 | PetscCall(MatDenseCUDARestoreArray(ctx->A,&d_pv)); |
| 86 | PetscFunctionReturn(PETSC_SUCCESS); | ||
| 87 | } | ||
| 88 | |||
| 89 | ✗ | PetscErrorCode BVMultInPlaceHermitianTranspose_Mat_CUDA(BV V,Mat Q,PetscInt s,PetscInt e) | |
| 90 | { | ||
| 91 | ✗ | BV_MAT *ctx = (BV_MAT*)V->data; | |
| 92 | ✗ | PetscScalar *d_pv; | |
| 93 | ✗ | const PetscScalar *d_q; | |
| 94 | ✗ | PetscInt ldq; | |
| 95 | |||
| 96 | ✗ | PetscFunctionBegin; | |
| 97 | ✗ | if (s>=e || !V->n) PetscFunctionReturn(PETSC_SUCCESS); | |
| 98 | ✗ | PetscCall(MatDenseGetLDA(Q,&ldq)); | |
| 99 | ✗ | PetscCall(MatDenseCUDAGetArray(ctx->A,&d_pv)); | |
| 100 | ✗ | PetscCall(BV_MatDenseCUDAGetArrayRead(V,Q,&d_q)); | |
| 101 | ✗ | PetscCall(BVMultInPlace_BLAS_CUDA(V,V->n,V->k-V->l,s-V->l,e-V->l,d_pv+(V->nc+V->l)*V->ld,V->ld,d_q+V->l*ldq+V->l,ldq,PETSC_TRUE)); | |
| 102 | ✗ | PetscCall(BV_MatDenseCUDARestoreArrayRead(V,Q,&d_q)); | |
| 103 | ✗ | PetscCall(MatDenseCUDARestoreArray(ctx->A,&d_pv)); | |
| 104 | PetscFunctionReturn(PETSC_SUCCESS); | ||
| 105 | } | ||
| 106 | |||
| 107 | 4939 | PetscErrorCode BVDot_Mat_CUDA(BV X,BV Y,Mat M) | |
| 108 | { | ||
| 109 | 4939 | BV_MAT *x = (BV_MAT*)X->data,*y = (BV_MAT*)Y->data; | |
| 110 | 4939 | const PetscScalar *d_px,*d_py; | |
| 111 | 4939 | PetscScalar *pm; | |
| 112 | 4939 | PetscInt ldm; | |
| 113 | |||
| 114 | 4939 | PetscFunctionBegin; | |
| 115 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4939 | PetscCall(MatDenseGetLDA(M,&ldm)); |
| 116 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4939 | PetscCall(MatDenseCUDAGetArrayRead(x->A,&d_px)); |
| 117 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4939 | PetscCall(MatDenseCUDAGetArrayRead(y->A,&d_py)); |
| 118 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4939 | PetscCall(MatDenseGetArrayWrite(M,&pm)); |
| 119 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4939 | PetscCall(BVDot_BLAS_CUDA(X,Y->k-Y->l,X->k-X->l,X->n,d_py+(Y->nc+Y->l)*Y->ld,Y->ld,d_px+(X->nc+X->l)*X->ld,X->ld,pm+X->l*ldm+Y->l,ldm,x->mpi)); |
| 120 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4939 | PetscCall(MatDenseRestoreArrayWrite(M,&pm)); |
| 121 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4939 | PetscCall(MatDenseCUDARestoreArrayRead(x->A,&d_px)); |
| 122 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
4939 | PetscCall(MatDenseCUDARestoreArrayRead(y->A,&d_py)); |
| 123 | PetscFunctionReturn(PETSC_SUCCESS); | ||
| 124 | } | ||
| 125 | |||
| 126 | 37952 | PetscErrorCode BVDotVec_Mat_CUDA(BV X,Vec y,PetscScalar *q) | |
| 127 | { | ||
| 128 | 37952 | BV_MAT *x = (BV_MAT*)X->data; | |
| 129 | 37952 | const PetscScalar *d_px,*d_py; | |
| 130 | 37952 | Vec z = y; | |
| 131 | |||
| 132 | 37952 | PetscFunctionBegin; | |
| 133 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
|
37952 | if (PetscUnlikely(X->matrix)) { |
| 134 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
12414 | PetscCall(BV_IPMatMult(X,y)); |
| 135 | 12414 | z = X->Bx; | |
| 136 | } | ||
| 137 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
37952 | PetscCall(MatDenseCUDAGetArrayRead(x->A,&d_px)); |
| 138 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
37952 | PetscCall(VecCUDAGetArrayRead(z,&d_py)); |
| 139 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
37952 | PetscCall(BVDotVec_BLAS_CUDA(X,X->n,X->k-X->l,d_px+(X->nc+X->l)*X->ld,X->ld,d_py,q,x->mpi)); |
| 140 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
37952 | PetscCall(VecCUDARestoreArrayRead(z,&d_py)); |
| 141 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
37952 | PetscCall(MatDenseCUDARestoreArrayRead(x->A,&d_px)); |
| 142 | PetscFunctionReturn(PETSC_SUCCESS); | ||
| 143 | } | ||
| 144 | |||
| 145 | 1462 | PetscErrorCode BVDotVec_Local_Mat_CUDA(BV X,Vec y,PetscScalar *m) | |
| 146 | { | ||
| 147 | 1462 | BV_MAT *x = (BV_MAT*)X->data; | |
| 148 | 1462 | const PetscScalar *d_px,*d_py; | |
| 149 | 1462 | Vec z = y; | |
| 150 | |||
| 151 | 1462 | PetscFunctionBegin; | |
| 152 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
1462 | if (PetscUnlikely(X->matrix)) { |
| 153 | ✗ | PetscCall(BV_IPMatMult(X,y)); | |
| 154 | ✗ | z = X->Bx; | |
| 155 | } | ||
| 156 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
1462 | PetscCall(MatDenseCUDAGetArrayRead(x->A,&d_px)); |
| 157 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
1462 | PetscCall(VecCUDAGetArrayRead(z,&d_py)); |
| 158 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
1462 | PetscCall(BVDotVec_BLAS_CUDA(X,X->n,X->k-X->l,d_px+(X->nc+X->l)*X->ld,X->ld,d_py,m,PETSC_FALSE)); |
| 159 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
1462 | PetscCall(VecCUDARestoreArrayRead(z,&d_py)); |
| 160 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
1462 | PetscCall(MatDenseCUDARestoreArrayRead(x->A,&d_px)); |
| 161 | PetscFunctionReturn(PETSC_SUCCESS); | ||
| 162 | } | ||
| 163 | |||
| 164 | 24716 | PetscErrorCode BVScale_Mat_CUDA(BV bv,PetscInt j,PetscScalar alpha) | |
| 165 | { | ||
| 166 | 24716 | BV_MAT *ctx = (BV_MAT*)bv->data; | |
| 167 | 24716 | PetscScalar *d_array,*d_A; | |
| 168 | 24716 | PetscInt n=0; | |
| 169 | |||
| 170 | 24716 | PetscFunctionBegin; | |
| 171 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
24716 | if (!bv->n) PetscFunctionReturn(PETSC_SUCCESS); |
| 172 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
24716 | PetscCall(MatDenseCUDAGetArray(ctx->A,&d_array)); |
| 173 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
|
24716 | if (PetscUnlikely(j<0)) { |
| 174 | 220 | d_A = d_array+(bv->nc+bv->l)*bv->ld; | |
| 175 | 220 | n = (bv->k-bv->l)*bv->ld; | |
| 176 | } else { | ||
| 177 | 24496 | d_A = d_array+(bv->nc+j)*bv->ld; | |
| 178 | 24496 | n = bv->n; | |
| 179 | } | ||
| 180 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
24716 | PetscCall(BVScale_BLAS_CUDA(bv,n,d_A,alpha)); |
| 181 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
24716 | PetscCall(MatDenseCUDARestoreArray(ctx->A,&d_array)); |
| 182 | PetscFunctionReturn(PETSC_SUCCESS); | ||
| 183 | } | ||
| 184 | |||
| 185 | 1226 | PetscErrorCode BVNorm_Mat_CUDA(BV bv,PetscInt j,NormType type,PetscReal *val) | |
| 186 | { | ||
| 187 | 1226 | BV_MAT *ctx = (BV_MAT*)bv->data; | |
| 188 | 1226 | const PetscScalar *array,*d_array,*d_A; | |
| 189 | 1226 | PetscInt n=0; | |
| 190 | |||
| 191 | 1226 | PetscFunctionBegin; | |
| 192 |
6/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
|
1226 | if (!ctx->mpi && ((j<0 && type==NORM_FROBENIUS && bv->ld==bv->n) || (j>=0 && type==NORM_2))) { |
| 193 | /* compute on GPU with cuBLAS - TODO: include the MPI case here */ | ||
| 194 | 1094 | *val = 0.0; | |
| 195 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
1094 | if (!bv->n) PetscFunctionReturn(PETSC_SUCCESS); |
| 196 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
1094 | PetscCall(MatDenseCUDAGetArrayRead(ctx->A,&d_array)); |
| 197 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
|
1094 | if (PetscUnlikely(j<0)) { |
| 198 | 20 | d_A = d_array+(bv->nc+bv->l)*bv->ld; | |
| 199 | 20 | n = (bv->k-bv->l)*bv->ld; | |
| 200 | } else { | ||
| 201 | 1074 | d_A = d_array+(bv->nc+j)*bv->ld; | |
| 202 | 1074 | n = bv->n; | |
| 203 | } | ||
| 204 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
1094 | PetscCall(BVNorm_BLAS_CUDA(bv,n,d_A,val)); |
| 205 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
1094 | PetscCall(MatDenseCUDARestoreArrayRead(ctx->A,&d_array)); |
| 206 | } else { | ||
| 207 | /* compute on CPU */ | ||
| 208 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
132 | PetscCall(MatDenseGetArrayRead(ctx->A,&array)); |
| 209 |
3/4✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
|
132 | 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)); |
| 210 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
40 | else PetscCall(BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->ld,bv->ld,type,val,ctx->mpi)); |
| 211 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
132 | PetscCall(MatDenseRestoreArrayRead(ctx->A,&array)); |
| 212 | } | ||
| 213 | PetscFunctionReturn(PETSC_SUCCESS); | ||
| 214 | } | ||
| 215 | |||
| 216 | 1380 | PetscErrorCode BVNorm_Local_Mat_CUDA(BV bv,PetscInt j,NormType type,PetscReal *val) | |
| 217 | { | ||
| 218 | 1380 | BV_MAT *ctx = (BV_MAT*)bv->data; | |
| 219 | 1380 | const PetscScalar *array,*d_array,*d_A; | |
| 220 | 1380 | PetscInt n=0; | |
| 221 | |||
| 222 | 1380 | PetscFunctionBegin; | |
| 223 |
2/6✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
|
1380 | if ((j<0 && type==NORM_FROBENIUS && bv->ld==bv->n) || (j>=0 && type==NORM_2)) { |
| 224 | /* compute on GPU with cuBLAS */ | ||
| 225 | 1380 | *val = 0.0; | |
| 226 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
1380 | if (!bv->n) PetscFunctionReturn(PETSC_SUCCESS); |
| 227 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
1380 | PetscCall(MatDenseCUDAGetArrayRead(ctx->A,&d_array)); |
| 228 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
1380 | if (PetscUnlikely(j<0)) { |
| 229 | ✗ | d_A = d_array+(bv->nc+bv->l)*bv->ld; | |
| 230 | ✗ | n = (bv->k-bv->l)*bv->ld; | |
| 231 | } else { | ||
| 232 | 1380 | d_A = d_array+(bv->nc+j)*bv->ld; | |
| 233 | 1380 | n = bv->n; | |
| 234 | } | ||
| 235 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
1380 | PetscCall(BVNorm_BLAS_CUDA(bv,n,d_A,val)); |
| 236 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
1380 | PetscCall(MatDenseCUDARestoreArrayRead(ctx->A,&d_array)); |
| 237 | } else { | ||
| 238 | /* compute on CPU */ | ||
| 239 | ✗ | PetscCall(MatDenseGetArrayRead(ctx->A,&array)); | |
| 240 | ✗ | 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)); | |
| 241 | ✗ | else PetscCall(BVNorm_LAPACK_Private(bv,bv->n,1,array+(bv->nc+j)*bv->ld,bv->ld,type,val,PETSC_FALSE)); | |
| 242 | ✗ | PetscCall(MatDenseRestoreArrayRead(ctx->A,&array)); | |
| 243 | } | ||
| 244 | PetscFunctionReturn(PETSC_SUCCESS); | ||
| 245 | } | ||
| 246 | |||
| 247 | 511 | PetscErrorCode BVNormalize_Mat_CUDA(BV bv,PetscScalar *eigi) | |
| 248 | { | ||
| 249 | 511 | BV_MAT *ctx = (BV_MAT*)bv->data; | |
| 250 | 511 | PetscScalar *array,*d_array,*wi=NULL; | |
| 251 | |||
| 252 | 511 | PetscFunctionBegin; | |
| 253 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
|
511 | if (eigi) wi = eigi+bv->l; |
| 254 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
|
511 | if (!ctx->mpi) { |
| 255 | /* compute on GPU with cuBLAS - TODO: include the MPI case here */ | ||
| 256 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
493 | if (!bv->n) PetscFunctionReturn(PETSC_SUCCESS); |
| 257 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
493 | PetscCall(MatDenseCUDAGetArray(ctx->A,&d_array)); |
| 258 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
493 | PetscCall(BVNormalize_BLAS_CUDA(bv,bv->n,bv->k-bv->l,d_array+(bv->nc+bv->l)*bv->ld,bv->ld,wi)); |
| 259 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
493 | PetscCall(MatDenseCUDARestoreArray(ctx->A,&d_array)); |
| 260 | } else { | ||
| 261 | /* compute on CPU */ | ||
| 262 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
18 | PetscCall(MatDenseGetArray(ctx->A,&array)); |
| 263 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
18 | PetscCall(BVNormalize_LAPACK_Private(bv,bv->n,bv->k-bv->l,array+(bv->nc+bv->l)*bv->ld,bv->ld,wi,ctx->mpi)); |
| 264 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
18 | PetscCall(MatDenseRestoreArray(ctx->A,&array)); |
| 265 | } | ||
| 266 | PetscFunctionReturn(PETSC_SUCCESS); | ||
| 267 | } | ||
| 268 | |||
| 269 | 1061 | PetscErrorCode BVMatMult_Mat_CUDA(BV V,Mat A,BV W) | |
| 270 | { | ||
| 271 | 1061 | BV_MAT *v = (BV_MAT*)V->data,*w = (BV_MAT*)W->data; | |
| 272 | 1061 | Mat Vmat,Wmat; | |
| 273 | 1061 | const PetscScalar *d_pv; | |
| 274 | 1061 | PetscScalar *d_pw; | |
| 275 | 1061 | PetscInt j; | |
| 276 | |||
| 277 | 1061 | PetscFunctionBegin; | |
| 278 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
|
1061 | if (V->vmm) { |
| 279 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
1055 | PetscCall(BVGetMat(V,&Vmat)); |
| 280 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
1055 | PetscCall(BVGetMat(W,&Wmat)); |
| 281 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
1055 | PetscCall(MatProductCreateWithMat(A,Vmat,NULL,Wmat)); |
| 282 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
1055 | PetscCall(MatProductSetType(Wmat,MATPRODUCT_AB)); |
| 283 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
1055 | PetscCall(MatProductSetFromOptions(Wmat)); |
| 284 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
1055 | PetscCall(MatProductSymbolic(Wmat)); |
| 285 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
1055 | PetscCall(MatProductNumeric(Wmat)); |
| 286 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
1055 | PetscCall(MatProductClear(Wmat)); |
| 287 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
1055 | PetscCall(BVRestoreMat(V,&Vmat)); |
| 288 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
1055 | PetscCall(BVRestoreMat(W,&Wmat)); |
| 289 | } else { | ||
| 290 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
6 | PetscCall(MatDenseCUDAGetArrayRead(v->A,&d_pv)); |
| 291 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
6 | PetscCall(MatDenseCUDAGetArrayWrite(w->A,&d_pw)); |
| 292 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
|
52 | for (j=0;j<V->k-V->l;j++) { |
| 293 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
46 | PetscCall(VecCUDAPlaceArray(V->cv[1],(PetscScalar *)d_pv+(V->nc+V->l+j)*V->ld)); |
| 294 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
46 | PetscCall(VecCUDAPlaceArray(W->cv[1],d_pw+(W->nc+W->l+j)*W->ld)); |
| 295 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
46 | PetscCall(MatMult(A,V->cv[1],W->cv[1])); |
| 296 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
46 | PetscCall(VecCUDAResetArray(V->cv[1])); |
| 297 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
46 | PetscCall(VecCUDAResetArray(W->cv[1])); |
| 298 | } | ||
| 299 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
6 | PetscCall(MatDenseCUDARestoreArrayRead(v->A,&d_pv)); |
| 300 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
6 | PetscCall(MatDenseCUDARestoreArrayWrite(w->A,&d_pw)); |
| 301 | } | ||
| 302 | PetscFunctionReturn(PETSC_SUCCESS); | ||
| 303 | } | ||
| 304 | |||
| 305 | 1285 | PetscErrorCode BVCopy_Mat_CUDA(BV V,BV W) | |
| 306 | { | ||
| 307 | 1285 | BV_MAT *v = (BV_MAT*)V->data,*w = (BV_MAT*)W->data; | |
| 308 | 1285 | const PetscScalar *d_pv; | |
| 309 | 1285 | PetscScalar *d_pw; | |
| 310 | |||
| 311 | 1285 | PetscFunctionBegin; | |
| 312 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
1285 | PetscCall(MatDenseCUDAGetArrayRead(v->A,&d_pv)); |
| 313 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
1285 | PetscCall(MatDenseCUDAGetArray(w->A,&d_pw)); |
| 314 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
1285 | PetscCallCUDA(cudaMemcpy2D(d_pw+(W->nc+W->l)*W->ld,W->ld*sizeof(PetscScalar),d_pv+(V->nc+V->l)*V->ld,V->ld*sizeof(PetscScalar),V->n*sizeof(PetscScalar),V->k-V->l,cudaMemcpyDeviceToDevice)); |
| 315 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
1285 | PetscCall(MatDenseCUDARestoreArrayRead(v->A,&d_pv)); |
| 316 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
1285 | PetscCall(MatDenseCUDARestoreArray(w->A,&d_pw)); |
| 317 | PetscFunctionReturn(PETSC_SUCCESS); | ||
| 318 | } | ||
| 319 | |||
| 320 | 665 | PetscErrorCode BVCopyColumn_Mat_CUDA(BV V,PetscInt j,PetscInt i) | |
| 321 | { | ||
| 322 | 665 | BV_MAT *v = (BV_MAT*)V->data; | |
| 323 | 665 | PetscScalar *d_pv; | |
| 324 | |||
| 325 | 665 | PetscFunctionBegin; | |
| 326 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
665 | PetscCall(MatDenseCUDAGetArray(v->A,&d_pv)); |
| 327 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
665 | PetscCallCUDA(cudaMemcpy(d_pv+(V->nc+i)*V->ld,d_pv+(V->nc+j)*V->ld,V->n*sizeof(PetscScalar),cudaMemcpyDeviceToDevice)); |
| 328 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
665 | PetscCall(MatDenseCUDARestoreArray(v->A,&d_pv)); |
| 329 | PetscFunctionReturn(PETSC_SUCCESS); | ||
| 330 | } | ||
| 331 | |||
| 332 | 145146 | PetscErrorCode BVGetColumn_Mat_CUDA(BV bv,PetscInt j,Vec*) | |
| 333 | { | ||
| 334 | 145146 | BV_MAT *ctx = (BV_MAT*)bv->data; | |
| 335 | 145146 | PetscScalar *d_pv; | |
| 336 | 145146 | PetscInt l; | |
| 337 | |||
| 338 | 145146 | PetscFunctionBegin; | |
| 339 |
3/4✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
|
145146 | l = BVAvailableVec; |
| 340 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
145146 | PetscCall(MatDenseCUDAGetArray(ctx->A,&d_pv)); |
| 341 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
145146 | PetscCall(VecCUDAPlaceArray(bv->cv[l],d_pv+(bv->nc+j)*bv->ld)); |
| 342 | PetscFunctionReturn(PETSC_SUCCESS); | ||
| 343 | } | ||
| 344 | |||
| 345 | 145146 | PetscErrorCode BVRestoreColumn_Mat_CUDA(BV bv,PetscInt j,Vec*) | |
| 346 | { | ||
| 347 | 145146 | BV_MAT *ctx = (BV_MAT*)bv->data; | |
| 348 | 145146 | PetscInt l; | |
| 349 | |||
| 350 | 145146 | PetscFunctionBegin; | |
| 351 | 145146 | l = (j==bv->ci[0])? 0: 1; | |
| 352 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
145146 | PetscCall(VecCUDAResetArray(bv->cv[l])); |
| 353 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
145146 | PetscCall(MatDenseCUDARestoreArray(ctx->A,NULL)); |
| 354 | PetscFunctionReturn(PETSC_SUCCESS); | ||
| 355 | } | ||
| 356 | |||
| 357 | 72 | PetscErrorCode BVRestoreSplit_Mat_CUDA(BV bv,BV *L,BV *R) | |
| 358 | { | ||
| 359 | 72 | Mat A; | |
| 360 | 72 | const PetscScalar *d_pv; | |
| 361 | 72 | PetscObjectState lstate,rstate; | |
| 362 | 72 | PetscBool change=PETSC_FALSE; | |
| 363 | |||
| 364 | 72 | PetscFunctionBegin; | |
| 365 | /* force sync flag to PETSC_CUDA_BOTH */ | ||
| 366 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
|
72 | if (L) { |
| 367 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
64 | PetscCall(PetscObjectStateGet((PetscObject)*L,&lstate)); |
| 368 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
|
64 | if (lstate != bv->lstate) { |
| 369 | 8 | A = ((BV_MAT*)bv->L->data)->A; | |
| 370 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
8 | PetscCall(MatDenseCUDAGetArrayRead(A,&d_pv)); |
| 371 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
8 | PetscCall(MatDenseCUDARestoreArrayRead(A,&d_pv)); |
| 372 | change = PETSC_TRUE; | ||
| 373 | } | ||
| 374 | } | ||
| 375 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
|
72 | if (R) { |
| 376 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
16 | PetscCall(PetscObjectStateGet((PetscObject)*R,&rstate)); |
| 377 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
|
16 | if (rstate != bv->rstate) { |
| 378 | 8 | A = ((BV_MAT*)bv->R->data)->A; | |
| 379 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
8 | PetscCall(MatDenseCUDAGetArrayRead(A,&d_pv)); |
| 380 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
8 | PetscCall(MatDenseCUDARestoreArrayRead(A,&d_pv)); |
| 381 | change = PETSC_TRUE; | ||
| 382 | } | ||
| 383 | } | ||
| 384 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
|
64 | if (change) { |
| 385 | 16 | A = ((BV_MAT*)bv->data)->A; | |
| 386 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
16 | PetscCall(MatDenseCUDAGetArray(A,(PetscScalar **)&d_pv)); |
| 387 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
16 | PetscCall(MatDenseCUDARestoreArray(A,(PetscScalar **)&d_pv)); |
| 388 | } | ||
| 389 | PetscFunctionReturn(PETSC_SUCCESS); | ||
| 390 | } | ||
| 391 | |||
| 392 | 36 | PetscErrorCode BVRestoreSplitRows_Mat_CUDA(BV bv,IS,IS,BV *U,BV *L) | |
| 393 | { | ||
| 394 | 36 | Mat A; | |
| 395 | 36 | const PetscScalar *d_pv; | |
| 396 | 36 | PetscObjectState lstate,rstate; | |
| 397 | 36 | PetscBool change=PETSC_FALSE; | |
| 398 | |||
| 399 | 36 | PetscFunctionBegin; | |
| 400 | /* force sync flag to PETSC_CUDA_BOTH */ | ||
| 401 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
36 | if (U) { |
| 402 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
36 | PetscCall(PetscObjectStateGet((PetscObject)*U,&rstate)); |
| 403 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
36 | if (rstate != bv->rstate) { |
| 404 | 36 | A = ((BV_MAT*)bv->R->data)->A; | |
| 405 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
36 | PetscCall(MatDenseCUDAGetArrayRead(A,&d_pv)); |
| 406 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
36 | PetscCall(MatDenseCUDARestoreArrayRead(A,&d_pv)); |
| 407 | change = PETSC_TRUE; | ||
| 408 | } | ||
| 409 | } | ||
| 410 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
36 | if (L) { |
| 411 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
36 | PetscCall(PetscObjectStateGet((PetscObject)*L,&lstate)); |
| 412 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
36 | if (lstate != bv->lstate) { |
| 413 | 36 | A = ((BV_MAT*)bv->L->data)->A; | |
| 414 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
36 | PetscCall(MatDenseCUDAGetArrayRead(A,&d_pv)); |
| 415 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
36 | PetscCall(MatDenseCUDARestoreArrayRead(A,&d_pv)); |
| 416 | change = PETSC_TRUE; | ||
| 417 | } | ||
| 418 | } | ||
| 419 | ✗ | if (change) { | |
| 420 | 36 | A = ((BV_MAT*)bv->data)->A; | |
| 421 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
36 | PetscCall(MatDenseCUDAGetArray(A,(PetscScalar **)&d_pv)); |
| 422 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
36 | PetscCall(MatDenseCUDARestoreArray(A,(PetscScalar **)&d_pv)); |
| 423 | } | ||
| 424 | PetscFunctionReturn(PETSC_SUCCESS); | ||
| 425 | } | ||
| 426 | |||
| 427 | 2393 | PetscErrorCode BVGetMat_Mat_CUDA(BV bv,Mat *A) | |
| 428 | { | ||
| 429 | 2393 | BV_MAT *ctx = (BV_MAT*)bv->data; | |
| 430 | 2393 | PetscScalar *vv,*aa; | |
| 431 | 2393 | PetscBool create=PETSC_FALSE; | |
| 432 | 2393 | PetscInt m,cols; | |
| 433 | |||
| 434 | 2393 | PetscFunctionBegin; | |
| 435 | 2393 | m = bv->k-bv->l; | |
| 436 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
|
2393 | if (!bv->Aget) create=PETSC_TRUE; |
| 437 | else { | ||
| 438 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
2061 | PetscCall(MatDenseCUDAGetArray(bv->Aget,&aa)); |
| 439 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
2061 | PetscCheck(!aa,PetscObjectComm((PetscObject)bv),PETSC_ERR_ARG_WRONGSTATE,"BVGetMat already called on this BV"); |
| 440 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
2061 | PetscCall(MatGetSize(bv->Aget,NULL,&cols)); |
| 441 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
|
2061 | if (cols!=m) { |
| 442 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
187 | PetscCall(MatDestroy(&bv->Aget)); |
| 443 | create=PETSC_TRUE; | ||
| 444 | } | ||
| 445 | } | ||
| 446 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
2393 | PetscCall(MatDenseCUDAGetArray(ctx->A,&vv)); |
| 447 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
|
2393 | if (create) { |
| 448 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
519 | PetscCall(MatCreateDenseFromVecType(PetscObjectComm((PetscObject)bv),bv->vtype,bv->n,PETSC_DECIDE,bv->N,m,bv->ld,vv,&bv->Aget)); /* pass a pointer to avoid allocation of storage */ |
| 449 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
519 | PetscCall(MatDenseCUDAReplaceArray(bv->Aget,NULL)); /* replace with a null pointer, the value after BVRestoreMat */ |
| 450 | } | ||
| 451 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
2393 | PetscCall(MatDenseCUDAPlaceArray(bv->Aget,vv+(bv->nc+bv->l)*bv->ld)); /* set the actual pointer */ |
| 452 | 2393 | *A = bv->Aget; | |
| 453 | 2393 | PetscFunctionReturn(PETSC_SUCCESS); | |
| 454 | } | ||
| 455 | |||
| 456 | 2393 | PetscErrorCode BVRestoreMat_Mat_CUDA(BV bv,Mat *A) | |
| 457 | { | ||
| 458 | 2393 | BV_MAT *ctx = (BV_MAT*)bv->data; | |
| 459 | 2393 | PetscScalar *vv,*aa; | |
| 460 | |||
| 461 | 2393 | PetscFunctionBegin; | |
| 462 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
2393 | PetscCall(MatDenseCUDAGetArray(bv->Aget,&aa)); |
| 463 | 2393 | vv = aa-(bv->nc+bv->l)*bv->ld; | |
| 464 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
2393 | PetscCall(MatDenseCUDAResetArray(bv->Aget)); |
| 465 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
2393 | PetscCall(MatDenseCUDARestoreArray(ctx->A,&vv)); |
| 466 | 2393 | *A = NULL; | |
| 467 | 2393 | PetscFunctionReturn(PETSC_SUCCESS); | |
| 468 | } | ||
| 469 |