Actual source code: meshmgsnes.c
1: #include <petscmg.h> /*I "petscmg.h" I*/
2: #include <petscdmmg.h> /*I "petscdmmg.h" I*/
3: #ifdef PETSC_HAVE_SIEVE
4: #include <petscmesh.h> /*I "petscmesh.h" I*/
5: #include <Selection.hh>
6: #endif
8: /* Just to set iterations */
9: #include include/private/snesimpl.h
11: #ifdef PETSC_HAVE_SIEVE
12: PetscErrorCode DMMGFormFunctionMesh(SNES snes, Vec X, Vec F, void *ptr);
14: #if 0
15: PetscErrorCode CreateNullSpace(DMMG dmmg, Vec *nulls) {
16: Mesh mesh = (Mesh) dmmg->dm;
17: Vec nS = nulls[0];
18: SectionReal nullSpace;
22: MeshGetSectionReal(mesh, "nullSpace", &nullSpace);
23: {
24: ALE::Obj<ALE::Mesh> m;
25: ALE::Obj<ALE::Mesh::real_section_type> s;
27: MeshGetMesh(mesh, m);
28: SectionRealGetSection(nullSpace, s);
29: ALE::Obj<ALE::Discretization> disc = m->getDiscretization("p");
30: const int dim = m->getDimension();
32: for(int d = 0; d <= dim; ++d) {
33: const int numDof = disc->getNumDof(d);
35: if (numDof) {
36: const ALE::Obj<ALE::Mesh::label_sequence>& stratum = m->depthStratum(d);
37: const ALE::Mesh::label_sequence::iterator end = stratum->end();
38: double *values = new double[numDof];
40: for(ALE::Mesh::label_sequence::iterator p_iter = stratum->begin(); p_iter != end; ++p_iter) {
41: for(int i = 0; i < numDof; ++i) values[i] = 1.0;
42: s->updatePoint(*p_iter, values);
43: }
44: }
45: }
46: }
47: SectionRealToVec(nullSpace, mesh, SCATTER_FORWARD, nS);
48: std::cout << "Null space:" << std::endl;
49: VecView(nS, PETSC_VIEWER_STDOUT_SELF);
50: SectionRealDestroy(nullSpace);
51: return(0);
52: }
53: #endif
55: /* Nonlinear relaxation on all the equations with an initial guess in x */
59: PetscErrorCode Relax_Mesh(DMMG *dmmg, Mesh mesh, MatSORType flag, int its, Vec X, Vec B)
60: {
61: SectionReal sectionX, sectionB, cellX;
62: Mesh smallMesh;
63: DMMG *smallDmmg;
64: DALocalFunction1 func;
65: DALocalFunction1 jac;
66: ALE::Obj<ALE::Mesh> m;
67: ALE::Obj<ALE::Mesh::real_section_type> sX;
68: ALE::Obj<ALE::Mesh::real_section_type> sB;
69: PetscTruth fasDebug;
70: PetscErrorCode ierr;
73: PetscOptionsHasName(dmmg[0]->prefix, "-dmmg_fas_debug", &fasDebug);
74: if (fasDebug) {PetscPrintf(dmmg[0]->comm, " FAS mesh relaxation\n");}
75: if (its <= 0) SETERRQ1(PETSC_ERR_ARG_WRONG, "Relaxation requires global its %D positive", its);
76: MeshCreate(PETSC_COMM_SELF, &smallMesh);
77: DMMGCreate(PETSC_COMM_SELF, -1, PETSC_NULL, &smallDmmg);
78: //DMMGSetMatType(smallDmmg, MATSEQDENSE);
79: DMMGSetPrefix(smallDmmg, "fas_");
80: DMMGSetUser(smallDmmg, 0, DMMGGetUser(dmmg, 0));
81: DMMGGetSNESLocal(dmmg, &func, &jac);
82: MeshGetMesh(mesh, m);
83: MeshGetSectionReal(mesh, "default", §ionX);
84: SectionRealToVec(sectionX, mesh, SCATTER_REVERSE, X);
85: SectionRealGetSection(sectionX, sX);
86: MeshGetSectionReal(mesh, "constant", §ionB);
87: SectionRealToVec(sectionB, mesh, SCATTER_REVERSE, B);
88: SectionRealGetSection(sectionB, sB);
89: SectionRealCreate(PETSC_COMM_SELF, &cellX);
90: const ALE::Obj<ALE::Mesh::sieve_type>& sieve = m->getSieve();
91: const ALE::Obj<ALE::Mesh::label_sequence>& cells = m->heightStratum(0);
92: const int depth = m->depth();
93: const ALE::Obj<ALE::Mesh::label_type>& marker = m->getLabel("marker");
94: //const int cellDof = m->sizeWithBC(sX, *cells->begin());
96: ALE::Obj<ALE::Mesh::names_type> fields = m->getDiscretizations();
97: std::map<std::string, ALE::Obj<ALE::Discretization> > sDiscs;
99: for(ALE::Mesh::names_type::iterator f_iter = fields->begin(); f_iter != fields->end(); ++f_iter) {
100: const ALE::Obj<ALE::Discretization>& disc = m->getDiscretization(*f_iter);
101: ALE::Obj<ALE::Discretization> sDisc = new ALE::Discretization(disc->comm(), disc->debug());
103: sDisc->setQuadratureSize(disc->getQuadratureSize());
104: sDisc->setQuadraturePoints(disc->getQuadraturePoints());
105: sDisc->setQuadratureWeights(disc->getQuadratureWeights());
106: sDisc->setBasisSize(disc->getBasisSize());
107: sDisc->setBasis(disc->getBasis());
108: sDisc->setBasisDerivatives(disc->getBasisDerivatives());
109: for(int d = 0; d <= m->getDimension(); ++d) {
110: sDisc->setNumDof(d, disc->getNumDof(d));
111: sDisc->setDofClass(d, disc->getDofClass(d));
112: }
113: if (disc->getBoundaryConditions()->size()) {
114: if (fasDebug) {std::cout << "Adding BC for field " << *f_iter << std::endl;}
115: ALE::Obj<ALE::BoundaryCondition> sBC = new ALE::BoundaryCondition(disc->comm(), disc->debug());
116: sBC->setLabelName("marker");
117: sBC->setMarker(1);
118: sBC->setFunction(PETSC_NULL);
119: sBC->setDualIntegrator(PETSC_NULL);
120: sDisc->setBoundaryCondition(sBC);
121: }
122: sDiscs[*f_iter] = sDisc;
123: }
124: while(its--) {
125: if (fasDebug) {PetscPrintf(dmmg[0]->comm, " forward sweep %d\n", its);}
126: if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
127: // Loop over all cells
128: // This is an overlapping block SOR, but it is easier and seems more natural than doing each unknown
129: for(ALE::Mesh::label_sequence::iterator c_iter = cells->begin(); c_iter != cells->end(); ++c_iter) {
130: ALE::Obj<ALE::Mesh::sieve_type::supportSet> cellBlock = sieve->nSupport(sieve->nCone(*c_iter, depth), depth);
131: ALE::Obj<ALE::Mesh> sm = ALE::Selection<ALE::Mesh>::submesh(m, cellBlock);
132: ALE::Obj<ALE::Mesh::real_section_type> ssX = sm->getRealSection("default");
133: const ALE::Obj<ALE::Mesh::label_type>& cellMarker = sm->createLabel("marker");
135: if (fasDebug) {PetscPrintf(dmmg[0]->comm, " forward sweep cell %d\n", *c_iter);}
136: SectionRealSetSection(cellX, ssX);
137: // Assign BC to mesh
138: for(ALE::Mesh::sieve_type::supportSet::iterator b_iter = cellBlock->begin(); b_iter != cellBlock->end(); ++b_iter) {
139: const ALE::Obj<ALE::Mesh::coneArray> closure = ALE::SieveAlg<ALE::Mesh>::closure(m, *b_iter);
140: const ALE::Mesh::coneArray::iterator end = closure->end();
141: const bool isCell = *b_iter == *c_iter;
143: for(ALE::Mesh::coneArray::iterator cl_iter = closure->begin(); cl_iter != end; ++cl_iter) {
144: if (isCell) {
145: sm->setValue(cellMarker, *cl_iter, m->getValue(marker, *cl_iter));
146: } else {
147: if (sm->height(*cl_iter) == 0) {
148: sm->setValue(cellMarker, *cl_iter, 2);
149: } else if (sm->getValue(cellMarker, *cl_iter, -1) < 0) {
150: sm->setValue(cellMarker, *cl_iter, 1);
151: }
152: }
153: }
154: }
155: for(std::map<std::string, ALE::Obj<ALE::Discretization> >::iterator d_iter = sDiscs.begin(); d_iter != sDiscs.end(); ++d_iter) {
156: sm->setDiscretization(d_iter->first, d_iter->second);
157: }
158: // Create field
159: sm->setupField(ssX, 2, true);
160: // Setup constant
161: sm->setRealSection("constant", sB);
162: // Setup DMMG
163: MeshSetMesh(smallMesh, sm);
164: DMMGSetDM(smallDmmg, (DM) smallMesh);
165: DMMGSetSNESLocal(smallDmmg, func, jac, 0, 0);
166: // TODO: Construct null space, if necessary
167: //DMMGSetNullSpace(smallDmmg, PETSC_FALSE, 1, CreateNullSpace);
168: //ALE::Obj<ALE::Mesh::real_section_type> nullSpace = sm->getRealSection("nullSpace");
169: //sm->setupField(nullSpace, 2, true);
170: // Fill in intial guess with BC values
171: for(ALE::Mesh::sieve_type::supportSet::iterator b_iter = cellBlock->begin(); b_iter != cellBlock->end(); ++b_iter) {
172: sm->updateAll(ssX, *b_iter, m->restrictNew(sX, *b_iter));
173: }
174: if (fasDebug) {
175: sX->view("Initial solution guess");
176: ssX->view("Cell solution guess");
177: }
178: // Solve
179: DMMGSolve(smallDmmg);
180: // Update global solution with local solution
181: SectionRealToVec(cellX, smallMesh, SCATTER_REVERSE, DMMGGetx(smallDmmg));
182: m->updateAll(sX, *c_iter, sm->restrictNew(ssX, *c_iter));
183: if (fasDebug) {
184: ssX->view("Cell solution final");
185: sX->view("Final solution");
186: }
187: }
188: }
189: if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
190: }
191: }
192: sB->zero();
193: SectionRealToVec(sectionX, mesh, SCATTER_FORWARD, X);
194: SectionRealDestroy(sectionX);
195: SectionRealDestroy(sectionB);
196: SectionRealDestroy(cellX);
197: DMMGDestroy(smallDmmg);
198: MeshDestroy(smallMesh);
199: return(0);
200: }
203: /*
204: This is alpha FAS code.
206: R is the usual multigrid restriction (e.g. the tranpose of piecewise linear interpolation)
207: Q is either a scaled injection or the usual R
208: */
211: PetscErrorCode DMMGSolveFAS_Mesh(DMMG *dmmg, PetscInt level)
212: {
213: SNES snes = dmmg[level]->snes;
214: PetscReal norm;
215: PetscInt i, j, k;
216: PetscTruth fasDebug;
220: PetscOptionsHasName(dmmg[0]->prefix, "-dmmg_fas_debug", &fasDebug);
221: VecSet(dmmg[level]->r, 0.0);
222: /* for(j = 1; j <= level; ++j) { */
223: /* if (!dmmg[j]->inject) { */
224: /* DMGetInjection(dmmg[j-1]->dm, dmmg[j]->dm, &dmmg[j]->inject); */
225: /* } */
226: /* } */
228: for(i = 0, snes->iter = 1; i < 100; ++i, ++snes->iter) {
229: PetscPrintf(dmmg[0]->comm, "FAS iteration %d\n", i);
230: for(j = level; j > 0; j--) {
231: if (dmmg[j]->monitorall) {PetscPrintf(dmmg[0]->comm, " FAS level %d\n", j);}
232: /* Relax on fine mesh to obtain x^{new}_{fine}, residual^{new}_{fine} = F_{fine}(x^{new}_{fine}) \approx 0 */
233: Relax_Mesh(dmmg, (Mesh) dmmg[j]->dm, SOR_SYMMETRIC_SWEEP, dmmg[j]->presmooth, dmmg[j]->x, dmmg[j]->r);
234: DMMGFormFunctionMesh(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);
236: /* residual^{old}_fine} - residual^{new}_{fine} = F(x^{old}_{fine}) - residual^{new}_{fine} */
237: VecAYPX(dmmg[j]->w,-1.0,dmmg[j]->r);
239: if (j == level || dmmg[j]->monitorall) {
240: /* norm( residual_fine - f(x_fine) ) */
241: VecNorm(dmmg[j]->w,NORM_2,&norm);
242: if (dmmg[j]->monitorall) {
243: for (k=0; k<level-j+1; k++) {PetscPrintf(dmmg[j]->comm," ");}
244: PetscPrintf(dmmg[j]->comm,"FAS lvl %d function norm %G\n",j,norm);
245: }
246: if (j == level) {
247: if (norm < dmmg[level]->abstol) goto theend;
248: if (i == 0) {
249: dmmg[level]->rrtol = norm*dmmg[level]->rtol;
250: } else {
251: if (norm < dmmg[level]->rrtol) goto theend;
252: }
253: }
254: }
256: /* residual^{new}_{coarse} = R*(residual^{old}_fine} - residual^{new}_{fine}) */
257: MatRestrict(dmmg[j]->R, dmmg[j]->w, dmmg[j-1]->r);
258:
259: /* F_{coarse}(R*x^{new}_{fine}) */
260: MatRestrict(dmmg[j]->R, dmmg[j]->x, dmmg[j-1]->x);
261: /* VecScatterBegin(dmmg[j]->inject,dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD); */
262: /* VecScatterEnd(dmmg[j]->inject,dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD); */
263: DMMGFormFunctionMesh(0,dmmg[j-1]->x,dmmg[j-1]->w,dmmg[j-1]);
265: /* residual_coarse = F_{coarse}(R*x_{fine}) + R*(residual^{old}_fine} - residual^{new}_{fine}) */
266: VecAYPX(dmmg[j-1]->r,1.0,dmmg[j-1]->w);
268: /* save R*x^{new}_{fine} into b (needed when interpolating compute x back up) */
269: VecCopy(dmmg[j-1]->x,dmmg[j-1]->b);
270: }
272: if (dmmg[0]->monitorall) {
273: for (k=0; k<level+1; k++) {PetscPrintf(dmmg[0]->comm," ");}
274: PetscPrintf(dmmg[0]->comm, "FAS coarse grid\n");
275: }
276: if (level == 0) {
277: DMMGFormFunctionMesh(0,dmmg[0]->x,dmmg[0]->w,dmmg[0]);
278: VecAYPX(dmmg[j]->w,-1.0,dmmg[j]->r);
279: VecNorm(dmmg[0]->w,NORM_2,&norm);
280: if (norm < dmmg[level]->abstol) goto theend;
281: if (i == 0) {
282: dmmg[level]->rrtol = norm*dmmg[level]->rtol;
283: }
284: }
285: Relax_Mesh(dmmg, (Mesh) dmmg[0]->dm, SOR_SYMMETRIC_SWEEP, dmmg[0]->coarsesmooth, dmmg[0]->x, dmmg[0]->r);
286: if (level == 0 || dmmg[0]->monitorall) {
287: DMMGFormFunctionMesh(0,dmmg[0]->x,dmmg[0]->w,dmmg[0]);
288: if (fasDebug) {
289: SectionReal residual;
291: MeshGetSectionReal((Mesh) dmmg[0]->dm, "default", &residual);
292: SectionRealView(residual, PETSC_VIEWER_STDOUT_WORLD);
293: SectionRealDestroy(residual);
294: }
295: VecAXPY(dmmg[0]->w,-1.0,dmmg[0]->r);
296: VecNorm(dmmg[0]->w,NORM_2,&norm);
297: for (k=0; k<level+1; k++) {PetscPrintf(dmmg[0]->comm," ");}
298: PetscPrintf(dmmg[0]->comm,"FAS coarse grid function norm %G\n",norm);
299: if (level == 0) {
300: if (norm < dmmg[level]->abstol) goto theend;
301: if (norm < dmmg[level]->rrtol) goto theend;
302: }
303: }
305: for (j=1; j<=level; j++) {
306: PetscPrintf(dmmg[0]->comm, " FAS level %d\n", j);
307: /* x^{new}_{coarse} - R*x^{new}_{fine} */
308: VecAXPY(dmmg[j-1]->x,-1.0,dmmg[j-1]->b);
309: /* x_fine = x_fine + R'*(x^{new}_{coarse} - R*x^{new}_{fine}) */
310: MatInterpolateAdd(dmmg[j]->R, dmmg[j-1]->x, dmmg[j]->x, dmmg[j]->x);
312: if (dmmg[j]->monitorall) {
313: /* norm( F(x_fine) - residual_fine ) */
314: DMMGFormFunctionMesh(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);
315: VecAXPY(dmmg[j]->w,-1.0,dmmg[j]->r);
316: VecNorm(dmmg[j]->w,NORM_2,&norm);
317: for (k=0; k<level-j+1; k++) {PetscPrintf(dmmg[j]->comm," ");}
318: PetscPrintf(dmmg[j]->comm,"FAS lvl %d function norm before postsmooth %G\n",j,norm);
319: }
321: /* Relax residual_fine - F(x_fine) = 0 */
322: for (k=0; k<dmmg[j]->postsmooth; k++) {
323: Relax_Mesh(dmmg, (Mesh) dmmg[j]->dm, SOR_SYMMETRIC_SWEEP, 1, dmmg[j]->x, dmmg[j]->r);
324: }
326: if ((j == level) || dmmg[j]->monitorall) {
327: /* norm( F(x_fine) - residual_fine ) */
328: DMMGFormFunctionMesh(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);
329: VecAXPY(dmmg[j]->w,-1.0,dmmg[j]->r);
330: VecNorm(dmmg[j]->w,NORM_2,&norm);
331: for (k=0; k<level-j+1; k++) {PetscPrintf(dmmg[j]->comm," ");}
332: PetscPrintf(dmmg[j]->comm,"FAS lvl %d function norm %G\n",j,norm);
333: if (j == level) {
334: if (norm < dmmg[level]->abstol) goto theend;
335: if (norm < dmmg[level]->rrtol) goto theend;
336: }
337: }
338: }
340: if (dmmg[level]->monitor){
341: DMMGFormFunctionMesh(0,dmmg[level]->x,dmmg[level]->w,dmmg[level]);
342: VecNorm(dmmg[level]->w,NORM_2,&norm);
343: PetscPrintf(dmmg[level]->comm,"%D FAS function norm %G\n",i+1,norm);
344: }
345: }
346: theend:
347: return(0);
348: }
349: #endif