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
1376 slepc 1
/*
2
       This file implements a wrapper to the BLOPEX solver
1199 slepc 3
 
1376 slepc 4
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1672 slepc 5
   SLEPc - Scalable Library for Eigenvalue Problem Computations
2575 eromero 6
   Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain
1376 slepc 7
 
1672 slepc 8
   This file is part of SLEPc.
9
 
10
   SLEPc is free software: you can redistribute it and/or modify it under  the
11
   terms of version 3 of the GNU Lesser General Public License as published by
12
   the Free Software Foundation.
13
 
14
   SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
15
   WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
16
   FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
17
   more details.
18
 
19
   You  should have received a copy of the GNU Lesser General  Public  License
20
   along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
1376 slepc 21
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1199 slepc 22
*/
1376 slepc 23
 
2729 jroman 24
#include <slepc-private/stimpl.h>       /*I "slepcst.h" I*/
25
#include <slepc-private/epsimpl.h>      /*I "slepceps.h" I*/
1499 slepc 26
#include "slepc-interface.h"
2423 jroman 27
#include <blopex_lobpcg.h>
28
#include <blopex_interpreter.h>
29
#include <blopex_multivector.h>
30
#include <blopex_temp_multivector.h>
1199 slepc 31
 
1947 jroman 32
PetscErrorCode EPSSolve_BLOPEX(EPS);
33
 
1199 slepc 34
typedef struct {
35
  lobpcg_Tolerance           tol;
36
  lobpcg_BLASLAPACKFunctions blap_fn;
37
  mv_MultiVectorPtr          eigenvectors;
1499 slepc 38
  mv_MultiVectorPtr          Y;
1199 slepc 39
  mv_InterfaceInterpreter    ii;
40
  KSP                        ksp;
2435 jroman 41
  Vec                        w;
1199 slepc 42
} EPS_BLOPEX;
43
 
44
#undef __FUNCT__  
45
#define __FUNCT__ "Precond_FnSingleVector"
46
static void Precond_FnSingleVector(void *data,void *x,void *y)
47
{
2317 jroman 48
  PetscErrorCode ierr;
49
  EPS            eps = (EPS)data;
50
  EPS_BLOPEX     *blopex = (EPS_BLOPEX*)eps->data;
51
  PetscInt       lits;
1199 slepc 52
 
53
  PetscFunctionBegin;
2730 jroman 54
  ierr = KSPSolve(blopex->ksp,(Vec)x,(Vec)y); CHKERRABORT(((PetscObject)eps)->comm,ierr);
55
  ierr = KSPGetIterationNumber(blopex->ksp,&lits); CHKERRABORT(((PetscObject)eps)->comm,ierr);
1499 slepc 56
  eps->OP->lineariterations+= lits;
1199 slepc 57
  PetscFunctionReturnVoid();
58
}
59
 
60
#undef __FUNCT__  
61
#define __FUNCT__ "Precond_FnMultiVector"
62
static void Precond_FnMultiVector(void *data,void *x,void *y)
63
{
2317 jroman 64
  EPS        eps = (EPS)data;
65
  EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data;
1199 slepc 66
 
67
  PetscFunctionBegin;
68
  blopex->ii.Eval(Precond_FnSingleVector,data,x,y);
69
  PetscFunctionReturnVoid();
70
}
71
 
72
#undef __FUNCT__  
73
#define __FUNCT__ "OperatorASingleVector"
74
static void OperatorASingleVector(void *data,void *x,void *y)
75
{
2316 jroman 76
  PetscErrorCode ierr;
77
  EPS            eps = (EPS)data;
2435 jroman 78
  EPS_BLOPEX     *blopex = (EPS_BLOPEX*)eps->data;
79
  Mat            A,B;
1499 slepc 80
 
1199 slepc 81
  PetscFunctionBegin;
2730 jroman 82
  ierr = STGetOperators(eps->OP,&A,&B);CHKERRABORT(((PetscObject)eps)->comm,ierr);
83
  ierr = MatMult(A,(Vec)x,(Vec)y);CHKERRABORT(((PetscObject)eps)->comm,ierr);
2435 jroman 84
  if (eps->OP->sigma != 0.0) {
2730 jroman 85
    if (B) { ierr = MatMult(B,(Vec)x,blopex->w);CHKERRABORT(((PetscObject)eps)->comm,ierr); }
86
    else { ierr = VecCopy((Vec)x,blopex->w);CHKERRABORT(((PetscObject)eps)->comm,ierr); }
87
    ierr = VecAXPY((Vec)y,-eps->OP->sigma,blopex->w);CHKERRABORT(((PetscObject)eps)->comm,ierr);
2435 jroman 88
  }
1199 slepc 89
  PetscFunctionReturnVoid();
90
}
91
 
92
#undef __FUNCT__  
93
#define __FUNCT__ "OperatorAMultiVector"
94
static void OperatorAMultiVector(void *data,void *x,void *y)
95
{
2316 jroman 96
  EPS        eps = (EPS)data;
97
  EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data;
1199 slepc 98
 
99
  PetscFunctionBegin;
100
  blopex->ii.Eval(OperatorASingleVector,data,x,y);
101
  PetscFunctionReturnVoid();
102
}
103
 
104
#undef __FUNCT__  
1499 slepc 105
#define __FUNCT__ "OperatorBSingleVector"
106
static void OperatorBSingleVector(void *data,void *x,void *y)
107
{
2316 jroman 108
  PetscErrorCode ierr;
109
  EPS            eps = (EPS)data;
110
  Mat            B;
1499 slepc 111
 
112
  PetscFunctionBegin;
2730 jroman 113
  ierr = STGetOperators(eps->OP,PETSC_NULL,&B);CHKERRABORT(((PetscObject)eps)->comm,ierr);
114
  ierr = MatMult(B,(Vec)x,(Vec)y);CHKERRABORT(((PetscObject)eps)->comm,ierr);
1499 slepc 115
  PetscFunctionReturnVoid();
116
}
117
 
118
#undef __FUNCT__  
119
#define __FUNCT__ "OperatorBMultiVector"
120
static void OperatorBMultiVector(void *data,void *x,void *y)
121
{
2316 jroman 122
  EPS        eps = (EPS)data;
123
  EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data;
1499 slepc 124
 
125
  PetscFunctionBegin;
126
  blopex->ii.Eval(OperatorBSingleVector,data,x,y);
127
  PetscFunctionReturnVoid();
128
}
129
 
130
#undef __FUNCT__  
1199 slepc 131
#define __FUNCT__ "EPSSetUp_BLOPEX"
132
PetscErrorCode EPSSetUp_BLOPEX(EPS eps)
133
{
2739 jroman 134
#if defined(PETSC_MISSING_LAPACK_POTRF) || defined(PETSC_MISSING_LAPACK_SYGV)
135
  PetscFunctionBegin;
136
  SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRF/SYGV - Lapack routine is unavailable.");
137
#else
2317 jroman 138
  PetscErrorCode ierr;
2372 jroman 139
  PetscInt       i;
2317 jroman 140
  EPS_BLOPEX     *blopex = (EPS_BLOPEX *)eps->data;
2461 jroman 141
  PetscBool      isPrecond;
1199 slepc 142
 
143
  PetscFunctionBegin;
2762 jroman 144
  if (!eps->ishermitian) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"blopex only works for hermitian problems");
1942 jroman 145
  if (!eps->which) eps->which = EPS_SMALLEST_REAL;
2762 jroman 146
  if (eps->which!=EPS_SMALLEST_REAL) SETERRQ(((PetscObject)eps)->comm,1,"Wrong value of eps->which");
1499 slepc 147
 
2104 eromero 148
  /* Change the default sigma to inf if necessary */
149
  if (eps->which == EPS_LARGEST_MAGNITUDE || eps->which == EPS_LARGEST_REAL ||
150
      eps->which == EPS_LARGEST_IMAGINARY) {
2331 jroman 151
    ierr = STSetDefaultShift(eps->OP,3e300);CHKERRQ(ierr);
2104 eromero 152
  }
153
 
2330 jroman 154
  ierr = STSetUp(eps->OP);CHKERRQ(ierr);
2823 jroman 155
  ierr = PetscObjectTypeCompare((PetscObject)eps->OP,STPRECOND,&isPrecond);CHKERRQ(ierr);
2331 jroman 156
  if (!isPrecond) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"blopex only works with STPRECOND");
157
  ierr = STGetKSP(eps->OP,&blopex->ksp);CHKERRQ(ierr);
1199 slepc 158
 
1928 jroman 159
  eps->ncv = eps->nev = PetscMin(eps->nev,eps->n);
2499 jroman 160
  if (eps->mpd) { ierr = PetscInfo(eps,"Warning: parameter mpd ignored\n");CHKERRQ(ierr); }
1928 jroman 161
  if (!eps->max_it) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
1499 slepc 162
 
163
  ierr = EPSAllocateSolution(eps);CHKERRQ(ierr);
164
  ierr = EPSDefaultGetWork(eps,1);CHKERRQ(ierr);
1199 slepc 165
 
2621 eromero 166
  blopex->tol.absolute = eps->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL:eps->tol;
1199 slepc 167
  blopex->tol.relative = 1e-50;
168
 
1499 slepc 169
  SLEPCSetupInterpreter(&blopex->ii);
170
  blopex->eigenvectors = mv_MultiVectorCreateFromSampleVector(&blopex->ii,eps->ncv,eps->V);
2372 jroman 171
  for (i=0;i<eps->ncv;i++) { ierr = PetscObjectReference((PetscObject)eps->V[i]);CHKERRQ(ierr); }
1199 slepc 172
  mv_MultiVectorSetRandom(blopex->eigenvectors,1234);
2435 jroman 173
  ierr = VecDuplicate(eps->V[0],&blopex->w);CHKERRQ(ierr);
1499 slepc 174
 
175
  if (eps->nds > 0) {
2331 jroman 176
    blopex->Y = mv_MultiVectorCreateFromSampleVector(&blopex->ii,eps->nds,eps->DS);
2372 jroman 177
    for (i=0;i<eps->nds;i++) { ierr = PetscObjectReference((PetscObject)eps->DS[i]);CHKERRQ(ierr); }
1499 slepc 178
  } else
179
    blopex->Y = PETSC_NULL;
1974 eromero 180
 
2320 jroman 181
#if defined(PETSC_USE_COMPLEX)
1974 eromero 182
  blopex->blap_fn.zpotrf = PETSC_zpotrf_interface;
183
  blopex->blap_fn.zhegv = PETSC_zsygv_interface;
184
#else
1199 slepc 185
  blopex->blap_fn.dpotrf = PETSC_dpotrf_interface;
186
  blopex->blap_fn.dsygv = PETSC_dsygv_interface;
1974 eromero 187
#endif
1199 slepc 188
 
2499 jroman 189
  if (eps->extraction) { ierr = PetscInfo(eps,"Warning: extraction type ignored\n");CHKERRQ(ierr); }
1426 slepc 190
 
1947 jroman 191
  /* dispatch solve method */
2214 jroman 192
  if (eps->leftvecs) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Left vectors not supported in this solver");
1947 jroman 193
  eps->ops->solve = EPSSolve_BLOPEX;
1199 slepc 194
  PetscFunctionReturn(0);
2739 jroman 195
#endif
1199 slepc 196
}
197
 
198
#undef __FUNCT__  
199
#define __FUNCT__ "EPSSolve_BLOPEX"
200
PetscErrorCode EPSSolve_BLOPEX(EPS eps)
201
{
2725 jroman 202
  EPS_BLOPEX     *blopex = (EPS_BLOPEX *)eps->data;
203
  int            i,j,info,its,nconv;
204
  double         *residhist=PETSC_NULL;
2587 jroman 205
  PetscErrorCode ierr;
2725 jroman 206
#if defined(PETSC_USE_COMPLEX)
207
  komplex        *lambdahist=PETSC_NULL;
208
#else
209
  double         *lambdahist=PETSC_NULL;
210
#endif
1199 slepc 211
 
212
  PetscFunctionBegin;
2587 jroman 213
  if (eps->numbermonitors>0) {
2725 jroman 214
#if defined(PETSC_USE_COMPLEX)
215
    ierr = PetscMalloc(eps->ncv*(eps->max_it+1)*sizeof(komplex),&lambdahist);CHKERRQ(ierr);
216
#else
2587 jroman 217
    ierr = PetscMalloc(eps->ncv*(eps->max_it+1)*sizeof(double),&lambdahist);CHKERRQ(ierr);
2725 jroman 218
#endif
2587 jroman 219
    ierr = PetscMalloc(eps->ncv*(eps->max_it+1)*sizeof(double),&residhist);CHKERRQ(ierr);
220
  }
221
 
2372 jroman 222
#if defined(PETSC_USE_COMPLEX)
223
  info = lobpcg_solve_complex(blopex->eigenvectors,eps,OperatorAMultiVector,
2316 jroman 224
        eps->isgeneralized?eps:PETSC_NULL,eps->isgeneralized?OperatorBMultiVector:PETSC_NULL,
1499 slepc 225
        eps,Precond_FnMultiVector,blopex->Y,
1524 slepc 226
        blopex->blap_fn,blopex->tol,eps->max_it,0,&its,
2587 jroman 227
        (komplex*)eps->eigr,lambdahist,eps->ncv,eps->errest,residhist,eps->ncv);
2372 jroman 228
#else
229
  info = lobpcg_solve_double(blopex->eigenvectors,eps,OperatorAMultiVector,
230
        eps->isgeneralized?eps:PETSC_NULL,eps->isgeneralized?OperatorBMultiVector:PETSC_NULL,
231
        eps,Precond_FnMultiVector,blopex->Y,
232
        blopex->blap_fn,blopex->tol,eps->max_it,0,&its,
2587 jroman 233
        eps->eigr,lambdahist,eps->ncv,eps->errest,residhist,eps->ncv);
2372 jroman 234
#endif
2214 jroman 235
  if (info>0) SETERRQ1(((PetscObject)eps)->comm,PETSC_ERR_LIB,"Error in blopex (code=%d)",info);
1199 slepc 236
 
2587 jroman 237
  if (eps->numbermonitors>0) {
238
    for (i=0;i<its;i++) {
239
      nconv = 0;
240
      for (j=0;j<eps->ncv;j++) { if (residhist[j+i*eps->ncv]>eps->tol) break; else nconv++; }
2725 jroman 241
      ierr = EPSMonitor(eps,i,nconv,(PetscScalar*)lambdahist+i*eps->ncv,eps->eigi,residhist+i*eps->ncv,eps->ncv);CHKERRQ(ierr);
2587 jroman 242
    }
243
    ierr = PetscFree(lambdahist);CHKERRQ(ierr);
244
    ierr = PetscFree(residhist);CHKERRQ(ierr);
245
  }
246
 
1509 slepc 247
  eps->its = its;
1199 slepc 248
  eps->nconv = eps->ncv;
2435 jroman 249
  if (eps->OP->sigma != 0.0) {
250
    for (i=0;i<eps->nconv;i++) eps->eigr[i]+=eps->OP->sigma;
251
  }
1199 slepc 252
  if (info==-1) eps->reason = EPS_DIVERGED_ITS;
253
  else eps->reason = EPS_CONVERGED_TOL;
254
  PetscFunctionReturn(0);
255
}
256
 
257
#undef __FUNCT__  
2348 jroman 258
#define __FUNCT__ "EPSReset_BLOPEX"
259
PetscErrorCode EPSReset_BLOPEX(EPS eps)
260
{
261
  PetscErrorCode ierr;
262
  EPS_BLOPEX     *blopex = (EPS_BLOPEX *)eps->data;
263
 
264
  PetscFunctionBegin;
265
  mv_MultiVectorDestroy(blopex->eigenvectors);
266
  mv_MultiVectorDestroy(blopex->Y);
2435 jroman 267
  ierr = VecDestroy(&blopex->w);CHKERRQ(ierr);
2348 jroman 268
  ierr = EPSReset_Default(eps);CHKERRQ(ierr);
269
  PetscFunctionReturn(0);
270
}
271
 
272
#undef __FUNCT__  
1199 slepc 273
#define __FUNCT__ "EPSDestroy_BLOPEX"
274
PetscErrorCode EPSDestroy_BLOPEX(EPS eps)
275
{
276
  PetscErrorCode ierr;
277
 
278
  PetscFunctionBegin;
279
  LOBPCG_DestroyRandomContext();
280
  ierr = PetscFree(eps->data);CHKERRQ(ierr);
281
  PetscFunctionReturn(0);
282
}
283
 
2448 eromero 284
#undef __FUNCT__  
285
#define __FUNCT__ "EPSSetFromOptions_BLOPEX"
286
PetscErrorCode EPSSetFromOptions_BLOPEX(EPS eps)
287
{
288
  PetscErrorCode  ierr;
289
 
290
  PetscFunctionBegin;
291
  ierr = PetscOptionsHead("EPS BLOPEX Options");CHKERRQ(ierr);
292
  LOBPCG_SetFromOptionsRandomContext();
293
  ierr = PetscOptionsTail();CHKERRQ(ierr);
294
  PetscFunctionReturn(0);
295
}
296
 
1199 slepc 297
EXTERN_C_BEGIN
298
#undef __FUNCT__  
299
#define __FUNCT__ "EPSCreate_BLOPEX"
300
PetscErrorCode EPSCreate_BLOPEX(EPS eps)
301
{
302
  PetscErrorCode ierr;
303
 
304
  PetscFunctionBegin;
2372 jroman 305
  ierr = PetscNewLog(eps,EPS_BLOPEX,&eps->data);CHKERRQ(ierr);
1199 slepc 306
  eps->ops->setup                = EPSSetUp_BLOPEX;
2448 eromero 307
  eps->ops->setfromoptions       = EPSSetFromOptions_BLOPEX;
1199 slepc 308
  eps->ops->destroy              = EPSDestroy_BLOPEX;
2348 jroman 309
  eps->ops->reset                = EPSReset_BLOPEX;
1199 slepc 310
  eps->ops->backtransform        = EPSBackTransform_Default;
311
  eps->ops->computevectors       = EPSComputeVectors_Default;
2348 jroman 312
  ierr = STSetType(eps->OP,STPRECOND);CHKERRQ(ierr);
313
  ierr = STPrecondSetKSPHasMat(eps->OP,PETSC_TRUE);CHKERRQ(ierr);
2695 jroman 314
  LOBPCG_InitRandomContext(((PetscObject)eps)->comm,eps->rand);
1199 slepc 315
  PetscFunctionReturn(0);
316
}
317
EXTERN_C_END