Skip to content

Commit 82e0c4d

Browse files
committed
Added random vertices distribution
1 parent d793726 commit 82e0c4d

File tree

6 files changed

+70
-14
lines changed

6 files changed

+70
-14
lines changed
Lines changed: 6 additions & 9 deletions
Original file line numberDiff line numberDiff line change
@@ -1,10 +1,7 @@
11
run,mesh_res,error_u,error_p
2-
0,2,23.34221,0.01392
3-
1,3,14.88356,0.00793
4-
2,4,11.04429,0.00485
5-
3,5,8.90114,0.00337
6-
4,6,6.49004,0.00260
7-
5,7,4.87041,0.00207
8-
6,8,3.99474,0.00165
9-
7,9,3.24332,0.00136
10-
8,10,2.65430,0.00114
2+
0,2,26.96934,0.01516
3+
1,3,18.26207,0.00982
4+
2,4,14.36768,0.00532
5+
3,5,11.16971,0.00438
6+
4,6,9.12436,0.00348
7+
5,7,7.85373,0.00257
-1.92 KB
Loading

examples/Hdiv-mixed/include/setup-dm.h

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -12,5 +12,6 @@
1212
// ---------------------------------------------------------------------------
1313
PetscErrorCode CreateDM(MPI_Comm comm, VecType vec_type, DM *dm);
1414
PetscErrorCode PerturbVerticesSmooth(DM dm);
15+
PetscErrorCode PerturbVerticesRandom(DM dm);
1516

1617
#endif // setupdm_h

examples/Hdiv-mixed/main.c

Lines changed: 7 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -30,6 +30,7 @@
3030
// (boundary is not working)
3131
// ./main -pc_type svd -problem darcy2d -dm_plex_dim 2 -dm_plex_box_faces 4,4 -bc_pressure 1
3232
// ./main -pc_type svd -problem darcy2d -dm_plex_dim 2 -dm_plex_box_faces 4,4 -bc_pressure 1,2,3,4
33+
// ./main -pc_type svd -problem darcy3d -dm_plex_dim 3 -dm_plex_box_faces 4,4,4 -view_solution -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 1,0.1,1
3334
const char help[] = "Solve H(div)-mixed problem using PETSc and libCEED\n";
3435

3536
#include "main.h"
@@ -83,8 +84,8 @@ int main(int argc, char **argv) {
8384
// ---------------------------------------------------------------------------
8485
// -- Initialize backend
8586
Ceed ceed;
86-
//CeedInit("/cpu/self/ref/serial", &ceed);
87-
CeedInit(app_ctx->ceed_resource, &ceed);
87+
CeedInit("/cpu/self/ref/serial", &ceed);
88+
//CeedInit(app_ctx->ceed_resource, &ceed);
8889
CeedMemType mem_type_backend;
8990
CeedGetPreferredMemType(ceed, &mem_type_backend);
9091

@@ -118,9 +119,11 @@ int main(int argc, char **argv) {
118119
PetscCall( CreateDM(comm, vec_type, &dm_H1) );
119120
// TODO: add mesh option
120121
// perturb to have smooth random mesh
121-
PetscCall( PerturbVerticesSmooth(dm) );
122-
PetscCall( PerturbVerticesSmooth(dm_H1) );
122+
//PetscCall( PerturbVerticesSmooth(dm) );
123+
//PetscCall( PerturbVerticesSmooth(dm_H1) );
123124

125+
PetscCall(PerturbVerticesRandom(dm) );
126+
PetscCall(PerturbVerticesRandom(dm_H1) );
124127
// ---------------------------------------------------------------------------
125128
// Process command line options
126129
// ---------------------------------------------------------------------------

examples/Hdiv-mixed/problems/darcy2d.c

Lines changed: 0 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -91,7 +91,6 @@ PetscErrorCode Hdiv_DARCY2D(Ceed ceed, ProblemData problem_data, DM dm,
9191
PetscCall( DMGetBoundingBox(dm, domain_min, domain_max) );
9292
for (PetscInt i=0; i<2; i++) domain_size[i] = domain_max[i] - domain_min[i];
9393

94-
printf("lx %f, ly %f \n",domain_size[0], domain_size[1]);
9594
darcy_ctx->kappa = kappa;
9695
darcy_ctx->rho_a0 = rho_a0;
9796
darcy_ctx->g = g;

examples/Hdiv-mixed/src/setup-dm.c

Lines changed: 56 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -68,3 +68,59 @@ PetscErrorCode PerturbVerticesSmooth(DM dm) {
6868
PetscCall( VecRestoreArray(coordinates,&coords) );
6969
PetscFunctionReturn(0);
7070
}
71+
72+
73+
PetscErrorCode PerturbVerticesRandom(DM dm) {
74+
75+
PetscFunctionBegin;
76+
Vec coordinates;
77+
PetscSection coordSection;
78+
PetscScalar *coords;
79+
PetscInt v,vStart,vEnd,offset, dim;
80+
PetscReal x, y, z;
81+
82+
PetscCall( DMGetDimension(dm,&dim) );
83+
PetscCall( DMGetCoordinateSection(dm, &coordSection) );
84+
PetscCall( DMGetCoordinatesLocal(dm, &coordinates) );
85+
PetscCall( DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd) );
86+
PetscCall( VecGetArray(coordinates,&coords) );
87+
PetscInt c_end, c_start, num_elem;
88+
PetscCall( DMPlexGetHeightStratum(dm, 0, &c_start, &c_end) );
89+
num_elem = c_end - c_start;
90+
91+
for(v=vStart; v<vEnd; v++) {
92+
PetscCall( PetscSectionGetOffset(coordSection,v,&offset) );
93+
if(dim==2) {
94+
PetscScalar nx = sqrt(num_elem);
95+
PetscReal domain_min[2], domain_max[2], domain_size[2];
96+
PetscCall( DMGetBoundingBox(dm, domain_min, domain_max) );
97+
for (PetscInt i=0; i<2; i++) domain_size[i] = domain_max[i] - domain_min[i];
98+
PetscReal hx = domain_size[0]/nx, hy = domain_size[1]/nx;
99+
x = coords[offset]; y = coords[offset+1];
100+
// perturb randomly O(h*sqrt(2)/3)
101+
PetscReal rx = ((PetscReal)rand())/((PetscReal)RAND_MAX)*(hx*0.471404);
102+
PetscReal ry = ((PetscReal)rand())/((PetscReal)RAND_MAX)*(hy*0.471404);
103+
PetscReal t = ((PetscReal)rand())/((PetscReal)RAND_MAX)*PETSC_PI;
104+
coords[offset ] = x + rx*PetscCosReal(t);
105+
coords[offset+1] = y + ry*PetscSinReal(t);
106+
} else {
107+
PetscScalar nx = cbrt(num_elem);
108+
PetscReal domain_min[3], domain_max[3], domain_size[3];
109+
PetscCall( DMGetBoundingBox(dm, domain_min, domain_max) );
110+
for (PetscInt i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i];
111+
PetscReal hx = domain_size[0]/nx, hy = domain_size[1]/nx,
112+
hz = domain_size[2]/nx;
113+
x = coords[offset]; y = coords[offset+1], z = coords[offset+2];
114+
// This is because 'boundary' is broken in 3D
115+
PetscReal rx = ((PetscReal)rand())/((PetscReal)RAND_MAX)*(hx*0.471404);
116+
PetscReal ry = ((PetscReal)rand())/((PetscReal)RAND_MAX)*(hy*0.471404);
117+
PetscReal rz = ((PetscReal)rand())/((PetscReal)RAND_MAX)*(hz*0.471404);
118+
PetscReal t = ((PetscReal)rand())/((PetscReal)RAND_MAX)*PETSC_PI;
119+
coords[offset ] = x + rx*PetscCosReal(t);
120+
coords[offset+1] = y + ry*PetscCosReal(t);
121+
coords[offset+2] = z + rz*PetscSinReal(t);
122+
}
123+
}
124+
PetscCall( VecRestoreArray(coordinates,&coords) );
125+
PetscFunctionReturn(0);
126+
}

0 commit comments

Comments
 (0)