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 7f0b74427SPierre Jolivet 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 25d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 26d71ae5a4SJacob Faibussowitsch { 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 34327415f7SBarry Smith PetscFunctionBeginUser; 35c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, nullptr, help)); 369566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 379566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 3808401ef6SPierre Jolivet PetscCheck(size <= 2, PETSC_COMM_WORLD, PETSC_ERR_SUP, "This example is for <=2 procs"); 39c4762a1bSJed Brown 409566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-f0", file[0], sizeof(file[0]), &flg)); 4128b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Must indicate dof indices file for rank 0 with -f0 option"); 42c4762a1bSJed Brown if (size == 2) { 439566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-f1", file[1], sizeof(file[1]), &flg)); 4428b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Must indicate dof indices file for rank 1 with -f1 option"); 45c4762a1bSJed Brown } 46c4762a1bSJed Brown 47c4762a1bSJed Brown if (size == 1) { 48c4762a1bSJed Brown first_local_index = 0; 49c4762a1bSJed Brown last_local_index = 26559; 50c4762a1bSJed Brown } else { 51c4762a1bSJed Brown if (rank == 0) { 52c4762a1bSJed Brown first_local_index = 0; 53c4762a1bSJed Brown last_local_index = 13911; 54c4762a1bSJed Brown } else { 55c4762a1bSJed Brown first_local_index = 13911; 56c4762a1bSJed Brown last_local_index = 26559; 57c4762a1bSJed Brown } 58c4762a1bSJed Brown } 59c4762a1bSJed Brown 60c4762a1bSJed Brown // Read element dof indices from files 61c4762a1bSJed Brown std::vector<std::vector<std::vector<PetscInt>>> elem_dof_indices(size); 62c4762a1bSJed Brown for (PetscInt proc_id = 0; proc_id < size; ++proc_id) { 63c4762a1bSJed Brown std::string line; 64c4762a1bSJed Brown std::ifstream dof_file(file[proc_id]); 65966bd95aSPierre Jolivet PetscCheck(dof_file.good(), PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Could not open file %s", file[proc_id]); 66c4762a1bSJed Brown while (std::getline(dof_file, line)) { 67c4762a1bSJed Brown std::vector<PetscInt> dof_indices; 68c4762a1bSJed Brown std::stringstream sstream(line); 69c4762a1bSJed Brown std::string token; 70ad540459SPierre Jolivet while (std::getline(sstream, token, ' ')) dof_indices.push_back(std::atoi(token.c_str())); 71c4762a1bSJed Brown elem_dof_indices[proc_id].push_back(dof_indices); 72c4762a1bSJed Brown } 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) 93*92f616caSPierre Jolivet for (unsigned int j = 0; j < dof_indices.size(); ++j) sparsity[dof_indices[i]].insert((unsigned int)dof_indices[j]); 94c4762a1bSJed Brown } 95c4762a1bSJed Brown 96c4762a1bSJed Brown // Determine the local nonzeros on this processor 97c4762a1bSJed Brown const unsigned int n_local_dofs = last_local_index - first_local_index; 98c4762a1bSJed Brown std::vector<PetscInt> n_nz(n_local_dofs); 99c4762a1bSJed Brown std::vector<PetscInt> n_oz(n_local_dofs); 100c4762a1bSJed Brown 101c4762a1bSJed Brown for (unsigned int i = 0; i < n_local_dofs; ++i) { 102c4762a1bSJed Brown for (std::set<unsigned int>::iterator iter = sparsity[i + first_local_index].begin(); iter != sparsity[i + first_local_index].end(); iter++) { 103c4762a1bSJed Brown unsigned int dof = *iter; 104c4762a1bSJed Brown if ((dof >= first_local_index) && (dof < last_local_index)) n_nz[i]++; 105c4762a1bSJed Brown else n_oz[i]++; 106c4762a1bSJed Brown } 107c4762a1bSJed Brown } 108c4762a1bSJed Brown 1094cf0e950SBarry Smith // Debugging: print number of on/off-diagonal nonzeros 110c4762a1bSJed Brown // for (unsigned int i=0; i<n_nz.size(); ++i) 111c4762a1bSJed Brown // std::cout << n_nz[i] << " "; 112c4762a1bSJed Brown // std::cout << std::endl; 113c4762a1bSJed Brown 114c4762a1bSJed Brown // for (unsigned int i=0; i<n_oz.size(); ++i) 115c4762a1bSJed Brown // std::cout << n_oz[i] << " "; 116c4762a1bSJed Brown // std::cout << std::endl; 117c4762a1bSJed Brown 118c4762a1bSJed Brown // Compute and print max number of on- and off-diagonal nonzeros. 119c4762a1bSJed Brown PetscInt n_nz_max = *std::max_element(n_nz.begin(), n_nz.end()); 120c4762a1bSJed Brown PetscInt n_oz_max = *std::max_element(n_oz.begin(), n_oz.end()); 121c4762a1bSJed Brown 1229566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "Max on-diagonal non-zeros: = %" PetscInt_FMT "\n", n_nz_max)); 1239566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "Max off-diagonal non-zeros: = %" PetscInt_FMT "\n", n_oz_max)); 1249566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 125c4762a1bSJed Brown 126c4762a1bSJed Brown // Initialize the matrix similar to what we do in the PetscMatrix 127c4762a1bSJed Brown // ctor and init() routines. 128c4762a1bSJed Brown Mat mat; 1299566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &mat)); 1309566063dSJacob Faibussowitsch PetscCall(MatSetSizes(mat, n_local_dofs, n_local_dofs, n_dofs, n_dofs)); 1319566063dSJacob Faibussowitsch PetscCall(MatSetBlockSize(mat, 1)); 1329566063dSJacob Faibussowitsch PetscCall(MatSetType(mat, MATAIJ)); // Automatically chooses seqaij or mpiaij 1339566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(mat, 0, &n_nz[0])); 1349566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(mat, 0, &n_nz[0], 0, &n_oz[0])); 1359566063dSJacob Faibussowitsch PetscCall(MatSetOption(mat, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE)); 136c4762a1bSJed Brown 137c4762a1bSJed Brown // Local "element" loop 138c4762a1bSJed Brown for (unsigned int k = 0; k < elem_dof_indices[rank].size(); k++) { 139c4762a1bSJed Brown std::vector<PetscInt> &dof_indices = elem_dof_indices[rank][k]; 140c4762a1bSJed Brown // DenseMatrix< Number > zero_mat( dof_indices.size(), dof_indices.size()); 141c4762a1bSJed Brown // B.add_matrix( zero_mat, dof_indices); 142c4762a1bSJed Brown std::vector<PetscScalar> ones(dof_indices.size() * dof_indices.size(), 1.); 1439566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat, dof_indices.size(), &dof_indices[0], dof_indices.size(), &dof_indices[0], &ones[0], ADD_VALUES)); 144c4762a1bSJed Brown } 145c4762a1bSJed Brown 146c4762a1bSJed Brown // Matrix assembly 1479566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY)); 1489566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY)); 149c4762a1bSJed Brown 150c4762a1bSJed Brown // Clean up 1519566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat)); 1529566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 153b122ec5aSJacob Faibussowitsch return 0; 154c4762a1bSJed Brown } 155c4762a1bSJed Brown /*TEST 156c4762a1bSJed Brown build: 157dfd57a17SPierre Jolivet requires: !defined(PETSC_HAVE_SUN_CXX) 158c4762a1bSJed Brown 159c4762a1bSJed Brown test: 160c4762a1bSJed Brown nsize: 2 161c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/meshes/moose_dof_indices_np_2_proc_0.txt -f1 ${DATAFILESPATH}/meshes/moose_dof_indices_np_2_proc_1.txt 162c4762a1bSJed Brown # Skip the test for Sun C++ compiler because of its warnings/errors in petscmat.h 163c4762a1bSJed Brown requires: datafilespath 164c4762a1bSJed Brown 165c4762a1bSJed Brown TEST*/ 166