NavierStokesCahnHilliard.cc 14.2 KB
Newer Older
1
2
3
4
#include "NavierStokesCahnHilliard.h"
// #include "Views.h"
#include "SignedDistFunctors.h"
#include "PhaseFieldConvert.h"
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
using namespace AMDiS;

NavierStokesCahnHilliard::NavierStokesCahnHilliard(const std::string &name_, bool createProblem) :
  super(name_, createProblem),
  useMobility(false),
  doubleWell(0),
  gamma(1.0),
  eps(0.1),
  minusEps(-0.1),
  epsInv(10.0),
  minusEpsInv(-10.0),
  epsSqr(0.01),
  minusEpsSqr(-0.01),
  multiPhase(NULL),
  densityPhase(NULL),
  viscosityPhase(NULL),
  viscosity1(1.0),
  viscosity2(1.0),
  density1(1.0),
  density2(1.0),    
  refFunction(NULL),
  refinement(NULL),
  sigma(0.0),
  surfaceTension(0.0)
  {
    dow = Global::getGeo(WORLD);
  Initfile::get(name + "->viscosity1", viscosity1); // viscosity of fluid 1
  Initfile::get(name + "->viscosity2", viscosity2); // viscosity of fluid 2
  Initfile::get(name + "->density1", density1); // density of fluid 1
  Initfile::get(name + "->density2", density2); // density of fluid 2

  Initfile::get(name + "->non-linear term", nonLinTerm); 
  Initfile::get(name + "->theta", theta); 
  Initfile::get("adapt->timestep", tau); 
    
    
    Parameters::get(name + "->epsilon", eps); // interface width
    Parameters::get(name + "->sigma", sigma); 
    double a_fac;
    Parameters::get(name + "->a_factor", a_fac);
    a=1.0/sigma*sqr(eps)*a_fac* (std::max(density1,density2)/tau);
    b=1.0;
    ab=a*b;   
    surfaceTension = -sigma/eps * 3.0/2.0/sqrt(2.0) * a ;

  theta1 = 1.0-theta;
  minusTheta1 = -theta1;

  force.set(0.0);
  Initfile::get(name + "->force", force);
  
  // parameters for CH
  Parameters::get(name + "->use mobility", useMobility); // mobility
  Parameters::get(name + "->gamma", gamma); // mobility

  
  // type of double well: 0= [0,1], 1= [-1,1]
  Parameters::get(name + "->double-well type", doubleWell); 
  
  // transformation of the parameters
  minusEps = -eps;
  minus1 = -1.0;
  epsInv = 1.0/eps;
  minusEpsInv = -epsInv;
  epsSqr = sqr(eps);
  minusEpsSqr = -epsSqr*b;
  
  
}


NavierStokesCahnHilliard::~NavierStokesCahnHilliard() 
{ FUNCNAME("NavierStokesCahnHilliard::~NavierStokesCahnHilliard()");   
  delete densityPhase;
  delete viscosityPhase;
}


void NavierStokesCahnHilliard::setTime(AdaptInfo* adaptInfo) 
  {
    super::setTime(adaptInfo);
    delta = gamma*adaptInfo->getTimestep();
  }
  
  


void NavierStokesCahnHilliard::initData()
{ FUNCNAME("NavierStokes_TH_MultiPhase::initTimeInterface()");

  #ifdef HAVE_PARALLEL_DOMAIN_AMDIS
    string initFileStr = name + "->space->solver", solverType = "";
    Parameters::get(initFileStr, solverType);
    if (solverType == "petsc-nsch") {
      PetscSolverNSCH* solver = dynamic_cast<PetscSolverNSCH*>(prob->getSolver()); 
      if (solver)
      {
	solver->setChData(&eps, &delta);
	solver->setStokesData( getInvTau(), prob->getSolution(), &viscosity1, &viscosity2, &density1, &density2);
      }
    }
    #endif

  densityPhase = new DOFVector<double>(prob->getFeSpace(0), "densityPhase");
  viscosityPhase = new DOFVector<double>(prob->getFeSpace(0), "viscosityPhase");

  densityPhase->set(density1);
  viscosityPhase->set(viscosity1);
  
   if (velocity == NULL)
    velocity = new DOFVector<WorldVector<double> >(getFeSpace(0), "velocity");
 
  //fileWriter = new FileVectorWriter(name + "->velocity->output", getFeSpace()->getMesh(), velocity);

  super::initData();
};


void NavierStokesCahnHilliard::closeTimestep(AdaptInfo *adaptInfo)
{ super::closeTimestep(adaptInfo);

}

void NavierStokesCahnHilliard::initTimestep(AdaptInfo *adaptInfo)
{ FUNCNAME("NavierStokes_TH_MultiPhase::initTimestep()");
  
  super::initTimestep(adaptInfo);
  refinement->refine(2);
  LinearInterpolation1 dLI(density1, density2);
  LinearInterpolation1 vLI(viscosity1, viscosity2);
  transformDOF(prob->getSolution()->getDOFVector(dow+2), densityPhase, &dLI);
  transformDOF(prob->getSolution()->getDOFVector(dow+2), viscosityPhase, &vLI);
};



void NavierStokesCahnHilliard::solveInitialProblem(AdaptInfo *adaptInfo)
  {
    // meshFunction for klocal refinement around the interface of the phasefield-DOFVector
    refFunction = new PhaseFieldRefinement(prob->getMesh());

    if (getDoubleWellType() == 0) {
      refinement = new RefinementLevelDOF(
	prob->getFeSpace(),
	refFunction,
	prob->getSolution()->getDOFVector(dow+2)); // phaseField-DOFVector in [0,1]
    } else {
      refinement = new RefinementLevelDOF(
	prob->getFeSpace(),
	refFunction,
	new PhaseDOFView<double>(prob->getSolution()->getDOFVector(dow+2))); // phaseField-DOFVector in [-1,1]
    }

    // initial refinement
    refinement->refine(0);

    for (int i = 0; i < 3; i++) {
      solveInitialProblem2(adaptInfo); 	// update phaseField-DOFVector
      refinement->refine((i < 4 ? 4 : 10));	// do k steps of local refinement
    }

    // solve all initial problems of the problems added to the CouplingTimeInterface    
  }



void NavierStokesCahnHilliard::solveInitialProblem2(AdaptInfo *adaptInfo) 
{ FUNCNAME("NavierStokesCahnHilliard::solveInitialProblem()");

  Flag initFlag = initDataFromFile(adaptInfo);

  if (!initFlag.isSet(DATA_ADOPTED)) {
    int initialInterface = 0;
    Initfile::get(name + "->initial interface", initialInterface);
    double initialEps = eps;
    Initfile::get(name + "->initial epsilon", initialEps);

    if (initialInterface == 0) {
      /// horizontale Linie
      double a= 0.0, dir= -1.0;
      Initfile::get(name + "->line->pos", a);
      Initfile::get(name + "->line->direction", dir);
      prob->getSolution()->getDOFVector(1+3)->interpol(new Plane(a, dir));
    }
    else if (initialInterface == 1) {
      /// schraege Linie
      double theta = m_pi/4.0;
      prob->getSolution()->getDOFVector(1+3)->interpol(new PlaneRotation(0.0, theta, 1.0));
      transformDOFInterpolation(prob->getSolution()->getDOFVector(1+3),new PlaneRotation(0.0, -theta, -1.0), new AMDiS::Min<double>);
    }
    else if (initialInterface == 2) {
      /// Ellipse
      double a= 1.0, b= 1.0;
      Initfile::get(name + "->ellipse->a", a);
      Initfile::get(name + "->ellipse->b", b);
      prob->getSolution()->getDOFVector(1+3)->interpol(new Ellipse(a,b));
    }
    else if (initialInterface == 3) {
      /// zwei horizontale Linien
      double a= 0.75, b= 0.375;
      Initfile::get(name + "->lines->pos1", a);
      Initfile::get(name + "->lines->pos2", b);
      prob->getSolution()->getDOFVector(1+3)->interpol(new Plane(a, -1.0));
      transformDOFInterpolation(prob->getSolution()->getDOFVector(1+3),new Plane(b, 1.0), new AMDiS::Max<double>);
    }
    else if (initialInterface == 4) {
      /// Kreis
      double radius= 1.0, x=0, y=0;
      Initfile::get(name + "->circle->radius", radius);
      Initfile::get(name + "->circle->center_x", x);
      Initfile::get(name + "->circle->center_y", y);
      WorldVector<double> w;
      w[0]=x; w[1]=y+adaptInfo->getTime();
      prob->getSolution()->getDOFVector(1+3)->interpol(new Circle(radius, w));
    } else if (initialInterface == 5) {
      /// Rechteck
      double width = 0.5;
      double height = 0.3;
      WorldVector<double> center; center.set(0.5);
      Initfile::get(name + "->rectangle->width", width);
      Initfile::get(name + "->rectangle->height", height);
      Initfile::get(name + "->rectangle->center", center);
      prob->getSolution()->getDOFVector(1+3)->interpol(new Rectangle(width, height, center));
    }

    /// create phase-field from signed-dist-function
    if (doubleWell == 0) {
      transformDOF(prob->getSolution()->getDOFVector(1+3),
        new SignedDistToPhaseField(initialEps));
    } else {
      transformDOF(prob->getSolution()->getDOFVector(1+3),
        new SignedDistToCh(initialEps));
    }
  }
}


void NavierStokesCahnHilliard::fillOperators()
{ FUNCNAME("NavierStokesCahnHilliard::fillOperators()");
  MSG("NavierStokesCahnHilliard::fillOperators()\n");

  const FiniteElemSpace* feSpace = prob->getFeSpace(0);
 
  // c
  Operator *opChMnew = new Operator(feSpace,feSpace);
  opChMnew->addZeroOrderTerm(new Simple_ZOT);
  Operator *opChMold = new Operator(feSpace,feSpace);
  opChMold->addZeroOrderTerm(new VecAtQP_ZOT(prob->getSolution()->getDOFVector(1+3)));
  // -nabla*(grad(c))
  Operator *opChL = new Operator(feSpace,feSpace);
  opChL->addSecondOrderTerm(new Simple_SOT);
  
  // div(M(c)grad(mu)), with M(c)=gamma/4*(c^2-1)^2
  Operator *opChLM = new Operator(feSpace,feSpace);
  opChLM->addSecondOrderTerm(new Simple_SOT(gamma*ab));
  
  // -2*c_old^3 + 3/2*c_old^2
  Operator *opChMPowExpl = new Operator(feSpace,feSpace);
  opChMPowExpl->addZeroOrderTerm(new VecAtQP_ZOT(
    prob->getSolution()->getDOFVector(1+3),
    new Pow3Functor(-2.0)));
  if (doubleWell == 0) {
    opChMPowExpl->addZeroOrderTerm(new VecAtQP_ZOT(
      prob->getSolution()->getDOFVector(1+3),
      new Pow2Functor(3.0/2.0)));
  }

  // -3*c_old^2 * c
  Operator *opChMPowImpl = new Operator(feSpace,feSpace);
  opChMPowImpl->addZeroOrderTerm(new VecAtQP_ZOT(
    prob->getSolution()->getDOFVector(1+3),
    new Pow2Functor(-3.0)));
  if (doubleWell == 0) {
    opChMPowImpl->addZeroOrderTerm(new VecAtQP_ZOT(
      prob->getSolution()->getDOFVector(1+3),
      new AMDiS::Factor<double>(3.0)));
    opChMPowImpl->addZeroOrderTerm(new Simple_ZOT(-0.5));
  } else {
    opChMPowImpl->addZeroOrderTerm(new Simple_ZOT(1.0));
  }

  // mu + eps^2*laplace(c) + c - 3*(c_old^2)*c = -2*c_old^3 [+ BC]
  // ----------------------------------------------------------------------
  prob->addMatrixOperator(*opChMnew,0+3,0+3, &ab);              /// < phi*mu , psi >
  
  prob->addMatrixOperator(*opChMPowImpl,0+3,1+3, &b);          /// < -3*phi*c*c_old^2 , psi >
  prob->addMatrixOperator(*opChL,0+3,1+3, &minusEpsSqr);   /// < -eps^2*phi*grad(c) , grad(psi) >
  // . . . vectorOperators . . . . . . . . . . . . . . .
  prob->addVectorOperator(*opChMPowExpl,0+3, &b);            /// < -2*phi*c_old^3 , psi >

//   setAssembleMatrixOnlyOnce_butTimestepChange(0,1);
  
  // dt(c) = laplace(mu) - u*grad(c)
  // -----------------------------------
  prob->addMatrixOperator(*opChMnew,1+3,1+3, &b); /// < phi*c/tau , psi >
  prob->addMatrixOperator(*opChLM,1+3,0+3, getTau());                /// < phi*grad(mu) , grad(psi) >
  // . . . vectorOperators . . . . . . . . . . . . . . .
  prob->addVectorOperator(*opChMold,1+3, &b );   /// < phi*c^old/tau , psi >
  
  /**/
  
  
  /// Navier-Stokes part
    
  WorldVector<DOFVector<double>* > vel;
  for (unsigned i = 0; i < dow; i++){
    vel[i]= prob->getSolution()->getDOFVector(i+2); 
  }


  // fill operators for prob
  for (unsigned i = 0; i < dow; ++i) {
    
    /// < (1/tau)*rho*u'_i , psi >
    Operator *opTime = new Operator(getFeSpace(i), getFeSpace(i));
    if (density1==density2)
      opTime->addTerm(new Simple_ZOT(density1));  
    else
      opTime->addTerm(new VecAtQP_ZOT(densityPhase, NULL));
    opTime->setUhOld(prob->getSolution()->getDOFVector(i));
    prob->addMatrixOperator(*opTime, i, i, getInvTau(), getInvTau());    
    prob->addVectorOperator(*opTime, i, getInvTau(), getInvTau());    
    
    
 
    /// < u^old*grad(u_i^old) , psi >
/*    Operator *opUGradU0 = new Operator(getFeSpace(i),getFeSpace(i));
    if (density1==density2)
        opUGradU0->addTerm(new WorldVec_FOT(vel, -1.0), GRD_PHI);
    else
        opUGradU0->addTerm(new WorldVecPhase_FOT(densityPhase, vel, -1.0), GRD_PHI);
    opUGradU0->setUhOld(prob->getSolution()->getDOFVector(i));
    if (nonLinTerm == 0) {
      prob->addVectorOperator(*opUGradU0, i);
    } else {
      prob->addVectorOperator(*opUGradU0, i, &theta1, &theta1);
    }

    if (nonLinTerm == 1) {
      /// < u'*grad(u_i^old) , psi >
      for (unsigned j = 2; j < 2+dow; ++j) {
        Operator *opUGradU1 = new Operator(getFeSpace(i),getFeSpace(i));
	if (density1==density2)
	  opUGradU1->addTerm(new PartialDerivative_ZOT(prob->getSolution()->getDOFVector(i), j-2));
	else
	  opUGradU1->addTerm(new VecAndPartialDerivative_ZOT(     densityPhase,    prob->getSolution()->getDOFVector(i), j-2));
        prob->addMatrixOperator(*opUGradU1, i, j, &theta, &theta);
      }
    } else if (nonLinTerm == 2) {
      /// < u^old*grad(u'_i) , psi >
*/   for(unsigned j = 2; j < 2+dow; ++j) {
        Operator *opUGradU2 = new Operator(getFeSpace(i),getFeSpace(i));
        opUGradU2->addTerm(new Vec2ProductPartial_FOT(   densityPhase,    prob->getSolution()->getDOFVector(j-2), j-2), GRD_PHI);
        prob->addMatrixOperator(*opUGradU2, i, i);
      }
    
  /**/
    /// Diffusion-Operator (including Stress-Tensor for space-dependent viscosity
    Operator *opLaplaceUi = new Operator(getFeSpace(i), getFeSpace(i));
    opLaplaceUi->addTerm(new VecAtQP_SOT(viscosityPhase, NULL));
    prob->addMatrixOperator(*opLaplaceUi, i, i);  
  
    /// < p , d_i(psi) >
    Operator *opGradP = new Operator(getFeSpace(i),getFeSpace(2));
    opGradP->addTerm(new PartialDerivative_FOT(i, -1.0), GRD_PSI);
    prob->addMatrixOperator(*opGradP, i, 2);
    
    /// external force, i.e. gravitational force
    if (force[i]*force[i] > DBL_TOL) {
      Operator *opForce = new Operator(getFeSpace(i), getFeSpace(i));
      opForce->addZeroOrderTerm(new VecAtQP_ZOT(densityPhase, new Force(force[i])));
      prob->addVectorOperator(*opForce, i);
    }
  
    /// < d_i(u'_i) , psi >
    Operator *opDivU = new Operator(getFeSpace(2),getFeSpace(i));
    opDivU->addTerm(new PartialDerivative_FOT(i), GRD_PHI);
    prob->addMatrixOperator(*opDivU, 2, i);
    
    
    
     ///coupling Operators 
      Operator *opNuGradC = new Operator(getFeSpace(i), getFeSpace(dow+1));
      opNuGradC->addTerm(new PartialDerivative_ZOT(prob->getSolution()->getDOFVector(dow+2), i));
      prob->addMatrixOperator(opNuGradC, i, dow+1, &surfaceTension, &surfaceTension);
    
      Operator *opVGradC = new Operator(getFeSpace(dow+2), getFeSpace(i));
      opVGradC->addTerm(new PartialDerivative_ZOT(prob->getSolution()->getDOFVector(dow+2), i, b));
      prob->addMatrixOperator(opVGradC, dow+2, i,  getTau());
      /**/
     
  }
  
  /**/
};


void NavierStokesCahnHilliard::addLaplaceTerm(int i)
{ FUNCNAME("NavierStokes_TH_MultiPhase::addLaplaceTerm()");

  /// < alpha*[grad(u)+grad(u)^t] , grad(psi) >
  if (viscosity1!=viscosity2) {
    for (unsigned j = 0; j < dow; ++j) {
      Operator *opLaplaceUi1 = new Operator(getFeSpace(i), getFeSpace(j));
      opLaplaceUi1->addTerm(new VecAtQP_IJ_SOT(  viscosityPhase, NULL, j, i));      
      prob->addMatrixOperator(*opLaplaceUi1, i, j);      
    }
  }/**/
  
  /// < alpha*grad(u'_i) , grad(psi) >
  
 
};