Actual source code: ex10.c
petsc-3.5.4 2015-05-23
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,"-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: TSSetFromOptions(ts);
156: MyMonitorSetUp(ts);
158: InitialConditions(da,C);
160: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
161: Solve the ODE system
162: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
163: TSSolve(ts,C);
165: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
166: Free work space.
167: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
168: VecDestroy(&C);
169: TSDestroy(&ts);
170: DMDestroy(&da);
171: PetscFinalize();
172: return(0);
173: }
175: /*
176: cHeV is "trick" to allow easy accessing of the values in the HeV portion of the Concentrations.
177: cHeV[i] points to the beginning of each row of HeV[] with V indexing starting a 1.
179: */
182: PetscErrorCode cHeVCreate(PetscReal ***cHeV)
183: {
187: PetscMalloc(MHeV*sizeof(PetscScalar),cHeV);
188: (*cHeV)--;
189: return(0);
190: }
194: PetscErrorCode cHeVInitialize(const PetscScalar *start,PetscReal **cHeV)
195: {
196: PetscInt i;
199: cHeV[1] = ((PetscScalar*) start) - 1 + NHe + NV + NI;
200: for (i=1; i<MHeV; i++) {
201: cHeV[i+1] = cHeV[i] + NHeV[i];
202: }
203: return(0);
204: }
208: PetscErrorCode cHeVDestroy(PetscReal **cHeV)
209: {
213: cHeV++;
214: PetscFree(cHeV);
215: return(0);
216: }
218: /* ------------------------------------------------------------------- */
221: PetscErrorCode InitialConditions(DM da,Vec C)
222: {
224: PetscInt i,I,He,V,xs,xm,Mx,cnt = 0;
225: Concentrations *c;
226: PetscReal hx,x,**cHeV;
227: char string[16];
230: 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);
231: hx = 1.0/(PetscReal)(Mx-1);
233: /* Name each of the concentrations */
234: for (He=1; He<NHe+1; He++) {
235: PetscSNPrintf(string,16,"%d-He",He);
236: DMDASetFieldName(da,cnt++,string);
237: }
238: for (V=1; V<NV+1; V++) {
239: PetscSNPrintf(string,16,"%d-V",V);
240: DMDASetFieldName(da,cnt++,string);
241: }
242: for (I=1; I<NI+1; I++) {
243: PetscSNPrintf(string,16,"%d-I",I);
244: DMDASetFieldName(da,cnt++,string);
245: }
246: for (He=1; He<MHeV+1; He++) {
247: for (V=1; V<NHeV[He]+1; V++) {
248: PetscSNPrintf(string,16,"%d-He-%d-V",He,V);
249: DMDASetFieldName(da,cnt++,string);
250: }
251: }
253: /*
254: Get pointer to vector data
255: */
256: DMDAVecGetArray(da,C,&c);
257: /* Shift the c pointer to allow accessing with index of 1, instead of 0 */
258: c = (Concentrations*)(((PetscScalar*)c)-1);
260: /*
261: Get local grid boundaries
262: */
263: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
265: /*
266: Compute function over the locally owned part of the grid
267: */
268: cHeVCreate(&cHeV);
269: for (i=xs; i<xs+xm; i++) {
270: x = i*hx;
271: for (He=1; He<NHe+1; He++) c[i].He[He] = 0.0;
272: for (V=1; V<NV+1; V++) c[i].V[V] = 1.0;
273: for (I=1; I <NI+1; I++) c[i].I[I] = 1.0;
274: cHeVInitialize(&c[i].He[1],cHeV);
275: for (V=1; V<MHeV+1; V++) {
276: for (He=1; He<NHeV[V]+1; He++) cHeV[V][He] = 0.0;
277: }
278: }
279: cHeVDestroy(cHeV);
281: /*
282: Restore vectors
283: */
284: c = (Concentrations*)(((PetscScalar*)c)+1);
285: DMDAVecRestoreArray(da,C,&c);
286: return(0);
287: }
289: /* ------------------------------------------------------------------- */
292: /*
293: RHSFunction - Evaluates nonlinear function that defines the ODE
295: Input Parameters:
296: . ts - the TS context
297: . U - input vector
298: . ptr - optional user-defined context
300: Output Parameter:
301: . F - function values
302: */
303: PetscErrorCode RHSFunction(TS ts,PetscReal ftime,Vec C,Vec F,void *ptr)
304: {
305: AppCtx *ctx = (AppCtx*) ptr;
306: DM da;
308: PetscInt xi,Mx,xs,xm,He,he,V,v,I,i;
309: PetscReal hx,sx,x,**cHeV,**fHeV;
310: Concentrations *c,*f;
311: Vec localC;
314: TSGetDM(ts,&da);
315: DMGetLocalVector(da,&localC);
316: 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);
317: hx = 8.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
318: cHeVCreate(&cHeV);
319: cHeVCreate(&fHeV);
321: /*
322: Scatter ghost points to local vector,using the 2-step process
323: DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
324: By placing code between these two statements, computations can be
325: done while messages are in transition.
326: */
327: DMGlobalToLocalBegin(da,C,INSERT_VALUES,localC);
328: DMGlobalToLocalEnd(da,C,INSERT_VALUES,localC);
330: VecSet(F,0.0);
332: /*
333: Get pointers to vector data
334: */
335: DMDAVecGetArray(da,localC,&c);
336: /* Shift the c pointer to allow accessing with index of 1, instead of 0 */
337: c = (Concentrations*)(((PetscScalar*)c)-1);
338: DMDAVecGetArray(da,F,&f);
339: f = (Concentrations*)(((PetscScalar*)f)-1);
341: /*
342: Get local grid boundaries
343: */
344: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
346: /*
347: Loop over grid points computing ODE terms for each grid point
348: */
349: for (xi=xs; xi<xs+xm; xi++) {
350: x = xi*hx;
352: /* -------------------------------------------------------------
353: ---- Compute diffusion over the locally owned part of the grid
354: */
355: /* He clusters larger than 5 do not diffuse -- are immobile */
356: for (He=1; He<PetscMin(NHe+1,6); He++) {
357: f[xi].He[He] += ctx->HeDiffusion[He]*(-2.0*c[xi].He[He] + c[xi-1].He[He] + c[xi+1].He[He])*sx;
358: }
360: /* V and I clusters ONLY of size 1 diffuse */
361: f[xi].V[1] += ctx->VDiffusion[1]*(-2.0*c[xi].V[1] + c[xi-1].V[1] + c[xi+1].V[1])*sx;
362: f[xi].I[1] += ctx->IDiffusion[1]*(-2.0*c[xi].I[1] + c[xi-1].I[1] + c[xi+1].I[1])*sx;
364: /* Mixed He - V clusters are immobile */
366: /* ----------------------------------------------------------------
367: ---- Compute forcing that produces He of cluster size 1
368: Crude cubic approximation of graph from Tibo's notes
369: */
370: f[xi].He[1] += ctx->forcingScale*PetscMax(0.0,0.0006*x*x*x - 0.0087*x*x + 0.0300*x);
372: cHeVInitialize(&c[xi].He[1],cHeV);
373: cHeVInitialize(&f[xi].He[1],fHeV);
375: /* -------------------------------------------------------------------------
376: ---- Compute dissociation terms that removes an item from a cluster
377: I assume dissociation means losing only a single item from a cluster
378: I cannot tell from the notes if clusters can break up into any sub-size.
379: */
380: /* He[He] -> He[He-1] + He[1] */
381: for (He=2; He<NHe+1; He++) {
382: f[xi].He[He-1] += ctx->dissociationScale*c[xi].He[He];
383: f[xi].He[1] += ctx->dissociationScale*c[xi].He[He];
384: f[xi].He[He] -= ctx->dissociationScale*c[xi].He[He];
385: }
387: /* V[V] -> V[V-1] + V[1] */
388: for (V=2; V<NV+1; V++) {
389: f[xi].V[V-1] += ctx->dissociationScale*c[xi].V[V];
390: f[xi].V[1] += ctx->dissociationScale*c[xi].V[V];
391: f[xi].V[V] -= ctx->dissociationScale*c[xi].V[V];
392: }
394: /* I[I] -> I[I-1] + I[1] */
395: for (I=2; I<NI+1; I++) {
396: f[xi].I[I-1] += ctx->dissociationScale*c[xi].I[I];
397: f[xi].I[1] += ctx->dissociationScale*c[xi].I[I];
398: f[xi].I[I] -= ctx->dissociationScale*c[xi].I[I];
399: }
401: /* He[He]-V[1] -> He[He] + V[1] */
402: for (He=1; He<NHeV[1]+1; He++) {
403: f[xi].He[He] += 1000*ctx->dissociationScale*cHeV[1][He];
404: f[xi].V[1] += 1000*ctx->dissociationScale*cHeV[1][He];
405: fHeV[1][He] -= 1000*ctx->dissociationScale*cHeV[1][He];
406: }
408: /* He[1]-V[V] -> He[1] + V[V] */
409: for (V=2; V<MHeV+1; V++) {
410: f[xi].He[1] += 1000*ctx->dissociationScale*cHeV[V][1];
411: f[xi].V[V] += 1000*ctx->dissociationScale*cHeV[V][1];
412: fHeV[V][1] -= 1000*ctx->dissociationScale*cHeV[V][1];
413: }
415: /* He[He]-V[V] -> He[He-1]-V[V] + He[1] */
416: for (V=2; V<MHeV+1; V++) {
417: for (He=2; He<NHeV[V]+1; He++) {
418: f[xi].He[1] += 1000*ctx->dissociationScale*cHeV[V][He];
419: fHeV[V][He-1] += 1000*ctx->dissociationScale*cHeV[V][He];
420: fHeV[V][He] -= 1000*ctx->dissociationScale*cHeV[V][He];
421: }
422: }
424: /* He[He]-V[V] -> He[He]-V[V-1] + V[1] */
425: for (V=2; V<MHeV+1; V++) {
426: for (He=2; He<NHeV[V-1]+1; He++) {
427: f[xi].V[1] += 1000*ctx->dissociationScale*cHeV[V][He];
428: fHeV[V-1][He] += 1000*ctx->dissociationScale*cHeV[V][He];
429: fHeV[V][He] -= 1000*ctx->dissociationScale*cHeV[V][He];
430: }
431: }
433: /* He[He]-V[V] -> He[He]-V[V+1] + I[1] */
434: for (V=1; V<MHeV; V++) {
435: for (He=1; He<NHeV[V]+1; He++) {
436: fHeV[V+1][He] += 1000*ctx->dissociationScale*cHeV[V][He];
437: f[xi].I[1] += 1000*ctx->dissociationScale*cHeV[V][He];
438: fHeV[V][He] -= 1000*ctx->dissociationScale*cHeV[V][He];
439: }
440: }
442: /* ----------------------------------------------------------------
443: ---- Compute reaction terms that can create a cluster of given size
444: */
445: /* He[He] + He[he] -> He[He+he] */
446: for (He=2; He<NHe+1; He++) {
447: /* compute all pairs of clusters of smaller size that can combine to create a cluster of size He,
448: remove the upper half since they are symmetric to the lower half of the pairs. For example
449: when He = 5 (cluster size 5) the pairs are
450: 1 4
451: 2 2
452: 3 2 these last two are not needed in the sum since they repeat from above
453: 4 1 this is why he < (He/2) + 1 */
454: for (he=1; he<(He/2)+1; he++) {
455: f[xi].He[He] += ctx->reactionScale*c[xi].He[he]*c[xi].He[He-he];
457: /* remove the two clusters that merged to form the larger cluster */
458: f[xi].He[he] -= ctx->reactionScale*c[xi].He[he]*c[xi].He[He-he];
459: f[xi].He[He-he] -= ctx->reactionScale*c[xi].He[he]*c[xi].He[He-he];
460: }
461: }
463: /* V[V] + V[v] -> V[V+v] */
464: for (V=2; V<NV+1; V++) {
465: for (v=1; v<(V/2)+1; v++) {
466: f[xi].V[V] += ctx->reactionScale*c[xi].V[v]*c[xi].V[V-v];
467: f[xi].V[v] -= ctx->reactionScale*c[xi].V[v]*c[xi].V[V-v];
468: f[xi].V[V-v] -= ctx->reactionScale*c[xi].V[v]*c[xi].V[V-v];
469: }
470: }
472: /* I[I] + I[i] -> I[I+i] */
473: for (I=2; I<NI+1; I++) {
474: for (i=1; i<(I/2)+1; i++) {
475: f[xi].I[I] += ctx->reactionScale*c[xi].I[i]*c[xi].I[I-i];
476: f[xi].I[i] -= ctx->reactionScale*c[xi].I[i]*c[xi].I[I-i];
477: f[xi].I[I-i] -= ctx->reactionScale*c[xi].I[i]*c[xi].I[I-i];
478: }
479: }
481: /* He[1] + V[1] -> He[1]-V[1] */
482: fHeV[1][1] += 1000*ctx->reactionScale*c[xi].He[1]*c[xi].V[1];
483: f[xi].He[1] -= 1000*ctx->reactionScale*c[xi].He[1]*c[xi].V[1];
484: f[xi].V[1] -= 1000*ctx->reactionScale*c[xi].He[1]*c[xi].V[1];
486: /* He[He]-V[V] + He[he] -> He[He+he]-V[V] */
487: for (V=1; V<MHeV+1; V++) {
488: for (He=1; He<NHeV[V]; He++) {
489: for (he=1; he+He<NHeV[V]+1; he++) {
490: fHeV[V][He+he] += ctx->reactionScale*cHeV[V][He]*c[xi].He[he];
491: f[xi].He[he] -= ctx->reactionScale*cHeV[V][He]*c[xi].He[he];
492: fHeV[V][He] -= ctx->reactionScale*cHeV[V][He]*c[xi].He[he];
493: }
494: }
495: }
497: /* He[He]-V[V] + V[1] -> He[He][V+1] */
498: for (V=1; V<MHeV; V++) {
499: for (He=1; He<NHeV[V+1]; He++) {
500: fHeV[V+1][He] += ctx->reactionScale*cHeV[V][He]*c[xi].V[1];
501: /* remove the two clusters that merged to form the larger cluster */
502: f[xi].V[1] -= ctx->reactionScale*cHeV[V][He]*c[xi].V[1];
503: fHeV[V][He] -= ctx->reactionScale*cHeV[V][He]*c[xi].V[1];
504: }
505: }
507: /* He[He]-V[V] + He[he]-V[v] -> He[He+he][V+v] */
508: /* Currently the reaction rates for this are zero */
511: /* V[V] + I[I] -> V[V-I] if V > I else I[I-V] */
512: for (V=1; V<NV+1; V++) {
513: for (I=1; I<PetscMin(V,NI); I++) {
514: f[xi].V[V-I] += ctx->reactionScale*c[xi].V[V]*c[xi].I[I];
515: f[xi].V[V] -= ctx->reactionScale*c[xi].V[V]*c[xi].I[I];
516: f[xi].I[I] -= ctx->reactionScale*c[xi].V[V]*c[xi].I[I];
517: }
518: for (I=V+1; I<NI+1; I++) {
519: f[xi].I[I-V] += ctx->reactionScale*c[xi].V[V]*c[xi].I[I];
520: f[xi].V[V] -= ctx->reactionScale*c[xi].V[V]*c[xi].I[I];
521: f[xi].I[I] -= ctx->reactionScale*c[xi].V[V]*c[xi].I[I];
522: }
523: }
524: }
526: /*
527: Restore vectors
528: */
529: c = (Concentrations*)(((PetscScalar*)c)+1);
530: DMDAVecRestoreArray(da,localC,&c);
531: f = (Concentrations*)(((PetscScalar*)f)+1);
532: DMDAVecRestoreArray(da,F,&f);
533: DMRestoreLocalVector(da,&localC);
534: cHeVDestroy(cHeV);
535: cHeVDestroy(fHeV);
536: return(0);
537: }
541: /*
542: Compute the Jacobian entries based on IFuction() and insert them into the matrix
543: */
544: PetscErrorCode RHSJacobian(TS ts,PetscReal ftime,Vec C,Mat A,Mat J,void *ptr)
545: {
546: AppCtx *ctx = (AppCtx*) ptr;
547: DM da;
548: PetscErrorCode ierr;
549: PetscInt xi,Mx,xs,xm,He,he,V,v,I,i;
550: PetscInt row[3],col[3];
551: PetscReal hx,sx,x,val[6];
552: const Concentrations *c,*f;
553: Vec localC;
554: const PetscReal *rowstart,*colstart;
555: const PetscReal **cHeV,**fHeV;
556: static PetscBool initialized = PETSC_FALSE;
559: cHeVCreate((PetscScalar***)&cHeV);
560: cHeVCreate((PetscScalar***)&fHeV);
561: MatZeroEntries(J);
562: TSGetDM(ts,&da);
563: DMGetLocalVector(da,&localC);
564: 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);
565: hx = 8.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
567: DMGlobalToLocalBegin(da,C,INSERT_VALUES,localC);
568: DMGlobalToLocalEnd(da,C,INSERT_VALUES,localC);
570: /*
571: The f[] is dummy, values are never set into it. It is only used to determine the
572: local row for the entries in the Jacobian
573: */
574: DMDAVecGetArray(da,localC,&c);
575: /* Shift the c pointer to allow accessing with index of 1, instead of 0 */
576: c = (Concentrations*)(((PetscScalar*)c)-1);
577: DMDAVecGetArray(da,C,&f);
578: f = (Concentrations*)(((PetscScalar*)f)-1);
580: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
582: rowstart = &f[xs].He[1] - DOF;
583: colstart = &c[xs-1].He[1];
585: if (!initialized) {
586: /*
587: Loop over grid points computing Jacobian terms for each grid point
588: */
589: for (xi=xs; xi<xs+xm; xi++) {
590: x = xi*hx;
591:
592: cHeVInitialize(&c[xi].He[1],(PetscScalar**)cHeV);
593: cHeVInitialize(&f[xi].He[1],(PetscScalar**)fHeV);
594:
595: /* -------------------------------------------------------------
596: ---- Compute diffusion over the locally owned part of the grid
597: */
598: /* He clusters larger than 5 do not diffuse -- are immobile */
599: for (He=1; He<PetscMin(NHe+1,6); He++) {
600: row[0] = &f[xi].He[He] - rowstart;
601: col[0] = &c[xi-1].He[He] - colstart;
602: col[1] = &c[xi].He[He] - colstart;
603: col[2] = &c[xi+1].He[He] - colstart;
604: val[0] = ctx->HeDiffusion[He]*sx;
605: val[1] = -2.0*ctx->HeDiffusion[He]*sx;
606: val[2] = ctx->HeDiffusion[He]*sx;
607: MatSetValuesLocal(J,1,row,3,col,val,ADD_VALUES);
608: }
610: /* V and I clusters ONLY of size 1 diffuse */
611: row[0] = &f[xi].V[1] - rowstart;
612: col[0] = &c[xi-1].V[1] - colstart;
613: col[1] = &c[xi].V[1] - colstart;
614: col[2] = &c[xi+1].V[1] - colstart;
615: val[0] = ctx->VDiffusion[1]*sx;
616: val[1] = -2.0*ctx->VDiffusion[1]*sx;
617: val[2] = ctx->VDiffusion[1]*sx;
618: MatSetValuesLocal(J,1,row,3,col,val,ADD_VALUES);
619:
620: row[0] = &f[xi].I[1] - rowstart;
621: col[0] = &c[xi-1].I[1] - colstart;
622: col[1] = &c[xi].I[1] - colstart;
623: col[2] = &c[xi+1].I[1] - colstart;
624: val[0] = ctx->IDiffusion[1]*sx;
625: val[1] = -2.0*ctx->IDiffusion[1]*sx;
626: val[2] = ctx->IDiffusion[1]*sx;
627: MatSetValuesLocal(J,1,row,3,col,val,ADD_VALUES);
628:
629: /* Mixed He - V clusters are immobile */
630:
631: /* -------------------------------------------------------------------------
632: ---- Compute dissociation terms that removes an item from a cluster
633: I assume dissociation means losing only a single item from a cluster
634: I cannot tell from the notes if clusters can break up into any sub-size.
635: */
636:
637: /* He[He] -> He[He-1] + He[1] */
638: for (He=2; He<NHe+1; He++) {
639: row[0] = &f[xi].He[He-1] - rowstart;
640: row[1] = &f[xi].He[1] - rowstart;
641: row[2] = &f[xi].He[He] - rowstart;
642: col[0] = &c[xi].He[He] - colstart;
643: val[0] = ctx->dissociationScale;
644: val[1] = ctx->dissociationScale;
645: val[2] = -ctx->dissociationScale;
646: MatSetValuesLocal(J,3,row,1,col,val,ADD_VALUES);
647: }
648:
649: /* V[V] -> V[V-1] + V[1] */
650: for (V=2; V<NV+1; V++) {
651: row[0] = &f[xi].V[V-1] - rowstart;
652: row[1] = &f[xi].V[1] - rowstart;
653: row[2] = &f[xi].V[V] - rowstart;
654: col[0] = &c[xi].V[V] - colstart;
655: val[0] = ctx->dissociationScale;
656: val[1] = ctx->dissociationScale;
657: val[2] = -ctx->dissociationScale;
658: MatSetValuesLocal(J,3,row,1,col,val,ADD_VALUES);
659: }
660:
661: /* I[I] -> I[I-1] + I[1] */
662: for (I=2; I<NI+1; I++) {
663: row[0] = &f[xi].I[I-1] - rowstart;
664: row[1] = &f[xi].I[1] - rowstart;
665: row[2] = &f[xi].I[I] - rowstart;
666: col[0] = &c[xi].I[I] - colstart;
667: val[0] = ctx->dissociationScale;
668: val[1] = ctx->dissociationScale;
669: val[2] = -ctx->dissociationScale;
670: MatSetValuesLocal(J,3,row,1,col,val,ADD_VALUES);
671: }
672:
673: /* He[He]-V[1] -> He[He] + V[1] */
674: for (He=1; He<NHeV[1]+1; He++) {
675: row[0] = &f[xi].He[He] - rowstart;
676: row[1] = &f[xi].V[1] - rowstart;
677: row[2] = &fHeV[1][He] - rowstart;
678: col[0] = &cHeV[1][He] - colstart;
679: val[0] = 1000*ctx->dissociationScale;
680: val[1] = 1000*ctx->dissociationScale;
681: val[2] = -1000*ctx->dissociationScale;
682: MatSetValuesLocal(J,3,row,1,col,val,ADD_VALUES);
683: }
684:
685: /* He[1]-V[V] -> He[1] + V[V] */
686: for (V=2; V<MHeV+1; V++) {
687: row[0] = &f[xi].He[1] - rowstart;
688: row[1] = &f[xi].V[V] - rowstart;
689: row[2] = &fHeV[V][1] - rowstart;
690: col[0] = &cHeV[V][1] - colstart;
691: val[0] = 1000*ctx->dissociationScale;
692: val[1] = 1000*ctx->dissociationScale;
693: val[2] = -1000*ctx->dissociationScale;
694: MatSetValuesLocal(J,3,row,1,col,val,ADD_VALUES);
695: }
696:
697: /* He[He]-V[V] -> He[He-1]-V[V] + He[1] */
698: for (V=2; V<MHeV+1; V++) {
699: for (He=2; He<NHeV[V]+1; He++) {
700: row[0] = &f[xi].He[1] - rowstart;
701: row[1] = &fHeV[V][He-1] - rowstart;
702: row[2] = &fHeV[V][He] - rowstart;
703: col[0] = &cHeV[V][He] - colstart;
704: val[0] = 1000*ctx->dissociationScale;
705: val[1] = 1000*ctx->dissociationScale;
706: val[2] = -1000*ctx->dissociationScale;
707: MatSetValuesLocal(J,3,row,1,col,val,ADD_VALUES);
708: }
709: }
710:
711: /* He[He]-V[V] -> He[He]-V[V-1] + V[1] */
712: for (V=2; V<MHeV+1; V++) {
713: for (He=2; He<NHeV[V-1]+1; He++) {
714: row[0] = &f[xi].V[1] - rowstart;
715: row[1] = &fHeV[V-1][He] - rowstart;
716: row[2] = &fHeV[V][He] - rowstart;
717: col[0] = &cHeV[V][He] - colstart;
718: val[0] = 1000*ctx->dissociationScale;
719: val[1] = 1000*ctx->dissociationScale;
720: val[2] = -1000*ctx->dissociationScale;
721: MatSetValuesLocal(J,3,row,1,col,val,ADD_VALUES);
722: }
723: }
724:
725: /* He[He]-V[V] -> He[He]-V[V+1] + I[1] */
726: for (V=1; V<MHeV; V++) {
727: for (He=1; He<NHeV[V]+1; He++) {
728: row[0] = &fHeV[V+1][He] - rowstart;
729: row[1] = &f[xi].I[1] - rowstart;
730: row[2] = &fHeV[V][He] - rowstart;
731: col[0] = &cHeV[V][He] - colstart;
732: val[0] = 1000*ctx->dissociationScale;
733: val[1] = 1000*ctx->dissociationScale;
734: val[2] = -1000*ctx->dissociationScale;
735: MatSetValuesLocal(J,3,row,1,col,val,ADD_VALUES);
736: }
737: }
738: }
739: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
740: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
741: MatSetOption(J,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);
742: MatStoreValues(J);
743: MatSetFromOptions(J);
744: initialized = PETSC_TRUE;
745: } else {
746: MatRetrieveValues(J);
747: }
749: /*
750: Loop over grid points computing Jacobian terms for each grid point for reaction terms
751: */
752: for (xi=xs; xi<xs+xm; xi++) {
753: x = xi*hx;
754: cHeVInitialize(&c[xi].He[1],(PetscScalar**)cHeV);
755: cHeVInitialize(&f[xi].He[1],(PetscScalar**)fHeV);
756: /* ----------------------------------------------------------------
757: ---- Compute reaction terms that can create a cluster of given size
758: */
759: /* He[He] + He[he] -> He[He+he] */
760: for (He=2; He<NHe+1; He++) {
761: /* compute all pairs of clusters of smaller size that can combine to create a cluster of size He,
762: remove the upper half since they are symmetric to the lower half of the pairs. For example
763: when He = 5 (cluster size 5) the pairs are
764: 1 4
765: 2 2
766: 3 2 these last two are not needed in the sum since they repeat from above
767: 4 1 this is why he < (He/2) + 1 */
768: for (he=1; he<(He/2)+1; he++) {
769: row[0] = &f[xi].He[He] - rowstart;
770: row[1] = &f[xi].He[he] - rowstart;
771: row[2] = &f[xi].He[He-he] - rowstart;
772: col[0] = &c[xi].He[he] - colstart;
773: col[1] = &c[xi].He[He-he] - colstart;
774: val[0] = ctx->reactionScale*c[xi].He[He-he];
775: val[1] = ctx->reactionScale*c[xi].He[he];
776: val[2] = -ctx->reactionScale*c[xi].He[He-he];
777: val[3] = -ctx->reactionScale*c[xi].He[he];
778: val[4] = -ctx->reactionScale*c[xi].He[He-he];
779: val[5] = -ctx->reactionScale*c[xi].He[he];
780: MatSetValuesLocal(J,3,row,2,col,val,ADD_VALUES);
781: }
782: }
784: /* V[V] + V[v] -> V[V+v] */
785: for (V=2; V<NV+1; V++) {
786: for (v=1; v<(V/2)+1; v++) {
787: row[0] = &f[xi].V[V] - rowstart;
788: row[1] = &f[xi].V[v] - rowstart;
789: row[2] = &f[xi].V[V-v] - rowstart;
790: col[0] = &c[xi].V[v] - colstart;
791: col[1] = &c[xi].V[V-v] - colstart;
792: val[0] = ctx->reactionScale*c[xi].V[V-v];
793: val[1] = ctx->reactionScale*c[xi].V[v];
794: val[2] = -ctx->reactionScale*c[xi].V[V-v];
795: val[3] = -ctx->reactionScale*c[xi].V[v];
796: val[4] = -ctx->reactionScale*c[xi].V[V-v];
797: val[5] = -ctx->reactionScale*c[xi].V[v];
798: MatSetValuesLocal(J,3,row,2,col,val,ADD_VALUES);
799: }
800: }
802: /* I[I] + I[i] -> I[I+i] */
803: for (I=2; I<NI+1; I++) {
804: for (i=1; i<(I/2)+1; i++) {
805: row[0] = &f[xi].I[I] - rowstart;
806: row[1] = &f[xi].I[i] - rowstart;
807: row[2] = &f[xi].I[I-i] - rowstart;
808: col[0] = &c[xi].I[i] - colstart;
809: col[1] = &c[xi].I[I-i] - colstart;
810: val[0] = ctx->reactionScale*c[xi].I[I-i];
811: val[1] = ctx->reactionScale*c[xi].I[i];
812: val[2] = -ctx->reactionScale*c[xi].I[I-i];
813: val[3] = -ctx->reactionScale*c[xi].I[i];
814: val[4] = -ctx->reactionScale*c[xi].I[I-i];
815: val[5] = -ctx->reactionScale*c[xi].I[i];
816: MatSetValuesLocal(J,3,row,2,col,val,ADD_VALUES);
817: }
818: }
820: /* He[1] + V[1] -> He[1]-V[1] */
821: row[0] = &fHeV[1][1] - rowstart;
822: row[1] = &f[xi].He[1] - rowstart;
823: row[2] = &f[xi].V[1] - rowstart;
824: col[0] = &c[xi].He[1] - colstart;
825: col[1] = &c[xi].V[1] - colstart;
826: val[0] = 1000*ctx->reactionScale*c[xi].V[1];
827: val[1] = 1000*ctx->reactionScale*c[xi].He[1];
828: val[2] = -1000*ctx->reactionScale*c[xi].V[1];
829: val[3] = -1000*ctx->reactionScale*c[xi].He[1];
830: val[4] = -1000*ctx->reactionScale*c[xi].V[1];
831: val[5] = -1000*ctx->reactionScale*c[xi].He[1];
832: MatSetValuesLocal(J,3,row,2,col,val,ADD_VALUES);
834: /* He[He]-V[V] + He[he] -> He[He+he]-V[V] */
835: for (V=1; V<MHeV+1; V++) {
836: for (He=1; He<NHeV[V]; He++) {
837: for (he=1; he+He<NHeV[V]+1; he++) {
838: row[0] = &fHeV[V][He+he] - rowstart;
839: row[1] = &f[xi].He[he] - rowstart;
840: row[2] = &fHeV[V][He] - rowstart;
841: col[0] = &c[xi].He[he] - colstart;
842: col[1] = &cHeV[V][He] - colstart;
843: val[0] = ctx->reactionScale*cHeV[V][He];
844: val[1] = ctx->reactionScale*c[xi].He[he];
845: val[2] = -ctx->reactionScale*cHeV[V][He];
846: val[3] = -ctx->reactionScale*c[xi].He[he];
847: val[4] = -ctx->reactionScale*cHeV[V][He];
848: val[5] = -ctx->reactionScale*c[xi].He[he];
849: MatSetValuesLocal(J,3,row,2,col,val,ADD_VALUES);
850: }
851: }
852: }
854: /* He[He]-V[V] + V[1] -> He[He][V+1] */
855: for (V=1; V<MHeV; V++) {
856: for (He=1; He<NHeV[V+1]; He++) {
857: row[0] = &fHeV[V+1][He] - rowstart;
858: row[1] = &f[xi].V[1] - rowstart;
859: row[2] = &fHeV[V][He] - rowstart;
860: col[0] = &c[xi].V[1] - colstart;
861: col[1] = &cHeV[V][He] - colstart;
862: val[0] = ctx->reactionScale*cHeV[V][He];
863: val[1] = ctx->reactionScale*c[xi].V[1];
864: val[2] = -ctx->reactionScale*cHeV[V][He];
865: val[3] = -ctx->reactionScale*c[xi].V[1];
866: val[4] = -ctx->reactionScale*cHeV[V][He];
867: val[5] = -ctx->reactionScale*c[xi].V[1];
868: MatSetValuesLocal(J,3,row,2,col,val,ADD_VALUES);
869: }
870: }
872: /* He[He]-V[V] + He[he]-V[v] -> He[He+he][V+v] */
873: /* Currently the reaction rates for this are zero */
876: /* V[V] + I[I] -> V[V-I] if V > I else I[I-V] */
877: for (V=1; V<NV+1; V++) {
878: for (I=1; I<PetscMin(V,NI); I++) {
879: row[0] = &f[xi].V[V-I] - rowstart;
880: row[1] = &f[xi].V[V] - rowstart;
881: row[2] = &f[xi].I[I] - rowstart;
882: col[0] = &c[xi].V[V] - colstart;
883: col[1] = &c[xi].I[I] - colstart;
884: val[0] = ctx->reactionScale*c[xi].I[I];
885: val[1] = ctx->reactionScale*c[xi].V[V];
886: val[2] = -ctx->reactionScale*c[xi].I[I];
887: val[3] = -ctx->reactionScale*c[xi].V[V];
888: val[4] = -ctx->reactionScale*c[xi].I[I];
889: val[5] = -ctx->reactionScale*c[xi].V[V];
890: MatSetValuesLocal(J,3,row,2,col,val,ADD_VALUES);
891: }
892: for (I=V+1; I<NI+1; I++) {
893: row[0] = &f[xi].I[I-V] - rowstart;
894: row[1] = &f[xi].V[V] - rowstart;
895: row[2] = &f[xi].I[I] - rowstart;
896: col[0] = &c[xi].V[V] - colstart;
897: col[1] = &c[xi].I[I] - colstart;
898: val[0] = ctx->reactionScale*c[xi].I[I];
899: val[1] = ctx->reactionScale*c[xi].V[V];
900: val[2] = -ctx->reactionScale*c[xi].I[I];
901: val[3] = -ctx->reactionScale*c[xi].V[V];
902: val[4] = -ctx->reactionScale*c[xi].I[I];
903: val[5] = -ctx->reactionScale*c[xi].V[V];
904: MatSetValuesLocal(J,3,row,2,col,val,ADD_VALUES);
905: }
906: }
907: }
909: /*
910: Restore vectors
911: */
912: c = (Concentrations*)(((PetscScalar*)c)+1);
913: DMDAVecRestoreArray(da,localC,&c);
914: f = (Concentrations*)(((PetscScalar*)f)+1);
915: DMDAVecRestoreArray(da,C,&f);
916: DMRestoreLocalVector(da,&localC);
917: cHeVDestroy((PetscScalar**)cHeV);
918: cHeVDestroy((PetscScalar**)fHeV);
920: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
921: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
922: if (A != J) {
923: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
924: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
925: }
926: return(0);
927: }
931: /*
932: Determines the nonzero structure within the diagonal blocks of the Jacobian that represent coupling resulting from reactions and
933: dissasociations of the clusters
934: */
935: PetscErrorCode GetDfill(PetscInt *dfill, void *ptr)
936: {
937: PetscInt He,he,V,v,I,i,j,k,rows[3],cols[2];
938: Concentrations *c;
939: PetscScalar *idxstart,**cHeV;
942: /* ensure fill for the diagonal of matrix */
943: for (i=0; i<(DOF); i++) {
944: dfill[i*DOF + i] = 1;
945: }
947: /*
948: c is never used except for computing offsets between variables which are used to fill the non-zero
949: structure of the matrix
950: */
951: PetscMalloc(sizeof(Concentrations),&c);
952: c = (Concentrations*)(((PetscScalar*)c)-1);
953: cHeVCreate(&cHeV);
954: cHeVInitialize(&c->He[1],cHeV);
955: idxstart = (PetscScalar*)&c->He[1];
957: /* -------------------------------------------------------------------------
958: ---- Compute dissociation terms that removes an item from a cluster
959: I assume dissociation means losing only a single item from a cluster
960: I cannot tell from the notes if clusters can break up into any sub-size.
961: */
962: /* He[He] -> He[He-1] + He[1] */
963: for (He=2; He<NHe+1; He++) {
964: rows[0] = &c->He[He-1] - idxstart;
965: rows[1] = &c->He[1] - idxstart;
966: rows[2] = &c->He[He] - idxstart;
967: cols[0] = &c->He[He] - idxstart;
968: for (j=0; j<3; j++) {
969: dfill[rows[j]*DOF + cols[0]] = 1;
970: }
971: }
973: /* V[V] -> V[V-1] + V[1] */
974: for (V=2; V<NV+1; V++) {
975: rows[0] = &c->V[V] - idxstart;
976: rows[1] = &c->V[1] - idxstart;
977: rows[2] = &c->V[V-1] - idxstart;
978: cols[0] = &c->V[V] - idxstart;
979: for (j=0; j<3; j++) {
980: dfill[rows[j]*DOF + cols[0]] = 1;
981: }
982: }
983:
984: /* I[I] -> I[I-1] + I[1] */
985: for (I=2; I<NI+1; I++) {
986: rows[0] = &c->I[I] - idxstart;
987: rows[1] = &c->I[1] - idxstart;
988: rows[2] = &c->I[I-1] - idxstart;
989: cols[0] = &c->I[I] - idxstart;
990: for (j=0; j<3; j++) {
991: dfill[rows[j]*DOF + cols[0]] = 1;
992: }
993: }
994:
995: /* He[He]-V[1] -> He[He] + V[1] */
996: for (He=1; He<NHeV[1]+1; He++) {
997: rows[0] = &c->He[He] - idxstart;
998: rows[1] = &c->V[1] - idxstart;
999: rows[2] = &cHeV[1][He] - idxstart;
1000: cols[0] = &cHeV[1][He] - idxstart;
1001: for (j=0; j<3; j++) {
1002: dfill[rows[j]*DOF + cols[0]] = 1;
1003: }
1004: }
1005:
1006: /* He[1]-V[V] -> He[1] + V[V] */
1007: for (V=2; V<MHeV+1; V++) {
1008: rows[0] = &c->He[1] - idxstart;
1009: rows[1] = &c->V[V] - idxstart;
1010: rows[2] = &cHeV[V][1] - idxstart;
1011: cols[0] = &cHeV[V][1] - idxstart;
1012: for (j=0; j<3; j++) {
1013: dfill[rows[j]*DOF + cols[0]] = 1;
1014: }
1015: }
1016:
1017: /* He[He]-V[V] -> He[He-1]-V[V] + He[1] */
1018: for (V=2; V<MHeV+1; V++) {
1019: for (He=2; He<NHeV[V]+1; He++) {
1020: rows[0] = &c->He[1] - idxstart;
1021: rows[1] = &cHeV[V][He] - idxstart;
1022: rows[2] = &cHeV[V][He-1] - idxstart;
1023: cols[0] = &cHeV[V][He] - idxstart;
1024: for (j=0; j<3; j++) {
1025: dfill[rows[j]*DOF + cols[0]] = 1;
1026: }
1027: }
1028: }
1029:
1030: /* He[He]-V[V] -> He[He]-V[V-1] + V[1] */
1031: for (V=2; V<MHeV+1; V++) {
1032: for (He=2; He<NHeV[V-1]+1; He++) {
1033: rows[0] = &c->V[1] - idxstart;
1034: rows[1] = &cHeV[V][He] - idxstart;
1035: rows[2] = &cHeV[V-1][He] - idxstart;
1036: cols[0] = &cHeV[V][He] - idxstart;
1037: for (j=0; j<3; j++) {
1038: dfill[rows[j]*DOF + cols[0]] = 1;
1039: }
1040: }
1041: }
1042:
1043: /* He[He]-V[V] -> He[He]-V[V+1] + I[1] */
1044: for (V=1; V<MHeV; V++) {
1045: for (He=1; He<NHeV[V]+1; He++) {
1046: rows[0] = &c->I[1] - idxstart;
1047: rows[1] = &cHeV[V+1][He] - idxstart;
1048: rows[2] = &cHeV[V][He] - idxstart;
1049: cols[0] = &cHeV[V][He] - idxstart;
1050: for (j=0; j<3; j++) {
1051: dfill[rows[j]*DOF + cols[0]] = 1;
1052: }
1053: }
1054: }
1056: /* These are the reaction terms in the diagonal block */
1057: for (He=2; He<NHe+1; He++) {
1058: for (he=1; he<(He/2)+1; he++) {
1059: rows[0] = &c->He[He] - idxstart;
1060: rows[1] = &c->He[he] - idxstart;
1061: rows[2] = &c->He[He-he] - idxstart;
1062: cols[0] = &c->He[he] - idxstart;
1063: cols[1] = &c->He[He-he] - idxstart;
1064: for (j=0; j<3; j++) {
1065: for (k=0; k<2; k++) {
1066: dfill[rows[j]*DOF + cols[k]] = 1;
1067: }
1068: }
1069: }
1070: }
1072: /* V[V] + V[v] -> V[V+v] */
1073: for (V=2; V<NV+1; V++) {
1074: for (v=1; v<(V/2)+1; v++) {
1075: rows[0] = &c->V[V] - idxstart;
1076: rows[1] = &c->V[v] - idxstart;
1077: rows[2] = &c->V[V-v] - idxstart;
1078: cols[0] = &c->V[v] - idxstart;
1079: cols[1] = &c->V[V-v] - idxstart;
1080: for (j=0; j<3; j++) {
1081: for (k=0; k<2; k++) {
1082: dfill[rows[j]*DOF + cols[k]] = 1;
1083: }
1084: }
1085: }
1086: }
1087:
1088: /* I[I] + I[i] -> I[I+i] */
1089: for (I=2; I<NI+1; I++) {
1090: for (i=1; i<(I/2)+1; i++) {
1091: rows[0] = &c->I[I] - idxstart;
1092: rows[1] = &c->I[i] - idxstart;
1093: rows[2] = &c->I[I-i] - idxstart;
1094: cols[0] = &c->I[i] - idxstart;
1095: cols[1] = &c->I[I-i] - idxstart;
1096: for (j=0; j<3; j++) {
1097: for (k=0; k<2; k++) {
1098: dfill[rows[j]*DOF + cols[k]] = 1;
1099: }
1100: }
1101: }
1102: }
1103:
1104: /* He[1] + V[1] -> He[1]-V[1] */
1105: rows[0] = &cHeV[1][1] - idxstart;
1106: rows[1] = &c->He[1] - idxstart;
1107: rows[2] = &c->V[1] - idxstart;
1108: cols[0] = &c->He[1] - idxstart;
1109: cols[1] = &c->V[1] - idxstart;
1110: for (j=0; j<3; j++) {
1111: for (k=0; k<2; k++) {
1112: dfill[rows[j]*DOF + cols[k]] = 1;
1113: }
1114: }
1115:
1116: /* He[He]-V[V] + He[he] -> He[He+he]-V[V] */
1117: for (V=1; V<MHeV+1; V++) {
1118: for (He=1; He<NHeV[V]; He++) {
1119: for (he=1; he+He<NHeV[V]+1; he++) {
1120: rows[0] = &cHeV[V][He+he] - idxstart;
1121: rows[1] = &c->He[he] - idxstart;
1122: rows[2] = &cHeV[V][He] - idxstart;
1123: cols[0] = &cHeV[V][He] - idxstart;
1124: cols[1] = &c->He[he] - idxstart;
1125: for (j=0; j<3; j++) {
1126: for (k=0; k<2; k++) {
1127: dfill[rows[j]*DOF + cols[k]] = 1;
1128: }
1129: }
1130: }
1131: }
1132: }
1133: /* He[He]-V[V] + V[1] -> He[He][V+1] */
1134: for (V=1; V<MHeV; V++) {
1135: for (He=1; He<NHeV[V+1]; He++) {
1136: rows[0] = &cHeV[V+1][He] - idxstart;
1137: rows[1] = &c->V[1] - idxstart;
1138: rows[2] = &cHeV[V][He] - idxstart;
1139: cols[0] = &cHeV[V][He] - idxstart;
1140: cols[1] = &c->V[1] - idxstart;
1141: for (j=0; j<3; j++) {
1142: for (k=0; k<2; k++) {
1143: dfill[rows[j]*DOF + cols[k]] = 1;
1144: }
1145: }
1146: }
1147: }
1149: /* He[He]-V[V] + He[he]-V[v] -> He[He+he][V+v] */
1150: /* Currently the reaction rates for this are zero */
1151:
1152: /* V[V] + I[I] -> V[V-I] if V > I else I[I-V] */
1153: for (V=1; V<NV+1; V++) {
1154: for (I=1; I<PetscMin(V,NI); I++) {
1155: rows[0] = &c->V[V-I] - idxstart;
1156: rows[1] = &c->V[V] - idxstart;
1157: rows[2] = &c->I[I] - idxstart;
1158: cols[0] = &c->V[V] - idxstart;
1159: cols[1] = &c->I[I] - idxstart;
1160: for (j=0; j<3; j++) {
1161: for (k=0; k<2; k++) {
1162: dfill[rows[j]*DOF + cols[k]] = 1;
1163: }
1164: }
1165: }
1166: for (I=V+1; I<NI+1; I++) {
1167: rows[0] = &c->I[I-V] - idxstart;
1168: rows[1] = &c->V[V] - idxstart;
1169: rows[2] = &c->I[I] - idxstart;
1170: cols[0] = &c->V[V] - idxstart;
1171: cols[1] = &c->I[I] - idxstart;
1172: for (j=0; j<3; j++) {
1173: for (k=0; k<2; k++) {
1174: dfill[rows[j]*DOF + cols[k]] = 1;
1175: }
1176: }
1177: }
1178: }
1180: c = (Concentrations*)(((PetscScalar*)c)+1);
1181: cHeVDestroy(cHeV);
1182: PetscFree(c);
1183: return(0);
1184: }
1185: /* ------------------------------------------------------------------- */
1187: typedef struct {
1188: DM Heda,Vda,HeVda; /* defines the 2d layout of the He subvector */
1189: Vec He,V,HeV;
1190: VecScatter Hescatter,Vscatter,HeVscatter;
1191: PetscViewer Heviewer,Vviewer,HeVviewer;
1192: } MyMonitorCtx;
1196: /*
1197: Display He as a function of space and cluster size for each time step
1198: */
1199: PetscErrorCode MyMonitorMonitor(TS ts,PetscInt timestep,PetscReal time,Vec solution, void *ictx)
1200: {
1201: MyMonitorCtx *ctx = (MyMonitorCtx*)ictx;
1205: VecScatterBegin(ctx->Hescatter,solution,ctx->He,INSERT_VALUES,SCATTER_FORWARD);
1206: VecScatterEnd(ctx->Hescatter,solution,ctx->He,INSERT_VALUES,SCATTER_FORWARD);
1207: VecView(ctx->He,ctx->Heviewer);
1209: VecScatterBegin(ctx->Vscatter,solution,ctx->V,INSERT_VALUES,SCATTER_FORWARD);
1210: VecScatterEnd(ctx->Vscatter,solution,ctx->V,INSERT_VALUES,SCATTER_FORWARD);
1211: VecView(ctx->V,ctx->Vviewer);
1213: VecScatterBegin(ctx->HeVscatter,solution,ctx->HeV,INSERT_VALUES,SCATTER_FORWARD);
1214: VecScatterEnd(ctx->HeVscatter,solution,ctx->HeV,INSERT_VALUES,SCATTER_FORWARD);
1215: VecView(ctx->HeV,ctx->HeVviewer);
1216: return(0);
1217: }
1221: /*
1222: Frees all data structures associated with the monitor
1223: */
1224: PetscErrorCode MyMonitorDestroy(void **ictx)
1225: {
1226: MyMonitorCtx **ctx = (MyMonitorCtx**)ictx;
1230: VecScatterDestroy(&(*ctx)->Hescatter);
1231: VecDestroy(&(*ctx)->He);
1232: DMDestroy(&(*ctx)->Heda);
1233: PetscViewerDestroy(&(*ctx)->Heviewer);
1235: VecScatterDestroy(&(*ctx)->Vscatter);
1236: VecDestroy(&(*ctx)->V);
1237: DMDestroy(&(*ctx)->Vda);
1238: PetscViewerDestroy(&(*ctx)->Vviewer);
1240: VecScatterDestroy(&(*ctx)->HeVscatter);
1241: VecDestroy(&(*ctx)->HeV);
1242: DMDestroy(&(*ctx)->HeVda);
1243: PetscViewerDestroy(&(*ctx)->HeVviewer);
1244: PetscFree(*ctx);
1245: return(0);
1246: }
1250: /*
1251: Sets up a monitor that will display He as a function of space and cluster size for each time step
1252: */
1253: PetscErrorCode MyMonitorSetUp(TS ts)
1254: {
1255: DM da;
1257: PetscInt xi,xs,xm,*idx,M,xj,cnt = 0;
1258: const PetscInt *lx;
1259: Vec C;
1260: MyMonitorCtx *ctx;
1261: PetscBool flg;
1262: IS is;
1263: char ycoor[32];
1264: PetscReal valuebounds[4] = {0, 1.2, 0, 1.2};
1267: PetscOptionsHasName(NULL,"-mymonitor",&flg);
1268: if (!flg) return(0);
1270: TSGetDM(ts,&da);
1271: PetscNew(MyMonitorCtx,&ctx);
1272: PetscNew(&ctx);
1273: PetscViewerDrawOpen(PetscObjectComm((PetscObject)da),NULL,"",PETSC_DECIDE,PETSC_DECIDE,600,400,&ctx->viewer);
1275: /* setup visualization for He */
1276: PetscViewerDrawOpen(PetscObjectComm((PetscObject)da),NULL,"",PETSC_DECIDE,PETSC_DECIDE,600,400,&ctx->Heviewer);
1277: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
1278: 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);
1279: DMDAGetOwnershipRanges(da,&lx,NULL,NULL);
1280: DMDACreate2d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,M,NHe,PETSC_DETERMINE,1,1,1,lx,NULL,&ctx->Heda);
1281: DMDASetFieldName(ctx->Heda,0,"He");
1282: DMDASetCoordinateName(ctx->Heda,0,"X coordinate direction");
1283: PetscSNPrintf(ycoor,32,"%D ... Cluster size ... 1",NHe);
1284: DMDASetCoordinateName(ctx->Heda,1,ycoor);
1285: DMCreateGlobalVector(ctx->Heda,&ctx->He);
1286: PetscMalloc(NHe*xm*sizeof(PetscInt),&idx);
1287: cnt = 0;
1288: for (xj=0; xj<NHe; xj++) {
1289: for (xi=xs; xi<xs+xm; xi++) {
1290: idx[cnt++] = DOF*xi + xj;
1291: }
1292: }
1293: ISCreateGeneral(PetscObjectComm((PetscObject)ts),NHe*xm,idx,PETSC_OWN_POINTER,&is);
1294: TSGetSolution(ts,&C);
1295: VecScatterCreate(C,is,ctx->He,NULL,&ctx->Hescatter);
1296: ISDestroy(&is);
1297: /* sets the bounds on the contour plot values so the colors mean the same thing for different timesteps */
1298: PetscViewerDrawSetBounds(ctx->Heviewer,2,valuebounds);
1300: /* setup visualization for V */
1301: PetscViewerDrawOpen(PetscObjectComm((PetscObject)da),NULL,"",PETSC_DECIDE,PETSC_DECIDE,600,400,&ctx->Vviewer);
1302: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
1303: 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);
1304: DMDAGetOwnershipRanges(da,&lx,NULL,NULL);
1305: DMDACreate2d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,M,NV,PETSC_DETERMINE,1,1,1,lx,NULL,&ctx->Vda);
1306: DMDASetFieldName(ctx->Vda,0,"V");
1307: DMDASetCoordinateName(ctx->Vda,0,"X coordinate direction");
1308: PetscSNPrintf(ycoor,32,"%D ... Cluster size ... 1",NV);
1309: DMDASetCoordinateName(ctx->Vda,1,ycoor);
1310: DMCreateGlobalVector(ctx->Vda,&ctx->V);
1311: PetscMalloc(NV*xm*sizeof(PetscInt),&idx);
1312: DMCreateGlobalVector(ctx->da,&ctx->He);
1313: PetscMalloc1(2*N*xm,&idx);
1314: cnt = 0;
1315: for (xj=0; xj<NV; xj++) {
1316: for (xi=xs; xi<xs+xm; xi++) {
1317: idx[cnt++] = NHe + DOF*xi + xj;
1318: }
1319: }
1320: ISCreateGeneral(PetscObjectComm((PetscObject)ts),NV*xm,idx,PETSC_OWN_POINTER,&is);
1321: TSGetSolution(ts,&C);
1322: VecScatterCreate(C,is,ctx->V,NULL,&ctx->Vscatter);
1323: ISDestroy(&is);
1324: /* sets the bounds on the contour plot values so the colors mean the same thing for different timesteps */
1325: PetscViewerDrawSetBounds(ctx->Vviewer,2,valuebounds);
1327: /* setup visualization for HeV[1][] */
1328: PetscViewerDrawOpen(PetscObjectComm((PetscObject)da),NULL,"",PETSC_DECIDE,PETSC_DECIDE,600,400,&ctx->HeVviewer);
1329: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
1330: 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);
1331: DMDAGetOwnershipRanges(da,&lx,NULL,NULL);
1332: 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);
1333: DMDASetFieldName(ctx->HeVda,0,"HeV[1][]");
1334: DMDASetCoordinateName(ctx->HeVda,0,"X coordinate direction");
1335: PetscSNPrintf(ycoor,32,"%D ... Cluster size ... 1",NHeV[1]);
1336: DMDASetCoordinateName(ctx->HeVda,1,ycoor);
1337: DMCreateGlobalVector(ctx->HeVda,&ctx->HeV);
1338: PetscMalloc(NHeV[1]*xm*sizeof(PetscInt),&idx);
1339: cnt = 0;
1340: for (xj=0; xj<NHeV[1]; xj++) {
1341: for (xi=xs; xi<xs+xm; xi++) {
1342: idx[cnt++] = NHe + NV + NI + DOF*xi + xj;
1343: }
1344: }
1345: ISCreateGeneral(PetscObjectComm((PetscObject)ts),NHeV[1]*xm,idx,PETSC_OWN_POINTER,&is);
1346: TSGetSolution(ts,&C);
1347: VecScatterCreate(C,is,ctx->HeV,NULL,&ctx->HeVscatter);
1348: ISDestroy(&is);
1349: /* sets the bounds on the contour plot values so the colors mean the same thing for different timesteps */
1350: PetscViewerDrawSetBounds(ctx->HeVviewer,2,valuebounds);
1352: TSMonitorSet(ts,MyMonitorMonitor,ctx,MyMonitorDestroy);
1353: return(0);
1354: }
1358: PetscErrorCode MyLoadData(MPI_Comm comm,const char *filename)
1359: {
1361: FILE *fp;
1362: char buff[256];
1363: PetscInt He,V,I,lc = 0;
1364: char Hebindstr[32],Vbindstr[32],Ibindstr[32],trapbindstr[32],*sharp;
1365: PetscReal Hebind,Vbind,Ibind,trapbind;
1368: PetscFOpen(comm,filename,"r",&fp);
1369: PetscSynchronizedFGets(comm,fp,256,buff);
1370: while (buff[0]) {
1371: PetscStrchr(buff,'#',&sharp);
1372: if (!sharp) {
1373: sscanf(buff,"%d %d %d %s %s %s %s",&He,&V,&I,Hebindstr,Vbindstr,Ibindstr,trapbindstr);
1374: Hebind = strtod(Hebindstr,NULL);
1375: Vbind = strtod(Vbindstr,NULL);
1376: Ibind = strtod(Ibindstr,NULL);
1377: trapbind = strtod(trapbindstr,NULL);
1378: if (V <= NV) {
1379: if (He > NHe && V == 0 && I == 0) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Recompile with correct NHe %d %d",He,NHe);
1380: if (He == 0 && V > NV && I == 0) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Recompile with correct V %d %d",V,NV);
1381: if (He == 0 && V == 0 && I > NI) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Recompile with correct NI %d %d",I,NI);
1382: if (lc++ > DOF) SETERRQ4(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Recompile with correct NHe %d NV %d NI %d MNHeV %",NHe,NV,NI,MNHeV);
1383: if (He > 0 && V > 0) { /* assumes the He are sorted in increasing order */
1384: NHeV[V] = He;
1385: }
1386: }
1387: }
1388: PetscSynchronizedFGets(comm,fp,256,buff);
1389: }
1390: 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);
1391: return(0);
1392: }