File: chts.c

package info (click to toggle)
illuminator 0.10.0-4
  • links: PTS
  • area: main
  • in suites: etch, etch-m68k
  • size: 2,488 kB
  • ctags: 671
  • sloc: sh: 8,903; ansic: 6,285; makefile: 252
file content (639 lines) | stat: -rw-r--r-- 21,705 bytes parent folder | download
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
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
/***************************************
  $Header: /cvsroot/petscgraphics/chts.c,v 1.35 2006/01/30 04:53:58 hazelsct Exp $

  This is the Cahn Hilliard timestepping code.  It is provided here as an
  example usage of the Illuminator Distributed Visualization Library.
***************************************/

static char help[] = "Solves a nonlinear system in parallel using PETSc timestepping routines.\n\
  \n\
The 2D or 3D transient Cahn-Hilliard phase field equation is solved on a 1x1\n\
 square or 1x1x1 cube.\n\
The model is governed by the following parameters:\n\
  -twodee : obvious (default is 3-D)\n\
  -random : random initial condition (default is a box)\n\
  -kappa <kap>, where <kap> = diffusivity\n\
  -epsilon <eps>, where <eps> = interface thickness (default 1/mx)\n\
  -gamma <gam>, where <gam> = interfacial tension (default 1)\n\
  -mparam <m>, where <m> = free energy parameter, bet 0 & 1/2 for 0 stable\n\
Model parameters alpha and beta are set from epsilon and gamma according to:\n\
  alpha = gam*eps        beta=gam/eps\n\
Low-level options:\n\
  -mx <xg>, where <xg> = number of grid points in the x-direction\n\
  -my <yg>, where <yg> = number of grid points in the y-direction\n\
  -mz <zg>, where <zg> = number of grid points in the z-direction\n\
  -printg : print grid information\n\
Graphics of the contours of C are drawn by default at each timestep:\n\
  -no_contours : do not draw contour plots of solution\n\
Parallelism can be invoked based on the DA construct:\n\
  -Nx <npx>, where <npx> = number of processors in the x-direction\n\
  -Ny <npy>, where <npy> = number of processors in the y-direction\n\
  -Nz <npz>, where <npz> = number of processors in the z-direction\n\n";


/* Including cahnhill.h includes the necessary PETSc header files */
#include <stdlib.h>
#include "cahnhill.h"
#include "illuminator.h"


/* Set #if to 1 to debug, 0 otherwise */
#undef DPRINTF
#ifdef DEBUG
#define DPRINTF(fmt, args...) ierr = PetscPrintf (PETSC_COMM_WORLD, "%s: " fmt, __FUNCT__, args); CHKERRQ(ierr)
#else
#define DPRINTF(fmt, args...)
#endif


/* Routines given below in this file */
int FormInitialCondition(AppCtx*,Vec);
int UserLevelEnd(AppCtx*,Vec);
int InitializeProblem(AppCtx*,Vec*);

ISurface Surf;
IDisplay Disp;

/*++++++++++++++++++++++++++++++++++++++
  Wrapper for
  +latex+{\tt ch\_residual\_vector\_2d()} in {\tt cahnhill.c}.
  +html+ <tt>ch_residual_vector_2d()</tt> in <tt>cahnhill.c</tt>.

  int ch_ts_residual_vector_2d Usual return: zero or error.

  TS thets Timestepping context, ignored here.

  PetscScalar time Current time, also ignored.

  Vec X Current solution vector.

  Vec F Function vector to return.

  void *ptr User data pointer.
  ++++++++++++++++++++++++++++++++++++++*/

int ch_ts_residual_vector_2d
(TS thets, PetscScalar time, Vec X, Vec F, void *ptr)
{ return ch_residual_vector_2d (X, F, ptr); }


/*++++++++++++++++++++++++++++++++++++++
  Wrapper for
  +latex+{\tt ch\_residual\_vector\_3d()} in {\tt cahnhill.c}.
  +html+ <tt>ch_residual_vector_3d()</tt> in <tt>cahnhill.c</tt>.

  int ch_ts_residual_vector_3d Usual return: zero or error.

  TS thets Timestepping context, ignored here.

  PetscScalar time Current time, also ignored.

  Vec X Current solution vector.

  Vec F Function vector to return.

  void *ptr User data pointer.
  ++++++++++++++++++++++++++++++++++++++*/

int ch_ts_residual_vector_3d
(TS thets, PetscScalar time, Vec X, Vec F, void *ptr)
{ return ch_residual_vector_3d (X, F, ptr); }


/*++++++++++++++++++++++++++++++++++++++
  Wrapper for
  +latex+{\tt ch\_jacobian\_2d()} in {\tt cahnhill.c}.
  +html+ <tt>ch_jacobian_2d()</tt> in <tt>cahnhill.c</tt>.

  int ch_ts_jacobian_2d Usual return: zero or error.

  TS thets Timestepping context, ignored here.

  PetscScalar time Current time, also ignored.

  Vec X Current solution vector.

  Mat *A Place to put the new Jacobian.

  Mat *B Place to put the new conditioning matrix.

  MatStructure *flag Flag describing the volatility of the structure.

  void *ptr User data pointer.
  ++++++++++++++++++++++++++++++++++++++*/

int ch_ts_jacobian_2d (TS thets, PetscScalar time, Vec X, Mat *A, Mat *B,
		       MatStructure *flag, void *ptr) {
  return ch_jacobian_2d (X, A, B, flag, ptr); }


/*++++++++++++++++++++++++++++++++++++++
  Wrapper for
  +latex+{\tt ch\_jacobian\_3d()} in {\tt cahnhill.c}.
  +html+ <tt>ch_jacobian_3d()</tt> in <tt>cahnhill.c</tt>.

  int ch_ts_jacobian_3d Usual return: zero or error.

  TS thets Timestepping context, ignored here.

  PetscScalar time Current time, also ignored.

  Vec X Current solution vector.

  Mat *A Place to put the new Jacobian.

  Mat *B Place to put the new conditioning matrix.

  MatStructure *flag Flag describing the volatility of the structure.

  void *ptr User data pointer.
  ++++++++++++++++++++++++++++++++++++++*/

int ch_ts_jacobian_3d (TS thets, PetscScalar time, Vec X, Mat *A, Mat *B,
		       MatStructure *flag, void *ptr) {
  return ch_jacobian_3d (X, A, B, flag, ptr); }


#undef __FUNCT__
#define __FUNCT__ "ch_ts_monitor"

/*++++++++++++++++++++++++++++++++++++++
  Monitor routine which displays the current state, using
  +latex+{\tt Illuminator}'s {\tt geomview}
  +html+ <tt>Illuminator</tt>'s <tt>geomview</tt>
  front-end (unless -no_contours is used); and also saves it using
  +latex+{\tt IlluMultiSave()}
  +html+ <tt>IlluMultiSave()</tt>
  if -save_data is specified.

  int ch_ts_monitor Usual return: zero or error.

  TS thets Timestepping context, ignored here.

  int stepno Current time step number.

  PetscScalar time Current time.

  Vec X Vector of current solved field values.

  void *ptr User data pointer.
  ++++++++++++++++++++++++++++++++++++++*/

int ch_ts_monitor (TS thets, int stepno, PetscScalar time, Vec X, void *ptr)
{
  AppCtx *user;
  int temp, ierr;
  PetscReal xmin, xmax;
  PetscScalar minmax[6] = { 0.,1., 0.,1., 0.,1. };
  /* Example colors and isoquant values to pass to DATriangulate() */
  /* PetscScalar colors[16] = { 1,0,0,.5, 1,1,0,.5, 0,1,0,.5, 0,0,1,.5 };
     PetscScalar isovals[4] = { .2, .4, .6, .8 }; */

  user = (AppCtx *)ptr;

  ierr = VecMin (X, &temp, &xmin); CHKERRQ (ierr);
  ierr = VecMax (X, &temp, &xmax); CHKERRQ (ierr);
  PetscPrintf (user->comm,"Step %d, Time %g, C values from %g to %g\n",
	   stepno, time, xmin, xmax);

  if (!user->no_contours)
    {
      if (user->threedee)
	{
#ifdef GEOMVIEW
	  DPRINTF ("Starting triangulation\n",0);
	  ierr = DATriangulate (Surf, user->da, X, user->chvar, minmax,
				PETSC_DECIDE, PETSC_NULL, PETSC_NULL,
				PETSC_FALSE, PETSC_FALSE, PETSC_FALSE);
	  CHKERRQ (ierr);
	  DPRINTF ("Done, sending to Geomview\n",0);
	  ierr = GeomviewDisplayTriangulation
	    (user->comm, Surf, Disp, minmax, "Cahn-Hilliard", PETSC_TRUE);
	  CHKERRQ (ierr);
#endif

	  ierr = SurfClear (Surf); CHKERRQ (ierr);
	  DPRINTF ("Done.\n",0);
	}
      else
	{
	  ierr = VecView (X,user->theviewer); CHKERRQ (ierr);
	}
    }

  if (user->save_data)
    {
      PetscReal deltat;
      field_plot_type chtype=FIELD_SCALAR;
      char filename [100], *paramdata [4],
	*paramnames [4] = { "time", "timestep", "gamma", "kappa" };

      for (ierr=0; ierr<4; ierr++)
	paramdata[ierr] = (char *) malloc (50*sizeof(char));
      snprintf (filename, 99, "chtsout.time%.3d", stepno);
      TSGetTimeStep (thets, &deltat);
      snprintf (paramdata[0], 49, "%g", time);
      snprintf (paramdata[1], 49, "%g", deltat);
      snprintf (paramdata[2], 49, "%g", user->gamma);
      snprintf (paramdata[3], 49, "%g", user->kappa);
      DPRINTF ("Storing data using IlluMultiSave()\n",0);
      ierr = IlluMultiSave
	(PETSC_COMM_WORLD, user->da, X, filename, 1.,1.,1., &chtype, 4,
	 paramnames, paramdata, COMPRESS_INT_NONE | COMPRESS_GZIP_FAST);
      CHKERRQ (ierr);
    }

  DPRINTF ("Completed timestep monitor tasks.\n",0);

  return 0;
}


#undef __FUNCT__
#define __FUNCT__ "main"

/*++++++++++++++++++++++++++++++++++++++
  The usual main function.

  int main Returns 0 or error.

  int argc Number of args.

  char **argv The args.
  ++++++++++++++++++++++++++++++++++++++*/

int main (int argc, char **argv)
{
  TS            thets;               /* the timestepping solver */
  Vec           x;                   /* solution vector */
  AppCtx        user;                /* user-defined work context */
  int           dim, ierr;
  PetscDraw     draw;
  PetscTruth    matfree;
  PetscReal     xmin, xmax;
  int           temp;
  PetscScalar   ftime, time_total_max = 100.0; /* default max total time */
  int           steps = 0, time_steps_max = 5; /* default max timesteps */

  PetscInitialize (&argc, &argv, (char *)0, help);
  ierr = MPI_Comm_rank (PETSC_COMM_WORLD, &user.rank); CHKERRQ (ierr);
  ierr = MPI_Comm_size (PETSC_COMM_WORLD, &user.size); CHKERRQ (ierr);
  user.comm = PETSC_COMM_WORLD;
  user.mc = 1; user.chvar = 0; /* Sets number of variables, which is conc */

  /* Create user context, set problem data, create vector data structures.
     Also, set the initial condition */

  DPRINTF ("About to initialize problem\n",0);
  ierr = InitializeProblem (&user, &x); CHKERRQ (ierr);
  if (user.load_data > -1)
    steps = user.load_data;
  ierr = VecDuplicate (user.localX, &user.localF); CHKERRQ (ierr);

  /* Set up displays to show graph of the solution */

  if (!user.no_contours)
    {
      if (user.threedee)
	{
	  ierr = SurfCreate (&Surf); CHKERRQ (ierr);
#ifdef GEOMVIEW
	  ierr = GeomviewBegin (PETSC_COMM_WORLD, &Disp); CHKERRQ (ierr);
#endif
	}
      else
	{
	  ierr = PetscViewerDrawOpen (PETSC_COMM_WORLD, 0, "", PETSC_DECIDE,
				      PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE,
				      &user.theviewer); CHKERRQ (ierr);
	  ierr = PetscViewerDrawGetDraw (user.theviewer, 0, &draw);
	  CHKERRQ (ierr);
	  ierr = PetscDrawSetDoubleBuffer (draw); CHKERRQ (ierr);
	}
    }

  /* Create timestepping solver context */
  DPRINTF ("Creating timestepping context\n",0);
  ierr = TSCreate (PETSC_COMM_WORLD, &thets); CHKERRQ (ierr);
  ierr = TSSetProblemType (thets, TS_NONLINEAR); CHKERRQ (ierr);
  ierr = VecGetSize (x, &dim); CHKERRQ (ierr);
  ierr = PetscPrintf (user.comm, "global size = %d, kappa = %g, epsilon = %g, "
		      "gamma = %g, mparam = %g\n", dim, user.kappa,
		      user.epsilon, user.gamma, user.mparam); CHKERRQ (ierr);

  /* Set function evaluation routine and monitor */
  DPRINTF ("Setting RHS function\n",0);
  if (user.threedee)
    ierr = TSSetRHSFunction (thets, ch_ts_residual_vector_3d, &user);
  else
    ierr = TSSetRHSFunction (thets, ch_ts_residual_vector_2d, &user);
  CHKERRQ(ierr);
  ierr = TSSetMonitor (thets, ch_ts_monitor, &user, PETSC_NULL); CHKERRQ(ierr);

  /* This condition prevents the construction of the Jacobian if we're
     running matrix-free. */
  ierr = PetscOptionsHasName (PETSC_NULL, "-snes_mf", &matfree); CHKERRQ(ierr);

  if (!matfree)
    {
      /* Set up the Jacobian using cahnhill.c subroutines */
      DPRINTF ("Using analytical Jacobian from cahnhill.c\n",0);
      if (user.threedee) {
	ierr = ch_jacobian_alpha_3d (&user); CHKERRQ (ierr);
	ierr = TSSetRHSJacobian (thets, user.J, user.J, ch_ts_jacobian_3d,
				 &user); CHKERRQ (ierr); }
      else {
	ierr = ch_jacobian_alpha_2d (&user); CHKERRQ (ierr);
	ierr = TSSetRHSJacobian (thets, user.J, user.J, ch_ts_jacobian_2d,
				 &user); CHKERRQ (ierr);}
      /*}*/
    }

  /* Set solution vector and initial timestep (currently a fraction of the
     explicit diffusion stability criterion */
  ierr = TSSetInitialTimeStep (thets, 0.0, 0.1/(user.mx-1)/(user.mx-1));
  CHKERRQ (ierr);
  ierr = TSSetSolution (thets, x); CHKERRQ (ierr);

  /* Customize timestepping solver, default to Crank-Nicholson */
  ierr = TSSetDuration (thets, time_steps_max, time_total_max); CHKERRQ (ierr);
  ierr = TSSetType (thets, TS_CRANK_NICHOLSON); CHKERRQ (ierr);
  ierr = TSSetFromOptions (thets); CHKERRQ (ierr);

  /* Run the solver */
  DPRINTF ("Running the solver\n",0);
  ierr = TSStep (thets, &steps, &ftime); CHKERRQ (ierr);

  /* Final clean-up */
  DPRINTF ("Cleaning up\n",0);
  if (!user.no_contours)
    {
      if (user.threedee)
	{
#ifdef GEOMVIEW
	  ierr = GeomviewEnd (PETSC_COMM_WORLD, Disp); CHKERRQ (ierr);
#endif
	  ierr = ISurfDestroy (Surf); CHKERRQ (ierr);
	}
      else
	{
	  ierr = PetscViewerDestroy (user.theviewer); CHKERRQ (ierr);
	}
    }
  ierr = VecDestroy (user.localX); CHKERRQ (ierr);
  ierr = VecDestroy (x); CHKERRQ (ierr);
  ierr = VecDestroy (user.localF); CHKERRQ (ierr);
  ierr = TSDestroy (thets); CHKERRQ (ierr);  
  ierr = PetscFree (user.label); CHKERRQ (ierr);
  ierr = DADestroy (user.da); CHKERRQ (ierr);

  PetscFinalize ();
  printf ("Game over man!\n");
  return 0;
}


#undef __FUNCT__
#define __FUNCT__ "FormInitialCondition"

/*++++++++++++++++++++++++++++++++++++++
  Like it says, put together the initial condition.

  int FormInitialCondition Returns zero or error.

  AppCtx *user The user context structure.

  Vec X Vector in which to place the initial condition.
  ++++++++++++++++++++++++++++++++++++++*/

int FormInitialCondition (AppCtx *user, Vec X)
{
  int     i,j,k, row, mc, chvar, mx,my,mz, ierr, xs,ys,zs, xm,ym,zm;
  int     gxm,gym,gzm, gxs,gys,gzs;
  PetscScalar  mparam;
  PetscScalar  *x;
  Vec     localX = user->localX;
  PetscRandom rand;

  mc = user->mc;
  chvar = user->chvar;
  mx = user->mx; my = user->my; mz = user->mz;
  mparam = user->mparam;

  /* Get a pointer to vector data.
       - For default PETSc vectors, VecGetArray() returns a pointer to
         the data array.  Otherwise, the routine is implementation dependent.
       - You MUST call VecRestoreArray() when you no longer need access to
         the array. */
  if (user->random)
    {
      DPRINTF ("Setting initial condition as random numbers\n",0);
      ierr = PetscRandomCreate (user->comm, &rand);
      CHKERRQ (ierr);
      ierr = PetscRandomSetInterval (rand, 0.35, 0.65); CHKERRQ (ierr);
      ierr = VecSetRandom (X, rand); CHKERRQ (ierr);
      ierr = PetscRandomDestroy (rand); CHKERRQ (ierr);
    }
  else if (user->load_data > -1)
    {
      int usermetacount;
      char basename [1000], **usermetanames, **usermetadata;
      sprintf (basename, "chtsout.time%.3d", user->load_data);
      DPRINTF ("Loading data for timestep %d, basename %s\n", user->load_data,
	       basename);
      IlluMultiRead (PETSC_COMM_WORLD, user->da, X, basename, &usermetacount,
		     &usermetanames, &usermetadata);
      /* TODO: free these */
      for (i=0; i<usermetacount; i++)
	PetscPrintf (PETSC_COMM_WORLD, "Parameter \"%s\" = \"%s\"\n",
		     usermetanames [i], usermetadata [i]);
    }
  else
    {
      DPRINTF ("Getting local array\n",0);
      ierr = VecGetArray(localX,&x); CHKERRQ (ierr);

      /* Get local grid boundaries (for 2-dimensional DA):
	 xs, ys   - starting grid indices (no ghost points)
	 xm, ym   - widths of local grid (no ghost points)
	 gxs, gys - starting grid indices (including ghost points)
	 gxm, gym - widths of local grid (including ghost points) */
      DPRINTF ("Getting corners and ghost corners\n",0);
      ierr = DAGetCorners (user->da, &xs,&ys,&zs, &xm,&ym,&zm); CHKERRQ (ierr);
      ierr = DAGetGhostCorners (user->da, &gxs,&gys,&gzs, &gxm,&gym,&gzm);
      CHKERRQ (ierr);

      /* Set up 2-D so it works */
      if (!user->threedee) { zs = gzs = 0; zm = 1; mz = 10; }

      /* Compute initial condition over the locally owned part of the grid
	 Initial condition is motionless fluid and equilibrium temperature */
      DPRINTF ("Looping to set initial condition\n",0);
      for (k=zs; k<zs+zm; k++)
	for (j=ys; j<ys+ym; j++)
	  for (i=xs; i<xs+xm; i++)
	    {
	      row = i - gxs + (j - gys)*gxm + (k-gzs)*gxm*gym;
	      /* x[C(row)] = (i>=mx*0.43921) ? 1.0 : 0.0; */
	      x[C(row)]     = ((i<mx/3 || i>2*mx/3) && (j<my/3 || j>2*my/3) &&
			       (k<mz/3 || k>2*mz/3)) ? 1.0 : 0.0;
	      /* x[C(row)] = (i<mx/2 && j<my/2 && k<mz/2) ? 1.0 : 0.0; */
	      /* x[V(row)]     = 0.0;
		 x[Omega(row)] = 0.0;
		 x[Temp(row)]  = (mparam>0)*(PetscScalar)(i)/(PetscScalar)(mx-1); */
	    }

      /* Restore vector */
      DPRINTF ("Restoring array to local vector\n",0);
      ierr = VecRestoreArray (localX,&x); CHKERRQ (ierr);

      /* Insert values into global vector */
      DPRINTF ("Inserting local vector values into global vector\n",0);
      ierr = DALocalToGlobal (user->da,localX,INSERT_VALUES,X); CHKERRQ (ierr);
    }

  return 0;
}


#undef __FUNCT__
#define __FUNCT__ "InitializeProblem"

/*++++++++++++++++++++++++++++++++++++++
  This takes the gory details of initialization out of the way (importing
  parameters into the user context, etc.).

  int InitializeProblem Returns zero or error.

  AppCtx *user The user context to fill.

  Vec *xvec Vector into which to put the initial condition.
  ++++++++++++++++++++++++++++++++++++++*/

int InitializeProblem (AppCtx *user, Vec *xvec)
{
  int        Nx,Ny,Nz;  /* number of processors in x-, y- and z- directions */
  int        xs,xm,ys,ym,zs,zm, Nlocal,ierr;
  Vec        xv;
  PetscTruth twodee;

  /* Initialize problem parameters */
  DPRINTF ("Initializing user->parameters\n",0);
  user->mx = 20;
  user->my = 20; 
  user->mz = 20; 
  ierr = PetscOptionsGetInt(PETSC_NULL, "-mx", &user->mx, PETSC_NULL);
  CHKERRQ (ierr);
  ierr = PetscOptionsGetInt(PETSC_NULL, "-my", &user->my, PETSC_NULL);
  CHKERRQ (ierr);
  ierr = PetscOptionsGetInt(PETSC_NULL, "-mz", &user->mz, PETSC_NULL);
  CHKERRQ (ierr);
  /* No. of components in the unknown vector and auxiliary vector */
  user->mc = 1;
  /* Problem parameters (kappa, epsilon, gamma, and mparam) */
  user->kappa = 1.0;
  ierr = PetscOptionsGetScalar(PETSC_NULL, "-kappa", &user->kappa, PETSC_NULL);
  CHKERRQ (ierr);
  user->epsilon = 1.0/(user->mx-1);
  ierr = PetscOptionsGetScalar (PETSC_NULL, "-epsilon", &user->epsilon,
				PETSC_NULL); CHKERRQ (ierr);
  user->gamma = 1.0;
  ierr = PetscOptionsGetScalar(PETSC_NULL, "-gamma", &user->gamma, PETSC_NULL);
  CHKERRQ (ierr);
  user->mparam = 0.0;
  ierr = PetscOptionsGetScalar (PETSC_NULL, "-mparam", &user->mparam,
				PETSC_NULL); CHKERRQ (ierr);
  ierr = PetscOptionsHasName (PETSC_NULL, "-twodee", &twodee);
  user->threedee = !twodee;
  CHKERRQ (ierr);
  ierr = PetscOptionsHasName (PETSC_NULL, "-printv", &user->print_vecs);
  CHKERRQ (ierr);
  ierr = PetscOptionsHasName (PETSC_NULL, "-printg", &user->print_grid);
  CHKERRQ (ierr);
  ierr = PetscOptionsHasName (PETSC_NULL, "-no_contours", &user->no_contours);
  CHKERRQ (ierr);
  ierr = PetscOptionsHasName (PETSC_NULL, "-random", &user->random);
  CHKERRQ (ierr);
  ierr = PetscOptionsHasName (PETSC_NULL, "-save_data", &user->save_data);
  CHKERRQ (ierr);
  user->load_data = -1;
  ierr = PetscOptionsGetInt (PETSC_NULL, "-load_data", &user->load_data,
			     PETSC_NULL); CHKERRQ (ierr);

  /* Create distributed array (DA) to manage parallel grid and vectors
     for principal unknowns (x) and governing residuals (f)
     Note the stencil width of 2 for this 4th-order equation. */
  DPRINTF ("Creating distributed array\n",0);
  Nx = PETSC_DECIDE;
  Ny = PETSC_DECIDE;
  Nz = PETSC_DECIDE;
  ierr = PetscOptionsGetInt (PETSC_NULL, "-Nx", &Nx, PETSC_NULL);CHKERRQ(ierr);
  ierr = PetscOptionsGetInt (PETSC_NULL, "-Ny", &Ny, PETSC_NULL);CHKERRQ(ierr);
  ierr = PetscOptionsGetInt (PETSC_NULL, "-Nz", &Nz, PETSC_NULL);CHKERRQ(ierr);
  if (user->threedee)
    {
      DPRINTF ("Three Dee!\n",0);
      user->period = DA_XYZPERIODIC;
      ierr = DACreate3d (PETSC_COMM_WORLD, user->period, DA_STENCIL_BOX,
			 user->mx, user->my, user->mz, Nx, Ny, Nz, user->mc, 2,
			 PETSC_NULL, PETSC_NULL, PETSC_NULL, &user->da);
      CHKERRQ (ierr);
    }
  else
    {
      user->period = DA_XYPERIODIC;
      ierr = DACreate2d (PETSC_COMM_WORLD, user->period, DA_STENCIL_BOX,
			 user->mx, user->my, Nx, Ny, user->mc, 2,
			 PETSC_NULL, PETSC_NULL, &user->da); CHKERRQ (ierr);
      user->mz = Nz = 1;
    }

  /* Extract global vector from DA */
  DPRINTF ("Extracting global vector from DA...\n",0);
  ierr = DACreateGlobalVector(user->da,&xv); CHKERRQ (ierr);

  /* Label PDE components */
  DPRINTF ("Labeling PDE components\n",0);
  ierr = PetscMalloc (user->mc * sizeof(char*), &(user->label));CHKERRQ (ierr);
  user->label[0] = "Concentration (C)";
  /* user->label[1] = "Velocity (V)";
     user->label[2] = "Omega";
     user->label[3] = "Temperature"; */
  ierr = DASetFieldName (user->da, 0, user->label[0]); CHKERRQ(ierr);

  user->x_old = 0;

  /* Get local vector */
  DPRINTF ("Getting local vector\n",0);
  ierr = DACreateLocalVector (user->da, &user->localX); CHKERRQ (ierr);

  /* Print grid info */
  if (user->print_grid)
    {
      DPRINTF ("Printing grid information\n",0);
      ierr = DAView(user->da,PETSC_VIEWER_STDOUT_SELF); CHKERRQ (ierr);
      ierr = DAGetCorners(user->da,&xs,&ys,&zs,&xm,&ym,&zm); CHKERRQ (ierr);
      if (!user->threedee) {
	zs = 0; zm = 1; }
      ierr = PetscPrintf(PETSC_COMM_WORLD,
			 "global grid: %d X %d X %d with %d components per"
			 " node ==> global vector dimension %d\n",
			 user->mx, user->my, user->mz, user->mc,
			 user->mc*user->mx*user->my*user->mz);
      CHKERRQ (ierr);
      fflush(stdout);
      ierr = VecGetLocalSize (xv, &Nlocal); CHKERRQ (ierr);
      ierr = PetscSynchronizedPrintf
	(PETSC_COMM_WORLD,"[%d] local grid %d X %d X %d with %d components"
	 " per node ==> local vector dimension %d\n",
	 user->rank, xm, ym, zm, user->mc, Nlocal); CHKERRQ (ierr);
      ierr = PetscSynchronizedFlush (PETSC_COMM_WORLD); CHKERRQ (ierr);
  }  

  /* Compute initial condition */
  DPRINTF ("Forming inital condition\n",0);
  ierr = FormInitialCondition (user, xv); CHKERRQ (ierr);

  *xvec = xv;
  return 0;
}