Mercurial > fem-fenics-eugenio
view src/SubDomain.cc @ 268:61830a4f9ab9
Improve formatting
author | Eugenio Gianniti <eugenio.gianniti@mail.polimi.it> |
---|---|
date | Thu, 14 Aug 2014 12:26:55 +0200 |
parents | 53039ac90368 |
children | f4d6ae912a08 |
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 <dolfin.h> #include "subdomain.h" DEFUN_DLD (SubDomain, args, nargout, "-*- texinfo -*-\n\ @deftypefn {Function File} {[@var{subdomain}]} = \ SubDomain (@var{inside}, @var{boundary})\n\ Initialize a subdomain with a function handle, @var{inside}, and a binary \ flag. When the latter is true, @var{subdomain} will only contain boundary \ entities. @var{subdomain} can be used to mark a meshfunction. \n\ @seealso{mark, MeshFunction}\n\ @end deftypefn") { int nargin = args.length (); octave_value retval; if (nargin < 2 || nargin > 2 || nargout > 1) { print_usage (); } else { if (args(0).is_function_handle () && args(1).is_bool_scalar ()) { if (! subdomain_type_loaded) { subdomain::register_type (); subdomain_type_loaded = true; mlock (); } octave_fcn_handle * pfh = args(0).fcn_handle_value (); bool on_boundary = args(1).bool_value (); retval = new subdomain (*pfh, on_boundary); } else { error ("invalid input arguments"); } } return retval; }