static char help[] = "Solves the incompressible, variable viscosity stokes equation in 3d using Q1Q1 elements, \n\ stabilized with Bochev's polynomial projection method. Note that implementation here assumes \n\ all boundaries are free-slip, i.e. zero normal flow and zero tangential stress \n\ -mx : number elements in x-direction \n\ -my : number elements in y-direction \n\ -mz : number elements in z-direction \n\ -stokes_ksp_monitor_blocks : active monitor for each component u,v,w,p \n\ -model : defines viscosity and forcing function. Choose either: 0 (isoviscous), 1 (manufactured-broken), 2 (solcx-2d), 3 (sinker), 4 (subdomain jumps) \n"; /* Contributed by Dave May */ #include #include #define PROFILE_TIMING #define ASSEMBLE_LOWER_TRIANGULAR #define NSD 3 /* number of spatial dimensions */ #define NODES_PER_EL 8 /* nodes per element */ #define U_DOFS 3 /* degrees of freedom per velocity node */ #define P_DOFS 1 /* degrees of freedom per pressure node */ #define GAUSS_POINTS 8 /* Gauss point based evaluation */ typedef struct { PetscScalar gp_coords[NSD*GAUSS_POINTS]; PetscScalar eta[GAUSS_POINTS]; PetscScalar fx[GAUSS_POINTS]; PetscScalar fy[GAUSS_POINTS]; PetscScalar fz[GAUSS_POINTS]; PetscScalar hc[GAUSS_POINTS]; } GaussPointCoefficients; typedef struct { PetscScalar u_dof; PetscScalar v_dof; PetscScalar w_dof; PetscScalar p_dof; } StokesDOF; typedef struct _p_CellProperties *CellProperties; struct _p_CellProperties { PetscInt ncells; PetscInt mx,my,mz; PetscInt sex,sey,sez; GaussPointCoefficients *gpc; }; /* elements */ PetscErrorCode CellPropertiesCreate(DM da_stokes,CellProperties *C) { PetscErrorCode ierr; CellProperties cells; PetscInt mx,my,mz,sex,sey,sez; PetscFunctionBeginUser; ierr = PetscNew(&cells);CHKERRQ(ierr); ierr = DMDAGetElementsCorners(da_stokes,&sex,&sey,&sez);CHKERRQ(ierr); ierr = DMDAGetElementsSizes(da_stokes,&mx,&my,&mz);CHKERRQ(ierr); cells->mx = mx; cells->my = my; cells->mz = mz; cells->ncells = mx * my * mz; cells->sex = sex; cells->sey = sey; cells->sez = sez; ierr = PetscMalloc1(mx*my*mz,&cells->gpc);CHKERRQ(ierr); *C = cells; PetscFunctionReturn(0); } PetscErrorCode CellPropertiesDestroy(CellProperties *C) { PetscErrorCode ierr; CellProperties cells; PetscFunctionBeginUser; if (!C) PetscFunctionReturn(0); cells = *C; ierr = PetscFree(cells->gpc);CHKERRQ(ierr); ierr = PetscFree(cells);CHKERRQ(ierr); *C = NULL; PetscFunctionReturn(0); } PetscErrorCode CellPropertiesGetCell(CellProperties C,PetscInt II,PetscInt J,PetscInt K,GaussPointCoefficients **G) { PetscFunctionBeginUser; *G = &C->gpc[(II-C->sex) + (J-C->sey)*C->mx + (K-C->sez)*C->mx*C->my]; PetscFunctionReturn(0); } /* FEM routines */ /* Element: Local basis function ordering 1-----2 | | | | 0-----3 */ static void ShapeFunctionQ13D_Evaluate(PetscScalar _xi[],PetscScalar Ni[]) { PetscReal xi = PetscRealPart(_xi[0]); PetscReal eta = PetscRealPart(_xi[1]); PetscReal zeta = PetscRealPart(_xi[2]); Ni[0] = 0.125 * (1.0 - xi) * (1.0 - eta) * (1.0 - zeta); Ni[1] = 0.125 * (1.0 - xi) * (1.0 + eta) * (1.0 - zeta); Ni[2] = 0.125 * (1.0 + xi) * (1.0 + eta) * (1.0 - zeta); Ni[3] = 0.125 * (1.0 + xi) * (1.0 - eta) * (1.0 - zeta); Ni[4] = 0.125 * (1.0 - xi) * (1.0 - eta) * (1.0 + zeta); Ni[5] = 0.125 * (1.0 - xi) * (1.0 + eta) * (1.0 + zeta); Ni[6] = 0.125 * (1.0 + xi) * (1.0 + eta) * (1.0 + zeta); Ni[7] = 0.125 * (1.0 + xi) * (1.0 - eta) * (1.0 + zeta); } static void ShapeFunctionQ13D_Evaluate_dxi(PetscScalar _xi[],PetscScalar GNi[][NODES_PER_EL]) { PetscReal xi = PetscRealPart(_xi[0]); PetscReal eta = PetscRealPart(_xi[1]); PetscReal zeta = PetscRealPart(_xi[2]); /* xi deriv */ GNi[0][0] = -0.125 * (1.0 - eta) * (1.0 - zeta); GNi[0][1] = -0.125 * (1.0 + eta) * (1.0 - zeta); GNi[0][2] = 0.125 * (1.0 + eta) * (1.0 - zeta); GNi[0][3] = 0.125 * (1.0 - eta) * (1.0 - zeta); GNi[0][4] = -0.125 * (1.0 - eta) * (1.0 + zeta); GNi[0][5] = -0.125 * (1.0 + eta) * (1.0 + zeta); GNi[0][6] = 0.125 * (1.0 + eta) * (1.0 + zeta); GNi[0][7] = 0.125 * (1.0 - eta) * (1.0 + zeta); /* eta deriv */ GNi[1][0] = -0.125 * (1.0 - xi) * (1.0 - zeta); GNi[1][1] = 0.125 * (1.0 - xi) * (1.0 - zeta); GNi[1][2] = 0.125 * (1.0 + xi) * (1.0 - zeta); GNi[1][3] = -0.125 * (1.0 + xi) * (1.0 - zeta); GNi[1][4] = -0.125 * (1.0 - xi) * (1.0 + zeta); GNi[1][5] = 0.125 * (1.0 - xi) * (1.0 + zeta); GNi[1][6] = 0.125 * (1.0 + xi) * (1.0 + zeta); GNi[1][7] = -0.125 * (1.0 + xi) * (1.0 + zeta); /* zeta deriv */ GNi[2][0] = -0.125 * (1.0 - xi) * (1.0 - eta); GNi[2][1] = -0.125 * (1.0 - xi) * (1.0 + eta); GNi[2][2] = -0.125 * (1.0 + xi) * (1.0 + eta); GNi[2][3] = -0.125 * (1.0 + xi) * (1.0 - eta); GNi[2][4] = 0.125 * (1.0 - xi) * (1.0 - eta); GNi[2][5] = 0.125 * (1.0 - xi) * (1.0 + eta); GNi[2][6] = 0.125 * (1.0 + xi) * (1.0 + eta); GNi[2][7] = 0.125 * (1.0 + xi) * (1.0 - eta); } static void matrix_inverse_3x3(PetscScalar A[3][3],PetscScalar B[3][3]) { PetscScalar t4, t6, t8, t10, t12, t14, t17; t4 = A[2][0] * A[0][1]; t6 = A[2][0] * A[0][2]; t8 = A[1][0] * A[0][1]; t10 = A[1][0] * A[0][2]; t12 = A[0][0] * A[1][1]; t14 = A[0][0] * A[1][2]; t17 = 0.1e1 / (t4 * A[1][2] - t6 * A[1][1] - t8 * A[2][2] + t10 * A[2][1] + t12 * A[2][2] - t14 * A[2][1]); B[0][0] = (A[1][1] * A[2][2] - A[1][2] * A[2][1]) * t17; B[0][1] = -(A[0][1] * A[2][2] - A[0][2] * A[2][1]) * t17; B[0][2] = (A[0][1] * A[1][2] - A[0][2] * A[1][1]) * t17; B[1][0] = -(-A[2][0] * A[1][2] + A[1][0] * A[2][2]) * t17; B[1][1] = (-t6 + A[0][0] * A[2][2]) * t17; B[1][2] = -(-t10 + t14) * t17; B[2][0] = (-A[2][0] * A[1][1] + A[1][0] * A[2][1]) * t17; B[2][1] = -(-t4 + A[0][0] * A[2][1]) * t17; B[2][2] = (-t8 + t12) * t17; } static void ShapeFunctionQ13D_Evaluate_dx(PetscScalar GNi[][NODES_PER_EL],PetscScalar GNx[][NODES_PER_EL],PetscScalar coords[],PetscScalar *det_J) { PetscScalar J00,J01,J02,J10,J11,J12,J20,J21,J22; PetscInt n; PetscScalar iJ[3][3],JJ[3][3]; J00 = J01 = J02 = 0.0; J10 = J11 = J12 = 0.0; J20 = J21 = J22 = 0.0; for (n=0; neta; /* initialise element stiffness matrix */ ierr = PetscMemzero(Ae,sizeof(Ae));CHKERRQ(ierr); ierr = PetscMemzero(Ge,sizeof(Ge));CHKERRQ(ierr); ierr = PetscMemzero(De,sizeof(De));CHKERRQ(ierr); ierr = PetscMemzero(Ce,sizeof(Ce));CHKERRQ(ierr); /* form element stiffness matrix */ FormStressOperatorQ13D(Ae,el_coords,prop_eta); FormGradientOperatorQ13D(Ge,el_coords); /*#if defined(ASSEMBLE_LOWER_TRIANGULAR)*/ FormDivergenceOperatorQ13D(De,el_coords); /*#endif*/ FormStabilisationOperatorQ13D(Ce,el_coords,prop_eta); /* insert element matrix into global matrix */ ierr = DMDAGetElementEqnums3D_up(u_eqn,p_eqn,ei,ej,ek);CHKERRQ(ierr); for (n=0; neta; /* initialise element stiffness matrix */ ierr = PetscMemzero(Ae,sizeof(Ae));CHKERRQ(ierr); ierr = PetscMemzero(Ge,sizeof(Ge));CHKERRQ(ierr); ierr = PetscMemzero(De,sizeof(De));CHKERRQ(ierr); ierr = PetscMemzero(Ce,sizeof(Ce));CHKERRQ(ierr); /* form element stiffness matrix */ FormStressOperatorQ13D(Ae,el_coords,prop_eta); FormGradientOperatorQ13D(Ge,el_coords); /* FormDivergenceOperatorQ13D(De,el_coords); */ FormScaledMassMatrixOperatorQ13D(Ce,el_coords,prop_eta); /* insert element matrix into global matrix */ ierr = DMDAGetElementEqnums3D_up(u_eqn,p_eqn,ei,ej,ek);CHKERRQ(ierr); for (n=0; nfx; prop_fy = props->fy; prop_fz = props->fz; prop_hc = props->hc; /* initialise element stiffness matrix */ ierr = PetscMemzero(Fe,sizeof(Fe));CHKERRQ(ierr); ierr = PetscMemzero(He,sizeof(He));CHKERRQ(ierr); /* form element stiffness matrix */ FormMomentumRhsQ13D(Fe,el_coords,prop_fx,prop_fy,prop_fz); FormContinuityRhsQ13D(He,el_coords,prop_hc); /* insert element matrix into global matrix */ ierr = DMDAGetElementEqnums3D_up(u_eqn,p_eqn,ei,ej,ek);CHKERRQ(ierr); for (n=0; n\n"); /* coords */ ierr = DMDAGetGhostCorners(da,&si,&sj,&sk,&nx,&ny,&nz);CHKERRQ(ierr); N = nx * ny * nz; PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"\n",byte_order); PetscFPrintf(PETSC_COMM_SELF,vtk_fp," \n",si,si+nx-1,sj,sj+ny-1,sk,sk+nz-1); PetscFPrintf(PETSC_COMM_SELF,vtk_fp," \n",si,si+nx-1,sj,sj+ny-1,sk,sk+nz-1); memory_offset = 0; PetscFPrintf(PETSC_COMM_SELF,vtk_fp," \n"); PetscFPrintf(PETSC_COMM_SELF,vtk_fp," \n"); /* copy coordinates */ ierr = DMGetCoordinateDM(da,&cda);CHKERRQ(ierr); ierr = DMGetCoordinatesLocal(da,&coords);CHKERRQ(ierr); PetscFPrintf(PETSC_COMM_SELF,vtk_fp," \n",memory_offset); memory_offset = memory_offset + sizeof(PetscInt) + sizeof(PetscScalar)*N*3; PetscFPrintf(PETSC_COMM_SELF,vtk_fp," \n"); PetscFPrintf(PETSC_COMM_SELF,vtk_fp," \n"); for (f=0; f\n", field_name,memory_offset); memory_offset = memory_offset + sizeof(PetscInt) + sizeof(PetscScalar)*N; } PetscFPrintf(PETSC_COMM_SELF,vtk_fp," \n"); PetscFPrintf(PETSC_COMM_SELF,vtk_fp," \n"); PetscFPrintf(PETSC_COMM_SELF,vtk_fp," \n"); ierr = PetscMalloc1(N,&buffer);CHKERRQ(ierr); ierr = DMGetLocalVector(da,&l_FIELD);CHKERRQ(ierr); ierr = DMGlobalToLocalBegin(da, FIELD,INSERT_VALUES,l_FIELD);CHKERRQ(ierr); ierr = DMGlobalToLocalEnd(da,FIELD,INSERT_VALUES,l_FIELD);CHKERRQ(ierr); ierr = VecGetArray(l_FIELD,&_L_FIELD);CHKERRQ(ierr); PetscFPrintf(PETSC_COMM_SELF,vtk_fp," \n"); PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"_"); /* write coordinates */ { int length = sizeof(PetscScalar)*N*3; PetscScalar *allcoords; fwrite(&length,sizeof(int),1,vtk_fp); ierr = VecGetArray(coords,&allcoords);CHKERRQ(ierr); fwrite(allcoords,sizeof(PetscScalar),3*N,vtk_fp); ierr = VecRestoreArray(coords,&allcoords);CHKERRQ(ierr); } /* write fields */ for (f=0; f\n"); PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"\n"); ierr = PetscFree(buffer);CHKERRQ(ierr); ierr = VecRestoreArray(l_FIELD,&_L_FIELD);CHKERRQ(ierr); ierr = DMRestoreLocalVector(da,&l_FIELD);CHKERRQ(ierr); if (vtk_fp) { fclose(vtk_fp); vtk_fp = NULL; } PetscFunctionReturn(0); } PetscErrorCode DAViewVTK_write_PieceExtend(FILE *vtk_fp,PetscInt indent_level,DM da,const char local_file_prefix[]) { PetscMPIInt size,rank; MPI_Comm comm; const PetscInt *lx,*ly,*lz; PetscInt M,N,P,pM,pN,pP,sum,*olx,*oly,*olz; PetscInt *osx,*osy,*osz,*oex,*oey,*oez; PetscInt i,j,k,II,stencil; PetscErrorCode ierr; PetscFunctionBeginUser; /* create file name */ PetscObjectGetComm((PetscObject)da,&comm); MPI_Comm_size(comm,&size); MPI_Comm_rank(comm,&rank); ierr = DMDAGetInfo(da,0,&M,&N,&P,&pM,&pN,&pP,0,&stencil,0,0,0,0);CHKERRQ(ierr); ierr = DMDAGetOwnershipRanges(da,&lx,&ly,&lz);CHKERRQ(ierr); /* generate start,end list */ ierr = PetscMalloc1(pM+1,&olx);CHKERRQ(ierr); ierr = PetscMalloc1(pN+1,&oly);CHKERRQ(ierr); ierr = PetscMalloc1(pP+1,&olz);CHKERRQ(ierr); sum = 0; for (i=0; iM) oex[i]=M; } for (i=0; iM)oey[i]=N; } for (i=0; iP) oez[i]=P; } for (k=0; k\n", osx[i],oex[i]-1, osy[j],oey[j]-1, osz[k],oez[k]-1,name); } } } ierr = PetscFree(olx);CHKERRQ(ierr); ierr = PetscFree(oly);CHKERRQ(ierr); ierr = PetscFree(olz);CHKERRQ(ierr); ierr = PetscFree(osx);CHKERRQ(ierr); ierr = PetscFree(osy);CHKERRQ(ierr); ierr = PetscFree(osz);CHKERRQ(ierr); ierr = PetscFree(oex);CHKERRQ(ierr); ierr = PetscFree(oey);CHKERRQ(ierr); ierr = PetscFree(oez);CHKERRQ(ierr); PetscFunctionReturn(0); } PetscErrorCode DAView_3DVTK_PStructuredGrid(DM da,const char file_prefix[],const char local_file_prefix[]) { MPI_Comm comm; PetscMPIInt size,rank; char vtk_filename[PETSC_MAX_PATH_LEN]; FILE *vtk_fp = NULL; const char *byte_order = PetscBinaryBigEndian() ? "BigEndian" : "LittleEndian"; PetscInt M,N,P,si,sj,sk,nx,ny,nz; PetscInt i,dofs; PetscErrorCode ierr; PetscFunctionBeginUser; /* only master generates this file */ PetscObjectGetComm((PetscObject)da,&comm); MPI_Comm_size(comm,&size); MPI_Comm_rank(comm,&rank); if (rank != 0) PetscFunctionReturn(0); /* create file name */ ierr = PetscSNPrintf(vtk_filename,sizeof(vtk_filename),"%s.pvts",file_prefix);CHKERRQ(ierr); vtk_fp = fopen(vtk_filename,"w"); if (!vtk_fp) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SYS,"Cannot open file = %s \n",vtk_filename); /* (VTK) generate pvts header */ PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"\n"); PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"\n",byte_order); /* define size of the nodal mesh based on the cell DM */ ierr = DMDAGetInfo(da,0,&M,&N,&P,0,0,0,&dofs,0,0,0,0,0);CHKERRQ(ierr); ierr = DMDAGetGhostCorners(da,&si,&sj,&sk,&nx,&ny,&nz);CHKERRQ(ierr); PetscFPrintf(PETSC_COMM_SELF,vtk_fp," \n",0,M-1,0,N-1,0,P-1); /* note overlap = 1 for Q1 */ /* DUMP THE CELL REFERENCES */ PetscFPrintf(PETSC_COMM_SELF,vtk_fp," \n"); PetscFPrintf(PETSC_COMM_SELF,vtk_fp," \n"); PetscFPrintf(PETSC_COMM_SELF,vtk_fp," \n"); PetscFPrintf(PETSC_COMM_SELF,vtk_fp," \n",NSD); PetscFPrintf(PETSC_COMM_SELF,vtk_fp," \n"); PetscFPrintf(PETSC_COMM_SELF,vtk_fp," \n"); for (i=0; i\n",fieldname); } PetscFPrintf(PETSC_COMM_SELF,vtk_fp," \n"); /* write out the parallel information */ ierr = DAViewVTK_write_PieceExtend(vtk_fp,2,da,local_file_prefix);CHKERRQ(ierr); /* close the file */ PetscFPrintf(PETSC_COMM_SELF,vtk_fp," \n"); PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"\n"); if (vtk_fp) { fclose(vtk_fp); vtk_fp = NULL; } PetscFunctionReturn(0); } PetscErrorCode DAView3DPVTS(DM da, Vec x,const char NAME[]) { char vts_filename[PETSC_MAX_PATH_LEN]; char pvts_filename[PETSC_MAX_PATH_LEN]; PetscErrorCode ierr; PetscFunctionBeginUser; ierr = PetscSNPrintf(vts_filename,sizeof(vts_filename),"%s-mesh",NAME);CHKERRQ(ierr); ierr = DAView_3DVTK_StructuredGrid_appended(da,x,vts_filename);CHKERRQ(ierr); ierr = PetscSNPrintf(pvts_filename,sizeof(pvts_filename),"%s-mesh",NAME);CHKERRQ(ierr); ierr = DAView_3DVTK_PStructuredGrid(da,pvts_filename,vts_filename);CHKERRQ(ierr); PetscFunctionReturn(0); } PetscErrorCode KSPMonitorStokesBlocks(KSP ksp,PetscInt n,PetscReal rnorm,void *dummy) { PetscErrorCode ierr; PetscReal norms[4]; Vec Br,v,w; Mat A; PetscFunctionBeginUser; ierr = KSPGetOperators(ksp,&A,NULL);CHKERRQ(ierr); ierr = MatCreateVecs(A,&w,&v);CHKERRQ(ierr); ierr = KSPBuildResidual(ksp,v,w,&Br);CHKERRQ(ierr); ierr = VecStrideNorm(Br,0,NORM_2,&norms[0]);CHKERRQ(ierr); ierr = VecStrideNorm(Br,1,NORM_2,&norms[1]);CHKERRQ(ierr); ierr = VecStrideNorm(Br,2,NORM_2,&norms[2]);CHKERRQ(ierr); ierr = VecStrideNorm(Br,3,NORM_2,&norms[3]);CHKERRQ(ierr); ierr = VecDestroy(&v);CHKERRQ(ierr); ierr = VecDestroy(&w);CHKERRQ(ierr); ierr = PetscPrintf(PETSC_COMM_WORLD,"%3D KSP Component U,V,W,P residual norm [ %1.12e, %1.12e, %1.12e, %1.12e ]\n",n,norms[0],norms[1],norms[2],norms[3]);CHKERRQ(ierr); PetscFunctionReturn(0); } static PetscErrorCode PCMGSetupViaCoarsen(PC pc,DM da_fine) { PetscInt nlevels,k; PETSC_UNUSED PetscInt finest; DM *da_list,*daclist; Mat R; PetscErrorCode ierr; PetscFunctionBeginUser; nlevels = 1; PetscOptionsGetInt(NULL,NULL,"-levels",&nlevels,0); ierr = PetscMalloc1(nlevels,&da_list);CHKERRQ(ierr); for (k=0; kgp_coords[NSD*p] = gp_x; cell->gp_coords[NSD*p+1] = gp_y; cell->gp_coords[NSD*p+2] = gp_z; } } } } ierr = PetscOptionsGetInt(NULL,NULL,"-model",&model_definition,NULL);CHKERRQ(ierr); switch (model_definition) { case 0: /* isoviscous */ for (ek = sez; ek < sez+Kmz; ek++) { for (ej = sey; ej < sey+Jmy; ej++) { for (ei = sex; ei < sex+Imx; ei++) { GaussPointCoefficients *cell; ierr = CellPropertiesGetCell(cell_properties,ei,ej,ek,&cell);CHKERRQ(ierr); for (p = 0; p < GAUSS_POINTS; p++) { PetscReal coord_x = PetscRealPart(cell->gp_coords[NSD*p]); PetscReal coord_y = PetscRealPart(cell->gp_coords[NSD*p+1]); PetscReal coord_z = PetscRealPart(cell->gp_coords[NSD*p+2]); cell->eta[p] = 1.0; cell->fx[p] = 0.0*coord_x; cell->fy[p] = -PetscSinReal(2.2*PETSC_PI*coord_y)*PetscCosReal(1.0*PETSC_PI*coord_x); cell->fz[p] = 0.0*coord_z; cell->hc[p] = 0.0; } } } } break; case 1: /* manufactured */ for (ek = sez; ek < sez+Kmz; ek++) { for (ej = sey; ej < sey+Jmy; ej++) { for (ei = sex; ei < sex+Imx; ei++) { PetscReal eta,Fm[NSD],Fc,pos[NSD]; GaussPointCoefficients *cell; ierr = CellPropertiesGetCell(cell_properties,ei,ej,ek,&cell);CHKERRQ(ierr); for (p = 0; p < GAUSS_POINTS; p++) { PetscReal coord_x = PetscRealPart(cell->gp_coords[NSD*p]); PetscReal coord_y = PetscRealPart(cell->gp_coords[NSD*p+1]); PetscReal coord_z = PetscRealPart(cell->gp_coords[NSD*p+2]); pos[0] = coord_x; pos[1] = coord_y; pos[2] = coord_z; evaluate_MS_FrankKamentski(pos,NULL,NULL,&eta,Fm,&Fc); cell->eta[p] = eta; cell->fx[p] = Fm[0]; cell->fy[p] = Fm[1]; cell->fz[p] = Fm[2]; cell->hc[p] = Fc; } } } } break; case 2: /* solcx */ for (ek = sez; ek < sez+Kmz; ek++) { for (ej = sey; ej < sey+Jmy; ej++) { for (ei = sex; ei < sex+Imx; ei++) { GaussPointCoefficients *cell; ierr = CellPropertiesGetCell(cell_properties,ei,ej,ek,&cell);CHKERRQ(ierr); for (p = 0; p < GAUSS_POINTS; p++) { PetscReal coord_x = PetscRealPart(cell->gp_coords[NSD*p]); PetscReal coord_y = PetscRealPart(cell->gp_coords[NSD*p+1]); PetscReal coord_z = PetscRealPart(cell->gp_coords[NSD*p+2]); cell->eta[p] = 1.0; cell->fx[p] = 0.0; cell->fy[p] = -PetscSinReal(3.0*PETSC_PI*coord_y)*PetscCosReal(1.0*PETSC_PI*coord_x); cell->fz[p] = 0.0*coord_z; cell->hc[p] = 0.0; } } } } break; case 3: /* sinker */ for (ek = sez; ek < sez+Kmz; ek++) { for (ej = sey; ej < sey+Jmy; ej++) { for (ei = sex; ei < sex+Imx; ei++) { GaussPointCoefficients *cell; ierr = CellPropertiesGetCell(cell_properties,ei,ej,ek,&cell);CHKERRQ(ierr); for (p = 0; p < GAUSS_POINTS; p++) { PetscReal xp = PetscRealPart(cell->gp_coords[NSD*p]); PetscReal yp = PetscRealPart(cell->gp_coords[NSD*p+1]); PetscReal zp = PetscRealPart(cell->gp_coords[NSD*p+2]); cell->eta[p] = 1.0e-2; cell->fx[p] = 0.0; cell->fy[p] = 0.0; cell->fz[p] = 0.0; cell->hc[p] = 0.0; if ((PetscAbs(xp-0.5) < 0.2) && (PetscAbs(yp-0.5) < 0.2) && (PetscAbs(zp-0.5) < 0.2)) { cell->eta[p] = 1.0; cell->fz[p] = 1.0; } } } } } break; case 4: /* subdomain jumps */ for (ek = sez; ek < sez+Kmz; ek++) { for (ej = sey; ej < sey+Jmy; ej++) { for (ei = sex; ei < sex+Imx; ei++) { PetscReal opts_mag,opts_eta0; PetscInt px,py,pz; PetscBool jump; PetscMPIInt rr; GaussPointCoefficients *cell; opts_mag = 1.0; opts_eta0 = 1.e-2; ierr = PetscOptionsGetReal(NULL,NULL,"-jump_eta0",&opts_eta0,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetReal(NULL,NULL,"-jump_magnitude",&opts_mag,NULL);CHKERRQ(ierr); ierr = DMDAGetInfo(da_Stokes,NULL,NULL,NULL,NULL,&px,&py,&pz,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); rr = PetscGlobalRank%(px*py); if (px%2) jump = (PetscBool)(rr%2); else jump = (PetscBool)((rr/px)%2 ? rr%2 : !(rr%2)); rr = PetscGlobalRank/(px*py); if (rr%2) jump = (PetscBool)!jump; ierr = CellPropertiesGetCell(cell_properties,ei,ej,ek,&cell);CHKERRQ(ierr); for (p = 0; p < GAUSS_POINTS; p++) { PetscReal xp = PetscRealPart(cell->gp_coords[NSD*p]); PetscReal yp = PetscRealPart(cell->gp_coords[NSD*p+1]); cell->eta[p] = jump ? PetscPowReal(10.0,opts_mag) : opts_eta0; cell->fx[p] = 0.0; cell->fy[p] = -PetscSinReal(2.2*PETSC_PI*yp)*PetscCosReal(1.0*PETSC_PI*xp); cell->fz[p] = 0.0; cell->hc[p] = 0.0; } } } } break; default: SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"No default model is supported. Choose either -model {0,1,2,3}"); break; } ierr = DMDAVecRestoreArray(vel_cda,vel_coords,&_vel_coords);CHKERRQ(ierr); /* Generate a matrix with the correct non-zero pattern of type AIJ. This will work in parallel and serial */ ierr = DMCreateMatrix(da_Stokes,&A);CHKERRQ(ierr); ierr = DMCreateMatrix(da_Stokes,&B);CHKERRQ(ierr); ierr = DMCreateGlobalVector(da_Stokes,&X);CHKERRQ(ierr); ierr = DMCreateGlobalVector(da_Stokes,&f);CHKERRQ(ierr); /* assemble A11 */ ierr = MatZeroEntries(A);CHKERRQ(ierr); ierr = MatZeroEntries(B);CHKERRQ(ierr); ierr = VecZeroEntries(f);CHKERRQ(ierr); ierr = AssembleA_Stokes(A,da_Stokes,cell_properties);CHKERRQ(ierr); ierr = MatViewFromOptions(A,NULL,"-amat_view");CHKERRQ(ierr); ierr = AssembleA_PCStokes(B,da_Stokes,cell_properties);CHKERRQ(ierr); ierr = MatViewFromOptions(B,NULL,"-bmat_view");CHKERRQ(ierr); /* build force vector */ ierr = AssembleF_Stokes(f,da_Stokes,cell_properties);CHKERRQ(ierr); /* SOLVE */ ierr = KSPCreate(PETSC_COMM_WORLD,&ksp_S);CHKERRQ(ierr); ierr = KSPSetOptionsPrefix(ksp_S,"stokes_"); /* stokes */ CHKERRQ(ierr); ierr = KSPSetOperators(ksp_S,A,B);CHKERRQ(ierr); ierr = KSPSetFromOptions(ksp_S);CHKERRQ(ierr); { PC pc; const PetscInt ufields[] = {0,1,2},pfields[] = {3}; ierr = KSPGetPC(ksp_S,&pc);CHKERRQ(ierr); ierr = PCFieldSplitSetBlockSize(pc,4);CHKERRQ(ierr); ierr = PCFieldSplitSetFields(pc,"u",3,ufields,ufields);CHKERRQ(ierr); ierr = PCFieldSplitSetFields(pc,"p",1,pfields,pfields);CHKERRQ(ierr); } { PC pc; PetscBool same = PETSC_FALSE; ierr = KSPGetPC(ksp_S,&pc);CHKERRQ(ierr); ierr = PetscObjectTypeCompare((PetscObject)pc,PCMG,&same);CHKERRQ(ierr); if (same) { ierr = PCMGSetupViaCoarsen(pc,da_Stokes);CHKERRQ(ierr); } } { PC pc; PetscBool same = PETSC_FALSE; ierr = KSPGetPC(ksp_S,&pc);CHKERRQ(ierr); ierr = PetscObjectTypeCompare((PetscObject)pc,PCBDDC,&same);CHKERRQ(ierr); if (same) { ierr = KSPSetOperators(ksp_S,A,A);CHKERRQ(ierr); } } { PetscBool stokes_monitor = PETSC_FALSE; ierr = PetscOptionsGetBool(NULL,NULL,"-stokes_ksp_monitor_blocks",&stokes_monitor,0);CHKERRQ(ierr); if (stokes_monitor) { ierr = KSPMonitorSet(ksp_S,KSPMonitorStokesBlocks,NULL,NULL);CHKERRQ(ierr); } } if (resolve) { /* Test changing matrix data structure and resolve */ ierr = VecDuplicate(X,&X1);CHKERRQ(ierr); } ierr = KSPSolve(ksp_S,f,X);CHKERRQ(ierr); if (resolve) { Mat C; ierr = MatDuplicate(A,MAT_COPY_VALUES,&C);CHKERRQ(ierr); ierr = KSPSetOperators(ksp_S,C,C);CHKERRQ(ierr); ierr = KSPSolve(ksp_S,f,X1);CHKERRQ(ierr); ierr = MatDestroy(&C);CHKERRQ(ierr); ierr = VecDestroy(&X1);CHKERRQ(ierr); } ierr = PetscOptionsGetBool(NULL,NULL,"-write_pvts",&write_output,NULL);CHKERRQ(ierr); if (write_output) {ierr = DAView3DPVTS(da_Stokes,X,"up");CHKERRQ(ierr);} { PetscBool flg = PETSC_FALSE; char filename[PETSC_MAX_PATH_LEN]; ierr = PetscOptionsGetString(NULL,NULL,"-write_binary",filename,sizeof(filename),&flg);CHKERRQ(ierr); if (flg) { PetscViewer viewer; /* ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename[0]?filename:"ex42-binaryoutput",FILE_MODE_WRITE,&viewer);CHKERRQ(ierr); */ ierr = PetscViewerVTKOpen(PETSC_COMM_WORLD,"ex42.vts",FILE_MODE_WRITE,&viewer);CHKERRQ(ierr); ierr = VecView(X,viewer);CHKERRQ(ierr); ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); } } ierr = KSPGetIterationNumber(ksp_S,&its);CHKERRQ(ierr); /* verify */ if (model_definition == 1) { DM da_Stokes_analytic; Vec X_analytic; ierr = DMDACreateManufacturedSolution(mx,my,mz,&da_Stokes_analytic,&X_analytic);CHKERRQ(ierr); if (write_output) { ierr = DAView3DPVTS(da_Stokes_analytic,X_analytic,"ms");CHKERRQ(ierr); } ierr = DMDAIntegrateErrors3D(da_Stokes_analytic,X,X_analytic);CHKERRQ(ierr); if (write_output) { ierr = DAView3DPVTS(da_Stokes,X,"up2");CHKERRQ(ierr); } ierr = DMDestroy(&da_Stokes_analytic);CHKERRQ(ierr); ierr = VecDestroy(&X_analytic);CHKERRQ(ierr); } ierr = KSPDestroy(&ksp_S);CHKERRQ(ierr); ierr = VecDestroy(&X);CHKERRQ(ierr); ierr = VecDestroy(&f);CHKERRQ(ierr); ierr = MatDestroy(&A);CHKERRQ(ierr); ierr = MatDestroy(&B);CHKERRQ(ierr); ierr = CellPropertiesDestroy(&cell_properties);CHKERRQ(ierr); ierr = DMDestroy(&da_Stokes);CHKERRQ(ierr); PetscFunctionReturn(0); } int main(int argc,char **args) { PetscErrorCode ierr; PetscInt mx,my,mz; ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; mx = my = mz = 10; ierr = PetscOptionsGetInt(NULL,NULL,"-mx",&mx,NULL);CHKERRQ(ierr); my = mx; mz = mx; ierr = PetscOptionsGetInt(NULL,NULL,"-my",&my,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(NULL,NULL,"-mz",&mz,NULL);CHKERRQ(ierr); ierr = solve_stokes_3d_coupled(mx,my,mz);CHKERRQ(ierr); ierr = PetscFinalize(); return ierr; } /*TEST test: args: -stokes_ksp_monitor_short -stokes_ksp_converged_reason -stokes_pc_type lu test: suffix: 2 nsize: 3 args: -stokes_ksp_monitor_short -stokes_ksp_converged_reason -stokes_pc_type redundant -stokes_redundant_pc_type lu test: suffix: bddc_stokes nsize: 6 args: -mx 5 -my 4 -mz 3 -stokes_ksp_monitor_short -stokes_ksp_converged_reason -stokes_pc_type bddc -dm_mat_type is -stokes_pc_bddc_dirichlet_pc_type svd -stokes_pc_bddc_neumann_pc_type svd -stokes_pc_bddc_coarse_redundant_pc_type svd test: suffix: bddc_stokes_deluxe nsize: 6 args: -mx 5 -my 4 -mz 3 -stokes_ksp_monitor_short -stokes_ksp_converged_reason -stokes_pc_type bddc -dm_mat_type is -stokes_pc_bddc_dirichlet_pc_type svd -stokes_pc_bddc_neumann_pc_type svd -stokes_pc_bddc_coarse_redundant_pc_type svd -stokes_pc_bddc_use_deluxe_scaling -stokes_sub_schurs_posdef 0 -stokes_sub_schurs_symmetric -stokes_sub_schurs_mat_solver_type petsc test: requires: !single suffix: bddc_stokes_subdomainjump_deluxe nsize: 9 args: -model 4 -jump_magnitude 4 -mx 6 -my 6 -mz 2 -stokes_ksp_monitor_short -stokes_ksp_converged_reason -stokes_pc_type bddc -dm_mat_type is -stokes_pc_bddc_use_deluxe_scaling -stokes_sub_schurs_posdef 0 -stokes_sub_schurs_symmetric -stokes_sub_schurs_mat_solver_type petsc -stokes_pc_bddc_schur_layers 1 test: requires: !single suffix: 3 args: -stokes_ksp_converged_reason -stokes_pc_type fieldsplit -resolve test: suffix: tut_1 nsize: 4 requires: !single args: -stokes_ksp_monitor test: suffix: tut_2 nsize: 4 requires: !single args: -stokes_ksp_monitor -stokes_pc_type fieldsplit -stokes_pc_fieldsplit_type schur test: suffix: tut_3 nsize: 4 requires: !single args: -mx 20 -stokes_ksp_monitor -stokes_pc_type fieldsplit -stokes_pc_fieldsplit_type schur TEST*/