Line data Source code
1 : /*
2 : * Copyright (c) 2010-2015: G-CSC, Goethe University Frankfurt
3 : * Author: Andreas Vogel
4 : *
5 : * This file is part of UG4.
6 : *
7 : * UG4 is free software: you can redistribute it and/or modify it under the
8 : * terms of the GNU Lesser General Public License version 3 (as published by the
9 : * Free Software Foundation) with the following additional attribution
10 : * requirements (according to LGPL/GPL v3 §7):
11 : *
12 : * (1) The following notice must be displayed in the Appropriate Legal Notices
13 : * of covered and combined works: "Based on UG4 (www.ug4.org/license)".
14 : *
15 : * (2) The following notice must be displayed at a prominent place in the
16 : * terminal output of covered works: "Based on UG4 (www.ug4.org/license)".
17 : *
18 : * (3) The following bibliography is recommended for citation and must be
19 : * preserved in all covered files:
20 : * "Reiter, S., Vogel, A., Heppner, I., Rupp, M., and Wittum, G. A massively
21 : * parallel geometric multigrid solver on hierarchically distributed grids.
22 : * Computing and visualization in science 16, 4 (2013), 151-164"
23 : * "Vogel, A., Reiter, S., Rupp, M., Nägel, A., and Wittum, G. UG4 -- a novel
24 : * flexible software system for simulating pde based models on high performance
25 : * computers. Computing and visualization in science 16, 4 (2013), 165-179"
26 : *
27 : * This program is distributed in the hope that it will be useful,
28 : * but WITHOUT ANY WARRANTY; without even the implied warranty of
29 : * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30 : * GNU Lesser General Public License for more details.
31 : */
32 :
33 : #ifndef __H__UG__LIB_DISC__ASSEMBLE__
34 : #define __H__UG__LIB_DISC__ASSEMBLE__
35 :
36 : #include "lib_grid/tools/selector_grid.h"
37 : #include "lib_grid/tools/grid_level.h"
38 : #include "lib_disc/spatial_disc/ass_tuner.h"
39 :
40 : namespace ug{
41 :
42 : //predeclaration
43 : template <typename TAlgebra>
44 : class IConstraint;
45 :
46 :
47 : /**
48 : * \brief Assemblings.
49 : *
50 : * Interfaces to assemble problems.
51 : *
52 : * \defgroup lib_disc_assemble Assembling
53 : * \ingroup lib_discretization
54 : */
55 :
56 : /// \ingroup lib_disc_assemble
57 : /// @{
58 :
59 : //////////////////////////////////////////////////////////////////////////
60 : /// Interface providing Jacobian and Defect of a discretization.
61 : /**
62 : * Interface to generate Jacobi-Matrices and Defect-Vectors for a nonlinear
63 : * problem and to compute Matrix and Right-Hand Side for a linear problem.
64 : *
65 : * The Interface can be used directly to compute Jacobian and Defect
66 : * (resp. Mass-Stiffness-matrix_type and right-hand side) in case of nonlinear
67 : * (resp. linear) problem. Furthermore, the interface is used in the
68 : * NewtonSolver and MultiGridSolver. The NewtonSolver uses the functions
69 : * assemble_defect+assemble_jacobian or assemble_jacobian_defect.
70 : * The MultiGridSolver uses the function assemble_jacobian in order to
71 : * assemble coarse grid Matrices.
72 : *
73 : * To give an idea of the usage, an implementation should obey this rules:
74 : *
75 : * - Time dependent nonlinear problem \f$ \partial_t u + A(u) = f \f$.
76 : * Using a l-step time stepping, the defect will be
77 : * \f[
78 : * d(u^k) = \sum_{i=0}^{l-1} s_{m,i} M(u^{k-i})
79 : * + s_{a,i} \{A(u^{k-i}) - f\}
80 : * \f]
81 : * and assemble_linear will return false.
82 : *
83 : * - Time dependent linear problem \f$ \partial_t u + A*u = f \f$.
84 : * Using a l-step time stepping, the defect will be
85 : * \f[
86 : * d(u^k) = \sum_{i=0}^{l-1} s_{m,i} M*u^{k-i} + s_{a,i} \{A*u^{k-i} - f \}
87 : * J(u^k) = s_{m,0}*M + s_{a,0} A
88 : * \f]
89 : * and assemble_linear will compute \f$ A = s_{m,0}*M + s_{a,0} A \f$ and
90 : * \f$ b = - \{ \sum_{i=1}^{l-1} s_{m,i} M*u^{k-i}
91 : * + s_{a,i} \{A*u^{k-i} - f\} \} + s_{a,0} * f \f$
92 : *
93 : * - Stationary Non-linear Problem \f$ A(u) = f \f$. Then, the defect will be
94 : * \f[
95 : * d(u) = A(u) - f
96 : * \f]
97 : * and assemble_linear will return false.
98 : *
99 : * - Stationary linear Problem \f$ A*u = f \f$. Then, the defect will be
100 : * \f[
101 : * d(u) = A*u - f
102 : * J(u) = A
103 : * \f]
104 : * and assemble_linear will compute \f$ A = A \f$ and \f$ b = f \f$.
105 : *
106 : * \tparam TAlgebra Algebra type
107 : */
108 : template <typename TAlgebra>
109 : class IAssemble
110 : {
111 : public:
112 : /// Algebra type
113 : typedef TAlgebra algebra_type;
114 :
115 : /// Type of algebra matrix
116 : typedef typename TAlgebra::matrix_type matrix_type;
117 :
118 : /// Type of algebra vector
119 : typedef typename TAlgebra::vector_type vector_type;
120 :
121 : public:
122 : /// assembles Jacobian (or Approximation of Jacobian)
123 : /**
124 : * Assembles Jacobian at a given iterate u.
125 : *
126 : * \param[out] J Jacobian J(u) matrix to be filled
127 : * \param[in] u Current iterate
128 : * \param[in] gl Grid Level
129 : */
130 : virtual void assemble_jacobian(matrix_type& J, const vector_type& u, const GridLevel& gl) = 0;
131 0 : void assemble_jacobian(matrix_type& J, const vector_type& u)
132 0 : {assemble_jacobian(J,u,GridLevel());}
133 :
134 : /// assembles Defect
135 : /**
136 : * Assembles Defect at a given Solution u.
137 : *
138 : * \param[out] d Defect d(u) to be filled
139 : * \param[in] u Current iterate
140 : * \param[in] gl Grid Level
141 : */
142 : virtual void assemble_defect(vector_type& d, const vector_type& u, const GridLevel& gl) = 0;
143 0 : void assemble_defect(vector_type& d, const vector_type& u)
144 0 : {assemble_defect(d,u, GridLevel());}
145 :
146 : /// Assembles Matrix and Right-Hand-Side for a linear problem
147 : /**
148 : * Assembles matrix_type and Right-Hand-Side for a linear problem
149 : *
150 : * \param[out] A Mass-/Stiffness- Matrix
151 : * \param[out] b Right-Hand-Side
152 : * \param[in] gl Grid Level
153 : */
154 : virtual void assemble_linear(matrix_type& A, vector_type& b, const GridLevel& gl) = 0;
155 0 : void assemble_linear(matrix_type& A, vector_type& b)
156 0 : {assemble_linear(A,b, GridLevel());}
157 :
158 : /// assembles rhs
159 : virtual void assemble_rhs(vector_type& rhs, const vector_type& u, const GridLevel& gl) = 0;
160 0 : virtual void assemble_rhs(vector_type& rhs, const vector_type& u)
161 0 : {assemble_rhs(rhs, u, GridLevel());}
162 :
163 : /// Assembles Right-Hand-Side for a linear problem
164 : /**
165 : * Assembles Right-Hand-Side for a linear problem
166 : *
167 : * \param[out] b Right-Hand-Side
168 : * \param[in] gl Grid Level
169 : */
170 : virtual void assemble_rhs(vector_type& b, const GridLevel& gl) = 0;
171 0 : void assemble_rhs(vector_type& b)
172 0 : {assemble_rhs(b, GridLevel());}
173 :
174 : /// sets dirichlet values in solution vector
175 : /**
176 : * Sets dirichlet values of the NumericalSolution u when components
177 : * are dirichlet
178 : *
179 : * \param[out] u Numerical Solution
180 : * \param[in] gl Grid Level
181 : */
182 : virtual void adjust_solution(vector_type& u, const GridLevel& gl) = 0;
183 0 : void adjust_solution(vector_type& u)
184 0 : {adjust_solution(u, GridLevel());}
185 :
186 : /// assembles mass matrix
187 0 : virtual void assemble_mass_matrix(matrix_type& M, const vector_type& u, const GridLevel& gl)
188 0 : {UG_THROW("IAssemble: assemble_mass_matrix not implemented.");}
189 0 : void assemble_mass_matrix(matrix_type& M, const vector_type& u)
190 0 : {assemble_mass_matrix(M,u,GridLevel());}
191 :
192 : /// assembles stiffness matrix
193 0 : virtual void assemble_stiffness_matrix(matrix_type& A, const vector_type& u, const GridLevel& gl)
194 0 : {UG_THROW("IAssemble: assemble_stiffness_matrix not implemented.");}
195 0 : void assemble_stiffness_matrix(matrix_type& A, const vector_type& u)
196 0 : {assemble_stiffness_matrix(A,u,GridLevel());}
197 :
198 : /// \{
199 : virtual SmartPtr<AssemblingTuner<TAlgebra> > ass_tuner() = 0;
200 : virtual ConstSmartPtr<AssemblingTuner<TAlgebra> > ass_tuner() const = 0;
201 : /// \}
202 :
203 : /// returns the number of constraints
204 : virtual size_t num_constraints() const = 0;
205 :
206 : /// returns the i'th constraint
207 : virtual SmartPtr<IConstraint<TAlgebra> > constraint(size_t i) = 0;
208 :
209 : /// Virtual Destructor
210 : virtual ~IAssemble(){};
211 :
212 : };
213 :
214 : /// @}
215 :
216 : }; // name space ug
217 :
218 : #endif /* __H__UG__LIB_DISC__ASSEMBLE__ */
|