<div dir="ltr"><div class="gmail_extra"><div class="gmail_quote">On Mon, Sep 22, 2014 at 9:15 PM, Jean-Arthur Louis Olive <span dir="ltr"><<a href="mailto:jaolive@mit.edu" target="_blank">jaolive@mit.edu</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex">



<div style="word-wrap:break-word">
Hi all,
<div>I am using PETSc (dev version) to solve the Stokes + temperature equations. My DM has fields (vx, vy, p, T).</div></div></blockquote><div><br></div><div>I have finally had time to look at this. I have tried to reproduce this setup in a PETSc example. Here is SNES ex19:</div><div><br></div><div>cd src/snes/examples/tutorials</div><div>make ex19</div><div>./ex19 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_block_size 4 -pc_fieldsplit_type SCHUR -pc_fieldsplit_0_fields 0,1,2 -pc_fieldsplit_1_fields 3 -fieldsplit_1_pc_type lu -fieldsplit_0_pc_type fieldsplit -fieldsplit_0_pc_fieldsplit_block_size 3 -fieldsplit_0_pc_fieldsplit_0_fields 0,1 -fieldsplit_0_pc_fieldsplit_1_fields 2 -fieldsplit_0_pc_fieldsplit_type schur -fieldsplit_0_fieldsplit_0_pc_type lu -fieldsplit_0_fieldsplit_1_pc_type lu -snes_monitor_short -ksp_monitor_short -snes_view<br></div><div><br></div><div><div>lid velocity = 0.0625, prandtl # = 1, grashof # = 1</div><div>  0 SNES Function norm 0.239155 </div><div>    0 KSP Residual norm 0.239155 </div><div>    1 KSP Residual norm 8.25786e-07 </div><div>  1 SNES Function norm 6.82106e-05 </div><div>    0 KSP Residual norm 6.82106e-05 </div><div>    1 KSP Residual norm 1.478e-11 </div><div>  2 SNES Function norm 1.533e-11 </div><div>SNES Object: 1 MPI processes</div><div>  type: newtonls</div><div>  maximum iterations=50, maximum function evaluations=10000</div><div>  tolerances: relative=1e-08, absolute=1e-50, solution=1e-08</div><div>  total number of linear solver iterations=2</div><div>  total number of function evaluations=3</div><div>  SNESLineSearch Object:   1 MPI processes</div><div>    type: bt</div><div>      interpolation: cubic</div><div>      alpha=1.000000e-04</div><div>    maxstep=1.000000e+08, minlambda=1.000000e-12</div><div>    tolerances: relative=1.000000e-08, absolute=1.000000e-15, lambda=1.000000e-08</div><div>    maximum iterations=40</div><div>  KSP Object:   1 MPI processes</div><div>    type: fgmres</div><div>      GMRES: restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement</div><div>      GMRES: happy breakdown tolerance 1e-30</div><div>    maximum iterations=10000, initial guess is zero</div><div>    tolerances:  relative=1e-05, absolute=1e-50, divergence=10000</div><div>    right preconditioning</div><div>    using UNPRECONDITIONED norm type for convergence test</div><div>  PC Object:   1 MPI processes</div><div>    type: fieldsplit</div><div>      FieldSplit with Schur preconditioner, blocksize = 4, factorization FULL</div><div>      Preconditioner for the Schur complement formed from A11</div><div>      Split info:</div><div>      Split number 0 Fields  0, 1, 2</div><div>      Split number 1 Fields  3</div><div>      KSP solver for A00 block</div><div>        KSP Object:        (fieldsplit_0_)         1 MPI processes</div><div>          type: gmres</div><div>            GMRES: restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement</div><div>            GMRES: happy breakdown tolerance 1e-30</div><div>          maximum iterations=10000, initial guess is zero</div><div>          tolerances:  relative=1e-05, absolute=1e-50, divergence=10000</div><div>          left preconditioning</div><div>          using PRECONDITIONED norm type for convergence test</div><div>        PC Object:        (fieldsplit_0_)         1 MPI processes</div><div>          type: fieldsplit</div><div>            FieldSplit with Schur preconditioner, blocksize = 3, factorization FULL</div><div>            Preconditioner for the Schur complement formed from A11</div><div>            Split info:</div><div>            Split number 0 Fields  0, 1</div><div>            Split number 1 Fields  2</div><div>            KSP solver for A00 block</div><div>              KSP Object:              (fieldsplit_0_fieldsplit_0_)               1 MPI processes</div><div>                type: gmres</div><div>                  GMRES: restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement</div><div>                  GMRES: happy breakdown tolerance 1e-30</div><div>                maximum iterations=10000, initial guess is zero</div><div>                tolerances:  relative=1e-05, absolute=1e-50, divergence=10000</div><div>                left preconditioning</div><div>                using PRECONDITIONED norm type for convergence test</div><div>              PC Object:              (fieldsplit_0_fieldsplit_0_)               1 MPI processes</div><div>                type: lu</div><div>                  LU: out-of-place factorization</div><div>                  tolerance for zero pivot 2.22045e-14</div><div>                  matrix ordering: nd</div><div>                  factor fill ratio given 5, needed 1.875</div><div>                    Factored matrix follows:</div><div>                      Mat Object:                       1 MPI processes</div><div>                        type: seqaij</div><div>                        rows=32, cols=32</div><div>                        package used to perform factorization: petsc</div><div>                        total: nonzeros=480, allocated nonzeros=480</div><div>                        total number of mallocs used during MatSetValues calls =0</div><div>                          using I-node routines: found 12 nodes, limit used is 5</div><div>                linear system matrix = precond matrix:</div><div>                Mat Object:                (fieldsplit_0_fieldsplit_0_)                 1 MPI processes</div><div>                  type: seqaij</div><div>                  rows=32, cols=32</div><div>                  total: nonzeros=256, allocated nonzeros=256</div><div>                  total number of mallocs used during MatSetValues calls =0</div><div>                    using I-node routines: found 16 nodes, limit used is 5</div><div>            KSP solver for S = A11 - A10 inv(A00) A01 </div><div>              KSP Object:              (fieldsplit_0_fieldsplit_1_)               1 MPI processes</div><div>                type: gmres</div><div>                  GMRES: restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement</div><div>                  GMRES: happy breakdown tolerance 1e-30</div><div>                maximum iterations=10000, initial guess is zero</div><div>                tolerances:  relative=1e-05, absolute=1e-50, divergence=10000</div><div>                left preconditioning</div><div>                using PRECONDITIONED norm type for convergence test</div><div>              PC Object:              (fieldsplit_0_fieldsplit_1_)               1 MPI processes</div><div>                type: lu</div><div>                  LU: out-of-place factorization</div><div>                  tolerance for zero pivot 2.22045e-14</div><div>                  matrix ordering: nd</div><div>                  factor fill ratio given 5, needed 1.875</div><div>                    Factored matrix follows:</div><div>                      Mat Object:                       1 MPI processes</div><div>                        type: seqaij</div><div>                        rows=16, cols=16</div><div>                        package used to perform factorization: petsc</div><div>                        total: nonzeros=120, allocated nonzeros=120</div><div>                        total number of mallocs used during MatSetValues calls =0</div><div>                          using I-node routines: found 12 nodes, limit used is 5</div><div>                linear system matrix followed by preconditioner matrix:</div><div>                Mat Object:                (fieldsplit_0_fieldsplit_1_)                 1 MPI processes</div><div>                  type: schurcomplement</div><div>                  rows=16, cols=16</div><div>                    Schur complement A11 - A10 inv(A00) A01</div><div>                    A11</div><div>                      Mat Object:                      (fieldsplit_0_fieldsplit_1_)                       1 MPI processes</div><div>                        type: seqaij</div><div>                        rows=16, cols=16</div><div>                        total: nonzeros=64, allocated nonzeros=64</div><div>                        total number of mallocs used during MatSetValues calls =0</div><div>                          not using I-node routines</div><div>                    A10</div><div>                      Mat Object:                       1 MPI processes</div><div>                        type: seqaij</div><div>                        rows=16, cols=32</div><div>                        total: nonzeros=128, allocated nonzeros=128</div><div>                        total number of mallocs used during MatSetValues calls =0</div><div>                          not using I-node routines</div><div>                    KSP of A00</div><div>                      KSP Object:                      (fieldsplit_0_fieldsplit_0_)                       1 MPI processes</div><div>                        type: gmres</div><div>                          GMRES: restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement</div><div>                          GMRES: happy breakdown tolerance 1e-30</div><div>                        maximum iterations=10000, initial guess is zero</div><div>                        tolerances:  relative=1e-05, absolute=1e-50, divergence=10000</div><div>                        left preconditioning</div><div>                        using PRECONDITIONED norm type for convergence test</div><div>                      PC Object:                      (fieldsplit_0_fieldsplit_0_)                       1 MPI processes</div><div>                        type: lu</div><div>                          LU: out-of-place factorization</div><div>                          tolerance for zero pivot 2.22045e-14</div><div>                          matrix ordering: nd</div><div>                          factor fill ratio given 5, needed 1.875</div><div>                            Factored matrix follows:</div><div>                              Mat Object:                               1 MPI processes</div><div>                                type: seqaij</div><div>                                rows=32, cols=32</div><div>                                package used to perform factorization: petsc</div><div>                                total: nonzeros=480, allocated nonzeros=480</div><div>                                total number of mallocs used during MatSetValues calls =0</div><div>                                  using I-node routines: found 12 nodes, limit used is 5</div><div>                        linear system matrix = precond matrix:</div><div>                        Mat Object:                        (fieldsplit_0_fieldsplit_0_)                         1 MPI processes</div><div>                          type: seqaij</div><div>                          rows=32, cols=32</div><div>                          total: nonzeros=256, allocated nonzeros=256</div><div>                          total number of mallocs used during MatSetValues calls =0</div><div>                            using I-node routines: found 16 nodes, limit used is 5</div><div>                    A01</div><div>                      Mat Object:                       1 MPI processes</div><div>                        type: seqaij</div><div>                        rows=32, cols=16</div><div>                        total: nonzeros=128, allocated nonzeros=128</div><div>                        total number of mallocs used during MatSetValues calls =0</div><div>                          using I-node routines: found 16 nodes, limit used is 5</div><div>                Mat Object:                (fieldsplit_0_fieldsplit_1_)                 1 MPI processes</div><div>                  type: seqaij</div><div>                  rows=16, cols=16</div><div>                  total: nonzeros=64, allocated nonzeros=64</div><div>                  total number of mallocs used during MatSetValues calls =0</div><div>                    not using I-node routines</div><div>          linear system matrix = precond matrix:</div><div>          Mat Object:          (fieldsplit_0_)           1 MPI processes</div><div>            type: seqaij</div><div>            rows=48, cols=48</div><div>            total: nonzeros=576, allocated nonzeros=576</div><div>            total number of mallocs used during MatSetValues calls =0</div><div>              using I-node routines: found 16 nodes, limit used is 5</div><div>      KSP solver for S = A11 - A10 inv(A00) A01 </div><div>        KSP Object:        (fieldsplit_1_)         1 MPI processes</div><div>          type: gmres</div><div>            GMRES: restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement</div><div>            GMRES: happy breakdown tolerance 1e-30</div><div>          maximum iterations=10000, initial guess is zero</div><div>          tolerances:  relative=1e-05, absolute=1e-50, divergence=10000</div><div>          left preconditioning</div><div>          using PRECONDITIONED norm type for convergence test</div><div>        PC Object:        (fieldsplit_1_)         1 MPI processes</div><div>          type: lu</div><div>            LU: out-of-place factorization</div><div>            tolerance for zero pivot 2.22045e-14</div><div>            matrix ordering: nd</div><div>            factor fill ratio given 5, needed 1.875</div><div>              Factored matrix follows:</div><div>                Mat Object:                 1 MPI processes</div><div>                  type: seqaij</div><div>                  rows=16, cols=16</div><div>                  package used to perform factorization: petsc</div><div>                  total: nonzeros=120, allocated nonzeros=120</div><div>                  total number of mallocs used during MatSetValues calls =0</div><div>                    using I-node routines: found 12 nodes, limit used is 5</div><div>          linear system matrix followed by preconditioner matrix:</div><div>          Mat Object:          (fieldsplit_1_)           1 MPI processes</div><div>            type: schurcomplement</div><div>            rows=16, cols=16</div><div>              Schur complement A11 - A10 inv(A00) A01</div><div>              A11</div><div>                Mat Object:                (fieldsplit_1_)                 1 MPI processes</div><div>                  type: seqaij</div><div>                  rows=16, cols=16</div><div>                  total: nonzeros=64, allocated nonzeros=64</div><div>                  total number of mallocs used during MatSetValues calls =0</div><div>                    not using I-node routines</div><div>              A10</div><div>                Mat Object:                 1 MPI processes</div><div>                  type: seqaij</div><div>                  rows=16, cols=48</div><div>                  total: nonzeros=192, allocated nonzeros=192</div><div>                  total number of mallocs used during MatSetValues calls =0</div><div>                    not using I-node routines</div><div>              KSP of A00</div><div>                KSP Object:                (fieldsplit_0_)                 1 MPI processes</div><div>                  type: gmres</div><div>                    GMRES: restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement</div><div>                    GMRES: happy breakdown tolerance 1e-30</div><div>                  maximum iterations=10000, initial guess is zero</div><div>                  tolerances:  relative=1e-05, absolute=1e-50, divergence=10000</div><div>                  left preconditioning</div><div>                  using PRECONDITIONED norm type for convergence test</div><div>                PC Object:                (fieldsplit_0_)                 1 MPI processes</div><div>                  type: fieldsplit</div><div>                    FieldSplit with Schur preconditioner, blocksize = 3, factorization FULL</div><div>                    Preconditioner for the Schur complement formed from A11</div><div>                    Split info:</div><div>                    Split number 0 Fields  0, 1</div><div>                    Split number 1 Fields  2</div><div>                    KSP solver for A00 block</div><div>                      KSP Object:                      (fieldsplit_0_fieldsplit_0_)                       1 MPI processes</div><div>                        type: gmres</div><div>                          GMRES: restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement</div><div>                          GMRES: happy breakdown tolerance 1e-30</div><div>                        maximum iterations=10000, initial guess is zero</div><div>                        tolerances:  relative=1e-05, absolute=1e-50, divergence=10000</div><div>                        left preconditioning</div><div>                        using PRECONDITIONED norm type for convergence test</div><div>                      PC Object:                      (fieldsplit_0_fieldsplit_0_)                       1 MPI processes</div><div>                        type: lu</div><div>                          LU: out-of-place factorization</div><div>                          tolerance for zero pivot 2.22045e-14</div><div>                          matrix ordering: nd</div><div>                          factor fill ratio given 5, needed 1.875</div><div>                            Factored matrix follows:</div><div>                              Mat Object:                               1 MPI processes</div><div>                                type: seqaij</div><div>                                rows=32, cols=32</div><div>                                package used to perform factorization: petsc</div><div>                                total: nonzeros=480, allocated nonzeros=480</div><div>                                total number of mallocs used during MatSetValues calls =0</div><div>                                  using I-node routines: found 12 nodes, limit used is 5</div><div>                        linear system matrix = precond matrix:</div><div>                        Mat Object:                        (fieldsplit_0_fieldsplit_0_)                         1 MPI processes</div><div>                          type: seqaij</div><div>                          rows=32, cols=32</div><div>                          total: nonzeros=256, allocated nonzeros=256</div><div>                          total number of mallocs used during MatSetValues calls =0</div><div>                            using I-node routines: found 16 nodes, limit used is 5</div><div>                    KSP solver for S = A11 - A10 inv(A00) A01 </div><div>                      KSP Object:                      (fieldsplit_0_fieldsplit_1_)                       1 MPI processes</div><div>                        type: gmres</div><div>                          GMRES: restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement</div><div>                          GMRES: happy breakdown tolerance 1e-30</div><div>                        maximum iterations=10000, initial guess is zero</div><div>                        tolerances:  relative=1e-05, absolute=1e-50, divergence=10000</div><div>                        left preconditioning</div><div>                        using PRECONDITIONED norm type for convergence test</div><div>                      PC Object:                      (fieldsplit_0_fieldsplit_1_)                       1 MPI processes</div><div>                        type: lu</div><div>                          LU: out-of-place factorization</div><div>                          tolerance for zero pivot 2.22045e-14</div><div>                          matrix ordering: nd</div><div>                          factor fill ratio given 5, needed 1.875</div><div>                            Factored matrix follows:</div><div>                              Mat Object:                               1 MPI processes</div><div>                                type: seqaij</div><div>                                rows=16, cols=16</div><div>                                package used to perform factorization: petsc</div><div>                                total: nonzeros=120, allocated nonzeros=120</div><div>                                total number of mallocs used during MatSetValues calls =0</div><div>                                  using I-node routines: found 12 nodes, limit used is 5</div><div>                        linear system matrix followed by preconditioner matrix:</div><div>                        Mat Object:                        (fieldsplit_0_fieldsplit_1_)                         1 MPI processes</div><div>                          type: schurcomplement</div><div>                          rows=16, cols=16</div><div>                            Schur complement A11 - A10 inv(A00) A01</div><div>                            A11</div><div>                              Mat Object:                              (fieldsplit_0_fieldsplit_1_)                               1 MPI processes</div><div>                                type: seqaij</div><div>                                rows=16, cols=16</div><div>                                total: nonzeros=64, allocated nonzeros=64</div><div>                                total number of mallocs used during MatSetValues calls =0</div><div>                                  not using I-node routines</div><div>                            A10</div><div>                              Mat Object:                               1 MPI processes</div><div>                                type: seqaij</div><div>                                rows=16, cols=32</div><div>                                total: nonzeros=128, allocated nonzeros=128</div><div>                                total number of mallocs used during MatSetValues calls =0</div><div>                                  not using I-node routines</div><div>                            KSP of A00</div><div>                              KSP Object:                              (fieldsplit_0_fieldsplit_0_)                               1 MPI processes</div><div>                                type: gmres</div><div>                                  GMRES: restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement</div><div>                                  GMRES: happy breakdown tolerance 1e-30</div><div>                                maximum iterations=10000, initial guess is zero</div><div>                                tolerances:  relative=1e-05, absolute=1e-50, divergence=10000</div><div>                                left preconditioning</div><div>                                using PRECONDITIONED norm type for convergence test</div><div>                              PC Object:                              (fieldsplit_0_fieldsplit_0_)                               1 MPI processes</div><div>                                type: lu</div><div>                                  LU: out-of-place factorization</div><div>                                  tolerance for zero pivot 2.22045e-14</div><div>                                  matrix ordering: nd</div><div>                                  factor fill ratio given 5, needed 1.875</div><div>                                    Factored matrix follows:</div><div>                                      Mat Object:                                       1 MPI processes</div><div>                                        type: seqaij</div><div>                                        rows=32, cols=32</div><div>                                        package used to perform factorization: petsc</div><div>                                        total: nonzeros=480, allocated nonzeros=480</div><div>                                        total number of mallocs used during MatSetValues calls =0</div><div>                                          using I-node routines: found 12 nodes, limit used is 5</div><div>                                linear system matrix = precond matrix:</div><div>                                Mat Object:                                (fieldsplit_0_fieldsplit_0_)                                 1 MPI processes</div><div>                                  type: seqaij</div><div>                                  rows=32, cols=32</div><div>                                  total: nonzeros=256, allocated nonzeros=256</div><div>                                  total number of mallocs used during MatSetValues calls =0</div><div>                                    using I-node routines: found 16 nodes, limit used is 5</div><div>                            A01</div><div>                              Mat Object:                               1 MPI processes</div><div>                                type: seqaij</div><div>                                rows=32, cols=16</div><div>                                total: nonzeros=128, allocated nonzeros=128</div><div>                                total number of mallocs used during MatSetValues calls =0</div><div>                                  using I-node routines: found 16 nodes, limit used is 5</div><div>                        Mat Object:                        (fieldsplit_0_fieldsplit_1_)                         1 MPI processes</div><div>                          type: seqaij</div><div>                          rows=16, cols=16</div><div>                          total: nonzeros=64, allocated nonzeros=64</div><div>                          total number of mallocs used during MatSetValues calls =0</div><div>                            not using I-node routines</div><div>                  linear system matrix = precond matrix:</div><div>                  Mat Object:                  (fieldsplit_0_)                   1 MPI processes</div><div>                    type: seqaij</div><div>                    rows=48, cols=48</div><div>                    total: nonzeros=576, allocated nonzeros=576</div><div>                    total number of mallocs used during MatSetValues calls =0</div><div>                      using I-node routines: found 16 nodes, limit used is 5</div><div>              A01</div><div>                Mat Object:                 1 MPI processes</div><div>                  type: seqaij</div><div>                  rows=48, cols=16</div><div>                  total: nonzeros=192, allocated nonzeros=192</div><div>                  total number of mallocs used during MatSetValues calls =0</div><div>                    using I-node routines: found 16 nodes, limit used is 5</div><div>          Mat Object:          (fieldsplit_1_)           1 MPI processes</div><div>            type: seqaij</div><div>            rows=16, cols=16</div><div>            total: nonzeros=64, allocated nonzeros=64</div><div>            total number of mallocs used during MatSetValues calls =0</div><div>              not using I-node routines</div><div>    linear system matrix = precond matrix:</div><div>    Mat Object:     1 MPI processes</div><div>      type: seqaij</div><div>      rows=64, cols=64, bs=4</div><div>      total: nonzeros=1024, allocated nonzeros=1024</div><div>      total number of mallocs used during MatSetValues calls =0</div><div>        using I-node routines: found 16 nodes, limit used is 5</div><div>Number of SNES iterations = 2</div></div><div><br></div><div>I cannot replicate your failure here. I am running on the 'next' branch here. What are you using?</div><div>Also, are you using the DMDA for data layout?</div><div><br></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex"><div style="word-wrap:break-word">
<div>I would like to use nested fieldsplits to separate the T part from the Stokes part, and apply a Schur complement approach to the Stokes block.</div>
<div>Unfortunately, I keep getting this error message:</div>
<div><span style="color:rgb(41,249,20);font-family:'Andale Mono';background-color:rgb(0,0,0)">[1]PETSC ERROR: DMCreateFieldDecomposition() line 1274 in /home/jolive/petsc/src/dm/interface/dm.c Decomposition defined only after DMSetUp</span></div>
<div><br>
</div>
<div>Here are the command line options I tried:</div>
<div><br>
</div>
<div>
<div><span style="font-family:Menlo;font-size:11px">-snes_type ksponly \</span></div>
<div>
<div style="margin:0px;font-size:11px;font-family:Menlo">-ksp_type fgmres \</div>
<div style="margin:0px;font-size:11px;font-family:Menlo"><span style="color:rgb(0,132,0)"># define 2 fields: [vx vy p] and [T] </span></div>
<div style="margin:0px;font-size:11px;font-family:Menlo">-pc_type fieldsplit -pc_fieldsplit_<span style="color:rgb(39,42,216)">0</span>_fields <span style="color:rgb(39,42,216)">0</span>,<span style="color:rgb(39,42,216)">1</span>,<span style="color:rgb(39,42,216)">2</span> -pc_fieldsplit_<span style="color:rgb(39,42,216)">1</span>_fields <span style="color:rgb(39,42,216)">3</span> \</div>
<div style="margin:0px;font-size:11px;font-family:Menlo"><span style="color:rgb(0,132,0)"># split [vx vy p] into 2 fields: [vx vy] and [p] </span></div>
<div style="margin:0px;font-size:11px;font-family:Menlo">-fieldsplit_<span style="color:rgb(39,42,216)">0</span>_pc_type fieldsplit \</div>
<div style="margin:0px;font-size:11px;font-family:Menlo">
<div style="margin:0px">-pc_fieldsplit_<span style="color:rgb(39,42,216)">0</span>_fieldsplit_<span style="color:rgb(39,42,216)">0</span>_fields <span style="color:rgb(39,42,216)">0</span>,<span style="color:rgb(39,42,216)">1</span> -pc_fieldsplit_<span style="color:rgb(39,42,216)">0</span>_fieldsplit_<span style="color:rgb(39,42,216)">1</span>_fields <span style="color:rgb(39,42,216)">2</span> \</div></div></div></div></div></blockquote><div><br></div><div>Note that the 2 options above are wrong. It should be -fieldsplit_0_pc_fieldsplit_0_fields 0,1</div><div><br></div><div>  Thanks,</div><div><br></div><div>     Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex"><div style="word-wrap:break-word"><div><div><div style="margin:0px;font-size:11px;font-family:Menlo">
<div style="margin:0px"><span style="color:rgb(0,132,0)"># apply schur complement to [vx vy p]</span></div>
</div>
<div style="margin:0px;font-size:11px;font-family:Menlo">-fieldsplit_<span style="color:rgb(39,42,216)">0</span>_pc_fieldsplit_type schur \</div>
<div style="margin:0px;font-size:11px;font-family:Menlo">-fieldsplit_<span style="color:rgb(39,42,216)">0</span>_pc_fieldsplit_schur_factorization_type upper \</div>
<div style="margin:0px;font-size:11px;font-family:Menlo"><br>
</div>
<div style="margin:0px;font-size:11px;font-family:Menlo"><span style="color:rgb(0,132,0)"># solve everything with lu, just for testing</span></div>
<div style="margin:0px;font-size:11px;font-family:Menlo">-fieldsplit_<span style="color:rgb(39,42,216)">0</span>_fieldsplit_<span style="color:rgb(39,42,216)">0</span>_ksp_type preonly \</div>
<div style="margin:0px;font-size:11px;font-family:Menlo">-fieldsplit_<span style="color:rgb(39,42,216)">0</span>_fieldsplit_<span style="color:rgb(39,42,216)">0</span>_pc_type lu -fieldsplit_<span style="color:rgb(39,42,216)">0</span>_fieldsplit_<span style="color:rgb(39,42,216)">0</span>_pc_factor_mat_solver_package
 superlu_dist \</div>
<div style="margin:0px;font-size:11px;font-family:Menlo">-fieldsplit_<span style="color:rgb(39,42,216)">0</span>_fieldsplit_<span style="color:rgb(39,42,216)">1</span>_ksp_type preonly \</div>
<div style="margin:0px;font-size:11px;font-family:Menlo">-fieldsplit_<span style="color:rgb(39,42,216)">0</span>_fieldsplit_<span style="color:rgb(39,42,216)">1</span>_pc_type lu -fieldsplit_<span style="color:rgb(39,42,216)">0</span>_fieldsplit_<span style="color:rgb(39,42,216)">1</span>_pc_factor_mat_solver_package
 superlu_dist \</div>
<div style="margin:0px;font-size:11px;font-family:Menlo">-fieldsplit_<span style="color:rgb(39,42,216)">1</span>_ksp_type preonly \</div>
<div style="margin:0px;font-size:11px;font-family:Menlo">-fieldsplit_<span style="color:rgb(39,42,216)">1</span>_pc_type lu -fieldsplit_<span style="color:rgb(39,42,216)">1</span>_pc_factor_mat_solver_package superlu_dist \</div>
</div>
</div>
<div style="margin:0px;font-size:11px;font-family:Menlo"><br>
</div>
<div style="margin:0px;font-size:11px;font-family:Menlo"><span style="font-family:Helvetica;font-size:12px">Any idea what could be causing this?</span></div>
<div style="margin:0px;font-size:11px;font-family:Menlo"><span style="font-family:Helvetica;font-size:12px">Thanks a lot,</span></div>
<div style="margin:0px;font-size:11px;font-family:Menlo"><span style="font-family:Helvetica;font-size:12px">Arthur</span></div>
</div>

</blockquote></div><br><br clear="all"><div><br></div>-- <br>What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br>-- Norbert Wiener
</div></div>