Actual source code: dmproject.c

petsc-3.6.4 2016-04-12
Report Typos and Errors
  2: #include <petsc/private/petscimpl.h>
  3: #include <petscdm.h>     /*I "petscdm.h" I*/
  4: #include <petscdmplex.h> /*I "petscdmplex.h" I*/
  5: #include <petscksp.h>    /*I "petscksp.h" I*/

  7: typedef struct _projectConstraintsCtx
  8: {
  9:   DM  dm;
 10:   Vec mask;
 11: }
 12: projectConstraintsCtx;

 16: PetscErrorCode MatMult_GlobalToLocalNormal(Mat CtC, Vec x, Vec y)
 17: {
 18:   DM                    dm;
 19:   Vec                   local, mask;
 20:   projectConstraintsCtx *ctx;
 21:   PetscErrorCode        ierr;

 24:   MatShellGetContext(CtC,&ctx);
 25:   dm   = ctx->dm;
 26:   mask = ctx->mask;
 27:   DMGetLocalVector(dm,&local);
 28:   DMGlobalToLocalBegin(dm,x,INSERT_VALUES,local);
 29:   DMGlobalToLocalEnd(dm,x,INSERT_VALUES,local);
 30:   if (mask) {VecPointwiseMult(local,mask,local);}
 31:   VecSet(y,0.);
 32:   DMLocalToGlobalBegin(dm,local,ADD_VALUES,y);
 33:   DMLocalToGlobalEnd(dm,local,ADD_VALUES,y);
 34:   DMRestoreLocalVector(dm,&local);
 35:   return(0);
 36: }

 40: /*@
 41:   DMGlobalToLocalSolve - Solve for the global vector that is mapped to a given local vector by DMGlobalToLocalBegin()/DMGlobalToLocalEnd() with mode
 42:   = INSERT_VALUES.  It is assumed that the sum of all the local vector sizes is greater than or equal to the global vector size, so the solution is
 43:   a least-squares solution.  It is also assumed that DMLocalToGlobalBegin()/DMLocalToGlobalEnd() with mode = ADD_VALUES is the adjoint of the
 44:   global-to-local map, so that the least-squares solution may be found by the normal equations.

 46:   collective

 48:   Input Parameters:
 49: + dm - The DM object
 50: . x - The local vector
 51: - y - The global vector: the input value of globalVec is used as an initial guess

 53:   Output Parameters:
 54: . y - The least-squares solution

 56:   Level: advanced

 58:   Note: If the DM is of type DMPLEX, then y is the solution of L' * D * L * y = L' * D * x, where D is a diagonal mask that is 1 for every point in
 59:   the union of the closures of the local cells and 0 otherwise.  This difference is only relevant if there are anchor points that are not in the
 60:   closure of any local cell (see DMPlexGetAnchors()/DMPlexSetAnchors()).

 62: .seealso: DMGlobalToLocalBegin(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMLocalToGlobalEnd(), DMPlexGetAnchors(), DMPlexSetAnchors()
 63: @*/
 64: PetscErrorCode DMGlobalToLocalSolve(DM dm, Vec x, Vec y)
 65: {
 66:   Mat                   CtC;
 67:   PetscInt              n, N, cStart, cEnd, c;
 68:   PetscBool             isPlex;
 69:   KSP                   ksp;
 70:   PC                    pc;
 71:   Vec                   global, mask=NULL;
 72:   projectConstraintsCtx ctx;

 76:   PetscObjectTypeCompare((PetscObject)dm,DMPLEX,&isPlex);
 77:   if (isPlex) {
 78:     /* mark points in the closure */
 79:     DMGetLocalVector(dm,&mask);
 80:     VecSet(mask,0.0);
 81:     DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);
 82:     if (cEnd > cStart) {
 83:       PetscScalar *ones;
 84:       PetscInt numValues, i;

 86:       DMPlexVecGetClosure(dm,NULL,mask,cStart,&numValues,NULL);
 87:       PetscMalloc1(numValues,&ones);
 88:       for (i = 0; i < numValues; i++) {
 89:         ones[i] = 1.;
 90:       }
 91:       for (c = cStart; c < cEnd; c++) {
 92:         DMPlexVecSetClosure(dm,NULL,mask,c,ones,INSERT_VALUES);
 93:       }
 94:       PetscFree(ones);
 95:     }
 96:   }
 97:   ctx.dm   = dm;
 98:   ctx.mask = mask;
 99:   VecGetSize(y,&N);
100:   VecGetLocalSize(y,&n);
101:   MatCreate(PetscObjectComm((PetscObject)dm),&CtC);
102:   MatSetSizes(CtC,n,n,N,N);
103:   MatSetType(CtC,MATSHELL);
104:   MatSetUp(CtC);
105:   MatShellSetContext(CtC,&ctx);
106:   MatShellSetOperation(CtC,MATOP_MULT,(void(*)(void))MatMult_GlobalToLocalNormal);
107:   KSPCreate(PetscObjectComm((PetscObject)dm),&ksp);
108:   KSPSetOperators(ksp,CtC,CtC);
109:   KSPSetType(ksp,KSPCG);
110:   KSPGetPC(ksp,&pc);
111:   PCSetType(pc,PCNONE);
112:   KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
113:   KSPSetUp(ksp);
114:   DMGetGlobalVector(dm,&global);
115:   VecSet(global,0.);
116:   if (mask) {VecPointwiseMult(x,mask,x);}
117:   DMLocalToGlobalBegin(dm,x,ADD_VALUES,global);
118:   DMLocalToGlobalEnd(dm,x,ADD_VALUES,global);
119:   KSPSetFromOptions(ksp);
120:   KSPSolve(ksp,global,y);
121:   DMRestoreGlobalVector(dm,&global);
122:   /* clean up */
123:   if (isPlex) {
124:     DMRestoreLocalVector(dm,&mask);
125:   }
126:   KSPDestroy(&ksp);
127:   MatDestroy(&CtC);
128:   return(0);
129: }

133: /*@C
134:   DMPlexProjectField - This projects the given function of the fields into the function space provided.

136:   Input Parameters:
137: + dm      - The DM
138: . U       - The input field vector
139: . funcs   - The functions to evaluate, one per field
140: - mode    - The insertion mode for values

142:   Output Parameter:
143: . X       - The output vector

145:   Level: developer

147: .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff()
148: @*/
149: PetscErrorCode DMPlexProjectField(DM dm, Vec U, void (**funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal [], PetscScalar []), InsertMode mode, Vec X)
150: {
151:   Vec            localX, localU;

156:   DMGetLocalVector(dm, &localX);
157:   DMGetLocalVector(dm, &localU);
158:   DMGlobalToLocalBegin(dm, U, INSERT_VALUES, localU);
159:   DMGlobalToLocalEnd(dm, U, INSERT_VALUES, localU);
160:   DMPlexProjectFieldLocal(dm, localU, funcs, mode, localX);
161:   DMLocalToGlobalBegin(dm, localX, mode, X);
162:   DMLocalToGlobalEnd(dm, localX, mode, X);
163:   if (mode == INSERT_VALUES || mode == INSERT_ALL_VALUES || mode == INSERT_BC_VALUES) {
164:     Mat cMat;

166:     DMGetDefaultConstraints(dm, NULL, &cMat);
167:     if (cMat) {
168:       DMGlobalToLocalSolve(dm, localX, X);
169:     }
170:   }
171:   DMRestoreLocalVector(dm, &localX);
172:   DMRestoreLocalVector(dm, &localU);
173:   return(0);
174: }