view src/interpolate.cc @ 219:a13b7d744b86

Added Python-like interpolate
author Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
date Sat, 17 May 2014 17:16:58 +0200
parents 8a3361bfa434
children 470586565dc7
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 "coefficient.h"
#include "function.h"
#include "functionspace.h"

DEFUN_DLD (interpolate, args, nargout, "-*- texinfo -*-\n\
@deftypefn {Function File} @var{interp} = \
interpolate (@var{Function}, @var{v})\n\
interpolate (@var{v}, @var{FunctionSpace})\n\
Interpolate a function on a FunctionSpace.\n\
This function can be used in two ways:\n\
@itemize @bullet\n\
@item To interpolate a Function, Expression or Constant, \
@var{v}, on the same FunctionSpace of a @var{Function}\n\
@item To interpolate a Function, Expression or Constant, \
@var{v}, on a given @var{FunctionSpace}\n\
@end itemize\n\
@seealso{Function, Expression, Constant, FunctionSpace}\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 (! coefficient_type_loaded)
        {
          coefficient::register_type ();
          coefficient_type_loaded = true;
          mlock ();
        }

      if (! functionspace_type_loaded)
        {
          functionspace::register_type ();
          functionspace_type_loaded = true;
          mlock ();
        }

      if (args(1).type_id () == functionspace::static_type_id ())
        {
          const functionspace & u1 =
            static_cast<const functionspace&> (args(1).get_rep ());

          if (! error_state)
            {
              boost::shared_ptr<dolfin::Function>
                output (new dolfin::Function (u1.get_pfsp ()));

              if (args(0).type_id () == function::static_type_id ())
                {
                  const function & u0 =
                    static_cast<const function&> (args(0).get_rep ());

                  if (! error_state)
                    {
                      output->interpolate (u0.get_fun ());
                      std::string name = u0.get_str ();
                      retval = new function (name, output);
                    }
                }
              else if (args(0).type_id () == coefficient::static_type_id ())
                {
                  const coefficient & u0 =
                    static_cast<const coefficient&> (args(0).get_rep ());

                  if (! error_state)
                    {
                      output->interpolate (* u0.get_expr ());
                      std::string name = u0.get_str ();
                      retval = new function (name, output);
                    }
                }
              else
                error ("interpolate: invalid arguments");
            }
        }
      else if (args(0).type_id () == function::static_type_id ())
        {
          const function & u0 =
            static_cast<const function&> (args(0).get_rep ());

          if (! error_state)
            {
              boost::shared_ptr<dolfin::Function>
                output (new dolfin::Function (u0.get_fun ()));

              if (args(1).type_id () == function::static_type_id ())
                {
                  const function & u1 =
                    static_cast<const function&> (args(1).get_rep ());

                  if (! error_state)
                    {
                      output->interpolate (u1.get_fun ());
                      std::string name = u1.get_str ();
                      retval = new function (name, output);
                    }
                }
              else if (args(1).type_id () == coefficient::static_type_id ())
                {
                  const coefficient & u1 =
                    static_cast<const coefficient&> (args(1).get_rep ());

                  if (! error_state)
                    {
                      output->interpolate (* u1.get_expr ());
                      std::string name = u1.get_str ();
                      retval = new function (name, output);
                    }
                }
              else
                error ("interpolate: invalid arguments");
            }
        }
      else
        error ("interpolate: invalid arguments");
    }
  return retval;
}