[petsc-users] odd SNES behavior

Mark F. Adams mark.adams at columbia.edu
Wed Mar 6 14:03:02 CST 2013


Good catch Peter,

I seem to be able to fix this by "redoing" the PC each solve (see below).  I see this when '-pc_gamg_reuse_interpolation false':

MatMult             5552 1.0 2.4250e+02 1.0 1.16e+09 1.0 2.9e+04 1.8e+03 7.5e+03 95 27 66 63 55  95 27 66 63 55    19
MatMultAdd           501 1.0 1.8848e-01 1.2 9.73e+07 1.0 3.2e+03 8.9e+02 0.0e+00  0  2  7  3  0   0  2  7  3  0  2055

and this with '-pc_gamg_reuse_interpolation true':

MatMult MF          1940 1.0 2.4999e+02 1.0 5.68e+08 1.0 0.0e+00 0.0e+00 7.7e+03 96 12  0  0 81  96 12  0  0 81     9
MatMult             6038 1.0 2.5118e+02 1.0 1.33e+09 1.0 3.3e+04 1.9e+03 7.7e+03 96 28 79 79 81  96 28 79 79 81    21

The offending GAMG code is appended.  It looks like the code in the first 'else' clause is causing the problem.  Do you see what is causing this problem?  I am marching up the grids redoing the Galerkin coarse grids.

Mark

redo_mesh_setup = (PetscBool)(!pc_gamg->reuse_prol);

  if (pc_gamg->setup_count++ > 0) {
    if (redo_mesh_setup) {
      /* reset everything */
      ierr            = PCReset_MG(pc);CHKERRQ(ierr);
      pc->setupcalled = 0;
    } else {
      PC_MG_Levels **mglevels = mg->levels;
      /* just do Galerkin grids */
      Mat          B,dA,dB;

     if (!pc->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCSetUp() has not been called yet");
      if (pc_gamg->Nlevels > 1) {
        /* currently only handle case where mat and pmat are the same on coarser levels */
        ierr = KSPGetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,&dA,&dB,NULL);CHKERRQ(ierr);
        /* (re)set to get dirty flag */
        ierr = KSPSetOperators(mglevels[pc_gamg->Nlevels-1]->smoothd,dA,dB,SAME_NONZERO_PATTERN);CHKERRQ(ierr);

        for (level=pc_gamg->Nlevels-2; level>=0; level--) {
          /* the first time through the matrix structure has changed from repartitioning */
          if (pc_gamg->setup_count==2 && (pc_gamg->repart || level==0)) {
            ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
            ierr = MatDestroy(&mglevels[level]->A);CHKERRQ(ierr);

            mglevels[level]->A = B;
          } else {
            ierr = KSPGetOperators(mglevels[level]->smoothd,NULL,&B,NULL);CHKERRQ(ierr);
            ierr = MatPtAP(dB,mglevels[level+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
          }
          ierr = KSPSetOperators(mglevels[level]->smoothd,B,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
          dB   = B;
        }
      }

      ierr = PCSetUp_MG(pc);CHKERRQ(ierr);

      /* PCSetUp_MG seems to insists on setting this to GMRES */
      ierr = KSPSetType(mglevels[0]->smoothd, KSPPREONLY);CHKERRQ(ierr);
      PetscFunctionReturn(0);
    }
  }



More information about the petsc-users mailing list