[Trilinos-Users] Possible Bug (or feature??) in importing non local elements

rrossi at cimne.upc.edu rrossi at cimne.upc.edu
Sun Jan 16 04:39:35 MST 2011



Dear All, and in particular dear Alan, 

i finally managed to import non
local elements. 

I attach the file i used to test. 

Unfortunately the
thing only works if at line 149 i use "Add". A different (and wrong) result
is obtained if i use "Insert". To my understanding in this case it should
NOT change! 

Is this a bug or a feature? 

i post my test here. I don't
know who should decide on this, but i would work it up a bit more and
prepare it as a "ex25" example for didasko, as it took me a while to figure
out how it works... 

greetings and thank you all 

Riccardo
-------------- next part --------------
An HTML attachment was scrubbed...
URL: https://software.sandia.gov/pipermail/trilinos-users/attachments/20110116/4316e637/attachment.html 
-------------- next part --------------
// @HEADER
// ***********************************************************************
// 
//                      Didasko Tutorial Package
//                 Copyright (2005) Sandia Corporation
// 
// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
// license for use of this work by or on behalf of the U.S. Government.
// 
// This library is free software; you can 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.
//  
// This library is distributed in the hope that it will be useful, but
// WITHOUT ANY WARRANTY; without even the implied warranty of
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
// Lesser General Public License for more details.
//  
// You should have received a copy of the GNU Lesser General Public
// License along with this library; if not, write to the Free Software
// Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307
// USA
// Questions? Contact Michael A. Heroux (maherou at sandia.gov) 
// 
// ***********************************************************************
// @HEADER

// Show the usage of RowMap, ColumnMap, RangeMap and DomainMap
// This code must be run with two processes

#include "Didasko_ConfigDefs.h"
//#if defined(HAVE_DIDASKO_EPETRA)

#include "Epetra_ConfigDefs.h"
//#ifdef HAVE_MPI
#include "mpi.h"
#include "Epetra_MpiComm.h"
#include "Epetra_Vector.h"
#include "Epetra_Map.h"
#include "Epetra_FECrsGraph.h"
#include "Epetra_FECrsMatrix.h"
#include "Epetra_SerialDenseMatrix.h"
#include "Epetra_Import.h"

int main(int argc, char *argv[])
{

  MPI_Init(&argc, &argv);
  Epetra_MpiComm Comm(MPI_COMM_WORLD);

  int MyPID = Comm.MyPID();
  int NumProc = Comm.NumProc();

  if( NumProc != 2 ) {
    MPI_Finalize();

    return(EXIT_SUCCESS);
  }

  // I define two maps, with and without overlapping elements:
  int non_overlap_nelements = 0;
  int overlap_nelements = 0;
  int ierr;
  
  
  switch( MyPID ) {
  case 0:
    non_overlap_nelements = 2;
    overlap_nelements = 3;
    break;

  case 1:
    non_overlap_nelements = 2;
    overlap_nelements = 3;
    break;
  }

  int* non_overlapping_indices = new int[non_overlap_nelements];
  int* overlapping_indices = new int[overlap_nelements];

  //i assing the indices so that proc0 owns 0 1 and proc2 owns 2 3
  //1 and 2 are overlapping in the overlappign map
  switch( MyPID ) {
  case 0:
    non_overlapping_indices[0] = 0; non_overlapping_indices[1] = 1; 
    overlapping_indices[0] = 0; overlapping_indices[1] = 1; overlapping_indices[2] = 2;
    break;
  case 1:
    non_overlapping_indices[0] = 2; non_overlapping_indices[1] = 3; 
    overlapping_indices[0] = 1; overlapping_indices[1] = 2; overlapping_indices[2] = 3;
    break;
  }
 
  Epetra_Map NonOverlappingMap(-1, non_overlap_nelements, non_overlapping_indices, 0, Comm);  
  Epetra_Map OverlappingMap(-1, overlap_nelements, overlapping_indices, 0, Comm);  
  
  //create the matrix graph
  int guess_row_size = 3;
  Epetra_FECrsGraph Agraph(Copy, NonOverlappingMap, guess_row_size);
  
  //fill the graph
  int* assembling_temp = new int[2];
  switch( MyPID ) {
  case 0:
    assembling_temp[0] = 0; assembling_temp[1] = 1;
    ierr = Agraph.InsertGlobalIndices(2,assembling_temp,2,assembling_temp);
    assembling_temp[0] = 1; assembling_temp[1] = 2;
    ierr = Agraph.InsertGlobalIndices(2,assembling_temp,2,assembling_temp);
    break;
  case 1:
    assembling_temp[0] = 2; assembling_temp[1] = 3;
    ierr = Agraph.InsertGlobalIndices(2,assembling_temp,2,assembling_temp);
    break;
  } 
  delete [] assembling_temp;

  int graph_assemble_ierr = Agraph.GlobalAssemble();
  
  //create a FE matrix (using the precomputed graph) and set it to zero
  Epetra_FECrsMatrix A(Copy,Agraph);
  A.PutScalar(0.0);

  //now do some simple FE assembly of 3 elements of the type (1 -1; -1 1) 
  Epetra_IntSerialDenseVector indices(2);
  Epetra_SerialDenseMatrix values(2, 2);
  values(0,0) = 1.0; values(0,1) = -1.0; values(1,0) = -1.0; values(1,1) = 1.0;

  switch( MyPID ) {
    case 0:
    indices(0) = 0; indices(1) = 1;
    ierr = A.SumIntoGlobalValues(indices, values);
    indices(0) = 1; indices(1) = 2;
    ierr = A.SumIntoGlobalValues(indices, values);
    break;
  case 1:
    indices(0) = 2; indices(1) = 3;
    ierr = A.SumIntoGlobalValues(indices, values);
    break;
  } 
  A.GlobalAssemble();

  cout << "***********************************" << std::endl;
  
  cout << A;
  
  Epetra_Import my_importer(OverlappingMap, NonOverlappingMap);
  Epetra_CrsMatrix B(Copy,OverlappingMap,3);
  ierr = B.Import(A,my_importer,Add); 
  if(ierr != 0) cout << "Epetra failure found";

  cout << "***********************************" << std::endl;
  cout << B;

  MPI_Finalize();


  return(EXIT_SUCCESS);
  
}



More information about the Trilinos-Users mailing list