1 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2 // vi: set et ts=4 sw=2 sts=2:
3 #include <config.h>
4 
5 #include <vector>
6 #include <cmath>
7 
8 #include <dune/common/bitsetvector.hh>
9 #include <dune/common/indices.hh>
10 #include <dune/common/transpose.hh>
11 
12 #include <dune/geometry/quadraturerules.hh>
13 
14 #include <dune/grid/yaspgrid.hh>
15 #include <dune/grid/io/file/vtk/subsamplingvtkwriter.hh>
16 
17 #include <dune/istl/matrix.hh>
18 #include <dune/istl/bcrsmatrix.hh>
19 #include <dune/istl/matrixindexset.hh>
20 #include <dune/istl/preconditioners.hh>
21 #include <dune/istl/solvers.hh>
22 
23 #include <dune/functions/functionspacebases/interpolate.hh>
24 #include <dune/functions/functionspacebases/lagrangebasis.hh>
25 #include <dune/functions/functionspacebases/raviartthomasbasis.hh>
26 #include <dune/functions/backends/istlvectorbackend.hh>
27 #include <dune/functions/functionspacebases/compositebasis.hh>
28 #include <dune/functions/functionspacebases/subspacebasis.hh>
29 #include <dune/functions/gridfunctions/discreteglobalbasisfunction.hh>
30 #include <dune/functions/gridfunctions/gridviewfunction.hh>
31 
32 #define DIM2 // Use a two-dimensional test, otherwise three-dimensional
33 
34 using namespace Dune;
35 
36 // Compute the stiffness matrix for a single element
37 template <class LocalView, class MatrixType>
getLocalMatrix(const LocalView & localView,MatrixType & elementMatrix)38 void getLocalMatrix(const LocalView& localView,
39                     MatrixType& elementMatrix)
40 {
41   // Get the grid element from the local FE basis view (use test space)
42   typedef typename LocalView::Element Element;
43   const Element& element = localView.element();
44 
45   const int dim = Element::dimension;
46   auto geometry = element.geometry();
47 
48   // Set all matrix entries to zero
49   elementMatrix.setSize(localView.size(), localView.size());
50   elementMatrix = 0;      // fills the entire matrix with zeroes
51 
52   // Get set of shape functions for this element
53   using namespace Dune::Indices;
54   const auto& fluxLocalFiniteElement     = localView.tree().child(_0).finiteElement();
55   const auto& pressureLocalFiniteElement = localView.tree().child(_1).finiteElement();
56 
57   // Get a quadrature rule
58   int order = std::max(2*(dim*fluxLocalFiniteElement.localBasis().order()-1), 2*(dim*pressureLocalFiniteElement.localBasis().order()));
59   const auto& quad = QuadratureRules<double, dim>::rule(element.type(), order);
60 
61   // Loop over all quadrature points
62   for (const auto& quadPoint : quad)
63   {
64     // Position of the current quadrature point in the reference element
65     const auto quadPos = quadPoint.position();
66 
67     // The transposed inverse Jacobian of the map from the reference element to the element
68     const auto jacInvTrans = geometry.jacobianInverseTransposed(quadPos);
69 
70     // The multiplicative factor in the integral transformation formula
71     const auto integrationElement = geometry.integrationElement(quadPos);
72 
73     ///////////////////////////////////////////////////////////////////////////
74     // Shape functions - flux
75     ///////////////////////////////////////////////////////////////////////////
76 
77     // Values of the flux shape functions on the current element
78     std::vector<FieldVector<double,dim> > fluxValues(fluxLocalFiniteElement.size());
79     fluxLocalFiniteElement.localBasis().evaluateFunction(quadPos, fluxValues);
80 
81     // Gradients of the flux shape function gradients on the reference element
82     std::vector<FieldMatrix<double,dim,dim> > fluxReferenceJacobians(fluxLocalFiniteElement.size());
83     fluxLocalFiniteElement.localBasis().evaluateJacobian(quadPos, fluxReferenceJacobians);
84 
85     // Helper function to compute the trace of a matrix
86     auto trace = [](const auto& matrix) {
87       double r=0;
88       for (size_t j=0; j<matrix.N(); j++)
89         r += matrix[j][j];
90       return r;
91     };
92 
93     // Domain transformation of Jacobians and computation of div = trace(Jacobian)
94     std::vector<double> fluxDivergence(fluxValues.size(), 0.0);
95     for (size_t i=0; i<fluxReferenceJacobians.size(); i++)
96       fluxDivergence[i] = trace(fluxReferenceJacobians[i] * transpose(jacInvTrans));
97 
98     ///////////////////////////////////////////////////////////////////////////
99     // Shape functions - pressure
100     ///////////////////////////////////////////////////////////////////////////
101 
102     // Values of the pressure shape functions
103     std::vector<FieldVector<double,1> > pressureValues(pressureLocalFiniteElement.size());
104     pressureLocalFiniteElement.localBasis().evaluateFunction(quadPos, pressureValues);
105 
106     /////////////////////////////
107     // Flux--flux coupling
108     /////////////////////////////
109 
110     for (size_t i=0; i<fluxLocalFiniteElement.size(); i++)
111     {
112       size_t row = localView.tree().child(_0).localIndex(i);
113       for (size_t j=0; j<fluxLocalFiniteElement.size(); j++)
114       {
115         size_t col = localView.tree().child(_0).localIndex(j);
116         elementMatrix[row][col] += (fluxValues[i] * fluxValues[j]) * quadPoint.weight() * integrationElement;
117       }
118     }
119 
120     /////////////////////////////
121     // Flux--pressure coupling
122     /////////////////////////////
123 
124     for (size_t i=0; i<fluxLocalFiniteElement.size(); i++)
125     {
126       size_t fluxIndex     = localView.tree().child(_0).localIndex(i);
127       for (size_t j=0; j<pressureLocalFiniteElement.size(); j++)
128       {
129         size_t pressureIndex = localView.tree().child(_1).localIndex(j);
130 
131         // Pre-compute matrix contribution
132         double tmp = - (fluxDivergence[i] * pressureValues[j]) * quadPoint.weight() * integrationElement;
133 
134         elementMatrix[fluxIndex][pressureIndex] += tmp;
135         elementMatrix[pressureIndex][fluxIndex] += tmp;
136       }
137     }
138   }
139 }
140 
141 
142 // Compute the source term for a single element
143 template <class LocalView, class LocalVolumeTerm>
getVolumeTerm(const LocalView & localView,BlockVector<FieldVector<double,1>> & localRhs,LocalVolumeTerm && localVolumeTerm)144 void getVolumeTerm( const LocalView& localView,
145                     BlockVector<FieldVector<double,1> >& localRhs,
146                     LocalVolumeTerm&& localVolumeTerm)
147 {
148   // Get the grid element from the local FE basis view
149   typedef typename LocalView::Element Element;
150   const Element& element = localView.element();
151 
152   const int dim = Element::dimension;
153 
154   // Set all entries to zero
155   localRhs.resize(localView.size());
156   localRhs = 0;
157 
158   // Get set of shape functions for this element
159   using namespace Dune::Indices;
160   const auto& fluxLocalFiniteElement     = localView.tree().child(_0).finiteElement();
161   const auto& pressureLocalFiniteElement = localView.tree().child(_1).finiteElement();
162 
163   // A quadrature rule
164   int order = std::max(dim*fluxLocalFiniteElement.localBasis().order(), dim*pressureLocalFiniteElement.localBasis().order());
165   const auto& quad = QuadratureRules<double, dim>::rule(element.type(), order);
166 
167   // Loop over all quadrature points
168   for (const auto& quadPoint : quad)
169   {
170     // Position of the current quadrature point in the reference element
171     const auto& quadPos = quadPoint.position();
172 
173     // The multiplicative factor in the integral transformation formula
174     const double integrationElement = element.geometry().integrationElement(quadPos);
175 
176     double functionValue = localVolumeTerm(quadPos);
177 
178     // Evaluate all shape function values at this point
179     std::vector<FieldVector<double,1> > pressureValues;
180     pressureLocalFiniteElement.localBasis().evaluateFunction(quadPos, pressureValues);
181 
182     // Actually compute the vector entries
183     for (size_t j=0; j<pressureLocalFiniteElement.size(); j++)
184     {
185       size_t pressureIndex = localView.tree().child(_1).localIndex(j);
186       localRhs[pressureIndex] += - pressureValues[j] * functionValue * quadPoint.weight() * integrationElement;
187     }
188   }
189 }
190 
191 // Get the occupation pattern of the stiffness matrix
192 template <class Basis>
getOccupationPattern(const Basis & basis,std::array<std::array<MatrixIndexSet,2>,2> & nb)193 void getOccupationPattern(const Basis& basis,
194                           std::array<std::array<MatrixIndexSet,2>,2>& nb)
195 {
196   // Set sizes of the 2x2 submatrices
197   for (size_t i=0; i<2; ++i)
198     for (size_t j=0; j<2; j++)
199       nb[i][j].resize(basis.size({i}), basis.size({j}));
200 
201   // A view on the FE basis on a single element
202   auto localView = basis.localView();
203 
204   // Loop over all leaf elements
205   for(const auto& element : elements(basis.gridView()))
206   {
207     // Bind the local FE basis view to the current element
208     localView.bind(element);
209 
210     // Add element stiffness matrix onto global stiffness matrix
211     for (size_t i=0; i<localView.size(); i++)
212     {
213       // The global index of the i-th local degree of freedom of the element 'e'
214       auto row = localView.index(i);
215 
216       for (size_t j=0; j<localView.size(); ++j)
217       {
218         // The global index set of the j-th local degree of freedom of the element 'e'
219         auto col = localView.index(j);
220 
221         nb[row[0]][col[0]].add(row[1], col[1]);
222       }
223     }
224   }
225 }
226 
227 /** \brief Assemble the divergence stiffness matrix on the given grid view */
228 template <class Basis, class MatrixType>
assembleMixedPoissonMatrix(const Basis & basis,MatrixType & matrix)229 void assembleMixedPoissonMatrix(const Basis& basis,
230                                 MatrixType& matrix)
231 {
232   // Get the grid view from the finite element basis (use of test space arbitrary)
233   typedef typename Basis::GridView GridView;
234   GridView gridView = basis.gridView();
235 
236   // MatrixIndexSets store the occupation pattern of a sparse matrix.
237   // They are not particularly efficient, but simple to use.
238   std::array<std::array<MatrixIndexSet, 2>, 2> occupationPattern;
239   getOccupationPattern(basis, occupationPattern);
240   // ... and give it the occupation pattern we want.
241   matrix.setSize(2,2);
242   for (size_t i=0; i<2; i++)
243     for (size_t j=0; j<2; j++)
244       occupationPattern[i][j].exportIdx(matrix[i][j]);
245 
246   // Set all entries to zero
247   matrix = 0;
248 
249   // A view on the FE basis on a single element
250   auto localView = basis.localView();
251 
252   // A loop over all elements of the grid
253   for(const auto& element : elements(gridView))
254   {
255     // Bind the local FE basis view to the current element
256     localView.bind(element);
257 
258     // Now let's get the element stiffness matrix
259     // A dense matrix is used for the element stiffness matrix
260     Matrix<FieldMatrix<double,1,1> > elementMatrix;
261     getLocalMatrix(localView, elementMatrix);
262 
263     // Add element stiffness matrix onto the global stiffness matrix
264     for (size_t i=0; i<elementMatrix.N(); i++)
265     {
266       // The global index of the i-th local degree of freedom of the element 'e'
267       auto row = localView.index(i);
268 
269       for (size_t j=0; j<elementMatrix.M(); j++ )
270       {
271         // The global index of the j-th local degree of freedom of the element 'e'
272         auto col = localView.index(j);
273         matrix[row[0]][col[0]][row[1]][col[1]] += elementMatrix[i][j];
274       }
275     }
276   }
277 }
278 
279 /** \brief Assemble the divergence stiffness matrix on the given grid view */
280 template <class Basis, class VectorType, class VolumeTerm>
assembleMixedPoissonRhs(const Basis & basis,VectorType & rhs,VolumeTerm && volumeTerm)281 void assembleMixedPoissonRhs(const Basis& basis,
282                              VectorType& rhs,
283                              VolumeTerm&& volumeTerm)
284 {
285   // Get the grid view from the finite element basis (use of test space arbitrary)
286   typedef typename Basis::GridView GridView;
287   GridView gridView = basis.gridView();
288 
289   auto localVolumeTerm = localFunction(Functions::makeGridViewFunction(volumeTerm, gridView));
290 
291   // set rhs to correct length -- the total number of basis vectors in the basis
292   using Functions::istlVectorBackend;
293   istlVectorBackend(rhs).resize(basis);
294 
295   // Set all entries to zero
296   rhs = 0;
297 
298   // A view on the FE basis on a single element
299   auto localView = basis.localView();
300 
301   // A loop over all elements of the grid
302   for(const auto& element : elements(gridView))
303   {
304     // Bind the local FE basis view to the current element
305     localView.bind(element);
306 
307     // Now get the local contribution to the right-hand side vector
308     BlockVector<FieldVector<double,1> > localRhs;
309     localVolumeTerm.bind(element);
310     getVolumeTerm(localView, localRhs, localVolumeTerm);
311 
312     for (size_t i=0; i<localRhs.size(); i++)
313     {
314       // The global index of the i-th vertex of the element 'e'
315       auto row = localView.index(i);
316       rhs[row[0]][row[1]] += localRhs[i];
317     }
318   }
319 }
320 
main(int argc,char * argv[])321 int main (int argc, char *argv[])
322 {
323   // Set up MPI, if available
324   MPIHelper::instance(argc, argv);
325 
326   ///////////////////////////////////
327   //   Generate the grid
328   ///////////////////////////////////
329 
330 #ifdef DIM2
331   const int dim = 2;
332   std::array<int,dim> elements = {{50, 50}};
333 #else
334   const int dim = 3;
335   std::array<int,dim> elements = {{10, 10, 10}};
336 #endif
337   typedef YaspGrid<dim> GridType;
338   FieldVector<double,dim> l(1);
339   GridType grid(l,elements);
340 
341   typedef GridType::LeafGridView GridView;
342   GridView gridView = grid.leafGridView();
343 
344   /////////////////////////////////////////////////////////
345   //   Choose a finite element space
346   /////////////////////////////////////////////////////////
347 
348   using namespace Functions::BasisFactory;
349 
350   const int k = 0;
351 
352   auto basis = makeBasis(
353     gridView,
354     composite(
355       raviartThomas<k>(),
356       lagrange<k>()
357     ));
358 
359   using namespace Indices;
360   auto fluxBasis = Functions::subspaceBasis(basis, _0);
361   auto pressureBasis = Functions::subspaceBasis(basis, _1);
362 
363 
364   /////////////////////////////////////////////////////////
365   //   Stiffness matrix and right hand side vector
366   /////////////////////////////////////////////////////////
367 
368   typedef Matrix<BCRSMatrix<FieldMatrix<double,1,1> > > MatrixType;
369   typedef BlockVector<BlockVector<FieldVector<double,1> > > VectorType;
370 
371   using Functions::istlVectorBackend;
372 
373   MatrixType stiffnessMatrix;
374   VectorType rhs;
375 
376   /////////////////////////////////////////////////////////
377   //  Assemble the system
378   /////////////////////////////////////////////////////////
379 
380   using Domain = GridType::template Codim<0>::Geometry::GlobalCoordinate;
381 
382   assembleMixedPoissonMatrix(basis, stiffnessMatrix);
383 
384   auto rightHandSide = [] (const Domain& x) { return 10; };
385   assembleMixedPoissonRhs(basis, rhs, rightHandSide);
386 
387   auto pi = std::acos(-1.0);
388   auto topFluxBC   = [&pi] (const Domain& x) { return -0.05 * (1. - x[0]) * std::sin(2.*pi*x[0]); };
389   auto lowerFluxBC = [&pi] (const Domain& x) { return  0.05 * (1. - x[0]) * std::sin(2.*pi*x[0]); };
390 
391   using namespace Dune::Indices;
392   using BitVectorType = BlockVector<BlockVector<FieldVector<char,1> > >;
393 
394   BitVectorType isTopBoundary;
395   BitVectorType isLowerBoundary;
396 
397   // Use a messy way of defining a boundary. Using the same way as for PQk elements, only the RT0 basis DOF are triggered,
398   // as they suffice to interpolate a constant function on an edge. In particular RT1, BDM1 DOF are not triggered,
399   // which are linear (but non-constant) on edges.
400 
401   // Mark top boundary.
402   auto topBoundaryIndicator = [&l] (Domain x)
403   {
404 #ifdef DIM2
405     double isBoundary = x[dim-1] > l[dim-1] - 1e-8 ? x[0]: 0.0;
406 #else
407     double isBoundary = x[dim-1] > l[dim-1] - 1e-8 ? x[0] * x[1]: 0.0;
408 #endif
409     return isBoundary;
410   };
411 
412   // Mark lower boundary.
413   auto lowerBoundaryIndicator = [] (Domain x)
414   {
415 #ifdef DIM2
416     double isBoundary = x[dim-1] < 1e-8 ? x[0]: 0.0;
417 #else
418     double isBoundary = x[dim-1] < 1e-8 ? x[0] * x[1]: 0.0;
419 #endif
420     return isBoundary;
421   };
422 
423   VectorType isTopBoundaryTmp, isLowerBoundaryTmp;
424 
425   // Use double-valued interpolation and transfer to char-valued vectors.
426   interpolate(fluxBasis, istlVectorBackend(isTopBoundaryTmp),   topBoundaryIndicator);
427   interpolate(fluxBasis, istlVectorBackend(isLowerBoundaryTmp), lowerBoundaryIndicator);
428   istlVectorBackend(isTopBoundary).resize(basis);
429   istlVectorBackend(isLowerBoundary).resize(basis);
430   isTopBoundary = 0;
431   isLowerBoundary = 0;
432   for (size_t i=0; i<isTopBoundaryTmp[0].size(); i++)
433   {
434     isTopBoundary[  0][i] = isTopBoundaryTmp[  0][i]!=0 ? 1: 0;
435     isLowerBoundary[0][i] = isLowerBoundaryTmp[0][i]!=0 ? 1: 0;
436   }
437 
438 
439   interpolate(fluxBasis, istlVectorBackend(rhs), topFluxBC,   istlVectorBackend(isTopBoundary));
440   interpolate(fluxBasis,istlVectorBackend(rhs), lowerFluxBC, istlVectorBackend(isLowerBoundary));
441 
442   ////////////////////////////////////////////
443   //   Modify Dirichlet rows
444   ////////////////////////////////////////////
445 
446   // loop over the matrix rows
447   for (size_t i=0; i<stiffnessMatrix[0][0].N(); i++)
448   {
449     if (isTopBoundary[0][i] or isLowerBoundary[0][i])
450     {
451       // Lower right matrix block
452       auto cIt    = stiffnessMatrix[0][0][i].begin();
453       auto cEndIt = stiffnessMatrix[0][0][i].end();
454       // loop over nonzero matrix entries in current row
455       for (; cIt!=cEndIt; ++cIt)
456         *cIt = (i==cIt.index()) ? 1. : 0.;
457 
458       // Lower left matrix block
459       for (auto&& entry: stiffnessMatrix[0][1][i])
460         entry = 0.0;
461     }
462 
463   }
464 
465   ////////////////////////////
466   //   Compute solution
467   ////////////////////////////
468 
469   // Start from the rhs vector; that way the Dirichlet entries are already correct
470   VectorType x = rhs;
471 
472   // Technicality:  turn the matrix into a linear operator
473   MatrixAdapter<MatrixType,VectorType,VectorType> op(stiffnessMatrix);
474 
475   // Fancy (but only) way to not have a preconditioner at all
476   Richardson<VectorType,VectorType> preconditioner(1.0);
477 
478   // Preconditioned GMRES / BiCGSTAB solver
479   //RestartedGMResSolver<VectorType> solver (op, preconditioner, 1e-6, 1000, 10000, 2);
480   BiCGSTABSolver<VectorType> solver(op, preconditioner, 1e-6, 4000, 2);
481 
482   // Object storing some statistics about the solving process
483   InverseOperatorResult statistics;
484 
485   // Solve!
486   solver.apply(x, rhs, statistics);
487 
488   ////////////////////////////////////////////////////////////////////////////
489   //  Make a discrete function from the FE basis and the coefficient vector
490   ////////////////////////////////////////////////////////////////////////////
491 
492   using FluxRange = FieldVector<double,dim>;
493   using PressureRange = double;
494 
495   auto fluxFunction = Functions::makeDiscreteGlobalBasisFunction<FluxRange>(fluxBasis, istlVectorBackend(x));
496   auto pressureFunction = Functions::makeDiscreteGlobalBasisFunction<PressureRange>(pressureBasis, istlVectorBackend(x));
497 
498   //////////////////////////////////////////////////////////////////////////////////////////////
499   //  Write result to VTK file
500   //////////////////////////////////////////////////////////////////////////////////////////////
501   SubsamplingVTKWriter<GridView> vtkWriter(gridView, Dune::refinementLevels(0));
502   vtkWriter.addVertexData(fluxFunction, VTK::FieldInfo("flux", VTK::FieldInfo::Type::vector, dim));
503   if (k==0)
504     vtkWriter.addCellData(pressureFunction, VTK::FieldInfo("pressure", VTK::FieldInfo::Type::scalar, 1));
505   else
506     vtkWriter.addVertexData(pressureFunction, VTK::FieldInfo("pressure", VTK::FieldInfo::Type::scalar, 1));
507   vtkWriter.write("poisson-mfem-result");
508 }
509