changeset 253:5e9b5bbdc56b

Support both DOLFIN 1.3.0 and 1.4.0 * src/dolfin_compat.h: use a macro to set the correct shared_ptr (std or boost)
author Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
date Tue, 29 Jul 2014 18:05:56 +0200
parents 7f33554e439a
children 2b51546a28f7
files inst/private/generate_fs.m inst/private/generate_fun.m inst/private/generate_lhs.m inst/private/generate_rhs.m inst/private/get_vars.m.in src/DirichletBC.cc src/Function.cc src/Makefile.in src/Mesh.cc src/SubSpace.cc src/assemble.cc src/assemble_system.cc src/boundarycondition.h src/coefficient.h src/dolfin_compat.h src/feval.cc src/form.h src/function.h src/functionspace.h src/interpolate.cc src/mesh.h src/plot_func.cc src/plot_mesh.cc src/save.cc
diffstat 24 files changed, 132 insertions(+), 64 deletions(-) [+]
line wrap: on
line diff
--- a/inst/private/generate_fs.m	Tue Jul 29 17:11:59 2014 +0200
+++ b/inst/private/generate_fs.m	Tue Jul 29 18:05:56 2014 +0200
@@ -23,6 +23,7 @@
 #include <fem-fenics/functionspace.h>\n\
 #include <fem-fenics/mesh.h>\n\
 #include ""@@UFL_NAME@@.h""\n\
+#include <fem-fenics/dolfin_compat.h>\n\
 \n\
 DEFUN_DLD (@@UFL_NAME@@_FunctionSpace, args, , ""initialize a fs from a mesh declared with fem_init_mesh"")\n\
 {\n\
@@ -45,7 +46,7 @@
         {\n\
           const mesh & msho = static_cast<const mesh&> (args(0).get_rep ());\n\
           const dolfin::Mesh & mshd = msho.get_msh ();\n\
-          boost::shared_ptr <const dolfin::FunctionSpace> g (new @@UFL_NAME@@::FunctionSpace (mshd));\n\
+          SHARED_PTR <const dolfin::FunctionSpace> g (new @@UFL_NAME@@::FunctionSpace (mshd));\n\
 \n\
           if (! functionspace_type_loaded)\n\
             {\n\
--- a/inst/private/generate_fun.m	Tue Jul 29 17:11:59 2014 +0200
+++ b/inst/private/generate_fun.m	Tue Jul 29 18:05:56 2014 +0200
@@ -25,6 +25,7 @@
 #include <fem-fenics/form.h>\n\
 #include <fem-fenics/coefficient.h>\n\
 #include <fem-fenics/function.h>\n\
+#include <fem-fenics/dolfin_compat.h>\n\
 \n\
 DEFUN_DLD (@@UFL_NAME@@_Functional, args, , ""A = fem_rhs_@@UFL_NAME@@ (FUNCTIONAL SPACE, COEFF)"")\n\
 {\n\
@@ -82,7 +83,7 @@
                         = static_cast <const coefficient&> (args(i).get_rep ());\n\
 \n\
                       std::size_t n = M.coefficient_number (cf.get_str ());\n\
-                      const boost::shared_ptr<const expression> & pexp = cf.get_expr ();\n\
+                      const SHARED_PTR <const expression> & pexp = cf.get_expr ();\n\
                       M.set_coefficient (n, pexp);\n\
                       ++nc;\n\
                     }\n\
@@ -93,7 +94,7 @@
                         = static_cast <const function&> (args(i).get_rep ());\n\
 \n\
                       std::size_t n = M.coefficient_number (fun.get_str ());\n\
-                      const boost::shared_ptr<const dolfin::Function> & pfun = fun.get_pfun ();\n\
+                      const SHARED_PTR <const dolfin::Function> & pfun = fun.get_pfun ();\n\
                       M.set_coefficient (n, pfun);\n\
                       ++nc;\n\
                     }\n\
--- a/inst/private/generate_lhs.m	Tue Jul 29 17:11:59 2014 +0200
+++ b/inst/private/generate_lhs.m	Tue Jul 29 18:05:56 2014 +0200
@@ -25,6 +25,7 @@
 #include <fem-fenics/coefficient.h>\n\
 #include <fem-fenics/function.h>\n\
 #include <fem-fenics/functionspace.h>\n\
+#include <fem-fenics/dolfin_compat.h>\n\
 \n\
 DEFUN_DLD (@@UFL_NAME@@_LinearForm, args, , "" b = fem_lhs_@@UFL_NAME@@ (FUNCTIONAL SPACE, COEFF)"")\n\
 {\n\
@@ -82,7 +83,7 @@
                         = static_cast <const coefficient&> (args(i).get_rep ());\n\
 \n\
                       std::size_t n = L.coefficient_number (cf.get_str ());\n\
-                      const boost::shared_ptr<const expression> & pexp = cf.get_expr ();\n\
+                      const SHARED_PTR <const expression> & pexp = cf.get_expr ();\n\
                       L.set_coefficient (n, pexp);\n\
                       ++nc;\n\
                     }\n\
@@ -93,7 +94,7 @@
                         = static_cast <const function&> (args(i).get_rep ());\n\
 \n\
                       std::size_t n = L.coefficient_number (fun.get_str ());\n\
-                      const boost::shared_ptr<const dolfin::Function> & pfun = fun.get_pfun ();\n\
+                      const SHARED_PTR <const dolfin::Function> & pfun = fun.get_pfun ();\n\
                       L.set_coefficient (n, pfun);\n\
                       ++nc;\n\
                     }\n\
--- a/inst/private/generate_rhs.m	Tue Jul 29 17:11:59 2014 +0200
+++ b/inst/private/generate_rhs.m	Tue Jul 29 18:05:56 2014 +0200
@@ -25,6 +25,7 @@
 #include <fem-fenics/coefficient.h>\n\
 #include <fem-fenics/function.h>\n\
 #include <fem-fenics/functionspace.h>\n\
+#include <fem-fenics/dolfin_compat.h>\n\
 \n\
 DEFUN_DLD (@@UFL_NAME@@_BilinearForm, args, , ""A = fem_rhs_@@UFL_NAME@@ (FUNCTIONAL SPACE, COEFF)"")\n\
 {\n\
@@ -86,7 +87,7 @@
                         = static_cast <const coefficient&> (args(i).get_rep ());\n\
 \n\
                       std::size_t n = a.coefficient_number (cf.get_str ());\n\
-                      const boost::shared_ptr<const expression> & pexp = cf.get_expr ();\n\
+                      const SHARED_PTR <const expression> & pexp = cf.get_expr ();\n\
                       a.set_coefficient (n, pexp);\n\
                       ++nc;\n\
                     }\n\
@@ -97,7 +98,7 @@
                         = static_cast <const function&> (args(i).get_rep ());\n\
 \n\
                       std::size_t n = a.coefficient_number (fun.get_str ());\n\
-                      const boost::shared_ptr<const dolfin::Function> & pfun = fun.get_pfun ();\n\
+                      const SHARED_PTR <const dolfin::Function> & pfun = fun.get_pfun ();\n\
                       a.set_coefficient (n, pfun);\n\
                       ++nc;\n\
                     }\n\
--- a/inst/private/get_vars.m.in	Tue Jul 29 17:11:59 2014 +0200
+++ b/inst/private/get_vars.m.in	Tue Jul 29 18:05:56 2014 +0200
@@ -25,6 +25,9 @@
   switch (select)
     case "CPPFLAGS"
       var = "@DENSE_CPPFLAGS@ @DOLFIN_CPPFLAGS@ @EIGEN_CPPFLAGS@";
+      if (! isempty (strfind (var, 'DOLFIN_VERSION="1.4.0"')))
+        var = [var, " -DLATEST_DOLFIN"];
+      endif
     case "LIBS"
       tokens = strsplit (strtrim ("@EIGEN_LIBS@ @DOLFIN_LIBS@"));
       var = strjoin (tokens, '" "');
@@ -32,4 +35,4 @@
     otherwise
       error ("get_vars: the requested variable is not available.");
   endswitch
-endfunction
\ No newline at end of file
+endfunction
--- a/src/DirichletBC.cc	Tue Jul 29 17:11:59 2014 +0200
+++ b/src/DirichletBC.cc	Tue Jul 29 18:05:56 2014 +0200
@@ -18,6 +18,7 @@
 #include "boundarycondition.h"
 #include "functionspace.h"
 #include "expression.h"
+#include "dolfin_compat.h"
 
 DEFUN_DLD (DirichletBC, args, ,
 "-*- texinfo -*-\n\
@@ -72,7 +73,7 @@
 
           if (!error_state)
             {
-              const boost::shared_ptr <const dolfin::FunctionSpace>
+              const SHARED_PTR <const dolfin::FunctionSpace>
                 & V (fspo.get_pfsp ());
 
               octave_value_list b (3, 1);
@@ -86,7 +87,7 @@
               else
                 pf = new expression (*fh);
 
-              boost::shared_ptr<const expression> f (pf);
+              SHARED_PTR <const expression> f (pf);
               boundarycondition * pbc = new boundarycondition ();
 
               for (octave_idx_type i = 0;
--- a/src/Function.cc	Tue Jul 29 17:11:59 2014 +0200
+++ b/src/Function.cc	Tue Jul 29 18:05:56 2014 +0200
@@ -18,6 +18,7 @@
 #include <dolfin.h>
 #include "functionspace.h"
 #include "function.h"
+#include "dolfin_compat.h"
 
 DEFUN_DLD (Function, args, , "-*- texinfo -*-\n\
 @deftypefn {Function File} {[@var{func}]} = \
@@ -81,7 +82,7 @@
 
           if (!error_state)
             {
-              const boost::shared_ptr<const dolfin::FunctionSpace>
+              const SHARED_PTR <const dolfin::FunctionSpace>
                 & V = fspo.get_pfsp ();
 
               if (V->dim () != myu.length ())
@@ -89,14 +90,18 @@
                        and of the vector must agree");
               else
                 {
+#ifdef LATEST_DOLFIN
+                  dolfin::Vector du (MPI_COMM_WORLD, myu.length ());
+#else
                   dolfin::Vector du(myu.length ());
+#endif
                   for (std::size_t i = 0; i < myu.length (); ++i)
                     du.setitem (i, myu(i));
 
-                  boost::shared_ptr<dolfin::Vector> 
+                  SHARED_PTR <dolfin::Vector> 
                     uu (new dolfin::Vector(du));
 
-                  boost::shared_ptr <const dolfin::Function>
+                  SHARED_PTR <const dolfin::Function>
                     u (new dolfin::Function(V, uu));
 
                   retval = new function (str, u);
@@ -112,7 +117,7 @@
 
           if (!error_state)
             {
-              const boost::shared_ptr<const dolfin::Function>
+              const SHARED_PTR <const dolfin::Function>
                 & f = fspo.get_pfun ();
 
               if (f->value_rank () < 1)
@@ -124,7 +129,7 @@
                     error ("Function: index out of bounds");
                  else
                    {
-                      boost::shared_ptr <const dolfin::Function>
+                      SHARED_PTR <const dolfin::Function>
                         u (new dolfin::Function((*f)[idx - 1]));
 
                       retval = new function (str, u);
--- a/src/Makefile.in	Tue Jul 29 17:11:59 2014 +0200
+++ b/src/Makefile.in	Tue Jul 29 18:05:56 2014 +0200
@@ -53,7 +53,7 @@
 Mesh.oct: Mesh.o
 	CPPFLAGS=$(CPPFLAGS) $(MKOCTFILE) Mesh.o -o $@ $(LIBS)
 
-Mesh.o: Mesh.cc mesh.h
+Mesh.o: Mesh.cc mesh.h dolfin_compat.h
 	CPPFLAGS=$(CPPFLAGS) $(MKOCTFILE) -c Mesh.cc -o $@
 
 fem_get_mesh.oct: fem_get_mesh.o
@@ -65,7 +65,7 @@
 DirichletBC.oct: DirichletBC.o
 	CPPFLAGS=$(CPPFLAGS) $(MKOCTFILE) DirichletBC.o -o $@ $(LIBS)
 
-DirichletBC.o: DirichletBC.cc functionspace.h boundarycondition.h expression.h
+DirichletBC.o: DirichletBC.cc functionspace.h boundarycondition.h expression.h dolfin_compat.h
 	CPPFLAGS=$(CPPFLAGS) $(MKOCTFILE) -c DirichletBC.cc -o $@
 
 Expression.oct: Expression.o
@@ -77,27 +77,27 @@
 Function.oct: Function.o
 	CPPFLAGS=$(CPPFLAGS) $(MKOCTFILE) Function.o -o $@ $(LIBS)
 
-Function.o: Function.cc function.h
+Function.o: Function.cc function.h dolfin_compat.h
 	CPPFLAGS=$(CPPFLAGS) $(MKOCTFILE) -c Function.cc -o $@
 
 assemble.oct: assemble.o libfemfenics_factories.a
 	CPPFLAGS=$(CPPFLAGS) $(MKOCTFILE) assemble.o -o $@ \
     libfemfenics_factories.a $(LIBS)
 
-assemble.o: assemble.cc form.h boundarycondition.h femfenics_factory.h
+assemble.o: assemble.cc form.h boundarycondition.h femfenics_factory.h dolfin_compat.h
 	CPPFLAGS=$(CPPFLAGS) $(MKOCTFILE) -c assemble.cc -o $@
 
 assemble_system.oct: assemble_system.o libfemfenics_factories.a
 	CPPFLAGS=$(CPPFLAGS) $(MKOCTFILE) assemble_system.o -o $@ \
     libfemfenics_factories.a $(LIBS)
 
-assemble_system.o: assemble_system.cc form.h boundarycondition.h femfenics_factory.h
+assemble_system.o: assemble_system.cc form.h boundarycondition.h femfenics_factory.h dolfin_compat.h
 	CPPFLAGS=$(CPPFLAGS) $(MKOCTFILE) -c assemble_system.cc -o $@
 
 save.oct: save.o mkfunction
 	CPPFLAGS=$(CPPFLAGS) $(MKOCTFILE) save.o -o ./@function/save.oct $(LIBS)
 
-save.o: save.cc 
+save.o: save.cc dolfin_compat.h
 	CPPFLAGS=$(CPPFLAGS) $(MKOCTFILE) -c save.cc -o $@
 
 mkfunction:
@@ -109,7 +109,7 @@
 plot_mesh.oct: plot_mesh.o mkmesh
 	CPPFLAGS=$(CPPFLAGS) $(MKOCTFILE) plot_mesh.o -o ./@mesh/plot.oct $(LIBS)
 
-plot_mesh.o: plot_mesh.cc Plot_2d.h mesh.h Plot_3d.h
+plot_mesh.o: plot_mesh.cc Plot_2d.h mesh.h Plot_3d.h dolfin_compat.h
 	CPPFLAGS=$(CPPFLAGS) $(MKOCTFILE) -c plot_mesh.cc -o $@
 
 Plot_2d.h: Plot_2d.ufl
@@ -121,22 +121,22 @@
 plot_func.oct: plot_func.o mkfunction
 	CPPFLAGS=$(CPPFLAGS) $(MKOCTFILE) plot_func.o -o ./@function/plot.oct $(LIBS)
 
-plot_func.o: plot_func.cc Plot_2d.h mesh.h Plot_3d.h
+plot_func.o: plot_func.cc Plot_2d.h mesh.h Plot_3d.h dolfin_compat.h
 	CPPFLAGS=$(CPPFLAGS) $(MKOCTFILE) -c plot_func.cc -o $@
 
-SubSpace.oct: SubSpace.cc functionspace.h
+SubSpace.oct: SubSpace.cc functionspace.h dolfin_compat.h
 	CPPFLAGS=$(CPPFLAGS) $(MKOCTFILE) -I. SubSpace.cc $(LIBS)
 
 feval.oct: feval.o mkfunction
 	CPPFLAGS=$(CPPFLAGS) $(MKOCTFILE) feval.o -o ./@function/feval.oct $(LIBS)
 
-feval.o: feval.cc function.h
+feval.o: feval.cc function.h dolfin_compat.h
 	CPPFLAGS=$(CPPFLAGS) $(MKOCTFILE) -c feval.cc -o $@
 
 interpolate.oct: interpolate.o
 	CPPFLAGS=$(CPPFLAGS) $(MKOCTFILE) interpolate.o -o $@ $(LIBS)
 
-interpolate.o: interpolate.cc function.h
+interpolate.o: interpolate.cc function.h dolfin_compat.h
 	CPPFLAGS=$(CPPFLAGS) $(MKOCTFILE) -c interpolate.cc -o $@
 
 femfenics_factory.o: femfenics_factory.cc femfenics_factory.h femfenics_base_factory.h\
--- a/src/Mesh.cc	Tue Jul 29 17:11:59 2014 +0200
+++ b/src/Mesh.cc	Tue Jul 29 18:05:56 2014 +0200
@@ -16,6 +16,7 @@
 */
 
 #include "mesh.h"
+#include "dolfin_compat.h"
 
 DEFUN_DLD (Mesh, args, ,"-*- texinfo -*-\n\
 @deftypefn {Function File} {[@var{mesh_out}]} = \
@@ -90,7 +91,7 @@
   else
     {
       dolfin::MeshEditor editor;
-      boost::shared_ptr<dolfin::Mesh> msh (new dolfin::Mesh ());
+      SHARED_PTR <dolfin::Mesh> msh (new dolfin::Mesh ());
       editor.open (*msh, D, D);
       editor.init_vertices (p.cols ());
       editor.init_cells (t.cols ());
--- a/src/SubSpace.cc	Tue Jul 29 17:11:59 2014 +0200
+++ b/src/SubSpace.cc	Tue Jul 29 18:05:56 2014 +0200
@@ -15,6 +15,7 @@
  this program; if not, see <http://www.gnu.org/licenses/>.
 */
 #include "functionspace.h"
+#include "dolfin_compat.h"
 
 DEFUN_DLD (SubSpace, args, , "-*- texinfo -*-\n\
 @deftypefn {Function File} {[@var{V1}]} = \
@@ -59,7 +60,7 @@
 
               else
                 {
-                  boost::shared_ptr <const dolfin::FunctionSpace>
+                  SHARED_PTR <const dolfin::FunctionSpace>
                     g (new dolfin::SubSpace (V, idx - 1));
 
                   retval = new functionspace (g);
--- a/src/assemble.cc	Tue Jul 29 17:11:59 2014 +0200
+++ b/src/assemble.cc	Tue Jul 29 18:05:56 2014 +0200
@@ -18,6 +18,7 @@
 #include "form.h"
 #include "boundarycondition.h"
 #include "femfenics_factory.h"
+#include "dolfin_compat.h"
 
 DEFUN_DLD (assemble, args, nargout, 
 "-*- texinfo -*-\n\
@@ -95,7 +96,7 @@
                                 (args(i).get_rep ());
 
                               const 
-                                std::vector<boost::shared_ptr
+                                std::vector<SHARED_PTR
                                             <const dolfin::DirichletBC> > 
                                 & pbc = bc.get_bc ();
                               
@@ -125,7 +126,7 @@
                                 = static_cast<const boundarycondition&> 
                                 (args(i).get_rep ());
 
-                              const std::vector<boost::shared_ptr 
+                              const std::vector<SHARED_PTR 
                                                 <const dolfin::DirichletBC> >
                                 & pbc = bc.get_bc ();
 
@@ -175,7 +176,11 @@
                       dolfin::Vector A;
                       dolfin::assemble (A, a);
 
+#ifdef LATEST_DOLFIN
+                      dolfin::Vector x (MPI_COMM_WORLD, myx.length ());
+#else
                       dolfin::Vector x (myx.length ());
+#endif
                       for (std::size_t i = 0; i < myx.length (); ++i)
                         x.setitem (i, myx.xelem (i));
 
@@ -188,7 +193,7 @@
                                 = static_cast<const boundarycondition&> 
                                 (args(i).get_rep ());
 
-                              const std::vector<boost::shared_ptr 
+                              const std::vector<SHARED_PTR 
                                                 <const dolfin::DirichletBC> > 
                                 & pbc = bc.get_bc ();
 
--- a/src/assemble_system.cc	Tue Jul 29 17:11:59 2014 +0200
+++ b/src/assemble_system.cc	Tue Jul 29 18:05:56 2014 +0200
@@ -18,6 +18,7 @@
 #include "form.h"
 #include "boundarycondition.h"
 #include "femfenics_factory.h"
+#include "dolfin_compat.h"
 
 DEFUN_DLD (assemble_system, args, nargout,
 "-*- texinfo -*-\n\
@@ -96,7 +97,7 @@
                                 = static_cast<const boundarycondition&> 
                                 (args(i).get_rep ());
 
-                              const std::vector<boost::shared_ptr 
+                              const std::vector<SHARED_PTR 
                                                 <const dolfin::DirichletBC> > 
                               & pbc = bc.get_bc ();
 
@@ -146,7 +147,11 @@
                       dolfin::assemble (A, a);
                       dolfin::Vector B;
                       dolfin::assemble (B, b);
+#ifdef LATEST_DOLFIN
+                      dolfin::Vector x (MPI_COMM_WORLD, myx.length ());
+#else
                       dolfin::Vector x (myx.length ());
+#endif
 
                       for (std::size_t i = 0; i < myx.length (); ++i)
                         x.setitem (i, myx.xelem (i));
@@ -160,7 +165,7 @@
                                 = static_cast<const boundarycondition&> 
                                   (args(i).get_rep ());
 
-                              const std::vector<boost::shared_ptr 
+                              const std::vector<SHARED_PTR 
                                                <const dolfin::DirichletBC> >
                               & pbc = bc.get_bc ();
 
--- a/src/boundarycondition.h	Tue Jul 29 17:11:59 2014 +0200
+++ b/src/boundarycondition.h	Tue Jul 29 18:05:56 2014 +0200
@@ -22,6 +22,7 @@
 #include <vector>
 #include <dolfin.h>
 #include <octave/oct.h>
+#include "dolfin_compat.h"
 
 class boundarycondition : public octave_base_value
 {
@@ -44,24 +45,24 @@
   is_defined (void) const
     { return true; }
 
-  const std::vector< boost::shared_ptr 
+  const std::vector< SHARED_PTR 
                    <const dolfin::DirichletBC> > &
   get_bc (void) const 
     { return bcu; }
 
   void
-  add_bc (const boost::shared_ptr<const dolfin::FunctionSpace> & V,
-          boost::shared_ptr<const dolfin::GenericFunction> f,
+  add_bc (const SHARED_PTR <const dolfin::FunctionSpace> & V,
+          SHARED_PTR <const dolfin::GenericFunction> f,
           std::size_t n)
     {
-      boost::shared_ptr<const dolfin::DirichletBC> 
+      SHARED_PTR <const dolfin::DirichletBC> 
         p (new dolfin::DirichletBC (V, f, n));
       bcu.push_back(p);
     }
 
  private:
 
-  std::vector<boost::shared_ptr 
+  std::vector<SHARED_PTR 
              <const dolfin::DirichletBC> > bcu;
 
   DECLARE_OCTAVE_ALLOCATOR;
--- a/src/coefficient.h	Tue Jul 29 17:11:59 2014 +0200
+++ b/src/coefficient.h	Tue Jul 29 18:05:56 2014 +0200
@@ -19,6 +19,7 @@
 #define _COEFFICIENT_OCTAVE_
 
 #include "expression.h"
+#include "dolfin_compat.h"
 
 class coefficient : public octave_base_value
 {
@@ -48,7 +49,7 @@
   is_defined (void) const 
     { return true; }
 
-  const boost::shared_ptr <const expression> & 
+  const SHARED_PTR <const expression> & 
   get_expr (void) const
     { return expr; }
 
@@ -59,7 +60,7 @@
  private:
 
   std::string str;
-  boost::shared_ptr <const expression> expr;
+  SHARED_PTR <const expression> expr;
 
   DECLARE_OCTAVE_ALLOCATOR;
   DECLARE_OV_TYPEID_FUNCTIONS_AND_DATA;
--- /dev/null	Thu Jan 01 00:00:00 1970 +0000
+++ b/src/dolfin_compat.h	Tue Jul 29 18:05:56 2014 +0200
@@ -0,0 +1,28 @@
+/*
+ Copyright (C) 2014 Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
+
+ This program is free software; you can redistribute it and/or modify it under
+ the terms of the GNU General Public License as published by the Free Software
+ Foundation; either version 3 of the License, or (at your option) any later
+ version.
+
+ This program is distributed in the hope that it will be useful, but WITHOUT
+ ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
+ FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
+ details.
+
+ You should have received a copy of the GNU General Public License along with
+ this program; if not, see <http://www.gnu.org/licenses/>.
+*/
+
+#ifndef __FEMFENICS_DOLFIN_COMPAT__
+#define __FEMFENICS_DOLFIN_COMPAT__
+
+#define SHARED_PTR boost::shared_ptr
+
+#ifdef LATEST_DOLFIN
+#undef SHARED_PTR
+#define SHARED_PTR std::shared_ptr
+#endif
+
+#endif
--- a/src/feval.cc	Tue Jul 29 17:11:59 2014 +0200
+++ b/src/feval.cc	Tue Jul 29 18:05:56 2014 +0200
@@ -19,6 +19,7 @@
 #include "function.h"
 #include <stdexcept>
 #include <vector>
+#include "dolfin_compat.h"
 
 DEFUN_DLD (feval, args, nargout, "-*- texinfo -*-\n\
 @deftypefn {Function File} {[@var{fx}, @var{fy}, @var{fz}]} = \
@@ -53,7 +54,7 @@
 
       if (!error_state)
         {
-          const boost::shared_ptr<const dolfin::Function> 
+          const SHARED_PTR <const dolfin::Function> 
             & f = fspo.get_pfun ();
 
           octave_idx_type pdim = f->geometric_dimension ();
--- a/src/form.h	Tue Jul 29 17:11:59 2014 +0200
+++ b/src/form.h	Tue Jul 29 18:05:56 2014 +0200
@@ -22,6 +22,7 @@
 #include <vector>
 #include <dolfin.h>
 #include <octave/oct.h>
+#include "dolfin_compat.h"
 
 class form : public octave_base_value
 {
@@ -34,7 +35,7 @@
   form (const dolfin::Form _frm)
     : octave_base_value (), frm (new dolfin::Form (_frm)) {}
 
-  form (boost::shared_ptr <const dolfin::Form> _frm)
+  form (SHARED_PTR <const dolfin::Form> _frm)
     : octave_base_value (), frm (_frm) {}
 
   void
@@ -55,13 +56,13 @@
   get_form (void) const
     { return (*frm); }
 
-  const boost::shared_ptr <const dolfin::Form> & 
+  const SHARED_PTR <const dolfin::Form> & 
   get_pform (void) const
     { return frm; }
 
  private:
 
-  boost::shared_ptr <const dolfin::Form> frm;
+  SHARED_PTR <const dolfin::Form> frm;
 
   DECLARE_OCTAVE_ALLOCATOR;
   DECLARE_OV_TYPEID_FUNCTIONS_AND_DATA;
--- a/src/function.h	Tue Jul 29 17:11:59 2014 +0200
+++ b/src/function.h	Tue Jul 29 18:05:56 2014 +0200
@@ -23,6 +23,7 @@
 #include <dolfin.h>
 #include <octave/oct.h>
 #include <octave/parse.h>
+#include "dolfin_compat.h"
 
 class function : public octave_base_value
 {
@@ -33,7 +34,7 @@
     : octave_base_value (), fun () {}
 
   function (std::string & _str, 
-            boost::shared_ptr <const dolfin::Function> _fun)
+            SHARED_PTR <const dolfin::Function> _fun)
     : octave_base_value (), str(_str), fun (_fun) {}
 
   function (function const& _func)
@@ -54,7 +55,7 @@
   get_fun (void) const
     { return (*fun); }
 
-  const boost::shared_ptr <const dolfin::Function> &
+  const SHARED_PTR <const dolfin::Function> &
   get_pfun (void) const
     { return fun; }
 
@@ -62,7 +63,7 @@
   set_fun (dolfin::Function & _fun)
     {
       dolfin::Function * p = new dolfin::Function (_fun);
-      fun = boost::shared_ptr<const dolfin::Function> (p);
+      fun = SHARED_PTR <const dolfin::Function> (p);
     }
 
   const std::string & 
@@ -121,7 +122,7 @@
  private:
 
   std::string str;
-  boost::shared_ptr <const dolfin::Function> fun;
+  SHARED_PTR <const dolfin::Function> fun;
 
   DECLARE_OCTAVE_ALLOCATOR;
   DECLARE_OV_TYPEID_FUNCTIONS_AND_DATA;
--- a/src/functionspace.h	Tue Jul 29 17:11:59 2014 +0200
+++ b/src/functionspace.h	Tue Jul 29 18:05:56 2014 +0200
@@ -22,6 +22,7 @@
 #include <vector>
 #include <dolfin.h>
 #include <octave/oct.h>
+#include "dolfin_compat.h"
 
 class functionspace : public octave_base_value
 {
@@ -31,7 +32,7 @@
   functionspace () 
     : octave_base_value (), fsp () {}
 
-  functionspace (boost::shared_ptr <const dolfin::FunctionSpace> _fsp)
+  functionspace (SHARED_PTR <const dolfin::FunctionSpace> _fsp)
     : octave_base_value (), fsp (_fsp) {}
 
   void
@@ -48,7 +49,7 @@
   get_fsp (void) const
     { return (*fsp); }
 
-  const boost::shared_ptr <const dolfin::FunctionSpace> &
+  const SHARED_PTR <const dolfin::FunctionSpace> &
   get_pfsp (void) const
     { return fsp; }
 
@@ -56,12 +57,12 @@
   set_fsp (dolfin::FunctionSpace & _fsp)
     {
       dolfin::FunctionSpace * p = new dolfin::FunctionSpace (_fsp);
-      fsp = boost::shared_ptr<const dolfin::FunctionSpace> (p);
+      fsp = SHARED_PTR <const dolfin::FunctionSpace> (p);
     }
 
  private:
 
-  boost::shared_ptr <const dolfin::FunctionSpace> fsp;
+  SHARED_PTR <const dolfin::FunctionSpace> fsp;
 
   DECLARE_OCTAVE_ALLOCATOR;
   DECLARE_OV_TYPEID_FUNCTIONS_AND_DATA;
--- a/src/interpolate.cc	Tue Jul 29 17:11:59 2014 +0200
+++ b/src/interpolate.cc	Tue Jul 29 18:05:56 2014 +0200
@@ -19,6 +19,7 @@
 #include "function.h"
 #include "functionspace.h"
 #include <stdexcept>
+#include "dolfin_compat.h"
 
 DEFUN_DLD (interpolate, args, nargout, "-*- texinfo -*-\n\
 @deftypefn {Function File} @var{interp} = \
@@ -77,7 +78,7 @@
 
           if (! error_state)
             {
-              boost::shared_ptr<dolfin::Function>
+              SHARED_PTR <dolfin::Function>
                 output (new dolfin::Function (u1.get_pfsp ()));
 
               if (args(0+offset).type_id () == function::static_type_id ())
@@ -132,7 +133,7 @@
 
           if (! error_state)
             {
-              boost::shared_ptr<dolfin::Function>
+              SHARED_PTR <dolfin::Function>
                 output (new dolfin::Function (u0.get_fun ()));
 
               if (args(0+offset).type_id () == function::static_type_id ())
--- a/src/mesh.h	Tue Jul 29 17:11:59 2014 +0200
+++ b/src/mesh.h	Tue Jul 29 18:05:56 2014 +0200
@@ -21,6 +21,7 @@
 #include <dolfin.h>
 #include <octave/oct.h>
 #include <octave/oct-map.h>
+#include "dolfin_compat.h"
 
 class mesh : public octave_base_value
 {
@@ -57,7 +58,7 @@
   const dolfin::Mesh & get_msh (void) const
   { return *pmsh; }
 
-  const boost::shared_ptr<const dolfin::Mesh> & 
+  const SHARED_PTR <const dolfin::Mesh> & 
   get_pmsh (void) const
     { return pmsh; }
 
@@ -66,7 +67,7 @@
 
  private:
 
-  boost::shared_ptr<const dolfin::Mesh> pmsh;
+  SHARED_PTR <const dolfin::Mesh> pmsh;
 
   DECLARE_OCTAVE_ALLOCATOR;
   DECLARE_OV_TYPEID_FUNCTIONS_AND_DATA;
--- a/src/plot_func.cc	Tue Jul 29 17:11:59 2014 +0200
+++ b/src/plot_func.cc	Tue Jul 29 18:05:56 2014 +0200
@@ -16,6 +16,7 @@
 */
 
 #include "function.h"
+#include "dolfin_compat.h"
 
 DEFUN_DLD (plot, args, , "-*- texinfo -*-\n\
 @deftypefn {Function File} \
@@ -46,7 +47,7 @@
 
           if (!error_state)
             {
-              const boost::shared_ptr<const dolfin::Function>
+              const SHARED_PTR <const dolfin::Function>
                 & u = uo.get_pfun ();
 
               dolfin::plot (*u);
--- a/src/plot_mesh.cc	Tue Jul 29 17:11:59 2014 +0200
+++ b/src/plot_mesh.cc	Tue Jul 29 18:05:56 2014 +0200
@@ -19,6 +19,7 @@
 #include "mesh.h"
 #include "Plot_2d.h"
 #include "Plot_3d.h"
+#include "dolfin_compat.h"
 
 DEFUN_DLD (plot, args, , "-*- texinfo -*-\n\
 @deftypefn {Function File} \
@@ -52,7 +53,7 @@
             {
               if (!error_state)
                 {
-                  const boost::shared_ptr<const dolfin::Mesh>
+                  const SHARED_PTR <const dolfin::Mesh>
                     & mshd = msho.get_pmsh ();
                   dolfin::plot (*mshd);
                   dolfin::interactive ();
@@ -67,18 +68,18 @@
               if (!error_state)
                 {
                   const dolfin::Mesh & mshd = msho.get_msh ();
-                  boost::shared_ptr <const dolfin::FunctionSpace> V;
+                  SHARED_PTR <const dolfin::FunctionSpace> V;
                   uint D = mshd.topology ().dim ();
                   if (D == 2)
                     {
-                      boost::shared_ptr <const dolfin::FunctionSpace>
+                      SHARED_PTR <const dolfin::FunctionSpace>
                         Vt (new Plot_2d::FunctionSpace (mshd));
                       V = Vt;
                     }
 
                   else if(D == 3)
                     {
-                      boost::shared_ptr <const dolfin::FunctionSpace>
+                      SHARED_PTR <const dolfin::FunctionSpace>
                         Vt (new Plot_3d::FunctionSpace (mshd));
                       V = Vt;
                     }
@@ -92,14 +93,18 @@
                             and the size of the vector must agree");
                   else
                     {
+#ifdef LATEST_DOLFIN
+                      dolfin::Vector du (MPI_COMM_WORLD, myu.length ());
+#else
                       dolfin::Vector du(myu.length ());
+#endif
                       for (std::size_t i = 0; i < myu.length (); ++i)
                         du.setitem (i, myu(i));
 
-                      boost::shared_ptr<dolfin::Vector>
+                      SHARED_PTR <dolfin::Vector>
                         uu (new dolfin::Vector(du));
 
-                      boost::shared_ptr <const dolfin::Function>
+                      SHARED_PTR <const dolfin::Function>
                         u (new dolfin::Function(V, uu));
 
                       dolfin::plot (*u);
--- a/src/save.cc	Tue Jul 29 17:11:59 2014 +0200
+++ b/src/save.cc	Tue Jul 29 18:05:56 2014 +0200
@@ -15,6 +15,7 @@
  this program; if not, see <http://www.gnu.org/licenses/>.
 */  
 #include "function.h"
+#include "dolfin_compat.h"
 
 DEFUN_DLD (save, args, , "-*- texinfo -*-\n\
 @deftypefn {Function File} \
@@ -51,7 +52,7 @@
 
           if (!error_state)
             {
-              const boost::shared_ptr<const dolfin::Function>
+              const SHARED_PTR <const dolfin::Function>
                 & u = uo.get_pfun ();
               str += ".pvd";
               dolfin::File file (str);