Mercurial > fem-fenics-eugenio
comparison 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 |
comparison
equal
deleted
inserted
replaced
217:5292e0614efc | 218:8a3361bfa434 |
---|---|
1 /* | |
2 Copyright (C) 2014 Eugenio Gianniti <eugenio.gianniti@mail.polimi.it> | |
3 | |
4 This program is free software; you can redistribute it and/or modify it under | |
5 the terms of the GNU General Public License as published by the Free Software | |
6 Foundation; either version 3 of the License, or (at your option) any later | |
7 version. | |
8 | |
9 This program is distributed in the hope that it will be useful, but WITHOUT | |
10 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or | |
11 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more | |
12 details. | |
13 | |
14 You should have received a copy of the GNU General Public License along with | |
15 this program; if not, see <http://www.gnu.org/licenses/>. | |
16 */ | |
17 | |
18 #include "function.h" | |
19 | |
20 DEFUN_DLD (interpolate, args, nargout, "-*- texinfo -*-\n\ | |
21 @deftypefn {Function File} @var{interpolated} =\ | |
22 interpolate (@var{f}, @var{g})\n\ | |
23 Interpolate Function @var{g} on @var{f}'s space. \n\ | |
24 @seealso{Function}\n\ | |
25 @end deftypefn") | |
26 { | |
27 | |
28 int nargin = args.length (); | |
29 octave_value retval; | |
30 | |
31 if (nargin < 2 || nargin > 2 || nargout > 1) | |
32 print_usage (); | |
33 else | |
34 { | |
35 if (! function_type_loaded) | |
36 { | |
37 function::register_type (); | |
38 function_type_loaded = true; | |
39 mlock (); | |
40 } | |
41 | |
42 if (args(0).type_id () == function::static_type_id () && | |
43 args(1).type_id () == function::static_type_id ()) | |
44 { | |
45 const function & u0 = | |
46 static_cast<const function&> (args(0).get_rep ()); | |
47 const function & u1 = | |
48 static_cast<const function&> (args(1).get_rep ()); | |
49 | |
50 if (!error_state) | |
51 { | |
52 boost::shared_ptr<dolfin::Function> | |
53 output (new dolfin::Function (u0.get_fun ())); | |
54 const dolfin::Function & | |
55 input = u1.get_fun (); | |
56 | |
57 output->interpolate (input); | |
58 std::string name = u1.get_str (); | |
59 | |
60 retval = new function (name, output); | |
61 } | |
62 } | |
63 } | |
64 return retval; | |
65 } |