[Trilinos-Users] [EXTERNAL] Re: Ifpack2: GlobalMPISession Error

Hoemmen, Mark mhoemme at sandia.gov
Mon Aug 1 00:49:55 EDT 2016


See Tpetra tutorial:

https://trilinos.org/docs/dev/packages/tpetra/doc/html/index.html

(scroll to bottom of page).  In particular, Lessons 02 and 05 will explain how to create different kinds of Maps and use them to redistribute data.  Linear solvers want to solve with uniquely owned Maps (see Lesson 02).  If you want to redistribute the solution vector from the uniquely owned Map back to the overlapping distribution, create and use an Import object (see Lesson 05).

mfh



On 7/31/16, 8:46 PM, "V M Krushnarao Kotteda" <kvmkrao at gmail.com<mailto:kvmkrao at gmail.com>> wrote:

Hi All,

I am using the code to solve the Tri-diagonal system of equations resulting
from two overlapping sub-domains. Each sub-domain contains 10 nodes/equations
including two overlapping nodes.

Proc # 0 (sub-domain # 1)
lglob[1] = 1;            lglob[2] = 2;           lglob[3] = 3;           lglob[4] = 4;             lglob[5] = 5;
lglob[6] = 10;          lglob[7] = 11;         lglob[8] = 12;          lglob[9] = 13;           lglob[10] = 14;
Proc # 1 (sub-domain # 2)
lglob[1] = 5;            lglob[2] = 6;           lglob[3] = 7;            lglob[4] = 8;            lglob[5] = 9;
lglob[6] = 14;          lglob[7] = 15;          lglob[8] = 16;          lglob[9] = 17;          lglob[10] = 18;
lglob array gives the global node number corresponding to a local node.
lglob[] represents the over-lapping nodes.

The following code arranges the linear system of equations on each processor and
solves the equations for 18 unknowns.
The local map (Tpetra) has only 9 nodes on a processor. Therefore, I get
solution at only 9 nodes. The global equation number corresponds to the
local node is different from values shown above. Please see the output
from the code.

Is there a way to get the solution x[[lglob[i]] at all 10 nodes in a sub-domain.
code:
int main ()
{

.......
  const global_size_t numGlobalElements =  18;

  RCP<const map_type> map (new map_type (numGlobalElements, 0, comm));
  const size_t numMyElements =  10;

  ArrayRCP<size_t> myGlobalElements = arcp<size_t> (numMyElements);


   for(int i = 0; i < numMyElements; ++i) {
   myGlobalElements[i] = glob[i]-1;
   }

  sparse_mat_type A (map, 0);

  for (size_t i = 0; i < numMyElements; ++i) {
   .........
      A.insertGlobalValues (myGlobalElements[i],
                            tuple<GO> (myGlobalElements[i]-1, myGlobalElements[i], myGlobalElements[i]+1),
                            tuple<ST> (negOne, two, negOne));
   ........
  }

A.fillComplete ();

multivector_type b (map, 1);

b.putScalar (STS::one ());

multivector_type x (map, 1, true);
x.putScalar (STS::zero ());

  try {
    Ifpack2::Test::solve (A, b, x, 1, 1, 100, 1000, 1.0e-8, false, "RILUK");
  }
  catch (std::exception& e) {
    lclSuccess = 0;
    errStrm << e.what ();
  }

   ArrayRCP<ST> view;
   int size = x.getLocalLength ();
   Array<ST> copy1(numMyElements);
   view = x.get1dViewNonConst();
   x.get1dCopy(copy1(),numMyElements);

   for(int i=0; i < numMyElements; i ++) {
      std::cout << myRank  << " " << i <<"  " <<copy1[i] << std::endl;
   }
.....
}
Output :
 Tpetra::CrsMatrix (Kokkos refactor):
  Template parameters:
   Scalar: double
   LocalOrdinal: int
   GlobalOrdinal: int
   Node: Kokkos::Compat::KokkosDeviceWrapperNode<Kokkos::Serial, Kokkos::HostSpace>
  isFillComplete: true
  Global dimensions: [18, 18]
  Global number of entries: 52
  Global number of diagonal entries: 18
  Global max number of entries in a row: 3

  Row Map:

   "Tpetra::Map":
    Template parameters:
     LocalOrdinal: int
     GlobalOrdinal: int
     Node: Kokkos::Compat::KokkosDeviceWrapperNode<Kokkos::Serial, Kokkos::HostSpace>
    Global number of entries: 18
    Minimum global index: 0
    Maximum global index: 17
    Index base: 0
    Number of processes: 2
    Uniform: true
    Contiguous: true
    Distributed: true
    Process 0 of 2:
     My number of entries: 9
     My minimum global index: 0
     My maximum global index: 8
     My global indices: [0, 1, 2, 3, 4, 5, 6, 7, 8]
    Process 1 of 2:
     My number of entries: 9
     My minimum global index: 9
     My maximum global index: 17
     My global indices: [9, 10, 11, 12, 13, 14, 15, 16, 17]
  Column Map:
   "Tpetra::Map":
    Template parameters:
     LocalOrdinal: int
     GlobalOrdinal: int
     Node: Kokkos::Compat::KokkosDeviceWrapperNode<Kokkos::Serial, Kokkos::HostSpace>
    Global number of entries: 20
    Minimum global index: 0
    Maximum global index: 17
    Index base: 0
    Number of processes: 2
    Uniform: false
    Contiguous: false
    Distributed: true
    Process 0 of 2:
     My number of entries: 10
     My minimum global index: 0
     My maximum global index: 9
     My global indices: [0, 1, 2, 3, 4, 5, 6, 7, 8, 9]
    Process 1 of 2:
     My number of entries: 10
     My minimum global index: 8
     My maximum global index: 17
     My global indices: [9, 10, 11, 12, 13, 14, 15, 16, 17, 8]
  Domain Map: same as row Map
  Range Map: same as domain Map
  Process rank: 0
   Number of allocated entries: 26
   Number of entries: 26
   Number of diagonal entries: 9
   Max number of entries per row: 3
  Process rank: 1
   Number of allocated entries: 26
   Number of entries: 26
   Number of diagonal entries: 9
   Max number of entries per row: 3
      Proc Rank   Global Row  Num Entries(Index,Value)
              0            0            2 (0, 2)  (1, -1)
              0            1            3 (0, -1)  (1, 2)  (2, -1)
              0            2            3 (1, -1)  (2, 2)  (3, -1)
              0            3            3 (2, -1)  (3, 2)  (4, -1)
              0            4            3 (3, -2)  (4, 4)  (5, -2)
              0            5            3 (4, -1)  (5, 2)  (6, -1)
              0            6            3 (5, -1)  (6, 2)  (7, -1)
              0            7            3 (6, -1)  (7, 2)  (8, -1)
              0            8            3 (7, -1)  (8, 2)  (9, -1)
      Proc Rank   Global Row  Num Entries(Index,Value)
              1            9            3 (9, 2)  (10, -1)  (8, -1)
              1           10            3 (9, -1)  (10, 2)  (11, -1)
              1           11            3 (10, -1)  (11, 2)  (12, -1)
              1           12            3 (11, -1)  (12, 2)  (13, -1)
              1           13            3 (12, -2)  (13, 4)  (14, -2)
              1           14            3 (13, -1)  (14, 2)  (15, -1)
              1           15            3 (14, -1)  (15, 2)  (16, -1)
              1           16            3 (15, -1)  (16, 2)  (17, -1)
              1           17            2 (16, -1)  (17, 2)

 "Tpetra::MultiVector":
  Template parameters:
   Scalar: double
   LocalOrdinal: int
   GlobalOrdinal: int
   Node: Serial/Wrapper
  Global number of rows: 18
  Number of columns: 1
  Process 0 of 2:
   Local number of rows: 9

   Column stride: 9
   Values:
   [8.5; 16; 22.5; 28; 32.5; 36.5; 39.5; 41.5; 42.5]
  Process 1 of 2:
   Local number of rows: 9
   Column stride: 9
   Values:
   [42.5; 41.5; 39.5; 36.5; 32.5; 28; 22.5; 16; 8.5]

The Belos solve took 4 iterations to converge.
It achieved a tolerance of: 0

   Proc #      node#   x[node#]
           0           1   8.5000000000000053
           0           2   16.000000000000011
           0           3   22.500000000000014
           0           4   28.000000000000021
           0           5   32.500000000000028
           0           6   36.500000000000036
           0           7   39.500000000000036
           0           8   41.500000000000043
           0           9   42.500000000000043
           0          10   0.0000000000000000
           1           1   42.500000000000043
           1           2   41.500000000000043
           1           3   39.500000000000036
           1           4   36.500000000000036
           1           5   32.500000000000036
           1           6   28.000000000000036
           1           7   22.500000000000025
           1           8   16.000000000000014
           1           9   8.5000000000000089
           1          10   0.0000000000000000

Thanks.
On Tue, Jul 26, 2016 at 11:58 AM, Hoemmen, Mark <mhoemme at sandia.gov<mailto:mhoemme at sandia.gov>> wrote:
On 7/26/16, 11:04 AM, "V M Krushnarao Kotteda" <kvmkrao at gmail.com<mailto:kvmkrao at gmail.com>> wrote:
> Ifpack2 provides various solvers and pre-conditioners.
> Please let me know if there is a way to choose proper pre-conditioner and solver for the given system of equations?

There are too many different kinds of linear systems for us to try to choose a preconditioner and solver automatically.  I can only recommend reading standard textbooks on numerical linear algebra.

  Thank you, Mark.

mfh


Best regards,
V M Krushnarao Kotteda
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <https://trilinos.org/pipermail/trilinos-users/attachments/20160801/9bffa980/attachment.html>


More information about the Trilinos-Users mailing list