|
Galeri
Development
|
File galeri/example/Maps.cpp
// @HEADER // ************************************************************************ // // Galeri: Finite Element and Matrix Generation Package // Copyright (2006) ETHZ/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. // // 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 about Galeri? Contact Marzio Sala (marzio.sala _AT_ gmail.com) // // ************************************************************************ // @HEADER #include "Galeri_Maps.h" #ifdef HAVE_MPI #include "Epetra_MpiComm.h" #include "mpi.h" #else #include "Epetra_SerialComm.h" #endif #include "Teuchos_ParameterList.hpp" using namespace Galeri; int main(int argc, char* argv[]) { #ifdef HAVE_MPI MPI_Init(&argc, &argv); Epetra_MpiComm Comm(MPI_COMM_WORLD); #else Epetra_SerialComm Comm; #endif // Creates an Epetra_Map corresponding to a 2D Cartesian grid // on the unit square. For parallel runs, the nodes are divided into // strips, so that the total number of subdomains is Comm.NumProc() x 1. // Pointer to the object to be created Epetra_Map* Map = 0; // Type of the object string MapType = "Cartesian2D"; // Container for parameters Teuchos::ParameterList GaleriList; GaleriList.set("nx", 2 * Comm.NumProc()); GaleriList.set("ny", 2); GaleriList.set("mx", Comm.NumProc()); GaleriList.set("my", 1); try { // Creation of the map #ifndef GALERI_TEST_USE_LONGLONG_GO Map = CreateMap("Cartesian2D", Comm, GaleriList); #else Map = CreateMap64("Cartesian2D", Comm, GaleriList); #endif // print out the map cout << *Map; // To created object must be free'd using delete delete Map; } catch (Exception& rhs) { if (Comm.MyPID() == 0) { cerr << "Caught exception: "; rhs.Print(); } } #ifdef HAVE_MPI MPI_Finalize(); #endif return(EXIT_SUCCESS); }
1.7.6.1