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
1=== modified file 'aclocal.m4'
2--- aclocal.m4 2012-06-15 12:17:04 +0000
3+++ aclocal.m4 2012-07-21 21:13:21 +0000
4@@ -497,7 +497,7 @@
5 PetscInt i,j,II,JJ,m,n,its
6 PetscInt Istart,Iend,ione
7 PetscErrorCode ierr
8-#if PETSC_VERSION_MINOR==2
9+#if PETSC_VERSION_MINOR>=2
10 PetscBool flg
11 #else
12 PetscTruth flg
13
14=== modified file 'configure'
15--- configure 2012-07-20 07:26:35 +0000
16+++ configure 2012-07-21 21:13:21 +0000
17@@ -13002,7 +13002,7 @@
18 PetscInt i,j,II,JJ,m,n,its
19 PetscInt Istart,Iend,ione
20 PetscErrorCode ierr
21-#if PETSC_VERSION_MINOR==2
22+#if PETSC_VERSION_MINOR>=2
23 PetscBool flg
24 #else
25 PetscTruth flg
26
27=== modified file 'configure.in'
28=== modified file 'femtools/Multigrid.F90'
29--- femtools/Multigrid.F90 2011-10-26 14:37:06 +0000
30+++ femtools/Multigrid.F90 2012-07-21 21:13:21 +0000
31@@ -40,9 +40,15 @@
32 #include "finclude/petscsys.h"
33 #endif
34 #endif
35-#if PETSC_VERSION_MINOR==2
36+#if PETSC_VERSION_MINOR>=2
37 #define KSP_NORM_NO KSP_NORM_NONE
38 #endif
39+#if PETSC_VERSION_MINOR>=3
40+#define MatCreateSeqAIJ myMatCreateSeqAIJ
41+#define MatCreateMPIAIJ myMatCreateMPIAIJ
42+#define MatCreateSeqBAIJ myMatCreateSeqBAIJ
43+#define MatCreateMPIBAIJ myMatCreateMPIBAIJ
44+#endif
45
46 !! Some parameters that change the behaviour of
47 !! the smoothed aggregation method. All of
48@@ -687,13 +693,17 @@
49 PC:: pc
50 PetscErrorCode:: ierr
51
52+#if PETSC_VERSION_MINOR>=3
53+ call KSPSetType(ksp, KSPCHEBYSHEV, ierr)
54+#else
55 call KSPSetType(ksp, KSPCHEBYCHEV, ierr)
56+#endif
57 call KSPSetOperators(ksp, matrix, matrix, SAME_PRECONDITIONER, ierr)
58 call KSPSetTolerances(ksp, PETSC_DEFAULT_DOUBLE_PRECISION, &
59 PETSC_DEFAULT_DOUBLE_PRECISION, PETSC_DEFAULT_DOUBLE_PRECISION, &
60 iterations, ierr)
61-#ifdef DOUBLEP
62- call KSPChebychevSetEigenvalues(ksp, emax, emin, ierr)
63+#if PETSC_VERSION_MINOR>=3
64+ call KSPChebyshevSetEigenvalues(ksp, emax, emin, ierr)
65 #else
66 call KSPChebychevSetEigenvalues(ksp, emax, emin, ierr)
67 #endif
68@@ -710,7 +720,7 @@
69 integer, intent(out):: maxlevels, coarsesize
70 integer, intent(out):: nosmd, nosmu, clustersize
71
72-#if PETSC_VERSION_MINOR==2
73+#if PETSC_VERSION_MINOR>=2
74 PetscBool flag
75 #else
76 PetscTruth flag
77@@ -811,7 +821,7 @@
78
79 !
80 call VecCopy(diag, sqrt_diag, ierr)
81-#if PETSC_VERSION_MINOR==2
82+#if PETSC_VERSION_MINOR>=3
83 call VecSqrtAbs(sqrt_diag, ierr)
84 #else
85 call VecSqrt(sqrt_diag, ierr)
86
87=== modified file 'femtools/Petsc_Tools.F90'
88--- femtools/Petsc_Tools.F90 2012-05-14 08:49:05 +0000
89+++ femtools/Petsc_Tools.F90 2012-07-21 21:13:21 +0000
90@@ -139,7 +139,13 @@
91 public addup_global_assembly
92 ! for petsc_numbering:
93 public incref, decref, addref
94-
95+#if PETSC_VERSION_MINOR>=3
96+#define MatCreateSeqAIJ myMatCreateSeqAIJ
97+#define MatCreateMPIAIJ myMatCreateMPIAIJ
98+#define MatCreateSeqBAIJ myMatCreateSeqBAIJ
99+#define MatCreateMPIBAIJ myMatCreateMPIBAIJ
100+ public MatCreateSeqAIJ, MatCreateMPIAIJ, MatCreateSeqBAIJ, MatCreateMPIBAIJ
101+#endif
102 contains
103
104 ! Note about definitions in this module:
105@@ -1467,7 +1473,69 @@
106 PetscErrorCode :: ierr
107 call PetscInitialize(PETSC_NULL_CHARACTER, ierr); CHKERRQ(ierr);
108 end subroutine Initialize_Petsc
109-
110+
111+! In petsc-3.3 the MatCreate[B]{Seq|MPI}() routines have changed to MatCreate[B]Aij
112+! and MatSetup always needs to be called
113+#if PETSC_VERSION_MINOR>=3
114+ subroutine MatCreateSeqAIJ(MPI_Comm, nrows, ncols, &
115+ nz, nnz, M, ierr)
116+ integer, intent(in):: MPI_Comm
117+ PetscInt, intent(in):: nrows, ncols, nz
118+ PetscInt, dimension(:), intent(in):: nnz
119+ Mat, intent(out):: M
120+ PetscErrorCode, intent(out):: ierr
121+
122+ call MatCreateAij(MPI_Comm, nrows, ncols, nrows, ncols, &
123+ nz, nnz, 0, PETSC_NULL_INTEGER, M, ierr)
124+ call MatSetup(M, ierr)
125+
126+ end subroutine MatCreateSeqAIJ
127+
128+ subroutine MatCreateMPIAIJ(MPI_Comm, nprows, npcols, &
129+ nrows, ncols, &
130+ dnz, dnnz, onz, onnz, M, ierr)
131+ integer, intent(in):: MPI_Comm
132+ PetscInt, intent(in):: nprows, npcols,nrows, ncols, dnz, onz
133+ PetscInt, dimension(:), intent(in):: dnnz, onnz
134+ Mat, intent(out):: M
135+ PetscErrorCode, intent(out):: ierr
136+
137+ call MatCreateAij(MPI_Comm, nprows, npcols, nrows, ncols, &
138+ dnz, dnnz, onz, onnz, M, ierr)
139+ call MatSetup(M, ierr)
140+
141+ end subroutine MatCreateMPIAIJ
142+
143+ subroutine MatCreateSeqBAIJ(MPI_Comm, bs, nrows, ncols, &
144+ nz, nnz, M, ierr)
145+ integer, intent(in):: MPI_Comm
146+ PetscInt, intent(in):: bs, nrows, ncols, nz
147+ PetscInt, dimension(:), intent(in):: nnz
148+ Mat, intent(out):: M
149+ PetscErrorCode, intent(out):: ierr
150+
151+ call MatCreateBAij(MPI_Comm, bs, nrows, ncols, nrows, ncols, &
152+ nz, nnz, 0, PETSC_NULL_INTEGER, M, ierr)
153+ call MatSetup(M, ierr)
154+
155+ end subroutine MatCreateSeqBAIJ
156+
157+ subroutine MatCreateMPIBAIJ(MPI_Comm, bs, nprows, npcols, &
158+ nrows, ncols, &
159+ dnz, dnnz, onz, onnz, M, ierr)
160+ integer, intent(in):: MPI_Comm
161+ PetscInt, intent(in):: bs, nprows, npcols,nrows, ncols, dnz, onz
162+ PetscInt, dimension(:), intent(in):: dnnz, onnz
163+ Mat, intent(out):: M
164+ PetscErrorCode, intent(out):: ierr
165+
166+ call MatCreateBAij(MPI_Comm, bs, nprows, npcols, nrows, ncols, &
167+ dnz, dnnz, onz, onnz, M, ierr)
168+ call MatSetup(M, ierr)
169+
170+ end subroutine MatCreateMPIBAIJ
171+#endif
172+
173 #include "Reference_count_petsc_numbering_type.F90"
174 end module Petsc_Tools
175
176
177=== modified file 'femtools/Solvers.F90'
178--- femtools/Solvers.F90 2012-04-16 16:41:59 +0000
179+++ femtools/Solvers.F90 2012-07-21 21:13:21 +0000
180@@ -1700,7 +1700,7 @@
181 FLAbort("Need petsc_numbering for monitor")
182 end if
183 call petsc_monitor_setup(petsc_numbering, max_its)
184- call KSPMonitorSet(ksp,MyKSPMonitor,PETSC_NULL_INTEGER, &
185+ call KSPMonitorSet(ksp,MyKSPMonitor,PETSC_NULL_OBJECT, &
186 & PETSC_NULL_FUNCTION,ierr)
187 end if
188
189@@ -1842,7 +1842,7 @@
190 PCType:: pctype
191 PetscReal:: rtol, atol, dtol
192 PetscInt:: maxits
193-#if PETSC_VERSION_MINOR==2
194+#if PETSC_VERSION_MINOR>=2
195 PetscBool:: flag
196 #else
197 PetscTruth:: flag
198@@ -2165,9 +2165,6 @@
199 PC:: pc
200 Vec:: dummy_vec, r, rhs
201
202- ! Stop warnings.
203- ierr = dummy
204-
205 ! Build the solution vector
206 call VecZeroEntries(petsc_monitor_x,ierr)
207 ! don't pass PETSC_NULL_OBJECT instead of dummy_vec, as petsc
208
209=== modified file 'femtools/Sparse_Tools_Petsc.F90'
210--- femtools/Sparse_Tools_Petsc.F90 2012-04-17 08:52:34 +0000
211+++ femtools/Sparse_Tools_Petsc.F90 2012-07-21 21:13:21 +0000
212@@ -73,6 +73,12 @@
213 #include "finclude/petscis.h"
214 #endif
215 #endif
216+#if PETSC_VERSION_MINOR>=3
217+#define MatCreateSeqAIJ myMatCreateSeqAIJ
218+#define MatCreateMPIAIJ myMatCreateMPIAIJ
219+#define MatCreateSeqBAIJ myMatCreateSeqBAIJ
220+#define MatCreateMPIBAIJ myMatCreateMPIBAIJ
221+#endif
222
223 private
224
225
226=== modified file 'femtools/tests/test_multigrid.F90'
227--- femtools/tests/test_multigrid.F90 2011-02-28 16:03:39 +0000
228+++ femtools/tests/test_multigrid.F90 2012-07-21 21:13:21 +0000
229@@ -12,6 +12,7 @@
230 use solvers
231 use fields
232 use parallel_tools
233+#include "petscversion.h"
234 #ifdef HAVE_PETSC
235 #ifdef HAVE_PETSC_MODULES
236 use petsc
237
238=== modified file 'tools/petsc_readnsolve.F90'
239--- tools/petsc_readnsolve.F90 2011-11-17 11:59:39 +0000
240+++ tools/petsc_readnsolve.F90 2012-07-21 21:13:21 +0000
241@@ -70,7 +70,7 @@
242 #endif
243 #endif
244 ! hack around PetscTruth->PetscBool change in petsc 3.2
245-#if PETSC_VERSION_MINOR==2
246+#if PETSC_VERSION_MINOR>=2
247 #define PetscTruth PetscBool
248 #endif
249 ! options read from command-line (-prns_... options)
250@@ -93,7 +93,7 @@
251 ! read in the matrix equation and init. guess:
252 call PetscViewerBinaryOpen(MPI_COMM_FEMTOOLS, trim(filename), &
253 FILE_MODE_READ, viewer, ierr)
254-#if PETSC_VERSION_MINOR==2
255+#if PETSC_VERSION_MINOR>=2
256 ! in petsc 3.2 MatLoad and VecLoad do no longer create a new vector
257 call MatCreate(MPI_COMM_FEMTOOLS, matrix, ierr)
258 call VecCreate(MPI_COMM_FEMTOOLS, rhs, ierr)
259@@ -753,7 +753,7 @@
260 ncomponents=size(petsc_numbering%gnn2unn, 2)
261 allocate(unns(1:n*ncomponents))
262 unns=reshape( petsc_numbering%gnn2unn(1:n,:), (/ n*ncomponents /))
263-#if PETSC_VERSION_MINOR==2
264+#if PETSC_VERSION_MINOR>=2
265 ! for petsc 3.2 we have an extra PetscCopyMode argument
266 call ISCreateGeneral(MPI_COMM_FEMTOOLS, &
267 size(unns), unns, PETSC_COPY_VALUES, row_indexset, ierr)
268
269=== modified file 'tools/petsc_readnsolve_main.cpp'
270--- tools/petsc_readnsolve_main.cpp 2011-10-27 09:56:35 +0000
271+++ tools/petsc_readnsolve_main.cpp 2012-07-21 21:13:21 +0000
272@@ -40,7 +40,7 @@
273 char flml_extension[]=".flml";
274 char *flml_file=NULL;
275 PetscErrorCode ierr;
276-#if PETSC_VERSION_MINOR==2
277+#if PETSC_VERSION_MINOR>=2
278 PetscBool flg;
279 #else
280 PetscTruth flg;
281
282=== modified file 'tools/test_pressure_solve.F90'
283--- tools/test_pressure_solve.F90 2011-10-26 14:37:06 +0000
284+++ tools/test_pressure_solve.F90 2012-07-21 21:13:21 +0000
285@@ -320,7 +320,7 @@
286 #include "finclude/petscpc.h"
287 #endif
288
289-#if PETSC_VERSION_MINOR==2
290+#if PETSC_VERSION_MINOR>=2
291 PetscBool:: flag
292 #else
293 PetscTruth:: flag