256 const PetscReal V0 = 1.0;
257 const PetscReal rho = 1.0;
258 const PetscReal p0 = 0.0;
261 const PetscReal nu = (simCtx->
ren > 0) ? (1.0 / simCtx->
ren) : 0.0;
263 const PetscReal k = 1.0;
264 const PetscReal t = simCtx->
ti;
266 LOG_ALLOW(
GLOBAL,
LOG_TRACE,
"TGV Setup: t = %.4f, V0* = %.4f, rho* = %.4f, k = %.4f, p0* = %4.f, nu = %.6f.\n",simCtx->
ti,V0,rho,k,p0,nu);
268 const PetscReal vel_decay = exp(-2.0 * nu * k * k * t);
269 const PetscReal prs_decay = exp(-4.0 * nu * k * k * t);
271 PetscFunctionBeginUser;
273 for (PetscInt bi = 0; bi < simCtx->
block_number; bi++) {
274 UserCtx* user = &user_finest[bi];
275 DMDALocalInfo info = user->
info;
276 PetscInt xs = info.xs, xe = info.xs + info.xm;
277 PetscInt ys = info.ys, ye = info.ys + info.ym;
278 PetscInt zs = info.zs, ze = info.zs + info.zm;
279 PetscInt mx = info.mx, my = info.my, mz = info.mz;
282 const Cmpnts ***cent, ***cent_x, ***cent_y, ***cent_z;
286 PetscInt lxs = (xs == 0) ? xs + 1 : xs, lxe = (xe == mx) ? xe - 1 : xe;
287 PetscInt lys = (ys == 0) ? ys + 1 : ys, lye = (ye == my) ? ye - 1 : ye;
288 PetscInt lzs = (zs == 0) ? zs + 1 : zs, lze = (ze == mz) ? ze - 1 : ze;
291 ierr = DMDAVecGetArray(user->
fda, user->
Ucat, &ucat); CHKERRQ(ierr);
292 ierr = DMDAVecGetArray(user->
da, user->
P, &p); CHKERRQ(ierr);
293 ierr = DMDAVecGetArray(user->
fda, user->
Bcs.
Ubcs, &ubcs); CHKERRQ(ierr);
294 ierr = DMDAVecGetArrayRead(user->
fda, user->
Cent, ¢); CHKERRQ(ierr);
295 ierr = DMDAVecGetArrayRead(user->
fda, user->
Centx, ¢_x); CHKERRQ(ierr);
296 ierr = DMDAVecGetArrayRead(user->
fda, user->
Centy, ¢_y); CHKERRQ(ierr);
297 ierr = DMDAVecGetArrayRead(user->
fda, user->
Centz, ¢_z); CHKERRQ(ierr);
300 for (PetscInt k_cell = lzs; k_cell < lze; k_cell++) {
301 for (PetscInt j_cell = lys; j_cell < lye; j_cell++) {
302 for (PetscInt i_cell = lxs; i_cell < lxe; i_cell++) {
303 const PetscReal cx = cent[k_cell][j_cell][i_cell].
x, cy = cent[k_cell][j_cell][i_cell].
y, cz = cent[k_cell][j_cell][i_cell].
z;
304 ucat[k_cell][j_cell][i_cell].
x = V0 * sin(k*cx) * cos(k*cy) * cos(k*cz) * vel_decay;
305 ucat[k_cell][j_cell][i_cell].
y = -V0 * cos(k*cx) * sin(k*cy) * cos(k*cz) * vel_decay;
306 ucat[k_cell][j_cell][i_cell].
z = 0.0;
313 for (PetscInt k_cell = lzs; k_cell < lze; k_cell++) {
314 for (PetscInt j_cell = lys; j_cell < lye; j_cell++) {
315 for (PetscInt i_cell = lxs; i_cell < lxe; i_cell++) {
316 const PetscReal cx = cent[k_cell][j_cell][i_cell].
x, cy = cent[k_cell][j_cell][i_cell].
y;
317 p[k_cell][j_cell][i_cell] = p0 + (rho * V0 * V0 / 4.0) * (cos(2*k*cx) + cos(2*k*cy)) * prs_decay;
324 if (xs == 0)
for (PetscInt k=zs; k<ze; k++)
for (PetscInt j=ys; j<ye; j++) {
325 const PetscReal fcx=cent_x[k][j][xs].
x, fcy=cent_x[k][j][xs].
y, fcz=cent_x[k][j][xs].
z;
326 ubcs[k][j][xs].
x = V0*sin(k*fcx)*cos(k*fcy)*cos(k*fcz)*vel_decay; ubcs[k][j][xs].
y = -V0*cos(k*fcx)*sin(k*fcy)*cos(k*fcz)*vel_decay; ubcs[k][j][xs].
z = 0.0;
328 if (xe == mx)
for (PetscInt k=zs; k<ze; k++)
for (PetscInt j=ys; j<ye; j++) {
329 const PetscReal fcx=cent_x[k][j][xe-1].
x, fcy=cent_x[k][j][xe-1].
y, fcz=cent_x[k][j][xe-1].
z;
330 ubcs[k][j][xe-1].
x = V0*sin(k*fcx)*cos(k*fcy)*cos(k*fcz)*vel_decay; ubcs[k][j][xe-1].
y = -V0*cos(k*fcx)*sin(k*fcy)*cos(k*fcz)*vel_decay; ubcs[k][j][xe-1].
z = 0.0;
332 if (ys == 0)
for (PetscInt k=zs; k<ze; k++)
for (PetscInt i=xs; i<xe; i++) {
333 const PetscReal fcx=cent_y[k][ys][i].
x, fcy=cent_y[k][ys][i].
y, fcz=cent_y[k][ys][i].
z;
334 ubcs[k][ys][i].
x = V0*sin(k*fcx)*cos(k*fcy)*cos(k*fcz)*vel_decay; ubcs[k][ys][i].
y = -V0*cos(k*fcx)*sin(k*fcy)*cos(k*fcz)*vel_decay; ubcs[k][ys][i].
z = 0.0;
336 if (ye == my)
for (PetscInt k=zs; k<ze; k++)
for (PetscInt i=xs; i<xe; i++) {
337 const PetscReal fcx=cent_y[k][ye-1][i].
x, fcy=cent_y[k][ye-1][i].
y, fcz=cent_y[k][ye-1][i].
z;
338 ubcs[k][ye-1][i].
x = V0*sin(k*fcx)*cos(k*fcy)*cos(k*fcz)*vel_decay; ubcs[k][ye-1][i].
y = -V0*cos(k*fcx)*sin(k*fcy)*cos(k*fcz)*vel_decay; ubcs[k][ye-1][i].
z = 0.0;
340 if (zs == 0)
for (PetscInt j=ys; j<ye; j++)
for (PetscInt i=xs; i<xe; i++) {
341 const PetscReal fcx=cent_z[zs][j][i].
x, fcy=cent_z[zs][j][i].
y, fcz=cent_z[zs][j][i].
z;
342 ubcs[zs][j][i].
x = V0*sin(k*fcx)*cos(k*fcy)*cos(k*fcz)*vel_decay; ubcs[zs][j][i].
y = -V0*cos(k*fcx)*sin(k*fcy)*cos(k*fcz)*vel_decay; ubcs[zs][j][i].
z = 0.0;
344 if (ze == mz)
for (PetscInt j=ys; j<ye; j++)
for (PetscInt i=xs; i<xe; i++) {
345 const PetscReal fcx=cent_z[ze-1][j][i].
x, fcy=cent_z[ze-1][j][i].
y, fcz=cent_z[ze-1][j][i].
z;
346 ubcs[ze-1][j][i].
x = V0*sin(k*fcx)*cos(k*fcy)*cos(k*fcz)*vel_decay; ubcs[ze-1][j][i].
y = -V0*cos(k*fcx)*sin(k*fcy)*cos(k*fcz)*vel_decay; ubcs[ze-1][j][i].
z = 0.0;
350 if (xs == 0)
for (PetscInt k=lzs; k<lze; k++)
for (PetscInt j=lys; j<lye; j++) p[k][j][xs] = p[k][j][xs+1];
351 if (xe == mx)
for (PetscInt k=lzs; k<lze; k++)
for (PetscInt j=lys; j<lye; j++) p[k][j][xe-1] = p[k][j][xe-2];
353 if (ys == 0)
for (PetscInt k=lzs; k<lze; k++)
for (PetscInt i=lxs; i<lxe; i++) p[k][ys][i] = p[k][ys+1][i];
354 if (ye == my)
for (PetscInt k=lzs; k<lze; k++)
for (PetscInt i=lxs; i<lxe; i++) p[k][ye-1][i] = p[k][ye-2][i];
356 if (zs == 0)
for (PetscInt j=lys; j<lye; j++)
for (PetscInt i=lxs; i<lxe; i++) p[zs][j][i] = p[zs+1][j][i];
357 if (ze == mz)
for (PetscInt j=lys; j<lye; j++)
for (PetscInt i=lxs; i<lxe; i++) p[ze-1][j][i] = p[ze-2][j][i];
360 ierr = DMDAVecRestoreArray(user->
fda, user->
Ucat, &ucat); CHKERRQ(ierr);
361 ierr = DMDAVecRestoreArray(user->
da, user->
P, &p); CHKERRQ(ierr);
362 ierr = DMDAVecRestoreArray(user->
fda, user->
Bcs.
Ubcs, &ubcs); CHKERRQ(ierr);
363 ierr = DMDAVecRestoreArrayRead(user->
fda, user->
Cent, ¢); CHKERRQ(ierr);
364 ierr = DMDAVecRestoreArrayRead(user->
fda, user->
Centx, ¢_x); CHKERRQ(ierr);
365 ierr = DMDAVecRestoreArrayRead(user->
fda, user->
Centy, ¢_y); CHKERRQ(ierr);
366 ierr = DMDAVecRestoreArrayRead(user->
fda, user->
Centz, ¢_z); CHKERRQ(ierr);
382 PetscFunctionReturn(0);