navierstokes.cc 9.18 KB
Newer Older
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
#include "AMDiS.h"

using namespace std;
using namespace AMDiS;

// ===========================================================================
// ===== function definitions ================================================
// ===========================================================================

class Project : public AbstractFunction<double, WorldVector<double> >
{
public:
  MEMORY_MANAGED(Project);

  /** \brief
   * Constructor
   */
  Project(int i) : comp(i) {};

  /** \brief
   * Implementation of AbstractFunction::operator().
   */
  const double& operator()(const WorldVector<double>& x) const
  {
    static double result;
    result = x[comp];
    return result;  
  };

protected:
  /** \brief
   * coordinate plane to project at
   */
  int comp;
};


/** \brief
 * RHS function
 */
class ConstantFct : public AbstractFunction<double, WorldVector<double> >
{
public:
  MEMORY_MANAGED(ConstantFct);

  ConstantFct(double c) 
    : AbstractFunction<double, WorldVector<double> >(0) ,
      constant(c)
  {};

  /** \brief
   * Implementation of AbstractFunction::operator().
   */
  const double& operator()(const WorldVector<double>& x) const 
  {
    static double result = 0.0;
    result = constant;
    return result;
  };

protected:
  double constant;
};

class F : public AbstractFunction<double, WorldVector<double> >
{
public:
  MEMORY_MANAGED(F);

  F() 
    : AbstractFunction<double, WorldVector<double> >(0)
  {};

  /** \brief
   * Implementation of AbstractFunction::operator().
   */
  const double& operator()(const WorldVector<double>& x) const 
  {
    static double result = 0.0;
    result = (x[0] == 0.0) || (x[0] == 1.0) ? 0.0 : 1.0;
    return result;
  };

protected:
  double constant;
};



// ===========================================================================
// ===== instationary problem ================================================
// ===========================================================================

/** \brief
 * Instationary problem
 */
class Navierstokes : public ProblemInstatVec
{
public:
  MEMORY_MANAGED(Navierstokes);

  /** \brief
   *  Constructor
   */
  Navierstokes(ProblemVec *navierstokesSpace) 
    : ProblemInstatVec("navierstokes", navierstokesSpace)
  {
    int i, dow = Global::getGeo(WORLD);
    TEST_EXIT(dow == navierstokesSpace->getNumComponents() - 1)
      ("number of components != dow + 1\n");
    boundaryFct.resize(dow+1, NULL);

    epsilon = -1.0;
    GET_PARAMETER(0, "epsilon", "%f", &epsilon);
    TEST_EXIT(epsilon != -1.0)("invalid epsilon\n");
  };

  // ===== initial problem methods =====================================

  /** \brief
   * Used by \ref problemInitial to solve the system of the initial problem
   */
  void solve(AdaptInfo *adaptInfo) 
  {
    int i, size = problemStat->getNumComponents()-1;
    for(i = 0; i < size; i++) {
      //problemStat->getMesh(i)->dofCompress();
      problemStat->getSolution()->getDOFVector(i)->interpol(boundaryFct[i]);
    }
  };

  /** \brief
   * Used by \ref problemInitial to do error estimation for the initial
   * problem.
   */
  void estimate(AdaptInfo *adaptInfo) 
  {
    int i, size = problemStat->getNumComponents()-1;
    double errMax, errSum;

    for(i = 0; i < size; i++) {
      errSum = 
	Error<double>::L2Err(*boundaryFct[i],
			     *(problemStat->getSolution()->getDOFVector(i)), 
			     0, &errMax, false);

      adaptInfo->setEstSum(errSum, i);
    }
  };

  /** \brief
   * Sets \ref boundaryFct;
   */
  void setBoundaryFct(AbstractFunction<double, WorldVector<double> > *fct, int i) 
  {
    boundaryFct[i] = fct;
  };

  /** \brief
   * Set the time in all needed functions! Called by adaption loop.
   */
  void setTime(AdaptInfo *adaptInfo) 
  {
    time = adaptInfo->getTime();
    tau1 = 1.0 / adaptInfo->getTimestep();
  };
  
  /** \brief
   * Returns address of \ref time.
   */
  double *getTimePtr() {
    return &time;
  };

  /** \brief
   * Returns address of \ref tau1.
   */
  double *getTau1Ptr() {
    return &tau1;
  };

  /** \brief
   * Implementation of ProblemInstatBase::closeTimestep()
   * pressure correction step
   */
  void closeTimestep(AdaptInfo *adaptInfo) {
    problemStat->writeFiles(adaptInfo, true);

    WAIT;

    int i, numComponents = problemStat->getNumComponents();
    int dow = Global::getGeo(WORLD);

    double timestep = adaptInfo->getTimestep();

    static DOFVector<double> pInterpol(problemStat->getFESpace(0), "p interpolated");
    static DOFVector<WorldVector<double> > gradPInterpol(problemStat->getFESpace(0), "gradient p interpolated");
    DOFVector<double> *p = problemStat->getSolution()->getDOFVector(numComponents-1);

    pInterpol.interpol(p, timestep * epsilon);
    pInterpol.getRecoveryGradient(&gradPInterpol);

    for(i = 0; i < dow; i++) {
      DOFVector<double>::Iterator 
	vIt(problemStat->getSolution()->getDOFVector(i), USED_DOFS);
      DOFVector<WorldVector<double> >::Iterator 
	gradIt(&gradPInterpol, USED_DOFS);

      for(vIt.reset(), gradIt.reset(); !vIt.end(); ++vIt, ++gradIt) {
	*vIt -= (*gradIt)[i];
      }
    }

    problemStat->writeFiles(adaptInfo, true);
  };

private:
  /** \brief
   * boundary function
   */
  std::vector<AbstractFunction<double, WorldVector<double> >*> boundaryFct;

  /** \brief
   * time
   */
  double time;

  /** \brief
   * 1.0/timestep
   */
  double tau1;

  /** \brief
   */
  double epsilon;
};

// ===========================================================================
// ===== main program ========================================================
// ===========================================================================

int main(int argc, char** argv)
{
  FUNCNAME("main");

  // ===== check for init file =====
  TEST_EXIT(argc > 1)("usage: navierstokes initfile\n");

  // ===== init parameters =====
  Parameters::init(false, argv[1]);
  Parameters::readArgv(argc, argv);

  // read nu (1/Re)
  double nu = -1.0;
  GET_PARAMETER(0, "nu", "%f", &nu);
  TEST_EXIT(nu != -1.0)("invalid nu\n");

  // ===== create and init stationary problem =====
  ProblemVec *navierstokesSpace = NEW ProblemVec("navierstokes->space");
  navierstokesSpace->initialize(INIT_ALL);

  int numComponents = navierstokesSpace->getNumComponents();

  // ===== create instationary problem =====
  Navierstokes *navierstokes = NEW Navierstokes(navierstokesSpace);
  navierstokes->initialize(INIT_ALL);

  // ===== create adapt Info =====
  AdaptInfo *adaptInfo = NEW AdaptInfo("navierstokes->adapt",
				       numComponents);

  // create initial adapt
  AdaptInfo *adaptInfoInitial = NEW AdaptInfo("navierstokes->initial->adapt",
					      numComponents);

  // create instationary adapt
  AdaptInstationary *adaptInstat =
    NEW AdaptInstationary("navierstokes->adapt",
			  navierstokesSpace,
			  adaptInfo,
			  navierstokes,
			  adaptInfoInitial);

  // ===== create boundary functions =====
  
  double *tau1 = navierstokes->getTau1Ptr();
  int i;

  // === boundary functions ===
  ConstantFct *boundaryFct = NEW ConstantFct(0.0);

  for(i = 0; i < numComponents-1; i++) {
    navierstokes->setBoundaryFct(boundaryFct, i);
    navierstokesSpace->addDirichletBC(1, i, boundaryFct);
  }
  navierstokesSpace->addDirichletBC(2, 0, NEW F());
  navierstokes->setBoundaryFct(boundaryFct, numComponents-1);
  

  // ===== create operators =====

  // ===== velocity =====
  for(i = 0; i < numComponents-1; i++) {
    // create laplace
    Operator *laplaceV = NEW Operator(Operator::MATRIX_OPERATOR,
				      navierstokesSpace->getFESpace(i),
				      navierstokesSpace->getFESpace(i));

    laplaceV->addSecondOrderTerm(NEW FactorLaplace_SOT(nu));
    navierstokesSpace->addMatrixOperator(laplaceV, i, i);

    // create zero order operator
    Operator *dtv = NEW Operator(Operator::MATRIX_OPERATOR |
				 Operator::VECTOR_OPERATOR,
				 navierstokesSpace->getFESpace(i),
				 navierstokesSpace->getFESpace(i));

    dtv->setUhOld(navierstokes->getOldSolution()->getDOFVector(i));
    dtv->addZeroOrderTerm(NEW Simple_ZOT);

    navierstokesSpace->addMatrixOperator(dtv, i, i, tau1, tau1);
    navierstokesSpace->addVectorOperator(dtv, i, tau1, tau1);

    
    Operator *gradP = NEW Operator(Operator::VECTOR_OPERATOR,
				   navierstokesSpace->getFESpace(i));
    gradP->addZeroOrderTerm(NEW FctGradient_ZOT(navierstokes->getOldSolution()->
						getDOFVector(numComponents-1),
						NEW Project(i)));
    double minusOne = 1.0;
    navierstokesSpace->addVectorOperator(gradP, i, &minusOne, &minusOne);
  }

  // ===== pressure =====
  Operator *laplaceP = NEW Operator(Operator::MATRIX_OPERATOR,
				    navierstokesSpace->getFESpace(numComponents-1),
				    navierstokesSpace->getFESpace(numComponents-1));

  laplaceP->addSecondOrderTerm(NEW Laplace_SOT());

  navierstokesSpace->addMatrixOperator(laplaceP, 
				       numComponents-1, 
				       numComponents-1);
  
  WorldVector<double> b;

  for(i = 0; i < numComponents-1; i++) {
    Operator *divV = NEW Operator(Operator::MATRIX_OPERATOR,
				  navierstokesSpace->getFESpace(numComponents-1),
				  navierstokesSpace->getFESpace(i));

    b.set(0.0);
    b[i] = 1.0;

    divV->addFirstOrderTerm(NEW Vector_FOT(b), GRD_PHI);

    navierstokesSpace->addMatrixOperator(divV, numComponents-1, i, tau1, tau1);    
  }

  // ===== start adaption loop =====
  int errorCode = adaptInstat->adapt();
  
  return errorCode;
}