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