<html>
  <head>
    <meta content="text/html; charset=ISO-8859-1"
      http-equiv="Content-Type">
  </head>
  <body text="#000000" bgcolor="#FFFFFF">
    <br>
    <div class="moz-cite-prefix">On 07/17/2012 11:03 AM, Hong Zhang
      wrote:<br>
    </div>
    <blockquote
cite="mid:CAGCphBv-S3yMWrr30XoN7As+G1Ry4ohTJzPmtr8+StP-fDiJww@mail.gmail.com"
      type="cite">Michele :<br>
      <div class="gmail_quote">
        <blockquote class="gmail_quote" style="margin:0 0 0
          .8ex;border-left:1px #ccc solid;padding-left:1ex">
          <div text="#000000" bgcolor="#FFFFFF"><font face="Ubuntu"><br>
              I have some problems with the block jacobi preconditioner.<br>
              I am solving a 3D  Poisson equation with periodic BCs,
              discretized by using finite differences (7-points
              stencil).<br>
              Thus the problem is singular and the nullspace has to be
              removed.<br>
            </font></div>
        </blockquote>
        <div><br>
        </div>
        <div>For  <span class="Apple-style-span"
            style="font-family:Ubuntu">Poisson equations, multigrid
            precondition should be the method of</span></div>
        <div><span class="Apple-style-span" style="font-family:Ubuntu">choice.
            <br>
          </span></div>
      </div>
    </blockquote>
    Thank you for the suggestion. I do not have any experience with
    multigrid, but I will try.<br>
    <blockquote
cite="mid:CAGCphBv-S3yMWrr30XoN7As+G1Ry4ohTJzPmtr8+StP-fDiJww@mail.gmail.com"
      type="cite">
      <div class="gmail_quote">
        <div><span class="Apple-style-span" style="font-family:Ubuntu"><br>
          </span></div>
        <blockquote class="gmail_quote" style="margin:0 0 0
          .8ex;border-left:1px #ccc solid;padding-left:1ex">
          <div text="#000000" bgcolor="#FFFFFF"><font face="Ubuntu"> If
              I solve with the PCG method + JACOBI preconditioner the
              results are fine.<br>
              If I use PCG + Block Jacobi preconditioner + ICC on each
              block the results are fine on the majority of the
              processors,<br>
              but on few of them the error is very large.  <br>
            </font></div>
        </blockquote>
        <div> </div>
        <div>How do you know "<span class="Apple-style-span"
            style="font-family:Ubuntu"> few of them"? <br>
          </span></div>
      </div>
    </blockquote>
     Basically the solution is not correct on some grid points, say 6
    grid nodes out of 64^3. The 6 grid nodes with problems belongs to 2
    of the 128 processors<br>
    I am using.<br>
    <blockquote
cite="mid:CAGCphBv-S3yMWrr30XoN7As+G1Ry4ohTJzPmtr8+StP-fDiJww@mail.gmail.com"
      type="cite">
      <div class="gmail_quote">
        <div><font class="Apple-style-span" face="Ubuntu"><br>
          </font></div>
        <blockquote class="gmail_quote" style="margin:0 0 0
          .8ex;border-left:1px #ccc solid;padding-left:1ex">
          <div text="#000000" bgcolor="#FFFFFF"><font face="Ubuntu"> Do
              you have any idea/suggestions on how to fix this problem?
              <br>
              This is the fragment of code I am using ( petsc 3.1 and
              Fortran 90):<br>
            </font></div>
        </blockquote>
        <div> </div>
        <div>Please update to petsc-3.3. petsc-3.1 is too old.</div>
      </div>
    </blockquote>
        I would do that but the version installed on the platform
    (Intrepid at ALCF) I am working on is 3.1-p2.<br>
    <br>
    <blockquote
cite="mid:CAGCphBv-S3yMWrr30XoN7As+G1Ry4ohTJzPmtr8+StP-fDiJww@mail.gmail.com"
      type="cite">
      <div class="gmail_quote">
        <blockquote class="gmail_quote" style="margin:0 0 0
          .8ex;border-left:1px #ccc solid;padding-left:1ex">
          <div text="#000000" bgcolor="#FFFFFF">
            <font face="Ubuntu"> <br>
                  PetscErrorCode                petsc_err     <br>
                  Mat                                         A<br>
                  PC                                           pc, subpc<br>
                  KSP                                         ksp<br>
                  KSP                                         subksp(1)<br>
                  :<br>
                  :<br>
                  :<br>
                  call KSPCreate(PETSC_COMM_WORLD,ksp,petsc_err)<br>
                  call
              KSPSetOperators(ksp,A,A,DIFFERENT_NONZERO_PATTERN,petsc_err)<br>
            </font></div>
        </blockquote>
        <div> </div>
        <div>call KSPSetType(ksp,KSPCG, ) !the default type is gmres. I
          guess you want CG</div>
        <div><br>
        </div>
        <blockquote class="gmail_quote" style="margin:0 0 0
          .8ex;border-left:1px #ccc solid;padding-left:1ex">
          <div text="#000000" bgcolor="#FFFFFF"><font face="Ubuntu">    
              call KSPGetPC(ksp,pc,petsc_err)<br>
                  call PCSetType(pc,PCBJACOBI,petsc_err)<br>
              !    call KSPSetUp(ksp,petsc_err)    call this at the end<br>
                  <br>
                  ! KSP context for each single block<br>
                  call
              PCBJacobiGetSubKSP(pc,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,subksp(1),petsc_err)

              <br>
                  call KSPGetPC(subksp(1),subpc,petsc_err)<br>
                  call PCSetType(subpc,PCICC,petsc_err)<br>
            </font></div>
        </blockquote>
        <div> </div>
        <blockquote class="gmail_quote" style="margin:0 0 0
          .8ex;border-left:1px #ccc solid;padding-left:1ex">
          <div text="#000000" bgcolor="#FFFFFF">
            <font face="Ubuntu">     call
              KSPSetType(subksp(1),KSPPREONLY petsc_err)<br>
            </font></div>
        </blockquote>
        <div> </div>
        <blockquote class="gmail_quote" style="margin:0 0 0
          .8ex;border-left:1px #ccc solid;padding-left:1ex">
          <div text="#000000" bgcolor="#FFFFFF">
            <font face="Ubuntu">     call KSPSetTolerances(subksp(1),tol
              ,PETSC_DEFAULT_DOUBLE_PRECISION,&<br>
                       &
              PETSC_DEFAULT_DOUBLE_PRECISION,PETSC_DEFAULT_INTEGER,petsc_err)<br>
              <br>
                 ! Remove nullspace from the singular system (Check
              PETSC_NULL)<br>
                  call
MatNullSpaceCreate(MPI_COMM_WORLD,PETSC_TRUE,0,PETSC_NULL,nullspace,petsc_err)<br>
                  call KSPSetNullSpace(ksp, nullspace, petsc_err)<br>
                  call MatNullSpaceRemove(nullspace, b,
              PETSC_NULL,petsc_err)    <br>
              <br>
                  call KSPSolve(ksp,b,x,petsc_err)<br>
            </font></div>
        </blockquote>
        <div><br>
        </div>
        <div>I modified your code slightly. All these options can be
          provided at runtime:</div>
        <div>'-ksp_type cg -pc_type bjacobi -sub_pc_type icc'</div>
        <div><br>
        </div>
        <div>Hong</div>
        <blockquote class="gmail_quote" style="margin:0 0 0
          .8ex;border-left:1px #ccc solid;padding-left:1ex">
          <div text="#000000" bgcolor="#FFFFFF"><font face="Ubuntu"> <br>
              <br>
              <br>
              <br>
              <br>
              <br>
            </font>
            <div class="im"><font face="Ubuntu"><br>
                <br>
                <br>
              </font>
              <div>On 07/13/2012 12:14 PM, Hong Zhang wrote:<br>
              </div>
            </div>
            <div>
              <div class="h5">
                <blockquote type="cite">Michele :
                  <div class="gmail_quote">
                    <blockquote class="gmail_quote" style="margin:0 0 0
                      .8ex;border-left:1px #ccc solid;padding-left:1ex">
                      <div text="#000000" bgcolor="#FFFFFF"><font
                          face="Ubuntu"> <br>
                          I need to use the ICC factorization as
                          preconditioner, but I noticed that no parallel
                          version is supported.<br>
                          Is that correct?<br>
                        </font></div>
                    </blockquote>
                    <div>Correct.</div>
                    <div> </div>
                    <blockquote class="gmail_quote" style="margin:0 0 0
                      .8ex;border-left:1px #ccc solid;padding-left:1ex">
                      <div text="#000000" bgcolor="#FFFFFF"> <font
                          face="Ubuntu"> If so, is there a work around,
                          like building  the preconditioner "by hand" by
                          using PETSc functions?<br>
                        </font></div>
                    </blockquote>
                    <div>You may try block jacobi with icc in the blocks
                       '-ksp_type cg -pc_type bjacobi -sub_pc_type icc'</div>
                    <div><br>
                    </div>
                    <div>Hong </div>
                    <blockquote class="gmail_quote" style="margin:0 0 0
                      .8ex;border-left:1px #ccc solid;padding-left:1ex">
                      <div text="#000000" bgcolor="#FFFFFF"> </div>
                    </blockquote>
                  </div>
                  <br>
                </blockquote>
                <br>
                <br>
              </div>
            </div>
          </div>
        </blockquote>
      </div>
      <br>
    </blockquote>
    <br>
    <br>
  </body>
</html>