Skip to content
GitLab
Menu
Projects
Groups
Snippets
/
Help
Help
Support
Community forum
Keyboard shortcuts
?
Submit feedback
Contribute to GitLab
Sign in
Toggle navigation
Menu
Open sidebar
Backofen, Rainer
amdis
Commits
06e1522f
Commit
06e1522f
authored
Jun 25, 2008
by
Thomas Witkowski
Browse files
* Parallel file i/o
parent
35c63e58
Changes
19
Show whitespace changes
Inline
Side-by-side
AMDiS/src/AdaptInstationary.cc
View file @
06e1522f
...
@@ -204,7 +204,21 @@ namespace AMDiS {
...
@@ -204,7 +204,21 @@ namespace AMDiS {
iterationTimestamp_
=
time
(
NULL
);
iterationTimestamp_
=
time
(
NULL
);
problemTime_
->
initTimestep
(
adaptInfo_
);
problemTime_
->
initTimestep
(
adaptInfo_
);
#ifdef _OPENMP
#pragma omp parallel sections
{
#pragma omp section
problemTime_
->
startDelayedTimestepCalculation
();
#pragma omp section
oneTimestep
();
oneTimestep
();
}
#else
problemTime_
->
startDelayedTimestepCalculation
();
oneTimestep
();
#endif
problemTime_
->
closeTimestep
(
adaptInfo_
);
problemTime_
->
closeTimestep
(
adaptInfo_
);
if
(
breakWhenStable
&&
(
adaptInfo_
->
getSolverIterations
()
==
0
))
{
if
(
breakWhenStable
&&
(
adaptInfo_
->
getSolverIterations
()
==
0
))
{
...
@@ -220,6 +234,8 @@ namespace AMDiS {
...
@@ -220,6 +234,8 @@ namespace AMDiS {
}
}
}
}
problemTime_
->
startDelayedTimestepCalculation
();
return
errorCode
;
return
errorCode
;
}
}
...
...
AMDiS/src/DataCollector.cc
View file @
06e1522f
...
@@ -56,6 +56,19 @@ namespace AMDiS {
...
@@ -56,6 +56,19 @@ namespace AMDiS {
DELETE
dofCoords_
;
DELETE
dofCoords_
;
}
}
void
DataCollector
::
fillAllData
()
{
if
(
!
elementDataCollected_
)
{
startCollectingElementData
();
}
if
(
!
periodicDataCollected_
)
{
startCollectingPeriodicData
();
}
if
(
!
valueDataCollected_
)
{
startCollectingValueData
();
}
}
int
DataCollector
::
startCollectingElementData
()
int
DataCollector
::
startCollectingElementData
()
{
{
FUNCNAME
(
"DataCollector::startCollectingElementData()"
);
FUNCNAME
(
"DataCollector::startCollectingElementData()"
);
...
@@ -427,7 +440,7 @@ namespace AMDiS {
...
@@ -427,7 +440,7 @@ namespace AMDiS {
::
std
::
list
<
PeriodicInfo
>*
DataCollector
::
getPeriodicInfos
()
::
std
::
list
<
PeriodicInfo
>*
DataCollector
::
getPeriodicInfos
()
{
{
if
(
!
periodicDataCollected_
)
{
if
(
!
periodicDataCollected_
)
{
startCollectingPeriodicData
();
startCollectingPeriodicData
();
}
}
...
@@ -463,7 +476,7 @@ namespace AMDiS {
...
@@ -463,7 +476,7 @@ namespace AMDiS {
int
DataCollector
::
getNumberConnections
()
int
DataCollector
::
getNumberConnections
()
{
{
if
(
!
periodicDataCollected_
)
{
if
(
!
periodicDataCollected_
)
{
startCollectingPeriodicData
();
startCollectingPeriodicData
();
}
}
...
...
AMDiS/src/DataCollector.h
View file @
06e1522f
...
@@ -61,6 +61,11 @@ namespace AMDiS {
...
@@ -61,6 +61,11 @@ namespace AMDiS {
~
DataCollector
();
~
DataCollector
();
/** \brief
* Fills the DataCollector with all possible datas.
*/
void
fillAllData
();
/** \brief
/** \brief
* Returns list with element information.
* Returns list with element information.
*/
*/
...
...
AMDiS/src/DuneSolver.cc
View file @
06e1522f
...
@@ -71,6 +71,13 @@ namespace AMDiS {
...
@@ -71,6 +71,13 @@ namespace AMDiS {
{
{
FUNCNAME
(
"DuneSolver::solveSystem()"
);
FUNCNAME
(
"DuneSolver::solveSystem()"
);
// Calculate residual reduction which should be achieved by the dune solver.
*
p
=
*
x
;
*
p
*=
-
1.0
;
mv
->
matVec
(
NoTranspose
,
*
p
,
*
r
);
*
r
+=
*
b
;
double
reduction
=
this
->
tolerance
/
norm
(
r
);
typedef
Dune
::
Matrix
<
DuneMatrix
>
SystemDuneMatrix
;
typedef
Dune
::
Matrix
<
DuneMatrix
>
SystemDuneMatrix
;
typedef
Dune
::
BlockVector
<
DuneVector
>
SystemDuneVector
;
typedef
Dune
::
BlockVector
<
DuneVector
>
SystemDuneVector
;
...
@@ -117,7 +124,7 @@ namespace AMDiS {
...
@@ -117,7 +124,7 @@ namespace AMDiS {
Dune
::
InverseOperator
<
SystemDuneVector
,
SystemDuneVector
>
*
duneSolver
=
Dune
::
InverseOperator
<
SystemDuneVector
,
SystemDuneVector
>
*
duneSolver
=
getSolver
<
SystemDuneMatrix
,
SystemDuneVector
>
(
&
op
,
dunePreconditioner
);
getSolver
<
SystemDuneMatrix
,
SystemDuneVector
>
(
&
op
,
dunePreconditioner
);
duneSolver
->
apply
(
duneVecX
,
duneVecB
,
opRes
);
duneSolver
->
apply
(
duneVecX
,
duneVecB
,
reduction
,
opRes
);
for
(
int
i
=
0
;
i
<
nComponents
;
i
++
)
{
for
(
int
i
=
0
;
i
<
nComponents
;
i
++
)
{
DuneVector
::
Iterator
duneVecIt
;
DuneVector
::
Iterator
duneVecIt
;
...
@@ -129,6 +136,16 @@ namespace AMDiS {
...
@@ -129,6 +136,16 @@ namespace AMDiS {
}
}
}
}
// Calculate and print the residual.
*
p
=
*
x
;
*
p
*=
-
1.0
;
mv
->
matVec
(
NoTranspose
,
*
p
,
*
r
);
*
r
+=
*
b
;
this
->
residual
=
norm
(
r
);
MSG
(
"Residual: %e
\n
"
,
this
->
residual
);
delete
dunePreconditioner
;
delete
dunePreconditioner
;
delete
duneSolver
;
delete
duneSolver
;
...
...
AMDiS/src/DuneSolver.h
View file @
06e1522f
...
@@ -123,12 +123,20 @@ namespace AMDiS {
...
@@ -123,12 +123,20 @@ namespace AMDiS {
/** \brief
/** \brief
* Implements OEMSolver<VectorType>::init().
* Implements OEMSolver<VectorType>::init().
*/
*/
void
init
()
{};
void
init
()
{
p
=
this
->
vectorCreator
->
create
();
r
=
this
->
vectorCreator
->
create
();
};
/** \brief
/** \brief
* Implements OEMSolver<VectorType>::exit().
* Implements OEMSolver<VectorType>::exit().
*/
*/
void
exit
()
{};
void
exit
()
{
this
->
vectorCreator
->
free
(
p
);
this
->
vectorCreator
->
free
(
r
);
};
/** \brief
/** \brief
* Creates a dune preconditioner using the informations in the used init file.
* Creates a dune preconditioner using the informations in the used init file.
...
@@ -168,6 +176,12 @@ namespace AMDiS {
...
@@ -168,6 +176,12 @@ namespace AMDiS {
* to be allocated before.
* to be allocated before.
*/
*/
void
mapDOFVector
(
DOFVector
<
double
>
*
dofVector
,
DuneVector
*
duneVector
);
void
mapDOFVector
(
DOFVector
<
double
>
*
dofVector
,
DuneVector
*
duneVector
);
private:
/** \brief
* These vectors are justed to calculate the final residual of the solution.
*/
VectorType
*
r
,
*
p
;
};
};
}
}
...
...
AMDiS/src/DuneSolver.hh
View file @
06e1522f
AMDiS/src/ElementFileWriter.h
View file @
06e1522f
...
@@ -30,6 +30,8 @@ namespace AMDiS {
...
@@ -30,6 +30,8 @@ namespace AMDiS {
Flag
traverseFlag
=
Mesh
::
CALL_LEAF_EL
,
Flag
traverseFlag
=
Mesh
::
CALL_LEAF_EL
,
bool
(
*
writeElem
)(
ElInfo
*
)
=
NULL
);
bool
(
*
writeElem
)(
ElInfo
*
)
=
NULL
);
void
writeDelayedFiles
()
{};
protected:
protected:
/**
/**
* Writes element data in tecplot format.
* Writes element data in tecplot format.
...
...
AMDiS/src/FileWriter.cc
View file @
06e1522f
...
@@ -5,7 +5,6 @@
...
@@ -5,7 +5,6 @@
#include
"ValueWriter.h"
#include
"ValueWriter.h"
#include
"MacroWriter.h"
#include
"MacroWriter.h"
#include
"VtkWriter.h"
#include
"VtkWriter.h"
#include
"DataCollector.h"
#include
"FiniteElemSpace.h"
#include
"FiniteElemSpace.h"
#include
"AdaptInfo.h"
#include
"AdaptInfo.h"
#include
"Flag.h"
#include
"Flag.h"
...
@@ -28,6 +27,8 @@ namespace AMDiS {
...
@@ -28,6 +27,8 @@ namespace AMDiS {
solutionVecs_
.
resize
(
1
);
solutionVecs_
.
resize
(
1
);
solutionVecs_
[
0
]
=
vec
;
solutionVecs_
[
0
]
=
vec
;
dataCollectors_
.
resize
(
1
);
}
}
...
@@ -48,6 +49,8 @@ namespace AMDiS {
...
@@ -48,6 +49,8 @@ namespace AMDiS {
feSpace
=
vecs
[
0
]
->
getFESpace
();
feSpace
=
vecs
[
0
]
->
getFESpace
();
solutionVecs_
=
vecs
;
solutionVecs_
=
vecs
;
dataCollectors_
.
resize
(
vecs
.
size
());
}
}
...
@@ -79,6 +82,8 @@ namespace AMDiS {
...
@@ -79,6 +82,8 @@ namespace AMDiS {
}
}
feSpace
=
vec
->
getFESpace
();
feSpace
=
vec
->
getFESpace
();
dataCollectors_
.
resize
(
nTmpSolutions_
);
}
}
...
@@ -91,6 +96,12 @@ namespace AMDiS {
...
@@ -91,6 +96,12 @@ namespace AMDiS {
DELETE
solutionVecs_
[
i
];
DELETE
solutionVecs_
[
i
];
}
}
}
}
for
(
int
i
=
0
;
i
<
static_cast
<
int
>
(
dataCollectors_
.
size
());
i
++
)
{
if
(
dataCollectors_
[
i
])
{
DELETE
dataCollectors_
[
i
];
}
}
}
}
...
@@ -111,6 +122,9 @@ namespace AMDiS {
...
@@ -111,6 +122,9 @@ namespace AMDiS {
indexDecimals
=
3
;
indexDecimals
=
3
;
tsModulo
=
1
;
tsModulo
=
1
;
nTmpSolutions_
=
0
;
nTmpSolutions_
=
0
;
delayWriting_
=
0
;
writingIsDelayed_
=
false
;
delayedFilename_
=
""
;
paraViewAnimationFrames_
.
resize
(
0
);
paraViewAnimationFrames_
.
resize
(
0
);
readParameters
();
readParameters
();
...
@@ -135,6 +149,7 @@ namespace AMDiS {
...
@@ -135,6 +149,7 @@ namespace AMDiS {
GET_PARAMETER
(
0
,
name
+
"->index length"
,
"%d"
,
&
indexLength
);
GET_PARAMETER
(
0
,
name
+
"->index length"
,
"%d"
,
&
indexLength
);
GET_PARAMETER
(
0
,
name
+
"->index decimals"
,
"%d"
,
&
indexDecimals
);
GET_PARAMETER
(
0
,
name
+
"->index decimals"
,
"%d"
,
&
indexDecimals
);
GET_PARAMETER
(
0
,
name
+
"->write every i-th timestep"
,
"%d"
,
&
tsModulo
);
GET_PARAMETER
(
0
,
name
+
"->write every i-th timestep"
,
"%d"
,
&
tsModulo
);
GET_PARAMETER
(
0
,
name
+
"->delay"
,
"%d"
,
&
delayWriting_
);
}
}
void
FileWriter
::
writeFiles
(
AdaptInfo
*
adaptInfo
,
void
FileWriter
::
writeFiles
(
AdaptInfo
*
adaptInfo
,
...
@@ -148,6 +163,22 @@ namespace AMDiS {
...
@@ -148,6 +163,22 @@ namespace AMDiS {
if
((
adaptInfo
->
getTimestepNumber
()
%
tsModulo
!=
0
)
&&
!
force
)
if
((
adaptInfo
->
getTimestepNumber
()
%
tsModulo
!=
0
)
&&
!
force
)
return
;
return
;
if
(
writingIsDelayed_
)
{
ERROR_EXIT
(
"This should not happen!
\n
"
);
}
if
(
writeElem
)
{
for
(
int
i
=
0
;
i
<
static_cast
<
int
>
(
dataCollectors_
.
size
());
i
++
)
{
dataCollectors_
[
i
]
=
NEW
DataCollector
(
feSpace
,
solutionVecs_
[
i
],
level
,
flag
,
writeElem
);
}
}
else
{
for
(
int
i
=
0
;
i
<
static_cast
<
int
>
(
dataCollectors_
.
size
());
i
++
)
{
dataCollectors_
[
i
]
=
NEW
DataCollector
(
feSpace
,
solutionVecs_
[
i
],
traverseLevel
,
flag
|
traverseFlag
,
writeElement
);
}
}
::
std
::
string
fn
=
filename
;
::
std
::
string
fn
=
filename
;
if
(
appendIndex
)
{
if
(
appendIndex
)
{
...
@@ -156,30 +187,28 @@ namespace AMDiS {
...
@@ -156,30 +187,28 @@ namespace AMDiS {
TEST_EXIT
(
indexDecimals
<
indexLength
)(
"index length <= index decimals
\n
"
);
TEST_EXIT
(
indexDecimals
<
indexLength
)(
"index length <= index decimals
\n
"
);
char
formatStr
[
9
];
char
formatStr
[
9
];
sprintf
(
formatStr
,
"%%0%d.%df"
,
indexLength
,
indexDecimals
);
char
timeStr
[
20
];
char
timeStr
[
20
];
sprintf
(
formatStr
,
"%%0%d.%df"
,
indexLength
,
indexDecimals
);
sprintf
(
timeStr
,
formatStr
,
adaptInfo
?
adaptInfo
->
getTime
()
:
0.0
);
sprintf
(
timeStr
,
formatStr
,
adaptInfo
?
adaptInfo
->
getTime
()
:
0.0
);
fn
+=
timeStr
;
fn
+=
timeStr
;
}
}
::
std
::
vector
<
DataCollector
*
>
dc
(
solutionVecs_
.
size
());
if
(
writeElem
)
{
if
(
delayWriting_
)
{
for
(
int
i
=
0
;
i
<
static_cast
<
int
>
(
dc
.
size
());
i
++
)
{
if
(
writeTecPlotFormat
||
writeAMDiSFormat
||
writePeriodicFormat
)
{
dc
[
i
]
=
NEW
DataCollector
(
feSpace
,
solutionVecs_
[
i
],
ERROR_EXIT
(
"Delay writing only supported for ParaView file format!
\n
"
);
level
,
flag
,
writeElem
);
}
}
else
{
for
(
int
i
=
0
;
i
<
static_cast
<
int
>
(
dc
.
size
());
i
++
)
{
dc
[
i
]
=
NEW
DataCollector
(
feSpace
,
solutionVecs_
[
i
],
traverseLevel
,
flag
|
traverseFlag
,
writeElement
);
}
}
for
(
int
i
=
0
;
i
<
static_cast
<
int
>
(
dataCollectors_
.
size
());
i
++
)
{
dataCollectors_
[
i
]
->
fillAllData
();
}
}
writingIsDelayed_
=
true
;
delayedFilename_
=
fn
;
return
;
}
if
(
writeTecPlotFormat
)
{
if
(
writeTecPlotFormat
)
{
TecPlotWriter
<
DOFVector
<
double
>
>::
writeValues
(
solutionVecs_
[
0
],
TecPlotWriter
<
DOFVector
<
double
>
>::
writeValues
(
solutionVecs_
[
0
],
...
@@ -191,42 +220,67 @@ namespace AMDiS {
...
@@ -191,42 +220,67 @@ namespace AMDiS {
if
(
writeAMDiSFormat
)
{
if
(
writeAMDiSFormat
)
{
TEST_EXIT
(
mesh
)(
"no mesh
\n
"
);
TEST_EXIT
(
mesh
)(
"no mesh
\n
"
);
MacroWriter
::
writeMacro
(
d
c
[
0
],
MacroWriter
::
writeMacro
(
d
ataCollectors_
[
0
],
const_cast
<
char
*>
(
(
fn
+
amdisMeshExt
).
c_str
()),
const_cast
<
char
*>
((
fn
+
amdisMeshExt
).
c_str
()),
adaptInfo
?
adaptInfo
->
getTime
()
:
0.0
);
adaptInfo
?
adaptInfo
->
getTime
()
:
0.0
);
MSG
(
"macro file written to %s
\n
"
,
(
fn
+
amdisMeshExt
).
c_str
());
MSG
(
"macro file written to %s
\n
"
,
(
fn
+
amdisMeshExt
).
c_str
());
ValueWriter
::
writeValues
(
d
c
[
0
],
ValueWriter
::
writeValues
(
d
ataCollectors_
[
0
],
(
fn
+
amdisDataExt
).
c_str
(),
(
fn
+
amdisDataExt
).
c_str
(),
adaptInfo
?
adaptInfo
->
getTime
()
:
0.0
);
adaptInfo
?
adaptInfo
->
getTime
()
:
0.0
);
MSG
(
"value file written to %s
\n
"
,
(
fn
+
amdisDataExt
).
c_str
());
MSG
(
"value file written to %s
\n
"
,
(
fn
+
amdisDataExt
).
c_str
());
}
}
if
(
writePeriodicFormat
)
{
if
(
writePeriodicFormat
)
{
MacroWriter
::
writePeriodicFile
(
d
c
[
0
],
MacroWriter
::
writePeriodicFile
(
d
ataCollectors_
[
0
],
(
fn
+
periodicFileExt
).
c_str
());
(
fn
+
periodicFileExt
).
c_str
());
MSG
(
"periodic file written to %s
\n
"
,
(
fn
+
periodicFileExt
).
c_str
());
MSG
(
"periodic file written to %s
\n
"
,
(
fn
+
periodicFileExt
).
c_str
());
}
}
if
(
writeParaViewFormat
)
{
if
(
writeParaViewFormat
)
{
VtkWriter
vtkWriter
(
&
dc
);
VtkWriter
vtkWriter
(
&
dataCollectors_
);
vtkWriter
.
writeFile
(
const_cast
<
char
*>
((
fn
+
paraViewFileExt
).
c_str
()));
vtkWriter
.
writeFile
(
const_cast
<
char
*>
(
(
fn
+
paraViewFileExt
).
c_str
()));
MSG
(
"ParaView file written to %s
\n
"
,
(
fn
+
paraViewFileExt
).
c_str
());
MSG
(
"ParaView file written to %s
\n
"
,
(
fn
+
paraViewFileExt
).
c_str
());
}
}
if
(
writeParaViewAnimation
)
{
if
(
writeParaViewAnimation
)
{
VtkWriter
vtkWriter
(
&
dc
);
VtkWriter
vtkWriter
(
&
dataCollectors_
);
vtkWriter
.
updateAnimationFile
(
fn
+
paraViewFileExt
,
vtkWriter
.
updateAnimationFile
(
fn
+
paraViewFileExt
,
&
paraViewAnimationFrames_
,
&
paraViewAnimationFrames_
,
const_cast
<
char
*>
(
(
filename
+
".pvd"
).
c_str
()));
const_cast
<
char
*>
((
filename
+
".pvd"
).
c_str
()));
}
for
(
int
i
=
0
;
i
<
static_cast
<
int
>
(
dataCollectors_
.
size
());
i
++
)
{
DELETE
dataCollectors_
[
i
];
}
}
}
void
FileWriter
::
writeDelayedFiles
()
{
if
(
!
writingIsDelayed_
)
{
return
;
}
for
(
int
i
=
0
;
i
<
static_cast
<
int
>
(
dc
.
size
());
i
++
)
{
if
(
writeParaViewFormat
)
{
DELETE
dc
[
i
];
VtkWriter
vtkWriter
(
&
dataCollectors_
);
vtkWriter
.
writeFile
(
const_cast
<
char
*>
((
delayedFilename_
+
paraViewFileExt
).
c_str
()));
}
}
if
(
writeParaViewAnimation
)
{
VtkWriter
vtkWriter
(
&
dataCollectors_
);
vtkWriter
.
updateAnimationFile
(
delayedFilename_
+
paraViewFileExt
,
&
paraViewAnimationFrames_
,
const_cast
<
char
*>
((
filename
+
".pvd"
).
c_str
()));
}
for
(
int
i
=
0
;
i
<
static_cast
<
int
>
(
dataCollectors_
.
size
());
i
++
)
{
DELETE
dataCollectors_
[
i
];
}
writingIsDelayed_
=
false
;
delayedFilename_
=
""
;
}
}
}
}
AMDiS/src/FileWriter.h
View file @
06e1522f
...
@@ -31,6 +31,7 @@
...
@@ -31,6 +31,7 @@
#include
"MemoryManager.h"
#include
"MemoryManager.h"
#include
"MatrixVector.h"
#include
"MatrixVector.h"
#include
"Mesh.h"
#include
"Mesh.h"
#include
"DataCollector.h"
namespace
AMDiS
{
namespace
AMDiS
{
...
@@ -69,6 +70,8 @@ namespace AMDiS {
...
@@ -69,6 +70,8 @@ namespace AMDiS {
Flag
traverseFlag
=
Mesh
::
CALL_LEAF_EL
,
Flag
traverseFlag
=
Mesh
::
CALL_LEAF_EL
,
bool
(
*
writeElem
)(
ElInfo
*
)
=
NULL
)
=
0
;
bool
(
*
writeElem
)(
ElInfo
*
)
=
NULL
)
=
0
;
virtual
void
writeDelayedFiles
()
=
0
;
void
setTraverseProperties
(
int
level
,
void
setTraverseProperties
(
int
level
,
Flag
flag
,
Flag
flag
,
bool
(
*
writeElem
)(
ElInfo
*
))
bool
(
*
writeElem
)(
ElInfo
*
))
...
@@ -150,6 +153,8 @@ namespace AMDiS {
...
@@ -150,6 +153,8 @@ namespace AMDiS {
Flag
traverseFlag
=
Mesh
::
CALL_LEAF_EL
,
Flag
traverseFlag
=
Mesh
::
CALL_LEAF_EL
,
bool
(
*
writeElem
)(
ElInfo
*
)
=
NULL
);
bool
(
*
writeElem
)(
ElInfo
*
)
=
NULL
);
virtual
void
writeDelayedFiles
();
protected:
protected:
/** \brief
/** \brief
* Initialization of the filewriter.
* Initialization of the filewriter.
...
@@ -272,6 +277,28 @@ namespace AMDiS {
...
@@ -272,6 +277,28 @@ namespace AMDiS {
* stored in solutionVecs_ must be deleted in the destructor.
* stored in solutionVecs_ must be deleted in the destructor.
*/
*/
int
nTmpSolutions_
;
int
nTmpSolutions_
;
/** \brief
* Containers, which store the data to be written;
*/
::
std
::
vector
<
DataCollector
*
>
dataCollectors_
;
/** \brief
* If set to 1, the FileWriter will delay the file writing to the future, where
* it can be executed in parallel with some other independent calculations.
*/
int
delayWriting_
;
/** \brief
* If set to true, the FileWriter was filled with data, but the Files are currently
* not written to disk.
*/
bool
writingIsDelayed_
;
/** \brief
* Here the filename for the file, which should be written to the next, is stored.
*/
::
std
::
string
delayedFilename_
;
};
};
}
}
...
...
AMDiS/src/ProblemInstat.cc
View file @
06e1522f
...
@@ -95,7 +95,7 @@ namespace AMDiS {
...
@@ -95,7 +95,7 @@ namespace AMDiS {
oldSolution
=
NEW
DOFVector
<
double
>
(
problemStat
->
getFESpace
(),
name
+
"->uOld"
);
oldSolution
=
NEW
DOFVector
<
double
>
(
problemStat
->
getFESpace
(),
name
+
"->uOld"
);
oldSolution
->
refineInterpol
(
true
);
oldSolution
->
refineInterpol
(
true
);
oldSolution
->
setCoarsenOperation
(
COARSE_INTERPOL
);
oldSolution
->
setCoarsenOperation
(
COARSE_INTERPOL
);
if
(
problemStat
->
getEstimator
())
{
if
(
problemStat
->
getEstimator
())
{
dynamic_cast
<
Estimator
*>
(
problemStat
->
getEstimator
())
dynamic_cast
<
Estimator
*>
(
problemStat
->
getEstimator
())
->
addUhOldToSystem
(
0
,
oldSolution
);
->
addUhOldToSystem
(
0
,
oldSolution
);
}
}
...
@@ -184,4 +184,14 @@ namespace AMDiS {
...
@@ -184,4 +184,14 @@ namespace AMDiS {
oldSolution
->
copy
(
*
(
problemStat
->
getSolution
()));
oldSolution
->
copy
(
*
(
problemStat
->
getSolution
()));
}
}
void
ProblemInstatScal
::
startDelayedTimestepCalculation
()
{
problemStat
->
writeDelayedFiles
();
}
void
ProblemInstatVec
::
startDelayedTimestepCalculation
()
{
problemStat
->
writeDelayedFiles
();
}
}
}
AMDiS/src/ProblemInstat.h
View file @
06e1522f
...
@@ -121,6 +121,11 @@ namespace AMDiS {
...
@@ -121,6 +121,11 @@ namespace AMDiS {
*/
*/
virtual
void
solveInitialProblem
(
AdaptInfo
*
adaptInfo
);
virtual
void
solveInitialProblem
(
AdaptInfo
*
adaptInfo
);