1 #include <petscconvest.h> /*I "petscconvest.h" I*/ 2 #include <petscts.h> 3 #include <petscdmplex.h> 4 5 #include <petsc/private/petscconvestimpl.h> 6 7 static PetscErrorCode PetscConvEstSetTS_Private(PetscConvEst ce, PetscObject solver) 8 { 9 PetscClassId id; 10 PetscErrorCode ierr; 11 12 PetscFunctionBegin; 13 ierr = PetscObjectGetClassId(ce->solver, &id);CHKERRQ(ierr); 14 PetscCheckFalse(id != TS_CLASSID,PetscObjectComm((PetscObject) ce), PETSC_ERR_ARG_WRONG, "Solver was not a TS"); 15 ierr = TSGetDM((TS) ce->solver, &ce->idm);CHKERRQ(ierr); 16 PetscFunctionReturn(0); 17 } 18 19 static PetscErrorCode PetscConvEstInitGuessTS_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u) 20 { 21 PetscErrorCode ierr; 22 23 PetscFunctionBegin; 24 ierr = TSComputeInitialCondition((TS) ce->solver, u);CHKERRQ(ierr); 25 PetscFunctionReturn(0); 26 } 27 28 static PetscErrorCode PetscConvEstComputeErrorTS_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[]) 29 { 30 TS ts = (TS) ce->solver; 31 PetscErrorCode (*exactError)(TS, Vec, Vec); 32 PetscErrorCode ierr; 33 34 PetscFunctionBegin; 35 ierr = TSGetComputeExactError(ts, &exactError);CHKERRQ(ierr); 36 if (exactError) { 37 Vec e; 38 PetscInt f; 39 40 ierr = VecDuplicate(u, &e);CHKERRQ(ierr); 41 ierr = TSComputeExactError(ts, u, e);CHKERRQ(ierr); 42 ierr = VecNorm(e, NORM_2, errors);CHKERRQ(ierr); 43 for (f = 1; f < ce->Nf; ++f) errors[f] = errors[0]; 44 ierr = VecDestroy(&e);CHKERRQ(ierr); 45 } else { 46 PetscReal t; 47 48 ierr = TSGetSolveTime(ts, &t);CHKERRQ(ierr); 49 ierr = DMComputeL2FieldDiff(dm, t, ce->exactSol, ce->ctxs, u, errors);CHKERRQ(ierr); 50 } 51 PetscFunctionReturn(0); 52 } 53 54 static PetscErrorCode PetscConvEstGetConvRateTS_Temporal_Private(PetscConvEst ce, PetscReal alpha[]) 55 { 56 TS ts = (TS) ce->solver; 57 Vec u; 58 PetscReal *dt, *x, *y, slope, intercept; 59 PetscInt Ns, oNs, Nf = ce->Nf, f, Nr = ce->Nr, r; 60 PetscErrorCode ierr; 61 62 PetscFunctionBegin; 63 ierr = TSGetSolution(ts, &u);CHKERRQ(ierr); 64 ierr = PetscMalloc1(Nr+1, &dt);CHKERRQ(ierr); 65 ierr = TSGetTimeStep(ts, &dt[0]);CHKERRQ(ierr); 66 ierr = TSGetMaxSteps(ts, &oNs);CHKERRQ(ierr); 67 Ns = oNs; 68 for (r = 0; r <= Nr; ++r) { 69 if (r > 0) { 70 dt[r] = dt[r-1]/ce->r; 71 Ns = PetscCeilReal(Ns*ce->r); 72 } 73 ierr = TSSetTime(ts, 0.0);CHKERRQ(ierr); 74 ierr = TSSetStepNumber(ts, 0);CHKERRQ(ierr); 75 ierr = TSSetTimeStep(ts, dt[r]);CHKERRQ(ierr); 76 ierr = TSSetMaxSteps(ts, Ns);CHKERRQ(ierr); 77 ierr = PetscConvEstComputeInitialGuess(ce, r, NULL, u);CHKERRQ(ierr); 78 ierr = TSSolve(ts, u);CHKERRQ(ierr); 79 ierr = PetscLogEventBegin(ce->event, ce, 0, 0, 0);CHKERRQ(ierr); 80 ierr = PetscConvEstComputeError(ce, r, ce->idm, u, &ce->errors[r*Nf]);CHKERRQ(ierr); 81 ierr = PetscLogEventEnd(ce->event, ce, 0, 0, 0);CHKERRQ(ierr); 82 for (f = 0; f < Nf; ++f) { 83 ce->dofs[r*Nf+f] = 1.0/dt[r]; 84 ierr = PetscLogEventSetDof(ce->event, f, ce->dofs[r*Nf+f]);CHKERRQ(ierr); 85 ierr = PetscLogEventSetError(ce->event, f, ce->errors[r*Nf+f]);CHKERRQ(ierr); 86 } 87 /* Monitor */ 88 ierr = PetscConvEstMonitorDefault(ce, r);CHKERRQ(ierr); 89 } 90 /* Fit convergence rate */ 91 if (Nr) { 92 ierr = PetscMalloc2(Nr+1, &x, Nr+1, &y);CHKERRQ(ierr); 93 for (f = 0; f < Nf; ++f) { 94 for (r = 0; r <= Nr; ++r) { 95 x[r] = PetscLog10Real(dt[r]); 96 y[r] = PetscLog10Real(ce->errors[r*Nf+f]); 97 } 98 ierr = PetscLinearRegression(Nr+1, x, y, &slope, &intercept);CHKERRQ(ierr); 99 /* Since lg err = s lg dt + b */ 100 alpha[f] = slope; 101 } 102 ierr = PetscFree2(x, y);CHKERRQ(ierr); 103 } 104 /* Reset solver */ 105 ierr = TSSetConvergedReason(ts, TS_CONVERGED_ITERATING);CHKERRQ(ierr); 106 ierr = TSSetTime(ts, 0.0);CHKERRQ(ierr); 107 ierr = TSSetStepNumber(ts, 0);CHKERRQ(ierr); 108 ierr = TSSetTimeStep(ts, dt[0]);CHKERRQ(ierr); 109 ierr = TSSetMaxSteps(ts, oNs);CHKERRQ(ierr); 110 ierr = PetscConvEstComputeInitialGuess(ce, 0, NULL, u);CHKERRQ(ierr); 111 ierr = PetscFree(dt);CHKERRQ(ierr); 112 PetscFunctionReturn(0); 113 } 114 115 static PetscErrorCode PetscConvEstGetConvRateTS_Spatial_Private(PetscConvEst ce, PetscReal alpha[]) 116 { 117 TS ts = (TS) ce->solver; 118 Vec uInitial; 119 DM *dm; 120 PetscObject disc; 121 PetscReal *x, *y, slope, intercept; 122 PetscInt Nr = ce->Nr, r, Nf = ce->Nf, f, dim, oldlevel, oldnlev; 123 void *ctx; 124 PetscErrorCode ierr; 125 126 PetscFunctionBegin; 127 PetscCheckFalse(ce->r != 2.0,PetscObjectComm((PetscObject) ce), PETSC_ERR_SUP, "Only refinement factor 2 is currently supported (not %g)", (double) ce->r); 128 ierr = DMGetDimension(ce->idm, &dim);CHKERRQ(ierr); 129 ierr = DMGetApplicationContext(ce->idm, &ctx);CHKERRQ(ierr); 130 ierr = DMPlexSetRefinementUniform(ce->idm, PETSC_TRUE);CHKERRQ(ierr); 131 ierr = DMGetRefineLevel(ce->idm, &oldlevel);CHKERRQ(ierr); 132 ierr = PetscMalloc1((Nr+1), &dm);CHKERRQ(ierr); 133 ierr = TSGetSolution(ts, &uInitial);CHKERRQ(ierr); 134 /* Loop over meshes */ 135 dm[0] = ce->idm; 136 for (r = 0; r <= Nr; ++r) { 137 Vec u; 138 #if defined(PETSC_USE_LOG) 139 PetscLogStage stage; 140 #endif 141 char stageName[PETSC_MAX_PATH_LEN]; 142 const char *dmname, *uname; 143 144 ierr = PetscSNPrintf(stageName, PETSC_MAX_PATH_LEN-1, "ConvEst Refinement Level %D", r);CHKERRQ(ierr); 145 #if defined(PETSC_USE_LOG) 146 ierr = PetscLogStageGetId(stageName, &stage);CHKERRQ(ierr); 147 if (stage < 0) {ierr = PetscLogStageRegister(stageName, &stage);CHKERRQ(ierr);} 148 #endif 149 ierr = PetscLogStagePush(stage);CHKERRQ(ierr); 150 if (r > 0) { 151 if (!ce->noRefine) { 152 ierr = DMRefine(dm[r-1], MPI_COMM_NULL, &dm[r]);CHKERRQ(ierr); 153 ierr = DMSetCoarseDM(dm[r], dm[r-1]);CHKERRQ(ierr); 154 } else { 155 DM cdm, rcdm; 156 157 ierr = DMClone(dm[r-1], &dm[r]);CHKERRQ(ierr); 158 ierr = DMCopyDisc(dm[r-1], dm[r]);CHKERRQ(ierr); 159 ierr = DMGetCoordinateDM(dm[r-1], &cdm);CHKERRQ(ierr); 160 ierr = DMGetCoordinateDM(dm[r], &rcdm);CHKERRQ(ierr); 161 ierr = DMCopyDisc(cdm, rcdm);CHKERRQ(ierr); 162 } 163 ierr = DMCopyTransform(ce->idm, dm[r]);CHKERRQ(ierr); 164 ierr = PetscObjectGetName((PetscObject) dm[r-1], &dmname);CHKERRQ(ierr); 165 ierr = PetscObjectSetName((PetscObject) dm[r], dmname);CHKERRQ(ierr); 166 for (f = 0; f <= Nf; ++f) { 167 PetscErrorCode (*nspconstr)(DM, PetscInt, PetscInt, MatNullSpace *); 168 169 ierr = DMGetNullSpaceConstructor(dm[r-1], f, &nspconstr);CHKERRQ(ierr); 170 ierr = DMSetNullSpaceConstructor(dm[r], f, nspconstr);CHKERRQ(ierr); 171 } 172 } 173 ierr = DMViewFromOptions(dm[r], NULL, "-conv_dm_view");CHKERRQ(ierr); 174 /* Create solution */ 175 ierr = DMCreateGlobalVector(dm[r], &u);CHKERRQ(ierr); 176 ierr = DMGetField(dm[r], 0, NULL, &disc);CHKERRQ(ierr); 177 ierr = PetscObjectGetName(disc, &uname);CHKERRQ(ierr); 178 ierr = PetscObjectSetName((PetscObject) u, uname);CHKERRQ(ierr); 179 /* Setup solver */ 180 ierr = TSReset(ts);CHKERRQ(ierr); 181 ierr = TSSetDM(ts, dm[r]);CHKERRQ(ierr); 182 ierr = DMTSSetBoundaryLocal(dm[r], DMPlexTSComputeBoundary, ctx);CHKERRQ(ierr); 183 ierr = DMTSSetIFunctionLocal(dm[r], DMPlexTSComputeIFunctionFEM, ctx);CHKERRQ(ierr); 184 ierr = DMTSSetIJacobianLocal(dm[r], DMPlexTSComputeIJacobianFEM, ctx);CHKERRQ(ierr); 185 ierr = TSSetTime(ts, 0.0);CHKERRQ(ierr); 186 ierr = TSSetStepNumber(ts, 0);CHKERRQ(ierr); 187 ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 188 /* Create initial guess */ 189 ierr = PetscConvEstComputeInitialGuess(ce, r, dm[r], u);CHKERRQ(ierr); 190 ierr = TSSolve(ts, u);CHKERRQ(ierr); 191 ierr = PetscLogEventBegin(ce->event, ce, 0, 0, 0);CHKERRQ(ierr); 192 ierr = PetscConvEstComputeError(ce, r, dm[r], u, &ce->errors[r*Nf]);CHKERRQ(ierr); 193 ierr = PetscLogEventEnd(ce->event, ce, 0, 0, 0);CHKERRQ(ierr); 194 for (f = 0; f < Nf; ++f) { 195 PetscSection s, fs; 196 PetscInt lsize; 197 198 /* Could use DMGetOutputDM() to add in Dirichlet dofs */ 199 ierr = DMGetLocalSection(dm[r], &s);CHKERRQ(ierr); 200 ierr = PetscSectionGetField(s, f, &fs);CHKERRQ(ierr); 201 ierr = PetscSectionGetConstrainedStorageSize(fs, &lsize);CHKERRQ(ierr); 202 ierr = MPI_Allreduce(&lsize, &ce->dofs[r*Nf+f], 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) ts));CHKERRMPI(ierr); 203 ierr = PetscLogEventSetDof(ce->event, f, ce->dofs[r*Nf+f]);CHKERRQ(ierr); 204 ierr = PetscLogEventSetError(ce->event, f, ce->errors[r*Nf+f]);CHKERRQ(ierr); 205 } 206 /* Monitor */ 207 ierr = PetscConvEstMonitorDefault(ce, r);CHKERRQ(ierr); 208 if (!r) { 209 /* PCReset() does not wipe out the level structure */ 210 SNES snes; 211 KSP ksp; 212 PC pc; 213 214 ierr = TSGetSNES(ts, &snes);CHKERRQ(ierr); 215 ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr); 216 ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 217 ierr = PCMGGetLevels(pc, &oldnlev);CHKERRQ(ierr); 218 } 219 /* Cleanup */ 220 ierr = VecDestroy(&u);CHKERRQ(ierr); 221 ierr = PetscLogStagePop();CHKERRQ(ierr); 222 } 223 for (r = 1; r <= Nr; ++r) { 224 ierr = DMDestroy(&dm[r]);CHKERRQ(ierr); 225 } 226 /* Fit convergence rate */ 227 ierr = PetscMalloc2(Nr+1, &x, Nr+1, &y);CHKERRQ(ierr); 228 for (f = 0; f < Nf; ++f) { 229 for (r = 0; r <= Nr; ++r) { 230 x[r] = PetscLog10Real(ce->dofs[r*Nf+f]); 231 y[r] = PetscLog10Real(ce->errors[r*Nf+f]); 232 } 233 ierr = PetscLinearRegression(Nr+1, x, y, &slope, &intercept);CHKERRQ(ierr); 234 /* Since h^{-dim} = N, lg err = s lg N + b = -s dim lg h + b */ 235 alpha[f] = -slope * dim; 236 } 237 ierr = PetscFree2(x, y);CHKERRQ(ierr); 238 ierr = PetscFree(dm);CHKERRQ(ierr); 239 /* Restore solver */ 240 ierr = TSReset(ts);CHKERRQ(ierr); 241 { 242 /* PCReset() does not wipe out the level structure */ 243 SNES snes; 244 KSP ksp; 245 PC pc; 246 247 ierr = TSGetSNES(ts, &snes);CHKERRQ(ierr); 248 ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr); 249 ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 250 ierr = PCMGSetLevels(pc, oldnlev, NULL);CHKERRQ(ierr); 251 ierr = DMSetRefineLevel(ce->idm, oldlevel);CHKERRQ(ierr); /* The damn DMCoarsen() calls in PCMG can reset this */ 252 } 253 ierr = TSSetDM(ts, ce->idm);CHKERRQ(ierr); 254 ierr = DMTSSetBoundaryLocal(ce->idm, DMPlexTSComputeBoundary, ctx);CHKERRQ(ierr); 255 ierr = DMTSSetIFunctionLocal(ce->idm, DMPlexTSComputeIFunctionFEM, ctx);CHKERRQ(ierr); 256 ierr = DMTSSetIJacobianLocal(ce->idm, DMPlexTSComputeIJacobianFEM, ctx);CHKERRQ(ierr); 257 ierr = TSSetConvergedReason(ts, TS_CONVERGED_ITERATING);CHKERRQ(ierr); 258 ierr = TSSetTime(ts, 0.0);CHKERRQ(ierr); 259 ierr = TSSetStepNumber(ts, 0);CHKERRQ(ierr); 260 ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 261 ierr = TSSetSolution(ts, uInitial);CHKERRQ(ierr); 262 ierr = PetscConvEstComputeInitialGuess(ce, 0, NULL, uInitial);CHKERRQ(ierr); 263 PetscFunctionReturn(0); 264 } 265 266 PetscErrorCode PetscConvEstUseTS(PetscConvEst ce, PetscBool checkTemporal) 267 { 268 PetscFunctionBegin; 269 PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1); 270 ce->ops->setsolver = PetscConvEstSetTS_Private; 271 ce->ops->initguess = PetscConvEstInitGuessTS_Private; 272 ce->ops->computeerror = PetscConvEstComputeErrorTS_Private; 273 if (checkTemporal) { 274 ce->ops->getconvrate = PetscConvEstGetConvRateTS_Temporal_Private; 275 } else { 276 ce->ops->getconvrate = PetscConvEstGetConvRateTS_Spatial_Private; 277 } 278 PetscFunctionReturn(0); 279 } 280