Source

PetIGA / demo / NavierStokesVMS.c

Lisandro Dalcin f3d7eb8 

Lisandro Dalcin 05825b3 


Lisandro Dalcin f3d7eb8 


Nathan Collier 5de65db 
Lisandro Dalcin f3d7eb8 



Lisandro Dalcin 8192bd7 
Lisandro Dalcin f3d7eb8 







































































Lisandro Dalcin 1dac882 


Lisandro Dalcin f3d7eb8 

Lisandro Dalcin 05825b3 

Lisandro Dalcin f3d7eb8 


Lisandro Dalcin 055cf57 
Lisandro Dalcin f3d7eb8 

Lisandro Dalcin 38c7bce 



Lisandro Dalcin f3d7eb8 














Lisandro Dalcin 8192bd7 
Lisandro Dalcin 1798b28 
Lisandro Dalcin f3d7eb8 
Lisandro Dalcin 8192bd7 
Lisandro Dalcin f3d7eb8 
















































Lisandro Dalcin 05825b3 

Lisandro Dalcin f3d7eb8 




Lisandro Dalcin 1dac882 


Lisandro Dalcin f3d7eb8 

Lisandro Dalcin 05825b3 

Lisandro Dalcin f3d7eb8 



Lisandro Dalcin 38c7bce 
Lisandro Dalcin f3d7eb8 



Lisandro Dalcin 8192bd7 
Lisandro Dalcin 9bb9f2d 
Lisandro Dalcin f3d7eb8 
Lisandro Dalcin 8192bd7 
Lisandro Dalcin f3d7eb8 
























































Lisandro Dalcin 05825b3 

Lisandro Dalcin f3d7eb8 










Nathan Collier 7e031e9 


























Lisandro Dalcin f3d7eb8 

Nathan Collier 7e031e9 










Lisandro Dalcin f3d7eb8 




Nathan Collier 5de65db 

Lisandro Dalcin f3d7eb8 

Nathan Collier 5de65db 




Lisandro Dalcin f3d7eb8 









Lisandro Dalcin 05825b3 


Lisandro Dalcin f3d7eb8 










Lisandro Dalcin 05825b3 
Lisandro Dalcin 6b70e16 

Lisandro Dalcin f3d7eb8 


Lisandro Dalcin 6b70e16 
Lisandro Dalcin 1798b28 





Lisandro Dalcin f3d7eb8 
Lisandro Dalcin 6b70e16 






Lisandro Dalcin 05825b3 
Lisandro Dalcin 6b70e16 






Lisandro Dalcin f3d7eb8 




Nathan Collier 5de65db 
Lisandro Dalcin f3d7eb8 


Lisandro Dalcin 073d86a 
Lisandro Dalcin f3d7eb8 
Lisandro Dalcin e1d8929 

Lisandro Dalcin f3d7eb8 

Lisandro Dalcin 073d86a 
Lisandro Dalcin e1d8929 

Lisandro Dalcin f3d7eb8 

Lisandro Dalcin 073d86a 
Lisandro Dalcin f3d7eb8 
Lisandro Dalcin e1d8929 

Lisandro Dalcin f3d7eb8 









Lisandro Dalcin 66a3a54 


Lisandro Dalcin f3d7eb8 











Lisandro Dalcin 1798b28 
Lisandro Dalcin f3d7eb8 





Lisandro Dalcin 16f254a 
Lisandro Dalcin 1798b28 
Lisandro Dalcin ce59474 
Lisandro Dalcin 16f254a 
Lisandro Dalcin ce59474 
Lisandro Dalcin f3d7eb8 



Lisandro Dalcin 05825b3 
Lisandro Dalcin f3d7eb8 

  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
443
444
#include "petiga.h"

PetscLogEvent NS_Residual = 0;
PetscLogEvent NS_Tangent  = 0;

typedef struct {
  PetscReal nu;
  PetscScalar fx,fy,fz;
  IGA iga;
} AppCtx;

#undef __FUNCT__
#define __FUNCT__ "Tau"
void Tau(PetscReal J[3][3],
         PetscReal dt,PetscScalar u[],PetscReal nu,
         PetscScalar *tauM,PetscScalar *tauC)
{
  PetscReal C_I = 1.0/12.0;

  PetscInt i,j,k;

  PetscReal G[3][3] = {{0}};
  for (i=0;i<3;i++)
  for (j=0;j<3;j++)
  for (k=0;k<3;k++)
    G[i][j] += J[i][k]*J[j][k];
  PetscReal g[3] = {0};
  for (i=0;i<3;i++)
  for (j=0;j<3;j++)
    g[i] += J[i][j];

  PetscReal G_G = 0;
  for (i=0;i<3;i++)
  for (j=0;j<3;j++)
    G_G += G[i][j]*G[i][j];
  PetscReal g_g = 0;
  for (i=0;i<3;i++)
    g_g += g[i]*g[i];

  PetscScalar u_G_u = 0;
  for (i=0;i<3;i++)
  for (j=0;j<3;j++)
    u_G_u += u[i]*G[i][j]*u[j];

  // Eqn 63
  *tauM = 4/(dt*dt) + u_G_u + C_I * nu*nu * G_G;
  *tauM = 1/sqrt(*tauM);
  // Eqn 64
  *tauC = (*tauM) * g_g;
  *tauC = 1/(*tauC);
}

#undef __FUNCT__
#define __FUNCT__ "FineScale"
void FineScale(const AppCtx *user,PetscScalar tauM,PetscScalar tauC,
               PetscScalar ux,
               PetscScalar ux_t,
               PetscScalar ux_x ,PetscScalar ux_y, PetscScalar ux_z,
               PetscScalar ux_xx,PetscScalar ux_yy,PetscScalar ux_zz,
               PetscScalar uy,
               PetscScalar uy_t,
               PetscScalar uy_x ,PetscScalar uy_y, PetscScalar uy_z,
               PetscScalar uy_xx,PetscScalar uy_yy,PetscScalar uy_zz,
               PetscScalar uz,
               PetscScalar uz_t,
               PetscScalar uz_x ,PetscScalar uz_y, PetscScalar uz_z,
               PetscScalar uz_xx,PetscScalar uz_yy,PetscScalar uz_zz,
               PetscScalar p,PetscScalar p_x,PetscScalar p_y,PetscScalar p_z,
               PetscScalar *ux_s,PetscScalar *uy_s,PetscScalar *uz_s,PetscScalar *p_s)
{
  // Eqn 61
  (*ux_s) = ux_t + (ux*ux_x + uy*ux_y + uz*ux_z) + p_x - user->nu*(ux_xx + ux_yy + ux_zz) - user->fx;
  (*uy_s) = uy_t + (ux*uy_x + uy*uy_y + uz*uy_z) + p_y - user->nu*(uy_xx + uy_yy + uy_zz) - user->fy;
  (*uz_s) = uz_t + (ux*uz_x + uy*uz_y + uz*uz_z) + p_z - user->nu*(uz_xx + uz_yy + uz_zz) - user->fz;
  // Eqn 62
  (*p_s)  = ux_x + uy_y + uz_z;
  // Eqn 58
  (*ux_s) *= -tauM;
  (*uy_s) *= -tauM;
  (*uz_s) *= -tauM;
  (*p_s)  *= -tauC;
}


#undef  __FUNCT__
#define __FUNCT__ "Residual"
PetscErrorCode Residual(IGAPoint pnt,PetscReal dt,
                        PetscReal shift,const PetscScalar *V,
                        PetscReal t,const PetscScalar *U,
                        PetscScalar *Re,void *ctx)
{
  //PetscLogEventBegin(NS_Residual,0,0,0,0);

  AppCtx *user = (AppCtx *)ctx;
  PetscReal nu = user->nu;

  PetscScalar u_t[4],u[4];
  PetscScalar grad_u[4][3];
  PetscScalar der2_u[4][3][3];
  IGAPointFormValue(pnt,V,&u_t[0]);
  IGAPointFormValue(pnt,U,&u[0]);
  IGAPointFormGrad (pnt,U,&grad_u[0][0]);
  IGAPointFormHess (pnt,U,&der2_u[0][0][0]);

  PetscScalar ux=u[0],ux_t=u_t[0];
  PetscScalar uy=u[1],uy_t=u_t[1];
  PetscScalar uz=u[2],uz_t=u_t[2];
  PetscScalar p =u[3];

  PetscScalar ux_x=grad_u[0][0],ux_y=grad_u[0][1],ux_z=grad_u[0][2];
  PetscScalar uy_x=grad_u[1][0],uy_y=grad_u[1][1],uy_z=grad_u[1][2];
  PetscScalar uz_x=grad_u[2][0],uz_y=grad_u[2][1],uz_z=grad_u[2][2];
  PetscScalar p_x =grad_u[3][0],p_y =grad_u[3][1],p_z =grad_u[3][2];

  PetscScalar ux_xx=der2_u[0][0][0],ux_yy=der2_u[0][1][1],ux_zz=der2_u[0][2][2];
  PetscScalar uy_xx=der2_u[1][0][0],uy_yy=der2_u[1][1][1],uy_zz=der2_u[1][2][2];
  PetscScalar uz_xx=der2_u[2][0][0],uz_yy=der2_u[2][1][1],uz_zz=der2_u[2][2][2];

  PetscReal InvGradMap[3][3];
  IGAPointFormGradMap(pnt,NULL,&InvGradMap[0][0]);
  PetscScalar tauM,tauC;
  Tau(InvGradMap,dt,u,nu,&tauM,&tauC);
  PetscScalar ux_s,uy_s,uz_s,p_s;
  FineScale(user,tauM,tauC,
            ux,ux_t,ux_x,ux_y,ux_z,ux_xx,ux_yy,ux_zz,
            uy,uy_t,uy_x,uy_y,uy_z,uy_xx,uy_yy,uy_zz,
            uz,uz_t,uz_x,uz_y,uz_z,uz_xx,uz_yy,uz_zz,
            p,p_x,p_y,p_z,
            &ux_s,&uy_s,&uz_s,&p_s);

  PetscReal  *N0 = pnt->shape[0];
  PetscReal (*N1)[3] = (PetscReal (*)[3]) pnt->shape[1];

  PetscScalar (*R)[4] = (PetscScalar (*)[4])Re;
  PetscInt a,nen=pnt->nen;
  for (a=0; a<nen; a++) {
    PetscReal Na    = N0[a];
    PetscReal Na_x  = N1[a][0];
    PetscReal Na_y  = N1[a][1];
    PetscReal Na_z  = N1[a][2];
    /* ----- */
    PetscScalar Rux,Ruy,Ruz,Rp;
    // -L(W)
    Rux = -Na*user->fx;
    Ruy = -Na*user->fy;
    Ruz = -Na*user->fz;
    Rp  = 0.0;
    // += B_1(W,U)
    Rux += Na*ux_t - Na_x*p + nu*( Na_x*( ux_x + ux_x ) + Na_y*( ux_y + uy_x ) + Na_z*( ux_z + uz_x ) );
    Ruy += Na*uy_t - Na_y*p + nu*( Na_x*( uy_x + ux_y ) + Na_y*( uy_y + uy_y ) + Na_z*( uy_z + uz_y ) );
    Ruz += Na*uz_t - Na_z*p + nu*( Na_x*( uz_x + ux_z ) + Na_y*( uz_y + uy_z ) + Na_z*( uz_z + uz_z ) );
    Rp  += Na*( ux_x + uy_y + uz_z );
    // += Btilde_1(W,U')
    Rux += - ( Na_x*p_s );
    Ruy += - ( Na_y*p_s );
    Ruz += - ( Na_z*p_s );
    Rp  += - ( Na_x*ux_s + Na_y*uy_s + Na_z*uz_s );
    // += B_2(W,U,U+U') [after integration by parts, div(u+u')=0]
    Rux += + Na * ( (ux+ux_s)*ux_x + (uy+uy_s)*ux_y + (uz+uz_s)*ux_z );
    Ruy += + Na * ( (ux+ux_s)*uy_x + (uy+uy_s)*uy_y + (uz+uz_s)*uy_z );
    Ruz += + Na * ( (ux+ux_s)*uz_x + (uy+uy_s)*uz_y + (uz+uz_s)*uz_z );
    // += B_2(W,U',U+U')
    Rux += - ( Na_x*ux_s*(ux+ux_s) + Na_y*ux_s*(uy+uy_s) + Na_z*ux_s*(uz+uz_s) );
    Ruy += - ( Na_x*uy_s*(ux+ux_s) + Na_y*uy_s*(uy+uy_s) + Na_z*uy_s*(uz+uz_s) );
    Ruz += - ( Na_x*uz_s*(ux+ux_s) + Na_y*uz_s*(uy+uy_s) + Na_z*uz_s*(uz+uz_s) );
    /* ----- */
    R[a][0] = Rux;
    R[a][1] = Ruy;
    R[a][2] = Ruz;
    R[a][3] = Rp;
  }

  //PetscLogEventEnd(NS_Residual,0,0,0,0);
  return 0;
}

#undef  __FUNCT__
#define __FUNCT__ "Tangent"
PetscErrorCode Tangent(IGAPoint pnt,PetscReal dt,
                       PetscReal shift,const PetscScalar *V,
                       PetscReal t,const PetscScalar *U,
                       PetscScalar *Ke,void *ctx)
{
  //PetscLogEventBegin(NS_Tangent,0,0,0,0);

  AppCtx *user = (AppCtx *)ctx;
  PetscReal nu = user->nu;

  PetscScalar u[4];
  IGAPointFormValue(pnt,U,&u[0]);
  PetscScalar ux=u[0];
  PetscScalar uy=u[1];
  PetscScalar uz=u[2];

  PetscReal InvGradMap[3][3];
  IGAPointFormGradMap(pnt,0,&InvGradMap[0][0]);
  PetscScalar tauM,tauC;
  Tau(InvGradMap,dt,u,nu,&tauM,&tauC);

  PetscReal *N0 = pnt->shape[0];
  PetscReal (*N1)[3] = (PetscReal (*)[3]) pnt->shape[1];

  PetscInt a,b,nen=pnt->nen;
  PetscScalar (*K)[4][nen][4] = (PetscScalar (*)[4][nen][4])Ke;
  for (a=0; a<nen; a++) {
    PetscReal Na    = N0[a];
    PetscReal Na_x  = N1[a][0];
    PetscReal Na_y  = N1[a][1];
    PetscReal Na_z  = N1[a][2];
    for (b=0; b<nen; b++) {
      PetscReal Nb    = N0[b];
      PetscReal Nb_x  = N1[b][0];
      PetscReal Nb_y  = N1[b][1];
      PetscReal Nb_z  = N1[b][2];
      /* ----- */
      PetscInt    i,j;
      PetscScalar T[4][4];
      PetscScalar Tii =
        (+ shift * Na * Nb
         + Na * (ux * Nb_x + uy * Nb_y + uz * Nb_z)
         + nu * (Na_x * Nb_x + Na_y * Nb_y + Na_z * Nb_z)

         + tauM * (ux * Na_x + uy * Na_y + uz * Na_z) *
         /**/     (shift * Nb + (ux * Nb_x + uy * Nb_y + uz * Nb_z))
         );
      T[0][0] = /*+ Na * Nb * ux_x*/  +  nu * Na_x * Nb_x  +  tauC * Na_x * Nb_x;
      T[0][1] = /*+ Na * Nb * ux_y*/  +  nu * Na_y * Nb_x  +  tauC * Na_x * Nb_y;
      T[0][2] = /*+ Na * Nb * ux_z*/  +  nu * Na_z * Nb_x  +  tauC * Na_x * Nb_z;
      //
      T[1][0] = /*+ Na * Nb * uy_x*/  +  nu * Na_x * Nb_y  +  tauC * Na_y * Nb_x;
      T[1][1] = /*+ Na * Nb * uy_y*/  +  nu * Na_y * Nb_y  +  tauC * Na_y * Nb_y;
      T[1][2] = /*+ Na * Nb * uy_z*/  +  nu * Na_z * Nb_y  +  tauC * Na_y * Nb_z;
      //
      T[2][0] = /*+ Na * Nb * uz_x*/  +  nu * Na_x * Nb_z  +  tauC * Na_z * Nb_x;
      T[2][1] = /*+ Na * Nb * uz_y*/  +  nu * Na_y * Nb_z  +  tauC * Na_z * Nb_y;
      T[2][2] = /*+ Na * Nb * uz_z*/  +  nu * Na_z * Nb_z  +  tauC * Na_z * Nb_z;
      T[0][0] += Tii;
      T[1][1] += Tii;
      T[2][2] += Tii;
      // G as in Eq. (104)
      T[0][3] = - Na_x * Nb  +  tauM * (ux * Na_x + uy * Na_y + uz * Na_z) * Nb_x;
      T[1][3] = - Na_y * Nb  +  tauM * (ux * Na_x + uy * Na_y + uz * Na_z) * Nb_y;
      T[2][3] = - Na_z * Nb  +  tauM * (ux * Na_x + uy * Na_y + uz * Na_z) * Nb_z;
      // D as in Eq. (106)
      T[3][0] = + Na * Nb_x  +  tauM * Na_x * (shift * Nb + (ux * Nb_x + uy * Nb_y + uz * Nb_z));
      T[3][1] = + Na * Nb_y  +  tauM * Na_y * (shift * Nb + (ux * Nb_x + uy * Nb_y + uz * Nb_z));
      T[3][2] = + Na * Nb_z  +  tauM * Na_z * (shift * Nb + (ux * Nb_x + uy * Nb_y + uz * Nb_z));
      // L as in Eq. (108)
      T[3][3] = + tauM * (Na_x * Nb_x + Na_y * Nb_y + Na_z * Nb_z);
      /* ----- */
      for (i=0;i<4;i++)
        for (j=0;j<4;j++)
          K[a][i][b][j] += T[i][j];
    }
  }

  //PetscLogEventEnd(NS_Tangent,0,0,0,0);
  return 0;
}

#undef __FUNCT__
#define __FUNCT__ "FormInitialCondition"
PetscErrorCode FormInitialCondition(AppCtx *user,IGA iga,const char datafile[],PetscReal t,Vec U)
{
  PetscErrorCode ierr;
  PetscFunctionBegin;
  ierr = VecZeroEntries(U);CHKERRQ(ierr);

  if (datafile[0] != 0) { /* initial condition from datafile */
    MPI_Comm comm;
    PetscViewer viewer;
    ierr = PetscObjectGetComm((PetscObject)U,&comm);CHKERRQ(ierr);
    ierr = PetscViewerBinaryOpen(comm,datafile,FILE_MODE_READ,&viewer);CHKERRQ(ierr);
    ierr = VecLoad(U,viewer);CHKERRQ(ierr);
    ierr = PetscViewerDestroy(&viewer);
  } else { 
    DM da;
    ierr = IGACreateNodeDM(iga,4,&da);CHKERRQ(ierr);
    DMDALocalInfo info;
    ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
    PetscScalar ****u;
    ierr = DMDAVecGetArrayDOF(da,U,&u);CHKERRQ(ierr);
    
    PetscScalar H    = 2;//user->Ly;
    PetscScalar visc = user->nu;
    PetscScalar dpdx = user->fx;
    PetscInt  i,j,k;
    PetscReal jmax = (info.my-1);
    for(k=info.zs;k<info.zs+info.zm;k++){
      for(j=info.ys;j<info.ys+info.ym;j++){
	for(i=info.xs;i<info.xs+info.xm;i++){
	  PetscReal   y = (j/jmax)*H;
	  PetscScalar ux = 1/(2*visc) * dpdx * y*(H-y);
	  u[k][j][i][0] = ux;
	}
      }
    }

    ierr = DMDAVecRestoreArrayDOF(da,U,&u);CHKERRQ(ierr);
    ierr = DMDestroy(&da);CHKERRQ(ierr);

    PetscReal v; Vec pert;
    ierr = VecDuplicate(U,&pert);CHKERRQ(ierr);
    ierr = VecSetRandom(pert,NULL);CHKERRQ(ierr);
    ierr = VecMax(U,NULL,&v);CHKERRQ(ierr);
    ierr = VecAXPY(U,0.05*v,pert);CHKERRQ(ierr);
    ierr = VecStrideSet(U,3,0.0);CHKERRQ(ierr);
    ierr = VecDestroy(&pert);CHKERRQ(ierr);
  }
  PetscFunctionReturn(0);
}

#undef __FUNCT__
#define __FUNCT__ "OutputMonitor"
PetscErrorCode OutputMonitor(TS ts,PetscInt it_number,PetscReal c_time,Vec U,void *mctx)
{
  PetscFunctionBegin;
  PetscErrorCode ierr;
  AppCtx *user = (AppCtx *)mctx;
  char           filename[256];
  sprintf(filename,"./nsvms%d.dat",it_number);
  ierr = IGAWriteVec(user->iga,U,filename);CHKERRQ(ierr);
  PetscFunctionReturn(0);
}

#undef __FUNCT__
#define __FUNCT__ "main"
int main(int argc, char *argv[]) {

  PetscErrorCode  ierr;
  ierr = PetscInitialize(&argc,&argv,0,0);CHKERRQ(ierr);

  ierr = PetscLogEventRegister("NS_Residual",IGA_CLASSID,&NS_Residual);CHKERRQ(ierr);
  ierr = PetscLogEventRegister("NS_Tangent", IGA_CLASSID,&NS_Tangent );CHKERRQ(ierr);

  AppCtx user;
  user.nu = 1.47200e-4;
  user.fx = 3.37204e-3;
  user.fy = 0.0;
  user.fz = 0.0;

  PetscReal PI = 4.0*atan(1.0);
  PetscReal Lx = 2.0*PI;
  PetscReal Ly = 2.0;
  PetscReal Lz = 2.0/3.0*PI;

  PetscInt N[3] = {16,16,16}, nN = 3;
  PetscInt p[3] = { 2, 2, 2}, np = 3;
  PetscInt C[3] = {-1,-1,-1}, nC = 3;
  PetscBool output = PETSC_FALSE;
  PetscBool monitor = PETSC_FALSE;
  char initial[PETSC_MAX_PATH_LEN] = {0};
  ierr = PetscOptionsBegin(PETSC_COMM_WORLD,"","NavierStokesVMS Options","IGA");CHKERRQ(ierr);
  ierr = PetscOptionsIntArray("-N","number of elements",     __FILE__,N,&nN,NULL);CHKERRQ(ierr);
  ierr = PetscOptionsIntArray("-p","polynomial order",       __FILE__,p,&np,NULL);CHKERRQ(ierr);
  ierr = PetscOptionsIntArray("-C","global continuity order",__FILE__,C,&nC,NULL);CHKERRQ(ierr);
  ierr = PetscOptionsString("-ns_initial","Load initial solution from file",__FILE__,initial,initial,sizeof(initial),NULL);CHKERRQ(ierr);
  ierr = PetscOptionsBool  ("-ns_output","Enable output files",__FILE__,output,&output,NULL);CHKERRQ(ierr);
  ierr = PetscOptionsBool  ("-ns_monitor","Compute and show statistics of solution",__FILE__,monitor,&monitor,NULL);CHKERRQ(ierr);
  ierr = PetscOptionsEnd();CHKERRQ(ierr);
  if (nN == 1) N[2] = N[1] = N[0]; if (nN == 2) N[2] = N[0];
  if (np == 1) p[2] = p[1] = p[0]; if (np == 2) p[2] = p[0];
  if (nC == 1) C[2] = C[1] = C[0]; if (nC == 2) C[2] = C[0];
  if (C[0] == -1) C[0] = p[0]-1;
  if (C[1] == -1) C[1] = p[1]-1;
  if (C[2] == -1) C[2] = p[2]-1;
  PetscInt i;
  for (i=0; i<3; i++) {
    if (p[i] < 2 || C[i] < 1) /* Problem requires a p>=2 C1 basis */
      SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,
              "Problem requires minimum of p = 2 and C = 1");
    if (p[i] <= C[i])         /* Check C < p */
      SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,
              "Discretization inconsistent: polynomial order must be greater than degree of continuity");
  }

  IGA iga;
  ierr = IGACreate(PETSC_COMM_WORLD,&iga);CHKERRQ(ierr);
  ierr = IGASetDim(iga,3);CHKERRQ(ierr);
  ierr = IGASetDof(iga,4);CHKERRQ(ierr);
  user.iga = iga;

  IGAAxis axis0;
  ierr = IGAGetAxis(iga,0,&axis0);CHKERRQ(ierr);
  ierr = IGAAxisSetDegree(axis0,p[0]);CHKERRQ(ierr);
  ierr = IGAAxisSetPeriodic(axis0,PETSC_TRUE);CHKERRQ(ierr);
  //ierr = IGAAxisInitUniform(axis0,N,-0.5*Lx,0.5*Lx);CHKERRQ(ierr);
  ierr = IGAAxisInitUniform(axis0,N[0],0.0,Lx,C[0]);CHKERRQ(ierr);
  IGAAxis axis1;
  ierr = IGAGetAxis(iga,1,&axis1);CHKERRQ(ierr);
  ierr = IGAAxisSetDegree(axis1,p[1]);CHKERRQ(ierr);
  //ierr = IGAAxisInitUniform(axis1,N,-0.5*Ly,0.5*Ly);CHKERRQ(ierr);
  ierr = IGAAxisInitUniform(axis1,N[1],0.0,Ly,C[1]);CHKERRQ(ierr);
  IGAAxis axis2;
  ierr = IGAGetAxis(iga,2,&axis2);CHKERRQ(ierr);
  ierr = IGAAxisSetDegree(axis2,p[2]);CHKERRQ(ierr);
  ierr = IGAAxisSetPeriodic(axis2,PETSC_TRUE);CHKERRQ(ierr);
  //ierr = IGAAxisInitUniform(axis2,N,-0.5*Lz,0.5*Lz);CHKERRQ(ierr);
  ierr = IGAAxisInitUniform(axis2,N[2],0.0,Lz,C[2]);CHKERRQ(ierr);

  IGABoundary bnd;
  PetscInt dir=1,side,field;
  for (side=0;side<2;side++) {
    ierr = IGAGetBoundary(iga,dir,side,&bnd);CHKERRQ(ierr);
    for (field=0;field<3;field++) {
      ierr = IGABoundarySetValue(bnd,field,0.0);CHKERRQ(ierr);
    }
  }

  ierr = IGASetFromOptions(iga);CHKERRQ(ierr);
  ierr = IGASetUp(iga);CHKERRQ(ierr);

  ierr = IGASetUserIFunction(iga,Residual,&user);CHKERRQ(ierr);
  ierr = IGASetUserIJacobian(iga,Tangent,&user);CHKERRQ(ierr);

  TS ts;
  ierr = IGACreateTS(iga,&ts);CHKERRQ(ierr);
  ierr = TSSetDuration(ts,1000000,1000.0);CHKERRQ(ierr);
  ierr = TSSetTimeStep(ts,1.0e-2);CHKERRQ(ierr);

  ierr = TSSetType(ts,TSALPHA);CHKERRQ(ierr);
  ierr = TSAlphaSetRadius(ts,0.5);CHKERRQ(ierr);

  if (output) {
    ierr = TSMonitorSet(ts,OutputMonitor,&user,NULL);CHKERRQ(ierr);
  }
  ierr = TSSetFromOptions(ts);CHKERRQ(ierr);

  PetscReal t=0; Vec U;
  ierr = IGACreateVec(iga,&U);CHKERRQ(ierr);
  ierr = FormInitialCondition(&user,iga,initial,t,U);CHKERRQ(ierr);
#if PETSC_VERSION_LE(3,3,0)
  ierr = TSSolve(ts,U,NULL);CHKERRQ(ierr);
#else
  ierr = TSSolve(ts,U);CHKERRQ(ierr);
#endif

  ierr = VecDestroy(&U);CHKERRQ(ierr);
  ierr = TSDestroy(&ts);CHKERRQ(ierr);
  ierr = IGADestroy(&iga);CHKERRQ(ierr);

  ierr = PetscFinalize();CHKERRQ(ierr);
  return 0;
}
Tip: Filter by directory path e.g. /media app.js to search for public/media/app.js.
Tip: Use camelCasing e.g. ProjME to search for ProjectModifiedEvent.java.
Tip: Filter by extension type e.g. /repo .js to search for all .js files in the /repo directory.
Tip: Separate your search with spaces e.g. /ssh pom.xml to search for src/ssh/pom.xml.
Tip: Use ↑ and ↓ arrow keys to navigate and return to view the file.
Tip: You can also navigate files with Ctrl+j (next) and Ctrl+k (previous) and view the file with Ctrl+o.
Tip: You can also navigate files with Alt+j (next) and Alt+k (previous) and view the file with Alt+o.