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