changeset 123:31fb4c16dcb0

Polymorphism for the feavl() applied to objevt of type function.
author gedeone-octave <marcovass89@hotmail.it>
date Mon, 02 Sep 2013 23:00:05 +0200
parents 6d4dc85d7f64
children 2191111a1cad
files post_install.m src/Makefile.in src/fem_eval.cc src/feval.cc
diffstat 4 files changed, 77 insertions(+), 68 deletions(-) [+]
line wrap: on
line diff
--- a/post_install.m	Sun Sep 01 23:39:51 2013 +0200
+++ b/post_install.m	Mon Sep 02 23:00:05 2013 +0200
@@ -1,11 +1,20 @@
 
 function post_install (desc)
   private = fullfile (desc.dir, "private");
-  mkdir(private);
+  %mkdir(private);
   [status, msg] = mkdir (private);
   if (status != 1)
     error ("couldn't create private directory: %s", msg);
   endif
 
   movefile ('./src/*.h', private, 'f');
+
+  func = fullfile (desc.dir, "@function");
+  %mkdir(func);
+  [status, msg] = mkdir (func);
+  if (status != 1)
+    error ("couldn't create @function directory: %s", msg);
+  endif
+
+  movefile ('./src/feval.oct', func, 'f');
 endfunction
--- a/src/Makefile.in	Sun Sep 01 23:39:51 2013 +0200
+++ b/src/Makefile.in	Mon Sep 02 23:00:05 2013 +0200
@@ -14,7 +14,7 @@
              Plot_2d.oct \
              Plot_3d.oct \
              SubSpace.oct \
-             fem_eval.oct \
+             feval.oct \
 
 
 LIBS += -ldolfin
@@ -93,8 +93,8 @@
 SubSpace.oct: SubSpace.cc functionspace.h
 	$(MKOCTFILE) $(CPPFLAGS) -I. SubSpace.cc $(LDFLAGS) $(LIBS)
 
-fem_eval.oct: fem_eval.cc function.h
-	$(MKOCTFILE) $(CPPFLAGS) -I. fem_eval.cc $(LDFLAGS) $(LIBS)
+feval.oct: feval.cc function.h
+	$(MKOCTFILE) $(CPPFLAGS) -I. feval.cc $(LDFLAGS) $(LIBS)
 
 Plot_3d.h: Plot_3d.ufl
 	$(FFC) -l dolfin Plot_3d.ufl
--- a/src/fem_eval.cc	Sun Sep 01 23:39:51 2013 +0200
+++ /dev/null	Thu Jan 01 00:00:00 1970 +0000
@@ -1,64 +0,0 @@
-/*
- Copyright (C) 2013 Marco Vassallo
-
- 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 2 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"
-
-DEFUN_DLD (fem_eval, args, , "-*- texinfo -*-\n\
-@deftypefn {Function File} {[@var{value}]} = \
-fem_eval (@var{function_name}, @var{Coordinate})\n\
-@end deftypefn")
-{
-
-  int nargin = args.length ();
-  octave_value retval=0;
-  
-  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 & fspo =
-            static_cast<const function&> (args(0).get_rep ());
-          Array<double> point= args(1).array_value ();
-
-          if (!error_state)
-            {
-              const boost::shared_ptr<const dolfin::Function> & f = fspo.get_pfun ();
-              dim_vector dims;
-              dims.resize (2);
-              dims(0) = 1;
-              dims(1) = f->value_dimension (0);
-              Array<double> res (dims);
-              dolfin::Array<double> x(point.length (), point.fortran_vec ());
-              dolfin::Array<double> values(res.length (), res.fortran_vec ());
-
-              f->eval (values, x);
-
-              retval = octave_value (res);
-            }
-        }
-    }
-  return retval;
-}
--- /dev/null	Thu Jan 01 00:00:00 1970 +0000
+++ b/src/feval.cc	Mon Sep 02 23:00:05 2013 +0200
@@ -0,0 +1,64 @@
+/*
+ Copyright (C) 2013 Marco Vassallo
+
+ 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 2 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"
+
+DEFUN_DLD (feval, args, , "-*- texinfo -*-\n\
+@deftypefn {Function File} {[@var{value}]} = \
+feval (@var{function_name}, @var{Coordinate})\n\
+@end deftypefn")
+{
+
+  int nargin = args.length ();
+  octave_value retval=0;
+  
+  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 & fspo =
+            static_cast<const function&> (args(0).get_rep ());
+          Array<double> point= args(1).array_value ();
+
+          if (!error_state)
+            {
+              const boost::shared_ptr<const dolfin::Function> & f = fspo.get_pfun ();
+              dim_vector dims;
+              dims.resize (2);
+              dims(0) = 1;
+              dims(1) = f->value_dimension (0);
+              Array<double> res (dims);
+              dolfin::Array<double> x(point.length (), point.fortran_vec ());
+              dolfin::Array<double> values(res.length (), res.fortran_vec ());
+
+              f->eval (values, x);
+
+              retval = octave_value (res);
+            }
+        }
+    }
+  return retval;
+}