#include #include #include #include "../src/dm/impls/swarm/data_bucket.h" PetscErrorCode private_DMSwarmSetPointCoordinatesCellwise_PLEX(DM,DM,PetscInt,PetscReal*xi); static PetscErrorCode private_PetscFECreateDefault_scalar_pk1(DM dm, PetscInt dim, PetscBool isSimplex, PetscInt qorder, PetscFE *fem) { const PetscInt Nc = 1; PetscQuadrature q, fq; DM K; PetscSpace P; PetscDualSpace Q; PetscInt order, quadPointsPerEdge; PetscBool tensor = isSimplex ? PETSC_FALSE : PETSC_TRUE; PetscFunctionBegin; /* Create space */ PetscCall(PetscSpaceCreate(PetscObjectComm((PetscObject) dm), &P)); /* PetscCall(PetscObjectSetOptionsPrefix((PetscObject) P, prefix)); */ PetscCall(PetscSpacePolynomialSetTensor(P, tensor)); /* PetscCall(PetscSpaceSetFromOptions(P)); */ PetscCall(PetscSpaceSetType(P,PETSCSPACEPOLYNOMIAL)); PetscCall(PetscSpaceSetDegree(P,1,PETSC_DETERMINE)); PetscCall(PetscSpaceSetNumComponents(P, Nc)); PetscCall(PetscSpaceSetNumVariables(P, dim)); PetscCall(PetscSpaceSetUp(P)); PetscCall(PetscSpaceGetDegree(P, &order, NULL)); PetscCall(PetscSpacePolynomialGetTensor(P, &tensor)); /* Create dual space */ PetscCall(PetscDualSpaceCreate(PetscObjectComm((PetscObject) dm), &Q)); PetscCall(PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE)); /* PetscCall(PetscObjectSetOptionsPrefix((PetscObject) Q, prefix)); */ PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, DMPolytopeTypeSimpleShape(dim, isSimplex), &K)); PetscCall(PetscDualSpaceSetDM(Q, K)); PetscCall(DMDestroy(&K)); PetscCall(PetscDualSpaceSetNumComponents(Q, Nc)); PetscCall(PetscDualSpaceSetOrder(Q, order)); PetscCall(PetscDualSpaceLagrangeSetTensor(Q, tensor)); /* PetscCall(PetscDualSpaceSetFromOptions(Q)); */ PetscCall(PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE)); PetscCall(PetscDualSpaceSetUp(Q)); /* Create element */ PetscCall(PetscFECreate(PetscObjectComm((PetscObject) dm), fem)); /* PetscCall(PetscObjectSetOptionsPrefix((PetscObject) *fem, prefix)); */ /* PetscCall(PetscFESetFromOptions(*fem)); */ PetscCall(PetscFESetType(*fem,PETSCFEBASIC)); PetscCall(PetscFESetBasisSpace(*fem, P)); PetscCall(PetscFESetDualSpace(*fem, Q)); PetscCall(PetscFESetNumComponents(*fem, Nc)); PetscCall(PetscFESetUp(*fem)); PetscCall(PetscSpaceDestroy(&P)); PetscCall(PetscDualSpaceDestroy(&Q)); /* Create quadrature (with specified order if given) */ qorder = qorder >= 0 ? qorder : order; quadPointsPerEdge = PetscMax(qorder + 1,1); if (isSimplex) { PetscCall(PetscDTStroudConicalQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0, &q)); PetscCall(PetscDTStroudConicalQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq)); } else { PetscCall(PetscDTGaussTensorQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0, &q)); PetscCall(PetscDTGaussTensorQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq)); } PetscCall(PetscFESetQuadrature(*fem, q)); PetscCall(PetscFESetFaceQuadrature(*fem, fq)); PetscCall(PetscQuadratureDestroy(&q)); PetscCall(PetscQuadratureDestroy(&fq)); PetscFunctionReturn(0); } PetscErrorCode subdivide_triangle(PetscReal v1[2],PetscReal v2[2],PetscReal v3[2],PetscInt depth,PetscInt max,PetscReal xi[],PetscInt *np) { PetscReal v12[2],v23[2],v31[2]; PetscInt i; PetscFunctionBegin; if (depth == max) { PetscReal cx[2]; cx[0] = (v1[0] + v2[0] + v3[0])/3.0; cx[1] = (v1[1] + v2[1] + v3[1])/3.0; xi[2*(*np)+0] = cx[0]; xi[2*(*np)+1] = cx[1]; (*np)++; PetscFunctionReturn(0); } /* calculate midpoints of each side */ for (i = 0; i < 2; i++) { v12[i] = (v1[i]+v2[i])/2.0; v23[i] = (v2[i]+v3[i])/2.0; v31[i] = (v3[i]+v1[i])/2.0; } /* recursively subdivide new triangles */ PetscCall(subdivide_triangle(v1,v12,v31,depth+1,max,xi,np)); PetscCall(subdivide_triangle(v2,v23,v12,depth+1,max,xi,np)); PetscCall(subdivide_triangle(v3,v31,v23,depth+1,max,xi,np)); PetscCall(subdivide_triangle(v12,v23,v31,depth+1,max,xi,np)); PetscFunctionReturn(0); } PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX2D_SubDivide(DM dm,DM dmc,PetscInt nsub) { const PetscInt dim = 2; PetscInt q,npoints_q,e,nel,npe,pcnt,ps,pe,d,k,depth; PetscReal *xi; PetscReal **basis; Vec coorlocal; PetscSection coordSection; PetscScalar *elcoor = NULL; PetscReal *swarm_coor; PetscInt *swarm_cellid; PetscReal v1[2],v2[2],v3[2]; PetscFunctionBegin; npoints_q = 1; for (d=0; dcell, 1->edge, 2->vert */ PetscCall(DMPlexGetHeightStratum(dmc,0,&ps,&pe)); nel = pe - ps; PetscCall(DMSwarmSetLocalSizes(dm,npoints_q*nel,-1)); PetscCall(DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor)); PetscCall(DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid)); PetscCall(DMGetCoordinatesLocal(dmc,&coorlocal)); PetscCall(DMGetCoordinateSection(dmc,&coordSection)); pcnt = 0; for (e=0; ecell, 1->edge, 2->vert */ PetscCall(DMPlexGetHeightStratum(dmc,0,&ps,&pe)); nel = pe - ps; PetscCall(DMSwarmSetLocalSizes(dm,npoints_q*nel,-1)); PetscCall(DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor)); PetscCall(DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid)); PetscCall(DMGetCoordinatesLocal(dmc,&coorlocal)); PetscCall(DMGetCoordinateSection(dmc,&coordSection)); pcnt = 0; for (e=0; eT[0][q*nbasis + k] * PetscRealPart(elcoor[dim*k+d]); } } swarm_cellid[pcnt] = e; pcnt++; } PetscCall(DMPlexVecRestoreClosure(dmc,coordSection,coorlocal,ps+e,NULL,&elcoor)); } PetscCall(DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid)); PetscCall(DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor)); PetscCall(PetscFEDestroy(&fe)); PetscFunctionReturn(0); } PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX2D_Regular(DM dm,DM dmc,PetscInt npoints) { PetscInt dim; PetscInt ii,jj,q,npoints_q,e,nel,npe,pcnt,ps,pe,d,k,nfaces; PetscReal *xi,ds,ds2; PetscReal **basis; Vec coorlocal; PetscSection coordSection; PetscScalar *elcoor = NULL; PetscReal *swarm_coor; PetscInt *swarm_cellid; PetscBool is_simplex; PetscFunctionBegin; PetscCall(DMGetDimension(dmc,&dim)); PetscCheckFalse(dim != 2,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only 2D is supported"); is_simplex = PETSC_FALSE; PetscCall(DMPlexGetHeightStratum(dmc,0,&ps,&pe)); PetscCall(DMPlexGetConeSize(dmc, ps, &nfaces)); if (nfaces == (dim+1)) { is_simplex = PETSC_TRUE; } PetscCheck(is_simplex,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only the simplex is supported"); PetscCall(PetscMalloc1(dim*npoints*npoints,&xi)); pcnt = 0; ds = 1.0/((PetscReal)(npoints-1)); ds2 = 1.0/((PetscReal)(npoints)); for (jj = 0; jjcell, 1->edge, 2->vert */ PetscCall(DMPlexGetHeightStratum(dmc,0,&ps,&pe)); nel = pe - ps; PetscCall(DMSwarmSetLocalSizes(dm,npoints_q*nel,-1)); PetscCall(DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor)); PetscCall(DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid)); PetscCall(DMGetCoordinatesLocal(dmc,&coorlocal)); PetscCall(DMGetCoordinateSection(dmc,&coordSection)); pcnt = 0; for (e=0; e mapped to element %D xi(%+1.4e,%+1.4e)\n",p,point.x,point.y,e,xi_p[0],xi_p[1]); } */ PetscCall(ComputeLocalCoordinateAffine2d(coor_p,elcoor,xi_p,&dJ)); /* PetscPrintf(PETSC_COMM_SELF,"[p=%D] x(%+1.4e,%+1.4e) -> mapped to element %D xi(%+1.4e,%+1.4e)\n",p,point.x,point.y,e,xi_p[0],xi_p[1]); */ /* point_located = PETSC_TRUE; if (xi_p[0] < 0.0) { if (xi_p[0] > -PLEX_C_EPS) { xi_p[0] = 0.0; } else { point_located = PETSC_FALSE; } } if (xi_p[1] < 0.0) { if (xi_p[1] > -PLEX_C_EPS) { xi_p[1] = 0.0; } else { point_located = PETSC_FALSE; } } if (xi_p[1] > (1.0-xi_p[0])) { if ((xi_p[1] - 1.0 + xi_p[0]) < PLEX_C_EPS) { xi_p[1] = 1.0 - xi_p[0]; } else { point_located = PETSC_FALSE; } } if (!point_located) { PetscPrintf(PETSC_COMM_SELF,"[Error] xi,eta = %+1.8e, %+1.8e\n",xi_p[0],xi_p[1]); PetscPrintf(PETSC_COMM_SELF,"[Error] Failed to locate point (%1.8e,%1.8e) in local mesh (cell %D) with triangle coords (%1.8e,%1.8e) : (%1.8e,%1.8e) : (%1.8e,%1.8e)\n",coor_p[0],coor_p[1],e,elcoor[0],elcoor[1],elcoor[2],elcoor[3],elcoor[4],elcoor[5]); } PetscCheck(point_located,PETSC_COMM_SELF,PETSC_ERR_SUP,"Failed to locate point (%1.8e,%1.8e) in local mesh (cell %D)",coor_p[0],coor_p[1],e); */ Ni[0] = 1.0 - xi_p[0] - xi_p[1]; Ni[1] = xi_p[0]; Ni[2] = xi_p[1]; point_located = PETSC_TRUE; for (k=0; k<3; k++) { if (PetscRealPart(Ni[k]) < -PLEX_C_EPS) point_located = PETSC_FALSE; if (PetscRealPart(Ni[k]) > (1.0+PLEX_C_EPS)) point_located = PETSC_FALSE; } if (!point_located) { PetscCall(PetscPrintf(PETSC_COMM_SELF,"[Error] xi,eta = %+1.8e, %+1.8e\n",(double)xi_p[0],(double)xi_p[1])); PetscCall(PetscPrintf(PETSC_COMM_SELF,"[Error] Failed to locate point (%1.8e,%1.8e) in local mesh (cell %D) with triangle coords (%1.8e,%1.8e) : (%1.8e,%1.8e) : (%1.8e,%1.8e)\n",(double)coor_p[0],(double)coor_p[1],e,(double)PetscRealPart(elcoor[0]),(double)PetscRealPart(elcoor[1]),(double)PetscRealPart(elcoor[2]),(double)PetscRealPart(elcoor[3]),(double)PetscRealPart(elcoor[4]),(double)PetscRealPart(elcoor[5]))); } PetscCheck(point_located,PETSC_COMM_SELF,PETSC_ERR_SUP,"Failed to locate point (%1.8e,%1.8e) in local mesh (cell %D)",(double)coor_p[0],(double)coor_p[1],e); for (k=0; k<3; k++) { Ni[k] = Ni[k] * dJ; elfield[k] = Ni[k] * swarm_field[p]; } PetscCall(DMPlexVecRestoreClosure(dm,coordSection,coor_l,e,NULL,&elcoor)); PetscCall(DMPlexVecSetClosure(dm, NULL,v_field_l, e, elfield, ADD_VALUES)); PetscCall(DMPlexVecSetClosure(dm, NULL,denom_l, e, Ni, ADD_VALUES)); } PetscCall(DMSwarmRestoreField(swarm,DMSwarmPICField_cellid,NULL,NULL,(void**)&mpfield_cell)); PetscCall(DMSwarmRestoreField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&mpfield_coor)); PetscCall(DMLocalToGlobalBegin(dm,v_field_l,ADD_VALUES,v_field)); PetscCall(DMLocalToGlobalEnd(dm,v_field_l,ADD_VALUES,v_field)); PetscCall(DMLocalToGlobalBegin(dm,denom_l,ADD_VALUES,denom)); PetscCall(DMLocalToGlobalEnd(dm,denom_l,ADD_VALUES,denom)); PetscCall(VecPointwiseDivide(v_field,v_field,denom)); PetscCall(DMRestoreLocalVector(dm,&v_field_l)); PetscCall(DMRestoreLocalVector(dm,&denom_l)); PetscCall(DMRestoreGlobalVector(dm,&denom)); PetscFunctionReturn(0); } PetscErrorCode private_DMSwarmProjectFields_PLEX(DM swarm,DM celldm,PetscInt project_type,PetscInt nfields,DMSwarmDataField dfield[],Vec vecs[]) { PetscInt f,dim; PetscFunctionBegin; PetscCall(DMGetDimension(swarm,&dim)); switch (dim) { case 2: for (f=0; f 0.0,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Points do not fail inside the simplex domain"); } } else if (is_tensorcell) { for (p=0; p 1.0,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Points do not fail inside the tensor domain [-1,1]^d"); } } } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only support for d-simplex and d-tensorcell"); PetscCall(PetscQuadratureCreate(PetscObjectComm((PetscObject)dm),&quadrature)); PetscCall(PetscMalloc1(npoints*dim,&xiq)); PetscCall(PetscArraycpy(xiq,xi,npoints*dim)); PetscCall(PetscQuadratureSetData(quadrature,dim,1,npoints,(const PetscReal*)xiq,NULL)); PetscCall(private_PetscFECreateDefault_scalar_pk1(dmc, dim, is_simplex, 0, &fe)); PetscCall(PetscFESetQuadrature(fe,quadrature)); PetscCall(PetscFEGetDimension(fe,&nbasis)); PetscCall(PetscFEGetCellTabulation(fe, 1, &T)); /* for each cell, interpolate coordaintes and insert the interpolated points coordinates into swarm */ /* 0->cell, 1->edge, 2->vert */ PetscCall(DMPlexGetHeightStratum(dmc,0,&ps,&pe)); nel = pe - ps; PetscCall(DMSwarmSetLocalSizes(dm,npoints*nel,-1)); PetscCall(DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor)); PetscCall(DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid)); PetscCall(DMGetCoordinatesLocal(dmc,&coorlocal)); PetscCall(DMGetCoordinateSection(dmc,&coordSection)); pcnt = 0; for (e=0; eT[0][p*nbasis + k] * PetscRealPart(elcoor[dim*k+d]); } } swarm_cellid[pcnt] = e; pcnt++; } PetscCall(DMPlexVecRestoreClosure(dmc,coordSection,coorlocal,ps+e,NULL,&elcoor)); } PetscCall(DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid)); PetscCall(DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor)); PetscCall(PetscQuadratureDestroy(&quadrature)); PetscCall(PetscFEDestroy(&fe)); PetscFunctionReturn(0); }