Actual source code: meshmgsnes.c

  1: #include <petscmg.h>      /*I      "petscmg.h"    I*/
  2: #include <petscdmmg.h>    /*I      "petscdmmg.h"  I*/
  3: #include <petscmesh.h>    /*I      "petscmesh.h"  I*/
  4: #include <Selection.hh>

  6: /* Just to set iterations */
 7:  #include private/snesimpl.h

  9: PetscErrorCode DMMGFormFunctionMesh(SNES snes, Vec X, Vec F, void *ptr);

 11: #if 0
 12: PetscErrorCode CreateNullSpace(DMMG dmmg, Vec *nulls) {
 13:   Mesh           mesh = (Mesh) dmmg->dm;
 14:   Vec            nS   = nulls[0];
 15:   SectionReal    nullSpace;

 19:   MeshGetSectionReal(mesh, "nullSpace", &nullSpace);
 20:   {
 21:     ALE::Obj<PETSC_MESH_TYPE> m;
 22:     ALE::Obj<PETSC_MESH_TYPE::real_section_type> s;

 24:     MeshGetMesh(mesh, m);
 25:     SectionRealGetSection(nullSpace, s);
 26:     ALE::Obj<ALE::Discretization> disc = m->getDiscretization("p");
 27:     const int dim = m->getDimension();

 29:     for(int d = 0; d <= dim; ++d) {
 30:       const int numDof = disc->getNumDof(d);

 32:       if (numDof) {
 33:         const ALE::Obj<PETSC_MESH_TYPE::label_sequence>& stratum = m->depthStratum(d);
 34:         const PETSC_MESH_TYPE::label_sequence::iterator  end     = stratum->end();
 35:         double                                    *values  = new double[numDof];

 37:         for(PETSC_MESH_TYPE::label_sequence::iterator p_iter = stratum->begin(); p_iter != end; ++p_iter) {
 38:           for(int i = 0; i < numDof; ++i) values[i] = 1.0;
 39:           s->updatePoint(*p_iter, values);
 40:         }
 41:       }
 42:     }
 43:   }
 44:   SectionRealToVec(nullSpace, mesh, SCATTER_FORWARD, nS);
 45:   std::cout << "Null space:" << std::endl;
 46:   VecView(nS, PETSC_VIEWER_STDOUT_SELF);
 47:   SectionRealDestroy(nullSpace);
 48:   return(0);
 49: }
 50: #endif

 52: /* Nonlinear relaxation on all the equations with an initial guess in x */
 56: PetscErrorCode  Relax_Mesh(DMMG *dmmg, Mesh mesh, MatSORType flag, int its, Vec X, Vec B)
 57: {
 58:   SectionReal      sectionX, sectionB, cellX;
 59:   Mesh             smallMesh;
 60:   DMMG            *smallDmmg;
 61:   DALocalFunction1 func;
 62:   DALocalFunction1 jac;
 63:   ALE::Obj<PETSC_MESH_TYPE> m;
 64:   ALE::Obj<PETSC_MESH_TYPE::real_section_type> sX;
 65:   ALE::Obj<PETSC_MESH_TYPE::real_section_type> sB;
 66:   PetscTruth       fasDebug;
 67:   PetscErrorCode   ierr;

 70:   PetscOptionsHasName(dmmg[0]->prefix, "-dmmg_fas_debug", &fasDebug);
 71:   if (fasDebug) {PetscPrintf(dmmg[0]->comm, "  FAS mesh relaxation\n");}
 72:   if (its <= 0) SETERRQ1(PETSC_ERR_ARG_WRONG, "Relaxation requires global its %D positive", its);
 73:   MeshCreate(PETSC_COMM_SELF, &smallMesh);
 74:   DMMGCreate(PETSC_COMM_SELF, -1, PETSC_NULL, &smallDmmg);
 75:   //DMMGSetMatType(smallDmmg, MATSEQDENSE);
 76:   DMMGSetOptionsPrefix(smallDmmg, "fas_");
 77:   DMMGSetUser(smallDmmg, 0, DMMGGetUser(dmmg, 0));
 78:   DMMGGetSNESLocal(dmmg, &func, &jac);
 79:   MeshGetMesh(mesh, m);
 80:   MeshGetSectionReal(mesh, "default", &sectionX);
 81:   SectionRealToVec(sectionX, mesh, SCATTER_REVERSE, X);
 82:   SectionRealGetSection(sectionX, sX);
 83:   MeshGetSectionReal(mesh, "constant", &sectionB);
 84:   SectionRealToVec(sectionB, mesh, SCATTER_REVERSE, B);
 85:   SectionRealGetSection(sectionB, sB);
 86:   SectionRealCreate(PETSC_COMM_SELF, &cellX);
 87:   //const ALE::Obj<PETSC_MESH_TYPE::sieve_type>&     sieve   = m->getSieve();
 88:   //const ALE::Obj<PETSC_MESH_TYPE::label_sequence>& cells   = m->heightStratum(0);
 89:   //const int                                  depth   = m->depth();
 90:   //const ALE::Obj<PETSC_MESH_TYPE::label_type>&     marker  = m->getLabel("marker");
 91:   //const int                                  cellDof = m->sizeWithBC(sX, *cells->begin());

 93: #ifdef PETSC_OPT_SIEVE
 94:   SETERRQ(PETSC_ERR_SUP, "I am being lazy, bug me.");
 95: #else
 96:   ALE::Obj<PETSC_MESH_TYPE::names_type> fields = m->getDiscretizations();
 97:   std::map<std::string, ALE::Obj<ALE::Discretization> > sDiscs;

 99:   for(PETSC_MESH_TYPE::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(PETSC_MESH_TYPE::label_sequence::iterator c_iter = cells->begin(); c_iter != cells->end(); ++c_iter) {
130:         ALE::Obj<PETSC_MESH_TYPE::sieve_type::supportSet> cellBlock  = sieve->nSupport(sieve->nCone(*c_iter, depth), depth);
131:         ALE::Obj<PETSC_MESH_TYPE>                         sm         = ALE::Selection<PETSC_MESH_TYPE>::submesh(m, cellBlock);
132:         ALE::Obj<PETSC_MESH_TYPE::real_section_type>      ssX        = sm->getRealSection("default");
133:         const ALE::Obj<PETSC_MESH_TYPE::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(PETSC_MESH_TYPE::sieve_type::supportSet::iterator b_iter = cellBlock->begin(); b_iter != cellBlock->end(); ++b_iter) {
139:           const ALE::Obj<PETSC_MESH_TYPE::coneArray> closure = ALE::SieveAlg<PETSC_MESH_TYPE>::closure(m, *b_iter);
140:           const PETSC_MESH_TYPE::coneArray::iterator end     = closure->end();
141:           const bool                           isCell  = *b_iter == *c_iter;

143:           for(PETSC_MESH_TYPE::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:         DMMGSetFromOptions(smallDmmg);
167:         // TODO: Construct null space, if necessary
168:         //DMMGSetNullSpace(smallDmmg, PETSC_FALSE, 1, CreateNullSpace);
169:         //ALE::Obj<PETSC_MESH_TYPE::real_section_type> nullSpace = sm->getRealSection("nullSpace");
170:         //sm->setupField(nullSpace, 2, true);
171:         // Fill in intial guess with BC values
172:         for(PETSC_MESH_TYPE::sieve_type::supportSet::iterator b_iter = cellBlock->begin(); b_iter != cellBlock->end(); ++b_iter) {
173:           sm->updateAll(ssX, *b_iter, m->restrictNew(sX, *b_iter));
174:         }
175:         if (fasDebug) {
176:           sX->view("Initial solution guess");
177:           ssX->view("Cell solution guess");
178:         }
179:         // Solve
180:         DMMGSolve(smallDmmg);
181:         // Update global solution with local solution
182:         SectionRealToVec(cellX, smallMesh, SCATTER_REVERSE, DMMGGetx(smallDmmg));
183:         m->updateAll(sX, *c_iter, sm->restrictNew(ssX, *c_iter));
184:         if (fasDebug) {
185:           ssX->view("Cell solution final");
186:           sX->view("Final solution");
187:         }
188:       }
189:     }
190:     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
191:     }
192:   }
193: #endif
194:   sB->zero();
195:   SectionRealToVec(sectionX, mesh, SCATTER_FORWARD, X);
196:   SectionRealDestroy(sectionX);
197:   SectionRealDestroy(sectionB);
198:   SectionRealDestroy(cellX);
199:   DMMGDestroy(smallDmmg);
200:   MeshDestroy(smallMesh);
201:   return(0);
202: }

205: /*
206:  This is alpha FAS code.

208:   R is the usual multigrid restriction (e.g. the tranpose of piecewise linear interpolation)
209:   Q is either a scaled injection or the usual R
210: */
213: PetscErrorCode DMMGSolveFAS_Mesh(DMMG *dmmg, PetscInt level)
214: {
215:   SNES           snes = dmmg[level]->snes;
216:   PetscReal      norm;
217:   PetscInt       i, j, k;
218:   PetscTruth     fasDebug;

222:   PetscOptionsHasName(dmmg[0]->prefix, "-dmmg_fas_debug", &fasDebug);
223:   VecSet(dmmg[level]->r, 0.0);
224: /*   for(j = 1; j <= level; ++j) { */
225: /*     if (!dmmg[j]->inject) { */
226: /*       DMGetInjection(dmmg[j-1]->dm, dmmg[j]->dm, &dmmg[j]->inject); */
227: /*     } */
228: /*   } */

230:   for(i = 0, snes->iter = 1; i < 100; ++i, ++snes->iter) {
231:     PetscPrintf(dmmg[0]->comm, "FAS iteration %d\n", i);
232:     for(j = level; j > 0; j--) {
233:       if (dmmg[j]->monitorall) {PetscPrintf(dmmg[0]->comm, "  FAS level %d\n", j);}
234:       /* Relax on fine mesh to obtain x^{new}_{fine}, residual^{new}_{fine} = F_{fine}(x^{new}_{fine}) \approx 0 */
235:       Relax_Mesh(dmmg, (Mesh) dmmg[j]->dm, SOR_SYMMETRIC_SWEEP, dmmg[j]->presmooth, dmmg[j]->x, dmmg[j]->r);
236:       DMMGFormFunctionMesh(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);

238:       /* residual^{old}_fine} - residual^{new}_{fine} = F(x^{old}_{fine}) - residual^{new}_{fine} */
239:       VecAYPX(dmmg[j]->w,-1.0,dmmg[j]->r);

241:       if (j == level || dmmg[j]->monitorall) {
242:         /* norm( residual_fine - f(x_fine) ) */
243:         VecNorm(dmmg[j]->w,NORM_2,&norm);
244:         if (dmmg[j]->monitorall) {
245:           for (k=0; k<level-j+1; k++) {PetscPrintf(dmmg[j]->comm,"  ");}
246:           PetscPrintf(dmmg[j]->comm,"FAS lvl %d function norm %G\n",j,norm);
247:         }
248:         if (j == level) {
249:           if (norm < dmmg[level]->abstol) goto theend;
250:           if (i == 0) {
251:             dmmg[level]->rrtol = norm*dmmg[level]->rtol;
252:           } else {
253:             if (norm < dmmg[level]->rrtol) goto theend;
254:           }
255:         }
256:       }

258:       /* residual^{new}_{coarse} = R*(residual^{old}_fine} - residual^{new}_{fine}) */
259:       MatRestrict(dmmg[j]->R, dmmg[j]->w, dmmg[j-1]->r);
260: 
261:       /* F_{coarse}(R*x^{new}_{fine}) */
262:       MatRestrict(dmmg[j]->R, dmmg[j]->x, dmmg[j-1]->x);
263: /*       VecScatterBegin(dmmg[j]->inject,dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD); */
264: /*       VecScatterEnd(dmmg[j]->inject,dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD); */
265:       DMMGFormFunctionMesh(0,dmmg[j-1]->x,dmmg[j-1]->w,dmmg[j-1]);

267:       /* residual_coarse = F_{coarse}(R*x_{fine}) + R*(residual^{old}_fine} - residual^{new}_{fine}) */
268:       VecAYPX(dmmg[j-1]->r,1.0,dmmg[j-1]->w);

270:       /* save R*x^{new}_{fine} into b (needed when interpolating compute x back up) */
271:       VecCopy(dmmg[j-1]->x,dmmg[j-1]->b);
272:     }

274:     if (dmmg[0]->monitorall) {
275:       for (k=0; k<level+1; k++) {PetscPrintf(dmmg[0]->comm,"  ");}
276:       PetscPrintf(dmmg[0]->comm, "FAS coarse grid\n");
277:     }
278:     if (level == 0) {
279:       DMMGFormFunctionMesh(0,dmmg[0]->x,dmmg[0]->w,dmmg[0]);
280:       VecAYPX(dmmg[j]->w,-1.0,dmmg[j]->r);
281:       VecNorm(dmmg[0]->w,NORM_2,&norm);
282:       if (norm < dmmg[level]->abstol) goto theend;
283:       if (i == 0) {
284:         dmmg[level]->rrtol = norm*dmmg[level]->rtol;
285:       }
286:     }
287:     Relax_Mesh(dmmg, (Mesh) dmmg[0]->dm, SOR_SYMMETRIC_SWEEP, dmmg[0]->coarsesmooth, dmmg[0]->x, dmmg[0]->r);
288:     if (level == 0 || dmmg[0]->monitorall) {
289:       DMMGFormFunctionMesh(0,dmmg[0]->x,dmmg[0]->w,dmmg[0]);
290:       if (fasDebug) {
291:         SectionReal residual;

293:         MeshGetSectionReal((Mesh) dmmg[0]->dm, "default", &residual);
294:         SectionRealView(residual, PETSC_VIEWER_STDOUT_WORLD);
295:         SectionRealDestroy(residual);
296:       }
297:       VecAXPY(dmmg[0]->w,-1.0,dmmg[0]->r);
298:       VecNorm(dmmg[0]->w,NORM_2,&norm);
299:       for (k=0; k<level+1; k++) {PetscPrintf(dmmg[0]->comm,"  ");}
300:       PetscPrintf(dmmg[0]->comm,"FAS coarse grid function norm %G\n",norm);
301:       if (level == 0) {
302:         if (norm < dmmg[level]->abstol) goto theend;
303:         if (norm < dmmg[level]->rrtol)  goto theend;
304:       }
305:     }

307:     for (j=1; j<=level; j++) {
308:       PetscPrintf(dmmg[0]->comm, "  FAS level %d\n", j);
309:       /* x^{new}_{coarse} - R*x^{new}_{fine} */
310:       VecAXPY(dmmg[j-1]->x,-1.0,dmmg[j-1]->b);
311:       /* x_fine = x_fine + R'*(x^{new}_{coarse} - R*x^{new}_{fine}) */
312:       MatInterpolateAdd(dmmg[j]->R, dmmg[j-1]->x, dmmg[j]->x, dmmg[j]->x);

314:       if (dmmg[j]->monitorall) {
315:         /* norm( F(x_fine) - residual_fine ) */
316:         DMMGFormFunctionMesh(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);
317:         VecAXPY(dmmg[j]->w,-1.0,dmmg[j]->r);
318:         VecNorm(dmmg[j]->w,NORM_2,&norm);
319:         for (k=0; k<level-j+1; k++) {PetscPrintf(dmmg[j]->comm,"  ");}
320:         PetscPrintf(dmmg[j]->comm,"FAS lvl %d function norm before postsmooth %G\n",j,norm);
321:       }

323:       /* Relax residual_fine - F(x_fine)  = 0 */
324:       for (k=0; k<dmmg[j]->postsmooth; k++) {
325:         Relax_Mesh(dmmg, (Mesh) dmmg[j]->dm, SOR_SYMMETRIC_SWEEP, 1, dmmg[j]->x, dmmg[j]->r);
326:       }

328:       if ((j == level) || dmmg[j]->monitorall) {
329:         /* norm( F(x_fine) - residual_fine ) */
330:         DMMGFormFunctionMesh(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);
331:         VecAXPY(dmmg[j]->w,-1.0,dmmg[j]->r);
332:         VecNorm(dmmg[j]->w,NORM_2,&norm);
333:         for (k=0; k<level-j+1; k++) {PetscPrintf(dmmg[j]->comm,"  ");}
334:         PetscPrintf(dmmg[j]->comm,"FAS lvl %d function norm %G\n",j,norm);
335:         if (j == level) {
336:           if (norm < dmmg[level]->abstol) goto theend;
337:           if (norm < dmmg[level]->rrtol) goto theend;
338:         }
339:       }
340:     }

342:     if (dmmg[level]->monitor){
343:       DMMGFormFunctionMesh(0,dmmg[level]->x,dmmg[level]->w,dmmg[level]);
344:       VecNorm(dmmg[level]->w,NORM_2,&norm);
345:       PetscPrintf(dmmg[level]->comm,"%D FAS function norm %G\n",i+1,norm);
346:     }
347:   }
348:   theend:
349:   return(0);
350: }