Mercurial > fem-fenics-eugenio
view src/subdomain.h @ 267:53039ac90368
Mark meshfunction using subdomain
* src/SubDomain.cc: return subdomain to use for marking
* src/MeshFunction.cc: return a user defined meshfunction
* src/mark.cc: with the information from subdomain, mark a meshfunction
author | Eugenio Gianniti <eugenio.gianniti@mail.polimi.it> |
---|---|
date | Tue, 12 Aug 2014 15:42:50 +0200 |
parents | |
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/>. */ #ifndef _SUBDOMAIN_OCTAVE_ #define _SUBDOMAIN_OCTAVE_ #include "subdomain_rep.h" #include "dolfin_compat.h" class subdomain : public octave_base_value { public: subdomain () : octave_base_value () {} subdomain (octave_fcn_handle const & fh, bool on_boundary) : octave_base_value (), rep (new subdomain_rep (fh, on_boundary)) {} ~subdomain (void) {} void print (std::ostream& os, bool pr_as_read_syntax = false) const { os << "SubDomain" << std::endl; } bool is_defined (void) const { return true; } SHARED_PTR <subdomain_rep const> const & get_psd (void) const { return rep; } private: SHARED_PTR <subdomain_rep const> rep; DECLARE_OCTAVE_ALLOCATOR; DECLARE_OV_TYPEID_FUNCTIONS_AND_DATA; }; static bool subdomain_type_loaded = false; DEFINE_OCTAVE_ALLOCATOR (subdomain); DEFINE_OV_TYPEID_FUNCTIONS_AND_DATA (subdomain, "subdomain", "subdomain"); #endif