Commits

Jed Brown  committed 270f95d

Handle naming with missing prefixes for VecView_Nest and MatView_Nest

Hg-commit: 2cbe4b66b29dadccdf9e4ac68bbc3c5ff41ea175

  • Participants
  • Parent commits 2a4a745

Comments (0)

Files changed (4)

File src/ksp/ksp/examples/tests/output/ex11_1.out

     0 KSP Residual norm 18874.5 
     Linear solve converged due to CONVERGED_RTOL iterations 30
     Linear solve converged due to CONVERGED_RTOL iterations 30
-    1 KSP Residual norm 5980.61 
+    1 KSP Residual norm 5980.6 
     Linear solve converged due to CONVERGED_RTOL iterations 34
     Linear solve converged due to CONVERGED_RTOL iterations 30
-    2 KSP Residual norm 4874.27 
+    2 KSP Residual norm 4874.25 
     Linear solve converged due to CONVERGED_RTOL iterations 33
     Linear solve converged due to CONVERGED_RTOL iterations 30
-    3 KSP Residual norm 735.443 
+    3 KSP Residual norm 735.445 
     Linear solve converged due to CONVERGED_RTOL iterations 30
     Linear solve converged due to CONVERGED_RTOL iterations 30
     4 KSP Residual norm 186.933 
     5 KSP Residual norm 49.2586 
     Linear solve converged due to CONVERGED_RTOL iterations 33
     Linear solve converged due to CONVERGED_RTOL iterations 30
-    6 KSP Residual norm 4.93617 
+    6 KSP Residual norm 4.93616 
     Linear solve converged due to CONVERGED_RTOL iterations 33
     Linear solve converged due to CONVERGED_RTOL iterations 30
-    7 KSP Residual norm 0.721076 
+    7 KSP Residual norm 0.721077 
     Linear solve converged due to CONVERGED_RTOL iterations 31
     Linear solve converged due to CONVERGED_RTOL iterations 30
     8 KSP Residual norm 0.104728 
-  3 KSP Residual norm 2.072e-10 
+  3 KSP Residual norm 2.073e-10 
+    Linear solve converged due to CONVERGED_RTOL iterations 34
     Linear solve converged due to CONVERGED_RTOL iterations 34
-    Linear solve converged due to CONVERGED_RTOL iterations 29
     Residual norms for fc_fieldsplit_pressure_ solve.
-    0 KSP Residual norm 5061.58 
-    Linear solve converged due to CONVERGED_RTOL iterations 32
+    0 KSP Residual norm 5052.31 
+    Linear solve converged due to CONVERGED_RTOL iterations 34
     Linear solve converged due to CONVERGED_RTOL iterations 33
-    1 KSP Residual norm 3963.48 
+    1 KSP Residual norm 3867.45 
     Linear solve converged due to CONVERGED_RTOL iterations 30
     Linear solve converged due to CONVERGED_RTOL iterations 30
-    2 KSP Residual norm 1285.3 
+    2 KSP Residual norm 1289.57 
     Linear solve converged due to CONVERGED_RTOL iterations 33
     Linear solve converged due to CONVERGED_RTOL iterations 30
-    3 KSP Residual norm 744.933 
+    3 KSP Residual norm 754.745 
     Linear solve converged due to CONVERGED_RTOL iterations 33
     Linear solve converged due to CONVERGED_RTOL iterations 30
-    4 KSP Residual norm 109.99 
+    4 KSP Residual norm 107.047 
     Linear solve converged due to CONVERGED_RTOL iterations 35
     Linear solve converged due to CONVERGED_RTOL iterations 30
-    5 KSP Residual norm 22.2872 
+    5 KSP Residual norm 22.3026 
     Linear solve converged due to CONVERGED_RTOL iterations 33
     Linear solve converged due to CONVERGED_RTOL iterations 30
-    6 KSP Residual norm 1.70277 
+    6 KSP Residual norm 1.71227 
     Linear solve converged due to CONVERGED_RTOL iterations 34
     Linear solve converged due to CONVERGED_RTOL iterations 30
-    7 KSP Residual norm 0.21094 
+    7 KSP Residual norm 0.20216 
     Linear solve converged due to CONVERGED_RTOL iterations 34
     Linear solve converged due to CONVERGED_RTOL iterations 30
-    8 KSP Residual norm 0.0373685 
+    8 KSP Residual norm 0.0518672 
+    Linear solve converged due to CONVERGED_RTOL iterations 31
+    Linear solve converged due to CONVERGED_RTOL iterations 30
+    9 KSP Residual norm 0.0102098 
   4 KSP Residual norm < 1.e-11
 KSP Object:(fc_) 1 MPI processes
   type: fgmres
                 total number of mallocs used during MatSetValues calls =0
                     block size is 1
         linear system matrix = precond matrix:
-        Matrix Object:         1 MPI processes        
+        Matrix Object:        (a11_)         1 MPI processes        
           type: seqaij
           rows=2046, cols=2046
           total: nonzeros=34952, allocated nonzeros=34952
           rows=1024, cols=1024
             Schur complement A11 - A10 inv(A00) A01
             A11
-              Matrix Object:               1 MPI processes              
+              Matrix Object:              (a22_)               1 MPI processes              
                 type: seqaij
                 rows=1024, cols=1024
                 total: nonzeros=1024, allocated nonzeros=1024
                         total number of mallocs used during MatSetValues calls =0
                             block size is 1
                 linear system matrix = precond matrix:
-                Matrix Object:                 1 MPI processes                
+                Matrix Object:                (a11_)                 1 MPI processes                
                   type: seqaij
                   rows=2046, cols=2046
                   total: nonzeros=34952, allocated nonzeros=34952
                 total: nonzeros=7936, allocated nonzeros=7936
                 total number of mallocs used during MatSetValues calls =0
                   using I-node routines: found 1085 nodes, limit used is 5
-        Matrix Object:         1 MPI processes        
+        Matrix Object:        (a22_)         1 MPI processes        
           type: seqaij
           rows=1024, cols=1024
           total: nonzeros=1024, allocated nonzeros=1024
       Matrix object: 
         type=nest, rows=2, cols=2 
         MatNest structure: 
-        (0,0) : name="(null)", type=seqaij, rows=2046, cols=2046 
-        (0,1) : name="(null)", type=seqaij, rows=2046, cols=1024 
-        (1,0) : name="(null)", type=seqaij, rows=1024, cols=2046 
-        (1,1) : name="(null)", type=seqaij, rows=1024, cols=1024 
+        (0,0) : prefix="a11_", type=seqaij, rows=2046, cols=2046 
+        (0,1) : type=seqaij, rows=2046, cols=1024 
+        (1,0) : type=seqaij, rows=1024, cols=2046 
+        (1,1) : prefix="a22_", type=seqaij, rows=1024, cols=1024 
 -- vector vector values --
   Min(u)  = -0.017280 [loc=2035]
   Max(u)  = 0.018577 [loc=2024]
   Min(p)  = -6.039639 [loc=844]
   Max(p)  = 4.836615 [loc=16]
   Norm(p) = 94.125988 
-  Sum(p)  = -0.000000 
+  Sum(p)  = 0.000000 
 -- Full vector values --
   Min(u,p)  = -6.039639 [loc=2890]
   Max(u,p)  = 4.836615 [loc=2062]

File src/ksp/ksp/examples/tests/output/ex11_2.out

     1 KSP Residual norm 0.542365 
     Linear solve converged due to CONVERGED_RTOL iterations 23
     Linear solve converged due to CONVERGED_RTOL iterations 15
-    2 KSP Residual norm 0.219772 
+    2 KSP Residual norm 0.219773 
     Linear solve converged due to CONVERGED_RTOL iterations 28
     Linear solve converged due to CONVERGED_RTOL iterations 18
     3 KSP Residual norm 0.0566767 
     4 KSP Residual norm 0.0143373 
     Linear solve converged due to CONVERGED_RTOL iterations 27
     Linear solve converged due to CONVERGED_RTOL iterations 12
-    5 KSP Residual norm 0.0021517 
+    5 KSP Residual norm 0.00215169 
     Linear solve converged due to CONVERGED_RTOL iterations 21
     Linear solve converged due to CONVERGED_RTOL iterations 9
     6 KSP Residual norm 0.000294117 
     Linear solve converged due to CONVERGED_RTOL iterations 21
     Linear solve converged due to CONVERGED_RTOL iterations 10
     8 KSP Residual norm 3.68237e-06 
-  2 KSP Residual norm 8.15484e-07 
+  2 KSP Residual norm 9.083e-07 
     Residual norms for fc_fieldsplit_pressure_ solve.
-    0 KSP Residual norm 0.209052 
-    Linear solve converged due to CONVERGED_RTOL iterations 14
-    Linear solve converged due to CONVERGED_RTOL iterations 6
-    1 KSP Residual norm 0.0334117 
-    Linear solve converged due to CONVERGED_RTOL iterations 12
-    Linear solve converged due to CONVERGED_RTOL iterations 4
-    2 KSP Residual norm 0.0185075 
-    Linear solve converged due to CONVERGED_RTOL iterations 12
-    Linear solve converged due to CONVERGED_RTOL iterations 4
-    3 KSP Residual norm 0.0116486 
-    Linear solve converged due to CONVERGED_RTOL iterations 16
-    Linear solve converged due to CONVERGED_RTOL iterations 12
-    4 KSP Residual norm 0.00823444 
+    0 KSP Residual norm 0.188118 
+    Linear solve converged due to CONVERGED_RTOL iterations 18
+    Linear solve converged due to CONVERGED_RTOL iterations 7
+    1 KSP Residual norm 0.0711197 
+    Linear solve converged due to CONVERGED_RTOL iterations 20
     Linear solve converged due to CONVERGED_RTOL iterations 12
-    Linear solve converged due to CONVERGED_RTOL iterations 4
-    5 KSP Residual norm 0.00549928 
+    2 KSP Residual norm 0.0413885 
+    Linear solve converged due to CONVERGED_RTOL iterations 18
     Linear solve converged due to CONVERGED_RTOL iterations 10
-    Linear solve converged due to CONVERGED_RTOL iterations 4
-    6 KSP Residual norm 0.00232252 
+    3 KSP Residual norm 0.025695 
+    Linear solve converged due to CONVERGED_RTOL iterations 17
     Linear solve converged due to CONVERGED_RTOL iterations 10
+    4 KSP Residual norm 0.0149315 
+    Linear solve converged due to CONVERGED_RTOL iterations 18
+    Linear solve converged due to CONVERGED_RTOL iterations 7
+    5 KSP Residual norm 0.00778653 
+    Linear solve converged due to CONVERGED_RTOL iterations 16
     Linear solve converged due to CONVERGED_RTOL iterations 4
-    7 KSP Residual norm 0.00109302 
+    6 KSP Residual norm 0.00501512 
     Linear solve converged due to CONVERGED_RTOL iterations 9
     Linear solve converged due to CONVERGED_RTOL iterations 4
-    8 KSP Residual norm 0.000457367 
-    Linear solve converged due to CONVERGED_RTOL iterations 12
+    7 KSP Residual norm 0.00392468 
     Linear solve converged due to CONVERGED_RTOL iterations 4
-    9 KSP Residual norm 0.000393394 
-    Linear solve converged due to CONVERGED_RTOL iterations 9
     Linear solve converged due to CONVERGED_RTOL iterations 4
-   10 KSP Residual norm 0.000164953 
+    8 KSP Residual norm 0.00287529 
+    Linear solve converged due to CONVERGED_RTOL iterations 4
+    Linear solve converged due to CONVERGED_RTOL iterations 4
+    9 KSP Residual norm 0.00144937 
     Linear solve converged due to CONVERGED_RTOL iterations 9
+    Linear solve converged due to CONVERGED_RTOL iterations 6
+   10 KSP Residual norm 0.000304755 
+    Linear solve converged due to CONVERGED_RTOL iterations 5
     Linear solve converged due to CONVERGED_RTOL iterations 4
-   11 KSP Residual norm 8.01212e-05 
+   11 KSP Residual norm 8.25974e-05 
     Linear solve converged due to CONVERGED_RTOL iterations 4
+    Linear solve converged due to CONVERGED_RTOL iterations 3
+   12 KSP Residual norm 4.04453e-05 
+    Linear solve converged due to CONVERGED_RTOL iterations 7
+    Linear solve converged due to CONVERGED_RTOL iterations 4
+   13 KSP Residual norm 2.58232e-05 
     Linear solve converged due to CONVERGED_RTOL iterations 4
-   12 KSP Residual norm 4.63833e-05 
     Linear solve converged due to CONVERGED_RTOL iterations 4
+   14 KSP Residual norm 1.47445e-05 
     Linear solve converged due to CONVERGED_RTOL iterations 4
-   13 KSP Residual norm 3.09627e-05 
     Linear solve converged due to CONVERGED_RTOL iterations 4
+   15 KSP Residual norm 1.15697e-05 
     Linear solve converged due to CONVERGED_RTOL iterations 4
-   14 KSP Residual norm 1.92262e-05 
     Linear solve converged due to CONVERGED_RTOL iterations 4
+   16 KSP Residual norm 9.7534e-06 
     Linear solve converged due to CONVERGED_RTOL iterations 4
-   15 KSP Residual norm 1.16141e-05 
     Linear solve converged due to CONVERGED_RTOL iterations 4
+   17 KSP Residual norm 8.06467e-06 
+    Linear solve converged due to CONVERGED_RTOL iterations 3
     Linear solve converged due to CONVERGED_RTOL iterations 4
-   16 KSP Residual norm 5.33869e-06 
+   18 KSP Residual norm 7.37128e-06 
+    Linear solve converged due to CONVERGED_RTOL iterations 3
+    Linear solve converged due to CONVERGED_RTOL iterations 3
+   19 KSP Residual norm 6.50594e-06 
+    Linear solve converged due to CONVERGED_RTOL iterations 3
     Linear solve converged due to CONVERGED_RTOL iterations 4
+   20 KSP Residual norm 5.17356e-06 
+    Linear solve converged due to CONVERGED_RTOL iterations 3
     Linear solve converged due to CONVERGED_RTOL iterations 4
-   17 KSP Residual norm 2.72645e-06 
+   21 KSP Residual norm 3.4868e-06 
     Linear solve converged due to CONVERGED_RTOL iterations 4
     Linear solve converged due to CONVERGED_RTOL iterations 4
-   18 KSP Residual norm 1.09577e-06 
-  3 KSP Residual norm 1.49647e-08 
+   22 KSP Residual norm 2.08233e-06 
+    Linear solve converged due to CONVERGED_RTOL iterations 5
+    Linear solve converged due to CONVERGED_RTOL iterations 4
+   23 KSP Residual norm 9.20065e-07 
+  3 KSP Residual norm 2.00482e-08 
     Residual norms for fc_fieldsplit_pressure_ solve.
-    0 KSP Residual norm 0.999176 
+    0 KSP Residual norm 0.999577 
     Linear solve converged due to CONVERGED_RTOL iterations 21
-    Linear solve converged due to CONVERGED_RTOL iterations 7
-    1 KSP Residual norm 0.199486 
-    Linear solve converged due to CONVERGED_RTOL iterations 27
-    Linear solve converged due to CONVERGED_RTOL iterations 19
-    2 KSP Residual norm 0.125199 
-    Linear solve converged due to CONVERGED_RTOL iterations 24
-    Linear solve converged due to CONVERGED_RTOL iterations 18
-    3 KSP Residual norm 0.0911001 
-    Linear solve converged due to CONVERGED_RTOL iterations 31
     Linear solve converged due to CONVERGED_RTOL iterations 15
-    4 KSP Residual norm 0.0250606 
-    Linear solve converged due to CONVERGED_RTOL iterations 21
-    Linear solve converged due to CONVERGED_RTOL iterations 12
-    5 KSP Residual norm 0.00393383 
-    Linear solve converged due to CONVERGED_RTOL iterations 21
-    Linear solve converged due to CONVERGED_RTOL iterations 9
-    6 KSP Residual norm 0.00065741 
-    Linear solve converged due to CONVERGED_RTOL iterations 20
-    Linear solve converged due to CONVERGED_RTOL iterations 12
-    7 KSP Residual norm 0.000122053 
-    Linear solve converged due to CONVERGED_RTOL iterations 21
-    Linear solve converged due to CONVERGED_RTOL iterations 10
-    8 KSP Residual norm 2.04244e-05 
+    1 KSP Residual norm 0.845249 
     Linear solve converged due to CONVERGED_RTOL iterations 25
-    Linear solve converged due to CONVERGED_RTOL iterations 10
-    9 KSP Residual norm 2.17897e-06 
-  4 KSP Residual norm 9.456e-11 
-    Residual norms for fc_fieldsplit_pressure_ solve.
-    0 KSP Residual norm 0.999983 
-    Linear solve converged due to CONVERGED_RTOL iterations 21
     Linear solve converged due to CONVERGED_RTOL iterations 18
-    1 KSP Residual norm 0.120117 
-    Linear solve converged due to CONVERGED_RTOL iterations 26
-    Linear solve converged due to CONVERGED_RTOL iterations 15
-    2 KSP Residual norm 0.0913282 
-    Linear solve converged due to CONVERGED_RTOL iterations 23
-    Linear solve converged due to CONVERGED_RTOL iterations 15
-    3 KSP Residual norm 0.0163947 
+    2 KSP Residual norm 0.160848 
     Linear solve converged due to CONVERGED_RTOL iterations 27
     Linear solve converged due to CONVERGED_RTOL iterations 15
-    4 KSP Residual norm 0.0031777 
+    3 KSP Residual norm 0.0868043 
     Linear solve converged due to CONVERGED_RTOL iterations 27
-    Linear solve converged due to CONVERGED_RTOL iterations 12
-    5 KSP Residual norm 0.00039371 
+    Linear solve converged due to CONVERGED_RTOL iterations 15
+    4 KSP Residual norm 0.0192392 
     Linear solve converged due to CONVERGED_RTOL iterations 21
+    Linear solve converged due to CONVERGED_RTOL iterations 12
+    5 KSP Residual norm 0.00315774 
+    Linear solve converged due to CONVERGED_RTOL iterations 28
     Linear solve converged due to CONVERGED_RTOL iterations 9
-    6 KSP Residual norm 3.84559e-05 
-    Linear solve converged due to CONVERGED_RTOL iterations 22
-    Linear solve converged due to CONVERGED_RTOL iterations 10
-    7 KSP Residual norm 7.60788e-06 
-  5 KSP Residual norm < 1.e-11
+    6 KSP Residual norm 0.000505015 
+    Linear solve converged due to CONVERGED_RTOL iterations 20
+    Linear solve converged due to CONVERGED_RTOL iterations 9
+    7 KSP Residual norm 0.000112515 
+    Linear solve converged due to CONVERGED_RTOL iterations 23
+    Linear solve converged due to CONVERGED_RTOL iterations 9
+    8 KSP Residual norm 1.39695e-05 
+    Linear solve converged due to CONVERGED_RTOL iterations 20
+    Linear solve converged due to CONVERGED_RTOL iterations 9
+    9 KSP Residual norm 1.49078e-06 
+  4 KSP Residual norm < 1.e-11
 KSP Object:(fc_) 4 MPI processes
   type: fgmres
     GMRES: restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement
           maximum iterations=10000, initial guess is zero
           tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
           left preconditioning
-          using PRECONDITIONED norm type for convergence test
+          using NONE norm type for convergence test
         PC Object:        (fc_fieldsplit_velocity_sub_)         1 MPI processes        
           type: cholesky
             Cholesky: out-of-place factorization
                   total number of mallocs used during MatSetValues calls =0
                       block size is 1
           linear system matrix = precond matrix:
-          Matrix Object:           1 MPI processes          
+          Matrix Object:          (a11_)           1 MPI processes          
             type: seqaij
             rows=512, cols=512
             total: nonzeros=8140, allocated nonzeros=8140
             total number of mallocs used during MatSetValues calls =0
               using I-node routines: found 279 nodes, limit used is 5
         linear system matrix = precond matrix:
-        Matrix Object:         4 MPI processes        
+        Matrix Object:        (a11_)         4 MPI processes        
           type: mpiaij
           rows=2046, cols=2046
           total: nonzeros=34952, allocated nonzeros=34952
               maximum iterations=10000, initial guess is zero
               tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
               left preconditioning
-              using PRECONDITIONED norm type for convergence test
+              using NONE norm type for convergence test
             PC Object:            (fc_fieldsplit_pressure_lsc_sub_)             1 MPI processes            
               type: icc
                 0 levels of fill
           rows=1024, cols=1024
             Schur complement A11 - A10 inv(A00) A01
             A11
-              Matrix Object:               4 MPI processes              
+              Matrix Object:              (a22_)               4 MPI processes              
                 type: mpiaij
                 rows=1024, cols=1024
                 total: nonzeros=1024, allocated nonzeros=1024
                   maximum iterations=10000, initial guess is zero
                   tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
                   left preconditioning
-                  using PRECONDITIONED norm type for convergence test
+                  using NONE norm type for convergence test
                 PC Object:                (fc_fieldsplit_velocity_sub_)                 1 MPI processes                
                   type: cholesky
                     Cholesky: out-of-place factorization
                           total number of mallocs used during MatSetValues calls =0
                               block size is 1
                   linear system matrix = precond matrix:
-                  Matrix Object:                   1 MPI processes                  
+                  Matrix Object:                  (a11_)                   1 MPI processes                  
                     type: seqaij
                     rows=512, cols=512
                     total: nonzeros=8140, allocated nonzeros=8140
                     total number of mallocs used during MatSetValues calls =0
                       using I-node routines: found 279 nodes, limit used is 5
                 linear system matrix = precond matrix:
-                Matrix Object:                 4 MPI processes                
+                Matrix Object:                (a11_)                 4 MPI processes                
                   type: mpiaij
                   rows=2046, cols=2046
                   total: nonzeros=34952, allocated nonzeros=34952
                 total: nonzeros=7936, allocated nonzeros=7936
                 total number of mallocs used during MatSetValues calls =0
                   using I-node (on process 0) routines: found 279 nodes, limit used is 5
-        Matrix Object:         4 MPI processes        
+        Matrix Object:        (a22_)         4 MPI processes        
           type: mpiaij
           rows=1024, cols=1024
           total: nonzeros=1024, allocated nonzeros=1024
       Matrix object: 
         type=nest, rows=2, cols=2 
         MatNest structure: 
-        (0,0) : name="(null)", type=mpiaij, rows=2046, cols=2046 
-        (0,1) : name="(null)", type=mpiaij, rows=2046, cols=1024 
-        (1,0) : name="(null)", type=mpiaij, rows=1024, cols=2046 
-        (1,1) : name="(null)", type=mpiaij, rows=1024, cols=1024 
+        (0,0) : prefix="a11_", type=mpiaij, rows=2046, cols=2046 
+        (0,1) : type=mpiaij, rows=2046, cols=1024 
+        (1,0) : type=mpiaij, rows=1024, cols=2046 
+        (1,1) : prefix="a22_", type=mpiaij, rows=1024, cols=1024 
 -- vector vector values --
   Min(u)  = -0.017280 [loc=2035]
   Max(u)  = 0.018577 [loc=2024]

File src/mat/impls/nest/matnest.c

     for (i=0; i<bA->nr; i++) {
       for (j=0; j<bA->nc; j++) {
         const MatType type;
-        const char *name;
+        char name[256] = "",prefix[256] = "";
         PetscInt NR,NC;
         PetscBool isNest = PETSC_FALSE;
 
         }
         ierr = MatGetSize(bA->m[i][j],&NR,&NC);CHKERRQ(ierr);
         ierr = MatGetType( bA->m[i][j], &type );CHKERRQ(ierr);
-        name = ((PetscObject)bA->m[i][j])->prefix;
+        if (((PetscObject)bA->m[i][j])->name) {ierr = PetscSNPrintf(name,sizeof name,"name=\"%s\", ",((PetscObject)bA->m[i][j])->name);CHKERRQ(ierr);}
+        if (((PetscObject)bA->m[i][j])->prefix) {ierr = PetscSNPrintf(prefix,sizeof prefix,"prefix=\"%s\", ",((PetscObject)bA->m[i][j])->prefix);CHKERRQ(ierr);}
         ierr = PetscTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);CHKERRQ(ierr);
 
-        PetscViewerASCIIPrintf( viewer, "(%D,%D) : name=\"%s\", type=%s, rows=%D, cols=%D \n",i,j,name,type,NR,NC);
+        ierr = PetscViewerASCIIPrintf(viewer,"(%D,%D) : %s%stype=%s, rows=%D, cols=%D \n",i,j,name,prefix,type,NR,NC);CHKERRQ(ierr);
 
         if (isNest) {
-          PetscViewerASCIIPushTab( viewer );  /* push1 */
-          ierr = MatView( bA->m[i][j], viewer );CHKERRQ(ierr);
-          PetscViewerASCIIPopTab(viewer);    /* pop1 */
+          ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);  /* push1 */
+          ierr = MatView(bA->m[i][j],viewer);CHKERRQ(ierr);
+          ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);    /* pop1 */
         }
       }
     }

File src/vec/vec/impls/nest/vecnest.c

     ierr = PetscViewerASCIIPrintf(viewer,"VecNest  structure: \n");CHKERRQ(ierr);
     for (i=0; i<bx->nb; i++) {
       const VecType type;
-      const char *name;
+      char name[256] = "",prefix[256] = "";
       PetscInt NR;
 
       ierr = VecGetSize(bx->v[i],&NR);CHKERRQ(ierr);
       ierr = VecGetType(bx->v[i],&type);CHKERRQ(ierr);
-      name = ((PetscObject)bx->v[i])->prefix;
+      if (((PetscObject)bx->v[i])->name) {ierr = PetscSNPrintf(name,sizeof name,"name=\"%s\", ",((PetscObject)bx->v[i])->name);CHKERRQ(ierr);}
+      if (((PetscObject)bx->v[i])->prefix) {ierr = PetscSNPrintf(prefix,sizeof prefix,"prefix=\"%s\", ",((PetscObject)bx->v[i])->prefix);CHKERRQ(ierr);}
 
-      ierr = PetscViewerASCIIPrintf(viewer,"(%D) : name=\"%s\", type=%s, rows=%D \n",i,name,type,NR);CHKERRQ(ierr);
+      ierr = PetscViewerASCIIPrintf(viewer,"(%D) : %s%stype=%s, rows=%D \n",i,name,prefix,type,NR);CHKERRQ(ierr);
 
       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);             /* push1 */
       ierr = VecView(bx->v[i],viewer);CHKERRQ(ierr);