112 DMBoundaryType xperiod = (simCtx->
i_periodic) ? DM_BOUNDARY_PERIODIC : DM_BOUNDARY_NONE;
113 DMBoundaryType yperiod = (simCtx->
j_periodic) ? DM_BOUNDARY_PERIODIC : DM_BOUNDARY_NONE;
114 DMBoundaryType zperiod = (simCtx->
k_periodic) ? DM_BOUNDARY_PERIODIC : DM_BOUNDARY_NONE;
117 PetscInt *lx = NULL, *ly = NULL, *lz = NULL;
120 PetscFunctionBeginUser;
128 DMDAGetInfo(coarse_user->
da, NULL, NULL, NULL, NULL, &m, &n, &p, NULL, NULL, NULL, NULL, NULL, NULL);
132 PetscInt *lx_contrib, *ly_contrib, *lz_contrib;
133 ierr = PetscMalloc3(m, &lx_contrib, n, &ly_contrib, p, &lz_contrib); CHKERRQ(ierr);
134 ierr = PetscMemzero(lx_contrib, m *
sizeof(PetscInt)); CHKERRQ(ierr);
135 ierr = PetscMemzero(ly_contrib, n *
sizeof(PetscInt)); CHKERRQ(ierr);
136 ierr = PetscMemzero(lz_contrib, p *
sizeof(PetscInt)); CHKERRQ(ierr);
139 DMDAGetLocalInfo(coarse_user->
da, &info);
140 PetscInt xs = info.xs, xe = info.xs + info.xm, mx = info.mx;
141 PetscInt ys = info.ys, ye = info.ys + info.ym, my = info.my;
142 PetscInt zs = info.zs, ze = info.zs + info.zm, mz = info.mz;
145 ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank); CHKERRQ(ierr);
146 PetscInt proc_i = rank % m;
147 PetscInt proc_j = (rank / m) % n;
148 PetscInt proc_k = rank / (m * n);
151 if (user->
isc) lx_contrib[proc_i] = (xe - xs);
153 if (m == 1) lx_contrib[0] = user->
IM + 1;
154 else if (xs == 0) lx_contrib[0] = 2 * xe - 1;
155 else if (xe == mx) lx_contrib[proc_i] = user->
IM + 1 - (2 * xs - 1);
156 else lx_contrib[proc_i] = (xe - xs) * 2;
160 if (user->
jsc) ly_contrib[proc_j] = (ye - ys);
162 if (n == 1) ly_contrib[0] = user->
JM + 1;
163 else if (ys == 0) ly_contrib[0] = 2 * ye - 1;
164 else if (ye == my) ly_contrib[proc_j] = user->
JM + 1 - (2 * ys - 1);
165 else ly_contrib[proc_j] = (ye - ys) * 2;
169 if (user->
ksc) lz_contrib[proc_k] = (ze - zs);
171 if (p == 1) lz_contrib[0] = user->
KM + 1;
172 else if (zs == 0) lz_contrib[0] = 2 * ze - 1;
173 else if (ze == mz) lz_contrib[proc_k] = user->
KM + 1 - (2 * zs - 1);
174 else lz_contrib[proc_k] = (ze - zs) * 2;
176 LOG_ALLOW_SYNC(
LOCAL,
LOG_VERBOSE,
"Rank %d: Calculated this rank's node contribution to fine grid: lx=%d, ly=%d, lz=%d\n", simCtx->
rank, lx_contrib[proc_i], ly_contrib[proc_j], lz_contrib[proc_k]);
179 ierr = PetscMalloc3(m, &lx, n, &ly, p, &lz); CHKERRQ(ierr);
180 ierr = MPI_Allreduce(lx_contrib, lx, m, MPIU_INT, MPI_MAX, PETSC_COMM_WORLD); CHKERRQ(ierr);
181 ierr = MPI_Allreduce(ly_contrib, ly, n, MPIU_INT, MPI_MAX, PETSC_COMM_WORLD); CHKERRQ(ierr);
182 ierr = MPI_Allreduce(lz_contrib, lz, p, MPIU_INT, MPI_MAX, PETSC_COMM_WORLD); CHKERRQ(ierr);
184 ierr = PetscFree3(lx_contrib, ly_contrib, lz_contrib); CHKERRQ(ierr);
199 m = n = p = PETSC_DECIDE;
207 ierr = DMDACreate3d(PETSC_COMM_WORLD, xperiod, yperiod, zperiod, DMDA_STENCIL_BOX,
208 user->
IM + 1, user->
JM + 1, user->
KM + 1,
210 1, stencil_width, lx, ly, lz, &user->
da); CHKERRQ(ierr);
213 ierr = PetscFree3(lx, ly, lz); CHKERRQ(ierr);
217 ierr = DMSetUp(user->
da); CHKERRQ(ierr);
218 ierr = DMGetCoordinateDM(user->
da, &user->
fda); CHKERRQ(ierr);
219 ierr = DMDASetUniformCoordinates(user->
da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0); CHKERRQ(ierr);
220 ierr = DMDAGetLocalInfo(user->
da, &user->
info); CHKERRQ(ierr);
225 PetscFunctionReturn(0);
375 PetscFunctionBeginUser;
385 FILE *grid_file_handle = NULL;
387 if (simCtx->
rank == 0) {
388 grid_file_handle = fopen(simCtx->
grid_file,
"r");
389 if (!grid_file_handle) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN,
"Cannot open file: %s", simCtx->
grid_file);
394 char dummy_buffer[2048];
395 for (PetscInt s = 0; s < headerLines; ++s) {
396 if (!fgets(dummy_buffer,
sizeof(dummy_buffer), grid_file_handle)) {
397 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_READ,
"Unexpected EOF while skipping grid header");
408 if (simCtx->
rank == 0) {
409 fclose(grid_file_handle);
417 ierr = DMGetCoordinates(user->
da, &gCoor); CHKERRQ(ierr);
418 ierr = DMGetCoordinatesLocal(user->
da, &lCoor); CHKERRQ(ierr);
421 ierr = DMLocalToGlobalBegin(user->
fda, lCoor, INSERT_VALUES, gCoor); CHKERRQ(ierr);
422 ierr = DMLocalToGlobalEnd(user->
fda, lCoor, INSERT_VALUES, gCoor); CHKERRQ(ierr);
424 ierr = DMGlobalToLocalBegin(user->
fda, gCoor, INSERT_VALUES, lCoor); CHKERRQ(ierr);
425 ierr = DMGlobalToLocalEnd(user->
fda, gCoor, INSERT_VALUES, lCoor); CHKERRQ(ierr);
429 PetscFunctionReturn(0);
526 PetscMPIInt rank = simCtx->
rank;
527 PetscInt block_index = user->
_this;
528 PetscInt IM = user->
IM, JM = user->
JM, KM = user->
KM;
532 PetscReal *gc = NULL;
534 PetscFunctionBeginUser;
539 simCtx->
rank, block_index);
543 PetscInt n_nodes = (IM) * (JM) * (KM);
544 ierr = PetscMalloc1(3 * n_nodes, &gc); CHKERRQ(ierr);
548 if (!fd) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN,
"Recieved a NULL file handle.\n");
551 for (PetscInt k = 0; k < KM; k++) {
552 for (PetscInt j = 0; j < JM; j++) {
553 for (PetscInt i = 0; i < IM; i++) {
554 PetscInt base_index = 3 * ((k * (JM) + j) * (IM) + i);
555 if (fscanf(fd,
"%le %le %le\n", &gc[base_index], &gc[base_index + 1], &gc[base_index + 2]) != 3) {
556 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_READ,
"Error reading coordinates for node (i,j,k)=(%d,%d,%d) in block %d", i, j, k, block_index);
565 ierr = MPI_Bcast(gc, 3 * n_nodes, MPIU_REAL, 0, PETSC_COMM_WORLD); CHKERRQ(ierr);
568 ierr = DMDAGetLocalInfo(user->
da, &info); CHKERRQ(ierr);
569 ierr = DMGetCoordinatesLocal(user->
da, &lCoor); CHKERRQ(ierr);
570 ierr = DMDAVecGetArray(user->
fda, lCoor, &coor); CHKERRQ(ierr);
572 for (PetscInt k = info.zs; k < info.zs + info.zm; k++) {
573 for (PetscInt j = info.ys; j < info.ys + info.ym; j++) {
574 for (PetscInt i = info.xs; i < info.xs + info.xm; i++) {
575 if(k< KM && j < JM && i < IM){
576 PetscInt base_idx = 3 * ((k * (JM) + j) * (IM) + i);
577 coor[k][j][i].
x = gc[base_idx];
578 coor[k][j][i].
y = gc[base_idx + 1];
579 coor[k][j][i].
z = gc[base_idx + 2];
586 ierr = DMDAVecRestoreArray(user->
fda, lCoor, &coor); CHKERRQ(ierr);
587 ierr = PetscFree(gc); CHKERRQ(ierr);
591 PetscFunctionReturn(0);
603 Vec c_lCoor, f_lCoor;
606 DMDALocalInfo c_info;
609 PetscFunctionBeginUser;
616 ierr = DMDAGetLocalInfo(coarse_user->
da, &c_info); CHKERRQ(ierr);
618 ierr = DMGetCoordinatesLocal(coarse_user->
da, &c_lCoor); CHKERRQ(ierr);
619 ierr = DMGetCoordinatesLocal(fine_user->
da, &f_lCoor); CHKERRQ(ierr);
621 ierr = DMDAVecGetArray(coarse_user->
fda, c_lCoor, &c_coor); CHKERRQ(ierr);
622 ierr = DMDAVecGetArrayRead(fine_user->
fda, f_lCoor, &f_coor); CHKERRQ(ierr);
625 PetscInt xs = c_info.xs, xe = c_info.xs + c_info.xm;
626 PetscInt ys = c_info.ys, ye = c_info.ys + c_info.ym;
627 PetscInt zs = c_info.zs, ze = c_info.zs + c_info.zm;
630 PetscInt mx = c_info.mx, my = c_info.my, mz = c_info.mz;
639 for (PetscInt k = zs; k < ze; k++) {
640 for (PetscInt j = ys; j < ye; j++) {
641 for (PetscInt i = xs; i < xe; i++) {
644 ih = coarse_user->
isc ? i : 2 * i;
645 jh = coarse_user->
jsc ? j : 2 * j;
646 kh = coarse_user->
ksc ? k : 2 * k;
650 c_coor[k][j][i] = f_coor[kh][jh][ih];
655 ierr = DMDAVecRestoreArray(coarse_user->
fda, c_lCoor, &c_coor); CHKERRQ(ierr);
656 ierr = DMDAVecRestoreArrayRead(fine_user->
fda, f_lCoor, &f_coor); CHKERRQ(ierr);
660 ierr = DMGetCoordinates(coarse_user->
da, &c_gCoor); CHKERRQ(ierr);
661 ierr = DMLocalToGlobalBegin(coarse_user->
fda, c_lCoor, INSERT_VALUES, c_gCoor); CHKERRQ(ierr);
662 ierr = DMLocalToGlobalEnd(coarse_user->
fda, c_lCoor, INSERT_VALUES, c_gCoor); CHKERRQ(ierr);
663 ierr = DMGlobalToLocalBegin(coarse_user->
fda, c_gCoor, INSERT_VALUES, c_lCoor); CHKERRQ(ierr);
664 ierr = DMGlobalToLocalEnd(coarse_user->
fda, c_gCoor, INSERT_VALUES, c_lCoor); CHKERRQ(ierr);
668 PetscFunctionReturn(0);
1030 Cmpnts *face_center, PetscReal *face_area)
1032 PetscErrorCode ierr;
1038 PetscReal local_sum[3] = {0.0, 0.0, 0.0};
1039 PetscReal localAreaSum = 0.0;
1040 PetscCount local_n_points = 0;
1045 PetscReal global_sum[3] = {0.0, 0.0, 0.0};
1046 PetscReal globalAreaSum = 0.0;
1047 PetscCount global_n_points = 0;
1055 PetscInt xs = info.xs, xe = info.xs + info.xm;
1056 PetscInt ys = info.ys, ye = info.ys + info.ym;
1057 PetscInt zs = info.zs, ze = info.zs + info.zm;
1060 PetscInt mx = info.mx, my = info.my, mz = info.mz;
1061 PetscInt IM = user->
IM;
1062 PetscInt JM = user->
JM;
1063 PetscInt KM = user->
KM;
1069 PetscInt lxs = xs;
if(xs == 0) lxs = xs + 1;
1070 PetscInt lxe = xe;
if(xe == mx) lxe = xe - 1;
1071 PetscInt lys = ys;
if(ys == 0) lys = ys + 1;
1072 PetscInt lye = ye;
if(ye == my) lye = ye - 1;
1073 PetscInt lzs = zs;
if(zs == 0) lzs = zs + 1;
1074 PetscInt lze = ze;
if(ze == mz) lze = ze - 1;
1080 PetscInt i_max = (xe == mx) ? mx - 1 : xe;
1081 PetscInt j_max = (ye == my) ? my - 1 : ye;
1082 PetscInt k_max = (ze == mz) ? mz - 1 : ze;
1089 Cmpnts ***csi, ***eta, ***zet;
1092 PetscFunctionBeginUser;
1098 PetscBool owns_face = PETSC_FALSE;
1104 ierr = DMGetCoordinatesLocal(user->
da, &lCoor); CHKERRQ(ierr);
1105 ierr = DMDAVecGetArrayRead(user->
fda, lCoor, &coor); CHKERRQ(ierr);
1106 ierr = DMDAVecGetArrayRead(user->
da, user->
lNvert, &nvert); CHKERRQ(ierr);
1107 ierr = DMDAVecGetArrayRead(user->
fda, user->
lCsi, &csi); CHKERRQ(ierr);
1108 ierr = DMDAVecGetArrayRead(user->
fda, user->
lEta, &eta); CHKERRQ(ierr);
1109 ierr = DMDAVecGetArrayRead(user->
fda, user->
lZet, &zet); CHKERRQ(ierr);
1126 for (PetscInt k = zs; k < k_max; k++) {
1127 for (PetscInt j = ys; j < j_max; j++) {
1129 local_sum[0] += coor[k][j][i].
x;
1130 local_sum[1] += coor[k][j][i].
y;
1131 local_sum[2] += coor[k][j][i].
z;
1139 for (PetscInt k = lzs; k < lze; k++) {
1140 for (PetscInt j = lys; j < lye; j++) {
1144 if (nvert[k][j][i+1] < 0.1) {
1147 localAreaSum += sqrt(csi[k][j][i].x * csi[k][j][i].x +
1148 csi[k][j][i].y * csi[k][j][i].y +
1149 csi[k][j][i].z * csi[k][j][i].z);
1161 PetscInt i = mx - 2;
1164 for (PetscInt k = zs; k < k_max; k++) {
1165 for (PetscInt j = ys; j < j_max; j++) {
1166 local_sum[0] += coor[k][j][i].
x;
1167 local_sum[1] += coor[k][j][i].
y;
1168 local_sum[2] += coor[k][j][i].
z;
1174 for (PetscInt k = lzs; k < lze; k++) {
1175 for (PetscInt j = lys; j < lye; j++) {
1179 if (nvert[k][j][i] < 0.1) {
1181 localAreaSum += sqrt(csi[k][j][i].x * csi[k][j][i].x +
1182 csi[k][j][i].y * csi[k][j][i].y +
1183 csi[k][j][i].z * csi[k][j][i].z);
1199 for (PetscInt k = zs; k < k_max; k++) {
1200 for (PetscInt i = xs; i < i_max; i++) {
1201 local_sum[0] += coor[k][j][i].
x;
1202 local_sum[1] += coor[k][j][i].
y;
1203 local_sum[2] += coor[k][j][i].
z;
1210 for (PetscInt k = lzs; k < lze; k++) {
1211 for (PetscInt i = lxs; i < lxe; i++) {
1213 if (nvert[k][j+1][i] < 0.1) {
1215 localAreaSum += sqrt(eta[k][j][i].x * eta[k][j][i].x +
1216 eta[k][j][i].y * eta[k][j][i].y +
1217 eta[k][j][i].z * eta[k][j][i].z);
1229 PetscInt j = my - 2;
1232 for (PetscInt k = zs; k < k_max; k++) {
1233 for (PetscInt i = xs; i < i_max; i++) {
1234 local_sum[0] += coor[k][j][i].
x;
1235 local_sum[1] += coor[k][j][i].
y;
1236 local_sum[2] += coor[k][j][i].
z;
1242 for (PetscInt k = lzs; k < lze; k++) {
1243 for (PetscInt i = lxs; i < lxe; i++) {
1245 if (nvert[k][j][i] < 0.1) {
1247 localAreaSum += sqrt(eta[k][j][i].x * eta[k][j][i].x +
1248 eta[k][j][i].y * eta[k][j][i].y +
1249 eta[k][j][i].z * eta[k][j][i].z);
1265 for (PetscInt j = ys; j < j_max; j++) {
1266 for (PetscInt i = xs; i < i_max; i++) {
1267 local_sum[0] += coor[k][j][i].
x;
1268 local_sum[1] += coor[k][j][i].
y;
1269 local_sum[2] += coor[k][j][i].
z;
1276 for (PetscInt j = lys; j < lye; j++) {
1277 for (PetscInt i = lxs; i < lxe; i++) {
1280 if (nvert[k+1][j][i] < 0.1) {
1282 localAreaSum += sqrt(zet[k][j][i].x * zet[k][j][i].x +
1283 zet[k][j][i].y * zet[k][j][i].y +
1284 zet[k][j][i].z * zet[k][j][i].z);
1296 PetscInt k = mz - 2;
1300 for (PetscInt j = ys; j < j_max; j++) {
1301 for (PetscInt i = xs; i < i_max; i++) {
1302 local_sum[0] += coor[k][j][i].
x;
1303 local_sum[1] += coor[k][j][i].
y;
1304 local_sum[2] += coor[k][j][i].
z;
1311 for (PetscInt j = lys; j < lye; j++) {
1312 for (PetscInt i = lxs; i < lxe; i++) {
1315 if (nvert[k][j][i] < 0.1) {
1317 localAreaSum += sqrt(zet[k][j][i].x * zet[k][j][i].x +
1318 zet[k][j][i].y * zet[k][j][i].y +
1319 zet[k][j][i].z * zet[k][j][i].z);
1327 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE,
1328 "Unknown face_id %d in CalculateFaceCenterAndArea", face_id);
1334 ierr = DMDAVecRestoreArrayRead(user->
fda, lCoor, &coor); CHKERRQ(ierr);
1335 ierr = DMDAVecRestoreArrayRead(user->
da, user->
lNvert, &nvert); CHKERRQ(ierr);
1336 ierr = DMDAVecRestoreArrayRead(user->
fda, user->
lCsi, &csi); CHKERRQ(ierr);
1337 ierr = DMDAVecRestoreArrayRead(user->
fda, user->
lEta, &eta); CHKERRQ(ierr);
1338 ierr = DMDAVecRestoreArrayRead(user->
fda, user->
lZet, &zet); CHKERRQ(ierr);
1344 ierr = MPI_Allreduce(local_sum, global_sum, 3, MPI_DOUBLE, MPI_SUM,
1345 PETSC_COMM_WORLD); CHKERRQ(ierr);
1348 ierr = MPI_Allreduce(&local_n_points, &global_n_points, 1, MPI_COUNT, MPI_SUM,
1349 PETSC_COMM_WORLD); CHKERRQ(ierr);
1352 ierr = MPI_Allreduce(&localAreaSum, &globalAreaSum, 1, MPI_DOUBLE, MPI_SUM,
1353 PETSC_COMM_WORLD); CHKERRQ(ierr);
1358 if (global_n_points > 0) {
1359 face_center->
x = global_sum[0] / global_n_points;
1360 face_center->
y = global_sum[1] / global_n_points;
1361 face_center->
z = global_sum[2] / global_n_points;
1363 "Calculated center for Face %s: (x=%.4f, y=%.4f, z=%.4f) from %lld nodes\n",
1365 face_center->
x, face_center->
y, face_center->
z,
1366 (
long long)global_n_points);
1370 "WARNING: Face %s identified but no grid points found. Center not calculated.\n",
1372 face_center->
x = face_center->
y = face_center->
z = 0.0;
1378 *face_area = globalAreaSum;
1380 "Calculated area for Face %s: Area=%.6f\n",
1383 PetscFunctionReturn(0);