10e2ec84fSDave May #include <petscdm.h> 20e2ec84fSDave May #include <petscdmda.h> 30e2ec84fSDave May #include <petscdmswarm.h> 41f52046fSDave May #include <petsc/private/dmswarmimpl.h> 51f52046fSDave May #include "data_bucket.h" 60e2ec84fSDave May 70e2ec84fSDave May PetscErrorCode private_DMSwarmCreateCellLocalCoords_DA_Q1_Regular(PetscInt dim,PetscInt np[],PetscInt *_npoints,PetscReal **_xi) 80e2ec84fSDave May { 90e2ec84fSDave May PetscErrorCode ierr; 100e2ec84fSDave May PetscReal *xi; 11f9368e01SSatish Balay PetscInt d,npoints=0,cnt; 120e2ec84fSDave May PetscReal ds[] = {0.0,0.0,0.0}; 130e2ec84fSDave May PetscInt ii,jj,kk; 140e2ec84fSDave May 150e2ec84fSDave May PetscFunctionBegin; 160e2ec84fSDave May switch (dim) { 170e2ec84fSDave May case 1: 180e2ec84fSDave May npoints = np[0]; 190e2ec84fSDave May break; 200e2ec84fSDave May case 2: 210e2ec84fSDave May npoints = np[0]*np[1]; 220e2ec84fSDave May break; 230e2ec84fSDave May case 3: 240e2ec84fSDave May npoints = np[0]*np[1]*np[2]; 250e2ec84fSDave May break; 260e2ec84fSDave May } 270e2ec84fSDave May for (d=0; d<dim; d++) { 280e2ec84fSDave May ds[d] = 2.0 / ((PetscReal)np[d]); 290e2ec84fSDave May } 300e2ec84fSDave May 310e2ec84fSDave May ierr = PetscMalloc1(dim*npoints,&xi);CHKERRQ(ierr); 320e2ec84fSDave May 330e2ec84fSDave May switch (dim) { 340e2ec84fSDave May case 1: 350e2ec84fSDave May cnt = 0; 360e2ec84fSDave May for (ii=0; ii<np[0]; ii++) { 370e2ec84fSDave May xi[dim*cnt+0] = -1.0 + 0.5*ds[d] + ii*ds[0]; 380e2ec84fSDave May cnt++; 390e2ec84fSDave May } 400e2ec84fSDave May break; 410e2ec84fSDave May 420e2ec84fSDave May case 2: 430e2ec84fSDave May cnt = 0; 440e2ec84fSDave May for (jj=0; jj<np[1]; jj++) { 450e2ec84fSDave May for (ii=0; ii<np[0]; ii++) { 460e2ec84fSDave May xi[dim*cnt+0] = -1.0 + 0.5*ds[0] + ii*ds[0]; 470e2ec84fSDave May xi[dim*cnt+1] = -1.0 + 0.5*ds[1] + jj*ds[1]; 480e2ec84fSDave May cnt++; 490e2ec84fSDave May } 500e2ec84fSDave May } 510e2ec84fSDave May break; 520e2ec84fSDave May 530e2ec84fSDave May case 3: 540e2ec84fSDave May cnt = 0; 550e2ec84fSDave May for (kk=0; kk<np[2]; kk++) { 560e2ec84fSDave May for (jj=0; jj<np[1]; jj++) { 570e2ec84fSDave May for (ii=0; ii<np[0]; ii++) { 580e2ec84fSDave May xi[dim*cnt+0] = -1.0 + 0.5*ds[0] + ii*ds[0]; 590e2ec84fSDave May xi[dim*cnt+1] = -1.0 + 0.5*ds[1] + jj*ds[1]; 600e2ec84fSDave May xi[dim*cnt+2] = -1.0 + 0.5*ds[2] + kk*ds[2]; 610e2ec84fSDave May cnt++; 620e2ec84fSDave May } 630e2ec84fSDave May } 640e2ec84fSDave May } 650e2ec84fSDave May break; 660e2ec84fSDave May } 670e2ec84fSDave May 680e2ec84fSDave May *_npoints = npoints; 690e2ec84fSDave May *_xi = xi; 700e2ec84fSDave May 710e2ec84fSDave May PetscFunctionReturn(0); 720e2ec84fSDave May } 730e2ec84fSDave May 740e2ec84fSDave May PetscErrorCode private_DMSwarmCreateCellLocalCoords_DA_Q1_Gauss(PetscInt dim,PetscInt np_1d,PetscInt *_npoints,PetscReal **_xi) 750e2ec84fSDave May { 760e2ec84fSDave May PetscErrorCode ierr; 770e2ec84fSDave May PetscQuadrature quadrature; 780e2ec84fSDave May const PetscReal *quadrature_xi; 790e2ec84fSDave May PetscReal *xi; 800e2ec84fSDave May PetscInt d,q,npoints_q; 810e2ec84fSDave May 820e2ec84fSDave May PetscFunctionBegin; 8320b7b1ceSDave May ierr = PetscDTGaussTensorQuadrature(dim,1,np_1d,-1.0,1.0,&quadrature);CHKERRQ(ierr); 8420b7b1ceSDave May ierr = PetscQuadratureGetData(quadrature,NULL,NULL,&npoints_q,&quadrature_xi,NULL);CHKERRQ(ierr); 850e2ec84fSDave May ierr = PetscMalloc1(dim*npoints_q,&xi);CHKERRQ(ierr); 860e2ec84fSDave May for (q=0; q<npoints_q; q++) { 870e2ec84fSDave May for (d=0; d<dim; d++) { 880e2ec84fSDave May xi[dim*q+d] = quadrature_xi[dim*q+d]; 890e2ec84fSDave May } 900e2ec84fSDave May } 910e2ec84fSDave May ierr = PetscQuadratureDestroy(&quadrature);CHKERRQ(ierr); 920e2ec84fSDave May 930e2ec84fSDave May *_npoints = npoints_q; 940e2ec84fSDave May *_xi = xi; 950e2ec84fSDave May 960e2ec84fSDave May PetscFunctionReturn(0); 970e2ec84fSDave May } 980e2ec84fSDave May 99d3393ee1SDave May PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_DA_Q1(DM dm,DM dmc,PetscInt npoints,DMSwarmPICLayoutType layout) 1000e2ec84fSDave May { 1010e2ec84fSDave May PetscErrorCode ierr; 1020e2ec84fSDave May PetscInt dim,npoints_q; 1030e2ec84fSDave May PetscInt nel,npe,e,q,k,d; 1040e2ec84fSDave May const PetscInt *element_list; 1050e2ec84fSDave May PetscReal **basis; 1060e2ec84fSDave May PetscReal *xi; 1070e2ec84fSDave May Vec coor; 1080e2ec84fSDave May const PetscScalar *_coor; 1090e2ec84fSDave May PetscReal *elcoor; 1100e2ec84fSDave May PetscReal *swarm_coor; 1110e2ec84fSDave May PetscInt *swarm_cellid; 1120e2ec84fSDave May PetscInt pcnt; 1130e2ec84fSDave May 1140e2ec84fSDave May PetscFunctionBegin; 1150e2ec84fSDave May ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 1160e2ec84fSDave May switch (layout) { 1170e2ec84fSDave May case DMSWARMPIC_LAYOUT_REGULAR: 1180e2ec84fSDave May { 1190e2ec84fSDave May PetscInt np_dir[3]; 1200e2ec84fSDave May np_dir[0] = np_dir[1] = np_dir[2] = npoints; 1210e2ec84fSDave May ierr = private_DMSwarmCreateCellLocalCoords_DA_Q1_Regular(dim,np_dir,&npoints_q,&xi);CHKERRQ(ierr); 1220e2ec84fSDave May } 1230e2ec84fSDave May break; 1240e2ec84fSDave May case DMSWARMPIC_LAYOUT_GAUSS: 1250e2ec84fSDave May ierr = private_DMSwarmCreateCellLocalCoords_DA_Q1_Gauss(dim,npoints,&npoints_q,&xi);CHKERRQ(ierr); 1260e2ec84fSDave May break; 127d3393ee1SDave May 128d3393ee1SDave May case DMSWARMPIC_LAYOUT_SUBDIVISION: 129d3393ee1SDave May { 130d3393ee1SDave May PetscInt s,nsub; 131d3393ee1SDave May PetscInt np_dir[3]; 132d3393ee1SDave May nsub = npoints; 133d3393ee1SDave May np_dir[0] = 1; 134d3393ee1SDave May for (s=0; s<nsub; s++) { 135d3393ee1SDave May np_dir[0] *= 2; 136d3393ee1SDave May } 137d3393ee1SDave May np_dir[1] = np_dir[0]; 138d3393ee1SDave May np_dir[2] = np_dir[0]; 139d3393ee1SDave May ierr = private_DMSwarmCreateCellLocalCoords_DA_Q1_Regular(dim,np_dir,&npoints_q,&xi);CHKERRQ(ierr); 140d3393ee1SDave May } 141d3393ee1SDave May break; 1420e2ec84fSDave May default: 1430e2ec84fSDave May SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"A valid DMSwarmPIC layout must be provided"); 1440e2ec84fSDave May break; 1450e2ec84fSDave May } 1460e2ec84fSDave May 1470e2ec84fSDave May ierr = DMDAGetElements(dmc,&nel,&npe,&element_list);CHKERRQ(ierr); 1480e2ec84fSDave May 1490e2ec84fSDave May ierr = PetscMalloc1(dim*npe,&elcoor);CHKERRQ(ierr); 1500e2ec84fSDave May ierr = PetscMalloc1(npoints_q,&basis);CHKERRQ(ierr); 1510e2ec84fSDave May for (q=0; q<npoints_q; q++) { 1520e2ec84fSDave May ierr = PetscMalloc1(npe,&basis[q]);CHKERRQ(ierr); 1530e2ec84fSDave May 1540e2ec84fSDave May switch (dim) { 1550e2ec84fSDave May case 1: 1560e2ec84fSDave May basis[q][0] = 0.5*(1.0 - xi[dim*q+0]); 1570e2ec84fSDave May basis[q][1] = 0.5*(1.0 + xi[dim*q+0]); 1580e2ec84fSDave May break; 1590e2ec84fSDave May case 2: 1600e2ec84fSDave May basis[q][0] = 0.25*(1.0 - xi[dim*q+0])*(1.0 - xi[dim*q+1]); 1610e2ec84fSDave May basis[q][1] = 0.25*(1.0 + xi[dim*q+0])*(1.0 - xi[dim*q+1]); 1620e2ec84fSDave May basis[q][2] = 0.25*(1.0 + xi[dim*q+0])*(1.0 + xi[dim*q+1]); 1630e2ec84fSDave May basis[q][3] = 0.25*(1.0 - xi[dim*q+0])*(1.0 + xi[dim*q+1]); 1640e2ec84fSDave May break; 1650e2ec84fSDave May 1660e2ec84fSDave May case 3: 1670e2ec84fSDave May basis[q][0] = 0.125*(1.0 - xi[dim*q+0])*(1.0 - xi[dim*q+1])*(1.0 - xi[dim*q+2]); 1680e2ec84fSDave May basis[q][1] = 0.125*(1.0 + xi[dim*q+0])*(1.0 - xi[dim*q+1])*(1.0 - xi[dim*q+2]); 1690e2ec84fSDave May basis[q][2] = 0.125*(1.0 + xi[dim*q+0])*(1.0 + xi[dim*q+1])*(1.0 - xi[dim*q+2]); 1700e2ec84fSDave May basis[q][3] = 0.125*(1.0 - xi[dim*q+0])*(1.0 + xi[dim*q+1])*(1.0 - xi[dim*q+2]); 1710e2ec84fSDave May basis[q][4] = 0.125*(1.0 - xi[dim*q+0])*(1.0 - xi[dim*q+1])*(1.0 + xi[dim*q+2]); 1720e2ec84fSDave May basis[q][5] = 0.125*(1.0 + xi[dim*q+0])*(1.0 - xi[dim*q+1])*(1.0 + xi[dim*q+2]); 1730e2ec84fSDave May basis[q][6] = 0.125*(1.0 + xi[dim*q+0])*(1.0 + xi[dim*q+1])*(1.0 + xi[dim*q+2]); 1740e2ec84fSDave May basis[q][7] = 0.125*(1.0 - xi[dim*q+0])*(1.0 + xi[dim*q+1])*(1.0 + xi[dim*q+2]); 1750e2ec84fSDave May break; 1760e2ec84fSDave May } 1770e2ec84fSDave May } 1780e2ec84fSDave May 1790e2ec84fSDave May ierr = DMSwarmSetLocalSizes(dm,npoints_q*nel,-1);CHKERRQ(ierr); 1800e2ec84fSDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 1810e2ec84fSDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 1820e2ec84fSDave May 1830e2ec84fSDave May ierr = DMGetCoordinatesLocal(dmc,&coor);CHKERRQ(ierr); 1840e2ec84fSDave May ierr = VecGetArrayRead(coor,&_coor);CHKERRQ(ierr); 1850e2ec84fSDave May pcnt = 0; 1860e2ec84fSDave May for (e=0; e<nel; e++) { 1870e2ec84fSDave May const PetscInt *element = &element_list[npe*e]; 1880e2ec84fSDave May 1890e2ec84fSDave May for (k=0; k<npe; k++) { 1900e2ec84fSDave May for (d=0; d<dim; d++) { 1910ca99263SDave May elcoor[dim*k+d] = PetscRealPart(_coor[ dim*element[k] + d ]); 1920e2ec84fSDave May } 1930e2ec84fSDave May } 1940e2ec84fSDave May 1950e2ec84fSDave May for (q=0; q<npoints_q; q++) { 1960e2ec84fSDave May for (d=0; d<dim; d++) { 1970e2ec84fSDave May swarm_coor[dim*pcnt+d] = 0.0; 1980e2ec84fSDave May } 1990e2ec84fSDave May for (k=0; k<npe; k++) { 2000e2ec84fSDave May for (d=0; d<dim; d++) { 2010e2ec84fSDave May swarm_coor[dim*pcnt+d] += basis[q][k] * elcoor[dim*k+d]; 2020e2ec84fSDave May } 2030e2ec84fSDave May } 2040e2ec84fSDave May swarm_cellid[pcnt] = e; 2050e2ec84fSDave May pcnt++; 2060e2ec84fSDave May } 2070e2ec84fSDave May } 2080e2ec84fSDave May ierr = VecRestoreArrayRead(coor,&_coor);CHKERRQ(ierr); 2090e2ec84fSDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 2100e2ec84fSDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 2110e2ec84fSDave May ierr = DMDARestoreElements(dmc,&nel,&npe,&element_list);CHKERRQ(ierr); 2120e2ec84fSDave May 2130e2ec84fSDave May ierr = PetscFree(xi);CHKERRQ(ierr); 2140e2ec84fSDave May ierr = PetscFree(elcoor);CHKERRQ(ierr); 2150e2ec84fSDave May for (q=0; q<npoints_q; q++) { 2160e2ec84fSDave May ierr = PetscFree(basis[q]);CHKERRQ(ierr); 2170e2ec84fSDave May } 2180e2ec84fSDave May ierr = PetscFree(basis);CHKERRQ(ierr); 2190e2ec84fSDave May 2200e2ec84fSDave May PetscFunctionReturn(0); 2210e2ec84fSDave May } 2220e2ec84fSDave May 2230e2ec84fSDave May PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_DA(DM dm,DM celldm,DMSwarmPICLayoutType layout,PetscInt layout_param) 2240e2ec84fSDave May { 2250e2ec84fSDave May PetscErrorCode ierr; 2260e2ec84fSDave May DMDAElementType etype; 2270e2ec84fSDave May PetscInt dim; 2280e2ec84fSDave May 2290e2ec84fSDave May PetscFunctionBegin; 2300e2ec84fSDave May ierr = DMDAGetElementType(celldm,&etype);CHKERRQ(ierr); 2310e2ec84fSDave May ierr = DMGetDimension(celldm,&dim);CHKERRQ(ierr); 2320e2ec84fSDave May switch (etype) { 2330e2ec84fSDave May case DMDA_ELEMENT_P1: 2340e2ec84fSDave May SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"DA support is not currently available for DMDA_ELEMENT_P1"); 2350e2ec84fSDave May break; 2360e2ec84fSDave May case DMDA_ELEMENT_Q1: 2370e2ec84fSDave May if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Support only available for dim = 2, 3"); 2380e2ec84fSDave May ierr = private_DMSwarmInsertPointsUsingCellDM_DA_Q1(dm,celldm,layout_param,layout);CHKERRQ(ierr); 2390e2ec84fSDave May break; 2400e2ec84fSDave May } 2410e2ec84fSDave May PetscFunctionReturn(0); 2420e2ec84fSDave May } 2431f52046fSDave May 2441f52046fSDave May PetscErrorCode DMSwarmProjectField_ApproxQ1_DA_2D(DM swarm,PetscReal *swarm_field,DM dm,Vec v_field) 2451f52046fSDave May { 2461f52046fSDave May PetscErrorCode ierr; 2471f52046fSDave May Vec v_field_l,denom_l,coor_l,denom; 248a5f152d1SDave May PetscScalar *_field_l,*_denom_l; 2491f52046fSDave May PetscInt k,p,e,npoints,nel,npe; 2501f52046fSDave May PetscInt *mpfield_cell; 2511f52046fSDave May PetscReal *mpfield_coor; 2521f52046fSDave May const PetscInt *element_list; 2531f52046fSDave May const PetscInt *element; 254a5f152d1SDave May PetscScalar xi_p[2],Ni[4]; 255a5f152d1SDave May const PetscScalar *_coor; 2561f52046fSDave May 2571f52046fSDave May PetscFunctionBegin; 2581f52046fSDave May ierr = VecZeroEntries(v_field);CHKERRQ(ierr); 2591f52046fSDave May 2601f52046fSDave May ierr = DMGetLocalVector(dm,&v_field_l);CHKERRQ(ierr); 2611f52046fSDave May ierr = DMGetGlobalVector(dm,&denom);CHKERRQ(ierr); 2621f52046fSDave May ierr = DMGetLocalVector(dm,&denom_l);CHKERRQ(ierr); 2631f52046fSDave May ierr = VecZeroEntries(v_field_l);CHKERRQ(ierr); 2641f52046fSDave May ierr = VecZeroEntries(denom);CHKERRQ(ierr); 2651f52046fSDave May ierr = VecZeroEntries(denom_l);CHKERRQ(ierr); 2661f52046fSDave May 2671f52046fSDave May ierr = VecGetArray(v_field_l,&_field_l);CHKERRQ(ierr); 2681f52046fSDave May ierr = VecGetArray(denom_l,&_denom_l);CHKERRQ(ierr); 2691f52046fSDave May 2701f52046fSDave May ierr = DMGetCoordinatesLocal(dm,&coor_l);CHKERRQ(ierr); 2711f52046fSDave May ierr = VecGetArrayRead(coor_l,&_coor);CHKERRQ(ierr); 2721f52046fSDave May 2731f52046fSDave May ierr = DMDAGetElements(dm,&nel,&npe,&element_list);CHKERRQ(ierr); 2741f52046fSDave May ierr = DMSwarmGetLocalSize(swarm,&npoints);CHKERRQ(ierr); 2751f52046fSDave May ierr = DMSwarmGetField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&mpfield_coor);CHKERRQ(ierr); 2761f52046fSDave May ierr = DMSwarmGetField(swarm,DMSwarmPICField_cellid,NULL,NULL,(void**)&mpfield_cell);CHKERRQ(ierr); 2771f52046fSDave May 2781f52046fSDave May for (p=0; p<npoints; p++) { 2791f52046fSDave May PetscReal *coor_p; 280a5f152d1SDave May const PetscScalar *x0; 281a5f152d1SDave May const PetscScalar *x2; 282a5f152d1SDave May PetscScalar dx[2]; 2831f52046fSDave May 2841f52046fSDave May e = mpfield_cell[p]; 2851f52046fSDave May coor_p = &mpfield_coor[2*p]; 2861f52046fSDave May element = &element_list[npe*e]; 2871f52046fSDave May 2881f52046fSDave May /* compute local coordinates: (xp-x0)/dx = (xip+1)/2 */ 2891f52046fSDave May x0 = &_coor[2*element[0]]; 2901f52046fSDave May x2 = &_coor[2*element[2]]; 2911f52046fSDave May 2921f52046fSDave May dx[0] = x2[0] - x0[0]; 2931f52046fSDave May dx[1] = x2[1] - x0[1]; 2941f52046fSDave May 2951f52046fSDave May xi_p[0] = 2.0 * (coor_p[0] - x0[0])/dx[0] - 1.0; 2961f52046fSDave May xi_p[1] = 2.0 * (coor_p[1] - x0[1])/dx[1] - 1.0; 2971f52046fSDave May 2981f52046fSDave May /* evaluate basis functions */ 2991f52046fSDave May Ni[0] = 0.25*(1.0 - xi_p[0])*(1.0 - xi_p[1]); 3001f52046fSDave May Ni[1] = 0.25*(1.0 + xi_p[0])*(1.0 - xi_p[1]); 3011f52046fSDave May Ni[2] = 0.25*(1.0 + xi_p[0])*(1.0 + xi_p[1]); 3021f52046fSDave May Ni[3] = 0.25*(1.0 - xi_p[0])*(1.0 + xi_p[1]); 3031f52046fSDave May 3041f52046fSDave May for (k=0; k<npe; k++) { 3051f52046fSDave May _field_l[ element[k] ] += Ni[k] * swarm_field[p]; 3061f52046fSDave May _denom_l[ element[k] ] += Ni[k]; 3071f52046fSDave May } 3081f52046fSDave May } 3091f52046fSDave May 3101f52046fSDave May ierr = DMSwarmRestoreField(swarm,DMSwarmPICField_cellid,NULL,NULL,(void**)&mpfield_cell);CHKERRQ(ierr); 3111f52046fSDave May ierr = DMSwarmRestoreField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&mpfield_coor);CHKERRQ(ierr); 3121f52046fSDave May ierr = DMDARestoreElements(dm,&nel,&npe,&element_list);CHKERRQ(ierr); 3131f52046fSDave May ierr = VecRestoreArrayRead(coor_l,&_coor);CHKERRQ(ierr); 3141f52046fSDave May ierr = VecRestoreArray(v_field_l,&_field_l);CHKERRQ(ierr); 3151f52046fSDave May ierr = VecRestoreArray(denom_l,&_denom_l);CHKERRQ(ierr); 3161f52046fSDave May 3171f52046fSDave May ierr = DMLocalToGlobalBegin(dm,v_field_l,ADD_VALUES,v_field);CHKERRQ(ierr); 3181f52046fSDave May ierr = DMLocalToGlobalEnd(dm,v_field_l,ADD_VALUES,v_field);CHKERRQ(ierr); 3191f52046fSDave May ierr = DMLocalToGlobalBegin(dm,denom_l,ADD_VALUES,denom);CHKERRQ(ierr); 3201f52046fSDave May ierr = DMLocalToGlobalEnd(dm,denom_l,ADD_VALUES,denom);CHKERRQ(ierr); 3211f52046fSDave May 3221f52046fSDave May ierr = VecPointwiseDivide(v_field,v_field,denom);CHKERRQ(ierr); 3231f52046fSDave May 3241f52046fSDave May ierr = DMRestoreLocalVector(dm,&v_field_l);CHKERRQ(ierr); 3251f52046fSDave May ierr = DMRestoreLocalVector(dm,&denom_l);CHKERRQ(ierr); 3261f52046fSDave May ierr = DMRestoreGlobalVector(dm,&denom);CHKERRQ(ierr); 3271f52046fSDave May 3281f52046fSDave May PetscFunctionReturn(0); 3291f52046fSDave May } 3301f52046fSDave May 331*77048351SPatrick Sanan PetscErrorCode private_DMSwarmProjectFields_DA(DM swarm,DM celldm,PetscInt project_type,PetscInt nfields,DMSwarmDataField dfield[],Vec vecs[]) 3321f52046fSDave May { 3331f52046fSDave May PetscErrorCode ierr; 3341f52046fSDave May PetscInt f,dim; 3351f52046fSDave May DMDAElementType etype; 3361f52046fSDave May 3371f52046fSDave May PetscFunctionBegin; 3381f52046fSDave May ierr = DMDAGetElementType(celldm,&etype);CHKERRQ(ierr); 3391f52046fSDave May if (etype == DMDA_ELEMENT_P1) SETERRQ(PetscObjectComm((PetscObject)swarm),PETSC_ERR_SUP,"Only Q1 DMDA supported"); 3401f52046fSDave May 3411f52046fSDave May ierr = DMGetDimension(swarm,&dim);CHKERRQ(ierr); 3421f52046fSDave May switch (dim) { 3431f52046fSDave May case 2: 3441f52046fSDave May for (f=0; f<nfields; f++) { 3451f52046fSDave May PetscReal *swarm_field; 3461f52046fSDave May 347*77048351SPatrick Sanan ierr = DMSwarmDataFieldGetEntries(dfield[f],(void**)&swarm_field);CHKERRQ(ierr); 3481f52046fSDave May ierr = DMSwarmProjectField_ApproxQ1_DA_2D(swarm,swarm_field,celldm,vecs[f]);CHKERRQ(ierr); 3491f52046fSDave May } 3501f52046fSDave May break; 3511f52046fSDave May case 3: 3521f52046fSDave May SETERRQ(PetscObjectComm((PetscObject)swarm),PETSC_ERR_SUP,"No support for 3D"); 3531f52046fSDave May break; 3541f52046fSDave May default: 3551f52046fSDave May break; 3561f52046fSDave May } 3571f52046fSDave May PetscFunctionReturn(0); 3581f52046fSDave May } 359