Index: /issm/trunk-jpl/src/c/solutionsequences/solutionsequence_schurcg.cpp
===================================================================
--- /issm/trunk-jpl/src/c/solutionsequences/solutionsequence_schurcg.cpp	(revision 23924)
+++ /issm/trunk-jpl/src/c/solutionsequences/solutionsequence_schurcg.cpp	(revision 23925)
@@ -88,10 +88,20 @@
 	 * and I the Schur preconditioner (stored here, because the space was allocated either way) 
 	 *         */
+	#if _PETSC_MINOR_>10
+	MatCreateSubMatrix(Kff,isv,isv,MAT_INITIAL_MATRIX,&A);
+	MatCreateSubMatrix(Kff,isv,isp,MAT_INITIAL_MATRIX,&B);
+	MatCreateSubMatrix(Kff,isp,isv,MAT_INITIAL_MATRIX,&BT);
+	#else
 	MatGetSubMatrix(Kff,isv,isv,MAT_INITIAL_MATRIX,&A);
 	MatGetSubMatrix(Kff,isv,isp,MAT_INITIAL_MATRIX,&B);
 	MatGetSubMatrix(Kff,isp,isv,MAT_INITIAL_MATRIX,&BT);
+	#endif
 	
 	/* Extract preconditioner matrix on the pressure space*/
+	#if _PETSC_MINOR_>10
+	MatCreateSubMatrix(Kff,isp,isp,MAT_INITIAL_MATRIX,&IP);
+	#else
 	MatGetSubMatrix(Kff,isp,isp,MAT_INITIAL_MATRIX,&IP);
+	#endif
 
 
@@ -522,14 +532,20 @@
 	}
 
-  /* Note: SchurCG also constructs the Schur preconditioner and stores it in the free block of Kff */
-  /*			[A    B]
+  /* Note: SchurCG also constructs the Schur preconditioner and stores it in the free block of Kff 
+   *			[A    B]
 	* Kff =  |      |
 	*			[B^T IP]
-  /* To calculate the residual, only the necessary blocks need to be extracted */
-
-		/*Extract A, B, B^T */
-		MatGetSubMatrix(Kff->pmatrix->matrix,isv,isv,MAT_INITIAL_MATRIX,&A);
-		MatGetSubMatrix(Kff->pmatrix->matrix,isv,isp,MAT_INITIAL_MATRIX,&B);
-		MatGetSubMatrix(Kff->pmatrix->matrix,isp,isv,MAT_INITIAL_MATRIX,&BT);
+   * To calculate the residual, only the necessary blocks need to be extracted */
+
+	/*Extract A, B, B^T */
+	#if _PETSC_MINOR_>10
+	MatCreateSubMatrix(Kff->pmatrix->matrix,isv,isv,MAT_INITIAL_MATRIX,&A);
+	MatCreateSubMatrix(Kff->pmatrix->matrix,isv,isp,MAT_INITIAL_MATRIX,&B);
+	MatCreateSubMatrix(Kff->pmatrix->matrix,isp,isv,MAT_INITIAL_MATRIX,&BT);
+	#else
+	MatGetSubMatrix(Kff->pmatrix->matrix,isv,isv,MAT_INITIAL_MATRIX,&A);
+	MatGetSubMatrix(Kff->pmatrix->matrix,isv,isp,MAT_INITIAL_MATRIX,&B);
+	MatGetSubMatrix(Kff->pmatrix->matrix,isp,isv,MAT_INITIAL_MATRIX,&BT);
+	#endif
 	
 		/*no. of free nodes in velocity/pressure space*/
@@ -688,6 +704,4 @@
 	Vector<IssmDouble>* df  = NULL;
 	Vector<IssmDouble>* ys  = NULL;
-
-
 
 	/*parameters:*/
@@ -707,5 +721,4 @@
 	int  count=0;
 	bool converged=false;
-
 	
 	/*Start non-linear iteration using input velocity: */
@@ -716,5 +729,4 @@
 	InputUpdateFromConstantx(femmodel,converged,ConvergedEnum);
 	InputUpdateFromSolutionx(femmodel,ug);
-
 	
 	for(;;){
@@ -724,6 +736,4 @@
 		delete ug;
 
-
-		
 		/*Get stiffness matrix and Load vector*/
 		SystemMatricesx(&Kff,&Kfs,&pf,&df,NULL,femmodel);
