[Trilinos-Users] [EXTERNAL] Tpetra::MatrixMarket: Write and Read MultiVector

Sander Schaffner ssander at student.ethz.ch
Mon Feb 9 07:21:29 MST 2015


Hello again,

the new version helpt in the sense that it reads now the file. But still 
it behaves like described before:

[1 2 3 ... 63] -> [0 1 64 65 2 3 ...]
[64 65 66 ... 127] -> [32 33 96 97 34 35 ...]

I wrote a small test program: First I create a 64x2 MV on 32 cores and 
write into the first MV  [1 2 3 ... 63] and into the second [64 65 66 
... 127]. Then I write out the map and the MV. After that I read both 
files to create a new MV to finally write it out again. Since I changed 
nothing inbetween and only read and write a valid file i expect two 
identical files. But this is not the case.

Can someone tell me if I'm doing something wrong in my Code, or is it 
maybe a bug?

Thanks

Sander

Output first File (test_64_2):

%%MatrixMarket matrix array real general
64 2
0
1
64
65
2
3
66
67
4
5
68
69
6
7
...

Second (test_64):

%%MatrixMarket matrix array real general
64 2
0
1
32
33
64
65
96
97
2
3
34
35
66
67

Am 06.02.2015 um 19:04 schrieb Andrey Prokopenko:
> Hi Sander,
>
> I believe that was fixed (see 
> https://software.sandia.gov/bugzilla/show_bug.cgi?id=6264).
> Which version of Trilinos do you use?
>
> -Andrey
> On 02/06/2015 02:25 AM, Sander Schaffner wrote:
>> Hi
>>
>> I'm trying to write a MultiVector with 
>> Tpetra::MatrixMarket::Writer.writeDenseFile. The goal is to read the 
>> same MultiVector afterwards with 
>> Tpetra::MatrixMarket::Reader.readDenseFile. But there are some 
>> reasons why this doesn't work:
>>
>> If i generate a MultiVector with two Vectors (dimension 128) where 
>> the first holds the values 0 and the second 1 I get this output on 32 
>> cores:
>>
>> %%MatrixMarket matrix array real general
>> 128 2
>> 0
>> 0
>> 0
>> 0
>>
>> 1
>> 1
>> 1
>> 1
>>
>> 0
>> 0
>> 0
>> 0
>>
>> 1
>> 1
>> 1
>> 1
>>
>> and so on. So I get first the local part of processor 0 for both 
>> Vectors and so on.
>>
>> The reader does not accept this file. It complains about the spaces. 
>> If I manually get ridd of them it reads the file. But then it handles 
>> the first half of the values as Vector 1 and the second part as 
>> Vector 2. So in the End we have the following:
>>
>> [0 0 0 0 0 ... 0] -> [0 0 0 0 1 1 1 1 ...]
>> [1 1 1 1 ... 1] -> [0 0 0 0 1 1 1 1 ...]
>>
>> Is there a way of directly read in a MatrixMarket file created by the 
>> Writer? What other options are there to write/read MultiVectors?
>>
>> Any help is appreciated.
>>
>> Sander
>> _______________________________________________
>> Trilinos-Users mailing list
>> Trilinos-Users at software.sandia.gov
>> https://software.sandia.gov/mailman/listinfo/trilinos-users
>



Am 06.02.2015 um 19:04 schrieb Andrey Prokopenko:
> Hi Sander,
>
> I believe that was fixed (see 
> https://software.sandia.gov/bugzilla/show_bug.cgi?id=6264).
> Which version of Trilinos do you use?
>
> -Andrey
> On 02/06/2015 02:25 AM, Sander Schaffner wrote:
>> Hi
>>
>> I'm trying to write a MultiVector with 
>> Tpetra::MatrixMarket::Writer.writeDenseFile. The goal is to read the 
>> same MultiVector afterwards with 
>> Tpetra::MatrixMarket::Reader.readDenseFile. But there are some 
>> reasons why this doesn't work:
>>
>> If i generate a MultiVector with two Vectors (dimension 128) where 
>> the first holds the values 0 and the second 1 I get this output on 32 
>> cores:
>>
>> %%MatrixMarket matrix array real general
>> 128 2
>> 0
>> 0
>> 0
>> 0
>>
>> 1
>> 1
>> 1
>> 1
>>
>> 0
>> 0
>> 0
>> 0
>>
>> 1
>> 1
>> 1
>> 1
>>
>> and so on. So I get first the local part of processor 0 for both 
>> Vectors and so on.
>>
>> The reader does not accept this file. It complains about the spaces. 
>> If I manually get ridd of them it reads the file. But then it handles 
>> the first half of the values as Vector 1 and the second part as 
>> Vector 2. So in the End we have the following:
>>
>> [0 0 0 0 0 ... 0] -> [0 0 0 0 1 1 1 1 ...]
>> [1 1 1 1 ... 1] -> [0 0 0 0 1 1 1 1 ...]
>>
>> Is there a way of directly read in a MatrixMarket file created by the 
>> Writer? What other options are there to write/read MultiVectors?
>>
>> Any help is appreciated.
>>
>> Sander
>> _______________________________________________
>> Trilinos-Users mailing list
>> Trilinos-Users at software.sandia.gov
>> https://software.sandia.gov/mailman/listinfo/trilinos-users
>

-------------- next part --------------
/*
//@HEADER
// ***********************************************************************
//
//     EpetraExt: Epetra Extended - Linear Algebra Services Package
//                 Copyright (2011) Sandia Corporation
//
// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
// the U.S. Government retains certain rights in this software.
//
// Redistribution and use in source and binary forms, with or without
// modification, are permitted provided that the following conditions are
// met:
//
// 1. Redistributions of source code must retain the above copyright
// notice, this list of conditions and the following disclaimer.
//
// 2. Redistributions in binary form must reproduce the above copyright
// notice, this list of conditions and the following disclaimer in the
// documentation and/or other materials provided with the distribution.
//
// 3. Neither the name of the Corporation nor the names of the
// contributors may be used to endorse or promote products derived from
// this software without specific prior written permission.
//
// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
//
// Questions? Contact Michael A. Heroux (maherou at sandia.gov)
//
// ***********************************************************************
//@HEADER
*/

//////
// Std
//////

#include <fstream> 
#include <string>

//////////
// Teuchos
//////////

#include <Teuchos_Array.hpp>
#include <Teuchos_ScalarTraits.hpp>
#include <Teuchos_OrdinalTraits.hpp>
#include <Teuchos_RCP.hpp>
#include <Teuchos_GlobalMPISession.hpp>
#include <Teuchos_oblackholestream.hpp>
#include <Teuchos_TimeMonitor.hpp>
#include <Teuchos_ArrayRCPDecl.hpp> // debug: vector

// Include header for Teuchos serial dense matrix
#include "Teuchos_SerialDenseMatrix.hpp"

#include "Teuchos_ArrayViewDecl.hpp"

#include "Teuchos_ParameterList.hpp"

/////////
// Tpetra
/////////

#include "Tpetra_DefaultPlatform.hpp"
#include "Tpetra_Version.hpp"
#include "Tpetra_Map.hpp"
#include "Tpetra_MultiVector.hpp"
#include "Tpetra_Vector.hpp"
#include "Tpetra_CrsMatrix.hpp"

#include <MatrixMarket_Tpetra.hpp>

// Include header for sparse matrix operations
#include <TpetraExt_MatrixMatrix_def.hpp>

using namespace std;

  //
  // Specify types
  //                                   
  typedef double                                                                Scalar;
  typedef int                                                                   Integer;
  typedef int                                                                   Local_Ordinal;
  typedef int                                                                   Global_Ordinal;
  typedef Teuchos::ScalarTraits<Scalar>                                         SCT; 
  typedef Tpetra::DefaultPlatform::DefaultPlatformType                          Platform; 
  typedef Tpetra::DefaultPlatform::DefaultPlatformType::NodeType                Node;
  typedef Tpetra::CrsMatrix<Scalar,Local_Ordinal,Global_Ordinal,Node>           CrsMatrix;
  typedef Tpetra::Vector<Scalar,Local_Ordinal,Global_Ordinal,Node>              Vector;
  typedef Tpetra::MultiVector<Scalar,Local_Ordinal,Global_Ordinal,Node>         MV;
  typedef Tpetra::Operator<Scalar,Local_Ordinal,Global_Ordinal,Node>            OP;
  typedef Tpetra::Map<Local_Ordinal, Global_Ordinal>                            Map;
  typedef Tpetra::CrsMatrix<Scalar, Local_Ordinal, Global_Ordinal>              Matrix;
  typedef Tpetra::MatrixMarket::Writer<CrsMatrix>                               Writer;
  typedef Tpetra::MatrixMarket::Reader<CrsMatrix>                               Reader;

int main (int argc, char **argv)
{
  // global_size_t: Tpetra defines this unsigned integer type big
  // enough to hold any global dimension or amount of data.
  using Tpetra::global_size_t;
  using Teuchos::Array;
  using Teuchos::ArrayView;
  using Teuchos::ArrayRCP;
  using Teuchos::arcp;
  using Teuchos::RCP;
  using Teuchos::rcp;
  using Teuchos::tuple;
  using std::cerr;
  using std::cout;
  using std::endl;
  using Teuchos::Time;
  using Teuchos::TimeMonitor;

  // initialize MPI. We need the rawComm for the class HDF5!
  MPI_Init (&argc, &argv);
  const MPI_Comm rawComm = MPI_COMM_WORLD;

  // initialize the MPI-Communicator for Tpetra
  RCP<const Teuchos::MpiComm<int> > Comm = rcp (new Teuchos::MpiComm<int> (rawComm));

  const size_t myRank = Comm->getRank();

  // Make an output stream (for verbose output) that only prints on
  // Proc 0 of the communicator.
  Teuchos::oblackholestream blackHole;
  std::ostream& out = (myRank == 0) ? std::cout : blackHole;
  std::ostream& err = (myRank == 0) ? std::cerr : blackHole;

  RCP<const Map> RangeMap = rcp (new Map (64, 0, Comm));

  out << endl << "*) Create MultiVector!" << endl << endl;

  RCP<MV> mv = Teuchos::rcp( new MV(RangeMap, 2) );

  for(Global_Ordinal i=mv->getMap()->getMinGlobalIndex(); i<mv->getMap()->getMaxGlobalIndex()+1; i++){
    for(Integer j=0; j<2; j++){
      mv->replaceGlobalValue(i,j,j*64+i);
    }
  }

  out << endl << "*) Write MultiVector!" << endl << endl;

  std::string name = "test_64";
  Writer::writeDenseFile(name,mv);
  name += "_map";
  RCP<const Map> mv_map = mv->getMap();
  Writer::writeMapFile(name,*mv_map);

  out << endl << "*) Read MultiVector!" << endl << endl;

  Platform &platform = Tpetra::DefaultPlatform::getDefaultPlatform();
  RCP<Node>                      node = platform.getNode(); 
  RCP<const Map> read_mv_map = Reader::readMapFile("test_64_map", Comm, node, false, false);    

  out << endl << "*) Write MultiVector again!" << endl << endl;

  RCP<MV> read_mv = Teuchos::rcp( new MV(read_mv_map, 2) );
  read_mv = Reader::readDenseFile("test_64", Comm, node, read_mv_map,false,false);

  name = "test_64_2";
  Writer::writeDenseFile(name,read_mv);

  out << endl << "*) End of programm" << endl << endl;

#ifdef HAVE_MPI
  MPI_Finalize();
#endif

  return(EXIT_SUCCESS);
}


More information about the Trilinos-Users mailing list