Actual source code: ex10.c
petsc-3.4.2 2013-07-02
1: static char help[] = "Solves C_t = -D*C_xx + F(C) + R(C) + D(C) from Brian Wirth's SciDAC project.\n";
3: /*
4: C_t = -D*C_xx + F(C) + R(C) + D(C) from Brian Wirth's SciDAC project.
6: D*C_xx - diffusion of He[1-5] and V[1] and I[1]
7: F(C) - forcing function; He being created.
8: R(C) - reaction terms (clusters combining)
9: D(C) - dissociation terms (cluster breaking up)
11: Sample Options:
12: -ts_monitor_draw_solution -- plot the solution for each concentration as a function of x each in a separate 1d graph
13: -draw_fields_by_name 1-He-2-V,1-He -- only plot the solution for these two concentrations
14: -mymonitor -- plot the concentrations of He and V as a function of x and cluster size (2d contour plot)
15: -da_refine <n=1,2,...> -- run on a finer grid
16: -ts_max_steps maxsteps -- maximum number of time-steps to take
17: -ts_final_time time -- maximum time to compute to
19: Rules for maximum number of He allowed for V in cluster
22: */
24: #include <petscdmda.h>
25: #include <petscts.h>
27: /* Hard wire the number of cluster sizes for He, V, and I */
28: #define N 15
30: /*
31: Define all the concentrations (there is one of these unions at each grid point)
33: He[He] represents the clusters of pure Helium of size He
34: V[V] the Vacencies of size V,
35: I[I] represents the clusters of Interstials of size I, and
36: HeV[He][V] the mixed Helium-Vacancy clusters of size He and V
38: The variables He, V, I are always used to index into the concentrations of He, V, and I respectively
39: Note that unlike in traditional C code the indices for He[], V[] and I[] run from 1 to N, NOT 0 to N-1
40: (the use of the union below "tricks" the C compiler to allow the indices to start at 1.)
42: */
43: typedef struct {
44: PetscScalar He[N];
45: PetscScalar V[N];
46: union {
47: PetscScalar I[N];
48: PetscScalar HeV[N+1][N]; /* actual size is N by N, the N+1 is there only to "trick" the compiler to have indices start at 1.*/
49: };
50: } Concentrations;
52: /*
53: Holds problem specific options and data
54: */
55: typedef struct {
56: PetscBool noreactions; /* run without the reaction terms */
57: PetscBool nodissociations; /* run without the dissociation terms */
58: PetscScalar HeDiffusion[6];
59: PetscScalar VDiffusion[2];
60: PetscScalar IDiffusion[2];
61: PetscScalar forcingScale;
62: PetscScalar reactionScale;
63: PetscScalar dissociationScale;
64: } AppCtx;
66: extern PetscErrorCode IFunction(TS,PetscReal,Vec,Vec,Vec,void*);
67: extern PetscErrorCode InitialConditions(DM,Vec);
68: extern PetscErrorCode MyMonitorSetUp(TS);
69: extern PetscErrorCode GetDfill(PetscInt*,void*);
73: int main(int argc,char **argv)
74: {
75: TS ts; /* nonlinear solver */
76: Vec C; /* solution */
78: DM da; /* manages the grid data */
79: AppCtx ctx; /* holds problem specific paramters */
80: PetscInt He,dof = 3*N+N*N,*ofill,*dfill;
81:
83: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
84: Initialize program
85: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
86: PetscInitialize(&argc,&argv,(char*)0,help);
89: ctx.noreactions = PETSC_FALSE;
90: ctx.nodissociations = PETSC_FALSE;
92: PetscOptionsHasName(NULL,"-noreactions",&ctx.noreactions);
93: PetscOptionsHasName(NULL,"-nodissociations",&ctx.nodissociations);
95: ctx.HeDiffusion[1] = 1000*2.95e-4; /* From Tibo's notes times 1,000 */
96: ctx.HeDiffusion[2] = 1000*3.24e-4;
97: ctx.HeDiffusion[3] = 1000*2.26e-4;
98: ctx.HeDiffusion[4] = 1000*1.68e-4;
99: ctx.HeDiffusion[5] = 1000*5.20e-5;
100: ctx.VDiffusion[1] = 1000*2.71e-3;
101: ctx.IDiffusion[1] = 1000*2.13e-4;
102: ctx.forcingScale = 100.; /* made up numbers */
103: ctx.reactionScale = .001;
104: ctx.dissociationScale = .0001;
105: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106: Create distributed array (DMDA) to manage parallel grid and vectors
107: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
108: DMDACreate1d(PETSC_COMM_WORLD, DMDA_BOUNDARY_MIRROR,-8,dof,1,NULL,&da);
110: /* The only spatial coupling in the Jacobian (diffusion) is for the first 5 He, the first V, and the first I.
111: The ofill (thought of as a dof by dof 2d (row-oriented) array represents the nonzero coupling between degrees
112: of freedom at one point with degrees of freedom on the adjacent point to the left or right. A 1 at i,j in the
113: ofill array indicates that the degree of freedom i at a point is coupled to degree of freedom j at the
114: adjacent point. In this case ofill has only a few diagonal entries since the only spatial coupling is regular diffusion. */
115: PetscMalloc(dof*dof*sizeof(PetscInt),&ofill);
116: PetscMalloc(dof*dof*sizeof(PetscInt),&dfill);
117: PetscMemzero(ofill,dof*dof*sizeof(PetscInt));
118: PetscMemzero(dfill,dof*dof*sizeof(PetscInt));
120: for (He=0; He<PetscMin(N,5); He++) ofill[He*dof + He] = 1;
121: ofill[N*dof + N] = ofill[2*N*dof + 2*N] = 1;
123: DMDASetBlockFills(da,NULL,ofill);
124: PetscFree(ofill);
125: GetDfill(dfill,&ctx);
126: PetscFree(dfill);
128:
130: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
131: Extract global vector from DMDA to hold solution
132: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
133: DMCreateGlobalVector(da,&C);
135: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
136: Create timestepping solver context
137: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
138: TSCreate(PETSC_COMM_WORLD,&ts);
139: TSSetType(ts,TSARKIMEX);
140: TSSetDM(ts,da);
141: TSSetProblemType(ts,TS_NONLINEAR);
142: TSSetIFunction(ts,NULL,IFunction,&ctx);
143: TSSetSolution(ts,C);
145: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
146: Set solver options
147: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
148: TSSetInitialTimeStep(ts,0.0,.001);
149: TSSetDuration(ts,100,50.0);
150: TSSetFromOptions(ts);
151: MyMonitorSetUp(ts);
153: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
154: Set initial conditions
155: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
156: InitialConditions(da,C);
158: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
159: Solve the ODE system
160: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
161: TSSolve(ts,C);
163: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
164: Free work space.
165: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
166: VecDestroy(&C);
167: TSDestroy(&ts);
168: DMDestroy(&da);
169: PetscFinalize();
170: return(0);
171: }
173: /* ------------------------------------------------------------------- */
176: PetscErrorCode InitialConditions(DM da,Vec C)
177: {
179: PetscInt i,I,He,V,xs,xm,Mx,cnt = 0;
180: Concentrations *c;
181: PetscReal hx,x;
182: char string[16];
185: DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
186: hx = 1.0/(PetscReal)(Mx-1);
188: /* Name each of the concentrations */
189: for (He=1; He<N+1; He++) {
190: PetscSNPrintf(string,16,"%d-He",He);
191: DMDASetFieldName(da,cnt++,string);
192: }
193: for (V=1; V<N+1; V++) {
194: PetscSNPrintf(string,16,"%d-V",V);
195: DMDASetFieldName(da,cnt++,string);
196: }
197: for (I=1; I<N+1; I++) {
198: PetscSNPrintf(string,16,"%d-I",I);
199: DMDASetFieldName(da,cnt++,string);
200: }
201: for (He=1; He<N+1; He++) {
202: for (V=1; V<N+1; V++) {
203: PetscSNPrintf(string,16,"%d-He-%d-V",He,V);
204: DMDASetFieldName(da,cnt++,string);
205: }
206: }
208: /*
209: Get pointer to vector data
210: */
211: DMDAVecGetArray(da,C,&c);
212: /* Shift the c pointer to allow accessing with index of 1, instead of 0 */
213: c = (Concentrations*)(((PetscScalar*)c)-1);
215: /*
216: Get local grid boundaries
217: */
218: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
220: /*
221: Compute function over the locally owned part of the grid
222: */
223: for (i=xs; i<xs+xm; i++) {
224: x = i*hx;
225: for (He=1; He<N+1; He++) c[i].He[He] = 0.0;
226: for (V=1; V<N+1; V++) c[i].V[V] = 1.0;
227: for (I=1; I<N+1; I++) c[i].I[I] = 1.0;
228: for (He=1; He<N+1; He++) {
229: for (V=1; V<N+1; V++) c[i].HeV[He][V] = 0.0;
230: }
231: }
233: /*
234: Restore vectors
235: */
236: c = (Concentrations*)(((PetscScalar*)c)+1);
237: DMDAVecRestoreArray(da,C,&c);
238: return(0);
239: }
241: /* ------------------------------------------------------------------- */
244: /*
245: IFunction - Evaluates nonlinear function that defines the ODE
247: Input Parameters:
248: . ts - the TS context
249: . U - input vector
250: . ptr - optional user-defined context
252: Output Parameter:
253: . F - function values
254: */
255: PetscErrorCode IFunction(TS ts,PetscReal ftime,Vec C,Vec Cdot,Vec F,void *ptr)
256: {
257: AppCtx *ctx = (AppCtx*) ptr;
258: DM da;
260: PetscInt xi,Mx,xs,xm,He,he,V,v,I,i;
261: PetscReal hx,sx,x;
262: Concentrations *c,*f;
263: Vec localC;
266: TSGetDM(ts,&da);
267: DMGetLocalVector(da,&localC);
268: DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
269: hx = 8.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
271: /*
272: F = Cdot + all the diffusion and reaction terms added below
273: */
274: VecCopy(Cdot,F);
276: /*
277: Scatter ghost points to local vector,using the 2-step process
278: DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
279: By placing code between these two statements, computations can be
280: done while messages are in transition.
281: */
282: DMGlobalToLocalBegin(da,C,INSERT_VALUES,localC);
283: DMGlobalToLocalEnd(da,C,INSERT_VALUES,localC);
285: /*
286: Get pointers to vector data
287: */
288: DMDAVecGetArray(da,localC,&c);
289: /* Shift the c pointer to allow accessing with index of 1, instead of 0 */
290: c = (Concentrations*)(((PetscScalar*)c)-1);
291: DMDAVecGetArray(da,F,&f);
292: f = (Concentrations*)(((PetscScalar*)f)-1);
294: /*
295: Get local grid boundaries
296: */
297: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
299: /*
300: Loop over grid points computing ODE terms for each grid point
301: */
302: for (xi=xs; xi<xs+xm; xi++) {
303: x = xi*hx;
305: /* -------------------------------------------------------------
306: ---- Compute diffusion over the locally owned part of the grid
307: */
308: /* He clusters larger than 5 do not diffuse -- are immobile */
309: for (He=1; He<PetscMin(N+1,6); He++) {
310: f[xi].He[He] -= ctx->HeDiffusion[He]*(-2.0*c[xi].He[He] + c[xi-1].He[He] + c[xi+1].He[He])*sx;
311: }
313: /* V and I clusters ONLY of size 1 diffuse */
314: f[xi].V[1] -= ctx->VDiffusion[1]*(-2.0*c[xi].V[1] + c[xi-1].V[1] + c[xi+1].V[1])*sx;
315: f[xi].I[1] -= ctx->IDiffusion[1]*(-2.0*c[xi].I[1] + c[xi-1].I[1] + c[xi+1].I[1])*sx;
317: /* Mixed He - V clusters are immobile */
319: /* ----------------------------------------------------------------
320: ---- Compute forcing that produces He of cluster size 1
321: Crude cubic approximation of graph from Tibo's notes
322: */
323: f[xi].He[1] -= ctx->forcingScale*PetscMax(0.0,0.0006*x*x*x - 0.0087*x*x + 0.0300*x);
324: /* Are V or I produced? */
326: if (ctx->noreactions) continue;
327: /* ----------------------------------------------------------------
328: ---- Compute reaction terms that can create a cluster of given size
329: */
330: /* He[He] + He[he] -> He[He+he] */
331: for (He=2; He<N+1; He++) {
332: /* compute all pairs of clusters of smaller size that can combine to create a cluster of size He,
333: remove the upper half since they are symmetric to the lower half of the pairs. For example
334: when He = 5 (cluster size 5) the pairs are
335: 1 4
336: 2 2
337: 3 2 these last two are not needed in the sum since they repeat from above
338: 4 1 this is why he < (He/2) + 1 */
339: for (he=1; he<(He/2)+1; he++) {
340: f[xi].He[He] -= ctx->reactionScale*c[xi].He[he]*c[xi].He[He-he];
342: /* remove the two clusters that merged to form the larger cluster */
343: f[xi].He[he] += ctx->reactionScale*c[xi].He[he]*c[xi].He[He-he];
344: f[xi].He[He-he] += ctx->reactionScale*c[xi].He[he]*c[xi].He[He-he];
345: }
346: }
347: /* V[V] + V[v] -> V[V+v] */
348: for (V=2; V<N+1; V++) {
349: for (v=1; v<(V/2)+1; v++) {
350: f[xi].V[V] -= ctx->reactionScale*c[xi].V[v]*c[xi].V[V-v];
351: /* remove the clusters that merged to form the larger cluster */
352: f[xi].V[v] += ctx->reactionScale*c[xi].V[v]*c[xi].V[V-v];
353: f[xi].V[V-v] += ctx->reactionScale*c[xi].V[v]*c[xi].V[V-v];
354: }
355: }
356: /* I[I] + I[i] -> I[I+i] */
357: for (I=2; I<N+1; I++) {
358: for (i=1; i<(I/2)+1; i++) {
359: f[xi].I[I] -= ctx->reactionScale*c[xi].I[i]*c[xi].I[I-i];
360: /* remove the clusters that merged to form the larger cluster */
361: f[xi].I[i] += ctx->reactionScale*c[xi].I[i]*c[xi].I[I-i];
362: f[xi].I[I-i] += ctx->reactionScale*c[xi].I[i]*c[xi].I[I-i];
363: }
364: }
365: /* He[1] + V[1] -> He[1]-V[1] */
366: f[xi].HeV[1][1] -= 1000*ctx->reactionScale*c[xi].He[1]*c[xi].V[1];
367: /* remove the He and V that merged to form the He-V cluster */
368: f[xi].He[1] += 1000*ctx->reactionScale*c[xi].He[1]*c[xi].V[1];
369: f[xi].V[1] += 1000*ctx->reactionScale*c[xi].He[1]*c[xi].V[1];
370: /* He[He]-V[V] + He[he] -> He[He+he]-V[V] */
371: for (He=1; He<N; He++) {
372: for (V=1; V<N+1; V++) {
373: for (he=1; he<N-He+1; he++) {
374: f[xi].HeV[He+he][V] -= ctx->reactionScale*c[xi].HeV[He][V]*c[xi].He[he];
375: /* remove the two clusters that merged to form the larger cluster */
376: f[xi].He[he] += ctx->reactionScale*c[xi].HeV[He][V]*c[xi].He[he];
377: f[xi].HeV[He][V] += ctx->reactionScale*c[xi].HeV[He][V]*c[xi].He[he];
378: }
379: }
380: }
381: /* He[He]-V[V] + V[v] -> He[He][V+v] */
382: for (He=1; He<N+1; He++) {
383: for (V=1; V<N; V++) {
384: for (v=1; v<N-V+1; v++) {
385: f[xi].HeV[He][V+v] -= ctx->reactionScale*c[xi].HeV[He][V]*c[xi].V[v];
386: /* remove the two clusters that merged to form the larger cluster */
387: f[xi].V[v] += ctx->reactionScale*c[xi].HeV[He][V]*c[xi].V[v];
388: f[xi].HeV[He][V] += ctx->reactionScale*c[xi].HeV[He][V]*c[xi].V[v];
389: }
390: }
391: }
392: /* He[He]-V[V] + He[he]-V[v] -> He[He+he][V+v] */
393: /* Currently the reaction rates for this are zero */
394: for (He=1; He<N; He++) {
395: for (V=1; V<N; V++) {
396: for (he=1; he<N-He+1; he++) {
397: for (v=1; v<N-V+1; v++) {
398: f[xi].HeV[He+he][V+v] -= 0.0*c[xi].HeV[He][V]*c[xi].HeV[he][v];
399: /* remove the two clusters that merged to form the larger cluster */
400: f[xi].HeV[he][V] += 0.0*c[xi].HeV[He][V]*c[xi].HeV[he][v];
401: f[xi].HeV[He][V] += 0.0*c[xi].HeV[He][V]*c[xi].HeV[he][v];
402: }
403: }
404: }
405: }
406: /* V[V] + I[I] -> V[V-I] if V > I else I[I-V] */
407: /* What should the correct reaction rate should be? */
408: for (V=1; V<N+1; V++) {
409: for (I=1; I<V; I++) {
410: f[xi].V[V-I] -= ctx->reactionScale*c[xi].V[V]*c[xi].I[I];
411: f[xi].V[V] += ctx->reactionScale*c[xi].V[V]*c[xi].I[I];
412: f[xi].I[I] += ctx->reactionScale*c[xi].V[V]*c[xi].I[I];
413: }
414: for (I=V+1; I<N+1; I++) {
415: f[xi].I[I-V] -= ctx->reactionScale*c[xi].V[V]*c[xi].I[I];
416: f[xi].V[V] += ctx->reactionScale*c[xi].V[V]*c[xi].I[I];
417: f[xi].I[I] += ctx->reactionScale*c[xi].V[V]*c[xi].I[I];
418: }
419: }
423: if (ctx->nodissociations) continue;
424: /* -------------------------------------------------------------------------
425: ---- Compute dissociation terms that removes an item from a cluster
426: I assume dissociation means losing only a single item from a cluster
427: I cannot tell from the notes if clusters can break up into any sub-size.
428: */
429: /* He[He] -> He[He-1] + He[1] */
430: for (He=2; He<N+1; He++) {
431: f[xi].He[He-1] -= ctx->dissociationScale*c[xi].He[He];
432: f[xi].He[1] -= ctx->dissociationScale*c[xi].He[He];
433: f[xi].He[He] += ctx->dissociationScale*c[xi].He[He];
434: }
435: /* V[V] -> V[V-1] + V[1] */
436: for (V=2; V<N+1; V++) {
437: f[xi].V[V-1] -= ctx->dissociationScale*c[xi].V[V];
438: f[xi].V[1] -= ctx->dissociationScale*c[xi].V[V];
439: f[xi].V[V] += ctx->dissociationScale*c[xi].V[V];
440: }
441: /* I[I] -> I[I-1] + I[1] */
442: for (I=2; I<N+1; I++) {
443: f[xi].I[I-1] -= ctx->dissociationScale*c[xi].I[I];
444: f[xi].I[1] -= ctx->dissociationScale*c[xi].I[I];
445: f[xi].I[I] += ctx->dissociationScale*c[xi].I[I];
446: }
447: /* He[1]-V[1] -> He[1] + V[1] */
448: f[xi].He[1] -= 1000*ctx->reactionScale*c[xi].HeV[1][1];
449: f[xi].V[1] -= 1000*ctx->reactionScale*c[xi].HeV[1][1];
450: f[xi].HeV[1][1] += 1000*ctx->reactionScale*c[xi].HeV[1][1];
451: /* He[He]-V[1] -> He[He] + V[1] */
452: for (He=2; He<N+1; He++) {
453: f[xi].He[He] -= 1000*ctx->reactionScale*c[xi].HeV[He][1];
454: f[xi].V[1] -= 1000*ctx->reactionScale*c[xi].HeV[He][1];
455: f[xi].HeV[He][1] += 1000*ctx->reactionScale*c[xi].HeV[He][1];
456: }
457: /* He[1]-V[V] -> He[1] + V[V] */
458: for (V=2; V<N+1; V++) {
459: f[xi].He[1] -= 1000*ctx->reactionScale*c[xi].HeV[1][V];
460: f[xi].V[V] -= 1000*ctx->reactionScale*c[xi].HeV[1][V];
461: f[xi].HeV[1][V] += 1000*ctx->reactionScale*c[xi].HeV[1][V];
462: }
463: /* He[He]-V[V] -> He[He-1]-V[V] + He[1] */
464: for (He=2; He<N+1; He++) {
465: for (V=2; V<N+1; V++) {
466: f[xi].He[1] -= 1000*ctx->reactionScale*c[xi].HeV[He][V];
467: f[xi].HeV[He-1][V] -= 1000*ctx->reactionScale*c[xi].HeV[He][V];
468: f[xi].HeV[He][V] += 1000*ctx->reactionScale*c[xi].HeV[He][V];
469: }
470: }
471: /* He[He]-V[V] -> He[He]-V[V-1] + V[1] */
472: for (He=2; He<N+1; He++) {
473: for (V=2; V<N+1; V++) {
474: f[xi].V[1] -= 1000*ctx->reactionScale*c[xi].HeV[He][V];
475: f[xi].HeV[He][V-1] -= 1000*ctx->reactionScale*c[xi].HeV[He][V];
476: f[xi].HeV[He][V] += 1000*ctx->reactionScale*c[xi].HeV[He][V];
477: }
478: }
479: /* He[He]-V[V] -> He[He]-V[V+1] + I[1] */
480: /* Again, what is the reasonable dissociation rate? */
481: for (He=1; He<N+1; He++) {
482: for (V=1; V<N; V++) {
483: f[xi].HeV[He][V+1] -= 1000*ctx->reactionScale*c[xi].HeV[He][V];
484: f[xi].I[1] -= 1000*ctx->reactionScale*c[xi].HeV[He][V];
485: f[xi].HeV[He][V] += 1000*ctx->reactionScale*c[xi].HeV[He][V];
486: }
487: }
489: }
491: /*
492: Restore vectors
493: */
494: c = (Concentrations*)(((PetscScalar*)c)+1);
495: DMDAVecRestoreArray(da,localC,&c);
496: f = (Concentrations*)(((PetscScalar*)f)+1);
497: DMDAVecRestoreArray(da,F,&f);
498: DMRestoreLocalVector(da,&localC);
499: return(0);
500: }
505: PetscErrorCode GetDfill(PetscInt *dfill, void *ptr)
506: {
507: AppCtx *ctx = (AppCtx*) ptr;
508: PetscInt He,he,V,v,I,i,dof = 3*N+N*N,reactants[3],row,col1,col2,j;
510: if (!ctx->noreactions) {
511:
512: for (He=2; He<N+1; He++) {
513: /* compute all pairs of clusters of smaller size that can combine to create a cluster of size He,
514: remove the upper half since they are symmetric to the lower half of the pairs. For example
515: when He = 5 (cluster size 5) the pairs are
516: 1 4
517: 2 2
518: 3 2 these last two are not needed in the sum since they repeat from above
519: 4 1 this is why he < (He/2) + 1 */
520: for (he=1; he<(He/2)+1; he++) {
521: reactants[0] = he, reactants[1] = He-he, reactants[2] = He;
522: for (j=0; j<3; j++) {
523: row = reactants[j], col1 = reactants[0], col2 = reactants[1];
524: dfill[(row-1)*dof + col1 - 1] = 1, dfill[(row-1)*dof + col2 - 1] = 1;
525: }
526: }
527: }
528: /* V[V] + V[v] -> V[V+v] */
529: for (V=2; V<N+1; V++) {
530: for (v=1; v<(V/2)+1; v++) {
531: reactants[0] = N+v, reactants[1] = N+V-v, reactants[2] = N+V;
532: for (j=0; j<3; j++) {
533: row = reactants[j], col1 = reactants[0], col2 = reactants[1];
534: dfill[(row-1)*dof + col1 - 1] = 1, dfill[(row-1)*dof + col2 - 1] = 1;
535: }
536: }
537: }
539: /* I[I] + I[i] -> I[I+i] */
540: for (I=2; I<N+1; I++) {
541: for (i=1; i<(I/2)+1; i++) {
542: reactants[0] = 2*N+i, reactants[1] = 2*N+I-i, reactants[2] = 2*N+I;
543: for (j=0; j<3; j++) {
544: row = reactants[j], col1 = reactants[0], col2 = reactants[1];
545: dfill[(row-1)*dof + col1 - 1] = 1, dfill[(row-1)*dof + col2 - 1] = 1;
546: }
547: }
548: }
549:
550: /* He[1] + V[1] -> He[1]-V[1] */
551: reactants[0] = 1, reactants[1] = N+1, reactants[2] = 3*N+1;
552: for (j=0; j<3; j++) {
553: row = reactants[j], col1 = reactants[0], col2 = reactants[1];
554: dfill[(row-1)*dof + col1 - 1] = 1, dfill[(row-1)*dof + col2 - 1] = 1;
555: }
557: /* He[He]-V[V] + He[he] -> He[He+he]-V[V] */
558: for (He=1; He<N; He++) {
559: for (V=1; V<N+1; V++) {
560: for (he=1; he<N-He+1; he++) {
561: reactants[0] = 3*N + (He-1)*N + V, reactants[1] = he, reactants[2] = 3*N+(He+he-1)*N+V;
562: for (j=0; j<3; j++) {
563: row = reactants[j], col1 = reactants[0], col2 = reactants[1];
564: dfill[(row-1)*dof + col1 - 1] = 1, dfill[(row-1)*dof + col2 - 1] = 1;
565: }
566: }
567: }
568: }
569: /* He[He]-V[V] + V[v] -> He[He][V+v] */
570: for (He=1; He<N+1; He++) {
571: for (V=1; V<N; V++) {
572: for (v=1; v<N-V+1; v++) {
573: reactants[0] = 3*N+(He-1)*N+V, reactants[1] = N+v, reactants[2] = 3*N+(He-1)*N+V+v;
574: for (j=0; j<3; j++) {
575: row = reactants[j], col1 = reactants[0], col2 = reactants[1];
576: dfill[(row-1)*dof + col1 - 1] = 1, dfill[(row-1)*dof + col2 - 1] = 1;
577: }
578: }
579: }
580: }
582: /* He[He]-V[V] + He[he]-V[v] -> He[He+he][V+v] */
583: /* Currently the reaction rates for this are zero */
584: for (He=1; He<N; He++) {
585: for (V=1; V<N; V++) {
586: for (he=1; he<N-He+1; he++) {
587: for (v=1; v<N-V+1; v++) {
588: reactants[0] = 3*N+(He-1)*N+V, reactants[1] = 3*N+(he-1)*N+V, reactants[2] = 3*N + (He+he-1)*N + V+v;
589: for (j=0; j<3; j++) {
590: row = reactants[j], col1 = reactants[0], col2 = reactants[1];
591: dfill[(row-1)*dof + col1 - 1] = 1, dfill[(row-1)*dof + col2 - 1] = 1;
592: }
593: }
594: }
595: }
596: }
597: /* V[V] + I[I] -> V[V-I] if V > I else I[I-V] */
598: /* What should the correct reaction rate should be? */
599: for (V=1; V<N+1; V++) {
600: for (I=1; I<V; I++) {
601: reactants[0] = N+V, reactants[1] = 2*N+I, reactants[2] = N+V-I;
602: for (j=0; j<3; j++) {
603: row = reactants[j], col1 = reactants[0], col2 = reactants[1];
604: dfill[(row-1)*dof + col1 - 1] = 1, dfill[(row-1)*dof + col2 - 1] = 1;
605: }
606: }
607: for (I=V+1; I<N+1; I++) {
608: reactants[0] = N+V, reactants[1] = 2*N+I, reactants[2] = 2*N+I-V;
609: for (j=0; j<3; j++) {
610: row = reactants[j], col1 = reactants[0], col2 = reactants[1];
611: dfill[(row-1)*dof + col1 - 1] = 1, dfill[(row-1)*dof + col2 - 1] = 1;
612: }
613: }
614: }
615: }
616: /* -------------------------------------------------------------------------
617: ---- Compute dissociation terms that removes an item from a cluster
618: I assume dissociation means losing only a single item from a cluster
619: I cannot tell from the notes if clusters can break up into any sub-size.
620: */
621: if (!ctx->nodissociations) {
622: /* He[He] -> He[He-1] + He[1] */
623: for (He=2; He<N+1; He++) {
624: reactants[0] = He, reactants[1] = He-1, reactants[2] = 1;
626: for (j=0; j<3; j++) {
627: row = reactants[j], col1 = reactants[0];
628: dfill[(row-1)*dof + col1 - 1] = 1;
629: }
630: }
631: /* V[V] -> V[V-1] + V[1] */
632: for (V=2; V<N+1; V++) {
633: reactants[0] = N+V, reactants[1] = N+V-1, reactants[2] = N+1;
635: for (j=0; j<3; j++) {
636: row = reactants[j], col1 = reactants[0];
637: dfill[(row-1)*dof + col1 - 1] = 1;
638: }
639: }
641: /* I[I] -> I[I-1] + I[1] */
642: for (I=2; I<N+1; I++) {
643: reactants[0] = 2*N+I, reactants[1] = 2*N+I-1, reactants[2] = 2*N+1;
645: for (j=0; j<3; j++) {
646: row = reactants[j], col1 = reactants[0];
647: dfill[(row-1)*dof + col1 - 1] = 1;
648: }
649: }
651: /* He[1]-V[1] -> He[1] + V[1] */
652: reactants[0] = 3*N+1, reactants[1] = 1, reactants[2] = N+1;
654: for (j=0; j<3; j++) {
655: row = reactants[j], col1 = reactants[0];
656: dfill[(row-1)*dof + col1 - 1] = 1;
657: }
658:
659: /* He[He]-V[1] -> He[He] + V[1] */
660: for (He=2; He<N+1; He++) {
661: reactants[0] = 3*N+(He-1)*N+1, reactants[1] = He, reactants[2] = N+1;
663: for (j=0; j<3; j++) {
664: row = reactants[j], col1 = reactants[0];
665: dfill[(row-1)*dof + col1 - 1] = 1;
666: }
667: }
669: /* He[1]-V[V] -> He[1] + V[V] */
670: for (V=2; V<N+1; V++) {
671: reactants[0] = 3*N+V, reactants[1] = 1, reactants[2] = N+V;
673: for (j=0; j<3; j++) {
674: row = reactants[j], col1 = reactants[0];
675: dfill[(row-1)*dof + col1 - 1] = 1;
676: }
677: }
679: /* He[He]-V[V] -> He[He-1]-V[V] + He[1] */
680: for (He=2; He<N+1; He++) {
681: for (V=2; V<N+1; V++) {
682: reactants[0] = 3*N+(He-1)*N+V, reactants[1] = 3*N+(He-2)*N+V, reactants[2] = 1;
684: for (j=0; j<3; j++) {
685: row = reactants[j], col1 = reactants[0];
686: dfill[(row-1)*dof + col1 - 1] = 1;
687: }
688: }
689: }
690:
691: /* He[He]-V[V] -> He[He]-V[V-1] + V[1] */
692: for (He=2; He<N+1; He++) {
693: for (V=2; V<N+1; V++) {
694: reactants[0] = 3*N+(He-1)*N+V, reactants[1] = 3*N+(He-1)*N+V-1, reactants[2] = N+1;
696: for (j=0; j<3; j++) {
697: row = reactants[j], col1 = reactants[0];
698: dfill[(row-1)*dof + col1 - 1] = 1;
699: }
700: }
701: }
703: /* He[He]-V[V] -> He[He]-V[V+1] + I[1] */
704: /* Again, what is the reasonable dissociation rate? */
705: for (He=1; He<N+1; He++) {
706: for (V=1; V<N; V++) {
707: reactants[0] = 3*N+(He-1)*N+V, reactants[1] = 3*N+(He-1)*N+V+1, reactants[2] = 2*N+1;
709: for (j=0; j<3; j++) {
710: row = reactants[j], col1 = reactants[0];
711: dfill[(row-1)*dof + col1 - 1] = 1;
712: }
713: }
714: }
715: }
716: }
717: /* ------------------------------------------------------------------- */
719: typedef struct {
720: DM da; /* defines the 2d layout of the He subvector */
721: Vec He;
722: VecScatter scatter;
723: PetscViewer viewer;
724: } MyMonitorCtx;
728: /*
729: Display He and V as a function of space and cluster size for each time step
730: */
731: PetscErrorCode MyMonitorMonitor(TS ts,PetscInt timestep,PetscReal time,Vec solution, void *ictx)
732: {
733: MyMonitorCtx *ctx = (MyMonitorCtx*)ictx;
737: VecScatterBegin(ctx->scatter,solution,ctx->He,INSERT_VALUES,SCATTER_FORWARD);
738: VecScatterEnd(ctx->scatter,solution,ctx->He,INSERT_VALUES,SCATTER_FORWARD);
739: VecView(ctx->He,ctx->viewer);
740: return(0);
741: }
745: /*
746: Frees all data structures associated with the monitor
747: */
748: PetscErrorCode MyMonitorDestroy(void **ictx)
749: {
750: MyMonitorCtx **ctx = (MyMonitorCtx**)ictx;
754: VecScatterDestroy(&(*ctx)->scatter);
755: VecDestroy(&(*ctx)->He);
756: DMDestroy(&(*ctx)->da);
757: PetscViewerDestroy(&(*ctx)->viewer);
758: PetscFree(*ctx);
759: return(0);
760: }
764: /*
765: Sets up a monitor that will display He as a function of space and cluster size for each time step
766: */
767: PetscErrorCode MyMonitorSetUp(TS ts)
768: {
769: DM da;
771: PetscInt xi,xs,xm,*idx,M,xj,cnt = 0,dof = 3*N + N*N;
772: const PetscInt *lx;
773: Vec C;
774: MyMonitorCtx *ctx;
775: PetscBool flg;
776: IS is;
777: char ycoor[32];
778: PetscReal valuebounds[4] = {0, 1.2, 0, 1.2};
781: PetscOptionsHasName(NULL,"-mymonitor",&flg);
782: if (!flg) return(0);
784: TSGetDM(ts,&da);
785: PetscNew(MyMonitorCtx,&ctx);
786: PetscViewerDrawOpen(PetscObjectComm((PetscObject)da),NULL,"",PETSC_DECIDE,PETSC_DECIDE,600,400,&ctx->viewer);
788: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
789: DMDAGetInfo(da,PETSC_IGNORE,&M,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
790: DMDAGetOwnershipRanges(da,&lx,NULL,NULL);
791: DMDACreate2d(PetscObjectComm((PetscObject)da),DMDA_BOUNDARY_NONE,DMDA_BOUNDARY_NONE,DMDA_STENCIL_STAR,M,N,PETSC_DETERMINE,1,2,1,lx,NULL,&ctx->da);
792: DMDASetFieldName(ctx->da,0,"He");
793: DMDASetFieldName(ctx->da,1,"V");
794: DMDASetCoordinateName(ctx->da,0,"X coordinate direction");
795: PetscSNPrintf(ycoor,32,"%D ... Cluster size ... 1",N);
796: DMDASetCoordinateName(ctx->da,1,ycoor);
798: DMCreateGlobalVector(ctx->da,&ctx->He);
799: PetscMalloc(2*N*xm*sizeof(PetscInt),&idx);
800: cnt = 0;
801: for (xj=0; xj<N; xj++) {
802: for (xi=xs; xi<xs+xm; xi++) {
803: idx[cnt++] = dof*xi + xj;
804: idx[cnt++] = dof*xi + xj + N;
805: }
806: }
807: ISCreateGeneral(PetscObjectComm((PetscObject)ts),2*N*xm,idx,PETSC_OWN_POINTER,&is);
808: TSGetSolution(ts,&C);
809: VecScatterCreate(C,is,ctx->He,NULL,&ctx->scatter);
810: ISDestroy(&is);
812: /* sets the bounds on the contour plot values so the colors mean the same thing for different timesteps */
813: PetscViewerDrawSetBounds(ctx->viewer,2,valuebounds);
815: TSMonitorSet(ts,MyMonitorMonitor,ctx,MyMonitorDestroy);
816: return(0);
817: }