1 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2 // vi: set et ts=4 sw=2 sts=2:
3 
4 #include <config.h>
5 
6 #include <set>
7 #include <map>
8 #include <memory>
9 
10 #include <dune/grid/uggrid.hh>
11 
12 #if ModelP and DUNE_UGGRID_HAVE_PPIFCONTEXT
13 #  include <dune/uggrid/parallel/ppif/ppifcontext.hh>
14 #endif
15 
16 /** \todo Remove the following two includes once getAllSubfaces... is gone */
17 #include <list>
18 #include <iterator>
19 #include <dune/common/stdstreams.hh>
20 #include <dune/grid/common/mcmgmapper.hh>
21 
22 namespace Dune {
23 
24 //***********************************************************************
25 //
26 // --UGGrid
27 // --Grid
28 //
29 //***********************************************************************
30 
31 template<> int UGGrid<2>::numOfUGGrids = 0;
32 template<> int UGGrid<3>::numOfUGGrids = 0;
33 
34 
35 template <int dim>
UGGrid(UGCollectiveCommunication comm)36 UGGrid < dim >::UGGrid(UGCollectiveCommunication comm)
37   : multigrid_(nullptr),
38     ccobj_(comm),
39     leafIndexSet_(*this),
40     idSet_(*this),
41     refinementType_(LOCAL),
42     closureType_(GREEN),
43     someElementHasBeenMarkedForRefinement_(false),
44     someElementHasBeenMarkedForCoarsening_(false),
45     numBoundarySegments_(0)
46 {
47   // If no UGGrid object exists yet start up UG for 2d and 3d
48   if ( (UGGrid<2>::numOfUGGrids+UGGrid<3>::numOfUGGrids)==0) {
49 
50     /* Why do we start UG for 2d and 3d?  Why not just the one for the grid dimension
51      * that we actually need?  The reason is a certain error condition that you get
52      * when using 2d and 3d UGGrids together. This was FlySpray task 887.
53      *
54      * UG has an internal directory-tree-like structure for storing different types of
55      * miscellaneous data.  It is called the 'environment heap'. The method InitUg
56      * creates a directory 'BVP' there, which holds, among other things, a 'problem name'
57      * (a string) for each existing grid.  When you call InitUg twice (i.e., once for 2d
58      * and once for 3d), the second call will not notice that a directory 'BVP' already exists.
59      * This is possible because these directories have a name and a 'type', (an integer
60      * of which I have not really understood what you need it for), and the new 'BVP' directory
61      * has a different type. The new one will appear in the list before the old one.
62      * However, some methods that look for entries in the environment heap only look for the name,
63      * and hence only find the first 'BVP' entry, which only contains the information about
64      * the second grid.  Surprisingly, this doesn't appear to be a problem, at least the dune-grid
65      * unit test for UG successfully tests a 2d and a 3d grid side by side.  However when trying
66      * to delete grids the method DisposeMultiGrid tries to remove a grid's information from the
67      * BVP directory.  However the information of the first grid is in the second BVP directory,
68      * which is effectively hidden by the first one.  DisposeMultiGrid doesn't find the
69      * grid information and reports an error.
70      *
71      * If both calls to InitUg are directly in a row, then all information will be stored
72      * in the first instance of 'BVP', and there is no error.
73      */
74     int argc = 1;
75     char* arg = strdup("dune.exe");
76     char** argv = &arg;
77 
78     if (UG_NS<2>::InitUg(&argc, &argv))
79       DUNE_THROW(GridError, "UG" << dim << "d::InitUg() returned an error code!");
80 
81     if (UG_NS<3>::InitUg(&argc, &argv))
82       DUNE_THROW(GridError, "UG" << dim << "d::InitUg() returned an error code!");
83 
84     free(arg);
85   }
86 
87   // Create a dummy problem
88   typename UG_NS<dim>::CoeffProcPtr coeffs[1] = {nullptr};
89   typename UG_NS<dim>::UserProcPtr upp[1] = {nullptr};
90 
91   // Create unique problem name
92   std::stringstream numberAsAscii;
93   numberAsAscii << numOfUGGrids;
94   name_ = "DuneUGGrid_" + std::string((dim==2) ? "2" : "3") + std::string("d_") + numberAsAscii.str();
95 
96   std::string problemName = name_ + "_Problem";
97 
98   if (UG_NS<dim>::CreateBoundaryValueProblem(problemName.c_str(), 1,coeffs,1,upp) == nullptr)
99     DUNE_THROW(GridError, "UG" << dim << "d::CreateBoundaryValueProblem() returned an error code!");
100 
101   numOfUGGrids++;
102 
103   dverb << "UGGrid<" << dim << "> with name " << name_ << " created!" << std::endl;
104 
105 }
106 
107 
108 template < int dim >
~UGGrid()109 UGGrid < dim >::~UGGrid() noexcept(false)
110 {
111   // Delete the UG multigrid if there is one (== createEnd() has been called)
112   if (multigrid_) {
113     // Set UG's currBVP variable to the BVP corresponding to this
114     // grid.  This is necessary if we have more than one UGGrid in use.
115     // DisposeMultiGrid will crash if we don't do this
116     UG_NS<dim>::Set_Current_BVP(multigrid_->theBVP);
117     if (UG_NS<dim>::DisposeMultiGrid(multigrid_) != 0)
118       DUNE_THROW(GridError, "UG" << dim << "d::DisposeMultiGrid returned error code!");
119   }
120 
121   // DisposeMultiGrid cleans up the BVP as well.  But if there was no
122   // multigrid we have to take care of the BVP ourselves.
123   std::string problemName = name_ + "_Problem";
124   void** BVP = UG_NS<dim>::BVP_GetByName(problemName.c_str());
125 
126   if (BVP)
127     if (UG_NS<dim>::BVP_Dispose(BVP))
128       DUNE_THROW(GridError, "Couldn't dispose of UG boundary value problem!");
129 
130 
131 
132   numOfUGGrids--;
133 
134   // Shut down UG if this was the last existing UGGrid object.
135   // Since we started both 2d and 3d versions of UG even if we don't need both,
136   // we have to shut them down both, too.
137   if (UGGrid<2>::numOfUGGrids + UGGrid<3>::numOfUGGrids == 0) {
138     UG_NS<2>::ExitUg();
139     UG_NS<3>::ExitUg();
140   }
141 }
142 
143 template < int dim >
maxLevel() const144 int UGGrid < dim >::maxLevel() const
145 {
146   if (!multigrid_)
147     DUNE_THROW(GridError, "The grid has not been properly initialized!");
148 
149   return multigrid_->topLevel;
150 }
151 
152 
153 template < int dim >
size(int level,int codim) const154 int UGGrid < dim >::size (int level, int codim) const
155 {
156   return levelIndexSet(level).size(codim);
157 }
158 
159 
160 template < int dim >
mark(int refCount,const typename Traits::template Codim<0>::Entity & e)161 bool UGGrid < dim >::mark(int refCount,
162                           const typename Traits::template Codim<0>::Entity& e )
163 {
164   typename UG_NS<dim>::Element* target = e.impl().target_;
165 
166   // No refinement requested
167   if (refCount==0) {
168     if (UG_NS<dim>::MarkForRefinement(target,
169                                       UG_NS<dim>::NO_REFINEMENT,      // unset
170                                       0)      // Irrelevant if refinement rule is not BLUE
171         ) DUNE_THROW(GridError, "UG" << dim << "d::MarkForRefinement returned error code!");
172     return true;
173   }
174 
175   // Check whether element can be marked for refinement
176   if (!EstimateHere(target))
177     return false;
178 
179   if (refCount==1) {
180     if (UG_NS<dim>::MarkForRefinement(target,
181                                       UG_NS<dim>::RED,      // red refinement rule
182                                       0)      // Irrelevant if refinement rule is not BLUE
183         ) DUNE_THROW(GridError, "UG" << dim << "d::MarkForRefinement returned error code!");
184 
185     someElementHasBeenMarkedForRefinement_ = true;
186     return true;
187   } else if (refCount==-1) {
188 
189     if (UG_NS<dim>::MarkForRefinement(target,
190                                       UG_NS<dim>::COARSE,      // coarsen the element
191                                       0)      // Irrelevant if refinement rule is not BLUE
192         ) DUNE_THROW(GridError, "UG" << dim << "d::MarkForRefinement returned error code!");
193 
194     someElementHasBeenMarkedForCoarsening_ = true;
195     return true;
196   } else
197     DUNE_THROW(GridError, "UGGrid only supports refCount values -1, 0, and 1 for mark()!");
198 
199 }
200 
201 template < int dim >
mark(const typename Traits::template Codim<0>::Entity & e,typename UG_NS<dim>::RefinementRule rule,int side)202 bool UGGrid < dim >::mark(const typename Traits::template Codim<0>::Entity& e,
203                           typename UG_NS<dim>::RefinementRule rule,
204                           int side)
205 {
206   typename UG_NS<dim>::Element* target = e.impl().target_;
207 
208   if (!UG_NS<dim>::isLeaf(target))
209     return false;
210 
211   someElementHasBeenMarkedForRefinement_ = true;
212 
213   return UG_NS<dim>::MarkForRefinement(target, rule, side);
214 
215 }
216 
217 template <int dim>
getMark(const typename Traits::template Codim<0>::Entity & e) const218 int UGGrid<dim>::getMark(const typename Traits::template Codim<0>::Entity& e) const
219 {
220   typename UG_NS<dim>::Element* target = e.impl().target_;
221 
222   // Return -1 if element is marked for coarsening
223   if (UG_NS<dim>::ReadCW(target,UG_NS<dim>::COARSEN_CE))
224     return -1;
225 
226   // If the element is not regular, it's mark is actually stored on its regular
227   // ancestor.  That's what we look for here.
228   if (dim==2)
229     target = (typename UG_NS<dim>::Element*) UG::D2::ELEMENT_TO_MARK((UG::D2::element*)target);
230   else
231     target = (typename UG_NS<dim>::Element*) UG::D3::ELEMENT_TO_MARK((UG::D3::element*)target);
232 
233   // Return 0 if element is not marked at all
234   if (UG_NS<dim>::ReadCW(target,UG_NS<dim>::MARK_CE)==UG_NS<dim>::NO_REFINEMENT)
235     return 0;
236 
237   // Else return 1
238   return 1;
239 }
240 
241 template <int dim>
preAdapt()242 bool UGGrid <dim>::preAdapt()
243 {
244   if( closureType_ == GREEN )
245   {
246     // when conforming refinement is enabled
247     // green closure has to be removed although not
248     // marked for coarsening
249     return someElementHasBeenMarkedForCoarsening_ ||
250            someElementHasBeenMarkedForRefinement_ ;
251   }
252   else
253   {
254     // in non-conforming meshes only elements marked for coarsening
255     // will be coarsened (hopefully)
256     return someElementHasBeenMarkedForCoarsening_;
257   }
258 }
259 
260 template < int dim >
adapt()261 bool UGGrid < dim >::adapt()
262 {
263   assert(multigrid_);
264 
265   // Set UG's currBVP variable to the BVP corresponding to this
266   // grid.  This is necessary if we have more than one UGGrid in use.
267   UG_NS<dim>::Set_Current_BVP(multigrid_->theBVP);
268 
269   int mode = UG_NS<dim>::GM_REFINE_TRULY_LOCAL;
270 
271   if (refinementType_==COPY)
272     mode = mode | UG_NS<dim>::GM_COPY_ALL;
273 
274   if (closureType_==NONE)
275     mode = mode | UG_NS<dim>::GM_REFINE_NOT_CLOSED;
276 
277   // I don't really know what this means
278   int seq = UG_NS<dim>::GM_REFINE_PARALLEL;
279 
280   // Skip test whether we have enough memory available
281   int mgtest = UG_NS<dim>::GM_REFINE_NOHEAPTEST;
282 
283   int rv = AdaptMultiGrid(multigrid_,mode,seq,mgtest);
284 
285   if (rv!=0)
286     DUNE_THROW(GridError, "UG::adapt() returned with error code " << rv);
287 
288   // Renumber everything
289   setIndices(false, nullptr);
290 
291   // Return true iff the grid hierarchy changed
292   //return !(bool)multigrid_->status;
293 
294   // grid has changed if at least one element was marked for refinement
295   return someElementHasBeenMarkedForRefinement_;
296 }
297 
298 template <int dim>
postAdapt()299 void UGGrid <dim>::postAdapt()
300 {
301   for (int i=0; i<=maxLevel(); i++)
302     for (const auto& element : elements(this->levelGridView(i)))
303       UG_NS<dim>::WriteCW(element.impl().target_, UG_NS<dim>::NEWEL_CE, 0);
304 
305   // reset marker flags
306   someElementHasBeenMarkedForRefinement_ = false;
307   someElementHasBeenMarkedForCoarsening_ = false;
308 }
309 
310 template < int dim >
globalRefine(int n)311 void UGGrid < dim >::globalRefine(int n)
312 {
313   for (int i=0; i<n; i++) {
314 
315     // mark all entities for grid refinement
316     for (const auto& element : elements(this->leafGridView()))
317       mark(1, element);
318 
319     this->preAdapt();
320     adapt();
321 
322   }
323 
324   this->postAdapt();
325 
326 }
327 
328 template <int dim>
getChildrenOfSubface(const typename Traits::template Codim<0>::Entity & e,int elementSide,int maxl,std::vector<typename Traits::template Codim<0>::Entity> & childElements,std::vector<unsigned char> & childElementSides) const329 void UGGrid<dim>::getChildrenOfSubface(const typename Traits::template Codim<0>::Entity & e,
330                                        int elementSide,
331                                        int maxl,
332                                        std::vector<typename Traits::template Codim<0>::Entity>& childElements,
333                                        std::vector<unsigned char>& childElementSides) const
334 {
335 
336   typedef std::pair<typename UG_NS<dim>::Element*,int> ListEntryType;
337 
338   std::list<ListEntryType> list;
339 
340   // //////////////////////////////////////////////////////////////////////
341   //   Change the input face number from Dune numbering to UG numbering
342   // //////////////////////////////////////////////////////////////////////
343 
344   elementSide = UGGridRenumberer<dim>::facesDUNEtoUG(elementSide, e.type());
345 
346   // ///////////////
347   //   init list
348   // ///////////////
349   if (!e.isLeaf()   // Get_Sons_of_ElementSide returns GM_FATAL when called for a leaf !?!
350       && e.level() < maxl) {
351 
352     typename UG_NS<dim>::Element* theElement = e.impl().target_;
353     list.emplace_back(theElement, elementSide);
354 
355   }
356 
357   // //////////////////////////////////////////////////
358   //   Traverse and collect all children of the side
359   // //////////////////////////////////////////////////
360 
361   for (const auto& f : list)
362   {
363     typename UG_NS<dim>::Element* theElement = f.first;
364     int side                                  = f.second;
365 
366     UG::INT Sons_of_Side = 0;
367     typename UG_NS<dim>::Element* SonList[UG_NS<dim>::MAX_SONS];
368     UG::INT SonSides[UG_NS<dim>::MAX_SONS];
369 
370     if (UG_NS<dim>::myLevel(theElement) < maxl) {
371 
372       int rv = Get_Sons_of_ElementSide(theElement,
373                                        side,          // Input element side number
374                                        &Sons_of_Side, // Number of topological sons of the element side
375                                        SonList,       // Output elements
376                                        SonSides,      // Output element side numbers
377                                        true,
378                                        true);
379 
380       if (rv != 0)
381         DUNE_THROW(GridError, "Get_Sons_of_ElementSide returned with error value " << rv);
382 
383       for (int i=0; i<Sons_of_Side; i++)
384         list.emplace_back(SonList[i], SonSides[i]);
385 
386     }
387   }
388 
389   // Remove the element itself. We are only interested in the children.
390   list.pop_front();
391 
392   // //////////////////////////////
393   //   Extract result from stack
394   // //////////////////////////////
395 
396   childElements.clear();
397   childElements.reserve(std::distance(list.begin(), list.end()));
398   childElementSides.resize(childElements.size());
399 
400   int i=0;
401   for (const auto &f : list)
402   {
403     // Set element
404     typedef typename Traits::template Codim< 0 >::Entity Entity;
405     childElements.push_back( Entity( UGGridEntity< 0, dim, const UGGrid< dim > >( f.first, this ) ) );
406 
407     int side = f.second;
408 
409     // Dune numbers the faces of several elements differently than UG.
410     // The following switch does the transformation
411     childElementSides[i] = UGGridRenumberer<dim>::facesUGtoDUNE(side, childElements[i].type());
412     ++i;
413   }
414 
415 }
416 
417 template < int dim >
loadBalance(int minlevel)418 bool UGGrid < dim >::loadBalance(int minlevel)
419 {
420   // Do nothing if we are on a single process
421   if (comm().size()==1)
422     return true;
423 
424   std::stringstream levelarg;
425   levelarg << minlevel;
426   UG_NS<dim>::lbs(levelarg.str().c_str(), multigrid_);
427 
428   // Renumber everything.
429   // Note: this must not be called when on a single process, because it renumbers the zero-level
430   // elements and vertices.
431   setIndices(true, nullptr);
432 
433   return true;
434 }
435 
436 template < int dim >
loadBalance(const std::vector<Rank> & targetProcessors,unsigned int fromLevel)437 bool UGGrid < dim >::loadBalance(const std::vector<Rank>& targetProcessors, unsigned int fromLevel)
438 {
439   // Do nothing if we are on a single process
440   if (comm().size()==1)
441     return true;
442 
443 #ifdef ModelP  // The Partition field only exists if ModelP is set
444   if (int(targetProcessors.size()) != this->leafGridView().size(0))
445     DUNE_THROW(Exception, "targetProcessors argument does not have the correct size");
446 
447   // Get unique consecutive index across different element types
448   typedef MultipleCodimMultipleGeomTypeMapper<typename Base::LeafGridView> ElementMapper;
449   ElementMapper elementMapper(this->leafGridView(), mcmgElementLayout());
450 
451   // Loop over all elements of all level, in decreasing level number.
452   // If the element is a leaf, take its target rank from the input targetProcessors array.
453   // If it is not, assign it to the processor most of its children are assigned to.
454   for (int i=maxLevel(); i>=0; i--) {
455     for (const auto& element : elements(this->levelGridView(i), Partitions::interior)) {
456 
457       if (element.isLeaf()) {
458 
459         int targetRank = targetProcessors[elementMapper.index(element)];
460 
461         // sanity check
462         if (targetRank >= comm().size())
463           DUNE_THROW(GridError, "Requesting target processor " << targetRank <<
464                      ", but only " << comm().size() << " processors are available.");
465 
466         UG_NS<dim>::Partition(element.impl().target_) = targetRank;
467       } else {
468 
469         std::map<Rank,unsigned int> rank;
470         Rank mostFrequentRank = 0;    // which rank occurred most often?
471         unsigned int mostFrequentCount = 0;   // how often did it occur?
472 
473         // Loop over all children and collect the ranks they are assigned to
474         for (const auto& child : descendantElements(element, element.level() + 1)) {
475 
476           auto childRank = UG_NS<dim>::Partition(child.impl().target_);
477 
478           if (rank.find(childRank) == rank.end())
479             rank[childRank] = 1;
480           else
481             rank[childRank]++;
482 
483           if (rank[childRank] > mostFrequentCount) {
484             mostFrequentRank = childRank;
485             mostFrequentCount = rank[childRank];
486           }
487 
488         }
489 
490         // Assign rank that occurred most often
491         UG_NS<dim>::Partition(element.impl().target_) = mostFrequentRank;
492 
493       }
494     }
495   }
496 #endif
497 
498   int errCode = UG_NS<dim>::TransferGridFromLevel(multigrid_, fromLevel);
499 
500   if (errCode)
501     DUNE_THROW(GridError, "UG" << dim << "d::TransferGridFromLevel returned error code " << errCode);
502 
503   // Renumber everything.
504   // Note: this must not be called when on a single process, because it renumbers the zero-level
505   // elements and vertices.
506   setIndices(true, nullptr);
507 
508   return true;
509 }
510 
511 #ifdef ModelP
512 template <int dim>
findDDDInterfaces(InterfaceType iftype,int codim) const513 std::vector<typename UG_NS<dim>::DDD_IF> UGGrid<dim>::findDDDInterfaces(InterfaceType iftype,
514                                                                         int codim) const
515 {
516   std::vector<typename UG_NS<dim>::DDD_IF> dddIfaces;
517 
518   // UGGrid does not have overlap or front entities
519   if (iftype == Overlap_OverlapFront_Interface || iftype == Overlap_All_Interface)
520     return dddIfaces;
521 
522   if (codim == 0)
523   {
524     switch (iftype) {
525     case InteriorBorder_InteriorBorder_Interface :
526       // do not communicate anything: Elements cannot be in
527       // the interior of two processes at the same time
528       break;
529     case InteriorBorder_All_Interface :
530       dddIfaces.push_back(UG_NS<dim>::ElementVHIF(multigrid_->dddContext()));
531       break;
532     case All_All_Interface :
533       dddIfaces.push_back(UG_NS<dim>::ElementSymmVHIF(multigrid_->dddContext()));
534       break;
535     default :
536       DUNE_THROW(GridError,
537                  "Element communication not supported for "
538                  "interfaces of type  "
539                  << iftype);
540     }
541   }
542   else if (codim == dim)
543   {
544     switch (iftype)
545     {
546     case InteriorBorder_InteriorBorder_Interface :
547       dddIfaces.push_back(UG_NS<dim>::BorderNodeSymmIF(multigrid_->dddContext()));
548       break;
549     case InteriorBorder_All_Interface :
550       dddIfaces.push_back(UG_NS<dim>::NodeInteriorBorderAllIF(multigrid_->dddContext()));
551       break;
552     case All_All_Interface :
553       dddIfaces.push_back(UG_NS<dim>::NodeAllIF(multigrid_->dddContext()));
554       break;
555     default :
556       DUNE_THROW(GridError,
557                  "Node communication not supported for "
558                  "interfaces of type  "
559                  << iftype);
560     }
561   }
562   else if (codim == dim-1)
563   {
564     switch (iftype)
565     {
566     case InteriorBorder_InteriorBorder_Interface :
567       dddIfaces.push_back(UG_NS<dim>::BorderEdgeSymmIF(multigrid_->dddContext()));
568       break;
569     case InteriorBorder_All_Interface :
570       dddIfaces.push_back(UG_NS<dim>::EdgeVHIF(multigrid_->dddContext()));
571       break;
572     case All_All_Interface :
573       dddIfaces.push_back(UG_NS<dim>::EdgeSymmVHIF(multigrid_->dddContext()));
574       break;
575     default :
576       DUNE_THROW(GridError,
577                  "Edge communication not supported for "
578                  "interfaces of type  "
579                  << iftype);
580     }
581   }
582   else if (codim == 1)
583   {
584     switch (iftype)
585     {
586     case InteriorBorder_InteriorBorder_Interface :
587       dddIfaces.push_back(UG_NS<dim>::BorderVectorSymmIF(multigrid_->dddContext()));
588       break;
589     case InteriorBorder_All_Interface :
590       dddIfaces.push_back(UG_NS<dim>::FacetInteriorBorderAllIF(multigrid_->dddContext()));
591       break;
592     case All_All_Interface :
593       dddIfaces.push_back(UG_NS<dim>::FacetAllAllIF(multigrid_->dddContext()));
594       break;
595     default :
596       DUNE_THROW(GridError,
597                  "Face communication not supported for "
598                  "interfaces of type  "
599                  << iftype);
600     }
601   }
602   else
603   {
604     DUNE_THROW(GridError,
605                "Communication codimension must be between 0 and " << dim << "!");
606   }
607 
608   return dddIfaces;
609 };
610 #endif // ModelP
611 
612 
613 template < int dim >
setPosition(const typename Traits::template Codim<dim>::Entity & e,const FieldVector<double,dim> & pos)614 void UGGrid < dim >::setPosition(const typename Traits::template Codim<dim>::Entity& e,
615                                  const FieldVector<double, dim>& pos)
616 {
617   typename UG_NS<dim>::Node* target = e.impl().target_;
618 
619   for (int i=0; i<dim; i++)
620     target->myvertex->iv.x[i] = pos[i];
621 }
622 
623 template <int dim>
saveState(const std::string & filename) const624 void UGGrid<dim>::saveState(const std::string& filename) const
625 {
626   const char* type = "asc";
627   const char* comment = "written by DUNE";
628 
629   if (dim==2)
630     UG::D2::SaveMultiGrid((UG::D2::multigrid*)multigrid_,
631                           filename.c_str(),
632                           type,
633                           comment,
634                           0,      // autosave
635                           0       // rename
636                           );
637   else
638     UG::D3::SaveMultiGrid((UG::D3::multigrid*)multigrid_,
639                           filename.c_str(),
640                           type,
641                           comment,
642                           0,      // autosave
643                           0       // rename
644                           );
645 }
646 
647 
648 template <int dim>
loadState(const std::string & filename)649 void UGGrid<dim>::loadState(const std::string& filename)
650 {
651   const char* type = "asc";
652   std::string problemName = name_ + "_Problem";
653   std::string formatName = "DuneFormat2d";
654 
655   if (dim==2) {
656     std::string formatName = "DuneFormat2d";
657     multigrid_ = (typename UG_NS<dim>::MultiGrid*) UG::D2::LoadMultiGrid(
658       name_.c_str(),
659       filename.c_str(),
660       type,
661       problemName.c_str(),
662       formatName.c_str(),
663       0,    // dummy heap size
664       true, //force,
665       true, //optimizedIO,
666       false //autosave
667 #if ModelP and DUNE_UGGRID_HAVE_PPIFCONTEXT
668       , std::make_shared<PPIF::PPIFContext>(comm())
669 #endif
670       );
671   } else {
672     std::string formatName = "DuneFormat3d";
673     multigrid_ = (typename UG_NS<dim>::MultiGrid*) UG::D3::LoadMultiGrid(
674       name_.c_str(),
675       filename.c_str(),
676       type,
677       problemName.c_str(),
678       formatName.c_str(),
679       0,    // dummy heap size
680       true, //force,
681       true, //optimizedIO,
682       false //autosave
683 #if ModelP and DUNE_UGGRID_HAVE_PPIFCONTEXT
684       , std::make_shared<PPIF::PPIFContext>(comm())
685 #endif
686       );
687   }
688 
689   if (multigrid_==nullptr)
690     DUNE_THROW(GridError, "In loadState()");
691 }
692 
693 template < int dim >
setIndices(bool setLevelZero,std::vector<unsigned int> * nodePermutation)694 void UGGrid < dim >::setIndices(bool setLevelZero,
695                                       std::vector<unsigned int>* nodePermutation)
696 {
697   // Create new level index sets if necessary
698   for (int i=levelIndexSets_.size(); i<=maxLevel(); i++)
699     levelIndexSets_.push_back(std::make_shared<UGGridLevelIndexSet<const UGGrid<dim> > >());
700 
701   // Update the zero level LevelIndexSet.  It is updated only once, at the time
702   // of creation of the coarse grid.  After that it is not touched anymore.
703   if (setLevelZero)
704     levelIndexSets_[0]->update(*this, 0, nodePermutation);
705 
706   // Update the remaining level index sets
707   for (int i=1; i<=maxLevel(); i++)
708     if (levelIndexSets_[i])
709       levelIndexSets_[i]->update(*this, i);
710 
711   leafIndexSet_.update(nodePermutation);
712 
713   // id sets don't need updating
714 }
715 
716 // /////////////////////////////////////////////////////////////////////////////////
717 //   Explicit instantiation of the dimensions that are actually supported by UG.
718 //   g++-4.0 wants them to be _after_ the method implementations.
719 // /////////////////////////////////////////////////////////////////////////////////
720 
721 template class UGGrid<2>;
722 template class UGGrid<3>;
723 
724 } /* namespace Dune */
725