| Line | Branch | Exec | Source |
|---|---|---|---|
| 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 | #pragma once | ||
| 12 | |||
| 13 | #include <slepcsys.h> | ||
| 14 | #include <petsc/private/petscimpl.h> | ||
| 15 | |||
| 16 | /* SUBMANSEC = Sys */ | ||
| 17 | |||
| 18 | SLEPC_INTERN PetscBool SlepcBeganPetsc; | ||
| 19 | |||
| 20 | /* SlepcSwap - swap two variables a,b of the same type using a temporary variable t */ | ||
| 21 | #define SlepcSwap(a,b,t) do {t=a;a=b;b=t;} while (0) | ||
| 22 | |||
| 23 | /*MC | ||
| 24 | SlepcHeaderCreate - Creates a SLEPc object | ||
| 25 | |||
| 26 | Input Parameters: | ||
| 27 | + classid - the classid associated with this object | ||
| 28 | . class_name - string name of class; should be static | ||
| 29 | . descr - string containing short description; should be static | ||
| 30 | . mansec - string indicating section in manual pages; should be static | ||
| 31 | . comm - the MPI Communicator | ||
| 32 | . destroy - the destroy routine for this object | ||
| 33 | - view - the view routine for this object | ||
| 34 | |||
| 35 | Output Parameter: | ||
| 36 | . h - the newly created object | ||
| 37 | |||
| 38 | Note: | ||
| 39 | This is equivalent to `PetscHeaderCreate()` but makes sure that `SlepcInitialize()` | ||
| 40 | has been called. | ||
| 41 | |||
| 42 | Level: developer | ||
| 43 | M*/ | ||
| 44 | #define SlepcHeaderCreate(h,classid,class_name,descr,mansec,comm,destroy,view) \ | ||
| 45 | ((PetscErrorCode)((!SlepcInitializeCalled && \ | ||
| 46 | PetscError(comm,__LINE__,PETSC_FUNCTION_NAME,__FILE__,PETSC_ERR_ORDER,PETSC_ERROR_INITIAL, \ | ||
| 47 | "Must call SlepcInitialize instead of PetscInitialize to use SLEPc classes")) || \ | ||
| 48 | PetscHeaderCreate(h,classid,class_name,descr,mansec,comm,destroy,view))) | ||
| 49 | |||
| 50 | /* context for structured eigenproblem matrices created via MatCreateXXX */ | ||
| 51 | struct _n_SlepcMatStruct { | ||
| 52 | PetscInt cookie; /* identify which structured matrix */ | ||
| 53 | PetscScalar s; /* in BSE sign of the bottom part of the vector */ | ||
| 54 | }; | ||
| 55 | typedef struct _n_SlepcMatStruct* SlepcMatStruct; | ||
| 56 | |||
| 57 | #define SLEPC_MAT_STRUCT_BSE 88101 | ||
| 58 | #define SLEPC_MAT_STRUCT_HAMILT 88102 | ||
| 59 | #define SLEPC_MAT_STRUCT_LREP 88103 | ||
| 60 | |||
| 61 | /* | ||
| 62 | SlepcCheckMatStruct - Check that a given Mat is a structured matrix of the wanted type. | ||
| 63 | |||
| 64 | Returns true/false in flg if it is given, otherwise yields an error if the check fails. | ||
| 65 | If cookie==0 it will check for any type. | ||
| 66 | */ | ||
| 67 | 8345 | static inline PetscErrorCode SlepcCheckMatStruct(Mat A,PetscInt cookie,PetscBool *flg) | |
| 68 | { | ||
| 69 | 8345 | PetscContainer container; | |
| 70 | 8345 | SlepcMatStruct mctx; | |
| 71 | |||
| 72 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
8345 | PetscFunctionBegin; |
| 73 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 20 times.
|
8345 | if (flg) *flg = PETSC_FALSE; |
| 74 |
4/6✓ Branch 0 taken 6 times.
✓ Branch 1 taken 24 times.
✓ Branch 2 taken 6 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 6 times.
|
8345 | PetscCall(PetscObjectQuery((PetscObject)A,"SlepcMatStruct",(PetscObject*)&container)); |
| 75 |
10/16✓ Branch 0 taken 10 times.
✓ Branch 1 taken 20 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
|
8345 | if (flg && !container) PetscFunctionReturn(PETSC_SUCCESS); |
| 76 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 20 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
8266 | PetscCheck(container,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"The Mat is not a structured matrix"); |
| 77 |
1/2✓ Branch 0 taken 20 times.
✗ Branch 1 not taken.
|
8266 | if (cookie) { |
| 78 |
4/6✓ Branch 0 taken 4 times.
✓ Branch 1 taken 16 times.
✓ Branch 2 taken 4 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 4 times.
|
8266 | PetscCall(PetscContainerGetPointer(container,(void**)&mctx)); |
| 79 |
1/18✗ Branch 0 not taken.
✓ Branch 1 taken 20 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✗ Branch 17 not taken.
|
8266 | if (flg && (!mctx || mctx->cookie!=cookie)) PetscFunctionReturn(PETSC_SUCCESS); |
| 80 |
2/6✓ Branch 0 taken 20 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 20 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
8266 | PetscCheck(mctx && mctx->cookie==cookie,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"The type of structured matrix is different from the expected one"); |
| 81 | } | ||
| 82 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 20 times.
|
8266 | if (flg) *flg = PETSC_TRUE; |
| 83 |
6/12✓ Branch 0 taken 4 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 4 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 4 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 4 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 4 times.
|
1394 | PetscFunctionReturn(PETSC_SUCCESS); |
| 84 | } | ||
| 85 | |||
| 86 | /* | ||
| 87 | SlepcCheckMatLREPReduced - Check whether a given LREP Mat is of reduced form or not. | ||
| 88 | |||
| 89 | True if the matrix has reduced form [0 K; M 0], otherwise it has the form [A B; -B -A]. | ||
| 90 | */ | ||
| 91 | 7572 | static inline PetscErrorCode SlepcCheckMatLREPReduced(Mat A,PetscBool *flg) | |
| 92 | { | ||
| 93 | 7572 | Mat A00; | |
| 94 | |||
| 95 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
7572 | PetscFunctionBegin; |
| 96 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
7572 | PetscAssertPointer(flg,2); |
| 97 | 7572 | *flg = PETSC_FALSE; | |
| 98 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
7572 | PetscCall(SlepcCheckMatStruct(A,SLEPC_MAT_STRUCT_LREP,NULL)); /* error out if not LREP */ |
| 99 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
7572 | PetscCall(MatNestGetSubMat(A,0,0,&A00)); |
| 100 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
7572 | if (A00==NULL) *flg = PETSC_TRUE; |
| 101 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
1284 | PetscFunctionReturn(PETSC_SUCCESS); |
| 102 | } | ||
| 103 | |||
| 104 | /* | ||
| 105 | SlepcPrintEigenvalueASCII - Print an eigenvalue on an ASCII viewer. | ||
| 106 | */ | ||
| 107 | 55256 | static inline PetscErrorCode SlepcPrintEigenvalueASCII(PetscViewer viewer,PetscScalar eigr,PetscScalar eigi) | |
| 108 | { | ||
| 109 | 55256 | PetscReal re,im; | |
| 110 | |||
| 111 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
55256 | PetscFunctionBegin; |
| 112 | #if defined(PETSC_USE_COMPLEX) | ||
| 113 | 27591 | re = PetscRealPart(eigr); | |
| 114 | 27591 | im = PetscImaginaryPart(eigr); | |
| 115 | 27591 | (void)eigi; | |
| 116 | #else | ||
| 117 | 27665 | re = eigr; | |
| 118 | 27665 | im = eigi; | |
| 119 | #endif | ||
| 120 | /* print zero instead of tiny value */ | ||
| 121 |
8/8✓ Branch 0 taken 25 times.
✓ Branch 1 taken 30 times.
✓ Branch 2 taken 20 times.
✓ Branch 3 taken 25 times.
✓ Branch 4 taken 25 times.
✓ Branch 5 taken 25 times.
✓ Branch 6 taken 23 times.
✓ Branch 7 taken 22 times.
|
55256 | if (PetscAbs(im) && PetscAbs(re)/PetscAbs(im)<PETSC_SMALL) re = 0.0; |
| 122 |
8/8✓ Branch 0 taken 30 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 25 times.
✓ Branch 3 taken 30 times.
✓ Branch 4 taken 20 times.
✓ Branch 5 taken 30 times.
✓ Branch 6 taken 26 times.
✓ Branch 7 taken 29 times.
|
54870 | if (PetscAbs(re) && PetscAbs(im)/PetscAbs(re)<PETSC_SMALL) im = 0.0; |
| 123 | /* print as real if imaginary part is zero */ | ||
| 124 |
6/8✓ Branch 0 taken 25 times.
✓ Branch 1 taken 9 times.
✓ Branch 2 taken 5 times.
✓ Branch 3 taken 20 times.
✓ Branch 4 taken 5 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 5 times.
|
19762 | if (im!=(PetscReal)0.0) PetscCall(PetscViewerASCIIPrintf(viewer,"%.5f%+.5fi",(double)re,(double)im)); |
| 125 |
4/6✓ Branch 0 taken 6 times.
✓ Branch 1 taken 24 times.
✓ Branch 2 taken 6 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 6 times.
|
50779 | else PetscCall(PetscViewerASCIIPrintf(viewer,"%.5f",(double)re)); |
| 126 |
6/12✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 6 times.
✓ Branch 4 taken 6 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 6 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 6 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 6 times.
|
10855 | PetscFunctionReturn(PETSC_SUCCESS); |
| 127 | } | ||
| 128 | |||
| 129 | /* | ||
| 130 | SlepcViewEigenvector - Outputs an eigenvector xr,xi to a viewer. | ||
| 131 | In complex scalars only xr is written. | ||
| 132 | The name of xr,xi is set before writing, based on the label, the index, and the name of obj. | ||
| 133 | */ | ||
| 134 | 135 | static inline PetscErrorCode SlepcViewEigenvector(PetscViewer viewer,Vec xr,Vec xi,const char *label,PetscInt index,PetscObject obj) | |
| 135 | { | ||
| 136 | 135 | size_t count; | |
| 137 | 135 | char vname[30]; | |
| 138 | 135 | const char *pname; | |
| 139 | |||
| 140 |
1/2✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
|
135 | PetscFunctionBegin; |
| 141 |
4/6✓ Branch 0 taken 6 times.
✓ Branch 1 taken 24 times.
✓ Branch 2 taken 6 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 6 times.
|
135 | PetscCall(PetscObjectGetName(obj,&pname)); |
| 142 |
4/6✓ Branch 0 taken 6 times.
✓ Branch 1 taken 24 times.
✓ Branch 2 taken 6 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 6 times.
|
135 | PetscCall(PetscSNPrintfCount(vname,sizeof(vname),"%s%s",&count,label,PetscDefined(USE_COMPLEX)?"":"r")); |
| 143 | 135 | count--; | |
| 144 |
4/6✓ Branch 0 taken 6 times.
✓ Branch 1 taken 24 times.
✓ Branch 2 taken 6 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 6 times.
|
135 | PetscCall(PetscSNPrintf(vname+count,sizeof(vname)-count,"%" PetscInt_FMT "_%s",index,pname)); |
| 145 |
4/6✓ Branch 0 taken 6 times.
✓ Branch 1 taken 24 times.
✓ Branch 2 taken 6 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 6 times.
|
135 | PetscCall(PetscObjectSetName((PetscObject)xr,vname)); |
| 146 |
4/6✓ Branch 0 taken 6 times.
✓ Branch 1 taken 24 times.
✓ Branch 2 taken 6 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 6 times.
|
135 | PetscCall(VecView(xr,viewer)); |
| 147 | 83 | if (!PetscDefined(USE_COMPLEX)) { | |
| 148 | 70 | vname[count-1] = 'i'; | |
| 149 |
4/6✓ Branch 0 taken 3 times.
✓ Branch 1 taken 12 times.
✓ Branch 2 taken 3 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 3 times.
|
70 | PetscCall(PetscObjectSetName((PetscObject)xi,vname)); |
| 150 |
4/6✓ Branch 0 taken 3 times.
✓ Branch 1 taken 12 times.
✓ Branch 2 taken 3 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 3 times.
|
83 | PetscCall(VecView(xi,viewer)); |
| 151 | } | ||
| 152 |
6/12✓ Branch 0 taken 6 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 6 times.
✓ Branch 4 taken 6 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 6 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 6 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 6 times.
|
27 | PetscFunctionReturn(PETSC_SUCCESS); |
| 153 | } | ||
| 154 | |||
| 155 | /* Macros for strings with different value in real and complex */ | ||
| 156 | #if defined(PETSC_USE_COMPLEX) | ||
| 157 | #define SLEPC_STRING_HERMITIAN "hermitian" | ||
| 158 | #else | ||
| 159 | #define SLEPC_STRING_HERMITIAN "symmetric" | ||
| 160 | #endif | ||
| 161 | |||
| 162 | /* Private functions that are shared by several classes */ | ||
| 163 | SLEPC_SINGLE_LIBRARY_INTERN PetscErrorCode SlepcBasisReference_Private(PetscInt,Vec*,PetscInt*,Vec**); | ||
| 164 | SLEPC_SINGLE_LIBRARY_INTERN PetscErrorCode SlepcBasisDestroy_Private(PetscInt*,Vec**); | ||
| 165 | SLEPC_SINGLE_LIBRARY_INTERN PetscErrorCode SlepcMonitorMakeKey_Internal(const char[],PetscViewerType,PetscViewerFormat,char[]); | ||
| 166 | SLEPC_SINGLE_LIBRARY_INTERN PetscErrorCode PetscViewerAndFormatCreate_Internal(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**); | ||
| 167 | |||
| 168 | SLEPC_INTERN PetscErrorCode SlepcCitationsInitialize(void); | ||
| 169 | SLEPC_INTERN PetscErrorCode SlepcInitialize_DynamicLibraries(void); | ||
| 170 | SLEPC_INTERN PetscErrorCode SlepcInitialize_Packages(void); | ||
| 171 | |||
| 172 | /* Macro to check a sequential Mat (including GPU) */ | ||
| 173 | #if !defined(PETSC_USE_DEBUG) | ||
| 174 | #define SlepcMatCheckSeq(h) do {(void)(h);} while (0) | ||
| 175 | #else | ||
| 176 | #if defined(PETSC_HAVE_CUDA) | ||
| 177 | #define SlepcMatCheckSeq(h) do { PetscCheckTypeNames((h),MATSEQDENSE,MATSEQDENSECUDA); } while (0) | ||
| 178 | #elif defined(PETSC_HAVE_HIP) | ||
| 179 | #define SlepcMatCheckSeq(h) do { PetscCheckTypeNames((h),MATSEQDENSE,MATSEQDENSEHIP); } while (0) | ||
| 180 | #else | ||
| 181 | #define SlepcMatCheckSeq(h) do { PetscCheckTypeName((h),MATSEQDENSE); } while (0) | ||
| 182 | #endif | ||
| 183 | #endif | ||
| 184 | |||
| 185 | /* Definitions needed to work with GPU kernels */ | ||
| 186 | #if defined(PETSC_HAVE_CUPM) | ||
| 187 | #include <petscdevice_cupm.h> | ||
| 188 | |||
| 189 | #define X_AXIS 0 | ||
| 190 | #define Y_AXIS 1 | ||
| 191 | |||
| 192 | #define SLEPC_TILE_SIZE_X 32 | ||
| 193 | #define SLEPC_BLOCK_SIZE_X 128 | ||
| 194 | #define SLEPC_TILE_SIZE_Y 32 | ||
| 195 | #define SLEPC_BLOCK_SIZE_Y 128 | ||
| 196 | |||
| 197 | 3392 | static inline PetscErrorCode SlepcKernelSetGrid1D(PetscInt rows,dim3 *dimGrid,dim3 *dimBlock,PetscInt *dimGrid_xcount) | |
| 198 | { | ||
| 199 | 3392 | int card; | |
| 200 | #if defined(PETSC_HAVE_CUDA) | ||
| 201 | 3392 | struct cudaDeviceProp devprop; | |
| 202 | #elif defined(PETSC_HAVE_HIP) | ||
| 203 | hipDeviceProp_t devprop; | ||
| 204 | #endif | ||
| 205 | |||
| 206 | 3392 | PetscFunctionBegin; | |
| 207 | #if defined(PETSC_HAVE_CUDA) | ||
| 208 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
3392 | PetscCallCUDA(cudaGetDevice(&card)); |
| 209 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
3392 | PetscCallCUDA(cudaGetDeviceProperties(&devprop,card)); |
| 210 | #elif defined(PETSC_HAVE_HIP) | ||
| 211 | PetscCallHIP(hipGetDevice(&card)); | ||
| 212 | PetscCallHIP(hipGetDeviceProperties(&devprop,card)); | ||
| 213 | #endif | ||
| 214 | 3392 | *dimGrid_xcount = 1; | |
| 215 | |||
| 216 | /* X axis */ | ||
| 217 | 3392 | dimGrid->x = 1; | |
| 218 | 3392 | dimBlock->x = SLEPC_BLOCK_SIZE_X; | |
| 219 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 4 times.
|
3392 | if (rows>SLEPC_BLOCK_SIZE_X) dimGrid->x = (rows+SLEPC_BLOCK_SIZE_X-1)/SLEPC_BLOCK_SIZE_X; |
| 220 | 3392 | else dimBlock->x = rows; | |
| 221 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 4 times.
|
3392 | if (dimGrid->x>(unsigned)devprop.maxGridSize[X_AXIS]) { |
| 222 | ✗ | *dimGrid_xcount = (dimGrid->x+(devprop.maxGridSize[X_AXIS]-1))/devprop.maxGridSize[X_AXIS]; | |
| 223 | ✗ | dimGrid->x = devprop.maxGridSize[X_AXIS]; | |
| 224 | } | ||
| 225 | PetscFunctionReturn(PETSC_SUCCESS); | ||
| 226 | } | ||
| 227 | |||
| 228 | 72 | static inline PetscErrorCode SlepcKernelSetGrid2DTiles(PetscInt rows,PetscInt cols,dim3 *dimGrid,dim3 *dimBlock,PetscInt *dimGrid_xcount,PetscInt *dimGrid_ycount) | |
| 229 | { | ||
| 230 | 72 | int card; | |
| 231 | #if defined(PETSC_HAVE_CUDA) | ||
| 232 | 72 | struct cudaDeviceProp devprop; | |
| 233 | #elif defined(PETSC_HAVE_HIP) | ||
| 234 | hipDeviceProp_t devprop; | ||
| 235 | #endif | ||
| 236 | |||
| 237 | 72 | PetscFunctionBegin; | |
| 238 | #if defined(PETSC_HAVE_CUDA) | ||
| 239 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
72 | PetscCallCUDA(cudaGetDevice(&card)); |
| 240 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
72 | PetscCallCUDA(cudaGetDeviceProperties(&devprop,card)); |
| 241 | #elif defined(PETSC_HAVE_HIP) | ||
| 242 | PetscCallHIP(hipGetDevice(&card)); | ||
| 243 | PetscCallHIP(hipGetDeviceProperties(&devprop,card)); | ||
| 244 | #endif | ||
| 245 | 72 | *dimGrid_xcount = *dimGrid_ycount = 1; | |
| 246 | |||
| 247 | /* X axis */ | ||
| 248 | 72 | dimGrid->x = 1; | |
| 249 | 72 | dimBlock->x = SLEPC_BLOCK_SIZE_X; | |
| 250 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
72 | if (rows>SLEPC_BLOCK_SIZE_X*SLEPC_TILE_SIZE_X) dimGrid->x = (rows+SLEPC_BLOCK_SIZE_X*SLEPC_TILE_SIZE_X-1)/(SLEPC_BLOCK_SIZE_X*SLEPC_TILE_SIZE_X); |
| 251 | 72 | else dimBlock->x = (rows+SLEPC_TILE_SIZE_X-1)/SLEPC_TILE_SIZE_X; | |
| 252 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
72 | if (dimGrid->x>(unsigned)devprop.maxGridSize[X_AXIS]) { |
| 253 | ✗ | *dimGrid_xcount = (dimGrid->x+(devprop.maxGridSize[X_AXIS]-1))/devprop.maxGridSize[X_AXIS]; | |
| 254 | ✗ | dimGrid->x = devprop.maxGridSize[X_AXIS]; | |
| 255 | } | ||
| 256 | |||
| 257 | /* Y axis */ | ||
| 258 | 72 | dimGrid->y = 1; | |
| 259 | 72 | dimBlock->y = SLEPC_BLOCK_SIZE_Y; | |
| 260 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
72 | if (cols>SLEPC_BLOCK_SIZE_Y*SLEPC_TILE_SIZE_Y) dimGrid->y = (cols+SLEPC_BLOCK_SIZE_Y*SLEPC_TILE_SIZE_Y-1)/(SLEPC_BLOCK_SIZE_Y*SLEPC_TILE_SIZE_Y); |
| 261 | 72 | else dimBlock->y = (cols+SLEPC_TILE_SIZE_Y-1)/SLEPC_TILE_SIZE_Y; | |
| 262 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
72 | if (dimGrid->y>(unsigned)devprop.maxGridSize[Y_AXIS]) { |
| 263 | ✗ | *dimGrid_ycount = (dimGrid->y+(devprop.maxGridSize[Y_AXIS]-1))/devprop.maxGridSize[Y_AXIS]; | |
| 264 | ✗ | dimGrid->y = devprop.maxGridSize[Y_AXIS]; | |
| 265 | } | ||
| 266 | PetscFunctionReturn(PETSC_SUCCESS); | ||
| 267 | } | ||
| 268 | #undef X_AXIS | ||
| 269 | #undef Y_AXIS | ||
| 270 | #endif | ||
| 271 |