view src/DirichletBC.cc @ 268:61830a4f9ab9

Improve formatting
author Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
date Thu, 14 Aug 2014 12:26:55 +0200
parents ab35a8b0deef
children
line wrap: on
line source

/*
 Copyright (C) 2013 Marco Vassallo <gedeone-octave@users.sourceforge.net>

 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 "boundarycondition.h"
#include "functionspace.h"
#include "expression.h"
#include "meshfunction.h"
#include "dolfin_compat.h"

DEFUN_DLD (DirichletBC, args, ,
"-*- texinfo -*-\n\
@deftypefn {Function File} {[@var{bc}]} = \
DirichletBC (@var{FunctionSpace}, [@var{Boundaries},] @var{Boundary_Label},\
             @var{Function_handle}) \n\
Specify essential boundary condition on a specific side.\n\
The input parameters are\n\
@itemize @bullet \n\
@item @var{FunctionSpace} is the fem-fenics space where \
we want to apply the BC\n\
@item @var{Function_handle} is a function handle which contains the expression \
that we want to apply as a BC. If we have a Vector field, we can just use a\n\
vector of function handles: \
@var{Function handle} = ['''@'''(x, y) f1, '''@'''(x, y) f2, ...]\n\
@item @var{Boundary_Label} is an Array which contains the label(s) of the \
side(s) where the BC has to be applied.\n\
@item @var{Boundaries} is an optional MeshFunction storing the markers \
set on mesh facets. In parallel execution you must supply this argument.\n\
@end itemize\n\
The output @var{bc} is an object which contains the boundary conditions\n\
@seealso{Mesh, FunctionSpace}\n\
@end deftypefn")
{
  int nargin = args.length ();
  octave_value retval;

  if (nargin < 3 || nargin > 4)
    { print_usage (); }
  else
    {

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

      if (! boundarycondition_type_loaded)
        {
          boundarycondition::register_type ();
          boundarycondition_type_loaded = true;
          mlock ();
        }

      if (! meshfunction_type_loaded)
        {
          meshfunction::register_type ();
          meshfunction_type_loaded = true;
          mlock ();
        }

      if (args(0).type_id () == functionspace::static_type_id ())
        {
          const functionspace & fspo =
            static_cast<const functionspace &> (args(0).get_rep ());
          octave_fcn_handle * fh = args(1).fcn_handle_value ();

          if (nargin == 3)
            {
              Array<octave_idx_type> side = args(2).array_value ();

              if (!error_state)
                {
                  const SHARED_PTR <const dolfin::FunctionSpace>
                    & V (fspo.get_pfsp ());

                  octave_value_list b (3, 1);
                  octave_value_list tmp = feval (fh->function_value (), b);
                  Array<double> res = tmp(0).array_value ();
                  std::size_t l = res.length ();

                  expression * pf;
                  if (l > 1)
                    { pf = new expression (*fh, l); }
                  else
                    { pf = new expression (*fh); }

                  SHARED_PTR <const expression> f (pf);
                  boundarycondition * pbc = new boundarycondition ();

                  for (octave_idx_type i = 0;
                       i < side.length (); ++i)
                    { pbc->add_bc (V, f, side(i)); }
                  retval = octave_value (pbc);
                }
            }
          else
            {
              Array<octave_idx_type> side = args(3).array_value ();
              meshfunction const & mf =
                static_cast <meshfunction const &> (args(2).get_rep ());

              if (! error_state)
                {
                  const SHARED_PTR <const dolfin::FunctionSpace>
                    & V (fspo.get_pfsp ());

                  octave_value_list b (3, 1);
                  octave_value_list tmp = feval (fh->function_value (), b);
                  Array<double> res = tmp(0).array_value ();
                  std::size_t l = res.length ();

                  expression * pf;
                  if (l > 1)
                    { pf = new expression (*fh, l); }
                  else
                    { pf = new expression (*fh); }

                  SHARED_PTR <const expression> f (pf);
                  boundarycondition * pbc = new boundarycondition ();

                  SHARED_PTR <dolfin::MeshFunction <std::size_t> const> const &
                    subdomains = mf.get_pmf ();

                  for (octave_idx_type i = 0;
                       i < side.length (); ++i)
                    { pbc->add_bc (V, f, subdomains, side (i)); }
                  retval = octave_value (pbc);
                }
            }
        }
    }
  return retval;
}