Subversion Repositories slepc-dev

Rev

Go to most recent revision | Details | Compare with Previous | Last modification | View Log | RSS feed

Rev Author Line No. Line
6 dsic.upv.es!jroman 1
 
183 dsic.upv.es!antodo 2
#include "slepc.h" /*I "slepc.h" I*/
341 dsic.upv.es!antodo 3
#include "slepceps.h"
4
#include "slepcst.h"
6 dsic.upv.es!jroman 5
 
6
#undef __FUNCT__  
7
#define __FUNCT__ "SlepcPrintVersion"
8
/*
9
   SlepcPrintVersion - Prints SLEPc version info.
10
 
11
   Collective on MPI_Comm
12
*/
476 dsic.upv.es!antodo 13
PetscErrorCode SlepcPrintVersion(MPI_Comm comm)
6 dsic.upv.es!jroman 14
{
15
  int  info = 0;
16
 
17
  PetscFunctionBegin;
18
 
19
  info = (*PetscHelpPrintf)(comm,"--------------------------------------------\
20
------------------------------\n"); CHKERRQ(info);
21
  info = (*PetscHelpPrintf)(comm,"\t   %s\n",SLEPC_VERSION_NUMBER); CHKERRQ(info);
22
  info = (*PetscHelpPrintf)(comm,"%s",SLEPC_AUTHOR_INFO); CHKERRQ(info);
761 dsic.upv.es!jroman 23
  info = (*PetscHelpPrintf)(comm,"See docs/manual.html for help. \n"); CHKERRQ(info);
6 dsic.upv.es!jroman 24
#if !defined(PARCH_win32)
25
  info = (*PetscHelpPrintf)(comm,"SLEPc libraries linked from %s\n",SLEPC_LIB_DIR); CHKERRQ(info);
26
#endif
27
  info = (*PetscHelpPrintf)(comm,"--------------------------------------------\
28
------------------------------\n"); CHKERRQ(info);
29
 
30
  PetscFunctionReturn(info);
31
}
32
 
33
#undef __FUNCT__  
34
#define __FUNCT__ "SlepcPrintHelpIntro"
35
/*
36
   SlepcPrintHelpIntro - Prints introductory SLEPc help info.
37
 
38
   Collective on MPI_Comm
39
*/
476 dsic.upv.es!antodo 40
PetscErrorCode SlepcPrintHelpIntro(MPI_Comm comm)
6 dsic.upv.es!jroman 41
{
42
  int  info = 0;
43
 
44
  PetscFunctionBegin;
45
 
46
  info = (*PetscHelpPrintf)(comm,"--------------------------------------------\
47
------------------------------\n"); CHKERRQ(info);
48
  info = (*PetscHelpPrintf)(comm,"SLEPc help information includes that for the PETSc libraries, which provide\n"); CHKERRQ(info);
49
  info = (*PetscHelpPrintf)(comm,"low-level system infrastructure and linear algebra tools.\n"); CHKERRQ(info);
50
  info = (*PetscHelpPrintf)(comm,"--------------------------------------------\
51
------------------------------\n"); CHKERRQ(info);
52
 
53
  PetscFunctionReturn(info);
54
}
55
 
56
#undef __FUNCT__  
57
#define __FUNCT__ "SlepcRegisterEvents"
780 dsic.upv.es!jroman 58
PetscEvent EPS_SetUp, EPS_Solve, ST_SetUp, ST_Apply, ST_ApplyB, ST_ApplyNoB, ST_ApplyTranspose, EPS_Orthogonalize, ST_InnerProduct, EPS_ReverseProjection;
6 dsic.upv.es!jroman 59
 
60
/*
61
   SlepcRegisterEvents - Registers SLEPc events for use in performance logging.
62
*/
476 dsic.upv.es!antodo 63
PetscErrorCode SlepcRegisterEvents()
6 dsic.upv.es!jroman 64
{
65
  int  info = 0;
66
 
67
  PetscFunctionBegin;
68
 
69
  info = PetscLogEventRegister(&EPS_SetUp,"EPSSetUp",PETSC_NULL); CHKERRQ(info);
70
  info = PetscLogEventRegister(&EPS_Solve,"EPSSolve",PETSC_NULL); CHKERRQ(info);
71
  info = PetscLogEventRegister(&ST_SetUp,"STSetUp",PETSC_NULL); CHKERRQ(info);
72
  info = PetscLogEventRegister(&ST_Apply,"STApply",PETSC_NULL); CHKERRQ(info);
73
  info = PetscLogEventRegister(&ST_ApplyB,"STApplyB",PETSC_NULL); CHKERRQ(info);
74
  info = PetscLogEventRegister(&ST_ApplyNoB,"STApplyNoB",PETSC_NULL); CHKERRQ(info);
780 dsic.upv.es!jroman 75
  info = PetscLogEventRegister(&ST_ApplyTranspose,"STApplyTranspose",PETSC_NULL); CHKERRQ(info);
756 dsic.upv.es!antodo 76
  info = PetscLogEventRegister(&EPS_Orthogonalize,"EPSOrthogonalize",PETSC_NULL); CHKERRQ(info);
382 dsic.upv.es!antodo 77
  info = PetscLogEventRegister(&ST_InnerProduct,"STInnerProduct",PETSC_NULL); CHKERRQ(info);
756 dsic.upv.es!antodo 78
  info = PetscLogEventRegister(&EPS_ReverseProjection,"EPSReverseProjection",PETSC_NULL); CHKERRQ(info);
6 dsic.upv.es!jroman 79
 
80
  PetscFunctionReturn(info);
81
}
82
 
83
/* ------------------------Nasty global variables -------------------------------*/
84
/*
85
   Indicates whether SLEPc started PETSc, or whether it was
86
   already started before SLEPc was initialized.
87
*/
88
PetscTruth  SlepcBeganPetsc = PETSC_FALSE;
89
PetscTruth  SlepcInitializeCalled = PETSC_FALSE;
476 dsic.upv.es!antodo 90
PetscCookie EPS_COOKIE = 0;
91
PetscCookie ST_COOKIE = 0;
6 dsic.upv.es!jroman 92
 
93
#if defined(PETSC_USE_DYNAMIC_LIBRARIES)
94
extern PetscDLLibraryList DLLibrariesLoaded;
95
#endif
96
 
97
#undef __FUNCT__  
98
#define __FUNCT__ "SlepcInitialize"
99
/*@C
100
   SlepcInitialize - Initializes the SLEPc library. SlepcInitialize() calls
101
   PetscInitialize() if that has not been called yet, so this routine should
102
   always be called near the beginning of your program.
103
 
104
   Collective on MPI_COMM_WORLD or PETSC_COMM_WORLD if it has been set
105
 
106
   Input Parameters:
107
+  argc - count of number of command line arguments
108
.  args - the command line arguments
109
.  file - [optional] PETSc database file, defaults to ~username/.petscrc
110
          (use PETSC_NULL for default)
111
-  help - [optional] Help message to print, use PETSC_NULL for no message
112
 
113
   Fortran Note:
114
   Fortran syntax is very similar to that of PetscInitialize()
115
 
116
   Level: beginner
117
 
118
.seealso: SlepcInitializeFortran(), SlepcFinalize(), PetscInitialize()
119
@*/
476 dsic.upv.es!antodo 120
PetscErrorCode SlepcInitialize(int *argc,char ***args,char file[],const char help[])
6 dsic.upv.es!jroman 121
{
476 dsic.upv.es!antodo 122
  PetscErrorCode ierr,info=0;
373 dsic.upv.es!antodo 123
#if defined(PETSC_USE_DYNAMIC_LIBRARIES)
476 dsic.upv.es!antodo 124
  char           libs[PETSC_MAX_PATH_LEN],dlib[PETSC_MAX_PATH_LEN];
125
  PetscTruth     found;
373 dsic.upv.es!antodo 126
#endif
6 dsic.upv.es!jroman 127
 
128
  PetscFunctionBegin;
129
 
130
  if (SlepcInitializeCalled==PETSC_TRUE) {
131
    PetscFunctionReturn(0);
132
  }
133
 
134
#if !defined(PARCH_t3d)
135
  info = PetscSetHelpVersionFunctions(SlepcPrintHelpIntro,SlepcPrintVersion);CHKERRQ(info);
136
#endif
137
 
138
  if (!PetscInitializeCalled) {
139
    info = PetscInitialize(argc,args,file,help);CHKERRQ(info);
140
    SlepcBeganPetsc = PETSC_TRUE;
141
  }
142
 
143
  EPS_COOKIE = 0;
144
  ierr = PetscLogClassRegister(&EPS_COOKIE,"Eigenproblem Solver");CHKERRQ(ierr);
145
  ST_COOKIE = 0;
146
  ierr = PetscLogClassRegister(&ST_COOKIE,"Spectral Transform");CHKERRQ(ierr);
147
 
148
  /*
149
      Load the dynamic libraries
150
  */
151
 
152
#if defined(PETSC_USE_DYNAMIC_LIBRARIES)
153
  ierr = PetscStrcpy(libs,SLEPC_LIB_DIR);CHKERRQ(ierr);
154
  ierr = PetscStrcat(libs,"/libslepc");CHKERRQ(ierr);
155
  ierr = PetscDLLibraryRetrieve(PETSC_COMM_WORLD,libs,dlib,1024,&found);CHKERRQ(ierr);
156
  if (found) {
157
    ierr = PetscDLLibraryAppend(PETSC_COMM_WORLD,&DLLibrariesLoaded,libs);CHKERRQ(ierr);
158
  } else {
159
    SETERRQ1(1,"Unable to locate SLEPc dynamic library %s \n You cannot move the dynamic libraries!\n or remove USE_DYNAMIC_LIBRARIES from ${PETSC_DIR}/bmake/$PETSC_ARCH/petscconf.h\n and rebuild libraries before moving",libs);
160
  }
161
#else
162
 
163
  ierr = EPSRegisterAll(PETSC_NULL);CHKERRQ(ierr);
164
  ierr = STRegisterAll(PETSC_NULL);CHKERRQ(ierr);
165
 
166
#endif
167
 
168
  /*
169
      Register SLEPc events
170
  */
171
  info = SlepcRegisterEvents();CHKERRQ(info);
172
  SlepcInitializeCalled = PETSC_TRUE;
784 dsic.upv.es!antodo 173
  PetscLogInfo((0,"SlepcInitialize: SLEPc successfully started\n"));
6 dsic.upv.es!jroman 174
  PetscFunctionReturn(info);
175
}
176
 
177
#undef __FUNCT__  
178
#define __FUNCT__ "SlepcFinalize"
179
/*@
180
   SlepcFinalize - Checks for options to be called at the conclusion
181
   of the SLEPc program and calls PetscFinalize().
182
 
183
   Collective on PETSC_COMM_WORLD
184
 
185
   Level: beginner
186
 
187
.seealso: SlepcInitialize(), PetscFinalize()
188
@*/
476 dsic.upv.es!antodo 189
PetscErrorCode SlepcFinalize(void)
6 dsic.upv.es!jroman 190
{
572 dsic.upv.es!antodo 191
  PetscErrorCode info=0;
6 dsic.upv.es!jroman 192
 
193
  PetscFunctionBegin;
784 dsic.upv.es!antodo 194
  PetscLogInfo((0,"SlepcFinalize: SLEPc successfully ended!\n"));
6 dsic.upv.es!jroman 195
 
196
  if (SlepcBeganPetsc) {
197
    info = PetscFinalize();CHKERRQ(info);
198
  }
199
 
200
  SlepcInitializeCalled = PETSC_FALSE;
201
 
202
  PetscFunctionReturn(info);
203
}
204