LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
fix_dof.h
1#ifndef _LF_FIXDOF_H
2#define _LF_FIXDOF_H
3/***************************************************************************
4 * LehrFEM++ - A simple C++ finite element libray for teaching
5 * Developed from 2018 at the Seminar of Applied Mathematics of ETH Zurich,
6 * lead developers Dr. R. Casagrande and Prof. R. Hiptmair
7 ***************************************************************************/
8
17#include "coomatrix.h"
18
19namespace lf::assemble {
20
86template <typename SCALAR, typename SELECTOR, typename RHSVECTOR>
87void FixFlaggedSolutionComponents(SELECTOR &&selectvals, COOMatrix<SCALAR> &A,
88 RHSVECTOR &b) {
89 const lf::assemble::size_type N(A.cols());
90 LF_ASSERT_MSG(A.rows() == N, "Matrix must be square!");
91 LF_ASSERT_MSG(N == b.size(),
92 "Mismatch N = " << N << " <-> b.size() = " << b.size());
93 {
94 // Multiply sparse matrix with the vector of fixed components and subtract
95 // result from right hand side
96 // To skip irrelevant components of fixed_vec rely on a temporary vector
97 Eigen::Matrix<SCALAR, Eigen::Dynamic, 1> tmp_vec(N);
98 for (lf::assemble::gdof_idx_t k = 0; k < N; ++k) {
99 const auto selval{selectvals(k)};
100 if (selval.first) {
101 tmp_vec[k] = selval.second;
102 } else {
103 tmp_vec[k] = SCALAR();
104 }
105 }
106 A.MatVecMult(-1.0, tmp_vec, b);
107 }
108 // Set vector components of right-hand-side vector for prescribed values
109 for (lf::assemble::gdof_idx_t k = 0; k < N; ++k) {
110 const auto selval{selectvals(k)};
111 if (selval.first) {
112 b[k] = selval.second;
113 }
114 }
115 // Set rows and columns of the sparse matrix corresponding to the fixed
116 // solution components to zero
117 A.setZero([&selectvals](gdof_idx_t i, gdof_idx_t j) {
118 return (selectvals(i).first || selectvals(j).first);
119 });
120 // Old implementation, demonstrating what is going on
121 // lf::assemble::COOMatrix<double>::TripletVec::iterator new_last =
122 // std::remove_if(
123 // A.triplets().begin(), A.triplets().end(),
124 // [&fixed_comp_flags](
125 // typename lf::assemble::COOMatrix<double>::Triplet &triplet) {
126 // return (fixed_comp_flags[triplet.row()] ||
127 // fixed_comp_flags[triplet.col()]);
128 // });
129 // // Adjust size of triplet vector
130 // A.triplets().erase(new_last, A.triplets().end());
131 // Add Unit diagonal entries corrresponding to fixed components
132 for (lf::assemble::gdof_idx_t dofnum = 0; dofnum < N; ++dofnum) {
133 const auto selval{selectvals(dofnum)};
134 if (selval.first) {
135 A.AddToEntry(dofnum, dofnum, 1.0);
136 }
137 }
138}
139
183template <typename SCALAR, typename SELECTOR, typename RHSVECTOR>
184void FixFlaggedSolutionCompAlt(SELECTOR &&selectvals, COOMatrix<SCALAR> &A,
185 RHSVECTOR &b) {
186 const lf::assemble::size_type N(A.cols());
187 LF_ASSERT_MSG(A.rows() == N, "Matrix must be square!");
188 LF_ASSERT_MSG(N == b.size(),
189 "Mismatch N = " << N << " <-> b.size() = " << b.size());
190
191 // Set vector components of right-hand-side vector for prescribed values
192 for (lf::assemble::gdof_idx_t k = 0; k < N; ++k) {
193 const auto selval{selectvals(k)};
194 if (selval.first) {
195 b[k] = selval.second;
196 }
197 }
198 // Set rows and columns of the sparse matrix corresponding to the fixed
199 // solution components to zero
200 A.setZero([&selectvals](gdof_idx_t i, gdof_idx_t /*unused*/) {
201 return (selectvals(i).first);
202 });
203 // Old implementation showing the algorithm:
204 // lf::assemble::COOMatrix<double>::TripletVec::iterator new_last =
205 // std::remove_if(
206 // A.triplets().begin(), A.triplets().end(),
207 // [&fixed_comp_flags](
208 // typename lf::assemble::COOMatrix<double>::Triplet &triplet) {
209 // return (fixed_comp_flags[triplet.row()]);
210 // });
211 // // Adjust size of triplet vector
212 // A.triplets().erase(new_last, A.triplets().end());
213 // Add Unit diagonal entries corrresponding to fixed components
214 for (lf::assemble::gdof_idx_t dofnum = 0; dofnum < N; ++dofnum) {
215 if (selectvals(dofnum).first) {
216 A.AddToEntry(dofnum, dofnum, 1.0);
217 }
218 }
219}
220
224template <typename SCALAR>
226 std::vector<std::pair<lf::assemble::gdof_idx_t, SCALAR>>;
227
250template <typename SCALAR, typename RHSVECTOR>
252 const fixed_components_t<SCALAR> &fixed_components, COOMatrix<SCALAR> &A,
253 RHSVECTOR &b) {
254 const lf::assemble::size_type N(A.cols());
255 LF_ASSERT_MSG(A.rows() == N, "Matrix must be square!");
256 LF_ASSERT_MSG(N == b.size(),
257 "Mismatch N = " << N << " <-> b.size() = " << b.size());
258 // Set up a vector containing the negative fixed solution components and a
259 // flag vector, whose entry is `true` if the corresponding
260 // vector component is meant to be fixed
261 Eigen::Matrix<SCALAR, Eigen::Dynamic, 1> fixed_vec(b.size());
262 fixed_vec.setZero();
263 std::vector<bool> fixed_comp_flags(N, false);
264 for (typename fixed_components_t<SCALAR>::const_reference idx_val_pair :
265 fixed_components) {
266 LF_ASSERT_MSG(idx_val_pair.first < N,
267 "Index " << idx_val_pair.first << " >= N = " << N);
268 fixed_vec[idx_val_pair.first] += idx_val_pair.second;
269 fixed_comp_flags[idx_val_pair.first] = true;
270 }
271 // FixFlaggedSolutionComponents<SCALAR>(fixed_comp_flags, fixed_vec, A, b);
272 FixFlaggedSolutionCompAlt<SCALAR>(
273 [&fixed_comp_flags,
274 &fixed_vec](lf::assemble::gdof_idx_t i) -> std::pair<bool, double> {
275 LF_ASSERT_MSG((i < fixed_comp_flags.size()) && (i < fixed_vec.size()),
276 "Illegal index " << i);
277 return std::make_pair(fixed_comp_flags[i], fixed_vec[i]);
278 },
279 A, b);
280}
281
282} // namespace lf::assemble
283
284#endif
A temporary data structure for matrix in COO format.
Definition: coomatrix.h:52
void setZero()
Erase all entries of the matrix.
Definition: coomatrix.h:98
Index cols() const
return number of column
Definition: coomatrix.h:72
Eigen::Matrix< SCALAR, Eigen::Dynamic, 1 > MatVecMult(SCALAR alpha, const VECTOR &vec) const
Computes the product of a (scaled) vector with the matrix in COO format.
Definition: coomatrix.h:236
Index rows() const
return number of rows
Definition: coomatrix.h:70
void AddToEntry(gdof_idx_t i, gdof_idx_t j, SCALAR increment)
Add a value to the specified entry.
Definition: coomatrix.h:87
D.o.f. index mapping and assembly facilities.
Definition: assemble.h:30
lf::base::size_type size_type
void FixSolutionComponentsLse(const fixed_components_t< SCALAR > &fixed_components, COOMatrix< SCALAR > &A, RHSVECTOR &b)
manipulate a square linear system of equations with a coefficient matrix in COO format so that some s...
Definition: fix_dof.h:251
void FixFlaggedSolutionComponents(SELECTOR &&selectvals, COOMatrix< SCALAR > &A, RHSVECTOR &b)
enforce prescribed solution components
Definition: fix_dof.h:87
Eigen::Index gdof_idx_t
void FixFlaggedSolutionCompAlt(SELECTOR &&selectvals, COOMatrix< SCALAR > &A, RHSVECTOR &b)
Setting unknowns of a sparse linear system of equations to fixed values.
Definition: fix_dof.h:184
std::vector< std::pair< lf::assemble::gdof_idx_t, SCALAR > > fixed_components_t
Information about fixed solution components.
Definition: fix_dof.h:226