/** 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 dUnitsView(dUnits un,dViewer viewer) { dBool iascii; dErr err; dFunctionBegin; dValidHeader(un,dUNITS_CLASSID,1); if (!viewer) {err = PetscViewerASCIIGetStdout(((dObject)un)->comm,&viewer);dCHK(err);} dValidHeader(viewer,PETSC_VIEWER_CLASSID,2); PetscCheckSameComm(un,1,viewer,2); err = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);dCHK(err); if (iascii) { err = PetscObjectPrintClassNamePrefixType((PetscObject)un,viewer,"Units Manager");dCHK(err); err = PetscViewerASCIIPushTab(viewer);dCHK(err); for (dInt i=0; i<un->nalloc && un->list[i]; i++) { dUnit u = un->list[i]; err = PetscViewerASCIIPrintf(viewer,"%-12s: 1 internal unit = %10.4e %s (%s) = %10.4e %s\n", dUnitQuantityName(u),dUnitDimensionalize(u,1.0),dUnitName(u),dUnitShortName(u),dUnitDimensionalizeSI(u,1.0),dUnitSIName(u));dCHK(err); err = PetscViewerASCIIPrintf(viewer,"%-12s 1 %s = %10.4e %s\n","",dUnitShortName(u),dUnitDimensionalizeSI(u,dUnitNonDimensionalize(u,1.0)),dUnitSIName(u));dCHK(err); } err = PetscViewerASCIIPopTab(viewer);dCHK(err); } else dERROR(((dObject)un)->comm,PETSC_ERR_SUP,"Viewer type %s not supported",((PetscObject)viewer)->type_name); dFunctionReturn(0); }
static dErr VecStateSync_Private(Vec x,Vec y) { dInt xstate,ystate; dErr err; dFunctionBegin; dValidHeader(x,VEC_CLASSID,1); dValidHeader(y,VEC_CLASSID,2); err = PetscObjectStateQuery((dObject)x,&xstate);dCHK(err); err = PetscObjectStateQuery((dObject)y,&ystate);dCHK(err); err = PetscObjectSetState((dObject)x,dMaxInt(xstate,ystate));dCHK(err); err = PetscObjectSetState((dObject)y,dMaxInt(xstate,ystate));dCHK(err); dFunctionReturn(0); }
dErr dFSRotationView(dFSRotation rot,PetscViewer viewer) { dErr err; dFunctionBegin; dValidHeader(rot,dFSROT_CLASSID,1); if (!viewer) { err = PetscViewerASCIIGetStdout(((dObject)rot)->comm,&viewer);dCHK(err); } dValidHeader(viewer,PETSC_VIEWER_CLASSID,2); PetscCheckSameComm(rot,1,viewer,2); dERROR(PETSC_COMM_SELF,1,"not implemented"); dFunctionReturn(0); }
/** 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); }
dErr dUnitsSetFromOptions(dUnits un) { dErr err; dFunctionBegin; dValidHeader(un,dUNITS_CLASSID,1); err = PetscOptionsBegin(((PetscObject)un)->comm,((PetscObject)un)->prefix,"Units manager","dUnits");dCHK(err); for (dUnitsBaseType btype = 0; btype < dUNITS_MAX; btype++) { char opt[256],help[256],uspec[256]; dReal commonpersi = 1.0,scale = 1.0; dBool flg; err = PetscSNPrintf(opt,sizeof opt,"-units_%s",dUnitsBaseTypes[btype]);dCHK(err); err = PetscSNPrintf(uspec,sizeof uspec,"%s:%s:%f:%f",dUnitsBaseNamesSI[btype],dUnitsBaseNamesShortSI[btype],commonpersi,scale); err = PetscSNPrintf(help,sizeof help,"Common name:short name:one common unit of %s expressed in %s:common units per non-dimensionalized",dUnitsBaseTypes[btype],dUnitsBaseNamesSI[btype]);dCHK(err); err = PetscOptionsString(opt,help,"dUnitsSetBase",uspec,uspec,sizeof uspec,&flg);dCHK(err); if (flg) { char *longname,*shortname,*buf1,*buf2; longname = uspec; if (!(shortname = strchr(longname,':'))) dERROR(((dObject)un)->comm,PETSC_ERR_USER,"The field specification '%s' is ':' delimited",opt); *shortname++ = 0; if (!(buf1 = strchr(shortname,':'))) dERROR(((dObject)un)->comm,PETSC_ERR_USER,"The field specification for '%s' needs four arguments, but only two given",longname); *buf1++ = 0; if (!(buf2 = strchr(buf1,':'))) dERROR(((dObject)un)->comm,PETSC_ERR_USER,"The field specification for '%s' needs four arguments, but only three given",longname); *buf2++ = 0; if (sscanf(buf1,"%lf",&commonpersi) != 1) dERROR(((dObject)un)->comm,PETSC_ERR_USER,"Size of common unit '%s' could not be parsed from '%s'",longname,buf1); if (sscanf(buf2,"%lf",&scale) != 1) dERROR(((dObject)un)->comm,PETSC_ERR_USER,"Scale for common unit '%s' could not be parsed from '%s'",longname,buf2); err = dUnitsSetBase(un,btype,longname,shortname,commonpersi,scale,NULL);dCHK(err); } } err = PetscOptionsEnd();dCHK(err); 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); }
/** Apply rotation to global vector * * @note does nothing if rotation is NULL **/ dErr dFSRotationApply(dFSRotation rot,Vec g,dFSRotateMode rmode,dFSHomogeneousMode hmode) { dErr err; Vec gc,lf; dFunctionBegin; if (!rot) dFunctionReturn(0); dValidHeader(rot,dFSROT_CLASSID,1); dValidHeader(g,VEC_CLASSID,2); if (rot->ops->apply) { err = rot->ops->apply(rot,g,rmode,hmode);dCHK(err); } else { err = VecDohpGetClosure(g,&gc);dCHK(err); err = VecGhostGetLocalForm(gc,&lf);dCHK(err); err = dFSRotationApplyLocal(rot,lf,rmode,hmode);dCHK(err); /* \todo only rotate the global portion */ err = VecGhostRestoreLocalForm(gc,&lf);dCHK(err); err = VecDohpRestoreClosure(g,&gc);dCHK(err); } dFunctionReturn(0); }
dErr dFSRotationDestroy(dFSRotation *rot) { dErr err; dFunctionBegin; if (!*rot) dFunctionReturn(0); dValidHeader(*rot,dFSROT_CLASSID,1); err = ISDestroy(&(*rot)->is);dCHK(err); err = VecDestroy(&(*rot)->strong);dCHK(err); err = dFree2((*rot)->rmat,(*rot)->nstrong);dCHK(err); err = PetscHeaderDestroy(rot);dCHK(err); dFunctionReturn(0); }
/** Get the number of subelements and number of vertices of subelements. **/ dErr dFSGetSubElementMesh(dFS fs,dInt nelems,dInt nconn,dEntTopology topo[],dInt off[],dInt ind[]) { dErr err; dFunctionBegin; dValidHeader(fs,DM_CLASSID,1); dValidIntPointer(topo,4); dValidIntPointer(off,5); dValidIntPointer(ind,6); if (!fs->ops->getsubelementmesh) dERROR(PETSC_COMM_SELF,1,"not implemented"); err = (*fs->ops->getsubelementmesh)(fs,nelems,nconn,topo,off,ind); dCHK(err); dFunctionReturn(0); }
/** Get the number of subelements and number of vertices of subelements. * * @note the number of vertices is the same as the number of local nodes in closure vertor. **/ dErr dFSGetSubElementMeshSize(dFS fs,dInt *nelems,dInt *nverts,dInt *nconn) { dErr err; dFunctionBegin; dValidHeader(fs,DM_CLASSID,1); dValidIntPointer(nelems,2); dValidIntPointer(nverts,3); dValidIntPointer(nconn,4); if (!fs->ops->getsubelementmeshsize) dERROR(PETSC_COMM_SELF,1,"not implemented"); err = (*fs->ops->getsubelementmeshsize)(fs,nelems,nverts,nconn); dCHK(err); 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 VecDohpZeroEntries(Vec v) { dErr err; dBool isdohp; Vec c; dFunctionBegin; dValidHeader(v,VEC_CLASSID,1); err = PetscTypeCompare((dObject)v,VECDOHP,&isdohp);dCHK(err); if (!isdohp) dERROR(PETSC_COMM_SELF,PETSC_ERR_SUP,"Vector type %s",((dObject)v)->type_name); err = VecDohpGetClosure(v,&c);dCHK(err); err = VecZeroEntries(c);dCHK(err); err = VecDohpRestoreClosure(v,&c);dCHK(err); 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); }
// Generates a mesh of a brick using run-time parameters. // The new mesh populates the given root set. // This should be converted to have a useful programmatic API. dErr dMeshGenerateBlock(dMesh dmesh,dMeshESH root,PetscBool do_geom) { const char pTagName[]="OWNING_PART", pSetName[]="PARALLEL_PARTITION"; PetscBool assoc_with_brick=0,do_color_bdy=0,do_material = 1,do_uniform = 1,do_global_number = 0,do_global_id = 1; PetscBool do_partition = 1,do_pressure = 0,do_faces = 1,do_edges = 1; dReal rotate_y = 0; dInt verbose = 0; iMesh_Instance mesh; iBase_EntityHandle *entbuf; iBase_EntitySetHandle facesets[6]; iBase_TagHandle pTag; MeshListEH v=MLZ,e=MLZ,f=MLZ,r=MLZ,c=MLZ; MeshListReal x=MLZ; MeshListInt s=MLZ,part=MLZ; dInt *face[6],facecount[6]={0}; int err,i,j,k,m,n,p,M,N,P,I,J,K,order=iBase_INTERLEAVED; Box box; PetscViewer viewer; dFunctionBegin; dValidHeader(dmesh,dMESH_CLASSID,1); err = PetscViewerASCIIGetStdout(((PetscObject)dmesh)->comm,&viewer);dCHK(err); err = PetscOptionsBegin(((PetscObject)dmesh)->comm,((PetscObject)dmesh)->prefix,"dMeshGenerate Block: generate cartesian meshes",NULL);dCHK(err); { char boxstr[256] = "-1:1,-1:1,-1:1",mnp[256] = "5,5,5",MNP[256] = "2,2,2"; err = PetscOptionsInt("-dmeshgen_block_verbose","verbosity of output","none",verbose,&verbose,NULL);dCHK(err); if (do_geom) { err = PetscOptionsBool("-dmeshgen_block_assoc_with_brick","associate boundaries with brick","none",assoc_with_brick,&assoc_with_brick,NULL);dCHK(err); } err = PetscOptionsBool("-dmeshgen_block_color_bdy","color boundary sets","none",do_color_bdy,&do_color_bdy,NULL);dCHK(err); err = PetscOptionsBool("-dmeshgen_block_material","create material sets","none",do_material,&do_material,NULL);dCHK(err); err = PetscOptionsBool("-dmeshgen_block_uniform","create uniform sets","none",do_uniform,&do_uniform,NULL);dCHK(err); err = PetscOptionsBool("-dmeshgen_block_global_number","create global_number tags","none",do_global_number,&do_global_number,NULL);dCHK(err); err = PetscOptionsBool("-dmeshgen_block_global_id","create GLOBAL_ID tags","none",do_global_id,&do_global_id,NULL);dCHK(err); err = PetscOptionsBool("-dmeshgen_block_partition","create partition sets","none",do_partition,&do_partition,NULL);dCHK(err); err = PetscOptionsBool("-dmeshgen_block_pressure","create pressure sets","none",do_pressure,&do_pressure,NULL);dCHK(err); err = PetscOptionsBool("-dmeshgen_block_faces","create face entities","none",do_faces,&do_faces,NULL);dCHK(err); err = PetscOptionsBool("-dmeshgen_block_edges","create face entities","none",do_edges,&do_edges,NULL);dCHK(err); err = PetscOptionsReal("-dmeshgen_block_rotate_y","rotate domain by given angle (degrees) around y axis)","none",rotate_y,&rotate_y,NULL);dCHK(err); rotate_y *= PETSC_PI/180.; err = PetscOptionsString("-dmeshgen_block_box","box x0:x1,y0:y1,z0:z1","none",boxstr,boxstr,sizeof(boxstr),NULL);dCHK(err); err = PetscOptionsString("-dmeshgen_block_mnp","number of points m,n,p","none",mnp,mnp,sizeof(mnp),NULL);dCHK(err); err = PetscOptionsString("-dmeshgen_block_procs_mnp","number of procs M,N,P","none",MNP,MNP,sizeof(MNP),NULL);dCHK(err); i = sscanf(boxstr,"%lf:%lf,%lf:%lf,%lf:%lf",&box.x0,&box.x1,&box.y0,&box.y1,&box.z0,&box.z1); if (i != 6) dERROR(PETSC_COMM_SELF,1,"Failed to parse bounding box."); i = sscanf(mnp,"%d,%d,%d",&m,&n,&p); if (i != 3) dERROR(PETSC_COMM_SELF,1,"Failed to parse size."); i = sscanf(MNP,"%d,%d,%d",&M,&N,&P); if (i != 3) dERROR(PETSC_COMM_SELF,1,"Failed to parse partition size."); } err = PetscOptionsEnd(); err = dMeshGetInstance(dmesh,&mesh);dCHK(err); /* Allocate buffers */ err = dMallocA(m*n*p*3,&entbuf);dCHK(err); /* More than enough to hold all entities of any given type */ for (i=0; i<6; i++) { int n2max = dSqrInt(dMaxInt(m,dMaxInt(n,p))); err = dMallocA(2*n2max,&face[i]);dCHK(err); iMesh_createEntSet(mesh,0,&facesets[i],&err);dICHK(mesh,err); } /* Create vertices */ x.a = x.s = m*n*p*3; x.v = malloc(x.a*sizeof(double)); for (i=0; i<m; i++) { for (j=0; j<n; j++) { for (k=0; k<p; k++) { dReal X,Y,Z; I = (i*n+j)*p+k; if (i==0) AddToFace(face,facecount,3,I); else if (i==m-1) AddToFace(face,facecount,1,I); else if (j==0) AddToFace(face,facecount,0,I); else if (j==n-1) AddToFace(face,facecount,2,I); else if (k==0) AddToFace(face,facecount,4,I); else if (k==p-1) AddToFace(face,facecount,5,I); X = box.x0 + (box.x1-box.x0)*(1.*i/(m-1)); Y = box.y0 + (box.y1-box.y0)*(1.*j/(n-1)); Z = box.z0 + (box.z1-box.z0)*(1.*k/(p-1)); x.v[3*I+0] = cos(rotate_y) * X - sin(rotate_y) * Z; x.v[3*I+1] = Y; x.v[3*I+2] = sin(rotate_y) * X + cos(rotate_y) * Z; } } } iMesh_createVtxArr(mesh,m*n*p,order,x.v,x.s,&v.v,&v.a,&v.s,&err);dICHK(mesh,err); err = CommitToFaceSets(mesh,v.v,face,facecount,facesets,entbuf); MeshListFree(x); /* Create regions */ c.a = c.s = (m-1)*(n-1)*(p-1)*8; c.v = malloc(c.a*sizeof(iBase_EntityHandle)); /* connectivity */ I=0; for (i=0; i<m-1; i++) { for (j=0; j<n-1; j++) { for (k=0; k<p-1; k++) { c.v[I++] = v.v[((i+0)*n+(j+0))*p+(k+0)]; c.v[I++] = v.v[((i+1)*n+(j+0))*p+(k+0)]; c.v[I++] = v.v[((i+1)*n+(j+1))*p+(k+0)]; c.v[I++] = v.v[((i+0)*n+(j+1))*p+(k+0)]; c.v[I++] = v.v[((i+0)*n+(j+0))*p+(k+1)]; c.v[I++] = v.v[((i+1)*n+(j+0))*p+(k+1)]; c.v[I++] = v.v[((i+1)*n+(j+1))*p+(k+1)]; c.v[I++] = v.v[((i+0)*n+(j+1))*p+(k+1)]; } } } if (I != c.s) dERROR(PETSC_COMM_SELF,1,"Wrong number of regions."); iMesh_createEntArr(mesh,iMesh_HEXAHEDRON,c.v,c.s,&r.v,&r.a,&r.s,&s.v,&s.a,&s.s,&err);dICHK(mesh,err); if (r.s != (m-1)*(n-1)*(p-1)) dERROR(PETSC_COMM_SELF,1,"Wrong number of regions created."); if (verbose > 0) {err = PetscViewerASCIIPrintf(viewer,"region size %d, status size %d\n",r.s,s.s);dCHK(err);} if (do_global_number) {err = doGlobalNumber(mesh,root);dCHK(err);} if (do_global_id) {err = doGlobalID(mesh,root);dCHK(err);} if (do_partition) { /* Partition tags */ /* Create partition. */ part.a = part.s = r.s; part.v = malloc(part.a*sizeof(int)); for (i=0; i<m-1; i++) { for (j=0; j<n-1; j++) { for (k=0; k<p-1; k++) { I = i*M/(m-1); J = j*N/(n-1); K = k*P/(p-1); part.v[(i*(n-1)+j)*(p-1)+k] = (I*N+J)*P+K; } } } /* MATERIAL_SET is a special name associated with all iMesh instances * If we are using a different name, we can assume it is not special. */ if (strcmp(pTagName,"MATERIAL_SET")) { iMesh_createTag(mesh,pTagName,1,iBase_INTEGER,&pTag,&err,sizeof(pTagName));dICHK(mesh,err); } else { iMesh_getTagHandle(mesh,"MATERIAL_SET",&pTag,&err,sizeof("MATERIAL_SET"));dICHK(mesh,err); } iMesh_setIntArrData(mesh,r.v,r.s,pTag,part.v,part.s,&err);dICHK(mesh,err); MeshListFree(part); } if (do_partition) /* Partition sets */ { int ii,jj,kk; iBase_EntitySetHandle partset; iBase_EntityHandle *entp; /* reuse some stuff, set up the a partition set */ iMesh_createTag(mesh,pSetName,1,iBase_INTEGER,&pTag,&err,sizeof(pSetName));dICHK(mesh,err); for (i=0; i<M; i++) { for (j=0; j<N; j++) { for (k=0; k<P; k++) { iMesh_createEntSet(mesh,0,&partset,&err);dICHK(mesh,err); entp = entbuf; for (ii=i*(m-1)/M; ii<(i+1)*(m-1)/M; ii++) { for (jj=j*(n-1)/N; jj<(j+1)*(n-1)/N; jj++) { for (kk=k*(p-1)/P; kk<(k+1)*(p-1)/P; kk++) { *entp++ = r.v[(ii*(n-1)+jj)*(p-1)+kk]; } } } if (verbose > 0) {err = PetscViewerASCIIPrintf(viewer,"part[%d (%d,%d,%d)] has %d regions\n",(i*N+j)*P+k,i,j,k,(int)(entp-entbuf));dCHK(err);} iMesh_addEntArrToSet(mesh,entbuf,(int)(entp-entbuf),partset,&err);dICHK(mesh,err); iMesh_setEntSetIntData(mesh,partset,pTag,(i*N+j)*P+k,&err);dICHK(mesh,err); } } } } MeshListFree(r); MeshListFree(s); MeshListFree(c); if (do_faces) { /* Create faces */ c.a = c.s = 4*((m-1)*(n-1)*p + (m-1)*n*(p-1) + m*(n-1)*(p-1)); c.v = malloc(c.a*sizeof(iBase_EntityHandle)); I = 0; for (i=0; i<m-1; i++) { /* Faces with normal pointing in positive z direction */ for (j=0; j<n-1; j++) { for (k=0; k<p; k++) { if (k==0) AddToFace(face,facecount,4,I/4); if (k==p-1) AddToFace(face,facecount,5,I/4); c.v[I++] = v.v[((i+0)*n+(j+0))*p+k]; c.v[I++] = v.v[((i+1)*n+(j+0))*p+k]; c.v[I++] = v.v[((i+1)*n+(j+1))*p+k]; c.v[I++] = v.v[((i+0)*n+(j+1))*p+k]; } } } for (i=0; i<m-1; i++) { /* Faces with normal pointing in negative y direction */ for (j=0; j<n; j++) { for (k=0; k<p-1; k++) { if (j==0) AddToFace(face,facecount,0,I/4); if (j==n-1) AddToFace(face,facecount,2,I/4); c.v[I++] = v.v[((i+0)*n+j)*p+(k+0)]; c.v[I++] = v.v[((i+1)*n+j)*p+(k+0)]; c.v[I++] = v.v[((i+1)*n+j)*p+(k+1)]; c.v[I++] = v.v[((i+0)*n+j)*p+(k+1)]; } } } for (i=0; i<m; i++) { /* Faces with normal pointing in positive x direction */ for (j=0; j<n-1; j++) { for (k=0; k<p-1; k++) { if (i==0) AddToFace(face,facecount,3,I/4); if (i==m-1) AddToFace(face,facecount,1,I/4); c.v[I++] = v.v[(i*n+(j+0))*p+(k+0)]; c.v[I++] = v.v[(i*n+(j+1))*p+(k+0)]; c.v[I++] = v.v[(i*n+(j+1))*p+(k+1)]; c.v[I++] = v.v[(i*n+(j+0))*p+(k+1)]; } } } if (I != c.s) dERROR(PETSC_COMM_SELF,1, "Wrong number of faces."); iMesh_createEntArr(mesh,iMesh_QUADRILATERAL,c.v,c.s,&f.v,&f.a,&f.s,&s.v,&s.a,&s.s,&err);dICHK(mesh,err); err = CommitToFaceSets(mesh,f.v,face,facecount,facesets,entbuf);dCHK(err); if (verbose > 0) {err = PetscViewerASCIIPrintf(viewer,"face size %d, status size %d\n",f.s,s.s);dCHK(err);} MeshListFree(f); MeshListFree(s); MeshListFree(c); } if (do_edges) { /* Create edges */ c.a = c.s = 2*(m*n*(p-1) + m*(n-1)*p + (m-1)*n*p); c.v = malloc(c.a*sizeof(iBase_EntityHandle)); I = 0; for (i=0; i<m; i++) { for (j=0; j<n; j++) { for (k=0; k<p-1; k++) { if (i==0) AddToFace(face,facecount,0,I/2); else if (i==m-1) AddToFace(face,facecount,2,I/2); else if (j==0) AddToFace(face,facecount,3,I/2); else if (j==n-1) AddToFace(face,facecount,1,I/2); c.v[I++] = v.v[(i*n+j)*p+(k+0)]; c.v[I++] = v.v[(i*n+j)*p+(k+1)]; } } } for (i=0; i<m; i++) { for (j=0; j<n-1; j++) { for (k=0; k<p; k++) { if (i==0) AddToFace(face,facecount,0,I/2); else if (i==m-1) AddToFace(face,facecount,2,I/2); else if (k==0) AddToFace(face,facecount,4,I/2); else if (k==p-1) AddToFace(face,facecount,5,I/2); c.v[I++] = v.v[(i*n+(j+0))*p+k]; c.v[I++] = v.v[(i*n+(j+1))*p+k]; } } } for (i=0; i<m-1; i++) { for (j=0; j<n; j++) { for (k=0; k<p; k++) { if (j==0) AddToFace(face,facecount,3,I/2); else if (j==n-1) AddToFace(face,facecount,1,I/2); else if (k==0) AddToFace(face,facecount,4,I/2); else if (k==p-1) AddToFace(face,facecount,5,I/2); c.v[I++] = v.v[((i+0)*n+j)*p+k]; c.v[I++] = v.v[((i+1)*n+j)*p+k]; } } } if (I != c.s) dERROR(PETSC_COMM_SELF,1, "Wrong number of edges."); iMesh_createEntArr(mesh,iMesh_LINE_SEGMENT,c.v,c.s, &e.v,&e.a,&e.s, &s.v,&s.a,&s.s,&err);dICHK(mesh,err); err = CommitToFaceSets(mesh,e.v,face,facecount,facesets,entbuf);dCHK(err); if (verbose > 0) {err = PetscViewerASCIIPrintf(viewer,"edge size %d, status size %d\n",e.s,s.s);dCHK(err);} MeshListFree(e); MeshListFree(s); MeshListFree(c); } /* We are done with the master vertex record. */ MeshListFree(v); /* Create boundary sets, these are not related to geometry here */ { dMeshESH wallset,topset,bottomset,senseSet; iBase_TagHandle bdyTag,senseTag; iMesh_getTagHandle(mesh,"NEUMANN_SET",&bdyTag,&err,sizeof("NEUMANN_SET"));dICHK(mesh,err); iMesh_createTag(mesh,"SENSE",1,iBase_INTEGER,&senseTag,&err,sizeof "SENSE");dICHK(mesh,err); iMesh_createEntSet(mesh,0,&wallset,&err);dICHK(mesh,err); iMesh_createEntSet(mesh,0,&topset,&err);dICHK(mesh,err); iMesh_createEntSet(mesh,0,&bottomset,&err);dICHK(mesh,err); iMesh_setEntSetIntData(mesh,wallset,bdyTag,100,&err);dICHK(mesh,err); iMesh_setEntSetIntData(mesh,topset,bdyTag,200,&err);dICHK(mesh,err); iMesh_setEntSetIntData(mesh,bottomset,bdyTag,300,&err);dICHK(mesh,err); for (i=0; i<4; i++) {iMesh_addEntSet(mesh,facesets[i],wallset,&err);dICHK(mesh,err);} iMesh_addEntSet(mesh,facesets[5],topset,&err);dICHK(mesh,err); iMesh_addEntSet(mesh,facesets[4],bottomset,&err);dICHK(mesh,err); /* Deal with SENSE on the walls */ iMesh_createEntSet(mesh,0,&senseSet,&err);dICHK(mesh,err); iMesh_addEntSet(mesh,facesets[2],senseSet,&err);dICHK(mesh,err); iMesh_addEntSet(mesh,facesets[3],senseSet,&err);dICHK(mesh,err); iMesh_setEntSetIntData(mesh,senseSet,senseTag,-1,&err);dICHK(mesh,err); iMesh_addEntSet(mesh,senseSet,wallset,&err);dICHK(mesh,err); /* Deal with SENSE on the bottom */ iMesh_createEntSet(mesh,0,&senseSet,&err);dICHK(mesh,err); iMesh_addEntSet(mesh,facesets[4],senseSet,&err);dICHK(mesh,err); iMesh_setEntSetIntData(mesh,senseSet,senseTag,-1,&err);dICHK(mesh,err); iMesh_addEntSet(mesh,senseSet,bottomset,&err);dICHK(mesh,err); for (i=0; i<6; i++) {err = dFree(face[i]);} err = dFree(entbuf);dCHK(err); } if (do_material) {err = doMaterial(mesh,root);dCHK(err);} /* Add a real valued tag over the vertices. */ if (do_pressure) { static const char *myTagName = "pressure"; iBase_TagHandle myTag; double *myData; iMesh_getEntities(mesh,root,iBase_VERTEX,iMesh_POINT,&v.v,&v.a,&v.s,&err);dICHK(mesh,err); iMesh_createTag(mesh,myTagName,1,iBase_DOUBLE,&myTag,&err,(int)strlen(myTagName));dICHK(mesh,err); err = PetscMalloc(v.s*sizeof(double),&myData);dCHK(err); for (i=0; i<v.s; i++) { myData[i] = 1.0 * i; } iMesh_setDblArrData(mesh,v.v,v.s,myTag,myData,v.s,&err);dICHK(mesh,err); err = PetscFree(myData);dCHK(err); MeshListFree(v); } if (do_uniform) {err = createUniformTags(mesh,root);dCHK(err);} if (do_geom) #ifndef dHAVE_ITAPS_REL dERROR(((dObject)dmesh)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Dohp has not been configured with support for geometry"); #else { const char geom_options[] = ";ENGINE=OCC;"; const char rel_options[] = ""; iGeom_Instance geom; iRel_Instance assoc; iRel_PairHandle pair; iBase_EntityHandle brick; iGeom_newGeom(geom_options,&geom,&err,sizeof geom_options);dIGCHK(geom,err); iRel_create(rel_options,&assoc,&err,sizeof rel_options);dIRCHK(assoc,err); iRel_createPair(assoc,geom,0,iRel_IGEOM_IFACE,iRel_ACTIVE,mesh,1,iRel_IMESH_IFACE,iRel_ACTIVE,&pair,&err);dIGCHK(assoc,err); iGeom_createBrick(geom,box.x1-box.x0,box.y1-box.y0,box.z1-box.z0,&brick,&err);dIGCHK(geom,err); iGeom_moveEnt(geom,brick,0.5*(box.x0+box.x1),0.5*(box.y0+box.y1),0.5*(box.z0+box.z1),&err);dIGCHK(geom,err); if (verbose > 0) {err = BoundingBoxView(geom,brick,"brick",viewer);dCHK(err);} { iBase_EntityHandle gface[6],*gface_p=gface; int gface_a=6,gface_s; iGeom_getEntAdj(geom,brick,2,&gface_p,&gface_a,&gface_s,&err);dIGCHK(geom,err); for (i=0; i<6; i++) { char name[20]; sprintf(name,"face_%d",i); err = BoundingBoxView(geom,gface[i],name,viewer);dCHK(err); } if (assoc_with_brick) { for (i=0; i<6; i++) { iRel_setEntSetRelation(assoc,pair,brick,facesets[i],&err);dIRCHK(assoc,err); } } else { /* Set associations. With the current Lasso implementation, these will not be saved */ iRel_setEntSetRelation(assoc,pair,gface[0],facesets[3],&err);dIRCHK(assoc,err); iRel_setEntSetRelation(assoc,pair,gface[1],facesets[1],&err);dIRCHK(assoc,err); iRel_setEntSetRelation(assoc,pair,gface[2],facesets[0],&err);dIRCHK(assoc,err); iRel_setEntSetRelation(assoc,pair,gface[3],facesets[2],&err);dIRCHK(assoc,err); iRel_setEntSetRelation(assoc,pair,gface[4],facesets[4],&err);dIRCHK(assoc,err); iRel_setEntSetRelation(assoc,pair,gface[5],facesets[5],&err);dIRCHK(assoc,err); } } { dMeshTag meshGlobalIDTag,meshGeomDimTag,geomGlobalIDTag; /* Manually set association tags, these are set so that the associations above can be inferred. */ iMesh_getTagHandle(mesh,"GLOBAL_ID",&meshGlobalIDTag,&err,sizeof "GLOBAL_ID");dICHK(mesh,err); iMesh_getTagHandle(mesh,"GEOM_DIMENSION",&meshGeomDimTag,&err,sizeof "GEOM_DIMENSION");dICHK(mesh,err); iGeom_getTagHandle(geom,"GLOBAL_ID",&geomGlobalIDTag,&err,sizeof "GLOBAL_ID");dIGCHK(geom,err); for (i=0; i<6; i++) { iBase_EntityHandle gface; int gid,gdim; iRel_getSetEntRelation(assoc,pair,facesets[i],1,&gface,&err);dIRCHK(assoc,err); iGeom_getEntType(geom,gface,&gdim,&err);dIGCHK(geom,err); if (gdim != 2) dERROR(PETSC_COMM_SELF,1,"Geometric dimension is %d, expected 2",gdim); iGeom_getIntData(geom,gface,geomGlobalIDTag,&gid,&err);dIGCHK(geom,err); iMesh_setEntSetIntData(mesh,facesets[i],meshGeomDimTag,2,&err);dICHK(mesh,err); /* If the following line is disabled, Lasso will pick up the wrong relations, but at least they will still be with * surfaces. Wouldn't it be better to not find relations? */ iMesh_setEntSetIntData(mesh,facesets[i],meshGlobalIDTag,gid,&err);dICHK(mesh,err); } } err = dMeshSetGeometryRelation(dmesh,geom,assoc);dCHK(err); } #endif 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); }
/** Apply rotation to local vector * * @note does nothing if rotation is NULL **/ dErr dFSRotationApplyLocal(dFSRotation rot,Vec l,dFSRotateMode rmode,dFSHomogeneousMode hmode) { dErr err; dFunctionBegin; if (!rot) dFunctionReturn(0); dValidHeader(rot,dFSROT_CLASSID,1); dValidHeader(l,VEC_CLASSID,2); if (rot->ops->applylocal) { err = rot->ops->applylocal(rot,l,rmode,hmode);dCHK(err); } else { const dInt *idx,bs = rot->bs; const dScalar *s; dScalar *v,*strong; dScalar tmp[16]; if (bs > 16) dERROR(PETSC_COMM_SELF,1,"large block size"); err = ISGetIndices(rot->is,&idx);dCHK(err); err = VecGetArray(l,&v);dCHK(err); if (hmode == dFS_INHOMOGENEOUS) { err = VecGetArray(rot->strong,&strong);dCHK(err); s = strong; } else { s = 0; } for (dInt i=0; i<rot->n; i++) { const dInt ii = idx[i]; const dReal *rmat = &rot->rmat[ii*bs*bs]; for (dInt j=0; j<bs; j++) { tmp[j] = v[ii*bs+j]; v[ii*bs+j] = 0; } switch (rmode) { case dFS_ROTATE_FORWARD: for (dInt j=0; j<bs; j++) { for (dInt k=0; k<bs; k++) { v[ii*bs+k] += rmat[k*bs+j] * tmp[j]; } } switch (hmode) { case dFS_HOMOGENEOUS: for (dInt j=0; j<s[i]; j++) v[ii*bs+j] = 0; break; case dFS_INHOMOGENEOUS: for (dInt j=0; j<s[i]; j++) v[ii*bs+j] = *s++; break; default: dERROR(PETSC_COMM_SELF,1,"Invalid homogeneous mode"); } break; case dFS_ROTATE_REVERSE: switch (hmode) { case dFS_HOMOGENEOUS: for (dInt j=0; j<s[i]; j++) tmp[j] = 0; break; case dFS_INHOMOGENEOUS: for (dInt j=0; j<s[i]; j++) tmp[j] = *s++; break; default: dERROR(PETSC_COMM_SELF,1,"Invalid homogeneous mode"); } for (dInt j=0; j<bs; j++) { for (dInt k=0; k<bs; k++) { v[ii*bs+k] += rmat[j*bs+k] * tmp[k]; } } break; default: dERROR(PETSC_COMM_SELF,1,"Invalid rotate mode"); } } err = ISRestoreIndices(rot->is,&idx);dCHK(err); err = VecRestoreArray(l,&v);dCHK(err); if (hmode == dFS_HOMOGENEOUS) { dInt n; err = VecGetSize(rot->strong,&n);dCHK(err); if (s-strong != n) dERROR(PETSC_COMM_SELF,1,"should not happen"); err = VecRestoreArray(rot->strong,&strong);dCHK(err); } } dFunctionReturn(0); }