navier_stokes.cc 3.27 KB
Newer Older
1
2
3
4
5
#include <iostream>
#include <ctime>
#include <cmath>

#include <dune/amdis/AMDiS.hpp>
6
#include <dune/amdis/AdaptInstationary.hpp>
7
#include <dune/amdis/ProblemStat.hpp>
8
#include <dune/amdis/ProblemInstat.hpp>
9
10
#include <dune/amdis/Operators.hpp>
#include <dune/amdis/assembler/StokesOperator.hpp>
11

12
13
#ifdef DOW
#undef DOW
14
#endif
15
#define DOW AMDIS_DOW
16

17
18
19
20
#ifndef STRATEGY
#define STRATEGY 1
#endif

21
22
using namespace AMDiS;

23
24
25
26
27
struct NavierStokesBasis
{
  using Grid = Dune::AlbertaGrid<AMDIS_DIM,AMDIS_DOW>;
  using GlobalBasis = typename TaylorHoodBasis<Grid::LeafGridView>::GlobalBasis;
};
28

29
30
using StokesProblem = ProblemStat<NavierStokesBasis>;
using StokesProblemInstat = ProblemInstat<NavierStokesBasis>;
31

32
33
int main(int argc, char** argv)
{
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
  AMDiS::init(argc, argv);

  StokesProblem prob("stokes");
  prob.initialize(INIT_ALL);

  StokesProblemInstat probInstat("stokes", prob);
  probInstat.initialize(INIT_UH_OLD);

  double viscosity = 1.0;
  double density = 1.0;
  double vel = 1.0;
  Parameters::get("stokes->viscosity", viscosity);
  Parameters::get("stokes->density", density);
  Parameters::get("stokes->boundary velocity", vel);

  AdaptInfo adaptInfo("adapt", prob.getNumComponents());

  // tree-paths for components
  auto _v = 0_c;
  auto _p = 1_c;

  // <1/tau * u, v>
56
  auto opTime = makeOperator(tag::testvec_trialvec{}, density);
57
58
59
  prob.addMatrixOperator(opTime, _v, _v);

  // <1/tau * u^old, v>
60
  auto opTimeOld = makeOperator(tag::testvec{}, density * prob.getSolution(_v));
61
62
63
  prob.addVectorOperator(opTimeOld, _v);


64
65
66
  // sum_i <grad(u_i),grad(v_i)> + <p, div(v)> + <div(u), q>
  auto opStokes = makeOperator(tag::stokes{}, viscosity);
  prob.addMatrixOperator(opStokes, treepath(), treepath());
67

68
69
70
  // <(u * nabla)u_i^old, v_i>
  auto opNonlin1 = makeOperator(tag::testvec_trialvec{}, density * trans(gradientAtQP(prob.getSolution(_v))));
  prob.addMatrixOperator(opNonlin1, _v, _v);
71

72
  for (std::size_t i = 0; i < DOW; ++i) {
73
    // <(u^old * nabla)u_i, v_i>
74
    auto opNonlin2 = makeOperator(tag::test_gradtrial{}, density * prob.getSolution(_v));
75
76
    prob.addMatrixOperator(opNonlin2, treepath(_v,i), treepath(_v,i));
  }
77

78
79
80
81
  // <(u^old * grad(u_i^old)), v_i>
  auto opNonlin3 = makeOperator(tag::testvec{}, trans(gradientAtQP(prob.getSolution(_v))) * prob.getSolution(_v));
  prob.addVectorOperator(opNonlin3, _v);

82
83
84
  // define boundary regions
  auto left = [](auto const& x) { return x[0] < 1.e-8; };
  auto not_left = [](auto const& x) { return x[0] > 1.0 - 1.e-8 || x[1] < 1.e-8 || x[1] > 1.0 - 1.e-8; };
85

86
87
88
89
90
  // define boundary values
  auto parabolic_y = [](auto const& x) -> Dune::FieldVector<double,2>
  {
    return {0.0, x[1]*(1.0 - x[1])};
  };
91

92
93
94
95
  auto zero = [](auto const& x) -> Dune::FieldVector<double,2>
  {
    return {0.0, 0.0};
  };
96

97
98
99
  // set boundary conditions for velocity
  prob.addDirichletBC(left, _v, _v, parabolic_y);
  prob.addDirichletBC(not_left, _v, _v, zero);
100

101
102
103
  // set point constraint for pressure
  prob.addDirichletBC([](auto const& x) { return x[0] < 1.e-8 && x[1] < 1.e-8; }, _p, _p, 0.0);

104
105
106
  // set initial conditions
  prob.getSolution(_v).interpolate(parabolic_y);
  prob.getSolution(_p).interpolate(0.0);
107
108
109
110
111
112
113
114
115

  AdaptInstationary adapt("adapt", prob, adaptInfo, probInstat, adaptInfo);
  adapt.adapt();

  // output solution
  prob.writeFiles(adaptInfo);

  AMDiS::finalize();
  return 0;
116
}