1*0e2ec84fSDave May 2*0e2ec84fSDave May #include <petsc.h> 3*0e2ec84fSDave May #include <petscdm.h> 4*0e2ec84fSDave May #include <petscdmplex.h> 5*0e2ec84fSDave May #include <petscdmswarm.h> 6*0e2ec84fSDave May 7*0e2ec84fSDave May 8*0e2ec84fSDave May PetscErrorCode subdivide_triangle(PetscReal v1[2],PetscReal v2[2],PetscReal v3[3],PetscInt depth,PetscInt max,PetscReal xi[],PetscInt *np) 9*0e2ec84fSDave May { 10*0e2ec84fSDave May PetscReal v12[2],v23[2],v31[2]; 11*0e2ec84fSDave May PetscInt i; 12*0e2ec84fSDave May PetscErrorCode ierr; 13*0e2ec84fSDave May 14*0e2ec84fSDave May PetscFunctionBegin; 15*0e2ec84fSDave May if (depth == max) { 16*0e2ec84fSDave May PetscReal cx[2]; 17*0e2ec84fSDave May 18*0e2ec84fSDave May //printf("centroid \n"); 19*0e2ec84fSDave May cx[0] = (v1[0] + v2[0] + v3[0])/3.0; 20*0e2ec84fSDave May cx[1] = (v1[1] + v2[1] + v3[1])/3.0; 21*0e2ec84fSDave May 22*0e2ec84fSDave May xi[2*(*np)+0] = cx[0]; 23*0e2ec84fSDave May xi[2*(*np)+1] = cx[1]; 24*0e2ec84fSDave May (*np)++; 25*0e2ec84fSDave May PetscFunctionReturn(0); 26*0e2ec84fSDave May } 27*0e2ec84fSDave May 28*0e2ec84fSDave May /* calculate midpoints of each side */ 29*0e2ec84fSDave May for (i = 0; i < 2; i++) { 30*0e2ec84fSDave May v12[i] = (v1[i]+v2[i])/2.0; 31*0e2ec84fSDave May v23[i] = (v2[i]+v3[i])/2.0; 32*0e2ec84fSDave May v31[i] = (v3[i]+v1[i])/2.0; 33*0e2ec84fSDave May } 34*0e2ec84fSDave May 35*0e2ec84fSDave May /* recursively subdivide new triangles */ 36*0e2ec84fSDave May ierr = subdivide_triangle(v1,v12,v31,depth+1,max,xi,np);CHKERRQ(ierr); 37*0e2ec84fSDave May ierr = subdivide_triangle(v2,v23,v12,depth+1,max,xi,np);CHKERRQ(ierr); 38*0e2ec84fSDave May ierr = subdivide_triangle(v3,v31,v23,depth+1,max,xi,np);CHKERRQ(ierr); 39*0e2ec84fSDave May ierr = subdivide_triangle(v12,v23,v31,depth+1,max,xi,np);CHKERRQ(ierr); 40*0e2ec84fSDave May PetscFunctionReturn(0); 41*0e2ec84fSDave May } 42*0e2ec84fSDave May 43*0e2ec84fSDave May 44*0e2ec84fSDave May PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX2D_SubDivide(DM dm,DM dmc,PetscInt nsub) 45*0e2ec84fSDave May { 46*0e2ec84fSDave May PetscErrorCode ierr; 47*0e2ec84fSDave May const PetscInt dim = 2; 48*0e2ec84fSDave May PetscInt q,npoints_q,e,nel,npe,pcnt,ps,pe,d,k,depth; 49*0e2ec84fSDave May PetscReal *xi; 50*0e2ec84fSDave May PetscReal **basis; 51*0e2ec84fSDave May Vec coorlocal; 52*0e2ec84fSDave May PetscSection coordSection; 53*0e2ec84fSDave May PetscScalar *elcoor = NULL; 54*0e2ec84fSDave May PetscReal *swarm_coor; 55*0e2ec84fSDave May PetscInt *swarm_cellid; 56*0e2ec84fSDave May PetscReal v1[2],v2[2],v3[2]; 57*0e2ec84fSDave May 58*0e2ec84fSDave May PetscFunctionBegin; 59*0e2ec84fSDave May 60*0e2ec84fSDave May npoints_q = 1; 61*0e2ec84fSDave May for (d=0; d<nsub; d++) { npoints_q *= 4; } 62*0e2ec84fSDave May ierr = PetscMalloc1(dim*npoints_q,&xi);CHKERRQ(ierr); 63*0e2ec84fSDave May 64*0e2ec84fSDave May v1[0] = 0.0; v1[1] = 0.0; 65*0e2ec84fSDave May v2[0] = 1.0; v2[1] = 0.0; 66*0e2ec84fSDave May v3[0] = 0.0; v3[1] = 1.0; 67*0e2ec84fSDave May depth = 0; 68*0e2ec84fSDave May pcnt = 0; 69*0e2ec84fSDave May ierr = subdivide_triangle(v1,v2,v3,depth,nsub,xi,&pcnt);CHKERRQ(ierr); 70*0e2ec84fSDave May 71*0e2ec84fSDave May npe = 3; /* nodes per element (triangle) */ 72*0e2ec84fSDave May ierr = PetscMalloc1(npoints_q,&basis);CHKERRQ(ierr); 73*0e2ec84fSDave May for (q=0; q<npoints_q; q++) { 74*0e2ec84fSDave May ierr = PetscMalloc1(npe,&basis[q]);CHKERRQ(ierr); 75*0e2ec84fSDave May 76*0e2ec84fSDave May basis[q][0] = 1.0 - xi[dim*q+0] - xi[dim*q+1]; 77*0e2ec84fSDave May basis[q][1] = xi[dim*q+0]; 78*0e2ec84fSDave May basis[q][2] = xi[dim*q+1]; 79*0e2ec84fSDave May } 80*0e2ec84fSDave May 81*0e2ec84fSDave May /* 0->cell, 1->edge, 2->vert */ 82*0e2ec84fSDave May ierr = DMPlexGetHeightStratum(dmc,0,&ps,&pe);CHKERRQ(ierr); 83*0e2ec84fSDave May nel = pe - ps; 84*0e2ec84fSDave May 85*0e2ec84fSDave May ierr = DMSwarmSetLocalSizes(dm,npoints_q*nel,-1);CHKERRQ(ierr); 86*0e2ec84fSDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 87*0e2ec84fSDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 88*0e2ec84fSDave May 89*0e2ec84fSDave May ierr = DMGetCoordinatesLocal(dmc,&coorlocal);CHKERRQ(ierr); 90*0e2ec84fSDave May ierr = DMGetCoordinateSection(dmc,&coordSection);CHKERRQ(ierr); 91*0e2ec84fSDave May 92*0e2ec84fSDave May pcnt = 0; 93*0e2ec84fSDave May for (e=0; e<nel; e++) { 94*0e2ec84fSDave May ierr = DMPlexVecGetClosure(dmc,coordSection,coorlocal,e,NULL,&elcoor);CHKERRQ(ierr); 95*0e2ec84fSDave May 96*0e2ec84fSDave May for (q=0; q<npoints_q; q++) { 97*0e2ec84fSDave May for (d=0; d<dim; d++) { 98*0e2ec84fSDave May swarm_coor[dim*pcnt+d] = 0.0; 99*0e2ec84fSDave May for (k=0; k<npe; k++) { 100*0e2ec84fSDave May swarm_coor[dim*pcnt+d] += basis[q][k] * elcoor[dim*k+d]; 101*0e2ec84fSDave May } 102*0e2ec84fSDave May } 103*0e2ec84fSDave May swarm_cellid[pcnt] = e; 104*0e2ec84fSDave May pcnt++; 105*0e2ec84fSDave May } 106*0e2ec84fSDave May ierr = DMPlexVecRestoreClosure(dmc,coordSection,coorlocal,e,NULL,&elcoor);CHKERRQ(ierr); 107*0e2ec84fSDave May } 108*0e2ec84fSDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 109*0e2ec84fSDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 110*0e2ec84fSDave May 111*0e2ec84fSDave May ierr = PetscFree(xi);CHKERRQ(ierr); 112*0e2ec84fSDave May for (q=0; q<npoints_q; q++) { 113*0e2ec84fSDave May ierr = PetscFree(basis[q]);CHKERRQ(ierr); 114*0e2ec84fSDave May } 115*0e2ec84fSDave May ierr = PetscFree(basis);CHKERRQ(ierr); 116*0e2ec84fSDave May 117*0e2ec84fSDave May PetscFunctionReturn(0); 118*0e2ec84fSDave May } 119*0e2ec84fSDave May 120*0e2ec84fSDave May PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX2D_Regular(DM dm,DM dmc,PetscInt npoints) 121*0e2ec84fSDave May { 122*0e2ec84fSDave May PetscErrorCode ierr; 123*0e2ec84fSDave May const PetscInt dim = 2; 124*0e2ec84fSDave May PetscInt ii,jj,q,npoints_q,e,nel,npe,pcnt,ps,pe,d,k; 125*0e2ec84fSDave May PetscReal *xi,ds,ds2; 126*0e2ec84fSDave May PetscReal **basis; 127*0e2ec84fSDave May Vec coorlocal; 128*0e2ec84fSDave May PetscSection coordSection; 129*0e2ec84fSDave May PetscScalar *elcoor = NULL; 130*0e2ec84fSDave May PetscReal *swarm_coor; 131*0e2ec84fSDave May PetscInt *swarm_cellid; 132*0e2ec84fSDave May 133*0e2ec84fSDave May PetscFunctionBegin; 134*0e2ec84fSDave May ierr = PetscMalloc1(dim*npoints*npoints,&xi);CHKERRQ(ierr); 135*0e2ec84fSDave May pcnt = 0; 136*0e2ec84fSDave May ds = 1.0/((PetscReal)(npoints-1)); 137*0e2ec84fSDave May ds2 = 1.0/((PetscReal)(npoints)); 138*0e2ec84fSDave May for (jj = 0; jj<npoints; jj++) { 139*0e2ec84fSDave May for (ii=0; ii<npoints-jj; ii++) { 140*0e2ec84fSDave May xi[dim*pcnt+0] = ii * ds; 141*0e2ec84fSDave May xi[dim*pcnt+1] = jj * ds; 142*0e2ec84fSDave May 143*0e2ec84fSDave May xi[dim*pcnt+0] *= (1.0 - 1.2*ds2); 144*0e2ec84fSDave May xi[dim*pcnt+1] *= (1.0 - 1.2*ds2); 145*0e2ec84fSDave May 146*0e2ec84fSDave May xi[dim*pcnt+0] += 0.35*ds2; 147*0e2ec84fSDave May xi[dim*pcnt+1] += 0.35*ds2; 148*0e2ec84fSDave May 149*0e2ec84fSDave May pcnt++; 150*0e2ec84fSDave May } 151*0e2ec84fSDave May } 152*0e2ec84fSDave May npoints_q = pcnt; 153*0e2ec84fSDave May 154*0e2ec84fSDave May npe = 3; /* nodes per element (triangle) */ 155*0e2ec84fSDave May ierr = PetscMalloc1(npoints_q,&basis);CHKERRQ(ierr); 156*0e2ec84fSDave May for (q=0; q<npoints_q; q++) { 157*0e2ec84fSDave May ierr = PetscMalloc1(npe,&basis[q]);CHKERRQ(ierr); 158*0e2ec84fSDave May 159*0e2ec84fSDave May basis[q][0] = 1.0 - xi[dim*q+0] - xi[dim*q+1]; 160*0e2ec84fSDave May basis[q][1] = xi[dim*q+0]; 161*0e2ec84fSDave May basis[q][2] = xi[dim*q+1]; 162*0e2ec84fSDave May } 163*0e2ec84fSDave May 164*0e2ec84fSDave May /* 0->cell, 1->edge, 2->vert */ 165*0e2ec84fSDave May ierr = DMPlexGetHeightStratum(dmc,0,&ps,&pe);CHKERRQ(ierr); 166*0e2ec84fSDave May nel = pe - ps; 167*0e2ec84fSDave May 168*0e2ec84fSDave May ierr = DMSwarmSetLocalSizes(dm,npoints_q*nel,-1);CHKERRQ(ierr); 169*0e2ec84fSDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 170*0e2ec84fSDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 171*0e2ec84fSDave May 172*0e2ec84fSDave May ierr = DMGetCoordinatesLocal(dmc,&coorlocal);CHKERRQ(ierr); 173*0e2ec84fSDave May ierr = DMGetCoordinateSection(dmc,&coordSection);CHKERRQ(ierr); 174*0e2ec84fSDave May 175*0e2ec84fSDave May pcnt = 0; 176*0e2ec84fSDave May for (e=0; e<nel; e++) { 177*0e2ec84fSDave May ierr = DMPlexVecGetClosure(dmc,coordSection,coorlocal,e,NULL,&elcoor);CHKERRQ(ierr); 178*0e2ec84fSDave May 179*0e2ec84fSDave May for (q=0; q<npoints_q; q++) { 180*0e2ec84fSDave May for (d=0; d<dim; d++) { 181*0e2ec84fSDave May swarm_coor[dim*pcnt+d] = 0.0; 182*0e2ec84fSDave May for (k=0; k<npe; k++) { 183*0e2ec84fSDave May swarm_coor[dim*pcnt+d] += basis[q][k] * elcoor[dim*k+d]; 184*0e2ec84fSDave May } 185*0e2ec84fSDave May } 186*0e2ec84fSDave May swarm_cellid[pcnt] = e; 187*0e2ec84fSDave May pcnt++; 188*0e2ec84fSDave May } 189*0e2ec84fSDave May ierr = DMPlexVecRestoreClosure(dmc,coordSection,coorlocal,e,NULL,&elcoor);CHKERRQ(ierr); 190*0e2ec84fSDave May } 191*0e2ec84fSDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 192*0e2ec84fSDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 193*0e2ec84fSDave May 194*0e2ec84fSDave May ierr = PetscFree(xi);CHKERRQ(ierr); 195*0e2ec84fSDave May for (q=0; q<npoints_q; q++) { 196*0e2ec84fSDave May ierr = PetscFree(basis[q]);CHKERRQ(ierr); 197*0e2ec84fSDave May } 198*0e2ec84fSDave May ierr = PetscFree(basis);CHKERRQ(ierr); 199*0e2ec84fSDave May 200*0e2ec84fSDave May PetscFunctionReturn(0); 201*0e2ec84fSDave May } 202*0e2ec84fSDave May 203*0e2ec84fSDave May PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX(DM dm,DM celldm,DMSwarmPICLayoutType layout,PetscInt layout_param) 204*0e2ec84fSDave May { 205*0e2ec84fSDave May PetscErrorCode ierr; 206*0e2ec84fSDave May PetscInt dim; 207*0e2ec84fSDave May 208*0e2ec84fSDave May PetscFunctionBegin; 209*0e2ec84fSDave May ierr = DMGetDimension(celldm,&dim);CHKERRQ(ierr); 210*0e2ec84fSDave May if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"No 3D support for PLEX"); 211*0e2ec84fSDave May switch (layout) { 212*0e2ec84fSDave May case DMSWARMPIC_LAYOUT_REGULAR: 213*0e2ec84fSDave May ierr = private_DMSwarmInsertPointsUsingCellDM_PLEX2D_Regular(dm,celldm,layout_param);CHKERRQ(ierr); 214*0e2ec84fSDave May break; 215*0e2ec84fSDave May case DMSWARMPIC_LAYOUT_GAUSS: 216*0e2ec84fSDave May SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Gauss layout not supported for PLEX"); 217*0e2ec84fSDave May break; 218*0e2ec84fSDave May case DMSWARMPIC_LAYOUT_SUBDIVISION: 219*0e2ec84fSDave May ierr = private_DMSwarmInsertPointsUsingCellDM_PLEX2D_SubDivide(dm,celldm,layout_param);CHKERRQ(ierr); 220*0e2ec84fSDave May break; 221*0e2ec84fSDave May } 222*0e2ec84fSDave May PetscFunctionReturn(0); 223*0e2ec84fSDave May } 224