view src/interpolate.cc @ 253:5e9b5bbdc56b

Support both DOLFIN 1.3.0 and 1.4.0 * src/dolfin_compat.h: use a macro to set the correct shared_ptr (std or boost)
author Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
date Tue, 29 Jul 2014 18:05:56 +0200
parents 958a0e0e8102
children 61830a4f9ab9
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"
#include <stdexcept>
#include "dolfin_compat.h"

DEFUN_DLD (interpolate, args, nargout, "-*- texinfo -*-\n\
@deftypefn {Function File} @var{interp} = \
interpolate (@var{name}, @var{u}, @var{v})\n\
Interpolate a function on a FunctionSpace.\n\
@var{u} is the Function, Expression or Constant \
to interpolate. @var{v} may be a FunctionSpace \
or a Function. In the latter case @var{u} is interpolated \
on the same FunctionSpace where @var{v} is defined.\n\
Note that @var{name} is optional: if not provided, \
it will default to the one assigned to @var{u}.\n\
@seealso{Function, Expression, Constant, FunctionSpace}\n\
@end deftypefn")
{

  int nargin = args.length ();
  octave_value retval;
  
  if (nargin < 2 || nargin > 3 || 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 ();
        }

      octave_idx_type offset = 0;
      std::string name;
      if (args(0).is_string ())
        {
          name = args(0).string_value ();
          offset = 1;
        }

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

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

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

                  if (! error_state)
                    {
                      try { output->interpolate (u0.get_fun ()); }
                      catch (std::runtime_error &)
                        {
                          error ("unable to interpolate on this function space");
                        }
                      if (! error_state)
                        {
                          if (name.empty ())
                            name = u0.get_str ();
                          retval = new function (name, output);
                        }
                    }
                }
              else if (args(0+offset).type_id () ==
                       coefficient::static_type_id ())
                {
                  const coefficient & u0 =
                    static_cast<const coefficient&> (args(0+offset).get_rep ());

                  if (! error_state)
                    {
                      try { output->interpolate (* u0.get_expr ()); }
                      catch (std::runtime_error &)
                        {
                          error ("unable to interpolate on this function space");
                        }
                      if (! error_state)
                        {
                          if (name.empty ())
                            name = u0.get_str ();
                          retval = new function (name, output);
                        }
                    }
                }
              else
                error ("interpolate: invalid arguments");
            }
        }
      else if (args(1+offset).type_id () == function::static_type_id ())
        {
          const function & u0 =
            static_cast<const function&> (args(1+offset).get_rep ());

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

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

                  if (! error_state)
                    {
                      try { output->interpolate (u1.get_fun ()); }
                      catch (std::runtime_error &)
                        {
                          error ("unable to interpolate on this function space");
                        }
                      if (! error_state)
                        {
                          if (name.empty ())
                            name = u1.get_str ();
                          retval = new function (name, output);
                        }
                    }
                }
              else if (args(0+offset).type_id () ==
                       coefficient::static_type_id ())
                {
                  const coefficient & u1 =
                    static_cast<const coefficient&> (args(0+offset).get_rep ());

                  if (! error_state)
                    {
                      try { output->interpolate (* u1.get_expr ()); }
                      catch (std::runtime_error &)
                        {
                          error ("unable to interpolate on this function space");
                        }
                      if (! error_state)
                        {
                          if (name.empty ())
                            name = u1.get_str ();
                          retval = new function (name, output);
                        }
                    }
                }
              else
                error ("interpolate: invalid arguments");
            }
        }
      else
        error ("interpolate: invalid arguments");
    }
  return retval;
}