Actual source code: ex10.c
petsc-3.7.5 2017-01-01
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: */
21: #include <petscdm.h>
22: #include <petscdmda.h>
23: #include <petscts.h>
25: /* Hard wire the number of cluster sizes for He, V, and I, and He-V */
26: #define NHe 9
27: #define NV 10 /* 50 */
28: #define NI 2
29: #define MHeV 10 /* 50 */ /* maximum V size in He-V */
30: PetscInt NHeV[MHeV+1]; /* maximum He size in an He-V with given V */
31: #define MNHeV 451 /* 6778 */
32: #define DOF (NHe + NV + NI + MNHeV)
34: /*
35: Define all the concentrations (there is one of these structs at each grid point)
37: He[He] represents the clusters of pure Helium of size He
38: V[V] the Vacencies of size V,
39: I[I] represents the clusters of Interstials of size I, and
40: HeV[He][V] the mixed Helium-Vacancy clusters of size He and V
42: The variables He, V, I are always used to index into the concentrations of He, V, and I respectively
43: Note that unlike in traditional C code the indices for He[], V[] and I[] run from 1 to N, NOT 0 to N-1
45: */
46: typedef struct {
47: PetscScalar He[NHe];
48: PetscScalar V[NV];
49: PetscScalar I[NI];
50: PetscScalar HeV[MNHeV];
51: } Concentrations;
55: /*
56: Holds problem specific options and data
57: */
58: typedef struct {
59: PetscScalar HeDiffusion[6];
60: PetscScalar VDiffusion[2];
61: PetscScalar IDiffusion[2];
62: PetscScalar forcingScale;
63: PetscScalar reactionScale;
64: PetscScalar dissociationScale;
65: } AppCtx;
67: extern PetscErrorCode RHSFunction(TS,PetscReal,Vec,Vec,void*);
68: extern PetscErrorCode RHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*);
69: extern PetscErrorCode InitialConditions(DM,Vec);
70: extern PetscErrorCode MyMonitorSetUp(TS);
71: extern PetscErrorCode GetDfill(PetscInt*,void*);
72: extern PetscErrorCode MyLoadData(MPI_Comm,const char*);
76: int main(int argc,char **argv)
77: {
78: TS ts; /* nonlinear solver */
79: Vec C; /* solution */
81: DM da; /* manages the grid data */
82: AppCtx ctx; /* holds problem specific paramters */
83: PetscInt He,*ofill,*dfill;
84: char filename[PETSC_MAX_PATH_LEN];
85: PetscBool flg;
87: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
88: Initialize program
89: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
90: PetscInitialize(&argc,&argv,(char*)0,help);
93: PetscOptionsGetString(NULL,NULL,"-file",filename,PETSC_MAX_PATH_LEN,&flg);
94: if (flg) {
95: MyLoadData(PETSC_COMM_WORLD,filename);
96: }
99: ctx.HeDiffusion[1] = 1000*2.95e-4; /* From Tibo's notes times 1,000 */
100: ctx.HeDiffusion[2] = 1000*3.24e-4;
101: ctx.HeDiffusion[3] = 1000*2.26e-4;
102: ctx.HeDiffusion[4] = 1000*1.68e-4;
103: ctx.HeDiffusion[5] = 1000*5.20e-5;
104: ctx.VDiffusion[1] = 1000*2.71e-3;
105: ctx.IDiffusion[1] = 1000*2.13e-4;
106: ctx.forcingScale = 100.; /* made up numbers */
107: ctx.reactionScale = .001;
108: ctx.dissociationScale = .0001;
110: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
111: Create distributed array (DMDA) to manage parallel grid and vectors
112: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
113: DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_MIRROR,-2,DOF,1,NULL,&da);
115: /* The only spatial coupling in the Jacobian (diffusion) is for the first 5 He, the first V, and the first I.
116: The ofill (thought of as a DOF by DOF 2d (row-oriented) array) represents the nonzero coupling between degrees
117: 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
118: ofill array indicates that the degree of freedom i at a point is coupled to degree of freedom j at the
119: adjacent point. In this case ofill has only a few diagonal entries since the only spatial coupling is regular diffusion. */
120: PetscMalloc1(dof*dof,&ofill);
121: PetscMalloc1(dof*dof,&dfill);
122: PetscMemzero(ofill,dof*dof*sizeof(PetscInt));
123: PetscMemzero(dfill,dof*dof*sizeof(PetscInt));
125: /*
126: dfil (thought of as a DOF by DOF 2d (row-oriented) array) repesents the nonzero coupling between degrees of
127: freedom within a single grid point, i.e. the reaction and dissassociation interactions. */
128: PetscMalloc(DOF*DOF*sizeof(PetscInt),&dfill);
129: PetscMemzero(dfill,DOF*DOF*sizeof(PetscInt));
130: GetDfill(dfill,&ctx);
131: DMDASetBlockFills(da,dfill,ofill);
132: PetscFree(ofill);
133: PetscFree(dfill);
135: /* Extract global vector to hold solution */
136: DMCreateGlobalVector(da,&C);
138: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
139: Create timestepping solver context
140: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
141: TSCreate(PETSC_COMM_WORLD,&ts);
142: TSSetType(ts,TSARKIMEX);
143: TSARKIMEXSetFullyImplicit(ts,PETSC_TRUE);
144: TSSetDM(ts,da);
145: TSSetProblemType(ts,TS_NONLINEAR);
146: TSSetRHSFunction(ts,NULL,RHSFunction,&ctx);
147: TSSetRHSJacobian(ts,NULL,NULL,RHSJacobian,&ctx);
148: TSSetSolution(ts,C);
150: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
151: Set solver options
152: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
153: TSSetInitialTimeStep(ts,0.0,.001);
154: TSSetDuration(ts,100,50.0);
155: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
156: TSSetFromOptions(ts);
157: MyMonitorSetUp(ts);
159: InitialConditions(da,C);
161: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
162: Solve the ODE system
163: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
164: TSSolve(ts,C);
166: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
167: Free work space.
168: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
169: VecDestroy(&C);
170: TSDestroy(&ts);
171: DMDestroy(&da);
172: PetscFinalize();
173: return(0);
174: }
176: /*
177: cHeV is "trick" to allow easy accessing of the values in the HeV portion of the Concentrations.
178: cHeV[i] points to the beginning of each row of HeV[] with V indexing starting a 1.
180: */
183: PetscErrorCode cHeVCreate(PetscReal ***cHeV)
184: {
188: PetscMalloc(MHeV*sizeof(PetscScalar),cHeV);
189: (*cHeV)--;
190: return(0);
191: }
195: PetscErrorCode cHeVInitialize(const PetscScalar *start,PetscReal **cHeV)
196: {
197: PetscInt i;
200: cHeV[1] = ((PetscScalar*) start) - 1 + NHe + NV + NI;
201: for (i=1; i<MHeV; i++) {
202: cHeV[i+1] = cHeV[i] + NHeV[i];
203: }
204: return(0);
205: }
209: PetscErrorCode cHeVDestroy(PetscReal **cHeV)
210: {
214: cHeV++;
215: PetscFree(cHeV);
216: return(0);
217: }
219: /* ------------------------------------------------------------------- */
222: PetscErrorCode InitialConditions(DM da,Vec C)
223: {
225: PetscInt i,I,He,V,xs,xm,Mx,cnt = 0;
226: Concentrations *c;
227: PetscReal hx,x,**cHeV;
228: char string[16];
231: 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);
232: hx = 1.0/(PetscReal)(Mx-1);
234: /* Name each of the concentrations */
235: for (He=1; He<NHe+1; He++) {
236: PetscSNPrintf(string,16,"%d-He",He);
237: DMDASetFieldName(da,cnt++,string);
238: }
239: for (V=1; V<NV+1; V++) {
240: PetscSNPrintf(string,16,"%d-V",V);
241: DMDASetFieldName(da,cnt++,string);
242: }
243: for (I=1; I<NI+1; I++) {
244: PetscSNPrintf(string,16,"%d-I",I);
245: DMDASetFieldName(da,cnt++,string);
246: }
247: for (He=1; He<MHeV+1; He++) {
248: for (V=1; V<NHeV[He]+1; V++) {
249: PetscSNPrintf(string,16,"%d-He-%d-V",He,V);
250: DMDASetFieldName(da,cnt++,string);
251: }
252: }
254: /*
255: Get pointer to vector data
256: */
257: DMDAVecGetArrayRead(da,C,&c);
258: /* Shift the c pointer to allow accessing with index of 1, instead of 0 */
259: c = (Concentrations*)(((PetscScalar*)c)-1);
261: /*
262: Get local grid boundaries
263: */
264: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
266: /*
267: Compute function over the locally owned part of the grid
268: */
269: cHeVCreate(&cHeV);
270: for (i=xs; i<xs+xm; i++) {
271: x = i*hx;
272: for (He=1; He<NHe+1; He++) c[i].He[He] = 0.0;
273: for (V=1; V<NV+1; V++) c[i].V[V] = 1.0;
274: for (I=1; I <NI+1; I++) c[i].I[I] = 1.0;
275: cHeVInitialize(&c[i].He[1],cHeV);
276: for (V=1; V<MHeV+1; V++) {
277: for (He=1; He<NHeV[V]+1; He++) cHeV[V][He] = 0.0;
278: }
279: }
280: cHeVDestroy(cHeV);
282: /*
283: Restore vectors
284: */
285: c = (Concentrations*)(((PetscScalar*)c)+1);
286: DMDAVecRestoreArrayRead(da,C,&c);
287: return(0);
288: }
290: /* ------------------------------------------------------------------- */
293: /*
294: RHSFunction - Evaluates nonlinear function that defines the ODE
296: Input Parameters:
297: . ts - the TS context
298: . U - input vector
299: . ptr - optional user-defined context
301: Output Parameter:
302: . F - function values
303: */
304: PetscErrorCode RHSFunction(TS ts,PetscReal ftime,Vec C,Vec F,void *ptr)
305: {
306: AppCtx *ctx = (AppCtx*) ptr;
307: DM da;
309: PetscInt xi,Mx,xs,xm,He,he,V,v,I,i;
310: PetscReal hx,sx,x,**cHeV,**fHeV;
311: Concentrations *c,*f;
312: Vec localC;
315: TSGetDM(ts,&da);
316: DMGetLocalVector(da,&localC);
317: 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);
318: hx = 8.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
319: cHeVCreate(&cHeV);
320: cHeVCreate(&fHeV);
322: /*
323: Scatter ghost points to local vector,using the 2-step process
324: DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
325: By placing code between these two statements, computations can be
326: done while messages are in transition.
327: */
328: DMGlobalToLocalBegin(da,C,INSERT_VALUES,localC);
329: DMGlobalToLocalEnd(da,C,INSERT_VALUES,localC);
331: VecSet(F,0.0);
333: /*
334: Get pointers to vector data
335: */
336: DMDAVecGetArrayRead(da,localC,&c);
337: /* Shift the c pointer to allow accessing with index of 1, instead of 0 */
338: c = (Concentrations*)(((PetscScalar*)c)-1);
339: DMDAVecGetArray(da,F,&f);
340: f = (Concentrations*)(((PetscScalar*)f)-1);
342: /*
343: Get local grid boundaries
344: */
345: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
347: /*
348: Loop over grid points computing ODE terms for each grid point
349: */
350: for (xi=xs; xi<xs+xm; xi++) {
351: x = xi*hx;
353: /* -------------------------------------------------------------
354: ---- Compute diffusion over the locally owned part of the grid
355: */
356: /* He clusters larger than 5 do not diffuse -- are immobile */
357: for (He=1; He<PetscMin(NHe+1,6); He++) {
358: f[xi].He[He] += ctx->HeDiffusion[He]*(-2.0*c[xi].He[He] + c[xi-1].He[He] + c[xi+1].He[He])*sx;
359: }
361: /* V and I clusters ONLY of size 1 diffuse */
362: f[xi].V[1] += ctx->VDiffusion[1]*(-2.0*c[xi].V[1] + c[xi-1].V[1] + c[xi+1].V[1])*sx;
363: f[xi].I[1] += ctx->IDiffusion[1]*(-2.0*c[xi].I[1] + c[xi-1].I[1] + c[xi+1].I[1])*sx;
365: /* Mixed He - V clusters are immobile */
367: /* ----------------------------------------------------------------
368: ---- Compute forcing that produces He of cluster size 1
369: Crude cubic approximation of graph from Tibo's notes
370: */
371: f[xi].He[1] += ctx->forcingScale*PetscMax(0.0,0.0006*x*x*x - 0.0087*x*x + 0.0300*x);
373: cHeVInitialize(&c[xi].He[1],cHeV);
374: cHeVInitialize(&f[xi].He[1],fHeV);
376: /* -------------------------------------------------------------------------
377: ---- Compute dissociation terms that removes an item from a cluster
378: I assume dissociation means losing only a single item from a cluster
379: I cannot tell from the notes if clusters can break up into any sub-size.
380: */
381: /* He[He] -> He[He-1] + He[1] */
382: for (He=2; He<NHe+1; He++) {
383: f[xi].He[He-1] += ctx->dissociationScale*c[xi].He[He];
384: f[xi].He[1] += ctx->dissociationScale*c[xi].He[He];
385: f[xi].He[He] -= ctx->dissociationScale*c[xi].He[He];
386: }
388: /* V[V] -> V[V-1] + V[1] */
389: for (V=2; V<NV+1; V++) {
390: f[xi].V[V-1] += ctx->dissociationScale*c[xi].V[V];
391: f[xi].V[1] += ctx->dissociationScale*c[xi].V[V];
392: f[xi].V[V] -= ctx->dissociationScale*c[xi].V[V];
393: }
395: /* I[I] -> I[I-1] + I[1] */
396: for (I=2; I<NI+1; I++) {
397: f[xi].I[I-1] += ctx->dissociationScale*c[xi].I[I];
398: f[xi].I[1] += ctx->dissociationScale*c[xi].I[I];
399: f[xi].I[I] -= ctx->dissociationScale*c[xi].I[I];
400: }
402: /* He[He]-V[1] -> He[He] + V[1] */
403: for (He=1; He<NHeV[1]+1; He++) {
404: f[xi].He[He] += 1000*ctx->dissociationScale*cHeV[1][He];
405: f[xi].V[1] += 1000*ctx->dissociationScale*cHeV[1][He];
406: fHeV[1][He] -= 1000*ctx->dissociationScale*cHeV[1][He];
407: }
409: /* He[1]-V[V] -> He[1] + V[V] */
410: for (V=2; V<MHeV+1; V++) {
411: f[xi].He[1] += 1000*ctx->dissociationScale*cHeV[V][1];
412: f[xi].V[V] += 1000*ctx->dissociationScale*cHeV[V][1];
413: fHeV[V][1] -= 1000*ctx->dissociationScale*cHeV[V][1];
414: }
416: /* He[He]-V[V] -> He[He-1]-V[V] + He[1] */
417: for (V=2; V<MHeV+1; V++) {
418: for (He=2; He<NHeV[V]+1; He++) {
419: f[xi].He[1] += 1000*ctx->dissociationScale*cHeV[V][He];
420: fHeV[V][He-1] += 1000*ctx->dissociationScale*cHeV[V][He];
421: fHeV[V][He] -= 1000*ctx->dissociationScale*cHeV[V][He];
422: }
423: }
425: /* He[He]-V[V] -> He[He]-V[V-1] + V[1] */
426: for (V=2; V<MHeV+1; V++) {
427: for (He=2; He<NHeV[V-1]+1; He++) {
428: f[xi].V[1] += 1000*ctx->dissociationScale*cHeV[V][He];
429: fHeV[V-1][He] += 1000*ctx->dissociationScale*cHeV[V][He];
430: fHeV[V][He] -= 1000*ctx->dissociationScale*cHeV[V][He];
431: }
432: }
434: /* He[He]-V[V] -> He[He]-V[V+1] + I[1] */
435: for (V=1; V<MHeV; V++) {
436: for (He=1; He<NHeV[V]+1; He++) {
437: fHeV[V+1][He] += 1000*ctx->dissociationScale*cHeV[V][He];
438: f[xi].I[1] += 1000*ctx->dissociationScale*cHeV[V][He];
439: fHeV[V][He] -= 1000*ctx->dissociationScale*cHeV[V][He];
440: }
441: }
443: /* ----------------------------------------------------------------
444: ---- Compute reaction terms that can create a cluster of given size
445: */
446: /* He[He] + He[he] -> He[He+he] */
447: for (He=2; He<NHe+1; He++) {
448: /* compute all pairs of clusters of smaller size that can combine to create a cluster of size He,
449: remove the upper half since they are symmetric to the lower half of the pairs. For example
450: when He = 5 (cluster size 5) the pairs are
451: 1 4
452: 2 2
453: 3 2 these last two are not needed in the sum since they repeat from above
454: 4 1 this is why he < (He/2) + 1 */
455: for (he=1; he<(He/2)+1; he++) {
456: f[xi].He[He] += ctx->reactionScale*c[xi].He[he]*c[xi].He[He-he];
458: /* remove the two clusters that merged to form the larger cluster */
459: f[xi].He[he] -= ctx->reactionScale*c[xi].He[he]*c[xi].He[He-he];
460: f[xi].He[He-he] -= ctx->reactionScale*c[xi].He[he]*c[xi].He[He-he];
461: }
462: }
464: /* V[V] + V[v] -> V[V+v] */
465: for (V=2; V<NV+1; V++) {
466: for (v=1; v<(V/2)+1; v++) {
467: f[xi].V[V] += ctx->reactionScale*c[xi].V[v]*c[xi].V[V-v];
468: f[xi].V[v] -= ctx->reactionScale*c[xi].V[v]*c[xi].V[V-v];
469: f[xi].V[V-v] -= ctx->reactionScale*c[xi].V[v]*c[xi].V[V-v];
470: }
471: }
473: /* I[I] + I[i] -> I[I+i] */
474: for (I=2; I<NI+1; I++) {
475: for (i=1; i<(I/2)+1; i++) {
476: f[xi].I[I] += ctx->reactionScale*c[xi].I[i]*c[xi].I[I-i];
477: f[xi].I[i] -= ctx->reactionScale*c[xi].I[i]*c[xi].I[I-i];
478: f[xi].I[I-i] -= ctx->reactionScale*c[xi].I[i]*c[xi].I[I-i];
479: }
480: }
482: /* He[1] + V[1] -> He[1]-V[1] */
483: fHeV[1][1] += 1000*ctx->reactionScale*c[xi].He[1]*c[xi].V[1];
484: f[xi].He[1] -= 1000*ctx->reactionScale*c[xi].He[1]*c[xi].V[1];
485: f[xi].V[1] -= 1000*ctx->reactionScale*c[xi].He[1]*c[xi].V[1];
487: /* He[He]-V[V] + He[he] -> He[He+he]-V[V] */
488: for (V=1; V<MHeV+1; V++) {
489: for (He=1; He<NHeV[V]; He++) {
490: for (he=1; he+He<NHeV[V]+1; he++) {
491: fHeV[V][He+he] += ctx->reactionScale*cHeV[V][He]*c[xi].He[he];
492: f[xi].He[he] -= ctx->reactionScale*cHeV[V][He]*c[xi].He[he];
493: fHeV[V][He] -= ctx->reactionScale*cHeV[V][He]*c[xi].He[he];
494: }
495: }
496: }
498: /* He[He]-V[V] + V[1] -> He[He][V+1] */
499: for (V=1; V<MHeV; V++) {
500: for (He=1; He<NHeV[V+1]; He++) {
501: fHeV[V+1][He] += ctx->reactionScale*cHeV[V][He]*c[xi].V[1];
502: /* remove the two clusters that merged to form the larger cluster */
503: f[xi].V[1] -= ctx->reactionScale*cHeV[V][He]*c[xi].V[1];
504: fHeV[V][He] -= ctx->reactionScale*cHeV[V][He]*c[xi].V[1];
505: }
506: }
508: /* He[He]-V[V] + He[he]-V[v] -> He[He+he][V+v] */
509: /* Currently the reaction rates for this are zero */
512: /* V[V] + I[I] -> V[V-I] if V > I else I[I-V] */
513: for (V=1; V<NV+1; V++) {
514: for (I=1; I<PetscMin(V,NI); I++) {
515: f[xi].V[V-I] += ctx->reactionScale*c[xi].V[V]*c[xi].I[I];
516: f[xi].V[V] -= ctx->reactionScale*c[xi].V[V]*c[xi].I[I];
517: f[xi].I[I] -= ctx->reactionScale*c[xi].V[V]*c[xi].I[I];
518: }
519: for (I=V+1; I<NI+1; I++) {
520: f[xi].I[I-V] += ctx->reactionScale*c[xi].V[V]*c[xi].I[I];
521: f[xi].V[V] -= ctx->reactionScale*c[xi].V[V]*c[xi].I[I];
522: f[xi].I[I] -= ctx->reactionScale*c[xi].V[V]*c[xi].I[I];
523: }
524: }
525: }
527: /*
528: Restore vectors
529: */
530: c = (Concentrations*)(((PetscScalar*)c)+1);
531: DMDAVecRestoreArray(da,localC,&c);
532: f = (Concentrations*)(((PetscScalar*)f)+1);
533: DMDAVecRestoreArray(da,F,&f);
534: DMRestoreLocalVector(da,&localC);
535: cHeVDestroy(cHeV);
536: cHeVDestroy(fHeV);
537: return(0);
538: }
542: /*
543: Compute the Jacobian entries based on IFuction() and insert them into the matrix
544: */
545: PetscErrorCode RHSJacobian(TS ts,PetscReal ftime,Vec C,Mat A,Mat J,void *ptr)
546: {
547: AppCtx *ctx = (AppCtx*) ptr;
548: DM da;
549: PetscErrorCode ierr;
550: PetscInt xi,Mx,xs,xm,He,he,V,v,I,i;
551: PetscInt row[3],col[3];
552: PetscReal hx,sx,x,val[6];
553: const Concentrations *c,*f;
554: Vec localC;
555: const PetscReal *rowstart,*colstart;
556: const PetscReal **cHeV,**fHeV;
557: static PetscBool initialized = PETSC_FALSE;
560: cHeVCreate((PetscScalar***)&cHeV);
561: cHeVCreate((PetscScalar***)&fHeV);
562: MatZeroEntries(J);
563: TSGetDM(ts,&da);
564: DMGetLocalVector(da,&localC);
565: 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);
566: hx = 8.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
568: DMGlobalToLocalBegin(da,C,INSERT_VALUES,localC);
569: DMGlobalToLocalEnd(da,C,INSERT_VALUES,localC);
571: /*
572: The f[] is dummy, values are never set into it. It is only used to determine the
573: local row for the entries in the Jacobian
574: */
575: DMDAVecGetArray(da,localC,&c);
576: /* Shift the c pointer to allow accessing with index of 1, instead of 0 */
577: c = (Concentrations*)(((PetscScalar*)c)-1);
578: DMDAVecGetArray(da,C,&f);
579: f = (Concentrations*)(((PetscScalar*)f)-1);
581: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
583: rowstart = &f[xs].He[1] - DOF;
584: colstart = &c[xs-1].He[1];
586: if (!initialized) {
587: /*
588: Loop over grid points computing Jacobian terms for each grid point
589: */
590: for (xi=xs; xi<xs+xm; xi++) {
591: x = xi*hx;
592:
593: cHeVInitialize(&c[xi].He[1],(PetscScalar**)cHeV);
594: cHeVInitialize(&f[xi].He[1],(PetscScalar**)fHeV);
595:
596: /* -------------------------------------------------------------
597: ---- Compute diffusion over the locally owned part of the grid
598: */
599: /* He clusters larger than 5 do not diffuse -- are immobile */
600: for (He=1; He<PetscMin(NHe+1,6); He++) {
601: row[0] = &f[xi].He[He] - rowstart;
602: col[0] = &c[xi-1].He[He] - colstart;
603: col[1] = &c[xi].He[He] - colstart;
604: col[2] = &c[xi+1].He[He] - colstart;
605: val[0] = ctx->HeDiffusion[He]*sx;
606: val[1] = -2.0*ctx->HeDiffusion[He]*sx;
607: val[2] = ctx->HeDiffusion[He]*sx;
608: MatSetValuesLocal(J,1,row,3,col,val,ADD_VALUES);
609: }
611: /* V and I clusters ONLY of size 1 diffuse */
612: row[0] = &f[xi].V[1] - rowstart;
613: col[0] = &c[xi-1].V[1] - colstart;
614: col[1] = &c[xi].V[1] - colstart;
615: col[2] = &c[xi+1].V[1] - colstart;
616: val[0] = ctx->VDiffusion[1]*sx;
617: val[1] = -2.0*ctx->VDiffusion[1]*sx;
618: val[2] = ctx->VDiffusion[1]*sx;
619: MatSetValuesLocal(J,1,row,3,col,val,ADD_VALUES);
620:
621: row[0] = &f[xi].I[1] - rowstart;
622: col[0] = &c[xi-1].I[1] - colstart;
623: col[1] = &c[xi].I[1] - colstart;
624: col[2] = &c[xi+1].I[1] - colstart;
625: val[0] = ctx->IDiffusion[1]*sx;
626: val[1] = -2.0*ctx->IDiffusion[1]*sx;
627: val[2] = ctx->IDiffusion[1]*sx;
628: MatSetValuesLocal(J,1,row,3,col,val,ADD_VALUES);
629:
630: /* Mixed He - V clusters are immobile */
631:
632: /* -------------------------------------------------------------------------
633: ---- Compute dissociation terms that removes an item from a cluster
634: I assume dissociation means losing only a single item from a cluster
635: I cannot tell from the notes if clusters can break up into any sub-size.
636: */
637:
638: /* He[He] -> He[He-1] + He[1] */
639: for (He=2; He<NHe+1; He++) {
640: row[0] = &f[xi].He[He-1] - rowstart;
641: row[1] = &f[xi].He[1] - rowstart;
642: row[2] = &f[xi].He[He] - rowstart;
643: col[0] = &c[xi].He[He] - colstart;
644: val[0] = ctx->dissociationScale;
645: val[1] = ctx->dissociationScale;
646: val[2] = -ctx->dissociationScale;
647: MatSetValuesLocal(J,3,row,1,col,val,ADD_VALUES);
648: }
649:
650: /* V[V] -> V[V-1] + V[1] */
651: for (V=2; V<NV+1; V++) {
652: row[0] = &f[xi].V[V-1] - rowstart;
653: row[1] = &f[xi].V[1] - rowstart;
654: row[2] = &f[xi].V[V] - rowstart;
655: col[0] = &c[xi].V[V] - colstart;
656: val[0] = ctx->dissociationScale;
657: val[1] = ctx->dissociationScale;
658: val[2] = -ctx->dissociationScale;
659: MatSetValuesLocal(J,3,row,1,col,val,ADD_VALUES);
660: }
661:
662: /* I[I] -> I[I-1] + I[1] */
663: for (I=2; I<NI+1; I++) {
664: row[0] = &f[xi].I[I-1] - rowstart;
665: row[1] = &f[xi].I[1] - rowstart;
666: row[2] = &f[xi].I[I] - rowstart;
667: col[0] = &c[xi].I[I] - colstart;
668: val[0] = ctx->dissociationScale;
669: val[1] = ctx->dissociationScale;
670: val[2] = -ctx->dissociationScale;
671: MatSetValuesLocal(J,3,row,1,col,val,ADD_VALUES);
672: }
673:
674: /* He[He]-V[1] -> He[He] + V[1] */
675: for (He=1; He<NHeV[1]+1; He++) {
676: row[0] = &f[xi].He[He] - rowstart;
677: row[1] = &f[xi].V[1] - rowstart;
678: row[2] = &fHeV[1][He] - rowstart;
679: col[0] = &cHeV[1][He] - colstart;
680: val[0] = 1000*ctx->dissociationScale;
681: val[1] = 1000*ctx->dissociationScale;
682: val[2] = -1000*ctx->dissociationScale;
683: MatSetValuesLocal(J,3,row,1,col,val,ADD_VALUES);
684: }
685:
686: /* He[1]-V[V] -> He[1] + V[V] */
687: for (V=2; V<MHeV+1; V++) {
688: row[0] = &f[xi].He[1] - rowstart;
689: row[1] = &f[xi].V[V] - rowstart;
690: row[2] = &fHeV[V][1] - rowstart;
691: col[0] = &cHeV[V][1] - colstart;
692: val[0] = 1000*ctx->dissociationScale;
693: val[1] = 1000*ctx->dissociationScale;
694: val[2] = -1000*ctx->dissociationScale;
695: MatSetValuesLocal(J,3,row,1,col,val,ADD_VALUES);
696: }
697:
698: /* He[He]-V[V] -> He[He-1]-V[V] + He[1] */
699: for (V=2; V<MHeV+1; V++) {
700: for (He=2; He<NHeV[V]+1; He++) {
701: row[0] = &f[xi].He[1] - rowstart;
702: row[1] = &fHeV[V][He-1] - rowstart;
703: row[2] = &fHeV[V][He] - rowstart;
704: col[0] = &cHeV[V][He] - colstart;
705: val[0] = 1000*ctx->dissociationScale;
706: val[1] = 1000*ctx->dissociationScale;
707: val[2] = -1000*ctx->dissociationScale;
708: MatSetValuesLocal(J,3,row,1,col,val,ADD_VALUES);
709: }
710: }
711:
712: /* He[He]-V[V] -> He[He]-V[V-1] + V[1] */
713: for (V=2; V<MHeV+1; V++) {
714: for (He=2; He<NHeV[V-1]+1; He++) {
715: row[0] = &f[xi].V[1] - rowstart;
716: row[1] = &fHeV[V-1][He] - rowstart;
717: row[2] = &fHeV[V][He] - rowstart;
718: col[0] = &cHeV[V][He] - colstart;
719: val[0] = 1000*ctx->dissociationScale;
720: val[1] = 1000*ctx->dissociationScale;
721: val[2] = -1000*ctx->dissociationScale;
722: MatSetValuesLocal(J,3,row,1,col,val,ADD_VALUES);
723: }
724: }
725:
726: /* He[He]-V[V] -> He[He]-V[V+1] + I[1] */
727: for (V=1; V<MHeV; V++) {
728: for (He=1; He<NHeV[V]+1; He++) {
729: row[0] = &fHeV[V+1][He] - rowstart;
730: row[1] = &f[xi].I[1] - rowstart;
731: row[2] = &fHeV[V][He] - rowstart;
732: col[0] = &cHeV[V][He] - colstart;
733: val[0] = 1000*ctx->dissociationScale;
734: val[1] = 1000*ctx->dissociationScale;
735: val[2] = -1000*ctx->dissociationScale;
736: MatSetValuesLocal(J,3,row,1,col,val,ADD_VALUES);
737: }
738: }
739: }
740: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
741: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
742: MatSetOption(J,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);
743: MatStoreValues(J);
744: MatSetFromOptions(J);
745: initialized = PETSC_TRUE;
746: } else {
747: MatRetrieveValues(J);
748: }
750: /*
751: Loop over grid points computing Jacobian terms for each grid point for reaction terms
752: */
753: for (xi=xs; xi<xs+xm; xi++) {
754: x = xi*hx;
755: cHeVInitialize(&c[xi].He[1],(PetscScalar**)cHeV);
756: cHeVInitialize(&f[xi].He[1],(PetscScalar**)fHeV);
757: /* ----------------------------------------------------------------
758: ---- Compute reaction terms that can create a cluster of given size
759: */
760: /* He[He] + He[he] -> He[He+he] */
761: for (He=2; He<NHe+1; He++) {
762: /* compute all pairs of clusters of smaller size that can combine to create a cluster of size He,
763: remove the upper half since they are symmetric to the lower half of the pairs. For example
764: when He = 5 (cluster size 5) the pairs are
765: 1 4
766: 2 2
767: 3 2 these last two are not needed in the sum since they repeat from above
768: 4 1 this is why he < (He/2) + 1 */
769: for (he=1; he<(He/2)+1; he++) {
770: row[0] = &f[xi].He[He] - rowstart;
771: row[1] = &f[xi].He[he] - rowstart;
772: row[2] = &f[xi].He[He-he] - rowstart;
773: col[0] = &c[xi].He[he] - colstart;
774: col[1] = &c[xi].He[He-he] - colstart;
775: val[0] = ctx->reactionScale*c[xi].He[He-he];
776: val[1] = ctx->reactionScale*c[xi].He[he];
777: val[2] = -ctx->reactionScale*c[xi].He[He-he];
778: val[3] = -ctx->reactionScale*c[xi].He[he];
779: val[4] = -ctx->reactionScale*c[xi].He[He-he];
780: val[5] = -ctx->reactionScale*c[xi].He[he];
781: MatSetValuesLocal(J,3,row,2,col,val,ADD_VALUES);
782: }
783: }
785: /* V[V] + V[v] -> V[V+v] */
786: for (V=2; V<NV+1; V++) {
787: for (v=1; v<(V/2)+1; v++) {
788: row[0] = &f[xi].V[V] - rowstart;
789: row[1] = &f[xi].V[v] - rowstart;
790: row[2] = &f[xi].V[V-v] - rowstart;
791: col[0] = &c[xi].V[v] - colstart;
792: col[1] = &c[xi].V[V-v] - colstart;
793: val[0] = ctx->reactionScale*c[xi].V[V-v];
794: val[1] = ctx->reactionScale*c[xi].V[v];
795: val[2] = -ctx->reactionScale*c[xi].V[V-v];
796: val[3] = -ctx->reactionScale*c[xi].V[v];
797: val[4] = -ctx->reactionScale*c[xi].V[V-v];
798: val[5] = -ctx->reactionScale*c[xi].V[v];
799: MatSetValuesLocal(J,3,row,2,col,val,ADD_VALUES);
800: }
801: }
803: /* I[I] + I[i] -> I[I+i] */
804: for (I=2; I<NI+1; I++) {
805: for (i=1; i<(I/2)+1; i++) {
806: row[0] = &f[xi].I[I] - rowstart;
807: row[1] = &f[xi].I[i] - rowstart;
808: row[2] = &f[xi].I[I-i] - rowstart;
809: col[0] = &c[xi].I[i] - colstart;
810: col[1] = &c[xi].I[I-i] - colstart;
811: val[0] = ctx->reactionScale*c[xi].I[I-i];
812: val[1] = ctx->reactionScale*c[xi].I[i];
813: val[2] = -ctx->reactionScale*c[xi].I[I-i];
814: val[3] = -ctx->reactionScale*c[xi].I[i];
815: val[4] = -ctx->reactionScale*c[xi].I[I-i];
816: val[5] = -ctx->reactionScale*c[xi].I[i];
817: MatSetValuesLocal(J,3,row,2,col,val,ADD_VALUES);
818: }
819: }
821: /* He[1] + V[1] -> He[1]-V[1] */
822: row[0] = &fHeV[1][1] - rowstart;
823: row[1] = &f[xi].He[1] - rowstart;
824: row[2] = &f[xi].V[1] - rowstart;
825: col[0] = &c[xi].He[1] - colstart;
826: col[1] = &c[xi].V[1] - colstart;
827: val[0] = 1000*ctx->reactionScale*c[xi].V[1];
828: val[1] = 1000*ctx->reactionScale*c[xi].He[1];
829: val[2] = -1000*ctx->reactionScale*c[xi].V[1];
830: val[3] = -1000*ctx->reactionScale*c[xi].He[1];
831: val[4] = -1000*ctx->reactionScale*c[xi].V[1];
832: val[5] = -1000*ctx->reactionScale*c[xi].He[1];
833: MatSetValuesLocal(J,3,row,2,col,val,ADD_VALUES);
835: /* He[He]-V[V] + He[he] -> He[He+he]-V[V] */
836: for (V=1; V<MHeV+1; V++) {
837: for (He=1; He<NHeV[V]; He++) {
838: for (he=1; he+He<NHeV[V]+1; he++) {
839: row[0] = &fHeV[V][He+he] - rowstart;
840: row[1] = &f[xi].He[he] - rowstart;
841: row[2] = &fHeV[V][He] - rowstart;
842: col[0] = &c[xi].He[he] - colstart;
843: col[1] = &cHeV[V][He] - colstart;
844: val[0] = ctx->reactionScale*cHeV[V][He];
845: val[1] = ctx->reactionScale*c[xi].He[he];
846: val[2] = -ctx->reactionScale*cHeV[V][He];
847: val[3] = -ctx->reactionScale*c[xi].He[he];
848: val[4] = -ctx->reactionScale*cHeV[V][He];
849: val[5] = -ctx->reactionScale*c[xi].He[he];
850: MatSetValuesLocal(J,3,row,2,col,val,ADD_VALUES);
851: }
852: }
853: }
855: /* He[He]-V[V] + V[1] -> He[He][V+1] */
856: for (V=1; V<MHeV; V++) {
857: for (He=1; He<NHeV[V+1]; He++) {
858: row[0] = &fHeV[V+1][He] - rowstart;
859: row[1] = &f[xi].V[1] - rowstart;
860: row[2] = &fHeV[V][He] - rowstart;
861: col[0] = &c[xi].V[1] - colstart;
862: col[1] = &cHeV[V][He] - colstart;
863: val[0] = ctx->reactionScale*cHeV[V][He];
864: val[1] = ctx->reactionScale*c[xi].V[1];
865: val[2] = -ctx->reactionScale*cHeV[V][He];
866: val[3] = -ctx->reactionScale*c[xi].V[1];
867: val[4] = -ctx->reactionScale*cHeV[V][He];
868: val[5] = -ctx->reactionScale*c[xi].V[1];
869: MatSetValuesLocal(J,3,row,2,col,val,ADD_VALUES);
870: }
871: }
873: /* He[He]-V[V] + He[he]-V[v] -> He[He+he][V+v] */
874: /* Currently the reaction rates for this are zero */
877: /* V[V] + I[I] -> V[V-I] if V > I else I[I-V] */
878: for (V=1; V<NV+1; V++) {
879: for (I=1; I<PetscMin(V,NI); I++) {
880: row[0] = &f[xi].V[V-I] - rowstart;
881: row[1] = &f[xi].V[V] - rowstart;
882: row[2] = &f[xi].I[I] - rowstart;
883: col[0] = &c[xi].V[V] - colstart;
884: col[1] = &c[xi].I[I] - colstart;
885: val[0] = ctx->reactionScale*c[xi].I[I];
886: val[1] = ctx->reactionScale*c[xi].V[V];
887: val[2] = -ctx->reactionScale*c[xi].I[I];
888: val[3] = -ctx->reactionScale*c[xi].V[V];
889: val[4] = -ctx->reactionScale*c[xi].I[I];
890: val[5] = -ctx->reactionScale*c[xi].V[V];
891: MatSetValuesLocal(J,3,row,2,col,val,ADD_VALUES);
892: }
893: for (I=V+1; I<NI+1; I++) {
894: row[0] = &f[xi].I[I-V] - rowstart;
895: row[1] = &f[xi].V[V] - rowstart;
896: row[2] = &f[xi].I[I] - rowstart;
897: col[0] = &c[xi].V[V] - colstart;
898: col[1] = &c[xi].I[I] - colstart;
899: val[0] = ctx->reactionScale*c[xi].I[I];
900: val[1] = ctx->reactionScale*c[xi].V[V];
901: val[2] = -ctx->reactionScale*c[xi].I[I];
902: val[3] = -ctx->reactionScale*c[xi].V[V];
903: val[4] = -ctx->reactionScale*c[xi].I[I];
904: val[5] = -ctx->reactionScale*c[xi].V[V];
905: MatSetValuesLocal(J,3,row,2,col,val,ADD_VALUES);
906: }
907: }
908: }
910: /*
911: Restore vectors
912: */
913: c = (Concentrations*)(((PetscScalar*)c)+1);
914: DMDAVecRestoreArray(da,localC,&c);
915: f = (Concentrations*)(((PetscScalar*)f)+1);
916: DMDAVecRestoreArray(da,C,&f);
917: DMRestoreLocalVector(da,&localC);
918: cHeVDestroy((PetscScalar**)cHeV);
919: cHeVDestroy((PetscScalar**)fHeV);
921: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
922: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
923: if (A != J) {
924: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
925: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
926: }
927: return(0);
928: }
932: /*
933: Determines the nonzero structure within the diagonal blocks of the Jacobian that represent coupling resulting from reactions and
934: dissasociations of the clusters
935: */
936: PetscErrorCode GetDfill(PetscInt *dfill, void *ptr)
937: {
938: PetscInt He,he,V,v,I,i,j,k,rows[3],cols[2];
939: Concentrations *c;
940: PetscScalar *idxstart,**cHeV;
943: /* ensure fill for the diagonal of matrix */
944: for (i=0; i<(DOF); i++) {
945: dfill[i*DOF + i] = 1;
946: }
948: /*
949: c is never used except for computing offsets between variables which are used to fill the non-zero
950: structure of the matrix
951: */
952: PetscMalloc(sizeof(Concentrations),&c);
953: c = (Concentrations*)(((PetscScalar*)c)-1);
954: cHeVCreate(&cHeV);
955: cHeVInitialize(&c->He[1],cHeV);
956: idxstart = (PetscScalar*)&c->He[1];
958: /* -------------------------------------------------------------------------
959: ---- Compute dissociation terms that removes an item from a cluster
960: I assume dissociation means losing only a single item from a cluster
961: I cannot tell from the notes if clusters can break up into any sub-size.
962: */
963: /* He[He] -> He[He-1] + He[1] */
964: for (He=2; He<NHe+1; He++) {
965: rows[0] = &c->He[He-1] - idxstart;
966: rows[1] = &c->He[1] - idxstart;
967: rows[2] = &c->He[He] - idxstart;
968: cols[0] = &c->He[He] - idxstart;
969: for (j=0; j<3; j++) {
970: dfill[rows[j]*DOF + cols[0]] = 1;
971: }
972: }
974: /* V[V] -> V[V-1] + V[1] */
975: for (V=2; V<NV+1; V++) {
976: rows[0] = &c->V[V] - idxstart;
977: rows[1] = &c->V[1] - idxstart;
978: rows[2] = &c->V[V-1] - idxstart;
979: cols[0] = &c->V[V] - idxstart;
980: for (j=0; j<3; j++) {
981: dfill[rows[j]*DOF + cols[0]] = 1;
982: }
983: }
984:
985: /* I[I] -> I[I-1] + I[1] */
986: for (I=2; I<NI+1; I++) {
987: rows[0] = &c->I[I] - idxstart;
988: rows[1] = &c->I[1] - idxstart;
989: rows[2] = &c->I[I-1] - idxstart;
990: cols[0] = &c->I[I] - idxstart;
991: for (j=0; j<3; j++) {
992: dfill[rows[j]*DOF + cols[0]] = 1;
993: }
994: }
995:
996: /* He[He]-V[1] -> He[He] + V[1] */
997: for (He=1; He<NHeV[1]+1; He++) {
998: rows[0] = &c->He[He] - idxstart;
999: rows[1] = &c->V[1] - idxstart;
1000: rows[2] = &cHeV[1][He] - idxstart;
1001: cols[0] = &cHeV[1][He] - idxstart;
1002: for (j=0; j<3; j++) {
1003: dfill[rows[j]*DOF + cols[0]] = 1;
1004: }
1005: }
1006:
1007: /* He[1]-V[V] -> He[1] + V[V] */
1008: for (V=2; V<MHeV+1; V++) {
1009: rows[0] = &c->He[1] - idxstart;
1010: rows[1] = &c->V[V] - idxstart;
1011: rows[2] = &cHeV[V][1] - idxstart;
1012: cols[0] = &cHeV[V][1] - idxstart;
1013: for (j=0; j<3; j++) {
1014: dfill[rows[j]*DOF + cols[0]] = 1;
1015: }
1016: }
1017:
1018: /* He[He]-V[V] -> He[He-1]-V[V] + He[1] */
1019: for (V=2; V<MHeV+1; V++) {
1020: for (He=2; He<NHeV[V]+1; He++) {
1021: rows[0] = &c->He[1] - idxstart;
1022: rows[1] = &cHeV[V][He] - idxstart;
1023: rows[2] = &cHeV[V][He-1] - idxstart;
1024: cols[0] = &cHeV[V][He] - idxstart;
1025: for (j=0; j<3; j++) {
1026: dfill[rows[j]*DOF + cols[0]] = 1;
1027: }
1028: }
1029: }
1030:
1031: /* He[He]-V[V] -> He[He]-V[V-1] + V[1] */
1032: for (V=2; V<MHeV+1; V++) {
1033: for (He=2; He<NHeV[V-1]+1; He++) {
1034: rows[0] = &c->V[1] - idxstart;
1035: rows[1] = &cHeV[V][He] - idxstart;
1036: rows[2] = &cHeV[V-1][He] - idxstart;
1037: cols[0] = &cHeV[V][He] - idxstart;
1038: for (j=0; j<3; j++) {
1039: dfill[rows[j]*DOF + cols[0]] = 1;
1040: }
1041: }
1042: }
1043:
1044: /* He[He]-V[V] -> He[He]-V[V+1] + I[1] */
1045: for (V=1; V<MHeV; V++) {
1046: for (He=1; He<NHeV[V]+1; He++) {
1047: rows[0] = &c->I[1] - idxstart;
1048: rows[1] = &cHeV[V+1][He] - idxstart;
1049: rows[2] = &cHeV[V][He] - idxstart;
1050: cols[0] = &cHeV[V][He] - idxstart;
1051: for (j=0; j<3; j++) {
1052: dfill[rows[j]*DOF + cols[0]] = 1;
1053: }
1054: }
1055: }
1057: /* These are the reaction terms in the diagonal block */
1058: for (He=2; He<NHe+1; He++) {
1059: for (he=1; he<(He/2)+1; he++) {
1060: rows[0] = &c->He[He] - idxstart;
1061: rows[1] = &c->He[he] - idxstart;
1062: rows[2] = &c->He[He-he] - idxstart;
1063: cols[0] = &c->He[he] - idxstart;
1064: cols[1] = &c->He[He-he] - idxstart;
1065: for (j=0; j<3; j++) {
1066: for (k=0; k<2; k++) {
1067: dfill[rows[j]*DOF + cols[k]] = 1;
1068: }
1069: }
1070: }
1071: }
1073: /* V[V] + V[v] -> V[V+v] */
1074: for (V=2; V<NV+1; V++) {
1075: for (v=1; v<(V/2)+1; v++) {
1076: rows[0] = &c->V[V] - idxstart;
1077: rows[1] = &c->V[v] - idxstart;
1078: rows[2] = &c->V[V-v] - idxstart;
1079: cols[0] = &c->V[v] - idxstart;
1080: cols[1] = &c->V[V-v] - idxstart;
1081: for (j=0; j<3; j++) {
1082: for (k=0; k<2; k++) {
1083: dfill[rows[j]*DOF + cols[k]] = 1;
1084: }
1085: }
1086: }
1087: }
1088:
1089: /* I[I] + I[i] -> I[I+i] */
1090: for (I=2; I<NI+1; I++) {
1091: for (i=1; i<(I/2)+1; i++) {
1092: rows[0] = &c->I[I] - idxstart;
1093: rows[1] = &c->I[i] - idxstart;
1094: rows[2] = &c->I[I-i] - idxstart;
1095: cols[0] = &c->I[i] - idxstart;
1096: cols[1] = &c->I[I-i] - idxstart;
1097: for (j=0; j<3; j++) {
1098: for (k=0; k<2; k++) {
1099: dfill[rows[j]*DOF + cols[k]] = 1;
1100: }
1101: }
1102: }
1103: }
1104:
1105: /* He[1] + V[1] -> He[1]-V[1] */
1106: rows[0] = &cHeV[1][1] - idxstart;
1107: rows[1] = &c->He[1] - idxstart;
1108: rows[2] = &c->V[1] - idxstart;
1109: cols[0] = &c->He[1] - idxstart;
1110: cols[1] = &c->V[1] - idxstart;
1111: for (j=0; j<3; j++) {
1112: for (k=0; k<2; k++) {
1113: dfill[rows[j]*DOF + cols[k]] = 1;
1114: }
1115: }
1116:
1117: /* He[He]-V[V] + He[he] -> He[He+he]-V[V] */
1118: for (V=1; V<MHeV+1; V++) {
1119: for (He=1; He<NHeV[V]; He++) {
1120: for (he=1; he+He<NHeV[V]+1; he++) {
1121: rows[0] = &cHeV[V][He+he] - idxstart;
1122: rows[1] = &c->He[he] - idxstart;
1123: rows[2] = &cHeV[V][He] - idxstart;
1124: cols[0] = &cHeV[V][He] - idxstart;
1125: cols[1] = &c->He[he] - idxstart;
1126: for (j=0; j<3; j++) {
1127: for (k=0; k<2; k++) {
1128: dfill[rows[j]*DOF + cols[k]] = 1;
1129: }
1130: }
1131: }
1132: }
1133: }
1134: /* He[He]-V[V] + V[1] -> He[He][V+1] */
1135: for (V=1; V<MHeV; V++) {
1136: for (He=1; He<NHeV[V+1]; He++) {
1137: rows[0] = &cHeV[V+1][He] - idxstart;
1138: rows[1] = &c->V[1] - idxstart;
1139: rows[2] = &cHeV[V][He] - idxstart;
1140: cols[0] = &cHeV[V][He] - idxstart;
1141: cols[1] = &c->V[1] - idxstart;
1142: for (j=0; j<3; j++) {
1143: for (k=0; k<2; k++) {
1144: dfill[rows[j]*DOF + cols[k]] = 1;
1145: }
1146: }
1147: }
1148: }
1150: /* He[He]-V[V] + He[he]-V[v] -> He[He+he][V+v] */
1151: /* Currently the reaction rates for this are zero */
1152:
1153: /* V[V] + I[I] -> V[V-I] if V > I else I[I-V] */
1154: for (V=1; V<NV+1; V++) {
1155: for (I=1; I<PetscMin(V,NI); I++) {
1156: rows[0] = &c->V[V-I] - idxstart;
1157: rows[1] = &c->V[V] - idxstart;
1158: rows[2] = &c->I[I] - idxstart;
1159: cols[0] = &c->V[V] - idxstart;
1160: cols[1] = &c->I[I] - idxstart;
1161: for (j=0; j<3; j++) {
1162: for (k=0; k<2; k++) {
1163: dfill[rows[j]*DOF + cols[k]] = 1;
1164: }
1165: }
1166: }
1167: for (I=V+1; I<NI+1; I++) {
1168: rows[0] = &c->I[I-V] - idxstart;
1169: rows[1] = &c->V[V] - idxstart;
1170: rows[2] = &c->I[I] - idxstart;
1171: cols[0] = &c->V[V] - idxstart;
1172: cols[1] = &c->I[I] - idxstart;
1173: for (j=0; j<3; j++) {
1174: for (k=0; k<2; k++) {
1175: dfill[rows[j]*DOF + cols[k]] = 1;
1176: }
1177: }
1178: }
1179: }
1181: c = (Concentrations*)(((PetscScalar*)c)+1);
1182: cHeVDestroy(cHeV);
1183: PetscFree(c);
1184: return(0);
1185: }
1186: /* ------------------------------------------------------------------- */
1188: typedef struct {
1189: DM Heda,Vda,HeVda; /* defines the 2d layout of the He subvector */
1190: Vec He,V,HeV;
1191: VecScatter Hescatter,Vscatter,HeVscatter;
1192: PetscViewer Heviewer,Vviewer,HeVviewer;
1193: } MyMonitorCtx;
1197: /*
1198: Display He as a function of space and cluster size for each time step
1199: */
1200: PetscErrorCode MyMonitorMonitor(TS ts,PetscInt timestep,PetscReal time,Vec solution, void *ictx)
1201: {
1202: MyMonitorCtx *ctx = (MyMonitorCtx*)ictx;
1206: VecScatterBegin(ctx->Hescatter,solution,ctx->He,INSERT_VALUES,SCATTER_FORWARD);
1207: VecScatterEnd(ctx->Hescatter,solution,ctx->He,INSERT_VALUES,SCATTER_FORWARD);
1208: VecView(ctx->He,ctx->Heviewer);
1210: VecScatterBegin(ctx->Vscatter,solution,ctx->V,INSERT_VALUES,SCATTER_FORWARD);
1211: VecScatterEnd(ctx->Vscatter,solution,ctx->V,INSERT_VALUES,SCATTER_FORWARD);
1212: VecView(ctx->V,ctx->Vviewer);
1214: VecScatterBegin(ctx->HeVscatter,solution,ctx->HeV,INSERT_VALUES,SCATTER_FORWARD);
1215: VecScatterEnd(ctx->HeVscatter,solution,ctx->HeV,INSERT_VALUES,SCATTER_FORWARD);
1216: VecView(ctx->HeV,ctx->HeVviewer);
1217: return(0);
1218: }
1222: /*
1223: Frees all data structures associated with the monitor
1224: */
1225: PetscErrorCode MyMonitorDestroy(void **ictx)
1226: {
1227: MyMonitorCtx **ctx = (MyMonitorCtx**)ictx;
1231: VecScatterDestroy(&(*ctx)->Hescatter);
1232: VecDestroy(&(*ctx)->He);
1233: DMDestroy(&(*ctx)->Heda);
1234: PetscViewerDestroy(&(*ctx)->Heviewer);
1236: VecScatterDestroy(&(*ctx)->Vscatter);
1237: VecDestroy(&(*ctx)->V);
1238: DMDestroy(&(*ctx)->Vda);
1239: PetscViewerDestroy(&(*ctx)->Vviewer);
1241: VecScatterDestroy(&(*ctx)->HeVscatter);
1242: VecDestroy(&(*ctx)->HeV);
1243: DMDestroy(&(*ctx)->HeVda);
1244: PetscViewerDestroy(&(*ctx)->HeVviewer);
1245: PetscFree(*ctx);
1246: return(0);
1247: }
1251: /*
1252: Sets up a monitor that will display He as a function of space and cluster size for each time step
1253: */
1254: PetscErrorCode MyMonitorSetUp(TS ts)
1255: {
1256: DM da;
1258: PetscInt xi,xs,xm,*idx,M,xj,cnt = 0;
1259: const PetscInt *lx;
1260: Vec C;
1261: MyMonitorCtx *ctx;
1262: PetscBool flg;
1263: IS is;
1264: char ycoor[32];
1265: PetscReal valuebounds[4] = {0, 1.2, 0, 1.2};
1268: PetscOptionsHasName(NULL,"-mymonitor",&flg);
1269: if (!flg) return(0);
1271: TSGetDM(ts,&da);
1272: PetscNew(MyMonitorCtx,&ctx);
1273: PetscNew(&ctx);
1274: PetscViewerDrawOpen(PetscObjectComm((PetscObject)da),NULL,"",PETSC_DECIDE,PETSC_DECIDE,600,400,&ctx->viewer);
1276: /* setup visualization for He */
1277: PetscViewerDrawOpen(PetscObjectComm((PetscObject)da),NULL,"",PETSC_DECIDE,PETSC_DECIDE,600,400,&ctx->Heviewer);
1278: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
1279: 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);
1280: DMDAGetOwnershipRanges(da,&lx,NULL,NULL);
1281: DMDACreate2d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,M,NHe,PETSC_DETERMINE,1,1,1,lx,NULL,&ctx->Heda);
1282: DMDASetFieldName(ctx->Heda,0,"He");
1283: DMDASetCoordinateName(ctx->Heda,0,"X coordinate direction");
1284: PetscSNPrintf(ycoor,32,"%D ... Cluster size ... 1",NHe);
1285: DMDASetCoordinateName(ctx->Heda,1,ycoor);
1286: DMCreateGlobalVector(ctx->Heda,&ctx->He);
1287: PetscMalloc(NHe*xm*sizeof(PetscInt),&idx);
1288: cnt = 0;
1289: for (xj=0; xj<NHe; xj++) {
1290: for (xi=xs; xi<xs+xm; xi++) {
1291: idx[cnt++] = DOF*xi + xj;
1292: }
1293: }
1294: ISCreateGeneral(PetscObjectComm((PetscObject)ts),NHe*xm,idx,PETSC_OWN_POINTER,&is);
1295: TSGetSolution(ts,&C);
1296: VecScatterCreate(C,is,ctx->He,NULL,&ctx->Hescatter);
1297: ISDestroy(&is);
1298: /* sets the bounds on the contour plot values so the colors mean the same thing for different timesteps */
1299: PetscViewerDrawSetBounds(ctx->Heviewer,2,valuebounds);
1301: /* setup visualization for V */
1302: PetscViewerDrawOpen(PetscObjectComm((PetscObject)da),NULL,"",PETSC_DECIDE,PETSC_DECIDE,600,400,&ctx->Vviewer);
1303: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
1304: 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);
1305: DMDAGetOwnershipRanges(da,&lx,NULL,NULL);
1306: DMDACreate2d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,M,NV,PETSC_DETERMINE,1,1,1,lx,NULL,&ctx->Vda);
1307: DMDASetFieldName(ctx->Vda,0,"V");
1308: DMDASetCoordinateName(ctx->Vda,0,"X coordinate direction");
1309: PetscSNPrintf(ycoor,32,"%D ... Cluster size ... 1",NV);
1310: DMDASetCoordinateName(ctx->Vda,1,ycoor);
1311: DMCreateGlobalVector(ctx->Vda,&ctx->V);
1312: PetscMalloc(NV*xm*sizeof(PetscInt),&idx);
1313: DMCreateGlobalVector(ctx->da,&ctx->He);
1314: PetscMalloc1(2*N*xm,&idx);
1315: cnt = 0;
1316: for (xj=0; xj<NV; xj++) {
1317: for (xi=xs; xi<xs+xm; xi++) {
1318: idx[cnt++] = NHe + DOF*xi + xj;
1319: }
1320: }
1321: ISCreateGeneral(PetscObjectComm((PetscObject)ts),NV*xm,idx,PETSC_OWN_POINTER,&is);
1322: TSGetSolution(ts,&C);
1323: VecScatterCreate(C,is,ctx->V,NULL,&ctx->Vscatter);
1324: ISDestroy(&is);
1325: /* sets the bounds on the contour plot values so the colors mean the same thing for different timesteps */
1326: PetscViewerDrawSetBounds(ctx->Vviewer,2,valuebounds);
1328: /* setup visualization for HeV[1][] */
1329: PetscViewerDrawOpen(PetscObjectComm((PetscObject)da),NULL,"",PETSC_DECIDE,PETSC_DECIDE,600,400,&ctx->HeVviewer);
1330: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
1331: 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);
1332: DMDAGetOwnershipRanges(da,&lx,NULL,NULL);
1333: DMDACreate2d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,M,NHeV[1],PETSC_DETERMINE,1,1,1,lx,NULL,&ctx->HeVda);
1334: DMDASetFieldName(ctx->HeVda,0,"HeV[1][]");
1335: DMDASetCoordinateName(ctx->HeVda,0,"X coordinate direction");
1336: PetscSNPrintf(ycoor,32,"%D ... Cluster size ... 1",NHeV[1]);
1337: DMDASetCoordinateName(ctx->HeVda,1,ycoor);
1338: DMCreateGlobalVector(ctx->HeVda,&ctx->HeV);
1339: PetscMalloc(NHeV[1]*xm*sizeof(PetscInt),&idx);
1340: cnt = 0;
1341: for (xj=0; xj<NHeV[1]; xj++) {
1342: for (xi=xs; xi<xs+xm; xi++) {
1343: idx[cnt++] = NHe + NV + NI + DOF*xi + xj;
1344: }
1345: }
1346: ISCreateGeneral(PetscObjectComm((PetscObject)ts),NHeV[1]*xm,idx,PETSC_OWN_POINTER,&is);
1347: TSGetSolution(ts,&C);
1348: VecScatterCreate(C,is,ctx->HeV,NULL,&ctx->HeVscatter);
1349: ISDestroy(&is);
1350: /* sets the bounds on the contour plot values so the colors mean the same thing for different timesteps */
1351: PetscViewerDrawSetBounds(ctx->HeVviewer,2,valuebounds);
1353: TSMonitorSet(ts,MyMonitorMonitor,ctx,MyMonitorDestroy);
1354: return(0);
1355: }
1359: PetscErrorCode MyLoadData(MPI_Comm comm,const char *filename)
1360: {
1362: FILE *fp;
1363: char buff[256];
1364: PetscInt He,V,I,lc = 0;
1365: char Hebindstr[32],Vbindstr[32],Ibindstr[32],trapbindstr[32],*sharp;
1366: PetscReal Hebind,Vbind,Ibind,trapbind;
1369: PetscFOpen(comm,filename,"r",&fp);
1370: PetscSynchronizedFGets(comm,fp,256,buff);
1371: while (buff[0]) {
1372: PetscStrchr(buff,'#',&sharp);
1373: if (!sharp) {
1374: sscanf(buff,"%d %d %d %s %s %s %s",&He,&V,&I,Hebindstr,Vbindstr,Ibindstr,trapbindstr);
1375: Hebind = strtod(Hebindstr,NULL);
1376: Vbind = strtod(Vbindstr,NULL);
1377: Ibind = strtod(Ibindstr,NULL);
1378: trapbind = strtod(trapbindstr,NULL);
1379: if (V <= NV) {
1380: if (He > NHe && V == 0 && I == 0) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Recompile with correct NHe %d %d",He,NHe);
1381: if (He == 0 && V > NV && I == 0) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Recompile with correct V %d %d",V,NV);
1382: if (He == 0 && V == 0 && I > NI) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Recompile with correct NI %d %d",I,NI);
1383: if (lc++ > DOF) SETERRQ4(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Recompile with correct NHe %d NV %d NI %d MNHeV %",NHe,NV,NI,MNHeV);
1384: if (He > 0 && V > 0) { /* assumes the He are sorted in increasing order */
1385: NHeV[V] = He;
1386: }
1387: }
1388: }
1389: PetscSynchronizedFGets(comm,fp,256,buff);
1390: }
1391: if (lc != DOF) SETERRQ5(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Recompile with correct NHe %d NV %d NI %d MNHeV %d Actual DOF %d",NHe,NV,NI,MNHeV,lc);
1392: return(0);
1393: }