Merge lp:~ldeo-magma/dolfin/devfixes into lp:dolfin/1.0.x

Proposed by Cian Wilson
Status: Merged
Merged at revision: 6351
Proposed branch: lp:~ldeo-magma/dolfin/devfixes
Merge into: lp:dolfin/1.0.x
Diff against target: 382 lines (+190/-128)
3 files modified
dolfin/common/defines.cpp (+170/-0)
dolfin/common/defines.h (+16/-124)
site-packages/dolfin_utils/meshconvert.py (+4/-4)
To merge this branch: bzr merge lp:~ldeo-magma/dolfin/devfixes
Reviewer Review Type Date Requested Status
Registry Administrators Pending
Review via email: mp+79484@code.launchpad.net

Description of the change

The following are so trivial they may not be worth your while merging from this branch but I thought I'd flag them anyway.

meshconvert.py still produces xml headers with 'meshfunction' in them, which produces a warning at runtime. Changing to 'mesh_function'.

defines.h was producing multiple definitions of its various subroutines during linking (I have multiple libraries which include dolfin.h). I've split it up into definitions and implementation now and this fixes my problem but don't know if you have a preferred solution.

This all appears to work for me but hasn't been thoroughly tested.

To post a comment you must log in.
Revision history for this message
Johan Hake (johan-hake) wrote :

Cian!

Thanks for the fixes. They all look good. Have you looked at the webpage for
contributing code?

http://fenicsproject.org/contributing/contributing_code.html

There are a copyright consent you need to sign (if you have not done it).
Basically giving the consent that FEniCS can be distributed with LGPL3.0 or
later.

Johan

On Saturday October 15 2011 17:47:29 Cian Wilson wrote:
> Cian Wilson has proposed merging lp:~ldeo-magma/dolfin/devfixes into
> lp:dolfin.
>
> Requested reviews:
> DOLFIN Core Team (dolfin-core)
>
> For more details, see:
> https://code.launchpad.net/~ldeo-magma/dolfin/devfixes/+merge/79484
>
> The following are so trivial they may not be worth your while merging from
> this branch but I thought I'd flag them anyway.
>
> meshconvert.py still produces xml headers with 'meshfunction' in them,
> which produces a warning at runtime. Changing to 'mesh_function'.
>
> defines.h was producing multiple definitions of its various subroutines
> during linking (I have multiple libraries which include dolfin.h). I've
> split it up into definitions and implementation now and this fixes my
> problem but don't know if you have a preferred solution.
>
> This all appears to work for me but hasn't been thoroughly tested.

Revision history for this message
Cian Wilson (cwilson) wrote :

Hi Johan,

I haven't signed the consent form no but am familiar with the procedure. Will get on it straight away - slow bit might be getting my employer to do it, though they've done it in the past.

Cheers,
Cian

Revision history for this message
Cian Wilson (cwilson) wrote :

Johan,

Forms are in the mail.

Cheers,
Cian

Revision history for this message
Johan Hake (johan-hake) wrote :

On Tuesday October 18 2011 12:48:41 Cian Wilson wrote:
> Johan,
>
> Forms are in the mail.

Execellent!

Johan

> Cheers,
> Cian

Preview Diff

[H/L] Next/Prev Comment, [J/K] Next/Prev File, [N/P] Next/Prev Hunk
1=== added file 'dolfin/common/defines.cpp'
2--- dolfin/common/defines.cpp 1970-01-01 00:00:00 +0000
3+++ dolfin/common/defines.cpp 2011-10-16 00:46:23 +0000
4@@ -0,0 +1,170 @@
5+// Copyright (C) 2009-2011 Johan Hake
6+//
7+// This file is part of DOLFIN.
8+//
9+// DOLFIN is free software: you can redistribute it and/or modify
10+// it under the terms of the GNU Lesser General Public License as published by
11+// the Free Software Foundation, either version 3 of the License, or
12+// (at your option) any later version.
13+//
14+// DOLFIN is distributed in the hope that it will be useful,
15+// but WITHOUT ANY WARRANTY; without even the implied warranty of
16+// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17+// GNU Lesser General Public License for more details.
18+//
19+// You should have received a copy of the GNU Lesser General Public License
20+// along with DOLFIN. If not, see <http://www.gnu.org/licenses/>.
21+//
22+// First added: 2011-10-15
23+// Last changed: 2011-10-15
24+
25+#include <dolfin/common/defines.h>
26+
27+// Return true if DOLFIN is compiled with OpenMP
28+bool has_openmp()
29+{
30+#ifdef HAS_OPENMP
31+ return true;
32+#else
33+ return false;
34+#endif
35+}
36+
37+// Return true if DOLFIN is compiled with MPI
38+bool has_mpi()
39+{
40+#ifdef HAS_MPI
41+ return true;
42+#else
43+ return false;
44+#endif
45+}
46+
47+// Return true if DOLFIN is compiled with SLEPc
48+bool has_slepc()
49+{
50+#ifdef HAS_SLEPC
51+ return true;
52+#else
53+ return false;
54+#endif
55+}
56+
57+// Return true if DOLFIN is compiled with Trilinos
58+bool has_trilinos()
59+{
60+#ifdef HAS_TRILINOS
61+ return true;
62+#else
63+ return false;
64+#endif
65+}
66+
67+// Return true if DOLFIN is compiled with Scotch
68+bool has_scotch()
69+{
70+#ifdef HAS_SCOTCH
71+ return true;
72+#else
73+ return false;
74+#endif
75+}
76+
77+// Return true if DOLFIN is compiled with CGAL
78+bool has_cgal()
79+{
80+#ifdef HAS_CGAL
81+ return true;
82+#else
83+ return false;
84+#endif
85+}
86+
87+// Return true if DOLFIN is compiled with Umfpack
88+bool has_umfpack()
89+{
90+#ifdef HAS_UMFPACK
91+ return true;
92+#else
93+ return false;
94+#endif
95+}
96+
97+// Return true if DOLFIN is compiled with Cholmod
98+bool has_cholmod()
99+{
100+#ifdef HAS_CHOLMOD
101+ return true;
102+#else
103+ return false;
104+#endif
105+}
106+
107+// Return true if DOLFIN is compiled with parmetis
108+bool has_parmetis()
109+{
110+#ifdef HAS_PARMETIS
111+ return true;
112+#else
113+ return false;
114+#endif
115+}
116+
117+// Return true if DOLFIN is compiled with GMP
118+bool has_gmp()
119+{
120+#ifdef HAS_GMP
121+ return true;
122+#else
123+ return false;
124+#endif
125+}
126+
127+// Return true if DOLFIN is compiled with ZLIB
128+bool has_zlib()
129+{
130+#ifdef HAS_ZLIB
131+ return true;
132+#else
133+ return false;
134+#endif
135+}
136+
137+// Return true if a specific linear algebra backend is supported
138+bool has_linear_algebra_backend(std::string backend)
139+{
140+ if (backend == "uBLAS")
141+ {
142+ return true;
143+ }
144+ else if (backend == "PETSc")
145+ {
146+#ifdef HAS_PETSC
147+ return true;
148+#else
149+ return false;
150+#endif
151+ }
152+ else if (backend == "Epetra")
153+ {
154+#ifdef HAS_TRILINOS
155+ return true;
156+#else
157+ return false;
158+#endif
159+ }
160+ else if (backend == "MTL4")
161+ {
162+#ifdef HAS_MTL4
163+ return true;
164+#else
165+ return false;
166+#endif
167+ }
168+ else if (backend == "STL")
169+ {
170+ return true;
171+ }
172+ return false;
173+}
174+
175
176=== modified file 'dolfin/common/defines.h'
177--- dolfin/common/defines.h 2011-10-10 00:06:19 +0000
178+++ dolfin/common/defines.h 2011-10-16 00:46:23 +0000
179@@ -19,153 +19,45 @@
180 // First added: 2009-09-03
181 // Last changed: 2011-10-09
182
183+#ifndef __DOLFIN_DEFINES_H
184+#define __DOLFIN_DEFINES_H
185+
186 #include <string>
187
188 // Return true if DOLFIN is compiled with OpenMP
189-bool has_openmp()
190-{
191-#ifdef HAS_OPENMP
192- return true;
193-#else
194- return false;
195-#endif
196-}
197+bool has_openmp();
198
199 // Return true if DOLFIN is compiled with MPI
200-bool has_mpi()
201-{
202-#ifdef HAS_MPI
203- return true;
204-#else
205- return false;
206-#endif
207-}
208+bool has_mpi();
209
210 // Return true if DOLFIN is compiled with SLEPc
211-bool has_slepc()
212-{
213-#ifdef HAS_SLEPC
214- return true;
215-#else
216- return false;
217-#endif
218-}
219+bool has_slepc();
220
221 // Return true if DOLFIN is compiled with Trilinos
222-bool has_trilinos()
223-{
224-#ifdef HAS_TRILINOS
225- return true;
226-#else
227- return false;
228-#endif
229-}
230+bool has_trilinos();
231
232 // Return true if DOLFIN is compiled with Scotch
233-bool has_scotch()
234-{
235-#ifdef HAS_SCOTCH
236- return true;
237-#else
238- return false;
239-#endif
240-}
241+bool has_scotch();
242
243 // Return true if DOLFIN is compiled with CGAL
244-bool has_cgal()
245-{
246-#ifdef HAS_CGAL
247- return true;
248-#else
249- return false;
250-#endif
251-}
252+bool has_cgal();
253
254 // Return true if DOLFIN is compiled with Umfpack
255-bool has_umfpack()
256-{
257-#ifdef HAS_UMFPACK
258- return true;
259-#else
260- return false;
261-#endif
262-}
263+bool has_umfpack();
264
265 // Return true if DOLFIN is compiled with Cholmod
266-bool has_cholmod()
267-{
268-#ifdef HAS_CHOLMOD
269- return true;
270-#else
271- return false;
272-#endif
273-}
274+bool has_cholmod();
275
276 // Return true if DOLFIN is compiled with parmetis
277-bool has_parmetis()
278-{
279-#ifdef HAS_PARMETIS
280- return true;
281-#else
282- return false;
283-#endif
284-}
285+bool has_parmetis();
286
287 // Return true if DOLFIN is compiled with GMP
288-bool has_gmp()
289-{
290-#ifdef HAS_GMP
291- return true;
292-#else
293- return false;
294-#endif
295-}
296+bool has_gmp();
297
298 // Return true if DOLFIN is compiled with ZLIB
299-bool has_zlib()
300-{
301-#ifdef HAS_ZLIB
302- return true;
303-#else
304- return false;
305-#endif
306-}
307+bool has_zlib();
308
309 // Return true if a specific linear algebra backend is supported
310-bool has_linear_algebra_backend(std::string backend)
311-{
312- if (backend == "uBLAS")
313- {
314- return true;
315- }
316- else if (backend == "PETSc")
317- {
318-#ifdef HAS_PETSC
319- return true;
320-#else
321- return false;
322-#endif
323- }
324- else if (backend == "Epetra")
325- {
326-#ifdef HAS_TRILINOS
327- return true;
328-#else
329- return false;
330-#endif
331- }
332- else if (backend == "MTL4")
333- {
334-#ifdef HAS_MTL4
335- return true;
336-#else
337- return false;
338-#endif
339- }
340- else if (backend == "STL")
341- {
342- return true;
343- }
344- return false;
345-}
346+bool has_linear_algebra_backend(std::string backend);
347
348+#endif
349
350=== modified file 'site-packages/dolfin_utils/meshconvert.py'
351--- site-packages/dolfin_utils/meshconvert.py 2011-07-11 15:55:33 +0000
352+++ site-packages/dolfin_utils/meshconvert.py 2011-10-16 00:46:23 +0000
353@@ -742,7 +742,7 @@
354 def write_header_meshfunction(ofile, dimensions, size):
355 header = """<?xml version="1.0" encoding="UTF-8"?>
356 <dolfin xmlns:dolfin="http://www.fenics.org/dolfin/">
357- <meshfunction type="uint" dim="%d" size="%d">
358+ <mesh_function type="uint" dim="%d" size="%d">
359 """ % (dimensions, size)
360 ofile.write(header)
361
362@@ -751,7 +751,7 @@
363 """ % (index, value))
364
365 def write_footer_meshfunction(ofile):
366- ofile.write(""" </meshfunction>
367+ ofile.write(""" </mesh_function>
368 </dolfin>""")
369
370 def diffpack2xml(ifilename, ofilename):
371@@ -763,9 +763,9 @@
372 meshfunction_header = """\
373 <?xml version="1.0" encoding="UTF-8"?>\n
374 <dolfin xmlns:dolfin="http://www.fenics.org/dolfin/">
375- <meshfunction type="uint" dim="%d" size="%d">\n"""
376+ <mesh_function type="uint" dim="%d" size="%d">\n"""
377 meshfunction_entity = " <entity index=\"%d\" value=\"%d\"/>\n"
378- meshfunction_footer = " </meshfunction>\n</dolfin>"
379+ meshfunction_footer = " </mesh_function>\n</dolfin>"
380
381 # Open files
382 ifile = open(ifilename, "r")

Subscribers

People subscribed via source and target branches