dalcinl / PetIGA

PETSc-based framework for isogeometric analysis (IGA) methods

Geek Repo:Geek Repo

Github PK Tool:Github PK Tool

Possible issue with IGALocateElement in older version / no IGALocateElement for new PetIGA

shaunakshende opened this issue · comments

I am using the new version of PetIGA posted here, but I need to use IGALocateElement for integration of an embedded manifold. As I have been using a uniform discretization in the past, the element index in each direction can be found through (PetscInt) ParametricCoordinate/knot_parametric_coordinate (so it is floored). However, the previous implementation of IGALocateElement does not have this behavior, because points on the right edges of elements are assigned an index corresponding to the next element. I think this is incorrect, as if I assign ranks to my foreground particles in this way, probe cannot find them. I think

if(pnt[i] > U[j] && pnt[i] <= U[j+1]) ID[i] = e;

has to be changed to

if((pnt[i] >= U[j] && pnt[i] < U[j+1]) || sqrt((pnt[i]-1.0)*(pnt[i]-1.0)) ) ID[i] = e;

in this function to be consistent with the IGA partition. This came up for me when I switched to a non-uniform mesh and needed to actually search for the element.

#undef  __FUNCT__
#define __FUNCT__ "IGALocateElement_1"
PetscBool IGALocateElement_1(IGA iga,PetscReal *pnt,IGAElement element, AppCtx *user, ParticleManager &manager)
{
  PetscFunctionBegin;
  user->debug = PETSC_FALSE;
  PetscMPIInt rank;
  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
  PetscInt i,j,e,m,dim=iga->dim;
  PetscInt *ID = element->ID;
  PetscScalar *U;

  element->nen  = 1;
  element->nval = 0;
  element->nvec = 0;
  element->nmat = 0;
  for(i=0;i<dim;i++){
    element->nen *= (iga->axis[i]->p+1);
    m = iga->axis[i]->m;
    U = iga->axis[i]->U;
    e = -1;
    ID[i] = 0;
    /* find which nonzero span this point is located in */
    for(j=0;j<m;j++){
      if(U[j+1]-U[j]>1.0e-13) e += 1;
      if((pnt[i] >= U[j] && pnt[i] < U[j+1]) || sqrt((pnt[i]-1.0)*(pnt[i]-1.0)) ){ID[i] = e;
      break;}
    }
    /* reject if the element is not in this partition */
    if(ID[i] < iga->elem_start[i] || ID[i] >= iga->elem_start[i]+iga->elem_width[i]){
    PetscPrintf(PETSC_COMM_SELF, "Error: Particle on domain but not located by IGALocateElement! This Probably means that the rank of the point \n and the rank of the Subdomain assigned to rank are not consistent when this function is being called! Check pointRank. Exiting...\n");
    user->debug=PETSC_TRUE;
    int CalculatedRank = manager.pointRank(pnt, user);
    PetscPrintf(PETSC_COMM_SELF, "Rank is %d elemX = %d elemY = %d and elemZ = %d. Calculated Rank is %d!\n", rank, ID[0], ID[1], ID[2], CalculatedRank);
    MPI_Barrier(PETSC_COMM_WORLD);
    exit(1);}
  }

  element->index = 0;
{
    PetscErrorCode ierr;
    ierr = IGAElementBuildClosure(element);CHKERRCONTINUE(ierr);
    if (PetscUnlikely(ierr)) PetscFunctionReturn(PETSC_FALSE);
    ierr = IGAElementBuildFix(element);CHKERRCONTINUE(ierr);
    if (PetscUnlikely(ierr)) PetscFunctionReturn(PETSC_FALSE);
  }

  return PETSC_TRUE;
}

I removed IGALocateElement() simply because it was broken. That happened long ago, and I do not remember the details by now. The needs I had at the time for IGALocateElement() were replaced with IGAProbe and related routines. If you want to reinstate the old locate element routine (or a refinement of it) for your own needs, you are free to do so, but the sqrt call you are using looks a bit suspicious.

The version I posted fixed it for me, but yes there was a typo I meant sqrt((pnt[i]-1.0)*(pnt[i]-1.0))<1.0e-13 (points on the right edge are in element index N-1). Just wondering, is there another routine associated with probe that can build the element around a point like IGALocateElement?