1 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2 // vi: set et ts=4 sw=2 sts=2:
3 
4 #ifndef DUNE_UGGRID_FACTORY_HH
5 #define DUNE_UGGRID_FACTORY_HH
6 
7 /** \file
8     \brief The specialization of the generic GridFactory for UGGrid
9     \author Oliver Sander
10  */
11 
12 #include <array>
13 #include <memory>
14 #include <vector>
15 
16 #include <dune/common/fvector.hh>
17 
18 #include <dune/grid/common/boundarysegment.hh>
19 #include <dune/grid/common/gridfactory.hh>
20 #include <dune/grid/uggrid.hh>
21 
22 namespace Dune {
23 
24 
25   /** \brief Specialization of the generic GridFactory for UGGrid
26 
27       \ingroup GridFactory
28 
29      <p>
30      If you want to write a routine that reads a grid from some
31      file into a Dune UGGrid object you have to know how to use the UGGrid
32      grid factory.  In the following we
33      assume that you have a grid in some file format and an
34      empty UGGrid object, created by one of its constructors.
35      Hence, your file reader method signature may look like this:
36      </p>
37 
38      <pre>
39      UGGrid&lt;3&gt;* readMyFileFormat(const std::string& filename)
40      </pre>
41 
42      Now, in order to create a valid UGGrid object do the
43      following steps:
44 
45      <h2> 1) Create a GridFactory Object </h2>
46      <p>
47      Get a new GridFactory object by calling
48      <pre>
49      GridFactory<UGGrid<dim> > factory;
50      </pre>
51 
52      <h2> 2)  Enter the Vertices </h2>
53 
54      <p>
55      Insert the grid vertices by calling
56      </p>
57 
58      <pre>
59      factory.insertVertex(const FieldVector&lt;double,dimworld&gt;& position);
60      </pre>
61 
62      <p>
63      for each vertex.  The order of insertion determines the level- and leaf indices
64      of your level 0 vertices.
65      </p>
66 
67 
68      <h2> 3) Enter the elements </h2>
69 
70      <p>
71      For each element call
72      </p>
73 
74      <pre>
75      factory.insertElement(Dune::GeometryType type, const std::vector&lt;int&gt;& vertices);
76      </pre>
77 
78      <p>
79      The parameters are
80      </p>
81 
82      <ul>
83      <li> <b>type</b> - The element type.  UG supports the types <i>simplex</i> and
84      <i>cube</i> in 2d, and <i>simplex, cube, prism</i>, and <i>pyramid</i> in 3d.
85      <li> <b>vertices</b> - The Ids of the vertices of this element.</li>
86      </ul>
87 
88      <p>
89      The numbering of the vertices of each element is expected to follow the DUNE conventions.
90      Refer to the page on reference elements for the details.
91 
92      <h2> 4) Parametrized Domains </h2>
93 
94      <p>
95      UGGrid supports parametrized domains.  That means that you can provide a
96      smooth description of your grid boundary.  The actual grid will always
97      be piecewise linear; however, as you refine, the grid will approach your
98      prescribed boundary.  You don't have to do this.  If your
99      coarse grid boundary describes your domain well read on at Section 5.
100      </p>
101 
102      <p>
103      In order to create curved boundary segments, for each segment you have to write
104      a class which implements the correct geometry.  These classes are then handed
105      over to the UGGrid object.  Boundary segment implementations must be derived
106      from
107      <pre>
108      template &lt;int dimworld&gt; Dune::BoundarySegment
109      </pre>
110      This is an abstract base class which requires you to overload the method
111 
112      <pre>
113      virtual FieldVector&lt; double, dimworld &gt; operator() (const FieldVector&lt; double, dimworld-1 &gt; &local)
114      </pre>
115 
116      <p>
117      This methods must compute the world coordinates from the local ones on the
118      boundary segment.  Give these classes to your grid factory by calling
119      </p>
120      <pre>
121      factory.insertBoundarySegment(const std::vector&lt;int&gt;& vertices,
122                              const BoundarySegment&lt;dimworld&gt; *boundarySegment);
123      </pre>
124 
125      <p>
126      Control over the allocated objects is taken from you, and the grid object
127      will take care of their destruction.
128      </p>
129 
130      <h2> 5) Finish construction </h2>
131 
132      <p>
133      To finish off the construction of the UGGrid object call
134      </p>
135 
136      <pre>
137      std::unique_ptr<UGGrid<dim> > grid = factory.createGrid();
138      </pre>
139 
140      <p>
141      This time it is you who gets full responsibility for the allocated object.
142      </p>
143 
144      <h2> Loading a Grid on a Parallel Machine </h2>
145      <p>
146      If you're working on a parallel machine, and you want to set up a
147      parallel grid, proceed as described only on the rank-0 process.
148      On the other processes just create a GridFactory and call createGrid()
149      to obtain the grid object. This will create the grid on the master process
150      and set up UG correctly on all other process.  Call <tt>loadBalance()</tt>
151      to actually distribute the grid.
152      </p>
153 
154      <p>\warning To use a parametrized boundary on a parallel machine you need
155      to hand over the boundary segments to the grid factory on <b>all</b> processes.
156      This behavior violates the Dune grid interface specification and will be
157      corrected in the future.
158      </p>
159    */
160   template <int dimworld>
161   class GridFactory<UGGrid<dimworld> > : public GridFactoryInterface<UGGrid<dimworld> > {
162 
163     /** \brief Type used by the grid for coordinates */
164     typedef typename UGGrid<dimworld>::ctype ctype;
165 
166     // UGGrid only in 2d and 3d
167     static_assert(dimworld==2 || dimworld == 3, "UGGrid only in 2d and 3d");
168 
169   public:
170 
171     /** \brief Default constructor */
172     GridFactory();
173 
174     /** \brief Constructor for a given grid object
175 
176        If you already have your grid object constructed you can
177        hand it over using this constructor.
178 
179        If you construct your factory class using this constructor
180        the pointer handed over to you by the method createGrid() is
181        the one you supplied here.
182      */
183     GridFactory(UGGrid<dimworld>* grid);
184 
185     /** \brief Destructor */
186     ~GridFactory();
187 
188     /** \brief Insert a vertex into the coarse grid */
189     virtual void insertVertex(const FieldVector<ctype,dimworld>& pos);
190 
191     /** \brief Insert an element into the coarse grid
192         \param type The GeometryType of the new element
193         \param vertices The vertices of the new element, using the DUNE numbering
194      */
195     virtual void insertElement(const GeometryType& type,
196                                const std::vector<unsigned int>& vertices);
197 
198     /** \brief Method to insert a boundary segment into a coarse grid
199 
200        Using this method is optional.  It only influences the ordering of the segments
201 
202         \param vertices The indices of the vertices of the segment
203      */
204     void insertBoundarySegment(const std::vector<unsigned int>& vertices);
205 
206     /** \brief Method to insert an arbitrarily shaped boundary segment into a coarse grid
207         \param vertices The indices of the vertices of the segment
208         \param boundarySegment Class implementing the geometry of the boundary segment.
209      */
210     void insertBoundarySegment(const std::vector<unsigned int>& vertices,
211                                const std::shared_ptr<BoundarySegment<dimworld> > &boundarySegment);
212 
213 
214     /** \brief Finalize grid creation and hand over the grid
215 
216        The receiver takes responsibility of the memory allocated for the grid
217      */
218     virtual std::unique_ptr<UGGrid<dimworld>> createGrid();
219 
220     static const int dimension = UGGrid<dimworld>::dimension;
221 
222     template< int codim >
223     struct Codim
224     {
225       typedef typename UGGrid<dimworld>::template Codim< codim >::Entity Entity;
226     };
227 
228     /** \brief Return the number of the element in the order of insertion into the factory
229      *
230      * For UGGrid elements this number is the same as the element level index
231      */
232     virtual unsigned int
insertionIndex(const typename Codim<0>::Entity & entity) const233     insertionIndex ( const typename Codim< 0 >::Entity &entity ) const
234     {
235       return UG_NS<dimension>::levelIndex(entity.impl().target_);
236     }
237 
238     /** \brief Return the number of the vertex in the order of insertion into the factory
239      *
240      * For UGGrid vertices this number is the same as the vertex level index
241      */
242     virtual unsigned int
insertionIndex(const typename Codim<dimension>::Entity & entity) const243     insertionIndex ( const typename Codim< dimension >::Entity &entity ) const
244     {
245       return UG_NS<dimension>::levelIndex(entity.impl().target_);
246     }
247 
248     /** \brief Return the number of the intersection in the order of insertion into the factory
249      *
250      * For UGGrid intersections this number is the same as the boundary segment index
251      */
252     virtual unsigned int
insertionIndex(const typename UGGrid<dimworld>::LeafIntersection & intersection) const253     insertionIndex ( const typename UGGrid<dimworld>::LeafIntersection &intersection ) const
254     {
255       return intersection.boundarySegmentIndex();
256     }
257 
258     /** \brief Return true if the intersection has been explictily insterted into the factory */
259     virtual bool
wasInserted(const typename UGGrid<dimworld>::LeafIntersection & intersection) const260     wasInserted ( const typename UGGrid<dimworld>::LeafIntersection &intersection ) const
261     {
262       return (insertionIndex( intersection ) < boundarySegmentVertices_.size());
263     }
264 
265     using Communication = typename UGGrid<dimworld>::CollectiveCommunication;
266 
267     /** \brief Return the Communication used by the grid factory
268      *
269      * Use the Communication available from the grid.
270      */
comm() const271     Communication comm() const
272     {
273         return grid_->comm();
274     }
275 
276   private:
277 
278     // Initialize the grid structure in UG
279     void createBegin();
280 
281     // Pointer to the grid being built
282     UGGrid<dimworld>* grid_;
283 
284     // True if the factory allocated the grid itself, false if the
285     // grid was handed over from the outside
286     bool factoryOwnsGrid_;
287 
288     /** \brief Buffer for the vertices of each explicitly given boundary segment */
289     std::vector<std::array<int, dimworld*2-2> > boundarySegmentVertices_;
290 
291     /** \brief While inserting the elements this array records the number of
292         vertices of each element. */
293     std::vector<unsigned char> elementTypes_;
294 
295     /** \brief While inserting the elements this array records the vertices
296         of the elements. */
297     std::vector<unsigned int> elementVertices_;
298 
299     /** \brief Buffer the vertices until createend() is called */
300     std::vector<FieldVector<double, dimworld> > vertexPositions_;
301 
302   };
303 
304 }
305 
306 #endif
307