Merge lp:~jobh/cbc.block/fenics-1.0-updates into lp:~cbc.block/cbc.block/kent

Proposed by Joachim Haga
Status: Merged
Merged at revision: 91
Proposed branch: lp:~jobh/cbc.block/fenics-1.0-updates
Merge into: lp:~cbc.block/cbc.block/kent
Diff against target: 140 lines (+35/-23)
5 files modified
block/algebraic/trilinos/Epetra.py (+6/-13)
demo/biot.py (+1/-1)
demo/navierstokes.py (+4/-2)
demo/regression-test.sh (+19/-5)
demo/stokes.py (+5/-2)
To merge this branch: bzr merge lp:~jobh/cbc.block/fenics-1.0-updates
Reviewer Review Type Date Requested Status
Kent-Andre Mardal Pending
Review via email: mp+81329@code.launchpad.net

Description of the change

Changes here and there, to support fenics-1.0 and trilinos-8.4. (And a few other small things)

To post a comment you must log in.
91. By Kent-Andre Mardal

merge with jobh

Preview Diff

[H/L] Next/Prev Comment, [J/K] Next/Prev File, [N/P] Next/Prev Hunk
1=== modified file 'block/algebraic/trilinos/Epetra.py'
2--- block/algebraic/trilinos/Epetra.py 2011-08-09 10:59:22 +0000
3+++ block/algebraic/trilinos/Epetra.py 2011-11-04 21:09:23 +0000
4@@ -25,17 +25,10 @@
5 except AttributeError:
6 return NotImplemented
7
8- try:
9- x = self.create_vec()
10- if len(x) != len(b):
11- raise RuntimeError, \
12- 'incompatible dimensions for %s matvec, %d != %d'%(self.__class__.__name__,len(x),len(b))
13- except TypeError:
14- # FIXME: This implies an unnecessary vector copy, and can be
15- # removed when dolfin.EpetraVector is changed to accept BlockMap
16- # instead of just Map (patch posted at https://bugs.launchpad.net/dolfin/+bug/735785)
17- from dolfin import Vector
18- x = Vector(b)
19+ x = self.create_vec()
20+ if len(x) != len(b):
21+ raise RuntimeError, \
22+ 'incompatible dimensions for %s matvec, %d != %d'%(self.__class__.__name__,len(x),len(b))
23
24 x.down_cast().vec().Multiply(1.0, self.v, b_vec, 0.0)
25 return x
26@@ -291,14 +284,14 @@
27 from dolfin import EpetraMatrix
28 rowmap = vec.down_cast().vec().Map()
29 graph = Epetra.CrsGraph(Epetra.Copy, rowmap, 1)
30- indices = numpy.array([0])
31+ indices = numpy.array([0], dtype=numpy.int32)
32 for row in rowmap.MyGlobalElements():
33 indices[0] = row
34 graph.InsertGlobalIndices(row, indices)
35 graph.FillComplete()
36
37 matrix = EpetraMatrix(graph)
38- indices = numpy.array(rowmap.MyGlobalElements(), dtype='I')
39+ indices = numpy.array(rowmap.MyGlobalElements(), dtype=numpy.uint32)
40 if val == 0:
41 matrix.zero(indices)
42 else:
43
44=== modified file 'demo/biot.py'
45--- demo/biot.py 2011-05-08 04:38:28 +0000
46+++ demo/biot.py 2011-11-04 21:09:23 +0000
47@@ -108,7 +108,7 @@
48 h = mesh.hmin()
49 fluid_source_domain = compile_subdomains('{min}<x[0] && x[0]<{max} && {min}<x[1] && x[1]<{max}'
50 .format(min=c-h, max=c+h))
51-topload_source = Expression("-sin(2*t*pi)*sin(x[0]*pi/2)/3")
52+topload_source = Expression("-sin(2*t*pi)*sin(x[0]*pi/2)/3", t=0)
53
54 bc_u_bedrock = DirichletBC(V, [0]*dim, boundary.bottom)
55 bc_u_topload = DirichletBC(V.sub(dim-1), topload_source, boundary.top)
56
57=== modified file 'demo/navierstokes.py'
58--- demo/navierstokes.py 2011-08-09 11:22:44 +0000
59+++ demo/navierstokes.py 2011-11-04 21:09:23 +0000
60@@ -11,6 +11,7 @@
61 # PyTrilinos.
62 import scipy
63 import PyTrilinos
64+import os
65
66 from dolfin import *
67 from block import *
68@@ -20,8 +21,9 @@
69 dolfin.set_log_level(15)
70
71 # Load mesh and subdomains
72-mesh = Mesh("dolfin-2.xml.gz")
73-sub_domains = MeshFunction("uint", mesh, "subdomains.xml.gz")
74+path = os.path.join(os.path.dirname(__file__), '')
75+mesh = Mesh(path+"dolfin-2.xml.gz")
76+sub_domains = MeshFunction("uint", mesh, path+"subdomains.xml.gz")
77
78 # Define function spaces
79 V = VectorFunctionSpace(mesh, "CG", 2)
80
81=== modified file 'demo/regression-test.sh'
82--- demo/regression-test.sh 2011-08-09 11:22:44 +0000
83+++ demo/regression-test.sh 2011-11-04 21:09:23 +0000
84@@ -1,13 +1,27 @@
85-#!/bin/sh
86+#!/bin/bash
87+
88+# Script that runs all cbc.block demos and exits if anything goes wrong. If
89+# this script finishes successfully, we are reasonably confident that there are
90+# no api-related breakages (and other types of error will often lead to a lack
91+# of convergence, so in practice they may be apparent too).
92+
93 set -e
94 export DOLFIN_NOPLOT=1
95
96-if which parallel >/dev/null 2>/dev/null; then
97- parallel -j +0 --halt-on-error=2 -n 1 python ::: *.py
98+cd ${0%/*}
99+demos=$(find . -name \*.py)
100+
101+if which parallel &>/dev/null; then
102+ parallel -j +0 --halt-on-error=2 -v -n 1 python ::: $demos
103 else
104- echo *.py | xargs -n 1 python
105+ for demo in $demos; do
106+ echo python $demo
107+ python $demo
108+ done
109 fi
110
111-# hodge8 is the only one that runs in parallel, since it doesn't use Dirichlet BCs.
112+# Only demos that don't use symmetric Dirichlet BCs can run in parallel
113 mpirun -np 3 python fenics-book/hodge.py N=4
114 mpirun -np 3 python parallelmixedpoisson.py
115+
116+ps -o etime,cputime $$
117
118=== modified file 'demo/stokes.py'
119--- demo/stokes.py 2011-04-26 11:51:16 +0000
120+++ demo/stokes.py 2011-11-04 21:09:23 +0000
121@@ -28,11 +28,14 @@
122 from block.iterative import *
123 from block.algebraic.trilinos import *
124
125+import os
126+
127 dolfin.set_log_level(15)
128
129 # Load mesh and subdomains
130-mesh = Mesh("dolfin-2.xml.gz")
131-sub_domains = MeshFunction("uint", mesh, "subdomains.xml.gz")
132+path = os.path.join(os.path.dirname(__file__), '')
133+mesh = Mesh(path+"dolfin-2.xml.gz")
134+sub_domains = MeshFunction("uint", mesh, path+"subdomains.xml.gz")
135
136 # Define function spaces
137 V = VectorFunctionSpace(mesh, "CG", 2)
138
139=== modified file 'demo/subdomains.xml.gz'
140Binary files demo/subdomains.xml.gz 2011-04-26 11:12:12 +0000 and demo/subdomains.xml.gz 2011-11-04 21:09:23 +0000 differ

Subscribers

People subscribed via source and target branches