1*0e2ec84fSDave May 2*0e2ec84fSDave May #include <petsc.h> 3*0e2ec84fSDave May #include <petscdm.h> 4*0e2ec84fSDave May #include <petscdmda.h> 5*0e2ec84fSDave May #include <petscdmswarm.h> 6*0e2ec84fSDave May 7*0e2ec84fSDave May PetscErrorCode private_DMSwarmCreateCellLocalCoords_DA_Q1_Regular(PetscInt dim,PetscInt np[],PetscInt *_npoints,PetscReal **_xi) 8*0e2ec84fSDave May { 9*0e2ec84fSDave May PetscErrorCode ierr; 10*0e2ec84fSDave May PetscReal *xi; 11*0e2ec84fSDave May PetscInt d,npoints,cnt; 12*0e2ec84fSDave May PetscReal ds[] = {0.0,0.0,0.0}; 13*0e2ec84fSDave May PetscInt ii,jj,kk; 14*0e2ec84fSDave May 15*0e2ec84fSDave May PetscFunctionBegin; 16*0e2ec84fSDave May switch (dim) { 17*0e2ec84fSDave May case 1: 18*0e2ec84fSDave May npoints = np[0]; 19*0e2ec84fSDave May break; 20*0e2ec84fSDave May case 2: 21*0e2ec84fSDave May npoints = np[0]*np[1]; 22*0e2ec84fSDave May break; 23*0e2ec84fSDave May case 3: 24*0e2ec84fSDave May npoints = np[0]*np[1]*np[2]; 25*0e2ec84fSDave May break; 26*0e2ec84fSDave May } 27*0e2ec84fSDave May for (d=0; d<dim; d++) { 28*0e2ec84fSDave May ds[d] = 2.0 / ((PetscReal)np[d]); 29*0e2ec84fSDave May } 30*0e2ec84fSDave May 31*0e2ec84fSDave May ierr = PetscMalloc1(dim*npoints,&xi);CHKERRQ(ierr); 32*0e2ec84fSDave May 33*0e2ec84fSDave May switch (dim) { 34*0e2ec84fSDave May case 1: 35*0e2ec84fSDave May cnt = 0; 36*0e2ec84fSDave May for (ii=0; ii<np[0]; ii++) { 37*0e2ec84fSDave May xi[dim*cnt+0] = -1.0 + 0.5*ds[d] + ii*ds[0]; 38*0e2ec84fSDave May cnt++; 39*0e2ec84fSDave May } 40*0e2ec84fSDave May break; 41*0e2ec84fSDave May 42*0e2ec84fSDave May case 2: 43*0e2ec84fSDave May cnt = 0; 44*0e2ec84fSDave May for (jj=0; jj<np[1]; jj++) { 45*0e2ec84fSDave May for (ii=0; ii<np[0]; ii++) { 46*0e2ec84fSDave May xi[dim*cnt+0] = -1.0 + 0.5*ds[0] + ii*ds[0]; 47*0e2ec84fSDave May xi[dim*cnt+1] = -1.0 + 0.5*ds[1] + jj*ds[1]; 48*0e2ec84fSDave May cnt++; 49*0e2ec84fSDave May } 50*0e2ec84fSDave May } 51*0e2ec84fSDave May break; 52*0e2ec84fSDave May 53*0e2ec84fSDave May case 3: 54*0e2ec84fSDave May cnt = 0; 55*0e2ec84fSDave May for (kk=0; kk<np[2]; kk++) { 56*0e2ec84fSDave May for (jj=0; jj<np[1]; jj++) { 57*0e2ec84fSDave May for (ii=0; ii<np[0]; ii++) { 58*0e2ec84fSDave May xi[dim*cnt+0] = -1.0 + 0.5*ds[0] + ii*ds[0]; 59*0e2ec84fSDave May xi[dim*cnt+1] = -1.0 + 0.5*ds[1] + jj*ds[1]; 60*0e2ec84fSDave May xi[dim*cnt+2] = -1.0 + 0.5*ds[2] + kk*ds[2]; 61*0e2ec84fSDave May cnt++; 62*0e2ec84fSDave May } 63*0e2ec84fSDave May } 64*0e2ec84fSDave May } 65*0e2ec84fSDave May break; 66*0e2ec84fSDave May } 67*0e2ec84fSDave May 68*0e2ec84fSDave May *_npoints = npoints; 69*0e2ec84fSDave May *_xi = xi; 70*0e2ec84fSDave May 71*0e2ec84fSDave May PetscFunctionReturn(0); 72*0e2ec84fSDave May } 73*0e2ec84fSDave May 74*0e2ec84fSDave May PetscErrorCode private_DMSwarmCreateCellLocalCoords_DA_Q1_Gauss(PetscInt dim,PetscInt np_1d,PetscInt *_npoints,PetscReal **_xi) 75*0e2ec84fSDave May { 76*0e2ec84fSDave May PetscErrorCode ierr; 77*0e2ec84fSDave May PetscQuadrature quadrature; 78*0e2ec84fSDave May const PetscReal *quadrature_xi; 79*0e2ec84fSDave May PetscReal *xi; 80*0e2ec84fSDave May PetscInt d,q,npoints_q; 81*0e2ec84fSDave May 82*0e2ec84fSDave May PetscFunctionBegin; 83*0e2ec84fSDave May ierr = PetscDTGaussTensorQuadrature(dim,np_1d,-1.0,1.0,&quadrature);CHKERRQ(ierr); 84*0e2ec84fSDave May ierr = PetscQuadratureGetData(quadrature,NULL,&npoints_q,&quadrature_xi,NULL);CHKERRQ(ierr); 85*0e2ec84fSDave May ierr = PetscMalloc1(dim*npoints_q,&xi);CHKERRQ(ierr); 86*0e2ec84fSDave May for (q=0; q<npoints_q; q++) { 87*0e2ec84fSDave May for (d=0; d<dim; d++) { 88*0e2ec84fSDave May xi[dim*q+d] = quadrature_xi[dim*q+d]; 89*0e2ec84fSDave May } 90*0e2ec84fSDave May } 91*0e2ec84fSDave May ierr = PetscQuadratureDestroy(&quadrature);CHKERRQ(ierr); 92*0e2ec84fSDave May 93*0e2ec84fSDave May *_npoints = npoints_q; 94*0e2ec84fSDave May *_xi = xi; 95*0e2ec84fSDave May 96*0e2ec84fSDave May PetscFunctionReturn(0); 97*0e2ec84fSDave May } 98*0e2ec84fSDave May 99*0e2ec84fSDave May PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_DA_Q1(DM dm,DM dmc,PetscInt npoints,PetscInt layout) 100*0e2ec84fSDave May { 101*0e2ec84fSDave May PetscErrorCode ierr; 102*0e2ec84fSDave May PetscInt dim,npoints_q; 103*0e2ec84fSDave May PetscInt nel,npe,e,q,k,d; 104*0e2ec84fSDave May const PetscInt *element_list; 105*0e2ec84fSDave May PetscReal **basis; 106*0e2ec84fSDave May PetscReal *xi; 107*0e2ec84fSDave May Vec coor; 108*0e2ec84fSDave May const PetscScalar *_coor; 109*0e2ec84fSDave May PetscReal *elcoor; 110*0e2ec84fSDave May PetscReal *swarm_coor; 111*0e2ec84fSDave May PetscInt *swarm_cellid; 112*0e2ec84fSDave May PetscInt pcnt; 113*0e2ec84fSDave May 114*0e2ec84fSDave May PetscFunctionBegin; 115*0e2ec84fSDave May ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 116*0e2ec84fSDave May switch (layout) { 117*0e2ec84fSDave May case DMSWARMPIC_LAYOUT_REGULAR: 118*0e2ec84fSDave May { 119*0e2ec84fSDave May PetscInt np_dir[3]; 120*0e2ec84fSDave May np_dir[0] = np_dir[1] = np_dir[2] = npoints; 121*0e2ec84fSDave May ierr = private_DMSwarmCreateCellLocalCoords_DA_Q1_Regular(dim,np_dir,&npoints_q,&xi);CHKERRQ(ierr); 122*0e2ec84fSDave May } 123*0e2ec84fSDave May break; 124*0e2ec84fSDave May case DMSWARMPIC_LAYOUT_GAUSS: 125*0e2ec84fSDave May ierr = private_DMSwarmCreateCellLocalCoords_DA_Q1_Gauss(dim,npoints,&npoints_q,&xi);CHKERRQ(ierr); 126*0e2ec84fSDave May break; 127*0e2ec84fSDave May default: 128*0e2ec84fSDave May SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"A valid DMSwarmPIC layout must be provided"); 129*0e2ec84fSDave May break; 130*0e2ec84fSDave May } 131*0e2ec84fSDave May 132*0e2ec84fSDave May ierr = DMDAGetElements(dmc,&nel,&npe,&element_list);CHKERRQ(ierr); 133*0e2ec84fSDave May 134*0e2ec84fSDave May ierr = PetscMalloc1(dim*npe,&elcoor);CHKERRQ(ierr); 135*0e2ec84fSDave May ierr = PetscMalloc1(npoints_q,&basis);CHKERRQ(ierr); 136*0e2ec84fSDave May for (q=0; q<npoints_q; q++) { 137*0e2ec84fSDave May ierr = PetscMalloc1(npe,&basis[q]);CHKERRQ(ierr); 138*0e2ec84fSDave May 139*0e2ec84fSDave May switch (dim) { 140*0e2ec84fSDave May case 1: 141*0e2ec84fSDave May basis[q][0] = 0.5*(1.0 - xi[dim*q+0]); 142*0e2ec84fSDave May basis[q][1] = 0.5*(1.0 + xi[dim*q+0]); 143*0e2ec84fSDave May break; 144*0e2ec84fSDave May case 2: 145*0e2ec84fSDave May basis[q][0] = 0.25*(1.0 - xi[dim*q+0])*(1.0 - xi[dim*q+1]); 146*0e2ec84fSDave May basis[q][1] = 0.25*(1.0 + xi[dim*q+0])*(1.0 - xi[dim*q+1]); 147*0e2ec84fSDave May basis[q][2] = 0.25*(1.0 + xi[dim*q+0])*(1.0 + xi[dim*q+1]); 148*0e2ec84fSDave May basis[q][3] = 0.25*(1.0 - xi[dim*q+0])*(1.0 + xi[dim*q+1]); 149*0e2ec84fSDave May break; 150*0e2ec84fSDave May 151*0e2ec84fSDave May case 3: 152*0e2ec84fSDave May basis[q][0] = 0.125*(1.0 - xi[dim*q+0])*(1.0 - xi[dim*q+1])*(1.0 - xi[dim*q+2]); 153*0e2ec84fSDave May basis[q][1] = 0.125*(1.0 + xi[dim*q+0])*(1.0 - xi[dim*q+1])*(1.0 - xi[dim*q+2]); 154*0e2ec84fSDave May basis[q][2] = 0.125*(1.0 + xi[dim*q+0])*(1.0 + xi[dim*q+1])*(1.0 - xi[dim*q+2]); 155*0e2ec84fSDave May basis[q][3] = 0.125*(1.0 - xi[dim*q+0])*(1.0 + xi[dim*q+1])*(1.0 - xi[dim*q+2]); 156*0e2ec84fSDave May basis[q][4] = 0.125*(1.0 - xi[dim*q+0])*(1.0 - xi[dim*q+1])*(1.0 + xi[dim*q+2]); 157*0e2ec84fSDave May basis[q][5] = 0.125*(1.0 + xi[dim*q+0])*(1.0 - xi[dim*q+1])*(1.0 + xi[dim*q+2]); 158*0e2ec84fSDave May basis[q][6] = 0.125*(1.0 + xi[dim*q+0])*(1.0 + xi[dim*q+1])*(1.0 + xi[dim*q+2]); 159*0e2ec84fSDave May basis[q][7] = 0.125*(1.0 - xi[dim*q+0])*(1.0 + xi[dim*q+1])*(1.0 + xi[dim*q+2]); 160*0e2ec84fSDave May break; 161*0e2ec84fSDave May } 162*0e2ec84fSDave May } 163*0e2ec84fSDave May 164*0e2ec84fSDave May ierr = DMSwarmSetLocalSizes(dm,npoints_q*nel,-1);CHKERRQ(ierr); 165*0e2ec84fSDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 166*0e2ec84fSDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 167*0e2ec84fSDave May 168*0e2ec84fSDave May ierr = DMGetCoordinatesLocal(dmc,&coor);CHKERRQ(ierr); 169*0e2ec84fSDave May ierr = VecGetArrayRead(coor,&_coor);CHKERRQ(ierr); 170*0e2ec84fSDave May pcnt = 0; 171*0e2ec84fSDave May for (e=0; e<nel; e++) { 172*0e2ec84fSDave May const PetscInt *element = &element_list[npe*e]; 173*0e2ec84fSDave May 174*0e2ec84fSDave May for (k=0; k<npe; k++) { 175*0e2ec84fSDave May for (d=0; d<dim; d++) { 176*0e2ec84fSDave May elcoor[dim*k+d] = _coor[ dim*element[k] + d ]; 177*0e2ec84fSDave May } 178*0e2ec84fSDave May } 179*0e2ec84fSDave May 180*0e2ec84fSDave May for (q=0; q<npoints_q; q++) { 181*0e2ec84fSDave May for (d=0; d<dim; d++) { 182*0e2ec84fSDave May swarm_coor[dim*pcnt+d] = 0.0; 183*0e2ec84fSDave May } 184*0e2ec84fSDave May for (k=0; k<npe; k++) { 185*0e2ec84fSDave May for (d=0; d<dim; d++) { 186*0e2ec84fSDave May swarm_coor[dim*pcnt+d] += basis[q][k] * elcoor[dim*k+d]; 187*0e2ec84fSDave May } 188*0e2ec84fSDave May } 189*0e2ec84fSDave May swarm_cellid[pcnt] = e; 190*0e2ec84fSDave May pcnt++; 191*0e2ec84fSDave May } 192*0e2ec84fSDave May } 193*0e2ec84fSDave May ierr = VecRestoreArrayRead(coor,&_coor);CHKERRQ(ierr); 194*0e2ec84fSDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 195*0e2ec84fSDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 196*0e2ec84fSDave May ierr = DMDARestoreElements(dmc,&nel,&npe,&element_list);CHKERRQ(ierr); 197*0e2ec84fSDave May 198*0e2ec84fSDave May ierr = PetscFree(xi);CHKERRQ(ierr); 199*0e2ec84fSDave May ierr = PetscFree(elcoor);CHKERRQ(ierr); 200*0e2ec84fSDave May for (q=0; q<npoints_q; q++) { 201*0e2ec84fSDave May ierr = PetscFree(basis[q]);CHKERRQ(ierr); 202*0e2ec84fSDave May } 203*0e2ec84fSDave May ierr = PetscFree(basis);CHKERRQ(ierr); 204*0e2ec84fSDave May 205*0e2ec84fSDave May PetscFunctionReturn(0); 206*0e2ec84fSDave May } 207*0e2ec84fSDave May 208*0e2ec84fSDave May PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_DA(DM dm,DM celldm,DMSwarmPICLayoutType layout,PetscInt layout_param) 209*0e2ec84fSDave May { 210*0e2ec84fSDave May PetscErrorCode ierr; 211*0e2ec84fSDave May DMDAElementType etype; 212*0e2ec84fSDave May PetscInt dim; 213*0e2ec84fSDave May 214*0e2ec84fSDave May PetscFunctionBegin; 215*0e2ec84fSDave May ierr = DMDAGetElementType(celldm,&etype);CHKERRQ(ierr); 216*0e2ec84fSDave May ierr = DMGetDimension(celldm,&dim);CHKERRQ(ierr); 217*0e2ec84fSDave May switch (etype) { 218*0e2ec84fSDave May case DMDA_ELEMENT_P1: 219*0e2ec84fSDave May SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"DA support is not currently available for DMDA_ELEMENT_P1"); 220*0e2ec84fSDave May break; 221*0e2ec84fSDave May case DMDA_ELEMENT_Q1: 222*0e2ec84fSDave May if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Support only available for dim = 2, 3"); 223*0e2ec84fSDave May ierr = private_DMSwarmInsertPointsUsingCellDM_DA_Q1(dm,celldm,layout_param,layout);CHKERRQ(ierr); 224*0e2ec84fSDave May break; 225*0e2ec84fSDave May } 226*0e2ec84fSDave May PetscFunctionReturn(0); 227*0e2ec84fSDave May } 228