Actual source code: stagintern.c

petsc-3.13.6 2020-09-29
Report Typos and Errors
  1: /* DMStag dimension-independent internal functions. If added to the public API,
  2:    these would move to stagutils.c */

  4:  #include <petsc/private/dmstagimpl.h>

  6: /* Note: this is an internal function but we provide a man page in case it's made public */
  7: /*@C
  8:   DMStagDuplicateWithoutSetup - duplicate a DMStag object without setting it up

 10:   Collective

 12:   Input Parameters:
 13: + dm - The original DM object
 14: - comm - the MPI communicator for the new DM (MPI_COMM_NULL to use the same communicator as dm)

 16:   Output Parameter:
 17: . newdm  - The new DM object

 19:   Developer Notes:
 20:   Copies over all of the state for a DMStag object, except that which is
 21:   populated during DMSetUp().  This function is used within (all) other
 22:   functions that require an un-setup clone, which is common when duplicating,
 23:   coarsening, refining, or creating compatible DMs with different fields.  For
 24:   this reason it also accepts an MPI communicator as an argument (though note
 25:   that at the time of this writing, implementations of DMCoarsen and DMRefine
 26:   don't usually seem to respect their "comm" arguments). This function could be
 27:   pushed up to the general DM API (and perhaps given a different name).

 29:   Level: developer

 31:   .seealso: DMClone(), DMStagCreateCompatibleDMStag(), DMCoarsen(), DMRefine()
 32: @*/
 33: PetscErrorCode DMStagDuplicateWithoutSetup(DM dm, MPI_Comm comm, DM *newdm)
 34: {
 35:   PetscErrorCode  ierr;
 36:   DM_Stag * const stag  = (DM_Stag*)dm->data;
 37:   DM_Stag         *newstag;
 38:   PetscInt        dim;
 39:   MPI_Comm        newcomm;

 42:   newcomm = (comm == MPI_COMM_NULL) ? PetscObjectComm((PetscObject)dm) : comm;
 43:   DMCreate(newcomm,newdm);
 44:   DMGetDimension(dm,&dim);
 45:   DMSetDimension(*newdm,dim);

 47:   /* Call routine to define all data required for setup */
 48:   DMStagInitialize(stag->boundaryType[0],stag->boundaryType[1],stag->boundaryType[2],stag->N[0],stag->N[1],stag->N[2],stag->nRanks[0],stag->nRanks[1],stag->nRanks[2],stag->dof[0],stag->dof[1],stag->dof[2],stag->dof[3],stag->stencilType,stag->stencilWidth,stag->l[0],stag->l[1],stag->l[2],*newdm);

 50:   /* Copy all data unrelated to setup */
 51:   newstag = (DM_Stag*)(*newdm)->data;
 52:   PetscStrallocpy(stag->coordinateDMType,(char**)&newstag->coordinateDMType);
 53:   return(0);
 54: }

 56: /* Populate data created after DMCreate_Stag() is called, which is used by DMSetUp_Stag(),
 57:    such as the grid dimensions and dof information. Arguments are ignored for dimensions
 58:    less than three. */
 59: PetscErrorCode DMStagInitialize(DMBoundaryType bndx,DMBoundaryType bndy,DMBoundaryType bndz,PetscInt M,PetscInt N,PetscInt P,PetscInt m,PetscInt n,PetscInt p,PetscInt dof0,PetscInt dof1,PetscInt dof2,PetscInt dof3,DMStagStencilType stencilType,PetscInt stencilWidth,const PetscInt lx[],const PetscInt ly[],const PetscInt lz[],DM dm)
 60: {

 64:   DMSetType(dm,DMSTAG);
 65:   DMStagSetBoundaryTypes(dm,bndx,bndy,bndz);
 66:   DMStagSetGlobalSizes(dm,M,N,P);
 67:   DMStagSetNumRanks(dm,m,n,p);
 68:   DMStagSetStencilType(dm,stencilType);
 69:   DMStagSetStencilWidth(dm,stencilWidth);
 70:   DMStagSetDOF(dm,dof0,dof1,dof2,dof3);
 71:   DMStagSetOwnershipRanges(dm,lx,ly,lz);
 72:   return(0);
 73: }