Mercurial > fem-fenics-eugenio
view src/MeshFunction.cc @ 264:41b76530fe5e
Do not force extension for imported MeshFunction
To be consistent with Mesh
author | Eugenio Gianniti <eugenio.gianniti@mail.polimi.it> |
---|---|
date | Thu, 07 Aug 2014 20:34:50 +0200 |
parents | a61fc34334ca |
children | 53039ac90368 |
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 "mesh.h" #include "meshfunction.h" #include "dolfin_compat.h" DEFUN_DLD (MeshFunction, args, nargout, "-*- texinfo -*-\n\ @deftypefn {Function File} {[@var{meshfunc}]} = \ MeshFunction (@var{Mesh}, @var{filename})\n\ @deftypefnx {Function File} {[@var{meshfunc}]} = \ MeshFunction (@var{type}, @var{Mesh}, @var{filename})\n\ Initialize a meshfunction with the values contained in a file.\n\ The output @var{meshfunc} is an object which contains a representation of the \ meshfunction in @var{filename} which can be used to mark subdomains or \ facets where Dirichlet boundary conditions are to be applied.\n\ Optionally, you can set a label identifying the @var{type} of markers, \ in order to use @var{meshfunc} to mark subdomains in BilinearForm or \ LinearForm. Three labels are recognized: ""dx"" for cell markers, \ ""ds"" for exterior facets and ""dS"" for interior facets.\n\ @seealso{DirichletBC, save, BilinearForm, LinearForm}\n\ @end deftypefn") { int nargin = args.length (); octave_value retval; if (nargin < 2 || nargin > 3 || nargout > 1) print_usage (); else { if (! mesh_type_loaded) { mesh::register_type (); mesh_type_loaded = true; mlock (); } if (! meshfunction_type_loaded) { meshfunction::register_type (); meshfunction_type_loaded = true; mlock (); } octave_idx_type offset = 0; std::string label; if (args(0).is_string ()) { label = args(0).string_value (); offset = 1; } if (args(0 + offset).type_id () == mesh::static_type_id () && args(1 + offset).is_string ()) { std::string filename = args(1 + offset).string_value (); mesh const & msh_arg = static_cast<mesh const &> (args(0 + offset).get_rep ()); if (!error_state) { SHARED_PTR <dolfin::Mesh const> const & pmsh = msh_arg.get_pmsh (); try { retval = new meshfunction (label, pmsh, filename); } catch (std::runtime_error &) { error ("error reading file"); } } } else error ("invalid input arguments"); } return retval; }