/usr/include/dolfin/fem/SystemAssembler.h is in libdolfin1.0-dev 1.0.0-1.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 | // Copyright (C) 2008-2009 Kent-Andre Mardal and Garth N. Wells
//
// This file is part of DOLFIN.
//
// DOLFIN is free software: you can redistribute it and/or modify
// it under the terms of the GNU Lesser General Public License as published by
// the Free Software Foundation, either version 3 of the License, or
// (at your option) any later version.
//
// DOLFIN 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 Lesser General Public License for more details.
//
// You should have received a copy of the GNU Lesser General Public License
// along with DOLFIN. If not, see <http://www.gnu.org/licenses/>.
//
// Modified by Anders Logg 2008-2011
//
// First added: 2009-06-22
// Last changed: 2011-01-25
#ifndef __SYSTEM_ASSEMBLER_H
#define __SYSTEM_ASSEMBLER_H
#include <map>
#include <vector>
#include <dolfin/common/types.h>
#include "DirichletBC.h"
namespace dolfin
{
// Forward declarations
class GenericMatrix;
class GenericTensor;
class GenericVector;
class Form;
class Mesh;
class SubDomain;
class UFC;
class Cell;
class Facet;
class Function;
template<typename T> class MeshFunction;
/// This class provides implements an assembler for systems
/// of the form Ax = b. It differs from the default DOLFIN
/// assembler in that it assembles both A and b and the same
/// time (leading to better performance) and in that it applies
/// boundary conditions at the time of assembly.
class SystemAssembler
{
public:
/// Assemble system (A, b)
static void assemble(GenericMatrix& A,
GenericVector& b,
const Form& a,
const Form& L,
bool reset_sparsity=true,
bool add_values=false,
bool finalize_tensor=true);
/// Assemble system (A, b) and apply Dirichlet boundary condition
static void assemble(GenericMatrix& A,
GenericVector& b,
const Form& a,
const Form& L,
const DirichletBC& bc,
bool reset_sparsity=true,
bool add_values=true,
bool finalize_tensor=true);
/// Assemble system (A, b) and apply Dirichlet boundary conditions
static void assemble(GenericMatrix& A,
GenericVector& b,
const Form& a,
const Form& L,
const std::vector<const DirichletBC*>& bcs,
bool reset_sparsity=true,
bool add_values=false,
bool finalize_tensor=true);
/// Assemble system (A, b) and apply Dirichlet boundary conditions
static void assemble(GenericMatrix& A,
GenericVector& b,
const Form& a,
const Form& L,
const std::vector<const DirichletBC*>& bcs,
const MeshFunction<uint>* cell_domains,
const MeshFunction<uint>* exterior_facet_domains,
const MeshFunction<uint>* interior_facet_domains,
const GenericVector* x0,
bool reset_sparsity=true,
bool add_values=false,
bool finalize_tensor=true);
private:
class Scratch;
static void compute_tensor_on_one_interior_facet(const Form& a,
UFC& ufc,
const Cell& cell1,
const Cell& cell2,
const Facet& facet,
const MeshFunction<uint>* exterior_facet_domains);
static void cell_wise_assembly(GenericMatrix& A, GenericVector& b,
const Form& a, const Form& L,
UFC& A_ufc, UFC& b_ufc, Scratch& data,
const DirichletBC::Map& boundary_values,
const MeshFunction<uint>* cell_domains,
const MeshFunction<uint>* exterior_facet_domains);
static void facet_wise_assembly(GenericMatrix& A, GenericVector& b,
const Form& a, const Form& L,
UFC& A_ufc, UFC& b_ufc, Scratch& data,
const DirichletBC::Map& boundary_values,
const MeshFunction<uint>* cell_domains,
const MeshFunction<uint>* exterior_facet_domains,
const MeshFunction<uint>* interior_facet_domains);
static void assemble_interior_facet(GenericMatrix& A, GenericVector& b,
UFC& A_ufc, UFC& b_ufc,
const Form& a, const Form& L,
const Cell& cell0, const Cell& cell1,
const Facet& facet,
Scratch& data,
const DirichletBC::Map& boundary_values);
static void assemble_exterior_facet(GenericMatrix& A, GenericVector& b,
UFC& A_ufc, UFC& b_ufc,
const Form& a,
const Form& L,
const Cell& cell, const Facet& facet,
Scratch& data,
const DirichletBC::Map& boundary_values);
static void apply_bc(double* A, double* b,
const DirichletBC::Map& boundary_values,
const std::vector<const std::vector<uint>* >& global_dofs);
// Class to hold temporary data
class Scratch
{
public:
Scratch(const Form& a, const Form& L);
~Scratch();
void zero_cell();
std::vector<double> Ae;
std::vector<double> be;
};
};
}
#endif
|