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
=== modified file 'block/algebraic/trilinos/Epetra.py'
--- block/algebraic/trilinos/Epetra.py 2011-08-09 10:59:22 +0000
+++ block/algebraic/trilinos/Epetra.py 2011-11-04 21:09:23 +0000
@@ -25,17 +25,10 @@
25 except AttributeError:25 except AttributeError:
26 return NotImplemented26 return NotImplemented
2727
28 try:28 x = self.create_vec()
29 x = self.create_vec()29 if len(x) != len(b):
30 if len(x) != len(b):30 raise RuntimeError, \
31 raise RuntimeError, \31 'incompatible dimensions for %s matvec, %d != %d'%(self.__class__.__name__,len(x),len(b))
32 'incompatible dimensions for %s matvec, %d != %d'%(self.__class__.__name__,len(x),len(b))
33 except TypeError:
34 # FIXME: This implies an unnecessary vector copy, and can be
35 # removed when dolfin.EpetraVector is changed to accept BlockMap
36 # instead of just Map (patch posted at https://bugs.launchpad.net/dolfin/+bug/735785)
37 from dolfin import Vector
38 x = Vector(b)
3932
40 x.down_cast().vec().Multiply(1.0, self.v, b_vec, 0.0)33 x.down_cast().vec().Multiply(1.0, self.v, b_vec, 0.0)
41 return x34 return x
@@ -291,14 +284,14 @@
291 from dolfin import EpetraMatrix284 from dolfin import EpetraMatrix
292 rowmap = vec.down_cast().vec().Map()285 rowmap = vec.down_cast().vec().Map()
293 graph = Epetra.CrsGraph(Epetra.Copy, rowmap, 1)286 graph = Epetra.CrsGraph(Epetra.Copy, rowmap, 1)
294 indices = numpy.array([0])287 indices = numpy.array([0], dtype=numpy.int32)
295 for row in rowmap.MyGlobalElements():288 for row in rowmap.MyGlobalElements():
296 indices[0] = row289 indices[0] = row
297 graph.InsertGlobalIndices(row, indices)290 graph.InsertGlobalIndices(row, indices)
298 graph.FillComplete()291 graph.FillComplete()
299292
300 matrix = EpetraMatrix(graph)293 matrix = EpetraMatrix(graph)
301 indices = numpy.array(rowmap.MyGlobalElements(), dtype='I')294 indices = numpy.array(rowmap.MyGlobalElements(), dtype=numpy.uint32)
302 if val == 0:295 if val == 0:
303 matrix.zero(indices)296 matrix.zero(indices)
304 else:297 else:
305298
=== modified file 'demo/biot.py'
--- demo/biot.py 2011-05-08 04:38:28 +0000
+++ demo/biot.py 2011-11-04 21:09:23 +0000
@@ -108,7 +108,7 @@
108h = mesh.hmin()108h = mesh.hmin()
109fluid_source_domain = compile_subdomains('{min}<x[0] && x[0]<{max} && {min}<x[1] && x[1]<{max}'109fluid_source_domain = compile_subdomains('{min}<x[0] && x[0]<{max} && {min}<x[1] && x[1]<{max}'
110 .format(min=c-h, max=c+h))110 .format(min=c-h, max=c+h))
111topload_source = Expression("-sin(2*t*pi)*sin(x[0]*pi/2)/3")111topload_source = Expression("-sin(2*t*pi)*sin(x[0]*pi/2)/3", t=0)
112112
113bc_u_bedrock = DirichletBC(V, [0]*dim, boundary.bottom)113bc_u_bedrock = DirichletBC(V, [0]*dim, boundary.bottom)
114bc_u_topload = DirichletBC(V.sub(dim-1), topload_source, boundary.top)114bc_u_topload = DirichletBC(V.sub(dim-1), topload_source, boundary.top)
115115
=== modified file 'demo/navierstokes.py'
--- demo/navierstokes.py 2011-08-09 11:22:44 +0000
+++ demo/navierstokes.py 2011-11-04 21:09:23 +0000
@@ -11,6 +11,7 @@
11# PyTrilinos.11# PyTrilinos.
12import scipy12import scipy
13import PyTrilinos13import PyTrilinos
14import os
1415
15from dolfin import *16from dolfin import *
16from block import *17from block import *
@@ -20,8 +21,9 @@
20dolfin.set_log_level(15)21dolfin.set_log_level(15)
2122
22# Load mesh and subdomains23# Load mesh and subdomains
23mesh = Mesh("dolfin-2.xml.gz")24path = os.path.join(os.path.dirname(__file__), '')
24sub_domains = MeshFunction("uint", mesh, "subdomains.xml.gz")25mesh = Mesh(path+"dolfin-2.xml.gz")
26sub_domains = MeshFunction("uint", mesh, path+"subdomains.xml.gz")
2527
26# Define function spaces28# Define function spaces
27V = VectorFunctionSpace(mesh, "CG", 2)29V = VectorFunctionSpace(mesh, "CG", 2)
2830
=== modified file 'demo/regression-test.sh'
--- demo/regression-test.sh 2011-08-09 11:22:44 +0000
+++ demo/regression-test.sh 2011-11-04 21:09:23 +0000
@@ -1,13 +1,27 @@
1#!/bin/sh1#!/bin/bash
2
3# Script that runs all cbc.block demos and exits if anything goes wrong. If
4# this script finishes successfully, we are reasonably confident that there are
5# no api-related breakages (and other types of error will often lead to a lack
6# of convergence, so in practice they may be apparent too).
7
2set -e8set -e
3export DOLFIN_NOPLOT=19export DOLFIN_NOPLOT=1
410
5if which parallel >/dev/null 2>/dev/null; then11cd ${0%/*}
6 parallel -j +0 --halt-on-error=2 -n 1 python ::: *.py12demos=$(find . -name \*.py)
13
14if which parallel &>/dev/null; then
15 parallel -j +0 --halt-on-error=2 -v -n 1 python ::: $demos
7else16else
8 echo *.py | xargs -n 1 python17 for demo in $demos; do
18 echo python $demo
19 python $demo
20 done
9fi21fi
1022
11# hodge8 is the only one that runs in parallel, since it doesn't use Dirichlet BCs.23# Only demos that don't use symmetric Dirichlet BCs can run in parallel
12mpirun -np 3 python fenics-book/hodge.py N=424mpirun -np 3 python fenics-book/hodge.py N=4
13mpirun -np 3 python parallelmixedpoisson.py25mpirun -np 3 python parallelmixedpoisson.py
26
27ps -o etime,cputime $$
1428
=== modified file 'demo/stokes.py'
--- demo/stokes.py 2011-04-26 11:51:16 +0000
+++ demo/stokes.py 2011-11-04 21:09:23 +0000
@@ -28,11 +28,14 @@
28from block.iterative import *28from block.iterative import *
29from block.algebraic.trilinos import *29from block.algebraic.trilinos import *
3030
31import os
32
31dolfin.set_log_level(15)33dolfin.set_log_level(15)
3234
33# Load mesh and subdomains35# Load mesh and subdomains
34mesh = Mesh("dolfin-2.xml.gz")36path = os.path.join(os.path.dirname(__file__), '')
35sub_domains = MeshFunction("uint", mesh, "subdomains.xml.gz")37mesh = Mesh(path+"dolfin-2.xml.gz")
38sub_domains = MeshFunction("uint", mesh, path+"subdomains.xml.gz")
3639
37# Define function spaces40# Define function spaces
38V = VectorFunctionSpace(mesh, "CG", 2)41V = VectorFunctionSpace(mesh, "CG", 2)
3942
=== modified file 'demo/subdomains.xml.gz'
40Binary files demo/subdomains.xml.gz 2011-04-26 11:12:12 +0000 and demo/subdomains.xml.gz 2011-11-04 21:09:23 +0000 differ43Binary 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