vecellipt.cc 1.55 KB
Newer Older
1
2
3
4
#ifdef HAVE_CONFIG_H
#include "config.h"
#endif

5
6
#include <iostream>

7
#include <amdis/AMDiS.hpp>
8
#include <amdis/LocalOperators.hpp>
9
#include <amdis/ProblemStat.hpp>
10
11
12
13

using namespace AMDiS;

// 1 component with polynomial degree 1
14
15
//using Grid = Dune::AlbertaGrid<AMDIS_DIM, AMDIS_DOW>;
using ElliptParam   = YaspGridBasis<AMDIS_DIM, 2,2>;
16
17
18
19
20
21
22
23
24
25
26
using ElliptProblem = ProblemStat<ElliptParam>;

int main(int argc, char** argv)
{
  AMDiS::init(argc, argv);

  ElliptProblem prob("ellipt");
  prob.initialize(INIT_ALL);

  AdaptInfo adaptInfo("adapt");

27
  auto opL = makeOperator(tag::gradtest_gradtrial{}, 1.0);
28
29
  prob.addMatrixOperator(opL, 1, 1);

30
31
32
  auto opM = makeOperator(tag::test_trial{}, 1.0);
  prob.addMatrixOperator(std::ref(opM), 0, 0);
  prob.addMatrixOperator(std::ref(opM), 1, 0);
33

34
  auto opForce = makeOperator(tag::test{}, [](auto const& x) { return -1.0; }, 0);
35
36
37
38
39
40
41
42
  prob.addVectorOperator(opForce, 0);


  // set boundary condition
  auto predicate = [](auto const& x){ return x[0] < 1.e-8 || x[1] < 1.e-8; }; // define boundary
  auto dbcValues = [](auto const& x){ return 0.0; }; // set value
  prob.addDirichletBC(predicate, 1, 1, dbcValues);

43
  prob.buildAfterAdapt(adaptInfo, Flag(0));
44

45
#ifdef DEBUG_MTL
46
47
48
  // write matrix to file
  if (Parameters::get<int>("elliptMesh->global refinements").value_or(0) < 4) {
    mtl::io::matrix_market_ostream out("matrix.mtx");
Praetorius, Simon's avatar
Praetorius, Simon committed
49
    out << prob.systemMatrix().matrix();
50

Praetorius, Simon's avatar
Praetorius, Simon committed
51
    std::cout << prob.systemMatrix().matrix() << '\n';
52
  }
53
#endif
54
55
56
57
58
59
60

  prob.solve(adaptInfo);
  prob.writeFiles(adaptInfo, true);

  AMDiS::finalize();
  return 0;
}