DOFMappingTest.cpp 3.64 KB
Newer Older
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
#include <amdis/AMDiS.hpp>

#include <cassert>
#include <iostream>
#include <sstream>

#include <dune/common/parallel/mpihelper.hh>
#include <dune/common/parallel/indexset.hh>
#include <dune/common/parallel/plocalindex.hh>
#include <dune/common/parallel/remoteindices.hh>

#include <amdis/linearalgebra/DOFMapping.hpp>

#if HAVE_PETSC
#include <petscvec.h>
#endif

#if HAVE_PMTL
#include <boost/numeric/mtl/mtl.hpp>
#endif

namespace Attribute {
  enum Type { owner = 1, overlap = 2, copy = 3 };
}

struct GlobalIndex : std::pair<int,int>
{
  GlobalIndex(int i = 0, int j = 0) : std::pair<int,int>(i,j) {}
};

struct LexicographicOrdering
{
  LexicographicOrdering(int Nx, int Ny, int shiftx, int shifty)
    : Nx_(Nx)
    , Ny_(Ny)
    , shiftx_(shiftx)
    , shifty_(shifty)
  {}

  int operator()(int i, int j) const
  {
    return Nx_*(i - shifty_) + (j - shiftx_);
  }

  int Nx_,Ny_,shiftx_,shifty_;
};


std::ostream& operator<<(std::ostream& out, GlobalIndex const& p)
{
  out << "(" << p.first << ", " << p.second << ")";
  return out;
}

template <class PIS>
struct Communication
{
  using ParallelIndexSet = PIS;
  using RemoteIndices = Dune::RemoteIndices<PIS>;

  Communication()
    : pis_()
    , rIndices_(pis_, pis_, Dune::MPIHelper::getCommunicator())
  {}

  ParallelIndexSet& indexSet() { return pis_; }
  ParallelIndexSet const& indexSet() const { return pis_; }

  RemoteIndices& remoteIndices() { return rIndices_; }
  RemoteIndices const& remoteIndices() const { return rIndices_; }

  auto communicator() const { return Dune::MPIHelper::getCommunicator(); }

private:
  ParallelIndexSet pis_;
  RemoteIndices rIndices_;
};

int main(int argc, char** argv)
{
  AMDiS::Environment env(argc, argv);
  MPI_Comm comm = MPI_COMM_WORLD;

  using GI = GlobalIndex;
  using LI = Dune::ParallelLocalIndex<Attribute::Type>;
  using PIS = Dune::ParallelIndexSet<GI,LI>;

  int r = env.mpiRank();
  int s = env.mpiSize();
  assert( s == 2 );

  int N = 3;
  int M = 2;

  auto index = LexicographicOrdering{N+1, M+1, 0, r};

  Communication<PIS> c;

  PIS& pis = c.indexSet();
  pis.beginResize();
  int i0 = r;
  int i1 = i0 + M + 1;
  for (int i = i0; i < i1; ++i) {
    for (int j = 0; j < N+1; ++j) {
      bool overlap_region = (i >= 1 && i <= 2);
      pis.add(GlobalIndex{i,j}, LI(index(i,j), !overlap_region || r == 0 ? Attribute::owner : Attribute::overlap, true));
    }
  }
  pis.endResize();

  auto& remoteIndices = c.remoteIndices();
  remoteIndices.template rebuild<true>();

  AMDiS::DofMapping<PIS,std::size_t> mapping(c);
  for (int i = i0; i < i1; ++i) {
    for (int j = 0; j < N+1; ++j) {
      std::cout << "[" << r << "]    (" << i << "," << j << ") = " << index(i,j) << " => " << mapping.global(index(i,j)) << "\n";
    }
  }

#if HAVE_PETSC
  Vec v;
  VecCreateMPI(comm, mapping.localSize(), mapping.globalSize(), &v);
  VecSet(v,0.0);

  for (int i = i0; i < i1; ++i) {
    for (int j = 0; j < N+1; ++j) {
      PetscInt row = mapping.global(index(i,j));
      VecSetValue(v, row, PetscScalar(index(i,j)), ADD_VALUES);
    }
  }
  VecAssemblyBegin(v);
  VecAssemblyEnd(v);

  VecView(v,PETSC_VIEWER_STDOUT_WORLD);
  VecDestroy(&v);
#endif

#if HAVE_PMTL
  mtl::par::block_distribution dist(mapping.globalStarts(), c.communicator());

  using vector_type = mtl::vector::distributed<mtl::dense_vector<double> >;
  vector_type vec(mapping.globalSize(), dist);

  { mtl::vector::inserter<vector_type, mtl::update_plus<double>> ins(vec);

    for (int i = i0; i < i1; ++i) {
      for (int j = 0; j < N+1; ++j) {
        PetscInt row = mapping.global(index(i,j));
        ins[row] << index(i,j);
      }
    }
  }
  mtl::par::rout << "My local part of v is " << local(vec) << '\n';
#endif
}