1c4762a1bSJed Brown static char help[] = "3D, tri-quadratic hexahedra (Q1), displacement finite element formulation\n\ 2c4762a1bSJed Brown of linear elasticity. E=1.0, nu=1/3.\n\ 3c4762a1bSJed Brown Unit cube domain with Dirichlet boundary\n\n"; 4c4762a1bSJed Brown 5c4762a1bSJed Brown #include <petscdmplex.h> 6c4762a1bSJed Brown #include <petscsnes.h> 7c4762a1bSJed Brown #include <petscds.h> 8c4762a1bSJed Brown #include <petscdmforest.h> 9c4762a1bSJed Brown 10c4762a1bSJed Brown static PetscReal s_soft_alpha=1.e-3; 11c4762a1bSJed Brown static PetscReal s_mu=0.4; 12c4762a1bSJed Brown static PetscReal s_lambda=0.4; 13c4762a1bSJed Brown 14c4762a1bSJed Brown static void f0_bd_u_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 15c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 16c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 17c4762a1bSJed Brown PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 18c4762a1bSJed Brown { 19c4762a1bSJed Brown f0[0] = 1; /* x direction pull */ 20c4762a1bSJed Brown f0[1] = -x[2]; /* add a twist around x-axis */ 21c4762a1bSJed Brown f0[2] = x[1]; 22c4762a1bSJed Brown } 23c4762a1bSJed Brown 24c4762a1bSJed Brown static void f1_bd_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 25c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 26c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 27c4762a1bSJed Brown PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 28c4762a1bSJed Brown { 29c4762a1bSJed Brown const PetscInt Ncomp = dim; 30c4762a1bSJed Brown PetscInt comp, d; 31c4762a1bSJed Brown for (comp = 0; comp < Ncomp; ++comp) { 32c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 33c4762a1bSJed Brown f1[comp*dim+d] = 0.0; 34c4762a1bSJed Brown } 35c4762a1bSJed Brown } 36c4762a1bSJed Brown } 37c4762a1bSJed Brown 38c4762a1bSJed Brown /* gradU[comp*dim+d] = {u_x, u_y} or {u_x, u_y, u_z} */ 39c4762a1bSJed Brown static void f1_u_3d_alpha(PetscInt dim, PetscInt Nf, PetscInt NfAux, 40c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 41c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 42c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 43c4762a1bSJed Brown { 44c4762a1bSJed Brown PetscReal trace,mu=s_mu,lambda=s_lambda,rad; 45c4762a1bSJed Brown PetscInt i,j; 46c4762a1bSJed Brown for (i=0,rad=0.;i<dim;i++) { 47c4762a1bSJed Brown PetscReal t=x[i]; 48c4762a1bSJed Brown rad += t*t; 49c4762a1bSJed Brown } 50c4762a1bSJed Brown rad = PetscSqrtReal(rad); 51c4762a1bSJed Brown if (rad>0.25) { 52c4762a1bSJed Brown mu *= s_soft_alpha; 53c4762a1bSJed Brown lambda *= s_soft_alpha; /* we could keep the bulk the same like rubberish */ 54c4762a1bSJed Brown } 55c4762a1bSJed Brown for (i=0,trace=0; i < dim; ++i) { 56c4762a1bSJed Brown trace += PetscRealPart(u_x[i*dim+i]); 57c4762a1bSJed Brown } 58c4762a1bSJed Brown for (i=0; i < dim; ++i) { 59c4762a1bSJed Brown for (j=0; j < dim; ++j) { 60c4762a1bSJed Brown f1[i*dim+j] = mu*(u_x[i*dim+j]+u_x[j*dim+i]); 61c4762a1bSJed Brown } 62c4762a1bSJed Brown f1[i*dim+i] += lambda*trace; 63c4762a1bSJed Brown } 64c4762a1bSJed Brown } 65c4762a1bSJed Brown 66c4762a1bSJed Brown /* gradU[comp*dim+d] = {u_x, u_y} or {u_x, u_y, u_z} */ 67c4762a1bSJed Brown static void f1_u_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 68c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 69c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 70c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 71c4762a1bSJed Brown { 72c4762a1bSJed Brown PetscReal trace,mu=s_mu,lambda=s_lambda; 73c4762a1bSJed Brown PetscInt i,j; 74c4762a1bSJed Brown for (i=0,trace=0; i < dim; ++i) { 75c4762a1bSJed Brown trace += PetscRealPart(u_x[i*dim+i]); 76c4762a1bSJed Brown } 77c4762a1bSJed Brown for (i=0; i < dim; ++i) { 78c4762a1bSJed Brown for (j=0; j < dim; ++j) { 79c4762a1bSJed Brown f1[i*dim+j] = mu*(u_x[i*dim+j]+u_x[j*dim+i]); 80c4762a1bSJed Brown } 81c4762a1bSJed Brown f1[i*dim+i] += lambda*trace; 82c4762a1bSJed Brown } 83c4762a1bSJed Brown } 84c4762a1bSJed Brown 85c4762a1bSJed Brown /* 3D elasticity */ 86c4762a1bSJed Brown #define IDX(ii,jj,kk,ll) (27*ii+9*jj+3*kk+ll) 87c4762a1bSJed Brown 88c4762a1bSJed Brown void g3_uu_3d_private( PetscScalar g3[], const PetscReal mu, const PetscReal lambda) 89c4762a1bSJed Brown { 90c4762a1bSJed Brown if (1) { 91c4762a1bSJed Brown g3[0] += lambda; 92c4762a1bSJed Brown g3[0] += mu; 93c4762a1bSJed Brown g3[0] += mu; 94c4762a1bSJed Brown g3[4] += lambda; 95c4762a1bSJed Brown g3[8] += lambda; 96c4762a1bSJed Brown g3[10] += mu; 97c4762a1bSJed Brown g3[12] += mu; 98c4762a1bSJed Brown g3[20] += mu; 99c4762a1bSJed Brown g3[24] += mu; 100c4762a1bSJed Brown g3[28] += mu; 101c4762a1bSJed Brown g3[30] += mu; 102c4762a1bSJed Brown g3[36] += lambda; 103c4762a1bSJed Brown g3[40] += lambda; 104c4762a1bSJed Brown g3[40] += mu; 105c4762a1bSJed Brown g3[40] += mu; 106c4762a1bSJed Brown g3[44] += lambda; 107c4762a1bSJed Brown g3[50] += mu; 108c4762a1bSJed Brown g3[52] += mu; 109c4762a1bSJed Brown g3[56] += mu; 110c4762a1bSJed Brown g3[60] += mu; 111c4762a1bSJed Brown g3[68] += mu; 112c4762a1bSJed Brown g3[70] += mu; 113c4762a1bSJed Brown g3[72] += lambda; 114c4762a1bSJed Brown g3[76] += lambda; 115c4762a1bSJed Brown g3[80] += lambda; 116c4762a1bSJed Brown g3[80] += mu; 117c4762a1bSJed Brown g3[80] += mu; 118c4762a1bSJed Brown } else { 119c4762a1bSJed Brown int i,j,k,l; 120c4762a1bSJed Brown static int cc=-1; 121c4762a1bSJed Brown cc++; 122c4762a1bSJed Brown for (i = 0; i < 3; ++i) { 123c4762a1bSJed Brown for (j = 0; j < 3; ++j) { 124c4762a1bSJed Brown for (k = 0; k < 3; ++k) { 125c4762a1bSJed Brown for (l = 0; l < 3; ++l) { 126c4762a1bSJed Brown if (k==l && i==j) g3[IDX(i,j,k,l)] += lambda; 127c4762a1bSJed Brown if (i==k && j==l) g3[IDX(i,j,k,l)] += mu; 128c4762a1bSJed Brown if (i==l && j==k) g3[IDX(i,j,k,l)] += mu; 129c4762a1bSJed Brown if (k==l && i==j && !cc) (void) PetscPrintf(PETSC_COMM_WORLD,"g3[%d] += lambda;\n",IDX(i,j,k,l)); 130c4762a1bSJed Brown if (i==k && j==l && !cc) (void) PetscPrintf(PETSC_COMM_WORLD,"g3[%d] += mu;\n",IDX(i,j,k,l)); 131c4762a1bSJed Brown if (i==l && j==k && !cc) (void) PetscPrintf(PETSC_COMM_WORLD,"g3[%d] += mu;\n",IDX(i,j,k,l)); 132c4762a1bSJed Brown } 133c4762a1bSJed Brown } 134c4762a1bSJed Brown } 135c4762a1bSJed Brown } 136c4762a1bSJed Brown } 137c4762a1bSJed Brown } 138c4762a1bSJed Brown 139c4762a1bSJed Brown static void g3_uu_3d_alpha(PetscInt dim, PetscInt Nf, PetscInt NfAux, 140c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 141c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 142c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 143c4762a1bSJed Brown { 144c4762a1bSJed Brown PetscReal mu=s_mu, lambda=s_lambda,rad; 145c4762a1bSJed Brown PetscInt i; 146c4762a1bSJed Brown for (i=0,rad=0.;i<dim;i++) { 147c4762a1bSJed Brown PetscReal t=x[i]; 148c4762a1bSJed Brown rad += t*t; 149c4762a1bSJed Brown } 150c4762a1bSJed Brown rad = PetscSqrtReal(rad); 151c4762a1bSJed Brown if (rad>0.25) { 152c4762a1bSJed Brown mu *= s_soft_alpha; 153c4762a1bSJed Brown lambda *= s_soft_alpha; /* we could keep the bulk the same like rubberish */ 154c4762a1bSJed Brown } 155c4762a1bSJed Brown g3_uu_3d_private(g3,mu,lambda); 156c4762a1bSJed Brown } 157c4762a1bSJed Brown 158c4762a1bSJed Brown static void g3_uu_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 159c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 160c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 161c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 162c4762a1bSJed Brown { 163c4762a1bSJed Brown g3_uu_3d_private(g3,s_mu,s_lambda); 164c4762a1bSJed Brown } 165c4762a1bSJed Brown 166c4762a1bSJed Brown static void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 167c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 168c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 169c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 170c4762a1bSJed Brown { 171c4762a1bSJed Brown const PetscInt Ncomp = dim; 172c4762a1bSJed Brown PetscInt comp; 173c4762a1bSJed Brown 174c4762a1bSJed Brown for (comp = 0; comp < Ncomp; ++comp) f0[comp] = 0.0; 175c4762a1bSJed Brown } 176c4762a1bSJed Brown 177c4762a1bSJed Brown /* PI_i (x_i^4 - x_i^2) */ 178c4762a1bSJed Brown static void f0_u_x4(PetscInt dim, PetscInt Nf, PetscInt NfAux, 179c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 180c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 181c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 182c4762a1bSJed Brown { 183c4762a1bSJed Brown const PetscInt Ncomp = dim; 184c4762a1bSJed Brown PetscInt comp,i; 185c4762a1bSJed Brown 186c4762a1bSJed Brown for (comp = 0; comp < Ncomp; ++comp) { 187c4762a1bSJed Brown f0[comp] = 1e5; 188c4762a1bSJed Brown for (i = 0; i < Ncomp; ++i) { 189c4762a1bSJed Brown f0[comp] *= /* (comp+1)* */(x[i]*x[i]*x[i]*x[i] - x[i]*x[i]); /* assumes (0,1]^D domain */ 190c4762a1bSJed Brown } 191c4762a1bSJed Brown } 192c4762a1bSJed Brown } 193c4762a1bSJed Brown 194c4762a1bSJed Brown PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 195c4762a1bSJed Brown { 196c4762a1bSJed Brown const PetscInt Ncomp = dim; 197c4762a1bSJed Brown PetscInt comp; 198c4762a1bSJed Brown 199c4762a1bSJed Brown for (comp = 0; comp < Ncomp; ++comp) u[comp] = 0; 200c4762a1bSJed Brown return 0; 201c4762a1bSJed Brown } 202c4762a1bSJed Brown 203c4762a1bSJed Brown int main(int argc,char **args) 204c4762a1bSJed Brown { 205c4762a1bSJed Brown Mat Amat; 206c4762a1bSJed Brown PetscErrorCode ierr; 207c4762a1bSJed Brown SNES snes; 208c4762a1bSJed Brown KSP ksp; 209c4762a1bSJed Brown MPI_Comm comm; 210c4762a1bSJed Brown PetscMPIInt rank; 211956f8c0dSBarry Smith #if defined(PETSC_USE_LOG) 212c4762a1bSJed Brown PetscLogStage stage[17]; 213956f8c0dSBarry Smith #endif 214c4762a1bSJed Brown PetscBool test_nonzero_cols = PETSC_FALSE,use_nearnullspace = PETSC_TRUE,attach_nearnullspace = PETSC_FALSE; 215c4762a1bSJed Brown Vec xx,bb; 216c4762a1bSJed Brown PetscInt iter,i,N,dim = 3,cells[3] = {1,1,1},max_conv_its,local_sizes[7],run_type = 1; 217c4762a1bSJed Brown DM dm,distdm,basedm; 218c4762a1bSJed Brown PetscBool flg; 219c4762a1bSJed Brown char convType[256]; 220c4762a1bSJed Brown PetscReal Lx,mdisp[10],err[10]; 221c4762a1bSJed Brown const char * const options[10] = {"-ex56_dm_refine 0", 222c4762a1bSJed Brown "-ex56_dm_refine 1", 223c4762a1bSJed Brown "-ex56_dm_refine 2", 224c4762a1bSJed Brown "-ex56_dm_refine 3", 225c4762a1bSJed Brown "-ex56_dm_refine 4", 226c4762a1bSJed Brown "-ex56_dm_refine 5", 227c4762a1bSJed Brown "-ex56_dm_refine 6", 228c4762a1bSJed Brown "-ex56_dm_refine 7", 229c4762a1bSJed Brown "-ex56_dm_refine 8", 230c4762a1bSJed Brown "-ex56_dm_refine 9"}; 231c4762a1bSJed Brown PetscFunctionBeginUser; 232c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 233c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 2345f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm, &rank)); 235c4762a1bSJed Brown /* options */ 236c4762a1bSJed Brown ierr = PetscOptionsBegin(comm,NULL,"3D bilinear Q1 elasticity options","");CHKERRQ(ierr); 237c4762a1bSJed Brown { 238c4762a1bSJed Brown i = 3; 2395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsIntArray("-cells", "Number of (flux tube) processor in each dimension", "ex56.c", cells, &i, NULL)); 240c4762a1bSJed Brown 241c4762a1bSJed Brown Lx = 1.; /* or ne for rod */ 242c4762a1bSJed Brown max_conv_its = 3; 2435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-max_conv_its","Number of iterations in convergence study","",max_conv_its,&max_conv_its,NULL)); 2442c71b3e2SJacob Faibussowitsch PetscCheckFalse(max_conv_its<=0 || max_conv_its>7,PETSC_COMM_WORLD, PETSC_ERR_USER, "Bad number of iterations for convergence test (%D)",max_conv_its); 2455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-lx","Length of domain","",Lx,&Lx,NULL)); 2465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-alpha","material coefficient inside circle","",s_soft_alpha,&s_soft_alpha,NULL)); 2475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-test_nonzero_cols","nonzero test","",test_nonzero_cols,&test_nonzero_cols,NULL)); 2485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-use_mat_nearnullspace","MatNearNullSpace API test","",use_nearnullspace,&use_nearnullspace,NULL)); 2495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-attach_mat_nearnullspace","MatNearNullSpace API test (via MatSetNearNullSpace)","",attach_nearnullspace,&attach_nearnullspace,NULL)); 2505f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-run_type","0: twisting load on cantalever, 1: 3rd order accurate convergence test","",run_type,&run_type,NULL)); 251c4762a1bSJed Brown } 252c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 2535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStageRegister("Mesh Setup", &stage[16])); 254c4762a1bSJed Brown for (iter=0 ; iter<max_conv_its ; iter++) { 255c4762a1bSJed Brown char str[] = "Solve 0"; 256c4762a1bSJed Brown str[6] += iter; 2575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStageRegister(str, &stage[iter])); 258c4762a1bSJed Brown } 259c4762a1bSJed Brown /* create DM, Plex calls DMSetup */ 2605f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePush(stage[16])); 2615f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, cells, NULL, NULL, NULL, PETSC_TRUE, &dm)); 262c4762a1bSJed Brown { 263c4762a1bSJed Brown DMLabel label; 264c4762a1bSJed Brown IS is; 2655f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLabel(dm, "boundary")); 2665f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "boundary", &label)); 2675f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexMarkBoundaryFaces(dm, 1, label)); 268c4762a1bSJed Brown if (!run_type) { 2695f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetStratumIS(dm, "boundary", 1, &is)); 2705f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLabel(dm,"Faces")); 271c4762a1bSJed Brown if (is) { 272c4762a1bSJed Brown PetscInt d, f, Nf; 273c4762a1bSJed Brown const PetscInt *faces; 274c4762a1bSJed Brown PetscInt csize; 275c4762a1bSJed Brown PetscSection cs; 276c4762a1bSJed Brown Vec coordinates ; 277c4762a1bSJed Brown DM cdm; 2785f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(is, &Nf)); 2795f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(is, &faces)); 2805f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinatesLocal(dm, &coordinates)); 2815f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDM(dm, &cdm)); 2825f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalSection(cdm, &cs)); 283c4762a1bSJed Brown /* Check for each boundary face if any component of its centroid is either 0.0 or 1.0 */ 284c4762a1bSJed Brown for (f = 0; f < Nf; ++f) { 285c4762a1bSJed Brown PetscReal faceCoord; 286c4762a1bSJed Brown PetscInt b,v; 287c4762a1bSJed Brown PetscScalar *coords = NULL; 288c4762a1bSJed Brown PetscInt Nv; 2895f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexVecGetClosure(cdm, cs, coordinates, faces[f], &csize, &coords)); 290c4762a1bSJed Brown Nv = csize/dim; /* Calculate mean coordinate vector */ 291c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 292c4762a1bSJed Brown faceCoord = 0.0; 293c4762a1bSJed Brown for (v = 0; v < Nv; ++v) faceCoord += PetscRealPart(coords[v*dim+d]); 294c4762a1bSJed Brown faceCoord /= Nv; 295c4762a1bSJed Brown for (b = 0; b < 2; ++b) { 296c4762a1bSJed Brown if (PetscAbs(faceCoord - b) < PETSC_SMALL) { /* domain have not been set yet, still [0,1]^3 */ 2975f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetLabelValue(dm, "Faces", faces[f], d*2+b+1)); 298c4762a1bSJed Brown } 299c4762a1bSJed Brown } 300c4762a1bSJed Brown } 3015f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexVecRestoreClosure(cdm, cs, coordinates, faces[f], &csize, &coords)); 302c4762a1bSJed Brown } 3035f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(is, &faces)); 304c4762a1bSJed Brown } 3055f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&is)); 3065f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "Faces", &label)); 3075f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexLabelComplete(dm, label)); 308c4762a1bSJed Brown } 309c4762a1bSJed Brown } 310c4762a1bSJed Brown { 311c4762a1bSJed Brown PetscInt dimEmbed, i; 312c4762a1bSJed Brown PetscInt nCoords; 313c4762a1bSJed Brown PetscScalar *coords,bounds[] = {0,1,-.5,.5,-.5,.5,}; /* x_min,x_max,y_min,y_max */ 314c4762a1bSJed Brown Vec coordinates; 315c4762a1bSJed Brown bounds[1] = Lx; 316c4762a1bSJed Brown if (run_type==1) { 317c4762a1bSJed Brown for (i = 0; i < 2*dim; i++) bounds[i] = (i%2) ? 1 : 0; 318c4762a1bSJed Brown } 3195f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinatesLocal(dm,&coordinates)); 3205f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDim(dm,&dimEmbed)); 3212c71b3e2SJacob Faibussowitsch PetscCheckFalse(dimEmbed != dim,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"dimEmbed != dim %D",dimEmbed); 3225f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(coordinates,&nCoords)); 3232c71b3e2SJacob Faibussowitsch PetscCheckFalse(nCoords % dimEmbed,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Coordinate vector the wrong size"); 3245f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(coordinates,&coords)); 325c4762a1bSJed Brown for (i = 0; i < nCoords; i += dimEmbed) { 326c4762a1bSJed Brown PetscInt j; 327c4762a1bSJed Brown PetscScalar *coord = &coords[i]; 328c4762a1bSJed Brown for (j = 0; j < dimEmbed; j++) { 329c4762a1bSJed Brown coord[j] = bounds[2 * j] + coord[j] * (bounds[2 * j + 1] - bounds[2 * j]); 330c4762a1bSJed Brown } 331c4762a1bSJed Brown } 3325f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(coordinates,&coords)); 3335f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetCoordinatesLocal(dm,coordinates)); 334c4762a1bSJed Brown } 335c4762a1bSJed Brown 336c4762a1bSJed Brown /* convert to p4est, and distribute */ 337c4762a1bSJed Brown 338c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");CHKERRQ(ierr); 3395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsFList("-dm_type","Convert DMPlex to another format (should not be Plex!)","ex56.c",DMList,DMPLEX,convType,256,&flg)); 3401e1ea65dSPierre Jolivet ierr = PetscOptionsEnd();CHKERRQ(ierr); 341c4762a1bSJed Brown if (flg) { 342c4762a1bSJed Brown DM newdm; 3435f80ce2aSJacob Faibussowitsch CHKERRQ(DMConvert(dm,convType,&newdm)); 344c4762a1bSJed Brown if (newdm) { 345c4762a1bSJed Brown const char *prefix; 346c4762a1bSJed Brown PetscBool isForest; 3475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetOptionsPrefix((PetscObject)dm,&prefix)); 3485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject)newdm,prefix)); 3495f80ce2aSJacob Faibussowitsch CHKERRQ(DMIsForest(newdm,&isForest)); 350*28b400f6SJacob Faibussowitsch PetscCheck(isForest,PETSC_COMM_WORLD, PETSC_ERR_USER, "Converted to non Forest?"); 3515f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 352c4762a1bSJed Brown dm = newdm; 353c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_USER, "Convert failed?"); 354c4762a1bSJed Brown } else { 355c4762a1bSJed Brown PetscPartitioner part; 356c4762a1bSJed Brown /* Plex Distribute mesh over processes */ 3575f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetPartitioner(dm,&part)); 3585f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPartitionerSetFromOptions(part)); 3595f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexDistribute(dm, 0, NULL, &distdm)); 360c4762a1bSJed Brown if (distdm) { 361c4762a1bSJed Brown const char *prefix; 3625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetOptionsPrefix((PetscObject)dm,&prefix)); 3635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject)distdm,prefix)); 3645f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 365c4762a1bSJed Brown dm = distdm; 366c4762a1bSJed Brown } 367c4762a1bSJed Brown } 3685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePop()); 369c4762a1bSJed Brown basedm = dm; dm = NULL; 370c4762a1bSJed Brown 371c4762a1bSJed Brown for (iter=0 ; iter<max_conv_its ; iter++) { 3725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePush(stage[16])); 373c4762a1bSJed Brown /* make new DM */ 3745f80ce2aSJacob Faibussowitsch CHKERRQ(DMClone(basedm, &dm)); 3755f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) dm, "ex56_")); 3765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName( (PetscObject)dm,"Mesh")); 377c4762a1bSJed Brown if (max_conv_its > 1) { 3780e75e42fSMark /* If max_conv_its == 1, then we are not doing a convergence study. */ 3795f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInsertString(NULL,options[iter])); 380c4762a1bSJed Brown } 3815f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(dm)); /* refinement done here in Plex, p4est */ 382c4762a1bSJed Brown /* snes */ 3835f80ce2aSJacob Faibussowitsch CHKERRQ(SNESCreate(comm, &snes)); 3845f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetDM(snes, dm)); 385c4762a1bSJed Brown /* fem */ 386c4762a1bSJed Brown { 387c4762a1bSJed Brown const PetscInt Ncomp = dim; 388c4762a1bSJed Brown const PetscInt components[] = {0,1,2}; 389c4762a1bSJed Brown const PetscInt Nfid = 1, Npid = 1; 390c4762a1bSJed Brown const PetscInt fid[] = {1}; /* The fixed faces (x=0) */ 391c4762a1bSJed Brown const PetscInt pid[] = {2}; /* The faces with loading (x=L_x) */ 392c4762a1bSJed Brown PetscFE fe; 393c4762a1bSJed Brown PetscDS prob; 39445480ffeSMatthew G. Knepley DMLabel label; 395c4762a1bSJed Brown DM cdm = dm; 396c4762a1bSJed Brown 3975f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, dim, PETSC_FALSE, NULL, PETSC_DECIDE, &fe)); /* elasticity */ 3985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) fe, "deformation")); 399c4762a1bSJed Brown /* FEM prob */ 4005f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetField(dm, 0, NULL, (PetscObject) fe)); 4015f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateDS(dm)); 4025f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDS(dm, &prob)); 403c4762a1bSJed Brown /* setup problem */ 404c4762a1bSJed Brown if (run_type==1) { 4055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu_3d)); 4065f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(prob, 0, f0_u_x4, f1_u_3d)); 407c4762a1bSJed Brown } else { 40845480ffeSMatthew G. Knepley PetscWeakForm wf; 40945480ffeSMatthew G. Knepley PetscInt bd, i; 41045480ffeSMatthew G. Knepley 4115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu_3d_alpha)); 4125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(prob, 0, f0_u, f1_u_3d_alpha)); 41345480ffeSMatthew G. Knepley 4145f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "Faces", &label)); 4155f80ce2aSJacob Faibussowitsch CHKERRQ(DMAddBoundary(dm, DM_BC_NATURAL, "traction", label, Npid, pid, 0, Ncomp, components, NULL, NULL, NULL, &bd)); 4165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSGetBoundary(prob, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 4175f80ce2aSJacob Faibussowitsch for (i = 0; i < Npid; ++i) CHKERRQ(PetscWeakFormSetIndexBdResidual(wf, label, pid[i], 0, 0, 0, f0_bd_u_3d, 0, f1_bd_u)); 418c4762a1bSJed Brown } 419c4762a1bSJed Brown /* bcs */ 420c4762a1bSJed Brown if (run_type==1) { 421c4762a1bSJed Brown PetscInt id = 1; 4225f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "boundary", &label)); 4235f80ce2aSJacob Faibussowitsch CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) zero, NULL, NULL, NULL)); 424c4762a1bSJed Brown } else { 4255f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "Faces", &label)); 4265f80ce2aSJacob Faibussowitsch CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "fixed", label, Nfid, fid, 0, Ncomp, components, (void (*)(void)) zero, NULL, NULL, NULL)); 427c4762a1bSJed Brown } 428c4762a1bSJed Brown while (cdm) { 4295f80ce2aSJacob Faibussowitsch CHKERRQ(DMCopyDisc(dm, cdm)); 4305f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoarseDM(cdm, &cdm)); 431c4762a1bSJed Brown } 4325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEDestroy(&fe)); 433c4762a1bSJed Brown } 434c4762a1bSJed Brown /* vecs & mat */ 4355f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(dm,&xx)); 4365f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(xx, &bb)); 4375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) bb, "b")); 4385f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) xx, "u")); 4395f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(dm, &Amat)); 4405f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(Amat,MAT_SYMMETRIC,PETSC_TRUE)); /* Some matrix kernels can take advantage of symmetry if we set this. */ 4415f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(Amat,MAT_SYMMETRY_ETERNAL,PETSC_TRUE)); /* Inform PETSc that Amat is always symmetric, so info set above isn't lost. */ 4425f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetBlockSize(Amat,3)); 4435f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(Amat,MAT_SPD,PETSC_TRUE)); 4445f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(bb,&N)); 445c4762a1bSJed Brown local_sizes[iter] = N; 4465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(snes,"%D global equations, %D vertices\n",N,N/dim)); 447c4762a1bSJed Brown if ((use_nearnullspace || attach_nearnullspace) && N/dim > 1) { 448c4762a1bSJed Brown /* Set up the near null space (a.k.a. rigid body modes) that will be used by the multigrid preconditioner */ 449c4762a1bSJed Brown DM subdm; 450c4762a1bSJed Brown MatNullSpace nearNullSpace; 451c4762a1bSJed Brown PetscInt fields = 0; 452c4762a1bSJed Brown PetscObject deformation; 4535f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateSubDM(dm, 1, &fields, NULL, &subdm)); 4545f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateRigidBody(subdm, 0, &nearNullSpace)); 4555f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetField(dm, 0, NULL, &deformation)); 4565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectCompose(deformation, "nearnullspace", (PetscObject) nearNullSpace)); 4575f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&subdm)); 458c4762a1bSJed Brown if (attach_nearnullspace) { 4595f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetNearNullSpace(Amat,nearNullSpace)); 460c4762a1bSJed Brown } 4615f80ce2aSJacob Faibussowitsch CHKERRQ(MatNullSpaceDestroy(&nearNullSpace)); /* created by DM and destroyed by Mat */ 462c4762a1bSJed Brown } 4635f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSetSNESLocalFEM(dm,NULL,NULL,NULL)); 4645f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetJacobian(snes, Amat, Amat, NULL, NULL)); 4655f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFromOptions(snes)); 4665f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(dm)); 4675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePop()); 4685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePush(stage[16])); 469c4762a1bSJed Brown /* ksp */ 4705f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetKSP(snes, &ksp)); 4715f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetComputeSingularValues(ksp,PETSC_TRUE)); 472c4762a1bSJed Brown /* test BCs */ 4735f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(xx)); 474c4762a1bSJed Brown if (test_nonzero_cols) { 475dd400576SPatrick Sanan if (rank == 0) { 4765f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValue(xx,0,1.0,INSERT_VALUES)); 477c4762a1bSJed Brown } 4785f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(xx)); 4795f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(xx)); 480c4762a1bSJed Brown } 4815f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(bb)); 4825f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(bb,&i)); 483c4762a1bSJed Brown local_sizes[iter] = i; 4845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(snes,"%D equations in vector, %D vertices\n",i,i/dim)); 4855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePop()); 486c4762a1bSJed Brown /* solve */ 4875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePush(stage[iter])); 4885f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSolve(snes, bb, xx)); 4895f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePop()); 4905f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(xx,NORM_INFINITY,&mdisp[iter])); 4915f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(dm, NULL, "-dm_view")); 492c4762a1bSJed Brown { 493c4762a1bSJed Brown PetscViewer viewer = NULL; 494c4762a1bSJed Brown PetscViewerFormat fmt; 4955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetViewer(comm,NULL,"ex56_","-vec_view",&viewer,&fmt,&flg)); 496c4762a1bSJed Brown if (flg) { 4975f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPushFormat(viewer,fmt)); 4985f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(xx,viewer)); 4995f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(bb,viewer)); 5005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPopFormat(viewer)); 501c4762a1bSJed Brown } 5025f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&viewer)); 503c4762a1bSJed Brown } 504c4762a1bSJed Brown /* Free work space */ 5055f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 5065f80ce2aSJacob Faibussowitsch CHKERRQ(SNESDestroy(&snes)); 5075f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&xx)); 5085f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&bb)); 5095f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Amat)); 510c4762a1bSJed Brown } 5115f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&basedm)); 512c4762a1bSJed Brown if (run_type==1) err[0] = 59.975208 - mdisp[0]; /* error with what I think is the exact solution */ 513c4762a1bSJed Brown else err[0] = 171.038 - mdisp[0]; 514c4762a1bSJed Brown for (iter=1 ; iter<max_conv_its ; iter++) { 515c4762a1bSJed Brown if (run_type==1) err[iter] = 59.975208 - mdisp[iter]; 516c4762a1bSJed Brown else err[iter] = 171.038 - mdisp[iter]; 517c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"[%d] %D) N=%12D, max displ=%9.7e, disp diff=%9.2e, error=%4.3e, rate=%3.2g\n",rank,iter,local_sizes[iter],(double)mdisp[iter], 518c4762a1bSJed Brown (double)(mdisp[iter]-mdisp[iter-1]),(double)err[iter],(double)(PetscLogReal(err[iter-1]/err[iter])/PetscLogReal(2.)));CHKERRQ(ierr); 519c4762a1bSJed Brown } 520c4762a1bSJed Brown 521c4762a1bSJed Brown ierr = PetscFinalize(); 522c4762a1bSJed Brown return ierr; 523c4762a1bSJed Brown } 524c4762a1bSJed Brown 525c4762a1bSJed Brown /*TEST 526c4762a1bSJed Brown 527c4762a1bSJed Brown test: 528c4762a1bSJed Brown suffix: 0 529c4762a1bSJed Brown nsize: 4 530c4762a1bSJed Brown requires: !single 53173f7197eSJed Brown args: -cells 2,2,1 -max_conv_its 2 -petscspace_degree 2 -snes_max_it 2 -ksp_max_it 100 -ksp_type cg -ksp_rtol 1.e-10 -ksp_norm_type unpreconditioned -snes_rtol 1.e-10 -pc_type gamg -pc_gamg_esteig_ksp_max_it 10 -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 -pc_gamg_coarse_eq_limit 10 -pc_gamg_reuse_interpolation true -pc_gamg_square_graph 1 -pc_gamg_threshold 0.05 -pc_gamg_threshold_scale .0 -ksp_converged_reason -snes_monitor_short -ksp_monitor_short -snes_converged_reason -use_mat_nearnullspace true -mg_levels_ksp_max_it 2 -mg_levels_ksp_type chebyshev -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.1 -mg_levels_pc_type jacobi -petscpartitioner_type simple -matptap_via scalable -ex56_dm_view 532c4762a1bSJed Brown timeoutfactor: 2 533c4762a1bSJed Brown 534c4762a1bSJed Brown # HYPRE PtAP broken with complex numbers 535c4762a1bSJed Brown test: 536c4762a1bSJed Brown suffix: hypre 537263f2b91SStefano Zampini requires: hypre !single !complex !defined(PETSC_HAVE_HYPRE_DEVICE) 538c4762a1bSJed Brown nsize: 4 5390e75e42fSMark args: -cells 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -petscspace_degree 2 -ksp_type cg -ksp_monitor_short -ksp_rtol 1.e-8 -pc_type hypre -pc_hypre_type boomeramg -pc_hypre_boomeramg_no_CF true -pc_hypre_boomeramg_agg_nl 1 -pc_hypre_boomeramg_coarsen_type HMIS -pc_hypre_boomeramg_interp_type ext+i -ksp_converged_reason -use_mat_nearnullspace true -petscpartitioner_type simple 540c4762a1bSJed Brown 541c4762a1bSJed Brown test: 542c4762a1bSJed Brown suffix: ml 543c4762a1bSJed Brown requires: ml !single 544c4762a1bSJed Brown nsize: 4 5450e75e42fSMark args: -cells 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -petscspace_degree 2 -ksp_type cg -ksp_monitor_short -ksp_converged_reason -ksp_rtol 1.e-8 -pc_type ml -mg_levels_ksp_type chebyshev -mg_levels_ksp_max_it 3 -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 -mg_levels_pc_type sor -petscpartitioner_type simple -use_mat_nearnullspace 546c4762a1bSJed Brown 547c4762a1bSJed Brown test: 548c4762a1bSJed Brown suffix: hpddm 549dfd57a17SPierre Jolivet requires: hpddm slepc !single defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES) 550c4762a1bSJed Brown nsize: 4 5510e75e42fSMark args: -cells 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -petscspace_degree 2 -ksp_type fgmres -ksp_monitor_short -ksp_converged_reason -ksp_rtol 1.e-8 -pc_type hpddm -petscpartitioner_type simple -pc_hpddm_levels_1_sub_pc_type lu -pc_hpddm_levels_1_eps_nev 6 -pc_hpddm_coarse_p 1 -pc_hpddm_coarse_pc_type svd 552c4762a1bSJed Brown 553c4762a1bSJed Brown test: 55463b77682SMark Adams suffix: repart 555c4762a1bSJed Brown nsize: 4 556c4762a1bSJed Brown requires: parmetis !single 55773f7197eSJed Brown args: -cells 8,2,2 -max_conv_its 1 -petscspace_degree 2 -snes_max_it 4 -ksp_max_it 100 -ksp_type cg -ksp_rtol 1.e-2 -ksp_norm_type unpreconditioned -snes_rtol 1.e-3 -pc_type gamg -pc_gamg_esteig_ksp_max_it 10 -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 -pc_gamg_square_graph 1 -pc_gamg_threshold 0.05 -pc_gamg_threshold_scale .0 -use_mat_nearnullspace true -mg_levels_ksp_max_it 2 -mg_levels_ksp_type chebyshev -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 -mg_levels_pc_type jacobi -pc_gamg_mat_partitioning_type parmetis -pc_gamg_repartition true -snes_converged_reason -pc_gamg_process_eq_limit 20 -pc_gamg_coarse_eq_limit 10 -ksp_converged_reason -snes_converged_reason -pc_gamg_reuse_interpolation true 558c4762a1bSJed Brown 559c4762a1bSJed Brown test: 560c4762a1bSJed Brown suffix: bddc 561c4762a1bSJed Brown nsize: 4 562c4762a1bSJed Brown requires: !single 5630e75e42fSMark args: -cells 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -petscspace_degree 2 -ksp_type cg -ksp_monitor_short -ksp_rtol 1.e-8 -ksp_converged_reason -petscpartitioner_type simple -ex56_dm_mat_type is -matis_localmat_type {{sbaij baij aij}} -pc_type bddc 564c4762a1bSJed Brown 565c4762a1bSJed Brown testset: 566c4762a1bSJed Brown nsize: 4 567c4762a1bSJed Brown requires: !single 5680e75e42fSMark args: -cells 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -petscspace_degree 2 -ksp_type cg -ksp_monitor_short -ksp_rtol 1.e-10 -ksp_converged_reason -petscpartitioner_type simple -ex56_dm_mat_type is -matis_localmat_type aij -pc_type bddc -attach_mat_nearnullspace {{0 1}separate output} 569c4762a1bSJed Brown test: 570c4762a1bSJed Brown suffix: bddc_approx_gamg 57173f7197eSJed Brown args: -pc_bddc_switch_static -prefix_push pc_bddc_dirichlet_ -approximate -pc_type gamg -pc_gamg_esteig_ksp_max_it 10 -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 -pc_gamg_reuse_interpolation true -pc_gamg_square_graph 1 -pc_gamg_threshold 0.05 -pc_gamg_threshold_scale .0 -mg_levels_ksp_max_it 1 -mg_levels_ksp_type chebyshev -prefix_pop -prefix_push pc_bddc_neumann_ -approximate -pc_type gamg -pc_gamg_esteig_ksp_max_it 10 -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 -pc_gamg_coarse_eq_limit 10 -pc_gamg_reuse_interpolation true -pc_gamg_square_graph 1 -pc_gamg_threshold 0.05 -pc_gamg_threshold_scale .0 -mg_levels_ksp_max_it 1 -mg_levels_ksp_type chebyshev -prefix_pop 572c4762a1bSJed Brown # HYPRE PtAP broken with complex numbers 573c4762a1bSJed Brown test: 574263f2b91SStefano Zampini requires: hypre !complex !defined(PETSC_HAVE_HYPRE_DEVICE) 575c4762a1bSJed Brown suffix: bddc_approx_hypre 576c4762a1bSJed Brown args: -pc_bddc_switch_static -prefix_push pc_bddc_dirichlet_ -pc_type hypre -pc_hypre_boomeramg_no_CF true -pc_hypre_boomeramg_strong_threshold 0.75 -pc_hypre_boomeramg_agg_nl 1 -pc_hypre_boomeramg_coarsen_type HMIS -pc_hypre_boomeramg_interp_type ext+i -prefix_pop -prefix_push pc_bddc_neumann_ -pc_type hypre -pc_hypre_boomeramg_no_CF true -pc_hypre_boomeramg_strong_threshold 0.75 -pc_hypre_boomeramg_agg_nl 1 -pc_hypre_boomeramg_coarsen_type HMIS -pc_hypre_boomeramg_interp_type ext+i -prefix_pop 577c4762a1bSJed Brown test: 578c4762a1bSJed Brown requires: ml 579c4762a1bSJed Brown suffix: bddc_approx_ml 5800e75e42fSMark args: -pc_bddc_switch_static -prefix_push pc_bddc_dirichlet_ -approximate -pc_type ml -mg_levels_ksp_max_it 1 -mg_levels_ksp_type chebyshev -prefix_pop -prefix_push pc_bddc_neumann_ -approximate -pc_type ml -mg_levels_ksp_max_it 1 -mg_levels_ksp_type chebyshev -prefix_pop 581c4762a1bSJed Brown 582c4762a1bSJed Brown test: 583c4762a1bSJed Brown suffix: fetidp 584c4762a1bSJed Brown nsize: 4 585c4762a1bSJed Brown requires: !single 5860e75e42fSMark args: -cells 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -petscspace_degree 2 -ksp_type fetidp -fetidp_ksp_type cg -ksp_monitor_short -ksp_rtol 1.e-8 -ksp_converged_reason -petscpartitioner_type simple -ex56_dm_mat_type is -matis_localmat_type {{sbaij baij aij}} 587c4762a1bSJed Brown 588c4762a1bSJed Brown test: 589c4762a1bSJed Brown suffix: bddc_elast 590c4762a1bSJed Brown nsize: 4 591c4762a1bSJed Brown requires: !single 5920e75e42fSMark args: -cells 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -petscspace_degree 2 -ksp_type cg -ksp_monitor_short -ksp_rtol 1.e-8 -ksp_converged_reason -petscpartitioner_type simple -ex56_dm_mat_type is -matis_localmat_type sbaij -pc_type bddc -pc_bddc_monolithic -attach_mat_nearnullspace 593c4762a1bSJed Brown 594c4762a1bSJed Brown test: 595c4762a1bSJed Brown suffix: fetidp_elast 596c4762a1bSJed Brown nsize: 4 597c4762a1bSJed Brown requires: !single 5980e75e42fSMark args: -cells 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -petscspace_degree 2 -ksp_type fetidp -fetidp_ksp_type cg -ksp_monitor_short -ksp_rtol 1.e-8 -ksp_converged_reason -petscpartitioner_type simple -ex56_dm_mat_type is -matis_localmat_type sbaij -fetidp_bddc_pc_bddc_monolithic -attach_mat_nearnullspace 599c4762a1bSJed Brown 60035990778SJunchao Zhang testset: 60135990778SJunchao Zhang nsize: 4 60235990778SJunchao Zhang requires: !single 60373f7197eSJed Brown args: -cells 2,2,1 -max_conv_its 2 -petscspace_degree 2 -snes_max_it 2 -ksp_max_it 100 -ksp_type cg -ksp_rtol 1.e-10 -ksp_norm_type unpreconditioned -snes_rtol 1.e-10 -pc_type gamg -pc_gamg_esteig_ksp_max_it 10 -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 -pc_gamg_coarse_eq_limit 10 -pc_gamg_reuse_interpolation true -pc_gamg_square_graph 1 -pc_gamg_threshold 0.05 -pc_gamg_threshold_scale .0 -use_mat_nearnullspace true -mg_levels_ksp_max_it 2 -mg_levels_ksp_type chebyshev -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 -mg_levels_pc_type jacobi -ksp_monitor_short -ksp_converged_reason -snes_converged_reason -snes_monitor_short -ex56_dm_view -petscpartitioner_type simple -pc_gamg_process_eq_limit 20 60435990778SJunchao Zhang output_file: output/ex56_cuda.out 60535990778SJunchao Zhang 606c4762a1bSJed Brown test: 607c4762a1bSJed Brown suffix: cuda 60835990778SJunchao Zhang requires: cuda 60935990778SJunchao Zhang args: -ex56_dm_mat_type aijcusparse -ex56_dm_vec_type cuda 6106cb74228SMark Adams 6116cb74228SMark Adams test: 6126cb74228SMark Adams suffix: viennacl 61335990778SJunchao Zhang requires: viennacl 61435990778SJunchao Zhang args: -ex56_dm_mat_type aijviennacl -ex56_dm_vec_type viennacl 6156cb74228SMark Adams 61635990778SJunchao Zhang test: 61735990778SJunchao Zhang suffix: kokkos 6183078479eSJunchao Zhang requires: !sycl kokkos_kernels 61935990778SJunchao Zhang args: -ex56_dm_mat_type aijkokkos -ex56_dm_vec_type kokkos 620dea3b165SRichard Tran Mills # Don't run AIJMKL caes with complex scalars because of convergence issues. 621dea3b165SRichard Tran Mills # Note that we need to test both single and multiple MPI rank cases, because these use different sparse MKL routines to implement the PtAP operation. 622a8736bf8SRichard Tran Mills test: 623a8736bf8SRichard Tran Mills suffix: seqaijmkl 624a8736bf8SRichard Tran Mills nsize: 1 625dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MKL_SPARSE_OPTIMIZE) !single !complex 626a8736bf8SRichard Tran Mills args: -cells 2,2,1 -max_conv_its 2 -petscspace_degree 2 -snes_max_it 2 -ksp_max_it 100 -ksp_type cg -ksp_rtol 1.e-11 -ksp_norm_type unpreconditioned -snes_rtol 1.e-10 -pc_type gamg -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 -pc_gamg_coarse_eq_limit 1000 -pc_gamg_reuse_interpolation true -pc_gamg_square_graph 1 -pc_gamg_threshold 0.05 -pc_gamg_threshold_scale .0 -ksp_converged_reason -snes_monitor_short -ksp_monitor_short -snes_converged_reason -use_mat_nearnullspace true -mg_levels_ksp_max_it 1 -mg_levels_ksp_type chebyshev -pc_gamg_esteig_ksp_type cg -pc_gamg_esteig_ksp_max_it 10 -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.1 -mg_levels_pc_type jacobi -petscpartitioner_type simple -mat_block_size 3 -ex56_dm_view -run_type 1 -mat_seqaij_type seqaijmkl 627a8736bf8SRichard Tran Mills timeoutfactor: 2 628a8736bf8SRichard Tran Mills 629dea3b165SRichard Tran Mills test: 630dea3b165SRichard Tran Mills suffix: mpiaijmkl 631dea3b165SRichard Tran Mills nsize: 2 632dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MKL_SPARSE_OPTIMIZE) !single !complex 633dea3b165SRichard Tran Mills args: -cells 2,2,1 -max_conv_its 2 -petscspace_degree 2 -snes_max_it 2 -ksp_max_it 100 -ksp_type cg -ksp_rtol 1.e-11 -ksp_norm_type unpreconditioned -snes_rtol 1.e-10 -pc_type gamg -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 -pc_gamg_coarse_eq_limit 1000 -pc_gamg_reuse_interpolation true -pc_gamg_square_graph 1 -pc_gamg_threshold 0.05 -pc_gamg_threshold_scale .0 -ksp_converged_reason -snes_monitor_short -ksp_monitor_short -snes_converged_reason -use_mat_nearnullspace true -mg_levels_ksp_max_it 1 -mg_levels_ksp_type chebyshev -pc_gamg_esteig_ksp_type cg -pc_gamg_esteig_ksp_max_it 10 -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.1 -mg_levels_pc_type jacobi -petscpartitioner_type simple -mat_block_size 3 -ex56_dm_view -run_type 1 -mat_seqaij_type seqaijmkl 634dea3b165SRichard Tran Mills timeoutfactor: 2 635dea3b165SRichard Tran Mills 636c4762a1bSJed Brown TEST*/ 637