Mercurial > fem-fenics-eugenio
view src/save_mf.cc @ 268:61830a4f9ab9
Improve formatting
author | Eugenio Gianniti <eugenio.gianniti@mail.polimi.it> |
---|---|
date | Thu, 14 Aug 2014 12:26:55 +0200 |
parents | 0b9dc516ba36 |
children |
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 "meshfunction.h" DEFUN_DLD (save, args, nargout, "-*- texinfo -*-\n\ @deftypefn {Function File} \ save (@var{MeshFunction}, @var{Name})\n\ Save a meshfunction in XDMF format.\n\ The input parameters are\n\ @itemize @bullet \n\ @item @var{MeshFunction} is the meshfunction you want to save\n\ @item @var{Name} is a string for the output name\n\ @end itemize\n\ The output is a .xdmf file.\n\ @seealso{MeshFunction}\n\ @end deftypefn") { int nargin = args.length (); octave_value retval; if (nargin < 2 || nargin > 2 || nargout > 1) { print_usage (); } else { if (! meshfunction_type_loaded) { meshfunction::register_type (); meshfunction_type_loaded = true; mlock (); } if (args(0).type_id () == meshfunction::static_type_id () && args(1).is_string ()) { meshfunction const & mf_arg = static_cast<meshfunction const &> (args(0).get_rep ()); std::string str = args(1).string_value (); if (!error_state) { dolfin::MeshFunction <std::size_t> const & mf = mf_arg.get_mf (); str += ".xdmf"; try { dolfin::File file (str); file << mf; } catch (std::runtime_error &) { error ("error saving meshfunction"); } retval = 0; } } else { error ("invalid input arguments"); } } return retval; }