1 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 3 #include <petsc/private/snesimpl.h> 4 #include <petscds.h> 5 #include <petscfv.h> 6 7 static PetscErrorCode DMTSConvertPlex(DM dm, DM *plex, PetscBool copy) 8 { 9 PetscBool isPlex; 10 PetscErrorCode ierr; 11 12 PetscFunctionBegin; 13 ierr = PetscObjectTypeCompare((PetscObject) dm, DMPLEX, &isPlex);CHKERRQ(ierr); 14 if (isPlex) { 15 *plex = dm; 16 ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr); 17 } else { 18 ierr = PetscObjectQuery((PetscObject) dm, "dm_plex", (PetscObject *) plex);CHKERRQ(ierr); 19 if (!*plex) { 20 ierr = DMConvert(dm,DMPLEX,plex);CHKERRQ(ierr); 21 ierr = PetscObjectCompose((PetscObject) dm, "dm_plex", (PetscObject) *plex);CHKERRQ(ierr); 22 if (copy) { 23 ierr = DMCopyDMTS(dm, *plex);CHKERRQ(ierr); 24 ierr = DMCopyDMSNES(dm, *plex);CHKERRQ(ierr); 25 ierr = DMCopyAuxiliaryVec(dm, *plex);CHKERRQ(ierr); 26 } 27 } else { 28 ierr = PetscObjectReference((PetscObject) *plex);CHKERRQ(ierr); 29 } 30 } 31 PetscFunctionReturn(0); 32 } 33 34 /*@ 35 DMPlexTSComputeRHSFunctionFVM - Form the local forcing F from the local input X using pointwise functions specified by the user 36 37 Input Parameters: 38 + dm - The mesh 39 . t - The time 40 . locX - Local solution 41 - user - The user context 42 43 Output Parameter: 44 . F - Global output vector 45 46 Level: developer 47 48 .seealso: DMPlexComputeJacobianActionFEM() 49 @*/ 50 PetscErrorCode DMPlexTSComputeRHSFunctionFVM(DM dm, PetscReal time, Vec locX, Vec F, void *user) 51 { 52 Vec locF; 53 IS cellIS; 54 DM plex; 55 PetscInt depth; 56 PetscHashFormKey key = {NULL, 0, 0}; 57 PetscErrorCode ierr; 58 59 PetscFunctionBegin; 60 ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr); 61 ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr); 62 ierr = DMGetStratumIS(plex, "dim", depth, &cellIS);CHKERRQ(ierr); 63 if (!cellIS) { 64 ierr = DMGetStratumIS(plex, "depth", depth, &cellIS);CHKERRQ(ierr); 65 } 66 ierr = DMGetLocalVector(plex, &locF);CHKERRQ(ierr); 67 ierr = VecZeroEntries(locF);CHKERRQ(ierr); 68 ierr = DMPlexComputeResidual_Internal(plex, key, cellIS, time, locX, NULL, time, locF, user);CHKERRQ(ierr); 69 ierr = DMLocalToGlobalBegin(plex, locF, ADD_VALUES, F);CHKERRQ(ierr); 70 ierr = DMLocalToGlobalEnd(plex, locF, ADD_VALUES, F);CHKERRQ(ierr); 71 ierr = DMRestoreLocalVector(plex, &locF);CHKERRQ(ierr); 72 ierr = ISDestroy(&cellIS);CHKERRQ(ierr); 73 ierr = DMDestroy(&plex);CHKERRQ(ierr); 74 PetscFunctionReturn(0); 75 } 76 77 /*@ 78 DMPlexTSComputeBoundary - Insert the essential boundary values for the local input X and/or its time derivative X_t using pointwise functions specified by the user 79 80 Input Parameters: 81 + dm - The mesh 82 . t - The time 83 . locX - Local solution 84 . locX_t - Local solution time derivative, or NULL 85 - user - The user context 86 87 Level: developer 88 89 .seealso: DMPlexComputeJacobianActionFEM() 90 @*/ 91 PetscErrorCode DMPlexTSComputeBoundary(DM dm, PetscReal time, Vec locX, Vec locX_t, void *user) 92 { 93 DM plex; 94 Vec faceGeometryFVM = NULL; 95 PetscInt Nf, f; 96 PetscErrorCode ierr; 97 98 PetscFunctionBegin; 99 ierr = DMTSConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr); 100 ierr = DMGetNumFields(plex, &Nf);CHKERRQ(ierr); 101 if (!locX_t) { 102 /* This is the RHS part */ 103 for (f = 0; f < Nf; f++) { 104 PetscObject obj; 105 PetscClassId id; 106 107 ierr = DMGetField(plex, f, NULL, &obj);CHKERRQ(ierr); 108 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 109 if (id == PETSCFV_CLASSID) { 110 ierr = DMPlexGetGeometryFVM(plex, &faceGeometryFVM, NULL, NULL);CHKERRQ(ierr); 111 break; 112 } 113 } 114 } 115 ierr = DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, time, faceGeometryFVM, NULL, NULL);CHKERRQ(ierr); 116 ierr = DMPlexInsertTimeDerivativeBoundaryValues(plex, PETSC_TRUE, locX_t, time, faceGeometryFVM, NULL, NULL);CHKERRQ(ierr); 117 ierr = DMDestroy(&plex);CHKERRQ(ierr); 118 PetscFunctionReturn(0); 119 } 120 121 /*@ 122 DMPlexTSComputeIFunctionFEM - Form the local residual F from the local input X using pointwise functions specified by the user 123 124 Input Parameters: 125 + dm - The mesh 126 . t - The time 127 . locX - Local solution 128 . locX_t - Local solution time derivative, or NULL 129 - user - The user context 130 131 Output Parameter: 132 . locF - Local output vector 133 134 Level: developer 135 136 .seealso: DMPlexComputeJacobianActionFEM() 137 @*/ 138 PetscErrorCode DMPlexTSComputeIFunctionFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, Vec locF, void *user) 139 { 140 DM plex; 141 IS allcellIS; 142 PetscInt Nds, s; 143 PetscErrorCode ierr; 144 145 PetscFunctionBegin; 146 ierr = DMTSConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr); 147 ierr = DMPlexGetAllCells_Internal(plex, &allcellIS);CHKERRQ(ierr); 148 ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr); 149 for (s = 0; s < Nds; ++s) { 150 PetscDS ds; 151 IS cellIS; 152 PetscHashFormKey key; 153 154 ierr = DMGetRegionNumDS(dm, s, &key.label, NULL, &ds);CHKERRQ(ierr); 155 key.value = 0; 156 key.field = 0; 157 if (!key.label) { 158 ierr = PetscObjectReference((PetscObject) allcellIS);CHKERRQ(ierr); 159 cellIS = allcellIS; 160 } else { 161 IS pointIS; 162 163 key.value = 1; 164 ierr = DMLabelGetStratumIS(key.label, key.value, &pointIS);CHKERRQ(ierr); 165 ierr = ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS);CHKERRQ(ierr); 166 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 167 } 168 ierr = DMPlexComputeResidual_Internal(plex, key, cellIS, time, locX, locX_t, time, locF, user);CHKERRQ(ierr); 169 ierr = ISDestroy(&cellIS);CHKERRQ(ierr); 170 } 171 ierr = ISDestroy(&allcellIS);CHKERRQ(ierr); 172 ierr = DMDestroy(&plex);CHKERRQ(ierr); 173 PetscFunctionReturn(0); 174 } 175 176 /*@ 177 DMPlexTSComputeIJacobianFEM - Form the local Jacobian J from the local input X using pointwise functions specified by the user 178 179 Input Parameters: 180 + dm - The mesh 181 . t - The time 182 . locX - Local solution 183 . locX_t - Local solution time derivative, or NULL 184 . X_tshift - The multiplicative parameter for dF/du_t 185 - user - The user context 186 187 Output Parameter: 188 . locF - Local output vector 189 190 Level: developer 191 192 .seealso: DMPlexComputeJacobianActionFEM() 193 @*/ 194 PetscErrorCode DMPlexTSComputeIJacobianFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, PetscReal X_tShift, Mat Jac, Mat JacP, void *user) 195 { 196 DM plex; 197 IS allcellIS; 198 PetscBool hasJac, hasPrec; 199 PetscInt Nds, s; 200 PetscErrorCode ierr; 201 202 PetscFunctionBegin; 203 ierr = DMTSConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr); 204 ierr = DMPlexGetAllCells_Internal(plex, &allcellIS);CHKERRQ(ierr); 205 ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr); 206 for (s = 0; s < Nds; ++s) { 207 PetscDS ds; 208 IS cellIS; 209 PetscHashFormKey key; 210 211 ierr = DMGetRegionNumDS(dm, s, &key.label, NULL, &ds);CHKERRQ(ierr); 212 key.value = 0; 213 key.field = 0; 214 if (!key.label) { 215 ierr = PetscObjectReference((PetscObject) allcellIS);CHKERRQ(ierr); 216 cellIS = allcellIS; 217 } else { 218 IS pointIS; 219 220 key.value = 1; 221 ierr = DMLabelGetStratumIS(key.label, key.value, &pointIS);CHKERRQ(ierr); 222 ierr = ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS);CHKERRQ(ierr); 223 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 224 } 225 if (!s) { 226 ierr = PetscDSHasJacobian(ds, &hasJac);CHKERRQ(ierr); 227 ierr = PetscDSHasJacobianPreconditioner(ds, &hasPrec);CHKERRQ(ierr); 228 if (hasJac && hasPrec) {ierr = MatZeroEntries(Jac);CHKERRQ(ierr);} 229 ierr = MatZeroEntries(JacP);CHKERRQ(ierr); 230 } 231 ierr = DMPlexComputeJacobian_Internal(plex, key, cellIS, time, X_tShift, locX, locX_t, Jac, JacP, user);CHKERRQ(ierr); 232 ierr = ISDestroy(&cellIS);CHKERRQ(ierr); 233 } 234 ierr = ISDestroy(&allcellIS);CHKERRQ(ierr); 235 ierr = DMDestroy(&plex);CHKERRQ(ierr); 236 PetscFunctionReturn(0); 237 } 238 239 /*@C 240 DMTSCheckResidual - Check the residual of the exact solution 241 242 Input Parameters: 243 + ts - the TS object 244 . dm - the DM 245 . t - the time 246 . u - a DM vector 247 . u_t - a DM vector 248 - tol - A tolerance for the check, or -1 to print the results instead 249 250 Output Parameters: 251 . residual - The residual norm of the exact solution, or NULL 252 253 Level: developer 254 255 .seealso: DNTSCheckFromOptions(), DMTSCheckJacobian(), DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckJacobian() 256 @*/ 257 PetscErrorCode DMTSCheckResidual(TS ts, DM dm, PetscReal t, Vec u, Vec u_t, PetscReal tol, PetscReal *residual) 258 { 259 MPI_Comm comm; 260 Vec r; 261 PetscReal res; 262 PetscErrorCode ierr; 263 264 PetscFunctionBegin; 265 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 266 PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 267 PetscValidHeaderSpecific(u, VEC_CLASSID, 3); 268 if (residual) PetscValidRealPointer(residual, 5); 269 ierr = PetscObjectGetComm((PetscObject) ts, &comm);CHKERRQ(ierr); 270 ierr = DMComputeExactSolution(dm, t, u, u_t);CHKERRQ(ierr); 271 ierr = VecDuplicate(u, &r);CHKERRQ(ierr); 272 ierr = TSComputeIFunction(ts, t, u, u_t, r, PETSC_FALSE);CHKERRQ(ierr); 273 ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr); 274 if (tol >= 0.0) { 275 if (res > tol) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "L_2 Residual %g exceeds tolerance %g", (double) res, (double) tol); 276 } else if (residual) { 277 *residual = res; 278 } else { 279 ierr = PetscPrintf(comm, "L_2 Residual: %g\n", (double)res);CHKERRQ(ierr); 280 ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr); 281 ierr = PetscObjectCompose((PetscObject) r, "__Vec_bc_zero__", (PetscObject) dm);CHKERRQ(ierr); 282 ierr = PetscObjectSetName((PetscObject) r, "Initial Residual");CHKERRQ(ierr); 283 ierr = PetscObjectSetOptionsPrefix((PetscObject)r,"res_");CHKERRQ(ierr); 284 ierr = VecViewFromOptions(r, NULL, "-vec_view");CHKERRQ(ierr); 285 ierr = PetscObjectCompose((PetscObject) r, "__Vec_bc_zero__", NULL);CHKERRQ(ierr); 286 } 287 ierr = VecDestroy(&r);CHKERRQ(ierr); 288 PetscFunctionReturn(0); 289 } 290 291 /*@C 292 DMTSCheckJacobian - Check the Jacobian of the exact solution against the residual using the Taylor Test 293 294 Input Parameters: 295 + ts - the TS object 296 . dm - the DM 297 . t - the time 298 . u - a DM vector 299 . u_t - a DM vector 300 - tol - A tolerance for the check, or -1 to print the results instead 301 302 Output Parameters: 303 + isLinear - Flag indicaing that the function looks linear, or NULL 304 - convRate - The rate of convergence of the linear model, or NULL 305 306 Level: developer 307 308 .seealso: DNTSCheckFromOptions(), DMTSCheckResidual(), DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckResidual() 309 @*/ 310 PetscErrorCode DMTSCheckJacobian(TS ts, DM dm, PetscReal t, Vec u, Vec u_t, PetscReal tol, PetscBool *isLinear, PetscReal *convRate) 311 { 312 MPI_Comm comm; 313 PetscDS ds; 314 Mat J, M; 315 MatNullSpace nullspace; 316 PetscReal dt, shift, slope, intercept; 317 PetscBool hasJac, hasPrec, isLin = PETSC_FALSE; 318 PetscErrorCode ierr; 319 320 PetscFunctionBegin; 321 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 322 PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 323 PetscValidHeaderSpecific(u, VEC_CLASSID, 3); 324 if (isLinear) PetscValidBoolPointer(isLinear, 5); 325 if (convRate) PetscValidRealPointer(convRate, 5); 326 ierr = PetscObjectGetComm((PetscObject) ts, &comm);CHKERRQ(ierr); 327 ierr = DMComputeExactSolution(dm, t, u, u_t);CHKERRQ(ierr); 328 /* Create and view matrices */ 329 ierr = TSGetTimeStep(ts, &dt);CHKERRQ(ierr); 330 shift = 1.0/dt; 331 ierr = DMCreateMatrix(dm, &J);CHKERRQ(ierr); 332 ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 333 ierr = PetscDSHasJacobian(ds, &hasJac);CHKERRQ(ierr); 334 ierr = PetscDSHasJacobianPreconditioner(ds, &hasPrec);CHKERRQ(ierr); 335 if (hasJac && hasPrec) { 336 ierr = DMCreateMatrix(dm, &M);CHKERRQ(ierr); 337 ierr = TSComputeIJacobian(ts, t, u, u_t, shift, J, M, PETSC_FALSE);CHKERRQ(ierr); 338 ierr = PetscObjectSetName((PetscObject) M, "Preconditioning Matrix");CHKERRQ(ierr); 339 ierr = PetscObjectSetOptionsPrefix((PetscObject) M, "jacpre_");CHKERRQ(ierr); 340 ierr = MatViewFromOptions(M, NULL, "-mat_view");CHKERRQ(ierr); 341 ierr = MatDestroy(&M);CHKERRQ(ierr); 342 } else { 343 ierr = TSComputeIJacobian(ts, t, u, u_t, shift, J, J, PETSC_FALSE);CHKERRQ(ierr); 344 } 345 ierr = PetscObjectSetName((PetscObject) J, "Jacobian");CHKERRQ(ierr); 346 ierr = PetscObjectSetOptionsPrefix((PetscObject) J, "jac_");CHKERRQ(ierr); 347 ierr = MatViewFromOptions(J, NULL, "-mat_view");CHKERRQ(ierr); 348 /* Check nullspace */ 349 ierr = MatGetNullSpace(J, &nullspace);CHKERRQ(ierr); 350 if (nullspace) { 351 PetscBool isNull; 352 ierr = MatNullSpaceTest(nullspace, J, &isNull);CHKERRQ(ierr); 353 if (!isNull) SETERRQ(comm, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid."); 354 } 355 /* Taylor test */ 356 { 357 PetscRandom rand; 358 Vec du, uhat, uhat_t, r, rhat, df; 359 PetscReal h; 360 PetscReal *es, *hs, *errors; 361 PetscReal hMax = 1.0, hMin = 1e-6, hMult = 0.1; 362 PetscInt Nv, v; 363 364 /* Choose a perturbation direction */ 365 ierr = PetscRandomCreate(comm, &rand);CHKERRQ(ierr); 366 ierr = VecDuplicate(u, &du);CHKERRQ(ierr); 367 ierr = VecSetRandom(du, rand);CHKERRQ(ierr); 368 ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr); 369 ierr = VecDuplicate(u, &df);CHKERRQ(ierr); 370 ierr = MatMult(J, du, df);CHKERRQ(ierr); 371 /* Evaluate residual at u, F(u), save in vector r */ 372 ierr = VecDuplicate(u, &r);CHKERRQ(ierr); 373 ierr = TSComputeIFunction(ts, t, u, u_t, r, PETSC_FALSE);CHKERRQ(ierr); 374 /* Look at the convergence of our Taylor approximation as we approach u */ 375 for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv); 376 ierr = PetscCalloc3(Nv, &es, Nv, &hs, Nv, &errors);CHKERRQ(ierr); 377 ierr = VecDuplicate(u, &uhat);CHKERRQ(ierr); 378 ierr = VecDuplicate(u, &uhat_t);CHKERRQ(ierr); 379 ierr = VecDuplicate(u, &rhat);CHKERRQ(ierr); 380 for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv) { 381 ierr = VecWAXPY(uhat, h, du, u);CHKERRQ(ierr); 382 ierr = VecWAXPY(uhat_t, h*shift, du, u_t);CHKERRQ(ierr); 383 /* F(\hat u, \hat u_t) \approx F(u, u_t) + J(u, u_t) (uhat - u) + J_t(u, u_t) (uhat_t - u_t) = F(u) + h * J(u) du + h * shift * J_t(u) du = F(u) + h F' du */ 384 ierr = TSComputeIFunction(ts, t, uhat, uhat_t, rhat, PETSC_FALSE);CHKERRQ(ierr); 385 ierr = VecAXPBYPCZ(rhat, -1.0, -h, 1.0, r, df);CHKERRQ(ierr); 386 ierr = VecNorm(rhat, NORM_2, &errors[Nv]);CHKERRQ(ierr); 387 388 es[Nv] = PetscLog10Real(errors[Nv]); 389 hs[Nv] = PetscLog10Real(h); 390 } 391 ierr = VecDestroy(&uhat);CHKERRQ(ierr); 392 ierr = VecDestroy(&uhat_t);CHKERRQ(ierr); 393 ierr = VecDestroy(&rhat);CHKERRQ(ierr); 394 ierr = VecDestroy(&df);CHKERRQ(ierr); 395 ierr = VecDestroy(&r);CHKERRQ(ierr); 396 ierr = VecDestroy(&du);CHKERRQ(ierr); 397 for (v = 0; v < Nv; ++v) { 398 if ((tol >= 0) && (errors[v] > tol)) break; 399 else if (errors[v] > PETSC_SMALL) break; 400 } 401 if (v == Nv) isLin = PETSC_TRUE; 402 ierr = PetscLinearRegression(Nv, hs, es, &slope, &intercept);CHKERRQ(ierr); 403 ierr = PetscFree3(es, hs, errors);CHKERRQ(ierr); 404 /* Slope should be about 2 */ 405 if (tol >= 0) { 406 if (!isLin && PetscAbsReal(2 - slope) > tol) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Taylor approximation convergence rate should be 2, not %0.2f", (double) slope); 407 } else if (isLinear || convRate) { 408 if (isLinear) *isLinear = isLin; 409 if (convRate) *convRate = slope; 410 } else { 411 if (!isLin) {ierr = PetscPrintf(comm, "Taylor approximation converging at order %3.2f\n", (double) slope);CHKERRQ(ierr);} 412 else {ierr = PetscPrintf(comm, "Function appears to be linear\n");CHKERRQ(ierr);} 413 } 414 } 415 ierr = MatDestroy(&J);CHKERRQ(ierr); 416 PetscFunctionReturn(0); 417 } 418 419 /*@C 420 DMTSCheckFromOptions - Check the residual and Jacobian functions using the exact solution by outputting some diagnostic information 421 422 Input Parameters: 423 + ts - the TS object 424 - u - representative TS vector 425 426 Note: The user must call PetscDSSetExactSolution() beforehand 427 428 Level: developer 429 @*/ 430 PetscErrorCode DMTSCheckFromOptions(TS ts, Vec u) 431 { 432 DM dm; 433 SNES snes; 434 Vec sol, u_t; 435 PetscReal t; 436 PetscBool check; 437 PetscErrorCode ierr; 438 439 PetscFunctionBegin; 440 ierr = PetscOptionsHasName(((PetscObject)ts)->options,((PetscObject)ts)->prefix, "-dmts_check", &check);CHKERRQ(ierr); 441 if (!check) PetscFunctionReturn(0); 442 ierr = VecDuplicate(u, &sol);CHKERRQ(ierr); 443 ierr = VecCopy(u, sol);CHKERRQ(ierr); 444 ierr = TSSetSolution(ts, u);CHKERRQ(ierr); 445 ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 446 ierr = TSSetUp(ts);CHKERRQ(ierr); 447 ierr = TSGetSNES(ts, &snes);CHKERRQ(ierr); 448 ierr = SNESSetSolution(snes, u);CHKERRQ(ierr); 449 450 ierr = TSGetTime(ts, &t);CHKERRQ(ierr); 451 ierr = DMSNESCheckDiscretization(snes, dm, t, sol, -1.0, NULL);CHKERRQ(ierr); 452 ierr = DMGetGlobalVector(dm, &u_t);CHKERRQ(ierr); 453 ierr = DMTSCheckResidual(ts, dm, t, sol, u_t, -1.0, NULL);CHKERRQ(ierr); 454 ierr = DMTSCheckJacobian(ts, dm, t, sol, u_t, -1.0, NULL, NULL);CHKERRQ(ierr); 455 ierr = DMRestoreGlobalVector(dm, &u_t);CHKERRQ(ierr); 456 457 ierr = VecDestroy(&sol);CHKERRQ(ierr); 458 PetscFunctionReturn(0); 459 } 460