Computes the Right-Hand Side (RHS) of the momentum equations.
This function calculates the contribution of the convective and diffusive terms. It is called by the momentum solvers (e.g., RungeKutta).
Computes the Right-Hand Side (RHS) of the momentum equations.
This function orchestrates the calculation of the spatial discretization of the convective and diffusive terms. It calls specialized helper functions (Convection
, Viscous
) to compute these terms and then combines them with the pressure gradient to form the final RHS vector.
This is a minimally-edited version of the legacy kernel. It operates on a single UserCtx (one grid block) and retrieves all global parameters (Re, rans, les, etc.) from the master SimCtx.
1232{
1233 PetscErrorCode ierr;
1235 DM da = user->
da, fda = user->
fda;
1236 DMDALocalInfo info = user->
info;
1237 PetscInt i,j,k;
1238 PetscReal dpdc, dpde,dpdz;
1239
1240 PetscInt xs = info.xs, xe = xs + info.xm, mx = info.mx;
1241 PetscInt ys = info.ys, ye = ys + info.ym, my = info.my;
1242 PetscInt zs = info.zs, ze = zs + info.zm, mz = info.mz;
1243 PetscInt lxs = (xs==0) ? xs+1 : xs;
1244 PetscInt lys = (ys==0) ? ys+1 : ys;
1245 PetscInt lzs = (zs==0) ? zs+1 : zs;
1246 PetscInt lxe = (xe==mx) ? xe-1 : xe;
1247 PetscInt lye = (ye==my) ? ye-1 : ye;
1248 PetscInt lze = (ze==mz) ? ze-1 : ze;
1249 PetscInt mz_end = (user->
bctype[5] == 8) ? mz - 2 : mz - 3;
1250
1251
1252 Cmpnts ***csi, ***eta, ***zet, ***icsi, ***ieta, ***izet, ***jcsi, ***jeta, ***jzet, ***kcsi, ***keta, ***kzet;
1253 PetscReal ***p, ***iaj, ***jaj, ***kaj, ***aj, ***nvert, ***nvert_o;
1254 Cmpnts ***rhs, ***rc, ***rct,***ucont;
1255
1256
1257 Vec Conv, Visc, Rc, Rct;
1258
1259 PetscFunctionBeginUser;
1262
1263
1264 ierr = DMDAVecGetArrayRead(fda, user->
lCsi, &csi); CHKERRQ(ierr);
1265 ierr = DMDAVecGetArrayRead(fda, user->
lEta, &eta); CHKERRQ(ierr);
1266 ierr = DMDAVecGetArrayRead(fda, user->
lZet, &zet); CHKERRQ(ierr);
1267 ierr = DMDAVecGetArrayRead(da, user->
lAj, &aj); CHKERRQ(ierr);
1268 ierr = DMDAVecGetArrayRead(fda, user->
lICsi, &icsi); CHKERRQ(ierr);
1269 ierr = DMDAVecGetArrayRead(fda, user->
lIEta, &ieta); CHKERRQ(ierr);
1270 ierr = DMDAVecGetArrayRead(fda, user->
lIZet, &izet); CHKERRQ(ierr);
1271 ierr = DMDAVecGetArrayRead(fda, user->
lJCsi, &jcsi); CHKERRQ(ierr);
1272 ierr = DMDAVecGetArrayRead(fda, user->
lJEta, &jeta); CHKERRQ(ierr);
1273 ierr = DMDAVecGetArrayRead(fda, user->
lJZet, &jzet); CHKERRQ(ierr);
1274 ierr = DMDAVecGetArrayRead(fda, user->
lKCsi, &kcsi); CHKERRQ(ierr);
1275 ierr = DMDAVecGetArrayRead(fda, user->
lKEta, &keta); CHKERRQ(ierr);
1276 ierr = DMDAVecGetArrayRead(fda, user->
lKZet, &kzet); CHKERRQ(ierr);
1277 ierr = DMDAVecGetArrayRead(da, user->
lIAj, &iaj); CHKERRQ(ierr);
1278 ierr = DMDAVecGetArrayRead(da, user->
lJAj, &jaj); CHKERRQ(ierr);
1279 ierr = DMDAVecGetArrayRead(da, user->
lKAj, &kaj); CHKERRQ(ierr);
1280 ierr = DMDAVecGetArrayRead(da, user->
lP, &p); CHKERRQ(ierr);
1281 ierr = DMDAVecGetArrayRead(da, user->
lNvert, &nvert); CHKERRQ(ierr);
1282 ierr = DMDAVecGetArray(fda, Rhs, &rhs); CHKERRQ(ierr);
1283
1284
1285 ierr = VecDuplicate(user->
lUcont, &Rc); CHKERRQ(ierr);
1286 ierr = VecDuplicate(Rc, &Rct); CHKERRQ(ierr);
1287 ierr = VecDuplicate(Rct, &Conv); CHKERRQ(ierr);
1288 ierr = VecDuplicate(Rct, &Visc); CHKERRQ(ierr);
1289
1290
1291
1292
1293
1294
1296
1297
1300
1301 } else {
1303 }
1304
1305
1307 ierr = VecSet(Visc, 0.0); CHKERRQ(ierr);
1308 } else {
1311 }
1312
1313
1314 ierr = VecWAXPY(Rc, -1.0, Conv, Visc); CHKERRQ(ierr);
1315
1316
1318 ierr = DMDAVecGetArray(fda, Rct, &rct); CHKERRQ(ierr);
1319 ierr = DMDAVecGetArray(fda, Rc, &rc); CHKERRQ(ierr);
1320
1321 for (k = lzs; k < lze; k++) {
1322 for (j = lys; j < lye; j++) {
1323 for (i = lxs; i < lxe; i++) {
1324 rct[k][j][i].
x = aj[k][j][i] *
1325 (0.5 * (csi[k][j][i].
x + csi[k][j][i-1].
x) * rc[k][j][i].x +
1326 0.5 * (csi[k][j][i].y + csi[k][j][i-1].y) * rc[k][j][i].
y +
1327 0.5 * (csi[k][j][i].
z + csi[k][j][i-1].
z) * rc[k][j][i].z);
1328 rct[k][j][i].
y = aj[k][j][i] *
1329 (0.5 * (eta[k][j][i].
x + eta[k][j-1][i].
x) * rc[k][j][i].x +
1330 0.5 * (eta[k][j][i].y + eta[k][j-1][i].y) * rc[k][j][i].
y +
1331 0.5 * (eta[k][j][i].
z + eta[k][j-1][i].
z) * rc[k][j][i].z);
1332 rct[k][j][i].
z = aj[k][j][i] *
1333 (0.5 * (zet[k][j][i].
x + zet[k-1][j][i].
x) * rc[k][j][i].x +
1334 0.5 * (zet[k][j][i].y + zet[k-1][j][i].y) * rc[k][j][i].
y +
1335 0.5 * (zet[k][j][i].
z + zet[k-1][j][i].
z) * rc[k][j][i].z);
1336 }
1337 }
1338 }
1339 ierr = DMDAVecRestoreArray(fda, Rct, &rct); CHKERRQ(ierr);
1340 ierr = DMDAVecRestoreArray(fda, Rc, &rc); CHKERRQ(ierr);
1341
1342 PetscBarrier(NULL);
1343
1344 DMLocalToLocalBegin(fda, Rct, INSERT_VALUES, Rct);
1345 DMLocalToLocalEnd(fda, Rct, INSERT_VALUES, Rct);
1346
1347 DMDAVecGetArray(fda, Rct, &rct);
1348
1350 if (xs==0){
1351 i=xs;
1352 for (k=lzs; k<lze; k++) {
1353 for (j=lys; j<lye; j++) {
1354 rct[k][j][i].
x=rct[k][j][i-2].
x;
1355 }
1356 }
1357 }
1358
1359 if (xe==mx){
1360 i=mx-1;
1361 for (k=lzs; k<lze; k++) {
1362 for (j=lys; j<lye; j++) {
1363 rct[k][j][i].
x=rct[k][j][i+2].
x;
1364 }
1365 }
1366 }
1367 }
1368
1370 if (ys==0){
1371 j=ys;
1372 for (k=lzs; k<lze; k++) {
1373 for (i=lxs; i<lxe; i++) {
1374 rct[k][j][i].
y=rct[k][j-2][i].
y;
1375 }
1376 }
1377 }
1378
1379 if (ye==my){
1380 j=my-1;
1381 for (k=lzs; k<lze; k++) {
1382 for (i=lxs; i<lxe; i++) {
1383 rct[k][j][i].
y=rct[k][j+2][i].
y;
1384 }
1385 }
1386 }
1387 }
1389 if (zs==0){
1390 k=zs;
1391 for (j=lys; j<lye; j++) {
1392 for (i=lxs; i<lxe; i++) {
1393 rct[k][j][i].
z=rct[k-2][j][i].
z;
1394 }
1395 }
1396 }
1397 if (ze==mz){
1398 k=mz-1;
1399 for (j=lys; j<lye; j++) {
1400 for (i=lxs; i<lxe; i++) {
1401 rct[k][j][i].
z=rct[k+2][j][i].
z;
1402 }
1403 }
1404 }
1405 }
1406
1407
1408
1409
1411
1412 ierr = DMDAVecRestoreArray(fda, Rct, &rct); CHKERRQ(ierr);
1413
1414 ierr = DMLocalToLocalBegin(fda, Rct, INSERT_VALUES, Rct); CHKERRQ(ierr);
1415 ierr = DMLocalToLocalEnd(fda, Rct, INSERT_VALUES, Rct); CHKERRQ(ierr);
1416
1417 ierr = DMDAVecGetArray(fda, Rct, &rct); CHKERRQ(ierr);
1418
1419 for (k = lzs; k < lze; k++) {
1420 for (j = lys; j < lye; j++) {
1421 for (i = lxs; i < lxe; i++) {
1422 PetscReal dpdc, dpde, dpdz;
1423 dpdc = p[k][j][i+1] - p[k][j][i];
1424
1425 if ((j==my-2 && user->
bctype[2]!=7)|| nvert[k][j+1][i]+nvert[k][j+1][i+1] > 0.1) {
1426 if (nvert[k][j-1][i] + nvert[k][j-1][i+1] < 0.1 && (j!=1 || (j==1 && user->
bctype[2]==7))) {
1427 dpde = (p[k][j ][i] + p[k][j ][i+1] -
1428 p[k][j-1][i] - p[k][j-1][i+1]) * 0.5;
1429 }
1430 }
1431 else if ((j==my-2 || j==1) && user->
bctype[2]==7 && nvert[k][j+1][i]+nvert[k][j+1][i+1] > 0.1) {
1432 if (nvert[k][j-1][i] + nvert[k][j-1][i+1] < 0.1) {
1433 dpde = (p[k][j ][i] + p[k][j ][i+1] -
1434 p[k][j-1][i] - p[k][j-1][i+1]) * 0.5;
1435 }
1436 }
1437 else if ((j == 1 && user->
bctype[2]!=7) || nvert[k][j-1][i] + nvert[k][j-1][i+1] > 0.1) {
1438 if (nvert[k][j+1][i] + nvert[k][j+1][i+1] < 0.1) {
1439 dpde = (p[k][j+1][i] + p[k][j+1][i+1] -
1440 p[k][j ][i] - p[k][j ][i+1]) * 0.5;
1441 }
1442 }
1443 else if ((j == 1 || j==my-2) && user->
bctype[2]==7 && nvert[k][j-1][i] + nvert[k][j-1][i+1] > 0.1) {
1444 if (nvert[k][j+1][i] + nvert[k][j+1][i+1] < 0.1) {
1445 dpde = (p[k][j+1][i] + p[k][j+1][i+1] -
1446 p[k][j ][i] - p[k][j ][i+1]) * 0.5;
1447 }
1448 }
1449 else {
1450 dpde = (p[k][j+1][i] + p[k][j+1][i+1] -
1451 p[k][j-1][i] - p[k][j-1][i+1]) * 0.25;
1452 }
1453
1454 if ((k == mz-2 && user->
bctype[4]!=7) || nvert[k+1][j][i] + nvert[k+1][j][i+1] > 0.1) {
1455 if (nvert[k-1][j][i] + nvert[k-1][j][i+1] < 0.1 && (k!=1 || (k==1 && user->
bctype[4]==7))) {
1456 dpdz = (p[k ][j][i] + p[k ][j][i+1] -
1457 p[k-1][j][i] - p[k-1][j][i+1]) * 0.5;
1458 }
1459 }
1460 else if ((k == mz-2 || k==1) && user->
bctype[4]==7 && nvert[k+1][j][i] + nvert[k+1][j][i+1] > 0.1) {
1461 if (nvert[k-1][j][i] + nvert[k-1][j][i+1] < 0.1) {
1462 dpdz = (p[k ][j][i] + p[k ][j][i+1] -
1463 p[k-1][j][i] - p[k-1][j][i+1]) * 0.5;
1464 }
1465 }
1466 else if ((k == 1 && user->
bctype[4]!=7)|| nvert[k-1][j][i] + nvert[k-1][j][i+1] > 0.1) {
1467 if (nvert[k+1][j][i] + nvert[k+1][j][i+1] < 0.1) {
1468 dpdz = (p[k+1][j][i] + p[k+1][j][i+1] -
1469 p[k ][j][i] - p[k ][j][i+1]) * 0.5;
1470 }
1471 }
1472 else if ((k == 1 || k==mz-2) && user->
bctype[4]==7 && nvert[k-1][j][i] + nvert[k-1][j][i+1] > 0.1) {
1473 if (nvert[k+1][j][i] + nvert[k+1][j][i+1] < 0.1) {
1474 dpdz = (p[k+1][j][i] + p[k+1][j][i+1] -
1475 p[k ][j][i] - p[k ][j][i+1]) * 0.5;
1476 }
1477 }
1478 else {
1479 dpdz = (p[k+1][j][i] + p[k+1][j][i+1] -
1480 p[k-1][j][i] - p[k-1][j][i+1]) * 0.25;
1481 }
1482
1483 rhs[k][j][i].
x =0.5 * (rct[k][j][i].
x + rct[k][j][i+1].
x);
1484
1485
1487 (dpdc * (icsi[k][j][i].
x * icsi[k][j][i].
x +
1488 icsi[k][j][i].
y * icsi[k][j][i].
y +
1489 icsi[k][j][i].
z * icsi[k][j][i].
z)+
1490 dpde * (ieta[k][j][i].x * icsi[k][j][i].x +
1491 ieta[k][j][i].y * icsi[k][j][i].y +
1492 ieta[k][j][i].z * icsi[k][j][i].z)+
1493 dpdz * (izet[k][j][i].
x * icsi[k][j][i].
x +
1494 izet[k][j][i].
y * icsi[k][j][i].
y +
1495 izet[k][j][i].
z * icsi[k][j][i].
z)) * iaj[k][j][i];
1496
1497 if ((i == mx-2 && user->
bctype[0]!=7) || nvert[k][j][i+1] + nvert[k][j+1][i+1] > 0.1) {
1498 if (nvert[k][j][i-1] + nvert[k][j+1][i-1] < 0.1 && (i!=1 || (i==1 && user->
bctype[0]==7))) {
1499 dpdc = (p[k][j][i ] + p[k][j+1][i ] -
1500 p[k][j][i-1] - p[k][j+1][i-1]) * 0.5;
1501 }
1502 }
1503 else if ((i == mx-2 || i==1) && user->
bctype[0]==7 && nvert[k][j][i+1] + nvert[k][j+1][i+1] > 0.1) {
1504 if (nvert[k][j][i-1] + nvert[k][j+1][i-1] < 0.1) {
1505 dpdc = (p[k][j][i ] + p[k][j+1][i ] -
1506 p[k][j][i-1] - p[k][j+1][i-1]) * 0.5;
1507 }
1508 }
1509 else if ((i == 1 && user->
bctype[0]!=7)|| nvert[k][j][i-1] + nvert[k][j+1][i-1] > 0.1) {
1510 if (nvert[k][j][i+1] + nvert[k][j+1][i+1] < 0.1) {
1511 dpdc = (p[k][j][i+1] + p[k][j+1][i+1] -
1512 p[k][j][i ] - p[k][j+1][i ]) * 0.5;
1513 }
1514 }
1515 else if ((i == 1 || i==mx-2) && user->
bctype[0]==7 && nvert[k][j][i-1] + nvert[k][j+1][i-1] > 0.1) {
1516 if (nvert[k][j][i+1] + nvert[k][j+1][i+1] < 0.1) {
1517 dpdc = (p[k][j][i+1] + p[k][j+1][i+1] -
1518 p[k][j][i ] - p[k][j+1][i ]) * 0.5;
1519 }
1520 }
1521 else {
1522 dpdc = (p[k][j][i+1] + p[k][j+1][i+1] -
1523 p[k][j][i-1] - p[k][j+1][i-1]) * 0.25;
1524 }
1525
1526 dpde = p[k][j+1][i] - p[k][j][i];
1527
1528 if ((k == mz-2 && user->
bctype[4]!=7)|| nvert[k+1][j][i] + nvert[k+1][j+1][i] > 0.1) {
1529 if (nvert[k-1][j][i] + nvert[k-1][j+1][i] < 0.1 && (k!=1 || (k==1 && user->
bctype[4]==7))) {
1530 dpdz = (p[k ][j][i] + p[k ][j+1][i] -
1531 p[k-1][j][i] - p[k-1][j+1][i]) * 0.5;
1532 }
1533 }
1534 else if ((k == mz-2 || k==1) && user->
bctype[4]==7 && nvert[k+1][j][i] + nvert[k+1][j+1][i] > 0.1) {
1535 if (nvert[k-1][j][i] + nvert[k-1][j+1][i] < 0.1) {
1536 dpdz = (p[k ][j][i] + p[k ][j+1][i] -
1537 p[k-1][j][i] - p[k-1][j+1][i]) * 0.5;
1538 }
1539 }
1540 else if ((k == 1 && user->
bctype[4]!=7)|| nvert[k-1][j][i] + nvert[k-1][j+1][i] > 0.1) {
1541 if (nvert[k+1][j][i] + nvert[k+1][j+1][i] < 0.1) {
1542 dpdz = (p[k+1][j][i] + p[k+1][j+1][i] -
1543 p[k ][j][i] - p[k ][j+1][i]) * 0.5;
1544 }
1545 }
1546 else if ((k == 1 || k==mz-2) && user->
bctype[4]==7 && nvert[k-1][j][i] + nvert[k-1][j+1][i] > 0.1) {
1547 if (nvert[k+1][j][i] + nvert[k+1][j+1][i] < 0.1) {
1548 dpdz = (p[k+1][j][i] + p[k+1][j+1][i] -
1549 p[k ][j][i] - p[k ][j+1][i]) * 0.5;
1550 }
1551 }
1552 else {
1553 dpdz = (p[k+1][j][i] + p[k+1][j+1][i] -
1554 p[k-1][j][i] - p[k-1][j+1][i]) * 0.25;
1555 }
1556
1557 rhs[k][j][i].
y =0.5 * (rct[k][j][i].
y + rct[k][j+1][i].
y);
1558
1559
1561 (dpdc * (jcsi[k][j][i].
x * jeta[k][j][i].
x +
1562 jcsi[k][j][i].
y * jeta[k][j][i].
y +
1563 jcsi[k][j][i].
z * jeta[k][j][i].
z) +
1564 dpde * (jeta[k][j][i].x * jeta[k][j][i].x +
1565 jeta[k][j][i].y * jeta[k][j][i].y +
1566 jeta[k][j][i].z * jeta[k][j][i].z) +
1567 dpdz * (jzet[k][j][i].
x * jeta[k][j][i].
x +
1568 jzet[k][j][i].
y * jeta[k][j][i].
y +
1569 jzet[k][j][i].
z * jeta[k][j][i].
z)) * jaj[k][j][i];
1570
1571 if ((i == mx-2 && user->
bctype[0]!=7)|| nvert[k][j][i+1] + nvert[k+1][j][i+1] > 0.1) {
1572 if (nvert[k][j][i-1] + nvert[k+1][j][i-1] < 0.1 && (i!=1 || (i==1 && user->
bctype[0]==7))) {
1573 dpdc = (p[k][j][i ] + p[k+1][j][i ] -
1574 p[k][j][i-1] - p[k+1][j][i-1]) * 0.5;
1575 }
1576 }
1577 else if ((i == mx-2 || i==1) && user->
bctype[0]==7 && nvert[k][j][i+1] + nvert[k+1][j][i+1] > 0.1) {
1578 if (nvert[k][j][i-1] + nvert[k+1][j][i-1] < 0.1) {
1579 dpdc = (p[k][j][i ] + p[k+1][j][i ] -
1580 p[k][j][i-1] - p[k+1][j][i-1]) * 0.5;
1581 }
1582 }
1583 else if ((i == 1 && user->
bctype[0]!=7)|| nvert[k][j][i-1] + nvert[k+1][j][i-1] > 0.1) {
1584 if (nvert[k][j][i+1] + nvert[k+1][j][i+1] < 0.1) {
1585 dpdc = (p[k][j][i+1] + p[k+1][j][i+1] -
1586 p[k][j][i ] - p[k+1][j][i ]) * 0.5;
1587 }
1588 }
1589 else if ((i == 1 || i==mx-2) && user->
bctype[0]==7 && nvert[k][j][i-1] + nvert[k+1][j][i-1] > 0.1) {
1590 if (nvert[k][j][i+1] + nvert[k+1][j][i+1] < 0.1) {
1591 dpdc = (p[k][j][i+1] + p[k+1][j][i+1] -
1592 p[k][j][i ] - p[k+1][j][i ]) * 0.5;
1593 }
1594 }
1595 else {
1596 dpdc = (p[k][j][i+1] + p[k+1][j][i+1] -
1597 p[k][j][i-1] - p[k+1][j][i-1]) * 0.25;
1598 }
1599
1600 if ((j == my-2 && user->
bctype[2]!=7)|| nvert[k][j+1][i] + nvert[k+1][j+1][i] > 0.1) {
1601 if (nvert[k][j-1][i] + nvert[k+1][j-1][i] < 0.1 && (j!=1 || (j==1 && user->
bctype[2]==7))) {
1602 dpde = (p[k][j ][i] + p[k+1][j ][i] -
1603 p[k][j-1][i] - p[k+1][j-1][i]) * 0.5;
1604 }
1605 }
1606 else if ((j == my-2 || j==1) && user->
bctype[2]==7 && nvert[k][j+1][i] + nvert[k+1][j+1][i] > 0.1) {
1607 if (nvert[k][j-1][i] + nvert[k+1][j-1][i] < 0.1) {
1608 dpde = (p[k][j ][i] + p[k+1][j ][i] -
1609 p[k][j-1][i] - p[k+1][j-1][i]) * 0.5;
1610 }
1611 }
1612 else if ((j == 1 && user->
bctype[2]!=7)|| nvert[k][j-1][i] + nvert[k+1][j-1][i] > 0.1) {
1613 if (nvert[k][j+1][i] + nvert[k+1][j+1][i] < 0.1) {
1614 dpde = (p[k][j+1][i] + p[k+1][j+1][i] -
1615 p[k][j ][i] - p[k+1][j ][i]) * 0.5;
1616 }
1617 }
1618 else if ((j == 1 || j==my-2) && user->
bctype[2]==7 && nvert[k][j-1][i] + nvert[k+1][j-1][i] > 0.1) {
1619 if (nvert[k][j+1][i] + nvert[k+1][j+1][i] < 0.1) {
1620 dpde = (p[k][j+1][i] + p[k+1][j+1][i] -
1621 p[k][j ][i] - p[k+1][j ][i]) * 0.5;
1622 }
1623 }
1624 else {
1625 dpde = (p[k][j+1][i] + p[k+1][j+1][i] -
1626 p[k][j-1][i] - p[k+1][j-1][i]) * 0.25;
1627 }
1628
1629 dpdz = (p[k+1][j][i] - p[k][j][i]);
1630
1631 rhs[k][j][i].
z =0.5 * (rct[k][j][i].
z + rct[k+1][j][i].
z);
1632
1634 (dpdc * (kcsi[k][j][i].
x * kzet[k][j][i].
x +
1635 kcsi[k][j][i].
y * kzet[k][j][i].
y +
1636 kcsi[k][j][i].
z * kzet[k][j][i].
z) +
1637 dpde * (keta[k][j][i].x * kzet[k][j][i].x +
1638 keta[k][j][i].y * kzet[k][j][i].y +
1639 keta[k][j][i].z * kzet[k][j][i].z) +
1640 dpdz * (kzet[k][j][i].
x * kzet[k][j][i].
x +
1641 kzet[k][j][i].
y * kzet[k][j][i].
y +
1642 kzet[k][j][i].
z * kzet[k][j][i].
z)) * kaj[k][j][i];
1643
1644 }
1645 }
1646 }
1647
1648
1649
1650
1651
1652 if (user->
bctype[0]==7 && xs==0){
1653 for (k=lzs; k<lze; k++) {
1654 for (j=lys; j<lye; j++) {
1655 i=xs;
1656
1657 dpdc = p[k][j][i+1] - p[k][j][i];
1658
1659 if ((j==my-2 && user->
bctype[2]!=7)|| nvert[k][j+1][i]+nvert[k][j+1][i+1] > 0.1) {
1660 if (nvert[k][j-1][i] + nvert[k][j-1][i+1] < 0.1 && (j!=1 || (j==1 && user->
bctype[2]==7))) {
1661 dpde = (p[k][j ][i] + p[k][j ][i+1] -
1662 p[k][j-1][i] - p[k][j-1][i+1]) * 0.5;
1663 }
1664 }
1665 else if ((j==my-2 || j==1) && user->
bctype[2]==7 && nvert[k][j+1][i]+nvert[k][j+1][i+1] > 0.1) {
1666 if (nvert[k][j-1][i] + nvert[k][j-1][i+1] < 0.1) {
1667 dpde = (p[k][j ][i] + p[k][j ][i+1] -
1668 p[k][j-1][i] - p[k][j-1][i+1]) * 0.5;
1669 }
1670 }
1671 else if ((j == 1 && user->
bctype[2]!=7) || nvert[k][j-1][i] + nvert[k][j-1][i+1] > 0.1) {
1672 if (nvert[k][j+1][i] + nvert[k][j+1][i+1] < 0.1) {
1673 dpde = (p[k][j+1][i] + p[k][j+1][i+1] -
1674 p[k][j ][i] - p[k][j ][i+1]) * 0.5;
1675 }
1676 }
1677 else if ((j == 1 || j==my-2) && user->
bctype[2]==7 && nvert[k][j-1][i] + nvert[k][j-1][i+1] > 0.1) {
1678 if (nvert[k][j+1][i] + nvert[k][j+1][i+1] < 0.1) {
1679 dpde = (p[k][j+1][i] + p[k][j+1][i+1] -
1680 p[k][j ][i] - p[k][j ][i+1]) * 0.5;
1681 }
1682 }
1683 else {
1684 dpde = (p[k][j+1][i] + p[k][j+1][i+1] -
1685 p[k][j-1][i] - p[k][j-1][i+1]) * 0.25;
1686 }
1687
1688 if ((k == mz-2 && user->
bctype[4]!=7) || nvert[k+1][j][i] + nvert[k+1][j][i+1] > 0.1) {
1689 if (nvert[k-1][j][i] + nvert[k-1][j][i+1] < 0.1 && (k!=1 || (k==1 && user->
bctype[4]==7))) {
1690 dpdz = (p[k ][j][i] + p[k ][j][i+1] -
1691 p[k-1][j][i] - p[k-1][j][i+1]) * 0.5;
1692 }
1693 }
1694 else if ((k == mz-2 || k==1) && user->
bctype[4]==7 && nvert[k+1][j][i] + nvert[k+1][j][i+1] > 0.1) {
1695 if (nvert[k-1][j][i] + nvert[k-1][j][i+1] < 0.1) {
1696 dpdz = (p[k ][j][i] + p[k ][j][i+1] -
1697 p[k-1][j][i] - p[k-1][j][i+1]) * 0.5;
1698 }
1699 }
1700 else if ((k == 1 && user->
bctype[4]!=7)|| nvert[k-1][j][i] + nvert[k-1][j][i+1] > 0.1) {
1701 if (nvert[k+1][j][i] + nvert[k+1][j][i+1] < 0.1) {
1702 dpdz = (p[k+1][j][i] + p[k+1][j][i+1] -
1703 p[k ][j][i] - p[k ][j][i+1]) * 0.5;
1704 }
1705 }
1706 else if ((k == 1 || k==mz-2) && user->
bctype[4]==7 && nvert[k-1][j][i] + nvert[k-1][j][i+1] > 0.1) {
1707 if (nvert[k+1][j][i] + nvert[k+1][j][i+1] < 0.1) {
1708 dpdz = (p[k+1][j][i] + p[k+1][j][i+1] -
1709 p[k ][j][i] - p[k ][j][i+1]) * 0.5;
1710 }
1711 }
1712 else {
1713 dpdz = (p[k+1][j][i] + p[k+1][j][i+1] -
1714 p[k-1][j][i] - p[k-1][j][i+1]) * 0.25;
1715 }
1716
1717 rhs[k][j][i].
x =0.5 * (rct[k][j][i].
x + rct[k][j][i+1].
x);
1719 (dpdc * (icsi[k][j][i].
x * icsi[k][j][i].
x +
1720 icsi[k][j][i].
y * icsi[k][j][i].
y +
1721 icsi[k][j][i].
z * icsi[k][j][i].
z)+
1722 dpde * (ieta[k][j][i].x * icsi[k][j][i].x +
1723 ieta[k][j][i].y * icsi[k][j][i].y +
1724 ieta[k][j][i].z * icsi[k][j][i].z)+
1725 dpdz * (izet[k][j][i].
x * icsi[k][j][i].
x +
1726 izet[k][j][i].
y * icsi[k][j][i].
y +
1727 izet[k][j][i].
z * icsi[k][j][i].
z)) * iaj[k][j][i];
1728 }
1729 }
1730 }
1731
1732
1733 if (user->
bctype[2]==7 && ys==0){
1734 for (k=lzs; k<lze; k++) {
1735 for (i=lxs; i<lxe; i++) {
1736
1737 j=ys;
1738
1739 if ((i == mx-2 && user->
bctype[0]!=7) || nvert[k][j][i+1] + nvert[k][j+1][i+1] > 0.1) {
1740 if (nvert[k][j][i-1] + nvert[k][j+1][i-1] < 0.1 && (i!=1 || (i==1 && user->
bctype[0]==7))) {
1741 dpdc = (p[k][j][i ] + p[k][j+1][i ] -
1742 p[k][j][i-1] - p[k][j+1][i-1]) * 0.5;
1743 }
1744 }
1745 else if ((i == mx-2 || i==1) && user->
bctype[0]==7 && nvert[k][j][i+1] + nvert[k][j+1][i+1] > 0.1) {
1746 if (nvert[k][j][i-1] + nvert[k][j+1][i-1] < 0.1) {
1747 dpdc = (p[k][j][i ] + p[k][j+1][i ] -
1748 p[k][j][i-1] - p[k][j+1][i-1]) * 0.5;
1749 }
1750 }
1751 else if ((i == 1 && user->
bctype[0]!=7)|| nvert[k][j][i-1] + nvert[k][j+1][i-1] > 0.1) {
1752 if (nvert[k][j][i+1] + nvert[k][j+1][i+1] < 0.1) {
1753 dpdc = (p[k][j][i+1] + p[k][j+1][i+1] -
1754 p[k][j][i ] - p[k][j+1][i ]) * 0.5;
1755 }
1756 }
1757 else if ((i == 1 || i==mx-2) && user->
bctype[0]==7 && nvert[k][j][i-1] + nvert[k][j+1][i-1] > 0.1) {
1758 if (nvert[k][j][i+1] + nvert[k][j+1][i+1] < 0.1) {
1759 dpdc = (p[k][j][i+1] + p[k][j+1][i+1] -
1760 p[k][j][i ] - p[k][j+1][i ]) * 0.5;
1761 }
1762 }
1763 else {
1764 dpdc = (p[k][j][i+1] + p[k][j+1][i+1] -
1765 p[k][j][i-1] - p[k][j+1][i-1]) * 0.25;
1766 }
1767
1768 dpde = p[k][j+1][i] - p[k][j][i];
1769
1770 if ((k == mz-2 && user->
bctype[4]!=7)|| nvert[k+1][j][i] + nvert[k+1][j+1][i] > 0.1) {
1771 if (nvert[k-1][j][i] + nvert[k-1][j+1][i] < 0.1 && (k!=1 || (k==1 && user->
bctype[4]==7))) {
1772 dpdz = (p[k ][j][i] + p[k ][j+1][i] -
1773 p[k-1][j][i] - p[k-1][j+1][i]) * 0.5;
1774 }
1775 }
1776 else if ((k == mz-2 || k==1 ) && user->
bctype[4]==7 && nvert[k+1][j][i] + nvert[k+1][j+1][i] > 0.1) {
1777 if (nvert[k-1][j][i] + nvert[k-1][j+1][i] < 0.1) {
1778 dpdz = (p[k ][j][i] + p[k ][j+1][i] -
1779 p[k-1][j][i] - p[k-1][j+1][i]) * 0.5;
1780 }
1781 }
1782 else if ((k == 1 && user->
bctype[4]!=7)|| nvert[k-1][j][i] + nvert[k-1][j+1][i] > 0.1) {
1783 if (nvert[k+1][j][i] + nvert[k+1][j+1][i] < 0.1) {
1784 dpdz = (p[k+1][j][i] + p[k+1][j+1][i] -
1785 p[k ][j][i] - p[k ][j+1][i]) * 0.5;
1786 }
1787 }
1788 else if ((k == 1 || k==mz-2) && user->
bctype[4]==7 && nvert[k-1][j][i] + nvert[k-1][j+1][i] > 0.1) {
1789 if (nvert[k+1][j][i] + nvert[k+1][j+1][i] < 0.1) {
1790 dpdz = (p[k+1][j][i] + p[k+1][j+1][i] -
1791 p[k ][j][i] - p[k ][j+1][i]) * 0.5;
1792 }
1793 }
1794 else {
1795 dpdz = (p[k+1][j][i] + p[k+1][j+1][i] -
1796 p[k-1][j][i] - p[k-1][j+1][i]) * 0.25;
1797 }
1798
1799 rhs[k][j][i].
y =0.5 * (rct[k][j][i].
y + rct[k][j+1][i].
y);
1800
1802 (dpdc * (jcsi[k][j][i].
x * jeta[k][j][i].
x +
1803 jcsi[k][j][i].
y * jeta[k][j][i].
y +
1804 jcsi[k][j][i].
z * jeta[k][j][i].
z)+
1805 dpde * (jeta[k][j][i].x * jeta[k][j][i].x +
1806 jeta[k][j][i].y * jeta[k][j][i].y +
1807 jeta[k][j][i].z * jeta[k][j][i].z)+
1808 dpdz * (jzet[k][j][i].
x * jeta[k][j][i].
x +
1809 jzet[k][j][i].
y * jeta[k][j][i].
y +
1810 jzet[k][j][i].
z * jeta[k][j][i].
z)) * jaj[k][j][i];
1811
1812 }
1813 }
1814 }
1815
1816
1817 if (user->
bctype[4]==7&& zs==0){
1818 for (j=lys; j<lye; j++) {
1819 for (i=lxs; i<lxe; i++) {
1820
1821 k=zs;
1822
1823 if ((i == mx-2 && user->
bctype[0]!=7)|| nvert[k][j][i+1] + nvert[k+1][j][i+1] > 0.1) {
1824 if (nvert[k][j][i-1] + nvert[k+1][j][i-1] < 0.1 && (i!=1 || (i==1 && user->
bctype[0]==7))) {
1825 dpdc = (p[k][j][i ] + p[k+1][j][i ] -
1826 p[k][j][i-1] - p[k+1][j][i-1]) * 0.5;
1827 }
1828 }
1829 else if ((i == mx-2 || i==1) && user->
bctype[0]==7 && nvert[k][j][i+1] + nvert[k+1][j][i+1] > 0.1) {
1830 if (nvert[k][j][i-1] + nvert[k+1][j][i-1] < 0.1) {
1831 dpdc = (p[k][j][i ] + p[k+1][j][i ] -
1832 p[k][j][i-1] - p[k+1][j][i-1]) * 0.5;
1833 }
1834 }
1835 else if ((i == 1 && user->
bctype[0]!=7)|| nvert[k][j][i-1] + nvert[k+1][j][i-1] > 0.1) {
1836 if (nvert[k][j][i+1] + nvert[k+1][j][i+1] < 0.1) {
1837 dpdc = (p[k][j][i+1] + p[k+1][j][i+1] -
1838 p[k][j][i ] - p[k+1][j][i ]) * 0.5;
1839 }
1840 }
1841 else if ((i == 1 || i==mx-2) && user->
bctype[0]==7 && nvert[k][j][i-1] + nvert[k+1][j][i-1] > 0.1) {
1842 if (nvert[k][j][i+1] + nvert[k+1][j][i+1] < 0.1) {
1843 dpdc = (p[k][j][i+1] + p[k+1][j][i+1] -
1844 p[k][j][i ] - p[k+1][j][i ]) * 0.5;
1845 }
1846 }
1847 else {
1848 dpdc = (p[k][j][i+1] + p[k+1][j][i+1] -
1849 p[k][j][i-1] - p[k+1][j][i-1]) * 0.25;
1850 }
1851
1852 if ((j == my-2 && user->
bctype[2]!=7)|| nvert[k][j+1][i] + nvert[k+1][j+1][i] > 0.1) {
1853 if (nvert[k][j-1][i] + nvert[k+1][j-1][i] < 0.1 && (j!=1 || (j==1 && user->
bctype[2]==7))) {
1854 dpde = (p[k][j ][i] + p[k+1][j ][i] -
1855 p[k][j-1][i] - p[k+1][j-1][i]) * 0.5;
1856 }
1857 }
1858 else if ((j == my-2 || j==1) && user->
bctype[2]==7 && nvert[k][j+1][i] + nvert[k+1][j+1][i] > 0.1) {
1859 if (nvert[k][j-1][i] + nvert[k+1][j-1][i] < 0.1) {
1860 dpde = (p[k][j ][i] + p[k+1][j ][i] -
1861 p[k][j-1][i] - p[k+1][j-1][i]) * 0.5;
1862 }
1863 }
1864 else if ((j == 1 && user->
bctype[2]!=7)|| nvert[k][j-1][i] + nvert[k+1][j-1][i] > 0.1) {
1865 if (nvert[k][j+1][i] + nvert[k+1][j+1][i] < 0.1) {
1866 dpde = (p[k][j+1][i] + p[k+1][j+1][i] -
1867 p[k][j ][i] - p[k+1][j ][i]) * 0.5;
1868 }
1869 }
1870 else if ((j == 1 || j==my-2) && user->
bctype[2]==7 && nvert[k][j-1][i] + nvert[k+1][j-1][i] > 0.1) {
1871 if (nvert[k][j+1][i] + nvert[k+1][j+1][i] < 0.1) {
1872 dpde = (p[k][j+1][i] + p[k+1][j+1][i] -
1873 p[k][j ][i] - p[k+1][j ][i]) * 0.5;
1874 }
1875 }
1876 else {
1877 dpde = (p[k][j+1][i] + p[k+1][j+1][i] -
1878 p[k][j-1][i] - p[k+1][j-1][i]) * 0.25;
1879 }
1880
1881 dpdz = (p[k+1][j][i] - p[k][j][i]);
1882
1883 rhs[k][j][i].
z =0.5 * (rct[k][j][i].
z + rct[k+1][j][i].
z);
1884
1886 (dpdc * (kcsi[k][j][i].
x * kzet[k][j][i].
x +
1887 kcsi[k][j][i].
y * kzet[k][j][i].
y +
1888 kcsi[k][j][i].
z * kzet[k][j][i].
z)+
1889 dpde * (keta[k][j][i].x * kzet[k][j][i].x +
1890 keta[k][j][i].y * kzet[k][j][i].y +
1891 keta[k][j][i].z * kzet[k][j][i].z)+
1892 dpdz * (kzet[k][j][i].
x * kzet[k][j][i].
x +
1893 kzet[k][j][i].
y * kzet[k][j][i].
y +
1894 kzet[k][j][i].
z * kzet[k][j][i].
z)) * kaj[k][j][i];
1895
1896 }
1897 }
1898 }
1899
1900 ierr = DMDAVecRestoreArray(fda, Rct, &rct); CHKERRQ(ierr);
1901
1903 PetscInt TwoD = simCtx->
TwoD;
1904
1906
1907
1908 for (k=lzs; k<lze; k++) {
1909 for (j=lys; j<lye; j++) {
1910 for (i=lxs; i<lxe; i++) {
1911 if (TwoD==1)
1913 else if (TwoD==2)
1915 else if (TwoD==3)
1917
1918 if (nvert[k][j][i]>0.1) {
1922 }
1923 if (nvert[k][j][i+1]>0.1) {
1925 }
1926 if (nvert[k][j+1][i]>0.1) {
1928 }
1929 if (nvert[k+1][j][i]>0.1) {
1931 }
1932 }
1933 }
1934 }
1936
1937
1938
1939
1940
1941 DMDAVecRestoreArray(fda, Rhs, &rhs);
1943
1944 DMDAVecRestoreArray(fda, user->
lCsi, &csi);
1945 DMDAVecRestoreArray(fda, user->
lEta, &eta);
1946 DMDAVecRestoreArray(fda, user->
lZet, &zet);
1947 DMDAVecRestoreArray(da, user->
lAj, &aj);
1949
1950 DMDAVecRestoreArray(fda, user->
lICsi, &icsi);
1951 DMDAVecRestoreArray(fda, user->
lIEta, &ieta);
1952 DMDAVecRestoreArray(fda, user->
lIZet, &izet);
1953 DMDAVecRestoreArray(da, user->
lIAj, &iaj);
1955
1956 DMDAVecRestoreArray(fda, user->
lJCsi, &jcsi);
1957 DMDAVecRestoreArray(fda, user->
lJEta, &jeta);
1958 DMDAVecRestoreArray(fda, user->
lJZet, &jzet);
1959 DMDAVecRestoreArray(da, user->
lJAj, &jaj);
1961
1962 DMDAVecRestoreArray(fda, user->
lKCsi, &kcsi);
1963 DMDAVecRestoreArray(fda, user->
lKEta, &keta);
1964 DMDAVecRestoreArray(fda, user->
lKZet, &kzet);
1965 DMDAVecRestoreArray(da, user->
lKAj, &kaj);
1967
1968 DMDAVecRestoreArray(da, user->
lP, &p);
1970
1971 DMDAVecRestoreArray(da, user->
lNvert, &nvert);
1973
1975
1976
1977 ierr = VecDestroy(&Conv); CHKERRQ(ierr);
1978 ierr = VecDestroy(&Visc); CHKERRQ(ierr);
1979 ierr = VecDestroy(&Rc); CHKERRQ(ierr);
1980 ierr = VecDestroy(&Rct); CHKERRQ(ierr);
1982
1985 PetscFunctionReturn(0);
1986}
#define LOCAL
Logging scope definitions for controlling message output.
PetscErrorCode Viscous(UserCtx *user, Vec Ucont, Vec Ucat, Vec Visc)
PetscErrorCode Convection(UserCtx *user, Vec Ucont, Vec Ucat, Vec Conv)
PetscErrorCode Contra2Cart(UserCtx *user)
Reconstructs Cartesian velocity (Ucat) at cell centers from contravariant velocity (Ucont) defined on...