1c4762a1bSJed Brown static char help[] = "A test for MatAssembly that heavily relies on PetscSortIntWithArrayPair\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown The characteristic of the array (about 4M in length) to sort in this test is that it has 5c4762a1bSJed Brown many duplicated values that already clustered together (around 95 duplicates per unique integer). 6c4762a1bSJed Brown 7c4762a1bSJed Brown It was gotten from a petsc performance bug report from the Moose project. One can use 8c4762a1bSJed Brown it for future performance study. 9c4762a1bSJed Brown 10c4762a1bSJed Brown Contributed-by: Fande Kong <fdkong.jd@gmail.com>, John Peterson <jwpeterson@gmail.com> 11c4762a1bSJed Brown */ 12c4762a1bSJed Brown 13c4762a1bSJed Brown // PETSc includes 14c4762a1bSJed Brown #include <petscmat.h> 15c4762a1bSJed Brown 16c4762a1bSJed Brown // C++ includes 17c4762a1bSJed Brown #include <iostream> 18c4762a1bSJed Brown #include <fstream> 19c4762a1bSJed Brown #include <sstream> 20c4762a1bSJed Brown #include <algorithm> 21c4762a1bSJed Brown #include <vector> 22c4762a1bSJed Brown #include <string> 23c4762a1bSJed Brown #include <set> 24c4762a1bSJed Brown 25c4762a1bSJed Brown int main (int argc, char** argv) 26c4762a1bSJed Brown { 27c4762a1bSJed Brown PetscMPIInt size, rank; 28c4762a1bSJed Brown char file[2][PETSC_MAX_PATH_LEN]; 29c4762a1bSJed Brown PetscBool flg; 30c4762a1bSJed Brown const unsigned int n_dofs = 26559; 31c4762a1bSJed Brown unsigned int first_local_index; 32c4762a1bSJed Brown unsigned int last_local_index; 33c4762a1bSJed Brown 34*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 355f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 365f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 372c71b3e2SJacob Faibussowitsch PetscCheckFalse(size > 2,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This example is for <=2 procs"); 38c4762a1bSJed Brown 395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f0",file[0],sizeof(file[0]),&flg)); 4028b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate dof indices file for rank 0 with -f0 option"); 41c4762a1bSJed Brown if (size == 2) { 425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f1",file[1],sizeof(file[1]),&flg)); 4328b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate dof indices file for rank 1 with -f1 option"); 44c4762a1bSJed Brown } 45c4762a1bSJed Brown 46c4762a1bSJed Brown if (size == 1) { 47c4762a1bSJed Brown first_local_index = 0; 48c4762a1bSJed Brown last_local_index = 26559; 49c4762a1bSJed Brown } else { 50c4762a1bSJed Brown if (rank == 0) { 51c4762a1bSJed Brown first_local_index = 0; 52c4762a1bSJed Brown last_local_index = 13911; 53c4762a1bSJed Brown } else { 54c4762a1bSJed Brown first_local_index = 13911; 55c4762a1bSJed Brown last_local_index = 26559; 56c4762a1bSJed Brown } 57c4762a1bSJed Brown } 58c4762a1bSJed Brown 59c4762a1bSJed Brown // Read element dof indices from files 60c4762a1bSJed Brown std::vector<std::vector<std::vector<PetscInt> > > elem_dof_indices(size); 61c4762a1bSJed Brown for (PetscInt proc_id = 0; proc_id < size; ++proc_id) { 62c4762a1bSJed Brown std::string line; 63c4762a1bSJed Brown std::ifstream dof_file(file[proc_id]); 64c4762a1bSJed Brown if (dof_file.good()) { 65c4762a1bSJed Brown while (std::getline (dof_file,line)) { 66c4762a1bSJed Brown std::vector<PetscInt> dof_indices; 67c4762a1bSJed Brown std::stringstream sstream(line); 68c4762a1bSJed Brown std::string token; 69c4762a1bSJed Brown while (std::getline(sstream, token, ' ')) {dof_indices.push_back(std::atoi(token.c_str()));} 70c4762a1bSJed Brown elem_dof_indices[proc_id].push_back(dof_indices); 71c4762a1bSJed Brown } 7298921bdaSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Could not open file %s",file[proc_id]); 73c4762a1bSJed Brown } 74c4762a1bSJed Brown 75c4762a1bSJed Brown // Debugging: Verify we read in elem_dof_indices correctly 76c4762a1bSJed Brown // for (unsigned int i=0; i<elem_dof_indices.size(); ++i) 77c4762a1bSJed Brown // { 78c4762a1bSJed Brown // for (unsigned int j=0; j<elem_dof_indices[i].size(); ++j) 79c4762a1bSJed Brown // { 80c4762a1bSJed Brown // for (unsigned int k=0; k<elem_dof_indices[i][j].size(); ++k) 81c4762a1bSJed Brown // std::cout << elem_dof_indices[i][j][k] << " "; 82c4762a1bSJed Brown // std::cout << std::endl; 83c4762a1bSJed Brown // } 84c4762a1bSJed Brown // std::cout << std::endl; 85c4762a1bSJed Brown // } 86c4762a1bSJed Brown 87c4762a1bSJed Brown // Set up the (global) sparsity pattern 88c4762a1bSJed Brown std::vector< std::set< unsigned int > > sparsity(n_dofs); 89c4762a1bSJed Brown for (PetscInt proc_id = 0; proc_id < size; ++proc_id) 90c4762a1bSJed Brown for (unsigned int k = 0; k < elem_dof_indices[proc_id].size(); k++) { 91c4762a1bSJed Brown std::vector<PetscInt>& dof_indices = elem_dof_indices[proc_id][k]; 92c4762a1bSJed Brown for (unsigned int i = 0; i < dof_indices.size(); ++i) 93c4762a1bSJed Brown for (unsigned int j = 0; j < dof_indices.size(); ++j) 94c4762a1bSJed Brown sparsity[dof_indices[i]].insert(dof_indices[j]); 95c4762a1bSJed Brown } 96c4762a1bSJed Brown 97c4762a1bSJed Brown // Determine the local nonzeros on this processor 98c4762a1bSJed Brown const unsigned int n_local_dofs = last_local_index - first_local_index; 99c4762a1bSJed Brown std::vector<PetscInt> n_nz(n_local_dofs); 100c4762a1bSJed Brown std::vector<PetscInt> n_oz(n_local_dofs); 101c4762a1bSJed Brown 102c4762a1bSJed Brown for (unsigned int i = 0; i < n_local_dofs; ++i) { 103c4762a1bSJed Brown for (std::set<unsigned int>::iterator iter = sparsity[i+first_local_index].begin(); iter != sparsity[i+first_local_index].end(); iter++) { 104c4762a1bSJed Brown unsigned int dof = *iter; 105c4762a1bSJed Brown if ((dof >= first_local_index) && (dof < last_local_index)) n_nz[i]++; 106c4762a1bSJed Brown else n_oz[i]++; 107c4762a1bSJed Brown } 108c4762a1bSJed Brown } 109c4762a1bSJed Brown 110c4762a1bSJed Brown // Debugging: print number of on/off diagonal nonzeros 111c4762a1bSJed Brown // for (unsigned int i=0; i<n_nz.size(); ++i) 112c4762a1bSJed Brown // std::cout << n_nz[i] << " "; 113c4762a1bSJed Brown // std::cout << std::endl; 114c4762a1bSJed Brown 115c4762a1bSJed Brown // for (unsigned int i=0; i<n_oz.size(); ++i) 116c4762a1bSJed Brown // std::cout << n_oz[i] << " "; 117c4762a1bSJed Brown // std::cout << std::endl; 118c4762a1bSJed Brown 119c4762a1bSJed Brown // Compute and print max number of on- and off-diagonal nonzeros. 120c4762a1bSJed Brown PetscInt n_nz_max = *std::max_element(n_nz.begin(), n_nz.end()); 121c4762a1bSJed Brown PetscInt n_oz_max = *std::max_element(n_oz.begin(), n_oz.end()); 122c4762a1bSJed Brown 1235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "Max on-diagonal non-zeros: = %" PetscInt_FMT "\n", n_nz_max)); 1245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "Max off-diagonal non-zeros: = %" PetscInt_FMT "\n", n_oz_max)); 1255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT)); 126c4762a1bSJed Brown 127c4762a1bSJed Brown // Initialize the matrix similar to what we do in the PetscMatrix 128c4762a1bSJed Brown // ctor and init() routines. 129c4762a1bSJed Brown Mat mat; 1305f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD, &mat)); 1315f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(mat, n_local_dofs, n_local_dofs, n_dofs, n_dofs)); 1325f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetBlockSize(mat, 1)); 1335f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(mat, MATAIJ)); // Automatically chooses seqaij or mpiaij 1345f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJSetPreallocation(mat, 0, &n_nz[0])); 1355f80ce2aSJacob Faibussowitsch CHKERRQ(MatMPIAIJSetPreallocation(mat, 0, &n_nz[0], 0, &n_oz[0])); 1365f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(mat, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE)); 137c4762a1bSJed Brown 138c4762a1bSJed Brown // Local "element" loop 139c4762a1bSJed Brown for (unsigned int k = 0; k < elem_dof_indices[rank].size(); k++) { 140c4762a1bSJed Brown std::vector<PetscInt>& dof_indices = elem_dof_indices[rank][k]; 141c4762a1bSJed Brown // DenseMatrix< Number > zero_mat( dof_indices.size(), dof_indices.size()); 142c4762a1bSJed Brown // B.add_matrix( zero_mat, dof_indices); 143c4762a1bSJed Brown std::vector<PetscScalar> ones(dof_indices.size() * dof_indices.size(), 1.); 1445f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(mat, dof_indices.size(), &dof_indices[0], dof_indices.size(), &dof_indices[0], &ones[0], ADD_VALUES)); 145c4762a1bSJed Brown } 146c4762a1bSJed Brown 147c4762a1bSJed Brown // Matrix assembly 1485f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY)); 1495f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY)); 150c4762a1bSJed Brown 151c4762a1bSJed Brown // Clean up 1525f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&mat)); 153*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 154*b122ec5aSJacob Faibussowitsch return 0; 155c4762a1bSJed Brown } 156c4762a1bSJed Brown /*TEST 157c4762a1bSJed Brown build: 158dfd57a17SPierre Jolivet requires: !defined(PETSC_HAVE_SUN_CXX) 159c4762a1bSJed Brown 160c4762a1bSJed Brown test: 161c4762a1bSJed Brown nsize: 2 162c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/meshes/moose_dof_indices_np_2_proc_0.txt -f1 ${DATAFILESPATH}/meshes/moose_dof_indices_np_2_proc_1.txt 163c4762a1bSJed Brown # Skip the test for Sun C++ compiler because of its warnings/errors in petscmat.h 164c4762a1bSJed Brown requires: datafilespath 165c4762a1bSJed Brown 166c4762a1bSJed Brown TEST*/ 167