changeset 259:598c5e9e0a9e

Add utilities for meshfunction * src/save_mf.cc: save meshfunction to .xml * src/MeshFunction.cc: import meshfunction from .xml * src/meshfunction.h: add print method
author Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
date Wed, 06 Aug 2014 16:48:25 +0200
parents ab35a8b0deef
children 1e2a9be8083a
files INDEX src/FILES src/Makefile.in src/MeshFunction.cc src/meshfunction.h src/save.cc src/save_func.cc src/save_mf.cc
diffstat 8 files changed, 253 insertions(+), 71 deletions(-) [+]
line wrap: on
line diff
--- a/INDEX	Fri Aug 01 18:59:18 2014 +0200
+++ b/INDEX	Wed Aug 06 16:48:25 2014 +0200
@@ -14,6 +14,7 @@
   Constant
   Expression
   Function
+  MeshFunction
   DirichletBC
   interpolate
 Definition of the abstract Variational problem
@@ -31,6 +32,7 @@
   @function/plot
   @mesh/plot
   @function/feval
+  @meshfunction/save
 MPI
   barrier
   is_master_node
--- a/src/FILES	Fri Aug 01 18:59:18 2014 +0200
+++ b/src/FILES	Wed Aug 06 16:48:25 2014 +0200
@@ -2,3 +2,4 @@
 *.m
 ./@mesh/
 ./@function/
+./@meshfunction/
--- a/src/Makefile.in	Fri Aug 01 18:59:18 2014 +0200
+++ b/src/Makefile.in	Wed Aug 06 16:48:25 2014 +0200
@@ -33,7 +33,8 @@
              DirichletBC.oct \
              Expression.oct \
              Function.oct \
-             save.oct \
+             save_func.oct \
+             save_mf.oct \
              assemble.oct \
              assemble_system.oct \
              plot_func.oct \
@@ -42,7 +43,8 @@
              feval.oct \
              interpolate.oct \
              is_master_node.oct \
-             barrier.oct
+             barrier.oct \
+             MeshFunction.oct
 
 
 all: $(OCTFILES)
@@ -59,6 +61,12 @@
 Mesh.o: Mesh.cc mesh.h dolfin_compat.h meshfunction.h
 	CPPFLAGS=$(CPPFLAGS) $(MKOCTFILE) -c Mesh.cc -o $@
 
+MeshFunction.oct: MeshFunction.o
+	CPPFLAGS=$(CPPFLAGS) $(MKOCTFILE) $< -o $@ $(LIBS)
+
+MeshFunction.o: MeshFunction.cc mesh.h dolfin_compat.h meshfunction.h
+	CPPFLAGS=$(CPPFLAGS) $(MKOCTFILE) -c $< -o $@
+
 fem_get_mesh.oct: fem_get_mesh.o
 	CPPFLAGS=$(CPPFLAGS) $(MKOCTFILE) fem_get_mesh.o -o $@ $(LIBS)
 
@@ -97,15 +105,24 @@
 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_func.oct: save_func.o mkfunction
+	CPPFLAGS=$(CPPFLAGS) $(MKOCTFILE) $< -o ./@function/save.oct $(LIBS)
+
+save_func.o: save_func.cc function.h dolfin_compat.h
+	CPPFLAGS=$(CPPFLAGS) $(MKOCTFILE) -c $< -o $@
 
-save.o: save.cc dolfin_compat.h
-	CPPFLAGS=$(CPPFLAGS) $(MKOCTFILE) -c save.cc -o $@
+save_mf.oct: save_mf.o mkmeshfunction
+	CPPFLAGS=$(CPPFLAGS) $(MKOCTFILE) $< -o ./@meshfunction/save.oct $(LIBS)
+
+save_mf.o: save_mf.cc meshfunction.h
+	CPPFLAGS=$(CPPFLAGS) $(MKOCTFILE) -c $< -o $@
 
 mkfunction:
 	 mkdir -p @function
 
+mkmeshfunction:
+	 mkdir -p @meshfunction
+
 mkmesh:
 	 mkdir -p @mesh
 
--- /dev/null	Thu Jan 01 00:00:00 1970 +0000
+++ b/src/MeshFunction.cc	Wed Aug 06 16:48:25 2014 +0200
@@ -0,0 +1,79 @@
+/*
+ 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/>.
+*/
+
+#include <dolfin.h>
+#include "mesh.h"
+#include "meshfunction.h"
+#include "dolfin_compat.h"
+
+DEFUN_DLD (MeshFunction, args, nargout, "-*- texinfo -*-\n\
+@deftypefn {Function File} {[@var{meshfunc}]} = \
+MeshFunction (@var{Mesh}, @var{filename})\n\
+Initialize a meshfunction with the values contained in a XDMF file.\n\
+The output @var{meshfunc} is an object which contains a representation of the \
+meshfunction in @var{filename} which can be used to mark subdomains or \
+facets where Dirichlet boundary conditions are to be applied.\n\
+@seealso{DirichletBC, save}\n\
+@end deftypefn")
+{
+
+  int nargin = args.length ();
+  octave_value retval;
+  
+  if (nargin < 2 || nargin > 2 || nargout > 1)
+    print_usage ();
+  else
+    {
+      if (! mesh_type_loaded)
+        {
+          mesh::register_type ();
+          mesh_type_loaded = true;
+          mlock ();
+        }
+
+      if (! meshfunction_type_loaded)
+        {
+          meshfunction::register_type ();
+          meshfunction_type_loaded = true;
+          mlock ();
+        }
+
+      if (args(0).type_id () == mesh::static_type_id ()
+          && args(1).is_string ())
+        {
+          std::string filename = args(1).string_value ();
+          mesh const & msh_arg =
+            static_cast<mesh const &> (args(0).get_rep ());
+
+          if (!error_state)
+            {
+              SHARED_PTR <dolfin::Mesh const> const &
+                pmsh = msh_arg.get_pmsh ();
+              filename += ".xdmf";
+
+              try
+                { retval = new meshfunction (pmsh, filename); }
+              catch (std::runtime_error &)
+                { error ("error reading file"); }
+            }
+        }
+      else
+        error ("invalid input arguments");
+    }
+
+  return retval;
+}
--- a/src/meshfunction.h	Fri Aug 01 18:59:18 2014 +0200
+++ b/src/meshfunction.h	Wed Aug 06 16:48:25 2014 +0200
@@ -32,10 +32,19 @@
   meshfunction (dolfin::MeshFunction <std::size_t> const & _mf)
     : octave_base_value (), pmf (new dolfin::MeshFunction <std::size_t> (_mf)) {}
 
+  meshfunction (SHARED_PTR <dolfin::Mesh const> mesh,
+                std::string const & filename)
+    : octave_base_value (),
+      pmf (new dolfin::MeshFunction <std::size_t> (mesh, filename)) {}
+
   bool
   is_defined (void) const
     { return true; }
 
+  void
+  print (std::ostream& os, bool pr_as_read_syntax = false) const
+    { os << "MeshFunction <std::size_t>: " << get_pmf ()->str (true) << std::endl; }
+
   dolfin::MeshFunction <std::size_t> const &
   get_mf (void) const
     { return *pmf; }
--- a/src/save.cc	Fri Aug 01 18:59:18 2014 +0200
+++ /dev/null	Thu Jan 01 00:00:00 1970 +0000
@@ -1,65 +0,0 @@
-/*
- Copyright (C) 2013 Marco Vassallo <gedeone-octave@users.sourceforge.net>
-
- 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/>.
-*/  
-#include "function.h"
-#include "dolfin_compat.h"
-
-DEFUN_DLD (save, args, , "-*- texinfo -*-\n\
-@deftypefn {Function File} \
-fem_save (@var{Function}, @var{Name})\n\
-Save a function in vtu format.\n\
-The input parameters are\n\
-@itemize @bullet \n\
-@item @var{Function} is the function that you want to save\n\
-@item @var{Name} is a string for the output name\n\
-@end itemize\n\
-The output is a file in format .vtu\n\
-@seealso{plot, Function}\n\
-@end deftypefn")
-{
-
-  int nargin = args.length ();
-  octave_value retval;
-  
-  if (nargin < 2 || nargin > 2)
-    print_usage ();
-  else
-    {
-      if (! function_type_loaded)
-        {
-          function::register_type ();
-          function_type_loaded = true;
-          mlock ();
-        }
-      if (args(0).type_id () == function::static_type_id ())
-        {
-          const function & uo =
-            static_cast<const function&> (args(0).get_rep ());
-          std::string str = args(1).string_value ();
-
-          if (!error_state)
-            {
-              const SHARED_PTR <const dolfin::Function>
-                & u = uo.get_pfun ();
-              str += ".pvd";
-              dolfin::File file (str);
-              file << (*u);
-              retval = 0;
-            }
-        }
-    }
-  return retval;
-}
--- /dev/null	Thu Jan 01 00:00:00 1970 +0000
+++ b/src/save_func.cc	Wed Aug 06 16:48:25 2014 +0200
@@ -0,0 +1,65 @@
+/*
+ Copyright (C) 2013 Marco Vassallo <gedeone-octave@users.sourceforge.net>
+
+ 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/>.
+*/  
+#include "function.h"
+#include "dolfin_compat.h"
+
+DEFUN_DLD (save, args, , "-*- texinfo -*-\n\
+@deftypefn {Function File} \
+save (@var{Function}, @var{Name})\n\
+Save a function in vtu format.\n\
+The input parameters are\n\
+@itemize @bullet \n\
+@item @var{Function} is the function that you want to save\n\
+@item @var{Name} is a string for the output name\n\
+@end itemize\n\
+The output is a file in format .vtu\n\
+@seealso{plot, Function}\n\
+@end deftypefn")
+{
+
+  int nargin = args.length ();
+  octave_value retval;
+  
+  if (nargin < 2 || nargin > 2)
+    print_usage ();
+  else
+    {
+      if (! function_type_loaded)
+        {
+          function::register_type ();
+          function_type_loaded = true;
+          mlock ();
+        }
+      if (args(0).type_id () == function::static_type_id ())
+        {
+          const function & uo =
+            static_cast<const function&> (args(0).get_rep ());
+          std::string str = args(1).string_value ();
+
+          if (!error_state)
+            {
+              const SHARED_PTR <const dolfin::Function>
+                & u = uo.get_pfun ();
+              str += ".pvd";
+              dolfin::File file (str);
+              file << (*u);
+              retval = 0;
+            }
+        }
+    }
+  return retval;
+}
--- /dev/null	Thu Jan 01 00:00:00 1970 +0000
+++ b/src/save_mf.cc	Wed Aug 06 16:48:25 2014 +0200
@@ -0,0 +1,74 @@
+/*
+ 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/>.
+*/
+
+#include "meshfunction.h"
+
+DEFUN_DLD (save, args, nargout, "-*- texinfo -*-\n\
+@deftypefn {Function File} \
+save (@var{MeshFunction}, @var{Name})\n\
+Save a meshfunction in XDMF format.\n\
+The input parameters are\n\
+@itemize @bullet \n\
+@item @var{MeshFunction} is the meshfunction that you want to save\n\
+@item @var{Name} is a string for the output name\n\
+@end itemize\n\
+The output is a .xdmf file.\n\
+@seealso{plot, MeshFunction}\n\
+@end deftypefn")
+{
+  int nargin = args.length ();
+  octave_value retval;
+
+  if (nargin < 2 || nargin > 2 || nargout > 1)
+    print_usage ();
+  else
+    {
+      if (! meshfunction_type_loaded)
+        {
+          meshfunction::register_type ();
+          meshfunction_type_loaded = true;
+          mlock ();
+        }
+
+      if (args(0).type_id () == meshfunction::static_type_id () &&
+          args(1).is_string ())
+        {
+          meshfunction const & mf_arg =
+            static_cast<meshfunction const &> (args(0).get_rep ());
+          std::string str = args(1).string_value ();
+
+          if (!error_state)
+            {
+              dolfin::MeshFunction <std::size_t> const &
+                mf = mf_arg.get_mf ();
+              str += ".xdmf";
+              try
+                {
+                  dolfin::File file (str);
+                  file << mf;
+                }
+              catch (std::runtime_error &)
+                { error ("error saving meshfunction"); }
+              retval = 0;
+            }
+        }
+      else
+        error ("invalid input arguments");
+    }
+
+  return retval;
+}