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"
756 dsic.upv.es!antodo 58
PetscEvent EPS_SetUp, EPS_Solve, ST_SetUp, ST_Apply, ST_ApplyB, ST_ApplyNoB, 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);
756 dsic.upv.es!antodo 75
  info = PetscLogEventRegister(&EPS_Orthogonalize,"EPSOrthogonalize",PETSC_NULL); CHKERRQ(info);
382 dsic.upv.es!antodo 76
  info = PetscLogEventRegister(&ST_InnerProduct,"STInnerProduct",PETSC_NULL); CHKERRQ(info);
756 dsic.upv.es!antodo 77
  info = PetscLogEventRegister(&EPS_ReverseProjection,"EPSReverseProjection",PETSC_NULL); CHKERRQ(info);
6 dsic.upv.es!jroman 78
 
79
  PetscFunctionReturn(info);
80
}
81
 
82
/* ------------------------Nasty global variables -------------------------------*/
83
/*
84
   Indicates whether SLEPc started PETSc, or whether it was
85
   already started before SLEPc was initialized.
86
*/
87
PetscTruth  SlepcBeganPetsc = PETSC_FALSE;
88
PetscTruth  SlepcInitializeCalled = PETSC_FALSE;
476 dsic.upv.es!antodo 89
PetscCookie EPS_COOKIE = 0;
90
PetscCookie ST_COOKIE = 0;
6 dsic.upv.es!jroman 91
 
92
#if defined(PETSC_USE_DYNAMIC_LIBRARIES)
93
extern PetscDLLibraryList DLLibrariesLoaded;
94
#endif
95
 
96
#undef __FUNCT__  
97
#define __FUNCT__ "SlepcInitialize"
98
/*@C
99
   SlepcInitialize - Initializes the SLEPc library. SlepcInitialize() calls
100
   PetscInitialize() if that has not been called yet, so this routine should
101
   always be called near the beginning of your program.
102
 
103
   Collective on MPI_COMM_WORLD or PETSC_COMM_WORLD if it has been set
104
 
105
   Input Parameters:
106
+  argc - count of number of command line arguments
107
.  args - the command line arguments
108
.  file - [optional] PETSc database file, defaults to ~username/.petscrc
109
          (use PETSC_NULL for default)
110
-  help - [optional] Help message to print, use PETSC_NULL for no message
111
 
112
   Fortran Note:
113
   Fortran syntax is very similar to that of PetscInitialize()
114
 
115
   Level: beginner
116
 
117
.seealso: SlepcInitializeFortran(), SlepcFinalize(), PetscInitialize()
118
@*/
476 dsic.upv.es!antodo 119
PetscErrorCode SlepcInitialize(int *argc,char ***args,char file[],const char help[])
6 dsic.upv.es!jroman 120
{
476 dsic.upv.es!antodo 121
  PetscErrorCode ierr,info=0;
373 dsic.upv.es!antodo 122
#if defined(PETSC_USE_DYNAMIC_LIBRARIES)
476 dsic.upv.es!antodo 123
  char           libs[PETSC_MAX_PATH_LEN],dlib[PETSC_MAX_PATH_LEN];
124
  PetscTruth     found;
373 dsic.upv.es!antodo 125
#endif
6 dsic.upv.es!jroman 126
 
127
  PetscFunctionBegin;
128
 
129
  if (SlepcInitializeCalled==PETSC_TRUE) {
130
    PetscFunctionReturn(0);
131
  }
132
 
133
#if !defined(PARCH_t3d)
134
  info = PetscSetHelpVersionFunctions(SlepcPrintHelpIntro,SlepcPrintVersion);CHKERRQ(info);
135
#endif
136
 
137
  if (!PetscInitializeCalled) {
138
    info = PetscInitialize(argc,args,file,help);CHKERRQ(info);
139
    SlepcBeganPetsc = PETSC_TRUE;
140
  }
141
 
142
  EPS_COOKIE = 0;
143
  ierr = PetscLogClassRegister(&EPS_COOKIE,"Eigenproblem Solver");CHKERRQ(ierr);
144
  ST_COOKIE = 0;
145
  ierr = PetscLogClassRegister(&ST_COOKIE,"Spectral Transform");CHKERRQ(ierr);
146
 
147
  /*
148
      Load the dynamic libraries
149
  */
150
 
151
#if defined(PETSC_USE_DYNAMIC_LIBRARIES)
152
  ierr = PetscStrcpy(libs,SLEPC_LIB_DIR);CHKERRQ(ierr);
153
  ierr = PetscStrcat(libs,"/libslepc");CHKERRQ(ierr);
154
  ierr = PetscDLLibraryRetrieve(PETSC_COMM_WORLD,libs,dlib,1024,&found);CHKERRQ(ierr);
155
  if (found) {
156
    ierr = PetscDLLibraryAppend(PETSC_COMM_WORLD,&DLLibrariesLoaded,libs);CHKERRQ(ierr);
157
  } else {
158
    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);
159
  }
160
#else
161
 
162
  ierr = EPSRegisterAll(PETSC_NULL);CHKERRQ(ierr);
163
  ierr = STRegisterAll(PETSC_NULL);CHKERRQ(ierr);
164
 
165
#endif
166
 
167
  /*
168
      Register SLEPc events
169
  */
170
  info = SlepcRegisterEvents();CHKERRQ(info);
171
  SlepcInitializeCalled = PETSC_TRUE;
172
  PetscLogInfo(0,"SlepcInitialize: SLEPc successfully started\n");
173
  PetscFunctionReturn(info);
174
}
175
 
176
#undef __FUNCT__  
177
#define __FUNCT__ "SlepcFinalize"
178
/*@
179
   SlepcFinalize - Checks for options to be called at the conclusion
180
   of the SLEPc program and calls PetscFinalize().
181
 
182
   Collective on PETSC_COMM_WORLD
183
 
184
   Level: beginner
185
 
186
.seealso: SlepcInitialize(), PetscFinalize()
187
@*/
476 dsic.upv.es!antodo 188
PetscErrorCode SlepcFinalize(void)
6 dsic.upv.es!jroman 189
{
572 dsic.upv.es!antodo 190
  PetscErrorCode info=0;
6 dsic.upv.es!jroman 191
 
192
  PetscFunctionBegin;
193
  PetscLogInfo(0,"SlepcFinalize: SLEPc successfully ended!\n");
194
 
195
  if (SlepcBeganPetsc) {
196
    info = PetscFinalize();CHKERRQ(info);
197
  }
198
 
199
  SlepcInitializeCalled = PETSC_FALSE;
200
 
201
  PetscFunctionReturn(info);
202
}
203