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 (261 lines) | stat: -rw-r--r-- 9,412 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
// ---------------------------------------------------------------------
//
// Copyright (C) 2009 - 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/config.h>

#ifdef DEAL_II_WITH_P4EST

#include <deal.II/lac/vector.h>
#include <deal.II/lac/block_vector.h>
#include <deal.II/lac/la_parallel_vector.h>
#include <deal.II/lac/la_parallel_block_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/trilinos_parallel_block_vector.h>

#include <deal.II/distributed/solution_transfer.h>
#include <deal.II/distributed/tria.h>
#include <deal.II/dofs/dof_tools.h>
#include <deal.II/dofs/dof_accessor.h>
#include <deal.II/grid/tria_accessor.h>
#include <deal.II/grid/tria_iterator.h>

#include <functional>

DEAL_II_NAMESPACE_OPEN

namespace parallel
{
  namespace distributed
  {

    template <int dim, typename VectorType, typename DoFHandlerType>
    SolutionTransfer<dim, VectorType, DoFHandlerType>::SolutionTransfer (const DoFHandlerType &dof)
      :
      dof_handler(&dof, typeid(*this).name()),
      offset (numbers::invalid_unsigned_int)
    {
      Assert ((dynamic_cast<const parallel::distributed::Triangulation<dim,DoFHandlerType::space_dimension>*>
               (&dof_handler->get_triangulation()) != nullptr),
              ExcMessage("parallel::distributed::SolutionTransfer requires a parallel::distributed::Triangulation object."));
    }



    template <int dim, typename VectorType, typename DoFHandlerType>
    void
    SolutionTransfer<dim, VectorType, DoFHandlerType>::
    prepare_for_coarsening_and_refinement (const std::vector<const VectorType *> &all_in)
    {
      input_vectors = all_in;
      register_data_attach( get_data_size() * input_vectors.size() );
    }



    template <int dim, typename VectorType, typename DoFHandlerType>
    void
    SolutionTransfer<dim, VectorType, DoFHandlerType>::register_data_attach (const std::size_t size)
    {
      Assert(size > 0, ExcMessage("Please transfer at least one vector!"));

//TODO: casting away constness is bad
      parallel::distributed::Triangulation<dim,DoFHandlerType::space_dimension> *tria
        = (dynamic_cast<parallel::distributed::Triangulation<dim,DoFHandlerType::space_dimension>*>
           (const_cast<dealii::Triangulation<dim,DoFHandlerType::space_dimension>*>
            (&dof_handler->get_triangulation())));
      Assert (tria != nullptr, ExcInternalError());

      offset
        = tria->register_data_attach(size,
                                     std::bind(&SolutionTransfer<dim, VectorType,
                                               DoFHandlerType>::pack_callback,
                                               this,
                                               std::placeholders::_1,
                                               std::placeholders::_2,
                                               std::placeholders::_3));

    }



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



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



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



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



    template <int dim, typename VectorType, typename DoFHandlerType>
    void
    SolutionTransfer<dim, VectorType, DoFHandlerType>::deserialize (std::vector<VectorType *> &all_in)
    {
      register_data_attach( get_data_size() * all_in.size() );

      // this makes interpolate() happy
      input_vectors.resize(all_in.size());

      interpolate(all_in);
    }


    template <int dim, typename VectorType, typename DoFHandlerType>
    void
    SolutionTransfer<dim, VectorType, DoFHandlerType>::interpolate (std::vector<VectorType *> &all_out)
    {
      Assert(input_vectors.size()==all_out.size(),
             ExcDimensionMismatch(input_vectors.size(), all_out.size()) );

//TODO: casting away constness is bad
      parallel::distributed::Triangulation<dim,DoFHandlerType::space_dimension> *tria
        = (dynamic_cast<parallel::distributed::Triangulation<dim,DoFHandlerType::space_dimension>*>
           (const_cast<dealii::Triangulation<dim,DoFHandlerType::space_dimension>*>
            (&dof_handler->get_triangulation())));
      Assert (tria != nullptr, ExcInternalError());

      tria->notify_ready_to_unpack(offset,
                                   std::bind(&SolutionTransfer<dim, VectorType,
                                             DoFHandlerType>::unpack_callback,
                                             this,
                                             std::placeholders::_1,
                                             std::placeholders::_2,
                                             std::placeholders::_3,
                                             std::ref(all_out)));


      for (typename std::vector<VectorType *>::iterator it=all_out.begin();
           it !=all_out.end();
           ++it)
        (*it)->compress(::dealii::VectorOperation::insert);

      input_vectors.clear();
    }



    template <int dim, typename VectorType, typename DoFHandlerType>
    void
    SolutionTransfer<dim, VectorType, DoFHandlerType>::interpolate (VectorType &out)
    {
      std::vector<VectorType *> all_out(1, &out);
      interpolate(all_out);
    }



    template <int dim, typename VectorType, typename DoFHandlerType>
    unsigned int
    SolutionTransfer<dim, VectorType, DoFHandlerType>::get_data_size() const
    {
      return sizeof(typename VectorType::value_type)* DoFTools::max_dofs_per_cell(*dof_handler);
    }


    template <int dim, typename VectorType, typename DoFHandlerType>
    void
    SolutionTransfer<dim, VectorType, DoFHandlerType>::
    pack_callback(const typename Triangulation<dim,DoFHandlerType::space_dimension>::cell_iterator &cell_,
                  const typename Triangulation<dim,DoFHandlerType::space_dimension>::CellStatus /*status*/,
                  void *data)
    {
      typename VectorType::value_type *data_store = reinterpret_cast<typename VectorType::value_type *>(data);

      typename DoFHandlerType::cell_iterator cell(*cell_, dof_handler);

      const unsigned int dofs_per_cell=cell->get_fe().dofs_per_cell;
      ::dealii::Vector<typename VectorType::value_type> dofvalues(dofs_per_cell);
      for (typename std::vector<const VectorType *>::iterator it=input_vectors.begin();
           it !=input_vectors.end();
           ++it)
        {
          cell->get_interpolated_dof_values(*(*it), dofvalues);
          std::memcpy(data_store, &dofvalues(0), sizeof(typename VectorType::value_type)*dofs_per_cell);
          data_store += dofs_per_cell;
        }
    }


    template <int dim, typename VectorType, typename DoFHandlerType>
    void
    SolutionTransfer<dim, VectorType, DoFHandlerType>::unpack_callback
    (const typename Triangulation<dim,DoFHandlerType::space_dimension>::cell_iterator &cell_,
     const typename Triangulation<dim,DoFHandlerType::space_dimension>::CellStatus    /*status*/,
     const void                                           *data,
     std::vector<VectorType *>                            &all_out)
    {
      typename DoFHandlerType::cell_iterator
      cell(*cell_, dof_handler);

      const unsigned int dofs_per_cell=cell->get_fe().dofs_per_cell;
      ::dealii::Vector<typename VectorType::value_type> dofvalues(dofs_per_cell);
      const typename VectorType::value_type *data_store = reinterpret_cast<const typename VectorType::value_type *>(data);

      for (typename std::vector<VectorType *>::iterator it = all_out.begin();
           it != all_out.end();
           ++it)
        {
          std::memcpy(&dofvalues(0), data_store, sizeof(typename VectorType::value_type)*dofs_per_cell);
          cell->set_dof_values_by_interpolation(dofvalues, *(*it));
          data_store += dofs_per_cell;
        }
    }


  }
}


// explicit instantiations
#include "solution_transfer.inst"

DEAL_II_NAMESPACE_CLOSE

#endif