Epetra Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
test/Bugs_LL/Bug_6079_DistObject_CombineMode_flags_LL/cxx_main.cpp
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 // $Id$
3 //
4 // Copyright (C) 2014 by the deal.II authors
5 //
6 // This file is part of the deal.II library.
7 //
8 // The deal.II library is free software; you can use it, redistribute
9 // it, and/or modify it under the terms of the GNU Lesser General
10 // Public License as published by the Free Software Foundation; either
11 // version 2.1 of the License, or (at your option) any later version.
12 // The full text of the license can be found in the file LICENSE at
13 // the top level of the deal.II distribution.
14 //
15 // ---------------------------------------------------------------------
16 
17 
18 
19 // check correct behavior of sadd of Trilinos vectors
20 // if they have different Epetra maps
21 
22 #include <iostream>
23 #include <vector>
24 #include <Epetra_Import.h>
25 #include <Epetra_Map.h>
26 #include <Epetra_FEVector.h>
27 #ifdef HAVE_MPI
28 #include <Epetra_MpiComm.h>
29 #include <mpi.h>
30 #else
31 #include <Epetra_SerialComm.h>
32 #endif
33 
34 void test ()
35 {
36  using std::cout;
37  using std::endl;
38 
39  int n_proc;
40 #ifdef HAVE_MPI
41  MPI_Comm_size(MPI_COMM_WORLD, &n_proc);
42 #else
43  n_proc = 1;
44 #endif
45  int my_id;
46 #ifdef HAVE_MPI
47  MPI_Comm_rank(MPI_COMM_WORLD, &my_id);
48 #else
49  my_id = 0;
50 #endif
51 
52  //All processes should own 10 entries
53  const int entries_per_process = 10;
54 
55  const long long begin_index = ((long long)my_id)*entries_per_process;
56  const long long end_index = ((long long)(my_id+1))*entries_per_process;
57 
58  const long long local_begin = std::max(0LL, begin_index-entries_per_process/2);
59  const long long local_end = entries_per_process*n_proc;
60 
61  //create Epetra maps
62  std::vector<unsigned long long> ghosted_indices;
63  ghosted_indices.reserve(local_end-local_begin);
64  for (long long i = local_begin; i< local_end; ++i)
65  ghosted_indices.push_back(i);
66  Epetra_Map map_ghosted
67  (-1LL,
68  local_end-local_begin,
69  reinterpret_cast<long long*>(&ghosted_indices[0]),
70  0,
71 #ifdef HAVE_MPI
72  Epetra_MpiComm(MPI_COMM_WORLD));
73 #else
75 #endif
76 
77  std::vector<unsigned long long> distributed_indices;
78  distributed_indices.reserve(entries_per_process*n_proc);
79  for (long long i = begin_index; i< end_index; ++i)
80  distributed_indices.push_back(i);
81  Epetra_Map map_distributed
82  (entries_per_process*n_proc,
83  entries_per_process,
84  reinterpret_cast<long long*>(&distributed_indices[0]),
85  0,
86 #ifdef HAVE_MPI
87  Epetra_MpiComm(MPI_COMM_WORLD));
88 #else
90 #endif
91 
92  Epetra_FEVector v_ghosted(map_ghosted);
93  Epetra_FEVector v_distributed(map_distributed);
94 
95  v_distributed.PutScalar(2.);
96  v_ghosted.PutScalar(1.);
97 
98  Epetra_Import data_exchange (v_distributed.Map(), v_ghosted.Map());
99  int ierr = v_distributed.Import(v_ghosted, data_exchange, Epetra_AddLocalAlso);
100 
101  // Do a reduction over processes to ensure that the Import succeeded everywhere.
102  {
103  int lclErr = (ierr == 0) ? 0 : 1;
104  int gblErr = 0;
105  const int count = 1;
106 
107  const Epetra_Comm& comm = map_ghosted.Comm ();
108  (void) comm.MinAll (&lclErr, &gblErr, count);
109  if (comm.MyPID () == 0 && gblErr != 0) {
110  // Match the test's expected failure string.
111  cout << "tests FAILED: Import failed (returned nonzero error code) "
112  "on at least one process." << endl;
113  }
114  }
115 
116  cout << "Distributed:" << endl;
117  for (long long i=begin_index; i<end_index; ++i)
118  {
119  const int trilinos_i = v_distributed.Map().LID(i);
120  double value = v_distributed[0][trilinos_i];
121  cout << "proc " << my_id << " " << i << ": " << value << endl;
122  if (value != 3) {
123  cout << "tests FAILED: value = " << value << endl;
124  }
125  }
126 }
127 
128 
129 
130 int main (int argc, char **argv)
131 {
132 #ifdef HAVE_MPI
133  MPI_Init(&argc, &argv);
134 #endif
135  test();
136 #ifdef HAVE_MPI
137  MPI_Finalize();
138 #endif
139 }
140 
141 
Epetra_Map: A class for partitioning vectors and matrices.
Definition: Epetra_Map.h:119
int PutScalar(double ScalarConstant)
Initialize all values in a multi-vector with constant value.
virtual int MinAll(double *PartialMins, double *GlobalMins, int Count) const =0
Epetra_Comm Global Min function.
Epetra_MpiComm: The Epetra MPI Communication Class.
Epetra_Import: This class builds an import object for efficient importing of off-processor elements...
Definition: Epetra_Import.h:63
virtual int MyPID() const =0
Return my process ID.
Epetra_Comm: The Epetra Communication Abstract Base Class.
Definition: Epetra_Comm.h:73
int LID(int GID) const
Returns local ID of global ID, return -1 if not found on this processor.
const Epetra_Comm & Comm() const
Access function for Epetra_Comm communicator.
Epetra_SerialComm: The Epetra Serial Communication Class.
Epetra Finite-Element Vector.
int main(int argc, char *argv[])
int Import(const Epetra_SrcDistObject &A, const Epetra_Import &Importer, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor=0)
Imports an Epetra_DistObject using the Epetra_Import object.
const Epetra_BlockMap & Map() const
Returns the address of the Epetra_BlockMap for this multi-vector.