DOLFIN
DOLFIN C++ interface
Toggle main menu visibility
Loading...
Searching...
No Matches
dolfin
adaptivity
adapt.h
1
// Copyright (C) 2010-2011 Anders Logg, Marie Rognes and Garth N. Wells
2
//
3
// This file is part of DOLFIN.
4
//
5
// DOLFIN is free software: you can redistribute it and/or modify
6
// it under the terms of the GNU Lesser General Public License as published by
7
// the Free Software Foundation, either version 3 of the License, or
8
// (at your option) any later version.
9
//
10
// DOLFIN is distributed in the hope that it will be useful,
11
// but WITHOUT ANY WARRANTY; without even the implied warranty of
12
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13
// GNU Lesser General Public License for more details.
14
//
15
// You should have received a copy of the GNU Lesser General Public License
16
// along with DOLFIN. If not, see <http://www.gnu.org/licenses/>.
17
//
18
// First added: 2010-02-10
19
// Last changed: 2011-10-11
20
//
21
// This file defines free functions for refinement/adaption of meshes,
22
// function spaces, functions etc.
23
24
#ifndef __DOLFIN_ADAPT_H
25
#define __DOLFIN_ADAPT_H
26
27
#include <vector>
28
29
namespace
dolfin
30
{
31
32
// Forward declarations
33
class
DirichletBC
;
34
class
ErrorControl
;
35
class
Form
;
36
class
Function
;
37
class
FunctionSpace
;
38
class
GenericFunction
;
39
class
LinearVariationalProblem
;
40
class
Mesh
;
41
template
<
typename
T>
class
MeshFunction
;
42
class
NonlinearVariationalProblem
;
43
44
//--- Refinement of meshes ---
45
50
std::shared_ptr<Mesh>
adapt
(
const
Mesh
& mesh);
51
57
std::shared_ptr<Mesh>
adapt
(
const
Mesh
& mesh,
58
const
MeshFunction<bool>
& cell_markers);
59
60
//--- Refinement of function spaces ---
61
65
std::shared_ptr<FunctionSpace>
adapt
(
const
FunctionSpace
& space);
66
71
std::shared_ptr<FunctionSpace>
adapt
(
const
FunctionSpace
& space,
72
const
MeshFunction<bool>
& cell_markers);
73
78
std::shared_ptr<FunctionSpace>
adapt
(
const
FunctionSpace
& space,
79
std::shared_ptr<const Mesh> adapted_mesh);
80
81
//--- Refinement of functions ---
82
96
std::shared_ptr<Function>
adapt
(
const
Function
& function,
97
std::shared_ptr<const Mesh> adapted_mesh,
98
bool
interpolate=
true
);
99
109
std::shared_ptr<GenericFunction>
110
adapt
(std::shared_ptr<const GenericFunction> function,
111
std::shared_ptr<const Mesh> adapted_mesh);
112
114
std::shared_ptr<MeshFunction<std::size_t>>
115
adapt
(
const
MeshFunction<std::size_t>
& mesh_function,
116
std::shared_ptr<const Mesh> adapted_mesh);
117
118
//--- Refinement of boundary conditions ---
119
121
std::shared_ptr<DirichletBC>
adapt
(
const
DirichletBC
& bc,
122
std::shared_ptr<const Mesh> adapted_mesh,
123
const
FunctionSpace
& S);
124
126
void
adapt_markers
(std::vector<std::size_t>& refined_markers,
127
const
Mesh
& adapted_mesh,
128
const
std::vector<std::size_t>& markers,
129
const
Mesh
& mesh);
130
131
//--- Refinement of forms ---
132
146
std::shared_ptr<Form>
adapt
(
const
Form
& form,
147
std::shared_ptr<const Mesh> adapted_mesh,
148
bool
adapt_coefficients=
true
);
149
150
//--- Refinement of variational problems ---
151
153
std::shared_ptr<LinearVariationalProblem>
154
adapt
(
const
LinearVariationalProblem
& problem,
155
std::shared_ptr<const Mesh> adapted_mesh);
156
158
std::shared_ptr<NonlinearVariationalProblem>
159
adapt
(
const
NonlinearVariationalProblem
& problem,
160
std::shared_ptr<const Mesh> adapted_mesh);
161
175
std::shared_ptr<ErrorControl>
adapt
(
const
ErrorControl
& ec,
176
std::shared_ptr<const Mesh> adapted_mesh,
177
bool
adapt_coefficients=
true
);
178
179
180
}
181
182
#endif
dolfin::DirichletBC
Interface for setting (strong) Dirichlet boundary conditions.
Definition
DirichletBC.h:125
dolfin::ErrorControl
(Goal-oriented) Error Control class.
Definition
ErrorControl.h:51
dolfin::Form
Base class for UFC code generated by FFC for DOLFIN with option -l.
Definition
Form.h:86
dolfin::FunctionSpace
Definition
FunctionSpace.h:54
dolfin::Function
Definition
Function.h:66
dolfin::GenericFunction
Definition
GenericFunction.h:54
dolfin::LinearVariationalProblem
Definition
LinearVariationalProblem.h:43
dolfin::MeshFunction
Definition
MeshFunction.h:59
dolfin::Mesh
Definition
Mesh.h:83
dolfin::NonlinearVariationalProblem
Definition
NonlinearVariationalProblem.h:47
dolfin
Definition
adapt.h:30
dolfin::adapt
std::shared_ptr< Mesh > adapt(const Mesh &mesh)
Definition
adapt.cpp:56
dolfin::adapt_markers
void adapt_markers(std::vector< std::size_t > &refined_markers, const Mesh &adapted_mesh, const std::vector< std::size_t > &markers, const Mesh &mesh)
Helper function for refinement of boundary conditions.
Definition
adapt.cpp:590
Generated by
1.17.0