Subversion Repositories slepc-dev

Rev

Details | Compare with Previous | Last modification | View Log | RSS feed

Rev Author Line No. Line
1376 slepc 1
/*
1187 slepc 2
       This file implements a wrapper to the PRIMME library
1376 slepc 3
 
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
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1187 slepc 22
*/
23
 
2283 jroman 24
#include <petscsys.h>
25
#include <private/epsimpl.h>    /*I "slepceps.h" I*/
2654 jroman 26
#include <private/stimpl.h>
1187 slepc 27
 
1947 jroman 28
PetscErrorCode EPSSolve_PRIMME(EPS);
29
 
1187 slepc 30
EXTERN_C_BEGIN
2283 jroman 31
#include <primme.h>
1187 slepc 32
EXTERN_C_END
33
 
34
typedef struct {
35
  primme_params primme;           /* param struc */
36
  primme_preset_method method;    /* primme method */
2718 jroman 37
  Mat       A;                    /* problem matrix */
38
  EPS       eps;                  /* EPS current context */
39
  KSP       ksp;                  /* linear solver and preconditioner */
40
  Vec       x,y;                  /* auxiliary vectors */
41
  PetscReal target;               /* a copy of eps's target */
1187 slepc 42
} EPS_PRIMME;
43
 
1239 slepc 44
EPSPRIMMEMethod methodN[] = {
1940 jroman 45
  EPS_PRIMME_DYNAMIC,
46
  EPS_PRIMME_DEFAULT_MIN_TIME,
47
  EPS_PRIMME_DEFAULT_MIN_MATVECS,
48
  EPS_PRIMME_ARNOLDI,
49
  EPS_PRIMME_GD,
50
  EPS_PRIMME_GD_PLUSK,
51
  EPS_PRIMME_GD_OLSEN_PLUSK,
52
  EPS_PRIMME_JD_OLSEN_PLUSK,
53
  EPS_PRIMME_RQI,
54
  EPS_PRIMME_JDQR,
55
  EPS_PRIMME_JDQMR,
56
  EPS_PRIMME_JDQMR_ETOL,
57
  EPS_PRIMME_SUBSPACE_ITERATION,
58
  EPS_PRIMME_LOBPCG_ORTHOBASIS,
59
  EPS_PRIMME_LOBPCG_ORTHOBASISW
1239 slepc 60
};
1187 slepc 61
 
2331 jroman 62
static void multMatvec_PRIMME(void *in,void *out,int *blockSize,primme_params *primme);
63
static void applyPreconditioner_PRIMME(void *in,void *out,int *blockSize,struct primme_params *primme);
1239 slepc 64
 
2331 jroman 65
void par_GlobalSumDouble(void *sendBuf,void *recvBuf,int *count,primme_params *primme)
2317 jroman 66
{
1472 slepc 67
  PetscErrorCode ierr;
2331 jroman 68
  ierr = MPI_Allreduce((double*)sendBuf,(double*)recvBuf,*count,MPI_DOUBLE,MPI_SUM,((PetscObject)(primme->commInfo))->comm);CHKERRABORT(((PetscObject)(primme->commInfo))->comm,ierr);
1187 slepc 69
}
70
 
71
#undef __FUNCT__  
72
#define __FUNCT__ "EPSSetUp_PRIMME"
73
PetscErrorCode EPSSetUp_PRIMME(EPS eps)
74
{
75
  PetscErrorCode ierr;
2331 jroman 76
  PetscMPIInt    numProcs,procID;
2718 jroman 77
  EPS_PRIMME     *ops = (EPS_PRIMME*)eps->data;
78
  primme_params  *primme = &ops->primme;
2216 jroman 79
  PetscBool      t;
1187 slepc 80
 
81
  PetscFunctionBegin;
1472 slepc 82
  ierr = MPI_Comm_size(((PetscObject)eps)->comm,&numProcs);CHKERRQ(ierr);
83
  ierr = MPI_Comm_rank(((PetscObject)eps)->comm,&procID);CHKERRQ(ierr);
1187 slepc 84
 
1212 slepc 85
  /* Check some constraints and set some default values */
1928 jroman 86
  if (!eps->max_it) eps->max_it = PetscMax(1000,eps->n);
2331 jroman 87
  ierr = STGetOperators(eps->OP,&ops->A,PETSC_NULL);
2214 jroman 88
  if (!ops->A) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_WRONGSTATE,"The problem matrix has to be specified first");
1187 slepc 89
  if (!eps->ishermitian)
2214 jroman 90
    SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"PRIMME is only available for Hermitian problems");
1189 slepc 91
  if (eps->isgeneralized)
2214 jroman 92
    SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"PRIMME is not available for generalized problems");
2104 eromero 93
  if (!eps->which) eps->which = EPS_LARGEST_REAL;
94
 
95
  /* Change the default sigma to inf if necessary */
96
  if (eps->which == EPS_LARGEST_MAGNITUDE || eps->which == EPS_LARGEST_REAL ||
97
      eps->which == EPS_LARGEST_IMAGINARY) {
2331 jroman 98
    ierr = STSetDefaultShift(eps->OP,3e300);CHKERRQ(ierr);
2104 eromero 99
  }
100
 
2485 eromero 101
  /* Avoid setting the automatic shift when a target is set */
102
  ierr = STSetDefaultShift(eps->OP,0.0);CHKERRQ(ierr);
103
 
2330 jroman 104
  ierr = STSetUp(eps->OP);CHKERRQ(ierr);
2331 jroman 105
  ierr = PetscTypeCompare((PetscObject)eps->OP,STPRECOND,&t);CHKERRQ(ierr);
106
  if (!t) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"PRIMME only works with STPRECOND");
1187 slepc 107
 
108
  /* Transfer SLEPc options to PRIMME options */
1928 jroman 109
  primme->n = eps->n;
110
  primme->nLocal = eps->nloc;
1187 slepc 111
  primme->numEvals = eps->nev;
1189 slepc 112
  primme->matrix = ops;
1212 slepc 113
  primme->commInfo = eps;
1187 slepc 114
  primme->maxMatvecs = eps->max_it;
2621 eromero 115
  primme->eps = eps->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL:eps->tol;
1187 slepc 116
  primme->numProcs = numProcs;
117
  primme->procID = procID;
1700 eromero 118
  primme->printLevel = 0;
2001 eromero 119
  primme->correctionParams.precondition = 1;
1187 slepc 120
 
1942 jroman 121
  if (!eps->which) eps->which = EPS_LARGEST_REAL;
1187 slepc 122
  switch(eps->which) {
1203 slepc 123
    case EPS_LARGEST_REAL:
1187 slepc 124
      primme->target = primme_largest;
125
      break;
1203 slepc 126
    case EPS_SMALLEST_REAL:
1187 slepc 127
      primme->target = primme_smallest;
128
      break;
2485 eromero 129
    case EPS_TARGET_MAGNITUDE:
130
    case EPS_TARGET_REAL:
131
      primme->target = primme_closest_abs;
132
      primme->numTargetShifts = 1;
2718 jroman 133
      ops->target = PetscRealPart(eps->target);
134
      primme->targetShifts = &ops->target;
2485 eromero 135
      break;
1187 slepc 136
    default:
2485 eromero 137
      SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"'which' value does not supported by PRIMME");
1187 slepc 138
      break;  
139
  }
1207 slepc 140
 
2331 jroman 141
  if (primme_set_method(ops->method,primme) < 0)
2214 jroman 142
    SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"PRIMME method not valid");
1187 slepc 143
 
1207 slepc 144
  /* If user sets ncv, maxBasisSize is modified. If not, ncv is set as maxBasisSize */
145
  if (eps->ncv) primme->maxBasisSize = eps->ncv;
146
  else eps->ncv = primme->maxBasisSize;
147
  if (eps->ncv < eps->nev+primme->maxBlockSize)  
2214 jroman 148
    SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"PRIMME needs ncv >= nev+maxBlockSize");
2499 jroman 149
  if (eps->mpd) { ierr = PetscInfo(eps,"Warning: parameter mpd ignored\n");CHKERRQ(ierr); }
1207 slepc 150
 
2499 jroman 151
  if (eps->extraction) { ierr = PetscInfo(eps,"Warning: extraction type ignored\n");CHKERRQ(ierr); }
1426 slepc 152
 
1187 slepc 153
  /* Set workspace */
1596 slepc 154
  ierr = EPSAllocateSolution(eps);CHKERRQ(ierr);
1207 slepc 155
 
2001 eromero 156
  /* Setup the preconditioner */
157
  ops->eps = eps;
1207 slepc 158
  if (primme->correctionParams.precondition) {
2331 jroman 159
    ierr = STGetKSP(eps->OP,&ops->ksp);CHKERRQ(ierr);
160
    ierr = PetscTypeCompare((PetscObject)ops->ksp,KSPPREONLY,&t);CHKERRQ(ierr);
161
    if (!t) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"PRIMME only works with KSPPREONLY");
1189 slepc 162
    primme->preconditioner = PETSC_NULL;
1187 slepc 163
    primme->applyPreconditioner = applyPreconditioner_PRIMME;
164
  }
165
 
1212 slepc 166
  /* Prepare auxiliary vectors */
1928 jroman 167
  ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,eps->nloc,eps->n,PETSC_NULL,&ops->x);CHKERRQ(ierr);
168
  ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,eps->nloc,eps->n,PETSC_NULL,&ops->y);CHKERRQ(ierr);
1189 slepc 169
 
1947 jroman 170
  /* dispatch solve method */
2214 jroman 171
  if (eps->leftvecs) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Left vectors not supported in this solver");
1947 jroman 172
  eps->ops->solve = EPSSolve_PRIMME;
1187 slepc 173
  PetscFunctionReturn(0);
174
}
175
 
176
#undef __FUNCT__  
177
#define __FUNCT__ "EPSSolve_PRIMME"
178
PetscErrorCode EPSSolve_PRIMME(EPS eps)
179
{
180
  PetscErrorCode ierr;
1201 slepc 181
  EPS_PRIMME     *ops = (EPS_PRIMME *)eps->data;
182
  PetscScalar    *a;
2320 jroman 183
#if defined(PETSC_USE_COMPLEX)
1509 slepc 184
  PetscInt       i;
1201 slepc 185
  PetscReal      *evals;
1187 slepc 186
#endif
187
 
188
  PetscFunctionBegin;
1240 slepc 189
  /* Reset some parameters left from previous runs */
190
  ops->primme.aNorm    = 0.0;
1933 jroman 191
  ops->primme.initSize = eps->nini;
1240 slepc 192
  ops->primme.iseed[0] = -1;
193
 
1187 slepc 194
  /* Call PRIMME solver */
2331 jroman 195
  ierr = VecGetArray(eps->V[0],&a);CHKERRQ(ierr);
2320 jroman 196
#if !defined(PETSC_USE_COMPLEX)
2331 jroman 197
  ierr = dprimme(eps->eigr,a,eps->errest,&ops->primme);
1187 slepc 198
#else
1212 slepc 199
  /* PRIMME returns real eigenvalues, but SLEPc works with complex ones */
2330 jroman 200
  ierr = PetscMalloc(eps->ncv*sizeof(PetscReal),&evals);CHKERRQ(ierr);
2331 jroman 201
  ierr = zprimme(evals,(Complex_Z*)a,eps->errest,&ops->primme);
1201 slepc 202
  for (i=0;i<eps->ncv;i++)
203
    eps->eigr[i] = evals[i];
2330 jroman 204
  ierr = PetscFree(evals);CHKERRQ(ierr);
1187 slepc 205
#endif
2331 jroman 206
  ierr = VecRestoreArray(eps->V[0],&a);CHKERRQ(ierr);
1187 slepc 207
 
208
  switch(ierr) {
209
    case 0: /* Successful */
210
      break;
211
 
212
    case -1:
2214 jroman 213
      SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"PRIMME: Failed to open output file");
1187 slepc 214
      break;
215
 
216
    case -2:
2214 jroman 217
      SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"PRIMME: Insufficient integer or real workspace allocated");
1187 slepc 218
      break;
219
 
220
    case -3:
2214 jroman 221
      SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"PRIMME: main_iter encountered a problem");
1187 slepc 222
      break;
223
 
224
    default:
2214 jroman 225
      SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"PRIMME: some parameters wrong configured");
1187 slepc 226
      break;
227
  }
228
 
1240 slepc 229
  eps->nconv = ops->primme.initSize >= 0 ? ops->primme.initSize : 0;
230
  eps->reason = eps->ncv >= eps->nev ? EPS_CONVERGED_TOL : EPS_DIVERGED_ITS;
1236 slepc 231
  eps->its = ops->primme.stats.numOuterIterations;
232
  eps->OP->applys = ops->primme.stats.numMatvecs;
1187 slepc 233
  PetscFunctionReturn(0);
234
}
235
 
1201 slepc 236
#undef __FUNCT__  
237
#define __FUNCT__ "multMatvec_PRIMME"
2331 jroman 238
static void multMatvec_PRIMME(void *in,void *out,int *blockSize,primme_params *primme)
1187 slepc 239
{
1212 slepc 240
  PetscErrorCode ierr;
2331 jroman 241
  PetscInt       i,N = primme->n;
1212 slepc 242
  EPS_PRIMME     *ops = (EPS_PRIMME *)primme->matrix;
2331 jroman 243
  Vec            x = ops->x,y = ops->y;
1212 slepc 244
  Mat            A = ops->A;
1187 slepc 245
 
1201 slepc 246
  PetscFunctionBegin;
1187 slepc 247
  for (i=0;i<*blockSize;i++) {
248
    /* build vectors using 'in' an 'out' workspace */
2331 jroman 249
    ierr = VecPlaceArray(x,(PetscScalar*)in+N*i);CHKERRABORT(PETSC_COMM_WORLD,ierr);
250
    ierr = VecPlaceArray(y,(PetscScalar*)out+N*i);CHKERRABORT(PETSC_COMM_WORLD,ierr);
1187 slepc 251
 
2331 jroman 252
    ierr = MatMult(A,x,y);CHKERRABORT(PETSC_COMM_WORLD,ierr);
1187 slepc 253
 
2331 jroman 254
    ierr = VecResetArray(x);CHKERRABORT(PETSC_COMM_WORLD,ierr);
255
    ierr = VecResetArray(y);CHKERRABORT(PETSC_COMM_WORLD,ierr);
1187 slepc 256
  }
1201 slepc 257
  PetscFunctionReturnVoid();
1187 slepc 258
}
259
 
1201 slepc 260
#undef __FUNCT__  
261
#define __FUNCT__ "applyPreconditioner_PRIMME"
2331 jroman 262
static void applyPreconditioner_PRIMME(void *in,void *out,int *blockSize,struct primme_params *primme)
1187 slepc 263
{
1212 slepc 264
  PetscErrorCode ierr;
2331 jroman 265
  PetscInt       i,N = primme->n,lits;
1212 slepc 266
  EPS_PRIMME     *ops = (EPS_PRIMME *)primme->matrix;
2331 jroman 267
  Vec            x = ops->x,y = ops->y;
1189 slepc 268
 
1201 slepc 269
  PetscFunctionBegin;
1187 slepc 270
  for (i=0;i<*blockSize;i++) {
271
    /* build vectors using 'in' an 'out' workspace */
2331 jroman 272
    ierr = VecPlaceArray(x,(PetscScalar*)in+N*i);CHKERRABORT(PETSC_COMM_WORLD,ierr);
273
    ierr = VecPlaceArray(y,(PetscScalar*)out+N*i);CHKERRABORT(PETSC_COMM_WORLD,ierr);
1187 slepc 274
 
2331 jroman 275
    ierr = KSPSolve(ops->ksp,x,y);CHKERRABORT(PETSC_COMM_WORLD,ierr);
276
    ierr = KSPGetIterationNumber(ops->ksp,&lits);CHKERRABORT(PETSC_COMM_WORLD,ierr);
2001 eromero 277
    ops->eps->OP->lineariterations+= lits;
1187 slepc 278
 
2331 jroman 279
    ierr = VecResetArray(x);CHKERRABORT(PETSC_COMM_WORLD,ierr);
280
    ierr = VecResetArray(y);CHKERRABORT(PETSC_COMM_WORLD,ierr);
1187 slepc 281
  }
1201 slepc 282
  PetscFunctionReturnVoid();
1187 slepc 283
}
284
 
285
#undef __FUNCT__  
2348 jroman 286
#define __FUNCT__ "EPSReset_PRIMME"
2395 jroman 287
PetscErrorCode EPSReset_PRIMME(EPS eps)
1187 slepc 288
{
289
  PetscErrorCode ierr;
2317 jroman 290
  EPS_PRIMME     *ops = (EPS_PRIMME *)eps->data;
1187 slepc 291
 
292
  PetscFunctionBegin;
293
  primme_Free(&ops->primme);
2305 jroman 294
  ierr = VecDestroy(&ops->x);CHKERRQ(ierr);
295
  ierr = VecDestroy(&ops->y);CHKERRQ(ierr);
2348 jroman 296
  ierr = EPSFreeSolution(eps);CHKERRQ(ierr);
297
  PetscFunctionReturn(0);
298
}
299
 
300
#undef __FUNCT__  
301
#define __FUNCT__ "EPSDestroy_PRIMME"
302
PetscErrorCode EPSDestroy_PRIMME(EPS eps)
303
{
304
  PetscErrorCode ierr;
305
 
306
  PetscFunctionBegin;
1218 slepc 307
  ierr = PetscFree(eps->data);CHKERRQ(ierr);
1925 jroman 308
  ierr = PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPRIMMESetBlockSize_C","",PETSC_NULL);CHKERRQ(ierr);
309
  ierr = PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPRIMMESetMethod_C","",PETSC_NULL);CHKERRQ(ierr);
310
  ierr = PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPRIMMEGetBlockSize_C","",PETSC_NULL);CHKERRQ(ierr);
311
  ierr = PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPRIMMEGetMethod_C","",PETSC_NULL);CHKERRQ(ierr);
1187 slepc 312
  PetscFunctionReturn(0);
313
}
314
 
315
#undef __FUNCT__  
316
#define __FUNCT__ "EPSView_PRIMME"
317
PetscErrorCode EPSView_PRIMME(EPS eps,PetscViewer viewer)
318
{
2216 jroman 319
  PetscErrorCode  ierr;
320
  PetscBool       isascii;
321
  primme_params   *primme = &((EPS_PRIMME *)eps->data)->primme;
1239 slepc 322
  EPSPRIMMEMethod methodn;
2482 jroman 323
  PetscMPIInt     rank;
1239 slepc 324
 
1187 slepc 325
  PetscFunctionBegin;
2215 jroman 326
  ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
1189 slepc 327
  if (!isascii) {
2214 jroman 328
    SETERRQ1(((PetscObject)eps)->comm,1,"Viewer type %s not supported for EPSPRIMME",((PetscObject)viewer)->type_name);
1189 slepc 329
  }
330
 
2365 jroman 331
  ierr = PetscViewerASCIIPrintf(viewer,"  PRIMME: block size=%d\n",primme->maxBlockSize);CHKERRQ(ierr);
2322 jroman 332
  ierr = EPSPRIMMEGetMethod(eps,&methodn);CHKERRQ(ierr);
2365 jroman 333
  ierr = PetscViewerASCIIPrintf(viewer,"  PRIMME: solver method: %s\n",EPSPRIMMEMethods[methodn]);CHKERRQ(ierr);
1189 slepc 334
 
1700 eromero 335
  /* Display PRIMME params */
2482 jroman 336
  ierr = MPI_Comm_rank(((PetscObject)eps)->comm,&rank);CHKERRQ(ierr);
337
  if (!rank) primme_display_params(*primme);
1187 slepc 338
  PetscFunctionReturn(0);
339
}
340
 
341
#undef __FUNCT__  
342
#define __FUNCT__ "EPSSetFromOptions_PRIMME"
343
PetscErrorCode EPSSetFromOptions_PRIMME(EPS eps)
344
{
2322 jroman 345
  PetscErrorCode  ierr;
346
  EPS_PRIMME      *ctx = (EPS_PRIMME *)eps->data;
347
  PetscInt        bs;
348
  EPSPRIMMEMethod meth;
349
  PetscBool       flg;
1187 slepc 350
 
351
  PetscFunctionBegin;
2384 jroman 352
  ierr = PetscOptionsHead("EPS PRIMME Options");CHKERRQ(ierr);
2322 jroman 353
  ierr = PetscOptionsInt("-eps_primme_block_size","Maximum block size","EPSPRIMMESetBlockSize",ctx->primme.maxBlockSize,&bs,&flg);CHKERRQ(ierr);
354
  if (flg) { ierr = EPSPRIMMESetBlockSize(eps,bs);CHKERRQ(ierr); }
2330 jroman 355
  ierr = PetscOptionsEnum("-eps_primme_method","Method for solving the eigenproblem","EPSPRIMMESetMethod",EPSPRIMMEMethods,(PetscEnum)ctx->method,(PetscEnum*)&meth,&flg);CHKERRQ(ierr);
2322 jroman 356
  if (flg) { ierr = EPSPRIMMESetMethod(eps,meth);CHKERRQ(ierr); }
2384 jroman 357
  ierr = PetscOptionsTail();CHKERRQ(ierr);
1187 slepc 358
  PetscFunctionReturn(0);
359
}
360
 
361
EXTERN_C_BEGIN
362
#undef __FUNCT__  
363
#define __FUNCT__ "EPSPRIMMESetBlockSize_PRIMME"
1509 slepc 364
PetscErrorCode EPSPRIMMESetBlockSize_PRIMME(EPS eps,PetscInt bs)
1187 slepc 365
{
2331 jroman 366
  EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
1187 slepc 367
 
368
  PetscFunctionBegin;
369
  if (bs == PETSC_DEFAULT) ops->primme.maxBlockSize = 1;
370
  else if (bs <= 0) {
2331 jroman 371
    SETERRQ(((PetscObject)eps)->comm,1,"PRIMME: wrong block size");
1187 slepc 372
  } else ops->primme.maxBlockSize = bs;
373
  PetscFunctionReturn(0);
374
}
375
EXTERN_C_END
376
 
377
#undef __FUNCT__  
378
#define __FUNCT__ "EPSPRIMMESetBlockSize"
1201 slepc 379
/*@
2317 jroman 380
   EPSPRIMMESetBlockSize - The maximum block size the code will try to use.
381
   The user should set
382
   this based on the architecture specifics of the target computer,
383
   as well as any a priori knowledge of multiplicities. The code does
384
   NOT require BlockSize > 1 to find multiple eigenvalues.  For some
385
   methods, keeping BlockSize = 1 yields the best overall performance.
1187 slepc 386
 
387
   Collective on EPS
388
 
389
   Input Parameters:
390
+  eps - the eigenproblem solver context
391
-  bs - block size
392
 
393
   Options Database Key:
1212 slepc 394
.  -eps_primme_block_size - Sets the max allowed block size value
1187 slepc 395
 
1212 slepc 396
   Notes:
397
   If the block size is not set, the value established by primme_initialize
398
   is used.
1187 slepc 399
 
400
   Level: advanced
401
.seealso: EPSPRIMMEGetBlockSize()
402
@*/
1509 slepc 403
PetscErrorCode EPSPRIMMESetBlockSize(EPS eps,PetscInt bs)
1187 slepc 404
{
2221 jroman 405
  PetscErrorCode ierr;
1187 slepc 406
 
407
  PetscFunctionBegin;
2213 jroman 408
  PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
2326 jroman 409
  PetscValidLogicalCollectiveInt(eps,bs,2);
2221 jroman 410
  ierr = PetscTryMethod(eps,"EPSPRIMMESetBlockSize_C",(EPS,PetscInt),(eps,bs));CHKERRQ(ierr);
1187 slepc 411
  PetscFunctionReturn(0);
412
}
413
 
414
EXTERN_C_BEGIN
415
#undef __FUNCT__  
416
#define __FUNCT__ "EPSPRIMMEGetBlockSize_PRIMME"
1509 slepc 417
PetscErrorCode EPSPRIMMEGetBlockSize_PRIMME(EPS eps,PetscInt *bs)
1187 slepc 418
{
2331 jroman 419
  EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
1187 slepc 420
 
421
  PetscFunctionBegin;
422
  if (bs) *bs = ops->primme.maxBlockSize;
423
  PetscFunctionReturn(0);
424
}
425
EXTERN_C_END
426
 
427
#undef __FUNCT__  
428
#define __FUNCT__ "EPSPRIMMEGetBlockSize"
1201 slepc 429
/*@
2317 jroman 430
   EPSPRIMMEGetBlockSize - Get the maximum block size the code will try to use.
1187 slepc 431
 
2317 jroman 432
   Collective on EPS
433
 
1187 slepc 434
   Input Parameters:
1206 slepc 435
.  eps - the eigenproblem solver context
1187 slepc 436
 
437
   Output Parameters:  
1206 slepc 438
.  bs - returned block size
1187 slepc 439
 
440
   Level: advanced
441
.seealso: EPSPRIMMESetBlockSize()
442
@*/
1509 slepc 443
PetscErrorCode EPSPRIMMEGetBlockSize(EPS eps,PetscInt *bs)
1187 slepc 444
{
2221 jroman 445
  PetscErrorCode ierr;
1187 slepc 446
 
447
  PetscFunctionBegin;
2213 jroman 448
  PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
2221 jroman 449
  ierr = PetscTryMethod(eps,"EPSPRIMMEGetBlockSize_C",(EPS,PetscInt*),(eps,bs));CHKERRQ(ierr);
1187 slepc 450
  PetscFunctionReturn(0);
451
}
452
 
453
EXTERN_C_BEGIN
454
#undef __FUNCT__  
455
#define __FUNCT__ "EPSPRIMMESetMethod_PRIMME"
2326 jroman 456
PetscErrorCode EPSPRIMMESetMethod_PRIMME(EPS eps,EPSPRIMMEMethod method)
1187 slepc 457
{
2331 jroman 458
  EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
1187 slepc 459
 
460
  PetscFunctionBegin;
1446 slepc 461
  if (method == PETSC_DEFAULT) ops->method = DEFAULT_MIN_TIME;
1187 slepc 462
  else ops->method = (primme_preset_method)method;
463
  PetscFunctionReturn(0);
464
}
465
EXTERN_C_END
466
 
467
#undef __FUNCT__  
468
#define __FUNCT__ "EPSPRIMMESetMethod"
1201 slepc 469
/*@
1187 slepc 470
   EPSPRIMMESetMethod - Sets the method for the PRIMME library.
471
 
472
   Collective on EPS
473
 
474
   Input Parameters:
1206 slepc 475
+  eps - the eigenproblem solver context
1212 slepc 476
-  method - method that will be used by PRIMME. It must be one of:
1940 jroman 477
    EPS_PRIMME_DYNAMIC, EPS_PRIMME_DEFAULT_MIN_TIME(EPS_PRIMME_JDQMR_ETOL),
478
    EPS_PRIMME_DEFAULT_MIN_MATVECS(EPS_PRIMME_GD_OLSEN_PLUSK), EPS_PRIMME_ARNOLDI,
479
    EPS_PRIMME_GD, EPS_PRIMME_GD_PLUSK, EPS_PRIMME_GD_OLSEN_PLUSK,
480
    EPS_PRIMME_JD_OLSEN_PLUSK, EPS_PRIMME_RQI, EPS_PRIMME_JDQR, EPS_PRIMME_JDQMR,
481
    EPS_PRIMME_JDQMR_ETOL, EPS_PRIMME_SUBSPACE_ITERATION,
482
    EPS_PRIMME_LOBPCG_ORTHOBASIS, EPS_PRIMME_LOBPCG_ORTHOBASISW
1187 slepc 483
 
484
   Options Database Key:
1206 slepc 485
.  -eps_primme_set_method - Sets the method for the PRIMME library.
1187 slepc 486
 
1212 slepc 487
   Note:
1940 jroman 488
   If not set, the method defaults to EPS_PRIMME_DEFAULT_MIN_TIME.
1187 slepc 489
 
490
   Level: advanced
1364 slepc 491
 
492
.seealso: EPSPRIMMEGetMethod(), EPSPRIMMEMethod
1187 slepc 493
@*/
2326 jroman 494
PetscErrorCode EPSPRIMMESetMethod(EPS eps,EPSPRIMMEMethod method)
1187 slepc 495
{
2221 jroman 496
  PetscErrorCode ierr;
1187 slepc 497
 
498
  PetscFunctionBegin;
2213 jroman 499
  PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
2326 jroman 500
  PetscValidLogicalCollectiveEnum(eps,method,2);
2221 jroman 501
  ierr = PetscTryMethod(eps,"EPSPRIMMESetMethod_C",(EPS,EPSPRIMMEMethod),(eps,method));CHKERRQ(ierr);
1187 slepc 502
  PetscFunctionReturn(0);
503
}
504
 
505
EXTERN_C_BEGIN
506
#undef __FUNCT__  
507
#define __FUNCT__ "EPSPRIMMEGetMethod_PRIMME"
2326 jroman 508
PetscErrorCode EPSPRIMMEGetMethod_PRIMME(EPS eps,EPSPRIMMEMethod *method)
1187 slepc 509
{
2331 jroman 510
  EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
1187 slepc 511
 
512
  PetscFunctionBegin;
2317 jroman 513
  if (method) *method = (EPSPRIMMEMethod)ops->method;
1187 slepc 514
  PetscFunctionReturn(0);
515
}
516
EXTERN_C_END
517
 
518
#undef __FUNCT__  
519
#define __FUNCT__ "EPSPRIMMEGetMethod"
520
/*@C
521
    EPSPRIMMEGetMethod - Gets the method for the PRIMME library.
522
 
523
    Mon Collective on EPS
524
 
525
   Input Parameters:
1206 slepc 526
.  eps - the eigenproblem solver context
1187 slepc 527
 
528
   Output Parameters:
1212 slepc 529
.  method - method that will be used by PRIMME. It must be one of:
1940 jroman 530
    EPS_PRIMME_DYNAMIC, EPS_PRIMME_DEFAULT_MIN_TIME(EPS_PRIMME_JDQMR_ETOL),
531
    EPS_PRIMME_DEFAULT_MIN_MATVECS(EPS_PRIMME_GD_OLSEN_PLUSK), EPS_PRIMME_ARNOLDI,
532
    EPS_PRIMME_GD, EPS_PRIMME_GD_PLUSK, EPS_PRIMME_GD_OLSEN_PLUSK,
533
    EPS_PRIMME_JD_OLSEN_PLUSK, EPS_PRIMME_RQI, EPS_PRIMME_JDQR, EPS_PRIMME_JDQMR,
534
    EPS_PRIMME_JDQMR_ETOL, EPS_PRIMME_SUBSPACE_ITERATION,
535
    EPS_PRIMME_LOBPCG_ORTHOBASIS, EPS_PRIMME_LOBPCG_ORTHOBASISW
1187 slepc 536
 
537
    Level: advanced
1364 slepc 538
 
539
.seealso: EPSPRIMMESetMethod(), EPSPRIMMEMethod
1187 slepc 540
@*/
1190 slepc 541
PetscErrorCode EPSPRIMMEGetMethod(EPS eps, EPSPRIMMEMethod *method)
1187 slepc 542
{
2221 jroman 543
  PetscErrorCode ierr;
1187 slepc 544
 
545
  PetscFunctionBegin;
2213 jroman 546
  PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
2221 jroman 547
  ierr = PetscTryMethod(eps,"EPSPRIMMEGetMethod_C",(EPS,EPSPRIMMEMethod*),(eps,method));CHKERRQ(ierr);
1187 slepc 548
  PetscFunctionReturn(0);
549
}
550
 
551
EXTERN_C_BEGIN
552
#undef __FUNCT__  
553
#define __FUNCT__ "EPSCreate_PRIMME"
554
PetscErrorCode EPSCreate_PRIMME(EPS eps)
555
{
556
  PetscErrorCode ierr;
557
  EPS_PRIMME     *primme;
558
 
559
  PetscFunctionBegin;
2331 jroman 560
  ierr = STSetType(eps->OP,STPRECOND);CHKERRQ(ierr);
561
  ierr = STPrecondSetKSPHasMat(eps->OP,PETSC_TRUE);CHKERRQ(ierr);
2001 eromero 562
 
2329 jroman 563
  ierr = PetscNewLog(eps,EPS_PRIMME,&primme);CHKERRQ(ierr);
2331 jroman 564
  eps->data                      = (void*)primme;
1187 slepc 565
  eps->ops->setup                = EPSSetUp_PRIMME;
566
  eps->ops->setfromoptions       = EPSSetFromOptions_PRIMME;
567
  eps->ops->destroy              = EPSDestroy_PRIMME;
2348 jroman 568
  eps->ops->reset                = EPSReset_PRIMME;
1187 slepc 569
  eps->ops->view                 = EPSView_PRIMME;
570
  eps->ops->backtransform        = EPSBackTransform_Default;
571
  eps->ops->computevectors       = EPSComputeVectors_Default;
572
 
573
  primme_initialize(&primme->primme);
574
  primme->primme.matrixMatvec = multMatvec_PRIMME;
1239 slepc 575
  primme->primme.globalSumDouble = par_GlobalSumDouble;
1940 jroman 576
  primme->method = (primme_preset_method)EPS_PRIMME_DEFAULT_MIN_TIME;
1942 jroman 577
 
1187 slepc 578
  ierr = PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPRIMMESetBlockSize_C","EPSPRIMMESetBlockSize_PRIMME",EPSPRIMMESetBlockSize_PRIMME);CHKERRQ(ierr);
579
  ierr = PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPRIMMESetMethod_C","EPSPRIMMESetMethod_PRIMME",EPSPRIMMESetMethod_PRIMME);CHKERRQ(ierr);
580
  ierr = PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPRIMMEGetBlockSize_C","EPSPRIMMEGetBlockSize_PRIMME",EPSPRIMMEGetBlockSize_PRIMME);CHKERRQ(ierr);
581
  ierr = PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPRIMMEGetMethod_C","EPSPRIMMEGetMethod_PRIMME",EPSPRIMMEGetMethod_PRIMME);CHKERRQ(ierr);
582
  PetscFunctionReturn(0);
583
}
584
EXTERN_C_END