/** Create a cache for Dirichlet part of closure vector, and scatter from global closure to Dirichlet cache. @arg[in] gvec Global vector @arg[out] dcache New vector to hold the Dirichlet values @arg[out] dscat Scatter from global closure to \a dcache @note This could be local but it doesn't cost anything to make it global. **/ dErr VecDohpCreateDirichletCache(Vec gvec,Vec *dcache,VecScatter *dscat) { MPI_Comm comm; dErr err; dBool isdohp; IS from; Vec gc; dInt n,nc,crstart; dFunctionBegin; dValidHeader(gvec,VEC_CLASSID,1); dValidPointer(dcache,2); dValidPointer(dscat,3); err = PetscTypeCompare((PetscObject)gvec,VECDOHP,&isdohp);dCHK(err); if (!isdohp) dERROR(PETSC_COMM_SELF,PETSC_ERR_SUP,"Vec type %s",((PetscObject)gvec)->type_name); err = PetscObjectGetComm((PetscObject)gvec,&comm);dCHK(err); err = VecGetLocalSize(gvec,&n);dCHK(err); err = VecDohpGetClosure(gvec,&gc);dCHK(err); err = VecGetLocalSize(gc,&nc);dCHK(err); err = VecGetOwnershipRange(gc,&crstart,NULL);dCHK(err); err = VecCreateMPI(comm,nc-n,PETSC_DECIDE,dcache);dCHK(err); err = ISCreateStride(comm,nc-n,crstart+n,1,&from);dCHK(err); err = VecScatterCreate(gc,from,*dcache,NULL,dscat);dCHK(err); err = VecDohpRestoreClosure(gvec,&gc);dCHK(err); err = ISDestroy(&from);dCHK(err); /* \todo deal with rotations */ dFunctionReturn(0); }
// Logically collective dErr dUnitsCreateUnit(dUnits un,const char *type,const char *longname,const char *shortname,dInt n,const dReal expon[],dUnit *newunit) { dErr err; dUnit unit; dFunctionBegin; dValidHeader(un,dUNITS_CLASSID,1); if (n < 1 || n > dUNITS_MAX) dERROR(((dObject)un)->comm,PETSC_ERR_ARG_OUTOFRANGE,"The number of exponents %D must be positive, but no larger than %D",n,(dInt)dUNITS_MAX); dValidRealPointer(expon,5); dValidPointer(newunit,6); err = dUnitsGetEmptyUnit_Private(un,&unit);dCHK(err); err = PetscStrallocpy(type,&unit->quantity);dCHK(err); err = dUnitsAssignName(un,dUnitName,longname,n,expon,&unit->longname);dCHK(err); err = dUnitsAssignName(un,dUnitShortName,shortname,n,expon,&unit->shortname);dCHK(err); err = dUnitsAssignName(un,dUnitSIName,NULL,n,expon,&unit->siname);dCHK(err); unit->toSI = 1.0; unit->toCommon = 1.0; for (dInt i=0; i<n; i++) { dUnit base; err = dUnitsGetBase(un,i,&base);dCHK(err); unit->toCommon *= PetscPowScalar(dUnitDimensionalize(base,1.0),expon[i]); unit->toSI *= PetscPowScalar(dUnitDimensionalizeSI(base,1.0),expon[i]); unit->expon[i] = expon[i]; } *newunit = unit; dFunctionReturn(0); }
/** Set rotation for certain nodes in a function space * * @param fs the function space * @param is index set of nodes to rotate, sequential, with respect blocks of local vector * @param rot Rotation matrices at all nodes in \a is. Should have length \c bs*bs*size(is). * @param ns number of dofs to enforce strongly at each node (every entry must have 0<=ns[i]<=bs) * @param v Vector of values for strongly enforced dofs * * @example Consider 2D flow over a bed with known melt rates. Suppose the local velocity vector is * * [u0x,u0y; u1x,u1y; u2x,u2y; u3x,u3y | u4x,u4y] * * (4 owned blocks, one ghosted block) and nodes 1,4 are on the slip boundary with normal and tangent vectors n1,t1,n4,t4 * and melt rates r1,r4. To enforce the melt rate strongly, use * * \a is = [1,4] * \a rot = [n10,n11,t10,t11, n40,n41,t40,t41] * \a ns = [1,1] * \a v = [r1,r4] * * The rotated vector will become (. = \cdot) * * [u0x,u0y; u1.n1,u1.t1; u2x,u2y; u3x,u3y | u4.n4,u4.t4] * * and strongly enforcing melt rate produces the global vector * * [u0x,u0y; r1,u1.t1; u2x,u2y; u3x,u3y | r4,u4.t4] . * * This is what the solver sees, the Jacobian will always have rows and columns of the identity corresponding to the * strongly enforced components (2,8 of the local vector) and the residual will always be 0 in these components. Hence * the Newton step v will always be of the form * * [v0x,v0y; 0,v1y; v2x,v2y; v3x,v3y | 0,v4y] . **/ dErr dFSRotationCreate(dFS fs,IS is,dReal rmat[],dInt ns[],Vec v,dFSRotation *inrot) { dFSRotation rot; dInt bs,n; dErr err; dFunctionBegin; dValidHeader(fs,DM_CLASSID,1); dValidHeader(is,IS_CLASSID,2); dValidRealPointer(rmat,3); dValidIntPointer(ns,4); dValidHeader(v,VEC_CLASSID,5); dValidPointer(inrot,6); *inrot = 0; err = PetscHeaderCreate(rot,_p_dFSRotation,struct _dFSRotationOps,dFSROT_CLASSID,"dFSRotation","Local function space rotation","FS",PETSC_COMM_SELF,dFSRotationDestroy,dFSRotationView);dCHK(err); err = dFSGetBlockSize(fs,&bs);dCHK(err); rot->bs = bs; err = ISGetSize(is,&n);dCHK(err); rot->n = n; err = PetscObjectReference((PetscObject)is);dCHK(err); rot->is = is; err = PetscObjectReference((PetscObject)v);dCHK(err); rot->strong = v; for (dInt i=0; i<n; i++) { if (ns[i] < 0 || bs < ns[i]) dERROR(PETSC_COMM_SELF,1,"Number of strong dofs must be between 0 and bs=%d (inclusive)",bs); /* \todo Check that every rmat is orthogonal */ } err = dMallocA2(n*bs*bs,&rot->rmat,n,&rot->nstrong);dCHK(err); err = dMemcpy(rot->rmat,rmat,n*bs*bs*sizeof rmat[0]);dCHK(err); err = dMemcpy(rot->nstrong,ns,n*sizeof ns[0]);dCHK(err); *inrot = rot; dFunctionReturn(0); }
dErr VecCreateDohp(MPI_Comm comm,dInt bs,dInt n,dInt nc,dInt nghosts,const dInt ghosts[],Vec *v) { Vec_MPI *vmpi; Vec vc,vg; dScalar *a; dErr err; dFunctionBegin; dValidPointer(v,7); *v = 0; err = VecCreateGhostBlock(comm,bs,nc*bs,PETSC_DECIDE,nghosts,ghosts,&vc);dCHK(err); err = VecGetArray(vc,&a);dCHK(err); err = VecCreateMPIWithArray(comm,n*bs,PETSC_DECIDE,a,&vg);dCHK(err); err = VecRestoreArray(vc,&a);dCHK(err); err = VecSetBlockSize(vg,bs);dCHK(err); vmpi = vg->data; if (vmpi->localrep) dERROR(PETSC_COMM_SELF,1,"Vector has localrep, expected no localrep"); vmpi->localrep = vc; /* subvert this field to mean closed rep */ /* Since we subvect .localrep, VecDestroy_MPI will automatically destroy the closed form */ vg->ops->duplicate = VecDuplicate_Dohp; //vg->ops->destroy = VecDestroy_Dohp; /* It might be useful to set the (block) LocalToGlobal mapping here, but in the use case I have in mind, the user is * always working with the closed form anyway (in function evaluation). The \e matrix does need a customized * LocalToGlobal mapping. */ err = PetscObjectChangeTypeName((dObject)vg,VECDOHP);dCHK(err); *v = vg; dFunctionReturn(0); }
/** Get coordinates for every node in closure (every subelement vertex) * * @note This function cannot be implemented for all \a dFS types. For most purposes, users should * \a dFSGetGeometryVectorExpanded and evaluate (element by element) on the nodes of their choice * (with a self-quadrature). * * @param fs Function space * @param inx the new vector with block size 3 and the same number of blocks as the closure vector **/ dErr dFSGetNodalCoordinatesGlobal(dFS fs,Vec *inx) { dErr err; Vec Expanded,Ones,X,Count,Xclosure,Countclosure; dFS fs3; dFunctionBegin; dValidHeader(fs,DM_CLASSID,1); dValidPointer(inx,2); *inx = 0; err = dFSGetNodalCoordinateFS(fs,&fs3); dCHK(err); err = dFSGetNodalCoordinatesExpanded(fs,&Expanded); dCHK(err); if (!fs->nodalcoord.global) { err = dFSCreateGlobalVector(fs3,&fs->nodalcoord.global); dCHK(err); } X = fs->nodalcoord.global; /* Count the number of occurances of each node in the closure. */ err = VecDuplicate(Expanded,&Ones); dCHK(err); err = VecDuplicate(X,&Count); dCHK(err); err = VecDohpZeroEntries(Count); dCHK(err); err = VecSet(Ones,1.); dCHK(err); err = dFSExpandedToGlobal(fs3,Ones,Count,dFS_INHOMOGENEOUS,ADD_VALUES); dCHK(err); err = VecDestroy(&Ones); dCHK(err); err = VecDohpZeroEntries(X); dCHK(err); err = dFSExpandedToGlobal(fs3,Expanded,X,dFS_INHOMOGENEOUS,ADD_VALUES); dCHK(err); err = VecDohpGetClosure(X,&Xclosure); dCHK(err); err = VecDohpGetClosure(Count,&Countclosure); dCHK(err); err = VecPointwiseDivide(Xclosure,Xclosure,Countclosure); dCHK(err); err = VecDohpRestoreClosure(X,&Xclosure); dCHK(err); err = VecDohpRestoreClosure(Count,&Countclosure); dCHK(err); err = VecDestroy(&Count); dCHK(err); *inx = X; dFunctionReturn(0); }
/** Get displacements of every node in expanded vector. * * @param fs Function space * @param ingeom the new geometry vector with block size 3 and same number of blocks as points in the expanded vector. * * @note This is not suitable for use as a function space, * * @note It is important to get geometry associated with boundary sets because it will frequently be projected against * the boundary. **/ dErr dFSGetGeometryVectorExpanded(dFS fs,Vec *ingeom) { dErr err; dFunctionBegin; dValidHeader(fs,DM_CLASSID,1); dValidPointer(ingeom,2); *ingeom = 0; if (!fs->geometry.expanded) { err = dFSCreateGeometryFromMesh_Private(fs); dCHK(err); } *ingeom = fs->geometry.expanded; dFunctionReturn(0); }
dErr VecDohpRestoreClosure(Vec v,Vec *c) { dErr err; dBool isdohp; dFunctionBegin; dValidHeader(v,VEC_CLASSID,1); dValidPointer(c,2); err = PetscTypeCompare((dObject)v,VECDOHP,&isdohp);dCHK(err); if (!isdohp) dERROR(PETSC_COMM_SELF,1,"Vector type %s does not have closure",((dObject)v)->type_name); if (*c != ((Vec_MPI*)v->data)->localrep) dERROR(PETSC_COMM_SELF,1,"attempting to restore incorrect closure"); err = VecStateSync_Private(v,*c);dCHK(err); err = PetscObjectDereference((dObject)*c);dCHK(err); *c = NULL; dFunctionReturn(0); }
/** Get the closed form of a Dohp vector * * @note Dohp vectors are basically just MPI vectors, the only difference is that instead of a local form, we have a * closed form. We subvert .localrep to mean the closed form. * **/ dErr VecDohpGetClosure(Vec v,Vec *c) { Vec_MPI *vmpi; dBool isdohp; dErr err; dFunctionBegin; dValidHeader(v,VEC_CLASSID,1); dValidPointer(c,2); err = PetscTypeCompare((dObject)v,VECDOHP,&isdohp);dCHK(err); if (!isdohp) dERROR(PETSC_COMM_SELF,1,"Vector type %s does not have closure",((dObject)v)->type_name); vmpi = v->data; if (!vmpi->localrep) dERROR(PETSC_COMM_SELF,1,"Vector has no closure"); *c = vmpi->localrep; err = VecStateSync_Private(v,*c);dCHK(err); err = PetscObjectReference((dObject)*c);dCHK(err); dFunctionReturn(0); }
static dErr VecDuplicate_Dohp(Vec x,Vec *iny) { Vec y,xc,yc; Vec_MPI *ympi; dScalar *a; dErr err; dFunctionBegin; dValidHeader(x,VEC_CLASSID,1); dValidPointer(iny,2); *iny = 0; err = VecDohpGetClosure(x,&xc);dCHK(err); err = VecDuplicate(xc,&yc);dCHK(err); err = VecDohpRestoreClosure(x,&xc);dCHK(err); /* The rest is mostly the same as VecDuplicate_MPI, but we can't call that because it allocates memory. * Unfortunately, this is fragile if the VecMPI implementation changes. I think this part of PETSc is quite stable and * I will be sufficiently involved to notice changes here. Famous last words. */ err = VecCreate(((dObject)x)->comm,&y);dCHK(err); err = PetscLayoutReference(x->map,&y->map);dCHK(err); err = VecGetArray(yc,&a);dCHK(err); err = VecCreate_MPI_Private(y,PETSC_FALSE,0,a);dCHK(err); err = VecRestoreArray(yc,&a);dCHK(err); ympi = y->data; err = dMemcpy(y->ops,x->ops,sizeof(struct _VecOps));dCHK(err); ympi->localrep = yc; /* subverting .localrep to mean closed form */ y->stash.donotstash = x->stash.donotstash; y->stash.ignorenegidx = x->stash.ignorenegidx; err = PetscOListDuplicate(((dObject)x)->olist,&((dObject)y)->olist);dCHK(err); err = PetscFListDuplicate(((dObject)x)->qlist,&((dObject)y)->qlist);dCHK(err); y->map->bs = x->map->bs; y->bstash.bs = x->bstash.bs; err = PetscObjectChangeTypeName((dObject)y,VECDOHP);dCHK(err); *iny = y; dFunctionReturn(0); }
/** dFSRedimension - Create a new function space with the same layout, but different number of dofs per node * */ dErr dFSRedimension(dFS fs,dInt bs,dFSClosureMode mode,dFS *infs) { dErr err; dFS rfs; dJacobi jac; dFunctionBegin; dValidHeader(fs,DM_CLASSID,1); if (bs < 1) dERROR(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Block size must be at least 1, was %D",bs); dValidPointer(infs,4); *infs = NULL; if (!fs->spacebuilt) dERROR(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Cannot redimension a space that has not been built"); if (!fs->mesh || !fs->set.active) dERROR(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Space has been built, but does not have a mesh (should not happen)"); err = dFSCreate(((dObject)fs)->comm,&rfs); dCHK(err); err = dFSSetBlockSize(rfs,bs); dCHK(err); err = dFSSetMesh(rfs,fs->mesh,fs->set.active); dCHK(err); err = dFSGetJacobi(fs,&jac); dCHK(err); err = dFSSetDegree(rfs,jac,fs->tag.degree); dCHK(err); switch (mode) { case dFS_CLOSURE: dERROR(PETSC_COMM_SELF,PETSC_ERR_SUP,"Probably not what you want"); /* because ordering would be different, but there is nothing to do for this choice. */ break; case dFS_INTERIOR: { dMeshESH *bsets; dInt nsets; err = dMeshGetNumSubsets(fs->mesh,fs->set.boundaries,0,&nsets); dCHK(err); err = dMallocA(nsets,&bsets); dCHK(err); err = dMeshGetSubsets(fs->mesh,fs->set.boundaries,0,bsets,nsets,NULL); dCHK(err); for (dInt i=0; i<nsets; i++) { dFSBStatus bstat; err = dMeshTagSGetData(fs->mesh,fs->tag.bstatus,&bsets[i],1,&bstat,sizeof(bstat),dDATA_BYTE); dCHK(err); err = dFSRegisterBoundarySet(rfs,bsets[i],bstat,NULL,NULL); dCHK(err); /* TODO: handle dFSConstraintCtx */ } err = dFree(bsets); dCHK(err); } break; default: dERROR(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for dFSClosureMode %s",dFSClosureModes[mode]); } err = dFSSetType(rfs,((dObject)fs)->type_name); dCHK(err); err = dFSSetOrderingType(rfs,fs->orderingtype); dCHK(err); rfs->set.ordered = fs->set.ordered; { // The FS has the layout, ordering, and boundary status tags set so we are ready to build the function space. dMesh mesh; dMeshAdjacency meshadj; err = dFSGetMesh(rfs,&mesh); dCHK(err); err = dMeshGetAdjacency(mesh,rfs->set.ordered,&meshadj); dCHK(err); err = dFSPopulatePartitionedSets_Private(rfs,meshadj); dCHK(err); err = dFSBuildSpaceWithOrderedSet_Private(rfs,meshadj); dCHK(err); err = dMeshRestoreAdjacency(mesh,rfs->set.ordered,&meshadj); dCHK(err); } *infs = rfs; dFunctionReturn(0); }