Ice Sheet System Model  4.18
Code documentation
Functions
SurfaceAreax.cpp File Reference
#include "./SurfaceAreax.h"
#include "../../shared/shared.h"
#include "../../toolkits/toolkits.h"
#include "../InputUpdateFromConstantx/InputUpdateFromConstantx.h"

Go to the source code of this file.

Functions

void SurfaceAreax (IssmDouble *pS, FemModel *femmodel)
 

Function Documentation

◆ SurfaceAreax()

void SurfaceAreax ( IssmDouble pS,
FemModel femmodel 
)

Definition at line 11 of file SurfaceAreax.cpp.

11  {
12 
13  /*Intermediary*/
14  Element* element=NULL;
15 
16  /*output: */
17  IssmDouble S = 0.;
18  IssmDouble S_sum;
19 
20  /*Compute gradients: */
21  for(int i=0;i<femmodel->elements->Size();i++){
22  element=xDynamicCast<Element*>(femmodel->elements->GetObjectByOffset(i));
23  S+=element->SurfaceArea();
24  }
25 
26  /*Sum all J from all cpus of the cluster:*/
29  S=S_sum;
30 
31  /*add surface area to element inputs:*/
33 
34  /*Assign output pointers: */
35  if(pS) *pS=S;
36 }
DataSet::Size
int Size()
Definition: DataSet.cpp:399
IssmDouble
double IssmDouble
Definition: types.h:37
Element::SurfaceArea
virtual IssmDouble SurfaceArea(void)=0
ISSM_MPI_SUM
#define ISSM_MPI_SUM
Definition: issmmpi.h:134
IssmComm::GetComm
static ISSM_MPI_Comm GetComm(void)
Definition: IssmComm.cpp:30
Element
Definition: Element.h:41
ISSM_MPI_DOUBLE
#define ISSM_MPI_DOUBLE
Definition: issmmpi.h:125
SurfaceAreaEnum
@ SurfaceAreaEnum
Definition: EnumDefinitions.h:820
InputUpdateFromConstantx
void InputUpdateFromConstantx(FemModel *femmodel, bool constant, int name)
Definition: InputUpdateFromConstantx.cpp:10
FemModel::elements
Elements * elements
Definition: FemModel.h:44
ISSM_MPI_Bcast
int ISSM_MPI_Bcast(void *buffer, int count, ISSM_MPI_Datatype datatype, int root, ISSM_MPI_Comm comm)
Definition: issmmpi.cpp:162
DataSet::GetObjectByOffset
Object * GetObjectByOffset(int offset)
Definition: DataSet.cpp:334
ISSM_MPI_Reduce
int ISSM_MPI_Reduce(void *sendbuf, void *recvbuf, int count, ISSM_MPI_Datatype datatype, ISSM_MPI_Op op, int root, ISSM_MPI_Comm comm)
Definition: issmmpi.cpp:373
femmodel
FemModel * femmodel
Definition: esmfbinders.cpp:16