Actual source code: neprefine.c

slepc-3.6.1 2015-09-03
Report Typos and Errors
  1: /*
  2:    Newton refinement for NEP, simple version.

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.

 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.

 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.

 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/>.
 21:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 22: */

 24: #include <slepc/private/nepimpl.h>
 25: #include <slepcblaslapack.h>

 27: #define NREF_MAXIT 10

 29: typedef struct {
 30:   PetscSubcomm  subc;
 31:   VecScatter    *scatter_id;
 32:   Mat           *A;
 33:   Vec           vg,v;
 34:   FN            *fn;
 35: } NEPSimpNRefctx;

 39: static PetscErrorCode NEPSimpleNRefSetUp(NEP nep,NEPSimpNRefctx **ctx_)
 40: {
 42:   PetscInt       i,si,j,n0,m0,nloc,*idx1,*idx2;
 43:   IS             is1,is2;
 44:   NEPSimpNRefctx *ctx;
 45:   Vec            v;

 48:   PetscMalloc1(1,ctx_);
 49:   ctx = *ctx_;
 50:   if (nep->npart==1) {
 51:     ctx->subc = NULL;
 52:     ctx->scatter_id = NULL;
 53:     ctx->A = nep->A;
 54:     ctx->fn = nep->f;
 55:   } else {
 56:     PetscMalloc2(nep->nt,&ctx->A,nep->npart,&ctx->scatter_id);

 58:     /* Split in subcomunicators */
 59:     PetscSubcommCreate(PetscObjectComm((PetscObject)nep),&ctx->subc);
 60:     PetscSubcommSetNumber(ctx->subc,nep->npart);
 61:     PetscSubcommSetType(ctx->subc,PETSC_SUBCOMM_CONTIGUOUS);
 62:     PetscLogObjectMemory((PetscObject)nep,sizeof(PetscSubcomm));

 64:     /* Duplicate matrices */
 65:     for (i=0;i<nep->nt;i++) {
 66:       MatCreateRedundantMatrix(nep->A[i],0,PetscSubcommChild(ctx->subc),MAT_INITIAL_MATRIX,&ctx->A[i]);
 67:     }
 68:     MatCreateVecs(ctx->A[0],&ctx->v,NULL);

 70:     /* Duplicate FNs */
 71:     PetscMalloc1(nep->nt,&ctx->fn);
 72:     for (i=0;i<nep->nt;i++) {
 73:       FNDuplicate(nep->f[i],PetscSubcommChild(ctx->subc),&ctx->fn[i]);
 74:     }

 76:     /* Create scatters for sending vectors to each subcommucator */
 77:     BVGetColumn(nep->V,0,&v);
 78:     VecGetOwnershipRange(v,&n0,&m0);
 79:     BVRestoreColumn(nep->V,0,&v);
 80:     VecGetLocalSize(ctx->v,&nloc);
 81:     PetscMalloc2(m0-n0,&idx1,m0-n0,&idx2);
 82:     VecCreateMPI(PetscObjectComm((PetscObject)nep),nloc,PETSC_DECIDE,&ctx->vg);
 83:     for (si=0;si<nep->npart;si++) {
 84:       j = 0;
 85:       for (i=n0;i<m0;i++) {
 86:         idx1[j]   = i;
 87:         idx2[j++] = i+nep->n*si;
 88:       }
 89:       ISCreateGeneral(PetscObjectComm((PetscObject)nep),(m0-n0),idx1,PETSC_COPY_VALUES,&is1);
 90:       ISCreateGeneral(PetscObjectComm((PetscObject)nep),(m0-n0),idx2,PETSC_COPY_VALUES,&is2);
 91:       BVGetColumn(nep->V,0,&v);
 92:       VecScatterCreate(v,is1,ctx->vg,is2,&ctx->scatter_id[si]);
 93:       BVRestoreColumn(nep->V,0,&v);
 94:       ISDestroy(&is1);
 95:       ISDestroy(&is2);
 96:     }
 97:     PetscFree2(idx1,idx2);
 98:   }
 99:   return(0);  
100: }

102: /*
103:   Gather Eigenpair idx from subcommunicator with color sc
104: */
107: PetscErrorCode NEPSimpleNRefGatherEigenpair(NEP nep,NEPSimpNRefctx *ctx,PetscInt sc,PetscInt idx)
108: {
110:   PetscMPIInt    nproc,p;
111:   MPI_Comm       comm=((PetscObject)nep)->comm;
112:   Vec            v;
113:   PetscScalar    *array;

116:   /* The eigenvalue information is in the last process of the 
117:      subcommunicator sc. p is its mapping in the general comm */
118:   MPI_Comm_size(comm,&nproc);
119:   p = (nproc/nep->npart)*(sc+1)+PetscMin(nproc%nep->npart,sc+1)-1;
120:   MPI_Bcast(&nep->eigr[idx],1,MPIU_SCALAR,p,comm);

122:   if (nep->npart>1) {
123:     /* Gather nep->V[idx] from the subcommuniator sc */
124:     BVGetColumn(nep->V,idx,&v);
125:     if (ctx->subc->color==sc) {
126:       VecGetArray(ctx->v,&array);
127:       VecPlaceArray(ctx->vg,array);
128:     }
129:     VecScatterBegin(ctx->scatter_id[sc],ctx->vg,v,INSERT_VALUES,SCATTER_REVERSE);
130:     VecScatterEnd(ctx->scatter_id[sc],ctx->vg,v,INSERT_VALUES,SCATTER_REVERSE);
131:     if (ctx->subc->color==sc) {
132:       VecResetArray(ctx->vg);
133:       VecRestoreArray(ctx->v,&array);
134:     }
135:     BVRestoreColumn(nep->V,idx,&v);
136:   }
137:   return(0);
138: }

142: PetscErrorCode NEPSimpleNRefScatterEigenvector(NEP nep,NEPSimpNRefctx *ctx,PetscInt sc,PetscInt idx)
143: {
145:   Vec            v;
146:   PetscScalar    *array;
147:   
149:   if (nep->npart>1) {
150:     BVGetColumn(nep->V,idx,&v);
151:     if (ctx->subc->color==sc) {
152:       VecGetArray(ctx->v,&array);
153:       VecPlaceArray(ctx->vg,array);
154:     }
155:     VecScatterBegin(ctx->scatter_id[sc],v,ctx->vg,INSERT_VALUES,SCATTER_FORWARD);
156:     VecScatterEnd(ctx->scatter_id[sc],v,ctx->vg,INSERT_VALUES,SCATTER_FORWARD);
157:     if (ctx->subc->color==sc) {
158:       VecResetArray(ctx->vg);
159:       VecRestoreArray(ctx->v,&array);
160:     }
161:     BVRestoreColumn(nep->V,idx,&v);
162:   }
163:   return(0);
164: }

168: static PetscErrorCode NEPSimpleNRefSetUpSystem(NEP nep,NEPSimpNRefctx *ctx,Mat *A,PetscInt idx,Mat *M,Mat *T,PetscBool ini,Vec *t,Vec v)
169: {
170:   PetscErrorCode    ierr;
171:   PetscInt          i,st,ml,m0,m1,mg;
172:   PetscInt          *dnz,*onz,ncols,*cols2,*nnz,nt=nep->nt;
173:   PetscScalar       zero=0.0,*coeffs;
174:   PetscMPIInt       rank,size;
175:   MPI_Comm          comm;
176:   const PetscInt    *cols;
177:   const PetscScalar *vals,*array;
178:   Vec               w=t[1],q=t[0];

181:   comm = PetscObjectComm((PetscObject)A[0]);
182:   PetscMalloc1(nt,&coeffs);
183:   MPI_Comm_rank(comm,&rank);
184:   MPI_Comm_size(comm,&size);
185:   if (ini) {
186:     MatDuplicate(A[0],MAT_COPY_VALUES,T);
187:   } else {
188:     MatCopy(A[0],*T,SUBSET_NONZERO_PATTERN);
189:   }
190:   for (i=0;i<nt;i++) {
191:     FNEvaluateFunction(ctx->fn[i],nep->eigr[idx],coeffs+i);
192:   }
193:   if (coeffs[0]!=1.0) {
194:     MatScale(*T,coeffs[0]);
195:   }
196:   for (i=1;i<nt;i++) {
197:     MatAXPY(*T,coeffs[i],A[i],(ini)?DIFFERENT_NONZERO_PATTERN:SUBSET_NONZERO_PATTERN);
198:   }
199:   MatGetSize(*T,&mg,NULL);
200:   MatGetOwnershipRange(*T,&m0,&m1);
201:   if (ini) {
202:     MatCreate(comm,M);
203:     MatGetLocalSize(*T,&ml,NULL);
204:     if (rank==size-1) ml++;
205:     MatSetSizes(*M,ml,ml,mg+1,mg+1);
206:     MatSetFromOptions(*M);
207:     MatSetUp(*M);
208:     /* Preallocate M */
209:     if (size>1) {
210:       MatPreallocateInitialize(comm,ml,ml,dnz,onz);
211:       for (i=m0;i<m1;i++) {
212:         MatGetRow(*T,i,&ncols,&cols,NULL);
213:         MatPreallocateSet(i,ncols,cols,dnz,onz);
214:         MatPreallocateSet(i,1,&mg,dnz,onz);
215:         MatRestoreRow(*T,i,&ncols,&cols,NULL);
216:       }
217:       if (rank==size-1) {
218:         PetscCalloc1(mg+1,&cols2);
219:         for (i=0;i<mg+1;i++) cols2[i]=i;
220:         MatPreallocateSet(m1,mg+1,cols2,dnz,onz);
221:         PetscFree(cols2);
222:       }
223:       MatMPIAIJSetPreallocation(*M,0,dnz,0,onz);
224:       MatPreallocateFinalize(dnz,onz);
225:     } else {
226:       PetscCalloc1(mg+1,&nnz);
227:       for (i=0;i<mg;i++) {
228:         MatGetRow(*T,i,&ncols,NULL,NULL);
229:         nnz[i] = ncols+1;
230:         MatRestoreRow(*T,i,&ncols,NULL,NULL);
231:       }
232:       nnz[mg] = mg+1;
233:       MatSeqAIJSetPreallocation(*M,0,nnz);
234:       PetscFree(nnz);
235:     }
236:   }
237:   for (i=0;i<nt;i++) {
238:     FNEvaluateDerivative(ctx->fn[i],nep->eigr[idx],coeffs+i);
239:   }
240:   st = 0;
241:   for (i=0;i<nt && PetscAbsScalar(coeffs[i])==0.0;i++) st++;
242:   MatMult(A[st],v,w);
243:   if (coeffs[st]!=1.0) {
244:     VecScale(w,coeffs[st]);
245:   }
246:   for (i=st+1;i<nt;i++) {
247:     MatMult(A[i],v,q);
248:     VecAXPY(w,coeffs[i],q);
249:   }
250:   /* Set values */
251:   PetscMalloc1(m1-m0,&cols2);
252:   for (i=0;i<m1-m0;i++) cols2[i]=m0+i;
253:   VecGetArrayRead(w,&array);
254:   for (i=m0;i<m1;i++) {
255:     MatGetRow(*T,i,&ncols,&cols,&vals);
256:     MatSetValues(*M,1,&i,ncols,cols,vals,INSERT_VALUES);
257:     MatRestoreRow(*T,i,&ncols,&cols,&vals);
258:     MatSetValues(*M,1,&i,1,&mg,array+i-m0,INSERT_VALUES);
259:   }
260:   VecRestoreArrayRead(w,&array);
261:   VecConjugate(v);
262:   VecGetArrayRead(v,&array);
263:   MatSetValues(*M,1,&mg,m1-m0,cols2,array,INSERT_VALUES);
264:   MatSetValues(*M,1,&mg,1,&mg,&zero,INSERT_VALUES);
265:   VecRestoreArrayRead(v,&array);
266:   VecConjugate(v);
267:   MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);
268:   MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);  
269:   PetscFree(cols2);
270:   PetscFree(coeffs);
271:   return(0);
272: }

276: PetscErrorCode NEPNewtonRefinementSimple(NEP nep,PetscInt *maxits,PetscReal *tol,PetscInt k)
277: {
278:   PetscErrorCode    ierr;
279:   PetscInt          i,n,its,idx=0,*idx_sc,*its_sc,color;
280:   PetscMPIInt       rank,size;
281:   KSP               ksp;
282:   Mat               M=NULL,T=NULL;
283:   MPI_Comm          comm;
284:   Vec               r,v,dv,rr=NULL,dvv=NULL,t[2];
285:   const PetscScalar *array;
286:   PetscScalar       *array2;
287:   PetscReal         norm,error;
288:   PetscBool         ini=PETSC_TRUE,sc_pend,solved=PETSC_FALSE;
289:   NEPSimpNRefctx    *ctx;

292:   PetscLogEventBegin(NEP_Refine,nep,0,0,0);
293:   NEPSimpleNRefSetUp(nep,&ctx);
294:   its = (maxits)?*maxits:NREF_MAXIT;
295:   comm = (nep->npart==1)?PetscObjectComm((PetscObject)nep):PetscSubcommChild(ctx->subc);
296:   KSPCreate(comm,&ksp);
297:   if (nep->npart==1) {
298:     BVGetColumn(nep->V,0,&v);
299:   } else v = ctx->v;
300:   VecDuplicate(v,&r);
301:   VecDuplicate(v,&dv);
302:   VecDuplicate(v,&t[0]);
303:   VecDuplicate(v,&t[1]);
304:   if (nep->npart==1) { BVRestoreColumn(nep->V,0,&v); }
305:   MPI_Comm_size(comm,&size);
306:   MPI_Comm_rank(comm,&rank);
307:   VecGetLocalSize(r,&n);
308:   PetscMalloc2(nep->npart,&idx_sc,nep->npart,&its_sc);
309:   for (i=0;i<nep->npart;i++) its_sc[i] = 0;
310:   color = (nep->npart==1)?0:ctx->subc->color;
311:    
312:   /* Loop performing iterative refinements */
313:   while (!solved) {
314:     for (i=0;i<nep->npart;i++) {
315:       sc_pend = PETSC_TRUE;
316:       if (its_sc[i]==0) {
317:         idx_sc[i] = idx++;
318:         if (idx_sc[i]>=k) {
319:           sc_pend = PETSC_FALSE;
320:         } else {
321:           NEPSimpleNRefScatterEigenvector(nep,ctx,i,idx_sc[i]);
322:         }
323:       }  else { /* Gather Eigenpair from subcommunicator i */
324:         NEPSimpleNRefGatherEigenpair(nep,ctx,i,idx_sc[i]);
325:       }
326:       while (sc_pend) {
327:         if (tol) {
328:           NEPComputeError(nep,idx_sc[i],NEP_ERROR_RELATIVE,&error);
329:         }
330:         if (error<=*tol || its_sc[i]>=its) {
331:           idx_sc[i] = idx++;
332:           its_sc[i] = 0;
333:           if (idx_sc[i]<k) { NEPSimpleNRefScatterEigenvector(nep,ctx,i,idx_sc[i]); }
334:         } else {
335:           sc_pend = PETSC_FALSE;
336:           its_sc[i]++;
337:         }
338:         if (idx_sc[i]>=k) sc_pend = PETSC_FALSE;
339:       }
340:     }
341:     solved = PETSC_TRUE;
342:     for (i=0;i<nep->npart&&solved;i++) solved = (idx_sc[i]<k)?PETSC_FALSE:PETSC_TRUE;
343:     if (idx_sc[color]<k) {
344: #if !defined(PETSC_USE_COMPLEX)
345:       if (nep->eigi[idx_sc[color]]!=0.0) SETERRQ(PetscObjectComm((PetscObject)nep),1,"Simple Refinement not implemented in real scalar for complex eigenvalues");
346: #endif
347:       if (nep->npart==1) {
348:         BVGetColumn(nep->V,idx_sc[color],&v);
349:       } else v = ctx->v; 
350:       NEPSimpleNRefSetUpSystem(nep,ctx,ctx->A,idx_sc[color],&M,&T,ini,t,v);
351:       KSPSetOperators(ksp,M,M);
352:       if (ini) {
353:         KSPSetFromOptions(ksp);
354:         MatCreateVecs(M,&dvv,NULL);
355:         VecDuplicate(dvv,&rr);
356:         ini = PETSC_FALSE;
357:       }
358:       MatMult(T,v,r);
359:       VecGetArrayRead(r,&array);
360:       if (rank==size-1) {
361:         VecGetArray(rr,&array2);
362:         PetscMemcpy(array2,array,n*sizeof(PetscScalar));
363:         array2[n] = 0.0;
364:         VecRestoreArray(rr,&array2);
365:       } else {
366:         VecPlaceArray(rr,array);
367:       }
368:       KSPSolve(ksp,rr,dvv);
369:       if (rank != size-1) {
370:         VecResetArray(rr);
371:       }
372:       VecRestoreArrayRead(r,&array);
373:       VecGetArrayRead(dvv,&array);
374:       VecPlaceArray(dv,array);
375:       VecAXPY(v,-1.0,dv);
376:       VecNorm(v,NORM_2,&norm);
377:       VecScale(v,1.0/norm);
378:       VecResetArray(dv);
379:       if (rank==size-1) nep->eigr[idx_sc[color]] -= array[n];
380:       VecRestoreArrayRead(dvv,&array);
381:       if (nep->npart==1) { BVRestoreColumn(nep->V,idx_sc[color],&v); } 
382:     }
383:   }
384:   KSPDestroy(&ksp);
385:   MatDestroy(&M);
386:   MatDestroy(&T);
387:   VecDestroy(&t[0]);
388:   VecDestroy(&t[1]);
389:   VecDestroy(&dv);
390:   VecDestroy(&dvv);
391:   VecDestroy(&r);
392:   VecDestroy(&rr);
393:   PetscFree2(idx_sc,its_sc);
394:   if (nep->npart>1) {
395:     VecDestroy(&ctx->vg);
396:     VecDestroy(&ctx->v);
397:     PetscSubcommDestroy(&ctx->subc);
398:     for (i=0;i<nep->nt;i++) {
399:       MatDestroy(&ctx->A[i]);
400:     }
401:     for (i=0;i<nep->npart;i++) {
402:       VecScatterDestroy(&ctx->scatter_id[i]);
403:     }
404:     PetscFree2(ctx->A,ctx->scatter_id);
405:     for (i=0;i<nep->nt;i++) { FNDestroy(&ctx->fn[i]); }
406:     PetscFree(ctx->fn);
407:   }
408:   PetscFree(ctx);
409:   PetscLogEventEnd(NEP_Refine,nep,0,0,0);
410:   return(0);
411: }