Skip to main content
  • Home
  • Development
  • Documentation
  • Donate
  • Operational login
  • Browse the archive

swh logo
SoftwareHeritage
Software
Heritage
Archive
Features
  • Search

  • Downloads

  • Save code now

  • Add forge now

  • Help

https://github.com/geodynamics/citcoms
26 October 2024, 18:30:49 UTC
  • Code
  • Branches (31)
  • Releases (16)
  • Visits
    • Branches
    • Releases
    • HEAD
    • refs/heads/CitcomS
    • refs/heads/CitcomS-2_x
    • refs/heads/compressible
    • refs/heads/cxx
    • refs/heads/eheien
    • refs/heads/eheien_dev
    • refs/heads/gui-launcher-branch
    • refs/heads/kommu/aspect-citcom-benchmarks
    • refs/heads/master
    • refs/heads/python-removal
    • refs/heads/rajesh-petsc
    • refs/heads/rajesh-petsc-schur
    • refs/heads/rajesh_dev
    • refs/heads/v2
    • refs/heads/v2.2
    • refs/heads/v3.0
    • refs/heads/v3.1
    • refs/remotes/svn/CitcomS
    • refs/remotes/svn/CitcomS-2_x
    • refs/remotes/svn/compressible
    • refs/remotes/svn/cxx
    • refs/remotes/svn/eheien
    • refs/remotes/svn/eheien_dev
    • refs/remotes/svn/gui-launcher-branch
    • refs/remotes/svn/trunk
    • refs/remotes/svn/v2
    • refs/remotes/svn/v2.2
    • refs/remotes/svn/v3.0
    • refs/remotes/svn/v3.1
    • refs/tags/v3.3.0
    • refs/tags/v3.3.1
    • v3.2.0
    • v3.1.2
    • v3.1.1
    • v3.1.0
    • v3.0.2
    • v3.0.1
    • v3.0.0b
    • v3.0.0
    • v2.2.2
    • v2.2.1
    • v2.2.0
    • v2.1.0
    • v2.0.2
    • v2.0.1
    • pre-2.0
    • 3.2.0
  • 2ec1c85
  • /
  • lib
  • /
  • Global_operations.c
Raw File Download
Take a new snapshot of a software origin

If the archived software origin currently browsed is not synchronized with its upstream version (for instance when new commits have been issued), you can explicitly request Software Heritage to take a new snapshot of it.

Use the form below to proceed. Once a request has been submitted and accepted, it will be processed as soon as possible. You can then check its processing state by visiting this dedicated page.
swh spinner

Processing "take a new snapshot" request ...

Permalinks

To reference or cite the objects present in the Software Heritage archive, permalinks based on SoftWare Hash IDentifiers (SWHIDs) must be used.
Select below a type of object currently browsed in order to display its associated SWHID and permalink.

  • content
  • directory
  • revision
  • snapshot
origin badgecontent badge Iframe embedding
swh:1:cnt:ecf0eb4ecf9e35311da8efc13459ae58dc2ea291
origin badgedirectory badge Iframe embedding
swh:1:dir:6555f79617665d181317e44231cba813dccd61e7
origin badgerevision badge
swh:1:rev:4eaf0dce8350a336dd310a788f96f236db7f0d7f
origin badgesnapshot badge
swh:1:snp:6144a327b54924075059de040a51fd832e5692cc
Citations

This interface enables to generate software citations, provided that the root directory of browsed objects contains a citation.cff or codemeta.json file.
Select below a type of object currently browsed in order to generate citations for them.

  • content
  • directory
  • revision
  • snapshot
Generate software citation in BibTex format (requires biblatex-software package)
Generating citation ...
Generate software citation in BibTex format (requires biblatex-software package)
Generating citation ...
Generate software citation in BibTex format (requires biblatex-software package)
Generating citation ...
Generate software citation in BibTex format (requires biblatex-software package)
Generating citation ...
Tip revision: 4eaf0dce8350a336dd310a788f96f236db7f0d7f authored by Eh Tan on 08 September 2008, 19:18:16 UTC
Updated ChangeLog to r12826
Tip revision: 4eaf0dc
Global_operations.c
/*
 *~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
 *
 *<LicenseText>
 *
 * CitcomS by Louis Moresi, Shijie Zhong, Lijie Han, Eh Tan,
 * Clint Conrad, Michael Gurnis, and Eun-seo Choi.
 * Copyright (C) 1994-2005, California Institute of Technology.
 *
 * This program is free software; you can redistribute it and/or modify
 * it under the terms of the GNU General Public License as published by
 * the Free Software Foundation; either version 2 of the License, or
 * (at your option) any later version.
 *
 * This program is distributed in the hope that it will be useful,
 * but WITHOUT ANY WARRANTY; without even the implied warranty of
 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
 * GNU General Public License for more details.
 *
 * You should have received a copy of the GNU General Public License
 * along with this program; if not, write to the Free Software
 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
 *
 *</LicenseText>
 *
 *~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
 */
#include <mpi.h>

#include <math.h>
#include <sys/types.h>
#include "element_definitions.h"
#include "global_defs.h"

/* ===============================================
   strips horizontal average from nodal field X.
   Assumes orthogonal mesh, otherwise, horizontals
   aren't & another method is required.
   =============================================== */

void remove_horiz_ave(E,X,H,store_or_not)
     struct All_variables *E;
     double **X, *H;
     int store_or_not;

{
    int m,i,j,k,n,nox,noz,noy;
    void return_horiz_ave();

    const int dims = E->mesh.nsd;

    noy = E->lmesh.noy;
    noz = E->lmesh.noz;
    nox = E->lmesh.nox;

    return_horiz_ave(E,X,H);

  for(m=1;m<=E->sphere.caps_per_proc;m++)
    for(k=1;k<=noy;k++)
      for(j=1;j<=nox;j++)
	for(i=1;i<=noz;i++) {
            n = i+(j-1)*noz+(k-1)*noz*nox;
            X[m][n] -= H[i];
	}

   return;
}


void remove_horiz_ave2(struct All_variables *E, double **X)
{
    double *H;

    H = (double *)malloc( (E->lmesh.noz+1)*sizeof(double));
    remove_horiz_ave(E, X, H, 0);
    free ((void *) H);
}


void return_horiz_ave(E,X,H)
     struct All_variables *E;
     double **X, *H;
{
  const int dims = E->mesh.nsd;
  int m,i,j,k,d,nint,noz,nox,noy,el,elz,elx,ely,j1,j2,i1,i2,k1,k2,nproc;
  int top,lnode[5], sizeofH, noz2,iroot;
  double *Have,*temp,aa[5];
  struct Shape_function1 M;
  struct Shape_function1_dA dGamma;
  void get_global_1d_shape_fn();

  sizeofH = (2*E->lmesh.noz+2)*sizeof(double);

  Have = (double *)malloc(sizeofH);
  temp = (double *)malloc(sizeofH);

  noz = E->lmesh.noz;
  noy = E->lmesh.noy;
  elz = E->lmesh.elz;
  elx = E->lmesh.elx;
  ely = E->lmesh.ely;
  noz2 = 2*noz;

  for (i=1;i<=elz;i++)  {
    temp[i] = temp[i+noz] = 0.0;
    temp[i+1] = temp[i+1+noz] = 0.0;
    top = 0;
    if (i==elz) top = 1;
    for (m=1;m<=E->sphere.caps_per_proc;m++)
      for (k=1;k<=ely;k++)
        for (j=1;j<=elx;j++)     {
          el = i + (j-1)*elz + (k-1)*elx*elz;
          get_global_1d_shape_fn(E,el,&M,&dGamma,top,m);

          lnode[1] = E->ien[m][el].node[1];
          lnode[2] = E->ien[m][el].node[2];
          lnode[3] = E->ien[m][el].node[3];
          lnode[4] = E->ien[m][el].node[4];

          for(nint=1;nint<=onedvpoints[E->mesh.nsd];nint++)   {
            for(d=1;d<=onedvpoints[E->mesh.nsd];d++)
              temp[i] += X[m][lnode[d]] * E->M.vpt[GMVINDEX(d,nint)]
                          * dGamma.vpt[GMVGAMMA(0,nint)];
            temp[i+noz] += dGamma.vpt[GMVGAMMA(0,nint)];
            }

          if (i==elz)  {
            lnode[1] = E->ien[m][el].node[5];
            lnode[2] = E->ien[m][el].node[6];
            lnode[3] = E->ien[m][el].node[7];
            lnode[4] = E->ien[m][el].node[8];

            for(nint=1;nint<=onedvpoints[E->mesh.nsd];nint++)   {
              for(d=1;d<=onedvpoints[E->mesh.nsd];d++)
                temp[i+1] += X[m][lnode[d]] * E->M.vpt[GMVINDEX(d,nint)]
                          * dGamma.vpt[GMVGAMMA(1,nint)];
              temp[i+1+noz] += dGamma.vpt[GMVGAMMA(1,nint)];
              }

            }   /* end of if i==elz    */
          }   /* end of j  and k, and m  */
     }        /* Done for i */

  MPI_Allreduce(temp,Have,noz2+1,MPI_DOUBLE,MPI_SUM,E->parallel.horizontal_comm);

  for (i=1;i<=noz;i++) {
    if(Have[i+noz] != 0.0)
       H[i] = Have[i]/Have[i+noz];
    }
 /* if (E->parallel.me==0)
    for(i=1;i<=noz;i++)
      fprintf(stderr,"area %d %d %g\n",E->parallel.me,i,Have[i+noz]);
*/
  free ((void *) Have);
  free ((void *) temp);

  return;
  }

void return_horiz_ave_f(E,X,H)
     struct All_variables *E;
     float **X, *H;
{
  const int dims = E->mesh.nsd;
  int m,i,j,k,d,nint,noz,nox,noy,el,elz,elx,ely,j1,j2,i1,i2,k1,k2,nproc;
  int top,lnode[5], sizeofH, noz2,iroot;
  float *Have,*temp,aa[5];
  struct Shape_function1 M;
  struct Shape_function1_dA dGamma;
  void get_global_1d_shape_fn();

  sizeofH = (2*E->lmesh.noz+2)*sizeof(float);

  Have = (float *)malloc(sizeofH);
  temp = (float *)malloc(sizeofH);

  noz = E->lmesh.noz;
  noy = E->lmesh.noy;
  elz = E->lmesh.elz;
  elx = E->lmesh.elx;
  ely = E->lmesh.ely;
  noz2 = 2*noz;

  for (i=1;i<=elz;i++)  {
    temp[i] = temp[i+noz] = 0.0;
    temp[i+1] = temp[i+1+noz] = 0.0;
    top = 0;
    if (i==elz) top = 1;
    for (m=1;m<=E->sphere.caps_per_proc;m++)
      for (k=1;k<=ely;k++)
        for (j=1;j<=elx;j++)     {
          el = i + (j-1)*elz + (k-1)*elx*elz;
          get_global_1d_shape_fn(E,el,&M,&dGamma,top,m);

          lnode[1] = E->ien[m][el].node[1];
          lnode[2] = E->ien[m][el].node[2];
          lnode[3] = E->ien[m][el].node[3];
          lnode[4] = E->ien[m][el].node[4];

          for(nint=1;nint<=onedvpoints[E->mesh.nsd];nint++)   {
            for(d=1;d<=onedvpoints[E->mesh.nsd];d++)
              temp[i] += X[m][lnode[d]] * E->M.vpt[GMVINDEX(d,nint)]
                          * dGamma.vpt[GMVGAMMA(0,nint)];
            temp[i+noz] += dGamma.vpt[GMVGAMMA(0,nint)];
            }

          if (i==elz)  {
            lnode[1] = E->ien[m][el].node[5];
            lnode[2] = E->ien[m][el].node[6];
            lnode[3] = E->ien[m][el].node[7];
            lnode[4] = E->ien[m][el].node[8];

            for(nint=1;nint<=onedvpoints[E->mesh.nsd];nint++)   {
              for(d=1;d<=onedvpoints[E->mesh.nsd];d++)
                temp[i+1] += X[m][lnode[d]] * E->M.vpt[GMVINDEX(d,nint)]
                          * dGamma.vpt[GMVGAMMA(1,nint)];
              temp[i+1+noz] += dGamma.vpt[GMVGAMMA(1,nint)];
              }

            }   /* end of if i==elz    */
          }   /* end of j  and k, and m  */
     }        /* Done for i */

  MPI_Allreduce(temp,Have,noz2+1,MPI_FLOAT,MPI_SUM,E->parallel.horizontal_comm);

  for (i=1;i<=noz;i++) {
    if(Have[i+noz] != 0.0)
       H[i] = Have[i]/Have[i+noz];
    }
 /* if (E->parallel.me==0)
    for(i=1;i<=noz;i++)
      fprintf(stderr,"area %d %d %g\n",E->parallel.me,i,Have[i+noz]);
*/
  free ((void *) Have);
  free ((void *) temp);

  return;
  }


/******* RETURN ELEMENTWISE HORIZ AVE ********************************/
/*                                                                   */
/* This function is similar to return_horiz_ave in the citcom code   */
/* however here, elemental horizontal averages are given rather than */
/* nodal averages. Also note, here is average per element            */

void return_elementwise_horiz_ave(E,X,H)
     struct All_variables *E;
     double **X, *H;
{

  int m,i,j,k,d,noz,noy,el,elz,elx,ely,nproc;
  int sizeofH;
  int elz2;
  double *Have,*temp;

  sizeofH = (2*E->lmesh.elz+2)*sizeof(double);

  Have = (double *)malloc(sizeofH);
  temp = (double *)malloc(sizeofH);

  noz = E->lmesh.noz;
  noy = E->lmesh.noy;
  elz = E->lmesh.elz;
  elx = E->lmesh.elx;
  ely = E->lmesh.ely;
  elz2 = 2*elz;

  for (i=0;i<=(elz*2+1);i++)
  {
    temp[i]=0.0;
  }

  for (i=1;i<=elz;i++)
  {
    for (m=1;m<=E->sphere.caps_per_proc;m++)
    {
      for (k=1;k<=ely;k++)
      {
        for (j=1;j<=elx;j++)
        {
          el = i + (j-1)*elz + (k-1)*elx*elz;
          temp[i] += X[m][el]*E->ECO[E->mesh.levmax][m][el].area;
          temp[i+elz] += E->ECO[E->mesh.levmax][m][el].area;
        }
      }
    }
  }



/* determine which processors should get the message from me for
               computing the layer averages */

  MPI_Allreduce(temp,Have,elz2+1,MPI_DOUBLE,MPI_SUM,E->parallel.horizontal_comm);

  for (i=1;i<=elz;i++) {
    if(Have[i+elz] != 0.0)
       H[i] = Have[i]/Have[i+elz];
    }


  free ((void *) Have);
  free ((void *) temp);

  return;
}

float return_bulk_value(E,Z,average)
     struct All_variables *E;
     float **Z;
     int average;

{
    void get_global_shape_fn();
    void float_global_operation();

    double rtf[4][9];

    int n,i,j,k,el,m;
    float volume,integral,volume1,integral1;

    struct Shape_function GN;
    struct Shape_function_dx GNx;
    struct Shape_function_dA dOmega;

    const int vpts = vpoints[E->mesh.nsd];
    const int ends = enodes[E->mesh.nsd];
    const int sphere_key=1;

    volume1=0.0;
    integral1=0.0;

    for (m=1;m<=E->sphere.caps_per_proc;m++)
       for (el=1;el<=E->lmesh.nel;el++)  {

	  get_global_shape_fn(E,el,&GN,&GNx,&dOmega,0,sphere_key,rtf,E->mesh.levmax,m);

	  for(j=1;j<=vpts;j++)
	    for(i=1;i<=ends;i++) {
		n = E->ien[m][el].node[i];
		volume1 += E->N.vpt[GNVINDEX(i,j)] * dOmega.vpt[j];
		integral1 += Z[m][n] * E->N.vpt[GNVINDEX(i,j)] * dOmega.vpt[j];
                }

          }


    MPI_Allreduce(&volume1  ,&volume  ,1,MPI_FLOAT,MPI_SUM,E->parallel.world);
    MPI_Allreduce(&integral1,&integral,1,MPI_FLOAT,MPI_SUM,E->parallel.world);

    if(average && volume != 0.0)
 	   integral /= volume;

    return((float)integral);
}

/************ RETURN BULK VALUE_D *****************************************/
/*                                                                        */
/* Same as return_bulk_value but allowing double instead of float.        */
/* I think when integer average =1, volume average is returned.           */
/*         when integer average =0, integral is returned.           */


double return_bulk_value_d(E,Z,average)
     struct All_variables *E;
     double **Z;
     int average;

{
    void get_global_shape_fn();

    double rtf[4][9];
    int n,i,j,el,m;
    double volume,integral,volume1,integral1;

    struct Shape_function GN;
    struct Shape_function_dx GNx;
    struct Shape_function_dA dOmega;

    const int vpts = vpoints[E->mesh.nsd];
    const int ends = enodes[E->mesh.nsd];
    const int sphere_key=1;

    volume1=0.0;
    integral1=0.0;

    for (m=1;m<=E->sphere.caps_per_proc;m++)
       for (el=1;el<=E->lmesh.nel;el++)  {

	  get_global_shape_fn(E,el,&GN,&GNx,&dOmega,0,sphere_key,rtf,E->mesh.levmax,m);
          for(j=1;j<=vpts;j++)
            for(i=1;i<=ends;i++) {
                n = E->ien[m][el].node[i];
                volume1 += E->N.vpt[GNVINDEX(i,j)] * dOmega.vpt[j];
                integral1 += Z[m][n] * E->N.vpt[GNVINDEX(i,j)] * dOmega.vpt[j];
            }

       }


    MPI_Allreduce(&volume1  ,&volume  ,1,MPI_DOUBLE,MPI_SUM,E->parallel.world);
    MPI_Allreduce(&integral1,&integral,1,MPI_DOUBLE,MPI_SUM,E->parallel.world);

    if(average && volume != 0.0)
           integral /= volume;

    return((double)integral);
}

/* ================================================== */
float find_max_horizontal(E,Tmax)
struct All_variables *E;
float Tmax;
{
 float ttmax;

 MPI_Allreduce(&Tmax,&ttmax,1,MPI_FLOAT,MPI_MAX,E->parallel.horizontal_comm);

 return(ttmax);
 }

/* ================================================== */
void sum_across_surface(E,data,total)
struct All_variables *E;
float *data;
int total;
{
 int j,d;
 float *temp;

 temp = (float *)malloc((total+1)*sizeof(float));
 MPI_Allreduce(data,temp,total,MPI_FLOAT,MPI_SUM,E->parallel.horizontal_comm);

 for (j=0;j<total;j++) {
   data[j] = temp[j];
 }

 free((void *)temp);

 return;
}

/* ================================================== */
/* ================================================== */

/* ================================================== */
void sum_across_surf_sph1(E,sphc,sphs)
struct All_variables *E;
float *sphc,*sphs;
{
 int jumpp,total,j,d;
 float *sphcs,*temp;

 temp = (float *) malloc((E->sphere.hindice*2)*sizeof(float));
 sphcs = (float *) malloc((E->sphere.hindice*2)*sizeof(float));

 /* pack */
 jumpp = E->sphere.hindice;
 total = E->sphere.hindice*2;
 for (j=0;j<E->sphere.hindice;j++)   {
   sphcs[j] = sphc[j];
   sphcs[j+jumpp] = sphs[j];
 }

 /* sum across processors in horizontal direction */
 MPI_Allreduce(sphcs,temp,total,MPI_FLOAT,MPI_SUM,E->parallel.horizontal_comm);

 /* unpack */
 for (j=0;j<E->sphere.hindice;j++)   {
   sphc[j] = temp[j];
   sphs[j] = temp[j+jumpp];
 }

 free((void *)temp);
 free((void *)sphcs);

 return;
}

/* ================================================== */


float global_fvdot(E,A,B,lev)
   struct All_variables *E;
   float **A,**B;
   int lev;

{
  int m,i,neq;
  float prod, temp,temp1;

  neq=E->lmesh.NEQ[lev];

  temp = 0.0;
  temp1 = 0.0;
  prod = 0.0;
  for (m=1;m<=E->sphere.caps_per_proc;m++)  {
    neq=E->lmesh.NEQ[lev];
    temp1 = 0.0;
    for (i=0;i<neq;i++)
      temp += A[m][i]*B[m][i];

    for (i=1;i<=E->parallel.Skip_neq[lev][m];i++)
       temp1 += A[m][E->parallel.Skip_id[lev][m][i]]*B[m][E->parallel.Skip_id[lev][m][i]];

    temp -= temp1;

    }

  MPI_Allreduce(&temp, &prod,1,MPI_FLOAT,MPI_SUM,E->parallel.world);

  return (prod);
}


double kineticE_radial(E,A,lev)
   struct All_variables *E;
   double **A;
   int lev;

{
  int m,i,neq;
  double prod, temp,temp1;

    temp = 0.0;
    prod = 0.0;

  for (m=1;m<=E->sphere.caps_per_proc;m++)  {
    neq=E->lmesh.NEQ[lev];
    temp1 = 0.0;
    for (i=0;i<neq;i++)
      if ((i+1)%3==0)
        temp += A[m][i]*A[m][i];

    for (i=1;i<=E->parallel.Skip_neq[lev][m];i++)
      if ((E->parallel.Skip_id[lev][m][i]+1)%3==0)
        temp1 += A[m][E->parallel.Skip_id[lev][m][i]]*A[m][E->parallel.Skip_id[lev][m][i]];

    temp -= temp1;

    }

  MPI_Allreduce(&temp, &prod,1,MPI_DOUBLE,MPI_SUM,E->parallel.world);

  return (prod);
}

double global_vdot(E,A,B,lev)
   struct All_variables *E;
   double **A,**B;
   int lev;

{
  int m,i,neq;
  double prod, temp,temp1;

    temp = 0.0;
    prod = 0.0;

  for (m=1;m<=E->sphere.caps_per_proc;m++)  {
    neq=E->lmesh.NEQ[lev];
    temp1 = 0.0;
    for (i=0;i<neq;i++)
      temp += A[m][i]*B[m][i];

    for (i=1;i<=E->parallel.Skip_neq[lev][m];i++)
       temp1 += A[m][E->parallel.Skip_id[lev][m][i]]*B[m][E->parallel.Skip_id[lev][m][i]];

    temp -= temp1;

    }

  MPI_Allreduce(&temp, &prod,1,MPI_DOUBLE,MPI_SUM,E->parallel.world);

  return (prod);
}


double global_pdot(E,A,B,lev)
   struct All_variables *E;
   double **A,**B;
   int lev;

{
  int i,m,npno;
  double prod, temp;

  npno=E->lmesh.NPNO[lev];

  temp = 0.0;
  prod = 0.0;
  for (m=1;m<=E->sphere.caps_per_proc;m++)  {
    npno=E->lmesh.NPNO[lev];
    for (i=1;i<=npno;i++)
      temp += A[m][i]*B[m][i];
    }

  MPI_Allreduce(&temp, &prod,1,MPI_DOUBLE,MPI_SUM,E->parallel.world);

  return (prod);
  }


double global_tdot_d(E,A,B,lev)
   struct All_variables *E;
   double **A,**B;
   int lev;

{
  int i,nno,m;
  double prod, temp;

  nno=E->lmesh.NNO[lev];

  temp = 0.0;
  prod = 0.0;
  for (m=1;m<=E->sphere.caps_per_proc;m++)  {
    nno=E->lmesh.NNO[lev];
    for (i=1;i<=nno;i++)
    if (!(E->NODE[lev][m][i] & SKIP))
      temp += A[m][i];
    }

  MPI_Allreduce(&temp, &prod,1,MPI_DOUBLE,MPI_SUM,E->parallel.world);

  return (prod);
  }

float global_tdot(E,A,B,lev)
   struct All_variables *E;
   float **A,**B;
   int lev;

{
  int i,nno,m;
  float prod, temp;


  temp = 0.0;
  prod = 0.0;
  for (m=1;m<=E->sphere.caps_per_proc;m++)  {
    nno=E->lmesh.NNO[lev];
    for (i=1;i<=nno;i++)
      if (!(E->NODE[lev][m][i] & SKIP))
        temp += A[m][i]*B[m][i];
    }

  MPI_Allreduce(&temp, &prod,1,MPI_FLOAT,MPI_SUM,E->parallel.world);

  return (prod);
  }


float global_fmin(E,a)
   struct All_variables *E;
   float a;
{
  float temp;
  MPI_Allreduce(&a, &temp,1,MPI_FLOAT,MPI_MIN,E->parallel.world);
  return (temp);
  }

double global_dmax(E,a)
   struct All_variables *E;
   double a;
{
  double temp;
  MPI_Allreduce(&a, &temp,1,MPI_DOUBLE,MPI_MAX,E->parallel.world);
  return (temp);
  }


float global_fmax(E,a)
   struct All_variables *E;
   float a;
{
  float temp;
  MPI_Allreduce(&a, &temp,1,MPI_FLOAT,MPI_MAX,E->parallel.world);
  return (temp);
  }

double Tmaxd(E,T)
  struct All_variables *E;
  double **T;
{
  double global_dmax(),temp,temp1;
  int i,m;

  temp = -10.0;
  for (m=1;m<=E->sphere.caps_per_proc;m++)
    for(i=1;i<=E->lmesh.nno;i++)
      temp = max(T[m][i],temp);

  temp1 = global_dmax(E,temp);
  return (temp1);
  }


float Tmax(E,T)
  struct All_variables *E;
  float **T;
{
  float global_fmax(),temp,temp1;
  int i,m;

  temp = -10.0;
  for (m=1;m<=E->sphere.caps_per_proc;m++)
    for(i=1;i<=E->lmesh.nno;i++)
      temp = max(T[m][i],temp);

  temp1 = global_fmax(E,temp);
  return (temp1);
  }


double  vnorm_nonnewt(E,dU,U,lev)
  struct All_variables *E;
  double **dU,**U;
  int lev;
{
 double temp1,temp2,dtemp,temp;
 int a,e,i,m,node;
 const int dims = E->mesh.nsd;
 const int ends = enodes[dims];
 const int nel=E->lmesh.nel;

 dtemp=0.0;
 temp=0.0;
for (m=1;m<=E->sphere.caps_per_proc;m++)
  for (e=1;e<=nel;e++)
   /*if (E->mat[m][e]==1)*/
     for (i=1;i<=dims;i++)
       for (a=1;a<=ends;a++) {
	 node = E->IEN[lev][m][e].node[a];
         dtemp += dU[m][ E->ID[lev][m][node].doff[i] ]*
                  dU[m][ E->ID[lev][m][node].doff[i] ];
         temp += U[m][ E->ID[lev][m][node].doff[i] ]*
                 U[m][ E->ID[lev][m][node].doff[i] ];
         }


  MPI_Allreduce(&dtemp, &temp2,1,MPI_DOUBLE,MPI_SUM,E->parallel.world);
  MPI_Allreduce(&temp, &temp1,1,MPI_DOUBLE,MPI_SUM,E->parallel.world);

  temp1 = sqrt(temp2/temp1);

  return (temp1);
}


void sum_across_depth_sph1(E,sphc,sphs)
     struct All_variables *E;
     float *sphc,*sphs;
{
    int jumpp,total,j;

    float *sphcs,*temp;

    if (E->parallel.nprocz > 1)  {
	total = E->sphere.hindice*2;
	temp = (float *) malloc(total*sizeof(float));
	sphcs = (float *) malloc(total*sizeof(float));

	/* pack sphc[] and sphs[] into sphcs[] */
	jumpp = E->sphere.hindice;
	for (j=0;j<E->sphere.hindice;j++)   {
	    sphcs[j] = sphc[j];
	    sphcs[j+jumpp] = sphs[j];
	}

	/* sum across processors in z direction */
	MPI_Allreduce(sphcs, temp, total, MPI_FLOAT, MPI_SUM,
		      E->parallel.vertical_comm);

	/* unpack */
	for (j=0;j<E->sphere.hindice;j++)   {
	    sphc[j] = temp[j];
	    sphs[j] = temp[j+jumpp];
	}

	free(temp);
	free(sphcs);
    }


    return;
}


/* ================================================== */
/* ================================================== */
void broadcast_vertical(struct All_variables *E,
                        float *sphc, float *sphs,
                        int root)
{
    int jumpp, total, j;
    float *temp;

    if(E->parallel.nprocz == 1) return;

    jumpp = E->sphere.hindice;
    total = E->sphere.hindice*2;
    temp = (float *) malloc(total*sizeof(float));

    if (E->parallel.me_loc[3] == root) {
        /* pack */
        for (j=0; j<E->sphere.hindice; j++)   {
            temp[j] = sphc[j];
            temp[j+jumpp] = sphs[j];
        }
    }

    MPI_Bcast(temp, total, MPI_FLOAT, root, E->parallel.vertical_comm);

    if (E->parallel.me_loc[3] != root) {
        /* unpack */
        for (j=0; j<E->sphere.hindice; j++)   {
            sphc[j] = temp[j];
            sphs[j] = temp[j+jumpp];
        }
    }

    free((void *)temp);

    return;
}


/*
 * remove rigid body rotation from the velocity
 */

void remove_rigid_rot(struct All_variables *E)
{
    void velo_from_element_d();
    double myatan();
    double wx, wy, wz, v_theta, v_phi;
    double vx[9], vy[9], vz[9];
    double r, t, f;

    double exyz[4], fxyz[4];

    int m, e, i, k, j, node;

    const int nno = E->lmesh.nno;
    const int ends = ENODES3D;
    const int ppts = PPOINTS3D;
    const int vpts = VPOINTS3D;
    const int sphere_key = 1;
    double VV[4][9];
    double rot, fr, tr;

    /* Note: no need to weight in rho(r) here. */
    double moment_of_inertia = (8.0*M_PI/15.0)*
        (pow(E->sphere.ro,(double)5.0) - pow(E->sphere.ri,(double)5.0));

    /* compute and add angular momentum components */

    exyz[1] = exyz[2] = exyz[3] = 0.0;
    fxyz[1] = fxyz[2] = fxyz[3] = 0.0;


    for (m=1;m<=E->sphere.caps_per_proc;m++) {

        for (e=1;e<=E->lmesh.nel;e++) {

            t = E->eco[m][e].centre[1];
            f = E->eco[m][e].centre[2];
            r = E->eco[m][e].centre[3];

            velo_from_element_d(E,VV,m,e,sphere_key);

            for (j=1;j<=ppts;j++)   {
                vx[j] = 0.0;
                vy[j] = 0.0;
            }

            for (j=1;j<=ppts;j++)   {
                for (i=1;i<=ends;i++)   {
                    vx[j] += VV[1][i]*E->N.ppt[GNPINDEX(i,j)];
                    vy[j] += VV[2][i]*E->N.ppt[GNPINDEX(i,j)];
                }
            }

            wx = -r*vy[1];
            wy = r*vx[1];

            exyz[1] += (wx*cos(t)*cos(f)-wy*sin(f)) * E->eco[m][e].area;
            exyz[2] += (wx*cos(t)*sin(f)+wy*cos(f)) * E->eco[m][e].area;
            exyz[3] -= (wx*sin(t)                 ) * E->eco[m][e].area;
        }
    } /* end cap */

    MPI_Allreduce(exyz,fxyz,4,MPI_DOUBLE,MPI_SUM,E->parallel.world);

    fxyz[1] = fxyz[1] / moment_of_inertia;
    fxyz[2] = fxyz[2] / moment_of_inertia;
    fxyz[3] = fxyz[3] / moment_of_inertia;


    rot = sqrt(fxyz[1]*fxyz[1] + fxyz[2]*fxyz[2] + fxyz[3]*fxyz[3]);
    fr = myatan(fxyz[2], fxyz[1]);
    tr = acos(fxyz[3] / rot);

    if (E->parallel.me==0) {
        fprintf(E->fp,"Rigid rotation: rot=%e tr=%e fr=%e\n",rot,tr*180/M_PI,fr*180/M_PI);
        fprintf(stderr,"Rigid rotation: rot=%e tr=%e fr=%e\n",rot,tr*180/M_PI,fr*180/M_PI);
    }


    /* remove rigid rotation */

    for (m=1;m<=E->sphere.caps_per_proc;m++)  {
        for (node=1;node<=nno;node++)   {

            v_theta = E->sx[m][3][node] * rot * sin(tr)
                * sin(fr - E->sx[m][2][node]);
            v_phi = E->sx[m][3][node] * rot
                * ( sin(E->sx[m][1][node]) * cos(tr)
                    - cos(E->sx[m][1][node]) * sin(tr)
                    * cos(fr-E->sx[m][2][node]) );

            E->sphere.cap[m].V[1][node] -= v_theta;
            E->sphere.cap[m].V[2][node] -= v_phi;

        }
    }


    return;

}

back to top

Software Heritage — Copyright (C) 2015–2025, The Software Heritage developers. License: GNU AGPLv3+.
The source code of Software Heritage itself is available on our development forge.
The source code files archived by Software Heritage are available under their own copyright and licenses.
Terms of use: Archive access, API— Contact— JavaScript license information— Web API