DOLFIN
DOLFIN C++ interface
PETScDMCollection.h
1 // Copyright (C) 2016 Patrick E. Farrell and Garth N. Wells
2 //
3 // This file is part of DOLFIN.
4 //
5 // DOLFIN is free software: you can redistribute it and/or modify
6 // it under the terms of the GNU Lesser General Public License as published by
7 // the Free Software Foundation, either version 3 of the License, or
8 // (at your option) any later version.
9 //
10 // DOLFIN is distributed in the hope that it will be useful,
11 // but WITHOUT ANY WARRANTY; without even the implied warranty of
12 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 // GNU Lesser General Public License for more details.
14 //
15 // You should have received a copy of the GNU Lesser General Public License
16 // along with DOLFIN. If not, see <http://www.gnu.org/licenses/>.
17 
18 #ifndef __PETSC_DM_COLLECTION_H
19 #define __PETSC_DM_COLLECTION_H
20 
21 #ifdef HAS_PETSC
22 
23 #include <memory>
24 #include <vector>
25 #include <petscdm.h>
26 #include <petscvec.h>
27 #include <dolfin/common/MPI.h>
28 #include <dolfin/la/PETScMatrix.h>
29 #include <dolfin/log/log.h>
30 
31 namespace dolfin
32 {
33 
34  class FunctionSpace;
35  class BoundingBoxTree;
36 
42 
44  {
45  public:
46 
50  PETScDMCollection(std::vector<std::shared_ptr<const FunctionSpace>> function_spaces);
51 
54 
58  DM get_dm(int i);
59 
61  void check_ref_count() const;
62 
64  void reset(int i);
65 
68  static std::shared_ptr<PETScMatrix>
69  create_transfer_matrix(const FunctionSpace& coarse_space,
70  const FunctionSpace& fine_space);
71 
72  private:
73 
74  // Find the nearest cells to points which lie outside the domain
75  static void find_exterior_points(MPI_Comm mpi_comm,
76  std::shared_ptr<const BoundingBoxTree> treec,
77  int dim, int data_size,
78  const std::vector<double>& send_points,
79  const std::vector<int>& send_indices,
80  std::vector<int>& indices,
81  std::vector<std::size_t>& cell_ids,
82  std::vector<double>& points);
83 
84 
85  // Pointers to functions that are used in PETSc DM call-backs
86  static PetscErrorCode create_global_vector(DM dm, Vec* vec);
87  static PetscErrorCode create_interpolation(DM dmc, DM dmf, Mat *mat,
88  Vec *vec);
89  static PetscErrorCode coarsen(DM dmf, MPI_Comm comm, DM* dmc);
90  static PetscErrorCode refine(DM dmc, MPI_Comm comm, DM* dmf);
91 
92 
93  // The FunctionSpaces associated with each level, starting with
94  // the coarest space
95  std::vector<std::shared_ptr<const FunctionSpace>> _spaces;
96 
97  // The PETSc DM objects
98  std::vector<DM> _dms;
99 
100  };
101 
102 }
103 
104 #endif
105 
106 #endif
static std::shared_ptr< PETScMatrix > create_transfer_matrix(const FunctionSpace &coarse_space, const FunctionSpace &fine_space)
Definition: PETScDMCollection.cpp:212
DM get_dm(int i)
Definition: PETScDMCollection.cpp:186
Definition: adapt.h:29
Definition: FunctionSpace.h:53
Definition: PETScDMCollection.h:43
void reset(int i)
Debugging use - to be removed.
Definition: PETScDMCollection.cpp:203
~PETScDMCollection()
Destructor.
Definition: PETScDMCollection.cpp:176
Definition: PETScObject.h:33
PETScDMCollection(std::vector< std::shared_ptr< const FunctionSpace >> function_spaces)
Definition: PETScDMCollection.cpp:135
void check_ref_count() const
These are test/debugging functions that will be removed.
Definition: PETScDMCollection.cpp:193