vecheat.cc 10.6 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
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
#include "AMDiS.h"

using namespace std;
using namespace AMDiS;

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

/** \brief
 * Dirichlet boundary function
 */
class G : public AbstractFunction<double, WorldVector<double> >,
	  public TimedObject
{
public:
  MEMORY_MANAGED(G);

  /** \brief
   * Implementation of AbstractFunction::operator().
   */
  const double& operator()(const WorldVector<double>& x) const
  {
    static double result;
    result = sin(M_PI*(*timePtr)) * exp(-10.0*(x*x));
    return result;  
  };
};

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

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

  /** \brief
   * Implementation of AbstractFunction::operator().
   */
  const double& operator()(const WorldVector<double>& x) const 
  {
    static double result;
    int dim = x.getSize();
    double r2 = (x*x);
    double ux = sin(M_PI * (*timePtr)) * exp(-10.0*r2);
    double ut = M_PI * cos(M_PI*(*timePtr)) * exp(-10.0*r2);
    result = ut -(400.0*r2 - 20.0*dim)*ux;
    return result;
  };
};

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

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

  /** \brief
   *  Constructor
   */
  Vecheat(ProblemVec *vecheatSpace) 
    : ProblemInstatVec("vecheat", vecheatSpace)
  {
    // init theta scheme
    theta = -1.0;
    GET_PARAMETER(0, name + "->theta", "%f", &theta);
    TEST_EXIT(theta >= 0)("theta not set!\n");
    if (theta == 0.0) {
      WARNING("You are using the explicit Euler scheme\n");
      WARNING("Use a sufficiently small time step size!!!\n");
    }
    MSG("theta = %f\n", theta);
    theta1 = theta - 1;
  };

  // ===== ProblemInstatBase methods ===================================

  /** \brief
   * set the time in all needed functions!
   */
  void setTime(AdaptInfo *adaptInfo) 
  {
    rhsTime = adaptInfo->getTime() - (1-theta) * adaptInfo->getTimestep();
    boundaryTime = adaptInfo->getTime();
    tau1 = 1.0 / adaptInfo->getTimestep();    
  };

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

  /** \brief
   * Used by \ref problemInitial to solve the system of the initial problem
   */
  void solveInitialProblem(AdaptInfo *adaptInfo) 
  {
    int i, size = problemStat->getNumComponents();

    boundaryTime = rhsTime = 0.0;
    for(i = 0; i < size; i++) {
      problemStat->getMesh(i)->dofCompress();
      problemStat->getSolution()->getDOFVector(i)->interpol(boundaryFct);
    }
  };

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

    boundaryTime = 0.0;

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

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

  /** \brief
   * Used by \ref problemInitial to build before refinement.
   */
  void buildBeforeRefineInitial(Flag) 
  {};

  /** \brief
   * Used by \ref problemInitial to build before coarsening.
   */
  void buildBeforeCoarsenInitial(Flag) 
  {};

  /** \brief
   * Used by \ref problemInitial to build after coarsening.
   */
  void buildAfterCoarsenInitial(Flag) 
  {};

  /** \brief
   * Used by \ref problemInitial to mark elements.
   */
  Flag markElementsInitial() 
  { 
    return 0; 
  };

  /** \brief
   * Used by \ref problemInitial
   */
  Flag refineMeshInitial() 
  { 
    return 0; 
  };

  /** \brief
   * Used by \ref problemInitial
   */
  Flag coarsenMeshInitial() 
  { 
    return 0; 
  };

  /** \brief
   * Used by \ref problemInitial
   */
  void beginIterationInitial(int ) 
  {};

  /** \brief
   * Used by \ref problemInitial
   */
  void endIterationInitial(int ) 
  {};

  // ===== setting methods ===============================================

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

  // ===== getting methods ===============================================

  /** \brief
   * Returns pointer to \ref theta.
   */
  double *getThetaPtr() 
  { 
    return &theta; 
  };

  /** \brief
   * Returns pointer to \ref theta1.
   */
  double *getTheta1Ptr() 
  { 
    return &theta1; 
  };

  /** \brief
   * Returns pointer to \ref tau1
   */
  double *getTau1Ptr() 
  { 
    return &tau1; 
  };

  /** \brief
   * Returns pointer to \ref rhsTime.
   */
  double *getRHSTimePtr() 
  { 
    return &rhsTime; 
  };

  /** \brief
   * Returns pointer to \ref theta1.
   */
  double *getBoundaryTimePtr() 
  { 
    return &boundaryTime; 
  };

  /** \brief
   * Returns \ref boundaryFct;
   */
  AbstractFunction<double, WorldVector<double> > *getBoundaryFct() 
  {
    return boundaryFct;
  };

protected:
  /** \brief
   * Used for theta scheme.
   */
  double theta;

  /** \brief
   * theta - 1
   */
  double theta1;

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

  /** \brief
   * time for right hand side functions.
   */
  double rhsTime;

  /** \brief
   * time for boundary functions.
   */
  double boundaryTime;

  /** \brief
   * Pointer to boundary function. Needed for initial problem.
   */
  AbstractFunction<double, WorldVector<double> > *boundaryFct;
};

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

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

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

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

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

  // ===== create instationary problem =====
  Vecheat *vecheat = NEW Vecheat(vecheatSpace);
  vecheat->initialize(INIT_ALL);

  // create adapt info
  AdaptInfo *adaptInfo = 
    NEW AdaptInfo("vecheat->adapt",
		  vecheatSpace->getNumComponents());

  // create initial adapt info
  AdaptInfo *adaptInfoInitial = 
    NEW AdaptInfo("vecheat->initial->adapt", 
		  vecheatSpace->getNumComponents());

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

  double fac = (*(vecheat->getThetaPtr())) != 0 ? 1.0 : 1.0e-3;
  int nRefine = 0, dim = vecheatSpace->getMesh(0)->getDim();;
  GET_PARAMETER(0, vecheatSpace->getMesh(0)->getName() 
		+ "->global refinements", "%d", &nRefine);
  if((*(vecheat->getThetaPtr())) == 0.5) {
    (*(adaptInfo->getTimestepPtr())) *= 
      fac * pow(2.0, ((double) -nRefine)/dim);
  } else {
    (*(adaptInfo->getTimestepPtr())) *= 
      fac * pow(2.0, -nRefine);
  }


  // ===== create boundary functions =====
  G *boundaryFct = NEW G;
  boundaryFct->setTimePtr(vecheat->getBoundaryTimePtr());
  vecheat->setBoundaryFct(boundaryFct);

  vecheatSpace->addDirichletBC(DIRICHLET, 0, boundaryFct);
  vecheatSpace->addDirichletBC(DIRICHLET, 1, boundaryFct);

  // ===== create rhs functions =====
  F *rhsFct0 = NEW F(vecheatSpace->getFESpace(0)->getBasisFcts()->getDegree());
  rhsFct0->setTimePtr(vecheat->getRHSTimePtr());
  F *rhsFct1 = NEW F(vecheatSpace->getFESpace(1)->getBasisFcts()->getDegree());
  rhsFct1->setTimePtr(vecheat->getRHSTimePtr());

  // ===== create operators =====
  double one = 1.0;
  double zero = 0.0;

  // create laplace
  Operator *A00 = NEW Operator(Operator::MATRIX_OPERATOR | 
			       Operator::VECTOR_OPERATOR,
			       vecheatSpace->getFESpace(0),
			       vecheatSpace->getFESpace(0));

  A00->addSecondOrderTerm(new Laplace_SOT);
  A00->setUhOld(vecheat->getOldSolution()->getDOFVector(0));

  Operator *A11 = NEW Operator(Operator::MATRIX_OPERATOR | 
			       Operator::VECTOR_OPERATOR,
			       vecheatSpace->getFESpace(1),
			       vecheatSpace->getFESpace(1));

  A11->addSecondOrderTerm(new Laplace_SOT);
  A11->setUhOld(vecheat->getOldSolution()->getDOFVector(1));

  if((*(vecheat->getThetaPtr())) != 0.0) {
    vecheatSpace->addMatrixOperator(A00, 0, 0, 
				      vecheat->getThetaPtr(), 
				      &one);
    vecheatSpace->addMatrixOperator(A11, 1, 1, 
				      vecheat->getThetaPtr(), 
				      &one);
  }

  if((*(vecheat->getTheta1Ptr())) != 0.0) {
    vecheatSpace->addVectorOperator(A00, 0, 
				      vecheat->getTheta1Ptr(),
				      &zero);
    vecheatSpace->addVectorOperator(A11, 1, 
				      vecheat->getTheta1Ptr(),
				      &zero);
  }

  // create zero order operator
  Operator *C00 = NEW Operator(Operator::MATRIX_OPERATOR | 
			       Operator::VECTOR_OPERATOR,
			       vecheatSpace->getFESpace(0),
			       vecheatSpace->getFESpace(0));

  C00->addZeroOrderTerm(NEW Simple_ZOT);

  C00->setUhOld(vecheat->getOldSolution()->getDOFVector(0));

  vecheatSpace->addMatrixOperator(C00, 0, 0, 
				  vecheat->getTau1Ptr(), 
				  vecheat->getTau1Ptr());

  vecheatSpace->addVectorOperator(C00, 0,
				  vecheat->getTau1Ptr(), 
				  vecheat->getTau1Ptr());

  Operator *C11 = NEW Operator(Operator::MATRIX_OPERATOR | 
			       Operator::VECTOR_OPERATOR,
			       vecheatSpace->getFESpace(1),
			       vecheatSpace->getFESpace(1));

  C11->addZeroOrderTerm(NEW Simple_ZOT);

  C11->setUhOld(vecheat->getOldSolution()->getDOFVector(1));

  vecheatSpace->addMatrixOperator(C11, 1, 1, 
				    vecheat->getTau1Ptr(), 
				    vecheat->getTau1Ptr());

  vecheatSpace->addVectorOperator(C11, 1,
				    vecheat->getTau1Ptr(), 
				    vecheat->getTau1Ptr());

  // create RHS operator
  Operator *F0 = NEW Operator(Operator::VECTOR_OPERATOR,
			      vecheatSpace->getFESpace(0));

  F0->addZeroOrderTerm(NEW CoordsAtQP_ZOT(rhsFct0));

  vecheatSpace->addVectorOperator(F0, 0);

  Operator *F1 = NEW Operator(Operator::VECTOR_OPERATOR,
			      vecheatSpace->getFESpace(1));

  F1->addZeroOrderTerm(NEW CoordsAtQP_ZOT(rhsFct1));

  vecheatSpace->addVectorOperator(F1, 1);

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

  return errorCode;
}