Actual source code: asa.c

  1: #define PETSCKSP_DLL

  3: /*  -------------------------------------------------------------------- 

  5:       Contributed by Arvid Bessen, Columbia University, June 2007
  6:      
  7:      This file implements a ASA preconditioner in PETSc as part of PC.

  9:      The adaptive smoothed aggregation algorithm is described in the paper
 10:      "Adaptive Smoothed Aggregation (ASA)", M. Brezina, R. Falgout, S. MacLachlan,
 11:      T. Manteuffel, S. McCormick, and J. Ruge, SIAM Journal on Scientific Computing,
 12:      SISC Volume 25 Issue 6, Pages 1896-1920.

 14:      For an example usage of this preconditioner, see, e.g.
 15:      $PETSC_DIR/src/ksp/ksp/examples/tutorials/ex38.c ex39.c
 16:      and other files in that directory.

 18:      This code is still somewhat experimental. A number of improvements would be
 19:      - keep vectors allocated on each level, instead of destroying them
 20:        (see mainly PCApplyVcycleOnLevel_ASA)
 21:      - in PCCreateTransferOp_ASA we get all of the submatrices at once, this could
 22:        be optimized by differentiating between local and global matrices
 23:      - the code does not handle it gracefully if there is just one level
 24:      - if relaxation is sufficient, exit of PCInitializationStage_ASA is not
 25:        completely clean
 26:      - default values could be more reasonable, especially for parallel solves,
 27:        where we need a parallel LU or similar
 28:      - the richardson scaling parameter is somewhat special, should be treated in a
 29:        good default way
 30:      - a number of parameters for smoother (sor_omega, etc.) that we store explicitly
 31:        could be kept in the respective smoothers themselves
 32:      - some parameters have to be set via command line options, there are no direct
 33:        function calls available
 34:      - numerous other stuff

 36:      Example runs in parallel would be with parameters like
 37:      mpiexec ./program -pc_asa_coarse_mat_type aijmumps -pc_asa_direct_solver 200
 38:      -pc_asa_max_cand_vecs 4 -pc_asa_mu_initial 50 -pc_asa_richardson_scale 1.0
 39:      -pc_asa_rq_improve 0.9 -asa_smoother_pc_type asm -asa_smoother_sub_pc_type sor

 41:     -------------------------------------------------------------------- */

 43: /*
 44:   This defines the data structures for the smoothed aggregation procedure
 45: */
 46:  #include src/ksp/pc/impls/asa/asa.h

 48: /*
 49:   We need the QR algorithm from LAPACK
 50: */
 51:  #include petscblaslapack.h

 53: /* -------------------------------------------------------------------------- */

 55: /* Event logging */

 57: PetscEvent PC_InitializationStage_ASA, PC_GeneralSetupStage_ASA;
 58: PetscEvent PC_CreateTransferOp_ASA, PC_CreateVcycle_ASA;
 59: PetscTruth asa_events_registered = PETSC_FALSE;


 64: /*@C
 65:     PCASASetDM - Sets the coarse grid information for the grids

 67:     Collective on PC

 69:     Input Parameter:
 70: +   pc - the context
 71: -   dm - the DA or ADDA or VecPack object

 73:     Level: advanced

 75: @*/
 76: PetscErrorCode  PCASASetDM(PC pc,DM dm)
 77: {
 78:   PetscErrorCode ierr,(*f)(PC,DM);

 82:   PetscObjectQueryFunction((PetscObject)pc,"PCASASetDM_C",(void (**)(void))&f);
 83:   if (f) {
 84:     (*f)(pc,dm);
 85:   }
 86:   return(0);
 87: }

 91: PetscErrorCode  PCASASetDM_ASA(PC pc, DM dm)
 92: {
 94:   PC_ASA         *asa = (PC_ASA *) pc->data;

 97:   PetscObjectReference((PetscObject)dm);
 98:   asa->dm = dm;
 99:   return(0);
100: }

104: /*@C
105:     PCASASetTolerances - Sets the convergence thresholds for ASA algorithm

107:     Collective on PC

109:     Input Parameter:
110: +   pc - the context
111: .   rtol - the relative convergence tolerance
112:     (relative decrease in the residual norm)
113: .   abstol - the absolute convergence tolerance 
114:     (absolute size of the residual norm)
115: .   dtol - the divergence tolerance
116:     (amount residual can increase before KSPDefaultConverged()
117:     concludes that the method is diverging)
118: -   maxits - maximum number of iterations to use

120:     Notes:
121:     Use PETSC_DEFAULT to retain the default value of any of the tolerances.

123:     Level: advanced
124: @*/
125: PetscErrorCode  PCASASetTolerances(PC pc, PetscReal rtol, PetscReal abstol,PetscReal dtol, PetscInt maxits)
126: {
127:   PetscErrorCode ierr,(*f)(PC,PetscReal,PetscReal,PetscReal,PetscInt);

131:   PetscObjectQueryFunction((PetscObject)pc,"PCASASetTolerances_C",(void (**)(void))&f);
132:   if (f) {
133:     (*f)(pc,rtol,abstol,dtol,maxits);
134:   }
135:   return(0);
136: }

140: PetscErrorCode  PCASASetTolerances_ASA(PC pc, PetscReal rtol, PetscReal abstol,PetscReal dtol, PetscInt maxits)
141: {
142:   PC_ASA         *asa = (PC_ASA *) pc->data;

146:   if (rtol != PETSC_DEFAULT)   asa->rtol   = rtol;
147:   if (abstol != PETSC_DEFAULT)   asa->abstol   = abstol;
148:   if (dtol != PETSC_DEFAULT)   asa->divtol = dtol;
149:   if (maxits != PETSC_DEFAULT) asa->max_it = maxits;
150:   return(0);
151: }

155: /*
156:    PCCreateLevel_ASA - Creates one level for the ASA algorithm

158:    Input Parameters:
159: +  level - current level
160: .  comm - MPI communicator object
161: .  next - pointer to next level
162: .  prev - pointer to previous level
163: .  ksptype - the KSP type for the smoothers on this level
164: -  pctype - the PC type for the smoothers on this level

166:    Output Parameters:
167: .  new_asa_lev - the newly created level

169: .keywords: ASA, create, levels, multigrid
170: */
171: PetscErrorCode  PCCreateLevel_ASA(PC_ASA_level **new_asa_lev, int level,MPI_Comm comm, PC_ASA_level *prev,
172:                                                     PC_ASA_level *next,KSPType ksptype, PCType pctype)
173: {
175:   PC_ASA_level   *asa_lev;
176: 
178:   PetscMalloc(sizeof(PC_ASA_level), &asa_lev);

180:   asa_lev->level = level;
181:   asa_lev->size  = 0;

183:   asa_lev->A = 0;
184:   asa_lev->B = 0;
185:   asa_lev->x = 0;
186:   asa_lev->b = 0;
187:   asa_lev->r = 0;
188: 
189:   asa_lev->dm           = 0;
190:   asa_lev->aggnum       = 0;
191:   asa_lev->agg          = 0;
192:   asa_lev->loc_agg_dofs = 0;
193:   asa_lev->agg_corr     = 0;
194:   asa_lev->bridge_corr  = 0;
195: 
196:   asa_lev->P = 0;
197:   asa_lev->Pt = 0;
198:   asa_lev->smP = 0;
199:   asa_lev->smPt = 0;

201:   asa_lev->comm = comm;

203:   asa_lev->smoothd = 0;
204:   asa_lev->smoothu = 0;

206:   asa_lev->prev = prev;
207:   asa_lev->next = next;
208: 
209:   *new_asa_lev = asa_lev;
210:   return(0);
211: }

215: PetscErrorCode SafeMatDestroy(Mat *m)
216: {
217:   PetscErrorCode 0;

220:   if (m && *m) {MatDestroy(*m); *m=0;}
221:   PetscFunctionReturn(ierr);
222: }

226: PetscErrorCode SafeVecDestroy(Vec *v)
227: {
228:   PetscErrorCode 0;

231:   if (v && *v) {VecDestroy(*v); *v=0;}
232:   PetscFunctionReturn(ierr);
233: }

237: PetscErrorCode PrintResNorm(Mat A, Vec x, Vec b, Vec r)
238: {
240:   PetscTruth     destroyr = PETSC_FALSE;
241:   PetscReal      resnorm;
242:   MPI_Comm       Acomm;

245:   if (!r) {
246:     MatGetVecs(A, PETSC_NULL, &r);
247:     destroyr = PETSC_TRUE;
248:   }
249:   MatMult(A, x, r);
250:   VecAYPX(r, -1.0, b);
251:   VecNorm(r, NORM_2, &resnorm);
252:   PetscObjectGetComm((PetscObject) A, &Acomm);
253:   PetscPrintf(Acomm, "Residual norm is %f.\n", resnorm);

255:   if (destroyr) {
256:     VecDestroy(r);
257:   }
258: 
259:   return(0);
260: }

264: PetscErrorCode PrintEnergyNormOfDiff(Mat A, Vec x, Vec y)
265: {
267:   Vec            vecdiff, Avecdiff;
268:   PetscScalar    dotprod;
269:   PetscReal      dotabs;
270:   MPI_Comm       Acomm;
272:   VecDuplicate(x, &vecdiff);
273:   VecWAXPY(vecdiff, -1.0, x, y);
274:   MatGetVecs(A, PETSC_NULL, &Avecdiff);
275:   MatMult(A, vecdiff, Avecdiff);
276:   VecDot(vecdiff, Avecdiff, &dotprod);
277:   dotabs = PetscAbsScalar(dotprod);
278:   PetscObjectGetComm((PetscObject) A, &Acomm);
279:   PetscPrintf(Acomm, "Energy norm %f.\n", dotabs);
280:   VecDestroy(vecdiff);
281:   VecDestroy(Avecdiff);
282:   return(0);
283: }

285: /* -------------------------------------------------------------------------- */
286: /*
287:    PCDestroyLevel_ASA - Destroys one level of the ASA preconditioner

289:    Input Parameter:
290: .  asa_lev - pointer to level that should be destroyed

292: */
295: PetscErrorCode PCDestroyLevel_ASA(PC_ASA_level *asa_lev)
296: {

300:   SafeMatDestroy(&(asa_lev->A));
301:   SafeMatDestroy(&(asa_lev->B));
302:   SafeVecDestroy(&(asa_lev->x));
303:   SafeVecDestroy(&(asa_lev->b));
304:   SafeVecDestroy(&(asa_lev->r));

306:   if (asa_lev->dm) {DMDestroy(asa_lev->dm);}

308:   SafeMatDestroy(&(asa_lev->agg));
309:   PetscFree(asa_lev->loc_agg_dofs);
310:   SafeMatDestroy(&(asa_lev->agg_corr));
311:   SafeMatDestroy(&(asa_lev->bridge_corr));

313:   SafeMatDestroy(&(asa_lev->P));
314:   SafeMatDestroy(&(asa_lev->Pt));
315:   SafeMatDestroy(&(asa_lev->smP));
316:   SafeMatDestroy(&(asa_lev->smPt));

318:   if (asa_lev->smoothd != asa_lev->smoothu) {
319:     if (asa_lev->smoothd) {KSPDestroy(asa_lev->smoothd);}
320:   }
321:   if (asa_lev->smoothu) {KSPDestroy(asa_lev->smoothu);}

323:   PetscFree(asa_lev);
324:   return(0);
325: }

327: /* -------------------------------------------------------------------------- */
328: /*
329:    PCComputeSpectralRadius_ASA - Computes the spectral radius of asa_lev->A
330:    and stores it it asa_lev->spec_rad

332:    Input Parameters:
333: .  asa_lev - the level we are treating

335:    Compute spectral radius with  sqrt(||A||_1 ||A||_inf) >= ||A||_2 >= rho(A) 

337: */
340: PetscErrorCode PCComputeSpectralRadius_ASA(PC_ASA_level *asa_lev)
341: {
343:   PetscReal      norm_1, norm_inf;

346:   MatNorm(asa_lev->A, NORM_1, &norm_1);
347:   MatNorm(asa_lev->A, NORM_INFINITY, &norm_inf);
348:   asa_lev->spec_rad = sqrt(norm_1*norm_inf);
349:   return(0);
350: }

354: PetscErrorCode PCSetRichardsonScale_ASA(KSP ksp, PetscReal spec_rad, PetscReal richardson_scale) {
356:   PC             pc;
357:   PetscTruth     flg;
358:   PetscReal      spec_rad_inv;

361:   KSPSetInitialGuessNonzero(ksp, PETSC_TRUE);
362:   if (richardson_scale != PETSC_DECIDE) {
363:     KSPRichardsonSetScale(ksp, richardson_scale);
364:   } else {
365:     KSPGetPC(ksp, &pc);
366:     PetscTypeCompare((PetscObject)(pc), PCNONE, &flg);
367:     if (flg) {
368:       /* WORK: this is just an educated guess. Any number between 0 and 2/rho(A)
369:          should do. asa_lev->spec_rad has to be an upper bound on rho(A). */
370:       spec_rad_inv = 1.0/spec_rad;
371:       KSPRichardsonSetScale(ksp, spec_rad_inv);
372:     } else {
373:       SETERRQ(PETSC_ERR_SUP, "Unknown PC type for smoother. Please specify scaling factor with -pc_asa_richardson_scale\n");
374:     }
375:   }
376:   return(0);
377: }

381: PetscErrorCode PCSetSORomega_ASA(PC pc, PetscReal sor_omega)
382: {

386:   PCSORSetSymmetric(pc, SOR_SYMMETRIC_SWEEP);
387:   if (sor_omega != PETSC_DECIDE) {
388:     PCSORSetOmega(pc, sor_omega);
389:   }
390:   return(0);
391: }


394: /* -------------------------------------------------------------------------- */
395: /*
396:    PCSetupSmoothersOnLevel_ASA - Creates the smoothers of the level.
397:    We assume that asa_lev->A and asa_lev->spec_rad are correctly computed

399:    Input Parameters:
400: +  asa - the data structure for the ASA preconditioner
401: .  asa_lev - the level we are treating
402: -  maxits - maximum number of iterations to use
403: */
406: PetscErrorCode PCSetupSmoothersOnLevel_ASA(PC_ASA *asa, PC_ASA_level *asa_lev, PetscInt maxits)
407: {
408:   PetscErrorCode    ierr;
409:   PetscTruth        flg;
410:   PC                pc;

413:   /* destroy old smoothers */
414:   if (asa_lev->smoothu && asa_lev->smoothu != asa_lev->smoothd) {
415:     KSPDestroy(asa_lev->smoothu);
416:   }
417:   asa_lev->smoothu = 0;
418:   if (asa_lev->smoothd) {
419:     KSPDestroy(asa_lev->smoothd);
420:   }
421:   asa_lev->smoothd = 0;
422:   /* create smoothers */
423:   KSPCreate(asa_lev->comm,&asa_lev->smoothd);
424:   KSPSetType(asa_lev->smoothd, asa->ksptype_smooth);
425:   KSPGetPC(asa_lev->smoothd,&pc);
426:   PCSetType(pc,asa->pctype_smooth);

428:   /* set up problems for smoothers */
429:   KSPSetOperators(asa_lev->smoothd, asa_lev->A, asa_lev->A, DIFFERENT_NONZERO_PATTERN);
430:   KSPSetTolerances(asa_lev->smoothd, asa->smoother_rtol, asa->smoother_abstol, asa->smoother_dtol, maxits);
431:   PetscTypeCompare((PetscObject)(asa_lev->smoothd), KSPRICHARDSON, &flg);
432:   if (flg) {
433:     /* special parameters for certain smoothers */
434:     KSPSetInitialGuessNonzero(asa_lev->smoothd, PETSC_TRUE);
435:     KSPGetPC(asa_lev->smoothd, &pc);
436:     PetscTypeCompare((PetscObject)pc, PCSOR, &flg);
437:     if (flg) {
438:       PCSetSORomega_ASA(pc, asa->sor_omega);
439:     } else {
440:       /* just set asa->richardson_scale to get some very basic smoother */
441:       PCSetRichardsonScale_ASA(asa_lev->smoothd, asa_lev->spec_rad, asa->richardson_scale);
442:     }
443:     /* this would be the place to add support for other preconditioners */
444:   }
445:   KSPSetOptionsPrefix(asa_lev->smoothd, "asa_smoother_");
446:   KSPSetFromOptions(asa_lev->smoothd);
447:   /* set smoothu equal to smoothd, this could change later */
448:   asa_lev->smoothu = asa_lev->smoothd;
449:   return(0);
450: }

452: /* -------------------------------------------------------------------------- */
453: /*
454:    PCSetupDirectSolversOnLevel_ASA - Creates the direct solvers on the coarsest level.
455:    We assume that asa_lev->A and asa_lev->spec_rad are correctly computed

457:    Input Parameters:
458: +  asa - the data structure for the ASA preconditioner
459: .  asa_lev - the level we are treating
460: -  maxits - maximum number of iterations to use
461: */
464: PetscErrorCode PCSetupDirectSolversOnLevel_ASA(PC_ASA *asa, PC_ASA_level *asa_lev, PetscInt maxits)
465: {
466:   PetscErrorCode    ierr;
467:   PetscTruth        flg;
468:   PetscInt          comm_size;
469:   PC                pc;

472:   if (asa_lev->smoothu && asa_lev->smoothu != asa_lev->smoothd) {
473:     KSPDestroy(asa_lev->smoothu);
474:   }
475:   asa_lev->smoothu = 0;
476:   if (asa_lev->smoothd) {
477:     KSPDestroy(asa_lev->smoothd);
478:     asa_lev->smoothd = 0;
479:   }
480:   PetscStrcmp(asa->ksptype_direct, KSPPREONLY, &flg);
481:   if (flg) {
482:     PetscStrcmp(asa->pctype_direct, PCLU, &flg);
483:     if (flg) {
484:       MPI_Comm_size(asa_lev->comm, &comm_size);
485:       if (comm_size > 1) {
486:         /* the LU PC will call MatSolve, we may have to set the correct type for the matrix
487:            to have support for this in parallel */
488:         MatConvert(asa_lev->A, asa->coarse_mat_type, MAT_REUSE_MATRIX, &(asa_lev->A));
489:       }
490:     }
491:   }
492:   /* create new solvers */
493:   KSPCreate(asa_lev->comm,&asa_lev->smoothd);
494:   KSPSetType(asa_lev->smoothd, asa->ksptype_direct);
495:   KSPGetPC(asa_lev->smoothd,&pc);
496:   PCSetType(pc,asa->pctype_direct);
497:   /* set up problems for direct solvers */
498:   KSPSetOperators(asa_lev->smoothd, asa_lev->A, asa_lev->A, DIFFERENT_NONZERO_PATTERN);
499:   KSPSetTolerances(asa_lev->smoothd, asa->direct_rtol, asa->direct_abstol, asa->direct_dtol, maxits);
500:   /* user can set any option by using -pc_asa_direct_xxx */
501:   KSPSetOptionsPrefix(asa_lev->smoothd, "asa_coarse_");
502:   KSPSetFromOptions(asa_lev->smoothd);
503:   /* set smoothu equal to 0, not used */
504:   asa_lev->smoothu = 0;
505:   return(0);
506: }

508: /* -------------------------------------------------------------------------- */
509: /*
510:    PCCreateAggregates_ASA - Creates the aggregates

512:    Input Parameters:
513: .  asa_lev - the level for which we should create the projection matrix

515: */
518: PetscErrorCode PCCreateAggregates_ASA(PC_ASA_level *asa_lev)
519: {
520:   PetscInt          m,n, m_loc,n_loc;
521:   PetscInt          m_loc_s, m_loc_e;
522:   const PetscScalar one = 1.0;
523:   PetscErrorCode    ierr;

526:   /* Create nodal aggregates A_i^l */
527:   /* we use the DM grid information for that */
528:   if (asa_lev->dm) {
529:     /* coarsen DM and get the restriction matrix */
530:     DMCoarsen(asa_lev->dm, PETSC_NULL, &(asa_lev->next->dm));
531:     DMGetAggregates(asa_lev->next->dm, asa_lev->dm, &(asa_lev->agg));
532:     MatGetSize(asa_lev->agg, &m, &n);
533:     MatGetLocalSize(asa_lev->agg, &m_loc, &n_loc);
534:     if (n!=asa_lev->size) SETERRQ(PETSC_ERR_ARG_SIZ,"DM interpolation matrix has incorrect size!\n");
535:     asa_lev->next->size = m;
536:     asa_lev->aggnum     = m;
537:     /* create the correlators, right now just identity matrices */
538:     MatCreateMPIAIJ(asa_lev->comm, n_loc, n_loc, n, n, 1, PETSC_NULL, 1, PETSC_NULL,&(asa_lev->agg_corr));
539:     MatGetOwnershipRange(asa_lev->agg_corr, &m_loc_s, &m_loc_e);
540:     for (m=m_loc_s; m<m_loc_e; m++) {
541:       MatSetValues(asa_lev->agg_corr, 1, &m, 1, &m, &one, INSERT_VALUES);
542:     }
543:     MatAssemblyBegin(asa_lev->agg_corr, MAT_FINAL_ASSEMBLY);
544:     MatAssemblyEnd(asa_lev->agg_corr, MAT_FINAL_ASSEMBLY);
545: /*     MatShift(asa_lev->agg_corr, 1.0); */
546:   } else {
547:     /* somehow define the aggregates without knowing the geometry */
548:     /* future WORK */
549:     SETERRQ(PETSC_ERR_SUP, "Currently pure algebraic coarsening is not supported!");
550:   }
551:   return(0);
552: }

554: /* -------------------------------------------------------------------------- */
555: /*
556:    PCCreateTransferOp_ASA - Creates the transfer operator P_{l+1}^l for current level

558:    Input Parameters:
559: +  asa_lev - the level for which should create the transfer operator
560: -  construct_bridge - true, if we should construct a bridge operator, false for normal prolongator

562:    If we add a second, third, ... candidate vector (i.e. more than one column in B), we
563:    have to relate the additional dimensions to the original aggregates. This is done through
564:    the "aggregate correlators" agg_corr and bridge_corr.
565:    The aggregate that is used in the construction is then given by
566:    asa_lev->agg * asa_lev->agg_corr
567:    for the regular prolongator construction and
568:    asa_lev->agg * asa_lev->bridge_corr
569:    for the bridging prolongator constructions.
570: */
573: PetscErrorCode PCCreateTransferOp_ASA(PC_ASA_level *asa_lev, PetscTruth construct_bridge)
574: {

577:   const PetscReal Ca = 1e-3;
578:   PetscReal      cutoff;
579:   PetscInt       nodes_on_lev;

581:   Mat            logical_agg;
582:   PetscInt       mat_agg_loc_start, mat_agg_loc_end, mat_agg_loc_size;
583:   PetscInt       a;
584:   const PetscInt *agg = 0;
585:   PetscInt       **agg_arr = 0;

587:   IS             *idxm_is_B_arr = 0;
588:   PetscInt       *idxn_B = 0;
589:   IS             idxn_is_B, *idxn_is_B_arr = 0;

591:   Mat            *b_submat_arr = 0;

593:   PetscScalar    *b_submat = 0, *b_submat_tp = 0;
594:   PetscInt       *idxm = 0, *idxn = 0;
595:   PetscInt       cand_vecs_num;
596:   PetscInt       *cand_vec_length = 0;
597:   PetscInt       max_cand_vec_length = 0;
598:   PetscScalar    **b_orth_arr = 0;

600:   PetscInt       i,j;

602:   PetscScalar    *tau = 0, *work = 0;
603:   PetscBLASInt   info;

605:   PetscInt       max_cand_vecs_to_add;
606:   PetscInt       *new_loc_agg_dofs = 0;

608:   PetscInt       total_loc_cols = 0;
609:   PetscReal      norm;

611:   PetscInt       a_loc_m, a_loc_n;
612:   PetscInt       mat_loc_col_start, mat_loc_col_end, mat_loc_col_size;
613:   PetscInt       loc_agg_dofs_sum;
614:   PetscInt       row, col;
615:   PetscScalar    val;
616:   PetscInt       comm_size, comm_rank;
617:   PetscInt       *loc_cols = 0;


622:   MatGetSize(asa_lev->B, &nodes_on_lev, PETSC_NULL);

624:   /* If we add another candidate vector, we want to be able to judge, how much the new candidate
625:      improves our current projection operators and whether it is worth adding it.
626:      This is the precomputation necessary for implementing Notes (4.1) to (4.7).
627:      We require that all candidate vectors x stored in B are normalized such that
628:      <A x, x> = 1 and we thus do not have to compute this.
629:      For each aggregate A we can now test condition (4.5) and (4.6) by computing
630:      || quantity to check ||_{A}^2 <= cutoff * card(A)/N_l */
631:   cutoff = Ca/(asa_lev->spec_rad);

633:   /* compute logical aggregates by using the correlators */
634:   if (construct_bridge) {
635:     /* construct bridging operator */
636:     MatMatMult(asa_lev->agg, asa_lev->bridge_corr, MAT_INITIAL_MATRIX, 1.0, &logical_agg);
637:   } else {
638:     /* construct "regular" prolongator */
639:     MatMatMult(asa_lev->agg, asa_lev->agg_corr, MAT_INITIAL_MATRIX, 1.0, &logical_agg);
640:   }

642:   /* destroy correlator matrices for next level, these will be rebuilt in this routine */
643:   if (asa_lev->next) {
644:     SafeMatDestroy(&(asa_lev->next->agg_corr));
645:     SafeMatDestroy(&(asa_lev->next->bridge_corr));
646:   }

648:   /* find out the correct local row indices */
649:   MatGetOwnershipRange(logical_agg, &mat_agg_loc_start, &mat_agg_loc_end);
650:   mat_agg_loc_size = mat_agg_loc_end-mat_agg_loc_start;
651: 
652:   cand_vecs_num = asa_lev->cand_vecs;

654:   /* construct column indices idxn_B for reading from B */
655:   PetscMalloc(sizeof(PetscInt)*(cand_vecs_num), &idxn_B);
656:   for (i=0; i<cand_vecs_num; i++) {
657:     idxn_B[i] = i;
658:   }
659:   ISCreateGeneral(asa_lev->comm, asa_lev->cand_vecs, idxn_B, &idxn_is_B);
660:   PetscFree(idxn_B);
661:   PetscMalloc(sizeof(IS)*mat_agg_loc_size, &idxn_is_B_arr);
662:   for (a=0; a<mat_agg_loc_size; a++) {
663:     idxn_is_B_arr[a] = idxn_is_B;
664:   }
665:   /* allocate storage for row indices idxm_B */
666:   PetscMalloc(sizeof(IS)*mat_agg_loc_size, &idxm_is_B_arr);

668:   /* Storage for the orthogonalized  submatrices of B and their sizes */
669:   PetscMalloc(sizeof(PetscInt)*mat_agg_loc_size, &cand_vec_length);
670:   PetscMalloc(sizeof(PetscScalar*)*mat_agg_loc_size, &b_orth_arr);
671:   /* Storage for the information about each aggregate */
672:   PetscMalloc(sizeof(PetscInt*)*mat_agg_loc_size, &agg_arr);
673:   /* Storage for the number of candidate vectors that are orthonormal and used in each submatrix */
674:   PetscMalloc(sizeof(PetscInt)*mat_agg_loc_size, &new_loc_agg_dofs);

676:   /* loop over local aggregates */
677:   for (a=0; a<mat_agg_loc_size; a++) {
678:        /* get info about current aggregate, this gives the rows we have to get from B */
679:        MatGetRow(logical_agg, a+mat_agg_loc_start, &cand_vec_length[a], &agg, 0);
680:        /* copy aggregate information */
681:        PetscMalloc(sizeof(PetscInt)*cand_vec_length[a], &(agg_arr[a]));
682:        PetscMemcpy(agg_arr[a], agg, sizeof(PetscInt)*cand_vec_length[a]);
683:        /* restore row */
684:        MatRestoreRow(logical_agg, a+mat_agg_loc_start, &cand_vec_length[a], &agg, 0);
685: 
686:        /* create index sets */
687:        ISCreateGeneral(PETSC_COMM_SELF, cand_vec_length[a], agg_arr[a], &(idxm_is_B_arr[a]));
688:        /* maximum candidate vector length */
689:        if (cand_vec_length[a] > max_cand_vec_length) { max_cand_vec_length = cand_vec_length[a]; }
690:   }
691:   /* destroy logical_agg, no longer needed */
692:   SafeMatDestroy(&logical_agg);

694:   /* get the entries for aggregate from B */
695:   MatGetSubMatrices(asa_lev->B, mat_agg_loc_size, idxm_is_B_arr, idxn_is_B_arr, MAT_INITIAL_MATRIX, &b_submat_arr);
696: 
697:   /* clean up all the index sets */
698:   for (a=0; a<mat_agg_loc_size; a++) { ISDestroy(idxm_is_B_arr[a]); }
699:   PetscFree(idxm_is_B_arr);
700:   ISDestroy(idxn_is_B);
701:   PetscFree(idxn_is_B_arr);
702: 
703:   /* storage for the values from each submatrix */
704:   PetscMalloc(sizeof(PetscScalar)*max_cand_vec_length*cand_vecs_num, &b_submat);
705:   PetscMalloc(sizeof(PetscScalar)*max_cand_vec_length*cand_vecs_num, &b_submat_tp);
706:   PetscMalloc(sizeof(PetscInt)*max_cand_vec_length, &idxm);
707:   for (i=0; i<max_cand_vec_length; i++) { idxm[i] = i; }
708:   PetscMalloc(sizeof(PetscInt)*cand_vecs_num, &idxn);
709:   for (i=0; i<cand_vecs_num; i++) { idxn[i] = i; }
710:   /* work storage for QR algorithm */
711:   PetscMalloc(sizeof(PetscScalar)*max_cand_vec_length, &tau);
712:   PetscMalloc(sizeof(PetscScalar)*cand_vecs_num, &work);

714:   /* orthogonalize all submatrices and store them in b_orth_arr */
715:   for (a=0; a<mat_agg_loc_size; a++) {
716:        /* Get the entries for aggregate from B. This is row ordered (although internally
717:           it is column ordered and we will waste some energy transposing it).
718:           WORK: use something like MatGetArray(b_submat_arr[a], &b_submat) but be really
719:           careful about all the different matrix types */
720:        MatGetValues(b_submat_arr[a], cand_vec_length[a], idxm, cand_vecs_num, idxn, b_submat);

722:        if (construct_bridge) {
723:          /* if we are constructing a bridging restriction/interpolation operator, we have
724:             to use the same number of dofs as in our previous construction */
725:          max_cand_vecs_to_add = asa_lev->loc_agg_dofs[a];
726:        } else {
727:          /* for a normal restriction/interpolation operator, we should make sure that we
728:             do not create linear dependence by accident */
729:          max_cand_vecs_to_add = PetscMin(cand_vec_length[a], cand_vecs_num);
730:        }

732:        /* We use LAPACK to compute the QR decomposition of b_submat. For LAPACK we have to
733:           transpose the matrix. We might throw out some column vectors during this process.
734:           We are keeping count of the number of column vectors that we use (and therefore the
735:           number of dofs on the lower level) in new_loc_agg_dofs[a]. */
736:        new_loc_agg_dofs[a] = 0;
737:        for (j=0; j<max_cand_vecs_to_add; j++) {
738:          /* check for condition (4.5) */
739:          norm = 0.0;
740:          for (i=0; i<cand_vec_length[a]; i++) {
741:            norm += PetscRealPart(b_submat[i*cand_vecs_num+j])*PetscRealPart(b_submat[i*cand_vecs_num+j])
742:              + PetscImaginaryPart(b_submat[i*cand_vecs_num+j])*PetscImaginaryPart(b_submat[i*cand_vecs_num+j]);
743:          }
744:          /* only add candidate vector if bigger than cutoff or first candidate */
745:          if ((!j) || (norm > cutoff*((PetscReal) cand_vec_length[a])/((PetscReal) nodes_on_lev))) {
746:            /* passed criterion (4.5), we have not implemented criterion (4.6) yet */
747:            for (i=0; i<cand_vec_length[a]; i++) {
748:              b_submat_tp[new_loc_agg_dofs[a]*cand_vec_length[a]+i] = b_submat[i*cand_vecs_num+j];
749:            }
750:            new_loc_agg_dofs[a]++;
751:          }
752:          /* #ifdef PCASA_VERBOSE */
753:          else {
754:            PetscPrintf(asa_lev->comm, "Cutoff criteria invoked\n");
755:          }
756:          /* #endif */
757:        }

759:        CHKMEMQ;
760:        /* orthogonalize b_submat_tp using the QR algorithm from LAPACK */
761:        LAPACKgeqrf_(cand_vec_length+a, new_loc_agg_dofs+a, b_submat_tp, cand_vec_length+a, tau, work, new_loc_agg_dofs+a, &info);
762:        if (info) SETERRQ(PETSC_ERR_LIB, "LAPACKgeqrf_ LAPACK routine failed");
763: #if !defined(PETSC_MISSING_LAPACK_ORGQR) 
764:        LAPACKungqr_(cand_vec_length+a, new_loc_agg_dofs+a, new_loc_agg_dofs+a, b_submat_tp, cand_vec_length+a, tau, work, new_loc_agg_dofs+a, &info);
765: #else
766:        SETERRQ(PETSC_ERR_SUP,"ORGQR - Lapack routine is unavailable\nIf linking with ESSL you MUST also link with full LAPACK, for example\nuse config/configure.py with --with-blas-lib=libessl.a --with-lapack-lib=/usr/local/lib/liblapack.a'");
767: #endif
768:        if (info) SETERRQ(PETSC_ERR_LIB, "LAPACKungqr_ LAPACK routine failed");

770:        /* Transpose b_submat_tp and store it in b_orth_arr[a]. If we are constructing a
771:           bridging restriction/interpolation operator, we could end up with less dofs than
772:           we previously had. We fill those up with zeros. */
773:        if (!construct_bridge) {
774:          PetscMalloc(sizeof(PetscScalar)*cand_vec_length[a]*new_loc_agg_dofs[a], b_orth_arr+a);
775:          for (j=0; j<new_loc_agg_dofs[a]; j++) {
776:            for (i=0; i<cand_vec_length[a]; i++) {
777:              b_orth_arr[a][i*new_loc_agg_dofs[a]+j] = b_submat_tp[j*cand_vec_length[a]+i];
778:            }
779:          }
780:        } else {
781:          /* bridge, might have to fill up */
782:          PetscMalloc(sizeof(PetscScalar)*cand_vec_length[a]*max_cand_vecs_to_add, b_orth_arr+a);
783:          for (j=0; j<new_loc_agg_dofs[a]; j++) {
784:            for (i=0; i<cand_vec_length[a]; i++) {
785:              b_orth_arr[a][i*max_cand_vecs_to_add+j] = b_submat_tp[j*cand_vec_length[a]+i];
786:            }
787:          }
788:          for (j=new_loc_agg_dofs[a]; j<max_cand_vecs_to_add; j++) {
789:            for (i=0; i<cand_vec_length[a]; i++) {
790:              b_orth_arr[a][i*max_cand_vecs_to_add+j] = 0.0;
791:            }
792:          }
793:          new_loc_agg_dofs[a] = max_cand_vecs_to_add;
794:        }
795:        /* the number of columns in asa_lev->P that are local to this process */
796:        total_loc_cols += new_loc_agg_dofs[a];
797:   } /* end of loop over local aggregates */

799:   /* destroy the submatrices, also frees all allocated space */
800:   MatDestroyMatrices(mat_agg_loc_size, &b_submat_arr);
801:   /* destroy all other workspace */
802:   PetscFree(b_submat);
803:   PetscFree(b_submat_tp);
804:   PetscFree(idxm);
805:   PetscFree(idxn);
806:   PetscFree(tau);
807:   PetscFree(work);

809:   /* destroy old matrix P, Pt */
810:   SafeMatDestroy(&(asa_lev->P));
811:   SafeMatDestroy(&(asa_lev->Pt));

813:   MatGetLocalSize(asa_lev->A, &a_loc_m, &a_loc_n);

815:   /* determine local range */
816:   MPI_Comm_size(asa_lev->comm, &comm_size);
817:   MPI_Comm_rank(asa_lev->comm, &comm_rank);
818:   PetscMalloc(comm_size*sizeof(PetscInt), &loc_cols);
819:   MPI_Allgather(&total_loc_cols, 1, MPI_INT, loc_cols, 1, MPI_INT, asa_lev->comm);
820:   mat_loc_col_start = 0;
821:   for (i=0;i<comm_rank;i++) {
822:     mat_loc_col_start += loc_cols[i];
823:   }
824:   mat_loc_col_end = mat_loc_col_start + loc_cols[i];
825:   mat_loc_col_size = mat_loc_col_end-mat_loc_col_start;
826:   if (mat_loc_col_size != total_loc_cols) SETERRQ(PETSC_ERR_COR, "Local size does not match matrix size");
827:   PetscFree(loc_cols);

829:   /* we now have enough information to create asa_lev->P */
830:   MatCreateMPIAIJ(asa_lev->comm, a_loc_n,  total_loc_cols, asa_lev->size, PETSC_DETERMINE,
831:                          cand_vecs_num, PETSC_NULL, cand_vecs_num, PETSC_NULL, &(asa_lev->P));
832:   /* create asa_lev->Pt */
833:   MatCreateMPIAIJ(asa_lev->comm, total_loc_cols, a_loc_n, PETSC_DETERMINE, asa_lev->size,
834:                          max_cand_vec_length, PETSC_NULL, max_cand_vec_length, PETSC_NULL, &(asa_lev->Pt));
835:   if (asa_lev->next) {
836:     /* create correlator for aggregates of next level */
837:     MatCreateMPIAIJ(asa_lev->comm, mat_agg_loc_size, total_loc_cols, PETSC_DETERMINE, PETSC_DETERMINE,
838:                            cand_vecs_num, PETSC_NULL, cand_vecs_num, PETSC_NULL, &(asa_lev->next->agg_corr));
839:     /* create asa_lev->next->bridge_corr matrix */
840:     MatCreateMPIAIJ(asa_lev->comm, mat_agg_loc_size, total_loc_cols, PETSC_DETERMINE, PETSC_DETERMINE,
841:                            cand_vecs_num, PETSC_NULL, cand_vecs_num, PETSC_NULL, &(asa_lev->next->bridge_corr));
842:   }

844:   /* this is my own hack, but it should give the columns that we should write to */
845:   MatGetOwnershipRangeColumn(asa_lev->P, &mat_loc_col_start, &mat_loc_col_end);
846:   mat_loc_col_size = mat_loc_col_end-mat_loc_col_start;
847:   if (mat_loc_col_size != total_loc_cols) SETERRQ(PETSC_ERR_ARG_SIZ, "The number of local columns in asa_lev->P assigned to this processor does not match the local vector size");

849:   loc_agg_dofs_sum = 0;
850:   /* construct P, Pt, agg_corr, bridge_corr */
851:   for (a=0; a<mat_agg_loc_size; a++) {
852:     /* store b_orth_arr[a] in P */
853:     for (i=0; i<cand_vec_length[a]; i++) {
854:       row = agg_arr[a][i];
855:       for (j=0; j<new_loc_agg_dofs[a]; j++) {
856:         col = mat_loc_col_start + loc_agg_dofs_sum + j;
857:         val = b_orth_arr[a][i*new_loc_agg_dofs[a] + j];
858:         MatSetValues(asa_lev->P, 1, &row, 1, &col, &val, INSERT_VALUES);
859:         val = PetscConj(val);
860:         MatSetValues(asa_lev->Pt, 1, &col, 1, &row, &val, INSERT_VALUES);
861:       }
862:     }

864:     /* compute aggregate correlation matrices */
865:     if (asa_lev->next) {
866:       row = a+mat_agg_loc_start;
867:       for (i=0; i<new_loc_agg_dofs[a]; i++) {
868:         col = mat_loc_col_start + loc_agg_dofs_sum + i;
869:         val = 1.0;
870:         MatSetValues(asa_lev->next->agg_corr, 1, &row, 1, &col, &val, INSERT_VALUES);
871:         /* for the bridge operator we leave out the newest candidates, i.e.
872:            we set bridge_corr to 1.0 for all columns up to asa_lev->loc_agg_dofs[a] and to
873:            0.0 between asa_lev->loc_agg_dofs[a] and new_loc_agg_dofs[a] */
874:         if (!(asa_lev->loc_agg_dofs && (i >= asa_lev->loc_agg_dofs[a]))) {
875:           MatSetValues(asa_lev->next->bridge_corr, 1, &row, 1, &col, &val, INSERT_VALUES);
876:         }
877:       }
878:     }

880:     /* move to next entry point col */
881:     loc_agg_dofs_sum += new_loc_agg_dofs[a];
882:   } /* end of loop over local aggregates */

884:   MatAssemblyBegin(asa_lev->P,MAT_FINAL_ASSEMBLY);
885:   MatAssemblyEnd(asa_lev->P,MAT_FINAL_ASSEMBLY);
886:   MatAssemblyBegin(asa_lev->Pt,MAT_FINAL_ASSEMBLY);
887:   MatAssemblyEnd(asa_lev->Pt,MAT_FINAL_ASSEMBLY);
888:   if (asa_lev->next) {
889:     MatAssemblyBegin(asa_lev->next->agg_corr,MAT_FINAL_ASSEMBLY);
890:     MatAssemblyEnd(asa_lev->next->agg_corr,MAT_FINAL_ASSEMBLY);
891:     MatAssemblyBegin(asa_lev->next->bridge_corr,MAT_FINAL_ASSEMBLY);
892:     MatAssemblyEnd(asa_lev->next->bridge_corr,MAT_FINAL_ASSEMBLY);
893:   }

895:   /* if we are not constructing a bridging operator, switch asa_lev->loc_agg_dofs
896:      and new_loc_agg_dofs */
897:   if (construct_bridge) {
898:     PetscFree(new_loc_agg_dofs);
899:   } else {
900:     if (asa_lev->loc_agg_dofs) {
901:       PetscFree(asa_lev->loc_agg_dofs);
902:     }
903:     asa_lev->loc_agg_dofs = new_loc_agg_dofs;
904:   }

906:   /* clean up */
907:   for (a=0; a<mat_agg_loc_size; a++) {
908:     PetscFree(b_orth_arr[a]);
909:     PetscFree(agg_arr[a]);
910:   }
911:   PetscFree(cand_vec_length);
912:   PetscFree(b_orth_arr);
913:   PetscFree(agg_arr);

916:   return(0);
917: }

919: /* -------------------------------------------------------------------------- */
920: /*
921:    PCSmoothProlongator_ASA - Computes the smoothed prolongators I and It on the level

923:    Input Parameters:
924: .  asa_lev - the level for which the smoothed prolongator is constructed
925: */
928: PetscErrorCode PCSmoothProlongator_ASA(PC_ASA_level *asa_lev)
929: {

933:   SafeMatDestroy(&(asa_lev->smP));
934:   SafeMatDestroy(&(asa_lev->smPt));
935:   /* compute prolongator I_{l+1}^l = S_l P_{l+1}^l */
936:   /* step 1: compute I_{l+1}^l = A_l P_{l+1}^l */
937:   MatMatMult(asa_lev->A, asa_lev->P, MAT_INITIAL_MATRIX, 1, &(asa_lev->smP));
938:   MatMatMult(asa_lev->Pt, asa_lev->A, MAT_INITIAL_MATRIX, 1, &(asa_lev->smPt));
939:   /* step 2: shift and scale to get I_{l+1}^l = P_{l+1}^l - 4/(3/rho) A_l P_{l+1}^l */
940:   MatAYPX(asa_lev->smP, -4./(3.*(asa_lev->spec_rad)), asa_lev->P, SUBSET_NONZERO_PATTERN);
941:   MatAYPX(asa_lev->smPt, -4./(3.*(asa_lev->spec_rad)), asa_lev->Pt, SUBSET_NONZERO_PATTERN);

943:   return(0);
944: }


947: /* -------------------------------------------------------------------------- */
948: /*
949:    PCCreateVcycle_ASA - Creates the V-cycle, when aggregates are already defined

951:    Input Parameters:
952: .  asa - the preconditioner context
953: */
956: PetscErrorCode PCCreateVcycle_ASA(PC_ASA *asa)
957: {
959:   PC_ASA_level   *asa_lev, *asa_next_lev;
960:   Mat            AI;


965:   if (!asa) SETERRQ(PETSC_ERR_ARG_NULL, "asa pointer is NULL");
966:   if (!(asa->levellist)) SETERRQ(PETSC_ERR_ARG_NULL, "no levels found");
967:   asa_lev = asa->levellist;
968:   PCComputeSpectralRadius_ASA(asa_lev);
969:   PCSetupSmoothersOnLevel_ASA(asa, asa_lev, asa->nu);

971:   while(asa_lev->next) {
972:     asa_next_lev = asa_lev->next;
973:     /* (a) aggregates are already constructed */

975:     /* (b) construct B_{l+1} and P_{l+1}^l using (2.11) */
976:     /* construct P_{l+1}^l */
977:     PCCreateTransferOp_ASA(asa_lev, PETSC_FALSE);

979:     /* construct B_{l+1} */
980:     SafeMatDestroy(&(asa_next_lev->B));
981:     MatMatMult(asa_lev->Pt, asa_lev->B, MAT_INITIAL_MATRIX, 1, &(asa_next_lev->B));
982:     asa_next_lev->cand_vecs = asa_lev->cand_vecs;

984:     /* (c) construct smoothed prolongator */
985:     PCSmoothProlongator_ASA(asa_lev);
986: 
987:     /* (d) construct coarse matrix */
988:     /* Define coarse matrix A_{l+1} = (I_{l+1}^l)^T A_l I_{l+1}^l */
989:     SafeMatDestroy(&(asa_next_lev->A));
990:        MatMatMult(asa_lev->A, asa_lev->smP, MAT_INITIAL_MATRIX, 1.0, &AI);
991:      MatMatMult(asa_lev->smPt, AI, MAT_INITIAL_MATRIX, 1.0, &(asa_next_lev->A));
992:      SafeMatDestroy(&AI);
993:     /*     MatPtAP(asa_lev->A, asa_lev->smP, MAT_INITIAL_MATRIX, 1, &(asa_next_lev->A)); */
994:     MatGetSize(asa_next_lev->A, PETSC_NULL, &(asa_next_lev->size));
995:     PCComputeSpectralRadius_ASA(asa_next_lev);
996:     PCSetupSmoothersOnLevel_ASA(asa, asa_next_lev, asa->nu);
997:     /* create corresponding vectors x_{l+1}, b_{l+1}, r_{l+1} */
998:     SafeVecDestroy(&(asa_next_lev->x));
999:     SafeVecDestroy(&(asa_next_lev->b));
1000:     SafeVecDestroy(&(asa_next_lev->r));
1001:     MatGetVecs(asa_next_lev->A, &(asa_next_lev->x), &(asa_next_lev->b));
1002:     MatGetVecs(asa_next_lev->A, PETSC_NULL, &(asa_next_lev->r));

1004:     /* go to next level */
1005:     asa_lev = asa_lev->next;
1006:   } /* end of while loop over the levels */
1007:   /* asa_lev now points to the coarsest level, set up direct solver there */
1008:   PCComputeSpectralRadius_ASA(asa_lev);
1009:   PCSetupDirectSolversOnLevel_ASA(asa, asa_lev, asa->nu);

1012:   return(0);
1013: }

1015: /* -------------------------------------------------------------------------- */
1016: /*
1017:    PCAddCandidateToB_ASA - Inserts a candidate vector in B

1019:    Input Parameters:
1020: +  B - the matrix to insert into
1021: .  col_idx - the column we should insert to
1022: .  x - the vector to insert
1023: -  A - system matrix

1025:    Function will insert normalized x into B, such that <A x, x> = 1
1026:    (x itself is not changed). If B is projected down then this property
1027:    is kept. If <A_l x_l, x_l> = 1 and the next level is defined by
1028:    x_{l+1} = Pt x_l  and  A_{l+1} = Pt A_l P then
1029:    <A_{l+1} x_{l+1}, x_l> = <Pt A_l P Pt x_l, Pt x_l>
1030:    = <A_l P Pt x_l, P Pt x_l> = <A_l x_l, x_l> = 1
1031:    because of the definition of P in (2.11).
1032: */
1035: PetscErrorCode PCAddCandidateToB_ASA(Mat B, PetscInt col_idx, Vec x, Mat A)
1036: {
1038:   Vec            Ax;
1039:   PetscScalar    dotprod;
1040:   PetscReal      norm;
1041:   PetscInt       i, loc_start, loc_end;
1042:   PetscScalar    val, *vecarray;

1045:   MatGetVecs(A, PETSC_NULL, &Ax);
1046:   MatMult(A, x, Ax);
1047:   VecDot(Ax, x, &dotprod);
1048:   norm = PetscAbsScalar(PetscSqrtScalar(PetscAbsScalar(dotprod))); /* there has to be a better way */
1049:   VecGetOwnershipRange(x, &loc_start, &loc_end);
1050:   VecGetArray(x, &vecarray);
1051:   for (i=loc_start; i<loc_end; i++) {
1052:     val = vecarray[i-loc_start]/norm;
1053:     MatSetValues(B, 1, &i, 1, &col_idx, &val, INSERT_VALUES);
1054:   }
1055:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
1056:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
1057:   VecRestoreArray(x, &vecarray);
1058:   VecDestroy(Ax);
1059:   return(0);
1060: }

1062: /* -------------------------------------------------------------------------- */
1063: /*
1064: -  x - a starting guess for a hard to approximate vector, if PETSC_NULL, will be generated
1065: */
1068: PetscErrorCode PCInitializationStage_ASA(PC_ASA *asa, Vec x)
1069: {
1071:   PetscInt       l;
1072:   PC_ASA_level   *asa_lev, *asa_next_lev;
1073:   PetscRandom    rctx;     /* random number generator context */

1075:   Vec            ax;
1076:   PetscScalar    tmp;
1077:   PetscReal      prevnorm, norm;

1079:   PetscTruth     skip_steps_f_i = PETSC_FALSE;
1080:   PetscTruth     sufficiently_coarsened = PETSC_FALSE;

1082:   PetscInt       vec_size, vec_loc_size;
1083:   PetscInt       loc_vec_low, loc_vec_high;
1084:   PetscInt       i,j;

1086: /*   Vec            xhat = 0; */

1088:   Mat            AI;

1090:   Vec            cand_vec, cand_vec_new;
1091:   PetscTruth     isrichardson;
1092:   PC             coarse_pc;

1096:   l=1;
1097:   /* create first level */
1098:   PCCreateLevel_ASA(&(asa->levellist), l, asa->comm, 0, 0, asa->ksptype_smooth, asa->pctype_smooth);
1099:   asa_lev = asa->levellist;

1101:   /* Set matrix */
1102:   asa_lev->A = asa->A;
1103:   MatGetSize(asa_lev->A, &i, &j);
1104:   asa_lev->size = i;
1105:   PCComputeSpectralRadius_ASA(asa_lev);
1106:   PCSetupSmoothersOnLevel_ASA(asa, asa_lev, asa->mu_initial);

1108:   /* Set DM */
1109:   asa_lev->dm = asa->dm;
1110:   PetscObjectReference((PetscObject)asa->dm);

1112:   PetscPrintf(asa_lev->comm, "Initialization stage\n");

1114:   if (x) {
1115:     /* use starting guess */
1116:     SafeVecDestroy(&(asa_lev->x));
1117:     VecDuplicate(x, &(asa_lev->x));
1118:     VecCopy(x, asa_lev->x);
1119:   } else {
1120:     /* select random starting vector */
1121:     SafeVecDestroy(&(asa_lev->x));
1122:     MatGetVecs(asa_lev->A, &(asa_lev->x), 0);
1123:     PetscRandomCreate(asa_lev->comm,&rctx);
1124:     PetscRandomSetFromOptions(rctx);
1125:     VecSetRandom(asa_lev->x, rctx);
1126:     PetscRandomDestroy(rctx);
1127:   }

1129:   /* create right hand side */
1130:   SafeVecDestroy(&(asa_lev->b));
1131:   MatGetVecs(asa_lev->A, &(asa_lev->b), 0);
1132:   VecSet(asa_lev->b, 0.0);

1134:   /* relax and check whether that's enough already */
1135:   /* compute old norm */
1136:   MatGetVecs(asa_lev->A, 0, &ax);
1137:   MatMult(asa_lev->A, asa_lev->x, ax);
1138:   VecDot(asa_lev->x, ax, &tmp);
1139:   prevnorm = PetscAbsScalar(tmp);
1140:   PetscPrintf(asa_lev->comm, "Residual norm of starting guess: %f\n", prevnorm);

1142:   /* apply mu_initial relaxations */
1143:   KSPSolve(asa_lev->smoothd, asa_lev->b, asa_lev->x);
1144:   /* compute new norm */
1145:   MatMult(asa_lev->A, asa_lev->x, ax);
1146:   VecDot(asa_lev->x, ax, &tmp);
1147:   norm = PetscAbsScalar(tmp);
1148:   SafeVecDestroy(&(ax));
1149:   PetscPrintf(asa_lev->comm, "Residual norm of relaxation after %g %d relaxations: %g %g\n", asa->epsilon,asa->mu_initial, norm,prevnorm);

1151:   /* Check if it already converges by itself */
1152:   if (norm/prevnorm <= PetscAbsScalar(PetscPowScalar(asa->epsilon, asa->mu_initial))) {
1153:     /* converges by relaxation alone */
1154:     SETERRQ(PETSC_ERR_SUP, "Relaxation should be sufficient to treat this problem. "
1155:             "Use relaxation or decrease epsilon with -pc_asa_epsilon");
1156:   } else {
1157:     /* set the number of relaxations to asa->mu from asa->mu_initial */
1158:     PCSetupSmoothersOnLevel_ASA(asa, asa_lev, asa->mu);

1160:     /* Let's do some multigrid ! */
1161:     sufficiently_coarsened = PETSC_FALSE;

1163:     /* do the whole initialization stage loop */
1164:     while (!sufficiently_coarsened) {
1165:       PetscPrintf(asa_lev->comm, "Initialization stage: creating level %d\n", asa_lev->level+1);

1167:       /* (a) Set candidate matrix B_l = x_l */
1168:       /* get the correct vector sizes and data */
1169:       VecGetSize(asa_lev->x, &vec_size);
1170:       VecGetOwnershipRange(asa_lev->x, &loc_vec_low, &loc_vec_high);
1171:       vec_loc_size = loc_vec_high - loc_vec_low;

1173:       /* create matrix for candidates */
1174:       MatCreateMPIDense(asa_lev->comm, vec_loc_size, PETSC_DECIDE, vec_size, asa->max_cand_vecs, PETSC_NULL, &(asa_lev->B));
1175:       /* set the first column */
1176:       PCAddCandidateToB_ASA(asa_lev->B, 0, asa_lev->x, asa_lev->A);
1177:       asa_lev->cand_vecs = 1;

1179:       /* create next level */
1180:       PCCreateLevel_ASA(&(asa_lev->next), asa_lev->level+1,  asa_lev->comm, asa_lev, PETSC_NULL, asa->ksptype_smooth, asa->pctype_smooth);
1181:       asa_next_lev = asa_lev->next;

1183:       /* (b) Create nodal aggregates A_i^l */
1184:       PCCreateAggregates_ASA(asa_lev);
1185: 
1186:       /* (c) Define tentatative prolongator P_{l+1}^l and candidate matrix B_{l+1}
1187:              using P_{l+1}^l B_{l+1} = B_l and (P_{l+1}^l)^T P_{l+1}^l = I */
1188:       PCCreateTransferOp_ASA(asa_lev, PETSC_FALSE);

1190:       /* future WORK: set correct fill ratios for all the operations below */
1191:       MatMatMult(asa_lev->Pt, asa_lev->B, MAT_INITIAL_MATRIX, 1, &(asa_next_lev->B));
1192:       asa_next_lev->cand_vecs = asa_lev->cand_vecs;

1194:       /* (d) Define prolongator I_{l+1}^l = S_l P_{l+1}^l */
1195:       PCSmoothProlongator_ASA(asa_lev);

1197:       /* (e) Define coarse matrix A_{l+1} = (I_{l+1}^l)^T A_l I_{l+1}^l */
1198:             MatMatMult(asa_lev->A, asa_lev->smP, MAT_INITIAL_MATRIX, 1.0, &AI);
1199:       MatMatMult(asa_lev->smPt, AI, MAT_INITIAL_MATRIX, 1.0, &(asa_next_lev->A));
1200:       SafeMatDestroy(&AI);
1201:       /*      MatPtAP(asa_lev->A, asa_lev->smP, MAT_INITIAL_MATRIX, 1, &(asa_next_lev->A)); */
1202:       MatGetSize(asa_next_lev->A, PETSC_NULL, &(asa_next_lev->size));
1203:       PCComputeSpectralRadius_ASA(asa_next_lev);
1204:       PCSetupSmoothersOnLevel_ASA(asa, asa_next_lev, asa->mu);

1206:       /* coarse enough for direct solver? */
1207:       MatGetSize(asa_next_lev->A, &i, &j);
1208:       if (PetscMax(i,j) <= asa->direct_solver) {
1209:         PetscPrintf(asa_lev->comm, "Level %d can be treated directly.\n"
1210:                            "Algorithm will use %d levels.\n", asa_next_lev->level,
1211:                            asa_next_lev->level);
1212:         break; /* go to step 5 */
1213:       }

1215:       if (skip_steps_f_i == PETSC_FALSE) {
1216:         /* (f) Set x_{l+1} = B_{l+1}, we just compute it again */
1217:         SafeVecDestroy(&(asa_next_lev->x));
1218:         MatGetVecs(asa_lev->P, &(asa_next_lev->x), 0);
1219:         MatMult(asa_lev->Pt, asa_lev->x, asa_next_lev->x);

1221: /*         /\* (g) Make copy \hat{x}_{l+1} = x_{l+1} *\/ */
1222: /*         VecDuplicate(asa_next_lev->x, &xhat); */
1223: /*         VecCopy(asa_next_lev->x, xhat); */
1224: 
1225:         /* Create b_{l+1} */
1226:         SafeVecDestroy(&(asa_next_lev->b));
1227:         MatGetVecs(asa_next_lev->A, &(asa_next_lev->b), 0);
1228:         VecSet(asa_next_lev->b, 0.0);

1230:         /* (h) Relax mu times on A_{l+1} x = 0 */
1231:         /* compute old norm */
1232:         MatGetVecs(asa_next_lev->A, 0, &ax);
1233:         MatMult(asa_next_lev->A, asa_next_lev->x, ax);
1234:         VecDot(asa_next_lev->x, ax, &tmp);
1235:         prevnorm = PetscAbsScalar(tmp);
1236:         PetscPrintf(asa_next_lev->comm, "Residual norm of starting guess on level %d: %f\n", asa_next_lev->level, prevnorm);
1237:         /* apply mu relaxations: WORK, make sure that mu is set correctly */
1238:         KSPSolve(asa_next_lev->smoothd, asa_next_lev->b, asa_next_lev->x);
1239:         /* compute new norm */
1240:         MatMult(asa_next_lev->A, asa_next_lev->x, ax);
1241:         VecDot(asa_next_lev->x, ax, &tmp);
1242:         norm = PetscAbsScalar(tmp);
1243:         SafeVecDestroy(&(ax));
1244:         PetscPrintf(asa_next_lev->comm, "Residual norm after Richardson iteration  on level %d: %f\n", asa_next_lev->level, norm);
1245:         /* (i) Check if it already converges by itself */
1246:         if (norm/prevnorm <= PetscAbsScalar(PetscPowScalar(asa->epsilon, asa->mu))) {
1247:           /* relaxation reduces error sufficiently */
1248:           skip_steps_f_i = PETSC_TRUE;
1249:         }
1250:       }
1251:       /* (j) go to next coarser level */
1252:       l++;
1253:       asa_lev = asa_next_lev;
1254:     }
1255:     /* Step 5. */
1256:     asa->levels = asa_next_lev->level; /* WORK: correct? */

1258:     /* Set up direct solvers on coarsest level */
1259:     if (asa_next_lev->smoothd != asa_next_lev->smoothu) {
1260:       if (asa_next_lev->smoothu) { KSPDestroy(asa_next_lev->smoothu); }
1261:     }
1262:     KSPSetType(asa_next_lev->smoothd, asa->ksptype_direct);
1263:     PetscTypeCompare((PetscObject)(asa_next_lev->smoothd), KSPRICHARDSON, &isrichardson);
1264:     if (isrichardson) {
1265:       KSPSetInitialGuessNonzero(asa_next_lev->smoothd, PETSC_TRUE);
1266:     } else {
1267:       KSPSetInitialGuessNonzero(asa_next_lev->smoothd, PETSC_FALSE);
1268:     }
1269:     KSPGetPC(asa_next_lev->smoothd, &coarse_pc);
1270:     PCSetType(coarse_pc, asa->pctype_direct);
1271:     asa_next_lev->smoothu = asa_next_lev->smoothd;
1272:     PCSetupDirectSolversOnLevel_ASA(asa, asa_next_lev, asa->nu);

1274:     /* update finest-level candidate matrix B_1 = I_2^1 I_3^2 ... I_{L-1}^{L-2} x_{L-1} */
1275:     if (!asa_lev->prev) {
1276:       /* just one relaxation level */
1277:       VecDuplicate(asa_lev->x, &cand_vec);
1278:       VecCopy(asa_lev->x, cand_vec);
1279:     } else {
1280:       /* interpolate up the chain */
1281:       cand_vec = asa_lev->x;
1282:       asa_lev->x = 0;
1283:       while(asa_lev->prev) {
1284:         /* interpolate to higher level */
1285:         MatGetVecs(asa_lev->prev->smP, 0, &cand_vec_new);
1286:         MatMult(asa_lev->prev->smP, cand_vec, cand_vec_new);
1287:         SafeVecDestroy(&(cand_vec));
1288:         cand_vec = cand_vec_new;
1289: 
1290:         /* destroy all working vectors on the way */
1291:         SafeVecDestroy(&(asa_lev->x));
1292:         SafeVecDestroy(&(asa_lev->b));

1294:         /* move to next higher level */
1295:         asa_lev = asa_lev->prev;
1296:       }
1297:     }
1298:     /* set the first column of B1 */
1299:     PCAddCandidateToB_ASA(asa_lev->B, 0, cand_vec, asa_lev->A);
1300:     SafeVecDestroy(&(cand_vec));

1302:     /* Step 6. Create V-cycle */
1303:     PCCreateVcycle_ASA(asa);
1304:   }
1306:   return(0);
1307: }

1309: /* -------------------------------------------------------------------------- */
1310: /*
1311:    PCApplyVcycleOnLevel_ASA - Applies current V-cycle

1313:    Input Parameters:
1314: +  asa_lev - the current level we should recurse on
1315: -  gamma - the number of recursive cycles we should run

1317: */
1320: PetscErrorCode PCApplyVcycleOnLevel_ASA(PC_ASA_level *asa_lev, PetscInt gamma)
1321: {
1323:   PC_ASA_level   *asa_next_lev;
1324:   PetscInt       g;

1327:   if (!asa_lev) SETERRQ(PETSC_ERR_ARG_NULL, "Level is empty in PCApplyVcycleOnLevel_ASA");
1328:   asa_next_lev = asa_lev->next;

1330:   if (asa_next_lev) {
1331:     /* 1. Presmoothing */
1332:     KSPSolve(asa_lev->smoothd, asa_lev->b, asa_lev->x);
1333:     /* 2. Coarse grid corrections */
1334: /*     MatGetVecs(asa_lev->A, 0, &tmp); */
1335: /*     MatGetVecs(asa_lev->smP, &(asa_next_lev->b), 0); */
1336: /*     MatGetVecs(asa_next_lev->A, &(asa_next_lev->x), 0); */
1337:     for (g=0; g<gamma; g++) {
1338:       /* (a) get coarsened b_{l+1} = (I_{l+1}^l)^T (b_l - A_l x_l) */
1339:       MatMult(asa_lev->A, asa_lev->x, asa_lev->r);
1340:       VecAYPX(asa_lev->r, -1.0, asa_lev->b);
1341:       MatMult(asa_lev->smPt, asa_lev->r, asa_next_lev->b);

1343:       /* (b) Set x_{l+1} = 0 and recurse */
1344:       VecSet(asa_next_lev->x, 0.0);
1345:       PCApplyVcycleOnLevel_ASA(asa_next_lev, gamma);

1347:       /* (c) correct solution x_l = x_l + I_{l+1}^l x_{l+1} */
1348:       MatMultAdd(asa_lev->smP, asa_next_lev->x, asa_lev->x, asa_lev->x);
1349:     }
1350: /*     SafeVecDestroy(&(asa_lev->r)); */
1351: /*     /\* discard x_{l+1}, b_{l+1} *\/ */
1352: /*     SafeVecDestroy(&(asa_next_lev->x)); */
1353: /*     SafeVecDestroy(&(asa_next_lev->b)); */
1354: 
1355:     /* 3. Postsmoothing */
1356:     KSPSolve(asa_lev->smoothu, asa_lev->b, asa_lev->x);
1357:   } else {
1358:     /* Base case: solve directly */
1359:     KSPSolve(asa_lev->smoothd, asa_lev->b, asa_lev->x);
1360:   }
1361:   return(0);
1362: }


1365: /* -------------------------------------------------------------------------- */
1366: /*
1367:    PCGeneralSetupStage_ASA - Applies the ASA preconditioner to a vector. Algorithm
1368:                              4 from the ASA paper

1370:    Input Parameters:
1371: +  asa - the data structure for the ASA algorithm
1372: -  cand - a possible candidate vector, if PETSC_NULL, will be constructed randomly

1374:    Output Parameters:
1375: .  cand_added - PETSC_TRUE, if new candidate vector added, PETSC_FALSE otherwise
1376: */
1379: PetscErrorCode PCGeneralSetupStage_ASA(PC_ASA *asa, Vec cand, PetscTruth *cand_added)
1380: {
1382:   PC_ASA_level   *asa_lev, *asa_next_lev;

1384:   PetscRandom    rctx;     /* random number generator context */
1385:   PetscReal      r;
1386:   PetscScalar    rs;
1387:   PetscTruth     nd_fast;

1389:   Vec            ax;
1390:   PetscScalar    tmp;
1391:   PetscReal      norm, prevnorm = 0.0;
1392:   PetscInt       c;

1394:   PetscInt       loc_vec_low, loc_vec_high;
1395:   PetscInt       i;

1397:   PetscTruth     skip_steps_d_j = PETSC_FALSE;

1399:   PetscInt       *idxm, *idxn;
1400:   PetscScalar    *v;

1402:   Mat            AI;

1404:   Vec            cand_vec, cand_vec_new;

1407:   *cand_added = PETSC_FALSE;
1408: 
1409:   asa_lev = asa->levellist;
1410:   if (asa_lev == 0) SETERRQ(PETSC_ERR_ARG_NULL, "No levels found in PCGeneralSetupStage_ASA");
1411:   asa_next_lev = asa_lev->next;
1412:   if (asa_next_lev == 0) SETERRQ(PETSC_ERR_ARG_NULL, "Just one level, not implemented yet");
1413: 
1414:   PetscPrintf(asa_lev->comm, "General setup stage\n");


1418:   /* 1. If max. dof per node on level 2 equals K, stop */
1419:   if (asa_next_lev->cand_vecs >= asa->max_dof_lev_2) {
1420:     PetscPrintf(PETSC_COMM_WORLD,
1421:                        "Maximum dof on level 2 reached: %d\n"
1422:                        "Consider increasing this limit by setting it with -pc_asa_max_dof_lev_2\n",
1423:                        asa->max_dof_lev_2);
1424:     return(0);
1425:   }

1427:   /* 2. Create copy of B_1 (skipped, we just replace the last column in step 8.) */
1428: 
1429:   if (!cand) {
1430:     /* 3. Select a random x_1 */
1431:     SafeVecDestroy(&(asa_lev->x));
1432:     MatGetVecs(asa_lev->A, &(asa_lev->x), 0);
1433:     PetscRandomCreate(asa_lev->comm,&rctx);
1434:     PetscRandomSetFromOptions(rctx);
1435:     VecGetOwnershipRange(asa_lev->x, &loc_vec_low, &loc_vec_high);
1436:     for (i=loc_vec_low; i<loc_vec_high; i++) {
1437:       PetscRandomGetValueReal(rctx, &r);
1438:       rs = r;
1439:       VecSetValues(asa_lev->x, 1, &i, &rs, INSERT_VALUES);
1440:     }
1441:     VecAssemblyBegin(asa_lev->x);
1442:     VecAssemblyEnd(asa_lev->x);
1443:     PetscRandomDestroy(rctx);
1444:   } else {
1445:     SafeVecDestroy(&(asa_lev->x));
1446:     VecDuplicate(cand, &(asa_lev->x));
1447:     VecCopy(cand, asa_lev->x);
1448:   }

1450:   /* create right hand side */
1451:   SafeVecDestroy(&(asa_lev->b));
1452:   MatGetVecs(asa_lev->A, &(asa_lev->b), 0);
1453:   VecSet(asa_lev->b, 0.0);
1454: 
1455:   /* Apply mu iterations of current V-cycle */
1456:   nd_fast = PETSC_FALSE;
1457:   MatGetVecs(asa_lev->A, 0, &ax);
1458:   for (c=0; c<asa->mu; c++) {
1459:     PCApplyVcycleOnLevel_ASA(asa_lev, asa->gamma);
1460: 
1461:     MatMult(asa_lev->A, asa_lev->x, ax);
1462:     VecDot(asa_lev->x, ax, &tmp);
1463:     norm = PetscAbsScalar(tmp);
1464:     if (c>0) {
1465:       if (norm/prevnorm < asa->epsilon) {
1466:         nd_fast = PETSC_TRUE;
1467:         break;
1468:       }
1469:     }
1470:     prevnorm = norm;
1471:   }
1472:   SafeVecDestroy(&(ax));

1474:   /* 4. If energy norm decreases sufficiently fast, then stop */
1475:   if (nd_fast) {
1476:     PetscPrintf(asa_lev->comm, "nd_fast is true\n");
1477:     return(0);
1478:   }

1480:   /* 5. Update B_1, by adding new column x_1 */
1481:   if (asa_lev->cand_vecs >= asa->max_cand_vecs) {
1482:     SETERRQ(PETSC_ERR_MEM, "Number of candidate vectors will exceed allocated storage space");
1483:   } else {
1484:     PetscPrintf(asa_lev->comm, "Adding candidate vector %d\n", asa_lev->cand_vecs+1);
1485:   }
1486:   PCAddCandidateToB_ASA(asa_lev->B, asa_lev->cand_vecs, asa_lev->x, asa_lev->A);
1487:   *cand_added = PETSC_TRUE;
1488:   asa_lev->cand_vecs++;

1490:   /* 6. loop over levels */
1491:   while(asa_next_lev && asa_next_lev->next) {
1492:     PetscPrintf(asa_lev->comm, "General setup stage: processing level %d\n", asa_next_lev->level);
1493:     /* (a) define B_{l+1} and P_{l+1}^L */
1494:     /* construct P_{l+1}^l */
1495:     PCCreateTransferOp_ASA(asa_lev, PETSC_FALSE);

1497:     /* construct B_{l+1} */
1498:     SafeMatDestroy(&(asa_next_lev->B));
1499:     MatMatMult(asa_lev->Pt, asa_lev->B, MAT_INITIAL_MATRIX, 1.0, &(asa_next_lev->B));
1500:     /* do not increase asa_next_lev->cand_vecs until step (j) */
1501: 
1502:     /* (b) construct prolongator I_{l+1}^l = S_l P_{l+1}^l */
1503:     PCSmoothProlongator_ASA(asa_lev);
1504: 
1505:     /* (c) construct coarse matrix A_{l+1} = (I_{l+1}^l)^T A_l I_{l+1}^l */
1506:     SafeMatDestroy(&(asa_next_lev->A));
1507:        MatMatMult(asa_lev->A, asa_lev->smP, MAT_INITIAL_MATRIX, 1.0, &AI);
1508:     MatMatMult(asa_lev->smPt, AI, MAT_INITIAL_MATRIX, 1.0, &(asa_next_lev->A));
1509:     SafeMatDestroy(&AI);
1510:                                  /* MatPtAP(asa_lev->A, asa_lev->smP, MAT_INITIAL_MATRIX, 1, &(asa_next_lev->A)); */
1511:     MatGetSize(asa_next_lev->A, PETSC_NULL, &(asa_next_lev->size));
1512:     PCComputeSpectralRadius_ASA(asa_next_lev);
1513:     PCSetupSmoothersOnLevel_ASA(asa, asa_next_lev, asa->mu);

1515:     if (! skip_steps_d_j) {
1516:       /* (d) get vector x_{l+1} from last column in B_{l+1} */
1517:       SafeVecDestroy(&(asa_next_lev->x));
1518:       MatGetVecs(asa_next_lev->B, 0, &(asa_next_lev->x));

1520:       VecGetOwnershipRange(asa_next_lev->x, &loc_vec_low, &loc_vec_high);
1521:       PetscMalloc(sizeof(PetscInt)*(loc_vec_high-loc_vec_low), &idxm);
1522:       for (i=loc_vec_low; i<loc_vec_high; i++)
1523:         idxm[i-loc_vec_low] = i;
1524:       PetscMalloc(sizeof(PetscInt)*1, &idxn);
1525:       idxn[0] = asa_next_lev->cand_vecs;

1527:       PetscMalloc(sizeof(PetscScalar)*(loc_vec_high-loc_vec_low), &v);
1528:       MatGetValues(asa_next_lev->B, loc_vec_high-loc_vec_low, idxm, 1, idxn, v);

1530:       VecSetValues(asa_next_lev->x, loc_vec_high-loc_vec_low, idxm, v, INSERT_VALUES);
1531:       VecAssemblyBegin(asa_next_lev->x);
1532:       VecAssemblyEnd(asa_next_lev->x);

1534:       PetscFree(v);
1535:       PetscFree(idxm);
1536:       PetscFree(idxn);
1537: 
1538:       /* (e) create bridge transfer operator P_{l+2}^{l+1}, by using the previously
1539:          computed candidates */
1540:       PCCreateTransferOp_ASA(asa_next_lev, PETSC_TRUE);

1542:       /* (f) construct bridging prolongator I_{l+2}^{l+1} = S_{l+1} P_{l+2}^{l+1} */
1543:       PCSmoothProlongator_ASA(asa_next_lev);

1545:       /* (g) compute <A_{l+1} x_{l+1}, x_{l+1}> and save it */
1546:       MatGetVecs(asa_next_lev->A, 0, &ax);
1547:       MatMult(asa_next_lev->A, asa_next_lev->x, ax);
1548:       VecDot(asa_next_lev->x, ax, &tmp);
1549:       prevnorm = PetscAbsScalar(tmp);
1550:       SafeVecDestroy(&(ax));

1552:       /* (h) apply mu iterations of current V-cycle */
1553:       /* set asa_next_lev->b */
1554:       SafeVecDestroy(&(asa_next_lev->b));
1555:       SafeVecDestroy(&(asa_next_lev->r));
1556:       MatGetVecs(asa_next_lev->A, &(asa_next_lev->b), &(asa_next_lev->r));
1557:       VecSet(asa_next_lev->b, 0.0);
1558:       /* apply V-cycle */
1559:       for (c=0; c<asa->mu; c++) {
1560:         PCApplyVcycleOnLevel_ASA(asa_next_lev, asa->gamma);
1561:       }

1563:       /* (i) check convergence */
1564:       /* compute <A_{l+1} x_{l+1}, x_{l+1}> and save it */
1565:       MatGetVecs(asa_next_lev->A, 0, &ax);
1566:       MatMult(asa_next_lev->A, asa_next_lev->x, ax);
1567:       VecDot(asa_next_lev->x, ax, &tmp);
1568:       norm = PetscAbsScalar(tmp);
1569:       SafeVecDestroy(&(ax));

1571:       if (norm/prevnorm <= PetscAbsScalar(PetscPowScalar(asa->epsilon, asa->mu))) skip_steps_d_j = PETSC_TRUE;
1572: 
1573:       /* (j) update candidate B_{l+1} */
1574:       PCAddCandidateToB_ASA(asa_next_lev->B, asa_next_lev->cand_vecs, asa_next_lev->x, asa_next_lev->A);
1575:       asa_next_lev->cand_vecs++;
1576:     }
1577:     /* go to next level */
1578:     asa_lev = asa_lev->next;
1579:     asa_next_lev = asa_next_lev->next;
1580:   }

1582:   /* 7. update the fine-level candidate */
1583:   if (! asa_lev->prev) {
1584:     /* just one coarsening level */
1585:     VecDuplicate(asa_lev->x, &cand_vec);
1586:     VecCopy(asa_lev->x, cand_vec);
1587:   } else {
1588:     cand_vec = asa_lev->x;
1589:     asa_lev->x = 0;
1590:     while(asa_lev->prev) {
1591:       /* interpolate to higher level */
1592:       MatGetVecs(asa_lev->prev->smP, 0, &cand_vec_new);
1593:       MatMult(asa_lev->prev->smP, cand_vec, cand_vec_new);
1594:       SafeVecDestroy(&(cand_vec));
1595:       cand_vec = cand_vec_new;

1597:       /* destroy all working vectors on the way */
1598:       SafeVecDestroy(&(asa_lev->x));
1599:       SafeVecDestroy(&(asa_lev->b));

1601:       /* move to next higher level */
1602:       asa_lev = asa_lev->prev;
1603:     }
1604:   }
1605:   /* 8. update B_1 by setting the last column of B_1 */
1606:   PCAddCandidateToB_ASA(asa_lev->B, asa_lev->cand_vecs-1, cand_vec, asa_lev->A);
1607:   SafeVecDestroy(&(cand_vec));

1609:   /* 9. create V-cycle */
1610:   PCCreateVcycle_ASA(asa);
1611: 
1613:   return(0);
1614: }

1616: /* -------------------------------------------------------------------------- */
1617: /*
1618:    PCConstructMultigrid_ASA - creates the multigrid preconditionier, this is a fairly
1619:    involved process, which runs extensive testing to compute good candidate vectors

1621:    Input Parameters:
1622: .  pc - the preconditioner context

1624:  */
1627: PetscErrorCode PCConstructMultigrid_ASA(PC pc)
1628: {
1630:   PC_ASA         *asa = (PC_ASA*)pc->data;
1631:   PC_ASA_level   *asa_lev;
1632:   PetscInt       i, ls, le;
1633:   PetscScalar    *d;
1634:   PetscTruth     zeroflag = PETSC_FALSE;
1635:   PetscReal      rnorm, rnorm_start;
1636:   PetscReal      rq, rq_prev;
1637:   PetscScalar    rq_nom, rq_denom;
1638:   PetscTruth     cand_added;
1639:   PetscRandom    rctx;


1643:   /* check if we should scale with diagonal */
1644:   if (asa->scale_diag) {
1645:     /* Get diagonal scaling factors */
1646:     MatGetVecs(pc->pmat,&(asa->invsqrtdiag),0);
1647:     MatGetDiagonal(pc->pmat,asa->invsqrtdiag);
1648:     /* compute (inverse) sqrt of diagonal */
1649:     VecGetOwnershipRange(asa->invsqrtdiag, &ls, &le);
1650:     VecGetArray(asa->invsqrtdiag, &d);
1651:     for (i=0; i<le-ls; i++) {
1652:       if (d[i] == 0.0) {
1653:         d[i]     = 1.0;
1654:         zeroflag = PETSC_TRUE;
1655:       } else {
1656:         d[i] = 1./sqrt(PetscAbsScalar(d[i]));
1657:       }
1658:     }
1659:     VecRestoreArray(asa->invsqrtdiag,&d);
1660:     VecAssemblyBegin(asa->invsqrtdiag);
1661:     VecAssemblyEnd(asa->invsqrtdiag);
1662:     if (zeroflag) {
1663:       PetscInfo(pc,"Zero detected in diagonal of matrix, using 1 at those locations\n");
1664:     }
1665: 
1666:     /* scale the matrix and store it: D^{-1/2} A D^{-1/2} */
1667:     MatDuplicate(pc->pmat, MAT_COPY_VALUES, &(asa->A)); /* probably inefficient */
1668:     MatDiagonalScale(asa->A, asa->invsqrtdiag, asa->invsqrtdiag);
1669:   } else {
1670:     /* don't scale */
1671:     asa->A = pc->pmat;
1672:   }
1673:   /* Initialization stage */
1674:   PCInitializationStage_ASA(asa, PETSC_NULL);
1675: 
1676:   /* get first level */
1677:   asa_lev = asa->levellist;

1679:   PetscRandomCreate(asa->comm,&rctx);
1680:   PetscRandomSetFromOptions(rctx);
1681:   VecSetRandom(asa_lev->x,rctx);

1683:   /* compute starting residual */
1684:   SafeVecDestroy(&(asa_lev->r));
1685:   MatGetVecs(asa_lev->A, PETSC_NULL, &(asa_lev->r));
1686:   MatMult(asa_lev->A, asa_lev->x, asa_lev->r);
1687:   /* starting residual norm */
1688:   VecNorm(asa_lev->r, NORM_2, &rnorm_start);
1689:   /* compute Rayleigh quotients */
1690:   VecDot(asa_lev->x, asa_lev->r, &rq_nom);
1691:   VecDot(asa_lev->x, asa_lev->x, &rq_denom);
1692:   rq_prev = PetscAbsScalar(rq_nom / rq_denom);

1694:   /* check if we have to add more candidates */
1695:   for (i=0; i<asa->max_it; i++) {
1696:     if (asa_lev->cand_vecs >= asa->max_cand_vecs) {
1697:       /* reached limit for candidate vectors */
1698:       break;
1699:     }
1700:     /* apply V-cycle */
1701:     PCApplyVcycleOnLevel_ASA(asa_lev, asa->gamma);
1702:     /* check convergence */
1703:     MatMult(asa_lev->A, asa_lev->x, asa_lev->r);
1704:     VecNorm(asa_lev->r, NORM_2, &rnorm);
1705:     PetscPrintf(asa->comm, "After %d iterations residual norm is %f\n", i+1, rnorm);
1706:     if (rnorm < rnorm_start*(asa->rtol) || rnorm < asa->abstol) {
1707:       /* convergence */
1708:       break;
1709:     }
1710:     /* compute new Rayleigh quotient */
1711:     VecDot(asa_lev->x, asa_lev->r, &rq_nom);
1712:     VecDot(asa_lev->x, asa_lev->x, &rq_denom);
1713:     rq = PetscAbsScalar(rq_nom / rq_denom);
1714:     PetscPrintf(asa->comm, "After %d iterations Rayleigh quotient of residual is %f\n", i+1, rq);
1715:     /* test Rayleigh quotient decrease and add more candidate vectors if necessary */
1716:     if (i && (rq > asa->rq_improve*rq_prev)) {
1717:       /* improve interpolation by adding another candidate vector */
1718:       PCGeneralSetupStage_ASA(asa, asa_lev->r, &cand_added);
1719:       if (!cand_added) {
1720:         /* either too many candidates for storage or cycle is already effective */
1721:         PetscPrintf(asa->comm, "either too many candidates for storage or cycle is already effective\n");
1722:         break;
1723:       }
1724:       VecSetRandom(asa_lev->x, rctx);
1725:       rq_prev = rq*10000.; /* give the new V-cycle some grace period */
1726:     } else {
1727:       rq_prev = rq;
1728:     }
1729:   }

1731:   SafeVecDestroy(&(asa_lev->x));
1732:   SafeVecDestroy(&(asa_lev->b));
1733:   PetscRandomDestroy(rctx);
1734:   asa->multigrid_constructed = PETSC_TRUE;
1735:   return(0);
1736: }

1738: /* -------------------------------------------------------------------------- */
1739: /*
1740:    PCApply_ASA - Applies the ASA preconditioner to a vector.

1742:    Input Parameters:
1743: .  pc - the preconditioner context
1744: .  x - input vector

1746:    Output Parameter:
1747: .  y - output vector

1749:    Application Interface Routine: PCApply()
1750:  */
1753: PetscErrorCode PCApply_ASA(PC pc,Vec x,Vec y)
1754: {
1755:   PC_ASA         *asa = (PC_ASA*)pc->data;
1756:   PC_ASA_level   *asa_lev;

1761:   if (!asa->multigrid_constructed) {
1762:     PCConstructMultigrid_ASA(pc);
1763:   }

1765:   /* get first level */
1766:   asa_lev = asa->levellist;

1768:   /* set the right hand side */
1769:   VecDuplicate(x, &(asa->b));
1770:   VecCopy(x, asa->b);
1771:   /* set starting vector */
1772:   SafeVecDestroy(&(asa->x));
1773:   MatGetVecs(asa->A, &(asa->x), PETSC_NULL);
1774:   VecSet(asa->x, 0.0);
1775: 
1776:   /* set vectors */
1777:   asa_lev->x = asa->x;
1778:   asa_lev->b = asa->b;

1780:   PCApplyVcycleOnLevel_ASA(asa_lev, asa->gamma);
1781: 
1782:   /* Return solution */
1783:   VecCopy(asa->x, y);

1785:   /* delete working vectors */
1786:   SafeVecDestroy(&(asa->x));
1787:   SafeVecDestroy(&(asa->b));
1788:   asa_lev->x = PETSC_NULL;
1789:   asa_lev->b = PETSC_NULL;

1791:   return(0);
1792: }

1794: /* -------------------------------------------------------------------------- */
1795: /*
1796:    PCApplyRichardson_ASA - Applies the ASA iteration to solve a linear system

1798:    Input Parameters:
1799: .  pc - the preconditioner context
1800: .  b - the right hand side

1802:    Output Parameter:
1803: .  x - output vector

1805:   DOES NOT WORK!!!!!

1807:  */
1810: PetscErrorCode PCApplyRichardson_ASA(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its)
1811: {
1812:   PC_ASA         *asa = (PC_ASA*)pc->data;
1813:   PC_ASA_level   *asa_lev;
1814:   PetscInt       i;
1815:   PetscReal      rnorm, rnorm_start;

1820:   if (! asa->multigrid_constructed) {
1821:     PCConstructMultigrid_ASA(pc);
1822:   }

1824:   /* get first level */
1825:   asa_lev = asa->levellist;

1827:   /* set the right hand side */
1828:   VecDuplicate(b, &(asa->b));
1829:   if (asa->scale_diag) {
1830:     VecPointwiseMult(asa->b, asa->invsqrtdiag, b);
1831:   } else {
1832:     VecCopy(b, asa->b);
1833:   }
1834:   /* set starting vector */
1835:   VecDuplicate(x, &(asa->x));
1836:   VecCopy(x, asa->x);
1837: 
1838:   /* compute starting residual */
1839:   SafeVecDestroy(&(asa->r));
1840:   MatGetVecs(asa->A, &(asa->r), PETSC_NULL);
1841:   MatMult(asa->A, asa->x, asa->r);
1842:   VecAYPX(asa->r, -1.0, asa->b);
1843:   /* starting residual norm */
1844:   VecNorm(asa->r, NORM_2, &rnorm_start);

1846:   /* set vectors */
1847:   asa_lev->x = asa->x;
1848:   asa_lev->b = asa->b;

1850:   /* **************** Full algorithm loop *********************************** */
1851:   for (i=0; i<its; i++) {
1852:     /* apply V-cycle */
1853:     PCApplyVcycleOnLevel_ASA(asa_lev, asa->gamma);
1854:     /* check convergence */
1855:     MatMult(asa->A, asa->x, asa->r);
1856:     VecAYPX(asa->r, -1.0, asa->b);
1857:     VecNorm(asa->r, NORM_2, &rnorm);
1858:     PetscPrintf(asa->comm, "After %d iterations residual norm is %f\n", i+1, rnorm);
1859:     if (rnorm < rnorm_start*(rtol) || rnorm < asa->abstol) {
1860:       /* convergence */
1861:       break;
1862:     }
1863:     if (rnorm > rnorm_start*(dtol)) {
1864:       /* divergence */
1865:       break;
1866:     }
1867:   }
1868: 
1869:   /* Return solution */
1870:   if (asa->scale_diag) {
1871:     VecPointwiseMult(x, asa->x, asa->invsqrtdiag);
1872:   } else {
1873:     VecCopy(x, asa->x);
1874:   }

1876:   /* delete working vectors */
1877:   SafeVecDestroy(&(asa->x));
1878:   SafeVecDestroy(&(asa->b));
1879:   SafeVecDestroy(&(asa->r));
1880:   asa_lev->x = PETSC_NULL;
1881:   asa_lev->b = PETSC_NULL;
1882:   return(0);
1883: }

1885: /* -------------------------------------------------------------------------- */
1886: /*
1887:    PCDestroy_ASA - Destroys the private context for the ASA preconditioner
1888:    that was created with PCCreate_ASA().

1890:    Input Parameter:
1891: .  pc - the preconditioner context

1893:    Application Interface Routine: PCDestroy()
1894: */
1897: static PetscErrorCode PCDestroy_ASA(PC pc)
1898: {
1899:   PC_ASA         *asa;
1900:   PC_ASA_level   *asa_lev;
1901:   PC_ASA_level   *asa_next_level;

1906:   asa = (PC_ASA*)pc->data;
1907:   asa_lev = asa->levellist;

1909:   /* Delete top level data */
1910:   PetscFree(asa->ksptype_smooth);
1911:   PetscFree(asa->pctype_smooth);
1912:   PetscFree(asa->ksptype_direct);
1913:   PetscFree(asa->pctype_direct);
1914:   PetscFree(asa->coarse_mat_type);

1916:   /* this is destroyed by the levels below  */
1917: /*   SafeMatDestroy(&(asa->A)); */
1918:   SafeVecDestroy(&(asa->invsqrtdiag));
1919:   SafeVecDestroy(&(asa->b));
1920:   SafeVecDestroy(&(asa->x));
1921:   SafeVecDestroy(&(asa->r));

1923:   if (asa->dm) {DMDestroy(asa->dm);}

1925:   /* Destroy each of the levels */
1926:   while(asa_lev) {
1927:     asa_next_level = asa_lev->next;
1928:     PCDestroyLevel_ASA(asa_lev);
1929:     asa_lev = asa_next_level;
1930:   }

1932:   PetscFree(asa);
1933:   return(0);
1934: }

1938: static PetscErrorCode PCSetFromOptions_ASA(PC pc)
1939: {
1940:   PC_ASA         *asa = (PC_ASA*)pc->data;
1941:   PetscTruth     flg;
1943:   char           type[20];


1948:   PetscOptionsHead("ASA options");
1949:   /* convergence parameters */
1950:   PetscOptionsInt("-pc_asa_nu","Number of cycles to run smoother","No manual page yet",asa->nu,&(asa->nu),&flg);
1951:   PetscOptionsInt("-pc_asa_gamma","Number of cycles to run coarse grid correction","No manual page yet",asa->gamma,&(asa->gamma),&flg);
1952:   PetscOptionsReal("-pc_asa_epsilon","Tolerance for the relaxation method","No manual page yet",asa->epsilon,&(asa->epsilon),&flg);
1953:   PetscOptionsInt("-pc_asa_mu","Number of cycles to relax in setup stages","No manual page yet",asa->mu,&(asa->mu),&flg);
1954:   PetscOptionsInt("-pc_asa_mu_initial","Number of cycles to relax for generating first candidate vector","No manual page yet",asa->mu_initial,&(asa->mu_initial),&flg);
1955:   PetscOptionsInt("-pc_asa_direct_solver","For which matrix size should we use the direct solver?","No manual page yet",asa->direct_solver,&(asa->direct_solver),&flg);
1956:   PetscOptionsTruth("-pc_asa_scale_diag","Should we scale the matrix with the inverse of its diagonal?","No manual page yet",asa->scale_diag,&(asa->scale_diag),&flg);
1957:   /* type of smoother used */
1958:   PetscOptionsList("-pc_asa_smoother_ksp_type","The type of KSP to be used in the smoothers","No manual page yet",KSPList,asa->ksptype_smooth,type,20,&flg);
1959:   if (flg) {
1960:     PetscFree(asa->ksptype_smooth);
1961:     PetscStrallocpy(type,&(asa->ksptype_smooth));
1962:   }
1963:   PetscOptionsList("-pc_asa_smoother_pc_type","The type of PC to be used in the smoothers","No manual page yet",PCList,asa->pctype_smooth,type,20,&flg);
1964:   if (flg) {
1965:     PetscFree(asa->pctype_smooth);
1966:     PetscStrallocpy(type,&(asa->pctype_smooth));
1967:   }
1968:   PetscOptionsList("-pc_asa_direct_ksp_type","The type of KSP to be used in the direct solver","No manual page yet",KSPList,asa->ksptype_direct,type,20,&flg);
1969:   if (flg) {
1970:     PetscFree(asa->ksptype_direct);
1971:     PetscStrallocpy(type,&(asa->ksptype_direct));
1972:   }
1973:   PetscOptionsList("-pc_asa_direct_pc_type","The type of PC to be used in the direct solver","No manual page yet",PCList,asa->pctype_direct,type,20,&flg);
1974:   if (flg) {
1975:     PetscFree(asa->pctype_direct);
1976:     PetscStrallocpy(type,&(asa->pctype_direct));
1977:   }
1978:   /* options specific for certain smoothers */
1979:   PetscOptionsReal("-pc_asa_richardson_scale","Scaling parameter for preconditioning in relaxation, if smoothing KSP is Richardson","No manual page yet",asa->richardson_scale,&(asa->richardson_scale),&flg);
1980:   PetscOptionsReal("-pc_asa_sor_omega","Scaling parameter for preconditioning in relaxation, if smoothing KSP is Richardson","No manual page yet",asa->sor_omega,&(asa->sor_omega),&flg);
1981:   /* options for direct solver */
1982:   PetscOptionsString("-pc_asa_coarse_mat_type","The coarse level matrix type (e.g. SuperLU, MUMPS, ...)","No manual page yet",asa->coarse_mat_type, type,20,&flg);
1983:   if (flg) {
1984:     PetscFree(asa->coarse_mat_type);
1985:     PetscStrallocpy(type,&(asa->coarse_mat_type));
1986:   }
1987:   /* storage allocation parameters */
1988:   PetscOptionsInt("-pc_asa_max_cand_vecs","Maximum number of candidate vectors","No manual page yet",asa->max_cand_vecs,&(asa->max_cand_vecs),&flg);
1989:   PetscOptionsInt("-pc_asa_max_dof_lev_2","The maximum number of degrees of freedom per node on level 2 (K in paper)","No manual page yet",asa->max_dof_lev_2,&(asa->max_dof_lev_2),&flg);
1990:   /* construction parameters */
1991:   PetscOptionsReal("-pc_asa_rq_improve","Threshold in RQ improvement for adding another candidate","No manual page yet",asa->rq_improve,&(asa->rq_improve),&flg);
1992:   PetscOptionsTail();
1993:   return(0);
1994: }

1998: static PetscErrorCode PCView_ASA(PC pc,PetscViewer viewer)
1999: {
2000:   PC_ASA          *asa = (PC_ASA*)pc->data;
2002:   PetscTruth     iascii;
2003:   PC_ASA_level   *asa_lev = asa->levellist;

2006:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
2007:   if (iascii) {
2008:     PetscViewerASCIIPrintf(viewer,"  ASA:\n");
2009:     asa_lev = asa->levellist;
2010:     while (asa_lev) {
2011:       if (!asa_lev->next) {
2012:         PetscViewerASCIIPrintf(viewer,"Coarse gride solver -- level %D -------------------------------\n",0);
2013:       } else {
2014:         PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level ? -------------------------------\n");
2015:       }
2016:       PetscViewerASCIIPushTab(viewer);
2017:       KSPView(asa_lev->smoothd,viewer);
2018:       PetscViewerASCIIPopTab(viewer);
2019:       if (asa_lev->next && asa_lev->smoothd == asa_lev->smoothu) {
2020:         PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");
2021:       } else if (asa_lev->next){
2022:         PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level ? -------------------------------\n");
2023:         PetscViewerASCIIPushTab(viewer);
2024:         KSPView(asa_lev->smoothu,viewer);
2025:         PetscViewerASCIIPopTab(viewer);
2026:       }
2027:       asa_lev = asa_lev->next;
2028:     }
2029:   } else {
2030:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCASA",((PetscObject)viewer)->type_name);
2031:   }
2032:   return(0);
2033: }

2035: /* -------------------------------------------------------------------------- */
2036: /*
2037:    PCCreate_ASA - Creates a ASA preconditioner context, PC_ASA, 
2038:    and sets this as the private data within the generic preconditioning 
2039:    context, PC, that was created within PCCreate().

2041:    Input Parameter:
2042: .  pc - the preconditioner context

2044:    Application Interface Routine: PCCreate()
2045: */
2049: PetscErrorCode  PCCreate_ASA(PC pc)
2050: {
2052:   PC_ASA         *asa;


2057:   /*
2058:       Set the pointers for the functions that are provided above.
2059:       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
2060:       are called, they will automatically call these functions.  Note we
2061:       choose not to provide a couple of these functions since they are
2062:       not needed.
2063:   */
2064:   pc->ops->apply               = PCApply_ASA;
2065:   /*  pc->ops->applytranspose      = PCApply_ASA;*/
2066:   pc->ops->applyrichardson     = PCApplyRichardson_ASA;
2067:   pc->ops->setup               = 0;
2068:   pc->ops->destroy             = PCDestroy_ASA;
2069:   pc->ops->setfromoptions      = PCSetFromOptions_ASA;
2070:   pc->ops->view                = PCView_ASA;

2072:   /* Set the data to pointer to 0 */
2073:   pc->data                = (void*)0;

2075:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASASetDM_C","PCASASetDM_ASA",PCASASetDM_ASA);
2076:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASASetTolerances_C","PCASASetTolerances_ASA",PCASASetTolerances_ASA);

2078:   /* register events */
2079:   if (! asa_events_registered) {
2084:     asa_events_registered = PETSC_TRUE;
2085:   }

2087:   /* Create new PC_ASA object */
2088:   PetscNewLog(pc,PC_ASA,&asa);
2089:   pc->data = (void*)asa;

2091:   /* WORK: find some better initial values  */
2092:   asa->nu             = 3;
2093:   asa->gamma          = 1;
2094:   asa->epsilon        = 1e-4;
2095:   asa->mu             = 3;
2096:   asa->mu_initial     = 20;
2097:   asa->direct_solver  = 100;
2098:   asa->scale_diag     = PETSC_TRUE;
2099:   PetscStrallocpy(KSPRICHARDSON, (char **) &(asa->ksptype_smooth));
2100:   PetscStrallocpy(PCSOR, (char **) &(asa->pctype_smooth));
2101:   asa->smoother_rtol    = 1e-10;
2102:   asa->smoother_abstol  = 1e-20;
2103:   asa->smoother_dtol    = PETSC_DEFAULT;
2104:   PetscStrallocpy(KSPPREONLY, (char **) &(asa->ksptype_direct));
2105:   PetscStrallocpy(PCREDUNDANT, (char **) &(asa->pctype_direct));
2106:   asa->direct_rtol      = 1e-10;
2107:   asa->direct_abstol    = 1e-20;
2108:   asa->direct_dtol      = PETSC_DEFAULT;
2109:   asa->richardson_scale = PETSC_DECIDE;
2110:   asa->sor_omega        = PETSC_DECIDE;
2111:   PetscStrallocpy(MATSAME, (char **) &(asa->coarse_mat_type));

2113:   asa->max_cand_vecs    = 4;
2114:   asa->max_dof_lev_2    = 640; /* I don't think this parameter really matters, 640 should be enough for everyone! */

2116:   asa->multigrid_constructed = PETSC_FALSE;

2118:   asa->rtol       = 1e-10;
2119:   asa->abstol     = 1e-15;
2120:   asa->divtol     = 1e5;
2121:   asa->max_it     = 10000;
2122:   asa->rq_improve = 0.9;
2123: 
2124:   asa->A           = 0;
2125:   asa->invsqrtdiag = 0;
2126:   asa->b           = 0;
2127:   asa->x           = 0;
2128:   asa->r           = 0;

2130:   asa->dm = 0;
2131: 
2132:   asa->levels    = 0;
2133:   asa->levellist = 0;

2135:   asa->comm = ((PetscObject)pc)->comm;
2136:   return(0);
2137: }