Actual source code: slepcimpl.h
slepc-3.19.0 2023-03-31
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
11: #if !defined(SLEPCIMPL_H)
12: #define SLEPCIMPL_H
14: #include <slepcsys.h>
15: #include <petsc/private/petscimpl.h>
17: /* SUBMANSEC = sys */
19: SLEPC_INTERN PetscBool SlepcBeganPetsc;
21: /*@C
22: SlepcHeaderCreate - Creates a SLEPc object
24: Input Parameters:
25: + classid - the classid associated with this object
26: . class_name - string name of class; should be static
27: . descr - string containing short description; should be static
28: . mansec - string indicating section in manual pages; should be static
29: . comm - the MPI Communicator
30: . destroy - the destroy routine for this object
31: - view - the view routine for this object
33: Output Parameter:
34: . h - the newly created object
36: Note:
37: This is equivalent to PetscHeaderCreate but makes sure that SlepcInitialize
38: has been called.
40: Level: developer
41: @*/
42: #define SlepcHeaderCreate(h,classid,class_name,descr,mansec,comm,destroy,view) \
43: ((PetscErrorCode)((!SlepcInitializeCalled && \
44: PetscError(comm,__LINE__,PETSC_FUNCTION_NAME,__FILE__,PETSC_ERR_ORDER,PETSC_ERROR_INITIAL, \
45: "Must call SlepcInitialize instead of PetscInitialize to use SLEPc classes")) || \
46: PetscHeaderCreate(h,classid,class_name,descr,mansec,comm,destroy,view)))
48: /* context for monitors of type XXXMonitorConverged */
49: struct _n_SlepcConvMon {
50: void *ctx;
51: PetscInt oldnconv; /* previous value of nconv */
52: };
54: /*
55: SlepcPrintEigenvalueASCII - Print an eigenvalue on an ASCII viewer.
56: */
57: static inline PetscErrorCode SlepcPrintEigenvalueASCII(PetscViewer viewer,PetscScalar eigr,PetscScalar eigi)
58: {
59: PetscReal re,im;
61: PetscFunctionBegin;
62: #if defined(PETSC_USE_COMPLEX)
63: re = PetscRealPart(eigr);
64: im = PetscImaginaryPart(eigr);
65: #else
66: re = eigr;
67: im = eigi;
68: #endif
69: /* print zero instead of tiny value */
70: if (PetscAbs(im) && PetscAbs(re)/PetscAbs(im)<PETSC_SMALL) re = 0.0;
71: if (PetscAbs(re) && PetscAbs(im)/PetscAbs(re)<PETSC_SMALL) im = 0.0;
72: /* print as real if imaginary part is zero */
73: if (im!=0.0) PetscCall(PetscViewerASCIIPrintf(viewer,"%.5f%+.5fi",(double)re,(double)im));
74: else PetscCall(PetscViewerASCIIPrintf(viewer,"%.5f",(double)re));
75: PetscFunctionReturn(PETSC_SUCCESS);
76: }
78: /*
79: SlepcViewEigenvector - Outputs an eigenvector xr,xi to a viewer.
80: In complex scalars only xr is written.
81: The name of xr,xi is set before writing, based on the label, the index, and the name of obj.
82: */
83: static inline PetscErrorCode SlepcViewEigenvector(PetscViewer viewer,Vec xr,Vec xi,const char *label,PetscInt index,PetscObject obj)
84: {
85: size_t count;
86: char vname[30];
87: const char *pname;
89: PetscFunctionBegin;
90: PetscCall(PetscObjectGetName(obj,&pname));
91: PetscCall(PetscSNPrintfCount(vname,sizeof(vname),"%s%s",&count,label,PetscDefined(USE_COMPLEX)?"":"r"));
92: count--;
93: PetscCall(PetscSNPrintf(vname+count,sizeof(vname)-count,"%" PetscInt_FMT "_%s",index,pname));
94: PetscCall(PetscObjectSetName((PetscObject)xr,vname));
95: PetscCall(VecView(xr,viewer));
96: #if !defined(PETSC_USE_COMPLEX)
97: vname[count-1] = 'i';
98: PetscCall(PetscObjectSetName((PetscObject)xi,vname));
99: PetscCall(VecView(xi,viewer));
100: #endif
101: PetscFunctionReturn(PETSC_SUCCESS);
102: }
104: /* Macros for strings with different value in real and complex */
105: #if defined(PETSC_USE_COMPLEX)
106: #define SLEPC_STRING_HERMITIAN "hermitian"
107: #else
108: #define SLEPC_STRING_HERMITIAN "symmetric"
109: #endif
111: /* Private functions that are shared by several classes */
112: SLEPC_EXTERN PetscErrorCode SlepcBasisReference_Private(PetscInt,Vec*,PetscInt*,Vec**);
113: SLEPC_EXTERN PetscErrorCode SlepcBasisDestroy_Private(PetscInt*,Vec**);
114: SLEPC_EXTERN PetscErrorCode SlepcMonitorMakeKey_Internal(const char[],PetscViewerType,PetscViewerFormat,char[]);
115: SLEPC_EXTERN PetscErrorCode PetscViewerAndFormatCreate_Internal(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**);
117: SLEPC_INTERN PetscErrorCode SlepcCitationsInitialize(void);
118: SLEPC_INTERN PetscErrorCode SlepcInitialize_DynamicLibraries(void);
119: SLEPC_INTERN PetscErrorCode SlepcInitialize_Packages(void);
121: /* Definitions needed to work with CUDA kernels */
122: #if defined(PETSC_HAVE_CUDA)
123: #include <petscdevice_cuda.h>
125: #define X_AXIS 0
126: #define Y_AXIS 1
128: #define SLEPC_TILE_SIZE_X 32
129: #define SLEPC_BLOCK_SIZE_X 128
130: #define SLEPC_TILE_SIZE_Y 32
131: #define SLEPC_BLOCK_SIZE_Y 128
133: static inline PetscErrorCode SlepcKernelSetGrid1D(PetscInt rows,dim3 *dimGrid,dim3 *dimBlock,PetscInt *dimGrid_xcount)
134: {
135: int card;
136: struct cudaDeviceProp devprop;
138: PetscFunctionBegin;
139: PetscCallCUDA(cudaGetDevice(&card));
140: PetscCallCUDA(cudaGetDeviceProperties(&devprop,card));
141: *dimGrid_xcount = 1;
143: /* X axis */
144: dimGrid->x = 1;
145: dimBlock->x = SLEPC_BLOCK_SIZE_X;
146: if (rows>SLEPC_BLOCK_SIZE_X) dimGrid->x = (rows+SLEPC_BLOCK_SIZE_X-1)/SLEPC_BLOCK_SIZE_X;
147: else dimBlock->x = rows;
148: if (dimGrid->x>(unsigned)devprop.maxGridSize[X_AXIS]) {
149: *dimGrid_xcount = (dimGrid->x+(devprop.maxGridSize[X_AXIS]-1))/devprop.maxGridSize[X_AXIS];
150: dimGrid->x = devprop.maxGridSize[X_AXIS];
151: }
152: PetscFunctionReturn(PETSC_SUCCESS);
153: }
155: static inline PetscErrorCode SlepcKernelSetGrid2DTiles(PetscInt rows,PetscInt cols,dim3 *dimGrid,dim3 *dimBlock,PetscInt *dimGrid_xcount,PetscInt *dimGrid_ycount)
156: {
157: int card;
158: struct cudaDeviceProp devprop;
160: PetscFunctionBegin;
161: PetscCallCUDA(cudaGetDevice(&card));
162: PetscCallCUDA(cudaGetDeviceProperties(&devprop,card));
163: *dimGrid_xcount = *dimGrid_ycount = 1;
165: /* X axis */
166: dimGrid->x = 1;
167: dimBlock->x = SLEPC_BLOCK_SIZE_X;
168: 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);
169: else dimBlock->x = (rows+SLEPC_TILE_SIZE_X-1)/SLEPC_TILE_SIZE_X;
170: if (dimGrid->x>(unsigned)devprop.maxGridSize[X_AXIS]) {
171: *dimGrid_xcount = (dimGrid->x+(devprop.maxGridSize[X_AXIS]-1))/devprop.maxGridSize[X_AXIS];
172: dimGrid->x = devprop.maxGridSize[X_AXIS];
173: }
175: /* Y axis */
176: dimGrid->y = 1;
177: dimBlock->y = SLEPC_BLOCK_SIZE_Y;
178: 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);
179: else dimBlock->y = (cols+SLEPC_TILE_SIZE_Y-1)/SLEPC_TILE_SIZE_Y;
180: if (dimGrid->y>(unsigned)devprop.maxGridSize[Y_AXIS]) {
181: *dimGrid_ycount = (dimGrid->y+(devprop.maxGridSize[Y_AXIS]-1))/devprop.maxGridSize[Y_AXIS];
182: dimGrid->y = devprop.maxGridSize[Y_AXIS];
183: }
184: PetscFunctionReturn(PETSC_SUCCESS);
185: }
186: #undef X_AXIS
187: #undef Y_AXIS
188: #endif
190: #endif