Mercurial > fem-fenics-eugenio
view src/interpolate.cc @ 218:8a3361bfa434
interpolate function added
author | Eugenio Gianniti <eugenio.gianniti@mail.polimi.it> |
---|---|
date | Mon, 12 May 2014 21:38:35 +0200 |
parents | |
children | a13b7d744b86 |
line wrap: on
line source
/* 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 "function.h" DEFUN_DLD (interpolate, args, nargout, "-*- texinfo -*-\n\ @deftypefn {Function File} @var{interpolated} =\ interpolate (@var{f}, @var{g})\n\ Interpolate Function @var{g} on @var{f}'s space. \n\ @seealso{Function}\n\ @end deftypefn") { int nargin = args.length (); octave_value retval; if (nargin < 2 || nargin > 2 || nargout > 1) print_usage (); else { if (! function_type_loaded) { function::register_type (); function_type_loaded = true; mlock (); } if (args(0).type_id () == function::static_type_id () && args(1).type_id () == function::static_type_id ()) { const function & u0 = static_cast<const function&> (args(0).get_rep ()); const function & u1 = static_cast<const function&> (args(1).get_rep ()); if (!error_state) { boost::shared_ptr<dolfin::Function> output (new dolfin::Function (u0.get_fun ())); const dolfin::Function & input = u1.get_fun (); output->interpolate (input); std::string name = u1.get_str (); retval = new function (name, output); } } } return retval; }