linear_solvers.html 19.2 KB
Newer Older
Praetorius, Simon's avatar
Praetorius, Simon committed
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 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801
<!DOCTYPE html>
<html>
  <head>
    <title>AMDiS - Adaptive Multi-Dimensional Simulations</title>
    <meta http-equiv="Content-Type" content="text/html; charset=UTF-8"/>
    <style type="text/css">
      @import url(https://fonts.googleapis.com/css?family=Yanone+Kaffeesatz);
      @import url(https://fonts.googleapis.com/css?family=Raleway);
      @import url(https://fonts.googleapis.com/css?family=Ubuntu);
      @import url(https://fonts.googleapis.com/css?family=Droid+Serif:400,700,400italic);
      @import url(https://fonts.googleapis.com/css?family=Ubuntu+Mono:400,700,400italic);

      table tr:nth-child(even) {
          background-color: #eee;
      }
      table tr:nth-child(odd) {
          background-color: #fff;
      }
      table tr td {
        padding: 10px;
        border-left: 1px solid #aaa;
      }
    </style>
    <!--<link rel="stylesheet" type="text/css" href="style_display.css" />-->
    <link rel="stylesheet" type="text/css" href="style_print.css" />
  </head>
  <body>
    <textarea id="source">

class: center, middle

# Session 9
## Linear solvers and preconditioners

---

# Agenda

### Wednesday
- Boundary conditions and Composite FEM
- MTL4 - a linear algebra library
- **Linear solvers and preconditioners**
- Talk by Sebastian Aland

---

# Assemble - Solve

General solution process:
- Assemble linear system, results in
  ```
  Matrix<DOFMatrix*> systemMatrix;
  ```
  either `NULL` or the corresponding block-matrix
--

- SystemMatrix is transformed into `solverMatrix`
  ```
  SolverMatrix<Matrix<DOFMatrix*> > solverMatrix;
  ```
---

# Assemble - Solve

General solution process:
- Assemble linear system, results in `systemMatrix`
- SystemMatrix is transformed into `solverMatrix`
- When a method calls `solverMatrix.getMatrix()` a global matrix is constructed:
  ```
  BlockMapper mapper(solverMatrix);
  using MatrixMapper = MatMap<const SolverMatrix< Matrix<DOFMatrix*> >,
                              BlockMapper>;
  MatrixMapper matMap(solverMatrix, mapper);

  globalMatrix.change_dim(mapper.getNumRows(), mapper.getNumCols());
  set_to_zero(globalMatrix);
  globalMatrix << matMap;
  ```
  using a mapper that assigns each value inside a block to the global matrix.
--

- Thus, every solver works on **global assembled matrix**.

---

# Sequential/Parallel solvers

- In sequential AMDiS, the linear algebra backend MTL4 is used by default
- In parallel AMDiS or if the configuration flag `ENABLE_SEQ_PETSC=ON` the
  libeary PETSc is used as linear algebra backend.

### Each backend provides linear solvers and general preconditioners

<table width="100%">
  <thead>
  <tr>
    <th></th><th>MTL4</th><th>PETSc</th>
  </tr>
  </thead>
  <tbody>
  <tr>
    <th>Solvers</th>
    <td>cg, cgs, bicg, bicgstab, bicgstab2, qmr, tfqmr, bicgstab(l), gmres, fgmres, idr_s, minres, gcr, preonly, umfpack, hypre</td>
    <td>richardson, cg, cgs, fcg, bicg, gmres, fgmres, minres, tcqmr, tfqmr, bicgstab, fbicgstab, bicgstab(l), preonly, umfpack, mumps, superlu, ...</td>
  </tr>
  <tr>
    <th>Preconditioners</th>
    <td>diagonal, identity, mass-lumping, ilu(0), ic(0), solver</td>
    <td>jacobi, none, sor, lu, bjacobi, ilu, icc, ksp, ...</td>
  </tr>
  </tbody>
</table>

---

# Configure a solver for a problem

```matlab
prob->solver:                  solver_name
prob->solver->backend:         mtl    % [mtl, bmtl, petsc, p_mtl, p_petsc]

prob->solver->max iteration:   1000   % default parameters are shown
prob->solver->tolerance:       1.e-8
prob->solver->relative tolerance: 0
prob->solver->print cycle:     100
prob->solver->info:            2
prob->solver->break if tolerance not reached: true

prob->solver->left precon:     no
prob->solver->right precon:    no
```

- If multiple linear algebra **backend**s are installed, the one to use can be chosen.
- **tolerance** reached, if `\(\|b-Ax\| < TOL\)`
- **relative tolerance** reached, if `\(\frac{\|b-Ax\|}{\|b\|} < RTOL\)`
- **print cycle**: print residuum every i'th solver iteration
- **info** greater than 0 means, output of residuum, tolerance, iterations, ...
- Stop the simulation of tolerance are not reached, when **`break if tolerance not reached`** is set to `true`.

---

# Configure a solver for a problem

For parallel/PETSc solvers, some more parameters can/must be set
```matlab
prob->solver->petsc prefix: name_ % must be set for each problem

prob->solver->ksp: ...            % all command line petsc parameters

% set initial guess to zero (only parallel petsc)
prob->solver->use zero start vector: false
```

- Every PETSc solver should have a unique **prefix**
- Using the `solver->ksp` parameter, **command-line arguments** can be pushed to petsc, e.g.
  ```
  prob->solver->ksp: -name_ksp_type cg -name_pc_type jacobi
  ```
  Here, the petsc prefix is used infront of every argument.
- If **`info`** > 10, the flag `-ksp_monitor` is set, if `info` > 20, the flag `ksp_monitor_true_residual` is set.

---

# Write your own solver

Every solver derives from the abstract class `LinearSolverInterface`:

```
class LinearSolverInterface
{
  /* ... */
  protected:
    /// main methods that all solvers must implement
    virtual int solveLinearSystem(SolverMatrix<Matrix<DOFMatrix*>> const& A,
                                  SystemVector& x,
                                  SystemVector& b,
                                  bool createMatrixData,
                                  bool storeMatrixData) = 0;
};
```
where
- **createMatrixData** controls whether the globalMatrix is created from the solverMatrix
- **storeMatrixData** set whether the matrix can be reused in subsequent calls, e.g.
  necessary for reusing a matrix-factorization for directo solvers

---

# Example: Richardson iteration

\\[
x^{n+1} = x^{n} - \alpha B(A x^{n} - b)
\\]
where `\(B\)` is a preconditioner and `\(\alpha\)` a dumping factor. For simplicity of the example, we set `\(B=I\)`.

```
class RichardsonSolver : public LinearSolverInterface
{
  public:
    RichardsonSolver(std::string name)
      : LinearSolverInterface(name)
      , alpha(0.1)
    {
      Parameters::get(name + "->alpha", alpha);
    }

    // ...

  private:
    double alpha;
};
```

---

# Example: Richardson iteration

\\[
x^{n+1} = x^{n} - \alpha B(A x^{n} - b)
\\]
where `\(B\)` is a preconditioner and `\(\alpha\)` a dumping factor.


```
class RichardsonSolver : public LinearSolverInterface
{
  // ...

  protected:
    virtual int solveLinearSystem(SolverMatrix<Matrix<DOFMatrix*>> const& A,
                                  SystemVector& x,
                                  SystemVector& b,
                                  bool /* createMatrixData */,
                                  bool /* storeMatrixData */) override
    {
      SystemVector r(b);

      auto const& mat = *A.getOriginalMat();
      for (size_t i = 0; i < max_iter; ++i) {
        mv(mat, x, r);
        r -= b;
        axpy(-alpha, r, x);
      }
    }
};
```

---

# Example: Richardson iteration
## Register new solver to AMDiS

```
class RichardsonSolver : public LinearSolverInterface
{
public:
  struct Creator : LinearSolverCreator
  {
    LinearSolverInterface* create()
    {
      return new RichardsonSolver(this->name);
    }
  };

  // ...
};

int main(...) {
  AMDiS::init(...);
  CreatorMap<LinearSolverInterface>::addCreator("mtl_richard",
    new RichardsonSolver::Creator);
}
```
- register for specific backend (here: `mtl`)
- can now be selected in init-file...

---

# Example: Richardson iteration
## Register new solver to AMDiS

```matlab
prob->solver:                  richard
prob->solver->backend:         mtl
prob->solver->tolerance:       1.e-6
prob->solver->alpha:           0.1
prob->solver->max iteration:   100000
```

---

# Example2: Add ITL solver

Sometimes it is easier to work with MTL data structures directly:

```
struct richardson
{
  richardson(std::string name)
  {
    Parameters::get(name + "->alpha", alpha);
  }

  template <class Matrix, class Vector, class Precon, class Precon2, class Iter>
  int operator()(Matrix const& A, Vector& x, Vector const& b,
                 Precon const& P, Precon2 const&, Iter& iter)
  {
    Vector r(b - A*x);
    double residual = two_norm(r);
    for (; !iter.finished(residual); ++iter) {
      x += alpha * Vector(solve(P, r));
      r = b - A*x;
      residual = two_norm(r);
    }
    return iter;
  }

private:
  double alpha = 0.1;
};

using RichardsonSolver = ITL_Solver<richardson>;
```

---

# Write your own preconditioner

Structure of preconditioner is backend depended. Here: MTL preconditioner.

Implement a preconditioner interface
```
template <class MatrixType, class VectorType>
struct ITL_PreconditionerBase : public PreconditionerInterface
{
  virtual void init(SolverMatrix<Matrix<DOFMatrix*>> const& A,
                    MatrixType const& globalMatrix) = 0;
  virtual void exit() {}

  // P x = b
  virtual void solve(VectorType const& b, VectorType& x) const = 0;
  // P^T x = b
  virtual void adjoint_solve(VectorType const& b, VectorType& x) const = 0;
};
```

---

# Example: SOR preconditioner
\\[
x^{k+1}\_i = (1-\omega) x^{k}\_i + \frac{\omega}{a\_{ii}}\left(b\_i - \sum\_{j=0}^{i-1} a\_{ij} x^{k+1}\_j - \sum\_{j=i+1}^{n} a\_{ij} x^{k}\_j\right)
\\]

```
class SORPrecon
  : public ITL_PreconditionerBase<MTLTypes::MTLMatrix, MTLTypes::MTLVector>
{
  virtual void init(SolverMatrix<Matrix<DOFMatrix*>> const& A,
                    MatrixType const& globalMatrix) override
  {
    mat = &globalMatrix;
    diag = diagonal(globalMatrix);

    for (auto const& d : diag)
      assert( std::abs(d) > 1.e-10 );
  }

  // ...

  MTLTypes::MTLMatrix const* mat = NULL;
  MTLTypes::MTLVector diag;
};
```

---

# Example: SOR preconditioner
\\[
x^{k+1}\_i = (1-\omega) x^{k}\_i + \frac{\omega}{a\_{ii}}\left(b\_i - \sum\_{j=0}^{i-1} a\_{ij} x^{k+1}\_j - \sum\_{j=i+1}^{n} a\_{ij} x^{k}\_j\right)
\\]

```
class SORPrecon
  : public ITL_PreconditionerBase<MTLTypes::MTLMatrix, MTLTypes::MTLVector>
{
  // ...

  virtual void solve(VectorType const& b, VectorType& x) const override
  {
    auto col = col_map(*mat);
    auto value = const_value_map(*mat);

    x.change_dim(num_rows(b));
    for (auto i : rows_of(*mat)) {
      double delta = 0.0;
      for (auto j : nz_of(i)) // traverse non-zeros in row i
        delta += value(j) * x[col(j)];

      // update solution component
      delta = (b[i.value()] - delta) * omega / diag[i.value()];
      x[i.value()] += delta;
    }
  }
  double omega = 0.8;
};

```

---

# Example: SOR preconditioner
## Register new preconditioner to AMDiS

```
using PreconBase = ITL_PreconditionerBase<MTLTypes::MTLMatrix,
                                          MTLTypes::MTLVector>;

class SORPrecon : public PreconBase
{
public:
  struct Creator : CreatorInterfaceName<PreconBase>
  {
    PreconBase* create()
    {
      return new SORPrecon(this->name);
    }
  };

  // ...
};

int main(...) {
  AMDiS::init(...);
  CreatorMap<PreconBase>::addCreator("sor", new SORPrecon::Creator);
}
```
- can now be selected in init-file...

---

# Example: SOR preconditioner
## Register new preconditioner to AMDiS

```matlab
prob->solver:                  fgmres
prob->solver->tolerance:       1.e-6
prob->solver->left precon:     sor
```

---

# Block-Preconditioner

Consider system of equations. E.g. Cahn-Hilliard equation
\\[
  \partial\_t c = \nabla\cdot(M(c)\nabla\mu) \\\\
  \mu = -\epsilon^2\Delta c + f'(c)
\\]
Results in block system
\\[
  \begin{pmatrix} M & \tau L \\\\ -\epsilon^2 L - J & M \end{pmatrix}
\\]
with `\(M_{ij} = (\theta_i, \theta_j)_\Omega\)`, `\(L_{ij} = (\nabla\theta_i, \nabla\theta_j)_\Omega\)` and `\(J\)` a discretization of the jacobian of `\(f'(c)\)`.

--

### Block-diagonal preconditioner

\\[
P = \begin{pmatrix} M & 0 \\\\ 0 & M \end{pmatrix}
\\]

---

# Block-Preconditioner
### Block-diagonal preconditioner

\\[
P = \begin{pmatrix} M & 0 \\\\ 0 & M \end{pmatrix}
\\]
Invert block system, by applying a **krylov subspace method** to the blocks `\(M\)`.

--

```
class BlockDiagonalPrecon
  : public BlockPreconditioner<MTLTypes::MTLMatrix>
{
  using Super = BlockPreconditioner<MTLTypes::MTLMatrix>;
  using MatrixType = MTLTypes::MTLMatrix;
  using VectorType = MTLTypes::MTLVector;

  BlockDiagonalPrecon(std::string name)
  {
    Super::createSubSolver(name + "->subsolver", solver, runner, ...);
  }

  // ...

  LinearSolverInterface* solver = NULL;
  RunnerBase<MatrixType, VectorType>* runner = NULL;
};
```


---

# Block-Preconditioner
### Block-diagonal preconditioner

```
class BlockDiagonalPrecon
  : public BlockPreconditioner<MTLTypes::MTLMatrix>
{
  // ...

  virtual void init(SolverMatrix<Matrix<DOFMatrix*> > const& A,
                    MatrixType const& globalMatrix) override
  {
    Super::init(A, globalMatrix);
    runner->init(A, A.getSubMatrix(0,0));
  }

  virtual void exit() override
  {
    runner->exit();
  }

  // ...
};
```

---

# Block-Preconditioner
### Block-diagonal preconditioner

```
class BlockDiagonalPrecon
  : public BlockPreconditioner<MTLTypes::MTLMatrix>
{
  // ...

  virtual void solve(VectorType const& b, VectorType& x) const override
  {
    x.change_dim(num_rows(b));

    const VectorType b0(b[getRowRange(0)]);
    const VectorType b1(b[getRowRange(1)]);

    VectorType x0(x[getColRange(0)]);
    VectorType x1(x[getColRange(1)]);

    runner->solve(getSubMatrix(0,0), x0, b0);
    runner->solve(getSubMatrix(1,1), x1, b1);
  }
};
```

---

# Block-Preconditioner
### Block-diagonal preconditioner

Configure the block preconditioner in the init-file!

```
int main(...) {
  AMDiS::init(...);
  CreatorMap<PreconBase>::addCreator("mtl_block_diag",
    new BlockDiagonalPrecon::Creator);
}
```

```matlab
prob->solver:        fgmres
prob->right precon:  block_diag
prob->right precon->subsolver:                     cg
prob->right precon->subsolver->max iterations:     100
prob->right precon->subsolver->relative tolerance: 1.e-4
prob->right precon->subsolver->left precon:        diag
```

- When iterative methods used in preconditioner: Use flexible krylov-subspace methods, e.g. `fgmres`

---

# Matrix-Shells
## (Advanced)

Instead of assembling a full matrix, only **assemble individual operators** and implement composed
**matrix-vector product**. Many krylov-subpsace methods based one matrix-vector products.

### Example: Cahn-Hilliard
\\[
  \begin{pmatrix} M & \tau L \\\\ -\epsilon^2 L - J & M \end{pmatrix}\begin{pmatrix}x\_0 \\\\ x\_1\end{pmatrix}
\\]
leads to
\\[
  M x\_0 + \tau L x\_1 \\\\ -\epsilon^2 L x\_0 - J x\_0 + M x\_1
\\]
- Only three matrices necessary: `\(M, L, J\)`,
- Only `\(J\)` changes over time (when mesh is fixed).

---

# Matrix-Shells

MTL provides matrix-shell interface, based on multiplication method:

```
class MatrixShell : BlockMTLMatrix
{
public:
  template <class VectorIn, class VectorOut, class Assign>
  void mult(VectorIn const& b, VectorOut& x, Assign) const
  {
    // perform matrix-vector multiplication
  }

  template <class VectorIn>
  mtl::vector::mat_cvec_multiplier<MatrixShell, VectorIn>
  operator*(VectorIn const& v) const
  {
    return mtl::vector::mat_cvec_multiplier<MatrixShell, VectorIn>(*this, v);
  }

  friend size_t size(MatrixShell const& A) { return (A.n) * (A.m); }
  friend size_t num_rows(MatrixShell const& A) { return A.n; }
  friend size_t num_cols(MatrixShell const& A) { return A.m; }

private:
  size_t n, m;
};
```

---

# Matrix-Shells

MTL provides matrix-shell interface, based on multiplication method:

```
class MatrixShell;

namespace mtl
{
  template <> struct Collection<MatrixShell>
  {
    using value_type = double;
    using size_type  = size_t;
  };

  namespace ashape
  {
    template <> struct ashape_aux<MatrixShell>
    {
      using type = nonscal;
    };
  } // end namespace ashape
} // end namespace mtl
```

---

# Matrix-Shells

MTL provides matrix-shell interface, based on multiplication method:

\\[
  y\_0 := M x\_0 + \tau L x\_1 \\\\ y\_1 := -\epsilon^2 L x\_0 - J x\_0 + M x\_1
\\]
```
template <class VectorIn, class VectorOut, class Assign>
void MatrixShell::mult(VectorIn const& x, VectorOut& y, Assign) const
{
  VectorOut y0(y[getColRange(0)]);
  VectorOut y1(y[getColRange(1)]);

  const VectorIn x0(x[getRowRange(0)]);
  const VectorIn x1(x[getRowRange(1)]);

  // y0 = M*x0 + tau*L*x1
  Assign::first_update(y0, M * x0 );
  Assign::update(      y0, L * (tau*x1) );

  // y1 = -eps^3*L*x0 - J*x0 + M*x1
  Assign::first_update(y1, L * ((-eps*eps) * x0) );
  Assign::update(      y1, J * (-x0) );
  Assign::update(      y1, M * x1 );
}
```

---

# Matrix-Shells

Assemble necessary operators

```
void fillOperators()
{
  DOFVector<double>* c = getSolution(0);
  const FiniteElemSpace* feSpace = c->getFeSpace();

  // M := (phi, psi)
  Operator *mass_op = new Operator(feSpace);
  addZOT(mass_op, 1.0);
  addMatrixOperator(*mass_op, 0,0);
  setAssembleMatrixOnlyOnce(0,0,true);

  // L := (grad(phi), grad(psi))
  Operator *laplace_op = new Operator(feSpace);
  addSOT(laplace_op, 1.0);
  addMatrixOperator(*laplace_op, 1,1);
  setAssembleMatrixOnlyOnce(1,1,true);

  // J := (-2*c^2 * phi, psi)
  Operator *rho3_op2 = new Operator(feSpace);
  addZOT(rho3_op2, 1.0 - 3.0 * pow<2>(valueOf(c)));
  addMatrixOperator(*rho3_op2, 1,0);
  setAssembleMatrixOnlyOnce(1,0,false);

  // ...
```

---

# Matrix-Shells

Assemble necessary operators

```
  // ...

  // --- vector operators ----------------------------------------------------

  /// (c, psi)
  Operator *mass_op_rhs = new Operator(feSpace);
  addZOT(mass_op_rhs, valueOf(c));
  addVectorOperator(*mass_op_rhs, 0);

  /// (-c^3, psi)
  Operator *rho3_op_rhs = new Operator(feSpace);
  addZOT(rho3_op_rhs, -2.0 * pow<3>(valueOf(c)));
  addVectorOperator(*rho3_op_rhs, 1);
}
```

--

### Register matrix as solver-backend

```
int main(...) {
  AMDiS::init(...);
  initCreatorMap<MatrixShell>("cahn_hilliard");
  // ...
}
```

---

# Matrix-Shells

Assemble necessary operators

### Register matrix as solver-backend

```
int main(...) {
  AMDiS::init(...);
  initCreatorMap<MatrixShell>("cahn_hilliard");
  // ...
}
```

Init-file:
```
prob->solver:             fgmres
prob->solver->backend:    cahn_hilliard
```
- Can be combined with block-preconditioners
- Example: `extensions/demo/cahn_hilliard/src/cahnHilliard_shell.cc`

---

class: center, middle
# Exercise

---

# Exercise
Incorporate the new symmetric Dirichlet boundary condition into AMDiS:

1. Create a new class that derives from `DirichletBC` or `BoundaryCondition`.
2. Override the methods `fillBoundaryCondition(DOFMatrix*,...)` and `fillBoundaryCondition(DOFVectorBase<double>*,...)`
3. set `bool isDirichlet() { return false; }` to enforce your manual implementation.
4. Create a new class that derives from `ProblemStat` and override the method `addDirichletBC` to add your new boundary condition implementation.

    </textarea>
    <script src="lib/remark.js" type="text/javascript"></script>
802
    <script type="text/javascript" src="MathJax/MathJax.js?config=TeX-AMS_HTML"></script>
Praetorius, Simon's avatar
Praetorius, Simon committed
803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825

    <script type="text/javascript">
      var slideshow = remark.create({
        ratio: "4:3",
        highlightLanguage: "cpp"
      });

      // Setup MathJax
      MathJax.Hub.Config({
          tex2jax: {
          skipTags: ['script', 'noscript', 'style', 'textarea', 'pre']
          }
      });
      MathJax.Hub.Queue(function() {
          $(MathJax.Hub.getAllJax()).map(function(index, elem) {
              return(elem.SourceElement());
          }).parent().addClass('has-jax');
      });

      MathJax.Hub.Configured();
    </script>
  </body>
</html>