Merge lp:~g-gorman/fluidity/fluidity-petsc-3.3 into lp:fluidity

Proposed by Gerard
Status: Merged
Merged at revision: 4082
Proposed branch: lp:~g-gorman/fluidity/fluidity-petsc-3.3
Merge into: lp:fluidity
Diff against target: 293 lines (+101/-19)
10 files modified
aclocal.m4 (+1/-1)
configure (+1/-1)
femtools/Multigrid.F90 (+15/-5)
femtools/Petsc_Tools.F90 (+70/-2)
femtools/Solvers.F90 (+2/-5)
femtools/Sparse_Tools_Petsc.F90 (+6/-0)
femtools/tests/test_multigrid.F90 (+1/-0)
tools/petsc_readnsolve.F90 (+3/-3)
tools/petsc_readnsolve_main.cpp (+1/-1)
tools/test_pressure_solve.F90 (+1/-1)
To merge this branch: bzr merge lp:~g-gorman/fluidity/fluidity-petsc-3.3
Reviewer Review Type Date Requested Status
Stephan Kramer Approve
Review via email: mp+116130@code.launchpad.net

Description of the change

This is a branch of fluidity-petsc-dev specifically intended to support PETSc 3.3.

To post a comment you must log in.
Revision history for this message
Stephan Kramer (s-kramer) wrote :

Looks good. This can go once a buildbot queue has been set up and returned green.

review: Approve
Revision history for this message
Lawrence Mitchell (wence) wrote :

femtools/tests/test_multigrid.F90 needs a fix too:

#include "petscversion.h" needs to happen before the use statements too

In general, I don't like all this sprinkling of petsc version knowledge everywhere in the source tree. Can we do something like this:

file_that_needs_petsc.F90:

module foo

  use bar
#include "petsc_modules_maybe.h"
  implicit none
#include "petsc_headers_maybe.h"

  ...
  Write code that is petsc version-agnostic
end module foo

petsc_modules_maybe.h is something like

#ifdef HAVE_PETSC
#ifdef HAVE_PETSC_MODULES
   use petsc
#include "petscversion.h"
#if PETSC_VERSION_MINOR == 0
  use ...
#endif
#endif
#endif

and so forth.

Can we then get away without sprinkling all the code with VecSqrt/VecSqrtAbs etc...?

Revision history for this message
James Robert Percival (j-percival) wrote :

While I remember and before this gets merged, I think femtools/tests/test_petsc_csr_matrix.F90 has slipped through without being changed for support. Post petsc-3.2 the first call to error handlers must be an MPI_Comm communicator. I have the code which does this but it isn't petsc-3.1 safe at the moment.

Revision history for this message
Stephan Kramer (s-kramer) wrote :

Just to recap for those not present at the dev-meeting: we decided to go ahead with this merge as is, so that people depending on petsc 3.3 don't have to wait. A cleanup of the petsc-version specific hackery is still very much appreciated and shall be done in someone's copious free time. Also we discussed possibly dropping support for petsc 3.1 and petsc 3.2, and roll our own petsc 3.3 packages on a ppa as long as Ubuntu's catching up. I'll try to friendly nudge the Debian package maintainer and ask for an upgrade of their package (Debian unstable is on 3.2), so we could simply back(trans?)-port it to Ubuntu.

James: would you mind committing that fix to the unit-test with another ugly "#if PETSC_VERSION_MINOR>=2"...#else ... to make it work with 3.1. Or commit what you have and I'll add the hacks? This stuff will go away as soon as we drop support for 3.1

Revision history for this message
AMCG Buildbot (amcg-buildbot) wrote :

Let me know if you get any movement from Debian and if so I'll push this to the PPA. If not I suspect it's not too tricky to roll our own.

Revision history for this message
Jon Hill (jon-hill) wrote :

Just a note that HECToR currently supports:

petsc/3.0.0.10
petsc/3.1.09(default)
petsc/3.2.01
petsc/3.2.02
petsc-complex/3.0.0.10
petsc-complex/3.1.09(default)
petsc-complex/3.2.01
petsc-complex/3.2.02

So let me know in plenty of time if we drop 3.1 and 3.2 - it'll try and switch the default petsc we use to 3.2.02 soon.

Revision history for this message
Rhodri Davies (rhodri-davies) wrote :

Following on from John's note - I have asked the HECToR team to install a module for PETSC 3.3. This was around 6 weeks back. There is still no module available, but have just given them another poke for an update.

Revision history for this message
James Robert Percival (j-percival) wrote :

Stephan, just to confirm, I don't seem to be able upload to this branch since it's private, so I'm emailing my attempt at a protected version of the test supporting both versions directly to Gerard and you.

Revision history for this message
Stephan Kramer (s-kramer) wrote :

This branch has now been moved to ~/fluidity-core/fluidity/fluidity-petsc-3.3 so other people can commit.

Revision history for this message
Rhodri Davies (rhodri-davies) wrote :

Note that this branch should no longer be considered for merge. The newer version under fluidity-core should be. I will add a further comment there.

Preview Diff

[H/L] Next/Prev Comment, [J/K] Next/Prev File, [N/P] Next/Prev Hunk
=== modified file 'aclocal.m4'
--- aclocal.m4 2012-06-15 12:17:04 +0000
+++ aclocal.m4 2012-07-21 21:13:21 +0000
@@ -497,7 +497,7 @@
497 PetscInt i,j,II,JJ,m,n,its497 PetscInt i,j,II,JJ,m,n,its
498 PetscInt Istart,Iend,ione498 PetscInt Istart,Iend,ione
499 PetscErrorCode ierr499 PetscErrorCode ierr
500#if PETSC_VERSION_MINOR==2500#if PETSC_VERSION_MINOR>=2
501 PetscBool flg501 PetscBool flg
502#else502#else
503 PetscTruth flg503 PetscTruth flg
504504
=== modified file 'configure'
--- configure 2012-07-20 07:26:35 +0000
+++ configure 2012-07-21 21:13:21 +0000
@@ -13002,7 +13002,7 @@
13002 PetscInt i,j,II,JJ,m,n,its13002 PetscInt i,j,II,JJ,m,n,its
13003 PetscInt Istart,Iend,ione13003 PetscInt Istart,Iend,ione
13004 PetscErrorCode ierr13004 PetscErrorCode ierr
13005#if PETSC_VERSION_MINOR==213005#if PETSC_VERSION_MINOR>=2
13006 PetscBool flg13006 PetscBool flg
13007#else13007#else
13008 PetscTruth flg13008 PetscTruth flg
1300913009
=== modified file 'configure.in'
=== modified file 'femtools/Multigrid.F90'
--- femtools/Multigrid.F90 2011-10-26 14:37:06 +0000
+++ femtools/Multigrid.F90 2012-07-21 21:13:21 +0000
@@ -40,9 +40,15 @@
40#include "finclude/petscsys.h"40#include "finclude/petscsys.h"
41#endif41#endif
42#endif42#endif
43#if PETSC_VERSION_MINOR==243#if PETSC_VERSION_MINOR>=2
44#define KSP_NORM_NO KSP_NORM_NONE44#define KSP_NORM_NO KSP_NORM_NONE
45#endif45#endif
46#if PETSC_VERSION_MINOR>=3
47#define MatCreateSeqAIJ myMatCreateSeqAIJ
48#define MatCreateMPIAIJ myMatCreateMPIAIJ
49#define MatCreateSeqBAIJ myMatCreateSeqBAIJ
50#define MatCreateMPIBAIJ myMatCreateMPIBAIJ
51#endif
4652
47!! Some parameters that change the behaviour of 53!! Some parameters that change the behaviour of
48!! the smoothed aggregation method. All of54!! the smoothed aggregation method. All of
@@ -687,13 +693,17 @@
687 PC:: pc693 PC:: pc
688 PetscErrorCode:: ierr694 PetscErrorCode:: ierr
689 695
696#if PETSC_VERSION_MINOR>=3
697 call KSPSetType(ksp, KSPCHEBYSHEV, ierr)
698#else
690 call KSPSetType(ksp, KSPCHEBYCHEV, ierr)699 call KSPSetType(ksp, KSPCHEBYCHEV, ierr)
700#endif
691 call KSPSetOperators(ksp, matrix, matrix, SAME_PRECONDITIONER, ierr)701 call KSPSetOperators(ksp, matrix, matrix, SAME_PRECONDITIONER, ierr)
692 call KSPSetTolerances(ksp, PETSC_DEFAULT_DOUBLE_PRECISION, &702 call KSPSetTolerances(ksp, PETSC_DEFAULT_DOUBLE_PRECISION, &
693 PETSC_DEFAULT_DOUBLE_PRECISION, PETSC_DEFAULT_DOUBLE_PRECISION, &703 PETSC_DEFAULT_DOUBLE_PRECISION, PETSC_DEFAULT_DOUBLE_PRECISION, &
694 iterations, ierr)704 iterations, ierr)
695#ifdef DOUBLEP705#if PETSC_VERSION_MINOR>=3
696 call KSPChebychevSetEigenvalues(ksp, emax, emin, ierr)706 call KSPChebyshevSetEigenvalues(ksp, emax, emin, ierr)
697#else707#else
698 call KSPChebychevSetEigenvalues(ksp, emax, emin, ierr)708 call KSPChebychevSetEigenvalues(ksp, emax, emin, ierr)
699#endif709#endif
@@ -710,7 +720,7 @@
710integer, intent(out):: maxlevels, coarsesize720integer, intent(out):: maxlevels, coarsesize
711integer, intent(out):: nosmd, nosmu, clustersize721integer, intent(out):: nosmd, nosmu, clustersize
712722
713#if PETSC_VERSION_MINOR==2723#if PETSC_VERSION_MINOR>=2
714 PetscBool flag724 PetscBool flag
715#else725#else
716 PetscTruth flag726 PetscTruth flag
@@ -811,7 +821,7 @@
811 821
812 !822 !
813 call VecCopy(diag, sqrt_diag, ierr)823 call VecCopy(diag, sqrt_diag, ierr)
814#if PETSC_VERSION_MINOR==2824#if PETSC_VERSION_MINOR>=3
815 call VecSqrtAbs(sqrt_diag, ierr)825 call VecSqrtAbs(sqrt_diag, ierr)
816#else826#else
817 call VecSqrt(sqrt_diag, ierr)827 call VecSqrt(sqrt_diag, ierr)
818828
=== modified file 'femtools/Petsc_Tools.F90'
--- femtools/Petsc_Tools.F90 2012-05-14 08:49:05 +0000
+++ femtools/Petsc_Tools.F90 2012-07-21 21:13:21 +0000
@@ -139,7 +139,13 @@
139 public addup_global_assembly139 public addup_global_assembly
140 ! for petsc_numbering:140 ! for petsc_numbering:
141 public incref, decref, addref141 public incref, decref, addref
142 142#if PETSC_VERSION_MINOR>=3
143#define MatCreateSeqAIJ myMatCreateSeqAIJ
144#define MatCreateMPIAIJ myMatCreateMPIAIJ
145#define MatCreateSeqBAIJ myMatCreateSeqBAIJ
146#define MatCreateMPIBAIJ myMatCreateMPIBAIJ
147 public MatCreateSeqAIJ, MatCreateMPIAIJ, MatCreateSeqBAIJ, MatCreateMPIBAIJ
148#endif
143contains149contains
144150
145 ! Note about definitions in this module:151 ! Note about definitions in this module:
@@ -1467,7 +1473,69 @@
1467 PetscErrorCode :: ierr1473 PetscErrorCode :: ierr
1468 call PetscInitialize(PETSC_NULL_CHARACTER, ierr); CHKERRQ(ierr);1474 call PetscInitialize(PETSC_NULL_CHARACTER, ierr); CHKERRQ(ierr);
1469 end subroutine Initialize_Petsc1475 end subroutine Initialize_Petsc
1470 1476
1477! In petsc-3.3 the MatCreate[B]{Seq|MPI}() routines have changed to MatCreate[B]Aij
1478! and MatSetup always needs to be called
1479#if PETSC_VERSION_MINOR>=3
1480 subroutine MatCreateSeqAIJ(MPI_Comm, nrows, ncols, &
1481 nz, nnz, M, ierr)
1482 integer, intent(in):: MPI_Comm
1483 PetscInt, intent(in):: nrows, ncols, nz
1484 PetscInt, dimension(:), intent(in):: nnz
1485 Mat, intent(out):: M
1486 PetscErrorCode, intent(out):: ierr
1487
1488 call MatCreateAij(MPI_Comm, nrows, ncols, nrows, ncols, &
1489 nz, nnz, 0, PETSC_NULL_INTEGER, M, ierr)
1490 call MatSetup(M, ierr)
1491
1492 end subroutine MatCreateSeqAIJ
1493
1494 subroutine MatCreateMPIAIJ(MPI_Comm, nprows, npcols, &
1495 nrows, ncols, &
1496 dnz, dnnz, onz, onnz, M, ierr)
1497 integer, intent(in):: MPI_Comm
1498 PetscInt, intent(in):: nprows, npcols,nrows, ncols, dnz, onz
1499 PetscInt, dimension(:), intent(in):: dnnz, onnz
1500 Mat, intent(out):: M
1501 PetscErrorCode, intent(out):: ierr
1502
1503 call MatCreateAij(MPI_Comm, nprows, npcols, nrows, ncols, &
1504 dnz, dnnz, onz, onnz, M, ierr)
1505 call MatSetup(M, ierr)
1506
1507 end subroutine MatCreateMPIAIJ
1508
1509 subroutine MatCreateSeqBAIJ(MPI_Comm, bs, nrows, ncols, &
1510 nz, nnz, M, ierr)
1511 integer, intent(in):: MPI_Comm
1512 PetscInt, intent(in):: bs, nrows, ncols, nz
1513 PetscInt, dimension(:), intent(in):: nnz
1514 Mat, intent(out):: M
1515 PetscErrorCode, intent(out):: ierr
1516
1517 call MatCreateBAij(MPI_Comm, bs, nrows, ncols, nrows, ncols, &
1518 nz, nnz, 0, PETSC_NULL_INTEGER, M, ierr)
1519 call MatSetup(M, ierr)
1520
1521 end subroutine MatCreateSeqBAIJ
1522
1523 subroutine MatCreateMPIBAIJ(MPI_Comm, bs, nprows, npcols, &
1524 nrows, ncols, &
1525 dnz, dnnz, onz, onnz, M, ierr)
1526 integer, intent(in):: MPI_Comm
1527 PetscInt, intent(in):: bs, nprows, npcols,nrows, ncols, dnz, onz
1528 PetscInt, dimension(:), intent(in):: dnnz, onnz
1529 Mat, intent(out):: M
1530 PetscErrorCode, intent(out):: ierr
1531
1532 call MatCreateBAij(MPI_Comm, bs, nprows, npcols, nrows, ncols, &
1533 dnz, dnnz, onz, onnz, M, ierr)
1534 call MatSetup(M, ierr)
1535
1536 end subroutine MatCreateMPIBAIJ
1537#endif
1538
1471#include "Reference_count_petsc_numbering_type.F90"1539#include "Reference_count_petsc_numbering_type.F90"
1472end module Petsc_Tools1540end module Petsc_Tools
14731541
14741542
=== modified file 'femtools/Solvers.F90'
--- femtools/Solvers.F90 2012-04-16 16:41:59 +0000
+++ femtools/Solvers.F90 2012-07-21 21:13:21 +0000
@@ -1700,7 +1700,7 @@
1700 FLAbort("Need petsc_numbering for monitor")1700 FLAbort("Need petsc_numbering for monitor")
1701 end if1701 end if
1702 call petsc_monitor_setup(petsc_numbering, max_its)1702 call petsc_monitor_setup(petsc_numbering, max_its)
1703 call KSPMonitorSet(ksp,MyKSPMonitor,PETSC_NULL_INTEGER, &1703 call KSPMonitorSet(ksp,MyKSPMonitor,PETSC_NULL_OBJECT, &
1704 & PETSC_NULL_FUNCTION,ierr)1704 & PETSC_NULL_FUNCTION,ierr)
1705 end if1705 end if
17061706
@@ -1842,7 +1842,7 @@
1842 PCType:: pctype1842 PCType:: pctype
1843 PetscReal:: rtol, atol, dtol1843 PetscReal:: rtol, atol, dtol
1844 PetscInt:: maxits1844 PetscInt:: maxits
1845#if PETSC_VERSION_MINOR==21845#if PETSC_VERSION_MINOR>=2
1846 PetscBool:: flag1846 PetscBool:: flag
1847#else1847#else
1848 PetscTruth:: flag1848 PetscTruth:: flag
@@ -2165,9 +2165,6 @@
2165 PC:: pc2165 PC:: pc
2166 Vec:: dummy_vec, r, rhs2166 Vec:: dummy_vec, r, rhs
21672167
2168 ! Stop warnings.
2169 ierr = dummy
2170
2171 ! Build the solution vector 2168 ! Build the solution vector
2172 call VecZeroEntries(petsc_monitor_x,ierr)2169 call VecZeroEntries(petsc_monitor_x,ierr)
2173 ! don't pass PETSC_NULL_OBJECT instead of dummy_vec, as petsc2170 ! don't pass PETSC_NULL_OBJECT instead of dummy_vec, as petsc
21742171
=== modified file 'femtools/Sparse_Tools_Petsc.F90'
--- femtools/Sparse_Tools_Petsc.F90 2012-04-17 08:52:34 +0000
+++ femtools/Sparse_Tools_Petsc.F90 2012-07-21 21:13:21 +0000
@@ -73,6 +73,12 @@
73#include "finclude/petscis.h"73#include "finclude/petscis.h"
74#endif74#endif
75#endif75#endif
76#if PETSC_VERSION_MINOR>=3
77#define MatCreateSeqAIJ myMatCreateSeqAIJ
78#define MatCreateMPIAIJ myMatCreateMPIAIJ
79#define MatCreateSeqBAIJ myMatCreateSeqBAIJ
80#define MatCreateMPIBAIJ myMatCreateMPIBAIJ
81#endif
7682
77 private83 private
78 84
7985
=== modified file 'femtools/tests/test_multigrid.F90'
--- femtools/tests/test_multigrid.F90 2011-02-28 16:03:39 +0000
+++ femtools/tests/test_multigrid.F90 2012-07-21 21:13:21 +0000
@@ -12,6 +12,7 @@
12 use solvers12 use solvers
13 use fields13 use fields
14 use parallel_tools14 use parallel_tools
15#include "petscversion.h"
15#ifdef HAVE_PETSC16#ifdef HAVE_PETSC
16#ifdef HAVE_PETSC_MODULES17#ifdef HAVE_PETSC_MODULES
17 use petsc 18 use petsc
1819
=== modified file 'tools/petsc_readnsolve.F90'
--- tools/petsc_readnsolve.F90 2011-11-17 11:59:39 +0000
+++ tools/petsc_readnsolve.F90 2012-07-21 21:13:21 +0000
@@ -70,7 +70,7 @@
70#endif70#endif
71#endif71#endif
72! hack around PetscTruth->PetscBool change in petsc 3.272! hack around PetscTruth->PetscBool change in petsc 3.2
73#if PETSC_VERSION_MINOR==273#if PETSC_VERSION_MINOR>=2
74#define PetscTruth PetscBool74#define PetscTruth PetscBool
75#endif75#endif
76 ! options read from command-line (-prns_... options)76 ! options read from command-line (-prns_... options)
@@ -93,7 +93,7 @@
93 ! read in the matrix equation and init. guess:93 ! read in the matrix equation and init. guess:
94 call PetscViewerBinaryOpen(MPI_COMM_FEMTOOLS, trim(filename), &94 call PetscViewerBinaryOpen(MPI_COMM_FEMTOOLS, trim(filename), &
95 FILE_MODE_READ, viewer, ierr)95 FILE_MODE_READ, viewer, ierr)
96#if PETSC_VERSION_MINOR==296#if PETSC_VERSION_MINOR>=2
97 ! in petsc 3.2 MatLoad and VecLoad do no longer create a new vector97 ! in petsc 3.2 MatLoad and VecLoad do no longer create a new vector
98 call MatCreate(MPI_COMM_FEMTOOLS, matrix, ierr)98 call MatCreate(MPI_COMM_FEMTOOLS, matrix, ierr)
99 call VecCreate(MPI_COMM_FEMTOOLS, rhs, ierr)99 call VecCreate(MPI_COMM_FEMTOOLS, rhs, ierr)
@@ -753,7 +753,7 @@
753 ncomponents=size(petsc_numbering%gnn2unn, 2)753 ncomponents=size(petsc_numbering%gnn2unn, 2)
754 allocate(unns(1:n*ncomponents))754 allocate(unns(1:n*ncomponents))
755 unns=reshape( petsc_numbering%gnn2unn(1:n,:), (/ n*ncomponents /))755 unns=reshape( petsc_numbering%gnn2unn(1:n,:), (/ n*ncomponents /))
756#if PETSC_VERSION_MINOR==2756#if PETSC_VERSION_MINOR>=2
757 ! for petsc 3.2 we have an extra PetscCopyMode argument757 ! for petsc 3.2 we have an extra PetscCopyMode argument
758 call ISCreateGeneral(MPI_COMM_FEMTOOLS, &758 call ISCreateGeneral(MPI_COMM_FEMTOOLS, &
759 size(unns), unns, PETSC_COPY_VALUES, row_indexset, ierr)759 size(unns), unns, PETSC_COPY_VALUES, row_indexset, ierr)
760760
=== modified file 'tools/petsc_readnsolve_main.cpp'
--- tools/petsc_readnsolve_main.cpp 2011-10-27 09:56:35 +0000
+++ tools/petsc_readnsolve_main.cpp 2012-07-21 21:13:21 +0000
@@ -40,7 +40,7 @@
40 char flml_extension[]=".flml";40 char flml_extension[]=".flml";
41 char *flml_file=NULL;41 char *flml_file=NULL;
42 PetscErrorCode ierr;42 PetscErrorCode ierr;
43#if PETSC_VERSION_MINOR==243#if PETSC_VERSION_MINOR>=2
44 PetscBool flg;44 PetscBool flg;
45#else45#else
46 PetscTruth flg;46 PetscTruth flg;
4747
=== modified file 'tools/test_pressure_solve.F90'
--- tools/test_pressure_solve.F90 2011-10-26 14:37:06 +0000
+++ tools/test_pressure_solve.F90 2012-07-21 21:13:21 +0000
@@ -320,7 +320,7 @@
320#include "finclude/petscpc.h"320#include "finclude/petscpc.h"
321#endif321#endif
322322
323#if PETSC_VERSION_MINOR==2323#if PETSC_VERSION_MINOR>=2
324 PetscBool:: flag324 PetscBool:: flag
325#else325#else
326 PetscTruth:: flag326 PetscTruth:: flag