Skip to content
GitLab
Explore
Sign in
Register
Primary navigation
Search or go to…
Project
gmsh
Manage
Activity
Members
Labels
Plan
Issues
Issue boards
Milestones
Wiki
Code
Merge requests
Repository
Branches
Commits
Tags
Repository graph
Compare revisions
Snippets
Build
Pipelines
Jobs
Pipeline schedules
Artifacts
Deploy
Releases
Model registry
Operate
Environments
Monitor
Incidents
Analyze
Value stream analytics
Contributor analytics
CI/CD analytics
Repository analytics
Model experiments
Help
Help
Support
GitLab documentation
Compare GitLab plans
GitLab community forum
Contribute to GitLab
Provide feedback
Terms and privacy
Keyboard shortcuts
?
Snippets
Groups
Projects
Show more breadcrumbs
Larry Price
gmsh
Commits
75a67f77
Commit
75a67f77
authored
Mar 16, 2005
by
Jean-François Remacle
Browse files
Options
Downloads
Patches
Plain Diff
*** empty log message ***
parent
3813c5e1
No related branches found
No related tags found
No related merge requests found
Changes
4
Show whitespace changes
Inline
Side-by-side
Showing
4 changed files
Common/AdaptiveViews.cpp
+96
-598
96 additions, 598 deletions
Common/AdaptiveViews.cpp
Common/AdaptiveViews.h
+59
-22
59 additions, 22 deletions
Common/AdaptiveViews.h
Common/GmshMatrix.h
+1
-0
1 addition, 0 deletions
Common/GmshMatrix.h
Plugin/Levelset.cpp
+9
-9
9 additions, 9 deletions
Plugin/Levelset.cpp
with
165 additions
and
629 deletions
Common/AdaptiveViews.cpp
+
96
−
598
View file @
75a67f77
...
@@ -33,10 +33,14 @@
...
@@ -33,10 +33,14 @@
void
computeShapeFunctions
(
Double_Matrix
*
coeffs
,
Double_Matrix
*
eexps
,
double
u
,
double
v
,
double
w
,
double
*
sf
);
void
computeShapeFunctions
(
Double_Matrix
*
coeffs
,
Double_Matrix
*
eexps
,
double
u
,
double
v
,
double
w
,
double
*
sf
);
std
::
set
<
adapt_point
>
adapt_point
::
all_points
;
std
::
set
<
adapt_point
>
adapt_point
::
all_points
;
std
::
list
<
adapt_triangle
*>
adapt_triangle
::
all_triangles
;
std
::
list
<
adapt_triangle
*>
adapt_triangle
::
all_elems
;
std
::
list
<
adapt_tet
*>
adapt_tet
::
all_tets
;
std
::
list
<
adapt_tet
*>
adapt_tet
::
all_elems
;
std
::
list
<
adapt_quad
*>
adapt_quad
::
all_quads
;
std
::
list
<
adapt_quad
*>
adapt_quad
::
all_elems
;
std
::
list
<
adapt_hex
*>
adapt_hex
::
all_hexes
;
std
::
list
<
adapt_hex
*>
adapt_hex
::
all_elems
;
int
adapt_triangle
::
nbNod
=
3
;
int
adapt_tet
::
nbNod
=
4
;
int
adapt_quad
::
nbNod
=
4
;
int
adapt_hex
::
nbNod
=
8
;
adapt_point
*
adapt_point
::
New
(
double
x
,
double
y
,
double
z
,
Double_Matrix
*
coeffs
,
Double_Matrix
*
eexps
)
adapt_point
*
adapt_point
::
New
(
double
x
,
double
y
,
double
z
,
Double_Matrix
*
coeffs
,
Double_Matrix
*
eexps
)
{
{
...
@@ -54,60 +58,13 @@ adapt_point * adapt_point::New ( double x, double y, double z, Double_Matrix *co
...
@@ -54,60 +58,13 @@ adapt_point * adapt_point::New ( double x, double y, double z, Double_Matrix *co
else
else
return
(
adapt_point
*
)
&
(
*
it
);
return
(
adapt_point
*
)
&
(
*
it
);
}
}
void
adapt_triangle
::
clean
()
{
std
::
list
<
adapt_triangle
*>::
iterator
it
=
all_triangles
.
begin
();
std
::
list
<
adapt_triangle
*>::
iterator
ite
=
all_triangles
.
end
();
for
(;
it
!=
ite
;
++
it
)
{
delete
*
it
;
}
all_triangles
.
clear
();
adapt_point
::
all_points
.
clear
();
}
void
adapt_quad
::
clean
()
{
std
::
list
<
adapt_quad
*>::
iterator
it
=
all_quads
.
begin
();
std
::
list
<
adapt_quad
*>::
iterator
ite
=
all_quads
.
end
();
for
(;
it
!=
ite
;
++
it
)
{
delete
*
it
;
}
all_quads
.
clear
();
adapt_point
::
all_points
.
clear
();
}
void
adapt_tet
::
clean
()
{
std
::
list
<
adapt_tet
*>::
iterator
it
=
all_tets
.
begin
();
std
::
list
<
adapt_tet
*>::
iterator
ite
=
all_tets
.
end
();
for
(;
it
!=
ite
;
++
it
)
{
delete
*
it
;
}
all_tets
.
clear
();
adapt_point
::
all_points
.
clear
();
}
void
adapt_hex
::
clean
()
{
std
::
list
<
adapt_hex
*>::
iterator
it
=
all_hexes
.
begin
();
std
::
list
<
adapt_hex
*>::
iterator
ite
=
all_hexes
.
end
();
for
(;
it
!=
ite
;
++
it
)
{
delete
*
it
;
}
all_hexes
.
clear
();
adapt_point
::
all_points
.
clear
();
}
void
adapt_triangle
::
Create
(
int
maxlevel
,
Double_Matrix
*
coeffs
,
Double_Matrix
*
eexps
)
void
adapt_triangle
::
Create
(
int
maxlevel
,
Double_Matrix
*
coeffs
,
Double_Matrix
*
eexps
)
{
{
printf
(
"creating the sub-elements
\n
"
);
printf
(
"creating the sub-elements
\n
"
);
int
level
=
0
;
int
level
=
0
;
clean
();
clean
Element
<
adapt_triangle
>
();
adapt_point
*
p1
=
adapt_point
::
New
(
0
,
0
,
0
,
coeffs
,
eexps
);
adapt_point
*
p1
=
adapt_point
::
New
(
0
,
0
,
0
,
coeffs
,
eexps
);
adapt_point
*
p2
=
adapt_point
::
New
(
0
,
1
,
0
,
coeffs
,
eexps
);
adapt_point
*
p2
=
adapt_point
::
New
(
0
,
1
,
0
,
coeffs
,
eexps
);
adapt_point
*
p3
=
adapt_point
::
New
(
1
,
0
,
0
,
coeffs
,
eexps
);
adapt_point
*
p3
=
adapt_point
::
New
(
1
,
0
,
0
,
coeffs
,
eexps
);
...
@@ -119,7 +76,7 @@ void adapt_quad::Create (int maxlevel, Double_Matrix *coeffs, Double_Matrix *eex
...
@@ -119,7 +76,7 @@ void adapt_quad::Create (int maxlevel, Double_Matrix *coeffs, Double_Matrix *eex
{
{
printf
(
"creating the sub-elements
\n
"
);
printf
(
"creating the sub-elements
\n
"
);
int
level
=
0
;
int
level
=
0
;
clean
();
clean
Element
<
adapt_quad
>
();
adapt_point
*
p1
=
adapt_point
::
New
(
-
1
,
-
1
,
0
,
coeffs
,
eexps
);
adapt_point
*
p1
=
adapt_point
::
New
(
-
1
,
-
1
,
0
,
coeffs
,
eexps
);
adapt_point
*
p2
=
adapt_point
::
New
(
-
1
,
1
,
0
,
coeffs
,
eexps
);
adapt_point
*
p2
=
adapt_point
::
New
(
-
1
,
1
,
0
,
coeffs
,
eexps
);
adapt_point
*
p3
=
adapt_point
::
New
(
1
,
1
,
0
,
coeffs
,
eexps
);
adapt_point
*
p3
=
adapt_point
::
New
(
1
,
1
,
0
,
coeffs
,
eexps
);
...
@@ -132,7 +89,7 @@ void adapt_tet::Create (int maxlevel, Double_Matrix *coeffs, Double_Matrix *eexp
...
@@ -132,7 +89,7 @@ void adapt_tet::Create (int maxlevel, Double_Matrix *coeffs, Double_Matrix *eexp
{
{
Msg
(
INFO
,
"creating sub-elements"
);
Msg
(
INFO
,
"creating sub-elements"
);
int
level
=
0
;
int
level
=
0
;
clean
();
clean
Element
<
adapt_tet
>
();
adapt_point
*
p1
=
adapt_point
::
New
(
0
,
0
,
0
,
coeffs
,
eexps
);
adapt_point
*
p1
=
adapt_point
::
New
(
0
,
0
,
0
,
coeffs
,
eexps
);
adapt_point
*
p2
=
adapt_point
::
New
(
0
,
1
,
0
,
coeffs
,
eexps
);
adapt_point
*
p2
=
adapt_point
::
New
(
0
,
1
,
0
,
coeffs
,
eexps
);
adapt_point
*
p3
=
adapt_point
::
New
(
1
,
0
,
0
,
coeffs
,
eexps
);
adapt_point
*
p3
=
adapt_point
::
New
(
1
,
0
,
0
,
coeffs
,
eexps
);
...
@@ -145,7 +102,7 @@ void adapt_hex::Create (int maxlevel, Double_Matrix *coeffs, Double_Matrix *eexp
...
@@ -145,7 +102,7 @@ void adapt_hex::Create (int maxlevel, Double_Matrix *coeffs, Double_Matrix *eexp
{
{
Msg
(
INFO
,
"creating sub-elements"
);
Msg
(
INFO
,
"creating sub-elements"
);
int
level
=
0
;
int
level
=
0
;
clean
();
clean
Element
<
adapt_hex
>
();
adapt_point
*
p1
=
adapt_point
::
New
(
-
1
,
-
1
,
-
1
,
coeffs
,
eexps
);
adapt_point
*
p1
=
adapt_point
::
New
(
-
1
,
-
1
,
-
1
,
coeffs
,
eexps
);
adapt_point
*
p2
=
adapt_point
::
New
(
-
1
,
1
,
-
1
,
coeffs
,
eexps
);
adapt_point
*
p2
=
adapt_point
::
New
(
-
1
,
1
,
-
1
,
coeffs
,
eexps
);
adapt_point
*
p3
=
adapt_point
::
New
(
1
,
1
,
-
1
,
coeffs
,
eexps
);
adapt_point
*
p3
=
adapt_point
::
New
(
1
,
1
,
-
1
,
coeffs
,
eexps
);
...
@@ -160,7 +117,7 @@ void adapt_hex::Create (int maxlevel, Double_Matrix *coeffs, Double_Matrix *eexp
...
@@ -160,7 +117,7 @@ void adapt_hex::Create (int maxlevel, Double_Matrix *coeffs, Double_Matrix *eexp
void
adapt_triangle
::
Recur_Create
(
adapt_triangle
*
t
,
int
maxlevel
,
int
level
,
Double_Matrix
*
coeffs
,
Double_Matrix
*
eexps
)
void
adapt_triangle
::
Recur_Create
(
adapt_triangle
*
t
,
int
maxlevel
,
int
level
,
Double_Matrix
*
coeffs
,
Double_Matrix
*
eexps
)
{
{
all_
triang
les
.
push_back
(
t
);
all_
e
le
m
s
.
push_back
(
t
);
if
(
level
++
>=
maxlevel
)
if
(
level
++
>=
maxlevel
)
return
;
return
;
...
@@ -197,7 +154,7 @@ void adapt_triangle::Recur_Create (adapt_triangle *t, int maxlevel, int level ,
...
@@ -197,7 +154,7 @@ void adapt_triangle::Recur_Create (adapt_triangle *t, int maxlevel, int level ,
void
adapt_quad
::
Recur_Create
(
adapt_quad
*
q
,
int
maxlevel
,
int
level
,
Double_Matrix
*
coeffs
,
Double_Matrix
*
eexps
)
void
adapt_quad
::
Recur_Create
(
adapt_quad
*
q
,
int
maxlevel
,
int
level
,
Double_Matrix
*
coeffs
,
Double_Matrix
*
eexps
)
{
{
all_
quad
s
.
push_back
(
q
);
all_
elem
s
.
push_back
(
q
);
if
(
level
++
>=
maxlevel
)
if
(
level
++
>=
maxlevel
)
return
;
return
;
...
@@ -237,7 +194,7 @@ void adapt_quad::Recur_Create (adapt_quad *q, int maxlevel, int level , Double_M
...
@@ -237,7 +194,7 @@ void adapt_quad::Recur_Create (adapt_quad *q, int maxlevel, int level , Double_M
void
adapt_tet
::
Recur_Create
(
adapt_tet
*
t
,
int
maxlevel
,
int
level
,
Double_Matrix
*
coeffs
,
Double_Matrix
*
eexps
)
void
adapt_tet
::
Recur_Create
(
adapt_tet
*
t
,
int
maxlevel
,
int
level
,
Double_Matrix
*
coeffs
,
Double_Matrix
*
eexps
)
{
{
all_
tet
s
.
push_back
(
t
);
all_
elem
s
.
push_back
(
t
);
if
(
level
++
>=
maxlevel
)
if
(
level
++
>=
maxlevel
)
return
;
return
;
...
@@ -277,7 +234,7 @@ void adapt_tet::Recur_Create (adapt_tet *t, int maxlevel, int level , Double_Mat
...
@@ -277,7 +234,7 @@ void adapt_tet::Recur_Create (adapt_tet *t, int maxlevel, int level , Double_Mat
void
adapt_hex
::
Recur_Create
(
adapt_hex
*
h
,
int
maxlevel
,
int
level
,
Double_Matrix
*
coeffs
,
Double_Matrix
*
eexps
)
void
adapt_hex
::
Recur_Create
(
adapt_hex
*
h
,
int
maxlevel
,
int
level
,
Double_Matrix
*
coeffs
,
Double_Matrix
*
eexps
)
{
{
all_
hexe
s
.
push_back
(
h
);
all_
elem
s
.
push_back
(
h
);
if
(
level
++
>=
maxlevel
)
if
(
level
++
>=
maxlevel
)
return
;
return
;
...
@@ -338,24 +295,24 @@ void adapt_hex::Recur_Create (adapt_hex *h, int maxlevel, int level , Double_Mat
...
@@ -338,24 +295,24 @@ void adapt_hex::Recur_Create (adapt_hex *h, int maxlevel, int level , Double_Mat
void
adapt_triangle
::
Error
(
double
AVG
,
double
tol
)
void
adapt_triangle
::
Error
(
double
AVG
,
double
tol
)
{
{
adapt_triangle
*
t
=
*
all_
triang
les
.
begin
();
adapt_triangle
*
t
=
*
all_
e
le
m
s
.
begin
();
Recur_Error
(
t
,
AVG
,
tol
);
Recur_Error
(
t
,
AVG
,
tol
);
}
}
void
adapt_quad
::
Error
(
double
AVG
,
double
tol
)
void
adapt_quad
::
Error
(
double
AVG
,
double
tol
)
{
{
adapt_quad
*
q
=
*
all_
quad
s
.
begin
();
adapt_quad
*
q
=
*
all_
elem
s
.
begin
();
Recur_Error
(
q
,
AVG
,
tol
);
Recur_Error
(
q
,
AVG
,
tol
);
}
}
void
adapt_tet
::
Error
(
double
AVG
,
double
tol
)
void
adapt_tet
::
Error
(
double
AVG
,
double
tol
)
{
{
adapt_tet
*
t
=
*
all_
tet
s
.
begin
();
adapt_tet
*
t
=
*
all_
elem
s
.
begin
();
Recur_Error
(
t
,
AVG
,
tol
);
Recur_Error
(
t
,
AVG
,
tol
);
}
}
void
adapt_hex
::
Error
(
double
AVG
,
double
tol
)
void
adapt_hex
::
Error
(
double
AVG
,
double
tol
)
{
{
adapt_hex
*
h
=
*
all_
hexe
s
.
begin
();
adapt_hex
*
h
=
*
all_
elem
s
.
begin
();
Recur_Error
(
h
,
AVG
,
tol
);
Recur_Error
(
h
,
AVG
,
tol
);
}
}
...
@@ -572,22 +529,28 @@ void adapt_hex::Recur_Error ( adapt_hex *h, double AVG, double tol )
...
@@ -572,22 +529,28 @@ void adapt_hex::Recur_Error ( adapt_hex *h, double AVG, double tol )
static
double
t0
,
t1
,
t2
,
t3
;
static
double
t0
,
t1
,
t2
,
t3
;
void
Adaptive_Post_View
::
zoomTriangle
(
Post_View
*
view
,
template
<
class
ELEM
>
void
Adaptive_Post_View
::
zoomElement
(
Post_View
*
view
,
int
ielem
,
int
ielem
,
int
level
,
int
level
,
GMSH_Post_Plugin
*
plug
)
GMSH_Post_Plugin
*
plug
,
List_T
*
theList
,
int
*
counter
)
{
{
std
::
set
<
adapt_point
>::
iterator
it
=
adapt_point
::
all_points
.
begin
();
std
::
set
<
adapt_point
>::
iterator
ite
=
adapt_point
::
all_points
.
end
();
const
int
nbNod
=
ELEM
::
nbNod
;
typename
std
::
set
<
adapt_point
>::
iterator
it
=
adapt_point
::
all_points
.
begin
();
typename
std
::
set
<
adapt_point
>::
iterator
ite
=
adapt_point
::
all_points
.
end
();
double
c0
=
Cpu
();
double
c0
=
Cpu
();
const
int
N
=
_coefs
->
size1
();
const
int
N
=
_coefs
->
size1
();
Double_Vector
val
(
N
),
res
(
adapt_point
::
all_points
.
size
());
Double_Vector
val
(
N
),
res
(
adapt_point
::
all_points
.
size
());
Double_Matrix
xyz
(
3
,
3
);
Double_Matrix
xyz
(
nbNod
,
3
);
Double_Matrix
XYZ
(
adapt_point
::
all_points
.
size
(),
3
);
Double_Matrix
XYZ
(
adapt_point
::
all_points
.
size
(),
3
);
for
(
int
k
=
0
;
k
<
3
;
++
k
)
for
(
int
k
=
0
;
k
<
nbNod
;
++
k
)
{
{
xyz
(
k
,
0
)
=
(
*
_STposX
)
(
ielem
,
k
);
xyz
(
k
,
0
)
=
(
*
_STposX
)
(
ielem
,
k
);
xyz
(
k
,
1
)
=
(
*
_STposY
)
(
ielem
,
k
);
xyz
(
k
,
1
)
=
(
*
_STposY
)
(
ielem
,
k
);
...
@@ -618,8 +581,8 @@ void Adaptive_Post_View:: zoomTriangle (Post_View * view ,
...
@@ -618,8 +581,8 @@ void Adaptive_Post_View:: zoomTriangle (Post_View * view ,
double
c2
=
Cpu
();
double
c2
=
Cpu
();
std
::
list
<
adapt_triangle
*>::
iterator
itt
=
adapt_triangle
::
all_triang
les
.
begin
();
typename
std
::
list
<
ELEM
*>::
iterator
itt
=
ELEM
::
all_e
le
m
s
.
begin
();
std
::
list
<
adapt_triangle
*>::
iterator
itte
=
adapt_triangle
::
all_triang
les
.
end
();
typename
std
::
list
<
ELEM
*>::
iterator
itte
=
ELEM
::
all_e
le
m
s
.
end
();
for
(
;
itt
!=
itte
;
itt
++
)
for
(
;
itt
!=
itte
;
itt
++
)
{
{
...
@@ -629,133 +592,27 @@ void Adaptive_Post_View:: zoomTriangle (Post_View * view ,
...
@@ -629,133 +592,27 @@ void Adaptive_Post_View:: zoomTriangle (Post_View * view ,
if
(
!
plug
||
tol
!=
0.0
)
if
(
!
plug
||
tol
!=
0.0
)
{
{
adapt_triangle
::
Error
(
max
-
min
,
tol
);
ELEM
::
Error
(
max
-
min
,
tol
);
}
}
if
(
plug
)
if
(
plug
)
plug
->
assign_specific_visibility
();
plug
->
assign_specific_visibility
();
double
c3
=
Cpu
();
double
c3
=
Cpu
();
itt
=
adapt_triangle
::
all_triangles
.
begin
();
itt
=
ELEM
::
all_elems
.
begin
();
for
(
;
itt
!=
itte
;
itt
++
)
{
if
((
*
itt
)
->
visible
)
{
adapt_point
*
p1
=
(
*
itt
)
->
p
[
0
];
adapt_point
*
p2
=
(
*
itt
)
->
p
[
1
];
adapt_point
*
p3
=
(
*
itt
)
->
p
[
2
];
List_Add
(
view
->
ST
,
&
p1
->
X
);
List_Add
(
view
->
ST
,
&
p2
->
X
);
List_Add
(
view
->
ST
,
&
p3
->
X
);
List_Add
(
view
->
ST
,
&
p1
->
Y
);
List_Add
(
view
->
ST
,
&
p2
->
Y
);
List_Add
(
view
->
ST
,
&
p3
->
Y
);
List_Add
(
view
->
ST
,
&
p1
->
Z
);
List_Add
(
view
->
ST
,
&
p2
->
Z
);
List_Add
(
view
->
ST
,
&
p3
->
Z
);
List_Add
(
view
->
ST
,
&
p1
->
val
);
List_Add
(
view
->
ST
,
&
p2
->
val
);
List_Add
(
view
->
ST
,
&
p3
->
val
);
view
->
NbST
++
;
}
}
double
c4
=
Cpu
();
t0
+=
c1
-
c0
;
t1
+=
c2
-
c1
;
t2
+=
c3
-
c2
;
t3
+=
c4
-
c3
;
}
void
Adaptive_Post_View
::
zoomQuad
(
Post_View
*
view
,
int
ielem
,
int
level
,
GMSH_Post_Plugin
*
plug
)
{
std
::
set
<
adapt_point
>::
iterator
it
=
adapt_point
::
all_points
.
begin
();
std
::
set
<
adapt_point
>::
iterator
ite
=
adapt_point
::
all_points
.
end
();
double
c0
=
Cpu
();
const
int
N
=
_coefs
->
size1
();
Double_Vector
val
(
N
),
res
(
adapt_point
::
all_points
.
size
());
Double_Matrix
xyz
(
4
,
3
);
Double_Matrix
XYZ
(
adapt_point
::
all_points
.
size
(),
3
);
for
(
int
k
=
0
;
k
<
4
;
++
k
)
{
xyz
(
k
,
0
)
=
(
*
_STposX
)
(
ielem
,
k
);
xyz
(
k
,
1
)
=
(
*
_STposY
)
(
ielem
,
k
);
xyz
(
k
,
2
)
=
(
*
_STposZ
)
(
ielem
,
k
);
}
for
(
int
k
=
0
;
k
<
N
;
++
k
)
{
val
(
k
)
=
(
*
_STval
)(
ielem
,
k
);
}
_Interpolate
->
mult
(
val
,
res
);
_Geometry
->
mult
(
xyz
,
XYZ
);
double
c1
=
Cpu
();
int
kk
=
0
;
for
(
;
it
!=
ite
;
++
it
)
{
adapt_point
*
p
=
(
adapt_point
*
)
&
(
*
it
);
p
->
val
=
res
(
kk
);
p
->
X
=
XYZ
(
kk
,
0
);
p
->
Y
=
XYZ
(
kk
,
1
);
p
->
Z
=
XYZ
(
kk
,
2
);
if
(
min
>
p
->
val
)
min
=
p
->
val
;
if
(
max
<
p
->
val
)
max
=
p
->
val
;
kk
++
;
}
double
c2
=
Cpu
();
std
::
list
<
adapt_quad
*>::
iterator
itt
=
adapt_quad
::
all_quads
.
begin
();
std
::
list
<
adapt_quad
*>::
iterator
itte
=
adapt_quad
::
all_quads
.
end
();
for
(
;
itt
!=
itte
;
itt
++
)
{
(
*
itt
)
->
visible
=
false
;
}
if
(
!
plug
||
tol
!=
0.0
)
adapt_point
**
p
;
{
adapt_quad
::
Error
(
max
-
min
,
tol
);
}
if
(
plug
)
plug
->
assign_specific_visibility
();
double
c3
=
Cpu
();
itt
=
adapt_quad
::
all_quads
.
begin
();
for
(
;
itt
!=
itte
;
itt
++
)
for
(
;
itt
!=
itte
;
itt
++
)
{
{
if
((
*
itt
)
->
visible
)
if
((
*
itt
)
->
visible
)
{
{
adapt_point
*
p1
=
(
*
itt
)
->
p
[
0
];
p
=
(
*
itt
)
->
p
;
adapt_point
*
p2
=
(
*
itt
)
->
p
[
1
];
for
(
int
k
=
0
;
k
<
nbNod
;
++
k
)
List_Add
(
theList
,
&
p
[
k
]
->
X
);
adapt_point
*
p3
=
(
*
itt
)
->
p
[
2
];
for
(
int
k
=
0
;
k
<
nbNod
;
++
k
)
List_Add
(
theList
,
&
p
[
k
]
->
Y
);
adapt_point
*
p4
=
(
*
itt
)
->
p
[
3
];
for
(
int
k
=
0
;
k
<
nbNod
;
++
k
)
List_Add
(
theList
,
&
p
[
k
]
->
Z
);
List_Add
(
view
->
SQ
,
&
p1
->
X
);
for
(
int
k
=
0
;
k
<
nbNod
;
++
k
)
List_Add
(
theList
,
&
p
[
k
]
->
val
);
List_Add
(
view
->
SQ
,
&
p2
->
X
);
(
*
counter
)
++
;
List_Add
(
view
->
SQ
,
&
p3
->
X
);
List_Add
(
view
->
SQ
,
&
p4
->
X
);
List_Add
(
view
->
SQ
,
&
p1
->
Y
);
List_Add
(
view
->
SQ
,
&
p2
->
Y
);
List_Add
(
view
->
SQ
,
&
p3
->
Y
);
List_Add
(
view
->
SQ
,
&
p4
->
Y
);
List_Add
(
view
->
SQ
,
&
p1
->
Z
);
List_Add
(
view
->
SQ
,
&
p2
->
Z
);
List_Add
(
view
->
SQ
,
&
p3
->
Z
);
List_Add
(
view
->
SQ
,
&
p4
->
Z
);
List_Add
(
view
->
SQ
,
&
p1
->
val
);
List_Add
(
view
->
SQ
,
&
p2
->
val
);
List_Add
(
view
->
SQ
,
&
p3
->
val
);
List_Add
(
view
->
SQ
,
&
p4
->
val
);
view
->
NbSQ
++
;
}
}
}
}
double
c4
=
Cpu
();
double
c4
=
Cpu
();
...
@@ -764,253 +621,43 @@ void Adaptive_Post_View:: zoomQuad (Post_View * view ,
...
@@ -764,253 +621,43 @@ void Adaptive_Post_View:: zoomQuad (Post_View * view ,
t1
+=
c2
-
c1
;
t1
+=
c2
-
c1
;
t2
+=
c3
-
c2
;
t2
+=
c3
-
c2
;
t3
+=
c4
-
c3
;
t3
+=
c4
-
c3
;
}
void
Adaptive_Post_View
::
zoomTet
(
Post_View
*
view
,
int
ielem
,
int
level
,
GMSH_Post_Plugin
*
plug
,
Double_Vector
&
va2l
,
Double_Vector
&
re2s
,
Double_Matrix
&
XY2Z
)
{
std
::
set
<
adapt_point
>::
iterator
it
=
adapt_point
::
all_points
.
begin
();
std
::
set
<
adapt_point
>::
iterator
ite
=
adapt_point
::
all_points
.
end
();
double
c0
=
Cpu
();
Double_Matrix
xyz
(
4
,
3
);
for
(
int
k
=
0
;
k
<
4
;
++
k
)
{
xyz
(
k
,
0
)
=
(
*
_STposX
)
(
ielem
,
k
);
xyz
(
k
,
1
)
=
(
*
_STposY
)
(
ielem
,
k
);
xyz
(
k
,
2
)
=
(
*
_STposZ
)
(
ielem
,
k
);
}
const
int
N
=
_coefs
->
size1
();
Double_Vector
val
(
N
),
res
(
adapt_point
::
all_points
.
size
());
Double_Matrix
XYZ
(
adapt_point
::
all_points
.
size
(),
3
);
for
(
int
k
=
0
;
k
<
N
;
++
k
)
{
val
(
k
)
=
(
*
_STval
)(
ielem
,
k
);
}
}
_Interpolate
->
mult
(
val
,
res
);
_Geometry
->
mult
(
xyz
,
XYZ
);
double
c1
=
Cpu
();
int
kk
=
0
;
/*
for
(
;
it
!=
ite
;
++
it
)
We first do the adaptive stuff at level 2 and will only process
{
elements that have reached the maximal recursion level
adapt_point
*
p
=
(
adapt_point
*
)
&
(
*
it
);
p
->
val
=
res
(
kk
);
p
->
X
=
XYZ
(
kk
,
0
);
p
->
Y
=
XYZ
(
kk
,
1
);
p
->
Z
=
XYZ
(
kk
,
2
);
if
(
min
>
p
->
val
)
min
=
p
->
val
;
if
(
max
<
p
->
val
)
max
=
p
->
val
;
kk
++
;
}
double
c2
=
Cpu
();
We compute first the matrix at max recursion level (those should be stored
on disk once in the GMSHPLUGINSHOME directory, i'll do that at some point).
std
::
list
<
adapt_tet
*>::
iterator
itt
=
adapt_tet
::
all_tets
.
begin
();
std
::
list
<
adapt_tet
*>::
iterator
itte
=
adapt_tet
::
all_tets
.
end
();
for
(
;
itt
!=
itte
;
itt
++
)
{
(
*
itt
)
->
visible
=
false
;
}
*/
if
(
!
plug
||
tol
!=
0.0
)
{
adapt_tet
::
Error
(
max
-
min
,
tol
);
}
if
(
plug
)
plug
->
assign_specific_visibility
();
double
c3
=
Cpu
();
void
Adaptive_Post_View
::
setAdaptiveResolutionLevel
(
Post_View
*
view
,
int
level
,
GMSH_Post_Plugin
*
plug
)
itt
=
adapt_tet
::
all_tets
.
begin
();
for
(
;
itt
!=
itte
;
itt
++
)
{
if
((
*
itt
)
->
visible
)
{
{
adapt_point
*
p1
=
(
*
itt
)
->
p
[
0
];
setAdaptiveResolutionLevel_TEMPL
<
adapt_triangle
>
(
view
,
level
,
plug
,
&
(
view
->
ST
),
&
(
view
->
NbST
))
;
adapt_point
*
p2
=
(
*
itt
)
->
p
[
1
];
setAdaptiveResolutionLevel_TEMPL
<
adapt_quad
>
(
view
,
level
,
plug
,
&
(
view
->
SQ
),
&
(
view
->
NbSQ
))
;
adapt_point
*
p3
=
(
*
itt
)
->
p
[
2
];
setAdaptiveResolutionLevel_TEMPL
<
adapt_hex
>
(
view
,
level
,
plug
,
&
(
view
->
SH
),
&
(
view
->
NbSH
))
;
adapt_point
*
p4
=
(
*
itt
)
->
p
[
3
];
setAdaptiveResolutionLevel_TEMPL
<
adapt_tet
>
(
view
,
level
,
plug
,
&
(
view
->
SS
),
&
(
view
->
NbSS
))
;
List_Add
(
view
->
SS
,
&
p1
->
X
);
List_Add
(
view
->
SS
,
&
p2
->
X
);
List_Add
(
view
->
SS
,
&
p3
->
X
);
List_Add
(
view
->
SS
,
&
p4
->
X
);
List_Add
(
view
->
SS
,
&
p1
->
Y
);
List_Add
(
view
->
SS
,
&
p2
->
Y
);
List_Add
(
view
->
SS
,
&
p3
->
Y
);
List_Add
(
view
->
SS
,
&
p4
->
Y
);
List_Add
(
view
->
SS
,
&
p1
->
Z
);
List_Add
(
view
->
SS
,
&
p2
->
Z
);
List_Add
(
view
->
SS
,
&
p3
->
Z
);
List_Add
(
view
->
SS
,
&
p4
->
Z
);
List_Add
(
view
->
SS
,
&
p1
->
val
);
List_Add
(
view
->
SS
,
&
p2
->
val
);
List_Add
(
view
->
SS
,
&
p3
->
val
);
List_Add
(
view
->
SS
,
&
p4
->
val
);
view
->
NbSS
++
;
}
}
}
double
c4
=
Cpu
();
t0
+=
c1
-
c0
;
t1
+=
c2
-
c1
;
t2
+=
c3
-
c2
;
t3
+=
c4
-
c3
;
}
void
Adaptive_Post_View
::
zoomHex
(
Post_View
*
view
,
int
ielem
,
int
level
,
GMSH_Post_Plugin
*
plug
,
Double_Vector
&
va2l
,
Double_Vector
&
re2s
,
Double_Matrix
&
XY2Z
)
{
std
::
set
<
adapt_point
>::
iterator
it
=
adapt_point
::
all_points
.
begin
();
std
::
set
<
adapt_point
>::
iterator
ite
=
adapt_point
::
all_points
.
end
();
double
c0
=
Cpu
();
template
<
class
ELEM
>
void
Adaptive_Post_View
::
setAdaptiveResolutionLevel_TEMPL
(
Post_View
*
view
,
int
level
,
GMSH_Post_Plugin
*
plug
,
List_T
**
myList
,
int
*
counter
)
Double_Matrix
xyz
(
8
,
3
);
for
(
int
k
=
0
;
k
<
8
;
++
k
)
{
{
xyz
(
k
,
0
)
=
(
*
_STposX
)
(
ielem
,
k
);
xyz
(
k
,
1
)
=
(
*
_STposY
)
(
ielem
,
k
);
xyz
(
k
,
2
)
=
(
*
_STposZ
)
(
ielem
,
k
);
}
const
int
N
=
_coefs
->
size1
();
const
int
N
=
_coefs
->
size1
();
Double_Vector
val
(
N
),
res
(
adapt_point
::
all_points
.
size
());
Double_Matrix
XYZ
(
adapt_point
::
all_points
.
size
(),
3
);
for
(
int
k
=
0
;
k
<
N
;
++
k
)
{
val
(
k
)
=
(
*
_STval
)(
ielem
,
k
);
}
_Interpolate
->
mult
(
val
,
res
);
_Geometry
->
mult
(
xyz
,
XYZ
);
double
c1
=
Cpu
();
int
kk
=
0
;
for
(
;
it
!=
ite
;
++
it
)
{
adapt_point
*
p
=
(
adapt_point
*
)
&
(
*
it
);
p
->
val
=
res
(
kk
);
p
->
X
=
XYZ
(
kk
,
0
);
p
->
Y
=
XYZ
(
kk
,
1
);
p
->
Z
=
XYZ
(
kk
,
2
);
if
(
min
>
p
->
val
)
min
=
p
->
val
;
if
(
max
<
p
->
val
)
max
=
p
->
val
;
kk
++
;
}
double
c2
=
Cpu
();
std
::
list
<
adapt_hex
*>::
iterator
itt
=
adapt_hex
::
all_hexes
.
begin
();
std
::
list
<
adapt_hex
*>::
iterator
itte
=
adapt_hex
::
all_hexes
.
end
();
for
(
;
itt
!=
itte
;
itt
++
)
{
(
*
itt
)
->
visible
=
false
;
}
if
(
!
plug
||
tol
!=
0.0
)
{
adapt_hex
::
Error
(
max
-
min
,
tol
);
}
if
(
plug
)
plug
->
assign_specific_visibility
();
double
c3
=
Cpu
();
itt
=
adapt_hex
::
all_hexes
.
begin
();
for
(
;
itt
!=
itte
;
itt
++
)
{
if
((
*
itt
)
->
visible
)
{
adapt_point
*
p1
=
(
*
itt
)
->
p
[
0
];
adapt_point
*
p2
=
(
*
itt
)
->
p
[
1
];
adapt_point
*
p3
=
(
*
itt
)
->
p
[
2
];
adapt_point
*
p4
=
(
*
itt
)
->
p
[
3
];
adapt_point
*
p5
=
(
*
itt
)
->
p
[
4
];
adapt_point
*
p6
=
(
*
itt
)
->
p
[
5
];
adapt_point
*
p7
=
(
*
itt
)
->
p
[
6
];
adapt_point
*
p8
=
(
*
itt
)
->
p
[
7
];
List_Add
(
view
->
SH
,
&
p1
->
X
);
List_Add
(
view
->
SH
,
&
p2
->
X
);
List_Add
(
view
->
SH
,
&
p3
->
X
);
List_Add
(
view
->
SH
,
&
p4
->
X
);
List_Add
(
view
->
SH
,
&
p5
->
X
);
List_Add
(
view
->
SH
,
&
p6
->
X
);
List_Add
(
view
->
SH
,
&
p7
->
X
);
List_Add
(
view
->
SH
,
&
p8
->
X
);
List_Add
(
view
->
SH
,
&
p1
->
Y
);
List_Add
(
view
->
SH
,
&
p2
->
Y
);
List_Add
(
view
->
SH
,
&
p3
->
Y
);
List_Add
(
view
->
SH
,
&
p4
->
Y
);
List_Add
(
view
->
SH
,
&
p5
->
Y
);
List_Add
(
view
->
SH
,
&
p6
->
Y
);
List_Add
(
view
->
SH
,
&
p7
->
Y
);
List_Add
(
view
->
SH
,
&
p8
->
Y
);
List_Add
(
view
->
SH
,
&
p1
->
Z
);
List_Add
(
view
->
SH
,
&
p2
->
Z
);
List_Add
(
view
->
SH
,
&
p3
->
Z
);
List_Add
(
view
->
SH
,
&
p4
->
Z
);
List_Add
(
view
->
SH
,
&
p5
->
Z
);
List_Add
(
view
->
SH
,
&
p6
->
Z
);
List_Add
(
view
->
SH
,
&
p7
->
Z
);
List_Add
(
view
->
SH
,
&
p8
->
Z
);
List_Add
(
view
->
SH
,
&
p1
->
val
);
List_Add
(
view
->
SH
,
&
p2
->
val
);
List_Add
(
view
->
SH
,
&
p3
->
val
);
List_Add
(
view
->
SH
,
&
p4
->
val
);
List_Add
(
view
->
SH
,
&
p5
->
val
);
List_Add
(
view
->
SH
,
&
p6
->
val
);
List_Add
(
view
->
SH
,
&
p7
->
val
);
List_Add
(
view
->
SH
,
&
p8
->
val
);
view
->
NbSH
++
;
}
}
double
c4
=
Cpu
();
t0
+=
c1
-
c0
;
printf
(
"coucou %d %d
\n
"
,
*
counter
,
ELEM
::
nbNod
);
t1
+=
c2
-
c1
;
t2
+=
c3
-
c2
;
t3
+=
c4
-
c3
;
}
void
Adaptive_Post_View
::
setAdaptiveResolutionLevel
(
Post_View
*
view
,
int
level
,
GMSH_Post_Plugin
*
plug
)
{
const
int
N
=
_coefs
->
size1
();
if
(
presentTol
==
tol
&&
presentZoomLevel
==
level
&&
!
plug
)
return
;
if
(
presentTol
==
tol
&&
presentZoomLevel
==
level
&&
!
plug
)
return
;
if
(
view
->
NbST
)
if
(
*
counter
)
{
{
adapt_triangle
::
Create
(
level
,
_coefs
,
_eexps
);
double
sf
[
ELEM
::
nbNod
];
ELEM
::
Create
(
level
,
_coefs
,
_eexps
);
std
::
set
<
adapt_point
>::
iterator
it
=
adapt_point
::
all_points
.
begin
();
std
::
set
<
adapt_point
>::
iterator
it
=
adapt_point
::
all_points
.
begin
();
std
::
set
<
adapt_point
>::
iterator
ite
=
adapt_point
::
all_points
.
end
();
std
::
set
<
adapt_point
>::
iterator
ite
=
adapt_point
::
all_points
.
end
();
...
@@ -1019,7 +666,7 @@ void Adaptive_Post_View:: setAdaptiveResolutionLevel (Post_View * view , int lev
...
@@ -1019,7 +666,7 @@ void Adaptive_Post_View:: setAdaptiveResolutionLevel (Post_View * view , int lev
if
(
_Geometry
)
if
(
_Geometry
)
delete
_Geometry
;
delete
_Geometry
;
_Interpolate
=
new
Double_Matrix
(
adapt_point
::
all_points
.
size
(),
N
);
_Interpolate
=
new
Double_Matrix
(
adapt_point
::
all_points
.
size
(),
N
);
_Geometry
=
new
Double_Matrix
(
adapt_point
::
all_points
.
size
(),
3
);
_Geometry
=
new
Double_Matrix
(
adapt_point
::
all_points
.
size
(),
ELEM
::
nbNod
);
int
kk
=
0
;
int
kk
=
0
;
for
(
;
it
!=
ite
;
++
it
)
for
(
;
it
!=
ite
;
++
it
)
...
@@ -1029,182 +676,30 @@ void Adaptive_Post_View:: setAdaptiveResolutionLevel (Post_View * view , int lev
...
@@ -1029,182 +676,30 @@ void Adaptive_Post_View:: setAdaptiveResolutionLevel (Post_View * view , int lev
{
{
(
*
_Interpolate
)(
kk
,
k
)
=
p
->
shape_functions
[
k
];
(
*
_Interpolate
)(
kk
,
k
)
=
p
->
shape_functions
[
k
];
}
}
(
*
_Geometry
)(
kk
,
0
)
=
(
1.
-
p
->
x
-
p
->
y
);
ELEM
::
GSF
(
p
->
x
,
p
->
y
,
p
->
z
,
sf
);
(
*
_Geometry
)(
kk
,
1
)
=
p
->
x
;
for
(
int
k
=
0
;
k
<
ELEM
::
nbNod
;
k
++
)(
*
_Geometry
)(
kk
,
k
)
=
sf
[
k
];
(
*
_Geometry
)(
kk
,
2
)
=
p
->
y
;
kk
++
;
kk
++
;
}
}
List_Delete
(
view
->
ST
);
view
->
ST
=
0
;
List_Delete
(
*
myList
);
*
myList
=
0
;
view
->
NbST
=
0
;
*
counter
=
0
;
/// for now, that's all we do, 1 TS
/// for now, that's all we do, 1 TS
view
->
NbTimeStep
=
1
;
view
->
NbTimeStep
=
1
;
int
nbelm
=
_STposX
->
size1
();
int
nbelm
=
_STposX
->
size1
();
view
->
ST
=
List_Create
(
nbelm
*
4
,
nbelm
,
sizeof
(
double
));
*
myList
=
List_Create
(
nbelm
*
4
,
nbelm
,
sizeof
(
double
));
t0
=
t1
=
t2
=
t3
=
0
;
t0
=
t1
=
t2
=
t3
=
0
;
for
(
int
i
=
0
;
i
<
nbelm
;
++
i
)
for
(
int
i
=
0
;
i
<
nbelm
;
++
i
)
{
{
zoom
Triangle
(
view
,
i
,
level
,
plug
);
zoom
Element
<
ELEM
>
(
view
,
i
,
level
,
plug
,
*
myList
,
counter
);
}
}
printf
(
"finished %g %g %g %g
\n
"
,
t0
,
t1
,
t2
,
t3
);
printf
(
"finished %g %g %g %g
\n
"
,
t0
,
t1
,
t2
,
t3
);
}
else
if
(
view
->
NbSS
)
{
adapt_tet
::
Create
(
level
,
_coefs
,
_eexps
);
std
::
set
<
adapt_point
>::
iterator
it
=
adapt_point
::
all_points
.
begin
();
std
::
set
<
adapt_point
>::
iterator
ite
=
adapt_point
::
all_points
.
end
();
if
(
_Interpolate
)
delete
_Interpolate
;
if
(
_Geometry
)
delete
_Geometry
;
_Interpolate
=
new
Double_Matrix
(
adapt_point
::
all_points
.
size
(),
N
);
_Geometry
=
new
Double_Matrix
(
adapt_point
::
all_points
.
size
(),
4
);
int
kk
=
0
;
for
(
;
it
!=
ite
;
++
it
)
{
adapt_point
*
p
=
(
adapt_point
*
)
&
(
*
it
);
for
(
int
k
=
0
;
k
<
N
;
++
k
)
{
(
*
_Interpolate
)(
kk
,
k
)
=
p
->
shape_functions
[
k
];
}
(
*
_Geometry
)(
kk
,
0
)
=
(
1.
-
p
->
x
-
p
->
y
-
p
->
z
);
(
*
_Geometry
)(
kk
,
1
)
=
p
->
x
;
(
*
_Geometry
)(
kk
,
2
)
=
p
->
y
;
(
*
_Geometry
)(
kk
,
3
)
=
p
->
z
;
kk
++
;
}
List_Delete
(
view
->
SS
);
view
->
SS
=
0
;
view
->
NbSS
=
0
;
/// for now, that's all we do, 1 TS
view
->
NbTimeStep
=
1
;
int
nbelm
=
_STposX
->
size1
();
view
->
SS
=
List_Create
(
nbelm
*
4
,
nbelm
,
sizeof
(
double
));
t0
=
t1
=
t2
=
t3
=
0
;
Double_Vector
val
(
N
),
res
(
adapt_point
::
all_points
.
size
());
Double_Matrix
XYZ
(
adapt_point
::
all_points
.
size
(),
3
);
for
(
int
i
=
0
;
i
<
nbelm
;
++
i
)
{
zoomTet
(
view
,
i
,
level
,
plug
,
val
,
res
,
XYZ
);
}
printf
(
"finished B %g %g %g %g
\n
"
,
t0
,
t1
,
t2
,
t3
);
}
else
if
(
view
->
NbSH
)
{
adapt_hex
::
Create
(
level
,
_coefs
,
_eexps
);
std
::
set
<
adapt_point
>::
iterator
it
=
adapt_point
::
all_points
.
begin
();
std
::
set
<
adapt_point
>::
iterator
ite
=
adapt_point
::
all_points
.
end
();
if
(
_Interpolate
)
delete
_Interpolate
;
if
(
_Geometry
)
delete
_Geometry
;
_Interpolate
=
new
Double_Matrix
(
adapt_point
::
all_points
.
size
(),
N
);
_Geometry
=
new
Double_Matrix
(
adapt_point
::
all_points
.
size
(),
8
);
int
kk
=
0
;
for
(
;
it
!=
ite
;
++
it
)
{
adapt_point
*
p
=
(
adapt_point
*
)
&
(
*
it
);
for
(
int
k
=
0
;
k
<
N
;
++
k
)
{
(
*
_Interpolate
)(
kk
,
k
)
=
p
->
shape_functions
[
k
];
}
(
*
_Geometry
)(
kk
,
0
)
=
(
1.
-
p
->
x
)
*
(
1.
-
p
->
y
)
*
(
1.
-
p
->
z
)
*
.125
;
(
*
_Geometry
)(
kk
,
1
)
=
(
1.
+
p
->
x
)
*
(
1.
-
p
->
y
)
*
(
1.
-
p
->
z
)
*
.125
;
(
*
_Geometry
)(
kk
,
2
)
=
(
1.
+
p
->
x
)
*
(
1.
+
p
->
y
)
*
(
1.
-
p
->
z
)
*
.125
;
(
*
_Geometry
)(
kk
,
3
)
=
(
1.
-
p
->
x
)
*
(
1.
+
p
->
y
)
*
(
1.
-
p
->
z
)
*
.125
;
(
*
_Geometry
)(
kk
,
4
)
=
(
1.
-
p
->
x
)
*
(
1.
-
p
->
y
)
*
(
1.
+
p
->
z
)
*
.125
;
(
*
_Geometry
)(
kk
,
5
)
=
(
1.
+
p
->
x
)
*
(
1.
-
p
->
y
)
*
(
1.
+
p
->
z
)
*
.125
;
(
*
_Geometry
)(
kk
,
6
)
=
(
1.
+
p
->
x
)
*
(
1.
+
p
->
y
)
*
(
1.
+
p
->
z
)
*
.125
;
(
*
_Geometry
)(
kk
,
7
)
=
(
1.
-
p
->
x
)
*
(
1.
+
p
->
y
)
*
(
1.
+
p
->
z
)
*
.125
;
kk
++
;
}
List_Delete
(
view
->
SH
);
view
->
SH
=
0
;
view
->
NbSH
=
0
;
/// for now, that's all we do, 1 TS
view
->
NbTimeStep
=
1
;
int
nbelm
=
_STposX
->
size1
();
view
->
SH
=
List_Create
(
nbelm
*
4
,
nbelm
,
sizeof
(
double
));
t0
=
t1
=
t2
=
t3
=
0
;
Double_Vector
val
(
N
),
res
(
adapt_point
::
all_points
.
size
());
Double_Matrix
XYZ
(
adapt_point
::
all_points
.
size
(),
3
);
for
(
int
i
=
0
;
i
<
nbelm
;
++
i
)
{
zoomHex
(
view
,
i
,
level
,
plug
,
val
,
res
,
XYZ
);
}
printf
(
"finished B %g %g %g %g
\n
"
,
t0
,
t1
,
t2
,
t3
);
}
else
if
(
view
->
NbSQ
)
{
adapt_quad
::
Create
(
level
,
_coefs
,
_eexps
);
std
::
set
<
adapt_point
>::
iterator
it
=
adapt_point
::
all_points
.
begin
();
std
::
set
<
adapt_point
>::
iterator
ite
=
adapt_point
::
all_points
.
end
();
if
(
_Interpolate
)
delete
_Interpolate
;
if
(
_Geometry
)
delete
_Geometry
;
_Interpolate
=
new
Double_Matrix
(
adapt_point
::
all_points
.
size
(),
N
);
_Geometry
=
new
Double_Matrix
(
adapt_point
::
all_points
.
size
(),
4
);
int
kk
=
0
;
for
(
;
it
!=
ite
;
++
it
)
{
adapt_point
*
p
=
(
adapt_point
*
)
&
(
*
it
);
for
(
int
k
=
0
;
k
<
N
;
++
k
)
{
(
*
_Interpolate
)(
kk
,
k
)
=
p
->
shape_functions
[
k
];
}
(
*
_Geometry
)(
kk
,
0
)
=
(
1.
-
p
->
x
)
*
(
1.
-
p
->
y
)
*
.25
;
(
*
_Geometry
)(
kk
,
1
)
=
(
1.
+
p
->
x
)
*
(
1.
-
p
->
y
)
*
.25
;
(
*
_Geometry
)(
kk
,
2
)
=
(
1.
+
p
->
x
)
*
(
1.
+
p
->
y
)
*
.25
;
(
*
_Geometry
)(
kk
,
3
)
=
(
1.
-
p
->
x
)
*
(
1.
+
p
->
y
)
*
.25
;
kk
++
;
}
List_Delete
(
view
->
SQ
);
view
->
SQ
=
0
;
view
->
NbSQ
=
0
;
/// for now, that's all we do, 1 TS
view
->
NbTimeStep
=
1
;
int
nbelm
=
_STposX
->
size1
();
view
->
SQ
=
List_Create
(
nbelm
*
4
,
nbelm
,
sizeof
(
double
));
t0
=
t1
=
t2
=
t3
=
0
;
for
(
int
i
=
0
;
i
<
nbelm
;
++
i
)
{
zoomQuad
(
view
,
i
,
level
,
plug
);
}
printf
(
"finished Q %g %g %g %g
\n
"
,
t0
,
t1
,
t2
,
t3
);
}
else
return
;
view
->
Changed
=
1
;
view
->
Changed
=
1
;
presentZoomLevel
=
level
;
presentZoomLevel
=
level
;
presentTol
=
tol
;
presentTol
=
tol
;
}
}
}
void
computeShapeFunctions
(
Double_Matrix
*
coeffs
,
Double_Matrix
*
eexps
,
double
u
,
double
v
,
double
w
,
double
*
sf
)
void
computeShapeFunctions
(
Double_Matrix
*
coeffs
,
Double_Matrix
*
eexps
,
double
u
,
double
v
,
double
w
,
double
*
sf
)
{
{
...
@@ -1335,5 +830,8 @@ Adaptive_Post_View::~Adaptive_Post_View()
...
@@ -1335,5 +830,8 @@ Adaptive_Post_View::~Adaptive_Post_View()
delete
_STval
;
delete
_STval
;
if
(
_Interpolate
)
delete
_Interpolate
;
if
(
_Interpolate
)
delete
_Interpolate
;
if
(
_Geometry
)
delete
_Geometry
;
if
(
_Geometry
)
delete
_Geometry
;
adapt_triangle
::
clean
();
cleanElement
<
adapt_triangle
>
();
cleanElement
<
adapt_tet
>
();
cleanElement
<
adapt_hex
>
();
cleanElement
<
adapt_quad
>
();
}
}
This diff is collapsed.
Click to expand it.
Common/AdaptiveViews.h
+
59
−
22
View file @
75a67f77
...
@@ -73,11 +73,16 @@ public:
...
@@ -73,11 +73,16 @@ public:
{
{
return
(
p
[
0
]
->
val
+
p
[
1
]
->
val
+
p
[
2
]
->
val
)
/
3.
;
return
(
p
[
0
]
->
val
+
p
[
1
]
->
val
+
p
[
2
]
->
val
)
/
3.
;
}
}
inline
static
void
GSF
(
const
double
u
,
const
double
v
,
double
w
,
double
sf
[])
{
sf
[
0
]
=
1
-
u
-
v
;
sf
[
1
]
=
u
;
sf
[
2
]
=
v
;
}
void
print
()
void
print
()
{
{
printf
(
"p1 %g %g p2 %g %g p3 %g %g
\n
"
,
p
[
0
]
->
x
,
p
[
0
]
->
y
,
p
[
1
]
->
x
,
p
[
1
]
->
y
,
p
[
2
]
->
x
,
p
[
2
]
->
y
);
printf
(
"p1 %g %g p2 %g %g p3 %g %g
\n
"
,
p
[
0
]
->
x
,
p
[
0
]
->
y
,
p
[
1
]
->
x
,
p
[
1
]
->
y
,
p
[
2
]
->
x
,
p
[
2
]
->
y
);
}
}
static
void
clean
();
static
void
Create
(
int
maxlevel
,
Double_Matrix
*
coeffs
,
Double_Matrix
*
eexps
)
;
static
void
Create
(
int
maxlevel
,
Double_Matrix
*
coeffs
,
Double_Matrix
*
eexps
)
;
static
void
Recur_Create
(
adapt_triangle
*
t
,
int
maxlevel
,
int
level
,
Double_Matrix
*
coeffs
,
Double_Matrix
*
eexps
);
static
void
Recur_Create
(
adapt_triangle
*
t
,
int
maxlevel
,
int
level
,
Double_Matrix
*
coeffs
,
Double_Matrix
*
eexps
);
static
void
Error
(
double
AVG
,
double
tol
);
static
void
Error
(
double
AVG
,
double
tol
);
...
@@ -85,7 +90,8 @@ public:
...
@@ -85,7 +90,8 @@ public:
bool
visible
;
bool
visible
;
adapt_point
*
p
[
3
];
adapt_point
*
p
[
3
];
adapt_triangle
*
t
[
4
];
adapt_triangle
*
t
[
4
];
static
std
::
list
<
adapt_triangle
*>
all_triangles
;
static
std
::
list
<
adapt_triangle
*>
all_elems
;
static
int
nbNod
;
};
};
class
adapt_quad
class
adapt_quad
...
@@ -105,11 +111,17 @@ public:
...
@@ -105,11 +111,17 @@ public:
{
{
return
(
p
[
0
]
->
val
+
p
[
1
]
->
val
+
p
[
2
]
->
val
+
p
[
3
]
->
val
)
/
4.
;
return
(
p
[
0
]
->
val
+
p
[
1
]
->
val
+
p
[
2
]
->
val
+
p
[
3
]
->
val
)
/
4.
;
}
}
inline
static
void
GSF
(
const
double
u
,
const
double
v
,
double
w
,
double
sf
[])
{
sf
[
0
]
=
0.25
*
(
1
-
u
)
*
(
1
-
v
);
sf
[
1
]
=
0.25
*
(
1
+
u
)
*
(
1
-
v
);
sf
[
2
]
=
0.25
*
(
1
+
u
)
*
(
1
+
v
);
sf
[
3
]
=
0.25
*
(
1
-
u
)
*
(
1
+
v
);
}
void
print
()
void
print
()
{
{
printf
(
"p1 %g %g p2 %g %g p3 %g %g
\n
"
,
p
[
0
]
->
x
,
p
[
0
]
->
y
,
p
[
1
]
->
x
,
p
[
1
]
->
y
,
p
[
2
]
->
x
,
p
[
2
]
->
y
);
printf
(
"p1 %g %g p2 %g %g p3 %g %g
\n
"
,
p
[
0
]
->
x
,
p
[
0
]
->
y
,
p
[
1
]
->
x
,
p
[
1
]
->
y
,
p
[
2
]
->
x
,
p
[
2
]
->
y
);
}
}
static
void
clean
();
static
void
Create
(
int
maxlevel
,
Double_Matrix
*
coeffs
,
Double_Matrix
*
eexps
)
;
static
void
Create
(
int
maxlevel
,
Double_Matrix
*
coeffs
,
Double_Matrix
*
eexps
)
;
static
void
Recur_Create
(
adapt_quad
*
q
,
int
maxlevel
,
int
level
,
Double_Matrix
*
coeffs
,
Double_Matrix
*
eexps
);
static
void
Recur_Create
(
adapt_quad
*
q
,
int
maxlevel
,
int
level
,
Double_Matrix
*
coeffs
,
Double_Matrix
*
eexps
);
static
void
Error
(
double
AVG
,
double
tol
);
static
void
Error
(
double
AVG
,
double
tol
);
...
@@ -117,7 +129,8 @@ public:
...
@@ -117,7 +129,8 @@ public:
bool
visible
;
bool
visible
;
adapt_point
*
p
[
4
];
adapt_point
*
p
[
4
];
adapt_quad
*
q
[
4
];
adapt_quad
*
q
[
4
];
static
std
::
list
<
adapt_quad
*>
all_quads
;
static
std
::
list
<
adapt_quad
*>
all_elems
;
static
int
nbNod
;
};
};
class
adapt_tet
class
adapt_tet
...
@@ -134,6 +147,13 @@ public:
...
@@ -134,6 +147,13 @@ public:
t
[
4
]
=
t
[
5
]
=
t
[
6
]
=
t
[
7
]
=
0
;
t
[
4
]
=
t
[
5
]
=
t
[
6
]
=
t
[
7
]
=
0
;
}
}
inline
static
void
GSF
(
const
double
u
,
const
double
v
,
double
w
,
double
sf
[])
{
sf
[
0
]
=
1
-
u
-
v
-
w
;
sf
[
1
]
=
u
;
sf
[
2
]
=
v
;
sf
[
3
]
=
w
;
}
inline
double
V
()
const
inline
double
V
()
const
{
{
return
(
p
[
0
]
->
val
+
p
[
1
]
->
val
+
p
[
2
]
->
val
+
p
[
3
]
->
val
)
/
4.
;
return
(
p
[
0
]
->
val
+
p
[
1
]
->
val
+
p
[
2
]
->
val
+
p
[
3
]
->
val
)
/
4.
;
...
@@ -142,7 +162,6 @@ public:
...
@@ -142,7 +162,6 @@ public:
{
{
printf
(
"p1 %g %g p2 %g %g p3 %g %g
\n
"
,
p
[
0
]
->
x
,
p
[
0
]
->
y
,
p
[
1
]
->
x
,
p
[
1
]
->
y
,
p
[
2
]
->
x
,
p
[
2
]
->
y
);
printf
(
"p1 %g %g p2 %g %g p3 %g %g
\n
"
,
p
[
0
]
->
x
,
p
[
0
]
->
y
,
p
[
1
]
->
x
,
p
[
1
]
->
y
,
p
[
2
]
->
x
,
p
[
2
]
->
y
);
}
}
static
void
clean
();
static
void
Create
(
int
maxlevel
,
Double_Matrix
*
coeffs
,
Double_Matrix
*
eexps
)
;
static
void
Create
(
int
maxlevel
,
Double_Matrix
*
coeffs
,
Double_Matrix
*
eexps
)
;
static
void
Recur_Create
(
adapt_tet
*
t
,
int
maxlevel
,
int
level
,
Double_Matrix
*
coeffs
,
Double_Matrix
*
eexps
);
static
void
Recur_Create
(
adapt_tet
*
t
,
int
maxlevel
,
int
level
,
Double_Matrix
*
coeffs
,
Double_Matrix
*
eexps
);
static
void
Error
(
double
AVG
,
double
tol
);
static
void
Error
(
double
AVG
,
double
tol
);
...
@@ -150,7 +169,8 @@ public:
...
@@ -150,7 +169,8 @@ public:
bool
visible
;
bool
visible
;
adapt_point
*
p
[
4
];
adapt_point
*
p
[
4
];
adapt_tet
*
t
[
8
];
adapt_tet
*
t
[
8
];
static
std
::
list
<
adapt_tet
*>
all_tets
;
static
std
::
list
<
adapt_tet
*>
all_elems
;
static
int
nbNod
;
};
};
class
adapt_hex
class
adapt_hex
...
@@ -171,6 +191,17 @@ public:
...
@@ -171,6 +191,17 @@ public:
h
[
4
]
=
h
[
5
]
=
h
[
6
]
=
h
[
7
]
=
0
;
h
[
4
]
=
h
[
5
]
=
h
[
6
]
=
h
[
7
]
=
0
;
}
}
inline
static
void
GSF
(
const
double
u
,
const
double
v
,
double
w
,
double
sf
[])
{
sf
[
0
]
=
0.125
*
(
1
-
u
)
*
(
1
-
v
)
*
(
1
-
w
);
sf
[
1
]
=
0.125
*
(
1
+
u
)
*
(
1
-
v
)
*
(
1
-
w
);
sf
[
2
]
=
0.125
*
(
1
+
u
)
*
(
1
+
v
)
*
(
1
-
w
);
sf
[
3
]
=
0.125
*
(
1
-
u
)
*
(
1
+
v
)
*
(
1
-
w
);
sf
[
4
]
=
0.125
*
(
1
-
u
)
*
(
1
-
v
)
*
(
1
+
w
);
sf
[
5
]
=
0.125
*
(
1
+
u
)
*
(
1
-
v
)
*
(
1
+
w
);
sf
[
6
]
=
0.125
*
(
1
+
u
)
*
(
1
+
v
)
*
(
1
+
w
);
sf
[
7
]
=
0.125
*
(
1
-
u
)
*
(
1
+
v
)
*
(
1
+
w
);
}
inline
double
V
()
const
inline
double
V
()
const
{
{
return
(
p
[
0
]
->
val
+
p
[
1
]
->
val
+
p
[
2
]
->
val
+
p
[
3
]
->
val
+
p
[
4
]
->
val
+
p
[
5
]
->
val
+
p
[
6
]
->
val
+
p
[
7
]
->
val
)
/
8.
;
return
(
p
[
0
]
->
val
+
p
[
1
]
->
val
+
p
[
2
]
->
val
+
p
[
3
]
->
val
+
p
[
4
]
->
val
+
p
[
5
]
->
val
+
p
[
6
]
->
val
+
p
[
7
]
->
val
)
/
8.
;
...
@@ -179,7 +210,6 @@ public:
...
@@ -179,7 +210,6 @@ public:
{
{
printf
(
"p1 %g %g p2 %g %g p3 %g %g
\n
"
,
p
[
0
]
->
x
,
p
[
0
]
->
y
,
p
[
1
]
->
x
,
p
[
1
]
->
y
,
p
[
2
]
->
x
,
p
[
2
]
->
y
);
printf
(
"p1 %g %g p2 %g %g p3 %g %g
\n
"
,
p
[
0
]
->
x
,
p
[
0
]
->
y
,
p
[
1
]
->
x
,
p
[
1
]
->
y
,
p
[
2
]
->
x
,
p
[
2
]
->
y
);
}
}
static
void
clean
();
static
void
Create
(
int
maxlevel
,
Double_Matrix
*
coeffs
,
Double_Matrix
*
eexps
)
;
static
void
Create
(
int
maxlevel
,
Double_Matrix
*
coeffs
,
Double_Matrix
*
eexps
)
;
static
void
Recur_Create
(
adapt_hex
*
h
,
int
maxlevel
,
int
level
,
Double_Matrix
*
coeffs
,
Double_Matrix
*
eexps
);
static
void
Recur_Create
(
adapt_hex
*
h
,
int
maxlevel
,
int
level
,
Double_Matrix
*
coeffs
,
Double_Matrix
*
eexps
);
static
void
Error
(
double
AVG
,
double
tol
);
static
void
Error
(
double
AVG
,
double
tol
);
...
@@ -187,7 +217,8 @@ public:
...
@@ -187,7 +217,8 @@ public:
bool
visible
;
bool
visible
;
adapt_point
*
p
[
8
];
adapt_point
*
p
[
8
];
adapt_hex
*
h
[
8
];
adapt_hex
*
h
[
8
];
static
std
::
list
<
adapt_hex
*>
all_hexes
;
static
std
::
list
<
adapt_hex
*>
all_elems
;
static
int
nbNod
;
};
};
class
Adaptive_Post_View
class
Adaptive_Post_View
...
@@ -212,24 +243,30 @@ public:
...
@@ -212,24 +243,30 @@ public:
{
{
setAdaptiveResolutionLevel
(
view
,
level
);
setAdaptiveResolutionLevel
(
view
,
level
);
}
}
template
<
class
ELEM
>
void
setAdaptiveResolutionLevel_TEMPL
(
Post_View
*
view
,
int
level
,
GMSH_Post_Plugin
*
plug
,
List_T
**
myList
,
int
*
counter
);
void
setAdaptiveResolutionLevel
(
Post_View
*
view
,
int
levelmax
,
GMSH_Post_Plugin
*
plug
=
0
);
void
setAdaptiveResolutionLevel
(
Post_View
*
view
,
int
levelmax
,
GMSH_Post_Plugin
*
plug
=
0
);
void
initWithLowResolution
(
Post_View
*
view
);
void
initWithLowResolution
(
Post_View
*
view
);
void
setTolerance
(
const
double
eps
)
{
tol
=
eps
;}
void
setTolerance
(
const
double
eps
)
{
tol
=
eps
;}
double
getTolerance
()
const
{
return
tol
;}
double
getTolerance
()
const
{
return
tol
;}
void
zoomQuad
(
Post_View
*
view
,
template
<
class
ELEM
>
int
ielem
,
int
level
,
GMSH_Post_Plugin
*
plug
);
void
zoomElement
(
Post_View
*
view
,
void
zoomTriangle
(
Post_View
*
view
,
int
ielem
,
int
level
,
GMSH_Post_Plugin
*
plug
,
List_T
*
theList
,
int
*
counter
);
int
ielem
,
int
level
,
GMSH_Post_Plugin
*
plug
);
void
zoomTet
(
Post_View
*
view
,
int
ielem
,
int
level
,
GMSH_Post_Plugin
*
plug
,
Double_Vector
&
val
,
Double_Vector
&
res
,
Double_Matrix
&
XYZ
);
void
zoomHex
(
Post_View
*
view
,
int
ielem
,
int
level
,
GMSH_Post_Plugin
*
plug
,
Double_Vector
&
val
,
Double_Vector
&
res
,
Double_Matrix
&
XYZ
);
};
};
template
<
class
ELEM
>
void
cleanElement
()
{
typename
std
::
list
<
ELEM
*>::
iterator
it
=
ELEM
::
all_elems
.
begin
();
typename
std
::
list
<
ELEM
*>::
iterator
ite
=
ELEM
::
all_elems
.
end
();
for
(;
it
!=
ite
;
++
it
)
{
delete
*
it
;
}
ELEM
::
all_elems
.
clear
();
adapt_point
::
all_points
.
clear
();
}
#endif
#endif
This diff is collapsed.
Click to expand it.
Common/GmshMatrix.h
+
1
−
0
View file @
75a67f77
...
@@ -160,6 +160,7 @@ public:
...
@@ -160,6 +160,7 @@ public:
{
{
gsl_blas_dgemm
(
CblasNoTrans
,
CblasNoTrans
,
1.0
,
data
,
x
.
data
,
1.0
,
b
.
data
);
gsl_blas_dgemm
(
CblasNoTrans
,
CblasNoTrans
,
1.0
,
data
,
x
.
data
,
1.0
,
b
.
data
);
}
}
inline
void
mult
(
const
GSL_Vector
&
x
,
GSL_Vector
&
b
)
inline
void
mult
(
const
GSL_Vector
&
x
,
GSL_Vector
&
b
)
{
{
gsl_blas_dgemv
(
CblasNoTrans
,
1.0
,
data
,
x
.
data
,
1.0
,
b
.
data
);
gsl_blas_dgemv
(
CblasNoTrans
,
1.0
,
data
,
x
.
data
,
1.0
,
b
.
data
);
...
...
...
...
This diff is collapsed.
Click to expand it.
Plugin/Levelset.cpp
+
9
−
9
View file @
75a67f77
// $Id: Levelset.cpp,v 1.2
4
2005-0
1-03 04:09:27 geuzain
e Exp $
// $Id: Levelset.cpp,v 1.2
5
2005-0
3-16 16:44:38 remacl
e Exp $
//
//
// Copyright (C) 1997-2005 C. Geuzaine, J.-F. Remacle
// Copyright (C) 1997-2005 C. Geuzaine, J.-F. Remacle
//
//
...
@@ -916,24 +916,24 @@ static bool recur_sign_change (adapt_quad *q, double val, const GMSH_LevelsetPlu
...
@@ -916,24 +916,24 @@ static bool recur_sign_change (adapt_quad *q, double val, const GMSH_LevelsetPlu
void
GMSH_LevelsetPlugin
::
assign_specific_visibility
()
const
void
GMSH_LevelsetPlugin
::
assign_specific_visibility
()
const
{
{
if
(
adapt_triangle
::
all_
triang
les
.
size
())
if
(
adapt_triangle
::
all_
e
le
m
s
.
size
())
{
{
adapt_triangle
*
t
=
*
adapt_triangle
::
all_
triang
les
.
begin
();
adapt_triangle
*
t
=
*
adapt_triangle
::
all_
e
le
m
s
.
begin
();
t
->
visible
=
!
recur_sign_change
(
t
,
_valueView
,
this
);
t
->
visible
=
!
recur_sign_change
(
t
,
_valueView
,
this
);
}
}
if
(
adapt_tet
::
all_
tet
s
.
size
())
if
(
adapt_tet
::
all_
elem
s
.
size
())
{
{
adapt_tet
*
te
=
*
adapt_tet
::
all_
tet
s
.
begin
();
adapt_tet
*
te
=
*
adapt_tet
::
all_
elem
s
.
begin
();
te
->
visible
=
!
recur_sign_change
(
te
,
_valueView
,
this
);
te
->
visible
=
!
recur_sign_change
(
te
,
_valueView
,
this
);
}
}
if
(
adapt_quad
::
all_
quad
s
.
size
())
if
(
adapt_quad
::
all_
elem
s
.
size
())
{
{
adapt_quad
*
qe
=
*
adapt_quad
::
all_
quad
s
.
begin
();
adapt_quad
*
qe
=
*
adapt_quad
::
all_
elem
s
.
begin
();
qe
->
visible
=
!
recur_sign_change
(
qe
,
_valueView
,
this
);
qe
->
visible
=
!
recur_sign_change
(
qe
,
_valueView
,
this
);
}
}
if
(
adapt_hex
::
all_
hexe
s
.
size
())
if
(
adapt_hex
::
all_
elem
s
.
size
())
{
{
adapt_hex
*
he
=
*
adapt_hex
::
all_
hexe
s
.
begin
();
adapt_hex
*
he
=
*
adapt_hex
::
all_
elem
s
.
begin
();
he
->
visible
=
!
recur_sign_change
(
he
,
_valueView
,
this
);
he
->
visible
=
!
recur_sign_change
(
he
,
_valueView
,
this
);
}
}
}
}
This diff is collapsed.
Click to expand it.
Preview
0%
Loading
Try again
or
attach a new file
.
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Save comment
Cancel
Please
register
or
sign in
to comment