Actual source code: slepcimpl.h
slepc-3.18.1 2022-11-02
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: ((!SlepcInitializeCalled && \
44: PetscError(comm,__LINE__,PETSC_FUNCTION_NAME,__FILE__,1,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: #if defined(PETSC_USE_COMPLEX)
62: re = PetscRealPart(eigr);
63: im = PetscImaginaryPart(eigr);
64: #else
65: re = eigr;
66: im = eigi;
67: #endif
68: /* print zero instead of tiny value */
69: if (PetscAbs(im) && PetscAbs(re)/PetscAbs(im)<PETSC_SMALL) re = 0.0;
70: if (PetscAbs(re) && PetscAbs(im)/PetscAbs(re)<PETSC_SMALL) im = 0.0;
71: /* print as real if imaginary part is zero */
72: if (im!=0.0) PetscViewerASCIIPrintf(viewer,"%.5f%+.5fi",(double)re,(double)im);
73: else PetscViewerASCIIPrintf(viewer,"%.5f",(double)re);
74: return 0;
75: }
77: /*
78: SlepcViewEigenvector - Outputs an eigenvector xr,xi to a viewer.
79: In complex scalars only xr is written.
80: The name of xr,xi is set before writing, based on the label, the index, and the name of obj.
81: */
82: static inline PetscErrorCode SlepcViewEigenvector(PetscViewer viewer,Vec xr,Vec xi,const char *label,PetscInt index,PetscObject obj)
83: {
84: size_t count;
85: char vname[30];
86: const char *pname;
88: PetscObjectGetName(obj,&pname);
89: PetscSNPrintfCount(vname,sizeof(vname),"%s%s",&count,label,PetscDefined(USE_COMPLEX)?"":"r");
90: count--;
91: PetscSNPrintf(vname+count,sizeof(vname)-count,"%" PetscInt_FMT "_%s",index,pname);
92: PetscObjectSetName((PetscObject)xr,vname);
93: VecView(xr,viewer);
94: #if !defined(PETSC_USE_COMPLEX)
95: vname[count-1] = 'i';
96: PetscObjectSetName((PetscObject)xi,vname);
97: VecView(xi,viewer);
98: #endif
99: return 0;
100: }
102: /* Macros for strings with different value in real and complex */
103: #if defined(PETSC_USE_COMPLEX)
104: #define SLEPC_STRING_HERMITIAN "hermitian"
105: #else
106: #define SLEPC_STRING_HERMITIAN "symmetric"
107: #endif
109: /* Private functions that are shared by several classes */
110: SLEPC_EXTERN PetscErrorCode SlepcBasisReference_Private(PetscInt,Vec*,PetscInt*,Vec**);
111: SLEPC_EXTERN PetscErrorCode SlepcBasisDestroy_Private(PetscInt*,Vec**);
112: SLEPC_EXTERN PetscErrorCode SlepcMonitorMakeKey_Internal(const char[],PetscViewerType,PetscViewerFormat,char[]);
113: SLEPC_EXTERN PetscErrorCode PetscViewerAndFormatCreate_Internal(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**);
115: SLEPC_INTERN PetscErrorCode SlepcCitationsInitialize(void);
116: SLEPC_INTERN PetscErrorCode SlepcInitialize_DynamicLibraries(void);
117: SLEPC_INTERN PetscErrorCode SlepcInitialize_Packages(void);
119: /* Definitions needed to work with CUDA kernels */
120: #if defined(PETSC_HAVE_CUDA)
121: #include <petscdevice_cuda.h>
123: #define X_AXIS 0
124: #define Y_AXIS 1
126: #define SLEPC_TILE_SIZE_X 32
127: #define SLEPC_BLOCK_SIZE_X 128
128: #define SLEPC_TILE_SIZE_Y 32
129: #define SLEPC_BLOCK_SIZE_Y 128
131: static inline PetscErrorCode SlepcKernelSetGrid1D(PetscInt rows,dim3 *dimGrid,dim3 *dimBlock,PetscInt *dimGrid_xcount)
132: {
133: int card;
134: struct cudaDeviceProp devprop;
136: cudaGetDevice(&card);
137: cudaGetDeviceProperties(&devprop,card);
138: *dimGrid_xcount = 1;
140: /* X axis */
141: dimGrid->x = 1;
142: dimBlock->x = SLEPC_BLOCK_SIZE_X;
143: if (rows>SLEPC_BLOCK_SIZE_X) dimGrid->x = (rows+SLEPC_BLOCK_SIZE_X-1)/SLEPC_BLOCK_SIZE_X;
144: else dimBlock->x = rows;
145: if (dimGrid->x>(unsigned)devprop.maxGridSize[X_AXIS]) {
146: *dimGrid_xcount = (dimGrid->x+(devprop.maxGridSize[X_AXIS]-1))/devprop.maxGridSize[X_AXIS];
147: dimGrid->x = devprop.maxGridSize[X_AXIS];
148: }
149: return 0;
150: }
152: static inline PetscErrorCode SlepcKernelSetGrid2DTiles(PetscInt rows,PetscInt cols,dim3 *dimGrid,dim3 *dimBlock,PetscInt *dimGrid_xcount,PetscInt *dimGrid_ycount)
153: {
154: int card;
155: struct cudaDeviceProp devprop;
157: cudaGetDevice(&card);
158: cudaGetDeviceProperties(&devprop,card);
159: *dimGrid_xcount = *dimGrid_ycount = 1;
161: /* X axis */
162: dimGrid->x = 1;
163: dimBlock->x = SLEPC_BLOCK_SIZE_X;
164: 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);
165: else dimBlock->x = (rows+SLEPC_TILE_SIZE_X-1)/SLEPC_TILE_SIZE_X;
166: if (dimGrid->x>(unsigned)devprop.maxGridSize[X_AXIS]) {
167: *dimGrid_xcount = (dimGrid->x+(devprop.maxGridSize[X_AXIS]-1))/devprop.maxGridSize[X_AXIS];
168: dimGrid->x = devprop.maxGridSize[X_AXIS];
169: }
171: /* Y axis */
172: dimGrid->y = 1;
173: dimBlock->y = SLEPC_BLOCK_SIZE_Y;
174: 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);
175: else dimBlock->y = (cols+SLEPC_TILE_SIZE_Y-1)/SLEPC_TILE_SIZE_Y;
176: if (dimGrid->y>(unsigned)devprop.maxGridSize[Y_AXIS]) {
177: *dimGrid_ycount = (dimGrid->y+(devprop.maxGridSize[Y_AXIS]-1))/devprop.maxGridSize[Y_AXIS];
178: dimGrid->y = devprop.maxGridSize[Y_AXIS];
179: }
180: return 0;
181: }
182: #undef X_AXIS
183: #undef Y_AXIS
184: #endif
186: #endif