Ice Sheet System Model  4.18
Code documentation
FemModel.h
Go to the documentation of this file.
1 /*
2  * FemModel.h:
3  */
4 
5 #ifndef _FEMMODEL_H_
6 #define _FEMMODEL_H_
7 
8 /*Headers:*/
9 /*{{{*/
10 #include "../toolkits/toolkits.h"
11 class DataSet;
12 class Parameters;
13 class Inputs2;
14 class Nodes;
15 class Vertices;
16 class Results;
17 class Constraints;
18 class Loads;
19 class Materials;
20 class SealevelMasks;
21 class Profiler;
22 class Elements;
23 #if defined(_HAVE_NEOPZ_) && !defined(_HAVE_AD_)
25 #endif
26 #if defined(_HAVE_BAMG_) && !defined(_HAVE_AD_)
27 #include "./AmrBamg.h"
28 #endif
29 /*}}}*/
30 
31 class FemModel {
32 
33  /*no private members, as we need access to these datasets quite often!:*/
34 
35  public:
36 
37  int analysis_counter; //counter into analysis_type_list
38  int *analysis_type_list; //list of analyses this femmodel is going to carry out
39  int nummodels;
41 
42  Profiler* profiler; //keep time, cpu and mem statistics while we are running.
43 
44  Elements *elements; //elements (one set for all analyses)
45  Materials *materials; //one set of materials, for each element
46  Parameters *parameters; //one set of parameters, independent of the analysis_type
47  Inputs2 *inputs2; //one set of inputs, independent of the analysis_type
48  Results *results; //results that cannot be fit into the elements
49  Vertices *vertices; //one set of vertices
50 
51  /*Analysis dependent datasets*/
58 
59  //FIXME: do we want only one class and have virtual functions? or keep 2 classes, at least rename AdaptiveMeshRefinement -> AmrNeopz
60  #if defined(_HAVE_NEOPZ_) && !defined(_HAVE_AD_)
61  AdaptiveMeshRefinement *amr; //adaptive mesh refinement object. It keeps coarse mesh and execute refinement process
62  #endif
63 
64  #if defined(_HAVE_BAMG_) && !defined(_HAVE_AD_)
65  AmrBamg *amrbamg; //adaptive mesh refinement object. It keeps coarse mesh and execute refinement process
66  #endif
67 
68  /*constructors, destructors: */
69  FemModel(int argc,char** argv,ISSM_MPI_Comm comm_init,bool trace=false);
70  FemModel(char* rootpath, char* inputfilename, char* outputfilename, char* toolkitsfilename, char* lockfilename, char* restartfilename, ISSM_MPI_Comm incomm, int solution_type,IssmPDouble* X);
71  ~FemModel();
72 
73  /*Methods:*/
74  int AnalysisIndex(int);
75  void CheckPoint(void);
76  void CleanUp(void);
77  FemModel* copy();
78  void Echo();
79  int GetElementsWidth(){return 3;};//just tria elements in this first version
80  void InitFromFiles(char* rootpath, char* inputfilename, char* outputfilename, char* petscfilename, char* lockfilename, char* restartfilename, const int solution_type,bool trace,IssmPDouble* X=NULL);
81  void InitFromFids(char* rootpath, FILE* IOMODEL, FILE* toolkitsoptionsfid, int in_solution_type, bool trace, IssmPDouble* X=NULL);
82  void Marshall(char** pmarshalled_data, int* pmarshalled_data_size, int marshall_direction);
83  void Restart(void);
84  void SetCurrentConfiguration(int configuration_type);
85  void SetCurrentConfiguration(int configuration_type,int analysis_type);
86  int Size(void);
87  void SolutionAnalysesList(int** panalyses,int* pnumanalyses,IoModel* iomodel,int solutiontype);
88  void Solve(void);
89 
90  /*Modules*/
92  void CalvingRateVonmisesx();
93  void CalvingRateLevermannx();
94  void CalvingFluxLevelsetx();
96  void DeviatoricStressx();
97  void Divergencex(IssmDouble* pdiv);
98  void ElementOperationx(void (Element::*function)(void));
99  void ElementResponsex(IssmDouble* presponse,int response_enum);
100  void FloatingAreax(IssmDouble* pV, bool scaled);
103  void GetLocalVectorWithClonesVertices(IssmDouble** plocal_vector,Vector<IssmDouble> *vector);
104  void SyncLocalVectorWithClonesVertices(IssmDouble* local_vector);
106  void GetLocalVectorWithClonesNodes(IssmDouble** plocal_vector,Vector<IssmDouble> *vector);
107  void GroundedAreax(IssmDouble* pV, bool scaled);
108  void IcefrontAreax();
109  void IcefrontMassFluxx(IssmDouble* presponse, bool scaled);
110  void IcefrontMassFluxLevelsetx(IssmDouble* presponse, bool scaled);
111  void GroundinglineMassFluxx(IssmDouble* presponse, bool scaled);
112  void IceMassx(IssmDouble* pV, bool scaled);
113  void IceVolumex(IssmDouble* pV, bool scaled);
114  void IceVolumeAboveFloatationx(IssmDouble* pV, bool scaled);
115  void InputMakeDiscontinuous(int enum_in);
116  void MassFluxx(IssmDouble* presponse);
117  void MaxAbsVxx(IssmDouble* presponse);
118  void MaxAbsVyx(IssmDouble* presponse);
119  void MaxAbsVzx(IssmDouble* presponse);
120  void MaxDivergencex(IssmDouble* pdiv);
121  void MaxVelx(IssmDouble* presponse);
122  void MaxVxx(IssmDouble* presponse);
123  void MaxVyx(IssmDouble* presponse);
124  void MaxVzx(IssmDouble* presponse);
125  void MinVelx(IssmDouble* presponse);
126  void MinVxx(IssmDouble* presponse);
127  void MinVyx(IssmDouble* presponse);
128  void MinVzx(IssmDouble* presponse);
129  void DistanceToFieldValue(int fieldenum,IssmDouble fieldvalue,int distanceenum);
130  void ResetLevelset();
131  void StrainRateparallelx();
133  void StrainRateeffectivex();
134  void StressIntensityFactorx();
136  void TotalCalvingFluxLevelsetx(IssmDouble* pGbmb, bool scaled);
137  void TotalCalvingMeltingFluxLevelsetx(IssmDouble* pGbmb, bool scaled);
138  void TotalFloatingBmbx(IssmDouble* pFbmb, bool scaled);
139  void TotalGroundedBmbx(IssmDouble* pGbmb, bool scaled);
140  void TotalSmbx(IssmDouble* pSmb, bool scaled);
141  #ifdef _HAVE_DAKOTA_
142  void DakotaResponsesx(double* d_responses,char** responses_descriptors,int numresponsedescriptors,int d_numresponses);
143  #endif
144  void CostFunctionx(IssmDouble* pJ,IssmDouble** pJlist,int* pn);
145  void OutputControlsx(Results **presults);
146  void RequestedDependentsx(void);
147  void RequestedOutputsx(Results **presults,char** requested_outputs, int numoutputs,bool save_results=true);
148  void RequestedOutputsx(Results **presults,int* requested_outputs, int numoutputs,bool save_results=true);
149  void Responsex(IssmDouble* presponse,int response_descriptor_enum);
150  void Responsex(IssmDouble* presponse,const char* response_descriptor);
151  void SurfaceAbsMisfitx( IssmDouble* pJ);
152  void OmegaAbsGradientx( IssmDouble* pJ);
153  void EtaDiffx( IssmDouble* pJ);
154  void ThicknessAverage();
156  void ThicknessPositivex(IssmDouble* pJ);
157  #ifdef _HAVE_GIA_
158  void Deflection(Vector<IssmDouble>* wg,Vector<IssmDouble>* dwgdt, IssmDouble* x, IssmDouble* y);
159  #endif
160  #ifdef _HAVE_ESA_
161  void EsaGeodetic2D(Vector<IssmDouble>* pUp, Vector<IssmDouble>* pNorth, Vector<IssmDouble>* pEast, Vector<IssmDouble>* pX, Vector<IssmDouble>* pY, IssmDouble* xx, IssmDouble* yy);
162  void EsaGeodetic3D(Vector<IssmDouble>* pUp, Vector<IssmDouble>* pNorth, Vector<IssmDouble>* pEast, IssmDouble* latitude, IssmDouble* longitude, IssmDouble* radius, IssmDouble* xx, IssmDouble* yy, IssmDouble* zz);
163  #endif
164  #ifdef _HAVE_SEALEVELRISE_
165  void SealevelriseEustatic(Vector<IssmDouble>* pSgi, IssmDouble* poceanarea, IssmDouble* peustatic, SealevelMasks* masks);
166  void SealevelriseNonEustatic(Vector<IssmDouble>* pSgo, Vector<IssmDouble>* pSg_old, SealevelMasks* masks,bool verboseconvolution);
167  void SealevelriseRotationalFeedback(Vector<IssmDouble>* pRSLgo_rot, Vector<IssmDouble>* pRSLg_old, IssmDouble* pIxz, IssmDouble* pIyz, IssmDouble* pIzz, SealevelMasks* masks);
168  void SealevelriseElastic(Vector<IssmDouble>* pUp, Vector<IssmDouble>* pNorth, Vector<IssmDouble>* pEast, Vector<IssmDouble>* pSg_old, SealevelMasks* masks);
169  IssmDouble SealevelriseOceanAverage(Vector<IssmDouble>* Sg,SealevelMasks* masks, IssmDouble oceanarea);
170  #endif
171  void HydrologyEPLupdateDomainx(IssmDouble* pEplcount);
172  void HydrologyIDSupdateDomainx(IssmDouble* pIDScount);
173  void TimeAdaptx(IssmDouble* pdt);
177  void InitTransientInputx(int* transientinput_enum,int numoutputs);
178  void StackTransientInputx(int* input_enum,int* transientinput_enum,IssmDouble hydrotime,int numoutputs);
179  void AverageTransientInputx(int* transientinput_enum,int* averagedinput_enum,IssmDouble init_time,IssmDouble end_time,int numoutputs,int averaging_method);
180  void UpdateConstraintsx(void);
181  int UpdateVertexPositionsx(void);
182 
183  #ifdef _HAVE_JAVASCRIPT_
184  FemModel(IssmDouble* buffer, int buffersize, char* toolkits, char* solution, char* modelname,ISSM_MPI_Comm incomm, bool trace=false);
185  void CleanUpJs(char** poutput, size_t* psize);
186  void InitFromBuffers(char* buffer, int buffersize, char* toolkits, int solution_type,bool trace,IssmPDouble* X=NULL);
187  #endif
188 
189  /*AMR*/
190  #if !defined(_HAVE_AD_)
191  void ReMesh(void);
192  void BedrockFromMismipPlus(void);
193  void AdjustBaseThicknessAndMask(void);
194  void GetMesh(Vertices* femmodel_vertices,Elements* femmodel_elements,IssmDouble** px, IssmDouble** py, int** pelementslist);
195  void GetMesh(int** elementslist, IssmDouble** x, IssmDouble** y, int* numberofvertices, int* numberofelements);
196  void SetMesh(int** elementslist, IssmDouble** x, IssmDouble** y, int* numberofvertices, int* numberofelements);
197  void GetMeshOnPartition(Vertices* femmodel_vertices,Elements* femmodel_elements,IssmDouble** px, IssmDouble** py, IssmDouble** pz, int** pelementslist,int** psidtoindex);
198  void CreateElements(int newnumberofelements,int elementswidth,int* newelementslist,bool* my_elements,Elements* elements);
199  void CreateMaterials(int newnumberofelements,bool* my_elements,Materials* materials);
200  void CreateConstraints(Vertices* newfemmodel_vertices,int analysis_enum,Constraints* newfemmodel_constraints);
201  void GetInputs(int* pnumP0inputs,IssmDouble** pP0inputs,int** pP0input_enums,int** pP0input_interp,int* pnumP1inputs,IssmDouble** pP1inputs,int** pP1input_enums,int** pP1input_interp);
202  void InterpolateInputs(Vertices* newfemmodel_vertices,Elements* newfemmodel_elements,Inputs2* new_inputs);
203  void UpdateElements(int newnumberofelements,int* newelementslist,bool* my_elements,int analysis_counter,Elements* newelements);
204  void WriteMeshInResults(void);
206  void SmoothedDeviatoricStressTensor(IssmDouble** ptauxx,IssmDouble** ptauyy,IssmDouble** ptauxy); //nodal values, just for SSA-P1: TauXX, TauYY, TauXY
207  void ZZErrorEstimator(IssmDouble** pelementerror);
208  void SmoothedGradThickness(IssmDouble** pdHdx,IssmDouble** pdHdy);
209  void ThicknessZZErrorEstimator(IssmDouble** pelementerror);
210  void MeanGroundedIceLevelSet(IssmDouble** pmasklevelset);
212  void GetZeroLevelSetPoints(IssmDouble** pzerolevelset_points,int &numberofpoints,int levelset_type);
213  #endif
214 
215  #if defined(_HAVE_BAMG_) && !defined(_HAVE_AD_)
216  void ReMeshBamg(int* pnewnumberofvertices,int* pnewnumberofelements,IssmDouble** pnewx,IssmDouble** pnewy,IssmDouble** pnewz,int** pnewelementslist);
217  void InitializeAdaptiveRefinementBamg(void);
218  void GethmaxVerticesFromZeroLevelSetDistance(IssmDouble* hmaxvertices,int levelset_type);
219  void GethmaxVerticesFromEstimators(IssmDouble* hmaxvertices,int errorestimator_type);
220  void GetVerticeDistanceToZeroLevelSet(IssmDouble** pverticedistance,int leveset_type);
221  #endif
222 
223  #if defined(_HAVE_NEOPZ_) && !defined(_HAVE_AD_)
224  void ReMeshNeopz(int* pnewnumberofvertices,int* pnewnumberofelements,IssmDouble** pnewx,IssmDouble** pnewy,IssmDouble** pnewz,int** pnewelementslist);
225  void InitializeAdaptiveRefinementNeopz(void);
226  void GetElementDistanceToZeroLevelSet(IssmDouble** pelementdistance,int levelset_type);
227  void SetRefPatterns(void);
228  #endif
229 };
230 
231 #endif
FemModel::TotalCalvingFluxLevelsetx
void TotalCalvingFluxLevelsetx(IssmDouble *pGbmb, bool scaled)
Definition: FemModel.cpp:2929
FemModel::MaxAbsVzx
void MaxAbsVzx(IssmDouble *presponse)
Definition: FemModel.cpp:1811
FemModel::AnalysisIndex
int AnalysisIndex(int)
Definition: FemModel.cpp:220
Vertices
Declaration of Vertices class.
Definition: Vertices.h:15
FemModel::Solve
void Solve(void)
Definition: FemModel.cpp:883
AmrBamg.h
FemModel::CreateElements
void CreateElements(int newnumberofelements, int elementswidth, int *newelementslist, bool *my_elements, Elements *elements)
Definition: FemModel.cpp:3675
FemModel::ElementOperationx
void ElementOperationx(void(Element::*function)(void))
Definition: FemModel.cpp:1243
FemModel::InitFromFiles
void InitFromFiles(char *rootpath, char *inputfilename, char *outputfilename, char *petscfilename, char *lockfilename, char *restartfilename, const int solution_type, bool trace, IssmPDouble *X=NULL)
Definition: FemModel.cpp:393
FemModel::vertices
Vertices * vertices
Definition: FemModel.h:49
IssmDouble
double IssmDouble
Definition: types.h:37
FemModel::MaxVelx
void MaxVelx(IssmDouble *presponse)
Definition: FemModel.cpp:1855
Nodes
Declaration of Nodes class.
Definition: Nodes.h:19
FemModel::MaxAbsVyx
void MaxAbsVyx(IssmDouble *presponse)
Definition: FemModel.cpp:1786
FemModel::GetLocalVectorWithClonesGset
void GetLocalVectorWithClonesGset(IssmDouble **plocal_ug, Vector< IssmDouble > *ug)
Definition: FemModel.cpp:1337
AdaptiveMeshRefinement
Definition: AdaptiveMeshRefinement.h:12
FemModel::InputMakeDiscontinuous
void InputMakeDiscontinuous(int enum_in)
Definition: FemModel.cpp:1611
FemModel::MinVyx
void MinVyx(IssmDouble *presponse)
Definition: FemModel.cpp:2005
FemModel::InitTransientInputx
void InitTransientInputx(int *transientinput_enum, int numoutputs)
Definition: FemModel.cpp:5202
FemModel::nummodels
int nummodels
Definition: FemModel.h:39
FemModel::BalancethicknessMisfitx
void BalancethicknessMisfitx(IssmDouble *pV)
Definition: FemModel.cpp:962
Parameters
Declaration of Parameters class.
Definition: Parameters.h:18
FemModel::EtaDiffx
void EtaDiffx(IssmDouble *pJ)
Definition: FemModel.cpp:2112
FemModel::parameters
Parameters * parameters
Definition: FemModel.h:46
Constraints
Declaration of Constraints class.
Definition: Constraints.h:13
FemModel::Size
int Size(void)
Definition: FemModel.cpp:638
FemModel::StrainRateeffectivex
void StrainRateeffectivex()
Definition: FemModel.cpp:2646
FemModel::CostFunctionx
void CostFunctionx(IssmDouble *pJ, IssmDouble **pJlist, int *pn)
Definition: FemModel.cpp:1062
FemModel::GetMeshOnPartition
void GetMeshOnPartition(Vertices *femmodel_vertices, Elements *femmodel_elements, IssmDouble **px, IssmDouble **py, IssmDouble **pz, int **pelementslist, int **psidtoindex)
Definition: FemModel.cpp:3837
Elements
Declaration of Elements class.
Definition: Elements.h:17
FemModel::SolutionAnalysesList
void SolutionAnalysesList(int **panalyses, int *pnumanalyses, IoModel *iomodel, int solutiontype)
Definition: FemModel.cpp:646
FemModel::results
Results * results
Definition: FemModel.h:48
FemModel::CalvingMeltingFluxLevelsetx
void CalvingMeltingFluxLevelsetx()
Definition: FemModel.cpp:1054
FemModel::Echo
void Echo()
Definition: FemModel.cpp:382
FemModel::TotalCalvingMeltingFluxLevelsetx
void TotalCalvingMeltingFluxLevelsetx(IssmDouble *pGbmb, bool scaled)
Definition: FemModel.cpp:2945
FemModel::~FemModel
~FemModel()
Definition: FemModel.cpp:172
FemModel::AdjustBaseThicknessAndMask
void AdjustBaseThicknessAndMask(void)
Definition: FemModel.cpp:3314
FemModel::ThicknessAverage
void ThicknessAverage()
Definition: FemModel.cpp:2774
FemModel::Marshall
void Marshall(char **pmarshalled_data, int *pmarshalled_data_size, int marshall_direction)
Definition: FemModel.cpp:464
FemModel::OutputControlsx
void OutputControlsx(Results **presults)
Definition: FemModel.cpp:2171
FemModel::CalvingFluxLevelsetx
void CalvingFluxLevelsetx()
Definition: FemModel.cpp:1046
FemModel::TimeAdaptx
void TimeAdaptx(IssmDouble *pdt)
Definition: FemModel.cpp:2895
FemModel::constraints_list
Constraints ** constraints_list
Definition: FemModel.h:53
FemModel::SmoothedDeviatoricStressTensor
void SmoothedDeviatoricStressTensor(IssmDouble **ptauxx, IssmDouble **ptauyy, IssmDouble **ptauxy)
Definition: FemModel.cpp:4018
FemModel::TotalSmbx
void TotalSmbx(IssmDouble *pSmb, bool scaled)
Definition: FemModel.cpp:2993
FemModel::RequestedDependentsx
void RequestedDependentsx(void)
Definition: FemModel.cpp:2220
FemModel::GetLocalVectorWithClonesVertices
void GetLocalVectorWithClonesVertices(IssmDouble **plocal_vector, Vector< IssmDouble > *vector)
Definition: FemModel.cpp:1342
FemModel::InitFromFids
void InitFromFids(char *rootpath, FILE *IOMODEL, FILE *toolkitsoptionsfid, int in_solution_type, bool trace, IssmPDouble *X=NULL)
Definition: FemModel.cpp:424
FemModel::MeanGroundedIceLevelSet
void MeanGroundedIceLevelSet(IssmDouble **pmasklevelset)
Definition: FemModel.cpp:4330
FemModel::Divergencex
void Divergencex(IssmDouble *pdiv)
Definition: FemModel.cpp:1227
Element
Definition: Element.h:41
FemModel::CalvingRateLevermannx
void CalvingRateLevermannx()
Definition: FemModel.cpp:1038
FemModel::OmegaAbsGradientx
void OmegaAbsGradientx(IssmDouble *pJ)
Definition: FemModel.cpp:2055
FemModel::MaxAbsVxx
void MaxAbsVxx(IssmDouble *presponse)
Definition: FemModel.cpp:1761
FemModel::SyncLocalVectorWithClonesVertices
void SyncLocalVectorWithClonesVertices(IssmDouble *local_vector)
Definition: FemModel.cpp:1366
FemModel::HydrologyIDSupdateDomainx
void HydrologyIDSupdateDomainx(IssmDouble *pIDScount)
Definition: FemModel.cpp:5072
FemModel::HydrologyEPLupdateDomainx
void HydrologyEPLupdateDomainx(IssmDouble *pEplcount)
Definition: FemModel.cpp:4965
Materials
Declaration of Materials class.
Definition: Materials.h:16
FemModel::StressIntensityFactorx
void StressIntensityFactorx()
Definition: FemModel.cpp:2654
FemModel::UpdateConstraintsx
void UpdateConstraintsx(void)
Definition: FemModel.cpp:3027
FemModel::IceVolumex
void IceVolumex(IssmDouble *pV, bool scaled)
Definition: FemModel.cpp:1687
FemModel::nodes
Nodes * nodes
Definition: FemModel.h:56
FemModel::MinVelx
void MinVelx(IssmDouble *presponse)
Definition: FemModel.cpp:1955
FemModel::analysis_counter
int analysis_counter
Definition: FemModel.h:37
FemModel::GetInputLocalMinMaxOnNodesx
void GetInputLocalMinMaxOnNodesx(IssmDouble **pmin, IssmDouble **pmax, IssmDouble *ug)
Definition: FemModel.cpp:1308
Profiler
Definition: Profiler.h:33
FemModel::materials
Materials * materials
Definition: FemModel.h:45
FemModel::DeviatoricStressx
void DeviatoricStressx()
Definition: FemModel.cpp:1095
FemModel::MinVxx
void MinVxx(IssmDouble *presponse)
Definition: FemModel.cpp:1980
FemModel::analysis_type_list
int * analysis_type_list
Definition: FemModel.h:38
FemModel::GetElementsWidth
int GetElementsWidth()
Definition: FemModel.h:79
FemModel::CalvingRateVonmisesx
void CalvingRateVonmisesx()
Definition: FemModel.cpp:1030
FemModel::TotalFloatingBmbx
void TotalFloatingBmbx(IssmDouble *pFbmb, bool scaled)
Definition: FemModel.cpp:2961
FemModel::ZZErrorEstimator
void ZZErrorEstimator(IssmDouble **pelementerror)
Definition: FemModel.cpp:4106
FemModel::AverageTransientInputx
void AverageTransientInputx(int *transientinput_enum, int *averagedinput_enum, IssmDouble init_time, IssmDouble end_time, int numoutputs, int averaging_method)
Definition: FemModel.cpp:5246
FemModel::MaxVyx
void MaxVyx(IssmDouble *presponse)
Definition: FemModel.cpp:1905
FemModel::MaxVxx
void MaxVxx(IssmDouble *presponse)
Definition: FemModel.cpp:1880
FemModel::inputs2
Inputs2 * inputs2
Definition: FemModel.h:47
FemModel::ElementResponsex
void ElementResponsex(IssmDouble *presponse, int response_enum)
Definition: FemModel.cpp:1252
FemModel::constraints
Constraints * constraints
Definition: FemModel.h:52
FemModel::UpdateConstraintsExtrudeFromBasex
void UpdateConstraintsExtrudeFromBasex()
Definition: FemModel.cpp:3009
FemModel::GroundinglineMassFluxx
void GroundinglineMassFluxx(IssmDouble *presponse, bool scaled)
Definition: FemModel.cpp:1627
FemModel::GetLocalVectorWithClonesNodes
void GetLocalVectorWithClonesNodes(IssmDouble **plocal_vector, Vector< IssmDouble > *vector)
Definition: FemModel.cpp:1472
Inputs2
Declaration of Inputs class.
Definition: Inputs2.h:23
FemModel::loads
Loads * loads
Definition: FemModel.h:54
FemModel::UpdateConstraintsL2ProjectionEPLx
void UpdateConstraintsL2ProjectionEPLx(IssmDouble *pL2count)
Definition: FemModel.cpp:5160
FemModel::nodes_list
Nodes ** nodes_list
Definition: FemModel.h:57
FemModel::BedrockFromMismipPlus
void BedrockFromMismipPlus(void)
Definition: FemModel.cpp:3283
FemModel::solution_type
int solution_type
Definition: FemModel.h:40
FemModel::StackTransientInputx
void StackTransientInputx(int *input_enum, int *transientinput_enum, IssmDouble hydrotime, int numoutputs)
Definition: FemModel.cpp:5213
FemModel::elements
Elements * elements
Definition: FemModel.h:44
FemModel::Responsex
void Responsex(IssmDouble *presponse, int response_descriptor_enum)
Definition: FemModel.cpp:2558
FemModel::RignotMeltParameterizationx
void RignotMeltParameterizationx()
Definition: FemModel.cpp:2622
FemModel::TotalGroundedBmbx
void TotalGroundedBmbx(IssmDouble *pGbmb, bool scaled)
Definition: FemModel.cpp:2977
FemModel::GetZeroLevelSetPoints
void GetZeroLevelSetPoints(IssmDouble **pzerolevelset_points, int &numberofpoints, int levelset_type)
Definition: FemModel.cpp:4400
FemModel::IceMassx
void IceMassx(IssmDouble *pV, bool scaled)
Definition: FemModel.cpp:1655
FemModel::ThicknessPositivex
void ThicknessPositivex(IssmDouble *pJ)
Definition: FemModel.cpp:2837
FemModel
Definition: FemModel.h:31
FemModel::IcefrontMassFluxx
void IcefrontMassFluxx(IssmDouble *presponse, bool scaled)
Definition: FemModel.cpp:1579
FemModel::MinVzx
void MinVzx(IssmDouble *presponse)
Definition: FemModel.cpp:2030
FemModel::StrainRateperpendicularx
void StrainRateperpendicularx()
Definition: FemModel.cpp:2638
FemModel::IcefrontAreax
void IcefrontAreax()
Definition: FemModel.cpp:1546
FemModel::SyncLocalVectorWithClonesVerticesAdd
void SyncLocalVectorWithClonesVerticesAdd(IssmDouble *local_vector)
Definition: FemModel.cpp:1405
FemModel::IceVolumeAboveFloatationx
void IceVolumeAboveFloatationx(IssmDouble *pV, bool scaled)
Definition: FemModel.cpp:1671
FemModel::ResetLevelset
void ResetLevelset()
Definition: FemModel.cpp:2545
Loads
Declaration of Loads class.
Definition: Loads.h:16
FemModel::CheckPoint
void CheckPoint(void)
Definition: FemModel.cpp:238
FemModel::CreateMaterials
void CreateMaterials(int newnumberofelements, bool *my_elements, Materials *materials)
Definition: FemModel.cpp:3721
FemModel::GetInputs
void GetInputs(int *pnumP0inputs, IssmDouble **pP0inputs, int **pP0input_enums, int **pP0input_interp, int *pnumP1inputs, IssmDouble **pP1inputs, int **pP1input_enums, int **pP1input_interp)
Definition: FemModel.cpp:3368
FemModel::UpdateElements
void UpdateElements(int newnumberofelements, int *newelementslist, bool *my_elements, int analysis_counter, Elements *newelements)
Definition: FemModel.cpp:3993
FemModel::MassFluxx
void MassFluxx(IssmDouble *presponse)
Definition: FemModel.cpp:1703
FemModel::Restart
void Restart(void)
Definition: FemModel.cpp:549
FemModel::CreateConstraints
void CreateConstraints(Vertices *newfemmodel_vertices, int analysis_enum, Constraints *newfemmodel_constraints)
Definition: FemModel.cpp:3894
FemModel::MaxVzx
void MaxVzx(IssmDouble *presponse)
Definition: FemModel.cpp:1930
ISSM_MPI_Comm
int ISSM_MPI_Comm
Definition: issmmpi.h:118
FemModel::SurfaceAbsMisfitx
void SurfaceAbsMisfitx(IssmDouble *pJ)
Definition: FemModel.cpp:2663
FemModel::RequestedOutputsx
void RequestedOutputsx(Results **presults, char **requested_outputs, int numoutputs, bool save_results=true)
Definition: FemModel.cpp:2267
FemModel::SetCurrentConfiguration
void SetCurrentConfiguration(int configuration_type)
Definition: FemModel.cpp:634
FemModel::GetElementCenterCoordinates
void GetElementCenterCoordinates(IssmDouble **pxc, IssmDouble **pyc)
Definition: FemModel.cpp:4355
FemModel::CleanUp
void CleanUp(void)
Definition: FemModel.cpp:279
FemModel::copy
FemModel * copy()
Definition: FemModel.cpp:320
FemModel::FloatingAreax
void FloatingAreax(IssmDouble *pV, bool scaled)
Definition: FemModel.cpp:1292
FemModel::InterpolateInputs
void InterpolateInputs(Vertices *newfemmodel_vertices, Elements *newfemmodel_elements, Inputs2 *new_inputs)
Definition: FemModel.cpp:3499
FemModel::StrainRateparallelx
void StrainRateparallelx()
Definition: FemModel.cpp:2630
FemModel::IcefrontMassFluxLevelsetx
void IcefrontMassFluxLevelsetx(IssmDouble *presponse, bool scaled)
Definition: FemModel.cpp:1595
FemModel::ThicknessAbsGradientx
void ThicknessAbsGradientx(IssmDouble *pJ)
Definition: FemModel.cpp:2718
IoModel
Definition: IoModel.h:48
FemModel::ThicknessZZErrorEstimator
void ThicknessZZErrorEstimator(IssmDouble **pelementerror)
Definition: FemModel.cpp:4260
FemModel::loads_list
Loads ** loads_list
Definition: FemModel.h:55
IssmPDouble
IssmDouble IssmPDouble
Definition: types.h:38
DataSet
Declaration of DataSet class.
Definition: DataSet.h:14
FemModel::SetMesh
void SetMesh(int **elementslist, IssmDouble **x, IssmDouble **y, int *numberofvertices, int *numberofelements)
Definition: FemModel.cpp:3819
FemModel::DistanceToFieldValue
void DistanceToFieldValue(int fieldenum, IssmDouble fieldvalue, int distanceenum)
Definition: FemModel.cpp:1103
FemModel::SmoothedGradThickness
void SmoothedGradThickness(IssmDouble **pdHdx, IssmDouble **pdHdy)
Definition: FemModel.cpp:4180
FemModel::WriteMeshInResults
void WriteMeshInResults(void)
Definition: FemModel.cpp:3610
FemModel::profiler
Profiler * profiler
Definition: FemModel.h:42
Vector< IssmDouble >
AmrBamg
Definition: AmrBamg.h:12
Results
Declaration of Results class.
Definition: Results.h:14
FemModel::ReMesh
void ReMesh(void)
Definition: FemModel.cpp:3101
FemModel::UpdateConstraintsExtrudeFromTopx
void UpdateConstraintsExtrudeFromTopx()
Definition: FemModel.cpp:3018
FemModel::MaxDivergencex
void MaxDivergencex(IssmDouble *pdiv)
Definition: FemModel.cpp:1836
FemModel::GetMesh
void GetMesh(Vertices *femmodel_vertices, Elements *femmodel_elements, IssmDouble **px, IssmDouble **py, int **pelementslist)
Definition: FemModel.cpp:3731
FemModel::UpdateVertexPositionsx
int UpdateVertexPositionsx(void)
Definition: FemModel.cpp:3055
SealevelMasks
Definition: SealevelMasks.h:10
FemModel::WriteErrorEstimatorsInResults
void WriteErrorEstimatorsInResults(void)
Definition: FemModel.cpp:3642
FemModel::FemModel
FemModel(int argc, char **argv, ISSM_MPI_Comm comm_init, bool trace=false)
Definition: FemModel.cpp:56
AdaptiveMeshRefinement.h
FemModel::GroundedAreax
void GroundedAreax(IssmDouble *pV, bool scaled)
Definition: FemModel.cpp:1530