1 // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at 2 // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights 3 // reserved. See files LICENSE and NOTICE for details. 4 // 5 // This file is part of CEED, a collection of benchmarks, miniapps, software 6 // libraries and APIs for efficient high-order finite element and spectral 7 // element discretizations for exascale applications. For more information and 8 // source code availability see http://github.com/ceed. 9 // 10 // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 11 // a collaborative effort of two U.S. Department of Energy organizations (Office 12 // of Science and the National Nuclear Security Administration) responsible for 13 // the planning and preparation of a capable exascale ecosystem, including 14 // software, applications, hardware, advanced system engineering and early 15 // testbed platforms, in support of the nation's exascale computing imperative. 16 17 // libCEED + PETSc Example: Navier-Stokes 18 // 19 // This example demonstrates a simple usage of libCEED with PETSc to solve a 20 // Navier-Stokes problem. 21 // 22 // The code is intentionally "raw", using only low-level communication 23 // primitives. 24 // 25 // Build with: 26 // 27 // make [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>] navierstokes 28 // 29 // Sample runs: 30 // 31 // ./navierstokes -ceed /cpu/self -problem density_current -degree 1 32 // ./navierstokes -ceed /gpu/occa -problem advection -degree 1 33 // 34 //TESTARGS -ceed {ceed_resource} -test -degree 1 35 36 /// @file 37 /// Navier-Stokes example using PETSc 38 39 const char help[] = "Solve Navier-Stokes using PETSc and libCEED\n"; 40 41 #include <petscts.h> 42 #include <petscdmplex.h> 43 #include <ceed.h> 44 #include <stdbool.h> 45 #include <petscsys.h> 46 #include "common.h" 47 #include "advection.h" 48 #include "advection2d.h" 49 #include "densitycurrent.h" 50 51 // Problem Options 52 typedef enum { 53 NS_DENSITY_CURRENT = 0, 54 NS_ADVECTION = 1, 55 NS_ADVECTION2D = 2, 56 } problemType; 57 static const char *const problemTypes[] = { 58 "density_current", 59 "advection", 60 "advection2d", 61 "problemType","NS_",0 62 }; 63 64 typedef enum { 65 STAB_NONE = 0, 66 STAB_SU = 1, // Streamline Upwind 67 STAB_SUPG = 2, // Streamline Upwind Petrov-Galerkin 68 } StabilizationType; 69 static const char *const StabilizationTypes[] = { 70 "NONE", 71 "SU", 72 "SUPG", 73 "StabilizationType", "STAB_", NULL 74 }; 75 76 // Problem specific data 77 typedef struct { 78 CeedInt dim, qdatasize; 79 CeedQFunctionUser setup, ics, apply_rhs, apply_ifunction; 80 PetscErrorCode (*bc)(PetscInt, PetscReal, const PetscReal[], PetscInt, 81 PetscScalar[], void *); 82 const char *setup_loc, *ics_loc, *apply_rhs_loc, *apply_ifunction_loc; 83 const bool non_zero_time; 84 } problemData; 85 86 problemData problemOptions[] = { 87 [NS_DENSITY_CURRENT] = { 88 .dim = 3, 89 .qdatasize = 10, 90 .setup = Setup, 91 .setup_loc = Setup_loc, 92 .ics = ICsDC, 93 .ics_loc = ICsDC_loc, 94 .apply_rhs = DC, 95 .apply_rhs_loc = DC_loc, 96 .apply_ifunction = IFunction_DC, 97 .apply_ifunction_loc = IFunction_DC_loc, 98 .bc = Exact_DC, 99 .non_zero_time = false, 100 }, 101 [NS_ADVECTION] = { 102 .dim = 3, 103 .qdatasize = 10, 104 .setup = Setup, 105 .setup_loc = Setup_loc, 106 .ics = ICsAdvection, 107 .ics_loc = ICsAdvection_loc, 108 .apply_rhs = Advection, 109 .apply_rhs_loc = Advection_loc, 110 .apply_ifunction = IFunction_Advection, 111 .apply_ifunction_loc = IFunction_Advection_loc, 112 .bc = Exact_Advection, 113 .non_zero_time = false, 114 }, 115 [NS_ADVECTION2D] = { 116 .dim = 2, 117 .qdatasize = 5, 118 .setup = Setup2d, 119 .setup_loc = Setup2d_loc, 120 .ics = ICsAdvection2d, 121 .ics_loc = ICsAdvection2d_loc, 122 .apply_rhs = Advection2d, 123 .apply_rhs_loc = Advection2d_loc, 124 .apply_ifunction = IFunction_Advection2d, 125 .apply_ifunction_loc = IFunction_Advection2d_loc, 126 .bc = Exact_Advection2d, 127 .non_zero_time = true, 128 }, 129 }; 130 131 // PETSc user data 132 typedef struct User_ *User; 133 typedef struct Units_ *Units; 134 135 struct User_ { 136 MPI_Comm comm; 137 PetscInt outputfreq; 138 DM dm; 139 DM dmviz; 140 Mat interpviz; 141 Ceed ceed; 142 Units units; 143 CeedVector qceed, qdotceed, gceed; 144 CeedOperator op_rhs, op_ifunction; 145 Vec M; 146 char outputfolder[PETSC_MAX_PATH_LEN]; 147 PetscInt contsteps; 148 }; 149 150 struct Units_ { 151 // fundamental units 152 PetscScalar meter; 153 PetscScalar kilogram; 154 PetscScalar second; 155 PetscScalar Kelvin; 156 // derived units 157 PetscScalar Pascal; 158 PetscScalar JperkgK; 159 PetscScalar mpersquareds; 160 PetscScalar WpermK; 161 PetscScalar kgpercubicm; 162 PetscScalar kgpersquaredms; 163 PetscScalar Joulepercubicm; 164 }; 165 166 typedef struct SimpleBC_ *SimpleBC; 167 struct SimpleBC_ { 168 PetscInt nwall, nslip[3]; 169 PetscInt walls[6], slips[3][6]; 170 }; 171 172 // Essential BC dofs are encoded in closure indices as -(i+1). 173 static PetscInt Involute(PetscInt i) { 174 return i >= 0 ? i : -(i+1); 175 } 176 177 // Utility function to create local CEED restriction 178 static PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt P, 179 CeedElemRestriction *Erestrict) { 180 181 PetscSection section; 182 PetscInt c, cStart, cEnd, Nelem, Ndof, *erestrict, eoffset, nfields, dim; 183 PetscErrorCode ierr; 184 Vec Uloc; 185 186 PetscFunctionBeginUser; 187 ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr); 188 ierr = DMGetLocalSection(dm,§ion); CHKERRQ(ierr); 189 ierr = PetscSectionGetNumFields(section, &nfields); CHKERRQ(ierr); 190 PetscInt ncomp[nfields], fieldoff[nfields+1]; 191 fieldoff[0] = 0; 192 for (PetscInt f=0; f<nfields; f++) { 193 ierr = PetscSectionGetFieldComponents(section, f, &ncomp[f]); CHKERRQ(ierr); 194 fieldoff[f+1] = fieldoff[f] + ncomp[f]; 195 } 196 197 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd); CHKERRQ(ierr); 198 Nelem = cEnd - cStart; 199 ierr = PetscMalloc1(Nelem*PetscPowInt(P, dim), &erestrict); CHKERRQ(ierr); 200 for (c=cStart,eoffset=0; c<cEnd; c++) { 201 PetscInt numindices, *indices, nnodes; 202 ierr = DMPlexGetClosureIndices(dm, section, section, c, &numindices, 203 &indices, NULL); CHKERRQ(ierr); 204 if (numindices % fieldoff[nfields]) SETERRQ1(PETSC_COMM_SELF, 205 PETSC_ERR_ARG_INCOMP, "Number of closure indices not compatible with Cell %D", 206 c); 207 nnodes = numindices / fieldoff[nfields]; 208 for (PetscInt i=0; i<nnodes; i++) { 209 // Check that indices are blocked by node and thus can be coalesced as a single field with 210 // fieldoff[nfields] = sum(ncomp) components. 211 for (PetscInt f=0; f<nfields; f++) { 212 for (PetscInt j=0; j<ncomp[f]; j++) { 213 if (Involute(indices[fieldoff[f]*nnodes + i*ncomp[f] + j]) 214 != Involute(indices[i*ncomp[0]]) + fieldoff[f] + j) 215 SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP, 216 "Cell %D closure indices not interlaced for node %D field %D component %D", 217 c, i, f, j); 218 } 219 } 220 // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode. 221 PetscInt loc = Involute(indices[i*ncomp[0]]); 222 erestrict[eoffset++] = loc / fieldoff[nfields]; 223 } 224 ierr = DMPlexRestoreClosureIndices(dm, section, section, c, &numindices, 225 &indices, NULL); CHKERRQ(ierr); 226 } 227 if (eoffset != Nelem*PetscPowInt(P, dim)) SETERRQ3(PETSC_COMM_SELF, 228 PETSC_ERR_LIB, "ElemRestriction of size (%D,%D) initialized %D nodes", Nelem, 229 PetscPowInt(P, dim),eoffset); 230 ierr = DMGetLocalVector(dm, &Uloc); CHKERRQ(ierr); 231 ierr = VecGetLocalSize(Uloc, &Ndof); CHKERRQ(ierr); 232 ierr = DMRestoreLocalVector(dm, &Uloc); CHKERRQ(ierr); 233 CeedElemRestrictionCreate(ceed, CEED_INTERLACED, Nelem, PetscPowInt(P, dim), 234 Ndof/fieldoff[nfields], fieldoff[nfields], 235 CEED_MEM_HOST, CEED_COPY_VALUES, erestrict, Erestrict); 236 ierr = PetscFree(erestrict); CHKERRQ(ierr); 237 PetscFunctionReturn(0); 238 } 239 240 static int CreateVectorFromPetscVec(Ceed ceed, Vec p, CeedVector *v) { 241 PetscErrorCode ierr; 242 PetscInt m; 243 244 PetscFunctionBeginUser; 245 ierr = VecGetLocalSize(p, &m); CHKERRQ(ierr); 246 ierr = CeedVectorCreate(ceed, m, v); CHKERRQ(ierr); 247 PetscFunctionReturn(0); 248 } 249 250 static int VectorPlacePetscVec(CeedVector c, Vec p) { 251 PetscErrorCode ierr; 252 PetscInt mceed,mpetsc; 253 PetscScalar *a; 254 255 PetscFunctionBeginUser; 256 ierr = CeedVectorGetLength(c, &mceed); CHKERRQ(ierr); 257 ierr = VecGetLocalSize(p, &mpetsc); CHKERRQ(ierr); 258 if (mceed != mpetsc) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP, 259 "Cannot place PETSc Vec of length %D in CeedVector of length %D", 260 mpetsc, mceed); 261 ierr = VecGetArray(p, &a); CHKERRQ(ierr); 262 CeedVectorSetArray(c, CEED_MEM_HOST, CEED_USE_POINTER, a); 263 PetscFunctionReturn(0); 264 } 265 266 static PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm, 267 PetscBool insertEssential, Vec Qloc, PetscReal time, Vec faceGeomFVM, 268 Vec cellGeomFVM, Vec gradFVM) { 269 PetscErrorCode ierr; 270 Vec Qbc; 271 272 PetscFunctionBegin; 273 ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr); 274 ierr = VecAXPY(Qloc, 1., Qbc); CHKERRQ(ierr); 275 ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr); 276 PetscFunctionReturn(0); 277 } 278 279 // This is the RHS of the ODE, given as u_t = G(t,u) 280 // This function takes in a state vector Q and writes into G 281 static PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *userData) { 282 PetscErrorCode ierr; 283 User user = *(User *)userData; 284 PetscScalar *q, *g; 285 Vec Qloc, Gloc; 286 287 // Global-to-local 288 PetscFunctionBeginUser; 289 ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 290 ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr); 291 ierr = VecZeroEntries(Qloc); CHKERRQ(ierr); 292 ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr); 293 ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0, 294 NULL, NULL, NULL); CHKERRQ(ierr); 295 ierr = VecZeroEntries(Gloc); CHKERRQ(ierr); 296 297 // Ceed Vectors 298 ierr = VecGetArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr); 299 ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr); 300 CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER, q); 301 CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g); 302 303 // Apply CEED operator 304 CeedOperatorApply(user->op_rhs, user->qceed, user->gceed, 305 CEED_REQUEST_IMMEDIATE); 306 307 // Restore vectors 308 ierr = VecRestoreArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr); 309 ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr); 310 311 ierr = VecZeroEntries(G); CHKERRQ(ierr); 312 ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr); 313 314 // Inverse of the lumped mass matrix 315 ierr = VecPointwiseMult(G, G, user->M); // M is Minv 316 CHKERRQ(ierr); 317 318 ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 319 ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr); 320 PetscFunctionReturn(0); 321 } 322 323 static PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Qdot, Vec G, 324 void *userData) { 325 PetscErrorCode ierr; 326 User user = *(User *)userData; 327 const PetscScalar *q, *qdot; 328 PetscScalar *g; 329 Vec Qloc, Qdotloc, Gloc; 330 331 // Global-to-local 332 PetscFunctionBeginUser; 333 ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 334 ierr = DMGetLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr); 335 ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr); 336 ierr = VecZeroEntries(Qloc); CHKERRQ(ierr); 337 ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr); 338 ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0, 339 NULL, NULL, NULL); CHKERRQ(ierr); 340 ierr = VecZeroEntries(Qdotloc); CHKERRQ(ierr); 341 ierr = DMGlobalToLocal(user->dm, Qdot, INSERT_VALUES, Qdotloc); CHKERRQ(ierr); 342 ierr = VecZeroEntries(Gloc); CHKERRQ(ierr); 343 344 // Ceed Vectors 345 ierr = VecGetArrayRead(Qloc, &q); CHKERRQ(ierr); 346 ierr = VecGetArrayRead(Qdotloc, &qdot); CHKERRQ(ierr); 347 ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr); 348 CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER, 349 (PetscScalar *)q); 350 CeedVectorSetArray(user->qdotceed, CEED_MEM_HOST, CEED_USE_POINTER, 351 (PetscScalar *)qdot); 352 CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g); 353 354 // Apply CEED operator 355 CeedOperatorApply(user->op_ifunction, user->qceed, user->gceed, 356 CEED_REQUEST_IMMEDIATE); 357 358 // Restore vectors 359 ierr = VecRestoreArrayRead(Qloc, &q); CHKERRQ(ierr); 360 ierr = VecRestoreArrayRead(Qdotloc, &qdot); CHKERRQ(ierr); 361 ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr); 362 363 ierr = VecZeroEntries(G); CHKERRQ(ierr); 364 ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr); 365 366 ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 367 ierr = DMRestoreLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr); 368 ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr); 369 PetscFunctionReturn(0); 370 } 371 372 // User provided TS Monitor 373 static PetscErrorCode TSMonitor_NS(TS ts, PetscInt stepno, PetscReal time, 374 Vec Q, void *ctx) { 375 User user = ctx; 376 Vec Qloc; 377 char filepath[PETSC_MAX_PATH_LEN]; 378 PetscViewer viewer; 379 PetscErrorCode ierr; 380 381 // Set up output 382 PetscFunctionBeginUser; 383 // Print every 'outputfreq' steps 384 if (stepno % user->outputfreq != 0) 385 PetscFunctionReturn(0); 386 ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 387 ierr = PetscObjectSetName((PetscObject)Qloc, "StateVec"); CHKERRQ(ierr); 388 ierr = VecZeroEntries(Qloc); CHKERRQ(ierr); 389 ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr); 390 391 // Output 392 ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-%03D.vtu", 393 user->outputfolder, stepno + user->contsteps); 394 CHKERRQ(ierr); 395 ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), filepath, 396 FILE_MODE_WRITE, &viewer); CHKERRQ(ierr); 397 ierr = VecView(Qloc, viewer); CHKERRQ(ierr); 398 if (user->dmviz) { 399 Vec Qrefined, Qrefined_loc; 400 char filepath_refined[PETSC_MAX_PATH_LEN]; 401 PetscViewer viewer_refined; 402 403 ierr = DMGetGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr); 404 ierr = DMGetLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr); 405 ierr = PetscObjectSetName((PetscObject)Qrefined_loc, "Refined"); 406 CHKERRQ(ierr); 407 ierr = MatInterpolate(user->interpviz, Q, Qrefined); CHKERRQ(ierr); 408 ierr = VecZeroEntries(Qrefined_loc); CHKERRQ(ierr); 409 ierr = DMGlobalToLocal(user->dmviz, Qrefined, INSERT_VALUES, Qrefined_loc); 410 CHKERRQ(ierr); 411 ierr = PetscSNPrintf(filepath_refined, sizeof filepath_refined, 412 "%s/nsrefined-%03D.vtu", 413 user->outputfolder, stepno + user->contsteps); 414 CHKERRQ(ierr); 415 ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Qrefined), 416 filepath_refined, 417 FILE_MODE_WRITE, &viewer_refined); CHKERRQ(ierr); 418 ierr = VecView(Qrefined_loc, viewer_refined); CHKERRQ(ierr); 419 ierr = DMRestoreLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr); 420 ierr = DMRestoreGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr); 421 ierr = PetscViewerDestroy(&viewer_refined); CHKERRQ(ierr); 422 } 423 ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 424 ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 425 426 // Save data in a binary file for continuation of simulations 427 ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin", 428 user->outputfolder); CHKERRQ(ierr); 429 ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer); 430 CHKERRQ(ierr); 431 ierr = VecView(Q, viewer); CHKERRQ(ierr); 432 ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 433 434 // Save time stamp 435 // Dimensionalize time back 436 time /= user->units->second; 437 ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin", 438 user->outputfolder); CHKERRQ(ierr); 439 ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer); 440 CHKERRQ(ierr); 441 #if PETSC_VERSION_GE(3,13,0) 442 ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL); 443 #else 444 ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL, true); 445 #endif 446 CHKERRQ(ierr); 447 ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 448 449 PetscFunctionReturn(0); 450 } 451 452 static PetscErrorCode ICs_FixMultiplicity(CeedOperator op_ics, 453 CeedVector xcorners, CeedVector q0ceed, DM dm, Vec Qloc, Vec Q, 454 CeedElemRestriction restrictq, SetupContext ctxSetup, CeedScalar time) { 455 PetscErrorCode ierr; 456 CeedVector multlvec; 457 Vec Multiplicity, MultiplicityLoc; 458 459 ctxSetup->time = time; 460 ierr = VecZeroEntries(Qloc); CHKERRQ(ierr); 461 ierr = VectorPlacePetscVec(q0ceed, Qloc); CHKERRQ(ierr); 462 CeedOperatorApply(op_ics, xcorners, q0ceed, CEED_REQUEST_IMMEDIATE); 463 ierr = VecZeroEntries(Q); CHKERRQ(ierr); 464 ierr = DMLocalToGlobal(dm, Qloc, ADD_VALUES, Q); CHKERRQ(ierr); 465 466 // Fix multiplicity for output of ICs 467 ierr = DMGetLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr); 468 CeedElemRestrictionCreateVector(restrictq, &multlvec, NULL); 469 ierr = VectorPlacePetscVec(multlvec, MultiplicityLoc); CHKERRQ(ierr); 470 CeedElemRestrictionGetMultiplicity(restrictq, multlvec); 471 CeedVectorDestroy(&multlvec); 472 ierr = DMGetGlobalVector(dm, &Multiplicity); CHKERRQ(ierr); 473 ierr = VecZeroEntries(Multiplicity); CHKERRQ(ierr); 474 ierr = DMLocalToGlobal(dm, MultiplicityLoc, ADD_VALUES, Multiplicity); 475 CHKERRQ(ierr); 476 ierr = VecPointwiseDivide(Q, Q, Multiplicity); CHKERRQ(ierr); 477 ierr = VecPointwiseDivide(Qloc, Qloc, MultiplicityLoc); CHKERRQ(ierr); 478 ierr = DMRestoreLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr); 479 ierr = DMRestoreGlobalVector(dm, &Multiplicity); CHKERRQ(ierr); 480 481 PetscFunctionReturn(0); 482 } 483 484 static PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm, 485 CeedElemRestriction restrictq, CeedBasis basisq, 486 CeedElemRestriction restrictqdi, CeedVector qdata, Vec M) { 487 PetscErrorCode ierr; 488 CeedQFunction qf_mass; 489 CeedOperator op_mass; 490 CeedVector mceed; 491 Vec Mloc; 492 CeedInt ncompq, qdatasize; 493 494 PetscFunctionBeginUser; 495 CeedElemRestrictionGetNumComponents(restrictq, &ncompq); 496 CeedElemRestrictionGetNumComponents(restrictqdi, &qdatasize); 497 // Create the Q-function that defines the action of the mass operator 498 CeedQFunctionCreateInterior(ceed, 1, Mass, Mass_loc, &qf_mass); 499 CeedQFunctionAddInput(qf_mass, "q", ncompq, CEED_EVAL_INTERP); 500 CeedQFunctionAddInput(qf_mass, "qdata", qdatasize, CEED_EVAL_NONE); 501 CeedQFunctionAddOutput(qf_mass, "v", ncompq, CEED_EVAL_INTERP); 502 503 // Create the mass operator 504 CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass); 505 CeedOperatorSetField(op_mass, "q", restrictq, basisq, CEED_VECTOR_ACTIVE); 506 CeedOperatorSetField(op_mass, "qdata", restrictqdi, 507 CEED_BASIS_COLLOCATED, qdata); 508 CeedOperatorSetField(op_mass, "v", restrictq, basisq, CEED_VECTOR_ACTIVE); 509 510 ierr = DMGetLocalVector(dm, &Mloc); CHKERRQ(ierr); 511 ierr = VecZeroEntries(Mloc); CHKERRQ(ierr); 512 CeedElemRestrictionCreateVector(restrictq, &mceed, NULL); 513 ierr = VectorPlacePetscVec(mceed, Mloc); CHKERRQ(ierr); 514 515 { 516 // Compute a lumped mass matrix 517 CeedVector onesvec; 518 CeedElemRestrictionCreateVector(restrictq, &onesvec, NULL); 519 CeedVectorSetValue(onesvec, 1.0); 520 CeedOperatorApply(op_mass, onesvec, mceed, CEED_REQUEST_IMMEDIATE); 521 CeedVectorDestroy(&onesvec); 522 CeedOperatorDestroy(&op_mass); 523 CeedVectorDestroy(&mceed); 524 } 525 CeedQFunctionDestroy(&qf_mass); 526 527 ierr = VecZeroEntries(M); CHKERRQ(ierr); 528 ierr = DMLocalToGlobal(dm, Mloc, ADD_VALUES, M); CHKERRQ(ierr); 529 ierr = DMRestoreLocalVector(dm, &Mloc); CHKERRQ(ierr); 530 531 // Invert diagonally lumped mass vector for RHS function 532 ierr = VecReciprocal(M); CHKERRQ(ierr); 533 PetscFunctionReturn(0); 534 } 535 536 static PetscErrorCode SetUpDM(DM dm, problemData *problem, PetscInt degree, 537 SimpleBC bc, void *ctxSetup) { 538 PetscErrorCode ierr; 539 540 PetscFunctionBeginUser; 541 { 542 // Configure the finite element space and boundary conditions 543 PetscFE fe; 544 PetscInt ncompq = 5; 545 ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, ncompq, 546 PETSC_FALSE, degree, PETSC_DECIDE, 547 &fe); 548 ierr = PetscObjectSetName((PetscObject)fe, "Q"); CHKERRQ(ierr); 549 ierr = DMAddField(dm,NULL,(PetscObject)fe); CHKERRQ(ierr); 550 ierr = DMCreateDS(dm); CHKERRQ(ierr); 551 // Wall boundary conditions are zero energy density and zero flux for 552 // velocity in advection/advection2d, and zero velocity and zero flux 553 // for mass density and energy density in density_current 554 { 555 if (problem->bc == Exact_Advection || problem->bc == Exact_Advection2d) { 556 PetscInt comps[1] = {4}; 557 ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "Face Sets", 0, 558 1, comps, (void(*)(void))problem->bc, 559 bc->nwall, bc->walls, ctxSetup); CHKERRQ(ierr); 560 } else if (problem->bc == Exact_DC) { 561 PetscInt comps[3] = {1, 2, 3}; 562 ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "Face Sets", 0, 563 3, comps, (void(*)(void))problem->bc, 564 bc->nwall, bc->walls, ctxSetup); CHKERRQ(ierr); 565 } else { 566 MPI_Comm comm; 567 ierr = PetscObjectGetComm((PetscObject)dm, &comm); CHKERRQ(ierr); 568 SETERRQ(comm, PETSC_ERR_ARG_NULL, 569 "Undefined boundary conditions for this problem"); 570 } 571 } 572 { 573 PetscInt comps[1] = {1}; 574 ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", "Face Sets", 0, 575 1, comps, (void(*)(void))NULL, bc->nslip[0], 576 bc->slips[0], ctxSetup); CHKERRQ(ierr); 577 comps[0] = 2; 578 ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", "Face Sets", 0, 579 1, comps, (void(*)(void))NULL, bc->nslip[1], 580 bc->slips[1], ctxSetup); CHKERRQ(ierr); 581 comps[0] = 3; 582 ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", "Face Sets", 0, 583 1, comps, (void(*)(void))NULL, bc->nslip[2], 584 bc->slips[2], ctxSetup); CHKERRQ(ierr); 585 } 586 ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL); 587 CHKERRQ(ierr); 588 ierr = PetscFEDestroy(&fe); CHKERRQ(ierr); 589 } 590 { 591 // Empty name for conserved field (because there is only one field) 592 PetscSection section; 593 ierr = DMGetLocalSection(dm, §ion); CHKERRQ(ierr); 594 ierr = PetscSectionSetFieldName(section, 0, ""); CHKERRQ(ierr); 595 ierr = PetscSectionSetComponentName(section, 0, 0, "Density"); 596 CHKERRQ(ierr); 597 ierr = PetscSectionSetComponentName(section, 0, 1, "MomentumX"); 598 CHKERRQ(ierr); 599 ierr = PetscSectionSetComponentName(section, 0, 2, "MomentumY"); 600 CHKERRQ(ierr); 601 ierr = PetscSectionSetComponentName(section, 0, 3, "MomentumZ"); 602 CHKERRQ(ierr); 603 ierr = PetscSectionSetComponentName(section, 0, 4, "EnergyDensity"); 604 CHKERRQ(ierr); 605 } 606 PetscFunctionReturn(0); 607 } 608 609 int main(int argc, char **argv) { 610 PetscInt ierr; 611 MPI_Comm comm; 612 DM dm, dmcoord, dmviz; 613 Mat interpviz; 614 TS ts; 615 TSAdapt adapt; 616 User user; 617 Units units; 618 char ceedresource[4096] = "/cpu/self"; 619 PetscInt cStart, cEnd, localNelem, lnodes, steps; 620 const PetscInt ncompq = 5; 621 PetscMPIInt rank; 622 PetscScalar ftime; 623 Vec Q, Qloc, Xloc; 624 Ceed ceed; 625 CeedInt numP, numQ; 626 CeedVector xcorners, qdata, q0ceed; 627 CeedBasis basisx, basisxc, basisq; 628 CeedElemRestriction restrictx, restrictxcoord, restrictq, restrictqdi; 629 CeedQFunction qf_setup, qf_ics, qf_rhs, qf_ifunction; 630 CeedOperator op_setup, op_ics; 631 CeedScalar Rd; 632 PetscScalar WpermK, Pascal, JperkgK, mpersquareds, kgpercubicm, 633 kgpersquaredms, Joulepercubicm; 634 problemType problemChoice; 635 problemData *problem = NULL; 636 StabilizationType stab; 637 PetscBool test, implicit; 638 PetscInt viz_refine = 0; 639 struct SimpleBC_ bc = { 640 .nwall = 6, 641 .walls = {1,2,3,4,5,6}, 642 }; 643 double start, cpu_time_used; 644 645 // Create the libCEED contexts 646 PetscScalar meter = 1e-2; // 1 meter in scaled length units 647 PetscScalar second = 1e-2; // 1 second in scaled time units 648 PetscScalar kilogram = 1e-6; // 1 kilogram in scaled mass units 649 PetscScalar Kelvin = 1; // 1 Kelvin in scaled temperature units 650 CeedScalar theta0 = 300.; // K 651 CeedScalar thetaC = -15.; // K 652 CeedScalar P0 = 1.e5; // Pa 653 CeedScalar N = 0.01; // 1/s 654 CeedScalar cv = 717.; // J/(kg K) 655 CeedScalar cp = 1004.; // J/(kg K) 656 CeedScalar g = 9.81; // m/s^2 657 CeedScalar lambda = -2./3.; // - 658 CeedScalar mu = 75.; // Pa s, dynamic viscosity 659 // mu = 75 is not physical for air, but is good for numerical stability 660 CeedScalar k = 0.02638; // W/(m K) 661 CeedScalar CtauS = 0.; // dimensionless 662 CeedScalar strong_form = 0.; // [0,1] 663 PetscScalar lx = 8000.; // m 664 PetscScalar ly = 8000.; // m 665 PetscScalar lz = 4000.; // m 666 CeedScalar rc = 1000.; // m (Radius of bubble) 667 PetscScalar resx = 1000.; // m (resolution in x) 668 PetscScalar resy = 1000.; // m (resolution in y) 669 PetscScalar resz = 1000.; // m (resolution in z) 670 PetscInt outputfreq = 10; // - 671 PetscInt contsteps = 0; // - 672 PetscInt degree = 1; // - 673 PetscInt qextra = 2; // - 674 PetscReal center[3], dc_axis[3] = {0, 0, 0}; 675 676 ierr = PetscInitialize(&argc, &argv, NULL, help); 677 if (ierr) return ierr; 678 679 // Allocate PETSc context 680 ierr = PetscCalloc1(1, &user); CHKERRQ(ierr); 681 ierr = PetscMalloc1(1, &units); CHKERRQ(ierr); 682 683 // Parse command line options 684 comm = PETSC_COMM_WORLD; 685 ierr = PetscOptionsBegin(comm, NULL, "Navier-Stokes in PETSc with libCEED", 686 NULL); CHKERRQ(ierr); 687 ierr = PetscOptionsString("-ceed", "CEED resource specifier", 688 NULL, ceedresource, ceedresource, 689 sizeof(ceedresource), NULL); CHKERRQ(ierr); 690 ierr = PetscOptionsBool("-test", "Run in test mode", 691 NULL, test=PETSC_FALSE, &test, NULL); CHKERRQ(ierr); 692 problemChoice = NS_DENSITY_CURRENT; 693 ierr = PetscOptionsEnum("-problem", "Problem to solve", NULL, 694 problemTypes, (PetscEnum)problemChoice, 695 (PetscEnum *)&problemChoice, NULL); CHKERRQ(ierr); 696 problem = &problemOptions[problemChoice]; 697 ierr = PetscOptionsEnum("-stab", "Stabilization method", NULL, 698 StabilizationTypes, (PetscEnum)(stab = STAB_NONE), 699 (PetscEnum *)&stab, NULL); CHKERRQ(ierr); 700 ierr = PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", 701 NULL, implicit=PETSC_FALSE, &implicit, NULL); 702 CHKERRQ(ierr); 703 if (!implicit && stab != STAB_NONE) { 704 ierr = PetscPrintf(comm, "Warning! Use -stab only with -implicit\n"); 705 CHKERRQ(ierr); 706 } 707 { 708 PetscInt len; 709 PetscBool flg; 710 ierr = PetscOptionsIntArray("-bc_wall", 711 "Use wall boundary conditions on this list of faces", 712 NULL, bc.walls, 713 (len = sizeof(bc.walls) / sizeof(bc.walls[0]), 714 &len), &flg); CHKERRQ(ierr); 715 if (flg) bc.nwall = len; 716 for (PetscInt j=0; j<3; j++) { 717 const char *flags[3] = {"-bc_slip_x", "-bc_slip_y", "-bc_slip_z"}; 718 ierr = PetscOptionsIntArray(flags[j], 719 "Use slip boundary conditions on this list of faces", 720 NULL, bc.slips[j], 721 (len = sizeof(bc.slips[j]) / sizeof(bc.slips[j][0]), 722 &len), &flg); 723 CHKERRQ(ierr); 724 if (flg) bc.nslip[j] = len; 725 } 726 } 727 ierr = PetscOptionsInt("-viz_refine", 728 "Regular refinement levels for visualization", 729 NULL, viz_refine, &viz_refine, NULL); 730 CHKERRQ(ierr); 731 ierr = PetscOptionsScalar("-units_meter", "1 meter in scaled length units", 732 NULL, meter, &meter, NULL); CHKERRQ(ierr); 733 meter = fabs(meter); 734 ierr = PetscOptionsScalar("-units_second","1 second in scaled time units", 735 NULL, second, &second, NULL); CHKERRQ(ierr); 736 second = fabs(second); 737 ierr = PetscOptionsScalar("-units_kilogram","1 kilogram in scaled mass units", 738 NULL, kilogram, &kilogram, NULL); CHKERRQ(ierr); 739 kilogram = fabs(kilogram); 740 ierr = PetscOptionsScalar("-units_Kelvin", 741 "1 Kelvin in scaled temperature units", 742 NULL, Kelvin, &Kelvin, NULL); CHKERRQ(ierr); 743 Kelvin = fabs(Kelvin); 744 ierr = PetscOptionsScalar("-theta0", "Reference potential temperature", 745 NULL, theta0, &theta0, NULL); CHKERRQ(ierr); 746 ierr = PetscOptionsScalar("-thetaC", "Perturbation of potential temperature", 747 NULL, thetaC, &thetaC, NULL); CHKERRQ(ierr); 748 ierr = PetscOptionsScalar("-P0", "Atmospheric pressure", 749 NULL, P0, &P0, NULL); CHKERRQ(ierr); 750 ierr = PetscOptionsScalar("-N", "Brunt-Vaisala frequency", 751 NULL, N, &N, NULL); CHKERRQ(ierr); 752 ierr = PetscOptionsScalar("-cv", "Heat capacity at constant volume", 753 NULL, cv, &cv, NULL); CHKERRQ(ierr); 754 ierr = PetscOptionsScalar("-cp", "Heat capacity at constant pressure", 755 NULL, cp, &cp, NULL); CHKERRQ(ierr); 756 ierr = PetscOptionsScalar("-g", "Gravitational acceleration", 757 NULL, g, &g, NULL); CHKERRQ(ierr); 758 ierr = PetscOptionsScalar("-lambda", 759 "Stokes hypothesis second viscosity coefficient", 760 NULL, lambda, &lambda, NULL); CHKERRQ(ierr); 761 ierr = PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", 762 NULL, mu, &mu, NULL); CHKERRQ(ierr); 763 ierr = PetscOptionsScalar("-k", "Thermal conductivity", 764 NULL, k, &k, NULL); CHKERRQ(ierr); 765 ierr = PetscOptionsScalar("-CtauS", 766 "Scale coefficient for tau (nondimensional)", 767 NULL, CtauS, &CtauS, NULL); CHKERRQ(ierr); 768 if (stab == STAB_NONE && CtauS != 0) { 769 ierr = PetscPrintf(comm, 770 "Warning! Use -CtauS only with -stab su or -stab supg\n"); 771 CHKERRQ(ierr); 772 } 773 ierr = PetscOptionsScalar("-strong_form", 774 "Strong (1) or weak/integrated by parts (0) advection residual", 775 NULL, strong_form, &strong_form, NULL); 776 CHKERRQ(ierr); 777 ierr = PetscOptionsScalar("-lx", "Length scale in x direction", 778 NULL, lx, &lx, NULL); CHKERRQ(ierr); 779 ierr = PetscOptionsScalar("-ly", "Length scale in y direction", 780 NULL, ly, &ly, NULL); CHKERRQ(ierr); 781 ierr = PetscOptionsScalar("-lz", "Length scale in z direction", 782 NULL, lz, &lz, NULL); CHKERRQ(ierr); 783 ierr = PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble", 784 NULL, rc, &rc, NULL); CHKERRQ(ierr); 785 ierr = PetscOptionsScalar("-resx","Target resolution in x", 786 NULL, resx, &resx, NULL); CHKERRQ(ierr); 787 ierr = PetscOptionsScalar("-resy","Target resolution in y", 788 NULL, resy, &resy, NULL); CHKERRQ(ierr); 789 ierr = PetscOptionsScalar("-resz","Target resolution in z", 790 NULL, resz, &resz, NULL); CHKERRQ(ierr); 791 PetscInt n = problem->dim; 792 center[0] = 0.5 * lx; 793 center[1] = 0.5 * ly; 794 center[2] = 0.5 * lz; 795 ierr = PetscOptionsRealArray("-center", "Location of bubble center", 796 NULL, center, &n, NULL); CHKERRQ(ierr); 797 n = problem->dim; 798 ierr = PetscOptionsRealArray("-dc_axis", 799 "Axis of density current cylindrical anomaly, or {0,0,0} for spherically symmetric", 800 NULL, dc_axis, &n, NULL); CHKERRQ(ierr); 801 { 802 PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) + 803 PetscSqr(dc_axis[1]) + PetscSqr(dc_axis[2])); 804 if (norm > 0) { 805 for (int i=0; i<3; i++) dc_axis[i] /= norm; 806 } 807 } 808 ierr = PetscOptionsInt("-output_freq", 809 "Frequency of output, in number of steps", 810 NULL, outputfreq, &outputfreq, NULL); CHKERRQ(ierr); 811 ierr = PetscOptionsInt("-continue", "Continue from previous solution", 812 NULL, contsteps, &contsteps, NULL); CHKERRQ(ierr); 813 ierr = PetscOptionsInt("-degree", "Polynomial degree of finite elements", 814 NULL, degree, °ree, NULL); CHKERRQ(ierr); 815 ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points", 816 NULL, qextra, &qextra, NULL); CHKERRQ(ierr); 817 PetscStrncpy(user->outputfolder, ".", 2); 818 ierr = PetscOptionsString("-of", "Output folder", 819 NULL, user->outputfolder, user->outputfolder, 820 sizeof(user->outputfolder), NULL); CHKERRQ(ierr); 821 ierr = PetscOptionsEnd(); CHKERRQ(ierr); 822 823 // Define derived units 824 Pascal = kilogram / (meter * PetscSqr(second)); 825 JperkgK = PetscSqr(meter) / (PetscSqr(second) * Kelvin); 826 mpersquareds = meter / PetscSqr(second); 827 WpermK = kilogram * meter / (pow(second,3) * Kelvin); 828 kgpercubicm = kilogram / pow(meter,3); 829 kgpersquaredms = kilogram / (PetscSqr(meter) * second); 830 Joulepercubicm = kilogram / (meter * PetscSqr(second)); 831 832 // Scale variables to desired units 833 theta0 *= Kelvin; 834 thetaC *= Kelvin; 835 P0 *= Pascal; 836 N *= (1./second); 837 cv *= JperkgK; 838 cp *= JperkgK; 839 Rd = cp - cv; 840 g *= mpersquareds; 841 mu *= Pascal * second; 842 k *= WpermK; 843 lx = fabs(lx) * meter; 844 ly = fabs(ly) * meter; 845 lz = fabs(lz) * meter; 846 rc = fabs(rc) * meter; 847 resx = fabs(resx) * meter; 848 resy = fabs(resy) * meter; 849 resz = fabs(resz) * meter; 850 for (int i=0; i<3; i++) center[i] *= meter; 851 852 const CeedInt dim = problem->dim, ncompx = problem->dim, 853 qdatasize = problem->qdatasize; 854 // Set up the libCEED context 855 struct SetupContext_ ctxSetup = { 856 .theta0 = theta0, 857 .thetaC = thetaC, 858 .P0 = P0, 859 .N = N, 860 .cv = cv, 861 .cp = cp, 862 .Rd = Rd, 863 .g = g, 864 .rc = rc, 865 .lx = lx, 866 .ly = ly, 867 .lz = lz, 868 .center[0] = center[0], 869 .center[1] = center[1], 870 .center[2] = center[2], 871 .dc_axis[0] = dc_axis[0], 872 .dc_axis[1] = dc_axis[1], 873 .dc_axis[2] = dc_axis[2], 874 .time = 0, 875 }; 876 877 // Create the mesh 878 { 879 const PetscReal scale[3] = {lx, ly, lz}; 880 ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, NULL, NULL, scale, 881 NULL, PETSC_TRUE, &dm); 882 CHKERRQ(ierr); 883 } 884 885 // Distribute the mesh over processes 886 { 887 DM dmDist = NULL; 888 PetscPartitioner part; 889 890 ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr); 891 ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr); 892 ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr); 893 if (dmDist) { 894 ierr = DMDestroy(&dm); CHKERRQ(ierr); 895 dm = dmDist; 896 } 897 } 898 ierr = DMViewFromOptions(dm, NULL, "-dm_view"); CHKERRQ(ierr); 899 900 // Setup DM 901 ierr = DMLocalizeCoordinates(dm); CHKERRQ(ierr); 902 ierr = DMSetFromOptions(dm); CHKERRQ(ierr); 903 ierr = SetUpDM(dm, problem, degree, &bc, &ctxSetup); CHKERRQ(ierr); 904 905 // Print FEM space information 906 if (!test) { 907 ierr = PetscPrintf(PETSC_COMM_WORLD, 908 "Degree of FEM space: %D\n", 909 degree); CHKERRQ(ierr); 910 } 911 912 // Refine DM for high-order viz 913 dmviz = NULL; 914 interpviz = NULL; 915 if (viz_refine) { 916 DM dmhierarchy[viz_refine+1]; 917 918 ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr); 919 dmhierarchy[0] = dm; 920 for (PetscInt i = 0, d = degree; i < viz_refine; i++) { 921 Mat interp_next; 922 923 ierr = DMRefine(dmhierarchy[i], MPI_COMM_NULL, &dmhierarchy[i+1]); 924 CHKERRQ(ierr); 925 ierr = DMSetCoarseDM(dmhierarchy[i+1], dmhierarchy[i]); CHKERRQ(ierr); 926 d = (d + 1) / 2; 927 if (i + 1 == viz_refine) d = 1; 928 ierr = SetUpDM(dmhierarchy[i+1], problem, d, &bc, &ctxSetup); CHKERRQ(ierr); 929 ierr = DMCreateInterpolation(dmhierarchy[i], dmhierarchy[i+1], 930 &interp_next, NULL); CHKERRQ(ierr); 931 if (!i) interpviz = interp_next; 932 else { 933 Mat C; 934 ierr = MatMatMult(interp_next, interpviz, MAT_INITIAL_MATRIX, 935 PETSC_DECIDE, &C); CHKERRQ(ierr); 936 ierr = MatDestroy(&interp_next); CHKERRQ(ierr); 937 ierr = MatDestroy(&interpviz); CHKERRQ(ierr); 938 interpviz = C; 939 } 940 } 941 for (PetscInt i=1; i<viz_refine; i++) { 942 ierr = DMDestroy(&dmhierarchy[i]); CHKERRQ(ierr); 943 } 944 dmviz = dmhierarchy[viz_refine]; 945 } 946 ierr = DMCreateGlobalVector(dm, &Q); CHKERRQ(ierr); 947 ierr = DMGetLocalVector(dm, &Qloc); CHKERRQ(ierr); 948 ierr = VecGetSize(Qloc, &lnodes); CHKERRQ(ierr); 949 lnodes /= ncompq; 950 951 { 952 // Print grid information 953 CeedInt gdofs, odofs; 954 int comm_size; 955 char box_faces_str[PETSC_MAX_PATH_LEN] = "NONE"; 956 ierr = VecGetSize(Q, &gdofs); CHKERRQ(ierr); 957 ierr = VecGetLocalSize(Q, &odofs); CHKERRQ(ierr); 958 ierr = MPI_Comm_size(comm, &comm_size); CHKERRQ(ierr); 959 ierr = PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str, 960 sizeof(box_faces_str), NULL); CHKERRQ(ierr); 961 if (!test) { 962 ierr = PetscPrintf(comm, "Global FEM dofs: %D (%D owned) on %d rank(s)\n", 963 gdofs, odofs, comm_size); CHKERRQ(ierr); 964 ierr = PetscPrintf(comm, "Local FEM nodes: %D\n", lnodes); CHKERRQ(ierr); 965 ierr = PetscPrintf(comm, "dm_plex_box_faces: %s\n", box_faces_str); 966 CHKERRQ(ierr); 967 } 968 969 } 970 971 // Set up global mass vector 972 ierr = VecDuplicate(Q, &user->M); CHKERRQ(ierr); 973 974 // Set up CEED 975 // CEED Bases 976 CeedInit(ceedresource, &ceed); 977 numP = degree + 1; 978 numQ = numP + qextra; 979 CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompq, numP, numQ, CEED_GAUSS, 980 &basisq); 981 CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numQ, CEED_GAUSS, 982 &basisx); 983 CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numP, 984 CEED_GAUSS_LOBATTO, &basisxc); 985 986 ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr); 987 ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL); 988 CHKERRQ(ierr); 989 990 // CEED Restrictions 991 ierr = CreateRestrictionFromPlex(ceed, dm, degree+1, &restrictq); 992 CHKERRQ(ierr); 993 ierr = CreateRestrictionFromPlex(ceed, dmcoord, 2, &restrictx); CHKERRQ(ierr); 994 DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd); CHKERRQ(ierr); 995 localNelem = cEnd - cStart; 996 CeedInt numQdim = CeedIntPow(numQ, dim); 997 CeedElemRestrictionCreateStrided(ceed, localNelem, numQdim, 998 localNelem*numQdim, qdatasize, 999 CEED_STRIDES_BACKEND, &restrictqdi); 1000 CeedElemRestrictionCreateStrided(ceed, localNelem, PetscPowInt(numP, dim), 1001 localNelem*PetscPowInt(numP, dim), ncompx, 1002 CEED_STRIDES_BACKEND, &restrictxcoord); 1003 1004 ierr = DMGetCoordinatesLocal(dm, &Xloc); CHKERRQ(ierr); 1005 ierr = CreateVectorFromPetscVec(ceed, Xloc, &xcorners); CHKERRQ(ierr); 1006 1007 // Create the CEED vectors that will be needed in setup 1008 CeedInt Nqpts; 1009 CeedBasisGetNumQuadraturePoints(basisq, &Nqpts); 1010 CeedVectorCreate(ceed, qdatasize*localNelem*Nqpts, &qdata); 1011 CeedElemRestrictionCreateVector(restrictq, &q0ceed, NULL); 1012 1013 // Create the Q-function that builds the quadrature data for the NS operator 1014 CeedQFunctionCreateInterior(ceed, 1, problem->setup, problem->setup_loc, 1015 &qf_setup); 1016 CeedQFunctionAddInput(qf_setup, "dx", ncompx*dim, CEED_EVAL_GRAD); 1017 CeedQFunctionAddInput(qf_setup, "weight", 1, CEED_EVAL_WEIGHT); 1018 CeedQFunctionAddOutput(qf_setup, "qdata", qdatasize, CEED_EVAL_NONE); 1019 1020 // Create the Q-function that sets the ICs of the operator 1021 CeedQFunctionCreateInterior(ceed, 1, problem->ics, problem->ics_loc, &qf_ics); 1022 CeedQFunctionAddInput(qf_ics, "x", ncompx, CEED_EVAL_INTERP); 1023 CeedQFunctionAddOutput(qf_ics, "q0", ncompq, CEED_EVAL_NONE); 1024 1025 qf_rhs = NULL; 1026 if (problem->apply_rhs) { // Create the Q-function that defines the action of the RHS operator 1027 CeedQFunctionCreateInterior(ceed, 1, problem->apply_rhs, 1028 problem->apply_rhs_loc, &qf_rhs); 1029 CeedQFunctionAddInput(qf_rhs, "q", ncompq, CEED_EVAL_INTERP); 1030 CeedQFunctionAddInput(qf_rhs, "dq", ncompq*dim, CEED_EVAL_GRAD); 1031 CeedQFunctionAddInput(qf_rhs, "qdata", qdatasize, CEED_EVAL_NONE); 1032 CeedQFunctionAddInput(qf_rhs, "x", ncompx, CEED_EVAL_INTERP); 1033 CeedQFunctionAddOutput(qf_rhs, "v", ncompq, CEED_EVAL_INTERP); 1034 CeedQFunctionAddOutput(qf_rhs, "dv", ncompq*dim, CEED_EVAL_GRAD); 1035 } 1036 1037 qf_ifunction = NULL; 1038 if (problem->apply_ifunction) { // Create the Q-function that defines the action of the IFunction 1039 CeedQFunctionCreateInterior(ceed, 1, problem->apply_ifunction, 1040 problem->apply_ifunction_loc, &qf_ifunction); 1041 CeedQFunctionAddInput(qf_ifunction, "q", ncompq, CEED_EVAL_INTERP); 1042 CeedQFunctionAddInput(qf_ifunction, "dq", ncompq*dim, CEED_EVAL_GRAD); 1043 CeedQFunctionAddInput(qf_ifunction, "qdot", ncompq, CEED_EVAL_INTERP); 1044 CeedQFunctionAddInput(qf_ifunction, "qdata", qdatasize, CEED_EVAL_NONE); 1045 CeedQFunctionAddInput(qf_ifunction, "x", ncompx, CEED_EVAL_INTERP); 1046 CeedQFunctionAddOutput(qf_ifunction, "v", ncompq, CEED_EVAL_INTERP); 1047 CeedQFunctionAddOutput(qf_ifunction, "dv", ncompq*dim, CEED_EVAL_GRAD); 1048 } 1049 1050 // Create the operator that builds the quadrature data for the NS operator 1051 CeedOperatorCreate(ceed, qf_setup, NULL, NULL, &op_setup); 1052 CeedOperatorSetField(op_setup, "dx", restrictx, basisx, CEED_VECTOR_ACTIVE); 1053 CeedOperatorSetField(op_setup, "weight", CEED_ELEMRESTRICTION_NONE, 1054 basisx, CEED_VECTOR_NONE); 1055 CeedOperatorSetField(op_setup, "qdata", restrictqdi, 1056 CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 1057 1058 // Create the operator that sets the ICs 1059 CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics); 1060 CeedOperatorSetField(op_ics, "x", restrictx, basisxc, CEED_VECTOR_ACTIVE); 1061 CeedOperatorSetField(op_ics, "q0", restrictq, 1062 CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 1063 1064 CeedElemRestrictionCreateVector(restrictq, &user->qceed, NULL); 1065 CeedElemRestrictionCreateVector(restrictq, &user->qdotceed, NULL); 1066 CeedElemRestrictionCreateVector(restrictq, &user->gceed, NULL); 1067 1068 if (qf_rhs) { // Create the RHS physics operator 1069 CeedOperator op; 1070 CeedOperatorCreate(ceed, qf_rhs, NULL, NULL, &op); 1071 CeedOperatorSetField(op, "q", restrictq, basisq, CEED_VECTOR_ACTIVE); 1072 CeedOperatorSetField(op, "dq", restrictq, basisq, CEED_VECTOR_ACTIVE); 1073 CeedOperatorSetField(op, "qdata", restrictqdi, 1074 CEED_BASIS_COLLOCATED, qdata); 1075 CeedOperatorSetField(op, "x", restrictx, basisx, xcorners); 1076 CeedOperatorSetField(op, "v", restrictq, basisq, CEED_VECTOR_ACTIVE); 1077 CeedOperatorSetField(op, "dv", restrictq, basisq, CEED_VECTOR_ACTIVE); 1078 user->op_rhs = op; 1079 } 1080 1081 if (qf_ifunction) { // Create the IFunction operator 1082 CeedOperator op; 1083 CeedOperatorCreate(ceed, qf_ifunction, NULL, NULL, &op); 1084 CeedOperatorSetField(op, "q", restrictq, basisq, CEED_VECTOR_ACTIVE); 1085 CeedOperatorSetField(op, "dq", restrictq, basisq, CEED_VECTOR_ACTIVE); 1086 CeedOperatorSetField(op, "qdot", restrictq, basisq, user->qdotceed); 1087 CeedOperatorSetField(op, "qdata", restrictqdi, 1088 CEED_BASIS_COLLOCATED, qdata); 1089 CeedOperatorSetField(op, "x", restrictx, basisx, xcorners); 1090 CeedOperatorSetField(op, "v", restrictq, basisq, CEED_VECTOR_ACTIVE); 1091 CeedOperatorSetField(op, "dv", restrictq, basisq, CEED_VECTOR_ACTIVE); 1092 user->op_ifunction = op; 1093 } 1094 1095 CeedQFunctionSetContext(qf_ics, &ctxSetup, sizeof ctxSetup); 1096 CeedScalar ctxNS[8] = {lambda, mu, k, cv, cp, g, Rd}; 1097 struct Advection2dContext_ ctxAdvection2d = { 1098 .CtauS = CtauS, 1099 .strong_form = strong_form, 1100 .stabilization = stab, 1101 }; 1102 switch (problemChoice) { 1103 case NS_DENSITY_CURRENT: 1104 if (qf_rhs) CeedQFunctionSetContext(qf_rhs, &ctxNS, sizeof ctxNS); 1105 if (qf_ifunction) CeedQFunctionSetContext(qf_ifunction, &ctxNS, 1106 sizeof ctxNS); 1107 break; 1108 case NS_ADVECTION: 1109 case NS_ADVECTION2D: 1110 if (qf_rhs) CeedQFunctionSetContext(qf_rhs, &ctxAdvection2d, 1111 sizeof ctxAdvection2d); 1112 if (qf_ifunction) CeedQFunctionSetContext(qf_ifunction, &ctxAdvection2d, 1113 sizeof ctxAdvection2d); 1114 } 1115 1116 // Set up PETSc context 1117 // Set up units structure 1118 units->meter = meter; 1119 units->kilogram = kilogram; 1120 units->second = second; 1121 units->Kelvin = Kelvin; 1122 units->Pascal = Pascal; 1123 units->JperkgK = JperkgK; 1124 units->mpersquareds = mpersquareds; 1125 units->WpermK = WpermK; 1126 units->kgpercubicm = kgpercubicm; 1127 units->kgpersquaredms = kgpersquaredms; 1128 units->Joulepercubicm = Joulepercubicm; 1129 1130 // Set up user structure 1131 user->comm = comm; 1132 user->outputfreq = outputfreq; 1133 user->contsteps = contsteps; 1134 user->units = units; 1135 user->dm = dm; 1136 user->dmviz = dmviz; 1137 user->interpviz = interpviz; 1138 user->ceed = ceed; 1139 1140 // Calculate qdata and ICs 1141 // Set up state global and local vectors 1142 ierr = VecZeroEntries(Q); CHKERRQ(ierr); 1143 1144 ierr = VectorPlacePetscVec(q0ceed, Qloc); CHKERRQ(ierr); 1145 1146 // Apply Setup Ceed Operators 1147 ierr = VectorPlacePetscVec(xcorners, Xloc); CHKERRQ(ierr); 1148 CeedOperatorApply(op_setup, xcorners, qdata, CEED_REQUEST_IMMEDIATE); 1149 ierr = ComputeLumpedMassMatrix(ceed, dm, restrictq, basisq, restrictqdi, qdata, 1150 user->M); CHKERRQ(ierr); 1151 1152 ierr = ICs_FixMultiplicity(op_ics, xcorners, q0ceed, dm, Qloc, Q, restrictq, 1153 &ctxSetup, 0.0); 1154 if (1) { // Record boundary values from initial condition and override DMPlexInsertBoundaryValues() 1155 // We use this for the main simulation DM because the reference DMPlexInsertBoundaryValues() is very slow. If we 1156 // disable this, we should still get the same results due to the problem->bc function, but with potentially much 1157 // slower execution. 1158 Vec Qbc; 1159 ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr); 1160 ierr = VecCopy(Qloc, Qbc); CHKERRQ(ierr); 1161 ierr = VecZeroEntries(Qloc); CHKERRQ(ierr); 1162 ierr = DMGlobalToLocal(dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr); 1163 ierr = VecAXPY(Qbc, -1., Qloc); CHKERRQ(ierr); 1164 ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr); 1165 ierr = PetscObjectComposeFunction((PetscObject)dm, 1166 "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_NS); 1167 CHKERRQ(ierr); 1168 } 1169 1170 MPI_Comm_rank(comm, &rank); 1171 if (!rank) {ierr = PetscMkdir(user->outputfolder); CHKERRQ(ierr);} 1172 // Gather initial Q values 1173 // In case of continuation of simulation, set up initial values from binary file 1174 if (contsteps) { // continue from existent solution 1175 PetscViewer viewer; 1176 char filepath[PETSC_MAX_PATH_LEN]; 1177 // Read input 1178 ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin", 1179 user->outputfolder); 1180 CHKERRQ(ierr); 1181 ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer); 1182 CHKERRQ(ierr); 1183 ierr = VecLoad(Q, viewer); CHKERRQ(ierr); 1184 ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 1185 } 1186 ierr = DMRestoreLocalVector(dm, &Qloc); CHKERRQ(ierr); 1187 1188 // Create and setup TS 1189 ierr = TSCreate(comm, &ts); CHKERRQ(ierr); 1190 ierr = TSSetDM(ts, dm); CHKERRQ(ierr); 1191 if (implicit) { 1192 ierr = TSSetType(ts, TSBDF); CHKERRQ(ierr); 1193 if (user->op_ifunction) { 1194 ierr = TSSetIFunction(ts, NULL, IFunction_NS, &user); CHKERRQ(ierr); 1195 } else { // Implicit integrators can fall back to using an RHSFunction 1196 ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr); 1197 } 1198 } else { 1199 if (!user->op_rhs) SETERRQ(comm,PETSC_ERR_ARG_NULL, 1200 "Problem does not provide RHSFunction"); 1201 ierr = TSSetType(ts, TSRK); CHKERRQ(ierr); 1202 ierr = TSRKSetType(ts, TSRK5F); CHKERRQ(ierr); 1203 ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr); 1204 } 1205 ierr = TSSetMaxTime(ts, 500. * units->second); CHKERRQ(ierr); 1206 ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER); CHKERRQ(ierr); 1207 ierr = TSSetTimeStep(ts, 1.e-2 * units->second); CHKERRQ(ierr); 1208 if (test) {ierr = TSSetMaxSteps(ts, 1); CHKERRQ(ierr);} 1209 ierr = TSGetAdapt(ts, &adapt); CHKERRQ(ierr); 1210 ierr = TSAdaptSetStepLimits(adapt, 1211 1.e-12 * units->second, 1212 1.e2 * units->second); CHKERRQ(ierr); 1213 ierr = TSSetFromOptions(ts); CHKERRQ(ierr); 1214 if (!contsteps) { // print initial condition 1215 if (!test) { 1216 ierr = TSMonitor_NS(ts, 0, 0., Q, user); CHKERRQ(ierr); 1217 } 1218 } else { // continue from time of last output 1219 PetscReal time; 1220 PetscInt count; 1221 PetscViewer viewer; 1222 char filepath[PETSC_MAX_PATH_LEN]; 1223 ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin", 1224 user->outputfolder); CHKERRQ(ierr); 1225 ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer); 1226 CHKERRQ(ierr); 1227 ierr = PetscViewerBinaryRead(viewer, &time, 1, &count, PETSC_REAL); 1228 CHKERRQ(ierr); 1229 ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 1230 ierr = TSSetTime(ts, time * user->units->second); CHKERRQ(ierr); 1231 } 1232 if (!test) { 1233 ierr = TSMonitorSet(ts, TSMonitor_NS, user, NULL); CHKERRQ(ierr); 1234 } 1235 1236 // Solve 1237 start = MPI_Wtime(); 1238 ierr = PetscBarrier((PetscObject)ts); CHKERRQ(ierr); 1239 ierr = TSSolve(ts, Q); CHKERRQ(ierr); 1240 cpu_time_used = MPI_Wtime() - start; 1241 ierr = TSGetSolveTime(ts, &ftime); CHKERRQ(ierr); 1242 ierr = MPI_Allreduce(MPI_IN_PLACE, &cpu_time_used, 1, MPI_DOUBLE, MPI_MIN, 1243 comm); CHKERRQ(ierr); 1244 if (!test) { 1245 ierr = PetscPrintf(PETSC_COMM_WORLD, 1246 "Time taken for solution (sec): %g\n", 1247 (double)cpu_time_used); CHKERRQ(ierr); 1248 } 1249 1250 // Get error 1251 if (problem->non_zero_time && !test) { 1252 Vec Qexact, Qexactloc; 1253 PetscReal norm; 1254 ierr = DMCreateGlobalVector(dm, &Qexact); CHKERRQ(ierr); 1255 ierr = DMGetLocalVector(dm, &Qexactloc); CHKERRQ(ierr); 1256 ierr = VecGetSize(Qexactloc, &lnodes); CHKERRQ(ierr); 1257 1258 ierr = ICs_FixMultiplicity(op_ics, xcorners, q0ceed, dm, Qexactloc, Qexact, 1259 restrictq, &ctxSetup, ftime); CHKERRQ(ierr); 1260 1261 ierr = VecAXPY(Q, -1.0, Qexact); CHKERRQ(ierr); 1262 ierr = VecNorm(Q, NORM_MAX, &norm); CHKERRQ(ierr); 1263 CeedVectorDestroy(&q0ceed); 1264 ierr = PetscPrintf(PETSC_COMM_WORLD, 1265 "Max Error: %g\n", 1266 (double)norm); CHKERRQ(ierr); 1267 // Clean up vectors 1268 ierr = DMRestoreLocalVector(dm, &Qexactloc); CHKERRQ(ierr); 1269 ierr = VecDestroy(&Qexact); CHKERRQ(ierr); 1270 } 1271 1272 // Output Statistics 1273 ierr = TSGetStepNumber(ts,&steps); CHKERRQ(ierr); 1274 if (!test) { 1275 ierr = PetscPrintf(PETSC_COMM_WORLD, 1276 "Time integrator took %D time steps to reach final time %g\n", 1277 steps, (double)ftime); CHKERRQ(ierr); 1278 } 1279 1280 // Clean up libCEED 1281 CeedVectorDestroy(&qdata); 1282 CeedVectorDestroy(&user->qceed); 1283 CeedVectorDestroy(&user->qdotceed); 1284 CeedVectorDestroy(&user->gceed); 1285 CeedVectorDestroy(&xcorners); 1286 CeedBasisDestroy(&basisq); 1287 CeedBasisDestroy(&basisx); 1288 CeedBasisDestroy(&basisxc); 1289 CeedElemRestrictionDestroy(&restrictq); 1290 CeedElemRestrictionDestroy(&restrictx); 1291 CeedElemRestrictionDestroy(&restrictqdi); 1292 CeedElemRestrictionDestroy(&restrictxcoord); 1293 CeedQFunctionDestroy(&qf_setup); 1294 CeedQFunctionDestroy(&qf_ics); 1295 CeedQFunctionDestroy(&qf_rhs); 1296 CeedQFunctionDestroy(&qf_ifunction); 1297 CeedOperatorDestroy(&op_setup); 1298 CeedOperatorDestroy(&op_ics); 1299 CeedOperatorDestroy(&user->op_rhs); 1300 CeedOperatorDestroy(&user->op_ifunction); 1301 CeedDestroy(&ceed); 1302 1303 // Clean up PETSc 1304 ierr = VecDestroy(&Q); CHKERRQ(ierr); 1305 ierr = VecDestroy(&user->M); CHKERRQ(ierr); 1306 ierr = MatDestroy(&interpviz); CHKERRQ(ierr); 1307 ierr = DMDestroy(&dmviz); CHKERRQ(ierr); 1308 ierr = TSDestroy(&ts); CHKERRQ(ierr); 1309 ierr = DMDestroy(&dm); CHKERRQ(ierr); 1310 ierr = PetscFree(units); CHKERRQ(ierr); 1311 ierr = PetscFree(user); CHKERRQ(ierr); 1312 return PetscFinalize(); 1313 } 1314 1315