File: solution_transfer.cc

package info (click to toggle)
deal.ii 9.0.1-1
  • links: PTS, VCS
  • area: main
  • in suites: buster
  • size: 181,876 kB
  • sloc: cpp: 265,739; ansic: 52,054; python: 1,507; perl: 645; sh: 506; xml: 437; makefile: 73
file content (590 lines) | stat: -rw-r--r-- 22,566 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
// ---------------------------------------------------------------------
//
// Copyright (C) 1999 - 2017 by the deal.II authors
//
// This file is part of the deal.II library.
//
// The deal.II library is free software; you can use it, redistribute
// it, and/or modify it under the terms of the GNU Lesser General
// Public License as published by the Free Software Foundation; either
// version 2.1 of the License, or (at your option) any later version.
// The full text of the license can be found in the file LICENSE at
// the top level of the deal.II distribution.
//
// ---------------------------------------------------------------------

#include <deal.II/base/memory_consumption.h>
#include <deal.II/grid/tria.h>
#include <deal.II/grid/tria_accessor.h>
#include <deal.II/grid/tria_iterator.h>
#include <deal.II/distributed/tria.h>
#include <deal.II/dofs/dof_handler.h>
#include <deal.II/dofs/dof_accessor.h>
#include <deal.II/fe/fe.h>
#include <deal.II/lac/vector_element_access.h>
#include <deal.II/lac/vector.h>
#include <deal.II/lac/la_vector.h>
#include <deal.II/lac/la_parallel_vector.h>
#include <deal.II/lac/petsc_parallel_vector.h>
#include <deal.II/lac/petsc_parallel_block_vector.h>
#include <deal.II/lac/trilinos_vector.h>
#include <deal.II/lac/block_vector.h>
#include <deal.II/lac/la_parallel_block_vector.h>
#include <deal.II/lac/trilinos_parallel_block_vector.h>
#include <deal.II/numerics/solution_transfer.h>

DEAL_II_NAMESPACE_OPEN

template <int dim, typename VectorType, typename DoFHandlerType>
SolutionTransfer<dim, VectorType, DoFHandlerType>::
SolutionTransfer(const DoFHandlerType &dof)
  :
  dof_handler(&dof, typeid(*this).name()),
  n_dofs_old(0),
  prepared_for(none)
{
  Assert ((dynamic_cast<const parallel::distributed::Triangulation<DoFHandlerType::dimension, DoFHandlerType::space_dimension>*>
           (&dof_handler->get_triangulation())
           == nullptr),
          ExcMessage ("You are calling the dealii::SolutionTransfer class "
                      "with a DoF handler that is built on a "
                      "parallel::distributed::Triangulation. This will not "
                      "work for parallel computations. You probably want to "
                      "use the parallel::distributed::SolutionTransfer class."));
}



template <int dim, typename VectorType, typename DoFHandlerType>
SolutionTransfer<dim, VectorType, DoFHandlerType>::~SolutionTransfer()
{
  clear ();
}



template <int dim, typename VectorType, typename DoFHandlerType>
void SolutionTransfer<dim, VectorType, DoFHandlerType>::clear ()
{
  indices_on_cell.clear();
  dof_values_on_cell.clear();
  cell_map.clear();

  prepared_for=none;
}



template <int dim, typename VectorType, typename DoFHandlerType>
void SolutionTransfer<dim, VectorType, DoFHandlerType>::prepare_for_pure_refinement()
{
  Assert(prepared_for!=pure_refinement, ExcAlreadyPrepForRef());
  Assert(prepared_for!=coarsening_and_refinement,
         ExcAlreadyPrepForCoarseAndRef());

  clear();

  const unsigned int n_active_cells = dof_handler->get_triangulation().n_active_cells();
  n_dofs_old=dof_handler->n_dofs();

  // efficient reallocation of indices_on_cell
  std::vector<std::vector<types::global_dof_index> > (n_active_cells)
  .swap(indices_on_cell);

  typename DoFHandlerType::active_cell_iterator cell = dof_handler->begin_active(),
                                                endc = dof_handler->end();

  for (unsigned int i=0; cell!=endc; ++cell, ++i)
    {
      indices_on_cell[i].resize(cell->get_fe().dofs_per_cell);
      // on each cell store the indices of the
      // dofs. after refining we get the values
      // on the children by taking these
      // indices, getting the respective values
      // out of the data vectors and prolonging
      // them to the children
      cell->get_dof_indices(indices_on_cell[i]);
      cell_map[std::make_pair(cell->level(),cell->index())]
        = Pointerstruct(&indices_on_cell[i], cell->active_fe_index());
    }
  prepared_for=pure_refinement;
}



template <int dim, typename VectorType, typename DoFHandlerType>
void
SolutionTransfer<dim, VectorType, DoFHandlerType>::refine_interpolate
(const VectorType &in,
 VectorType       &out) const
{
  Assert(prepared_for==pure_refinement, ExcNotPrepared());
  Assert(in.size()==n_dofs_old, ExcDimensionMismatch(in.size(),n_dofs_old));
  Assert(out.size()==dof_handler->n_dofs(),
         ExcDimensionMismatch(out.size(),dof_handler->n_dofs()));
  Assert(&in != &out,
         ExcMessage ("Vectors cannot be used as input and output"
                     " at the same time!"));

  Vector<typename VectorType::value_type> local_values(0);

  typename DoFHandlerType::cell_iterator cell = dof_handler->begin(),
                                         endc = dof_handler->end();

  typename std::map<std::pair<unsigned int, unsigned int>, Pointerstruct>::const_iterator
  pointerstruct,
  cell_map_end=cell_map.end();

  for (; cell!=endc; ++cell)
    {
      pointerstruct=cell_map.find(std::make_pair(cell->level(),cell->index()));

      if (pointerstruct!=cell_map_end)
        // this cell was refined or not
        // touched at all, so we can get
        // the new values by just setting
        // or interpolating to the children,
        // which is both done by one
        // function
        {
          const unsigned int this_fe_index = pointerstruct->second.active_fe_index;
          const unsigned int dofs_per_cell=cell->get_dof_handler().get_fe(this_fe_index).dofs_per_cell;
          local_values.reinit(dofs_per_cell, true);

          // make sure that the size of the stored indices is the same as
          // dofs_per_cell. since we store the desired fe_index, we know
          // what this size should be
          Assert(dofs_per_cell==(*pointerstruct->second.indices_ptr).size(),
                 ExcInternalError());
          for (unsigned int i=0; i<dofs_per_cell; ++i)
            local_values(i)=internal::ElementAccess<VectorType>::get(
                              in,(*pointerstruct->second.indices_ptr)[i]);
          cell->set_dof_values_by_interpolation(local_values, out,
                                                this_fe_index);
        }
    }
}



namespace internal
{
  /**
   * Generate a table that contains
   * interpolation matrices between
   * each combination of finite
   * elements used in a DoFHandler of
   * some kind. Since not all
   * elements can be interpolated
   * onto each other, the table may
   * contain empty matrices for those
   * combinations of elements for
   * which no such interpolation is
   * implemented.
   */
  template <typename DoFHandlerType>
  void extract_interpolation_matrices (const DoFHandlerType &,
                                       dealii::Table<2,FullMatrix<double> > &)
  {}

  template <int dim, int spacedim>
  void extract_interpolation_matrices (const dealii::hp::DoFHandler<dim,spacedim> &dof,
                                       dealii::Table<2,FullMatrix<double> > &matrices)
  {
    const dealii::hp::FECollection<dim,spacedim> &fe = dof.get_fe_collection();
    matrices.reinit (fe.size(), fe.size());
    for (unsigned int i=0; i<fe.size(); ++i)
      for (unsigned int j=0; j<fe.size(); ++j)
        if (i != j)
          {
            matrices(i,j).reinit (fe[i].dofs_per_cell, fe[j].dofs_per_cell);

            // see if we can get the interpolation matrices for this
            // combination of elements. if not, reset the matrix sizes to zero
            // to indicate that this particular combination isn't
            // supported. this isn't an outright error right away since we may
            // never need to actually interpolate between these two elements
            // on actual cells; we simply have to trigger an error if someone
            // actually tries
            try
              {
                fe[i].get_interpolation_matrix (fe[j], matrices(i,j));
              }
            catch (const typename FiniteElement<dim,spacedim>::ExcInterpolationNotImplemented &)
              {
                matrices(i,j).reinit (0,0);
              }
          }
  }


  template <int dim, int spacedim>
  void restriction_additive (const FiniteElement<dim,spacedim> &,
                             std::vector<std::vector<bool> > &)
  {}

  template <int dim, int spacedim>
  void restriction_additive (const dealii::hp::FECollection<dim,spacedim> &fe,
                             std::vector<std::vector<bool> > &restriction_is_additive)
  {
    restriction_is_additive.resize (fe.size());
    for (unsigned int f=0; f<fe.size(); ++f)
      {
        restriction_is_additive[f].resize (fe[f].dofs_per_cell);
        for (unsigned int i=0; i<fe[f].dofs_per_cell; ++i)
          restriction_is_additive[f][i] = fe[f].restriction_is_additive(i);
      }
  }
}



template <int dim, typename VectorType, typename DoFHandlerType>
void
SolutionTransfer<dim, VectorType, DoFHandlerType>::
prepare_for_coarsening_and_refinement(const std::vector<VectorType> &all_in)
{
  Assert (prepared_for!=pure_refinement, ExcAlreadyPrepForRef());
  Assert (prepared_for!=coarsening_and_refinement,
          ExcAlreadyPrepForCoarseAndRef());

  const unsigned int in_size=all_in.size();
  Assert(in_size!=0,
         ExcMessage("The array of input vectors you pass to this "
                    "function has no elements. This is not useful."));

  clear();

  const unsigned int n_active_cells = dof_handler->get_triangulation().n_active_cells();
  (void)n_active_cells;
  n_dofs_old = dof_handler->n_dofs();

  for (unsigned int i=0; i<in_size; ++i)
    {
      Assert(all_in[i].size()==n_dofs_old,
             ExcDimensionMismatch(all_in[i].size(),n_dofs_old));
    }

  // first count the number
  // of cells that will be coarsened
  // and that'll stay or be refined
  unsigned int n_cells_to_coarsen=0;
  unsigned int n_cells_to_stay_or_refine=0;
  for (typename DoFHandlerType::active_cell_iterator act_cell = dof_handler->begin_active();
       act_cell!=dof_handler->end(); ++act_cell)
    {
      if (act_cell->coarsen_flag_set())
        ++n_cells_to_coarsen;
      else
        ++n_cells_to_stay_or_refine;
    }
  Assert((n_cells_to_coarsen+n_cells_to_stay_or_refine)==n_active_cells,
         ExcInternalError());

  unsigned int n_coarsen_fathers=0;
  for (typename DoFHandlerType::cell_iterator cell=dof_handler->begin();
       cell!=dof_handler->end(); ++cell)
    if (!cell->active() && cell->child(0)->coarsen_flag_set())
      ++n_coarsen_fathers;
  Assert(n_cells_to_coarsen>=2*n_coarsen_fathers, ExcInternalError());

  // allocate the needed memory. initialize
  // the following arrays in an efficient
  // way, without copying much
  std::vector<std::vector<types::global_dof_index> >(n_cells_to_stay_or_refine)
  .swap(indices_on_cell);

  std::vector<std::vector<Vector<typename VectorType::value_type> > >
  (n_coarsen_fathers,
   std::vector<Vector<typename VectorType::value_type> > (in_size))
  .swap(dof_values_on_cell);

  Table<2,FullMatrix<double> > interpolation_hp;
  std::vector<std::vector<bool> > restriction_is_additive;

  internal::extract_interpolation_matrices (*dof_handler, interpolation_hp);
  internal::restriction_additive (dof_handler->get_fe_collection(), restriction_is_additive);

  // we need counters for
  // the 'to_stay_or_refine' cells 'n_sr' and
  // the 'coarsen_fathers' cells 'n_cf',
  unsigned int n_sr=0, n_cf=0;
  for (typename DoFHandlerType::cell_iterator cell=dof_handler->begin();
       cell!=dof_handler->end(); ++cell)
    {
      // CASE 1: active cell that remains as it is
      if (cell->active() && !cell->coarsen_flag_set())
        {
          const unsigned int dofs_per_cell=cell->get_fe().dofs_per_cell;
          indices_on_cell[n_sr].resize(dofs_per_cell);
          // cell will not be coarsened,
          // so we get away by storing the
          // dof indices and later
          // interpolating to the children
          cell->get_dof_indices(indices_on_cell[n_sr]);
          cell_map[std::make_pair(cell->level(), cell->index())]
            = Pointerstruct(&indices_on_cell[n_sr], cell->active_fe_index());
          ++n_sr;
        }

      // CASE 2: cell is inactive but will become active
      else if (cell->has_children() && cell->child(0)->coarsen_flag_set())
        {
          // note that if one child has the
          // coarsen flag, then all should
          // have if Tria::prepare_* has
          // worked correctly
          for (unsigned int i=1; i<cell->n_children(); ++i)
            Assert(cell->child(i)->coarsen_flag_set(),
                   ExcMessage("It looks like you didn't call "
                              "Triangulation::prepare_coarsening_and_refinement before "
                              "calling the current function. This can't work."));

          // we will need to interpolate from the children of this cell
          // to the current one. in the hp context, this also means
          // we need to figure out which finite element space to interpolate
          // to since that is not implied by the global FE as in the non-hp
          // case.
          bool different_fe_on_children = false;
          for (unsigned int child=1; child<cell->n_children(); ++child)
            if (cell->child(child)->active_fe_index()
                != cell->child(0)->active_fe_index())
              {
                different_fe_on_children = true;
                break;
              }

          // take FE index from the child with most
          // degrees of freedom locally
          unsigned int most_general_child = 0;
          if (different_fe_on_children == true)
            for (unsigned int child=1; child<cell->n_children(); ++child)
              if (cell->child(child)->get_fe().dofs_per_cell >
                  cell->child(most_general_child)->get_fe().dofs_per_cell)
                most_general_child = child;
          const unsigned int target_fe_index = cell->child(most_general_child)->active_fe_index();

          const unsigned int dofs_per_cell=cell->get_dof_handler().get_fe(target_fe_index).dofs_per_cell;

          std::vector<Vector<typename VectorType::value_type> >(in_size,
                                                                Vector<typename VectorType::value_type>(dofs_per_cell))
          .swap(dof_values_on_cell[n_cf]);


          // store the data of each of the input vectors. get this data
          // as interpolated onto a finite element space that encompasses
          // that of all the children. note that cell->get_interpolated_dof_values
          // already does all of the interpolations between spaces
          for (unsigned int j=0; j<in_size; ++j)
            cell->get_interpolated_dof_values(all_in[j],
                                              dof_values_on_cell[n_cf][j],
                                              target_fe_index);
          cell_map[std::make_pair(cell->level(), cell->index())]
            = Pointerstruct(&dof_values_on_cell[n_cf], target_fe_index);
          ++n_cf;
        }
    }
  Assert(n_sr==n_cells_to_stay_or_refine, ExcInternalError());
  Assert(n_cf==n_coarsen_fathers, ExcInternalError());

  prepared_for=coarsening_and_refinement;
}



template <int dim, typename VectorType, typename DoFHandlerType>
void
SolutionTransfer<dim, VectorType, DoFHandlerType>::prepare_for_coarsening_and_refinement
(const VectorType &in)
{
  std::vector<VectorType> all_in=std::vector<VectorType>(1, in);
  prepare_for_coarsening_and_refinement(all_in);
}



template <int dim, typename VectorType, typename DoFHandlerType>
void SolutionTransfer<dim, VectorType, DoFHandlerType>::
interpolate (const std::vector<VectorType> &all_in,
             std::vector<VectorType>       &all_out) const
{
  Assert(prepared_for==coarsening_and_refinement, ExcNotPrepared());
  const unsigned int size=all_in.size();
  Assert(all_out.size()==size, ExcDimensionMismatch(all_out.size(), size));
  for (unsigned int i=0; i<size; ++i)
    Assert (all_in[i].size() == n_dofs_old,
            ExcDimensionMismatch(all_in[i].size(), n_dofs_old));
  for (unsigned int i=0; i<all_out.size(); ++i)
    Assert (all_out[i].size() == dof_handler->n_dofs(),
            ExcDimensionMismatch(all_out[i].size(), dof_handler->n_dofs()));
  for (unsigned int i=0; i<size; ++i)
    for (unsigned int j=0; j<size; ++j)
      Assert(&all_in[i] != &all_out[j],
             ExcMessage ("Vectors cannot be used as input and output"
                         " at the same time!"));

  Vector<typename VectorType::value_type> local_values;
  std::vector<types::global_dof_index> dofs;

  typename std::map<std::pair<unsigned int, unsigned int>, Pointerstruct>::const_iterator
  pointerstruct,
  cell_map_end=cell_map.end();

  Table<2,FullMatrix<double> > interpolation_hp;
  internal::extract_interpolation_matrices (*dof_handler, interpolation_hp);
  Vector<typename VectorType::value_type> tmp, tmp2;

  typename DoFHandlerType::cell_iterator cell = dof_handler->begin(),
                                         endc = dof_handler->end();
  for (; cell!=endc; ++cell)
    {
      pointerstruct=cell_map.find(std::make_pair(cell->level(),cell->index()));

      if (pointerstruct!=cell_map_end)
        {
          const std::vector<types::global_dof_index> *const indexptr
            =pointerstruct->second.indices_ptr;

          const std::vector<Vector<typename VectorType::value_type> > *const valuesptr
            =pointerstruct->second.dof_values_ptr;

          // cell stayed as it was or was refined
          if (indexptr)
            {
              Assert (valuesptr == nullptr,
                      ExcInternalError());

              const unsigned int old_fe_index =
                pointerstruct->second.active_fe_index;

              // get the values of
              // each of the input
              // data vectors on this
              // cell and prolong it
              // to its children
              unsigned int in_size = indexptr->size();
              for (unsigned int j=0; j<size; ++j)
                {
                  tmp.reinit (in_size, true);
                  for (unsigned int i=0; i<in_size; ++i)
                    tmp(i) = internal::ElementAccess<VectorType>::get(all_in[j],
                                                                      (*indexptr)[i]);

                  cell->set_dof_values_by_interpolation (tmp, all_out[j],
                                                         old_fe_index);
                }
            }
          else if (valuesptr)
            // the children of this cell were
            // deleted
            {
              Assert (!cell->has_children(), ExcInternalError());
              Assert (indexptr == nullptr,
                      ExcInternalError());

              const unsigned int dofs_per_cell = cell->get_fe().dofs_per_cell;
              dofs.resize(dofs_per_cell);
              // get the local
              // indices
              cell->get_dof_indices(dofs);

              // distribute the
              // stored data to the
              // new vectors
              for (unsigned int j=0; j<size; ++j)
                {
                  // make sure that the size of
                  // the stored indices is the
                  // same as
                  // dofs_per_cell. this is
                  // kind of a test if we use
                  // the same fe in the hp
                  // case. to really do that
                  // test we would have to
                  // store the fe_index of all
                  // cells
                  const Vector<typename VectorType::value_type> *data = nullptr;
                  const unsigned int active_fe_index = cell->active_fe_index();
                  if (active_fe_index != pointerstruct->second.active_fe_index)
                    {
                      const unsigned int old_index = pointerstruct->second.active_fe_index;
                      tmp.reinit (dofs_per_cell, true);
                      AssertDimension ((*valuesptr)[j].size(),
                                       interpolation_hp(active_fe_index,old_index).n());
                      AssertDimension (tmp.size(),
                                       interpolation_hp(active_fe_index,old_index).m());
                      interpolation_hp(active_fe_index,old_index).vmult (tmp, (*valuesptr)[j]);
                      data = &tmp;
                    }
                  else
                    data = &(*valuesptr)[j];


                  for (unsigned int i=0; i<dofs_per_cell; ++i)
                    internal::ElementAccess<VectorType>::set((*data)(i), dofs[i],
                                                             all_out[j]);
                }
            }
          // undefined status
          else
            Assert(false, ExcInternalError());
        }
    }
}



template <int dim, typename VectorType, typename DoFHandlerType>
void SolutionTransfer<dim, VectorType, DoFHandlerType>::interpolate
(const VectorType &in,
 VectorType       &out) const
{
  Assert (in.size()==n_dofs_old,
          ExcDimensionMismatch(in.size(), n_dofs_old));
  Assert (out.size()==dof_handler->n_dofs(),
          ExcDimensionMismatch(out.size(), dof_handler->n_dofs()));

  std::vector<VectorType> all_in(1);
  all_in[0] = in;
  std::vector<VectorType> all_out(1);
  all_out[0] = out;
  interpolate(all_in,
              all_out);
  out=all_out[0];
}



template <int dim, typename VectorType, typename DoFHandlerType>
std::size_t
SolutionTransfer<dim, VectorType, DoFHandlerType>::memory_consumption () const
{
  // at the moment we do not include the memory
  // consumption of the cell_map as we have no
  // real idea about memory consumption of a
  // std::map
  return (MemoryConsumption::memory_consumption (dof_handler) +
          MemoryConsumption::memory_consumption (n_dofs_old) +
          sizeof (prepared_for) +
          MemoryConsumption::memory_consumption (indices_on_cell) +
          MemoryConsumption::memory_consumption (dof_values_on_cell));
}



template <int dim, typename VectorType, typename DoFHandlerType>
std::size_t
SolutionTransfer<dim, VectorType, DoFHandlerType>::Pointerstruct::memory_consumption () const
{
  return sizeof(*this);
}


/*-------------- Explicit Instantiations -------------------------------*/
#define SPLIT_INSTANTIATIONS_COUNT 4
#ifndef SPLIT_INSTANTIATIONS_INDEX
#define SPLIT_INSTANTIATIONS_INDEX 0
#endif
#include "solution_transfer.inst"

DEAL_II_NAMESPACE_CLOSE