Actual source code: ex31.c
2: static char help[] = "Model multi-physics solver\n\n";
4: /*
5: A model "multi-physics" solver based on the Vincent Mousseau's reactor core pilot code.
7: There a three grids:
9: ---------------------- ---------------
10: | | / / /
11: | | / / /
12: | | / / /
13: | | / / /
14: | | / / /
15: ---------------------- ---------------
17: A 1d grid along the left edge, a 2d grid in the middle and another 2d grid on the right.
19: */
21: #include petscdmmg.h
29: int main(int argc,char **argv)
30: {
31: DMMG *dmmg; /* multilevel grid structure */
33: MPI_Comm comm;
34: DA da;
35: VecPack pack;
37: PetscInitialize(&argc,&argv,(char *)0,help);
38: comm = PETSC_COMM_WORLD;
41: PreLoadBegin(PETSC_TRUE,"SetUp");
43: /*
44: Create the VecPack object to manage the three grids/physics.
45: We only support a 1d decomposition along the y direction (since one of the grids is 1d).
47: */
48: VecPackCreate(comm,&pack);
49: DACreate1d(comm,DA_NONPERIODIC,-6,1,1,0,&da);
50: VecPackAddDA(pack,da);
51: DADestroy(da);
52: DACreate2d(comm,DA_NONPERIODIC,DA_STENCIL_STAR,-6,-6,1,PETSC_DETERMINE,1,1,0,0,&da);
53: VecPackAddDA(pack,da);
54: DADestroy(da);
55: DACreate2d(comm,DA_NONPERIODIC,DA_STENCIL_STAR,-6,-6,1,PETSC_DETERMINE,1,1,0,0,&da);
56: VecPackAddDA(pack,da);
57: DADestroy(da);
58:
59: /*
60: Create the solver object and attach the grid/physics info
61: */
62: DMMGCreate(comm,2,0,&dmmg);
63: DMMGSetDM(dmmg,(DM)pack);
64: VecPackDestroy(pack);
65: CHKMEMQ;
68: DMMGSetInitialGuessLocal(dmmg,(PetscErrorCode (*)(void)) FormInitialGuessLocal);
69: CHKMEMQ;
70: DMMGSetSNESLocal(dmmg,(PetscErrorCode (*)(void))0,0,0,0);
71: CHKMEMQ;
73: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
74: Solve the nonlinear system
75: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
77: PreLoadStage("Solve");
78: DMMGSolve(dmmg);
80: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
81: Free work space. All PETSc objects should be destroyed when they
82: are no longer needed.
83: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
85: DMMGDestroy(dmmg);
86: PreLoadEnd();
88: PetscFinalize();
89: return 0;
90: }
92: /* ------------------------------------------------------------------- */
97: /*
98: FormInitialGuessLocal - Forms initial approximation.
100: Input Parameters:
101: user - user-defined application context
102: X - vector
104: Output Parameter:
105: X - vector
106: */
107: PetscErrorCode FormInitialGuessLocal(DALocalInfo *info1,PetscScalar *x1,DALocalInfo *info2,PetscScalar **x2,DALocalInfo *info3,PetscScalar **x3)
108: {
109: PetscInt i;
110: PetscReal hx,dhx;
113: dhx = (PetscReal)(info1->mx-1);
114: hx = 1.0/dhx;
116: for (i=info1->xs; i<info1->xs+info1->xm; i++) {
117: x1[i] = 22.3;
118: }
119: return(0);
120: }
121: