10e2ec84fSDave May #define PETSCDM_DLL 20e2ec84fSDave May #include <petsc/private/dmswarmimpl.h> /*I "petscdmswarm.h" I*/ 30e2ec84fSDave May #include <petscsf.h> 4b799feefSDave May #include <petscdmda.h> 5b799feefSDave May #include <petscdmplex.h> 635b38c60SMatthew G. Knepley #include <petscdt.h> 7279f676cSBarry Smith #include "../src/dm/impls/swarm/data_bucket.h" 80e2ec84fSDave May 935b38c60SMatthew G. Knepley #include <petsc/private/petscfeimpl.h> /* For CoordinatesRefToReal() */ 1035b38c60SMatthew G. Knepley 110e2ec84fSDave May /* 125627991aSBarry Smith Error checking to ensure the swarm type is correct and that a cell DM has been set 130e2ec84fSDave May */ 140e2ec84fSDave May #define DMSWARMPICVALID(dm) \ 150e2ec84fSDave May { \ 160e2ec84fSDave May DM_Swarm *_swarm = (DM_Swarm*)(dm)->data; \ 175f80ce2aSJacob Faibussowitsch PetscCheck(_swarm->swarm_type == DMSWARM_PIC,PetscObjectComm((PetscObject)(dm)),PETSC_ERR_SUP,"Valid only for DMSwarm-PIC. You must call DMSwarmSetType(dm,DMSWARM_PIC)"); \ 180e2ec84fSDave May else \ 19*28b400f6SJacob Faibussowitsch PetscCheck(_swarm->dmcell,PetscObjectComm((PetscObject)(dm)),PETSC_ERR_SUP,"Valid only for DMSwarmPIC if the cell DM is set. You must call DMSwarmSetCellDM(dm,celldm)"); \ 200e2ec84fSDave May } 210e2ec84fSDave May 220e2ec84fSDave May /* Coordinate insertition/addition API */ 230e2ec84fSDave May /*@C 240e2ec84fSDave May DMSwarmSetPointsUniformCoordinates - Set point coordinates in a DMSwarm on a regular (ijk) grid 250e2ec84fSDave May 26d083f849SBarry Smith Collective on dm 270e2ec84fSDave May 280e2ec84fSDave May Input parameters: 290e2ec84fSDave May + dm - the DMSwarm 300e2ec84fSDave May . min - minimum coordinate values in the x, y, z directions (array of length dim) 310e2ec84fSDave May . max - maximum coordinate values in the x, y, z directions (array of length dim) 320e2ec84fSDave May . npoints - number of points in each spatial direction (array of length dim) 330e2ec84fSDave May - mode - indicates whether to append points to the swarm (ADD_VALUES), or over-ride existing points (INSERT_VALUES) 340e2ec84fSDave May 350e2ec84fSDave May Level: beginner 360e2ec84fSDave May 370e2ec84fSDave May Notes: 380e2ec84fSDave May When using mode = INSERT_VALUES, this method will reset the number of particles in the DMSwarm 390e2ec84fSDave May to be npoints[0]*npoints[1] (2D) or npoints[0]*npoints[1]*npoints[2] (3D). When using mode = ADD_VALUES, 400e2ec84fSDave May new points will be appended to any already existing in the DMSwarm 410e2ec84fSDave May 420e2ec84fSDave May .seealso: DMSwarmSetType(), DMSwarmSetCellDM(), DMSwarmType 430e2ec84fSDave May @*/ 440e2ec84fSDave May PETSC_EXTERN PetscErrorCode DMSwarmSetPointsUniformCoordinates(DM dm,PetscReal min[],PetscReal max[],PetscInt npoints[],InsertMode mode) 450e2ec84fSDave May { 460e2ec84fSDave May PetscReal gmin[] = {PETSC_MAX_REAL ,PETSC_MAX_REAL, PETSC_MAX_REAL}; 470e2ec84fSDave May PetscReal gmax[] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL}; 480e2ec84fSDave May PetscInt i,j,k,N,bs,b,n_estimate,n_curr,n_new_est,p,n_found; 490e2ec84fSDave May Vec coorlocal; 500e2ec84fSDave May const PetscScalar *_coor; 510e2ec84fSDave May DM celldm; 520e2ec84fSDave May PetscReal dx[3]; 532e3d3761SDave May PetscInt _npoints[] = { 0, 0, 1 }; 540e2ec84fSDave May Vec pos; 550e2ec84fSDave May PetscScalar *_pos; 560e2ec84fSDave May PetscReal *swarm_coor; 570e2ec84fSDave May PetscInt *swarm_cellid; 580e2ec84fSDave May PetscSF sfcell = NULL; 590e2ec84fSDave May const PetscSFNode *LA_sfcell; 600e2ec84fSDave May 610e2ec84fSDave May PetscFunctionBegin; 620e2ec84fSDave May DMSWARMPICVALID(dm); 635f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetCellDM(dm,&celldm)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinatesLocal(celldm,&coorlocal)); 655f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(coorlocal,&N)); 665f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetBlockSize(coorlocal,&bs)); 670e2ec84fSDave May N = N / bs; 685f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(coorlocal,&_coor)); 690e2ec84fSDave May for (i=0; i<N; i++) { 700e2ec84fSDave May for (b=0; b<bs; b++) { 71a5f152d1SDave May gmin[b] = PetscMin(gmin[b],PetscRealPart(_coor[bs*i+b])); 72a5f152d1SDave May gmax[b] = PetscMax(gmax[b],PetscRealPart(_coor[bs*i+b])); 730e2ec84fSDave May } 740e2ec84fSDave May } 755f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(coorlocal,&_coor)); 760e2ec84fSDave May 770e2ec84fSDave May for (b=0; b<bs; b++) { 7871844bb8SDave May if (npoints[b] > 1) { 790e2ec84fSDave May dx[b] = (max[b] - min[b])/((PetscReal)(npoints[b]-1)); 8071844bb8SDave May } else { 8171844bb8SDave May dx[b] = 0.0; 8271844bb8SDave May } 832e3d3761SDave May _npoints[b] = npoints[b]; 840e2ec84fSDave May } 850e2ec84fSDave May 860e2ec84fSDave May /* determine number of points living in the bounding box */ 870e2ec84fSDave May n_estimate = 0; 882e3d3761SDave May for (k=0; k<_npoints[2]; k++) { 892e3d3761SDave May for (j=0; j<_npoints[1]; j++) { 902e3d3761SDave May for (i=0; i<_npoints[0]; i++) { 910e2ec84fSDave May PetscReal xp[] = {0.0,0.0,0.0}; 920e2ec84fSDave May PetscInt ijk[3]; 930e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 940e2ec84fSDave May 950e2ec84fSDave May ijk[0] = i; 960e2ec84fSDave May ijk[1] = j; 970e2ec84fSDave May ijk[2] = k; 980e2ec84fSDave May for (b=0; b<bs; b++) { 990e2ec84fSDave May xp[b] = min[b] + ijk[b] * dx[b]; 1000e2ec84fSDave May } 1010e2ec84fSDave May for (b=0; b<bs; b++) { 1020e2ec84fSDave May if (xp[b] < gmin[b]) { point_inside = PETSC_FALSE; } 1030e2ec84fSDave May if (xp[b] > gmax[b]) { point_inside = PETSC_FALSE; } 1040e2ec84fSDave May } 1050e2ec84fSDave May if (point_inside) { n_estimate++; } 1060e2ec84fSDave May } 1070e2ec84fSDave May } 1080e2ec84fSDave May } 1090e2ec84fSDave May 1100e2ec84fSDave May /* create candidate list */ 1115f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_SELF,&pos)); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(pos,bs*n_estimate,PETSC_DECIDE)); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetBlockSize(pos,bs)); 1145f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(pos)); 1155f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(pos,&_pos)); 1160e2ec84fSDave May 1170e2ec84fSDave May n_estimate = 0; 1182e3d3761SDave May for (k=0; k<_npoints[2]; k++) { 1192e3d3761SDave May for (j=0; j<_npoints[1]; j++) { 1202e3d3761SDave May for (i=0; i<_npoints[0]; i++) { 1210e2ec84fSDave May PetscReal xp[] = {0.0,0.0,0.0}; 1220e2ec84fSDave May PetscInt ijk[3]; 1230e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 1240e2ec84fSDave May 1250e2ec84fSDave May ijk[0] = i; 1260e2ec84fSDave May ijk[1] = j; 1270e2ec84fSDave May ijk[2] = k; 1280e2ec84fSDave May for (b=0; b<bs; b++) { 1290e2ec84fSDave May xp[b] = min[b] + ijk[b] * dx[b]; 1300e2ec84fSDave May } 1310e2ec84fSDave May for (b=0; b<bs; b++) { 1320e2ec84fSDave May if (xp[b] < gmin[b]) { point_inside = PETSC_FALSE; } 1330e2ec84fSDave May if (xp[b] > gmax[b]) { point_inside = PETSC_FALSE; } 1340e2ec84fSDave May } 1350e2ec84fSDave May if (point_inside) { 1360e2ec84fSDave May for (b=0; b<bs; b++) { 1370e2ec84fSDave May _pos[bs*n_estimate+b] = xp[b]; 1380e2ec84fSDave May } 1390e2ec84fSDave May n_estimate++; 1400e2ec84fSDave May } 1410e2ec84fSDave May } 1420e2ec84fSDave May } 1430e2ec84fSDave May } 1445f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(pos,&_pos)); 1450e2ec84fSDave May 1460e2ec84fSDave May /* locate points */ 1475f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocatePoints(celldm,pos,DM_POINTLOCATION_NONE,&sfcell)); 1485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell)); 1490e2ec84fSDave May n_found = 0; 1500e2ec84fSDave May for (p=0; p<n_estimate; p++) { 1510e2ec84fSDave May if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { 1520e2ec84fSDave May n_found++; 1530e2ec84fSDave May } 1540e2ec84fSDave May } 1550e2ec84fSDave May 1560e2ec84fSDave May /* adjust size */ 1570e2ec84fSDave May if (mode == ADD_VALUES) { 1585f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetLocalSize(dm,&n_curr)); 1590e2ec84fSDave May n_new_est = n_curr + n_found; 1605f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSetLocalSizes(dm,n_new_est,-1)); 1610e2ec84fSDave May } 1620e2ec84fSDave May if (mode == INSERT_VALUES) { 1630e2ec84fSDave May n_curr = 0; 1640e2ec84fSDave May n_new_est = n_found; 1655f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSetLocalSizes(dm,n_new_est,-1)); 1660e2ec84fSDave May } 1670e2ec84fSDave May 1680e2ec84fSDave May /* initialize new coords, cell owners, pid */ 1695f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(pos,&_coor)); 1705f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor)); 1715f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid)); 1720e2ec84fSDave May n_found = 0; 1730e2ec84fSDave May for (p=0; p<n_estimate; p++) { 1740e2ec84fSDave May if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { 1750e2ec84fSDave May for (b=0; b<bs; b++) { 1760ca99263SDave May swarm_coor[bs*(n_curr + n_found) + b] = PetscRealPart(_coor[bs*p+b]); 1770e2ec84fSDave May } 1780e2ec84fSDave May swarm_cellid[n_curr + n_found] = LA_sfcell[p].index; 1790e2ec84fSDave May n_found++; 1800e2ec84fSDave May } 1810e2ec84fSDave May } 1825f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid)); 1835f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor)); 1845f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(pos,&_coor)); 1850e2ec84fSDave May 1865f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFDestroy(&sfcell)); 1875f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&pos)); 1880e2ec84fSDave May PetscFunctionReturn(0); 1890e2ec84fSDave May } 1900e2ec84fSDave May 1910e2ec84fSDave May /*@C 1920e2ec84fSDave May DMSwarmSetPointCoordinates - Set point coordinates in a DMSwarm from a user defined list 1930e2ec84fSDave May 194d083f849SBarry Smith Collective on dm 1950e2ec84fSDave May 1960e2ec84fSDave May Input parameters: 1970e2ec84fSDave May + dm - the DMSwarm 1980e2ec84fSDave May . npoints - the number of points to insert 1990e2ec84fSDave May . coor - the coordinate values 2000e2ec84fSDave May . redundant - if set to PETSC_TRUE, it is assumed that npoints and coor[] are only valid on rank 0 and should be broadcast to other ranks 2010e2ec84fSDave May - mode - indicates whether to append points to the swarm (ADD_VALUES), or over-ride existing points (INSERT_VALUES) 2020e2ec84fSDave May 2030e2ec84fSDave May Level: beginner 2040e2ec84fSDave May 2050e2ec84fSDave May Notes: 2060e2ec84fSDave May If the user has specified redundant = PETSC_FALSE, the cell DM will attempt to locate the coordinates provided by coor[] within 2070e2ec84fSDave May its sub-domain. If they any values within coor[] are not located in the sub-domain, they will be ignored and will not get 2080e2ec84fSDave May added to the DMSwarm. 2090e2ec84fSDave May 2100e2ec84fSDave May .seealso: DMSwarmSetType(), DMSwarmSetCellDM(), DMSwarmType, DMSwarmSetPointsUniformCoordinates() 2110e2ec84fSDave May @*/ 2120e2ec84fSDave May PETSC_EXTERN PetscErrorCode DMSwarmSetPointCoordinates(DM dm,PetscInt npoints,PetscReal coor[],PetscBool redundant,InsertMode mode) 2130e2ec84fSDave May { 2140e2ec84fSDave May PetscReal gmin[] = {PETSC_MAX_REAL ,PETSC_MAX_REAL, PETSC_MAX_REAL}; 2150e2ec84fSDave May PetscReal gmax[] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL}; 2160e2ec84fSDave May PetscInt i,N,bs,b,n_estimate,n_curr,n_new_est,p,n_found; 2170e2ec84fSDave May Vec coorlocal; 2180e2ec84fSDave May const PetscScalar *_coor; 2190e2ec84fSDave May DM celldm; 2200e2ec84fSDave May Vec pos; 2210e2ec84fSDave May PetscScalar *_pos; 2220e2ec84fSDave May PetscReal *swarm_coor; 2230e2ec84fSDave May PetscInt *swarm_cellid; 2240e2ec84fSDave May PetscSF sfcell = NULL; 2250e2ec84fSDave May const PetscSFNode *LA_sfcell; 2260e2ec84fSDave May PetscReal *my_coor; 2270e2ec84fSDave May PetscInt my_npoints; 2280e2ec84fSDave May PetscMPIInt rank; 2290e2ec84fSDave May MPI_Comm comm; 2300e2ec84fSDave May 2310e2ec84fSDave May PetscFunctionBegin; 2320e2ec84fSDave May DMSWARMPICVALID(dm); 2335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject)dm,&comm)); 2345f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm,&rank)); 2350e2ec84fSDave May 2365f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetCellDM(dm,&celldm)); 2375f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinatesLocal(celldm,&coorlocal)); 2385f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(coorlocal,&N)); 2395f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetBlockSize(coorlocal,&bs)); 2400e2ec84fSDave May N = N / bs; 2415f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(coorlocal,&_coor)); 2420e2ec84fSDave May for (i=0; i<N; i++) { 2430e2ec84fSDave May for (b=0; b<bs; b++) { 244a5f152d1SDave May gmin[b] = PetscMin(gmin[b],PetscRealPart(_coor[bs*i+b])); 245a5f152d1SDave May gmax[b] = PetscMax(gmax[b],PetscRealPart(_coor[bs*i+b])); 2460e2ec84fSDave May } 2470e2ec84fSDave May } 2485f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(coorlocal,&_coor)); 2490e2ec84fSDave May 2500e2ec84fSDave May /* broadcast points from rank 0 if requested */ 2510e2ec84fSDave May if (redundant) { 2520e2ec84fSDave May my_npoints = npoints; 2535f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Bcast(&my_npoints,1,MPIU_INT,0,comm)); 2540e2ec84fSDave May 2550e2ec84fSDave May if (rank > 0) { /* allocate space */ 2565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(bs*my_npoints,&my_coor)); 2570e2ec84fSDave May } else { 2580e2ec84fSDave May my_coor = coor; 2590e2ec84fSDave May } 2605f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Bcast(my_coor,bs*my_npoints,MPIU_REAL,0,comm)); 2610e2ec84fSDave May } else { 2620e2ec84fSDave May my_npoints = npoints; 2630e2ec84fSDave May my_coor = coor; 2640e2ec84fSDave May } 2650e2ec84fSDave May 2660e2ec84fSDave May /* determine the number of points living in the bounding box */ 2670e2ec84fSDave May n_estimate = 0; 2680e2ec84fSDave May for (i=0; i<my_npoints; i++) { 2690e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 2700e2ec84fSDave May 2710e2ec84fSDave May for (b=0; b<bs; b++) { 2720e2ec84fSDave May if (my_coor[bs*i+b] < gmin[b]) { point_inside = PETSC_FALSE; } 2730e2ec84fSDave May if (my_coor[bs*i+b] > gmax[b]) { point_inside = PETSC_FALSE; } 2740e2ec84fSDave May } 2750e2ec84fSDave May if (point_inside) { n_estimate++; } 2760e2ec84fSDave May } 2770e2ec84fSDave May 2780e2ec84fSDave May /* create candidate list */ 2795f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_SELF,&pos)); 2805f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(pos,bs*n_estimate,PETSC_DECIDE)); 2815f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetBlockSize(pos,bs)); 2825f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(pos)); 2835f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(pos,&_pos)); 2840e2ec84fSDave May 2850e2ec84fSDave May n_estimate = 0; 2860e2ec84fSDave May for (i=0; i<my_npoints; i++) { 2870e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 2880e2ec84fSDave May 2890e2ec84fSDave May for (b=0; b<bs; b++) { 2900e2ec84fSDave May if (my_coor[bs*i+b] < gmin[b]) { point_inside = PETSC_FALSE; } 2910e2ec84fSDave May if (my_coor[bs*i+b] > gmax[b]) { point_inside = PETSC_FALSE; } 2920e2ec84fSDave May } 2930e2ec84fSDave May if (point_inside) { 2940e2ec84fSDave May for (b=0; b<bs; b++) { 2950e2ec84fSDave May _pos[bs*n_estimate+b] = my_coor[bs*i+b]; 2960e2ec84fSDave May } 2970e2ec84fSDave May n_estimate++; 2980e2ec84fSDave May } 2990e2ec84fSDave May } 3005f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(pos,&_pos)); 3010e2ec84fSDave May 3020e2ec84fSDave May /* locate points */ 3035f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocatePoints(celldm,pos,DM_POINTLOCATION_NONE,&sfcell)); 3040e2ec84fSDave May 3055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell)); 3060e2ec84fSDave May n_found = 0; 3070e2ec84fSDave May for (p=0; p<n_estimate; p++) { 3080e2ec84fSDave May if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { 3090e2ec84fSDave May n_found++; 3100e2ec84fSDave May } 3110e2ec84fSDave May } 3120e2ec84fSDave May 3130e2ec84fSDave May /* adjust size */ 3140e2ec84fSDave May if (mode == ADD_VALUES) { 3155f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetLocalSize(dm,&n_curr)); 3160e2ec84fSDave May n_new_est = n_curr + n_found; 3175f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSetLocalSizes(dm,n_new_est,-1)); 3180e2ec84fSDave May } 3190e2ec84fSDave May if (mode == INSERT_VALUES) { 3200e2ec84fSDave May n_curr = 0; 3210e2ec84fSDave May n_new_est = n_found; 3225f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSetLocalSizes(dm,n_new_est,-1)); 3230e2ec84fSDave May } 3240e2ec84fSDave May 3250e2ec84fSDave May /* initialize new coords, cell owners, pid */ 3265f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(pos,&_coor)); 3275f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor)); 3285f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid)); 3290e2ec84fSDave May n_found = 0; 3300e2ec84fSDave May for (p=0; p<n_estimate; p++) { 3310e2ec84fSDave May if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { 3320e2ec84fSDave May for (b=0; b<bs; b++) { 3330ca99263SDave May swarm_coor[bs*(n_curr + n_found) + b] = PetscRealPart(_coor[bs*p+b]); 3340e2ec84fSDave May } 3350e2ec84fSDave May swarm_cellid[n_curr + n_found] = LA_sfcell[p].index; 3360e2ec84fSDave May n_found++; 3370e2ec84fSDave May } 3380e2ec84fSDave May } 3395f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid)); 3405f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor)); 3415f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(pos,&_coor)); 3420e2ec84fSDave May 3430e2ec84fSDave May if (redundant) { 3440e2ec84fSDave May if (rank > 0) { 3455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(my_coor)); 3460e2ec84fSDave May } 3470e2ec84fSDave May } 3485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFDestroy(&sfcell)); 3495f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&pos)); 3500e2ec84fSDave May PetscFunctionReturn(0); 3510e2ec84fSDave May } 3520e2ec84fSDave May 3530e2ec84fSDave May extern PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_DA(DM,DM,DMSwarmPICLayoutType,PetscInt); 3540e2ec84fSDave May extern PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX(DM,DM,DMSwarmPICLayoutType,PetscInt); 3550e2ec84fSDave May 3560e2ec84fSDave May /*@C 3570e2ec84fSDave May DMSwarmInsertPointsUsingCellDM - Insert point coordinates within each cell 3580e2ec84fSDave May 3590e2ec84fSDave May Not collective 3600e2ec84fSDave May 3610e2ec84fSDave May Input parameters: 3620e2ec84fSDave May + dm - the DMSwarm 3630e2ec84fSDave May . layout_type - method used to fill each cell with the cell DM 3640e2ec84fSDave May - fill_param - parameter controlling how many points per cell are added (the meaning of this parameter is dependent on the layout type) 3650e2ec84fSDave May 3660e2ec84fSDave May Level: beginner 3670e2ec84fSDave May 3680e2ec84fSDave May Notes: 369e7af74c8SDave May 370e7af74c8SDave May The insert method will reset any previous defined points within the DMSwarm. 371e7af74c8SDave May 372e7af74c8SDave May When using a DMDA both 2D and 3D are supported for all layout types provided you are using DMDA_ELEMENT_Q1. 373e7af74c8SDave May 374e7af74c8SDave May When using a DMPLEX the following case are supported: 375ea3d7275SDave May (i) DMSWARMPIC_LAYOUT_REGULAR: 2D (triangle), 376ea3d7275SDave May (ii) DMSWARMPIC_LAYOUT_GAUSS: 2D and 3D provided the cell is a tri/tet or a quad/hex, 377ea3d7275SDave May (iii) DMSWARMPIC_LAYOUT_SUBDIVISION: 2D and 3D for quad/hex and 2D tri. 3780e2ec84fSDave May 3790e2ec84fSDave May .seealso: DMSwarmPICLayoutType, DMSwarmSetType(), DMSwarmSetCellDM(), DMSwarmType 3800e2ec84fSDave May @*/ 3810e2ec84fSDave May PETSC_EXTERN PetscErrorCode DMSwarmInsertPointsUsingCellDM(DM dm,DMSwarmPICLayoutType layout_type,PetscInt fill_param) 3820e2ec84fSDave May { 3830e2ec84fSDave May DM celldm; 3840e2ec84fSDave May PetscBool isDA,isPLEX; 3850e2ec84fSDave May 3860e2ec84fSDave May PetscFunctionBegin; 3870e2ec84fSDave May DMSWARMPICVALID(dm); 3885f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetCellDM(dm,&celldm)); 3895f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)celldm,DMDA,&isDA)); 3905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)celldm,DMPLEX,&isPLEX)); 3910e2ec84fSDave May if (isDA) { 3925f80ce2aSJacob Faibussowitsch CHKERRQ(private_DMSwarmInsertPointsUsingCellDM_DA(dm,celldm,layout_type,fill_param)); 3930e2ec84fSDave May } else if (isPLEX) { 3945f80ce2aSJacob Faibussowitsch CHKERRQ(private_DMSwarmInsertPointsUsingCellDM_PLEX(dm,celldm,layout_type,fill_param)); 3950e2ec84fSDave May } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only supported for cell DMs of type DMDA and DMPLEX"); 3960e2ec84fSDave May PetscFunctionReturn(0); 3970e2ec84fSDave May } 3980e2ec84fSDave May 399431382f2SDave May extern PetscErrorCode private_DMSwarmSetPointCoordinatesCellwise_PLEX(DM,DM,PetscInt,PetscReal*); 400431382f2SDave May 401431382f2SDave May /*@C 402431382f2SDave May DMSwarmSetPointCoordinatesCellwise - Insert point coordinates (defined over the reference cell) within each cell 403431382f2SDave May 404431382f2SDave May Not collective 405431382f2SDave May 406431382f2SDave May Input parameters: 407431382f2SDave May + dm - the DMSwarm 408431382f2SDave May . celldm - the cell DM 409431382f2SDave May . npoints - the number of points to insert in each cell 410431382f2SDave May - xi - the coordinates (defined in the local coordinate system for each cell) to insert 411431382f2SDave May 412431382f2SDave May Level: beginner 413431382f2SDave May 414431382f2SDave May Notes: 415431382f2SDave May The method will reset any previous defined points within the DMSwarm. 416ea3d7275SDave May Only supported for DMPLEX. If you are using a DMDA it is recommended to either use 417e7af74c8SDave May DMSwarmInsertPointsUsingCellDM(), or extract and set the coordinates yourself the following code 418e7af74c8SDave May 419e7af74c8SDave May $ PetscReal *coor; 420e7af74c8SDave May $ DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&coor); 421e7af74c8SDave May $ // user code to define the coordinates here 422e7af74c8SDave May $ DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&coor); 423431382f2SDave May 424431382f2SDave May .seealso: DMSwarmSetCellDM(), DMSwarmInsertPointsUsingCellDM() 425431382f2SDave May @*/ 426431382f2SDave May PETSC_EXTERN PetscErrorCode DMSwarmSetPointCoordinatesCellwise(DM dm,PetscInt npoints,PetscReal xi[]) 4270e2ec84fSDave May { 428431382f2SDave May DM celldm; 429431382f2SDave May PetscBool isDA,isPLEX; 430431382f2SDave May 4310e2ec84fSDave May PetscFunctionBegin; 432431382f2SDave May DMSWARMPICVALID(dm); 4335f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetCellDM(dm,&celldm)); 4345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)celldm,DMDA,&isDA)); 4355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)celldm,DMPLEX,&isPLEX)); 436*28b400f6SJacob Faibussowitsch PetscCheck(!isDA,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only supported for cell DMs of type DMPLEX. Recommended you use DMSwarmInsertPointsUsingCellDM()"); 437431382f2SDave May else if (isPLEX) { 4385f80ce2aSJacob Faibussowitsch CHKERRQ(private_DMSwarmSetPointCoordinatesCellwise_PLEX(dm,celldm,npoints,xi)); 439431382f2SDave May } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only supported for cell DMs of type DMDA and DMPLEX"); 4400e2ec84fSDave May PetscFunctionReturn(0); 4410e2ec84fSDave May } 442431382f2SDave May 4430e2ec84fSDave May /* Field projection API */ 44477048351SPatrick Sanan extern PetscErrorCode private_DMSwarmProjectFields_DA(DM swarm,DM celldm,PetscInt project_type,PetscInt nfields,DMSwarmDataField dfield[],Vec vecs[]); 44577048351SPatrick Sanan extern PetscErrorCode private_DMSwarmProjectFields_PLEX(DM swarm,DM celldm,PetscInt project_type,PetscInt nfields,DMSwarmDataField dfield[],Vec vecs[]); 44694f7d2dcSDave May 44794f7d2dcSDave May /*@C 44894f7d2dcSDave May DMSwarmProjectFields - Project a set of swarm fields onto the cell DM 44994f7d2dcSDave May 450d083f849SBarry Smith Collective on dm 45194f7d2dcSDave May 45294f7d2dcSDave May Input parameters: 45394f7d2dcSDave May + dm - the DMSwarm 45494f7d2dcSDave May . nfields - the number of swarm fields to project 45594f7d2dcSDave May . fieldnames - the textual names of the swarm fields to project 45694f7d2dcSDave May . fields - an array of Vec's of length nfields 45794f7d2dcSDave May - reuse - flag indicating whether the array and contents of fields should be re-used or internally allocated 45894f7d2dcSDave May 45994f7d2dcSDave May Currently, the only available projection method consists of 46094f7d2dcSDave May phi_i = \sum_{p=0}^{np} N_i(x_p) phi_p dJ / \sum_{p=0}^{np} N_i(x_p) dJ 46194f7d2dcSDave May where phi_p is the swarm field at point p, 46294f7d2dcSDave May N_i() is the cell DM basis function at vertex i, 46394f7d2dcSDave May dJ is the determinant of the cell Jacobian and 46494f7d2dcSDave May phi_i is the projected vertex value of the field phi. 46594f7d2dcSDave May 46694f7d2dcSDave May Level: beginner 46794f7d2dcSDave May 46894f7d2dcSDave May Notes: 469e7af74c8SDave May 470e7af74c8SDave May If reuse = PETSC_FALSE, this function will allocate the array of Vec's, and each individual Vec. 47194f7d2dcSDave May The user is responsible for destroying both the array and the individual Vec objects. 472e7af74c8SDave May 473e7af74c8SDave May Only swarm fields registered with data type = PETSC_REAL can be projected onto the cell DM. 474e7af74c8SDave May 475e7af74c8SDave May Only swarm fields of block size = 1 can currently be projected. 476e7af74c8SDave May 477e7af74c8SDave May The only projection methods currently only support the DA (2D) and PLEX (triangles 2D). 47894f7d2dcSDave May 47994f7d2dcSDave May .seealso: DMSwarmSetType(), DMSwarmSetCellDM(), DMSwarmType 48094f7d2dcSDave May @*/ 48194f7d2dcSDave May PETSC_EXTERN PetscErrorCode DMSwarmProjectFields(DM dm,PetscInt nfields,const char *fieldnames[],Vec **fields,PetscBool reuse) 4820e2ec84fSDave May { 48394f7d2dcSDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 48477048351SPatrick Sanan DMSwarmDataField *gfield; 48594f7d2dcSDave May DM celldm; 48694f7d2dcSDave May PetscBool isDA,isPLEX; 48794f7d2dcSDave May Vec *vecs; 48894f7d2dcSDave May PetscInt f,nvecs; 48994f7d2dcSDave May PetscInt project_type = 0; 49094f7d2dcSDave May 4910e2ec84fSDave May PetscFunctionBegin; 49294f7d2dcSDave May DMSWARMPICVALID(dm); 4935f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetCellDM(dm,&celldm)); 4945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nfields,&gfield)); 49594f7d2dcSDave May nvecs = 0; 49694f7d2dcSDave May for (f=0; f<nfields; f++) { 4975f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldnames[f],&gfield[f])); 4985f80ce2aSJacob Faibussowitsch PetscCheck(gfield[f]->petsc_type == PETSC_REAL,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Projection only valid for fields using a data type = PETSC_REAL"); 4995f80ce2aSJacob Faibussowitsch PetscCheck(gfield[f]->bs == 1,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Projection only valid for fields with block size = 1"); 50094f7d2dcSDave May nvecs += gfield[f]->bs; 50194f7d2dcSDave May } 50294f7d2dcSDave May if (!reuse) { 5035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nvecs,&vecs)); 50494f7d2dcSDave May for (f=0; f<nvecs; f++) { 5055f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(celldm,&vecs[f])); 5065f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)vecs[f],gfield[f]->name)); 50794f7d2dcSDave May } 50894f7d2dcSDave May } else { 50994f7d2dcSDave May vecs = *fields; 51094f7d2dcSDave May } 51194f7d2dcSDave May 5125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)celldm,DMDA,&isDA)); 5135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)celldm,DMPLEX,&isPLEX)); 51494f7d2dcSDave May if (isDA) { 5155f80ce2aSJacob Faibussowitsch CHKERRQ(private_DMSwarmProjectFields_DA(dm,celldm,project_type,nfields,gfield,vecs)); 51694f7d2dcSDave May } else if (isPLEX) { 5175f80ce2aSJacob Faibussowitsch CHKERRQ(private_DMSwarmProjectFields_PLEX(dm,celldm,project_type,nfields,gfield,vecs)); 51894f7d2dcSDave May } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only supported for cell DMs of type DMDA and DMPLEX"); 51994f7d2dcSDave May 5205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(gfield)); 5215f80ce2aSJacob Faibussowitsch if (!reuse) *fields = vecs; 5220e2ec84fSDave May PetscFunctionReturn(0); 5230e2ec84fSDave May } 5240e2ec84fSDave May 5250e2ec84fSDave May /*@C 526b799feefSDave May DMSwarmCreatePointPerCellCount - Count the number of points within all cells in the cell DM 5270e2ec84fSDave May 5280e2ec84fSDave May Not collective 5290e2ec84fSDave May 5300e2ec84fSDave May Input parameter: 5310e2ec84fSDave May . dm - the DMSwarm 5320e2ec84fSDave May 5330e2ec84fSDave May Output parameters: 5340e2ec84fSDave May + ncells - the number of cells in the cell DM (optional argument, pass NULL to ignore) 535b799feefSDave May - count - array of length ncells containing the number of points per cell 5360e2ec84fSDave May 5370e2ec84fSDave May Level: beginner 5380e2ec84fSDave May 5390e2ec84fSDave May Notes: 5400e2ec84fSDave May The array count is allocated internally and must be free'd by the user. 5410e2ec84fSDave May 5420e2ec84fSDave May .seealso: DMSwarmSetType(), DMSwarmSetCellDM(), DMSwarmType 5430e2ec84fSDave May @*/ 5440e2ec84fSDave May PETSC_EXTERN PetscErrorCode DMSwarmCreatePointPerCellCount(DM dm,PetscInt *ncells,PetscInt **count) 5450e2ec84fSDave May { 546b799feefSDave May PetscBool isvalid; 547b799feefSDave May PetscInt nel; 548b799feefSDave May PetscInt *sum; 549b799feefSDave May 5500e2ec84fSDave May PetscFunctionBegin; 5515f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSortGetIsValid(dm,&isvalid)); 552b799feefSDave May nel = 0; 553b799feefSDave May if (isvalid) { 554b799feefSDave May PetscInt e; 555b799feefSDave May 5565f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSortGetSizes(dm,&nel,NULL)); 557b799feefSDave May 5585f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nel,&sum)); 559b799feefSDave May for (e=0; e<nel; e++) { 5605f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSortGetNumberOfPointsPerCell(dm,e,&sum[e])); 561b799feefSDave May } 562b799feefSDave May } else { 563b799feefSDave May DM celldm; 564b799feefSDave May PetscBool isda,isplex,isshell; 565b799feefSDave May PetscInt p,npoints; 566b799feefSDave May PetscInt *swarm_cellid; 567b799feefSDave May 568b799feefSDave May /* get the number of cells */ 5695f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetCellDM(dm,&celldm)); 5705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)celldm,DMDA,&isda)); 5715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)celldm,DMPLEX,&isplex)); 5725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)celldm,DMSHELL,&isshell)); 573b799feefSDave May if (isda) { 574b799feefSDave May PetscInt _nel,_npe; 575b799feefSDave May const PetscInt *_element; 576b799feefSDave May 5775f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetElements(celldm,&_nel,&_npe,&_element)); 578b799feefSDave May nel = _nel; 5795f80ce2aSJacob Faibussowitsch CHKERRQ(DMDARestoreElements(celldm,&_nel,&_npe,&_element)); 580b799feefSDave May } else if (isplex) { 581b799feefSDave May PetscInt ps,pe; 582b799feefSDave May 5835f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(celldm,0,&ps,&pe)); 584b799feefSDave May nel = pe - ps; 585b799feefSDave May } else if (isshell) { 586b799feefSDave May PetscErrorCode (*method_DMShellGetNumberOfCells)(DM,PetscInt*); 587b799feefSDave May 5885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectQueryFunction((PetscObject)celldm,"DMGetNumberOfCells_C",&method_DMShellGetNumberOfCells)); 589b799feefSDave May if (method_DMShellGetNumberOfCells) { 5905f80ce2aSJacob Faibussowitsch CHKERRQ(method_DMShellGetNumberOfCells(celldm,&nel)); 591b799feefSDave May } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Cannot determine the number of cells for the DMSHELL object. User must provide a method via PetscObjectComposeFunction( (PetscObject)shelldm, \"DMGetNumberOfCells_C\", your_function_to_compute_number_of_cells);"); 592b799feefSDave May } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Cannot determine the number of cells for a DM not of type DA, PLEX or SHELL"); 593b799feefSDave May 5945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nel,&sum)); 5955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArrayzero(sum,nel)); 5965f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetLocalSize(dm,&npoints)); 5975f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid)); 598b799feefSDave May for (p=0; p<npoints; p++) { 599b799feefSDave May if (swarm_cellid[p] != DMLOCATEPOINT_POINT_NOT_FOUND) { 600b799feefSDave May sum[ swarm_cellid[p] ]++; 601b799feefSDave May } 602b799feefSDave May } 6035f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid)); 604b799feefSDave May } 605b799feefSDave May if (ncells) { *ncells = nel; } 606b799feefSDave May *count = sum; 6070e2ec84fSDave May PetscFunctionReturn(0); 6080e2ec84fSDave May } 60935b38c60SMatthew G. Knepley 61035b38c60SMatthew G. Knepley /*@ 61135b38c60SMatthew G. Knepley DMSwarmGetNumSpecies - Get the number of particle species 61235b38c60SMatthew G. Knepley 61335b38c60SMatthew G. Knepley Not collective 61435b38c60SMatthew G. Knepley 61535b38c60SMatthew G. Knepley Input parameter: 61635b38c60SMatthew G. Knepley . dm - the DMSwarm 61735b38c60SMatthew G. Knepley 61835b38c60SMatthew G. Knepley Output parameters: 61935b38c60SMatthew G. Knepley . Ns - the number of species 62035b38c60SMatthew G. Knepley 62135b38c60SMatthew G. Knepley Level: intermediate 62235b38c60SMatthew G. Knepley 62335b38c60SMatthew G. Knepley .seealso: DMSwarmSetNumSpecies(), DMSwarmSetType(), DMSwarmType 62435b38c60SMatthew G. Knepley @*/ 62535b38c60SMatthew G. Knepley PetscErrorCode DMSwarmGetNumSpecies(DM sw, PetscInt *Ns) 62635b38c60SMatthew G. Knepley { 62735b38c60SMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *) sw->data; 62835b38c60SMatthew G. Knepley 62935b38c60SMatthew G. Knepley PetscFunctionBegin; 63035b38c60SMatthew G. Knepley *Ns = swarm->Ns; 63135b38c60SMatthew G. Knepley PetscFunctionReturn(0); 63235b38c60SMatthew G. Knepley } 63335b38c60SMatthew G. Knepley 63435b38c60SMatthew G. Knepley /*@ 63535b38c60SMatthew G. Knepley DMSwarmSetNumSpecies - Set the number of particle species 63635b38c60SMatthew G. Knepley 63735b38c60SMatthew G. Knepley Not collective 63835b38c60SMatthew G. Knepley 63935b38c60SMatthew G. Knepley Input parameter: 64035b38c60SMatthew G. Knepley + dm - the DMSwarm 64135b38c60SMatthew G. Knepley - Ns - the number of species 64235b38c60SMatthew G. Knepley 64335b38c60SMatthew G. Knepley Level: intermediate 64435b38c60SMatthew G. Knepley 64535b38c60SMatthew G. Knepley .seealso: DMSwarmGetNumSpecies(), DMSwarmSetType(), DMSwarmType 64635b38c60SMatthew G. Knepley @*/ 64735b38c60SMatthew G. Knepley PetscErrorCode DMSwarmSetNumSpecies(DM sw, PetscInt Ns) 64835b38c60SMatthew G. Knepley { 64935b38c60SMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *) sw->data; 65035b38c60SMatthew G. Knepley 65135b38c60SMatthew G. Knepley PetscFunctionBegin; 65235b38c60SMatthew G. Knepley swarm->Ns = Ns; 65335b38c60SMatthew G. Knepley PetscFunctionReturn(0); 65435b38c60SMatthew G. Knepley } 65535b38c60SMatthew G. Knepley 65635b38c60SMatthew G. Knepley /*@C 65735b38c60SMatthew G. Knepley DMSwarmComputeLocalSize - Compute the local number and distribution of particles based upon a density function 65835b38c60SMatthew G. Knepley 65935b38c60SMatthew G. Knepley Not collective 66035b38c60SMatthew G. Knepley 66135b38c60SMatthew G. Knepley Input Parameters: 66235b38c60SMatthew G. Knepley + sw - The DMSwarm 66335b38c60SMatthew G. Knepley . N - The target number of particles 66435b38c60SMatthew G. Knepley - density - The density field for the particle layout, normalized to unity 66535b38c60SMatthew G. Knepley 66635b38c60SMatthew G. Knepley Note: One particle will be created for each species. 66735b38c60SMatthew G. Knepley 66835b38c60SMatthew G. Knepley Level: advanced 66935b38c60SMatthew G. Knepley 67035b38c60SMatthew G. Knepley .seealso: DMSwarmComputeLocalSizeFromOptions() 67135b38c60SMatthew G. Knepley @*/ 67235b38c60SMatthew G. Knepley PetscErrorCode DMSwarmComputeLocalSize(DM sw, PetscInt N, PetscProbFunc density) 67335b38c60SMatthew G. Knepley { 67435b38c60SMatthew G. Knepley DM dm; 67535b38c60SMatthew G. Knepley PetscQuadrature quad; 67635b38c60SMatthew G. Knepley const PetscReal *xq, *wq; 67735b38c60SMatthew G. Knepley PetscInt *npc, *cellid; 67835b38c60SMatthew G. Knepley PetscReal xi0[3], scale[1] = {.01}; 67935b38c60SMatthew G. Knepley PetscInt Ns, cStart, cEnd, c, dim, d, Nq, q, Np = 0, p; 68035b38c60SMatthew G. Knepley PetscBool simplex; 68135b38c60SMatthew G. Knepley 68235b38c60SMatthew G. Knepley PetscFunctionBegin; 6835f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetNumSpecies(sw, &Ns)); 6845f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetCellDM(sw, &dm)); 6855f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 6865f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 6875f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexIsSimplex(dm, &simplex)); 6885f80ce2aSJacob Faibussowitsch if (simplex) CHKERRQ(PetscDTStroudConicalQuadrature(dim, 1, 5, -1.0, 1.0, &quad)); 6895f80ce2aSJacob Faibussowitsch else CHKERRQ(PetscDTGaussTensorQuadrature(dim, 1, 5, -1.0, 1.0, &quad)); 6905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscQuadratureGetData(quad, NULL, NULL, &Nq, &xq, &wq)); 6915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(cEnd-cStart, &npc)); 69235b38c60SMatthew G. Knepley /* Integrate the density function to get the number of particles in each cell */ 69335b38c60SMatthew G. Knepley for (d = 0; d < dim; ++d) xi0[d] = -1.0; 69435b38c60SMatthew G. Knepley for (c = 0; c < cEnd-cStart; ++c) { 69535b38c60SMatthew G. Knepley const PetscInt cell = c + cStart; 69635b38c60SMatthew G. Knepley PetscReal v0[3], J[9], invJ[9], detJ; 69735b38c60SMatthew G. Knepley PetscReal n_int = 0.; 69835b38c60SMatthew G. Knepley 6995f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexComputeCellGeometryFEM(dm, cell, NULL, v0, J, invJ, &detJ)); 70035b38c60SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 70135b38c60SMatthew G. Knepley PetscReal xr[3], den[3]; 70235b38c60SMatthew G. Knepley 70335b38c60SMatthew G. Knepley CoordinatesRefToReal(dim, dim, xi0, v0, J, &xq[q*dim], xr); 7045f80ce2aSJacob Faibussowitsch CHKERRQ(density(xr, scale, den)); 70535b38c60SMatthew G. Knepley n_int += N*den[0]*wq[q]; 70635b38c60SMatthew G. Knepley } 70735b38c60SMatthew G. Knepley npc[c] = (PetscInt) n_int; 70835b38c60SMatthew G. Knepley npc[c] *= Ns; 70935b38c60SMatthew G. Knepley Np += npc[c]; 71035b38c60SMatthew G. Knepley } 7115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscQuadratureDestroy(&quad)); 7125f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSetLocalSizes(sw, Np, 0)); 71335b38c60SMatthew G. Knepley 7145f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid)); 71535b38c60SMatthew G. Knepley for (c = 0, p = 0; c < cEnd-cStart; ++c) { 71635b38c60SMatthew G. Knepley for (q = 0; q < npc[c]; ++q, ++p) cellid[p] = c; 71735b38c60SMatthew G. Knepley } 7185f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid)); 7195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(npc)); 72035b38c60SMatthew G. Knepley PetscFunctionReturn(0); 72135b38c60SMatthew G. Knepley } 72235b38c60SMatthew G. Knepley 72335b38c60SMatthew G. Knepley /*@ 72435b38c60SMatthew G. Knepley DMSwarmComputeLocalSizeFromOptions - Compute the local number and distribution of particles based upon a density function determined by options 72535b38c60SMatthew G. Knepley 72635b38c60SMatthew G. Knepley Not collective 72735b38c60SMatthew G. Knepley 72835b38c60SMatthew G. Knepley Input Parameters: 72935b38c60SMatthew G. Knepley , sw - The DMSwarm 73035b38c60SMatthew G. Knepley 73135b38c60SMatthew G. Knepley Level: advanced 73235b38c60SMatthew G. Knepley 73335b38c60SMatthew G. Knepley .seealso: DMSwarmComputeLocalSize() 73435b38c60SMatthew G. Knepley @*/ 73535b38c60SMatthew G. Knepley PetscErrorCode DMSwarmComputeLocalSizeFromOptions(DM sw) 73635b38c60SMatthew G. Knepley { 73735b38c60SMatthew G. Knepley DTProbDensityType den = DTPROB_DENSITY_CONSTANT; 73835b38c60SMatthew G. Knepley PetscProbFunc pdf; 73935b38c60SMatthew G. Knepley PetscInt N, Ns, dim; 74035b38c60SMatthew G. Knepley PetscBool flg; 74135b38c60SMatthew G. Knepley const char *prefix; 74235b38c60SMatthew G. Knepley PetscErrorCode ierr; 74335b38c60SMatthew G. Knepley 74435b38c60SMatthew G. Knepley PetscFunctionBegin; 74535b38c60SMatthew G. Knepley ierr = PetscOptionsBegin(PetscObjectComm((PetscObject) sw), "", "DMSwarm Options", "DMSWARM");CHKERRQ(ierr); 7465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-dm_swarm_num_particles", "The target number of particles", "", N, &N, NULL)); 7475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-dm_swarm_num_species", "The number of species", "DMSwarmSetNumSpecies", Ns, &Ns, &flg)); 7485f80ce2aSJacob Faibussowitsch if (flg) CHKERRQ(DMSwarmSetNumSpecies(sw, Ns)); 7495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsEnum("-dm_swarm_density", "Method to compute particle density <constant, gaussian>", "", DTProbDensityTypes, (PetscEnum) den, (PetscEnum *) &den, NULL)); 75035b38c60SMatthew G. Knepley ierr = PetscOptionsEnd();CHKERRQ(ierr); 75135b38c60SMatthew G. Knepley 7525f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(sw, &dim)); 7535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetOptionsPrefix((PetscObject) sw, &prefix)); 7545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscProbCreateFromOptions(dim, prefix, "-dm_swarm_coordinate_density", &pdf, NULL, NULL)); 7555f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmComputeLocalSize(sw, N, pdf)); 75635b38c60SMatthew G. Knepley PetscFunctionReturn(0); 75735b38c60SMatthew G. Knepley } 75835b38c60SMatthew G. Knepley 75935b38c60SMatthew G. Knepley /*@ 76035b38c60SMatthew G. Knepley DMSwarmInitializeCoordinates - Determine the initial coordinates of particles for a PIC method 76135b38c60SMatthew G. Knepley 76235b38c60SMatthew G. Knepley Not collective 76335b38c60SMatthew G. Knepley 76435b38c60SMatthew G. Knepley Input Parameters: 76535b38c60SMatthew G. Knepley , sw - The DMSwarm 76635b38c60SMatthew G. Knepley 76735b38c60SMatthew G. Knepley Note: Currently, we randomly place particles in their assigned cell 76835b38c60SMatthew G. Knepley 76935b38c60SMatthew G. Knepley Level: advanced 77035b38c60SMatthew G. Knepley 77135b38c60SMatthew G. Knepley .seealso: DMSwarmComputeLocalSize(), DMSwarmInitializeVelocities() 77235b38c60SMatthew G. Knepley @*/ 77335b38c60SMatthew G. Knepley PetscErrorCode DMSwarmInitializeCoordinates(DM sw) 77435b38c60SMatthew G. Knepley { 77535b38c60SMatthew G. Knepley DM dm; 77635b38c60SMatthew G. Knepley PetscRandom rnd; 77735b38c60SMatthew G. Knepley PetscScalar *weight; 77835b38c60SMatthew G. Knepley PetscReal *x, xi0[3]; 77935b38c60SMatthew G. Knepley PetscInt *species; 78035b38c60SMatthew G. Knepley PetscBool removePoints = PETSC_TRUE; 78135b38c60SMatthew G. Knepley PetscDataType dtype; 78235b38c60SMatthew G. Knepley PetscInt Ns, cStart, cEnd, c, dim, d, s, bs; 78335b38c60SMatthew G. Knepley 78435b38c60SMatthew G. Knepley PetscFunctionBeginUser; 7855f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetNumSpecies(sw, &Ns)); 7865f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetCellDM(sw, &dm)); 7875f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 7885f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 78935b38c60SMatthew G. Knepley 79035b38c60SMatthew G. Knepley /* Set particle position randomly in cell, set weights to 1 */ 7915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PetscObjectComm((PetscObject) dm), &rnd)); 7925f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetInterval(rnd, -1.0, 1.0)); 7935f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rnd)); 7945f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **) &weight)); 7955f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **) &x)); 7965f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(sw, "species", NULL, NULL, (void **) &species)); 7975f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSortGetAccess(sw)); 79835b38c60SMatthew G. Knepley for (d = 0; d < dim; ++d) xi0[d] = -1.0; 79935b38c60SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 80035b38c60SMatthew G. Knepley PetscReal v0[3], J[9], invJ[9], detJ; 80135b38c60SMatthew G. Knepley PetscInt *pidx, Npc, q; 80235b38c60SMatthew G. Knepley 8035f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSortGetPointsPerCell(sw, c, &Npc, &pidx)); 8045f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); 80535b38c60SMatthew G. Knepley for (q = 0; q < Npc; ++q) { 80635b38c60SMatthew G. Knepley const PetscInt p = pidx[q]; 80735b38c60SMatthew G. Knepley PetscReal xref[3]; 80835b38c60SMatthew G. Knepley 8095f80ce2aSJacob Faibussowitsch for (d = 0; d < dim; ++d) CHKERRQ(PetscRandomGetValueReal(rnd, &xref[d])); 81035b38c60SMatthew G. Knepley CoordinatesRefToReal(dim, dim, xi0, v0, J, xref, &x[p*dim]); 81135b38c60SMatthew G. Knepley 81235b38c60SMatthew G. Knepley weight[p] = 1.0; 81335b38c60SMatthew G. Knepley for (s = 0; s < Ns; ++s) species[p] = p % Ns; 81435b38c60SMatthew G. Knepley } 8155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(pidx)); 81635b38c60SMatthew G. Knepley } 8175f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rnd)); 8185f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSortRestoreAccess(sw)); 8195f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **) &weight)); 8205f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **) &x)); 8215f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **) &species)); 8225f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmMigrate(sw, removePoints)); 8235f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalizeCoordinates(sw)); 82435b38c60SMatthew G. Knepley PetscFunctionReturn(0); 82535b38c60SMatthew G. Knepley } 82635b38c60SMatthew G. Knepley 82735b38c60SMatthew G. Knepley /*@C 82835b38c60SMatthew G. Knepley DMSwarmInitializeVelocities - Set the initial velocities of particles using a distribution. 82935b38c60SMatthew G. Knepley 83035b38c60SMatthew G. Knepley Collective on dm 83135b38c60SMatthew G. Knepley 83235b38c60SMatthew G. Knepley Input Parameters: 83335b38c60SMatthew G. Knepley + sw - The DMSwarm object 83435b38c60SMatthew G. Knepley . sampler - A function which uniformly samples the velocity PDF 83535b38c60SMatthew G. Knepley - v0 - The velocity scale for nondimensionalization for each species 83635b38c60SMatthew G. Knepley 83735b38c60SMatthew G. Knepley Level: advanced 83835b38c60SMatthew G. Knepley 83935b38c60SMatthew G. Knepley .seealso: DMSwarmComputeLocalSize(), DMSwarmInitializeCoordinates(), DMSwarmInitializeVelocitiesFromOptions() 84035b38c60SMatthew G. Knepley @*/ 84135b38c60SMatthew G. Knepley PetscErrorCode DMSwarmInitializeVelocities(DM sw, PetscProbFunc sampler, const PetscReal v0[]) 84235b38c60SMatthew G. Knepley { 84335b38c60SMatthew G. Knepley PetscRandom rnd; 84435b38c60SMatthew G. Knepley PetscReal *v; 84535b38c60SMatthew G. Knepley PetscInt *species; 84635b38c60SMatthew G. Knepley PetscInt dim, Np, p; 84735b38c60SMatthew G. Knepley 84835b38c60SMatthew G. Knepley PetscFunctionBegin; 8495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PetscObjectComm((PetscObject) sw), &rnd)); 8505f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetInterval(rnd, 0, 1.)); 8515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rnd)); 85235b38c60SMatthew G. Knepley 8535f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(sw, &dim)); 8545f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetLocalSize(sw, &Np)); 8555f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **) &v)); 8565f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(sw, "species", NULL, NULL, (void **) &species)); 85735b38c60SMatthew G. Knepley for (p = 0; p < Np; ++p) { 85835b38c60SMatthew G. Knepley PetscInt s = species[p], d; 85935b38c60SMatthew G. Knepley PetscReal a[3], vel[3]; 86035b38c60SMatthew G. Knepley 8615f80ce2aSJacob Faibussowitsch for (d = 0; d < dim; ++d) CHKERRQ(PetscRandomGetValueReal(rnd, &a[d])); 8625f80ce2aSJacob Faibussowitsch CHKERRQ(sampler(a, NULL, vel)); 86335b38c60SMatthew G. Knepley for (d = 0; d < dim; ++d) {v[p*dim+d] = (v0[s] / v0[0]) * vel[d];} 86435b38c60SMatthew G. Knepley } 8655f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **) &v)); 8665f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **) &species)); 8675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rnd)); 86835b38c60SMatthew G. Knepley PetscFunctionReturn(0); 86935b38c60SMatthew G. Knepley } 87035b38c60SMatthew G. Knepley 87135b38c60SMatthew G. Knepley /*@ 87235b38c60SMatthew G. Knepley DMSwarmInitializeVelocitiesFromOptions - Set the initial velocities of particles using a distribution determined from options. 87335b38c60SMatthew G. Knepley 87435b38c60SMatthew G. Knepley Collective on dm 87535b38c60SMatthew G. Knepley 87635b38c60SMatthew G. Knepley Input Parameters: 87735b38c60SMatthew G. Knepley + sw - The DMSwarm object 87835b38c60SMatthew G. Knepley - v0 - The velocity scale for nondimensionalization for each species 87935b38c60SMatthew G. Knepley 88035b38c60SMatthew G. Knepley Level: advanced 88135b38c60SMatthew G. Knepley 88235b38c60SMatthew G. Knepley .seealso: DMSwarmComputeLocalSize(), DMSwarmInitializeCoordinates(), DMSwarmInitializeVelocities() 88335b38c60SMatthew G. Knepley @*/ 88435b38c60SMatthew G. Knepley PetscErrorCode DMSwarmInitializeVelocitiesFromOptions(DM sw, const PetscReal v0[]) 88535b38c60SMatthew G. Knepley { 88635b38c60SMatthew G. Knepley PetscProbFunc sampler; 88735b38c60SMatthew G. Knepley PetscInt dim; 88835b38c60SMatthew G. Knepley const char *prefix; 88935b38c60SMatthew G. Knepley 89035b38c60SMatthew G. Knepley PetscFunctionBegin; 8915f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(sw, &dim)); 8925f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetOptionsPrefix((PetscObject) sw, &prefix)); 8935f80ce2aSJacob Faibussowitsch CHKERRQ(PetscProbCreateFromOptions(dim, prefix, "-dm_swarm_velocity_density", NULL, NULL, &sampler)); 8945f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmInitializeVelocities(sw, sampler, v0)); 89535b38c60SMatthew G. Knepley PetscFunctionReturn(0); 89635b38c60SMatthew G. Knepley } 897