1735d7f90SBarry Smith static char help[] = "Demonstrates the use of the COO interface to PETSc matrices for finite element computations\n\n"; 2735d7f90SBarry Smith 3735d7f90SBarry Smith /* 4735d7f90SBarry Smith The COO interface for PETSc matrices provides a convenient way to provide finite element element stiffness matrices to PETSc matrix that should work 5735d7f90SBarry Smith well on both CPUs and GPUs. It is an alternative to using MatSetValues() 6735d7f90SBarry Smith 7735d7f90SBarry Smith This example is intended for people who are NOT using DMPLEX or libCEED or any other higher-level infrastructure for finite elements; 8735d7f90SBarry Smith it is only to demonstrate the concepts in a simple way for those people who are interested and for those people who are using PETSc for 9735d7f90SBarry Smith linear algebra solvers but are managing their own finite element process. 10735d7f90SBarry Smith 11735d7f90SBarry Smith Please do NOT use this example as a starting point to writing your own finite element code from scratch! 12735d7f90SBarry Smith 13735d7f90SBarry Smith Each element in this example has three vertices; hence the the usage below needs to be adjusted for elements of a different number of vertices. 14735d7f90SBarry Smith */ 15735d7f90SBarry Smith 16735d7f90SBarry Smith #include <petscmat.h> 17735d7f90SBarry Smith #include "ex18.h" 18735d7f90SBarry Smith 19735d7f90SBarry Smith static PetscErrorCode CreateFEStruct(FEStruct *fe) 20735d7f90SBarry Smith { 21735d7f90SBarry Smith PetscErrorCode ierr; 22735d7f90SBarry Smith 23735d7f90SBarry Smith PetscFunctionBeginUser; 24735d7f90SBarry Smith fe->Nv = 5; 25735d7f90SBarry Smith fe->Ne = 3; 26735d7f90SBarry Smith ierr = PetscMalloc1(3*fe->Ne,&fe->vertices);CHKERRQ(ierr); 27735d7f90SBarry Smith /* the three vertices associated with each element in order of element */ 28735d7f90SBarry Smith fe->vertices[0 + 0] = 0; 29735d7f90SBarry Smith fe->vertices[0 + 1] = 1; 30735d7f90SBarry Smith fe->vertices[0 + 2] = 2; 31735d7f90SBarry Smith fe->vertices[3 + 0] = 2; 32735d7f90SBarry Smith fe->vertices[3 + 1] = 1; 33735d7f90SBarry Smith fe->vertices[3 + 2] = 3; 34735d7f90SBarry Smith fe->vertices[6 + 0] = 2; 35735d7f90SBarry Smith fe->vertices[6 + 1] = 4; 36735d7f90SBarry Smith fe->vertices[6 + 2] = 3; 37735d7f90SBarry Smith fe->n = 5; 38735d7f90SBarry Smith PetscFunctionReturn(0); 39735d7f90SBarry Smith } 40735d7f90SBarry Smith 41735d7f90SBarry Smith static PetscErrorCode DestroyFEStruct(FEStruct *fe) 42735d7f90SBarry Smith { 43735d7f90SBarry Smith PetscErrorCode ierr; 44735d7f90SBarry Smith 45735d7f90SBarry Smith PetscFunctionBeginUser; 46735d7f90SBarry Smith ierr = PetscFree(fe->vertices);CHKERRQ(ierr); 47735d7f90SBarry Smith ierr = PetscFree(fe->coo);CHKERRQ(ierr); 48735d7f90SBarry Smith PetscFunctionReturn(0); 49735d7f90SBarry Smith } 50735d7f90SBarry Smith 51735d7f90SBarry Smith static PetscErrorCode CreateMatrix(FEStruct *fe,Mat *A) 52735d7f90SBarry Smith { 53735d7f90SBarry Smith PetscErrorCode ierr; 54735d7f90SBarry Smith PetscInt *oor,*ooc,cnt = 0; 55735d7f90SBarry Smith 56735d7f90SBarry Smith PetscFunctionBeginUser; 57735d7f90SBarry Smith ierr = MatCreate(PETSC_COMM_WORLD,A);CHKERRQ(ierr); 58735d7f90SBarry Smith ierr = MatSetSizes(*A,fe->n,fe->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 59735d7f90SBarry Smith ierr = MatSetFromOptions(*A);CHKERRQ(ierr); 60735d7f90SBarry Smith 61735d7f90SBarry Smith /* determine for each entry in each element stiffness matrix the global row and colum */ 62735d7f90SBarry Smith /* since the element is triangular with piecewise linear basis functions there are three degrees of freedom per element, one for each vertex */ 63735d7f90SBarry Smith ierr = PetscMalloc2(3*3*fe->Ne,&oor,3*3*fe->Ne,&ooc);CHKERRQ(ierr); 64735d7f90SBarry Smith for (PetscInt e=0; e<fe->Ne; e++) { 65735d7f90SBarry Smith for (PetscInt vi=0; vi<3; vi++) { 66735d7f90SBarry Smith for (PetscInt vj=0; vj<3; vj++) { 67735d7f90SBarry Smith oor[cnt] = fe->vertices[3*e+vi]; 68735d7f90SBarry Smith ooc[cnt++] = fe->vertices[3*e+vj]; 69735d7f90SBarry Smith } 70735d7f90SBarry Smith } 71735d7f90SBarry Smith } 72735d7f90SBarry Smith ierr = MatSetPreallocationCOO(*A,3*3*fe->Ne,oor,ooc);CHKERRQ(ierr); 73735d7f90SBarry Smith ierr = PetscFree2(oor,ooc);CHKERRQ(ierr); 74735d7f90SBarry Smith 75735d7f90SBarry Smith /* determine the offset into the COO value array the offset of each element stiffness; there are 9 = 3*3 entries for each element stiffness */ 76735d7f90SBarry Smith /* for lists of elements with different numbers of degrees of freedom assocated with each element the offsets will not be uniform */ 77735d7f90SBarry Smith ierr = PetscMalloc1(fe->Ne,&fe->coo);CHKERRQ(ierr); 78735d7f90SBarry Smith fe->coo[0] = 0; 79735d7f90SBarry Smith for (PetscInt e=1; e<fe->Ne; e++) { 80735d7f90SBarry Smith fe->coo[e] = fe->coo[e-1] + 3*3; 81735d7f90SBarry Smith } 82735d7f90SBarry Smith PetscFunctionReturn(0); 83735d7f90SBarry Smith } 84735d7f90SBarry Smith 85735d7f90SBarry Smith static PetscErrorCode FillMatrixCPU(FEStruct *fe,Mat A) 86735d7f90SBarry Smith { 87735d7f90SBarry Smith PetscErrorCode ierr; 88735d7f90SBarry Smith PetscScalar s[9]; 89735d7f90SBarry Smith 90735d7f90SBarry Smith PetscFunctionBeginUser; 91735d7f90SBarry Smith /* simulation of traditional PETSc CPU based finite assembly process */ 92735d7f90SBarry Smith for (PetscInt e=0; e<fe->Ne; e++) { 93735d7f90SBarry Smith for (PetscInt vi=0; vi<3; vi++) { 94735d7f90SBarry Smith for (PetscInt vj=0; vj<3; vj++) { 95735d7f90SBarry Smith s[3*vi+vj] = vi+2*vj; 96735d7f90SBarry Smith } 97735d7f90SBarry Smith } 98735d7f90SBarry Smith ierr = MatSetValues(A,3,fe->vertices + 3*e,3, fe->vertices + 3*e,s,ADD_VALUES);CHKERRQ(ierr); 99735d7f90SBarry Smith } 100735d7f90SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 101735d7f90SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 102735d7f90SBarry Smith PetscFunctionReturn(0); 103735d7f90SBarry Smith } 104735d7f90SBarry Smith 105735d7f90SBarry Smith /* 106735d7f90SBarry Smith Shows an example of tracking element offsets explicitly, which allows for 107735d7f90SBarry Smith mixed-topology meshes and combining both volume and surface parts into the weak form. 108735d7f90SBarry Smith */ 109735d7f90SBarry Smith static PetscErrorCode FillMatrixCPUCOO(FEStruct *fe,Mat A) 110735d7f90SBarry Smith { 111735d7f90SBarry Smith PetscErrorCode ierr; 112735d7f90SBarry Smith PetscScalar *v,*s; 113735d7f90SBarry Smith 114735d7f90SBarry Smith PetscFunctionBeginUser; 115735d7f90SBarry Smith /* simulation of CPU based finite assembly process with COO */ 116735d7f90SBarry Smith ierr = PetscMalloc1(3*3*fe->Ne,&v);CHKERRQ(ierr); 117735d7f90SBarry Smith for (PetscInt e=0; e<fe->Ne; e++) { 118735d7f90SBarry Smith s = v + fe->coo[e]; /* point to location in COO of current element stiffness */ 119735d7f90SBarry Smith for (PetscInt vi=0; vi<3; vi++) { 120735d7f90SBarry Smith for (PetscInt vj=0; vj<3; vj++) { 121735d7f90SBarry Smith s[3*vi+vj] = vi+2*vj; 122735d7f90SBarry Smith } 123735d7f90SBarry Smith } 124735d7f90SBarry Smith } 125735d7f90SBarry Smith ierr = MatSetValuesCOO(A,v,ADD_VALUES);CHKERRQ(ierr); 126735d7f90SBarry Smith ierr = PetscFree(v);CHKERRQ(ierr); 127735d7f90SBarry Smith PetscFunctionReturn(0); 128735d7f90SBarry Smith } 129735d7f90SBarry Smith 130735d7f90SBarry Smith /* 131735d7f90SBarry Smith Uses a multi-dimensional indexing technique that works for homogeneous meshes 132735d7f90SBarry Smith such as single-topology with volume integral only. 133735d7f90SBarry Smith */ 134735d7f90SBarry Smith static PetscErrorCode FillMatrixCPUCOO3d(FEStruct *fe,Mat A) 135735d7f90SBarry Smith { 136735d7f90SBarry Smith PetscErrorCode ierr; 137735d7f90SBarry Smith PetscScalar (*s)[3][3]; 138735d7f90SBarry Smith 139735d7f90SBarry Smith PetscFunctionBeginUser; 140735d7f90SBarry Smith /* simulation of CPU based finite assembly process with COO */ 141735d7f90SBarry Smith ierr = PetscMalloc1(fe->Ne,&s);CHKERRQ(ierr); 142735d7f90SBarry Smith for (PetscInt e=0; e<fe->Ne; e++) { 143735d7f90SBarry Smith for (PetscInt vi=0; vi<3; vi++) { 144735d7f90SBarry Smith for (PetscInt vj=0; vj<3; vj++) { 145735d7f90SBarry Smith s[e][vi][vj] = vi+2*vj; 146735d7f90SBarry Smith } 147735d7f90SBarry Smith } 148735d7f90SBarry Smith } 149735d7f90SBarry Smith ierr = MatSetValuesCOO(A,(PetscScalar*)s,INSERT_VALUES);CHKERRQ(ierr); 150735d7f90SBarry Smith ierr = PetscFree(s);CHKERRQ(ierr); 151735d7f90SBarry Smith PetscFunctionReturn(0); 152735d7f90SBarry Smith } 153735d7f90SBarry Smith 154735d7f90SBarry Smith int main(int argc, char **args) 155735d7f90SBarry Smith { 156735d7f90SBarry Smith Mat A; 157735d7f90SBarry Smith PetscErrorCode ierr; 158735d7f90SBarry Smith FEStruct fe; 159735d7f90SBarry Smith PetscMPIInt size; 160735d7f90SBarry Smith PetscBool is_kokkos,is_cuda; 161735d7f90SBarry Smith 162735d7f90SBarry Smith ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 163735d7f90SBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 164*54c59aa7SJacob Faibussowitsch PetscCheck(size <= 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Demonstration is only for sequential runs"); 165735d7f90SBarry Smith 166735d7f90SBarry Smith ierr = CreateFEStruct(&fe);CHKERRQ(ierr); 167735d7f90SBarry Smith ierr = CreateMatrix(&fe,&A);CHKERRQ(ierr); 168735d7f90SBarry Smith 169735d7f90SBarry Smith ierr = FillMatrixCPU(&fe,A);CHKERRQ(ierr); 170735d7f90SBarry Smith ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 171735d7f90SBarry Smith 172735d7f90SBarry Smith ierr = MatZeroEntries(A);CHKERRQ(ierr); 173735d7f90SBarry Smith ierr = FillMatrixCPUCOO(&fe,A);CHKERRQ(ierr); 174735d7f90SBarry Smith ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 175735d7f90SBarry Smith 176735d7f90SBarry Smith ierr = MatZeroEntries(A);CHKERRQ(ierr); 177735d7f90SBarry Smith ierr = FillMatrixCPUCOO3d(&fe,A);CHKERRQ(ierr); 178735d7f90SBarry Smith ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 179735d7f90SBarry Smith 180735d7f90SBarry Smith ierr = MatZeroEntries(A);CHKERRQ(ierr); 181735d7f90SBarry Smith ierr = PetscObjectBaseTypeCompare((PetscObject)A,MATSEQAIJKOKKOS,&is_kokkos);CHKERRQ(ierr); 182735d7f90SBarry Smith ierr = PetscObjectBaseTypeCompare((PetscObject)A,MATSEQAIJCUSPARSE,&is_cuda);CHKERRQ(ierr); 183735d7f90SBarry Smith #if defined(PETSC_HAVE_KOKKOS) 184735d7f90SBarry Smith if (is_kokkos) {ierr = FillMatrixKokkosCOO(&fe,A);CHKERRQ(ierr);} 185735d7f90SBarry Smith #endif 186735d7f90SBarry Smith #if defined(PETSC_HAVE_CUDA) 187735d7f90SBarry Smith if (is_cuda) {ierr = FillMatrixCUDACOO(&fe,A);CHKERRQ(ierr);} 188735d7f90SBarry Smith #endif 189735d7f90SBarry Smith ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 190735d7f90SBarry Smith 191735d7f90SBarry Smith ierr = MatDestroy(&A);CHKERRQ(ierr); 192735d7f90SBarry Smith ierr = DestroyFEStruct(&fe);CHKERRQ(ierr); 193735d7f90SBarry Smith ierr = PetscFinalize(); 194735d7f90SBarry Smith return ierr; 195735d7f90SBarry Smith } 196735d7f90SBarry Smith 197735d7f90SBarry Smith /*TEST 198735d7f90SBarry Smith build: 199735d7f90SBarry Smith requires: cuda kokkos_kernels 200735d7f90SBarry Smith depends: ex18cu.cu ex18kok.kokkos.cxx 201735d7f90SBarry Smith 202735d7f90SBarry Smith testset: 203735d7f90SBarry Smith filter: grep -v "type" 204735d7f90SBarry Smith output_file: output/ex18_1.out 205735d7f90SBarry Smith 206735d7f90SBarry Smith test: 207735d7f90SBarry Smith suffix: kok 208735d7f90SBarry Smith requires: kokkos_kernels 209735d7f90SBarry Smith args: -mat_type aijkokkos 210735d7f90SBarry Smith 211735d7f90SBarry Smith test: 212735d7f90SBarry Smith suffix: cuda 213735d7f90SBarry Smith requires: cuda 214735d7f90SBarry Smith args: -mat_type aijcusparse 215735d7f90SBarry Smith 216735d7f90SBarry Smith TEST*/ 217