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
2044 antodo 1
/*                      
2
 
3
   Q-Arnoldi method for quadratic eigenproblems.
4
 
5
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
6
   SLEPc - Scalable Library for Eigenvalue Problem Computations
2116 eromero 7
   Copyright (c) 2002-2010, Universidad Politecnica de Valencia, Spain
2044 antodo 8
 
9
   This file is part of SLEPc.
10
 
11
   SLEPc is free software: you can redistribute it and/or modify it under  the
12
   terms of version 3 of the GNU Lesser General Public License as published by
13
   the Free Software Foundation.
14
 
15
   SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
16
   WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
17
   FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
18
   more details.
19
 
20
   You  should have received a copy of the GNU Lesser General  Public  License
21
   along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
22
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
23
*/
24
 
2284 jroman 25
#include <private/qepimpl.h>         /*I "slepcqep.h" I*/
26
#include <petscblaslapack.h>
2044 antodo 27
 
28
typedef struct {
29
  KSP ksp;
30
} QEP_QARNOLDI;
31
 
32
#undef __FUNCT__  
33
#define __FUNCT__ "QEPSetUp_QARNOLDI"
34
PetscErrorCode QEPSetUp_QARNOLDI(QEP qep)
35
{
36
  PetscErrorCode ierr;
37
  QEP_QARNOLDI   *ctx = (QEP_QARNOLDI *)qep->data;
38
 
39
  PetscFunctionBegin;
40
  if (qep->ncv) { /* ncv set */
2214 jroman 41
    if (qep->ncv<qep->nev) SETERRQ(((PetscObject)qep)->comm,1,"The value of ncv must be at least nev");
2044 antodo 42
  }
43
  else if (qep->mpd) { /* mpd set */
44
    qep->ncv = PetscMin(qep->n,qep->nev+qep->mpd);
45
  }
46
  else { /* neither set: defaults depend on nev being small or large */
47
    if (qep->nev<500) qep->ncv = PetscMin(qep->n,PetscMax(2*qep->nev,qep->nev+15));
48
    else { qep->mpd = 500; qep->ncv = PetscMin(qep->n,qep->nev+qep->mpd); }
49
  }
50
  if (!qep->mpd) qep->mpd = qep->ncv;
2214 jroman 51
  if (qep->ncv>qep->nev+qep->mpd) SETERRQ(((PetscObject)qep)->comm,1,"The value of ncv must not be larger than nev+mpd");
2044 antodo 52
  if (!qep->max_it) qep->max_it = PetscMax(100,2*qep->n/qep->ncv);
53
  if (!qep->which) qep->which = QEP_LARGEST_MAGNITUDE;
54
  if (qep->problem_type != QEP_GENERAL)
2214 jroman 55
    SETERRQ(((PetscObject)qep)->comm,1,"Wrong value of qep->problem_type");
2044 antodo 56
 
57
  ierr = PetscFree(qep->T);CHKERRQ(ierr);
58
  ierr = PetscMalloc(qep->ncv*qep->ncv*sizeof(PetscScalar),&qep->T);CHKERRQ(ierr);
59
  ierr = QEPDefaultGetWork(qep,4);CHKERRQ(ierr);
60
 
2051 jroman 61
  ierr = KSPSetOperators(ctx->ksp,qep->M,qep->M,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
2044 antodo 62
  ierr = KSPSetUp(ctx->ksp);CHKERRQ(ierr);
63
  PetscFunctionReturn(0);
64
}
65
 
66
#undef __FUNCT__  
67
#define __FUNCT__ "QEPQArnoldiCGS"
68
/*
2066 jroman 69
  Compute a step of Classical Gram-Schmidt orthogonalization
2044 antodo 70
*/
71
PetscErrorCode QEPQArnoldiCGS(QEP qep,PetscScalar *H,PetscBLASInt ldh,PetscScalar *h,PetscBLASInt j,Vec *V,Vec t,Vec v,Vec w,PetscReal *onorm,PetscReal *norm,PetscScalar *work)
72
{
73
  PetscErrorCode ierr;
74
  PetscBLASInt   ione = 1, j_1 = j+1;
2060 jroman 75
  PetscReal      x, y;
76
  PetscScalar    dot, one = 1.0, zero = 0.0;
2044 antodo 77
 
78
  PetscFunctionBegin;
79
  /* compute norm of v and w */
80
  if (onorm) {
81
    ierr = VecNorm(v,NORM_2,&x);CHKERRQ(ierr);
82
    ierr = VecNorm(w,NORM_2,&y);CHKERRQ(ierr);
2060 jroman 83
    *onorm = sqrt(x*x+y*y);
2044 antodo 84
  }
85
 
86
  /* orthogonalize: compute h */
87
  ierr = VecMDot(v,j_1,V,h);CHKERRQ(ierr);
88
  ierr = VecMDot(w,j_1,V,work);CHKERRQ(ierr);
89
  if (j>0)
2067 jroman 90
    BLASgemv_("C",&j_1,&j,&one,H,&ldh,work,&ione,&one,h,&ione);
2060 jroman 91
  ierr = VecDot(t,w,&dot);CHKERRQ(ierr);
92
  h[j] += dot;
2044 antodo 93
 
94
  /* orthogonalize: update v and w */
95
  ierr = SlepcVecMAXPBY(v,1.0,-1.0,j_1,h,V);CHKERRQ(ierr);
96
  if (j>0) {
97
    BLASgemv_("N",&j_1,&j,&one,H,&ldh,h,&ione,&zero,work,&ione);
98
    ierr = SlepcVecMAXPBY(w,1.0,-1.0,j_1,work,V);CHKERRQ(ierr);
99
  }
100
  ierr = VecAXPY(w,-h[j],t);CHKERRQ(ierr);
101
 
102
  /* compute norm of v and w */
103
  if (norm) {
104
    ierr = VecNorm(v,NORM_2,&x);CHKERRQ(ierr);
105
    ierr = VecNorm(w,NORM_2,&y);CHKERRQ(ierr);
2060 jroman 106
    *norm = sqrt(x*x+y*y);CHKERRQ(ierr);
2044 antodo 107
  }
108
  PetscFunctionReturn(0);
109
}
110
 
111
#undef __FUNCT__  
112
#define __FUNCT__ "QEPQArnoldi"
113
/*
2066 jroman 114
  Compute a run of Q-Arnoldi iterations
2044 antodo 115
*/
2216 jroman 116
PetscErrorCode QEPQArnoldi(QEP qep,PetscScalar *H,PetscInt ldh,Vec *V,PetscInt k,PetscInt *M,Vec v,Vec w,PetscReal *beta,PetscBool *breakdown,PetscScalar *work)
2044 antodo 117
{
118
  PetscErrorCode ierr;
119
  PetscInt       i,j,l,m = *M;
120
  QEP_QARNOLDI   *ctx = (QEP_QARNOLDI *)qep->data;
121
  Vec            t = qep->work[2], u = qep->work[3];
122
  IPOrthogonalizationRefinementType refinement;
123
  PetscReal      norm,onorm,eta;
124
  PetscScalar    *c = work + m;
125
 
126
  PetscFunctionBegin;
127
  ierr = IPGetOrthogonalization(qep->ip,PETSC_NULL,&refinement,&eta);CHKERRQ(ierr);
128
  ierr = VecCopy(v,qep->V[k]);CHKERRQ(ierr);
129
 
130
  for (j=k;j<m;j++) {
131
    /* apply operator */
132
    ierr = VecCopy(w,t);CHKERRQ(ierr);
2051 jroman 133
    ierr = MatMult(qep->K,v,u);CHKERRQ(ierr);
2050 antodo 134
    ierr = MatMult(qep->C,t,w);CHKERRQ(ierr);
135
    ierr = VecAXPY(u,qep->sfactor,w);CHKERRQ(ierr);
2044 antodo 136
    ierr = KSPSolve(ctx->ksp,u,w);CHKERRQ(ierr);
2051 jroman 137
    ierr = VecScale(w,-1.0/(qep->sfactor*qep->sfactor));CHKERRQ(ierr);
2044 antodo 138
    ierr = VecCopy(t,v);CHKERRQ(ierr);
139
 
140
    /* orthogonalize */
141
    switch (refinement) {
142
      case IP_ORTH_REFINE_NEVER:
143
        ierr = QEPQArnoldiCGS(qep,H,ldh,H+ldh*j,j,V,t,v,w,PETSC_NULL,&norm,work);CHKERRQ(ierr);
144
        *breakdown = PETSC_FALSE;
145
        break;
146
      case IP_ORTH_REFINE_ALWAYS:
147
        ierr = QEPQArnoldiCGS(qep,H,ldh,H+ldh*j,j,V,t,v,w,PETSC_NULL,PETSC_NULL,work);CHKERRQ(ierr);
148
        ierr = QEPQArnoldiCGS(qep,H,ldh,c,j,V,t,v,w,&onorm,&norm,work);CHKERRQ(ierr);
149
        for (i=0;i<j;i++) H[ldh*j+i] += c[i];
150
        if (norm < eta * onorm) *breakdown = PETSC_TRUE;
151
        else *breakdown = PETSC_FALSE;
152
        break;
153
      case IP_ORTH_REFINE_IFNEEDED:
154
        ierr = QEPQArnoldiCGS(qep,H,ldh,H+ldh*j,j,V,t,v,w,&onorm,&norm,work);CHKERRQ(ierr);
155
        /* ||q|| < eta ||h|| */
156
        l = 1;
157
        while (l<3 && norm < eta * onorm) {
158
          l++;
159
          onorm = norm;
160
          ierr = QEPQArnoldiCGS(qep,H,ldh,c,j,V,t,v,w,PETSC_NULL,&norm,work);CHKERRQ(ierr);
161
          for (i=0;i<j;i++) H[ldh*j+i] += c[i];
162
        }
163
        if (norm < eta * onorm) *breakdown = PETSC_TRUE;
164
        else *breakdown = PETSC_FALSE;
165
        break;
2214 jroman 166
      default: SETERRQ(((PetscObject)qep)->comm,1,"Wrong value of ip->orth_ref");
2044 antodo 167
    }
168
    ierr = VecScale(v,1.0/norm);CHKERRQ(ierr);
169
    ierr = VecScale(w,1.0/norm);CHKERRQ(ierr);
170
 
171
    if (j<m-1) {
172
      H[j+1+ldh*j] = norm;
173
      ierr = VecCopy(v,V[j+1]);CHKERRQ(ierr);
174
    }
175
  }
176
  *beta = norm;
177
  PetscFunctionReturn(0);
178
}
179
 
180
 
181
#undef __FUNCT__  
182
#define __FUNCT__ "QEPProjectedKSNonsym"
183
/*
184
   QEPProjectedKSNonsym - Solves the projected eigenproblem in the Krylov-Schur
185
   method (non-symmetric case).
186
 
187
   On input:
188
     l is the number of vectors kept in previous restart (0 means first restart)
189
     S is the projected matrix (leading dimension is lds)
190
 
191
   On output:
192
     S has (real) Schur form with diagonal blocks sorted appropriately
193
     Q contains the corresponding Schur vectors (order n, leading dimension n)
194
*/
195
PetscErrorCode QEPProjectedKSNonsym(QEP qep,PetscInt l,PetscScalar *S,PetscInt lds,PetscScalar *Q,PetscInt n)
196
{
197
  PetscErrorCode ierr;
198
  PetscInt       i;
199
 
200
  PetscFunctionBegin;
201
  if (l==0) {
202
    ierr = PetscMemzero(Q,n*n*sizeof(PetscScalar));CHKERRQ(ierr);
203
    for (i=0;i<n;i++)
204
      Q[i*(n+1)] = 1.0;
205
  } else {
206
    /* Reduce S to Hessenberg form, S <- Q S Q' */
207
    ierr = EPSDenseHessenberg(n,qep->nconv,S,lds,Q);CHKERRQ(ierr);
208
  }
209
  /* Reduce S to (quasi-)triangular form, S <- Q S Q' */
210
  ierr = EPSDenseSchur(n,qep->nconv,S,lds,Q,qep->eigr,qep->eigi);CHKERRQ(ierr);
211
  /* Sort the remaining columns of the Schur form */
212
  ierr = QEPSortDenseSchur(qep,n,qep->nconv,S,lds,Q,qep->eigr,qep->eigi);CHKERRQ(ierr);    
213
  PetscFunctionReturn(0);
214
}
215
 
216
#undef __FUNCT__  
217
#define __FUNCT__ "QEPSolve_QARNOLDI"
218
PetscErrorCode QEPSolve_QARNOLDI(QEP qep)
219
{
220
  PetscErrorCode ierr;
221
  PetscInt       i,j,k,l,lwork,nv;
222
  Vec            v=qep->work[0],w=qep->work[1];
2060 jroman 223
  PetscScalar    *S=qep->T,*Q,*work;
224
  PetscReal      beta,norm,x,y;
2216 jroman 225
  PetscBool      breakdown;
2044 antodo 226
 
227
  PetscFunctionBegin;
228
  ierr = PetscMemzero(S,qep->ncv*qep->ncv*sizeof(PetscScalar));CHKERRQ(ierr);
229
  ierr = PetscMalloc(qep->ncv*qep->ncv*sizeof(PetscScalar),&Q);CHKERRQ(ierr);
2066 jroman 230
  lwork = 7*qep->ncv;
2044 antodo 231
  ierr = PetscMalloc(lwork*sizeof(PetscScalar),&work);CHKERRQ(ierr);
232
 
233
  /* Get the starting Arnoldi vector */
234
  if (qep->nini>0) {
235
    ierr = VecCopy(qep->V[0],v);CHKERRQ(ierr);
236
  } else {
237
    ierr = SlepcVecSetRandom(v,qep->rand);CHKERRQ(ierr);
238
  }
239
  /* w is always a random vector */
240
  ierr = SlepcVecSetRandom(w,qep->rand);CHKERRQ(ierr);
241
  ierr = VecNorm(v,NORM_2,&x);CHKERRQ(ierr);
242
  ierr = VecNorm(w,NORM_2,&y);CHKERRQ(ierr);
2060 jroman 243
  norm = sqrt(x*x+y*y);CHKERRQ(ierr);
2066 jroman 244
  ierr = VecScale(v,1.0/norm);CHKERRQ(ierr);
245
  ierr = VecScale(w,1.0/norm);CHKERRQ(ierr);
2044 antodo 246
 
247
  /* Restart loop */
248
  l = 0;
2112 eromero 249
  while (qep->reason == QEP_CONVERGED_ITERATING) {
2044 antodo 250
    qep->its++;
251
 
252
    /* Compute an nv-step Arnoldi factorization */
253
    nv = PetscMin(qep->nconv+qep->mpd,qep->ncv);
254
    ierr = QEPQArnoldi(qep,S,qep->ncv,qep->V,qep->nconv+l,&nv,v,w,&beta,&breakdown,work);CHKERRQ(ierr);
255
 
256
    /* Solve projected problem */
257
    ierr = QEPProjectedKSNonsym(qep,l,S,qep->ncv,Q,nv);CHKERRQ(ierr);
258
 
2066 jroman 259
    /* Check convergence */
260
    ierr = QEPKrylovConvergence(qep,qep->nconv,nv-qep->nconv,S,qep->ncv,Q,nv,beta,&k,work);CHKERRQ(ierr);
2112 eromero 261
    if (qep->its >= qep->max_it) qep->reason = QEP_DIVERGED_ITS;
262
    if (k >= qep->nev) qep->reason = QEP_CONVERGED_TOL;
2044 antodo 263
 
264
    /* Update l */
2112 eromero 265
    if (qep->reason != QEP_CONVERGED_ITERATING || breakdown) l = 0;
2044 antodo 266
    else {
267
      l = (nv-k)/2;
268
#if !defined(PETSC_USE_COMPLEX)
269
      if (S[(k+l-1)*(qep->ncv+1)+1] != 0.0) {
270
        if (k+l<nv-1) l = l+1;
271
        else l = l-1;
272
      }
273
#endif
274
    }
275
 
2112 eromero 276
    if (qep->reason == QEP_CONVERGED_ITERATING) {
2044 antodo 277
      if (breakdown) {
278
        /* Stop if breakdown */
279
        PetscInfo2(qep,"Breakdown Quadratic Arnoldi method (it=%i norm=%g)\n",qep->its,beta);
2112 eromero 280
        qep->reason = QEP_DIVERGED_BREAKDOWN;
2044 antodo 281
      } else {
282
        /* Prepare the Rayleigh quotient for restart */
283
        for (i=k;i<k+l;i++) {
284
          S[i*qep->ncv+k+l] = Q[(i+1)*nv-1]*beta;
285
        }
286
      }
287
    }
288
    /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
289
    ierr = SlepcUpdateVectors(nv,qep->V,qep->nconv,k+l,Q,nv,PETSC_FALSE);CHKERRQ(ierr);
290
 
291
    qep->nconv = k;
2313 jroman 292
    ierr = QEPMonitor(qep,qep->its,qep->nconv,qep->eigr,qep->eigi,qep->errest,nv);CHKERRQ(ierr);
2044 antodo 293
  }
294
 
295
  for (j=0;j<qep->nconv;j++) {
2051 jroman 296
    qep->eigr[j] *= qep->sfactor;
297
    qep->eigi[j] *= qep->sfactor;
2044 antodo 298
  }
299
 
300
  /* Compute eigenvectors */
2049 antodo 301
  if (qep->nconv > 0) {
302
    ierr = QEPComputeVectors_Schur(qep);
303
  }
2044 antodo 304
 
305
  ierr = PetscFree(Q);CHKERRQ(ierr);
306
  ierr = PetscFree(work);CHKERRQ(ierr);
307
  PetscFunctionReturn(0);
308
}
309
 
310
#undef __FUNCT__  
311
#define __FUNCT__ "QEPSetFromOptions_QARNOLDI"
312
PetscErrorCode QEPSetFromOptions_QARNOLDI(QEP qep)
313
{
314
  PetscErrorCode ierr;
315
  QEP_QARNOLDI   *ctx = (QEP_QARNOLDI *)qep->data;
316
 
317
  PetscFunctionBegin;
318
  ierr = KSPSetFromOptions(ctx->ksp);CHKERRQ(ierr);
319
  PetscFunctionReturn(0);
320
}
321
 
322
#undef __FUNCT__  
323
#define __FUNCT__ "QEPView_QARNOLDI"
324
PetscErrorCode QEPView_QARNOLDI(QEP qep,PetscViewer viewer)
325
{
326
  PetscErrorCode ierr;
327
  QEP_QARNOLDI   *ctx = (QEP_QARNOLDI *)qep->data;
328
 
329
  PetscFunctionBegin;
330
  ierr = KSPView(ctx->ksp,viewer);CHKERRQ(ierr);
331
  PetscFunctionReturn(0);
332
}
333
 
334
#undef __FUNCT__  
335
#define __FUNCT__ "QEPDestroy_QARNOLDI"
336
PetscErrorCode QEPDestroy_QARNOLDI(QEP qep)
337
{
338
  PetscErrorCode ierr;
339
  QEP_QARNOLDI   *ctx = (QEP_QARNOLDI *)qep->data;
340
 
341
  PetscFunctionBegin;
2305 jroman 342
  ierr = KSPDestroy(&ctx->ksp);CHKERRQ(ierr);
2044 antodo 343
  ierr = QEPDestroy_Default(qep);CHKERRQ(ierr);
344
  PetscFunctionReturn(0);
345
}
346
 
347
EXTERN_C_BEGIN
348
#undef __FUNCT__  
349
#define __FUNCT__ "QEPCreate_QARNOLDI"
350
PetscErrorCode QEPCreate_QARNOLDI(QEP qep)
351
{
352
  PetscErrorCode ierr;
353
  QEP_QARNOLDI   *ctx;
354
 
355
  PetscFunctionBegin;
356
  ierr = PetscNew(QEP_QARNOLDI,&ctx);CHKERRQ(ierr);
357
  PetscLogObjectMemory(qep,sizeof(QEP_QARNOLDI));
358
  qep->data                      = ctx;
359
  qep->ops->solve                = QEPSolve_QARNOLDI;
360
  qep->ops->setup                = QEPSetUp_QARNOLDI;
361
  qep->ops->setfromoptions       = QEPSetFromOptions_QARNOLDI;
362
  qep->ops->destroy              = QEPDestroy_QARNOLDI;
363
  qep->ops->view                 = QEPView_QARNOLDI;
364
 
365
  ierr = KSPCreate(((PetscObject)qep)->comm,&ctx->ksp);CHKERRQ(ierr);
366
  ierr = KSPSetOptionsPrefix(ctx->ksp,((PetscObject)qep)->prefix);CHKERRQ(ierr);
367
  ierr = KSPAppendOptionsPrefix(ctx->ksp,"qep_");CHKERRQ(ierr);
368
  ierr = PetscObjectIncrementTabLevel((PetscObject)ctx->ksp,(PetscObject)qep,1);CHKERRQ(ierr);  
369
  PetscLogObjectParent(qep,ctx->ksp);
370
  PetscFunctionReturn(0);
371
}
372
EXTERN_C_END
373