Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
source
meshworker
mesh_worker.cc
1
// ---------------------------------------------------------------------
2
//
3
// Copyright (C) 2006 - 2018 by the deal.II authors
4
//
5
// This file is part of the deal.II library.
6
//
7
// The deal.II library is free software; you can use it, redistribute
8
// it, and/or modify it under the terms of the GNU Lesser General
9
// Public License as published by the Free Software Foundation; either
10
// version 2.1 of the License, or (at your option) any later version.
11
// The full text of the license can be found in the file LICENSE.md at
12
// the top level directory of deal.II.
13
//
14
// ---------------------------------------------------------------------
15
16
17
#include <deal.II/lac/block_indices.h>
18
19
#include <deal.II/meshworker/local_integrator.h>
20
#include <deal.II/meshworker/local_results.h>
21
22
DEAL_II_NAMESPACE_OPEN
23
24
namespace
MeshWorker
25
{
26
template
<
typename
number>
27
void
28
LocalResults<number>::reinit
(
const
BlockIndices
&bi)
29
{
30
for
(
unsigned
int
i = 0; i < J.size(); ++i)
31
J[i] = 0.;
32
for
(
unsigned
int
i = 0; i < R.size(); ++i)
33
R[i].reinit(bi);
34
for
(
unsigned
int
i = 0; i < M1.size(); ++i)
35
M1[i].matrix.reinit(bi.
block_size
(M1[i].row),
36
bi.
block_size
(M1[i].column));
37
for
(
unsigned
int
i = 0; i < M2.size(); ++i)
38
M2[i].matrix.reinit(bi.
block_size
(M2[i].row),
39
bi.
block_size
(M2[i].column));
40
quadrature_data.reset_values();
41
}
42
43
44
template
<
typename
number>
45
std::size_t
46
LocalResults<number>::memory_consumption
()
const
47
{
48
std::size_t mem =
sizeof
(*this) +
MemoryConsumption::memory_consumption
(J) +
49
MemoryConsumption::memory_consumption
(R) +
50
MemoryConsumption::memory_consumption
(M1) +
51
MemoryConsumption::memory_consumption
(M2) +
52
MemoryConsumption::memory_consumption
(quadrature_data);
53
return
mem;
54
}
55
56
57
template
class
LocalResults<float>
;
58
template
class
LocalResults<double>
;
59
60
template
<
int
dim,
int
spacedim,
typename
number>
61
LocalIntegrator<dim, spacedim, number>::LocalIntegrator
()
62
: use_cell(true)
63
, use_boundary(true)
64
, use_face(true)
65
{}
66
67
68
template
<
int
dim,
int
spacedim,
typename
number>
69
LocalIntegrator<dim, spacedim, number>::LocalIntegrator
(
bool
c,
70
bool
b,
71
bool
f)
72
: use_cell(c)
73
, use_boundary(b)
74
, use_face(f)
75
{}
76
77
78
79
template
<
int
dim,
int
spacedim,
typename
number>
80
void
81
LocalIntegrator<dim, spacedim, number>::cell
(
82
DoFInfo<dim, spacedim, number>
&,
83
IntegrationInfo<dim, spacedim>
&)
const
84
{
85
Assert
(
false
, ExcPureFunction());
86
}
87
88
89
template
<
int
dim,
int
spacedim,
typename
number>
90
void
91
LocalIntegrator<dim, spacedim, number>::boundary
(
92
DoFInfo<dim, spacedim, number>
&,
93
IntegrationInfo<dim, spacedim>
&)
const
94
{
95
Assert
(
false
, ExcPureFunction());
96
}
97
98
99
template
<
int
dim,
int
spacedim,
typename
number>
100
void
101
LocalIntegrator<dim, spacedim, number>::face
(
102
DoFInfo<dim, spacedim, number>
&,
103
DoFInfo<dim, spacedim, number>
&,
104
IntegrationInfo<dim, spacedim>
&,
105
IntegrationInfo<dim, spacedim>
&)
const
106
{
107
Assert
(
false
, ExcPureFunction());
108
}
109
110
111
template
class
LocalIntegrator<1, 1, float>
;
112
template
class
LocalIntegrator<1, 1, double>
;
113
template
class
LocalIntegrator<1, 2, float>
;
114
template
class
LocalIntegrator<1, 2, double>
;
115
template
class
LocalIntegrator<1, 3, float>
;
116
template
class
LocalIntegrator<1, 3, double>
;
117
template
class
LocalIntegrator<2, 2, float>
;
118
template
class
LocalIntegrator<2, 2, double>
;
119
template
class
LocalIntegrator<2, 3, float>
;
120
template
class
LocalIntegrator<2, 3, double>
;
121
template
class
LocalIntegrator<3, 3, float>
;
122
template
class
LocalIntegrator<3, 3, double>
;
123
}
// namespace MeshWorker
124
125
126
DEAL_II_NAMESPACE_CLOSE
MeshWorker::LocalIntegrator::boundary
virtual void boundary(DoFInfo< dim, spacedim, number > &dinfo, IntegrationInfo< dim, spacedim > &info) const
Definition:
mesh_worker.cc:91
MeshWorker
Definition:
assemble_flags.h:30
BlockIndices::block_size
size_type block_size(const unsigned int i) const
Definition:
block_indices.h:376
MeshWorker::IntegrationInfo
Definition:
integration_info.h:78
MeshWorker::LocalIntegrator::cell
virtual void cell(DoFInfo< dim, spacedim, number > &dinfo, IntegrationInfo< dim, spacedim > &info) const
Definition:
mesh_worker.cc:81
MeshWorker::LocalIntegrator
Definition:
local_integrator.h:55
MeshWorker::LocalResults
Definition:
local_results.h:211
MemoryConsumption::memory_consumption
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
Definition:
memory_consumption.h:274
MeshWorker::LocalIntegrator::face
virtual void face(DoFInfo< dim, spacedim, number > &dinfo1, DoFInfo< dim, spacedim, number > &dinfo2, IntegrationInfo< dim, spacedim > &info1, IntegrationInfo< dim, spacedim > &info2) const
Definition:
mesh_worker.cc:101
Assert
#define Assert(cond, exc)
Definition:
exceptions.h:1407
MeshWorker::LocalResults::reinit
void reinit(const BlockIndices &local_sizes)
Definition:
mesh_worker.cc:28
MeshWorker::LocalResults::memory_consumption
std::size_t memory_consumption() const
Definition:
mesh_worker.cc:46
MeshWorker::DoFInfo
Definition:
dof_info.h:74
MeshWorker::LocalIntegrator::LocalIntegrator
LocalIntegrator()
Definition:
mesh_worker.cc:61
BlockIndices
Definition:
block_indices.h:60
Generated by
1.8.16